LORENE
Lorene::Tensor Class Reference

Tensor handling. More...

#include <tensor.h>

Inheritance diagram for Lorene::Tensor:
Lorene::Scalar Lorene::Tensor_sym Lorene::Vector Lorene::Sym_tensor Lorene::Vector_divfree Lorene::Sym_tensor_trans Lorene::Sym_tensor_tt

Public Member Functions

 Tensor (const Map &map, int val, const Itbl &tipe, const Base_vect &triad_i)
 Standard constructor. More...
 
 Tensor (const Map &map, int val, const Itbl &tipe, const Base_vect *triad_i)
 Standard constructor with the triad passed as a pointer. More...
 
 Tensor (const Map &map, int val, int tipe, const Base_vect &triad_i)
 Standard constructor when all the indices are of the same type. More...
 
 Tensor (const Tensor &)
 Copy constructor. More...
 
 Tensor (const Map &map, const Base_vect &triad_i, FILE *fich)
 Constructor from a file (see sauve(FILE*) ). More...
 
virtual ~Tensor ()
 Destructor. More...
 
virtual void set_etat_nondef ()
 Sets the logical state of all components to ETATNONDEF
(undefined state). More...
 
virtual void set_etat_zero ()
 Sets the logical state of all components to ETATZERO
(zero state). More...
 
virtual void set_etat_qcq ()
 Sets the logical state of all components to ETATQCQ
(ordinary state). More...
 
virtual void allocate_all ()
 Performs the memory allocation of all the elements, down to the double arrays of the Tbl s. More...
 
virtual void change_triad (const Base_vect &new_triad)
 Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly. More...
 
void set_triad (const Base_vect &new_triad)
 Assigns a new vectorial basis (triad) of decomposition. More...
 
virtual void operator= (const Tensor &)
 Assignment to a Tensor. More...
 
Scalarset (const Itbl &ind)
 Returns the value of a component (read/write version). More...
 
Scalarset (int i1, int i2)
 Returns the value of a component for a tensor of valence 2 (read/write version). More...
 
Scalarset (int i1, int i2, int i3)
 Returns the value of a component for a tensor of valence 3 (read/write version). More...
 
Scalarset (int i1, int i2, int i3, int i4)
 Returns the value of a component for a tensor of valence 4 (read/write version). More...
 
void annule_domain (int l)
 Sets the Tensor to zero in a given domain. More...
 
virtual void annule (int l_min, int l_max)
 Sets the Tensor to zero in several domains. More...
 
void annule_extern_cn (int l_0, int deg)
 Performs a smooth (C^n) transition in a given domain to zero. More...
 
virtual void std_spectral_base ()
 Sets the standard spectal bases of decomposition for each component. More...
 
virtual void std_spectral_base_odd ()
 Sets the standard odd spectal bases of decomposition for each component. More...
 
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). More...
 
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). More...
 
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 ). More...
 
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 ). More...
 
virtual void exponential_filter_ylm_phi (int lzmin, int lzmax, int p_r, int p_tet, int p_phi, double alpha=-16.)
 Applies exponential filters to all components (see Scalar::exponential_filter_ylm_phi ). More...
 
const Tensorderive_cov (const Metric &gam) const
 Returns the covariant derivative of this with respect to some metric $\gamma$. More...
 
const Tensorderive_con (const Metric &gam) const
 Returns the "contravariant" derivative of this with respect to some metric $\gamma$, by raising the last index of the covariant derivative (cf. More...
 
const Tensordivergence (const Metric &gam) const
 Computes the divergence of this with respect to some metric $\gamma$. More...
 
Tensor derive_lie (const Vector &v) const
 Computes the Lie derivative of this with respect to some vector field v. More...
 
Tensor up (int ind, const Metric &gam) const
 Computes a new tensor by raising an index of *this. More...
 
Tensor down (int ind, const Metric &gam) const
 Computes a new tensor by lowering an index of *this. More...
 
Tensor up_down (const Metric &gam) const
 Computes a new tensor by raising or lowering all the indices of *this . More...
 
Tensor trace (int ind1, int ind2) const
 Trace on two different type indices. More...
 
Tensor trace (int ind1, int ind2, const Metric &gam) const
 Trace with respect to a given metric. More...
 
Scalar trace () const
 Trace on two different type indices for a valence 2 tensor. More...
 
Scalar trace (const Metric &gam) const
 Trace with respect to a given metric for a valence 2 tensor. More...
 
virtual int position (const Itbl &ind) const
 Returns the position in the array cmp of a component given by its indices. More...
 
virtual Itbl indices (int pos) const
 Returns the indices of a component given by its position in the array cmp . More...
 
const Mapget_mp () const
 Returns the mapping. More...
 
const Base_vectget_triad () const
 Returns the vectorial basis (triad) on which the components are defined. More...
 
int get_valence () const
 Returns the valence. More...
 
int get_n_comp () const
 Returns the number of stored components. More...
 
int get_index_type (int i) const
 Gives the type (covariant or contravariant) of the index number i . More...
 
Itbl get_index_type () const
 Returns the types of all the indices. More...
 
int & set_index_type (int i)
 Sets the type of the index number i . More...
 
Itblset_index_type ()
 Sets the types of all the indices. More...
 
const Scalaroperator() (const Itbl &ind) const
 Returns the value of a component (read-only version). More...
 
const Scalaroperator() (int i1, int i2) const
 Returns the value of a component for a tensor of valence 2 (read-only version). More...
 
const Scalaroperator() (int i1, int i2, int i3) const
 Returns the value of a component for a tensor of valence 3 (read-only version). More...
 
const Scalaroperator() (int i1, int i2, int i3, int i4) const
 Returns the value of a component for a tensor of valence 4 (read-only version). More...
 
void operator+= (const Tensor &)
 += Tensor More...
 
void operator-= (const Tensor &)
 -= Tensor More...
 
virtual void sauve (FILE *) const
 Save in a binary file. More...
 
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. More...
 

Protected Member Functions

 Tensor (const Map &map)
 Constructor for a scalar field: to be used only by the derived class Scalar . More...
 
 Tensor (const Map &map, int val, const Itbl &tipe, int n_comp_i, const Base_vect &triad_i)
 Constructor to be used by derived classes, with symmetries among the components. More...
 
 Tensor (const Map &map, int val, int tipe, int n_comp_i, const Base_vect &triad_i)
 Constructor used by derived classes, with symmetries among the components, when all the indices are of the same type. More...
 
virtual void del_deriv () const
 Deletes the derived quantities. More...
 
void set_der_0x0 () const
 Sets the pointers on derived quantities to 0x0. More...
 
virtual void del_derive_met (int) const
 Logical destructor of the derivatives depending on the i-th element of met_depend . More...
 
void set_der_met_0x0 (int) const
 Sets all the i-th components of met_depend , p_derive_cov , etc... More...
 
void set_dependance (const Metric &) const
 To be used to describe the fact that the derivatives members have been calculated with met . More...
 
int get_place_met (const Metric &) const
 Returns the position of the pointer on metre in the array met_depend . More...
 
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 ). More...
 

Protected Attributes

const Map *const mp
 Mapping on which the numerical values at the grid points are defined. More...
 
int valence
 Valence of the tensor (0 = scalar, 1 = vector, etc...) More...
 
const Base_vecttriad
 Vectorial basis (triad) with respect to which the tensor components are defined. More...
 
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. More...
 
int n_comp
 Number of stored components, depending on the symmetry. More...
 
Scalar ** cmp
 Array of size n_comp of pointers onto the components. More...
 
const Metricmet_depend [N_MET_MAX]
 Array on the Metric 's which were used to compute derived quantities, like p_derive_cov , etc... More...
 
Tensorp_derive_cov [N_MET_MAX]
 Array of pointers on the covariant derivatives of this with respect to various metrics. More...
 
Tensorp_derive_con [N_MET_MAX]
 Array of pointers on the contravariant derivatives of this with respect to various metrics. More...
 
Tensorp_divergence [N_MET_MAX]
 Array of pointers on the divergence of this with respect to various metrics. More...
 

Friends

class Scalar
 
class Vector
 
class Sym_tensor
 
class Tensor_sym
 
class Metric
 
ostream & operator<< (ostream &, const Tensor &)
 
Scalar operator+ (const Tensor &, const Scalar &)
 Tensor + Scalar. The Tensor must be of valence 0. More...
 
Scalar operator+ (const Scalar &, const Tensor &)
 Scalar + Tensor. The Tensor must be of valence 0. More...
 
Scalar operator- (const Tensor &, const Scalar &)
 Tensor - Scalar. The Tensor must be of valence 0. More...
 
Scalar operator- (const Scalar &, const Tensor &)
 Scalar - Tensor. The Tensor must be of valence 0. More...
 
Tensor operator* (const Tensor &, const Tensor &)
 Tensorial product. More...
 
Tensor_sym operator* (const Tensor &, const Tensor_sym &)
 Tensorial product with symmetries. More...
 
Tensor_sym operator* (const Tensor_sym &, const Tensor &)
 Tensorial product with symmetries. More...
 
Tensor_sym operator* (const Tensor_sym &, const Tensor_sym &)
 Tensorial product of two symmetric tensors. More...
 

Detailed Description

Tensor handling.

()

This class has been devised to replace Tenseur and Cmp (the latter via the derived class Scalar ).

The Tensor class is intended to store the components of a tensorial field with respect to a specific basis (triad).
All this is 3D meaning that the indices go from 1 to 3.

Definition at line 294 of file tensor.h.

Constructor & Destructor Documentation

◆ Tensor() [1/8]

Lorene::Tensor::Tensor ( const Map map,
int  val,
const Itbl tipe,
const Base_vect triad_i 
)

Standard constructor.

Parameters
mapthe mapping
valvalence of the tensor
tipe1-D 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, with the following storage convention:
  • tipe(0) : type of the first index
  • tipe(1) : type of the second index
  • and so on...
triad_ivectorial basis (triad) with respect to which the tensor components are defined

Definition at line 211 of file tensor.C.

References cmp, Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), n_comp, set_der_0x0(), and valence.

