LORENE
Lorene::Map_radial Class Referenceabstract

Base class for pure radial mappings. More...

#include <map.h>

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

Public Member Functions

virtual ~Map_radial ()
 Destructor. More...
 
virtual void operator= (const Map_af &)=0
 Assignment to an affine mapping. More...
 
virtual void sauve (FILE *) const
 Save in a file. More...
 
virtual double val_r_jk (int l, double xi, int j, int k) const =0
 Returns the value of the radial coordinate r for a given $\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 =0
 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 =0
 Comparison operator (egality) More...
 
virtual void reevaluate (const Map *mp_prev, int nzet, Cmp &uu) const
 Recomputes the values of a Cmp at the collocation points after a change in the mapping. More...
 
virtual void reevaluate_symy (const Map *mp_prev, int nzet, Cmp &uu) const
 Recomputes the values of a Cmp at the collocation points after a change in the mapping. More...
 
virtual void reevaluate (const Map *mp_prev, int nzet, Scalar &uu) const
 Recomputes the values of a Scalar at the collocation points after a change in the mapping. More...
 
virtual void reevaluate_symy (const Map *mp_prev, int nzet, Scalar &uu) const
 Recomputes the values of a Scalar at the collocation points after a change in the mapping. More...
 
virtual void mult_r (Scalar &uu) const
 Multiplication by r of a Scalar, the dzpuis
of uu is not changed. More...
 
virtual void mult_r (Cmp &ci) const
 Multiplication by r of a Cmp. More...
 
virtual void mult_r_zec (Scalar &) const
 Multiplication by r (in the compactified external domain only) of a Scalar. More...
 
virtual void mult_rsint (Scalar &) const
 Multiplication by $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...
 
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...
 
void set_ori (double xa0, double ya0, double za0)
 Sets a new origin. More...
 
void set_rot_phi (double phi0)
 Sets a new rotation angle. More...
 
virtual void homothetie (double lambda)=0
 Sets a new radial scale. More...
 
virtual void resize (int l, double lambda)=0
 Rescales the outer boundary of one domain. More...
 
virtual void adapt (const Cmp &ent, const Param &par, int nbr=0)=0
 Adaptation of the mapping to a given scalar field. More...
 
virtual void dsdxi (const Cmp &ci, Cmp &resu) const =0
 Computes $\partial/ \partial \xi$ 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 Cmp &ci, Cmp &resu) const =0
 Computes $\partial/ \partial r$ of a Cmp . More...
 
virtual void dsdr (const Scalar &uu, Scalar &resu) const =0
 Computes $\partial/ \partial r$ of a Scalar . More...
 
virtual void srdsdt (const Cmp &ci, Cmp &resu) const =0
 Computes $1/r \partial/ \partial \theta$ of a Cmp . 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 Cmp &ci, Cmp &resu) const =0
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Cmp . 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 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 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 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_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 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

 Map_radial (const Mg3d &mgrid)
 Constructor from a grid (protected to make Map_radial an abstract class) More...
 
 Map_radial (const Map_radial &mp)
 Copy constructor. More...
 
 Map_radial (const Mg3d &, FILE *)
 Constructor from a file (see sauve(FILE* ) ) More...
 
virtual void reset_coord ()
 Resets all the member Coords. More...
 

Protected Attributes

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

Detailed Description

Base class for pure radial mappings.

()

A pure radial mapping is a mapping of the type $r=R(\xi, \theta', \phi')$, $\theta=\theta'$, $\phi=\phi'$. The class Map_radial is an abstract class. Effective implementations of radial mapping are performed by the derived class Map_af and Map_et .

Definition at line 1551 of file map.h.

Constructor & Destructor Documentation

◆ Map_radial() [1/3]

Lorene::Map_radial::Map_radial ( const Mg3d mgrid)
protected

Constructor from a grid (protected to make Map_radial an abstract class)

Definition at line 92 of file map_radial.C.

◆ Map_radial() [2/3]

Lorene::Map_radial::Map_radial ( const Map_radial mp)
protected

Copy constructor.

Definition at line 99 of file map_radial.C.

◆ Map_radial() [3/3]

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

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

Definition at line 106 of file map_radial.C.

◆ ~Map_radial()

Lorene::Map_radial::~Map_radial ( )
virtual

Destructor.

Definition at line 113 of file map_radial.C.

Member Function Documentation

◆ adapt()

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

Adaptation of the mapping to a given scalar field.

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

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

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

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
virtual

Cmp version

Implements Lorene::Map.

Definition at line 179 of file map_radial_comp_rtp.C.

References 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
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

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
virtual

Cmp version

Implements Lorene::Map.

Definition at line 68 of file map_radial_comp_rtp.C.

References 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
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

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
virtual

Cmp version

Implements Lorene::Map.

Definition at line 124 of file map_radial_comp_rtp.C.

References 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
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

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
virtual

Cmp version

Implements Lorene::Map.

Definition at line 71 of file map_radial_comp_xyz.C.

References 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
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

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
virtual

Cmp version

Implements Lorene::Map.

Definition at line 129 of file map_radial_comp_xyz.C.

References 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
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

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
virtual

Cmp version

Implements Lorene::Map.

Definition at line 187 of file map_radial_comp_xyz.C.

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

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ dec2_dzpuis()

void Lorene::Map_radial::dec2_dzpuis ( Scalar ci) const
virtual

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

Implements Lorene::Map.

Definition at line 751 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ dec_dzpuis()

void Lorene::Map_radial::dec_dzpuis ( Scalar ci) const
virtual

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

Implements Lorene::Map.

Definition at line 646 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ div_cost()

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

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
virtual

Division by r of a Scalar.

Implements Lorene::Map.

Definition at line 517 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ div_r_zec()

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

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

Implements Lorene::Map.

Definition at line 588 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ div_rsint()

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

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
virtual

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
virtual

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

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

Function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara .

It constructs the sets of parameters used for each scalar Poisson equation from the one for the vectorial one.

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ dsdr() [1/2]

virtual void Lorene::Map::dsdr ( const Cmp ci,
Cmp resu 
) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ dsdr() [2/2]

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ dsdradial()

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ dsdt()

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

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

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

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

◆ dsdxi() [1/2]

virtual void Lorene::Map::dsdxi ( const Cmp ci,
Cmp resu 
) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ dsdxi() [2/2]

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

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

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

Sets a new radial scale.

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

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

◆ inc2_dzpuis()

void Lorene::Map_radial::inc2_dzpuis ( Scalar ci) const
virtual

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

Implements Lorene::Map.

Definition at line 802 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ inc_dzpuis()

void Lorene::Map_radial::inc_dzpuis ( Scalar ci) const
virtual

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

Implements Lorene::Map.

Definition at line 699 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ integrale()

virtual Tbl* Lorene::Map::integrale ( const Cmp ) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ lapang()

virtual void Lorene::Map::lapang ( const Scalar uu,
Scalar lap 
) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ laplacien() [1/2]

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ laplacien() [2/2]

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

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

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

◆ mp_angu()

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

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

Valid only for the class Map_af.

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

◆ mult_cost()

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

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
virtual

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

Implements Lorene::Map.

Definition at line 147 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ mult_r() [2/2]

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

Multiplication by r of a Cmp.

In the CED, there is only a decrement of dzpuis

Implements Lorene::Map.

Definition at line 222 of file map_radial_r_manip.C.

References Lorene::Cmp::get_etat().

◆ mult_r_zec()

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

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

Implements Lorene::Map.

Definition at line 299 of file map_radial_r_manip.C.

References Lorene::Scalar::get_etat().

◆ mult_rsint()

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

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
virtual

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

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

Assignment to an affine mapping.

Implements Lorene::Map.

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

◆ operator==()

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

Comparison operator (egality)

Implements Lorene::Map.

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

◆ poisson()

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ poisson2d()

virtual void Lorene::Map::poisson2d ( const Cmp source_mat,
const Cmp source_quad,
Param par,
Cmp uu 
) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ poisson_angu()

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

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

Implements Lorene::Map.

Definition at line 456 of file map_radial_poisson_cpt.C.

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

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

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

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ poisson_interne()

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ poisson_tau()

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ primr()

virtual void Lorene::Map::primr ( const Scalar uu,
Scalar resu,
bool  null_infty 
) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ reevaluate() [1/2]

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

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
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 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
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 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
virtual

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

Case where the Scalar is symmetric with respect to the plane y=0.

Parameters
mp_prev[input] Previous value of the mapping.
nzet[input] Number of domains where the computation must be done: the computation is performed for the domains of index $0\le {\tt l} \le {\tt nzet-1}$; uu is set to zero in the other domains.
uu[input/output] input : Scalar previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Scalar at the grid points defined by *this.

Implements Lorene::Map.

Definition at line 199 of file map_radial_reeval_symy.C.

References Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.

◆ reset_coord()

void Lorene::Map_radial::reset_coord ( )
protectedvirtual

Resets all the member Coords.

Reimplemented from Lorene::Map.

Reimplemented in Lorene::Map_et.

Definition at line 129 of file map_radial.C.

References d2rdtdx, d2rdx2, Lorene::Coord::del_t(), drdt, dxdr, lapr_tp, Lorene::Map::reset_coord(), sr2d2rdt2, sr2drdt, sr2stdrdp, srdrdt, srstdrdp, sstd2rdpdx, stdrdp, and xsr.

◆ resize()

virtual void Lorene::Map::resize ( int  l,
double  lambda 
)
pure virtualinherited

Rescales the outer boundary of one domain.

The inner boundary is unchanged. The inner boundary of the next domain is changed to match the new outer boundary.

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ sauve()

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

Save in a file.

Reimplemented from Lorene::Map.

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

Definition at line 119 of file map_radial.C.

References Lorene::Map::sauve().

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

virtual void Lorene::Map::srdsdt ( const Cmp ci,
Cmp resu 
) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ srdsdt() [2/2]

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ srstdsdp() [1/2]

virtual void Lorene::Map::srstdsdp ( const Cmp ci,
Cmp resu 
) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ srstdsdp() [2/2]

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ stdsdp()

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

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

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

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

◆ val_lx() [1/2]

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

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

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ val_lx_jk()

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

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

Parameters
rr[input] value of r
j[input] index of the collocation point in $\theta$
k[input] index of the collocation point in $\phi$
par[input/output] parameters to control the accuracy of the computation
l[output] value of the domain index
xi[output] value of $\xi$

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

◆ val_r()

virtual double Lorene::Map::val_r ( int  l,
double  xi,
double  theta,
double  pphi 
) const
pure virtualinherited

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_eps, Lorene::Map_star, Lorene::Map_log, Lorene::Map_et, and Lorene::Map_af.

◆ val_r_jk()

virtual double Lorene::Map_radial::val_r_jk ( int  l,
double  xi,
int  j,
int  k 
) const
pure 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)$

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

Member Data Documentation

◆ 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

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

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

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

$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

$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

$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

$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

$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

$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

$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

$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

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

$\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: