|
LORENE
|
Base class for pure radial mappings. More...
#include <map.h>
Public Member Functions | |
| virtual | ~Map_radial () |
| Destructor. | |
| virtual void | operator= (const Map_af &)=0 |
| Assignment to an affine mapping. | |
| virtual void | sauve (FILE *) const |
| Save in a file. | |
| virtual double | val_r_jk (int l, double xi, int j, int k) const =0 |
Returns the value of the radial coordinate r for a given ![]() ![]() | |
| virtual void | val_lx_jk (double rr, int j, int k, const Param &par, int &l, double &xi) const =0 |
Computes the domain index l and the value of ![]() ![]() | |
| virtual bool | operator== (const Map &) const =0 |
| Comparison operator (egality). | |
| 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_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 (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, 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 ![]() | |
Public Attributes | |
| 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 |
![]() ![]() | |
Protected Member Functions | |
| Map_radial (const Mg3d &mgrid) | |
Constructor from a grid (protected to make Map_radial an abstract class). | |
| Map_radial (const Map_radial &mp) | |
| Copy constructor. | |
| Map_radial (const Mg3d &, FILE *) | |
Constructor from a file (see sauve(FILE* ) ). | |
| virtual void | reset_coord () |
Resets all the member Coords. | |
Base class for pure radial mappings.
()
A pure radial mapping is a mapping of the type 


Map_radial is an abstract class. Effective implementations of radial mapping are performed by the derived class Map_af and Map_et .
|
protected |
Constructor from a grid (protected to make Map_radial an abstract class).
Definition at line 89 of file map_radial.C.
References Lorene::Map().
|
protected |
Copy constructor.
Definition at line 96 of file map_radial.C.
References Lorene::Map(), and Map_radial().
|
protected |
Constructor from a file (see sauve(FILE* ) ).
Definition at line 103 of file map_radial.C.
References Lorene::Map().
|
virtual |
Destructor.
Definition at line 110 of file map_radial.C.
|
virtual |
Cmp version
Definition at line 176 of file map_radial_comp_rtp.C.
References comp_p_from_cartesian().
|
virtual |
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] ![]() |
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().
|
virtual |
Cmp version
Definition at line 65 of file map_radial_comp_rtp.C.
References comp_r_from_cartesian().
|
virtual |
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 |
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().
|
virtual |
Cmp version
Definition at line 121 of file map_radial_comp_rtp.C.
References comp_t_from_cartesian().
|
virtual |
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] ![]() |
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().
|
virtual |
Cmp version
Definition at line 68 of file map_radial_comp_xyz.C.
References comp_x_from_spherical().
|
virtual |
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 |
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().
|
virtual |
Cmp version
Definition at line 126 of file map_radial_comp_xyz.C.
References comp_y_from_spherical().
|
virtual |
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 |
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().
|
virtual |
Cmp version
Definition at line 184 of file map_radial_comp_xyz.C.
References comp_z_from_spherical().
|
virtual |
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 |
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().
|
virtual |
Decreases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the
compactified external domain (CED).
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::Valeur::mult2_xm1_zec(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and xsr.
|
virtual |
Decreases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED).
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::Valeur::mult_xm1_zec(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and xsr.
|
virtual |
Division by 
Scalar.
Definition at line 85 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Valeur::scost(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
|
virtual |
Division by r of a Scalar.
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::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sx(), and xsr.
|
virtual |
Division by r (in the compactified external domain only) of a Scalar.
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::Scalar::get_spectral_va(), Lorene::Valeur::mult_xm1_zec(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and xsr.
|
virtual |
Division by 
Scalar.
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::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::ssint(), Lorene::Valeur::sx(), and xsr.
|
virtual |
Division by 
Scalar.
Definition at line 133 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Valeur::ssint().
|
virtual |
Division by 
Scalar.
Definition at line 158 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Valeur::mult_ct(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Valeur::ssint().
|
virtual |
Increases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the
compactified external domain (CED).
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::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and xsr.
|
virtual |
Increases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the
compactified external domain (CED).
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::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and xsr.
|
virtual |
Multiplication by 
Scalar.
Definition at line 61 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Valeur::mult_ct(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
|
virtual |
Multiplication by r of a Cmp.
In the CED, there is only a decrement of dzpuis
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::Valeur::mult_x(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::va, and xsr.
|
virtual |
Multiplication by r of a Scalar, the dzpuis
of uu is not changed.
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::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and xsr.
|
virtual |
Multiplication by r (in the compactified external domain only) of a Scalar.
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::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and xsr.
|
virtual |
Multiplication by 
Scalar.
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::Valeur::mult_st(), Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and xsr.
|
virtual |
Multiplication by 
Scalar.
Definition at line 109 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Valeur::mult_st(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
|
pure virtual |
Assignment to an affine mapping.
Implemented in Lorene::Map_af, Lorene::Map_et, and Lorene::Map_log.
|
pure virtual |
Comparison operator (egality).
Implemented in Lorene::Map_af, Lorene::Map_et, and Lorene::Map_log.
References Lorene::Map().
|
virtual |
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 ![]() ![]() |
Definition at line 155 of file map_radial_poisson_cpt.C.
References Lorene::Cmp::annule(), Lorene::bvect_spher, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Valeur::d2sdx2(), Lorene::diffrel(), Lorene::Cmp::dsdr(), Lorene::Valeur::dsdx(), 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::Tenseur::get_triad(), Lorene::Valeur::lapang(), Lorene::Cmp::laplacien(), Lorene::max(), 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().
|
virtual |
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 ![]() ![]() |
Definition at line 453 of file map_radial_poisson_cpt.C.
References Lorene::Cmp::annule(), Lorene::Mtbl::annule_hard(), Lorene::Tbl::annule_hard(), Lorene::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::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Cmp::get_mp(), Lorene::Tenseur::get_mp(), Lorene::Tenseur::get_triad(), Lorene::Cmp::laplacien(), Lorene::Map_af::laplacien(), Lorene::max(), Lorene::min(), poisson_compact(), Lorene::Mtbl::set(), Lorene::Tbl::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::Valeur::ylm(), and Lorene::Valeur::ylm_i().
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. |
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::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Map(), Map_radial(), Lorene::r, Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::t, Lorene::Cmp::va, val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk().
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. |
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::Valeur::coef(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Map(), Map_radial(), Lorene::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, val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk().
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. |
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::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Map(), Map_radial(), Lorene::r, Lorene::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::t, Lorene::Cmp::va, val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk_symy().
|
virtual |
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. |
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::Valeur::coef(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Map(), Map_radial(), Lorene::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, val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk_symy().
|
protectedvirtual |
|
virtual |
Save in a file.
Reimplemented in Lorene::Map_af, Lorene::Map_et, and Lorene::Map_log.
Definition at line 116 of file map_radial.C.
|
pure 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 | [input/output] parameters to control the accuracy of the computation |
| l | [output] value of the domain index |
| xi | [output] value of ![]() |
Implemented in Lorene::Map_af, Lorene::Map_et, and Lorene::Map_log.
|
pure virtual |
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 ![]() |

Implemented in Lorene::Map_af, Lorene::Map_et, and Lorene::Map_log.
| Coord Lorene::Map_radial::d2rdtdx |
| Coord Lorene::Map_radial::d2rdx2 |
| Coord Lorene::Map_radial::drdt |
| Coord Lorene::Map_radial::dxdr |
| Coord Lorene::Map_radial::lapr_tp |
| Coord Lorene::Map_radial::sr2d2rdt2 |
| Coord Lorene::Map_radial::sr2drdt |
| Coord Lorene::Map_radial::sr2stdrdp |
| Coord Lorene::Map_radial::srdrdt |
| Coord Lorene::Map_radial::srstdrdp |
| Coord Lorene::Map_radial::sstd2rdpdx |
| Coord Lorene::Map_radial::stdrdp |
| Coord Lorene::Map_radial::xsr |