LORENE

Base class for coordinate mappings. More...

#include <map.h>

Inheritance diagram for Lorene::Map:
Lorene::Map_radial Lorene::Map_af Lorene::Map_eps Lorene::Map_et Lorene::Map_log Lorene::Map_star

Public Member Functions

virtual ~Map ()
 Destructor. More...
 
virtual void sauve (FILE *) const
 Save in a file. 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...
 
virtual const Map_afmp_angu (int) const =0
 Returns the "angular" mapping for the outside of domain l_zone. More...
 
void convert_absolute (double xx, double yy, double zz, double &rr, double &theta, double &pphi) const
 Determines the coordinates $(r,\theta,\phi)$ corresponding to given absolute Cartesian coordinates (X,Y,Z). More...
 
virtual double val_r (int l, double xi, double theta, double pphi) const =0
 Returns the value of the radial coordinate r for a given $(\xi, \theta', \phi')$ in a given domain. More...
 
virtual void val_lx (double rr, double theta, double pphi, int &l, double &xi) const =0
 Computes the domain index l and the value of $\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 =0
 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 =0
 Comparison operator (egality) More...
 
void set_ori (double xa0, double ya0, double za0)
 Sets a new origin. More...
 
void set_rot_phi (double phi0)
 Sets a new rotation angle. More...
 
virtual void homothetie (double lambda)=0
 Sets a new radial scale. More...
 
virtual void resize (int l, double lambda)=0
 Rescales the outer boundary of one domain. More...
 
virtual void operator= (const Map_af &)=0
 Assignment to an affine mapping. More...
 
virtual void adapt (const Cmp &ent, const Param &par, int nbr=0)=0
 Adaptation of the mapping to a given scalar field. More...
 
virtual void reevaluate (const Map *mp_prev, int nzet, Cmp &uu) const =0
 Recomputes the values of a Cmp at the collocation points after a change in the mapping. More...
 
virtual void reevaluate_symy (const Map *mp_prev, int nzet, Cmp &uu) const =0
 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 =0
 Recomputes the values of a Scalar at the collocation points after a change in the mapping. More...
 
virtual void reevaluate_symy (const Map *mp_prev, int nzet, Scalar &uu) const =0
 Recomputes the values of a Scalar at the collocation points after a change in the mapping. More...
 
virtual void dsdxi (const Cmp &ci, Cmp &resu) const =0
 Computes $\partial/ \partial \xi$ of a Cmp . More...
 
virtual void dsdr (const Cmp &ci, Cmp &resu) const =0
 Computes $\partial/ \partial r$ of a Cmp . More...
 
virtual void srdsdt (const Cmp &ci, Cmp &resu) const =0
 Computes $1/r \partial/ \partial \theta$ of a Cmp . More...
 
virtual void srstdsdp (const Cmp &ci, Cmp &resu) const =0
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Cmp . More...
 
virtual void dsdxi (const Scalar &uu, Scalar &resu) const =0
 Computes $\partial/ \partial xi$ of a Scalar . More...
 
virtual void dsdr (const Scalar &uu, Scalar &resu) const =0
 Computes $\partial/ \partial r$ of a Scalar . More...
 
virtual void dsdradial (const Scalar &uu, Scalar &resu) const =0
 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 =0
 Computes $1/r \partial/ \partial \theta$ of a Scalar . More...
 
virtual void srstdsdp (const Scalar &uu, Scalar &resu) const =0
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Scalar . More...
 
virtual void dsdt (const Scalar &uu, Scalar &resu) const =0
 Computes $\partial/ \partial \theta$ of a Scalar . More...
 
virtual void stdsdp (const Scalar &uu, Scalar &resu) const =0
 Computes $1/\sin\theta \partial/ \partial \varphi$ of a Scalar . More...
 
virtual void laplacien (const Scalar &uu, int zec_mult_r, Scalar &lap) const =0
 Computes the Laplacian of a scalar field. More...
 
virtual void laplacien (const Cmp &uu, int zec_mult_r, Cmp &lap) const =0
 Computes the Laplacian of a scalar field (Cmp version). More...
 
virtual void lapang (const Scalar &uu, Scalar &lap) const =0
 Computes the angular Laplacian of a scalar field. More...
 
virtual void primr (const Scalar &uu, Scalar &resu, bool null_infty) const =0
 Computes the radial primitive which vanishes for $r\to \infty$. More...
 
