LORENE
Lorene::Map_eps Class Reference
Inheritance diagram for Lorene::Map_eps:
Lorene::Map_radial Lorene::Map

Public Member Functions

 Map_eps (const Mg3d &mgrille, const double *r_limits)
 Standard Constructor. More...
 
 Map_eps (const Mg3d &mgrille, double a, double b, double c)
 Standard Constructor with Tbl. More...
 
 Map_eps (const Map_eps &)
 Copy constructor. More...
 
 Map_eps (const Mg3d &, const string &)
 Constructor from a formatted file. More...
 
 Map_eps (const Mg3d &, FILE *)
 Constructor from a file (see sauve(FILE*) ) More...
 
virtual ~Map_eps ()
 Destructor. More...
 
virtual const Map_afmp_angu (int) const
 Returns the "angular" mapping for the outside of domain l_zone. More...
 
virtual void operator= (const Map_af &)
 Assignment to an affine mapping. More...
 
virtual void operator= (const Map_eps &)
 Assignment to another Map_eps. More...
 
const Valeurget_alpha () const
 Returns the reference on the Tbl alpha. More...
 
const Valeurget_beta () const
 
double get_aa () const
 
double get_bb () const
 
double get_cc () const
 
void set_alpha (const Tbl &alpha0, int l)
 Modifies the value of $\alpha$ in domain no. l. More...
 
void set_beta (const Tbl &beta0, int l)
 
void set_alpha (const Valeur &alpha0)
 
void set_beta (const Valeur &beta0)
 
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 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 bool operator== (const Map &) const
 Comparison operator (egality) More...
 
virtual void dsdr (const Scalar &ci, Scalar &resu) const
 Computes $\partial/ \partial r$ of a Scalar. More...
 
virtual void srstdsdp (const Scalar &, Scalar &) const
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Scalar. More...
 
virtual void srdsdt (const Scalar &, Scalar &) const
 Computes $1/r \partial/ \partial \theta$ of a Scalar. More...
 
virtual void dsdt (const Scalar &, Scalar &) const
 Computes $\partial/ \partial \theta$ of a Scalar . More...
 
virtual void stdsdp (const Scalar &, Scalar &) const
 Computes $1/\sin\theta \partial/ \partial \varphi$ of a Scalar . More...
 
virtual void sauve (FILE *) const
 Save in a file. 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 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 ()
 Assignment of the building functions to the member Coords. More...
 
virtual ostream & operator>> (ostream &) const
 Operator >> More...
 
virtual void homothetie (double)
 Sets a new radial scale. More...
 
virtual void resize (int, double)
 < Not implemented More...
 
virtual void adapt (const Cmp &, const Param &, int)
 < Not implemented More...
 
virtual void dsdr (const Cmp &, Cmp &) const
 < Not implemented More...
 
virtual void dsdxi (const Cmp &, Cmp &) const
 < Not implemented More...
 
virtual void dsdxi (const Scalar &, Scalar &) const
 < Not implemented More...
 
virtual void dsdradial (const Scalar &uu, Scalar &resu) const
 < Not implemented More...
 
virtual void srdsdt (const Cmp &, Cmp &) const
 < Not implemented More...
 
virtual void srstdsdp (const Cmp &, Cmp &) const
 < Not implemented More...
 
virtual void laplacien (const Scalar &, int, Scalar &) const
 < Not implemented More...
 
virtual void laplacien (const Cmp &, int, Cmp &) const
 < Not implemented More...
 
virtual void lapang (const Scalar &, Scalar &) const
 < Not implemented More...
 
virtual void primr (const Scalar &, Scalar &, bool) const
 < Not implemented More...
 
virtual Tblintegrale (const Cmp &) const
 < Not implemented More...
 
virtual void poisson (const Cmp &, Param &, Cmp &) const
 < Not implemented More...
 
virtual void poisson_tau (const Cmp &, Param &, Cmp &) const
 < Not implemented More...
 
virtual void poisson_falloff (const Cmp &, Param &, Cmp &, int) const
 < Not implemented More...
 
