#include <MarmotMaterialHypoElasticAD.h>
Public Member Functions | |
void | setCharacteristicElementLength (double length) |
virtual void | computeStress (double *stress, double *dStressDDStrain, const double *FOld, const double *FNew, const double *timeOld, const double dT, double &pNewDT) override |
virtual void | computeStress (autodiff::dual *stress, const autodiff::dual *dStrain, const double *timeOld, const double dT, double &pNewDT)=0 |
virtual void | computePlaneStress (double *stress2D, double *dStress_dStrain2D, const double *dStrain2D, const double *timeOld, const double dT, double &pNewDT) |
virtual void | computeUniaxialStress (double *stress1D, double *dStress_dStrain1D, const double *dStrain, const double *timeOld, const double dT, double &pNewDT) |
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) |
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 () |
Public Attributes | |
double | characteristicElementLength |
Characteristic element length. More... | |
Public Attributes inherited from MarmotMaterial | |
const int | materialNumber |
Additional Inherited Members | |
Protected Attributes inherited from MarmotMaterial | |
const double * | materialProperties |
const int | nMaterialProperties |
double * | stateVars |
int | nStateVars |
Derived abstract base class for elastic materials expressed purely in rate form.
In general, the nominal stress rate tensor \( \sigRate \) can be written as a function of the nominal stress tensor \( \sig \), the stretching rate tensor \( \epsRate \) and the time \( t \).
\[ \displaystyle \sigRate = f( \sig, \epsRate, t, ...) \]
In course of numerical time integration, this relation will be formulated incrementally as
\[ \displaystyle \Delta \sig = f ( \sig_n, \Delta\eps, \Delta t, t_n, ...) \]
with
\[ \displaystyle \Delta\eps = \epsRate\, \Delta t \]
and the algorithmic tangent
\[ \displaystyle \frac{d \sig }{d \eps } = \frac{d \Delta \sig }{d \Delta \eps } \]
This formulation is compatible with an Abaqus interface.
void MarmotMaterialHypoElasticAD::setCharacteristicElementLength | ( | double | length | ) |
Set the characteristic element length at the considered quadrature point. It is needed for the regularization of materials with softening behavior based on the mesh-adjusted softening modulus.
[in] | length | characteristic length; will be assigned to characteristicElementLength |
|
overridevirtual |
For a given deformation gradient at the old and the current time, compute the Cauchy stress and the algorithmic tangent \(\frac{\partial\boldsymbol{\sigma}^{(n+1)}}{\partial\boldsymbol{F}^{(n+1)}}\).
[in,out] | stress | Cauchy stress |
[in,out] | dSdE | Algorithmic tangent representing the derivative of the Cauchy stress tensor with respect to the deformation gradient \(\boldsymbol{F}\) |
[in] | FOld | Deformation gradient at the old (pseudo-)time |
[in] | FNew | Deformation gradient at the current (pseudo-)time |
[in] | timeOld | Old (pseudo-)time |
[in] | dt | (Pseudo-)time increment from the old (pseudo-)time to the current (pseudo-)time |
[in,out] | pNewDT | Suggestion for a new time increment |
Implements MarmotMaterialMechanical.
|
pure virtual |
For a given linearized strain increment \(\Delta\boldsymbol{\varepsilon}\) at the old and the current time, compute the Cauchy stress and the algorithmic tangent \(\frac{\partial\boldsymbol{\sigma}^{(n+1)}}{\partial\boldsymbol{\varepsilon}^{(n+1)}}\).
[in,out] | stress | Cauchy stress |
[in,out] | dStressDDstrain | Algorithmic tangent representing the derivative of the Cauchy stress tensor with respect to the linearized strain |
[in] | dStrain | linearized strain increment |
[in] | timeOld | Old (pseudo-)time |
[in] | dt | (Pseudo-)time increment from the old (pseudo-)time to the current (pseudo-)time |
[in,out] | pNewDT | Suggestion for a new time increment |
|
virtual |
|
virtual |
void MarmotMaterialMechanical::computePlaneStress |
Plane stress implementation of computeStress.
void MarmotMaterialMechanical::computeUniaxialStress |
Uniaxial stress implementation of computeStress.
double MarmotMaterialHypoElasticAD::characteristicElementLength |
Characteristic element length.