virtual void mult_r (Scalar &uu) const =0
 Multiplication by r of a Scalar , the dzpuis
of uu is not changed. More...
 
virtual void mult_r (Cmp &ci) const =0
 Multiplication by r of a Cmp . More...
 
virtual void mult_r_zec (Scalar &) const =0
 Multiplication by r (in the compactified external domain only) of a Scalar. More...
 
virtual void mult_rsint (Scalar &) const =0
 Multiplication by $r\sin\theta$ of a Scalar. More...
 
virtual void div_rsint (Scalar &) const =0
 Division by $r\sin\theta$ of a Scalar. More...
 
virtual void div_r (Scalar &) const =0
 Division by r of a Scalar. More...
 
virtual void div_r_zec (Scalar &) const =0
 Division by r (in the compactified external domain only) of a Scalar. More...
 
virtual void mult_cost (Scalar &) const =0
 Multiplication by $\cos\theta$ of a Scalar. More...
 
virtual void div_cost (Scalar &) const =0
 Division by $\cos\theta$ of a Scalar. More...
 
virtual void mult_sint (Scalar &) const =0
 Multiplication by $\sin\theta$ of a Scalar. More...
 
virtual void div_sint (Scalar &) const =0
 Division by $\sin\theta$ of a Scalar. More...
 
virtual void div_tant (Scalar &) const =0
 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 =0
 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 =0
 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 =0
 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 =0
 Cmp version More...
 
virtual void comp_z_from_spherical (const Scalar &v_r, const Scalar &v_theta, Scalar &v_z) const =0
 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 =0
 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 =0
 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 =0
 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 =0
 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 =0
 Cmp version More...
 
virtual void comp_p_from_cartesian (const Scalar &v_x, const Scalar &v_y, Scalar &v_p) const =0
 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 =0
 Cmp version More...
 
virtual void dec_dzpuis (Scalar &) const =0
 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 =0
 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 =0
 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 =0
 Increases by 2 the value of dzpuis of a Scalar
and changes accordingly its values in the
compactified external domain (CED). More...
 
virtual Tblintegrale (const Cmp &) const =0
 Computes the integral over all space of a Cmp . More...
 
virtual void poisson (const Cmp &source, Param &par, Cmp &uu) const =0
 Computes the solution of a scalar Poisson equation. More...
 
virtual void poisson_tau (const Cmp &source, Param &par, Cmp &uu) const =0
 Computes the solution of a scalar Poisson equationwith a Tau method. More...
 
virtual void poisson_falloff (const Cmp &source, Param &par, Cmp &uu, int k_falloff) const =0
 
virtual void poisson_ylm (const Cmp &source, Param &par, Cmp &pot, int nylm, double *intvec) const =0
 
virtual void poisson_regular (const Cmp &source, int k_div, int nzet, double unsgam1, Param &par, Cmp &uu, Cmp &uu_regu, Cmp &uu_div, Tenseur &duu_div, Cmp &source_regu, Cmp &source_div) const =0
 Computes the solution of a scalar Poisson equation. More...
 
