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 122 of file StressDivergenceRZTensors.C.

Referenced by computeQpJacobian(), and computeQpOffDiagJacobian().

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

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

239 {
240  // calculate volume averaged value of shape function derivative
241  _avg_grad_test.resize(_test.size());
242  for (_i = 0; _i < _test.size(); ++_i)
243  {
244  _avg_grad_test[_i].resize(2);
245  _avg_grad_test[_i][_component] = 0.0;
246  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
247  {
248  if (_component == 0)
249  _avg_grad_test[_i][_component] +=
250  (_grad_test[_i][_qp](_component) + _test[_i][_qp] / _q_point[_qp](0)) * _JxW[_qp] *
251  _coord[_qp];
252  else
253  _avg_grad_test[_i][_component] += _grad_test[_i][_qp](_component) * _JxW[_qp] * _coord[_qp];
254  }
255  _avg_grad_test[_i][_component] /= _current_elem_volume;
256  }
257 }
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 335 of file StressDivergenceTensors.C.

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

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

161 {
163  {
166  }
167 
169  {
170  _finite_deform_Jacobian_mult.resize(_qrule->n_points());
171 
172  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
174 
176  }
177  else
178  Kernel::computeJacobian();
179 }
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 182 of file StressDivergenceTensors.C.

183 {
185  {
188  }
189 
191  {
192  _finite_deform_Jacobian_mult.resize(_qrule->n_points());
193 
194  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
196 
198  }
199  else
200  Kernel::computeOffDiagJacobian(jvar);
201 }
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 80 of file StressDivergenceRZTensors.C.

81 {
83 }
Real calculateJacobian(unsigned int ivar, unsigned int jvar)

◆ computeQpOffDiagJacobian()

Real StressDivergenceRZTensors::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtual

Reimplemented from StressDivergenceTensors.

Definition at line 86 of file StressDivergenceRZTensors.C.

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

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

124 {
125  DenseVector<Number> & re = _assembly.residualBlock(_var.number());
126  _local_re.resize(re.size());
127  _local_re.zero();
128 
131 
132  precalculateResidual();
133  for (_i = 0; _i < _test.size(); ++_i)
134  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
135  _local_re(_i) += _JxW[_qp] * _coord[_qp] * computeQpResidual();
136 
137  re += _local_re;
138 
139  if (_has_save_in)
140  {
141  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
142  for (const auto & var : _save_in)
143  var->sys().solution().add_vector(_local_re, var->dofIndices());
144  }
145 }
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 
45  if (getBlockCoordSystem() == Moose::COORD_RZ && _fe_problem.getAxisymmetricRadialCoord() != 0)
46  mooseError("rz_coord_axis=Y is the only supported option for StressDivergenceRZTensors");
47 }

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 81 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 71 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: