12 #include "MatrixTools.h"
14 #include "MooseException.h"
15 #include "RotationMatrix.h"
17 #include "libmesh/utility.h"
28 params.addClassDescription(
"Base class for multi-surface finite-strain plasticity");
29 params.addRangeCheckedParam<
unsigned int>(
"max_NR_iterations",
31 "max_NR_iterations>0",
32 "Maximum number of Newton-Raphson iterations allowed");
33 params.addRequiredParam<Real>(
"ep_plastic_tolerance",
34 "The Newton-Raphson process is only deemed "
35 "converged if the plastic strain increment "
36 "constraints have L2 norm less than this.");
37 params.addRangeCheckedParam<Real>(
40 "min_stepsize>0 & min_stepsize<=1",
41 "If ordinary Newton-Raphson + line-search fails, then the applied strain increment is "
42 "subdivided, and the return-map is tried again. This parameter is the minimum fraction of "
43 "applied strain increment that may be applied before the algorithm gives up entirely");
44 params.addRangeCheckedParam<Real>(
"max_stepsize_for_dumb",
46 "max_stepsize_for_dumb>0 & max_stepsize_for_dumb<=1",
47 "If your deactivation_scheme is 'something_to_dumb', then "
48 "'dumb' will only be used if the stepsize falls below this "
49 "value. This parameter is useful because the 'dumb' scheme is "
50 "computationally expensive");
51 MooseEnum deactivation_scheme(
"optimized safe dumb optimized_to_safe safe_to_dumb "
52 "optimized_to_safe_to_dumb optimized_to_dumb",
54 params.addParam<MooseEnum>(
55 "deactivation_scheme",
57 "Scheme by which constraints are deactivated. (NOTE: This is irrelevant if there is only "
58 "one yield surface.) safe: return to the yield surface and then deactivate constraints with "
59 "negative plasticity multipliers. optimized: deactivate a constraint as soon as its "
60 "plasticity multiplier becomes negative. dumb: iteratively try all combinations of active "
61 "constraints until the solution is found. You may specify fall-back options. Eg "
62 "optimized_to_safe: first use 'optimized', and if that fails, try the return with 'safe'.");
63 params.addParam<RealVectorValue>(
64 "transverse_direction",
65 "If this parameter is provided, before the return-map algorithm is "
66 "called a rotation is performed so that the 'z' axis in the new "
67 "frame lies along the transverse_direction in the original frame. "
68 "After returning, the inverse rotation is performed. The "
69 "transverse_direction will itself rotate with large strains. This "
70 "is so that transversely-isotropic plasticity models may be easily "
71 "defined in the frame where the isotropy holds in the x-y plane.");
72 params.addParam<
bool>(
"ignore_failures",
74 "The return-map algorithm will return with the best admissible "
75 "stresses and internal parameters that it can, even if they don't "
76 "fully correspond to the applied strain increment. To speed "
77 "computations, this flag can be set to true, the max_NR_iterations "
78 "set small, and the min_stepsize large.");
79 MooseEnum tangent_operator(
"elastic linear nonlinear",
"nonlinear");
80 params.addParam<MooseEnum>(
"tangent_operator",
82 "Type of tangent operator to return. 'elastic': return the "
83 "elasticity tensor. 'linear': return the consistent tangent operator "
84 "that is correct for plasticity with yield functions linear in "
85 "stress. 'nonlinear': return the full, general consistent tangent "
86 "operator. The calculations assume the hardening potentials are "
87 "independent of stress and hardening parameters.");
88 params.addParam<
bool>(
"perform_finite_strain_rotations",
90 "Tensors are correctly rotated in "
91 "finite-strain simulations. For "
92 "optimal performance you can set "
93 "this to 'false' if you are only "
94 "ever using small strains");
95 params.addClassDescription(
"Material for multi-surface finite-strain plasticity");
102 _max_iter(getParam<unsigned int>(
"max_NR_iterations")),
103 _min_stepsize(getParam<Real>(
"min_stepsize")),
104 _max_stepsize_for_dumb(getParam<Real>(
"max_stepsize_for_dumb")),
105 _ignore_failures(getParam<bool>(
"ignore_failures")),
109 _epp_tol(getParam<Real>(
"ep_plastic_tolerance")),
117 _n_supplied(parameters.isParamValid(
"transverse_direction")),
118 _n_input(_n_supplied ? getParam<RealVectorValue>(
"transverse_direction") : RealVectorValue()),
119 _rot(RealTensorValue()),
121 _perform_finite_strain_rotations(getParam<bool>(
"perform_finite_strain_rotations")),
123 _elasticity_tensor_name(_base_name +
"elasticity_tensor"),
124 _elasticity_tensor(getMaterialPropertyByName<
RankFourTensor>(_elasticity_tensor_name)),
125 _plastic_strain(declareProperty<
RankTwoTensor>(
"plastic_strain")),
126 _plastic_strain_old(getMaterialPropertyOld<
RankTwoTensor>(
"plastic_strain")),
127 _intnl(declareProperty<std::vector<Real>>(
"plastic_internal_parameter")),
128 _intnl_old(getMaterialPropertyOld<std::vector<Real>>(
"plastic_internal_parameter")),
129 _yf(declareProperty<std::vector<Real>>(
"plastic_yield_function")),
130 _iter(declareProperty<Real>(
"plastic_NR_iterations")),
132 _linesearch_needed(declareProperty<Real>(
"plastic_linesearch_needed")),
136 _ld_encountered(declareProperty<Real>(
137 "plastic_linear_dependence_encountered")),
139 _constraints_added(declareProperty<Real>(
"plastic_constraints_added")),
143 _n(declareProperty<RealVectorValue>(
"plastic_transverse_direction")),
144 _n_old(getMaterialPropertyOld<RealVectorValue>(
"plastic_transverse_direction")),
146 _strain_increment(getMaterialPropertyByName<
RankTwoTensor>(_base_name +
"strain_increment")),
147 _total_strain_old(getMaterialPropertyOldByName<
RankTwoTensor>(_base_name +
"total_strain")),
149 getMaterialPropertyByName<
RankTwoTensor>(_base_name +
"rotation_increment")),
151 _stress_old(getMaterialPropertyOld<
RankTwoTensor>(_base_name +
"stress")),
152 _elastic_strain_old(getMaterialPropertyOld<
RankTwoTensor>(_base_name +
"elastic_strain")),
157 hasMaterialProperty<
RankFourTensor>(
"elastic_flexural_rigidity_tensor")),
158 _curvature(_cosserat ? &getMaterialPropertyByName<
RankTwoTensor>(
"curvature") : nullptr),
159 _elastic_flexural_rigidity_tensor(
160 _cosserat ? &getMaterialPropertyByName<
RankFourTensor>(
"elastic_flexural_rigidity_tensor")
162 _couple_stress(_cosserat ? &declareProperty<
RankTwoTensor>(
"couple_stress") : nullptr),
163 _couple_stress_old(_cosserat ? &getMaterialPropertyOld<
RankTwoTensor>(
"couple_stress")
165 _Jacobian_mult_couple(_cosserat ? &declareProperty<
RankFourTensor>(
"couple_Jacobian_mult")
174 mooseError(
"ComputeMultiPlasticityStress: ep_plastic_tolerance must be positive");
181 "ComputeMultiPlasticityStress: transverse_direction vector must not have zero length");
211 (*_couple_stress)[_qp].zero();
216 mooseError(
"Finite-differencing completed. Exiting with no error");
237 mooseError(
"Finite-differencing completed. Exiting with no error");
242 unsigned int number_iterations;
243 bool linesearch_needed =
false;
244 bool ld_encountered =
false;
245 bool constraints_added =
false;
285 (*_couple_stress)[_qp] = (*_elastic_flexural_rigidity_tensor)[_qp] *
_my_curvature;
291 _iter[_qp] = 1.0 * number_iterations;
316 return tens.rotated(
_rot);
325 _rot = RotationMatrix::rotVecToZ(
_n[_qp]);
355 (*_Jacobian_mult_couple)[_qp].rotate(
_rot);
357 (*_couple_stress)[_qp].rotate(
_rot);
361 for (
unsigned int i = 0; i < LIBMESH_DIM; ++i)
364 for (
unsigned int j = 0; j < LIBMESH_DIM; ++j)
373 const std::vector<Real> & intnl_old,
374 std::vector<Real> & intnl,
375 std::vector<Real> & pm,
376 std::vector<Real> & cumulative_pm,
381 std::vector<Real> & yf,
382 unsigned int & iterations,
389 unsigned num_plastic_returns;
394 unsigned custom_model = 0;
395 bool successful_return =
returnMapAll(stress_old + E_ijkl * strain_increment,
410 if (num_plastic_returns == 0)
417 plastic_strain = plastic_strain_old;
418 if (successful_return && final_step)
421 consistent_tangent_operator = E_ijkl;
424 consistent_tangent_operator =
427 return successful_return;
429 else if (num_plastic_returns == 1 && successful_return)
434 plastic_strain = plastic_strain_old + delta_dp;
440 consistent_tangent_operator = E_ijkl;
443 std::vector<Real> custom_model_pm;
444 for (
unsigned surface = 0; surface <
_f[custom_model]->numberSurfaces(); ++surface)
446 consistent_tangent_operator =
447 _f[custom_model]->consistentTangentOperator(stress_old + E_ijkl * strain_increment,
448 intnl_old[custom_model],
458 consistent_tangent_operator =
464 mooseError(
"ComputeMultiPlasticityStress::quickStep should not get here!");
470 const std::vector<Real> & intnl_old,
471 std::vector<Real> & intnl,
476 std::vector<Real> & yf,
477 unsigned int & iterations,
478 bool & linesearch_needed,
479 bool & ld_encountered,
480 bool & constraints_added,
498 bool return_successful =
false;
501 unsigned int iter = 0;
503 Real step_size = 1.0;
504 Real time_simulated = 0.0;
511 for (
unsigned model = 0; model <
_num_models; ++model)
512 intnl_good[model] = intnl_old[model];
523 unsigned int num_consecutive_successes = 0;
527 return_successful =
returnMap(stress_good,
541 time_simulated + step_size >= 1,
542 consistent_tangent_operator,
546 if (return_successful)
548 num_consecutive_successes += 1;
549 time_simulated += step_size;
551 if (time_simulated < 1.0)
554 stress_good = stress;
555 plastic_strain_good = plastic_strain;
556 for (
unsigned model = 0; model <
_num_models; ++model)
557 intnl_good[model] = intnl[model];
558 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
559 yf_good[surface] = yf[surface];
560 if (num_consecutive_successes >= 2)
563 step_size = std::min(step_size, 1.0 - time_simulated);
568 num_consecutive_successes = 0;
569 stress = stress_good;
570 plastic_strain = plastic_strain_good;
571 for (
unsigned model = 0; model <
_num_models; ++model)
572 intnl[model] = intnl_good[model];
574 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
575 yf[surface] = yf_good[surface];
576 dep = step_size * this_strain_increment;
580 if (!return_successful)
584 stress = stress_good;
585 plastic_strain = plastic_strain_good;
586 for (
unsigned model = 0; model <
_num_models; ++model)
587 intnl[model] = intnl_good[model];
588 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
589 yf[surface] = yf_good[surface];
593 Moose::out <<
"After reducing the stepsize to " << step_size
594 <<
" with original strain increment with L2norm " << this_strain_increment.L2norm()
595 <<
" the returnMap algorithm failed\n";
602 for (
unsigned model = 0; model <
_num_models; ++model)
607 mooseError(
"Exiting\n");
611 return return_successful;
617 const std::vector<Real> & intnl_old,
618 std::vector<Real> & intnl,
623 std::vector<Real> & f,
625 bool can_revert_to_dumb,
626 bool & linesearch_needed,
627 bool & ld_encountered,
628 bool & constraints_added,
631 std::vector<Real> & cumulative_pm)
637 std::vector<Real> pm;
640 bool successful_return =
quickStep(stress_old,
652 consistent_tangent_operator,
656 if (successful_return)
657 return successful_return;
687 std::vector<Real> ic;
698 std::vector<bool> act;
712 bool single_step_success =
true;
718 std::vector<bool> initial_act;
726 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
727 initial_act[surface] = act[surface];
734 int dumb_iteration = 0;
735 std::vector<unsigned int> dumb_order;
752 std::set<unsigned int> actives_tried;
759 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
761 nr_res2 += 0.5 * Utility::pow<2>(f[surface] /
_f[
modelNumber(surface)]->_f_tol);
763 successful_return =
false;
765 bool still_finding_solution =
true;
766 while (still_finding_solution)
768 single_step_success =
true;
769 unsigned int local_iter = 0;
772 while (nr_res2 > 0.5 && local_iter++ <
_max_iter && single_step_success)
788 bool nr_good = (nr_res2 <= 0.5 && local_iter <=
_max_iter && single_step_success);
802 bool change_scheme =
false;
803 bool increment_dumb =
false;
805 if (!change_scheme && deact_scheme ==
dumb)
808 still_finding_solution = (change_scheme || increment_dumb);
823 if (!still_finding_solution)
826 successful_return =
false;
831 bool kt_good =
false;
838 if (deact_scheme !=
dumb)
843 still_finding_solution =
845 if (!still_finding_solution)
851 still_finding_solution =
true;
865 bool increment_dumb =
false;
867 still_finding_solution = increment_dumb;
872 if (!still_finding_solution)
876 successful_return =
false;
882 bool admissible =
false;
883 if (nr_good && kt_good)
886 std::vector<Real> all_f;
902 constraints_added =
true;
905 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
907 (!act[surface] && (all_f[surface] >
_f[
modelNumber(surface)]->_f_tol)))
908 act_plus[surface] =
true;
912 constraints_added =
true;
913 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
914 act[surface] = act_plus[surface];
917 add_constraints =
false;
920 bool change_scheme =
false;
921 bool increment_dumb =
false;
923 if (!add_constraints)
926 if (!add_constraints && !change_scheme && deact_scheme ==
dumb)
929 still_finding_solution = (add_constraints || change_scheme || increment_dumb);
944 if (!still_finding_solution)
947 successful_return =
false;
953 successful_return = (nr_good && admissible && kt_good);
954 if (successful_return)
957 if (still_finding_solution)
959 stress = initial_stress;
961 for (
unsigned model = 0; model <
_num_models; ++model)
962 intnl[model] = intnl_old[model];
968 successful_return =
false;
980 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
984 nr_res2 += 0.5 * Utility::pow<2>(f[ind] /
_f[
modelNumber(surface)]->_f_tol);
991 if (successful_return)
993 plastic_strain += delta_dp;
995 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
996 cumulative_pm[surface] += pm[surface];
999 consistent_tangent_operator =
1011 return successful_return;
1016 const std::vector<Real> & all_f)
1018 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1019 if (!act[surface] && (all_f[surface] >
_f[
modelNumber(surface)]->_f_tol))
1026 bool can_revert_to_dumb)
1051 bool can_revert_to_dumb,
1053 const std::vector<Real> & intnl_old,
1055 std::vector<bool> & act,
1056 int & dumb_iteration,
1057 std::vector<unsigned int> & dumb_order)
1059 if (current_deactivation_scheme ==
optimized &&
1063 current_deactivation_scheme =
safe;
1064 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1065 act[surface] = initial_act[surface];
1067 else if ((current_deactivation_scheme ==
safe &&
1070 can_revert_to_dumb) ||
1072 can_revert_to_dumb))
1074 current_deactivation_scheme =
dumb;
1084 const std::vector<Real> & intnl_old,
1085 std::vector<Real> & intnl,
1086 std::vector<Real> & pm,
1089 std::vector<Real> & f,
1091 std::vector<Real> & ic,
1092 std::vector<bool> & active,
1094 bool & linesearch_needed,
1095 bool & ld_encountered)
1097 bool successful_step;
1099 Real nr_res2_before_step = nr_res2;
1101 std::vector<Real> intnl_before_step;
1102 std::vector<Real> pm_before_step;
1108 stress_before_step = stress;
1110 for (
unsigned model = 0; model <
_num_models; ++model)
1111 intnl_before_step[model] = intnl[model];
1113 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1114 pm_before_step[surface] = pm[surface];
1115 delta_dp_before_step = delta_dp;
1122 std::vector<Real> dpm;
1129 std::vector<bool> deact_ld;
1137 bool constraints_changing =
true;
1138 bool reinstated_actives;
1139 while (constraints_changing)
1142 nrStep(stress, intnl_old, intnl, pm, E_inv, delta_dp, dstress, dpm, dintnl, active, deact_ld);
1144 for (
unsigned surface = 0; surface < deact_ld.size(); ++surface)
1145 if (deact_ld[surface])
1147 ld_encountered =
true;
1170 if (!successful_step)
1172 return successful_step;
1175 constraints_changing =
false;
1178 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1179 if (active[surface] && pm[surface] < 0.0)
1180 constraints_changing =
true;
1183 if (constraints_changing)
1185 stress = stress_before_step;
1186 delta_dp = delta_dp_before_step;
1187 nr_res2 = nr_res2_before_step;
1188 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1190 if (active[surface] && pm[surface] < 0.0)
1193 active[surface] =
false;
1194 pm_before_step[surface] = 0.0;
1199 pm[surface] = pm_before_step[surface];
1204 successful_step =
false;
1205 return successful_step;
1215 if (reinstated_actives)
1217 bool completely_converged =
true;
1218 if (successful_step && nr_res2 < 0.5)
1232 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1233 if (active[surface])
1235 completely_converged =
false;
1238 completely_converged =
false;
1240 if (!completely_converged)
1241 nr_res2 =
residual2(pm, f, epp, ic, active, deact_ld);
1244 return successful_step;
1249 std::vector<bool> & deactivated_due_to_ld)
1251 bool reinstated_actives =
false;
1252 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1253 if (deactivated_due_to_ld[surface])
1254 reinstated_actives =
true;
1257 return reinstated_actives;
1263 unsigned num_active = 0;
1264 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1265 if (active[surface])
1272 const std::vector<Real> & intnl,
1273 std::vector<Real> & all_f)
1275 std::vector<bool> act;
1280 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1289 const std::vector<Real> & pm,
1290 const std::vector<bool> & active)
1293 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1295 if (active[surface])
1298 if (pm[surface] != 0)
1301 else if (pm[surface] != 0)
1302 mooseError(
"Crash due to plastic multiplier not being zero. This occurred because of poor "
1305 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1306 if (pm[surface] < 0)
1314 const std::vector<Real> & pm,
1315 std::vector<bool> & active)
1317 bool turned_off =
false;
1321 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1323 if (active[surface])
1326 if (pm[surface] != 0)
1329 active[surface] =
false;
1332 else if (pm[surface] != 0)
1333 mooseError(
"Crash due to plastic multiplier not being zero. This occurred because of poor "
1340 int surface_to_turn_off = -1;
1342 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1343 if (pm[surface] < min_pm)
1345 min_pm = pm[surface];
1346 surface_to_turn_off = surface;
1348 if (surface_to_turn_off >= 0)
1349 active[surface_to_turn_off] =
false;
1355 const std::vector<Real> & f,
1357 const std::vector<Real> & ic,
1358 const std::vector<bool> & active,
1359 const std::vector<bool> & deactivated_due_to_ld)
1364 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1365 if (active[surface])
1367 if (!deactivated_due_to_ld[surface])
1369 if (!(pm[surface] == 0 && f[ind] <= 0))
1370 nr_res2 += 0.5 * Utility::pow<2>(f[ind] /
_f[
modelNumber(surface)]->_f_tol);
1372 else if (deactivated_due_to_ld[surface] && f[ind] > 0)
1373 nr_res2 += 0.5 * Utility::pow<2>(f[ind] /
_f[
modelNumber(surface)]->_f_tol);
1377 nr_res2 += 0.5 * Utility::pow<2>(epp.L2norm() /
_epp_tol);
1380 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1381 active_not_deact[surface] = (active[surface] && !deactivated_due_to_ld[surface]);
1383 for (
unsigned model = 0; model <
_num_models; ++model)
1385 nr_res2 += 0.5 * Utility::pow<2>(ic[ind++] /
_f[model]->_ic_tol);
1393 const std::vector<Real> & intnl_old,
1394 std::vector<Real> & intnl,
1395 std::vector<Real> & pm,
1399 const std::vector<Real> & dpm,
1400 const std::vector<Real> & dintnl,
1401 std::vector<Real> & f,
1403 std::vector<Real> & ic,
1404 const std::vector<bool> & active,
1405 const std::vector<bool> & deactivated_due_to_ld,
1406 bool & linesearch_needed)
1419 Real slope = -2 * nr_res2;
1426 std::vector<Real> ls_pm;
1427 ls_pm.resize(pm.size());
1433 std::vector<Real> ls_intnl;
1434 ls_intnl.resize(intnl.size());
1440 std::vector<RankTwoTensor> r;
1445 for (
unsigned alpha = 0; alpha < pm.size(); ++alpha)
1446 ls_pm[alpha] = pm[alpha] + dpm[alpha] * lam;
1447 ls_delta_dp = delta_dp - E_inv * dstress * lam;
1448 for (
unsigned a = 0; a < intnl.size(); ++a)
1449 ls_intnl[a] = intnl[a] + dintnl[a] * lam;
1450 ls_stress = stress + dstress * lam;
1453 calculateConstraints(ls_stress, intnl_old, ls_intnl, ls_pm, ls_delta_dp, f, r, epp, ic, active);
1456 nr_res2 =
residual2(ls_pm, f, epp, ic, active, deactivated_due_to_ld);
1458 if (nr_res2 < f0 + 1E-4 * lam * slope)
1460 else if (lam < lam_min)
1464 for (
unsigned alpha = 0; alpha < pm.size(); ++alpha)
1465 ls_pm[alpha] = pm[alpha];
1466 ls_delta_dp = delta_dp;
1467 for (
unsigned a = 0; a < intnl.size(); ++a)
1468 ls_intnl[a] = intnl[a];
1471 ls_stress, intnl_old, ls_intnl, ls_pm, ls_delta_dp, f, r, epp, ic, active);
1472 nr_res2 =
residual2(ls_pm, f, epp, ic, active, deactivated_due_to_ld);
1475 else if (lam == 1.0)
1478 tmp_lam = -slope / 2.0 / (nr_res2 - f0 - slope);
1483 Real rhs1 = nr_res2 - f0 - lam * slope;
1484 Real rhs2 = f2 - f0 - lam2 * slope;
1485 Real a = (rhs1 / Utility::pow<2>(lam) - rhs2 / Utility::pow<2>(lam2)) / (lam - lam2);
1487 (-lam2 * rhs1 / Utility::pow<2>(lam) + lam * rhs2 / Utility::pow<2>(lam2)) / (lam - lam2);
1489 tmp_lam = -slope / (2.0 * b);
1492 Real disc = Utility::pow<2>(b) - 3 * a * slope;
1494 tmp_lam = 0.5 * lam;
1496 tmp_lam = (-b + std::sqrt(disc)) / (3.0 * a);
1498 tmp_lam = -slope / (b + std::sqrt(disc));
1500 if (tmp_lam > 0.5 * lam)
1501 tmp_lam = 0.5 * lam;
1505 lam = std::max(tmp_lam, 0.1 * lam);
1509 linesearch_needed =
true;
1513 for (
unsigned alpha = 0; alpha < pm.size(); ++alpha)
1514 pm[alpha] = ls_pm[alpha];
1515 delta_dp = ls_delta_dp;
1516 for (
unsigned a = 0; a < intnl.size(); ++a)
1517 intnl[a] = ls_intnl[a];
1525 const std::vector<Real> & intnl,
1526 std::vector<unsigned int> & dumb_order)
1528 if (dumb_order.size() != 0)
1531 std::vector<bool> act;
1534 std::vector<Real> f;
1536 std::vector<RankTwoTensor> df_dstress;
1539 typedef std::pair<Real, unsigned> pair_for_sorting;
1541 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1543 dist[surface].first = f[surface] / df_dstress[surface].L2norm();
1544 dist[surface].second = surface;
1546 std::sort(dist.begin(), dist.end());
1556 const std::vector<unsigned int> & dumb_order,
1557 std::vector<bool> & act)
1559 dumb_iteration += 1;
1560 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1561 act[dumb_order[surface]] =
1577 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1579 num += (1 << surface);
1586 const std::vector<Real> & intnl,
1588 const std::vector<Real> & pm_this_step,
1589 const std::vector<Real> & cumulative_pm)
1599 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1600 act_at_some_step[surface] = (cumulative_pm[surface] > 0);
1607 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1608 act_vary[surface] = (pm_this_step[surface] > 0);
1610 std::vector<RankTwoTensor> df_dstress;
1612 std::vector<Real> df_dintnl;
1614 std::vector<RankTwoTensor> r;
1616 std::vector<RankFourTensor> dr_dstress_at_some_step;
1618 std::vector<RankTwoTensor> dr_dintnl_at_some_step;
1620 std::vector<Real> h;
1629 std::vector<RankTwoTensor> r_minus_stuff;
1631 for (
unsigned surface1 = 0; surface1 <
_num_surfaces; ++surface1)
1632 if (act_vary[surface1])
1634 r_minus_stuff.push_back(r[ind1]);
1636 for (
unsigned surface2 = 0; surface2 <
_num_surfaces; ++surface2)
1637 if (act_at_some_step[surface2])
1641 r_minus_stuff.back() -=
1642 cumulative_pm[surface2] * dr_dintnl_at_some_step[ind2] * h[ind1];
1649 unsigned int num_currently_active = 0;
1650 for (
unsigned surface = 0; surface <
_num_surfaces; ++surface)
1651 if (act_vary[surface])
1652 num_currently_active += 1;
1660 std::vector<PetscScalar> zzz;
1661 zzz.assign(num_currently_active * num_currently_active, 0.0);
1665 for (
unsigned surface1 = 0; surface1 <
_num_surfaces; ++surface1)
1666 if (act_vary[surface1])
1669 for (
unsigned surface2 = 0; surface2 <
_num_surfaces; ++surface2)
1670 if (act_vary[surface2])
1672 r2 = df_dstress[ind1] * (E_ijkl * r_minus_stuff[ind2]);
1673 zzz[ind1 * num_currently_active + ind2] += r2(0, 0) + r2(1, 1) + r2(2, 2);
1675 zzz[ind1 * num_currently_active + ind2] += df_dintnl[ind1] * h[ind2];
1681 if (num_currently_active > 0)
1686 MatrixTools::inverse(zzz, num_currently_active);
1688 catch (
const MooseException & e)
1697 for (
unsigned surface1 = 0; surface1 <
_num_surfaces; ++surface1)
1698 if (act_vary[surface1])
1702 for (
unsigned surface2 = 0; surface2 <
_num_surfaces; ++surface2)
1703 if (act_vary[surface2])
1706 for (
unsigned i = 0; i < 3; i++)
1707 for (
unsigned j = 0; j < 3; j++)
1708 for (
unsigned k = 0; k < 3; k++)
1709 for (
unsigned l = 0; l < 3; l++)
1710 strain_coeff(i, j, k, l) -=
1711 part1(i, j) * part2(k, l) * zzz[ind1 * num_currently_active + ind2];
1718 return strain_coeff;
1720 RankFourTensor stress_coeff(RankFourTensor::initIdentitySymmetricFour);
1724 for (
unsigned surface1 = 0; surface1 <
_num_surfaces; ++surface1)
1725 if (act_at_some_step[surface1])
1727 part3 += cumulative_pm[surface1] * E_ijkl * dr_dstress_at_some_step[ind1];
1731 stress_coeff += part3;
1733 part3 = part3.transposeMajor();
1737 for (
unsigned surface1 = 0; surface1 <
_num_surfaces; ++surface1)
1738 if (act_vary[surface1])
1742 for (
unsigned surface2 = 0; surface2 <
_num_surfaces; ++surface2)
1743 if (act_vary[surface2])
1746 for (
unsigned i = 0; i < 3; i++)
1747 for (
unsigned j = 0; j < 3; j++)
1748 for (
unsigned k = 0; k < 3; k++)
1749 for (
unsigned l = 0; l < 3; l++)
1750 stress_coeff(i, j, k, l) -=
1751 part1(i, j) * part2(k, l) * zzz[ind1 * num_currently_active + ind2];
1768 s_inv = stress_coeff.invSymm();
1770 catch (
const MooseException & e)
1772 return strain_coeff;
1776 return s_inv * strain_coeff;