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
Classes | Typedefs | Functions | Variables
DDR_yangmills

Implementation of the lowest order DDR scheme for the Yang-Mills problem. More...

Collaboration diagram for DDR_yangmills:

Classes

struct  HArDCore3D::YangMillsNorms
 Structure to store norm components. More...
 
struct  HArDCore3D::YangMills
 Assemble a YangMills problem. More...
 
struct  HArDCore3D::AssembleSystems
 Struct to store the systems that need assembling. More...
 
struct  HArDCore3D::AssembledSystems
 Struct to store the assembled systems. More...
 

Typedefs

typedef Eigen::SparseMatrix< doubleHArDCore3D::YangMills::SystemMatrixType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::YangMills::ForcingTermType
 Base functions.
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::YangMills::ElectricFieldType
 
typedef std::function< Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::YangMills::MagneticFieldType
 
typedef std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> HArDCore3D::YangMills::TForcingTermType
 Base functions with time component.
 
typedef std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> HArDCore3D::YangMills::TElectricFieldType
 
typedef std::function< Eigen::Vector3d(const double &, const Eigen::Vector3d &)> HArDCore3D::YangMills::TMagneticFieldType
 
typedef std::vector< ForcingTermTypeHArDCore3D::YangMills::LAForcingTermType
 Lie algebra valued functions (vector of dim Lie algebra)
 
typedef std::vector< ElectricFieldTypeHArDCore3D::YangMills::LAElectricFieldType
 
typedef std::vector< MagneticFieldTypeHArDCore3D::YangMills::LAMagneticFieldType
 
typedef std::vector< TForcingTermTypeHArDCore3D::YangMills::TLAForcingTermType
 Lie algebra valued functions with time component.
 
typedef std::vector< TElectricFieldTypeHArDCore3D::YangMills::TLAElectricFieldType
 
typedef std::vector< TMagneticFieldTypeHArDCore3D::YangMills::TLAMagneticFieldType
 

Functions

 HArDCore3D::YangMillsNorms::YangMillsNorms (double norm_E, double norm_A, double norm_lambda)
 Constructor.
 
 HArDCore3D::YangMills::YangMills (const DDRCore &ddrcore, const LieAlgebra &liealgebra, bool use_threads, std::ostream &output=std::cout)
 Constructor.
 
void HArDCore3D::YangMills::assembleLinearSystem (double dt)
 Assemble the global system

 
void HArDCore3D::YangMills::setSystemVector (const Eigen::VectorXd &interp_f, const Eigen::VectorXd &interp_dE, const Eigen::VectorXd &interp_A, const Eigen::VectorXd &E_old, const Eigen::VectorXd &A_old, double dt, double theta, double nonlinear_coeff)
 Sets system vector for the nonlinear problem.
 
void HArDCore3D::YangMills::assembleSystemNewton (const Eigen::VectorXd &E_i, const Eigen::VectorXd &A_i, const Eigen::VectorXd &Elambdac_k, const Eigen::VectorXd &sysVec, double dt, double theta, double nonlinear_coeff)
 Assembles the system for Newton iterations.
 
double HArDCore3D::YangMills::stoppingCrit (const Eigen::VectorXd &v, const Eigen::VectorXd &u)
 Stops once changes become small.
 
double HArDCore3D::YangMills::stoppingCrit2 (const Eigen::VectorXd &Elambda_k, const Eigen::VectorXd &E_i, const Eigen::VectorXd &A_i, const Eigen::VectorXd &sysVec, double dt, double theta, double nonlinear_coeff)
 Stops once close enough to system vector for nonlinear problem.
 
void HArDCore3D::YangMills::addBoundaryConditions (const LAMagneticFieldType &curl_A, double dt)
 Adds boundary condition for chosen solution to the system vector.
 
Eigen::VectorXd HArDCore3D::YangMills::computeInitialConditions (const Eigen::MatrixXd &Eh, const Eigen::MatrixXd &Ah, const double nonlinear_coeff, const size_t solver)
 Computes the projected initial conditions that satisfy the constraint.
 
Eigen::MatrixXd HArDCore3D::YangMills::epsBkt_v (size_t iT, boost::multi_array< double, 3 > &ebkt_T, const Eigen::VectorXd &v) const
 Calculates the matrix of *[v,.]^div in element index iT.
 
Eigen::MatrixXd HArDCore3D::YangMills::L2v_Bkt (size_t iT, boost::multi_array< double, 3 > &intPciPcjPgk, const Eigen::VectorXd &v, const size_t &entry) const
 Wrapper to plug vector into L2 integral product: int (Pcurl 1,[Pcurl 2, Pgrad 3])
 
Eigen::MatrixXd HArDCore3D::YangMills::L2v_epsBkt (size_t iT, boost::multi_array< double, 3 > &ebkt_T, const Eigen::VectorXd &v_T, const Eigen::MatrixXd &L2prod) const
 Calculates the bracket inside the L2prod with v.
 
size_t HArDCore3D::YangMills::dimension () const
 Returns the global problem dimension.
 
size_t HArDCore3D::YangMills::sizeSystem () const
 Returns the size of the system.
 
const LieAlgebraHArDCore3D::YangMills::lieAlg () const
 Returns the Lie algebra.
 
const LAXGradHArDCore3D::YangMills::laXGrad () const
 Returns the space LAXGrad.
 
const LAXCurlHArDCore3D::YangMills::laXCurl () const
 Returns the space LAXCurl.
 
const LAXDivHArDCore3D::YangMills::laXDiv () const
 Returns the space LAXCurl.
 
const XGradHArDCore3D::YangMills::xGrad () const
 Returns the space XDiv.
 
const XCurlHArDCore3D::YangMills::xCurl () const
 Returns the space XDiv.
 
const XDivHArDCore3D::YangMills::xDiv () const
 Returns the space XDiv.
 
const SystemMatrixTypeHArDCore3D::YangMills::systemMatrix () const
 Returns the linear system matrix.
 
SystemMatrixTypeHArDCore3D::YangMills::systemMatrix ()
 Returns the linear system matrix.
 
