LORENE
param_elliptic.h
1 /*
2  * Definition of Lorene class Param_elliptic
3  *
4  */
5 
6 /*
7  * Copyright (c) 2003 Philippe Grandclement
8  *
9  * This file is part of LORENE.
10  *
11  * LORENE is free software; you can redistribute it and/or modify
12  * it under the terms of the GNU General Public License version 2
13  * as published by the Free Software Foundation.
14  *
15  * LORENE is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18  * GNU General Public License for more details.
19  *
20  * You should have received a copy of the GNU General Public License
21  * along with LORENE; if not, write to the Free Software
22  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
23  *
24  */
25 
26 #ifndef __PARAM_ELLIPTIC_H_
27 #define __PARAM_ELLIPTIC_H_
28 
29 /*
30  * $Id: param_elliptic.h,v 1.21 2014/10/13 08:52:36 j_novak Exp $
31  * $Log: param_elliptic.h,v $
32  * Revision 1.21 2014/10/13 08:52:36 j_novak
33  * Lorene classes and functions now belong to the namespace Lorene.
34  *
35  * Revision 1.20 2007/05/06 10:48:08 p_grandclement
36  * Modification of a few operators for the vorton project
37  *
38  * Revision 1.19 2007/04/24 09:04:11 p_grandclement
39  * Addition of an operator for the vortons
40  *
41  * Revision 1.18 2007/01/16 15:05:59 n_vasset
42  * New constructor (taking a Scalar in mono-domain angular grid for
43  * boundary) for function sol_elliptic_boundary
44  *
45  * Revision 1.17 2005/11/30 11:09:03 p_grandclement
46  * Changes for the Bin_ns_bh project
47  *
48  * Revision 1.16 2005/08/26 14:02:38 p_grandclement
49  * Modification of the elliptic solver that matches with an oscillatory exterior solution
50  * small correction in Poisson tau also...
51  *
52  * Revision 1.15 2005/06/09 07:56:25 f_limousin
53  * Implement a new function sol_elliptic_boundary() and
54  * Vector::poisson_boundary(...) which solve the vectorial poisson
55  * equation (method 6) with an inner boundary condition.
56  *
57  * Revision 1.14 2005/02/15 15:43:16 j_novak
58  * First version of the block inversion for the vector Poisson equation (method 6).
59  *
60  * Revision 1.13 2004/12/23 16:30:14 j_novak
61  * New files and class for the solution of the rr component of the tensor Poisson
62  * equation.
63  *
64  * Revision 1.12 2004/08/24 09:14:40 p_grandclement
65  * Addition of some new operators, like Poisson in 2d... It now requieres the
66  * GSL library to work.
67  *
68  * Also, the way a variable change is stored by a Param_elliptic is changed and
69  * no longer uses Change_var but rather 2 Scalars. The codes using that feature
70  * will requiere some modification. (It should concern only the ones about monopoles)
71  *
72  * Revision 1.11 2004/06/22 08:49:57 p_grandclement
73  * Addition of everything needed for using the logarithmic mapping
74  *
75  * Revision 1.10 2004/06/14 15:23:07 j_novak
76  * Modif. comments.
77  *
78  * Revision 1.9 2004/06/14 15:07:10 j_novak
79  * New methods for the construction of the elliptic operator appearing in
80  * the vector Poisson equation (acting on eta).
81  *
82  * Revision 1.8 2004/05/14 08:51:00 p_grandclement
83  * *** empty log message ***
84  *
85  * Revision 1.7 2004/05/10 15:28:21 j_novak
86  * First version of functions for the solution of the r-component of the
87  * vector Poisson equation.
88  *
89  * Revision 1.6 2004/03/23 14:54:45 j_novak
90  * More documentation
91  *
92  * Revision 1.5 2004/03/17 15:58:47 p_grandclement
93  * Slight modification of sol_elliptic_no_zec
94  *
95  * Revision 1.4 2004/03/05 09:18:48 p_grandclement
96  * Addition of operator sec_order_r2
97  *
98  * Revision 1.3 2004/02/11 09:47:44 p_grandclement
99  * Addition of a new elliptic solver, matching with the homogeneous solution
100  * at the outer shell and not solving in the external domain (more details
101  * coming soon ; check your local Lorene dealer...)
102  *
103  * Revision 1.2 2004/01/28 16:46:22 p_grandclement
104  * Addition of the sol_elliptic_fixe_der_zero stuff
105  *
106  * Revision 1.1 2003/12/11 14:57:00 p_grandclement
107  * I had forgotten the .h (sorry folks...)
108  *
109  * Revision 1.2 2001/12/11 06:44:41 e_gourgoulhon
110  *
111  * $Header: /cvsroot/Lorene/C++/Include/param_elliptic.h,v 1.21 2014/10/13 08:52:36 j_novak Exp $
112  *
113  */
114 
115 #include "map.h"
116 #include "ope_elementary.h"
117 #include "scalar.h"
118 
119 #define MAP_AFF 0
120 #define MAP_LOG 1
121 
122 namespace Lorene {
136 
137  protected:
138  int type_map ;
139  const Map_af* mp_af ;
140  const Map_log* mp_log ;
141 
143 
146 
147  mutable Itbl done_F ;
148  mutable Itbl done_G ;
149  mutable Tbl val_F_plus ;
150  mutable Tbl val_F_minus ;
151  mutable Tbl val_dF_plus ;
152  mutable Tbl val_dF_minus ;
153  mutable Tbl val_G_plus ;
154  mutable Tbl val_G_minus ;
155  mutable Tbl val_dG_plus ;
156  mutable Tbl val_dG_minus ;
157 
158  private:
159  void compute_val_F(int, int, int) const ;
160  void compute_val_G(int) const ;
161 
162  public:
173  Param_elliptic (const Scalar&) ;
174  ~Param_elliptic() ;
175 
177  const Map_radial& get_mp() const ;
178  double get_alpha (int) const ;
179  double get_beta (int) const ;
180  int get_type (int) const ;
181 
182  public:
191  void set_helmholtz_minus (int zone, double mas, Scalar& so) ;
198  void set_poisson_pseudo_1d (Scalar& so) ;
206  void set_helmholtz_minus_pseudo_1d (int zone, double mas, Scalar& so) ;
207 
216  void set_helmholtz_plus (int zone, double mas, Scalar& so) ;
217 
221  void set_poisson_2d (const Scalar &, bool indic = false) ;
222 
229  void set_helmholtz_minus_2d (int zone, double mas, const Scalar&) ;
230 
240  void set_sec_order_r2 (int zone, double a, double b, double c) ;
241 
251  void set_sec_order (int zone, double a, double b, double c) ;
252 
259  void set_ope_vorton (int zone, Scalar& so) ;
269  void set_poisson_vect_r(int zone, bool only_l_zero = false) ;
270 
286  void set_poisson_vect_eta(int zone) ;
287 
294  void set_poisson_tens_rr(int zone) ;
295 
299  void inc_l_quant (int zone) ;
300 
304  void set_variable_F (const Scalar&) ;
305 
309  void set_variable_G (const Scalar&) ;
310 
315  double F_plus (int zone, int k, int j) const ;
316 
321  double F_minus (int zone, int k, int j) const ;
322 
328  double dF_plus (int zone, int k, int j) const ;
329 
335  double dF_minus (int zone, int k, int j) const ;
336 
341  double G_plus (int zone) const ;
342 
347  double G_minus (int zone) const ;
348 
354  double dG_plus (int zone) const ;
355 
361  double dG_minus (int zone) const ;
362 
363  // A lot of friend functions... Possiblement pas toutes utiles...
364  friend Mtbl_cf elliptic_solver (const Param_elliptic&, const Mtbl_cf&) ;
365  friend Mtbl_cf elliptic_solver_boundary (const Param_elliptic&, const Mtbl_cf&, const Mtbl_cf& bound, double fact_dir, double fact_neu ) ;
366  friend Mtbl_cf elliptic_solver_no_zec
367  (const Param_elliptic&, const Mtbl_cf&, double) ;
368  friend Mtbl_cf elliptic_solver_only_zec
369  (const Param_elliptic&, const Mtbl_cf&, double) ;
370  friend Mtbl_cf elliptic_solver_sin_zec
371  (const Param_elliptic&, const Mtbl_cf&, double*, double*) ;
372  friend Mtbl_cf elliptic_solver_fixe_der_zero
373  (double, const Param_elliptic&, const Mtbl_cf&) ;
374 
375  friend void Map_af::sol_elliptic(Param_elliptic&, const Scalar&, Scalar&) const ;
376  friend void Map_af::sol_elliptic_boundary(Param_elliptic&, const Scalar&, Scalar&, const Mtbl_cf& ,
377  double , double ) const ;
378  friend void Map_af::sol_elliptic_boundary(Param_elliptic&, const Scalar&, Scalar&, const Scalar& ,
379  double , double ) const ;
380 
381 
382  friend void Map_af::sol_elliptic_no_zec(Param_elliptic&, const Scalar&, Scalar&, double) const ;
383  friend void Map_af::sol_elliptic_only_zec(Param_elliptic&, const Scalar&, Scalar&, double) const ;
384  friend void Map_af::sol_elliptic_sin_zec(Param_elliptic&, const Scalar&, Scalar&, double*, double*) const ;
385  friend void Map_af::sol_elliptic_fixe_der_zero(double, Param_elliptic&, const Scalar&, Scalar&) const ;
386 
387  friend void Map_af::sol_elliptic_2d(Param_elliptic&, const Scalar&, Scalar&) const ;
388  friend void Map_af::sol_elliptic_pseudo_1d(Param_elliptic&, const Scalar&, Scalar&) const ;
389 
390  friend void Map_log::sol_elliptic(Param_elliptic&, const Scalar&, Scalar&) const ;
391  friend void Map_log::sol_elliptic_boundary(Param_elliptic&, const Scalar&, Scalar&, const Mtbl_cf&,
392  double, double) const ;
393 
394  friend void Map_log::sol_elliptic_boundary(Param_elliptic&, const Scalar&, Scalar&, const Scalar&,
395  double, double) const ;
396 
397 
398  friend void Map_log::sol_elliptic_no_zec(Param_elliptic&, const Scalar&, Scalar&, double) const ;
399 } ;
400 
401 }
402 #endif
void set_poisson_vect_eta(int zone)
Sets the operator to be a regular elliptic operator to solve for the component of the vector Poisson...
void set_helmholtz_minus(int zone, double mas, Scalar &so)
Set the operator to in one domain (not in the nucleus).
void set_variable_F(const Scalar &)
Changes the variable function F.
double dF_plus(int zone, int k, int j) const
Returns the value of the radial derivative of F, for a given angular point, at the outer boundary of ...
double dG_plus(int zone) const
Returns the value of the radial derivative of G, for a given angular point, at the outer boundary of ...
Lorene prototypes.
Definition: app_hor.h:67
void sol_elliptic_only_zec(Param_elliptic &params, const Scalar &so, Scalar &uu, double val) const
General elliptic solver.
void set_ope_vorton(int zone, Scalar &so)
Set the operator to in one domain (not implemented in the nucleus).
Tensor field of valence 0 (or component of a tensorial field).
Definition: scalar.h:393
void sol_elliptic_2d(Param_elliptic &, const Scalar &, Scalar &) const
General elliptic solver in a 2D case.
Basic integer array class.
Definition: itbl.h:122
Itbl done_G
Stores what has been computed for G.
Param_elliptic(const Scalar &)
Standard constructor from a Scalar.
void compute_val_F(int, int, int) const
Computes the various values of F.
Tbl val_F_minus
Values of F at the inner boundaries of the various domains.
const Map_af * mp_af
The mapping, if affine.
double G_plus(int zone) const
Returns the value of G, for a given angular point, at the outer boundary of the domain zone ;...
Tbl val_G_plus
Values of G at the outer boundaries of the various domains.
Scalar var_F
Additive variable change function.
Logarithmic radial mapping.
Definition: map.h:3616
double dF_minus(int zone, int k, int j) const
Returns the value of the radial derivative of F, for a given angular point, at the inner boundary of ...
void sol_elliptic(Param_elliptic &params, const Scalar &so, Scalar &uu) const
General elliptic solver.
Tbl val_G_minus
Values of G at the inner boundaries of the various domains.
int type_map
Type of mapping either MAP_AFF or MAP_LOG.
Tbl val_F_plus
Values of F at the outer boundaries of the various domains.
void sol_elliptic_no_zec(Param_elliptic &params, const Scalar &so, Scalar &uu, double val) const
General elliptic solver.
double F_minus(int zone, int k, int j) const
Returns the value of F, for a given angular point, at the inner boundary of the domain zone ;...
void set_variable_G(const Scalar &)
Changes the variable function G.
Base class for pure radial mappings.
Definition: map.h:1557
void set_helmholtz_minus_2d(int zone, double mas, const Scalar &)
Set the 2D Helmholtz operator (with minus sign)
const Map_radial & get_mp() const
Returns the mapping.
void set_poisson_tens_rr(int zone)
Sets the operator to in all domains, for only.
Scalar var_G
Multiplicative variable change that must be sphericaly symetric !
double G_minus(int zone) const
Returns the value of G, for a given angular point, at the inner boundary of the domain zone ;...
void set_helmholtz_plus(int zone, double mas, Scalar &so)
Set the operator to in one domain (only in the shells).
This class contains the parameters needed to call the general elliptic solver.
Tbl val_dG_minus
Values of the derivative of G at the inner boundaries of the various domains.
~Param_elliptic()
Destructor.
Ope_elementary ** operateurs
Array on the elementary operators.
const Map_log * mp_log
The mapping if log type.
void sol_elliptic_no_zec(Param_elliptic &params, const Scalar &so, Scalar &uu, double) const
General elliptic solver.
void sol_elliptic(Param_elliptic &params, const Scalar &so, Scalar &uu) const
General elliptic solver.
Tbl val_dG_plus
Values of the derivative of G at the outer boundaries of the various domains.
void set_poisson_2d(const Scalar &, bool indic=false)
Set everything to do a 2d-Poisson, with or without l=0 (not put by default...)
double dG_minus(int zone) const
Returns the value of the radial derivative of G, for a given angular point, at the inner boundary of ...
Affine radial mapping.
Definition: map.h:2048
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 ...
void inc_l_quant(int zone)
Increases the quantum number l in the domain zone .
Basic class for elementary elliptic operators.
void set_poisson_pseudo_1d(Scalar &so)
Set the operator to everywhere but in the compactified domain.
Tbl val_dF_minus
Values of the derivative of F at the inner boundaries of the various domains.
Coefficients storage for the multi-domain spectral method.
Definition: mtbl_cf.h:196
Itbl done_F
Stores what has been computed for F.
void set_sec_order_r2(int zone, double a, double b, double c)
Set the operator to in one domain (only in the shells).
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.
Basic array class.
Definition: tbl.h:164
void compute_val_G(int) const
Computes the various values of G.
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.
Tbl val_dF_plus
Values of the derivative of F at the outer boundaries of the various domains.
void sol_elliptic_sin_zec(Param_elliptic &params, const Scalar &so, Scalar &uu, double *coefs, double *) const
General elliptic solver.
void set_helmholtz_minus_pseudo_1d(int zone, double mas, Scalar &so)
Set the operator to in one domain.
double F_plus(int zone, int k, int j) const
Returns the value of F, for a given angular point, at the outer boundary of the domain zone ;...
void set_poisson_vect_r(int zone, bool only_l_zero=false)
Sets the operator to in all domains, for ; and to in all domains otherwise.
void set_sec_order(int zone, double a, double b, double c)
Set the operator to in one domain (only in the shells).