Line data Source code
1 : //* This file is part of the MOOSE framework
2 : //* https://mooseframework.inl.gov
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 :
10 : // MOOSE includes
11 : #include "DomainIntegralAction.h"
12 : #include "ExecFlagRegistry.h"
13 : #include "Factory.h"
14 : #include "FEProblem.h"
15 : #include "Moose.h"
16 : #include "Parser.h"
17 : #include "CrackFrontDefinition.h"
18 : #include "MooseMesh.h"
19 : #include "Conversion.h"
20 :
21 : #include "libmesh/string_to_enum.h"
22 :
23 : registerMooseAction("SolidMechanicsApp", DomainIntegralAction, "add_user_object");
24 :
25 : registerMooseAction("SolidMechanicsApp", DomainIntegralAction, "add_aux_variable");
26 :
27 : registerMooseAction("SolidMechanicsApp", DomainIntegralAction, "add_aux_kernel");
28 :
29 : registerMooseAction("SolidMechanicsApp", DomainIntegralAction, "add_postprocessor");
30 :
31 : registerMooseAction("SolidMechanicsApp", DomainIntegralAction, "add_material");
32 :
33 : InputParameters
34 2156 : DomainIntegralAction::validParams()
35 : {
36 2156 : InputParameters params = Action::validParams();
37 2156 : addCrackFrontDefinitionParams(params);
38 : MultiMooseEnum integral_vec(
39 : "JIntegral CIntegral KFromJIntegral InteractionIntegralKI InteractionIntegralKII "
40 2156 : "InteractionIntegralKIII InteractionIntegralT");
41 2156 : params.addClassDescription(
42 : "Creates the MOOSE objects needed to compute fraction domain integrals");
43 4312 : params.addRequiredParam<MultiMooseEnum>("integrals",
44 : integral_vec,
45 2156 : "Domain integrals to calculate. Choices are: " +
46 2156 : integral_vec.getRawNames());
47 4312 : params.addParam<std::vector<BoundaryName>>(
48 : "boundary", {}, "Boundary containing the crack front points");
49 4312 : params.addParam<std::vector<Point>>("crack_front_points", "Set of points to define crack front");
50 4312 : params.addParam<std::string>(
51 : "order", "FIRST", "Specifies the order of the FE shape function to use for q AuxVariables");
52 4312 : params.addParam<std::string>(
53 : "family", "LAGRANGE", "Specifies the family of FE shape functions to use for q AuxVariables");
54 4312 : params.addParam<std::vector<Real>>("radius_inner", "Inner radius for volume integral domain");
55 4312 : params.addParam<std::vector<Real>>("radius_outer", "Outer radius for volume integral domain");
56 4312 : params.addParam<unsigned int>("ring_first",
57 : "The first ring of elements for volume integral domain");
58 4312 : params.addParam<unsigned int>("ring_last",
59 : "The last ring of elements for volume integral domain");
60 4312 : params.addParam<std::vector<VariableName>>(
61 : "output_variable", "Variable values to be reported along the crack front");
62 4312 : params.addParam<Real>("poissons_ratio", "Poisson's ratio");
63 4312 : params.addParam<Real>("youngs_modulus", "Young's modulus");
64 4312 : params.addParam<std::vector<SubdomainName>>(
65 : "block", {}, "The block ids where integrals are defined");
66 4312 : params.addParam<std::vector<VariableName>>(
67 : "displacements",
68 : {},
69 : "The displacements appropriate for the simulation geometry and coordinate system");
70 4312 : params.addParam<VariableName>("temperature", "", "The temperature");
71 4312 : params.addParam<MaterialPropertyName>(
72 : "functionally_graded_youngs_modulus_crack_dir_gradient",
73 : "Gradient of the spatially varying Young's modulus provided in "
74 : "'functionally_graded_youngs_modulus' in the direction of crack extension.");
75 4312 : params.addParam<MaterialPropertyName>(
76 : "functionally_graded_youngs_modulus",
77 : "Spatially varying elasticity modulus variable. This input is required when "
78 : "using the functionally graded material capability.");
79 4312 : params.addCoupledVar("additional_eigenstrain_00",
80 : "Optional additional eigenstrain variable that will be accounted for in the "
81 : "interaction integral (component 00 or XX).");
82 4312 : params.addCoupledVar("additional_eigenstrain_01",
83 : "Optional additional eigenstrain variable that will be accounted for in the "
84 : "interaction integral (component 01 or XY).");
85 4312 : params.addCoupledVar("additional_eigenstrain_11",
86 : "Optional additional eigenstrain variable that will be accounted for in the "
87 : "interaction integral (component 11 or YY).");
88 4312 : params.addCoupledVar("additional_eigenstrain_22",
89 : "Optional additional eigenstrain variable that will be accounted for in the "
90 : "interaction integral (component 22 or ZZ).");
91 4312 : params.addCoupledVar("additional_eigenstrain_02",
92 : "Optional additional eigenstrain variable that will be accounted for in the "
93 : "interaction integral (component 02 or XZ).");
94 4312 : params.addCoupledVar("additional_eigenstrain_12",
95 : "Optional additional eigenstrain variable that will be accounted for in the "
96 : "interaction integral (component 12 or XZ).");
97 4312 : params.addParamNamesToGroup(
98 : "additional_eigenstrain_00 additional_eigenstrain_01 additional_eigenstrain_11 "
99 : "additional_eigenstrain_22 additional_eigenstrain_02 additional_eigenstrain_12",
100 : "Generic eigenstrains for the computation of the interaction integral.");
101 4312 : MooseEnum position_type("Angle Distance", "Distance");
102 4312 : params.addParam<MooseEnum>(
103 : "position_type",
104 : position_type,
105 2156 : "The method used to calculate position along crack front. Options are: " +
106 2156 : position_type.getRawNames());
107 4312 : MooseEnum q_function_type("Geometry Topology", "Geometry");
108 4312 : params.addParam<MooseEnum>("q_function_type",
109 : q_function_type,
110 2156 : "The method used to define the integration domain. Options are: " +
111 2156 : q_function_type.getRawNames());
112 4312 : params.addParam<bool>(
113 : "equivalent_k",
114 4312 : false,
115 : "Calculate an equivalent K from KI, KII and KIII, assuming self-similar crack growth.");
116 4312 : params.addParam<bool>("output_q", true, "Output q");
117 4312 : params.addRequiredParam<bool>(
118 : "incremental", "Flag to indicate whether an incremental or total model is being used.");
119 4312 : params.addParam<std::vector<MaterialPropertyName>>(
120 : "eigenstrain_names", "List of eigenstrains applied in the strain calculation");
121 6468 : params.addDeprecatedParam<bool>("convert_J_to_K",
122 4312 : false,
123 : "Convert J-integral to stress intensity factor K.",
124 : "This input parameter is deprecated and will be removed soon. "
125 : "Use 'integrals = KFromJIntegral' to request output of the "
126 : "conversion from the J-integral to stress intensity factors");
127 4312 : params.addParam<std::vector<MaterialName>>(
128 : "inelastic_models",
129 : "The material objects to use to calculate the strain energy rate density.");
130 4312 : params.addParam<MaterialPropertyName>("eigenstrain_gradient",
131 : "Material defining gradient of eigenstrain tensor");
132 4312 : params.addParam<MaterialPropertyName>("body_force", "Material defining body force");
133 4312 : params.addParam<bool>("use_automatic_differentiation",
134 4312 : false,
135 : "Flag to use automatic differentiation (AD) objects when possible");
136 4312 : params.addParam<bool>("output_vpp",
137 4312 : true,
138 : "Flag to control the vector postprocessor outputs. Select false to "
139 : "suppress the redundant csv files for each time step and ring");
140 2156 : return params;
141 2156 : }
142 :
143 2156 : DomainIntegralAction::DomainIntegralAction(const InputParameters & params)
144 : : Action(params),
145 4312 : _boundary_names(getParam<std::vector<BoundaryName>>("boundary")),
146 4312 : _closed_loop(getParam<bool>("closed_loop")),
147 2156 : _use_crack_front_points_provider(false),
148 4312 : _order(getParam<std::string>("order")),
149 4312 : _family(getParam<std::string>("family")),
150 4312 : _direction_method_moose_enum(getParam<MooseEnum>("crack_direction_method")),
151 4312 : _end_direction_method_moose_enum(getParam<MooseEnum>("crack_end_direction_method")),
152 4312 : _have_crack_direction_vector(isParamValid("crack_direction_vector")),
153 1752 : _crack_direction_vector(
154 2156 : _have_crack_direction_vector ? getParam<RealVectorValue>("crack_direction_vector") : 0.0),
155 4312 : _have_crack_direction_vector_end_1(isParamValid("crack_direction_vector_end_1")),
156 320 : _crack_direction_vector_end_1(_have_crack_direction_vector_end_1
157 2156 : ? getParam<RealVectorValue>("crack_direction_vector_end_1")
158 : : 0.0),
159 4312 : _have_crack_direction_vector_end_2(isParamValid("crack_direction_vector_end_2")),
160 320 : _crack_direction_vector_end_2(_have_crack_direction_vector_end_2
161 4312 : ? getParam<RealVectorValue>("crack_direction_vector_end_2")
162 : : 0.0),
163 4312 : _treat_as_2d(getParam<bool>("2d")),
164 4312 : _axis_2d(getParam<unsigned int>("axis_2d")),
165 4312 : _has_symmetry_plane(isParamValid("symmetry_plane")),
166 3628 : _symmetry_plane(_has_symmetry_plane ? getParam<unsigned int>("symmetry_plane")
167 : : std::numeric_limits<unsigned int>::max()),
168 4312 : _position_type(getParam<MooseEnum>("position_type")),
169 4312 : _q_function_type(getParam<MooseEnum>("q_function_type")),
170 4312 : _get_equivalent_k(getParam<bool>("equivalent_k")),
171 2156 : _use_displaced_mesh(false),
172 4312 : _output_q(getParam<bool>("output_q")),
173 4312 : _incremental(getParam<bool>("incremental")),
174 8624 : _convert_J_to_K(isParamValid("convert_J_to_K") ? getParam<bool>("convert_J_to_K") : false),
175 2156 : _fgm_crack(false),
176 6468 : _use_ad(getParam<bool>("use_automatic_differentiation"))
177 : {
178 :
179 4312 : if (isParamValid("functionally_graded_youngs_modulus_crack_dir_gradient") !=
180 4312 : isParamValid("functionally_graded_youngs_modulus"))
181 2 : paramError("functionally_graded_youngs_modulus_crack_dir_gradient",
182 : "You have selected to compute the interaction integral for a crack in FGM. That "
183 : "selection requires the user to provide a spatially varying elasticity modulus that "
184 : "defines the transition of material properties (i.e. "
185 : "'functionally_graded_youngs_modulus') and its "
186 : "spatial derivative in the crack direction (i.e. "
187 : "'functionally_graded_youngs_modulus_crack_dir_gradient').");
188 :
189 6462 : if (isParamValid("functionally_graded_youngs_modulus_crack_dir_gradient") &&
190 2322 : isParamValid("functionally_graded_youngs_modulus"))
191 : {
192 84 : _fgm_crack = true;
193 : _functionally_graded_youngs_modulus_crack_dir_gradient =
194 84 : getParam<MaterialPropertyName>("functionally_graded_youngs_modulus_crack_dir_gradient");
195 : _functionally_graded_youngs_modulus =
196 168 : getParam<MaterialPropertyName>("functionally_graded_youngs_modulus");
197 : }
198 :
199 2154 : if (_q_function_type == GEOMETRY)
200 : {
201 7776 : if (isParamValid("radius_inner") && isParamValid("radius_outer"))
202 : {
203 3888 : _radius_inner = getParam<std::vector<Real>>("radius_inner");
204 3888 : _radius_outer = getParam<std::vector<Real>>("radius_outer");
205 : }
206 : else
207 0 : mooseError("DomainIntegral error: must set radius_inner and radius_outer.");
208 9048 : for (unsigned int i = 0; i < _radius_inner.size(); ++i)
209 7104 : _ring_vec.push_back(i + 1);
210 : }
211 210 : else if (_q_function_type == TOPOLOGY)
212 : {
213 840 : if (isParamValid("ring_first") && isParamValid("ring_last"))
214 : {
215 420 : _ring_first = getParam<unsigned int>("ring_first");
216 420 : _ring_last = getParam<unsigned int>("ring_last");
217 : }
218 : else
219 0 : mooseError(
220 : "DomainIntegral error: must set ring_first and ring_last if q_function_type = Topology.");
221 924 : for (unsigned int i = _ring_first; i <= _ring_last; ++i)
222 714 : _ring_vec.push_back(i);
223 : }
224 : else
225 0 : paramError("q_function_type", "DomainIntegral error: invalid q_function_type.");
226 :
227 4308 : if (isParamValid("crack_front_points"))
228 612 : _crack_front_points = getParam<std::vector<Point>>("crack_front_points");
229 :
230 4308 : if (isParamValid("crack_front_points_provider"))
231 : {
232 0 : if (!isParamValid("number_points_from_provider"))
233 0 : paramError("number_points_from_provider",
234 : "DomainIntegral error: when crack_front_points_provider is used, "
235 : "number_points_from_provider must be provided.");
236 0 : _use_crack_front_points_provider = true;
237 0 : _crack_front_points_provider = getParam<UserObjectName>("crack_front_points_provider");
238 : }
239 4308 : else if (isParamValid("number_points_from_provider"))
240 0 : paramError("crack_front_points_provider",
241 : "DomainIntegral error: number_points_from_provider is provided but "
242 : "crack_front_points_provider cannot be found.");
243 4308 : if (isParamValid("crack_mouth_boundary"))
244 648 : _crack_mouth_boundary_names = getParam<std::vector<BoundaryName>>("crack_mouth_boundary");
245 4308 : if (isParamValid("intersecting_boundary"))
246 828 : _intersecting_boundary_names = getParam<std::vector<BoundaryName>>("intersecting_boundary");
247 2154 : if (_radius_inner.size() != _radius_outer.size())
248 0 : mooseError("Number of entries in 'radius_inner' and 'radius_outer' must match.");
249 :
250 : bool youngs_modulus_set(false);
251 : bool poissons_ratio_set(false);
252 :
253 : // All domain integral types block restrict the objects created by this action.
254 4308 : _blocks = getParam<std::vector<SubdomainName>>("block");
255 :
256 4308 : MultiMooseEnum integral_moose_enums = getParam<MultiMooseEnum>("integrals");
257 5554 : for (unsigned int i = 0; i < integral_moose_enums.size(); ++i)
258 : {
259 10200 : _displacements = getParam<std::vector<VariableName>>("displacements");
260 :
261 3400 : if (_displacements.size() < 2)
262 0 : paramError(
263 : "displacements",
264 : "DomainIntegral error: The size of the displacements vector should at least be 2.");
265 :
266 5248 : if (integral_moose_enums[i] != "JIntegral" && integral_moose_enums[i] != "CIntegral" &&
267 1848 : integral_moose_enums[i] != "KFromJIntegral")
268 : {
269 : // Check that parameters required for interaction integrals are defined
270 7224 : if (!(isParamValid("poissons_ratio")) || !(isParamValid("youngs_modulus")))
271 0 : mooseError(
272 : "DomainIntegral error: must set Poisson's ratio and Young's modulus for integral: ",
273 0 : integral_moose_enums[i]);
274 :
275 3612 : _poissons_ratio = getParam<Real>("poissons_ratio");
276 : poissons_ratio_set = true;
277 3612 : _youngs_modulus = getParam<Real>("youngs_modulus");
278 : youngs_modulus_set = true;
279 : }
280 :
281 3400 : _integrals.insert(INTEGRAL(int(integral_moose_enums.get(i))));
282 : }
283 :
284 2154 : if ((_integrals.count(J_INTEGRAL) != 0 && _integrals.count(C_INTEGRAL) != 0) ||
285 2154 : (_integrals.count(J_INTEGRAL) != 0 && _integrals.count(K_FROM_J_INTEGRAL) != 0) ||
286 2154 : (_integrals.count(C_INTEGRAL) != 0 && _integrals.count(K_FROM_J_INTEGRAL) != 0))
287 0 : paramError("integrals",
288 : "JIntegral, CIntegral, and KFromJIntegral options are mutually exclusive");
289 :
290 : // Acommodate deprecated parameter convert_J_to_K
291 2154 : if (_convert_J_to_K && _integrals.count(K_FROM_J_INTEGRAL) != 0)
292 : {
293 0 : _integrals.insert(K_FROM_J_INTEGRAL);
294 0 : _integrals.erase(J_INTEGRAL);
295 : }
296 :
297 4308 : if (isParamValid("temperature"))
298 4308 : _temp = getParam<VariableName>("temperature");
299 :
300 2916 : if (_temp != "" && !isParamValid("eigenstrain_names"))
301 0 : paramError(
302 : "eigenstrain_names",
303 : "DomainIntegral error: must provide `eigenstrain_names` when temperature is coupled.");
304 :
305 2154 : if (_get_equivalent_k && (_integrals.count(INTERACTION_INTEGRAL_KI) == 0 ||
306 390 : _integrals.count(INTERACTION_INTEGRAL_KII) == 0 ||
307 390 : _integrals.count(INTERACTION_INTEGRAL_KIII) == 0))
308 0 : paramError("integrals",
309 : "DomainIntegral error: must calculate KI, KII and KIII to get equivalent K.");
310 :
311 4308 : if (isParamValid("output_variable"))
312 : {
313 630 : _output_variables = getParam<std::vector<VariableName>>("output_variable");
314 210 : if (_crack_front_points.size() > 0)
315 0 : paramError("output_variables",
316 : "'output_variables' not yet supported with 'crack_front_points'");
317 : }
318 :
319 2154 : if (_integrals.count(K_FROM_J_INTEGRAL) != 0)
320 : {
321 168 : if (!isParamValid("youngs_modulus") || !isParamValid("poissons_ratio"))
322 0 : mooseError("DomainIntegral error: must set Young's modulus and Poisson's ratio "
323 : "if K_FROM_J_INTEGRAL is selected.");
324 42 : if (!youngs_modulus_set)
325 0 : _youngs_modulus = getParam<Real>("youngs_modulus");
326 42 : if (!poissons_ratio_set)
327 0 : _poissons_ratio = getParam<Real>("poissons_ratio");
328 : }
329 :
330 2154 : if (_integrals.count(J_INTEGRAL) != 0 || _integrals.count(C_INTEGRAL) != 0 ||
331 602 : _integrals.count(K_FROM_J_INTEGRAL) != 0)
332 : {
333 3188 : if (isParamValid("eigenstrain_gradient"))
334 2 : paramError("eigenstrain_gradient",
335 : "'eigenstrain_gradient' cannot be specified when the computed integrals include "
336 : "JIntegral, CIntegral, or KFromJIntegral");
337 3184 : if (isParamValid("body_force"))
338 2 : paramError("body_force",
339 : "'body_force' cannot be specified when the computed integrals include JIntegral, "
340 : "CIntegral, or KFromJIntegral");
341 : }
342 4384 : if (isParamValid("eigenstrain_gradient") && (_temp != "" || isParamValid("eigenstrain_names")))
343 2 : paramError("eigenstrain_gradient",
344 : "'eigenstrain_gradient' cannot be specified together with 'temperature' or "
345 : "'eigenstrain_names'. These are for separate, mutually exclusive systems for "
346 : "including the effect of eigenstrains");
347 2148 : }
348 :
349 10680 : DomainIntegralAction::~DomainIntegralAction() {}
350 :
351 : void
352 2142 : DomainIntegralAction::act()
353 : {
354 2142 : const std::string uo_name("crackFrontDefinition");
355 2142 : const std::string ak_base_name("q");
356 2142 : const std::string av_base_name("q");
357 2142 : const unsigned int num_crack_front_points = calcNumCrackFrontPoints();
358 2142 : const std::string aux_stress_base_name("aux_stress");
359 2142 : const std::string aux_grad_disp_base_name("aux_grad_disp");
360 :
361 : // checking if built with xfem and setting flags for vpps used by xfem
362 : std::vector<std::string> xfem_exec_flags;
363 2142 : for (const auto & item :
364 59976 : moose::internal::ExecFlagRegistry::getExecFlagRegistry().getFlags().items())
365 : {
366 55692 : if (item.name() == "XFEM_MARK")
367 : {
368 : // TODO the xfem_flags should be something like this:
369 : // xfem_exec_flags = {item, EXEC_NONLINEAR, EXEC_TIMESTEP_END};
370 : // but the item=XFEM_MARK flag causes xfem tests to diverge
371 0 : xfem_exec_flags = {EXEC_INITIAL, EXEC_NONLINEAR, EXEC_TIMESTEP_END};
372 : }
373 : }
374 :
375 2142 : std::string ad_prepend = "";
376 2142 : if (_use_ad)
377 : ad_prepend = "AD";
378 :
379 2142 : if (_current_task == "add_user_object")
380 : {
381 358 : const std::string uo_type_name("CrackFrontDefinition");
382 :
383 358 : InputParameters params = _factory.getValidParams(uo_type_name);
384 358 : if (_use_crack_front_points_provider && !xfem_exec_flags.empty())
385 : {
386 : // The CrackFrontDefinition updates the vpps and MUST execute before them
387 0 : params.set<int>("execution_order_group") = -1;
388 0 : params.set<ExecFlagEnum>("execute_on") = xfem_exec_flags;
389 : }
390 : else
391 1432 : params.set<ExecFlagEnum>("execute_on") = {EXEC_INITIAL, EXEC_TIMESTEP_END};
392 :
393 358 : params.set<MooseEnum>("crack_direction_method") = _direction_method_moose_enum;
394 358 : params.set<MooseEnum>("crack_end_direction_method") = _end_direction_method_moose_enum;
395 358 : if (_have_crack_direction_vector)
396 291 : params.set<RealVectorValue>("crack_direction_vector") = _crack_direction_vector;
397 358 : if (_have_crack_direction_vector_end_1)
398 53 : params.set<RealVectorValue>("crack_direction_vector_end_1") = _crack_direction_vector_end_1;
399 358 : if (_have_crack_direction_vector_end_2)
400 53 : params.set<RealVectorValue>("crack_direction_vector_end_2") = _crack_direction_vector_end_2;
401 358 : if (_crack_mouth_boundary_names.size() != 0)
402 72 : params.set<std::vector<BoundaryName>>("crack_mouth_boundary") = _crack_mouth_boundary_names;
403 358 : if (_intersecting_boundary_names.size() != 0)
404 92 : params.set<std::vector<BoundaryName>>("intersecting_boundary") = _intersecting_boundary_names;
405 358 : params.set<bool>("2d") = _treat_as_2d;
406 358 : params.set<unsigned int>("axis_2d") = _axis_2d;
407 358 : if (_has_symmetry_plane)
408 244 : params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
409 358 : if (_boundary_names.size() != 0)
410 648 : params.set<std::vector<BoundaryName>>("boundary") = _boundary_names;
411 358 : if (_crack_front_points.size() != 0)
412 68 : params.set<std::vector<Point>>("crack_front_points") = _crack_front_points;
413 358 : if (_use_crack_front_points_provider)
414 0 : params.applyParameters(parameters(),
415 : {"crack_front_points_provider, number_points_from_provider"});
416 358 : if (_closed_loop)
417 0 : params.set<bool>("closed_loop") = _closed_loop;
418 358 : params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
419 358 : if (_integrals.count(INTERACTION_INTEGRAL_T) != 0)
420 : {
421 42 : params.set<VariableName>("disp_x") = _displacements[0];
422 42 : params.set<VariableName>("disp_y") = _displacements[1];
423 21 : if (_displacements.size() == 3)
424 28 : params.set<VariableName>("disp_z") = _displacements[2];
425 21 : params.set<bool>("t_stress") = true;
426 : }
427 :
428 : unsigned int nrings = 0;
429 358 : if (_q_function_type == TOPOLOGY)
430 : {
431 35 : params.set<bool>("q_function_rings") = true;
432 35 : params.set<unsigned int>("last_ring") = _ring_last;
433 35 : params.set<unsigned int>("first_ring") = _ring_first;
434 35 : nrings = _ring_last - _ring_first + 1;
435 : }
436 323 : else if (_q_function_type == GEOMETRY)
437 : {
438 323 : params.set<std::vector<Real>>("j_integral_radius_inner") = _radius_inner;
439 646 : params.set<std::vector<Real>>("j_integral_radius_outer") = _radius_outer;
440 323 : nrings = _ring_vec.size();
441 : }
442 :
443 358 : params.set<unsigned int>("nrings") = nrings;
444 358 : params.set<MooseEnum>("q_function_type") = _q_function_type;
445 :
446 358 : _problem->addUserObject(uo_type_name, uo_name, params);
447 358 : }
448 1784 : else if (_current_task == "add_aux_variable" && _output_q)
449 : {
450 921 : for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
451 : {
452 : std::string aux_var_type;
453 719 : if (_family == "LAGRANGE")
454 : aux_var_type = "MooseVariable";
455 0 : else if (_family == "MONOMIAL")
456 : aux_var_type = "MooseVariableConstMonomial";
457 0 : else if (_family == "SCALAR")
458 : aux_var_type = "MooseVariableScalar";
459 : else
460 0 : mooseError("Unsupported finite element family in, " + name() +
461 : ". Please use LAGRANGE, MONOMIAL, or SCALAR");
462 :
463 719 : auto params = _factory.getValidParams(aux_var_type);
464 719 : params.set<MooseEnum>("order") = _order;
465 719 : params.set<MooseEnum>("family") = _family;
466 :
467 719 : if (_treat_as_2d && _use_crack_front_points_provider == false)
468 : {
469 516 : std::ostringstream av_name_stream;
470 516 : av_name_stream << av_base_name << "_" << _ring_vec[ring_index];
471 516 : _problem->addAuxVariable(aux_var_type, av_name_stream.str(), params);
472 516 : }
473 : else
474 : {
475 1280 : for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
476 : {
477 1077 : std::ostringstream av_name_stream;
478 1077 : av_name_stream << av_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
479 1077 : _problem->addAuxVariable(aux_var_type, av_name_stream.str(), params);
480 1077 : }
481 : }
482 719 : }
483 : }
484 :
485 1582 : else if (_current_task == "add_aux_kernel" && _output_q)
486 : {
487 : std::string ak_type_name;
488 : unsigned int nrings = 0;
489 202 : if (_q_function_type == GEOMETRY)
490 : {
491 : ak_type_name = "DomainIntegralQFunction";
492 188 : nrings = _ring_vec.size();
493 : }
494 14 : else if (_q_function_type == TOPOLOGY)
495 : {
496 : ak_type_name = "DomainIntegralTopologicalQFunction";
497 14 : nrings = _ring_last - _ring_first + 1;
498 : }
499 :
500 202 : InputParameters params = _factory.getValidParams(ak_type_name);
501 808 : params.set<ExecFlagEnum>("execute_on") = {EXEC_INITIAL, EXEC_TIMESTEP_END};
502 404 : params.set<UserObjectName>("crack_front_definition") = uo_name;
503 202 : params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
504 :
505 921 : for (unsigned int ring_index = 0; ring_index < nrings; ++ring_index)
506 : {
507 719 : if (_q_function_type == GEOMETRY)
508 : {
509 670 : params.set<Real>("j_integral_radius_inner") = _radius_inner[ring_index];
510 670 : params.set<Real>("j_integral_radius_outer") = _radius_outer[ring_index];
511 : }
512 49 : else if (_q_function_type == TOPOLOGY)
513 : {
514 49 : params.set<unsigned int>("ring_index") = _ring_first + ring_index;
515 : }
516 :
517 719 : if (_treat_as_2d && _use_crack_front_points_provider == false)
518 : {
519 516 : std::ostringstream ak_name_stream;
520 516 : ak_name_stream << ak_base_name << "_" << _ring_vec[ring_index];
521 516 : std::ostringstream av_name_stream;
522 516 : av_name_stream << av_base_name << "_" << _ring_vec[ring_index];
523 1032 : params.set<AuxVariableName>("variable") = av_name_stream.str();
524 516 : _problem->addAuxKernel(ak_type_name, ak_name_stream.str(), params);
525 516 : }
526 : else
527 : {
528 1280 : for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
529 : {
530 1077 : std::ostringstream ak_name_stream;
531 1077 : ak_name_stream << ak_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
532 1077 : std::ostringstream av_name_stream;
533 2154 : av_name_stream << av_base_name << "_" << cfp_index + 1 << "_" << _ring_vec[ring_index];
534 2154 : params.set<AuxVariableName>("variable") = av_name_stream.str();
535 1077 : params.set<unsigned int>("crack_front_point_index") = cfp_index;
536 1077 : _problem->addAuxKernel(ak_type_name, ak_name_stream.str(), params);
537 1077 : }
538 : }
539 : }
540 202 : }
541 :
542 1380 : else if (_current_task == "add_postprocessor")
543 : {
544 919 : for (std::set<INTEGRAL>::iterator sit = _integrals.begin(); sit != _integrals.end(); ++sit)
545 : {
546 : std::string pp_base_name;
547 563 : switch (*sit)
548 : {
549 : case J_INTEGRAL:
550 : pp_base_name = "J";
551 : break;
552 :
553 : case C_INTEGRAL:
554 : pp_base_name = "C";
555 : break;
556 :
557 : case K_FROM_J_INTEGRAL:
558 : pp_base_name = "K";
559 : break;
560 :
561 : case INTERACTION_INTEGRAL_KI:
562 : pp_base_name = "II_KI";
563 : break;
564 :
565 : case INTERACTION_INTEGRAL_KII:
566 : pp_base_name = "II_KII";
567 : break;
568 :
569 : case INTERACTION_INTEGRAL_KIII:
570 : pp_base_name = "II_KIII";
571 : break;
572 :
573 : case INTERACTION_INTEGRAL_T:
574 : pp_base_name = "II_T";
575 : break;
576 : }
577 563 : const std::string pp_type_name("VectorPostprocessorComponent");
578 1126 : InputParameters params = _factory.getValidParams(pp_type_name);
579 2544 : for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
580 : {
581 1981 : if (_treat_as_2d && _use_crack_front_points_provider == false)
582 : {
583 2804 : params.set<VectorPostprocessorName>("vectorpostprocessor") =
584 2804 : pp_base_name + "_2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
585 2804 : std::string pp_name = pp_base_name + +"_" + Moose::stringify(_ring_vec[ring_index]);
586 1402 : params.set<unsigned int>("index") = 0;
587 2804 : params.set<std::string>("vector_name") =
588 4206 : pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
589 1402 : _problem->addPostprocessor(pp_type_name, pp_name, params);
590 1402 : }
591 : else
592 : {
593 3120 : for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
594 : {
595 5082 : params.set<VectorPostprocessorName>("vectorpostprocessor") =
596 5082 : pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
597 5082 : std::string pp_name = pp_base_name + "_" + Moose::stringify(cfp_index + 1) + "_" +
598 5082 : Moose::stringify(_ring_vec[ring_index]);
599 2541 : params.set<unsigned int>("index") = cfp_index;
600 5082 : params.set<std::string>("vector_name") =
601 7623 : pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
602 2541 : _problem->addPostprocessor(pp_type_name, pp_name, params);
603 : }
604 : }
605 : }
606 : }
607 :
608 356 : if (_get_equivalent_k)
609 : {
610 65 : std::string pp_base_name("Keq");
611 65 : const std::string pp_type_name("VectorPostprocessorComponent");
612 65 : InputParameters params = _factory.getValidParams(pp_type_name);
613 285 : for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
614 : {
615 220 : if (_treat_as_2d && _use_crack_front_points_provider == false)
616 : {
617 300 : params.set<VectorPostprocessorName>("vectorpostprocessor") =
618 300 : pp_base_name + "_2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
619 300 : std::string pp_name = pp_base_name + +"_" + Moose::stringify(_ring_vec[ring_index]);
620 150 : params.set<unsigned int>("index") = 0;
621 300 : params.set<std::string>("vector_name") =
622 450 : pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
623 150 : _problem->addPostprocessor(pp_type_name, pp_name, params);
624 150 : }
625 : else
626 : {
627 280 : for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
628 : {
629 420 : params.set<VectorPostprocessorName>("vectorpostprocessor") =
630 420 : pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
631 420 : std::string pp_name = pp_base_name + "_" + Moose::stringify(cfp_index + 1) + "_" +
632 420 : Moose::stringify(_ring_vec[ring_index]);
633 210 : params.set<unsigned int>("index") = cfp_index;
634 420 : params.set<std::string>("vector_name") =
635 630 : pp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
636 210 : _problem->addPostprocessor(pp_type_name, pp_name, params);
637 : }
638 : }
639 : }
640 65 : }
641 :
642 391 : for (unsigned int i = 0; i < _output_variables.size(); ++i)
643 : {
644 35 : const std::string ov_base_name(_output_variables[i]);
645 35 : const std::string pp_type_name("CrackFrontData");
646 35 : InputParameters params = _factory.getValidParams(pp_type_name);
647 70 : params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
648 70 : params.set<UserObjectName>("crack_front_definition") = uo_name;
649 35 : if (_treat_as_2d && _use_crack_front_points_provider == false)
650 : {
651 0 : std::ostringstream pp_name_stream;
652 0 : pp_name_stream << ov_base_name << "_crack";
653 0 : params.set<VariableName>("variable") = _output_variables[i];
654 0 : _problem->addPostprocessor(pp_type_name, pp_name_stream.str(), params);
655 0 : }
656 : else
657 : {
658 140 : for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
659 : {
660 105 : std::ostringstream pp_name_stream;
661 105 : pp_name_stream << ov_base_name << "_crack_" << cfp_index + 1;
662 105 : params.set<VariableName>("variable") = _output_variables[i];
663 105 : params.set<unsigned int>("crack_front_point_index") = cfp_index;
664 105 : _problem->addPostprocessor(pp_type_name, pp_name_stream.str(), params);
665 105 : }
666 : }
667 35 : }
668 : }
669 :
670 1024 : else if (_current_task == "add_vector_postprocessor")
671 : {
672 356 : if (_integrals.count(J_INTEGRAL) != 0 || _integrals.count(C_INTEGRAL) != 0 ||
673 363 : _integrals.count(K_FROM_J_INTEGRAL) != 0)
674 : {
675 : std::string vpp_base_name;
676 263 : std::string jintegral_selection = "JIntegral";
677 :
678 263 : if (_integrals.count(J_INTEGRAL) != 0)
679 : {
680 : vpp_base_name = "J";
681 : jintegral_selection = "JIntegral";
682 : }
683 21 : else if (_integrals.count(K_FROM_J_INTEGRAL) != 0)
684 : {
685 : vpp_base_name = "K";
686 : jintegral_selection = "KFromJIntegral";
687 : }
688 14 : else if (_integrals.count(C_INTEGRAL) != 0)
689 : {
690 : vpp_base_name = "C";
691 : jintegral_selection = "CIntegral";
692 : }
693 :
694 263 : if (_treat_as_2d && _use_crack_front_points_provider == false)
695 : vpp_base_name += "_2DVPP";
696 :
697 263 : const std::string vpp_type_name("JIntegral");
698 263 : InputParameters params = _factory.getValidParams(vpp_type_name);
699 526 : if (!getParam<bool>("output_vpp"))
700 42 : params.set<std::vector<OutputName>>("outputs") = {"none"};
701 :
702 526 : params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
703 526 : params.set<UserObjectName>("crack_front_definition") = uo_name;
704 263 : params.set<std::vector<SubdomainName>>("block") = {_blocks};
705 263 : params.set<MooseEnum>("position_type") = _position_type;
706 :
707 263 : if (_integrals.count(K_FROM_J_INTEGRAL) != 0)
708 : {
709 7 : params.set<Real>("youngs_modulus") = _youngs_modulus;
710 7 : params.set<Real>("poissons_ratio") = _poissons_ratio;
711 : }
712 :
713 263 : if (_has_symmetry_plane)
714 214 : params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
715 :
716 : // Select the integral type to be computed in JIntegral vector postprocessor
717 263 : params.set<MooseEnum>("integral") = jintegral_selection;
718 :
719 263 : params.set<std::vector<VariableName>>("displacements") = _displacements;
720 263 : params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
721 1220 : for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
722 : {
723 957 : params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
724 1914 : params.set<MooseEnum>("q_function_type") = _q_function_type;
725 :
726 1914 : std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
727 957 : _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
728 : }
729 263 : }
730 :
731 356 : if (_integrals.count(INTERACTION_INTEGRAL_KI) != 0 ||
732 207 : _integrals.count(INTERACTION_INTEGRAL_KII) != 0 ||
733 207 : _integrals.count(INTERACTION_INTEGRAL_KIII) != 0 ||
734 207 : _integrals.count(INTERACTION_INTEGRAL_T) != 0)
735 : {
736 149 : if (_has_symmetry_plane && (_integrals.count(INTERACTION_INTEGRAL_KII) != 0 ||
737 233 : _integrals.count(INTERACTION_INTEGRAL_KIII) != 0))
738 0 : paramError("symmetry_plane",
739 : "In DomainIntegral, symmetry_plane option cannot be used with mode-II or "
740 : "mode-III interaction integral");
741 :
742 : std::string vpp_base_name;
743 149 : std::string vpp_type_name(ad_prepend + "InteractionIntegral");
744 :
745 149 : InputParameters params = _factory.getValidParams(vpp_type_name);
746 298 : if (!getParam<bool>("output_vpp"))
747 0 : params.set<std::vector<OutputName>>("outputs") = {"none"};
748 :
749 149 : if (_use_crack_front_points_provider && !xfem_exec_flags.empty())
750 0 : params.set<ExecFlagEnum>("execute_on") = xfem_exec_flags;
751 : else
752 447 : params.set<ExecFlagEnum>("execute_on") = {EXEC_TIMESTEP_END};
753 :
754 319 : if (isParamValid("additional_eigenstrain_00") && isParamValid("additional_eigenstrain_01") &&
755 184 : isParamValid("additional_eigenstrain_11") && isParamValid("additional_eigenstrain_22"))
756 : {
757 14 : params.set<CoupledName>("additional_eigenstrain_00") =
758 14 : getParam<CoupledName>("additional_eigenstrain_00");
759 14 : params.set<CoupledName>("additional_eigenstrain_01") =
760 14 : getParam<CoupledName>("additional_eigenstrain_01");
761 14 : params.set<CoupledName>("additional_eigenstrain_11") =
762 14 : getParam<CoupledName>("additional_eigenstrain_11");
763 14 : params.set<CoupledName>("additional_eigenstrain_22") =
764 21 : getParam<CoupledName>("additional_eigenstrain_22");
765 : }
766 :
767 298 : params.set<UserObjectName>("crack_front_definition") = uo_name;
768 149 : params.set<bool>("use_displaced_mesh") = _use_displaced_mesh;
769 149 : params.set<std::vector<SubdomainName>>("block") = {_blocks};
770 :
771 149 : if (_has_symmetry_plane)
772 84 : params.set<unsigned int>("symmetry_plane") = _symmetry_plane;
773 :
774 149 : if (_fgm_crack)
775 : {
776 14 : params.set<MaterialPropertyName>(
777 : "functionally_graded_youngs_modulus_crack_dir_gradient") = {
778 : _functionally_graded_youngs_modulus_crack_dir_gradient};
779 28 : params.set<MaterialPropertyName>("functionally_graded_youngs_modulus") = {
780 : _functionally_graded_youngs_modulus};
781 : }
782 :
783 149 : params.set<Real>("poissons_ratio") = _poissons_ratio;
784 149 : params.set<Real>("youngs_modulus") = _youngs_modulus;
785 149 : params.set<std::vector<VariableName>>("displacements") = _displacements;
786 149 : if (_temp != "")
787 42 : params.set<std::vector<VariableName>>("temperature") = {_temp};
788 :
789 298 : if (parameters().isParamValid("eigenstrain_gradient"))
790 7 : params.set<MaterialPropertyName>("eigenstrain_gradient") =
791 14 : parameters().get<MaterialPropertyName>("eigenstrain_gradient");
792 298 : if (parameters().isParamValid("body_force"))
793 7 : params.set<MaterialPropertyName>("body_force") =
794 14 : parameters().get<MaterialPropertyName>("body_force");
795 :
796 505 : for (std::set<INTEGRAL>::iterator sit = _integrals.begin(); sit != _integrals.end(); ++sit)
797 : {
798 356 : switch (*sit)
799 : {
800 42 : case J_INTEGRAL:
801 42 : continue;
802 :
803 7 : case C_INTEGRAL:
804 7 : continue;
805 :
806 7 : case K_FROM_J_INTEGRAL:
807 7 : continue;
808 :
809 : case INTERACTION_INTEGRAL_KI:
810 : vpp_base_name = "II_KI";
811 149 : params.set<Real>("K_factor") =
812 149 : 0.5 * _youngs_modulus / (1.0 - std::pow(_poissons_ratio, 2.0));
813 298 : params.set<MooseEnum>("sif_mode") = "KI";
814 149 : break;
815 :
816 : case INTERACTION_INTEGRAL_KII:
817 : vpp_base_name = "II_KII";
818 65 : params.set<Real>("K_factor") =
819 65 : 0.5 * _youngs_modulus / (1.0 - std::pow(_poissons_ratio, 2.0));
820 130 : params.set<MooseEnum>("sif_mode") = "KII";
821 65 : break;
822 :
823 : case INTERACTION_INTEGRAL_KIII:
824 : vpp_base_name = "II_KIII";
825 65 : params.set<Real>("K_factor") = 0.5 * _youngs_modulus / (1.0 + _poissons_ratio);
826 130 : params.set<MooseEnum>("sif_mode") = "KIII";
827 65 : break;
828 :
829 : case INTERACTION_INTEGRAL_T:
830 : vpp_base_name = "II_T";
831 21 : params.set<Real>("K_factor") = _youngs_modulus / (1 - std::pow(_poissons_ratio, 2));
832 42 : params.set<MooseEnum>("sif_mode") = "T";
833 21 : break;
834 : }
835 300 : if (_treat_as_2d && _use_crack_front_points_provider == false)
836 : vpp_base_name += "_2DVPP";
837 1324 : for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
838 : {
839 1024 : params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
840 2048 : params.set<MooseEnum>("q_function_type") = _q_function_type;
841 :
842 2048 : std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
843 1024 : _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
844 : }
845 : }
846 149 : }
847 :
848 356 : if (_get_equivalent_k)
849 : {
850 65 : std::string vpp_base_name("Keq");
851 65 : if (_treat_as_2d && _use_crack_front_points_provider == false)
852 : vpp_base_name += "_2DVPP";
853 65 : const std::string vpp_type_name("MixedModeEquivalentK");
854 65 : InputParameters params = _factory.getValidParams(vpp_type_name);
855 65 : params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
856 65 : params.set<Real>("poissons_ratio") = _poissons_ratio;
857 :
858 285 : for (unsigned int ring_index = 0; ring_index < _ring_vec.size(); ++ring_index)
859 : {
860 220 : std::string ki_name = "II_KI_";
861 220 : std::string kii_name = "II_KII_";
862 220 : std::string kiii_name = "II_KIII_";
863 220 : params.set<unsigned int>("ring_index") = _ring_vec[ring_index];
864 220 : if (_treat_as_2d && _use_crack_front_points_provider == false)
865 : {
866 300 : params.set<VectorPostprocessorName>("KI_vectorpostprocessor") =
867 300 : ki_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
868 300 : params.set<VectorPostprocessorName>("KII_vectorpostprocessor") =
869 300 : kii_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
870 300 : params.set<VectorPostprocessorName>("KIII_vectorpostprocessor") =
871 300 : kiii_name + "2DVPP_" + Moose::stringify(_ring_vec[ring_index]);
872 : }
873 : else
874 : {
875 140 : params.set<VectorPostprocessorName>("KI_vectorpostprocessor") =
876 140 : ki_name + Moose::stringify(_ring_vec[ring_index]);
877 140 : params.set<VectorPostprocessorName>("KII_vectorpostprocessor") =
878 140 : kii_name + Moose::stringify(_ring_vec[ring_index]);
879 140 : params.set<VectorPostprocessorName>("KIII_vectorpostprocessor") =
880 140 : kiii_name + Moose::stringify(_ring_vec[ring_index]);
881 : }
882 440 : params.set<std::string>("KI_vector_name") =
883 660 : ki_name + Moose::stringify(_ring_vec[ring_index]);
884 440 : params.set<std::string>("KII_vector_name") =
885 660 : kii_name + Moose::stringify(_ring_vec[ring_index]);
886 440 : params.set<std::string>("KIII_vector_name") =
887 660 : kiii_name + Moose::stringify(_ring_vec[ring_index]);
888 440 : std::string vpp_name = vpp_base_name + "_" + Moose::stringify(_ring_vec[ring_index]);
889 220 : _problem->addVectorPostprocessor(vpp_type_name, vpp_name, params);
890 : }
891 65 : }
892 :
893 356 : if (!_treat_as_2d || _use_crack_front_points_provider == true)
894 : {
895 172 : for (unsigned int i = 0; i < _output_variables.size(); ++i)
896 : {
897 35 : const std::string vpp_type_name("VectorOfPostprocessors");
898 35 : InputParameters params = _factory.getValidParams(vpp_type_name);
899 35 : params.set<ExecFlagEnum>("execute_on") = EXEC_TIMESTEP_END;
900 35 : std::ostringstream vpp_name_stream;
901 35 : vpp_name_stream << _output_variables[i] << "_crack";
902 : std::vector<PostprocessorName> postprocessor_names;
903 140 : for (unsigned int cfp_index = 0; cfp_index < num_crack_front_points; ++cfp_index)
904 : {
905 105 : std::ostringstream pp_name_stream;
906 210 : pp_name_stream << vpp_name_stream.str() << "_" << cfp_index + 1;
907 105 : postprocessor_names.push_back(pp_name_stream.str());
908 105 : }
909 35 : params.set<std::vector<PostprocessorName>>("postprocessors") = postprocessor_names;
910 35 : _problem->addVectorPostprocessor(vpp_type_name, vpp_name_stream.str(), params);
911 35 : }
912 : }
913 : }
914 :
915 668 : else if (_current_task == "add_material")
916 : {
917 358 : if (_temp != "")
918 : {
919 : std::string mater_name;
920 42 : const std::string mater_type_name("ThermalFractureIntegral");
921 : mater_name = "ThermalFractureIntegral";
922 :
923 42 : InputParameters params = _factory.getValidParams(mater_type_name);
924 84 : params.set<std::vector<MaterialPropertyName>>("eigenstrain_names") =
925 126 : getParam<std::vector<MaterialPropertyName>>("eigenstrain_names");
926 126 : params.set<std::vector<VariableName>>("temperature") = {_temp};
927 42 : params.set<std::vector<SubdomainName>>("block") = {_blocks};
928 42 : _problem->addMaterial(mater_type_name, mater_name, params);
929 42 : }
930 1074 : MultiMooseEnum integral_moose_enums = getParam<MultiMooseEnum>("integrals");
931 : bool have_j_integral = false;
932 : bool have_c_integral = false;
933 :
934 923 : for (auto ime : integral_moose_enums)
935 : {
936 1186 : if (ime == "JIntegral" || ime == "CIntegral" || ime == "KFromJIntegral" ||
937 537 : ime == "InteractionIntegralKI" || ime == "InteractionIntegralKII" ||
938 672 : ime == "InteractionIntegralKIII" || ime == "InteractionIntegralT")
939 : have_j_integral = true;
940 :
941 565 : if (ime == "CIntegral")
942 : have_c_integral = true;
943 : }
944 358 : if (have_j_integral)
945 : {
946 : std::string mater_name;
947 358 : const std::string mater_type_name(ad_prepend + "StrainEnergyDensity");
948 358 : mater_name = ad_prepend + "StrainEnergyDensity";
949 :
950 358 : InputParameters params = _factory.getValidParams(mater_type_name);
951 716 : _incremental = getParam<bool>("incremental");
952 358 : params.set<bool>("incremental") = _incremental;
953 358 : params.set<std::vector<SubdomainName>>("block") = {_blocks};
954 358 : _problem->addMaterial(mater_type_name, mater_name, params);
955 :
956 : {
957 : std::string mater_name;
958 358 : const std::string mater_type_name(ad_prepend + "EshelbyTensor");
959 358 : mater_name = ad_prepend + "EshelbyTensor";
960 :
961 358 : InputParameters params = _factory.getValidParams(mater_type_name);
962 716 : _displacements = getParam<std::vector<VariableName>>("displacements");
963 358 : params.set<std::vector<VariableName>>("displacements") = _displacements;
964 358 : params.set<std::vector<SubdomainName>>("block") = {_blocks};
965 :
966 358 : if (have_c_integral)
967 14 : params.set<bool>("compute_dissipation") = true;
968 :
969 358 : if (_temp != "")
970 126 : params.set<std::vector<VariableName>>("temperature") = {_temp};
971 :
972 358 : _problem->addMaterial(mater_type_name, mater_name, params);
973 356 : }
974 : // Strain energy rate density needed for C(t)/C* integral
975 356 : if (have_c_integral)
976 : {
977 : std::string mater_name;
978 14 : const std::string mater_type_name(ad_prepend + "StrainEnergyRateDensity");
979 14 : mater_name = ad_prepend + "StrainEnergyRateDensity";
980 :
981 14 : InputParameters params = _factory.getValidParams(mater_type_name);
982 14 : params.set<std::vector<SubdomainName>>("block") = {_blocks};
983 28 : params.set<std::vector<MaterialName>>("inelastic_models") =
984 28 : getParam<std::vector<MaterialName>>("inelastic_models");
985 :
986 14 : _problem->addMaterial(mater_type_name, mater_name, params);
987 14 : }
988 356 : }
989 356 : }
990 4989 : }
991 :
992 : unsigned int
993 2142 : DomainIntegralAction::calcNumCrackFrontPoints()
994 : {
995 : unsigned int num_points = 0;
996 2142 : if (_boundary_names.size() != 0)
997 : {
998 1938 : std::vector<BoundaryID> bids = _mesh->getBoundaryIDs(_boundary_names, true);
999 : std::set<unsigned int> nodes;
1000 :
1001 1938 : ConstBndNodeRange & bnd_nodes = *_mesh->getBoundaryNodeRange();
1002 207978 : for (ConstBndNodeRange::const_iterator nd = bnd_nodes.begin(); nd != bnd_nodes.end(); ++nd)
1003 : {
1004 206040 : const BndNode * bnode = *nd;
1005 206040 : BoundaryID boundary_id = bnode->_bnd_id;
1006 :
1007 407982 : for (unsigned int ibid = 0; ibid < bids.size(); ++ibid)
1008 : {
1009 206040 : if (boundary_id == bids[ibid])
1010 : {
1011 4098 : nodes.insert(bnode->_node->id());
1012 4098 : break;
1013 : }
1014 : }
1015 : }
1016 1938 : num_points = nodes.size();
1017 1938 : }
1018 204 : else if (_crack_front_points.size() != 0)
1019 204 : num_points = _crack_front_points.size();
1020 0 : else if (_use_crack_front_points_provider)
1021 0 : num_points = getParam<unsigned int>("number_points_from_provider");
1022 : else
1023 0 : mooseError("Must define either 'boundary' or 'crack_front_points'");
1024 2142 : return num_points;
1025 : }
1026 :
1027 : void
1028 6420 : DomainIntegralAction::addRelationshipManagers(Moose::RelationshipManagerType input_rm_type)
1029 : {
1030 6420 : if (_integrals.count(INTERACTION_INTEGRAL_T) != 0)
1031 : {
1032 378 : InputParameters params = _factory.getValidParams("CrackFrontDefinition");
1033 378 : addRelationshipManagers(input_rm_type, params);
1034 378 : }
1035 6420 : }
|