LORENE
map_af_poisson_angu.C
1 /*
2  * Resolution of the angular Poisson equation.
3  *
4  * (see file map.h for documentation)
5  *
6  */
7 
8 /*
9  * Copyright (c) 2003-2005 Eric Gourgoulhon & Jerome Novak
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 version 2
15  * as published by the Free Software Foundation.
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  * $Id: map_af_poisson_angu.C,v 1.6 2024/01/26 17:44:25 g_servignat Exp $
32  * $Log: map_af_poisson_angu.C,v $
33  * Revision 1.6 2024/01/26 17:44:25 g_servignat
34  * Updated the Pseudopolytrope_1D class to be consistent with the paper (i.e. with a GPP in the middle)
35  *
36  * Revision 1.5 2016/12/05 16:17:57 j_novak
37  * Suppression of some global variables (file names, loch, ...) to prevent redefinitions
38  *
39  * Revision 1.4 2014/10/13 08:53:03 j_novak
40  * Lorene classes and functions now belong to the namespace Lorene.
41  *
42  * Revision 1.3 2005/04/04 21:31:31 e_gourgoulhon
43  * Added argument lambda to method poisson_angu
44  * to deal with the generalized angular Poisson equation:
45  * Lap_ang u + lambda u = source.
46  *
47  * Revision 1.2 2003/10/16 08:49:23 j_novak
48  * Added a flag to decide wether the output is in the Ylm or in the standard base.
49  *
50  * Revision 1.1 2003/10/15 21:11:26 e_gourgoulhon
51  * Added method poisson_angu.
52  *
53  *
54  *
55  * $Header: /cvsroot/Lorene/C++/Source/Map/map_af_poisson_angu.C,v 1.6 2024/01/26 17:44:25 g_servignat Exp $
56  *
57  */
58 
59 // Lorene headers
60 #include "tensor.h"
61 #include "param.h"
62 
63 namespace Lorene {
64 void Map_af::poisson_angu(const Scalar& source, Param& par, Scalar& uu,
65  double lambda) const {
66 
67  assert(source.get_etat() != ETATNONDEF) ;
68 
69  assert(&(source.get_mp()) == this ) ;
70  assert(&(uu.get_mp()) == this ) ;
71 
72  // Spherical harmonic expansion of the source
73  // ------------------------------------------
74 
75  const Valeur& sourva = source.get_spectral_va() ;
76 
77  if (sourva.get_etat() == ETATZERO) {
78  uu.set_etat_zero() ;
79  return ;
80  }
81 
82  // Spectral coefficients of the source
83  assert(sourva.get_etat() == ETATQCQ) ;
84  sourva.coef() ;
85 
86  Valeur resu(mg) ;
87  resu = *(sourva.c_cf) ; // copy of the coefficients of the source
88 
89  resu.ylm() ; // spherical harmonic transform
90 
91  // Call to the Mtbl_cf version
92  // ---------------------------
93  (resu.c_cf)->poisson_angu(lambda) ;
94 
95  if (par.get_n_int() == 0) resu.ylm_i() ; // Back to standard bases
96  //in the case of no flag present
97  // in the Param
98 
99  // Final result returned as a Scalar
100  // ---------------------------------
101 
102  uu = resu ;
103 
104  uu.set_dzpuis( source.get_dzpuis() ) ; // dzpuis unchanged
105 }
106 }
Mtbl_cf * c_cf
Coefficients of the spectral expansion of the function.
Definition: valeur.h:312
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
void ylm()
Computes the coefficients of *this.
Definition: valeur_ylm.C:141
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
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:560
void set_dzpuis(int)
Modifies the dzpuis flag.
Definition: scalar.C:814
int get_etat() const
Returns the logical state.
Definition: valeur.h:760
int get_n_int() const
Returns the number of stored int 's addresses.
Definition: param.C:242
int get_dzpuis() const
Returns dzpuis.
Definition: scalar.h:563
Parameter storage.
Definition: param.h:125
const Mg3d * mg
Pointer on the multi-grid Mgd3 on which this is defined.
Definition: map.h:688
virtual void poisson_angu(const Scalar &source, Param &par, Scalar &uu, double lambda=0) const
Computes the solution of the generalized angular Poisson equation.
const Map & get_mp() const
Returns the mapping.
Definition: tensor.h:874
const Valeur & get_spectral_va() const
Returns va (read only version)
Definition: scalar.h:607