Line data Source code
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 :
10 : #include "ADComputeSmearedCrackingStress.h"
11 : #include "ElasticityTensorTools.h"
12 : #include "StressUpdateBase.h"
13 : #include "Conversion.h"
14 :
15 : registerMooseObject("TensorMechanicsApp", ADComputeSmearedCrackingStress);
16 :
17 : InputParameters
18 157 : ADComputeSmearedCrackingStress::validParams()
19 : {
20 157 : InputParameters params = ADComputeMultipleInelasticStress::validParams();
21 157 : params.addClassDescription(
22 : "Compute stress using a fixed smeared cracking model. Uses automatic differentiation");
23 314 : params.addRequiredParam<std::vector<MaterialName>>(
24 : "softening_models",
25 : "The material objects used to compute softening behavior for loading a crack."
26 : "Either 1 or 3 models must be specified. If a single model is specified, it is"
27 : "used for all directions. If 3 models are specified, they will be used for the"
28 : "3 crack directions in sequence");
29 314 : params.addRequiredCoupledVar(
30 : "cracking_stress",
31 : "The stress threshold beyond which cracking occurs. Negative values prevent cracking.");
32 314 : MultiMooseEnum direction("x y z");
33 314 : params.addParam<MultiMooseEnum>(
34 : "prescribed_crack_directions", direction, "Prescribed directions of first cracks");
35 314 : params.addParam<unsigned int>(
36 314 : "max_cracks", 3, "The maximum number of cracks allowed at a material point.");
37 314 : params.addRangeCheckedParam<Real>("cracking_neg_fraction",
38 : 0,
39 : "cracking_neg_fraction <= 1 & cracking_neg_fraction >= 0",
40 : "The fraction of the cracking strain at which "
41 : "a transition begins during decreasing "
42 : "strain to the original stiffness.");
43 314 : params.addParam<Real>(
44 : "max_stress_correction",
45 314 : 1.0,
46 : "Maximum permitted correction to the predicted stress as a ratio of the "
47 : "stress change to the predicted stress from the previous step's damage level. "
48 : "Values less than 1 will improve robustness, but not be as accurate.");
49 :
50 471 : params.addRangeCheckedParam<Real>(
51 : "shear_retention_factor",
52 314 : 0.0,
53 : "shear_retention_factor>=0 & shear_retention_factor<=1.0",
54 : "Fraction of original shear stiffness to be retained after cracking");
55 157 : params.set<std::vector<MaterialName>>("inelastic_models") = {};
56 :
57 314 : MooseEnum crackedElasticityType("DIAGONAL FULL", "DIAGONAL");
58 314 : crackedElasticityType.addDocumentation(
59 : "DIAGONAL", "Zero out terms coupling with directions orthogonal to a crack (legacy)");
60 314 : crackedElasticityType.addDocumentation(
61 : "FULL", "Consistently scale all entries based on damage (recommended)");
62 314 : params.addParam<MooseEnum>(
63 : "cracked_elasticity_type",
64 : crackedElasticityType,
65 : "Method to modify the local elasticity tensor to account for cracking");
66 :
67 157 : return params;
68 157 : }
69 :
70 118 : ADComputeSmearedCrackingStress::ADComputeSmearedCrackingStress(const InputParameters & parameters)
71 : : ADComputeMultipleInelasticStress(parameters),
72 118 : _cracking_stress(adCoupledValue("cracking_stress")),
73 236 : _max_cracks(getParam<unsigned int>("max_cracks")),
74 236 : _cracking_neg_fraction(getParam<Real>("cracking_neg_fraction")),
75 236 : _shear_retention_factor(getParam<Real>("shear_retention_factor")),
76 236 : _max_stress_correction(getParam<Real>("max_stress_correction")),
77 118 : _cracked_elasticity_type(
78 118 : getParam<MooseEnum>("cracked_elasticity_type").getEnum<CrackedElasticityType>()),
79 118 : _crack_damage(declareADProperty<RealVectorValue>(_base_name + "crack_damage")),
80 236 : _crack_damage_old(getMaterialPropertyOld<RealVectorValue>(_base_name + "crack_damage")),
81 118 : _crack_flags(declareADProperty<RealVectorValue>(_base_name + "crack_flags")),
82 118 : _crack_rotation(declareADProperty<RankTwoTensor>(_base_name + "crack_rotation")),
83 236 : _crack_rotation_old(getMaterialPropertyOld<RankTwoTensor>(_base_name + "crack_rotation")),
84 118 : _crack_initiation_strain(
85 118 : declareADProperty<RealVectorValue>(_base_name + "crack_initiation_strain")),
86 118 : _crack_initiation_strain_old(
87 118 : getMaterialPropertyOld<RealVectorValue>(_base_name + "crack_initiation_strain")),
88 118 : _crack_max_strain(declareADProperty<RealVectorValue>(_base_name + "crack_max_strain")),
89 472 : _crack_max_strain_old(getMaterialPropertyOld<RealVectorValue>(_base_name + "crack_max_strain"))
90 : {
91 : MultiMooseEnum prescribed_crack_directions =
92 236 : getParam<MultiMooseEnum>("prescribed_crack_directions");
93 118 : if (prescribed_crack_directions.size() > 0)
94 : {
95 36 : if (prescribed_crack_directions.size() > 3)
96 0 : mooseError("A maximum of three crack directions may be specified");
97 114 : for (unsigned int i = 0; i < prescribed_crack_directions.size(); ++i)
98 : {
99 138 : for (unsigned int j = 0; j < i; ++j)
100 60 : if (prescribed_crack_directions[i] == prescribed_crack_directions[j])
101 0 : mooseError("Entries in 'prescribed_crack_directions' cannot be repeated");
102 : _prescribed_crack_directions.push_back(
103 78 : static_cast<unsigned int>(prescribed_crack_directions.get(i)));
104 : }
105 :
106 : // Fill in the last remaining direction if 2 are specified
107 36 : if (_prescribed_crack_directions.size() == 2)
108 : {
109 6 : std::set<unsigned int> available_dirs = {0, 1, 2};
110 18 : for (auto dir : _prescribed_crack_directions)
111 12 : if (available_dirs.erase(dir) != 1)
112 0 : mooseError("Invalid prescribed crack direction:" + Moose::stringify(dir));
113 6 : if (available_dirs.size() != 1)
114 0 : mooseError("Error in finding remaining available crack direction");
115 6 : _prescribed_crack_directions.push_back(*available_dirs.begin());
116 : }
117 : }
118 236 : if (!isParamSetByUser("cracked_elasticity_type"))
119 25 : paramWarning(
120 : "cracked_elasticity_type",
121 : "Defaulting to the legacy option of 'DIAGONAL', but the 'FULL' option is preferred");
122 :
123 117 : _local_elastic_vector.resize(9);
124 117 : }
125 :
126 : void
127 448 : ADComputeSmearedCrackingStress::initQpStatefulProperties()
128 : {
129 448 : ADComputeMultipleInelasticStress::initQpStatefulProperties();
130 :
131 448 : _crack_damage[_qp] = 0.0;
132 :
133 448 : _crack_initiation_strain[_qp] = 0.0;
134 448 : _crack_max_strain[_qp](0) = 0.0;
135 :
136 448 : switch (_prescribed_crack_directions.size())
137 : {
138 288 : case 0:
139 : {
140 288 : _crack_rotation[_qp] = ADRankTwoTensor::Identity();
141 288 : break;
142 : }
143 64 : case 1:
144 : {
145 64 : _crack_rotation[_qp].zero();
146 64 : switch (_prescribed_crack_directions[0])
147 : {
148 32 : case 0:
149 : {
150 32 : _crack_rotation[_qp](0, 0) = 1.0;
151 32 : _crack_rotation[_qp](1, 1) = 1.0;
152 32 : _crack_rotation[_qp](2, 2) = 1.0;
153 32 : break;
154 : }
155 0 : case 1:
156 : {
157 0 : _crack_rotation[_qp](1, 0) = 1.0;
158 0 : _crack_rotation[_qp](0, 1) = 1.0;
159 0 : _crack_rotation[_qp](2, 2) = 1.0;
160 0 : break;
161 : }
162 32 : case 2:
163 : {
164 32 : _crack_rotation[_qp](2, 0) = 1.0;
165 32 : _crack_rotation[_qp](0, 1) = 1.0;
166 32 : _crack_rotation[_qp](1, 2) = 1.0;
167 32 : break;
168 : }
169 : }
170 : break;
171 : }
172 0 : case 2:
173 : {
174 0 : mooseError("Number of prescribed crack directions cannot be 2");
175 : break;
176 : }
177 : case 3:
178 : {
179 384 : for (unsigned int i = 0; i < _prescribed_crack_directions.size(); ++i)
180 : {
181 : ADRealVectorValue crack_dir_vec;
182 288 : crack_dir_vec(_prescribed_crack_directions[i]) = 1.0;
183 288 : _crack_rotation[_qp].fillColumn(i, crack_dir_vec);
184 : }
185 : }
186 : }
187 448 : }
188 :
189 : void
190 115 : ADComputeSmearedCrackingStress::initialSetup()
191 : {
192 115 : ADComputeMultipleInelasticStress::initialSetup();
193 :
194 230 : if (!hasGuaranteedMaterialProperty(_elasticity_tensor_name, Guarantee::ISOTROPIC))
195 0 : mooseError("ADComputeSmearedCrackingStress requires that the elasticity tensor be "
196 : "guaranteed isotropic");
197 :
198 345 : std::vector<MaterialName> soft_matls = getParam<std::vector<MaterialName>>("softening_models");
199 267 : for (auto soft_matl : soft_matls)
200 : {
201 : ADSmearedCrackSofteningBase * scsb =
202 152 : dynamic_cast<ADSmearedCrackSofteningBase *>(&getMaterialByName(soft_matl));
203 152 : if (scsb)
204 152 : _softening_models.push_back(scsb);
205 : else
206 0 : paramError("softening_models", "Model " + soft_matl + " is not a softening model");
207 : }
208 115 : if (_softening_models.size() == 1)
209 : {
210 : // Reuse the same model in all 3 directions
211 96 : _softening_models.push_back(_softening_models[0]);
212 96 : _softening_models.push_back(_softening_models[0]);
213 : }
214 19 : else if (_softening_models.size() != 3)
215 1 : paramError("softening_models", "Either 1 or 3 softening models must be specified");
216 114 : }
217 :
218 : void
219 415472 : ADComputeSmearedCrackingStress::computeQpStress()
220 : {
221 :
222 415472 : if (!previouslyCracked())
223 26688 : computeQpStressIntermediateConfiguration();
224 : else
225 : {
226 388784 : _elastic_strain[_qp] = _elastic_strain_old[_qp] + _strain_increment[_qp];
227 :
228 : // Propagate behavior from the (now inactive) inelastic models
229 388784 : _inelastic_strain[_qp] = _inelastic_strain_old[_qp];
230 388784 : for (auto model : _models)
231 : {
232 0 : model->setQp(_qp);
233 0 : model->propagateQpStatefulProperties();
234 : }
235 :
236 : // Since the other inelastic models are inactive, they will not limit the time step
237 388784 : _material_timestep_limit[_qp] = std::numeric_limits<Real>::max();
238 :
239 : // update _local_elasticity_tensor based on cracking state in previous time step
240 388784 : updateLocalElasticityTensor();
241 :
242 : // Calculate stress in intermediate configuration
243 388784 : _stress[_qp] = _local_elasticity_tensor * _elastic_strain[_qp];
244 : }
245 :
246 : // compute crack status and adjust stress
247 415472 : updateCrackingStateAndStress();
248 :
249 415472 : if (_perform_finite_strain_rotations)
250 : {
251 415472 : finiteStrainRotation();
252 415472 : _crack_rotation[_qp] = _rotation_increment[_qp] * _crack_rotation[_qp];
253 : }
254 415472 : }
255 :
256 : void
257 388784 : ADComputeSmearedCrackingStress::updateLocalElasticityTensor()
258 : {
259 : const ADReal youngs_modulus =
260 388784 : ElasticityTensorTools::getIsotropicYoungsModulus(_elasticity_tensor[_qp]);
261 :
262 : bool cracking_locally_active = false;
263 :
264 388784 : const ADReal cracking_stress = _cracking_stress[_qp];
265 :
266 388784 : if (cracking_stress > 0)
267 : {
268 388784 : ADRealVectorValue stiffness_ratio_local(1.0, 1.0, 1.0);
269 388784 : const ADRankTwoTensor & R = _crack_rotation_old[_qp];
270 388784 : ADRankTwoTensor ePrime(_elastic_strain_old[_qp]);
271 388784 : ePrime.rotate(R.transpose());
272 :
273 1555136 : for (unsigned int i = 0; i < 3; ++i)
274 : {
275 : // Update elasticity tensor based on crack status of the end of last time step
276 1166352 : if (_crack_damage_old[_qp](i) > 0.0)
277 : {
278 1118720 : if (_cracking_neg_fraction == 0.0 && MooseUtils::absoluteFuzzyLessThan(ePrime(i, i), 0.0))
279 34072 : stiffness_ratio_local(i) = 1.0;
280 0 : else if (_cracking_neg_fraction > 0.0 &&
281 1050576 : ePrime(i, i) < _crack_initiation_strain_old[_qp](i) * _cracking_neg_fraction &&
282 525288 : ePrime(i, i) > -_crack_initiation_strain_old[_qp](i) * _cracking_neg_fraction)
283 : {
284 0 : const ADReal etr = _cracking_neg_fraction * _crack_initiation_strain_old[_qp](i);
285 : const ADReal Eo = cracking_stress / _crack_initiation_strain_old[_qp](i);
286 0 : const ADReal Ec = Eo * (1.0 - _crack_damage_old[_qp](i));
287 0 : const ADReal a = (Ec - Eo) / (4 * etr);
288 0 : const ADReal b = (Ec + Eo) / 2;
289 : // Compute the ratio of the current transition stiffness to the original stiffness
290 0 : stiffness_ratio_local(i) = (2.0 * a * etr + b) / Eo;
291 : cracking_locally_active = true;
292 : }
293 : else
294 : {
295 525288 : stiffness_ratio_local(i) = (1.0 - _crack_damage_old[_qp](i));
296 : cracking_locally_active = true;
297 : }
298 : }
299 : }
300 :
301 388784 : if (cracking_locally_active)
302 : {
303 : // Update the elasticity tensor in the crack coordinate system
304 355192 : if (_cracked_elasticity_type == CrackedElasticityType::DIAGONAL)
305 : {
306 1168 : const bool c0_coupled = MooseUtils::absoluteFuzzyEqual(stiffness_ratio_local(0), 1.0);
307 1168 : const bool c1_coupled = MooseUtils::absoluteFuzzyEqual(stiffness_ratio_local(1), 1.0);
308 1168 : const bool c2_coupled = MooseUtils::absoluteFuzzyEqual(stiffness_ratio_local(2), 1.0);
309 :
310 2336 : const ADReal c01 = (c0_coupled && c1_coupled ? 1.0 : 0.0);
311 1856 : const ADReal c02 = (c0_coupled && c2_coupled ? 1.0 : 0.0);
312 1648 : const ADReal c12 = (c1_coupled && c2_coupled ? 1.0 : 0.0);
313 :
314 : const ADReal c01_shear_retention =
315 1168 : (c0_coupled && c1_coupled ? 1.0 : _shear_retention_factor);
316 : const ADReal c02_shear_retention =
317 1168 : (c0_coupled && c2_coupled ? 1.0 : _shear_retention_factor);
318 : const ADReal c12_shear_retention =
319 1168 : (c1_coupled && c2_coupled ? 1.0 : _shear_retention_factor);
320 :
321 1168 : _local_elastic_vector[0] = (c0_coupled ? _elasticity_tensor[_qp](0, 0, 0, 0)
322 1168 : : stiffness_ratio_local(0) * youngs_modulus);
323 2336 : _local_elastic_vector[1] = _elasticity_tensor[_qp](0, 0, 1, 1) * c01;
324 2336 : _local_elastic_vector[2] = _elasticity_tensor[_qp](0, 0, 2, 2) * c02;
325 1168 : _local_elastic_vector[3] = (c1_coupled ? _elasticity_tensor[_qp](1, 1, 1, 1)
326 1168 : : stiffness_ratio_local(1) * youngs_modulus);
327 2336 : _local_elastic_vector[4] = _elasticity_tensor[_qp](1, 1, 2, 2) * c12;
328 1168 : _local_elastic_vector[5] = (c2_coupled ? _elasticity_tensor[_qp](2, 2, 2, 2)
329 1168 : : stiffness_ratio_local(2) * youngs_modulus);
330 2336 : _local_elastic_vector[6] = _elasticity_tensor[_qp](1, 2, 1, 2) * c12_shear_retention;
331 2336 : _local_elastic_vector[7] = _elasticity_tensor[_qp](0, 2, 0, 2) * c02_shear_retention;
332 2336 : _local_elastic_vector[8] = _elasticity_tensor[_qp](0, 1, 0, 1) * c01_shear_retention;
333 : }
334 : else // _cracked_elasticity_type == CrackedElasticityType::FULL
335 : {
336 : const ADReal & c0 = stiffness_ratio_local(0);
337 : const ADReal & c1 = stiffness_ratio_local(1);
338 : const ADReal & c2 = stiffness_ratio_local(2);
339 :
340 : const ADReal c01 = c0 * c1;
341 : const ADReal c02 = c0 * c2;
342 : const ADReal c12 = c1 * c2;
343 :
344 354024 : const ADReal c01_shear_retention = std::max(c01, _shear_retention_factor);
345 354024 : const ADReal c02_shear_retention = std::max(c02, _shear_retention_factor);
346 354024 : const ADReal c12_shear_retention = std::max(c12, _shear_retention_factor);
347 :
348 708048 : _local_elastic_vector[0] = _elasticity_tensor[_qp](0, 0, 0, 0) * c0;
349 708048 : _local_elastic_vector[1] = _elasticity_tensor[_qp](0, 0, 1, 1) * c01;
350 708048 : _local_elastic_vector[2] = _elasticity_tensor[_qp](0, 0, 2, 2) * c02;
351 708048 : _local_elastic_vector[3] = _elasticity_tensor[_qp](1, 1, 1, 1) * c1;
352 708048 : _local_elastic_vector[4] = _elasticity_tensor[_qp](1, 1, 2, 2) * c12;
353 708048 : _local_elastic_vector[5] = _elasticity_tensor[_qp](2, 2, 2, 2) * c2;
354 708048 : _local_elastic_vector[6] = _elasticity_tensor[_qp](1, 2, 1, 2) * c12_shear_retention;
355 708048 : _local_elastic_vector[7] = _elasticity_tensor[_qp](0, 2, 0, 2) * c02_shear_retention;
356 708048 : _local_elastic_vector[8] = _elasticity_tensor[_qp](0, 1, 0, 1) * c01_shear_retention;
357 : }
358 :
359 : // Filling with 9 components is sufficient because these are the only nonzero entries
360 : // for isotropic or orthotropic materials.
361 355192 : _local_elasticity_tensor.fillFromInputVector(_local_elastic_vector,
362 : ADRankFourTensor::symmetric9);
363 :
364 : // Rotate the modified elasticity tensor back into global coordinates
365 355192 : _local_elasticity_tensor.rotate(R);
366 : }
367 : }
368 : if (!cracking_locally_active)
369 33592 : _local_elasticity_tensor = _elasticity_tensor[_qp];
370 388784 : }
371 :
372 : void
373 415472 : ADComputeSmearedCrackingStress::updateCrackingStateAndStress()
374 : {
375 : const ADReal youngs_modulus =
376 415472 : ElasticityTensorTools::getIsotropicYoungsModulus(_elasticity_tensor[_qp]);
377 :
378 415472 : ADReal cracking_stress = _cracking_stress[_qp];
379 :
380 415472 : if (cracking_stress > 0)
381 : {
382 : // Initializing crack states
383 415448 : _crack_rotation[_qp] = _crack_rotation_old[_qp];
384 :
385 1661792 : for (unsigned i = 0; i < 3; ++i)
386 : {
387 1246344 : _crack_max_strain[_qp](i) = _crack_max_strain_old[_qp](i);
388 1246344 : _crack_initiation_strain[_qp](i) = _crack_initiation_strain_old[_qp](i);
389 1246344 : _crack_damage[_qp](i) = _crack_damage_old[_qp](i);
390 : }
391 :
392 : // Compute crack orientations: updated _crack_rotation[_qp] based on current strain
393 : ADRealVectorValue strain_in_crack_dir;
394 415448 : computeCrackStrainAndOrientation(strain_in_crack_dir);
395 :
396 1661792 : for (unsigned i = 0; i < 3; ++i)
397 : {
398 1246344 : if (strain_in_crack_dir(i) > _crack_max_strain[_qp](i))
399 422443 : _crack_max_strain[_qp](i) = strain_in_crack_dir(i);
400 : }
401 :
402 : // Check for new cracks.
403 : // Rotate stress to cracked orientation.
404 415448 : const ADRankTwoTensor & R = _crack_rotation[_qp];
405 415448 : ADRankTwoTensor sigmaPrime(_stress[_qp]);
406 415448 : sigmaPrime.rotate(R.transpose()); // stress in crack coordinates
407 :
408 : unsigned int num_cracks = 0;
409 1661792 : for (unsigned int i = 0; i < 3; ++i)
410 : {
411 1246344 : if (_crack_damage_old[_qp](i) > 0.0)
412 559360 : ++num_cracks;
413 : }
414 :
415 : bool cracked(false);
416 : ADRealVectorValue sigma;
417 : mooseAssert(_softening_models.size() == 3, "Must have 3 softening models");
418 1661792 : for (unsigned int i = 0; i < 3; ++i)
419 : {
420 1246344 : sigma(i) = sigmaPrime(i, i);
421 :
422 1246344 : ADReal stiffness_ratio = 1.0 - _crack_damage[_qp](i);
423 :
424 1246344 : const bool pre_existing_crack = (_crack_damage_old[_qp](i) > 0.0);
425 : const bool met_stress_criterion = (sigma(i) > cracking_stress);
426 1246344 : const bool loading_existing_crack = (strain_in_crack_dir(i) >= _crack_max_strain[_qp](i));
427 1246344 : const bool allowed_to_crack = (pre_existing_crack || num_cracks < _max_cracks);
428 : bool new_crack = false;
429 :
430 1246344 : cracked |= pre_existing_crack;
431 :
432 : // Adjustments for newly created cracks
433 1246344 : if (met_stress_criterion && !pre_existing_crack && allowed_to_crack)
434 : {
435 : new_crack = true;
436 3944 : ++num_cracks;
437 :
438 : // Assume Poisson's ratio drops to zero for this direction. Stiffness is then Young's
439 : // modulus.
440 3944 : _crack_initiation_strain[_qp](i) = cracking_stress / youngs_modulus;
441 :
442 3944 : if (_crack_max_strain[_qp](i) < _crack_initiation_strain[_qp](i))
443 256 : _crack_max_strain[_qp](i) = _crack_initiation_strain[_qp](i);
444 : }
445 :
446 : // Update stress and stiffness ratio according to specified crack release model
447 1242400 : if (new_crack || (pre_existing_crack && loading_existing_crack))
448 : {
449 : cracked = true;
450 372122 : _softening_models[i]->computeCrackingRelease(sigma(i),
451 : stiffness_ratio,
452 : strain_in_crack_dir(i),
453 372122 : _crack_initiation_strain[_qp](i),
454 372122 : _crack_max_strain[_qp](i),
455 : cracking_stress,
456 : youngs_modulus);
457 744244 : _crack_damage[_qp](i) = 1.0 - stiffness_ratio;
458 : }
459 :
460 694118 : else if (cracked && _cracking_neg_fraction > 0 &&
461 1748444 : _crack_initiation_strain[_qp](i) * _cracking_neg_fraction > strain_in_crack_dir(i) &&
462 874222 : -_crack_initiation_strain[_qp](i) * _cracking_neg_fraction < strain_in_crack_dir(i))
463 : {
464 0 : const ADReal etr = _cracking_neg_fraction * _crack_initiation_strain[_qp](i);
465 0 : const ADReal Eo = cracking_stress / _crack_initiation_strain[_qp](i);
466 0 : const ADReal Ec = Eo * (1.0 - _crack_damage_old[_qp](i));
467 0 : const ADReal a = (Ec - Eo) / (4.0 * etr);
468 0 : const ADReal b = 0.5 * (Ec + Eo);
469 0 : const ADReal c = 0.25 * (Ec - Eo) * etr;
470 0 : sigma(i) = (a * strain_in_crack_dir(i) + b) * strain_in_crack_dir(i) + c;
471 : }
472 : }
473 :
474 415448 : if (cracked)
475 390016 : updateStressTensorForCracking(_stress[_qp], sigma);
476 : }
477 :
478 830944 : _crack_flags[_qp](0) = 1.0 - _crack_damage[_qp](2);
479 830944 : _crack_flags[_qp](1) = 1.0 - _crack_damage[_qp](1);
480 830944 : _crack_flags[_qp](2) = 1.0 - _crack_damage[_qp](0);
481 415472 : }
482 :
483 : void
484 415448 : ADComputeSmearedCrackingStress::computeCrackStrainAndOrientation(
485 : ADRealVectorValue & strain_in_crack_dir)
486 : {
487 : // The rotation tensor is ordered such that directions for pre-existing cracks appear first
488 : // in the list of columns. For example, if there is one existing crack, its direction is in the
489 : // first column in the rotation tensor.
490 415448 : const unsigned int num_known_dirs = getNumKnownCrackDirs();
491 :
492 415448 : if (num_known_dirs == 0)
493 : {
494 21896 : std::vector<ADReal> eigval(3, 0.0);
495 21896 : ADRankTwoTensor eigvec;
496 :
497 21896 : _elastic_strain[_qp].symmetricEigenvaluesEigenvectors(eigval, eigvec);
498 :
499 : // If the elastic strain is beyond the cracking strain, save the eigen vectors as
500 : // the rotation tensor. Reverse their order so that the third principal strain
501 : // (most tensile) will correspond to the first crack.
502 21896 : _crack_rotation[_qp].fillColumn(0, eigvec.column(2));
503 21896 : _crack_rotation[_qp].fillColumn(1, eigvec.column(1));
504 21896 : _crack_rotation[_qp].fillColumn(2, eigvec.column(0));
505 :
506 21896 : strain_in_crack_dir(0) = eigval[2];
507 21896 : strain_in_crack_dir(1) = eigval[1];
508 21896 : strain_in_crack_dir(2) = eigval[0];
509 : }
510 393552 : else if (num_known_dirs == 1)
511 : {
512 : // This is easily the most complicated case.
513 : // 1. Rotate the elastic strain to the orientation associated with the known
514 : // crack.
515 : // 2. Extract the lower 2x2 block into a separate tensor.
516 : // 3. Run the eigen solver on the result.
517 : // 4. Update the rotation tensor to reflect the effect of the 2 eigenvectors.
518 :
519 : // 1.
520 140112 : const ADRankTwoTensor & R = _crack_rotation[_qp];
521 140112 : RankTwoTensor ePrime(raw_value(_elastic_strain[_qp]));
522 280224 : ePrime.rotate(raw_value(R.transpose())); // elastic strain in crack coordinates
523 :
524 : // 2.
525 140112 : ColumnMajorMatrix e2x2(2, 2);
526 140112 : e2x2(0, 0) = ePrime(1, 1);
527 140112 : e2x2(1, 0) = ePrime(2, 1);
528 140112 : e2x2(0, 1) = ePrime(1, 2);
529 140112 : e2x2(1, 1) = ePrime(2, 2);
530 :
531 : // 3.
532 140112 : ColumnMajorMatrix e_val2x1(2, 1);
533 140112 : ColumnMajorMatrix e_vec2x2(2, 2);
534 140112 : e2x2.eigen(e_val2x1, e_vec2x2);
535 :
536 : // 4.
537 : ADRankTwoTensor eigvec(
538 140112 : 1.0, 0.0, 0.0, 0.0, e_vec2x2(0, 1), e_vec2x2(1, 1), 0.0, e_vec2x2(0, 0), e_vec2x2(1, 0));
539 :
540 140112 : _crack_rotation[_qp] = _crack_rotation_old[_qp] * eigvec; // Roe implementation
541 :
542 140112 : strain_in_crack_dir(0) = ePrime(0, 0);
543 140112 : strain_in_crack_dir(1) = e_val2x1(1, 0);
544 140112 : strain_in_crack_dir(2) = e_val2x1(0, 0);
545 : }
546 253440 : else if (num_known_dirs == 2 || num_known_dirs == 3)
547 : {
548 : // Rotate to cracked orientation and pick off the strains in the rotated
549 : // coordinate directions.
550 253440 : const ADRankTwoTensor & R = _crack_rotation[_qp];
551 253440 : ADRankTwoTensor ePrime(_elastic_strain[_qp]);
552 253440 : ePrime.rotate(R.transpose()); // elastic strain in crack coordinates
553 :
554 253440 : strain_in_crack_dir(0) = ePrime(0, 0);
555 253440 : strain_in_crack_dir(1) = ePrime(1, 1);
556 253440 : strain_in_crack_dir(2) = ePrime(2, 2);
557 : }
558 : else
559 0 : mooseError("Invalid number of known crack directions");
560 415448 : }
561 :
562 : unsigned int
563 415448 : ADComputeSmearedCrackingStress::getNumKnownCrackDirs() const
564 : {
565 : unsigned int num_known_dirs = 0;
566 1661792 : for (unsigned int i = 0; i < 3; ++i)
567 : {
568 1246344 : if (_crack_damage_old[_qp](i) > 0.0 || _prescribed_crack_directions.size() >= i + 1)
569 839200 : ++num_known_dirs;
570 : }
571 415448 : return num_known_dirs;
572 : }
573 :
574 : void
575 390016 : ADComputeSmearedCrackingStress::updateStressTensorForCracking(ADRankTwoTensor & tensor,
576 : const ADRealVectorValue & sigma)
577 : {
578 : // Get transformation matrix
579 390016 : const ADRankTwoTensor & R = _crack_rotation[_qp];
580 : // Rotate to crack frame
581 390016 : tensor.rotate(R.transpose());
582 :
583 : // Reset stress if cracked
584 1560064 : for (unsigned int i = 0; i < 3; ++i)
585 1170048 : if (_crack_damage[_qp](i) > 0.0)
586 : {
587 561848 : const ADReal stress_correction_ratio = (tensor(i, i) - sigma(i)) / tensor(i, i);
588 561848 : if (stress_correction_ratio > _max_stress_correction)
589 0 : tensor(i, i) *= (1.0 - _max_stress_correction);
590 561848 : else if (stress_correction_ratio < -_max_stress_correction)
591 0 : tensor(i, i) *= (1.0 + _max_stress_correction);
592 : else
593 561848 : tensor(i, i) = sigma(i);
594 : }
595 :
596 : // Rotate back to global frame
597 390016 : tensor.rotate(R);
598 390016 : }
599 :
600 : bool
601 415472 : ADComputeSmearedCrackingStress::previouslyCracked()
602 : {
603 599152 : for (unsigned int i = 0; i < 3; ++i)
604 572464 : if (_crack_damage_old[_qp](i) > 0.0)
605 : return true;
606 : return false;
607 : }
|