KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
List of all members
Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType > Class Template Reference

This object integrates the predictive stress using the plasticity theory by means of linear/exponential softening or hardening + softening evolution laws. More...

#include <generic_cl_integrator_plasticity.h>

Collaboration diagram for Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >:

Public Types

Enum's
enum class  HardeningCurveType {
  LinearSoftening = 0 , ExponentialSoftening = 1 , InitialHardeningExponentialSoftening = 2 , PerfectPlasticity = 3 ,
  CurveFittingHardening = 4 , LinearExponentialSoftening = 5 , CurveDefinedByPoints = 6
}
 

Public Member Functions

Life Cycle
 GenericConstitutiveLawIntegratorPlasticity ()
 Initialization constructor. More...
 
 GenericConstitutiveLawIntegratorPlasticity (GenericConstitutiveLawIntegratorPlasticity const &rOther)
 Copy constructor. More...
 
GenericConstitutiveLawIntegratorPlasticityoperator= (GenericConstitutiveLawIntegratorPlasticity const &rOther)
 Assignment operator. More...
 
virtual ~GenericConstitutiveLawIntegratorPlasticity ()
 Destructor. More...
 

Static Public Member Functions

Operations
static void IntegrateStressVector (array_1d< double, VoigtSize > &rPredictiveStressVector, Vector &rStrainVector, double &rUniaxialStress, double &rThreshold, double &rPlasticDenominator, array_1d< double, VoigtSize > &rFflux, array_1d< double, VoigtSize > &rGflux, double &rPlasticDissipation, array_1d< double, VoigtSize > &rPlasticStrainIncrement, Matrix &rConstitutiveMatrix, Vector &rPlasticStrain, ConstitutiveLaw::Parameters &rValues, const double CharacteristicLength)
 This method integrates the predictive stress vector with the CL using differents evolution laws using the backward euler scheme. More...
 
static double CalculatePlasticParameters (array_1d< double, VoigtSize > &rPredictiveStressVector, Vector &rStrainVector, double &rUniaxialStress, double &rThreshold, double &rPlasticDenominator, array_1d< double, VoigtSize > &rFflux, array_1d< double, VoigtSize > &rGflux, double &rPlasticDissipation, array_1d< double, VoigtSize > &rPlasticStrainIncrement, const Matrix &rConstitutiveMatrix, ConstitutiveLaw::Parameters &rValues, const double CharacteristicLength, const Vector &rPlasticStrain)
 This method calculates all the plastic parameters required for the integration of the PredictiveStressVector. More...
 
static void CalculateTangentMatrix (Matrix &rTangent, const Matrix &rElasticMatrix, const array_1d< double, VoigtSize > &rFFluxVector, const array_1d< double, VoigtSize > &rGFluxVector, const double Denominator)
 This method calculates the analytical tangent tensor. More...
 
static void CalculateFFluxVector (const array_1d< double, VoigtSize > &rPredictiveStressVector, const array_1d< double, VoigtSize > &rDeviator, const double J2, array_1d< double, VoigtSize > &rFFluxVector, ConstitutiveLaw::Parameters &rValues)
 This method calculates the derivative of the yield surface. More...
 
static void CalculateGFluxVector (const array_1d< double, VoigtSize > &rPredictiveStressVector, const array_1d< double, VoigtSize > &rDeviator, const double J2, array_1d< double, VoigtSize > &rGFluxVector, ConstitutiveLaw::Parameters &rValues)
 This method calculates the derivative of the plastic potential. More...
 
static void CalculateIndicatorsFactors (const array_1d< double, VoigtSize > &rPredictiveStressVector, double &rTensileIndicatorFactor, double &rCompressionIndicatorFactor)
 This method computes the tensile/compressive indicators. More...
 
static void CalculatePlasticDissipation (const array_1d< double, VoigtSize > &rPredictiveStressVector, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, const Vector &PlasticStrainInc, double &rPlasticDissipation, array_1d< double, VoigtSize > &rHCapa, ConstitutiveLaw::Parameters &rValues, const double CharacteristicLength)
 This method computes the plastic dissipation of the plasticity model. More...
 
