LORENE
Lorene::Vector_divfree Class Reference

Divergence-free vectors. More...

#include <vector.h>

Inheritance diagram for Lorene::Vector_divfree:
Lorene::Vector Lorene::Tensor

Public Member Functions

 Vector_divfree (const Map &map, const Base_vect &triad_i, const Metric &met)
 Standard constructor.
 
 Vector_divfree (const Vector_divfree &)
 Copy constructor.
 
 Vector_divfree (const Map &map, const Base_vect &triad_i, const Metric &met, FILE *fich)
 Constructor from a file (see Tensor::sauve(FILE*) ).
 
virtual ~Vector_divfree ()
 Destructor.
 
void operator= (const Vector_divfree &a)
 Assignment from another Vector_divfree.
 
virtual void operator= (const Vector &a)
 Assignment from a Vector.
 
virtual void operator= (const Tensor &a)
 Assignment from a Tensor.
 
void set_vr_mu (const Scalar &vr_i, const Scalar &mu_i)
 Sets the angular potentials $\mu$ (see member p_mu ), and the $V^r$ component of the vector.
 
void set_vr_eta_mu (const Scalar &vr_i, const Scalar &eta_i, const Scalar &mu_i)
 Defines the components through $V_r$, $\eta$ and $\mu$.
 
void set_A_mu (const Scalar &A_i, const Scalar &mu_i, const Param *par_bc)
 Defines the components through potentials $A$ and $\mu$.
 
virtual const Scalareta () const
 Gives the field $\eta$ such that the angular components $(V^\theta, V^\varphi)$ of the vector are written:
 
Vector_divfree poisson () const
 Computes the solution of a vectorial Poisson equation with *this $= \vec{V}$ as a source:
 
void update_etavr ()
 Computes the components $V^r$ and $\eta$ from the potential A and the divergence-free condition, according to :
 
Vector_divfree poisson (Param &par) const
 Computes the solution of a vectorial Poisson equation with *this $= \vec{V}$ as a source:
 
virtual void change_triad (const Base_vect &)
 Sets a new vectorial basis (triad) of decomposition and modifies the components accordingly.
 
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 Scalarpotential (const Metric &) const
 Returns the potential in the Helmholtz decomposition.
 
const Vector_divfreediv_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 ).
 
Scalarset (int)
 Read/write access to a component.
 
Scalarset (const Itbl &ind)
 Returns the value of a component (read/write version).
 
Scalarset (int i1, int i2)
 Returns the value of a component for a tensor of valence 2 (read/write version).
 
Scalarset (int i1, int i2, int i3)
 Returns the value of a component for a tensor of valence 3 (read/write version).
 
Scalarset (int i1, int i2, int i3, int i4)
 Returns the value of a component for a tensor of valence 4 (read/write version).
 
const Scalaroperator() (int) const
 Readonly access to a component.
 
const Scalaroperator() (const Itbl &ind) const
 Returns the value of a component (read-only version).
 
const Scalaroperator() (int i1, int i2) const
 Returns the value of a component for a tensor of valence 2 (read-only version).
 
const Scalaroperator() (int i1, int i2, int i3) const
 Returns the value of a component for a tensor of valence 3 (read-only version).
 
const Scalaroperator() (int i1, int i2, int i3, int i4) const
 Returns the value of a component for a tensor of valence 4 (read-only version).
 
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 Scalarmu () const
 Gives the field $\mu$ such that the angular components $(V^\theta, V^\varphi)$ of the vector are written:
 
virtual const ScalarA () const
 Gives the field $A$ defined by.
 
void update_vtvp ()
 Computes the components $V^\theta$ and $V^\varphi$ from the potential $\eta$ and $\mu$, according to:
 
