https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | List of all members
UpdatedLagrangianStressDivergenceBase< G > Class Template Reference

Enforce equilibrium with an updated Lagrangian formulation. More...

#include <UpdatedLagrangianStressDivergence.h>

Inheritance diagram for UpdatedLagrangianStressDivergenceBase< G >:
[legend]

Public Types

typedef std::vector< intJvarMap
 

Public Member Functions

 UpdatedLagrangianStressDivergenceBase (const InputParameters &parameters)
 
virtual void initialSetup () override
 
template<>
InputParameters validParams ()
 
template<>
void initialSetup ()
 
virtual void computeOffDiagJacobian (unsigned int jvar) override
 
unsigned int mapJvarToCvar (unsigned int jvar)
 
int mapJvarToCvar (unsigned int jvar, const JvarMap &jvar_map)
 
bool mapJvarToCvar (unsigned int jvar, unsigned int &cvar)
 
const JvarMapgetJvarMap ()
 
const JvarMapgetParameterJvarMap (std::string parameter_name)
 

Static Public Member Functions

static InputParameters baseParams ()
 
static InputParameters validParams ()
 

Protected Member Functions

virtual RankTwoTensor gradTest (unsigned int component) override
 
virtual RankTwoTensor gradTrial (unsigned int component) override
 
virtual void precalculateJacobianDisplacement (unsigned int component) override
 Prepare the average shape function gradients for stabilization. More...
 
virtual Real computeQpResidual () override
 
virtual Real computeQpJacobianDisplacement (unsigned int alpha, unsigned int beta) override
 
virtual Real computeQpJacobianTemperature (unsigned int cvar) override
 
virtual Real computeQpJacobianOutOfPlaneStrain () override
 
virtual void precalculateJacobian () override
 
virtual void precalculateOffDiagJacobian (unsigned int jvar) override
 
virtual Real computeQpJacobian () override
 
virtual Real computeQpOffDiagJacobian (unsigned int jvar) override
 

Protected Attributes

const MaterialProperty< RankTwoTensor > & _stress
 The Cauchy stress. More...
 
const MaterialProperty< RankFourTensor > & _material_jacobian
 
const bool _large_kinematics
 If true use large deformation kinematics. More...
 
const bool _stabilize_strain
 If true calculate the deformation gradient derivatives for F_bar. More...
 
const std::string _base_name
 Prepend to the material properties. More...
 
const unsigned int _alpha
 Which component of the vector residual this kernel is responsible for. More...
 
const unsigned int _ndisp
 Total number of displacements/size of residual vector. More...
 
std::vector< unsigned int_disp_nums
 The displacement numbers. More...
 
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
 
const MaterialProperty< RankTwoTensor > & _F_ust
 The unmodified deformation gradient. More...
 
const MaterialProperty< RankTwoTensor > & _F_avg
 The element-average deformation gradient. More...
 
const MaterialProperty< RankTwoTensor > & _f_inv
 The inverse increment deformation gradient. More...
 
const MaterialProperty< RankTwoTensor > & _F_inv
 The inverse deformation gradient. More...
 
const MaterialProperty< RankTwoTensor > & _F
 The actual (stabilized) deformation gradient. More...
 
const MooseVariable_temperature
 Temperature, if provided. This is used only to get the trial functions. More...
 
const MooseVariable_out_of_plane_strain
 Out-of-plane strain, if provided. More...
 
std::vector< std::vector< const MaterialProperty< RankTwoTensor > * > > _deigenstrain_dargs
 Eigenstrain derivatives wrt generate coupleds. More...
 
const unsigned int _n_args
 
Assembly_assembly_undisplaced
 
const VariablePhiGradient_grad_phi_undisplaced
 
const MooseArray< Real > & _JxW_undisplaced
 
const MooseArray< Real > & _coord_undisplaced
 
const MooseArray< Point > & _q_point_undisplaced
 

Private Member Functions

virtual RankTwoTensor gradTrialUnstabilized (unsigned int component)
 The unstabilized trial function gradient. More...
 
