HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
Classes | Typedefs | Functions | Variables
DDR_rmplate

Implementation of the DDR scheme for the Reissner-Mindlin plate problem. More...

Classes

struct  HArDCore2D::RMParameters
 Structure to store model data. More...
 
struct  HArDCore2D::RMNorms
 Structure to store component norms (for rotation, displacement, Kirchoff term and total energy) More...
 
struct  HArDCore2D::ReissnerMindlin
 Assemble a RM problem. More...
 

Typedefs

typedef Eigen::SparseMatrix< double > HArDCore2D::ReissnerMindlin::SystemMatrixType
 
typedef std::function< double(const RMParameters &, const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::ForcingTermType
 
typedef std::function< Eigen::Vector2d(const RMParameters &, const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::SolutionRotationType
 
typedef std::function< Eigen::Matrix2d(const RMParameters &, const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::GradientRotationType
 
typedef std::function< double(const RMParameters &, const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::SolutionDisplacementType
 

Functions

 HArDCore2D::RMParameters::RMParameters (const double thickness, const double young_modulus, const double poisson_ratio)
 Constructor. More...
 
 HArDCore2D::RMNorms::RMNorms (double norm_rotation, double norm_displacement, double norm_kirchoff)
 Constructor. More...
 
 HArDCore2D::ReissnerMindlin::ReissnerMindlin (const DDRCore &ddrcore, const RMParameters &para, const BoundaryConditions &BC_theta, const BoundaryConditions &BC_u, bool use_threads, std::ostream &output=std::cout)
 Constructor. More...
 
void HArDCore2D::ReissnerMindlin::assembleLinearSystem (const ForcingTermType &f, const SolutionRotationType &theta, const GradientRotationType &grad_theta, const SolutionDisplacementType &u)
 Assemble the global system
More...
 
size_t HArDCore2D::ReissnerMindlin::dimensionSpace () const
 Returns the dimension of the rotation + displacement space (with BC) More...
 
size_t HArDCore2D::ReissnerMindlin::nb_bdryDOFs () const
 Returns the nb of DOFs for BC. More...
 
size_t HArDCore2D::ReissnerMindlin::sizeSystem () const
 Returns the size of the system without BC. More...
 
const std::vector< std::pair< size_t, size_t > > & HArDCore2D::ReissnerMindlin::locUKN () const
 Returns the location of the unknowns among the DOFs. More...
 
std::vector< size_t > HArDCore2D::ReissnerMindlin::globalDOFIndices (const Cell &T) const
 Create the vector of DOF indices for cell T, which combines the DOFs for the spaces EXcurl and Xgrad. More...
 
const RMParametersHArDCore2D::ReissnerMindlin::para () const
 Returns the parameters. More...
 
const EXCurlHArDCore2D::ReissnerMindlin::exCurl () const
 Returns the space EXCurl. More...
 
const XGradHArDCore2D::ReissnerMindlin::xGrad () const
 Returns the space XGrad. More...
 
const SystemMatrixTypeHArDCore2D::ReissnerMindlin::systemMatrix () const
 Returns the linear system matrix. More...
 
SystemMatrixTypeHArDCore2D::ReissnerMindlin::systemMatrix ()
 Returns the linear system matrix. More...
 
const Eigen::VectorXd & HArDCore2D::ReissnerMindlin::systemVector () const
 Returns the linear system right-hand side vector. More...
 
Eigen::VectorXd & HArDCore2D::ReissnerMindlin::systemVector ()
 Returns the linear system right-hand side vector. More...
 
const SystemMatrixTypeHArDCore2D::ReissnerMindlin::bdryMatrix () const
 Returns the Matrix for BC. More...
 
const Eigen::VectorXd & HArDCore2D::ReissnerMindlin::bdryValues () const
 Returns the boundary values. More...
 
const double & HArDCore2D::ReissnerMindlin::stabilizationParameter () const
 Returns the stabilization parameter. More...
 
double & HArDCore2D::ReissnerMindlin::stabilizationParameter ()
 Returns the stabilization parameter. More...
 
RMNorms HArDCore2D::ReissnerMindlin::computeNorms (const Eigen::VectorXd &v) const
 Compute the discrete norms: rotation, displacement, Kirchoff term, and complete energy. More...
 
template<typename outValue , typename Fct >
std::function< outValue(const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::contractPara (const Fct &F) const
 Takes a function dependent on RMParameter and a position x, and returns a function depending only on x (using the parameters of this class) More...
 

Variables

double HArDCore2D::RMParameters::t
 
double HArDCore2D::RMParameters::E
 
double HArDCore2D::RMParameters::nu
 
double HArDCore2D::RMParameters::beta0
 
double HArDCore2D::RMParameters::beta1
 
double HArDCore2D::RMParameters::kappa
 
double HArDCore2D::RMNorms::rotation
 Norm of rotation. More...
 
double HArDCore2D::RMNorms::displacement
 Norm of displacement. More...
 
double HArDCore2D::RMNorms::kirchoff
 Norm of kirchoff term. More...
 
double HArDCore2D::RMNorms::energy
 Total energy. More...
 
static const double HArDCore2D::PI = boost::math::constants::pi<double>()
 
static ReissnerMindlin::SolutionRotationType HArDCore2D::constant_theta
 
static ReissnerMindlin::GradientRotationType HArDCore2D::constant_grad_theta
 
static ReissnerMindlin::SolutionDisplacementType HArDCore2D::constant_u
 
static ReissnerMindlin::ForcingTermType HArDCore2D::constant_f
 
static ReissnerMindlin::SolutionRotationType HArDCore2D::polynomial_theta
 
static ReissnerMindlin::GradientRotationType HArDCore2D::polynomial_grad_theta
 
static ReissnerMindlin::SolutionDisplacementType HArDCore2D::polynomial_u
 
static ReissnerMindlin::ForcingTermType HArDCore2D::polynomial_f
 
static std::function< double(const VectorRd &)> HArDCore2D::an_g = [](const VectorRd &x)->double { return sin(PI*x(0))*sin(PI*x(1)); }
 
static std::function< VectorRd(const VectorRd &)> HArDCore2D::an_GRAD_g
 
static std::function< MatrixRd(const VectorRd &)> HArDCore2D::an_HESS_g
 
static std::function< double(const VectorRd &)> HArDCore2D::an_LAPL_g = [](const VectorRd &x)->double { return -2*pow(PI,2)*sin(PI*x(0))*sin(PI*x(1)); }
 
static std::function< double(const VectorRd &)> HArDCore2D::an_V = [](const VectorRd &x)->double { return x(0) * exp(-x(0))*cos(x(1)); }
 
static std::function< VectorRd(const VectorRd &)> HArDCore2D::an_GRAD_V
 
static std::function< MatrixRd(const VectorRd &)> HArDCore2D::an_HESS_V
 
static std::function< double(const VectorRd &)> HArDCore2D::an_LAPL_V = [](const VectorRd &x)->double { return -2.* exp(-x(0))*cos(x(1)); }
 
static ReissnerMindlin::SolutionDisplacementType HArDCore2D::an_v
 
static ReissnerMindlin::SolutionRotationType HArDCore2D::an_GRAD_v
 
static ReissnerMindlin::GradientRotationType HArDCore2D::an_HESS_v
 
static ReissnerMindlin::SolutionDisplacementType HArDCore2D::an_LAPL_v
 
static ReissnerMindlin::SolutionRotationType HArDCore2D::analytical_theta = an_GRAD_v
 
static ReissnerMindlin::GradientRotationType HArDCore2D::analytical_grad_theta = an_HESS_v
 
static ReissnerMindlin::SolutionDisplacementType HArDCore2D::analytical_u
 
static ReissnerMindlin::ForcingTermType HArDCore2D::analytical_f
 
static ReissnerMindlin::SolutionRotationType HArDCore2D::ukn_theta
 
static ReissnerMindlin::GradientRotationType HArDCore2D::ukn_grad_theta
 
static ReissnerMindlin::SolutionDisplacementType HArDCore2D::ukn_u
 
static ReissnerMindlin::ForcingTermType HArDCore2D::ukn_f
 
static ReissnerMindlin::SolutionRotationType HArDCore2D::kir_theta
 
static ReissnerMindlin::GradientRotationType HArDCore2D::kir_grad_theta
 
static ReissnerMindlin::SolutionDisplacementType HArDCore2D::kir_u
 
static ReissnerMindlin::ForcingTermType HArDCore2D::kir_f
 

Detailed Description

Implementation of the DDR scheme for the Reissner-Mindlin plate problem.

Typedef Documentation

◆ ForcingTermType

typedef std::function<double(const RMParameters &, const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::ForcingTermType

◆ GradientRotationType

typedef std::function<Eigen::Matrix2d(const RMParameters &, const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::GradientRotationType

◆ SolutionDisplacementType

typedef std::function<double(const RMParameters &, const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::SolutionDisplacementType

◆ SolutionRotationType

typedef std::function<Eigen::Vector2d(const RMParameters &, const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::SolutionRotationType

◆ SystemMatrixType

typedef Eigen::SparseMatrix<double> HArDCore2D::ReissnerMindlin::SystemMatrixType

Function Documentation

◆ assembleLinearSystem()

void ReissnerMindlin::assembleLinearSystem ( const ForcingTermType f,
const SolutionRotationType theta,
const GradientRotationType grad_theta,
const SolutionDisplacementType u 
)

Assemble the global system

Parameters
fForcing term
thetaBoundary value
grad_thetaBoundary value
uBoundary value

◆ bdryMatrix()

const SystemMatrixType& HArDCore2D::ReissnerMindlin::bdryMatrix ( ) const
inline

Returns the Matrix for BC.

◆ bdryValues()

const Eigen::VectorXd& HArDCore2D::ReissnerMindlin::bdryValues ( ) const
inline

Returns the boundary values.

◆ computeNorms()

RMNorms ReissnerMindlin::computeNorms ( const Eigen::VectorXd &  v) const

Compute the discrete norms: rotation, displacement, Kirchoff term, and complete energy.

Parameters
vThe vector

◆ contractPara()

template<typename outValue , typename Fct >
std::function<outValue(const Eigen::Vector2d &)> HArDCore2D::ReissnerMindlin::contractPara ( const Fct &  F) const
inline

Takes a function dependent on RMParameter and a position x, and returns a function depending only on x (using the parameters of this class)

◆ dimensionSpace()

size_t HArDCore2D::ReissnerMindlin::dimensionSpace ( ) const
inline

Returns the dimension of the rotation + displacement space (with BC)

◆ exCurl()

const EXCurl& HArDCore2D::ReissnerMindlin::exCurl ( ) const
inline

Returns the space EXCurl.

◆ globalDOFIndices()

std::vector<size_t> HArDCore2D::ReissnerMindlin::globalDOFIndices ( const Cell &  T) const
inline

Create the vector of DOF indices for cell T, which combines the DOFs for the spaces EXcurl and Xgrad.

◆ locUKN()

const std::vector<std::pair<size_t,size_t> >& HArDCore2D::ReissnerMindlin::locUKN ( ) const
inline

Returns the location of the unknowns among the DOFs.

◆ nb_bdryDOFs()

size_t HArDCore2D::ReissnerMindlin::nb_bdryDOFs ( ) const
inline

Returns the nb of DOFs for BC.

◆ para()

const RMParameters& HArDCore2D::ReissnerMindlin::para ( ) const
inline

Returns the parameters.

◆ ReissnerMindlin()

ReissnerMindlin::ReissnerMindlin ( const DDRCore ddrcore,
const RMParameters para,
const BoundaryConditions BC_theta,
const BoundaryConditions BC_u,
bool  use_threads,
std::ostream &  output = std::cout 
)

Constructor.

Parameters
ddrcoreCore for the DDR space sequence
paraPhysical parameters
BC_thetaBoundary conditions for rotation
BC_uBoundary conditions for displacement
use_threadsTrue for parallel execution, false for sequential execution
outputOutput stream to print status messages

◆ RMNorms()

HArDCore2D::RMNorms::RMNorms ( double  norm_rotation,
double  norm_displacement,
double  norm_kirchoff 
)
inline

Constructor.

◆ RMParameters()

HArDCore2D::RMParameters::RMParameters ( const double  thickness,
const double  young_modulus,
const double  poisson_ratio 
)
inline

Constructor.

◆ sizeSystem()

size_t HArDCore2D::ReissnerMindlin::sizeSystem ( ) const
inline

Returns the size of the system without BC.

◆ stabilizationParameter() [1/2]

double& HArDCore2D::ReissnerMindlin::stabilizationParameter ( )
inline

Returns the stabilization parameter.

◆ stabilizationParameter() [2/2]

const double& HArDCore2D::ReissnerMindlin::stabilizationParameter ( ) const
inline

Returns the stabilization parameter.

◆ systemMatrix() [1/2]

SystemMatrixType& HArDCore2D::ReissnerMindlin::systemMatrix ( )
inline

Returns the linear system matrix.

◆ systemMatrix() [2/2]

const SystemMatrixType& HArDCore2D::ReissnerMindlin::systemMatrix ( ) const
inline

Returns the linear system matrix.

◆ systemVector() [1/2]

Eigen::VectorXd& HArDCore2D::ReissnerMindlin::systemVector ( )
inline

Returns the linear system right-hand side vector.

◆ systemVector() [2/2]

const Eigen::VectorXd& HArDCore2D::ReissnerMindlin::systemVector ( ) const
inline

Returns the linear system right-hand side vector.

◆ xGrad()

const XGrad& HArDCore2D::ReissnerMindlin::xGrad ( ) const
inline

Returns the space XGrad.

Variable Documentation

◆ an_g

std::function<double(const VectorRd &)> HArDCore2D::an_g = [](const VectorRd &x)->double { return sin(PI*x(0))*sin(PI*x(1)); }
static

◆ an_GRAD_g

std::function<VectorRd(const VectorRd &)> HArDCore2D::an_GRAD_g
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd( PI*cos(PI*x(0))*sin(PI*x(1)), PI*sin(PI*x(0))*cos(PI*x(1)) );
}
Eigen::Vector2d VectorRd
Definition: basis.hpp:55
static const double PI
Definition: ddr-rmplate.hpp:275
Eigen::Matrix< double, DIMENSION, 1 > VectorRd
Definition: Polytope2D.hpp:19

◆ an_GRAD_V

std::function<VectorRd(const VectorRd &)> HArDCore2D::an_GRAD_V
static
Initial value:
= [](const VectorRd & x) -> VectorRd {
return VectorRd((1.-x(0))*exp(-x(0))*cos(x(1)), -x(0)*exp(-x(0))*sin(x(1)));
}

◆ an_GRAD_v

ReissnerMindlin::SolutionRotationType HArDCore2D::an_GRAD_v
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> VectorRd {
return pow(para.t, 2) * an_GRAD_V(x/para.t) + an_GRAD_g(x);
}
static std::function< VectorRd(const VectorRd &)> an_GRAD_g
Definition: ddr-rmplate.hpp:359
double t
Definition: ddr-rmplate.hpp:55
static std::function< VectorRd(const VectorRd &)> an_GRAD_V
Definition: ddr-rmplate.hpp:378
Structure to store model data.
Definition: ddr-rmplate.hpp:39

◆ an_HESS_g

std::function<MatrixRd(const VectorRd &)> HArDCore2D::an_HESS_g
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd H = MatrixRd::Zero();
H.row(0) << -PI*PI*sin(PI*x(0))*sin(PI*x(1)), PI*PI*cos(PI*x(0))*cos(PI*x(1));
H.row(1) << PI*PI*cos(PI*x(0))*cos(PI*x(1)), -PI*PI*sin(PI*x(0))*sin(PI*x(1));
return H;
}
Eigen::Matrix2d MatrixRd
Definition: basis.hpp:54

◆ an_HESS_V

std::function<MatrixRd(const VectorRd &)> HArDCore2D::an_HESS_V
static
Initial value:
= [](const VectorRd & x) -> MatrixRd {
MatrixRd H = MatrixRd::Zero();
H.row(0) << (x(0)-2.)*exp(-x(0))*cos(x(1)), -(1.-x(0))*exp(-x(0))*sin(x(1));
H.row(1) << -(1.-x(0))*exp(-x(0))*sin(x(1)), -x(0)*exp(-x(0))*cos(x(1));
return H;
}

◆ an_HESS_v

ReissnerMindlin::GradientRotationType HArDCore2D::an_HESS_v
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> MatrixRd {
return para.t * an_HESS_V(x/para.t) + an_HESS_g(x);
}
static std::function< MatrixRd(const VectorRd &)> an_HESS_g
Definition: ddr-rmplate.hpp:364
static std::function< MatrixRd(const VectorRd &)> an_HESS_V
Definition: ddr-rmplate.hpp:383

◆ an_LAPL_g

std::function<double(const VectorRd &)> HArDCore2D::an_LAPL_g = [](const VectorRd &x)->double { return -2*pow(PI,2)*sin(PI*x(0))*sin(PI*x(1)); }
static

◆ an_LAPL_V

std::function<double(const VectorRd &)> HArDCore2D::an_LAPL_V = [](const VectorRd &x)->double { return -2.* exp(-x(0))*cos(x(1)); }
static

◆ an_LAPL_v

ReissnerMindlin::SolutionDisplacementType HArDCore2D::an_LAPL_v
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
return para.t * (an_LAPL_V)(x/para.t) + an_LAPL_g(x);
}
static std::function< double(const VectorRd &)> an_LAPL_g
Definition: ddr-rmplate.hpp:372
static std::function< double(const VectorRd &)> an_LAPL_V
Definition: ddr-rmplate.hpp:391

◆ an_V

std::function<double(const VectorRd &)> HArDCore2D::an_V = [](const VectorRd &x)->double { return x(0) * exp(-x(0))*cos(x(1)); }
static

◆ an_v

Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
return pow(para.t,3) * an_V(x/para.t) + an_g(x);
}
static std::function< double(const VectorRd &)> an_V
Definition: ddr-rmplate.hpp:375
static std::function< double(const VectorRd &)> an_g
Definition: ddr-rmplate.hpp:356

◆ analytical_f

ReissnerMindlin::ForcingTermType HArDCore2D::analytical_f
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
return (para.beta0+para.beta1) * 4 * pow(PI,4)*sin(PI*x(0))*sin(PI*x(1));
}
double beta1
Definition: ddr-rmplate.hpp:59
double beta0
Definition: ddr-rmplate.hpp:58

◆ analytical_grad_theta

ReissnerMindlin::GradientRotationType HArDCore2D::analytical_grad_theta = an_HESS_v
static

◆ analytical_theta

ReissnerMindlin::SolutionRotationType HArDCore2D::analytical_theta = an_GRAD_v
static

◆ analytical_u

ReissnerMindlin::SolutionDisplacementType HArDCore2D::analytical_u
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
return an_v(para, x) - pow(para.t,2) * ( (para.beta0+para.beta1)/para.kappa ) * an_LAPL_v(para, x);
}
static ReissnerMindlin::SolutionDisplacementType an_LAPL_v
Definition: ddr-rmplate.hpp:410
double kappa
Definition: ddr-rmplate.hpp:60
static ReissnerMindlin::SolutionDisplacementType an_v
Definition: ddr-rmplate.hpp:395

◆ beta0

double HArDCore2D::RMParameters::beta0

◆ beta1

double HArDCore2D::RMParameters::beta1

◆ constant_f

ReissnerMindlin::ForcingTermType HArDCore2D::constant_f
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
return 0.;
}

◆ constant_grad_theta

ReissnerMindlin::GradientRotationType HArDCore2D::constant_grad_theta
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> Eigen::Matrix2d {
return Eigen::Matrix2d::Zero();
}

◆ constant_theta

ReissnerMindlin::SolutionRotationType HArDCore2D::constant_theta
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> VectorRd {
return VectorRd(1., 2.);
}

◆ constant_u

ReissnerMindlin::SolutionDisplacementType HArDCore2D::constant_u
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
return 1.;
}

◆ displacement

double HArDCore2D::RMNorms::displacement

Norm of displacement.

◆ E

double HArDCore2D::RMParameters::E

◆ energy

double HArDCore2D::RMNorms::energy

Total energy.

◆ kappa

double HArDCore2D::RMParameters::kappa

◆ kir_f

ReissnerMindlin::ForcingTermType HArDCore2D::kir_f
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
return sin(PI*x(0))*sin(PI*x(1));
}

◆ kir_grad_theta

ReissnerMindlin::GradientRotationType HArDCore2D::kir_grad_theta
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> Eigen::Matrix2d {
double coef = 12.*(1-pow(para.nu,2))/ (para.E * 4.*pow(PI, 4));
Eigen::Matrix2d H = Eigen::Matrix2d::Zero();
H.row(0) << -PI*PI*sin(PI*x(0))*sin(PI*x(1)), PI*PI*cos(PI*x(0))*cos(PI*x(1));
H.row(1) << PI*PI*cos(PI*x(0))*cos(PI*x(1)), -PI*PI*sin(PI*x(0))*sin(PI*x(1));
return coef * H;
}
double E
Definition: ddr-rmplate.hpp:56
double nu
Definition: ddr-rmplate.hpp:57

◆ kir_theta

ReissnerMindlin::SolutionRotationType HArDCore2D::kir_theta
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> VectorRd {
double coef = 12.*(1-pow(para.nu,2))/ (para.E * 4.*pow(PI, 4));
return coef * PI * VectorRd(cos(PI*x(0))*sin(PI*x(1)), sin(PI*x(0))*cos(PI*x(1)));
}

◆ kir_u

Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
double coef = 12.*(1-pow(para.nu,2))/ (para.E * 4.*pow(PI, 4));
return coef * sin(PI*x(0))*sin(PI*x(1));
}

