#include <ikarus/finiteelements/mechanics/membranestrains.hh>
|
template<typename Geometry > |
auto | value (const Dune::FieldVector< double, 2 > &gpPos, const Geometry &geo, const auto &uFunction) const -> Eigen::Vector3< typename std::remove_cvref_t< decltype(uFunction)>::ctype > |
| Compute the strain vector at a given integration point. More...
|
|
template<typename Geometry , typename ScalarType > |
auto | derivative (const Dune::FieldVector< double, 2 > &gpPos, const Eigen::Matrix< ScalarType, 2, 3 > &jcur, const auto &dNAtGp, const Geometry &geo, const auto &uFunction, const auto &localBasis, const int node) const |
| Compute the strain-displacement matrix for a given node and integration point. More...
|
|
template<typename Geometry , typename ScalarType > |
auto | secondDerivative (const Dune::FieldVector< double, 2 > &gpPos, const auto &dNAtGp, const Geometry &geo, const auto &uFunction, const auto &localBasis, const Eigen::Vector3< ScalarType > &S, int I, int J) const |
| Compute the second derivative of the membrane strain for a given node pair and integration point. More...
|
|
◆ derivative()
template<typename Geometry , typename ScalarType >
auto Ikarus::DefaultMembraneStrain::derivative |
( |
const Dune::FieldVector< double, 2 > & |
gpPos, |
|
|
const Eigen::Matrix< ScalarType, 2, 3 > & |
jcur, |
|
|
const auto & |
dNAtGp, |
|
|
const Geometry & |
geo, |
|
|
const auto & |
uFunction, |
|
|
const auto & |
localBasis, |
|
|
const int |
node |
|
) |
| const |
|
inline |
- Parameters
-
gpPos | The position of the integration point. |
jcur | The Jacobian matrix. |
dNAtGp | The derivative of the shape functions at the integration point. |
geo | The geometry object of the finite element. |
uFunction | The function representing the displacement field. |
localBasis | The local basis object. |
node | The FE node index. |
- Template Parameters
-
Geometry | The type of the geometry object. |
ScalarType | The scalar type for the matrix elements. |
- Returns
- The strain-displacement matrix for the given node and integration point.
◆ secondDerivative()
template<typename Geometry , typename ScalarType >
auto Ikarus::DefaultMembraneStrain::secondDerivative |
( |
const Dune::FieldVector< double, 2 > & |
gpPos, |
|
|
const auto & |
dNAtGp, |
|
|
const Geometry & |
geo, |
|
|
const auto & |
uFunction, |
|
|
const auto & |
localBasis, |
|
|
const Eigen::Vector3< ScalarType > & |
S, |
|
|
int |
I, |
|
|
int |
J |
|
) |
| const |
|
inline |
This function computes the geometric tangent stiffness for a given node pair at a given integration point.
- Parameters
-
gpPos | The position of the integration point. |
dNAtGp | The derivative of the shape functions at the integration point. |
geo | The geometry object. |
uFunction | The function representing the displacement field. |
localBasis | The local basis object. |
S | The strain vector. |
I | The index of the first node. |
J | The index of the second node. |
- Template Parameters
-
Geometry | The type of the geometry object. |
ScalarType | The scalar type for the matrix elements. |
- Returns
- The second derivative of the membrane strain.
◆ value()
template<typename Geometry >
auto Ikarus::DefaultMembraneStrain::value |
( |
const Dune::FieldVector< double, 2 > & |
gpPos, |
|
|
const Geometry & |
geo, |
|
|
const auto & |
uFunction |
|
) |
| const -> Eigen::Vector3<typename std::remove_cvref_t<decltype(uFunction)>::ctype> |
|
inline |
- Parameters
-
gpPos | The position of the integration point. |
geo | The geometry object providing the transposed Jacobian. |
uFunction | The function representing the displacement field. |
- Template Parameters
-
Geometry | The type of the geometry object. |
- Returns
- The strain vector at the given integration point.
The documentation for this struct was generated from the following file: