www.mooseframework.org
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
 

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

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

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 30 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_var(_out_of_plane_strain_coupled ? coupled("out_of_plane_strain") : 0),
79  _out_of_plane_direction(getParam<MooseEnum>("out_of_plane_direction")),
80  _avg_grad_test(_test.size(), std::vector<Real>(3, 0.0)),
81  _avg_grad_phi(_phi.size(), std::vector<Real>(3, 0.0)),
82  _volumetric_locking_correction(getParam<bool>("volumetric_locking_correction"))
83 {
84  for (unsigned int i = 0; i < _ndisp; ++i)
85  _disp_var[i] = coupled("displacements", i);
86 
87  // Checking for consistency between mesh size and length of the provided displacements vector
88  if (_out_of_plane_direction != 2 && _ndisp != 3)
89  mooseError("For 2D simulations where the out-of-plane direction is x or y coordinate "
90  "directions the number of supplied displacements must be three.");
91  else if (_out_of_plane_direction == 2 && _ndisp != _mesh.dimension())
92  mooseError("The number of displacement variables supplied must match the mesh dimension");
93 
95  {
97  &getMaterialProperty<RankTwoTensor>(_base_name + "deformation_gradient");
99  &getMaterialPropertyOld<RankTwoTensor>(_base_name + "deformation_gradient");
100  _rotation_increment = &getMaterialProperty<RankTwoTensor>(_base_name + "rotation_increment");
101  }
102 
103  // Error if volumetic locking correction is turned on for 1D problems
105  mooseError("Volumetric locking correction should be set to false for 1-D problems.");
106 }
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.
unsigned int _ndisp
Coupled displacement variables.
std::vector< std::vector< Real > > _avg_grad_test
Gradient of test function averaged over the element. Used in volumetric locking correction calculatio...
const unsigned int _out_of_plane_strain_var
const MaterialProperty< RankTwoTensor > * _deformation_gradient_old
const unsigned int _out_of_plane_direction
const MaterialProperty< RankTwoTensor > * _rotation_increment
std::vector< unsigned int > _disp_var
ALEKernel(const InputParameters &parameters)
Definition: ALEKernel.C:24
const std::string name
Definition: Setup.h:22
const MaterialProperty< RankTwoTensor > * _deformation_gradient
const MaterialProperty< RankFourTensor > & _Jacobian_mult
const MaterialProperty< RankTwoTensor > & _stress
const MaterialProperty< RankTwoTensor > *const _deigenstrain_dT
d(strain)/d(temperature), if computed by ComputeThermalExpansionEigenstrain

Member Function Documentation

◆ computeAverageGradientPhi()

void StressDivergenceTensors::computeAverageGradientPhi ( )
protectedvirtual

Reimplemented in StressDivergenceRZTensors.

Definition at line 402 of file StressDivergenceTensors.C.

Referenced by computeJacobian(), and computeOffDiagJacobian().

403 {
404  // Calculate volume average derivatives for phi
405  _avg_grad_phi.resize(_phi.size());
406  for (_i = 0; _i < _phi.size(); ++_i)
407  {
408  _avg_grad_phi[_i].resize(3);
409  for (unsigned int component = 0; component < _mesh.dimension(); ++component)
410  {
411  _avg_grad_phi[_i][component] = 0.0;
412  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
413  _avg_grad_phi[_i][component] += _grad_phi[_i][_qp](component) * _JxW[_qp] * _coord[_qp];
414 
415  _avg_grad_phi[_i][component] /= _current_elem_volume;
416  }
417  }
418 }
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 StressDivergenceTensors::computeAverageGradientTest ( )
protectedvirtual

Reimplemented in StressDivergenceRZTensors.

Definition at line 386 of file StressDivergenceTensors.C.

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

387 {
388  // Calculate volume averaged value of shape function derivative
389  _avg_grad_test.resize(_test.size());
390  for (_i = 0; _i < _test.size(); ++_i)
391  {
392  _avg_grad_test[_i].resize(3);
393  _avg_grad_test[_i][_component] = 0.0;
394  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
395  _avg_grad_test[_i][_component] += _grad_test[_i][_qp](_component) * _JxW[_qp] * _coord[_qp];
396 
397  _avg_grad_test[_i][_component] /= _current_elem_volume;
398  }
399 }
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 ( )
protectedvirtual

Definition at line 329 of file StressDivergenceTensors.C.

Referenced by computeJacobian(), and 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 ( )
overridevirtual

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)
overridevirtual

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 StressDivergenceTensors::computeQpJacobian ( )
overrideprotectedvirtual

Reimplemented in StressDivergenceRSphericalTensors, StressDivergenceRZTensors, and DynamicStressDivergenceTensors.

Definition at line 198 of file StressDivergenceTensors.C.

Referenced by DynamicStressDivergenceTensors::computeQpJacobian().

199 {
202  _component,
203  _component,
204  _grad_test[_i][_qp],
205  _grad_phi_undisplaced[_j][_qp]);
206 
207  Real sum_C3x3 = _Jacobian_mult[_qp].sum3x3();
208  RealGradient sum_C3x1 = _Jacobian_mult[_qp].sum3x1();
209 
210  Real jacobian = 0.0;
211  // B^T_i * C * B_j
213  _Jacobian_mult[_qp], _component, _component, _grad_test[_i][_qp], _grad_phi[_j][_qp]);
214 
216  {
217  // jacobian = Bbar^T_i * C * Bbar_j where Bbar = B + Bvol
218  // 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
219 
220  // Bvol^T_i * C * Bvol_j
221  jacobian += sum_C3x3 * (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component)) *
222  (_avg_grad_phi[_j][_component] - _grad_phi[_j][_qp](_component)) / 9.0;
223 
224  // B^T_i * C * Bvol_j
225  jacobian += sum_C3x1(_component) * _grad_test[_i][_qp](_component) *
226  (_avg_grad_phi[_j][_component] - _grad_phi[_j][_qp](_component)) / 3.0;
227 
228  // Bvol^T_i * C * B_j
229  RankTwoTensor phi;
230  if (_component == 0)
231  {
232  phi(0, 0) = _grad_phi[_j][_qp](0);
233  phi(0, 1) = phi(1, 0) = _grad_phi[_j][_qp](1);
234  phi(0, 2) = phi(2, 0) = _grad_phi[_j][_qp](2);
235  }
236  else if (_component == 1)
237  {
238  phi(1, 1) = _grad_phi[_j][_qp](1);
239  phi(0, 1) = phi(1, 0) = _grad_phi[_j][_qp](0);
240  phi(1, 2) = phi(2, 1) = _grad_phi[_j][_qp](2);
241  }
242  else if (_component == 2)
243  {
244  phi(2, 2) = _grad_phi[_j][_qp](2);
245  phi(0, 2) = phi(2, 0) = _grad_phi[_j][_qp](0);
246  phi(1, 2) = phi(2, 1) = _grad_phi[_j][_qp](1);
247  }
248 
249  jacobian += (_Jacobian_mult[_qp] * phi).trace() *
250  (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component)) / 3.0;
251  }
252  return jacobian;
253 }
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...
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
const MaterialProperty< RankFourTensor > & _Jacobian_mult
const VariablePhiGradient & _grad_phi_undisplaced
Shape and test functions on the undisplaced mesh.
Definition: ALEKernel.h:39

◆ computeQpOffDiagJacobian()

Real StressDivergenceTensors::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtual

Reimplemented in StressDivergenceRSphericalTensors, StressDivergenceRZTensors, CosseratStressDivergenceTensors, and DynamicStressDivergenceTensors.

Definition at line 256 of file StressDivergenceTensors.C.

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

