https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Static Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Types | Private Member Functions | Private Attributes | List of all members
SolidMechanicsPlasticTensileMulti Class Reference

FiniteStrainTensileMulti implements rate-independent associative tensile failure with hardening/softening in the finite-strain framework, using planar (non-smoothed) surfaces. More...

#include <SolidMechanicsPlasticTensileMulti.h>

Inheritance diagram for SolidMechanicsPlasticTensileMulti:
[legend]

Public Types

typedef DataFileName DataFileParameterType
 

Public Member Functions

 SolidMechanicsPlasticTensileMulti (const InputParameters &parameters)
 
virtual unsigned int numberSurfaces () const override
 The number of yield surfaces for this plasticity model. More...
 
virtual void yieldFunctionV (const RankTwoTensor &stress, Real intnl, std::vector< Real > &f) const override
 Calculates the yield functions. More...
 
virtual void dyieldFunction_dstressV (const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &df_dstress) const override
 The derivative of yield functions with respect to stress. More...
 
virtual void dyieldFunction_dintnlV (const RankTwoTensor &stress, Real intnl, std::vector< Real > &df_dintnl) const override
 The derivative of yield functions with respect to the internal parameter. More...
 
virtual void flowPotentialV (const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &r) const override
 The flow potentials. More...
 
virtual void dflowPotential_dstressV (const RankTwoTensor &stress, Real intnl, std::vector< RankFourTensor > &dr_dstress) const override
 The derivative of the flow potential with respect to stress. More...
 
virtual void dflowPotential_dintnlV (const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &dr_dintnl) const override
 The derivative of the flow potential with respect to the internal parameter. More...
 
virtual void activeConstraints (const std::vector< Real > &f, const RankTwoTensor &stress, Real intnl, const RankFourTensor &Eijkl, std::vector< bool > &act, RankTwoTensor &returned_stress) const override
 The active yield surfaces, given a vector of yield functions. More...
 
virtual std::string modelName () const override
 
virtual bool useCustomReturnMap () const override
 Returns false. You will want to override this in your derived class if you write a custom returnMap function. More...
 
virtual bool useCustomCTO () const override
 Returns false. You will want to override this in your derived class if you write a custom consistent tangent operator function. More...
 
virtual bool returnMap (const RankTwoTensor &trial_stress, Real intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &returned_stress, Real &returned_intnl, std::vector< Real > &dpm, RankTwoTensor &delta_dp, std::vector< Real > &yf, bool &trial_stress_inadmissible) const override
 Performs a custom return-map. More...
 
virtual RankFourTensor consistentTangentOperator (const RankTwoTensor &trial_stress, Real intnl_old, const RankTwoTensor &stress, Real intnl, const RankFourTensor &E_ijkl, const std::vector< Real > &cumulative_pm) const override
 Calculates a custom consistent tangent operator. More...
 
void initialize ()
 
void execute ()
 
void finalize ()
 
virtual void hardPotentialV (const RankTwoTensor &stress, Real intnl, std::vector< Real > &h) const
 The hardening potential. More...
 
virtual void dhardPotential_dstressV (const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &dh_dstress) const
 The derivative of the hardening potential with respect to stress. More...
 
virtual void dhardPotential_dintnlV (const RankTwoTensor &stress, Real intnl, std::vector< Real > &dh_dintnl) const
 The derivative of the hardening potential with respect to the internal parameter. More...
 
bool KuhnTuckerSingleSurface (Real yf, Real dpm, Real dpm_tol) const
 Returns true if the Kuhn-Tucker conditions for the single surface are satisfied. More...
 
SubProblemgetSubProblem () const
 
bool shouldDuplicateInitialExecution () const
 
virtual Real spatialValue (const Point &) const
 
virtual const std::vector< Point > spatialPoints () const
 
void gatherSum (T &value)
 
void gatherMax (T &value)
 
void gatherMin (T &value)
 
void gatherProxyValueMax (T1 &proxy, T2 &value)
 
void gatherProxyValueMin (T1 &proxy, T2 &value)
 
void setPrimaryThreadCopy (UserObject *primary)
 
UserObjectprimaryThreadCopy ()
 
std::set< UserObjectName > getDependObjects () const
 
virtual bool needThreadedCopy () const
 
const std::set< std::string > & getRequestedItems () override
 
const std::set< std::string > & getSuppliedItems () override
 
unsigned int systemNumber () const
 
virtual bool enabled () const
 
std::shared_ptr< MooseObjectgetSharedPtr ()
 
std::shared_ptr< const MooseObjectgetSharedPtr () const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
virtual const std::string & name () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &nm) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
virtual void initialSetup ()
 
virtual void timestepSetup ()
 
virtual void jacobianSetup ()
 
virtual void residualSetup ()
 
virtual void customSetup (const ExecFlagType &)
 
const ExecFlagEnumgetExecuteOnEnum () const
 
UserObjectName getUserObjectName (const std::string &param_name) const
 
const T & getUserObject (const std::string &param_name, bool is_dependency=true) const
 
const T & getUserObjectByName (const UserObjectName &object_name, bool is_dependency=true) const
 
const UserObjectgetUserObjectBase (const std::string &param_name, bool is_dependency=true) const
 
const UserObjectgetUserObjectBaseByName (const UserObjectName &object_name, bool is_dependency=true) const
 
const std::vector< MooseVariableScalar *> & getCoupledMooseScalarVars ()
 
const std::set< TagID > & getScalarVariableCoupleableVectorTags () const
 
const std::set< TagID > & getScalarVariableCoupleableMatrixTags () const
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialProperty (const std::string &name, MaterialData &material_data, const unsigned int state=0)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialProperty (const std::string &name, const unsigned int state=0)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialProperty (const std::string &name, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialProperty (const std::string &name, MaterialData &material_data, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialProperty (const std::string &name, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialProperty (const std::string &name, const unsigned int state=0)
 
const ADMaterialProperty< T > & getADMaterialProperty (const std::string &name, MaterialData &material_data)
 
const ADMaterialProperty< T > & getADMaterialProperty (const std::string &name)
 
const ADMaterialProperty< T > & getADMaterialProperty (const std::string &name)
 
const MaterialProperty< T > & getMaterialPropertyOld (const std::string &name, MaterialData &material_data)
 
const MaterialProperty< T > & getMaterialPropertyOld (const std::string &name)
 
const MaterialProperty< T > & getMaterialPropertyOld (const std::string &name)
 
const MaterialProperty< T > & getMaterialPropertyOlder (const std::string &name, MaterialData &material_data)
 
const MaterialProperty< T > & getMaterialPropertyOlder (const std::string &name)
 
const MaterialProperty< T > & getMaterialPropertyOlder (const std::string &name)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialPropertyByName (const MaterialPropertyName &name, MaterialData &material_data, const unsigned int state)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialPropertyByName (const MaterialPropertyName &name, const unsigned int state=0)
 
const GenericMaterialProperty< T, is_ad > & getGenericMaterialPropertyByName (const MaterialPropertyName &name, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialPropertyByName (const MaterialPropertyName &name, MaterialData &material_data, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialPropertyByName (const MaterialPropertyName &name, const unsigned int state=0)
 
const MaterialProperty< T > & getMaterialPropertyByName (const MaterialPropertyName &name, const unsigned int state=0)
 
const ADMaterialProperty< T > & getADMaterialPropertyByName (const MaterialPropertyName &name, MaterialData &material_data)
 
const ADMaterialProperty< T > & getADMaterialPropertyByName (const MaterialPropertyName &name)
 
const ADMaterialProperty< T > & getADMaterialPropertyByName (const MaterialPropertyName &name)
 
const MaterialProperty< T > & getMaterialPropertyOldByName (const MaterialPropertyName &name, MaterialData &material_data)
 
const MaterialProperty< T > & getMaterialPropertyOldByName (const MaterialPropertyName &name)
 
const MaterialProperty< T > & getMaterialPropertyOldByName (const MaterialPropertyName &name)
 
const MaterialProperty< T > & getMaterialPropertyOlderByName (const MaterialPropertyName &name, MaterialData &material_data)
 
const MaterialProperty< T > & getMaterialPropertyOlderByName (const MaterialPropertyName &name)
 
const MaterialProperty< T > & getMaterialPropertyOlderByName (const MaterialPropertyName &name)
 
std::pair< const MaterialProperty< T > *, std::set< SubdomainID > > getBlockMaterialProperty (const MaterialPropertyName &name)
 
const GenericMaterialProperty< T, is_ad > & getGenericZeroMaterialProperty (const std::string &name)
 
const GenericMaterialProperty< T, is_ad > & getGenericZeroMaterialProperty ()
 
const GenericMaterialProperty< T, is_ad > & getGenericZeroMaterialPropertyByName (const std::string &prop_name)
 
const MaterialProperty< T > & getZeroMaterialProperty (Ts... args)
 
std::set< SubdomainIDgetMaterialPropertyBlocks (const std::string &name)
 
std::vector< SubdomainName > getMaterialPropertyBlockNames (const std::string &name)
 
std::set< BoundaryIDgetMaterialPropertyBoundaryIDs (const std::string &name)
 
std::vector< BoundaryName > getMaterialPropertyBoundaryNames (const std::string &name)
 
void checkBlockAndBoundaryCompatibility (std::shared_ptr< MaterialBase > discrete)
 
std::unordered_map< SubdomainID, std::vector< MaterialBase *> > buildRequiredMaterials (bool allow_stateful=true)
 
void statefulPropertiesAllowed (bool)
 
bool getMaterialPropertyCalled () const
 
virtual const std::unordered_set< unsigned int > & getMatPropDependencies () const
 
virtual void resolveOptionalProperties ()
 
const GenericMaterialProperty< T, is_ad > & getPossiblyConstantGenericMaterialPropertyByName (const MaterialPropertyName &prop_name, MaterialData &material_data, const unsigned int state)
 
bool isImplicit ()
 
Moose::StateArg determineState () const
 
virtual void threadJoin (const UserObject &) override
 
virtual void threadJoin (const UserObject &) override
 
virtual void subdomainSetup () override
 
virtual void subdomainSetup () override
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObject (const std::string &param_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
bool hasUserObjectByName (const UserObjectName &object_name) const
 
const GenericOptionalMaterialProperty< T, is_ad > & getGenericOptionalMaterialProperty (const std::string &name, const unsigned int state=0)
 
const GenericOptionalMaterialProperty< T, is_ad > & getGenericOptionalMaterialProperty (const std::string &name, const unsigned int state=0)
 
const OptionalMaterialProperty< T > & getOptionalMaterialProperty (const std::string &name, const unsigned int state=0)
 
const OptionalMaterialProperty< T > & getOptionalMaterialProperty (const std::string &name, const unsigned int state=0)
 
const OptionalADMaterialProperty< T > & getOptionalADMaterialProperty (const std::string &name)
 
const OptionalADMaterialProperty< T > & getOptionalADMaterialProperty (const std::string &name)
 
const OptionalMaterialProperty< T > & getOptionalMaterialPropertyOld (const std::string &name)
 
const OptionalMaterialProperty< T > & getOptionalMaterialPropertyOld (const std::string &name)
 
const OptionalMaterialProperty< T > & getOptionalMaterialPropertyOlder (const std::string &name)
 
const OptionalMaterialProperty< T > & getOptionalMaterialPropertyOlder (const std::string &name)
 
MaterialBasegetMaterial (const std::string &name)
 
MaterialBasegetMaterial (const std::string &name)
 
MaterialBasegetMaterialByName (const std::string &name, bool no_warn=false)
 
MaterialBasegetMaterialByName (const std::string &name, bool no_warn=false)
 
bool hasMaterialProperty (const std::string &name)
 
bool hasMaterialProperty (const std::string &name)
 
bool hasMaterialPropertyByName (const std::string &name)
 
bool hasMaterialPropertyByName (const std::string &name)
 
bool hasADMaterialProperty (const std::string &name)
 
bool hasADMaterialProperty (const std::string &name)
 
bool hasADMaterialPropertyByName (const std::string &name)
 
bool hasADMaterialPropertyByName (const std::string &name)
 
bool hasGenericMaterialProperty (const std::string &name)
 
bool hasGenericMaterialProperty (const std::string &name)
 
bool hasGenericMaterialPropertyByName (const std::string &name)
 
bool hasGenericMaterialPropertyByName (const std::string &name)
 
const FunctiongetFunction (const std::string &name) const
 
const FunctiongetFunctionByName (const FunctionName &name) const
 
bool hasFunction (const std::string &param_name) const
 
bool hasFunctionByName (const FunctionName &name) const
 
bool isDefaultPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessor (const std::string &param_name, const unsigned int index=0) const
 
bool hasPostprocessorByName (const PostprocessorName &name) const
 
std::size_t coupledPostprocessors (const std::string &param_name) const
 
const PostprocessorName & getPostprocessorName (const std::string &param_name, const unsigned int index=0) const
 
const VectorPostprocessorValuegetVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name, bool needs_broadcast) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const VectorPostprocessorValuegetVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name, bool needs_broadcast) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValue (const std::string &param_name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueOld (const std::string &param_name, const std::string &vector_name) const
 
const ScatterVectorPostprocessorValuegetScatterVectorPostprocessorValueOldByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
bool hasVectorPostprocessor (const std::string &param_name, const std::string &vector_name) const
 
bool hasVectorPostprocessor (const std::string &param_name) const
 
bool hasVectorPostprocessorByName (const VectorPostprocessorName &name, const std::string &vector_name) const
 
bool hasVectorPostprocessorByName (const VectorPostprocessorName &name) const
 
const VectorPostprocessorName & getVectorPostprocessorName (const std::string &param_name) const
 
T & getSampler (const std::string &name)
 
SamplergetSampler (const std::string &name)
 
T & getSamplerByName (const SamplerName &name)
 
SamplergetSamplerByName (const SamplerName &name)
 
virtual void meshChanged ()
 
virtual void meshDisplaced ()
 
PerfGraphperfGraph ()
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValue (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOld (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
const PostprocessorValuegetPostprocessorValueOlder (const std::string &param_name, const unsigned int index=0) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
virtual const PostprocessorValuegetPostprocessorValueByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOldByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
const PostprocessorValuegetPostprocessorValueOlderByName (const PostprocessorName &name) const
 
bool isVectorPostprocessorDistributed (const std::string &param_name) const
 
bool isVectorPostprocessorDistributed (const std::string &param_name) const
 
bool isVectorPostprocessorDistributedByName (const VectorPostprocessorName &name) const
 
bool isVectorPostprocessorDistributedByName (const VectorPostprocessorName &name) const
 
const DistributiongetDistribution (const std::string &name) const
 
const T & getDistribution (const std::string &name) const
 
const DistributiongetDistribution (const std::string &name) const
 
const T & getDistribution (const std::string &name) const
 
const DistributiongetDistributionByName (const DistributionName &name) const
 
const T & getDistributionByName (const std::string &name) const
 
const DistributiongetDistributionByName (const DistributionName &name) const
 
const T & getDistributionByName (const std::string &name) const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static InputParameters validParams ()
 
static void sort (typename std::vector< T > &vector)
 
static void sortDFS (typename std::vector< T > &vector)
 
static void cyclicDependencyError (CyclicDependencyException< T2 > &e, const std::string &header)
 

Public Attributes

const Real _f_tol
 Tolerance on yield function. More...
 
const Real _ic_tol
 Tolerance on internal constraint. More...
 
const ConsoleStream _console
 

Static Public Attributes

static constexpr PropertyValue::id_type default_property_id
 
static constexpr PropertyValue::id_type zero_property_id
 
static constexpr auto SYSTEM
 
static constexpr auto NAME
 

Protected Member Functions

virtual Real tensile_strength (const Real internal_param) const
 tensile strength as a function of residual value, rate, and internal_param More...
 
virtual Real dtensile_strength (const Real internal_param) const
 d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param More...
 
virtual Real yieldFunction (const RankTwoTensor &stress, Real intnl) const
 The following functions are what you should override when building single-plasticity models. More...
 
virtual RankTwoTensor dyieldFunction_dstress (const RankTwoTensor &stress, Real intnl) const
 The derivative of yield function with respect to stress. More...
 
virtual Real dyieldFunction_dintnl (const RankTwoTensor &stress, Real intnl) const
 The derivative of yield function with respect to the internal parameter. More...
 
virtual RankTwoTensor flowPotential (const RankTwoTensor &stress, Real intnl) const
 The flow potential. More...
 
virtual RankFourTensor dflowPotential_dstress (const RankTwoTensor &stress, Real intnl) const
 The derivative of the flow potential with respect to stress. More...
 
virtual RankTwoTensor dflowPotential_dintnl (const RankTwoTensor &stress, Real intnl) const
 The derivative of the flow potential with respect to the internal parameter. More...
 
virtual Real hardPotential (const RankTwoTensor &stress, Real intnl) const
 The hardening potential. More...
 
virtual RankTwoTensor dhardPotential_dstress (const RankTwoTensor &stress, Real intnl) const
 The derivative of the hardening potential with respect to stress. More...
 
virtual Real dhardPotential_dintnl (const RankTwoTensor &stress, Real intnl) const
 The derivative of the hardening potential with respect to the internal parameter. More...
 
virtual void addPostprocessorDependencyHelper (const PostprocessorName &name) const override
 
virtual void addVectorPostprocessorDependencyHelper (const VectorPostprocessorName &name) const override
 
virtual void addUserObjectDependencyHelper (const UserObject &uo) const override
 
void addReporterDependencyHelper (const ReporterName &reporter_name) override
 
const ReporterNamegetReporterName (const std::string &param_name) const
 
T & declareRestartableData (const std::string &data_name, Args &&... args)
 
ManagedValue< T > declareManagedRestartableDataWithContext (const std::string &data_name, void *context, Args &&... args)
 
const T & getRestartableData (const std::string &data_name) const
 
T & declareRestartableDataWithContext (const std::string &data_name, void *context, Args &&... args)
 
T & declareRecoverableData (const std::string &data_name, Args &&... args)
 
T & declareRestartableDataWithObjectName (const std::string &data_name, const std::string &object_name, Args &&... args)
 
T & declareRestartableDataWithObjectNameWithContext (const std::string &data_name, const std::string &object_name, void *context, Args &&... args)
 
std::string restartableName (const std::string &data_name) const
 
const T & getMeshProperty (const std::string &data_name, const std::string &prefix)
 
const T & getMeshProperty (const std::string &data_name)
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name) const
 
bool hasMeshProperty (const std::string &data_name) const
 
std::string meshPropertyName (const std::string &data_name) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level) const
 
PerfID registerTimedSection (const std::string &section_name, const unsigned int level, const std::string &live_message, const bool print_dots=true) const
 
std::string timedSectionName (const std::string &section_name) const
 
bool isCoupledScalar (const std::string &var_name, unsigned int i=0) const
 
unsigned int coupledScalarComponents (const std::string &var_name) const
 
unsigned int coupledScalar (const std::string &var_name, unsigned int comp=0) const
 
libMesh::Order coupledScalarOrder (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarValue (const std::string &var_name, unsigned int comp=0) const
 
const ADVariableValueadCoupledScalarValue (const std::string &var_name, unsigned int comp=0) const
 
const GenericVariableValue< is_ad > & coupledGenericScalarValue (const std::string &var_name, unsigned int comp=0) const
 
const GenericVariableValue< false > & coupledGenericScalarValue (const std::string &var_name, const unsigned int comp) const
 
const GenericVariableValue< true > & coupledGenericScalarValue (const std::string &var_name, const unsigned int comp) const
 
const VariableValuecoupledVectorTagScalarValue (const std::string &var_name, TagID tag, unsigned int comp=0) const
 
const VariableValuecoupledMatrixTagScalarValue (const std::string &var_name, TagID tag, unsigned int comp=0) const
 
const VariableValuecoupledScalarValueOld (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarValueOlder (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDot (const std::string &var_name, unsigned int comp=0) const
 
const ADVariableValueadCoupledScalarDot (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotDot (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotOld (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotDotOld (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotDu (const std::string &var_name, unsigned int comp=0) const
 
const VariableValuecoupledScalarDotDotDu (const std::string &var_name, unsigned int comp=0) const
 
const MooseVariableScalargetScalarVar (const std::string &var_name, unsigned int comp) const
 
virtual void checkMaterialProperty (const std::string &name, const unsigned int state)
 
void markMatPropRequested (const std::string &)
 
MaterialPropertyName getMaterialPropertyName (const std::string &name) const
 
void checkExecutionStage ()
 
const T & getReporterValue (const std::string &param_name, const std::size_t time_index=0)
 
const T & getReporterValue (const std::string &param_name, ReporterMode mode, const std::size_t time_index=0)
 
const T & getReporterValue (const std::string &param_name, const std::size_t time_index=0)
 
const T & getReporterValue (const std::string &param_name, ReporterMode mode, const std::size_t time_index=0)
 
const T & getReporterValueByName (const ReporterName &reporter_name, const std::size_t time_index=0)
 
const T & getReporterValueByName (const ReporterName &reporter_name, ReporterMode mode, const std::size_t time_index=0)
 
const T & getReporterValueByName (const ReporterName &reporter_name, const std::size_t time_index=0)
 
const T & getReporterValueByName (const ReporterName &reporter_name, ReporterMode mode, const std::size_t time_index=0)
 
bool hasReporterValue (const std::string &param_name) const
 
bool hasReporterValue (const std::string &param_name) const
 
bool hasReporterValue (const std::string &param_name) const
 
bool hasReporterValue (const std::string &param_name) const
 
bool hasReporterValueByName (const ReporterName &reporter_name) const
 
bool hasReporterValueByName (const ReporterName &reporter_name) const
 
bool hasReporterValueByName (const ReporterName &reporter_name) const
 
bool hasReporterValueByName (const ReporterName &reporter_name) const
 
const GenericMaterialProperty< T, is_ad > * defaultGenericMaterialProperty (const std::string &name)
 
const GenericMaterialProperty< T, is_ad > * defaultGenericMaterialProperty (const std::string &name)
 
const MaterialProperty< T > * defaultMaterialProperty (const std::string &name)
 
const MaterialProperty< T > * defaultMaterialProperty (const std::string &name)
 
const ADMaterialProperty< T > * defaultADMaterialProperty (const std::string &name)
 
const ADMaterialProperty< T > * defaultADMaterialProperty (const std::string &name)
 

Static Protected Member Functions

static std::string meshPropertyName (const std::string &data_name, const std::string &prefix)
 

Protected Attributes

SubProblem_subproblem
 
FEProblemBase_fe_problem
 
SystemBase_sys
 
const THREAD_ID _tid
 
Assembly_assembly
 
const Moose::CoordinateSystemType_coord_sys
 
const bool _duplicate_initial_execution
 
std::set< std::string > _depend_uo
 
const bool & _enabled
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
const ExecFlagEnum_execute_enum
 
const ExecFlagType_current_execute_flag
 
MooseApp_restartable_app
 
const std::string _restartable_system_name
 
const THREAD_ID _restartable_tid
 
const bool _restartable_read_only
 
FEProblemBase_mci_feproblem
 
FEProblemBase_mdi_feproblem
 
MooseApp_pg_moose_app
 
const std::string _prefix
 
FEProblemBase_sc_fe_problem
 
const THREAD_ID _sc_tid
 
const Real_real_zero
 
const VariableValue_scalar_zero
 
const Point & _point_zero
 
const InputParameters_mi_params
 
const std::string _mi_name
 
const MooseObjectName _mi_moose_object_name
 
FEProblemBase_mi_feproblem
 
SubProblem_mi_subproblem
 
const THREAD_ID _mi_tid
 
const Moose::MaterialDataType _material_data_type
 
MaterialData_material_data
 
bool _stateful_allowed
 
bool _get_material_property_called
 
std::vector< std::unique_ptr< PropertyValue > > _default_properties
 
std::unordered_set< unsigned int_material_property_dependencies
 
const MaterialPropertyName _get_suffix
 
const bool _use_interpolated_state
 
const InputParameters_ti_params
 
FEProblemBase_ti_feproblem
 
bool _is_implicit
 
Real_t
 
const Real_t_old
 
int_t_step
 
Real_dt
 
Real_dt_old
 
bool _is_transient
 
const Parallel::Communicator & _communicator
 

Static Protected Attributes

static const std::string _interpolated_old
 
static const std::string _interpolated_older
 

Private Types

enum  return_type { tip = 0, edge = 1, plane = 2 }
 

Private Member Functions

Real dot (const std::vector< Real > &a, const std::vector< Real > &b) const
 dot product of two 3-dimensional vectors More...
 
Real triple (const std::vector< Real > &a, const std::vector< Real > &b, const std::vector< Real > &c) const
 triple product of three 3-dimensional vectors More...
 
bool returnTip (const std::vector< Real > &eigvals, const std::vector< RealVectorValue > &n, std::vector< Real > &dpm, RankTwoTensor &returned_stress, Real intnl_old, Real initial_guess) const
 Tries to return-map to the Tensile tip. More...
 
bool returnEdge (const std::vector< Real > &eigvals, const std::vector< RealVectorValue > &n, std::vector< Real > &dpm, RankTwoTensor &returned_stress, Real intnl_old, Real initial_guess) const
 Tries to return-map to the Tensile edge. More...
 
bool returnPlane (const std::vector< Real > &eigvals, const std::vector< RealVectorValue > &n, std::vector< Real > &dpm, RankTwoTensor &returned_stress, Real intnl_old, Real initial_guess) const
 Tries to return-map to the Tensile plane The return value is true if the internal Newton-Raphson process has converged, otherwise it is false. More...
 
bool KuhnTuckerOK (const RankTwoTensor &returned_diagonal_stress, const std::vector< Real > &dpm, Real str, Real ep_plastic_tolerance) const
 Returns true if the Kuhn-Tucker conditions are satisfied. More...
 
virtual bool doReturnMap (const RankTwoTensor &trial_stress, Real intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &returned_stress, Real &returned_intnl, std::vector< Real > &dpm, RankTwoTensor &delta_dp, std::vector< Real > &yf, bool &trial_stress_inadmissible) const
 Just like returnMap, but a protected interface that definitely uses the algorithm, since returnMap itself does not use the algorithm if _use_returnMap=false. More...
 

Private Attributes

const SolidMechanicsHardeningModel_strength
 
const unsigned int _max_iters
 maximum iterations allowed in the custom return-map algorithm More...
 
const Real _shift
 yield function is shifted by this amount to avoid problems with stress-derivatives at equal eigenvalues More...
 
const bool _use_custom_returnMap
 Whether to use the custom return-map algorithm. More...
 
const bool _use_custom_cto
 Whether to use the custom consistent tangent operator calculation. More...
 

Detailed Description

FiniteStrainTensileMulti implements rate-independent associative tensile failure with hardening/softening in the finite-strain framework, using planar (non-smoothed) surfaces.

Definition at line 19 of file SolidMechanicsPlasticTensileMulti.h.

Member Enumeration Documentation

◆ return_type

Constructor & Destructor Documentation

◆ SolidMechanicsPlasticTensileMulti()

SolidMechanicsPlasticTensileMulti::SolidMechanicsPlasticTensileMulti ( const InputParameters parameters)

Definition at line 53 of file SolidMechanicsPlasticTensileMulti.C.

56  _strength(getUserObject<SolidMechanicsHardeningModel>("tensile_strength")),
57  _max_iters(getParam<unsigned int>("max_iterations")),
58  _shift(parameters.isParamValid("shift") ? getParam<Real>("shift") : _f_tol),
59  _use_custom_returnMap(getParam<bool>("use_custom_returnMap")),
60  _use_custom_cto(getParam<bool>("use_custom_cto"))
61 {
62  if (_shift < 0)
63  mooseError("Value of 'shift' in SolidMechanicsPlasticTensileMulti must not be negative\n");
64  if (_shift > _f_tol)
65  _console << "WARNING: value of 'shift' in SolidMechanicsPlasticTensileMulti is probably set "
66  "too high"
67  << std::endl;
68  if (LIBMESH_DIM != 3)
69  mooseError("SolidMechanicsPlasticTensileMulti is only defined for LIBMESH_DIM=3");
71 }
void seed(std::size_t i, unsigned int seed)
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.
SolidMechanicsPlasticModel(const InputParameters &parameters)
const SolidMechanicsHardeningModel & _strength
const Real _f_tol
Tolerance on yield function.
const bool _use_custom_cto
Whether to use the custom consistent tangent operator calculation.
const Real _shift
yield function is shifted by this amount to avoid problems with stress-derivatives at equal eigenvalu...
void mooseError(Args &&... args) const
const InputParameters & parameters() const
const ConsoleStream _console
const unsigned int _max_iters
maximum iterations allowed in the custom return-map algorithm
bool isParamValid(const std::string &name) const

Member Function Documentation

◆ activeConstraints()

void SolidMechanicsPlasticTensileMulti::activeConstraints ( const std::vector< Real > &  f,
const RankTwoTensor stress,
Real  intnl,
const RankFourTensor Eijkl,
std::vector< bool > &  act,
RankTwoTensor returned_stress 
) const
overridevirtual

The active yield surfaces, given a vector of yield functions.

This is used by FiniteStrainMultiPlasticity to determine the initial set of active constraints at the trial (stress, intnl) configuration. It is up to you (the coder) to determine how accurate you want the returned_stress to be. Currently it is only used by FiniteStrainMultiPlasticity to estimate a good starting value for the Newton-Rahson procedure, so currently it may not need to be super perfect.

Parameters
fvalues of the yield functions
stressstress tensor
intnlinternal parameter
Eijklelasticity tensor (stress = Eijkl*strain)
[out]actact[i] = true if the i_th yield function is active
[out]returned_stressApproximate value of the returned stress

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 163 of file SolidMechanicsPlasticTensileMulti.C.

169 {
170  act.assign(3, false);
171 
172  if (f[0] <= _f_tol && f[1] <= _f_tol && f[2] <= _f_tol)
173  {
174  returned_stress = stress;
175  return;
176  }
177 
178  Real returned_intnl;
179  std::vector<Real> dpm(3);
180  RankTwoTensor delta_dp;
181  std::vector<Real> yf(3);
182  bool trial_stress_inadmissible;
183  doReturnMap(stress,
184  intnl,
185  Eijkl,
186  0.0,
187  returned_stress,
188  returned_intnl,
189  dpm,
190  delta_dp,
191  yf,
192  trial_stress_inadmissible);
193 
194  for (unsigned i = 0; i < 3; ++i)
195  act[i] = (dpm[i] > 0);
196 }
Real f(Real x)
Test function for Brents method.
virtual bool doReturnMap(const RankTwoTensor &trial_stress, Real intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &returned_stress, Real &returned_intnl, std::vector< Real > &dpm, RankTwoTensor &delta_dp, std::vector< Real > &yf, bool &trial_stress_inadmissible) const
Just like returnMap, but a protected interface that definitely uses the algorithm, since returnMap itself does not use the algorithm if _use_returnMap=false.
const Real _f_tol
Tolerance on yield function.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ consistentTangentOperator()

RankFourTensor SolidMechanicsPlasticTensileMulti::consistentTangentOperator ( const RankTwoTensor trial_stress,
Real  intnl_old,
const RankTwoTensor stress,
Real  intnl,
const RankFourTensor E_ijkl,
const std::vector< Real > &  cumulative_pm 
) const
overridevirtual

Calculates a custom consistent tangent operator.

You may choose to over-ride this in your derived SolidMechanicsPlasticXXXX class.

(Note, if you over-ride returnMap, you will probably want to override consistentTangentOpertor too, otherwise it will default to E_ijkl.)

Parameters
stress_oldtrial stress before returning
intnl_oldinternal parameter before returning
stresscurrent returned stress state
intnlinternal parameter
E_ijklelasticity tensor
cumulative_pmthe cumulative plastic multipliers
Returns
the consistent tangent operator: E_ijkl if not over-ridden

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 617 of file SolidMechanicsPlasticTensileMulti.C.

624 {
625  if (!_use_custom_cto)
627  trial_stress, intnl_old, stress, intnl, E_ijkl, cumulative_pm);
628 
629  mooseAssert(cumulative_pm.size() == 3,
630  "SolidMechanicsPlasticTensileMulti size of cumulative_pm should be 3 but it is "
631  << cumulative_pm.size());
632 
633  if (cumulative_pm[2] <= 0) // All cumulative_pm are non-positive, so this is admissible
634  return E_ijkl;
635 
636  // Need the eigenvalues at the returned configuration
637  std::vector<Real> eigvals;
638  stress.symmetricEigenvalues(eigvals);
639 
640  // need to rotate to and from principal stress space
641  // using the eigenvectors of the trial configuration
642  // (not the returned configuration).
643  std::vector<Real> trial_eigvals;
644  RankTwoTensor trial_eigvecs;
645  trial_stress.symmetricEigenvaluesEigenvectors(trial_eigvals, trial_eigvecs);
646 
647  // The returnMap will have returned to the Tip, Edge or
648  // Plane. The consistentTangentOperator describes the
649  // change in stress for an arbitrary change in applied
650  // strain. I assume that the change in strain will not
651  // change the type of return (Tip remains Tip, Edge remains
652  // Edge, Plane remains Plane).
653  // I assume isotropic elasticity.
654  //
655  // The consistent tangent operator is a little different
656  // than cases where no rotation to principal stress space
657  // is made during the returnMap. Let S_ij be the stress
658  // in original coordinates, and s_ij be the stress in the
659  // principal stress coordinates, so that
660  // s_ij = diag(eigvals[0], eigvals[1], eigvals[2])
661  // We want dS_ij under an arbitrary change in strain (ep->ep+dep)
662  // dS = S(ep+dep) - S(ep)
663  // = R(ep+dep) s(ep+dep) R(ep+dep)^T - R(ep) s(ep) R(ep)^T
664  // Here R = the rotation to principal-stress space, ie
665  // R_ij = eigvecs[i][j] = i^th component of j^th eigenvector
666  // Expanding to first order in dep,
667  // dS = R(ep) (s(ep+dep) - s(ep)) R(ep)^T
668  // + dR/dep s(ep) R^T + R(ep) s(ep) dR^T/dep
669  // The first line is all that is usually calculated in the
670  // consistent tangent operator calculation, and is called
671  // cto below.
672  // The second line involves changes in the eigenvectors, and
673  // is called sec below.
674 
675  RankFourTensor cto;
676  const Real hard = dtensile_strength(intnl);
677  const Real la = E_ijkl(0, 0, 1, 1);
678  const Real mu = 0.5 * (E_ijkl(0, 0, 0, 0) - la);
679 
680  if (cumulative_pm[1] <= 0)
681  {
682  // only cumulative_pm[2] is positive, so this is return to the Plane
683  const Real denom = hard + la + 2 * mu;
684  const Real al = la * la / denom;
685  const Real be = la * (la + 2 * mu) / denom;
686  const Real ga = hard * (la + 2 * mu) / denom;
687  std::vector<Real> comps(9);
688  comps[0] = comps[4] = la + 2 * mu - al;
689  comps[1] = comps[3] = la - al;
690  comps[2] = comps[5] = comps[6] = comps[7] = la - be;
691  comps[8] = ga;
693  }
694  else if (cumulative_pm[0] <= 0)
695  {
696  // both cumulative_pm[2] and cumulative_pm[1] are positive, so Edge
697  const Real denom = 2 * hard + 2 * la + 2 * mu;
698  const Real al = hard * 2 * la / denom;
699  const Real be = hard * (2 * la + 2 * mu) / denom;
700  std::vector<Real> comps(9);
701  comps[0] = la + 2 * mu - 2 * la * la / denom;
702  comps[1] = comps[2] = al;
703  comps[3] = comps[6] = al;
704  comps[4] = comps[5] = comps[7] = comps[8] = be;
706  }
707  else
708  {
709  // all cumulative_pm are positive, so Tip
710  const Real denom = 3 * hard + 3 * la + 2 * mu;
711  std::vector<Real> comps(2);
712  comps[0] = hard * (3 * la + 2 * mu) / denom;
713  comps[1] = 0;
715  }
716 
717  cto.rotate(trial_eigvecs);
718 
719  // drdsig = change in eigenvectors under a small stress change
720  // drdsig(i,j,m,n) = dR(i,j)/dS_mn
721  // The formula below is fairly easily derived:
722  // S R = R s, so taking the variation
723  // dS R + S dR = dR s + R ds, and multiplying by R^T
724  // R^T dS R + R^T S dR = R^T dR s + ds .... (eqn 1)
725  // I demand that RR^T = 1 = R^T R, and also that
726  // (R+dR)(R+dR)^T = 1 = (R+dT)^T (R+dR), which means
727  // that dR = R*c, for some antisymmetric c, so Eqn1 reads
728  // R^T dS R + s c = c s + ds
729  // Grabbing the components of this gives ds/dS (already
730  // in RankTwoTensor), and c, which is:
731  // dR_ik/dS_mn = drdsig(i, k, m, n) = trial_eigvecs(m, b)*trial_eigvecs(n, k)*trial_eigvecs(i,
732  // b)/(trial_eigvals[k] - trial_eigvals[b]);
733  // (sum over b!=k).
734 
735  RankFourTensor drdsig;
736  for (unsigned k = 0; k < 3; ++k)
737  for (unsigned b = 0; b < 3; ++b)
738  {
739  if (b == k)
740  continue;
741  for (unsigned m = 0; m < 3; ++m)
742  for (unsigned n = 0; n < 3; ++n)
743  for (unsigned i = 0; i < 3; ++i)
744  drdsig(i, k, m, n) += trial_eigvecs(m, b) * trial_eigvecs(n, k) * trial_eigvecs(i, b) /
745  (trial_eigvals[k] - trial_eigvals[b]);
746  }
747 
748  // With diagla = diag(eigvals[0], eigvals[1], digvals[2])
749  // The following implements
750  // ans(i, j, a, b) += (drdsig(i, k, m, n)*trial_eigvecs(j, l)*diagla(k, l) + trial_eigvecs(i,
751  // k)*drdsig(j, l, m, n)*diagla(k, l))*E_ijkl(m, n, a, b);
752  // (sum over k, l, m and n)
753 
754  RankFourTensor ans;
755  for (unsigned i = 0; i < 3; ++i)
756  for (unsigned j = 0; j < 3; ++j)
757  for (unsigned a = 0; a < 3; ++a)
758  for (unsigned k = 0; k < 3; ++k)
759  for (unsigned m = 0; m < 3; ++m)
760  ans(i, j, a, a) += (drdsig(i, k, m, m) * trial_eigvecs(j, k) +
761  trial_eigvecs(i, k) * drdsig(j, k, m, m)) *
762  eigvals[k] * la; // E_ijkl(m, n, a, b) = la*(m==n)*(a==b);
763 
764  for (unsigned i = 0; i < 3; ++i)
765  for (unsigned j = 0; j < 3; ++j)
766  for (unsigned a = 0; a < 3; ++a)
767  for (unsigned b = 0; b < 3; ++b)
768  for (unsigned k = 0; k < 3; ++k)
769  {
770  ans(i, j, a, b) += (drdsig(i, k, a, b) * trial_eigvecs(j, k) +
771  trial_eigvecs(i, k) * drdsig(j, k, a, b)) *
772  eigvals[k] * mu; // E_ijkl(m, n, a, b) = mu*(m==a)*(n==b)
773  ans(i, j, a, b) += (drdsig(i, k, b, a) * trial_eigvecs(j, k) +
774  trial_eigvecs(i, k) * drdsig(j, k, b, a)) *
775  eigvals[k] * mu; // E_ijkl(m, n, a, b) = mu*(m==b)*(n==a)
776  }
777 
778  return cto + ans;
779 }
void symmetricEigenvalues(std::vector< Real > &eigvals) const
virtual Real dtensile_strength(const Real internal_param) const
d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param ...
virtual RankFourTensor consistentTangentOperator(const RankTwoTensor &trial_stress, Real intnl_old, const RankTwoTensor &stress, Real intnl, const RankFourTensor &E_ijkl, const std::vector< Real > &cumulative_pm) const
Calculates a custom consistent tangent operator.
void fillFromInputVector(const std::vector< T > &input, FillMethod fill_method)
static const std::string mu
Definition: NS.h:123
void rotate(const TypeTensor< T > &R)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const bool _use_custom_cto
Whether to use the custom consistent tangent operator calculation.
void symmetricEigenvaluesEigenvectors(std::vector< Real > &eigvals, RankTwoTensorTempl< Real > &eigvecs) const
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
static const std::string k
Definition: NS.h:130

◆ dflowPotential_dintnl()

RankTwoTensor SolidMechanicsPlasticModel::dflowPotential_dintnl ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The derivative of the flow potential with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
Returns
dr_dintnl(i, j) = dr(i, j)/dintnl

Reimplemented in SolidMechanicsPlasticMeanCapTC, SolidMechanicsPlasticDruckerPrager, SolidMechanicsPlasticJ2, SolidMechanicsPlasticWeakPlaneShear, SolidMechanicsPlasticWeakPlaneTensile, SolidMechanicsPlasticMohrCoulomb, SolidMechanicsPlasticTensile, SolidMechanicsPlasticMeanCap, SolidMechanicsPlasticSimpleTester, and SolidMechanicsPlasticWeakPlaneTensileN.

Definition at line 131 of file SolidMechanicsPlasticModel.C.

Referenced by SolidMechanicsPlasticModel::dflowPotential_dintnlV().

133 {
134  return RankTwoTensor();
135 }

◆ dflowPotential_dintnlV()

void SolidMechanicsPlasticTensileMulti::dflowPotential_dintnlV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankTwoTensor > &  dr_dintnl 
) const
overridevirtual

The derivative of the flow potential with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
[out]dr_dintnldr_dintnl[alpha](i, j) = dr[alpha](i, j)/dintnl

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 144 of file SolidMechanicsPlasticTensileMulti.C.

146 {
147  dr_dintnl.assign(3, RankTwoTensor());
148 }

◆ dflowPotential_dstress()

RankFourTensor SolidMechanicsPlasticModel::dflowPotential_dstress ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

◆ dflowPotential_dstressV()

void SolidMechanicsPlasticTensileMulti::dflowPotential_dstressV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankFourTensor > &  dr_dstress 
) const
overridevirtual

The derivative of the flow potential with respect to stress.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
[out]dr_dstressdr_dstress[alpha](i, j, k, l) = dr[alpha](i, j)/dstress(k, l)

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 137 of file SolidMechanicsPlasticTensileMulti.C.

139 {
140  stress.d2symmetricEigenvalues(dr_dstress);
141 }
void d2symmetricEigenvalues(std::vector< RankFourTensorTempl< Real >> &deriv) const

◆ dhardPotential_dintnl()

Real SolidMechanicsPlasticModel::dhardPotential_dintnl ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The derivative of the hardening potential with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlinternal parameter
Returns
the derivative

Reimplemented in SolidMechanicsPlasticMeanCapTC.

Definition at line 172 of file SolidMechanicsPlasticModel.C.

Referenced by SolidMechanicsPlasticModel::dhardPotential_dintnlV().

174 {
175  return 0.0;
176 }

◆ dhardPotential_dintnlV()

void SolidMechanicsPlasticModel::dhardPotential_dintnlV ( const RankTwoTensor stress,
Real  intnl,
std::vector< Real > &  dh_dintnl 
) const
virtualinherited

The derivative of the hardening potential with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlinternal parameter
[out]dh_dintnldh_dintnl[alpha] = dh[alpha]/dintnl

Definition at line 178 of file SolidMechanicsPlasticModel.C.

181 {
182  dh_dintnl.resize(numberSurfaces(), dhardPotential_dintnl(stress, intnl));
183 }
virtual Real dhardPotential_dintnl(const RankTwoTensor &stress, Real intnl) const
The derivative of the hardening potential with respect to the internal parameter. ...
virtual unsigned int numberSurfaces() const
The number of yield surfaces for this plasticity model.

◆ dhardPotential_dstress()

RankTwoTensor SolidMechanicsPlasticModel::dhardPotential_dstress ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The derivative of the hardening potential with respect to stress.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlinternal parameter
Returns
dh_dstress(i, j) = dh/dstress(i, j)

Reimplemented in SolidMechanicsPlasticMeanCapTC.

Definition at line 158 of file SolidMechanicsPlasticModel.C.

Referenced by SolidMechanicsPlasticModel::dhardPotential_dstressV().

160 {
161  return RankTwoTensor();
162 }

◆ dhardPotential_dstressV()

void SolidMechanicsPlasticModel::dhardPotential_dstressV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankTwoTensor > &  dh_dstress 
) const
virtualinherited

The derivative of the hardening potential with respect to stress.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlinternal parameter
[out]dh_dstressdh_dstress[alpha](i, j) = dh[alpha]/dstress(i, j)

Definition at line 164 of file SolidMechanicsPlasticModel.C.

167 {
168  dh_dstress.assign(numberSurfaces(), dhardPotential_dstress(stress, intnl));
169 }
virtual RankTwoTensor dhardPotential_dstress(const RankTwoTensor &stress, Real intnl) const
The derivative of the hardening potential with respect to stress.
virtual unsigned int numberSurfaces() const
The number of yield surfaces for this plasticity model.

◆ doReturnMap()

bool SolidMechanicsPlasticTensileMulti::doReturnMap ( const RankTwoTensor trial_stress,
Real  intnl_old,
const RankFourTensor E_ijkl,
Real  ep_plastic_tolerance,
RankTwoTensor returned_stress,
Real returned_intnl,
std::vector< Real > &  dpm,
RankTwoTensor delta_dp,
std::vector< Real > &  yf,
bool &  trial_stress_inadmissible 
) const
privatevirtual

Just like returnMap, but a protected interface that definitely uses the algorithm, since returnMap itself does not use the algorithm if _use_returnMap=false.

Definition at line 257 of file SolidMechanicsPlasticTensileMulti.C.

Referenced by activeConstraints(), and returnMap().

267 {
268  mooseAssert(dpm.size() == 3,
269  "SolidMechanicsPlasticTensileMulti size of dpm should be 3 but it is " << dpm.size());
270 
271  std::vector<Real> eigvals;
272  RankTwoTensor eigvecs;
273  trial_stress.symmetricEigenvaluesEigenvectors(eigvals, eigvecs);
274  eigvals[0] += _shift;
275  eigvals[2] -= _shift;
276 
277  Real str = tensile_strength(intnl_old);
278 
279  yf.resize(3);
280  yf[0] = eigvals[0] - str;
281  yf[1] = eigvals[1] - str;
282  yf[2] = eigvals[2] - str;
283 
284  if (yf[0] <= _f_tol && yf[1] <= _f_tol && yf[2] <= _f_tol)
285  {
286  // purely elastic (trial_stress, intnl_old)
287  trial_stress_inadmissible = false;
288  return true;
289  }
290 
291  trial_stress_inadmissible = true;
292  delta_dp.zero();
293  returned_stress.zero();
294 
295  // In the following i often assume that E_ijkl is
296  // for an isotropic situation. This reduces FLOPS
297  // substantially which is important since the returnMap
298  // is potentially the most compute-intensive function
299  // of a simulation.
300  // In many comments i write the general expression, and
301  // i hope that might guide future coders if they are
302  // generalising to a non-istropic E_ijkl
303 
304  // n[alpha] = E_ijkl*r[alpha]_kl expressed in principal stress space
305  // (alpha = 0, 1, 2, corresponding to the three surfaces)
306  // Note that in principal stress space, the flow
307  // directions are, expressed in 'vector' form,
308  // r[0] = (1,0,0), r[1] = (0,1,0), r[2] = (0,0,1).
309  // Similar for _n:
310  // so _n[0] = E_ij00*r[0], _n[1] = E_ij11*r[1], _n[2] = E_ij22*r[2]
311  // In the following I assume that the E_ijkl is
312  // for an isotropic situation.
313  // In the anisotropic situation, we couldn't express
314  // the flow directions as vectors in the same principal
315  // stress space as the stress: they'd be full rank-2 tensors
316  std::vector<RealVectorValue> n(3);
317  n[0](0) = E_ijkl(0, 0, 0, 0);
318  n[0](1) = E_ijkl(1, 1, 0, 0);
319  n[0](2) = E_ijkl(2, 2, 0, 0);
320  n[1](0) = E_ijkl(0, 0, 1, 1);
321  n[1](1) = E_ijkl(1, 1, 1, 1);
322  n[1](2) = E_ijkl(2, 2, 1, 1);
323  n[2](0) = E_ijkl(0, 0, 2, 2);
324  n[2](1) = E_ijkl(1, 1, 2, 2);
325  n[2](2) = E_ijkl(2, 2, 2, 2);
326 
327  // With non-zero Poisson's ratio and hardening
328  // it is not computationally cheap to know whether
329  // the trial stress will return to the tip, edge,
330  // or plane. The following is correct for zero
331  // Poisson's ratio and no hardening, and at least
332  // gives a not-completely-stupid guess in the
333  // more general case.
334  // trial_order[0] = type of return to try first
335  // trial_order[1] = type of return to try second
336  // trial_order[2] = type of return to try third
337  const unsigned int number_of_return_paths = 3;
338  std::vector<int> trial_order(number_of_return_paths);
339  if (yf[0] > _f_tol) // all the yield functions are positive, since eigvals are ordered eigvals[0]
340  // <= eigvals[1] <= eigvals[2]
341  {
342  trial_order[0] = tip;
343  trial_order[1] = edge;
344  trial_order[2] = plane;
345  }
346  else if (yf[1] > _f_tol) // two yield functions are positive
347  {
348  trial_order[0] = edge;
349  trial_order[1] = tip;
350  trial_order[2] = plane;
351  }
352  else
353  {
354  trial_order[0] = plane;
355  trial_order[1] = edge;
356  trial_order[2] = tip;
357  }
358 
359  unsigned trial;
360  bool nr_converged = false;
361  for (trial = 0; trial < number_of_return_paths; ++trial)
362  {
363  switch (trial_order[trial])
364  {
365  case tip:
366  nr_converged = returnTip(eigvals, n, dpm, returned_stress, intnl_old, 0);
367  break;
368  case edge:
369  nr_converged = returnEdge(eigvals, n, dpm, returned_stress, intnl_old, 0);
370  break;
371  case plane:
372  nr_converged = returnPlane(eigvals, n, dpm, returned_stress, intnl_old, 0);
373  break;
374  }
375 
376  str = tensile_strength(intnl_old + dpm[0] + dpm[1] + dpm[2]);
377 
378  if (nr_converged && KuhnTuckerOK(returned_stress, dpm, str, ep_plastic_tolerance))
379  break;
380  }
381 
382  if (trial == number_of_return_paths)
383  {
384  Moose::err << "Trial stress = \n";
385  trial_stress.print(Moose::err);
386  Moose::err << "Internal parameter = " << intnl_old << std::endl;
387  mooseError("SolidMechanicsPlasticTensileMulti: FAILURE! You probably need to implement a "
388  "line search\n");
389  // failure - must place yield function values at trial stress into yf
390  str = tensile_strength(intnl_old);
391  yf[0] = eigvals[0] - str;
392  yf[1] = eigvals[1] - str;
393  yf[2] = eigvals[2] - str;
394  return false;
395  }
396 
397  // success
398 
399  returned_intnl = intnl_old;
400  for (unsigned i = 0; i < 3; ++i)
401  {
402  yf[i] = returned_stress(i, i) - str;
403  delta_dp(i, i) = dpm[i];
404  returned_intnl += dpm[i];
405  }
406  returned_stress = eigvecs * returned_stress * (eigvecs.transpose());
407  delta_dp = eigvecs * delta_dp * (eigvecs.transpose());
408  return true;
409 }
void print(std::ostream &stm=Moose::out) const
bool returnPlane(const std::vector< Real > &eigvals, const std::vector< RealVectorValue > &n, std::vector< Real > &dpm, RankTwoTensor &returned_stress, Real intnl_old, Real initial_guess) const
Tries to return-map to the Tensile plane The return value is true if the internal Newton-Raphson proc...
bool KuhnTuckerOK(const RankTwoTensor &returned_diagonal_stress, const std::vector< Real > &dpm, Real str, Real ep_plastic_tolerance) const
Returns true if the Kuhn-Tucker conditions are satisfied.
bool returnTip(const std::vector< Real > &eigvals, const std::vector< RealVectorValue > &n, std::vector< Real > &dpm, RankTwoTensor &returned_stress, Real intnl_old, Real initial_guess) const
Tries to return-map to the Tensile tip.
const Real _f_tol
Tolerance on yield function.
RankTwoTensorTempl< Real > transpose() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void symmetricEigenvaluesEigenvectors(std::vector< Real > &eigvals, RankTwoTensorTempl< Real > &eigvecs) const
const Real _shift
yield function is shifted by this amount to avoid problems with stress-derivatives at equal eigenvalu...
bool returnEdge(const std::vector< Real > &eigvals, const std::vector< RealVectorValue > &n, std::vector< Real > &dpm, RankTwoTensor &returned_stress, Real intnl_old, Real initial_guess) const
Tries to return-map to the Tensile edge.
void mooseError(Args &&... args) const
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param

◆ dot()

Real SolidMechanicsPlasticTensileMulti::dot ( const std::vector< Real > &  a,
const std::vector< Real > &  b 
) const
private

dot product of two 3-dimensional vectors

Definition at line 199 of file SolidMechanicsPlasticTensileMulti.C.

201 {
202  return a[0] * b[0] + a[1] * b[1] + a[2] * b[2];
203 }

◆ dtensile_strength()

Real SolidMechanicsPlasticTensileMulti::dtensile_strength ( const Real  internal_param) const
protectedvirtual

d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param

Definition at line 157 of file SolidMechanicsPlasticTensileMulti.C.

Referenced by consistentTangentOperator(), dyieldFunction_dintnlV(), returnEdge(), returnPlane(), and returnTip().

158 {
159  return _strength.derivative(internal_param);
160 }
const SolidMechanicsHardeningModel & _strength
virtual Real derivative(Real intnl) const

◆ dyieldFunction_dintnl()

Real SolidMechanicsPlasticModel::dyieldFunction_dintnl ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The derivative of yield function with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
Returns
the derivative

Reimplemented in SolidMechanicsPlasticMeanCapTC, SolidMechanicsPlasticDruckerPrager, SolidMechanicsPlasticJ2, SolidMechanicsPlasticWeakPlaneShear, SolidMechanicsPlasticWeakPlaneTensile, SolidMechanicsPlasticMohrCoulomb, SolidMechanicsPlasticTensile, SolidMechanicsPlasticMeanCap, SolidMechanicsPlasticSimpleTester, and SolidMechanicsPlasticWeakPlaneTensileN.

Definition at line 90 of file SolidMechanicsPlasticModel.C.

Referenced by SolidMechanicsPlasticModel::dyieldFunction_dintnlV().

92 {
93  return 0.0;
94 }

◆ dyieldFunction_dintnlV()

void SolidMechanicsPlasticTensileMulti::dyieldFunction_dintnlV ( const RankTwoTensor stress,
Real  intnl,
std::vector< Real > &  df_dintnl 
) const
overridevirtual

The derivative of yield functions with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
[out]df_dintnldf_dintnl[alpha] = df[alpha]/dintnl

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 120 of file SolidMechanicsPlasticTensileMulti.C.

123 {
124  df_dintnl.assign(3, -dtensile_strength(intnl));
125 }
virtual Real dtensile_strength(const Real internal_param) const
d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param ...

◆ dyieldFunction_dstress()

RankTwoTensor SolidMechanicsPlasticModel::dyieldFunction_dstress ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

◆ dyieldFunction_dstressV()

void SolidMechanicsPlasticTensileMulti::dyieldFunction_dstressV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankTwoTensor > &  df_dstress 
) const
overridevirtual

The derivative of yield functions with respect to stress.

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
[out]df_dstressdf_dstress[alpha](i, j) = dyieldFunction[alpha]/dstress(i, j)

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 95 of file SolidMechanicsPlasticTensileMulti.C.

Referenced by flowPotentialV().

97 {
98  std::vector<Real> eigvals;
99  stress.dsymmetricEigenvalues(eigvals, df_dstress);
100 
101  if (eigvals[0] > eigvals[1] - 0.1 * _shift || eigvals[1] > eigvals[2] - 0.1 * _shift)
102  {
103  Real small_perturbation;
104  RankTwoTensor shifted_stress = stress;
105  while (eigvals[0] > eigvals[1] - 0.1 * _shift || eigvals[1] > eigvals[2] - 0.1 * _shift)
106  {
107  for (unsigned i = 0; i < 3; ++i)
108  for (unsigned j = 0; j <= i; ++j)
109  {
110  small_perturbation = 0.1 * _shift * 2 * (MooseRandom::rand() - 0.5);
111  shifted_stress(i, j) += small_perturbation;
112  shifted_stress(j, i) += small_perturbation;
113  }
114  shifted_stress.dsymmetricEigenvalues(eigvals, df_dstress);
115  }
116  }
117 }
void dsymmetricEigenvalues(std::vector< Real > &eigvals, std::vector< RankTwoTensorTempl< Real >> &deigvals) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real _shift
yield function is shifted by this amount to avoid problems with stress-derivatives at equal eigenvalu...
static Real rand()
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")

◆ execute()

void SolidMechanicsPlasticModel::execute ( )
virtualinherited

Implements GeneralUserObject.

Definition at line 45 of file SolidMechanicsPlasticModel.C.

46 {
47 }

◆ finalize()

void SolidMechanicsPlasticModel::finalize ( )
virtualinherited

Implements GeneralUserObject.

Definition at line 50 of file SolidMechanicsPlasticModel.C.

51 {
52 }

◆ flowPotential()

RankTwoTensor SolidMechanicsPlasticModel::flowPotential ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

◆ flowPotentialV()

void SolidMechanicsPlasticTensileMulti::flowPotentialV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankTwoTensor > &  r 
) const
overridevirtual

The flow potentials.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
[out]rr[alpha] is the flow potential for the "alpha" yield function

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 128 of file SolidMechanicsPlasticTensileMulti.C.

131 {
132  // This plasticity is associative so
133  dyieldFunction_dstressV(stress, intnl, r);
134 }
virtual void dyieldFunction_dstressV(const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &df_dstress) const override
The derivative of yield functions with respect to stress.

◆ hardPotential()

Real SolidMechanicsPlasticModel::hardPotential ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The hardening potential.

Parameters
stressthe stress at which to calculate the hardening potential
intnlinternal parameter
Returns
the hardening potential

Reimplemented in SolidMechanicsPlasticMeanCapTC.

Definition at line 145 of file SolidMechanicsPlasticModel.C.

Referenced by SolidMechanicsPlasticModel::hardPotentialV().

146 {
147  return -1.0;
148 }

◆ hardPotentialV()

void SolidMechanicsPlasticModel::hardPotentialV ( const RankTwoTensor stress,
Real  intnl,
std::vector< Real > &  h 
) const
virtualinherited

The hardening potential.

Parameters
stressthe stress at which to calculate the hardening potential
intnlinternal parameter
[out]hh[alpha] is the hardening potential for the "alpha" yield function

Definition at line 150 of file SolidMechanicsPlasticModel.C.

153 {
154  h.assign(numberSurfaces(), hardPotential(stress, intnl));
155 }
virtual Real hardPotential(const RankTwoTensor &stress, Real intnl) const
The hardening potential.
virtual unsigned int numberSurfaces() const
The number of yield surfaces for this plasticity model.

◆ initialize()

void SolidMechanicsPlasticModel::initialize ( )
virtualinherited

Implements GeneralUserObject.

Definition at line 40 of file SolidMechanicsPlasticModel.C.

41 {
42 }

◆ KuhnTuckerOK()

bool SolidMechanicsPlasticTensileMulti::KuhnTuckerOK ( const RankTwoTensor returned_diagonal_stress,
const std::vector< Real > &  dpm,
Real  str,
Real  ep_plastic_tolerance 
) const
private

Returns true if the Kuhn-Tucker conditions are satisfied.

Parameters
returned_diagonal_stressThe eigenvalues (sorted in ascending order as is standard in this Class) are stored in the diagonal components
dpmThe three plastic multipliers
strThe yield strength
ep_plastic_toleranceThe tolerance on the plastic strain (if dpm>-ep_plastic_tolerance then it is grouped as "non-negative" in the Kuhn-Tucker conditions).

Definition at line 604 of file SolidMechanicsPlasticTensileMulti.C.

Referenced by doReturnMap().

608 {
609  for (unsigned i = 0; i < 3; ++i)
611  returned_diagonal_stress(i, i) - str, dpm[i], ep_plastic_tolerance))
612  return false;
613  return true;
614 }
bool KuhnTuckerSingleSurface(Real yf, Real dpm, Real dpm_tol) const
Returns true if the Kuhn-Tucker conditions for the single surface are satisfied.

◆ KuhnTuckerSingleSurface()

bool SolidMechanicsPlasticModel::KuhnTuckerSingleSurface ( Real  yf,
Real  dpm,
Real  dpm_tol 
) const
inherited

Returns true if the Kuhn-Tucker conditions for the single surface are satisfied.

Parameters
yfYield function value
dpmplastic multiplier
dpm_toltolerance on plastic multiplier: viz dpm>-dpm_tol means "dpm is non-negative"

Definition at line 246 of file SolidMechanicsPlasticModel.C.

Referenced by SolidMechanicsPlasticMohrCoulombMulti::KuhnTuckerOK(), KuhnTuckerOK(), and SolidMechanicsPlasticModel::returnMap().

247 {
248  return (dpm == 0 && yf <= _f_tol) || (dpm > -dpm_tol && yf <= _f_tol && yf >= -_f_tol);
249 }
const Real _f_tol
Tolerance on yield function.

◆ modelName()

std::string SolidMechanicsPlasticTensileMulti::modelName ( ) const
overridevirtual

Implements SolidMechanicsPlasticModel.

Definition at line 215 of file SolidMechanicsPlasticTensileMulti.C.

216 {
217  return "TensileMulti";
218 }

◆ numberSurfaces()

unsigned int SolidMechanicsPlasticTensileMulti::numberSurfaces ( ) const
overridevirtual

The number of yield surfaces for this plasticity model.

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 74 of file SolidMechanicsPlasticTensileMulti.C.

75 {
76  return 3;
77 }

◆ returnEdge()

bool SolidMechanicsPlasticTensileMulti::returnEdge ( const std::vector< Real > &  eigvals,
const std::vector< RealVectorValue > &  n,
std::vector< Real > &  dpm,
RankTwoTensor returned_stress,
Real  intnl_old,
Real  initial_guess 
) const
private

Tries to return-map to the Tensile edge.

The return value is true if the internal Newton-Raphson process has converged, otherwise it is false

Parameters
eigvalsThe three stress eigenvalues, sorted in ascending order
nThe three return directions, n=E_ijkl*r. Note this algorithm assumes isotropic elasticity, so these are 3 vectors in principal stress space
dpmThe three plastic multipliers resulting from the return-map to the edge. This algorithm doesn't do Kuhn-Tucker checking, so these could be positive or negative or zero (dpm[0]=0 always for Edge return).
returned_stressThe returned stress. This will be diagonal, with the return-mapped eigenvalues in the diagonal positions, sorted in ascending order
intnl_oldThe internal parameter at stress=eigvals. This algorithm doesn't form the plastic strain, so you will have to use intnl=intnl_old+sum(dpm) if you need the new internal-parameter value at the returned point.
initial_guessA guess of dpm[1]+dpm[2]

Definition at line 489 of file SolidMechanicsPlasticTensileMulti.C.

Referenced by doReturnMap().

495 {
496  // work out the point to which we would return, "a". It is defined by
497  // f1 = 0 = f2, and the normality condition:
498  // (eigvals - a).(n1 x n2) = 0,
499  // where eigvals is the starting position
500  // (it is a vector in principal stress space).
501  // To get f1=0=f2, we need a = (a0, str, str), and a0 is found
502  // by expanding the normality condition to yield:
503  // a0 = (-str*n1xn2[1] - str*n1xn2[2] + edotn1xn2)/n1xn2[0];
504  // where edotn1xn2 = eigvals.(n1 x n2)
505  //
506  // We need to find the plastic multipliers, dpm, defined by
507  // eigvals - a = dpm[1]*n1 + dpm[2]*n2
508  // For the isotropic case, and defining eminusa = eigvals - a,
509  // the solution is easy:
510  // dpm[0] = 0;
511  // dpm[1] = (eminusa[1] - eminusa[0])/(n[1][1] - n[1][0]);
512  // dpm[2] = (eminusa[2] - eminusa[0])/(n[2][2] - n[2][0]);
513  //
514  // Now specialise to the isotropic case. Define
515  // x = dpm[1] + dpm[2] = (eigvals[1] + eigvals[2] - 2*str)/(n[0][0] + n[0][1])
516  // Notice that the RHS is a function of x, so we solve using
517  // Newton-Raphson starting with x=initial_guess
518  Real x = initial_guess;
519  const Real denom = n[0](0) + n[0](1);
520  Real str = tensile_strength(intnl_old + x);
521 
522  if (_strength.modelName().compare("Constant") != 0)
523  {
524  // Finding x is expensive. Therefore
525  // although x!=0 for Constant Hardening, solution
526  // for dpm above demonstrates that we don't
527  // actually need to know x to find the dpm for
528  // Constant Hardening.
529  //
530  // However, for nontrivial Hardening, the following
531  // is necessary
532  const Real eig = eigvals[1] + eigvals[2];
533  Real residual = denom * x - eig + 2 * str;
534  Real jacobian;
535  unsigned int iter = 0;
536  do
537  {
538  jacobian = denom + 2 * dtensile_strength(intnl_old + x);
539  x += -residual / jacobian;
540  if (iter > _max_iters) // not converging
541  return false;
542  str = tensile_strength(intnl_old + x);
543  residual = denom * x - eig + 2 * str;
544  iter++;
545  } while (residual * residual > _f_tol * _f_tol);
546  }
547 
548  dpm[0] = 0;
549  dpm[1] = ((eigvals[1] * n[0](0) - eigvals[2] * n[0](1)) / (n[0](0) - n[0](1)) - str) / denom;
550  dpm[2] = ((eigvals[2] * n[0](0) - eigvals[1] * n[0](1)) / (n[0](0) - n[0](1)) - str) / denom;
551 
552  returned_stress(0, 0) = eigvals[0] - n[0](1) * (dpm[1] + dpm[2]);
553  returned_stress(1, 1) = returned_stress(2, 2) = str;
554  return true;
555 }
virtual Real dtensile_strength(const Real internal_param) const
d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param ...
virtual std::string modelName() const =0
const std::vector< double > x
const SolidMechanicsHardeningModel & _strength
const Real _f_tol
Tolerance on yield function.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
const unsigned int _max_iters
maximum iterations allowed in the custom return-map algorithm

◆ returnMap()

bool SolidMechanicsPlasticTensileMulti::returnMap ( const RankTwoTensor trial_stress,
Real  intnl_old,
const RankFourTensor E_ijkl,
Real  ep_plastic_tolerance,
RankTwoTensor returned_stress,
Real returned_intnl,
std::vector< Real > &  dpm,
RankTwoTensor delta_dp,
std::vector< Real > &  yf,
bool &  trial_stress_inadmissible 
) const
overridevirtual

Performs a custom return-map.

You may choose to over-ride this in your derived SolidMechanicsPlasticXXXX class, and you may implement the return-map algorithm in any way that suits you. Eg, using a Newton-Raphson approach, or a radial-return, etc. This may also be used as a quick way of ascertaining whether (trial_stress, intnl_old) is in fact admissible.

For over-riding this function, please note the following.

(1) Denoting the return value of the function by "successful_return", the only possible output values should be: (A) trial_stress_inadmissible=false, successful_return=true. That is, (trial_stress, intnl_old) is in fact admissible (in the elastic domain). (B) trial_stress_inadmissible=true, successful_return=false. That is (trial_stress, intnl_old) is inadmissible (outside the yield surface), and you didn't return to the yield surface. (C) trial_stress_inadmissible=true, successful_return=true. That is (trial_stress, intnl_old) is inadmissible (outside the yield surface), but you did return to the yield surface. The default implementation only handles case (A) and (B): it does not attempt to do a return-map algorithm.

(2) you must correctly signal "successful_return" using the return value of this function. Don't assume the calling function will do Kuhn-Tucker checking and so forth!

(3) In cases (A) and (B) you needn't set returned_stress, returned_intnl, delta_dp, or dpm. This is for computational efficiency.

(4) In cases (A) and (B), you MUST place the yield function values at (trial_stress, intnl_old) into yf so the calling function can use this information optimally. You will have already calculated these yield function values, which can be quite expensive, and it's not very optimal for the calling function to have to re-calculate them.

(5) In case (C), you need to set: returned_stress (the returned value of stress) returned_intnl (the returned value of the internal variable) delta_dp (the change in plastic strain) dpm (the plastic multipliers needed to bring about the return) yf (yield function values at the returned configuration)

(Note, if you over-ride returnMap, you will probably want to override consistentTangentOpertor too, otherwise it will default to E_ijkl.)

Parameters
trial_stressThe trial stress
intnl_oldValue of the internal parameter
E_ijklElasticity tensor
ep_plastic_toleranceTolerance defined by the user for the plastic strain
[out]returned_stressIn case (C): lies on the yield surface after returning and produces the correct plastic strain (normality condition). Otherwise: not defined
[out]returned_intnlIn case (C): the value of the internal parameter after returning. Otherwise: not defined
[out]dpmIn case (C): the plastic multipliers needed to bring about the return. Otherwise: not defined
[out]delta_dpIn case (C): The change in plastic strain induced by the return process. Otherwise: not defined
[out]yfIn case (C): the yield function at (returned_stress, returned_intnl). Otherwise: the yield function at (trial_stress, intnl_old)
[out]trial_stress_inadmissibleShould be set to false if the trial_stress is admissible, and true if the trial_stress is inadmissible. This can be used by the calling prorgram
Returns
true if a successful return (or a return-map not needed), false if the trial_stress is inadmissible but the return process failed

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 221 of file SolidMechanicsPlasticTensileMulti.C.

231 {
233  return SolidMechanicsPlasticModel::returnMap(trial_stress,
234  intnl_old,
235  E_ijkl,
236  ep_plastic_tolerance,
237  returned_stress,
238  returned_intnl,
239  dpm,
240  delta_dp,
241  yf,
242  trial_stress_inadmissible);
243 
244  return doReturnMap(trial_stress,
245  intnl_old,
246  E_ijkl,
247  ep_plastic_tolerance,
248  returned_stress,
249  returned_intnl,
250  dpm,
251  delta_dp,
252  yf,
253  trial_stress_inadmissible);
254 }
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.
virtual bool returnMap(const RankTwoTensor &trial_stress, Real intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &returned_stress, Real &returned_intnl, std::vector< Real > &dpm, RankTwoTensor &delta_dp, std::vector< Real > &yf, bool &trial_stress_inadmissible) const
Performs a custom return-map.
virtual bool doReturnMap(const RankTwoTensor &trial_stress, Real intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &returned_stress, Real &returned_intnl, std::vector< Real > &dpm, RankTwoTensor &delta_dp, std::vector< Real > &yf, bool &trial_stress_inadmissible) const
Just like returnMap, but a protected interface that definitely uses the algorithm, since returnMap itself does not use the algorithm if _use_returnMap=false.

◆ returnPlane()

bool SolidMechanicsPlasticTensileMulti::returnPlane ( const std::vector< Real > &  eigvals,
const std::vector< RealVectorValue > &  n,
std::vector< Real > &  dpm,
RankTwoTensor returned_stress,
Real  intnl_old,
Real  initial_guess 
) const
private

Tries to return-map to the Tensile plane The return value is true if the internal Newton-Raphson process has converged, otherwise it is false.

Parameters
eigvalsThe three stress eigenvalues, sorted in ascending order
nThe three return directions, n=E_ijkl*r. Note this algorithm assumes isotropic elasticity, so these are 3 vectors in principal stress space
dpmThe three plastic multipliers resulting from the return-map to the plane. This algorithm doesn't do Kuhn-Tucker checking, so dpm[2] could be positive or negative or zero (dpm[0]=dpm[1]=0 always for Plane return).
returned_stressThe returned stress. This will be diagonal, with the return-mapped eigenvalues in the diagonal positions, sorted in ascending order
intnl_oldThe internal parameter at stress=eigvals. This algorithm doesn't form the plastic strain, so you will have to use intnl=intnl_old+sum(dpm) if you need the new internal-parameter value at the returned point.
initial_guessA guess of dpm[2]

Definition at line 558 of file SolidMechanicsPlasticTensileMulti.C.

Referenced by doReturnMap().

564 {
565  // the returned point, "a", is defined by f2=0 and
566  // a = p - dpm[2]*n2.
567  // This is a vector equation in
568  // principal stress space, and dpm[2] is the third
569  // plasticity multiplier (dpm[0]=0=dpm[1] for return
570  // to the plane) and "p" is the starting
571  // position (p=eigvals).
572  // (Kuhn-Tucker demands that dpm[2]>=0, but we leave checking
573  // that condition for later.)
574  // Since f2=0, we must have a[2]=tensile_strength,
575  // so we can just look at the [2] component of the
576  // equation, which yields
577  // n[2][2]*dpm[2] - eigvals[2] + str = 0
578  // For hardening, str=tensile_strength(intnl_old+dpm[2]),
579  // and we want to solve for dpm[2].
580  // Use Newton-Raphson with initial guess dpm[2] = initial_guess
581  dpm[2] = initial_guess;
582  Real residual = n[2](2) * dpm[2] - eigvals[2] + tensile_strength(intnl_old + dpm[2]);
583  Real jacobian;
584  unsigned int iter = 0;
585  do
586  {
587  jacobian = n[2](2) + dtensile_strength(intnl_old + dpm[2]);
588  dpm[2] += -residual / jacobian;
589  if (iter > _max_iters) // not converging
590  return false;
591  residual = n[2](2) * dpm[2] - eigvals[2] + tensile_strength(intnl_old + dpm[2]);
592  iter++;
593  } while (residual * residual > _f_tol * _f_tol);
594 
595  dpm[0] = 0;
596  dpm[1] = 0;
597  returned_stress(0, 0) = eigvals[0] - dpm[2] * n[2](0);
598  returned_stress(1, 1) = eigvals[1] - dpm[2] * n[2](1);
599  returned_stress(2, 2) = eigvals[2] - dpm[2] * n[2](2);
600  return true;
601 }
virtual Real dtensile_strength(const Real internal_param) const
d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param ...
const Real _f_tol
Tolerance on yield function.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
const unsigned int _max_iters
maximum iterations allowed in the custom return-map algorithm

◆ returnTip()

bool SolidMechanicsPlasticTensileMulti::returnTip ( const std::vector< Real > &  eigvals,
const std::vector< RealVectorValue > &  n,
std::vector< Real > &  dpm,
RankTwoTensor returned_stress,
Real  intnl_old,
Real  initial_guess 
) const
private

Tries to return-map to the Tensile tip.

The return value is true if the internal Newton-Raphson process has converged, otherwise it is false

Parameters
eigvalsThe three stress eigenvalues, sorted in ascending order
nThe three return directions, n=E_ijkl*r. Note this algorithm assumes isotropic elasticity, so these are 3 vectors in principal stress space
dpmThe three plastic multipliers resulting from the return-map to the tip. This algorithm doesn't do Kuhn-Tucker checking, so these could be positive or negative or zero
returned_stressThe returned stress. This will be diagonal, with the return-mapped eigenvalues in the diagonal positions, sorted in ascending order
intnl_oldThe internal parameter at stress=eigvals. This algorithm doesn't form the plastic strain, so you will have to use intnl=intnl_old+sum(dpm) if you need the new internal-parameter value at the returned point.
initial_guessA guess of dpm[0]+dpm[1]+dpm[2]

Definition at line 412 of file SolidMechanicsPlasticTensileMulti.C.

Referenced by doReturnMap().

418 {
419  // The returned point is defined by f0=f1=f2=0.
420  // that is, returned_stress = diag(str, str, str), where
421  // str = tensile_strength(intnl),
422  // where intnl = intnl_old + dpm[0] + dpm[1] + dpm[2]
423  // The 3 plastic multipliers, dpm, are defiend by the normality condition
424  // eigvals - str = dpm[0]*n[0] + dpm[1]*n[1] + dpm[2]*n[2]
425  // (Kuhn-Tucker demands that all dpm are non-negative, but we leave
426  // that checking for later.)
427  // This is a vector equation with solution (A):
428  // dpm[0] = triple(eigvals - str, n[1], n[2])/trip;
429  // dpm[1] = triple(eigvals - str, n[2], n[0])/trip;
430  // dpm[2] = triple(eigvals - str, n[0], n[1])/trip;
431  // where trip = triple(n[0], n[1], n[2]).
432  // By adding the three components of that solution together
433  // we can get an equation for x = dpm[0] + dpm[1] + dpm[2],
434  // and then our Newton-Raphson only involves one variable (x).
435  // In the following, i specialise to the isotropic situation.
436 
437  Real x = initial_guess;
438  const Real denom = (n[0](0) - n[0](1)) * (n[0](0) + 2 * n[0](1));
439  Real str = tensile_strength(intnl_old + x);
440 
441  if (_strength.modelName().compare("Constant") != 0)
442  {
443  // Finding x is expensive. Therefore
444  // although x!=0 for Constant Hardening, solution (A)
445  // demonstrates that we don't
446  // actually need to know x to find the dpm for
447  // Constant Hardening.
448  //
449  // However, for nontrivial Hardening, the following
450  // is necessary
451  const Real eig = eigvals[0] + eigvals[1] + eigvals[2];
452  const Real bul = (n[0](0) + 2 * n[0](1));
453 
454  // and finally, the equation we want to solve is:
455  // bul*x - eig + 3*str = 0
456  // where str=tensile_strength(intnl_old + x)
457  // and x = dpm[0] + dpm[1] + dpm[2]
458  // (Note this has units of stress, so using _f_tol as a convergence check is reasonable.)
459  // Use Netwon-Raphson with initial guess x = 0
460  Real residual = bul * x - eig + 3 * str;
461  Real jacobian;
462  unsigned int iter = 0;
463  do
464  {
465  jacobian = bul + 3 * dtensile_strength(intnl_old + x);
466  x += -residual / jacobian;
467  if (iter > _max_iters) // not converging
468  return false;
469  str = tensile_strength(intnl_old + x);
470  residual = bul * x - eig + 3 * str;
471  iter++;
472  } while (residual * residual > _f_tol * _f_tol);
473  }
474 
475  // The following is the solution (A) written above
476  // (dpm[0] = triple(eigvals - str, n[1], n[2])/trip, etc)
477  // in the isotropic situation
478  dpm[0] = (n[0](0) * (eigvals[0] - str) + n[0](1) * (eigvals[0] - eigvals[1] - eigvals[2] + str)) /
479  denom;
480  dpm[1] = (n[0](0) * (eigvals[1] - str) + n[0](1) * (eigvals[1] - eigvals[2] - eigvals[0] + str)) /
481  denom;
482  dpm[2] = (n[0](0) * (eigvals[2] - str) + n[0](1) * (eigvals[2] - eigvals[0] - eigvals[1] + str)) /
483  denom;
484  returned_stress(0, 0) = returned_stress(1, 1) = returned_stress(2, 2) = str;
485  return true;
486 }
virtual Real dtensile_strength(const Real internal_param) const
d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param ...
virtual std::string modelName() const =0
const std::vector< double > x
const SolidMechanicsHardeningModel & _strength
const Real _f_tol
Tolerance on yield function.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
const unsigned int _max_iters
maximum iterations allowed in the custom return-map algorithm

◆ tensile_strength()

Real SolidMechanicsPlasticTensileMulti::tensile_strength ( const Real  internal_param) const
protectedvirtual

tensile strength as a function of residual value, rate, and internal_param

Definition at line 151 of file SolidMechanicsPlasticTensileMulti.C.

Referenced by doReturnMap(), returnEdge(), returnPlane(), returnTip(), and yieldFunctionV().

152 {
153  return _strength.value(internal_param);
154 }
virtual Real value(Real intnl) const
const SolidMechanicsHardeningModel & _strength

◆ triple()

Real SolidMechanicsPlasticTensileMulti::triple ( const std::vector< Real > &  a,
const std::vector< Real > &  b,
const std::vector< Real > &  c 
) const
private

triple product of three 3-dimensional vectors

Definition at line 206 of file SolidMechanicsPlasticTensileMulti.C.

209 {
210  return a[0] * (b[1] * c[2] - b[2] * c[1]) - a[1] * (b[0] * c[2] - b[2] * c[0]) +
211  a[2] * (b[0] * c[1] - b[1] * c[0]);
212 }

◆ useCustomCTO()

bool SolidMechanicsPlasticTensileMulti::useCustomCTO ( ) const
overridevirtual

Returns false. You will want to override this in your derived class if you write a custom consistent tangent operator function.

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 788 of file SolidMechanicsPlasticTensileMulti.C.

789 {
790  return _use_custom_cto;
791 }
const bool _use_custom_cto
Whether to use the custom consistent tangent operator calculation.

◆ useCustomReturnMap()

bool SolidMechanicsPlasticTensileMulti::useCustomReturnMap ( ) const
overridevirtual

Returns false. You will want to override this in your derived class if you write a custom returnMap function.

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 782 of file SolidMechanicsPlasticTensileMulti.C.

783 {
784  return _use_custom_returnMap;
785 }
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.

◆ validParams()

InputParameters SolidMechanicsPlasticTensileMulti::validParams ( )
static

Definition at line 23 of file SolidMechanicsPlasticTensileMulti.C.

24 {
26  params.addClassDescription("Associative tensile plasticity with hardening/softening");
27  params.addRequiredParam<UserObjectName>(
28  "tensile_strength",
29  "A SolidMechanicsHardening UserObject that defines hardening of the tensile strength");
30  params.addParam<Real>("shift",
31  "Yield surface is shifted by this amount to avoid problems with "
32  "defining derivatives when eigenvalues are equal. If this is "
33  "larger than f_tol, a warning will be issued. Default = f_tol.");
34  params.addParam<unsigned int>("max_iterations",
35  10,
36  "Maximum number of Newton-Raphson iterations "
37  "allowed in the custom return-map algorithm. "
38  " For highly nonlinear hardening this may "
39  "need to be higher than 10.");
40  params.addParam<bool>("use_custom_returnMap",
41  true,
42  "Whether to use the custom returnMap "
43  "algorithm. Set to true if you are using "
44  "isotropic elasticity.");
45  params.addParam<bool>("use_custom_cto",
46  true,
47  "Whether to use the custom consistent tangent "
48  "operator computations. Set to true if you are "
49  "using isotropic elasticity.");
50  return params;
51 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void addRequiredParam(const std::string &name, const std::string &doc_string)
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)

◆ yieldFunction()

Real SolidMechanicsPlasticModel::yieldFunction ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The following functions are what you should override when building single-plasticity models.

The yield function

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
Returns
the yield function

Reimplemented in SolidMechanicsPlasticMeanCapTC, SolidMechanicsPlasticIsotropicSD, SolidMechanicsPlasticDruckerPrager, SolidMechanicsPlasticJ2, SolidMechanicsPlasticOrthotropic, SolidMechanicsPlasticWeakPlaneShear, SolidMechanicsPlasticWeakPlaneTensile, SolidMechanicsPlasticMohrCoulomb, SolidMechanicsPlasticTensile, SolidMechanicsPlasticDruckerPragerHyperbolic, SolidMechanicsPlasticMeanCap, SolidMechanicsPlasticSimpleTester, and SolidMechanicsPlasticWeakPlaneTensileN.

Definition at line 61 of file SolidMechanicsPlasticModel.C.

Referenced by SolidMechanicsPlasticModel::yieldFunctionV().

62 {
63  return 0.0;
64 }

◆ yieldFunctionV()

void SolidMechanicsPlasticTensileMulti::yieldFunctionV ( const RankTwoTensor stress,
Real  intnl,
std::vector< Real > &  f 
) const
overridevirtual

Calculates the yield functions.

Note that for single-surface plasticity you don't want to override this - override the private yieldFunction below

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
[out]fthe yield functions

Reimplemented from SolidMechanicsPlasticModel.

Definition at line 80 of file SolidMechanicsPlasticTensileMulti.C.

83 {
84  std::vector<Real> eigvals;
85  stress.symmetricEigenvalues(eigvals);
86  const Real str = tensile_strength(intnl);
87 
88  f.resize(3);
89  f[0] = eigvals[0] + _shift - str;
90  f[1] = eigvals[1] - str;
91  f[2] = eigvals[2] - _shift - str;
92 }
void symmetricEigenvalues(std::vector< Real > &eigvals) const
Real f(Real x)
Test function for Brents method.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real _shift
yield function is shifted by this amount to avoid problems with stress-derivatives at equal eigenvalu...
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param

Member Data Documentation

◆ _f_tol

const Real SolidMechanicsPlasticModel::_f_tol
inherited

◆ _ic_tol

const Real SolidMechanicsPlasticModel::_ic_tol
inherited

Tolerance on internal constraint.

Definition at line 173 of file SolidMechanicsPlasticModel.h.

◆ _max_iters

const unsigned int SolidMechanicsPlasticTensileMulti::_max_iters
private

maximum iterations allowed in the custom return-map algorithm

Definition at line 94 of file SolidMechanicsPlasticTensileMulti.h.

Referenced by returnEdge(), returnPlane(), and returnTip().

◆ _shift

const Real SolidMechanicsPlasticTensileMulti::_shift
private

yield function is shifted by this amount to avoid problems with stress-derivatives at equal eigenvalues

Definition at line 97 of file SolidMechanicsPlasticTensileMulti.h.

Referenced by doReturnMap(), dyieldFunction_dstressV(), SolidMechanicsPlasticTensileMulti(), and yieldFunctionV().

◆ _strength

const SolidMechanicsHardeningModel& SolidMechanicsPlasticTensileMulti::_strength
private

◆ _use_custom_cto

const bool SolidMechanicsPlasticTensileMulti::_use_custom_cto
private

Whether to use the custom consistent tangent operator calculation.

Definition at line 103 of file SolidMechanicsPlasticTensileMulti.h.

Referenced by consistentTangentOperator(), and useCustomCTO().

◆ _use_custom_returnMap

const bool SolidMechanicsPlasticTensileMulti::_use_custom_returnMap
private

Whether to use the custom return-map algorithm.

Definition at line 100 of file SolidMechanicsPlasticTensileMulti.h.

Referenced by returnMap(), and useCustomReturnMap().


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