LORENE
|
Affine radial mapping. More...
#include <map.h>
Public Member Functions | |
Map_af (const Mg3d &mgrille, const double *r_limits) | |
Standard Constructor. More... | |
Map_af (const Mg3d &mgrille, const Tbl &r_limits) | |
Standard Constructor with Tbl. More... | |
Map_af (const Map_af &) | |
Copy constructor. More... | |
Map_af (const Mg3d &, const string &) | |
Constructor from a formatted file. More... | |
Map_af (const Mg3d &, FILE *) | |
Constructor from a file (see sauve(FILE*) ) More... | |
Map_af (const Map &) | |
Constructor from a general mapping. More... | |
virtual | ~Map_af () |
Destructor. More... | |
virtual void | operator= (const Map_af &) |
Assignment to another affine mapping. More... | |
const double * | get_alpha () const |
Returns the pointer on the array alpha . More... | |
const double * | get_beta () const |
Returns the pointer on the array beta . More... | |
virtual const Map_af & | mp_angu (int) const |
Returns the "angular" mapping for the outside of domain l_zone . More... | |
virtual double | val_r (int l, double xi, double theta, double pphi) const |
Returns the value of the radial coordinate r for a given in a given domain. More... | |
virtual void | val_lx (double rr, double theta, double pphi, int &l, double &xi) const |
Computes the domain index l and the value of corresponding to a point given by its physical coordinates . More... | |
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 corresponding to a point given by its physical coordinates . More... | |
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 and a given collocation point in in a given domain. More... | |
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 corresponding to a point of arbitrary r but collocation values of . More... | |
virtual bool | operator== (const Map &) const |
Comparison operator (egality) More... | |
virtual void | sauve (FILE *) const |
Save in a file. More... | |
virtual void | homothetie (double lambda) |
Sets a new radial scale. More... | |
virtual void | resize (int l, double lambda) |
Rescales the outer boundary of one domain. More... | |
void | homothetie_interne (double lambda) |
Sets a new radial scale at the bondary between the nucleus and the first shell. More... | |
virtual void | adapt (const Cmp &ent, const Param &par, int nbr=0) |
Adaptation of the mapping to a given scalar field. More... | |
void | set_alpha (double alpha0, int l) |
Modifies the value of in domain no. l. More... | |
void | set_beta (double beta0, int l) |
Modifies the value of in domain no. l. More... | |
virtual void | dsdxi (const Cmp &ci, Cmp &resu) const |
Computes of a Cmp . More... | |
virtual void | dsdr (const Cmp &ci, Cmp &resu) const |
Computes of a Cmp . More... | |
virtual void | srdsdt (const Cmp &ci, Cmp &resu) const |
Computes of a Cmp . More... | |
virtual void | srstdsdp (const Cmp &ci, Cmp &resu) const |
Computes of a Cmp . More... | |
virtual void | dsdr (const Scalar &uu, Scalar &resu) const |
Computes of a Scalar . More... | |
virtual void | dsdxi (const Scalar &uu, Scalar &resu) const |
Computes of a Scalar . More... | |
virtual void | dsdradial (const Scalar &, Scalar &) const |
Computes of a Scalar . More... | |
virtual void | srdsdt (const Scalar &uu, Scalar &resu) const |
Computes of a Scalar . More... | |
virtual void | srstdsdp (const Scalar &uu, Scalar &resu) const |
Computes of a Scalar . More... | |
virtual void | dsdt (const Scalar &uu, Scalar &resu) const |
Computes of a Scalar . More... | |
virtual void | stdsdp (const Scalar &uu, Scalar &resu) const |
Computes of a Scalar . More... | |
virtual void | laplacien (const Scalar &uu, int zec_mult_r, Scalar &lap) const |
Computes the Laplacian of a scalar field. More... | |
virtual void | laplacien (const Cmp &uu, int zec_mult_r, Cmp &lap) const |
Computes the Laplacian of a scalar field (Cmp version). More... | |
virtual void | lapang (const Scalar &uu, Scalar &lap) const |
Computes the angular Laplacian of a scalar field. More... | |
virtual void | primr (const Scalar &uu, Scalar &resu, bool null_infty) const |
Computes the radial primitive which vanishes for . More... | |
virtual Tbl * | integrale (const Cmp &) const |
Computes the integral over all space of a Cmp . More... | |
virtual void | poisson (const Cmp &source, Param &par, Cmp &uu) const |
Computes the solution of a scalar Poisson equation. More... | |
virtual void | poisson_tau (const Cmp &source, Param &par, Cmp &uu) const |
Computes the solution of a scalar Poisson equation using a Tau method. More... | |
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. More... | |
virtual void | poisson_angu (const Scalar &source, Param &par, Scalar &uu, double lambda=0) const |
Computes the solution of the generalized angular Poisson equation. More... | |
virtual void | poisson_angu (const Cmp &source, Param &par, Cmp &uu, double lambda=0) const |
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 . More... | |
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. More... | |
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). More... | |
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. More... | |
double | integrale_surface (const Cmp &ci, double rayon) const |
Performs the surface integration of ci on the sphere of radius rayon . More... | |
double | integrale_surface (const Scalar &ci, double rayon) const |
Performs the surface integration of ci on the sphere of radius rayon . More... | |
double | integrale_surface_falloff (const Cmp &ci) const |
double | integrale_surface_infini (const Cmp &ci) const |
Performs the surface integration of ci at infinity. More... | |
double | integrale_surface_infini (const Scalar &ci) const |
Performs the surface integration of ci at infinity. More... | |
void | sol_elliptic (Param_elliptic ¶ms, const Scalar &so, Scalar &uu) const |
General elliptic solver. More... | |
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. More... | |
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. More... | |
void | sol_elliptic_no_zec (Param_elliptic ¶ms, const Scalar &so, Scalar &uu, double val) const |
General elliptic solver. More... | |
void | sol_elliptic_only_zec (Param_elliptic ¶ms, const Scalar &so, Scalar &uu, double val) const |
General elliptic solver. More... | |
void | sol_elliptic_sin_zec (Param_elliptic ¶ms, const Scalar &so, Scalar &uu, double *coefs, double *) const |
General elliptic solver. More... | |
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. More... | |
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. More... | |
void | sol_elliptic_2d (Param_elliptic &, const Scalar &, Scalar &) const |
General elliptic solver in a 2D case. More... | |
void | sol_elliptic_pseudo_1d (Param_elliptic &, const Scalar &, Scalar &) const |
General elliptic solver in a pseudo 1d case. More... | |
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. More... | |
Scalar | interpolate_from_map_star (const Scalar &f_s) const |
Interpolates from a Scalar defined on a Map_star and returns the new Scalar defined on *this . More... | |
virtual void | reevaluate (const Map *mp_prev, int nzet, Cmp &uu) const |
Recomputes the values of a Cmp at the collocation points after a change in the mapping. More... | |
virtual void | reevaluate (const Map *mp_prev, int nzet, Scalar &uu) const |
Recomputes the values of a Scalar at the collocation points after a change in the mapping. More... | |
virtual void | reevaluate_symy (const Map *mp_prev, int nzet, Cmp &uu) const |
Recomputes the values of a Cmp at the collocation points after a change in the mapping. More... | |
virtual void | reevaluate_symy (const Map *mp_prev, int nzet, Scalar &uu) const |
Recomputes the values of a Scalar at the collocation points after a change in the mapping. More... | |
virtual void | mult_r (Scalar &uu) const |
Multiplication by r of a Scalar , the dzpuis of uu is not changed. More... | |
virtual void | mult_r (Cmp &ci) const |
Multiplication by r of a Cmp . More... | |
virtual void | mult_r_zec (Scalar &) const |
Multiplication by r (in the compactified external domain only) of a Scalar . More... | |
virtual void | mult_rsint (Scalar &) const |
Multiplication by of a Scalar . More... | |
virtual void | div_rsint (Scalar &) const |
Division by of a Scalar . More... | |
virtual void | div_r (Scalar &) const |
Division by r of a Scalar . More... | |
virtual void | div_r_zec (Scalar &) const |
Division by r (in the compactified external domain only) of a Scalar . More... | |
virtual void | mult_cost (Scalar &) const |
Multiplication by of a Scalar . More... | |
virtual void | div_cost (Scalar &) const |
Division by of a Scalar . More... | |
virtual void | mult_sint (Scalar &) const |
Multiplication by of a Scalar . More... | |
virtual void | div_sint (Scalar &) const |
Division by of a Scalar . More... | |
virtual void | div_tant (Scalar &) const |
Division by of a Scalar . More... | |
virtual void | comp_x_from_spherical (const Scalar &v_r, const Scalar &v_theta, const Scalar &v_phi, Scalar &v_x) const |
Computes the Cartesian x component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher . More... | |
virtual void | comp_x_from_spherical (const Cmp &v_r, const Cmp &v_theta, const Cmp &v_phi, Cmp &v_x) const |
Cmp version More... | |
virtual void | comp_y_from_spherical (const Scalar &v_r, const Scalar &v_theta, const Scalar &v_phi, Scalar &v_y) const |
Computes the Cartesian y component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher . More... | |
virtual void | comp_y_from_spherical (const Cmp &v_r, const Cmp &v_theta, const Cmp &v_phi, Cmp &v_y) const |
Cmp version More... | |
virtual void | comp_z_from_spherical (const Scalar &v_r, const Scalar &v_theta, Scalar &v_z) const |
Computes the Cartesian z component (with respect to bvect_cart ) of a vector given by its spherical components with respect to bvect_spher . More... | |
virtual void | comp_z_from_spherical (const Cmp &v_r, const Cmp &v_theta, Cmp &v_z) const |
Cmp version More... | |
virtual void | comp_r_from_cartesian (const Scalar &v_x, const Scalar &v_y, const Scalar &v_z, Scalar &v_r) const |
Computes the Spherical r component (with respect to bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart . More... | |
virtual void | comp_r_from_cartesian (const Cmp &v_x, const Cmp &v_y, const Cmp &v_z, Cmp &v_r) const |
Cmp version More... | |
virtual void | comp_t_from_cartesian (const Scalar &v_x, const Scalar &v_y, const Scalar &v_z, Scalar &v_t) const |
Computes the Spherical component (with respect to bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart . More... | |
virtual void | comp_t_from_cartesian (const Cmp &v_x, const Cmp &v_y, const Cmp &v_z, Cmp &v_t) const |
Cmp version More... | |
virtual void | comp_p_from_cartesian (const Scalar &v_x, const Scalar &v_y, Scalar &v_p) const |
Computes the Spherical component (with respect to bvect_spher ) of a vector given by its cartesian components with respect to bvect_cart . More... | |
virtual void | comp_p_from_cartesian (const Cmp &v_x, const Cmp &v_y, Cmp &v_p) const |
Cmp version More... | |
virtual void | dec_dzpuis (Scalar &) const |
Decreases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED). More... | |
virtual void | dec2_dzpuis (Scalar &) const |
Decreases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED). More... | |
virtual void | inc_dzpuis (Scalar &) const |
Increases by 1 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED). More... | |
virtual void | inc2_dzpuis (Scalar &) const |
Increases by 2 the value of dzpuis of a Scalar and changes accordingly its values in the compactified external domain (CED). More... | |
virtual void | poisson_compact (const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const |
Resolution of the elliptic equation in the case where the stellar interior is covered by a single domain. More... | |
virtual void | poisson_compact (int nzet, const Cmp &source, const Cmp &aa, const Tenseur &bb, const Param &par, Cmp &psi) const |
Resolution of the elliptic equation in the case of a multidomain stellar interior. More... | |
const Mg3d * | get_mg () const |
Gives the Mg3d on which the mapping is defined. More... | |
double | get_ori_x () const |
Returns the x coordinate of the origin. More... | |
double | get_ori_y () const |
Returns the y coordinate of the origin. More... | |
double | get_ori_z () const |
Returns the z coordinate of the origin. More... | |
double | get_rot_phi () const |
Returns the angle between the x –axis and X –axis. More... | |
const Base_vect_spher & | get_bvect_spher () const |
Returns the orthonormal vectorial basis associated with the coordinates of the mapping. More... | |
const Base_vect_cart & | get_bvect_cart () const |
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e. More... | |
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 . More... | |
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 . More... | |
const Cmp & | cmp_zero () const |
Returns the null Cmp defined on *this . More... | |
void | convert_absolute (double xx, double yy, double zz, double &rr, double &theta, double &pphi) const |
Determines the coordinates corresponding to given absolute Cartesian coordinates (X,Y,Z). More... | |
void | set_ori (double xa0, double ya0, double za0) |
Sets a new origin. More... | |
void | set_rot_phi (double phi0) |
Sets a new rotation angle. More... | |
Public Attributes | |
Coord | xsr |
in the nucleus; \ 1/R in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | dxdr |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | drdt |
in the nucleus and in the non-compactified shells; \ in the compactified external domain (CED). More... | |
Coord | stdrdp |
in the nucleus and in the non-compactified shells; \ in the compactified external domain (CED). More... | |
Coord | srdrdt |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | srstdrdp |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | sr2drdt |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | sr2stdrdp |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | d2rdx2 |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | lapr_tp |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | d2rdtdx |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | sstd2rdpdx |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | sr2d2rdt2 |
in the nucleus and in the non-compactified shells; \ in the compactified outer domain. More... | |
Coord | r |
r coordinate centered on the grid More... | |
Coord | tet |
coordinate centered on the grid More... | |
Coord | phi |
coordinate centered on the grid More... | |
Coord | sint |
More... | |
Coord | cost |
More... | |
Coord | sinp |
More... | |
Coord | cosp |
More... | |
Coord | x |
x coordinate centered on the grid More... | |
Coord | y |
y coordinate centered on the grid More... | |
Coord | z |
z coordinate centered on the grid More... | |
Coord | xa |
Absolute x coordinate. More... | |
Coord | ya |
Absolute y coordinate. More... | |
Coord | za |
Absolute z coordinate. More... | |
Protected Member Functions | |
virtual void | reset_coord () |
Resets all the member Coords . More... | |
Protected Attributes | |
const Mg3d * | mg |
Pointer on the multi-grid Mgd3 on which this is defined. More... | |
double | ori_x |
Absolute coordinate x of the origin. More... | |
double | ori_y |
Absolute coordinate y of the origin. More... | |
double | ori_z |
Absolute coordinate z of the origin. More... | |
double | rot_phi |
Angle between the x –axis and X –axis. More... | |
Base_vect_spher | bvect_spher |
Orthonormal vectorial basis associated with the coordinates of the mapping. More... | |
Base_vect_cart | bvect_cart |
Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e. More... | |
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 . More... | |
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 . More... | |
Cmp * | p_cmp_zero |
The null Cmp. More... | |
Map_af * | p_mp_angu |
Pointer on the "angular" mapping. More... | |
Private Member Functions | |
void | set_coord () |
Assignment of the building functions to the member Coords . More... | |
virtual ostream & | operator>> (ostream &) const |
Operator >> More... | |
Private Attributes | |
double * | alpha |
Array (size: mg->nzone ) of the values of in each domain. More... | |
double * | beta |
Array (size: mg->nzone ) of the values of in each domain. More... | |
Friends | |
Mtbl * | map_af_fait_r (const Map *) |
Mtbl * | map_af_fait_tet (const Map *) |
Mtbl * | map_af_fait_phi (const Map *) |
Mtbl * | map_af_fait_sint (const Map *) |
Mtbl * | map_af_fait_cost (const Map *) |
Mtbl * | map_af_fait_sinp (const Map *) |
Mtbl * | map_af_fait_cosp (const Map *) |
Mtbl * | map_af_fait_x (const Map *) |
Mtbl * | map_af_fait_y (const Map *) |
Mtbl * | map_af_fait_z (const Map *) |
Mtbl * | map_af_fait_xa (const Map *) |
Mtbl * | map_af_fait_ya (const Map *) |
Mtbl * | map_af_fait_za (const Map *) |
Mtbl * | map_af_fait_xsr (const Map *) |
Mtbl * | map_af_fait_dxdr (const Map *) |
Mtbl * | map_af_fait_drdt (const Map *) |
Mtbl * | map_af_fait_stdrdp (const Map *) |
Mtbl * | map_af_fait_srdrdt (const Map *) |
Mtbl * | map_af_fait_srstdrdp (const Map *) |
Mtbl * | map_af_fait_sr2drdt (const Map *) |
Mtbl * | map_af_fait_sr2stdrdp (const Map *) |
Mtbl * | map_af_fait_d2rdx2 (const Map *) |
Mtbl * | map_af_fait_lapr_tp (const Map *) |
Mtbl * | map_af_fait_d2rdtdx (const Map *) |
Mtbl * | map_af_fait_sstd2rdpdx (const Map *) |
Mtbl * | map_af_fait_sr2d2rdt2 (const Map *) |
Affine radial mapping.
()
The affine radial mapping is the simplest one between the grid coordinates and the physical coordinates . It is defined by , and
Lorene::Map_af::Map_af | ( | const Mg3d & | mgrille, |
const double * | r_limits | ||
) |
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 216 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 261 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, and set_coord().
Lorene::Map_af::Map_af | ( | const Map_af & | mp | ) |
Copy constructor.
Definition at line 307 of file map_af.C.
References alpha, beta, Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and set_coord().
Lorene::Map_af::Map_af | ( | const Mg3d & | mgrille, |
const string & | filename | ||
) |
Constructor from a formatted file.
Definition at line 326 of file map_af.C.
References Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Map::mg, Lorene::search_file(), Lorene::Tbl::set(), and Lorene::Tbl::set_etat_qcq().
Lorene::Map_af::Map_af | ( | const Mg3d & | mgi, |
FILE * | fd | ||
) |
Constructor from a file (see sauve(FILE*)
)
Definition at line 431 of file map_af.C.
References alpha, beta, Lorene::fread_be(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and set_coord().
|
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 and to those of the Map_et
.
Definition at line 447 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 803 of file map_af.C.
References Lorene::c_est_pas_fait().
|
inlineinherited |
|
virtualinherited |
Computes the Spherical 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_p | [output] -component of the vector |
Implements Lorene::Map.
Definition at line 186 of file map_radial_comp_rtp.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 179 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_p_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 75 of file map_radial_comp_rtp.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 68 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_r_from_cartesian().
|
virtualinherited |
Computes the Spherical 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_t | [output] -component of the vector |
Implements Lorene::Map.
Definition at line 131 of file map_radial_comp_rtp.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 124 of file map_radial_comp_rtp.C.
References Lorene::Map_radial::comp_t_from_cartesian().
|
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] -component of the vector |
v_phi | [input] -component of the vector |
v_x | [output] x-component of the vector |
Implements Lorene::Map.
Definition at line 79 of file map_radial_comp_xyz.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 71 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_x_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] -component of the vector |
v_phi | [input] -component of the vector |
v_y | [output] y-component of the vector |
Implements Lorene::Map.
Definition at line 138 of file map_radial_comp_xyz.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 129 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_y_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] -component of the vector |
v_z | [output] z-component of the vector |
Implements Lorene::Map.
Definition at line 195 of file map_radial_comp_xyz.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Cmp
version
Implements Lorene::Map.
Definition at line 187 of file map_radial_comp_xyz.C.
References Lorene::Map_radial::comp_z_from_spherical().
|
inherited |
Determines the coordinates corresponding to given absolute Cartesian coordinates (X,Y,Z).
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 305 of file map.C.
References Lorene::Map::ori_x, Lorene::Map::ori_y, Lorene::Map::ori_z, Lorene::Map::rot_phi, and Lorene::sqrt().
|
virtual |
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 . The value of this variable gives the minimal value of (the shifted) 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 at time J+1 with boundary conditions defined by par.get_int(0) |
fJ | [input] solution at time J |
fJm1 | [input] solution at time J-1 |
source | [input] source of the d'Alembert equation . |
!!
Implements Lorene::Map.
Definition at line 125 of file map_af_dalembert.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Decreases by 2 the value of dzpuis
of a Scalar
and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 751 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Decreases by 1 the value of dzpuis
of a Scalar
and changes accordingly its values in the compactified external domain (CED).
Implements Lorene::Map.
Definition at line 646 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by of a Scalar
.
Implements Lorene::Map.
Definition at line 88 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by r of a Scalar
.
Implements Lorene::Map.
Definition at line 517 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by r (in the compactified external domain only) of a Scalar
.
Implements Lorene::Map.
Definition at line 588 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by of a Scalar
.
Implements Lorene::Map.
Definition at line 437 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by of a Scalar
.
Implements Lorene::Map.
Definition at line 136 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Division by of a Scalar
.
Implements Lorene::Map.
Definition at line 161 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
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 73 of file map_poisson_vect.C.
Computes of a 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 282 of file map_af_deriv.C.
References Lorene::Cmp::get_etat().
Computes 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.
uu | [input] field to consider |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 366 of file map_af_deriv.C.
References Lorene::Scalar::get_etat().
Computes 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.
uu | [input] field to consider |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 420 of file map_af_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a Scalar
.
uu | [input] scalar field |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 800 of file map_af_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a 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 139 of file map_af_deriv.C.
References Lorene::Cmp::get_etat().
Computes 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.
uu | [input] field to consider |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 223 of file map_af_deriv.C.
References Lorene::Scalar::get_etat().
|
inherited |
Returns the flat metric associated with the Cartesian coordinates and with components expressed in the triad bvect_cart
.
Definition at line 334 of file map.C.
References Lorene::Map::bvect_cart, and Lorene::Map::p_flat_met_cart.
|
inherited |
Returns the flat metric associated with the spherical coordinates and with components expressed in the triad bvect_spher
.
Definition at line 324 of file map.C.
References Lorene::Map::bvect_spher, and Lorene::Map::p_flat_met_spher.
const double * Lorene::Map_af::get_alpha | ( | ) | const |
const double * Lorene::Map_af::get_beta | ( | ) | const |
|
inlineinherited |
Returns the Cartesian basis associated with the coordinates (x,y,z) of the mapping, i.e.
the Cartesian coordinates related to by means of usual formulae.
Definition at line 809 of file map.h.
References Lorene::Map::bvect_cart.
|
inlineinherited |
Returns the orthonormal vectorial basis associated with the coordinates of the mapping.
Definition at line 801 of file map.h.
References Lorene::Map::bvect_spher.
|
inlineinherited |
Gives the Mg3d
on which the mapping is defined.
Definition at line 783 of file map.h.
References Lorene::Map::mg.
|
inlineinherited |
Returns the x coordinate of the origin.
Definition at line 786 of file map.h.
References Lorene::Map::ori_x.
|
inlineinherited |
Returns the y coordinate of the origin.
Definition at line 788 of file map.h.
References Lorene::Map::ori_y.
|
inlineinherited |
Returns the z coordinate of the origin.
Definition at line 790 of file map.h.
References Lorene::Map::ori_z.
|
inlineinherited |
Returns the angle between the x –axis and X –axis.
Definition at line 793 of file map.h.
References Lorene::Map::rot_phi.
|
virtual |
Sets a new radial scale.
lambda | [input] factor by which the value of r is to be multiplied |
Implements Lorene::Map.
Definition at line 667 of file map_af.C.
References Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
void Lorene::Map_af::homothetie_interne | ( | double | lambda | ) |
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 744 of file map_af.C.
References alpha, beta, and Lorene::Map_radial::reset_coord().
|
virtualinherited |
Increases by 2 the value of dzpuis
of a Scalar
and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 802 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Increases by 1 the value of dzpuis
of a Scalar
and changes accordingly its values in the
compactified external domain (CED).
Implements Lorene::Map.
Definition at line 699 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
Computes the integral over all space of a Cmp
.
The computed quantity is . The routine allocates a 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 84 of file map_af_integ.C.
References Lorene::Cmp::get_etat().
double Lorene::Map_af::integrale_surface | ( | const Cmp & | ci, |
double | rayon | ||
) | const |
Performs the surface integration of ci
on the sphere of radius rayon
.
Definition at line 99 of file map_af_integ_surf.C.
References Lorene::Cmp::get_etat().
double Lorene::Map_af::integrale_surface | ( | const Scalar & | ci, |
double | rayon | ||
) | const |
Performs the surface integration of ci
on the sphere of radius rayon
.
Definition at line 295 of file map_af_integ_surf.C.
References Lorene::Scalar::get_etat().
double Lorene::Map_af::integrale_surface_infini | ( | const Cmp & | ci | ) | const |
Performs the surface integration of ci
at infinity.
ci
must have dzpuis
=2.
Definition at line 214 of file map_af_integ_surf.C.
References Lorene::Cmp::get_etat().
double Lorene::Map_af::integrale_surface_infini | ( | const Scalar & | ci | ) | const |
Performs the surface integration of ci
at infinity.
ci
must have dzpuis
=2.
Definition at line 412 of file map_af_integ_surf.C.
References Lorene::Scalar::get_etat().
Interpolates from a Scalar
defined on a Map_star
and returns the new Scalar
defined on *this
.
f_s | [input] the field to be interpolated (must be defined on a Map_star ) |
Scalar
defined on *this
) Definition at line 179 of file map_interp_af_and_star.C.
References Lorene::Mg3d::get_type_t(), and Lorene::Map::mg.
Computes the angular Laplacian of a scalar field.
uu | [input] Scalar field u (represented as a Scalar ) the Laplacian of which is to be computed |
lap | [output] Angular Laplacian of u (see documentation of Scalar |
Implements Lorene::Map.
Definition at line 552 of file map_af_lap.C.
References Lorene::Scalar::get_etat().
Computes the Laplacian of a scalar field.
uu | [input] Scalar field u (represented as a Scalar ) the Laplacian of which is to be computed |
zec_mult_r | [input] Determines the quantity computed in the compactified external domain (CED) : \ zec_mult_r = 0 : \ zec_mult_r = 2 : \ zec_mult_r = 4 (default) : |
lap | [output] Laplacian of u |
Implements Lorene::Map.
Definition at line 182 of file map_af_lap.C.
References Lorene::Scalar::get_etat().
Computes the Laplacian of a scalar field (Cmp
version).
Implements Lorene::Map.
Definition at line 370 of file map_af_lap.C.
References Lorene::Cmp::get_etat().
|
virtual |
Returns the "angular" mapping for the outside of domain l_zone
.
Valid only for the class Map_af
.
Implements Lorene::Map.
Definition at line 786 of file map_af.C.
References Lorene::Mg3d::get_angu_1dom(), Lorene::Map::get_mg(), Map_af(), Lorene::Map::p_mp_angu, Lorene::Tbl::set(), Lorene::Tbl::set_etat_qcq(), and val_r_jk().
|
virtualinherited |
Multiplication by of a Scalar
.
Implements Lorene::Map.
Definition at line 64 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Multiplication by r of a Scalar
, the dzpuis
of uu
is not changed.
Implements Lorene::Map.
Definition at line 147 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Multiplication by r of a Cmp
.
In the CED, there is only a decrement of dzpuis
Implements Lorene::Map.
Definition at line 222 of file map_radial_r_manip.C.
References Lorene::Cmp::get_etat().
|
virtualinherited |
Multiplication by r (in the compactified external domain only) of a Scalar
.
Implements Lorene::Map.
Definition at line 299 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Multiplication by of a Scalar
.
Implements Lorene::Map.
Definition at line 358 of file map_radial_r_manip.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Multiplication by of a Scalar
.
Implements Lorene::Map.
Definition at line 112 of file map_radial_th_manip.C.
References Lorene::Scalar::get_etat().
|
virtual |
Assignment to another affine mapping.
Implements Lorene::Map_radial.
Definition at line 510 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().
|
virtual |
Comparison operator (egality)
Implements Lorene::Map_radial.
Definition at line 569 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.
|
privatevirtual |
Operator >>
Implements Lorene::Map.
Definition at line 633 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 of the Poisson equation . |
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 103 of file map_af_poisson.C.
References Lorene::Cmp::get_etat().
|
virtual |
Computes the solution of a 2-D Poisson equation.
The 2-D Poisson equation writes
source_mat | [input] Compactly supported part of the source of the 2-D Poisson equation (typically matter terms) |
source_quad | [input] Non-compactly supported part of the source of the 2-D Poisson equation (typically quadratic terms) |
par | [output] Parameter which contains the constant 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 Lorene::Map.
Definition at line 84 of file map_af_poisson2d.C.
References Lorene::Cmp::get_etat().
|
virtual |
Computes the solution of the generalized angular Poisson equation.
The generalized angular Poisson equation is , where .
source | [input] source of the equation . |
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 in the above equation (default value = 0) |
Implements Lorene::Map.
Definition at line 64 of file map_af_poisson_angu.C.
References Lorene::Scalar::get_etat().
|
virtualinherited |
Resolution of the elliptic equation in the case where the stellar interior is covered by a single domain.
source | [input] source 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 and is greater than 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 to start the iteration (if nothing is known a priori, psi must be set to zero); output: solution which satisfies . |
Implements Lorene::Map.
Definition at line 158 of file map_radial_poisson_cpt.C.
References Lorene::Cmp::get_etat().
|
virtualinherited |
Resolution of the elliptic equation in the case of a multidomain stellar interior.
nzet | [input] number of domains covering the stellar interior |
source | [input] source 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 which satisfies . |
Implements Lorene::Map.
Definition at line 456 of file map_radial_poisson_cpt.C.
References Lorene::Cmp::get_etat(), and Lorene::Map_radial::poisson_compact().
|
virtual |
Solver of the Poisson equation with boundary condition for the affine mapping case.
Implements Lorene::Map.
Definition at line 138 of file cmp_pde_frontiere.C.
References Lorene::Cmp::get_etat().
|
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 211 of file cmp_pde_frontiere.C.
References Lorene::Cmp::get_etat().
|
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 477 of file cmp_pde_frontiere.C.
References Lorene::Cmp::get_etat().
|
virtual |
Computes the solution of a scalar Poisson equation.
The regularized source is constructed and solved.
source | [input] source of the Poisson equation . |
k_div | [input] regularization degree of the procedure |
nzet | [input] number of domains covering the star |
unsgam1 | [input] parameter where 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 Lorene::Map.
Definition at line 112 of file map_af_poisson_regu.C.
References Lorene::Cmp::get_etat().
Computes the solution of a scalar Poisson equation using a Tau method.
source | [input] source of the Poisson equation . |
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 168 of file map_af_poisson.C.
References Lorene::Cmp::get_etat().
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 86 of file map_af_primr.C.
References Lorene::Scalar::get_etat(), MAX_BASE, R_CHEB, R_CHEBI, R_CHEBP, R_CHEBPIM_I, R_CHEBPIM_P, R_CHEBU, R_JACO02, R_LEG, R_LEGI, R_LEGP, and TRA_R.
|
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 64 of file map_radial_reevaluate.C.
References Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.
|
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 179 of file map_radial_reevaluate.C.
References Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.
|
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 65 of file map_radial_reeval_symy.C.
References Lorene::Cmp::get_dzpuis(), Lorene::Cmp::get_etat(), Lorene::Cmp::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.
|
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 199 of file map_radial_reeval_symy.C.
References Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Tensor::get_mp(), Lorene::Mg3d::get_nzone(), and Lorene::Map::mg.
|
protectedvirtualinherited |
Resets all the member Coords
.
Reimplemented from Lorene::Map.
Reimplemented in Lorene::Map_et.
Definition at line 129 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.
|
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.
l | [input] index of the domain |
lambda | [input] factor by which the value of defining the outer boundary of the domain is to be multiplied. |
Implements Lorene::Map.
Definition at line 690 of file map_af.C.
References Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
virtual |
Save in a file.
Reimplemented from Lorene::Map_radial.
Definition at line 619 of file map_af.C.
References alpha, beta, Lorene::fwrite_be(), Lorene::Mg3d::get_nzone(), Lorene::Map::mg, and Lorene::Map_radial::sauve().
void Lorene::Map_af::set_alpha | ( | double | alpha0, |
int | l | ||
) |
Modifies the value of in domain no. l.
Definition at line 760 of file map_af.C.
References alpha, and Lorene::Map_radial::reset_coord().
void Lorene::Map_af::set_beta | ( | double | beta0, |
int | l | ||
) |
Modifies the value of in domain no. l.
Definition at line 771 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 533 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.
|
inherited |
Sets a new origin.
Definition at line 256 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().
|
inherited |
Sets a new rotation angle.
Definition at line 266 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 103 of file map_af_elliptic.C.
References Lorene::Scalar::get_etat().
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 61 of file map_af_elliptic_2d.C.
References Lorene::Scalar::get_etat().
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 162 of file map_af_elliptic.C.
References Lorene::Scalar::get_etat().
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 225 of file map_af_elliptic.C.
References Lorene::Scalar::get_etat().
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 486 of file map_af_elliptic.C.
References Lorene::Scalar::get_etat().
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 315 of file map_af_elliptic.C.
References Lorene::Scalar::get_etat().
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 371 of file map_af_elliptic.C.
References Lorene::Scalar::get_etat().
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 62 of file map_af_elliptic_pseudo_1d.C.
References Lorene::Scalar::get_etat().
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 429 of file map_af_elliptic.C.
References Lorene::Scalar::get_etat().
Computes of a 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 478 of file map_af_deriv.C.
References Lorene::Cmp::get_etat().
Computes 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.
uu | [input] field to consider |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 572 of file map_af_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a 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 636 of file map_af_deriv.C.
References Lorene::Cmp::get_etat().
Computes 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.
uu | [input] field to consider |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 732 of file map_af_deriv.C.
References Lorene::Scalar::get_etat().
Computes of a Scalar
.
uu | [input] scalar field |
resu | [output] derivative of uu |
Implements Lorene::Map.
Definition at line 826 of file map_af_deriv.C.
References Lorene::Scalar::get_etat().
|
virtual |
Computes the domain index l and the value of corresponding to a point given by its physical coordinates .
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 131 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 corresponding to a point given by its physical coordinates .
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 195 of file map_af_radius.C.
References val_lx().
|
virtual |
Computes the domain index l and the value of corresponding to a point of arbitrary r but collocation values 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 218 of file map_af_radius.C.
References val_lx().
|
virtual |
Returns the value of the radial coordinate r for a given in a given domain.
l | [input] index of the domain |
xi | [input] value of |
theta | [input] value of |
pphi | [input] value of |
Implements Lorene::Map.
Definition at line 99 of file map_af_radius.C.
References Lorene::Mg3d::get_type_r(), and Lorene::Map::mg.
|
virtual |
Returns the value of the radial coordinate r for a given and a given collocation point in in a given domain.
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 208 of file map_af_radius.C.
References val_r().
|
private |
|
private |
|
protectedinherited |
|
protectedinherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
mutableprotectedinherited |
|
inherited |
|
inherited |
|
protectedinherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |
|
inherited |