Transverse symmetric tensors of rank 2. More...
#include <sym_tensor.h>
Public Member Functions | |
Sym_tensor_trans (const Map &map, const Base_vect &triad_i, const Metric &met) | |
Standard constructor. | |
Sym_tensor_trans (const Sym_tensor_trans &) | |
Copy constructor. | |
Sym_tensor_trans (const Map &map, const Base_vect &triad_i, const Metric &met, FILE *fich) | |
Constructor from a file (see Tensor::sauve(FILE*) ). | |
virtual | ~Sym_tensor_trans () |
Destructor. | |
const Metric & | get_met_div () const |
Returns the metric with respect to which the divergence and the trace are defined. | |
virtual void | operator= (const Sym_tensor_trans &a) |
Assignment to another Sym_tensor_trans . | |
virtual void | operator= (const Sym_tensor &a) |
Assignment to a Sym_tensor . | |
virtual void | operator= (const Tensor_sym &a) |
Assignment to a Tensor_sym . | |
virtual void | operator= (const Tensor &a) |
Assignment to a Tensor . | |
void | set_tt_trace (const Sym_tensor_tt &a, const Scalar &h, Param *par=0x0) |
Assigns the derived members p_tt and p_trace and updates the components accordingly. | |
const Scalar & | the_trace () const |
Returns the trace of the tensor with respect to metric *met_div . | |
const Sym_tensor_tt & | tt_part (Param *par=0x0) const |
Returns the transverse traceless part of the tensor, the trace being defined with respect to metric *met_div . | |
void | sol_Dirac_Abound (const Scalar &aaa, Scalar &tilde_mu, Scalar &x_new, Scalar bound_mu, const Param *par_bc) |
Same resolution as sol_Dirac_A, but with inner boundary conditions added. | |
void | sol_Dirac_A2 (const Scalar &aaa, Scalar &tilde_mu, Scalar &x_new, Scalar bound_mu, const Param *par_bc) |
Same resolution as sol_Dirac_Abound, but here the boundary conditions are the degenerate elliptic conditions encontered when solving the Kerr problem. | |
void | sol_Dirac_BC2 (const Scalar &bb, const Scalar &cc, const Scalar &hh, Scalar &hrr, Scalar &tilde_eta, Scalar &ww, Scalar bound_eta, double dir, double neum, double rhor, Param *par_bc, Param *par_mat) |
Same resolution as sol_Dirac_tilde_B, but with inner boundary conditions added. | |
void | sol_Dirac_BC3 (const Scalar &bb, const Scalar &hh, Scalar &hrr, Scalar &tilde_eta, Scalar &ww, Scalar bound_hrr, Scalar bound_eta, Param *par_bc, Param *par_mat) |
Same resolution as sol_Dirac_Abound, but here the boundary conditions are the degenerate elliptic conditions encontered when solving the Kerr problem. | |
void | sol_Dirac_l01_bound (const Scalar &hh, Scalar &hrr, Scalar &tilde_eta, Scalar &bound_hrr, Scalar &bound_eta, Param *par_mat) |
void | sol_Dirac_l01_2 (const Scalar &hh, Scalar &hrr, Scalar &tilde_eta, Param *par_mat) |
void | sol_elliptic_ABC (Sym_tensor &source, Scalar aaa, Scalar bbb, Scalar ccc) |
Finds spectral potentials A, B, C of solution of an tensorial TT elliptic equation, given the source. | |
void | trace_from_det_one (const Sym_tensor_tt &htt, double precis=1.e-14, int it_max=100) |
Assigns the derived member p_tt and computes the trace so that *this + the flat metric has a determinant equal to 1. | |
void | set_hrr_mu_det_one (const Scalar &hrr, const Scalar &mu_in, double precis=1.e-14, int it_max=100) |
Assigns the rr component and the derived member ![]() | |
void | set_tt_part_det_one (const Sym_tensor_tt &hijtt, const Scalar *h_prev=0x0, Param *par_mat=0x0, double precis=1.e-14, int it_max=100) |
Assignes the TT-part of the tensor. | |
void | set_AtBtt_det_one (const Scalar &a_in, const Scalar &tbtt_in, const Scalar *h_prev=0x0, Param *par_bc=0x0, Param *par_mat=0x0, double precis=1.e-14, int it_max=100) |
Assigns the derived member A and computes ![]() Sym_tensor::compute_tilde_B_tt() ). | |
void | set_AtB_trace (const Scalar &a_in, const Scalar &tb_in, const Scalar &trace, Param *par_bc=0x0, Param *par_mat=0x0) |
Assigns the derived members A , ![]() | |
Sym_tensor_trans | poisson (const Scalar *h_guess=0x0) const |
Computes the solution of a tensorial transverse Poisson equation with *this ![]()
In particular, it makes an iteration on the trace of the result, using | |
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. | |
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 ![]() p_eta , p_mu , p_www , p_xxx and p_ttt , fro, their values and ![]() ![]() | |
virtual void | exponential_filter_r (int lzmin, int lzmax, int p, double alpha=-16.) |
Applies exponential filters to all components (see Scalar::exponential_filter_r ). | |
virtual void | exponential_filter_ylm (int lzmin, int lzmax, int p, double alpha=-16.) |
Applies exponential filters to all components (see Scalar::exponential_filter_ylm ). | |
const Vector & | divergence (const Metric &) const |
Returns the divergence of this with respect to a Metric . | |
Sym_tensor | derive_lie (const Vector &v) const |
Computes the Lie derivative of this with respect to some vector field v . | |
const Sym_tensor_trans & | transverse (const Metric &gam, Param *par=0x0, int method_poisson=6) const |
Computes the transverse part ![]() | |
const Vector & | longit_pot (const Metric &gam, Param *par=0x0, int method_poisson=6) const |
Computes the vector potential ![]() transverse() above). | |
virtual const Scalar & | eta (Param *par=0x0) const |
Gives the field ![]() p_eta ). | |
const Scalar & | mu (Param *par=0x0) const |
Gives the field ![]() p_mu ). | |
const Scalar & | www () const |
Gives the field W (see member p_www ). | |
const Scalar & | xxx () const |
Gives the field X (see member p_xxx ). | |
const Scalar & | ttt () const |
Gives the field T (see member p_ttt ). | |
const Scalar & | compute_A (bool output_ylm=true, Param *par=0x0) const |
Gives the field A (see member p_aaa ). | |
const Scalar & | compute_tilde_B (bool output_ylm=true, Param *par=0x0) const |
Gives the field ![]() p_tilde_b ). | |
Scalar | compute_tilde_B_tt (bool output_ylm=true, Param *par=0x0) const |
Gives the field ![]() p_tilde_b ) associated with the TT-part of the Sym_tensor . | |
const Scalar & | compute_tilde_C (bool output_ylm=true, Param *par=0x0) const |
Gives the field ![]() p_tilde_c ). | |
int | sym_index1 () const |
Number of the first symmetric index (0<= id_sym1 < valence ). | |
int | sym_index2 () const |
Number of the second symmetric index (id_sym1 < id_sym2 < valence ). | |
virtual int | position (const Itbl &ind) const |
Returns the position in the array cmp of a component given by its indices. | |
virtual Itbl | indices (int pos) const |
Returns the indices of a component given by its position in the array cmp . | |
virtual void | sauve (FILE *) const |
Save in a binary file. | |
const Tensor_sym & | derive_cov (const Metric &gam) const |
Returns the covariant derivative of this with respect to some metric ![]() | |
const Tensor_sym & | derive_con (const Metric &gam) const |
Returns the "contravariant" derivative of this with respect to some metric ![]() | |
virtual void | set_etat_nondef () |
Sets the logical state of all components to ETATNONDEF (undefined state). | |
virtual void | set_etat_zero () |
Sets the logical state of all components to ETATZERO (zero state). | |
virtual void | set_etat_qcq () |
Sets the logical state of all components to ETATQCQ (ordinary state). | |
virtual void | allocate_all () |
Performs the memory allocation of all the elements, down to the double arrays of the Tbl s. | |
virtual void | change_triad (const Base_vect &new_triad) |
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly. | |
void | set_triad (const Base_vect &new_triad) |
Assigns a new vectorial basis (triad) of decomposition. | |
Scalar & | set (const Itbl &ind) |
Returns the value of a component (read/write version). | |
Scalar & | set (int i1, int i2) |
Returns the value of a component for a tensor of valence 2 (read/write version). | |
Scalar & | set (int i1, int i2, int i3) |
Returns the value of a component for a tensor of valence 3 (read/write version). | |
Scalar & | set (int i1, int i2, int i3, int i4) |
Returns the value of a component for a tensor of valence 4 (read/write version). | |
void | annule_domain (int l) |
Sets the Tensor to zero in a given domain. | |
virtual void | annule (int l_min, int l_max) |
Sets the Tensor to zero in several domains. | |
void | annule_extern_cn (int l_0, int deg) |
Performs a smooth (C^n) transition in a given domain to zero. | |
virtual void | std_spectral_base () |
Sets the standard spectal bases of decomposition for each component. | |
virtual void | std_spectral_base_odd () |
Sets the standard odd spectal bases of decomposition for each component. | |
virtual void | dec_dzpuis (int dec=1) |
Decreases by dec units the value of dzpuis and changes accordingly the values in the compactified external domain (CED). | |
virtual void | inc_dzpuis (int inc=1) |
Increases by inc units the value of dzpuis and changes accordingly the values in the compactified external domain (CED). | |
Tensor | up (int ind, const Metric &gam) const |
Computes a new tensor by raising an index of *this . | |
Tensor | down (int ind, const Metric &gam) const |
Computes a new tensor by lowering an index of *this . | |
Tensor | up_down (const Metric &gam) const |
Computes a new tensor by raising or lowering all the indices of *this . | |
Tensor | trace (int ind1, int ind2) const |
Trace on two different type indices. | |
Tensor | trace (int ind1, int ind2, const Metric &gam) const |
Trace with respect to a given metric. | |
Scalar | trace () const |
Trace on two different type indices for a valence 2 tensor. | |
Scalar | trace (const Metric &gam) const |
Trace with respect to a given metric for a valence 2 tensor. | |
const Map & | get_mp () const |
Returns the mapping. | |
const Base_vect * | get_triad () const |
Returns the vectorial basis (triad) on which the components are defined. | |
int | get_valence () const |
Returns the valence. | |
int | get_n_comp () const |
Returns the number of stored components. | |
int | get_index_type (int i) const |
Gives the type (covariant or contravariant) of the index number i . | |
Itbl | get_index_type () const |
Returns the types of all the indices. | |
int & | set_index_type (int i) |
Sets the type of the index number i . | |
Itbl & | set_index_type () |
Sets the types of all the indices. | |
const Scalar & | operator() (const Itbl &ind) const |
Returns the value of a component (read-only version). | |
const Scalar & | operator() (int i1, int i2) const |
Returns the value of a component for a tensor of valence 2 (read-only version). | |
const Scalar & | operator() (int i1, int i2, int i3) const |
Returns the value of a component for a tensor of valence 3 (read-only version). | |
const Scalar & | operator() (int i1, int i2, int i3, int i4) const |
Returns the value of a component for a tensor of valence 4 (read-only version). | |
void | operator+= (const Tensor &) |
+= Tensor | |
void | operator-= (const Tensor &) |
-= Tensor | |
virtual void | spectral_display (const char *comment=0x0, double threshold=1.e-7, int precision=4, ostream &ostr=cout) const |
Displays the spectral coefficients and the associated basis functions of each component. | |
Protected Member Functions | |
virtual void | del_deriv () const |
Deletes the derived quantities. | |
void | set_der_0x0 () const |
Sets the pointers on derived quantities to 0x0. | |
void | sol_Dirac_A (const Scalar &aaa, Scalar &tilde_mu, Scalar &xxx, const Param *par_bc=0x0) const |
Solves a system of two coupled first-order PDEs obtained from the divergence-free condition (Dirac gauge) and the requirement that the potential A (see Sym_tensor::p_aaa ) has a given value. | |
void | sol_Dirac_tilde_B (const Scalar &tilde_b, const Scalar &hh, Scalar &hrr, Scalar &tilde_eta, Scalar &www, Param *par_bc=0x0, Param *par_mat=0x0) const |
Solves a system of three coupled first-order PDEs obtained from divergence-free conditions (Dirac gauge) and the requirement that the potential ![]() Sym_tensor::p_tilde_b ) has a given value. | |
void | sol_Dirac_l01 (const Scalar &hh, Scalar &hrr, Scalar &tilde_eta, Param *par_mat) const |
Solves the same system as Sym_tensor_trans::sol_Dirac_tilde_B but only for ![]() | |
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. | |
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. | |
Scalar | get_tilde_B_from_TT_trace (const Scalar &tilde_B_tt_in, const Scalar &trace) const |
Computes ![]() Sym_tensor::p_tilde_b ) from its transverse-traceless part and the trace. | |
Sym_tensor * | inverse () const |
Returns a pointer on the inverse of the Sym_tensor (seen as a matrix). | |
void | set_dependance (const Metric &) const |
To be used to describe the fact that the derivatives members have been calculated with met . | |
int | get_place_met (const Metric &) const |
Returns the position of the pointer on metre in the array met_depend . | |
void | compute_derive_lie (const Vector &v, Tensor &resu) const |
Computes the Lie derivative of this with respect to some vector field v (protected method; the public interface is method derive_lie ). | |
Protected Attributes | |
const Metric *const | met_div |
Metric with respect to which the divergence and the trace are defined. | |
Scalar * | p_trace |
Trace with respect to the metric *met_div . | |
Sym_tensor_tt * | p_tt |
Traceless part with respect to the metric *met_div . | |
Sym_tensor_trans * | p_transverse [N_MET_MAX] |
Array of the transverse part ![]() | |
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 . | |
Scalar * | p_eta |
Field ![]() ![]()
| |
Scalar * | p_mu |
Field ![]() ![]()
| |
Scalar * | p_www |
Field W such that the components ![]() ![]()
| |
Scalar * | p_xxx |
Field X such that the components ![]() ![]()
| |
Scalar * | p_ttt |
Field T defined as ![]() | |
Scalar * | p_aaa |
Field A defined from X and ![]() Sym_tensor (only for ![]() | |
Scalar * | p_tilde_b |
Field ![]() ![]() Sym_tensor . | |
Scalar * | p_tilde_c |
Field ![]() ![]() Sym_tensor . | |
int | id_sym1 |
Number of the first symmetric index (0<= id_sym1 < valence ). | |
int | id_sym2 |
Number of the second symmetric index (id_sym1 < id_sym2 < valence ). | |
const Map *const | mp |
Mapping on which the numerical values at the grid points are defined. | |
int | valence |
Valence of the tensor (0 = scalar, 1 = vector, etc...). | |
const Base_vect * | triad |
Vectorial basis (triad) with respect to which the tensor components are defined. | |
Itbl | type_indice |
1D array of integers (class Itbl ) of size valence containing the type of each index: COV for a covariant one and CON for a contravariant one. | |
int | n_comp |
Number of stored components, depending on the symmetry. | |
Scalar ** | cmp |
Array of size n_comp of pointers onto the components. | |
const Metric * | met_depend [N_MET_MAX] |
Array on the Metric 's which were used to compute derived quantities, like p_derive_cov , etc. | |
Tensor * | p_derive_cov [N_MET_MAX] |
Array of pointers on the covariant derivatives of this with respect to various metrics. | |
Tensor * | p_derive_con [N_MET_MAX] |
Array of pointers on the contravariant derivatives of this with respect to various metrics. | |
Tensor * | p_divergence [N_MET_MAX] |
Array of pointers on the divergence of this with respect to various metrics. | |
Friends | |
class | Metric |
class | Sym_tensor |
class | Tensor_sym |
Tensor_sym | operator* (const Tensor &, const Tensor_sym &) |
Tensor_sym | operator* (const Tensor_sym &, const Tensor &) |
Tensor | operator* (const Tensor &, const Tensor &) |
Tensor_sym | operator* (const Tensor &, const Tensor_sym &) |
Tensor_sym | operator* (const Tensor_sym &, const Tensor &) |
Tensor_sym | operator* (const Tensor_sym &, const Tensor_sym &) |
Tensorial product of two symmetric tensors. | |
class | Scalar |
class | Vector |
ostream & | operator<< (ostream &, const Tensor &) |
Scalar | operator+ (const Tensor &, const Scalar &) |
Scalar | operator+ (const Scalar &, const Tensor &) |
Scalar | operator- (const Tensor &, const Scalar &) |
Scalar | operator- (const Scalar &, const Tensor &) |
Transverse symmetric tensors of rank 2.
()
This class is designed to store transverse (divergence-free) symmetric contravariant tensors of rank 2, with the component expressed in an orthonormal spherical basis .
Definition at line 604 of file sym_tensor.h.
Sym_tensor_trans::Sym_tensor_trans | ( | const Map & | map, | |
const Base_vect & | triad_i, | |||
const Metric & | met | |||
) |
Standard constructor.
map | the mapping | |
triad_i | vectorial basis (triad) with respect to which the tensor components are defined | |
met | the metric with respect to which the divergence is defined |
Definition at line 110 of file sym_tensor_trans.C.
References set_der_0x0().
Sym_tensor_trans::Sym_tensor_trans | ( | const Sym_tensor_trans & | source | ) |
Copy constructor.
Definition at line 121 of file sym_tensor_trans.C.
References p_trace, p_tt, and set_der_0x0().
Sym_tensor_trans::Sym_tensor_trans | ( | const Map & | map, | |
const Base_vect & | triad_i, | |||
const Metric & | met, | |||
FILE * | fich | |||
) |
Constructor from a file (see Tensor::sauve(FILE*)
).
map | the mapping | |
triad_i | vectorial basis (triad) with respect to which the tensor components are defined. It will be checked that it coincides with the basis saved in the file. | |
met | the metric with respect to which the divergence is defined | |
fich | file which has been used by the function sauve(FILE*) . |
Definition at line 135 of file sym_tensor_trans.C.
References set_der_0x0().
Sym_tensor_trans::~Sym_tensor_trans | ( | ) | [virtual] |
void Tensor::allocate_all | ( | ) | [virtual, inherited] |
Performs the memory allocation of all the elements, down to the double
arrays of the Tbl
s.
This function performs in fact recursive calls to set_etat_qcq()
on each element of the chain Scalar
-> Valeur
-> Mtbl
-> Tbl
.
Reimplemented in Scalar.
Definition at line 504 of file tensor.C.
References Scalar::allocate_all(), Tensor::cmp, Tensor::del_deriv(), and Tensor::n_comp.
void Tensor::annule | ( | int | l_min, | |
int | l_max | |||
) | [virtual, inherited] |
Sets the Tensor
to zero in several domains.
l_min | [input] The Tensor will be set (logically) to zero in the domains whose indices are in the range [l_min,l_max] . | |
l_max | [input] see the comments for l_min . |
Note that annule
(0,nz-1) , where nz
is the total number of domains, is equivalent to set_etat_zero()
.
Reimplemented in Scalar.
Definition at line 667 of file tensor.C.
References Scalar::annule(), Tensor::cmp, Tensor::del_deriv(), Map::get_mg(), Mg3d::get_nzone(), Tensor::mp, Tensor::n_comp, and Tensor::set_etat_zero().
void Tensor::annule_domain | ( | int | l | ) | [inherited] |
Sets the Tensor
to zero in a given domain.
l | [input] Index of the domain in which the Tensor will be set (logically) to zero. |
Definition at line 662 of file tensor.C.
References Tensor::annule().
void Tensor::annule_extern_cn | ( | int | l_0, | |
int | deg | |||
) | [inherited] |
Performs a smooth (C^n) transition in a given domain to zero.
l_0 | [input] in the domain of index l0 the tensor is multiplied by the right polynomial (of degree 2n+1), to ensure continuty of the function and its n first derivative at both ends of this domain. The tensor is unchanged in the domains l < l_0 and set to zero in domains l > l_0. | |
deg | [input] the degree n of smoothness of the transition. |
Definition at line 686 of file tensor.C.
References Scalar::allocate_all(), Scalar::annule(), Itbl::annule_hard(), Tensor::cmp, Tensor::del_deriv(), Map::get_mg(), Mg3d::get_nr(), Mg3d::get_nzone(), Mg3d::get_type_r(), Tensor::mp, Tensor::n_comp, pow(), Map::r, Tbl::set(), Itbl::set(), Scalar::set_domain(), Tbl::set_etat_qcq(), Scalar::std_spectral_base(), and Map::val_r().
void Tensor::change_triad | ( | const Base_vect & | new_triad | ) | [virtual, inherited] |
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Reimplemented in Scalar, and Vector.
Definition at line 73 of file tensor_change_triad.C.
References Map::comp_p_from_cartesian(), Map::comp_r_from_cartesian(), Map::comp_t_from_cartesian(), Map::comp_x_from_spherical(), Map::comp_y_from_spherical(), Map::comp_z_from_spherical(), Base_vect_cart::get_align(), Map::get_bvect_cart(), Map::get_bvect_spher(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nt(), Mg3d::get_nzone(), Tensor::mp, Tensor::set(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
const Scalar & Sym_tensor::compute_A | ( | bool | output_ylm = true , |
|
Param * | par = 0x0 | |||
) | const [inherited] |
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 312 of file sym_tensor_aux.C.
References Scalar::annule_l(), Scalar::div_r_dzpuis(), Scalar::div_tant(), Scalar::dsdt(), Scalar::get_dzpuis(), Tensor::mp, Tensor::operator()(), Sym_tensor::p_aaa, Map::poisson_angu(), Scalar::poisson_angu(), Scalar::set_spectral_va(), Scalar::stdsdp(), Tensor::triad, Sym_tensor::xxx(), Valeur::ylm(), and Valeur::ylm_i().
Computes the Lie derivative of this
with respect to some vector field v
(protected method; the public interface is method derive_lie
).
Definition at line 335 of file tensor_calculus.C.
References Tensor::cmp, contract(), Scalar::dec_dzpuis(), Tensor::derive_cov(), Map::flat_met_cart(), Map::flat_met_spher(), Scalar::get_dzpuis(), Tensor::get_n_comp(), Tensor::get_triad(), Tensor::indices(), Tensor::mp, Tensor::n_comp, Tensor::operator()(), Tensor::set(), Itbl::set(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
const Scalar & Sym_tensor::compute_tilde_B | ( | bool | output_ylm = true , |
|
Param * | par = 0x0 | |||
) | const [inherited] |
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 358 of file sym_tensor_aux.C.
References Scalar::annule_hard(), Mtbl_cf::annule_hard(), Valeur::c_cf, Scalar::div_r_dzpuis(), Scalar::div_tant(), Scalar::dsdr(), Scalar::dsdt(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), Tensor::mp, Tensor::operator()(), Sym_tensor::p_tilde_b, Map::poisson_angu(), Scalar::poisson_angu(), Mtbl_cf::set(), Scalar::set_dzpuis(), Valeur::set_etat_cf_qcq(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Scalar::stdsdp(), Tensor::triad, Sym_tensor::ttt(), Sym_tensor::www(), Valeur::ylm(), and Valeur::ylm_i().
Scalar Sym_tensor::compute_tilde_B_tt | ( | bool | output_ylm = true , |
|
Param * | par = 0x0 | |||
) | const [inherited] |
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 474 of file sym_tensor_aux.C.
References Scalar::annule_hard(), Valeur::c, Valeur::c_cf, Sym_tensor::compute_tilde_B(), Scalar::div_r_dzpuis(), Scalar::dsdr(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), Tensor::mp, Tensor::operator()(), Mtbl_cf::set(), Scalar::set_dzpuis(), Scalar::set_spectral_va(), Sym_tensor::ttt(), Valeur::ylm(), and Valeur::ylm_i().
const Scalar & Sym_tensor::compute_tilde_C | ( | bool | output_ylm = true , |
|
Param * | par = 0x0 | |||
) | const [inherited] |
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 592 of file sym_tensor_aux.C.
References Scalar::annule_hard(), Mtbl_cf::annule_hard(), Valeur::c_cf, Scalar::div_r_dzpuis(), Scalar::div_tant(), Scalar::dsdr(), Scalar::dsdt(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), Tensor::mp, Tensor::operator()(), Sym_tensor::p_tilde_c, Map::poisson_angu(), Scalar::poisson_angu(), Mtbl_cf::set(), Scalar::set_dzpuis(), Valeur::set_etat_cf_qcq(), Scalar::set_etat_qcq(), Scalar::set_etat_zero(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Scalar::stdsdp(), Tensor::triad, Sym_tensor::ttt(), Sym_tensor::www(), Valeur::ylm(), and Valeur::ylm_i().
void Tensor::dec_dzpuis | ( | int | dec = 1 |
) | [virtual, inherited] |
Decreases by dec
units the value of dzpuis
and changes accordingly the values in the compactified external domain (CED).
Reimplemented in Scalar.
Definition at line 804 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), and Tensor::n_comp.
void Sym_tensor_trans::del_deriv | ( | ) | const [protected, virtual] |
Deletes the derived quantities.
Reimplemented from Sym_tensor.
Reimplemented in Sym_tensor_tt.
Definition at line 160 of file sym_tensor_trans.C.
References p_trace, p_tt, and set_der_0x0().
void Sym_tensor::del_derive_met | ( | int | i | ) | const [protected, virtual, inherited] |
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 Tensor.
Definition at line 316 of file sym_tensor.C.
References Tensor::met_depend, Sym_tensor::p_longit_pot, Sym_tensor::p_transverse, and Sym_tensor::set_der_met_0x0().
const Tensor_sym & Tensor_sym::derive_con | ( | const Metric & | gam | ) | const [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 .
Reimplemented from Tensor.
Definition at line 200 of file tensor_sym_calculus.C.
const Tensor_sym & Tensor_sym::derive_cov | ( | const Metric & | gam | ) | const [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 Reimplemented from Tensor.
Definition at line 188 of file tensor_sym_calculus.C.
Sym_tensor Sym_tensor::derive_lie | ( | const Vector & | v | ) | const [inherited] |
Computes the Lie derivative of this
with respect to some vector field v
.
Reimplemented from Tensor_sym.
Definition at line 356 of file sym_tensor.C.
References Tensor::compute_derive_lie(), Tensor::mp, Tensor::triad, and Tensor::type_indice.
Returns the divergence of this
with respect to a Metric
.
The indices are assumed to be contravariant.
Reimplemented from Tensor.
Definition at line 345 of file sym_tensor.C.
Computes a new tensor by lowering an index of *this
.
ind | index to be lowered, with the following convention :
| |
gam | metric used to lower the index (contraction with the twice covariant form of the metric on the index ind ). |
Definition at line 261 of file tensor_calculus.C.
References contract(), Metric::cov(), Tensor::indices(), Tensor::mp, Tensor::n_comp, Tensor::set(), Itbl::set(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
Gives the field (see member
p_eta
).
Reimplemented in Sym_tensor_tt.
Definition at line 107 of file sym_tensor_aux.C.
References Scalar::div_tant(), Scalar::dsdt(), Scalar::get_dzpuis(), Tensor::mp, Scalar::mult_r_dzpuis(), Tensor::operator()(), Sym_tensor::p_eta, Map::poisson_angu(), Scalar::poisson_angu(), Scalar::stdsdp(), and Tensor::triad.
void Sym_tensor::exponential_filter_r | ( | int | lzmin, | |
int | lzmax, | |||
int | p, | |||
double | alpha = -16. | |||
) | [virtual, inherited] |
Applies exponential filters to all components (see Scalar::exponential_filter_r
).
Does a loop for Cartesian components, and works in terms of the rr-component, ,
,
W
, X
, T
for spherical components.
Reimplemented from Tensor.
Definition at line 442 of file sym_tensor.C.
References Tensor::cmp, Scalar::div_r(), Sym_tensor::eta(), Scalar::exponential_filter_r(), Map::get_bvect_cart(), Map::get_bvect_spher(), Base_vect::identify(), Tensor::mp, Sym_tensor::mu(), Tensor::n_comp, Tensor::operator()(), Sym_tensor::set_auxiliary(), Tensor::triad, Sym_tensor::ttt(), Sym_tensor::www(), and Sym_tensor::xxx().
void Sym_tensor::exponential_filter_ylm | ( | int | lzmin, | |
int | lzmax, | |||
int | p, | |||
double | alpha = -16. | |||
) | [virtual, inherited] |
Applies exponential filters to all components (see Scalar::exponential_filter_ylm
).
Does a loop for Cartesian components, and works in terms of the r-component, ,
,
W
, X
, T
for spherical components.
Reimplemented from Tensor.
Definition at line 467 of file sym_tensor.C.
References Tensor::cmp, Scalar::div_r(), Sym_tensor::eta(), Scalar::exponential_filter_ylm(), Map::get_bvect_cart(), Map::get_bvect_spher(), Base_vect::identify(), Tensor::mp, Sym_tensor::mu(), Tensor::n_comp, Tensor::operator()(), Sym_tensor::set_auxiliary(), Tensor::triad, Sym_tensor::ttt(), Sym_tensor::www(), and Sym_tensor::xxx().
Itbl Tensor::get_index_type | ( | ) | const [inline, inherited] |
Returns the types of all the indices.
Itbl
) of size valence
containing the type of each index, COV
for a covariant one and CON
for a contravariant one. Definition at line 892 of file tensor.h.
References Tensor::type_indice.
int Tensor::get_index_type | ( | int | i | ) | const [inline, inherited] |
Gives the type (covariant or contravariant) of the index number i
.
i
must be strictly lower than valence
and obey the following convention:
i
= 0 : first index i
= 1 : second index Definition at line 882 of file tensor.h.
References Tensor::type_indice.
const Metric& Sym_tensor_trans::get_met_div | ( | ) | const [inline] |
Returns the metric with respect to which the divergence and the trace are defined.
Definition at line 665 of file sym_tensor.h.
References met_div.
const Map& Tensor::get_mp | ( | ) | const [inline, inherited] |
int Tensor::get_n_comp | ( | ) | const [inline, inherited] |
Returns the number of stored components.
Definition at line 868 of file tensor.h.
References Tensor::n_comp.
int Tensor::get_place_met | ( | const Metric & | metre | ) | const [protected, inherited] |
Returns the position of the pointer on metre
in the array met_depend
.
Definition at line 439 of file tensor.C.
References Tensor::met_depend.
Scalar Sym_tensor::get_tilde_B_from_TT_trace | ( | const Scalar & | tilde_B_tt_in, | |
const Scalar & | trace | |||
) | const [protected, inherited] |
Computes (see
Sym_tensor::p_tilde_b
) from its transverse-traceless part and the trace.
Definition at line 527 of file sym_tensor_aux.C.
References Tbl::annule_hard(), Scalar::annule_hard(), Valeur::c, Valeur::c_cf, Scalar::div_r_dzpuis(), Scalar::dsdr(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Base_val::give_quant_numbers(), Tensor::mp, Base_val::mult_x(), Tensor::operator()(), Mtbl_cf::set(), Scalar::set_dzpuis(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Tensor::triad, and Valeur::ylm().
const Base_vect* Tensor::get_triad | ( | ) | const [inline, inherited] |
Returns the vectorial basis (triad) on which the components are defined.
Definition at line 862 of file tensor.h.
References Tensor::triad.
int Tensor::get_valence | ( | ) | const [inline, inherited] |
void Tensor::inc_dzpuis | ( | int | inc = 1 |
) | [virtual, inherited] |
Increases by inc
units the value of dzpuis
and changes accordingly the values in the compactified external domain (CED).
Reimplemented in Scalar.
Definition at line 812 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), and Tensor::n_comp.
Itbl Tensor_sym::indices | ( | int | pos | ) | const [virtual, inherited] |
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 Tensor.
Definition at line 306 of file tensor_sym.C.
References Tensor_sym::id_sym1, Tensor_sym::id_sym2, Tensor::n_comp, Itbl::set(), and Tensor::valence.
Sym_tensor * Sym_tensor::inverse | ( | ) | const [protected, inherited] |
Returns a pointer on the inverse of the Sym_tensor
(seen as a matrix).
Definition at line 368 of file sym_tensor.C.
References Tensor::mp, Tensor::operator()(), Tensor::set(), Sym_tensor::Sym_tensor(), Tensor::triad, and Tensor::type_indice.
const Vector & Sym_tensor::longit_pot | ( | const Metric & | gam, | |
Param * | par = 0x0 , |
|||
int | method_poisson = 6 | |||
) | const [inherited] |
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 139 of file sym_tensor_decomp.C.
References Tensor::dec_dzpuis(), Tensor_sym::derive_con(), diffrel(), Sym_tensor::divergence(), Map::get_mg(), Mg3d::get_nzone(), Tensor::get_place_met(), maxabs(), Tensor::mp, Sym_tensor::p_longit_pot, Vector::poisson(), and Tensor::set_dependance().
Gives the field (see member
p_mu
).
Definition at line 147 of file sym_tensor_aux.C.
References Scalar::div_tant(), Scalar::dsdt(), Scalar::get_dzpuis(), Tensor::mp, Scalar::mult_r_dzpuis(), Tensor::operator()(), Sym_tensor::p_mu, Map::poisson_angu(), Scalar::poisson_angu(), Scalar::stdsdp(), and Tensor::triad.
const Scalar & Tensor::operator() | ( | int | i1, | |
int | i2, | |||
int | i3, | |||
int | i4 | |||
) | const [inherited] |
Returns the value of a component for a tensor of valence 4 (read-only version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) | |
i3 | value of the third index (1, 2 or 3) | |
i4 | value of the fourth index (1, 2 or 3) |
(i1,i2,i3,i4) Definition at line 779 of file tensor.C.
References Tensor::cmp, Tensor::position(), Itbl::set(), and Tensor::valence.
const Scalar & Tensor::operator() | ( | int | i1, | |
int | i2, | |||
int | i3 | |||
) | const [inherited] |
Returns the value of a component for a tensor of valence 3 (read-only version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) | |
i3 | value of the third index (1, 2 or 3) |
(i1,i2,i3) Definition at line 767 of file tensor.C.
References Tensor::cmp, Tensor::position(), Itbl::set(), and Tensor::valence.
const Scalar & Tensor::operator() | ( | int | i1, | |
int | i2 | |||
) | const [inherited] |
Returns the value of a component for a tensor of valence 2 (read-only version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) |
(i1,i2) Definition at line 756 of file tensor.C.
References Tensor::cmp, Tensor::position(), Itbl::set(), and Tensor::valence.
Returns the value of a component (read-only version).
ind | 1-D Itbl of size valence containing the values of each index specifing the component, with the following storage convention:
|
ind
Definition at line 794 of file tensor.C.
References Tensor::cmp, Itbl::get_dim(), Itbl::get_ndim(), Tensor::position(), and Tensor::valence.
void Tensor::operator+= | ( | const Tensor & | t | ) | [inherited] |
+= Tensor
Reimplemented in Scalar.
Definition at line 567 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::indices(), Tensor::n_comp, Tensor::position(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
void Tensor::operator-= | ( | const Tensor & | t | ) | [inherited] |
-= Tensor
Reimplemented in Scalar.
Definition at line 583 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::indices(), Tensor::n_comp, Tensor::position(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
void Sym_tensor_trans::operator= | ( | const Tensor & | a | ) | [virtual] |
Assignment to a Tensor
.
Reimplemented from Sym_tensor.
Reimplemented in Sym_tensor_tt.
Definition at line 221 of file sym_tensor_trans.C.
References del_deriv(), and operator=().
void Sym_tensor_trans::operator= | ( | const Tensor_sym & | a | ) | [virtual] |
Assignment to a Tensor_sym
.
Reimplemented from Sym_tensor.
Reimplemented in Sym_tensor_tt.
Definition at line 209 of file sym_tensor_trans.C.
References del_deriv(), and operator=().
void Sym_tensor_trans::operator= | ( | const Sym_tensor & | a | ) | [virtual] |
Assignment to a Sym_tensor
.
Reimplemented in Sym_tensor_tt.
Definition at line 198 of file sym_tensor_trans.C.
References del_deriv(), and operator=().
void Sym_tensor_trans::operator= | ( | const Sym_tensor_trans & | a | ) | [virtual] |
Assignment to another Sym_tensor_trans
.
Reimplemented from Sym_tensor.
Reimplemented in Sym_tensor_tt.
Definition at line 182 of file sym_tensor_trans.C.
References del_deriv(), met_div, p_trace, and p_tt.
Sym_tensor_trans Sym_tensor_trans::poisson | ( | const Scalar * | h_guess = 0x0 |
) | const |
Computes the solution of a tensorial transverse Poisson equation with *this
as a source:
In particular, it makes an iteration on the trace of the result, using Sym_tensor::set_WX_det_one
.
h_guess | a pointer on a guess for the trace of the result; it is passed to Sym_tensor::set_WX_det_one . |
Definition at line 95 of file sym_tensor_trans_pde.C.
References Scalar::allocate_all(), Tensor::change_triad(), Tensor::dec_dzpuis(), Sym_tensor::divergence(), Map_af::get_alpha(), Map_af::get_beta(), Map::get_bvect_cart(), Map::get_bvect_spher(), Scalar::get_dzpuis(), Scalar::get_etat(), Map::get_mg(), Mg3d::get_non_axi(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Mg3d::get_type_r(), maxabs(), met_div, Tensor::mp, Tensor::operator()(), poisson(), set_AtBtt_det_one(), Scalar::set_dzpuis(), Scalar::set_etat_one(), Scalar::set_etat_zero(), Scalar::set_grid_point(), Scalar::set_spectral_base(), Tensor::triad, and Scalar::val_grid_point().
int Tensor_sym::position | ( | const Itbl & | ind | ) | const [virtual, inherited] |
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 Tensor.
Definition at line 241 of file tensor_sym.C.
References Itbl::get_dim(), Itbl::get_ndim(), Tensor_sym::id_sym1, Tensor_sym::id_sym2, Itbl::set(), and Tensor::valence.
void Tensor_sym::sauve | ( | FILE * | fd | ) | const [virtual, inherited] |
Save in a binary file.
Reimplemented from Tensor.
Definition at line 368 of file tensor_sym.C.
References fwrite_be(), Tensor_sym::id_sym1, and Tensor_sym::id_sym2.
Scalar & Tensor::set | ( | int | i1, | |
int | i2, | |||
int | i3, | |||
int | i4 | |||
) | [inherited] |
Returns the value of a component for a tensor of valence 4 (read/write version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) | |
i3 | value of the third index (1, 2 or 3) | |
i4 | value of the fourth index (1, 2 or 3) |
(i1,i2,i3,i4) Definition at line 633 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::position(), Itbl::set(), and Tensor::valence.
Scalar & Tensor::set | ( | int | i1, | |
int | i2, | |||
int | i3 | |||
) | [inherited] |
Returns the value of a component for a tensor of valence 3 (read/write version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) | |
i3 | value of the third index (1, 2 or 3) |
(i1,i2,i3) Definition at line 617 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::position(), Itbl::set(), and Tensor::valence.
Scalar & Tensor::set | ( | int | i1, | |
int | i2 | |||
) | [inherited] |
Returns the value of a component for a tensor of valence 2 (read/write version).
i1 | value of the first index (1, 2 or 3) | |
i2 | value of the second index (1, 2 or 3) |
(i1,i2) Definition at line 602 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::position(), Itbl::set(), and Tensor::valence.
Returns the value of a component (read/write version).
ind | 1-D Itbl of size valence containing the values of each index specifing the component, with the following storage convention:
|
ind
Definition at line 650 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Itbl::get_dim(), Itbl::get_ndim(), Tensor::position(), and Tensor::valence.
void Sym_tensor_trans::set_AtB_trace | ( | const Scalar & | a_in, | |
const Scalar & | tb_in, | |||
const Scalar & | trace, | |||
Param * | par_bc = 0x0 , |
|||
Param * | par_mat = 0x0 | |||
) |
Assigns the derived members A
, and the trace.
Other derived members are deduced from the divergence-free condition.
a_in | the A potential (see Sym_tensor::p_aaa ) | |
tb_in | the ![]() Sym_tensor::p_tilde_b ) | |
trace | the trace of the Sym_tensor . |
Definition at line 298 of file sym_tensor_trans_aux.C.
References Scalar::check_dzpuis(), Tensor::mp, Sym_tensor::p_aaa, Sym_tensor::p_tilde_b, Sym_tensor::set_auxiliary(), sol_Dirac_A(), sol_Dirac_tilde_B(), and Tensor::triad.
void Sym_tensor_trans::set_AtBtt_det_one | ( | const Scalar & | a_in, | |
const Scalar & | tbtt_in, | |||
const Scalar * | h_prev = 0x0 , |
|||
Param * | par_bc = 0x0 , |
|||
Param * | par_mat = 0x0 , |
|||
double | precis = 1.e-14 , |
|||
int | it_max = 100 | |||
) |
Assigns the derived member A
and computes from its TT-part (see
Sym_tensor::compute_tilde_B_tt()
).
Other derived members are deduced from the divergence-free condition. Finally, it computes the trace so that *this
+ the flat metric has a determinant equal to 1. It then updates the components accordingly. This function makes an iteration until the relative difference in the trace between two steps is lower than precis
.
a_in | the A potential (see Sym_tensor::p_aaa ) | |
tbtt_in | the TT-part of ![]() Sym_tensor::p_tilde_b and Sym_tensor::compute_tilde_B_tt() ) | |
h_prev | a pointer on a guess for the trace of *this ; if null, then the iteration starts from 0. | |
precis | relative difference in the trace computation to end the iteration. | |
it_max | maximal number of iterations. |
Definition at line 133 of file sym_tensor_trans_aux.C.
References abs(), Scalar::check_dzpuis(), Scalar::get_spectral_base(), Sym_tensor::get_tilde_B_from_TT_trace(), Tensor::inc_dzpuis(), max(), met_div, Tensor::mp, Sym_tensor::p_aaa, Sym_tensor::p_tilde_b, p_trace, p_tt, Sym_tensor::set_auxiliary(), Scalar::set_etat_zero(), Scalar::set_spectral_base(), sol_Dirac_A(), sol_Dirac_tilde_B(), and Tensor::triad.
void 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 | |||
) | [inherited] |
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 262 of file sym_tensor_aux.C.
References Scalar::check_dzpuis(), Sym_tensor::del_deriv(), Scalar::dsdt(), Scalar::get_dzpuis(), Scalar::lapang(), Scalar::mult_r_dzpuis(), Sym_tensor::p_eta, Sym_tensor::p_mu, Sym_tensor::p_ttt, Sym_tensor::p_www, Sym_tensor::p_xxx, Scalar::set_spectral_va(), Scalar::stdsdp(), Tensor::triad, and Valeur::ylm_i().
void Tensor::set_dependance | ( | const Metric & | met | ) | const [protected, inherited] |
To be used to describe the fact that the derivatives members have been calculated with met
.
First it sets a null element of met_depend
to &met
and puts this
in the list of the dependancies of met
.
Definition at line 449 of file tensor.C.
References Tensor::met_depend, and Metric::tensor_depend.
void Sym_tensor_trans::set_der_0x0 | ( | ) | const [protected] |
Sets the pointers on derived quantities to 0x0.
Reimplemented from Sym_tensor.
Reimplemented in Sym_tensor_tt.
Definition at line 170 of file sym_tensor_trans.C.
void Sym_tensor::set_der_met_0x0 | ( | int | i | ) | const [protected, inherited] |
Sets all the i-th components of met_depend
specific to the class Sym_tensor
(p_transverse
, etc.
..) to 0x0.
Reimplemented from Tensor.
Definition at line 331 of file sym_tensor.C.
References Sym_tensor::p_longit_pot, and Sym_tensor::p_transverse.
void Tensor::set_etat_nondef | ( | ) | [virtual, inherited] |
Sets the logical state of all components to ETATNONDEF
(undefined state).
Reimplemented in Scalar.
Definition at line 485 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::n_comp, and Scalar::set_etat_nondef().
void Tensor::set_etat_qcq | ( | ) | [virtual, inherited] |
Sets the logical state of all components to ETATQCQ
(ordinary state).
Reimplemented in Scalar.
Definition at line 477 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::n_comp, and Scalar::set_etat_qcq().
void Tensor::set_etat_zero | ( | ) | [virtual, inherited] |
Sets the logical state of all components to ETATZERO
(zero state).
Reimplemented in Scalar.
Definition at line 493 of file tensor.C.
References Tensor::cmp, Tensor::del_deriv(), Tensor::n_comp, and Scalar::set_etat_zero().
void Sym_tensor_trans::set_hrr_mu_det_one | ( | const Scalar & | hrr, | |
const Scalar & | mu_in, | |||
double | precis = 1.e-14 , |
|||
int | it_max = 100 | |||
) |
Assigns the rr component and the derived member .
Other derived members are deduced from the divergence-free condition. Finally, it computes T
(Sym_tensor::p_ttt
) so that *this
+ the flat metric has a determinant equal to 1. It then updates the components accordingly. This function makes an iteration until the relative difference in T
between two steps is lower than precis
.
hrr | the rr component of the tensor, | |
mu_in | the ![]() | |
precis | relative difference in the trace computation to end the iteration. | |
it_max | maximal number of iterations. |
Definition at line 113 of file sym_tensor_trans_aux.C.
References Scalar::check_dzpuis(), Tensor::dec_dzpuis(), Tensor::inc_dzpuis(), met_div, Tensor::mp, Sym_tensor::p_mu, Sym_tensor_tt::set_rr_mu(), trace_from_det_one(), and Tensor::triad.
Itbl& Tensor::set_index_type | ( | ) | [inline, inherited] |
Sets the types of all the indices.
Itbl
) of size valence
that can be modified (COV
for a covariant one and CON
for a contravariant one) Definition at line 914 of file tensor.h.
References Tensor::type_indice.
int& Tensor::set_index_type | ( | int | i | ) | [inline, inherited] |
Sets the type of the index number i
.
i
must be strictly lower than valence
and obey the following convention:
i
= 0 : first index i
= 1 : second index COV
for a covariant index, CON
for a contravariant one) Definition at line 905 of file tensor.h.
References Itbl::set(), and Tensor::type_indice.
void Sym_tensor::set_longit_trans | ( | const Vector & | v, | |
const Sym_tensor_trans & | a | |||
) | [inherited] |
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 84 of file sym_tensor_decomp.C.
References Tensor::dec_dzpuis(), Sym_tensor::del_deriv(), Tensor::get_index_type(), get_met_div(), Tensor::get_place_met(), Vector::ope_killing(), Sym_tensor::p_longit_pot, Sym_tensor::p_transverse, and Tensor::set_dependance().
void Tensor::set_triad | ( | const Base_vect & | new_triad | ) | [inherited] |
Assigns a new vectorial basis (triad) of decomposition.
NB: this function modifies only the member triad
and leave unchanged the components (member cmp
). In order to change them coherently with the new basis, the function change_triad
(const Base_vect& ) must be called instead.
Definition at line 515 of file tensor.C.
References Tensor::triad.
void Sym_tensor_trans::set_tt_part_det_one | ( | const Sym_tensor_tt & | hijtt, | |
const Scalar * | h_prev = 0x0 , |
|||
Param * | par_mat = 0x0 , |
|||
double | precis = 1.e-14 , |
|||
int | it_max = 100 | |||
) |
Assignes the TT-part of the tensor.
The trace is deduced from the divergence-free condition, through the Dirac system on , so that
*this
+ the flat metric has a determinant equal to 1. It then updates the components accordingly. This function makes an iteration until the relative difference in the trace between two steps is lower than precis
.
hijtt | the TT part for this . | |
h_prev | a pointer on a guess for the trace of *this ; if null, then the iteration starts from 0. | |
precis | relative difference in the trace computation to end the iteration. | |
it_max | maximal number of iterations. |
Definition at line 222 of file sym_tensor_trans_aux.C.
References abs(), Scalar::div_r(), Sym_tensor_tt::eta(), Scalar::get_spectral_base(), Sym_tensor::get_tilde_B_from_TT_trace(), Tensor::inc_dzpuis(), max(), Tensor::mp, Sym_tensor::mu(), p_trace, p_tt, Sym_tensor::set_auxiliary(), Scalar::set_etat_zero(), Scalar::set_spectral_base(), sol_Dirac_tilde_B(), Tensor::triad, Sym_tensor::www(), and Sym_tensor::xxx().
void Sym_tensor_trans::set_tt_trace | ( | const Sym_tensor_tt & | a, | |
const Scalar & | h, | |||
Param * | par = 0x0 | |||
) |
Assigns the derived members p_tt
and p_trace
and updates the components accordingly.
(see the documentation of these derived members for details)
Definition at line 231 of file sym_tensor_trans.C.
References Scalar::check_dzpuis(), Metric::con(), Tensor::dec_dzpuis(), del_deriv(), Tensor_sym::derive_con(), get_met_div(), met_div, Tensor::mp, p_trace, p_tt, and Scalar::poisson().
void Sym_tensor_trans::sol_Dirac_A | ( | const Scalar & | aaa, | |
Scalar & | tilde_mu, | |||
Scalar & | xxx, | |||
const Param * | par_bc = 0x0 | |||
) | const [protected] |
Solves a system of two coupled first-order PDEs obtained from the divergence-free condition (Dirac gauge) and the requirement that the potential A (see Sym_tensor::p_aaa
) has a given value.
The system reads:
Note that this is solved only for and that
(see
Sym_tensor::p_mu
).
aaa | [input] the source A | |
tilde_mu | [output] the solution ![]() | |
xxx | [output] the solution X | |
par_bc | [input] Param to control the boundary conditions |
Definition at line 75 of file sym_tensor_trans_dirac.C.
References Tensor::annule_domain(), Matrice::annule_hard(), Tbl::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::c, Valeur::c_cf, Mtbl_cf::dsdx(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Param::get_int(), Diff_id::get_matrice(), Diff_xdsdx::get_matrice(), Diff_sx::get_matrice(), Diff_dsdx::get_matrice(), Map::get_mg(), Param::get_n_int(), Param::get_n_tbl_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Param::get_tbl_mod(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Base_val::mult_x(), R_CHEBP, Mtbl_cf::set(), Tbl::set(), Matrice::set(), Valeur::set_etat_cf_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), Valeur::ylm(), and Valeur::ylm_i().
void Sym_tensor_trans::sol_Dirac_A2 | ( | const Scalar & | aaa, | |
Scalar & | tilde_mu, | |||
Scalar & | x_new, | |||
Scalar | bound_mu, | |||
const Param * | par_bc | |||
) |
Same resolution as sol_Dirac_Abound, but here the boundary conditions are the degenerate elliptic conditions encontered when solving the Kerr problem.
Definition at line 76 of file sym_tensor_trans_dirac_boundfree.C.
References Tensor::annule_domain(), Matrice::annule_hard(), Tbl::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::c, Valeur::c_cf, Mtbl_cf::dsdx(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Param::get_int(), Diff_sx::get_matrice(), Diff_id::get_matrice(), Diff_dsdx::get_matrice(), Diff_xdsdx::get_matrice(), Map::get_mg(), Param::get_n_int(), Param::get_n_tbl_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Param::get_tbl_mod(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Base_val::mult_x(), Tbl::set(), Matrice::set(), Mtbl_cf::set(), Valeur::set_etat_cf_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Tbl::set_etat_zero(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), Valeur::ylm(), and Valeur::ylm_i().
void Sym_tensor_trans::sol_Dirac_Abound | ( | const Scalar & | aaa, | |
Scalar & | tilde_mu, | |||
Scalar & | x_new, | |||
Scalar | bound_mu, | |||
const Param * | par_bc | |||
) |
Same resolution as sol_Dirac_A, but with inner boundary conditions added.
For now, only Robyn-type boundary conditions on can be imposed.
Definition at line 79 of file sym_tensor_trans_dirac_bound2.C.
References Tensor::annule_domain(), Matrice::annule_hard(), Tbl::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::c, Valeur::c_cf, Mtbl_cf::dsdx(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Param::get_int(), Diff_sx::get_matrice(), Diff_id::get_matrice(), Diff_dsdx::get_matrice(), Diff_xdsdx::get_matrice(), Map::get_mg(), Param::get_n_int(), Param::get_n_tbl_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Param::get_tbl_mod(), Mg3d::get_type_r(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Base_val::mult_x(), Tbl::set(), Matrice::set(), Mtbl_cf::set(), Valeur::set_etat_cf_qcq(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Tbl::set_etat_zero(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), Valeur::ylm(), and Valeur::ylm_i().
void Sym_tensor_trans::sol_Dirac_BC2 | ( | const Scalar & | bb, | |
const Scalar & | cc, | |||
const Scalar & | hh, | |||
Scalar & | hrr, | |||
Scalar & | tilde_eta, | |||
Scalar & | ww, | |||
Scalar | bound_eta, | |||
double | dir, | |||
double | neum, | |||
double | rhor, | |||
Param * | par_bc, | |||
Param * | par_mat | |||
) |
Same resolution as sol_Dirac_tilde_B, but with inner boundary conditions added.
The difference is here, one has to put B and C values in (and not only ). For now, only Robyn-type boundary conditions on
can be imposed.
Definition at line 578 of file sym_tensor_trans_dirac_bound2.C.
References Param::add_int_mod(), Param::add_itbl_mod(), Param::add_matrice_mod(), Param::add_tbl_mod(), Tensor::annule_domain(), Matrice::annule_hard(), Tbl::annule_hard(), Itbl::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Scalar::annule_l(), Valeur::c, Valeur::c_cf, Scalar::check_dzpuis(), Param::clean_all(), Valeur::coef_i(), Scalar::div_r_dzpuis(), Scalar::dsdr(), Mtbl_cf::dsdx(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Param::get_int(), Param::get_int_mod(), Param::get_itbl_mod(), Diff_id::get_matrice(), Diff_xdsdx::get_matrice(), Diff_sx::get_matrice(), Diff_dsdx::get_matrice(), Param::get_matrice_mod(), Map::get_mg(), Param::get_n_int(), Param::get_n_int_mod(), Param::get_n_itbl_mod(), Param::get_n_matrice_mod(), Param::get_n_tbl_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Param::get_tbl_mod(), Mg3d::get_type_p(), Mg3d::get_type_r(), Mg3d::get_type_t(), Base_val::give_lmax(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Scalar::mult_r_dzpuis(), R_CHEBP, Matrice::set(), Tbl::set(), Itbl::set(), Mtbl_cf::set(), Valeur::set_etat_cf_qcq(), Matrice::set_etat_qcq(), Scalar::set_etat_qcq(), Tbl::set_etat_qcq(), Itbl::set_etat_qcq(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), sol_Dirac_l01(), Scalar::std_spectral_base(), Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), Valeur::ylm(), and Valeur::ylm_i().
void Sym_tensor_trans::sol_Dirac_BC3 | ( | const Scalar & | bb, | |
const Scalar & | hh, | |||
Scalar & | hrr, | |||
Scalar & | tilde_eta, | |||
Scalar & | ww, | |||
Scalar | bound_hrr, | |||
Scalar | bound_eta, | |||
Param * | par_bc, | |||
Param * | par_mat | |||
) |
Same resolution as sol_Dirac_Abound, but here the boundary conditions are the degenerate elliptic conditions encontered when solving the Kerr problem.
Definition at line 596 of file sym_tensor_trans_dirac_boundfree.C.
References Param::add_int_mod(), Param::add_itbl_mod(), Param::add_matrice_mod(), Param::add_tbl_mod(), Tensor::annule_domain(), Matrice::annule_hard(), Tbl::annule_hard(), Itbl::annule_hard(), Scalar::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_l(), Valeur::c, Valeur::c_cf, Scalar::check_dzpuis(), Param::clean_all(), Scalar::div_r_dzpuis(), Scalar::dsdr(), Mtbl_cf::dsdx(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Param::get_int(), Param::get_int_mod(), Param::get_itbl_mod(), Diff_id::get_matrice(), Diff_xdsdx::get_matrice(), Diff_sx::get_matrice(), Diff_dsdx::get_matrice(), Param::get_matrice_mod(), Map::get_mg(), Param::get_n_int(), Param::get_n_int_mod(), Param::get_n_itbl_mod(), Param::get_n_matrice_mod(), Param::get_n_tbl_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Param::get_tbl_mod(), Mg3d::get_type_p(), Mg3d::get_type_r(), Mg3d::get_type_t(), Base_val::give_lmax(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Scalar::mult_r_dzpuis(), R_CHEBP, Mtbl_cf::set(), Matrice::set(), Tbl::set(), Itbl::set(), Valeur::set_etat_cf_qcq(), Matrice::set_etat_qcq(), Scalar::set_etat_qcq(), Tbl::set_etat_qcq(), Itbl::set_etat_qcq(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), Valeur::ylm(), and Valeur::ylm_i().
void Sym_tensor_trans::sol_Dirac_l01 | ( | const Scalar & | hh, | |
Scalar & | hrr, | |||
Scalar & | tilde_eta, | |||
Param * | par_mat | |||
) | const [protected] |
Solves the same system as Sym_tensor_trans::sol_Dirac_tilde_B
but only for .
In these particular cases, W =0 the system is simpler and homogeneous solutions are different.
Definition at line 1427 of file sym_tensor_trans_dirac.C.
References Param::add_matrice_mod(), Matrice::annule_hard(), Tbl::annule_hard(), Itbl::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_hard(), Valeur::c_cf, Scalar::div_r_dzpuis(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Diff_id::get_matrice(), Diff_xdsdx::get_matrice(), Diff_sx::get_matrice(), Diff_dsdx::get_matrice(), Param::get_matrice_mod(), Map::get_mg(), Param::get_n_matrice_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_base(), Scalar::get_spectral_va(), Mg3d::get_type_r(), Base_val::give_lmax(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Base_val::mult_x(), R_CHEBP, Mtbl_cf::set(), Tbl::set(), Itbl::set(), Matrice::set(), Tbl::set_etat_qcq(), Matrice::set_etat_qcq(), Matrice::set_lu(), Scalar::set_spectral_va(), Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), and Valeur::ylm().
void Sym_tensor_trans::sol_Dirac_tilde_B | ( | const Scalar & | tilde_b, | |
const Scalar & | hh, | |||
Scalar & | hrr, | |||
Scalar & | tilde_eta, | |||
Scalar & | www, | |||
Param * | par_bc = 0x0 , |
|||
Param * | par_mat = 0x0 | |||
) | const [protected] |
Solves a system of three coupled first-order PDEs obtained from divergence-free conditions (Dirac gauge) and the requirement that the potential (see
Sym_tensor::p_tilde_b
) has a given value.
The system reads:
Note that (for definitions, see derived members of
Sym_tensor
).
tilde_b | [input] the source ![]() | |
hh | [input] the trace of the tensor | |
hrr | [output] the rr component of the result | |
tilde_eta | [output] the solution ![]() | |
www | [output] the solution W | |
par_bc | [input] Param to control the boundary conditions | |
par_mat | [input/output] Param in which the operator matrix is stored. |
Definition at line 576 of file sym_tensor_trans_dirac.C.
References Param::add_int_mod(), Param::add_itbl_mod(), Param::add_matrice_mod(), Param::add_tbl_mod(), Tensor::annule_domain(), Matrice::annule_hard(), Tbl::annule_hard(), Itbl::annule_hard(), Scalar::annule_hard(), Mtbl_cf::annule_hard(), Scalar::annule_l(), Valeur::c, Valeur::c_cf, Scalar::check_dzpuis(), Param::clean_all(), Scalar::div_r_dzpuis(), Scalar::dsdr(), Mtbl_cf::dsdx(), Map_af::get_alpha(), Map_af::get_beta(), Scalar::get_etat(), Param::get_int(), Param::get_int_mod(), Param::get_itbl_mod(), Diff_id::get_matrice(), Diff_xdsdx::get_matrice(), Diff_sx::get_matrice(), Diff_dsdx::get_matrice(), Param::get_matrice_mod(), Map::get_mg(), Param::get_n_int(), Param::get_n_int_mod(), Param::get_n_itbl_mod(), Param::get_n_matrice_mod(), Param::get_n_tbl_mod(), Mg3d::get_np(), Mg3d::get_nr(), Mg3d::get_nt(), Mg3d::get_nzone(), Scalar::get_spectral_va(), Param::get_tbl_mod(), Mg3d::get_type_p(), Mg3d::get_type_r(), Mg3d::get_type_t(), Base_val::give_lmax(), Base_val::give_quant_numbers(), Matrice::inverse(), Tensor::mp, Scalar::mult_r(), Scalar::mult_r_dzpuis(), Base_val::mult_x(), R_CHEBP, Mtbl_cf::set(), Matrice::set(), Tbl::set(), Itbl::set(), Valeur::set_etat_cf_qcq(), Matrice::set_etat_qcq(), Scalar::set_etat_qcq(), Tbl::set_etat_qcq(), Itbl::set_etat_qcq(), Matrice::set_lu(), Scalar::set_spectral_base(), Scalar::set_spectral_va(), sol_Dirac_l01(), Mtbl_cf::val_in_bound_jk(), Mtbl_cf::val_out_bound_jk(), Valeur::ylm(), and Valeur::ylm_i().
void Sym_tensor_trans::sol_elliptic_ABC | ( | Sym_tensor & | source, | |
Scalar | aaa, | |||
Scalar | bbb, | |||
Scalar | ccc | |||
) |
Finds spectral potentials A, B, C of solution of an tensorial TT elliptic equation, given the source.
void Tensor::spectral_display | ( | const char * | comment = 0x0 , |
|
double | threshold = 1.e-7 , |
|||
int | precision = 4 , |
|||
ostream & | ostr = cout | |||
) | const [virtual, inherited] |
Displays the spectral coefficients and the associated basis functions of each component.
This function shows only the values greater than a given threshold.
comment | comment to be printed at top of the display (default: 0x0 = nothing printed) | |
threshold | [input] Value above which a coefficient is printed (default: 1.e-7) | |
precision | [input] Number of printed digits (default: 4) | |
ostr | [input] Output stream used for the printing (default: cout) |
Reimplemented in Scalar.
Definition at line 870 of file tensor.C.
References Tensor::cmp, Tensor::indices(), Tensor::n_comp, Scalar::spectral_display(), and Tensor::valence.
void Tensor::std_spectral_base | ( | ) | [virtual, inherited] |
Sets the standard spectal bases of decomposition for each component.
To be used only with valence
lower than or equal 2.
Reimplemented in Scalar, and Vector.
Definition at line 922 of file tensor.C.
References Tensor::cmp, Map::get_bvect_cart(), Map::get_bvect_spher(), Map::get_mg(), Base_vect::identify(), Tensor::indices(), Tensor::mp, Tensor::n_comp, Scalar::set_spectral_base(), Mg3d::std_base_vect_cart(), Mg3d::std_base_vect_spher(), Scalar::std_spectral_base(), Tensor::triad, and Tensor::valence.
void Tensor::std_spectral_base_odd | ( | ) | [virtual, inherited] |
Sets the standard odd spectal bases of decomposition for each component.
Currently only implemented for a scalar.
Reimplemented in Scalar.
Definition at line 978 of file tensor.C.
References Tensor::cmp, Scalar::std_spectral_base_odd(), and Tensor::valence.
int Tensor_sym::sym_index1 | ( | ) | const [inline, inherited] |
Number of the first symmetric index (0<=
id_sym1
< valence
).
Definition at line 1145 of file tensor.h.
References Tensor_sym::id_sym1.
int Tensor_sym::sym_index2 | ( | ) | const [inline, inherited] |
Number of the second symmetric index (id_sym1
< id_sym2
< valence
).
Definition at line 1150 of file tensor.h.
References Tensor_sym::id_sym2.
const Scalar & Sym_tensor_trans::the_trace | ( | ) | const |
Returns the trace of the tensor with respect to metric *met_div
.
Definition at line 266 of file sym_tensor_trans.C.
References met_div, p_trace, Tensor::trace(), and Tensor::type_indice.
Trace with respect to a given metric for a valence 2 tensor.
gam | metric used to raise or lower one of the indices, in order to take the trace |
Definition at line 193 of file tensor_calculus.C.
References Metric::con(), contract(), Metric::cov(), Tensor::trace(), Tensor::type_indice, and Tensor::valence.
Scalar Tensor::trace | ( | ) | const [inherited] |
Trace on two different type indices for a valence 2 tensor.
Definition at line 176 of file tensor_calculus.C.
References Tensor::mp, Tensor::operator()(), Scalar::set_etat_zero(), Tensor::type_indice, and Tensor::valence.
Trace with respect to a given metric.
ind1 | first index for the contraction, with the following convention :
| |
ind2 | second index for the contraction | |
gam | metric used to raise or lower ind1 in order that it has a opposite type than ind2 |
Definition at line 149 of file tensor_calculus.C.
References Metric::con(), contract(), Metric::cov(), Tensor::trace(), Tensor::type_indice, and Tensor::valence.
Tensor Tensor::trace | ( | int | ind1, | |
int | ind2 | |||
) | const [inherited] |
Trace on two different type indices.
ind1 | first index for the contraction, with the following convention :
| |
ind2 | second index for the contraction |
Definition at line 90 of file tensor_calculus.C.
References Tensor::cmp, Tensor::get_n_comp(), Tensor::indices(), Tensor::mp, Tensor::position(), Tensor::set(), Itbl::set(), Scalar::set_etat_zero(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
void Sym_tensor_trans::trace_from_det_one | ( | const Sym_tensor_tt & | htt, | |
double | precis = 1.e-14 , |
|||
int | it_max = 100 | |||
) |
Assigns the derived member p_tt
and computes the trace so that *this
+ the flat metric has a determinant equal to 1.
It then updates the components accordingly, with a dzpuis
= 2. This function makes an iteration until the relative difference in the trace between two steps is lower than precis
.
htt | the transverse traceless part; all components must have dzpuis = 2. | |
precis | relative difference in the trace computation to end the iteration. | |
it_max | maximal number of iterations. |
Definition at line 311 of file sym_tensor_trans.C.
References abs(), Scalar::check_dzpuis(), Tensor::cmp, Scalar::dec_dzpuis(), get_met_div(), Tensor::get_n_comp(), max(), met_div, Tensor::mp, Scalar::set_etat_zero(), and set_tt_trace().
const Sym_tensor_trans & Sym_tensor::transverse | ( | const Metric & | gam, | |
Param * | par = 0x0 , |
|||
int | method_poisson = 6 | |||
) | const [inherited] |
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 106 of file sym_tensor_decomp.C.
References Tensor::cmp, Tensor::get_place_met(), Tensor::inc_dzpuis(), Sym_tensor::longit_pot(), Tensor::mp, Tensor::n_comp, Vector::ope_killing(), Sym_tensor::p_transverse, Tensor::set_dependance(), Tensor::triad, and Tensor::type_indice.
const Sym_tensor_tt & Sym_tensor_trans::tt_part | ( | Param * | par = 0x0 |
) | const |
Returns the transverse traceless part of the tensor, the trace being defined with respect to metric *met_div
.
Definition at line 280 of file sym_tensor_trans.C.
References Metric::con(), Tensor::dec_dzpuis(), Tensor_sym::derive_con(), Scalar::derive_con(), Scalar::get_dzpuis(), Tensor::inc_dzpuis(), met_div, Tensor::mp, p_tt, Scalar::poisson(), the_trace(), and Tensor::triad.
const Scalar & Sym_tensor::ttt | ( | ) | const [inherited] |
Gives the field T (see member p_ttt
).
Definition at line 186 of file sym_tensor_aux.C.
References Sym_tensor::p_ttt, and 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 221 of file tensor_calculus.C.
References Metric::con(), contract(), Tensor::indices(), Tensor::mp, Tensor::n_comp, Tensor::set(), Itbl::set(), Tensor::triad, Tensor::type_indice, and Tensor::valence.
Computes a new tensor by raising or lowering all the indices of *this
.
gam | metric used to lower the contravariant indices and raising the covariant ones. |
Definition at line 301 of file tensor_calculus.C.
References Tensor::down(), Tensor::Tensor(), Tensor::type_indice, Tensor::up(), and Tensor::valence.
const Scalar & Sym_tensor::www | ( | ) | const [inherited] |
Gives the field W (see member p_www
).
Definition at line 205 of file sym_tensor_aux.C.
References Scalar::div_tant(), Scalar::dsdt(), Tensor::operator()(), Sym_tensor::p_www, Scalar::poisson_angu(), Scalar::stdsdp(), and Tensor::triad.
const Scalar & Sym_tensor::xxx | ( | ) | const [inherited] |
Gives the field X (see member p_xxx
).
Definition at line 236 of file sym_tensor_aux.C.
References Scalar::div_tant(), Scalar::dsdt(), Tensor::operator()(), Sym_tensor::p_xxx, Scalar::poisson_angu(), Scalar::stdsdp(), and Tensor::triad.
Tensor_sym operator* | ( | const Tensor_sym & | a, | |
const Tensor_sym & | b | |||
) | [friend, inherited] |
Tensorial product of two symmetric tensors.
NB: the output is an object of class Tensor_sym
, with the two symmetric indices corresponding to the symmetric indices of tensor a
. This means that the symmetries of tensor b
indices are not used in the storage, since there is currently no class in Lorene to manage tensors with more than two symmetric indices.
Definition at line 147 of file tensor_sym_calculus.C.
Scalar** Tensor::cmp [protected, inherited] |
int Tensor_sym::id_sym1 [protected, inherited] |
int Tensor_sym::id_sym2 [protected, inherited] |
const Metric* Tensor::met_depend[N_MET_MAX] [mutable, protected, inherited] |
const Metric* const Sym_tensor_trans::met_div [protected] |
Metric with respect to which the divergence and the trace are defined.
Definition at line 610 of file sym_tensor.h.
const Map* const Tensor::mp [protected, inherited] |
int Tensor::n_comp [protected, inherited] |
Scalar* Sym_tensor::p_aaa [mutable, protected, inherited] |
Field A defined from X and insensitive to the longitudinal part of the
Sym_tensor
(only for ).
Its definition reads
Definition at line 318 of file sym_tensor.h.
Tensor* Tensor::p_derive_con[N_MET_MAX] [mutable, protected, inherited] |
Array of pointers on the contravariant derivatives of this
with respect to various metrics.
See the comments of met_depend
. See also the comments of method derive_con()
for a precise definition of a "contravariant" derivative.
Tensor* Tensor::p_derive_cov[N_MET_MAX] [mutable, protected, inherited] |
Array of pointers on the covariant derivatives of this
with respect to various metrics.
See the comments of met_depend
. See also the comments of method derive_cov()
for the index convention of the covariant derivation.
Tensor* Tensor::p_divergence[N_MET_MAX] [mutable, protected, inherited] |
Array of pointers on the divergence of this
with respect to various metrics.
See the comments of met_depend
. See also the comments of method divergence()
for a precise definition of a the divergence with respect to a given metric.
Scalar* Sym_tensor::p_eta [mutable, protected, inherited] |
Field such that the components
of the tensor are written (has only meaning with spherical components!):
.
Definition at line 256 of file sym_tensor.h.
Vector* Sym_tensor::p_longit_pot[N_MET_MAX] [mutable, protected, inherited] |
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 242 of file sym_tensor.h.
Scalar* Sym_tensor::p_mu [mutable, protected, inherited] |
Field such that the components
of the tensor are written (has only meaning with spherical components!):
.
Definition at line 270 of file sym_tensor.h.
Scalar* Sym_tensor::p_tilde_b [mutable, protected, inherited] |
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 330 of file sym_tensor.h.
Scalar* Sym_tensor::p_tilde_c [mutable, protected, inherited] |
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 342 of file sym_tensor.h.
Scalar* Sym_tensor_trans::p_trace [mutable, protected] |
Trace with respect to the metric *met_div
.
Definition at line 613 of file sym_tensor.h.
Sym_tensor_trans* Sym_tensor::p_transverse[N_MET_MAX] [mutable, protected, inherited] |
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 235 of file sym_tensor.h.
Sym_tensor_tt* Sym_tensor_trans::p_tt [mutable, protected] |
Traceless part with respect to the metric *met_div
.
Definition at line 616 of file sym_tensor.h.
Scalar* Sym_tensor::p_ttt [mutable, protected, inherited] |
Field T defined as .
Definition at line 311 of file sym_tensor.h.
Scalar* Sym_tensor::p_www [mutable, protected, inherited] |
Field W such that the components and
of the tensor are written (has only meaning with spherical components!):
.
Definition at line 289 of file sym_tensor.h.
Scalar* Sym_tensor::p_xxx [mutable, protected, inherited] |
Field X such that the components and
of the tensor are written (has only meaning with spherical components!):
.
Definition at line 308 of file sym_tensor.h.
const Base_vect* Tensor::triad [protected, inherited] |
Itbl Tensor::type_indice [protected, inherited] |
int Tensor::valence [protected, inherited] |