static void CalculateEquivalentStressThreshold (const double PlasticDissipation, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, double &rEquivalentStressThreshold, double &rSlope, ConstitutiveLaw::Parameters &rValues, const double EquivalentPlasticStrain, const double CharacteristicLength)
 This method computes the uniaxial threshold that differentiates the elastic-plastic behaviour. More...
 
static void CalculateEquivalentStressThresholdHardeningCurveLinearSoftening (const double PlasticDissipation, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, double &rEquivalentStressThreshold, double &rSlope, ConstitutiveLaw::Parameters &rValues)
 This method computes the uniaxial threshold using a linear softening. More...
 
static void CalculateEquivalentStressThresholdHardeningCurveExponentialSoftening (const double PlasticDissipation, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, double &rEquivalentStressThreshold, double &rSlope, ConstitutiveLaw::Parameters &rValues, const double CharacteristicLength)
 This method computes the uniaxial threshold using a exponential softening. More...
 
static void CalculateEquivalentStressThresholdHardeningCurveInitialHardeningExponentialSoftening (const double PlasticDissipation, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, double &rEquivalentStressThreshold, double &rSlope, ConstitutiveLaw::Parameters &rValues)
 This method computes the uniaxial threshold using a hardening-softening law. More...
 
static void CalculateEquivalentStressThresholdHardeningCurvePerfectPlasticity (const double PlasticDissipation, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, double &rEquivalentStressThreshold, double &rSlope, ConstitutiveLaw::Parameters &rValues)
 This method computes the uniaxial threshold using a perfect plasticity law. More...
 
static void CalculateEquivalentStressThresholdCurveFittingHardening (const double PlasticDissipation, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, double &rEquivalentStressThreshold, double &rSlope, ConstitutiveLaw::Parameters &rValues, const double EquivalentPlasticStrain, const double CharacteristicLength)
 This method computes the uniaxial threshold using a perfect plasticity law. More...
 
static void CalculateEquivalentStressThresholdHardeningCurveLinearExponentialSoftening (const double PlasticDissipation, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, double &rEquivalentStressThreshold, double &rSlope, const double CharacteristicLength, ConstitutiveLaw::Parameters &rValues)
 This method computes the uniaxial threshold using a linear-exponential softening, which changes from one to the other through the platic_dissipation_limit. More...
 
static void CalculateEquivalentStressThresholdHardeningCurveDefinedByPoints (const double PlasticDissipation, const double TensileIndicatorFactor, const double CompressionIndicatorFactor, double &rEquivalentStressThreshold, double &rSlope, ConstitutiveLaw::Parameters &rValues, const double CharacteristicLength)
 This method computes the uniaxial threshold using a two region curve: More...
 
static void CalculateEquivalentPlasticStrain (const Vector &rStressVector, const double UniaxialStress, const Vector &rPlasticStrain, const double r0, ConstitutiveLaw::Parameters &rValues, double &rEquivalentPlasticStrain)
 This method returns the equivalent plastic strain. More...
 
static void GetInitialUniaxialThreshold (ConstitutiveLaw::Parameters &rValues, double &rThreshold)
 This method returns the initial uniaxial stress threshold. More...
 
static void CalculateHardeningParameter (const array_1d< double, VoigtSize > &rGFlux, const double SlopeThreshold, const array_1d< double, VoigtSize > &rHCapa, double &rHardeningParameter)
 This method computes hardening parameter needed for the algorithm. More...
 
static void CalculatePlasticDenominator (const array_1d< double, VoigtSize > &rFFlux, const array_1d< double, VoigtSize > &rGFlux, const Matrix &rConstitutiveMatrix, double &rHardeningParameter, double &rPlasticDenominator)
 This method computes the plastic denominator needed to compute the plastic consistency factor. More...
 
static int Check (const Properties &rMaterialProperties)
 This method defines in the CL integrator. More...
 

Type Definitions

typedef std::size_t IndexType
 Definition of index. More...
 
typedef TYieldSurfaceType YieldSurfaceType
 The type of yield surface. More...
 