virtual void poisson_ylm (const Cmp &, Param &, Cmp &, int, double *) const
 < Not implemented More...
 
virtual void poisson_regular (const Cmp &, int, int, double, Param &, Cmp &, Cmp &, Cmp &, Tenseur &, Cmp &, Cmp &) const
 < Not implemented More...
 
virtual void poisson_angu (const Scalar &, Param &, Scalar &, double=0) const
 < Not implemented More...
 
virtual void poisson_angu (const Cmp &, Param &, Cmp &, double=0) const
 < Not implemented More...
 
virtual Paramdonne_para_poisson_vect (Param &, int) const
 < Not implemented More...
 
virtual void poisson_frontiere (const Cmp &, const Valeur &, int, int, Cmp &, double=0., double=0.) const
 < Not implemented More...
 
virtual void poisson_frontiere_double (const Cmp &, const Valeur &, const Valeur &, int, Cmp &) const
 < Not implemented More...
 
virtual void poisson_interne (const Cmp &, const Valeur &, Param &, Cmp &) const
 < Not implemented More...
 
virtual void poisson2d (const Cmp &, const Cmp &, Param &, Cmp &) const
 < Not implemented More...
 
virtual void dalembert (Param &, Scalar &, const Scalar &, const Scalar &, const Scalar &) const
 < Not implemented More...
 

Private Attributes

double aa
 Array (size: mg->nzone*Nt*Np ) of the values of $\alpha$ in each domain. More...
 
double bb
 
double cc
 
Valeur alpha
 
Valeur beta
 

Friends

Mtblmap_eps_fait_r (const Map *)
 < Not implemented More...
 
Mtblmap_eps_fait_tet (const Map *)
 
Mtblmap_eps_fait_phi (const Map *)
 
Mtblmap_eps_fait_sint (const Map *)
 
Mtblmap_eps_fait_cost (const Map *)
 
Mtblmap_eps_fait_sinp (const Map *)
 
Mtblmap_eps_fait_cosp (const Map *)
 
Mtblmap_eps_fait_x (const Map *)
 
Mtblmap_eps_fait_y (const Map *)
 
Mtblmap_eps_fait_z (const Map *)
 
Mtblmap_eps_fait_xa (const Map *)
 
Mtblmap_eps_fait_ya (const Map *)
 
Mtblmap_eps_fait_za (const Map *)
 
Mtblmap_eps_fait_xsr (const Map *)
 
Mtblmap_eps_fait_dxdr (const Map *)
 
Mtblmap_eps_fait_drdt (const Map *)
 
Mtblmap_eps_fait_stdrdp (const Map *)
 
Mtblmap_eps_fait_srdrdt (const Map *)
 
Mtblmap_eps_fait_srstdrdp (const Map *)
 
Mtblmap_eps_fait_sr2drdt (const Map *)
 
Mtblmap_eps_fait_sr2stdrdp (const Map *)
 
Mtblmap_eps_fait_d2rdx2 (const Map *)
 
Mtblmap_eps_fait_lapr_tp (const Map *)
 
Mtblmap_eps_fait_d2rdtdx (const Map *)
 
Mtblmap_eps_fait_sstd2rdpdx (const Map *)
 
Mtblmap_eps_fait_sr2d2rdt2 (const Map *)
 

Detailed Description

Definition at line 4210 of file map.h.

Constructor & Destructor Documentation

◆ Map_eps() [1/5]

Lorene::Map_eps::Map_eps ( 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: Nt*Np ) (only 1 domain at the moment) of the value of r at the boundaries of the various domains :
  • r_limits[l,k,j] : outer boundary of the domain no. l in the angular direction j,k

Definition at line 52 of file map_eps.C.

References Lorene::c_est_pas_fait().

◆ Map_eps() [2/5]

Lorene::Map_eps::Map_eps ( const Mg3d mgrille,
double  a,
double  b,
double  c 
)

Standard Constructor with Tbl.

Parameters
mgrille[input] Multi-domain grid on which the mapping is defined
r_limits[input] Array (size: Nt*Np ) (only 1 domain at the moment) of the value of r at the boundaries of the various domains :
  • r_limits[l,k,j] : outer boundary of the domain no. l in the angular direction j,k