const Eigen::VectorXd & HArDCore3D::YangMills::systemVector () const
 Returns the linear system right-hand side vector.
 
Eigen::VectorXd & HArDCore3D::YangMills::systemVector ()
 Returns the linear system right-hand side vector.
 
const doubleHArDCore3D::YangMills::stabilizationParameter () const
 Returns the stabilization parameter.
 
doubleHArDCore3D::YangMills::stabilizationParameter ()
 Returns the stabilization parameter.
 
std::vector< YangMillsNormsHArDCore3D::YangMills::computeYangMillsNorms (const std::vector< Eigen::VectorXd > &list_dofs) const
 Compute the discrete L2 norms, for a family of Eigen::VectorXd representing the electric field and potential.
 
template<typename outValue , typename TFct >
std::function< outValue(const Eigen::Vector3d &)> HArDCore3D::YangMills::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< FctHArDCore3D::YangMills::contractParaLA (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.
 
double HArDCore3D::YangMills::computeResidual (const Eigen::VectorXd &x) const
 Calculates residual with current matrix and rhs from given x (m_A*x = m_b)
 
Eigen::VectorXd HArDCore3D::YangMills::computeConstraint (const Eigen::VectorXd &E, const Eigen::VectorXd &A, const double nonlinear_coeff)
 Computes the constraint [E, grad P'] + int <E, [A, P']> for the DOFs.
 
double HArDCore3D::YangMills::computeConstraintNorm (const Eigen::VectorXd &Ch, const size_t itersolver) const
 Solves a system to calculate the norm of the constraint (in the dual space of LaXgrad)
 
std::pair< double, doubleHArDCore3D::YangMills::computeConditionNum () const
 
template<typename outValue , typename Fct >
std::vector< FctHArDCore3D::sumLA (const std::vector< Fct > &LAF, const std::vector< Fct > &LAG, double lambda)
 Template to evaluate a vector of functions.
 
 HArDCore3D::AssembleSystems::AssembleSystems (std::vector< Eigen::MatrixXd > sys, std::vector< Eigen::VectorXd > vec)
 Constructor.
 
 HArDCore3D::AssembledSystems::AssembledSystems (std::vector< Eigen::SparseMatrix< double > > sys, std::vector< Eigen::VectorXd > vec)
 Constructor.
 

Variables

double HArDCore3D::YangMillsNorms::E
 Norm of electric field.
 
double HArDCore3D::YangMillsNorms::A
 Norm of potential.
 
double HArDCore3D::YangMillsNorms::lambda
 Norm of Lagrange multiplier.
 
static const double HArDCore3D::PI = boost::math::constants::pi<double>()
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_E1
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_E2
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_E3
 
static YangMills::TLAElectricFieldType HArDCore3D::trigonometric_E = {trigonometric_E1, trigonometric_E2, trigonometric_E3}
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_A1
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_A2
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_A3
 
static YangMills::TLAElectricFieldType HArDCore3D::trigonometric_A = {trigonometric_A1, trigonometric_A2, trigonometric_A3}
 
static YangMills::TMagneticFieldType HArDCore3D::trigonometric_B1_linear
 
static YangMills::TMagneticFieldType HArDCore3D::trigonometric_B2_linear
 
static YangMills::TMagneticFieldType HArDCore3D::trigonometric_B3_linear
 
static YangMills::TLAMagneticFieldType HArDCore3D::trigonometric_B_linear = {trigonometric_B1_linear, trigonometric_B2_linear, trigonometric_B3_linear}
 
static YangMills::TMagneticFieldType HArDCore3D::trigonometric_B1_nonlinear
 
static YangMills::TMagneticFieldType HArDCore3D::trigonometric_B2_nonlinear
 
static YangMills::TMagneticFieldType HArDCore3D::trigonometric_B3_nonlinear
 
static YangMills::TLAMagneticFieldType HArDCore3D::trigonometric_B_nonlinear = {trigonometric_B1_nonlinear, trigonometric_B2_nonlinear, trigonometric_B3_nonlinear}
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_dtE1
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_dtE2
 
static YangMills::TElectricFieldType HArDCore3D::trigonometric_dtE3
 
static YangMills::TLAElectricFieldType HArDCore3D::trigonometric_dtE = {trigonometric_dtE1, trigonometric_dtE2, trigonometric_dtE3}
 
static YangMills::TForcingTermType HArDCore3D::trigonometric_f1_linear
 
static YangMills::TForcingTermType HArDCore3D::trigonometric_f2_linear
 
static YangMills::TForcingTermType HArDCore3D::trigonometric_f3_linear
 
static YangMills::TLAForcingTermType HArDCore3D::trigonometric_f_linear = {trigonometric_f1_linear, trigonometric_f2_linear, trigonometric_f3_linear}
 
static YangMills::TForcingTermType HArDCore3D::trigonometric_f1_nonlinear
 
static YangMills::TForcingTermType HArDCore3D::trigonometric_f2_nonlinear
 
static YangMills::TForcingTermType HArDCore3D::trigonometric_f3_nonlinear
 
static YangMills::TLAForcingTermType HArDCore3D::trigonometric_f_nonlinear = {trigonometric_f1_nonlinear, trigonometric_f2_nonlinear, trigonometric_f3_nonlinear}
 
static YangMills::TElectricFieldType HArDCore3D::linear_E1
 
static YangMills::TElectricFieldType HArDCore3D::linear_E2
 
static YangMills::TElectricFieldType HArDCore3D::linear_E3
 
static YangMills::TLAElectricFieldType HArDCore3D::linear_E = {linear_E1, linear_E2, linear_E3}
 
static YangMills::TElectricFieldType HArDCore3D::linear_A1
 
static YangMills::TElectricFieldType HArDCore3D::linear_A2
 
static YangMills::TElectricFieldType HArDCore3D::linear_A3
 
static YangMills::TLAElectricFieldType HArDCore3D::linear_A = {linear_A1, linear_A2, linear_A3}
 
static YangMills::TMagneticFieldType HArDCore3D::linear_B1_linear
 
static YangMills::TMagneticFieldType HArDCore3D::linear_B2_linear
 
static YangMills::TMagneticFieldType HArDCore3D::linear_B3_linear
 
static YangMills::TLAMagneticFieldType HArDCore3D::linear_B_linear = {linear_B1_linear, linear_B2_linear, linear_B3_linear}
 
static YangMills::TMagneticFieldType HArDCore3D::linear_B1_nonlinear
 
static YangMills::TMagneticFieldType HArDCore3D::linear_B2_nonlinear
 
static YangMills::TMagneticFieldType HArDCore3D::linear_B3_nonlinear
 
static YangMills::TLAMagneticFieldType HArDCore3D::linear_B_nonlinear = {linear_B1_nonlinear, linear_B2_nonlinear, linear_B3_nonlinear}
 
static YangMills::TElectricFieldType HArDCore3D::linear_dtE1
 
static YangMills::TElectricFieldType HArDCore3D::linear_dtE2
 
static YangMills::TElectricFieldType HArDCore3D::linear_dtE3
 
static YangMills::TLAElectricFieldType HArDCore3D::linear_dtE = {linear_dtE1, linear_dtE2, linear_dtE3}
 
static YangMills::TForcingTermType HArDCore3D::linear_f1_linear
 
static YangMills::TForcingTermType HArDCore3D::linear_f2_linear
 
static YangMills::TForcingTermType HArDCore3D::linear_f3_linear
 
static YangMills::TLAForcingTermType HArDCore3D::linear_f_linear = {linear_f1_linear, linear_f2_linear, linear_f3_linear}
 
static YangMills::TForcingTermType HArDCore3D::linear_f1_nonlinear
 
static YangMills::TForcingTermType HArDCore3D::linear_f2_nonlinear
 
static YangMills::TForcingTermType HArDCore3D::linear_f3_nonlinear
 
static YangMills::TLAForcingTermType HArDCore3D::linear_f_nonlinear = {linear_f1_nonlinear, linear_f2_nonlinear, linear_f3_nonlinear}
 
static YangMills::TElectricFieldType HArDCore3D::const_E1
 
static YangMills::TElectricFieldType HArDCore3D::const_E2
 
static YangMills::TElectricFieldType HArDCore3D::const_E3
 
static YangMills::TLAElectricFieldType HArDCore3D::const_E = {const_E1, const_E2, const_E3}
 
static YangMills::TElectricFieldType HArDCore3D::const_A1
 
static YangMills::TElectricFieldType HArDCore3D::const_A2
 
static YangMills::TElectricFieldType HArDCore3D::const_A3
 
static YangMills::TLAElectricFieldType HArDCore3D::const_A = {const_A1, const_A2, const_A3}
 
static YangMills::TMagneticFieldType HArDCore3D::const_B1_linear
 
static YangMills::TMagneticFieldType HArDCore3D::const_B2_linear
 
static YangMills::TMagneticFieldType HArDCore3D::const_B3_linear
 
static YangMills::TLAMagneticFieldType HArDCore3D::const_B_linear = {const_B1_linear, const_B2_linear, const_B3_linear}
 
static YangMills::TMagneticFieldType HArDCore3D::const_B1_nonlinear
 
static YangMills::TMagneticFieldType HArDCore3D::const_B2_nonlinear
 
static YangMills::TMagneticFieldType HArDCore3D::const_B3_nonlinear
 
static YangMills::TLAMagneticFieldType HArDCore3D::const_B_nonlinear = {const_B1_nonlinear, const_B2_nonlinear, const_B3_nonlinear}
 
static YangMills::TElectricFieldType HArDCore3D::const_dtE1
 
static YangMills::TElectricFieldType HArDCore3D::const_dtE2
 
static YangMills::TElectricFieldType HArDCore3D::const_dtE3
 
static YangMills::TLAElectricFieldType HArDCore3D::const_dtE = {const_dtE1, const_dtE2, const_dtE3}
 
static YangMills::TForcingTermType HArDCore3D::const_f1_linear
 
static YangMills::TForcingTermType HArDCore3D::const_f2_linear
 
static YangMills::TForcingTermType HArDCore3D::const_f3_linear
 
static YangMills::TLAForcingTermType HArDCore3D::const_f_linear = {const_f1_linear, const_f2_linear, const_f3_linear}
 
static YangMills::TForcingTermType HArDCore3D::const_f1_nonlinear
 
static YangMills::TForcingTermType HArDCore3D::const_f2_nonlinear
 
static YangMills::TForcingTermType HArDCore3D::const_f3_nonlinear
 
static YangMills::TLAForcingTermType HArDCore3D::const_f_nonlinear = {const_f1_nonlinear, const_f2_nonlinear, const_f3_nonlinear}
 
std::vector< Eigen::MatrixXd > HArDCore3D::AssembleSystems::systems
 
std::vector< Eigen::VectorXd > HArDCore3D::AssembleSystems::vectors
 
std::vector< Eigen::SparseMatrix< double > > HArDCore3D::AssembledSystems::systems
 
std::vector< Eigen::VectorXd > HArDCore3D::AssembledSystems::vectors
 

Detailed Description

Implementation of the lowest order DDR scheme for the Yang-Mills problem.

Typedef Documentation

◆ ElectricFieldType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::YangMills::ElectricFieldType

◆ ForcingTermType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::YangMills::ForcingTermType

Base functions.

◆ LAElectricFieldType

◆ LAForcingTermType

Lie algebra valued functions (vector of dim Lie algebra)

◆ LAMagneticFieldType

◆ MagneticFieldType

typedef std::function<Eigen::Vector3d(const Eigen::Vector3d &)> HArDCore3D::YangMills::MagneticFieldType

◆ SystemMatrixType

◆ TElectricFieldType

typedef std::function<Eigen::Vector3d(const double &, const Eigen::Vector3d &)> HArDCore3D::YangMills::TElectricFieldType

◆ TForcingTermType

typedef std::function<Eigen::Vector3d(const double &, const Eigen::Vector3d &)> HArDCore3D::YangMills::TForcingTermType

Base functions with time component.

◆ TLAElectricFieldType

◆ TLAForcingTermType

Lie algebra valued functions with time component.

◆ TLAMagneticFieldType

◆ TMagneticFieldType

typedef std::function<Eigen::Vector3d(const double &, const Eigen::Vector3d &)> HArDCore3D::YangMills::TMagneticFieldType

Function Documentation

◆ addBoundaryConditions()

void YangMills::addBoundaryConditions ( const LAMagneticFieldType curl_A,
double  dt 
)

Adds boundary condition for chosen solution to the system vector.

Parameters
curl_ACurl of the potential
dtTime step

◆ AssembledSystems()

HArDCore3D::AssembledSystems::AssembledSystems ( std::vector< Eigen::SparseMatrix< double > >  sys,
std::vector< Eigen::VectorXd >  vec 
)
inline

Constructor.

◆ assembleLinearSystem()

void YangMills::assembleLinearSystem ( double  dt)

Assemble the global system

Parameters
dtTime step

◆ assembleSystemNewton()

void YangMills::assembleSystemNewton ( const Eigen::VectorXd &  E_i,
const Eigen::VectorXd &  A_i,
const Eigen::VectorXd &  Elambdac_k,
const Eigen::VectorXd &  sysVec,
double  dt,
double  theta,
double  nonlinear_coeff 
)

Assembles the system for Newton iterations.

Parameters
E_iElectric field at previous time
A_iPotential at previous time
Elambdac_kSolution position after previous Newton iteration
sysVecSystem vector of nonlinear problem
dtTime step
thetaTheta scheme
nonlinear_coeffScaling of the nonlinear terms (0 for Maxwell)

◆ AssembleSystems()

HArDCore3D::AssembleSystems::AssembleSystems ( std::vector< Eigen::MatrixXd >  sys,
std::vector< Eigen::VectorXd >  vec 
)
inline

Constructor.

◆ computeConditionNum()

std::pair< double, double > YangMills::computeConditionNum ( ) const

◆ computeConstraint()

Eigen::VectorXd YangMills::computeConstraint ( const Eigen::VectorXd &  E,
const Eigen::VectorXd &  A,
const double  nonlinear_coeff 
)

Computes the constraint [E, grad P'] + int <E, [A, P']> for the DOFs.

◆ computeConstraintNorm()

double YangMills::computeConstraintNorm ( const Eigen::VectorXd &  Ch,
const size_t  itersolver 
) const

Solves a system to calculate the norm of the constraint (in the dual space of LaXgrad)

◆ computeInitialConditions()

Eigen::VectorXd YangMills::computeInitialConditions ( const Eigen::MatrixXd &  Eh,
const Eigen::MatrixXd &  Ah,
const double  nonlinear_coeff,
const size_t  solver 
)

Computes the projected initial conditions that satisfy the constraint.

Parameters
EhInterpolate of the initial electric field
AhInterpolate of the initial potential field
nonlinear_coeffScaling of the nonlinear terms (0 for Maxwell)
solverSolver to use

◆ computeResidual()

double YangMills::computeResidual ( const Eigen::VectorXd &  x) const

Calculates residual with current matrix and rhs from given x (m_A*x = m_b)

◆ computeYangMillsNorms()

std::vector< YangMillsNorms > YangMills::computeYangMillsNorms ( const std::vector< Eigen::VectorXd > &  list_dofs) const

Compute the discrete L2 norms, for a family of Eigen::VectorXd representing the electric field and potential.

Parameters
list_dofsList of vectors representing the electric field and potential

◆ contractPara()

std::function< outValue(const Eigen::Vector3d &)> HArDCore3D::YangMills::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.

◆ contractParaLA()

std::vector< Fct > HArDCore3D::YangMills::contractParaLA ( 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.

◆ dimension()

size_t HArDCore3D::YangMills::dimension ( ) const
inline

Returns the global problem dimension.

◆ epsBkt_v()

Eigen::MatrixXd YangMills::epsBkt_v ( size_t  iT,
boost::multi_array< double, 3 > &  ebkt_T,
const Eigen::VectorXd &  v 
) const

Calculates the matrix of *[v,.]^div in element index iT.

Parameters
iTElement index
ebkt_TThe *[.,.]^div operator for cell values (doesn't include faces)
vVector to plug in

◆ L2v_Bkt()

Eigen::MatrixXd YangMills::L2v_Bkt ( size_t  iT,
boost::multi_array< double, 3 > &  intPciPcjPgk,
const Eigen::VectorXd &  v,
const size_t entry 
) const

Wrapper to plug vector into L2 integral product: int (Pcurl 1,[Pcurl 2, Pgrad 3])

Parameters
iTElement index
intPciPcjPgkL2 integral product (trilinear form)
vVector to plug in
entryPosition to put vector (1, 2, 3)

◆ L2v_epsBkt()

Eigen::MatrixXd YangMills::L2v_epsBkt ( size_t  iT,
boost::multi_array< double, 3 > &  ebkt_T,
const Eigen::VectorXd &  v_T,
const Eigen::MatrixXd &  L2prod 
) const

Calculates the bracket inside the L2prod with v.

Parameters
iTElement index
ebkt_TThe *[.,.]^div operator for cell values (doesn't include faces)
v_TLocal vector in T (because we often only have the local part of the vector e.g. (*[A,A]^div)_T )
L2prodL2 product to be combined with bracket

◆ laXCurl()

const LAXCurl & HArDCore3D::YangMills::laXCurl ( ) const
inline

Returns the space LAXCurl.

◆ laXDiv()

const LAXDiv & HArDCore3D::YangMills::laXDiv ( ) const
inline

Returns the space LAXCurl.

◆ laXGrad()

const LAXGrad & HArDCore3D::YangMills::laXGrad ( ) const
inline

Returns the space LAXGrad.

◆ lieAlg()

const LieAlgebra & HArDCore3D::YangMills::lieAlg ( ) const
inline

Returns the Lie algebra.

◆ setSystemVector()

void YangMills::setSystemVector ( const Eigen::VectorXd &  interp_f,
const Eigen::VectorXd &  interp_dE,
const Eigen::VectorXd &  interp_A,
const Eigen::VectorXd &  E_old,
const Eigen::VectorXd &  A_old,
double  dt,
double  theta,
double  nonlinear_coeff 
)

Sets system vector for the nonlinear problem.

Parameters
interp_fForcing term of first equation
interp_dEUsed in forcing term of second equation
interp_AUsed in forcing term of second equation
E_oldElectric field at previous time
A_oldPotential at previous time
dtTime step
thetaTheta scheme
nonlinear_coeffScaling of the nonlinear terms (0 for Maxwell)

◆ sizeSystem()

size_t HArDCore3D::YangMills::sizeSystem ( ) const
inline

Returns the size of the system.

◆ stabilizationParameter() [1/2]

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

Returns the stabilization parameter.

◆ stabilizationParameter() [2/2]

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

Returns the stabilization parameter.

◆ stoppingCrit()

double YangMills::stoppingCrit ( const Eigen::VectorXd &  v,
const Eigen::VectorXd &  u 
)

Stops once changes become small.

◆ stoppingCrit2()

double YangMills::stoppingCrit2 ( const Eigen::VectorXd &  Elambda_k,
const Eigen::VectorXd &  E_i,
const Eigen::VectorXd &  A_i,
const Eigen::VectorXd &  sysVec,
double  dt,
double  theta,
double  nonlinear_coeff 
)

Stops once close enough to system vector for nonlinear problem.

Parameters
Elambda_kSolution in Newton iterations
E_iElectric field at previous time
A_iPotential at previous time
sysVecSystem vector of nonlinear problem
dtTime step
thetaTheta scheme
nonlinear_coeffScaling of the nonlinear terms (0 for Maxwell)

◆ sumLA()

template<typename outValue , typename Fct >
std::vector< Fct > HArDCore3D::sumLA ( const std::vector< Fct > &  LAF,
const std::vector< Fct > &  LAG,
double  lambda 
)

Template to evaluate a vector of functions.

◆ systemMatrix() [1/2]

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

Returns the linear system matrix.

◆ systemMatrix() [2/2]

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

Returns the linear system matrix.

◆ systemVector() [1/2]

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

Returns the linear system right-hand side vector.

◆ systemVector() [2/2]

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

Returns the linear system right-hand side vector.

◆ xCurl()

const XCurl & HArDCore3D::YangMills::xCurl ( ) const
inline

Returns the space XDiv.

◆ xDiv()

const XDiv & HArDCore3D::YangMills::xDiv ( ) const
inline

Returns the space XDiv.

◆ xGrad()

const XGrad & HArDCore3D::YangMills::xGrad ( ) const
inline

Returns the space XDiv.

◆ YangMills()

YangMills::YangMills ( const DDRCore ddrcore,
const LieAlgebra liealgebra,
bool  use_threads,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
ddrcoreCore for the DDR space sequence
liealgebraLie algebra of the problem
use_threadsTrue for parallel execution, false for sequential execution
outputOutput stream to print status messages

◆ YangMillsNorms()

HArDCore3D::YangMillsNorms::YangMillsNorms ( double  norm_E,
double  norm_A,
double  norm_lambda 
)
inline

Constructor.

Variable Documentation

◆ A

double HArDCore3D::YangMillsNorms::A

Norm of potential.

◆ const_A

YangMills::TLAElectricFieldType HArDCore3D::const_A = {const_A1, const_A2, const_A3}
static

◆ const_A1

YangMills::TElectricFieldType HArDCore3D::const_A1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0.2, 0.4, 0.6);
}

◆ const_A2

YangMills::TElectricFieldType HArDCore3D::const_A2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0.8, 1.0, 1.2);
}

◆ const_A3

YangMills::TElectricFieldType HArDCore3D::const_A3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(1.4, 1.6, 1.8);
}

