|
LORENE
|
Affine radial mapping. More...
#include <map.h>
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_af & | mp_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 ![]() | |
| virtual void | val_lx (double rr, double theta, double pphi, int &l, double &xi) const |
Computes the domain index l and the value of ![]() ![]() | |
| 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 ![]() ![]() | |
| 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 ![]() ![]() | |
| 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 ![]() ![]() | |
| 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 ![]() | |
| void | set_beta (double beta0, int l) |
Modifies the value of ![]() | |
| virtual void | dsdxi (const Cmp &ci, Cmp &resu) const |
Computes ![]() Cmp. | |
| virtual void | dsdr (const Cmp &ci, Cmp &resu) const |
Computes ![]() Cmp. | |
| virtual void | srdsdt (const Cmp &ci, Cmp &resu) const |
Computes ![]() Cmp. | |
| virtual void | srstdsdp (const Cmp &ci, Cmp &resu) const |
Computes ![]() Cmp. | |
| virtual void | dsdr (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar. | |
| virtual void | dsdxi (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar. | |
| virtual void | dsdradial (const Scalar &, Scalar &) const |
Computes ![]() Scalar. | |
| virtual void | srdsdt (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar. | |
| virtual void | srstdsdp (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar. | |
| virtual void | dsdt (const Scalar &uu, Scalar &resu) const |
Computes ![]() Scalar. | |
| virtual void | stdsdp (const Scalar &uu, Scalar &resu) const |
Computes ![]() 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 ![]() | |
| virtual Tbl * | integrale (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 Param * | donne_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 ¶ms, const Scalar &so, Scalar &uu) const |
| General elliptic solver. | |
| void | sol_elliptic_boundary (Param_elliptic ¶ms, 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 ¶ms, 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 ¶ms, const Scalar &so, Scalar &uu, double val) const |
| General elliptic solver. | |
| void | sol_elliptic_only_zec (Param_elliptic ¶ms, const Scalar &so, Scalar &uu, double val) const |
| General elliptic solver. | |
| void | sol_elliptic_sin_zec (Param_elliptic ¶ms, const Scalar &so, Scalar &uu, double *coefs, double *) const |
| General elliptic solver. | |
| void | sol_elliptic_fixe_der_zero (double val, Param_elliptic ¶ms, 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 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_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 | 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 ![]() Scalar. | |
| virtual void | div_rsint (Scalar &) const |
Division by ![]() 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 ![]() Scalar. | |
| virtual void | div_cost (Scalar &) const |
Division by ![]() Scalar. | |
| virtual void | mult_sint (Scalar &) const |
Multiplication by ![]() Scalar. | |
| virtual void | div_sint (Scalar &) const |
Division by ![]() Scalar. | |
| virtual void | div_tant (Scalar &) const |
Division by ![]() 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 ![]() 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 ![]() 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 ![]() | |
| 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 ![]() | |
| const Mg3d * | get_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_spher & | get_bvect_spher () const |
Returns the orthonormal vectorial basis ![]() ![]() | |
| const Base_vect_cart & | get_bvect_cart () const |
Returns the Cartesian basis ![]() | |
| const Metric_flat & | flat_met_spher () const |
Returns the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher. | |
| const Metric_flat & | flat_met_cart () const |
Returns the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart. | |
| const Cmp & | cmp_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 ![]() | |
| 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 |
![]() ![]() | |
| Coord | dxdr |
![]() ![]() | |
| Coord | drdt |
![]() ![]() | |
| Coord | stdrdp |
![]() ![]() | |
| Coord | srdrdt |
![]() ![]() | |
| Coord | srstdrdp |
![]() ![]() | |
| Coord | sr2drdt |
![]() ![]() | |
| Coord | sr2stdrdp |
![]() ![]() | |
| Coord | d2rdx2 |
![]() ![]() | |
| Coord | lapr_tp |
![]() ![]() | |
| Coord | d2rdtdx |
![]() ![]() | |
| Coord | sstd2rdpdx |
![]() ![]() | |
| Coord | sr2d2rdt2 |
![]() ![]() | |
| Coord | r |
| r coordinate centered on the grid | |
| Coord | tet |
![]() | |
| Coord | phi |
![]() | |
| Coord | sint |
![]() | |
| Coord | cost |
![]() | |
| Coord | sinp |
![]() | |
| Coord | cosp |
![]() | |
| 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 Mg3d * | mg |
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 ![]() ![]() | |
| Base_vect_cart | bvect_cart |
Cartesian basis ![]() | |
| Metric_flat * | p_flat_met_spher |
Pointer onto the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher. | |
| Metric_flat * | p_flat_met_cart |
Pointer onto the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart. | |
| Cmp * | p_cmp_zero |
| The null Cmp. | |
| Map_af * | p_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 ![]() | |
| double * | beta |
Array (size: mg->nzone ) of the values of ![]() | |
Affine radial mapping.
()
The affine radial mapping is the simplest one between the grid coordinates 







Standard Constructor.
| 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 :
|
Definition at line 191 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and set_coord().
Standard Constructor with Tbl.
| 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 :
|
Definition at line 236 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and set_coord().
Copy constructor.
Definition at line 282 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and set_coord().
Constructor from a file (see sauve(FILE*) )
Definition at line 301 of file map_af.C.
References alpha, beta, Lorene::fread_be(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and set_coord().
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 

Map_et .
Definition at line 317 of file map_af.C.
References alpha, beta, get_alpha(), Lorene::Map_et::get_alpha(), get_beta(), Lorene::Map_et::get_beta(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_x(), Lorene::Map::get_ori_y(), Lorene::Map::get_ori_z(), Lorene::Map::get_rot_phi(), Lorene::Map::mg, set_coord(), Lorene::Map::set_ori(), and Lorene::Map::set_rot_phi().
|
virtual |
Adaptation of the mapping to a given scalar field.
Implements Lorene::Map.
Definition at line 673 of file map_af.C.
References Lorene::c_est_pas_fait().
|
virtualinherited |
Cmp version
Implements Lorene::Map.
Definition at line 176 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_p_from_cartesian().
|
virtualinherited |
Computes the Spherical 
bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart .
| v_x | [input] x-component of the vector |
| v_y | [input] y-component of the vector |
| v_p | [output] ![]() |
Implements Lorene::Map.
Definition at line 183 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_sp(), and Lorene::Scalar::set_dzpuis().
|
virtualinherited |
Cmp version
Implements Lorene::Map.
Definition at line 65 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_r_from_cartesian().
|
virtualinherited |
Computes the Spherical r component (with respect to bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart .
| 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 72 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().
|
virtualinherited |
Cmp version
Implements Lorene::Map.
Definition at line 121 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_t_from_cartesian().
|
virtualinherited |
Computes the Spherical 
bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart .
| 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] ![]() |
Implements Lorene::Map.
Definition at line 128 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().
|
virtualinherited |
Cmp version
Implements Lorene::Map.
Definition at line 68 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_x_from_spherical().
|
virtualinherited |
Computes the Cartesian x component (with respect to bvect_cart) of a vector given by its spherical components with respect to bvect_spher.
| v_r | [input] r -component of the vector |
| v_theta | [input] ![]() |
| v_phi | [input] ![]() |
| v_x | [output] x-component of the vector |
Implements Lorene::Map.
Definition at line 76 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().
|
virtualinherited |
Cmp version
Implements Lorene::Map.
Definition at line 126 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_y_from_spherical().
|
virtualinherited |
Computes the Cartesian y component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher .
| v_r | [input] r -component of the vector |
| v_theta | [input] ![]() |
| v_phi | [input] ![]() |
| v_y | [output] y-component of the vector |
Implements Lorene::Map.
Definition at line 135 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().
|
virtualinherited |
Cmp version
Implements Lorene::Map.
Definition at line 184 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_z_from_spherical().
|
virtualinherited |
Computes the Cartesian z component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher .
| v_r | [input] r -component of the vector |
| v_theta | [input] ![]() |
| v_z | [output] z-component of the vector |
Implements Lorene::Map.
Definition at line 192 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_ct(), Lorene::Valeur::mult_st(), and Lorene::Scalar::set_dzpuis().
|
inherited |
Determines the coordinates 
| 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 ![]() |
| pphi | [output] value of ![]() |
Definition at line 302 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 |
Performs one time-step integration of the d'Alembert scalar equation.
| 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 ![]() ![]() 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 ![]() par.get_int(0) |
| fJ | [input] solution ![]() |
| fJm1 | [input] solution ![]() |
| source | [input] source ![]() ![]() |
!!
Implements Lorene::Map.
Definition at line 122 of file map_af_dalembert.C.
References Lorene::Param::add_double_mod(), Lorene::Param::add_map(), Lorene::Param::add_tbl_mod(), alpha, Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Mtbl_cf::annule_hard(), Lorene::Scalar::annule_hard(), Lorene::Tbl::annule_hard(), Lorene::Valeur::base, beta, Lorene::Valeur::c, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Valeur::coef_i(), Lorene::Scalar::div_r(), Lorene::Scalar::domain(), Lorene::Scalar::dsdr(), Lorene::Base_val::dsdx(), Lorene::Param::get_double(), Lorene::Param::get_double_mod(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Param::get_map(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Param::get_n_double(), Lorene::Param::get_n_int(), Lorene::Param::get_n_int_mod(), Lorene::Param::get_n_tbl_mod(), Lorene::Param::get_n_tensor_mod(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_radial(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Param::get_tbl_mod(), Lorene::Param::get_tensor_mod(), Lorene::Mg3d::get_type_r(), Lorene::Base_val::give_quant_numbers(), Lorene::Scalar::laplacian(), Lorene::max(), Lorene::Map::mg, Lorene::Map::r, Lorene::Tbl::set(), Lorene::Mtbl_cf::set(), Lorene::Valeur::set_base(), Lorene::Valeur::set_base_t(), Lorene::Scalar::set_domain(), Lorene::Scalar::set_dzpuis(), Lorene::Valeur::set_etat_cf_qcq(), Lorene::Scalar::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_grid_point(), Lorene::Scalar::set_spectral_va(), Lorene::sqrt(), Lorene::Mg3d::std_base_scal(), Lorene::Scalar::std_spectral_base(), T_LEG_PP, Lorene::Scalar::val_grid_point(), Lorene::Mtbl_cf::val_out_bound_jk(), val_r(), Lorene::Base_val::ylm(), Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
Decreases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 748 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult2_xm1_zec(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Map_radial::xsr.
Decreases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED).
Implements Lorene::Map.
Definition at line 643 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Map_radial::xsr.
Division by 
Scalar.
Implements Lorene::Map.
Definition at line 85 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Valeur::scost(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
Division by r of a Scalar.
Implements Lorene::Map.
Definition at line 514 of file map_radial_r_manip.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sx(), and Lorene::Map_radial::xsr.
Division by r (in the compactified external domain only) of a Scalar.
Implements Lorene::Map.
Definition at line 585 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Map_radial::xsr.
Division by 
Scalar.
Implements Lorene::Map.
Definition at line 434 of file map_radial_r_manip.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::ssint(), Lorene::Valeur::sx(), and Lorene::Map_radial::xsr.
Division by 
Scalar.
Implements Lorene::Map.
Definition at line 133 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Valeur::ssint().
Division by 
Scalar.
Implements Lorene::Map.
Definition at line 158 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Valeur::mult_ct(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), and Lorene::Valeur::ssint().
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 Lorene::Map.
Definition at line 70 of file map_poisson_vect.C.
Computes 
Cmp.
Note that in the compactified external domain (CED), it computes 
| ci | [input] field to consider |
| resu | [output] derivative of ci |
Implements Lorene::Map.
Definition at line 279 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Valeur::dsdx(), Lorene::Map_radial::dxdr, Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Map::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_zero(), Lorene::Cmp::va, and Lorene::Map_radial::xsr.
Computes 
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.
| uu | [input] field to consider |
| resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 363 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::coef(), Lorene::Valeur::dsdx(), Lorene::Map_radial::dxdr, Lorene::Valeur::get_base(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_base(), and Lorene::Scalar::set_spectral_va().
Computes 
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.
| uu | [input] field to consider |
| resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 417 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::coef(), Lorene::Valeur::dsdx(), Lorene::Map_radial::dxdr, Lorene::Valeur::get_base(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_base(), and Lorene::Scalar::set_spectral_va().
Computes 
Scalar.
| uu | [input] scalar field |
| resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 797 of file map_af_deriv.C.
References Lorene::Valeur::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), and Lorene::Scalar::set_etat_zero().
Computes 
Cmp.
Note that in the compactified external domain (CED), it computes 
| ci | [input] field to consider |
| resu | [output] derivative of ci |
Implements Lorene::Map.
Definition at line 136 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Valeur::dsdx(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Map::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_zero(), Lorene::Cmp::va, and Lorene::Map_radial::xsr.
Computes 
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.
| uu | [input] field to consider |
| resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 220 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::coef(), Lorene::Valeur::dsdx(), Lorene::Valeur::get_base(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Map::r, Lorene::Valeur::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_base(), and Lorene::Scalar::set_spectral_va().
|
inherited |
Returns the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart.
Definition at line 331 of file map.C.
References Lorene::Map::bvect_cart, and Lorene::Map::p_flat_met_cart.
|
inherited |
Returns the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher.
Definition at line 321 of file map.C.
References Lorene::Map::bvect_spher, and Lorene::Map::p_flat_met_spher.
|
inlineinherited |
Returns the Cartesian basis 
the Cartesian coordinates related to 
Definition at line 791 of file map.h.
References Lorene::Map::bvect_cart.
|
inlineinherited |
Returns the orthonormal vectorial basis 

Definition at line 783 of file map.h.
References Lorene::Map::bvect_spher.
Gives the Mg3d on which the mapping is defined.
Definition at line 765 of file map.h.
References Lorene::Map::mg.
|
inlineinherited |
Returns the x coordinate of the origin.
Definition at line 768 of file map.h.
References Lorene::Map::ori_x.
|
inlineinherited |
Returns the y coordinate of the origin.
Definition at line 770 of file map.h.
References Lorene::Map::ori_y.
|
inlineinherited |
Returns the z coordinate of the origin.
Definition at line 772 of file map.h.
References Lorene::Map::ori_z.
|
inlineinherited |
Returns the angle between the x –axis and X –axis.
Definition at line 775 of file map.h.
References Lorene::Map::rot_phi.
Sets a new radial scale.
| lambda | [input] factor by which the value of r is to be multiplied |
Implements Lorene::Map.
Definition at line 537 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and Lorene::Map_radial::reset_coord().
Sets a new radial scale at the bondary between the nucleus and the first shell.
| lambda | [input] factor by which the value of r is to be multiplied |
Definition at line 614 of file map_af.C.
References alpha, beta, and Lorene::Map_radial::reset_coord().
Increases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 799 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Increases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 696 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Computes the integral over all space of a Cmp.
The computed quantity is 
Tbl (size: mg->nzone ) to store the result (partial integral) in each domain and returns a pointer to it.
Implements Lorene::Map.
Definition at line 81 of file map_af_integ.C.
References alpha, Lorene::Tbl::annule_hard(), Lorene::Mtbl_cf::base, beta, Lorene::Cmp::check_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Mtbl_cf::get_etat(), Lorene::Tbl::get_etat(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::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, Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl_cf::t, Lorene::Tbl::t, T_COS, T_COS_P, T_COSSIN_C, T_COSSIN_CP, and Lorene::Cmp::va.
Performs the surface integration of ci on the sphere of radius rayon .
Definition at line 90 of file map_af_integ_surf.C.
References Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Cmp::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Base_val::get_base_p(), Lorene::Base_val::get_base_r(), Lorene::Base_val::get_base_t(), Lorene::Cmp::get_etat(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), P_COSSIN, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, Lorene::Cmp::va, and val_lx().
Performs the surface integration of ci on the sphere of radius rayon .
Definition at line 277 of file map_af_integ_surf.C.
References Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Base_val::get_base_p(), Lorene::Base_val::get_base_r(), Lorene::Base_val::get_base_t(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Scalar::get_spectral_va(), P_COSSIN, P_COSSIN_P, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPI_I, R_CHEBPI_P, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, and val_lx().
Definition at line 56 of file map_af_integ_surf_falloff.C.
Performs the surface integration of ci at infinity.
ci must have dzpuis =2.
Definition at line 196 of file map_af_integ_surf.C.
References Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Cmp::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Base_val::get_base_p(), Lorene::Base_val::get_base_r(), Lorene::Base_val::get_base_t(), Lorene::Cmp::get_etat(), Lorene::Map::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), P_COSSIN, P_COSSIN_P, R_CHEBU, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, T_COSSIN_CP, and Lorene::Cmp::va.
Performs the surface integration of ci at infinity.
ci must have dzpuis =2.
Definition at line 385 of file map_af_integ_surf.C.
References Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Base_val::get_base_p(), Lorene::Base_val::get_base_r(), Lorene::Base_val::get_base_t(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), P_COSSIN, P_COSSIN_P, R_CHEBU, T_COS, T_COS_I, T_COS_P, T_COSSIN_C, T_COSSIN_CI, and T_COSSIN_CP.
Computes the angular Laplacian of a scalar field.
| uu | [input] Scalar field u (represented as a Scalar) the Laplacian ![]() |
| lap | [output] Angular Laplacian of u (see documentation of Scalar |
Implements Lorene::Map.
Definition at line 549 of file map_af_lap.C.
References Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_etat_zero().
Computes the Laplacian of a scalar field (Cmp version).
Implements Lorene::Map.
Definition at line 367 of file map_af_lap.C.
References Lorene::Cmp::annule(), Lorene::Valeur::base, Lorene::Cmp::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Valeur::coef_i(), Lorene::Valeur::d2sdx2(), Lorene::Map_radial::dxdr, Lorene::Cmp::get_etat(), Lorene::Map::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Valeur::lapang(), Lorene::Map::mg, Lorene::Valeur::mult2_xm1_zec(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Valeur::sx2(), Lorene::Cmp::va, Lorene::Map_radial::xsr, and Lorene::Valeur::ylm().
Computes the Laplacian of a scalar field.
| uu | [input] Scalar field u (represented as a Scalar ) the Laplacian ![]() |
| zec_mult_r | [input] Determines the quantity computed in the compactified external domain (CED) : \ zec_mult_r = 0 : ![]() ![]() ![]() |
| lap | [output] Laplacian of u |
Implements Lorene::Map.
Definition at line 179 of file map_af_lap.C.
References Lorene::Scalar::annule(), Lorene::Valeur::base, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Valeur::coef_i(), Lorene::Valeur::d2sdx2(), Lorene::Map_radial::dxdr, Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Valeur::lapang(), Lorene::Map::mg, Lorene::Valeur::mult2_xm1_zec(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sx2(), Lorene::Map_radial::xsr, and Lorene::Valeur::ylm().
Returns the "angular" mapping for the outside of domain l_zone.
Valid only for the class Map_af.
Implements Lorene::Map.
Definition at line 656 of file map_af.C.
References Lorene::Mg3d::get_angu_1dom(), Lorene::Map::get_mg(), Lorene::Map::p_mp_angu, Lorene::Tbl::set(), Lorene::Tbl::set_etat_qcq(), and val_r_jk().
Multiplication by 
Scalar.
Implements Lorene::Map.
Definition at line 61 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Valeur::mult_ct(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
Multiplication by r of a Cmp.
In the CED, there is only a decrement of dzpuis
Implements Lorene::Map.
Definition at line 219 of file map_radial_r_manip.C.
References Lorene::Cmp::annule(), Lorene::Valeur::base, Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_x(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::va, and Lorene::Map_radial::xsr.
Multiplication by r of a Scalar, the dzpuis
of uu is not changed.
Implements Lorene::Map.
Definition at line 144 of file map_radial_r_manip.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Multiplication by r (in the compactified external domain only) of a Scalar.
Implements Lorene::Map.
Definition at line 296 of file map_radial_r_manip.C.
References Lorene::Valeur::annule(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Multiplication by 
Scalar.
Implements Lorene::Map.
Definition at line 355 of file map_radial_r_manip.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_st(), Lorene::Valeur::mult_x(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::sxm1_zec(), and Lorene::Map_radial::xsr.
Multiplication by 
Scalar.
Implements Lorene::Map.
Definition at line 109 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat(), Lorene::Valeur::get_mg(), Lorene::Map::mg, Lorene::Valeur::mult_st(), Lorene::Scalar::set_etat_qcq(), and Lorene::Scalar::set_spectral_va().
Assignment to another affine mapping.
Implements Lorene::Map_radial.
Definition at line 380 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Map::mg, Lorene::Map::ori_x, Lorene::Map::ori_y, Lorene::Map::ori_z, Lorene::Map_radial::reset_coord(), Lorene::Map::rot_phi, Lorene::Map::set_ori(), and Lorene::Map::set_rot_phi().
Comparison operator (egality)
Implements Lorene::Map_radial.
Definition at line 439 of file map_af.C.
References alpha, beta, Lorene::Map::bvect_cart, Lorene::Map::bvect_spher, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Map::get_ori_x(), Lorene::Map::get_ori_y(), Lorene::Map::get_ori_z(), Lorene::Map::mg, Lorene::Map::ori_x, Lorene::Map::ori_y, and Lorene::Map::ori_z.
Operator >>
Implements Lorene::Map.
Definition at line 503 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, Lorene::Map::r, and val_r().
Computes the solution of a scalar Poisson equation.
| source | [input] source ![]() ![]() |
| par | [] not used by this Map_af version. |
| uu | [output] solution u with the boundary condition u =0 at spatial infinity. |
Implements Lorene::Map.
Definition at line 100 of file map_af_poisson.C.
References Lorene::Valeur::c_cf, Lorene::Cmp::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Cmp::dz_nonzero(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Cmp::get_mp(), Lorene::Map::mg, Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Cmp::va, and Lorene::Valeur::ylm().
|
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 \ .
\]](form_582.png)
| source_mat | [input] Compactly supported part of the source ![]() |
| source_quad | [input] Non-compactly supported part of the source ![]() |
| par | [output] Parameter which contains the constant ![]() 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 Lorene::Map.
Definition at line 81 of file map_af_poisson2d.C.
References Lorene::Cmp::allocate_all(), alpha, Lorene::Valeur::base, beta, Lorene::Valeur::c, Lorene::Cmp::check_dzpuis(), Lorene::Valeur::coef_i(), Lorene::Param::get_double_mod(), Lorene::Cmp::get_etat(), Lorene::Mtbl::get_etat(), Lorene::Tbl::get_etat(), Lorene::Valeur::get_etat(), Lorene::Mg3d::get_grille3d(), Lorene::Map::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Mg3d::get_type_t(), Lorene::Map::mg, Lorene::Cmp::set(), Lorene::Cmp::set_dzpuis(), Lorene::Mtbl::t, Lorene::Tbl::t, T_COS, T_COS_P, T_SIN, T_SIN_I, Lorene::Cmp::va, and Lorene::Grille3d::x.
|
virtual |
Computes the solution of the generalized angular Poisson equation.
The generalized angular Poisson equation is 

| source | [input] source ![]() ![]() |
| 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 ![]() |
Implements Lorene::Map.
Definition at line 58 of file map_af_poisson_angu.C.
References Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Tensor::get_mp(), Lorene::Param::get_n_int(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, poisson_angu(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
|
virtualinherited |
Resolution of the elliptic equation 
| source | [input] source ![]() |
| 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 ![]() ![]() par.get_double(0) \ par.get_double(1) : [input] relaxation parameter ![]() par.get_int_mod(0) : [output] number of iterations actually used to get the solution. |
| psi | [input/output]: input : previously computed value of ![]() psi must be set to zero); output: solution ![]() ![]() |
Implements Lorene::Map.
Definition at line 155 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(), Lorene::Map_radial::dxdr, Lorene::Param::get_double(), Lorene::Cmp::get_etat(), Lorene::Tenseur::get_etat(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Cmp::get_mp(), Lorene::Tenseur::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().
|
virtualinherited |
Resolution of the elliptic equation 
| nzet | [input] number of domains covering the stellar interior |
| source | [input] source ![]() |
| 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 ![]() ![]() |
Implements Lorene::Map.
Definition at line 453 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::Cmp::dsdr(), dsdr(), Lorene::Param::get_double(), Lorene::Cmp::get_etat(), Lorene::Tenseur::get_etat(), Lorene::Mg3d::get_grille3d(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Cmp::get_mp(), Lorene::Tenseur::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Tenseur::get_triad(), laplacien(), Lorene::Cmp::laplacien(), Lorene::max(), Lorene::Map::mg, Lorene::min(), Lorene::Map_radial::poisson_compact(), Lorene::Tbl::set(), Lorene::Mtbl::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().
|
virtual |
Implements Lorene::Map.
Definition at line 57 of file map_af_poisson_falloff.C.
|
virtual |
Solver of the Poisson equation with boundary condition for the affine mapping case.
Implements Lorene::Map.
Definition at line 135 of file cmp_pde_frontiere.C.
References Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Cmp::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Cmp::dz_nonzero(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Cmp::get_mp(), Lorene::Map::mg, Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Mtbl_cf::set_etat_zero(), Lorene::Cmp::va, and Lorene::Valeur::ylm().
|
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 Lorene::Map.
Definition at line 208 of file cmp_pde_frontiere.C.
References Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Mg3d::get_angu(), Lorene::Cmp::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Cmp::get_mp(), Lorene::Map::mg, Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Mtbl_cf::set_etat_zero(), Lorene::Cmp::va, and Lorene::Valeur::ylm().
|
virtual |
Computes the solution of a Poisson equation in the shell, imposing a boundary condition at the surface of the star.
| 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 Lorene::Map.
Definition at line 474 of file cmp_pde_frontiere.C.
References Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), Lorene::Mg3d::get_angu(), Lorene::Cmp::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Cmp::get_mp(), Lorene::Map::mg, Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Mtbl_cf::set_etat_zero(), Lorene::Cmp::va, and Lorene::Valeur::ylm().
|
virtual |
Computes the solution of a scalar Poisson equation.
The regularized source 
| source | [input] source ![]() ![]() |
| k_div | [input] regularization degree of the procedure |
| nzet | [input] number of domains covering the star |
| unsgam1 | [input] parameter ![]() ![]() |
| 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 Lorene::Map.
Definition at line 109 of file map_af_poisson_regu.C.
References alpha, Lorene::Cmp::annule(), Lorene::Base_val::b, Lorene::Mtbl_cf::base, Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Valeur::coef(), dsdt(), Lorene::Map_radial::dxdr, Lorene::Cmp::get_etat(), Lorene::Valeur::get_etat(), Lorene::Mg3d::get_grille3d(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Matrice::inverse(), Lorene::Map::mg, MSQ_P, MSQ_R, Lorene::pow(), R_CHEBPIM_I, Lorene::Tenseur::set(), Lorene::Tbl::set(), Lorene::Matrice::set(), Lorene::Mtbl_cf::set(), Lorene::Matrice::set_band(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_qcq(), Lorene::Matrice::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Tenseur::set_etat_qcq(), Lorene::Matrice::set_lu(), Lorene::Tenseur::set_triad(), Lorene::Mg3d::std_base_scal(), stdsdp(), T_LEG_P, Lorene::Cmp::va, Lorene::Grille3d::x, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
Computes the solution of a scalar Poisson equation using a Tau method.
| source | [input] source ![]() ![]() |
| par | [] not used by this Map_af version. |
| uu | [output] solution u with the boundary condition u =0 at spatial infinity. |
Implements Lorene::Map.
Definition at line 165 of file map_af_poisson.C.
References Lorene::Valeur::c_cf, Lorene::Cmp::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Cmp::dz_nonzero(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Cmp::get_mp(), Lorene::Map::mg, Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_qcq(), Lorene::Cmp::set_etat_zero(), Lorene::Cmp::va, and Lorene::Valeur::ylm().
|
virtual |
Implements Lorene::Map.
Definition at line 57 of file map_af_poisson_ylm.C.
Computes the radial primitive which vanishes for 
i.e. the function 
| 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 Lorene::Map.
Definition at line 74 of file map_af_primr.C.
References alpha, Lorene::Base_val::b, Lorene::Mtbl_cf::base, Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Base_val::get_b(), Lorene::Valeur::get_base(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), MAX_BASE, Lorene::Map::mg, MSQ_R, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, Lorene::Tbl::set(), Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Valeur::set_etat_cf_qcq(), Lorene::Mtbl_cf::set_etat_qcq(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Tbl::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Mtbl_cf::t, TRA_R, and Lorene::Mtbl_cf::val_out_bound_jk().
|
virtualinherited |
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
| 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 ![]() uu is set to zero in the other domains. |
| uu | [input/output] input : Cmp previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Cmp at the grid points defined by *this. |
Implements Lorene::Map.
Definition at line 58 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::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::t, Lorene::Cmp::va, Lorene::Map_radial::val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk().
|
virtualinherited |
Recomputes the values of a Scalar at the collocation points after a change in the mapping.
| 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 ![]() uu is set to zero in the other domains. |
| uu | [input/output] input : Scalar previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Scalar at the grid points defined by *this. |
Implements Lorene::Map.
Definition at line 173 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::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Mtbl::t, Lorene::Map_radial::val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk().
|
virtualinherited |
Recomputes the values of a Cmp at the collocation points after a change in the mapping.
Case where the Cmp is symmetric with respect to the plane y=0.
| 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 ![]() uu is set to zero in the other domains. |
| uu | [input/output] input : Cmp previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Cmp at the grid points defined by *this. |
Implements Lorene::Map.
Definition at line 59 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::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Mtbl::t, Lorene::Cmp::va, Lorene::Map_radial::val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk_symy().
|
virtualinherited |
Recomputes the values of a Scalar at the collocation points after a change in the mapping.
Case where the Scalar is symmetric with respect to the plane y=0.
| 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 ![]() uu is set to zero in the other domains. |
| uu | [input/output] input : Scalar previously computed on the mapping *mp_prev ; ouput : values of (logically) the same Scalar at the grid points defined by *this. |
Implements Lorene::Map.
Definition at line 193 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::Mtbl::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Mtbl::t, Lorene::Map_radial::val_lx_jk(), and Lorene::Mtbl_cf::val_point_jk_symy().
|
protectedvirtualinherited |
Resets all the member Coords.
Reimplemented from Lorene::Map.
Reimplemented in Lorene::Map_et.
Definition at line 126 of file map_radial.C.
References Lorene::Map_radial::d2rdtdx, Lorene::Map_radial::d2rdx2, Lorene::Coord::del_t(), Lorene::Map_radial::drdt, Lorene::Map_radial::dxdr, Lorene::Map_radial::lapr_tp, Lorene::Map::reset_coord(), Lorene::Map_radial::sr2d2rdt2, Lorene::Map_radial::sr2drdt, Lorene::Map_radial::sr2stdrdp, Lorene::Map_radial::srdrdt, Lorene::Map_radial::srstdrdp, Lorene::Map_radial::sstd2rdpdx, Lorene::Map_radial::stdrdp, and Lorene::Map_radial::xsr.
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.
| l | [input] index of the domain |
| lambda | [input] factor by which the value of ![]() |
Implements Lorene::Map.
Definition at line 560 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and Lorene::Map_radial::reset_coord().
Save in a file.
Reimplemented from Lorene::Map_radial.
Definition at line 489 of file map_af.C.
References alpha, beta, Lorene::fwrite_be(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and Lorene::Map_radial::sauve().
Modifies the value of 
Definition at line 630 of file map_af.C.
References alpha, and Lorene::Map_radial::reset_coord().
Modifies the value of 
Definition at line 641 of file map_af.C.
References beta, and Lorene::Map_radial::reset_coord().
|
private |
Assignment of the building functions to the member Coords.
Definition at line 403 of file map_af.C.
References Lorene::Map::cosp, Lorene::Map::cost, Lorene::Map_radial::d2rdtdx, Lorene::Map_radial::d2rdx2, Lorene::Map_radial::drdt, Lorene::Map_radial::dxdr, Lorene::Map_radial::lapr_tp, Lorene::Map::phi, Lorene::Map::r, Lorene::Coord::set(), Lorene::Map::sinp, Lorene::Map::sint, Lorene::Map_radial::sr2d2rdt2, Lorene::Map_radial::sr2drdt, Lorene::Map_radial::sr2stdrdp, Lorene::Map_radial::srdrdt, Lorene::Map_radial::srstdrdp, Lorene::Map_radial::sstd2rdpdx, Lorene::Map_radial::stdrdp, Lorene::Map::tet, Lorene::Map::x, Lorene::Map::xa, Lorene::Map_radial::xsr, Lorene::Map::y, Lorene::Map::ya, Lorene::Map::z, and Lorene::Map::za.
Sets a new origin.
Definition at line 253 of file map.C.
References Lorene::Map::bvect_spher, Lorene::Map::ori_x, Lorene::Map::ori_y, Lorene::Map::ori_z, Lorene::Map::reset_coord(), and Lorene::Base_vect_spher::set_ori().
Sets a new rotation angle.
Definition at line 263 of file map.C.
References Lorene::Map::bvect_cart, Lorene::Map::bvect_spher, Lorene::Map::reset_coord(), Lorene::Map::rot_phi, Lorene::Base_vect_cart::set_rot_phi(), and Lorene::Base_vect_spher::set_rot_phi().
| void Lorene::Map_af::sol_elliptic | ( | Param_elliptic & | params, |
| const Scalar & | so, | ||
| Scalar & | uu | ||
| ) | const |
General elliptic solver.
The field is zero at infinity.
| params | [input] : the operators and variables to be uses. |
| so | [input] : the source. |
| uu | [output] : the solution. |
Definition at line 100 of file map_af_elliptic.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| void Lorene::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.
| params | [input] : the operators and variables to be uses. |
| so | [input] : the source. |
| uu | [output] : the solution. |
Definition at line 58 of file map_af_elliptic_2d.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, and Lorene::Param_elliptic::var_G.
| void Lorene::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.
| 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 159 of file map_af_elliptic.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| void Lorene::Map_af::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.
Definition at line 222 of file map_af_elliptic.C.
References Lorene::Mtbl_cf::annule_hard(), Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Mg3d::get_angu_1dom(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Mtbl_cf::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| void Lorene::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.
| 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 483 of file map_af_elliptic.C.
References alpha, Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| void Lorene::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.
| 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 312 of file map_af_elliptic.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| void Lorene::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.
| 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 368 of file map_af_elliptic.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| void Lorene::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.
| params | [input] : the operators and variables to be uses. |
| so | [input] : the source. |
| uu | [output] : the solution. |
Definition at line 59 of file map_af_elliptic_pseudo_1d.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Valeur::val_propre_1d(), Lorene::Valeur::val_propre_1d_i(), Lorene::Param_elliptic::var_F, and Lorene::Param_elliptic::var_G.
| void Lorene::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 426 of file map_af_elliptic.C.
References Lorene::Valeur::c_cf, Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::get_etat(), Lorene::Valeur::get_etat(), Lorene::Map::get_mg(), Lorene::Valeur::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Scalar::set_spectral_va(), Lorene::Param_elliptic::var_F, Lorene::Param_elliptic::var_G, Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
Computes 
Cmp.
Note that in the compactified external domain (CED), it computes 
| ci | [input] field to consider |
| resu | [output] derivative of ci |
Implements Lorene::Map.
Definition at line 475 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Valeur::dsdt(), Lorene::Valeur::get_base(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Map::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set_base(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_zero(), Lorene::Valeur::sx(), Lorene::Cmp::va, and Lorene::Map_radial::xsr.
Computes 
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.
| uu | [input] field to consider |
| resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 569 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Valeur::coef(), Lorene::Valeur::dsdt(), Lorene::Valeur::get_base(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::set_base(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Valeur::sx(), and Lorene::Map_radial::xsr.
Computes 
Cmp.
Note that in the compactified external domain (CED), it computes 
| ci | [input] field to consider |
| resu | [output] derivative of ci |
Implements Lorene::Map.
Definition at line 633 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Valeur::dsdp(), Lorene::Valeur::get_base(), Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Map::get_mg(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::mult_xm1_zec(), Lorene::Valeur::set_base(), Lorene::Cmp::set_dzpuis(), Lorene::Cmp::set_etat_zero(), Lorene::Valeur::ssint(), Lorene::Valeur::sx(), Lorene::Cmp::va, and Lorene::Map_radial::xsr.
Computes 
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.
| uu | [input] field to consider |
| resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 729 of file map_af_deriv.C.
References Lorene::Valeur::annule(), Lorene::Valeur::coef(), Lorene::Valeur::dsdp(), Lorene::Valeur::get_base(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::Valeur::set_base(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), Lorene::Valeur::ssint(), Lorene::Valeur::sx(), and Lorene::Map_radial::xsr.
Computes 
Scalar.
| uu | [input] scalar field |
| resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 823 of file map_af_deriv.C.
References Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Tensor::get_mp(), Lorene::Scalar::get_spectral_va(), Lorene::Map::mg, Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_zero(), and Lorene::Valeur::stdsdp().
|
virtual |
Computes the domain index l and the value of 

| rr | [input] value of r |
| theta | [input] value of ![]() |
| pphi | [input] value of ![]() |
| par | [] unused by the Map_af version |
| l | [output] value of the domain index |
| xi | [output] value of ![]() |
Implements Lorene::Map.
Definition at line 192 of file map_af_radius.C.
References val_lx().
|
virtual |
Computes the domain index l and the value of 

| rr | [input] value of r |
| theta | [input] value of ![]() |
| pphi | [input] value of ![]() |
| l | [output] value of the domain index |
| xi | [output] value of ![]() |
Implements Lorene::Map.
Definition at line 128 of file map_af_radius.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
virtual |
Computes the domain index l and the value of 

| rr | [input] value of r |
| j | [input] index of the collocation point in ![]() |
| k | [input] index of the collocation point in ![]() |
| par | [] unused by the Map_af version |
| l | [output] value of the domain index |
| xi | [output] value of ![]() |
Implements Lorene::Map_radial.
Definition at line 215 of file map_af_radius.C.
References val_lx().
Returns the value of the radial coordinate r for a given 
| l | [input] index of the domain |
| xi | [input] value of ![]() |
| theta | [input] value of ![]() |
| pphi | [input] value of ![]() |

Implements Lorene::Map.
Definition at line 96 of file map_af_radius.C.
References alpha, beta, Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
Returns the value of the radial coordinate r for a given 

| l | [input] index of the domain |
| xi | [input] value of ![]() |
| j | [input] index of the collocation point in ![]() |
| k | [input] index of the collocation point in ![]() |

Implements Lorene::Map_radial.
Definition at line 205 of file map_af_radius.C.
References val_r().
Definition at line 492 of file map_af_fait.C.
Definition at line 424 of file map_af_fait.C.
Definition at line 825 of file map_af_fait.C.
Definition at line 789 of file map_af_fait.C.
Definition at line 681 of file map_af_fait.C.
Definition at line 617 of file map_af_fait.C.
Definition at line 807 of file map_af_fait.C.
Definition at line 218 of file map_af_fait.C.
Definition at line 116 of file map_af_fait.C.
Definition at line 458 of file map_af_fait.C.
Definition at line 390 of file map_af_fait.C.
Definition at line 861 of file map_af_fait.C.
Definition at line 753 of file map_af_fait.C.
Definition at line 771 of file map_af_fait.C.
Definition at line 717 of file map_af_fait.C.
Definition at line 735 of file map_af_fait.C.
Definition at line 843 of file map_af_fait.C.
Definition at line 699 of file map_af_fait.C.
Definition at line 180 of file map_af_fait.C.
Definition at line 256 of file map_af_fait.C.
Definition at line 310 of file map_af_fait.C.
Definition at line 532 of file map_af_fait.C.
Definition at line 274 of file map_af_fait.C.
Definition at line 340 of file map_af_fait.C.
Definition at line 292 of file map_af_fait.C.
Definition at line 370 of file map_af_fait.C.
|
private |
|
private |
|
protectedinherited |
|
protectedinherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
inherited |
|
inherited |
|
protectedinherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |