|
LORENE
|
Tensor field of valence 1. More...
#include <vector.h>
Public Member Functions | |
| Vector (const Map &map, int tipe, const Base_vect &triad_i) | |
| Standard constructor. | |
| Vector (const Map &map, int tipe, const Base_vect *triad_i) | |
| Standard constructor with the triad passed as a pointer. | |
| Vector (const Vector &a) | |
| Copy constructor. | |
| Vector (const Tensor &a) | |
Constructor from a Tensor . | |
| Vector (const Map &map, const Base_vect &triad_i, FILE *fich) | |
Constructor from a file (see Tensor::sauve(FILE*) ). | |
| virtual | ~Vector () |
| Destructor. | |
| virtual void | change_triad (const Base_vect &) |
| Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly. | |
| virtual void | operator= (const Vector &a) |
| Assignment from a Vector. | |
| virtual void | operator= (const Tensor &a) |
| Assignment from a Tensor. | |
| void | set_vr_eta_mu (const Scalar &vr_i, const Scalar &eta_i, const Scalar &mu_i) |
Defines the components through potentials ![]() ![]() p_eta and p_mu ), as well as the ![]() | |
| void | decompose_div (const Metric &) const |
Makes the Helmholtz decomposition (see documentation of p_potential ) of this with respect to a given Metric , only in the case of contravariant vectors. | |
| const Scalar & | potential (const Metric &) const |
| Returns the potential in the Helmholtz decomposition. | |
| const Vector_divfree & | div_free (const Metric &) const |
| Returns the div-free vector in the Helmholtz decomposition. | |
| 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 ). | |
| Scalar & | set (int) |
| Read/write access to a component. | |
| const Scalar & | operator() (int) const |
| Readonly access to a component. | |
| virtual int | position (const Itbl &idx) const |
Returns the position in the Scalar array cmp of a component given by its index. | |
| virtual Itbl | indices (int place) const |
Returns the index of a component given by its position in the Scalar array cmp . | |
| virtual void | std_spectral_base () |
| Sets the standard spectal bases of decomposition for each component. | |
| virtual void | pseudo_spectral_base () |
| Sets the standard spectal bases of decomposition for each component for a pseudo_vector. | |
| virtual const Scalar & | eta () const |
Gives the field ![]() ![]() | |
| virtual const Scalar & | mu () const |
Gives the field ![]() ![]() | |
| virtual const Scalar & | A () const |
Gives the field ![]() | |
| void | update_vtvp () |
Computes the components ![]() ![]() ![]() ![]() | |
| const Scalar & | divergence (const Metric &) const |
The divergence of this with respect to a Metric . | |
| const Vector_divfree | curl () const |
The curl of this with respect to a (flat) Metric . | |
| Vector | derive_lie (const Vector &v) const |
Computes the Lie derivative of this with respect to some vector field v. | |
| Sym_tensor | ope_killing (const Metric &gam) const |
| Computes the Killing operator associated with a given metric. | |
| Sym_tensor | ope_killing_conf (const Metric &gam) const |
| Computes the conformal Killing operator associated with a given metric. | |
| Vector | poisson (double lambda, int method=6) const |
Solves the vector Poisson equation with *this as a source. | |
| Vector | poisson (double lambda, const Metric_flat &met_f, int method=6) const |
Solves the vector Poisson equation with *this as a source. | |
| Vector | poisson (const double lambda, Param &par, int method=6) const |
Solves the vector Poisson equation with *this as a source and parameters controlling the solution. | |
| void | poisson_boundary (double lambda, const Mtbl_cf &limit_vr, const Mtbl_cf &limit_eta, const Mtbl_cf &limit_mu, int num_front, double fact_dir, double fact_neu, Vector &resu) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere. | |
| void | poisson_boundary2 (double lam, Vector &resu, Scalar boundvr, Scalar boundeta, Scalar boundmu, double dir_vr, double neum_vr, double dir_eta, double neum_eta, double dir_mu, double neum_mu) const |
| Alternative to previous poisson_boundary method for vectors ; this uses method 6 for vectorial solving, updated version (as in the poisson_vector_block routine). | |
| Vector | poisson_dirichlet (double lambda, const Valeur &limit_vr, const Valeur &limit_vt, const Valeur &limit_vp, int num_front) const |
| Vector | poisson_neumann (double lambda, const Valeur &limit_vr, const Valeur &limit_vt, const Valeur &limit_vp, int num_front) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere. | |
| Vector | poisson_robin (double lambda, const Valeur &limit_vr, const Valeur &limit_vt, const Valeur &limit_vp, double fact_dir, double fact_neu, int num_front) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere. | |
| double | flux (double radius, const Metric &met) const |
| Computes the flux of the vector accross a sphere r = const. | |
| void | poisson_block (double lambda, Vector &resu) const |
| void | visu_arrows (double xmin, double xmax, double ymin, double ymax, double zmin, double zmax, const char *title0=0x0, const char *filename0=0x0, bool start_dx=true, int nx=8, int ny=8, int nz=8) const |
| 3D visualization via OpenDX. | |
| void | visu_streamline (double xmin, double xmax, double ymin, double ymax, double zmin, double zmax, const char *title0=0x0, const char *filename0=0x0, bool start_dx=true, int nx=8, int ny=8, int nz=8) const |
| 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. | |
| 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_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). | |
| const Tensor & | derive_cov (const Metric &gam) const |
Returns the covariant derivative of this with respect to some metric ![]() | |
| const Tensor & | derive_con (const Metric &gam) const |
Returns the "contravariant" derivative of this with respect to some metric ![]() | |
| 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 | sauve (FILE *) const |
| Save in a binary file. | |
| 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) const |
Logical destructor of the derivatives depending on the i-th element of met_depend in the class Vector. | |
| void | set_der_met_0x0 (int) const |
Sets all the i-th components of met_depend in the class Vector (p_potential , etc...) to 0x0. | |
| 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 | |
| Scalar * | p_potential [N_MET_MAX] |
The potential ![]() ![]() | |
| Vector_divfree * | p_div_free [N_MET_MAX] |
The divergence-free vector ![]() ![]() | |
| Scalar * | p_eta |
Field ![]() ![]() | |
| Scalar * | p_mu |
Field ![]() ![]() | |
| Scalar * | p_A |
Field ![]() | |
| 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. | |
Standard constructor.
| map | the mapping |
| tipe | the type COV for a covariant vector (1-form) and CON for a contravariant one |
| triad_i | vectorial basis (triad) with respect to which the vector components are defined |
Definition at line 156 of file vector.C.
References Lorene::Map(), set_der_0x0(), and Lorene::Tensor::Tensor().
Standard constructor with the triad passed as a pointer.
| map | the mapping |
| tipe | the type COV for a covariant vector (1-form) and CON for a contravariant one |
| triad_i | pointer on the vectorial basis (triad) with respect to which the vector components are defined |
Definition at line 165 of file vector.C.
References Lorene::Map(), set_der_0x0(), and Lorene::Tensor::Tensor().
| Lorene::Vector::Vector | ( | const Vector & | a | ) |
Copy constructor.
Definition at line 173 of file vector.C.
References set_der_0x0(), Lorene::Tensor::Tensor(), Lorene::Tensor::valence, and Vector().
| Lorene::Vector::Vector | ( | const Tensor & | a | ) |
Constructor from a Tensor .
The Tensor must be of valence one.
Definition at line 184 of file vector.C.
References set_der_0x0(), Lorene::Tensor::Tensor(), and Lorene::Tensor::valence.
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. |
| fich | file which has been created by the function sauve(FILE*) . |
Definition at line 194 of file vector.C.
References Lorene::Map(), Lorene::Tensor::n_comp, set_der_0x0(), Lorene::Tensor::Tensor(), and Lorene::Tensor::valence.
|
virtual |
|
virtual |
Gives the field 
![\[ A = {\partial \eta \over \partial r} + { \eta \over r} - {V^r \over r}
\]](form_967.png)
Related to the curl, A is insensitive to the longitudinal part of the vector.
Definition at line 129 of file vector_etamu.C.
References Lorene::Tensor::cmp, Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::dsdr(), eta(), p_A, p_eta, Lorene::Scalar::set_dzpuis(), and Lorene::Tensor::triad.
|
virtual |
Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
Reimplemented from Lorene::Tensor.
Definition at line 75 of file vector_change_triad.C.
References Lorene::Tensor::cmp, Lorene::Base_vect_cart::get_align(), Lorene::Tensor::mp, set(), Lorene::Tensor::triad, and Vector().
| const Vector_divfree Lorene::Vector::curl | ( | ) | const |
The curl of this with respect to a (flat) Metric .
The Vector is assumed to be contravariant.
Definition at line 402 of file vector.C.
References Lorene::Tensor::cmp, Lorene::Scalar::div_r(), Lorene::Scalar::div_tant(), Lorene::Scalar::dsdr(), Lorene::Tensor::mp, Lorene::Scalar::mult_r(), pseudo_spectral_base(), set(), and Lorene::Tensor::triad.
| void Lorene::Vector::decompose_div | ( | const Metric & | metre | ) | const |
Makes the Helmholtz decomposition (see documentation of p_potential ) of this with respect to a given Metric , only in the case of contravariant vectors.
Definition at line 516 of file vector.C.
References Lorene::Scalar::annule(), Lorene::Tensor::annule_domain(), Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Tensor::cmp, Lorene::Valeur::coef(), Lorene::Tensor::dec_dzpuis(), divergence(), Lorene::Map_af::get_alpha(), Lorene::Map_af::get_beta(), Lorene::Scalar::get_etat(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Tensor::get_place_met(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Matrice::inverse(), Lorene::Tensor::mp, p_div_free, p_potential, Lorene::Scalar::poisson(), Lorene::pow(), R_CHEB, R_CHEBI, R_CHEBP, Lorene::Matrice::set(), Lorene::Mtbl_cf::set(), Lorene::Tbl::set(), Lorene::Matrice::set_band(), Lorene::Valeur::set_base_r(), Lorene::Tensor::set_dependance(), Lorene::Scalar::set_dzpuis(), Lorene::Matrice::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Matrice::set_lu(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::std_spectral_base(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, Vector(), Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
|
protectedvirtual |
Deletes the derived quantities.
Reimplemented from Lorene::Tensor.
Reimplemented in Lorene::Vector_divfree.
Definition at line 219 of file vector.C.
References Lorene::Tensor::del_deriv(), del_derive_met(), p_A, p_eta, p_mu, and set_der_0x0().
|
protectedvirtual |
Logical destructor of the derivatives depending on the i-th element of met_depend in the class Vector.
Reimplemented from Lorene::Tensor.
Definition at line 242 of file vector.C.
References Lorene::Tensor::del_derive_met(), Lorene::Tensor::met_depend, p_div_free, p_potential, and set_der_met_0x0().
Computes the Lie derivative of this with respect to some vector field v.
Definition at line 392 of file vector.C.
References Lorene::Tensor::compute_derive_lie(), Lorene::Tensor::mp, Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Vector().
| const Vector_divfree & Lorene::Vector::div_free | ( | const Metric & | metre | ) | const |
Returns the div-free vector in the Helmholtz decomposition.
It first makes the Helmholtz decomposition (see documentation of p_potential ) of this with respect to a given Metric and then returns 
Definition at line 504 of file vector.C.
References decompose_div(), Lorene::Tensor::get_place_met(), p_div_free, and Lorene::Tensor::set_dependance().
| const Scalar & Lorene::Vector::divergence | ( | const Metric & | metre | ) | const |
The divergence of this with respect to a Metric .
The Vector is assumed to be contravariant.
Definition at line 381 of file vector.C.
References Lorene::Tensor::divergence().
|
virtual |
Gives the field 

![\[ V^\theta = {\partial \eta \over \partial\theta} -
{1\over\sin\theta} {\partial \mu \over \partial\varphi}
\]](form_961.png)
![\[ V^\varphi = {1\over\sin\theta}
{\partial \eta \over \partial\varphi}
+ {\partial \mu \over \partial\theta}
\]](form_960.png)
Reimplemented in Lorene::Vector_divfree.
Definition at line 66 of file vector_etamu.C.
References Lorene::Tensor::cmp, Lorene::Scalar::div_tant(), p_eta, Lorene::Scalar::poisson_angu(), 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 r-component, 

Reimplemented from Lorene::Tensor.
Definition at line 850 of file vector.C.
References Lorene::Tensor::cmp, eta(), Lorene::Scalar::exponential_filter_r(), exponential_filter_r(), Lorene::Tensor::mp, mu(), Lorene::Tensor::n_comp, operator()(), set_vr_eta_mu(), and Lorene::Tensor::triad.
|
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, 

Reimplemented from Lorene::Tensor.
Definition at line 867 of file vector.C.
References Lorene::Tensor::cmp, eta(), Lorene::Scalar::exponential_filter_ylm(), exponential_filter_ylm(), Lorene::Tensor::mp, mu(), Lorene::Tensor::n_comp, operator()(), set_vr_eta_mu(), and Lorene::Tensor::triad.
| double Lorene::Vector::flux | ( | double | radius, |
| const Metric & | met ) const |
Computes the flux of the vector accross a sphere r = const.
| radius | radius of the sphere S on which the flux is to be taken; the center of S is assumed to be the center of the mapping (member mp). radius can take the value __infinity (to get the flux at spatial infinity). |
| met | metric ![]() |


*this and 

Definition at line 807 of file vector.C.
References change_triad(), Lorene::Map_af::integrale_surface(), Lorene::Map_af::integrale_surface_infini(), Lorene::Tensor::mp, Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Vector().
|
inlinevirtual |
Returns the index of a component given by its position in the Scalar array cmp .
Itbl ) of size 1 giving its value for the component located at the position place in the Scalar array cmp . The element of this Itbl Reimplemented from Lorene::Tensor.
|
virtual |
Gives the field 

![\[ V^\theta = {\partial \eta \over \partial\theta} -
{1\over\sin\theta} {\partial \mu \over \partial\varphi}
\]](form_966.png)
![\[ V^\varphi = {1\over\sin\theta}
{\partial \eta \over \partial\varphi}
+ {\partial \mu \over \partial\theta}
\]](form_960.png)
Definition at line 98 of file vector_etamu.C.
References Lorene::Tensor::cmp, Lorene::Scalar::div_tant(), p_mu, Lorene::Scalar::poisson_angu(), and Lorene::Tensor::triad.
| Sym_tensor Lorene::Vector::ope_killing | ( | const Metric & | gam | ) | const |
Computes the Killing operator associated with a given metric.
The Killing operator is defined by 

| gam | metric with respect to which the covariant derivative ![]() |
Definition at line 438 of file vector.C.
References Lorene::Tensor::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Tensor::mp, Lorene::Tensor::set(), Lorene::Tensor::triad, and Lorene::Tensor::type_indice.
| Sym_tensor Lorene::Vector::ope_killing_conf | ( | const Metric & | gam | ) | const |
Computes the conformal Killing operator associated with a given metric.
The conformal Killing operator is defined by 

| gam | metric ![]() ![]() |
Definition at line 462 of file vector.C.
References Lorene::Metric::con(), Lorene::Metric::cov(), Lorene::Tensor::derive_con(), Lorene::Tensor::derive_cov(), divergence(), Lorene::Tensor::mp, Lorene::Tensor::set(), Lorene::Tensor::trace(), Lorene::Tensor::triad, and Lorene::Tensor::type_indice.
| const Scalar & Lorene::Vector::operator() | ( | int | index | ) | const |
Readonly access to a component.
Definition at line 305 of file vector.C.
References Lorene::Tensor::cmp.
|
virtual |
Assignment from a Tensor.
Reimplemented from Lorene::Tensor.
Reimplemented in Lorene::Vector_divfree.
Definition at line 279 of file vector.C.
References Lorene::Tensor::cmp, del_deriv(), Lorene::Tensor::Tensor(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Lorene::Tensor::valence.
|
virtual |
Assignment from a Vector.
Reimplemented in Lorene::Vector_divfree.
Definition at line 267 of file vector.C.
References Lorene::Tensor::cmp, del_deriv(), Lorene::Tensor::triad, Lorene::Tensor::type_indice, and Vector().
Solves the vector Poisson equation with *this as a source and parameters controlling the solution.
The equatiopn solved is 
*this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential ), with a flat metric, deduced from the triad.
| lambda | [input] ![]() |
| par | [input/output] possible parameters |
| uu | [input/output] solution u with the boundary condition u =0 at spatial infinity. |
Definition at line 543 of file vector_poisson.C.
References Lorene::Param::add_cmp_mod(), Lorene::Param::add_double(), Lorene::Param::add_int(), Lorene::Param::add_int_mod(), Lorene::Cmp::annule_hard(), Lorene::Tenseur::change_triad(), Lorene::Tensor::cmp, Lorene::Tensor::dec_dzpuis(), Lorene::Scalar::derive_con(), div_free(), Lorene::Param::get_cmp_mod(), Lorene::Param::get_double(), Lorene::Param::get_int(), Lorene::Param::get_int_mod(), Lorene::Scalar::inc_dzpuis(), Lorene::Tensor::mp, Lorene::Scalar::poisson(), Lorene::Vector_divfree::poisson(), Lorene::Tenseur::poisson_vect(), potential(), Lorene::Tenseur::set(), set(), Lorene::Tenseur::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Tenseur::set_std_base(), Lorene::Tensor::triad, and Vector().
| Vector Lorene::Vector::poisson | ( | double | lambda, |
| const Metric_flat & | met_f, | ||
| int | method = 6 ) const |
Solves the vector Poisson equation with *this as a source.
The equation solved is 
*this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential ), with the flat metric met_f given in argument.
| lambda | [input] ![]() |
| met_f | [input] the flat metric for the Helmholtz decomposition. |
| method | [input] method used to solve the equation:
|

Definition at line 130 of file vector_poisson.C.
References Lorene::Scalar::annule_hard(), Lorene::Valeur::base, Lorene::Valeur::c, Lorene::Valeur::c_cf, Lorene::Tenseur::change_triad(), Lorene::Tensor::cmp, Lorene::Valeur::coef(), Lorene::Scalar::dec_dzpuis(), Lorene::Tensor::dec_dzpuis(), Lorene::Scalar::derive_con(), div_free(), Lorene::Scalar::div_r_dzpuis(), Lorene::Scalar::div_sint(), Lorene::Scalar::dsdr(), Lorene::Scalar::dsdt(), Lorene::Mtbl_cf::get_etat(), Lorene::Scalar::get_etat(), Lorene::Scalar::lapang(), Lorene::Tensor::mp, mu(), Lorene::Scalar::mult_r_dzpuis(), Lorene::Scalar::mult_sint(), Lorene::Scalar::poisson(), poisson(), Lorene::poisson_angu(), Lorene::Scalar::poisson_angu(), Lorene::Tenseur::poisson_vect(), Lorene::Tenseur::poisson_vect_oohara(), potential(), Lorene::Scalar::primr(), Lorene::Mtbl_cf::set(), Lorene::Tenseur::set(), set(), Lorene::Scalar::set_dzpuis(), Lorene::Scalar::set_etat_qcq(), Lorene::Tenseur::set_etat_qcq(), Lorene::Scalar::set_etat_zero(), Lorene::Tensor::set_etat_zero(), Lorene::Param_elliptic::set_poisson_vect_eta(), Lorene::Param_elliptic::set_poisson_vect_r(), Lorene::Scalar::set_spectral_va(), set_vr_eta_mu(), Lorene::Scalar::sol_elliptic(), Lorene::Tensor::triad, Vector(), Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| Vector Lorene::Vector::poisson | ( | double | lambda, |
| int | method = 6 ) const |
Solves the vector Poisson equation with *this as a source.
The equation solved is 
*this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential ), with a flat metric, deduced from the triad.
| lambda | [input] ![]() |
| method | [input] method used to solve the equation (see Vector::poisson(double, Metric_flat, int) for details). |

Definition at line 518 of file vector_poisson.C.
References Lorene::Tensor::mp, poisson(), Lorene::Tensor::triad, and Vector().
| void Lorene::Vector::poisson_block | ( | double | lambda, |
| Vector & | resu ) const |
Definition at line 75 of file vector_poisson_block.C.
| void Lorene::Vector::poisson_boundary | ( | double | lambda, |
| const Mtbl_cf & | limit_vr, | ||
| const Mtbl_cf & | limit_eta, | ||
| const Mtbl_cf & | limit_mu, | ||
| int | num_front, | ||
| double | fact_dir, | ||
| double | fact_neu, | ||
| Vector & | resu ) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere.
The equation solved is 
*this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential )
| lambda | [input] ![]() |
| resu | [output] the solution ![]() |
Definition at line 65 of file vector_poisson_boundary.C.
References Lorene::Mtbl_cf::annule_hard(), Lorene::Scalar::annule_hard(), Lorene::Tbl::annule_hard(), Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Tensor::cmp, eta(), Lorene::Map_af::get_alpha(), Lorene::Mg3d::get_angu(), Lorene::Map_af::get_beta(), Lorene::Scalar::get_etat(), Lorene::Diff_dsdx2::get_matrice(), Lorene::Diff_dsdx::get_matrice(), Lorene::Diff_x2dsdx2::get_matrice(), Lorene::Diff_xdsdx2::get_matrice(), Lorene::Diff_xdsdx::get_matrice(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Base_val::give_quant_numbers(), Lorene::Matrice::inverse(), Lorene::Tensor::mp, mu(), Lorene::pow(), Lorene::Matrice::set(), Lorene::Mtbl_cf::set(), Lorene::Tbl::set(), Lorene::Matrice::set_band(), Lorene::Valeur::set_etat_cf_qcq(), Lorene::Matrice::set_etat_qcq(), Lorene::Scalar::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Matrice::set_lu(), Lorene::Param_elliptic::set_poisson_vect_r(), Lorene::Scalar::set_spectral_base(), Lorene::Scalar::set_spectral_va(), set_vr_eta_mu(), Lorene::Scalar::sol_elliptic_boundary(), Lorene::Tbl::t, Lorene::Tensor::triad, Vector(), Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| void Lorene::Vector::poisson_boundary2 | ( | double | lam, |
| Vector & | resu, | ||
| Scalar | boundvr, | ||
| Scalar | boundeta, | ||
| Scalar | boundmu, | ||
| double | dir_vr, | ||
| double | neum_vr, | ||
| double | dir_eta, | ||
| double | neum_eta, | ||
| double | dir_mu, | ||
| double | neum_mu ) const |
Alternative to previous poisson_boundary method for vectors ; this uses method 6 for vectorial solving, updated version (as in the poisson_vector_block routine).
Boundary arguments are here required as scalar fields.
Definition at line 80 of file vector_poisson_boundary2.C.
References Lorene::Matrice::annule_hard(), Lorene::Mtbl_cf::annule_hard(), Lorene::Scalar::annule_hard(), Lorene::Tbl::annule_hard(), Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Tensor::cmp, eta(), Lorene::Map_af::get_alpha(), Lorene::Map_af::get_beta(), Lorene::Scalar::get_etat(), Lorene::Diff_dsdx2::get_matrice(), Lorene::Diff_dsdx::get_matrice(), Lorene::Diff_sx2::get_matrice(), Lorene::Diff_sxdsdx::get_matrice(), Lorene::Diff_x2dsdx2::get_matrice(), Lorene::Diff_xdsdx2::get_matrice(), Lorene::Diff_xdsdx::get_matrice(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Scalar::get_spectral_base(), Lorene::Scalar::get_spectral_va(), Lorene::Mg3d::get_type_r(), Lorene::Base_val::give_quant_numbers(), Lorene::Matrice::inverse(), Lorene::Tensor::mp, mu(), Lorene::Matrice::set(), Lorene::Mtbl_cf::set(), Lorene::Tbl::set(), Lorene::Scalar::set_dzpuis(), Lorene::Valeur::set_etat_cf_qcq(), Lorene::Matrice::set_etat_qcq(), Lorene::Scalar::set_etat_qcq(), Lorene::Tbl::set_etat_qcq(), Lorene::Matrice::set_lu(), Lorene::Param_elliptic::set_poisson_vect_r(), Lorene::Scalar::set_spectral_base(), Lorene::Scalar::set_spectral_va(), set_vr_eta_mu(), Lorene::Scalar::sol_elliptic_boundary(), Lorene::Scalar::std_spectral_base(), Lorene::Tbl::t, Lorene::Tensor::triad, Lorene::Mtbl_cf::val_in_bound_jk(), Lorene::Mtbl_cf::val_out_bound_jk(), Vector(), Lorene::Valeur::ylm(), and Lorene::Valeur::ylm_i().
| Vector Lorene::Vector::poisson_dirichlet | ( | double | lambda, |
| const Valeur & | limit_vr, | ||
| const Valeur & | limit_vt, | ||
| const Valeur & | limit_vp, | ||
| int | num_front ) const |
Definition at line 680 of file vector_poisson_boundary.C.
| Vector Lorene::Vector::poisson_neumann | ( | double | lambda, |
| const Valeur & | limit_vr, | ||
| const Valeur & | limit_vt, | ||
| const Valeur & | limit_vp, | ||
| int | num_front ) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere.
The equation solved is 
*this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential )
| lambda | [input] ![]() |
| resu | [output] the solution ![]() |
Definition at line 691 of file vector_poisson_boundary.C.
References Lorene::Tensor::mp, poisson_robin(), Lorene::Tensor::triad, and Vector().
| Vector Lorene::Vector::poisson_robin | ( | double | lambda, |
| const Valeur & | limit_vr, | ||
| const Valeur & | limit_vt, | ||
| const Valeur & | limit_vp, | ||
| double | fact_dir, | ||
| double | fact_neu, | ||
| int | num_front ) const |
Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere.
The equation solved is 
*this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential )
| lambda | [input] ![]() |
| resu | [output] the solution ![]() |
Definition at line 702 of file vector_poisson_boundary.C.
References Lorene::Scalar::annule_hard(), Lorene::Valeur::base, Lorene::Valeur::c_cf, Lorene::Tensor::cmp, Lorene::Valeur::coef(), Lorene::Scalar::div_tant(), Lorene::Scalar::dsdt(), Lorene::Valeur::get_base(), Lorene::Valeur::get_etat(), Lorene::Scalar::get_spectral_va(), Lorene::Tensor::mp, Lorene::norme(), Lorene::Scalar::poisson_angu(), poisson_boundary(), Lorene::Valeur::set(), Lorene::Valeur::set_base(), Lorene::Tensor::set_etat_zero(), Lorene::Scalar::set_grid_point(), Lorene::Scalar::set_spectral_va(), Lorene::Scalar::stdsdp(), Lorene::Tensor::triad, Lorene::Scalar::val_grid_point(), Vector(), and Lorene::Valeur::ylm().
|
inlinevirtual |
Returns the position in the Scalar array cmp of a component given by its index.
Scalar array cmp idx . idx must be a 1-D Itbl of size 1, the element of which must be one of the spatial indices 1, 2 or 3. Reimplemented from Lorene::Tensor.
Definition at line 392 of file vector.h.
References Lorene::Itbl::get_dim(), and Lorene::Itbl::get_ndim().
| const Scalar & Lorene::Vector::potential | ( | const Metric & | metre | ) | const |
Returns the potential in the Helmholtz decomposition.
It first makes the Helmholtz decomposition (see documentation of p_potential ) of this with respect to a given Metric and then returns 
Definition at line 492 of file vector.C.
References decompose_div(), Lorene::Tensor::get_place_met(), p_potential, and Lorene::Tensor::set_dependance().
|
virtual |
Sets the standard spectal bases of decomposition for each component for a pseudo_vector.
Definition at line 346 of file vector.C.
References Lorene::Tensor::cmp, Lorene::Tensor::mp, and Lorene::Tensor::triad.
| Scalar & Lorene::Vector::set | ( | int | index | ) |
Read/write access to a component.
Definition at line 296 of file vector.C.
References Lorene::Tensor::cmp, and del_deriv().
|
protected |
|
protected |
Sets all the i-th components of met_depend in the class Vector (p_potential , etc...) to 0x0.
Definition at line 258 of file vector.C.
References p_div_free, and p_potential.
| void Lorene::Vector::set_vr_eta_mu | ( | const Scalar & | vr_i, |
| const Scalar & | eta_i, | ||
| const Scalar & | mu_i ) |
Defines the components through potentials 

