www.mooseframework.org
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
StressDivergenceRZTensors Class Reference

StressDivergenceRZTensors is a modification of StressDivergenceTensors to accommodate the Axisymmetric material models that use cylindrical coordinates. More...

#include <StressDivergenceRZTensors.h>

Inheritance diagram for StressDivergenceRZTensors:
[legend]

Public Member Functions

 StressDivergenceRZTensors (const InputParameters &parameters)
 
virtual void computeJacobian () override
 
virtual void computeOffDiagJacobian (unsigned int jvar) override
 

Protected Member Functions

virtual void initialSetup () override
 
virtual Real computeQpResidual () override
 
virtual Real computeQpJacobian () override
 
virtual Real computeQpOffDiagJacobian (unsigned int jvar) override
 
virtual void computeAverageGradientTest () override
 
virtual void computeAverageGradientPhi () override
 
Real calculateJacobian (unsigned int ivar, unsigned int jvar)
 
virtual void computeResidual () override
 
virtual void computeFiniteDeformJacobian ()
 

Protected Attributes

std::string _base_name
 
bool _use_finite_deform_jacobian
 
const MaterialProperty< RankTwoTensor > & _stress
 
const MaterialProperty< RankFourTensor > & _Jacobian_mult
 
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
 
const MaterialProperty< RankTwoTensor > * _deformation_gradient
 
const MaterialProperty< RankTwoTensor > * _deformation_gradient_old
 
const MaterialProperty< RankTwoTensor > * _rotation_increment
 
const unsigned int _component
 
unsigned int _ndisp
 Coupled displacement variables. More...
 
std::vector< unsigned int > _disp_var
 
const bool _temp_coupled
 
const unsigned int _temp_var
 
const MaterialProperty< RankTwoTensor > *const _deigenstrain_dT
 d(strain)/d(temperature), if computed by ComputeThermalExpansionEigenstrain More...
 
std::vector< std::vector< Real > > _avg_grad_test
 Gradient of test function averaged over the element. Used in volumetric locking correction calculation. More...
 
std::vector< std::vector< Real > > _avg_grad_phi
 Gradient of phi function averaged over the element. Used in volumetric locking correction calculation. More...
 
bool _volumetric_locking_correction
 Flag for volumetric locking correction. More...
 
Assembly & _assembly_undisplaced
 undisplaced problem More...
 
MooseVariable & _var_undisplaced
 Reference to this Kernel's undisplaced MooseVariable object. More...
 
const VariablePhiGradient & _grad_phi_undisplaced
 Shape and test functions on the undisplaced mesh. More...
 
const VariableTestGradient & _grad_test_undisplaced
 

Detailed Description

StressDivergenceRZTensors is a modification of StressDivergenceTensors to accommodate the Axisymmetric material models that use cylindrical coordinates.

This kernel is for symmetrical loading only. The key modifications are a result of the circumferential stress' dependence on displacement in the axial direction. Reference: Cook et.al. Concepts and Applications of Finite Element Analysis, 4th Ed. 2002. p 510. Within this kernel, '_disp_x' refers to displacement in the radial direction, u_r, and '_disp_y' refers to displacement in the axial direction, u_z. The COORD_TYPE in the Problem block must be set to RZ.

Definition at line 29 of file StressDivergenceRZTensors.h.

Constructor & Destructor Documentation

StressDivergenceRZTensors::StressDivergenceRZTensors ( const InputParameters &  parameters)

Definition at line 29 of file StressDivergenceRZTensors.C.

30  : StressDivergenceTensors(parameters)
31 {
32 }
StressDivergenceTensors(const InputParameters &parameters)

Member Function Documentation

Real StressDivergenceRZTensors::calculateJacobian ( unsigned int  ivar,
unsigned int  jvar 
)
protected

Definition at line 115 of file StressDivergenceRZTensors.C.

Referenced by computeQpJacobian(), and computeQpOffDiagJacobian().

116 {
117  // B^T_i * C * B_j
118  RealGradient test, test_z, phi, phi_z;
119  Real first_term = 0.0;
120  if (ivar == 0) // Case grad_test for x, requires contributions from stress_xx, stress_xy, and
121  // stress_zz
122  {
123  test(0) = _grad_test[_i][_qp](0);
124  test(1) = _grad_test[_i][_qp](1);
125  test_z(2) = _test[_i][_qp] / _q_point[_qp](0);
126  }
127  else // Case grad_test for y
128  {
129  test(0) = _grad_test[_i][_qp](0);
130  test(1) = _grad_test[_i][_qp](1);
131  }
132 
133  if (jvar == 0)
134  {
135  phi(0) = _grad_phi[_j][_qp](0);
136  phi(1) = _grad_phi[_j][_qp](1);
137  phi_z(2) = _phi[_j][_qp] / _q_point[_qp](0);
138  }
139  else
140  {
141  phi(0) = _grad_phi[_j][_qp](0);
142  phi(1) = _grad_phi[_j][_qp](1);
143  }
144 
145  if (ivar == 0 &&
146  jvar == 0) // Case when both phi and test are functions of x and z; requires four terms
147  {
148  const Real first_sum = ElasticityTensorTools::elasticJacobian(
149  _Jacobian_mult[_qp], ivar, jvar, test, phi); // test_x and phi_x
150  const Real second_sum = ElasticityTensorTools::elasticJacobian(
151  _Jacobian_mult[_qp], 2, 2, test_z, phi_z); // test_z and phi_z
152  const Real mixed_sum1 = ElasticityTensorTools::elasticJacobian(
153  _Jacobian_mult[_qp], ivar, 2, test, phi_z); // test_x and phi_z
154  const Real mixed_sum2 = ElasticityTensorTools::elasticJacobian(
155  _Jacobian_mult[_qp], 2, jvar, test_z, phi); // test_z and phi_x
156 
157  first_term = first_sum + second_sum + mixed_sum1 + mixed_sum2;
158  }
159  else if (ivar == 0 && jvar == 1)
160  {
161  const Real first_sum = ElasticityTensorTools::elasticJacobian(
162  _Jacobian_mult[_qp], ivar, jvar, test, phi); // test_x and phi_y
163  const Real mixed_sum2 = ElasticityTensorTools::elasticJacobian(
164  _Jacobian_mult[_qp], 2, jvar, test_z, phi); // test_z and phi_y
165 
166  first_term = first_sum + mixed_sum2;
167  }
168  else if (ivar == 1 && jvar == 0)
169  {
170  const Real second_sum = ElasticityTensorTools::elasticJacobian(
171  _Jacobian_mult[_qp], ivar, jvar, test, phi); // test_y and phi_x
172  const Real mixed_sum1 = ElasticityTensorTools::elasticJacobian(
173  _Jacobian_mult[_qp], ivar, 2, test, phi_z); // test_y and phi_z
174 
175  first_term = second_sum + mixed_sum1;
176  }
177  else if (ivar == 1 && jvar == 1)
179  _Jacobian_mult[_qp], ivar, jvar, test, phi); // test_y and phi_y
180  else
181  mooseError("Invalid component in Jacobian Calculation");
182 
183  Real val = 0.0;
184  // volumetric locking correction
185  // K = Bbar^T_i * C * Bbar^T_j where Bbar = B + Bvol
186  // K = B^T_i * C * B_j + Bvol^T_i * C * Bvol_j + B^T_i * C * Bvol_j + Bvol^T_i * C * B_j
188  {
189  RealGradient new_test(2, 0.0);
190  RealGradient new_phi(2, 0.0);
191 
192  new_test(0) = _grad_test[_i][_qp](0) + _test[_i][_qp] / _q_point[_qp](0);
193  new_test(1) = _grad_test[_i][_qp](1);
194  new_phi(0) = _grad_phi[_j][_qp](0) + _phi[_j][_qp] / _q_point[_qp](0);
195  new_phi(1) = _grad_phi[_j][_qp](1);
196 
197  // Bvol^T_i * C * Bvol_j
198  val += _Jacobian_mult[_qp].sum3x3() * (_avg_grad_test[_i][ivar] - new_test(ivar)) *
199  (_avg_grad_phi[_j][jvar] - new_phi(jvar)) / 3.0;
200 
201  // B^T_i * C * Bvol_j
202  RealGradient sum_3x1 = _Jacobian_mult[_qp].sum3x1();
203  if (ivar == 0 && jvar == 0)
204  val += (sum_3x1(0) * test(0) + sum_3x1(2) * test_z(2)) * (_avg_grad_phi[_j][0] - new_phi(0));
205  else if (ivar == 0 && jvar == 1)
206  val += (sum_3x1(0) * test(0) + sum_3x1(2) * test_z(2)) * (_avg_grad_phi[_j][1] - new_phi(1));
207  else if (ivar == 1 && jvar == 0)
208  val += sum_3x1(1) * test(1) * (_avg_grad_phi[_j][0] - new_phi(0));
209  else
210  val += sum_3x1(1) * test(1) * (_avg_grad_phi[_j][1] - new_phi(1));
211 
212  // Bvol^T_i * C * B_j
213  // val = trace (C * B_j) *(avg_grad_test[_i][ivar] - new_test(ivar))
214  if (jvar == 0)
215  for (unsigned int i = 0; i < 3; ++i)
216  val +=
217  (_Jacobian_mult[_qp](i, i, 0, 0) * phi(0) + _Jacobian_mult[_qp](i, i, 0, 1) * phi(1) +
218  _Jacobian_mult[_qp](i, i, 2, 2) * phi_z(2)) *
219  (_avg_grad_test[_i][ivar] - new_test(ivar));
220  else if (jvar == 1)
221  for (unsigned int i = 0; i < 3; ++i)
222  val +=
223  (_Jacobian_mult[_qp](i, i, 0, 1) * phi(0) + _Jacobian_mult[_qp](i, i, 1, 1) * phi(1)) *
224  (_avg_grad_test[_i][ivar] - new_test(ivar));
225  }
226 
227  return val / 3.0 + first_term;
228 }
std::vector< std::vector< Real > > _avg_grad_phi
Gradient of phi function averaged over the element. Used in volumetric locking correction calculation...
bool _volumetric_locking_correction
Flag for volumetric locking correction.
std::vector< std::vector< Real > > _avg_grad_test
Gradient of test function averaged over the element. Used in volumetric locking correction calculatio...
Real elasticJacobian(const RankFourTensor &r4t, unsigned int i, unsigned int k, const RealGradient &grad_test, const RealGradient &grad_phi)
This is used for the standard kernel stress_ij*d(test)/dx_j, when varied wrt u_k Jacobian entry: d(st...
const MaterialProperty< RankFourTensor > & _Jacobian_mult
void StressDivergenceRZTensors::computeAverageGradientPhi ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 253 of file StressDivergenceRZTensors.C.

254 {
255  _avg_grad_phi.resize(_phi.size());
256  for (_i = 0; _i < _phi.size(); ++_i)
257  {
258  _avg_grad_phi[_i].resize(2);
259  for (unsigned int component = 0; component < 2; ++component)
260  {
261  _avg_grad_phi[_i][component] = 0.0;
262  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
263  {
264  if (component == 0)
265  _avg_grad_phi[_i][component] +=
266  (_grad_phi[_i][_qp](component) + _phi[_i][_qp] / _q_point[_qp](0)) * _JxW[_qp] *
267  _coord[_qp];
268  else
269  _avg_grad_phi[_i][component] += _grad_phi[_i][_qp](component) * _JxW[_qp] * _coord[_qp];
270  }
271  _avg_grad_phi[_i][component] /= _current_elem_volume;
272  }
273  }
274 }
std::vector< std::vector< Real > > _avg_grad_phi
Gradient of phi function averaged over the element. Used in volumetric locking correction calculation...
Real component(const SymmTensor &symm_tensor, unsigned int index)
void StressDivergenceRZTensors::computeAverageGradientTest ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 231 of file StressDivergenceRZTensors.C.

232 {
233  // calculate volume averaged value of shape function derivative
234  _avg_grad_test.resize(_test.size());
235  for (_i = 0; _i < _test.size(); ++_i)
236  {
237  _avg_grad_test[_i].resize(2);
238  _avg_grad_test[_i][_component] = 0.0;
239  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
240  {
241  if (_component == 0)
242  _avg_grad_test[_i][_component] +=
243  (_grad_test[_i][_qp](_component) + _test[_i][_qp] / _q_point[_qp](0)) * _JxW[_qp] *
244  _coord[_qp];
245  else
246  _avg_grad_test[_i][_component] += _grad_test[_i][_qp](_component) * _JxW[_qp] * _coord[_qp];
247  }
248  _avg_grad_test[_i][_component] /= _current_elem_volume;
249  }
250 }
std::vector< std::vector< Real > > _avg_grad_test
Gradient of test function averaged over the element. Used in volumetric locking correction calculatio...
void StressDivergenceTensors::computeFiniteDeformJacobian ( )
protectedvirtualinherited

Definition at line 301 of file StressDivergenceTensors.C.

Referenced by StressDivergenceTensors::computeJacobian(), and StressDivergenceTensors::computeOffDiagJacobian().

302 {
303  const RankTwoTensor I(RankTwoTensor::initIdentity);
304  const RankFourTensor II_ijkl = I.mixedProductIkJl(I);
305 
306  // Bring back to unrotated config
307  const RankTwoTensor unrotated_stress =
308  (*_rotation_increment)[_qp].transpose() * _stress[_qp] * (*_rotation_increment)[_qp];
309 
310  // Incremental deformation gradient Fhat
311  const RankTwoTensor Fhat =
312  (*_deformation_gradient)[_qp] * (*_deformation_gradient_old)[_qp].inverse();
313  const RankTwoTensor Fhatinv = Fhat.inverse();
314 
315  const RankTwoTensor rot_times_stress = (*_rotation_increment)[_qp] * unrotated_stress;
316  const RankFourTensor dstress_drot =
317  I.mixedProductIkJl(rot_times_stress) + I.mixedProductJkIl(rot_times_stress);
318  const RankFourTensor rot_rank_four =
319  (*_rotation_increment)[_qp].mixedProductIkJl((*_rotation_increment)[_qp]);
320  const RankFourTensor drot_dUhatinv = Fhat.mixedProductIkJl(I);
321 
322  const RankTwoTensor A = I - Fhatinv;
323 
324  // Ctilde = Chat^-1 - I
325  const RankTwoTensor Ctilde = A * A.transpose() - A - A.transpose();
326  const RankFourTensor dCtilde_dFhatinv =
327  -I.mixedProductIkJl(A) - I.mixedProductJkIl(A) + II_ijkl + I.mixedProductJkIl(I);
328 
329  // Second order approximation of Uhat - consistent with strain increment definition
330  // const RankTwoTensor Uhat = I - 0.5 * Ctilde - 3.0/8.0 * Ctilde * Ctilde;
331 
332  RankFourTensor dUhatinv_dCtilde =
333  0.5 * II_ijkl - 1.0 / 8.0 * (I.mixedProductIkJl(Ctilde) + Ctilde.mixedProductIkJl(I));
334  RankFourTensor drot_dFhatinv = drot_dUhatinv * dUhatinv_dCtilde * dCtilde_dFhatinv;
335 
336  drot_dFhatinv -= Fhat.mixedProductIkJl((*_rotation_increment)[_qp].transpose());
337  _finite_deform_Jacobian_mult[_qp] = dstress_drot * drot_dFhatinv;
338 
339  const RankFourTensor dstrain_increment_dCtilde =
340  -0.5 * II_ijkl + 0.25 * (I.mixedProductIkJl(Ctilde) + Ctilde.mixedProductIkJl(I));
342  rot_rank_four * _Jacobian_mult[_qp] * dstrain_increment_dCtilde * dCtilde_dFhatinv;
343  _finite_deform_Jacobian_mult[_qp] += Fhat.mixedProductJkIl(_stress[_qp]);
344 
345  const RankFourTensor dFhat_dFhatinv = -Fhat.mixedProductIkJl(Fhat.transpose());
346  const RankTwoTensor dJ_dFhatinv = dFhat_dFhatinv.innerProductTranspose(Fhat.ddet());
347 
348  // Component from Jacobian derivative
349  _finite_deform_Jacobian_mult[_qp] += _stress[_qp].outerProduct(dJ_dFhatinv);
350 
351  // Derivative of Fhatinv w.r.t. undisplaced coordinates
352  const RankTwoTensor Finv = (*_deformation_gradient)[_qp].inverse();
353  const RankFourTensor dFhatinv_dGradu = -Fhatinv.mixedProductIkJl(Finv.transpose());
354  _finite_deform_Jacobian_mult[_qp] = _finite_deform_Jacobian_mult[_qp] * dFhatinv_dGradu;
355 }
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
const MaterialProperty< RankFourTensor > & _Jacobian_mult
const MaterialProperty< RankTwoTensor > & _stress
void StressDivergenceTensors::computeJacobian ( )
overridevirtualinherited

Reimplemented from ALEKernel.

Definition at line 136 of file StressDivergenceTensors.C.

137 {
139  {
140  _finite_deform_Jacobian_mult.resize(_qrule->n_points());
141 
142  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
144 
146  }
147  else
148  {
150  {
153  }
154  Kernel::computeJacobian();
155  }
156 }
bool _volumetric_locking_correction
Flag for volumetric locking correction.
virtual void computeJacobian()
Definition: ALEKernel.C:33
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
virtual void computeFiniteDeformJacobian()
void StressDivergenceTensors::computeOffDiagJacobian ( unsigned int  jvar)
overridevirtualinherited

Reimplemented from ALEKernel.

Definition at line 159 of file StressDivergenceTensors.C.

160 {
162  {
163  _finite_deform_Jacobian_mult.resize(_qrule->n_points());
164 
165  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
167 
169  }
170  else
171  {
173  {
176  }
177  Kernel::computeOffDiagJacobian(jvar);
178  }
179 }
virtual void computeOffDiagJacobian(unsigned int jvar)
Definition: ALEKernel.C:40
bool _volumetric_locking_correction
Flag for volumetric locking correction.
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
virtual void computeFiniteDeformJacobian()
Real StressDivergenceRZTensors::computeQpJacobian ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 73 of file StressDivergenceRZTensors.C.

74 {
76 }
Real calculateJacobian(unsigned int ivar, unsigned int jvar)
Real StressDivergenceRZTensors::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 79 of file StressDivergenceRZTensors.C.

80 {
81  for (unsigned int i = 0; i < _ndisp; ++i)
82  {
83  if (jvar == _disp_var[i])
84  return calculateJacobian(_component, i);
85  }
86 
87  if (_temp_coupled && jvar == _temp_var)
88  {
89  Real jac = 0.0;
90  if (_component == 0)
91  {
92  for (unsigned k = 0; k < LIBMESH_DIM; ++k)
93  for (unsigned l = 0; l < LIBMESH_DIM; ++l)
94  jac -= (_grad_test[_i][_qp](0) * _Jacobian_mult[_qp](0, 0, k, l) +
95  _test[_i][_qp] / _q_point[_qp](0) * _Jacobian_mult[_qp](2, 2, k, l) +
96  _grad_test[_i][_qp](1) * _Jacobian_mult[_qp](0, 1, k, l)) *
97  (*_deigenstrain_dT)[_qp](k, l);
98  return jac * _phi[_j][_qp];
99  }
100  else if (_component == 1)
101  {
102  for (unsigned k = 0; k < LIBMESH_DIM; ++k)
103  for (unsigned l = 0; l < LIBMESH_DIM; ++l)
104  jac -= (_grad_test[_i][_qp](1) * _Jacobian_mult[_qp](1, 1, k, l) +
105  _grad_test[_i][_qp](0) * _Jacobian_mult[_qp](1, 0, k, l)) *
106  (*_deigenstrain_dT)[_qp](k, l);
107  return jac * _phi[_j][_qp];
108  }
109  }
110 
111  return 0.0;
112 }
unsigned int _ndisp
Coupled displacement variables.
std::vector< unsigned int > _disp_var
const MaterialProperty< RankFourTensor > & _Jacobian_mult
Real calculateJacobian(unsigned int ivar, unsigned int jvar)
Real StressDivergenceRZTensors::computeQpResidual ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 43 of file StressDivergenceRZTensors.C.

44 {
45  Real div = 0.0;
46  if (_component == 0)
47  {
48  div = _grad_test[_i][_qp](0) * _stress[_qp](0, 0) +
49  +(_test[_i][_qp] / _q_point[_qp](0)) * _stress[_qp](2, 2) +
50  +_grad_test[_i][_qp](1) * _stress[_qp](0, 1); // stress_{rz}
51 
52  // volumetric locking correction
54  div += (_avg_grad_test[_i][0] - _grad_test[_i][_qp](0) - _test[_i][_qp] / _q_point[_qp](0)) *
55  (_stress[_qp].trace()) / 3.0;
56  }
57  else if (_component == 1)
58  {
59  div = _grad_test[_i][_qp](1) * _stress[_qp](1, 1) +
60  +_grad_test[_i][_qp](0) * _stress[_qp](1, 0); // stress_{zr}
61 
62  // volumetric locking correction
64  div += (_avg_grad_test[_i][1] - _grad_test[_i][_qp](1)) * (_stress[_qp].trace()) / 3.0;
65  }
66  else
67  mooseError("Invalid component for this AxisymmetricRZ problem.");
68 
69  return div;
70 }
bool _volumetric_locking_correction
Flag for volumetric locking correction.
std::vector< std::vector< Real > > _avg_grad_test
Gradient of test function averaged over the element. Used in volumetric locking correction calculatio...
const MaterialProperty< RankTwoTensor > & _stress
void StressDivergenceTensors::computeResidual ( )
overrideprotectedvirtualinherited

Definition at line 99 of file StressDivergenceTensors.C.

100 {
101  DenseVector<Number> & re = _assembly.residualBlock(_var.number());
102  _local_re.resize(re.size());
103  _local_re.zero();
104 
107 
108  precalculateResidual();
109  for (_i = 0; _i < _test.size(); ++_i)
110  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
111  _local_re(_i) += _JxW[_qp] * _coord[_qp] * computeQpResidual();
112 
113  re += _local_re;
114 
115  if (_has_save_in)
116  {
117  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
118  for (const auto & var : _save_in)
119  var->sys().solution().add_vector(_local_re, var->dofIndices());
120  }
121 }
bool _volumetric_locking_correction
Flag for volumetric locking correction.
virtual Real computeQpResidual() override
void StressDivergenceRZTensors::initialSetup ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 35 of file StressDivergenceRZTensors.C.

36 {
37  if (getBlockCoordSystem() != Moose::COORD_RZ)
38  mooseError("The coordinate system in the Problem block must be set to RZ for axisymmetric "
39  "geometries.");
40 }

Member Data Documentation

Assembly& ALEKernel::_assembly_undisplaced
protectedinherited

undisplaced problem

Definition at line 29 of file ALEKernel.h.

std::vector<std::vector<Real> > StressDivergenceTensors::_avg_grad_phi
protectedinherited

Gradient of phi function averaged over the element. Used in volumetric locking correction calculation.

Definition at line 76 of file StressDivergenceTensors.h.

Referenced by calculateJacobian(), computeAverageGradientPhi(), StressDivergenceTensors::computeAverageGradientPhi(), StressDivergenceTensors::computeQpJacobian(), and StressDivergenceTensors::computeQpOffDiagJacobian().

std::vector<std::vector<Real> > StressDivergenceTensors::_avg_grad_test
protectedinherited
std::string StressDivergenceTensors::_base_name
protectedinherited
const unsigned int StressDivergenceTensors::_component
protectedinherited
const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_deformation_gradient
protectedinherited
const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_deformation_gradient_old
protectedinherited
const MaterialProperty<RankTwoTensor>* const StressDivergenceTensors::_deigenstrain_dT
protectedinherited

d(strain)/d(temperature), if computed by ComputeThermalExpansionEigenstrain

Definition at line 70 of file StressDivergenceTensors.h.

std::vector<unsigned int> StressDivergenceTensors::_disp_var
protectedinherited
std::vector<RankFourTensor> StressDivergenceTensors::_finite_deform_Jacobian_mult
protectedinherited
const VariablePhiGradient& ALEKernel::_grad_phi_undisplaced
protectedinherited

Shape and test functions on the undisplaced mesh.

Definition at line 35 of file ALEKernel.h.

Referenced by StressDivergenceTensors::computeQpJacobian(), and StressDivergenceTensors::computeQpOffDiagJacobian().

const VariableTestGradient& ALEKernel::_grad_test_undisplaced
protectedinherited

Definition at line 36 of file ALEKernel.h.

const MaterialProperty<RankFourTensor>& StressDivergenceTensors::_Jacobian_mult
protectedinherited
unsigned int StressDivergenceTensors::_ndisp
protectedinherited
const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_rotation_increment
protectedinherited
const MaterialProperty<RankTwoTensor>& StressDivergenceTensors::_stress
protectedinherited
const bool StressDivergenceTensors::_temp_coupled
protectedinherited
const unsigned int StressDivergenceTensors::_temp_var
protectedinherited
bool StressDivergenceTensors::_use_finite_deform_jacobian
protectedinherited
MooseVariable& ALEKernel::_var_undisplaced
protectedinherited

Reference to this Kernel's undisplaced MooseVariable object.

Definition at line 32 of file ALEKernel.h.

bool StressDivergenceTensors::_volumetric_locking_correction
protectedinherited

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