◆ const_B1_linear

YangMills::TMagneticFieldType HArDCore3D::const_B1_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_B1_nonlinear

YangMills::TMagneticFieldType HArDCore3D::const_B1_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-0.12, 0.24, -0.12);
}

◆ const_B2_linear

YangMills::TMagneticFieldType HArDCore3D::const_B2_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_B2_nonlinear

YangMills::TMagneticFieldType HArDCore3D::const_B2_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0.24, -0.48, 0.24);
}

◆ const_B3_linear

YangMills::TMagneticFieldType HArDCore3D::const_B3_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_B3_nonlinear

YangMills::TMagneticFieldType HArDCore3D::const_B3_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-0.12, 0.24, -0.12);
}

◆ const_B_linear

◆ const_B_nonlinear

◆ const_dtE

◆ const_dtE1

YangMills::TElectricFieldType HArDCore3D::const_dtE1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_dtE2

YangMills::TElectricFieldType HArDCore3D::const_dtE2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_dtE3

YangMills::TElectricFieldType HArDCore3D::const_dtE3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_E

YangMills::TLAElectricFieldType HArDCore3D::const_E = {const_E1, const_E2, const_E3}
static

◆ const_E1

YangMills::TElectricFieldType HArDCore3D::const_E1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ const_E2

YangMills::TElectricFieldType HArDCore3D::const_E2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ const_E3

YangMills::TElectricFieldType HArDCore3D::const_E3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d::Zero();
}

◆ const_f1_linear

YangMills::TForcingTermType HArDCore3D::const_f1_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_f1_nonlinear

YangMills::TForcingTermType HArDCore3D::const_f1_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(1.656, 0.144, -1.368);
}

◆ const_f2_linear

YangMills::TForcingTermType HArDCore3D::const_f2_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_f2_nonlinear

YangMills::TForcingTermType HArDCore3D::const_f2_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0.432, 0, -0.432);
}

◆ const_f3_linear

YangMills::TForcingTermType HArDCore3D::const_f3_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ const_f3_nonlinear

YangMills::TForcingTermType HArDCore3D::const_f3_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-0.792, -0.144, 0.504);
}

◆ const_f_linear

◆ const_f_nonlinear

◆ E

double HArDCore3D::YangMillsNorms::E

Norm of electric field.

◆ lambda

double HArDCore3D::YangMillsNorms::lambda

Norm of Lagrange multiplier.

◆ linear_A

YangMills::TLAElectricFieldType HArDCore3D::linear_A = {linear_A1, linear_A2, linear_A3}
static

◆ linear_A1

YangMills::TElectricFieldType HArDCore3D::linear_A1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(t*x(2), t*x(1), t*x(0));
}

