#include <ikarus/finiteelements/mechanics/assumedstress/asfunctions/pk2stress.hh>
|
template<typename GEO , typename AST > |
static auto | value (const GEO &geo, const Dune::FieldVector< double, GEO::mydimension > &gpPos, const AST &asFunction, const auto &beta) |
| Compute the stress vector at a given integration point or its index. More...
|
|
template<typename GEO , typename AST > |
static auto | firstDerivative (const GEO &geo, const auto &uFunction, const auto &localBasis, const auto &gpIndex, const Dune::FieldVector< double, GEO::mydimension > &gpPos, const AST &asFunction, const auto &beta, const int node=sNaN) |
| Compute the first derivative of the PK2 stress w.r.t beta for a given node and integration point. More...
|
|
template<typename GEO , typename ST , typename AST > |
static auto | secondDerivative (const GEO &geo, const auto &uFunction, const auto &localBasis, const auto &gpIndex, const Dune::FieldVector< double, GEO::mydimension > &gpPos, const Eigen::Vector< ST, GEO::mydimension *(GEO::mydimension+1)/2 > &S, const AST &asFunction, const auto &beta, const int I=sNaN, const int J=sNaN) |
| Compute the second derivative of the PK2 stress w.r.t beta for a given node and integration point. More...
|
|
static constexpr auto | name () |
| The name of the assumed stress type. More...
|
|
◆ firstDerivative()
template<typename GEO , typename AST >
static auto Ikarus::PS::PK2Stress::firstDerivative |
( |
const GEO & |
geo, |
|
|
const auto & |
uFunction, |
|
|
const auto & |
localBasis, |
|
|
const auto & |
gpIndex, |
|
|
const Dune::FieldVector< double, GEO::mydimension > & |
gpPos, |
|
|
const AST & |
asFunction, |
|
|
const auto & |
beta, |
|
|
const int |
node = sNaN |
|
) |
| |
|
inlinestatic |
- Parameters
-
geo | The geometry object of the finite element. |
uFunction | The function representing the displacement field. |
localBasis | The local basis object. |
gpPos | The position of the integration point. |
gpIndex | The index of the integration point. |
node | The FE node index (defaults to sNaN). |
asFunction | The AssumedStress function. |
beta | The coefficients of the PS function. |
- Template Parameters
-
GEO | The type of the geometry object. |
AST | The type of the Assumed Stress function function. |
- Returns
- The first derivative of the PK2 stress w.r.t beta.
◆ name()
static constexpr auto Ikarus::PS::PK2Stress::name |
( |
| ) |
|
|
inlinestaticconstexpr |
◆ secondDerivative()
template<typename GEO , typename ST , typename AST >
static auto Ikarus::PS::PK2Stress::secondDerivative |
( |
const GEO & |
geo, |
|
|
const auto & |
uFunction, |
|
|
const auto & |
localBasis, |
|
|
const auto & |
gpIndex, |
|
|
const Dune::FieldVector< double, GEO::mydimension > & |
gpPos, |
|
|
const Eigen::Vector< ST, GEO::mydimension *(GEO::mydimension+1)/2 > & |
S, |
|
|
const AST & |
asFunction, |
|
|
const auto & |
beta, |
|
|
const int |
I = sNaN , |
|
|
const int |
J = sNaN |
|
) |
| |
|
inlinestatic |
- Parameters
-
geo | The geometry object of the finite element. |
uFunction | The function representing the displacement field. |
localBasis | The local basis object. |
gpPos | The position of the integration point. |
gpIndex | The index of the integration point. |
I | The FE node index I (defaults to sNaN). |
J | The FE node index J (defaults to sNaN). |
S | The PK2 stress (in Voigt notation). |
asFunction | The AssumedStress function. |
beta | The coefficients of the PS function. |
- Template Parameters
-
GEO | The type of the geometry object. |
ST | The underlying scalar type. |
AST | The type of the Assumed Stress function. |
- Returns
- The second derivative of the PK2 stress w.r.t beta
◆ value()
template<typename GEO , typename AST >
static auto Ikarus::PS::PK2Stress::value |
( |
const GEO & |
geo, |
|
|
const Dune::FieldVector< double, GEO::mydimension > & |
gpPos, |
|
|
const AST & |
asFunction, |
|
|
const auto & |
beta |
|
) |
| |
|
inlinestatic |
- Parameters
-
geo | The geometry object providing the transposed Jacobian. |
gpPos | The position of the integration point. |
asFunction | The AssumedStress function. |
beta | The coefficients of the PS function. |
- Template Parameters
-
GEO | The type of the geometry object. |
AST | The type of the Assumed Stress function |
- Returns
- The PK2 stress vector at the given integration point.
The documentation for this struct was generated from the following file: