26 #ifndef __PARAM_ELLIPTIC_H_ 27 #define __PARAM_ELLIPTIC_H_ 116 #include "ope_elementary.h" 178 double get_alpha (
int)
const ;
179 double get_beta (
int)
const ;
180 int get_type (
int)
const ;
251 void set_sec_order (
int zone,
double a,
double b,
double c) ;
315 double F_plus (
int zone,
int k,
int j)
const ;
321 double F_minus (
int zone,
int k,
int j)
const ;
328 double dF_plus (
int zone,
int k,
int j)
const ;
335 double dF_minus (
int zone,
int k,
int j)
const ;
341 double G_plus (
int zone)
const ;
347 double G_minus (
int zone)
const ;
354 double dG_plus (
int zone)
const ;
366 friend Mtbl_cf elliptic_solver_no_zec
368 friend Mtbl_cf elliptic_solver_only_zec
370 friend Mtbl_cf elliptic_solver_sin_zec
372 friend Mtbl_cf elliptic_solver_fixe_der_zero
377 double ,
double )
const ;
379 double ,
double )
const ;
392 double,
double)
const ;
395 double,
double)
const ;
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 ...
void sol_elliptic_only_zec(Param_elliptic ¶ms, 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).
void sol_elliptic_2d(Param_elliptic &, const Scalar &, Scalar &) const
General elliptic solver in a 2D case.
Basic integer array class.
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.
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 ¶ms, 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 ¶ms, 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.
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 ¶ms, const Scalar &so, Scalar &uu, double) const
General elliptic solver.
void sol_elliptic(Param_elliptic ¶ms, 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 ...
void sol_elliptic_fixe_der_zero(double val, Param_elliptic ¶ms, 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.
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 ¶ms, const Scalar &so, Scalar &uu, const Mtbl_cf &bound, double fact_dir, double fact_neu) const
General elliptic solver including inner boundary conditions.
void compute_val_G(int) const
Computes the various values of G.
void sol_elliptic_boundary(Param_elliptic ¶ms, 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 ¶ms, 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).