https://mooseframework.inl.gov
AbaqusUMATStress.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
13 #include "DynamicLibraryLoader.h"
14 #include "ComputeFiniteStrain.h"
15 #include "RotationTensor.h"
17 
19 
24 {
25 public:
27 
29 
31  void initialSetup() override;
32 
34  void computeProperties() override;
35 
36 protected:
38  typedef void (*umat_t)(Real STRESS[],
39  Real STATEV[],
40  Real DDSDDE[],
41  Real * SSE,
42  Real * SPD,
43  Real * SCD,
44  Real * RPL,
45  Real DDSDDT[],
46  Real DRPLDE[],
47  Real * DRPLDT,
48  Real STRAN[],
49  Real DSTRAN[],
50  Real TIME[],
51  Real * DTIME,
52  Real * TEMP,
53  Real * DTEMP,
54  Real PREDEF[],
55  Real DPRED[],
56  char * CMNAME,
57  int * NDI,
58  int * NSHR,
59  int * NTENS,
60  int * NSTATV,
61  Real PROPS[],
62  int * NPROPS,
63  Real COORDS[],
64  Real DROT[],
65  Real * PNEWDT,
66  Real * CELENT,
67  Real DFGRD0[],
68  Real DFGRD1[],
69  int * NOEL,
70  unsigned int * NPT,
71  int * LAYER,
72  int * KSPT,
73  int * KSTEP,
74  int * KINC);
75 
76  // The plugin file name
77  FileName _plugin;
78 
79  // The plugin library wrapper
81 
82  // Function pointer to the dynamically loaded function
84 
91 
97 
100 
103 
106 
109 
112 
115 
117  char _aqCMNAME[80];
118 
120  int _aqNDI;
121 
123  int _aqNSHR;
124 
126  int _aqNTENS;
127 
130 
132  int _aqNOEL;
133 
135  unsigned int _aqNPT;
136 
138  int _aqLAYER;
139 
141  int _aqKSPT;
142 
144  int _aqKSTEP;
145 
148  int _aqKINC;
149 
150  // An array containing the solution-dependent state variables
151  std::vector<Real> _aqSTATEV;
152 
154  std::vector<Real> _aqDDSDDT;
155 
157  std::vector<Real> _aqDRPLDE;
158 
160  std::vector<Real> _aqSTRAN;
161 
163  std::vector<Real> _aqDFGRD0;
164 
166  std::vector<Real> _aqDFGRD1;
167 
169  std::vector<Real> _aqSTRESS;
170 
172  std::vector<Real> _aqDDSDDE;
173 
179  std::vector<Real> _aqDSTRAN;
180 
182  std::vector<Real> _aqPROPS;
183 
185  std::array<Real, 2> _aqTIME;
186 
188  std::array<Real, 3> _aqCOORDS;
189 
191  std::vector<Real> _aqDROT;
192 
195 
197  std::vector<Real> _aqPREDEF;
198 
200  std::vector<Real> _aqDPRED;
201 
202  void initQpStatefulProperties() override;
203  void computeQpStress() override;
204 
208 
211 
214 
217 
221 
224 
225  // Time step rotation increment
228 
229  // Coupled temperature field
232 
233  // Coupled user-defined field
234  const std::vector<const VariableValue *> _external_fields;
235  const std::vector<const VariableValue *> _external_fields_old;
236 
237  // External field names
238  std::vector<VariableName> _external_field_names;
239 
240  // Number of external fields provided by the user
241  const std::size_t _number_external_fields;
242 
243  const std::vector<MaterialPropertyName> _external_property_names;
244  const std::size_t _number_external_properties;
245  std::vector<const MaterialProperty<Real> *> _external_properties;
246  std::vector<const MaterialProperty<Real> *> _external_properties_old;
247 
250 
252  const bool _use_orientation;
256 
257 private:
260 
263 };
const MaterialProperty< std::vector< Real > > & _state_var_old
std::vector< Real > _aqSTATEV
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
const std::vector< MaterialPropertyName > _external_property_names
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.
std::vector< Real > _aqDPRED
Array of increments of predefined field variables.
const MaterialProperty< RankTwoTensor > & _total_rotation_old
DynamicLibraryLoader _library
Interface class for step user object.
std::array< Real, 3 > _aqCOORDS
An array containing the coordinates of this point.
MaterialProperty< std::vector< Real > > & _state_var
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
std::vector< Real > _aqDFGRD0
Array containing the deformation gradient at the beginning of the increment.
const std::size_t _number_external_fields
AbaqusUMATStress(const InputParameters &parameters)
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 ...
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
int _aqLAYER
Layer number (for composite shells and layered solids). (not supported)
std::vector< const MaterialProperty< Real > * > _external_properties
const RotationTensor _R
std::vector< VariableName > _external_field_names
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.
void initialSetup() override
check optional material properties for consistency
std::vector< Real > _aqDRPLDE
Variation of RPL with respect to the strain increments.
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
This is a RealTensor version of a rotation matrix It is instantiated with the Euler angles...
void(* umat_t)(Real STRESS[], Real STATEV[], Real DDSDDE[], Real *SSE, Real *SPD, Real *SCD, Real *RPL, Real DDSDDT[], Real DRPLDE[], Real *DRPLDT, Real STRAN[], Real DSTRAN[], Real TIME[], Real *DTIME, Real *TEMP, Real *DTEMP, Real PREDEF[], Real DPRED[], char *CMNAME, int *NDI, int *NSHR, int *NTENS, int *NSTATV, Real PROPS[], int *NPROPS, Real COORDS[], Real DROT[], Real *PNEWDT, Real *CELENT, Real DFGRD0[], Real DFGRD1[], int *NOEL, unsigned int *NPT, int *LAYER, int *KSPT, int *KSTEP, int *KINC)
function type for the external UMAT function
OutputTools< Real >::VariableValue VariableValue
static InputParameters validParams()
int _aqNOEL
Element number.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Real _aqSSE
specific elastic strain energy
const bool _use_orientation
Rotation information.
MaterialProperty< Real > & _creep_dissipation
Real _aqDTEMP
Increment of temperature.
class infix_ostream_iterator if void
MaterialProperty< Real > & _elastic_strain_energy
const InputParameters & parameters() const
Real _aqSCD
creep dissipation
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
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...
User object that provides analysis steps given user input.
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
Wrapper class to facilitate loading and lifetime management of dynamic libraries and obtaining pointe...
Real _aqSPD
plastic dissipation
const MaterialProperty< RankTwoTensor > & _stress_old
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.