Lorene::Map_radial Class Reference
[Mapping grid -> physical space (spherical coordinates)]

Base class for pure radial mappings. More...

#include <map.h>

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

List of all members.

Public Member Functions

virtual ~Map_radial ()
 Destructor.
virtual void operator= (const Map_af &)=0
 Assignment to an affine mapping.
virtual void sauve (FILE *) const
 Save in a file.
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.
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)$.
virtual bool operator== (const Map &) const =0
 Comparison operator (egality).
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.
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.
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.
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.
virtual void mult_r (Scalar &uu) const
 Multiplication by r of a Scalar, the dzpuis of uu is not changed.
virtual void mult_r (Cmp &ci) const
 Multiplication by r of a Cmp.
virtual void mult_r_zec (Scalar &) const
 Multiplication by r (in the compactified external domain only) of a Scalar.
virtual void mult_rsint (Scalar &) const
 Multiplication by $r\sin\theta$ of a Scalar.
virtual void div_rsint (Scalar &) const
 Division by $r\sin\theta$ of a Scalar.
virtual void div_r (Scalar &) const
 Division by r of a Scalar.
virtual void div_r_zec (Scalar &) const
 Division by r (in the compactified external domain only) of a Scalar.
virtual void mult_cost (Scalar &) const
 Multiplication by $\cos\theta$ of a Scalar.
virtual void div_cost (Scalar &) const
 Division by $\cos\theta$ of a Scalar.
virtual void mult_sint (Scalar &) const
 Multiplication by $\sin\theta$ of a Scalar.
virtual void div_sint (Scalar &) const
 Division by $\sin\theta$ of a Scalar.
virtual void div_tant (Scalar &) const
 Division by $\tan\theta$ of a Scalar.
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.
virtual void comp_x_from_spherical (const Cmp &v_r, const Cmp &v_theta, const Cmp &v_phi, Cmp &v_x) const
 Cmp version
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 .
virtual void comp_y_from_spherical (const Cmp &v_r, const Cmp &v_theta, const Cmp &v_phi, Cmp &v_y) const
 Cmp version
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 .
virtual void comp_z_from_spherical (const Cmp &v_r, const Cmp &v_theta, Cmp &v_z) const
 Cmp version
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 .
virtual void comp_r_from_cartesian (const Cmp &v_x, const Cmp &v_y, const Cmp &v_z, Cmp &v_r) const
 Cmp version
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 .
virtual void comp_t_from_cartesian (const Cmp &v_x, const Cmp &v_y, const Cmp &v_z, Cmp &v_t) const
 Cmp version
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 .
virtual void comp_p_from_cartesian (const Cmp &v_x, const Cmp &v_y, Cmp &v_p) const
 Cmp version
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).
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).
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).
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).
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.
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.
const Mg3dget_mg () const
 Gives the Mg3d on which the mapping is defined.
double get_ori_x () const
 Returns the x coordinate of the origin.
double get_ori_y () const
 Returns the y coordinate of the origin.
double get_ori_z () const
 Returns the z coordinate of the origin.
double get_rot_phi () const
 Returns the angle between the x --axis and X --axis.
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.
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.
const Metric_flatflat_met_spher () const
 Returns the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher.
const Metric_flatflat_met_cart () const
 Returns the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart.
const Cmpcmp_zero () const
 Returns the null Cmp defined on *this.
virtual const Map_afmp_angu (int) const =0
 Returns the "angular" mapping for the outside of domain l_zone.
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).
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.
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)$.
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)$.
virtual bool operator== (const Map &) const =0
 Comparison operator (egality).
void set_ori (double xa0, double ya0, double za0)
 Sets a new origin.
void set_rot_phi (double phi0)
 Sets a new rotation angle.
virtual void homothetie (double lambda)=0
 Sets a new radial scale.
virtual void resize (int l, double lambda)=0
 Rescales the outer boundary of one domain.
virtual void adapt (const Cmp &ent, const Param &par, int nbr=0)=0
 Adaptation of the mapping to a given scalar field.
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.
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.
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.
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.
virtual void dsdxi (const Cmp &ci, Cmp &resu) const =0
 Computes $\partial/ \partial \xi$ of a Cmp .
virtual void dsdxi (const Scalar &uu, Scalar &resu) const =0
 Computes $\partial/ \partial xi$ of a Scalar .
virtual void dsdr (const Cmp &ci, Cmp &resu) const =0
 Computes $\partial/ \partial r$ of a Cmp .
virtual void dsdr (const Scalar &uu, Scalar &resu) const =0
 Computes $\partial/ \partial r$ of a Scalar .
virtual void srdsdt (const Cmp &ci, Cmp &resu) const =0
 Computes $1/r \partial/ \partial \theta$ of a Cmp .
virtual void srdsdt (const Scalar &uu, Scalar &resu) const =0
 Computes $1/r \partial/ \partial \theta$ of a Scalar .
virtual void srstdsdp (const Cmp &ci, Cmp &resu) const =0
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Cmp .
virtual void srstdsdp (const Scalar &uu, Scalar &resu) const =0
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Scalar .
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.
virtual void dsdt (const Scalar &uu, Scalar &resu) const =0
 Computes $\partial/ \partial \theta$ of a Scalar .
virtual void stdsdp (const Scalar &uu, Scalar &resu) const =0
 Computes $1/\sin\theta \partial/ \partial \varphi$ of a Scalar .
virtual void laplacien (const Scalar &uu, int zec_mult_r, Scalar &lap) const =0
 Computes the Laplacian of a scalar field.
virtual void laplacien (const Cmp &uu, int zec_mult_r, Cmp &lap) const =0
 Computes the Laplacian of a scalar field (Cmp version).
virtual void lapang (const Scalar &uu, Scalar &lap) const =0
 Computes the angular Laplacian of a scalar field.
virtual void primr (const Scalar &uu, Scalar &resu, bool null_infty) const =0
 Computes the radial primitive which vanishes for $r\to \infty$.
virtual Tblintegrale (const Cmp &) const =0
 Computes the integral over all space of a Cmp .
virtual void poisson (const Cmp &source, Param &par, Cmp &uu) const =0
 Computes the solution of a scalar Poisson equation.
virtual void poisson_tau (const Cmp &source, Param &par, Cmp &uu) const =0
 Computes the solution of a scalar Poisson equationwith a Tau method.
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.
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.
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 .
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 .
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.
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.
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.

Public Attributes

Coord xsr
 $\xi/R$ in the nucleus; \ 1/R in the non-compactified shells; \ $(\xi-1)/U$ in the compactified outer domain.
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.
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).
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).
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.
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.
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.
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.
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.
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.
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.
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.
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.
Coord r
 r coordinate centered on the grid
Coord tet
 $\theta$ coordinate centered on the grid
Coord phi
 $\phi$ coordinate centered on the grid
Coord sint
 $\sin\theta$
Coord cost
 $\cos\theta$
Coord sinp
 $\sin\phi$
Coord cosp
 $\cos\phi$
Coord x
 x coordinate centered on the grid
Coord y
 y coordinate centered on the grid
Coord z
 z coordinate centered on the grid
Coord xa
 Absolute x coordinate.
Coord ya
 Absolute y coordinate.
Coord za
 Absolute z coordinate.

Protected Member Functions

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

Protected Attributes

const Mg3dmg
 Pointer on the multi-grid Mgd3 on which this is defined.
double ori_x
 Absolute coordinate x of the origin.
double ori_y
 Absolute coordinate y of the origin.
double ori_z
 Absolute coordinate z of the origin.
double rot_phi
 Angle between the x --axis and X --axis.
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.
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.
Metric_flatp_flat_met_spher
 Pointer onto the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher.
Metric_flatp_flat_met_cart
 Pointer onto the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart.
Cmpp_cmp_zero
 The null Cmp.
Map_afp_mp_angu
 Pointer on the "angular" mapping.

Friends

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

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 1539 of file map.h.


Constructor & Destructor Documentation

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.

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

Copy constructor.

Definition at line 99 of file map_radial.C.

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.

Lorene::Map_radial::~Map_radial (  )  [virtual]

Destructor.

Definition at line 113 of file map_radial.C.


Member Function Documentation

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

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_af, Lorene::Map_et, and Lorene::Map_log.

