Struct ae108::elements::ForceElement
template <std::size_t Dimension_, class ValueType_, class RealType_>
ClassList > ae108 > elements > ForceElement
A single-vertex element that applies a force at that vertex.
#include <ForceElement.h>
Public Types
| Type | Name |
|---|---|
| typedef ForceElement< Dimension_, ValueType_, RealType_ > | Element |
| typedef RealType_ | Energy |
| typedef tensor::Tensor< Value, Dimension > | Force |
| 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, NumVerts, DegreesOfFreedom > | NodalDisplacements The displacements per node. |
| typedef RealType_ | Real |
| typedef std::size_t | 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 RealType_ | Time |
| typedef ValueType_ | Value |
Public Attributes
| Type | Name |
|---|---|
| Force | force |
Public Static Attributes
| Type | Name |
|---|---|
| const Size | DegreesOfFreedom = Dimension\_Number of degrees of freedom. |
| const Size | Dimension = Dimension\_ |
| const Size | NumVerts = 1Number of element nodes / shape functions. |
Public Functions
| Type | Name |
|---|---|
| Energy | computeEnergy (const NodalDisplacements & u, const Time &) noexcept const Computes the energy for the given displacements. |
| Forces | computeForces (const NodalDisplacements &, const Time &) noexcept const Computes the forces for the given displacements. |
| StiffnessMatrix | computeStiffnessMatrix (const NodalDisplacements &, const Time &) noexcept const Computes the stiffness matrix for the given displacements. |
Public Types Documentation
typedef Element
using ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::Element = ForceElement<Dimension_, ValueType_, RealType_>;
typedef Energy
using ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::Energy = RealType_;
typedef Force
using ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::Force = tensor::Tensor<Value, Dimension>;
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::ForceElement< Dimension_, ValueType_, RealType_ >::Forces = NodalDisplacements;
typedef NodalDisplacements
The displacements per node.
using ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::NodalDisplacements = tensor::Tensor<Value, NumVerts, DegreesOfFreedom>;
typedef Real
using ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::Real = RealType_;
typedef Size
using ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::Size = std::size_t;
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::ForceElement< Dimension_, ValueType_, RealType_ >::StiffnessMatrix = Eigen::Matrix<Value, NumVerts * DegreesOfFreedom, NumVerts * DegreesOfFreedom, Eigen::RowMajor>;
typedef Time
using ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::Time = RealType_;
typedef Value
using ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::Value = ValueType_;
Public Attributes Documentation
variable force
Force ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::force;
Public Static Attributes Documentation
variable DegreesOfFreedom
Number of degrees of freedom.
const Size ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::DegreesOfFreedom;
variable Dimension
const Size ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::Dimension;
variable NumVerts
Number of element nodes / shape functions.
const Size ae108::elements::ForceElement< Dimension_, ValueType_, RealType_ >::NumVerts;
Public Functions Documentation
function computeEnergy
Computes the energy for the given displacements.
inline Energy ae108::elements::ForceElement::computeEnergy (
const NodalDisplacements & u,
const Time &
) noexcept const
function computeForces
Computes the forces for the given displacements.
inline Forces ae108::elements::ForceElement::computeForces (
const NodalDisplacements &,
const Time &
) noexcept const
function computeStiffnessMatrix
Computes the stiffness matrix for the given displacements.
inline StiffnessMatrix ae108::elements::ForceElement::computeStiffnessMatrix (
const NodalDisplacements &,
const Time &
) noexcept const
The documentation for this class was generated from the following file elements/src/include/ae108/elements/ForceElement.h