LCOV - code coverage report
Current view: top level - src/materials - CappedMohrCoulombStressUpdate.C (source / functions) Hit Total Coverage
Test: idaholab/moose solid_mechanics: f45d79 Lines: 421 445 94.6 %
Date: 2025-07-25 05:00:39 Functions: 13 14 92.9 %
Legend: Lines: hit not hit

          Line data    Source code
       1             : //* This file is part of the MOOSE framework
       2             : //* https://mooseframework.inl.gov
       3             : //*
       4             : //* All rights reserved, see COPYRIGHT for full restrictions
       5             : //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
       6             : //*
       7             : //* Licensed under LGPL 2.1, please see LICENSE for details
       8             : //* https://www.gnu.org/licenses/lgpl-2.1.html
       9             : 
      10             : #include "CappedMohrCoulombStressUpdate.h"
      11             : #include "libmesh/utility.h"
      12             : #include "ElasticityTensorTools.h"
      13             : 
      14             : registerMooseObject("SolidMechanicsApp", CappedMohrCoulombStressUpdate);
      15             : 
      16             : InputParameters
      17        1560 : CappedMohrCoulombStressUpdate::validParams()
      18             : {
      19        1560 :   InputParameters params = MultiParameterPlasticityStressUpdate::validParams();
      20        3120 :   params.addRequiredParam<UserObjectName>(
      21             :       "tensile_strength",
      22             :       "A SolidMechanicsHardening UserObject that defines hardening of the "
      23             :       "tensile strength.  In physical situations this is positive (and always "
      24             :       "must be greater than negative compressive-strength.");
      25        3120 :   params.addRequiredParam<UserObjectName>(
      26             :       "compressive_strength",
      27             :       "A SolidMechanicsHardening UserObject that defines hardening of the "
      28             :       "compressive strength.  In physical situations this is positive.");
      29        3120 :   params.addRequiredParam<UserObjectName>(
      30             :       "cohesion", "A SolidMechanicsHardening UserObject that defines hardening of the cohesion");
      31        3120 :   params.addRequiredParam<UserObjectName>("friction_angle",
      32             :                                           "A SolidMechanicsHardening UserObject "
      33             :                                           "that defines hardening of the "
      34             :                                           "friction angle (in radians)");
      35        3120 :   params.addRequiredParam<UserObjectName>(
      36             :       "dilation_angle",
      37             :       "A SolidMechanicsHardening UserObject that defines hardening of the "
      38             :       "dilation angle (in radians).  Unless you are quite confident, this should "
      39             :       "be set positive and not greater than the friction angle.");
      40        3120 :   params.addParam<bool>("perfect_guess",
      41        3120 :                         true,
      42             :                         "Provide a guess to the Newton-Raphson procedure "
      43             :                         "that is the result from perfect plasticity.  With "
      44             :                         "severe hardening/softening this may be "
      45             :                         "suboptimal.");
      46        1560 :   params.addClassDescription("Nonassociative, smoothed, Mohr-Coulomb plasticity capped with "
      47             :                              "tensile (Rankine) and compressive caps, with hardening/softening");
      48        1560 :   return params;
      49           0 : }
      50             : 
      51        1170 : CappedMohrCoulombStressUpdate::CappedMohrCoulombStressUpdate(const InputParameters & parameters)
      52             :   : MultiParameterPlasticityStressUpdate(parameters, 3, 12, 2),
      53        1170 :     _tensile_strength(getUserObject<SolidMechanicsHardeningModel>("tensile_strength")),
      54        1170 :     _compressive_strength(getUserObject<SolidMechanicsHardeningModel>("compressive_strength")),
      55        1170 :     _cohesion(getUserObject<SolidMechanicsHardeningModel>("cohesion")),
      56        1170 :     _phi(getUserObject<SolidMechanicsHardeningModel>("friction_angle")),
      57        1170 :     _psi(getUserObject<SolidMechanicsHardeningModel>("dilation_angle")),
      58        2340 :     _perfect_guess(getParam<bool>("perfect_guess")),
      59        1170 :     _poissons_ratio(0.0),
      60        1170 :     _shifter(_f_tol),
      61        2340 :     _eigvecs(RankTwoTensor())
      62             : {
      63        1170 :   if (_psi.value(0.0) <= 0.0 || _psi.value(0.0) > _phi.value(0.0))
      64           0 :     mooseWarning("Usually the Mohr-Coulomb dilation angle is positive and not greater than the "
      65             :                  "friction angle");
      66        1170 : }
      67             : 
      68             : void
      69      144064 : CappedMohrCoulombStressUpdate::computeStressParams(const RankTwoTensor & stress,
      70             :                                                    std::vector<Real> & stress_params) const
      71             : {
      72             :   // stress_params[0] = smallest eigenvalue, stress_params[2] = largest eigenvalue
      73      144064 :   stress.symmetricEigenvalues(stress_params);
      74      144064 : }
      75             : 
      76             : std::vector<RankTwoTensor>
      77      108560 : CappedMohrCoulombStressUpdate::dstress_param_dstress(const RankTwoTensor & stress) const
      78             : {
      79             :   std::vector<Real> sp;
      80             :   std::vector<RankTwoTensor> dsp;
      81      108560 :   stress.dsymmetricEigenvalues(sp, dsp);
      82      108560 :   return dsp;
      83             : }
      84             : 
      85             : std::vector<RankFourTensor>
      86           0 : CappedMohrCoulombStressUpdate::d2stress_param_dstress(const RankTwoTensor & stress) const
      87             : {
      88             :   std::vector<RankFourTensor> d2;
      89           0 :   stress.d2symmetricEigenvalues(d2);
      90           0 :   return d2;
      91             : }
      92             : 
      93             : void
      94       96016 : CappedMohrCoulombStressUpdate::preReturnMapV(const std::vector<Real> & /*trial_stress_params*/,
      95             :                                              const RankTwoTensor & stress_trial,
      96             :                                              const std::vector<Real> & /*intnl_old*/,
      97             :                                              const std::vector<Real> & /*yf*/,
      98             :                                              const RankFourTensor & Eijkl)
      99             : {
     100             :   std::vector<Real> eigvals;
     101       96016 :   stress_trial.symmetricEigenvaluesEigenvectors(eigvals, _eigvecs);
     102       96016 :   _poissons_ratio = ElasticityTensorTools::getIsotropicPoissonsRatio(Eijkl);
     103       96016 : }
     104             : 
     105             : void
     106       96016 : CappedMohrCoulombStressUpdate::setStressAfterReturnV(const RankTwoTensor & /*stress_trial*/,
     107             :                                                      const std::vector<Real> & stress_params,
     108             :                                                      Real /*gaE*/,
     109             :                                                      const std::vector<Real> & /*intnl*/,
     110             :                                                      const yieldAndFlow & /*smoothed_q*/,
     111             :                                                      const RankFourTensor & /*Eijkl*/,
     112             :                                                      RankTwoTensor & stress) const
     113             : {
     114             :   // form the diagonal stress
     115       96016 :   stress = RankTwoTensor(stress_params[0], stress_params[1], stress_params[2], 0.0, 0.0, 0.0);
     116             :   // rotate to the original frame
     117       96016 :   stress = _eigvecs * stress * (_eigvecs.transpose());
     118       96016 : }
     119             : 
     120             : void
     121      228560 : CappedMohrCoulombStressUpdate::yieldFunctionValuesV(const std::vector<Real> & stress_params,
     122             :                                                     const std::vector<Real> & intnl,
     123             :                                                     std::vector<Real> & yf) const
     124             : {
     125             :   // intnl[0] = shear, intnl[1] = tensile
     126      228560 :   const Real ts = _tensile_strength.value(intnl[1]);
     127      228560 :   const Real cs = _compressive_strength.value(intnl[1]);
     128      228560 :   const Real sinphi = std::sin(_phi.value(intnl[0]));
     129      228560 :   const Real cohcos = _cohesion.value(intnl[0]) * std::cos(_phi.value(intnl[0]));
     130      228560 :   const Real smax = stress_params[2]; // largest eigenvalue
     131      228560 :   const Real smid = stress_params[1];
     132      228560 :   const Real smin = stress_params[0]; // smallest eigenvalue
     133      228560 :   yf[0] = smax - ts;
     134      228560 :   yf[1] = smid - ts;
     135      228560 :   yf[2] = smin - ts;
     136      228560 :   yf[3] = -smin - cs;
     137      228560 :   yf[4] = -smid - cs;
     138      228560 :   yf[5] = -smax - cs;
     139      228560 :   yf[6] = 0.5 * (smax - smin) + 0.5 * (smax + smin) * sinphi - cohcos;
     140      228560 :   yf[7] = 0.5 * (smid - smin) + 0.5 * (smid + smin) * sinphi - cohcos;
     141      228560 :   yf[8] = 0.5 * (smax - smid) + 0.5 * (smax + smid) * sinphi - cohcos;
     142      228560 :   yf[9] = 0.5 * (smid - smax) + 0.5 * (smid + smax) * sinphi - cohcos;
     143      228560 :   yf[10] = 0.5 * (smin - smid) + 0.5 * (smin + smid) * sinphi - cohcos;
     144      228560 :   yf[11] = 0.5 * (smin - smax) + 0.5 * (smin + smax) * sinphi - cohcos;
     145      228560 : }
     146             : 
     147             : void
     148      820880 : CappedMohrCoulombStressUpdate::computeAllQV(const std::vector<Real> & stress_params,
     149             :                                             const std::vector<Real> & intnl,
     150             :                                             std::vector<yieldAndFlow> & all_q) const
     151             : {
     152      820880 :   const Real ts = _tensile_strength.value(intnl[1]);
     153      820880 :   const Real dts = _tensile_strength.derivative(intnl[1]);
     154      820880 :   const Real cs = _compressive_strength.value(intnl[1]);
     155      820880 :   const Real dcs = _compressive_strength.derivative(intnl[1]);
     156      820880 :   const Real sinphi = std::sin(_phi.value(intnl[0]));
     157      820880 :   const Real dsinphi = std::cos(_phi.value(intnl[0])) * _phi.derivative(intnl[0]);
     158      820880 :   const Real sinpsi = std::sin(_psi.value(intnl[0]));
     159      820880 :   const Real dsinpsi = std::cos(_psi.value(intnl[0])) * _psi.derivative(intnl[0]);
     160      820880 :   const Real cohcos = _cohesion.value(intnl[0]) * std::cos(_phi.value(intnl[0]));
     161             :   const Real dcohcos =
     162      820880 :       _cohesion.derivative(intnl[0]) * std::cos(_phi.value(intnl[0])) -
     163      820880 :       _cohesion.value(intnl[0]) * std::sin(_phi.value(intnl[0])) * _phi.derivative(intnl[0]);
     164      820880 :   const Real smax = stress_params[2]; // largest eigenvalue
     165      820880 :   const Real smid = stress_params[1];
     166      820880 :   const Real smin = stress_params[0]; // smallest eigenvalue
     167             : 
     168             :   // yield functions.  See comment in yieldFunctionValuesV
     169      820880 :   all_q[0].f = smax - ts;
     170      820880 :   all_q[1].f = smid - ts;
     171      820880 :   all_q[2].f = smin - ts;
     172      820880 :   all_q[3].f = -smin - cs;
     173      820880 :   all_q[4].f = -smid - cs;
     174      820880 :   all_q[5].f = -smax - cs;
     175      820880 :   all_q[6].f = 0.5 * (smax - smin) + 0.5 * (smax + smin) * sinphi - cohcos;
     176      820880 :   all_q[7].f = 0.5 * (smid - smin) + 0.5 * (smid + smin) * sinphi - cohcos;
     177      820880 :   all_q[8].f = 0.5 * (smax - smid) + 0.5 * (smax + smid) * sinphi - cohcos;
     178      820880 :   all_q[9].f = 0.5 * (smid - smax) + 0.5 * (smid + smax) * sinphi - cohcos;
     179      820880 :   all_q[10].f = 0.5 * (smin - smid) + 0.5 * (smin + smid) * sinphi - cohcos;
     180      820880 :   all_q[11].f = 0.5 * (smin - smax) + 0.5 * (smin + smax) * sinphi - cohcos;
     181             : 
     182             :   // d(yield function)/d(stress_params)
     183    10671440 :   for (unsigned yf = 0; yf < _num_yf; ++yf)
     184    39402240 :     for (unsigned a = 0; a < _num_sp; ++a)
     185    29551680 :       all_q[yf].df[a] = 0.0;
     186      820880 :   all_q[0].df[2] = 1.0;
     187      820880 :   all_q[1].df[1] = 1.0;
     188      820880 :   all_q[2].df[0] = 1.0;
     189      820880 :   all_q[3].df[0] = -1.0;
     190      820880 :   all_q[4].df[1] = -1.0;
     191      820880 :   all_q[5].df[2] = -1.0;
     192      820880 :   all_q[6].df[2] = 0.5 * (1.0 + sinphi);
     193      820880 :   all_q[6].df[0] = 0.5 * (-1.0 + sinphi);
     194      820880 :   all_q[7].df[1] = 0.5 * (1.0 + sinphi);
     195      820880 :   all_q[7].df[0] = 0.5 * (-1.0 + sinphi);
     196      820880 :   all_q[8].df[2] = 0.5 * (1.0 + sinphi);
     197      820880 :   all_q[8].df[1] = 0.5 * (-1.0 + sinphi);
     198      820880 :   all_q[9].df[1] = 0.5 * (1.0 + sinphi);
     199      820880 :   all_q[9].df[2] = 0.5 * (-1.0 + sinphi);
     200      820880 :   all_q[10].df[0] = 0.5 * (1.0 + sinphi);
     201      820880 :   all_q[10].df[1] = 0.5 * (-1.0 + sinphi);
     202      820880 :   all_q[11].df[0] = 0.5 * (1.0 + sinphi);
     203      820880 :   all_q[11].df[2] = 0.5 * (-1.0 + sinphi);
     204             : 
     205             :   // d(yield function)/d(intnl)
     206     5746160 :   for (unsigned i = 0; i < 6; ++i)
     207     4925280 :     all_q[i].df_di[0] = 0.0;
     208      820880 :   all_q[0].df_di[1] = all_q[1].df_di[1] = all_q[2].df_di[1] = -dts;
     209      820880 :   all_q[3].df_di[1] = all_q[4].df_di[1] = all_q[5].df_di[1] = -dcs;
     210     5746160 :   for (unsigned i = 6; i < 12; ++i)
     211     4925280 :     all_q[i].df_di[1] = 0.0;
     212      820880 :   all_q[6].df_di[0] = 0.5 * (smax + smin) * dsinphi - dcohcos;
     213      820880 :   all_q[7].df_di[0] = 0.5 * (smid + smin) * dsinphi - dcohcos;
     214      820880 :   all_q[8].df_di[0] = 0.5 * (smax + smid) * dsinphi - dcohcos;
     215      820880 :   all_q[9].df_di[0] = 0.5 * (smid + smax) * dsinphi - dcohcos;
     216      820880 :   all_q[10].df_di[0] = 0.5 * (smin + smid) * dsinphi - dcohcos;
     217      820880 :   all_q[11].df_di[0] = 0.5 * (smin + smax) * dsinphi - dcohcos;
     218             : 
     219             :   // the flow potential is just the yield function with phi->psi
     220             :   // d(flow potential)/d(stress_params)
     221     5746160 :   for (unsigned yf = 0; yf < 6; ++yf)
     222    19701120 :     for (unsigned a = 0; a < _num_sp; ++a)
     223    14775840 :       all_q[yf].dg[a] = all_q[yf].df[a];
     224      820880 :   all_q[6].dg[2] = all_q[7].dg[1] = all_q[8].dg[2] = all_q[9].dg[1] = all_q[10].dg[0] =
     225      820880 :       all_q[11].dg[0] = 0.5 * (1.0 + sinpsi);
     226      820880 :   all_q[6].dg[0] = all_q[7].dg[0] = all_q[8].dg[1] = all_q[9].dg[2] = all_q[10].dg[1] =
     227      820880 :       all_q[11].dg[2] = 0.5 * (-1.0 + sinpsi);
     228             : 
     229             :   // d(flow potential)/d(stress_params)/d(intnl)
     230    10671440 :   for (unsigned yf = 0; yf < _num_yf; ++yf)
     231    39402240 :     for (unsigned a = 0; a < _num_sp; ++a)
     232    88655040 :       for (unsigned i = 0; i < _num_intnl; ++i)
     233    59103360 :         all_q[yf].d2g_di[a][i] = 0.0;
     234      820880 :   all_q[6].d2g_di[2][0] = all_q[7].d2g_di[1][0] = all_q[8].d2g_di[2][0] = all_q[9].d2g_di[1][0] =
     235      820880 :       all_q[10].d2g_di[0][0] = all_q[11].d2g_di[0][0] = 0.5 * dsinpsi;
     236      820880 :   all_q[6].d2g_di[0][0] = all_q[7].d2g_di[0][0] = all_q[8].d2g_di[1][0] = all_q[9].d2g_di[2][0] =
     237      820880 :       all_q[10].d2g_di[1][0] = all_q[11].d2g_di[2][0] = 0.5 * dsinpsi;
     238             : 
     239             :   // d(flow potential)/d(stress_params)/d(stress_params)
     240    10671440 :   for (unsigned yf = 0; yf < _num_yf; ++yf)
     241    39402240 :     for (unsigned a = 0; a < _num_sp; ++a)
     242   118206720 :       for (unsigned b = 0; b < _num_sp; ++b)
     243    88655040 :         all_q[yf].d2g[a][b] = 0.0;
     244      820880 : }
     245             : 
     246             : void
     247       96016 : CappedMohrCoulombStressUpdate::setEffectiveElasticity(const RankFourTensor & Eijkl)
     248             : {
     249             :   // Eijkl is required to be isotropic, so we can use the
     250             :   // frame where stress is diagonal
     251      384064 :   for (unsigned a = 0; a < _num_sp; ++a)
     252     1152192 :     for (unsigned b = 0; b < _num_sp; ++b)
     253      864144 :       _Eij[a][b] = Eijkl(a, a, b, b);
     254       96016 :   _En = _Eij[2][2];
     255       96016 :   const Real denom = _Eij[0][0] * (_Eij[0][0] + _Eij[0][1]) - 2 * Utility::pow<2>(_Eij[0][1]);
     256      384064 :   for (unsigned a = 0; a < _num_sp; ++a)
     257             :   {
     258      288048 :     _Cij[a][a] = (_Eij[0][0] + _Eij[0][1]) / denom;
     259      576096 :     for (unsigned b = 0; b < a; ++b)
     260      288048 :       _Cij[a][b] = _Cij[b][a] = -_Eij[0][1] / denom;
     261             :   }
     262       96016 : }
     263             : 
     264             : void
     265      107632 : CappedMohrCoulombStressUpdate::initializeVarsV(const std::vector<Real> & trial_stress_params,
     266             :                                                const std::vector<Real> & intnl_old,
     267             :                                                std::vector<Real> & stress_params,
     268             :                                                Real & gaE,
     269             :                                                std::vector<Real> & intnl) const
     270             : {
     271      107632 :   if (!_perfect_guess)
     272             :   {
     273           0 :     for (unsigned i = 0; i < _num_sp; ++i)
     274           0 :       stress_params[i] = trial_stress_params[i];
     275           0 :     gaE = 0.0;
     276             :   }
     277             :   else
     278             :   {
     279      107632 :     const Real ts = _tensile_strength.value(intnl_old[1]);
     280      107632 :     const Real cs = _compressive_strength.value(intnl_old[1]);
     281      107632 :     const Real sinphi = std::sin(_phi.value(intnl_old[0]));
     282      107632 :     const Real cohcos = _cohesion.value(intnl_old[0]) * std::cos(_phi.value(intnl_old[0]));
     283             : 
     284      107632 :     const Real trial_tensile_yf = trial_stress_params[2] - ts;
     285      107632 :     const Real trial_compressive_yf = -trial_stress_params[0] - cs;
     286      107632 :     const Real trial_mc_yf = 0.5 * (trial_stress_params[2] - trial_stress_params[0]) +
     287      107632 :                              0.5 * (trial_stress_params[2] + trial_stress_params[0]) * sinphi -
     288      107632 :                              cohcos;
     289             : 
     290             :     /**
     291             :      * For CappedMohrCoulomb there are a number of possibilities for the
     292             :      * returned stress configuration:
     293             :      * - return to the Tensile yield surface to its plane
     294             :      * - return to the Tensile yield surface to its max=mid edge
     295             :      * - return to the Tensile yield surface to its tip
     296             :      * - return to the Compressive yield surface to its plane
     297             :      * - return to the Compressive yield surface to its max=mid edge
     298             :      * - return to the Compressive  yield surface to its tip
     299             :      * - return to the Mohr-Coulomb yield surface to its plane
     300             :      * - return to the Mohr-Coulomb yield surface to its max=mid edge
     301             :      * - return to the Mohr-Coulomb yield surface to its mid=min edge
     302             :      * - return to the Mohr-Coulomb yield surface to its tip
     303             :      * - return to the edge between Tensile and Mohr-Coulomb
     304             :      * - return to the edge between Tensile and Mohr-Coulomb, at max=mid point
     305             :      * - return to the edge between Tensile and Mohr-Coulomb, at mid=min point
     306             :      * - return to the edge between Compressive and Mohr-Coulomb
     307             :      * - return to the edge between Compressive and Mohr-Coulomb, at max=mid point
     308             :      * - return to the edge between Compressive and Mohr-Coulomb, at mid=min point
     309             :      * Which of these is more prelevant depends on the parameters
     310             :      * tensile strength, compressive strength, cohesion, etc.
     311             :      * I simply check each possibility until i find one that works.
     312             :      * _shifter is used to avoid equal eigenvalues
     313             :      */
     314             : 
     315             :     bool found_solution = false;
     316             : 
     317      107632 :     if (trial_tensile_yf <= _f_tol && trial_compressive_yf <= _f_tol && trial_mc_yf <= _f_tol)
     318             :     {
     319             :       // this is needed because although we know the smoothed yield function is
     320             :       // positive, the individual yield functions may not be
     321         768 :       for (unsigned i = 0; i < _num_sp; ++i)
     322         576 :         stress_params[i] = trial_stress_params[i];
     323         192 :       gaE = 0.0;
     324             :       found_solution = true;
     325             :     }
     326             : 
     327      107632 :     const bool tensile_possible = (ts < cohcos / sinphi); // tensile chops MC tip
     328             :     const bool mc_tip_possible = (cohcos / sinphi < ts);  // MC tip pokes through tensile
     329      107632 :     const bool mc_impossible = (0.5 * (ts + cs) + 0.5 * (ts - cs) * sinphi - cohcos <
     330      107632 :                                 _f_tol); // MC outside tensile and compressive
     331             : 
     332      107632 :     const Real sinpsi = std::sin(_psi.value(intnl_old[0]));
     333      107632 :     const Real halfplus = 0.5 + 0.5 * sinpsi;
     334      107632 :     const Real neghalfplus = -0.5 + 0.5 * sinpsi;
     335             : 
     336      107632 :     if (!found_solution && tensile_possible && trial_tensile_yf > _f_tol &&
     337        7856 :         (trial_compressive_yf <= _f_tol || (trial_compressive_yf > _f_tol && mc_impossible)))
     338             :     {
     339             :       // try pure tensile failure, return to the plane
     340             :       // This involves solving yf[0] = 0 and the three flow-direction equations
     341             :       // Don't try this if there is compressive failure, since returning to
     342             :       // the tensile yield surface will only make compressive failure worse
     343       33104 :       const Real ga = (trial_stress_params[2] - ts) / _Eij[2][2];
     344       33104 :       stress_params[2] = ts; // largest eigenvalue
     345       33104 :       stress_params[1] = trial_stress_params[1] - ga * _Eij[1][2];
     346       33104 :       stress_params[0] = trial_stress_params[0] - ga * _Eij[0][2];
     347             : 
     348             :       // if we have to return to the edge, or tip, do that
     349             :       Real dist_mod = 1.0;
     350       33104 :       const Real to_subtract1 = stress_params[1] - (ts - 0.5 * _shifter);
     351       33104 :       if (to_subtract1 > 0.0)
     352             :       {
     353       11744 :         dist_mod += Utility::pow<2>(to_subtract1 / (trial_stress_params[2] - ts));
     354       11744 :         stress_params[1] -= to_subtract1;
     355             :       }
     356       33104 :       const Real to_subtract0 = stress_params[0] - (ts - _shifter);
     357       33104 :       if (to_subtract0 > 0.0)
     358             :       {
     359        6368 :         dist_mod += Utility::pow<2>(to_subtract0 / (trial_stress_params[2] - ts));
     360        6368 :         stress_params[0] -= to_subtract0;
     361             :       }
     362       33104 :       if (mc_impossible) // might have to shift up to the compressive yield surface
     363             :       {
     364       26560 :         const Real to_add0 = -stress_params[0] - cs;
     365       26560 :         if (to_add0 > 0.0)
     366             :         {
     367           0 :           dist_mod += Utility::pow<2>(to_add0 / (trial_stress_params[2] - ts));
     368           0 :           stress_params[0] += to_add0;
     369             :         }
     370       26560 :         const Real to_add1 = -cs + 0.5 * _shifter - stress_params[1];
     371       26560 :         if (to_add1 > 0.0)
     372             :         {
     373           0 :           dist_mod += Utility::pow<2>(to_add1 / (trial_stress_params[2] - ts));
     374           0 :           stress_params[1] += to_add1;
     375             :         }
     376             :       }
     377             : 
     378       33104 :       const Real new_compressive_yf = -stress_params[0] - cs;
     379       33104 :       const Real new_mc_yf = 0.5 * (stress_params[2] - stress_params[0]) +
     380       33104 :                              0.5 * (stress_params[2] + stress_params[0]) * sinphi - cohcos;
     381       33104 :       if (new_mc_yf <= _f_tol && new_compressive_yf <= _f_tol)
     382             :       {
     383       29104 :         gaE = std::sqrt(dist_mod) * (trial_stress_params[2] - stress_params[2]);
     384             :         found_solution = true;
     385             :       }
     386             :     }
     387       78528 :     if (!found_solution && trial_compressive_yf > _f_tol &&
     388        7856 :         (trial_tensile_yf <= _f_tol || (trial_tensile_yf > _f_tol && mc_impossible)))
     389             :     {
     390             :       // try pure compressive failure
     391             :       // Don't try this if there is tensile failure, since returning to
     392             :       // the compressive yield surface will only make tensile failure worse
     393       31792 :       const Real ga = (trial_stress_params[0] + cs) / _Eij[0][0]; // this is negative
     394       31792 :       stress_params[0] = -cs;
     395       31792 :       stress_params[1] = trial_stress_params[1] - ga * _Eij[1][0];
     396       31792 :       stress_params[2] = trial_stress_params[2] - ga * _Eij[2][0];
     397             : 
     398             :       // if we have to return to the edge, or tip, do that
     399             :       Real dist_mod = 1.0;
     400       31792 :       const Real to_add1 = -cs + 0.5 * _shifter - stress_params[1];
     401       31792 :       if (to_add1 > 0.0)
     402             :       {
     403       11712 :         dist_mod += Utility::pow<2>(to_add1 / (trial_stress_params[0] + cs));
     404       11712 :         stress_params[1] += to_add1;
     405             :       }
     406       31792 :       const Real to_add2 = -cs + _shifter - stress_params[2];
     407       31792 :       if (to_add2 > 0.0)
     408             :       {
     409        6688 :         dist_mod += Utility::pow<2>(to_add2 / (trial_stress_params[0] + cs));
     410        6688 :         stress_params[2] += to_add2;
     411             :       }
     412       31792 :       if (mc_impossible) // might have to shift down to the tensile yield surface
     413             :       {
     414       23792 :         const Real to_subtract2 = stress_params[2] - ts;
     415       23792 :         if (to_subtract2 > 0.0)
     416             :         {
     417           0 :           dist_mod += Utility::pow<2>(to_subtract2 / (trial_stress_params[0] + cs));
     418           0 :           stress_params[2] -= to_subtract2;
     419             :         }
     420       23792 :         const Real to_subtract1 = stress_params[1] - (ts - 0.5 * _shifter);
     421       23792 :         if (to_subtract1 > 0.0)
     422             :         {
     423           0 :           dist_mod += Utility::pow<2>(to_subtract1 / (trial_stress_params[0] + cs));
     424           0 :           stress_params[1] -= to_subtract1;
     425             :         }
     426             :       }
     427             : 
     428       31792 :       const Real new_tensile_yf = stress_params[2] - ts;
     429       31792 :       const Real new_mc_yf = 0.5 * (stress_params[2] - stress_params[0]) +
     430       31792 :                              0.5 * (stress_params[2] + stress_params[0]) * sinphi - cohcos;
     431       31792 :       if (new_mc_yf <= _f_tol && new_tensile_yf <= _f_tol)
     432             :       {
     433       28528 :         gaE = std::sqrt(dist_mod) * (-trial_stress_params[0] + stress_params[0]);
     434             :         found_solution = true;
     435             :       }
     436             :     }
     437      107632 :     if (!found_solution && !mc_impossible && trial_mc_yf > _f_tol)
     438             :     {
     439             :       // try pure shear failure, return to the plane
     440             :       // This involves solving yf[6]=0 and the three flow-direction equations
     441       49808 :       const Real ga = ((trial_stress_params[2] - trial_stress_params[0]) +
     442       49808 :                        (trial_stress_params[2] + trial_stress_params[0]) * sinphi - 2.0 * cohcos) /
     443       49808 :                       (_Eij[2][2] - _Eij[0][2] + (_Eij[2][2] + _Eij[0][2]) * sinpsi * sinphi);
     444       49808 :       stress_params[2] =
     445       49808 :           trial_stress_params[2] - ga * (_Eij[2][2] * halfplus + _Eij[2][0] * neghalfplus);
     446       49808 :       stress_params[1] = trial_stress_params[1] - ga * _Eij[1][0] * sinpsi;
     447       49808 :       stress_params[0] =
     448       49808 :           trial_stress_params[0] - ga * (_Eij[0][0] * neghalfplus + _Eij[0][2] * halfplus);
     449       49808 :       const Real f7 = 0.5 * (stress_params[1] - stress_params[0]) +
     450       49808 :                       0.5 * (stress_params[1] + stress_params[0]) * sinphi - cohcos;
     451       49808 :       const Real f8 = 0.5 * (stress_params[2] - stress_params[1]) +
     452       49808 :                       0.5 * (stress_params[2] + stress_params[1]) * sinphi - cohcos;
     453       49808 :       const Real new_tensile_yf = stress_params[2] - ts;
     454       49808 :       const Real new_compressive_yf = -stress_params[0] - cs;
     455             : 
     456       49808 :       if (f7 <= _f_tol && f8 <= _f_tol && new_tensile_yf <= _f_tol && new_compressive_yf <= _f_tol)
     457             :       {
     458       14912 :         gaE = ((trial_stress_params[2] - trial_stress_params[0]) -
     459       14912 :                (stress_params[2] - stress_params[0])) /
     460       14912 :               (_Eij[2][2] - _Eij[0][2]) * _Eij[2][2];
     461             :         found_solution = true;
     462             :       }
     463             :     }
     464      107632 :     if (!found_solution && !mc_impossible && trial_mc_yf > _f_tol)
     465             :     {
     466             :       // Try return to the max=mid MC line.
     467             :       // To return to the max=mid line, we need to solve f6 = 0 = f7 and
     468             :       // the three flow-direction equations.  In the flow-direction equations
     469             :       // there are two plasticity multipliers, which i call ga6 and ga7,
     470             :       // corresponding to the amounts of strain normal to the f6 and f7
     471             :       // directions, respectively.
     472             :       // So:
     473             :       // Smax = Smax^trial - ga6 Emax,a dg6/dSa - ga7 Emax,a dg7/dSa
     474             :       //      = Smax^trial - ga6 cmax6 - ga7 cmax7  (with cmax6 and cmax7 evaluated below)
     475             :       // Smid = Smid^trial - ga6 Emid,a dg6/dSa - ga7 Emid,a dg7/dSa
     476             :       //      = Smid^trial - ga6 cmid6 - ga7 cmid7
     477             :       // Smin = Smin^trial - ga6 Emin,a dg6/dSa - ga7 Emin,a dg7/dSa
     478             :       //      = Smin^trial - ga6 cmin6 - ga7 cmin7
     479       34896 :       const Real cmax6 = _Eij[2][2] * halfplus + _Eij[2][0] * neghalfplus;
     480       34896 :       const Real cmax7 = _Eij[2][1] * halfplus + _Eij[2][0] * neghalfplus;
     481             :       // const Real cmid6 = _Eij[1][2] * halfplus + _Eij[1][0] * neghalfplus;
     482             :       // const Real cmid7 = _Eij[1][1] * halfplus + _Eij[1][0] * neghalfplus;
     483       34896 :       const Real cmin6 = _Eij[0][2] * halfplus + _Eij[0][0] * neghalfplus;
     484       34896 :       const Real cmin7 = _Eij[0][1] * halfplus + _Eij[0][0] * neghalfplus;
     485             :       // Substituting these into f6 = 0 yields
     486             :       // 0 = f6_trial - ga6 (0.5(cmax6 - cmin6) + 0.5(cmax6 + cmin6)sinphi) - ga7 (0.5(cmax6 -
     487             :       // cmin6) + 0.5(cmax6 + cmin6)sinphi) = f6_trial - ga6 c6 - ga7 c7, where
     488       34896 :       const Real c6 = 0.5 * (cmax6 - cmin6) + 0.5 * (cmax6 + cmin6) * sinphi;
     489       34896 :       const Real c7 = 0.5 * (cmax7 - cmin7) + 0.5 * (cmax7 + cmin7) * sinphi;
     490             :       // It isn't too hard to check that the other equation is
     491             :       // 0 = f7_trial - ga6 c7 - ga7 c6
     492             :       // These equations may be inverted to yield
     493       34896 :       if (c6 != c7)
     494             :       {
     495             :         const Real f6_trial = trial_mc_yf;
     496       34896 :         const Real f7_trial = 0.5 * (trial_stress_params[1] - trial_stress_params[0]) +
     497       34896 :                               0.5 * (trial_stress_params[1] + trial_stress_params[0]) * sinphi -
     498       34896 :                               cohcos;
     499       34896 :         const Real descr = Utility::pow<2>(c6) - Utility::pow<2>(c7);
     500       34896 :         Real ga6 = (c6 * f6_trial - c7 * f7_trial) / descr;
     501       34896 :         Real ga7 = (-c7 * f6_trial + c6 * f7_trial) / descr;
     502             :         // and finally
     503       34896 :         stress_params[2] = trial_stress_params[2] - ga6 * cmax6 - ga7 * cmax7;
     504       34896 :         stress_params[0] = trial_stress_params[0] - ga6 * cmin6 - ga7 * cmin7;
     505       34896 :         stress_params[1] = stress_params[2] - 0.5 * _shifter;
     506             : 
     507       34896 :         Real f8 = 0.5 * (stress_params[2] - stress_params[1]) +
     508       34896 :                   0.5 * (stress_params[2] + stress_params[1]) * sinphi - cohcos;
     509             : 
     510       34896 :         if (mc_tip_possible && f8 > _f_tol)
     511             :         {
     512        5568 :           stress_params[2] = cohcos / sinphi;
     513        5568 :           stress_params[1] = stress_params[2] - 0.5 * _shifter;
     514        5568 :           stress_params[0] = stress_params[2] - _shifter;
     515             :           f8 = 0.0;
     516             :           ga6 = 1.0;
     517             :           ga7 = 1.0;
     518             :         }
     519             : 
     520       34896 :         const Real new_tensile_yf = stress_params[2] - ts;
     521       34896 :         const Real new_compressive_yf = -stress_params[0] - cs;
     522             : 
     523       34896 :         if (f8 <= _f_tol && new_tensile_yf <= _f_tol && new_compressive_yf <= _f_tol &&
     524       19904 :             ga6 >= 0.0 && ga7 >= 0.0)
     525             :         {
     526        9696 :           gaE = ((trial_stress_params[2] - trial_stress_params[0]) -
     527        9696 :                  (stress_params[2] - stress_params[0])) /
     528        9696 :                 (_Eij[2][2] - _Eij[0][2]) * _Eij[2][2];
     529             :           found_solution = true;
     530             :         }
     531             :       }
     532             :     }
     533      107632 :     if (!found_solution && !mc_impossible && trial_mc_yf > _f_tol)
     534             :     {
     535             :       // Try return to the mid=min line.
     536             :       // To return to the mid=min line, we need to solve f6 = 0 = f8 and
     537             :       // the three flow-direction equations.  In the flow-direction equations
     538             :       // there are two plasticity multipliers, which i call ga6 and ga8,
     539             :       // corresponding to the amounts of strain normal to the f6 and f8
     540             :       // directions, respectively.
     541             :       // So:
     542             :       // Smax = Smax^trial - ga6 Emax,a dg6/dSa - ga8 Emax,a dg8/dSa
     543             :       //      = Smax^trial - ga6 cmax6 - ga8 cmax8  (with cmax6 and cmax8 evaluated below)
     544             :       // Smid = Smid^trial - ga6 Emid,a dg6/dSa - ga8 Emid,a dg8/dSa
     545             :       //      = Smid^trial - ga6 cmid6 - ga8 cmid8
     546             :       // Smin = Smin^trial - ga6 Emin,a dg6/dSa - ga8 Emin,a dg8/dSa
     547             :       //      = Smin^trial - ga6 cmin6 - ga8 cmin8
     548       25200 :       const Real cmax6 = _Eij[2][2] * halfplus + _Eij[2][0] * neghalfplus;
     549       25200 :       const Real cmax8 = _Eij[2][2] * halfplus + _Eij[2][1] * neghalfplus;
     550             :       // const Real cmid6 = _Eij[1][2] * halfplus + _Eij[1][0] * neghalfplus;
     551             :       // const Real cmid8 = _Eij[1][2] * halfplus + _Eij[1][1] * neghalfplus;
     552       25200 :       const Real cmin6 = _Eij[0][2] * halfplus + _Eij[0][0] * neghalfplus;
     553       25200 :       const Real cmin8 = _Eij[0][2] * halfplus + _Eij[0][1] * neghalfplus;
     554             :       // Substituting these into f6 = 0 yields
     555             :       // 0 = f6_trial - ga6 (0.5(cmax6 - cmin6) + 0.5(cmax6 + cmin6)sinphi) - ga8 (0.5(cmax6 -
     556             :       // cmin6) + 0.5(cmax6 + cmin6)sinphi) = f6_trial - ga6 c6 - ga8 c8, where
     557       25200 :       const Real c6 = 0.5 * (cmax6 - cmin6) + 0.5 * (cmax6 + cmin6) * sinphi;
     558       25200 :       const Real c8 = 0.5 * (cmax8 - cmin8) + 0.5 * (cmax8 + cmin8) * sinphi;
     559             :       // It isn't too hard to check that the other equation is
     560             :       // 0 = f8_trial - ga6 c8 - ga8 c6
     561             :       // These equations may be inverted to yield
     562       25200 :       if (c6 != c8)
     563             :       {
     564             :         const Real f6_trial = trial_mc_yf;
     565       25200 :         const Real f8_trial = 0.5 * (trial_stress_params[2] - trial_stress_params[1]) +
     566       25200 :                               0.5 * (trial_stress_params[2] + trial_stress_params[1]) * sinphi -
     567       25200 :                               cohcos;
     568       25200 :         const Real descr = Utility::pow<2>(c6) - Utility::pow<2>(c8);
     569       25200 :         Real ga6 = (c6 * f6_trial - c8 * f8_trial) / descr;
     570       25200 :         Real ga8 = (-c8 * f6_trial + c6 * f8_trial) / descr;
     571             :         // and finally
     572       25200 :         stress_params[2] = trial_stress_params[2] - ga6 * cmax6 - ga8 * cmax8;
     573       25200 :         stress_params[0] = trial_stress_params[0] - ga6 * cmin6 - ga8 * cmin8;
     574       25200 :         stress_params[1] = stress_params[0] + 0.5 * _shifter;
     575             : 
     576       25200 :         Real f7 = 0.5 * (stress_params[1] - stress_params[0]) +
     577       25200 :                   0.5 * (stress_params[1] + stress_params[0]) * sinphi - cohcos;
     578             : 
     579       25200 :         if (mc_tip_possible && f7 > _f_tol)
     580             :         {
     581           0 :           stress_params[2] = cohcos / sinphi;
     582           0 :           stress_params[1] = stress_params[2] - 0.5 * _shifter;
     583           0 :           stress_params[0] = stress_params[2] - _shifter;
     584             :           f7 = 0.0;
     585             :           ga6 = 1.0;
     586             :           ga8 = 1.0;
     587             :         }
     588             : 
     589       25200 :         const Real new_tensile_yf = stress_params[2] - ts;
     590       25200 :         const Real new_compressive_yf = -stress_params[0] - cs;
     591             : 
     592       25200 :         if (f7 <= _f_tol && new_tensile_yf <= _f_tol && new_compressive_yf <= _f_tol &&
     593       10416 :             ga6 >= 0.0 && ga8 >= 0.0)
     594             :         {
     595       10208 :           gaE = ((trial_stress_params[2] - trial_stress_params[0]) -
     596       10208 :                  (stress_params[2] - stress_params[0])) /
     597       10208 :                 (_Eij[2][2] - _Eij[0][2]) * _Eij[2][2];
     598             :           found_solution = true;
     599             :         }
     600             :       }
     601             :     }
     602      107632 :     if (!found_solution && !mc_impossible && tensile_possible && trial_tensile_yf > _f_tol)
     603             :     {
     604             :       // Return to the line where yf[0] = 0 = yf[6].
     605             :       // To return to this line, we need to solve f0 = 0 = f6 and
     606             :       // the three flow-direction equations.  In the flow-direction equations
     607             :       // there are two plasticity multipliers, which i call ga0 and ga6
     608             :       // corresponding to the amounts of strain normal to the f0 and f6
     609             :       // directions, respectively.
     610             :       // So:
     611             :       // Smax = Smax^trial - ga6 Emax,a dg6/dSa - ga0 Emax,a dg0/dSa
     612             :       //      = Smax^trial - ga6 cmax6 - ga0 cmax0  (with cmax6 and cmax0 evaluated below)
     613             :       // Smid = Smid^trial - ga6 Emid,a dg6/dSa - ga0 Emid,a dg0/dSa
     614             :       //      = Smid^trial - ga6 cmid6 - ga0 cmid0
     615             :       // Smin = Smin^trial - ga6 Emin,a dg6/dSa - ga0 Emin,a dg0/dSa
     616             :       //      = Smin^trial - ga6 cmin6 - ga0 cmin0
     617       11728 :       const Real cmax6 = _Eij[2][2] * halfplus + _Eij[2][0] * neghalfplus;
     618             :       const Real cmax0 = _Eij[2][2];
     619       11728 :       const Real cmid6 = _Eij[1][2] * halfplus + _Eij[1][0] * neghalfplus;
     620             :       const Real cmid0 = _Eij[1][2];
     621       11728 :       const Real cmin6 = _Eij[0][2] * halfplus + _Eij[0][0] * neghalfplus;
     622             :       const Real cmin0 = _Eij[0][2];
     623             :       // Substituting these into f6 = 0 yields
     624             :       // 0 = f6_trial - ga6 (0.5(cmax6 - cmin6) + 0.5(cmax6 + cmin6)sinphi) - ga0 (0.5(cmax0 -
     625             :       // cmin0) + 0.5(cmax0 + cmin0)sinphi) = f6_trial - ga6 c6 - ga0 c0, where
     626       11728 :       const Real c6 = 0.5 * (cmax6 - cmin6) + 0.5 * (cmax6 + cmin6) * sinphi;
     627       11728 :       const Real c0 = 0.5 * (cmax0 - cmin0) + 0.5 * (cmax0 + cmin0) * sinphi;
     628             :       // Substituting these into f0 = 0 yields
     629             :       // 0 = f0_trial - ga6 cmax6 - ga0 cmax0
     630             :       // These equations may be inverted to yield the following
     631       11728 :       const Real descr = c0 * cmax6 - c6 * cmax0;
     632       11728 :       if (descr != 0.0)
     633             :       {
     634       11728 :         const Real ga0 = (-c6 * trial_tensile_yf + cmax6 * trial_mc_yf) / descr;
     635       11728 :         const Real ga6 = (c0 * trial_tensile_yf - cmax0 * trial_mc_yf) / descr;
     636       11728 :         stress_params[2] = trial_stress_params[2] - ga6 * cmax6 - ga0 * cmax0;
     637       11728 :         stress_params[1] = trial_stress_params[1] - ga6 * cmid6 - ga0 * cmid0;
     638       11728 :         stress_params[0] = trial_stress_params[0] - ga6 * cmin6 - ga0 * cmin0;
     639             : 
     640             :         // enforce physicality (go to corners if necessary)
     641       11728 :         stress_params[0] =
     642       11728 :             std::min(stress_params[0],
     643       11728 :                      stress_params[2] - _shifter); // account for poor choice from user
     644             :         // if goto_corner then the max(min()) in the subsequent line will force the solution to lie
     645             :         // at the corner where max = mid = tensile.  This means the signs of ga0 and ga6 become
     646             :         // irrelevant in the check below
     647       11728 :         const bool goto_corner = (stress_params[1] >= stress_params[2] - 0.5 * _shifter);
     648       23456 :         stress_params[1] = std::max(std::min(stress_params[1], stress_params[2] - 0.5 * _shifter),
     649       11728 :                                     stress_params[0] + 0.5 * _shifter);
     650             : 
     651       11728 :         const Real new_compressive_yf = -stress_params[0] - cs;
     652       11728 :         if (new_compressive_yf <= _f_tol &&
     653       11728 :             (goto_corner || (ga0 >= 0.0 && ga6 >= 0.0))) // enforce ga>=0 unless going to a corner
     654             :         {
     655        4496 :           gaE = ((trial_stress_params[2] - trial_stress_params[0]) -
     656        4496 :                  (stress_params[2] - stress_params[0])) /
     657        4496 :                     (_Eij[2][2] - _Eij[0][2]) * _Eij[2][2] +
     658        4496 :                 (trial_stress_params[2] - stress_params[2]);
     659             :           found_solution = true;
     660             :         }
     661             :       }
     662             :     }
     663      107632 :     if (!found_solution && !mc_impossible)
     664             :     {
     665             :       // Return to the line where yf[3] = 0 = yf[6].
     666             :       // To return to this line, we need to solve f3 = 0 = f6 and
     667             :       // the three flow-direction equations.  In the flow-direction equations
     668             :       // there are two plasticity multipliers, which i call ga3 and ga6
     669             :       // corresponding to the amounts of strain normal to the f3 and f6
     670             :       // directions, respectively.
     671             :       // So:
     672             :       // Smax = Smax^trial - ga6 Emax,a dg6/dSa - ga3 Emax,a dg3/dSa
     673             :       //      = Smax^trial - ga6 cmax6 - ga3 cmax3  (with cmax6 and cmax3 evaluated below)
     674             :       // Smid = Smid^trial - ga6 Emid,a dg6/dSa - ga3 Emid,a dg3/dSa
     675             :       //      = Smid^trial - ga6 cmid6 - ga3 cmid3
     676             :       // Smin = Smin^trial - ga6 Emin,a dg6/dSa - ga3 Emin,a dg3/dSa
     677             :       //      = Smin^trial - ga6 cmin6 - ga3 cmin3
     678       10496 :       const Real cmax6 = _Eij[2][2] * halfplus + _Eij[2][0] * neghalfplus;
     679       10496 :       const Real cmax3 = -_Eij[2][0];
     680       10496 :       const Real cmid6 = _Eij[1][2] * halfplus + _Eij[1][0] * neghalfplus;
     681       10496 :       const Real cmid3 = -_Eij[1][0];
     682       10496 :       const Real cmin6 = _Eij[0][2] * halfplus + _Eij[0][0] * neghalfplus;
     683       10496 :       const Real cmin3 = -_Eij[0][0];
     684             :       // Substituting these into f6 = 0 yields
     685             :       // 0 = f6_trial - ga6 (0.5(cmax6 - cmin6) + 0.5(cmax6 + cmin6)sinphi) - ga3 (0.5(cmax3 -
     686             :       // cmin3) + 0.5(cmax3 + cmin3)sinphi) = f6_trial - ga6 c6 - ga3 c3, where
     687       10496 :       const Real c6 = 0.5 * (cmax6 - cmin6) + 0.5 * (cmax6 + cmin6) * sinphi;
     688       10496 :       const Real c3 = 0.5 * (cmax3 - cmin3) + 0.5 * (cmax3 + cmin3) * sinphi;
     689             :       // Substituting these into f3 = 0 yields
     690             :       // 0 = - f3_trial - ga6 cmin6 - ga3 cmin3
     691             :       // These equations may be inverted to yield the following
     692       10496 :       const Real descr = c3 * cmin6 - c6 * cmin3;
     693       10496 :       if (descr != 0.0)
     694             :       {
     695       10496 :         const Real ga3 = (c6 * trial_compressive_yf + cmin6 * trial_mc_yf) / descr;
     696       10496 :         const Real ga6 = (-c3 * trial_compressive_yf - cmin3 * trial_mc_yf) / descr;
     697       10496 :         stress_params[2] = trial_stress_params[2] - ga6 * cmax6 - ga3 * cmax3;
     698       10496 :         stress_params[1] = trial_stress_params[1] - ga6 * cmid6 - ga3 * cmid3;
     699       10496 :         stress_params[0] = trial_stress_params[0] - ga6 * cmin6 - ga3 * cmin3;
     700             : 
     701       10496 :         const Real new_tensile_yf = stress_params[2] - ts;
     702       10496 :         stress_params[2] =
     703       10496 :             std::max(stress_params[2],
     704       10496 :                      stress_params[0] + _shifter); // account for poor choice from user
     705       20992 :         stress_params[1] = std::max(std::min(stress_params[1], stress_params[2] - 0.5 * _shifter),
     706       10496 :                                     stress_params[0] + 0.5 * _shifter);
     707             : 
     708       10496 :         if (new_tensile_yf <= _f_tol && ga6 >= 0.0)
     709             :         {
     710       10496 :           gaE = ((trial_stress_params[2] - trial_stress_params[0]) -
     711       10496 :                  (stress_params[2] - stress_params[0])) /
     712       10496 :                     (_Eij[2][2] - _Eij[0][2]) * _Eij[2][2] +
     713       10496 :                 (-trial_stress_params[0] - stress_params[0]);
     714             : 
     715             :           found_solution = true;
     716             :         }
     717             :       }
     718             :     }
     719       97136 :     if (!found_solution)
     720             :     {
     721             :       // Cannot find an acceptable initialisation
     722           0 :       for (unsigned i = 0; i < _num_sp; ++i)
     723           0 :         stress_params[i] = trial_stress_params[i];
     724           0 :       gaE = 0.0;
     725           0 :       mooseWarning("CappedMohrCoulombStressUpdate cannot initialize from max = ",
     726             :                    stress_params[2],
     727             :                    " mid = ",
     728             :                    stress_params[1],
     729             :                    " min = ",
     730             :                    stress_params[0]);
     731             :     }
     732             :   }
     733      107632 :   setIntnlValuesV(trial_stress_params, stress_params, intnl_old, intnl);
     734      107632 : }
     735             : 
     736             : void
     737      928512 : CappedMohrCoulombStressUpdate::setIntnlValuesV(const std::vector<Real> & trial_stress_params,
     738             :                                                const std::vector<Real> & current_stress_params,
     739             :                                                const std::vector<Real> & intnl_old,
     740             :                                                std::vector<Real> & intnl) const
     741             : {
     742             :   // intnl[0] = shear, intnl[1] = tensile
     743      928512 :   const Real smax = current_stress_params[2];     // largest eigenvalue
     744      928512 :   const Real smin = current_stress_params[0];     // smallest eigenvalue
     745      928512 :   const Real trial_smax = trial_stress_params[2]; // largest eigenvalue
     746      928512 :   const Real trial_smin = trial_stress_params[0]; // smallest eigenvalue
     747      928512 :   const Real ga_shear = ((trial_smax - trial_smin) - (smax - smin)) / (_Eij[2][2] - _Eij[0][2]);
     748      928512 :   intnl[0] = intnl_old[0] + ga_shear;
     749      928512 :   const Real sinpsi = std::sin(_psi.value(intnl[0]));
     750      928512 :   const Real prefactor = (_Eij[2][2] + _Eij[0][2]) * sinpsi;
     751      928512 :   const Real shear_correction = prefactor * ga_shear;
     752      928512 :   const Real ga_tensile = (1.0 - _poissons_ratio) *
     753      928512 :                           ((trial_smax + trial_smin) - (smax + smin) - shear_correction) /
     754      928512 :                           _Eij[2][2];
     755      928512 :   intnl[1] = intnl_old[1] + ga_tensile;
     756      928512 : }
     757             : 
     758             : void
     759      460164 : CappedMohrCoulombStressUpdate::setIntnlDerivativesV(const std::vector<Real> & trial_stress_params,
     760             :                                                     const std::vector<Real> & current_stress_params,
     761             :                                                     const std::vector<Real> & intnl,
     762             :                                                     std::vector<std::vector<Real>> & dintnl) const
     763             : {
     764             :   // intnl[0] = shear, intnl[1] = tensile
     765      460164 :   const Real smax = current_stress_params[2];     // largest eigenvalue
     766      460164 :   const Real smin = current_stress_params[0];     // smallest eigenvalue
     767      460164 :   const Real trial_smax = trial_stress_params[2]; // largest eigenvalue
     768      460164 :   const Real trial_smin = trial_stress_params[0]; // smallest eigenvalue
     769      460164 :   const Real ga_shear = ((trial_smax - trial_smin) - (smax - smin)) / (_Eij[2][2] - _Eij[0][2]);
     770             :   const std::vector<Real> dga_shear = {
     771      460164 :       1.0 / (_Eij[2][2] - _Eij[0][2]), 0.0, -1.0 / (_Eij[2][2] - _Eij[0][2])};
     772             :   // intnl[0] = intnl_old[0] + ga_shear;
     773     1840656 :   for (std::size_t i = 0; i < _num_sp; ++i)
     774     1380492 :     dintnl[0][i] = dga_shear[i];
     775             : 
     776      460164 :   const Real sinpsi = std::sin(_psi.value(intnl[0]));
     777      460164 :   const Real dsinpsi_di0 = _psi.derivative(intnl[0]) * std::cos(_psi.value(intnl[0]));
     778             : 
     779      460164 :   const Real prefactor = (_Eij[2][2] + _Eij[0][2]) * sinpsi;
     780      460164 :   const Real dprefactor_di0 = (_Eij[2][2] + _Eij[0][2]) * dsinpsi_di0;
     781             :   // const Real shear_correction = prefactor * ga_shear;
     782      460164 :   std::vector<Real> dshear_correction(_num_sp);
     783     1840656 :   for (std::size_t i = 0; i < _num_sp; ++i)
     784     1380492 :     dshear_correction[i] = prefactor * dga_shear[i] + dprefactor_di0 * dintnl[0][i] * ga_shear;
     785             :   // const Real ga_tensile = (1 - _poissons_ratio) * ((trial_smax + trial_smin) - (smax + smin) -
     786             :   // shear_correction) /
     787             :   // _Eij[2][2];
     788             :   // intnl[1] = intnl_old[1] + ga_tensile;
     789     1840656 :   for (std::size_t i = 0; i < _num_sp; ++i)
     790     1380492 :     dintnl[1][i] = -(1.0 - _poissons_ratio) * dshear_correction[i] / _Eij[2][2];
     791      460164 :   dintnl[1][2] += -(1.0 - _poissons_ratio) / _Eij[2][2];
     792      460164 :   dintnl[1][0] += -(1.0 - _poissons_ratio) / _Eij[2][2];
     793      460164 : }
     794             : 
     795             : void
     796         928 : CappedMohrCoulombStressUpdate::consistentTangentOperatorV(
     797             :     const RankTwoTensor & stress_trial,
     798             :     const std::vector<Real> & trial_stress_params,
     799             :     const RankTwoTensor & /*stress*/,
     800             :     const std::vector<Real> & stress_params,
     801             :     Real /*gaE*/,
     802             :     const yieldAndFlow & /*smoothed_q*/,
     803             :     const RankFourTensor & elasticity_tensor,
     804             :     bool compute_full_tangent_operator,
     805             :     const std::vector<std::vector<Real>> & dvar_dtrial,
     806             :     RankFourTensor & cto)
     807             : {
     808         928 :   cto = elasticity_tensor;
     809         928 :   if (!compute_full_tangent_operator)
     810           0 :     return;
     811             : 
     812             :   // dvar_dtrial has been computed already, so
     813             :   // d(stress)/d(trial_stress) = d(eigvecs * stress_params * eigvecs.transpose())/d(trial_stress)
     814             :   // eigvecs is a rotation matrix, rot(i, j) = e_j(i) = i^th component of j^th eigenvector
     815             :   // d(rot_ij)/d(stress_kl) = d(e_j(i))/d(stress_kl)
     816             :   // = sum_a 0.5 * e_a(i) * (e_a(k)e_j(l) + e_a(l)e_j(k)) / (la_j - la_a)
     817             :   // = sum_a 0.5 * rot(i,a) * (rot(k,a)rot(l,j) + rot(l,a)*rot(k,j)) / (la_j - la_a)
     818         928 :   RankFourTensor drot_dstress;
     819        3712 :   for (unsigned i = 0; i < _tensor_dimensionality; ++i)
     820       11136 :     for (unsigned j = 0; j < _tensor_dimensionality; ++j)
     821       33408 :       for (unsigned k = 0; k < _tensor_dimensionality; ++k)
     822      100224 :         for (unsigned l = 0; l < _tensor_dimensionality; ++l)
     823      300672 :           for (unsigned a = 0; a < _num_sp; ++a)
     824             :           {
     825      225504 :             if (trial_stress_params[a] == trial_stress_params[j])
     826       75168 :               continue;
     827      150336 :             drot_dstress(i, j, k, l) +=
     828      150336 :                 0.5 * _eigvecs(i, a) *
     829      150336 :                 (_eigvecs(k, a) * _eigvecs(l, j) + _eigvecs(l, a) * _eigvecs(k, j)) /
     830      150336 :                 (trial_stress_params[j] - trial_stress_params[a]);
     831             :           }
     832             : 
     833         928 :   const RankTwoTensor eT = _eigvecs.transpose();
     834             : 
     835         928 :   RankFourTensor dstress_dtrial;
     836        3712 :   for (unsigned i = 0; i < _tensor_dimensionality; ++i)
     837       11136 :     for (unsigned j = 0; j < _tensor_dimensionality; ++j)
     838       33408 :       for (unsigned k = 0; k < _tensor_dimensionality; ++k)
     839      100224 :         for (unsigned l = 0; l < _tensor_dimensionality; ++l)
     840      300672 :           for (unsigned a = 0; a < _num_sp; ++a)
     841      225504 :             dstress_dtrial(i, j, k, l) +=
     842      225504 :                 drot_dstress(i, a, k, l) * stress_params[a] * eT(a, j) +
     843      225504 :                 _eigvecs(i, a) * stress_params[a] * drot_dstress(j, a, k, l);
     844             : 
     845         928 :   const std::vector<RankTwoTensor> dsp_trial = dstress_param_dstress(stress_trial);
     846        3712 :   for (unsigned i = 0; i < _tensor_dimensionality; ++i)
     847       11136 :     for (unsigned j = 0; j < _tensor_dimensionality; ++j)
     848       33408 :       for (unsigned k = 0; k < _tensor_dimensionality; ++k)
     849      100224 :         for (unsigned l = 0; l < _tensor_dimensionality; ++l)
     850      300672 :           for (unsigned a = 0; a < _num_sp; ++a)
     851      902016 :             for (unsigned b = 0; b < _num_sp; ++b)
     852      676512 :               dstress_dtrial(i, j, k, l) +=
     853      676512 :                   _eigvecs(i, a) * dvar_dtrial[a][b] * dsp_trial[b](k, l) * eT(a, j);
     854             : 
     855         928 :   cto = dstress_dtrial * elasticity_tensor;
     856             : }

Generated by: LCOV version 1.14