◆ kirchoff

double HArDCore2D::RMNorms::kirchoff

Norm of kirchoff term.

◆ nu

double HArDCore2D::RMParameters::nu

◆ PI

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

◆ polynomial_f

ReissnerMindlin::ForcingTermType HArDCore2D::polynomial_f
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
double val = 0;
val += 12*x(1)*(x(1)-1)*(5*x(0)*x(0)-5*x(0)+1) *
( 2*x(1)*x(1)*(x(1)-1)*(x(1)-1) + x(0)*(x(0)-1)*(5*x(1)*x(1)-5*x(1)+1) );
val += 12*x(0)*(x(0)-1)*(5*x(1)*x(1)-5*x(1)+1) *
( 2*x(0)*x(0)*(x(0)-1)*(x(0)-1) + x(1)*(x(1)-1)*(5*x(0)*x(0)-5*x(0)+1) );
return val * para.E / (12* (1.-para.nu*para.nu) );
}

◆ polynomial_grad_theta

ReissnerMindlin::GradientRotationType HArDCore2D::polynomial_grad_theta
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> Eigen::Matrix2d {
Eigen::Matrix2d G = Eigen::Matrix2d::Zero();
G(0,0) = 2.*pow(x(0),2)*pow(x(1),3)*pow(x(0)-1.0,2)*pow(x(1)-1.0,3)+2.*x(0)*pow(x(1),3)*(2.*x(0)-1.)*pow(x(0)-1.,2)*pow(x(1)-1.0,3)+pow(x(0),2)*pow(x(1),3)*(2.*x(0)-1.)*(2.*x(0)-2.)*pow(x(1)-1.,3);
G(0,1) = 3.*pow(x(0)*x(1),2)*(2.*x(0)-1.)*pow(x(0)-1.,2)*pow(x(1)-1.,3)+3.*pow(x(0),2)*pow(x(1),3)*(2.*x(0)-1.)*pow(x(0)-1.,2)*pow(x(1)-1.,2);
G(1,0) = 2.*pow(x(1),2)*pow(x(0),3)*pow(x(1)-1.0,2)*pow(x(0)-1.0,3)+2.*x(1)*pow(x(0),3)*(2.*x(1)-1.)*pow(x(1)-1.,2)*pow(x(0)-1.0,3)+pow(x(1),2)*pow(x(0),3)*(2.*x(1)-1.)*(2.*x(1)-2.)*pow(x(0)-1.,3);
G(1,1) = 3.*pow(x(1)*x(0),2)*(2.*x(1)-1.)*pow(x(1)-1.,2)*pow(x(0)-1.,3)+3.*pow(x(1),2)*pow(x(0),3)*(2.*x(1)-1.)*pow(x(1)-1.,2)*pow(x(0)-1.,2);
return G;
}

◆ polynomial_theta

ReissnerMindlin::SolutionRotationType HArDCore2D::polynomial_theta
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> VectorRd {
return VectorRd(
pow(x(1),3)*pow(x(1)-1,3)*pow(x(0),2)*pow(x(0)-1,2)*(2*x(0)-1),
pow(x(0),3)*pow(x(0)-1,3)*pow(x(1),2)*pow(x(1)-1,2)*(2*x(1)-1)
);
}

◆ polynomial_u

ReissnerMindlin::SolutionDisplacementType HArDCore2D::polynomial_u
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
double val = 0;
val += (1./3.)*pow(x(0),3)*pow(x(0)-1,3)*pow(x(1),3)*pow(x(1)-1,3);
val -= ( 2*pow(para.t,2) / (5*(1.-para.nu)) ) *
( pow(x(1),3)*pow(x(1)-1,3)*x(0)*(x(0)-1)*(5*x(0)*x(0)-5*x(0)+1)
+ pow(x(0),3)*pow(x(0)-1,3)*x(1)*(x(1)-1)*(5*x(1)*x(1)-5*x(1)+1) );
return val;
}

◆ rotation

double HArDCore2D::RMNorms::rotation

Norm of rotation.

◆ t

double HArDCore2D::RMParameters::t

◆ ukn_f

ReissnerMindlin::ForcingTermType HArDCore2D::ukn_f
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
double val = 1.;
if (x(0)<.5){
val = -1.;
}
return val;
}

◆ ukn_grad_theta

ReissnerMindlin::GradientRotationType HArDCore2D::ukn_grad_theta
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> Eigen::Matrix2d {
return Eigen::Matrix2d::Zero();
}

◆ ukn_theta

ReissnerMindlin::SolutionRotationType HArDCore2D::ukn_theta
static
Initial value:
= [](const RMParameters & para, const VectorRd & x) -> VectorRd {
return VectorRd(0., 0.);
}

◆ ukn_u

Initial value:
= [](const RMParameters & para, const VectorRd & x) -> double {
return x(0);
}