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 (MooseVariableFEBase &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...
 
const bool _out_of_plane_strain_coupled
 
const unsigned int _out_of_plane_strain_var
 
const unsigned int _out_of_plane_direction
 
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 32 of file StressDivergenceRZTensors.h.

Constructor & Destructor Documentation

◆ StressDivergenceRZTensors()

StressDivergenceRZTensors::StressDivergenceRZTensors ( const InputParameters &  parameters)

Definition at line 33 of file StressDivergenceRZTensors.C.

34  : StressDivergenceTensors(parameters)
35 {
36 }
StressDivergenceTensors(const InputParameters &parameters)

Member Function Documentation

◆ calculateJacobian()

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

Definition at line 119 of file StressDivergenceRZTensors.C.

Referenced by computeQpJacobian(), and computeQpOffDiagJacobian().

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

◆ computeAverageGradientPhi()

void StressDivergenceRZTensors::computeAverageGradientPhi ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 257 of file StressDivergenceRZTensors.C.

258 {
259  _avg_grad_phi.resize(_phi.size());
260  for (_i = 0; _i < _phi.size(); ++_i)
261  {
262  _avg_grad_phi[_i].resize(2);
263  for (unsigned int component = 0; component < 2; ++component)
264  {
265  _avg_grad_phi[_i][component] = 0.0;
266  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
267  {
268  if (component == 0)
269  _avg_grad_phi[_i][component] +=
270  (_grad_phi[_i][_qp](component) + _phi[_i][_qp] / _q_point[_qp](0)) * _JxW[_qp] *
271  _coord[_qp];
272  else
273  _avg_grad_phi[_i][component] += _grad_phi[_i][_qp](component) * _JxW[_qp] * _coord[_qp];
274  }
275  _avg_grad_phi[_i][component] /= _current_elem_volume;
276  }
277  }
278 }
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)

◆ computeAverageGradientTest()

void StressDivergenceRZTensors::computeAverageGradientTest ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 235 of file StressDivergenceRZTensors.C.

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

◆ computeFiniteDeformJacobian()

void StressDivergenceTensors::computeFiniteDeformJacobian ( )
protectedvirtualinherited

Definition at line 329 of file StressDivergenceTensors.C.

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

330 {
331  const RankTwoTensor I(RankTwoTensor::initIdentity);
332  const RankFourTensor II_ijkl = I.mixedProductIkJl(I);
333 
334  // Bring back to unrotated config
335  const RankTwoTensor unrotated_stress =
336  (*_rotation_increment)[_qp].transpose() * _stress[_qp] * (*_rotation_increment)[_qp];
337 
338  // Incremental deformation gradient Fhat
339  const RankTwoTensor Fhat =
340  (*_deformation_gradient)[_qp] * (*_deformation_gradient_old)[_qp].inverse();
341  const RankTwoTensor Fhatinv = Fhat.inverse();
342 
343  const RankTwoTensor rot_times_stress = (*_rotation_increment)[_qp] * unrotated_stress;
344  const RankFourTensor dstress_drot =
345  I.mixedProductIkJl(rot_times_stress) + I.mixedProductJkIl(rot_times_stress);
346  const RankFourTensor rot_rank_four =
347  (*_rotation_increment)[_qp].mixedProductIkJl((*_rotation_increment)[_qp]);
348  const RankFourTensor drot_dUhatinv = Fhat.mixedProductIkJl(I);
349 
350  const RankTwoTensor A = I - Fhatinv;
351 
352  // Ctilde = Chat^-1 - I
353  const RankTwoTensor Ctilde = A * A.transpose() - A - A.transpose();
354  const RankFourTensor dCtilde_dFhatinv =
355  -I.mixedProductIkJl(A) - I.mixedProductJkIl(A) + II_ijkl + I.mixedProductJkIl(I);
356 
357  // Second order approximation of Uhat - consistent with strain increment definition
358  // const RankTwoTensor Uhat = I - 0.5 * Ctilde - 3.0/8.0 * Ctilde * Ctilde;
359 
360  RankFourTensor dUhatinv_dCtilde =
361  0.5 * II_ijkl - 1.0 / 8.0 * (I.mixedProductIkJl(Ctilde) + Ctilde.mixedProductIkJl(I));
362  RankFourTensor drot_dFhatinv = drot_dUhatinv * dUhatinv_dCtilde * dCtilde_dFhatinv;
363 
364  drot_dFhatinv -= Fhat.mixedProductIkJl((*_rotation_increment)[_qp].transpose());
365  _finite_deform_Jacobian_mult[_qp] = dstress_drot * drot_dFhatinv;
366 
367  const RankFourTensor dstrain_increment_dCtilde =
368  -0.5 * II_ijkl + 0.25 * (I.mixedProductIkJl(Ctilde) + Ctilde.mixedProductIkJl(I));
370  rot_rank_four * _Jacobian_mult[_qp] * dstrain_increment_dCtilde * dCtilde_dFhatinv;
371  _finite_deform_Jacobian_mult[_qp] += Fhat.mixedProductJkIl(_stress[_qp]);
372 
373  const RankFourTensor dFhat_dFhatinv = -Fhat.mixedProductIkJl(Fhat.transpose());
374  const RankTwoTensor dJ_dFhatinv = dFhat_dFhatinv.innerProductTranspose(Fhat.ddet());
375 
376  // Component from Jacobian derivative
377  _finite_deform_Jacobian_mult[_qp] += _stress[_qp].outerProduct(dJ_dFhatinv);
378 
379  // Derivative of Fhatinv w.r.t. undisplaced coordinates
380  const RankTwoTensor Finv = (*_deformation_gradient)[_qp].inverse();
381  const RankFourTensor dFhatinv_dGradu = -Fhatinv.mixedProductIkJl(Finv.transpose());
382  _finite_deform_Jacobian_mult[_qp] = _finite_deform_Jacobian_mult[_qp] * dFhatinv_dGradu;
383 }
const MaterialProperty< RankTwoTensor > * _rotation_increment
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
const MaterialProperty< RankFourTensor > & _Jacobian_mult
const MaterialProperty< RankTwoTensor > & _stress

◆ computeJacobian()

void StressDivergenceTensors::computeJacobian ( )
overridevirtualinherited

Reimplemented from ALEKernel.

Definition at line 154 of file StressDivergenceTensors.C.

155 {
157  {
160  }
161 
163  {
164  _finite_deform_Jacobian_mult.resize(_qrule->n_points());
165 
166  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
168 
170  }
171  else
172  Kernel::computeJacobian();
173 }
bool _volumetric_locking_correction
Flag for volumetric locking correction.
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
virtual void computeFiniteDeformJacobian()
virtual void computeJacobian() override
Definition: ALEKernel.C:35

◆ computeOffDiagJacobian()

void StressDivergenceTensors::computeOffDiagJacobian ( MooseVariableFEBase &  jvar)
overridevirtualinherited

Reimplemented from ALEKernel.

Definition at line 176 of file StressDivergenceTensors.C.

177 {
179  {
182  }
183 
185  {
186  _finite_deform_Jacobian_mult.resize(_qrule->n_points());
187 
188  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
190 
192  }
193  else
194  Kernel::computeOffDiagJacobian(jvar);
195 }
bool _volumetric_locking_correction
Flag for volumetric locking correction.
virtual void computeOffDiagJacobian(MooseVariableFEBase &jvar) override
Definition: ALEKernel.C:42
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
virtual void computeFiniteDeformJacobian()

◆ computeQpJacobian()

Real StressDivergenceRZTensors::computeQpJacobian ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 77 of file StressDivergenceRZTensors.C.

78 {
80 }
Real calculateJacobian(unsigned int ivar, unsigned int jvar)

◆ computeQpOffDiagJacobian()

Real StressDivergenceRZTensors::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 83 of file StressDivergenceRZTensors.C.

84 {
85  for (unsigned int i = 0; i < _ndisp; ++i)
86  {
87  if (jvar == _disp_var[i])
88  return calculateJacobian(_component, i);
89  }
90 
91  if (_temp_coupled && jvar == _temp_var)
92  {
93  Real jac = 0.0;
94  if (_component == 0)
95  {
96  for (unsigned k = 0; k < LIBMESH_DIM; ++k)
97  for (unsigned l = 0; l < LIBMESH_DIM; ++l)
98  jac -= (_grad_test[_i][_qp](0) * _Jacobian_mult[_qp](0, 0, k, l) +
99  _test[_i][_qp] / _q_point[_qp](0) * _Jacobian_mult[_qp](2, 2, k, l) +
100  _grad_test[_i][_qp](1) * _Jacobian_mult[_qp](0, 1, k, l)) *
101  (*_deigenstrain_dT)[_qp](k, l);
102  return jac * _phi[_j][_qp];
103  }
104  else if (_component == 1)
105  {
106  for (unsigned k = 0; k < LIBMESH_DIM; ++k)
107  for (unsigned l = 0; l < LIBMESH_DIM; ++l)
108  jac -= (_grad_test[_i][_qp](1) * _Jacobian_mult[_qp](1, 1, k, l) +
109  _grad_test[_i][_qp](0) * _Jacobian_mult[_qp](1, 0, k, l)) *
110  (*_deigenstrain_dT)[_qp](k, l);
111  return jac * _phi[_j][_qp];
112  }
113  }
114 
115  return 0.0;
116 }
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)

◆ computeQpResidual()

Real StressDivergenceRZTensors::computeQpResidual ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 47 of file StressDivergenceRZTensors.C.

