www.mooseframework.org
ComputeMultiPlasticityStress.C
Go to the documentation of this file.
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 
12 #include "MatrixTools.h"
13 
14 #include "MooseException.h"
15 #include "RotationMatrix.h" // for rotVecToZ
16 
17 #include "libmesh/utility.h"
18 
20 
22 
23 InputParameters
25 {
26  InputParameters params = ComputeStressBase::validParams();
28  params.addClassDescription("Base class for multi-surface finite-strain plasticity");
29  params.addRangeCheckedParam<unsigned int>("max_NR_iterations",
30  20,
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>(
38  "min_stepsize",
39  0.01,
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",
45  0.01,
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",
53  "optimized");
54  params.addParam<MooseEnum>(
55  "deactivation_scheme",
56  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",
73  false,
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",
81  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",
89  true,
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");
96  return params;
97 }
98 
100  : ComputeStressBase(parameters),
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")),
106 
107  _tangent_operator_type((TangentOperatorEnum)(int)getParam<MooseEnum>("tangent_operator")),
108 
109  _epp_tol(getParam<Real>("ep_plastic_tolerance")),
110 
111  _dummy_pm(0),
112 
113  _cumulative_pm(0),
114 
115  _deactivation_scheme((DeactivationSchemeEnum)(int)getParam<MooseEnum>("deactivation_scheme")),
116 
117  _n_supplied(parameters.isParamValid("transverse_direction")),
118  _n_input(_n_supplied ? getParam<RealVectorValue>("transverse_direction") : RealVectorValue()),
119  _rot(RealTensorValue()),
120 
121  _perform_finite_strain_rotations(getParam<bool>("perform_finite_strain_rotations")),
122 
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")), // this is really an unsigned int, but
131  // for visualisation i convert it to Real
132  _linesearch_needed(declareProperty<Real>("plastic_linesearch_needed")), // this is really a
133  // boolean, but for
134  // visualisation i
135  // convert it to Real
136  _ld_encountered(declareProperty<Real>(
137  "plastic_linear_dependence_encountered")), // this is really a boolean, but for
138  // visualisation i convert it to Real
139  _constraints_added(declareProperty<Real>("plastic_constraints_added")), // this is really a
140  // boolean, but for
141  // visualisation i
142  // convert it to Real
143  _n(declareProperty<RealVectorValue>("plastic_transverse_direction")),
144  _n_old(getMaterialPropertyOld<RealVectorValue>("plastic_transverse_direction")),
145 
146  _strain_increment(getMaterialPropertyByName<RankTwoTensor>(_base_name + "strain_increment")),
147  _total_strain_old(getMaterialPropertyOldByName<RankTwoTensor>(_base_name + "total_strain")),
148  _rotation_increment(
149  getMaterialPropertyByName<RankTwoTensor>(_base_name + "rotation_increment")),
150 
151  _stress_old(getMaterialPropertyOld<RankTwoTensor>(_base_name + "stress")),
152  _elastic_strain_old(getMaterialPropertyOld<RankTwoTensor>(_base_name + "elastic_strain")),
153 
154  // TODO: This design does NOT work. It makes these materials construction order dependent and it
155  // disregards block restrictions.
156  _cosserat(hasMaterialProperty<RankTwoTensor>("curvature") &&
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")
161  : nullptr),
162  _couple_stress(_cosserat ? &declareProperty<RankTwoTensor>("couple_stress") : nullptr),
163  _couple_stress_old(_cosserat ? &getMaterialPropertyOld<RankTwoTensor>("couple_stress")
164  : nullptr),
165  _Jacobian_mult_couple(_cosserat ? &declareProperty<RankFourTensor>("couple_Jacobian_mult")
166  : nullptr),
167 
168  _my_elasticity_tensor(RankFourTensor()),
169  _my_strain_increment(RankTwoTensor()),
170  _my_flexural_rigidity_tensor(RankFourTensor()),
171  _my_curvature(RankTwoTensor())
172 {
173  if (_epp_tol <= 0)
174  mooseError("ComputeMultiPlasticityStress: ep_plastic_tolerance must be positive");
175 
176  if (_n_supplied)
177  {
178  // normalise the inputted transverse_direction
179  if (_n_input.norm() == 0)
180  mooseError(
181  "ComputeMultiPlasticityStress: transverse_direction vector must not have zero length");
182  else
183  _n_input /= _n_input.norm();
184  }
185 
186  if (_num_surfaces == 1)
188 }
189 
190 void
192 {
194 
195  _plastic_strain[_qp].zero();
196 
197  _intnl[_qp].assign(_num_models, 0);
198 
199  _yf[_qp].assign(_num_surfaces, 0);
200 
201  _dummy_pm.assign(_num_surfaces, 0);
202 
203  _iter[_qp] = 0.0; // this is really an unsigned int, but for visualisation i convert it to Real
204  _linesearch_needed[_qp] = 0;
205  _ld_encountered[_qp] = 0;
206  _constraints_added[_qp] = 0;
207 
208  _n[_qp] = _n_input;
209 
210  if (_cosserat)
211  (*_couple_stress)[_qp].zero();
212 
213  if (_fspb_debug == "jacobian")
214  {
216  mooseError("Finite-differencing completed. Exiting with no error");
217  }
218 }
219 
220 void
222 {
223  // the following "_my" variables can get rotated by preReturnMap and postReturnMap
226  if (_cosserat)
227  {
228  _my_flexural_rigidity_tensor = (*_elastic_flexural_rigidity_tensor)[_qp];
229  _my_curvature = (*_curvature)[_qp];
230  }
231 
232  if (_fspb_debug == "jacobian_and_linear_system")
233  {
234  // cannot do this at initQpStatefulProperties level since E_ijkl is not defined
235  checkJacobian(_elasticity_tensor[_qp].invSymm(), _intnl_old[_qp]);
236  checkSolution(_elasticity_tensor[_qp].invSymm());
237  mooseError("Finite-differencing completed. Exiting with no error");
238  }
239 
240  preReturnMap(); // do rotations to new frame if necessary
241 
242  unsigned int number_iterations;
243  bool linesearch_needed = false;
244  bool ld_encountered = false;
245  bool constraints_added = false;
246 
247  _cumulative_pm.assign(_num_surfaces, 0);
248  // try a "quick" return first - this can be purely elastic, or a customised plastic return defined
249  // by a TensorMechanicsPlasticXXXX UserObject
250  const bool found_solution = quickStep(rot(_stress_old[_qp]),
251  _stress[_qp],
252  _intnl_old[_qp],
253  _intnl[_qp],
254  _dummy_pm,
256  rot(_plastic_strain_old[_qp]),
257  _plastic_strain[_qp],
260  _yf[_qp],
261  number_iterations,
262  _Jacobian_mult[_qp],
264  true);
265 
266  // if not purely elastic or the customised stuff failed, do some plastic return
267  if (!found_solution)
269  _stress[_qp],
270  _intnl_old[_qp],
271  _intnl[_qp],
272  rot(_plastic_strain_old[_qp]),
273  _plastic_strain[_qp],
276  _yf[_qp],
277  number_iterations,
278  linesearch_needed,
279  ld_encountered,
280  constraints_added,
281  _Jacobian_mult[_qp]);
282 
283  if (_cosserat)
284  {
285  (*_couple_stress)[_qp] = (*_elastic_flexural_rigidity_tensor)[_qp] * _my_curvature;
286  (*_Jacobian_mult_couple)[_qp] = _my_flexural_rigidity_tensor;
287  }
288 
289  postReturnMap(); // rotate back from new frame if necessary
290 
291  _iter[_qp] = 1.0 * number_iterations;
292  _linesearch_needed[_qp] = linesearch_needed;
293  _ld_encountered[_qp] = ld_encountered;
294  _constraints_added[_qp] = constraints_added;
295 
296  // Update measures of strain
298  (_plastic_strain[_qp] - _plastic_strain_old[_qp]);
299 
300  // Rotate the tensors to the current configuration
302  {
303  _stress[_qp] = _rotation_increment[_qp] * _stress[_qp] * _rotation_increment[_qp].transpose();
304  _elastic_strain[_qp] =
305  _rotation_increment[_qp] * _elastic_strain[_qp] * _rotation_increment[_qp].transpose();
306  _plastic_strain[_qp] =
307  _rotation_increment[_qp] * _plastic_strain[_qp] * _rotation_increment[_qp].transpose();
308  }
309 }
310 
313 {
314  if (!_n_supplied)
315  return tens;
316  return tens.rotated(_rot);
317 }
318 
319 void
321 {
322  if (_n_supplied)
323  {
324  // this is a rotation matrix that will rotate _n to the "z" axis
325  _rot = RotationMatrix::rotVecToZ(_n[_qp]);
326 
327  // rotate the tensors to this frame
328  _my_elasticity_tensor.rotate(_rot);
329  _my_strain_increment.rotate(_rot);
330  if (_cosserat)
331  {
333  _my_curvature.rotate(_rot);
334  }
335  }
336 }
337 
338 void
340 {
341  if (_n_supplied)
342  {
343  // this is a rotation matrix that will rotate "z" axis back to _n
344  _rot = _rot.transpose();
345 
346  // rotate the tensors back to original frame where _n is correctly oriented
347  _my_elasticity_tensor.rotate(_rot);
348  _Jacobian_mult[_qp].rotate(_rot);
349  _my_strain_increment.rotate(_rot);
350  _stress[_qp].rotate(_rot);
351  _plastic_strain[_qp].rotate(_rot);
352  if (_cosserat)
353  {
355  (*_Jacobian_mult_couple)[_qp].rotate(_rot);
356  _my_curvature.rotate(_rot);
357  (*_couple_stress)[_qp].rotate(_rot);
358  }
359 
360  // Rotate n by _rotation_increment
361  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
362  {
363  _n[_qp](i) = 0;
364  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
365  _n[_qp](i) += _rotation_increment[_qp](i, j) * _n_old[_qp](j);
366  }
367  }
368 }
369 
370 bool
372  RankTwoTensor & stress,
373  const std::vector<Real> & intnl_old,
374  std::vector<Real> & intnl,
375  std::vector<Real> & pm,
376  std::vector<Real> & cumulative_pm,
377  const RankTwoTensor & plastic_strain_old,
378  RankTwoTensor & plastic_strain,
379  const RankFourTensor & E_ijkl,
380  const RankTwoTensor & strain_increment,
381  std::vector<Real> & yf,
382  unsigned int & iterations,
383  RankFourTensor & consistent_tangent_operator,
384  const quickStep_called_from_t called_from,
385  bool final_step)
386 {
387  iterations = 0;
388 
389  unsigned num_plastic_returns;
390  RankTwoTensor delta_dp;
391 
392  // the following does the customized returnMap algorithm
393  // for all the plastic models.
394  unsigned custom_model = 0;
395  bool successful_return = returnMapAll(stress_old + E_ijkl * strain_increment,
396  intnl_old,
397  E_ijkl,
398  _epp_tol,
399  stress,
400  intnl,
401  pm,
402  cumulative_pm,
403  delta_dp,
404  yf,
405  num_plastic_returns,
406  custom_model);
407 
408  // the following updates the plastic_strain, when necessary
409  // and calculates the consistent_tangent_operator, when necessary
410  if (num_plastic_returns == 0)
411  {
412  // if successful_return = true, then a purely elastic step
413  // if successful_return = false, then >=1 plastic model is in
414  // inadmissible zone and failed to return using its customized
415  // returnMap function.
416  // In either case:
417  plastic_strain = plastic_strain_old;
418  if (successful_return && final_step)
419  {
420  if (called_from == computeQpStress_function)
421  consistent_tangent_operator = E_ijkl;
422  else // cannot necessarily use E_ijkl since different plastic models may have been active
423  // during other substeps
424  consistent_tangent_operator =
425  consistentTangentOperator(stress, intnl, E_ijkl, pm, cumulative_pm);
426  }
427  return successful_return;
428  }
429  else if (num_plastic_returns == 1 && successful_return)
430  {
431  // one model has successfully completed its custom returnMap algorithm
432  // and the other models have signalled they are elastic at
433  // the trial stress
434  plastic_strain = plastic_strain_old + delta_dp;
435  if (final_step)
436  {
437  if (called_from == computeQpStress_function && _f[custom_model]->useCustomCTO())
438  {
440  consistent_tangent_operator = E_ijkl;
441  else
442  {
443  std::vector<Real> custom_model_pm;
444  for (unsigned surface = 0; surface < _f[custom_model]->numberSurfaces(); ++surface)
445  custom_model_pm.push_back(cumulative_pm[_surfaces_given_model[custom_model][surface]]);
446  consistent_tangent_operator =
447  _f[custom_model]->consistentTangentOperator(stress_old + E_ijkl * strain_increment,
448  intnl_old[custom_model],
449  stress,
450  intnl[custom_model],
451  E_ijkl,
452  custom_model_pm);
453  }
454  }
455  else // cannot necessarily use the custom consistentTangentOperator since different plastic
456  // models may have been active during other substeps or the custom model says not to use
457  // its custom CTO algorithm
458  consistent_tangent_operator =
459  consistentTangentOperator(stress, intnl, E_ijkl, pm, cumulative_pm);
460  }
461  return true;
462  }
463  else // presumably returnMapAll is incorrectly coded!
464  mooseError("ComputeMultiPlasticityStress::quickStep should not get here!");
465 }
466 
467 bool
469  RankTwoTensor & stress,
470  const std::vector<Real> & intnl_old,
471  std::vector<Real> & intnl,
472  const RankTwoTensor & plastic_strain_old,
473  RankTwoTensor & plastic_strain,
474  const RankFourTensor & E_ijkl,
475  const RankTwoTensor & strain_increment,
476  std::vector<Real> & yf,
477  unsigned int & iterations,
478  bool & linesearch_needed,
479  bool & ld_encountered,
480  bool & constraints_added,
481  RankFourTensor & consistent_tangent_operator)
482 {
498  bool return_successful = false;
499 
500  // total number of Newton-Raphson iterations used
501  unsigned int iter = 0;
502 
503  Real step_size = 1.0;
504  Real time_simulated = 0.0;
505 
506  // the "good" variables hold the latest admissible stress
507  // and internal parameters.
508  RankTwoTensor stress_good = stress_old;
509  RankTwoTensor plastic_strain_good = plastic_strain_old;
510  std::vector<Real> intnl_good(_num_models);
511  for (unsigned model = 0; model < _num_models; ++model)
512  intnl_good[model] = intnl_old[model];
513  std::vector<Real> yf_good(_num_surfaces);
514 
515  // Following is necessary because I want strain_increment to be "const"
516  // but I also want to be able to subdivide an initial_stress
517  RankTwoTensor this_strain_increment = strain_increment;
518 
519  RankTwoTensor dep = step_size * this_strain_increment;
520 
521  _cumulative_pm.assign(_num_surfaces, 0);
522 
523  unsigned int num_consecutive_successes = 0;
524  while (time_simulated < 1.0 && step_size >= _min_stepsize)
525  {
526  iter = 0;
527  return_successful = returnMap(stress_good,
528  stress,
529  intnl_good,
530  intnl,
531  plastic_strain_good,
532  plastic_strain,
533  E_ijkl,
534  dep,
535  yf,
536  iter,
537  step_size <= _max_stepsize_for_dumb,
538  linesearch_needed,
539  ld_encountered,
540  constraints_added,
541  time_simulated + step_size >= 1,
542  consistent_tangent_operator,
544  iterations += iter;
545 
546  if (return_successful)
547  {
548  num_consecutive_successes += 1;
549  time_simulated += step_size;
550 
551  if (time_simulated < 1.0) // this condition is just for optimization: if time_simulated=1 then
552  // the "good" quantities are no longer needed
553  {
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)
561  step_size *= 1.2;
562  }
563  step_size = std::min(step_size, 1.0 - time_simulated); // avoid overshoots
564  }
565  else
566  {
567  step_size *= 0.5;
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];
573  yf.resize(_num_surfaces); // might have excited with junk
574  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
575  yf[surface] = yf_good[surface];
576  dep = step_size * this_strain_increment;
577  }
578  }
579 
580  if (!return_successful)
581  {
582  if (_ignore_failures)
583  {
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];
590  }
591  else
592  {
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";
596 
597  _fspb_debug_stress = stress_good + E_ijkl * dep;
598  _fspb_debug_pm.assign(
600  1); // this is chosen arbitrarily - please change if a more suitable value occurs to you!
602  for (unsigned model = 0; model < _num_models; ++model)
603  _fspb_debug_intnl[model] = intnl_good[model];
607  mooseError("Exiting\n");
608  }
609  }
610 
611  return return_successful;
612 }
613 
614 bool
616  RankTwoTensor & stress,
617  const std::vector<Real> & intnl_old,
618  std::vector<Real> & intnl,
619  const RankTwoTensor & plastic_strain_old,
620  RankTwoTensor & plastic_strain,
621  const RankFourTensor & E_ijkl,
622  const RankTwoTensor & strain_increment,
623  std::vector<Real> & f,
624  unsigned int & iter,
625  bool can_revert_to_dumb,
626  bool & linesearch_needed,
627  bool & ld_encountered,
628  bool & constraints_added,
629  bool final_step,
630  RankFourTensor & consistent_tangent_operator,
631  std::vector<Real> & cumulative_pm)
632 {
633 
634  // The "consistency parameters" (plastic multipliers)
635  // Change in plastic strain in this timestep = pm*flowPotential
636  // Each pm must be non-negative
637  std::vector<Real> pm;
638  pm.assign(_num_surfaces, 0.0);
639 
640  bool successful_return = quickStep(stress_old,
641  stress,
642  intnl_old,
643  intnl,
644  pm,
645  cumulative_pm,
646  plastic_strain_old,
647  plastic_strain,
648  E_ijkl,
649  strain_increment,
650  f,
651  iter,
652  consistent_tangent_operator,
654  final_step);
655 
656  if (successful_return)
657  return successful_return;
658 
659  // Here we know that the trial stress and intnl_old
660  // is inadmissible, and we have to return from those values
661  // value to the yield surface. There are three
662  // types of constraints we have to satisfy, listed
663  // below, and calculated in calculateConstraints(...)
664  // f<=0, epp=0, ic=0 (up to tolerances), and these are
665  // guaranteed to hold if nr_res2<=0.5
666  //
667  // Kuhn-Tucker conditions must also be satisfied
668  // These are:
669  // pm>=0, which may not hold upon exit of the NR loops
670  // due to _deactivation_scheme!=optimized;
671  // pm*f=0 (up to tolerances), which may not hold upon exit
672  // of the NR loops if a constraint got deactivated
673  // due to linear dependence, and then f<0, and its pm>0
674 
675  // Plastic strain constraint, L2 norm must be zero (up to a tolerance)
676  RankTwoTensor epp;
677 
678  // Yield function constraint passed to this function as
679  // std::vector<Real> & f
680  // Each yield function must be <= 0 (up to tolerance)
681  // Note that only the constraints that are active will be
682  // contained in f until the final few lines of returnMap
683 
684  // Internal constraint(s), must be zero (up to a tolerance)
685  // Note that only the constraints that are active will be
686  // contained in ic.
687  std::vector<Real> ic;
688 
689  // Record the stress before Newton-Raphson in case of failure-and-restart
690  RankTwoTensor initial_stress = stress;
691 
692  iter = 0;
693 
694  // Initialize the set of active constraints
695  // At this stage, the active constraints are
696  // those that exceed their _f_tol
697  // active constraints.
698  std::vector<bool> act;
699  buildActiveConstraints(f, stress, intnl, E_ijkl, act);
700 
701  // Inverse of E_ijkl (assuming symmetric)
702  RankFourTensor E_inv = E_ijkl.invSymm();
703 
704  // convenience variable that holds the change in plastic strain incurred during the return
705  // delta_dp = plastic_strain - plastic_strain_old
706  // delta_dp = E^{-1}*(initial_stress - stress), where initial_stress = E*(strain -
707  // plastic_strain_old)
708  RankTwoTensor delta_dp = RankTwoTensor();
709 
710  // whether single step was successful (whether line search was successful, and whether turning off
711  // constraints was successful)
712  bool single_step_success = true;
713 
714  // deactivation scheme
716 
717  // For complicated deactivation schemes we have to record the initial active set
718  std::vector<bool> initial_act;
719  initial_act.resize(_num_surfaces);
723  {
724  // if "optimized" fails we can change the deactivation scheme to "safe", etc
725  deact_scheme = optimized;
726  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
727  initial_act[surface] = act[surface];
728  }
729 
731  deact_scheme = safe;
732 
733  // For "dumb" deactivation, the active set takes all combinations until a solution is found
734  int dumb_iteration = 0;
735  std::vector<unsigned int> dumb_order;
736 
737  if (_deactivation_scheme == dumb ||
738  (_deactivation_scheme == optimized_to_safe_to_dumb && can_revert_to_dumb) ||
739  (_deactivation_scheme == safe_to_dumb && can_revert_to_dumb) ||
740  (_deactivation_scheme == optimized_to_dumb && can_revert_to_dumb))
741  buildDumbOrder(stress, intnl, dumb_order);
742 
743  if (_deactivation_scheme == dumb)
744  {
745  incrementDumb(dumb_iteration, dumb_order, act);
746  yieldFunction(stress, intnl, act, f);
747  }
748 
749  // To avoid any re-trials of "act" combinations that
750  // we've already tried and rejected, i record the
751  // combinations in actives_tried
752  std::set<unsigned int> actives_tried;
753  actives_tried.insert(activeCombinationNumber(act));
754 
755  // The residual-squared that the line-search will reduce
756  // Later it will get contributions from epp and ic, but
757  // at present these are zero
758  Real nr_res2 = 0;
759  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
760  if (act[surface])
761  nr_res2 += 0.5 * Utility::pow<2>(f[surface] / _f[modelNumber(surface)]->_f_tol);
762 
763  successful_return = false;
764 
765  bool still_finding_solution = true;
766  while (still_finding_solution)
767  {
768  single_step_success = true;
769  unsigned int local_iter = 0;
770 
771  // The Newton-Raphson loops
772  while (nr_res2 > 0.5 && local_iter++ < _max_iter && single_step_success)
773  single_step_success = singleStep(nr_res2,
774  stress,
775  intnl_old,
776  intnl,
777  pm,
778  delta_dp,
779  E_inv,
780  f,
781  epp,
782  ic,
783  act,
784  deact_scheme,
785  linesearch_needed,
786  ld_encountered);
787 
788  bool nr_good = (nr_res2 <= 0.5 && local_iter <= _max_iter && single_step_success);
789 
790  iter += local_iter;
791 
792  // 'act' might have changed due to using deact_scheme = optimized, so
793  actives_tried.insert(activeCombinationNumber(act));
794 
795  if (!nr_good)
796  {
797  // failure of NR routine.
798  // We might be able to change the deactivation_scheme and
799  // then re-try the NR starting from the initial values
800  // Or, if deact_scheme == "dumb", we just increarse the
801  // dumb_iteration number and re-try
802  bool change_scheme = false;
803  bool increment_dumb = false;
804  change_scheme = canChangeScheme(deact_scheme, can_revert_to_dumb);
805  if (!change_scheme && deact_scheme == dumb)
806  increment_dumb = canIncrementDumb(dumb_iteration);
807 
808  still_finding_solution = (change_scheme || increment_dumb);
809 
810  if (change_scheme)
811  changeScheme(initial_act,
812  can_revert_to_dumb,
813  initial_stress,
814  intnl_old,
815  deact_scheme,
816  act,
817  dumb_iteration,
818  dumb_order);
819 
820  if (increment_dumb)
821  incrementDumb(dumb_iteration, dumb_order, act);
822 
823  if (!still_finding_solution)
824  {
825  // we cannot change the scheme, or have run out of "dumb" options
826  successful_return = false;
827  break;
828  }
829  }
830 
831  bool kt_good = false;
832  if (nr_good)
833  {
834  // check Kuhn-Tucker
835  kt_good = checkKuhnTucker(f, pm, act);
836  if (!kt_good)
837  {
838  if (deact_scheme != dumb)
839  {
840  applyKuhnTucker(f, pm, act);
841 
842  // true if we haven't tried this active set before
843  still_finding_solution =
844  (actives_tried.find(activeCombinationNumber(act)) == actives_tried.end());
845  if (!still_finding_solution)
846  {
847  // must have tried turning off the constraints already.
848  // so try changing the scheme
849  if (canChangeScheme(deact_scheme, can_revert_to_dumb))
850  {
851  still_finding_solution = true;
852  changeScheme(initial_act,
853  can_revert_to_dumb,
854  initial_stress,
855  intnl_old,
856  deact_scheme,
857  act,
858  dumb_iteration,
859  dumb_order);
860  }
861  }
862  }
863  else
864  {
865  bool increment_dumb = false;
866  increment_dumb = canIncrementDumb(dumb_iteration);
867  still_finding_solution = increment_dumb;
868  if (increment_dumb)
869  incrementDumb(dumb_iteration, dumb_order, act);
870  }
871 
872  if (!still_finding_solution)
873  {
874  // have tried turning off the constraints already,
875  // or have run out of "dumb" options
876  successful_return = false;
877  break;
878  }
879  }
880  }
881 
882  bool admissible = false;
883  if (nr_good && kt_good)
884  {
885  // check admissible
886  std::vector<Real> all_f;
887  if (_num_surfaces == 1)
888  admissible = true; // for a single surface if NR has exited successfully then (stress,
889  // intnl) must be admissible
890  else
891  admissible = checkAdmissible(stress, intnl, all_f);
892 
893  if (!admissible)
894  {
895  // Not admissible.
896  // We can try adding constraints back in
897  // We can try changing the deactivation scheme
898  // Or, if deact_scheme == dumb, just increase dumb_iteration
899  bool add_constraints = canAddConstraints(act, all_f);
900  if (add_constraints)
901  {
902  constraints_added = true;
903  std::vector<bool> act_plus(_num_surfaces,
904  false); // "act" with the positive constraints added in
905  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
906  if (act[surface] ||
907  (!act[surface] && (all_f[surface] > _f[modelNumber(surface)]->_f_tol)))
908  act_plus[surface] = true;
909  if (actives_tried.find(activeCombinationNumber(act_plus)) == actives_tried.end())
910  {
911  // haven't tried this combination of actives yet
912  constraints_added = true;
913  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
914  act[surface] = act_plus[surface];
915  }
916  else
917  add_constraints = false; // haven't managed to add a new combination
918  }
919 
920  bool change_scheme = false;
921  bool increment_dumb = false;
922 
923  if (!add_constraints)
924  change_scheme = canChangeScheme(deact_scheme, can_revert_to_dumb);
925 
926  if (!add_constraints && !change_scheme && deact_scheme == dumb)
927  increment_dumb = canIncrementDumb(dumb_iteration);
928 
929  still_finding_solution = (add_constraints || change_scheme || increment_dumb);
930 
931  if (change_scheme)
932  changeScheme(initial_act,
933  can_revert_to_dumb,
934  initial_stress,
935  intnl_old,
936  deact_scheme,
937  act,
938  dumb_iteration,
939  dumb_order);
940 
941  if (increment_dumb)
942  incrementDumb(dumb_iteration, dumb_order, act);
943 
944  if (!still_finding_solution)
945  {
946  // we cannot change the scheme, or have run out of "dumb" options
947  successful_return = false;
948  break;
949  }
950  }
951  }
952 
953  successful_return = (nr_good && admissible && kt_good);
954  if (successful_return)
955  break;
956 
957  if (still_finding_solution)
958  {
959  stress = initial_stress;
960  delta_dp = RankTwoTensor(); // back to zero change in plastic strain
961  for (unsigned model = 0; model < _num_models; ++model)
962  intnl[model] = intnl_old[model]; // back to old internal params
963  pm.assign(_num_surfaces, 0.0); // back to zero plastic multipliers
964 
965  unsigned num_active = numberActive(act);
966  if (num_active == 0)
967  {
968  successful_return = false;
969  break; // failure
970  }
971 
972  actives_tried.insert(activeCombinationNumber(act));
973 
974  // Since "act" set has changed, either by changing deact_scheme, or by KT failing, so need to
975  // re-calculate nr_res2
976  yieldFunction(stress, intnl, act, f);
977 
978  nr_res2 = 0;
979  unsigned ind = 0;
980  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
981  if (act[surface])
982  {
983  if (f[ind] > _f[modelNumber(surface)]->_f_tol)
984  nr_res2 += 0.5 * Utility::pow<2>(f[ind] / _f[modelNumber(surface)]->_f_tol);
985  ind++;
986  }
987  }
988  }
989 
990  // returned, with either success or failure
991  if (successful_return)
992  {
993  plastic_strain += delta_dp;
994 
995  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
996  cumulative_pm[surface] += pm[surface];
997 
998  if (final_step)
999  consistent_tangent_operator =
1000  consistentTangentOperator(stress, intnl, E_ijkl, pm, cumulative_pm);
1001 
1002  if (f.size() != _num_surfaces)
1003  {
1004  // at this stage f.size() = num_active, but we need to return with all the yield functions
1005  // evaluated, so:
1006  act.assign(_num_surfaces, true);
1007  yieldFunction(stress, intnl, act, f);
1008  }
1009  }
1010 
1011  return successful_return;
1012 }
1013 
1014 bool
1016  const std::vector<Real> & all_f)
1017 {
1018  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1019  if (!act[surface] && (all_f[surface] > _f[modelNumber(surface)]->_f_tol))
1020  return true;
1021  return false;
1022 }
1023 
1024 bool
1026  bool can_revert_to_dumb)
1027 {
1028  if (current_deactivation_scheme == optimized && _deactivation_scheme == optimized_to_safe)
1029  return true;
1030 
1031  if (current_deactivation_scheme == optimized && _deactivation_scheme == optimized_to_safe_to_dumb)
1032  return true;
1033 
1034  if (current_deactivation_scheme == safe && _deactivation_scheme == safe_to_dumb &&
1035  can_revert_to_dumb)
1036  return true;
1037 
1038  if (current_deactivation_scheme == safe && _deactivation_scheme == optimized_to_safe_to_dumb &&
1039  can_revert_to_dumb)
1040  return true;
1041 
1042  if (current_deactivation_scheme == optimized && _deactivation_scheme == optimized_to_dumb &&
1043  can_revert_to_dumb)
1044  return true;
1045 
1046  return false;
1047 }
1048 
1049 void
1050 ComputeMultiPlasticityStress::changeScheme(const std::vector<bool> & initial_act,
1051  bool can_revert_to_dumb,
1052  const RankTwoTensor & initial_stress,
1053  const std::vector<Real> & intnl_old,
1054  DeactivationSchemeEnum & current_deactivation_scheme,
1055  std::vector<bool> & act,
1056  int & dumb_iteration,
1057  std::vector<unsigned int> & dumb_order)
1058 {
1059  if (current_deactivation_scheme == optimized &&
1062  {
1063  current_deactivation_scheme = safe;
1064  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1065  act[surface] = initial_act[surface];
1066  }
1067  else if ((current_deactivation_scheme == safe &&
1070  can_revert_to_dumb) ||
1071  (current_deactivation_scheme == optimized && _deactivation_scheme == optimized_to_dumb &&
1072  can_revert_to_dumb))
1073  {
1074  current_deactivation_scheme = dumb;
1075  dumb_iteration = 0;
1076  buildDumbOrder(initial_stress, intnl_old, dumb_order);
1077  incrementDumb(dumb_iteration, dumb_order, act);
1078  }
1079 }
1080 
1081 bool
1083  RankTwoTensor & stress,
1084  const std::vector<Real> & intnl_old,
1085  std::vector<Real> & intnl,
1086  std::vector<Real> & pm,
1087  RankTwoTensor & delta_dp,
1088  const RankFourTensor & E_inv,
1089  std::vector<Real> & f,
1090  RankTwoTensor & epp,
1091  std::vector<Real> & ic,
1092  std::vector<bool> & active,
1093  DeactivationSchemeEnum deactivation_scheme,
1094  bool & linesearch_needed,
1095  bool & ld_encountered)
1096 {
1097  bool successful_step; // return value
1098 
1099  Real nr_res2_before_step = nr_res2;
1100  RankTwoTensor stress_before_step;
1101  std::vector<Real> intnl_before_step;
1102  std::vector<Real> pm_before_step;
1103  RankTwoTensor delta_dp_before_step;
1104 
1105  if (deactivation_scheme == optimized)
1106  {
1107  // we potentially use the "before_step" quantities, so record them here
1108  stress_before_step = stress;
1109  intnl_before_step.resize(_num_models);
1110  for (unsigned model = 0; model < _num_models; ++model)
1111  intnl_before_step[model] = intnl[model];
1112  pm_before_step.resize(_num_surfaces);
1113  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1114  pm_before_step[surface] = pm[surface];
1115  delta_dp_before_step = delta_dp;
1116  }
1117 
1118  // During the Newton-Raphson procedure, we'll be
1119  // changing the following parameters in order to
1120  // (attempt to) satisfy the constraints.
1121  RankTwoTensor dstress; // change in stress
1122  std::vector<Real> dpm; // change in plasticity multipliers ("consistency parameters"). For ALL
1123  // contraints (active and deactive)
1124  std::vector<Real>
1125  dintnl; // change in internal parameters. For ALL internal params (active and deactive)
1126 
1127  // The constraints that have been deactivated for this NR step
1128  // due to the flow directions being linearly dependent
1129  std::vector<bool> deact_ld;
1130  deact_ld.assign(_num_surfaces, false);
1131 
1132  /* After NR and linesearch, if _deactivation_scheme == "optimized", the
1133  * active plasticity multipliers are checked for non-negativity. If some
1134  * are negative then they are deactivated forever, and the NR step is
1135  * re-done starting from the _before_step quantities recorded above
1136  */
1137  bool constraints_changing = true;
1138  bool reinstated_actives;
1139  while (constraints_changing)
1140  {
1141  // calculate dstress, dpm and dintnl for one full Newton-Raphson step
1142  nrStep(stress, intnl_old, intnl, pm, E_inv, delta_dp, dstress, dpm, dintnl, active, deact_ld);
1143 
1144  for (unsigned surface = 0; surface < deact_ld.size(); ++surface)
1145  if (deact_ld[surface])
1146  {
1147  ld_encountered = true;
1148  break;
1149  }
1150 
1151  // perform a line search
1152  // The line-search will exit with updated values
1153  successful_step = lineSearch(nr_res2,
1154  stress,
1155  intnl_old,
1156  intnl,
1157  pm,
1158  E_inv,
1159  delta_dp,
1160  dstress,
1161  dpm,
1162  dintnl,
1163  f,
1164  epp,
1165  ic,
1166  active,
1167  deact_ld,
1168  linesearch_needed);
1169 
1170  if (!successful_step)
1171  // completely bomb out
1172  return successful_step;
1173 
1174  // See if any active constraints need to be removed, and the step re-done
1175  constraints_changing = false;
1176  if (deactivation_scheme == optimized)
1177  {
1178  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1179  if (active[surface] && pm[surface] < 0.0)
1180  constraints_changing = true;
1181  }
1182 
1183  if (constraints_changing)
1184  {
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)
1189  {
1190  if (active[surface] && pm[surface] < 0.0)
1191  {
1192  // turn off the constraint forever
1193  active[surface] = false;
1194  pm_before_step[surface] = 0.0;
1195  intnl_before_step[modelNumber(surface)] =
1196  intnl_old[modelNumber(surface)]; // don't want to muck-up hardening!
1197  }
1198  intnl[modelNumber(surface)] = intnl_before_step[modelNumber(surface)];
1199  pm[surface] = pm_before_step[surface];
1200  }
1201  if (numberActive(active) == 0)
1202  {
1203  // completely bomb out
1204  successful_step = false;
1205  return successful_step;
1206  }
1207  }
1208 
1209  // reinstate any active values that have been turned off due to linear-dependence
1210  reinstated_actives = reinstateLinearDependentConstraints(deact_ld);
1211  } // ends "constraints_changing" loop
1212 
1213  // if active constraints were reinstated then nr_res2 needs to be re-calculated so it is correct
1214  // upson returning
1215  if (reinstated_actives)
1216  {
1217  bool completely_converged = true;
1218  if (successful_step && nr_res2 < 0.5)
1219  {
1220  // Here we have converged to the correct solution if
1221  // all the yield functions are < 0. Excellent!
1222  //
1223  // This is quite tricky - perhaps i can refactor to make it more obvious.
1224  //
1225  // Because actives are now reinstated, the residual2
1226  // calculation below will give nr_res2 > 0.5, because it won't
1227  // realise that we only had to set the active-but-not-deactivated f = 0,
1228  // and not the entire active set. If we pass that nr_res2 back from
1229  // this function then the calling function will not realise we've converged!
1230  // Therefore, check for this case
1231  unsigned ind = 0;
1232  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1233  if (active[surface])
1234  if (f[ind++] > _f[modelNumber(surface)]->_f_tol)
1235  completely_converged = false;
1236  }
1237  else
1238  completely_converged = false;
1239 
1240  if (!completely_converged)
1241  nr_res2 = residual2(pm, f, epp, ic, active, deact_ld);
1242  }
1243 
1244  return successful_step;
1245 }
1246 
1247 bool
1249  std::vector<bool> & deactivated_due_to_ld)
1250 {
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;
1255 
1256  deactivated_due_to_ld.assign(_num_surfaces, false);
1257  return reinstated_actives;
1258 }
1259 
1260 unsigned
1261 ComputeMultiPlasticityStress::numberActive(const std::vector<bool> & active)
1262 {
1263  unsigned num_active = 0;
1264  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1265  if (active[surface])
1266  num_active++;
1267  return num_active;
1268 }
1269 
1270 bool
1272  const std::vector<Real> & intnl,
1273  std::vector<Real> & all_f)
1274 {
1275  std::vector<bool> act;
1276  act.assign(_num_surfaces, true);
1277 
1278  yieldFunction(stress, intnl, act, all_f);
1279 
1280  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1281  if (all_f[surface] > _f[modelNumber(surface)]->_f_tol)
1282  return false;
1283 
1284  return true;
1285 }
1286 
1287 bool
1289  const std::vector<Real> & pm,
1290  const std::vector<bool> & active)
1291 {
1292  unsigned ind = 0;
1293  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1294  {
1295  if (active[surface])
1296  {
1297  if (f[ind++] < -_f[modelNumber(surface)]->_f_tol)
1298  if (pm[surface] != 0)
1299  return false;
1300  }
1301  else if (pm[surface] != 0)
1302  mooseError("Crash due to plastic multiplier not being zero. This occurred because of poor "
1303  "coding!!");
1304  }
1305  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1306  if (pm[surface] < 0)
1307  return false;
1308 
1309  return true;
1310 }
1311 
1312 void
1314  const std::vector<Real> & pm,
1315  std::vector<bool> & active)
1316 {
1317  bool turned_off = false;
1318  unsigned ind = 0;
1319 
1320  // turn off all active surfaces that have f<0 and pm!=0
1321  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1322  {
1323  if (active[surface])
1324  {
1325  if (f[ind++] < -_f[modelNumber(surface)]->_f_tol)
1326  if (pm[surface] != 0)
1327  {
1328  turned_off = true;
1329  active[surface] = false;
1330  }
1331  }
1332  else if (pm[surface] != 0)
1333  mooseError("Crash due to plastic multiplier not being zero. This occurred because of poor "
1334  "coding!!");
1335  }
1336 
1337  // if didn't turn off anything yet, turn off surface with minimum pm
1338  if (!turned_off)
1339  {
1340  int surface_to_turn_off = -1;
1341  Real min_pm = 0;
1342  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1343  if (pm[surface] < min_pm)
1344  {
1345  min_pm = pm[surface];
1346  surface_to_turn_off = surface;
1347  }
1348  if (surface_to_turn_off >= 0)
1349  active[surface_to_turn_off] = false;
1350  }
1351 }
1352 
1353 Real
1354 ComputeMultiPlasticityStress::residual2(const std::vector<Real> & pm,
1355  const std::vector<Real> & f,
1356  const RankTwoTensor & epp,
1357  const std::vector<Real> & ic,
1358  const std::vector<bool> & active,
1359  const std::vector<bool> & deactivated_due_to_ld)
1360 {
1361  Real nr_res2 = 0;
1362  unsigned ind = 0;
1363 
1364  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1365  if (active[surface])
1366  {
1367  if (!deactivated_due_to_ld[surface])
1368  {
1369  if (!(pm[surface] == 0 && f[ind] <= 0))
1370  nr_res2 += 0.5 * Utility::pow<2>(f[ind] / _f[modelNumber(surface)]->_f_tol);
1371  }
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);
1374  ind++;
1375  }
1376 
1377  nr_res2 += 0.5 * Utility::pow<2>(epp.L2norm() / _epp_tol);
1378 
1379  std::vector<bool> active_not_deact(_num_surfaces);
1380  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1381  active_not_deact[surface] = (active[surface] && !deactivated_due_to_ld[surface]);
1382  ind = 0;
1383  for (unsigned model = 0; model < _num_models; ++model)
1384  if (anyActiveSurfaces(model, active_not_deact))
1385  nr_res2 += 0.5 * Utility::pow<2>(ic[ind++] / _f[model]->_ic_tol);
1386 
1387  return nr_res2;
1388 }
1389 
1390 bool
1392  RankTwoTensor & stress,
1393  const std::vector<Real> & intnl_old,
1394  std::vector<Real> & intnl,
1395  std::vector<Real> & pm,
1396  const RankFourTensor & E_inv,
1397  RankTwoTensor & delta_dp,
1398  const RankTwoTensor & dstress,
1399  const std::vector<Real> & dpm,
1400  const std::vector<Real> & dintnl,
1401  std::vector<Real> & f,
1402  RankTwoTensor & epp,
1403  std::vector<Real> & ic,
1404  const std::vector<bool> & active,
1405  const std::vector<bool> & deactivated_due_to_ld,
1406  bool & linesearch_needed)
1407 {
1408  // Line search algorithm straight out of "Numerical Recipes"
1409 
1410  bool success =
1411  true; // return value: will be false if linesearch couldn't reduce the residual-squared
1412 
1413  // Aim is to decrease residual2
1414 
1415  Real lam = 1.0; // the line-search parameter: 1.0 is a full Newton step
1416  Real lam_min =
1417  1E-10; // minimum value of lam allowed - perhaps this should be dynamically calculated?
1418  Real f0 = nr_res2; // initial value of residual2
1419  Real slope = -2 * nr_res2; // "Numerical Recipes" uses -b*A*x, in order to check for roundoff, but
1420  // i hope the nrStep would warn if there were problems.
1421  Real tmp_lam; // cached value of lam used in quadratic & cubic line search
1422  Real f2 = nr_res2; // cached value of f = residual2 used in the cubic in the line search
1423  Real lam2 = lam; // cached value of lam used in the cubic in the line search
1424 
1425  // pm during the line-search
1426  std::vector<Real> ls_pm;
1427  ls_pm.resize(pm.size());
1428 
1429  // delta_dp during the line-search
1430  RankTwoTensor ls_delta_dp;
1431 
1432  // internal parameter during the line-search
1433  std::vector<Real> ls_intnl;
1434  ls_intnl.resize(intnl.size());
1435 
1436  // stress during the line-search
1437  RankTwoTensor ls_stress;
1438 
1439  // flow directions (not used in line search, but calculateConstraints returns this parameter)
1440  std::vector<RankTwoTensor> r;
1441 
1442  while (true)
1443  {
1444  // update the variables using this line-search parameter
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;
1451 
1452  // calculate the new active yield functions, epp and active internal constraints
1453  calculateConstraints(ls_stress, intnl_old, ls_intnl, ls_pm, ls_delta_dp, f, r, epp, ic, active);
1454 
1455  // calculate the new residual-squared
1456  nr_res2 = residual2(ls_pm, f, epp, ic, active, deactivated_due_to_ld);
1457 
1458  if (nr_res2 < f0 + 1E-4 * lam * slope)
1459  break;
1460  else if (lam < lam_min)
1461  {
1462  success = false;
1463  // restore plastic multipliers, yield functions, etc to original values
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];
1469  ls_stress = stress;
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);
1473  break;
1474  }
1475  else if (lam == 1.0)
1476  {
1477  // model as a quadratic
1478  tmp_lam = -slope / 2.0 / (nr_res2 - f0 - slope);
1479  }
1480  else
1481  {
1482  // model as a cubic
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);
1486  Real b =
1487  (-lam2 * rhs1 / Utility::pow<2>(lam) + lam * rhs2 / Utility::pow<2>(lam2)) / (lam - lam2);
1488  if (a == 0)
1489  tmp_lam = -slope / (2.0 * b);
1490  else
1491  {
1492  Real disc = Utility::pow<2>(b) - 3 * a * slope;
1493  if (disc < 0)
1494  tmp_lam = 0.5 * lam;
1495  else if (b <= 0)
1496  tmp_lam = (-b + std::sqrt(disc)) / (3.0 * a);
1497  else
1498  tmp_lam = -slope / (b + std::sqrt(disc));
1499  }
1500  if (tmp_lam > 0.5 * lam)
1501  tmp_lam = 0.5 * lam;
1502  }
1503  lam2 = lam;
1504  f2 = nr_res2;
1505  lam = std::max(tmp_lam, 0.1 * lam);
1506  }
1507 
1508  if (lam < 1.0)
1509  linesearch_needed = true;
1510 
1511  // assign the quantities found in the line-search
1512  // back to the originals
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];
1518  stress = ls_stress;
1519 
1520  return success;
1521 }
1522 
1523 void
1525  const std::vector<Real> & intnl,
1526  std::vector<unsigned int> & dumb_order)
1527 {
1528  if (dumb_order.size() != 0)
1529  return;
1530 
1531  std::vector<bool> act;
1532  act.assign(_num_surfaces, true);
1533 
1534  std::vector<Real> f;
1535  yieldFunction(stress, intnl, act, f);
1536  std::vector<RankTwoTensor> df_dstress;
1537  dyieldFunction_dstress(stress, intnl, act, df_dstress);
1538 
1539  typedef std::pair<Real, unsigned> pair_for_sorting;
1540  std::vector<pair_for_sorting> dist(_num_surfaces);
1541  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1542  {
1543  dist[surface].first = f[surface] / df_dstress[surface].L2norm();
1544  dist[surface].second = surface;
1545  }
1546  std::sort(dist.begin(), dist.end()); // sorted in ascending order of f/df_dstress
1547 
1548  dumb_order.resize(_num_surfaces);
1549  for (unsigned i = 0; i < _num_surfaces; ++i)
1550  dumb_order[i] = dist[_num_surfaces - 1 - i].second;
1551  // now dumb_order[0] is the surface with the greatest f/df_dstress
1552 }
1553 
1554 void
1556  const std::vector<unsigned int> & dumb_order,
1557  std::vector<bool> & act)
1558 {
1559  dumb_iteration += 1;
1560  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1561  act[dumb_order[surface]] =
1562  (dumb_iteration &
1563  (1 << surface)); // returns true if the surface_th bit of dumb_iteration == 1
1564 }
1565 
1566 bool
1568 {
1569  // (1 << _num_surfaces) = 2^_num_surfaces
1570  return ((dumb_iteration + 1) < (1 << _num_surfaces));
1571 }
1572 
1573 unsigned int
1575 {
1576  unsigned num = 0;
1577  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1578  if (act[surface])
1579  num += (1 << surface); // (1 << x) = 2^x
1580 
1581  return num;
1582 }
1583 
1586  const std::vector<Real> & intnl,
1587  const RankFourTensor & E_ijkl,
1588  const std::vector<Real> & pm_this_step,
1589  const std::vector<Real> & cumulative_pm)
1590 {
1591 
1593  return E_ijkl;
1594 
1595  // Typically act_at_some_step = act, but it is possible
1596  // that when subdividing a strain increment, a surface
1597  // is only active for one sub-step
1598  std::vector<bool> act_at_some_step(_num_surfaces);
1599  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1600  act_at_some_step[surface] = (cumulative_pm[surface] > 0);
1601 
1602  // "act" might contain surfaces that are linearly dependent
1603  // with others. Only the plastic multipliers that are > 0
1604  // for this strain increment need to be varied to find
1605  // the consistent tangent operator
1606  std::vector<bool> act_vary(_num_surfaces);
1607  for (unsigned surface = 0; surface < _num_surfaces; ++surface)
1608  act_vary[surface] = (pm_this_step[surface] > 0);
1609 
1610  std::vector<RankTwoTensor> df_dstress;
1611  dyieldFunction_dstress(stress, intnl, act_vary, df_dstress);
1612  std::vector<Real> df_dintnl;
1613  dyieldFunction_dintnl(stress, intnl, act_vary, df_dintnl);
1614  std::vector<RankTwoTensor> r;
1615  flowPotential(stress, intnl, act_vary, r);
1616  std::vector<RankFourTensor> dr_dstress_at_some_step;
1617  dflowPotential_dstress(stress, intnl, act_at_some_step, dr_dstress_at_some_step);
1618  std::vector<RankTwoTensor> dr_dintnl_at_some_step;
1619  dflowPotential_dintnl(stress, intnl, act_at_some_step, dr_dintnl_at_some_step);
1620  std::vector<Real> h;
1621  hardPotential(stress, intnl, act_vary, h);
1622 
1623  unsigned ind1;
1624  unsigned ind2;
1625 
1626  // r_minus_stuff[alpha] = r[alpha] -
1627  // pm_cumulatve[gamma]*dr[gamma]_dintnl[a]_at_some_step*h[a][alpha], with alpha only being in
1628  // act_vary, but gamma being act_at_some_step
1629  std::vector<RankTwoTensor> r_minus_stuff;
1630  ind1 = 0;
1631  for (unsigned surface1 = 0; surface1 < _num_surfaces; ++surface1)
1632  if (act_vary[surface1])
1633  {
1634  r_minus_stuff.push_back(r[ind1]);
1635  ind2 = 0;
1636  for (unsigned surface2 = 0; surface2 < _num_surfaces; ++surface2)
1637  if (act_at_some_step[surface2])
1638  {
1639  if (modelNumber(surface1) == modelNumber(surface2))
1640  {
1641  r_minus_stuff.back() -=
1642  cumulative_pm[surface2] * dr_dintnl_at_some_step[ind2] * h[ind1];
1643  }
1644  ind2++;
1645  }
1646  ind1++;
1647  }
1648 
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;
1653 
1654  // zzz is a matrix in the form that can be easily
1655  // inverted by MatrixTools::inverse
1656  // Eg for num_currently_active = 3
1657  // (zzz[0] zzz[1] zzz[2])
1658  // (zzz[3] zzz[4] zzz[5])
1659  // (zzz[6] zzz[7] zzz[8])
1660  std::vector<PetscScalar> zzz;
1661  zzz.assign(num_currently_active * num_currently_active, 0.0);
1662 
1663  ind1 = 0;
1664  RankTwoTensor r2;
1665  for (unsigned surface1 = 0; surface1 < _num_surfaces; ++surface1)
1666  if (act_vary[surface1])
1667  {
1668  ind2 = 0;
1669  for (unsigned surface2 = 0; surface2 < _num_surfaces; ++surface2)
1670  if (act_vary[surface2])
1671  {
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);
1674  if (modelNumber(surface1) == modelNumber(surface2))
1675  zzz[ind1 * num_currently_active + ind2] += df_dintnl[ind1] * h[ind2];
1676  ind2++;
1677  }
1678  ind1++;
1679  }
1680 
1681  if (num_currently_active > 0)
1682  {
1683  // invert zzz, in place. if num_currently_active = 0 then zzz is not needed.
1684  try
1685  {
1686  MatrixTools::inverse(zzz, num_currently_active);
1687  }
1688  catch (const MooseException & e)
1689  {
1690  // in the very rare case of zzz being singular, just return the "elastic" tangent operator
1691  return E_ijkl;
1692  }
1693  }
1694 
1695  RankFourTensor strain_coeff = E_ijkl;
1696  ind1 = 0;
1697  for (unsigned surface1 = 0; surface1 < _num_surfaces; ++surface1)
1698  if (act_vary[surface1])
1699  {
1700  RankTwoTensor part1 = E_ijkl * r_minus_stuff[ind1];
1701  ind2 = 0;
1702  for (unsigned surface2 = 0; surface2 < _num_surfaces; ++surface2)
1703  if (act_vary[surface2])
1704  {
1705  RankTwoTensor part2 = E_ijkl * df_dstress[ind2];
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];
1712  ind2++;
1713  }
1714  ind1++;
1715  }
1716 
1718  return strain_coeff;
1719 
1720  RankFourTensor stress_coeff(RankFourTensor::initIdentitySymmetricFour);
1721 
1722  RankFourTensor part3;
1723  ind1 = 0;
1724  for (unsigned surface1 = 0; surface1 < _num_surfaces; ++surface1)
1725  if (act_at_some_step[surface1])
1726  {
1727  part3 += cumulative_pm[surface1] * E_ijkl * dr_dstress_at_some_step[ind1];
1728  ind1++;
1729  }
1730 
1731  stress_coeff += part3;
1732 
1733  part3 = part3.transposeMajor(); // this is because below i want df_dstress[ind2]*part3, and this
1734  // equals (part3.transposeMajor())*df_dstress[ind2]
1735 
1736  ind1 = 0;
1737  for (unsigned surface1 = 0; surface1 < _num_surfaces; ++surface1)
1738  if (act_vary[surface1])
1739  {
1740  RankTwoTensor part1 = E_ijkl * r_minus_stuff[ind1];
1741  ind2 = 0;
1742  for (unsigned surface2 = 0; surface2 < _num_surfaces; ++surface2)
1743  if (act_vary[surface2])
1744  {
1745  RankTwoTensor part2 = part3 * df_dstress[ind2];
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];
1752  ind2++;
1753  }
1754  ind1++;
1755  }
1756 
1757  // need to find the inverse of stress_coeff, but remember
1758  // stress_coeff does not have the symmetries commonly found
1759  // in tensor mechanics:
1760  // stress_coeff(i, j, k, l) = stress_coeff(j, i, k, l) = stress_coeff(i, j, l, k) !=
1761  // stress_coeff(k, l, i, j)
1762  // (note the final not-equals). We want s_inv, such that
1763  // s_inv(i, j, m, n)*stress_coeff(m, n, k, l) = (de_ik*de_jl + de_il*de_jk)/2
1764  // where de_ij = 1 if i=j, and 0 otherwise.
1765  RankFourTensor s_inv;
1766  try
1767  {
1768  s_inv = stress_coeff.invSymm();
1769  }
1770  catch (const MooseException & e)
1771  {
1772  return strain_coeff; // when stress_coeff is singular (perhaps for incompressible plasticity?)
1773  // return the "linear" tangent operator
1774  }
1775 
1776  return s_inv * strain_coeff;
1777 }
defineLegacyParams
defineLegacyParams(ComputeMultiPlasticityStress)
ComputeMultiPlasticityStress::_n
MaterialProperty< RealVectorValue > & _n
current value of transverse direction
Definition: ComputeMultiPlasticityStress.h:135
ComputeMultiPlasticityStress::safe
Definition: ComputeMultiPlasticityStress.h:82
MultiPlasticityRawComponentAssembler::dflowPotential_dintnl
virtual void dflowPotential_dintnl(const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankTwoTensor > &dr_dintnl)
The derivative of the active flow potentials with respect to the active internal parameters The UserO...
Definition: MultiPlasticityRawComponentAssembler.C:235
MultiPlasticityDebugger::validParams
static InputParameters validParams()
Definition: MultiPlasticityDebugger.C:17
ComputeMultiPlasticityStress::_linesearch_needed
MaterialProperty< Real > & _linesearch_needed
Whether a line-search was needed in the latest Newton-Raphson process (1 if true, 0 otherwise)
Definition: ComputeMultiPlasticityStress.h:126
ComputeMultiPlasticityStress::_dummy_pm
std::vector< Real > _dummy_pm
dummy "consistency parameters" (plastic multipliers) used in quickStep when called from computeQpStre...
Definition: ComputeMultiPlasticityStress.h:65
ComputeMultiPlasticityStress::elastic
Definition: ComputeMultiPlasticityStress.h:56
MultiPlasticityRawComponentAssembler::modelNumber
unsigned int modelNumber(unsigned int surface)
returns the model number, given the surface number
Definition: MultiPlasticityRawComponentAssembler.C:785
ComputeMultiPlasticityStress::_ld_encountered
MaterialProperty< Real > & _ld_encountered
Whether linear-dependence was encountered in the latest Newton-Raphson process (1 if true,...
Definition: ComputeMultiPlasticityStress.h:129
ComputeMultiPlasticityStress::changeScheme
void changeScheme(const std::vector< bool > &initial_act, bool can_revert_to_dumb, const RankTwoTensor &initial_stress, const std::vector< Real > &intnl_old, DeactivationSchemeEnum &current_deactivation_scheme, std::vector< bool > &act, int &dumb_iteration, std::vector< unsigned int > &dumb_order)
Definition: ComputeMultiPlasticityStress.C:1050
ComputeMultiPlasticityStress::checkKuhnTucker
virtual bool checkKuhnTucker(const std::vector< Real > &f, const std::vector< Real > &pm, const std::vector< bool > &active)
Checks Kuhn-Tucker conditions, and alters "active" if appropriate.
Definition: ComputeMultiPlasticityStress.C:1288
ComputeMultiPlasticityStress::lineSearch
virtual bool lineSearch(Real &nr_res2, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, std::vector< Real > &pm, const RankFourTensor &E_inv, RankTwoTensor &delta_dp, const RankTwoTensor &dstress, const std::vector< Real > &dpm, const std::vector< Real > &dintnl, std::vector< Real > &f, RankTwoTensor &epp, std::vector< Real > &ic, const std::vector< bool > &active, const std::vector< bool > &deactivated_due_to_ld, bool &linesearch_needed)
Performs a line search.
Definition: ComputeMultiPlasticityStress.C:1391
MultiPlasticityLinearSystem::nrStep
virtual void nrStep(const RankTwoTensor &stress, const std::vector< Real > &intnl_old, const std::vector< Real > &intnl, const std::vector< Real > &pm, const RankFourTensor &E_inv, const RankTwoTensor &delta_dp, RankTwoTensor &dstress, std::vector< Real > &dpm, std::vector< Real > &dintnl, const std::vector< bool > &active, std::vector< bool > &deactivated_due_to_ld)
Performs one Newton-Raphson step.
Definition: MultiPlasticityLinearSystem.C:615
ComputeMultiPlasticityStress::rot
RankTwoTensor rot(const RankTwoTensor &tens)
Definition: ComputeMultiPlasticityStress.C:312
ComputeMultiPlasticityStress::canIncrementDumb
bool canIncrementDumb(int dumb_iteration)
Definition: ComputeMultiPlasticityStress.C:1567
ComputeStressBase::_stress
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
Definition: ComputeStressBase.h:50
ComputeMultiPlasticityStress
ComputeMultiPlasticityStress performs the return-map algorithm and associated stress updates for plas...
Definition: ComputeMultiPlasticityStress.h:30
ComputeMultiPlasticityStress::_tangent_operator_type
enum ComputeMultiPlasticityStress::TangentOperatorEnum _tangent_operator_type
ComputeMultiPlasticityStress::returnMap
virtual bool returnMap(const RankTwoTensor &stress_old, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, const RankTwoTensor &plastic_strain_old, RankTwoTensor &plastic_strain, const RankFourTensor &E_ijkl, const RankTwoTensor &strain_increment, std::vector< Real > &f, unsigned int &iter, bool can_revert_to_dumb, bool &linesearch_needed, bool &ld_encountered, bool &constraints_added, bool final_step, RankFourTensor &consistent_tangent_operator, std::vector< Real > &cumulative_pm)
Implements the return map.
Definition: ComputeMultiPlasticityStress.C:615
MultiPlasticityRawComponentAssembler::flowPotential
virtual void flowPotential(const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankTwoTensor > &r)
The active flow potential(s) - one for each yield function.
Definition: MultiPlasticityRawComponentAssembler.C:180
ComputeMultiPlasticityStress::incrementDumb
virtual void incrementDumb(int &dumb_iteration, const std::vector< unsigned int > &dumb_order, std::vector< bool > &act)
Increments "dumb_iteration" by 1, and sets "act" appropriately (act[alpha] = true iff alpha_th bit of...
Definition: ComputeMultiPlasticityStress.C:1555
ComputeStressBase::_Jacobian_mult
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
Definition: ComputeStressBase.h:61
ComputeMultiPlasticityStress::DeactivationSchemeEnum
DeactivationSchemeEnum
Definition: ComputeMultiPlasticityStress.h:79
ComputeMultiPlasticityStress::_n_input
RealVectorValue _n_input
the supplied transverse direction vector
Definition: ComputeMultiPlasticityStress.h:94
ComputeMultiPlasticityStress::_max_stepsize_for_dumb
Real _max_stepsize_for_dumb
"dumb" deactivation will only be used if the stepsize falls below this quantity
Definition: ComputeMultiPlasticityStress.h:48
MultiPlasticityDebugger::_fspb_debug_intnl
std::vector< Real > _fspb_debug_intnl
Debug the Jacobian entires at these internal parameters.
Definition: MultiPlasticityDebugger.h:71
ComputeMultiPlasticityStress::_epp_tol
Real _epp_tol
Tolerance on the plastic strain increment ("direction") constraint.
Definition: ComputeMultiPlasticityStress.h:62
ComputeMultiPlasticityStress::_intnl_old
const MaterialProperty< std::vector< Real > > & _intnl_old
old values of internal parameters
Definition: ComputeMultiPlasticityStress.h:117
ComputeMultiPlasticityStress::checkAdmissible
virtual bool checkAdmissible(const RankTwoTensor &stress, const std::vector< Real > &intnl, std::vector< Real > &all_f)
Checks whether the yield functions are in the admissible region.
Definition: ComputeMultiPlasticityStress.C:1271
ComputeMultiPlasticityStress::_min_stepsize
Real _min_stepsize
Minimum fraction of applied strain that may be applied during adaptive stepsizing.
Definition: ComputeMultiPlasticityStress.h:45
ComputeMultiPlasticityStress::ComputeMultiPlasticityStress
ComputeMultiPlasticityStress(const InputParameters &parameters)
Definition: ComputeMultiPlasticityStress.C:99
MultiPlasticityLinearSystem::calculateConstraints
virtual void calculateConstraints(const RankTwoTensor &stress, const std::vector< Real > &intnl_old, const std::vector< Real > &intnl, const std::vector< Real > &pm, const RankTwoTensor &delta_dp, std::vector< Real > &f, std::vector< RankTwoTensor > &r, RankTwoTensor &epp, std::vector< Real > &ic, const std::vector< bool > &active)
The constraints.
Definition: MultiPlasticityLinearSystem.C:228
ComputeMultiPlasticityStress::residual2
virtual Real residual2(const std::vector< Real > &pm, const std::vector< Real > &f, const RankTwoTensor &epp, const std::vector< Real > &ic, const std::vector< bool > &active, const std::vector< bool > &deactivated_due_to_ld)
The residual-squared.
Definition: ComputeMultiPlasticityStress.C:1354
ComputeMultiPlasticityStress::dumb
Definition: ComputeMultiPlasticityStress.h:83
MultiPlasticityRawComponentAssembler::_num_surfaces
unsigned int _num_surfaces
Number of surfaces within the plastic models.
Definition: MultiPlasticityRawComponentAssembler.h:66
ComputeMultiPlasticityStress::_cumulative_pm
std::vector< Real > _cumulative_pm
the sum of the plastic multipliers over all the sub-steps.
Definition: ComputeMultiPlasticityStress.h:71
ComputeMultiPlasticityStress::_plastic_strain
MaterialProperty< RankTwoTensor > & _plastic_strain
plastic strain
Definition: ComputeMultiPlasticityStress.h:108
ComputeMultiPlasticityStress::computeQpStress
virtual void computeQpStress()
Compute the stress and store it in the _stress material property for the current quadrature point.
Definition: ComputeMultiPlasticityStress.C:221
MultiPlasticityRawComponentAssembler::hardPotential
virtual void hardPotential(const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< Real > &h)
The active hardening potentials (one for each internal parameter and for each yield function) by assu...
Definition: MultiPlasticityRawComponentAssembler.C:262
MultiPlasticityDebugger::_fspb_debug
MooseEnum _fspb_debug
none - don't do any debugging crash - currently inactive jacobian - check the jacobian entries jacobi...
Definition: MultiPlasticityDebugger.h:62
ComputeMultiPlasticityStress::_my_curvature
RankTwoTensor _my_curvature
Curvature that can be rotated by this class, and split into multiple increments (ie,...
Definition: ComputeMultiPlasticityStress.h:183
MultiPlasticityRawComponentAssembler::returnMapAll
bool returnMapAll(const RankTwoTensor &trial_stress, const std::vector< Real > &intnl_old, const RankFourTensor &E_ijkl, Real ep_plastic_tolerance, RankTwoTensor &stress, std::vector< Real > &intnl, std::vector< Real > &pm, std::vector< Real > &cumulative_pm, RankTwoTensor &delta_dp, std::vector< Real > &yf, unsigned &num_successful_plastic_returns, unsigned &custom_model)
Performs a returnMap for each plastic model using their inbuilt returnMap functions.
Definition: MultiPlasticityRawComponentAssembler.C:599
ComputeMultiPlasticityStress::returnMap_function
Definition: ComputeMultiPlasticityStress.h:448
ComputeMultiPlasticityStress::consistentTangentOperator
RankFourTensor consistentTangentOperator(const RankTwoTensor &stress, const std::vector< Real > &intnl, const RankFourTensor &E_ijkl, const std::vector< Real > &pm_this_step, const std::vector< Real > &cumulative_pm)
Computes the consistent tangent operator (another name for the jacobian = d(stress_rate)/d(strain_rat...
Definition: ComputeMultiPlasticityStress.C:1585
ComputeMultiPlasticityStress::quickStep
virtual bool quickStep(const RankTwoTensor &stress_old, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, std::vector< Real > &pm, std::vector< Real > &cumulative_pm, const RankTwoTensor &plastic_strain_old, RankTwoTensor &plastic_strain, const RankFourTensor &E_ijkl, const RankTwoTensor &strain_increment, std::vector< Real > &yf, unsigned int &iterations, RankFourTensor &consistent_tangent_operator, const quickStep_called_from_t called_from, bool final_step)
Attempts to find an admissible (stress, intnl) by using the customized return-map algorithms defined ...
Definition: ComputeMultiPlasticityStress.C:371
ComputeMultiPlasticityStress::_deactivation_scheme
enum ComputeMultiPlasticityStress::DeactivationSchemeEnum _deactivation_scheme
ComputeStressBase
ComputeStressBase is the base class for stress tensors.
Definition: ComputeStressBase.h:26
ComputeMultiPlasticityStress::_elasticity_tensor
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
Definition: ComputeMultiPlasticityStress.h:105
ComputeMultiPlasticityStress::_plastic_strain_old
const MaterialProperty< RankTwoTensor > & _plastic_strain_old
Old value of plastic strain.
Definition: ComputeMultiPlasticityStress.h:111
ComputeMultiPlasticityStress::plasticStep
virtual bool plasticStep(const RankTwoTensor &stress_old, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, const RankTwoTensor &plastic_strain_old, RankTwoTensor &plastic_strain, const RankFourTensor &E_ijkl, const RankTwoTensor &strain_increment, std::vector< Real > &yf, unsigned int &iterations, bool &linesearch_needed, bool &ld_encountered, bool &constraints_added, RankFourTensor &consistent_tangent_operator)
performs a plastic step
Definition: ComputeMultiPlasticityStress.C:468
ComputeMultiPlasticityStress::_rot
RealTensorValue _rot
rotation matrix that takes _n to (0, 0, 1)
Definition: ComputeMultiPlasticityStress.h:97
ComputeStressBase::validParams
static InputParameters validParams()
Definition: ComputeStressBase.C:17
MultiPlasticityDebugger::_fspb_debug_stress
RankTwoTensor _fspb_debug_stress
Debug the Jacobian entries at this stress.
Definition: MultiPlasticityDebugger.h:65
ComputeMultiPlasticityStress::buildDumbOrder
void buildDumbOrder(const RankTwoTensor &stress, const std::vector< Real > &intnl, std::vector< unsigned int > &dumb_order)
Builds the order which "dumb" activation will take.
Definition: ComputeMultiPlasticityStress.C:1524
ComputeMultiPlasticityStress::numberActive
virtual unsigned int numberActive(const std::vector< bool > &active)
counts the number of active constraints
Definition: ComputeMultiPlasticityStress.C:1261
MultiPlasticityRawComponentAssembler::_f
std::vector< const TensorMechanicsPlasticModel * > _f
User objects that define the yield functions, flow potentials, etc.
Definition: MultiPlasticityRawComponentAssembler.h:75
ComputeMultiPlasticityStress::_cosserat
bool _cosserat
whether Cosserat mechanics should be used
Definition: ComputeMultiPlasticityStress.h:156
ComputeMultiPlasticityStress::validParams
static InputParameters validParams()
Definition: ComputeMultiPlasticityStress.C:24
ComputeMultiPlasticityStress::_perform_finite_strain_rotations
bool _perform_finite_strain_rotations
whether to perform the rotations necessary in finite-strain simulations
Definition: ComputeMultiPlasticityStress.h:100
MultiPlasticityRawComponentAssembler::dyieldFunction_dintnl
virtual void dyieldFunction_dintnl(const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< Real > &df_dintnl)
The derivative of active yield function(s) with respect to their internal parameters (the user object...
Definition: MultiPlasticityRawComponentAssembler.C:153
ComputeStressBase::initQpStatefulProperties
virtual void initQpStatefulProperties() override
Definition: ComputeStressBase.C:43
MultiPlasticityRawComponentAssembler::buildActiveConstraints
virtual void buildActiveConstraints(const std::vector< Real > &f, const RankTwoTensor &stress, const std::vector< Real > &intnl, const RankFourTensor &Eijkl, std::vector< bool > &act)
Constructs a set of active constraints, given the yield functions, f.
Definition: MultiPlasticityRawComponentAssembler.C:344
ComputeMultiPlasticityStress::optimized_to_dumb
Definition: ComputeMultiPlasticityStress.h:87
ComputeMultiPlasticityStress::optimized_to_safe
Definition: ComputeMultiPlasticityStress.h:84
ComputeMultiPlasticityStress::reinstateLinearDependentConstraints
virtual bool reinstateLinearDependentConstraints(std::vector< bool > &deactivated_due_to_ld)
makes all deactivated_due_to_ld false, and if >0 of them were initially true, returns true
Definition: ComputeMultiPlasticityStress.C:1248
ComputeMultiPlasticityStress::quickStep_called_from_t
quickStep_called_from_t
The functions from which quickStep can be called.
Definition: ComputeMultiPlasticityStress.h:445
registerMooseObject
registerMooseObject("TensorMechanicsApp", ComputeMultiPlasticityStress)
ComputeMultiPlasticityStress::_rotation_increment
const MaterialProperty< RankTwoTensor > & _rotation_increment
Rotation increment (coming from ComputeIncrementalSmallStrain, for example)
Definition: ComputeMultiPlasticityStress.h:147
ComputeMultiPlasticityStress::singleStep
virtual bool singleStep(Real &nr_res2, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, std::vector< Real > &pm, RankTwoTensor &delta_dp, const RankFourTensor &E_inv, std::vector< Real > &f, RankTwoTensor &epp, std::vector< Real > &ic, std::vector< bool > &active, DeactivationSchemeEnum deactivation_scheme, bool &linesearch_needed, bool &ld_encountered)
Performs a single Newton-Raphson + linesearch step Constraints are deactivated and the step is re-don...
Definition: ComputeMultiPlasticityStress.C:1082
ComputeMultiPlasticityStress::_constraints_added
MaterialProperty< Real > & _constraints_added
Whether constraints were added in during the latest Newton-Raphson process (1 if true,...
Definition: ComputeMultiPlasticityStress.h:132
MultiPlasticityRawComponentAssembler::dyieldFunction_dstress
virtual void dyieldFunction_dstress(const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankTwoTensor > &df_dstress)
The derivative of the active yield function(s) with respect to stress.
Definition: MultiPlasticityRawComponentAssembler.C:125
ComputeMultiPlasticityStress::preReturnMap
virtual void preReturnMap()
Definition: ComputeMultiPlasticityStress.C:320
ComputeMultiPlasticityStress::optimized
Definition: ComputeMultiPlasticityStress.h:81
ComputeMultiPlasticityStress::computeQpStress_function
Definition: ComputeMultiPlasticityStress.h:447
ComputeMultiPlasticityStress::safe_to_dumb
Definition: ComputeMultiPlasticityStress.h:85
MultiPlasticityDebugger.h
ComputeMultiPlasticityStress::canAddConstraints
bool canAddConstraints(const std::vector< bool > &act, const std::vector< Real > &all_f)
Definition: ComputeMultiPlasticityStress.C:1015
MultiPlasticityDebugger::_fspb_debug_pm
std::vector< Real > _fspb_debug_pm
Debug the Jacobian entires at these plastic multipliers.
Definition: MultiPlasticityDebugger.h:68
RankTwoTensor
RankTwoTensorTempl< Real > RankTwoTensor
Definition: ACGrGrElasticDrivingForce.h:17
ComputeStressBase::_elastic_strain
MaterialProperty< RankTwoTensor > & _elastic_strain
Elastic strain material property.
Definition: ComputeStressBase.h:52
MultiPlasticityDebugger::checkDerivatives
void checkDerivatives()
Checks the derivatives, eg dyieldFunction_dstress by using finite difference approximations.
Definition: MultiPlasticityDebugger.C:85
ComputeMultiPlasticityStress::_strain_increment
const MaterialProperty< RankTwoTensor > & _strain_increment
strain increment (coming from ComputeIncrementalSmallStrain, for example)
Definition: ComputeMultiPlasticityStress.h:141
MultiPlasticityRawComponentAssembler::_num_models
unsigned int _num_models
Number of plastic models for this material.
Definition: MultiPlasticityRawComponentAssembler.h:57
MultiPlasticityDebugger
MultiPlasticityDebugger computes various finite-difference things to help developers remove bugs in t...
Definition: MultiPlasticityDebugger.h:24
MultiPlasticityRawComponentAssembler::anyActiveSurfaces
bool anyActiveSurfaces(int model, const std::vector< bool > &active)
returns true if any internal surfaces of the given model are active according to 'active'
Definition: MultiPlasticityRawComponentAssembler.C:791
ComputeMultiPlasticityStress::_n_old
const MaterialProperty< RealVectorValue > & _n_old
old value of transverse direction
Definition: ComputeMultiPlasticityStress.h:138
ComputeMultiPlasticityStress::_stress_old
const MaterialProperty< RankTwoTensor > & _stress_old
Old value of stress.
Definition: ComputeMultiPlasticityStress.h:150
ComputeMultiPlasticityStress.h
ComputeMultiPlasticityStress::initQpStatefulProperties
virtual void initQpStatefulProperties()
Definition: ComputeMultiPlasticityStress.C:191
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
ComputeMultiPlasticityStress::_intnl
MaterialProperty< std::vector< Real > > & _intnl
internal parameters
Definition: ComputeMultiPlasticityStress.h:114
ComputeMultiPlasticityStress::_max_iter
unsigned int _max_iter
Maximum number of Newton-Raphson iterations allowed.
Definition: ComputeMultiPlasticityStress.h:42
ComputeMultiPlasticityStress::canChangeScheme
bool canChangeScheme(DeactivationSchemeEnum current_deactivation_scheme, bool can_revert_to_dumb)
Definition: ComputeMultiPlasticityStress.C:1025
ComputeMultiPlasticityStress::activeCombinationNumber
unsigned int activeCombinationNumber(const std::vector< bool > &act)
Definition: ComputeMultiPlasticityStress.C:1574
MultiPlasticityDebugger::checkJacobian
void checkJacobian(const RankFourTensor &E_inv, const std::vector< Real > &intnl_old)
Checks the full Jacobian, which is just certain linear combinations of the dyieldFunction_dstress,...
Definition: MultiPlasticityDebugger.C:161
ComputeMultiPlasticityStress::_ignore_failures
bool _ignore_failures
Even if the returnMap fails, return the best values found for stress and internal parameters.
Definition: ComputeMultiPlasticityStress.h:51
MultiPlasticityDebugger::checkSolution
void checkSolution(const RankFourTensor &E_inv)
Checks that Ax does equal b in the NR procedure.
Definition: MultiPlasticityDebugger.C:367
ComputeMultiPlasticityStress::_my_flexural_rigidity_tensor
RankFourTensor _my_flexural_rigidity_tensor
Flexual rigidity tensor that can be rotated by this class (ie, its not const)
Definition: ComputeMultiPlasticityStress.h:180
ComputeMultiPlasticityStress::TangentOperatorEnum
TangentOperatorEnum
The type of tangent operator to return. tangent operator = d(stress_rate)/d(strain_rate).
Definition: ComputeMultiPlasticityStress.h:54
ComputeMultiPlasticityStress::_my_elasticity_tensor
RankFourTensor _my_elasticity_tensor
Elasticity tensor that can be rotated by this class (ie, its not const)
Definition: ComputeMultiPlasticityStress.h:174
MultiPlasticityRawComponentAssembler::yieldFunction
virtual void yieldFunction(const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< Real > &f)
The active yield function(s)
Definition: MultiPlasticityRawComponentAssembler.C:98
ComputeMultiPlasticityStress::_iter
MaterialProperty< Real > & _iter
Number of Newton-Raphson iterations used in the return-map.
Definition: ComputeMultiPlasticityStress.h:123
RankTwoTensorTempl< Real >
ComputeMultiPlasticityStress::linear
Definition: ComputeMultiPlasticityStress.h:57
ComputeMultiPlasticityStress::_my_strain_increment
RankTwoTensor _my_strain_increment
Strain increment that can be rotated by this class, and split into multiple increments (ie,...
Definition: ComputeMultiPlasticityStress.h:177
ComputeMultiPlasticityStress::applyKuhnTucker
virtual void applyKuhnTucker(const std::vector< Real > &f, const std::vector< Real > &pm, std::vector< bool > &active)
Checks Kuhn-Tucker conditions, and alters "active" if appropriate.
Definition: ComputeMultiPlasticityStress.C:1313
ComputeMultiPlasticityStress::_elastic_strain_old
const MaterialProperty< RankTwoTensor > & _elastic_strain_old
Old value of elastic strain.
Definition: ComputeMultiPlasticityStress.h:153
MultiPlasticityRawComponentAssembler::dflowPotential_dstress
virtual void dflowPotential_dstress(const RankTwoTensor &stress, const std::vector< Real > &intnl, const std::vector< bool > &active, std::vector< RankFourTensor > &dr_dstress)
The derivative of the active flow potential(s) with respect to stress.
Definition: MultiPlasticityRawComponentAssembler.C:207
ComputeMultiPlasticityStress::postReturnMap
virtual void postReturnMap()
Definition: ComputeMultiPlasticityStress.C:339
MultiPlasticityRawComponentAssembler::_surfaces_given_model
std::vector< std::vector< unsigned int > > _surfaces_given_model
_surfaces_given_model[model_number] = vector of surface numbers for this model
Definition: MultiPlasticityRawComponentAssembler.h:69
ComputeMultiPlasticityStress::optimized_to_safe_to_dumb
Definition: ComputeMultiPlasticityStress.h:86
ComputeMultiPlasticityStress::_n_supplied
bool _n_supplied
User supplied the transverse direction vector.
Definition: ComputeMultiPlasticityStress.h:91
ComputeMultiPlasticityStress::_yf
MaterialProperty< std::vector< Real > > & _yf
yield functions
Definition: ComputeMultiPlasticityStress.h:120