Divergence-free vectors. More...
#include <vector.h>
Public Member Functions | |
Vector_divfree (const Map &map, const Base_vect &triad_i, const Metric &met) | |
Standard constructor. | |
Vector_divfree (const Vector_divfree &) | |
Copy constructor. | |
Vector_divfree (const Map &map, const Base_vect &triad_i, const Metric &met, FILE *fich) | |
Constructor from a file (see Tensor::sauve(FILE*) ). | |
virtual | ~Vector_divfree () |
Destructor. | |
void | operator= (const Vector_divfree &a) |
Assignment from another Vector_divfree . | |
virtual void | operator= (const Vector &a) |
Assignment from a Vector . | |
virtual void | operator= (const Tensor &a) |
Assignment from a Tensor . | |
void | set_vr_mu (const Scalar &vr_i, const Scalar &mu_i) |
Sets the angular potentials ![]() p_mu ), and the ![]() | |
void | set_vr_eta_mu (const Scalar &vr_i, const Scalar &eta_i, const Scalar &mu_i) |
Defines the components through ![]() ![]() ![]() | |
void | set_A_mu (const Scalar &A_i, const Scalar &mu_i, const Param *par_bc) |
Defines the components through potentials ![]() ![]() | |
virtual const Scalar & | eta () const |
Gives the field ![]() ![]()
| |
Vector_divfree | poisson () const |
Computes the solution of a vectorial Poisson equation with *this ![]()
| |
void | update_etavr () |
Computes the components ![]() ![]()
| |
Vector_divfree | poisson (Param &par) const |
Computes the solution of a vectorial Poisson equation with *this ![]()
| |
virtual void | change_triad (const Base_vect &) |
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly. | |
void | decompose_div (const Metric &) const |
Makes the Helmholtz decomposition (see documentation of p_potential ) of this with respect to a given Metric , only in the case of contravariant vectors. | |
const Scalar & | potential (const Metric &) const |
Returns the potential in the Helmholtz decomposition. | |
const Vector_divfree & | div_free (const Metric &) const |
Returns the div-free vector in the Helmholtz decomposition. | |
virtual void | exponential_filter_r (int lzmin, int lzmax, int p, double alpha=-16.) |
Applies exponential filters to all components (see Scalar::exponential_filter_r ). | |
virtual void | exponential_filter_ylm (int lzmin, int lzmax, int p, double alpha=-16.) |
Applies exponential filters to all components (see Scalar::exponential_filter_ylm ). | |
Scalar & | set (int) |
Read/write access to a component. | |
Scalar & | set (const Itbl &ind) |
Returns the value of a component (read/write version). | |
Scalar & | set (int i1, int i2) |
Returns the value of a component for a tensor of valence 2 (read/write version). | |
Scalar & | set (int i1, int i2, int i3) |
Returns the value of a component for a tensor of valence 3 (read/write version). | |
Scalar & | set (int i1, int i2, int i3, int i4) |
Returns the value of a component for a tensor of valence 4 (read/write version). | |
const Scalar & | operator() (int) const |
Readonly access to a component. | |
const Scalar & | operator() (const Itbl &ind) const |
Returns the value of a component (read-only version). | |
const Scalar & | operator() (int i1, int i2) const |
Returns the value of a component for a tensor of valence 2 (read-only version). | |
const Scalar & | operator() (int i1, int i2, int i3) const |
Returns the value of a component for a tensor of valence 3 (read-only version). | |
const Scalar & | operator() (int i1, int i2, int i3, int i4) const |
Returns the value of a component for a tensor of valence 4 (read-only version). | |
virtual int | position (const Itbl &idx) const |
Returns the position in the Scalar array cmp of a component given by its index. | |
virtual Itbl | indices (int place) const |
Returns the index of a component given by its position in the Scalar array cmp . | |
virtual void | std_spectral_base () |
Sets the standard spectal bases of decomposition for each component. | |
virtual void | pseudo_spectral_base () |
Sets the standard spectal bases of decomposition for each component for a pseudo_vector. | |
virtual const Scalar & | mu () const |
Gives the field ![]() ![]()
| |
virtual const Scalar & | A () const |
Gives the field ![]()
Related to the curl, A is insensitive to the longitudinal part of the vector. | |
void | update_vtvp () |
Computes the components ![]() ![]() ![]() ![]()
| |
const Scalar & | divergence (const Metric &) const |
The divergence of this with respect to a Metric . | |
const Vector_divfree | curl () const |
The curl of this with respect to a (flat) Metric . | |
Vector | derive_lie (const Vector &v) const |
Computes the Lie derivative of this with respect to some vector field v . | |
Tensor | derive_lie (const Vector &v) const |
Computes the Lie derivative of this with respect to some vector field v . | |
Sym_tensor | ope_killing (const Metric &gam) const |
Computes the Killing operator associated with a given metric. | |
Sym_tensor | ope_killing_conf (const Metric &gam) const |
Computes the conformal Killing operator associated with a given metric. | |
Vector | poisson (double lambda, int method=6) const |
Solves the vector Poisson equation with *this as a source. | |
Vector | poisson (double lambda, const Metric_flat &met_f, int method=6) const |
Solves the vector Poisson equation with *this as a source. | |
Vector | poisson (const double lambda, Param &par, int method=6) const |
Solves the vector Poisson equation with *this as a source and parameters controlling the solution. | |
void | poisson_boundary (double lambda, const Mtbl_cf &limit_vr, const Mtbl_cf &limit_eta, const Mtbl_cf &limit_mu, int num_front, double fact_dir, double fact_neu, Vector &resu) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere. | |
void | poisson_boundary2 (double lam, Vector &resu, Scalar boundvr, Scalar boundeta, Scalar boundmu, double dir_vr, double neum_vr, double dir_eta, double neum_eta, double dir_mu, double neum_mu) const |
Alternative to previous poisson_boundary method for vectors ; this uses method 6 for vectorial solving, updated version (as in the poisson_vector_block routine). | |
Vector | poisson_dirichlet (double lambda, const Valeur &limit_vr, const Valeur &limit_vt, const Valeur &limit_vp, int num_front) const |
Vector | poisson_neumann (double lambda, const Valeur &limit_vr, const Valeur &limit_vt, const Valeur &limit_vp, int num_front) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere. | |
Vector | poisson_robin (double lambda, const Valeur &limit_vr, const Valeur &limit_vt, const Valeur &limit_vp, double fact_dir, double fact_neu, int num_front) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere. | |
double | flux (double radius, const Metric &met) const |
Computes the flux of the vector accross a sphere r = const. | |
void | poisson_block (double lambda, Vector &resu) const |
void | visu_arrows (double xmin, double xmax, double ymin, double ymax, double zmin, double zmax, const char *title0=0x0, const char *filename0=0x0, bool start_dx=true, int nx=8, int ny=8, int nz=8) const |
3D visualization via OpenDX. | |
void | visu_streamline (double xmin, double xmax, double ymin, double ymax, double zmin, double zmax, const char *title0=0x0, const char *filename0=0x0, bool start_dx=true, int nx=8, int ny=8, int nz=8) const |
virtual void | set_etat_nondef () |
Sets the logical state of all components to ETATNONDEF (undefined state). | |
virtual void | set_etat_zero () |
Sets the logical state of all components to ETATZERO (zero state). | |
virtual void | set_etat_qcq () |
Sets the logical state of all components to ETATQCQ (ordinary state). | |
virtual void | allocate_all () |
Performs the memory allocation of all the elements, down to the double arrays of the Tbl s. | |
void | set_triad (const Base_vect &new_triad) |
Assigns a new vectorial basis (triad) of decomposition. | |
void | annule_domain (int l) |
Sets the Tensor to zero in a given domain. | |
virtual void | annule (int l_min, int l_max) |
Sets the Tensor to zero in several domains. | |
void | annule_extern_cn (int l_0, int deg) |
Performs a smooth (C^n) transition in a given domain to zero. | |
virtual void | std_spectral_base_odd () |
Sets the standard odd spectal bases of decomposition for each component. | |
virtual void | dec_dzpuis (int dec=1) |
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified external domain (CED). | |
virtual void | inc_dzpuis (int inc=1) |
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified external domain (CED). | |
const Tensor & | derive_cov (const Metric &gam) const |
Returns the covariant derivative of this with respect to some metric ![]() | |
const Tensor & | derive_con (const Metric &gam) const |
Returns the "contravariant" derivative of this with respect to some metric ![]() | |
Tensor | up (int ind, const Metric &gam) const |
Computes a new tensor by raising an index of *this . | |
Tensor | down (int ind, const Metric &gam) const |
Computes a new tensor by lowering an index of *this . | |
Tensor | up_down (const Metric &gam) const |
Computes a new tensor by raising or lowering all the indices of *this . | |
Tensor | trace (int ind1, int ind2) const |
Trace on two different type indices. | |
Tensor | trace (int ind1, int ind2, const Metric &gam) const |
Trace with respect to a given metric. | |
Scalar | trace () const |
Trace on two different type indices for a valence 2 tensor. | |
Scalar | trace (const Metric &gam) const |
Trace with respect to a given metric for a valence 2 tensor. | |
const Map & | get_mp () const |
Returns the mapping. | |
const Base_vect * | get_triad () const |
Returns the vectorial basis (triad) on which the components are defined. | |
int | get_valence () const |
Returns the valence. | |
int | get_n_comp () const |
Returns the number of stored components. | |
int | get_index_type (int i) const |
Gives the type (covariant or contravariant) of the index number i . | |
Itbl | get_index_type () const |
Returns the types of all the indices. | |
int & | set_index_type (int i) |
Sets the type of the index number i . | |
Itbl & | set_index_type () |
Sets the types of all the indices. | |
void | operator+= (const Tensor &) |
+= Tensor | |
void | operator-= (const Tensor &) |
-= Tensor | |
virtual void | sauve (FILE *) const |
Save in a binary file. | |
virtual void | spectral_display (const char *comment=0x0, double threshold=1.e-7, int precision=4, ostream &ostr=cout) const |
Displays the spectral coefficients and the associated basis functions of each component. | |
Protected Member Functions | |
virtual void | del_deriv () const |
Deletes the derived quantities. | |
void | set_der_0x0 () const |
Sets the pointers on derived quantities to 0x0. | |
void | sol_Dirac_A (const Scalar &aaa, Scalar &eta, Scalar &vr, const Param *par_bc=0x0) const |
Solves a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value. | |
void | sol_Dirac_A_tau (const Scalar &aaa, Scalar &eta, Scalar &vr, const Param *par_bc=0x0) const |
Solves via a tau method a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value. | |
void | sol_Dirac_A_poisson (const Scalar &aaa, Scalar &eta, Scalar &vr, const Param *par_bc=0x0) const |
Solves via a poisson method a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value. | |
void | sol_Dirac_A_1z (const Scalar &aaa, Scalar &eta, Scalar &vr, const Param *par_bc=0x0) const |
Solves a one-domain system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value. | |
virtual void | del_derive_met (int) const |
Logical destructor of the derivatives depending on the i-th element of met_depend in the class Vector . | |
void | set_der_met_0x0 (int) const |
Sets all the i-th components of met_depend in the class Vector (p_potential , etc. | |
void | set_dependance (const Metric &) const |
To be used to describe the fact that the derivatives members have been calculated with met . | |
int | get_place_met (const Metric &) const |
Returns the position of the pointer on metre in the array met_depend . | |
void | compute_derive_lie (const Vector &v, Tensor &resu) const |
Computes the Lie derivative of this with respect to some vector field v (protected method; the public interface is method derive_lie ). | |
Protected Attributes | |
const Metric *const | met_div |
Metric with respect to which the divergence is defined. | |
Scalar * | p_potential [N_MET_MAX] |
The potential ![]() ![]() | |
Vector_divfree * | p_div_free [N_MET_MAX] |
The divergence-free vector ![]() ![]() | |
Scalar * | p_eta |
Field ![]() ![]()
| |
Scalar * | p_mu |
Field ![]() ![]()
| |
Scalar * | p_A |
Field ![]()
Insensitive to the longitudinal part of the vector, related to the curl. | |
const Map *const | mp |
Mapping on which the numerical values at the grid points are defined. | |
int | valence |
Valence of the tensor (0 = scalar, 1 = vector, etc...). | |
const Base_vect * | triad |
Vectorial basis (triad) with respect to which the tensor components are defined. | |
Itbl | type_indice |
1D array of integers (class Itbl ) of size valence containing the type of each index: COV for a covariant one and CON for a contravariant one. | |
int | n_comp |
Number of stored components, depending on the symmetry. | |
Scalar ** | cmp |
Array of size n_comp of pointers onto the components. | |
const Metric * | met_depend [N_MET_MAX] |
Array on the Metric 's which were used to compute derived quantities, like p_derive_cov , etc. | |
Tensor * | p_derive_cov [N_MET_MAX] |
Array of pointers on the covariant derivatives of this with respect to various metrics. | |
Tensor * | p_derive_con [N_MET_MAX] |
Array of pointers on the contravariant derivatives of this with respect to various metrics. | |
Tensor * | p_divergence [N_MET_MAX] |
Array of pointers on the divergence of this with respect to various metrics. | |
Friends | |
class | Vector |
class | Scalar |
class | Sym_tensor |
class | Tensor_sym |
class | Metric |
ostream & | operator<< (ostream &, const Tensor &) |
Scalar | operator+ (const Tensor &, const Scalar &) |
Scalar | operator+ (const Scalar &, const Tensor &) |
Scalar | operator- (const Tensor &, const Scalar &) |
Scalar | operator- (const Scalar &, const Tensor &) |
Tensor | operator* (const Tensor &, const Tensor &) |
Tensor_sym | operator* (const Tensor &, const Tensor_sym &) |
Tensor_sym | operator* (const Tensor_sym &, const Tensor &) |
Tensor_sym | operator* (const Tensor_sym &, const Tensor_sym &) |
Tensorial product of two symmetric tensors. |
Divergence-free vectors.
()
This class is designed to store divergence-free vectors, with the component expressed in a orthonormal spherical basis .
Definition at line 720 of file vector.h.
Standard constructor.
map | the mapping | |
triad_i | vectorial basis (triad) with respect to which the vector components are defined | |
met | the metric with respect to which the divergence is defined |
Definition at line 75 of file vector_divfree.C.
References set_der_0x0().
Vector_divfree::Vector_divfree | ( | const Vector_divfree & | source | ) |
Copy constructor.
Definition at line 87 of file vector_divfree.C.
References Vector::p_eta, Vector::p_mu, and set_der_0x0().
Vector_divfree::Vector_divfree | ( | const Map & | map, | |
const Base_vect & | triad_i, | |||
const Metric & | met, | |||
FILE * | fich | |||
) |
Constructor from a file (see Tensor::sauve(FILE*)
).
map | the mapping | |
triad_i | vectorial basis (triad) with respect to which the tensor components are defined. It will be checked that it coincides with the basis saved in the file. | |
met | the metric with respect to which the divergence is defined | |
fich | file which has been created by the function sauve(FILE*) . |
Definition at line 99 of file vector_divfree.C.
References set_der_0x0().
Vector_divfree::~Vector_divfree | ( | ) | [virtual] |
const Scalar & Vector::A | ( | ) | const [virtual, inherited] |
Gives the field defined by
Related to the curl, A is insensitive to the longitudinal part of the vector.
Definition at line 125 of file vector_etamu.C.
References Tensor::cmp, Scalar::div_r_dzpuis(), Scalar::dsdr(), Vector::eta(), Vector::p_A, Vector::p_eta, Scalar::set_dzpuis(), and Tensor::triad.
void Tensor::allocate_all | ( | ) | [virtual, inherited] |
Performs the memory allocation of all the elements, down to the double
arrays of the Tbl
s.
This function performs in fact recursive calls to set_etat_qcq()
on each element of the chain Scalar
-> Valeur
-> Mtbl
-> Tbl
.
Reimplemented in Scalar.
Definition at line 504 of file tensor.C.
References Scalar::allocate_all(), Tensor::cmp, Tensor::del_deriv(), and Tensor::n_comp.
void Tensor::annule | ( | int | l_min, | |
int | l_max | |||
) | [virtual, inherited] |
Sets the Tensor
to zero in several domains.
l_min | [input] The Tensor will be set (logically) to zero in the domains whose indices are in the range [l_min,l_max] . | |
l_max | [input] see the comments for l_min . |
Note that annule
(0,nz-1) , where nz
is the total number of domains, is equivalent to set_etat_zero()
.
Reimplemented in Scalar.
Definition at line 667 of file tensor.C.
References Scalar::annule(), Tensor::cmp, Tensor::del_deriv(), Map::get_mg(), Mg3d::get_nzone(), Tensor::mp, Tensor::n_comp, and Tensor::set_etat_zero().
void Tensor::annule_domain | ( | int | l | ) | [inherited] |
Sets the Tensor
to zero in a given domain.
l | [input] Index of the domain in which the Tensor will be set (logically) to zero. |
Definition at line 662 of file tensor.C.
References Tensor::annule().
void Tensor::annule_extern_cn | ( | int | l_0, | |
int | deg | |||
) | [inherited] |
Performs a smooth (C^n) transition in a given domain to zero.
l_0 | [input] in the domain of index l0 the tensor is multiplied by the right polynomial (of degree 2n+1), to ensure continuty of the function and its n first derivative at both ends of this domain. The tensor is unchanged in the domains l < l_0 and set to zero in domains l > l_0. | |
deg | [input] the degree n of smoothness of the transition. |
Definition at line 686 of file tensor.C.
References Scalar::allocate_all(), Scalar::annule(), Itbl::annule_hard(), Tensor::cmp, Tensor::del_deriv(), Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nzone(), Mg3d::get_type_r(), Tensor::mp, Tensor::n_comp, pow(), Map::r, Tbl::set(), Itbl::set(), Scalar::set_domain(), Tbl::set_etat_qcq(), Scalar::std_spectral_base(), and Map::val_r().
void Vector::change_triad | ( | const Base_vect & | new_triad | ) | [virtual, inherited] |
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Reimplemented from Tensor.
Definition at line 71 of file vector_change_triad.C.
References Tensor::cmp, Map::comp_p_from_cartesian(), Map::comp_r_from_cartesian(), Map::comp_t_from_cartesian(), Map::comp_x_from_spherical(), Map::comp_y_from_spherical(), Map::comp_z_from_spherical(), Base_vect_cart::get_align(), Map::get_bvect_cart(), Map::get_bvect_spher(), Map::get_mg(), Mg3d::get_nt(), Mg3d::get_nzone(), Tensor::mp, and Tensor::triad.
Computes the Lie derivative of this
with respect to some vector field v
(protected method; the public interface is method derive_lie
).
Definition at line 335 of file tensor_calculus.C.
References Tensor::cmp, contract(), Scalar::dec_dzpuis(), Tensor::derive_cov(), Map::flat_met_cart(), Map::flat_met_spher(), Scalar::get_dzpuis(), Tensor::get_n_comp(), Tensor::get_triad(), Tensor::indices(), Tensor::mp, Tensor::n_comp, Tensor::operator()(), Tensor::set(), Itbl::set(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
const Vector_divfree Vector::curl | ( | ) | const [inherited] |
The curl of this
with respect to a (flat) Metric
.
The Vector
is assumed to be contravariant.
Definition at line 398 of file vector.C.
References Tensor::cmp, Scalar::div_r(), Scalar::div_tant(), Scalar::dsdt(), Scalar::dsdx(), Scalar::dsdy(), Scalar::dsdz(), Map::flat_met_cart(), Map::flat_met_spher(), Map::get_bvect_cart(), Map::get_bvect_spher(), Base_vect::identify(), Tensor::mp, Scalar::mult_r(), Tensor::set(), Scalar::srstdsdp(), and Tensor::triad.
void Tensor::dec_dzpuis | ( | int | dec = 1 |
) | [virtual, inherited] |
Decreases by dec
units the value of dzpuis
and changes accordingly the values in the compactified external domain (CED).
Reimplemented in Scalar.
Definition at line 804 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), and Tensor::n_comp.
void Vector::decompose_div | ( | const Metric & | metre | ) | const [inherited] |
Makes the Helmholtz decomposition (see documentation of p_potential
) of this
with respect to a given Metric
, only in the case of contravariant vectors.
Definition at line 512 of file vector.C.
References Scalar::annule(), Tensor::annule_domain(), Valeur::base, Valeur::c_cf, Scalar::check_dzpuis(), Tensor::cmp, Valeur::coef(), Tensor::dec_dzpuis(), Scalar::derive_con(), Vector::divergence(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Tensor::get_place_met(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Matrice::inverse(), Tensor::mp, Vector::p_div_free, Vector::p_potential, Scalar::poisson(), pow(), Map::r, R_CHEB, R_CHEBI, R_CHEBP, Scalar::Scalar(), Mtbl_cf::set(), Tbl::set(), Matrice::set(), Matrice::set_band(), Valeur::set_base_r(), Tensor::set_dependance(), Scalar::set_dzpuis(), Matrice::set_etat_qcq(), Tbl::set_etat_qcq(), Scalar::set_etat_zero(), Matrice::set_lu(), Scalar::set_spectral_va(), Scalar::std_spectral_base(), Tensor::triad, Tensor::type_indice, Valeur::ylm(), and Valeur::ylm_i().
void Vector_divfree::del_deriv | ( | ) | const [protected, virtual] |
Deletes the derived quantities.
Reimplemented from Vector.
Definition at line 124 of file vector_divfree.C.
References set_der_0x0().
void Vector::del_derive_met | ( | int | j | ) | const [protected, virtual, inherited] |
Logical destructor of the derivatives depending on the i-th element of met_depend
in the class Vector
.
Reimplemented from Tensor.
Definition at line 238 of file vector.C.
References Tensor::met_depend, Vector::p_div_free, Vector::p_potential, and Vector::set_der_met_0x0().
Returns the "contravariant" derivative of this
with respect to some metric , by raising the last index of the covariant derivative (cf.
method derive_cov()
) with .
Reimplemented in Scalar, and Tensor_sym.
Definition at line 1010 of file tensor.C.
References Metric::con(), contract(), Tensor::derive_cov(), Tensor::get_index_type(), Tensor::get_place_met(), Tensor::mp, Tensor::p_derive_con, Itbl::set(), Tensor::set_dependance(), Tensor_sym::sym_index1(), Tensor_sym::sym_index2(), Tensor::Tensor(), Tensor::triad, and Tensor::valence.
Returns the covariant derivative of this
with respect to some metric .
denoting the tensor represented by
this
and its covariant derivative with respect to the metric
, the extra index (with respect to the indices of
) of
is chosen to be the last one. This convention agrees with that of MTW (see Eq. (10.17) of MTW). For instance, if
is a 1-form, whose components w.r.t. the triad
are
:
, then the covariant derivative of
is the bilinear form
whose components
are such that
gam | metric ![]() |
this
with respect to the connection Reimplemented in Scalar, and Tensor_sym.
Definition at line 998 of file tensor.C.
References Metric::connect(), Tensor::get_place_met(), Connection::p_derive_cov(), Tensor::p_derive_cov, and Tensor::set_dependance().
Computes the Lie derivative of this
with respect to some vector field v
.
Reimplemented in Scalar, Sym_tensor, and Tensor_sym.
Definition at line 477 of file tensor_calculus.C.
References Tensor::compute_derive_lie(), Tensor::mp, Tensor::triad, Tensor::type_indice, and Tensor::valence.
Computes the Lie derivative of this
with respect to some vector field v
.
Definition at line 388 of file vector.C.
References Tensor::compute_derive_lie(), Tensor::mp, Tensor::triad, and Tensor::type_indice.
const Vector_divfree & Vector::div_free | ( | const Metric & | metre | ) | const [inherited] |
Returns the div-free vector in the Helmholtz decomposition.
It first makes the Helmholtz decomposition (see documentation of p_potential
) of this
with respect to a given Metric
and then returns . Only in the case of contravariant vectors.
Definition at line 500 of file vector.C.
References Vector::decompose_div(), Tensor::get_place_met(), Vector::p_div_free, and Tensor::set_dependance().
Computes a new tensor by lowering an index of *this
.
ind | index to be lowered, with the following convention :
| |
gam | metric used to lower the index (contraction with the twice covariant form of the metric on the index ind ). |
Definition at line 261 of file tensor_calculus.C.
References contract(), Metric::cov(), Tensor::indices(), Tensor::mp, Tensor::n_comp, Tensor::set(), Itbl::set(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
const Scalar & Vector_divfree::eta | ( | ) | const [virtual] |
Gives the field such that the angular components
of the vector are written:
.
Reimplemented from Vector.
Definition at line 189 of file vector_divfree.C.
References Tensor::cmp, Scalar::dsdr(), Scalar::get_dzpuis(), Scalar::mult_r_dzpuis(), Vector::p_eta, Scalar::poisson_angu(), and Tensor::triad.
void Vector::exponential_filter_r | ( | int | lzmin, | |
int | lzmax, | |||
int | p, | |||
double | alpha = -16. | |||
) | [virtual, inherited] |
Applies exponential filters to all components (see Scalar::exponential_filter_r
).
Does a loop for Cartesian components, and works in terms of the r-component, and
for spherical components.
Reimplemented from Tensor.
Definition at line 846 of file vector.C.
References Tensor::cmp, Vector::eta(), Scalar::exponential_filter_r(), Map::get_bvect_cart(), Map::get_bvect_spher(), Base_vect::identify(), Tensor::mp, Vector::mu(), Tensor::n_comp, Vector::operator()(), Vector::set_vr_eta_mu(), and Tensor::triad.
void Vector::exponential_filter_ylm | ( | int | lzmin, | |
int | lzmax, | |||
int | p, | |||
double | alpha = -16. | |||
) | [virtual, inherited] |
Applies exponential filters to all components (see Scalar::exponential_filter_ylm
).
Does a loop for Cartesian components, and works in terms of the r-component, and
for spherical components.
Reimplemented from Tensor.
Definition at line 863 of file vector.C.
References Tensor::cmp, Vector::eta(), Scalar::exponential_filter_ylm(), Map::get_bvect_cart(), Map::get_bvect_spher(), Base_vect::identify(), Tensor::mp, Vector::mu(), Tensor::n_comp, Vector::operator()(), Vector::set_vr_eta_mu(), and Tensor::triad.
double Vector::flux | ( | double | radius, | |
const Metric & | met | |||
) | const [inherited] |
Computes the flux of the vector accross a sphere r = const.
radius | radius of the sphere S on which the flux is to be taken; the center of S is assumed to be the center of the mapping (member mp ). radius can take the value __infinity (to get the flux at spatial infinity). | |
met | metric ![]() |
*this
and Definition at line 803 of file vector.C.
References Vector::change_triad(), Map::get_bvect_spher(), Map_af::integrale_surface(), Map_af::integrale_surface_infini(), Tensor::mp, Tensor::triad, Tensor::type_indice, and Vector::Vector().
Itbl Tensor::get_index_type | ( | ) | const [inline, inherited] |
Returns the types of all the indices.
Itbl
) of size valence
containing the type of each index, COV
for a covariant one and CON
for a contravariant one. Definition at line 892 of file tensor.h.
References Tensor::type_indice.
int Tensor::get_index_type | ( | int | i | ) | const [inline, inherited] |
Gives the type (covariant or contravariant) of the index number i
.
i
must be strictly lower than valence
and obey the following convention:
i
= 0 : first index i
= 1 : second index Definition at line 882 of file tensor.h.
References Tensor::type_indice.
const Map& Tensor::get_mp | ( | ) | const [inline, inherited] |
int Tensor::get_n_comp | ( | ) | const [inline, inherited] |
Returns the number of stored components.
Definition at line 868 of file tensor.h.
References Tensor::n_comp.
int Tensor::get_place_met | ( | const Metric & | metre | ) | const [protected, inherited] |
Returns the position of the pointer on metre
in the array met_depend
.
Definition at line 439 of file tensor.C.
References Tensor::met_depend.
const Base_vect* Tensor::get_triad | ( | ) | const [inline, inherited] |
Returns the vectorial basis (triad) on which the components are defined.
Definition at line 862 of file tensor.h.
References Tensor::triad.
int Tensor::get_valence | ( | ) | const [inline, inherited] |
void Tensor::inc_dzpuis | ( | int | inc = 1 |
) | [virtual, inherited] |
Increases by inc
units the value of dzpuis
and changes accordingly the values in the compactified external domain (CED).
Reimplemented in Scalar.
Definition at line 812 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), and Tensor::n_comp.
virtual Itbl Vector::indices | ( | int | place | ) | const [inline, virtual, inherited] |
Returns the index of a component given by its position in the Scalar
array cmp
.
Itbl
) of size 1 giving its value for the component located at the position place
in the Scalar
array cmp
. The element of this Itbl
corresponds to a spatial index 1, 2 or 3. Reimplemented from Tensor.
const Scalar & Vector::mu | ( | ) | const [virtual, inherited] |
Gives the field such that the angular components
of the vector are written:
.
Definition at line 94 of file vector_etamu.C.
References Tensor::cmp, Scalar::div_tant(), Scalar::dsdt(), Vector::p_mu, Scalar::poisson_angu(), Scalar::stdsdp(), and Tensor::triad.
Sym_tensor Vector::ope_killing | ( | const Metric & | gam | ) | const [inherited] |
Computes the Killing operator associated with a given metric.
The Killing operator is defined by for a contravariant vector and by
for a covariant vector.
gam | metric with respect to which the covariant derivative ![]() |
Definition at line 434 of file vector.C.
References Tensor::derive_con(), Tensor::derive_cov(), Tensor::mp, Tensor::set(), Tensor::triad, and Tensor::type_indice.
Sym_tensor Vector::ope_killing_conf | ( | const Metric & | gam | ) | const [inherited] |
Computes the conformal Killing operator associated with a given metric.
The conformal Killing operator is defined by for a contravariant vector and by
for a covariant vector.
gam | metric ![]() ![]() |
Definition at line 458 of file vector.C.
References Metric::con(), Metric::cov(), Tensor::derive_con(), Tensor::derive_cov(), Vector::divergence(), Tensor::mp, Tensor::set(), Tensor::trace(), Tensor::triad, and Tensor::type_indice.
const Scalar & Tensor::operator() | ( | int | i1, | |
int | i2, | |||
int | i3, | |||
int | i4 | |||
) | const [inherited] |
Returns the value of a component for a tensor of valence 4 (read-only version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) | |
i3 | value of the third index (1, 2 or 3) | |
i4 | value of the fourth index (1, 2 or 3) |
(i1,i2,i3,i4) Definition at line 779 of file tensor.C.
References Tensor::cmp, Tensor::position(), Itbl::set(), and Tensor::valence.
const Scalar & Tensor::operator() | ( | int | i1, | |
int | i2, | |||
int | i3 | |||
) | const [inherited] |
Returns the value of a component for a tensor of valence 3 (read-only version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) | |
i3 | value of the third index (1, 2 or 3) |
(i1,i2,i3) Definition at line 767 of file tensor.C.
References Tensor::cmp, Tensor::position(), Itbl::set(), and Tensor::valence.
const Scalar & Tensor::operator() | ( | int | i1, | |
int | i2 | |||
) | const [inherited] |
Returns the value of a component for a tensor of valence 2 (read-only version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) |
(i1,i2) Definition at line 756 of file tensor.C.
References Tensor::cmp, Tensor::position(), Itbl::set(), and Tensor::valence.
Returns the value of a component (read-only version).
ind | 1-D Itbl of size valence containing the values of each index specifing the component, with the following storage convention:
|
ind
Definition at line 794 of file tensor.C.
References Tensor::cmp, Itbl::get_dim(), Itbl::get_ndim(), Tensor::position(), and Tensor::valence.
const Scalar & Vector::operator() | ( | int | index | ) | const [inherited] |
void Tensor::operator+= | ( | const Tensor & | t | ) | [inherited] |
+= Tensor
Reimplemented in Scalar.
Definition at line 567 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::indices(), Tensor::n_comp, Tensor::position(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
void Tensor::operator-= | ( | const Tensor & | t | ) | [inherited] |
-= Tensor
Reimplemented in Scalar.
Definition at line 583 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::indices(), Tensor::n_comp, Tensor::position(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
void Vector_divfree::operator= | ( | const Tensor & | a | ) | [virtual] |
Assignment from a Tensor
.
Reimplemented from Vector.
Definition at line 159 of file vector_divfree.C.
References del_deriv(), and operator=().
void Vector_divfree::operator= | ( | const Vector & | a | ) | [virtual] |
Assignment from a Vector
.
Definition at line 149 of file vector_divfree.C.
References del_deriv(), and operator=().
void Vector_divfree::operator= | ( | const Vector_divfree & | a | ) | [virtual] |
Assignment from another Vector_divfree
.
Reimplemented from Vector.
Definition at line 137 of file vector_divfree.C.
References del_deriv(), and met_div.
Solves the vector Poisson equation with *this
as a source and parameters controlling the solution.
The equatiopn solved is .
*this
must be given with dzpuis
= 4. It uses the Helmholtz decomposition (see documentation of p_potential
), with a flat metric, deduced from the triad.
lambda | [input] ![]() | |
par | [input/output] possible parameters | |
uu | [input/output] solution u with the boundary condition u =0 at spatial infinity. |
Definition at line 539 of file vector_poisson.C.
References Param::add_cmp_mod(), Param::add_double(), Param::add_int(), Param::add_int_mod(), Cmp::annule_hard(), Tensor::cmp, Tensor::dec_dzpuis(), Scalar::derive_con(), Vector::div_free(), Map::get_bvect_cart(), Map::get_bvect_spher(), Param::get_cmp_mod(), Param::get_double(), Param::get_int(), Param::get_int_mod(), Scalar::inc_dzpuis(), Tensor::mp, poisson(), Scalar::poisson(), Vector::potential(), Vector::set(), Tenseur::set(), Tenseur::set_etat_qcq(), Tensor::set_etat_qcq(), Scalar::set_etat_zero(), Tenseur::set_std_base(), and Tensor::triad.
Vector Vector::poisson | ( | double | lambda, | |
const Metric_flat & | met_f, | |||
int | method = 6 | |||
) | const [inherited] |
Solves the vector Poisson equation with *this
as a source.
The equation solved is .
*this
must be given with dzpuis
= 4. It uses the Helmholtz decomposition (see documentation of p_potential
), with the flat metric met_f
given in argument.
lambda | [input] ![]() | |
met_f | [input] the flat metric for the Helmholtz decomposition. | |
method | [input] method used to solve the equation:
|
Definition at line 126 of file vector_poisson.C.
References Scalar::annule_hard(), Valeur::base, Valeur::c, Valeur::c_cf, Tensor::cmp, Valeur::coef(), Scalar::dec_dzpuis(), Tensor::dec_dzpuis(), Scalar::derive_con(), Vector::div_free(), Scalar::div_r_dzpuis(), Scalar::div_sint(), Scalar::dsdr(), Scalar::dsdt(), Map::get_bvect_cart(), Map::get_bvect_spher(), Scalar::get_etat(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Mg3d::get_type_r(), Scalar::lapang(), Tensor::mp, Vector::mu(), Scalar::mult_r_dzpuis(), Scalar::mult_sint(), Vector::poisson(), Scalar::poisson(), Scalar::poisson_angu(), Vector::potential(), Scalar::primr(), Vector::set(), Mtbl_cf::set(), Scalar::set_dzpuis(), Tenseur::set_etat_qcq(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Tensor::set_etat_zero(), Param_elliptic::set_poisson_vect_eta(), Param_elliptic::set_poisson_vect_r(), Scalar::set_spectral_va(), Vector::set_vr_eta_mu(), Scalar::sol_elliptic(), Scalar::stdsdp(), Tensor::triad, Valeur::ylm(), and Valeur::ylm_i().
Vector Vector::poisson | ( | double | lambda, | |
int | method = 6 | |||
) | const [inherited] |
Solves the vector Poisson equation with *this
as a source.
The equation solved is .
*this
must be given with dzpuis
= 4. It uses the Helmholtz decomposition (see documentation of p_potential
), with a flat metric, deduced from the triad.
lambda | [input] ![]() | |
method | [input] method used to solve the equation (see Vector::poisson(double, Metric_flat, int) for details). |
Definition at line 514 of file vector_poisson.C.
References Map::flat_met_cart(), Map::flat_met_spher(), Tensor::mp, and Tensor::triad.
Vector_divfree Vector_divfree::poisson | ( | Param & | par | ) | const |
Computes the solution of a vectorial Poisson equation with *this
as a source:
.
Definition at line 65 of file vector_df_poisson.C.
References Param::add_cmp_mod(), Param::add_double(), Param::add_int(), Param::add_int_mod(), Tensor::cmp, Scalar::div_r(), Param::get_cmp_mod(), Param::get_double(), Param::get_int(), Param::get_int_mod(), met_div, Tensor::mp, Vector::mu(), Scalar::mult_r(), Scalar::poisson(), Scalar::set_etat_zero(), set_vr_mu(), and Tensor::triad.
Vector_divfree Vector_divfree::poisson | ( | ) | const |
Computes the solution of a vectorial Poisson equation with *this
as a source:
.
Definition at line 142 of file vector_df_poisson.C.
References Tbl::annule_hard(), Mtbl_cf::annule_hard(), Valeur::base, Valeur::c_cf, Tensor::cmp, eta(), Map_af::get_alpha(), Map_af::get_beta(), Diff_dsdx::get_matrice(), Diff_dsdx2::get_matrice(), Diff_xdsdx2::get_matrice(), Diff_xdsdx::get_matrice(), Diff_x2dsdx2::get_matrice(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), met_div, Tensor::mp, Vector::mu(), Scalar::poisson(), pow(), Mtbl_cf::set(), Tbl::set(), Matrice::set(), Matrice::set_band(), Valeur::set_etat_cf_qcq(), Scalar::set_etat_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Scalar::set_etat_zero(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), set_vr_eta_mu(), set_vr_mu(), Tbl::t, Tensor::triad, Valeur::ylm(), and Valeur::ylm_i().
void Vector::poisson_boundary | ( | double | lambda, | |
const Mtbl_cf & | limit_vr, | |||
const Mtbl_cf & | limit_eta, | |||
const Mtbl_cf & | limit_mu, | |||
int | num_front, | |||
double | fact_dir, | |||
double | fact_neu, | |||
Vector & | resu | |||
) | const [inherited] |
Solves the vector Poisson equation with *this
as a source with a boundary condition on the excised sphere.
The equation solved is .
*this
must be given with dzpuis
= 4. It uses the Helmholtz decomposition (see documentation of p_potential
)
lambda | [input] ![]() | |
resu | [output] the solution ![]() |
Definition at line 61 of file vector_poisson_boundary.C.
References Tbl::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::base, Valeur::c_cf, Tensor::cmp, Vector::eta(), Map_af::get_alpha(), Mg3d::get_angu(), Map_af::get_beta(), Scalar::get_etat(), Diff_dsdx::get_matrice(), Diff_xdsdx::get_matrice(), Diff_dsdx2::get_matrice(), Diff_xdsdx2::get_matrice(), Diff_x2dsdx2::get_matrice(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Vector::mu(), pow(), Mtbl_cf::set(), Tbl::set(), Matrice::set(), Matrice::set_band(), Valeur::set_etat_cf_qcq(), Scalar::set_etat_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Matrice::set_lu(), Param_elliptic::set_poisson_vect_r(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Vector::set_vr_eta_mu(), Scalar::sol_elliptic_boundary(), Tbl::t, Tensor::triad, Valeur::ylm(), and Valeur::ylm_i().
void Vector::poisson_boundary2 | ( | double | lam, | |
Vector & | resu, | |||
Scalar | boundvr, | |||
Scalar | boundeta, | |||
Scalar | boundmu, | |||
double | dir_vr, | |||
double | neum_vr, | |||
double | dir_eta, | |||
double | neum_eta, | |||
double | dir_mu, | |||
double | neum_mu | |||
) | const [inherited] |
Alternative to previous poisson_boundary method for vectors ; this uses method 6 for vectorial solving, updated version (as in the poisson_vector_block routine).
Boundary arguments are here required as scalar fields.
Definition at line 76 of file vector_poisson_boundary2.C.
References Matrice::annule_hard(), Tbl::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::base, Valeur::c_cf, Tensor::cmp, Vector::eta(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Diff_sx2::get_matrice(), Diff_sxdsdx::get_matrice(), Diff_dsdx::get_matrice(), Diff_xdsdx::get_matrice(), Diff_dsdx2::get_matrice(), Diff_xdsdx2::get_matrice(), Diff_x2dsdx2::get_matrice(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Vector::mu(), Mtbl_cf::set(), Tbl::set(), Matrice::set(), Scalar::set_dzpuis(), Valeur::set_etat_cf_qcq(), Scalar::set_etat_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Matrice::set_lu(), Param_elliptic::set_poisson_vect_r(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Vector::set_vr_eta_mu(), Scalar::sol_elliptic_boundary(), Scalar::std_spectral_base(), Tbl::t, Tensor::triad, Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), Valeur::ylm(), and Valeur::ylm_i().
Vector Vector::poisson_neumann | ( | double | lambda, | |
const Valeur & | limit_vr, | |||
const Valeur & | limit_vt, | |||
const Valeur & | limit_vp, | |||
int | num_front | |||
) | const [inherited] |
Solves the vector Poisson equation with *this
as a source with a boundary condition on the excised sphere.
The equation solved is .
*this
must be given with dzpuis
= 4. It uses the Helmholtz decomposition (see documentation of p_potential
)
lambda | [input] ![]() | |
resu | [output] the solution ![]() |
Definition at line 687 of file vector_poisson_boundary.C.
References Tensor::mp, Vector::poisson_robin(), and Tensor::triad.
Vector Vector::poisson_robin | ( | double | lambda, | |
const Valeur & | limit_vr, | |||
const Valeur & | limit_vt, | |||
const Valeur & | limit_vp, | |||
double | fact_dir, | |||
double | fact_neu, | |||
int | num_front | |||
) | const [inherited] |
Solves the vector Poisson equation with *this
as a source with a boundary condition on the excised sphere.
The equation solved is .
*this
must be given with dzpuis
= 4. It uses the Helmholtz decomposition (see documentation of p_potential
)
lambda | [input] ![]() | |
resu | [output] the solution ![]() |
Definition at line 698 of file vector_poisson_boundary.C.
References Scalar::annule_hard(), Valeur::base, Valeur::c_cf, Tensor::cmp, Valeur::coef(), Scalar::div_tant(), Scalar::dsdt(), Valeur::get_etat(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_nzone(), Tensor::mp, norme(), Scalar::poisson_angu(), Vector::poisson_boundary(), Tensor::set(), Valeur::set_base(), Tensor::set_etat_zero(), Scalar::set_grid_point(), Scalar::set_spectral_va(), Scalar::stdsdp(), Tensor::triad, Scalar::val_grid_point(), and Valeur::ylm().
virtual int Vector::position | ( | const Itbl & | idx | ) | const [inline, virtual, inherited] |
Returns the position in the Scalar
array cmp
of a component given by its index.
Scalar
array cmp
corresponding to the index given in idx
. idx
must be a 1-D Itbl
of size 1, the element of which must be one of the spatial indices 1, 2 or 3. Reimplemented from Tensor.
Definition at line 388 of file vector.h.
References Itbl::get_dim(), and Itbl::get_ndim().
Returns the potential in the Helmholtz decomposition.
It first makes the Helmholtz decomposition (see documentation of p_potential
) of this
with respect to a given Metric
and then returns . Only in the case of contravariant vectors.
Definition at line 488 of file vector.C.
References Vector::decompose_div(), Tensor::get_place_met(), Vector::p_potential, and Tensor::set_dependance().
void Vector::pseudo_spectral_base | ( | ) | [virtual, inherited] |
Sets the standard spectal bases of decomposition for each component for a pseudo_vector.
Definition at line 342 of file vector.C.
References Tensor::cmp, Map::get_bvect_cart(), Map::get_bvect_spher(), Map::get_mg(), Base_vect::identify(), Tensor::mp, Mg3d::pseudo_base_vect_cart(), Mg3d::pseudo_base_vect_spher(), Scalar::set_spectral_base(), and Tensor::triad.
void Tensor::sauve | ( | FILE * | fd | ) | const [virtual, inherited] |
Save in a binary file.
Reimplemented in Scalar, and Tensor_sym.
Definition at line 902 of file tensor.C.
References Tensor::cmp, fwrite_be(), Tensor::n_comp, Base_vect::sauve(), Itbl::sauve(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
Scalar & Tensor::set | ( | int | i1, | |
int | i2, | |||
int | i3, | |||
int | i4 | |||
) | [inherited] |
Returns the value of a component for a tensor of valence 4 (read/write version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) | |
i3 | value of the third index (1, 2 or 3) | |
i4 | value of the fourth index (1, 2 or 3) |
(i1,i2,i3,i4) Definition at line 633 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::position(), Itbl::set(), and Tensor::valence.
Scalar & Tensor::set | ( | int | i1, | |
int | i2, | |||
int | i3 | |||
) | [inherited] |
Returns the value of a component for a tensor of valence 3 (read/write version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) | |
i3 | value of the third index (1, 2 or 3) |
(i1,i2,i3) Definition at line 617 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::position(), Itbl::set(), and Tensor::valence.
Scalar & Tensor::set | ( | int | i1, | |
int | i2 | |||
) | [inherited] |
Returns the value of a component for a tensor of valence 2 (read/write version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) |
(i1,i2) Definition at line 602 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::position(), Itbl::set(), and Tensor::valence.
Returns the value of a component (read/write version).
ind | 1-D Itbl of size valence containing the values of each index specifing the component, with the following storage convention:
|
ind
Definition at line 650 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Itbl::get_dim(), Itbl::get_ndim(), Tensor::position(), and Tensor::valence.
Scalar & Vector::set | ( | int | index | ) | [inherited] |
Read/write access to a component.
Definition at line 292 of file vector.C.
References Tensor::cmp, and Vector::del_deriv().
Defines the components through potentials and
.
(see members p_A
and p_mu
),
A_i | [input] Angular potential ![]() | |
mu_i | [input] Angular potential ![]() |
Definition at line 119 of file vector_divfree_aux.C.
References Tensor::cmp, del_deriv(), Tensor::get_mp(), Tensor::mp, Vector::p_A, Vector::p_eta, Vector::p_mu, Scalar::set_dzpuis(), sol_Dirac_A(), Tensor::triad, and Vector::update_vtvp().
void Tensor::set_dependance | ( | const Metric & | met | ) | const [protected, inherited] |
To be used to describe the fact that the derivatives members have been calculated with met
.
First it sets a null element of met_depend
to &met
and puts this
in the list of the dependancies of met
.
Definition at line 449 of file tensor.C.
References Tensor::met_depend, and Metric::tensor_depend.
void Vector_divfree::set_der_0x0 | ( | ) | const [protected] |
Sets the pointers on derived quantities to 0x0.
Reimplemented from Vector.
Definition at line 131 of file vector_divfree.C.
void Vector::set_der_met_0x0 | ( | int | i | ) | const [protected, inherited] |
Sets all the i-th components of met_depend
in the class Vector
(p_potential
, etc.
..) to 0x0.
Reimplemented from Tensor.
Definition at line 254 of file vector.C.
References Vector::p_div_free, and Vector::p_potential.
void Tensor::set_etat_nondef | ( | ) | [virtual, inherited] |
Sets the logical state of all components to ETATNONDEF
(undefined state).
Reimplemented in Scalar.
Definition at line 485 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::n_comp, and Scalar::set_etat_nondef().
void Tensor::set_etat_qcq | ( | ) | [virtual, inherited] |
Sets the logical state of all components to ETATQCQ
(ordinary state).
Reimplemented in Scalar.
Definition at line 477 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::n_comp, and Scalar::set_etat_qcq().
void Tensor::set_etat_zero | ( | ) | [virtual, inherited] |
Sets the logical state of all components to ETATZERO
(zero state).
Reimplemented in Scalar.
Definition at line 493 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::n_comp, and Scalar::set_etat_zero().
Itbl& Tensor::set_index_type | ( | ) | [inline, inherited] |
Sets the types of all the indices.
Itbl
) of size valence
that can be modified (COV
for a covariant one and CON
for a contravariant one) Definition at line 914 of file tensor.h.
References Tensor::type_indice.
int& Tensor::set_index_type | ( | int | i | ) | [inline, inherited] |
Sets the type of the index number i
.
i
must be strictly lower than valence
and obey the following convention:
i
= 0 : first index i
= 1 : second index COV
for a covariant index, CON
for a contravariant one) Definition at line 905 of file tensor.h.
References Itbl::set(), and Tensor::type_indice.
void Tensor::set_triad | ( | const Base_vect & | new_triad | ) | [inherited] |
Assigns a new vectorial basis (triad) of decomposition.
NB: this function modifies only the member triad
and leave unchanged the components (member cmp
). In order to change them coherently with the new basis, the function change_triad
(const Base_vect& ) must be called instead.
Definition at line 515 of file tensor.C.
References Tensor::triad.
void Vector_divfree::set_vr_eta_mu | ( | const Scalar & | vr_i, | |
const Scalar & | eta_i, | |||
const Scalar & | mu_i | |||
) |
Defines the components through ,
and
.
(see members p_eta
and p_mu
),
vr_i | [input] r-component of the vector | |
eta_i | [input] Angular potential ![]() | |
mu_i | [input] Angular potential ![]() |
Reimplemented from Vector.
Definition at line 95 of file vector_divfree_aux.C.
References Tensor::cmp, Tensor::get_mp(), Vector::p_eta, Vector::p_mu, Tensor::triad, and Vector::update_vtvp().
Sets the angular potentials (see member
p_mu
), and the component of the vector.
The potential is then deduced from
by the divergence-free condition. The components
and
are updated consistently by a call to the method
update_vtvp()
.
vr_i | [input] component ![]() | |
mu_i | [input] angular potential ![]() |
Definition at line 170 of file vector_divfree.C.
References Tensor::cmp, del_deriv(), eta(), Vector::p_mu, Tensor::triad, and Vector::update_vtvp().
void Vector_divfree::sol_Dirac_A | ( | const Scalar & | aaa, | |
Scalar & | eta, | |||
Scalar & | vr, | |||
const Param * | par_bc = 0x0 | |||
) | const [protected] |
Solves a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.
The system reads :
Definition at line 73 of file vector_divfree_A.C.
References Tensor::annule_domain(), Matrice::annule_hard(), Tbl::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::c, Valeur::c_cf, Mtbl_cf::dsdx(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Param::get_int(), Diff_id::get_matrice(), Diff_xdsdx::get_matrice(), Diff_sx::get_matrice(), Diff_dsdx::get_matrice(), Map::get_mg(), Param::get_n_int(), Param::get_n_tbl_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Param::get_tbl_mod(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Base_val::mult_x(), R_CHEBP, Mtbl_cf::set(), Tbl::set(), Matrice::set(), Valeur::set_etat_cf_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), Valeur::ylm(), and Valeur::ylm_i().
void Vector_divfree::sol_Dirac_A_1z | ( | const Scalar & | aaa, | |
Scalar & | eta, | |||
Scalar & | vr, | |||
const Param * | par_bc = 0x0 | |||
) | const [protected] |
Solves a one-domain system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.
The system reads :
Definition at line 65 of file vector_divfree_A_1z.C.
References Tensor::annule_domain(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::c, Valeur::c_cf, Map_af::get_alpha(), Scalar::get_etat(), Param::get_int(), Diff_sx::get_matrice(), Diff_dsdx::get_matrice(), Map::get_mg(), Param::get_n_int(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Base_val::mult_x(), R_CHEBP, Mtbl_cf::set(), Tbl::set(), Matrice::set(), Valeur::set_etat_cf_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Valeur::ylm(), and Valeur::ylm_i().
void Vector_divfree::sol_Dirac_A_poisson | ( | const Scalar & | aaa, | |
Scalar & | eta, | |||
Scalar & | vr, | |||
const Param * | par_bc = 0x0 | |||
) | const [protected] |
Solves via a poisson method a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.
The system reads :
Definition at line 64 of file vector_divfree_A_poisson.C.
References Scalar::div_r(), Scalar::dsdr(), Scalar::lapang(), Scalar::mult_r_dzpuis(), and Scalar::poisson_tau().
void Vector_divfree::sol_Dirac_A_tau | ( | const Scalar & | aaa, | |
Scalar & | eta, | |||
Scalar & | vr, | |||
const Param * | par_bc = 0x0 | |||
) | const [protected] |
Solves via a tau method a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.
The system reads :
Definition at line 69 of file vector_divfree_A_tau.C.
References Tensor::annule_domain(), Tbl::annule_hard(), Matrice::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::c, Valeur::c_cf, Scalar::get_dzpuis(), Scalar::get_etat(), Param::get_int(), Diff_xdsdx::get_matrice(), Diff_id::get_matrice(), Diff_sx::get_matrice(), Diff_dsdx::get_matrice(), Map::get_mg(), Tensor::get_mp(), Param::get_n_int(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Base_val::mult_x(), R_CHEB, R_JACO02, Mtbl_cf::set(), Tbl::set(), Matrice::set(), Matrice::set_band(), Valeur::set_etat_cf_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Tbl::t, Valeur::ylm(), and Valeur::ylm_i().
void Tensor::spectral_display | ( | const char * | comment = 0x0 , |
|
double | threshold = 1.e-7 , |
|||
int | precision = 4 , |
|||
ostream & | ostr = cout | |||
) | const [virtual, inherited] |
Displays the spectral coefficients and the associated basis functions of each component.
This function shows only the values greater than a given threshold.
comment | comment to be printed at top of the display (default: 0x0 = nothing printed) | |
threshold | [input] Value above which a coefficient is printed (default: 1.e-7) | |
precision | [input] Number of printed digits (default: 4) | |
ostr | [input] Output stream used for the printing (default: cout) |
Reimplemented in Scalar.
Definition at line 870 of file tensor.C.
References Tensor::cmp, Tensor::indices(), Tensor::n_comp, Scalar::spectral_display(), and Tensor::valence.
void Vector::std_spectral_base | ( | ) | [virtual, inherited] |
Sets the standard spectal bases of decomposition for each component.
Reimplemented from Tensor.
Definition at line 312 of file vector.C.
References Tensor::cmp, Map::get_bvect_cart(), Map::get_bvect_spher(), Map::get_mg(), Base_vect::identify(), Tensor::mp, Scalar::set_spectral_base(), Mg3d::std_base_vect_cart(), Mg3d::std_base_vect_spher(), and Tensor::triad.
void Tensor::std_spectral_base_odd | ( | ) | [virtual, inherited] |
Sets the standard odd spectal bases of decomposition for each component.
Currently only implemented for a scalar.
Reimplemented in Scalar.
Definition at line 978 of file tensor.C.
References Tensor::cmp, Scalar::std_spectral_base_odd(), and Tensor::valence.
Trace with respect to a given metric for a valence 2 tensor.
gam | metric used to raise or lower one of the indices, in order to take the trace |
Definition at line 193 of file tensor_calculus.C.
References Metric::con(), contract(), Metric::cov(), Tensor::trace(), Tensor::type_indice, and Tensor::valence.
Scalar Tensor::trace | ( | ) | const [inherited] |
Trace on two different type indices for a valence 2 tensor.
Definition at line 176 of file tensor_calculus.C.
References Tensor::mp, Tensor::operator()(), Scalar::set_etat_zero(), Tensor::type_indice, and Tensor::valence.
Trace with respect to a given metric.
ind1 | first index for the contraction, with the following convention :
| |
ind2 | second index for the contraction | |
gam | metric used to raise or lower ind1 in order that it has a opposite type than ind2 |
Definition at line 149 of file tensor_calculus.C.
References Metric::con(), contract(), Metric::cov(), Tensor::trace(), Tensor::type_indice, and Tensor::valence.
Tensor Tensor::trace | ( | int | ind1, | |
int | ind2 | |||
) | const [inherited] |
Trace on two different type indices.
ind1 | first index for the contraction, with the following convention :
| |
ind2 | second index for the contraction |
Definition at line 90 of file tensor_calculus.C.
References Tensor::cmp, Tensor::get_n_comp(), Tensor::indices(), Tensor::mp, Tensor::position(), Tensor::set(), Itbl::set(), Scalar::set_etat_zero(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
Computes a new tensor by raising an index of *this
.
ind | index to be raised, with the following convention :
| |
gam | metric used to raise the index (contraction with the twice contravariant form of the metric on the index ind ). |
Definition at line 221 of file tensor_calculus.C.
References Metric::con(), contract(), Tensor::indices(), Tensor::mp, Tensor::n_comp, Tensor::set(), Itbl::set(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
Computes a new tensor by raising or lowering all the indices of *this
.
gam | metric used to lower the contravariant indices and raising the covariant ones. |
Definition at line 301 of file tensor_calculus.C.
References Tensor::down(), Tensor::Tensor(), Tensor::type_indice, Tensor::up(), and Tensor::valence.
void Vector_divfree::update_etavr | ( | ) |
Computes the components and
from the potential A and the divergence-free condition, according to :
.
Definition at line 72 of file vector_divfree_aux.C.
References Tensor::cmp, del_deriv(), Tensor::mp, Vector::p_A, Vector::p_eta, and sol_Dirac_A().
void Vector::update_vtvp | ( | ) | [inherited] |
Computes the components and
from the potential
and
, according to:
.
Definition at line 163 of file vector_etamu.C.
References Tensor::cmp, Vector::del_deriv(), Scalar::dsdt(), Vector::p_eta, Vector::p_mu, and Scalar::stdsdp().
void Vector::visu_arrows | ( | double | xmin, | |
double | xmax, | |||
double | ymin, | |||
double | ymax, | |||
double | zmin, | |||
double | zmax, | |||
const char * | title0 = 0x0 , |
|||
const char * | filename0 = 0x0 , |
|||
bool | start_dx = true , |
|||
int | nx = 8 , |
|||
int | ny = 8 , |
|||
int | nz = 8 | |||
) | const [inherited] |
3D visualization via OpenDX.
xmin | [input] defines with xmax the x range of the visualization box | |
xmax | [input] defines with xmin the x range of the visualization box | |
ymin | [input] defines with ymax the y range of the visualization box | |
ymax | [input] defines with ymin the y range of the visualization box | |
zmin | [input] defines with zmax the z range of the visualization box | |
zmax | [input] defines with zmin the z range of the visualization box | |
title | [input] title for the graph (for OpenDX legend) | |
filename | [input] name for the file which will be the input for OpenDX; the default 0x0 is transformed into "vector_arrows" | |
start_dx | [input] determines whether OpenDX must be launched (as a subprocess) to view the field; if set to false , only input files for future usage of OpenDX are created | |
nx | [input] number of points in the x direction (uniform sampling) | |
ny | [input] number of points in the y direction (uniform sampling) | |
nz | [input] number of points in the z direction (uniform sampling) |
Definition at line 58 of file vector_visu.C.
References Valeur::c_cf, Vector::change_triad(), Scalar::check_dzpuis(), Valeur::coef(), Map::convert_absolute(), Scalar::dec_dzpuis(), Map::get_bvect_cart(), Map::get_bvect_spher(), Scalar::get_dzpuis(), Tensor::mp, Vector::operator()(), Vector::set(), Tensor::triad, Map::val_lx(), Mtbl_cf::val_point(), and Vector::Vector().
Tensor_sym operator* | ( | const Tensor_sym & | a, | |
const Tensor_sym & | b | |||
) | [friend, inherited] |
Tensorial product of two symmetric tensors.
NB: the output is an object of class Tensor_sym
, with the two symmetric indices corresponding to the symmetric indices of tensor a
. This means that the symmetries of tensor b
indices are not used in the storage, since there is currently no class in Lorene to manage tensors with more than two symmetric indices.
Definition at line 147 of file tensor_sym_calculus.C.
Scalar** Tensor::cmp [protected, inherited] |
const Metric* Tensor::met_depend[N_MET_MAX] [mutable, protected, inherited] |
const Metric* const Vector_divfree::met_div [protected] |
const Map* const Tensor::mp [protected, inherited] |
int Tensor::n_comp [protected, inherited] |
Scalar* Vector::p_A [mutable, protected, inherited] |
Tensor* Tensor::p_derive_con[N_MET_MAX] [mutable, protected, inherited] |
Array of pointers on the contravariant derivatives of this
with respect to various metrics.
See the comments of met_depend
. See also the comments of method derive_con()
for a precise definition of a "contravariant" derivative.
Tensor* Tensor::p_derive_cov[N_MET_MAX] [mutable, protected, inherited] |
Array of pointers on the covariant derivatives of this
with respect to various metrics.
See the comments of met_depend
. See also the comments of method derive_cov()
for the index convention of the covariant derivation.
Vector_divfree* Vector::p_div_free[N_MET_MAX] [mutable, protected, inherited] |
Tensor* Tensor::p_divergence[N_MET_MAX] [mutable, protected, inherited] |
Array of pointers on the divergence of this
with respect to various metrics.
See the comments of met_depend
. See also the comments of method divergence()
for a precise definition of a the divergence with respect to a given metric.
Scalar* Vector::p_eta [mutable, protected, inherited] |
Scalar* Vector::p_mu [mutable, protected, inherited] |
Scalar* Vector::p_potential[N_MET_MAX] [mutable, protected, inherited] |
const Base_vect* Tensor::triad [protected, inherited] |
Itbl Tensor::type_indice [protected, inherited] |
int Tensor::valence [protected, inherited] |