48 {
49  Real div = 0.0;
50  if (_component == 0)
51  {
52  div = _grad_test[_i][_qp](0) * _stress[_qp](0, 0) +
53  +(_test[_i][_qp] / _q_point[_qp](0)) * _stress[_qp](2, 2) +
54  +_grad_test[_i][_qp](1) * _stress[_qp](0, 1); // stress_{rz}
55 
56  // volumetric locking correction
58  div += (_avg_grad_test[_i][0] - _grad_test[_i][_qp](0) - _test[_i][_qp] / _q_point[_qp](0)) *
59  (_stress[_qp].trace()) / 3.0;
60  }
61  else if (_component == 1)
62  {
63  div = _grad_test[_i][_qp](1) * _stress[_qp](1, 1) +
64  +_grad_test[_i][_qp](0) * _stress[_qp](1, 0); // stress_{zr}
65 
66  // volumetric locking correction
68  div += (_avg_grad_test[_i][1] - _grad_test[_i][_qp](1)) * (_stress[_qp].trace()) / 3.0;
69  }
70  else
71  mooseError("Invalid component for this AxisymmetricRZ problem.");
72 
73  return div;
74 }
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

◆ computeResidual()

void StressDivergenceTensors::computeResidual ( )
overrideprotectedvirtualinherited

Definition at line 117 of file StressDivergenceTensors.C.

118 {
119  DenseVector<Number> & re = _assembly.residualBlock(_var.number());
120  _local_re.resize(re.size());
121  _local_re.zero();
122 
125 
126  precalculateResidual();
127  for (_i = 0; _i < _test.size(); ++_i)
128  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
129  _local_re(_i) += _JxW[_qp] * _coord[_qp] * computeQpResidual();
130 
131  re += _local_re;
132 
133  if (_has_save_in)
134  {
135  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
136  for (const auto & var : _save_in)
137  var->sys().solution().add_vector(_local_re, var->dofIndices());
138  }
139 }
bool _volumetric_locking_correction
Flag for volumetric locking correction.
virtual Real computeQpResidual() override

◆ initialSetup()

void StressDivergenceRZTensors::initialSetup ( )
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 39 of file StressDivergenceRZTensors.C.

40 {
41  if (getBlockCoordSystem() != Moose::COORD_RZ)
42  mooseError("The coordinate system in the Problem block must be set to RZ for axisymmetric "
43  "geometries.");
44 }

Member Data Documentation

◆ _assembly_undisplaced

Assembly& ALEKernel::_assembly_undisplaced
protectedinherited

undisplaced problem

Definition at line 33 of file ALEKernel.h.

◆ _avg_grad_phi

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 83 of file StressDivergenceTensors.h.

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

◆ _avg_grad_test

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

◆ _base_name

std::string StressDivergenceTensors::_base_name
protectedinherited

◆ _component

const unsigned int StressDivergenceTensors::_component
protectedinherited

◆ _deformation_gradient

const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_deformation_gradient
protectedinherited

◆ _deformation_gradient_old

const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_deformation_gradient_old
protectedinherited

◆ _deigenstrain_dT

const MaterialProperty<RankTwoTensor>* const StressDivergenceTensors::_deigenstrain_dT
protectedinherited

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

Definition at line 73 of file StressDivergenceTensors.h.

◆ _disp_var

std::vector<unsigned int> StressDivergenceTensors::_disp_var
protectedinherited

◆ _finite_deform_Jacobian_mult

std::vector<RankFourTensor> StressDivergenceTensors::_finite_deform_Jacobian_mult
protectedinherited

◆ _grad_phi_undisplaced

const VariablePhiGradient& ALEKernel::_grad_phi_undisplaced
protectedinherited

Shape and test functions on the undisplaced mesh.

Definition at line 39 of file ALEKernel.h.

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

◆ _grad_test_undisplaced

const VariableTestGradient& ALEKernel::_grad_test_undisplaced
protectedinherited

Definition at line 40 of file ALEKernel.h.

◆ _Jacobian_mult

const MaterialProperty<RankFourTensor>& StressDivergenceTensors::_Jacobian_mult
protectedinherited

◆ _ndisp

unsigned int StressDivergenceTensors::_ndisp
protectedinherited

◆ _out_of_plane_direction

const unsigned int StressDivergenceTensors::_out_of_plane_direction
protectedinherited

◆ _out_of_plane_strain_coupled

const bool StressDivergenceTensors::_out_of_plane_strain_coupled
protectedinherited

◆ _out_of_plane_strain_var

const unsigned int StressDivergenceTensors::_out_of_plane_strain_var
protectedinherited

◆ _rotation_increment

const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_rotation_increment
protectedinherited

◆ _stress

const MaterialProperty<RankTwoTensor>& StressDivergenceTensors::_stress
protectedinherited

◆ _temp_coupled

const bool StressDivergenceTensors::_temp_coupled
protectedinherited

◆ _temp_var

const unsigned int StressDivergenceTensors::_temp_var
protectedinherited

◆ _use_finite_deform_jacobian

bool StressDivergenceTensors::_use_finite_deform_jacobian
protectedinherited

◆ _var_undisplaced

MooseVariable& ALEKernel::_var_undisplaced
protectedinherited

Reference to this Kernel's undisplaced MooseVariable object.

Definition at line 36 of file ALEKernel.h.

◆ _volumetric_locking_correction

bool StressDivergenceTensors::_volumetric_locking_correction
protectedinherited

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