Map_af Class Reference
[Mapping grid -> physical space (spherical coordinates)]

Affine radial mapping. More...

#include <map.h>

Inheritance diagram for Map_af:
Map_radial Map

List of all members.

Public Member Functions

 Map_af (const Mg3d &mgrille, const double *r_limits)
 Standard Constructor.
 Map_af (const Mg3d &mgrille, const Tbl &r_limits)
 Standard Constructor with Tbl.
 Map_af (const Map_af &)
 Copy constructor.
 Map_af (const Mg3d &, FILE *)
 Constructor from a file (see sauve(FILE*) ).
 Map_af (const Map &)
 Constructor from a general mapping.
virtual ~Map_af ()
 Destructor.
virtual void operator= (const Map_af &)
 Assignment to another affine mapping.
const double * get_alpha () const
 Returns the pointer on the array alpha.
const double * get_beta () const
 Returns the pointer on the array beta.
virtual const Map_afmp_angu (int) const
 Returns the "angular" mapping for the outside of domain l_zone.
virtual double val_r (int l, double xi, double theta, double pphi) const
 Returns the value of the radial coordinate r for a given $(\xi, \theta', \phi')$ in a given domain.
virtual void val_lx (double rr, double theta, double pphi, int &l, double &xi) const
 Computes the domain index l and the value of $\xi$ corresponding to a point given by its physical coordinates $(r, \theta, \phi)$.
virtual void val_lx (double rr, double theta, double pphi, const Param &par, int &l, double &xi) const
 Computes the domain index l and the value of $\xi$ corresponding to a point given by its physical coordinates $(r, \theta, \phi)$.
virtual double val_r_jk (int l, double xi, int j, int k) const
 Returns the value of the radial coordinate r for a given $\xi$ and a given collocation point in $(\theta', \phi')$ in a given domain.
virtual void val_lx_jk (double rr, int j, int k, const Param &par, int &l, double &xi) const
 Computes the domain index l and the value of $\xi$ corresponding to a point of arbitrary r but collocation values of $(\theta, \phi)$.
virtual bool operator== (const Map &) const
 Comparison operator (egality).
virtual void sauve (FILE *) const
 Save in a file.
virtual void homothetie (double lambda)
 Sets a new radial scale.
virtual void resize (int l, double lambda)
 Rescales the outer boundary of one domain.
void homothetie_interne (double lambda)
 Sets a new radial scale at the bondary between the nucleus and the first shell.
virtual void adapt (const Cmp &ent, const Param &par, int nbr=0)
 Adaptation of the mapping to a given scalar field.
void set_alpha (double alpha0, int l)
 Modifies the value of $\alpha$ in domain no. l.
void set_beta (double beta0, int l)
 Modifies the value of $\beta$ in domain no. l.
virtual void dsdxi (const Cmp &ci, Cmp &resu) const
 Computes $\partial/ \partial \xi$ of a Cmp.
virtual void dsdr (const Cmp &ci, Cmp &resu) const
 Computes $\partial/ \partial r$ of a Cmp.
virtual void srdsdt (const Cmp &ci, Cmp &resu) const
 Computes $1/r \partial/ \partial \theta$ of a Cmp.
virtual void srstdsdp (const Cmp &ci, Cmp &resu) const
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Cmp.
virtual void dsdr (const Scalar &uu, Scalar &resu) const
 Computes $\partial/ \partial r$ of a Scalar.
virtual void dsdxi (const Scalar &uu, Scalar &resu) const
 Computes $\partial/ \partial \xi$ of a Scalar.
virtual void dsdradial (const Scalar &, Scalar &) const
 Computes $\partial/ \partial r$ of a Scalar.
virtual void srdsdt (const Scalar &uu, Scalar &resu) const
 Computes $1/r \partial/ \partial \theta$ of a Scalar.
virtual void srstdsdp (const Scalar &uu, Scalar &resu) const
 Computes $1/(r\sin\theta) \partial/ \partial \phi$ of a Scalar.
virtual void dsdt (const Scalar &uu, Scalar &resu) const
 Computes $\partial/ \partial \theta$ of a Scalar.
virtual void stdsdp (const Scalar &uu, Scalar &resu) const
 Computes $1/\sin\theta \partial/ \partial \varphi$ of a Scalar.
virtual void laplacien (const Scalar &uu, int zec_mult_r, Scalar &lap) const
 Computes the Laplacian of a scalar field.
virtual void laplacien (const Cmp &uu, int zec_mult_r, Cmp &lap) const
 Computes the Laplacian of a scalar field (Cmp version).
virtual void lapang (const Scalar &uu, Scalar &lap) const
 Computes the angular Laplacian of a scalar field.
virtual void primr (const Scalar &uu, Scalar &resu, bool null_infty) const
 Computes the radial primitive which vanishes for $r\to \infty$.
virtual Tblintegrale (const Cmp &) const
 Computes the integral over all space of a Cmp.
virtual void poisson (const Cmp &source, Param &par, Cmp &uu) const
 Computes the solution of a scalar Poisson equation.
virtual void poisson_tau (const Cmp &source, Param &par, Cmp &uu) const
 Computes the solution of a scalar Poisson equation using a Tau method.
virtual void poisson_falloff (const Cmp &source, Param &par, Cmp &uu, int k_falloff) const
virtual void poisson_ylm (const Cmp &source, Param &par, Cmp &pot, int nylm, double *intvec) const
virtual void poisson_regular (const Cmp &source, int k_div, int nzet, double unsgam1, Param &par, Cmp &uu, Cmp &uu_regu, Cmp &uu_div, Tenseur &duu_div, Cmp &source_regu, Cmp &source_div) const
 Computes the solution of a scalar Poisson equation.
virtual void poisson_angu (const Scalar &source, Param &par, Scalar &uu, double lambda=0) const
 Computes the solution of the generalized angular Poisson equation.
virtual Paramdonne_para_poisson_vect (Param &par, int i) const
 Internal function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara .
virtual void poisson_frontiere (const Cmp &, const Valeur &, int, int, Cmp &, double=0., double=0.) const
 Solver of the Poisson equation with boundary condition for the affine mapping case.
virtual void poisson_frontiere_double (const Cmp &source, const Valeur &lim_func, const Valeur &lim_der, int num_zone, Cmp &pot) const
 Solver of the Poisson equation with boundary condition for the affine mapping case, cases with boundary conditions of both Dirichlet and Neumann type (no condition imposed at infinity).
virtual void poisson_interne (const Cmp &source, const Valeur &limite, Param &par, Cmp &pot) const
 Computes the solution of a Poisson equation in the shell, imposing a boundary condition at the surface of the star.
double integrale_surface (const Cmp &ci, double rayon) const
 Performs the surface integration of ci on the sphere of radius rayon .
double integrale_surface (const Scalar &ci, double rayon) const
 Performs the surface integration of ci on the sphere of radius rayon .
double integrale_surface_falloff (const Cmp &ci) const
double integrale_surface_infini (const Cmp &ci) const
 Performs the surface integration of ci at infinity.
double integrale_surface_infini (const Scalar &ci) const
 Performs the surface integration of ci at infinity.
void sol_elliptic (Param_elliptic &params, const Scalar &so, Scalar &uu) const
 General elliptic solver.
void sol_elliptic_boundary (Param_elliptic &params, const Scalar &so, Scalar &uu, const Mtbl_cf &bound, double fact_dir, double fact_neu) const
 General elliptic solver including inner boundary conditions.
void sol_elliptic_boundary (Param_elliptic &params, const Scalar &so, Scalar &uu, const Scalar &bound, double fact_dir, double fact_neu) const
 General elliptic solver including inner boundary conditions, with boundary given as a Scalar on mono-domain angular grid.
void sol_elliptic_no_zec (Param_elliptic &params, const Scalar &so, Scalar &uu, double val) const
 General elliptic solver.
void sol_elliptic_only_zec (Param_elliptic &params, const Scalar &so, Scalar &uu, double val) const
 General elliptic solver.
void sol_elliptic_sin_zec (Param_elliptic &params, const Scalar &so, Scalar &uu, double *coefs, double *) const
 General elliptic solver.
void sol_elliptic_fixe_der_zero (double val, Param_elliptic &params, const Scalar &so, Scalar &uu) const
 General elliptic solver fixing the derivative at the origin and relaxing the continuity of the first derivative at the boundary between the nucleus and the first shell.
virtual void poisson2d (const Cmp &source_mat, const Cmp &source_quad, Param &par, Cmp &uu) const
 Computes the solution of a 2-D Poisson equation.
void sol_elliptic_2d (Param_elliptic &, const Scalar &, Scalar &) const
 General elliptic solver in a 2D case.
void sol_elliptic_pseudo_1d (Param_elliptic &, const Scalar &, Scalar &) const
 General elliptic solver in a pseudo 1d case.