p_eta and p_mu ), as well as the 
| vr_i | [input] component ![]() |
| eta_i | [input] angular potential ![]() |
| mu_i | [input] angular potential ![]() |
Definition at line 189 of file vector_etamu.C.
References Lorene::Tensor::cmp, del_deriv(), Lorene::Tensor::get_mp(), p_eta, p_mu, Lorene::Tensor::triad, and update_vtvp().
|
virtual |
Sets the standard spectal bases of decomposition for each component.
Reimplemented from Lorene::Tensor.
Definition at line 316 of file vector.C.
References Lorene::Tensor::cmp, Lorene::Tensor::mp, and Lorene::Tensor::triad.
| void Lorene::Vector::update_vtvp | ( | ) |
Computes the components 



![\[ V^\theta = {\partial \eta \over \partial\theta} -
{1\over\sin\theta} {\partial \mu \over \partial\varphi}
\]](form_970.png)
![\[ V^\varphi = {1\over\sin\theta}
{\partial \eta \over \partial\varphi}
+ {\partial \mu \over \partial\theta}
\]](form_960.png)
Definition at line 167 of file vector_etamu.C.
References Lorene::Tensor::cmp, del_deriv(), p_eta, and p_mu.
| void Lorene::Vector::visu_arrows | ( | double | xmin, |
| double | xmax, | ||
| double | ymin, | ||
| double | ymax, | ||
| double | zmin, | ||
| double | zmax, | ||
| const char * | title0 = 0x0, | ||
| const char * | filename0 = 0x0, | ||
| bool | start_dx = true, | ||
| int | nx = 8, | ||
| int | ny = 8, | ||
| int | nz = 8 ) const |
3D visualization via OpenDX.
| xmin | [input] defines with xmax the x range of the visualization box |
| xmax | [input] defines with xmin the x range of the visualization box |
| ymin | [input] defines with ymax the y range of the visualization box |
| ymax | [input] defines with ymin the y range of the visualization box |
| zmin | [input] defines with zmax the z range of the visualization box |
| zmax | [input] defines with zmin the z range of the visualization box |
| title | [input] title for the graph (for OpenDX legend) |
| filename | [input] name for the file which will be the input for OpenDX; the default 0x0 is transformed into "vector_arrows" |
| start_dx | [input] determines whether OpenDX must be launched (as a subprocess) to view the field; if set to false , only input files for future usage of OpenDX are created |
| nx | [input] number of points in the x direction (uniform sampling) |
| ny | [input] number of points in the y direction (uniform sampling) |
| nz | [input] number of points in the z direction (uniform sampling) |
Definition at line 62 of file vector_visu.C.
References Lorene::Valeur::c_cf, change_triad(), Lorene::Scalar::check_dzpuis(), Lorene::Valeur::coef(), Lorene::Scalar::dec_dzpuis(), Lorene::Scalar::get_dzpuis(), Lorene::Tensor::mp, operator()(), set(), Lorene::Tensor::triad, Lorene::Mtbl_cf::val_point(), and Vector().
| void Lorene::Vector::visu_streamline | ( | double | xmin, |
| double | xmax, | ||
| double | ymin, | ||
| double | ymax, | ||
| double | zmin, | ||
| double | zmax, | ||
| const char * | title0 = 0x0, | ||
| const char * | filename0 = 0x0, | ||
| bool | start_dx = true, | ||
| int | nx = 8, | ||
| int | ny = 8, | ||
| int | nz = 8 ) const |
Definition at line 283 of file vector_visu.C.
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |
|
mutableprotected |