This class contains the parameters needed to call the general elliptic solver. More...
#include <param_elliptic.h>
Public Member Functions | |
Param_elliptic (const Scalar &) | |
Standard constructor from a Scalar . | |
~Param_elliptic () | |
Destructor. | |
const Map_radial & | get_mp () const |
Returns the mapping. | |
double | get_alpha (int) const |
double | get_beta (int) const |
int | get_type (int) const |
void | set_helmholtz_minus (int zone, double mas, Scalar &so) |
Set the operator to in one domain (not in the nucleus). | |
void | set_poisson_pseudo_1d (Scalar &so) |
Set the operator to everywhere but in the compactified domain. | |
void | set_helmholtz_minus_pseudo_1d (int zone, double mas, Scalar &so) |
Set the operator to in one domain. | |
void | set_helmholtz_plus (int zone, double mas, Scalar &so) |
Set the operator to in one domain (only in the shells). | |
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. | |
void | set_helmholtz_minus_2d (int zone, double mas, const Scalar &) |
Set the 2D Helmholtz operator (with minus sign). | |
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 | set_sec_order (int zone, double a, double b, double c) |
Set the operator to in one domain (only in the shells). | |
void | set_ope_vorton (int zone, Scalar &so) |
Set the operator to in one domain (not implemented in the nucleus). | |
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_poisson_vect_eta (int zone) |
Sets the operator to be a regular elliptic operator to solve for the component of the vector Poisson equation. | |
void | set_poisson_tens_rr (int zone) |
Sets the operator to in all domains, for only. | |
void | inc_l_quant (int zone) |
Increases the quantum number l in the domain zone . | |
void | set_variable_F (const Scalar &) |
Changes the variable function F. | |
void | set_variable_G (const Scalar &) |
Changes the variable function G. | |
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 ;. | |
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 ;. | |
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 the domain zone ;. | |
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 the domain zone ;. | |
double | G_plus (int zone) const |
Returns the value of G, for a given angular point, at the outer boundary of the domain zone ;. | |
double | G_minus (int zone) const |
Returns the value of G, for a given angular point, at the inner boundary of the domain zone ;. | |
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 the domain zone ;. | |
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 the domain zone ;. | |
Protected Attributes | |
int | type_map |
Type of mapping either MAP_AFF or MAP_LOG. | |
const Map_af * | mp_af |
The mapping, if affine. | |
const Map_log * | mp_log |
The mapping if log type. | |
Ope_elementary ** | operateurs |
Array on the elementary operators. | |
Scalar | var_F |
Additive variable change function. | |
Scalar | var_G |
Multiplicative variable change that must be sphericaly symetric ! | |
Itbl | done_F |
Stores what has been computed for F . | |
Itbl | done_G |
Stores what has been computed for G . | |
Tbl | val_F_plus |
Values of F at the outer boundaries of the various domains. | |
Tbl | val_F_minus |
Values of F at the inner boundaries of the various domains. | |
Tbl | val_dF_plus |
Values of the derivative of F at the outer boundaries of the various domains. | |
Tbl | val_dF_minus |
Values of the derivative of F at the inner boundaries of the various domains. | |
Tbl | val_G_plus |
Values of G at the outer boundaries of the various domains. | |
Tbl | val_G_minus |
Values of G at the inner boundaries of the various domains. | |
Tbl | val_dG_plus |
Values of the derivative of G at the outer boundaries of the various domains. | |
Tbl | val_dG_minus |
Values of the derivative of G at the inner boundaries of the various domains. | |
Private Member Functions | |
void | compute_val_F (int, int, int) const |
Computes the various values of F . | |
void | compute_val_G (int) const |
Computes the various values of G . | |
Friends | |
Mtbl_cf | elliptic_solver (const Param_elliptic &, const Mtbl_cf &) |
Mtbl_cf | elliptic_solver_boundary (const Param_elliptic &, const Mtbl_cf &, const Mtbl_cf &bound, double fact_dir, double fact_neu) |
Mtbl_cf | elliptic_solver_no_zec (const Param_elliptic &, const Mtbl_cf &, double) |
Mtbl_cf | elliptic_solver_only_zec (const Param_elliptic &, const Mtbl_cf &, double) |
Mtbl_cf | elliptic_solver_sin_zec (const Param_elliptic &, const Mtbl_cf &, double *, double *) |
Mtbl_cf | elliptic_solver_fixe_der_zero (double, const Param_elliptic &, const Mtbl_cf &) |
void | Map_af::sol_elliptic (Param_elliptic &, const Scalar &, Scalar &) const |
void | Map_af::sol_elliptic_boundary (Param_elliptic &, const Scalar &, Scalar &, const Mtbl_cf &, double, double) const |
void | Map_af::sol_elliptic_boundary (Param_elliptic &, const Scalar &, Scalar &, const Scalar &, double, double) const |
void | Map_af::sol_elliptic_no_zec (Param_elliptic &, const Scalar &, Scalar &, double) const |
void | Map_af::sol_elliptic_only_zec (Param_elliptic &, const Scalar &, Scalar &, double) const |
void | Map_af::sol_elliptic_sin_zec (Param_elliptic &, const Scalar &, Scalar &, double *, double *) const |
void | Map_af::sol_elliptic_fixe_der_zero (double, Param_elliptic &, const Scalar &, Scalar &) const |
void | Map_af::sol_elliptic_2d (Param_elliptic &, const Scalar &, Scalar &) const |
void | Map_af::sol_elliptic_pseudo_1d (Param_elliptic &, const Scalar &, Scalar &) const |
void | Map_log::sol_elliptic (Param_elliptic &, const Scalar &, Scalar &) const |
void | Map_log::sol_elliptic_boundary (Param_elliptic &, const Scalar &, Scalar &, const Mtbl_cf &, double, double) const |
void | Map_log::sol_elliptic_boundary (Param_elliptic &, const Scalar &, Scalar &, const Scalar &, double, double) const |
void | Map_log::sol_elliptic_no_zec (Param_elliptic &, const Scalar &, Scalar &, double) const |
This class contains the parameters needed to call the general elliptic solver.
For every domain and every spherical harmonics, it contains the appropriate operator of type Ope_elementary
and the appropriate variable given by a Change_var
. ()
This class is only defined on an affine mapping Map_af
or a logarithmic one Map_log
Definition at line 131 of file param_elliptic.h.
Param_elliptic::Param_elliptic | ( | const Scalar & | so | ) |
Standard constructor from a Scalar
.
so | [parameter] type of the source of the elliptic equation. The actual values are not used but *this will be constructed using the same number of points, domains and symetry than so . |
This constructor initializes everything to solve a Poisson equation with non variable changes from domains to another.
Definition at line 179 of file param_elliptic.C.
References Itbl::annule_hard(), Scalar::annule_hard(), Valeur::base, Valeur::coef(), done_F, done_G, Scalar::dz_nonzero(), Valeur::get_base(), Scalar::get_dzpuis(), Map::get_mg(), get_mp(), Tensor::get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Map_log::get_type(), Base_val::give_quant_numbers(), mp_af, mp_log, operateurs, Valeur::set_base(), Tbl::set_etat_qcq(), Scalar::set_etat_qcq(), Scalar::set_spectral_va(), Scalar::std_spectral_base(), type_map, val_dF_minus, val_dF_plus, val_dG_minus, val_dG_plus, val_F_minus, val_F_plus, val_G_minus, val_G_plus, var_F, var_G, and Valeur::ylm().
Param_elliptic::~Param_elliptic | ( | ) |
Destructor.
Definition at line 306 of file param_elliptic.C.
References Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_nzone(), and operateurs.
void Param_elliptic::compute_val_F | ( | int | zone, | |
int | k, | |||
int | j | |||
) | const [private] |
Computes the various values of F
.
Definition at line 117 of file param_elliptic_val_lim.C.
References Tbl::annule_hard(), Valeur::base, Valeur::c_cf, done_F, Mtbl_cf::get_etat(), Map::get_mg(), get_mp(), Mg3d::get_nr(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), Itbl::set(), Tbl::set(), Tbl::set_etat_qcq(), val_dF_minus, val_dF_plus, val_F_minus, val_F_plus, and var_F.
void Param_elliptic::compute_val_G | ( | int | zone | ) | const [private] |
Computes the various values of G
.
Definition at line 156 of file param_elliptic_val_lim.C.
References Tbl::annule_hard(), Valeur::base, Valeur::c_cf, done_G, Mtbl_cf::get_etat(), Map::get_mg(), get_mp(), Mg3d::get_nr(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), Itbl::set(), Tbl::set(), Tbl::set_etat_qcq(), val_dG_minus, val_dG_plus, val_G_minus, val_G_plus, and var_G.
double Param_elliptic::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 the domain zone
;.
Definition at line 75 of file param_elliptic_val_lim.C.
References compute_val_F(), done_F, and val_dF_minus.
double Param_elliptic::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 the domain zone
;.
Definition at line 67 of file param_elliptic_val_lim.C.
References compute_val_F(), done_F, and val_dF_plus.
double Param_elliptic::dG_minus | ( | int | zone | ) | const |
Returns the value of the radial derivative of G, for a given angular point, at the inner boundary of the domain zone
;.
Definition at line 108 of file param_elliptic_val_lim.C.
References compute_val_G(), done_G, and val_dG_minus.
double Param_elliptic::dG_plus | ( | int | zone | ) | const |
Returns the value of the radial derivative of G, for a given angular point, at the outer boundary of the domain zone
;.
Definition at line 100 of file param_elliptic_val_lim.C.
References compute_val_G(), done_G, and val_dG_plus.
double Param_elliptic::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
;.
Definition at line 59 of file param_elliptic_val_lim.C.
References compute_val_F(), done_F, and val_F_minus.
double Param_elliptic::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
;.
Definition at line 51 of file param_elliptic_val_lim.C.
References compute_val_F(), done_F, and val_F_plus.
double Param_elliptic::G_minus | ( | int | zone | ) | const |
Returns the value of G, for a given angular point, at the inner boundary of the domain zone
;.
Definition at line 92 of file param_elliptic_val_lim.C.
References compute_val_G(), done_G, and val_G_minus.
double Param_elliptic::G_plus | ( | int | zone | ) | const |
Returns the value of G, for a given angular point, at the outer boundary of the domain zone
;.
Definition at line 84 of file param_elliptic_val_lim.C.
References compute_val_G(), done_G, and val_G_plus.
const Map_radial & Param_elliptic::get_mp | ( | ) | const |
Returns the mapping.
Definition at line 114 of file param_elliptic.C.
void Param_elliptic::inc_l_quant | ( | int | zone | ) |
Increases the quantum number l in the domain zone
.
Definition at line 319 of file param_elliptic.C.
References Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_nzone(), Ope_elementary::inc_l_quant(), and operateurs.
void Param_elliptic::set_helmholtz_minus | ( | int | zone, | |
double | mas, | |||
Scalar & | so | |||
) |
Set the operator to in one domain (not in the nucleus).
zone | [input] : the domain. | |
mas | [input] : the masse . | |
so | [input] : the source (used only to get the right basis). |
Definition at line 339 of file param_elliptic.C.
References Valeur::base, Valeur::coef(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Map_log::get_type(), Base_val::give_quant_numbers(), mp_log, operateurs, Scalar::set_spectral_va(), type_map, and Valeur::ylm().
void Param_elliptic::set_helmholtz_minus_2d | ( | int | zone, | |
double | mas, | |||
const Scalar & | source | |||
) |
Set the 2D Helmholtz operator (with minus sign).
zone | [input] : the domain. | |
mas | [input] : the masse parameter. |
Definition at line 96 of file param_elliptic_2d.C.
References Valeur::base, Scalar::check_dzpuis(), Scalar::get_dzpuis(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), operateurs, and type_map.
void Param_elliptic::set_helmholtz_minus_pseudo_1d | ( | int | zone, | |
double | mas, | |||
Scalar & | so | |||
) |
Set the operator to in one domain.
zone | [input] : the domain. | |
mas | [input] : the masse . | |
so | [input] : the source (used only to get the right basis). |
Definition at line 91 of file param_elliptic_pseudo_1d.C.
References Valeur::base, Scalar::check_dzpuis(), Scalar::get_dzpuis(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), operateurs, and type_map.
void Param_elliptic::set_helmholtz_plus | ( | int | zone, | |
double | mas, | |||
Scalar & | so | |||
) |
Set the operator to in one domain (only in the shells).
zone | [input] : the domain. | |
mas | [input] : the masse . | |
so | [input] : the source (used only to get the right basis). |
Definition at line 407 of file param_elliptic.C.
References Valeur::base, Valeur::coef(), Ope_elementary::get_base_r(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Map_log::get_type(), Base_val::give_quant_numbers(), mp_log, operateurs, R_CHEBI, Scalar::set_spectral_va(), type_map, and Valeur::ylm().
void Param_elliptic::set_ope_vorton | ( | int | zone, | |
Scalar & | so | |||
) |
Set the operator to in one domain (not implemented in the nucleus).
zone | [input] : the domain. | |
so | [input] : the source (used only to get the right basis). |
Definition at line 701 of file param_elliptic.C.
References Valeur::base, Valeur::coef(), Scalar::get_dzpuis(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Map_log::get_type(), Base_val::give_quant_numbers(), mp_log, operateurs, Scalar::set_spectral_va(), type_map, and Valeur::ylm().
void Param_elliptic::set_poisson_2d | ( | const Scalar & | source, | |
bool | indic = false | |||
) |
Set everything to do a 2d-Poisson, with or without l=0 (not put by default.
..)
Definition at line 55 of file param_elliptic_2d.C.
References Valeur::base, Scalar::get_dzpuis(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), operateurs, and type_map.
void Param_elliptic::set_poisson_pseudo_1d | ( | Scalar & | so | ) |
Set the operator to everywhere but in the compactified domain.
so | [input] : the source (used only to get the right basis). |
Definition at line 55 of file param_elliptic_pseudo_1d.C.
References Valeur::base, Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), operateurs, and type_map.
void Param_elliptic::set_poisson_tens_rr | ( | int | zone | ) |
Sets the operator to in all domains, for only.
zone | [input] : the domain. |
Definition at line 556 of file param_elliptic.C.
References Ope_elementary::get_base_r(), Ope_poisson::get_dzpuis(), Ope_poisson::get_lquant(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), operateurs, and type_map.
void Param_elliptic::set_poisson_vect_eta | ( | int | zone | ) |
Sets the operator to be a regular elliptic operator to solve for the component of the vector Poisson equation.
The operator is (Poisson with the decrease of l by one unit) in all domains but the CED, for ; it is not defined for l = 0. In the CED, the operator is also the Laplace one, but with l increased by one unit: . This is intended to solve the equation for arising in the decomposition of the vector Poisson equation.
zone | [input] : the domain. |
Definition at line 529 of file param_elliptic.C.
References Ope_poisson::dec_l_quant(), Ope_poisson::get_lquant(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_nzone(), Mg3d::get_type_r(), Ope_poisson::inc_l_quant(), and operateurs.
void Param_elliptic::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.
zone | [input] : the domain. | |
only_l_zero | [input] : the operator in built only for l=0 |
Definition at line 483 of file param_elliptic.C.
References Ope_elementary::get_base_r(), Ope_poisson::get_dzpuis(), Ope_poisson::get_lquant(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), operateurs, and type_map.
void Param_elliptic::set_sec_order | ( | int | zone, | |
double | a, | |||
double | b, | |||
double | c | |||
) |
Set the operator to in one domain (only in the shells).
zone | [input] : the domain. | |
a | [input] : the parameter . | |
b | [input] : the parameter . | |
c | [input] : the parameter . |
Definition at line 665 of file param_elliptic.C.
References Ope_elementary::get_base_r(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Map_log::get_type(), mp_log, operateurs, R_CHEBI, and type_map.
void Param_elliptic::set_sec_order_r2 | ( | int | zone, | |
double | a, | |||
double | b, | |||
double | c | |||
) |
Set the operator to in one domain (only in the shells).
zone | [input] : the domain. | |
a | [input] : the parameter . | |
b | [input] : the parameter . | |
c | [input] : the parameter . |
Definition at line 598 of file param_elliptic.C.
References Ope_elementary::get_base_r(), Map::get_mg(), get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Map_log::get_type(), mp_log, operateurs, R_CHEBI, and type_map.
void Param_elliptic::set_variable_F | ( | const Scalar & | so | ) |
Changes the variable function F.
Definition at line 770 of file param_elliptic.C.
References Itbl::annule_hard(), done_F, Scalar::get_etat(), get_mp(), Tensor::get_mp(), and var_F.
void Param_elliptic::set_variable_G | ( | const Scalar & | so | ) |
Changes the variable function G.
Definition at line 779 of file param_elliptic.C.
References Itbl::annule_hard(), done_G, Scalar::get_etat(), get_mp(), Tensor::get_mp(), and var_G.
Itbl Param_elliptic::done_F [mutable, protected] |
Stores what has been computed for F
.
Definition at line 143 of file param_elliptic.h.
Itbl Param_elliptic::done_G [mutable, protected] |
Stores what has been computed for G
.
Definition at line 144 of file param_elliptic.h.
const Map_af* Param_elliptic::mp_af [protected] |
The mapping, if affine.
Definition at line 135 of file param_elliptic.h.
const Map_log* Param_elliptic::mp_log [protected] |
The mapping if log type.
Definition at line 136 of file param_elliptic.h.
Ope_elementary** Param_elliptic::operateurs [protected] |
Array on the elementary operators.
Definition at line 138 of file param_elliptic.h.
int Param_elliptic::type_map [protected] |
Type of mapping either MAP_AFF or MAP_LOG.
Definition at line 134 of file param_elliptic.h.
Tbl Param_elliptic::val_dF_minus [mutable, protected] |
Values of the derivative of F at the inner boundaries of the various domains.
Definition at line 148 of file param_elliptic.h.
Tbl Param_elliptic::val_dF_plus [mutable, protected] |
Values of the derivative of F at the outer boundaries of the various domains.
Definition at line 147 of file param_elliptic.h.
Tbl Param_elliptic::val_dG_minus [mutable, protected] |
Values of the derivative of G at the inner boundaries of the various domains.
Definition at line 152 of file param_elliptic.h.
Tbl Param_elliptic::val_dG_plus [mutable, protected] |
Values of the derivative of G at the outer boundaries of the various domains.
Definition at line 151 of file param_elliptic.h.
Tbl Param_elliptic::val_F_minus [mutable, protected] |
Values of F at the inner boundaries of the various domains.
Definition at line 146 of file param_elliptic.h.
Tbl Param_elliptic::val_F_plus [mutable, protected] |
Values of F at the outer boundaries of the various domains.
Definition at line 145 of file param_elliptic.h.
Tbl Param_elliptic::val_G_minus [mutable, protected] |
Values of G at the inner boundaries of the various domains.
Definition at line 150 of file param_elliptic.h.
Tbl Param_elliptic::val_G_plus [mutable, protected] |
Values of G at the outer boundaries of the various domains.
Definition at line 149 of file param_elliptic.h.
Scalar Param_elliptic::var_F [protected] |
Additive variable change function.
Definition at line 140 of file param_elliptic.h.
Scalar Param_elliptic::var_G [protected] |
Multiplicative variable change that must be sphericaly symetric !
Definition at line 141 of file param_elliptic.h.