virtual void dalembert (Param &par, Scalar &fJp1, const Scalar &fJ, const Scalar &fJm1, const Scalar &source) const
 Performs one time-step integration of the d'Alembert scalar equation.
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 (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 (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
 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
 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 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.
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).
void set_ori (double xa0, double ya0, double za0)
 Sets a new origin.
void set_rot_phi (double phi0)
 Sets a new rotation angle.

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

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.

Private Member Functions

void set_coord ()
 Assignment of the building functions to the member Coords.
virtual ostream & operator>> (ostream &) const
 Operator >>.

Private Attributes

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

Friends

Mtblmap_af_fait_r (const Map *)
Mtblmap_af_fait_tet (const Map *)
Mtblmap_af_fait_phi (const Map *)
Mtblmap_af_fait_sint (const Map *)
Mtblmap_af_fait_cost (const Map *)
Mtblmap_af_fait_sinp (const Map *)
Mtblmap_af_fait_cosp (const Map *)
Mtblmap_af_fait_x (const Map *)
Mtblmap_af_fait_y (const Map *)
Mtblmap_af_fait_z (const Map *)
Mtblmap_af_fait_xa (const Map *)
Mtblmap_af_fait_ya (const Map *)
Mtblmap_af_fait_za (const Map *)
Mtblmap_af_fait_xsr (const Map *)
Mtblmap_af_fait_dxdr (const Map *)
Mtblmap_af_fait_drdt (const Map *)
Mtblmap_af_fait_stdrdp (const Map *)
Mtblmap_af_fait_srdrdt (const Map *)
Mtblmap_af_fait_srstdrdp (const Map *)
Mtblmap_af_fait_sr2drdt (const Map *)
Mtblmap_af_fait_sr2stdrdp (const Map *)
Mtblmap_af_fait_d2rdx2 (const Map *)
Mtblmap_af_fait_lapr_tp (const Map *)
Mtblmap_af_fait_d2rdtdx (const Map *)
Mtblmap_af_fait_sstd2rdpdx (const Map *)
Mtblmap_af_fait_sr2d2rdt2 (const Map *)
ostream & operator<< (ostream &, const Map &)
 Operator <<.

Detailed Description

Affine radial mapping.

()

The affine radial mapping is the simplest one between the grid coordinates $(\xi, \theta', \phi')$ and the physical coordinates $(r, \theta, \phi)$. It is defined by $\theta=\theta'$, $\phi=\phi'$ and

Definition at line 2023 of file map.h.


Constructor & Destructor Documentation

Map_af::Map_af ( const Mg3d mgrille,
const double *  r_limits 
)

Standard Constructor.

Parameters:
mgrille [input] Multi-domain grid on which the mapping is defined
r_limits [input] Array (size: number of domains + 1) of the value of r at the boundaries of the various domains :

  • r_limits[l] : inner boundary of the domain no. l
  • r_limits[l+1] : outer boundary of the domain no. l

Definition at line 187 of file map_af.C.

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

Map_af::Map_af ( const Mg3d mgrille,
const Tbl r_limits 
)

Standard Constructor with Tbl.

Parameters:
mgrille [input] Multi-domain grid on which the mapping is defined
r_limits [input] Array (size: number of domains) of the value of r at the boundaries of the various domains :

  • r_limits[l] : inner boundary of the domain no. l
  • r_limits[l+1] : outer boundary of the domain no. l except in the last domain The last boundary is set to inifnity if the grid contains a compactified domain.

Definition at line 232 of file map_af.C.

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

Map_af::Map_af ( const Map_af mp  ) 

Copy constructor.

Definition at line 278 of file map_af.C.

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

Map_af::Map_af ( const Mg3d mgi,
FILE *  fd 
)

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

Definition at line 297 of file map_af.C.

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

Map_af::Map_af ( const Map mpi  )  [explicit]

Constructor from a general mapping.

If the input mapping belongs to the class Map_af , this constructor does the same job as the copy constructor Map_af(const Map_af& ) .

If the input mapping belongs to the class Map_et , this constructor sets in each domain, the values of $\alpha$ and $\beta$ to those of the Map_et .

Definition at line 313 of file map_af.C.

References alpha, beta, Map_et::get_alpha(), get_alpha(), Map_et::get_beta(), get_beta(), Mg3d::get_nzone(), Map::get_ori_x(), Map::get_ori_y(), Map::get_ori_z(), Map::get_rot_phi(), Map::mg, set_coord(), Map::set_ori(), and Map::set_rot_phi().

Map_af::~Map_af (  )  [virtual]

Destructor.

Definition at line 364 of file map_af.C.

References alpha, and beta.


Member Function Documentation

void Map_af::adapt ( const Cmp ent,
const Param par,
int  nbr = 0 
) [virtual]

Adaptation of the mapping to a given scalar field.

Implements Map.

Definition at line 669 of file map_af.C.

References c_est_pas_fait().

const Cmp& 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 803 of file map.h.

References Map::p_cmp_zero.

void Map_radial::comp_p_from_cartesian ( const Cmp v_x,
const Cmp v_y,
Cmp v_p 
) const [virtual, inherited]

Cmp version

Implements Map.

Definition at line 172 of file map_radial_comp_rtp.C.

References Map_radial::comp_p_from_cartesian().

void Map_radial::comp_p_from_cartesian ( const Scalar v_x,
const Scalar v_y,
Scalar v_p 
) const [virtual, inherited]

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

Definition at line 179 of file map_radial_comp_rtp.C.

References Scalar::check_dzpuis(), Scalar::dz_nonzero(), Scalar::get_dzpuis(), Scalar::get_etat(), Tensor::get_mp(), Scalar::get_spectral_va(), Valeur::mult_cp(), and Scalar::set_dzpuis().

void Map_radial::comp_r_from_cartesian ( const Cmp v_x,
const Cmp v_y,
const Cmp v_z,
Cmp v_r 
) const [virtual, inherited]

Cmp version

Implements Map.

Definition at line 61 of file map_radial_comp_rtp.C.

References Map_radial::comp_r_from_cartesian().

void Map_radial::comp_r_from_cartesian ( const Scalar v_x,
const Scalar v_y,
const Scalar v_z,
Scalar v_r 
) const [virtual, inherited]

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

Definition at line 68 of file map_radial_comp_rtp.C.

References Scalar::check_dzpuis(), Scalar::dz_nonzero(), Scalar::get_dzpuis(), Scalar::get_etat(), Tensor::get_mp(), Scalar::get_spectral_va(), Valeur::mult_cp(), Valeur::mult_ct(), Valeur::mult_sp(), Valeur::mult_st(), and Scalar::set_dzpuis().

void Map_radial::comp_t_from_cartesian ( const Cmp v_x,
const Cmp v_y,
const Cmp v_z,
Cmp v_t 
) const [virtual, inherited]

Cmp version

Implements Map.

Definition at line 117 of file map_radial_comp_rtp.C.

References Map_radial::comp_t_from_cartesian().

void Map_radial::comp_t_from_cartesian ( const Scalar v_x,
const Scalar v_y,
const Scalar v_z,
Scalar v_t 
) const [virtual, inherited]

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

Definition at line 124 of file map_radial_comp_rtp.C.

References Scalar::check_dzpuis(), Scalar::dz_nonzero(), Scalar::get_dzpuis(), Scalar::get_etat(), Tensor::get_mp(), Scalar::get_spectral_va(), Valeur::mult_cp(), Valeur::mult_ct(), Valeur::mult_sp(), Valeur::mult_st(), and Scalar::set_dzpuis().

void Map_radial::comp_x_from_spherical ( const Cmp v_r,
const Cmp v_theta,
const Cmp v_phi,
Cmp v_x 
) const [virtual, inherited]

Cmp version

Implements Map.

Definition at line 64 of file map_radial_comp_xyz.C.

References Map_radial::comp_x_from_spherical().

void Map_radial::comp_x_from_spherical ( const Scalar v_r,
const Scalar v_theta,
const Scalar v_phi,
Scalar v_x 
) const [virtual, inherited]

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

Definition at line 72 of file map_radial_comp_xyz.C.

References Scalar::check_dzpuis(), Scalar::dz_nonzero(), Scalar::get_dzpuis(), Scalar::get_etat(), Tensor::get_mp(), Scalar::get_spectral_va(), Valeur::mult_cp(), Valeur::mult_ct(), Valeur::mult_sp(), Valeur::mult_st(), and Scalar::set_dzpuis().

void Map_radial::comp_y_from_spherical ( const Cmp v_r,
const Cmp v_theta,
const Cmp v_phi,
Cmp v_y 
) const [virtual, inherited]

Cmp version

Implements Map.

Definition at line 122 of file map_radial_comp_xyz.C.

References Map_radial::comp_y_from_spherical().

void Map_radial::comp_y_from_spherical ( const Scalar v_r,
const Scalar v_theta,
const Scalar v_phi,
Scalar v_y 
) const [virtual, inherited]

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

Definition at line 131 of file map_radial_comp_xyz.C.

References Scalar::check_dzpuis(), Scalar::dz_nonzero(), Scalar::get_dzpuis(), Scalar::get_etat(), Tensor::get_mp(), Scalar::get_spectral_va(), Valeur::mult_cp(), Valeur::mult_ct(), Valeur::mult_sp(), Valeur::mult_st(), and Scalar::set_dzpuis().

void Map_radial::comp_z_from_spherical ( const Cmp v_r,
const Cmp v_theta,
Cmp v_z 
) const [virtual, inherited]

Cmp version

Implements Map.

Definition at line 180 of file map_radial_comp_xyz.C.

References Map_radial::comp_z_from_spherical().

void Map_radial::comp_z_from_spherical ( const Scalar v_r,
const Scalar v_theta,
Scalar v_z 
) const [virtual, inherited]

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

Definition at line 188 of file map_radial_comp_xyz.C.

References Scalar::check_dzpuis(), Scalar::dz_nonzero(), Scalar::get_dzpuis(), Scalar::get_etat(), Tensor::get_mp(), Scalar::get_spectral_va(), Valeur::mult_st(), and Scalar::set_dzpuis().

void 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 298 of file map.C.

References Map::ori_x, Map::ori_y, Map::ori_z, Map::rot_phi, and sqrt().

void Map_af::dalembert ( Param par,
Scalar fJp1,
const Scalar fJ,
const Scalar fJm1,
const Scalar source 
) const [virtual]

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

Parameters:
par [input/output] possible parameters to control the resolution of the d'Alembert equation: \ par.get_double(0) : [input] the time step dt ,\ par.get_int(0) : [input] the type of boundary conditions set at the outer boundary (0 : reflexion, 1 : Sommerfeld outgoing wave, valid only for l=0 components, 2 : Bayliss & Turkel outgoing wave, valid for l=0, 1, 2 components)\ par.get_int_mod(0) : [input/output] set to 0 at first call, is used as a working flag after (must not be modified after first call)\ par.get_int(1) : [input] (optional) if present, a shift of -1 is done in the multipolar spectrum in terms of $\ell$. The value of this variable gives the minimal value of (the shifted) $\ell$ for which the wave equation is solved.\ par.get_tensor_mod(0) : [input] (optional) if the wave equation is on a curved space-time, this is the potential in front of the Laplace operator. It has to be updated at every time-step (for a potential depending on time).\ Note: there are many other working objects attached to this Param , so one should not modify it.
fJp1 [output] solution $f^{J+1}$ at time J+1 with boundary conditions defined by par.get_int(0)
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$.

!!

Implements Map.

Definition at line 117 of file map_af_dalembert.C.

References Param::add_double_mod(), Param::add_map(), Param::add_tbl_mod(), alpha, Valeur::annule(), Tensor::annule_domain(), Tbl::annule_hard(), Scalar::annule_hard(), Mtbl_cf::annule_hard(), Valeur::base, beta, Valeur::c, Valeur::c_cf, Valeur::coef(), Valeur::coef_i(), Scalar::div_r(), Scalar::domain(), Scalar::dsdr(), Param::get_double(), Param::get_double_mod(), Scalar::get_dzpuis(), Valeur::get_etat(), Scalar::get_etat(), Param::get_int(), Param::get_int_mod(), Param::get_map(), Map::get_mg(), Tensor::get_mp(), Param::get_n_double(), Param::get_n_int(), Param::get_n_int_mod(), Param::get_n_tbl_mod(), Param::get_n_tensor_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Mg3d::get_radial(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Param::get_tbl_mod(), Param::get_tensor_mod(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Scalar::laplacian(), Map_af(), max(), Map::mg, Map::r, Tbl::set(), Mtbl_cf::set(), Valeur::set_base(), Valeur::set_base_t(), Scalar::set_domain(), Scalar::set_dzpuis(), Valeur::set_etat_cf_qcq(), Scalar::set_etat_qcq(), Tbl::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_grid_point(), Scalar::set_spectral_va(), sqrt(), Mg3d::std_base_scal(), Scalar::std_spectral_base(), T_LEG_PP, Scalar::val_grid_point(), Mtbl_cf::val_out_bound_jk(), val_r(), Valeur::ylm(), and Valeur::ylm_i().

void Map_radial::dec2_dzpuis ( Scalar ci  )  const [virtual, inherited]

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

Implements Map.

Definition at line 744 of file map_radial_r_manip.C.

References Valeur::annule(), Scalar::get_dzpuis(), Scalar::get_etat(), Valeur::get_mg(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_spectral_va(), and Map_radial::xsr.

void Map_radial::dec_dzpuis ( Scalar ci  )  const [virtual, inherited]

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

Implements Map.

Definition at line 639 of file map_radial_r_manip.C.

References Valeur::annule(), Scalar::get_dzpuis(), Scalar::get_etat(), Valeur::get_mg(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_spectral_va(), and Map_radial::xsr.

void Map_radial::div_cost ( Scalar ci  )  const [virtual, inherited]

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

Implements Map.

Definition at line 81 of file map_radial_th_manip.C.

References Scalar::get_etat(), Valeur::get_mg(), Map::mg, Valeur::scost(), Scalar::set_etat_qcq(), and Scalar::set_spectral_va().

void Map_radial::div_r ( Scalar ci  )  const [virtual, inherited]
void Map_radial::div_r_zec ( Scalar uu  )  const [virtual, inherited]
void Map_radial::div_rsint ( Scalar ci  )  const [virtual, inherited]
void Map_radial::div_sint ( Scalar ci  )  const [virtual, inherited]
void Map_radial::div_tant ( Scalar ci  )  const [virtual, inherited]
Param * Map_af::donne_para_poisson_vect ( Param par,
int  i 
) const [virtual]

Internal function intended to be used by Map::poisson_vect and Map::poisson_vect_oohara .

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

In the case of a Map_af the result is not used and the function only returns & par .

Implements Map.

Definition at line 66 of file map_poisson_vect.C.

void Map_af::dsdr ( const Scalar uu,
Scalar resu 
) const [virtual]

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

Note that in the compactified external domain (CED), the dzpuis flag of the output is 2 if the input has dzpuis = 0, and is increased by 1 in other cases.

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

Implements Map.

Definition at line 359 of file map_af_deriv.C.

References Valeur::annule(), Tensor::annule_domain(), Valeur::coef(), Scalar::dsdx(), Map_radial::dxdr, Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Map::mg, Valeur::mult_xm1_zec(), Valeur::set(), Scalar::set_dzpuis(), Scalar::set_etat_zero(), Scalar::set_spectral_base(), and Scalar::set_spectral_va().

void Map_af::dsdr ( const Cmp ci,
Cmp resu 
) const [virtual]

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

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

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

Implements Map.

Definition at line 275 of file map_af_deriv.C.

References Valeur::annule(), Valeur::base, Valeur::dsdx(), Map_radial::dxdr, Cmp::get_dzpuis(), Cmp::get_etat(), Map::get_mg(), Cmp::get_mp(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Valeur::mult_xm1_zec(), Coord::set(), Valeur::set(), Cmp::set_dzpuis(), Cmp::set_etat_zero(), Cmp::va, and Map_radial::xsr.

void Map_af::dsdradial ( const Scalar uu,
Scalar resu 
) const [virtual]

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

Note that in the compactified external domain (CED), the dzpuis flag of the output is 2 if the input has dzpuis = 0, and is increased by 1 in other cases.

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

Implements Map.

Definition at line 413 of file map_af_deriv.C.

References Valeur::annule(), Tensor::annule_domain(), Valeur::coef(), Scalar::dsdx(), Map_radial::dxdr, Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Map::mg, Valeur::mult_xm1_zec(), Valeur::set(), Scalar::set_dzpuis(), Scalar::set_etat_zero(), Scalar::set_spectral_base(), and Scalar::set_spectral_va().

void Map_af::dsdt ( const Scalar uu,
Scalar resu 
) const [virtual]

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

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

Implements Map.

Definition at line 793 of file map_af_deriv.C.

References Valeur::dsdt(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), and Scalar::set_etat_zero().

void Map_af::dsdxi ( const Scalar uu,
Scalar resu 
) const [virtual]

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

Note that in the compactified external domain (CED), the dzpuis flag of the output is 2 if the input has dzpuis = 0, and is increased by 1 in other cases.

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

Implements Map.

Definition at line 216 of file map_af_deriv.C.

References Valeur::annule(), Tensor::annule_domain(), Valeur::coef(), Scalar::dsdx(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Map::mg, Valeur::mult_xm1_zec(), Map::r, Valeur::set(), Scalar::set_dzpuis(), Scalar::set_etat_zero(), Scalar::set_spectral_base(), and Scalar::set_spectral_va().

void Map_af::dsdxi ( const Cmp ci,
Cmp resu 
) const [virtual]

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

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

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

Implements Map.

Definition at line 132 of file map_af_deriv.C.

References Valeur::annule(), Valeur::base, Valeur::dsdx(), Cmp::get_dzpuis(), Cmp::get_etat(), Map::get_mg(), Cmp::get_mp(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Valeur::mult_xm1_zec(), Coord::set(), Valeur::set(), Cmp::set_dzpuis(), Cmp::set_etat_zero(), Cmp::va, and Map_radial::xsr.

const Metric_flat & 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 327 of file map.C.

References Map::bvect_cart, and Map::p_flat_met_cart.

const Metric_flat & 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 317 of file map.C.

References Map::bvect_spher, and Map::p_flat_met_spher.

const double * Map_af::get_alpha (  )  const

Returns the pointer on the array alpha.

Definition at line 473 of file map_af.C.

References alpha.

const double * Map_af::get_beta (  )  const

Returns the pointer on the array beta.

Definition at line 477 of file map_af.C.

References beta.

const Base_vect_cart& 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 787 of file map.h.

References Map::bvect_cart.

const Base_vect_spher& 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 779 of file map.h.

References Map::bvect_spher.

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

Gives the Mg3d on which the mapping is defined.

Definition at line 761 of file map.h.

References Map::mg.

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

Returns the x coordinate of the origin.

Definition at line 764 of file map.h.

References Map::ori_x.

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

Returns the y coordinate of the origin.

Definition at line 766 of file map.h.

References Map::ori_y.

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

Returns the z coordinate of the origin.

Definition at line 768 of file map.h.

References Map::ori_z.

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

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

Definition at line 771 of file map.h.

References Map::rot_phi.

void Map_af::homothetie ( double  lambda  )  [virtual]

Sets a new radial scale.

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

Implements Map.

Definition at line 533 of file map_af.C.

References alpha, beta, Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, and Map_radial::reset_coord().

void Map_af::homothetie_interne ( double  lambda  ) 

Sets a new radial scale at the bondary between the nucleus and the first shell.

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

Definition at line 610 of file map_af.C.

References alpha, beta, and Map_radial::reset_coord().

void Map_radial::inc2_dzpuis ( Scalar ci  )  const [virtual, inherited]

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

Implements Map.

Definition at line 795 of file map_radial_r_manip.C.

References Valeur::annule(), Scalar::get_dzpuis(), Scalar::get_etat(), Valeur::get_mg(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_spectral_va(), and Map_radial::xsr.

void Map_radial::inc_dzpuis ( Scalar ci  )  const [virtual, inherited]

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

Implements Map.

Definition at line 692 of file map_radial_r_manip.C.

References Valeur::annule(), Scalar::get_dzpuis(), Scalar::get_etat(), Valeur::get_mg(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_spectral_va(), and Map_radial::xsr.

Tbl * Map_af::integrale ( const Cmp ci  )  const [virtual]

Computes the integral over all space of a Cmp.

The computed quantity is $\int u \, r^2 \sin\theta \, dr\, d\theta \, d\phi$. The routine allocates a Tbl (size: mg->nzone ) to store the result (partial integral) in each domain and returns a pointer to it.

Implements Map.

Definition at line 77 of file map_af_integ.C.

References alpha, Tbl::annule_hard(), Mtbl_cf::base, beta, Cmp::check_dzpuis(), Tbl::get_etat(), Mtbl_cf::get_etat(), Cmp::get_etat(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Map::mg, MSQ_P, MSQ_R, MSQ_T, P_COSSIN, P_COSSIN_P, R_CHEB, R_CHEBP, R_CHEBPI_P, R_CHEBPIM_P, R_CHEBU, R_JACO02, Tbl::set_etat_qcq(), Tbl::t, Mtbl_cf::t, T_COS, T_COS_P, T_COSSIN_C, T_COSSIN_CP, and Cmp::va.

double Map_af::integrale_surface ( const Scalar ci,
double  rayon 
) const
double Map_af::integrale_surface ( const Cmp ci,
double  rayon 
) const
double Map_af::integrale_surface_infini ( const Scalar ci  )  const
double Map_af::integrale_surface_infini ( const Cmp ci  )  const
void Map_af::lapang ( const Scalar uu,
Scalar lap 
) const [virtual]

Computes the angular Laplacian of a scalar field.

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

Implements Map.

Definition at line 544 of file map_af_lap.C.

References Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Valeur::lapang(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), and Valeur::ylm().

void Map_af::laplacien ( const Cmp uu,
int  zec_mult_r,
Cmp lap 
) const [virtual]
void Map_af::laplacien ( const Scalar uu,
int  zec_mult_r,
Scalar lap 
) const [virtual]

Computes the Laplacian of a scalar field.

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

Implements Map.

Definition at line 174 of file map_af_lap.C.

References Valeur::annule(), Scalar::annule(), Valeur::base, Scalar::check_dzpuis(), Valeur::coef(), Valeur::coef_i(), Valeur::d2sdx2(), Valeur::dsdx(), Map_radial::dxdr, Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Valeur::lapang(), Map::mg, Valeur::mult2_xm1_zec(), Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Valeur::sx(), Valeur::sx2(), Map_radial::xsr, and Valeur::ylm().

const Map_af & Map_af::mp_angu ( int  l_zone  )  const [virtual]

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

Valid only for the class Map_af.

Implements Map.

Definition at line 652 of file map_af.C.

References Mg3d::get_angu_1dom(), Map::get_mg(), Map_af(), Map::p_mp_angu, Tbl::set(), Tbl::set_etat_qcq(), and val_r_jk().

void Map_radial::mult_cost ( Scalar ci  )  const [virtual, inherited]

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

Implements Map.

Definition at line 57 of file map_radial_th_manip.C.

References Scalar::get_etat(), Valeur::get_mg(), Map::mg, Valeur::mult_ct(), Scalar::set_etat_qcq(), and Scalar::set_spectral_va().

void Map_radial::mult_r ( Cmp ci  )  const [virtual, inherited]

Multiplication by r of a Cmp.

In the CED, there is only a decrement of dzpuis

Implements Map.

Definition at line 215 of file map_radial_r_manip.C.

References Cmp::annule(), Valeur::base, Cmp::get_dzpuis(), Cmp::get_etat(), Valeur::get_mg(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Valeur::mult_x(), Cmp::set_dzpuis(), Cmp::va, and Map_radial::xsr.

void Map_radial::mult_r ( Scalar uu  )  const [virtual, inherited]
void Map_radial::mult_r_zec ( Scalar ci  )  const [virtual, inherited]

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

Implements Map.

Definition at line 292 of file map_radial_r_manip.C.

References Valeur::annule(), Valeur::base, Scalar::get_etat(), Valeur::get_mg(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Scalar::set_etat_qcq(), Scalar::set_spectral_va(), Valeur::sxm1_zec(), and Map_radial::xsr.

void Map_radial::mult_rsint ( Scalar ci  )  const [virtual, inherited]
void Map_radial::mult_sint ( Scalar ci  )  const [virtual, inherited]

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

Implements Map.

Definition at line 105 of file map_radial_th_manip.C.

References Scalar::get_etat(), Valeur::get_mg(), Map::mg, Valeur::mult_st(), Scalar::set_etat_qcq(), and Scalar::set_spectral_va().

void Map_af::operator= ( const Map_af mpi  )  [virtual]

Assignment to another affine mapping.

Definition at line 376 of file map_af.C.

References alpha, beta, Mg3d::get_nzone(), Map::mg, Map::ori_x, Map::ori_y, Map::ori_z, Map_radial::reset_coord(), Map::rot_phi, Map::set_ori(), and Map::set_rot_phi().

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

Comparison operator (egality).

bool Map_af::operator== ( const Map mpi  )  const [virtual]
ostream & Map_af::operator>> ( ostream &  ost  )  const [private, virtual]

Operator >>.

Implements Map.

Definition at line 499 of file map_af.C.

References alpha, beta, Mg3d::get_nr(), Mg3d::get_nzone(), Map::mg, Map::r, and val_r().

void Map_af::poisson ( const Cmp source,
Param par,
Cmp uu 
) const [virtual]

Computes the solution of a scalar Poisson equation.

Parameters:
source [input] source $\sigma$ of the Poisson equation $\Delta u = \sigma$.
par [] not used by this Map_af version.
uu [output] solution u with the boundary condition u =0 at spatial infinity.

Implements Map.

Definition at line 96 of file map_af_poisson.C.

References Valeur::c_cf, Cmp::check_dzpuis(), Valeur::coef(), Cmp::dz_nonzero(), Cmp::get_dzpuis(), Valeur::get_etat(), Cmp::get_etat(), Valeur::get_mg(), Map::get_mg(), Cmp::get_mp(), Map::mg, Cmp::set_dzpuis(), Cmp::set_etat_qcq(), Cmp::set_etat_zero(), and Cmp::va.

void Map_af::poisson2d ( const Cmp source_mat,
const Cmp source_quad,
Param par,
Cmp uu 
) const [virtual]

Computes the solution of a 2-D Poisson equation.

The 2-D Poisson equation writes

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

Parameters:
source_mat [input] Compactly supported part of the source $\sigma$ of the 2-D Poisson equation (typically matter terms)
source_quad [input] Non-compactly supported part of the source $\sigma$ of the 2-D Poisson equation (typically quadratic terms)
par [output] Parameter which contains the constant $\lambda$ used to fulfill the virial identity GRV2 : \ par.get_double_mod(0) : [output] constant lambda such that the source of the equation effectively solved is source_mat + lambda * source_quad .
uu [input/output] solution u with the boundary condition u =0 at spatial infinity.

Implements Map.

Definition at line 76 of file map_af_poisson2d.C.

References Cmp::allocate_all(), alpha, Valeur::base, beta, Valeur::c, Cmp::check_dzpuis(), Valeur::coef_i(), Param::get_double_mod(), Tbl::get_etat(), Mtbl::get_etat(), Valeur::get_etat(), Cmp::get_etat(), Mg3d::get_grille3d(), Map::get_mg(), Cmp::get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Mg3d::get_type_r(), Mg3d::get_type_t(), Map::mg, Cmp::set(), Cmp::set_dzpuis(), Tbl::t, Mtbl::t, T_COS, T_COS_P, T_SIN, T_SIN_I, Cmp::va, and Grille3d::x.

void Map_af::poisson_angu ( const Scalar source,
Param par,
Scalar uu,
double  lambda = 0 
) const [virtual]

Computes the solution of the generalized angular Poisson equation.

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

Parameters:
source [input] source $\sigma$ of the equation $\Delta_{\theta\varphi} u + \lambda u = \sigma$.
par [input/output] possible parameters to control the resolution of the Poisson equation. See the actual implementation in the derived class of Map for documentation.
uu [input/output] solution u
lambda [input] coefficient $\lambda$ in the above equation (default value = 0)

Implements Map.

Definition at line 54 of file map_af_poisson_angu.C.

References Valeur::c_cf, Valeur::coef(), Scalar::get_dzpuis(), Valeur::get_etat(), Scalar::get_etat(), Tensor::get_mp(), Param::get_n_int(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_zero(), Valeur::ylm(), and Valeur::ylm_i().

void Map_radial::poisson_compact ( int  nzet,
const Cmp source,
const Cmp aa,
const Tenseur bb,
const Param par,
Cmp psi 
) const [virtual, inherited]

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

Definition at line 449 of file map_radial_poisson_cpt.C.

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

void Map_radial::poisson_compact ( const Cmp source,
const Cmp aa,
const Tenseur bb,
const Param par,
Cmp psi 
) const [virtual, inherited]

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

Definition at line 151 of file map_radial_poisson_cpt.C.

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

void Map_af::poisson_frontiere ( const Cmp source,
const Valeur limite,
int  type_raccord,
int  num_front,
Cmp pot,
double  fact_dir = 0.,
double  fact_neu = 0. 
) const [virtual]
void Map_af::poisson_frontiere_double ( const Cmp source,
const Valeur lim_func,
const Valeur lim_der,
int  num_zone,
Cmp pot 
) const [virtual]

Solver of the Poisson equation with boundary condition for the affine mapping case, cases with boundary conditions of both Dirichlet and Neumann type (no condition imposed at infinity).

Implements Map.

Definition at line 204 of file cmp_pde_frontiere.C.

References Valeur::base, Valeur::c_cf, Valeur::coef(), Mg3d::get_angu(), Valeur::get_etat(), Cmp::get_etat(), Valeur::get_mg(), Map::get_mg(), Cmp::get_mp(), Map::mg, Cmp::set_etat_qcq(), Cmp::set_etat_zero(), Cmp::va, and Valeur::ylm().

void Map_af::poisson_interne ( const Cmp source,
const Valeur limite,
Param par,
Cmp pot 
) const [virtual]

Computes the solution of a Poisson equation in the shell, imposing a boundary condition at the surface of the star.

Parameters:
source [input] : source of the equation.
limite [input] : limite[num_front] contains the angular function being the boudary condition.
par [input] : parameters of the computation.
pot [output] : result.

Implements Map.

Definition at line 470 of file cmp_pde_frontiere.C.

References Valeur::base, Valeur::c_cf, Valeur::coef(), Mg3d::get_angu(), Valeur::get_etat(), Cmp::get_etat(), Valeur::get_mg(), Map::get_mg(), Cmp::get_mp(), Map::mg, Cmp::set_dzpuis(), Cmp::set_etat_qcq(), Cmp::set_etat_zero(), Cmp::va, and Valeur::ylm().

void Map_af::poisson_regular ( const Cmp source,
int  k_div,
int  nzet,
double  unsgam1,
Param par,
Cmp uu,
Cmp uu_regu,
Cmp uu_div,
Tenseur duu_div,
Cmp source_regu,
Cmp source_div 
) const [virtual]

Computes the solution of a scalar Poisson equation.

The regularized source $\sigma_{\rm regu} = \sigma - \sigma_{\rm div}$ is constructed and solved.

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 [] not used by this Map_af version.
uu [output] solution u with the boundary condition u =0 at spatial infinity.
uu_regu [output] solution of the regular part of the source.
uu_div [output] solution of the diverging part of the source.
duu_div [output] derivative of the diverging potential
source_regu [output] regularized source
source_div [output] diverging part of the source

Implements Map.

Definition at line 104 of file map_af_poisson_regu.C.

References alpha, Cmp::annule(), Base_val::b, Mtbl_cf::base, Valeur::base, Valeur::c_cf, Valeur::coef(), dsdt(), Map_radial::dxdr, Valeur::get_etat(), Cmp::get_etat(), Mg3d::get_grille3d(), Valeur::get_mg(), Map::get_mg(), Cmp::get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Matrice::inverse(), Map::mg, MSQ_P, MSQ_R, pow(), R_CHEBPIM_I, Tenseur::set(), Mtbl_cf::set(), Matrice::set(), Tbl::set(), Matrice::set_band(), Cmp::set_dzpuis(), Tenseur::set_etat_qcq(), Cmp::set_etat_qcq(), Matrice::set_etat_qcq(), Tbl::set_etat_qcq(), Matrice::set_lu(), Tenseur::set_triad(), Mg3d::std_base_scal(), stdsdp(), T_LEG_P, Cmp::va, Grille3d::x, and Valeur::ylm_i().

void Map_af::poisson_tau ( const Cmp source,
Param par,
Cmp uu 
) const [virtual]

Computes the solution of a scalar Poisson equation using a Tau method.

Parameters:
source [input] source $\sigma$ of the Poisson equation $\Delta u = \sigma$.
par [] not used by this Map_af version.
uu [output] solution u with the boundary condition u =0 at spatial infinity.

Implements Map.

Definition at line 161 of file map_af_poisson.C.

References Valeur::c_cf, Cmp::check_dzpuis(), Valeur::coef(), Cmp::dz_nonzero(), Cmp::get_dzpuis(), Valeur::get_etat(), Cmp::get_etat(), Valeur::get_mg(), Map::get_mg(), Cmp::get_mp(), Map::mg, Cmp::set_dzpuis(), Cmp::set_etat_qcq(), Cmp::set_etat_zero(), and Cmp::va.

void Map_af::primr ( const Scalar uu,
Scalar resu,
bool  null_infty 
) const [virtual]

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

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

Parameters:
uu [input] function f (must have a dzpuis = 2)
resu [input] function F
null_infty if true (default), the primitive is null at infinity (or on the grid boundary). F is null at the center otherwise

Implements Map.

Definition at line 70 of file map_af_primr.C.

References alpha, Base_val::b, Mtbl_cf::base, Valeur::c_cf, Scalar::check_dzpuis(), Valeur::coef(), Base_val::get_b(), Valeur::get_base(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), MAX_BASE, Map::mg, MSQ_R, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, Valeur::set(), Tbl::set(), Valeur::set_base(), Valeur::set_etat_cf_qcq(), Mtbl_cf::set_etat_qcq(), Scalar::set_etat_qcq(), Tbl::set_etat_zero(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Mtbl_cf::t, TRA_R, and Mtbl_cf::val_out_bound_jk().

virtual void 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 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 Map_radial::reevaluate ( const Map mp_prev,
int  nzet,
Scalar uu 
) const [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 Scalar at the grid points defined by *this.

Definition at line 169 of file map_radial_reevaluate.C.

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

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

Definition at line 54 of file map_radial_reevaluate.C.

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

virtual void 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 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 Map_radial::reevaluate_symy ( const Map mp_prev,
int  nzet,
Scalar uu 
) const [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 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 189 of file map_radial_reeval_symy.C.

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

void Map_radial::reevaluate_symy ( const Map mp_prev,
int  nzet,
Cmp uu 
) const [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 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 55 of file map_radial_reeval_symy.C.

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

void Map_radial::reset_coord (  )  [protected, virtual, inherited]
void Map_af::resize ( int  l,
double  lambda 
) [virtual]

Rescales the outer boundary of one domain.

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

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

Implements Map.

Definition at line 556 of file map_af.C.

References alpha, beta, Mg3d::get_type_r(), Map::mg, and Map_radial::reset_coord().

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

Save in a file.

Reimplemented from Map_radial.

Definition at line 485 of file map_af.C.

References alpha, beta, fwrite_be(), Mg3d::get_nzone(), and Map::mg.

void Map_af::set_alpha ( double  alpha0,
int  l 
)

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

Definition at line 626 of file map_af.C.

References alpha, and Map_radial::reset_coord().

void Map_af::set_beta ( double  beta0,
int  l 
)

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

Definition at line 637 of file map_af.C.

References beta, and Map_radial::reset_coord().

void Map_af::set_coord (  )  [private]
void Map::set_ori ( double  xa0,
double  ya0,
double  za0 
) [inherited]

Sets a new origin.

Definition at line 249 of file map.C.

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

void Map::set_rot_phi ( double  phi0  )  [inherited]
void Map_af::sol_elliptic ( Param_elliptic params,
const Scalar so,
Scalar uu 
) const

General elliptic solver.

The field is zero at infinity.

Parameters:
params [input] : the operators and variables to be uses.
so [input] : the source.
uu [output] : the solution.

Definition at line 96 of file map_af_elliptic.C.

References Scalar::check_dzpuis(), Valeur::coef(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Param_elliptic::var_F, Param_elliptic::var_G, Valeur::ylm(), and Valeur::ylm_i().

void Map_af::sol_elliptic_2d ( Param_elliptic ope_var,
const Scalar source,
Scalar pot 
) const

General elliptic solver in a 2D case.

The field is zero at infinity.

Parameters:
params [input] : the operators and variables to be uses.
so [input] : the source.
uu [output] : the solution.

Definition at line 53 of file map_af_elliptic_2d.C.

References Scalar::check_dzpuis(), Valeur::coef(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Param_elliptic::var_F, and Param_elliptic::var_G.

void Map_af::sol_elliptic_boundary ( Param_elliptic params,
const Scalar so,
Scalar uu,
const Scalar bound,
double  fact_dir,
double  fact_neu 
) const
void Map_af::sol_elliptic_boundary ( Param_elliptic params,
const Scalar so,
Scalar uu,
const Mtbl_cf bound,
double  fact_dir,
double  fact_neu 
) const

General elliptic solver including inner boundary conditions.

The field is zero at infinity.

Parameters:
params [input] : the operators and variables to be uses.
so [input] : the source.
uu [output] : the solution.
bound [input] : the boundary condition
fact_dir : 1 Dirchlet condition, 0 Neumann condition
fact_neu : 0 Dirchlet condition, 1 Neumann condition

Definition at line 155 of file map_af_elliptic.C.

References Scalar::check_dzpuis(), Valeur::coef(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Param_elliptic::var_F, Param_elliptic::var_G, Valeur::ylm(), and Valeur::ylm_i().

void Map_af::sol_elliptic_fixe_der_zero ( double  val,
Param_elliptic params,
const Scalar so,
Scalar uu 
) const

General elliptic solver fixing the derivative at the origin and relaxing the continuity of the first derivative at the boundary between the nucleus and the first shell.

Parameters:
val [input] : valeur of the derivative.
params [input] : the operators and variables to be uses.
so [input] : the source.
uu [output] : the solution.

Definition at line 479 of file map_af_elliptic.C.

References alpha, Scalar::check_dzpuis(), Valeur::coef(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Param_elliptic::var_F, Param_elliptic::var_G, Valeur::ylm(), and Valeur::ylm_i().

void Map_af::sol_elliptic_no_zec ( Param_elliptic params,
const Scalar so,
Scalar uu,
double  val 
) const

General elliptic solver.

The equation is not solved in the compactified domain.

Parameters:
params [input] : the operators and variables to be uses.
so [input] : the source.
uu [output] : the solution.
val [input] : value at the last shell.

Definition at line 308 of file map_af_elliptic.C.

References Scalar::check_dzpuis(), Valeur::coef(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Param_elliptic::var_F, Param_elliptic::var_G, Valeur::ylm(), and Valeur::ylm_i().

void Map_af::sol_elliptic_only_zec ( Param_elliptic params,
const Scalar so,
Scalar uu,
double  val 
) const

General elliptic solver.

The equation is solved only in the compactified domain.

Parameters:
params [input] : the operators and variables to be uses.
so [input] : the source.
uu [output] : the solution.
val [input] : value at the inner boundary.

Definition at line 364 of file map_af_elliptic.C.

References Scalar::check_dzpuis(), Valeur::coef(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Param_elliptic::var_F, Param_elliptic::var_G, Valeur::ylm(), and Valeur::ylm_i().

void Map_af::sol_elliptic_pseudo_1d ( Param_elliptic ope_var,
const Scalar source,
Scalar pot 
) const

General elliptic solver in a pseudo 1d case.

The field is zero at infinity.

Parameters:
params [input] : the operators and variables to be uses.
so [input] : the source.
uu [output] : the solution.

Definition at line 54 of file map_af_elliptic_pseudo_1d.C.

References Scalar::check_dzpuis(), Valeur::coef(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Valeur::val_propre_1d(), Valeur::val_propre_1d_i(), Param_elliptic::var_F, and Param_elliptic::var_G.

void Map_af::sol_elliptic_sin_zec ( Param_elliptic params,
const Scalar so,
Scalar uu,
double *  coefs,
double *  phases 
) const

General elliptic solver.

The equation is not solved in the compactified domain and the matching is done with an homogeneous solution.

Definition at line 422 of file map_af_elliptic.C.

References Scalar::check_dzpuis(), Valeur::coef(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_va(), Param_elliptic::var_F, Param_elliptic::var_G, Valeur::ylm(), and Valeur::ylm_i().

void Map_af::srdsdt ( const Scalar uu,
Scalar resu 
) const [virtual]

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

Note that in the compactified external domain (CED), the dzpuis flag of the output is 2 if the input has dzpuis = 0, and is increased by 1 in other cases.

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

Implements Map.

Definition at line 565 of file map_af_deriv.C.

References Valeur::annule(), Valeur::coef(), Valeur::get_base(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_zero(), Valeur::sx(), and Map_radial::xsr.

void Map_af::srdsdt ( const Cmp ci,
Cmp resu 
) const [virtual]

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

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

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

Implements Map.

Definition at line 471 of file map_af_deriv.C.

References Valeur::annule(), Valeur::dsdt(), Valeur::get_base(), Cmp::get_dzpuis(), Cmp::get_etat(), Map::get_mg(), Cmp::get_mp(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Valeur::set_base(), Cmp::set_dzpuis(), Cmp::set_etat_zero(), Valeur::sx(), Cmp::va, and Map_radial::xsr.

void Map_af::srstdsdp ( const Scalar uu,
Scalar resu 
) const [virtual]

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

Note that in the compactified external domain (CED), the dzpuis flag of the output is 2 if the input has dzpuis = 0, and is increased by 1 in other cases.

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

Implements Map.

Definition at line 725 of file map_af_deriv.C.

References Valeur::annule(), Valeur::coef(), Valeur::get_base(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_zero(), Valeur::ssint(), Valeur::sx(), and Map_radial::xsr.

void Map_af::srstdsdp ( const Cmp ci,
Cmp resu 
) const [virtual]

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

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

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

Implements Map.

Definition at line 629 of file map_af_deriv.C.

References Valeur::annule(), Valeur::dsdp(), Valeur::get_base(), Cmp::get_dzpuis(), Cmp::get_etat(), Map::get_mg(), Cmp::get_mp(), Mg3d::get_nzone(), Mg3d::get_type_r(), Map::mg, Valeur::set_base(), Cmp::set_dzpuis(), Cmp::set_etat_zero(), Valeur::ssint(), Valeur::sx(), Cmp::va, and Map_radial::xsr.

void Map_af::stdsdp ( const Scalar uu,
Scalar resu 
) const [virtual]

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

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

Implements Map.

Definition at line 819 of file map_af_deriv.C.

References Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Tensor::get_mp(), Scalar::get_spectral_va(), Map::mg, Scalar::set_dzpuis(), Scalar::set_etat_zero(), and Valeur::stdsdp().

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

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

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

Implements Map.

Definition at line 188 of file map_af_radius.C.

References val_lx().

void Map_af::val_lx ( double  rr,
double  theta,
double  pphi,
int &  l,
double &  xi 
) const [virtual]

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

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

Implements Map.

Definition at line 124 of file map_af_radius.C.

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

void Map_af::val_lx_jk ( double  rr,
int  j,
int  k,
const Param par,
int &  l,
double &  xi 
) const [virtual]

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

Parameters:
rr [input] value of r
j [input] index of the collocation point in $\theta$
k [input] index of the collocation point in $\phi$
par [] unused by the Map_af version
l [output] value of the domain index
xi [output] value of $\xi$

Implements Map_radial.

Definition at line 211 of file map_af_radius.C.

References val_lx().

double Map_af::val_r ( int  l,
double  xi,
double  theta,
double  pphi 
) const [virtual]

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

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

Implements Map.

Definition at line 92 of file map_af_radius.C.

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

double Map_af::val_r_jk ( int  l,
double  xi,
int  j,
int  k 
) const [virtual]

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

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

Implements Map_radial.

Definition at line 201 of file map_af_radius.C.

References val_r().


Friends And Related Function Documentation

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

Operator <<.


Member Data Documentation

double* Map_af::alpha [private]

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

Definition at line 2029 of file map.h.

double* Map_af::beta [private]

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

Definition at line 2031 of file map.h.

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

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

Coord Map::cosp [inherited]

$\cos\phi$

Definition at line 720 of file map.h.

Coord Map::cost [inherited]

$\cos\theta$

Definition at line 718 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 1636 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 1615 of file map.h.

Coord Map_radial::drdt [inherited]

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

Definition at line 1564 of file map.h.

Coord Map_radial::dxdr [inherited]

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

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

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

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

Definition at line 672 of file map.h.

double Map::ori_x [protected, inherited]

Absolute coordinate x of the origin.

Definition at line 674 of file map.h.

double Map::ori_y [protected, inherited]

Absolute coordinate y of the origin.

Definition at line 675 of file map.h.

double Map::ori_z [protected, inherited]

Absolute coordinate z of the origin.

Definition at line 676 of file map.h.

Cmp* 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 709 of file map.h.

Metric_flat* 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 703 of file map.h.

Metric_flat* 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 698 of file map.h.

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

Pointer on the "angular" mapping.

Definition at line 711 of file map.h.

Coord Map::phi [inherited]

$\phi$ coordinate centered on the grid

Definition at line 716 of file map.h.

Coord Map::r [inherited]

r coordinate centered on the grid

Definition at line 714 of file map.h.

double Map::rot_phi [protected, inherited]

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

Definition at line 677 of file map.h.

Coord Map::sinp [inherited]

$\sin\phi$

Definition at line 719 of file map.h.

Coord Map::sint [inherited]

$\sin\theta$

Definition at line 717 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 1653 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 1596 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 1604 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 1580 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 1588 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 1644 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 1572 of file map.h.

Coord Map::tet [inherited]

$\theta$ coordinate centered on the grid

Definition at line 715 of file map.h.

Coord Map::x [inherited]

x coordinate centered on the grid

Definition at line 722 of file map.h.

Coord Map::xa [inherited]

Absolute x coordinate.

Definition at line 726 of file map.h.

Coord Map_radial::xsr [inherited]

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

Definition at line 1545 of file map.h.

Coord Map::y [inherited]

y coordinate centered on the grid

Definition at line 723 of file map.h.

Coord Map::ya [inherited]

Absolute y coordinate.

Definition at line 727 of file map.h.

Coord Map::z [inherited]

z coordinate centered on the grid

Definition at line 724 of file map.h.

Coord Map::za [inherited]

Absolute z coordinate.

Definition at line 728 of file map.h.


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

Generated on 7 Oct 2014 for LORENE by  doxygen 1.6.1