◆ linear_A2

YangMills::TElectricFieldType HArDCore3D::linear_A2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(t*x(2), t*x(0), t*x(1));
}

◆ linear_A3

YangMills::TElectricFieldType HArDCore3D::linear_A3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(t*x(1) + t*x(2), t*x(0) + t*x(2), t*x(0) + t*x(1));
}

◆ linear_B1_linear

YangMills::TMagneticFieldType HArDCore3D::linear_B1_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0,
0,
0
);
}

◆ linear_B1_nonlinear

YangMills::TMagneticFieldType HArDCore3D::linear_B1_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2)),
1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1)),
-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))
);
}

◆ linear_B2_linear

YangMills::TMagneticFieldType HArDCore3D::linear_B2_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
t,
t,
t
);
}

◆ linear_B2_nonlinear

YangMills::TMagneticFieldType HArDCore3D::linear_B2_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)),
-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)),
1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2))
);
}

◆ linear_B3_linear

YangMills::TMagneticFieldType HArDCore3D::linear_B3_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0,
0,
0
);
}

◆ linear_B3_nonlinear

YangMills::TMagneticFieldType HArDCore3D::linear_B3_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2),
1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2),
1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)
);
}

◆ linear_B_linear

◆ linear_B_nonlinear

◆ linear_dtE

