LORENE
Lorene::Tenseur Class Reference

Tensor handling *** DEPRECATED : use class Tensor instead ***. More...

#include <tenseur.h>

Inheritance diagram for Lorene::Tenseur:
Lorene::Tenseur_sym

Public Member Functions

 Tenseur (const Map &map, const Metrique *met=0x0, double weight=0)
 Constructor for a scalar field. More...
 
 Tenseur (const Cmp &cmp, const Metrique *met=0x0, double weight=0)
 Constructor for a scalar field and from a Cmp . More...
 
 Tenseur (const Map &map, int val, const Itbl &tipe, const Base_vect &triad_i, const Metrique *met=0x0, double weight=0)
 Standard constructor. More...
 
 Tenseur (const Map &map, int val, const Itbl &tipe, const Base_vect *triad_i, const Metrique *met=0x0, double weight=0)
 Standard constructor with the triad passed as a pointer. More...
 
 Tenseur (const Map &map, int val, int tipe, const Base_vect &triad_i, const Metrique *met=0x0, double weight=0)
 Standard constructor when all the indices are of the same type. More...
 
 Tenseur (const Tenseur &)
 Copy constructor. More...
 
 Tenseur (const Tenseur_sym &)
 Constructor from a symmetric tensor. More...
 
 Tenseur (const Map &map, const Base_vect &triad_i, FILE *fich, const Metrique *met=0x0)
 Constructor from a file (see sauve(FILE*) ). More...
 
 Tenseur (const Map &map, FILE *fich, const Metrique *met=0x0)
 Constructor from a file for a scalar field (see sauve(FILE*) ). More...
 
virtual ~Tenseur ()
 Destructor. More...
 
void set_etat_nondef ()
 Sets the logical state to ETATNONDEF (undefined state). More...
 
void set_etat_zero ()
 Sets the logical state to ETATZERO (zero state). More...
 
void set_etat_qcq ()
 Sets the logical state to ETATQCQ (ordinary state). More...
 
void allocate_all ()
 Sets the logical state to ETATQCQ (ordinary state) and performs the memory allocation of all the elements, down to the double arrays of the Tbl s. More...
 
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...
 
void set_poids (double weight)
 Sets the weight for a tensor density. More...
 
void set_metric (const Metrique &met)
 Sets the pointer on the metric for a tensor density. More...
 
virtual void operator= (const Tenseur &tens)
 Assignment to another Tenseur. More...
 
void operator= (const Cmp &field)
 Assignment to a Cmp (scalar field only) More...
 
void operator= (double)
 Assignment to a double (scalar field only, except for zero) More...
 
void operator= (int)
 Assignment to a int (scalar field only, except for zero) More...
 
Cmpset ()
 Read/write for a scalar (see also operator=(const Cmp&) ). More...
 
Cmpset (int)
 Read/write for a vector. More...
 
Cmpset (int, int)
 Read/write for a tensor of valence 2. More...
 
Cmpset (int, int, int)
 Read/write for a tensor of valence 3. More...
 
Cmpset (const Itbl &)
 Read/write in the general case. More...
 
void annule (int l)
 Sets the Tenseur to zero in a given domain. More...
 
void annule (int l_min, int l_max)
 Sets the Tenseur to zero in several domains. More...
 
void set_std_base ()
 Set the standard spectal basis of decomposition for each component. More...
 
void dec_dzpuis ()
 dzpuis -= 1 ; More...
 
void inc_dzpuis ()
 dzpuis += 1 ; More...
 
void dec2_dzpuis ()
 dzpuis -= 2 ; More...
 
void inc2_dzpuis ()
 dzpuis += 2 ; More...
 
void mult_r_zec ()
 Multiplication by r in the external zone. More...
 
Tenseur inverse_poisson_vect (double lambda) const
 Compute $\Delta + \lambda \nabla\nabla$ of *this , *this being of valence 1. More...
 
virtual int donne_place (const Itbl &idx) const
 Returns the position in the Cmp 1-D array c of a component given by its indices. More...
 
virtual Itbl donne_indices (int place) const
 Returns the indices of a component given by its position in the Cmp 1-D array c . More...
 
const Mapget_mp () const
 Returns pointer on the mapping. More...
 
const Base_vectget_triad () const
 Returns the vectorial basis (triad) on which the components are defined. More...
 
int get_etat () const
 Returns the logical state. More...
 
int get_valence () const
 Returns the valence. More...
 
int get_n_comp () const
 Returns the number of components. More...
 
int get_type_indice (int i) const
 Returns the type of the index number i . More...
 
Itbl get_type_indice () const
 Returns the types of all the indices. More...
 
double get_poids () const
 Returns the weight. More...
 
const Metrique * get_metric () const
 Returns a pointer on the metric defining the conformal factor for tensor densities. More...
 
const Cmpoperator() () const
 Read only for a scalar. More...
 
const Cmpoperator() (int) const
 Read only for a vector. More...
 
const Cmpoperator() (int, int) const
 Read only for a tensor of valence 2. More...
 
const Cmpoperator() (int, int, int) const
 Read only for a tensor of valence 3. More...
 
const Cmpoperator() (const Itbl &) const
 Read only in the general case. More...
 
void sauve (FILE *) const
 Save in a file. More...
 
const Tenseurgradient () const
 Returns the gradient of *this (Cartesian coordinates) More...
 
const Tenseurgradient_spher () const
 Returns the gradient of *this (Spherical coordinates) (scalar field only). More...
 
const Tenseurderive_cov (const Metrique &met) const
 Returns the covariant derivative of *this , with respect to met . More...
 
const Tenseurderive_con (const Metrique &) const
 Returns the contravariant derivative of *this , with respect to met . More...
 
const Tenseurcarre_scal (const Metrique &) const
 Returns the scalar square of *this , with respect to met . More...
 
void poisson_vect (double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
 Solves the vectorial Poisson equation : $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$. More...
 
void poisson_vect_tau (double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
 
void poisson_vect_falloff (double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal, int *k_falloff) const
 
void poisson_vect_ylm (double lambda, Param &para, Tenseur &shift, Tenseur &vecteur, Tenseur &scalaire, int nylm, double *intvec) const
 
Tenseur poisson_vect (double lambda, Tenseur &vect, Tenseur &scal) const
 Solves the vectorial Poisson equation $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$. More...
 
Tenseur poisson_vect_tau (double lambda, Tenseur &vect, Tenseur &scal) const
 
Tenseur poisson_vect_falloff (double lambda, Tenseur &vect, Tenseur &scal, int *k_falloff) const
 
Tenseur poisson_vect_ylm (double lambda, Tenseur &vecteur, Tenseur &scalaire, int nylm, double *intvec) const
 
void poisson_vect_oohara (double lambda, Param &par, Tenseur &shift, Tenseur &scal) const
 Solves the vectorial Poisson equation $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$. More...
 
void poisson_vect_oohara_tau (double lambda, Param &par, Tenseur &shift, Tenseur &scal) const
 
Tenseur poisson_vect_oohara (double lambda, Tenseur &scal) const
 Solves the vectorial Poisson equation $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$. More...
 
Tenseur poisson_vect_oohara_tau (double lambda, Tenseur &scal) const
 
void poisson_vect_regu (int k_div, int nzet, double unsgam1, double lambda, Param &par, Tenseur &shift, Tenseur &vect, Tenseur &scal) const
 Solves the vectorial Poisson equation : $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$. More...
 

Protected Member Functions

bool verif () const
 Returns false for a tensor density without a defined metric. More...
 
void new_der_met ()
 Builds the arrays met_depend , p_derive_cov , p_derive_con and p_carre_scal and fills them with null pointers. More...
 
 Tenseur (const Map &map, int val, const Itbl &tipe, int n_comp, const Base_vect &triad_i, const Metrique *met=0x0, double weight=0)
 Constructor used by the derived classes. More...
 
 Tenseur (const Map &, int val, int tipe, int n_comp, const Base_vect &triad_i, const Metrique *met=0x0, double weight=0)
 Constructor used by the derived classes when all the indices are of the same type. More...
 
void del_t ()
 Logical destructor. More...
 
void del_derive_met (int i) const
 Logical destructor of the derivatives depending on the i-th element of *met_depend . More...
 
void del_derive () const
 Logical destructor of all the derivatives. More...
 
void set_der_met_0x0 (int i) const
 Sets the pointers of the derivatives depending on the i-th element of *met_depend to zero (as well as that i-th element). More...
 
void set_der_0x0 () const
 Sets the pointers of all the derivatives to zero. More...
 
virtual void fait_gradient () const
 Calculates, if needed, the gradient of *this . More...
 
void fait_gradient_spher () const
 Calculates, if needed, the gradient of *this in a spherical orthonormal basis (scalar field only). More...
 
virtual void fait_derive_cov (const Metrique &met, int i) const
 Calculates, if needed, the covariant derivative of *this , with respect to met . More...
 
virtual void fait_derive_con (const Metrique &, int i) const
 Calculates, if needed, the contravariant derivative of *this , with respect to met . More...
 
void fait_carre_scal (const Metrique &, int i) const
 Calculates, if needed, the scalar square of *this , with respect to met . More...
 
void set_dependance (const Metrique &met) const
 To be used to describe the fact that the derivatives members have been calculated with met . More...
 
int get_place_met (const Metrique &metre) const
 Returns the position of the pointer on metre in the array met_depend . More...
 

Protected Attributes

const Map *const mp
 Reference mapping. More...
 
int valence
 Valence. More...
 
const Base_vecttriad
 Vectorial basis (triad) with respect to which the tensor components are defined. More...
 
Itbl type_indice
 Array of size valence contening the type of each index, COV for a covariant one and CON for a contravariant one. More...
 
int n_comp
 Number of components, depending on the symmetry. More...
 
int etat
 Logical state ETATZERO , ETATQCQ or ETATNONDEF. More...
 
Cmp ** c
 The components. More...
 
double poids
 For tensor densities: the weight. More...
 
const Metrique * metric
 For tensor densities: the metric defining the conformal factor. More...
 
const Metrique ** met_depend
 Array of pointers on the Metrique 's used to calculate derivatives members. More...
 
Tenseurp_gradient
 Pointer on the gradient of *this . More...
 
Tenseurp_gradient_spher
 Pointer on the gradient of *this in a spherical orthonormal basis (scalar field only). More...
 
Tenseur ** p_derive_cov
 Array of pointers on the covariant derivatives of *this
with respect to the corresponding metric in *met_depend . More...
 
Tenseur ** p_derive_con
 Array of pointers on the contravariant derivatives of *this
with respect to the corresponding metric in *met_depend . More...
 
Tenseur ** p_carre_scal
 Array of pointers on the scalar squares of *this
with respect to the corresponding metric in *met_depend . More...
 

Friends

class Tenseur_sym
 
class Metrique
 
ostream & operator<< (ostream &, const Tenseur &)
 
Tenseur operator* (const Tenseur &, const Tenseur &)
 Tensorial product. More...
 
Tenseur operator% (const Tenseur &, const Tenseur &)
 Tensorial product with desaliasing. More...
 
Tenseur contract (const Tenseur &, int id1, int id2)
 Self contraction of two indices of a Tenseur . More...
 
Tenseur contract (const Tenseur &, int id1, const Tenseur &, int id2)
 Contraction of two Tenseur . More...
 
Tenseur contract_desal (const Tenseur &, int id1, const Tenseur &, int id2)
 
Tenseur flat_scalar_prod (const Tenseur &t1, const Tenseur &t2)
 Scalar product of two Tenseur when the metric is $\delta_{ij}$: performs the contraction of the last index of t1 with the first one of t2 , irrespective of the type of these indices. More...
 
Tenseur flat_scalar_prod_desal (const Tenseur &t1, const Tenseur &t2)
 Same as flat_scalar_prod but with desaliasing. More...
 
Tenseur manipule (const Tenseur &, const Metrique &, int idx)
 Raise or lower the index idx depending on its type, using the given Metrique . More...
 
Tenseur manipule (const Tenseur &, const Metrique &)
 Raise or lower all the indices, depending on their type, using the given Metrique . More...
 
Tenseur skxk (const Tenseur &)
 Contraction of the last index of (*this) with $x^k$ or $x_k$, depending on the type of S . More...
 
Tenseur lie_derive (const Tenseur &, const Tenseur &, const Metrique *)
 Lie Derivative of t with respect to x . More...
 

Detailed Description

Tensor handling *** DEPRECATED : use class Tensor instead ***.

()

This class is intended to store the components of a tensorial field in a specific basis. Tensor densities can also be stored. A tensor density $\tau^{i_1\ldots i_p}_{j_1\ldots j_q}$ is defined by: $ \tau^{i_1\ldots i_p}_{j_1\ldots j_q} = \gamma^{\frac{n}{2}} T^{i_1\ldots i_p}_{j_1\ldots j_q}$ where T is a q -covariant p -contravariant tensor and $\gamma$ is the determinant of the used 3-metric. n is called the weight of the tensor density.

All this is 3D meaning that the indices go from 0 to 2. Moreover, the components are described in orthonormal bases.

When first constructed, the memory for each component is not allocated.

Definition at line 304 of file tenseur.h.

Constructor & Destructor Documentation

◆ Tenseur() [1/11]

Lorene::Tenseur::Tenseur ( const Map map,
const Metrique *  met = 0x0,
double  weight = 0 
)
explicit

Constructor for a scalar field.

Definition at line 225 of file tenseur.C.

◆ Tenseur() [2/11]

Lorene::Tenseur::Tenseur ( const Cmp cmp,
const Metrique *  met = 0x0,
double  weight = 0 
)
explicit

Constructor for a scalar field and from a Cmp .

Definition at line 240 of file tenseur.C.

References Lorene::Cmp::get_etat().

◆ Tenseur() [3/11]

Lorene::Tenseur::Tenseur ( const Map map,
int  val,
const Itbl tipe,
const Base_vect triad_i,
const Metrique *  met = 0x0,
double  weight = 0 
)

Standard constructor.

Parameters
mapthe mapping
valvalence of the tensor
tipe1-D 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
metfor tensor densities only: a pointer on the metric defining the conformal factor
weightfor tensor densities: the weight

Definition at line 262 of file tenseur.C.

◆ Tenseur() [4/11]

Lorene::Tenseur::Tenseur ( const Map map,
int  val,
const Itbl tipe,
const Base_vect triad_i,
const Metrique *  met = 0x0,
double  weight = 0 
)

Standard constructor with the triad passed as a pointer.

Parameters
mapthe mapping
valvalence of the tensor
tipe1-D 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)
metfor tensor densities only: a pointer on the metric defining the conformal factor
weightfor tensor densities: the weight

Definition at line 284 of file tenseur.C.

◆ Tenseur() [5/11]

Lorene::Tenseur::Tenseur ( const Map map,
int  val,
int  tipe,
const Base_vect triad_i,
const Metrique *  met = 0x0,
double  weight = 0 
)

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

Parameters
mapthe mapping
valvalence of the tensor
tipethe type of the indices.
triad_ivectorial basis (triad) with respect to which the tensor components are defined.
metfor tensor densities only: a pointer on the metric defining the conformal factor
weightfor tensor densities: the weight

Definition at line 313 of file tenseur.C.

◆ Tenseur() [6/11]

Lorene::Tenseur::Tenseur ( const Tenseur source)

◆ Tenseur() [7/11]

Lorene::Tenseur::Tenseur ( const Tenseur_sym source)
explicit

◆ Tenseur() [8/11]

Lorene::Tenseur::Tenseur ( const Map map,
const Base_vect triad_i,
FILE *  fich,
const Metrique *  met = 0x0 
)

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*) .
metfor tensor densities only: a pointer on the metric defining the conformal factor

Definition at line 427 of file tenseur.C.

References Lorene::Base_vect::bvect_from_file(), c, etat, Lorene::fread_be(), n_comp, triad, and valence.

◆ Tenseur() [9/11]

Lorene::Tenseur::Tenseur ( const Map map,
FILE *  fich,
const Metrique *  met = 0x0 
)

Constructor from a file for a scalar field (see sauve(FILE*) ).

Parameters
mapthe mapping
fichfile which has been created by the function sauve(FILE*) .
metfor tensor densities only: a pointer on the metric defining the conformal factor

Definition at line 464 of file tenseur.C.

References c, etat, Lorene::fread_be(), n_comp, triad, and valence.

◆ Tenseur() [10/11]

Lorene::Tenseur::Tenseur ( const Map map,
int  val,
const Itbl tipe,
int  n_comp,
const Base_vect triad_i,
const Metrique *  met = 0x0,
double  weight = 0 
)
protected

Constructor used by the derived classes.

Parameters
mapthe mapping
valvalence of the tensor
tipe1-D 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_compthe number of components.
triad_ivectorial basis (triad) with respect to which the tensor components are defined
metfor tensor densities only: a pointer on the metric defining the conformal factor
weightfor tensor densities: the weight

Definition at line 501 of file tenseur.C.

◆ Tenseur() [11/11]

Lorene::Tenseur::Tenseur ( const Map map,
int  val,
int  tipe,
int  n_comp,
const Base_vect triad_i,
const Metrique *  met = 0x0,
double  weight = 0 
)
protected

Constructor used by the derived classes when all the indices are of the same type.

Parameters
mapthe mapping
valvalence of the tensor
tipethe type of the indices.
n_compthe number of components.
triad_ivectorial basis (triad) with respect to which the tensor components are defined
metfor tensor densities only: a pointer on the metric defining the conformal factor
weightfor tensor densities: the weight

Definition at line 525 of file tenseur.C.

◆ ~Tenseur()

Lorene::Tenseur::~Tenseur ( )
virtual

Destructor.

Definition at line 549 of file tenseur.C.

References c, del_t(), met_depend, p_carre_scal, p_derive_con, and p_derive_cov.

Member Function Documentation

◆ allocate_all()

void Lorene::Tenseur::allocate_all ( )

Sets the logical state to ETATQCQ (ordinary state) and 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 Tenseur -> Cmp -> Valeur -> Mtbl -> Tbl .

Definition at line 673 of file tenseur.C.

References Lorene::Cmp::allocate_all(), c, n_comp, and set_etat_qcq().

◆ annule() [1/2]

void Lorene::Tenseur::annule ( int  l)

Sets the Tenseur to zero in a given domain.

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

Definition at line 916 of file tenseur.C.

◆ annule() [2/2]

void Lorene::Tenseur::annule ( int  l_min,
int  l_max 
)

Sets the Tenseur to zero in several domains.

Parameters
l_min[input] The Tenseur 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 \c annule(0,nz-1), where \c nz  is the total number
of domains, is equivalent to set_etat_zero() .

Definition at line 921 of file tenseur.C.

References etat, Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), mp, and set_etat_zero().

◆ carre_scal()

const Tenseur & Lorene::Tenseur::carre_scal ( const Metrique &  metre) const

Returns the scalar square of *this , with respect to met .

Definition at line 1593 of file tenseur.C.

References fait_carre_scal(), get_place_met(), p_carre_scal, and set_dependance().

◆ change_triad()

void Lorene::Tenseur::change_triad ( const Base_vect new_triad)

Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.

Definition at line 684 of file tenseur.C.

References Lorene::Base_vect::change_basis().

◆ dec2_dzpuis()

void Lorene::Tenseur::dec2_dzpuis ( )

dzpuis -= 2 ;

Definition at line 1146 of file tenseur.C.

References etat.

◆ dec_dzpuis()

void Lorene::Tenseur::dec_dzpuis ( )

dzpuis -= 1 ;

Definition at line 1120 of file tenseur.C.

References etat.

◆ del_derive()

void Lorene::Tenseur::del_derive ( ) const
protected

Logical destructor of all the derivatives.

Definition at line 589 of file tenseur.C.

References del_derive_met(), p_gradient, p_gradient_spher, and set_der_0x0().

◆ del_derive_met()

void Lorene::Tenseur::del_derive_met ( int  i) const
protected

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

Definition at line 570 of file tenseur.C.

References met_depend, p_carre_scal, p_derive_con, p_derive_cov, and set_der_met_0x0().

◆ del_t()

void Lorene::Tenseur::del_t ( )
protected

Logical destructor.

Definition at line 561 of file tenseur.C.

References c, del_derive(), and n_comp.

◆ derive_con()

const Tenseur & Lorene::Tenseur::derive_con ( const Metrique &  metre) const

Returns the contravariant derivative of *this , with respect to met .

Definition at line 1584 of file tenseur.C.

References fait_derive_con(), get_place_met(), p_derive_con, and set_dependance().

◆ derive_cov()

const Tenseur & Lorene::Tenseur::derive_cov ( const Metrique &  met) const

Returns the covariant derivative of *this , with respect to met .

Definition at line 1570 of file tenseur.C.

References fait_derive_cov(), get_place_met(), gradient(), p_derive_cov, set_dependance(), and valence.

◆ donne_indices()

Itbl Lorene::Tenseur::donne_indices ( int  place) const
virtual

Returns the indices of a component given by its position in the Cmp 1-D array c .

Returns
1-D array of integers (Itbl ) of size valence giving the value of each index for the component located at the position place in the Cmp 1-D array c . Each element of this Itbl is 0, 1 or 2, which corresponds to spatial indices 1, 2 or 3 respectively. If (*this) is a scalar the function returns an undefined Itbl .

Reimplemented in Lorene::Tenseur_sym.

Definition at line 720 of file tenseur.C.

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

◆ donne_place()

int Lorene::Tenseur::donne_place ( const Itbl idx) const
virtual

Returns the position in the Cmp 1-D array c of a component given by its indices.

Returns
position in the Cmp 1-D array c
corresponding to the indices given in idx . idx must be a 1-D Itbl of size valence , each element of which must be 0, 1 or 2, corresponding to spatial indices 1, 2 or 3 respectively.

Reimplemented in Lorene::Tenseur_sym.

Definition at line 706 of file tenseur.C.

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

◆ fait_carre_scal()

void Lorene::Tenseur::fait_carre_scal ( const Metrique &  met,
int  i 
) const
protected

Calculates, if needed, the scalar square of *this , with respect to met .

The result is in *p_carre_scal[i]

Definition at line 1533 of file tenseur.C.

References contract, manipule, p_carre_scal, Tenseur(), and valence.

◆ fait_derive_con()

void Lorene::Tenseur::fait_derive_con ( const Metrique &  metre,
int  i 
) const
protectedvirtual

Calculates, if needed, the contravariant derivative of *this , with respect to met .

