www.mooseframework.org
StressDivergenceRZTensors.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "Assembly.h"
12 #include "ElasticityTensorTools.h"
13 #include "libmesh/quadrature.h"
14 
15 registerMooseObject("TensorMechanicsApp", StressDivergenceRZTensors);
16 
18 
19 InputParameters
21 {
22  InputParameters params = StressDivergenceTensors::validParams();
23  params.addClassDescription(
24  "Calculate stress divergence for an axisymmetric problem in cylindrical coordinates.");
25  params.addRequiredRangeCheckedParam<unsigned int>(
26  "component",
27  "component < 2",
28  "An integer corresponding to the direction the variable this kernel acts in. (0 "
29  "refers to the radial and 1 to the axial displacement.)");
30  params.set<bool>("use_displaced_mesh") = true;
31  return params;
32 }
33 
34 StressDivergenceRZTensors::StressDivergenceRZTensors(const InputParameters & parameters)
35  : StressDivergenceTensors(parameters)
36 {
37 }
38 
39 void
41 {
42  if (getBlockCoordSystem() != Moose::COORD_RZ)
43  mooseError("The coordinate system in the Problem block must be set to RZ for axisymmetric "
44  "geometries.");
45 
46  if (getBlockCoordSystem() == Moose::COORD_RZ && _fe_problem.getAxisymmetricRadialCoord() != 0)
47  mooseError("rz_coord_axis=Y is the only supported option for StressDivergenceRZTensors");
48 }
49 
50 Real
52 {
53  Real div = 0.0;
54  if (_component == 0)
55  {
56  div = _grad_test[_i][_qp](0) * _stress[_qp](0, 0) +
57  +(_test[_i][_qp] / _q_point[_qp](0)) * _stress[_qp](2, 2) +
58  +_grad_test[_i][_qp](1) * _stress[_qp](0, 1); // stress_{rz}
59 
60  // volumetric locking correction
62  div += (_avg_grad_test[_i][0] - _grad_test[_i][_qp](0) - _test[_i][_qp] / _q_point[_qp](0)) *
63  (_stress[_qp].trace()) / 3.0;
64  }
65  else if (_component == 1)
66  {
67  div = _grad_test[_i][_qp](1) * _stress[_qp](1, 1) +
68  +_grad_test[_i][_qp](0) * _stress[_qp](1, 0); // stress_{zr}
69 
70  // volumetric locking correction
72  div += (_avg_grad_test[_i][1] - _grad_test[_i][_qp](1)) * (_stress[_qp].trace()) / 3.0;
73  }
74  else
75  mooseError("Invalid component for this AxisymmetricRZ problem.");
76 
77  return div;
78 }
79 
80 Real
82 {
84 }
85 
86 Real
88 {
89  for (unsigned int i = 0; i < _ndisp; ++i)
90  {
91  if (jvar == _disp_var[i])
92  return calculateJacobian(_component, i);
93  }
94 
95  if (_temp_coupled && jvar == _temp_var)
96  {
97  Real jac = 0.0;
98  if (_component == 0)
99  {
100  for (unsigned k = 0; k < LIBMESH_DIM; ++k)
101  for (unsigned l = 0; l < LIBMESH_DIM; ++l)
102  jac -= (_grad_test[_i][_qp](0) * _Jacobian_mult[_qp](0, 0, k, l) +
103  _test[_i][_qp] / _q_point[_qp](0) * _Jacobian_mult[_qp](2, 2, k, l) +
104  _grad_test[_i][_qp](1) * _Jacobian_mult[_qp](0, 1, k, l)) *
105  (*_deigenstrain_dT)[_qp](k, l);
106  return jac * _phi[_j][_qp];
107  }
108  else if (_component == 1)
109  {
110  for (unsigned k = 0; k < LIBMESH_DIM; ++k)
111  for (unsigned l = 0; l < LIBMESH_DIM; ++l)
112  jac -= (_grad_test[_i][_qp](1) * _Jacobian_mult[_qp](1, 1, k, l) +
113  _grad_test[_i][_qp](0) * _Jacobian_mult[_qp](1, 0, k, l)) *
114  (*_deigenstrain_dT)[_qp](k, l);
115  return jac * _phi[_j][_qp];
116  }
117  }
118 
119  return 0.0;
120 }
121 
122 Real
123 StressDivergenceRZTensors::calculateJacobian(unsigned int ivar, unsigned int jvar)
124 {
125  // B^T_i * C * B_j
126  RealGradient test, test_z, phi, phi_z;
127  Real first_term = 0.0;
128  if (ivar == 0) // Case grad_test for x, requires contributions from stress_xx, stress_xy, and
129  // stress_zz
130  {
131  test(0) = _grad_test[_i][_qp](0);
132  test(1) = _grad_test[_i][_qp](1);
133  test_z(2) = _test[_i][_qp] / _q_point[_qp](0);
134  }
135  else // Case grad_test for y
136  {
137  test(0) = _grad_test[_i][_qp](0);
138  test(1) = _grad_test[_i][_qp](1);
139  }
140 
141  if (jvar == 0)
142  {
143  phi(0) = _grad_phi[_j][_qp](0);
144  phi(1) = _grad_phi[_j][_qp](1);
145  phi_z(2) = _phi[_j][_qp] / _q_point[_qp](0);
146  }
147  else
148  {
149  phi(0) = _grad_phi[_j][_qp](0);
150  phi(1) = _grad_phi[_j][_qp](1);
151  }
152 
153  if (ivar == 0 &&
154  jvar == 0) // Case when both phi and test are functions of x and z; requires four terms
155  {
156  const Real first_sum = ElasticityTensorTools::elasticJacobian(
157  _Jacobian_mult[_qp], ivar, jvar, test, phi); // test_x and phi_x
158  const Real second_sum = ElasticityTensorTools::elasticJacobian(
159  _Jacobian_mult[_qp], 2, 2, test_z, phi_z); // test_z and phi_z
160  const Real mixed_sum1 = ElasticityTensorTools::elasticJacobian(
161  _Jacobian_mult[_qp], ivar, 2, test, phi_z); // test_x and phi_z
162  const Real mixed_sum2 = ElasticityTensorTools::elasticJacobian(
163  _Jacobian_mult[_qp], 2, jvar, test_z, phi); // test_z and phi_x
164 
165  first_term = first_sum + second_sum + mixed_sum1 + mixed_sum2;
166  }
167  else if (ivar == 0 && jvar == 1)
168  {
169  const Real first_sum = ElasticityTensorTools::elasticJacobian(
170  _Jacobian_mult[_qp], ivar, jvar, test, phi); // test_x and phi_y
171  const Real mixed_sum2 = ElasticityTensorTools::elasticJacobian(
172  _Jacobian_mult[_qp], 2, jvar, test_z, phi); // test_z and phi_y
173 
174  first_term = first_sum + mixed_sum2;
175  }
176  else if (ivar == 1 && jvar == 0)
177  {
178  const Real second_sum = ElasticityTensorTools::elasticJacobian(
179  _Jacobian_mult[_qp], ivar, jvar, test, phi); // test_y and phi_x
180  const Real mixed_sum1 = ElasticityTensorTools::elasticJacobian(
181  _Jacobian_mult[_qp], ivar, 2, test, phi_z); // test_y and phi_z
182 
183  first_term = second_sum + mixed_sum1;
184  }
185  else if (ivar == 1 && jvar == 1)
187  _Jacobian_mult[_qp], ivar, jvar, test, phi); // test_y and phi_y
188  else
189  mooseError("Invalid component in Jacobian Calculation");
190 
191  Real val = 0.0;
192  // volumetric locking correction
193  // K = Bbar^T_i * C * Bbar^T_j where Bbar = B + Bvol
194  // 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
196  {
197  RealGradient new_test(2, 0.0);
198  RealGradient new_phi(2, 0.0);
199 
200  new_test(0) = _grad_test[_i][_qp](0) + _test[_i][_qp] / _q_point[_qp](0);
201  new_test(1) = _grad_test[_i][_qp](1);
202  new_phi(0) = _grad_phi[_j][_qp](0) + _phi[_j][_qp] / _q_point[_qp](0);
203  new_phi(1) = _grad_phi[_j][_qp](1);
204 
205  // Bvol^T_i * C * Bvol_j
206  val += _Jacobian_mult[_qp].sum3x3() * (_avg_grad_test[_i][ivar] - new_test(ivar)) *
207  (_avg_grad_phi[_j][jvar] - new_phi(jvar)) / 3.0;
208 
209  // B^T_i * C * Bvol_j
210  RealGradient sum_3x1 = _Jacobian_mult[_qp].sum3x1();
211  if (ivar == 0 && jvar == 0)
212  val += (sum_3x1(0) * test(0) + sum_3x1(2) * test_z(2)) * (_avg_grad_phi[_j][0] - new_phi(0));
213  else if (ivar == 0 && jvar == 1)
214  val += (sum_3x1(0) * test(0) + sum_3x1(2) * test_z(2)) * (_avg_grad_phi[_j][1] - new_phi(1));
215  else if (ivar == 1 && jvar == 0)
216  val += sum_3x1(1) * test(1) * (_avg_grad_phi[_j][0] - new_phi(0));
217  else
218  val += sum_3x1(1) * test(1) * (_avg_grad_phi[_j][1] - new_phi(1));
219 
220  // Bvol^T_i * C * B_j
221  // val = trace (C * B_j) *(avg_grad_test[_i][ivar] - new_test(ivar))
222  if (jvar == 0)
223  for (unsigned int i = 0; i < 3; ++i)
224  val +=
225  (_Jacobian_mult[_qp](i, i, 0, 0) * phi(0) + _Jacobian_mult[_qp](i, i, 0, 1) * phi(1) +
226  _Jacobian_mult[_qp](i, i, 2, 2) * phi_z(2)) *
227  (_avg_grad_test[_i][ivar] - new_test(ivar));
228  else if (jvar == 1)
229  for (unsigned int i = 0; i < 3; ++i)
230  val +=
231  (_Jacobian_mult[_qp](i, i, 0, 1) * phi(0) + _Jacobian_mult[_qp](i, i, 1, 1) * phi(1)) *
232  (_avg_grad_test[_i][ivar] - new_test(ivar));
233  }
234 
235  return val / 3.0 + first_term;
236 }
237 
238 void
240 {
241  // calculate volume averaged value of shape function derivative
242  _avg_grad_test.resize(_test.size());
243  for (_i = 0; _i < _test.size(); ++_i)
244  {
245  _avg_grad_test[_i].resize(2);
246  _avg_grad_test[_i][_component] = 0.0;
247  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
248  {
249  if (_component == 0)
250  _avg_grad_test[_i][_component] +=
251  (_grad_test[_i][_qp](_component) + _test[_i][_qp] / _q_point[_qp](0)) * _JxW[_qp] *
252  _coord[_qp];
253  else
254  _avg_grad_test[_i][_component] += _grad_test[_i][_qp](_component) * _JxW[_qp] * _coord[_qp];
255  }
256  _avg_grad_test[_i][_component] /= _current_elem_volume;
257  }
258 }
259 
260 void
262 {
263  _avg_grad_phi.resize(_phi.size());
264  for (_i = 0; _i < _phi.size(); ++_i)
265  {
266  _avg_grad_phi[_i].resize(2);
267  for (unsigned int component = 0; component < 2; ++component)
268  {
269  _avg_grad_phi[_i][component] = 0.0;
270  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
271  {
272  if (component == 0)
273  _avg_grad_phi[_i][component] +=
274  (_grad_phi[_i][_qp](component) + _phi[_i][_qp] / _q_point[_qp](0)) * _JxW[_qp] *
275  _coord[_qp];
276  else
277  _avg_grad_phi[_i][component] += _grad_phi[_i][_qp](component) * _JxW[_qp] * _coord[_qp];
278  }
279  _avg_grad_phi[_i][component] /= _current_elem_volume;
280  }
281  }
282 }
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
StressDivergenceRZTensors::computeQpResidual
virtual Real computeQpResidual() override
Definition: StressDivergenceRZTensors.C:51
StressDivergenceTensors
StressDivergenceTensors mostly copies from StressDivergence.
Definition: StressDivergenceTensors.h:27
StressDivergenceTensors::_component
const unsigned int _component
Definition: StressDivergenceTensors.h:62
libMesh::RealGradient
VectorValue< Real > RealGradient
Definition: GrainForceAndTorqueInterface.h:17
StressDivergenceTensors::_disp_var
std::vector< unsigned int > _disp_var
Definition: StressDivergenceTensors.h:66
StressDivergenceRZTensors::computeAverageGradientTest
virtual void computeAverageGradientTest() override
Definition: StressDivergenceRZTensors.C:239
StressDivergenceRZTensors.h
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
StressDivergenceRZTensors
StressDivergenceRZTensors is a modification of StressDivergenceTensors to accommodate the Axisymmetri...
Definition: StressDivergenceRZTensors.h:31
StressDivergenceRZTensors::calculateJacobian
Real calculateJacobian(unsigned int ivar, unsigned int jvar)
Definition: StressDivergenceRZTensors.C:123
StressDivergenceRZTensors::computeQpOffDiagJacobian
virtual Real computeQpOffDiagJacobian(unsigned int jvar) override
Definition: StressDivergenceRZTensors.C:87
StressDivergenceTensors::_volumetric_locking_correction
bool _volumetric_locking_correction
Flag for volumetric locking correction.
Definition: StressDivergenceTensors.h:89
StressDivergenceRZTensors::validParams
static InputParameters validParams()
Definition: StressDivergenceRZTensors.C:20
defineLegacyParams
defineLegacyParams(StressDivergenceRZTensors)
registerMooseObject
registerMooseObject("TensorMechanicsApp", StressDivergenceRZTensors)
StressDivergenceRZTensors::computeAverageGradientPhi
virtual void computeAverageGradientPhi() override
Definition: StressDivergenceRZTensors.C:261
ElasticityTensorTools.h
StressDivergenceTensors::validParams
static InputParameters validParams()
Definition: StressDivergenceTensors.C:26
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
StressDivergenceRZTensors::StressDivergenceRZTensors
StressDivergenceRZTensors(const InputParameters &parameters)
Definition: StressDivergenceRZTensors.C:34
StressDivergenceTensors::_stress
const MaterialProperty< RankTwoTensor > & _stress
Definition: StressDivergenceTensors.h:53
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::_ndisp
unsigned int _ndisp
Coupled displacement variables.
Definition: StressDivergenceTensors.h:65
StressDivergenceRZTensors::computeQpJacobian
virtual Real computeQpJacobian() override
Definition: StressDivergenceRZTensors.C:81
StressDivergenceRZTensors::initialSetup
virtual void initialSetup() override
Definition: StressDivergenceRZTensors.C:40
StressDivergenceTensors::_temp_var
const unsigned int _temp_var
Definition: StressDivergenceTensors.h:69