◆ linear_dtE1

YangMills::TElectricFieldType HArDCore3D::linear_dtE1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ linear_dtE2

YangMills::TElectricFieldType HArDCore3D::linear_dtE2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ linear_dtE3

YangMills::TElectricFieldType HArDCore3D::linear_dtE3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ linear_E

YangMills::TLAElectricFieldType HArDCore3D::linear_E = {linear_E1, linear_E2, linear_E3}
static

◆ linear_E1

YangMills::TElectricFieldType HArDCore3D::linear_E1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-x(2), -x(1), -x(0));
}

◆ linear_E2

YangMills::TElectricFieldType HArDCore3D::linear_E2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-x(2), -x(0), -x(1));
}

◆ linear_E3

YangMills::TElectricFieldType HArDCore3D::linear_E3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-x(1) - x(2), -x(0) - x(2), -x(0) - x(1));
}

◆ linear_f1_linear

YangMills::TForcingTermType HArDCore3D::linear_f1_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ linear_f1_nonlinear

YangMills::TForcingTermType HArDCore3D::linear_f1_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(1.0*pow(t, 2)*x(0) + 1.0*pow(t, 2)*x(1) - t*x(0)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) + t*x(1)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - 1.0*t*(t*x(0) + t*x(1)) - (t*x(0) + t*x(1))*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)) + t) + (t*x(0) + t*x(2))*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2)) + t), 1.0*pow(t, 2)*x(1) + 1.0*pow(t, 2)*x(2) - t*x(1)*(-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2)) + t*x(2)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - 1.0*t*(t*x(1) + t*x(2)) + (t*x(0) + t*x(1))*(1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)) + t) - (t*x(1) + t*x(2))*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2)) + t), 1.0*pow(t, 2)*x(0) + 1.0*pow(t, 2)*x(2) + t*x(0)*(-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2)) - t*x(2)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - 1.0*t*(t*x(0) + t*x(2)) - (t*x(0) + t*x(2))*(1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)) + t) + (t*x(1) + t*x(2))*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)) + t)
);
}

◆ linear_f2_linear

YangMills::TForcingTermType HArDCore3D::linear_f2_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ linear_f2_nonlinear

YangMills::TForcingTermType HArDCore3D::linear_f2_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-1.0*pow(t, 2)*x(0) - 1.0*pow(t, 2)*x(1) - t*x(0)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) + t*x(1)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) + 1.0*t*(t*x(0) + t*x(1)) - 1.0*t*(t*x(1) + t*x(2)) + (t*x(0) + t*x(1))*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1))) - (t*x(0) + t*x(2))*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))), -1.0*pow(t, 2)*x(0) - 1.0*pow(t, 2)*x(2) + t*x(0)*(-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2)) - t*x(2)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - (t*x(0) + t*x(1))*(1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2))) + (t*x(1) + t*x(2))*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))), -1.0*pow(t, 2)*x(1) - 1.0*pow(t, 2)*x(2) - t*x(1)*(-1.0*pow(t, 2)*pow(x(0), 2) + 1.0*pow(t, 2)*pow(x(1), 2)) + t*x(2)*(1.0*pow(t, 2)*x(0)*x(2) - 1.0*pow(t, 2)*x(1)*x(2)) - 1.0*t*(t*x(0) + t*x(1)) + 1.0*t*(t*x(1) + t*x(2)) + (t*x(0) + t*x(2))*(1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2))) - (t*x(1) + t*x(2))*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1)))
);
}

◆ linear_f3_linear

YangMills::TForcingTermType HArDCore3D::linear_f3_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(0, 0, 0);
}

◆ linear_f3_nonlinear

YangMills::TForcingTermType HArDCore3D::linear_f3_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(1.0*pow(t, 2)*x(0) - 1.0*pow(t, 2)*x(1) + 1.0*pow(t, 2)*x(2) + t*x(0)*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))) + t*x(0)*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)) + t) - t*x(1)*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1))) - t*x(1)*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2)) + t), 1.0*pow(t, 2)*x(2) - t*x(0)*(1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)) + t) + t*x(1)*(1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2))) - t*x(2)*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(2))) + t*x(2)*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(2)) + t), 2.0*pow(t, 2)*x(1) - 1.0*pow(t, 2)*x(2) - t*x(0)*(1.0*t*x(0)*(t*x(0) + t*x(1)) - 1.0*t*x(1)*(t*x(0) + t*x(2))) + t*x(1)*(1.0*t*x(0)*(t*x(0) + t*x(2)) - 1.0*t*x(1)*(t*x(0) + t*x(1)) + t) + t*x(2)*(1.0*t*x(1)*(t*x(1) + t*x(2)) - 1.0*t*x(2)*(t*x(0) + t*x(1))) - t*x(2)*(-1.0*t*x(0)*(t*x(1) + t*x(2)) + 1.0*t*x(2)*(t*x(0) + t*x(1)) + t)
);
}