The result is in *p_derive_con[i]

Reimplemented in Lorene::Tenseur_sym.

Definition at line 1515 of file tenseur.C.

References change_triad(), contract, derive_cov(), gradient(), p_derive_con, Tenseur(), and valence.

◆ fait_derive_cov()

void Lorene::Tenseur::fait_derive_cov ( const Metrique &  met,
int  i 
) const
protectedvirtual

Calculates, if needed, the covariant derivative of *this , with respect to met .

The result is in *p_derive_cov[i]

Reimplemented in Lorene::Tenseur_sym.

Definition at line 1445 of file tenseur.C.

References etat.

◆ fait_gradient()

void Lorene::Tenseur::fait_gradient ( ) const
protectedvirtual

Calculates, if needed, the gradient of *this .

The result is in *p_gradient

Reimplemented in Lorene::Tenseur_sym.

Definition at line 1361 of file tenseur.C.

References etat.

◆ fait_gradient_spher()

void Lorene::Tenseur::fait_gradient_spher ( ) const
protected

Calculates, if needed, the gradient of *this in a spherical orthonormal basis (scalar field only).

The result is in *p_gradient_spher

Definition at line 1408 of file tenseur.C.

References etat.

◆ get_etat()

int Lorene::Tenseur::get_etat ( ) const
inline

Returns the logical state.

Definition at line 710 of file tenseur.h.

References etat.

◆ get_metric()

const Metrique* Lorene::Tenseur::get_metric ( ) const
inline

Returns a pointer on the metric defining the conformal factor for tensor densities.

Otherwise (case of a pure tensor), it returns 0x0.

Definition at line 748 of file tenseur.h.

References metric.

◆ get_mp()

const Map* Lorene::Tenseur::get_mp ( ) const
inline

Returns pointer on the mapping.

Definition at line 702 of file tenseur.h.

References mp.

◆ get_n_comp()

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

Returns the number of components.

Definition at line 716 of file tenseur.h.

References n_comp.

◆ get_place_met()

int Lorene::Tenseur::get_place_met ( const Metrique &  metre) const
protected

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

Definition at line 614 of file tenseur.C.

References met_depend.

◆ get_poids()

double Lorene::Tenseur::get_poids ( ) const
inline

Returns the weight.

Definition at line 741 of file tenseur.h.

References poids.

◆ get_triad()

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

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

Definition at line 707 of file tenseur.h.

References triad.

◆ get_type_indice() [1/2]

int Lorene::Tenseur::get_type_indice ( int  i) const
inline

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

Definition at line 729 of file tenseur.h.

References type_indice.

◆ get_type_indice() [2/2]

Itbl Lorene::Tenseur::get_type_indice ( ) const
inline

Returns the types of all the indices.

Returns
1-D Itbl of size valence containing the type of each index, COV for a covariant one and CON
for a contravariant one.

Definition at line 738 of file tenseur.h.

References type_indice.

◆ get_valence()

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

Returns the valence.

Definition at line 713 of file tenseur.h.

References valence.

◆ gradient()

const Tenseur & Lorene::Tenseur::gradient ( ) const

Returns the gradient of *this (Cartesian coordinates)

Definition at line 1558 of file tenseur.C.

References fait_gradient(), and p_gradient.

◆ gradient_spher()

const Tenseur & Lorene::Tenseur::gradient_spher ( ) const

Returns the gradient of *this (Spherical coordinates) (scalar field only).

Definition at line 1564 of file tenseur.C.

References fait_gradient_spher(), and p_gradient_spher.

◆ inc2_dzpuis()

void Lorene::Tenseur::inc2_dzpuis ( )

dzpuis += 2 ;

Definition at line 1159 of file tenseur.C.

References etat.

◆ inc_dzpuis()

void Lorene::Tenseur::inc_dzpuis ( )

dzpuis += 1 ;

Definition at line 1133 of file tenseur.C.

References etat.

◆ inverse_poisson_vect()

Tenseur Lorene::Tenseur::inverse_poisson_vect ( double  lambda) const

Compute $\Delta + \lambda \nabla\nabla$ of *this , *this being of valence 1.

Definition at line 62 of file tenseur_inv_pois_vect.C.

References etat, and valence.

◆ mult_r_zec()

void Lorene::Tenseur::mult_r_zec ( )

Multiplication by r in the external zone.

Definition at line 1172 of file tenseur.C.

References etat.

◆ new_der_met()

void Lorene::Tenseur::new_der_met ( )
protected

Builds the arrays met_depend , p_derive_cov , p_derive_con and p_carre_scal and fills them with null pointers.

Definition at line 212 of file tenseur.C.

References met_depend, p_carre_scal, p_derive_con, p_derive_cov, and set_der_0x0().

◆ operator()() [1/5]

const Cmp & Lorene::Tenseur::operator() ( ) const

Read only for a scalar.

Definition at line 959 of file tenseur.C.

References etat, and valence.

◆ operator()() [2/5]

const Cmp & Lorene::Tenseur::operator() ( int  indice) const

Read only for a vector.

Definition at line 988 of file tenseur.C.

References etat, and valence.

◆ operator()() [3/5]

const Cmp & Lorene::Tenseur::operator() ( int  indice1,
int  indice2 
) const

Read only for a tensor of valence 2.

Definition at line 1016 of file tenseur.C.

References etat, and valence.

◆ operator()() [4/5]

const Cmp & Lorene::Tenseur::operator() ( int  indice1,
int  indice2,
int  indice3 
) const

Read only for a tensor of valence 3.

Definition at line 1050 of file tenseur.C.

References etat, and valence.

◆ operator()() [5/5]

const Cmp & Lorene::Tenseur::operator() ( const Itbl indices) const

Read only in the general case.

Definition at line 1087 of file tenseur.C.

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

◆ operator=() [1/4]

void Lorene::Tenseur::operator= ( const Tenseur tens)
virtual

Assignment to another Tenseur.

Reimplemented in Lorene::Tenseur_sym.

Definition at line 734 of file tenseur.C.

References etat, metric, poids, triad, type_indice, and valence.

◆ operator=() [2/4]

void Lorene::Tenseur::operator= ( const Cmp field)

Assignment to a Cmp (scalar field only)

Definition at line 777 of file tenseur.C.

References Lorene::Cmp::get_etat(), metric, poids, and valence.

◆ operator=() [3/4]

void Lorene::Tenseur::operator= ( double  x)

Assignment to a double (scalar field only, except for zero)

Definition at line 808 of file tenseur.C.

References c, metric, poids, set_etat_qcq(), set_etat_zero(), and valence.

◆ operator=() [4/4]

void Lorene::Tenseur::operator= ( int  x)

Assignment to a int (scalar field only, except for zero)

Definition at line 823 of file tenseur.C.

References c, metric, poids, set_etat_qcq(), set_etat_zero(), and valence.

◆ poisson_vect() [1/2]

void Lorene::Tenseur::poisson_vect ( double  lambda,
Param par,
Tenseur shift,
Tenseur vect,
Tenseur scal 
) const

Solves the vectorial Poisson equation : $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$.

with $\lambda \not= 1$.

*this must be given with dzpuis = 4.

It uses the Shibata scheme, where $N^i$ is given by :

\[ N^i = \frac{1}{2}\frac{\lambda+2}{\lambda+1}W^i-\frac{1}{2} \frac{\lambda}{\lambda+1}\left(\nabla^i\chi+\nabla^iW^kx_k\right) \]

with $\Delta W^i = S^i$ and $\Delta \chi = -x_kS^k$.

Parameters
lambda[input] $\lambda$.
par[input/output] see Map::donne_para_poisson_vect.
shift[input] solution $N^i$ at the previous step. Zero if nothing is known.
shift[output] solution at this step.
vect[input/output] the same thing than for shift but for $W^i$.
scal[input/output] the same thing than for shift but for $\chi$.

Definition at line 121 of file tenseur_pde.C.

References etat, get_type_indice(), get_valence(), type_indice, and valence.

◆ poisson_vect() [2/2]

Tenseur Lorene::Tenseur::poisson_vect ( double  lambda,
Tenseur vect,
Tenseur scal 
) const

Solves the vectorial Poisson equation $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$.

with $\lambda \not= 1$.

*this must be given with dzpuis = 4.

It uses the Shibata scheme, where $N^i$ is given by :

\[ N^i = \frac{1}{2}\frac{\lambda+2}{\lambda+1}W^i-\frac{1}{2} \frac{\lambda}{\lambda+1}\left(\nabla^i\chi+\nabla^iW^kx_k\right) \]

with $\Delta W^i = S^i$ and $\Delta \chi = -x_kS^k$.

This version is to be used only with an affine mapping.

Parameters
lambda[input] $\lambda$.
vect[input] $W^i$ at the previous step. Zero if nothing is known.
vect[output] $W^i$ at this step.
scal[input/output] the same thing than for shift but for $\chi$.
Returns
the solution $N^i$.

Definition at line 203 of file tenseur_pde.C.

References metric, mp, poids, poisson_vect(), set_etat_qcq(), triad, type_indice, and valence.

◆ poisson_vect_oohara() [1/2]

void Lorene::Tenseur::poisson_vect_oohara ( double  lambda,
Param par,
Tenseur shift,
Tenseur scal 
) const

Solves the vectorial Poisson equation $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$.

with $\lambda \not= 1$.

*this must be given with dzpuis = 3 or 4 and be continuous.

It uses the Oohara scheme, where $N^i$ is given by

\[ \Delta N^i = S^i-\lambda \nabla^i \chi \]

with $\chi$ solution of :

\[ \Delta \chi = \frac{1}{\lambda+1}\nabla_k S^k \]

Parameters
lambda[input] $\lambda$.
par[input/output] see Map::donne_para_poisson_vect.
shift[input] solution $N^i$ at the previous step. Zero if nothing is known.
shift[output] solution at this step.
scal[input/output] the same thing than for shift but for $\chi$.

Definition at line 221 of file tenseur_pde.C.

References etat, get_type_indice(), get_valence(), type_indice, and valence.

◆ poisson_vect_oohara() [2/2]

Tenseur Lorene::Tenseur::poisson_vect_oohara ( double  lambda,
Tenseur scal 
) const

Solves the vectorial Poisson equation $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$.

with $\lambda \not= 1$.

*this must be given with dzpuis = 3 or 4 and be continuous.

This version is to be used only with an affine mapping.

It uses the Oohara scheme, where $N^i$ is given by :

\[ \Delta N^i = S^i-\lambda \nabla^i \chi \]

with $\chi$ solution of :

\[ \Delta \chi = \frac{1}{\lambda+1}\nabla_k S^k \]

This version is to be used only with an affine mapping.

Parameters
lambda[input] $\lambda$.
scal[input] $\chi$ at the previous step. Zero if nothing is known.
scal[output] $\chi$ at this step.
Returns
the solution $N^i$.

Definition at line 293 of file tenseur_pde.C.

References metric, mp, poids, poisson_vect_oohara(), set_etat_qcq(), triad, type_indice, and valence.

◆ poisson_vect_regu()

void Lorene::Tenseur::poisson_vect_regu ( int  k_div,
int  nzet,
double  unsgam1,
double  lambda,
Param par,
Tenseur shift,
Tenseur vect,
Tenseur scal 
) const

Solves the vectorial Poisson equation : $\Delta N^i +\lambda \nabla^i \nabla_k N^k = S^i$.

