LORENE

Radial mapping of rather general form. More...

#include <map.h>

Inheritance diagram for Lorene::Map_et:
Lorene::Map_radial Lorene::Map

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 $F(\theta',\phi')$. More...
 
void set_gg (const Valeur &)
 Assigns a given value to the function $G(\theta',\phi')$. More...
 
virtual const Map_afmp_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 $\alpha$ in each domain) More...
 
const double * get_beta () const
 Returns a pointer on the array beta (values of $\beta$ in each domain) More...
 
const Valeurget_ff () const
 Returns a (constant) reference to the function $F(\theta',\phi')$. More...
 
const Valeurget_gg () const
 Returns a (constant) reference to the function $G(\theta',\phi')$. More...
 
virtual double val_r (int l, double xi, double theta, double pphi) const
 Returns the value of the radial coordinate r for a given $(\xi, \theta', \phi')$ in a given domain. More...
 
virtual void val_lx (double rr, double theta, double pphi, int &l, double &xi) const
 Computes the domain index l and the value of $\xi$ corresponding to a point given by its physical coordinates $(r, \theta, \phi)$. More...
 
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 $\xi$ corresponding to a point given by its physical coordinates $(r, \theta, \phi)$. More...
 
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 $\xi$ and a given collocation point in $(\theta', \phi')$ in a given domain. More...
 
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 $\xi$ corresponding to a point of arbitrary r but collocation values of $(\theta, \phi)$. More...
 
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 $\alpha$ in domain no. l. More...
 
void set_beta (double beta0, int l)
 Modifies the value of $\beta$ in domain no. l. More...
 
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 $\partial/ \partial \xi$ of a Cmp. More...
 
virtual void dsdr (const Cmp &ci, Cmp &resu) const
 Computes $\partial/ \partial r$ of a Cmp. More...
 
virtual void srdsdt (const Cmp &ci, Cmp &resu) const
 Computes $1/r \partial/ \partial \theta$ of a Cmp. More...
 
virtual void srstdsdp (const Cmp &ci, Cmp &resu) const
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Cmp. More...
 
virtual void dsdxi (const Scalar &uu, Scalar &resu) const
 Computes $\partial/ \partial \xi$ of a Scalar. More...
 
virtual void dsdr (const Scalar &uu, Scalar &resu) const
 Computes $\partial/ \partial r$ of a Scalar. More...
 
virtual void dsdradial (const Scalar &uu, Scalar &resu) const
 Computes $\partial/ \partial r$ of a Scalar if the description is affine and $\partial/ \partial \ln r$ if it is logarithmic. More...
 
virtual void srdsdt (const Scalar &uu, Scalar &resu) const
 Computes $1/r \partial/ \partial \theta$ of a Scalar. More...
 
virtual void srstdsdp (const Scalar &uu, Scalar &resu) const
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Scalar. More...
 
virtual void dsdt (const Scalar &uu, Scalar &resu) const
 Computes $\partial/ \partial \theta$ of a Scalar. More...
 
virtual void stdsdp (const Scalar &uu, Scalar &resu) const
 Computes $1/\sin\theta \partial/ \partial \varphi$ of a 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 $r\to \infty$. More...
 
