LORENE
map_et_poisson2d.C
1 /*
2  * Method of the class Map_et for the resolution of the 2-D Poisson
3  * equation.
4  *
5  * (see file map.h for documentation).
6  */
7 
8 /*
9  * Copyright (c) 2000-2001 Eric Gourgoulhon
10  *
11  * This file is part of LORENE.
12  *
13  * LORENE is free software; you can redistribute it and/or modify
14  * it under the terms of the GNU General Public License as published by
15  * the Free Software Foundation; either version 2 of the License, or
16  * (at your option) any later version.
17  *
18  * LORENE is distributed in the hope that it will be useful,
19  * but WITHOUT ANY WARRANTY; without even the implied warranty of
20  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
21  * GNU General Public License for more details.
22  *
23  * You should have received a copy of the GNU General Public License
24  * along with LORENE; if not, write to the Free Software
25  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
26  *
27  */
28 
29 
30 
31 
32 /*
33  * $Id: map_et_poisson2d.C,v 1.6 2026/03/04 15:47:04 j_novak Exp $
34  * $Log: map_et_poisson2d.C,v $
35  * Revision 1.6 2026/03/04 15:47:04 j_novak
36  * Added 'verbose' flag (true by default) to limit the output of poisson solvers.
37  *
38  * Revision 1.5 2016/12/05 16:17:58 j_novak
39  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
40  *
41  * Revision 1.4 2014/10/13 08:53:05 j_novak
42  * Lorene classes and functions now belong to the namespace Lorene.
43  *
44  * Revision 1.3 2014/10/06 15:13:13 j_novak
45  * Modified #include directives to use c++ syntax.
46  *
47  * Revision 1.2 2002/02/07 14:55:58 e_gourgoulhon
48  * Corrected a bug when the source is known only in the coefficient
49  * space.
50  *
51  * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
52  * LORENE
53  *
54  * Revision 2.4 2000/11/07 14:21:03 eric
55  * Correction d'une erreur dans le cas T_SIN_I (calcul de R(u)).
56  *
57  * Revision 2.3 2000/10/26 15:58:00 eric
58  * Correction cas T_COS_P : l'import de saff_q se fait par copie du Tbl.
59  *
60  * Revision 2.2 2000/10/12 15:37:43 eric
61  * Traitement des bases spectrales dans le cas T_COS_P.
62  *
63  * Revision 2.1 2000/10/11 15:15:43 eric
64  * 1ere version operationnelle.
65  *
66  * Revision 2.0 2000/10/09 13:47:17 eric
67  * *** empty log message ***
68  *
69  *
70  * $Header: /cvsroot/Lorene/C++/Source/Map/map_et_poisson2d.C,v 1.6 2026/03/04 15:47:04 j_novak Exp $
71  *
72  */
73 
74 // Headers C
75 #include <cmath>
76 
77 // Headers Lorene:
78 #include "map.h"
79 #include "cmp.h"
80 #include "param.h"
81 
82 //*****************************************************************************
83 
84 namespace Lorene {
85 
86 void Map_et::poisson2d(const Cmp& source_mat, const Cmp& source_quad,
87  Param& par, Cmp& uu, bool verbose) const {
88 
89  assert(source_mat.get_etat() != ETATNONDEF) ;
90  assert(source_quad.get_etat() != ETATNONDEF) ;
91  assert(source_mat.get_mp()->get_mg() == mg) ;
92  assert(source_quad.get_mp()->get_mg() == mg) ;
93  assert(uu.get_mp()->get_mg() == mg) ;
94 
95  assert( source_quad.check_dzpuis(4) ) ;
96 
97  double& lambda = par.get_double_mod(0) ;
98  int nz = mg->get_nzone() ;
99  int nzm1 = nz-1 ;
100 
101  // Special case of a vanishing source
102  // ----------------------------------
103 
104  if ( (source_mat.get_etat() == ETATZERO)
105  && (source_quad.get_etat() == ETATZERO) ) {
106 
107  uu = 0 ;
108  lambda = 1 ;
109  return ;
110  }
111 
112  int base_t = ((source_mat.va).base).get_base_t(0) ;
113 
114  switch (base_t) {
115 
116  //==================================================================
117  // case T_COS_P
118  //==================================================================
119 
120  case T_COS_P : {
121 
122  // Construction of a Map_af which coincides with *this on the equator
123  // ------------------------------------------------------------------
124 
125  double theta0 = M_PI / 2 ; // Equator
126  double phi0 = 0 ;
127 
128  Map_af mpaff(*this) ;
129 
130  for (int l=0 ; l<nz ; l++) {
131  double rmax = val_r(l, double(1), theta0, phi0) ;
132  switch ( mg->get_type_r(l) ) {
133  case RARE: {
134  double rmin = val_r(l, double(0), theta0, phi0) ;
135  mpaff.set_alpha(rmax - rmin, l) ;
136  mpaff.set_beta(rmin, l) ;
137  break ;
138  }
139 
140  case FIN: {
141  double rmin = val_r(l, double(-1), theta0, phi0) ;
142  mpaff.set_alpha( double(.5) * (rmax - rmin), l ) ;
143  mpaff.set_beta( double(.5) * (rmax + rmin), l) ;
144  break ;
145  }
146 
147  case UNSURR: {
148  double rmin = val_r(l, double(-1), theta0, phi0) ;
149  double umax = double(1) / rmin ;
150  double umin = double(1) / rmax ;
151  mpaff.set_alpha( double(.5) * (umin - umax), l) ;
152  mpaff.set_beta( double(.5) * (umin + umax), l) ;
153  break ;
154  }
155 
156  default: {
157  cerr << "Map_et::poisson2d: unknown type_r ! " << endl ;
158  abort () ;
159  break ;
160  }
161 
162  }
163  }
164 
165  // Importation of source_mat and source_quad of the affine mapping
166  // ---------------------------------------------------------------
167  Cmp saff_m(mpaff) ;
168  saff_m.import( nzm1, source_mat ) ;
169  (saff_m.va).set_base( (source_mat.va).base ) ;
170 
171  Cmp saff_q(mpaff) ;
172 
173  // In order to use Cmp::import with dzpuis != 0 :
174  Cmp set_q = source_quad ;
175  set_q.set_dzpuis(0) ; // dzpuis artificially set to 0
176 
177  saff_q.import( nzm1, set_q ) ;
178  (saff_q.va).set_base( (set_q.va).base ) ;
179 
180  // Copy in the external domain :
181  if ( (set_q.va).get_etat() == ETATQCQ) {
182  (set_q.va).coef_i() ; // the values in configuration space are required
183  assert( (set_q.va).c->get_etat() == ETATQCQ ) ;
184  assert( (saff_q.va).c->get_etat() == ETATQCQ ) ;
185  *( (saff_q.va).c->t[nzm1] ) = *( (set_q.va).c->t[nzm1] ) ;
186  }
187 
188  // the true dzpuis is restored :
189  saff_q.set_dzpuis( source_quad.get_dzpuis() ) ;
190 
191  // Resolution of the 2-D Poisson equation on the spherical domains
192  // ---------------------------------------------------------------
193 
194  Cmp uaff(mpaff) ;
195 
196  mpaff.poisson2d(saff_m, saff_q, par, uaff) ;
197 
198  // Importation of the solution on the Map_et mapping *this
199  // -------------------------------------------------------
200 
201  uu.import(uaff) ;
202 
203  uu.va.set_base( uaff.va.base ) ; // same spectral bases
204 
205  break ;
206  }
207 
208  //==================================================================
209  // case T_SIN_I
210  //==================================================================
211 
212  case T_SIN_I : {
213 
214  //-------------------------------
215  // Computation of the prefactor a ---> Cmp apre
216  //-------------------------------
217 
218  Mtbl unjj = 1 + srdrdt*srdrdt ;
219 
220  Mtbl apre1(*mg) ;
221  apre1.set_etat_qcq() ;
222  for (int l=0; l<nz; l++) {
223  *(apre1.t[l]) = alpha[l]*alpha[l] ;
224  }
225 
226  apre1 = apre1 * dxdr * dxdr * unjj ;
227 
228  Cmp apre(*this) ;
229  apre = apre1 ;
230 
231  Tbl amax0 = max(apre1) ; // maximum values in each domain
232 
233  // The maximum values of a in each domain are put in a Mtbl
234  Mtbl amax1(*mg) ;
235  amax1.set_etat_qcq() ;
236  for (int l=0; l<nz; l++) {
237  *(amax1.t[l]) = amax0(l) ;
238  }
239 
240  Cmp amax(*this) ;
241  amax = amax1 ;
242 
243 
244  //-------------------
245  // Initializations
246  //-------------------
247 
248  int nitermax = par.get_int() ;
249  int& niter = par.get_int_mod() ;
250  double lambda_relax = par.get_double() ;
251  double unmlambda_relax = 1. - lambda_relax ;
252  double precis = par.get_double(1) ;
253 
254  Cmp& ssj = par.get_cmp_mod() ;
255 
256  Cmp ssjm1 = ssj ;
257  Cmp ssjm2 = ssjm1 ;
258 
259  Cmp ssj_q(*this) ;
260  ssj_q = 0 ;
261 
262  Valeur& vuu = uu.va ;
263 
264  Valeur vuujm1(*mg) ;
265  if (uu.get_etat() == ETATZERO) {
266  vuujm1 = 1 ; // to take relative differences
267  vuujm1.set_base( vuu.base ) ;
268  }
269  else{
270  vuujm1 = vuu ;
271  }
272 
273  // Affine mapping for the Laplacian-tilde
274 
275  Map_af mpaff(*this) ;
276 
277  if (verbose)
278  cout << "Map_et::poisson2d : relat. diff. u^J <-> u^{J-1} : " << endl ;
279 
280 //==========================================================================
281 //==========================================================================
282 // Start of iteration
283 //==========================================================================
284 //==========================================================================
285 
286  Tbl tdiff(nz) ;
287  double diff ;
288  niter = 0 ;
289 
290  do {
291 
292  //====================================================================
293  // Computation of R(u) (the result is put in uu)
294  //====================================================================
295 
296 
297  //-----------------------
298  // First operations on uu
299  //-----------------------
300 
301  Valeur duudx = (uu.va).dsdx() ; // d/dx
302 
303  const Valeur& d2uudtdx = duudx.dsdt() ; // d^2/dxdtheta
304 
305  //-------------------
306  // 1/x d^2uu/dtheta^2
307  //-------------------
308 
309  Valeur sxlapang = uu.va ;
310 
311  sxlapang = sxlapang.d2sdt2() ;
312 
313  sxlapang = sxlapang.sx() ; // d^2(uu)/dth^2 /x in the nucleus
314  // d^2(uu)/dth^2 in the shells
315  // d^2(uu)/dth^2 /(x-1) in the ZEC
316 
317  //---------------------------------------------------------------
318  // Computation of
319  // [ (dR/dx)^{-1} ( A - 1 ) duu/dx + 1/R (B - 1) d^2uu/dth^2 ] / x
320  //
321  // with A = 1/(dRdx) R/(x+beta_l/alpha_l) unjj
322  // B = [1/(dRdx) R/(x+beta_l/alpha_l)]^2 unjj
323  //
324  // The result is put in uu (via vuu)
325  //---------------------------------------------------------------
326 
327  // Intermediate quantity jac which value is
328  // (dR/dx)^{-1} in the nucleus and the shells
329  // +(dU/dx)^{-1} in the ZEC
330 
331  Mtbl jac = dxdr ;
332  if (mg->get_type_r(nzm1) == UNSURR) {
333  jac.annule(nzm1, nzm1) ;
334  Mtbl jac_ext = dxdr ;
335  jac_ext.annule(0, nzm1-1) ;
336  jac_ext = - jac_ext ;
337  jac = jac + jac_ext ;
338  }
339 
340  uu.set_etat_qcq() ;
341 
342  Base_val sauve_base = duudx.base ;
343 
344  vuu = jac * ( rsxdxdr * unjj - 1.) * duudx
345  + ( rsxdxdr*rsxdxdr * unjj - 1.) * xsr * sxlapang ;
346 
347  vuu.set_base(sauve_base) ;
348 
349  vuu = vuu.sx() ;
350 
351  //---------------------------------------
352  // Computation of R(u)
353  //
354  // The result is put in uu (via vuu)
355  //---------------------------------------
356 
357 
358  sauve_base = vuu.base ;
359 
360  vuu = xsr * vuu
361  + 2. * dxdr * sr2drdt * d2uudtdx ;
362 
363  vuu += dxdr * ( sr2d2rdt2 + dxdr * (
364  dxdr* unjj * d2rdx2
365  - 2. * sr2drdt * d2rdtdx )
366  ) * duudx ;
367 
368  vuu.set_base(sauve_base) ;
369 
370  // Since the assignment is performed on vuu (uu.va), the treatment
371  // of uu.dzpuis must be performed by hand:
372 
373  uu.set_dzpuis(4) ;
374 
375  //====================================================================
376  // Computation of the effective source s^J of the ``affine''
377  // Poisson equation
378  //====================================================================
379 
380  ssj = lambda_relax * ssjm1 + unmlambda_relax * ssjm2 ;
381 
382  ssj = ( source_mat + source_quad + uu + (amax - apre) * ssj ) / amax ;
383 
384  (ssj.va).set_base((source_mat.va).base) ;
385 
386  //====================================================================
387  // Resolution of the ``affine'' Poisson equation
388  //====================================================================
389 
390  assert( uu.check_dzpuis( ssj.get_dzpuis() ) ) ;
391 
392  mpaff.poisson2d(ssj, ssj_q, par, uu) ;
393 
394  tdiff = diffrel(vuu, vuujm1) ;
395 
396  diff = max(tdiff) ;
397 
398  if (verbose) {
399  cout << " step " << niter << " : " ;
400  for (int l=0; l<nz; l++) {
401  cout << tdiff(l) << " " ;
402  }
403  cout << endl ;
404  }
405 
406  //=================================
407  // Updates for the next iteration
408  //=================================
409 
410  ssjm2 = ssjm1 ;
411  ssjm1 = ssj ;
412  vuujm1 = vuu ;
413 
414  niter++ ;
415 
416  } // End of iteration
417  while ( (diff > precis) && (niter < nitermax) ) ;
418 
419 //==========================================================================
420 //==========================================================================
421 // End of iteration
422 //==========================================================================
423 //==========================================================================
424 
425 
426  break ;
427  }
428 
429  default : {
430  cerr << "Map_et::poisson2d : unkown theta basis !" << endl ;
431  cerr << " basis : " << hex << base_t << endl ;
432  abort() ;
433  break ;
434  }
435  } // End of switch on base_t
436 }
437 
438 
439 }
const Map * get_mp() const
Returns the mapping.
Definition: cmp.h:904
Coord d2rdx2
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1694
Coord sr2d2rdt2
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1732
const Valeur & dsdt() const
Returns of *this.
Definition: valeur_dsdt.C:115
double & get_double_mod(int position=0) const
Returns the reference of a stored modifiable double .
Definition: param.C:501
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:449
Multi-domain array.
Definition: mtbl.h:118
double * alpha
Array (size: mg->nzone ) of the values of in each domain.
Definition: map.h:2876
Lorene prototypes.
Definition: app_hor.h:67
const Mg3d * get_mg() const
Gives the Mg3d on which the mapping is defined.
Definition: map.h:791
int get_etat() const
Returns the logical state.
Definition: cmp.h:902
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
Coord sr2drdt
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1675
void annule(int l_min, int l_max)
Sets the Mtbl to zero in some domains.
Definition: mtbl.C:332
const Valeur & sx() const
Returns (r -sampling = RARE ) \ Id (r sampling = FIN ) \ (r -sampling = UNSURR ) ...
Definition: valeur_sx.C:113
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
void set_base(const Base_val &)
Sets the bases for spectral expansions (member base )
Definition: valeur.C:813
Tbl diffrel(const Cmp &a, const Cmp &b)
Relative difference between two Cmp (norme version).
Definition: cmp_math.C:507
Cmp & get_cmp_mod(int position=0) const
Returns the reference of a modifiable Cmp stored in the list.
Definition: param.C:1052
void set_beta(double beta0, int l)
Modifies the value of in domain no. l.
Definition: map_af.C:771
const int & get_int(int position=0) const
Returns the reference of a int stored in the list.
Definition: param.C:295
#define T_COS_P
dev. cos seulement, harmoniques paires
Definition: type_parite.h:200
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.
Coord dxdr
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1635
Base_val base
Bases on which the spectral expansion is performed.
Definition: valeur.h:315
Parameter storage.
Definition: param.h:125
int get_nzone() const
Returns the number of domains.
Definition: grilles.h:467
int & get_int_mod(int position=0) const
Returns the reference of a modifiable int stored in the list.
Definition: param.C:433
Tbl max(const Cmp &)
Maximum values of a Cmp in each domain.
Definition: cmp_math.C:438
void set_alpha(double alpha0, int l)
Modifies the value of in domain no. l.
Definition: map_af.C:760
Coord xsr
in the nucleus; \ 1/R in the non-compactified shells; \ in the compactified outer domain...
Definition: map.h:1624
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: cmp.C:307
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: mtbl.C:302
Coord rsxdxdr
in the nucleus; \ in the shells; \ in the outermost compactified domain.
Definition: map.h:2952
const Mg3d * mg
Pointer on the multi-grid Mgd3 on which this is defined.
Definition: map.h:702
Bases of the spectral expansions.
Definition: base_val.h:325
Affine radial mapping.
Definition: map.h:2102
int get_dzpuis() const
Returns dzpuis.
Definition: cmp.h:906
bool check_dzpuis(int dzi) const
Returns false if the last domain is compactified and *this is not zero in this domain and dzpuis is n...
Definition: cmp.C:718
void set_dzpuis(int)
Set a value to dzpuis.
Definition: cmp.C:657
const double & get_double(int position=0) const
Returns the reference of a double stored in the list.
Definition: param.C:364
Basic array class.
Definition: tbl.h:164
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.
int get_type_r(int l) const
Returns the type of sampling in the radial direction in domain no.
Definition: grilles.h:493
void import(const Cmp &ci)
Assignment to another Cmp defined on a different mapping.
Definition: cmp_import.C:76
Tbl ** t
Array (size nzone ) of pointers on the Tbl &#39;s.
Definition: mtbl.h:132
Coord srdrdt
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1659
#define T_SIN_I
dev. sin seulement, harmoniques impaires
Definition: type_parite.h:206
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:467
const Valeur & d2sdt2() const
Returns of *this.
Coord d2rdtdx
in the nucleus and in the non-compactified shells; \ in the compactified outer domain.
Definition: map.h:1715