11 #include "RankFourTensor.h"
18 InputParameters params = emptyInputParameters();
19 MooseEnum specialIC(
"none rock joint",
"none");
20 params.addParam<MooseEnum>(
"specialIC",
22 "For certain combinations of plastic models, the set of active "
23 "constraints can be initialized optimally. 'none': no special "
24 "initialization is performed. For all other choices, the "
25 "plastic_models must be chosen to have the following types. 'rock': "
26 "'TensileMulti MohrCoulombMulti'. 'joint': 'WeakPlaneTensile "
28 params.addParam<std::vector<UserObjectName>>(
30 "List of names of user objects that define the plastic models that could "
31 "be active for this material. If no plastic_models are provided, only "
32 "elasticity will be used.");
33 params.addClassDescription(
"RawComponentAssembler class to calculate yield functions, etc, used "
34 "in multi-surface finite-strain plasticity");
39 const MooseObject * moose_object)
40 : UserObjectInterface(moose_object),
41 _params(moose_object->parameters()),
42 _num_models(_params.get<std::vector<UserObjectName>>(
"plastic_models").size()),
44 _specialIC(_params.get<MooseEnum>(
"specialIC"))
47 for (
unsigned model = 0; model <
_num_models; ++model)
48 _f[model] = &getUserObjectByName<TensorMechanicsPlasticModel>(
49 _params.get<std::vector<UserObjectName>>(
"plastic_models")[model]);
51 for (
unsigned model = 0; model <
_num_models; ++model)
56 unsigned int surface = 0;
57 for (
unsigned model = 0; model <
_num_models; ++model)
58 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
67 for (
unsigned model = 0; model <
_num_models; ++model)
70 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
78 mooseError(
"Choosing specialIC=rock, you must have plasticity models of type 'TensileMulti "
79 "MohrCoulombMulti'\n");
80 if (!(
_f[0]->modelName().compare(
"TensileMulti") == 0 &&
81 _f[1]->modelName().compare(
"MohrCoulombMulti") == 0))
82 mooseError(
"Choosing specialIC=rock, you must have plasticity models of type 'TensileMulti "
83 "MohrCoulombMulti'\n");
88 mooseError(
"Choosing specialIC=joint, you must have plasticity models of type "
89 "'WeakPlaneTensile WeakPlaneShear'\n");
90 if (!(
_f[0]->modelName().compare(
"WeakPlaneTensile") == 0 &&
91 _f[1]->modelName().compare(
"WeakPlaneShear") == 0))
92 mooseError(
"Choosing specialIC=joint, you must have plasticity models of type "
93 "'WeakPlaneTensile WeakPlaneShear'\n");
99 const std::vector<Real> & intnl,
100 const std::vector<bool> & active,
101 std::vector<Real> & f)
103 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
104 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
107 std::vector<unsigned int> active_surfaces_of_model;
108 std::vector<unsigned int>::iterator active_surface;
109 std::vector<Real> model_f;
110 for (
unsigned model = 0; model <
_num_models; ++model)
113 if (active_surfaces_of_model.size() > 0)
115 _f[model]->yieldFunctionV(stress, intnl[model], model_f);
116 for (active_surface = active_surfaces_of_model.begin();
117 active_surface != active_surfaces_of_model.end();
119 f.push_back(model_f[*active_surface]);
127 const std::vector<Real> & intnl,
128 const std::vector<bool> & active,
129 std::vector<RankTwoTensor> & df_dstress)
131 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
132 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
134 df_dstress.resize(0);
135 std::vector<unsigned int> active_surfaces_of_model;
136 std::vector<unsigned int>::iterator active_surface;
137 std::vector<RankTwoTensor> model_df_dstress;
138 for (
unsigned model = 0; model <
_num_models; ++model)
141 if (active_surfaces_of_model.size() > 0)
143 _f[model]->dyieldFunction_dstressV(stress, intnl[model], model_df_dstress);
144 for (active_surface = active_surfaces_of_model.begin();
145 active_surface != active_surfaces_of_model.end();
147 df_dstress.push_back(model_df_dstress[*active_surface]);
154 const std::vector<Real> & intnl,
155 const std::vector<bool> & active,
156 std::vector<Real> & df_dintnl)
158 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
159 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
162 std::vector<unsigned int> active_surfaces_of_model;
163 std::vector<unsigned int>::iterator active_surface;
164 std::vector<Real> model_df_dintnl;
165 for (
unsigned model = 0; model <
_num_models; ++model)
168 if (active_surfaces_of_model.size() > 0)
170 _f[model]->dyieldFunction_dintnlV(stress, intnl[model], model_df_dintnl);
171 for (active_surface = active_surfaces_of_model.begin();
172 active_surface != active_surfaces_of_model.end();
174 df_dintnl.push_back(model_df_dintnl[*active_surface]);
181 const std::vector<Real> & intnl,
182 const std::vector<bool> & active,
183 std::vector<RankTwoTensor> & r)
185 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
186 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
189 std::vector<unsigned int> active_surfaces_of_model;
190 std::vector<unsigned int>::iterator active_surface;
191 std::vector<RankTwoTensor> model_r;
192 for (
unsigned model = 0; model <
_num_models; ++model)
195 if (active_surfaces_of_model.size() > 0)
197 _f[model]->flowPotentialV(stress, intnl[model], model_r);
198 for (active_surface = active_surfaces_of_model.begin();
199 active_surface != active_surfaces_of_model.end();
201 r.push_back(model_r[*active_surface]);
209 const std::vector<Real> & intnl,
210 const std::vector<bool> & active,
211 std::vector<RankFourTensor> & dr_dstress)
213 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
214 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
216 dr_dstress.resize(0);
217 std::vector<unsigned int> active_surfaces_of_model;
218 std::vector<unsigned int>::iterator active_surface;
219 std::vector<RankFourTensor> model_dr_dstress;
220 for (
unsigned model = 0; model <
_num_models; ++model)
223 if (active_surfaces_of_model.size() > 0)
225 _f[model]->dflowPotential_dstressV(stress, intnl[model], model_dr_dstress);
226 for (active_surface = active_surfaces_of_model.begin();
227 active_surface != active_surfaces_of_model.end();
229 dr_dstress.push_back(model_dr_dstress[*active_surface]);
236 const std::vector<Real> & intnl,
237 const std::vector<bool> & active,
238 std::vector<RankTwoTensor> & dr_dintnl)
240 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
241 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
244 std::vector<unsigned int> active_surfaces_of_model;
245 std::vector<unsigned int>::iterator active_surface;
246 std::vector<RankTwoTensor> model_dr_dintnl;
247 for (
unsigned model = 0; model <
_num_models; ++model)
250 if (active_surfaces_of_model.size() > 0)
252 _f[model]->dflowPotential_dintnlV(stress, intnl[model], model_dr_dintnl);
253 for (active_surface = active_surfaces_of_model.begin();
254 active_surface != active_surfaces_of_model.end();
256 dr_dintnl.push_back(model_dr_dintnl[*active_surface]);
263 const std::vector<Real> & intnl,
264 const std::vector<bool> & active,
265 std::vector<Real> & h)
267 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
268 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
271 std::vector<unsigned int> active_surfaces_of_model;
272 std::vector<unsigned int>::iterator active_surface;
273 std::vector<Real> model_h;
274 for (
unsigned model = 0; model <
_num_models; ++model)
277 if (active_surfaces_of_model.size() > 0)
279 _f[model]->hardPotentialV(stress, intnl[model], model_h);
280 for (active_surface = active_surfaces_of_model.begin();
281 active_surface != active_surfaces_of_model.end();
283 h.push_back(model_h[*active_surface]);
291 const std::vector<Real> & intnl,
292 const std::vector<bool> & active,
293 std::vector<RankTwoTensor> & dh_dstress)
295 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
296 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
298 dh_dstress.resize(0);
299 std::vector<unsigned int> active_surfaces_of_model;
300 std::vector<unsigned int>::iterator active_surface;
301 std::vector<RankTwoTensor> model_dh_dstress;
302 for (
unsigned model = 0; model <
_num_models; ++model)
305 if (active_surfaces_of_model.size() > 0)
307 _f[model]->dhardPotential_dstressV(stress, intnl[model], model_dh_dstress);
308 for (active_surface = active_surfaces_of_model.begin();
309 active_surface != active_surfaces_of_model.end();
311 dh_dstress.push_back(model_dh_dstress[*active_surface]);
318 const std::vector<Real> & intnl,
319 const std::vector<bool> & active,
320 std::vector<Real> & dh_dintnl)
322 mooseAssert(intnl.size() ==
_num_models,
"Incorrect size of internal parameters");
323 mooseAssert(active.size() ==
_num_surfaces,
"Incorrect size of active");
326 std::vector<unsigned int> active_surfaces_of_model;
327 std::vector<unsigned int>::iterator active_surface;
328 std::vector<Real> model_dh_dintnl;
329 for (
unsigned model = 0; model <
_num_models; ++model)
332 if (active_surfaces_of_model.size() > 0)
334 _f[model]->dhardPotential_dintnlV(stress, intnl[model], model_dh_dintnl);
335 for (active_surface = active_surfaces_of_model.begin();
336 active_surface != active_surfaces_of_model.end();
338 dh_dintnl.push_back(model_dh_dintnl[*active_surface]);
346 const std::vector<Real> & intnl,
348 std::vector<bool> & act)
351 "buildActiveConstraints called with f.size = " << f.size() <<
" while there are "
354 "buildActiveConstraints called with intnl.size = "
355 << intnl.size() <<
" while there are " <<
_num_models <<
" models");
365 for (
unsigned model = 0; model <
_num_models; ++model)
367 std::vector<Real> model_f(0);
368 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
369 model_f.push_back(f[ind++]);
370 std::vector<bool> model_act;
372 _f[model]->activeConstraints(
373 model_f, stress, intnl[model], Eijkl, model_act, returned_stress);
374 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
375 act.push_back(model_act[model_surface]);
383 const std::vector<Real> & intnl,
385 std::vector<bool> & act)
387 act.assign(2,
false);
390 std::vector<bool> active_tensile;
391 std::vector<bool> active_shear;
392 std::vector<Real> f_single;
395 f_single.assign(1, 0);
397 _f[0]->activeConstraints(f_single, stress, intnl[0], Eijkl, active_tensile, returned_stress);
398 _f[1]->yieldFunctionV(returned_stress, intnl[1], f_single);
399 if (f_single[0] <=
_f[1]->_f_tol)
401 act[0] = active_tensile[0];
406 f_single.assign(1, 0);
408 _f[1]->activeConstraints(f_single, stress, intnl[1], Eijkl, active_shear, returned_stress);
409 _f[0]->yieldFunctionV(returned_stress, intnl[0], f_single);
410 if (f_single[0] <=
_f[0]->_f_tol)
412 act[1] = active_shear[0];
417 act[0] = act[1] =
true;
424 const std::vector<Real> & intnl,
426 std::vector<bool> & act)
428 act.assign(9,
false);
431 std::vector<bool> active_tensile;
432 std::vector<bool> active_MC;
433 std::vector<Real> f_single;
436 f_single.assign(3, 0);
440 _f[0]->activeConstraints(f_single, stress, intnl[0], Eijkl, active_tensile, returned_stress);
441 _f[1]->yieldFunctionV(returned_stress, intnl[1], f_single);
442 if (f_single[0] <=
_f[1]->_f_tol && f_single[1] <=
_f[1]->_f_tol &&
443 f_single[2] <=
_f[1]->_f_tol && f_single[3] <=
_f[1]->_f_tol &&
444 f_single[4] <=
_f[1]->_f_tol && f_single[5] <=
_f[1]->_f_tol)
446 act[0] = active_tensile[0];
447 act[1] = active_tensile[1];
448 act[2] = active_tensile[2];
453 f_single.assign(6, 0);
460 _f[1]->activeConstraints(f_single, stress, intnl[1], Eijkl, active_MC, returned_stress);
461 _f[0]->yieldFunctionV(returned_stress, intnl[0], f_single);
462 if (f_single[0] <=
_f[0]->_f_tol && f_single[1] <=
_f[0]->_f_tol && f_single[2] <=
_f[0]->_f_tol)
464 act[3] = active_MC[0];
465 act[4] = active_MC[1];
466 act[5] = active_MC[2];
467 act[6] = active_MC[3];
468 act[7] = active_MC[4];
469 act[8] = active_MC[5];
477 if (active_tensile[0] ==
false && active_tensile[1] ==
true && active_tensile[2] ==
true &&
478 active_MC[0] ==
true && active_MC[1] ==
true && active_MC[2] ==
false &&
479 active_MC[3] ==
true && active_MC[4] ==
false && active_MC[5] ==
false)
481 act[1] = act[2] = act[6] =
true;
485 if (active_tensile[0] ==
false && active_tensile[1] ==
true && active_tensile[2] ==
true &&
486 active_MC[0] ==
false && active_MC[1] ==
true && active_MC[2] ==
false &&
487 active_MC[3] ==
true && active_MC[4] ==
false && active_MC[5] ==
true)
489 act[1] = act[2] = act[6] =
true;
494 if (active_tensile[0] ==
false && active_tensile[1] ==
true && active_tensile[2] ==
true &&
495 active_MC[0] ==
false && active_MC[1] ==
true && active_MC[2] ==
false &&
496 active_MC[3] ==
true && active_MC[4] ==
false && active_MC[5] ==
false)
498 act[1] = act[2] = act[4] = act[6] =
true;
501 if (active_tensile[0] ==
false && active_tensile[1] ==
true && active_tensile[2] ==
true &&
502 active_MC[0] ==
false && active_MC[1] ==
false && active_MC[2] ==
false &&
503 active_MC[3] ==
true && active_MC[4] ==
false && active_MC[5] ==
true)
505 act[1] = act[2] = act[4] = act[6] =
true;
510 if (active_tensile[0] ==
false && active_tensile[1] ==
true && active_tensile[2] ==
true &&
511 active_MC[0] ==
false && active_MC[1] ==
false && active_MC[2] ==
false &&
512 active_MC[3] ==
true && active_MC[4] ==
false && active_MC[5] ==
false)
514 act[1] = act[2] = act[6] =
true;
519 if (active_tensile[0] ==
false && active_tensile[1] ==
false && active_tensile[2] ==
true &&
520 active_MC[0] ==
true && active_MC[1] ==
true && active_MC[2] ==
false &&
521 active_MC[3] ==
true && active_MC[4] ==
false && active_MC[5] ==
false)
523 act[2] = act[6] =
true;
528 if (active_tensile[0] ==
false && active_tensile[1] ==
false && active_tensile[2] ==
true &&
529 active_MC[0] ==
false && active_MC[1] ==
true && active_MC[2] ==
false &&
530 active_MC[3] ==
true && active_MC[4] ==
false && active_MC[5] ==
true)
532 act[2] = act[6] =
true;
538 if (active_tensile[0] ==
false && active_tensile[1] ==
false && active_tensile[2] ==
true &&
539 active_MC[0] ==
false && active_MC[1] ==
false && active_MC[2] ==
false &&
540 active_MC[3] ==
true && active_MC[4] ==
false && active_MC[5] ==
false)
542 act[1] = act[2] = act[6] =
true;
548 act[3] = active_MC[0];
549 act[4] = active_MC[1];
550 act[5] = active_MC[2];
551 act[6] = active_MC[3];
552 act[7] = active_MC[4];
553 act[8] = active_MC[5];
600 const std::vector<Real> & intnl_old,
602 Real ep_plastic_tolerance,
604 std::vector<Real> & intnl,
605 std::vector<Real> & pm,
606 std::vector<Real> & cumulative_pm,
608 std::vector<Real> & yf,
609 unsigned & num_successful_plastic_returns,
610 unsigned & custom_model)
613 "returnMapAll: Incorrect size of internal parameters");
614 mooseAssert(intnl.size() ==
_num_models,
"returnMapAll: Incorrect size of internal parameters");
615 mooseAssert(pm.size() ==
_num_surfaces,
"returnMapAll: Incorrect size of pm");
617 num_successful_plastic_returns = 0;
624 std::vector<Real> model_f;
626 std::vector<Real> model_pm;
627 bool trial_stress_inadmissible;
628 bool successful_return =
true;
629 unsigned the_single_plastic_model = 0;
630 bool using_custom_return_map =
true;
637 for (
unsigned model = 0; model <
_num_models; ++model)
639 if (using_custom_return_map)
641 model_pm.assign(
_f[model]->numberSurfaces(), 0.0);
642 bool model_returned =
_f[model]->returnMap(trial_stress,
645 ep_plastic_tolerance,
651 trial_stress_inadmissible);
652 if (!trial_stress_inadmissible)
656 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces();
658 yf.push_back(model_f[model_surface]);
660 else if (trial_stress_inadmissible && !model_returned)
667 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces();
669 yf.push_back(model_f[model_surface]);
672 using_custom_return_map =
false;
673 successful_return =
false;
681 if (trial_stress_inadmissible)
682 num_successful_plastic_returns++;
683 the_single_plastic_model = model;
684 stress = returned_stress;
701 _f[model]->yieldFunctionV(trial_stress, intnl_old[model], model_f);
702 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
703 yf.push_back(model_f[model_surface]);
707 if (num_successful_plastic_returns == 0)
713 stress = trial_stress;
714 for (
unsigned model = 0; model <
_num_models; ++model)
715 intnl[model] = intnl_old[model];
716 return successful_return;
723 std::vector<Real> yf_at_returned_stress(0);
724 bool all_admissible =
true;
725 for (
unsigned model = 0; model <
_num_models; ++model)
727 if (model == the_single_plastic_model)
730 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
731 yf_at_returned_stress.push_back(model_f[model_surface]);
734 _f[model]->yieldFunctionV(stress, intnl_old[model], model_f);
735 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
737 if (model_f[model_surface] >
_f[model]->_f_tol)
739 all_admissible =
false;
740 yf_at_returned_stress.push_back(model_f[model_surface]);
752 stress = trial_stress;
753 for (
unsigned model = 0; model <
_num_models; ++model)
754 intnl[model] = intnl_old[model];
756 for (
unsigned model = the_single_plastic_model; model <
_num_models; ++model)
758 _f[model]->yieldFunctionV(trial_stress, intnl[model], model_f);
759 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
760 yf.push_back(model_f[model_surface]);
762 num_successful_plastic_returns = 0;
770 for (
unsigned surface = 0; surface < yf_at_returned_stress.size(); ++surface)
771 yf.push_back(yf_at_returned_stress[surface]);
772 delta_dp = model_delta_dp;
773 for (
unsigned model_surface = 0; model_surface <
_f[the_single_plastic_model]->numberSurfaces();
777 model_pm[model_surface];
780 custom_model = the_single_plastic_model;
793 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
801 const std::vector<bool> & active,
802 std::vector<unsigned int> & active_surfaces)
804 active_surfaces.resize(0);
805 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
813 const std::vector<bool> & active,
814 std::vector<unsigned int> & active_surfaces_of_model)
816 active_surfaces_of_model.resize(0);
817 for (
unsigned model_surface = 0; model_surface <
_f[model]->numberSurfaces(); ++model_surface)
819 active_surfaces_of_model.push_back(model_surface);