LORENE
|
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 ![]() | |
void | set_poisson_pseudo_1d (Scalar &so) |
Set the operator to ![]() | |
void | set_helmholtz_minus_pseudo_1d (int zone, double mas, Scalar &so) |
Set the operator to ![]() | |
void | set_helmholtz_plus (int zone, double mas, Scalar &so) |
Set the operator to ![]() | |
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 ![]() | |
void | set_sec_order (int zone, double a, double b, double c) |
Set the operator to ![]() | |
void | set_ope_vorton (int zone, Scalar &so) |
Set the operator to ![]() | |
void | set_poisson_vect_r (int zone, bool only_l_zero=false) |
Sets the operator to ![]() ![]() ![]() | |
void | set_poisson_vect_eta (int zone) |
Sets the operator to be a regular elliptic operator to solve for the ![]() | |
void | set_poisson_tens_rr (int zone) |
Sets the operator to ![]() ![]() | |
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 . | |
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 135 of file param_elliptic.h.
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 183 of file param_elliptic.C.
References Lorene::Itbl::annule_hard(), Lorene::Scalar::annule_hard(), done_F, done_G, get_mp(), mp_af, mp_log, operateurs, Lorene::Valeur::set_base(), Lorene::Scalar::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::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, and var_G.
Lorene::Param_elliptic::~Param_elliptic | ( | ) |
Computes the various values of F
.
Definition at line 121 of file param_elliptic_val_lim.C.
References done_F, get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Itbl::set(), Lorene::Tbl::set(), val_dF_minus, val_dF_plus, val_F_minus, val_F_plus, and var_F.
Computes the various values of G
.
Definition at line 160 of file param_elliptic_val_lim.C.
References done_G, get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Itbl::set(), Lorene::Tbl::set(), val_dG_minus, val_dG_plus, val_G_minus, val_G_plus, and var_G.
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 79 of file param_elliptic_val_lim.C.
References compute_val_F(), done_F, and val_dF_minus.
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 71 of file param_elliptic_val_lim.C.
References compute_val_F(), done_F, and val_dF_plus.
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 112 of file param_elliptic_val_lim.C.
References compute_val_G(), done_G, and val_dG_minus.
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 104 of file param_elliptic_val_lim.C.
References compute_val_G(), done_G, and val_dG_plus.
Returns the value of F, for a given angular point, at the inner boundary of the domain zone
;.
Definition at line 63 of file param_elliptic_val_lim.C.
References compute_val_F(), done_F, and val_F_minus.
Returns the value of F, for a given angular point, at the outer boundary of the domain zone
;.
Definition at line 55 of file param_elliptic_val_lim.C.
References compute_val_F(), done_F, and val_F_plus.
Returns the value of G, for a given angular point, at the inner boundary of the domain zone
;.
Definition at line 96 of file param_elliptic_val_lim.C.
References compute_val_G(), done_G, and val_G_minus.
Returns the value of G, for a given angular point, at the outer boundary of the domain zone
;.
Definition at line 88 of file param_elliptic_val_lim.C.
References compute_val_G(), done_G, and val_G_plus.
Definition at line 134 of file param_elliptic.C.
Definition at line 150 of file param_elliptic.C.
const Map_radial & Lorene::Param_elliptic::get_mp | ( | ) | const |
Returns the mapping.
Definition at line 118 of file param_elliptic.C.
Definition at line 166 of file param_elliptic.C.
Increases the quantum number l in the domain zone
.
Definition at line 323 of file param_elliptic.C.
References get_mp(), Lorene::Ope_elementary::inc_l_quant(), and operateurs.
Set the operator to
zone | [input] : the domain. |
mas | [input] : the masse ![]() |
so | [input] : the source (used only to get the right basis). |
Definition at line 343 of file param_elliptic.C.
References get_mp(), mp_log, operateurs, and type_map.
Set the 2D Helmholtz operator (with minus sign)
zone | [input] : the domain. |
mas | [input] : the masse parameter. |
Definition at line 100 of file param_elliptic_2d.C.
References get_mp(), operateurs, and type_map.
Set the operator to
zone | [input] : the domain. |
mas | [input] : the masse ![]() |
so | [input] : the source (used only to get the right basis). |
Definition at line 95 of file param_elliptic_pseudo_1d.C.
References get_mp(), operateurs, and type_map.
Set the operator to
zone | [input] : the domain. |
mas | [input] : the masse ![]() |
so | [input] : the source (used only to get the right basis). |
Definition at line 411 of file param_elliptic.C.
References Lorene::Ope_elementary::get_base_r(), get_mp(), mp_log, operateurs, R_CHEBI, and type_map.
Set the operator to
zone | [input] : the domain. |
so | [input] : the source (used only to get the right basis). |
Definition at line 705 of file param_elliptic.C.
References get_mp(), mp_log, operateurs, and type_map.
Set everything to do a 2d-Poisson, with or without l=0 (not put by default...)
Definition at line 59 of file param_elliptic_2d.C.
References get_mp(), operateurs, and type_map.
Set the operator to
so | [input] : the source (used only to get the right basis). |
Definition at line 59 of file param_elliptic_pseudo_1d.C.
References get_mp(), operateurs, and type_map.
Sets the operator to
zone | [input] : the domain. |
Definition at line 560 of file param_elliptic.C.
References Lorene::Ope_elementary::get_base_r(), get_mp(), operateurs, and type_map.
Sets the operator to be a regular elliptic operator to solve for the
The operator is
zone | [input] : the domain. |
Definition at line 533 of file param_elliptic.C.
References get_mp(), and operateurs.
Sets the operator to
zone | [input] : the domain. |
only_l_zero | [input] : the operator in built only for l=0 |
Definition at line 487 of file param_elliptic.C.
References Lorene::Ope_elementary::get_base_r(), get_mp(), operateurs, and type_map.
Set the operator to
zone | [input] : the domain. |
a | [input] : the parameter ![]() |
b | [input] : the parameter ![]() |
c | [input] : the parameter ![]() |
Definition at line 669 of file param_elliptic.C.
References Lorene::Ope_elementary::get_base_r(), get_mp(), mp_log, operateurs, R_CHEBI, and type_map.
Set the operator to
zone | [input] : the domain. |
a | [input] : the parameter ![]() |
b | [input] : the parameter ![]() |
c | [input] : the parameter ![]() |
Definition at line 602 of file param_elliptic.C.
References Lorene::Ope_elementary::get_base_r(), get_mp(), mp_log, operateurs, R_CHEBI, and type_map.
Changes the variable function F.
Definition at line 774 of file param_elliptic.C.
References Lorene::Itbl::annule_hard(), done_F, get_mp(), and var_F.
Changes the variable function G.
Definition at line 783 of file param_elliptic.C.
References Lorene::Itbl::annule_hard(), done_G, get_mp(), and var_G.
|
friend |
Definition at line 84 of file sol_elliptic.C.
|
friend |
Definition at line 65 of file sol_elliptic_boundary.C.
|
friend |
Definition at line 46 of file sol_elliptic_fixe_der_zero.C.
|
friend |
Definition at line 73 of file sol_elliptic_no_zec.C.
|
friend |
Definition at line 66 of file sol_elliptic_only_zec.C.
|
friend |
Definition at line 76 of file sol_elliptic_sin_zec.C.
|
mutableprotected |
Stores what has been computed for F
.
Definition at line 147 of file param_elliptic.h.
|
mutableprotected |
Stores what has been computed for G
.
Definition at line 148 of file param_elliptic.h.
The mapping, if affine.
Definition at line 139 of file param_elliptic.h.
The mapping if log type.
Definition at line 140 of file param_elliptic.h.
|
protected |
Array on the elementary operators.
Definition at line 142 of file param_elliptic.h.
|
protected |
Type of mapping either MAP_AFF or MAP_LOG.
Definition at line 138 of file param_elliptic.h.
|
mutableprotected |
Values of the derivative of F at the inner boundaries of the various domains.
Definition at line 152 of file param_elliptic.h.
|
mutableprotected |
Values of the derivative of F at the outer boundaries of the various domains.
Definition at line 151 of file param_elliptic.h.
|
mutableprotected |
Values of the derivative of G at the inner boundaries of the various domains.
Definition at line 156 of file param_elliptic.h.
|
mutableprotected |
Values of the derivative of G at the outer boundaries of the various domains.
Definition at line 155 of file param_elliptic.h.
|
mutableprotected |
Values of F at the inner boundaries of the various domains.
Definition at line 150 of file param_elliptic.h.
|
mutableprotected |
Values of F at the outer boundaries of the various domains.
Definition at line 149 of file param_elliptic.h.
|
mutableprotected |
Values of G at the inner boundaries of the various domains.
Definition at line 154 of file param_elliptic.h.
|
mutableprotected |
Values of G at the outer boundaries of the various domains.
Definition at line 153 of file param_elliptic.h.
|
protected |
Additive variable change function.
Definition at line 144 of file param_elliptic.h.
|
protected |
Multiplicative variable change that must be sphericaly symetric !
Definition at line 145 of file param_elliptic.h.