typedef YieldSurfaceType::PlasticPotentialType PlasticPotentialType
 The type of plastic potential. More...
 
static constexpr double tolerance = std::numeric_limits<double>::epsilon()
 The machine precision tolerance. More...
 
static constexpr SizeType Dimension = YieldSurfaceType::Dimension
 The define the working dimension size, already defined in the yield surface. More...
 
static constexpr SizeType VoigtSize = YieldSurfaceType::VoigtSize
 The define the Voigt size, already defined in the yield surface. More...
 
 KRATOS_CLASS_POINTER_DEFINITION (GenericConstitutiveLawIntegratorPlasticity)
 Counted pointer of GenericConstitutiveLawIntegratorPlasticity. More...
 

Detailed Description

template<class TYieldSurfaceType>
class Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >

This object integrates the predictive stress using the plasticity theory by means of linear/exponential softening or hardening + softening evolution laws.

The definitions of these classes is completely static, the derivation is done in a static way

Template Parameters
TYieldSurfaceTypeThe yield surface considered The plasticity integrator requires the definition of the following properties:
  • SOFTENING_TYPE: The softening behaviour considered (linear, exponential,etc...)
  • HARDENING_CURVE: The type of considered hardening curve (linear, exponential, pure plastic, etc...)
  • MAXIMUM_STRESS: The maximum stress that defines the exponential hardening
  • MAXIMUM_STRESS_POSITION: The maximum stress position that defines the exponential hardening
  • FRACTURE_ENERGY: A fracture energy-based function is used to describe strength degradation in post-peak regime
  • YOUNG_MODULUS: It defines the relationship between stress (force per unit area) and strain (proportional deformation) in a material in the linear elasticity regime of a uniaxial deformation.
  • YIELD_STRESS: Yield stress is the amount of stress that an object needs to experience for it to be permanently deformed. Does not require to be defined simmetrically, one YIELD_STRESS_COMPRESSION and other YIELD_STRESS_TENSION can be defined for not symmetric cases
Author
Alejandro Cornejo & Lucia Barbu

Member Typedef Documentation

◆ IndexType

template<class TYieldSurfaceType >
typedef std::size_t Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::IndexType

Definition of index.

◆ PlasticPotentialType

template<class TYieldSurfaceType >
typedef YieldSurfaceType::PlasticPotentialType Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::PlasticPotentialType

The type of plastic potential.

◆ YieldSurfaceType

template<class TYieldSurfaceType >
typedef TYieldSurfaceType Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::YieldSurfaceType

The type of yield surface.

Member Enumeration Documentation

◆ HardeningCurveType

template<class TYieldSurfaceType >
enum Kratos::GenericConstitutiveLawIntegratorPlasticity::HardeningCurveType
strong
Enumerator
LinearSoftening 
ExponentialSoftening 
InitialHardeningExponentialSoftening 
PerfectPlasticity 
CurveFittingHardening 
LinearExponentialSoftening 
CurveDefinedByPoints 

Constructor & Destructor Documentation

◆ GenericConstitutiveLawIntegratorPlasticity() [1/2]

template<class TYieldSurfaceType >
Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::GenericConstitutiveLawIntegratorPlasticity ( )
inline

Initialization constructor.

◆ GenericConstitutiveLawIntegratorPlasticity() [2/2]

template<class TYieldSurfaceType >
Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::GenericConstitutiveLawIntegratorPlasticity ( GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType > const &  rOther)
inline

Copy constructor.

◆ ~GenericConstitutiveLawIntegratorPlasticity()

template<class TYieldSurfaceType >
virtual Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::~GenericConstitutiveLawIntegratorPlasticity ( )
inlinevirtual

Destructor.

Member Function Documentation

◆ CalculateEquivalentPlasticStrain()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentPlasticStrain ( const Vector rStressVector,
const double  UniaxialStress,
const Vector rPlasticStrain,
const double  r0,
ConstitutiveLaw::Parameters rValues,
double rEquivalentPlasticStrain 
)
inlinestatic

This method returns the equivalent plastic strain.

Parameters
rThresholdThe uniaxial stress threshold
rValuesParameters of the constitutive law
rStressVectorThe stress vector
r0The tensile indicator
rEquivalentPlasticStrainThe equivalent plastic strain
rPlasticStrainThe plastic strain vector

◆ CalculateEquivalentStressThreshold()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentStressThreshold ( const double  PlasticDissipation,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
double rEquivalentStressThreshold,
double rSlope,
ConstitutiveLaw::Parameters rValues,
const double  EquivalentPlasticStrain,
const double  CharacteristicLength 
)
inlinestatic

This method computes the uniaxial threshold that differentiates the elastic-plastic behaviour.

Parameters
PlasticDissipationThe internal variable of energy dissipation due to plasticity
TensileIndicatorFactorThe tensile indicator
CompressionIndicatorFactorThe compressive indicator
rEquivalentStressThresholdThe maximum uniaxial stress of the linear behaviour
rSlopeThe slope of the PlasticDiss-Threshold curve
rValuesParameters of the constitutive law
EquivalentPlasticStrainThe equivalent plastic strain
CharacteristicLengthCharacteristic length of the finite element

◆ CalculateEquivalentStressThresholdCurveFittingHardening()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentStressThresholdCurveFittingHardening ( const double  PlasticDissipation,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
double rEquivalentStressThreshold,
double rSlope,
ConstitutiveLaw::Parameters rValues,
const double  EquivalentPlasticStrain,
const double  CharacteristicLength 
)
inlinestatic

This method computes the uniaxial threshold using a perfect plasticity law.

Parameters
PlasticDissipationThe internal variable of energy dissipation due to plasticity
TensileIndicatorFactorThe tensile indicator
CompressionIndicatorFactorThe compressive indicator
rEquivalentStressThresholdThe maximum uniaxial stress of the linear behaviour
rSlopeThe slope of the PlasticDiss-Threshold curve
rValuesParameters of the constitutive law
EquivalentPlasticStrainThe Plastic Strain internal variable
CharacteristicLengthCharacteristic length of the finite element

◆ CalculateEquivalentStressThresholdHardeningCurveDefinedByPoints()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentStressThresholdHardeningCurveDefinedByPoints ( const double  PlasticDissipation,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
double rEquivalentStressThreshold,
double rSlope,
ConstitutiveLaw::Parameters rValues,
const double  CharacteristicLength 
)
inlinestatic

This method computes the uniaxial threshold using a two region curve:

  • point defined curve with linear interpolation followed by a
  • exponential curve.
    Parameters
    PlasticDissipationThe internal variable of energy dissipation due to plasticity
    TensileIndicatorFactorThe tensile indicator
    CompressionIndicatorFactorThe compressive indicator
    rEquivalentStressThresholdThe maximum uniaxial stress of the linear behaviour
    rSlopeThe slope of the PlasticDiss-Threshold curve
    rValuesParameters of the constitutive law
    CharacteristicLengthCharacteristic length of the finite element

◆ CalculateEquivalentStressThresholdHardeningCurveExponentialSoftening()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentStressThresholdHardeningCurveExponentialSoftening ( const double  PlasticDissipation,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
double rEquivalentStressThreshold,
double rSlope,
ConstitutiveLaw::Parameters rValues,
const double  CharacteristicLength 
)
inlinestatic

This method computes the uniaxial threshold using a exponential softening.

Parameters
PlasticDissipationThe internal variable of energy dissipation due to plasticity
TensileIndicatorFactorThe tensile indicator
CompressionIndicatorFactorThe compressive indicator
rEquivalentStressThresholdThe maximum uniaxial stress of the linear behaviour
rSlopeThe slope of the PlasticDiss-Threshold curve
rValuesParameters of the constitutive law
CharacteristicLengthCharacteristic length of the finite element

◆ CalculateEquivalentStressThresholdHardeningCurveInitialHardeningExponentialSoftening()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentStressThresholdHardeningCurveInitialHardeningExponentialSoftening ( const double  PlasticDissipation,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
double rEquivalentStressThreshold,
double rSlope,
ConstitutiveLaw::Parameters rValues 
)
inlinestatic

This method computes the uniaxial threshold using a hardening-softening law.

Parameters
PlasticDissipationThe internal variable of energy dissipation due to plasticity
TensileIndicatorFactorThe tensile indicator
CompressionIndicatorFactorThe compressive indicator
rEquivalentStressThresholdThe maximum uniaxial stress of the linear behaviour
rSlopeThe slope of the PlasticDiss-Threshold curve
rValuesParameters of the constitutive law

◆ CalculateEquivalentStressThresholdHardeningCurveLinearExponentialSoftening()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentStressThresholdHardeningCurveLinearExponentialSoftening ( const double  PlasticDissipation,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
double rEquivalentStressThreshold,
double rSlope,
const double  CharacteristicLength,
ConstitutiveLaw::Parameters rValues 
)
inlinestatic

This method computes the uniaxial threshold using a linear-exponential softening, which changes from one to the other through the platic_dissipation_limit.

Parameters
PlasticDissipationThe internal variable of energy dissipation due to plasticity
TensileIndicatorFactorThe tensile indicator
CompressionIndicatorFactorThe compressive indicator
rEquivalentStressThresholdThe maximum uniaxial stress of the linear behaviour
rSlopeThe slope of the PlasticDiss-Threshold curve
CharacteristicLengthCharacteristic length of the finite element
rValuesParameters of the constitutive law

◆ CalculateEquivalentStressThresholdHardeningCurveLinearSoftening()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentStressThresholdHardeningCurveLinearSoftening ( const double  PlasticDissipation,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
double rEquivalentStressThreshold,
double rSlope,
ConstitutiveLaw::Parameters rValues 
)
inlinestatic

This method computes the uniaxial threshold using a linear softening.

Parameters
PlasticDissipationThe internal variable of energy dissipation due to plasticity
TensileIndicatorFactorThe tensile indicator
CompressionIndicatorFactorThe compressive indicator
rEquivalentStressThresholdThe maximum uniaxial stress of the linear behaviour
rSlopeThe slope of the PlasticDiss-Threshold curve
rValuesParameters of the constitutive law

◆ CalculateEquivalentStressThresholdHardeningCurvePerfectPlasticity()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateEquivalentStressThresholdHardeningCurvePerfectPlasticity ( const double  PlasticDissipation,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
double rEquivalentStressThreshold,
double rSlope,
ConstitutiveLaw::Parameters rValues 
)
inlinestatic

This method computes the uniaxial threshold using a perfect plasticity law.

Parameters
PlasticDissipationThe internal variable of energy dissipation due to plasticity
TensileIndicatorFactorThe tensile indicator
CompressionIndicatorFactorThe compressive indicator
rEquivalentStressThresholdThe maximum uniaxial stress of the linear behaviour
rSlopeThe slope of the PlasticDiss-Threshold curve
rValuesParameters of the constitutive law

◆ CalculateFFluxVector()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateFFluxVector ( const array_1d< double, VoigtSize > &  rPredictiveStressVector,
const array_1d< double, VoigtSize > &  rDeviator,
const double  J2,
array_1d< double, VoigtSize > &  rFFluxVector,
ConstitutiveLaw::Parameters rValues 
)
inlinestatic

This method calculates the derivative of the yield surface.

Parameters
rPredictiveStressVectorThe predictive stress vector S = C:(E-Ep)
rDeviatorThe deviatoric part of the stress vector
J2The second invariant of the deviatoric part of the stress vector
rFFluxVectorThe derivative of the yield surface
rValuesParameters of the constitutive law

◆ CalculateGFluxVector()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateGFluxVector ( const array_1d< double, VoigtSize > &  rPredictiveStressVector,
const array_1d< double, VoigtSize > &  rDeviator,
const double  J2,
array_1d< double, VoigtSize > &  rGFluxVector,
ConstitutiveLaw::Parameters rValues 
)
inlinestatic

This method calculates the derivative of the plastic potential.

Parameters
rPredictiveStressVectorThe predictive stress vector S = C:(E-Ep)
rDeviatorThe deviatoric part of the stress vector
J2The second invariant of the deviatoric part of the stress vector
rGFluxVectorThe derivative of the yield surface
rValuesParameters of the constitutive law

◆ CalculateHardeningParameter()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateHardeningParameter ( const array_1d< double, VoigtSize > &  rGFlux,
const double  SlopeThreshold,
const array_1d< double, VoigtSize > &  rHCapa,
double rHardeningParameter 
)
inlinestatic

This method computes hardening parameter needed for the algorithm.

Parameters
rGFluxThe derivative of the plastic potential
SlopeThresholdThe slope of the PlasticDiss-Threshold curve
rHardeningParameterThe hardening parameter needed for the algorithm
rSlopeThe slope of the PlasticDiss-Threshold curve

◆ CalculateIndicatorsFactors()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateIndicatorsFactors ( const array_1d< double, VoigtSize > &  rPredictiveStressVector,
double rTensileIndicatorFactor,
double rCompressionIndicatorFactor 
)
inlinestatic

This method computes the tensile/compressive indicators.

Parameters
rPredictiveStressVectorThe predictive stress vector S = C:(E-Ep)
rTensileIndicatorFactorThe tensile indicator
rCompressionIndicatorFactorThe compressive indicator

◆ CalculatePlasticDenominator()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculatePlasticDenominator ( const array_1d< double, VoigtSize > &  rFFlux,
const array_1d< double, VoigtSize > &  rGFlux,
const Matrix rConstitutiveMatrix,
double rHardeningParameter,
double rPlasticDenominator 
)
inlinestatic

This method computes the plastic denominator needed to compute the plastic consistency factor.

Parameters
rFfluxThe derivative of the yield surface
rGfluxThe derivative of the plastic potential
rConstitutiveMatrixThe elastic constitutive matrix
rHardeningParameterThe hardening parameter needed for the algorithm
rPlasticDenominatorThe plasticity numerical value to obtain the pastic consistency factor

◆ CalculatePlasticDissipation()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculatePlasticDissipation ( const array_1d< double, VoigtSize > &  rPredictiveStressVector,
const double  TensileIndicatorFactor,
const double  CompressionIndicatorFactor,
const Vector PlasticStrainInc,
double rPlasticDissipation,
array_1d< double, VoigtSize > &  rHCapa,
ConstitutiveLaw::Parameters rValues,
const double  CharacteristicLength 
)
inlinestatic

This method computes the plastic dissipation of the plasticity model.

Parameters
rPredictiveStressVectorThe predictive stress vector S = C : (E-Ep)
TensileIndicatorFactorThe tensile indicator
CompressionIndicatorFactorThe compressive indicator
PlasticStrainIncrementThe increment of plastic strain of this time step
PlasticDissipationThe internal variable of energy dissipation due to plasticity
rHCapaThe slope of the PlasticDiss-Threshold curve
rValuesParameters of the constitutive law
CharacteristicLengthThe equivalent length of the FE

◆ CalculatePlasticParameters()

template<class TYieldSurfaceType >
static double Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculatePlasticParameters ( array_1d< double, VoigtSize > &  rPredictiveStressVector,
Vector rStrainVector,
double rUniaxialStress,
double rThreshold,
double rPlasticDenominator,
array_1d< double, VoigtSize > &  rFflux,
array_1d< double, VoigtSize > &  rGflux,
double rPlasticDissipation,
array_1d< double, VoigtSize > &  rPlasticStrainIncrement,
const Matrix rConstitutiveMatrix,
ConstitutiveLaw::Parameters rValues,
const double  CharacteristicLength,
const Vector rPlasticStrain 
)
inlinestatic

This method calculates all the plastic parameters required for the integration of the PredictiveStressVector.

Parameters
rPredictiveStressVectorThe predictive stress vector S = C:(E-Ep)
rStrainVectorThe equivalent strain vector of that integration point
rUniaxialStressThe equivalent uniaxial stress
rThresholdThe maximum uniaxial stress of the linear behaviour
rPlasticDenominatorThe plasticity numerical value to obtain the pastic consistency factor
rFfluxThe derivative of the yield surface
rGfluxThe derivative of the plastic potential
rPlasticDissipationThe internal variable of energy dissipation due to plasticity
rPlasticStrainIncrementThe increment of plastic strain of this time step
rConstitutiveMatrixThe elastic constitutive matrix
rValuesParameters of the constitutive law
CharacteristicLengthThe equivalent length of the FE
Returns
The threshold of plasticity

