https://mooseframework.inl.gov
ADComputeSmearedCrackingStress.C
Go to the documentation of this file.
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 
11 #include "ElasticityTensorTools.h"
12 #include "StressUpdateBase.h"
13 #include "Conversion.h"
14 
16 
19 {
21  params.addClassDescription(
22  "Compute stress using a fixed smeared cracking model. Uses automatic differentiation");
23  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  params.addRequiredCoupledVar(
30  "cracking_stress",
31  "The stress threshold beyond which cracking occurs. Negative values prevent cracking.");
32  MultiMooseEnum direction("x y z");
33  params.addParam<MultiMooseEnum>(
34  "prescribed_crack_directions", direction, "Prescribed directions of first cracks");
35  params.addParam<unsigned int>(
36  "max_cracks", 3, "The maximum number of cracks allowed at a material point.");
37  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  params.addParam<Real>(
44  "max_stress_correction",
45  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  params.addRangeCheckedParam<Real>(
51  "shear_retention_factor",
52  0.0,
53  "shear_retention_factor>=0 & shear_retention_factor<=1.0",
54  "Fraction of original shear stiffness to be retained after cracking");
55  params.set<std::vector<MaterialName>>("inelastic_models") = {};
56 
57  MooseEnum crackedElasticityType("DIAGONAL FULL", "DIAGONAL");
58  crackedElasticityType.addDocumentation(
59  "DIAGONAL", "Zero out terms coupling with directions orthogonal to a crack (legacy)");
60  crackedElasticityType.addDocumentation(
61  "FULL", "Consistently scale all entries based on damage (recommended)");
62  params.addParam<MooseEnum>(
63  "cracked_elasticity_type",
64  crackedElasticityType,
65  "Method to modify the local elasticity tensor to account for cracking");
66 
67  return params;
68 }
69 
72  _cracking_stress(adCoupledValue("cracking_stress")),
73  _max_cracks(getParam<unsigned int>("max_cracks")),
74  _cracking_neg_fraction(getParam<Real>("cracking_neg_fraction")),
75  _shear_retention_factor(getParam<Real>("shear_retention_factor")),
76  _max_stress_correction(getParam<Real>("max_stress_correction")),
77  _cracked_elasticity_type(
78  getParam<MooseEnum>("cracked_elasticity_type").getEnum<CrackedElasticityType>()),
79  _crack_damage(declareADProperty<RealVectorValue>(_base_name + "crack_damage")),
80  _crack_damage_old(getMaterialPropertyOld<RealVectorValue>(_base_name + "crack_damage")),
81  _crack_flags(declareADProperty<RealVectorValue>(_base_name + "crack_flags")),
82  _crack_rotation(declareADProperty<RankTwoTensor>(_base_name + "crack_rotation")),
83  _crack_rotation_old(getMaterialPropertyOld<RankTwoTensor>(_base_name + "crack_rotation")),
84  _crack_initiation_strain(
85  declareADProperty<RealVectorValue>(_base_name + "crack_initiation_strain")),
86  _crack_initiation_strain_old(
87  getMaterialPropertyOld<RealVectorValue>(_base_name + "crack_initiation_strain")),
88  _crack_max_strain(declareADProperty<RealVectorValue>(_base_name + "crack_max_strain")),
89  _crack_max_strain_old(getMaterialPropertyOld<RealVectorValue>(_base_name + "crack_max_strain"))
90 {
91  MultiMooseEnum prescribed_crack_directions =
92  getParam<MultiMooseEnum>("prescribed_crack_directions");
93  if (prescribed_crack_directions.size() > 0)
94  {
95  if (prescribed_crack_directions.size() > 3)
96  mooseError("A maximum of three crack directions may be specified");
97  for (unsigned int i = 0; i < prescribed_crack_directions.size(); ++i)
98  {
99  for (unsigned int j = 0; j < i; ++j)
100  if (prescribed_crack_directions[i] == prescribed_crack_directions[j])
101  mooseError("Entries in 'prescribed_crack_directions' cannot be repeated");
103  static_cast<unsigned int>(prescribed_crack_directions.get(i)));
104  }
105 
106  // Fill in the last remaining direction if 2 are specified
107  if (_prescribed_crack_directions.size() == 2)
108  {
109  std::set<unsigned int> available_dirs = {0, 1, 2};
110  for (auto dir : _prescribed_crack_directions)
111  if (available_dirs.erase(dir) != 1)
112  mooseError("Invalid prescribed crack direction:" + Moose::stringify(dir));
113  if (available_dirs.size() != 1)
114  mooseError("Error in finding remaining available crack direction");
115  _prescribed_crack_directions.push_back(*available_dirs.begin());
116  }
117  }
118  if (!isParamSetByUser("cracked_elasticity_type"))
119  paramWarning(
120  "cracked_elasticity_type",
121  "Defaulting to the legacy option of 'DIAGONAL', but the 'FULL' option is preferred");
122 
123  _local_elastic_vector.resize(9);
124 }
125 
126 void
128 {
130 
131  _crack_damage[_qp] = 0.0;
132 
134  _crack_max_strain[_qp](0) = 0.0;
135 
136  switch (_prescribed_crack_directions.size())
137  {
138  case 0:
139  {
141  break;
142  }
143  case 1:
144  {
145  _crack_rotation[_qp].zero();
146  switch (_prescribed_crack_directions[0])
147  {
148  case 0:
149  {
150  _crack_rotation[_qp](0, 0) = 1.0;
151  _crack_rotation[_qp](1, 1) = 1.0;
152  _crack_rotation[_qp](2, 2) = 1.0;
153  break;
154  }
155  case 1:
156  {
157  _crack_rotation[_qp](1, 0) = 1.0;
158  _crack_rotation[_qp](0, 1) = 1.0;
159  _crack_rotation[_qp](2, 2) = 1.0;
160  break;
161  }
162  case 2:
163  {
164  _crack_rotation[_qp](2, 0) = 1.0;
165  _crack_rotation[_qp](0, 1) = 1.0;
166  _crack_rotation[_qp](1, 2) = 1.0;
167  break;
168  }
169  }
170  break;
171  }
172  case 2:
173  {
174  mooseError("Number of prescribed crack directions cannot be 2");
175  break;
176  }
177  case 3:
178  {
179  for (unsigned int i = 0; i < _prescribed_crack_directions.size(); ++i)
180  {
181  ADRealVectorValue crack_dir_vec;
182  crack_dir_vec(_prescribed_crack_directions[i]) = 1.0;
183  _crack_rotation[_qp].fillColumn(i, crack_dir_vec);
184  }
185  }
186  }
187 }
188 
189 void
191 {
193 
195  mooseError("ADComputeSmearedCrackingStress requires that the elasticity tensor be "
196  "guaranteed isotropic");
197 
198  std::vector<MaterialName> soft_matls = getParam<std::vector<MaterialName>>("softening_models");
199  for (auto soft_matl : soft_matls)
200  {
202  dynamic_cast<ADSmearedCrackSofteningBase *>(&getMaterialByName(soft_matl));
203  if (scsb)
204  _softening_models.push_back(scsb);
205  else
206  paramError("softening_models", "Model " + soft_matl + " is not a softening model");
207  }
208  if (_softening_models.size() == 1)
209  {
210  // Reuse the same model in all 3 directions
211  _softening_models.push_back(_softening_models[0]);
212  _softening_models.push_back(_softening_models[0]);
213  }
214  else if (_softening_models.size() != 3)
215  paramError("softening_models", "Either 1 or 3 softening models must be specified");
216 }
217 
218 void
220 {
221 
222  if (!previouslyCracked())
224  else
225  {
227 
228  // Propagate behavior from the (now inactive) inelastic models
230  for (auto model : _models)
231  {
232  model->setQp(_qp);
233  model->propagateQpStatefulProperties();
234  }
235 
236  // Since the other inelastic models are inactive, they will not limit the time step
237  _material_timestep_limit[_qp] = std::numeric_limits<Real>::max();
238 
239  // update _local_elasticity_tensor based on cracking state in previous time step
241 
242  // Calculate stress in intermediate configuration
244  }
245 
246  // compute crack status and adjust stress
248 
250  {
253  }
254 }
255 
256 void
258 {
259  const ADReal youngs_modulus =
261 
262  bool cracking_locally_active = false;
263 
264  const ADReal cracking_stress = _cracking_stress[_qp];
265 
266  if (cracking_stress > 0)
267  {
268  ADRealVectorValue stiffness_ratio_local(1.0, 1.0, 1.0);
271  ePrime.rotate(R.transpose());
272 
273  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  if (_crack_damage_old[_qp](i) > 0.0)
277  {
278  if (_cracking_neg_fraction == 0.0 && MooseUtils::absoluteFuzzyLessThan(ePrime(i, i), 0.0))
279  stiffness_ratio_local(i) = 1.0;
280  else if (_cracking_neg_fraction > 0.0 &&
283  {
285  const ADReal Eo = cracking_stress / _crack_initiation_strain_old[_qp](i);
286  const ADReal Ec = Eo * (1.0 - _crack_damage_old[_qp](i));
287  const ADReal a = (Ec - Eo) / (4 * etr);
288  const ADReal b = (Ec + Eo) / 2;
289  // Compute the ratio of the current transition stiffness to the original stiffness
290  stiffness_ratio_local(i) = (2.0 * a * etr + b) / Eo;
291  cracking_locally_active = true;
292  }
293  else
294  {
295  stiffness_ratio_local(i) = (1.0 - _crack_damage_old[_qp](i));
296  cracking_locally_active = true;
297  }
298  }
299  }
300 
301  if (cracking_locally_active)
302  {
303  // Update the elasticity tensor in the crack coordinate system
305  {
306  const bool c0_coupled = MooseUtils::absoluteFuzzyEqual(stiffness_ratio_local(0), 1.0);
307  const bool c1_coupled = MooseUtils::absoluteFuzzyEqual(stiffness_ratio_local(1), 1.0);
308  const bool c2_coupled = MooseUtils::absoluteFuzzyEqual(stiffness_ratio_local(2), 1.0);
309 
310  const ADReal c01 = (c0_coupled && c1_coupled ? 1.0 : 0.0);
311  const ADReal c02 = (c0_coupled && c2_coupled ? 1.0 : 0.0);
312  const ADReal c12 = (c1_coupled && c2_coupled ? 1.0 : 0.0);
313 
314  const ADReal c01_shear_retention =
315  (c0_coupled && c1_coupled ? 1.0 : _shear_retention_factor);
316  const ADReal c02_shear_retention =
317  (c0_coupled && c2_coupled ? 1.0 : _shear_retention_factor);
318  const ADReal c12_shear_retention =
319  (c1_coupled && c2_coupled ? 1.0 : _shear_retention_factor);
320 
321  _local_elastic_vector[0] = (c0_coupled ? _elasticity_tensor[_qp](0, 0, 0, 0)
322  : stiffness_ratio_local(0) * youngs_modulus);
323  _local_elastic_vector[1] = _elasticity_tensor[_qp](0, 0, 1, 1) * c01;
324  _local_elastic_vector[2] = _elasticity_tensor[_qp](0, 0, 2, 2) * c02;
325  _local_elastic_vector[3] = (c1_coupled ? _elasticity_tensor[_qp](1, 1, 1, 1)
326  : stiffness_ratio_local(1) * youngs_modulus);
327  _local_elastic_vector[4] = _elasticity_tensor[_qp](1, 1, 2, 2) * c12;
328  _local_elastic_vector[5] = (c2_coupled ? _elasticity_tensor[_qp](2, 2, 2, 2)
329  : stiffness_ratio_local(2) * youngs_modulus);
330  _local_elastic_vector[6] = _elasticity_tensor[_qp](1, 2, 1, 2) * c12_shear_retention;
331  _local_elastic_vector[7] = _elasticity_tensor[_qp](0, 2, 0, 2) * c02_shear_retention;
332  _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  const ADReal c01_shear_retention = std::max(c01, _shear_retention_factor);
345  const ADReal c02_shear_retention = std::max(c02, _shear_retention_factor);
346  const ADReal c12_shear_retention = std::max(c12, _shear_retention_factor);
347 
348  _local_elastic_vector[0] = _elasticity_tensor[_qp](0, 0, 0, 0) * c0;
349  _local_elastic_vector[1] = _elasticity_tensor[_qp](0, 0, 1, 1) * c01;
350  _local_elastic_vector[2] = _elasticity_tensor[_qp](0, 0, 2, 2) * c02;
351  _local_elastic_vector[3] = _elasticity_tensor[_qp](1, 1, 1, 1) * c1;
352  _local_elastic_vector[4] = _elasticity_tensor[_qp](1, 1, 2, 2) * c12;
353  _local_elastic_vector[5] = _elasticity_tensor[_qp](2, 2, 2, 2) * c2;
354  _local_elastic_vector[6] = _elasticity_tensor[_qp](1, 2, 1, 2) * c12_shear_retention;
355  _local_elastic_vector[7] = _elasticity_tensor[_qp](0, 2, 0, 2) * c02_shear_retention;
356  _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.
363 
364  // Rotate the modified elasticity tensor back into global coordinates
366  }
367  }
368  if (!cracking_locally_active)
370 }
371 
372 void
374 {
375  const ADReal youngs_modulus =
377 
378  ADReal cracking_stress = _cracking_stress[_qp];
379 
380  if (cracking_stress > 0)
381  {
382  // Initializing crack states
384 
385  for (unsigned i = 0; i < 3; ++i)
386  {
390  }
391 
392  // Compute crack orientations: updated _crack_rotation[_qp] based on current strain
393  ADRealVectorValue strain_in_crack_dir;
394  computeCrackStrainAndOrientation(strain_in_crack_dir);
395 
396  for (unsigned i = 0; i < 3; ++i)
397  {
398  if (strain_in_crack_dir(i) > _crack_max_strain[_qp](i))
399  _crack_max_strain[_qp](i) = strain_in_crack_dir(i);
400  }
401 
402  // Check for new cracks.
403  // Rotate stress to cracked orientation.
405  ADRankTwoTensor sigmaPrime(_stress[_qp]);
406  sigmaPrime.rotate(R.transpose()); // stress in crack coordinates
407 
408  unsigned int num_cracks = 0;
409  for (unsigned int i = 0; i < 3; ++i)
410  {
411  if (_crack_damage_old[_qp](i) > 0.0)
412  ++num_cracks;
413  }
414 
415  bool cracked(false);
417  mooseAssert(_softening_models.size() == 3, "Must have 3 softening models");
418  for (unsigned int i = 0; i < 3; ++i)
419  {
420  sigma(i) = sigmaPrime(i, i);
421 
422  ADReal stiffness_ratio = 1.0 - _crack_damage[_qp](i);
423 
424  const bool pre_existing_crack = (_crack_damage_old[_qp](i) > 0.0);
425  const bool met_stress_criterion = (sigma(i) > cracking_stress);
426  const bool loading_existing_crack = (strain_in_crack_dir(i) >= _crack_max_strain[_qp](i));
427  const bool allowed_to_crack = (pre_existing_crack || num_cracks < _max_cracks);
428  bool new_crack = false;
429 
430  cracked |= pre_existing_crack;
431 
432  // Adjustments for newly created cracks
433  if (met_stress_criterion && !pre_existing_crack && allowed_to_crack)
434  {
435  new_crack = true;
436  ++num_cracks;
437 
438  // Assume Poisson's ratio drops to zero for this direction. Stiffness is then Young's
439  // modulus.
440  _crack_initiation_strain[_qp](i) = cracking_stress / youngs_modulus;
441 
444  }
445 
446  // Update stress and stiffness ratio according to specified crack release model
447  if (new_crack || (pre_existing_crack && loading_existing_crack))
448  {
449  cracked = true;
450  _softening_models[i]->computeCrackingRelease(sigma(i),
451  stiffness_ratio,
452  strain_in_crack_dir(i),
455  cracking_stress,
456  youngs_modulus);
457  _crack_damage[_qp](i) = 1.0 - stiffness_ratio;
458  }
459 
460  else if (cracked && _cracking_neg_fraction > 0 &&
461  _crack_initiation_strain[_qp](i) * _cracking_neg_fraction > strain_in_crack_dir(i) &&
462  -_crack_initiation_strain[_qp](i) * _cracking_neg_fraction < strain_in_crack_dir(i))
463  {
465  const ADReal Eo = cracking_stress / _crack_initiation_strain[_qp](i);
466  const ADReal Ec = Eo * (1.0 - _crack_damage_old[_qp](i));
467  const ADReal a = (Ec - Eo) / (4.0 * etr);
468  const ADReal b = 0.5 * (Ec + Eo);
469  const ADReal c = 0.25 * (Ec - Eo) * etr;
470  sigma(i) = (a * strain_in_crack_dir(i) + b) * strain_in_crack_dir(i) + c;
471  }
472  }
473 
474  if (cracked)
476  }
477 
478  _crack_flags[_qp](0) = 1.0 - _crack_damage[_qp](2);
479  _crack_flags[_qp](1) = 1.0 - _crack_damage[_qp](1);
480  _crack_flags[_qp](2) = 1.0 - _crack_damage[_qp](0);
481 }
482 
483 void
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  const unsigned int num_known_dirs = getNumKnownCrackDirs();
491 
492  if (num_known_dirs == 0)
493  {
494  std::vector<ADReal> eigval(3, 0.0);
495  ADRankTwoTensor eigvec;
496 
497  _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  _crack_rotation[_qp].fillColumn(0, eigvec.column(2));
503  _crack_rotation[_qp].fillColumn(1, eigvec.column(1));
504  _crack_rotation[_qp].fillColumn(2, eigvec.column(0));
505 
506  strain_in_crack_dir(0) = eigval[2];
507  strain_in_crack_dir(1) = eigval[1];
508  strain_in_crack_dir(2) = eigval[0];
509  }
510  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.
522  ePrime.rotate(raw_value(R.transpose())); // elastic strain in crack coordinates
523 
524  // 2.
525  ColumnMajorMatrix e2x2(2, 2);
526  e2x2(0, 0) = ePrime(1, 1);
527  e2x2(1, 0) = ePrime(2, 1);
528  e2x2(0, 1) = ePrime(1, 2);
529  e2x2(1, 1) = ePrime(2, 2);
530 
531  // 3.
532  ColumnMajorMatrix e_val2x1(2, 1);
533  ColumnMajorMatrix e_vec2x2(2, 2);
534  e2x2.eigen(e_val2x1, e_vec2x2);
535 
536  // 4.
537  ADRankTwoTensor eigvec(
538  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  _crack_rotation[_qp] = _crack_rotation_old[_qp] * eigvec; // Roe implementation
541 
542  strain_in_crack_dir(0) = ePrime(0, 0);
543  strain_in_crack_dir(1) = e_val2x1(1, 0);
544  strain_in_crack_dir(2) = e_val2x1(0, 0);
545  }
546  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.
552  ePrime.rotate(R.transpose()); // elastic strain in crack coordinates
553 
554  strain_in_crack_dir(0) = ePrime(0, 0);
555  strain_in_crack_dir(1) = ePrime(1, 1);
556  strain_in_crack_dir(2) = ePrime(2, 2);
557  }
558  else
559  mooseError("Invalid number of known crack directions");
560 }
561 
562 unsigned int
564 {
565  unsigned int num_known_dirs = 0;
566  for (unsigned int i = 0; i < 3; ++i)
567  {
568  if (_crack_damage_old[_qp](i) > 0.0 || _prescribed_crack_directions.size() >= i + 1)
569  ++num_known_dirs;
570  }
571  return num_known_dirs;
572 }
573 
574 void
576  const ADRealVectorValue & sigma)
577 {
578  // Get transformation matrix
580  // Rotate to crack frame
581  tensor.rotate(R.transpose());
582 
583  // Reset stress if cracked
584  for (unsigned int i = 0; i < 3; ++i)
585  if (_crack_damage[_qp](i) > 0.0)
586  {
587  const ADReal stress_correction_ratio = (tensor(i, i) - sigma(i)) / tensor(i, i);
588  if (stress_correction_ratio > _max_stress_correction)
589  tensor(i, i) *= (1.0 - _max_stress_correction);
590  else if (stress_correction_ratio < -_max_stress_correction)
591  tensor(i, i) *= (1.0 + _max_stress_correction);
592  else
593  tensor(i, i) = sigma(i);
594  }
595 
596  // Rotate back to global frame
597  tensor.rotate(R);
598 }
599 
600 bool
602 {
603  for (unsigned int i = 0; i < 3; ++i)
604  if (_crack_damage_old[_qp](i) > 0.0)
605  return true;
606  return false;
607 }
void updateStressTensorForCracking(ADRankTwoTensor &tensor, const ADRealVectorValue &sigma)
Updates the full stress tensor to account for the effect of cracking using the provided stresses in t...
const ADReal _cracking_neg_fraction
Defines transition to changed stiffness during unloading.
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
ADMaterialProperty< RealVectorValue > & _crack_initiation_strain
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
registerMooseObject("SolidMechanicsApp", ADComputeSmearedCrackingStress)
const MaterialProperty< RankTwoTensor > & _inelastic_strain_old
old value of inelastic strain
ADMaterialProperty< RankTwoTensor > & _inelastic_strain
The sum of the inelastic strains that come from the plastic models.
ADMaterialProperty< R2 > & _stress
The stress tensor to be calculated.
const MaterialProperty< RealVectorValue > & _crack_max_strain_old
virtual void computeQpStressIntermediateConfiguration()
Compute the stress for the current QP, but do not rotate tensors from the intermediate configuration ...
ADComputeMultipleInelasticStress computes the stress and a decomposition of the strain into elastic a...
const bool _perform_finite_strain_rotations
after updateQpState, rotate the stress, elastic_strain, and inelastic_strain using _rotation_incremen...
const ADReal _shear_retention_factor
Controls the amount of shear retained.
enum ADComputeSmearedCrackingStress::CrackedElasticityType _cracked_elasticity_type
std::vector< ADSmearedCrackSofteningBase * > _softening_models
The user-supplied list of softening models to be used in the 3 crack directions.
T & set(const std::string &name, bool quiet_mode=false)
auto raw_value(const Eigen::Map< T > &in)
unsigned int size() const
std::vector< unsigned int > _prescribed_crack_directions
User-prescribed cracking directions.
const unsigned int _max_cracks
Maximum number of cracks permitted at a material point.
void fillFromInputVector(const std::vector< T > &input, FillMethod fill_method)
static RankTwoTensorTempl Identity()
void eigen(ColumnMajorMatrixTempl< Real > &eval, ColumnMajorMatrixTempl< Real > &evec) const
void addRequiredParam(const std::string &name, const std::string &doc_string)
void rotate(const RankTwoTensorTempl< T > &R)
virtual unsigned int getNumKnownCrackDirs() const
Get the number of known crack directions.
ADMaterialProperty< RealVectorValue > & _crack_damage
const MaterialProperty< RankTwoTensor > & _crack_rotation_old
unsigned int _qp
const ADMaterialProperty< RankTwoTensor > & _rotation_increment
void updateLocalElasticityTensor()
Update the local elasticity tensor (_local_elasticity_tensor) due to the effects of cracking...
const MaterialProperty< RealVectorValue > & _crack_initiation_strain_old
ADMaterialProperty< RealVectorValue > & _crack_flags
Vector of values going from 1 to 0 as crack damage accumulates.
const ADMaterialProperty< R4 > & _elasticity_tensor
Elasticity tensor material property.
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
libMesh::VectorValue< T > column(const unsigned int i) const
bool absoluteFuzzyLessThan(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
ADMaterialProperty< RealVectorValue > & _crack_max_strain
ADSmearedCrackSofteningBase is the base class for a set of models that define the softening behavior ...
void paramError(const std::string &param, Args... args) const
ADComputeSmearedCrackingStress computes the stress for a finite strain material with smeared cracking...
std::string stringify(const T &t)
void rotate(const TypeTensor< T > &R)
const PertinentGeochemicalSystem model(database, {"H2O", "H+", "HCO3-", "O2(aq)", "Ca++", ">(s)FeOH", "radius_neg1", "radius_neg1.5"}, {"Calcite"}, {}, {"Calcite_asdf"}, {"CH4(aq)"}, {">(s)FeOCa+"}, "O2(aq)", "e-")
MaterialBase & getMaterialByName(const std::string &name, bool no_warn=false, bool no_dep=false)
std::vector< ADStressUpdateBase * > _models
The user supplied list of inelastic models to use in the simulation.
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
static const std::string R
Definition: NS.h:162
ADRankFourTensor _local_elasticity_tensor
Variables used by multiple methods within the calculation for a single material point.
bool isParamSetByUser(const std::string &nm) const
unsigned int get(unsigned int i) const
bool previouslyCracked()
Check to see whether there was cracking in any diretion in the previous time step.
ADComputeSmearedCrackingStress(const InputParameters &parameters)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addDocumentation(const std::string &name, const std::string &doc)
CrackedElasticityType
Enum defining the method used to adjust the elasticity tensor for cracking.
void mooseError(Args &&... args) const
virtual void updateCrackingStateAndStress()
Update all cracking-related state variables and the stress tensor due to cracking in all directions...
void addClassDescription(const std::string &doc_string)
ADMaterialProperty< RankTwoTensor > & _crack_rotation
std::vector< ADReal > _local_elastic_vector
Vector helper to update local elasticity tensor.
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const ADReal _max_stress_correction
Controls the maximum amount that the damaged elastic stress is corrected to folow the release model d...
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
ADMaterialProperty< R2 > & _elastic_strain
void paramWarning(const std::string &param, Args... args) const
bool hasGuaranteedMaterialProperty(const MaterialPropertyName &prop, Guarantee guarantee)
const MaterialProperty< RealVectorValue > & _crack_damage_old
T getIsotropicYoungsModulus(const RankFourTensorTempl< T > &elasticity_tensor)
Get the Young&#39;s modulus for an isotropic elasticity tensor param elasticity_tensor the tensor (must b...
const ADVariableValue & _cracking_stress
Input parameters for smeared crack models.
const MaterialProperty< R2 > & _elastic_strain_old
The old elastic strain is used to calculate the old stress in the case of variable elasticity tensors...
void ErrorVector unsigned int
void computeCrackStrainAndOrientation(ADRealVectorValue &strain_in_crack_dir)
Compute the crack strain in the crack coordinate system.
virtual void finiteStrainRotation()
Rotate _elastic_strain, _stress, and _inelastic_strain to the new configuration.