LORENE
map.h
1 /*
2  * Definition of Lorene classes Map
3  * Map_radial
4  * Map_af
5  * Map_et
6  * Map_log
7  *
8  */
9 
10 /*
11  * Copyright (c) 1999-2000 Jean-Alain Marck
12  * Copyright (c) 1999-2003 Eric Gourgoulhon
13  * Copyright (c) 1999-2001 Philippe Grandclement
14  * Copyright (c) 2000-2001 Jerome Novak
15  * Copyright (c) 2000-2001 Keisuke Taniguchi
16  *
17  * This file is part of LORENE.
18  *
19  * LORENE is free software; you can redistribute it and/or modify
20  * it under the terms of the GNU General Public License as published by
21  * the Free Software Foundation; either version 2 of the License, or
22  * (at your option) any later version.
23  *
24  * LORENE is distributed in the hope that it will be useful,
25  * but WITHOUT ANY WARRANTY; without even the implied warranty of
26  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
27  * GNU General Public License for more details.
28  *
29  * You should have received a copy of the GNU General Public License
30  * along with LORENE; if not, write to the Free Software
31  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
32  *
33  */
34 
35 
36 #ifndef __MAP_H_
37 #define __MAP_H_
38 
39 /*
40  * $Id: map.h,v 1.69 2026/03/04 15:47:04 j_novak Exp $
41  * $Log: map.h,v $
42  * Revision 1.69 2026/03/04 15:47:04 j_novak
43  * Added 'verbose' flag (true by default) to limit the output of poisson solvers.
44  *
45  * Revision 1.68 2025/03/28 09:12:06 j_novak
46  * New method to change the multi-grid, to be used for de-aliasing.
47  *
48  * Revision 1.67 2025/03/04 16:31:11 j_novak
49  * Updated method "integrale" to "Scalar" objects.
50  *
51  * Revision 1.66 2025/03/04 13:16:50 j_novak
52  * New complete versions of Map_af::poisson() and Map_et::poisson() for Scalar, not using Cmp.
53  *
54  * Revision 1.65 2024/07/03 13:51:02 j_novak
55  * Added methods to interpolate a Scalar between Map_af and Map_star efficiently.
56  *
57  * Revision 1.64 2024/06/25 13:56:24 j_novak
58  * Split map.h into map.h and map_star.h.
59  *
60  * Revision 1.63 2023/05/24 09:58:44 g_servignat
61  * Added Map_eps, an ellipsoidal radial mapping
62  *
63  * Revision 1.62 2022/06/20 14:21:55 g_servignat
64  * Implemented shells in Map_star
65  *
66  * Revision 1.61 2022/05/30 12:39:53 g_servignat
67  * New Map_star class for general adaptive domain to star boundaries
68  *
69  * Revision 1.60 2018/12/05 15:43:45 j_novak
70  * New Map_af constructor from a formatted file.
71  *
72  * Revision 1.59 2014/10/13 08:52:35 j_novak
73  * Lorene classes and functions now belong to the namespace Lorene.
74  *
75  * Revision 1.58 2014/10/06 15:09:40 j_novak
76  * Modified #include directives to use c++ syntax.
77  *
78  * Revision 1.57 2014/01/14 13:24:02 b_peres
79  * *** empty log message ***
80  *
81  * Revision 1.56 2012/01/24 14:58:54 j_novak
82  * Removed functions XXX_fait_xi()
83  *
84  * Revision 1.55 2012/01/17 15:34:20 j_penner
85  * *** empty log message ***
86  *
87  * Revision 1.54 2012/01/17 10:20:07 j_penner
88  * added a member cxi that allows for direct access to the computational coordinates in each domain
89  *
90  * Revision 1.53 2008/09/29 13:23:51 j_novak
91  * Implementation of the angular mapping associated with an affine
92  * mapping. Things must be improved to take into account the domain index.
93  *
94  * Revision 1.52 2007/10/16 21:52:10 e_gourgoulhon
95  * Added method poisson_compact for multi-domains.
96  *
97  * Revision 1.51 2007/05/15 12:44:18 p_grandclement
98  * Scalar version of reevaluate
99  *
100  * Revision 1.50 2007/05/06 10:48:08 p_grandclement
101  * Modification of a few operators for the vorton project
102  *
103  * Revision 1.49 2007/01/16 15:05:59 n_vasset
104  * New constructor (taking a Scalar in mono-domain angular grid for
105  * boundary) for function sol_elliptic_boundary
106  *
107  * Revision 1.48 2006/08/31 12:10:51 j_novak
108  * More comments for Map_af::avance_dalembert().
109  *
110  * Revision 1.47 2006/05/26 09:00:09 j_novak
111  * New members for multiplication or division by cos(theta).
112  *
113  * Revision 1.46 2006/04/25 07:21:54 p_grandclement
114  * Various changes for the NS_BH project
115  *
116  * Revision 1.45 2005/11/30 11:09:03 p_grandclement
117  * Changes for the Bin_ns_bh project
118  *
119  * Revision 1.44 2005/11/24 09:25:06 j_novak
120  * Added the Scalar version for the Laplacian
121  *
122  * Revision 1.43 2005/09/15 15:51:25 j_novak
123  * The "rotation" (change of triad) methods take now Scalars as default
124  * arguments.
125  *
126  * Revision 1.42 2005/08/26 14:02:38 p_grandclement
127  * Modification of the elliptic solver that matches with an oscillatory exterior solution
128  * small correction in Poisson tau also...
129  *
130  * Revision 1.41 2005/08/25 12:14:07 p_grandclement
131  * Addition of a new method to solve the scalar Poisson equation, based on a multi-domain Tau-method
132  *
133  * Revision 1.40 2005/06/09 07:56:24 f_limousin
134  * Implement a new function sol_elliptic_boundary() and
135  * Vector::poisson_boundary(...) which solve the vectorial poisson
136  * equation (method 6) with an inner boundary condition.
137  *
138  * Revision 1.39 2005/04/04 21:30:41 e_gourgoulhon
139  * Added argument lambda to method poisson_angu
140  * to treat the generalized angular Poisson equation:
141  * Lap_ang u + lambda u = source.
142  *
143  * Revision 1.38 2004/12/29 16:37:22 k_taniguchi
144  * Addition of some functions with the multipole falloff condition.
145  *
146  * Revision 1.37 2004/12/02 09:33:04 p_grandclement
147  * *** empty log message ***
148  *
149  * Revision 1.36 2004/11/30 20:42:05 k_taniguchi
150  * Addition of some functions with the falloff condition and a method
151  * to resize the external shell.
152  *
153  * Revision 1.35 2004/11/23 12:39:12 f_limousin
154  * Intoduce function poisson_dir_neu(...) to solve a scalar poisson
155  * equation with a mixed boundary condition (Dirichlet + Neumann).
156  *
157  * Revision 1.34 2004/10/11 15:08:59 j_novak
158  * The radial manipulation functions take Scalar as arguments, instead of Cmp.
159  * Added a conversion operator from Scalar to Cmp.
160  * The Cmp radial manipulation function make conversion to Scalar, call to the
161  * Map_radial version with a Scalar argument and back.
162  *
163  * Revision 1.33 2004/10/08 13:34:35 j_novak
164  * Scalar::div_r() does not need to pass through Cmp version anymore.
165  *
166  * Revision 1.32 2004/08/24 09:14:40 p_grandclement
167  * Addition of some new operators, like Poisson in 2d... It now requieres the
168  * GSL library to work.
169  *
170  * Also, the way a variable change is stored by a Param_elliptic is changed and
171  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
172  * will requiere some modification. (It should concern only the ones about monopoles)
173  *
174  * Revision 1.31 2004/07/27 08:24:26 j_novak
175  * Modif. comments
176  *
177  * Revision 1.30 2004/07/26 16:02:21 j_novak
178  * Added a flag to specify whether the primitive should be zero either at r=0
179  * or at r going to infinity.
180  *
181  * Revision 1.29 2004/06/22 08:49:57 p_grandclement
182  * Addition of everything needed for using the logarithmic mapping
183  *
184  * Revision 1.28 2004/06/14 15:23:53 e_gourgoulhon
185  * Added virtual function primr for computation of radial primitives.
186  *
187  * Revision 1.27 2004/03/31 11:22:23 f_limousin
188  * Methods Map_et::poisson_interne and Map_af::poisson_interne have been
189  * implemented to solve the continuity equation for strange stars.
190  *
191  * Revision 1.26 2004/03/22 13:12:41 j_novak
192  * Modification of comments to use doxygen instead of doc++
193  *
194  * Revision 1.24 2004/03/01 09:57:02 j_novak
195  * the wave equation is solved with Scalars. It now accepts a grid with a
196  * compactified external domain, which the solver ignores and where it copies
197  * the values of the field from one time-step to the next.
198  *
199  * Revision 1.23 2004/02/11 09:47:44 p_grandclement
200  * Addition of a new elliptic solver, matching with the homogeneous solution
201  * at the outer shell and not solving in the external domain (more details
202  * coming soon ; check your local Lorene dealer...)
203  *
204  * Revision 1.22 2004/01/29 08:50:01 p_grandclement
205  * Modification of Map::operator==(const Map&) and addition of the surface
206  * integrales using Scalar.
207  *
208  * Revision 1.21 2004/01/28 16:46:22 p_grandclement
209  * Addition of the sol_elliptic_fixe_der_zero stuff
210  *
211  * Revision 1.20 2004/01/28 10:35:52 j_novak
212  * Added new methods mult_r() for Scalars. These do not change the dzpuis flag.
213  *
214  * Revision 1.19 2004/01/27 09:33:46 j_novak
215  * New method Map_radial::div_r_zec
216  *
217  * Revision 1.18 2004/01/26 16:16:15 j_novak
218  * Methods of gradient for Scalar s. The input can have any dzpuis.
219  *
220  * Revision 1.17 2004/01/19 21:38:21 e_gourgoulhon
221  * Corrected sign error in comments of Map_radial::dxdr.
222  *
223  * Revision 1.16 2003/12/30 22:52:47 e_gourgoulhon
224  * Class Map: added methods flat_met_spher() and flat_met_cart() to get
225  * flat metric associated with the coordinates described by the mapping.
226  *
227  * Revision 1.15 2003/12/11 14:48:47 p_grandclement
228  * Addition of ALL (and that is a lot !) the files needed for the general elliptic solver ... UNDER DEVELOPEMENT...
229  *
230  * Revision 1.14 2003/11/06 14:43:37 e_gourgoulhon
231  * Gave a name to const arguments in certain method prototypes (e.g.
232  * constructors) to correct a bug of DOC++.
233  *
234  * Revision 1.13 2003/11/04 22:54:49 e_gourgoulhon
235  * Added new virtual methods mult_cost, mult_sint and div_sint.
236  *
237  * Revision 1.12 2003/10/16 08:49:21 j_novak
238  * Added a flag to decide wether the output is in the Ylm or in the standard base.
239  *
240  * Revision 1.11 2003/10/15 21:08:22 e_gourgoulhon
241  * Added method poisson_angu.
242  *
243  * Revision 1.10 2003/10/15 16:03:35 j_novak
244  * Added the angular Laplace operator for Scalar.
245  *
246  * Revision 1.9 2003/10/15 10:27:33 e_gourgoulhon
247  * Classes Map, Map_af and Map_et: added new methods dsdt, stdsdp and div_tant.
248  * Class Map_radial: added new Coord's : drdt and stdrdp.
249  *
250  * Revision 1.8 2003/06/20 14:14:53 f_limousin
251  * Add the operator== to compare two Cmp.
252  *
253  * Revision 1.7 2003/06/20 09:27:09 j_novak
254  * Modif commentaires.
255  *
256  * Revision 1.6 2002/10/16 14:36:29 j_novak
257  * Reorganization of #include instructions of standard C++, in order to
258  * use experimental version 3 of gcc.
259  *
260  * Revision 1.5 2002/09/13 09:17:33 j_novak
261  * Modif. commentaires
262  *
263  * Revision 1.4 2002/06/17 14:05:16 j_novak
264  * friend functions are now also declared outside the class definition
265  *
266  * Revision 1.3 2002/05/07 07:06:37 e_gourgoulhon
267  * Compatibily with xlC compiler on IBM SP2:
268  * added declaration of functions map_af_fait_* and map_et_fait_*
269  * outside the classes declarations.
270  *
271  * Revision 1.2 2002/01/15 15:53:06 p_grandclement
272  * I have had a constructor fot map_et using the equation of the surface
273  * of the star.
274  *
275  * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
276  * LORENE
277  *
278  * Revision 2.110 2001/10/29 15:31:55 novak
279  * Ajout de Map_radial::div_r
280  *
281  * Revision 2.109 2001/10/16 10:02:49 novak
282  * *** empty log message ***
283  *
284  * Revision 2.108 2001/07/19 14:01:00 novak
285  * new arguments for Map_af::dalembert
286  *
287  * Revision 2.107 2001/02/26 17:28:31 eric
288  * Ajout de la fonction virtuelle resize.
289  *
290  * Revision 2.106 2001/01/10 11:03:00 phil
291  * ajout de homothetie interne
292  *
293  * Revision 2.105 2001/01/02 10:51:55 phil
294  * ajout integrale de surface a l'infini
295  *
296  * Revision 2.104 2000/10/23 13:59:48 eric
297  * Map_et::adapt: changement des arguments (en autre, ajout de nz_search).
298  *
299  * Revision 2.103 2000/10/20 09:39:19 phil
300  * changement commentaires
301  *
302  * Revision 2.102 2000/10/19 14:33:23 novak
303  * corrige oubli pour Map_et?
304  *
305  * Revision 2.101 2000/10/19 14:11:20 novak
306  * Ajout des fonctions membres Map::dalembert et Map_af::dalembert
307  * (etat experimental)
308  *
309  * Revision 2.100 2000/10/09 13:46:39 eric
310  * Ajout de la fonction virtuelle poisson2d.
311  *
312  * Revision 2.99 2000/09/19 15:28:55 phil
313  * *** empty log message ***
314  *
315  * Revision 2.98 2000/09/19 15:24:19 phil
316  * ajout du passage de cartesienne en spheriques
317  *
318  * Revision 2.97 2000/09/19 13:05:38 phil
319  * ajout integrale_surface
320  *
321  * Revision 2.96 2000/09/11 15:54:03 eric
322  * Suppression des methodes deriv_x, deriv_y et deriv_z.
323  * Introduction des methodes comp_x_from_spherical, etc...
324  *
325  * Revision 2.95 2000/09/07 15:27:58 keisuke
326  * Add a new argument Cmp& uu in Map_af::poisson_regular and Map_et::poisson_regular.
327  *
328  * Revision 2.94 2000/09/04 15:30:56 keisuke
329  * Modify the arguments of Map_af::poisson_regular and Map_et::poisson_regular.
330  *
331  * Revision 2.93 2000/09/04 13:36:19 keisuke
332  * Modify the explanation for "uu_div" in Map_et::poisson_regular.
333  *
334  * Revision 2.92 2000/08/31 15:50:12 keisuke
335  * Modify Map_af::poisson_regular.
336  * Add Map_et::poisson_regular and Map::poisson_regular.
337  *
338  * Revision 2.91 2000/08/31 13:03:22 eric
339  * Ajout de la fonction virtuelle mult_rsint.
340  *
341  * Revision 2.90 2000/08/28 16:17:37 keisuke
342  * Add "int nzet" in the argumant of Map_af::poisson_regular.
343  *
344  * Revision 2.89 2000/08/18 11:10:12 eric
345  * Classe Map_et: ajout de l'operateur d'affectation a un autre Map_et.
346  *
347  * Revision 2.88 2000/08/11 08:50:18 keisuke
348  * Modif Map_af::poisson_regular
349  *
350  * Revision 2.87 2000/08/10 12:54:00 keisuke
351  * Ajout de Map_af::poisson_regular
352  *
353  * Revision 2.86 2000/07/20 14:21:07 eric
354  * Ajout de la fonction div_rsint.
355  *
356  * Revision 2.85 2000/05/25 13:54:41 eric
357  * Modif commentaires
358  *
359  * Revision 2.84 2000/05/22 14:38:51 phil
360  * ajout de inc_dzpuis et dec_dzpuis
361  *
362  * Revision 2.83 2000/04/27 15:18:54 phil
363  * *** empty log message ***
364  *
365  * Revision 2.82 2000/03/20 13:33:23 phil
366  * commentaires
367  *
368  * Revision 2.81 2000/03/17 17:32:48 phil
369  * *** empty log message ***
370  *
371  * Revision 2.80 2000/03/17 17:01:54 phil
372  * *** empty log message ***
373  *
374  * Revision 2.79 2000/03/17 16:58:48 phil
375  * ajout de poisson_frontiere
376  *
377  * Revision 2.78 2000/03/06 11:29:51 eric
378  * Ajout du membre reeavaluate_symy.
379  *
380  * Revision 2.77 2000/02/15 15:08:21 eric
381  * Changement du Param dans Map_et::adapt : fact_echelle est desormais
382  * passe en double_mod.
383  *
384  * Revision 2.76 2000/02/15 10:26:25 phil
385  * commentaire +
386  * suppression de poisson_vect et poisson_vect_oohara
387  *
388  * Revision 2.75 2000/02/11 13:37:43 eric
389  * Ajout de la fonction convert_absolute.
390  *
391  * Revision 2.74 2000/02/09 09:53:37 phil
392  * ajout de poisson_vect_oohara
393  *
394  * Revision 2.73 2000/01/26 13:07:02 eric
395  * Reprototypage complet des routines de derivation:
396  * le resultat est desormais suppose alloue a l'exterieur de la routine
397  * et est passe en argument (Cmp& resu), si bien que le prototypage
398  * complet devient:
399  * void DERIV(const Cmp& ci, Cmp& resu)
400  *
401  * Revision 2.72 2000/01/24 17:08:21 eric
402  * Class Map_af : suppression de la fonction convert.
403  * suppression du constructeur par convertion d'un Map_et.
404  * ajout du constructeur par conversion d'un Map.
405  *
406  * Revision 2.71 2000/01/24 16:41:43 eric
407  * Ajout de la fonction virtuelle operator=(const Map_af& ).
408  * Classe Map_af : ajout de la fonction convert(const Map& ).
409  *
410  * Revision 2.70 2000/01/21 12:48:34 phil
411  * changement prototypage de Map::poisson_vect
412  *
413  * Revision 2.69 2000/01/20 16:35:05 phil
414  * *** empty log message ***
415  *
416  * Revision 2.68 2000/01/20 15:44:42 phil
417  * *** empty log message ***
418  *
419  * Revision 2.67 2000/01/20 15:31:56 phil
420  * *** empty log message ***
421  *
422  * Revision 2.66 2000/01/20 14:18:06 phil
423  * *** empty log message ***
424  *
425  * Revision 2.65 2000/01/20 13:16:34 phil
426  * *** empty log message ***
427  *
428  * Revision 2.64 2000/01/20 12:51:24 phil
429  * *** empty log message ***
430  *
431  * Revision 2.63 2000/01/20 12:45:28 phil
432  * *** empty log message ***
433  *
434  * Revision 2.62 2000/01/20 12:40:27 phil
435  * *** empty log message ***
436  *
437  * Revision 2.61 2000/01/20 11:27:54 phil
438  * ajout de poisson_vect
439  *
440  * Revision 2.60 2000/01/13 15:31:55 eric
441  * Modif commentaires/
442  *
443  * Revision 2.59 2000/01/12 16:02:57 eric
444  * Modif commentaires poisson_compact.
445  *
446  * Revision 2.58 2000/01/12 12:54:23 eric
447  * Ajout du Cmp null, *p_cmp_zero, et de la methode associee cmp_zero().
448  *
449  * Revision 2.57 2000/01/10 13:27:43 eric
450  * Ajout des bases vectorielles associees aux coordonnees :
451  * membres bvect_spher et bvect_cart.
452  *
453  * Revision 2.56 2000/01/10 09:12:47 eric
454  * Reprototypage de poisson_compact : Valeur -> Cmp, Tenseur.
455  * Suppression de poisson_compact_boucle.
456  * poisson_compact est desormais implementee au niveau Map_radial.
457  *
458  * Revision 2.55 2000/01/04 15:23:11 eric
459  * Classe Map_radial : les data sont listees en premier
460  * Introduction de la fonction reevalutate.
461  *
462  * Revision 2.54 2000/01/03 13:30:32 eric
463  * Ajout de la fonction adapt.
464  *
465  * Revision 2.53 1999/12/22 17:09:52 eric
466  * Modif commentaires.
467  *
468  * Revision 2.52 1999/12/21 16:26:25 eric
469  * Ajout du constructeur par conversion Map_af::Map_af(const Map_et&).
470  * Ajout des fonctions Map_af::set_alpha et Map_af::set_beta.
471  *
472  * Revision 2.51 1999/12/21 13:01:29 eric
473  * Changement de prototype de la routine poisson : la solution est
474  * desormais passee en argument (et non plus en valeur de retour)
475  * pour permettre l'initialisation de methodes de resolution iterative.
476  *
477  * Revision 2.50 1999/12/21 10:12:09 eric
478  * Modif commentaires.
479  *
480  * Revision 2.49 1999/12/21 10:06:05 eric
481  * Ajout de l'argument Param& a poisson.
482  *
483  * Revision 2.48 1999/12/20 15:44:35 eric
484  * Modif commentaires.
485  *
486  * Revision 2.47 1999/12/20 10:47:45 eric
487  * Modif commentaires.
488  *
489  * Revision 2.46 1999/12/20 10:24:12 eric
490  * Ajout des fonctions de lecture des parametres de Map_et:
491  * get_alpha(), get_beta(), get_ff(), get_gg().
492  *
493  * Revision 2.45 1999/12/16 14:50:08 eric
494  * Modif commentaires.
495  *
496  * Revision 2.44 1999/12/16 14:17:54 eric
497  * Introduction de l'argument const Param& par dans val_lx et val_lx_jk.
498  * (en remplacement de l'argument Tbl& param).
499  *
500  * Revision 2.43 1999/12/09 10:45:24 eric
501  * Ajout de la fonction virtuelle integrale.
502  *
503  * Revision 2.42 1999/12/07 14:50:47 eric
504  * Changement ordre des arguments val_r, val_lx
505  * val_r_kj --> val_r_jk
506  * val_lx_kj -->val_lx_jk
507  *
508  * Revision 2.41 1999/12/06 16:45:20 eric
509  * Surcharge de val_lx avec la version sans param.
510  *
511  * Revision 2.40 1999/12/06 15:33:44 eric
512  * Ajout des fonctions val_r_kj et val_lx_kj.
513  *
514  * Revision 2.39 1999/12/06 13:11:54 eric
515  * Introduction des fonctions val_r, val_lx et homothetie.
516  *
517  * Revision 2.38 1999/12/02 14:28:22 eric
518  * Reprototypage de la fonction poisson: Valeur -> Cmp.
519  *
520  * Revision 2.37 1999/11/30 14:19:33 eric
521  * Reprototypage complet des fonctions membres mult_r, mult_r_zec,
522  * dec2_dzpuis et inc2_dzpuis : Valeur --> Cmp
523  *
524  * Revision 2.36 1999/11/29 13:17:57 eric
525  * Modif commentaires.
526  *
527  * Revision 2.35 1999/11/29 12:55:42 eric
528  * Changement prototype de la fonction laplacien : Valeur --> Cmp.
529  *
530  * Revision 2.34 1999/11/25 16:27:27 eric
531  * Reorganisation complete du calcul des derivees partielles.
532  *
533  * Revision 2.33 1999/11/24 16:31:17 eric
534  * Map_et: ajout des fonctions set_ff et set_gg.
535  *
536  * Revision 2.32 1999/11/24 14:31:48 eric
537  * Map_af: les membres alpha et beta deviennent prives.
538  * Map_af: introduction des fonctions get_alpha() et get_beta().
539  *
540  * Revision 2.31 1999/11/24 11:22:09 eric
541  * Map_et : Coords rendus publics
542  * Map_et : fonctions de constructions amies.
543  *
544  * Revision 2.30 1999/11/22 10:32:39 eric
545  * Introduction de la classe Map_et.
546  * Constructeurs de Map rendus protected.
547  * Fonction del_coord() rebaptisee reset_coord().
548  *
549  * Revision 2.29 1999/10/27 16:44:41 phil
550  * ajout de mult_r_zec
551  *
552  * Revision 2.28 1999/10/19 14:40:37 phil
553  * ajout de inc2_dzpuis()
554  *
555  * Revision 2.27 1999/10/15 14:12:20 eric
556  * *** empty log message ***
557  *
558  * Revision 2.26 1999/10/14 14:26:06 eric
559  * Depoussierage.
560  * Documentation.
561  *
562  * Revision 2.25 1999/10/11 11:16:29 phil
563  * changement prototypage de poisson_compact_boucle
564  *
565  * Revision 2.24 1999/10/11 10:48:51 phil
566  * changement de nom pour poisson a support compact
567  *
568  * Revision 2.23 1999/10/04 09:20:58 phil
569  * changement de prototypage de void Map_af::poisson_nul
570  *
571  * Revision 2.22 1999/09/30 18:38:32 phil
572  * *** empty log message ***
573  *
574  * Revision 2.21 1999/09/30 18:33:10 phil
575  * ajout de poisson_nul et poisson_nul_boucle
576  *
577  * Revision 2.20 1999/09/30 16:45:54 phil
578  * ajout de Map_af::poisson_nul(const Valeur&, int, int)
579  *
580  * Revision 2.19 1999/09/16 13:15:40 phil
581  * ajout de Valeur mult_r (const Valeur &)
582  *
583  * Revision 2.18 1999/09/15 10:42:11 phil
584  * ajout de Valeur dec2_dzpuis(const Valeur&)
585  *
586  * Revision 2.17 1999/09/14 13:45:45 phil
587  * suppression de la divergence
588  *
589  * Revision 2.16 1999/09/13 15:09:07 phil
590  * ajout de Map_af::divergence
591  *
592  * Revision 2.15 1999/09/13 13:52:23 phil
593  * ajout des derivations partielles par rapport a x,y et z.
594  *
595  * Revision 2.14 1999/09/07 14:35:20 phil
596  * ajout de la fonction Valeur** gradient(const Valeur&)
597  *
598  * Revision 2.13 1999/04/26 16:37:43 phil
599  * *** empty log message ***
600  *
601  * Revision 2.12 1999/04/26 16:33:28 phil
602  * *** empty log message ***
603  *
604  * Revision 2.11 1999/04/26 13:53:04 phil
605  * *** empty log message ***
606  *
607  * Revision 2.10 1999/04/26 13:51:19 phil
608  * ajout de Map_af::laplacien (2versions)
609  *
610  * Revision 2.9 1999/04/14 09:04:01 phil
611  * *** empty log message ***
612  *
613  * Revision 2.8 1999/04/14 08:53:27 phil
614  * *** empty log message ***
615  *
616  * Revision 2.7 1999/04/13 17:45:25 phil
617  * *** empty log message ***
618  *
619  * Revision 2.6 1999/04/13 17:02:41 phil
620  * ,
621  *
622  * Revision 2.5 1999/04/13 17:00:41 phil
623  * ajout de la resolution de poisson affine
624  *
625  * Revision 2.4 1999/03/04 13:10:53 eric
626  * Ajout des Coord representant les derivees du changement de variable
627  * dans la classe Map_radial.
628  *
629  * Revision 2.3 1999/03/01 17:00:38 eric
630  * *** empty log message ***
631  *
632  * Revision 2.2 1999/03/01 16:44:41 eric
633  * Operateurs << et >> sur les ostream.
634  * L'operateur >> est virtuel.
635  *
636  * Revision 2.1 1999/02/22 15:21:45 hyc
637  * *** empty log message ***
638  *
639  *
640  * Revision 2.0 1999/01/15 09:10:39 hyc
641  * *** empty log message ***
642  *
643  * $Header: /cvsroot/Lorene/C++/Include/map.h,v 1.69 2026/03/04 15:47:04 j_novak Exp $
644  *
645  */
646 
647 #include "coord.h"
648 #include "base_vect.h"
649 #include "valeur.h"
650 
651 namespace Lorene {
652 class Scalar ;
653 class Cmp ;
654 class Param ;
655 class Map_af ;
656 class Map_et ;
657 class Tenseur ;
658 class Param_elliptic ;
659 class Metric_flat ;
660 class Tbl ;
661 class Itbl ;
662 
663  //------------------------------------//
664  // class Map //
665  //------------------------------------//
666 
696 class Map {
697 
698  // Data :
699  // ----
700  protected:
702  const Mg3d* mg ;
703 
704  double ori_x ;
705  double ori_y ;
706  double ori_z ;
707  double rot_phi ;
708 
716 
724 
729 
734 
740 
741  mutable Map_af* p_mp_angu ;
742 
743  public:
744  Coord r ;
751 
752  Coord x ;
753  Coord y ;
754  Coord z ;
755 
759 
760 
761  // Constructors, destructor :
762  // ------------------------
763 
764  protected:
765  explicit Map(const Mg3d& ) ;
766  Map(const Map &) ;
767  Map(const Mg3d&, FILE* ) ;
768 
769  public:
770  virtual ~Map() ;
771 
772  // Memory management
773  // -----------------
774  protected:
775  virtual void reset_coord() ;
776 
777  // Outputs
778  // -------
779  private:
780  virtual ostream& operator>>(ostream &) const = 0 ;
781 
782  public:
783  virtual void sauve(FILE* ) const ;
784 
785 
786  // Extraction of information
787  // -------------------------
788 
789  public:
791  const Mg3d* get_mg() const {return mg; };
792 
794  double get_ori_x() const {return ori_x;} ;
796  double get_ori_y() const {return ori_y;} ;
798  double get_ori_z() const {return ori_z;} ;
799 
801  double get_rot_phi() const {return rot_phi;} ;
802 
809  const Base_vect_spher& get_bvect_spher() const {return bvect_spher;} ;
810 
817  const Base_vect_cart& get_bvect_cart() const {return bvect_cart;} ;
818 
822  const Metric_flat& flat_met_spher() const ;
823 
827  const Metric_flat& flat_met_cart() const ;
828 
833  const Cmp& cmp_zero() const {return *p_cmp_zero;} ;
834 
838  virtual const Map_af& mp_angu(int) const = 0 ;
839 
850  void convert_absolute(double xx, double yy, double zz,
851  double& rr, double& theta, double& pphi) const ;
852 
861  virtual double val_r(int l, double xi, double theta, double pphi)
862  const = 0 ;
863 
872  virtual void val_lx(double rr, double theta, double pphi,
873  int& l, double& xi) const = 0 ;
874 
887  virtual void val_lx(double rr, double theta, double pphi,
888  const Param& par, int& l, double& xi) const = 0 ;
889 
890 
892  virtual bool operator==(const Map& ) const = 0;
893 
894 
895 
896  // Modification of the origin, the orientation and the radial scale:
897  // ----------------------------------------------------------------
898  public:
899  void set_ori(double xa0, double ya0, double za0) ;
900  void set_rot_phi(double phi0) ;
901 
906  virtual void homothetie(double lambda) = 0 ;
907 
917  virtual void resize(int l, double lambda) = 0 ;
918 
919  // Modification of the mapping
920  // ---------------------------
921  public:
923  virtual void operator=(const Map_af& ) = 0 ;
924 
930  virtual void adapt(const Cmp& ent, const Param& par, int nbr=0) = 0 ;
931 
936  void set_new_grid(const Mg3d& new_mg) ;
937 
938 
939  // Values of a Cmp at the new grid points
940  // --------------------------------------
941 
954  virtual void reevaluate(const Map* mp_prev, int nzet,
955  Cmp& uu) const = 0 ;
956 
970  virtual void reevaluate_symy(const Map* mp_prev, int nzet,
971  Cmp& uu) const = 0 ;
972 
985  virtual void reevaluate(const Map* mp_prev, int nzet,
986  Scalar& uu) const = 0 ;
987 
1001  virtual void reevaluate_symy(const Map* mp_prev, int nzet,
1002  Scalar& uu) const = 0 ;
1003 
1004  // Differential operators:
1005  // ----------------------
1006  public:
1013  virtual void dsdxi(const Cmp& ci, Cmp& resu) const = 0 ;
1014 
1021  virtual void dsdr(const Cmp& ci, Cmp& resu) const = 0 ;
1022 
1029  virtual void srdsdt(const Cmp& ci, Cmp& resu) const = 0 ;
1030 
1038  virtual void srstdsdp(const Cmp& ci, Cmp& resu) const = 0 ;
1039 
1047  virtual void dsdxi(const Scalar& uu, Scalar& resu) const = 0 ;
1048 
1056  virtual void dsdr(const Scalar& uu, Scalar& resu) const = 0 ;
1057 
1066  virtual void dsdradial(const Scalar& uu, Scalar& resu) const = 0 ;
1067 
1075  virtual void srdsdt(const Scalar& uu, Scalar& resu) const = 0 ;
1076 
1084  virtual void srstdsdp(const Scalar& uu, Scalar& resu) const = 0 ;
1085 
1090  virtual void dsdt(const Scalar& uu, Scalar& resu) const = 0 ;
1091 
1096  virtual void stdsdp(const Scalar& uu, Scalar& resu) const = 0 ;
1097 
1108  virtual void laplacien(const Scalar& uu, int zec_mult_r,
1109  Scalar& lap) const = 0 ;
1110 
1112  virtual void laplacien(const Cmp& uu, int zec_mult_r,
1113  Cmp& lap) const = 0 ;
1114 
1121  virtual void lapang(const Scalar& uu, Scalar& lap) const = 0 ;
1122 
1123 
1134  virtual void primr(const Scalar& uu, Scalar& resu,
1135  bool null_infty) const = 0 ;
1136 
1137 
1138  // Various linear operators
1139  // ------------------------
1140  public:
1144  virtual void mult_r(Scalar& uu) const = 0 ;
1145 
1149  virtual void mult_r(Cmp& ci) const = 0 ;
1150 
1154  virtual void mult_r_zec(Scalar& ) const = 0 ;
1155 
1158  virtual void mult_rsint(Scalar& ) const = 0 ;
1159 
1162  virtual void div_rsint(Scalar& ) const = 0 ;
1163 
1166  virtual void div_r(Scalar& ) const = 0 ;
1167 
1171  virtual void div_r_zec(Scalar& ) const = 0 ;
1172 
1175  virtual void mult_cost(Scalar& ) const = 0 ;
1176 
1179  virtual void div_cost(Scalar& ) const = 0 ;
1180 
1183  virtual void mult_sint(Scalar& ) const = 0 ;
1184 
1187  virtual void div_sint(Scalar& ) const = 0 ;
1188 
1191  virtual void div_tant(Scalar& ) const = 0 ;
1192 
1202  virtual void comp_x_from_spherical(const Scalar& v_r, const Scalar& v_theta,
1203  const Scalar& v_phi, Scalar& v_x) const = 0 ;
1205  virtual void comp_x_from_spherical(const Cmp& v_r, const Cmp& v_theta,
1206  const Cmp& v_phi, Cmp& v_x) const = 0 ;
1207 
1217  virtual void comp_y_from_spherical(const Scalar& v_r, const Scalar& v_theta,
1218  const Scalar& v_phi, Scalar& v_y) const = 0 ;
1219 
1221  virtual void comp_y_from_spherical(const Cmp& v_r, const Cmp& v_theta,
1222  const Cmp& v_phi, Cmp& v_y) const = 0 ;
1223 
1232  virtual void comp_z_from_spherical(const Scalar& v_r, const Scalar& v_theta,
1233  Scalar& v_z) const = 0 ;
1234 
1236  virtual void comp_z_from_spherical(const Cmp& v_r, const Cmp& v_theta,
1237  Cmp& v_z) const = 0 ;
1238 
1248  virtual void comp_r_from_cartesian(const Scalar& v_x, const Scalar& v_y,
1249  const Scalar& v_z, Scalar& v_r) const = 0 ;
1251  virtual void comp_r_from_cartesian(const Cmp& v_x, const Cmp& v_y,
1252  const Cmp& v_z, Cmp& v_r) const = 0 ;
1253 
1263  virtual void comp_t_from_cartesian(const Scalar& v_x, const Scalar& v_y,
1264  const Scalar& v_z, Scalar& v_t) const = 0 ;
1265 
1267  virtual void comp_t_from_cartesian(const Cmp& v_x, const Cmp& v_y,
1268  const Cmp& v_z, Cmp& v_t) const = 0 ;
1269 
1278  virtual void comp_p_from_cartesian(const Scalar& v_x, const Scalar& v_y,
1279  Scalar& v_p) const = 0 ;
1280 
1282  virtual void comp_p_from_cartesian(const Cmp& v_x, const Cmp& v_y,
1283  Cmp& v_p) const = 0 ;
1284 
1289  virtual void dec_dzpuis(Scalar& ) const = 0 ;
1290 
1295  virtual void dec2_dzpuis(Scalar& ) const = 0 ;
1296 
1301  virtual void inc_dzpuis(Scalar& ) const = 0 ;
1302 
1307  virtual void inc2_dzpuis(Scalar& ) const = 0 ;
1308 
1316  virtual Tbl* integrale(const Scalar&) const = 0 ;
1317 
1325  virtual Tbl* integrale(const Cmp&) const = 0 ;
1326 
1327  // PDE resolution :
1328  // --------------
1329  public:
1342  virtual void poisson(const Cmp& source, Param& par, Cmp& uu, bool verbose=true) const = 0 ;
1343 
1354  virtual void poisson_tau(const Cmp& source, Param& par, Cmp& uu) const = 0 ;
1355 
1368  virtual void poisson(const Scalar& source, Param& par, Scalar& uu, bool verbose=true) const = 0 ;
1369 
1380  virtual void poisson_tau(const Scalar& source, Param& par, Scalar& uu) const = 0 ;
1381 
1382  virtual void poisson_falloff(const Cmp& source, Param& par, Cmp& uu,
1383  int k_falloff) const = 0 ;
1384 
1385  virtual void poisson_ylm(const Cmp& source, Param& par, Cmp& pot,
1386  int nylm, double* intvec) const = 0 ;
1387 
1410  virtual void poisson_regular(const Cmp& source, int k_div, int nzet,
1411  double unsgam1, Param& par, Cmp& uu,
1412  Cmp& uu_regu, Cmp& uu_div,
1413  Tenseur& duu_div, Cmp& source_regu,
1414  Cmp& source_div) const = 0 ;
1415 
1429  virtual void poisson_compact(const Cmp& source, const Cmp& aa,
1430  const Tenseur& bb, const Param& par,
1431  Cmp& psi) const = 0 ;
1432 
1447  virtual void poisson_compact(int nzet, const Cmp& source, const Cmp& aa,
1448  const Tenseur& bb, const Param& par,
1449  Cmp& psi) const = 0 ;
1450 
1468  virtual void poisson_angu(const Scalar& source, Param& par,
1469  Scalar& uu, double lambda=0) const = 0 ;
1470 
1471  virtual void poisson_angu(const Cmp& source, Param& par,
1472  Cmp& uu, double lambda=0) const = 0 ;
1473 
1474 
1475  public:
1499  virtual Param* donne_para_poisson_vect (Param& para, int i) const = 0;
1500 
1521  virtual void poisson_frontiere (const Cmp& source,const Valeur& limite,
1522  int raccord, int num_front, Cmp& pot,
1523  double = 0., double = 0.) const = 0 ;
1524 
1525  virtual void poisson_frontiere_double (const Cmp& source, const Valeur& lim_func,
1526  const Valeur& lim_der, int num_zone, Cmp& pot) const = 0 ;
1527 
1528 
1539  virtual void poisson_interne (const Cmp& source, const Valeur& limite,
1540  Param& par, Cmp& pot) const = 0 ;
1541 
1542 
1566  virtual void poisson2d(const Cmp& source_mat, const Cmp& source_quad,
1567  Param& par, Cmp& uu, bool verbose=true) const = 0 ;
1568 
1582  virtual void dalembert(Param& par, Scalar& fJp1, const Scalar& fJ,
1583  const Scalar& fJm1, const Scalar& source) const = 0 ;
1584 
1585  // Friend functions :
1586  // ----------------
1587  friend ostream& operator<<(ostream& , const Map& ) ;
1588 };
1589 ostream& operator<<(ostream& , const Map& ) ;
1590 
1591 
1592 
1593  //------------------------------------//
1594  // class Map_radial //
1595  //------------------------------------//
1596 
1597 
1598 
1611 class Map_radial : public Map {
1612 
1613  // Data :
1614  // ----
1615 
1616  // 0th order derivatives of the mapping
1617  // - - - - - - - - - - - - - - - - - -
1618  public:
1625 
1626  // 1st order derivatives of the mapping
1627  // - - - - - - - - - - - - - - - - - -
1628  public:
1636 
1644 
1652 
1660 
1668 
1676 
1684 
1685  // 2nd order derivatives of the mapping
1686  // - - - - - - - - - - - - - - - - - -
1687  public:
1695 
1707 
1708 
1716 
1724 
1725 
1733 
1734 
1735  // Constructors, destructor :
1736  // ------------------------
1737 
1738  protected:
1740  Map_radial(const Mg3d& mgrid ) ;
1741  Map_radial(const Map_radial& mp) ;
1742  Map_radial (const Mg3d&, FILE* ) ;
1743 
1744  public:
1745  virtual ~Map_radial() ;
1746 
1747  // Memory management
1748  // -----------------
1749  protected:
1750  virtual void reset_coord() ;
1751  // Modification of the mapping
1752  // ---------------------------
1753  public:
1755  virtual void operator=(const Map_af& ) = 0 ;
1756 
1757  // Outputs
1758  // -------
1759  public:
1760  virtual void sauve(FILE* ) const ;
1761 
1762  // Extraction of information
1763  // -------------------------
1773  virtual double val_r_jk(int l, double xi, int j, int k) const = 0 ;
1774 
1785  virtual void val_lx_jk(double rr, int j, int k, const Param& par,
1786  int& l, double& xi) const = 0 ;
1787 
1789  virtual bool operator==(const Map& ) const = 0;
1790 
1791  // Values of a Cmp at the new grid points
1792  // --------------------------------------
1805  virtual void reevaluate(const Map* mp_prev, int nzet, Cmp& uu) const ;
1806 
1820  virtual void reevaluate_symy(const Map* mp_prev, int nzet, Cmp& uu)
1821  const ;
1822 
1835  virtual void reevaluate(const Map* mp_prev, int nzet, Scalar& uu) const ;
1836 
1850  virtual void reevaluate_symy(const Map* mp_prev, int nzet, Scalar& uu)
1851  const ;
1852 
1853  // Various linear operators
1854  // ------------------------
1855  public:
1859  virtual void mult_r(Scalar& uu) const ;
1860 
1864  virtual void mult_r(Cmp& ci) const ;
1865 
1870  virtual void mult_r_zec(Scalar& ) const ;
1871 
1874  virtual void mult_rsint(Scalar& ) const ;
1875 
1878  virtual void div_rsint(Scalar& ) const ;
1879 
1882  virtual void div_r(Scalar& ) const ;
1883 
1887  virtual void div_r_zec(Scalar& ) const ;
1888 
1891  virtual void mult_cost(Scalar& ) const ;
1892 
1895  virtual void div_cost(Scalar& ) const ;
1896 
1899  virtual void mult_sint(Scalar& ) const ;
1900 
1903  virtual void div_sint(Scalar& ) const ;
1904 
1907  virtual void div_tant(Scalar& ) const ;
1908 
1918  virtual void comp_x_from_spherical(const Scalar& v_r, const Scalar& v_theta,
1919  const Scalar& v_phi, Scalar& v_x) const ;
1920 
1922  virtual void comp_x_from_spherical(const Cmp& v_r, const Cmp& v_theta,
1923  const Cmp& v_phi, Cmp& v_x) const ;
1924 
1934  virtual void comp_y_from_spherical(const Scalar& v_r, const Scalar& v_theta,
1935  const Scalar& v_phi, Scalar& v_y) const ;
1936 
1938  virtual void comp_y_from_spherical(const Cmp& v_r, const Cmp& v_theta,
1939  const Cmp& v_phi, Cmp& v_y) const ;
1940 
1949  virtual void comp_z_from_spherical(const Scalar& v_r, const Scalar& v_theta,
1950  Scalar& v_z) const ;
1951 
1953  virtual void comp_z_from_spherical(const Cmp& v_r, const Cmp& v_theta,
1954  Cmp& v_z) const ;
1955 
1965  virtual void comp_r_from_cartesian(const Scalar& v_x, const Scalar& v_y,
1966  const Scalar& v_z, Scalar& v_r) const ;
1967 
1969  virtual void comp_r_from_cartesian(const Cmp& v_x, const Cmp& v_y,
1970  const Cmp& v_z, Cmp& v_r) const ;
1971 
1981  virtual void comp_t_from_cartesian(const Scalar& v_x, const Scalar& v_y,
1982  const Scalar& v_z, Scalar& v_t) const ;
1983 
1985  virtual void comp_t_from_cartesian(const Cmp& v_x, const Cmp& v_y,
1986  const Cmp& v_z, Cmp& v_t) const ;
1987 
1996  virtual void comp_p_from_cartesian(const Scalar& v_x, const Scalar& v_y,
1997  Scalar& v_p) const ;
1998 
2000  virtual void comp_p_from_cartesian(const Cmp& v_x, const Cmp& v_y,
2001  Cmp& v_p) const ;
2002 
2008  virtual void dec_dzpuis(Scalar& ) const ;
2009 
2015  virtual void dec2_dzpuis(Scalar& ) const ;
2016 
2022  virtual void inc_dzpuis(Scalar& ) const ;
2023 
2029  virtual void inc2_dzpuis(Scalar& ) const ;
2030 
2031 
2032  // PDE resolution :
2033  // --------------
2034  public:
2057  virtual void poisson_compact(const Cmp& source, const Cmp& aa,
2058  const Tenseur& bb, const Param& par,
2059  Cmp& psi) const ;
2060 
2075  virtual void poisson_compact(int nzet, const Cmp& source, const Cmp& aa,
2076  const Tenseur& bb, const Param& par,
2077  Cmp& psi) const ;
2078 
2079 };
2080 
2081 
2082  //------------------------------------//
2083  // class Map_af //
2084  //------------------------------------//
2085 
2086 
2087 
2102 class Map_af : public Map_radial {
2103 
2104  // Data :
2105  // ----
2106  private:
2108  double* alpha ;
2110  double* beta ;
2111 
2112  // Constructors, destructor :
2113  // ------------------------
2114  public:
2126  Map_af(const Mg3d& mgrille, const double* r_limits) ;
2139  Map_af(const Mg3d& mgrille, const Tbl& r_limits) ;
2140 
2141  Map_af(const Map_af& ) ;
2142  Map_af(const Mg3d&, const string&) ;
2143  Map_af(const Mg3d&, FILE* ) ;
2144 
2156  explicit Map_af(const Map& ) ;
2157 
2158  virtual ~Map_af() ;
2159 
2160  // Assignment
2161  // ----------
2162  public:
2164  virtual void operator=(const Map_af& ) ;
2165 
2166  // Memory management
2167  // -----------------
2168  private:
2170  void set_coord() ;
2171 
2172  // Extraction of information
2173  // -------------------------
2174  public:
2176  const double* get_alpha() const ;
2177 
2179  const double* get_beta() const ;
2180 
2184  virtual const Map_af& mp_angu(int) const ;
2185 
2195  virtual double val_r(int l, double xi, double theta, double pphi) const ;
2196 
2206  virtual void val_lx(double rr, double theta, double pphi,
2207  int& l, double& xi) const ;
2208 
2218  virtual void val_lx(double rr, double theta, double pphi,
2219  const Param& par, int& l, double& xi) const ;
2220 
2230  virtual double val_r_jk(int l, double xi, int j, int k) const ;
2231 
2241  virtual void val_lx_jk(double rr, int j, int k, const Param& par,
2242  int& l, double& xi) const ;
2243 
2245  virtual bool operator==(const Map& ) const ;
2246 
2247 
2248  // Outputs
2249  // -------
2250  public:
2251  virtual void sauve(FILE* ) const ;
2252 
2253  private:
2254  virtual ostream& operator>>(ostream &) const ;
2255 
2256  // Modification of the mapping
2257  // ---------------------------
2258  public:
2263  virtual void homothetie(double lambda) ;
2264 
2274  virtual void resize(int l, double lambda) ;
2275 
2281  void homothetie_interne(double lambda) ;
2282 
2285  virtual void adapt(const Cmp& ent, const Param& par, int nbr=0) ;
2286 
2288  void set_alpha(double alpha0, int l) ;
2289 
2291  void set_beta(double beta0, int l) ;
2292 
2293  // Differential operators:
2294  // ----------------------
2295  public:
2302  virtual void dsdxi(const Cmp& ci, Cmp& resu) const ;
2303 
2310  virtual void dsdr(const Cmp& ci, Cmp& resu) const ;
2311 
2318  virtual void srdsdt(const Cmp& ci, Cmp& resu) const ;
2319 
2327  virtual void srstdsdp(const Cmp& ci, Cmp& resu) const ;
2328 
2336  virtual void dsdr(const Scalar& uu, Scalar& resu) const ;
2337 
2345  virtual void dsdxi(const Scalar& uu, Scalar& resu) const ;
2346 
2354  virtual void dsdradial(const Scalar&, Scalar&) const ;
2355 
2363  virtual void srdsdt(const Scalar& uu, Scalar& resu) const ;
2364 
2372  virtual void srstdsdp(const Scalar& uu, Scalar& resu) const ;
2373 
2378  virtual void dsdt(const Scalar& uu, Scalar& resu) const ;
2379 
2384  virtual void stdsdp(const Scalar& uu, Scalar& resu) const ;
2385 
2396  virtual void laplacien(const Scalar& uu, int zec_mult_r,
2397  Scalar& lap) const ;
2398 
2400  virtual void laplacien(const Cmp& uu, int zec_mult_r,
2401  Cmp& lap) const ;
2402 
2409  virtual void lapang(const Scalar& uu, Scalar& lap) const ;
2410 
2411 
2422  virtual void primr(const Scalar& uu, Scalar& resu,
2423  bool null_infty) const ;
2424 
2425 
2433  virtual Tbl* integrale(const Scalar&) const ;
2434 
2442  virtual Tbl* integrale(const Cmp&) const ;
2443 
2444 
2445  // PDE resolution :
2446  // --------------
2447  public:
2457  virtual void poisson(const Cmp& source, Param& par, Cmp& uu, bool verbose = true) const ;
2458 
2466  virtual void poisson_tau(const Cmp& source, Param& par, Cmp& uu) const ;
2467 
2477  virtual void poisson(const Scalar& source, Param& par, Scalar& uu, bool verbose=true) const ;
2478 
2486  virtual void poisson_tau(const Scalar& source, Param& par, Scalar& uu) const ;
2487 
2488  virtual void poisson_falloff(const Cmp& source, Param& par, Cmp& uu,
2489  int k_falloff) const ;
2490 
2491  virtual void poisson_ylm(const Cmp& source, Param& par, Cmp& pot,
2492  int nylm, double* intvec) const ;
2493 
2515  virtual void poisson_regular(const Cmp& source, int k_div, int nzet,
2516  double unsgam1, Param& par, Cmp& uu,
2517  Cmp& uu_regu, Cmp& uu_div,
2518  Tenseur& duu_div, Cmp& source_regu,
2519  Cmp& source_div) const ;
2520 
2538  virtual void poisson_angu(const Scalar& source, Param& par,
2539  Scalar& uu, double lambda=0) const ;
2540  virtual void poisson_angu(const Cmp& source, Param& par,
2541  Cmp& uu, double lambda=0) const ;
2542 
2552  virtual Param* donne_para_poisson_vect (Param& par, int i) const ;
2553 
2558  virtual void poisson_frontiere (const Cmp&, const Valeur&, int, int, Cmp&, double = 0., double = 0.) const ;
2559 
2565  virtual void poisson_frontiere_double (const Cmp& source, const Valeur& lim_func,
2566  const Valeur& lim_der, int num_zone, Cmp& pot) const ;
2567 
2578  virtual void poisson_interne (const Cmp& source, const Valeur& limite,
2579  Param& par, Cmp& pot) const ;
2580 
2585  double integrale_surface (const Cmp& ci, double rayon) const ;
2586 
2591  double integrale_surface (const Scalar& ci, double rayon) const ;
2592 
2593  double integrale_surface_falloff (const Cmp& ci) const ;
2594 
2599  double integrale_surface_infini (const Cmp& ci) const ;
2600 
2605  double integrale_surface_infini (const Scalar& ci) const ;
2606 
2614  void sol_elliptic (Param_elliptic& params,
2615  const Scalar& so, Scalar& uu) const ;
2616 
2617 
2629  void sol_elliptic_boundary (Param_elliptic& params,
2630  const Scalar& so, Scalar& uu, const Mtbl_cf& bound,
2631  double fact_dir, double fact_neu ) const ;
2632 
2637  void sol_elliptic_boundary (Param_elliptic& params,
2638  const Scalar& so, Scalar& uu, const Scalar& bound,
2639  double fact_dir, double fact_neu ) const ;
2640 
2650  void sol_elliptic_no_zec (Param_elliptic& params,
2651  const Scalar& so, Scalar& uu, double val) const ;
2652 
2662  void sol_elliptic_only_zec (Param_elliptic& params,
2663  const Scalar& so, Scalar& uu, double val) const ;
2664 
2670  void sol_elliptic_sin_zec (Param_elliptic& params,
2671  const Scalar& so, Scalar& uu,
2672  double* coefs, double*) const ;
2683  void sol_elliptic_fixe_der_zero (double val,
2684  Param_elliptic& params,
2685  const Scalar& so, Scalar& uu) const ;
2686 
2713  virtual void poisson2d(const Cmp& source_mat, const Cmp& source_quad,
2714  Param& par, Cmp& uu, bool verbose=true) const ;
2723  const Scalar&, Scalar&) const ;
2732  const Scalar&, Scalar&) const ;
2733 
2762  virtual void dalembert(Param& par, Scalar& fJp1, const Scalar& fJ,
2763  const Scalar& fJm1, const Scalar& source) const ;
2764 
2770  Scalar interpolate_from_map_star(const Scalar& f_s) const ;
2771 
2772  // Building functions for the Coord's
2773  // ----------------------------------
2774  friend Mtbl* map_af_fait_r(const Map* ) ;
2775  friend Mtbl* map_af_fait_tet(const Map* ) ;
2776  friend Mtbl* map_af_fait_phi(const Map* ) ;
2777  friend Mtbl* map_af_fait_sint(const Map* ) ;
2778  friend Mtbl* map_af_fait_cost(const Map* ) ;
2779  friend Mtbl* map_af_fait_sinp(const Map* ) ;
2780  friend Mtbl* map_af_fait_cosp(const Map* ) ;
2781 
2782  friend Mtbl* map_af_fait_x(const Map* ) ;
2783  friend Mtbl* map_af_fait_y(const Map* ) ;
2784  friend Mtbl* map_af_fait_z(const Map* ) ;
2785 
2786  friend Mtbl* map_af_fait_xa(const Map* ) ;
2787  friend Mtbl* map_af_fait_ya(const Map* ) ;
2788  friend Mtbl* map_af_fait_za(const Map* ) ;
2789 
2790  friend Mtbl* map_af_fait_xsr(const Map* ) ;
2791  friend Mtbl* map_af_fait_dxdr(const Map* ) ;
2792  friend Mtbl* map_af_fait_drdt(const Map* ) ;
2793  friend Mtbl* map_af_fait_stdrdp(const Map* ) ;
2794  friend Mtbl* map_af_fait_srdrdt(const Map* ) ;
2795  friend Mtbl* map_af_fait_srstdrdp(const Map* ) ;
2796  friend Mtbl* map_af_fait_sr2drdt(const Map* ) ;
2797  friend Mtbl* map_af_fait_sr2stdrdp(const Map* ) ;
2798  friend Mtbl* map_af_fait_d2rdx2(const Map* ) ;
2799  friend Mtbl* map_af_fait_lapr_tp(const Map* ) ;
2800  friend Mtbl* map_af_fait_d2rdtdx(const Map* ) ;
2801  friend Mtbl* map_af_fait_sstd2rdpdx(const Map* ) ;
2802  friend Mtbl* map_af_fait_sr2d2rdt2(const Map* ) ;
2803 
2804 };
2805 
2806  Mtbl* map_af_fait_r(const Map* ) ;
2807  Mtbl* map_af_fait_tet(const Map* ) ;
2808  Mtbl* map_af_fait_phi(const Map* ) ;
2809  Mtbl* map_af_fait_sint(const Map* ) ;
2810  Mtbl* map_af_fait_cost(const Map* ) ;
2811  Mtbl* map_af_fait_sinp(const Map* ) ;
2812  Mtbl* map_af_fait_cosp(const Map* ) ;
2813 
2814  Mtbl* map_af_fait_x(const Map* ) ;
2815  Mtbl* map_af_fait_y(const Map* ) ;
2816  Mtbl* map_af_fait_z(const Map* ) ;
2817 
2818  Mtbl* map_af_fait_xa(const Map* ) ;
2819  Mtbl* map_af_fait_ya(const Map* ) ;
2820  Mtbl* map_af_fait_za(const Map* ) ;
2821 
2822  Mtbl* map_af_fait_xsr(const Map* ) ;
2823  Mtbl* map_af_fait_dxdr(const Map* ) ;
2824  Mtbl* map_af_fait_drdt(const Map* ) ;
2825  Mtbl* map_af_fait_stdrdp(const Map* ) ;
2826  Mtbl* map_af_fait_srdrdt(const Map* ) ;
2827  Mtbl* map_af_fait_srstdrdp(const Map* ) ;
2828  Mtbl* map_af_fait_sr2drdt(const Map* ) ;
2829  Mtbl* map_af_fait_sr2stdrdp(const Map* ) ;
2830  Mtbl* map_af_fait_d2rdx2(const Map* ) ;
2831  Mtbl* map_af_fait_lapr_tp(const Map* ) ;
2832  Mtbl* map_af_fait_d2rdtdx(const Map* ) ;
2833  Mtbl* map_af_fait_sstd2rdpdx(const Map* ) ;
2834  Mtbl* map_af_fait_sr2d2rdt2(const Map* ) ;
2835 
2836 
2837 
2838 
2839  //------------------------------------//
2840  // class Map_et //
2841  //------------------------------------//
2842 
2843 
2844 
2870 class Map_et : public Map_radial {
2871 
2872  // Data :
2873  // ----
2874  private:
2876  double* alpha ;
2878  double* beta ;
2879 
2883  Tbl** aa ;
2884 
2888  Tbl** daa ;
2889 
2893  Tbl** ddaa ;
2894 
2897 
2900 
2905 
2910 
2914  Tbl** bb ;
2915 
2919  Tbl** dbb ;
2920 
2924  Tbl** ddbb ;
2925 
2928 
2931 
2938 
2945 
2946  public:
2953 
2960 
2961  // Constructors, destructor :
2962  // ------------------------
2963  public:
2975  Map_et(const Mg3d& mgrille, const double* r_limits) ;
2976 
2994  Map_et(const Mg3d& mgrille, const double* r_limits,const Tbl& tab);
2995  Map_et(const Map_et& ) ;
2996  Map_et(const Mg3d&, FILE* ) ;
2997 
2998  virtual ~Map_et() ;
2999 
3000  // Assignment
3001  // ----------
3002  public:
3004  virtual void operator=(const Map_et& mp) ;
3005 
3007  virtual void operator=(const Map_af& mpa) ;
3008 
3010  void set_ff(const Valeur& ) ;
3012  void set_gg(const Valeur& ) ;
3013 
3014  // Memory management
3015  // -----------------
3016  private:
3018  void set_coord() ;
3019  protected:
3021  virtual void reset_coord() ;
3022 
3023  private:
3025  void fait_poly() ;
3026 
3027  // Extraction of information
3028  // -------------------------
3029  public:
3033  virtual const Map_af& mp_angu(int) const ;
3034 
3038  const double* get_alpha() const ;
3039 
3043  const double* get_beta() const ;
3044 
3046  const Valeur& get_ff() const ;
3047 
3049  const Valeur& get_gg() const ;
3050 
3060  virtual double val_r(int l, double xi, double theta, double pphi) const ;
3061 
3071  virtual void val_lx(double rr, double theta, double pphi,
3072  int& l, double& xi) const ;
3073 
3092  virtual void val_lx(double rr, double theta, double pphi,
3093  const Param& par, int& l, double& xi) const ;
3094 
3096  virtual bool operator==(const Map& ) const ;
3097 
3107  virtual double val_r_jk(int l, double xi, int j, int k) const ;
3108 
3125  virtual void val_lx_jk(double rr, int j, int k, const Param& par,
3126  int& l, double& xi) const ;
3127 
3128 
3129 
3130  // Outputs
3131  // -------
3132  public:
3133  virtual void sauve(FILE* ) const ;
3134 
3135  private:
3136  virtual ostream& operator>>(ostream &) const ;
3137 
3138  // Modification of the radial scale
3139  // --------------------------------
3140  public:
3145  virtual void homothetie(double lambda) ;
3146 
3156  virtual void resize(int l, double lambda) ;
3157 
3164  void resize_extr(double lambda) ;
3165 
3167  void set_alpha(double alpha0, int l) ;
3168 
3170  void set_beta(double beta0, int l) ;
3171 
3172  // Modification of the mapping
3173  // ---------------------------
3218  virtual void adapt(const Cmp& ent, const Param& par, int nbr_filtre = 0) ;
3219 
3220  // Differential operators:
3221  // ----------------------
3222  public:
3229  virtual void dsdxi(const Cmp& ci, Cmp& resu) const ;
3230 
3237  virtual void dsdr(const Cmp& ci, Cmp& resu) const ;
3238 
3245  virtual void srdsdt(const Cmp& ci, Cmp& resu) const ;
3246 
3254  virtual void srstdsdp(const Cmp& ci, Cmp& resu) const ;
3255 
3263  virtual void dsdxi(const Scalar& uu, Scalar& resu) const ;
3264 
3272  virtual void dsdr(const Scalar& uu, Scalar& resu) const ;
3273 
3282  virtual void dsdradial(const Scalar& uu, Scalar& resu) const ;
3283 
3291  virtual void srdsdt(const Scalar& uu, Scalar& resu) const ;
3292 
3300  virtual void srstdsdp(const Scalar& uu, Scalar& resu) const ;
3301 
3306  virtual void dsdt(const Scalar& uu, Scalar& resu) const ;
3307 
3312  virtual void stdsdp(const Scalar& uu, Scalar& resu) const ;
3313 
3324  virtual void laplacien(const Scalar& uu, int zec_mult_r,
3325  Scalar& lap) const ;
3326 
3328  virtual void laplacien(const Cmp& uu, int zec_mult_r,
3329  Cmp& lap) const ;
3330 
3337  virtual void lapang(const Scalar& uu, Scalar& lap) const ;
3338 
3339 
3350  virtual void primr(const Scalar& uu, Scalar& resu,
3351  bool null_infty) const ;
3352 
3353 
3361  virtual Tbl* integrale(const Scalar&) const ;
3362 
3370  virtual Tbl* integrale(const Cmp&) const ;
3371 
3372 
3373  // PDE resolution :
3374  // --------------
3375  public:
3415  virtual void poisson(const Cmp& source, Param& par, Cmp& uu, bool verbose=true) const ;
3416 
3454  virtual void poisson_tau(const Cmp& source, Param& par, Cmp& uu) const ;
3455 
3495  virtual void poisson(const Scalar& source, Param& par, Scalar& uu, bool verbose=true) const ;
3496 
3534  virtual void poisson_tau(const Scalar& source, Param& par, Scalar& uu) const ;
3535 
3536  virtual void poisson_falloff(const Cmp& source, Param& par, Cmp& uu,
3537  int k_falloff) const ;
3538 
3539  virtual void poisson_ylm(const Cmp& source, Param& par, Cmp& uu,
3540  int nylm, double* intvec) const ;
3541 
3577  virtual void poisson_regular(const Cmp& source, int k_div, int nzet,
3578  double unsgam1, Param& par, Cmp& uu,
3579  Cmp& uu_regu, Cmp& uu_div,
3580  Tenseur& duu_div, Cmp& source_regu,
3581  Cmp& source_div) const ;
3582 
3600  virtual void poisson_angu(const Scalar& source, Param& par,
3601  Scalar& uu, double lambda=0) const ;
3602  virtual void poisson_angu(const Cmp& source, Param& par,
3603  Cmp& uu, double lambda=0) const ;
3604 
3628  virtual Param* donne_para_poisson_vect (Param& para, int i) const ;
3629 
3633  virtual void poisson_frontiere (const Cmp&, const Valeur&, int, int,
3634  Cmp&, double = 0., double = 0.) const ;
3635  virtual void poisson_frontiere_double (const Cmp& source,
3636  const Valeur& lim_func, const Valeur& lim_der,
3637  int num_zone, Cmp& pot) const ;
3638 
3649  virtual void poisson_interne (const Cmp& source, const Valeur& limite,
3650  Param& par, Cmp& pot) const ;
3651 
3652 
3691  virtual void poisson2d(const Cmp& source_mat, const Cmp& source_quad,
3692  Param& par, Cmp& uu, bool verbose=true) const ;
3693 
3697  virtual void dalembert(Param& par, Scalar& fJp1, const Scalar& fJ,
3698  const Scalar& fJm1, const Scalar& source) const ;
3699 
3700 
3701 
3702 
3703  // Building functions for the Coord's
3704  // ----------------------------------
3705  friend Mtbl* map_et_fait_r(const Map* ) ;
3706  friend Mtbl* map_et_fait_tet(const Map* ) ;
3707  friend Mtbl* map_et_fait_phi(const Map* ) ;
3708  friend Mtbl* map_et_fait_sint(const Map* ) ;
3709  friend Mtbl* map_et_fait_cost(const Map* ) ;
3710  friend Mtbl* map_et_fait_sinp(const Map* ) ;
3711  friend Mtbl* map_et_fait_cosp(const Map* ) ;
3712 
3713  friend Mtbl* map_et_fait_x(const Map* ) ;
3714  friend Mtbl* map_et_fait_y(const Map* ) ;
3715  friend Mtbl* map_et_fait_z(const Map* ) ;
3716 
3717  friend Mtbl* map_et_fait_xa(const Map* ) ;
3718  friend Mtbl* map_et_fait_ya(const Map* ) ;
3719  friend Mtbl* map_et_fait_za(const Map* ) ;
3720 
3721  friend Mtbl* map_et_fait_xsr(const Map* ) ;
3722  friend Mtbl* map_et_fait_dxdr(const Map* ) ;
3723  friend Mtbl* map_et_fait_drdt(const Map* ) ;
3724  friend Mtbl* map_et_fait_stdrdp(const Map* ) ;
3725  friend Mtbl* map_et_fait_srdrdt(const Map* ) ;
3726  friend Mtbl* map_et_fait_srstdrdp(const Map* ) ;
3727  friend Mtbl* map_et_fait_sr2drdt(const Map* ) ;
3728  friend Mtbl* map_et_fait_sr2stdrdp(const Map* ) ;
3729  friend Mtbl* map_et_fait_d2rdx2(const Map* ) ;
3730  friend Mtbl* map_et_fait_lapr_tp(const Map* ) ;
3731  friend Mtbl* map_et_fait_d2rdtdx(const Map* ) ;
3732  friend Mtbl* map_et_fait_sstd2rdpdx(const Map* ) ;
3733  friend Mtbl* map_et_fait_sr2d2rdt2(const Map* ) ;
3734 
3735  friend Mtbl* map_et_fait_rsxdxdr(const Map* ) ;
3736  friend Mtbl* map_et_fait_rsx2drdx(const Map* ) ;
3737 
3738 };
3739 
3740  Mtbl* map_et_fait_r(const Map* ) ;
3741  Mtbl* map_et_fait_tet(const Map* ) ;
3742  Mtbl* map_et_fait_phi(const Map* ) ;
3743  Mtbl* map_et_fait_sint(const Map* ) ;
3744  Mtbl* map_et_fait_cost(const Map* ) ;
3745  Mtbl* map_et_fait_sinp(const Map* ) ;
3746  Mtbl* map_et_fait_cosp(const Map* ) ;
3747 
3748  Mtbl* map_et_fait_x(const Map* ) ;
3749  Mtbl* map_et_fait_y(const Map* ) ;
3750  Mtbl* map_et_fait_z(const Map* ) ;
3751 
3752  Mtbl* map_et_fait_xa(const Map* ) ;
3753  Mtbl* map_et_fait_ya(const Map* ) ;
3754  Mtbl* map_et_fait_za(const Map* ) ;
3755 
3756  Mtbl* map_et_fait_xsr(const Map* ) ;
3757  Mtbl* map_et_fait_dxdr(const Map* ) ;
3758  Mtbl* map_et_fait_drdt(const Map* ) ;
3759  Mtbl* map_et_fait_stdrdp(const Map* ) ;
3760  Mtbl* map_et_fait_srdrdt(const Map* ) ;
3761  Mtbl* map_et_fait_srstdrdp(const Map* ) ;
3762  Mtbl* map_et_fait_sr2drdt(const Map* ) ;
3763  Mtbl* map_et_fait_sr2stdrdp(const Map* ) ;
3764  Mtbl* map_et_fait_d2rdx2(const Map* ) ;
3765  Mtbl* map_et_fait_lapr_tp(const Map* ) ;
3766  Mtbl* map_et_fait_d2rdtdx(const Map* ) ;
3767  Mtbl* map_et_fait_sstd2rdpdx(const Map* ) ;
3768  Mtbl* map_et_fait_sr2d2rdt2(const Map* ) ;
3769 
3770  Mtbl* map_et_fait_rsxdxdr(const Map* ) ;
3771  Mtbl* map_et_fait_rsx2drdx(const Map* ) ;
3772 
3773  //------------------------------------//
3774  // class Map_log //
3775  //------------------------------------//
3776 
3777 #define AFFINE 0
3778 #define LOG 1
3779 
3794 class Map_log : public Map_radial {
3795 
3796  // Data :
3797  // ----
3798  private:
3807 
3808  public:
3815 
3816  private:
3817  void set_coord() ;
3818 
3819  // Constructors, destructor :
3820  // ------------------------
3821  public:
3834  Map_log (const Mg3d& mgrille, const Tbl& r_limits, const Itbl& typevar) ;
3835 
3836 
3837  Map_log (const Map_log& ) ;
3838  Map_log (const Mg3d&, FILE* ) ;
3839 
3840  virtual ~Map_log() ;
3841 
3845  virtual const Map_af& mp_angu(int) const ;
3846 
3848  double get_alpha (int l) const {return alpha(l) ;} ;
3850  double get_beta (int l) const {return beta(l) ;} ;
3852  int get_type (int l) const {return type_var(l) ;} ;
3853 
3861  void sol_elliptic (Param_elliptic& params,
3862  const Scalar& so, Scalar& uu) const ;
3863 
3864 
3876  void sol_elliptic_boundary (Param_elliptic& params,
3877  const Scalar& so, Scalar& uu, const Mtbl_cf& bound,
3878  double fact_dir, double fact_neu ) const ;
3879 
3883  void sol_elliptic_boundary (Param_elliptic& params,
3884  const Scalar& so, Scalar& uu, const Scalar& bound,
3885  double fact_dir, double fact_neu ) const ;
3886 
3887 
3897  void sol_elliptic_no_zec (Param_elliptic& params,
3898  const Scalar& so, Scalar& uu, double) const ;
3899 
3900 
3901  virtual void sauve(FILE*) const ;
3902 
3904  virtual void operator=(const Map_af& mpa) ;
3905 
3906 
3907  virtual ostream& operator>> (ostream&) const ;
3908 
3918  virtual double val_r (int l, double xi, double theta, double pphi) const ;
3919 
3929  virtual void val_lx (double rr, double theta, double pphi, int& l, double& xi) const ;
3930 
3940  virtual void val_lx (double rr, double theta, double pphi, const Param& par, int& l, double& xi) const ;
3941 
3942 
3943  virtual bool operator== (const Map&) const ;
3944 
3954  virtual double val_r_jk (int l, double xi, int j, int k) const ;
3955 
3965  virtual void val_lx_jk (double rr, int j, int k, const Param& par, int& l, double& xi) const ;
3966 
3973  virtual void dsdr (const Scalar& ci, Scalar& resu) const ;
3974 
3981  virtual void dsdxi (const Scalar& ci, Scalar& resu) const ;
3982 
3991  virtual void dsdradial (const Scalar& uu, Scalar& resu) const ;
3992 
3993  virtual void homothetie (double) ;
3994  virtual void resize (int, double) ;
3995  virtual void adapt (const Cmp&, const Param&, int) ;
3996  virtual void dsdr (const Cmp&, Cmp&) const ;
3997  virtual void dsdxi (const Cmp&, Cmp&) const ;
3998  virtual void srdsdt (const Cmp&, Cmp&) const ;
3999  virtual void srstdsdp (const Cmp&, Cmp&) const ;
4000  virtual void srstdsdp (const Scalar&, Scalar&) const ;
4001  virtual void srdsdt (const Scalar&, Scalar&) const ;
4002  virtual void dsdt (const Scalar&, Scalar&) const ;
4003  virtual void stdsdp (const Scalar&, Scalar&) const ;
4004  virtual void laplacien (const Scalar&, int, Scalar&) const ;
4005  virtual void laplacien (const Cmp&, int, Cmp&) const ;
4006  virtual void lapang (const Scalar&, Scalar&) const ;
4007  virtual void primr(const Scalar&, Scalar&, bool) const ;
4008  virtual Tbl* integrale (const Scalar&) const ;
4009  virtual Tbl* integrale (const Cmp&) const ;
4010  virtual void poisson (const Cmp&, Param&, Cmp&, bool) const ;
4011  virtual void poisson (const Scalar&, Param&, Scalar&, bool) const ;
4012  virtual void poisson_tau (const Cmp&, Param&, Cmp&) const ;
4013  virtual void poisson_tau (const Scalar&, Param&, Scalar&) const ;
4014  virtual void poisson_falloff(const Cmp&, Param&, Cmp&, int) const ;
4015  virtual void poisson_ylm(const Cmp&, Param&, Cmp&, int, double*) const ;
4016  virtual void poisson_regular (const Cmp&, int, int, double, Param&, Cmp&, Cmp&, Cmp&,
4017  Tenseur&, Cmp&, Cmp&) const ;
4018  virtual void poisson_angu (const Scalar&, Param&, Scalar&, double=0) const ;
4019  virtual void poisson_angu (const Cmp&, Param&, Cmp&, double=0) const ;
4020  virtual Param* donne_para_poisson_vect (Param&, int) const ;
4021  virtual void poisson_frontiere (const Cmp&, const Valeur&, int, int, Cmp&, double = 0., double = 0.) const ;
4022  virtual void poisson_frontiere_double (const Cmp&, const Valeur&, const Valeur&, int, Cmp&) const ;
4023  virtual void poisson_interne (const Cmp&, const Valeur&, Param&, Cmp&) const ;
4024  virtual void poisson2d (const Cmp&, const Cmp&, Param&, Cmp&, bool=true) const ;
4025  virtual void dalembert (Param&, Scalar&, const Scalar&, const Scalar&, const Scalar&) const ;
4026 
4027 
4028  // Building functions for the Coord's
4029  // ----------------------------------
4030  friend Mtbl* map_log_fait_r(const Map* ) ;
4031  friend Mtbl* map_log_fait_tet(const Map* ) ;
4032  friend Mtbl* map_log_fait_phi(const Map* ) ;
4033  friend Mtbl* map_log_fait_sint(const Map* ) ;
4034  friend Mtbl* map_log_fait_cost(const Map* ) ;
4035  friend Mtbl* map_log_fait_sinp(const Map* ) ;
4036  friend Mtbl* map_log_fait_cosp(const Map* ) ;
4037 
4038  friend Mtbl* map_log_fait_x(const Map* ) ;
4039  friend Mtbl* map_log_fait_y(const Map* ) ;
4040  friend Mtbl* map_log_fait_z(const Map* ) ;
4041 
4042  friend Mtbl* map_log_fait_xa(const Map* ) ;
4043  friend Mtbl* map_log_fait_ya(const Map* ) ;
4044  friend Mtbl* map_log_fait_za(const Map* ) ;
4045 
4046  friend Mtbl* map_log_fait_xsr(const Map* ) ;
4047  friend Mtbl* map_log_fait_dxdr(const Map* ) ;
4048  friend Mtbl* map_log_fait_drdt(const Map* ) ;
4049  friend Mtbl* map_log_fait_stdrdp(const Map* ) ;
4050  friend Mtbl* map_log_fait_srdrdt(const Map* ) ;
4051  friend Mtbl* map_log_fait_srstdrdp(const Map* ) ;
4052  friend Mtbl* map_log_fait_sr2drdt(const Map* ) ;
4053  friend Mtbl* map_log_fait_sr2stdrdp(const Map* ) ;
4054  friend Mtbl* map_log_fait_d2rdx2(const Map* ) ;
4055  friend Mtbl* map_log_fait_lapr_tp(const Map* ) ;
4056  friend Mtbl* map_log_fait_d2rdtdx(const Map* ) ;
4057  friend Mtbl* map_log_fait_sstd2rdpdx(const Map* ) ;
4058  friend Mtbl* map_log_fait_sr2d2rdt2(const Map* ) ;
4059  friend Mtbl* map_log_fait_dxdlnr(const Map* ) ;
4060 
4061 };
4062 
4063 Mtbl* map_log_fait_r(const Map* ) ;
4064 Mtbl* map_log_fait_tet(const Map* ) ;
4065 Mtbl* map_log_fait_phi(const Map* ) ;
4066 Mtbl* map_log_fait_sint(const Map* ) ;
4067 Mtbl* map_log_fait_cost(const Map* ) ;
4068 Mtbl* map_log_fait_sinp(const Map* ) ;
4069 Mtbl* map_log_fait_cosp(const Map* ) ;
4070 
4071 Mtbl* map_log_fait_x(const Map* ) ;
4072 Mtbl* map_log_fait_y(const Map* ) ;
4073 Mtbl* map_log_fait_z(const Map* ) ;
4074 
4075 Mtbl* map_log_fait_xa(const Map* ) ;
4076 Mtbl* map_log_fait_ya(const Map* ) ;
4077 Mtbl* map_log_fait_za(const Map* ) ;
4078 
4079 Mtbl* map_log_fait_xsr(const Map* ) ;
4080 Mtbl* map_log_fait_dxdr(const Map* ) ;
4081 Mtbl* map_log_fait_drdt(const Map* ) ;
4082 Mtbl* map_log_fait_stdrdp(const Map* ) ;
4083 Mtbl* map_log_fait_srdrdt(const Map* ) ;
4084 Mtbl* map_log_fait_srstdrdp(const Map* ) ;
4085 Mtbl* map_log_fait_sr2drdt(const Map* ) ;
4086 Mtbl* map_log_fait_sr2stdrdp(const Map* ) ;
4087 Mtbl* map_log_fait_d2rdx2(const Map* ) ;
4088 Mtbl* map_log_fait_lapr_tp(const Map* ) ;
4089 Mtbl* map_log_fait_d2rdtdx(const Map* ) ;
4090 Mtbl* map_log_fait_sstd2rdpdx(const Map* ) ;
4091 Mtbl* map_log_fait_sr2d2rdt2(const Map* ) ;
4092 
4093 Mtbl* map_log_fait_dxdlnr (const Map*) ;
4094 
4095 }
4096 
4097  //------------------------//
4098  // Remaining Mappings //
4099  //------------------------//
4100 
4101 #include "map_star.h"
4102 
4103 #endif
Coord xa
Absolute x coordinate.
Definition: map.h:756
virtual double val_r(int l, double xi, double theta, double pphi) const
Returns the value of the radial coordinate r for a given in a given domain.
void set_coord()
Assignment of the building functions to the member Coords.
Definition: map_af.C:533
virtual void sauve(FILE *) const
Save in a file.
Definition: map_af.C:619
Tbl beta
Array (size: mg->nzone ) of the values of in each domain.
Definition: map.h:3802
virtual bool operator==(const Map &) const
Comparison operator (egality)
Definition: map_log.C:175
Coord d2rdx2
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1694
virtual Tbl * integrale(const Scalar &) const
< Not implemented
Tbl ** bb
Array (size: mg->nzone ) of Tbl which stores the values of in each domain.
Definition: map.h:2914
Coord sr2d2rdt2
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1732
Coord sr2stdrdp
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1683
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:449
virtual void stdsdp(const Scalar &, Scalar &) const
< Not implemented
Tbl zaasx
Values at the nr collocation points of in the outermost compactified domain.
Definition: map.h:2904
virtual void srstdsdp(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_af_deriv.C:636
const double * get_alpha() const
Returns the pointer on the array alpha.
Definition: map_af.C:607
virtual void primr(const Scalar &uu, Scalar &resu, bool null_infty) const
Computes the radial primitive which vanishes for .
Definition: map_af_primr.C:86
void fait_poly()
Construction of the polynomials and .
Definition: map_et.C:668
double integrale_surface(const Cmp &ci, double rayon) const
Performs the surface integration of ci on the sphere of radius rayon .
void set_gg(const Valeur &)
Assigns a given value to the function .
Definition: map_et.C:593
Radial mapping of rather general form.
Definition: map.h:2870
virtual const Map_af & mp_angu(int) const =0
Returns the "angular" mapping for the outside of domain l_zone.
virtual void primr(const Scalar &uu, Scalar &resu, bool null_infty) const
Computes the radial primitive which vanishes for .
Definition: map_et_integ.C:116
const double * get_alpha() const
Returns a pointer on the array alpha (values of in each domain)
Definition: map_et.C:1049
virtual void mult_sint(Scalar &) const =0
Multiplication by of a Scalar.
virtual void sauve(FILE *) const
Save in a file.
Definition: map.C:230
virtual void mult_cost(Scalar &) const
Multiplication by of a Scalar.
const Cmp & cmp_zero() const
Returns the null Cmp defined on *this.
Definition: map.h:833
const double * get_beta() const
Returns a pointer on the array beta (values of in each domain)
Definition: map_et.C:1053
virtual void dsdt(const Scalar &uu, Scalar &resu) const
Computes of a Scalar.
Definition: map_et_deriv.C:632
virtual void dsdxi(const Cmp &ci, Cmp &resu) const =0
Computes of a Cmp .
virtual void val_lx(double rr, double theta, double pphi, int &l, double &xi) const
Computes the domain index l and the value of corresponding to a point given by its physical coordina...
void set_rot_phi(double phi0)
Sets a new rotation angle.
Definition: map.C:269
virtual ostream & operator>>(ostream &) const =0
Operator >>
double integrale_surface_infini(const Cmp &ci) const
Performs the surface integration of ci at infinity.
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition: map_af.C:667
Multi-domain array.
Definition: mtbl.h:118
double * alpha
Array (size: mg->nzone ) of the values of in each domain.
Definition: map.h:2108
virtual void div_tant(Scalar &) const
Division by of a Scalar.
double get_ori_y() const
Returns the y coordinate of the origin.
Definition: map.h:796
const Base_vect_spher & get_bvect_spher() const
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:809
double * alpha
Array (size: mg->nzone ) of the values of in each domain.
Definition: map.h:2876
virtual void dsdt(const Scalar &uu, Scalar &resu) const =0
Computes of a Scalar .
virtual void val_lx_jk(double rr, int j, int k, const Param &par, int &l, double &xi) const
Computes the domain index l and the value of corresponding to a point of arbitrary r but collocation...
double get_alpha(int l) const
Returns in the domain l.
Definition: map.h:3848
Lorene prototypes.
Definition: app_hor.h:67
virtual void dsdradial(const Scalar &, Scalar &) const
Computes of a Scalar.
Definition: map_af_deriv.C:420
void sol_elliptic_only_zec(Param_elliptic &params, const Scalar &so, Scalar &uu, double val) const
General elliptic solver.
Cmp * p_cmp_zero
The null Cmp.
Definition: map.h:739
virtual void srdsdt(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_af_deriv.C:478
virtual void comp_p_from_cartesian(const Scalar &v_x, const Scalar &v_y, Scalar &v_p) const =0
Computes the Spherical component (with respect to bvect_spher ) of a vector given by its cartesian c...
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:791
virtual double val_r(int l, double xi, double theta, double pphi) const
Returns the value of the radial coordinate r for a given in a given domain.
Definition: map_et_radius.C:95
virtual void dsdr(const Cmp &ci, Cmp &resu) const =0
Computes of a Cmp .
Flat metric for tensor calculation.
Definition: metric.h:261
Coord sr2drdt
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1675
int get_type(int l) const
Returns the type of description in the domain l.
Definition: map.h:3852
Tbl bbsx
Values at the nr collocation points of in the nucleus.
Definition: map.h:2927
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:402
virtual void srstdsdp(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_et_deriv.C:488
virtual Param * donne_para_poisson_vect(Param &, int) const
< Not implemented
virtual void primr(const Scalar &uu, Scalar &resu, bool null_infty) const =0
Computes the radial primitive which vanishes for .
Base class for coordinate mappings.
Definition: map.h:696
double get_ori_x() const
Returns the x coordinate of the origin.
Definition: map.h:794
const Valeur & get_gg() const
Returns a (constant) reference to the function .
Definition: map_et.C:1061
void sol_elliptic_2d(Param_elliptic &, const Scalar &, Scalar &) const
General elliptic solver in a 2D case.
virtual void poisson_frontiere_double(const Cmp &source, const Valeur &lim_func, const Valeur &lim_der, int num_zone, Cmp &pot) const
Solver of the Poisson equation with boundary condition for the affine mapping case, cases with boundary conditions of both Dirichlet and Neumann type (no condition imposed at infinity).
virtual void div_tant(Scalar &) const =0
Division by of a Scalar.
virtual void dec_dzpuis(Scalar &) const
Decreases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactified...
virtual double val_r(int l, double xi, double theta, double pphi) const
Returns the value of the radial coordinate r for a given in a given domain.
Definition: map_af_radius.C:99
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)=0
Adaptation of the mapping to a given scalar field.
virtual void comp_y_from_spherical(const Scalar &v_r, const Scalar &v_theta, const Scalar &v_phi, Scalar &v_y) const
Computes the Cartesian y component (with respect to bvect_cart ) of a vector given by its spherical c...
virtual const Map_af & mp_angu(int) const
Returns the "angular" mapping for the outside of domain l_zone.
Definition: map_af.C:786
Basic integer array class.
Definition: itbl.h:122
virtual void poisson_regular(const Cmp &source, int k_div, int nzet, double unsgam1, Param &par, Cmp &uu, Cmp &uu_regu, Cmp &uu_div, Tenseur &duu_div, Cmp &source_regu, Cmp &source_div) const
Computes the solution of a scalar Poisson equation.
virtual void poisson_compact(const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const =0
Resolution of the elliptic equation in the case where the stellar interior is covered by a single do...
virtual void laplacien(const Scalar &, int, Scalar &) const
< Not implemented
virtual void dsdxi(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_et_deriv.C:98
virtual void inc_dzpuis(Scalar &) const =0
Increases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactifi...
virtual void dec2_dzpuis(Scalar &) const
Decreases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the compactifi...
virtual ostream & operator>>(ostream &) const
Operator >>
Definition: map_log.C:207
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
virtual void inc_dzpuis(Scalar &) const
Increases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactifi...
virtual void operator=(const Map_af &)
Assignment to another affine mapping.
Definition: map_af.C:510
virtual bool operator==(const Map &) const =0
Comparison operator (egality)
virtual void mult_rsint(Scalar &) const =0
Multiplication by of a Scalar.
Map_af * p_mp_angu
Pointer on the "angular" mapping.
Definition: map.h:741
virtual bool operator==(const Map &) const
Comparison operator (egality)
Definition: map_af.C:569
virtual void val_lx(double rr, double theta, double pphi, int &l, double &xi) const
Computes the domain index l and the value of corresponding to a point given by its physical coordina...
double get_rot_phi() const
Returns the angle between the x –axis and X –axis.
Definition: map.h:801
void set_new_grid(const Mg3d &new_mg)
Sets a new grid to a Map.
Definition: map.C:348
Tbl ** ddbb
Array (size: mg->nzone ) of Tbl which stores the values of in each domain.
Definition: map.h:2924
Itbl type_var
Array (size: mg->nzone ) of the type of variable in each domain.
Definition: map.h:3806
virtual void comp_x_from_spherical(const Scalar &v_r, const Scalar &v_theta, const Scalar &v_phi, Scalar &v_x) const =0
Computes the Cartesian x component (with respect to bvect_cart ) of a vector given by its spherical c...
virtual void srdsdt(const Cmp &ci, Cmp &resu) const =0
Computes of a Cmp .
const Metric_flat & flat_met_cart() const
Returns the flat metric associated with the Cartesian coordinates and with components expressed in th...
Definition: map.C:337
Tbl ** daa
Array (size: mg->nzone ) of Tbl which stores the values of in each domain.
Definition: map.h:2888
virtual void dalembert(Param &par, Scalar &fJp1, const Scalar &fJ, const Scalar &fJm1, const Scalar &source) const
Performs one time-step integration of the d&#39;Alembert scalar equation.
Coord srstdrdp
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1667
virtual Tbl * integrale(const Scalar &) const
Computes the integral over all space of a Scalar.
Definition: map_af_integ.C:89
virtual void div_r_zec(Scalar &) const =0
Division by r (in the compactified external domain only) of a Scalar.
virtual double val_r_jk(int l, double xi, int j, int k) const
Returns the value of the radial coordinate r for a given and a given collocation point in in a give...
virtual void poisson_tau(const Cmp &source, Param &par, Cmp &uu) const =0
Computes the solution of a scalar Poisson equationwith a Tau method (Cmp version).
Coord tet
coordinate centered on the grid
Definition: map.h:745
virtual void srstdsdp(const Cmp &, Cmp &) const
< Not implemented
virtual void stdsdp(const Scalar &uu, Scalar &resu) const
Computes of a Scalar.
Definition: map_af_deriv.C:826
virtual void dsdt(const Scalar &, Scalar &) const
< Not implemented
virtual ostream & operator>>(ostream &) const
Operator >>
Definition: map_af.C:633
void set_ori(double xa0, double ya0, double za0)
Sets a new origin.
Definition: map.C:259
virtual void stdsdp(const Scalar &uu, Scalar &resu) const
Computes of a Scalar.
Definition: map_et_deriv.C:677
Coord phi
coordinate centered on the grid
Definition: map.h:746
void set_beta(double beta0, int l)
Modifies the value of in domain no. l.
Definition: map_af.C:771
virtual void div_r(Scalar &) const
Division by r of a Scalar.
virtual void comp_z_from_spherical(const Scalar &v_r, const Scalar &v_theta, Scalar &v_z) const
Computes the Cartesian z component (with respect to bvect_cart ) of a vector given by its spherical c...
virtual void poisson2d(const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu, bool verbose=true) const =0
Computes the solution of a 2-D Poisson equation.
virtual void dsdxi(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_af_deriv.C:139
virtual void poisson_interne(const Cmp &source, const Valeur &limite, Param &par, Cmp &pot) const =0
Computes the solution of a Poisson equation in the shell, imposing a boundary condition at the surfac...
virtual Param * donne_para_poisson_vect(Param &para, int i) const
Internal function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara ...
Coord sint
Definition: map.h:747
virtual void dec2_dzpuis(Scalar &) const =0
Decreases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the compactifi...
virtual void comp_t_from_cartesian(const Scalar &v_x, const Scalar &v_y, const Scalar &v_z, Scalar &v_t) const
Computes the Spherical component (with respect to bvect_spher ) of a vector given by its cartesian c...
Coord dxdlnr
Same as dxdr if the domains where the description is affine and where it is logarithmic.
Definition: map.h:3814
virtual ~Map_af()
Destructor.
Definition: map_af.C:498
virtual void comp_r_from_cartesian(const Scalar &v_x, const Scalar &v_y, const Scalar &v_z, Scalar &v_r) const
Computes the Spherical r component (with respect to bvect_spher ) of a vector given by its cartesian ...
virtual void comp_t_from_cartesian(const Scalar &v_x, const Scalar &v_y, const Scalar &v_z, Scalar &v_t) const =0
Computes the Spherical component (with respect to bvect_spher ) of a vector given by its cartesian c...
virtual void poisson_frontiere(const Cmp &, const Valeur &, int, int, Cmp &, double=0., double=0.) const
Not yet implemented.
virtual double val_r(int l, double xi, double theta, double pphi) const =0
Returns the value of the radial coordinate r for a given in a given domain.
virtual void operator=(const Map_af &)=0
Assignment to an affine mapping.
virtual void dalembert(Param &par, Scalar &fJp1, const Scalar &fJ, const Scalar &fJm1, const Scalar &source) const =0
Performs one time-step integration of the d&#39;Alembert scalar equation.
Metric_flat * p_flat_met_spher
Pointer onto the flat metric associated with the spherical coordinates and with components expressed ...
Definition: map.h:728
Map_af(const Mg3d &mgrille, const double *r_limits)
Standard Constructor.
Definition: map_af.C:216
virtual void resize(int, double)
< Not implemented
virtual void comp_z_from_spherical(const Scalar &v_r, const Scalar &v_theta, Scalar &v_z) const =0
Computes the Cartesian z component (with respect to bvect_cart ) of a vector given by its spherical c...
const Valeur & get_ff() const
Returns a (constant) reference to the function .
Definition: map_et.C:1057
virtual void mult_r_zec(Scalar &) const
Multiplication by r (in the compactified external domain only) of a Scalar.
virtual void sauve(FILE *) const
Save in a file.
Definition: map_et.C:802
virtual void poisson(const Cmp &source, Param &par, Cmp &uu, bool verbose=true) const
Computes the solution of a scalar Poisson equation (Cmp version).
Tbl aasx2
Values at the nr collocation points of in the nucleus.
Definition: map.h:2899
virtual void poisson(const Cmp &source, Param &par, Cmp &uu, bool verbose=true) const =0
Computes the solution of a scalar Poisson equation (Cmp version).
Valeur ff
Values of the function at the nt*np angular collocation points in each domain.
Definition: map.h:2937
Coord stdrdp
in the nucleus and in the non-compactified shells; \ in the compactified external domain (CED)...
Definition: map.h:1651
virtual void div_rsint(Scalar &) const =0
Division by of a Scalar.
virtual void homothetie(double lambda)=0
Sets a new radial scale.
virtual void operator=(const Map_af &mpa)
Assignment to an affine mapping.
Definition: map_log.C:233
virtual void poisson2d(const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu, bool verbose=true) const
Computes the solution of a 2-D Poisson equation.
virtual void poisson_regular(const Cmp &source, int k_div, int nzet, double unsgam1, Param &par, Cmp &uu, Cmp &uu_regu, Cmp &uu_div, Tenseur &duu_div, Cmp &source_regu, Cmp &source_div) const
Computes the solution of a scalar Poisson equation.
virtual void reevaluate(const Map *mp_prev, int nzet, Cmp &uu) const
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
virtual void div_cost(Scalar &) const
Division by of a Scalar.
Logarithmic radial mapping.
Definition: map.h:3794
Coord dxdr
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1635
virtual void div_rsint(Scalar &) const
Division by of a Scalar.
void sol_elliptic(Param_elliptic &params, const Scalar &so, Scalar &uu) const
General elliptic solver.
double ori_y
Absolute coordinate y of the origin.
Definition: map.h:705
const double * get_beta() const
Returns the pointer on the array beta.
Definition: map_af.C:611
virtual void poisson_frontiere(const Cmp &, const Valeur &, int, int, Cmp &, double=0., double=0.) const
Solver of the Poisson equation with boundary condition for the affine mapping case.
double ori_z
Absolute coordinate z of the origin.
Definition: map.h:706
virtual Tbl * integrale(const Scalar &) const
Computes the integral over all space of a Scalar.
Definition: map_et_integ.C:82
void sol_elliptic_no_zec(Param_elliptic &params, const Scalar &so, Scalar &uu, double val) const
General elliptic solver.
virtual void mult_sint(Scalar &) const
Multiplication by of a Scalar.
virtual void laplacien(const Scalar &uu, int zec_mult_r, Scalar &lap) const
Computes the Laplacian of a scalar field.
Definition: map_et_lap.C:78
virtual void poisson(const Cmp &, Param &, Cmp &, bool) const
< Not implemented
Tbl alpha
Array (size: mg->nzone ) of the values of in each domain.
Definition: map.h:3800
virtual void resize(int l, double lambda)
Rescales the outer boundary of one domain.
Definition: map_et.C:951
Coord rsx2drdx
in the nucleus and the shells; \ in the outermost compactified domain.
Definition: map.h:2959
Parameter storage.
Definition: param.h:125
virtual void srstdsdp(const Cmp &ci, Cmp &resu) const =0
Computes of a Cmp .
Base class for pure radial mappings.
Definition: map.h:1611
virtual void dsdr(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_et_deriv.C:190
virtual void poisson_regular(const Cmp &source, int k_div, int nzet, double unsgam1, Param &par, Cmp &uu, Cmp &uu_regu, Cmp &uu_div, Tenseur &duu_div, Cmp &source_regu, Cmp &source_div) const =0
Computes the solution of a scalar Poisson equation.
double * beta
Array (size: mg->nzone ) of the values of in each domain.
Definition: map.h:2110
Coord sinp
Definition: map.h:749
Tbl ** aa
Array (size: mg->nzone ) of Tbl which stores the values of in each domain.
Definition: map.h:2883
void set_alpha(double alpha0, int l)
Modifies the value of in domain no. l.
Definition: map_et.C:447
virtual double val_r_jk(int l, double xi, int j, int k) const =0
Returns the value of the radial coordinate r for a given and a given collocation point in in a give...
virtual void adapt(const Cmp &, const Param &, int)
< Not implemented
virtual void homothetie(double lambda)
Sets a new radial scale.
Definition: map_et.C:928
virtual void reevaluate_symy(const Map *mp_prev, int nzet, Cmp &uu) const
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
Coord sstd2rdpdx
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1723
virtual void poisson(const Cmp &source, Param &par, Cmp &uu, bool verbose=true) const
Computes the solution of a scalar Poisson equation (Cmp version).
virtual void sauve(FILE *) const
Save in a file.
Definition: map_radial.C:119
void set_alpha(double alpha0, int l)
Modifies the value of in domain no. l.
Definition: map_af.C:760
virtual void val_lx(double rr, double theta, double pphi, int &l, double &xi) const
Computes the domain index l and the value of corresponding to a point given by its physical coordina...
virtual void poisson_interne(const Cmp &, const Valeur &, Param &, Cmp &) const
< Not implemented
virtual void mult_cost(Scalar &) const =0
Multiplication by of a Scalar.
Coord xsr
in the nucleus; \ 1/R in the non-compactified shells; \ in the compactified outer domain...
Definition: map.h:1624
Metric_flat * p_flat_met_cart
Pointer onto the flat metric associated with the Cartesian coordinates and with components expressed ...
Definition: map.h:733
virtual void poisson_interne(const Cmp &source, const Valeur &limite, Param &par, Cmp &pot) const
Computes the solution of a Poisson equation in the shell .
virtual Param * donne_para_poisson_vect(Param &para, int i) const =0
Function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara .
This class contains the parameters needed to call the general elliptic solver.
virtual void reset_coord()
Resets all the member Coords.
Definition: map.C:282
Coord drdt
in the nucleus and in the non-compactified shells; \ in the compactified external domain (CED)...
Definition: map.h:1643
virtual void dsdr(const Scalar &ci, Scalar &resu) const
Computes of a Scalar.
virtual void poisson_compact(const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const
Resolution of the elliptic equation in the case where the stellar interior is covered by a single do...
virtual void val_lx_jk(double rr, int j, int k, const Param &par, int &l, double &xi) const
Computes the domain index l and the value of corresponding to a point of arbitrary r but collocation...
virtual void lapang(const Scalar &uu, Scalar &lap) const
Computes the angular Laplacian of a scalar field.
Definition: map_af_lap.C:552
virtual void resize(int l, double lambda)=0
Rescales the outer boundary of one domain.
Active physical coordinates and mapping derivatives.
Definition: coord.h:90
Map_log(const Mg3d &mgrille, const Tbl &r_limits, const Itbl &typevar)
Standard Constructor.
Definition: map_log.C:70
void sol_elliptic_no_zec(Param_elliptic &params, const Scalar &so, Scalar &uu, double) const
General elliptic solver.
virtual void val_lx_jk(double rr, int j, int k, const Param &par, int &l, double &xi) const =0
Computes the domain index l and the value of corresponding to a point of arbitrary r but collocation...
Cartesian vectorial bases (triads).
Definition: base_vect.h:201
virtual void dsdradial(const Scalar &uu, Scalar &resu) const
Computes of a Scalar if the description is affine and if it is logarithmic.
Definition: map_et_deriv.C:277
void sol_elliptic(Param_elliptic &params, const Scalar &so, Scalar &uu) const
General elliptic solver.
virtual void primr(const Scalar &, Scalar &, bool) const
< Not implemented
Spherical orthonormal vectorial bases (triads).
Definition: base_vect.h:308
virtual void poisson_regular(const Cmp &, int, int, double, Param &, Cmp &, Cmp &, Cmp &, Tenseur &, Cmp &, Cmp &) const
< Not implemented
Coord rsxdxdr
in the nucleus; \ in the shells; \ in the outermost compactified domain.
Definition: map.h:2952
virtual void reset_coord()
Resets all the member Coords.
Definition: map_et.C:651
virtual void mult_r(Scalar &uu) const
Multiplication by r of a Scalar, the dzpuis of uu is not changed.
virtual void poisson_tau(const Cmp &source, Param &par, Cmp &uu) const
Computes the solution of a scalar Poisson equation using a Tau method (Cmp version).
friend Mtbl * map_log_fait_r(const Map *)
< Not implemented
Definition: map_log_fait.C:60
virtual void mult_r(Scalar &uu) const =0
Multiplication by r of a Scalar , the dzpuis of uu is not changed.
virtual void inc2_dzpuis(Scalar &) const =0
Increases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the compactifi...
virtual void poisson_ylm(const Cmp &, Param &, Cmp &, int, double *) const
< Not implemented
virtual void lapang(const Scalar &, Scalar &) const
< Not implemented
Coord ya
Absolute y coordinate.
Definition: map.h:757
virtual void poisson_angu(const Scalar &source, Param &par, Scalar &uu, double lambda=0) const
Computes the solution of the generalized angular Poisson equation.
virtual void poisson_falloff(const Cmp &, Param &, Cmp &, int) const
< Not implemented
const Mg3d * mg
Pointer on the multi-grid Mgd3 on which this is defined.
Definition: map.h:702
Multi-domain grid.
Definition: grilles.h:276
virtual void poisson_frontiere(const Cmp &source, const Valeur &limite, int raccord, int num_front, Cmp &pot, double=0., double=0.) const =0
Computes the solution of a Poisson equation from the domain num_front+1 .
virtual void adapt(const Cmp &ent, const Param &par, int nbr_filtre=0)
Adaptation of the mapping to a given scalar field.
Definition: map_et_adapt.C:111
virtual void dalembert(Param &par, Scalar &fJp1, const Scalar &fJ, const Scalar &fJm1, const Scalar &source) const
Not yet implemented.
double ori_x
Absolute coordinate x of the origin.
Definition: map.h:704
Affine radial mapping.
Definition: map.h:2102
Tbl zaasx2
Values at the nr collocation points of in the outermost compactified domain.
Definition: map.h:2909
virtual void dalembert(Param &, Scalar &, const Scalar &, const Scalar &, const Scalar &) const
< Not implemented
virtual void reevaluate_symy(const Map *mp_prev, int nzet, Cmp &uu) const =0
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
void sol_elliptic_fixe_der_zero(double val, Param_elliptic &params, const Scalar &so, Scalar &uu) const
General elliptic solver fixing the derivative at the origin and relaxing the continuity of the first ...
virtual Tbl * integrale(const Scalar &) const =0
Computes the integral over all space of a Scalar .
const Base_vect_cart & get_bvect_cart() const
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:817
virtual void poisson_angu(const Scalar &source, Param &par, Scalar &uu, double lambda=0) const
Computes the solution of the generalized angular Poisson equation.
virtual void div_r_zec(Scalar &) const
Division by r (in the compactified external domain only) of a Scalar.
Map_radial(const Mg3d &mgrid)
Constructor from a grid (protected to make Map_radial an abstract class)
Definition: map_radial.C:92
Valeur gg
Values of the function at the nt*np angular collocation points in each domain.
Definition: map.h:2944
Coord y
y coordinate centered on the grid
Definition: map.h:753
virtual void dsdradial(const Scalar &uu, Scalar &resu) const
Computes of a Scalar if the description is affine and if it is logarithmic.
virtual void inc2_dzpuis(Scalar &) const
Increases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the compactifi...
virtual void comp_x_from_spherical(const Scalar &v_r, const Scalar &v_theta, const Scalar &v_phi, Scalar &v_x) const
Computes the Cartesian x component (with respect to bvect_cart) of a vector given by its spherical co...
virtual void operator=(const Map_af &)=0
Assignment to an affine mapping.
virtual void homothetie(double)
Sets a new radial scale.
virtual void comp_r_from_cartesian(const Scalar &v_x, const Scalar &v_y, const Scalar &v_z, Scalar &v_r) const =0
Computes the Spherical r component (with respect to bvect_spher ) of a vector given by its cartesian ...
virtual void div_sint(Scalar &) const =0
Division by of a Scalar.
Coord za
Absolute z coordinate.
Definition: map.h:758
Coord lapr_tp
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1706
virtual const Map_af & mp_angu(int) const
Returns the "angular" mapping for the outside of domain l_zone.
Definition: map_et.C:1068
virtual ~Map_et()
Destructor.
Definition: map_et.C:509
double get_beta(int l) const
Returns in the domain l.
Definition: map.h:3850
virtual void reevaluate(const Map *mp_prev, int nzet, Cmp &uu) const =0
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
Coord cosp
Definition: map.h:750
double * beta
Array (size: mg->nzone ) of the values of in each domain.
Definition: map.h:2878
Coefficients storage for the multi-domain spectral method.
Definition: mtbl_cf.h:196
void resize_extr(double lambda)
Rescales the outer boundary of the outermost domain in the case of non-compactified external domain...
virtual void poisson_angu(const Scalar &, Param &, Scalar &, double=0) const
< Not implemented
virtual void dsdxi(const Scalar &ci, Scalar &resu) const
Computes of a Scalar.
Definition: map_log_deriv.C:57
virtual void operator=(const Map_et &mp)
Assignment to another Map_et.
Definition: map_et.C:536
virtual void poisson_interne(const Cmp &source, const Valeur &limite, Param &par, Cmp &pot) const
Computes the solution of a Poisson equation in the shell, imposing a boundary condition at the surfac...
Coord x
x coordinate centered on the grid
Definition: map.h:752
Base_vect_spher bvect_spher
Orthonormal vectorial basis associated with the coordinates of the mapping.
Definition: map.h:715
Scalar interpolate_from_map_star(const Scalar &f_s) const
Interpolates from a Scalar defined on a Map_star and returns the new Scalar defined on *this...
virtual void div_sint(Scalar &) const
Division by of a Scalar.
Map_et(const Mg3d &mgrille, const double *r_limits)
Standard Constructor.
Definition: map_et.C:155
double get_ori_z() const
Returns the z coordinate of the origin.
Definition: map.h:798
virtual void val_lx_jk(double rr, int j, int k, const Param &par, int &l, double &xi) const
Computes the domain index l and the value of corresponding to a point of arbitrary r but collocation...
virtual void srdsdt(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_et_deriv.C:340
virtual void reset_coord()
Resets all the member Coords.
Definition: map_radial.C:129
virtual void mult_r_zec(Scalar &) const =0
Multiplication by r (in the compactified external domain only) of a Scalar.
void sol_elliptic_pseudo_1d(Param_elliptic &, const Scalar &, Scalar &) const
General elliptic solver in a pseudo 1d case.
void sol_elliptic_boundary(Param_elliptic &params, const Scalar &so, Scalar &uu, const Mtbl_cf &bound, double fact_dir, double fact_neu) const
General elliptic solver including inner boundary conditions.
void set_ff(const Valeur &)
Assigns a given value to the function .
Definition: map_et.C:585
Basic array class.
Definition: tbl.h:164
virtual void mult_rsint(Scalar &) const
Multiplication by of a Scalar.
void set_coord()
Assignement of the building functions to the member Coords.
Definition: map_et.C:607
void homothetie_interne(double lambda)
Sets a new radial scale at the bondary between the nucleus and the first shell.
Definition: map_af.C:744
virtual void laplacien(const Scalar &uu, int zec_mult_r, Scalar &lap) const
Computes the Laplacian of a scalar field.
Definition: map_af_lap.C:182
virtual bool operator==(const Map &) const
Comparison operator (egality)
Definition: map_et.C:1007
virtual void poisson_frontiere(const Cmp &, const Valeur &, int, int, Cmp &, double=0., double=0.) const
< Not implemented
virtual void dsdradial(const Scalar &uu, Scalar &resu) const =0
Computes of a Scalar if the description is affine and if it is logarithmic.
virtual void poisson_tau(const Cmp &, Param &, Cmp &) const
< Not implemented
virtual void srdsdt(const Cmp &, Cmp &) const
< Not implemented
Tbl bbsx2
Values at the nr collocation points of in the nucleus.
Definition: map.h:2930
virtual void poisson2d(const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu, bool verbose=true) const
Computes the solution of a 2-D Poisson equation.
Map(const Mg3d &)
Constructor from a multi-domain 3D grid.
Definition: map.C:145
void set_beta(double beta0, int l)
Modifies the value of in domain no. l.
Definition: map_et.C:458
void sol_elliptic_boundary(Param_elliptic &params, const Scalar &so, Scalar &uu, const Mtbl_cf &bound, double fact_dir, double fact_neu) const
General elliptic solver including inner boundary conditions.
virtual ostream & operator>>(ostream &) const
Operator >>
Definition: map_et.C:821
virtual void dec_dzpuis(Scalar &) const =0
Decreases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactifi...
void convert_absolute(double xx, double yy, double zz, double &rr, double &theta, double &pphi) const
Determines the coordinates corresponding to given absolute Cartesian coordinates (X...
Definition: map.C:308
virtual void div_cost(Scalar &) const =0
Division by of a Scalar.
Coord srdrdt
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1659
virtual void poisson_frontiere_double(const Cmp &, const Valeur &, const Valeur &, int, Cmp &) const
< Not implemented
virtual ~Map_log()
Destructor.
Definition: map_log.C:163
virtual void adapt(const Cmp &ent, const Param &par, int nbr=0)
Adaptation of the mapping to a given scalar field.
Definition: map_af.C:803
virtual ~Map()
Destructor.
Definition: map.C:219
virtual void sauve(FILE *) const
Save in a file.
Definition: map_log.C:166
virtual void resize(int l, double lambda)
Rescales the outer boundary of one domain.
Definition: map_af.C:690
Coord z
z coordinate centered on the grid
Definition: map.h:754
virtual void div_r(Scalar &) const =0
Division by r of a Scalar.
double rot_phi
Angle between the x –axis and X –axis.
Definition: map.h:707
virtual void lapang(const Scalar &uu, Scalar &lap) const
Computes the angular Laplacian of a scalar field.
Definition: map_et_lap.C:286
void sol_elliptic_sin_zec(Param_elliptic &params, const Scalar &so, Scalar &uu, double *coefs, double *) const
General elliptic solver.
Tbl ** dbb
Array (size: mg->nzone ) of Tbl which stores the values of in each domain.
Definition: map.h:2919
virtual void comp_p_from_cartesian(const Scalar &v_x, const Scalar &v_y, Scalar &v_p) const
Computes the Spherical component (with respect to bvect_spher ) of a vector given by its cartesian c...
virtual void poisson_angu(const Scalar &source, Param &par, Scalar &uu, double lambda=0) const =0
Computes the solution of the generalized angular Poisson equation.
const Metric_flat & flat_met_spher() const
Returns the flat metric associated with the spherical coordinates and with components expressed in th...
Definition: map.C:327
virtual void dsdr(const Cmp &ci, Cmp &resu) const
Computes of a Cmp.
Definition: map_af_deriv.C:282
Tbl aasx
Values at the nr collocation points of in the nucleus.
Definition: map.h:2896
virtual void dsdt(const Scalar &uu, Scalar &resu) const
Computes of a Scalar.
Definition: map_af_deriv.C:800
virtual double val_r_jk(int l, double xi, int j, int k) const
< Comparison operator
Tensor handling *** DEPRECATED : use class Tensor instead ***.
Definition: tenseur.h:307
virtual void comp_y_from_spherical(const Scalar &v_r, const Scalar &v_theta, const Scalar &v_phi, Scalar &v_y) const =0
Computes the Cartesian y component (with respect to bvect_cart ) of a vector given by its spherical c...
Base_vect_cart bvect_cart
Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
Definition: map.h:723
virtual ~Map_radial()
Destructor.
Definition: map_radial.C:113
virtual void laplacien(const Scalar &uu, int zec_mult_r, Scalar &lap) const =0
Computes the Laplacian of a scalar field.
Tbl ** ddaa
Array (size: mg->nzone ) of Tbl which stores the values of in each domain.
Definition: map.h:2893
virtual Param * donne_para_poisson_vect(Param &par, int i) const
Internal function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara ...
Coord r
r coordinate centered on the grid
Definition: map.h:744
virtual void poisson2d(const Cmp &, const Cmp &, Param &, Cmp &, bool=true) const
< Not implemented
Coord d2rdtdx
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1715
virtual void lapang(const Scalar &uu, Scalar &lap) const =0
Computes the angular Laplacian of a scalar field.
friend ostream & operator<<(ostream &, const Map &)
Operator <<.
Definition: map.C:245
virtual void val_lx(double rr, double theta, double pphi, int &l, double &xi) const =0
Computes the domain index l and the value of corresponding to a point given by its physical coordina...
virtual void poisson_tau(const Cmp &source, Param &par, Cmp &uu) const
Computes the solution of a scalar Poisson equation with a Tau method (Cmp version).
virtual void stdsdp(const Scalar &uu, Scalar &resu) const =0
Computes of a Scalar .
virtual double val_r_jk(int l, double xi, int j, int k) const
Returns the value of the radial coordinate r for a given and a given collocation point in in a give...
virtual bool operator==(const Map &) const =0
Comparison operator (egality)
virtual const Map_af & mp_angu(int) const
Returns the "angular" mapping for the outside of domain l_zone.
Coord cost
Definition: map.h:748