virtual Tblintegrale (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 void poisson_angu (const Cmp &source, Param &par, Cmp &uu, double lambda=0) const
 
virtual Paramdonne_para_poisson_vect (Param &para, 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 $r\sin\theta$ of a Scalar. More...
 
virtual void div_rsint (Scalar &) const
 Division by $r\sin\theta$ 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 $\cos\theta$ of a Scalar. More...
 
virtual void div_cost (Scalar &) const
 Division by $\cos\theta$ of a Scalar. More...
 
virtual void mult_sint (Scalar &) const
 Multiplication by $\sin\theta$ of a Scalar. More...
 
virtual void div_sint (Scalar &) const
 Division by $\sin\theta$ of a Scalar. More...
 
virtual void div_tant (Scalar &) const
 Division by $\tan\theta$ 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 $\theta$ 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 $\phi$ 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 $ a \Delta\psi + {\bf b} \cdot \nabla \psi = \sigma$ 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 $ a \Delta\psi + {\bf b} \cdot \nabla \psi = \sigma$ in the case of a multidomain stellar interior. More...
 
const Mg3dget_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_spherget_bvect_spher () const
 Returns the orthonormal vectorial basis $(\partial/\partial r,1/r\partial/\partial \theta, 1/(r\sin\theta)\partial/\partial \phi)$ associated with the coordinates $(r, \theta, \phi)$ of the mapping. More...
 
const Base_vect_cartget_bvect_cart () const
 Returns the Cartesian basis $(\partial/\partial x,\partial/\partial y,\partial/\partial z)$ associated with the coordinates (x,y,z) of the mapping, i.e. More...
 
const Metric_flatflat_met_spher () const
 Returns the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher. More...
 
const Metric_flatflat_met_cart () const
 Returns the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart. More...
 
const Cmpcmp_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 $(r,\theta,\phi)$ corresponding to given absolute Cartesian coordinates (X,Y,Z). 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...
 

Public Attributes

Coord rsxdxdr
 $1/(\partial R/\partial \xi) R/\xi$ in the nucleus; \ $1/(\partial R/\partial \xi) R/(\xi + \beta/\alpha)$ in the shells; \ $1/(\partial U/\partial \xi) U/(\xi-1)$ in the outermost compactified domain. More...
 
Coord rsx2drdx
 $[ R/ (\alpha \xi + \beta) ]^2 (\partial R/\partial \xi) / \alpha$ in the nucleus and the shells; \ $\partial U/\partial \xi / \alpha$ in the outermost compactified domain. More...
 
Coord xsr
 $\xi/R$ in the nucleus; \ 1/R in the non-compactified shells; \ $(\xi-1)/U$ in the compactified outer domain. More...
 
Coord dxdr
 $1/(\partial R/\partial\xi) = \partial \xi /\partial r$ in the nucleus and in the non-compactified shells; \ $-1/(\partial U/\partial\xi) = - \partial \xi /\partial u$ in the compactified outer domain. More...
 
Coord drdt
 $\partial R/\partial\theta'$ in the nucleus and in the non-compactified shells; \ $-\partial U/\partial\theta'$ in the compactified external domain (CED). More...
 
Coord stdrdp
 ${1\over\sin\theta} \partial R/\partial\varphi'$ in the nucleus and in the non-compactified shells; \ $-{1\over\sin\theta}\partial U/\partial\varphi'$ in the compactified external domain (CED). More...
 
Coord srdrdt
 $1/R \times (\partial R/\partial\theta')$ in the nucleus and in the non-compactified shells; \ $-1/U \times (\partial U/\partial\theta)$ in the compactified outer domain. More...
 
Coord srstdrdp
 $1/(R\sin\theta) \times (\partial R/\partial\varphi')$ in the nucleus and in the non-compactified shells; \ $-1/(U\sin\theta) \times (\partial U/\partial\varphi')$ in the compactified outer domain. More...
 
Coord sr2drdt
 $1/R^2 \times (\partial R/\partial\theta')$ in the nucleus and in the non-compactified shells; \ $-1/U^2 \times (\partial U/\partial\theta')$ in the compactified outer domain. More...
 
Coord sr2stdrdp
 $1/(R^2\sin\theta) \times (\partial R/\partial\varphi')$ in the nucleus and in the non-compactified shells; \ $-1/(U^2\sin\theta) \times (\partial U/\partial\varphi')$ in the compactified outer domain. More...
 
Coord d2rdx2
 $\partial^2 R/\partial\xi^2$ in the nucleus and in the non-compactified shells; \ $-\partial^2 U/\partial\xi^2 $ in the compactified outer domain. More...
 
Coord lapr_tp
 $1/R^2 \times [ 1/\sin(\theta)\times \partial /\partial\theta' (\sin\theta \partial R /\partial\theta') + 1/\sin^2\theta \partial^2 R /\partial{\varphi'}^2] $ in the nucleus and in the non-compactified shells; \ $- 1/U^2 \times [ 1/\sin(\theta)\times \partial /\partial\theta' (\sin\theta \partial U /\partial\theta') + 1/\sin^2\theta \partial^2 U /\partial{\varphi'}^2] $ in the compactified outer domain. More...
 
Coord d2rdtdx
 $\partial^2 R/\partial\xi\partial\theta'$ in the nucleus and in the non-compactified shells; \ $-\partial^2 U/\partial\xi\partial\theta'$ in the compactified outer domain. More...
 
Coord sstd2rdpdx
 $1/\sin\theta \times \partial^2 R/\partial\xi\partial\varphi'$ in the nucleus and in the non-compactified shells; \ $-1/\sin\theta \times \partial^2 U/\partial\xi\partial\varphi' $ in the compactified outer domain. More...
 
Coord sr2d2rdt2
 $1/R^2 \partial^2 R/\partial{\theta'}^2$ in the nucleus and in the non-compactified shells; \ $-1/U^2 \partial^2 U/\partial{\theta'}^2$ in the compactified outer domain. More...
 
Coord r
 r coordinate centered on the grid More...
 
Coord tet
 $\theta$ coordinate centered on the grid More...
 
Coord phi
 $\phi$ coordinate centered on the grid More...
 
Coord sint
 $\sin\theta$ More...
 
Coord cost
 $\cos\theta$ More...
 
Coord sinp
 $\sin\phi$ More...
 
Coord cosp
 $\cos\phi$ 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

virtual void reset_coord ()
 Resets all the member Coords. More...
 

Protected Attributes

const Mg3dmg
 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 $(\partial/\partial r,1/r\partial/\partial \theta, 1/(r\sin\theta)\partial/\partial \phi)$ associated with the coordinates $(r, \theta, \phi)$ of the mapping. More...
 
Base_vect_cart bvect_cart
 Cartesian basis $(\partial/\partial x,\partial/\partial y,\partial/\partial z)$ associated with the coordinates (x,y,z) of the mapping, i.e. More...
 
Metric_flatp_flat_met_spher
 Pointer onto the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher. More...
 
Metric_flatp_flat_met_cart
 Pointer onto the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart. More...
 
Cmpp_cmp_zero
 The null Cmp. More...
 
Map_afp_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 $A(\xi)$ and $B(\xi)$. More...
 
virtual ostream & operator>> (ostream &) const
 Operator >> More...
 

Private Attributes

double * alpha
 Array (size: mg->nzone ) of the values of $\alpha$ in each domain. More...
 
double * beta
 Array (size: mg->nzone ) of the values of $\beta$ in each domain. More...
 
Tbl ** aa
 Array (size: mg->nzone ) of Tbl which stores the values of $A(\xi)$ in each domain. More...
 
Tbl ** daa
 Array (size: mg->nzone ) of Tbl which stores the values of $A'(\xi)$ in each domain. More...
 
Tbl ** ddaa
 Array (size: mg->nzone ) of Tbl which stores the values of $A''(\xi)$ in each domain. More...
 
Tbl aasx
 Values at the nr collocation points of $A(\xi)/\xi$ in the nucleus. More...
 
Tbl aasx2
 Values at the nr collocation points of $A(\xi)/\xi^2$ in the nucleus. More...
 
Tbl zaasx
 Values at the nr collocation points of $A(\xi)/(\xi-1)$ in the outermost compactified domain. More...
 
Tbl zaasx2
 Values at the nr collocation points of $A(\xi)/(\xi-1)^2$ in the outermost compactified domain. More...
 
Tbl ** bb
 Array (size: mg->nzone ) of Tbl which stores the values of $B(\xi)$ in each domain. More...
 
Tbl ** dbb
 Array (size: mg->nzone ) of Tbl which stores the values of $B'(\xi)$ in each domain. More...
 
Tbl ** ddbb
 Array (size: mg->nzone ) of Tbl which stores the values of $B''(\xi)$ in each domain. More...
 
Tbl bbsx
 Values at the nr collocation points of $B(\xi)/\xi$ in the nucleus. More...
 
Tbl bbsx2
 Values at the nr collocation points of $B(\xi)/\xi^2$ in the nucleus. More...
 
Valeur ff
 Values of the function $F(\theta', \phi')$ at the nt*np angular collocation points in each domain. More...
 
Valeur gg
 Values of the function $G(\theta', \phi')$ at the nt*np angular collocation points in each domain. More...
 

Friends

Mtblmap_et_fait_r (const Map *)
 
Mtblmap_et_fait_tet (const Map *)
 
Mtblmap_et_fait_phi (const Map *)
 
Mtblmap_et_fait_sint (const Map *)
 
Mtblmap_et_fait_cost (const Map *)
 
Mtblmap_et_fait_sinp (const Map *)
 
Mtblmap_et_fait_cosp (const Map *)
 
Mtblmap_et_fait_x (const Map *)
 
Mtblmap_et_fait_y (const Map *)
 
Mtblmap_et_fait_z (const Map *)
 
Mtblmap_et_fait_xa (const Map *)
 
Mtblmap_et_fait_ya (const Map *)
 
Mtblmap_et_fait_za (const Map *)
 
Mtblmap_et_fait_xsr (const Map *)
 
Mtblmap_et_fait_dxdr (const Map *)
 
Mtblmap_et_fait_drdt (const Map *)
 
Mtblmap_et_fait_stdrdp (const Map *)
 
Mtblmap_et_fait_srdrdt (const Map *)
 
Mtblmap_et_fait_srstdrdp (const Map *)
 
Mtblmap_et_fait_sr2drdt (const Map *)
 
Mtblmap_et_fait_sr2stdrdp (const Map *)
 
Mtblmap_et_fait_d2rdx2 (const Map *)
 
Mtblmap_et_fait_lapr_tp (const Map *)
 
Mtblmap_et_fait_d2rdtdx (const Map *)
 
Mtblmap_et_fait_sstd2rdpdx (const Map *)
 
Mtblmap_et_fait_sr2d2rdt2 (const Map *)
 
Mtblmap_et_fait_rsxdxdr (const Map *)
 
Mtblmap_et_fait_rsx2drdx (const Map *)
 

Detailed Description

Radial mapping of rather general form.

()

This mapping relates the grid coordinates $(\xi, \theta', \phi')$ and the physical coordinates $(r, \theta, \phi)$ in the following manner [see Bonazzola, Gourgoulhon & Marck, Phys. Rev. D 58 , 104020 (1998) for details]: $\theta=\theta'$, $\phi=\phi'$ and

  • $r = \alpha [\xi + A(\xi) F(\theta', \phi') + B(\xi) G(\theta', \phi')] + \beta$ in non-compactified domains,
  • $ u={1\over r} = \alpha [\xi + A(\xi) F(\theta', \phi')] + \beta$ in the (usually outermost) compactified domain, where $\alpha$ and $\beta$ are constant in each domain, $A(\xi)$ and $B(\xi)$ are constant polynomials defined by
  • $A(\xi) = 3 \xi^4 - 2 \xi^6$ and $B(\xi) = (5\xi^3 - 3\xi^5)/2$ in the nucleus (innermost domain which contains r =0);
  • $A(\xi) = (\xi^3 - 3\xi + 2)/4$ and $B(\xi) = (-\xi^3 + 3\xi +2)/4$ in the other domains. The functions $F(\theta', \phi')$ and $G(\theta', \phi')$ depend on the domain under consideration and define the boundaries of this domain.

Definition at line 2783 of file map.h.

Constructor & Destructor Documentation

◆ Map_et() [1/4]

Lorene::Map_et::Map_et ( const Mg3d mgrille,
const double *  r_limits 
)

Standard Constructor.

Parameters
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 :
  • r_limits[l] : inner boundary of the domain no. l
  • r_limits[l+1] : outer boundary of the domain no. l

Definition at line 155 of file map_et.C.

References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and set_coord().

◆ Map_et() [2/4]

Lorene::Map_et::Map_et ( const Mg3d mgrille,
const double *  r_limits,
const Tbl tab 
)

Constructor using the equation of the surface of the star.

Parameters
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 :
  • r_limits[l] : inner boundary of the domain no. l
  • r_limits[l+1] : outer boundary of the domain no. l The first value is not used.
tab[input] equation of the surface between the nucleus and the first shell in the form : ${\rm tab}(k,j) = r(\phi_k,\theta_j)$, where $\phi_k$ and $\theta_j$ are the values of the angular colocation points.

Definition at line 229 of file map_et.C.

References alpha, Lorene::Valeur::annule(), Lorene::Valeur::annule_hard(), beta, fait_poly(), ff, Lorene::Base_val::get_base_p(), 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, P_COSSIN_P, 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().

◆ Map_et() [3/4]

Lorene::Map_et::Map_et ( const Map_et mpi)

Copy constructor.

Definition at line 411 of file map_et.C.

References alpha, beta, fait_poly(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and set_coord().

◆ Map_et() [4/4]

Lorene::Map_et::Map_et ( const Mg3d mgi,
FILE *  fich 
)

Constructor from a file (see sauve(FILE*) )

Definition at line 472 of file map_et.C.

References alpha, beta, fait_poly(), Lorene::fread_be(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and set_coord().

◆ ~Map_et()

Lorene::Map_et::~Map_et ( )
virtual

Destructor.

Definition at line 509 of file map_et.C.

References aa, alpha, bb, beta, daa, dbb, ddaa, ddbb, Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.

Member Function Documentation

◆ adapt()

void Lorene::Map_et::adapt ( const Cmp ent,
const Param par,
int  nbr_filtre = 0 
)
virtual

Adaptation of the mapping to a given scalar field.

Computes the functions $F(\theta',\phi')$ and $G(\theta',\phi')$ as well as the factors $\alpha$ and $\beta$, so that the boundaries of some domains coincide with isosurfaces of the scalar field ent .

Parameters
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 $(\theta_*, \phi_*)$ [using the notations of Bonazzola, Gourgoulhon & Marck, Phys. Rev. D 58 , 104020 (1998)] defining an isosurface of ent \ par.get_int(5) : phi index of the collocation point $(\theta_*, \phi_*)$ [using the notations of Bonazzola, Gourgoulhon & Marck, Phys. Rev. D 58 , 104020 (1998)] defining an isosurface of 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 $\lambda$ and $\mu$ [using the notations of Bonazzola, Gourgoulhon & Marck, Phys. Rev. D 58 , 104020 (1998)] will be multiplied : 1 = regular mapping, 0 = contracting mapping \ 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 $\varphi$ set to zero.

Implements Lorene::Map.

Definition at line 111 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.

◆ cmp_zero()

const Cmp& Lorene::Map::cmp_zero ( ) const
inlineinherited

Returns the null Cmp defined on *this.

To be used by the Tenseur class when necessary to return a null Cmp .

Definition at line 825 of file map.h.

References Lorene::Map::p_cmp_zero.

◆ comp_p_from_cartesian() [1/2]

void Lorene::Map_radial::comp_p_from_cartesian ( const Scalar v_x,
const Scalar v_y,
Scalar v_p 
) const
virtualinherited

Computes the Spherical $\phi$ component (with respect to bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart .

Parameters
v_x[input] x-component of the vector
v_y[input] y-component of the vector
v_p[output] $\phi$-component of the vector

Implements Lorene::Map.

Definition at line 186 of file map_radial_comp_rtp.C.

References Lorene::Scalar::get_etat().

◆ comp_p_from_cartesian() [2/2]

void Lorene::Map_radial::comp_p_from_cartesian ( const Cmp v_x,
const Cmp v_y,
Cmp v_p 
) const
virtualinherited

Cmp version

Implements Lorene::Map.

Definition at line 179 of file map_radial_comp_rtp.C.

References Lorene::Map_radial::comp_p_from_cartesian().

◆ comp_r_from_cartesian() [1/2]

void Lorene::Map_radial::comp_r_from_cartesian ( const Scalar v_x,
const Scalar v_y,
const Scalar v_z,
Scalar v_r 
) const
virtualinherited

Computes the Spherical r component (with respect to bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart .

Parameters
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().

◆ comp_r_from_cartesian() [2/2]

void Lorene::Map_radial::comp_r_from_cartesian ( const Cmp v_x,
const Cmp v_y,
const Cmp v_z,
Cmp v_r 
) const
virtualinherited

Cmp version

Implements Lorene::Map.

Definition at line 68 of file map_radial_comp_rtp.C.

References Lorene::Map_radial::comp_r_from_cartesian().

◆ comp_t_from_cartesian() [1/2]

void Lorene::Map_radial::comp_t_from_cartesian ( const Scalar v_x,
const Scalar v_y,
const Scalar v_z,
Scalar v_t 
) const
virtualinherited

Computes the Spherical $\theta$ component (with respect to bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart .

Parameters
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] $\theta$-component of the vector

Implements Lorene::Map.

Definition at line 131 of file map_radial_comp_rtp.C.

References Lorene::Scalar::get_etat().

◆ comp_t_from_cartesian() [2/2]

void Lorene::Map_radial::comp_t_from_cartesian ( const Cmp v_x,
const Cmp v_y,
const Cmp v_z,
Cmp v_t 
) const
virtualinherited

Cmp version

Implements Lorene::Map.

Definition at line 124 of file map_radial_comp_rtp.C.

References Lorene::Map_radial::comp_t_from_cartesian().

◆ comp_x_from_spherical() [1/2]

void Lorene::Map_radial::comp_x_from_spherical ( const Scalar v_r,
const Scalar v_theta,
const Scalar v_phi,
Scalar v_x 
) const
virtualinherited

Computes the Cartesian x component (with respect to bvect_cart) of a vector given by its spherical components with respect to bvect_spher.

Parameters
v_r[input] r -component of the vector
v_theta[input] $\theta$-component of the vector
v_phi[input] $\phi$-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().

◆ comp_x_from_spherical() [2/2]

void Lorene::Map_radial::comp_x_from_spherical ( const Cmp v_r,
const Cmp v_theta,
const Cmp v_phi,
Cmp v_x 
) const
virtualinherited

Cmp version

Implements Lorene::Map.

Definition at line 71 of file map_radial_comp_xyz.C.

References Lorene::Map_radial::comp_x_from_spherical().

◆ comp_y_from_spherical() [1/2]

void Lorene::Map_radial::comp_y_from_spherical ( const Scalar v_r,
const Scalar v_theta,
const Scalar v_phi,
Scalar v_y 
) const
virtualinherited

Computes the Cartesian y component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher .

Parameters
v_r[input] r -component of the vector
v_theta[input] $\theta$-component of the vector
v_phi[input] $\phi$-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().

◆ comp_y_from_spherical() [2/2]

void Lorene::Map_radial::comp_y_from_spherical ( const Cmp v_r,
const Cmp v_theta,
const Cmp v_phi,
Cmp v_y 
) const
virtualinherited

Cmp version

Implements Lorene::Map.

Definition at line 129 of file map_radial_comp_xyz.C.

References Lorene::Map_radial::comp_y_from_spherical().

◆ comp_z_from_spherical() [1/2]

void Lorene::Map_radial::comp_z_from_spherical ( const Scalar v_r,
const Scalar v_theta,
Scalar v_z 
) const
virtualinherited

Computes the Cartesian z component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher .

Parameters
v_r[input] r -component of the vector
v_theta[input] $\theta$-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().

◆ comp_z_from_spherical() [2/2]

void Lorene::Map_radial::comp_z_from_spherical ( const Cmp v_r,
const Cmp v_theta,
Cmp v_z 
) const
virtualinherited

Cmp version

Implements Lorene::Map.

Definition at line 187 of file map_radial_comp_xyz.C.

References Lorene::Map_radial::comp_z_from_spherical().

◆ convert_absolute()

void Lorene::Map::convert_absolute ( double  xx,
double  yy,
double  zz,
double &  rr,
double &  theta,
double &  pphi 
) const
inherited

Determines the coordinates $(r,\theta,\phi)$ corresponding to given absolute Cartesian coordinates (X,Y,Z).

Parameters
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 $\theta$
pphi[output] value of $\phi$

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().

◆ dalembert()

void Lorene::Map_et::dalembert ( Param par,
Scalar fJp1,
const Scalar fJ,
const Scalar fJm1,
const Scalar source 
) const
virtual

Not yet implemented.

Implements Lorene::Map.

Definition at line 72 of file map_et_dalembert.C.

References Lorene::Scalar::get_etat().

◆ dec2_dzpuis()

void Lorene::Map_radial::dec2_dzpuis ( Scalar ci) const
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 751 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ dec_dzpuis()

void Lorene::Map_radial::dec_dzpuis ( Scalar ci) const
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 646 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ div_cost()

void Lorene::Map_radial::div_cost ( Scalar ci) const
virtualinherited

Division by $\cos\theta$ of a Scalar.

Implements Lorene::Map.

Definition at line 88 of file map_radial_th_manip.C.

References Lorene::Scalar::get_etat().

◆ div_r()

void Lorene::Map_radial::div_r ( Scalar ci) const
virtualinherited

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().

◆ div_r_zec()

void Lorene::Map_radial::div_r_zec ( Scalar uu) const
virtualinherited

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().

◆ div_rsint()

void Lorene::Map_radial::div_rsint ( Scalar ci) const
virtualinherited

Division by $r\sin\theta$ of a Scalar.

Implements Lorene::Map.

Definition at line 437 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ div_sint()

void Lorene::Map_radial::div_sint ( Scalar ci) const
virtualinherited

Division by $\sin\theta$ of a Scalar.

Implements Lorene::Map.

Definition at line 136 of file map_radial_th_manip.C.

References Lorene::Scalar::get_etat().

◆ div_tant()

void Lorene::Map_radial::div_tant ( Scalar ci) const
virtualinherited

Division by $\tan\theta$ of a Scalar.

Implements Lorene::Map.

Definition at line 161 of file map_radial_th_manip.C.

References Lorene::Scalar::get_etat().

◆ donne_para_poisson_vect()

Param * Lorene::Map_et::donne_para_poisson_vect ( Param para,
int  i 
) const
virtual

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.

Parameters
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).
Returns
the pointer on the parameter set used for solving the scalar Poisson equation labelled by i .

Implements Lorene::Map.

Definition at line 78 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().

◆ dsdr() [1/2]

void Lorene::Map_et::dsdr ( const Cmp ci,
Cmp resu 
) const
virtual

Computes $\partial/ \partial r$ of a Cmp.

Note that in the compactified external domain (CED), it computes $-\partial/ \partial u = r^2 \partial/ \partial r$.

Parameters
ci[input] field to consider
resu[output] derivative of ci

Implements Lorene::Map.

Definition at line 190 of file map_et_deriv.C.

References Lorene::Cmp::get_etat().

◆ dsdr() [2/2]

void Lorene::Map_et::dsdr ( const Scalar uu,
Scalar resu 
) const
virtual

Computes $\partial/ \partial r$ 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.

Parameters
uu[input] field to consider
resu[output] derivative of uu

Implements Lorene::Map.

Definition at line 218 of file map_et_deriv.C.

References Lorene::Scalar::get_etat().

◆ dsdradial()

void Lorene::Map_et::dsdradial ( const Scalar uu,
Scalar resu 
) const
virtual

Computes $\partial/ \partial r$ of a Scalar if the description is affine and $\partial/ \partial \ln r$ 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.

Parameters
uu[input] field to consider
resu[output] derivative of uu

Implements Lorene::Map.

Definition at line 277 of file map_et_deriv.C.

References Lorene::Scalar::get_etat().

◆ dsdt()

void Lorene::Map_et::dsdt ( const Scalar uu,
Scalar resu 
) const
virtual

Computes $\partial/ \partial \theta$ of a Scalar.

Parameters
uu[input] scalar field
resu[output] derivative of uu

Implements Lorene::Map.

Definition at line 632 of file map_et_deriv.C.

References Lorene::Scalar::get_etat().

◆ dsdxi() [1/2]

void Lorene::Map_et::dsdxi ( const Cmp ci,
Cmp resu 
) const
virtual

Computes $\partial/ \partial \xi$ of a Cmp.

Note that in the compactified external domain (CED), it computes $-\partial/ \partial u = \xi^2 \partial/ \partial \xi$.

Parameters
ci[input] field to consider
resu[output] derivative of ci

Implements Lorene::Map.

Definition at line 98 of file map_et_deriv.C.

References Lorene::Cmp::get_etat().

◆ dsdxi() [2/2]

void Lorene::Map_et::dsdxi ( const Scalar uu,
Scalar resu 
) const
virtual

Computes $\partial/ \partial \xi$ 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.

Parameters
uu[input] field to consider
resu[output] derivative of uu

Implements Lorene::Map.

Definition at line 126 of file map_et_deriv.C.

References Lorene::Scalar::get_etat().

◆ fait_poly()

void Lorene::Map_et::fait_poly ( )
private

Construction of the polynomials $A(\xi)$ and $B(\xi)$.

Definition at line 668 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.

◆ flat_met_cart()

const Metric_flat & Lorene::Map::flat_met_cart ( ) const
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.

◆ flat_met_spher()

const Metric_flat & Lorene::Map::flat_met_spher ( ) const
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.

◆ get_alpha()

const double * Lorene::Map_et::get_alpha ( ) const

Returns a pointer on the array alpha (values of $\alpha$ in each domain)

Definition at line 1049 of file map_et.C.

References alpha.

◆ get_beta()

const double * Lorene::Map_et::get_beta ( ) const

Returns a pointer on the array beta (values of $\beta$ in each domain)

Definition at line 1053 of file map_et.C.

References beta.

◆ get_bvect_cart()

const Base_vect_cart& Lorene::Map::get_bvect_cart ( ) const
inlineinherited

Returns the Cartesian basis $(\partial/\partial x,\partial/\partial y,\partial/\partial z)$ associated with the coordinates (x,y,z) of the mapping, i.e.

the Cartesian coordinates related to $(r, \theta, \phi)$ by means of usual formulae.

Definition at line 809 of file map.h.

References Lorene::Map::bvect_cart.

◆ get_bvect_spher()

const Base_vect_spher& Lorene::Map::get_bvect_spher ( ) const
inlineinherited

Returns the orthonormal vectorial basis $(\partial/\partial r,1/r\partial/\partial \theta, 1/(r\sin\theta)\partial/\partial \phi)$ associated with the coordinates $(r, \theta, \phi)$ of the mapping.

Definition at line 801 of file map.h.

References Lorene::Map::bvect_spher.

◆ get_ff()

const Valeur & Lorene::Map_et::get_ff ( ) const

Returns a (constant) reference to the function $F(\theta',\phi')$.

Definition at line 1057 of file map_et.C.

References ff.

◆ get_gg()

const Valeur & Lorene::Map_et::get_gg ( ) const

Returns a (constant) reference to the function $G(\theta',\phi')$.

Definition at line 1061 of file map_et.C.

References gg.

◆ get_mg()

const Mg3d* Lorene::Map::get_mg ( ) const
inlineinherited

Gives the Mg3d on which the mapping is defined.

Definition at line 783 of file map.h.

References Lorene::Map::mg.

◆ get_ori_x()

double Lorene::Map::get_ori_x ( ) const
inlineinherited

Returns the x coordinate of the origin.

Definition at line 786 of file map.h.

References Lorene::Map::ori_x.

◆ get_ori_y()

double Lorene::Map::get_ori_y ( ) const
inlineinherited

Returns the y coordinate of the origin.

Definition at line 788 of file map.h.

References Lorene::Map::ori_y.

◆ get_ori_z()

double Lorene::Map::get_ori_z ( ) const
inlineinherited

Returns the z coordinate of the origin.

Definition at line 790 of file map.h.

References Lorene::Map::ori_z.

◆ get_rot_phi()

double Lorene::Map::get_rot_phi ( ) const
inlineinherited

Returns the angle between the x –axis and X –axis.

Definition at line 793 of file map.h.

References Lorene::Map::rot_phi.

◆ homothetie()

void Lorene::Map_et::homothetie ( double  lambda)
virtual

Sets a new radial scale.

Parameters
lambda[input] factor by which the value of r is to be multiplied

Implements Lorene::Map.

Definition at line 928 of file map_et.C.

References Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.

◆ inc2_dzpuis()

void Lorene::Map_radial::inc2_dzpuis ( Scalar ci) const
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 802 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ inc_dzpuis()

void Lorene::Map_radial::inc_dzpuis ( Scalar ci) const
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 699 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ integrale()

Tbl * Lorene::Map_et::integrale ( const Cmp ci) const
virtual

Computes the integral over all space of a Cmp.

The computed quantity is $\int u \, r^2 \sin\theta \, dr\, d\theta \, d\phi$. 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 77 of file map_et_integ.C.

References Lorene::Cmp::get_etat().

◆ lapang()

void Lorene::Map_et::lapang ( const Scalar uu,
Scalar lap 
) const
virtual

Computes the angular Laplacian of a scalar field.

Parameters
uu[input] Scalar field u (represented as a Scalar) the Laplacian $\Delta u$ of which is to be computed
lap[output] Angular Laplacian of u (see documentation of Scalar

Implements Lorene::Map.

Definition at line 286 of file map_et_lap.C.

◆ laplacien() [1/2]

void Lorene::Map_et::laplacien ( const Scalar uu,
int  zec_mult_r,
Scalar lap 
) const
virtual

Computes the Laplacian of a scalar field.

Parameters
uu[input] Scalar field u (represented as a Scalar ) the Laplacian $\Delta u$ of which is to be computed
zec_mult_r[input] Determines the quantity computed in the compactified external domain (CED) : \ zec_mult_r = 0 : $\Delta u$ \ zec_mult_r = 2 : $r^2 \, \Delta u$ \ zec_mult_r = 4 (default) : $r^4 \, \Delta u$
lap[output] Laplacian of u

Implements Lorene::Map.

Definition at line 78 of file map_et_lap.C.

References Lorene::Scalar::get_etat().

◆ laplacien() [2/2]

void Lorene::Map_et::laplacien ( const Cmp uu,
int  zec_mult_r,
Cmp lap 
) const
virtual

Computes the Laplacian of a scalar field (Cmp version).

Implements Lorene::Map.

Definition at line 183 of file map_et_lap.C.

References Lorene::Cmp::get_etat().

◆ mp_angu()

const Map_af & Lorene::Map_et::mp_angu ( int  ) const
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 1068 of file map_et.C.

References Lorene::c_est_pas_fait(), and Lorene::Map::p_mp_angu.

◆ mult_cost()

void Lorene::Map_radial::mult_cost ( Scalar ci) const
virtualinherited

Multiplication by $\cos\theta$ of a Scalar.

Implements Lorene::Map.

Definition at line 64 of file map_radial_th_manip.C.

References Lorene::Scalar::get_etat().

◆ mult_r() [1/2]

void Lorene::Map_radial::mult_r ( Scalar uu) const
virtualinherited

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().

◆ mult_r() [2/2]

void Lorene::Map_radial::mult_r ( Cmp ci) const
virtualinherited

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().

◆ mult_r_zec()

void Lorene::Map_radial::mult_r_zec ( Scalar ci) const
virtualinherited

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().

◆ mult_rsint()

void Lorene::Map_radial::mult_rsint ( Scalar ci) const
virtualinherited

Multiplication by $r\sin\theta$ of a Scalar.

Implements Lorene::Map.

Definition at line 358 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ mult_sint()

void Lorene::Map_radial::mult_sint ( Scalar ci) const
virtualinherited

Multiplication by $\sin\theta$ of a Scalar.

Implements Lorene::Map.

Definition at line 112 of file map_radial_th_manip.C.

References Lorene::Scalar::get_etat().

◆ operator=() [1/2]

◆ operator=() [2/2]

◆ operator==()

◆ operator>>()

ostream & Lorene::Map_et::operator>> ( ostream &  ost) const
privatevirtual

◆ poisson()

void Lorene::Map_et::poisson ( const Cmp source,
Param par,
Cmp uu 
) const
virtual

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 $\Delta u = \sigma$ is re-written as $a \tilde\Delta u = \sigma + R(u)$, where $\tilde\Delta$ 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 $\tilde\Delta u^{J+1} = s^J$ where

\[ s^J = 1/a_l^{\rm max} \{ {\tt source} + R(u^J) + (a_l^{\rm max}-a) [ \lambda s^{J-1} + (1-\lambda) s^{J-2} ] \} \ , \]

with $a_l^{\rm max} := \max(a)$ in domain no. l and $\lambda$ is a relaxation parameter.

Parameters
source[input] source $\sigma$ of the Poisson equation
par[input/output] parameters for the iterative method: \ par.get_int(0) : [input] maximum number of iterations \ par.get_double(0) : [input] relaxation parameter $\lambda$ \ par.get_double(1) : [input] required precision: the iterative method is stopped as soon as the relative difference between $u^J$ and $u^{J-1}$ is greater than par.get_double(1) \ par.get_cmp_mod(0) : [input/output] input : Cmp containing $s^{J-1}$ (cf. the above equation) to start the iteration (if nothing is known a priori, this Cmp must be set to zero); output: value of $s^{J-1}$, to used in a next call to the routine \ 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 111 of file map_et_poisson.C.

References Lorene::Cmp::get_etat().

◆ poisson2d()

void Lorene::Map_et::poisson2d ( const Cmp source_mat,
const Cmp source_quad,
Param par,
Cmp uu 
) const
virtual

Computes the solution of a 2-D Poisson equation.

The 2-D Poisson equation writes

\[ {\partial^2 u\over\partial r^2} + {1\over r} {\partial u \over \partial r} + {1\over r^2} {\partial^2 u\over\partial \theta^2} = \sigma \ . \]

Parameters
source_mat[input] Compactly supported part of the source $\sigma$ of the 2-D Poisson equation (typically matter terms)
source_quad[input] Non-compactly supported part of the source $\sigma$ of the 2-D Poisson equation (typically quadratic terms)
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 $u^J$ and $u^{J-1}$ is greater than par.get_double(1) \ par.get_cmp_mod(0) : [input/output] input : Cmp containing $s^{J-1}$ to start the iteration (if nothing is known a priori, this Cmp must be set to zero); output: value of $s^{J-1}$, to used in a next call to the routine \ 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 83 of file map_et_poisson2d.C.

References Lorene::Cmp::get_etat().

◆ poisson_angu()

void Lorene::Map_et::poisson_angu ( const Scalar source,
Param par,
Scalar uu,
double  lambda = 0 
) const
virtual

Computes the solution of the generalized angular Poisson equation.

The generalized angular Poisson equation is $\Delta_{\theta\varphi} u + \lambda u = \sigma$, where $\Delta_{\theta\varphi} u := \frac{\partial^2 u} {\partial \theta^2} + \frac{1}{\tan \theta} \frac{\partial u} {\partial \theta} +\frac{1}{\sin^2 \theta}\frac{\partial^2 u} {\partial \varphi^2}$.

Parameters
source[input] source $\sigma$ of the equation $\Delta_{\theta\varphi} u + \lambda u = \sigma$.
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 $\lambda$ in the above equation (default value = 0)

Implements Lorene::Map.

Definition at line 602 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().

◆ poisson_compact() [1/2]

void Lorene::Map_radial::poisson_compact ( const Cmp source,
const Cmp aa,
const Tenseur bb,
const Param par,
Cmp psi 
) const
virtualinherited

Resolution of the elliptic equation $ a \Delta\psi + {\bf b} \cdot \nabla \psi = \sigma$ in the case where the stellar interior is covered by a single domain.

Parameters
source[input] source $\sigma$ 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 $\psi^J$ and $\psi^{J-1}$ is greater than par.get_double(0) \ par.get_double(1) : [input] relaxation parameter $\lambda$ \ par.get_int_mod(0) : [output] number of iterations actually used to get the solution.
psi[input/output]: input : previously computed value of $\psi$ to start the iteration (if nothing is known a priori, psi must be set to zero); output: solution $\psi$ which satisfies $\psi(0)=0$.

Implements Lorene::Map.

Definition at line 158 of file map_radial_poisson_cpt.C.

References Lorene::Cmp::get_etat().

◆ poisson_compact() [2/2]

void Lorene::Map_radial::poisson_compact ( int  nzet,
const Cmp source,
const Cmp aa,
const Tenseur bb,
const Param par,
Cmp psi 
) const
virtualinherited

Resolution of the elliptic equation $ a \Delta\psi + {\bf b} \cdot \nabla \psi = \sigma$ in the case of a multidomain stellar interior.

Parameters
nzet[input] number of domains covering the stellar interior
source[input] source $\sigma$ 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 $\psi$ which satisfies $\psi(0)=0$.

Implements Lorene::Map.

Definition at line 456 of file map_radial_poisson_cpt.C.

References Lorene::Cmp::get_etat(), and Lorene::Map_radial::poisson_compact().

◆ poisson_frontiere()

void Lorene::Map_et::poisson_frontiere ( const Cmp ,
const Valeur ,
int  ,
int  ,
Cmp ,
double  = 0.,
double  = 0. 
) const
virtual

Not yet implemented.

Implements Lorene::Map.

Definition at line 127 of file cmp_pde_frontiere.C.

◆ poisson_interne()

void Lorene::Map_et::poisson_interne ( const Cmp source,
const Valeur limite,
Param par,
Cmp pot 
) const
virtual

Computes the solution of a Poisson equation in the shell .

imposing a boundary condition at the surface of the star

Parameters
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 277 of file cmp_pde_frontiere.C.

References Lorene::Cmp::get_etat().

◆ poisson_regular()

void Lorene::Map_et::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
virtual

Computes the solution of a scalar Poisson equation.

The regularized source

Parameters
source[input] source $\sigma$ of the Poisson equation $\Delta u = \sigma$.
k_div[input] regularization degree of the procedure
nzet[input] number of domains covering the star
unsgam1[input] parameter $1/(\gamma-1)$ where $\gamma$ denotes the adiabatic index.
par[input/output] parameters for the iterative method: \ par.get_int(0) : [input] maximum number of iterations \ par.get_double(0) : [input] relaxation parameter $\lambda$ \ par.get_double(1) : [input] required precision: the iterative method is stopped as soon as the relative difference between $u^J$ and $u^{J-1}$ is greater than par.get_double(1) \ par.get_cmp_mod(0) : [input/output] input : Cmp containing $s^{J-1}$ (cf. the above equation) to start the iteration (if nothing is known a priori, this Cmp must be set to zero); output: value of $s^{J-1}$, to used in a next call to the routine \ 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 88 of file map_et_poisson_regu.C.

References Lorene::Cmp::get_etat().

◆ poisson_tau()

void Lorene::Map_et::poisson_tau ( const Cmp source,
Param par,
Cmp uu 
) const
virtual

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 $\Delta u = \sigma$ is re-written as $a \tilde\Delta u = \sigma + R(u)$, where $\tilde\Delta$ 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 $\tilde\Delta u^{J+1} = s^J$ where

\[ s^J = 1/a_l^{\rm max} \{ {\tt source} + R(u^J) + (a_l^{\rm max}-a) [ \lambda s^{J-1} + (1-\lambda) s^{J-2} ] \} \ , \]

with $a_l^{\rm max} := \max(a)$ in domain no. l and $\lambda$ is a relaxation parameter.

Parameters
source[input] source $\sigma$ of the Poisson equation
par[input/output] parameters for the iterative method: \ par.get_int(0) : [input] maximum number of iterations \ par.get_double(0) : [input] relaxation parameter $\lambda$ \ par.get_double(1) : [input] required precision: the iterative method is stopped as soon as the relative difference between $u^J$ and $u^{J-1}$ is greater than par.get_double(1) \ par.get_cmp_mod(0) : [input/output] input : Cmp containing $s^{J-1}$ (cf. the above equation) to start the iteration (if nothing is known a priori, this Cmp must be set to zero); output: value of $s^{J-1}$, to used in a next call to the routine \ 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 360 of file map_et_poisson.C.

References Lorene::Cmp::get_etat().

◆ primr()

void Lorene::Map_et::primr ( const Scalar uu,
Scalar resu,
bool  null_infty 
) const
virtual

Computes the radial primitive which vanishes for $r\to \infty$.

i.e. the function $ F(r,\theta,\varphi) = \int_r^\infty f(r',\theta,\varphi) \, dr' $

Parameters
uu[input] function f (must have a dzpuis = 2)
resu[input] function F
null_inftyif 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 113 of file map_et_integ.C.

◆ reevaluate() [1/2]

void Lorene::Map_radial::reevaluate ( const Map mp_prev,
int  nzet,
Cmp uu 
) const
virtualinherited

Recomputes the values of a Cmp at the collocation points after a change in the mapping.

Parameters
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 $0\le {\tt l} \le {\tt nzet-1}$; 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.

◆ reevaluate() [2/2]

void Lorene::Map_radial::reevaluate ( const Map mp_prev,
int  nzet,
Scalar uu 
) const
virtualinherited

Recomputes the values of a Scalar at the collocation points after a change in the mapping.

Parameters
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 $0\le {\tt l} \le {\tt nzet-1}$; 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.

◆ reevaluate_symy() [1/2]

void Lorene::Map_radial::reevaluate_symy ( const Map mp_prev,
int  nzet,
Cmp uu 
) const
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.

Parameters
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 $0\le {\tt l} \le {\tt nzet-1}$; 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.

◆ reevaluate_symy() [2/2]

void Lorene::Map_radial::reevaluate_symy ( const Map mp_prev,
int  nzet,
Scalar uu 
) const
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.

Parameters
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 $0\le {\tt l} \le {\tt nzet-1}$; 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.

◆ reset_coord()

void Lorene::Map_et::reset_coord ( )
protectedvirtual

Resets all the member Coords.

Reimplemented from Lorene::Map_radial.

Definition at line 651 of file map_et.C.

References Lorene::Coord::del_t(), Lorene::Map_radial::reset_coord(), rsx2drdx, and rsxdxdr.

◆ resize()

void Lorene::Map_et::resize ( int  l,
double  lambda 
)
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.

Parameters
l[input] index of the domain
lambda[input] factor by which the value of $R(\theta, \varphi)$ defining the outer boundary of the domain is to be multiplied.

Implements Lorene::Map.

Definition at line 951 of file map_et.C.

References Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.

◆ resize_extr()

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.

Parameters
lambda[input] factor by which the value of the radius of the outermost domain is to be multiplied.

Definition at line 58 of file map_et_resize_extr.C.

References Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.

◆ sauve()

void Lorene::Map_et::sauve ( FILE *  fich) const
virtual

Save in a file.

Reimplemented from Lorene::Map_radial.

Definition at line 802 of file map_et.C.

References alpha, beta, ff, Lorene::fwrite_be(), Lorene::Mg3d::get_nzone(), gg, Lorene::Map::mg, Lorene::Valeur::sauve(), and Lorene::Map_radial::sauve().

◆ set_alpha()

void Lorene::Map_et::set_alpha ( double  alpha0,
int  l 
)

Modifies the value of $\alpha$ in domain no. l.

Definition at line 447 of file map_et.C.

References alpha, and reset_coord().

◆ set_beta()

void Lorene::Map_et::set_beta ( double  beta0,
int  l 
)

Modifies the value of $\beta$ in domain no. l.

Definition at line 458 of file map_et.C.

References beta, and reset_coord().

◆ set_coord()

◆ set_ff()

void Lorene::Map_et::set_ff ( const Valeur ffi)

Assigns a given value to the function $F(\theta',\phi')$.

Definition at line 585 of file map_et.C.

References ff, and reset_coord().

◆ set_gg()

void Lorene::Map_et::set_gg ( const Valeur ggi)

Assigns a given value to the function $G(\theta',\phi')$.

Definition at line 593 of file map_et.C.

References gg, and reset_coord().

◆ set_ori()

void Lorene::Map::set_ori ( double  xa0,
double  ya0,
double  za0 
)
inherited

◆ set_rot_phi()

void Lorene::Map::set_rot_phi ( double  phi0)
inherited

◆ srdsdt() [1/2]

void Lorene::Map_et::srdsdt ( const Cmp ci,
Cmp resu 
) const
virtual

Computes $1/r \partial/ \partial \theta$ of a Cmp.

Note that in the compactified external domain (CED), it computes $1/u \partial/ \partial \theta = r \partial/ \partial \theta$.

Parameters
ci[input] field to consider
resu[output] derivative of ci

Implements Lorene::Map.

Definition at line 340 of file map_et_deriv.C.

References Lorene::Cmp::get_etat().

◆ srdsdt() [2/2]

void Lorene::Map_et::srdsdt ( const Scalar uu,
Scalar resu 
) const
virtual

Computes $1/r \partial/ \partial \theta$ 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.

Parameters
uu[input] field to consider
resu[output] derivative of uu

Implements Lorene::Map.

Definition at line 394 of file map_et_deriv.C.

References Lorene::Scalar::get_etat().

◆ srstdsdp() [1/2]

void Lorene::Map_et::srstdsdp ( const Cmp ci,
Cmp resu 
) const
virtual

Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Cmp.

Note that in the compactified external domain (CED), it computes $1/(u\sin\theta) \partial/ \partial \phi = r/\sin\theta \partial/ \partial \phi$.

Parameters
ci[input] field to consider
resu[output] derivative of ci

Implements Lorene::Map.

Definition at line 488 of file map_et_deriv.C.

References Lorene::Cmp::get_etat().

◆ srstdsdp() [2/2]

void Lorene::Map_et::srstdsdp ( const Scalar uu,
Scalar resu 
) const
virtual

Computes $1/(r\sin\theta) \partial/ \partial \phi$ 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.

Parameters
uu[input] field to consider
resu[output] derivative of uu

Implements Lorene::Map.

Definition at line 543 of file map_et_deriv.C.

References Lorene::Scalar::get_etat().

◆ stdsdp()

void Lorene::Map_et::stdsdp ( const Scalar uu,
Scalar resu 
) const
virtual

Computes $1/\sin\theta \partial/ \partial \varphi$ of a Scalar.

Parameters
uu[input] scalar field
resu[output] derivative of uu

Implements Lorene::Map.

Definition at line 677 of file map_et_deriv.C.

References Lorene::Scalar::get_etat().

◆ val_lx() [1/2]

void Lorene::Map_et::val_lx ( double  rr,
double  theta,
double  pphi,
int &  l,
double &  xi 
) const
virtual

Computes the domain index l and the value of $\xi$ corresponding to a point given by its physical coordinates $(r, \theta, \phi)$.

Parameters
rr[input] value of r
theta[input] value of $\theta$
pphi[input] value of $\phi$
l[output] value of the domain index
xi[output] value of $\xi$

Implements Lorene::Map.

Definition at line 149 of file map_et_radius.C.

References Lorene::Param::add_double(), Lorene::Param::add_int(), and Lorene::Param::add_int_mod().

◆ val_lx() [2/2]

void Lorene::Map_et::val_lx ( double  rr,
double  theta,
double  pphi,
const Param par,
int &  l,
double &  xi 
) const
virtual

Computes the domain index l and the value of $\xi$ corresponding to a point given by its physical coordinates $(r, \theta, \phi)$.

This version enables to pass some parameters to control the accuracy of the computation.

Parameters
rr[input] value of r
theta[input] value of $\theta$
pphi[input] value of $\phi$
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 $\xi$ \ par.get_int_mod(0) : [output] effective number of iterations used \ par.get_double(0) : [input] absolute precision in the secant method to locate $\xi$
l[output] value of the domain index
xi[output] value of $\xi$

Implements Lorene::Map.

Definition at line 172 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.

◆ val_lx_jk()

void Lorene::Map_et::val_lx_jk ( double  rr,
int  j,
int  k,
const Param par,
int &  l,
double &  xi 
) const
virtual

Computes the domain index l and the value of $\xi$ corresponding to a point of arbitrary r but collocation values of $(\theta, \phi)$.

Parameters
rr[input] value of r
j[input] index of the collocation point in $\theta$
k[input] index of the collocation point in $\phi$
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 $\xi$ \ par.get_int_mod(0) : [output] effective number of iterations used \ par.get_double(0) : [input] absolute precision in the secant method to locate $\xi$
l[output] value of the domain index
xi[output] value of $\xi$

Implements Lorene::Map_radial.

Definition at line 421 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.

◆ val_r()

double Lorene::Map_et::val_r ( int  l,
double  xi,
double  theta,
double  pphi 
) const
virtual

Returns the value of the radial coordinate r for a given $(\xi, \theta', \phi')$ in a given domain.

Parameters
l[input] index of the domain
xi[input] value of $\xi$
theta[input] value of $\theta'$
pphi[input] value of $\phi'$
Returns
value of $r=R_l(\xi, \theta', \phi')$

Implements Lorene::Map.

Definition at line 95 of file map_et_radius.C.

References ff, Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and Lorene::Valeur::val_point().

◆ val_r_jk()

double Lorene::Map_et::val_r_jk ( int  l,
double  xi,
int  j,
int  k 
) const
virtual

Returns the value of the radial coordinate r for a given $\xi$ and a given collocation point in $(\theta', \phi')$ in a given domain.

Parameters
l[input] index of the domain
xi[input] value of $\xi$
j[input] index of the collocation point in $\theta'$
k[input] index of the collocation point in $\phi'$
Returns
value of $r=R_l(\xi, {\theta'}_j, {\phi'}_k)$

Implements Lorene::Map_radial.

Definition at line 370 of file map_et_radius.C.

References ff, Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.

Member Data Documentation

◆ aa

Tbl** Lorene::Map_et::aa
private

Array (size: mg->nzone ) of Tbl which stores the values of $A(\xi)$ in each domain.

Definition at line 2796 of file map.h.

◆ aasx

Tbl Lorene::Map_et::aasx
private

Values at the nr collocation points of $A(\xi)/\xi$ in the nucleus.

Definition at line 2809 of file map.h.

◆ aasx2

Tbl Lorene::Map_et::aasx2
private

Values at the nr collocation points of $A(\xi)/\xi^2$ in the nucleus.

Definition at line 2812 of file map.h.

◆ alpha

double* Lorene::Map_et::alpha
private

Array (size: mg->nzone ) of the values of $\alpha$ in each domain.

Definition at line 2789 of file map.h.

◆ bb

Tbl** Lorene::Map_et::bb
private

Array (size: mg->nzone ) of Tbl which stores the values of $B(\xi)$ in each domain.

Definition at line 2827 of file map.h.

◆ bbsx

Tbl Lorene::Map_et::bbsx
private

Values at the nr collocation points of $B(\xi)/\xi$ in the nucleus.

Definition at line 2840 of file map.h.

◆ bbsx2

Tbl Lorene::Map_et::bbsx2
private

Values at the nr collocation points of $B(\xi)/\xi^2$ in the nucleus.

Definition at line 2843 of file map.h.

◆ beta

double* Lorene::Map_et::beta
private

Array (size: mg->nzone ) of the values of $\beta$ in each domain.

Definition at line 2791 of file map.h.

◆ bvect_cart

Base_vect_cart Lorene::Map::bvect_cart
protectedinherited

Cartesian basis $(\partial/\partial x,\partial/\partial y,\partial/\partial z)$ associated with the coordinates (x,y,z) of the mapping, i.e.

the Cartesian coordinates related to $(r, \theta, \phi)$ by means of usual formulae.

Definition at line 715 of file map.h.

◆ bvect_spher

Base_vect_spher Lorene::Map::bvect_spher
protectedinherited

Orthonormal vectorial basis $(\partial/\partial r,1/r\partial/\partial \theta, 1/(r\sin\theta)\partial/\partial \phi)$ associated with the coordinates $(r, \theta, \phi)$ of the mapping.

Definition at line 707 of file map.h.

◆ cosp

Coord Lorene::Map::cosp
inherited

$\cos\phi$

Definition at line 742 of file map.h.

◆ cost

Coord Lorene::Map::cost
inherited

$\cos\theta$

Definition at line 740 of file map.h.

◆ d2rdtdx

Coord Lorene::Map_radial::d2rdtdx
inherited

$\partial^2 R/\partial\xi\partial\theta'$ in the nucleus and in the non-compactified shells; \ $-\partial^2 U/\partial\xi\partial\theta'$ in the compactified outer domain.

Definition at line 1661 of file map.h.

◆ d2rdx2

Coord Lorene::Map_radial::d2rdx2
inherited

$\partial^2 R/\partial\xi^2$ in the nucleus and in the non-compactified shells; \ $-\partial^2 U/\partial\xi^2 $ in the compactified outer domain.

Definition at line 1640 of file map.h.

◆ daa

Tbl** Lorene::Map_et::daa
private

Array (size: mg->nzone ) of Tbl which stores the values of $A'(\xi)$ in each domain.

Definition at line 2801 of file map.h.

◆ dbb

Tbl** Lorene::Map_et::dbb
private

Array (size: mg->nzone ) of Tbl which stores the values of $B'(\xi)$ in each domain.

Definition at line 2832 of file map.h.

◆ ddaa

Tbl** Lorene::Map_et::ddaa
private

Array (size: mg->nzone ) of Tbl which stores the values of $A''(\xi)$ in each domain.

Definition at line 2806 of file map.h.

◆ ddbb

Tbl** Lorene::Map_et::ddbb
private

Array (size: mg->nzone ) of Tbl which stores the values of $B''(\xi)$ in each domain.

Definition at line 2837 of file map.h.

◆ drdt

Coord Lorene::Map_radial::drdt
inherited

$\partial R/\partial\theta'$ in the nucleus and in the non-compactified shells; \ $-\partial U/\partial\theta'$ in the compactified external domain (CED).

Definition at line 1589 of file map.h.

◆ dxdr

Coord Lorene::Map_radial::dxdr
inherited

$1/(\partial R/\partial\xi) = \partial \xi /\partial r$ in the nucleus and in the non-compactified shells; \ $-1/(\partial U/\partial\xi) = - \partial \xi /\partial u$ in the compactified outer domain.

Definition at line 1581 of file map.h.

◆ ff

Valeur Lorene::Map_et::ff
private

Values of the function $F(\theta', \phi')$ at the nt*np angular collocation points in each domain.

The Valeur ff is defined on the multi-grid mg->g_angu (cf. class Mg3d ).

Definition at line 2850 of file map.h.

◆ gg

Valeur Lorene::Map_et::gg
private

Values of the function $G(\theta', \phi')$ at the nt*np angular collocation points in each domain.

The Valeur gg is defined on the multi-grid mg->g_angu (cf. class Mg3d ).

Definition at line 2857 of file map.h.

◆ lapr_tp

Coord Lorene::Map_radial::lapr_tp
inherited

$1/R^2 \times [ 1/\sin(\theta)\times \partial /\partial\theta' (\sin\theta \partial R /\partial\theta') + 1/\sin^2\theta \partial^2 R /\partial{\varphi'}^2] $ in the nucleus and in the non-compactified shells; \ $- 1/U^2 \times [ 1/\sin(\theta)\times \partial /\partial\theta' (\sin\theta \partial U /\partial\theta') + 1/\sin^2\theta \partial^2 U /\partial{\varphi'}^2] $ in the compactified outer domain.

Definition at line 1652 of file map.h.

◆ mg

const Mg3d* Lorene::Map::mg
protectedinherited

Pointer on the multi-grid Mgd3 on which this is defined.

Definition at line 694 of file map.h.

◆ ori_x

double Lorene::Map::ori_x
protectedinherited

Absolute coordinate x of the origin.

Definition at line 696 of file map.h.

◆ ori_y

double Lorene::Map::ori_y
protectedinherited

Absolute coordinate y of the origin.

Definition at line 697 of file map.h.

◆ ori_z

double Lorene::Map::ori_z
protectedinherited

Absolute coordinate z of the origin.

Definition at line 698 of file map.h.

◆ p_cmp_zero

Cmp* Lorene::Map::p_cmp_zero
protectedinherited

The null Cmp.

To be used by the Tenseur class when necessary to return a null Cmp .

Definition at line 731 of file map.h.

◆ p_flat_met_cart

Metric_flat* Lorene::Map::p_flat_met_cart
mutableprotectedinherited

Pointer onto the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart.

Definition at line 725 of file map.h.

◆ p_flat_met_spher

Metric_flat* Lorene::Map::p_flat_met_spher
mutableprotectedinherited

Pointer onto the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher.

Definition at line 720 of file map.h.

◆ p_mp_angu

Map_af* Lorene::Map::p_mp_angu
mutableprotectedinherited

Pointer on the "angular" mapping.

Definition at line 733 of file map.h.

◆ phi

Coord Lorene::Map::phi
inherited

$\phi$ coordinate centered on the grid

Definition at line 738 of file map.h.

◆ r

Coord Lorene::Map::r
inherited

r coordinate centered on the grid

Definition at line 736 of file map.h.

◆ rot_phi

double Lorene::Map::rot_phi
protectedinherited

Angle between the x –axis and X –axis.

Definition at line 699 of file map.h.

◆ rsx2drdx

Coord Lorene::Map_et::rsx2drdx

$[ R/ (\alpha \xi + \beta) ]^2 (\partial R/\partial \xi) / \alpha$ in the nucleus and the shells; \ $\partial U/\partial \xi / \alpha$ in the outermost compactified domain.

Definition at line 2872 of file map.h.

◆ rsxdxdr

Coord Lorene::Map_et::rsxdxdr

$1/(\partial R/\partial \xi) R/\xi$ in the nucleus; \ $1/(\partial R/\partial \xi) R/(\xi + \beta/\alpha)$ in the shells; \ $1/(\partial U/\partial \xi) U/(\xi-1)$ in the outermost compactified domain.

Definition at line 2865 of file map.h.

◆ sinp

Coord Lorene::Map::sinp
inherited

$\sin\phi$

Definition at line 741 of file map.h.

◆ sint

Coord Lorene::Map::sint
inherited

$\sin\theta$

Definition at line 739 of file map.h.

◆ sr2d2rdt2

Coord Lorene::Map_radial::sr2d2rdt2
inherited

$1/R^2 \partial^2 R/\partial{\theta'}^2$ in the nucleus and in the non-compactified shells; \ $-1/U^2 \partial^2 U/\partial{\theta'}^2$ in the compactified outer domain.

Definition at line 1678 of file map.h.

◆ sr2drdt

Coord Lorene::Map_radial::sr2drdt
inherited

$1/R^2 \times (\partial R/\partial\theta')$ in the nucleus and in the non-compactified shells; \ $-1/U^2 \times (\partial U/\partial\theta')$ in the compactified outer domain.

Definition at line 1621 of file map.h.

◆ sr2stdrdp

Coord Lorene::Map_radial::sr2stdrdp
inherited

$1/(R^2\sin\theta) \times (\partial R/\partial\varphi')$ in the nucleus and in the non-compactified shells; \ $-1/(U^2\sin\theta) \times (\partial U/\partial\varphi')$ in the compactified outer domain.

Definition at line 1629 of file map.h.

◆ srdrdt

Coord Lorene::Map_radial::srdrdt
inherited

$1/R \times (\partial R/\partial\theta')$ in the nucleus and in the non-compactified shells; \ $-1/U \times (\partial U/\partial\theta)$ in the compactified outer domain.

Definition at line 1605 of file map.h.

◆ srstdrdp

Coord Lorene::Map_radial::srstdrdp
inherited

$1/(R\sin\theta) \times (\partial R/\partial\varphi')$ in the nucleus and in the non-compactified shells; \ $-1/(U\sin\theta) \times (\partial U/\partial\varphi')$ in the compactified outer domain.

Definition at line 1613 of file map.h.

◆ sstd2rdpdx

Coord Lorene::Map_radial::sstd2rdpdx
inherited

$1/\sin\theta \times \partial^2 R/\partial\xi\partial\varphi'$ in the nucleus and in the non-compactified shells; \ $-1/\sin\theta \times \partial^2 U/\partial\xi\partial\varphi' $ in the compactified outer domain.

Definition at line 1669 of file map.h.

◆ stdrdp

Coord Lorene::Map_radial::stdrdp
inherited

${1\over\sin\theta} \partial R/\partial\varphi'$ in the nucleus and in the non-compactified shells; \ $-{1\over\sin\theta}\partial U/\partial\varphi'$ in the compactified external domain (CED).

Definition at line 1597 of file map.h.

◆ tet

Coord Lorene::Map::tet
inherited

$\theta$ coordinate centered on the grid

Definition at line 737 of file map.h.

◆ x

Coord Lorene::Map::x
inherited

x coordinate centered on the grid

Definition at line 744 of file map.h.

◆ xa

Coord Lorene::Map::xa
inherited

Absolute x coordinate.

Definition at line 748 of file map.h.

◆ xsr

Coord Lorene::Map_radial::xsr
inherited

$\xi/R$ in the nucleus; \ 1/R in the non-compactified shells; \ $(\xi-1)/U$ in the compactified outer domain.

Definition at line 1570 of file map.h.

◆ y

Coord Lorene::Map::y
inherited

y coordinate centered on the grid

Definition at line 745 of file map.h.

◆ ya

Coord Lorene::Map::ya
inherited

Absolute y coordinate.

Definition at line 749 of file map.h.

◆ z

Coord Lorene::Map::z
inherited

z coordinate centered on the grid

Definition at line 746 of file map.h.

◆ za

Coord Lorene::Map::za
inherited

Absolute z coordinate.

Definition at line 750 of file map.h.

◆ zaasx

Tbl Lorene::Map_et::zaasx
private

Values at the nr collocation points of $A(\xi)/(\xi-1)$ in the outermost compactified domain.

Definition at line 2817 of file map.h.

◆ zaasx2

Tbl Lorene::Map_et::zaasx2
private

Values at the nr collocation points of $A(\xi)/(\xi-1)^2$ in the outermost compactified domain.

Definition at line 2822 of file map.h.


The documentation for this class was generated from the following files: