RankFourTensor
TensorMechanicsPlasticTensileMulti::consistentTangentOperator(const RankTwoTensor & trial_stress, const RankTwoTensor & stress, const Real & intnl,
                                                       const RankFourTensor & E_ijkl, const std::vector<Real> & cumulative_pm) const
{
  if (!_use_custom_cto)
    return TensorMechanicsPlasticModel::consistentTangentOperator(trial_stress, stress, intnl, E_ijkl, cumulative_pm);

  mooseAssert(cumulative_pm.size() == 3, "TensorMechanicsPlasticTensileMulti size of cumulative_pm should be 3 but it is " << cumulative_pm.size());

  if (cumulative_pm[2] <= 0) // All cumulative_pm are non-positive, so this is admissible
    return E_ijkl;


  // Need the eigenvalues at the returned configuration
  std::vector<Real> eigvals;
  stress.symmetricEigenvalues(eigvals);

  // need to rotate to and from principal stress space
  // using the eigenvectors of the trial configuration
  // (not the returned configuration).
  std::vector<Real> trial_eigvals;
  RankTwoTensor trial_eigvecs;
  trial_stress.symmetricEigenvaluesEigenvectors(trial_eigvals, trial_eigvecs);

  // The returnMap will have returned to the Tip, Edge or
  // Plane.  The consistentTangentOperator describes the
  // change in stress for an arbitrary change in applied
  // strain.  I assume that the change in strain will not
  // change the type of return (Tip remains Tip, Edge remains
  // Edge, Plane remains Plane).
  // I assume isotropic elasticity.
  //
  // The consistent tangent operator is a little different
  // than cases where no rotation to principal stress space
  // is made during the returnMap.  Let S_ij be the stress
  // in original coordinates, and s_ij be the stress in the
  // principal stress coordinates, so that
  // s_ij = diag(eigvals[0], eigvals[1], eigvals[2])
  // We want dS_ij under an arbitrary change in strain (ep->ep+dep)
  // dS = S(ep+dep) - S(ep)
  //    = R(ep+dep) s(ep+dep) R(ep+dep)^T - R(ep) s(ep) R(ep)^T
  // Here R = the rotation to principal-stress space, ie
  // R_ij = eigvecs[i][j] = i^th component of j^th eigenvector
  // Expanding to first order in dep,
  // dS = R(ep) (s(ep+dep) - s(ep)) R(ep)^T
  //      + dR/dep s(ep) R^T + R(ep) s(ep) dR^T/dep
  // The first line is all that is usually calculated in the
  // consistent tangent operator calculation, and is called
  // cto below.
  // The second line involves changes in the eigenvectors, and
  // is called sec below.

  RankFourTensor cto;
  Real hard = dtensile_strength(intnl);
  Real la = E_ijkl(0,0,1,1);
  Real mu = 0.5*(E_ijkl(0,0,0,0) - la);

  if (cumulative_pm[1] <= 0)
  {
    // only cumulative_pm[2] is positive, so this is return to the Plane
    Real denom = hard + la + 2*mu;
    Real al = la*la/denom;
    Real be = la*(la + 2*mu)/denom;
    Real ga = hard*(la + 2*mu)/denom;
    std::vector<Real> comps(9);
    comps[0] = comps[4] = la + 2*mu - al;
    comps[1] = comps[3] = la - al;
    comps[2] = comps[5] = comps[6] = comps[7] = la - be;
    comps[8] = ga;
    cto.fillFromInputVector(comps, RankFourTensor::principal);
  }
  else if (cumulative_pm[0] <= 0)
  {
    // both cumulative_pm[2] and cumulative_pm[1] are positive, so Edge
    Real denom = 2*hard + 2*la + 2*mu;
    Real al = hard*2*la/denom;
    Real be = hard*(2*la + 2*mu)/denom;
    std::vector<Real> comps(9);
    comps[0] = la + 2*mu - 2*la*la/denom;
    comps[1] = comps[2] = al;
    comps[3] = comps[6] = al;
    comps[4] = comps[5] = comps[7] = comps[8] = be;
    cto.fillFromInputVector(comps, RankFourTensor::principal);
  }
  else
  {
    // all cumulative_pm are positive, so Tip
    Real denom = 3*hard + 3*la + 2*mu;
    std::vector<Real> comps(2);
    comps[0] = hard*(3*la + 2*mu)/denom;
    comps[1] = 0;
    cto.fillFromInputVector(comps, RankFourTensor::symmetric_isotropic);
  }

  cto.rotate(trial_eigvecs);


  // drdsig = change in eigenvectors under a small stress change
  // drdsig(i,j,m,n) = dR(i,j)/dS_mn
  // The formula below is fairly easily derived:
  // S R = R s, so taking the variation
  // dS R + S dR = dR s + R ds, and multiplying by R^T
  // R^T dS R + R^T S dR = R^T dR s + ds .... (eqn 1)
  // I demand that RR^T = 1 = R^T R, and also that
  // (R+dR)(R+dR)^T = 1 = (R+dT)^T (R+dR), which means
  // that dR = R*c, for some antisymmetric c, so Eqn1 reads
  // R^T dS R + s c = c s + ds
  // Grabbing the components of this gives ds/dS (already
  // in RankTwoTensor), and c, which is:
  //   dR_ik/dS_mn = drdsig(i, k, m, n) = trial_eigvecs(m, b)*trial_eigvecs(n, k)*trial_eigvecs(i, b)/(trial_eigvals[k] - trial_eigvals[b]);
  // (sum over b!=k).

  RankFourTensor drdsig;
  for (unsigned k = 0 ; k < 3 ; ++k)
    for (unsigned b = 0 ; b < 3 ; ++b)
    {
      if (b == k)
        continue;
      for (unsigned m = 0 ; m < 3 ; ++m)
        for (unsigned n = 0 ; n < 3 ; ++n)
          for (unsigned i = 0 ; i < 3 ; ++i)
            drdsig(i, k, m, n) += trial_eigvecs(m, b)*trial_eigvecs(n, k)*trial_eigvecs(i, b)/(trial_eigvals[k] - trial_eigvals[b]);
    }



  // With diagla = diag(eigvals[0], eigvals[1], digvals[2])
  // The following implements
  // ans(i, j, a, b) += (drdsig(i, k, m, n)*trial_eigvecs(j, l)*diagla(k, l) + trial_eigvecs(i, k)*drdsig(j, l, m, n)*diagla(k, l))*E_ijkl(m, n, a, b);
  // (sum over k, l, m and n)

  RankFourTensor ans;
  for (unsigned i = 0 ; i < 3 ; ++i)
    for (unsigned j = 0 ; j < 3 ; ++j)
      for (unsigned a = 0 ; a < 3 ; ++a)
        for (unsigned k = 0 ; k < 3 ; ++k)
          for (unsigned m = 0 ; m < 3 ; ++m)
            ans(i, j, a, a) += (drdsig(i, k, m, m)*trial_eigvecs(j, k) + trial_eigvecs(i, k)*drdsig(j, k, m, m))*eigvals[k]*la;  //E_ijkl(m, n, a, b) = la*(m==n)*(a==b);

  for (unsigned i = 0 ; i < 3 ; ++i)
    for (unsigned j = 0 ; j < 3 ; ++j)
      for (unsigned a = 0 ; a < 3 ; ++a)
        for (unsigned b = 0 ; b < 3 ; ++b)
          for (unsigned k = 0 ; k < 3 ; ++k)
          {
            ans(i, j, a, b) += (drdsig(i, k, a, b)*trial_eigvecs(j, k) + trial_eigvecs(i, k)*drdsig(j, k, a, b))*eigvals[k]*mu;  //E_ijkl(m, n, a, b) = mu*(m==a)*(n==b)
            ans(i, j, a, b) += (drdsig(i, k, b, a)*trial_eigvecs(j, k) + trial_eigvecs(i, k)*drdsig(j, k, b, a))*eigvals[k]*mu;  //E_ijkl(m, n, a, b) = mu*(m==b)*(n==a)
          }


  return cto + ans;

}
bool
TensorMechanicsPlasticTensileMulti::doReturnMap(const RankTwoTensor & trial_stress, const Real & intnl_old, const RankFourTensor & E_ijkl,
                                                Real /*ep_plastic_tolerance*/, RankTwoTensor & returned_stress, Real & returned_intnl,
                                                std::vector<Real> & dpm, RankTwoTensor & delta_dp, std::vector<Real> & yf,
                                                bool & trial_stress_inadmissible) const
{
  mooseAssert(dpm.size() == 3, "TensorMechanicsPlasticTensileMulti size of dpm should be 3 but it is " << dpm.size());

  std::vector<Real> eigvals;
  RankTwoTensor eigvecs;
  trial_stress.symmetricEigenvaluesEigenvectors(eigvals, eigvecs);
  eigvals[0] += _shift;
  eigvals[2] -= _shift;

  Real str = tensile_strength(intnl_old);

  yf.resize(3);
  yf[0] = eigvals[0] - str;
  yf[1] = eigvals[1] - str;
  yf[2] = eigvals[2] - str;

  if (yf[0] <= _f_tol && yf[1] <= _f_tol && yf[2] <= _f_tol)
  {
    // purely elastic (trial_stress, intnl_old)
    trial_stress_inadmissible = false;
    return true;
  }

  trial_stress_inadmissible = true;
  delta_dp.zero();
  returned_stress.zero();

  // In the following i often assume that E_ijkl is
  // for an isotropic situation.  This reduces FLOPS
  // substantially which is important since the returnMap
  // is potentially the most compute-intensive function
  // of a simulation.
  // In many comments i write the general expression, and
  // i hope that might guide future coders if they are
  // generalising to a non-istropic E_ijkl

  // n[alpha] = E_ijkl*r[alpha]_kl expressed in principal stress space
  // (alpha = 0, 1, 2, corresponding to the three surfaces)
  // Note that in principal stress space, the flow
  // directions are, expressed in 'vector' form,
  // r[0] = (1,0,0), r[1] = (0,1,0), r[2] = (0,0,1).
  // Similar for _n:
  // so _n[0] = E_ij00*r[0], _n[1] = E_ij11*r[1], _n[2] = E_ij22*r[2]
  // In the following I assume that the E_ijkl is
  // for an isotropic situation.
  // In the anisotropic situation, we couldn't express
  // the flow directions as vectors in the same principal
  // stress space as the stress: they'd be full rank-2 tensors
  std::vector<std::vector<Real> > n(3);
  for (unsigned i = 0 ; i < 3 ; ++i)
    n[i].resize(3);
  n[0][0] = E_ijkl(0,0,0,0);
  n[0][1] = E_ijkl(1,1,0,0);
  n[0][2] = E_ijkl(2,2,0,0);
  n[1][0] = E_ijkl(0,0,1,1);
  n[1][1] = E_ijkl(1,1,1,1);
  n[1][2] = E_ijkl(2,2,1,1);
  n[2][0] = E_ijkl(0,0,2,2);
  n[2][1] = E_ijkl(1,1,2,2);
  n[2][2] = E_ijkl(2,2,2,2);


  // With non-zero Poisson's ratio and hardening
  // it is not computationally cheap to know whether
  // the trial stress will return to the tip, edge,
  // or plane.  The following is correct for zero
  // Poisson's ratio and no hardening, and at least
  // gives a not-completely-stupid guess in the
  // more general case.
  // trial_order[0] = type of return to try first
  // trial_order[1] = type of return to try second
  // trial_order[2] = type of return to try third
  std::vector<int> trial_order(3);
  if (yf[0] > 0) // all the yield functions are positive, since eigvals are ordered eigvals[0] <= eigvals[1] <= eigvals[2]
  {
    trial_order[0] = tip;
    trial_order[1] = edge;
    trial_order[2] = plane;
  }
  else if (yf[1] > 0)  // two yield functions are positive
  {
    trial_order[0] = edge;
    trial_order[1] = tip;
    trial_order[2] = plane;
  }
  else
  {
    trial_order[0] = plane;
    trial_order[1] = edge;
    trial_order[2] = tip;
  }

  unsigned trial;
  bool nr_converged;
  for (trial = 0 ; trial < 3 ; ++trial)
  {
    switch (trial_order[trial])
    {
      case tip:
        nr_converged = returnTip(eigvals, n, dpm, returned_stress, intnl_old, 0);
        break;
      case edge:
        nr_converged = returnEdge(eigvals, n, dpm, returned_stress, intnl_old, 0);
        break;
      case plane:
        nr_converged = returnPlane(eigvals, n, dpm, returned_stress, intnl_old, 0);
        break;
    }
    str = tensile_strength(intnl_old + dpm[0] + dpm[1] + dpm[2]);
    if (nr_converged && KuhnTuckerOK(returned_stress, dpm, str))
      break;
  }

  if (trial == 3)
  {
    Moose::err << "Trial stress = \n";
    trial_stress.print(Moose::err);
    Moose::err << "Internal parameter = " << intnl_old << "\n";
    mooseError("TensorMechanicsPlasticTensileMulti: FAILURE!  You probably need to implement a line search\n");
    // failure - must place yield function values at trial stress into yf
    str = tensile_strength(intnl_old);
    yf[0] = eigvals[0] - str;
    yf[1] = eigvals[1] - str;
    yf[2] = eigvals[2] - str;
    return false;
  }

  // success

  returned_intnl = intnl_old;
  for (unsigned i = 0 ; i < 3 ; ++i)
  {
    yf[i] = returned_stress(i, i) - str;
    delta_dp(i, i) = dpm[i];
    returned_intnl += dpm[i];
  }
  returned_stress = eigvecs*returned_stress*(eigvecs.transpose());
  delta_dp = eigvecs*delta_dp*(eigvecs.transpose());
  return true;
}
void
TensorMechanicsPlasticTensileMulti::activeConstraints(const std::vector<Real> & f, const RankTwoTensor & stress, const Real & intnl, const RankFourTensor & Eijkl, std::vector<bool> & act, RankTwoTensor & returned_stress) const
{
  act.assign(3, false);

  if (f[0] <= _f_tol && f[1] <= _f_tol && f[2] <= _f_tol)
  {
    returned_stress = stress;
    return;
  }

  returned_stress = RankTwoTensor();

  std::vector<Real> eigvals;
  RankTwoTensor eigvecs;
  stress.symmetricEigenvaluesEigenvectors(eigvals, eigvecs);
  eigvals[0] += _shift;
  eigvals[2] -= _shift;

  Real str = tensile_strength(intnl);
  std::vector<Real> v(3);
  v[0] = eigvals[0] - str;
  v[1] = eigvals[1] - str;
  v[2] = eigvals[2] - str;

  // these are the normals to the 3 yield surfaces
  std::vector<std::vector<Real> > n(3);
  n[0].resize(3);
  n[0][0] = 1 ; n[0][1] = 0 ; n[0][2] = 0;
  n[1].resize(3);
  n[1][0] = 0 ; n[1][1] = 1 ; n[1][2] = 0;
  n[2].resize(3);
  n[2][0] = 0 ; n[2][1] = 0 ; n[2][2] = 1;

  // the flow directions are these n multiplied by Eijkl.
  // I re-use the name "n" for the flow directions
  // In the following I assume that the Eijkl is
  // for an isotropic situation.  This is the most
  // common when using TensileMulti, and remember
  // that the returned_stress need not be perfect
  // anyway.
  // I divide by E(0,0,0,0) so the n remain of order 1
  Real ratio = Eijkl(1,1,0,0)/Eijkl(0,0,0,0);
  n[0][1] = n[0][2] = ratio;
  n[1][0] = n[1][2] = ratio;
  n[2][0] = n[2][1] = ratio;


  // 111 (tip)
  // For tip-return to satisfy Kuhn-Tucker, we need
  // v = alpha*n[0] + beta*n[1] * gamma*n[2]
  // with alpha, beta, and gamma all being non-negative (they are
  // the plasticity multipliers)
  Real denom = triple(n[0], n[1], n[2]);
  if (triple(v, n[0], n[1])/denom >= 0 && triple(v, n[1], n[2])/denom >= 0 && triple(v, n[2], n[0])/denom >= 0)
  {
    act[0] = act[1] = act[2] = true;
    returned_stress(0, 0) = returned_stress(1, 1) = returned_stress(2, 2) = str;
    returned_stress = eigvecs*returned_stress*(eigvecs.transpose());
    return;
  }

  // 011 (edge)
  std::vector<Real> n1xn2(3);
  n1xn2[0] = n[1][1]*n[2][2] - n[1][2]*n[2][1];
  n1xn2[1] = n[1][2]*n[2][0] - n[1][0]*n[2][2];
  n1xn2[2] = n[1][0]*n[2][1] - n[1][1]*n[2][0];
  // work out the point to which we would return, "a".  It is defined by
  // f1 = 0 = f2, and that (p - a).(n1 x n2) = 0, where "p" is the
  // starting position (p = eigvals).
  // In the following a = (a0, str, str)
  Real pdotn1xn2 = dot(eigvals, n1xn2);
  Real a0 = (-str*n1xn2[1] - str*n1xn2[2] + pdotn1xn2)/n1xn2[0];
  // we need p - a = alpha*n1 + beta*n2, where alpha and beta are non-negative
  // for Kuhn-Tucker to be satisfied
  std::vector<Real> pminusa(3);
  pminusa[0] = eigvals[0] - a0;
  pminusa[1] = v[1];
  pminusa[2] = v[2];
  if ((pminusa[2] - pminusa[0])/(1.0 - ratio) >= 0 && (pminusa[1] - pminusa[0])/(1.0 - ratio) >= 0)
  {
    returned_stress(0, 0) = a0;
    returned_stress(1, 1) = str;
    returned_stress(2, 2) = str;
    returned_stress = eigvecs*returned_stress*(eigvecs.transpose());
    act[1] = act[2] = true;
    return;
  }

  // 001 (plane)
  // the returned point, "a", is defined by f2=0 and
  // a = p - alpha*n2
  Real alpha = (eigvals[2] - str)/n[2][2];
  act[2] = true;
  returned_stress(0, 0) = eigvals[0] - alpha*n[2][0];
  returned_stress(1, 1) = eigvals[1] - alpha*n[2][1];
  returned_stress(2, 2) = str;
  returned_stress = eigvecs*returned_stress*(eigvecs.transpose());
  return;
}
bool
TensorMechanicsPlasticMohrCoulombMulti::doReturnMap(const RankTwoTensor & trial_stress, Real intnl_old, const RankFourTensor & E_ijkl,
                                                    Real ep_plastic_tolerance, RankTwoTensor & returned_stress, Real & returned_intnl,
                                                    std::vector<Real> & dpm, RankTwoTensor & delta_dp, std::vector<Real> & yf,
                                                    bool & trial_stress_inadmissible) const
{
  mooseAssert(dpm.size() == 6, "TensorMechanicsPlasticMohrCoulombMulti size of dpm should be 6 but it is " << dpm.size());

  std::vector<Real> eigvals;
  RankTwoTensor eigvecs;
  trial_stress.symmetricEigenvaluesEigenvectors(eigvals, eigvecs);
  eigvals[0] += _shift;
  eigvals[2] -= _shift;

  Real sinphi = std::sin(phi(intnl_old));
  Real cosphi = std::cos(phi(intnl_old));
  Real coh = cohesion(intnl_old);
  Real cohcos = coh*cosphi;

  yieldFunctionEigvals(eigvals[0], eigvals[1], eigvals[2], sinphi, cohcos, yf);

  if (yf[0] <= _f_tol && yf[1] <= _f_tol && yf[2] <= _f_tol && yf[3] <= _f_tol && yf[4] <= _f_tol && yf[5] <= _f_tol)
  {
    // purely elastic (trial_stress, intnl_old)
    trial_stress_inadmissible = false;
    return true;
  }

  trial_stress_inadmissible = true;
  delta_dp.zero();
  returned_stress = RankTwoTensor();

  // these are the normals to the 6 yield surfaces, which are const because of the assumption of no psi hardening
  std::vector<RealVectorValue> norm(6);
  const Real sinpsi = std::sin(psi(intnl_old));
  const Real oneminus = 0.5*(1 - sinpsi);
  const Real oneplus = 0.5*(1 + sinpsi);
  norm[0](0) = oneplus; norm[0](1) = -oneminus; norm[0](2) = 0;
  norm[1](0) = -oneminus; norm[1](1) = oneplus; norm[1](2) = 0;
  norm[2](0) = oneplus; norm[2](1) = 0; norm[2](2) = -oneminus;
  norm[3](0) = -oneminus; norm[3](1) = 0; norm[3](2) = oneplus;
  norm[4](0) = 0; norm[4](1) = oneplus; norm[4](2) = -oneminus;
  norm[5](0) = 0; norm[5](1) = -oneminus; norm[5](2) = oneplus;

  // the flow directions are these norm multiplied by Eijkl.
  // I call the flow directions "n".
  // In the following I assume that the Eijkl is
  // for an isotropic situation.  Then I don't have to
  // rotate to the principal-stress frame, and i don't
  // have to worry about strange off-diagonal things
  std::vector<RealVectorValue> n(6);
  for (unsigned ys = 0; ys < 6; ++ys)
    for (unsigned i = 0; i < 3; ++i)
      for (unsigned j = 0; j < 3; ++j)
        n[ys](i) += E_ijkl(i,i,j,j)*norm[ys](j);
  const Real mag_E = E_ijkl(0, 0, 0, 0);

  // With non-zero Poisson's ratio and hardening
  // it is not computationally cheap to know whether
  // the trial stress will return to the tip, edge,
  // or plane.  The following at least
  // gives a not-completely-stupid guess
  // trial_order[0] = type of return to try first
  // trial_order[1] = type of return to try second
  // trial_order[2] = type of return to try third
  // trial_order[3] = type of return to try fourth
  // trial_order[4] = type of return to try fifth
  // In the following the "binary" stuff indicates the
  // deactive (0) and active (1) surfaces, eg
  // 110100 means that surfaces 0, 1 and 3 are active
  // and 2, 4 and 5 are deactive
  const unsigned int number_of_return_paths = 5;
  std::vector<int> trial_order(number_of_return_paths);
  if (yf[1] > _f_tol && yf[3] > _f_tol && yf[5] > _f_tol)
  {
    trial_order[0] = tip110100;
    trial_order[1] = edge010100;
    trial_order[2] = plane000100;
    trial_order[3] = edge000101;
    trial_order[4] = tip010101;
  }
  else if (yf[1] <= _f_tol && yf[3] > _f_tol && yf[5] > _f_tol)
  {
    trial_order[0] = edge000101;
    trial_order[1] = plane000100;
    trial_order[2] = tip110100;
    trial_order[3] = tip010101;
    trial_order[4] = edge010100;
  }
  else if (yf[1] <= _f_tol && yf[3] > _f_tol && yf[5] <= _f_tol)
  {
    trial_order[0] = plane000100;
    trial_order[1] = edge000101;
    trial_order[2] = edge010100;
    trial_order[3] = tip110100;
    trial_order[4] = tip010101;
  }
  else
  {
    trial_order[0] = edge010100;
    trial_order[1] = plane000100;
    trial_order[2] = edge000101;
    trial_order[3] = tip110100;
    trial_order[4] = tip010101;
  }

  unsigned trial;
  bool nr_converged = false;
  bool kt_success = false;
  std::vector<RealVectorValue> ntip(3);
  std::vector<Real> dpmtip(3);

  for (trial = 0; trial < number_of_return_paths; ++trial)
  {
    switch (trial_order[trial])
    {
      case tip110100:
        for (unsigned int i = 0; i < 3; ++i)
        {
          ntip[0](i) = n[0](i);
          ntip[1](i) = n[1](i);
          ntip[2](i) = n[3](i);
        }
        kt_success = returnTip(eigvals, ntip, dpmtip, returned_stress, intnl_old, sinphi, cohcos, 0, nr_converged, ep_plastic_tolerance, yf);
        if (nr_converged && kt_success)
        {
          dpm[0] = dpmtip[0];
          dpm[1] = dpmtip[1];
          dpm[3] = dpmtip[2];
          dpm[2] = dpm[4] = dpm[5] = 0;
        }
        break;

      case tip010101:
        for (unsigned int i = 0; i < 3; ++i)
        {
          ntip[0](i) = n[1](i);
          ntip[1](i) = n[3](i);
          ntip[2](i) = n[5](i);
        }
        kt_success = returnTip(eigvals, ntip, dpmtip, returned_stress, intnl_old, sinphi, cohcos, 0, nr_converged, ep_plastic_tolerance, yf);
        if (nr_converged && kt_success)
        {
          dpm[1] = dpmtip[0];
          dpm[3] = dpmtip[1];
          dpm[5] = dpmtip[2];
          dpm[0] = dpm[2] = dpm[4] = 0;
        }
        break;

      case edge000101:
        kt_success = returnEdge000101(eigvals, n, dpm, returned_stress, intnl_old, sinphi, cohcos, 0, mag_E, nr_converged, ep_plastic_tolerance, yf);
        break;

      case edge010100:
        kt_success = returnEdge010100(eigvals, n, dpm, returned_stress, intnl_old, sinphi, cohcos, 0, mag_E, nr_converged, ep_plastic_tolerance, yf);
        break;

      case plane000100:
        kt_success = returnPlane(eigvals, n, dpm, returned_stress, intnl_old, sinphi, cohcos, 0, nr_converged, ep_plastic_tolerance, yf);
        break;
    }

    if (nr_converged && kt_success)
      break;
  }

  if (trial == number_of_return_paths)
  {
    sinphi = std::sin(phi(intnl_old));
    cosphi = std::cos(phi(intnl_old));
    coh = cohesion(intnl_old);
    cohcos = coh*cosphi;
    yieldFunctionEigvals(eigvals[0], eigvals[1], eigvals[2], sinphi, cohcos, yf);
    Moose::err << "Trial stress = \n";
    trial_stress.print(Moose::err);
    Moose::err << "which has eigenvalues = " << eigvals[0] << " " << eigvals[1] << " " << eigvals[2] << "\n";
    Moose::err << "and yield functions = " << yf[0] << " " << yf[1] << " " << yf[2] << " " << yf[3] << " " << yf[4] << " " << yf[5] << "\n";
    Moose::err << "Internal parameter = " << intnl_old << "\n";
    mooseError("TensorMechanicsPlasticMohrCoulombMulti: FAILURE!  You probably need to implement a line search if your hardening is too severe, or you need to tune your tolerances (eg, yield_function_tolerance should be a little smaller than (young modulus)*ep_plastic_tolerance).\n");
    return false;
  }

  // success

  returned_intnl = intnl_old;
  for (unsigned i = 0; i < 6; ++i)
    returned_intnl += dpm[i];
  for (unsigned i = 0; i < 6; ++i)
    for (unsigned j = 0; j < 3; ++j)
      delta_dp(j, j) += dpm[i]*norm[i](j);
  returned_stress = eigvecs*returned_stress*(eigvecs.transpose());
  delta_dp = eigvecs*delta_dp*(eigvecs.transpose());
  return true;
}
Ejemplo n.º 5
0
void
ComputeFiniteStrain::computeQpIncrements(RankTwoTensor & total_strain_increment, RankTwoTensor & rotation_increment)
{
  switch (_decomposition_method)
  {
    case DecompMethod::TaylorExpansion:
    {
      // inverse of _Fhat
      RankTwoTensor invFhat(_Fhat[_qp].inverse());

      // A = I - _Fhat^-1
      RankTwoTensor A(RankTwoTensor::initIdentity);
      A -= invFhat;

      // Cinv - I = A A^T - A - A^T;
      RankTwoTensor Cinv_I = A * A.transpose() - A - A.transpose();

      // strain rate D from Taylor expansion, Chat = (-1/2(Chat^-1 - I) + 1/4*(Chat^-1 - I)^2 + ...
      total_strain_increment = -Cinv_I * 0.5 + Cinv_I * Cinv_I * 0.25;

      const Real a[3] = {
        invFhat(1, 2) - invFhat(2, 1),
        invFhat(2, 0) - invFhat(0, 2),
        invFhat(0, 1) - invFhat(1, 0)
      };

      Real q = (a[0] * a[0] + a[1] * a[1] + a[2] * a[2]) / 4.0;
      Real trFhatinv_1 = invFhat.trace() - 1.0;
      const Real p = trFhatinv_1 * trFhatinv_1 / 4.0;

      // cos theta_a
      const Real C1 = std::sqrt(p + 3.0 * std::pow(p, 2.0) * (1.0 - (p + q)) / std::pow(p + q, 2.0) - 2.0 * std::pow(p, 3.0) * (1.0 - (p + q)) / std::pow(p + q, 3.0));

      Real C2;
      if (q > 0.01)
        // (1-cos theta_a)/4q
        C2 = (1.0 - C1) / (4.0 * q);
      else
        //alternate form for small q
        C2 = 0.125 + q * 0.03125 * (std::pow(p, 2.0) - 12.0 * (p - 1.0)) / std::pow(p, 2.0)
              + std::pow(q, 2.0) * (p - 2.0) * (std::pow(p, 2.0) - 10.0 * p + 32.0) / std::pow(p, 3.0)
              + std::pow(q, 3.0) * (1104.0 - 992.0 * p + 376.0 * std::pow(p, 2.0) - 72.0 * std::pow(p, 3.0) + 5.0 * std::pow(p, 4.0)) / (512.0 * std::pow(p, 4.0));

      const Real C3 = 0.5 * std::sqrt((p * q * (3.0 - q) + std::pow(p, 3.0) + std::pow(q, 2.0)) / std::pow(p + q, 3.0)); //sin theta_a/(2 sqrt(q))

      // Calculate incremental rotation. Note that this value is the transpose of that from Rashid, 93, so we transpose it before storing
      RankTwoTensor R_incr;
      R_incr.addIa(C1);
      for (unsigned int i = 0; i < 3; ++i)
        for (unsigned int j = 0; j < 3; ++j)
          R_incr(i,j) += C2 * a[i] * a[j];

      R_incr(0,1) += C3 * a[2];
      R_incr(0,2) -= C3 * a[1];
      R_incr(1,0) -= C3 * a[2];
      R_incr(1,2) += C3 * a[0];
      R_incr(2,0) += C3 * a[1];
      R_incr(2,1) -= C3 * a[0];

      rotation_increment = R_incr.transpose();
      break;
    }

    case DecompMethod::EigenSolution:
    {
      std::vector<Real> e_value(3);
      RankTwoTensor e_vector, N1, N2, N3;

      RankTwoTensor Chat = _Fhat[_qp].transpose() * _Fhat[_qp];
      Chat.symmetricEigenvaluesEigenvectors(e_value, e_vector);

      const Real lambda1 = std::sqrt(e_value[0]);
      const Real lambda2 = std::sqrt(e_value[1]);
      const Real lambda3 = std::sqrt(e_value[2]);

      N1.vectorOuterProduct(e_vector.column(0), e_vector.column(0));
      N2.vectorOuterProduct(e_vector.column(1), e_vector.column(1));
      N3.vectorOuterProduct(e_vector.column(2), e_vector.column(2));

      RankTwoTensor Uhat =  N1 * lambda1 + N2 * lambda2 + N3 * lambda3;
      RankTwoTensor invUhat(Uhat.inverse());

      rotation_increment = _Fhat[_qp] * invUhat;

      total_strain_increment = N1 * std::log(lambda1) + N2 * std::log(lambda2) + N3 * std::log(lambda3);
      break;
    }

    default:
      mooseError("ComputeFiniteStrain Error: Pass valid decomposition type: TaylorExpansion or EigenSolution.");
  }
}