13 #include "MooseMesh.h"
14 #include "RankTwoTensor.h"
15 #include "libmesh/string_to_enum.h"
16 #include "libmesh/quadrature.h"
17 #include "libmesh/utility.h"
22 return MooseEnum(
"Geometry Topology",
"Geometry");
28 return MooseEnum(
"KI KII KIII T",
"KI");
37 InputParameters params = validParams<ElementIntegralPostprocessor>();
38 params.addRequiredCoupledVar(
40 "The displacements appropriate for the simulation geometry and coordinate system");
41 params.addCoupledVar(
"temperature",
42 "The temperature (optional). Must be provided to correctly compute "
43 "stress intensity factors in models with thermal strain gradients.");
44 params.addRequiredParam<UserObjectName>(
"crack_front_definition",
45 "The CrackFrontDefinition user object name");
46 params.addParam<
unsigned int>(
47 "crack_front_point_index",
48 "The index of the point on the crack front corresponding to this q function");
49 params.addParam<Real>(
50 "K_factor",
"Conversion factor between interaction integral and stress intensity factor K");
51 params.addParam<
unsigned int>(
"symmetry_plane",
52 "Account for a symmetry plane passing through "
53 "the plane of the crack, normal to the specified "
54 "axis (0=x, 1=y, 2=z)");
55 params.addParam<Real>(
"poissons_ratio",
"Poisson's ratio for the material.");
56 params.addParam<Real>(
"youngs_modulus",
"Young's modulus of the material.");
57 params.set<
bool>(
"use_displaced_mesh") =
false;
58 params.addParam<
unsigned int>(
"ring_first",
59 "The first ring of elements for volume integral domain");
60 params.addParam<
unsigned int>(
"ring_index",
"Ring ID");
61 params.addParam<MooseEnum>(
"q_function_type",
63 "The method used to define the integration domain. Options are: " +
65 params.addRequiredParam<MooseEnum>(
"sif_mode",
67 "Stress intensity factor to calculate. Choices are: " +
74 : ElementIntegralPostprocessor(parameters),
75 _ndisp(coupledComponents(
"displacements")),
77 _has_crack_front_point_index(isParamValid(
"crack_front_point_index")),
78 _crack_front_point_index(
79 _has_crack_front_point_index ? getParam<unsigned int>(
"crack_front_point_index") : 0),
81 _stress(getMaterialPropertyByName<
SymmTensor>(
"stress")),
82 _strain(getMaterialPropertyByName<
SymmTensor>(
"elastic_strain")),
84 _has_temp(isCoupled(
"temperature")),
85 _grad_temp(_has_temp ? coupledGradient(
"temperature") : _grad_zero),
86 _current_instantaneous_thermal_expansion_coef(
87 hasMaterialProperty<Real>(
"current_instantaneous_thermal_expansion_coef")
88 ? &getMaterialProperty<Real>(
"current_instantaneous_thermal_expansion_coef")
90 _K_factor(getParam<Real>(
"K_factor")),
91 _has_symmetry_plane(isParamValid(
"symmetry_plane")),
92 _poissons_ratio(getParam<Real>(
"poissons_ratio")),
93 _youngs_modulus(getParam<Real>(
"youngs_modulus")),
94 _ring_index(getParam<unsigned int>(
"ring_index")),
95 _q_function_type(getParam<MooseEnum>(
"q_function_type").getEnum<
QMethod>()),
96 _sif_mode(getParam<MooseEnum>(
"sif_mode").getEnum<
SifMethod>())
99 mooseError(
"To include thermal strain term in interaction integral, must both couple "
100 "temperature in DomainIntegral block and compute thermal expansion property in "
101 "material model using compute_InteractionIntegral = true.");
108 if (
_ndisp != _mesh.dimension())
110 "The number of variables supplied in 'displacements' must match the mesh dimension.");
112 for (
unsigned int i = 0; i <
_ndisp; ++i)
113 _grad_disp[i] = &coupledGradient(
"displacements", i);
116 for (
unsigned i =
_ndisp; i < 3; ++i)
129 gatherSum(_integral_value);
142 RealVectorValue grad_q(0.0, 0.0, 0.0);
144 const std::vector<std::vector<Real>> & phi_curr_elem = *
_phi_curr_elem;
145 const std::vector<std::vector<RealGradient>> & dphi_curr_elem = *
_dphi_curr_elem;
147 for (
unsigned int i = 0; i < _current_elem->n_nodes(); ++i)
151 for (
unsigned int j = 0; j < _current_elem->dim(); ++j)
152 grad_q(j) += dphi_curr_elem[i][_qp](j) *
_q_curr_elem[i];
156 RealVectorValue crack_direction(0.0);
157 crack_direction(0) = 1.0;
160 Point p(_q_point[_qp]);
172 stress(0, 0) =
_stress[_qp].xx();
173 stress(0, 1) =
_stress[_qp].xy();
174 stress(0, 2) =
_stress[_qp].xz();
175 stress(1, 0) =
_stress[_qp].xy();
176 stress(1, 1) =
_stress[_qp].yy();
177 stress(1, 2) =
_stress[_qp].yz();
178 stress(2, 0) =
_stress[_qp].xz();
179 stress(2, 1) =
_stress[_qp].yz();
180 stress(2, 2) =
_stress[_qp].zz();
183 strain(0, 0) =
_strain[_qp].xx();
184 strain(0, 1) =
_strain[_qp].xy();
185 strain(0, 2) =
_strain[_qp].xz();
186 strain(1, 0) =
_strain[_qp].xy();
187 strain(1, 1) =
_strain[_qp].yy();
188 strain(1, 2) =
_strain[_qp].yz();
189 strain(2, 0) =
_strain[_qp].xz();
190 strain(2, 1) =
_strain[_qp].yz();
191 strain(2, 2) =
_strain[_qp].zz();
196 RealVectorValue grad_q_cf =
204 RealVectorValue grad_temp_cf =
208 dq(0, 0) = crack_direction(0) * grad_q_cf(0);
209 dq(0, 1) = crack_direction(0) * grad_q_cf(1);
210 dq(0, 2) = crack_direction(0) * grad_q_cf(2);
216 Real term1 = aux_du.doubleContraction(tmp1);
220 Real term2 = grad_disp_cf(0, 0) * tmp2(0, 0) + grad_disp_cf(1, 0) * tmp2(0, 1) +
221 grad_disp_cf(2, 0) * tmp2(0, 2);
224 Real term3 = dq(0, 0) * aux_stress.doubleContraction(strain_cf);
230 term4 = scalar_q * aux_stress.trace() * (*_current_instantaneous_thermal_expansion_coef)[_qp] *
233 Real q_avg_seg = 1.0;
242 Real eq = term1 + term2 - term3 + term4;
247 return eq / q_avg_seg;
256 FEType fe_type(Utility::string_to_enum<Order>(
"first"),
257 Utility::string_to_enum<FEFamily>(
"lagrange"));
258 const unsigned int dim = _current_elem->dim();
259 std::unique_ptr<FEBase> fe(FEBase::build(dim, fe_type));
260 fe->attach_quadrature_rule(const_cast<QBase *>(_qrule));
263 fe->reinit(_current_elem);
269 for (
unsigned int i = 0; i < _current_elem->n_nodes(); ++i)
271 const Node * this_node = _current_elem->node_ptr(i);
284 for (_qp = 0; _qp < _qrule->n_points(); _qp++)
292 RealVectorValue k(0.0);
302 Real tt2 = 3.0 *
_theta / 2.0;
303 Real st = std::sin(t);
304 Real ct = std::cos(t);
305 Real st2 = std::sin(t2);
306 Real ct2 = std::cos(t2);
307 Real stt2 = std::sin(tt2);
308 Real ctt2 = std::cos(tt2);
309 Real ct2sq = Utility::pow<2>(ct2);
310 Real ct2cu = Utility::pow<3>(ct2);
311 Real sqrt2PiR = std::sqrt(2.0 * libMesh::pi *
_r);
317 1.0 / sqrt2PiR * (k(0) * ct2 * (1.0 - st2 * stt2) - k(1) * st2 * (2.0 + ct2 * ctt2));
318 aux_stress(1, 1) = 1.0 / sqrt2PiR * (k(0) * ct2 * (1.0 + st2 * stt2) + k(1) * st2 * ct2 * ctt2);
319 aux_stress(0, 1) = 1.0 / sqrt2PiR * (k(0) * ct2 * st2 * ctt2 + k(1) * ct2 * (1.0 - st2 * stt2));
320 aux_stress(0, 2) = -1.0 / sqrt2PiR * k(2) * st2;
321 aux_stress(1, 2) = 1.0 / sqrt2PiR * k(2) * ct2;
325 aux_stress(2, 2) =
_poissons_ratio * (aux_stress(0, 0) + aux_stress(1, 1));
327 aux_stress(1, 0) = aux_stress(0, 1);
328 aux_stress(2, 0) = aux_stress(0, 2);
329 aux_stress(2, 1) = aux_stress(1, 2);
335 (ct * ct2 *
_kappa + ct * ct2 - 2.0 * ct * ct2cu + st * st2 *
_kappa +
336 st * st2 - 6.0 * st * st2 * ct2sq) +
338 (ct * st2 *
_kappa + ct * st2 + 2.0 * ct * st2 * ct2sq - st * ct2 *
_kappa +
339 3.0 * st * ct2 - 6.0 * st * ct2cu);
342 (ct * st2 *
_kappa + ct * st2 - 2.0 * ct * st2 * ct2sq - st * ct2 *
_kappa -
343 5.0 * st * ct2 + 6.0 * st * ct2cu) +
345 (-ct * ct2 *
_kappa + 3.0 * ct * ct2 - 2.0 * ct * ct2cu -
346 st * st2 *
_kappa + 3.0 * st * st2 - 6.0 * st * st2 * ct2sq);
348 grad_disp(0, 2) = k(2) / (
_shear_modulus * sqrt2PiR) * (st2 * ct - ct2 * st);
356 Real st = std::sin(t);
357 Real ct = std::cos(t);
358 Real stsq = Utility::pow<2>(st);
359 Real ctsq = Utility::pow<2>(ct);
360 Real ctcu = Utility::pow<3>(ct);
361 Real oneOverPiR = 1.0 / (libMesh::pi *
_r);
364 aux_stress(0, 0) = -oneOverPiR * ctcu;
365 aux_stress(0, 1) = -oneOverPiR * st * ctsq;
366 aux_stress(1, 0) = -oneOverPiR * st * ctsq;
367 aux_stress(1, 1) = -oneOverPiR * ct * stsq;