HArD::Core3D
Hybrid Arbitrary Degree::Core 3D - Library to implement 3D schemes with vertex, edge, face and cell polynomials as unknowns
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Public Attributes | List of all members
HArDCore3D::Einstein Struct Reference

Class to assemble and solve an Einstein problem. More...

#include <ecddr-einstein-dtb.hpp>

Public Types

typedef Eigen::SparseMatrix< doubleSystemMatrixType
 
typedef std::function< double(const Eigen::Vector3d &)> ZeroFormType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> OneFormType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> TwoFormType
 
typedef std::function< double(const Eigen::Vector3d &)> ThreeFormType
 
typedef TwoFormType ForcingTermType
 
typedef std::vector< ZeroFormTypeCollection0Forms
 
typedef std::vector< OneFormTypeCollection1Forms
 
typedef std::vector< TwoFormTypeCollection2Forms
 
typedef std::vector< ThreeFormTypeCollection3Forms
 
typedef std::function< double(const double t, const Eigen::Vector3d &)> TZeroFormType
 
typedef std::function< Eigen::Vector3d(const double t, const Eigen::Vector3d &)> TOneFormType
 
typedef std::function< Eigen::Vector3d(const double t, const Eigen::Vector3d &)> TTwoFormType
 
typedef std::function< double(const double t, const Eigen::Vector3d &)> TThreeFormType
 
typedef TTwoFormType TForcingTermType
 
typedef std::vector< ZeroFormTypeTCollection0Forms
 
typedef std::vector< TOneFormTypeTCollection1Forms
 
typedef std::vector< TTwoFormTypeTCollection2Forms
 
typedef std::vector< TThreeFormTypeTCollection3Forms
 
typedef std::function< Eigen::Vector3d(const Eigen::MatrixXd &E, const Eigen::MatrixXd &B)> CompositeHType
 
typedef std::vector< CompositeHTypeCollectionHforms
 
typedef std::function< Eigen::Vector3d(const Eigen::MatrixXd &D)> CompositeEType
 
typedef std::vector< CompositeETypeCollectionEforms
 
typedef std::function< Eigen::Vector3d(const Eigen::VectorXd &E0, const Eigen::VectorXd &D0, const Eigen::MatrixXd &H, const Eigen::MatrixXd &E, const Eigen::MatrixXd &D, const Eigen::MatrixXd &B)> CompositeUType
 
typedef std::vector< CompositeUTypeCollectionUforms
 
typedef Eigen::SparseMatrix< doubleSystemMatrixType
 
typedef std::function< double(const Eigen::Vector3d &)> ZeroFormType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> OneFormType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> TwoFormType
 
typedef std::function< double(const Eigen::Vector3d &)> ThreeFormType
 
typedef TwoFormType ForcingTermType
 
typedef std::vector< ZeroFormTypeCollection0Forms
 
typedef std::vector< OneFormTypeCollection1Forms
 
typedef std::vector< TwoFormTypeCollection2Forms
 
typedef std::vector< ThreeFormTypeCollection3Forms
 
typedef std::function< double(const double t, const Eigen::Vector3d &)> TZeroFormType
 
typedef std::function< Eigen::Vector3d(const double t, const Eigen::Vector3d &)> TOneFormType
 
typedef std::function< Eigen::Vector3d(const double t, const Eigen::Vector3d &)> TTwoFormType
 
typedef std::function< double(const double t, const Eigen::Vector3d &)> TThreeFormType
 
typedef TTwoFormType TForcingTermType
 
typedef std::vector< ZeroFormTypeTCollection0Forms
 
typedef std::vector< TOneFormTypeTCollection1Forms
 
typedef std::vector< TTwoFormTypeTCollection2Forms
 
typedef std::vector< TThreeFormTypeTCollection3Forms
 
typedef std::function< Eigen::Vector3d(const Eigen::MatrixXd &E, const Eigen::MatrixXd &B)> CompositeHType
 
typedef std::vector< CompositeHTypeCollectionHforms
 
typedef std::function< Eigen::Vector3d(const Eigen::MatrixXd &D)> CompositeEType
 
typedef std::vector< CompositeETypeCollectionEforms
 
typedef std::function< Eigen::Vector3d(const Eigen::VectorXd &E0, const Eigen::VectorXd &D0, const Eigen::MatrixXd &H, const Eigen::MatrixXd &E, const Eigen::MatrixXd &D, const Eigen::MatrixXd &B)> CompositeUType
 
typedef std::vector< CompositeUTypeCollectionUforms
 

Public Member Functions

 Einstein (const DDR_Spaces &ddrspaces, bool use_threads, double stab_par, std::ostream &output=std::cout)
 Constructor.
 
Eigen::VectorXd updateLapse (double dt, Eigen::VectorXd &lapse_n, Eigen::VectorXd &starDtB)
 
void precomputeL2 ()
 
void assembleLinearSystem (double dt, const ZeroFormType &lapse, const OneFormType &E0, const Eigen::VectorXd &starDtB)
 
void computeBCRHS (double dt, const ZeroFormType &lapse_n, const std::vector< TwoFormType > &D, const std::vector< OneFormType > &theta, const std::vector< TwoFormType > &B, const std::vector< TwoFormType > &dtD, const std::vector< OneFormType > &E, const std::vector< OneFormType > &H, const std::vector< TwoFormType > &dH)
 
void updateThetaRHS (const Eigen::VectorXd &starD)
 
double computeResidual (const Eigen::VectorXd &starDB) const
 Compute the residual.
 
Eigen::VectorXd computeConstraints (const ZeroFormType &lapse_n, const std::vector< OneFormType > &theta0, const std::vector< OneFormType > &thetan, const Eigen::VectorXd &starDB0, const Eigen::VectorXd &starDBn, LinearSolver< Einstein::SystemMatrixType > &solver) const
 Compute the constraints.
 
std::vector< EinsteinNormscomputeEinsteinNorms (const std::vector< Eigen::VectorXd > &list_dofs) const
 Compute the discrete L2 norms, for a family of Eigen::VectorXd.
 
std::pair< EinsteinNorms, EinsteinNormscomputeContinuousErrorsNorms (const Eigen::VectorXd &starDt, const std::vector< OneFormType > &starD, const std::vector< TwoFormType > &starTheta, const std::vector< TwoFormType > &starB) const
 
Eigen::MatrixXd basisChange (size_t k, const Eigen::MatrixXd &forms, const Eigen::MatrixXd &basis) const
 Change of basis for 2 forms.
 
Eigen::VectorXd contractIndex (size_t index, const Eigen::MatrixXd &forms) const
 Contract collection index with first form index (A_i)^{ji}.
 
Eigen::MatrixXd contractIndex (size_t index, const Eigen::MatrixXd &twoForm1, const Eigen::MatrixXd &twoForm2) const
 Contract collection index with first form index (A_i)^{ji}.
 
Eigen::Vector3d wedge11 (const Eigen::Vector3d &form1, const Eigen::Vector3d &form2) const
 wedge a 1-form and a 1-form to return 2-form as [A12, A13, A23]
 
Eigen::Matrix3d expand2Form (const Eigen::Vector3d &form) const
 expand a 2 form vector into a 3x3 matrix
 
size_t dimensionSpace () const
 Returns the global problem dimension.
 
const DDR_SpacesddrSpaces () const
 Returns the ECDDR sequence.
 
const SystemMatrixTypesystemMatrix () const
 Returns the linear system matrix.
 
SystemMatrixTypesystemMatrix ()
 Returns the linear system matrix.
 
const Eigen::VectorXd & systemVector () const
 Returns the linear system right-hand side vector.
 
Eigen::VectorXd & systemVector ()
 Returns the linear system right-hand side vector.
 
const SystemMatrixTypeL2_X1 () const
 Returns the L of the LLT decomposition of the L2 product of X1.
 
SystemMatrixTypeL2_X1 ()
 Returns the L of the LLT decomposition of the L2 product of X1.
 
const SystemMatrixTypeL2_X2 () const
 Returns the L of the LLT decomposition of the L2 product of X2.
 
SystemMatrixTypeL2_X2 ()
 Returns the L of the LLT decomposition of the L2 product of X2.
 
const SystemMatrixTypelhsTheta () const
 Returns the lhs of theta equation.
 
SystemMatrixTypelhsTheta ()
 Returns the lhs of theta equation.
 
const Eigen::VectorXd & rhsTheta () const
 Returns the rhs of theta equation.
 
Eigen::VectorXd & rhsTheta ()
 Returns the rhs of theta equation.
 
const doublestabilizationParameter () const
 Returns the stabilization parameter.
 
doublestabilizationParameter ()
 Returns the stabilization parameter.
 
doublemin_detTheta ()
 
template<typename outValue , typename TFct >
std::function< outValue(const Eigen::Vector3d &)> contractPara (const TFct &F, double t) const
 Takes a two parameter function and a value and returns a function with the first parameter fixed to value.
 
template<typename outValue , typename Fct , typename TFct >
std::vector< FctcontractParaVec (const std::vector< TFct > &TF, double t) const
 Takes a vector of two parameter functions and a value and returns a function with the first parameter fixed to value.
 
 Einstein (const DDR_Spaces &ddrspaces, bool use_threads, double stab_par, std::ostream &output=std::cout)
 Constructor.
 
Eigen::VectorXd updateLapse (double dt, Eigen::VectorXd &lapse_n, Eigen::VectorXd &starDtB)
 
void assembleLinearSystem (double dt, const ZeroFormType &lapse, const OneFormType &E0, const Eigen::VectorXd &starDt)
 
void computeBCRHS (double dt, const ZeroFormType &lapse_n, const std::vector< TwoFormType > &D, const std::vector< OneFormType > &theta, const std::vector< TwoFormType > &B, const std::vector< TwoFormType > &dtD, const std::vector< OneFormType > &E, const std::vector< OneFormType > &H, const std::vector< TwoFormType > &dH)
 
void updateThetaRHS (const Eigen::VectorXd &starD)
 
std::vector< EinsteinNormscomputeEinsteinNorms (const std::vector< Eigen::VectorXd > &list_dofs) const
 Compute the discrete L2 norms, for a family of Eigen::VectorXd.
 
std::pair< EinsteinNorms, EinsteinNormscomputeContinuousErrorsNorms (const Eigen::VectorXd &starDt, const std::vector< OneFormType > &starD, const std::vector< TwoFormType > &theta) const
 
Eigen::MatrixXd basisChange (size_t k, const Eigen::MatrixXd &forms, const Eigen::MatrixXd &basis) const
 Change of basis for 2 forms.
 
Eigen::MatrixXd contractIndex (size_t index, const Eigen::MatrixXd &twoForm1, const Eigen::MatrixXd &twoForm2) const
 Contract collection index with first form index (A_i)^{ji}.
 
Eigen::Vector3d wedge11 (const Eigen::Vector3d &form1, const Eigen::Vector3d &form2) const
 wedge a 1-form and a 1-form to return 2-form as [A12, A13, A23]
 
Eigen::Matrix3d expand2Form (const Eigen::Vector3d &form) const
 expand a 2 form vector into a 3x3 matrix
 
size_t dimensionSpace () const
 Returns the global problem dimension.
 
const DDR_SpacesddrSpaces () const
 Returns the ECDDR sequence.
 
const SystemMatrixTypesystemMatrix () const
 Returns the linear system matrix.
 
SystemMatrixTypesystemMatrix ()
 Returns the linear system matrix.
 
const Eigen::VectorXd & systemVector () const
 Returns the linear system right-hand side vector.
 
Eigen::VectorXd & systemVector ()
 Returns the linear system right-hand side vector.
 
const doublestabilizationParameter () const
 Returns the stabilization parameter.
 
doublestabilizationParameter ()
 Returns the stabilization parameter.
 
doublemin_detTheta ()
 
template<typename outValue , typename TFct >
std::function< outValue(const Eigen::Vector3d &)> contractPara (const TFct &F, double t) const
 Takes a two parameter function and a value and returns a function with the first parameter fixed to value.
 
template<typename outValue , typename Fct , typename TFct >
std::vector< FctcontractParaVec (const std::vector< TFct > &TF, double t) const
 Takes a vector of two parameter functions and a value and returns a function with the first parameter fixed to value.
 

Public Attributes

Einstein::CompositeHType composite_orth_H0
 
Einstein::CompositeHType composite_orth_H1
 
Einstein::CompositeHType composite_orth_H2
 
Einstein::CompositeHType composite_orth_H3
 
Einstein::CollectionHforms composite_orth_H = {composite_orth_H0, composite_orth_H1, composite_orth_H2, composite_orth_H3}
 
Einstein::CompositeEType composite_orth_E0
 
Einstein::CompositeEType composite_orth_E1
 
Einstein::CompositeEType composite_orth_E2
 
Einstein::CompositeEType composite_orth_E3
 
Einstein::CollectionEforms composite_orth_E = {composite_orth_E0, composite_orth_E1, composite_orth_E2, composite_orth_E3}
 
Einstein::CompositeUType composite_orth_U0
 
Einstein::CompositeUType composite_orth_U1
 
Einstein::CompositeUType composite_orth_U2
 
Einstein::CompositeUType composite_orth_U3
 
Einstein::CollectionUforms composite_orth_U = {composite_orth_U0, composite_orth_U1, composite_orth_U2, composite_orth_U3}
 

Detailed Description

Class to assemble and solve an Einstein problem.

Member Typedef Documentation

◆ Collection0Forms [1/2]

◆ Collection0Forms [2/2]

◆ Collection1Forms [1/2]

◆ Collection1Forms [2/2]

◆ Collection2Forms [1/2]

◆ Collection2Forms [2/2]

◆ Collection3Forms [1/2]

◆ Collection3Forms [2/2]

◆ CollectionEforms [1/2]

◆ CollectionEforms [2/2]

◆ CollectionHforms [1/2]

◆ CollectionHforms [2/2]

◆ CollectionUforms [1/2]

◆ CollectionUforms [2/2]

◆ CompositeEType [1/2]

typedef std::function<Eigen::Vector3d(const Eigen::MatrixXd & D)> HArDCore3D::Einstein::CompositeEType

◆ CompositeEType [2/2]

typedef std::function<Eigen::Vector3d(const Eigen::MatrixXd & D)> HArDCore3D::Einstein::CompositeEType

◆ CompositeHType [1/2]

typedef std::function<Eigen::Vector3d(const Eigen::MatrixXd & E, const Eigen::MatrixXd & B)> HArDCore3D::Einstein::CompositeHType

◆ CompositeHType [2/2]

typedef std::function<Eigen::Vector3d(const Eigen::MatrixXd & E, const Eigen::MatrixXd & B)> HArDCore3D::Einstein::CompositeHType

◆ CompositeUType [1/2]

typedef std::function<Eigen::Vector3d(const Eigen::VectorXd & E0, const Eigen::VectorXd & D0, const Eigen::MatrixXd & H, const Eigen::MatrixXd & E, const Eigen::MatrixXd & D, const Eigen::MatrixXd & B)> HArDCore3D::Einstein::CompositeUType

◆ CompositeUType [2/2]

typedef std::function<Eigen::Vector3d(const Eigen::VectorXd & E0, const Eigen::VectorXd & D0, const Eigen::MatrixXd & H, const Eigen::MatrixXd & E, const Eigen::MatrixXd & D, const Eigen::MatrixXd & B)> HArDCore3D::Einstein::CompositeUType

◆ ForcingTermType [1/2]

◆ ForcingTermType [2/2]

◆ OneFormType [1/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Einstein::OneFormType

◆ OneFormType [2/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Einstein::OneFormType

◆ SystemMatrixType [1/2]

◆ SystemMatrixType [2/2]

◆ TCollection0Forms [1/2]

◆ TCollection0Forms [2/2]

◆ TCollection1Forms [1/2]

◆ TCollection1Forms [2/2]

◆ TCollection2Forms [1/2]

◆ TCollection2Forms [2/2]

◆ TCollection3Forms [1/2]

◆ TCollection3Forms [2/2]

◆ TForcingTermType [1/2]

◆ TForcingTermType [2/2]

◆ ThreeFormType [1/2]

typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::Einstein::ThreeFormType

◆ ThreeFormType [2/2]

typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::Einstein::ThreeFormType

◆ TOneFormType [1/2]

typedef std::function<Eigen::Vector3d(const double t, const Eigen::Vector3d &)> HArDCore3D::Einstein::TOneFormType

◆ TOneFormType [2/2]

typedef std::function<Eigen::Vector3d(const double t, const Eigen::Vector3d &)> HArDCore3D::Einstein::TOneFormType

◆ TThreeFormType [1/2]

◆ TThreeFormType [2/2]

◆ TTwoFormType [1/2]

typedef std::function<Eigen::Vector3d(const double t, const Eigen::Vector3d &)> HArDCore3D::Einstein::TTwoFormType

◆ TTwoFormType [2/2]

typedef std::function<Eigen::Vector3d(const double t, const Eigen::Vector3d &)> HArDCore3D::Einstein::TTwoFormType

◆ TwoFormType [1/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Einstein::TwoFormType

◆ TwoFormType [2/2]

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::Einstein::TwoFormType

◆ TZeroFormType [1/2]

◆ TZeroFormType [2/2]

◆ ZeroFormType [1/2]

typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::Einstein::ZeroFormType

◆ ZeroFormType [2/2]

typedef std::function<double(const Eigen::Vector3d &)> HArDCore3D::Einstein::ZeroFormType

Constructor & Destructor Documentation

◆ Einstein() [1/2]

Einstein::Einstein ( const DDR_Spaces ddrspaces,
bool  use_threads,
double  stab_par,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
ddrspacesCore for the ECDDR space sequence
use_threadsTrue for parallel execution, false for sequential execution
outputOutput stream to print status messages

◆ Einstein() [2/2]

HArDCore3D::Einstein::Einstein ( const DDR_Spaces ddrspaces,
bool  use_threads,
double  stab_par,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
ddrspacesCore for the ECDDR space sequence
use_threadsTrue for parallel execution, false for sequential execution
outputOutput stream to print status messages

Member Function Documentation

◆ assembleLinearSystem() [1/2]

void HArDCore3D::Einstein::assembleLinearSystem ( double  dt,
const ZeroFormType lapse,
const OneFormType E0,
const Eigen::VectorXd &  starDt 
)

◆ assembleLinearSystem() [2/2]

void Einstein::assembleLinearSystem ( double  dt,
const ZeroFormType lapse,
const OneFormType E0,
const Eigen::VectorXd &  starDtB 
)

◆ basisChange() [1/2]

Eigen::MatrixXd Einstein::basisChange ( size_t  k,
const Eigen::MatrixXd &  forms,
const Eigen::MatrixXd &  basis 
) const

Change of basis for 2 forms.

Parameters
kform degree 1 or 2
formsForms either 1-form [A_1, A_2, A_3] or 2-form [B_12, B_13, B_23] in column vectors
basisBasis vectors by column

◆ basisChange() [2/2]

Eigen::MatrixXd HArDCore3D::Einstein::basisChange ( size_t  k,
const Eigen::MatrixXd &  forms,
const Eigen::MatrixXd &  basis 
) const

Change of basis for 2 forms.

Parameters
kform degree 1 or 2
formsForms either 1-form [A_1, A_2, A_3] or 2-form [B_12, B_13, B_23] in column vectors
basisBasis vectors by column

◆ computeBCRHS() [1/2]

void HArDCore3D::Einstein::computeBCRHS ( double  dt,
const ZeroFormType lapse_n,
const std::vector< TwoFormType > &  D,
const std::vector< OneFormType > &  theta,
const std::vector< TwoFormType > &  B,
const std::vector< TwoFormType > &  dtD,
const std::vector< OneFormType > &  E,
const std::vector< OneFormType > &  H,
const std::vector< TwoFormType > &  dH 
)

◆ computeBCRHS() [2/2]

void HArDCore3D::Einstein::computeBCRHS ( double  dt,
const ZeroFormType lapse_n,
const std::vector< TwoFormType > &  D,
const std::vector< OneFormType > &  theta,
const std::vector< TwoFormType > &  B,
const std::vector< TwoFormType > &  dtD,
const std::vector< OneFormType > &  E,
const std::vector< OneFormType > &  H,
const std::vector< TwoFormType > &  dH 
)

◆ computeConstraints()

Eigen::VectorXd Einstein::computeConstraints ( const ZeroFormType lapse_n,
const std::vector< OneFormType > &  theta0,
const std::vector< OneFormType > &  thetan,
const Eigen::VectorXd &  starDB0,
const Eigen::VectorXd &  starDBn,
LinearSolver< Einstein::SystemMatrixType > &  solver 
) const

Compute the constraints.

◆ computeContinuousErrorsNorms() [1/2]

std::pair< EinsteinNorms, EinsteinNorms > Einstein::computeContinuousErrorsNorms ( const Eigen::VectorXd &  starDt,
const std::vector< OneFormType > &  starD,
const std::vector< TwoFormType > &  starTheta,
const std::vector< TwoFormType > &  starB 
) const
Parameters
starDtList of vectors representing velocities and pressures

◆ computeContinuousErrorsNorms() [2/2]

std::pair< EinsteinNorms, EinsteinNorms > Einstein::computeContinuousErrorsNorms ( const Eigen::VectorXd &  starDt,
const std::vector< OneFormType > &  starD,
const std::vector< TwoFormType > &  theta 
) const
Parameters
starDtList of vectors representing velocities and pressures

◆ computeEinsteinNorms() [1/2]

std::vector< EinsteinNorms > Einstein::computeEinsteinNorms ( const std::vector< Eigen::VectorXd > &  list_dofs) const

Compute the discrete L2 norms, for a family of Eigen::VectorXd.

Parameters
list_dofsList of vectors representing velocities and pressures

◆ computeEinsteinNorms() [2/2]

std::vector< EinsteinNorms > HArDCore3D::Einstein::computeEinsteinNorms ( const std::vector< Eigen::VectorXd > &  list_dofs) const

Compute the discrete L2 norms, for a family of Eigen::VectorXd.

Parameters
list_dofsList of vectors representing velocities and pressures

◆ computeResidual()

double Einstein::computeResidual ( const Eigen::VectorXd &  starDB) const

Compute the residual.

◆ contractIndex() [1/3]

Eigen::VectorXd Einstein::contractIndex ( size_t  index,
const Eigen::MatrixXd &  forms 
) const

Contract collection index with first form index (A_i)^{ji}.

Contracts a vector of 2-forms with the first or second index: (A_i)ij (indices=12); (A_i)ji (indices=13)

Parameters
indexfirst or both
formsForms either 1-form [A_1, A_2, A_3] or 2-form [B_12, B_13, B_23] in column vectors

◆ contractIndex() [2/3]

Eigen::MatrixXd Einstein::contractIndex ( size_t  index,
const Eigen::MatrixXd &  twoForm1,
const Eigen::MatrixXd &  twoForm2 
) const

Contract collection index with first form index (A_i)^{ji}.

Parameters
indexfirst or both
twoForm1Forms either 1-form [A_1, A_2, A_3] or 2-form [B_12, B_13, B_23] in column vectors
twoForm2Forms either 1-form [A_1, A_2, A_3] or 2-form [B_12, B_13, B_23] in column vectors

◆ contractIndex() [3/3]

Eigen::MatrixXd HArDCore3D::Einstein::contractIndex ( size_t  index,
const Eigen::MatrixXd &  twoForm1,
const Eigen::MatrixXd &  twoForm2 
) const

Contract collection index with first form index (A_i)^{ji}.

Parameters
indexfirst or both
twoForm1Forms either 1-form [A_1, A_2, A_3] or 2-form [B_12, B_13, B_23] in column vectors
twoForm2Forms either 1-form [A_1, A_2, A_3] or 2-form [B_12, B_13, B_23] in column vectors

◆ contractPara() [1/2]

std::function< outValue(const Eigen::Vector3d &)> HArDCore3D::Einstein::contractPara ( const TFct F,
double  t 
) const
inline

Takes a two parameter function and a value and returns a function with the first parameter fixed to value.

◆ contractPara() [2/2]

std::function< outValue(const Eigen::Vector3d &)> HArDCore3D::Einstein::contractPara ( const TFct F,
double  t 
) const
inline

Takes a two parameter function and a value and returns a function with the first parameter fixed to value.

◆ contractParaVec() [1/2]

std::vector< Fct > HArDCore3D::Einstein::contractParaVec ( const std::vector< TFct > &  TF,
double  t 
) const
inline

Takes a vector of two parameter functions and a value and returns a function with the first parameter fixed to value.

◆ contractParaVec() [2/2]

std::vector< Fct > HArDCore3D::Einstein::contractParaVec ( const std::vector< TFct > &  TF,
double  t 
) const
inline

Takes a vector of two parameter functions and a value and returns a function with the first parameter fixed to value.

◆ ddrSpaces() [1/2]

const DDR_Spaces & HArDCore3D::Einstein::ddrSpaces ( ) const
inline

Returns the ECDDR sequence.

◆ ddrSpaces() [2/2]

const DDR_Spaces & HArDCore3D::Einstein::ddrSpaces ( ) const
inline

Returns the ECDDR sequence.

◆ dimensionSpace() [1/2]

size_t HArDCore3D::Einstein::dimensionSpace ( ) const
inline

Returns the global problem dimension.

◆ dimensionSpace() [2/2]

size_t HArDCore3D::Einstein::dimensionSpace ( ) const
inline

Returns the global problem dimension.

◆ expand2Form() [1/2]

Eigen::Matrix3d Einstein::expand2Form ( const Eigen::Vector3d &  form) const

expand a 2 form vector into a 3x3 matrix

Parameters
form2 form as [A_12, A_13, A_23]

◆ expand2Form() [2/2]

Eigen::Matrix3d HArDCore3D::Einstein::expand2Form ( const Eigen::Vector3d &  form) const

expand a 2 form vector into a 3x3 matrix

Parameters
form2 form as [A_12, A_13, A_23]

◆ L2_X1() [1/2]

SystemMatrixType & HArDCore3D::Einstein::L2_X1 ( )
inline

Returns the L of the LLT decomposition of the L2 product of X1.

◆ L2_X1() [2/2]

const SystemMatrixType & HArDCore3D::Einstein::L2_X1 ( ) const
inline

Returns the L of the LLT decomposition of the L2 product of X1.

◆ L2_X2() [1/2]

SystemMatrixType & HArDCore3D::Einstein::L2_X2 ( )
inline

Returns the L of the LLT decomposition of the L2 product of X2.

◆ L2_X2() [2/2]

const SystemMatrixType & HArDCore3D::Einstein::L2_X2 ( ) const
inline

Returns the L of the LLT decomposition of the L2 product of X2.

◆ lhsTheta() [1/2]

SystemMatrixType & HArDCore3D::Einstein::lhsTheta ( )
inline

Returns the lhs of theta equation.

◆ lhsTheta() [2/2]

const SystemMatrixType & HArDCore3D::Einstein::lhsTheta ( ) const
inline

Returns the lhs of theta equation.

◆ min_detTheta() [1/2]

double & HArDCore3D::Einstein::min_detTheta ( )
inline

◆ min_detTheta() [2/2]

double & HArDCore3D::Einstein::min_detTheta ( )
inline

◆ precomputeL2()

void Einstein::precomputeL2 ( )

◆ rhsTheta() [1/2]

Eigen::VectorXd & HArDCore3D::Einstein::rhsTheta ( )
inline

Returns the rhs of theta equation.

◆ rhsTheta() [2/2]

const Eigen::VectorXd & HArDCore3D::Einstein::rhsTheta ( ) const
inline

Returns the rhs of theta equation.

◆ stabilizationParameter() [1/4]

double & HArDCore3D::Einstein::stabilizationParameter ( )
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [2/4]

double & HArDCore3D::Einstein::stabilizationParameter ( )
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [3/4]

const double & HArDCore3D::Einstein::stabilizationParameter ( ) const
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [4/4]

const double & HArDCore3D::Einstein::stabilizationParameter ( ) const
inline

Returns the stabilization parameter.

◆ systemMatrix() [1/4]

SystemMatrixType & HArDCore3D::Einstein::systemMatrix ( )
inline

Returns the linear system matrix.

◆ systemMatrix() [2/4]

SystemMatrixType & HArDCore3D::Einstein::systemMatrix ( )
inline

Returns the linear system matrix.

◆ systemMatrix() [3/4]

const SystemMatrixType & HArDCore3D::Einstein::systemMatrix ( ) const
inline

Returns the linear system matrix.

◆ systemMatrix() [4/4]

const SystemMatrixType & HArDCore3D::Einstein::systemMatrix ( ) const
inline

Returns the linear system matrix.

◆ systemVector() [1/4]

Eigen::VectorXd & HArDCore3D::Einstein::systemVector ( )
inline

Returns the linear system right-hand side vector.

◆ systemVector() [2/4]

Eigen::VectorXd & HArDCore3D::Einstein::systemVector ( )
inline

Returns the linear system right-hand side vector.

◆ systemVector() [3/4]

const Eigen::VectorXd & HArDCore3D::Einstein::systemVector ( ) const
inline

Returns the linear system right-hand side vector.

◆ systemVector() [4/4]

const Eigen::VectorXd & HArDCore3D::Einstein::systemVector ( ) const
inline

Returns the linear system right-hand side vector.

◆ updateLapse() [1/2]

Eigen::VectorXd Einstein::updateLapse ( double  dt,
Eigen::VectorXd &  lapse_n,
Eigen::VectorXd &  starDtB 
)

Approximations of D, theta, B in the cell

◆ updateLapse() [2/2]

Eigen::VectorXd HArDCore3D::Einstein::updateLapse ( double  dt,
Eigen::VectorXd &  lapse_n,
Eigen::VectorXd &  starDtB 
)

◆ updateThetaRHS() [1/2]

void Einstein::updateThetaRHS ( const Eigen::VectorXd &  starD)

◆ updateThetaRHS() [2/2]

void HArDCore3D::Einstein::updateThetaRHS ( const Eigen::VectorXd &  starD)

◆ wedge11() [1/2]

Eigen::Vector3d Einstein::wedge11 ( const Eigen::Vector3d &  form1,
const Eigen::Vector3d &  form2 
) const

wedge a 1-form and a 1-form to return 2-form as [A12, A13, A23]

◆ wedge11() [2/2]

Eigen::Vector3d HArDCore3D::Einstein::wedge11 ( const Eigen::Vector3d &  form1,
const Eigen::Vector3d &  form2 
) const

wedge a 1-form and a 1-form to return 2-form as [A12, A13, A23]

Member Data Documentation

◆ composite_orth_E

◆ composite_orth_E0

Einstein::CompositeEType HArDCore3D::Einstein::composite_orth_E0
Initial value:
= [](const Eigen::Matrix3d & orth_D) -> Eigen::Vector3d {
return Eigen::Vector3d(
0,
0,
0
);
}
@ Matrix
Definition basis.hpp:67

◆ composite_orth_E1

Einstein::CompositeEType HArDCore3D::Einstein::composite_orth_E1
Initial value:
= [](const Eigen::Matrix3d & orth_D) -> Eigen::Vector3d {
Eigen::Matrix3d starD = Eigen::Matrix<double,3,3>{{0.,0.,1.},{0.,-1.,0.},{1.,0.,0.}} * orth_D;
return Eigen::Vector3d(
-0.5*starD(0,0)+0.5*(starD(1,1)+starD(2,2)),
-starD(0,1),
-starD(0,2)
);
}

◆ composite_orth_E2

Einstein::CompositeEType HArDCore3D::Einstein::composite_orth_E2
Initial value:
= [](const Eigen::Matrix3d & orth_D) -> Eigen::Vector3d {
Eigen::Matrix3d starD = Eigen::Matrix<double,3,3>{{0.,0.,1.},{0.,-1.,0.},{1.,0.,0.}} * orth_D;
return Eigen::Vector3d(
-starD(1,0),
-0.5*starD(1,1)+0.5*starD(0,0)+0.5*starD(2,2),
-starD(1,2)
);
}

◆ composite_orth_E3

Einstein::CompositeEType HArDCore3D::Einstein::composite_orth_E3
Initial value:
= [](const Eigen::Matrix3d & orth_D) -> Eigen::Vector3d {
Eigen::Matrix3d starD = Eigen::Matrix<double,3,3>{{0.,0.,1.},{0.,-1.,0.},{1.,0.,0.}} * orth_D;
return Eigen::Vector3d(
-starD(2,0),
-starD(2,1),
-0.5*starD(2,2)+0.5*starD(0,0)+0.5*starD(1,1)
);
}

◆ composite_orth_H

◆ composite_orth_H0

Einstein::CompositeHType HArDCore3D::Einstein::composite_orth_H0
Initial value:
= [](const Eigen::VectorXd & orth_E0, const Eigen::Matrix3d & orth_B) -> Eigen::Vector3d {
return Eigen::Vector3d(
0,
0,
0
);
}

◆ composite_orth_H1

Einstein::CompositeHType HArDCore3D::Einstein::composite_orth_H1
Initial value:
= [](const Eigen::VectorXd & orth_E0, const Eigen::Matrix3d & orth_B) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*(-orth_B(1,1)+orth_B(0,2))+0.5*orth_B(2,0),
orth_E0(2)+orth_B(2,1),
-orth_E0(1)+orth_B(2,2)
);
}

◆ composite_orth_H2

Einstein::CompositeHType HArDCore3D::Einstein::composite_orth_H2
Initial value:
= [](const Eigen::VectorXd & orth_E0, const Eigen::Matrix3d & orth_B) -> Eigen::Vector3d {
return Eigen::Vector3d(
-orth_E0(2)-orth_B(1,0),
-0.5*(orth_B(0,2)+orth_B(2,0))-0.5*orth_B(1,1),
orth_E0(0)-orth_B(1,2)
);
}

◆ composite_orth_H3

Einstein::CompositeHType HArDCore3D::Einstein::composite_orth_H3
Initial value:
= [](const Eigen::VectorXd & orth_E0, const Eigen::Matrix3d & orth_B) -> Eigen::Vector3d {
return Eigen::Vector3d(
orth_E0(1)+orth_B(0,0),
-orth_E0(0)+orth_B(0,1),
-0.5*(orth_B(2,0)-orth_B(1,1))+0.5*orth_B(0,2)
);
}

◆ composite_orth_U

◆ composite_orth_U0

Einstein::CompositeUType HArDCore3D::Einstein::composite_orth_U0
Initial value:
= [](const Eigen::Vector3d & E0, const Eigen::Vector3d & D0, const Eigen::Matrix3d & H, const Eigen::Matrix3d & E, const Eigen::Matrix3d & D, const Eigen::Matrix3d & B) -> Eigen::Vector3d {
return Eigen::Vector3d(
0,
0,
0
);
}

◆ composite_orth_U1

Einstein::CompositeUType HArDCore3D::Einstein::composite_orth_U1
Initial value:
= [](const Eigen::Vector3d & E0, const Eigen::Vector3d & D0, const Eigen::Matrix3d & H, const Eigen::Matrix3d & E, const Eigen::Matrix3d & D, const Eigen::Matrix3d & B) -> Eigen::Vector3d {
Eigen::Matrix<double,3,3> hstar{{0.,0.,1.},{0.,-1.,0.},{1.,0.,0.}};
Eigen::Vector3d starD0 = hstar * D0;
Eigen::Matrix3d starD = hstar * D;
Eigen::Matrix3d starB = hstar * B;
Eigen::MatrixXd starU = Eigen::MatrixXd::Identity(3,3)
* 0.5 * (E0.transpose()*starD0 + (E.transpose()*starD).trace()
- (H.transpose()*starB).trace())
- starD0*E0.transpose() - starD * E.transpose() + starB * H.transpose();
return (hstar*starU).col(0);
}

◆ composite_orth_U2

Einstein::CompositeUType HArDCore3D::Einstein::composite_orth_U2
Initial value:
= [](const Eigen::Vector3d & E0, const Eigen::Vector3d & D0, const Eigen::Matrix3d & H, const Eigen::Matrix3d & E, const Eigen::Matrix3d & D, const Eigen::Matrix3d & B) -> Eigen::Vector3d {
Eigen::Matrix<double,3,3> hstar{{0.,0.,1.},{0.,-1.,0.},{1.,0.,0.}};
Eigen::Vector3d starD0 = hstar * D0;
Eigen::Matrix3d starD = hstar * D;
Eigen::Matrix3d starB = hstar * B;
Eigen::MatrixXd starU = Eigen::MatrixXd::Identity(3,3)
* 0.5 * (E0.transpose()*starD0 + (E.transpose()*starD).trace()
- (H.transpose()*starB).trace())
- starD0*E0.transpose() - starD * E.transpose() + starB * H.transpose();
return (hstar*starU).col(1);
}

◆ composite_orth_U3

Einstein::CompositeUType HArDCore3D::Einstein::composite_orth_U3
Initial value:
= [](const Eigen::Vector3d & E0, const Eigen::Vector3d & D0, const Eigen::Matrix3d & H, const Eigen::Matrix3d & E, const Eigen::Matrix3d & D, const Eigen::Matrix3d & B) -> Eigen::Vector3d {
Eigen::Matrix<double,3,3> hstar{{0.,0.,1.},{0.,-1.,0.},{1.,0.,0.}};
Eigen::Vector3d starD0 = hstar * D0;
Eigen::Matrix3d starD = hstar * D;
Eigen::Matrix3d starB = hstar * B;
Eigen::MatrixXd starU = Eigen::MatrixXd::Identity(3,3)
* 0.5 * (E0.transpose()*starD0 + (E.transpose()*starD).trace()
- (H.transpose()*starB).trace())
- starD0*E0.transpose() - starD * E.transpose() + starB * H.transpose();
return (hstar*starU).col(2);
}

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