Public Member Functions | List of all members
MarmotMaterialGradientEnhancedHypoElastic Class Referenceabstract

#include <MarmotMaterialGradientEnhancedHypoElastic.h>

Inheritance diagram for MarmotMaterialGradientEnhancedHypoElastic:
[legend]

Public Member Functions

virtual void computeStress (double *stress, double &KLocal, double &nonLocalRadius, double *dStress_dDeformationGradient, double *dKLocal_dDeformationGradient, double *dStress_dK, const double *FOld, const double *FNew, const double KOld, const double dK, const double *timeOld, const double dT, double &pNewDT) override
 
virtual void computeStress (double *stress, double &K_local, double &nonLocalRadius, double *dStressDDStrain, double *dK_localDDStrain, double *dStressDK, const double *dStrain, double KOld, double dK, const double *timeOld, const double dT, double &pNewDT)=0
 
virtual void computePlaneStress (double *stress2D, double &KLocal2D, double &nonLocalRadius, double *dStress_DStrain2D, double *dKLocal_dStrain2D, double *dStress_dK2D, const double *dStrain2D, double KOld, double dK, const double *timeOld, const double dT, double &pNewDT)
 
virtual void computePlaneStress (double *stress, double &K_local, double &nonLocalRadius, double *dStressDDFNew, double *dK_localDDFNew, double *dStressDK, const double *FOld, const double *FNew, double KOld, double dK, const double *timeOld, const double dT, double &pNewDT)
 
- Public Member Functions inherited from MarmotMaterialGradientEnhancedMechanical
virtual void computePlaneStress (double *stress, double &K_local, double &nonLocalRadius, double *dStressDDFNew, double *dK_localDDFNew, double *dStressDK, const double *FOld, const double *FNew, double KOld, double dK, const double *timeOld, const double dT, double &pNewDT)
 
virtual void computeUniaxialStress (double *stress, double &K_local, double &nonLocalRadius, double *dStressDDFNew, double *dK_localDDFNew, double *dStressDK, const double *FOld, const double *FNew, double KOld, double dK, 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 gradient-enhanced hypo-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 \) the local variable \(\kappa\), the nonlocal variable \(\bar{\kappa}\), internal variables \(\boldsymbol{x}\) and the time \( t \)

\[ \displaystyle \sigRate = f( \sig, \kappa, \epsRate, \dot{\bar{\kappa}}, \boldsymbol{x}, t, ...) \]

In course of numerical time integration, this relation will be formulated incrementally as

\[ \displaystyle \Delta \sig = f ( \sig_n, \kappa_n, \Delta\eps, \Delta\bar{\kappa}, \boldsymbol{x}_n, \Delta t, t_n, ...) \]

with

\[ \displaystyle \Delta\eps = \epsRate\, \Delta t \]

\[ \displaystyle \Delta\bar{\kappa} = \dot{\bar{\kappa}}\, \Delta t \]

and the algorithmic tangents

\[ \displaystyle \frac{d \sig }{d \eps } = \frac{d \Delta \sig }{d \Delta \eps } \]

\[ \displaystyle \frac{d \sig }{d \bar{\kappa}} = \frac{d \Delta \sig }{d \Delta \bar{\kappa}} \]

\[ \displaystyle \frac{d \kappa }{d \eps }= \frac{d \Delta \kappa }{d \Delta \eps } \]

\[ \displaystyle \frac{d \kappa }{d \bar{\kappa}}= \frac{d \Delta \kappa }{d \Delta \bar{\kappa}} = l_{\text{nonlocal}} \]

Member Function Documentation

◆ computeStress() [1/2]

void MarmotMaterialGradientEnhancedHypoElastic::computeStress ( double *  stress,
double &  KLocal,
double &  nonLocalRadius,
double *  dStress_dDeformationGradient,
double *  dKLocal_dDeformationGradient,
double *  dStress_dK,
const double *  FOld,
const double *  FNew,
const double  KOld,
const double  dK,
const double *  timeOld,
const double  dT,
double &  pNewDT 
)
overridevirtual

For a given deformation gradient and a nonlocal variable at the old and the current time, compute the Cauchy stress and the local variable and the algorithmic tangents.

Parameters
[in,out]stressCauchy stress
[in,out]K_locallocal variable
[in,out]nonLocalRadiusnonlocal radius representing the tangent \(\frac{d \Delta \kappa }{d \Delta \bar{\kappa}}\)
[in,out]dStressDDDeformationGradientDerivative of the Cauchy stress tensor with respect to the deformation gradient \(\boldsymbol{F}\)
[in,out]dK_localDDeformationGradientDerivative of the local variable with respect to the deformation gradient \(\boldsymbol{F}\)
[in,out]dStressDKDerivative of the Cauchy stress tensor with respect to the nonlocal variable
[in]FOldDeformation gradient at the old (pseudo-)time
[in]FNewDeformation gradient at the current (pseudo-)time
[in]KOldLocal variable at the old (pseudo-)time
[in]dKIncrement of the local variable
[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 MarmotMaterialGradientEnhancedMechanical.

◆ computeStress() [2/2]

virtual void MarmotMaterialGradientEnhancedHypoElastic::computeStress ( double *  stress,
double &  K_local,
double &  nonLocalRadius,
double *  dStressDDStrain,
double *  dK_localDDStrain,
double *  dStressDK,
const double *  dStrain,
double  KOld,
double  dK,
const double *  timeOld,
const double  dT,
double &  pNewDT 
)
pure virtual

For a given deformation gradient and a nonlocal variable at the old and the current time, compute the Cauchy stress and the local variable and the algorithmic tangents.

Parameters
[in,out]stressCauchy stress
[in,out]K_locallocal variable
[in,out]nonLocalRadiusnonlocal radius representing the tangent \(\frac{d \Delta \kappa }{d \Delta \bar{\kappa}}\)
[in,out]dStressDDstrainAlgorithmic tangent representing the derivative of the Cauchy stress tensor with respect to the linearized strain
[in,out]dK_localDDStrainDerivative of the local variable with respect to the linearized strain
[in,out]dStressDKDerivative of the Cauchy stress tensor with respect to the nonlocal variable
[in]KOldLocal variable at the old (pseudo-)time
[in]dKIncrement of the local variable
[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

◆ computePlaneStress() [1/2]

void MarmotMaterialGradientEnhancedHypoElastic::computePlaneStress ( double *  stress2D,
double &  KLocal2D,
double &  nonLocalRadius,
double *  dStress_DStrain2D,
double *  dKLocal_dStrain2D,
double *  dStress_dK2D,
const double *  dStrain2D,
double  KOld,
double  dK,
const double *  timeOld,
const double  dT,
double &  pNewDT 
)
virtual

Plane stress implementation of computeStress.

◆ computePlaneStress() [2/2]

virtual void MarmotMaterialGradientEnhancedMechanical::computePlaneStress
inline

Plane stress implementation of computeStress.

Todo:
why is dStrain an in-out parameter?

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