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

StressDivergenceTensors mostly copies from StressDivergence. More...

#include <StressDivergenceTensors.h>

Inheritance diagram for StressDivergenceTensors:
[legend]

Public Member Functions

 StressDivergenceTensors (const InputParameters &parameters)
 
virtual void computeJacobian () override
 
virtual void computeOffDiagJacobian (MooseVariableFEBase &jvar) override
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

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

Protected Attributes

const 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 VariableValue * _out_of_plane_strain
 
const unsigned int _out_of_plane_strain_var
 
const unsigned int _out_of_plane_direction
 
const bool _use_displaced_mesh
 Whether this object is acting on the displaced mesh. 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

StressDivergenceTensors mostly copies from StressDivergence.

There are small changes to use RankFourTensor and RankTwoTensors instead of SymmElasticityTensors and SymmTensors. This is done to allow for more mathematical transparancy.

Definition at line 27 of file StressDivergenceTensors.h.

Constructor & Destructor Documentation

◆ StressDivergenceTensors()

StressDivergenceTensors::StressDivergenceTensors ( const InputParameters &  parameters)

Definition at line 62 of file StressDivergenceTensors.C.

63  : ALEKernel(parameters),
64  _base_name(isParamValid("base_name") ? getParam<std::string>("base_name") + "_" : ""),
65  _use_finite_deform_jacobian(getParam<bool>("use_finite_deform_jacobian")),
66  _stress(getMaterialPropertyByName<RankTwoTensor>(_base_name + "stress")),
67  _Jacobian_mult(getMaterialPropertyByName<RankFourTensor>(_base_name + "Jacobian_mult")),
68  _component(getParam<unsigned int>("component")),
69  _ndisp(coupledComponents("displacements")),
71  _temp_coupled(isCoupled("temperature")),
72  _temp_var(_temp_coupled ? coupled("temperature") : 0),
73  _deigenstrain_dT(_temp_coupled ? &getMaterialPropertyDerivative<RankTwoTensor>(
74  getParam<std::string>("thermal_eigenstrain_name"),
75  getVar("temperature", 0)->name())
76  : nullptr),
77  _out_of_plane_strain_coupled(isCoupled("out_of_plane_strain")),
78  _out_of_plane_strain(_out_of_plane_strain_coupled ? &coupledValue("out_of_plane_strain")
79  : nullptr),
80  _out_of_plane_strain_var(_out_of_plane_strain_coupled ? coupled("out_of_plane_strain") : 0),
81  _out_of_plane_direction(getParam<MooseEnum>("out_of_plane_direction")),
82  _use_displaced_mesh(getParam<bool>("use_displaced_mesh")),
83  _avg_grad_test(_test.size(), std::vector<Real>(3, 0.0)),
84  _avg_grad_phi(_phi.size(), std::vector<Real>(3, 0.0)),
85  _volumetric_locking_correction(getParam<bool>("volumetric_locking_correction"))
86 {
87  for (unsigned int i = 0; i < _ndisp; ++i)
88  _disp_var[i] = coupled("displacements", i);
89 
90  // Checking for consistency between mesh size and length of the provided displacements vector
91  if (_out_of_plane_direction != 2 && _ndisp != 3)
92  mooseError("For 2D simulations where the out-of-plane direction is x or y coordinate "
93  "directions the number of supplied displacements must be three.");
94  else if (_out_of_plane_direction == 2 && _ndisp != _mesh.dimension())
95  mooseError("The number of displacement variables supplied must match the mesh dimension");
96 
98  {
100  &getMaterialProperty<RankTwoTensor>(_base_name + "deformation_gradient");
102  &getMaterialPropertyOld<RankTwoTensor>(_base_name + "deformation_gradient");
103  _rotation_increment = &getMaterialProperty<RankTwoTensor>(_base_name + "rotation_increment");
104  }
105 
106  // Error if volumetric locking correction is turned on for 1D problems
108  mooseError("Volumetric locking correction should be set to false for 1-D problems.");
109 
110  // Generate warning when volumetric locking correction is used with second order elements
111  if (_mesh.hasSecondOrderElements() && _volumetric_locking_correction)
112  mooseWarning("Volumteric locking correction is not required for second order elements. Using "
113  "volumetric locking with second order elements could cause zigzag patterns in "
114  "stresses and strains.");
115 }

Member Function Documentation