Definition at line 59 of file map_eps.C.

References aa, Lorene::Valeur::annule_hard(), Lorene::Mg3d::get_grille3d(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and set_coord().

◆ Map_eps() [3/5]

Lorene::Map_eps::Map_eps ( const Map_eps mp)

Copy constructor.

Definition at line 115 of file map_eps.C.

References set_coord().

◆ Map_eps() [4/5]

Lorene::Map_eps::Map_eps ( const Mg3d ,
const string &   
)

Constructor from a formatted file.

◆ Map_eps() [5/5]

Lorene::Map_eps::Map_eps ( const Mg3d mgi,
FILE *  fd 
)

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

Definition at line 124 of file map_eps.C.

References set_coord().

◆ ~Map_eps()

Lorene::Map_eps::~Map_eps ( )
virtual

Destructor.

Definition at line 134 of file map_eps.C.

Member Function Documentation

◆ adapt()

void Lorene::Map_eps::adapt ( const Cmp ,
const Param ,
int   
)
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 51 of file map_eps_pas_fait.C.

◆ 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 819 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_eps::dalembert ( Param ,
Scalar ,
const Scalar ,
const Scalar ,
const Scalar  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 134 of file map_eps_pas_fait.C.

◆ 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_eps::donne_para_poisson_vect ( Param ,
int   
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 113 of file map_eps_pas_fait.C.

◆ dsdr() [1/2]

void Lorene::Map_eps::dsdr ( const Scalar ci,
Scalar resu 
) const
virtual

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

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 39 of file map_eps_deriv.C.

References Lorene::Scalar::get_etat().

◆ dsdr() [2/2]

void Lorene::Map_eps::dsdr ( const Cmp ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 55 of file map_eps_pas_fait.C.

◆ dsdradial()

void Lorene::Map_eps::dsdradial ( const Scalar uu,
Scalar resu 
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 154 of file map_eps_pas_fait.C.

◆ dsdt()

void Lorene::Map_eps::dsdt ( const Scalar ci,
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 231 of file map_eps_deriv.C.

References Lorene::Scalar::get_etat().

◆ dsdxi() [1/2]

void Lorene::Map_eps::dsdxi ( const Cmp ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 59 of file map_eps_pas_fait.C.

◆ dsdxi() [2/2]

void Lorene::Map_eps::dsdxi ( const Scalar ,
Scalar  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 63 of file map_eps_pas_fait.C.

◆ 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 Valeur & Lorene::Map_eps::get_alpha ( ) const

Returns the reference on the Tbl alpha.

Definition at line 243 of file map_eps.C.

◆ 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 803 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 795 of file map.h.

References Lorene::Map::bvect_spher.

◆ get_mg()

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

Gives the Mg3d on which the mapping is defined.

Definition at line 777 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 780 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 782 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 784 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 787 of file map.h.

References Lorene::Map::rot_phi.

◆ homothetie()

void Lorene::Map_eps::homothetie ( double  lambda)
privatevirtual

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 43 of file map_eps_pas_fait.C.

◆ 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_eps::integrale ( const Cmp ) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 87 of file map_eps_pas_fait.C.

◆ lapang()

void Lorene::Map_eps::lapang ( const Scalar ,
Scalar  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 83 of file map_eps_pas_fait.C.

◆ laplacien() [1/2]

void Lorene::Map_eps::laplacien ( const Scalar ,
int  ,
Scalar  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 79 of file map_eps_pas_fait.C.

◆ laplacien() [2/2]

void Lorene::Map_eps::laplacien ( const Cmp ,
int  ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 75 of file map_eps_pas_fait.C.

◆ mp_angu()

const Map_af & Lorene::Map_eps::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 138 of file map_eps_pas_fait.C.

◆ 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]

void Lorene::Map_eps::operator= ( const Map_af mpi)
virtual

Assignment to an affine mapping.

Implements Lorene::Map_radial.

Definition at line 161 of file map_eps.C.

References Lorene::c_est_pas_fait().

◆ operator=() [2/2]

void Lorene::Map_eps::operator= ( const Map_eps mpi)
virtual

◆ operator==()

◆ operator>>()

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

◆ poisson()

void Lorene::Map_eps::poisson ( const Cmp ,
Param ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 92 of file map_eps_pas_fait.C.

◆ poisson2d()

void Lorene::Map_eps::poisson2d ( const Cmp ,
const Cmp ,
Param ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 130 of file map_eps_pas_fait.C.

◆ poisson_angu() [1/2]

void Lorene::Map_eps::poisson_angu ( const Scalar ,
Param ,
Scalar ,
double  = 0 
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 105 of file map_eps_pas_fait.C.

◆ poisson_angu() [2/2]

void Lorene::Map_eps::poisson_angu ( const Cmp ,
Param ,
Cmp ,
double  = 0 
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 109 of file map_eps_pas_fait.C.

◆ 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_falloff()

void Lorene::Map_eps::poisson_falloff ( const Cmp ,
Param ,
Cmp ,
int   
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 146 of file map_eps_pas_fait.C.

◆ poisson_frontiere()

void Lorene::Map_eps::poisson_frontiere ( const Cmp ,
const Valeur ,
int  ,
int  ,
Cmp ,
double  = 0.,
double  = 0. 
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 118 of file map_eps_pas_fait.C.

◆ poisson_frontiere_double()

void Lorene::Map_eps::poisson_frontiere_double ( const Cmp ,
const Valeur ,
const Valeur ,
int  ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 122 of file map_eps_pas_fait.C.

◆ poisson_interne()

void Lorene::Map_eps::poisson_interne ( const Cmp ,
const Valeur ,
Param ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 126 of file map_eps_pas_fait.C.

◆ poisson_regular()

void Lorene::Map_eps::poisson_regular ( const Cmp ,
int  ,
int  ,
double  ,
Param ,
Cmp ,
Cmp ,
Cmp ,
Tenseur ,
Cmp ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 100 of file map_eps_pas_fait.C.

◆ poisson_tau()

void Lorene::Map_eps::poisson_tau ( const Cmp ,
Param ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 96 of file map_eps_pas_fait.C.

◆ poisson_ylm()

void Lorene::Map_eps::poisson_ylm ( const Cmp ,
Param ,
Cmp ,
int  ,
double *   
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 150 of file map_eps_pas_fait.C.

◆ primr()

void Lorene::Map_eps::primr ( const Scalar ,
Scalar ,
bool   
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 142 of file map_eps_pas_fait.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()

◆ resize()

void Lorene::Map_eps::resize ( int  ,
double   
)
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 47 of file map_eps_pas_fait.C.

◆ sauve()

void Lorene::Map_eps::sauve ( FILE *  fd) const
virtual

Save in a file.

Reimplemented from Lorene::Map_radial.

Definition at line 254 of file map_eps.C.

References Lorene::Valeur::sauve(), and Lorene::Map_radial::sauve().

◆ set_alpha()

void Lorene::Map_eps::set_alpha ( const Tbl alpha0,
int  l 
)

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

Definition at line 314 of file map_eps.C.

References Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Map::mg, Lorene::Map_radial::reset_coord(), and Lorene::Valeur::set().

◆ set_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_eps::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 271 of file map_eps_deriv.C.

References Lorene::Scalar::get_etat().

◆ srdsdt() [2/2]

void Lorene::Map_eps::srdsdt ( const Cmp ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 67 of file map_eps_pas_fait.C.

◆ srstdsdp() [1/2]

void Lorene::Map_eps::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 140 of file map_eps_deriv.C.

References Lorene::Scalar::get_etat().

◆ srstdsdp() [2/2]

void Lorene::Map_eps::srstdsdp ( const Cmp ,
Cmp  
) const
privatevirtual

< Not implemented

Implements Lorene::Map.

Definition at line 71 of file map_eps_pas_fait.C.

◆ stdsdp()

void Lorene::Map_eps::stdsdp ( const Scalar ci,
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 95 of file map_eps_deriv.C.

References Lorene::Scalar::get_etat().

◆ val_lx() [1/2]

void Lorene::Map_eps::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 82 of file map_eps_radius.C.

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

◆ val_lx() [2/2]

void Lorene::Map_eps::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)$.

Parameters
rr[input] value of r
theta[input] value of $\theta$
pphi[input] value of $\phi$
par[] unused by the Map_eps version
l[output] value of the domain index
xi[output] value of $\xi$

Implements Lorene::Map.

Definition at line 159 of file map_eps_radius.C.

References val_lx().

◆ val_lx_jk()

void Lorene::Map_eps::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[] unused by the Map_eps version
l[output] value of the domain index
xi[output] value of $\xi$

Implements Lorene::Map_radial.

Definition at line 182 of file map_eps_radius.C.

References Lorene::Map::phi, Lorene::Map::tet, and val_lx().

◆ val_r()

double Lorene::Map_eps::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 43 of file map_eps_radius.C.

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

◆ val_r_jk()

double Lorene::Map_eps::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 172 of file map_eps_radius.C.

References Lorene::Map::phi, Lorene::Map::tet, and val_r().

Friends And Related Function Documentation

◆ map_eps_fait_r

Mtbl* map_eps_fait_r ( const Map cvi)
friend

< Not implemented

Definition at line 43 of file map_eps_fait.C.

Member Data Documentation

◆ aa

double Lorene::Map_eps::aa
private

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

Definition at line 4217 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 709 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 701 of file map.h.

◆ cosp

Coord Lorene::Map::cosp
inherited

$\cos\phi$

Definition at line 736 of file map.h.

◆ cost

Coord Lorene::Map::cost
inherited

$\cos\theta$

Definition at line 734 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 1655 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 1634 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 1583 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 1575 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 1646 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 688 of file map.h.

◆ ori_x

double Lorene::Map::ori_x
protectedinherited

Absolute coordinate x of the origin.

Definition at line 690 of file map.h.

◆ ori_y

double Lorene::Map::ori_y
protectedinherited

Absolute coordinate y of the origin.

Definition at line 691 of file map.h.

◆ ori_z

double Lorene::Map::ori_z
protectedinherited

Absolute coordinate z of the origin.

Definition at line 692 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 725 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 719 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 714 of file map.h.

◆ p_mp_angu

Map_af* Lorene::Map::p_mp_angu
mutableprotectedinherited

Pointer on the "angular" mapping.

Definition at line 727 of file map.h.

◆ phi

Coord Lorene::Map::phi
inherited

$\phi$ coordinate centered on the grid

Definition at line 732 of file map.h.

◆ r

Coord Lorene::Map::r
inherited

r coordinate centered on the grid

Definition at line 730 of file map.h.

◆ rot_phi

double Lorene::Map::rot_phi
protectedinherited

Angle between the x –axis and X –axis.

Definition at line 693 of file map.h.

◆ sinp

Coord Lorene::Map::sinp
inherited

$\sin\phi$

Definition at line 735 of file map.h.

◆ sint

Coord Lorene::Map::sint
inherited

$\sin\theta$

Definition at line 733 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 1672 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 1615 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 1623 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 1599 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 1607 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 1663 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 1591 of file map.h.

◆ tet

Coord Lorene::Map::tet
inherited

$\theta$ coordinate centered on the grid

Definition at line 731 of file map.h.

◆ x

Coord Lorene::Map::x
inherited

x coordinate centered on the grid

Definition at line 738 of file map.h.

◆ xa

Coord Lorene::Map::xa
inherited

Absolute x coordinate.

Definition at line 742 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 1564 of file map.h.

◆ y

Coord Lorene::Map::y
inherited

y coordinate centered on the grid

Definition at line 739 of file map.h.

◆ ya

Coord Lorene::Map::ya
inherited

Absolute y coordinate.

Definition at line 743 of file map.h.

◆ z

Coord Lorene::Map::z
inherited

z coordinate centered on the grid

Definition at line 740 of file map.h.

◆ za

Coord Lorene::Map::za
inherited

Absolute z coordinate.

Definition at line 744 of file map.h.


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