◆ linear_f_linear

◆ linear_f_nonlinear

◆ PI

const double HArDCore3D::PI = boost::math::constants::pi<double>()
static

◆ systems [1/2]

std::vector<Eigen::MatrixXd> HArDCore3D::AssembleSystems::systems

◆ systems [2/2]

std::vector<Eigen::SparseMatrix<double> > HArDCore3D::AssembledSystems::systems

◆ trigonometric_A

◆ trigonometric_A1

YangMills::TElectricFieldType HArDCore3D::trigonometric_A1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)),
sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
-0.5*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))
);
}

◆ trigonometric_A2

YangMills::TElectricFieldType HArDCore3D::trigonometric_A2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)),
sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)),
-0.5*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1))
);
}

◆ trigonometric_A3

YangMills::TElectricFieldType HArDCore3D::trigonometric_A3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(-0.5*sin(t)*pow(sin(PI*x(1)), 2),
cos(t)*pow(cos(PI*x(2)), 2),
-0.5*sin(t)*pow(cos(PI*x(0)), 2)
);
}

◆ trigonometric_B1_linear

YangMills::TMagneticFieldType HArDCore3D::trigonometric_B1_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)),
0,
-1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2))
);
}

◆ trigonometric_B1_nonlinear

YangMills::TMagneticFieldType HArDCore3D::trigonometric_B1_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2),
-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)),
0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3)
);
}

◆ trigonometric_B2_linear

YangMills::TMagneticFieldType HArDCore3D::trigonometric_B2_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)),
0,
-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2))
);
}

◆ trigonometric_B2_nonlinear

YangMills::TMagneticFieldType HArDCore3D::trigonometric_B2_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2),
0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)),
-0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3)
);
}

◆ trigonometric_B3_linear

YangMills::TMagneticFieldType HArDCore3D::trigonometric_B3_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
2*PI*sin(PI*x(2))*cos(t)*cos(PI*x(2)),
-1.0*PI*sin(t)*sin(PI*x(0))*cos(PI*x(0)),
1.0*PI*sin(t)*sin(PI*x(1))*cos(PI*x(1))
);
}

◆ trigonometric_B3_nonlinear

YangMills::TMagneticFieldType HArDCore3D::trigonometric_B3_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0,
0,
0
);
}

◆ trigonometric_B_linear

◆ trigonometric_B_nonlinear

◆ trigonometric_dtE

◆ trigonometric_dtE1

YangMills::TElectricFieldType HArDCore3D::trigonometric_dtE1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)),
sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
-0.5*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))
);
}

◆ trigonometric_dtE2

YangMills::TElectricFieldType HArDCore3D::trigonometric_dtE2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)),
sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)),
-0.5*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1))
);
}

◆ trigonometric_dtE3

YangMills::TElectricFieldType HArDCore3D::trigonometric_dtE3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(t)*pow(sin(PI*x(1)), 2),
cos(t)*pow(cos(PI*x(2)), 2),
-0.5*sin(t)*pow(cos(PI*x(0)), 2)
);
}

◆ trigonometric_E

◆ trigonometric_E1

YangMills::TElectricFieldType HArDCore3D::trigonometric_E1
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)),
sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)),
-0.5*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1))
);
}

◆ trigonometric_E2

YangMills::TElectricFieldType HArDCore3D::trigonometric_E2
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0.5*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)),
-sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
0.5*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))
);
}

◆ trigonometric_E3

YangMills::TElectricFieldType HArDCore3D::trigonometric_E3
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0.5*pow(sin(PI*x(1)), 2)*cos(t),
sin(t)*pow(cos(PI*x(2)), 2),
0.5*cos(t)*pow(cos(PI*x(0)), 2)
);
}

◆ trigonometric_f1_linear

YangMills::TForcingTermType HArDCore3D::trigonometric_f1_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)) + 1.5*pow(PI, 2)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)),
-3.0*pow(PI, 2)*sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
-0.5*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)) + 1.5*pow(PI, 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))
);
}

◆ trigonometric_f1_nonlinear

YangMills::TForcingTermType HArDCore3D::trigonometric_f1_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0.5*(0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*pow(cos(PI*x(0)), 2) + (-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2)) - 0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3))*cos(t)*pow(cos(PI*x(2)), 2) + 0.75*PI*pow(sin(t), 2)*sin(PI*x(0))*sin(PI*x(2))*pow(cos(PI*x(0)), 2)*cos(PI*x(1)) - 2.25*PI*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) - 0.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(t)*pow(cos(PI*x(2)), 3),
0.5*(-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2)) - 0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3))*sin(t)*pow(sin(PI*x(1)), 2) - 0.5*(1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)) + 0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2))*sin(t)*pow(cos(PI*x(0)), 2) - 0.5*PI*pow(sin(t), 2)*sin(PI*x(0))*pow(sin(PI*x(1)), 3)*cos(PI*x(2)) - 0.5*PI*pow(sin(t), 2)*sin(PI*x(0))*sin(PI*x(1))*pow(cos(PI*x(1)), 2)*cos(PI*x(2)) - 0.5*PI*pow(sin(t), 2)*sin(PI*x(1))*sin(PI*x(2))*pow(cos(PI*x(0)), 3) + 2.0*PI*sin(t)*pow(sin(PI*x(2)), 2)*cos(t)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) - 1.0*PI*sin(t)*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 3),
-0.5*(0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*pow(sin(PI*x(1)), 2) - (1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)) + 0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2))*cos(t)*pow(cos(PI*x(2)), 2) - 1.0*PI*pow(sin(t), 2)*pow(sin(PI*x(0)), 2)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 0.25*PI*pow(sin(t), 2)*sin(PI*x(0))*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(1)) - 0.25*PI*pow(sin(t), 2)*pow(cos(PI*x(0)), 3)*cos(PI*x(1))*cos(PI*x(2)) + 1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0))*pow(cos(PI*x(2)), 2)
);
}