const Cmp& Lorene::Map::cmp_zero (  )  const [inline, inherited]

Returns the null Cmp defined on *this.

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

Definition at line 810 of file map.h.

References Lorene::Map::p_cmp_zero.

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

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::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), and Lorene::Scalar::set_dzpuis().

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

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::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().

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

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::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().

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

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::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().

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

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::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_cp(), Lorene::Valeur::mult_ct(), Lorene::Valeur::mult_sp(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().

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

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::check_dzpuis(), Lorene::Scalar::dz_nonzero(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().

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

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

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_af, Lorene::Map_et, and Lorene::Map_log.

void Lorene::Map_radial::dec2_dzpuis ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::dec_dzpuis ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::div_cost ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::div_r ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::div_r_zec ( Scalar uu  )  const [virtual]
void Lorene::Map_radial::div_rsint ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::div_sint ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::div_tant ( Scalar ci  )  const [virtual]
virtual Param* Lorene::Map::donne_para_poisson_vect ( Param para,
int  i 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual void Lorene::Map::dsdr ( const Scalar uu,
Scalar resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual void Lorene::Map::dsdr ( const Cmp ci,
Cmp resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual void Lorene::Map::dsdradial ( const Scalar uu,
Scalar resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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

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

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

virtual void Lorene::Map::dsdxi ( const Scalar uu,
Scalar resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual void Lorene::Map::dsdxi ( const Cmp ci,
Cmp resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

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.

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.

const Base_vect_cart& Lorene::Map::get_bvect_cart (  )  const [inline, inherited]

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 794 of file map.h.

References Lorene::Map::bvect_cart.

const Base_vect_spher& Lorene::Map::get_bvect_spher (  )  const [inline, inherited]

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 786 of file map.h.

References Lorene::Map::bvect_spher.

const Mg3d* Lorene::Map::get_mg (  )  const [inline, inherited]

Gives the Mg3d on which the mapping is defined.

Definition at line 768 of file map.h.

References Lorene::Map::mg.

double Lorene::Map::get_ori_x (  )  const [inline, inherited]

Returns the x coordinate of the origin.

Definition at line 771 of file map.h.

References Lorene::Map::ori_x.

double Lorene::Map::get_ori_y (  )  const [inline, inherited]

Returns the y coordinate of the origin.

Definition at line 773 of file map.h.

References Lorene::Map::ori_y.

double Lorene::Map::get_ori_z (  )  const [inline, inherited]

Returns the z coordinate of the origin.

Definition at line 775 of file map.h.

References Lorene::Map::ori_z.

double Lorene::Map::get_rot_phi (  )  const [inline, inherited]

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

Definition at line 778 of file map.h.

References Lorene::Map::rot_phi.

virtual void Lorene::Map::homothetie ( double  lambda  )  [pure virtual, inherited]

Sets a new radial scale.

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

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

void Lorene::Map_radial::inc2_dzpuis ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::inc_dzpuis ( Scalar ci  )  const [virtual]
virtual Tbl* Lorene::Map::integrale ( const Cmp  )  const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual void Lorene::Map::lapang ( const Scalar uu,
Scalar lap 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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

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

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

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual const Map_af& Lorene::Map::mp_angu ( int   )  const [pure virtual, inherited]

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

Valid only for the class Map_af.

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

void Lorene::Map_radial::mult_cost ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::mult_r ( Cmp ci  )  const [virtual]
void Lorene::Map_radial::mult_r ( Scalar uu  )  const [virtual]
void Lorene::Map_radial::mult_r_zec ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::mult_rsint ( Scalar ci  )  const [virtual]
void Lorene::Map_radial::mult_sint ( Scalar ci  )  const [virtual]
virtual void Lorene::Map_radial::operator= ( const Map_af  )  [pure virtual]

Assignment to an affine mapping.

Implements Lorene::Map.

Implemented in Lorene::Map_et, and Lorene::Map_log.

virtual bool Lorene::Map::operator== ( const Map  )  const [pure virtual, inherited]

Comparison operator (egality).

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

Comparison operator (egality).

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

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

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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_af, Lorene::Map_et, and Lorene::Map_log.

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::annule(), Lorene::Mtbl::annule_hard(), Lorene::Tbl::annule_hard(), Lorene::Map::bvect_spher, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::diffrel(), Lorene::Map_af::dsdr(), Lorene::Cmp::dsdr(), Lorene::Param::get_double(), Lorene::Tenseur::get_etat(), Lorene::Cmp::get_etat(), Lorene::Mg3d::get_grille3d(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Tenseur::get_mp(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Tenseur::get_triad(), Lorene::Cmp::laplacien(), Lorene::Map_af::laplacien(), Lorene::max(), Lorene::Map::mg, Lorene::min(), poisson_compact(), Lorene::Mtbl::set(), Lorene::Tbl::set(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Cmp::srdsdt(), Lorene::Cmp::srstdsdp(), Lorene::Valeur::std_base_scal(), Lorene::Tbl::t, Lorene::Cmp::va, Lorene::Grille3d::x, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().

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::annule(), Lorene::Map::bvect_spher, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Valeur::d2sdx2(), Lorene::diffrel(), Lorene::Cmp::dsdr(), Lorene::Valeur::dsdx(), dxdr, Lorene::Param::get_double(), Lorene::Tenseur::get_etat(), Lorene::Cmp::get_etat(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Tenseur::get_mp(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nzone(), Lorene::Tenseur::get_triad(), Lorene::Valeur::lapang(), Lorene::Cmp::laplacien(), Lorene::max(), Lorene::Map::mg, Lorene::min(), Lorene::Valeur::mult_x(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Cmp::srdsdt(), Lorene::Cmp::srstdsdp(), Lorene::Valeur::std_base_scal(), Lorene::Valeur::sx(), Lorene::Cmp::va, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().

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, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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_af, Lorene::Map_et, and Lorene::Map_log.

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, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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_infty if true (default), the primitive is null at infinity (or on the grid boundary). F is null at the center otherwise

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

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

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.
virtual void Lorene::Map::reevaluate ( const Map mp_prev,
int  nzet,
Cmp uu 
) const [pure virtual, inherited]

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

Definition at line 176 of file map_radial_reevaluate.C.

References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Scalar::annule(), Lorene::Coord::c, Lorene::Valeur::c, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Coord::fait(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, Lorene::Map::r, Lorene::Valeur::set_etat_c_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Mtbl::t, val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk().

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.

Definition at line 61 of file map_radial_reevaluate.C.

References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Cmp::annule(), Lorene::Coord::c, Lorene::Coord::fait(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, Lorene::Map::r, Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::set_etat_qcq(), Lorene::Mtbl::t, Lorene::Cmp::va, val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk().

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

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.
virtual void Lorene::Map::reevaluate_symy ( const Map mp_prev,
int  nzet,
Cmp uu 
) const [pure virtual, inherited]

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

Definition at line 196 of file map_radial_reeval_symy.C.

References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Scalar::annule(), Lorene::Coord::c, Lorene::Valeur::c, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Coord::fait(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_p(), Lorene::Map::mg, Lorene::Map::r, Lorene::Valeur::set_etat_c_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Mtbl::t, val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk_symy().

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.

Definition at line 62 of file map_radial_reeval_symy.C.

References Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Cmp::annule(), Lorene::Coord::c, Lorene::Coord::fait(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_p(), Lorene::Map::mg, Lorene::Map::r, Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::set_etat_qcq(), Lorene::Mtbl::t, Lorene::Cmp::va, val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk_symy().

void Lorene::Map_radial::reset_coord (  )  [protected, virtual]

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, sr2d2rdt2, sr2drdt, sr2stdrdp, srdrdt, srstdrdp, sstd2rdpdx, stdrdp, and xsr.

virtual void Lorene::Map::resize ( int  l,
double  lambda 
) [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

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

Save in a file.

Reimplemented from Lorene::Map.

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

Definition at line 119 of file map_radial.C.

void Lorene::Map::set_ori ( double  xa0,
double  ya0,
double  za0 
) [inherited]
void Lorene::Map::set_rot_phi ( double  phi0  )  [inherited]
virtual void Lorene::Map::srdsdt ( const Scalar uu,
Scalar resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual void Lorene::Map::srdsdt ( const Cmp ci,
Cmp resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual void Lorene::Map::srstdsdp ( const Scalar uu,
Scalar resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

virtual void Lorene::Map::srstdsdp ( const Cmp ci,
Cmp resu 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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

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

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

virtual void Lorene::Map::val_lx ( double  rr,
double  theta,
double  pphi,
const Param par,
int &  l,
double &  xi 
) const [pure virtual, inherited]

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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_af, Lorene::Map_et, and Lorene::Map_log.

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_af, Lorene::Map_et, and Lorene::Map_log.

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

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_af, Lorene::Map_et, and Lorene::Map_log.

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_af, Lorene::Map_et, and Lorene::Map_log.


Friends And Related Function Documentation

ostream& operator<< ( ostream &  ,
const Map  
) [friend, inherited]

Operator <<.


Member Data Documentation

Base_vect_cart Lorene::Map::bvect_cart [protected, inherited]

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 700 of file map.h.

Base_vect_spher Lorene::Map::bvect_spher [protected, inherited]

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 692 of file map.h.

Coord Lorene::Map::cosp [inherited]

$\cos\phi$

Definition at line 727 of file map.h.

Coord Lorene::Map::cost [inherited]

$\cos\theta$

Definition at line 725 of file map.h.

$\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 1643 of file map.h.

$\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 1622 of file map.h.

$\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 1571 of file map.h.

$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 1563 of file map.h.

$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 1634 of file map.h.

const Mg3d* Lorene::Map::mg [protected, inherited]

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

Definition at line 679 of file map.h.

double Lorene::Map::ori_x [protected, inherited]

Absolute coordinate x of the origin.

Definition at line 681 of file map.h.

double Lorene::Map::ori_y [protected, inherited]

Absolute coordinate y of the origin.

Definition at line 682 of file map.h.

double Lorene::Map::ori_z [protected, inherited]

Absolute coordinate z of the origin.

Definition at line 683 of file map.h.

Cmp* Lorene::Map::p_cmp_zero [protected, inherited]

The null Cmp.

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

Definition at line 716 of file map.h.

Metric_flat* Lorene::Map::p_flat_met_cart [mutable, protected, inherited]

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

Definition at line 710 of file map.h.

Metric_flat* Lorene::Map::p_flat_met_spher [mutable, protected, inherited]

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

Definition at line 705 of file map.h.

Map_af* Lorene::Map::p_mp_angu [mutable, protected, inherited]

Pointer on the "angular" mapping.

Definition at line 718 of file map.h.

Coord Lorene::Map::phi [inherited]

$\phi$ coordinate centered on the grid

Definition at line 723 of file map.h.

Coord Lorene::Map::r [inherited]

r coordinate centered on the grid

Definition at line 721 of file map.h.

double Lorene::Map::rot_phi [protected, inherited]

Angle between the x --axis and X --axis.

Definition at line 684 of file map.h.

Coord Lorene::Map::sinp [inherited]

$\sin\phi$

Definition at line 726 of file map.h.

Coord Lorene::Map::sint [inherited]

$\sin\theta$

Definition at line 724 of file map.h.

$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 1660 of file map.h.

$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 1603 of file map.h.

$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 1611 of file map.h.

$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 1587 of file map.h.

$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 1595 of file map.h.

$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 1651 of file map.h.

${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 1579 of file map.h.

Coord Lorene::Map::tet [inherited]

$\theta$ coordinate centered on the grid

Definition at line 722 of file map.h.

Coord Lorene::Map::x [inherited]

x coordinate centered on the grid

Definition at line 729 of file map.h.

Coord Lorene::Map::xa [inherited]

Absolute x coordinate.

Definition at line 733 of file map.h.

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

Definition at line 1552 of file map.h.

Coord Lorene::Map::y [inherited]

y coordinate centered on the grid

Definition at line 730 of file map.h.

Coord Lorene::Map::ya [inherited]

Absolute y coordinate.

Definition at line 734 of file map.h.

Coord Lorene::Map::z [inherited]

z coordinate centered on the grid

Definition at line 731 of file map.h.

Coord Lorene::Map::za [inherited]

Absolute z coordinate.

Definition at line 735 of file map.h.


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

Generated on 7 Dec 2019 for LORENE by  doxygen 1.6.1