16 #include "libmesh/int_range.h" 20 #define QUOTE(macro) stringifyName(macro) 30 "plugin",
"The path to the compiled dynamic library for the plugin you want to use");
32 "use_one_based_indexing",
33 "Parameter to control whether indexing for element and integration points as presented to " 34 "UMAT models is based on 1 (true) or 0 (false). This does not affect internal MOOSE " 35 "numbering. The option to use 0-based numbering is deprecated and will be removed soon.");
37 "constant_properties",
"Constant mechanical and thermal material properties (PROPS)");
39 "The number of state variables this UMAT is going to use");
40 params.
addCoupledVar(
"temperature", 0.0,
"Coupled temperature");
42 "The external fields that can be used in the UMAT subroutine");
43 params.
addParam<std::vector<MaterialPropertyName>>(
"external_properties", {},
"");
46 "Method to calculate the strain kinematics.");
50 "Whether or not this object should use the " 51 "displaced mesh for computing displacements and quantities based on the deformed state.");
53 "analysis_step_user_object",
54 "The AnalysisStepUserObject that provides times from simulation loading steps.");
57 "Euler angles that describe the orientation of the local material coordinate system.");
62 #error "The METHOD preprocessor symbol must be supplied by the build system." 67 _plugin(getParam<FileName>(
"plugin")),
68 _library(_plugin +
std::string(
"-") + QUOTE(METHOD) +
".plugin"),
69 _umat(_library.getFunction<umat_t>(
"umat_")),
70 _aqNSTATV(getParam<unsigned
int>(
"num_state_vars")),
72 _aqPROPS(getParam<
std::vector<
Real>>(
"constant_properties")),
73 _aqNPROPS(_aqPROPS.size()),
74 _stress_old(getMaterialPropertyOld<
RankTwoTensor>(_base_name +
"stress")),
75 _total_strain_old(getMaterialPropertyOld<
RankTwoTensor>(_base_name +
"total_strain")),
76 _strain_increment(getOptionalMaterialProperty<
RankTwoTensor>(_base_name +
"strain_increment")),
77 _jacobian_mult(declareProperty<
RankFourTensor>(_base_name +
"Jacobian_mult")),
78 _Fbar(getOptionalMaterialProperty<
RankTwoTensor>(_base_name +
"deformation_gradient")),
79 _Fbar_old(getOptionalMaterialPropertyOld<
RankTwoTensor>(_base_name +
"deformation_gradient")),
80 _state_var(declareProperty<
std::vector<
Real>>(_base_name +
"state_var")),
81 _state_var_old(getMaterialPropertyOld<
std::vector<
Real>>(_base_name +
"state_var")),
82 _elastic_strain_energy(declareProperty<
Real>(_base_name +
"elastic_strain_energy")),
83 _plastic_dissipation(declareProperty<
Real>(_base_name +
"plastic_dissipation")),
84 _creep_dissipation(declareProperty<
Real>(_base_name +
"creep_dissipation")),
85 _material_timestep(declareProperty<
Real>(_base_name +
"material_timestep_limit")),
87 getOptionalMaterialProperty<
RankTwoTensor>(_base_name +
"rotation_increment")),
88 _rotation_increment_old(
89 getOptionalMaterialPropertyOld<
RankTwoTensor>(_base_name +
"rotation_increment")),
90 _temperature(coupledValue(
"temperature")),
91 _temperature_old(coupledValueOld(
"temperature")),
92 _external_fields(isCoupled(
"external_fields") ? coupledValues(
"external_fields")
94 _external_fields_old(isCoupled(
"external_fields") ? coupledValuesOld(
"external_fields")
96 _number_external_fields(_external_fields.size()),
97 _external_property_names(getParam<std::vector<MaterialPropertyName>>(
"external_properties")),
98 _number_external_properties(_external_property_names.size()),
99 _external_properties(_number_external_properties),
100 _external_properties_old(_number_external_properties),
101 _use_one_based_indexing(getParam<bool>(
"use_one_based_indexing")),
102 _use_orientation(isParamValid(
"orientation")),
103 _R(_use_orientation ? getParam<RealVectorValue>(
"orientation") :
RealVectorValue(0.0)),
104 _total_rotation(declareProperty<RankTwoTensor>(
"total_rotation")),
105 _total_rotation_old(getMaterialPropertyOld<RankTwoTensor>(
"total_rotation")),
106 _decomposition_method(
107 getParam<MooseEnum>(
"decomposition_method").getEnum<ComputeFiniteStrain::DecompMethod>())
109 if (!_use_one_based_indexing)
111 "AbaqusUMATStress has transitioned to 1-based indexing in the element (NOEL) and " 112 "integration point (NPT) numbers to ensure maximum compatibility with legacy UMAT files. " 113 "Please ensure that any new UMAT plugins using these quantities are using the correct " 114 "indexing. 0-based indexing will be deprecated soon.");
117 for (std::size_t i = 0; i < _number_external_properties; ++i)
119 _external_properties[i] = &getMaterialProperty<Real>(_external_property_names[i]);
120 _external_properties_old[i] = &getMaterialPropertyOld<Real>(_external_property_names[i]);
128 _aqDDSDDT.resize(_aqNTENS);
129 _aqDRPLDE.resize(_aqNTENS);
130 _aqSTRAN.resize(_aqNTENS);
134 _aqSTRESS.resize(_aqNTENS);
135 _aqDDSDDE.resize(_aqNTENS * _aqNTENS);
136 _aqDSTRAN.resize(_aqNTENS);
137 _aqPREDEF.resize(_number_external_fields + _number_external_properties);
138 _aqDPRED.resize(_number_external_fields + _number_external_properties);
150 "': Incremental strain quantities are not available. You likely are using a total " 151 "strain formulation. Specify `incremental = true` in the tensor mechanics action, " 152 "or use ComputeIncrementalStrain in your input file.");
159 _step_user_object = &getUserObject<AnalysisStepUserObject>(
"analysis_step_user_object");
233 const Real * myDFGRD0 = &(FBar_old_fortran(0, 0));
234 const Real * myDFGRD1 = &(FBar_fortran(0, 0));
235 const Real * myDROT = &(DROT_fortran(0, 0));
273 static const std::array<Real, 6> strain_factor{{1, 1, 1, 2, 2, 2}};
276 static const std::array<std::pair<unsigned int, unsigned int>, 6>
component{
277 {{0, 0}, {1, 1}, {2, 2}, {0, 1}, {0, 2}, {1, 2}}};
284 _aqSTRAN[i] = total_strain_old(
a,
b) * strain_factor[i];
285 _aqDSTRAN[i] = strain_increment(
a,
b) * strain_factor[i];
297 _aqPNEWDT = std::numeric_limits<Real>::max();
393 const unsigned int ntens =
N * (
N + 1) / 2;
394 const int nskip =
N - 1;
408 :
_aqDDSDDE[(nskip + i +
j) * ntens +
k + nskip + l];
const MooseArray< Point > & _q_point
const MaterialProperty< std::vector< Real > > & _state_var_old
std::vector< Real > _aqSTATEV
FEProblemBase & _fe_problem
Real _aqDTIME
Time increment.
const OptionalMaterialProperty< RankTwoTensor > & _rotation_increment_old
ComputeGeneralStressBase is the direct base class for stress calculator materials that may leverage q...
const VariableValue & _temperature
std::vector< Real > _aqDSTRAN
Array of strain increments.
MaterialProperty< Real > & _material_timestep
recommended maximum timestep for this model under the current conditions
const std::vector< const VariableValue * > _external_fields_old
std::vector< Real > _aqSTRESS
Stress tensor (in: old stress, out: updated stress)
std::vector< const MaterialProperty< Real > * > _external_properties_old
Real _aqCELENT
Characteristic element length, which is a typical length of a line across an element for a first-orde...
MaterialProperty< RankFourTensor > & _jacobian_mult
Jacobian multiplier.
static const std::string component
std::vector< Real > _aqDPRED
Array of increments of predefined field variables.
const MaterialProperty< RankTwoTensor > & _total_rotation_old
std::array< Real, 3 > _aqCOORDS
An array containing the coordinates of this point.
virtual void computeProperties() override
MaterialProperty< std::vector< Real > > & _state_var
static constexpr std::size_t dim
std::vector< Real > _aqSTRAN
An array containing the total strains at the beginning of the increment.
int _aqKSTEP
The step number (as per Abaqus definition) can be set by the user.
std::vector< Real > _aqDDSDDE
Jacobian matrix of the model (out)
int _aqNDI
Number of direct stress components at this point.
int _aqNTENS
Size of the stress or strain component array (NDI + NSHR).
std::vector< Real > _aqPROPS
User-specified array of material constants associated with this user material.
std::vector< Real > _aqDDSDDT
Variation of the stress increments with respect to the temperature.
const std::size_t _number_external_properties
virtual const std::string & name() const
static RankTwoTensorTempl Identity()
std::vector< Real > _aqDFGRD0
Array containing the deformation gradient at the beginning of the increment.
const std::size_t _number_external_fields
virtual void resize(const std::size_t size) override final
void rotate(const RankTwoTensorTempl< Real > &R)
AbaqusUMATStress(const InputParameters ¶meters)
int _aqKINC
Increment number (_t_step).
Real _aqRPL
Volumetric heat generation per unit time at the end of the increment caused by mechanical working of ...
Real getStartTime(const unsigned int &step) const
const OptionalMaterialProperty< RankTwoTensor > & _rotation_increment
const ComputeFiniteStrain::DecompMethod _decomposition_method
Method being used to compute strain and rotation increments.
const MaterialProperty< RankTwoTensor > & _total_strain_old
virtual void initQpStatefulProperties() override
static InputParameters validParams()
int _aqLAYER
Layer number (for composite shells and layered solids). (not supported)
std::vector< const MaterialProperty< Real > * > _external_properties
int _aqKSPT
Section point number within the current layer. (not supported)
Real _aqDRPLDT
Variation of the volumetric heat generation (RPL) with respect to the temperature.
virtual void getAnalysisStepUserObject(const FEProblemBase &fe_problem, const AnalysisStepUserObject *&step_user_object, const std::string &name)
void initialSetup() override
check optional material properties for consistency
std::vector< Real > _aqDRPLDE
Variation of RPL with respect to the strain increments.
static MooseEnum decompositionType()
void mooseDeprecated(Args &&... args)
Real _aqPNEWDT
Ratio of suggested new time increment to the time increment being used (out)
int _aqNSHR
Number of engineering shear stress components at this point.
std::vector< Real > _aqDFGRD1
Array containing the deformation gradient at the end of the increment.
const bool _use_one_based_indexing
parameter to assist with the transition to 1-based indexing
const OptionalMaterialProperty< RankTwoTensor > & _strain_increment
unsigned int _aqNPT
Integration point number.
Real _aqTEMP
Temperature at the start of the increment.
Coupling material to use Abaqus UMAT models in MOOSE.
const std::vector< const VariableValue * > _external_fields
OutputTools< Real >::VariableValue VariableValue
static InputParameters validParams()
int _aqNOEL
Element number.
bool isParamSetByUser(const std::string &nm) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
registerMooseObject("SolidMechanicsApp", AbaqusUMATStress)
const bool _use_orientation
Rotation information.
MaterialProperty< Real > & _creep_dissipation
Real _aqDTEMP
Increment of temperature.
IntRange< T > make_range(T beg, T end)
void mooseError(Args &&... args) const
MaterialProperty< Real > & _elastic_strain_energy
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
std::array< Real, 2 > _aqTIME
Value of step time at the beginning of the current increment, total time at the beginning of the curr...
void initQpStatefulProperties() override
void computeProperties() override
perform per-element computation/initialization
std::vector< Real > _aqPREDEF
Array of interpolated values of predefined field variables at this point at the start of the incremen...
MaterialProperty< Real > & _plastic_dissipation
unsigned int getStep(const Real &time) const
int _aqNPROPS
User-defined number of material constants associated with this user material.
void computeQpStress() override
Compute the stress and store it in the _stress material property for the current quadrature point...
const OptionalMaterialProperty< RankTwoTensor > & _Fbar_old
const VariableValue & _temperature_old
int _aqNSTATV
Number of solution-dependent state variables that are associated with this material type...
MaterialProperty< RankTwoTensor > & _total_rotation
MooseUnits pow(const MooseUnits &, int)
static const std::string k
const MaterialProperty< RankTwoTensor > & _stress_old
void ErrorVector unsigned int
const Elem *const & _current_elem
const OptionalMaterialProperty< RankTwoTensor > & _Fbar
const AnalysisStepUserObject * _step_user_object
User object that determines step number.
char _aqCMNAME[80]
Model name buffer.
std::vector< Real > _aqDROT
Rotation increment matrix.