with $\lambda \not= 1$ by regularizing the source term.

*this must be given with dzpuis = 4.

It uses the Shibata scheme, where $N^i$ is given by :

\[ N^i = \frac{1}{2}\frac{\lambda+2}{\lambda+1}W^i-\frac{1}{2} \frac{\lambda}{\lambda+1}\left(\nabla^i\chi+\nabla^iW^kx_k\right) \]

with $\Delta W^i = S^i$ and $\Delta \chi = -x_kS^k$.

Parameters
k_div[input] regularization degree.
nzet[input] number of domains covering a star.
unsgam1[input] $1/(\gamma - 1)$.
lambda[input] $\lambda$.
par[input/output] see Map::donne_para_poisson_vect.
shift[input] solution $N^i$ at the previous step. Zero if nothing is known.
shift[output] solution at this step.
vect[input/output] the same thing than for shift but for $W^i$.
scal[input/output] the same thing than for shift but for $\chi$.

Definition at line 74 of file tenseur_pde_regu.C.

References etat, get_type_indice(), get_valence(), type_indice, and valence.

◆ sauve()

void Lorene::Tenseur::sauve ( FILE *  fd) const

Save in a file.

Definition at line 1341 of file tenseur.C.

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

◆ set() [1/5]

Cmp & Lorene::Tenseur::set ( )

Read/write for a scalar (see also operator=(const Cmp&) ).

Definition at line 840 of file tenseur.C.

References del_derive(), and etat.

◆ set() [2/5]

Cmp & Lorene::Tenseur::set ( int  ind)

Read/write for a vector.

Definition at line 850 of file tenseur.C.

References del_derive(), etat, and valence.

◆ set() [3/5]

Cmp & Lorene::Tenseur::set ( int  ind1,
int  ind2 
)

Read/write for a tensor of valence 2.

Definition at line 861 of file tenseur.C.

References del_derive(), etat, and valence.

◆ set() [4/5]

Cmp & Lorene::Tenseur::set ( int  ind1,
int  ind2,
int  ind3 
)

Read/write for a tensor of valence 3.

Definition at line 880 of file tenseur.C.

References del_derive(), etat, and valence.

◆ set() [5/5]

Cmp & Lorene::Tenseur::set ( const Itbl indices)

Read/write in the general case.

Definition at line 900 of file tenseur.C.

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

◆ set_dependance()

void Lorene::Tenseur::set_dependance ( const Metrique &  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 624 of file tenseur.C.

References met_depend.

◆ set_der_0x0()

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

Sets the pointers of all the derivatives to zero.

Definition at line 607 of file tenseur.C.

References p_gradient, p_gradient_spher, and set_der_met_0x0().

◆ set_der_met_0x0()

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

Sets the pointers of the derivatives depending on the i-th element of *met_depend to zero (as well as that i-th element).

Definition at line 599 of file tenseur.C.

References met_depend, p_carre_scal, p_derive_con, and p_derive_cov.

◆ set_etat_nondef()

void Lorene::Tenseur::set_etat_nondef ( )

Sets the logical state to ETATNONDEF (undefined state).

The components are not allocated.

Definition at line 666 of file tenseur.C.

References del_t(), and etat.

◆ set_etat_qcq()

void Lorene::Tenseur::set_etat_qcq ( )

Sets the logical state to ETATQCQ (ordinary state).

The components are now allocated and set to ETATNONDEF .

Definition at line 652 of file tenseur.C.

References c, del_derive(), etat, mp, and n_comp.

◆ set_etat_zero()

void Lorene::Tenseur::set_etat_zero ( )

Sets the logical state to ETATZERO (zero state).

The components are not allocated.

Definition at line 661 of file tenseur.C.

References del_t(), and etat.

◆ set_metric()

void Lorene::Tenseur::set_metric ( const Metrique &  met)

Sets the pointer on the metric for a tensor density.

Definition at line 701 of file tenseur.C.

References metric.

◆ set_poids()

void Lorene::Tenseur::set_poids ( double  weight)

Sets the weight for a tensor density.

Definition at line 696 of file tenseur.C.

References poids.

◆ set_std_base()

void Lorene::Tenseur::set_std_base ( )

Set the standard spectal basis of decomposition for each component.

To be used only with valence strictly lower than 3.

Definition at line 1186 of file tenseur.C.

References etat.

◆ set_triad()

void Lorene::Tenseur::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 c ). In order to change them coherently with the new basis, the function change_triad(const Base_vect&) must be called instead.

Definition at line 690 of file tenseur.C.

References triad.

◆ verif()

bool Lorene::Tenseur::verif ( ) const
protected

Returns false for a tensor density without a defined metric.

Definition at line 208 of file tenseur.C.

References metric, and poids.

Friends And Related Function Documentation

◆ contract [1/2]

Tenseur contract ( const Tenseur source,
int  id1,
int  id2 
)
friend

Self contraction of two indices of a Tenseur .

The two indices must be of different type, i.e. covariant and contravariant, or contravariant and covariant.

Parameters
id1[input] number of the first index for the contraction; id1 must be strictly lower than the valence of the tensor and obeys the following convention:
  • id1 = 0 : first index
  • id1 = 1 : second index
  • and so on...
id2[input] number of the second index for the contraction; id2 must be strictly lower than the valence of the tensor and obeys the following convention:
  • id2 = 0 : first index
  • id2 = 1 : second index
  • and so on...

Definition at line 282 of file tenseur_operateur.C.

◆ contract [2/2]

Tenseur contract ( const Tenseur t1,
int  id1,
const Tenseur t2,
int  id2 
)
friend

Contraction of two Tenseur .

The two indices must be of different type, i.e. covariant and contravariant, or contravariant and covariant.

