Public Member Functions | List of all members
MarmotMaterialHyperElastic Class Referenceabstract

#include <MarmotMaterialHyperElastic.h>

Inheritance diagram for MarmotMaterialHyperElastic:
[legend]

Public Member Functions

virtual void computeStress (double *S, double *dSdE, const double *FOld, const double *FNew, const double *timeOld, const double dT, double &pNewDT) override
 
virtual void computeStressPK2 (double *S, double *dSdE, const double *E, const double *timeOld, const double dT, double &pNewDT)=0
 
virtual void computePlaneStressPK2 (double *S2D, double *dSdE2D, const double *E2D, const double *timeOld, const double dT, double &pNewDT)
 
virtual void computeUniaxialStressPK2 (double *S1D, double *dSdE1D, const double *E1D, const double *timeOld, const double dT, double &pNewDT)
 
- Public Member Functions inherited from MarmotMaterialMechanical
virtual void computePlaneStress (double *stress2D, double *dStress_dF2DNew, const double *FOld2D, const double *FNew2D, const double *timeOld, const double dT, double &pNewDT)
 
virtual void computeUniaxialStress (double *stress1D, double *dStress1D_dF1DNew, const double *F1DOld, const double *F1DNew, const double *timeOld, const double dT, double &pNewDT)
 
 MarmotMaterial (const double *materialProperties, int nMaterialProperties, int materialNumber)
 
- Public Member Functions inherited from MarmotMaterial
 MarmotMaterial (const double *materialProperties, int nMaterialProperties, int materialNumber)
 
virtual ~MarmotMaterial ()
 
virtual int getNumberOfRequiredStateVars ()=0
 
virtual void assignStateVars (double *stateVars, int nStateVars)
 
virtual StateView getStateView (const std::string &stateName)=0
 
double * getAssignedStateVars ()
 
int getNumberOfAssignedStateVars ()
 
virtual void initializeYourself ()
 

Additional Inherited Members

- Public Attributes inherited from MarmotMaterial
const int materialNumber
 
- Protected Attributes inherited from MarmotMaterial
const double * materialProperties
 
const int nMaterialProperties
 
double * stateVars
 
int nStateVars
 

Detailed Description

Derived abstract base class for simple, purely hyperelastic materials to be used for finite elements based on the total lagrangian kinematic description (TL elements). The second Piola - Kirchhoff stress tensor \( S \) will be derived by

\[ \displaystyle S = \frac{\partial f(\boldsymbol{E},t )}{\partial \boldsymbol{E}} \]

with the Green - Lagrange strain tensor \( \boldsymbol{E} \)

\[ \displaystyle E = \frac{1}{2}\,\left(\boldsymbol{F}^T\cdot \boldsymbol{F} - \boldsymbol{I} \right) \]

as work conjugated measure and the variable \( \boldsymbol{F} \) denoting the deformation gradient. The algorithmic tangent will be calculated by

\[ \displaystyle \frac{d \boldsymbol{S}}{d \boldsymbol{E}} \]

Member Function Documentation

◆ computeStress()

void MarmotMaterialHyperElastic::computeStress ( double *  S,
double *  dSdE,
const double *  FOld,
const double *  FNew,
const double *  timeOld,
const double  dT,
double &  pNewDT 
)
overridevirtual

For a given deformation gradient at the old and the current time, compute the 2nd Piola-Kirchhoff stress and the algorithmic tangent \(\frac{\partial\boldsymbol{S}^{(n+1)}}{\partial\boldsymbol{E}^{(n+1)}}\).

Todo:
A default implementation is provided.
Parameters
[in,out]S2nd Piola-Kirchhoff stress
[in,out]dSdEAlgorithmic tangent representing the derivative of the 2nd Piola-Kirchhoff stress tensor with respect to the Green-Lagrange strain tensor \(\boldsymbol{E}\).
[in]FOldDeformation gradient at the old (pseudo-)time
[in]FNewDeformation gradient at the current (pseudo-)time
[in]timeOldOld (pseudo-)time
[in]dt(Pseudo-)time increment from the old (pseudo-)time to the current (pseudo-)time
[in,out]pNewDTSuggestion for a new time increment

Implements MarmotMaterialMechanical.

◆ computeStressPK2()

virtual void MarmotMaterialHyperElastic::computeStressPK2 ( double *  S,
double *  dSdE,
const double *  E,
const double *  timeOld,
const double  dT,
double &  pNewDT 
)
pure virtual

For a given Green-Lagrange strain, compute the 2nd Piola-Kirchhoff stress and the algorithmic tangent \(\frac{\partial\boldsymbol{S}^{(n+1)}}{\partial\boldsymbol{E}^{(n+1)}}\).

Todo:
Should we use function overloading in this case and simple also use computeStress for the function name?
Parameters
[in,out]S2nd Piola-Kirchhoff stress
[in,out]dSdEAlgorithmic tangent representing the derivative of the 2nd Piola-Kirchhoff stress tensor with respect to the Green-Lagrange strain tensor \(\boldsymbol{E}\).
[in]deltaEGreen-Lagrange strain increment
[in]timeOldOld (pseudo-)time
[in]dt(Pseudo-)time increment from the old (pseudo-)time to the current (pseudo-)time
[in,out]pNewDTSuggestion for a new time increment

◆ computePlaneStressPK2()

void MarmotMaterialHyperElastic::computePlaneStressPK2 ( double *  S2D,
double *  dSdE2D,
const double *  E2D,
const double *  timeOld,
const double  dT,
double &  pNewDT 
)
virtual

Plane stress implementation of computeStressPK2.

◆ computeUniaxialStressPK2()

void MarmotMaterialHyperElastic::computeUniaxialStressPK2 ( double *  S1D,
double *  dSdE1D,
const double *  E1D,
const double *  timeOld,
const double  dT,
double &  pNewDT 
)
virtual

Uniaxial stress implementation of computeStressPK2.


The documentation for this class was generated from the following files: