Struct ae108::elements::CoreElement
template <class MaterialModel_, class Integrator_, class ValueType_, class RealType_>
ClassList > ae108 > elements > CoreElement
An element that computes energy, forces, and stiffness matrix using a material model and an integrator.
#include <CoreElement.h>
Public Types
| Type | Name |
|---|---|
| typedef CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ > | Element |
| typedef Real | Energy |
| typedef NodalDisplacements | Forces The forces equal to $d_{ij} E$. Here, d_ij refers to the derivative with respect to jth degree of freedom of the ith node. |
| typedef Integrator_ | Integrator |
| typedef MaterialModel_ | MaterialModel |
| typedef tensor::Tensor< Value, NumVerts, DegreesOfFreedom > | NodalDisplacements The displacements per node. |
| typedef RealType_ | Real |
| typedef Eigen::Map< Eigen::Matrix< Value, NumVerts, NumVerts, Eigen::RowMajor >, 0, Eigen::Stride< Eigen::Index{NumVerts *DegreesOfFreedom *DegreesOfFreedom}, Eigen::Index{DegreesOfFreedom}> > | ResultSlice |
| typedef typename Integrator_::size_type | Size |
| typedef Eigen::Matrix< Value, NumVerts *DegreesOfFreedom, NumVerts *DegreesOfFreedom, Eigen::RowMajor > | StiffnessMatrix The stiffness matrix equal to $d_{ij} d_{kl} E$. Here, d_ij refers to the derivative with respect to jth degree of freedom of the ith node. |
| typedef std::conditional_t< DegreesOfFreedom==0||sizeof(typename MaterialModel::TangentMatrix)==DegreesOfFreedom *MaterialModel::dimension() *DegreesOfFreedom *MaterialModel::dimension() *sizeof(typename MaterialModel::value_type), Eigen::Map< const Eigen::Matrix< typename MaterialModel::value_type, MaterialModel::dimension(), MaterialModel::dimension(), Eigen::RowMajor >, 0, Eigen::Stride< Eigen::Index{MaterialModel::dimension() *Dimension}, Eigen::Index{1}> >, void > | TangentSlice |
| typedef Real | Time |
| typedef ValueType_ | Value |
Public Attributes
| Type | Name |
|---|---|
| Integrator | integrator The instance of the class used to integrate entities (e.g. energy). |
| MaterialModel | model The instance of the model used to compute entities at integration points (e.g. energy). |
Public Static Attributes
| Type | Name |
|---|---|
| const Size | DegreesOfFreedom = /* multi line expression */Number of degrees of freedom. |
| const Size | Dimension = Integrator\_::dimension() |
| const Size | NumVerts = Integrator\_::size()Number of element nodes / shape functions. |
Public Functions
| Type | Name |
|---|---|
| StiffnessMatrix | computeConsistentMassMatrix () noexcept const |
| Energy | computeEnergy (const NodalDisplacements & u, const Time & time) noexcept const Computes the energy for the given displacements. |
| Forces | computeForces (const NodalDisplacements & u, const Time & time) noexcept const Computes the forces for the given displacements. |
| StiffnessMatrix | computeLumpedMassMatrix () noexcept const |
| StiffnessMatrix | computeStiffnessMatrix (const NodalDisplacements & u, const Time & time) noexcept const Computes the stiffness matrix for the given displacements. |
Public Types Documentation
typedef Element
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Element = CoreElement<MaterialModel_, Integrator_, ValueType_, RealType_>;
typedef Energy
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Energy = Real;
typedef Forces
The forces equal to $d_{ij} E$. Here, d_ij refers to the derivative with respect to jth degree of freedom of the ith node.
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Forces = NodalDisplacements;
typedef Integrator
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Integrator = Integrator_;
typedef MaterialModel
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::MaterialModel = MaterialModel_;
typedef NodalDisplacements
The displacements per node.
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::NodalDisplacements = tensor::Tensor<Value, NumVerts, DegreesOfFreedom>;
typedef Real
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Real = RealType_;
typedef ResultSlice
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::ResultSlice = Eigen::Map<Eigen::Matrix<Value, NumVerts, NumVerts, Eigen::RowMajor>, 0, Eigen::Stride<Eigen::Index{NumVerts * DegreesOfFreedom * DegreesOfFreedom}, Eigen::Index{DegreesOfFreedom}> >;
typedef Size
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Size = typename Integrator_::size_type;
typedef StiffnessMatrix
The stiffness matrix equal to $d_{ij} d_{kl} E$. Here, d_ij refers to the derivative with respect to jth degree of freedom of the ith node.
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::StiffnessMatrix = Eigen::Matrix<Value, NumVerts * DegreesOfFreedom, NumVerts * DegreesOfFreedom, Eigen::RowMajor>;
typedef TangentSlice
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::TangentSlice = std::conditional_t< DegreesOfFreedom == 0 || sizeof(typename MaterialModel::TangentMatrix) == DegreesOfFreedom * MaterialModel::dimension() * DegreesOfFreedom * MaterialModel::dimension() * sizeof(typename MaterialModel::value_type), Eigen::Map< const Eigen::Matrix<typename MaterialModel::value_type, MaterialModel::dimension(), MaterialModel::dimension(), Eigen::RowMajor>, 0, Eigen::Stride<Eigen::Index{MaterialModel::dimension() * Dimension}, Eigen::Index{1}> >, void>;
typedef Time
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Time = Real;
typedef Value
using ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Value = ValueType_;
Public Attributes Documentation
variable integrator
The instance of the class used to integrate entities (e.g. energy).
Integrator ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::integrator;
variable model
The instance of the model used to compute entities at integration points (e.g. energy).
MaterialModel ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::model;
Public Static Attributes Documentation
variable DegreesOfFreedom
Number of degrees of freedom.
const Size ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::DegreesOfFreedom;
variable Dimension
const Size ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::Dimension;
variable NumVerts
Number of element nodes / shape functions.
const Size ae108::elements::CoreElement< MaterialModel_, Integrator_, ValueType_, RealType_ >::NumVerts;
Public Functions Documentation
function computeConsistentMassMatrix
inline StiffnessMatrix ae108::elements::CoreElement::computeConsistentMassMatrix () noexcept const
function computeEnergy
Computes the energy for the given displacements.
inline Energy ae108::elements::CoreElement::computeEnergy (
const NodalDisplacements & u,
const Time & time
) noexcept const
function computeForces
Computes the forces for the given displacements.
inline Forces ae108::elements::CoreElement::computeForces (
const NodalDisplacements & u,
const Time & time
) noexcept const
function computeLumpedMassMatrix
inline StiffnessMatrix ae108::elements::CoreElement::computeLumpedMassMatrix () noexcept const
function computeStiffnessMatrix
Computes the stiffness matrix for the given displacements.
inline StiffnessMatrix ae108::elements::CoreElement::computeStiffnessMatrix (
const NodalDisplacements & u,
const Time & time
) noexcept const
The documentation for this class was generated from the following file elements/src/include/ae108/elements/CoreElement.h