virtual void poisson_compact (const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const =0
 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 =0
 Resolution of the elliptic equation $ a \Delta\psi + {\bf b} \cdot \nabla \psi = \sigma$ in the case of a multidomain stellar interior. More...
 
virtual void poisson_angu (const Scalar &source, Param &par, Scalar &uu, double lambda=0) const =0
 Computes the solution of the generalized angular Poisson equation. More...
 
virtual void poisson_angu (const Cmp &source, Param &par, Cmp &uu, double lambda=0) const =0
 
virtual Paramdonne_para_poisson_vect (Param &para, int i) const =0
 Function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara . More...
 
virtual void poisson_frontiere (const Cmp &source, const Valeur &limite, int raccord, int num_front, Cmp &pot, double=0., double=0.) const =0
 Computes the solution of a Poisson equation from the domain num_front+1 . More...
 
virtual void poisson_frontiere_double (const Cmp &source, const Valeur &lim_func, const Valeur &lim_der, int num_zone, Cmp &pot) const =0
 
virtual void poisson_interne (const Cmp &source, const Valeur &limite, Param &par, Cmp &pot) const =0
 Computes the solution of a Poisson equation in the shell, imposing a boundary condition at the surface of the star. More...
 
virtual void poisson2d (const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu) const =0
 Computes the solution of a 2-D Poisson equation. More...
 
virtual void dalembert (Param &par, Scalar &fJp1, const Scalar &fJ, const Scalar &fJm1, const Scalar &source) const =0
 Performs one time-step integration of the d'Alembert scalar equation. More...
 

Public Attributes

Coord 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

 Map (const Mg3d &)
 Constructor from a multi-domain 3D grid. More...
 
 Map (const Map &)
 Copy constructor. More...
 
 Map (const Mg3d &, FILE *)
 Constructor from a file (see sauve(FILE* ) ) More...
 
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

virtual ostream & operator>> (ostream &) const =0
 Operator >> More...
 

Friends

ostream & operator<< (ostream &, const Map &)
 Operator <<. More...
 

Detailed Description

Base class for coordinate mappings.

()

This class is the basic class for describing the mapping between the grid coordinates $(\xi, \theta', \phi')$ and the physical coordinates $(r, \theta, \phi)$ [cf. Bonazzola, Gourgoulhon & Marck, Phys. Rev. D 58, 104020 (1998)]. The class Map is an abstract one: it cannot be instanciated. Specific implementation of coordinate mappings will be performed by derived classes of Map.

The class Map and its derived classes determine the methods for partial derivatives with respect to the physical coordinates, as well as resolution of basic partial differential equations (e.g. Poisson equations).

The mapping is defined with respect to some `‘absolute’' reference frame, whose Cartesian coordinates are denoted by (X,Y,Z). The coordinates (X, Y, Z) of center of the mapping (i.e. the point r =0) are given by the data members (ori_x,ori_y,ori_z). The Cartesian coordinate relative to the mapping (i.e. defined from $(r, \theta, \phi)$ by the usual formul $x=r\sin\theta\cos\phi, \ldots$) are denoted by (x,y,z). The planes (x,y) and (X,Y) are supposed to coincide, so that the relative orientation of the mapping with respect to the absolute reference frame is described by only one angle (the data member rot_phi).

Definition at line 688 of file map.h.

Constructor & Destructor Documentation

◆ Map() [1/3]

Lorene::Map::Map ( const Mg3d mgi)
explicitprotected

Constructor from a multi-domain 3D grid.

Definition at line 142 of file map.C.

References p_cmp_zero, and Lorene::Cmp::set_etat_zero().

◆ Map() [2/3]

Lorene::Map::Map ( const Map mp)
protected

Copy constructor.

Definition at line 160 of file map.C.

References p_cmp_zero, and Lorene::Cmp::set_etat_zero().

◆ Map() [3/3]

Lorene::Map::Map ( const Mg3d mgi,
FILE *  fd 
)
protected

◆ ~Map()

Lorene::Map::~Map ( )
virtual

Destructor.

Definition at line 216 of file map.C.

References p_cmp_zero, p_flat_met_cart, p_flat_met_spher, and p_mp_angu.

Member Function Documentation

◆ adapt()

virtual void Lorene::Map::adapt ( const Cmp ent,
const Param par,
int  nbr = 0 
)
pure virtual

Adaptation of the mapping to a given scalar field.

This is a virtual function: see the actual implementations in the derived classes for the meaning of the various parameters.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ cmp_zero()

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

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 p_cmp_zero.

◆ comp_p_from_cartesian() [1/2]

virtual void Lorene::Map::comp_p_from_cartesian ( const Scalar v_x,
const Scalar v_y,
Scalar v_p 
) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ comp_p_from_cartesian() [2/2]

virtual void Lorene::Map::comp_p_from_cartesian ( const Cmp v_x,
const Cmp v_y,
Cmp v_p 
) const
pure virtual

Cmp version

Implemented in Lorene::Map_radial.

◆ comp_r_from_cartesian() [1/2]

virtual void Lorene::Map::comp_r_from_cartesian ( const Scalar v_x,
const Scalar v_y,
const Scalar v_z,
Scalar v_r 
) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ comp_r_from_cartesian() [2/2]

virtual void Lorene::Map::comp_r_from_cartesian ( const Cmp v_x,
const Cmp v_y,
const Cmp v_z,
Cmp v_r 
) const
pure virtual

Cmp version

Implemented in Lorene::Map_radial.

◆ comp_t_from_cartesian() [1/2]

virtual void Lorene::Map::comp_t_from_cartesian ( const Scalar v_x,
const Scalar v_y,
const Scalar v_z,
Scalar v_t 
) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ comp_t_from_cartesian() [2/2]

virtual void Lorene::Map::comp_t_from_cartesian ( const Cmp v_x,
const Cmp v_y,
const Cmp v_z,
Cmp v_t 
) const
pure virtual

Cmp version

Implemented in Lorene::Map_radial.

◆ comp_x_from_spherical() [1/2]

virtual void Lorene::Map::comp_x_from_spherical ( const Scalar v_r,
const Scalar v_theta,
const Scalar v_phi,
Scalar v_x 
) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ comp_x_from_spherical() [2/2]

virtual void Lorene::Map::comp_x_from_spherical ( const Cmp v_r,
const Cmp v_theta,
const Cmp v_phi,
Cmp v_x 
) const
pure virtual

Cmp version

Implemented in Lorene::Map_radial.

◆ comp_y_from_spherical() [1/2]

virtual void Lorene::Map::comp_y_from_spherical ( const Scalar v_r,
const Scalar v_theta,
const Scalar v_phi,
Scalar v_y 
) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ comp_y_from_spherical() [2/2]

virtual void Lorene::Map::comp_y_from_spherical ( const Cmp v_r,
const Cmp v_theta,
const Cmp v_phi,
Cmp v_y 
) const
pure virtual

Cmp version

Implemented in Lorene::Map_radial.

◆ comp_z_from_spherical() [1/2]

virtual void Lorene::Map::comp_z_from_spherical ( const Scalar v_r,
const Scalar v_theta,
Scalar v_z 
) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ comp_z_from_spherical() [2/2]

virtual void Lorene::Map::comp_z_from_spherical ( const Cmp v_r,
const Cmp v_theta,
Cmp v_z 
) const
pure virtual

Cmp version

Implemented in Lorene::Map_radial.

◆ convert_absolute()

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

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 ori_x, ori_y, ori_z, rot_phi, and Lorene::sqrt().

◆ dalembert()

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

Performs one time-step integration of the d'Alembert scalar equation.

Parameters
par[input/output] possible parameters to control the resolution of the wave equation. See the actual implementation in the derived class of Map for documentation. Note that, at least, param must contain the time step as first double
parameter.
fJp1[output] solution $f^{J+1}$ at time J+1 with boundary conditions of outgoing radiation (not exact!)
fJ[input] solution $f^J$ at time J
fJm1[input] solution $f^{J-1}$ at time J-1
source[input] source $\sigma$ of the d'Alembert equation $\diamond u = \sigma$.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ dec2_dzpuis()

virtual void Lorene::Map::dec2_dzpuis ( Scalar ) const
pure virtual

Decreases by 2 the value of dzpuis of a Scalar
and changes accordingly its values in the
compactified external domain (CED).

Implemented in Lorene::Map_radial.

◆ dec_dzpuis()

virtual void Lorene::Map::dec_dzpuis ( Scalar ) const
pure virtual

Decreases by 1 the value of dzpuis of a Scalar
and changes accordingly its values in the
compactified external domain (CED).

Implemented in Lorene::Map_radial.

◆ div_cost()

virtual void Lorene::Map::div_cost ( Scalar ) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ div_r()

virtual void Lorene::Map::div_r ( Scalar ) const
pure virtual

Division by r of a Scalar.

Implemented in Lorene::Map_radial.

◆ div_r_zec()

virtual void Lorene::Map::div_r_zec ( Scalar ) const
pure virtual

Division by r (in the compactified external domain only) of a Scalar.

Implemented in Lorene::Map_radial.

◆ div_rsint()

virtual void Lorene::Map::div_rsint ( Scalar ) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ div_sint()

virtual void Lorene::Map::div_sint ( Scalar ) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ div_tant()

virtual void Lorene::Map::div_tant ( Scalar ) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ donne_para_poisson_vect()

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

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 .

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ dsdr() [1/2]

virtual void Lorene::Map::dsdr ( const Cmp ci,
Cmp resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ dsdr() [2/2]

virtual void Lorene::Map::dsdr ( const Scalar uu,
Scalar resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ dsdradial()

virtual void Lorene::Map::dsdradial ( const Scalar uu,
Scalar resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ dsdt()

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

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

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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ dsdxi() [1/2]

virtual void Lorene::Map::dsdxi ( const Cmp ci,
Cmp resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ dsdxi() [2/2]

virtual void Lorene::Map::dsdxi ( const Scalar uu,
Scalar resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ flat_met_cart()

const Metric_flat & Lorene::Map::flat_met_cart ( ) const

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 bvect_cart, and p_flat_met_cart.

◆ flat_met_spher()

const Metric_flat & Lorene::Map::flat_met_spher ( ) const

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 bvect_spher, and p_flat_met_spher.

◆ get_bvect_cart()

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

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 bvect_cart.

◆ get_bvect_spher()

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

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 bvect_spher.

◆ get_mg()

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

Gives the Mg3d on which the mapping is defined.

Definition at line 783 of file map.h.

References mg.

◆ get_ori_x()

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

Returns the x coordinate of the origin.

Definition at line 786 of file map.h.

References ori_x.

◆ get_ori_y()

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

Returns the y coordinate of the origin.

Definition at line 788 of file map.h.

References ori_y.

◆ get_ori_z()

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

Returns the z coordinate of the origin.

Definition at line 790 of file map.h.

References ori_z.

◆ get_rot_phi()

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

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

Definition at line 793 of file map.h.

References rot_phi.

◆ homothetie()

virtual void Lorene::Map::homothetie ( double  lambda)
pure virtual

Sets a new radial scale.

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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ inc2_dzpuis()

virtual void Lorene::Map::inc2_dzpuis ( Scalar ) const
pure virtual

Increases by 2 the value of dzpuis of a Scalar
and changes accordingly its values in the
compactified external domain (CED).

Implemented in Lorene::Map_radial.

◆ inc_dzpuis()

virtual void Lorene::Map::inc_dzpuis ( Scalar ) const
pure virtual

Increases by 1 the value of dzpuis of a Scalar
and changes accordingly its values in the
compactified external domain (CED).

Implemented in Lorene::Map_radial.

◆ integrale()

virtual Tbl* Lorene::Map::integrale ( const Cmp ) const
pure 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.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ lapang()

virtual void Lorene::Map::lapang ( const Scalar uu,
Scalar lap 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ laplacien() [1/2]

virtual void Lorene::Map::laplacien ( const Scalar uu,
int  zec_mult_r,
Scalar lap 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ laplacien() [2/2]

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

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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ mp_angu()

virtual const Map_af& Lorene::Map::mp_angu ( int  ) const
pure virtual

Returns the "angular" mapping for the outside of domain l_zone.

Valid only for the class Map_af.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ mult_cost()

virtual void Lorene::Map::mult_cost ( Scalar ) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ mult_r() [1/2]

virtual void Lorene::Map::mult_r ( Scalar uu) const
pure virtual

Multiplication by r of a Scalar , the dzpuis
of uu is not changed.

Implemented in Lorene::Map_radial.

◆ mult_r() [2/2]

virtual void Lorene::Map::mult_r ( Cmp ci) const
pure virtual

Multiplication by r of a Cmp .

In the CED, there is only a decrement of dzpuis

Implemented in Lorene::Map_radial.

◆ mult_r_zec()

virtual void Lorene::Map::mult_r_zec ( Scalar ) const
pure virtual

Multiplication by r (in the compactified external domain only) of a Scalar.

Implemented in Lorene::Map_radial.

◆ mult_rsint()

virtual void Lorene::Map::mult_rsint ( Scalar ) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ mult_sint()

virtual void Lorene::Map::mult_sint ( Scalar ) const
pure virtual

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

Implemented in Lorene::Map_radial.

◆ operator=()

virtual void Lorene::Map::operator= ( const Map_af )
pure virtual

Assignment to an affine mapping.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_radial, Lorene::Map_eps, and Lorene::Map_star.

◆ operator==()

virtual bool Lorene::Map::operator== ( const Map ) const
pure virtual

Comparison operator (egality)

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_radial, Lorene::Map_eps, and Lorene::Map_star.

◆ operator>>()

virtual ostream& Lorene::Map::operator>> ( ostream &  ) const
privatepure virtual

◆ poisson()

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

Computes the solution of a scalar Poisson equation.

Parameters
source[input] source $\sigma$ of the Poisson equation $\Delta 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 with the boundary condition u =0 at spatial infinity.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ poisson2d()

virtual void Lorene::Map::poisson2d ( const Cmp source_mat,
const Cmp source_quad,
Param par,
Cmp uu 
) const
pure 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] possible parameters to control the resolution of the Poisson equation. See the actual implementation in the derived class of Map for documentation.
uu[input/output] solution u with the boundary condition u =0 at spatial infinity.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ poisson_angu()

virtual void Lorene::Map::poisson_angu ( const Scalar source,
Param par,
Scalar uu,
double  lambda = 0 
) const
pure 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)

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ poisson_compact() [1/2]

virtual void Lorene::Map::poisson_compact ( const Cmp source,
const Cmp aa,
const Tenseur bb,
const Param par,
Cmp psi 
) const
pure virtual

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] 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$.

Implemented in Lorene::Map_radial.

◆ poisson_compact() [2/2]

virtual void Lorene::Map::poisson_compact ( int  nzet,
const Cmp source,
const Cmp aa,
const Tenseur bb,
const Param par,
Cmp psi 
) const
pure virtual

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$.

Implemented in Lorene::Map_radial.

◆ poisson_frontiere()

virtual void Lorene::Map::poisson_frontiere ( const Cmp source,
const Valeur limite,
int  raccord,
int  num_front,
Cmp pot,
double  = 0.,
double  = 0. 
) const
pure virtual

Computes the solution of a Poisson equation from the domain num_front+1 .

imposing a boundary condition at the boundary between the domains num_front and num_front+1 .

Parameters
source[input] : source of the equation.
limite[input] : limite[num_front] contains the angular function being the boudary condition.
raccord[input] : 1 for the Dirichlet problem and 2 for the Neumann one and 3 for Dirichlet-Neumann.
num_front[input] : index of the boudary at which the boundary condition has to be imposed.
pot[output] : result.
fact_dir[input] : Valeur by which we multiply the quantity we solve. (in the case of Dirichlet-Neumann boundary condition.)
fact_neu[input] : Valeur by which we multiply the radial derivative of the quantity we solve. (in the case of Dirichlet-Neumann boundary condition.)

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ poisson_interne()

virtual void Lorene::Map::poisson_interne ( const Cmp source,
const Valeur limite,
Param par,
Cmp pot 
) const
pure 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.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ poisson_regular()

virtual void Lorene::Map::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
pure 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] possible parameters to control the resolution of the Poisson equation. See the actual implementation in the derived class of Map for documentation.
uu[input/output] solution u with the boundary condition u =0 at spatial infinity.
uu_regu[output] solution of the regular part of the source.
uu_div[output] solution of the diverging part of the source.
duu_div[output] derivative of the diverging potential
source_regu[output] regularized source
source_div[output] diverging part of the source

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ poisson_tau()

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

Computes the solution of a scalar Poisson equationwith a Tau method.

Parameters
source[input] source $\sigma$ of the Poisson equation $\Delta 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 with the boundary condition u =0 at spatial infinity.

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ primr()

virtual void Lorene::Map::primr ( const Scalar uu,
Scalar resu,
bool  null_infty 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ reevaluate() [1/2]

virtual void Lorene::Map::reevaluate ( const Map mp_prev,
int  nzet,
Cmp uu 
) const
pure virtual

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.

Implemented in Lorene::Map_radial.

◆ reevaluate() [2/2]

virtual void Lorene::Map::reevaluate ( const Map mp_prev,
int  nzet,
Scalar uu 
) const
pure virtual

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 Scalae at the grid points defined by *this.

Implemented in Lorene::Map_radial.

◆ reevaluate_symy() [1/2]

virtual void Lorene::Map::reevaluate_symy ( const Map mp_prev,
int  nzet,
Cmp uu 
) const
pure virtual

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

Case where the Cmp is symmetric with respect 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.

Implemented in Lorene::Map_radial.

◆ reevaluate_symy() [2/2]

virtual void Lorene::Map::reevaluate_symy ( const Map mp_prev,
int  nzet,
Scalar uu 
) const
pure virtual

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

Case where the Scalar is symmetric with respect 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.

Implemented in Lorene::Map_radial.

◆ reset_coord()

void Lorene::Map::reset_coord ( )
protectedvirtual

Resets all the member Coords.

Reimplemented in Lorene::Map_et, and Lorene::Map_radial.

Definition at line 279 of file map.C.

References cosp, cost, Lorene::Coord::del_t(), phi, r, sinp, sint, tet, x, xa, y, ya, z, and za.

◆ resize()

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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ sauve()

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

◆ set_ori()

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

Sets a new origin.

Definition at line 256 of file map.C.

References bvect_spher, ori_x, ori_y, ori_z, reset_coord(), and Lorene::Base_vect_spher::set_ori().

◆ set_rot_phi()

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

Sets a new rotation angle.

Definition at line 266 of file map.C.

References bvect_cart, bvect_spher, reset_coord(), rot_phi, Lorene::Base_vect_cart::set_rot_phi(), and Lorene::Base_vect_spher::set_rot_phi().

◆ srdsdt() [1/2]

virtual void Lorene::Map::srdsdt ( const Cmp ci,
Cmp resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ srdsdt() [2/2]

virtual void Lorene::Map::srdsdt ( const Scalar uu,
Scalar resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ srstdsdp() [1/2]

virtual void Lorene::Map::srstdsdp ( const Cmp ci,
Cmp resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ srstdsdp() [2/2]

virtual void Lorene::Map::srstdsdp ( const Scalar uu,
Scalar resu 
) const
pure 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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ stdsdp()

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

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

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

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ val_lx() [1/2]

virtual void Lorene::Map::val_lx ( double  rr,
double  theta,
double  pphi,
int &  l,
double &  xi 
) const
pure 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$

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ val_lx() [2/2]

virtual void Lorene::Map::val_lx ( double  rr,
double  theta,
double  pphi,
const Param par,
int &  l,
double &  xi 
) const
pure 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
l[output] value of the domain index
xi[output] value of $\xi$

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

◆ val_r()

virtual double Lorene::Map::val_r ( int  l,
double  xi,
double  theta,
double  pphi 
) const
pure 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')$

Implemented in Lorene::Map_log, Lorene::Map_et, Lorene::Map_af, Lorene::Map_eps, and Lorene::Map_star.

Friends And Related Function Documentation

◆ operator<<

ostream& operator<< ( ostream &  o,
const Map cv 
)
friend

Operator <<.

Definition at line 242 of file map.C.

Member Data Documentation

◆ bvect_cart

Base_vect_cart Lorene::Map::bvect_cart
protected

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
protected

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

$\cos\phi$

Definition at line 742 of file map.h.

◆ cost

Coord Lorene::Map::cost

$\cos\theta$

Definition at line 740 of file map.h.

◆ mg

const Mg3d* Lorene::Map::mg
protected

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
protected

Absolute coordinate x of the origin.

Definition at line 696 of file map.h.

◆ ori_y

double Lorene::Map::ori_y
protected

Absolute coordinate y of the origin.

Definition at line 697 of file map.h.

◆ ori_z

double Lorene::Map::ori_z
protected

Absolute coordinate z of the origin.

Definition at line 698 of file map.h.

◆ p_cmp_zero

Cmp* Lorene::Map::p_cmp_zero
protected

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
mutableprotected

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
mutableprotected

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
mutableprotected

Pointer on the "angular" mapping.

Definition at line 733 of file map.h.

◆ phi

Coord Lorene::Map::phi

$\phi$ coordinate centered on the grid

Definition at line 738 of file map.h.

◆ r

Coord Lorene::Map::r

r coordinate centered on the grid

Definition at line 736 of file map.h.

◆ rot_phi

double Lorene::Map::rot_phi
protected

Angle between the x –axis and X –axis.

Definition at line 699 of file map.h.

◆ sinp

Coord Lorene::Map::sinp

$\sin\phi$

Definition at line 741 of file map.h.

◆ sint

Coord Lorene::Map::sint

$\sin\theta$

Definition at line 739 of file map.h.

◆ tet

Coord Lorene::Map::tet

$\theta$ coordinate centered on the grid

Definition at line 737 of file map.h.

◆ x

Coord Lorene::Map::x

x coordinate centered on the grid

Definition at line 744 of file map.h.

◆ xa

Coord Lorene::Map::xa

Absolute x coordinate.

Definition at line 748 of file map.h.

◆ y

Coord Lorene::Map::y

y coordinate centered on the grid

Definition at line 745 of file map.h.

◆ ya

Coord Lorene::Map::ya

Absolute y coordinate.

Definition at line 749 of file map.h.

◆ z

Coord Lorene::Map::z

z coordinate centered on the grid

Definition at line 746 of file map.h.

◆ za

Coord Lorene::Map::za

Absolute z coordinate.

Definition at line 750 of file map.h.


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