◆ Tensor() [2/8]

Lorene::Tensor::Tensor ( const Map map,
int  val,
const Itbl tipe,
const Base_vect triad_i 
)

Standard constructor with the triad passed as a pointer.

Parameters
mapthe mapping
valvalence of the tensor
tipe1-D 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, with the following storage convention:
  • tipe(0) : type of the first index
  • tipe(1) : type of the second index
  • and so on...
triad_ipointer on the vectorial basis (triad) with respect to which the tensor components are defined (can be set to 0x0 for a scalar field)

Definition at line 234 of file tensor.C.

References cmp, Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), n_comp, set_der_0x0(), and valence.

◆ Tensor() [3/8]

Lorene::Tensor::Tensor ( const Map map,
int  val,
int  tipe,
const Base_vect triad_i 
)

Standard constructor when all the indices are of the same type.

Parameters
mapthe mapping
valvalence of the tensor
tipethe type (COV or CON ) of the indices.
triad_ivectorial basis (triad) with respect to which the tensor components are defined.

Definition at line 260 of file tensor.C.

References cmp, n_comp, set_der_0x0(), type_indice, and valence.

◆ Tensor() [4/8]

Lorene::Tensor::Tensor ( const Tensor source)

Copy constructor.

Definition at line 280 of file tensor.C.

References cmp, indices(), n_comp, position(), Lorene::pow(), set_der_0x0(), and valence.

◆ Tensor() [5/8]

Lorene::Tensor::Tensor ( const Map map,
const Base_vect triad_i,
FILE *  fich 
)

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

Parameters
mapthe mapping
triad_ivectorial 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.
fichfile which has been created by the function sauve(FILE*) .

Definition at line 304 of file tensor.C.

References Lorene::Base_vect::bvect_from_file(), cmp, Lorene::fread_be(), Lorene::Map::get_mg(), mp, n_comp, set_der_0x0(), triad, and valence.

◆ Tensor() [6/8]

Lorene::Tensor::Tensor ( const Map map)
explicitprotected

Constructor for a scalar field: to be used only by the derived class Scalar .

Definition at line 331 of file tensor.C.

References cmp, n_comp, and set_der_0x0().

◆ Tensor() [7/8]

Lorene::Tensor::Tensor ( const Map map,
int  val,
const Itbl tipe,
int  n_comp_i,
const Base_vect triad_i 
)
protected

Constructor to be used by derived classes, with symmetries among the components.

The number of independent components is given as an argument (n_comp_i ), and not computed from the valence, as in the standard constructor.

Parameters
mapthe mapping
valvalence of the tensor
tipe1-D 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, with the following storage convention:
  • tipe(0) : type of the first index
  • tipe(1) : type of the second index
  • and so on...
n_comp_inumber of components to be stored
triad_ivectorial basis (triad) with respect to which the tensor components are defined

Definition at line 343 of file tensor.C.

References cmp, Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), n_comp, set_der_0x0(), and valence.

◆ Tensor() [8/8]

Lorene::Tensor::Tensor ( const Map map,
int  val,
int  tipe,
int  n_comp_i,
const Base_vect triad_i 
)
protected

Constructor used by derived classes, with symmetries among the components, when all the indices are of the same type.

The number of independent components is given as a argument (n_comp_i ), and not computed from the valence, as in the standard constructor.

Parameters
mapthe mapping
valvalence of the tensor
tipethe type of the indices.
n_comp_inumber of components to be stored
triad_ivectorial basis (triad) with respect to which the tensor components are defined

Definition at line 368 of file tensor.C.

References cmp, n_comp, set_der_0x0(), type_indice, and valence.

◆ ~Tensor()

Lorene::Tensor::~Tensor ( )
virtual

Destructor.

Definition at line 394 of file tensor.C.

References cmp, del_deriv(), and n_comp.

Member Function Documentation

◆ allocate_all()

void Lorene::Tensor::allocate_all ( )
virtual

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 Lorene::Scalar.

Definition at line 517 of file tensor.C.

References Lorene::Scalar::allocate_all(), cmp, del_deriv(), and n_comp.

◆ annule()

void Lorene::Tensor::annule ( int  l_min,
int  l_max 
)
virtual

Sets the Tensor to zero in several domains.

Parameters
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 Lorene::Scalar.

Definition at line 680 of file tensor.C.

References Lorene::Scalar::annule(), cmp, del_deriv(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), mp, n_comp, and set_etat_zero().

◆ annule_domain()

void Lorene::Tensor::annule_domain ( int  l)

Sets the Tensor to zero in a given domain.

Parameters
l[input] Index of the domain in which the Tensor will be set (logically) to zero.

Definition at line 675 of file tensor.C.

References annule().

◆ annule_extern_cn()

void Lorene::Tensor::annule_extern_cn ( int  l_0,
int  deg 
)

Performs a smooth (C^n) transition in a given domain to zero.

Parameters
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 699 of file tensor.C.

References Lorene::Map::get_mg(), Lorene::Mg3d::get_type_r(), and mp.

◆ change_triad()

◆ compute_derive_lie()

void Lorene::Tensor::compute_derive_lie ( const Vector v,
Tensor resu 
) const
protected

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 342 of file tensor_calculus.C.

References cmp, Lorene::contract(), Lorene::Scalar::dec_dzpuis(), derive_cov(), Lorene::Map::flat_met_cart(), Lorene::Map::flat_met_spher(), Lorene::Scalar::get_dzpuis(), get_n_comp(), get_triad(), indices(), mp, n_comp, operator()(), Lorene::Itbl::set(), set(), triad, type_indice, and valence.

◆ dec_dzpuis()

void Lorene::Tensor::dec_dzpuis ( int  dec = 1)
virtual

Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified external domain (CED).

Reimplemented in Lorene::Scalar.

Definition at line 817 of file tensor.C.

References cmp, del_deriv(), and n_comp.

◆ del_deriv()

void Lorene::Tensor::del_deriv ( ) const
protectedvirtual

Deletes the derived quantities.

Reimplemented in Lorene::Sym_tensor_tt, Lorene::Vector_divfree, Lorene::Sym_tensor_trans, Lorene::Scalar, Lorene::Sym_tensor, and Lorene::Vector.

Definition at line 407 of file tensor.C.

References del_derive_met(), and set_der_0x0().

◆ del_derive_met()

void Lorene::Tensor::del_derive_met ( int  j) const
protectedvirtual

Logical destructor of the derivatives depending on the i-th element of met_depend .

Reimplemented in Lorene::Sym_tensor, and Lorene::Vector.

Definition at line 423 of file tensor.C.

References met_depend, p_derive_con, p_derive_cov, p_divergence, set_der_met_0x0(), and Lorene::Metric::tensor_depend.

◆ derive_con()

const Tensor & Lorene::Tensor::derive_con ( const Metric gam) const

Returns the "contravariant" derivative of this with respect to some metric $\gamma$, by raising the last index of the covariant derivative (cf.

method derive_cov() ) with $\gamma$.

Definition at line 1023 of file tensor.C.

References Lorene::Metric::con(), Lorene::contract(), derive_cov(), get_index_type(), get_place_met(), mp, p_derive_con, Lorene::Itbl::set(), set_dependance(), Lorene::Tensor_sym::sym_index1(), Lorene::Tensor_sym::sym_index2(), Tensor(), triad, and valence.

◆ derive_cov()

const Tensor & Lorene::Tensor::derive_cov ( const Metric gam) const

Returns the covariant derivative of this with respect to some metric $\gamma$.

$T$ denoting the tensor represented by this and $\nabla T$ its covariant derivative with respect to the metric $\gamma$, the extra index (with respect to the indices of $T$) of $\nabla T$ is chosen to be the last one. This convention agrees with that of MTW (see Eq. (10.17) of MTW). For instance, if $T$ is a 1-form, whose components w.r.t. the triad $e^i$ are $T_i$: $T=T_i \; e^i$, then the covariant derivative of $T$ is the bilinear form $\nabla T$ whose components $\nabla_j T_i$ are such that

\[ \nabla T = \nabla_j T_i \; e^i \otimes e^j \]

Parameters
gammetric $\gamma$
Returns
covariant derivative $\nabla T$ of this with respect to the connection $\nabla$ associated with the metric $\gamma$

Definition at line 1011 of file tensor.C.

References Lorene::Metric::connect(), get_place_met(), Lorene::Connection::p_derive_cov(), p_derive_cov, and set_dependance().

◆ derive_lie()

Tensor Lorene::Tensor::derive_lie ( const Vector v) const

Computes the Lie derivative of this with respect to some vector field v.

Definition at line 484 of file tensor_calculus.C.

References compute_derive_lie(), mp, triad, type_indice, and valence.

◆ divergence()

const Tensor & Lorene::Tensor::divergence ( const Metric gam) const

Computes the divergence of this with respect to some metric $\gamma$.

The divergence is taken with respect of the last index of this which thus must be contravariant. For instance if the tensor $T$ represented by this is a twice contravariant tensor, whose components w.r.t. the triad $e_i$ are $T^{ij}$: $T = T^{ij} \; e_i \otimes e_j$, the divergence of $T$ w.r.t. $\gamma$ is the vector

\[ {\rm div}\, T = \nabla_k T^{ik} \; e_i \]

where $\nabla$ denotes the connection associated with the metric $\gamma$.

Parameters
gammetric $\gamma$
Returns
divergence of this with respect to $\gamma$.

Definition at line 1064 of file tensor.C.

References Lorene::Metric::connect(), get_place_met(), Lorene::Connection::p_divergence(), p_divergence, and set_dependance().

◆ down()

Tensor Lorene::Tensor::down ( int  ind,
const Metric gam 
) const

Computes a new tensor by lowering an index of *this.

Parameters
indindex to be lowered, with the following convention :
  • ind1 = 0 : first index of the tensor
  • ind1 = 1 : second index of the tensor
  • and so on... (ind must be of covariant type (CON )).
gammetric used to lower the index (contraction with the twice covariant form of the metric on the index ind ).

Definition at line 268 of file tensor_calculus.C.