const Scalardivergence (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.
 
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 Tensorderive_cov (const Metric &gam) const
 Returns the covariant derivative of this with respect to some metric $\gamma$.
 
const Tensorderive_con (const Metric &gam) const
 Returns the "contravariant" derivative of this with respect to some metric $\gamma$, by raising the last index of the covariant derivative (cf.
 
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 Mapget_mp () const
 Returns the mapping.
 
const Base_vectget_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.
 
intset_index_type (int i)
 Sets the type of the index number i .
 
Itblset_index_type ()
 Sets the types of all the indices.
 
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.
 
void sol_Dirac_A (const Scalar &aaa, Scalar &eta, Scalar &vr, const Param *par_bc=0x0) const
 Solves a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.
 
void sol_Dirac_A_tau (const Scalar &aaa, Scalar &eta, Scalar &vr, const Param *par_bc=0x0) const
 Solves via a tau method a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.
 
void sol_Dirac_A_poisson (const Scalar &aaa, Scalar &eta, Scalar &vr, const Param *par_bc=0x0) const
 Solves via a poisson method a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.
 
void sol_Dirac_A_1z (const Scalar &aaa, Scalar &eta, Scalar &vr, const Param *par_bc=0x0) const
 Solves a one-domain system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.
 
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

const Metric *const met_div
 Metric with respect to which the divergence is defined.
 
Scalarp_potential [N_MET_MAX]
 The potential $\phi$ giving the gradient part in the Helmholtz decomposition of any 3D vector $\vec{V}: \quad \vec{V} = 
\vec{\nabla} \phi + \vec{\nabla} \wedge \vec{\psi}$.
 
Vector_divfreep_div_free [N_MET_MAX]
 The divergence-free vector $\vec{W} =  \vec{\nabla} \wedge 
  \vec{\psi}$ of the Helmholtz decomposition of any 3D vector $\vec{V}: \quad \vec{V} = \vec{\nabla} \phi + \vec{\nabla} 
*\wedge \vec{\psi}$.
 
Scalarp_eta
 Field $\eta$ such that the angular components $(V^\theta, V^\varphi)$ of the vector are written:
 
Scalarp_mu
 Field $\mu$ such that the angular components $(V^\theta, V^\varphi)$ of the vector are written:
 
Scalarp_A
 Field $A$ defined by.
 
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_vecttriad
 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 Metricmet_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.
 
Tensorp_derive_cov [N_MET_MAX]
 Array of pointers on the covariant derivatives of this with respect to various metrics.
 
Tensorp_derive_con [N_MET_MAX]
 Array of pointers on the contravariant derivatives of this with respect to various metrics.
 
Tensorp_divergence [N_MET_MAX]
 Array of pointers on the divergence of this with respect to various metrics.
 

Detailed Description

Divergence-free vectors.

()

This class is designed to store divergence-free vectors, with the component expressed in a orthonormal spherical basis $(e_r,e_\theta,e_\varphi)$.

Definition at line 724 of file vector.h.

Constructor & Destructor Documentation

◆ Vector_divfree() [1/3]

Lorene::Vector_divfree::Vector_divfree ( const Map map,
const Base_vect triad_i,
const Metric met 
)

Standard constructor.

Parameters
mapthe mapping
triad_ivectorial basis (triad) with respect to which the vector components are defined
metthe metric with respect to which the divergence is defined

Definition at line 79 of file vector_divfree.C.

References set_der_0x0().

◆ Vector_divfree() [2/3]

Lorene::Vector_divfree::Vector_divfree ( const Vector_divfree source)

Copy constructor.

Definition at line 91 of file vector_divfree.C.

References Lorene::Vector::p_eta, Lorene::Vector::p_mu, and set_der_0x0().

◆ Vector_divfree() [3/3]

Lorene::Vector_divfree::Vector_divfree ( const Map map,
const Base_vect triad_i,
const Metric met,
FILE fich 
)

Constructor from a file (see Tensor::sauve(FILE*) ).

Parameters
mapthe mapping
triad_ivectorial basis (triad) with respect to which the tensor components are defined. It will be checked that it coincides with the basis saved in the file.
metthe metric with respect to which the divergence is defined
fichfile which has been created by the function sauve(FILE*) .

Definition at line 103 of file vector_divfree.C.

References set_der_0x0().

◆ ~Vector_divfree()

Lorene::Vector_divfree::~Vector_divfree ( )
virtual

Destructor.

Definition at line 117 of file vector_divfree.C.

References del_deriv().

Member Function Documentation

◆ A()

const Scalar & Lorene::Vector::A ( ) const
virtualinherited

Gives the field $A$ defined by.

\[ 
    A = {\partial \eta \over \partial r} + { \eta \over r} - {V^r \over r}
 \]

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::Vector::eta(), Lorene::Vector::p_A, Lorene::Vector::p_eta, and Lorene::Tensor::triad.

◆ change_triad()

void Lorene::Vector::change_triad ( const Base_vect new_triad)
virtualinherited

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::Tensor::mp, Lorene::Vector::set(), and Lorene::Tensor::triad.

◆ curl()

const Vector_divfree Lorene::Vector::curl ( ) const
inherited

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::dsdt(), Lorene::Scalar::dsdx(), Lorene::Scalar::dsdy(), Lorene::Scalar::dsdz(), Lorene::Tensor::mp, Lorene::Scalar::srstdsdp(), and Lorene::Tensor::triad.

◆ decompose_div()

◆ del_deriv()

void Lorene::Vector_divfree::del_deriv ( ) const
protectedvirtual

Deletes the derived quantities.

Reimplemented from Lorene::Vector.

Definition at line 128 of file vector_divfree.C.

References Lorene::Vector::del_deriv(), and set_der_0x0().

◆ del_derive_met()

void Lorene::Vector::del_derive_met ( int  j) const
protectedvirtualinherited

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, Lorene::Vector::p_div_free, Lorene::Vector::p_potential, and Lorene::Vector::set_der_met_0x0().

◆ derive_lie()

Vector Lorene::Vector::derive_lie ( const Vector v) const
inherited

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, and Lorene::Tensor::type_indice.

◆ div_free()

const Vector_divfree & Lorene::Vector::div_free ( const Metric metre) const
inherited

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 $\vec{W}$. Only in the case of contravariant vectors.

Definition at line 504 of file vector.C.

References Lorene::Vector::decompose_div(), Lorene::Tensor::get_place_met(), Lorene::Vector::p_div_free, and Lorene::Tensor::set_dependance().

◆ divergence()

const Scalar & Lorene::Vector::divergence ( const Metric metre) const
inherited

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().

◆ eta()

const Scalar & Lorene::Vector_divfree::eta ( ) const
virtual

Gives the field $\eta$ such that the angular components $(V^\theta, V^\varphi)$ of the vector are written:

\[
 V^\theta =  {\partial \eta \over \partial\theta} -
  {1\over\sin\theta} {\partial \mu \over \partial\varphi} 
\]

\[
 V^\varphi =  {1\over\sin\theta} 
            {\partial \eta \over \partial\varphi}
            + {\partial \mu \over \partial\theta} 
\]

Reimplemented from Lorene::Vector.

Definition at line 193 of file vector_divfree.C.

References Lorene::Tensor::cmp, Lorene::Scalar::dsdr(), Lorene::Scalar::get_dzpuis(), Lorene::Vector::p_eta, and Lorene::Tensor::triad.

◆ exponential_filter_r()

void Lorene::Vector::exponential_filter_r ( int  lzmin,
int  lzmax,
int  p,
double  alpha = -16. 
)
virtualinherited

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, $\eta$ and $\mu$ for spherical components.

Reimplemented from Lorene::Tensor.

Definition at line 850 of file vector.C.

References Lorene::Tensor::cmp, Lorene::Vector::eta(), Lorene::Vector::exponential_filter_r(), Lorene::Tensor::mp, Lorene::Vector::mu(), Lorene::Tensor::n_comp, Lorene::Vector::operator()(), Lorene::Vector::set_vr_eta_mu(), and Lorene::Tensor::triad.

◆ exponential_filter_ylm()

void Lorene::Vector::exponential_filter_ylm ( int  lzmin,
int  lzmax,
int  p,
double  alpha = -16. 
)
virtualinherited

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, $\eta$ and $\mu$ for spherical components.

Reimplemented from Lorene::Tensor.

Definition at line 867 of file vector.C.

References Lorene::Tensor::cmp, Lorene::Vector::eta(), Lorene::Vector::exponential_filter_ylm(), Lorene::Tensor::mp, Lorene::Vector::mu(), Lorene::Tensor::n_comp, Lorene::Vector::operator()(), Lorene::Vector::set_vr_eta_mu(), and Lorene::Tensor::triad.

◆ flux()

double Lorene::Vector::flux ( double  radius,
const Metric met 
) const
inherited

Computes the flux of the vector accross a sphere r = const.

Parameters
radiusradius 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).
metmetric $ \gamma $ giving the area element of the sphere
Returns
$ \oint_S V^i ds_i $, where $ V^i $ is the vector represented by *this and $ ds_i $ is the area element induced on S by $ \gamma $.

Definition at line 807 of file vector.C.

References Lorene::Map_af::integrale_surface(), Lorene::Map_af::integrale_surface_infini(), Lorene::Tensor::mp, Lorene::Tensor::triad, and Lorene::Tensor::type_indice.

◆ indices()

virtual Itbl Lorene::Vector::indices ( int  place) const
inlinevirtualinherited

Returns the index of a component given by its position in the Scalar array cmp .

Returns
the index is stored in an 1-D array (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
corresponds to a spatial index 1, 2 or 3.

Reimplemented from Lorene::Tensor.

Definition at line 411 of file vector.h.

◆ mu()

const Scalar & Lorene::Vector::mu ( ) const
virtualinherited

Gives the field $\mu$ such that the angular components $(V^\theta, V^\varphi)$ of the vector are written:

\[
 V^\theta =  {\partial \eta \over \partial\theta} -
  {1\over\sin\theta} {\partial \mu \over \partial\varphi}
\]

\[
 V^\varphi =  {1\over\sin\theta} 
            {\partial \eta \over \partial\varphi}
            + {\partial \mu \over \partial\theta} 
\]

Definition at line 98 of file vector_etamu.C.

References Lorene::Tensor::cmp, Lorene::Scalar::dsdt(), Lorene::Vector::p_mu, Lorene::Scalar::stdsdp(), and Lorene::Tensor::triad.

◆ ope_killing()

Sym_tensor Lorene::Vector::ope_killing ( const Metric gam) const
inherited

Computes the Killing operator associated with a given metric.

The Killing operator is defined by $ D^i V^j + D^j V^i $ for a contravariant vector and by $ D_i V_j + D_j V_i $ for a covariant vector.

Parameters
gammetric with respect to which the covariant derivative $ D_i $ is defined.

Definition at line 438 of file vector.C.

References Lorene::Tensor::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Tensor::mp, Lorene::Tensor::triad, and Lorene::Tensor::type_indice.

◆ ope_killing_conf()

Sym_tensor Lorene::Vector::ope_killing_conf ( const Metric gam) const
inherited

Computes the conformal Killing operator associated with a given metric.

The conformal Killing operator is defined by $ D^i V^j + D^j V^i - \frac{2}{3} D_k V^k \, \gamma^{ij} $ for a contravariant vector and by $ D_i V_j + D_j V_i - \frac{2}{3} D^k V_k \, \gamma_{ij}$ for a covariant vector.

Parameters
gammetric $\gamma_{ij}$ with respect to which the covariant derivative $ D_i $ is defined.

Definition at line 462 of file vector.C.

References Lorene::Metric::con(), Lorene::Metric::cov(), Lorene::Tensor::derive_con(), Lorene::Tensor::derive_cov(), Lorene::Vector::divergence(), Lorene::Tensor::mp, Lorene::Tensor::triad, and Lorene::Tensor::type_indice.

◆ operator()()

const Scalar & Lorene::Vector::operator() ( int  index) const
inherited

Readonly access to a component.

Definition at line 305 of file vector.C.

References Lorene::Tensor::cmp.

◆ operator=() [1/3]

void Lorene::Vector_divfree::operator= ( const Tensor a)
virtual

Assignment from a Tensor.

Reimplemented from Lorene::Vector.

Definition at line 163 of file vector_divfree.C.

References del_deriv(), and Lorene::Vector::operator=().

◆ operator=() [2/3]

void Lorene::Vector_divfree::operator= ( const Vector a)
virtual

Assignment from a Vector.

Reimplemented from Lorene::Vector.

Definition at line 153 of file vector_divfree.C.

References del_deriv(), and Lorene::Vector::operator=().

◆ operator=() [3/3]

void Lorene::Vector_divfree::operator= ( const Vector_divfree a)

Assignment from another Vector_divfree.

Definition at line 141 of file vector_divfree.C.

References del_deriv(), met_div, and Lorene::Vector::operator=().

◆ poisson() [1/5]

Vector_divfree Lorene::Vector_divfree::poisson ( ) const

Computes the solution of a vectorial Poisson equation with *this $= \vec{V}$ as a source:

\[
   \Delta \vec{W} = \vec{V}
\]

Returns
solution $\vec{W}$ of the above equation with the boundary condition $\vec{W}=0$ at spatial infinity.

Definition at line 146 of file vector_df_poisson.C.

References Lorene::Tensor::cmp, eta(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Base_val::give_quant_numbers(), met_div, Lorene::Tensor::mp, Lorene::Vector::mu(), Lorene::pow(), Lorene::Tbl::t, and Lorene::Tensor::triad.

◆ poisson() [2/5]

Vector Lorene::Vector::poisson ( const double  lambda,
Param par,
int  method = 6 
) const
inherited

Solves the vector Poisson equation with *this as a source and parameters controlling the solution.

The equatiopn solved is $\Delta N^i +\lambda \nabla^i 
\nabla_k N^k = S^i$. *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.

Parameters
lambda[input] $\lambda$.
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::Tensor::cmp, Lorene::Vector::div_free(), Lorene::Tensor::mp, Lorene::Vector::potential(), and Lorene::Tensor::triad.

◆ poisson() [3/5]

Vector Lorene::Vector::poisson ( double  lambda,
const Metric_flat met_f,
int  method = 6 
) const
inherited

Solves the vector Poisson equation with *this as a source.

The equation solved is $\Delta N^i +\lambda \nabla^i 
\nabla_k N^k = S^i$. *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.

Parameters
lambda[input] $\lambda$.
met_f[input] the flat metric for the Helmholtz decomposition.
method[input] method used to solve the equation:
  • 0 : It uses the Helmholtz decomposition (see documentation of p_potential ), with the flat metric met_f given in argument (the default).
  • 1 : It solves, first for the divergence (calculated using met_f ), then the r -component, the $\eta$ potential, and fianlly the $\mu$ potential (see documentation of Vector_div_free .
  • 2 : The sources is transformed to cartesian components and the equation is solved using Shibata method (see Granclement et al. JCPH 2001.
  • 6 : Solves for the r -component and $ \eta $ together in a system, and for the $ \mu $ potential (which decouples). The solution is then built from these fields through the method Vector::set_vr_eta_mu(). It is the default method.
Returns
the solution $N^i$.

Definition at line 130 of file vector_poisson.C.

References Lorene::Tensor::cmp, Lorene::Vector::div_free(), Lorene::Tensor::mp, Lorene::Vector::mu(), Lorene::Vector::poisson(), Lorene::Vector::potential(), Lorene::Scalar::stdsdp(), and Lorene::Tensor::triad.

◆ poisson() [4/5]

Vector Lorene::Vector::poisson ( double  lambda,
int  method = 6 
) const
inherited

Solves the vector Poisson equation with *this as a source.

The equation solved is $\Delta N^i +\lambda \nabla^i 
\nabla_k N^k = S^i$. *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.

Parameters
lambda[input] $\lambda$.
method[input] method used to solve the equation (see Vector::poisson(double, Metric_flat, int) for details).
Returns
the solution $N^i$.

Definition at line 518 of file vector_poisson.C.

References Lorene::Tensor::mp, Lorene::Vector::poisson(), and Lorene::Tensor::triad.

◆ poisson() [5/5]

Vector_divfree Lorene::Vector_divfree::poisson ( Param par) const

Computes the solution of a vectorial Poisson equation with *this $= \vec{V}$ as a source:

\[
   \Delta \vec{W} = \vec{V}
\]

Returns
solution $\vec{W}$ of the above equation with the boundary condition $\vec{W}=0$ at spatial infinity.

Definition at line 69 of file vector_df_poisson.C.

References Lorene::Tensor::cmp, Lorene::Scalar::div_r(), met_div, Lorene::Tensor::mp, Lorene::Vector::mu(), Lorene::Scalar::set_etat_zero(), and Lorene::Tensor::triad.

◆ poisson_block()

void Lorene::Vector::poisson_block ( double  lambda,
Vector resu 
) const
inherited

Definition at line 75 of file vector_poisson_block.C.

◆ poisson_boundary()

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
inherited

Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere.

The equation solved is $\Delta N^i +\lambda \nabla^i 
\nabla_k N^k = S^i$. *this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential )

Parameters
lambda[input] $\lambda$.
resu[output] the solution $N^i$.

Definition at line 65 of file vector_poisson_boundary.C.

References Lorene::Tensor::cmp, Lorene::Vector::eta(), Lorene::Mg3d::get_angu(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Vector::mu(), Lorene::pow(), Lorene::Tbl::t, and Lorene::Tensor::triad.

◆ poisson_boundary2()

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
inherited

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::Tensor::cmp, Lorene::Vector::eta(), Lorene::Diff_sx2::get_matrice(), Lorene::Mg3d::get_np(), Lorene::Mg3d::get_nr(), Lorene::Mg3d::get_nt(), Lorene::Mg3d::get_nzone(), Lorene::Mg3d::get_type_r(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Vector::mu(), Lorene::Tbl::t, and Lorene::Tensor::triad.

◆ poisson_dirichlet()

Vector Lorene::Vector::poisson_dirichlet ( double  lambda,
const Valeur limit_vr,
const Valeur limit_vt,
const Valeur limit_vp,
int  num_front 
) const
inherited

Definition at line 680 of file vector_poisson_boundary.C.

◆ poisson_neumann()

Vector Lorene::Vector::poisson_neumann ( double  lambda,
const Valeur limit_vr,
const Valeur limit_vt,
const Valeur limit_vp,
int  num_front 
) const
inherited

Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere.

The equation solved is $\Delta N^i +\lambda \nabla^i 
\nabla_k N^k = S^i$. *this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential )

Parameters
lambda[input] $\lambda$.
resu[output] the solution $N^i$.

Definition at line 691 of file vector_poisson_boundary.C.

References Lorene::Tensor::mp, Lorene::Vector::poisson_robin(), and Lorene::Tensor::triad.

◆ poisson_robin()

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
inherited

Solves the vector Poisson equation with *this as a source with a boundary condition on the excised sphere.

The equation solved is $\Delta N^i +\lambda \nabla^i 
\nabla_k N^k = S^i$. *this must be given with dzpuis = 4. It uses the Helmholtz decomposition (see documentation of p_potential )

Parameters
lambda[input] $\lambda$.
resu[output] the solution $N^i$.

Definition at line 702 of file vector_poisson_boundary.C.

References Lorene::Tensor::cmp, Lorene::Tensor::mp, Lorene::norme(), Lorene::Vector::poisson_boundary(), and Lorene::Tensor::triad.

◆ position()

virtual int Lorene::Vector::position ( const Itbl idx) const
inlinevirtualinherited

Returns the position in the Scalar array cmp of a component given by its index.


Returns
position in the Scalar array cmp
corresponding to the index given in 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().

◆ potential()

const Scalar & Lorene::Vector::potential ( const Metric metre) const
inherited

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 $\phi$. Only in the case of contravariant vectors.

Definition at line 492 of file vector.C.

References Lorene::Vector::decompose_div(), Lorene::Tensor::get_place_met(), Lorene::Vector::p_potential, and Lorene::Tensor::set_dependance().

◆ pseudo_spectral_base()

void Lorene::Vector::pseudo_spectral_base ( )
virtualinherited

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, Lorene::Scalar::set_spectral_base(), and Lorene::Tensor::triad.

◆ set()

Scalar & Lorene::Vector::set ( int  index)
inherited

Read/write access to a component.

Definition at line 296 of file vector.C.

References Lorene::Tensor::cmp, and Lorene::Vector::del_deriv().

◆ set_A_mu()

void Lorene::Vector_divfree::set_A_mu ( const Scalar A_i,
const Scalar mu_i,
const Param par_bc 
)

Defines the components through potentials $A$ and $\mu$.

(see members p_A and p_mu ),

Parameters
A_i[input] Angular potential $A$
mu_i[input] Angular potential $\mu$

Definition at line 123 of file vector_divfree_aux.C.

References Lorene::Tensor::cmp, del_deriv(), Lorene::Tensor::mp, Lorene::Vector::p_A, Lorene::Vector::p_eta, Lorene::Vector::p_mu, Lorene::Scalar::set_dzpuis(), sol_Dirac_A(), Lorene::Tensor::triad, and Lorene::Vector::update_vtvp().

◆ set_der_0x0()

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

Sets the pointers on derived quantities to 0x0.

Definition at line 135 of file vector_divfree.C.

◆ set_der_met_0x0()

void Lorene::Vector::set_der_met_0x0 ( int  i) const
protectedinherited

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 Lorene::Vector::p_div_free, and Lorene::Vector::p_potential.

◆ set_vr_eta_mu()

void Lorene::Vector_divfree::set_vr_eta_mu ( const Scalar vr_i,
const Scalar eta_i,
const Scalar mu_i 
)

Defines the components through $V_r$, $\eta$ and $\mu$.

(see members p_eta and p_mu ),

Parameters
vr_i[input] r-component of the vector
eta_i[input] Angular potential $\eta$
mu_i[input] Angular potential $\mu$

Definition at line 99 of file vector_divfree_aux.C.

References Lorene::Tensor::cmp, Lorene::Vector::p_eta, Lorene::Vector::p_mu, Lorene::Tensor::triad, and Lorene::Vector::update_vtvp().

◆ set_vr_mu()

void Lorene::Vector_divfree::set_vr_mu ( const Scalar vr_i,
const Scalar mu_i 
)

Sets the angular potentials $\mu$ (see member p_mu ), and the $V^r$ component of the vector.

The potential $\eta$ is then deduced from $V^r$ by the divergence-free condition. The components $V^\theta$ and $V^\varphi$ are updated consistently by a call to the method update_vtvp() .

Parameters
vr_i[input] component $V^r$ of the vector
mu_i[input] angular potential $\mu$

Definition at line 174 of file vector_divfree.C.

References Lorene::Tensor::cmp, del_deriv(), eta(), Lorene::Vector::p_mu, Lorene::Tensor::triad, and Lorene::Vector::update_vtvp().

◆ sol_Dirac_A()

void Lorene::Vector_divfree::sol_Dirac_A ( const Scalar aaa,
Scalar eta,
Scalar vr,
const Param par_bc = 0x0 
) const
protected

Solves a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.

The system reads :

\[ 
    {\partial \eta \over \partial r} + {\eta \over r} - {V^r \over r} = A
\]

\[ 
    {\partial V^r \over \partial r} + {2 V^r \over r} 
        + {1 \over r}\Delta_{\vartheta \varphi}\eta = 0
\]

Definition at line 77 of file vector_divfree_A.C.

References Lorene::Map_af::get_alpha(), Lorene::Map_af::get_beta(), Lorene::Map::get_mg(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Base_val::mult_x(), and R_CHEBP.

◆ sol_Dirac_A_1z()

void Lorene::Vector_divfree::sol_Dirac_A_1z ( const Scalar aaa,
Scalar eta,
Scalar vr,
const Param par_bc = 0x0 
) const
protected

Solves a one-domain system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.

The system reads :

\[ 
    {\partial \eta \over \partial r} + {\eta \over r} - {V^r \over r} = A
\]

\[ 
    {\partial V^r \over \partial r} + {2 V^r \over r} 
        + {1 \over r}\Delta_{\vartheta \varphi}\eta = 0
\]

Definition at line 69 of file vector_divfree_A_1z.C.

References Lorene::Map_af::get_alpha(), Lorene::Map::get_mg(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Base_val::mult_x(), and R_CHEBP.

◆ sol_Dirac_A_poisson()

void Lorene::Vector_divfree::sol_Dirac_A_poisson ( const Scalar aaa,
Scalar eta,
Scalar vr,
const Param par_bc = 0x0 
) const
protected

Solves via a poisson method a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.

The system reads :

\[ 
    {\partial \eta \over \partial r} + {\eta \over r} - {V^r \over r} = A
\]

\[ 
    {\partial V^r \over \partial r} + {2 V^r \over r} 
        + {1 \over r}\Delta_{\vartheta \varphi}\eta = 0
\]

Definition at line 68 of file vector_divfree_A_poisson.C.

◆ sol_Dirac_A_tau()

void Lorene::Vector_divfree::sol_Dirac_A_tau ( const Scalar aaa,
Scalar eta,
Scalar vr,
const Param par_bc = 0x0 
) const
protected

Solves via a tau method a system of two-coupled first-order PDEs obtained from the divergence-free condition and the requirement that the potential A has a given value.

The system reads :

\[ 
    {\partial \eta \over \partial r} + {\eta \over r} - {V^r \over r} = A
\]

\[ 
    {\partial V^r \over \partial r} + {2 V^r \over r} 
        + {1 \over r}\Delta_{\vartheta \varphi}\eta = 0
\]

Definition at line 73 of file vector_divfree_A_tau.C.

References Lorene::Map::get_mg(), Lorene::Base_val::give_quant_numbers(), Lorene::Tensor::mp, Lorene::Base_val::mult_x(), R_CHEB, and R_JACO02.

◆ std_spectral_base()

void Lorene::Vector::std_spectral_base ( )
virtualinherited

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, Lorene::Scalar::set_spectral_base(), and Lorene::Tensor::triad.

◆ update_etavr()

void Lorene::Vector_divfree::update_etavr ( )

Computes the components $V^r$ and $\eta$ from the potential A and the divergence-free condition, according to :

\[ 
     {\partial \eta \over \ partial r} + { \eta \over r} 
                    - {V^r \over r} = A
\]

\[
    {\partial V^r \over \partial r} + {2 V^r \over r} 
        + {1 \over r} \Delta_{\vartheta \varphi} \eta = 0
\]

Definition at line 76 of file vector_divfree_aux.C.

References Lorene::Tensor::cmp, Lorene::Vector::del_deriv(), Lorene::Tensor::mp, Lorene::Vector::p_A, Lorene::Vector::p_eta, and sol_Dirac_A().

◆ update_vtvp()

void Lorene::Vector::update_vtvp ( )
inherited

Computes the components $V^\theta$ and $V^\varphi$ from the potential $\eta$ and $\mu$, according to:

\[
 V^\theta =  {\partial \eta \over \partial\theta} - 
    {1\over\sin\theta} {\partial \mu \over \partial\varphi}
\]

\[
 V^\varphi =  {1\over\sin\theta} 
            {\partial \eta \over \partial\varphi}
            + {\partial \mu \over \partial\theta} 
\]

Definition at line 167 of file vector_etamu.C.

References Lorene::Tensor::cmp, Lorene::Vector::del_deriv(), Lorene::Scalar::dsdt(), Lorene::Vector::p_eta, Lorene::Vector::p_mu, and Lorene::Scalar::stdsdp().

◆ visu_arrows()

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
inherited

3D visualization via OpenDX.

Parameters
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::Tensor::mp, Lorene::Vector::operator()(), and Lorene::Tensor::triad.

◆ visu_streamline()

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
inherited

Definition at line 283 of file vector_visu.C.

Member Data Documentation

◆ met_div

const Metric* const Lorene::Vector_divfree::met_div
protected

Metric with respect to which the divergence is defined.

Definition at line 730 of file vector.h.

◆ p_A

Scalar* Lorene::Vector::p_A
mutableprotectedinherited

Field $A$ defined by.

\[ 
    A = {\partial \eta \over \ partial r} + { \eta \over r} - {V^r \over r}
 \]

Insensitive to the longitudinal part of the vector, related to the curl.

Definition at line 241 of file vector.h.

◆ p_div_free

Vector_divfree* Lorene::Vector::p_div_free[N_MET_MAX]
mutableprotectedinherited

The divergence-free vector $\vec{W} =  \vec{\nabla} \wedge 
  \vec{\psi}$ of the Helmholtz decomposition of any 3D vector $\vec{V}: \quad \vec{V} = \vec{\nabla} \phi + \vec{\nabla} 
*\wedge \vec{\psi}$.

Only in the case of contravariant vectors.

Definition at line 205 of file vector.h.

◆ p_eta

Scalar* Lorene::Vector::p_eta
mutableprotectedinherited

Field $\eta$ such that the angular components $(V^\theta, V^\varphi)$ of the vector are written:

\[
 V^\theta =   {\partial \eta \over \partial\theta} -
     {1\over\sin\theta} {\partial \mu \over \partial\varphi} 
\]

\[
 V^\varphi =  {1\over\sin\theta} 
            {\partial \eta \over \partial\varphi}
            + {\partial \mu \over \partial\theta} 
\]

Definition at line 219 of file vector.h.

◆ p_mu

Scalar* Lorene::Vector::p_mu
mutableprotectedinherited

Field $\mu$ such that the angular components $(V^\theta, V^\varphi)$ of the vector are written:

\[
 V^\theta =  {\partial \eta \over \partial\theta} -
  {1\over\sin\theta} {\partial \mu \over \partial\varphi} 
\]

\[
 V^\varphi =  {1\over\sin\theta} 
            {\partial \eta \over \partial\varphi}
            + {\partial \mu \over \partial\theta} 
\]

Definition at line 233 of file vector.h.

◆ p_potential

Scalar* Lorene::Vector::p_potential[N_MET_MAX]
mutableprotectedinherited

The potential $\phi$ giving the gradient part in the Helmholtz decomposition of any 3D vector $\vec{V}: \quad \vec{V} = 
\vec{\nabla} \phi + \vec{\nabla} \wedge \vec{\psi}$.

Only in the case of contravariant vectors.

Definition at line 198 of file vector.h.


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