◆ trigonometric_f2_linear

YangMills::TForcingTermType HArDCore3D::trigonometric_f2_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 1.5*pow(PI, 2)*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)),
-3.0*pow(PI, 2)*sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)) + sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)),
-0.5*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)) + 1.5*pow(PI, 2)*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1))
);
}

◆ trigonometric_f2_nonlinear

YangMills::TForcingTermType HArDCore3D::trigonometric_f2_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*(-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*pow(cos(PI*x(0)), 2) - (0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3) - 1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2)))*cos(t)*pow(cos(PI*x(2)), 2) - 0.75*PI*sin(t)*sin(PI*x(0))*sin(PI*x(2))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1)) + 2.25*PI*sin(t)*pow(sin(PI*x(1)), 2)*cos(t)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 0.5*PI*sin(PI*x(0))*sin(PI*x(1))*pow(cos(t), 2)*pow(cos(PI*x(2)), 3),
0.5*(-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2) + 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)))*sin(t)*pow(cos(PI*x(0)), 2) - 0.5*(0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3) - 1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2)))*sin(t)*pow(sin(PI*x(1)), 2) + 0.5*PI*sin(t)*sin(PI*x(0))*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(2)) + 0.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(t)*pow(cos(PI*x(1)), 2)*cos(PI*x(2)) + 0.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(t)*pow(cos(PI*x(0)), 3) - 2.0*PI*pow(sin(PI*x(2)), 2)*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 1.0*PI*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 3),
0.5*(-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*pow(sin(PI*x(1)), 2) + (-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2) + 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)))*cos(t)*pow(cos(PI*x(2)), 2) + 1.0*PI*sin(t)*pow(sin(PI*x(0)), 2)*cos(t)*cos(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) - 0.25*PI*sin(t)*sin(PI*x(0))*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(1)) + 0.25*PI*sin(t)*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(1))*cos(PI*x(2)) - 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*pow(cos(PI*x(2)), 2)
);
}

◆ trigonometric_f3_linear

YangMills::TForcingTermType HArDCore3D::trigonometric_f3_linear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
-0.5*sin(t)*pow(sin(PI*x(1)), 2) + 1.0*pow(PI, 2)*sin(t)*pow(sin(PI*x(1)), 2) - 1.0*pow(PI, 2)*sin(t)*pow(cos(PI*x(1)), 2),
2*pow(PI, 2)*pow(sin(PI*x(2)), 2)*cos(t) - 2*pow(PI, 2)*cos(t)*pow(cos(PI*x(2)), 2) + cos(t)*pow(cos(PI*x(2)), 2),
-1.0*pow(PI, 2)*sin(t)*pow(sin(PI*x(0)), 2) - 0.5*sin(t)*pow(cos(PI*x(0)), 2) + 1.0*pow(PI, 2)*sin(t)*pow(cos(PI*x(0)), 2)
);
}

◆ trigonometric_f3_nonlinear

YangMills::TForcingTermType HArDCore3D::trigonometric_f3_nonlinear
static
Initial value:
= [](const double t, const Eigen::Vector3d & x) -> Eigen::Vector3d {
return Eigen::Vector3d(
0.5*(-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)) - 0.5*(0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)))*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)) + (0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3) - 1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2)))*sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)) - (-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2)) - 0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3))*sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)),
-0.5*(-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2) + 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)))*sin(t)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)) + 0.5*(0.5*pow(sin(t), 2)*pow(sin(PI*x(1)), 3)*cos(PI*x(0))*cos(PI*x(2)) - 0.5*sin(t)*sin(PI*x(0))*cos(t)*cos(PI*x(1))*pow(cos(PI*x(2)), 3) - 1.5*PI*sin(PI*x(0))*sin(PI*x(1))*cos(t)*cos(PI*x(2)))*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) - 0.5*(-1.5*PI*sin(t)*sin(PI*x(0))*sin(PI*x(1))*cos(PI*x(2)) - 0.5*sin(t)*pow(sin(PI*x(1)), 3)*cos(t)*cos(PI*x(0))*cos(PI*x(2)) + 0.5*sin(PI*x(0))*pow(cos(t), 2)*cos(PI*x(1))*pow(cos(PI*x(2)), 3))*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)) + 0.5*(1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)) + 0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2))*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)),
-0.5*(-0.25*pow(sin(t), 2)*sin(PI*x(0))*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) + 0.25*pow(sin(t), 2)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(PI*x(0))*cos(PI*x(1)))*sin(t)*sin(PI*x(0))*cos(PI*x(1))*cos(PI*x(2)) + 0.5*(0.25*sin(t)*sin(PI*x(0))*cos(t)*pow(cos(PI*x(0)), 2)*cos(PI*x(1))*cos(PI*x(2)) - 0.25*sin(t)*pow(sin(PI*x(1)), 2)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1)))*sin(PI*x(0))*cos(t)*cos(PI*x(1))*cos(PI*x(2)) - (-0.5*pow(sin(t), 2)*sin(PI*x(1))*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) + 0.5*sin(t)*sin(PI*x(2))*cos(t)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2) + 1.5*PI*sin(PI*x(1))*sin(PI*x(2))*cos(t)*cos(PI*x(0)))*sin(t)*sin(PI*x(1))*cos(PI*x(0))*cos(PI*x(2)) + (1.5*PI*sin(t)*sin(PI*x(1))*sin(PI*x(2))*cos(PI*x(0)) + 0.5*sin(t)*sin(PI*x(1))*cos(t)*pow(cos(PI*x(0)), 3)*cos(PI*x(2)) - 0.5*sin(PI*x(2))*pow(cos(t), 2)*cos(PI*x(0))*cos(PI*x(1))*pow(cos(PI*x(2)), 2))*sin(PI*x(1))*cos(t)*cos(PI*x(0))*cos(PI*x(2)));
}

◆ trigonometric_f_linear

◆ trigonometric_f_nonlinear

◆ vectors [1/2]

std::vector<Eigen::VectorXd> HArDCore3D::AssembleSystems::vectors

◆ vectors [2/2]

std::vector<Eigen::VectorXd> HArDCore3D::AssembledSystems::vectors