LORENE
map_af_poisson.C
1 /*
2  * Method of the class Map_af for the resolution of the scalar Poisson
3  * equation
4  */
5 
6 /*
7  * Copyright (c) 1999-2001 Eric Gourgoulhon
8  * Copyright (c) 1999-2001 Philippe Grandclement
9  *
10  * This file is part of LORENE.
11  *
12  * LORENE is free software; you can redistribute it and/or modify
13  * it under the terms of the GNU General Public License as published by
14  * the Free Software Foundation; either version 2 of the License, or
15  * (at your option) any later version.
16  *
17  * LORENE is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20  * GNU General Public License for more details.
21  *
22  * You should have received a copy of the GNU General Public License
23  * along with LORENE; if not, write to the Free Software
24  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
25  *
26  */
27 
28 
29 
30 
31 /*
32  * $Id: map_af_poisson.C,v 1.10 2026/03/04 15:47:04 j_novak Exp $
33  * $Log: map_af_poisson.C,v $
34  * Revision 1.10 2026/03/04 15:47:04 j_novak
35  * Added 'verbose' flag (true by default) to limit the output of poisson solvers.
36  *
37  * Revision 1.9 2025/03/07 08:22:31 j_novak
38  * Removed unneeded display.
39  *
40  * Revision 1.8 2025/03/04 13:16:50 j_novak
41  * New complete versions of Map_af::poisson() and Map_et::poisson() for Scalar, not using Cmp.
42  *
43  * Revision 1.7 2016/12/05 16:17:57 j_novak
44  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
45  *
46  * Revision 1.6 2014/10/13 08:53:02 j_novak
47  * Lorene classes and functions now belong to the namespace Lorene.
48  *
49  * Revision 1.5 2005/08/25 12:14:09 p_grandclement
50  * Addition of a new method to solve the scalar Poisson equation, based on a multi-domain Tau-method
51  *
52  * Revision 1.4 2004/05/06 15:25:39 e_gourgoulhon
53  * The case dzpuis=5 with null value in CED is well treated now.
54  *
55  * Revision 1.3 2004/02/20 10:55:23 j_novak
56  * The versions dzpuis 5 -> 3 has been improved and polished. Should be
57  * operational now...
58  *
59  * Revision 1.2 2004/02/06 10:53:52 j_novak
60  * New dzpuis = 5 -> dzpuis = 3 case (not ready yet).
61  *
62  * Revision 1.1.1.1 2001/11/20 15:19:27 e_gourgoulhon
63  * LORENE
64  *
65  * Revision 1.9 2000/05/22 13:46:48 phil
66  * ajout du cas dzpuis = 3
67  *
68  * Revision 1.8 2000/02/09 14:44:24 eric
69  * Traitement de dzpuis ameliore.
70  *
71  * Revision 1.7 1999/12/22 16:37:10 eric
72  * Ajout de pot.set_dzpuis(0) a la fin.
73  *
74  * Revision 1.6 1999/12/22 15:11:03 eric
75  * Remplacement du test source.get_mp() == this par
76  * source.get_mp()->get_mg() == mg
77  * (idem pour pot),
78  * afin de permettre l'appel par Map_et::poisson.
79  *
80  * Revision 1.5 1999/12/21 13:02:37 eric
81  * Changement de prototype de la routine poisson : la solution est
82  * desormais passee en argument (et non plus en valeur de retour)
83  * pour s'adapter au prototype general de la fonction virtuelle
84  * Map::poisson.
85  *
86  * Revision 1.4 1999/12/21 10:06:29 eric
87  * Ajout de l'argument (muet) Param&.
88  *
89  * Revision 1.3 1999/12/07 16:48:50 phil
90  * On fait ylm_i avant de quitter
91  *
92  * Revision 1.2 1999/12/02 16:12:22 eric
93  * *** empty log message ***
94  *
95  * Revision 1.1 1999/12/02 14:30:07 eric
96  * Initial revision
97  *
98  *
99  * $Header: /cvsroot/Lorene/C++/Source/Map/map_af_poisson.C,v 1.10 2026/03/04 15:47:04 j_novak Exp $
100  *
101  */
102 
103 // Header Lorene:
104 #include "map.h"
105 #include "cmp.h"
106 #include "tensor.h"
107 
108 namespace Lorene {
109 Mtbl_cf sol_poisson(const Map_af&, const Mtbl_cf&, int, bool match = true) ;
110 Mtbl_cf sol_poisson_tau(const Map_af&, const Mtbl_cf&, int) ;
111 //*****************************************************************************
112 
113  void Map_af::poisson(const Scalar& source, Param& , Scalar& pot, bool) const {
114 
115  assert(source.get_etat() != ETATNONDEF) ;
116  assert(source.get_mp().get_mg() == mg) ;
117  assert(pot.get_mp().get_mg() == mg) ;
118 
119  assert( source.check_dzpuis(2) || source.check_dzpuis(4)
120  || source.check_dzpuis(3) || source.check_dzpuis(5) ) ;
121 
122  bool match = true ;
123 
124  int dzpuis ;
125 
126  if ( (source.dz_nonzero()) || (source.get_dzpuis() > 3)) { //##awkward??
127  dzpuis = source.get_dzpuis() ;
128  }
129  else{
130  dzpuis = 4 ;
131  }
132 
133  match = !(dzpuis == 5) ;
134 
135  // Spherical harmonic expansion of the source
136  // ------------------------------------------
137 
138  const Valeur& sourva = source.get_spectral_va() ;
139 
140  if (sourva.get_etat() == ETATZERO) {
141  pot.set_etat_zero() ;
142  return ;
143  }
144 
145  // Spectral coefficients of the source
146  assert(sourva.get_etat() == ETATQCQ) ;
147 
148  Valeur rho(sourva.get_mg()) ;
149  sourva.coef() ;
150  rho = *(sourva.c_cf) ; // copy of the coefficients of the source
151 
152  rho.ylm() ; // spherical harmonic transforms
153 
154  // Call to the Mtbl_cf version
155  // ---------------------------
156  Mtbl_cf resu = sol_poisson(*this, *(rho.c_cf), dzpuis, match) ;
157 
158  // Final result returned as a Scalar
159  // ----------------------------------
160 
161  pot.set_etat_zero() ; // to call Scalar::del_deriv().
162 
163  pot.set_etat_qcq() ;
164 
165  pot.set_spectral_va() = resu ;
166  pot.set_spectral_va().ylm_i() ; // On repasse en base standard.
167 
168  (dzpuis == 5) ? pot.set_dzpuis(3) : pot.set_dzpuis(0) ;
169 
170 }
171 
172 
173  //----------------------
174  // Tau version method
175  //---------------------
176 
177 
178 void Map_af::poisson_tau(const Scalar& source, Param& , Scalar& pot) const {
179 
180  assert(source.get_etat() != ETATNONDEF) ;
181  assert(source.get_mp().get_mg() == mg) ;
182  assert(pot.get_mp().get_mg() == mg) ;
183 
184  assert( source.check_dzpuis(2) || source.check_dzpuis(4)
185  || source.check_dzpuis(3)) ;
186 
187  int dzpuis ;
188 
189  if ( (source.dz_nonzero()) || (source.get_dzpuis() > 3)) { //##awkward??
190  dzpuis = source.get_dzpuis() ;
191  }
192  else{
193  dzpuis = 4 ;
194  }
195 
196  // Spherical harmonic expansion of the source
197  // ------------------------------------------
198 
199  const Valeur& sourva = source.get_spectral_va() ;
200 
201  if (sourva.get_etat() == ETATZERO) {
202  pot.set_etat_zero() ;
203  return ;
204  }
205 
206  // Spectral coefficients of the source
207  assert(sourva.get_etat() == ETATQCQ) ;
208 
209  Valeur rho(sourva.get_mg()) ;
210  sourva.coef() ;
211  rho = *(sourva.c_cf) ; // copy of the coefficients of the source
212 
213  rho.ylm() ; // spherical harmonic transforms
214 
215  // Call to the Mtbl_cf version
216  // ---------------------------
217 
218  Mtbl_cf resu = sol_poisson_tau(*this, *(rho.c_cf), dzpuis) ;
219 
220  // Final result returned as a Scalar
221  // ----------------------------------
222 
223  pot.set_etat_zero() ; // to call Scalar::del_deriv().
224 
225  pot.set_etat_qcq() ;
226 
227  pot.set_spectral_va() = resu ;
228  pot.set_spectral_va().ylm_i() ; // On repasse en base standard.
229  pot.set_dzpuis(0) ;
230 }
231 
232 //=============================================================================
233 
234  // ---------------------------
235  // Cmp versions
236  // ---------------------------
237 
238  void Map_af::poisson(const Cmp& source, Param& , Cmp& pot, bool) const {
239 
240  assert(source.get_etat() != ETATNONDEF) ;
241  assert(source.get_mp()->get_mg() == mg) ;
242  assert(pot.get_mp()->get_mg() == mg) ;
243 
244  assert( source.check_dzpuis(2) || source.check_dzpuis(4)
245  || source.check_dzpuis(3) || source.check_dzpuis(5) ) ;
246 
247  bool match = true ;
248 
249  int dzpuis ;
250 
251  if ( (source.dz_nonzero()) || (source.get_dzpuis() > 3)) { //##awkward??
252  dzpuis = source.get_dzpuis() ;
253  }
254  else{
255  dzpuis = 4 ;
256  }
257 
258  match = !(dzpuis == 5) ;
259 
260  // Spherical harmonic expansion of the source
261  // ------------------------------------------
262 
263  const Valeur& sourva = source.va ;
264 
265  if (sourva.get_etat() == ETATZERO) {
266  pot.set_etat_zero() ;
267  return ;
268  }
269 
270  // Spectral coefficients of the source
271  assert(sourva.get_etat() == ETATQCQ) ;
272 
273  Valeur rho(sourva.get_mg()) ;
274  sourva.coef() ;
275  rho = *(sourva.c_cf) ; // copy of the coefficients of the source
276 
277  rho.ylm() ; // spherical harmonic transforms
278 
279  // Call to the Mtbl_cf version
280  // ---------------------------
281  Mtbl_cf resu = sol_poisson(*this, *(rho.c_cf), dzpuis, match) ;
282 
283  // Final result returned as a Cmp
284  // ------------------------------
285 
286  pot.set_etat_zero() ; // to call Cmp::del_t().
287 
288  pot.set_etat_qcq() ;
289 
290  pot.va = resu ;
291  (pot.va).ylm_i() ; // On repasse en base standard.
292 
293  (dzpuis == 5) ? pot.set_dzpuis(3) : pot.set_dzpuis(0) ;
294 
295 }
296 
297 
298  //----------------------
299  // Tau version method
300  //---------------------
301 
302 
303 void Map_af::poisson_tau(const Cmp& source, Param& , Cmp& pot) const {
304 
305  assert(source.get_etat() != ETATNONDEF) ;
306  assert(source.get_mp()->get_mg() == mg) ;
307  assert(pot.get_mp()->get_mg() == mg) ;
308 
309  assert( source.check_dzpuis(2) || source.check_dzpuis(4)
310  || source.check_dzpuis(3)) ;
311 
312  int dzpuis ;
313 
314  if ( (source.dz_nonzero()) || (source.get_dzpuis() > 3)) { //##awkward??
315  dzpuis = source.get_dzpuis() ;
316  }
317  else{
318  dzpuis = 4 ;
319  }
320 
321  // Spherical harmonic expansion of the source
322  // ------------------------------------------
323 
324  const Valeur& sourva = source.va ;
325 
326  if (sourva.get_etat() == ETATZERO) {
327  pot.set_etat_zero() ;
328  return ;
329  }
330 
331  // Spectral coefficients of the source
332  assert(sourva.get_etat() == ETATQCQ) ;
333 
334  Valeur rho(sourva.get_mg()) ;
335  sourva.coef() ;
336  rho = *(sourva.c_cf) ; // copy of the coefficients of the source
337 
338  rho.ylm() ; // spherical harmonic transforms
339 
340  // Call to the Mtbl_cf version
341  // ---------------------------
342 
343  Mtbl_cf resu = sol_poisson_tau(*this, *(rho.c_cf), dzpuis) ;
344 
345  // Final result returned as a Cmp
346  // ------------------------------
347 
348  pot.set_etat_zero() ; // to call Cmp::del_t().
349 
350  pot.set_etat_qcq() ;
351 
352  pot.va = resu ;
353  (pot.va).ylm_i() ; // On repasse en base standard.
354  pot.set_dzpuis(0) ;
355 }
356 
357 }
const Map * get_mp() const
Returns the mapping.
Definition: cmp.h:904
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
Component of a tensorial field *** DEPRECATED : use class Scalar instead ***.
Definition: cmp.h:449
void ylm_i()
Inverse of ylm()
Definition: valeur_ylm_i.C:134
void coef() const
Computes the coeffcients of *this.
Definition: valeur_coef.C:151
virtual void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: scalar.C:330
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
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:402
Values and coefficients of a (real-value) function.
Definition: valeur.h:297
int get_etat() const
Returns the logical state ETATNONDEF (undefined), ETATZERO (null) or ETATQCQ (ordinary).
Definition: scalar.h:569
virtual void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: scalar.C:359
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:814
const Mg3d * get_mg() const
Returns the Mg3d on which the this is defined.
Definition: valeur.h:763
int get_etat() const
Returns the logical state.
Definition: valeur.h:760
void set_etat_zero()
Sets the logical state to ETATZERO (zero).
Definition: cmp.C:292
int get_dzpuis() const
Returns dzpuis.
Definition: scalar.h:572
bool dz_nonzero() const
Returns true if the last domain is compactified and *this is not zero in this domain.
Definition: cmp.C:663
Parameter storage.
Definition: param.h:125
virtual void poisson(const Cmp &source, Param &par, Cmp &uu, bool verbose=true) const
Computes the solution of a scalar Poisson equation (Cmp version).
void set_etat_qcq()
Sets the logical state to ETATQCQ (ordinary state).
Definition: cmp.C:307
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).
const Mg3d * mg
Pointer on the multi-grid Mgd3 on which this is defined.
Definition: map.h:702
bool dz_nonzero() const
Returns true if the last domain is compactified and *this is not zero in this domain.
Definition: scalar.C:820
int get_dzpuis() const
Returns dzpuis.
Definition: cmp.h:906
Coefficients storage for the multi-domain spectral method.
Definition: mtbl_cf.h:196
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
Valeur & set_spectral_va()
Returns va (read/write version)
Definition: scalar.h:619
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: scalar.C:879
const Map & get_mp() const
Returns the mapping.
Definition: tensor.h:902
Valeur va
The numerical value of the Cmp.
Definition: cmp.h:467
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:616