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