LORENE
|
Logarithmic radial mapping. More...
#include <map.h>
Public Member Functions | |
Map_log (const Mg3d &mgrille, const Tbl &r_limits, const Itbl &typevar) | |
Standard Constructor. | |
Map_log (const Map_log &) | |
Copy constructor. | |
Map_log (const Mg3d &, FILE *) | |
Constructor from a file (see sauve(FILE*) | |
virtual | ~Map_log () |
Destructor. | |
virtual const Map_af & | mp_angu (int) const |
Returns the "angular" mapping for the outside of domain l_zone . | |
double | get_alpha (int l) const |
Returns ![]() l . | |
double | get_beta (int l) const |
Returns ![]() l . | |
int | get_type (int l) const |
Returns the type of description in the domain l . | |
void | sol_elliptic (Param_elliptic ¶ms, const Scalar &so, Scalar &uu) const |
General elliptic solver. | |
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 | sol_elliptic_boundary (Param_elliptic ¶ms, const Scalar &so, Scalar &uu, const Scalar &bound, double fact_dir, double fact_neu) const |
General elliptic solver including inner boundary conditions, the bound being given as a Scalar on a mono-domain angular grid. | |
void | sol_elliptic_no_zec (Param_elliptic ¶ms, const Scalar &so, Scalar &uu, double) const |
General elliptic solver. | |
virtual void | sauve (FILE *) const |
Save in a file. | |
virtual void | operator= (const Map_af &mpa) |
Assignment to an affine mapping. | |
virtual ostream & | operator>> (ostream &) const |
Operator >> | |
virtual double | val_r (int l, double xi, double theta, double pphi) const |
Returns the value of the radial coordinate r for a given ![]() | |
virtual void | val_lx (double rr, double theta, double pphi, int &l, double &xi) const |
Computes the domain index l and the value of ![]() ![]() | |
virtual void | val_lx (double rr, double theta, double pphi, const Param &par, int &l, double &xi) const |
Computes the domain index l and the value of ![]() ![]() | |
virtual bool | operator== (const Map &) const |
Comparison operator (egality) | |
virtual double | val_r_jk (int l, double xi, int j, int k) const |
< Comparison operator | |
virtual void | val_lx_jk (double rr, int j, int k, const Param &par, int &l, double &xi) const |
Computes the domain index l and the value of ![]() ![]() | |
virtual void | dsdr (const Scalar &ci, Scalar &resu) const |
Computes ![]() Scalar . | |
virtual void | dsdxi (const Scalar &ci, Scalar &resu) const |
Computes ![]() Scalar . | |
virtual void | dsdradial (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar if the description is affine and ![]() | |
virtual void | homothetie (double) |
Sets a new radial scale. | |
virtual void | resize (int, double) |
< Not implemented | |
virtual void | adapt (const Cmp &, const Param &, int) |
< Not implemented | |
virtual void | dsdr (const Cmp &, Cmp &) const |
< Not implemented | |
virtual void | dsdxi (const Cmp &, Cmp &) const |
< Not implemented | |
virtual void | srdsdt (const Cmp &, Cmp &) const |
< Not implemented | |
virtual void | srstdsdp (const Cmp &, Cmp &) const |
< Not implemented | |
virtual void | srdsdt (const Scalar &, Scalar &) const |
< Not implemented | |
virtual void | srstdsdp (const Scalar &, Scalar &) const |
< Not implemented | |
virtual void | dsdt (const Scalar &, Scalar &) const |
< Not implemented | |
virtual void | stdsdp (const Scalar &, Scalar &) const |
< Not implemented | |
virtual void | laplacien (const Scalar &, int, Scalar &) const |
< Not implemented | |
virtual void | laplacien (const Cmp &, int, Cmp &) const |
< Not implemented | |
virtual void | lapang (const Scalar &, Scalar &) const |
< Not implemented | |
virtual void | primr (const Scalar &, Scalar &, bool) const |
< Not implemented | |
virtual Tbl * | integrale (const Cmp &) const |
< Not implemented | |
virtual void | poisson (const Cmp &, Param &, Cmp &) const |
< Not implemented | |
virtual void | poisson_tau (const Cmp &, Param &, Cmp &) const |
< Not implemented | |
virtual void | poisson_falloff (const Cmp &, Param &, Cmp &, int) const |
< Not implemented | |
virtual void | poisson_ylm (const Cmp &, Param &, Cmp &, int, double *) const |
< Not implemented | |
virtual void | poisson_regular (const Cmp &, int, int, double, Param &, Cmp &, Cmp &, Cmp &, Tenseur &, Cmp &, Cmp &) const |
< Not implemented | |
virtual void | poisson_angu (const Scalar &, Param &, Scalar &, double=0) const |
< Not implemented | |
virtual Param * | donne_para_poisson_vect (Param &, int) const |
< Not implemented | |
virtual void | poisson_frontiere (const Cmp &, const Valeur &, int, int, Cmp &, double=0., double=0.) const |
< Not implemented | |
virtual void | poisson_frontiere_double (const Cmp &, const Valeur &, const Valeur &, int, Cmp &) const |
< Not implemented | |
virtual void | poisson_interne (const Cmp &, const Valeur &, Param &, Cmp &) const |
< Not implemented | |
virtual void | poisson2d (const Cmp &, const Cmp &, Param &, Cmp &) const |
< Not implemented | |
virtual void | dalembert (Param &, Scalar &, const Scalar &, const Scalar &, const Scalar &) const |
< Not implemented | |
virtual void | reevaluate (const Map *mp_prev, int nzet, Cmp &uu) const |
Recomputes the values of a Cmp at the collocation points after a change in the mapping. | |
virtual void | reevaluate (const Map *mp_prev, int nzet, Scalar &uu) const |
Recomputes the values of a Scalar at the collocation points after a change in the mapping. | |
virtual void | reevaluate_symy (const Map *mp_prev, int nzet, Cmp &uu) const |
Recomputes the values of a Cmp at the collocation points after a change in the mapping. | |
virtual void | reevaluate_symy (const Map *mp_prev, int nzet, Scalar &uu) const |
Recomputes the values of a Scalar at the collocation points after a change in the mapping. | |
virtual void | mult_r (Scalar &uu) const |
Multiplication by r of a Scalar , the dzpuis of uu is not changed. | |
virtual void | mult_r (Cmp &ci) const |
Multiplication by r of a Cmp . | |
virtual void | mult_r_zec (Scalar &) const |
Multiplication by r (in the compactified external domain only) of a Scalar . | |
virtual void | mult_rsint (Scalar &) const |
Multiplication by ![]() Scalar . | |
virtual void | div_rsint (Scalar &) const |
Division by ![]() Scalar . | |
virtual void | div_r (Scalar &) const |
Division by r of a Scalar . | |
virtual void | div_r_zec (Scalar &) const |
Division by r (in the compactified external domain only) of a Scalar . | |
virtual void | mult_cost (Scalar &) const |
Multiplication by ![]() Scalar . | |
virtual void | div_cost (Scalar &) const |
Division by ![]() Scalar . | |
virtual void | mult_sint (Scalar &) const |
Multiplication by ![]() Scalar . | |
virtual void | div_sint (Scalar &) const |
Division by ![]() Scalar . | |
virtual void | div_tant (Scalar &) const |
Division by ![]() Scalar . | |
virtual void | comp_x_from_spherical (const Scalar &v_r, const Scalar &v_theta, const Scalar &v_phi, Scalar &v_x) const |
Computes the Cartesian x component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher . | |
virtual void | comp_x_from_spherical (const Cmp &v_r, const Cmp &v_theta, const Cmp &v_phi, Cmp &v_x) const |
Cmp version | |
virtual void | comp_y_from_spherical (const Scalar &v_r, const Scalar &v_theta, const Scalar &v_phi, Scalar &v_y) const |
Computes the Cartesian y component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher . | |
virtual void | comp_y_from_spherical (const Cmp &v_r, const Cmp &v_theta, const Cmp &v_phi, Cmp &v_y) const |
Cmp version | |
virtual void | comp_z_from_spherical (const Scalar &v_r, const Scalar &v_theta, Scalar &v_z) const |
Computes the Cartesian z component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher . | |
virtual void | comp_z_from_spherical (const Cmp &v_r, const Cmp &v_theta, Cmp &v_z) const |
Cmp version | |
virtual void | comp_r_from_cartesian (const Scalar &v_x, const Scalar &v_y, const Scalar &v_z, Scalar &v_r) const |
Computes the Spherical r component (with respect to bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart . | |
virtual void | comp_r_from_cartesian (const Cmp &v_x, const Cmp &v_y, const Cmp &v_z, Cmp &v_r) const |
Cmp version | |
virtual void | comp_t_from_cartesian (const Scalar &v_x, const Scalar &v_y, const Scalar &v_z, Scalar &v_t) const |
Computes the Spherical ![]() bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart . | |
virtual void | comp_t_from_cartesian (const Cmp &v_x, const Cmp &v_y, const Cmp &v_z, Cmp &v_t) const |
Cmp version | |
virtual void | comp_p_from_cartesian (const Scalar &v_x, const Scalar &v_y, Scalar &v_p) const |
Computes the Spherical ![]() bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart . | |
virtual void | comp_p_from_cartesian (const Cmp &v_x, const Cmp &v_y, Cmp &v_p) const |
Cmp version | |
virtual void | dec_dzpuis (Scalar &) const |
Decreases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED). | |
virtual void | dec2_dzpuis (Scalar &) const |
Decreases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED). | |
virtual void | inc_dzpuis (Scalar &) const |
Increases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED). | |
virtual void | inc2_dzpuis (Scalar &) const |
Increases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED). | |
virtual void | poisson_compact (const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const |
Resolution of the elliptic equation ![]() | |
virtual void | poisson_compact (int nzet, const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const |
Resolution of the elliptic equation ![]() | |
const Mg3d * | get_mg () const |
Gives the Mg3d on which the mapping is defined. | |
double | get_ori_x () const |
Returns the x coordinate of the origin. | |
double | get_ori_y () const |
Returns the y coordinate of the origin. | |
double | get_ori_z () const |
Returns the z coordinate of the origin. | |
double | get_rot_phi () const |
Returns the angle between the x –axis and X –axis. | |
const Base_vect_spher & | get_bvect_spher () const |
Returns the orthonormal vectorial basis ![]() ![]() | |
const Base_vect_cart & | get_bvect_cart () const |
Returns the Cartesian basis ![]() | |
const Metric_flat & | flat_met_spher () const |
Returns the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher . | |
const Metric_flat & | flat_met_cart () const |
Returns the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart . | |
const Cmp & | cmp_zero () const |
Returns the null Cmp defined on *this . | |
void | convert_absolute (double xx, double yy, double zz, double &rr, double &theta, double &pphi) const |
Determines the coordinates ![]() | |
void | set_ori (double xa0, double ya0, double za0) |
Sets a new origin. | |
void | set_rot_phi (double phi0) |
Sets a new rotation angle. | |
Public Attributes | |
Coord | dxdlnr |
Same as dxdr if the domains where the description is affine and ![]() | |
Coord | xsr |
![]() ![]() | |
Coord | dxdr |
![]() ![]() | |
Coord | drdt |
![]() ![]() | |
Coord | stdrdp |
![]() ![]() | |
Coord | srdrdt |
![]() ![]() | |
Coord | srstdrdp |
![]() ![]() | |
Coord | sr2drdt |
![]() ![]() | |
Coord | sr2stdrdp |
![]() ![]() | |
Coord | d2rdx2 |
![]() ![]() | |
Coord | lapr_tp |
![]() ![]() | |
Coord | d2rdtdx |
![]() ![]() | |
Coord | sstd2rdpdx |
![]() ![]() | |
Coord | sr2d2rdt2 |
![]() ![]() | |
Coord | r |
r coordinate centered on the grid | |
Coord | tet |
![]() | |
Coord | phi |
![]() | |
Coord | sint |
![]() | |
Coord | cost |
![]() | |
Coord | sinp |
![]() | |
Coord | cosp |
![]() | |
Coord | x |
x coordinate centered on the grid | |
Coord | y |
y coordinate centered on the grid | |
Coord | z |
z coordinate centered on the grid | |
Coord | xa |
Absolute x coordinate. | |
Coord | ya |
Absolute y coordinate. | |
Coord | za |
Absolute z coordinate. | |
Protected Member Functions | |
virtual void | reset_coord () |
Resets all the member Coords . | |
Protected Attributes | |
const Mg3d * | mg |
Pointer on the multi-grid Mgd3 on which this is defined | |
double | ori_x |
Absolute coordinate x of the origin. | |
double | ori_y |
Absolute coordinate y of the origin. | |
double | ori_z |
Absolute coordinate z of the origin. | |
double | rot_phi |
Angle between the x –axis and X –axis. | |
Base_vect_spher | bvect_spher |
Orthonormal vectorial basis ![]() ![]() | |
Base_vect_cart | bvect_cart |
Cartesian basis ![]() | |
Metric_flat * | p_flat_met_spher |
Pointer onto the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher . | |
Metric_flat * | p_flat_met_cart |
Pointer onto the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart . | |
Cmp * | p_cmp_zero |
The null Cmp. | |
Map_af * | p_mp_angu |
Pointer on the "angular" mapping. | |
Private Member Functions | |
void | set_coord () |
Private Attributes | |
Tbl | alpha |
Array (size: mg->nzone ) of the values of ![]() | |
Tbl | beta |
Array (size: mg->nzone ) of the values of ![]() | |
Itbl | type_var |
Array (size: mg->nzone ) of the type of variable in each domain. | |
Logarithmic radial mapping.
()
This mapping is a variation of the affine one.
In each domain the description can be either affine (cf. Map_af documentation) or logarithmic. In that case (implemented only in the shells) we have
Standard Constructor.
mgrille | [input] Multi-domain grid on which the mapping is defined |
r_limits | [input] Tbl (size: number of domains + 1) of the value of r at the boundaries of the various domains :
|
type_var | [input] Array (size: number of domains) defining the type f mapping in each domain. |
Definition at line 67 of file map_log.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::log(), Lorene::Map::mg, Lorene::Tbl::set(), Lorene::Tbl::set_etat_qcq(), and type_var.
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 176 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_p_from_cartesian().
|
virtualinherited |
Computes the Spherical bvect_spher
) of a vector given by its cartesian components with respect to bvect_cart
.
v_x | [input] x-component of the vector |
v_y | [input] y-component of the vector |
v_p | [output] ![]() |
Implements Lorene::Map.
Definition at line 183 of file map_radial_comp_rtp.C.
References Lorene::Scalar::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_sp(), and Lorene::Scalar::set_dzpuis().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 65 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_r_from_cartesian().
|
virtualinherited |
Computes the Spherical r component (with respect to bvect_spher
) of a vector given by its cartesian components with respect to bvect_cart
.
v_x | [input] x-component of the vector |
v_y | [input] y-component of the vector |
v_z | [input] z-component of the vector |
v_r | [output] r -component of the vector |
Implements Lorene::Map.
Definition at line 72 of file map_radial_comp_rtp.C.
References Lorene::Scalar::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 121 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_t_from_cartesian().
|
virtualinherited |
Computes the Spherical bvect_spher
) of a vector given by its cartesian components with respect to bvect_cart
.
v_x | [input] x-component of the vector |
v_y | [input] y-component of the vector |
v_z | [input] z-component of the vector |
v_t | [output] ![]() |
Implements Lorene::Map.
Definition at line 128 of file map_radial_comp_rtp.C.
References Lorene::Scalar::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 68 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_x_from_spherical().
|
virtualinherited |
Computes the Cartesian x component (with respect to bvect_cart
) of a vector given by its spherical components with respect to bvect_spher
.
v_r | [input] r -component of the vector |
v_theta | [input] ![]() |
v_phi | [input] ![]() |
v_x | [output] x-component of the vector |
Implements Lorene::Map.
Definition at line 76 of file map_radial_comp_xyz.C.
References Lorene::Scalar::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 126 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_y_from_spherical().
|
virtualinherited |
Computes the Cartesian y component (with respect to bvect_cart
) of a vector given by its spherical components with respect to bvect_spher
.
v_r | [input] r -component of the vector |
v_theta | [input] ![]() |
v_phi | [input] ![]() |
v_y | [output] y-component of the vector |
Implements Lorene::Map.
Definition at line 135 of file map_radial_comp_xyz.C.
References Lorene::Scalar::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 184 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_z_from_spherical().
|
virtualinherited |
Computes the Cartesian z component (with respect to bvect_cart
) of a vector given by its spherical components with respect to bvect_spher
.
v_r | [input] r -component of the vector |
v_theta | [input] ![]() |
v_z | [output] z-component of the vector |
Implements Lorene::Map.
Definition at line 192 of file map_radial_comp_xyz.C.
References Lorene::Scalar::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().
|
inherited |
Determines the coordinates
xx | [input] value of the coordinate x (absolute frame) |
yy | [input] value of the coordinate y (absolute frame) |
zz | [input] value of the coordinate z (absolute frame) |
rr | [output] value of r |
theta | [output] value of ![]() |
pphi | [output] value of ![]() |
Definition at line 302 of file map.C.
References Lorene::Map::ori_x, Lorene::Map::ori_y, Lorene::Map::ori_z, Lorene::Map::rot_phi, and Lorene::sqrt().
Decreases by 2 the value of dzpuis
of a Scalar
and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 748 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult2_xm1_zec(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Map_radial::xsr.
Decreases by 1 the value of dzpuis
of a Scalar
and changes accordingly its values in the compactified external domain (CED).
Implements Lorene::Map.
Definition at line 643 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Map_radial::xsr.
Division by Scalar
.
Implements Lorene::Map.
Definition at line 85 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Valeur::scost(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
Division by r of a Scalar
.
Implements Lorene::Map.
Definition at line 514 of file map_radial_r_manip.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sx(), and Lorene::Map_radial::xsr.
Division by r (in the compactified external domain only) of a Scalar
.
Implements Lorene::Map.
Definition at line 585 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Map_radial::xsr.
Division by Scalar
.
Implements Lorene::Map.
Definition at line 434 of file map_radial_r_manip.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::ssint(), Lorene::Valeur::sx(), and Lorene::Map_radial::xsr.
Division by Scalar
.
Implements Lorene::Map.
Definition at line 133 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Valeur::ssint().
Division by Scalar
.
Implements Lorene::Map.
Definition at line 158 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Valeur::mult_ct(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Valeur::ssint().
Computes Scalar
.
Note that in the compactified external domain (CED), it computes
ci | [input] field to consider |
resu | [output] derivative of ci |
Implements Lorene::Map.
Definition at line 108 of file map_log_deriv.C.
References Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::coef(), Lorene::Valeur::dsdx(), Lorene::Map_radial::dxdr, Lorene::Valeur::get_base(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_base(), and Lorene::Scalar::set_spectral_va().
Computes Scalar
if the description is affine and
Note that in the compactified external domain (CED), the dzpuis
flag of the output is 2 if the input has dzpuis
= 0, and is increased by 1 in other cases.
uu | [input] field to consider |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 162 of file map_log_deriv.C.
References Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::coef(), Lorene::Valeur::dsdx(), dxdlnr, Lorene::Valeur::get_base(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_base(), and Lorene::Scalar::set_spectral_va().
Computes Scalar
.
Note that in the compactified external domain (CED), it computes
ci | [input] field to consider |
resu | [output] derivative of ci |
Implements Lorene::Map.
Definition at line 54 of file map_log_deriv.C.
References Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::coef(), Lorene::Valeur::dsdx(), Lorene::Valeur::get_base(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_base(), and Lorene::Scalar::set_spectral_va().
|
inherited |
Returns the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart
.
Definition at line 331 of file map.C.
References Lorene::Map::bvect_cart, and Lorene::Map::p_flat_met_cart.
|
inherited |
Returns the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher
.
Definition at line 321 of file map.C.
References Lorene::Map::bvect_spher, and Lorene::Map::p_flat_met_spher.
|
inlineinherited |
Returns the Cartesian basis
the Cartesian coordinates related to
Definition at line 791 of file map.h.
References Lorene::Map::bvect_cart.
|
inlineinherited |
Returns the orthonormal vectorial basis
Definition at line 783 of file map.h.
References Lorene::Map::bvect_spher.
Gives the Mg3d
on which the mapping is defined.
Definition at line 765 of file map.h.
References Lorene::Map::mg.
|
inlineinherited |
Returns the x coordinate of the origin.
Definition at line 768 of file map.h.
References Lorene::Map::ori_x.
|
inlineinherited |
Returns the y coordinate of the origin.
Definition at line 770 of file map.h.
References Lorene::Map::ori_y.
|
inlineinherited |
Returns the z coordinate of the origin.
Definition at line 772 of file map.h.
References Lorene::Map::ori_z.
|
inlineinherited |
Returns the angle between the x –axis and X –axis.
Definition at line 775 of file map.h.
References Lorene::Map::rot_phi.
Sets a new radial scale.
lambda | [input] factor by which the value of r is to be multiplied |
Implements Lorene::Map.
Definition at line 93 of file map_log_pas_fait.C.
Increases by 2 the value of dzpuis
of a Scalar
and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 799 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Increases by 1 the value of dzpuis
of a Scalar
and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 696 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Returns the "angular" mapping for the outside of domain l_zone
.
Valid only for the class Map_af
.
Implements Lorene::Map.
Definition at line 196 of file map_log_pas_fait.C.
References Lorene::Map::p_mp_angu.
Multiplication by Scalar
.
Implements Lorene::Map.
Definition at line 61 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Valeur::mult_ct(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
Multiplication by r of a Cmp
.
In the CED, there is only a decrement of dzpuis
Implements Lorene::Map.
Definition at line 219 of file map_radial_r_manip.C.
References Lorene::Cmp::annule(), Lorene::Valeur::base, Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_x(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::va, and Lorene::Map_radial::xsr.
Multiplication by r of a Scalar
, the dzpuis
of uu
is not changed.
Implements Lorene::Map.
Definition at line 144 of file map_radial_r_manip.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Multiplication by r (in the compactified external domain only) of a Scalar
.
Implements Lorene::Map.
Definition at line 296 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Multiplication by Scalar
.
Implements Lorene::Map.
Definition at line 355 of file map_radial_r_manip.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_st(), Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Multiplication by Scalar
.
Implements Lorene::Map.
Definition at line 109 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Valeur::mult_st(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
Assignment to an affine mapping.
Implements Lorene::Map_radial.
Definition at line 230 of file map_log.C.
References alpha, beta, Lorene::Map_af::get_alpha(), Lorene::Map_af::get_beta(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_x(), Lorene::Map::get_ori_y(), Lorene::Map::get_ori_z(), Lorene::Map::get_rot_phi(), Lorene::Map::mg, Lorene::Map_radial::reset_coord(), Lorene::Tbl::set(), Lorene::Map::set_ori(), Lorene::Map::set_rot_phi(), and type_var.
Comparison operator (egality)
Implements Lorene::Map_radial.
Definition at line 172 of file map_log.C.
References alpha, beta, Lorene::Map::bvect_cart, Lorene::Map::bvect_spher, Lorene::diffrelmax(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Map::get_ori_x(), Lorene::Map::get_ori_y(), Lorene::Map::get_ori_z(), Lorene::Map::mg, Lorene::Map::ori_x, Lorene::Map::ori_y, Lorene::Map::ori_z, and type_var.
Operator >>
Implements Lorene::Map.
Definition at line 204 of file map_log.C.
References alpha, beta, Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, Lorene::Map::r, and type_var.
|
virtualinherited |
Resolution of the elliptic equation
source | [input] source ![]() |
aa | [input] factor a in the above equation |
bb | [input] vector b in the above equation |
par | [input/output] parameters of the iterative method of resolution : \ par.get_int(0) : [input] maximum number of iterations \ par.get_double(0) : [input] required precision: the iterative method is stopped as soon as the relative difference between ![]() ![]() par.get_double(0) \ par.get_double(1) : [input] relaxation parameter ![]() par.get_int_mod(0) : [output] number of iterations actually used to get the solution. |
psi | [input/output]: input : previously computed value of ![]() psi must be set to zero); output: solution ![]() ![]() |
Implements Lorene::Map.
Definition at line 155 of file map_radial_poisson_cpt.C.
References Lorene::Cmp::annule(), Lorene::Map::bvect_spher, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Valeur::d2sdx2(), Lorene::diffrel(), Lorene::Cmp::dsdr(), Lorene::Valeur::dsdx(), Lorene::Map_radial::dxdr, Lorene::Param::get_double(), Lorene::Cmp::get_etat(), Lorene::Tenseur::get_etat(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Cmp::get_mp(), Lorene::Tenseur::get_mp(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nzone(), Lorene::Tenseur::get_triad(), Lorene::Valeur::lapang(), Lorene::Cmp::laplacien(), Lorene::max(), Lorene::Map::mg, Lorene::min(), Lorene::Valeur::mult_x(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Cmp::srdsdt(), Lorene::Cmp::srstdsdp(), Lorene::Valeur::std_base_scal(), Lorene::Valeur::sx(), Lorene::Cmp::va, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
|
virtualinherited |
Resolution of the elliptic equation
nzet | [input] number of domains covering the stellar interior |
source | [input] source ![]() |
aa | [input] factor a in the above equation |
bb | [input] vector b in the above equation |
par | [input/output] possible parameters to control the resolution of the equation. See the actual implementation in the derived class of Map for documentation. |
psi | [input/output] solution ![]() ![]() |
Implements Lorene::Map.
Definition at line 453 of file map_radial_poisson_cpt.C.
References Lorene::Cmp::annule(), Lorene::Mtbl::annule_hard(), Lorene::Tbl::annule_hard(), Lorene::Map::bvect_spher, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::diffrel(), Lorene::Cmp::dsdr(), Lorene::Map_af::dsdr(), Lorene::Param::get_double(), Lorene::Cmp::get_etat(), Lorene::Tenseur::get_etat(), Lorene::Mg3d::get_grille3d(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Cmp::get_mp(), Lorene::Tenseur::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Tenseur::get_triad(), Lorene::Map_af::laplacien(), Lorene::Cmp::laplacien(), Lorene::max(), Lorene::Map::mg, Lorene::min(), Lorene::Map_radial::poisson_compact(), Lorene::Tbl::set(), Lorene::Mtbl::set(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Cmp::srdsdt(), Lorene::Cmp::srstdsdp(), Lorene::Valeur::std_base_scal(), Lorene::Tbl::t, Lorene::Cmp::va, Lorene::Grille3d::x, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
|
virtualinherited |
Recomputes the values of a Cmp
at the collocation points after a change in the mapping.
mp_prev | [input] Previous value of the mapping. |
nzet | [input] Number of domains where the computation must be done: the computation is performed for the domains of index ![]() uu is set to zero in the other domains. |
uu | [input/output] input : Cmp previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Cmp at the grid points defined by *this . |
Implements Lorene::Map.
Definition at line 58 of file map_radial_reevaluate.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Cmp::annule(), Lorene::Coord::c, Lorene::Coord::fait(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, Lorene::Map::r, Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::t, Lorene::Cmp::va, Lorene::Map_radial::val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk().
|
virtualinherited |
Recomputes the values of a Scalar
at the collocation points after a change in the mapping.
mp_prev | [input] Previous value of the mapping. |
nzet | [input] Number of domains where the computation must be done: the computation is performed for the domains of index ![]() uu is set to zero in the other domains. |
uu | [input/output] input : Scalar previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Scalar at the grid points defined by *this . |
Implements Lorene::Map.
Definition at line 173 of file map_radial_reevaluate.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Scalar::annule(), Lorene::Coord::c, Lorene::Valeur::c, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Coord::fait(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, Lorene::Map::r, Lorene::Valeur::set_etat_c_qcq(), Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Mtbl::t, Lorene::Map_radial::val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk().
|
virtualinherited |
Recomputes the values of a Cmp
at the collocation points after a change in the mapping.
Case where the Cmp
is symmetric with respect to the plane y=0.
mp_prev | [input] Previous value of the mapping. |
nzet | [input] Number of domains where the computation must be done: the computation is performed for the domains of index ![]() uu is set to zero in the other domains. |
uu | [input/output] input : Cmp previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Cmp at the grid points defined by *this . |
Implements Lorene::Map.
Definition at line 59 of file map_radial_reeval_symy.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Cmp::annule(), Lorene::Coord::c, Lorene::Coord::fait(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_p(), Lorene::Map::mg, Lorene::Map::r, Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::t, Lorene::Cmp::va, Lorene::Map_radial::val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk_symy().
|
virtualinherited |
Recomputes the values of a Scalar
at the collocation points after a change in the mapping.
Case where the Scalar
is symmetric with respect to the plane y=0.
mp_prev | [input] Previous value of the mapping. |
nzet | [input] Number of domains where the computation must be done: the computation is performed for the domains of index ![]() uu is set to zero in the other domains. |
uu | [input/output] input : Scalar previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Scalar at the grid points defined by *this . |
Implements Lorene::Map.
Definition at line 193 of file map_radial_reeval_symy.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Scalar::annule(), Lorene::Coord::c, Lorene::Valeur::c, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Coord::fait(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_p(), Lorene::Map::mg, Lorene::Map::r, Lorene::Valeur::set_etat_c_qcq(), Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Mtbl::t, Lorene::Map_radial::val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk_symy().
|
protectedvirtualinherited |
Resets all the member Coords
.
Reimplemented from Lorene::Map.
Reimplemented in Lorene::Map_et.
Definition at line 126 of file map_radial.C.
References Lorene::Map_radial::d2rdtdx, Lorene::Map_radial::d2rdx2, Lorene::Coord::del_t(), Lorene::Map_radial::drdt, Lorene::Map_radial::dxdr, Lorene::Map_radial::lapr_tp, Lorene::Map::reset_coord(), Lorene::Map_radial::sr2d2rdt2, Lorene::Map_radial::sr2drdt, Lorene::Map_radial::sr2stdrdp, Lorene::Map_radial::srdrdt, Lorene::Map_radial::srstdrdp, Lorene::Map_radial::sstd2rdpdx, Lorene::Map_radial::stdrdp, and Lorene::Map_radial::xsr.
Save in a file.
Reimplemented from Lorene::Map_radial.
Definition at line 163 of file map_log.C.
References alpha, beta, Lorene::Itbl::sauve(), Lorene::Map_radial::sauve(), Lorene::Tbl::sauve(), and type_var.
Sets a new origin.
Definition at line 253 of file map.C.
References Lorene::Map::bvect_spher, Lorene::Map::ori_x, Lorene::Map::ori_y, Lorene::Map::ori_z, Lorene::Map::reset_coord(), and Lorene::Base_vect_spher::set_ori().
Sets a new rotation angle.
Definition at line 263 of file map.C.
References Lorene::Map::bvect_cart, Lorene::Map::bvect_spher, Lorene::Map::reset_coord(), Lorene::Map::rot_phi, Lorene::Base_vect_cart::set_rot_phi(), and Lorene::Base_vect_spher::set_rot_phi().
void Lorene::Map_log::sol_elliptic | ( | Param_elliptic & | params, |
const Scalar & | so, | ||
Scalar & | uu | ||
) | const |
General elliptic solver.
The field is zero at infinity.
params | [input] : the operators and variables to be uses. |
so | [input] : the source. |
uu | [output] : the solution. |
Definition at line 75 of file map_log_elliptic.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
void Lorene::Map_log::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.
The field is zero at infinity.
params | [input] : the operators and variables to be uses. |
so | [input] : the source. |
uu | [output] : the solution. |
bound | [input] : the boundary condition |
fact_dir | : 1 Dirchlet condition, 0 Neumann condition |
fact_neu | : 0 Dirchlet condition, 1 Neumann condition |
Definition at line 132 of file map_log_elliptic.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
void Lorene::Map_log::sol_elliptic_boundary | ( | Param_elliptic & | params, |
const Scalar & | so, | ||
Scalar & | uu, | ||
const Scalar & | bound, | ||
double | fact_dir, | ||
double | fact_neu | ||
) | const |
General elliptic solver including inner boundary conditions, the bound being given as a Scalar on a mono-domain angular grid.
Definition at line 191 of file map_log_elliptic.C.
References Lorene::Mtbl_cf::annule_hard(), Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Mg3d::get_angu_1dom(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Mtbl_cf::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
void Lorene::Map_log::sol_elliptic_no_zec | ( | Param_elliptic & | params, |
const Scalar & | so, | ||
Scalar & | uu, | ||
double | val | ||
) | const |
General elliptic solver.
The equation is not solved in the compactified domain.
params | [input] : the operators and variables to be uses. |
so | [input] : the source. |
uu | [output] : the solution. |
val | [input] : value at the last shell. |
Definition at line 275 of file map_log_elliptic.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
|
virtual |
Computes the domain index l and the value of
rr | [input] value of r |
theta | [input] value of ![]() |
pphi | [input] value of ![]() |
par | [] unused by the Map_af version |
l | [output] value of the domain index |
xi | [output] value of ![]() |
Implements Lorene::Map.
Definition at line 204 of file map_log_radius.C.
References val_lx().
|
virtual |
Computes the domain index l and the value of
rr | [input] value of r |
theta | [input] value of ![]() |
pphi | [input] value of ![]() |
l | [output] value of the domain index |
xi | [output] value of ![]() |
Implements Lorene::Map.
Definition at line 118 of file map_log_radius.C.
References alpha, beta, Lorene::exp(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::log(), Lorene::Map::mg, and type_var.
|
virtual |
Computes the domain index l and the value of
rr | [input] value of r |
j | [input] index of the collocation point in ![]() |
k | [input] index of the collocation point in ![]() |
par | [] unused by the Map_af version |
l | [output] value of the domain index |
xi | [output] value of ![]() |
Implements Lorene::Map_radial.
Definition at line 227 of file map_log_radius.C.
References val_lx().
Returns the value of the radial coordinate r for a given
l | [input] index of the domain |
xi | [input] value of ![]() |
theta | [input] value of ![]() |
pphi | [input] value of ![]() |
Implements Lorene::Map.
Definition at line 61 of file map_log_radius.C.
References alpha, beta, Lorene::exp(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and type_var.
< Comparison operator
Returns the value of the radial coordinate r for a given
l | [input] index of the domain |
xi | [input] value of ![]() |
j | [input] index of the collocation point in ![]() |
k | [input] index of the collocation point in ![]() |
Implements Lorene::Map_radial.
Definition at line 217 of file map_log_radius.C.
References val_r().
Definition at line 467 of file map_log_fait.C.
Definition at line 399 of file map_log_fait.C.
Definition at line 923 of file map_log_fait.C.
Definition at line 822 of file map_log_fait.C.
Definition at line 714 of file map_log_fait.C.
Definition at line 977 of file map_log_fait.C.
Definition at line 616 of file map_log_fait.C.
Definition at line 905 of file map_log_fait.C.
Definition at line 193 of file map_log_fait.C.
< Not implemented
Definition at line 57 of file map_log_fait.C.
Definition at line 433 of file map_log_fait.C.
Definition at line 365 of file map_log_fait.C.
Definition at line 959 of file map_log_fait.C.
Definition at line 786 of file map_log_fait.C.
Definition at line 804 of file map_log_fait.C.
Definition at line 750 of file map_log_fait.C.
Definition at line 768 of file map_log_fait.C.
Definition at line 941 of file map_log_fait.C.
Definition at line 732 of file map_log_fait.C.
Definition at line 155 of file map_log_fait.C.
Definition at line 231 of file map_log_fait.C.
Definition at line 285 of file map_log_fait.C.
Definition at line 507 of file map_log_fait.C.
Definition at line 249 of file map_log_fait.C.
Definition at line 315 of file map_log_fait.C.
Definition at line 267 of file map_log_fait.C.
Definition at line 345 of file map_log_fait.C.
|
private |
|
private |
|
protectedinherited |
|
protectedinherited |
|
inherited |
|
inherited |
|
inherited |
Coord Lorene::Map_log::dxdlnr |
|
inherited |
|
inherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
inherited |
|
inherited |
|
protectedinherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
private |
|
inherited |
|
inherited |
|
inherited |
|
inherited |