|
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. | |
| Sym_tensor (const Map &map, int tipe, const Base_vect &triad_i) | |
| Standard constructor when both indices are of the same type. | |
| Sym_tensor (const Sym_tensor &a) | |
| Copy constructor. | |
| Sym_tensor (const Tensor &a) | |
Constructor from a Tensor . | |
| Sym_tensor (const Map &map, const Base_vect &triad_i, FILE *fich) | |
Constructor from a file (see sauve(FILE*) ). | |
| virtual | ~Sym_tensor () |
| Destructor. | |
| virtual void | operator= (const Sym_tensor &a) |
Assignment to another 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_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. | |
| 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...) to 0x0. | |
| 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 | |
| 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... The i-th element of this array is the Metric used to compute the i-th element of 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 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 223 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 142 of file sym_tensor.C.
References Lorene::Map(), set_der_0x0(), and Lorene::Tensor_sym::Tensor_sym().
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 152 of file sym_tensor.C.
References Lorene::Map(), set_der_0x0(), and Lorene::Tensor_sym::Tensor_sym().
| Lorene::Sym_tensor::Sym_tensor | ( | const Sym_tensor & | a | ) |
Copy constructor.
Definition at line 160 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(), set_der_0x0(), Sym_tensor(), and Lorene::Tensor_sym::Tensor_sym().
| 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 193 of file sym_tensor.C.
References Lorene::Tensor::cmp, Lorene::Tensor_sym::indices(), Lorene::Tensor::mp, Lorene::Tensor::n_comp, Lorene::Tensor::position(), set_der_0x0(), Lorene::Tensor_sym::Tensor_sym(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, 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 210 of file sym_tensor.C.
References Lorene::Map(), Lorene::Tensor::n_comp, set_der_0x0(), Lorene::Tensor_sym::Tensor_sym(), and Lorene::Tensor::valence.
|
virtual |
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 316 of file sym_tensor_aux.C.
References Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::div_tant(), Lorene::dsdr(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Tensor::mp, Lorene::Tensor::operator()(), p_aaa, Lorene::Scalar::poisson_angu(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, and xxx().
| const Scalar & Lorene::Sym_tensor::compute_tilde_B | ( | bool | output_ylm = true, |
| Param * | par = 0x0 ) const |
Gives the field 
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 362 of file sym_tensor_aux.C.
References Lorene::Scalar::annule_hard(), Lorene::Valeur::c_cf, Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::div_tant(), Lorene::Scalar::dsdr(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Tensor::operator()(), p_tilde_b, Lorene::Scalar::poisson_angu(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, ttt(), www(), and Lorene::Valeur::ylm().
Gives the field 
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 478 of file sym_tensor_aux.C.
References Lorene::Scalar::annule_hard(), Lorene::Valeur::c, Lorene::Valeur::c_cf, compute_tilde_B(), Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::dsdr(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Tensor::operator()(), Lorene::Mtbl_cf::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_spectral_va(), ttt(), Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| const Scalar & Lorene::Sym_tensor::compute_tilde_C | ( | bool | output_ylm = true, |
| Param * | par = 0x0 ) const |
Gives the field 
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 596 of file sym_tensor_aux.C.
References Lorene::Scalar::annule_hard(), Lorene::Valeur::c_cf, Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::div_tant(), Lorene::Scalar::dsdr(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Tensor::operator()(), p_tilde_c, Lorene::Scalar::poisson_angu(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, ttt(), www(), and Lorene::Valeur::ylm().
|
protectedvirtual |
Deletes the derived quantities.
Reimplemented from Lorene::Tensor.
Reimplemented in Lorene::Sym_tensor_trans, and Lorene::Sym_tensor_tt.
Definition at line 286 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 320 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().
| 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 360 of file sym_tensor.C.
References Lorene::Tensor::compute_derive_lie(), Lorene::Tensor::mp, Sym_tensor(), Lorene::Tensor::triad, and Lorene::Tensor::type_indice.
| const Vector & Lorene::Sym_tensor::divergence | ( | const Metric & | gam | ) | const |
Returns the divergence of this with respect to a Metric .
The indices are assumed to be contravariant.
Definition at line 349 of file sym_tensor.C.
References Lorene::Tensor::divergence().
Gives the field 
p_eta ).
Reimplemented in Lorene::Sym_tensor_tt.
Definition at line 111 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::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 446 of file sym_tensor.C.
References Lorene::Tensor::cmp, Lorene::Scalar::div_r(), eta(), Lorene::Scalar::exponential_filter_r(), exponential_filter_r(), 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 471 of file sym_tensor.C.
References Lorene::Tensor::cmp, Lorene::Scalar::div_r(), eta(), Lorene::Scalar::exponential_filter_ylm(), exponential_filter_ylm(), Lorene::Tensor::mp, mu(), Lorene::Tensor::n_comp, Lorene::Tensor::operator()(), set_auxiliary(), Lorene::Tensor::triad, ttt(), www(), and xxx().
|
protected |
Computes 
Sym_tensor::p_tilde_b ) from its transverse-traceless part and the trace.
Definition at line 531 of file sym_tensor_aux.C.
References Lorene::Scalar::annule_hard(), Lorene::Tbl::annule_hard(), Lorene::Valeur::c, Lorene::Valeur::c_cf, Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::dsdr(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::get_etat(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Base_val::mult_x(), Lorene::Tensor::operator()(), Lorene::Mtbl_cf::set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_spectral_base(), Lorene::Scalar::set_spectral_va(), Lorene::Tensor::triad, and Lorene::Valeur::ylm().
|
protected |
Returns a pointer on the inverse of the Sym_tensor
(seen as a matrix).
Definition at line 372 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 
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 143 of file sym_tensor_decomp.C.
References Lorene::Tensor::dec_dzpuis(), Lorene::Tensor::derive_con(), Lorene::diffrel(), divergence(), Lorene::Tensor::divergence(), Lorene::Vector::divergence(), Lorene::Tensor::get_place_met(), Lorene::maxabs(), Lorene::Tensor::mp, p_longit_pot, Lorene::Vector::poisson(), and Lorene::Tensor::set_dependance().
Gives the field 
p_mu ).
Definition at line 151 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::Scalar::poisson_angu(), Lorene::Scalar::stdsdp(), and Lorene::Tensor::triad.
|
virtual |
Assignment to another Sym_tensor.
Reimplemented in Lorene::Sym_tensor_trans, and Lorene::Sym_tensor_tt.
Definition at line 234 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, Lorene::Tensor::set_dependance(), and Sym_tensor().
|
virtual |
Assignment to a Tensor .
The symmetry is assumed but not checked.
Reimplemented from Lorene::Tensor_sym.
Reimplemented in Lorene::Sym_tensor_trans, and Lorene::Sym_tensor_tt.
Definition at line 275 of file sym_tensor.C.
References del_deriv(), and Lorene::Tensor_sym::operator=().
|
virtual |
Assignment to a Tensor_sym.
Reimplemented from Lorene::Tensor_sym.
Reimplemented in Lorene::Sym_tensor_trans, and Lorene::Sym_tensor_tt.
Definition at line 267 of file sym_tensor.C.
References del_deriv(), Lorene::Tensor_sym::operator=(), and Lorene::Tensor_sym::Tensor_sym().
| 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 
p_eta , p_mu , p_www, p_xxx and p_ttt , fro, their values and 

It updates the other components accordingly.
Definition at line 266 of file sym_tensor_aux.C.
References Lorene::Scalar::check_dzpuis(), del_deriv(), Lorene::Scalar::dsdt(), Lorene::Scalar::get_dzpuis(), Lorene::Scalar::lapang(), p_eta, p_mu, p_ttt, p_www, p_xxx, Lorene::Tensor::set(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, and Lorene::Valeur::ylm_i().
|
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 335 of file sym_tensor.C.
References p_longit_pot, and p_transverse.
| 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 88 of file sym_tensor_decomp.C.
References Lorene::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().
| const Sym_tensor_trans & Lorene::Sym_tensor::transverse | ( | const Metric & | gam, |
| Param * | par = 0x0, | ||
| int | method_poisson = 6 ) const |
Computes the transverse part 
Denoting *this by 
![\[ T^{ij} = {}^t T^{ij} + \nabla^i W^j + \nabla^j W^i
\qquad\mbox{with}\quad \nabla_j {}^t T^{ij} = 0
*\]](form_814.png)
where 


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 110 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(), Sym_tensor(), 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 190 of file sym_tensor_aux.C.
References p_ttt, and Lorene::Tensor::triad.
| const Scalar & Lorene::Sym_tensor::www | ( | ) | const |
Gives the field W (see member p_www ).
Definition at line 209 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 240 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.
|
friend |
Definition at line 588 of file sym_tensor.h.
|
mutableprotected |
Field A defined from X and 
Sym_tensor (only for 
Its definition reads
![\[A = \frac{\partial X}{\partial r} - \frac{\mu}{r^2}.
\]](form_827.png)
Definition at line 322 of file sym_tensor.h.
|
mutableprotected |
Field 

![\[ T^{r\theta} = {1\over r} \left( {\partial \eta \over \partial\theta} -
{1\over\sin\theta} {\partial \mu \over \partial\varphi} \right)
*\]](form_818.png)
![\[ T^{r\varphi} = {1\over r} \left( {1\over\sin\theta}
{\partial \eta \over \partial\varphi}
+ {\partial \mu \over \partial\theta} \right)
*\]](form_819.png)
Definition at line 260 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 246 of file sym_tensor.h.
|
mutableprotected |
Field 

![\[ T^{r\theta} = {1\over r} \left( {\partial \eta \over \partial\theta} -
{1\over\sin\theta} {\partial \mu \over \partial\varphi} \right)
*\]](form_820.png)
![\[ T^{r\varphi} = {1\over r} \left( {1\over\sin\theta}
{\partial \eta \over \partial\varphi}
+ {\partial \mu \over \partial\theta} \right)
*\]](form_819.png)
Definition at line 274 of file sym_tensor.h.
|
mutableprotected |
Field 

Sym_tensor.
It is defined for each multipolar momentum 
![\[\tilde{B} = (\ell + 2) \frac{\partial W}{\partial r} + \ell(\ell + 2)
\frac{W}{r} - \frac{2\eta}{r^2} + \frac{(\ell +2)T}{2r(\ell + 1)}
+ \frac{1}{2(\ell + 1)} \frac{\partial T}{\partial r} - \frac{h^{rr}}
{(\ell + 1)r}.
\]](form_830.png)
Definition at line 334 of file sym_tensor.h.
|
mutableprotected |
Field 

Sym_tensor.
It is defined for each multipolar momentum 
![\[\tilde{C} = - (\ell - 1) \frac{\partial W}{\partial r} + (\ell + 1)(\ell - 1)
\frac{W}{r} - \frac{2\eta}{r^2} + \frac{(\ell - 1)T}{2r\ell}
- \frac{1}{2 \ell } \frac{\partial T}{\partial r} - \frac{h^{rr}}
{\ell r}.
\]](form_832.png)
Definition at line 346 of file sym_tensor.h.
|
mutableprotected |
Array of the transverse part 
Denoting *this by 
![\[ T^{ij} = {}^t T^{ij} + \nabla^i W^j + \nabla^j W^i
\qquad\mbox{with}\quad \nabla_j {}^t T^{ij} = 0
*\]](form_814.png)
where 


p_longit_pot below)
Definition at line 239 of file sym_tensor.h.
|
mutableprotected |
Field T defined as 
Definition at line 315 of file sym_tensor.h.
|
mutableprotected |
Field W such that the components 

![\[ \frac{1}{2}\left(T^{\theta\theta} - T^{\varphi\varphi} \right)
= \frac{\partial^2 W}{\partial\theta^2} - \frac{1}{\tan
\theta} \frac{\partial W}{\partial \theta} - \frac{1}{\sin^2 \theta}
\frac{\partial^2 W}{\partial \varphi^2} - 2\frac{\partial}{\partial \theta}
\left( \frac{1}{\sin \theta} \frac{\partial X}{\partial \varphi} \right) ,
*\]](form_823.png)
![\[ T^{\theta\varphi} = \frac{\partial^2 X}{\partial\theta^2} - \frac{1}{\tan
\theta} \frac{\partial X}{\partial \theta} - \frac{1}{\sin^2 \theta}
\frac{\partial^2 X}{\partial \varphi^2} + 2\frac{\partial}{\partial \theta}
\left( \frac{1}{\sin \theta} \frac{\partial W}{\partial \varphi} \right) .
*\]](form_824.png)
Definition at line 293 of file sym_tensor.h.
|
mutableprotected |
Field X such that the components 

![\[ \frac{1}{2}\left(T^{\theta\theta} - T^{\varphi\varphi} \right)
= \frac{\partial^2 W}{\partial\theta^2} - \frac{1}{\tan
\theta} \frac{\partial W}{\partial \theta} - \frac{1}{\sin^2 \theta}
\frac{\partial^2 W}{\partial \varphi^2} - 2\frac{\partial}{\partial \theta}
\left( \frac{1}{\sin \theta} \frac{\partial X}{\partial \varphi} \right) ,
*\]](form_823.png)
![\[ T^{\theta\varphi} = \frac{\partial^2 X}{\partial\theta^2} - \frac{1}{\tan
\theta} \frac{\partial X}{\partial \theta} - \frac{1}{\sin^2 \theta}
\frac{\partial^2 X}{\partial \varphi^2} + 2\frac{\partial}{\partial \theta}
\left( \frac{1}{\sin \theta} \frac{\partial W}{\partial \varphi} \right) .
*\]](form_824.png)
Definition at line 312 of file sym_tensor.h.