◆ computeAverageGradientPhi()

void StressDivergenceTensors::computeAverageGradientPhi ( )
protectedvirtual

Reimplemented in StressDivergenceRZTensors.

Definition at line 423 of file StressDivergenceTensors.C.

424 {
425  // Calculate volume average derivatives for phi
426  _avg_grad_phi.resize(_phi.size());
427  for (_i = 0; _i < _phi.size(); ++_i)
428  {
429  _avg_grad_phi[_i].resize(3);
430  for (unsigned int component = 0; component < _mesh.dimension(); ++component)
431  {
432  _avg_grad_phi[_i][component] = 0.0;
433  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
434  _avg_grad_phi[_i][component] += _grad_phi[_i][_qp](component) * _JxW[_qp] * _coord[_qp];
435 
436  _avg_grad_phi[_i][component] /= _current_elem_volume;
437  }
438  }
439 }

Referenced by computeJacobian(), and computeOffDiagJacobian().

◆ computeAverageGradientTest()

void StressDivergenceTensors::computeAverageGradientTest ( )
protectedvirtual

Reimplemented in StressDivergenceRZTensors.

Definition at line 407 of file StressDivergenceTensors.C.

408 {
409  // Calculate volume averaged value of shape function derivative
410  _avg_grad_test.resize(_test.size());
411  for (_i = 0; _i < _test.size(); ++_i)
412  {
413  _avg_grad_test[_i].resize(3);
414  _avg_grad_test[_i][_component] = 0.0;
415  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
416  _avg_grad_test[_i][_component] += _grad_test[_i][_qp](_component) * _JxW[_qp] * _coord[_qp];
417 
418  _avg_grad_test[_i][_component] /= _current_elem_volume;
419  }
420 }

Referenced by computeJacobian(), computeOffDiagJacobian(), and computeResidual().

◆ computeFiniteDeformJacobian()

void StressDivergenceTensors::computeFiniteDeformJacobian ( )
protectedvirtual

Definition at line 350 of file StressDivergenceTensors.C.

351 {
352  const RankTwoTensor I(RankTwoTensor::initIdentity);
353  const RankFourTensor II_ijkl = I.mixedProductIkJl(I);
354 
355  // Bring back to unrotated config
356  const RankTwoTensor unrotated_stress =
357  (*_rotation_increment)[_qp].transpose() * _stress[_qp] * (*_rotation_increment)[_qp];
358 
359  // Incremental deformation gradient Fhat
360  const RankTwoTensor Fhat =
361  (*_deformation_gradient)[_qp] * (*_deformation_gradient_old)[_qp].inverse();
362  const RankTwoTensor Fhatinv = Fhat.inverse();
363 
364  const RankTwoTensor rot_times_stress = (*_rotation_increment)[_qp] * unrotated_stress;
365  const RankFourTensor dstress_drot =
366  I.mixedProductIkJl(rot_times_stress) + I.mixedProductJkIl(rot_times_stress);
367  const RankFourTensor rot_rank_four =
368  (*_rotation_increment)[_qp].mixedProductIkJl((*_rotation_increment)[_qp]);
369  const RankFourTensor drot_dUhatinv = Fhat.mixedProductIkJl(I);
370 
371  const RankTwoTensor A = I - Fhatinv;
372 
373  // Ctilde = Chat^-1 - I
374  const RankTwoTensor Ctilde = A * A.transpose() - A - A.transpose();
375  const RankFourTensor dCtilde_dFhatinv =
376  -I.mixedProductIkJl(A) - I.mixedProductJkIl(A) + II_ijkl + I.mixedProductJkIl(I);
377 
378  // Second order approximation of Uhat - consistent with strain increment definition
379  // const RankTwoTensor Uhat = I - 0.5 * Ctilde - 3.0/8.0 * Ctilde * Ctilde;
380 
381  RankFourTensor dUhatinv_dCtilde =
382  0.5 * II_ijkl - 1.0 / 8.0 * (I.mixedProductIkJl(Ctilde) + Ctilde.mixedProductIkJl(I));
383  RankFourTensor drot_dFhatinv = drot_dUhatinv * dUhatinv_dCtilde * dCtilde_dFhatinv;
384 
385  drot_dFhatinv -= Fhat.mixedProductIkJl((*_rotation_increment)[_qp].transpose());
386  _finite_deform_Jacobian_mult[_qp] = dstress_drot * drot_dFhatinv;
387 
388  const RankFourTensor dstrain_increment_dCtilde =
389  -0.5 * II_ijkl + 0.25 * (I.mixedProductIkJl(Ctilde) + Ctilde.mixedProductIkJl(I));
391  rot_rank_four * _Jacobian_mult[_qp] * dstrain_increment_dCtilde * dCtilde_dFhatinv;
392  _finite_deform_Jacobian_mult[_qp] += Fhat.mixedProductJkIl(_stress[_qp]);
393 
394  const RankFourTensor dFhat_dFhatinv = -Fhat.mixedProductIkJl(Fhat.transpose());
395  const RankTwoTensor dJ_dFhatinv = dFhat_dFhatinv.innerProductTranspose(Fhat.ddet());
396 
397  // Component from Jacobian derivative
398  _finite_deform_Jacobian_mult[_qp] += _stress[_qp].outerProduct(dJ_dFhatinv);
399 
400  // Derivative of Fhatinv w.r.t. undisplaced coordinates
401  const RankTwoTensor Finv = (*_deformation_gradient)[_qp].inverse();
402  const RankFourTensor dFhatinv_dGradu = -Fhatinv.mixedProductIkJl(Finv.transpose());
403  _finite_deform_Jacobian_mult[_qp] = _finite_deform_Jacobian_mult[_qp] * dFhatinv_dGradu;
404 }