Parameters
id1[input] number of the index of contraction for the first Tenseur ; id1 must be strictly lower than the valence of the tensor and obeys the following convention:
  • id1 = 0 : first index
  • id1 = 1 : second index
  • and so on...
id2[input] number of index of contraction for the second one; id2 must be strictly lower than the valence of the tensor and obeys the following convention:
  • id2 = 0 : first index
  • id2 = 1 : second index
  • and so on...

Definition at line 351 of file tenseur_operateur.C.

◆ flat_scalar_prod

Tenseur flat_scalar_prod ( const Tenseur t1,
const Tenseur t2 
)
friend

Scalar product of two Tenseur when the metric is $\delta_{ij}$: performs the contraction of the last index of t1 with the first one of t2 , irrespective of the type of these indices.

Definition at line 656 of file tenseur_operateur.C.

◆ flat_scalar_prod_desal

Tenseur flat_scalar_prod_desal ( const Tenseur t1,
const Tenseur t2 
)
friend

Same as flat_scalar_prod but with desaliasing.

Definition at line 738 of file tenseur_operateur.C.

◆ lie_derive

Tenseur lie_derive ( const Tenseur t,
const Tenseur x,
const Metrique *  met = 0x0 
)
friend

Lie Derivative of t with respect to x .

If no other argument is given, it uses partial derivatives with respect to cartesian coordinates to calculate the result (this is the default). Otherwise, it uses the covariant derivative associated to the metric given as last argument.

Definition at line 819 of file tenseur_operateur.C.

◆ manipule [1/2]

Tenseur manipule ( const Tenseur t1,
const Metrique &  met,
int  idx 
)
friend

Raise or lower the index idx depending on its type, using the given Metrique .

Definition at line 512 of file tenseur_operateur.C.

◆ manipule [2/2]

Tenseur manipule ( const Tenseur t1,
const Metrique &  met 
)
friend

Raise or lower all the indices, depending on their type, using the given Metrique .

Definition at line 565 of file tenseur_operateur.C.

◆ operator%

Tenseur operator% ( const Tenseur t1,
const Tenseur t2 
)
friend

Tensorial product with desaliasing.

Definition at line 202 of file tenseur_operateur.C.

◆ operator*

Tenseur operator* ( const Tenseur t1,
const Tenseur t2 
)
friend

Tensorial product.

Definition at line 122 of file tenseur_operateur.C.

◆ skxk

Tenseur skxk ( const Tenseur source)
friend

Contraction of the last index of (*this) with $x^k$ or $x_k$, depending on the type of S .

The calculation is performed to avoid singularities in the external zone. This is done only for a flat metric.

Definition at line 583 of file tenseur_operateur.C.

Member Data Documentation

◆ c

Cmp** Lorene::Tenseur::c
protected

The components.

Definition at line 325 of file tenseur.h.

◆ etat

int Lorene::Tenseur::etat
protected

Logical state ETATZERO , ETATQCQ or ETATNONDEF.

Definition at line 324 of file tenseur.h.

◆ met_depend

const Metrique** Lorene::Tenseur::met_depend
protected

Array of pointers on the Metrique 's used to calculate derivatives members.

If no such member has been calculated the pointer is set to zero. The array has the size N_MET_MAX abd the i-th corresponds to the i-th ones in p_derive_cov and p_derive_con.

Definition at line 340 of file tenseur.h.

◆ metric

const Metrique* Lorene::Tenseur::metric
protected

For tensor densities: the metric defining the conformal factor.

Definition at line 328 of file tenseur.h.

◆ mp

const Map* const Lorene::Tenseur::mp
protected

Reference mapping.

Definition at line 309 of file tenseur.h.

◆ n_comp

int Lorene::Tenseur::n_comp
protected

Number of components, depending on the symmetry.

Definition at line 323 of file tenseur.h.

◆ p_carre_scal

Tenseur** Lorene::Tenseur::p_carre_scal
protected

Array of pointers on the scalar squares of *this
with respect to the corresponding metric in *met_depend .

It is set to zero if it has not been calculated yet.

Definition at line 375 of file tenseur.h.

◆ p_derive_con

Tenseur** Lorene::Tenseur::p_derive_con
protected

Array of pointers on the contravariant derivatives of *this
with respect to the corresponding metric in *met_depend .

It is set to zero if it has not been calculated yet.

Definition at line 368 of file tenseur.h.

◆ p_derive_cov

Tenseur** Lorene::Tenseur::p_derive_cov
protected

Array of pointers on the covariant derivatives of *this
with respect to the corresponding metric in *met_depend .

It is set to zero if it has not been calculated yet, or for a scalar field.

Definition at line 361 of file tenseur.h.

◆ p_gradient

Tenseur* Lorene::Tenseur::p_gradient
mutableprotected

Pointer on the gradient of *this .

It is set to zero if it has not been calculated yet.

Definition at line 346 of file tenseur.h.

◆ p_gradient_spher

Tenseur* Lorene::Tenseur::p_gradient_spher
mutableprotected

Pointer on the gradient of *this in a spherical orthonormal basis (scalar field only).

It is set to zero if it has not been calculated yet.

Definition at line 353 of file tenseur.h.

◆ poids

double Lorene::Tenseur::poids
protected

For tensor densities: the weight.

Definition at line 326 of file tenseur.h.

◆ triad

const Base_vect* Lorene::Tenseur::triad
protected

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

Definition at line 315 of file tenseur.h.

◆ type_indice

Itbl Lorene::Tenseur::type_indice
protected

Array of size valence contening the type of each index, COV for a covariant one and CON for a contravariant one.

Definition at line 321 of file tenseur.h.

◆ valence

int Lorene::Tenseur::valence
protected

Valence.

Definition at line 310 of file tenseur.h.


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