257 {
258  // off-diagonal Jacobian with respect to a coupled displacement component
259  for (unsigned int coupled_component = 0; coupled_component < _ndisp; ++coupled_component)
260  if (jvar == _disp_var[coupled_component])
261  {
262  if (_out_of_plane_direction != 2)
263  {
264  if (coupled_component == _out_of_plane_direction)
265  continue;
266  }
267 
270  _component,
271  coupled_component,
272  _grad_test[_i][_qp],
273  _grad_phi_undisplaced[_j][_qp]);
274 
275  const Real sum_C3x3 = _Jacobian_mult[_qp].sum3x3();
276  const RealGradient sum_C3x1 = _Jacobian_mult[_qp].sum3x1();
277  Real jacobian = 0.0;
278 
279  // B^T_i * C * B_j
281  _component,
282  coupled_component,
283  _grad_test[_i][_qp],
284  _grad_phi[_j][_qp]);
285 
287  {
288  // jacobian = Bbar^T_i * C * Bbar_j where Bbar = B + Bvol
289  // jacobian = B^T_i * C * B_j + Bvol^T_i * C * Bvol_j + Bvol^T_i * C * B_j + B^T_i * C *
290  // Bvol_j
291 
292  // Bvol^T_i * C * Bvol_j
293  jacobian += sum_C3x3 * (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component)) *
294  (_avg_grad_phi[_j][coupled_component] - _grad_phi[_j][_qp](coupled_component)) /
295  9.0;
296 
297  // B^T_i * C * Bvol_j
298  jacobian += sum_C3x1(_component) * _grad_test[_i][_qp](_component) *
299  (_avg_grad_phi[_j][coupled_component] - _grad_phi[_j][_qp](coupled_component)) /
300  3.0;
301 
302  // Bvol^T_i * C * B_i
303  RankTwoTensor phi;
304  for (unsigned int i = 0; i < 3; ++i)
305  phi(coupled_component, i) = _grad_phi[_j][_qp](i);
306 
307  jacobian += (_Jacobian_mult[_qp] * phi).trace() *
308  (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component)) / 3.0;
309  }
310 
311  return jacobian;
312  }
313 
314  // off-diagonal Jacobian with respect to a coupled out_of_plane_strain variable
316  return _Jacobian_mult[_qp](
318  _grad_test[_i][_qp](_component) * _phi[_j][_qp];
319 
320  // off-diagonal Jacobian with respect to a coupled temperature variable
321  if (_temp_coupled && jvar == _temp_var)
322  return -((_Jacobian_mult[_qp] * (*_deigenstrain_dT)[_qp]) *
323  _grad_test[_i][_qp])(_component)*_phi[_j][_qp];
324 
325  return 0.0;
326 }
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.
unsigned int _ndisp
Coupled displacement variables.
std::vector< std::vector< Real > > _avg_grad_test
Gradient of test function averaged over the element. Used in volumetric locking correction calculatio...
const unsigned int _out_of_plane_strain_var
const unsigned int _out_of_plane_direction
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...
std::vector< unsigned int > _disp_var
std::vector< RankFourTensor > _finite_deform_Jacobian_mult
const MaterialProperty< RankFourTensor > & _Jacobian_mult
const VariablePhiGradient & _grad_phi_undisplaced
Shape and test functions on the undisplaced mesh.
Definition: ALEKernel.h:39

◆ computeQpResidual()

Real StressDivergenceTensors::computeQpResidual ( )
overrideprotectedvirtual

Reimplemented in StressDivergenceRSphericalTensors, StressDivergenceRZTensors, and DynamicStressDivergenceTensors.

Definition at line 142 of file StressDivergenceTensors.C.

Referenced by computeResidual().

143 {
144  Real residual = _stress[_qp].row(_component) * _grad_test[_i][_qp];
145  // volumetric locking correction
147  residual += _stress[_qp].trace() / 3.0 *
148  (_avg_grad_test[_i][_component] - _grad_test[_i][_qp](_component));
149 
150  return residual;
151 }
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 ( )
overrideprotectedvirtual

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 StressDivergenceTensors::initialSetup ( )
overrideprotectedvirtual

Reimplemented in StressDivergenceRSphericalTensors, and StressDivergenceRZTensors.

Definition at line 109 of file StressDivergenceTensors.C.

110 {
111  if (getBlockCoordSystem() != Moose::COORD_XYZ)
112  mooseError(
113  "The coordinate system in the Problem block must be set to XYZ for cartesian geometries.");
114 }

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
protected

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

Definition at line 83 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

std::string StressDivergenceTensors::_base_name
protected

Definition at line 51 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 58 of file StressDivergenceTensors.h.

Referenced by StressDivergenceTensors().

◆ _deformation_gradient_old

const MaterialProperty<RankTwoTensor>* StressDivergenceTensors::_deformation_gradient_old
protected

Definition at line 59 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 73 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 39 of file ALEKernel.h.

Referenced by computeQpJacobian(), and 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
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_coupled

const bool StressDivergenceTensors::_out_of_plane_strain_coupled
protected

Definition at line 75 of file StressDivergenceTensors.h.

Referenced by computeQpOffDiagJacobian().

◆ _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_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 36 of file ALEKernel.h.

◆ _volumetric_locking_correction

bool StressDivergenceTensors::_volumetric_locking_correction
protected

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