Referenced by computeJacobian(), and computeOffDiagJacobian().

◆ computeJacobian()

void StressDivergenceTensors::computeJacobian ( )
overridevirtual

Reimplemented from ALEKernel.

Definition at line 168 of file StressDivergenceTensors.C.

169 {
171  {
174  }
175 
177  {
178  _finite_deform_Jacobian_mult.resize(_qrule->n_points());
179 
180  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
182 
184  }
185  else
186  Kernel::computeJacobian();
187 }

◆ computeOffDiagJacobian()

void StressDivergenceTensors::computeOffDiagJacobian ( MooseVariableFEBase &  jvar)
overridevirtual

Reimplemented from ALEKernel.

Definition at line 190 of file StressDivergenceTensors.C.

191 {
193  {
196  }
197 
199  {
200  _finite_deform_Jacobian_mult.resize(_qrule->n_points());
201 
202  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
204 
206  }
207  else
208  Kernel::computeOffDiagJacobian(jvar);
209 }

◆ computeQpJacobian()

Real StressDivergenceTensors::computeQpJacobian ( )
overrideprotectedvirtual

Reimplemented in StressDivergenceRSphericalTensors, StressDivergenceRZTensors, and DynamicStressDivergenceTensors.

Definition at line 212 of file StressDivergenceTensors.C.

213 {
216  _component,
217  _component,
218  _grad_test[_i][_qp],
219  _grad_phi_undisplaced[_j][_qp]);
220 
221  Real sum_C3x3 = _Jacobian_mult[_qp].sum3x3();
222  RealGradient sum_C3x1 = _Jacobian_mult[_qp].sum3x1();
223 
224  Real jacobian = 0.0;
225  // B^T_i * C * B_j
227  _Jacobian_mult[_qp], _component, _component, _grad_test[_i][_qp], _grad_phi[_j][_qp]);
228 
230  {
231  // jacobian = Bbar^T_i * C * Bbar_j where Bbar = B + Bvol
232  // jacobian = B^T_i * C * B_j + Bvol^T_i * C * Bvol_j + Bvol^T_i * C * B_j + B^T_i * C * Bvol_j
233 
234  // Bvol^T_i * C * Bvol_j
235  jacobian += sum_C3x3 * (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component)) *
236  (_avg_grad_phi[_j][_component] - _grad_phi[_j][_qp](_component)) / 9.0;
237 
238  // B^T_i * C * Bvol_j
239  jacobian += sum_C3x1(_component) * _grad_test[_i][_qp](_component) *
240  (_avg_grad_phi[_j][_component] - _grad_phi[_j][_qp](_component)) / 3.0;
241 
242  // Bvol^T_i * C * B_j
243  RankTwoTensor phi;
244  if (_component == 0)
245  {
246  phi(0, 0) = _grad_phi[_j][_qp](0);
247  phi(0, 1) = phi(1, 0) = _grad_phi[_j][_qp](1);
248  phi(0, 2) = phi(2, 0) = _grad_phi[_j][_qp](2);
249  }
250  else if (_component == 1)
251  {
252  phi(1, 1) = _grad_phi[_j][_qp](1);
253  phi(0, 1) = phi(1, 0) = _grad_phi[_j][_qp](0);
254  phi(1, 2) = phi(2, 1) = _grad_phi[_j][_qp](2);
255  }
256  else if (_component == 2)
257  {
258  phi(2, 2) = _grad_phi[_j][_qp](2);
259  phi(0, 2) = phi(2, 0) = _grad_phi[_j][_qp](0);
260  phi(1, 2) = phi(2, 1) = _grad_phi[_j][_qp](1);
261  }
262 
263  jacobian += (_Jacobian_mult[_qp] * phi).trace() *
264  (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component)) / 3.0;
265  }
266 
268  {
269  const Real out_of_plane_thickness = std::exp((*_out_of_plane_strain)[_qp]);
270  jacobian *= out_of_plane_thickness;
271  }
272 
273  return jacobian;
274 }

Referenced by DynamicStressDivergenceTensors::computeQpJacobian().

◆ computeQpOffDiagJacobian()

Real StressDivergenceTensors::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtual

Reimplemented in StressDivergenceRSphericalTensors, StressDivergenceRZTensors, CosseratStressDivergenceTensors, and DynamicStressDivergenceTensors.

Definition at line 277 of file StressDivergenceTensors.C.

278 {
279  // off-diagonal Jacobian with respect to a coupled displacement component
280  for (unsigned int coupled_component = 0; coupled_component < _ndisp; ++coupled_component)
281  if (jvar == _disp_var[coupled_component])
282  {
283  if (_out_of_plane_direction != 2)
284  {
285  if (coupled_component == _out_of_plane_direction)
286  continue;
287  }
288 
291  _component,
292  coupled_component,
293  _grad_test[_i][_qp],
294  _grad_phi_undisplaced[_j][_qp]);
295 
296  const Real sum_C3x3 = _Jacobian_mult[_qp].sum3x3();
297  const RealGradient sum_C3x1 = _Jacobian_mult[_qp].sum3x1();
298  Real jacobian = 0.0;
299 
300  // B^T_i * C * B_j
302  _component,
303  coupled_component,
304  _grad_test[_i][_qp],
305  _grad_phi[_j][_qp]);
306 
308  {
309  // jacobian = Bbar^T_i * C * Bbar_j where Bbar = B + Bvol
310  // jacobian = B^T_i * C * B_j + Bvol^T_i * C * Bvol_j + Bvol^T_i * C * B_j + B^T_i * C *
311  // Bvol_j
312 
313  // Bvol^T_i * C * Bvol_j
314  jacobian += sum_C3x3 * (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component)) *
315  (_avg_grad_phi[_j][coupled_component] - _grad_phi[_j][_qp](coupled_component)) /
316  9.0;
317 
318  // B^T_i * C * Bvol_j
319  jacobian += sum_C3x1(_component) * _grad_test[_i][_qp](_component) *
320  (_avg_grad_phi[_j][coupled_component] - _grad_phi[_j][_qp](coupled_component)) /
321  3.0;
322 
323  // Bvol^T_i * C * B_i
324  RankTwoTensor phi;
325  for (unsigned int i = 0; i < 3; ++i)
326  phi(coupled_component, i) = _grad_phi[_j][_qp](i);
327 
328  jacobian += (_Jacobian_mult[_qp] * phi).trace() *
329  (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component)) / 3.0;
330  }
331 
332  return jacobian;
333  }
334 
335  // off-diagonal Jacobian with respect to a coupled out_of_plane_strain variable
337  return _Jacobian_mult[_qp](
339  _grad_test[_i][_qp](_component) * _phi[_j][_qp];
340 
341  // off-diagonal Jacobian with respect to a coupled temperature variable
342  if (_temp_coupled && jvar == _temp_var)
343  return -((_Jacobian_mult[_qp] * (*_deigenstrain_dT)[_qp]) *
344  _grad_test[_i][_qp])(_component)*_phi[_j][_qp];
345 
346  return 0.0;
347 }

Referenced by CosseratStressDivergenceTensors::computeQpOffDiagJacobian(), and DynamicStressDivergenceTensors::computeQpOffDiagJacobian().

◆ computeQpResidual()

Real StressDivergenceTensors::computeQpResidual ( )
overrideprotectedvirtual

Reimplemented in StressDivergenceRSphericalTensors, StressDivergenceRZTensors, and DynamicStressDivergenceTensors.

Definition at line 149 of file StressDivergenceTensors.C.

150 {
151 
152  Real residual = _stress[_qp].row(_component) * _grad_test[_i][_qp];
153  // volumetric locking correction
155  residual += _stress[_qp].trace() / 3.0 *
156  (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component));
157 
159  {
160  const Real out_of_plane_thickness = std::exp((*_out_of_plane_strain)[_qp]);
161  residual *= out_of_plane_thickness;
162  }
163 
164  return residual;
165 }

Referenced by computeResidual().

◆ computeResidual()

void StressDivergenceTensors::computeResidual ( )
overrideprotectedvirtual

Definition at line 126 of file StressDivergenceTensors.C.

127 {
128  prepareVectorTag(_assembly, _var.number());
129 
132 
133  precalculateResidual();
134  for (_i = 0; _i < _test.size(); ++_i)
135  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
136  _local_re(_i) += _JxW[_qp] * _coord[_qp] * computeQpResidual();
137 
138  accumulateTaggedLocalResidual();
139 
140  if (_has_save_in)
141  {
142  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
143  for (const auto & var : _save_in)
144  var->sys().solution().add_vector(_local_re, var->dofIndices());
145  }
146 }

◆ initialSetup()

void StressDivergenceTensors::initialSetup ( )
overrideprotectedvirtual

Reimplemented in StressDivergenceRSphericalTensors, and StressDivergenceRZTensors.

Definition at line 118 of file StressDivergenceTensors.C.

119 {
120  if (getBlockCoordSystem() != Moose::COORD_XYZ)
121  mooseError(
122  "The coordinate system in the Problem block must be set to XYZ for cartesian geometries.");
123 }

◆ validParams()

InputParameters StressDivergenceTensors::validParams ( )
static

Definition at line 26 of file StressDivergenceTensors.C.

27 {
28  InputParameters params = ALEKernel::validParams();
29  params.addClassDescription("Stress divergence kernel for the Cartesian coordinate system");
30  params.addRequiredParam<unsigned int>("component",
31  "An integer corresponding to the direction "
32  "the variable this kernel acts in. (0 for x, "
33  "1 for y, 2 for z)");
34  params.addRequiredCoupledVar("displacements",
35  "The string of displacements suitable for the problem statement");
36  params.addCoupledVar("temperature",
37  "The name of the temperature variable used in the "
38  "ComputeThermalExpansionEigenstrain. (Not required for "
39  "simulations without temperature coupling.)");
40  params.addParam<std::string>(
41  "thermal_eigenstrain_name",
42  "thermal_eigenstrain",
43  "The eigenstrain_name used in the ComputeThermalExpansionEigenstrain.");
44  params.addCoupledVar("out_of_plane_strain",
45  "The name of the out_of_plane_strain variable used in the "
46  "WeakPlaneStress kernel.");
47  MooseEnum out_of_plane_direction("x y z", "z");
48  params.addParam<MooseEnum>(
49  "out_of_plane_direction",
50  out_of_plane_direction,
51  "The direction of the out_of_plane_strain variable used in the WeakPlaneStress kernel.");
52  params.addParam<std::string>("base_name", "Material property base name");
53  params.set<bool>("use_displaced_mesh") = false;
54  params.addParam<bool>(
55  "use_finite_deform_jacobian", false, "Jacobian for corotational finite strain");
56  params.addParam<bool>("volumetric_locking_correction",
57  false,
58  "Set to false to turn off volumetric locking correction");
59  return params;
60 }

Referenced by DynamicStressDivergenceTensors::validParams(), CosseratStressDivergenceTensors::validParams(), StressDivergenceRZTensors::validParams(), and StressDivergenceRSphericalTensors::validParams().

Member Data Documentation

◆ _assembly_undisplaced

Assembly& ALEKernel::_assembly_undisplaced
protectedinherited

undisplaced problem

Definition at line 34 of file ALEKernel.h.

◆ _avg_grad_phi

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

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

Definition at line 86 of file StressDivergenceTensors.h.

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

◆ _avg_grad_test

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

◆ _base_name

const std::string StressDivergenceTensors::_base_name
protected

Definition at line 50 of file StressDivergenceTensors.h.

Referenced by StressDivergenceTensors().

◆ _component

const unsigned int StressDivergenceTensors::_component
protected

◆ _deformation_gradient

const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_deformation_gradient
protected

Definition at line 57 of file StressDivergenceTensors.h.

Referenced by StressDivergenceTensors().

◆ _deformation_gradient_old

const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_deformation_gradient_old
protected

Definition at line 58 of file StressDivergenceTensors.h.

Referenced by StressDivergenceTensors().

◆ _deigenstrain_dT

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

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

Definition at line 72 of file StressDivergenceTensors.h.

◆ _disp_var

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

◆ _finite_deform_Jacobian_mult

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

◆ _grad_phi_undisplaced

const VariablePhiGradient& ALEKernel::_grad_phi_undisplaced
protectedinherited

Shape and test functions on the undisplaced mesh.

Definition at line 40 of file ALEKernel.h.

Referenced by computeQpJacobian(), and computeQpOffDiagJacobian().

◆ _grad_test_undisplaced

const VariableTestGradient& ALEKernel::_grad_test_undisplaced
protectedinherited

Definition at line 41 of file ALEKernel.h.

◆ _Jacobian_mult

const MaterialProperty<RankFourTensor>& StressDivergenceTensors::_Jacobian_mult
protected

◆ _ndisp

unsigned int StressDivergenceTensors::_ndisp
protected

◆ _out_of_plane_direction

const unsigned int StressDivergenceTensors::_out_of_plane_direction
protected

Definition at line 77 of file StressDivergenceTensors.h.

Referenced by computeQpOffDiagJacobian(), and StressDivergenceTensors().

◆ _out_of_plane_strain

const VariableValue* StressDivergenceTensors::_out_of_plane_strain
protected

Definition at line 75 of file StressDivergenceTensors.h.

Referenced by computeQpJacobian(), and computeQpResidual().

◆ _out_of_plane_strain_coupled

const bool StressDivergenceTensors::_out_of_plane_strain_coupled
protected

◆ _out_of_plane_strain_var

const unsigned int StressDivergenceTensors::_out_of_plane_strain_var
protected

Definition at line 76 of file StressDivergenceTensors.h.

Referenced by computeQpOffDiagJacobian().

◆ _rotation_increment

const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_rotation_increment
protected

◆ _stress

const MaterialProperty<RankTwoTensor>& StressDivergenceTensors::_stress
protected

◆ _temp_coupled

const bool StressDivergenceTensors::_temp_coupled
protected

◆ _temp_var

const unsigned int StressDivergenceTensors::_temp_var
protected

◆ _use_displaced_mesh

const bool StressDivergenceTensors::_use_displaced_mesh
protected

Whether this object is acting on the displaced mesh.

Definition at line 80 of file StressDivergenceTensors.h.

Referenced by computeQpJacobian(), and computeQpResidual().

◆ _use_finite_deform_jacobian

bool StressDivergenceTensors::_use_finite_deform_jacobian
protected

◆ _var_undisplaced

MooseVariable& ALEKernel::_var_undisplaced
protectedinherited

Reference to this Kernel's undisplaced MooseVariable object.

Definition at line 37 of file ALEKernel.h.

◆ _volumetric_locking_correction

bool StressDivergenceTensors::_volumetric_locking_correction
protected

The documentation for this class was generated from the following files:
ElasticityTensorTools::elasticJacobian
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...
Definition: ElasticityTensorTools.C:21
ALEKernel::computeOffDiagJacobian
virtual void computeOffDiagJacobian(MooseVariableFEBase &jvar) override
Definition: ALEKernel.C:43
StressDivergenceTensors::_deigenstrain_dT
const MaterialProperty< RankTwoTensor > *const _deigenstrain_dT
d(strain)/d(temperature), if computed by ComputeThermalExpansionEigenstrain
Definition: StressDivergenceTensors.h:72
ALEKernel::validParams
static InputParameters validParams()
Definition: ALEKernel.C:18
StressDivergenceTensors::_out_of_plane_strain
const VariableValue * _out_of_plane_strain
Definition: StressDivergenceTensors.h:75
StressDivergenceTensors::_component
const unsigned int _component
Definition: StressDivergenceTensors.h:62
StressDivergenceTensors::computeAverageGradientPhi
virtual void computeAverageGradientPhi()
Definition: StressDivergenceTensors.C:423
libMesh::RealGradient
VectorValue< Real > RealGradient
Definition: GrainForceAndTorqueInterface.h:17
StressDivergenceTensors::_disp_var
std::vector< unsigned int > _disp_var
Definition: StressDivergenceTensors.h:66
StressDivergenceTensors::_rotation_increment
const MaterialProperty< RankTwoTensor > * _rotation_increment
Definition: StressDivergenceTensors.h:59
StressDivergenceTensors::_temp_coupled
const bool _temp_coupled
Definition: StressDivergenceTensors.h:68
StressDivergenceTensors::_avg_grad_phi
std::vector< std::vector< Real > > _avg_grad_phi
Gradient of phi function averaged over the element. Used in volumetric locking correction calculation...
Definition: StressDivergenceTensors.h:86
ALEKernel::computeJacobian
virtual void computeJacobian() override
Definition: ALEKernel.C:36
StressDivergenceTensors::computeFiniteDeformJacobian
virtual void computeFiniteDeformJacobian()
Definition: StressDivergenceTensors.C:350
StressDivergenceTensors::_use_displaced_mesh
const bool _use_displaced_mesh
Whether this object is acting on the displaced mesh.
Definition: StressDivergenceTensors.h:80
StressDivergenceTensors::_volumetric_locking_correction
bool _volumetric_locking_correction
Flag for volumetric locking correction.
Definition: StressDivergenceTensors.h:89
StressDivergenceTensors::_out_of_plane_direction
const unsigned int _out_of_plane_direction
Definition: StressDivergenceTensors.h:77
StressDivergenceTensors::_out_of_plane_strain_coupled
const bool _out_of_plane_strain_coupled
Definition: StressDivergenceTensors.h:74
StressDivergenceTensors::computeQpResidual
virtual Real computeQpResidual() override
Definition: StressDivergenceTensors.C:149
StressDivergenceTensors::_out_of_plane_strain_var
const unsigned int _out_of_plane_strain_var
Definition: StressDivergenceTensors.h:76
name
const std::string name
Definition: Setup.h:21
MaterialTensorCalculatorTools::component
Real component(const SymmTensor &symm_tensor, unsigned int index)
Definition: MaterialTensorCalculatorTools.C:16
StressDivergenceTensors::_Jacobian_mult
const MaterialProperty< RankFourTensor > & _Jacobian_mult
Definition: StressDivergenceTensors.h:54
StressDivergenceTensors::_finite_deform_Jacobian_mult
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
Definition: StressDivergenceTensors.h:56
StressDivergenceTensors::_stress
const MaterialProperty< RankTwoTensor > & _stress
Definition: StressDivergenceTensors.h:53
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
StressDivergenceTensors::_avg_grad_test
std::vector< std::vector< Real > > _avg_grad_test
Gradient of test function averaged over the element. Used in volumetric locking correction calculatio...
Definition: StressDivergenceTensors.h:83
StressDivergenceTensors::computeAverageGradientTest
virtual void computeAverageGradientTest()
Definition: StressDivergenceTensors.C:407
StressDivergenceTensors::_ndisp
unsigned int _ndisp
Coupled displacement variables.
Definition: StressDivergenceTensors.h:65
StressDivergenceTensors::_use_finite_deform_jacobian
bool _use_finite_deform_jacobian
Definition: StressDivergenceTensors.h:51
ALEKernel::_grad_phi_undisplaced
const VariablePhiGradient & _grad_phi_undisplaced
Shape and test functions on the undisplaced mesh.
Definition: ALEKernel.h:40
RankTwoTensorTempl< Real >
StressDivergenceTensors::_deformation_gradient
const MaterialProperty< RankTwoTensor > * _deformation_gradient
Definition: StressDivergenceTensors.h:57
StressDivergenceTensors::_temp_var
const unsigned int _temp_var
Definition: StressDivergenceTensors.h:69
StressDivergenceTensors::_deformation_gradient_old
const MaterialProperty< RankTwoTensor > * _deformation_gradient_old
Definition: StressDivergenceTensors.h:58
StressDivergenceTensors::_base_name
const std::string _base_name
Definition: StressDivergenceTensors.h:50
ALEKernel::ALEKernel
ALEKernel(const InputParameters &parameters)
Definition: ALEKernel.C:25