◆ CalculateTangentMatrix()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::CalculateTangentMatrix ( Matrix rTangent,
const Matrix rElasticMatrix,
const array_1d< double, VoigtSize > &  rFFluxVector,
const array_1d< double, VoigtSize > &  rGFluxVector,
const double  Denominator 
)
inlinestatic

This method calculates the analytical tangent tensor.

Parameters
rPredictiveStressVectorThe predictive stress vector S = C:(E-Ep)
rDeviatorThe deviatoric part of the stress vector
J2The second invariant of the deviatoric part of the stress vector
rFFluxVectorThe derivative of the yield surface
rValuesParameters of the constitutive law

◆ Check()

template<class TYieldSurfaceType >
static int Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::Check ( const Properties rMaterialProperties)
inlinestatic

This method defines in the CL integrator.

Returns
0 if OK, 1 otherwise

◆ GetInitialUniaxialThreshold()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::GetInitialUniaxialThreshold ( ConstitutiveLaw::Parameters rValues,
double rThreshold 
)
inlinestatic

This method returns the initial uniaxial stress threshold.

Parameters
rThresholdThe uniaxial stress threshold
rValuesParameters of the constitutive law

◆ IntegrateStressVector()

template<class TYieldSurfaceType >
static void Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::IntegrateStressVector ( array_1d< double, VoigtSize > &  rPredictiveStressVector,
Vector rStrainVector,
double rUniaxialStress,
double rThreshold,
double rPlasticDenominator,
array_1d< double, VoigtSize > &  rFflux,
array_1d< double, VoigtSize > &  rGflux,
double rPlasticDissipation,
array_1d< double, VoigtSize > &  rPlasticStrainIncrement,
Matrix rConstitutiveMatrix,
Vector rPlasticStrain,
ConstitutiveLaw::Parameters rValues,
const double  CharacteristicLength 
)
inlinestatic

This method integrates the predictive stress vector with the CL using differents evolution laws using the backward euler scheme.

Parameters
rPredictiveStressVectorThe predictive stress vector S = C:(E-Ep)
rStrainVectorThe equivalent strain vector of that integration point
rUniaxialStressThe equivalent uniaxial stress
rThresholdThe maximum uniaxial stress of the linear behaviour
rPlasticDenominatorThe plasticity numerical value to obtain the pastic consistency factor
rFfluxThe derivative of the yield surface
rGfluxThe derivative of the plastic potential
rPlasticDissipationThe internal variable of energy dissipation due to plasticity
rPlasticStrainIncrementThe increment of plastic strain of this time step
rConstitutiveMatrixThe elastic constitutive matrix
rPlasticStrainThe elastic constitutive matrix
rValuesParameters of the constitutive law
CharacteristicLengthThe equivalent length of the FE

◆ KRATOS_CLASS_POINTER_DEFINITION()

template<class TYieldSurfaceType >
Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::KRATOS_CLASS_POINTER_DEFINITION ( GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >  )

◆ operator=()

template<class TYieldSurfaceType >
GenericConstitutiveLawIntegratorPlasticity& Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::operator= ( GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType > const &  rOther)
inline

Assignment operator.

Member Data Documentation

◆ Dimension

template<class TYieldSurfaceType >
constexpr SizeType Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::Dimension = YieldSurfaceType::Dimension
staticconstexpr

The define the working dimension size, already defined in the yield surface.

◆ tolerance

template<class TYieldSurfaceType >
constexpr double Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::tolerance = std::numeric_limits<double>::epsilon()
staticconstexpr

The machine precision tolerance.

◆ VoigtSize

template<class TYieldSurfaceType >
constexpr SizeType Kratos::GenericConstitutiveLawIntegratorPlasticity< TYieldSurfaceType >::VoigtSize = YieldSurfaceType::VoigtSize
staticconstexpr

The define the Voigt size, already defined in the yield surface.


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