References Lorene::contract(), Lorene::Metric::cov(), indices(), mp, n_comp, Lorene::Itbl::set(), set(), triad, type_indice, and valence.

◆ exponential_filter_r()

void Lorene::Tensor::exponential_filter_r ( int  lzmin,
int  lzmax,
int  p,
double  alpha = -16. 
)
virtual

Applies exponential filters to all components (see Scalar::exponential_filter_r ).

Works only for Cartesian components.

Reimplemented in Lorene::Scalar, Lorene::Sym_tensor, and Lorene::Vector.

Definition at line 1075 of file tensor.C.

References cmp, Lorene::Map::get_bvect_cart(), Lorene::Base_vect::identify(), mp, n_comp, and triad.

◆ exponential_filter_ylm()

void Lorene::Tensor::exponential_filter_ylm ( int  lzmin,
int  lzmax,
int  p,
double  alpha = -16. 
)
virtual

Applies exponential filters to all components (see Scalar::exponential_filter_ylm ).

Works only for Cartesian components.

Reimplemented in Lorene::Scalar, Lorene::Sym_tensor, and Lorene::Vector.

Definition at line 1088 of file tensor.C.

References cmp, Lorene::Map::get_bvect_cart(), Lorene::Base_vect::identify(), mp, n_comp, and triad.

◆ exponential_filter_ylm_phi()

void Lorene::Tensor::exponential_filter_ylm_phi ( int  lzmin,
int  lzmax,
int  p_r,
int  p_tet,
int  p_phi,
double  alpha = -16. 
)
virtual

Applies exponential filters to all components (see Scalar::exponential_filter_ylm_phi ).

Works only for Cartesian components.

Reimplemented in Lorene::Scalar.

Definition at line 1101 of file tensor.C.

References cmp, Lorene::Map::get_bvect_cart(), Lorene::Base_vect::identify(), mp, n_comp, and triad.

◆ get_index_type() [1/2]

int Lorene::Tensor::get_index_type ( int  i) const
inline

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
  • and so on...
Returns
COV for a covariant index, CON for a contravariant one.

Definition at line 899 of file tensor.h.

References type_indice.

◆ get_index_type() [2/2]

Itbl Lorene::Tensor::get_index_type ( ) const
inline

Returns the types of all the indices.

Returns
1-D 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.

Definition at line 909 of file tensor.h.

References type_indice.

◆ get_mp()

const Map& Lorene::Tensor::get_mp ( ) const
inline

Returns the mapping.

Definition at line 874 of file tensor.h.

References mp.

◆ get_n_comp()

int Lorene::Tensor::get_n_comp ( ) const
inline

Returns the number of stored components.

Definition at line 885 of file tensor.h.

References n_comp.

◆ get_place_met()

int Lorene::Tensor::get_place_met ( const Metric metre) const
protected

Returns the position of the pointer on metre in the array met_depend .

Definition at line 452 of file tensor.C.

References met_depend.

◆ get_triad()

const Base_vect* Lorene::Tensor::get_triad ( ) const
inline

Returns the vectorial basis (triad) on which the components are defined.

Definition at line 879 of file tensor.h.

References triad.

◆ get_valence()

int Lorene::Tensor::get_valence ( ) const
inline

Returns the valence.

Definition at line 882 of file tensor.h.

References valence.

◆ inc_dzpuis()

void Lorene::Tensor::inc_dzpuis ( int  inc = 1)
virtual

Increases by inc units the value of dzpuis and changes accordingly the values in the compactified external domain (CED).

Reimplemented in Lorene::Scalar.

Definition at line 825 of file tensor.C.

References cmp, del_deriv(), and n_comp.

◆ indices()

Itbl Lorene::Tensor::indices ( int  pos) const
virtual

Returns the indices of a component given by its position in the array cmp .

Parameters
pos[input] position in the array cmp of the pointer to the Scalar representing a component
Returns
1-D array of integers (class Itbl ) of size valence giving the value of each index for the component located at the position pos in the array cmp, with the following storage convention:
  • Itbl(0) : value of the first index (1, 2 or 3)
  • Itbl(1) : value of the second index (1, 2 or 3)
  • and so on...

Reimplemented in Lorene::Tensor_sym, and Lorene::Vector.

Definition at line 548 of file tensor.C.

References n_comp, Lorene::Itbl::set(), and valence.

◆ operator()() [1/4]

const Scalar & Lorene::Tensor::operator() ( const Itbl ind) const

Returns the value of a component (read-only version).

Parameters
ind1-D Itbl of size valence containing the values of each index specifing the component, with the following storage convention:
  • ind(0) : value of the first index (1, 2 or 3)
  • ind(1) : value of the second index (1, 2 or 3)
  • and so on...
Returns
reference on the component specified by ind

Definition at line 807 of file tensor.C.

References cmp, Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), position(), and valence.

◆ operator()() [2/4]

const Scalar & Lorene::Tensor::operator() ( int  i1,
int  i2 
) const

Returns the value of a component for a tensor of valence 2 (read-only version).

Parameters
i1value of the first index (1, 2 or 3)
i2value of the second index (1, 2 or 3)
Returns
reference on the component specified by (i1,i2)

Definition at line 769 of file tensor.C.

References cmp, position(), Lorene::Itbl::set(), and valence.

◆ operator()() [3/4]

const Scalar & Lorene::Tensor::operator() ( int  i1,
int  i2,
int  i3 
) const

Returns the value of a component for a tensor of valence 3 (read-only version).

Parameters
i1value of the first index (1, 2 or 3)
i2value of the second index (1, 2 or 3)
i3value of the third index (1, 2 or 3)
Returns
reference on the component specified by (i1,i2,i3)

Definition at line 780 of file tensor.C.

References cmp, position(), Lorene::Itbl::set(), and valence.

◆ operator()() [4/4]

const Scalar & Lorene::Tensor::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).

Parameters
i1value of the first index (1, 2 or 3)
i2value of the second index (1, 2 or 3)
i3value of the third index (1, 2 or 3)
i4value of the fourth index (1, 2 or 3)
Returns
reference on the component specified by (i1,i2,i3,i4)

Definition at line 792 of file tensor.C.

References cmp, position(), Lorene::Itbl::set(), and valence.

◆ operator+=()

void Lorene::Tensor::operator+= ( const Tensor t)

+= Tensor

Definition at line 580 of file tensor.C.

References cmp, del_deriv(), indices(), n_comp, position(), triad, type_indice, and valence.

◆ operator-=()

void Lorene::Tensor::operator-= ( const Tensor t)

-= Tensor

Definition at line 596 of file tensor.C.

References cmp, del_deriv(), indices(), n_comp, position(), triad, type_indice, and valence.

◆ operator=()

void Lorene::Tensor::operator= ( const Tensor t)
virtual

◆ position()

int Lorene::Tensor::position ( const Itbl ind) const
virtual

Returns the position in the array cmp of a component given by its indices.

Parameters
ind[input] 1-D array of integers (class Itbl ) of size valence giving the values of each index specifing the component, with the following storage convention:
  • ind(0) : value of the first index (1, 2 or 3)
  • ind(1) : value of the second index (1, 2 or 3)
  • and so on...
Returns
position in the array cmp of the pointer to the Scalar containing the component specified by ind

Reimplemented in Lorene::Tensor_sym, and Lorene::Vector.

Definition at line 534 of file tensor.C.

References Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), and valence.

◆ sauve()

void Lorene::Tensor::sauve ( FILE *  fd) const
virtual

Save in a binary file.

Reimplemented in Lorene::Tensor_sym, and Lorene::Scalar.

Definition at line 915 of file tensor.C.

References cmp, Lorene::fwrite_be(), n_comp, Lorene::Base_vect::sauve(), Lorene::Itbl::sauve(), triad, type_indice, and valence.

◆ set() [1/4]

Scalar & Lorene::Tensor::set ( const Itbl ind)

Returns the value of a component (read/write version).

Parameters
ind1-D Itbl of size valence containing the values of each index specifing the component, with the following storage convention:
  • ind(0) : value of the first index (1, 2 or 3)
  • ind(1) : value of the second index (1, 2 or 3)
  • and so on...
Returns
modifiable reference on the component specified by ind

Definition at line 663 of file tensor.C.

References cmp, del_deriv(), Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), position(), and valence.

◆ set() [2/4]

Scalar & Lorene::Tensor::set ( int  i1,
int  i2 
)

Returns the value of a component for a tensor of valence 2 (read/write version).

Parameters
i1value of the first index (1, 2 or 3)
i2value of the second index (1, 2 or 3)
Returns
modifiable reference on the component specified by (i1,i2)

Definition at line 615 of file tensor.C.

References cmp, del_deriv(), position(), Lorene::Itbl::set(), and valence.

◆ set() [3/4]

Scalar & Lorene::Tensor::set ( int  i1,
int  i2,
int  i3 
)

Returns the value of a component for a tensor of valence 3 (read/write version).

Parameters
i1value of the first index (1, 2 or 3)
i2value of the second index (1, 2 or 3)
i3value of the third index (1, 2 or 3)
Returns
modifiable reference on the component specified by (i1,i2,i3)

Definition at line 630 of file tensor.C.

References cmp, del_deriv(), position(), Lorene::Itbl::set(), and valence.

◆ set() [4/4]

Scalar & Lorene::Tensor::set ( int  i1,
int  i2,
int  i3,
int  i4 
)

Returns the value of a component for a tensor of valence 4 (read/write version).

Parameters
i1value of the first index (1, 2 or 3)
i2value of the second index (1, 2 or 3)
i3value of the third index (1, 2 or 3)
i4value of the fourth index (1, 2 or 3)
Returns
modifiable reference on the component specified by (i1,i2,i3,i4)

Definition at line 646 of file tensor.C.

References cmp, del_deriv(), position(), Lorene::Itbl::set(), and valence.

◆ set_dependance()

void Lorene::Tensor::set_dependance ( const Metric met) const
protected

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 462 of file tensor.C.

References met_depend, and Lorene::Metric::tensor_depend.

◆ set_der_0x0()

void Lorene::Tensor::set_der_0x0 ( ) const
protected

Sets the pointers on derived quantities to 0x0.

Definition at line 416 of file tensor.C.

References set_der_met_0x0().

◆ set_der_met_0x0()

void Lorene::Tensor::set_der_met_0x0 ( int  i) const
protected

Sets all the i-th components of met_depend , p_derive_cov , etc...

to 0x0.

Definition at line 442 of file tensor.C.

References met_depend, p_derive_con, p_derive_cov, and p_divergence.

◆ set_etat_nondef()

void Lorene::Tensor::set_etat_nondef ( )
virtual

Sets the logical state of all components to ETATNONDEF
(undefined state).

Reimplemented in Lorene::Scalar.

Definition at line 498 of file tensor.C.

References cmp, del_deriv(), n_comp, and Lorene::Scalar::set_etat_nondef().

◆ set_etat_qcq()

void Lorene::Tensor::set_etat_qcq ( )
virtual

Sets the logical state of all components to ETATQCQ
(ordinary state).

Reimplemented in Lorene::Scalar.

Definition at line 490 of file tensor.C.

References cmp, del_deriv(), n_comp, and Lorene::Scalar::set_etat_qcq().

◆ set_etat_zero()

void Lorene::Tensor::set_etat_zero ( )
virtual

Sets the logical state of all components to ETATZERO
(zero state).

Reimplemented in Lorene::Scalar.

Definition at line 506 of file tensor.C.

References cmp, del_deriv(), n_comp, and Lorene::Scalar::set_etat_zero().

◆ set_index_type() [1/2]

int& Lorene::Tensor::set_index_type ( int  i)
inline

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
  • and so on...
Returns
reference on the type that can be modified (COV for a covariant index, CON for a contravariant one)

Definition at line 922 of file tensor.h.

References Lorene::Itbl::set(), and type_indice.

◆ set_index_type() [2/2]

Itbl& Lorene::Tensor::set_index_type ( )
inline

Sets the types of all the indices.

Returns
a reference on the 1-D array of integers (class Itbl ) of size valence that can be modified (COV for a covariant one and CON for a contravariant one)

Definition at line 931 of file tensor.h.

References type_indice.

◆ set_triad()

void Lorene::Tensor::set_triad ( const Base_vect new_triad)

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 528 of file tensor.C.

References triad.

◆ spectral_display()

void Lorene::Tensor::spectral_display ( const char *  comment = 0x0,
double  threshold = 1.e-7,
int  precision = 4,
ostream &  ostr = cout 
) const
virtual

Displays the spectral coefficients and the associated basis functions of each component.

This function shows only the values greater than a given threshold.

Parameters
commentcomment 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 Lorene::Scalar.

Definition at line 883 of file tensor.C.

References cmp, indices(), n_comp, Lorene::Scalar::spectral_display(), and valence.

◆ std_spectral_base()

void Lorene::Tensor::std_spectral_base ( )
virtual

Sets the standard spectal bases of decomposition for each component.

To be used only with valence lower than or equal 2.

Reimplemented in Lorene::Scalar, and Lorene::Vector.

Definition at line 935 of file tensor.C.

References cmp, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Base_vect::identify(), indices(), mp, n_comp, Lorene::Scalar::set_spectral_base(), Lorene::Mg3d::std_base_vect_cart(), Lorene::Mg3d::std_base_vect_spher(), Lorene::Scalar::std_spectral_base(), triad, and valence.

◆ std_spectral_base_odd()

void Lorene::Tensor::std_spectral_base_odd ( )
virtual

Sets the standard odd spectal bases of decomposition for each component.

Currently only implemented for a scalar.

Reimplemented in Lorene::Scalar.

Definition at line 991 of file tensor.C.

References cmp, Lorene::Scalar::std_spectral_base_odd(), and valence.

◆ trace() [1/4]

Tensor Lorene::Tensor::trace ( int  ind1,
int  ind2 
) const

Trace on two different type indices.

Parameters
ind1first index for the contraction, with the following convention :
  • ind1 = 0 : first index of the tensor
  • ind1 = 1 : second index of the tensor
  • and so on...
ind2second index for the contraction

Definition at line 97 of file tensor_calculus.C.

References cmp, get_n_comp(), indices(), mp, position(), Lorene::Itbl::set(), set(), Lorene::Scalar::set_etat_zero(), triad, type_indice, and valence.

◆ trace() [2/4]

Tensor Lorene::Tensor::trace ( int  ind1,
int  ind2,
const Metric gam 
) const

Trace with respect to a given metric.

Parameters
ind1first index for the contraction, with the following convention :
  • ind1 = 0 : first index of the tensor
  • ind1 = 1 : second index of the tensor
  • and so on...
ind2second index for the contraction
gammetric used to raise or lower ind1 in order that it has a opposite type than ind2

Definition at line 156 of file tensor_calculus.C.

References Lorene::Metric::con(), Lorene::contract(), Lorene::Metric::cov(), trace(), type_indice, and valence.

◆ trace() [3/4]

Scalar Lorene::Tensor::trace ( ) const

