LORENE
|
Base class for pure radial mappings. More...
#include <map.h>
Public Member Functions | |
virtual | ~Map_radial () |
Destructor. More... | |
virtual void | operator= (const Map_af &)=0 |
Assignment to an affine mapping. More... | |
virtual void | sauve (FILE *) const |
Save in a file. More... | |
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 and a given collocation point in in a given domain. More... | |
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 corresponding to a point of arbitrary r but collocation values of . More... | |
virtual bool | operator== (const Map &) const =0 |
Comparison operator (egality) 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_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 (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, 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 of a Scalar . More... | |
virtual void | div_rsint (Scalar &) const |
Division by of a 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 of a Scalar . More... | |
virtual void | div_cost (Scalar &) const |
Division by of a Scalar . More... | |
virtual void | mult_sint (Scalar &) const |
Multiplication by of a Scalar . More... | |
virtual void | div_sint (Scalar &) const |
Division by of a Scalar . More... | |
virtual void | div_tant (Scalar &) const |
Division by of a 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 component (with respect to 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 component (with respect to 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 in the case where the stellar interior is covered by a single domain. More... | |
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 in the case of a multidomain stellar interior. More... | |
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 associated with the coordinates of the mapping. More... | |
const Base_vect_cart & | get_bvect_cart () const |
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e. More... | |
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... | |
virtual const Map_af & | mp_angu (int) const =0 |
Returns the "angular" mapping for the outside of domain l_zone . More... | |
void | convert_absolute (double xx, double yy, double zz, double &rr, double &theta, double &pphi) const |
Determines the coordinates corresponding to given absolute Cartesian coordinates (X,Y,Z). More... | |
virtual double | val_r (int l, double xi, double theta, double pphi) const =0 |
Returns the value of the radial coordinate r for a given in a given domain. More... | |
virtual void | val_lx (double rr, double theta, double pphi, int &l, double &xi) const =0 |
Computes the domain index l and the value of corresponding to a point given by its physical coordinates . More... | |
virtual void | val_lx (double rr, double theta, double pphi, const Param &par, int &l, double &xi) const =0 |
Computes the domain index l and the value of corresponding to a point given by its physical coordinates . More... | |
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... | |
virtual void | homothetie (double lambda)=0 |
Sets a new radial scale. More... | |
virtual void | resize (int l, double lambda)=0 |
Rescales the outer boundary of one domain. More... | |
virtual void | adapt (const Cmp &ent, const Param &par, int nbr=0)=0 |
Adaptation of the mapping to a given scalar field. More... | |
virtual void | dsdxi (const Cmp &ci, Cmp &resu) const =0 |
Computes of a Cmp . More... | |
virtual void | dsdxi (const Scalar &uu, Scalar &resu) const =0 |
Computes of a Scalar . More... | |
virtual void | dsdr (const Cmp &ci, Cmp &resu) const =0 |
Computes of a Cmp . More... | |
virtual void | dsdr (const Scalar &uu, Scalar &resu) const =0 |
Computes of a Scalar . More... | |
virtual void | srdsdt (const Cmp &ci, Cmp &resu) const =0 |
Computes of a Cmp . More... | |
virtual void | srdsdt (const Scalar &uu, Scalar &resu) const =0 |
Computes of a Scalar . More... | |
virtual void | srstdsdp (const Cmp &ci, Cmp &resu) const =0 |
Computes of a Cmp . More... | |
virtual void | srstdsdp (const Scalar &uu, Scalar &resu) const =0 |
Computes of a Scalar . More... | |
virtual void | dsdradial (const Scalar &uu, Scalar &resu) const =0 |
Computes of a Scalar if the description is affine and if it is logarithmic. More... | |
virtual void | dsdt (const Scalar &uu, Scalar &resu) const =0 |
Computes of a Scalar . More... | |
virtual void | stdsdp (const Scalar &uu, Scalar &resu) const =0 |
Computes of a Scalar . More... | |
virtual void | laplacien (const Scalar &uu, int zec_mult_r, Scalar &lap) const =0 |
Computes the Laplacian of a scalar field. More... | |
virtual void | laplacien (const Cmp &uu, int zec_mult_r, Cmp &lap) const =0 |
Computes the Laplacian of a scalar field (Cmp version). More... | |
virtual void | lapang (const Scalar &uu, Scalar &lap) const =0 |
Computes the angular Laplacian of a scalar field. More... | |
virtual void | primr (const Scalar &uu, Scalar &resu, bool null_infty) const =0 |
Computes the radial primitive which vanishes for . More... | |
virtual Tbl * | integrale (const Cmp &) const =0 |
Computes the integral over all space of a Cmp . More... | |
virtual void | poisson (const Cmp &source, Param &par, Cmp &uu) const =0 |
Computes the solution of a scalar Poisson equation. More... | |
virtual void | poisson_tau (const Cmp &source, Param &par, Cmp &uu) const =0 |
Computes the solution of a scalar Poisson equationwith a Tau method. More... | |
virtual void | poisson_falloff (const Cmp &source, Param &par, Cmp &uu, int k_falloff) const =0 |
virtual void | poisson_ylm (const Cmp &source, Param &par, Cmp &pot, int nylm, double *intvec) const =0 |
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 =0 |
Computes the solution of a scalar Poisson equation. More... | |
virtual void | poisson_angu (const Scalar &source, Param &par, Scalar &uu, double lambda=0) const =0 |
Computes the solution of the generalized angular Poisson equation. More... | |
virtual void | poisson_angu (const Cmp &source, Param &par, Cmp &uu, double lambda=0) const =0 |
virtual Param * | donne_para_poisson_vect (Param ¶, int i) const =0 |
Function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara . More... | |
virtual void | poisson_frontiere (const Cmp &source, const Valeur &limite, int raccord, int num_front, Cmp &pot, double=0., double=0.) const =0 |
Computes the solution of a Poisson equation from the domain num_front+1 . More... | |
virtual void | poisson_frontiere_double (const Cmp &source, const Valeur &lim_func, const Valeur &lim_der, int num_zone, Cmp &pot) const =0 |
virtual void | poisson_interne (const Cmp &source, const Valeur &limite, Param &par, Cmp &pot) const =0 |
Computes the solution of a Poisson equation in the shell, imposing a boundary condition at the surface of the star. More... | |
virtual void | poisson2d (const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu) const =0 |
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 =0 |
Performs one time-step integration of the d'Alembert scalar equation. More... | |
Public Attributes | |
Coord | xsr |
in the nucleus; \ 1/R in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | dxdr |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | drdt |
in the nucleus and in the non-compactified shells; \ in the compactified external domain (CED). More... | |
Coord | stdrdp |
in the nucleus and in the non-compactified shells; \ in the compactified external domain (CED). More... | |
Coord | srdrdt |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | srstdrdp |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | sr2drdt |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | sr2stdrdp |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | d2rdx2 |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | lapr_tp |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | d2rdtdx |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | sstd2rdpdx |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | sr2d2rdt2 |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | r |
r coordinate centered on the grid More... | |
Coord | tet |
coordinate centered on the grid More... | |
Coord | phi |
coordinate centered on the grid More... | |
Coord | sint |
More... | |
Coord | cost |
More... | |
Coord | sinp |
More... | |
Coord | cosp |
More... | |
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 | |
Map_radial (const Mg3d &mgrid) | |
Constructor from a grid (protected to make Map_radial an abstract class) More... | |
Map_radial (const Map_radial &mp) | |
Copy constructor. More... | |
Map_radial (const Mg3d &, FILE *) | |
Constructor from a file (see sauve(FILE* ) ) More... | |
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 associated with the coordinates of the mapping. More... | |
Base_vect_cart | bvect_cart |
Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e. More... | |
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... | |
Base class for pure radial mappings.
()
A pure radial mapping is a mapping of the type , , . The class 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 92 of file map_radial.C.
|
protected |
Copy constructor.
Definition at line 99 of file map_radial.C.
|
protected |
Constructor from a file (see sauve(FILE* )
)
Definition at line 106 of file map_radial.C.
|
virtual |
Destructor.
Definition at line 113 of file map_radial.C.
|
pure virtualinherited |
Adaptation of the mapping to a given scalar field.
This is a virtual function: see the actual implementations in the derived classes for the meaning of the various parameters.
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
inlineinherited |
|
virtual |
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] -component of the vector |
Implements Lorene::Map.
Definition at line 186 of file map_radial_comp_rtp.C.
References Lorene::Scalar::get_etat().
|
virtual |
Cmp
version
Implements Lorene::Map.
Definition at line 179 of file map_radial_comp_rtp.C.
References comp_p_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 |
Implements Lorene::Map.
Definition at line 75 of file map_radial_comp_rtp.C.
References Lorene::Scalar::get_etat().
|
virtual |
Cmp
version
Implements Lorene::Map.
Definition at line 68 of file map_radial_comp_rtp.C.
References comp_r_from_cartesian().
|
virtual |
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] -component of the vector |
Implements Lorene::Map.
Definition at line 131 of file map_radial_comp_rtp.C.
References Lorene::Scalar::get_etat().
|
virtual |
Cmp
version
Implements Lorene::Map.
Definition at line 124 of file map_radial_comp_rtp.C.
References comp_t_from_cartesian().
|
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] -component of the vector |
v_phi | [input] -component of the vector |
v_x | [output] x-component of the vector |
Implements Lorene::Map.
Definition at line 79 of file map_radial_comp_xyz.C.
References Lorene::Scalar::get_etat().
|
virtual |
Cmp
version
Implements Lorene::Map.
Definition at line 71 of file map_radial_comp_xyz.C.
References comp_x_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] -component of the vector |
v_phi | [input] -component of the vector |
v_y | [output] y-component of the vector |
Implements Lorene::Map.
Definition at line 138 of file map_radial_comp_xyz.C.
References Lorene::Scalar::get_etat().
|
virtual |
Cmp
version
Implements Lorene::Map.
Definition at line 129 of file map_radial_comp_xyz.C.
References comp_y_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] -component of the vector |
v_z | [output] z-component of the vector |
Implements Lorene::Map.
Definition at line 195 of file map_radial_comp_xyz.C.
References Lorene::Scalar::get_etat().
|
virtual |
Cmp
version
Implements Lorene::Map.
Definition at line 187 of file map_radial_comp_xyz.C.
References comp_z_from_spherical().
|
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 305 of file map.C.
References Lorene::Map::ori_x, Lorene::Map::ori_y, Lorene::Map::ori_z, Lorene::Map::rot_phi, and Lorene::sqrt().
|
pure virtualinherited |
Performs one time-step integration of the d'Alembert scalar equation.
par | [input/output] possible parameters to control the resolution of the wave equation. See the actual implementation in the derived class of Map for documentation. Note that, at least, param must contain the time step as first double parameter. |
fJp1 | [output] solution at time J+1 with boundary conditions of outgoing radiation (not exact!) |
fJ | [input] solution at time J |
fJm1 | [input] solution at time J-1 |
source | [input] source of the d'Alembert equation . |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
virtual |
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 751 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
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 646 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Division by of a Scalar
.
Implements Lorene::Map.
Definition at line 88 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Division by r of a Scalar
.
Implements Lorene::Map.
Definition at line 517 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Division by r (in the compactified external domain only) of a Scalar
.
Implements Lorene::Map.
Definition at line 588 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Division by of a Scalar
.
Implements Lorene::Map.
Definition at line 437 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Division by of a Scalar
.
Implements Lorene::Map.
Definition at line 136 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Division by of a Scalar
.
Implements Lorene::Map.
Definition at line 161 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
pure virtualinherited |
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). |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
Computes of a Cmp
.
Note that in the compactified external domain (CED), it computes .
ci | [input] field to consider |
resu | [output] derivative of ci |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
Computes of a Scalar
.
uu | [input] scalar field |
resu | [output] derivative of uu |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
Computes of a Cmp
.
Note that in the compactified external domain (CED), it computes .
ci | [input] field to consider |
resu | [output] derivative of ci |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
inherited |
Returns the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart
.
Definition at line 334 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 324 of file map.C.
References Lorene::Map::bvect_spher, and Lorene::Map::p_flat_met_spher.
|
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 809 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 801 of file map.h.
References Lorene::Map::bvect_spher.
|
inlineinherited |
Gives the Mg3d
on which the mapping is defined.
Definition at line 783 of file map.h.
References Lorene::Map::mg.
|
inlineinherited |
Returns the x coordinate of the origin.
Definition at line 786 of file map.h.
References Lorene::Map::ori_x.
|
inlineinherited |
Returns the y coordinate of the origin.
Definition at line 788 of file map.h.
References Lorene::Map::ori_y.
|
inlineinherited |
Returns the z coordinate of the origin.
Definition at line 790 of file map.h.
References Lorene::Map::ori_z.
|
inlineinherited |
Returns the angle between the x –axis and X –axis.
Definition at line 793 of file map.h.
References Lorene::Map::rot_phi.
|
pure virtualinherited |
Sets a new radial scale.
lambda | [input] factor by which the value of r is to be multiplied |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
virtual |
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 802 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
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 699 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.
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
Computes the angular Laplacian of a scalar field.
uu | [input] Scalar field u (represented as a Scalar ) the Laplacian of which is to be computed |
lap | [output] Angular Laplacian of u (see documentation of Scalar |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
Computes the Laplacian of a scalar field.
uu | [input] Scalar field u (represented as a Scalar ) the Laplacian of which is to be computed |
zec_mult_r | [input] Determines the quantity computed in the compactified external domain (CED) : \ zec_mult_r = 0 : \ zec_mult_r = 2 : \ zec_mult_r = 4 (default) : |
lap | [output] Laplacian of u |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
Computes the Laplacian of a scalar field (Cmp
version).
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
Returns the "angular" mapping for the outside of domain l_zone
.
Valid only for the class Map_af
.
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
virtual |
Multiplication by of a Scalar
.
Implements Lorene::Map.
Definition at line 64 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Multiplication by r of a Scalar
, the dzpuis
of uu
is not changed.
Implements Lorene::Map.
Definition at line 147 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Multiplication by r of a Cmp
.
In the CED, there is only a decrement of dzpuis
Implements Lorene::Map.
Definition at line 222 of file map_radial_r_manip.C.
References Lorene::Cmp::get_etat().
|
virtual |
Multiplication by r (in the compactified external domain only) of a Scalar
.
Implements Lorene::Map.
Definition at line 299 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Multiplication by of a Scalar
.
Implements Lorene::Map.
Definition at line 358 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Multiplication by of a Scalar
.
Implements Lorene::Map.
Definition at line 112 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
pure virtual |
Assignment to an affine mapping.
Implements Lorene::Map.
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtual |
Comparison operator (egality)
Implements Lorene::Map.
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
Computes the solution of a scalar Poisson equation.
source | [input] source of the Poisson equation . |
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 with the boundary condition u =0 at spatial infinity. |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
Computes the solution of a 2-D Poisson equation.
The 2-D Poisson equation writes
source_mat | [input] Compactly supported part of the source of the 2-D Poisson equation (typically matter terms) |
source_quad | [input] Non-compactly supported part of the source of the 2-D Poisson equation (typically quadratic terms) |
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 with the boundary condition u =0 at spatial infinity. |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
Computes the solution of the generalized angular Poisson equation.
The generalized angular Poisson equation is , where .
source | [input] source of the equation . |
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 in the above equation (default value = 0) |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
virtual |
Resolution of the elliptic equation in the case where the stellar interior is covered by a single domain.
source | [input] source of the above equation |
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 and is greater than 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 to start the iteration (if nothing is known a priori, psi must be set to zero); output: solution which satisfies . |
Implements Lorene::Map.
Definition at line 158 of file map_radial_poisson_cpt.C.
References Lorene::Cmp::get_etat().
|
virtual |
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 of the above equation |
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 which satisfies . |
Implements Lorene::Map.
Definition at line 456 of file map_radial_poisson_cpt.C.
References Lorene::Cmp::get_etat(), and poisson_compact().
|
pure virtualinherited |
Computes the solution of a Poisson equation from the domain num_front+1
.
imposing a boundary condition at the boundary between the domains num_front
and num_front+1
.
source | [input] : source of the equation. |
limite | [input] : limite [num_front] contains the angular function being the boudary condition. |
raccord | [input] : 1 for the Dirichlet problem and 2 for the Neumann one and 3 for Dirichlet-Neumann. |
num_front | [input] : index of the boudary at which the boundary condition has to be imposed. |
pot | [output] : result. |
fact_dir | [input] : Valeur by which we multiply the quantity we solve. (in the case of Dirichlet-Neumann boundary condition.) |
fact_neu | [input] : Valeur by which we multiply the radial derivative of the quantity we solve. (in the case of Dirichlet-Neumann boundary condition.) |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
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. |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
Computes the solution of a scalar Poisson equation.
The regularized source
source | [input] source of the Poisson equation . |
k_div | [input] regularization degree of the procedure |
nzet | [input] number of domains covering the star |
unsgam1 | [input] parameter where denotes the adiabatic index. |
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 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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
Computes the solution of a scalar Poisson equationwith a Tau method.
source | [input] source of the Poisson equation . |
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 with the boundary condition u =0 at spatial infinity. |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
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 64 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.
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 179 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.
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 65 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.
|
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 . |
Implements Lorene::Map.
Definition at line 199 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.
Reimplemented in Lorene::Map_et.
Definition at line 129 of file map_radial.C.
References d2rdtdx, d2rdx2, Lorene::Coord::del_t(), drdt, dxdr, lapr_tp, Lorene::Map::reset_coord(), sr2d2rdt2, sr2drdt, sr2stdrdp, srdrdt, srstdrdp, sstd2rdpdx, stdrdp, and xsr.
|
pure virtualinherited |
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 defining the outer boundary of the domain is to be multiplied. |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
virtual |
Save in a file.
Reimplemented from Lorene::Map.
Reimplemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
Definition at line 119 of file map_radial.C.
References Lorene::Map::sauve().
|
inherited |
Sets a new origin.
Definition at line 256 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 266 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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
Computes of a Cmp
.
Note that in the compactified external domain (CED), it computes .
ci | [input] field to consider |
resu | [output] derivative of ci |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
Computes of a Scalar
.
uu | [input] scalar field |
resu | [output] derivative of uu |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
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 |
l | [output] value of the domain index |
xi | [output] value of |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure 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 |
l | [output] value of the domain index |
xi | [output] value of |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure virtualinherited |
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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
pure 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 |
Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.
|
protectedinherited |
|
protectedinherited |
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 |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
inherited |
|
inherited |
|
protectedinherited |
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 |
|
inherited |
|
inherited |
Coord Lorene::Map_radial::xsr |
|
inherited |
|
inherited |