LORENE
|
Radial mapping of rather general form. More...
#include <map.h>
Public Member Functions | |
Map_et (const Mg3d &mgrille, const double *r_limits) | |
Standard Constructor. More... | |
Map_et (const Mg3d &mgrille, const double *r_limits, const Tbl &tab) | |
Constructor using the equation of the surface of the star. More... | |
Map_et (const Map_et &) | |
Copy constructor. More... | |
Map_et (const Mg3d &, FILE *) | |
Constructor from a file (see sauve(FILE*) ) More... | |
virtual | ~Map_et () |
Destructor. More... | |
virtual void | operator= (const Map_et &mp) |
Assignment to another Map_et More... | |
virtual void | operator= (const Map_af &mpa) |
Assignment to an affine mapping. More... | |
void | set_ff (const Valeur &) |
Assigns a given value to the function ![]() | |
void | set_gg (const Valeur &) |
Assigns a given value to the function ![]() | |
virtual const Map_af & | mp_angu (int) const |
Returns the "angular" mapping for the outside of domain l_zone . More... | |
const double * | get_alpha () const |
Returns a pointer on the array alpha (values of ![]() | |
const double * | get_beta () const |
Returns a pointer on the array beta (values of ![]() | |
const Valeur & | get_ff () const |
Returns a (constant) reference to the function ![]() | |
const Valeur & | get_gg () const |
Returns a (constant) reference to the function ![]() | |
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) More... | |
virtual double | val_r_jk (int l, double xi, int j, int k) const |
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 |
Computes the domain index l and the value of ![]() ![]() | |
virtual void | sauve (FILE *) const |
Save in a file. More... | |
virtual void | homothetie (double lambda) |
Sets a new radial scale. More... | |
virtual void | resize (int l, double lambda) |
Rescales the outer boundary of one domain. More... | |
void | resize_extr (double lambda) |
Rescales the outer boundary of the outermost domain in the case of non-compactified external domain. More... | |
void | set_alpha (double alpha0, int l) |
Modifies the value of ![]() | |
void | set_beta (double beta0, int l) |
Modifies the value of ![]() | |
virtual void | adapt (const Cmp &ent, const Param &par, int nbr_filtre=0) |
Adaptation of the mapping to a given scalar field. More... | |
virtual void | dsdxi (const Cmp &ci, Cmp &resu) const |
Computes ![]() Cmp . More... | |
virtual void | dsdr (const Cmp &ci, Cmp &resu) const |
Computes ![]() Cmp . More... | |
virtual void | srdsdt (const Cmp &ci, Cmp &resu) const |
Computes ![]() Cmp . More... | |
virtual void | srstdsdp (const Cmp &ci, Cmp &resu) const |
Computes ![]() Cmp . More... | |
virtual void | dsdxi (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar . More... | |
virtual void | dsdr (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar . More... | |
virtual void | dsdradial (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar if the description is affine and ![]() | |
virtual void | srdsdt (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar . More... | |
virtual void | srstdsdp (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar . More... | |
virtual void | dsdt (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar . More... | |
virtual void | stdsdp (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar . More... | |
virtual void | laplacien (const Scalar &uu, int zec_mult_r, Scalar &lap) const |
Computes the Laplacian of a scalar field. More... | |
virtual void | laplacien (const Cmp &uu, int zec_mult_r, Cmp &lap) const |
Computes the Laplacian of a scalar field (Cmp version). More... | |
virtual void | lapang (const Scalar &uu, Scalar &lap) const |
Computes the angular Laplacian of a scalar field. More... | |
virtual void | primr (const Scalar &uu, Scalar &resu, bool null_infty) const |
Computes the radial primitive which vanishes for ![]() | |
virtual Tbl * | integrale (const Cmp &) const |
Computes the integral over all space of a Cmp . More... | |
virtual void | poisson (const Cmp &source, Param &par, Cmp &uu) const |
Computes the solution of a scalar Poisson equation. More... | |
virtual void | poisson_tau (const Cmp &source, Param &par, Cmp &uu) const |
Computes the solution of a scalar Poisson equation with a Tau method. More... | |
virtual void | poisson_falloff (const Cmp &source, Param &par, Cmp &uu, int k_falloff) const |
virtual void | poisson_ylm (const Cmp &source, Param &par, Cmp &uu, int nylm, double *intvec) const |
virtual void | poisson_regular (const Cmp &source, int k_div, int nzet, double unsgam1, Param &par, Cmp &uu, Cmp &uu_regu, Cmp &uu_div, Tenseur &duu_div, Cmp &source_regu, Cmp &source_div) const |
Computes the solution of a scalar Poisson equation. More... | |
virtual void | poisson_angu (const Scalar &source, Param &par, Scalar &uu, double lambda=0) const |
Computes the solution of the generalized angular Poisson equation. More... | |
virtual Param * | donne_para_poisson_vect (Param ¶, int i) const |
Internal function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara . More... | |
virtual void | poisson_frontiere (const Cmp &, const Valeur &, int, int, Cmp &, double=0., double=0.) const |
Not yet implemented. More... | |
virtual void | poisson_frontiere_double (const Cmp &source, const Valeur &lim_func, const Valeur &lim_der, int num_zone, Cmp &pot) const |
virtual void | poisson_interne (const Cmp &source, const Valeur &limite, Param &par, Cmp &pot) const |
Computes the solution of a Poisson equation in the shell . More... | |
virtual void | poisson2d (const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu) const |
Computes the solution of a 2-D Poisson equation. More... | |
virtual void | dalembert (Param &par, Scalar &fJp1, const Scalar &fJ, const Scalar &fJm1, const Scalar &source) const |
Not yet implemented. More... | |
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. More... | |
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. More... | |
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. More... | |
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. More... | |
virtual void | mult_r (Scalar &uu) const |
Multiplication by r of a Scalar , the dzpuis of uu is not changed. More... | |
virtual void | mult_r (Cmp &ci) const |
Multiplication by r of a Cmp . More... | |
virtual void | mult_r_zec (Scalar &) const |
Multiplication by r (in the compactified external domain only) of a Scalar . More... | |
virtual void | mult_rsint (Scalar &) const |
Multiplication by ![]() Scalar . More... | |
virtual void | div_rsint (Scalar &) const |
Division by ![]() Scalar . More... | |
virtual void | div_r (Scalar &) const |
Division by r of a Scalar . More... | |
virtual void | div_r_zec (Scalar &) const |
Division by r (in the compactified external domain only) of a Scalar . More... | |
virtual void | mult_cost (Scalar &) const |
Multiplication by ![]() Scalar . More... | |
virtual void | div_cost (Scalar &) const |
Division by ![]() Scalar . More... | |
virtual void | mult_sint (Scalar &) const |
Multiplication by ![]() Scalar . More... | |
virtual void | div_sint (Scalar &) const |
Division by ![]() Scalar . More... | |
virtual void | div_tant (Scalar &) const |
Division by ![]() Scalar . More... | |
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 . More... | |
virtual void | comp_x_from_spherical (const Cmp &v_r, const Cmp &v_theta, const Cmp &v_phi, Cmp &v_x) const |
Cmp version More... | |
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 . More... | |
virtual void | comp_y_from_spherical (const Cmp &v_r, const Cmp &v_theta, const Cmp &v_phi, Cmp &v_y) const |
Cmp version More... | |
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 . More... | |
virtual void | comp_z_from_spherical (const Cmp &v_r, const Cmp &v_theta, Cmp &v_z) const |
Cmp version More... | |
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 . More... | |
virtual void | comp_r_from_cartesian (const Cmp &v_x, const Cmp &v_y, const Cmp &v_z, Cmp &v_r) const |
Cmp version More... | |
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 . More... | |
virtual void | comp_t_from_cartesian (const Cmp &v_x, const Cmp &v_y, const Cmp &v_z, Cmp &v_t) const |
Cmp version More... | |
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 . More... | |
virtual void | comp_p_from_cartesian (const Cmp &v_x, const Cmp &v_y, Cmp &v_p) const |
Cmp version More... | |
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). More... | |
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). More... | |
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). More... | |
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). More... | |
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. More... | |
double | get_ori_x () const |
Returns the x coordinate of the origin. More... | |
double | get_ori_y () const |
Returns the y coordinate of the origin. More... | |
double | get_ori_z () const |
Returns the z coordinate of the origin. More... | |
double | get_rot_phi () const |
Returns the angle between the x –axis and X –axis. More... | |
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 . More... | |
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 . More... | |
const Cmp & | cmp_zero () const |
Returns the null Cmp defined on *this . More... | |
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. More... | |
void | set_rot_phi (double phi0) |
Sets a new rotation angle. More... | |
Public Attributes | |
Coord | rsxdxdr |
![]() ![]() ![]() | |
Coord | rsx2drdx |
![]() ![]() | |
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 More... | |
Coord | tet |
![]() | |
Coord | phi |
![]() | |
Coord | sint |
![]() | |
Coord | cost |
![]() | |
Coord | sinp |
![]() | |
Coord | cosp |
![]() | |
Coord | x |
x coordinate centered on the grid More... | |
Coord | y |
y coordinate centered on the grid More... | |
Coord | z |
z coordinate centered on the grid More... | |
Coord | xa |
Absolute x coordinate. More... | |
Coord | ya |
Absolute y coordinate. More... | |
Coord | za |
Absolute z coordinate. More... | |
Protected Member Functions | |
virtual void | reset_coord () |
Resets all the member Coords . More... | |
Protected Attributes | |
const Mg3d * | mg |
Pointer on the multi-grid Mgd3 on which this is defined More... | |
double | ori_x |
Absolute coordinate x of the origin. More... | |
double | ori_y |
Absolute coordinate y of the origin. More... | |
double | ori_z |
Absolute coordinate z of the origin. More... | |
double | rot_phi |
Angle between the x –axis and X –axis. More... | |
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 . More... | |
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 . More... | |
Cmp * | p_cmp_zero |
The null Cmp. More... | |
Map_af * | p_mp_angu |
Pointer on the "angular" mapping. More... | |
Private Member Functions | |
void | set_coord () |
Assignement of the building functions to the member Coords . More... | |
void | fait_poly () |
Construction of the polynomials ![]() ![]() | |
virtual ostream & | operator>> (ostream &) const |
Operator >> More... | |
Private Attributes | |
double * | alpha |
Array (size: mg->nzone ) of the values of ![]() | |
double * | beta |
Array (size: mg->nzone ) of the values of ![]() | |
Tbl ** | aa |
Array (size: mg->nzone ) of Tbl which stores the values of ![]() | |
Tbl ** | daa |
Array (size: mg->nzone ) of Tbl which stores the values of ![]() | |
Tbl ** | ddaa |
Array (size: mg->nzone ) of Tbl which stores the values of ![]() | |
Tbl | aasx |
Values at the nr collocation points of ![]() | |
Tbl | aasx2 |
Values at the nr collocation points of ![]() | |
Tbl | zaasx |
Values at the nr collocation points of ![]() | |
Tbl | zaasx2 |
Values at the nr collocation points of ![]() | |
Tbl ** | bb |
Array (size: mg->nzone ) of Tbl which stores the values of ![]() | |
Tbl ** | dbb |
Array (size: mg->nzone ) of Tbl which stores the values of ![]() | |
Tbl ** | ddbb |
Array (size: mg->nzone ) of Tbl which stores the values of ![]() | |
Tbl | bbsx |
Values at the nr collocation points of ![]() | |
Tbl | bbsx2 |
Values at the nr collocation points of ![]() | |
Valeur | ff |
Values of the function ![]() nt*np angular collocation points in each domain. More... | |
Valeur | gg |
Values of the function ![]() nt*np angular collocation points in each domain. More... | |
Friends | |
Mtbl * | map_et_fait_r (const Map *) |
Mtbl * | map_et_fait_tet (const Map *) |
Mtbl * | map_et_fait_phi (const Map *) |
Mtbl * | map_et_fait_sint (const Map *) |
Mtbl * | map_et_fait_cost (const Map *) |
Mtbl * | map_et_fait_sinp (const Map *) |
Mtbl * | map_et_fait_cosp (const Map *) |
Mtbl * | map_et_fait_x (const Map *) |
Mtbl * | map_et_fait_y (const Map *) |
Mtbl * | map_et_fait_z (const Map *) |
Mtbl * | map_et_fait_xa (const Map *) |
Mtbl * | map_et_fait_ya (const Map *) |
Mtbl * | map_et_fait_za (const Map *) |
Mtbl * | map_et_fait_xsr (const Map *) |
Mtbl * | map_et_fait_dxdr (const Map *) |
Mtbl * | map_et_fait_drdt (const Map *) |
Mtbl * | map_et_fait_stdrdp (const Map *) |
Mtbl * | map_et_fait_srdrdt (const Map *) |
Mtbl * | map_et_fait_srstdrdp (const Map *) |
Mtbl * | map_et_fait_sr2drdt (const Map *) |
Mtbl * | map_et_fait_sr2stdrdp (const Map *) |
Mtbl * | map_et_fait_d2rdx2 (const Map *) |
Mtbl * | map_et_fait_lapr_tp (const Map *) |
Mtbl * | map_et_fait_d2rdtdx (const Map *) |
Mtbl * | map_et_fait_sstd2rdpdx (const Map *) |
Mtbl * | map_et_fait_sr2d2rdt2 (const Map *) |
Mtbl * | map_et_fait_rsxdxdr (const Map *) |
Mtbl * | map_et_fait_rsx2drdx (const Map *) |
Radial mapping of rather general form.
()
This mapping relates the grid coordinates and the physical coordinates
in the following manner [see Bonazzola, Gourgoulhon & Marck, Phys. Rev. D 58 , 104020 (1998) for details]:
,
and
Lorene::Map_et::Map_et | ( | const Mg3d & | mgrille, |
const double * | r_limits | ||
) |
Standard Constructor.
mgrille | [input] Multi-domain grid on which the mapping is defined |
r_limits | [input] Array (size: number of domains + 1) of the value of r at the boundaries of the various domains :
|
Definition at line 149 of file map_et.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and set_coord().
Constructor using the equation of the surface of the star.
mgrille | [input] Multi-domain grid on which the mapping is defined It must contains at least one shell. |
r_limits | [input] Array (size: number of domains + 1) of the value of r at the boundaries of the various domains :
|
tab | [input] equation of the surface between the nucleus and the first shell in the form : ![]() ![]() ![]() |
Definition at line 223 of file map_et.C.
References alpha, Lorene::Valeur::annule(), Lorene::Valeur::annule_hard(), beta, fait_poly(), ff, get_alpha(), Lorene::Base_val::get_base_p(), get_beta(), Lorene::Tbl::get_dim(), Lorene::Tbl::get_ndim(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), gg, Lorene::max(), Lorene::min(), P_COSSIN, Lorene::Tbl::set(), Lorene::Valeur::set(), set_coord(), Lorene::Valeur::set_etat_c_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Valeur::std_base_scal(), and Lorene::Mg3d::std_base_scal().
Lorene::Map_et::Map_et | ( | const Map_et & | mpi | ) |
Copy constructor.
Definition at line 388 of file map_et.C.
References alpha, beta, fait_poly(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and set_coord().
Lorene::Map_et::Map_et | ( | const Mg3d & | mgi, |
FILE * | fich | ||
) |
Constructor from a file (see sauve(FILE*)
)
Definition at line 449 of file map_et.C.
References alpha, beta, fait_poly(), Lorene::fread_be(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and set_coord().
|
virtual |
Adaptation of the mapping to a given scalar field.
Computes the functions and
as well as the factors
and
, so that the boundaries of some domains coincide with isosurfaces of the scalar field
ent
.
ent | [input] scalar field, the isosurfaces of which are used to determine the mapping |
par | [input/output] parameters of the computation: \ par.get_int(0) : maximum number of iterations to locate zeros by the secant method \ par.get_int(1) : number of domains where the adjustment to the isosurfaces of ent is to be performed \ par.get_int(2) : number of domains nz_search where the isosurfaces will be searched : the routine scans the nz_search innermost domains, starting from the domain of index nz_search-1 . NB: the field ent must be continuous over these domains \ par.get_int(3) : 1 = performs the full computation, 0 = performs only the rescaling by the factor par.get_double_mod(0) \ par.get_int(4) : theta index of the collocation point ![]() ent \ par.get_int(5) : phi index of the collocation point ![]() ent \ par.get_int_mod(0) [output] : number of iterations actually used in the secant method \ par.get_double(0) : required absolute precision in the determination of zeros by the secant method \ par.get_double(1) : factor by which the values of ![]() ![]() par.get_double(2) : factor by which all the radial distances will be multiplied \ par.get_tbl(0) : array of values of the field ent to define the isosurfaces. |
nbr_filtre | [input] Number of the last coefficients in ![]() |
Implements Lorene::Map.
Definition at line 108 of file map_et_adapt.C.
References Lorene::Param::get_double(), Lorene::Tbl::get_etat(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Cmp::get_mp(), Lorene::Tbl::get_ndim(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Tbl::get_taille(), Lorene::Param::get_tbl(), and Lorene::Map::mg.
|
inlineinherited |
|
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 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_p | [output] ![]() |
Implements Lorene::Map.
Definition at line 183 of file map_radial_comp_rtp.C.
References Lorene::Scalar::get_etat().
|
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::get_etat().
|
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 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_t | [output] ![]() |
Implements Lorene::Map.
Definition at line 128 of file map_radial_comp_rtp.C.
References Lorene::Scalar::get_etat().
|
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::get_etat().
|
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::get_etat().
|
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::get_etat().
|
inherited |
Determines the coordinates corresponding to given absolute Cartesian coordinates (X,Y,Z).
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().
|
virtual |
Not yet implemented.
Implements Lorene::Map.
Definition at line 69 of file map_et_dalembert.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
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::Scalar::get_etat().
|
virtualinherited |
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::Scalar::get_etat().
|
virtualinherited |
Division by of a
Scalar
.
Implements Lorene::Map.
Definition at line 85 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by r of a Scalar
.
Implements Lorene::Map.
Definition at line 514 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
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::Scalar::get_etat().
|
virtualinherited |
Division by of a
Scalar
.
Implements Lorene::Map.
Definition at line 434 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by of a
Scalar
.
Implements Lorene::Map.
Definition at line 133 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by of a
Scalar
.
Implements Lorene::Map.
Definition at line 158 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
Internal function intended to be used by Map::poisson_vect
and Map::poisson_vect_oohara
.
It constructs the sets of parameters used for each scalar Poisson equation from the one for the vectorial one.
para | [input] : the Param used for the resolution of the vectorial Poisson equation : \ para.int() maximum number of iteration.\ para.double(0) relaxation parameter.\ para.double(1) required precision. \ para.tenseur_mod() source of the vectorial part at the previous step.\ para.cmp_mod() source of the scalar part at the previous step. |
i | [input] number of the scalar Poisson equation that is being solved (values from 0 to 2 for the componants of the vectorial part and 3 for the scalar one). |
Implements Lorene::Map.
Definition at line 75 of file map_poisson_vect.C.
References Lorene::Param::add_cmp_mod(), Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Param::get_cmp_mod(), Lorene::Param::get_double(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Param::get_tenseur_mod(), and Lorene::Tenseur::set().
Computes of a
Cmp
.
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 187 of file map_et_deriv.C.
References Lorene::Cmp::get_etat().
Computes of a
Scalar
.
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 215 of file map_et_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a
Scalar
if the description is affine and if it is logarithmic.
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 274 of file map_et_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a
Scalar
.
uu | [input] scalar field |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 629 of file map_et_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a
Cmp
.
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 95 of file map_et_deriv.C.
References Lorene::Cmp::get_etat().
Computes of a
Scalar
.
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 123 of file map_et_deriv.C.
References Lorene::Scalar::get_etat().
|
private |
Construction of the polynomials and
.
Definition at line 645 of file map_et.C.
References aa, bb, daa, dbb, ddaa, ddbb, Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
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.
const double * Lorene::Map_et::get_alpha | ( | ) | const |
const double * Lorene::Map_et::get_beta | ( | ) | const |
|
inlineinherited |
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
the Cartesian coordinates related to by means of usual formulae.
Definition at line 791 of file map.h.
References Lorene::Map::bvect_cart.
|
inlineinherited |
Returns the orthonormal vectorial basis associated with the coordinates
of the mapping.
Definition at line 783 of file map.h.
References Lorene::Map::bvect_spher.
const Valeur & Lorene::Map_et::get_ff | ( | ) | const |
const Valeur & Lorene::Map_et::get_gg | ( | ) | const |
|
inlineinherited |
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.
|
virtual |
Sets a new radial scale.
lambda | [input] factor by which the value of r is to be multiplied |
Implements Lorene::Map.
Definition at line 905 of file map_et.C.
References Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
virtualinherited |
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::Scalar::get_etat().
|
virtualinherited |
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::Scalar::get_etat().
Computes the integral over all space of a Cmp
.
The computed quantity is . The routine allocates a
Tbl
(size: mg->nzone
) to store the result (partial integral) in each domain and returns a pointer to it.
Implements Lorene::Map.
Definition at line 74 of file map_et_integ.C.
References Lorene::Cmp::get_etat().
Computes the angular Laplacian of a scalar field.
uu | [input] Scalar field u (represented as a Scalar ) the Laplacian ![]() |
lap | [output] Angular Laplacian of u (see documentation of Scalar |
Implements Lorene::Map.
Definition at line 283 of file map_et_lap.C.
Computes the Laplacian of a scalar field (Cmp
version).
Implements Lorene::Map.
Definition at line 180 of file map_et_lap.C.
References Lorene::Cmp::get_etat().
Computes the Laplacian of a scalar field.
uu | [input] Scalar field u (represented as a Scalar ) the Laplacian ![]() |
zec_mult_r | [input] Determines the quantity computed in the compactified external domain (CED) : \ zec_mult_r = 0 : ![]() ![]() ![]() |
lap | [output] Laplacian of u |
Implements Lorene::Map.
Definition at line 75 of file map_et_lap.C.
References Lorene::Scalar::get_etat().
|
virtual |
Returns the "angular" mapping for the outside of domain l_zone
.
Valid only for the class Map_af
.
Implements Lorene::Map.
Definition at line 1045 of file map_et.C.
References Lorene::c_est_pas_fait(), and Lorene::Map::p_mp_angu.
|
virtualinherited |
Multiplication by of a
Scalar
.
Implements Lorene::Map.
Definition at line 61 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
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::get_etat().
|
virtualinherited |
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::get_etat().
|
virtualinherited |
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::Scalar::get_etat().
|
virtualinherited |
Multiplication by of a
Scalar
.
Implements Lorene::Map.
Definition at line 355 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Multiplication by of a
Scalar
.
Implements Lorene::Map.
Definition at line 109 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Assignment to an affine mapping.
Implements Lorene::Map_radial.
Definition at line 539 of file map_et.C.
References alpha, beta, ff, 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(), gg, Lorene::Map::mg, reset_coord(), Lorene::Map::set_ori(), and Lorene::Map::set_rot_phi().
|
virtual |
Assignment to another Map_et
Definition at line 513 of file map_et.C.
References alpha, beta, ff, get_alpha(), 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(), gg, Lorene::Map::mg, reset_coord(), Lorene::Map::set_ori(), and Lorene::Map::set_rot_phi().
|
virtual |
Comparison operator (egality)
Implements Lorene::Map_radial.
Definition at line 984 of file map_et.C.
References alpha, beta, Lorene::Map::bvect_cart, Lorene::Map::bvect_spher, Lorene::diffrelmax(), ff, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_x(), Lorene::Map::get_ori_y(), Lorene::Map::get_ori_z(), gg, Lorene::max(), Lorene::Map::mg, Lorene::Map::ori_x, Lorene::Map::ori_y, and Lorene::Map::ori_z.
|
privatevirtual |
Operator >>
Implements Lorene::Map.
Definition at line 798 of file map_et.C.
References Lorene::Valeur::affiche_seuil(), alpha, beta, ff, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_p(), Lorene::Mg3d::get_type_t(), gg, Lorene::Map::mg, and val_r().
Computes the solution of a scalar Poisson equation.
Following the method explained in Sect. III.C of Bonazzola, Gourgoulhon & Marck, Phys. Rev. D 58 , 104020 (1998),
the Poisson equation is re-written as
, where
is the Laplacian in an affine mapping and R(u) contains the terms generated by the deviation of the mapping
*this
from spherical symmetry. This equation is solved by iterations. At each step J the equation effectively solved is where
with in domain no. l and
is a relaxation parameter.
source | [input] source ![]() |
par | [input/output] parameters for the iterative method: \ par.get_int(0) : [input] maximum number of iterations \ par.get_double(0) : [input] relaxation parameter ![]() par.get_double(1) : [input] required precision: the iterative method is stopped as soon as the relative difference between ![]() ![]() par.get_double(1) \ par.get_cmp_mod(0) : [input/output] input : Cmp containing ![]() Cmp must be set to zero); output: value of ![]() par.get_int_mod(0) : [output] number of iterations actually used to get the solution. |
uu | [input/output] input : previously computed value of u to start the iteration (term R(u) ) (if nothing is known a priori, uu must be set to zero); output: solution u with the boundary condition u =0 at spatial infinity. |
Implements Lorene::Map.
Definition at line 102 of file map_et_poisson.C.
References Lorene::Cmp::get_etat().
|
virtual |
Computes the solution of a 2-D Poisson equation.
The 2-D Poisson equation writes
source_mat | [input] Compactly supported part of the source ![]() |
source_quad | [input] Non-compactly supported part of the source ![]() |
par | [input/output] Parameters to control the resolution : \ par.get_double_mod(0) : [output] constant lambda such that the source of the equation effectively solved is source_mat + lambda * source_quad , in order to fulfill the virial identity GRV2. \ If the theta basis is T_SIN_I , the following arguments are required: \ par.get_int(0) : [input] maximum number of iterations \ par.get_double(0) : [input] relaxation parameter \ par.get_double(1) : [input] required precision: the iterative method is stopped as soon as the relative difference between ![]() ![]() par.get_double(1) \ par.get_cmp_mod(0) : [input/output] input : Cmp containing ![]() Cmp must be set to zero); output: value of ![]() par.get_int_mod(0) : [output] number of iterations actually used to get the solution. |
uu | [input/output] solution u with the boundary condition u =0 at spatial infinity. |
Implements Lorene::Map.
Definition at line 80 of file map_et_poisson2d.C.
References Lorene::Cmp::get_etat().
|
virtual |
Computes the solution of the generalized angular Poisson equation.
The generalized angular Poisson equation is , where
.
source | [input] source ![]() ![]() |
par | [input/output] possible parameters to control the resolution of the Poisson equation. See the actual implementation in the derived class of Map for documentation. |
uu | [input/output] solution u |
lambda | [input] coefficient ![]() |
Implements Lorene::Map.
Definition at line 593 of file map_et_poisson.C.
References Lorene::Param::get_cmp_mod(), Lorene::Param::get_double(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), and Lorene::Scalar::set_spectral_va().
|
virtualinherited |
Resolution of the elliptic equation in the case where the stellar interior is covered by a single domain.
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::get_etat().
|
virtualinherited |
Resolution of the elliptic equation in the case of a multidomain stellar interior.
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::get_etat(), and Lorene::Map_radial::poisson_compact().
|
virtual |
Computes the solution of a Poisson equation in the shell .
imposing a boundary condition at the surface of the star
source | [input] : source of the equation. |
limite | [input] : limite [num_front] contains the angular function being the boudary condition. |
par | [input] : parameters of the computation. |
pot | [output] : result. |
Implements Lorene::Map.
Definition at line 274 of file cmp_pde_frontiere.C.
References Lorene::Cmp::get_etat().
|
virtual |
Computes the solution of a scalar Poisson equation.
The regularized source
source | [input] source ![]() ![]() |
k_div | [input] regularization degree of the procedure |
nzet | [input] number of domains covering the star |
unsgam1 | [input] parameter ![]() ![]() |
par | [input/output] parameters for the iterative method: \ par.get_int(0) : [input] maximum number of iterations \ par.get_double(0) : [input] relaxation parameter ![]() par.get_double(1) : [input] required precision: the iterative method is stopped as soon as the relative difference between ![]() ![]() par.get_double(1) \ par.get_cmp_mod(0) : [input/output] input : Cmp containing ![]() Cmp must be set to zero); output: value of ![]() par.get_int_mod(0) : [output] number of iterations actually used to get the solution. |
uu | [input/output] input : previously computed value of u to start the iteration (term R(u) ) (if nothing is known a priori, uu must be set to zero); output: solution u with the boundary condition u =0 at spatial infinity. |
uu_regu | [output] solution of the regular part of the source. |
uu_div | [output] solution of the diverging part of the source. |
duu_div | [output] derivative of the diverging potential |
source_regu | [output] regularized source |
source_div | [output] diverging part of the source |
Implements Lorene::Map.
Definition at line 85 of file map_et_poisson_regu.C.
References Lorene::Cmp::get_etat().
Computes the solution of a scalar Poisson equation with a Tau method.
Following the method explained in Sect. III.C of Bonazzola, Gourgoulhon & Marck, Phys. Rev. D 58 , 104020 (1998),
the Poisson equation is re-written as
, where
is the Laplacian in an affine mapping and R(u) contains the terms generated by the deviation of the mapping
*this
from spherical symmetry. This equation is solved by iterations. At each step J the equation effectively solved is where
with in domain no. l and
is a relaxation parameter.
source | [input] source ![]() |
par | [input/output] parameters for the iterative method: \ par.get_int(0) : [input] maximum number of iterations \ par.get_double(0) : [input] relaxation parameter ![]() par.get_double(1) : [input] required precision: the iterative method is stopped as soon as the relative difference between ![]() ![]() par.get_double(1) \ par.get_cmp_mod(0) : [input/output] input : Cmp containing ![]() Cmp must be set to zero); output: value of ![]() par.get_int_mod(0) : [output] number of iterations actually used to get the solution. |
uu | [input/output] input : previously computed value of u to start the iteration (term R(u) ) (if nothing is known a priori, uu must be set to zero); output: solution u with the boundary condition u =0 at spatial infinity. |
Implements Lorene::Map.
Definition at line 351 of file map_et_poisson.C.
References Lorene::Cmp::get_etat().
Computes the radial primitive which vanishes for .
i.e. the function
uu | [input] function f (must have a dzpuis = 2) |
resu | [input] function F |
null_infty | if true (default), the primitive is null at infinity (or on the grid boundary). F is null at the center otherwise |
Implements Lorene::Map.
Definition at line 110 of file map_et_integ.C.
|
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::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.
|
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::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.
|
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::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.
|
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::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.
|
protectedvirtual |
Resets all the member Coords
.
Reimplemented from Lorene::Map_radial.
Definition at line 628 of file map_et.C.
References Lorene::Coord::del_t(), Lorene::Map_radial::reset_coord(), rsx2drdx, and rsxdxdr.
|
virtual |
Rescales the outer boundary of one domain.
The inner boundary is unchanged. The inner boundary of the next domain is changed to match the new outer boundary.
l | [input] index of the domain |
lambda | [input] factor by which the value of ![]() |
Implements Lorene::Map.
Definition at line 928 of file map_et.C.
References Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
void Lorene::Map_et::resize_extr | ( | double | lambda | ) |
Rescales the outer boundary of the outermost domain in the case of non-compactified external domain.
The inner boundary is unchanged.
lambda | [input] factor by which the value of the radius of the outermost domain is to be multiplied. |
Definition at line 55 of file map_et_resize_extr.C.
References Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
virtual |
Save in a file.
Reimplemented from Lorene::Map_radial.
Definition at line 779 of file map_et.C.
References alpha, beta, ff, Lorene::fwrite_be(), Lorene::Mg3d::get_nzone(), gg, Lorene::Map::mg, Lorene::Map_radial::sauve(), and Lorene::Valeur::sauve().
void Lorene::Map_et::set_alpha | ( | double | alpha0, |
int | l | ||
) |
Modifies the value of in domain no. l.
Definition at line 424 of file map_et.C.
References alpha, and reset_coord().
void Lorene::Map_et::set_beta | ( | double | beta0, |
int | l | ||
) |
Modifies the value of in domain no. l.
Definition at line 435 of file map_et.C.
References beta, and reset_coord().
|
private |
Assignement of the building functions to the member Coords
.
Definition at line 584 of file map_et.C.
References Lorene::Map::cosp, Lorene::Map::cost, Lorene::Map_radial::d2rdtdx, Lorene::Map_radial::d2rdx2, Lorene::Map_radial::drdt, Lorene::Map_radial::dxdr, Lorene::Map_radial::lapr_tp, Lorene::Map::phi, Lorene::Map::r, rsx2drdx, rsxdxdr, Lorene::Coord::set(), Lorene::Map::sinp, Lorene::Map::sint, 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, Lorene::Map::tet, Lorene::Map::x, Lorene::Map::xa, Lorene::Map_radial::xsr, Lorene::Map::y, Lorene::Map::ya, Lorene::Map::z, and Lorene::Map::za.
void Lorene::Map_et::set_ff | ( | const Valeur & | ffi | ) |
Assigns a given value to the function .
Definition at line 562 of file map_et.C.
References ff, and reset_coord().
void Lorene::Map_et::set_gg | ( | const Valeur & | ggi | ) |
Assigns a given value to the function .
Definition at line 570 of file map_et.C.
References gg, and reset_coord().
|
inherited |
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().
|
inherited |
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().
Computes of a
Cmp
.
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 337 of file map_et_deriv.C.
References Lorene::Cmp::get_etat().
Computes of a
Scalar
.
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 391 of file map_et_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a
Cmp
.
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 485 of file map_et_deriv.C.
References Lorene::Cmp::get_etat().
Computes of a
Scalar
.
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 540 of file map_et_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a
Scalar
.
uu | [input] scalar field |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 674 of file map_et_deriv.C.
References Lorene::Scalar::get_etat().
|
virtual |
Computes the domain index l and the value of corresponding to a point given by its physical coordinates
.
This version enables to pass some parameters to control the accuracy of the computation.
rr | [input] value of r |
theta | [input] value of ![]() |
pphi | [input] value of ![]() |
par | [input/output] parameters to control the accuracy of the computation: \ par.get_int(0) : [input] maximum number of iterations in the secant method to locate ![]() par.get_int_mod(0) : [output] effective number of iterations used \ par.get_double(0) : [input] absolute precision in the secant method to locate ![]() |
l | [output] value of the domain index |
xi | [output] value of ![]() |
Implements Lorene::Map.
Definition at line 169 of file map_et_radius.C.
References Lorene::Param::get_double(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
virtual |
Computes the domain index l and the value of corresponding to a point given by its physical coordinates
.
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 146 of file map_et_radius.C.
References Lorene::Param::add_double(), Lorene::Param::add_int(), and Lorene::Param::add_int_mod().
|
virtual |
Computes the domain index l and the value of corresponding to a point of arbitrary r but collocation values 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: \ par.get_int(0) : [input] maximum number of iterations in the secant method to locate ![]() par.get_int_mod(0) : [output] effective number of iterations used \ par.get_double(0) : [input] absolute precision in the secant method to locate ![]() |
l | [output] value of the domain index |
xi | [output] value of ![]() |
Implements Lorene::Map_radial.
Definition at line 418 of file map_et_radius.C.
References Lorene::Param::get_double(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
virtual |
Returns the value of the radial coordinate r for a given in a given domain.
l | [input] index of the domain |
xi | [input] value of ![]() |
theta | [input] value of ![]() |
pphi | [input] value of ![]() |
Implements Lorene::Map.
Definition at line 92 of file map_et_radius.C.
References ff, Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and Lorene::Valeur::val_point().
|
virtual |
Returns the value of the radial coordinate r for a given and a given collocation point in
in a given domain.
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 367 of file map_et_radius.C.
References ff, Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
private |
|
protectedinherited |
|
protectedinherited |
|
inherited |
|
inherited |
|
private |
|
private |
|
private |
|
private |
|
inherited |
|
inherited |
|
private |
|
private |
|
inherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
inherited |
|
inherited |
|
protectedinherited |
Coord Lorene::Map_et::rsx2drdx |
Coord Lorene::Map_et::rsxdxdr |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
private |
|
private |