Trace on two different type indices for a valence 2 tensor.

Definition at line 183 of file tensor_calculus.C.

References mp, operator()(), Lorene::Scalar::set_etat_zero(), type_indice, and valence.

◆ trace() [4/4]

Scalar Lorene::Tensor::trace ( const Metric gam) const

Trace with respect to a given metric for a valence 2 tensor.

Parameters
gammetric used to raise or lower one of the indices, in order to take the trace

Definition at line 200 of file tensor_calculus.C.

References Lorene::Metric::con(), Lorene::contract(), Lorene::Metric::cov(), trace(), type_indice, and valence.

◆ up()

Tensor Lorene::Tensor::up ( int  ind,
const Metric gam 
) const

Computes a new tensor by raising an index of *this.

Parameters
indindex to be raised, with the following convention :
  • ind1 = 0 : first index of the tensor
  • ind1 = 1 : second index of the tensor
  • and so on... (ind must be of covariant type (COV )).
gammetric used to raise the index (contraction with the twice contravariant form of the metric on the index ind ).

Definition at line 228 of file tensor_calculus.C.

References Lorene::Metric::con(), Lorene::contract(), indices(), mp, n_comp, Lorene::Itbl::set(), set(), triad, type_indice, and valence.

◆ up_down()

Tensor Lorene::Tensor::up_down ( const Metric gam) const

Computes a new tensor by raising or lowering all the indices of *this .

Parameters
gammetric used to lower the contravariant indices and raising the covariant ones.

Definition at line 308 of file tensor_calculus.C.

References down(), Tensor(), type_indice, up(), and valence.

Friends And Related Function Documentation

◆ operator* [1/4]

Tensor operator* ( const Tensor a,
const Tensor b 
)
friend

Tensorial product.

Definition at line 113 of file tensor_calculus_ext.C.

◆ operator* [2/4]

Tensor_sym operator* ( const Tensor a,
const Tensor_sym b 
)
friend

Tensorial product with symmetries.

Definition at line 111 of file tensor_sym_calculus.C.

◆ operator* [3/4]

Tensor_sym operator* ( const Tensor_sym a,
const Tensor b 
)
friend

Tensorial product with symmetries.

Definition at line 71 of file tensor_sym_calculus.C.

◆ operator* [4/4]

Tensor_sym operator* ( const Tensor_sym a,
const Tensor_sym b 
)
friend

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 154 of file tensor_sym_calculus.C.

◆ operator+ [1/2]

Scalar operator+ ( const Tensor a,
const Scalar b 
)
friend

Tensor + Scalar. The Tensor must be of valence 0.

Definition at line 123 of file tensor_arithm.C.

◆ operator+ [2/2]

Scalar operator+ ( const Scalar a,
const Tensor b 
)
friend

Scalar + Tensor. The Tensor must be of valence 0.

Definition at line 132 of file tensor_arithm.C.

◆ operator- [1/2]

Scalar operator- ( const Tensor a,
const Scalar b 
)
friend

Tensor - Scalar. The Tensor must be of valence 0.

Definition at line 167 of file tensor_arithm.C.

◆ operator- [2/2]

Scalar operator- ( const Scalar a,
const Tensor b 
)
friend

Scalar - Tensor. The Tensor must be of valence 0.

Definition at line 176 of file tensor_arithm.C.

Member Data Documentation

◆ cmp

Scalar** Lorene::Tensor::cmp
protected

Array of size n_comp of pointers onto the components.

Definition at line 321 of file tensor.h.

◆ met_depend

const Metric* Lorene::Tensor::met_depend[N_MET_MAX]
mutableprotected

Array on the Metric 's which were used to compute derived quantities, like p_derive_cov , etc...

The i-th element of this array is the Metric used to compute the i-th element of p_derive_cov , etc..

Definition at line 333 of file tensor.h.

◆ mp

const Map* const Lorene::Tensor::mp
protected

Mapping on which the numerical values at the grid points are defined.

Definition at line 301 of file tensor.h.

◆ n_comp

int Lorene::Tensor::n_comp
protected

Number of stored components, depending on the symmetry.

Definition at line 318 of file tensor.h.

◆ p_derive_con

Tensor* Lorene::Tensor::p_derive_con[N_MET_MAX]
mutableprotected

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.

Definition at line 349 of file tensor.h.

◆ p_derive_cov

Tensor* Lorene::Tensor::p_derive_cov[N_MET_MAX]
mutableprotected

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.

Definition at line 341 of file tensor.h.

◆ p_divergence

Tensor* Lorene::Tensor::p_divergence[N_MET_MAX]
mutableprotected

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.

Definition at line 357 of file tensor.h.

◆ triad

const Base_vect* Lorene::Tensor::triad
protected

Vectorial basis (triad) with respect to which the tensor components are defined.

Definition at line 309 of file tensor.h.

◆ type_indice

Itbl Lorene::Tensor::type_indice
protected

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.

Definition at line 316 of file tensor.h.

◆ valence

int Lorene::Tensor::valence
protected

Valence of the tensor (0 = scalar, 1 = vector, etc...)

Definition at line 304 of file tensor.h.


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