Struct ae108::elements::ElementBase
template <class Derived_, class SizeType_, class ValueType_, class RealType_, SizeType_ Size_, SizeType_ Dimension_, SizeType_ DegreesOfFreedom_>
ClassList > ae108 > elements > ElementBase
Public Types
Type | Name |
---|---|
typedef real_type | 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 tensor::Tensor< value_type, Size_, DegreesOfFreedom_ > | NodalDisplacements The displacements per node. |
typedef Eigen::Matrix< value_type, Size_ *DegreesOfFreedom_, Size_ *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 real_type | Time |
typedef RealType_ | real_type |
typedef SizeType_ | size_type |
typedef ValueType_ | value_type |
Public Functions
Type | Name |
---|---|
Energy | computeEnergy (const NodalDisplacements & displacements, const Time & time) const Computes the energy for the given displacements. |
Forces | computeForces (const NodalDisplacements & displacements, const Time & time) const Computes the forces for the given displacements. |
StiffnessMatrix | computeStiffnessMatrix (const NodalDisplacements & displacements, const Time & time) const Computes the stiffness matrix for the given displacements. |
Public Static Functions
Type | Name |
---|---|
constexpr size_type | degrees_of_freedom () noexcept Number of degrees of freedom. |
constexpr size_type | dimension () noexcept The dimension of physical space. |
constexpr size_type | size () noexcept Number of element nodes / shape functions. |
Public Types Documentation
typedef Energy
using ae108::elements::ElementBase< Derived_, SizeType_, ValueType_, RealType_, Size_, Dimension_, DegreesOfFreedom_ >::Energy = real_type;
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::ElementBase< Derived_, SizeType_, ValueType_, RealType_, Size_, Dimension_, DegreesOfFreedom_ >::Forces = NodalDisplacements;
typedef NodalDisplacements
The displacements per node.
using ae108::elements::ElementBase< Derived_, SizeType_, ValueType_, RealType_, Size_, Dimension_, DegreesOfFreedom_ >::NodalDisplacements = tensor::Tensor<value_type, Size_, DegreesOfFreedom_>;
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::ElementBase< Derived_, SizeType_, ValueType_, RealType_, Size_, Dimension_, DegreesOfFreedom_ >::StiffnessMatrix = Eigen::Matrix<value_type, Size_ * DegreesOfFreedom_, Size_ * DegreesOfFreedom_, Eigen::RowMajor>;
typedef Time
using ae108::elements::ElementBase< Derived_, SizeType_, ValueType_, RealType_, Size_, Dimension_, DegreesOfFreedom_ >::Time = real_type;
typedef real_type
using ae108::elements::ElementBase< Derived_, SizeType_, ValueType_, RealType_, Size_, Dimension_, DegreesOfFreedom_ >::real_type = RealType_;
typedef size_type
using ae108::elements::ElementBase< Derived_, SizeType_, ValueType_, RealType_, Size_, Dimension_, DegreesOfFreedom_ >::size_type = SizeType_;
typedef value_type
using ae108::elements::ElementBase< Derived_, SizeType_, ValueType_, RealType_, Size_, Dimension_, DegreesOfFreedom_ >::value_type = ValueType_;
Public Functions Documentation
function computeEnergy
Computes the energy for the given displacements.
inline Energy ae108::elements::ElementBase::computeEnergy (
const NodalDisplacements & displacements,
const Time & time
) const
Remark:
Delegates to the corresponding free function.
function computeForces
Computes the forces for the given displacements.
inline Forces ae108::elements::ElementBase::computeForces (
const NodalDisplacements & displacements,
const Time & time
) const
Remark:
Delegates to the corresponding free function.
function computeStiffnessMatrix
Computes the stiffness matrix for the given displacements.
inline StiffnessMatrix ae108::elements::ElementBase::computeStiffnessMatrix (
const NodalDisplacements & displacements,
const Time & time
) const
Remark:
Delegates to the corresponding free function.
Public Static Functions Documentation
function degrees_of_freedom
Number of degrees of freedom.
static inline constexpr size_type ae108::elements::ElementBase::degrees_of_freedom () noexcept
function dimension
The dimension of physical space.
static inline constexpr size_type ae108::elements::ElementBase::dimension () noexcept
function size
Number of element nodes / shape functions.
static inline constexpr size_type ae108::elements::ElementBase::size () noexcept
The documentation for this class was generated from the following file elements/src/include/ae108/elements/ElementBase.h