virtual RankTwoTensor gradTrialStabilized (unsigned int component)
 The stabilized trial function gradient. More...
 

Detailed Description

template<class G>
class UpdatedLagrangianStressDivergenceBase< G >

Enforce equilibrium with an updated Lagrangian formulation.

This class enforces equilibrium when used in conjunction with the corresponding strain calculator (CalculateStrainLagrangianKernel) and with either a stress calculator that provides the Cauchy stress ("stress") and the appropriate "cauchy_jacobian", which needs to be the derivative of the increment in Cauchy stress with respect to the increment in the spatial velocity gradient.

This kernel should be used with the new "ComputeLagrangianStressBase" stress update system and the "ComputeLagrangianStrain" system for strains.

use_displaced_mesh must be true for large deformation kinematics The kernel enforces this with an error

Definition at line 32 of file UpdatedLagrangianStressDivergence.h.

Constructor & Destructor Documentation

◆ UpdatedLagrangianStressDivergenceBase()

Definition at line 15 of file UpdatedLagrangianStressDivergence.C.

17  : LagrangianStressDivergenceBase(parameters),
18  _stress(getMaterialPropertyByName<RankTwoTensor>(_base_name + "cauchy_stress")),
19  _material_jacobian(getMaterialPropertyByName<RankFourTensor>(_base_name + "cauchy_jacobian")),
20 
21  // Assembly quantities in the reference frame for stabilization
22  _assembly_undisplaced(_fe_problem.assembly(_tid, this->_sys.number())),
27 {
28  // This kernel requires used_displaced_mesh to be true if large kinematics
29  // is on
30  if (_large_kinematics && (!getParam<bool>("use_displaced_mesh")))
31  mooseError("The UpdatedLagrangianStressDivergence kernels requires "
32  "used_displaced_mesh = true for large_kinematics = true");
33 
34  // Similarly, if large kinematics is off so should use_displaced_mesh
35  if (!_large_kinematics && (getParam<bool>("use_displaced_mesh")))
36  mooseError("The UpdatedLagrangianStressDivergence kernels requires "
37  "used_displaced_mesh = false for large_kinematics = false");
38 
39  // TODO: add weak plane stress support
41  mooseError("The UpdatedLagrangianStressDivergence kernels do not yet support the weak plane "
42  "stress formulation. Please use the TotalLagrangianStressDivergecen kernels.");
43 }
void mooseError(Args &&... args)
LagrangianStressDivergenceBase(const InputParameters &parameters)
const VariablePhiGradient & gradPhi() const
const MooseArray< Point > & qPoints() const
const MooseVariable * _out_of_plane_strain
Out-of-plane strain, if provided.
const bool _large_kinematics
If true use large deformation kinematics.
const MooseArray< Real > & coordTransformation() const
const MooseArray< Real > & JxW() const
const MaterialProperty< RankTwoTensor > & _stress
The Cauchy stress.
const MaterialProperty< RankFourTensor > & _material_jacobian
const std::string _base_name
Prepend to the material properties.

Member Function Documentation

◆ baseParams()

template<class G >
static InputParameters UpdatedLagrangianStressDivergenceBase< G >::baseParams ( )
inlinestatic

◆ computeQpJacobian()

Real LagrangianStressDivergenceBase::computeQpJacobian ( )
overrideprotectedvirtualinherited

Reimplemented in TotalLagrangianWeakPlaneStress.

Definition at line 115 of file LagrangianStressDivergenceBase.C.

116 {
118 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
virtual Real computeQpJacobianDisplacement(unsigned int alpha, unsigned int beta)=0

◆ computeQpJacobianDisplacement()

template<class G >
Real UpdatedLagrangianStressDivergenceBase< G >::computeQpJacobianDisplacement ( unsigned int  alpha,
unsigned int  beta 
)
overrideprotectedvirtual

Implements LagrangianStressDivergenceBase.

Definition at line 112 of file UpdatedLagrangianStressDivergence.C.

114 {
115  const auto grad_test = gradTest(alpha);
116  const auto grad_trial = gradTrial(beta);
117 
118  // J^{alpha beta} = J^{alpha beta}_material + J^{alpha beta}_geometric
119  // J^{alpha beta}_material = phi^alpha_{i, j} T_{ijkl} f^{-1}_{km} g^beta_{ml}
120  // J^{alpha beta}_geometric = sigma_{ij} (phi^alpha_{k, k} psi^beta_{i, j} -
121  // phi^alpha_{k, j} psi^beta_{i, k})
122 
123  // The material jacobian
124  Real J = grad_test.doubleContraction(_material_jacobian[_qp] * (_f_inv[_qp] * grad_trial));
125 
126  // The geometric jacobian
127  if (_large_kinematics)
128  {
129  // No stablization in the geometric jacobian
130  const auto grad_trial_ust = gradTrialUnstabilized(beta);
131  J += _stress[_qp].doubleContraction(grad_test) * grad_trial_ust.trace() -
132  _stress[_qp].doubleContraction(grad_test * grad_trial_ust);
133  }
134 
135  return J;
136 }
virtual RankTwoTensor gradTrial(unsigned int component) override
virtual RankTwoTensor gradTest(unsigned int component) override
const bool _large_kinematics
If true use large deformation kinematics.
const MaterialProperty< RankTwoTensor > & _stress
The Cauchy stress.
const MaterialProperty< RankFourTensor > & _material_jacobian
virtual RankTwoTensor gradTrialUnstabilized(unsigned int component)
The unstabilized trial function gradient.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string alpha
Definition: NS.h:134
const MaterialProperty< RankTwoTensor > & _f_inv
The inverse increment deformation gradient.

◆ computeQpJacobianOutOfPlaneStrain()

template<class G >
virtual Real UpdatedLagrangianStressDivergenceBase< G >::computeQpJacobianOutOfPlaneStrain ( )
inlineoverrideprotectedvirtual

Implements LagrangianStressDivergenceBase.

Definition at line 51 of file UpdatedLagrangianStressDivergence.h.

51 { return 0; }

◆ computeQpJacobianTemperature()

template<class G >
Real UpdatedLagrangianStressDivergenceBase< G >::computeQpJacobianTemperature ( unsigned int  cvar)
overrideprotectedvirtual

Implements LagrangianStressDivergenceBase.

Definition at line 140 of file UpdatedLagrangianStressDivergence.C.

141 {
142  // Multiple eigenstrains may depend on the same coupled var
143  RankTwoTensor total_deigen;
144  for (const auto deigen_darg : _deigenstrain_dargs[cvar])
145  total_deigen += (*deigen_darg)[_qp];
146 
148  RankFourTensor Csym = 0.5 * (C + C.transposeMajor().transposeIj().transposeMajor());
149 
150  return -(Csym * total_deigen).doubleContraction(gradTest(_alpha)) * _temperature->phi()[_j][_qp];
151 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
virtual RankTwoTensor gradTest(unsigned int component) override
const FieldVariablePhiValue & phi() const override
const MaterialProperty< RankFourTensor > & _material_jacobian
const MooseVariable * _temperature
Temperature, if provided. This is used only to get the trial functions.
std::vector< std::vector< const MaterialProperty< RankTwoTensor > * > > _deigenstrain_dargs
Eigenstrain derivatives wrt generate coupleds.
static const std::string C
Definition: NS.h:168

◆ computeQpOffDiagJacobian()

Real LagrangianStressDivergenceBase::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtualinherited

Reimplemented in TotalLagrangianWeakPlaneStress.

Definition at line 121 of file LagrangianStressDivergenceBase.C.

122 {
123  // Bail if jvar not coupled
124  if (getJvarMap()[jvar] < 0)
125  return 0.0;
126 
127  // Off diagonal terms for other displacements
128  for (auto beta : make_range(_ndisp))
129  if (jvar == _disp_nums[beta])
131 
132  // Off diagonal temperature term due to eigenstrain
133  if (_temperature && jvar == _temperature->number())
135 
136  // Off diagonal term due to weak plane stress
139 
140  return 0;
141 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
virtual Real computeQpJacobianDisplacement(unsigned int alpha, unsigned int beta)=0
unsigned int number() const
virtual Real computeQpJacobianOutOfPlaneStrain()=0
const MooseVariable * _out_of_plane_strain
Out-of-plane strain, if provided.
std::vector< unsigned int > _disp_nums
The displacement numbers.
const unsigned int _ndisp
Total number of displacements/size of residual vector.
const MooseVariable * _temperature
Temperature, if provided. This is used only to get the trial functions.
IntRange< T > make_range(T beg, T end)
virtual Real computeQpJacobianTemperature(unsigned int cvar)=0

◆ computeQpResidual()

template<class G >
Real UpdatedLagrangianStressDivergenceBase< G >::computeQpResidual ( )
overrideprotectedvirtual

Definition at line 105 of file UpdatedLagrangianStressDivergence.C.

106 {
108 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
virtual RankTwoTensor gradTest(unsigned int component) override
Real doubleContraction(const RankTwoTensorTempl< Real > &a) const
const MaterialProperty< RankTwoTensor > & _stress
The Cauchy stress.

◆ gradTest()

template<class G >
RankTwoTensor UpdatedLagrangianStressDivergenceBase< G >::gradTest ( unsigned int  component)
overrideprotectedvirtual

Implements LagrangianStressDivergenceBase.

Definition at line 47 of file UpdatedLagrangianStressDivergence.C.

48 {
49  // F-bar doesn't modify the test function
50  return G::gradOp(component, _grad_test[_i][_qp], _test[_i][_qp], _q_point[_qp]);
51 }
static const std::string component
Definition: NS.h:153

◆ gradTrial()

template<class G >
RankTwoTensor UpdatedLagrangianStressDivergenceBase< G >::gradTrial ( unsigned int  component)
overrideprotectedvirtual

Implements LagrangianStressDivergenceBase.

Definition at line 55 of file UpdatedLagrangianStressDivergence.C.

56 {
58 }
static const std::string component
Definition: NS.h:153
const bool _stabilize_strain
If true calculate the deformation gradient derivatives for F_bar.
virtual RankTwoTensor gradTrialStabilized(unsigned int component)
The stabilized trial function gradient.
virtual RankTwoTensor gradTrialUnstabilized(unsigned int component)
The unstabilized trial function gradient.

◆ gradTrialStabilized()

template<class G >
RankTwoTensor UpdatedLagrangianStressDivergenceBase< G >::gradTrialStabilized ( unsigned int  component)
privatevirtual

The stabilized trial function gradient.

Definition at line 70 of file UpdatedLagrangianStressDivergence.C.

71 {
72  // The base unstabilized trial function gradient
73  const auto Gb = G::gradOp(component, _grad_phi[_j][_qp], _phi[_j][_qp], _q_point[_qp]);
74  // The average trial function gradient
75  const auto Ga = _avg_grad_trial[component][_j];
76  // For updated Lagrangian, the modification is always linear
77  return Gb + (Ga.trace() - Gb.trace()) / 3.0 * RankTwoTensor::Identity();
78 }
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
static const std::string component
Definition: NS.h:153
static RankTwoTensorTempl Identity()

◆ gradTrialUnstabilized()

template<class G >
RankTwoTensor UpdatedLagrangianStressDivergenceBase< G >::gradTrialUnstabilized ( unsigned int  component)
privatevirtual

The unstabilized trial function gradient.

Definition at line 62 of file UpdatedLagrangianStressDivergence.C.

63 {
64  // Without F-bar stabilization, simply return the gradient of the trial functions
65  return G::gradOp(component, _grad_phi[_j][_qp], _phi[_j][_qp], _q_point[_qp]);
66 }
static const std::string component
Definition: NS.h:153

◆ initialSetup() [1/2]

template<class G >
virtual void UpdatedLagrangianStressDivergenceBase< G >::initialSetup ( )
overridevirtual

◆ initialSetup() [2/2]

Definition at line 89 of file UpdatedLagrangianStressDivergence.h.

90 {
91  if (getBlockCoordSystem() != Moose::COORD_XYZ)
92  mooseError("This kernel should only act in Cartesian coordinates.");
93 }
void mooseError(Args &&... args)

◆ precalculateJacobian()

void LagrangianStressDivergenceBase::precalculateJacobian ( )
overrideprotectedvirtualinherited

Reimplemented in TotalLagrangianWeakPlaneStress.

Definition at line 85 of file LagrangianStressDivergenceBase.C.

86 {
87  // Skip if we are not doing stabilization
88  if (!_stabilize_strain)
89  return;
90 
91  // We need the gradients of shape functions in the reference frame
92  _fe_problem.prepareShapes(_var.number(), _tid);
93  _avg_grad_trial[_alpha].resize(_phi.size());
95 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
virtual void precalculateJacobianDisplacement(unsigned int component)=0
Prepare the average shape function gradients for stabilization.
const bool _stabilize_strain
If true calculate the deformation gradient derivatives for F_bar.

◆ precalculateJacobianDisplacement()

template<class G >
void UpdatedLagrangianStressDivergenceBase< G >::precalculateJacobianDisplacement ( unsigned int  component)
overrideprotectedvirtual

Prepare the average shape function gradients for stabilization.

Implements LagrangianStressDivergenceBase.

Definition at line 82 of file UpdatedLagrangianStressDivergence.C.

83 {
84  // For updated Lagrangian, the averaging is taken on the reference frame. If large kinematics is
85  // used, the averaged gradients should be pushed forward to the current frame.
86  for (auto j : make_range(_phi.size()))
87  {
89  [this, component, j](unsigned int qp)
90  {
91  return G::gradOp(
93  },
97  // Push forward to the current frame.
98  // The average deformation gradient is the same at all qps.
99  _avg_grad_trial[component][j] *= _F_avg[0].inverse();
100  }
101 }
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
static const std::string component
Definition: NS.h:153
const MaterialProperty< RankTwoTensor > & _F_avg
The element-average deformation gradient.
const bool _large_kinematics
If true use large deformation kinematics.
IntRange< T > make_range(T beg, T end)
auto elementAverage(const Functor &f, const MooseArray< Real > &JxW, const MooseArray< Real > &coord)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")

◆ precalculateOffDiagJacobian()

void LagrangianStressDivergenceBase::precalculateOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtualinherited

Definition at line 98 of file LagrangianStressDivergenceBase.C.

99 {
100  // Skip if we are not doing stabilization
101  if (!_stabilize_strain)
102  return;
103 
104  for (auto beta : make_range(_ndisp))
105  if (jvar == _disp_nums[beta])
106  {
107  // We need the gradients of shape functions in the reference frame
108  _fe_problem.prepareShapes(jvar, _tid);
109  _avg_grad_trial[beta].resize(_phi.size());
111  }
112 }
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
virtual void precalculateJacobianDisplacement(unsigned int component)=0
Prepare the average shape function gradients for stabilization.
const bool _stabilize_strain
If true calculate the deformation gradient derivatives for F_bar.
std::vector< unsigned int > _disp_nums
The displacement numbers.
const unsigned int _ndisp
Total number of displacements/size of residual vector.
IntRange< T > make_range(T beg, T end)

◆ validParams() [1/2]

template<class G >
static InputParameters UpdatedLagrangianStressDivergenceBase< G >::validParams ( )
static

◆ validParams() [2/2]

Definition at line 79 of file UpdatedLagrangianStressDivergence.h.

80 {
82  params.addClassDescription(
83  "Enforce equilibrium with an updated Lagrangian formulation in Cartesian coordinates.");
84  return params;
85 }
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _alpha

const unsigned int LagrangianStressDivergenceBase::_alpha
protectedinherited

◆ _assembly_undisplaced

template<class G >
Assembly& UpdatedLagrangianStressDivergenceBase< G >::_assembly_undisplaced
protected

Definition at line 62 of file UpdatedLagrangianStressDivergence.h.

◆ _avg_grad_trial

std::vector<std::vector<RankTwoTensor> > LagrangianStressDivergenceBase::_avg_grad_trial
protectedinherited

◆ _base_name

const std::string LagrangianStressDivergenceBase::_base_name
protectedinherited

Prepend to the material properties.

Definition at line 69 of file LagrangianStressDivergenceBase.h.

◆ _coord_undisplaced

template<class G >
const MooseArray<Real>& UpdatedLagrangianStressDivergenceBase< G >::_coord_undisplaced
protected

Definition at line 65 of file UpdatedLagrangianStressDivergence.h.

◆ _deigenstrain_dargs

std::vector<std::vector<const MaterialProperty<RankTwoTensor> *> > LagrangianStressDivergenceBase::_deigenstrain_dargs
protectedinherited

Eigenstrain derivatives wrt generate coupleds.

Definition at line 107 of file LagrangianStressDivergenceBase.h.

Referenced by LagrangianStressDivergenceBase::LagrangianStressDivergenceBase().

◆ _disp_nums

std::vector<unsigned int> LagrangianStressDivergenceBase::_disp_nums
protectedinherited

◆ _F

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_F
protectedinherited

The actual (stabilized) deformation gradient.

Definition at line 98 of file LagrangianStressDivergenceBase.h.

◆ _F_avg

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_F_avg
protectedinherited

The element-average deformation gradient.

Definition at line 89 of file LagrangianStressDivergenceBase.h.

◆ _f_inv

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_f_inv
protectedinherited

The inverse increment deformation gradient.

Definition at line 92 of file LagrangianStressDivergenceBase.h.

◆ _F_inv

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_F_inv
protectedinherited

The inverse deformation gradient.

Definition at line 95 of file LagrangianStressDivergenceBase.h.

◆ _F_ust

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_F_ust
protectedinherited

The unmodified deformation gradient.

Definition at line 86 of file LagrangianStressDivergenceBase.h.

◆ _grad_phi_undisplaced

template<class G >
const VariablePhiGradient& UpdatedLagrangianStressDivergenceBase< G >::_grad_phi_undisplaced
protected

Definition at line 63 of file UpdatedLagrangianStressDivergence.h.

◆ _JxW_undisplaced

template<class G >
const MooseArray<Real>& UpdatedLagrangianStressDivergenceBase< G >::_JxW_undisplaced
protected

Definition at line 64 of file UpdatedLagrangianStressDivergence.h.

◆ _large_kinematics

const bool LagrangianStressDivergenceBase::_large_kinematics
protectedinherited

◆ _material_jacobian

template<class G >
const MaterialProperty<RankFourTensor>& UpdatedLagrangianStressDivergenceBase< G >::_material_jacobian
protected

Definition at line 58 of file UpdatedLagrangianStressDivergence.h.

◆ _ndisp

const unsigned int LagrangianStressDivergenceBase::_ndisp
protectedinherited

◆ _out_of_plane_strain

const MooseVariable* LagrangianStressDivergenceBase::_out_of_plane_strain
protectedinherited

◆ _q_point_undisplaced

template<class G >
const MooseArray<Point>& UpdatedLagrangianStressDivergenceBase< G >::_q_point_undisplaced
protected

Definition at line 66 of file UpdatedLagrangianStressDivergence.h.

◆ _stabilize_strain

const bool LagrangianStressDivergenceBase::_stabilize_strain
protectedinherited

If true calculate the deformation gradient derivatives for F_bar.

Definition at line 66 of file LagrangianStressDivergenceBase.h.

Referenced by LagrangianStressDivergenceBase::precalculateJacobian(), and LagrangianStressDivergenceBase::precalculateOffDiagJacobian().

◆ _stress

template<class G >
const MaterialProperty<RankTwoTensor>& UpdatedLagrangianStressDivergenceBase< G >::_stress
protected

The Cauchy stress.

Definition at line 54 of file UpdatedLagrangianStressDivergence.h.

◆ _temperature

const MooseVariable* LagrangianStressDivergenceBase::_temperature
protectedinherited

Temperature, if provided. This is used only to get the trial functions.

Definition at line 101 of file LagrangianStressDivergenceBase.h.

Referenced by LagrangianStressDivergenceBase::computeQpOffDiagJacobian().


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