HArD::Core2D
Hybrid Arbitrary Degree::Core 2D - Library to implement 2D schemes with edge and cell polynomials as unknowns
|
Implementation of the DDR scheme for the Kirchoff-Love plate problem. More...
Classes | |
struct | HArDCore2D::KirchhoffLove |
Assemble a Kirchhoff-Love plate problem. More... | |
Typedefs | |
typedef Eigen::SparseMatrix< double > | HArDCore2D::KirchhoffLove::SystemMatrixType |
typedef std::function< double(const Eigen::Vector2d &)> | HArDCore2D::KirchhoffLove::ForcingTermType |
typedef std::function< double(const Eigen::Vector2d &)> | HArDCore2D::KirchhoffLove::DeflectionType |
typedef std::function< Eigen::Vector2d(const Eigen::Vector2d &)> | HArDCore2D::KirchhoffLove::GradientDeflectionType |
typedef std::function< Eigen::Matrix2d(const Eigen::Vector2d &)> | HArDCore2D::KirchhoffLove::MomentTensorType |
typedef std::function< double(const Eigen::Vector2d &, const Edge &)> | HArDCore2D::KirchhoffLove::MomentTensorEdgeDerivativeType |
typedef XDivDiv::ConstitutiveLawType | HArDCore2D::KirchhoffLove::ConstitutiveLawType |
Functions | |
HArDCore2D::KirchhoffLove::KirchhoffLove (const PlatesCore &platescore, const ConstitutiveLawType &law, bool use_threads, std::ostream &output=std::cout) | |
Constructor. | |
size_t | HArDCore2D::KirchhoffLove::dimensionSpace () const |
Returns the dimension of the moment + deflection space. | |
size_t | HArDCore2D::KirchhoffLove::nbSCDOFs () const |
Returns the number of statically condensed DOFs (here, the cell moments DOFs) | |
size_t | HArDCore2D::KirchhoffLove::sizeSystem () const |
Returns the size of the statically condensed system. | |
const XDivDiv & | HArDCore2D::KirchhoffLove::xDivDiv () const |
Returns the space XDivDiv. | |
const GlobalDOFSpace | HArDCore2D::KirchhoffLove::polykm2Th () const |
Returns the space \(\mathbb{P}^{k-2}(\mathcal{T}_h)\). | |
const SystemMatrixType & | HArDCore2D::KirchhoffLove::systemMatrix () const |
Returns the linear system matrix. | |
SystemMatrixType & | HArDCore2D::KirchhoffLove::systemMatrix () |
Returns the linear system matrix. | |
const Eigen::VectorXd & | HArDCore2D::KirchhoffLove::systemVector () const |
Returns the linear system right-hand side vector. | |
Eigen::VectorXd & | HArDCore2D::KirchhoffLove::systemVector () |
Returns the linear system right-hand side vector. | |
const SystemMatrixType & | HArDCore2D::KirchhoffLove::scMatrix () const |
Returns the static condensation recovery operator. | |
Eigen::VectorXd & | HArDCore2D::KirchhoffLove::scVector () |
Returns the static condensation rhs. | |
const double & | HArDCore2D::KirchhoffLove::stabilizationParameter () const |
Returns the stabilization parameter. | |
double & | HArDCore2D::KirchhoffLove::stabilizationParameter () |
Returns the stabilization parameter. | |
void | HArDCore2D::KirchhoffLove::assembleLinearSystem (const ForcingTermType &f, const DeflectionType &u, const GradientDeflectionType &grad_u) |
Assemble the global system. | |
Eigen::VectorXd | HArDCore2D::KirchhoffLove::interpolateDeflection (const DeflectionType &u, int deg_quad=-1) const |
Interpolate deflection. | |
double | HArDCore2D::KirchhoffLove::computeNorm (const Eigen::VectorXd &v) const |
Compute the discrete norm. | |
Implementation of the DDR scheme for the Kirchoff-Love plate problem.
typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::KirchhoffLove::DeflectionType |
typedef std::function<double(const Eigen::Vector2d &)> HArDCore2D::KirchhoffLove::ForcingTermType |
typedef std::function<Eigen::Vector2d(const Eigen::Vector2d &)> HArDCore2D::KirchhoffLove::GradientDeflectionType |
typedef std::function<double(const Eigen::Vector2d &, const Edge &)> HArDCore2D::KirchhoffLove::MomentTensorEdgeDerivativeType |
typedef std::function<Eigen::Matrix2d(const Eigen::Vector2d &)> HArDCore2D::KirchhoffLove::MomentTensorType |
typedef Eigen::SparseMatrix<double> HArDCore2D::KirchhoffLove::SystemMatrixType |
void KirchhoffLove::assembleLinearSystem | ( | const ForcingTermType & | f, |
const DeflectionType & | u, | ||
const GradientDeflectionType & | grad_u | ||
) |
Assemble the global system.
f | Forcing term |
u | Deflection (for BC) |
grad_u | Gradient of deflection (for BC) |
double KirchhoffLove::computeNorm | ( | const Eigen::VectorXd & | v | ) | const |
Compute the discrete norm.
v | The vector |
|
inline |
Returns the dimension of the moment + deflection space.
Eigen::VectorXd KirchhoffLove::interpolateDeflection | ( | const DeflectionType & | u, |
int | deg_quad = -1 |
||
) | const |
Interpolate deflection.
u | The function to interpolate |
deg_quad | The degree of the quadrature rules (equal to \(2(k+1)\) with \(k\) equal to the degree of the complex by default) |
KirchhoffLove::KirchhoffLove | ( | const PlatesCore & | platescore, |
const ConstitutiveLawType & | law, | ||
bool | use_threads, | ||
std::ostream & | output = std::cout |
||
) |
Constructor.
platescore | Core for the DDR space sequence |
law | Constitutive law: law(sigma)+hess u = 0 |
use_threads | True for parallel execution, false for sequential execution |
output | Output stream to print status messages |
|
inline |
Returns the number of statically condensed DOFs (here, the cell moments DOFs)
|
inline |
Returns the space \(\mathbb{P}^{k-2}(\mathcal{T}_h)\).
|
inline |
Returns the static condensation recovery operator.
|
inline |
Returns the static condensation rhs.
|
inline |
Returns the size of the statically condensed system.
|
inline |
Returns the stabilization parameter.
|
inline |
Returns the stabilization parameter.
|
inline |
Returns the linear system matrix.
|
inline |
Returns the linear system matrix.
|
inline |
Returns the linear system right-hand side vector.
|
inline |
Returns the linear system right-hand side vector.
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |
|
static |