LORENE
|
Class intended to describe valence-2 symmetric tensors. More...
#include <sym_tensor.h>
Public Member Functions | |
Sym_tensor (const Map &map, const Itbl &tipe, const Base_vect &triad_i) | |
Standard constructor. More... | |
Sym_tensor (const Map &map, int tipe, const Base_vect &triad_i) | |
Standard constructor when both indices are of the same type. More... | |
Sym_tensor (const Sym_tensor &a) | |
Copy constructor. More... | |
Sym_tensor (const Tensor &a) | |
Constructor from a Tensor . More... | |
Sym_tensor (const Map &map, const Base_vect &triad_i, FILE *fich) | |
Constructor from a file (see sauve(FILE*) ). More... | |
virtual | ~Sym_tensor () |
Destructor. More... | |
virtual void | operator= (const Sym_tensor &a) |
Assignment to another Sym_tensor . More... | |
virtual void | operator= (const Tensor_sym &a) |
Assignment to a Tensor_sym . More... | |
virtual void | operator= (const Tensor &a) |
Assignment to a Tensor . More... | |
void | set_longit_trans (const Vector &v, const Sym_tensor_trans &a) |
Assigns the derived members p_longit_pot and p_transverse and updates the components accordingly. More... | |
void | set_auxiliary (const Scalar &trr, const Scalar &eta_over_r, const Scalar &mu_over_r, const Scalar &www, const Scalar &xxx, const Scalar &ttt) |
Assigns the component and the derived members p_eta , p_mu , p_www , p_xxx and p_ttt , fro, their values and , . 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... | |
const Vector & | divergence (const Metric &) const |
Returns the divergence of this with respect to a Metric . More... | |
Sym_tensor | derive_lie (const Vector &v) const |
Computes the Lie derivative of this with respect to some vector field v . More... | |
const Sym_tensor_trans & | transverse (const Metric &gam, Param *par=0x0, int method_poisson=6) const |
Computes the transverse part of the tensor with respect to a given metric, transverse meaning divergence-free with respect to that metric. More... | |
const Vector & | longit_pot (const Metric &gam, Param *par=0x0, int method_poisson=6) const |
Computes the vector potential of longitudinal part of the tensor (see documentation of method transverse() above). More... | |
virtual const Scalar & | eta (Param *par=0x0) const |
Gives the field (see member p_eta ). More... | |
const Scalar & | mu (Param *par=0x0) const |
Gives the field (see member p_mu ). More... | |
const Scalar & | www () const |
Gives the field W (see member p_www ). More... | |
const Scalar & | xxx () const |
Gives the field X (see member p_xxx ). More... | |
const Scalar & | ttt () const |
Gives the field T (see member p_ttt ). More... | |
const Scalar & | compute_A (bool output_ylm=true, Param *par=0x0) const |
Gives the field A (see member p_aaa ). More... | |
const Scalar & | compute_tilde_B (bool output_ylm=true, Param *par=0x0) const |
Gives the field (see member p_tilde_b ). More... | |
Scalar | compute_tilde_B_tt (bool output_ylm=true, Param *par=0x0) const |
Gives the field (see member p_tilde_b ) associated with the TT-part of the Sym_tensor . More... | |
const Scalar & | compute_tilde_C (bool output_ylm=true, Param *par=0x0) const |
Gives the field (see member p_tilde_c ). More... | |
int | sym_index1 () const |
Number of the first symmetric index (0<= id_sym1 < valence ) More... | |
int | sym_index2 () const |
Number of the second symmetric index (id_sym1 < id_sym2 < valence ) 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... | |
virtual void | sauve (FILE *) const |
Save in a binary file. More... | |
const Tensor_sym & | derive_cov (const Metric &gam) const |
Returns the covariant derivative of this with respect to some metric . More... | |
const Tensor_sym & | derive_con (const Metric &gam) const |
Returns the "contravariant" derivative of this with respect to some metric , by raising the last index of the covariant derivative (cf. 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... | |
Scalar & | set (const Itbl &ind) |
Returns the value of a component (read/write version). More... | |
Scalar & | set (int i1, int i2) |
Returns the value of a component for a tensor of valence 2 (read/write version). More... | |
Scalar & | set (int i1, int i2, int i3) |
Returns the value of a component for a tensor of valence 3 (read/write version). More... | |
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). 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_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... | |
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... | |
const Map & | get_mp () const |
Returns the mapping. More... | |
const Base_vect * | get_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... | |
Itbl & | set_index_type () |
Sets the types of all the indices. More... | |
const Scalar & | operator() (const Itbl &ind) const |
Returns the value of a component (read-only version). More... | |
const Scalar & | operator() (int i1, int i2) const |
Returns the value of a component for a tensor of valence 2 (read-only version). More... | |
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). More... | |
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). More... | |
void | operator+= (const Tensor &) |
+= Tensor More... | |
void | operator-= (const Tensor &) |
-= Tensor 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 | |
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 i) const |
Logical destructor of the derivatives depending on the i-th element of met_depend specific to the class Sym_tensor (p_transverse , etc...). More... | |
void | set_der_met_0x0 (int i) const |
Sets all the i-th components of met_depend specific to the class Sym_tensor (p_transverse , etc...) to 0x0. More... | |
Scalar | get_tilde_B_from_TT_trace (const Scalar &tilde_B_tt_in, const Scalar &trace) const |
Computes (see Sym_tensor::p_tilde_b ) from its transverse-traceless part and the trace. More... | |
Sym_tensor * | inverse () const |
Returns a pointer on the inverse of the Sym_tensor (seen as a matrix). 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 | |
Sym_tensor_trans * | p_transverse [N_MET_MAX] |
Array of the transverse part of the tensor with respect to various metrics, transverse meaning divergence-free with respect to a metric. More... | |
Vector * | p_longit_pot [N_MET_MAX] |
Array of the vector potential of the longitudinal part of the tensor with respect to various metrics (see documentation of member p_transverse . More... | |
Scalar * | p_eta |
Field such that the components of the tensor are written (has only meaning with spherical components!):
. More... | |
Scalar * | p_mu |
Field such that the components of the tensor are written (has only meaning with spherical components!):
. More... | |
Scalar * | p_www |
Field W such that the components and of the tensor are written (has only meaning with spherical components!):
. More... | |
Scalar * | p_xxx |
Field X such that the components and of the tensor are written (has only meaning with spherical components!):
. More... | |
Scalar * | p_ttt |
Field T defined as . More... | |
Scalar * | p_aaa |
Field A defined from X and insensitive to the longitudinal part of the Sym_tensor (only for ). More... | |
Scalar * | p_tilde_b |
Field defined from and h insensitive to the longitudinal part of the Sym_tensor . More... | |
Scalar * | p_tilde_c |
Field defined from and h insensitive to the longitudinal part of the Sym_tensor . More... | |
int | id_sym1 |
Number of the first symmetric index (0<= id_sym1 < valence ) More... | |
int | id_sym2 |
Number of the second symmetric index (id_sym1 < id_sym2 < valence ) More... | |
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_vect * | triad |
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 Metric * | met_depend [N_MET_MAX] |
Array on the Metric 's which were used to compute derived quantities, like p_derive_cov , etc... More... | |
Tensor * | p_derive_cov [N_MET_MAX] |
Array of pointers on the covariant derivatives of this with respect to various metrics. More... | |
Tensor * | p_derive_con [N_MET_MAX] |
Array of pointers on the contravariant derivatives of this with respect to various metrics. More... | |
Tensor * | p_divergence [N_MET_MAX] |
Array of pointers on the divergence of this with respect to various metrics. More... | |
Friends | |
class | Metric |
Class intended to describe valence-2 symmetric tensors.
The storage and the calculations are different and quicker than with an usual Tensor
.
The valence must be 2.()
Definition at line 226 of file sym_tensor.h.
Standard constructor.
map | the mapping |
tipe | 1-D array of integers (class Itbl ) of size 2 containing the type of each index, COV for a covariant one and CON for a contravariant one, with the following storage convention:
|
triad_i | vectorial basis (triad) with respect to which the tensor components are defined |
Definition at line 145 of file sym_tensor.C.
References set_der_0x0().
Standard constructor when both indices are of the same type.
map | the mapping |
tipe | the type of the indices. |
triad_i | vectorial basis (triad) with respect to which the tensor components are defined |
Definition at line 155 of file sym_tensor.C.
References set_der_0x0().
Lorene::Sym_tensor::Sym_tensor | ( | const Sym_tensor & | a | ) |
Copy constructor.
Definition at line 163 of file sym_tensor.C.
References Lorene::Tensor::get_place_met(), Lorene::Tensor::met_depend, p_eta, p_longit_pot, p_mu, p_transverse, p_www, p_xxx, Lorene::Tensor::set_dependance(), and set_der_0x0().
Lorene::Sym_tensor::Sym_tensor | ( | const Tensor & | a | ) |
Constructor from a Tensor
.
The symmetry of the input tensor is assumed but is not checked.
Definition at line 196 of file sym_tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor_sym::indices(), Lorene::Tensor::n_comp, Lorene::Tensor::position(), set_der_0x0(), and Lorene::Tensor::valence.
Constructor from a file (see 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. |
fich | file which has been used by the function sauve(FILE*) . |
Definition at line 213 of file sym_tensor.C.
References Lorene::Tensor::n_comp, set_der_0x0(), and Lorene::Tensor::valence.
|
virtual |
|
virtualinherited |
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(), Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), and Lorene::Tensor::n_comp.
|
virtualinherited |
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 Lorene::Scalar.
Definition at line 680 of file tensor.C.
References Lorene::Scalar::annule(), Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::mp, Lorene::Tensor::n_comp, and Lorene::Tensor::set_etat_zero().
|
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 675 of file tensor.C.
References Lorene::Tensor::annule().
|
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 699 of file tensor.C.
References Lorene::Map::get_mg(), Lorene::Mg3d::get_type_r(), and Lorene::Tensor::mp.
|
virtualinherited |
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Reimplemented in Lorene::Scalar, and Lorene::Vector.
Definition at line 80 of file tensor_change_triad.C.
References Lorene::Map::comp_p_from_cartesian(), Lorene::Map::comp_r_from_cartesian(), Lorene::Map::comp_t_from_cartesian(), Lorene::Map::comp_x_from_spherical(), Lorene::Map::comp_y_from_spherical(), Lorene::Map::comp_z_from_spherical(), Lorene::Base_vect_cart::get_align(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::mp, Lorene::Tensor::set(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Lorene::Tensor::valence.
Gives the field A (see member p_aaa
).
output_ylm | a flag to control the spectral decomposition base of the result: if true (default) the spherical harmonics base is used. |
Definition at line 319 of file sym_tensor_aux.C.
References Lorene::Scalar::annule_l(), Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::div_tant(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Tensor::mp, Lorene::Tensor::operator()(), p_aaa, Lorene::Map::poisson_angu(), Lorene::Scalar::poisson_angu(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, xxx(), Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
|
protectedinherited |
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 Lorene::Tensor::cmp, Lorene::contract(), Lorene::Scalar::dec_dzpuis(), Lorene::Tensor::derive_cov(), Lorene::Map::flat_met_cart(), Lorene::Map::flat_met_spher(), Lorene::Scalar::get_dzpuis(), Lorene::Tensor::get_n_comp(), Lorene::Tensor::get_triad(), Lorene::Tensor::indices(), Lorene::Tensor::mp, Lorene::Tensor::n_comp, Lorene::Tensor::operator()(), Lorene::Itbl::set(), Lorene::Tensor::set(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Lorene::Tensor::valence.
const Scalar & Lorene::Sym_tensor::compute_tilde_B | ( | bool | output_ylm = true , |
Param * | par = 0x0 |
||
) | const |
Gives the field (see member p_tilde_b
).
output_ylm | a flag to control the spectral decomposition base of the result: if true (default) the spherical harmonics base is used. |
Definition at line 365 of file sym_tensor_aux.C.
References Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::div_tant(), Lorene::Scalar::dsdr(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::mp, Lorene::Tensor::operator()(), p_tilde_b, Lorene::Map::poisson_angu(), Lorene::Scalar::poisson_angu(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, ttt(), www(), and Lorene::Valeur::ylm().
Gives the field (see member p_tilde_b
) associated with the TT-part of the Sym_tensor
.
output_ylm | a flag to control the spectral decomposition base of the result: if true (default) the spherical harmonics base is used. |
Definition at line 481 of file sym_tensor_aux.C.
References compute_tilde_B(), and Lorene::Scalar::get_etat().
const Scalar & Lorene::Sym_tensor::compute_tilde_C | ( | bool | output_ylm = true , |
Param * | par = 0x0 |
||
) | const |
Gives the field (see member p_tilde_c
).
output_ylm | a flag to control the spectral decomposition base of the result: if true (default) the spherical harmonics base is used. |
Definition at line 599 of file sym_tensor_aux.C.
References Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::div_tant(), Lorene::Scalar::dsdr(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::mp, Lorene::Tensor::operator()(), p_tilde_c, Lorene::Map::poisson_angu(), Lorene::Scalar::poisson_angu(), Lorene::Scalar::set_etat_qcq(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, ttt(), www(), and Lorene::Valeur::ylm().
|
virtualinherited |
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 Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), and Lorene::Tensor::n_comp.
|
protectedvirtual |
Deletes the derived quantities.
Reimplemented from Lorene::Tensor.
Reimplemented in Lorene::Sym_tensor_tt, and Lorene::Sym_tensor_trans.
Definition at line 289 of file sym_tensor.C.
References Lorene::Tensor::del_deriv(), del_derive_met(), p_aaa, p_eta, p_mu, p_tilde_b, p_tilde_c, p_ttt, p_www, p_xxx, and set_der_0x0().
|
protectedvirtual |
Logical destructor of the derivatives depending on the i-th element of met_depend
specific to the class Sym_tensor
(p_transverse
, etc...).
Reimplemented from Lorene::Tensor.
Definition at line 323 of file sym_tensor.C.
References Lorene::Tensor::del_derive_met(), Lorene::Tensor::met_depend, p_longit_pot, p_transverse, and set_der_met_0x0().
|
inherited |
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 .
Definition at line 207 of file tensor_sym_calculus.C.
References Lorene::Tensor::derive_con().
|
inherited |
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).
gam | metric |
this
with respect to the connection associated with the metric Definition at line 195 of file tensor_sym_calculus.C.
References Lorene::Tensor::derive_cov().
Sym_tensor Lorene::Sym_tensor::derive_lie | ( | const Vector & | v | ) | const |
Computes the Lie derivative of this
with respect to some vector field v
.
Definition at line 363 of file sym_tensor.C.
References Lorene::Tensor::compute_derive_lie(), Lorene::Tensor::mp, Lorene::Tensor::triad, and Lorene::Tensor::type_indice.
Returns the divergence of this
with respect to a Metric
.
The indices are assumed to be contravariant.
Definition at line 352 of file sym_tensor.C.
References Lorene::Tensor::divergence().
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 268 of file tensor_calculus.C.
References Lorene::contract(), Lorene::Metric::cov(), Lorene::Tensor::indices(), Lorene::Tensor::mp, Lorene::Tensor::n_comp, Lorene::Itbl::set(), Lorene::Tensor::set(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Lorene::Tensor::valence.
Gives the field (see member p_eta
).
Reimplemented in Lorene::Sym_tensor_tt.
Definition at line 114 of file sym_tensor_aux.C.
References Lorene::Scalar::div_tant(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Tensor::mp, Lorene::Scalar::mult_r_dzpuis(), Lorene::Tensor::operator()(), p_eta, Lorene::Map::poisson_angu(), Lorene::Scalar::poisson_angu(), Lorene::Scalar::stdsdp(), and Lorene::Tensor::triad.
|
virtual |
Applies exponential filters to all components (see Scalar::exponential_filter_r
).
Does a loop for Cartesian components, and works in terms of the rr-component, , , W
, X
, T
for spherical components.
Reimplemented from Lorene::Tensor.
Definition at line 449 of file sym_tensor.C.
References Lorene::Tensor::cmp, Lorene::Scalar::div_r(), eta(), Lorene::Scalar::exponential_filter_r(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Base_vect::identify(), Lorene::Tensor::mp, mu(), Lorene::Tensor::n_comp, Lorene::Tensor::operator()(), set_auxiliary(), Lorene::Tensor::triad, ttt(), www(), and xxx().
|
virtual |
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, , , W
, X
, T
for spherical components.
Reimplemented from Lorene::Tensor.
Definition at line 474 of file sym_tensor.C.
References Lorene::Tensor::cmp, Lorene::Scalar::div_r(), eta(), Lorene::Scalar::exponential_filter_ylm(), Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Base_vect::identify(), Lorene::Tensor::mp, mu(), Lorene::Tensor::n_comp, Lorene::Tensor::operator()(), set_auxiliary(), Lorene::Tensor::triad, ttt(), www(), and xxx().
|
virtualinherited |
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 Lorene::Tensor::cmp, Lorene::Map::get_bvect_cart(), Lorene::Base_vect::identify(), Lorene::Tensor::mp, Lorene::Tensor::n_comp, and Lorene::Tensor::triad.
|
inlineinherited |
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 899 of file tensor.h.
References Lorene::Tensor::type_indice.
|
inlineinherited |
Returns the types of all the indices.
Itbl
) of size valence
COV
for a covariant one and CON
Definition at line 909 of file tensor.h.
References Lorene::Tensor::type_indice.
|
inlineinherited |
|
inlineinherited |
Returns the number of stored components.
Definition at line 885 of file tensor.h.
References Lorene::Tensor::n_comp.
|
protectedinherited |
Returns the position of the pointer on metre
in the array met_depend
.
Definition at line 452 of file tensor.C.
References Lorene::Tensor::met_depend.
|
protected |
Computes (see Sym_tensor::p_tilde_b
) from its transverse-traceless part and the trace.
Definition at line 534 of file sym_tensor_aux.C.
References Lorene::Scalar::get_etat(), and Lorene::Tensor::triad.
|
inlineinherited |
Returns the vectorial basis (triad) on which the components are defined.
Definition at line 879 of file tensor.h.
References Lorene::Tensor::triad.
|
inlineinherited |
|
virtualinherited |
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 Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), and Lorene::Tensor::n_comp.
|
virtualinherited |
Returns the indices of a component given by its position in the array cmp
.
pos | [input] position in the array cmp of the pointer to the Scalar representing a component |
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) Reimplemented from Lorene::Tensor.
Definition at line 313 of file tensor_sym.C.
References Lorene::Tensor_sym::id_sym1, Lorene::Tensor_sym::id_sym2, Lorene::Tensor::n_comp, Lorene::Itbl::set(), and Lorene::Tensor::valence.
|
protected |
Returns a pointer on the inverse of the Sym_tensor
(seen as a matrix).
Definition at line 375 of file sym_tensor.C.
References Lorene::Tensor::mp, Lorene::Tensor::operator()(), Lorene::Tensor::set(), Sym_tensor(), Lorene::Tensor::triad, and Lorene::Tensor::type_indice.
const Vector & Lorene::Sym_tensor::longit_pot | ( | const Metric & | gam, |
Param * | par = 0x0 , |
||
int | method_poisson = 6 |
||
) | const |
Computes the vector potential of longitudinal part of the tensor (see documentation of method transverse()
above).
gam | metric with respect to the transverse decomposition is performed |
par | parameters for the vector Poisson equation |
method_poisson | type of method for solving the vector Poisson equation to get the longitudinal part (see method Vector::poisson ) |
Definition at line 149 of file sym_tensor_decomp.C.
References Lorene::Tensor::dec_dzpuis(), Lorene::Tensor::derive_con(), Lorene::Tensor_sym::derive_con(), Lorene::diffrel(), divergence(), Lorene::Tensor::divergence(), Lorene::Map::get_mg(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::get_place_met(), Lorene::maxabs(), Lorene::Tensor::mp, p_longit_pot, Lorene::Vector::poisson(), and Lorene::Tensor::set_dependance().
Gives the field (see member p_mu
).
Definition at line 154 of file sym_tensor_aux.C.
References Lorene::Scalar::div_tant(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Tensor::mp, Lorene::Scalar::mult_r_dzpuis(), Lorene::Tensor::operator()(), p_mu, Lorene::Map::poisson_angu(), Lorene::Scalar::poisson_angu(), Lorene::Scalar::stdsdp(), and Lorene::Tensor::triad.
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 807 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), Lorene::Tensor::position(), and Lorene::Tensor::valence.
|
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 769 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::position(), Lorene::Itbl::set(), and Lorene::Tensor::valence.
|
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 780 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::position(), Lorene::Itbl::set(), and Lorene::Tensor::valence.
|
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 792 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::position(), Lorene::Itbl::set(), and Lorene::Tensor::valence.
|
inherited |
+= Tensor
Definition at line 580 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Tensor::indices(), Lorene::Tensor::n_comp, Lorene::Tensor::position(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Lorene::Tensor::valence.
|
inherited |
-= Tensor
Definition at line 596 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Tensor::indices(), Lorene::Tensor::n_comp, Lorene::Tensor::position(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Lorene::Tensor::valence.
|
virtual |
Assignment to another Sym_tensor
.
Reimplemented in Lorene::Sym_tensor_tt, and Lorene::Sym_tensor_trans.
Definition at line 237 of file sym_tensor.C.
References del_deriv(), Lorene::Tensor::get_place_met(), Lorene::Tensor::met_depend, Lorene::Tensor_sym::operator=(), p_eta, p_longit_pot, p_mu, p_transverse, p_www, p_xxx, and Lorene::Tensor::set_dependance().
|
virtual |
Assignment to a Tensor_sym
.
Reimplemented from Lorene::Tensor_sym.
Reimplemented in Lorene::Sym_tensor_tt, and Lorene::Sym_tensor_trans.
Definition at line 270 of file sym_tensor.C.
References del_deriv(), and Lorene::Tensor_sym::operator=().
|
virtual |
Assignment to a Tensor
.
The symmetry is assumed but not checked.
Reimplemented from Lorene::Tensor_sym.
Reimplemented in Lorene::Sym_tensor_tt, and Lorene::Sym_tensor_trans.
Definition at line 278 of file sym_tensor.C.
References del_deriv(), and Lorene::Tensor_sym::operator=().
|
virtualinherited |
Returns the position in the array cmp
of a component given by its indices.
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:
|
cmp
of the pointer to the Scalar
containing the component specified by ind
Reimplemented from Lorene::Tensor.
Definition at line 248 of file tensor_sym.C.
References Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), Lorene::Tensor_sym::id_sym1, Lorene::Tensor_sym::id_sym2, Lorene::Itbl::set(), and Lorene::Tensor::valence.
|
virtualinherited |
Save in a binary file.
Reimplemented from Lorene::Tensor.
Definition at line 375 of file tensor_sym.C.
References Lorene::fwrite_be(), Lorene::Tensor_sym::id_sym1, Lorene::Tensor_sym::id_sym2, and Lorene::Tensor::sauve().
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 663 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Itbl::get_dim(), Lorene::Itbl::get_ndim(), Lorene::Tensor::position(), and Lorene::Tensor::valence.
|
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 615 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Tensor::position(), Lorene::Itbl::set(), and Lorene::Tensor::valence.
|
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 630 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Tensor::position(), Lorene::Itbl::set(), and Lorene::Tensor::valence.
|
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 646 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Tensor::position(), Lorene::Itbl::set(), and Lorene::Tensor::valence.
void Lorene::Sym_tensor::set_auxiliary | ( | const Scalar & | trr, |
const Scalar & | eta_over_r, | ||
const Scalar & | mu_over_r, | ||
const Scalar & | www, | ||
const Scalar & | xxx, | ||
const Scalar & | ttt | ||
) |
Assigns the component and the derived members p_eta
, p_mu
, p_www
, p_xxx
and p_ttt
, fro, their values and , .
It updates the other components accordingly.
Definition at line 269 of file sym_tensor_aux.C.
References Lorene::Scalar::check_dzpuis(), del_deriv(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::lapang(), Lorene::Scalar::mult_r_dzpuis(), p_eta, p_mu, p_ttt, p_www, p_xxx, Lorene::Scalar::set_spectral_va(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, and Lorene::Valeur::ylm_i().
|
protectedinherited |
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 Lorene::Tensor::met_depend, and Lorene::Metric::tensor_depend.
|
protected |
|
protected |
Sets all the i-th components of met_depend
specific to the class Sym_tensor
(p_transverse
, etc...) to 0x0.
Definition at line 338 of file sym_tensor.C.
References p_longit_pot, and p_transverse.
|
virtualinherited |
Sets the logical state of all components to ETATNONDEF
(undefined state).
Reimplemented in Lorene::Scalar.
Definition at line 498 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Tensor::n_comp, and Lorene::Scalar::set_etat_nondef().
|
virtualinherited |
Sets the logical state of all components to ETATQCQ
(ordinary state).
Reimplemented in Lorene::Scalar.
Definition at line 490 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Tensor::n_comp, and Lorene::Scalar::set_etat_qcq().
|
virtualinherited |
Sets the logical state of all components to ETATZERO
(zero state).
Reimplemented in Lorene::Scalar.
Definition at line 506 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::del_deriv(), Lorene::Tensor::n_comp, and Lorene::Scalar::set_etat_zero().
|
inlineinherited |
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 922 of file tensor.h.
References Lorene::Itbl::set(), and Lorene::Tensor::type_indice.
|
inlineinherited |
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 931 of file tensor.h.
References Lorene::Tensor::type_indice.
void Lorene::Sym_tensor::set_longit_trans | ( | const Vector & | v, |
const Sym_tensor_trans & | a | ||
) |
Assigns the derived members p_longit_pot
and p_transverse
and updates the components accordingly.
(see the documentation of these derived members for details)
Definition at line 94 of file sym_tensor_decomp.C.
References Lorene::Tensor::dec_dzpuis(), del_deriv(), Lorene::Tensor::get_index_type(), Lorene::Sym_tensor_trans::get_met_div(), Lorene::Tensor::get_place_met(), Lorene::Vector::ope_killing(), p_longit_pot, p_transverse, and Lorene::Tensor::set_dependance().
|
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 528 of file tensor.C.
References Lorene::Tensor::triad.
|
virtualinherited |
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 Lorene::Scalar.
Definition at line 883 of file tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor::indices(), Lorene::Tensor::n_comp, Lorene::Scalar::spectral_display(), and Lorene::Tensor::valence.
|
virtualinherited |
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 Lorene::Tensor::cmp, Lorene::Map::get_bvect_cart(), Lorene::Map::get_bvect_spher(), Lorene::Map::get_mg(), Lorene::Base_vect::identify(), Lorene::Tensor::indices(), Lorene::Tensor::mp, Lorene::Tensor::n_comp, Lorene::Scalar::set_spectral_base(), Lorene::Mg3d::std_base_vect_cart(), Lorene::Mg3d::std_base_vect_spher(), Lorene::Scalar::std_spectral_base(), Lorene::Tensor::triad, and Lorene::Tensor::valence.
|
virtualinherited |
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 Lorene::Tensor::cmp, Lorene::Scalar::std_spectral_base_odd(), and Lorene::Tensor::valence.
|
inlineinherited |
Number of the first symmetric index (0<=
id_sym1
< valence
)
Definition at line 1162 of file tensor.h.
References Lorene::Tensor_sym::id_sym1.
|
inlineinherited |
Number of the second symmetric index (id_sym1
< id_sym2
< valence
)
Definition at line 1167 of file tensor.h.
References Lorene::Tensor_sym::id_sym2.
|
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 97 of file tensor_calculus.C.
References Lorene::Tensor::cmp, Lorene::Tensor::get_n_comp(), Lorene::Tensor::indices(), Lorene::Tensor::mp, Lorene::Tensor::position(), Lorene::Itbl::set(), Lorene::Tensor::set(), Lorene::Scalar::set_etat_zero(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Lorene::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 156 of file tensor_calculus.C.
References Lorene::Metric::con(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Tensor::trace(), Lorene::Tensor::type_indice, and Lorene::Tensor::valence.
|
inherited |
Trace on two different type indices for a valence 2 tensor.
Definition at line 183 of file tensor_calculus.C.
References Lorene::Tensor::mp, Lorene::Tensor::operator()(), Lorene::Scalar::set_etat_zero(), Lorene::Tensor::type_indice, and Lorene::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 200 of file tensor_calculus.C.
References Lorene::Metric::con(), Lorene::contract(), Lorene::Metric::cov(), Lorene::Tensor::trace(), Lorene::Tensor::type_indice, and Lorene::Tensor::valence.
const Sym_tensor_trans & Lorene::Sym_tensor::transverse | ( | const Metric & | gam, |
Param * | par = 0x0 , |
||
int | method_poisson = 6 |
||
) | const |
Computes the transverse part of the tensor with respect to a given metric, transverse meaning divergence-free with respect to that metric.
Denoting *this
by , we then have
where denotes the covariant derivative with respect to the given metric and is the vector potential of the longitudinal part of (function longit_pot()
below)
gam | metric with respect to the transverse decomposition is performed |
par | parameters for the vector Poisson equation |
method_poisson | type of method for solving the vector Poisson equation to get the longitudinal part (see method Vector::poisson ) |
Definition at line 116 of file sym_tensor_decomp.C.
References Lorene::Tensor::cmp, Lorene::Tensor::get_place_met(), Lorene::Tensor::inc_dzpuis(), longit_pot(), Lorene::Tensor::mp, Lorene::Tensor::n_comp, Lorene::Vector::ope_killing(), p_transverse, Lorene::Tensor::set_dependance(), Lorene::Tensor::triad, and Lorene::Tensor::type_indice.
const Scalar & Lorene::Sym_tensor::ttt | ( | ) | const |
Gives the field T (see member p_ttt
).
Definition at line 193 of file sym_tensor_aux.C.
References p_ttt, and Lorene::Tensor::triad.
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 228 of file tensor_calculus.C.
References Lorene::Metric::con(), Lorene::contract(), Lorene::Tensor::indices(), Lorene::Tensor::mp, Lorene::Tensor::n_comp, Lorene::Itbl::set(), Lorene::Tensor::set(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Lorene::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 308 of file tensor_calculus.C.
References Lorene::Tensor::down(), Lorene::Tensor::Tensor(), Lorene::Tensor::type_indice, Lorene::Tensor::up(), and Lorene::Tensor::valence.
const Scalar & Lorene::Sym_tensor::www | ( | ) | const |
Gives the field W (see member p_www
).
Definition at line 212 of file sym_tensor_aux.C.
References Lorene::Scalar::div_tant(), Lorene::Scalar::dsdt(), Lorene::Tensor::operator()(), p_www, Lorene::Scalar::poisson_angu(), Lorene::Scalar::stdsdp(), and Lorene::Tensor::triad.
const Scalar & Lorene::Sym_tensor::xxx | ( | ) | const |
Gives the field X (see member p_xxx
).
Definition at line 243 of file sym_tensor_aux.C.
References Lorene::Scalar::div_tant(), Lorene::Scalar::dsdt(), Lorene::Tensor::operator()(), p_xxx, Lorene::Scalar::poisson_angu(), Lorene::Scalar::stdsdp(), and Lorene::Tensor::triad.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotectedinherited |
|
protectedinherited |
|
protectedinherited |
|
mutableprotected |
Field A defined from X and insensitive to the longitudinal part of the Sym_tensor
(only for ).
Its definition reads
Definition at line 325 of file sym_tensor.h.
|
mutableprotectedinherited |
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.
|
mutableprotectedinherited |
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.
|
mutableprotectedinherited |
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.
|
mutableprotected |
Field such that the components of the tensor are written (has only meaning with spherical components!):
.
Definition at line 263 of file sym_tensor.h.
|
mutableprotected |
Array of the vector potential of the longitudinal part of the tensor with respect to various metrics (see documentation of member p_transverse
.
Definition at line 249 of file sym_tensor.h.
|
mutableprotected |
Field such that the components of the tensor are written (has only meaning with spherical components!):
.
Definition at line 277 of file sym_tensor.h.
|
mutableprotected |
Field defined from and h insensitive to the longitudinal part of the Sym_tensor
.
It is defined for each multipolar momentum by
Definition at line 337 of file sym_tensor.h.
|
mutableprotected |
Field defined from and h insensitive to the longitudinal part of the Sym_tensor
.
It is defined for each multipolar momentum by
Definition at line 349 of file sym_tensor.h.
|
mutableprotected |
Array of the transverse part of the tensor with respect to various metrics, transverse meaning divergence-free with respect to a metric.
Denoting *this
by , we then have
where denotes the covariant derivative with respect to the given metric and is the vector potential of the longitudinal part of (member p_longit_pot
below)
Definition at line 242 of file sym_tensor.h.
|
mutableprotected |
Field T defined as .
Definition at line 318 of file sym_tensor.h.
|
mutableprotected |
Field W such that the components and of the tensor are written (has only meaning with spherical components!):
.
Definition at line 296 of file sym_tensor.h.
|
mutableprotected |
Field X such that the components and of the tensor are written (has only meaning with spherical components!):
.
Definition at line 315 of file sym_tensor.h.
|
protectedinherited |
|
protectedinherited |
|
protectedinherited |