www.mooseframework.org
Public Member Functions | Static Public Member Functions | Public Attributes | Protected Member Functions | Protected Attributes | List of all members
TensorMechanicsPlasticWeakPlaneShear Class Reference

Rate-independent associative weak-plane tensile failure with hardening/softening. More...

#include <TensorMechanicsPlasticWeakPlaneShear.h>

Inheritance diagram for TensorMechanicsPlasticWeakPlaneShear:
[legend]

Public Member Functions

 TensorMechanicsPlasticWeakPlaneShear (const InputParameters &parameters)
 
virtual void activeConstraints (const std::vector< Real > &f, const RankTwoTensor &stress, Real intnl, const RankFourTensor &Eijkl, std::vector< bool > &act, RankTwoTensor &returned_stress) const override
 The active yield surfaces, given a vector of yield functions. More...
 
virtual std::string modelName () const override
 
void initialize ()
 
void execute ()
 
void finalize ()
 
virtual unsigned int numberSurfaces () const
 The number of yield surfaces for this plasticity model. More...
 
virtual void yieldFunctionV (const RankTwoTensor &stress, Real intnl, std::vector< Real > &f) const
 Calculates the yield functions. More...
 
virtual void dyieldFunction_dstressV (const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &df_dstress) const
 The derivative of yield functions with respect to stress. More...
 
virtual void dyieldFunction_dintnlV (const RankTwoTensor &stress, Real intnl, std::vector< Real > &df_dintnl) const
 The derivative of yield functions with respect to the internal parameter. More...
 
virtual void flowPotentialV (const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &r) const
 The flow potentials. More...
 
virtual void dflowPotential_dstressV (const RankTwoTensor &stress, Real intnl, std::vector< RankFourTensor > &dr_dstress) const
 The derivative of the flow potential with respect to stress. More...
 
virtual void dflowPotential_dintnlV (const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &dr_dintnl) const
 The derivative of the flow potential with respect to the internal parameter. More...
 
virtual void hardPotentialV (const RankTwoTensor &stress, Real intnl, std::vector< Real > &h) const
 The hardening potential. More...
 
virtual void dhardPotential_dstressV (const RankTwoTensor &stress, Real intnl, std::vector< RankTwoTensor > &dh_dstress) const
 The derivative of the hardening potential with respect to stress. More...
 
virtual void dhardPotential_dintnlV (const RankTwoTensor &stress, Real intnl, std::vector< Real > &dh_dintnl) const
 The derivative of the hardening potential with respect to the internal parameter. More...
 
virtual bool useCustomReturnMap () const
 Returns false. You will want to override this in your derived class if you write a custom returnMap function. More...
 
virtual bool useCustomCTO () const
 Returns false. You will want to override this in your derived class if you write a custom consistent tangent operator function. More...
 
virtual bool returnMap (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
 Performs a custom return-map. More...
 
virtual RankFourTensor consistentTangentOperator (const RankTwoTensor &trial_stress, Real intnl_old, const RankTwoTensor &stress, Real intnl, const RankFourTensor &E_ijkl, const std::vector< Real > &cumulative_pm) const
 Calculates a custom consistent tangent operator. More...
 
bool KuhnTuckerSingleSurface (Real yf, Real dpm, Real dpm_tol) const
 Returns true if the Kuhn-Tucker conditions for the single surface are satisfied. More...
 

Static Public Member Functions

static InputParameters validParams ()
 

Public Attributes

const Real _f_tol
 Tolerance on yield function. More...
 
const Real _ic_tol
 Tolerance on internal constraint. More...
 

Protected Member Functions

Real yieldFunction (const RankTwoTensor &stress, Real intnl) const override
 The following functions are what you should override when building single-plasticity models. More...
 
RankTwoTensor dyieldFunction_dstress (const RankTwoTensor &stress, Real intnl) const override
 The derivative of yield function with respect to stress. More...
 
Real dyieldFunction_dintnl (const RankTwoTensor &stress, Real intnl) const override
 The derivative of yield function with respect to the internal parameter. More...
 
RankTwoTensor flowPotential (const RankTwoTensor &stress, Real intnl) const override
 The flow potential. More...
 
RankFourTensor dflowPotential_dstress (const RankTwoTensor &stress, Real intnl) const override
 The derivative of the flow potential with respect to stress. More...
 
RankTwoTensor dflowPotential_dintnl (const RankTwoTensor &stress, Real intnl) const override
 The derivative of the flow potential with respect to the internal parameter. More...
 
RankTwoTensor df_dsig (const RankTwoTensor &stress, Real _tan_phi_or_psi) const
 Function that's used in dyieldFunction_dstress and flowPotential. More...
 
virtual Real smooth (const RankTwoTensor &stress) const
 returns the 'a' parameter - see doco for _tip_scheme More...
 
virtual Real dsmooth (const RankTwoTensor &stress) const
 returns the da/dstress(2,2) - see doco for _tip_scheme More...
 
virtual Real d2smooth (const RankTwoTensor &stress) const
 returns the d^2a/dstress(2,2)^2 - see doco for _tip_scheme More...
 
virtual Real cohesion (const Real internal_param) const
 cohesion as a function of internal parameter More...
 
virtual Real dcohesion (const Real internal_param) const
 d(cohesion)/d(internal_param) More...
 
virtual Real tan_phi (const Real internal_param) const
 tan_phi as a function of internal parameter More...
 
virtual Real dtan_phi (const Real internal_param) const
 d(tan_phi)/d(internal_param); More...
 
virtual Real tan_psi (const Real internal_param) const
 tan_psi as a function of internal parameter More...
 
virtual Real dtan_psi (const Real internal_param) const
 d(tan_psi)/d(internal_param); More...
 
virtual Real hardPotential (const RankTwoTensor &stress, Real intnl) const
 The hardening potential. More...
 
virtual RankTwoTensor dhardPotential_dstress (const RankTwoTensor &stress, Real intnl) const
 The derivative of the hardening potential with respect to stress. More...
 
virtual Real dhardPotential_dintnl (const RankTwoTensor &stress, Real intnl) const
 The derivative of the hardening potential with respect to the internal parameter. More...
 

Protected Attributes

const TensorMechanicsHardeningModel_cohesion
 Hardening model for cohesion. More...
 
const TensorMechanicsHardeningModel_tan_phi
 Hardening model for tan(phi) More...
 
const TensorMechanicsHardeningModel_tan_psi
 Hardening model for tan(psi) More...
 
MooseEnum _tip_scheme
 The yield function is modified to f = sqrt(s_xz^2 + s_yz^2 + a) + s_zz*_tan_phi - _cohesion where "a" depends on the tip_scheme. More...
 
Real _small_smoother2
 smoothing parameter for the cone's tip - see doco for _tip_scheme More...
 
Real _cap_start
 smoothing parameter dictating when the 'cap' will start - see doco for _tip_scheme More...
 
Real _cap_rate
 dictates how quickly the 'cap' degenerates to a hemisphere - see doco for _tip_scheme More...
 

Detailed Description

Rate-independent associative weak-plane tensile failure with hardening/softening.

The cone's tip is smoothed.

Definition at line 24 of file TensorMechanicsPlasticWeakPlaneShear.h.

Constructor & Destructor Documentation

◆ TensorMechanicsPlasticWeakPlaneShear()

TensorMechanicsPlasticWeakPlaneShear::TensorMechanicsPlasticWeakPlaneShear ( const InputParameters &  parameters)

Definition at line 62 of file TensorMechanicsPlasticWeakPlaneShear.C.

64  : TensorMechanicsPlasticModel(parameters),
65  _cohesion(getUserObject<TensorMechanicsHardeningModel>("cohesion")),
66  _tan_phi(getUserObject<TensorMechanicsHardeningModel>("tan_friction_angle")),
67  _tan_psi(getUserObject<TensorMechanicsHardeningModel>("tan_dilation_angle")),
68  _tip_scheme(getParam<MooseEnum>("tip_scheme")),
69  _small_smoother2(Utility::pow<2>(getParam<Real>("smoother"))),
70  _cap_start(getParam<Real>("cap_start")),
71  _cap_rate(getParam<Real>("cap_rate"))
72 {
73  // With arbitary UserObjects, it is impossible to check everything, and
74  // I think this is the best I can do
75  if (tan_phi(0) < 0 || tan_psi(0) < 0)
76  mooseError("Weak-Plane-Shear friction and dilation angles must lie in [0, Pi/2]");
77  if (tan_phi(0) < tan_psi(0))
78  mooseError(
79  "Weak-Plane-Shear friction angle must not be less than Weak-Plane-Shear dilation angle");
80  if (cohesion(0) < 0)
81  mooseError("Weak-Plane-Shear cohesion must not be negative");
82 }

Member Function Documentation

◆ activeConstraints()

void TensorMechanicsPlasticWeakPlaneShear::activeConstraints ( const std::vector< Real > &  f,
const RankTwoTensor stress,
Real  intnl,
const RankFourTensor Eijkl,
std::vector< bool > &  act,
RankTwoTensor returned_stress 
) const
overridevirtual

The active yield surfaces, given a vector of yield functions.

This is used by FiniteStrainMultiPlasticity to determine the initial set of active constraints at the trial (stress, intnl) configuration. It is up to you (the coder) to determine how accurate you want the returned_stress to be. Currently it is only used by FiniteStrainMultiPlasticity to estimate a good starting value for the Newton-Rahson procedure, so currently it may not need to be super perfect.

Parameters
fvalues of the yield functions
stressstress tensor
intnlinternal parameter
Eijklelasticity tensor (stress = Eijkl*strain)
[out]actact[i] = true if the i_th yield function is active
[out]returned_stressApproximate value of the returned stress

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 278 of file TensorMechanicsPlasticWeakPlaneShear.C.

284 {
285  act.assign(1, false);
286  returned_stress = stress;
287 
288  if (f[0] <= _f_tol)
289  return;
290 
291  // in the following i will derive returned_stress for the case smoother=0
292 
293  Real tanpsi = tan_psi(intnl);
294  Real tanphi = tan_phi(intnl);
295 
296  // norm is the normal to the yield surface
297  // with f having psi (dilation angle) instead of phi:
298  // norm(0) = df/dsig(2,0) = df/dsig(0,2)
299  // norm(1) = df/dsig(2,1) = df/dsig(1,2)
300  // norm(2) = df/dsig(2,2)
301  std::vector<Real> norm(3, 0.0);
302  const Real tau = std::sqrt(Utility::pow<2>((stress(0, 2) + stress(2, 0)) / 2.0) +
303  Utility::pow<2>((stress(1, 2) + stress(2, 1)) / 2.0));
304  if (tau > 0.0)
305  {
306  norm[0] = 0.25 * (stress(0, 2) + stress(2, 0)) / tau;
307  norm[1] = 0.25 * (stress(1, 2) + stress(2, 1)) / tau;
308  }
309  else
310  {
311  returned_stress(2, 2) = cohesion(intnl) / tanphi;
312  act[0] = true;
313  return;
314  }
315  norm[2] = tanpsi;
316 
317  // to get the flow directions, we have to multiply norm by Eijkl.
318  // I assume that E(0,2,0,2) = E(1,2,1,2), and E(2,2,0,2) = 0 = E(0,2,1,2), etc
319  // with the usual symmetry. This makes finding the returned_stress
320  // much easier.
321  // returned_stress = stress - alpha*n
322  // where alpha is chosen so that f = 0
323  Real alpha = f[0] / (Eijkl(0, 2, 0, 2) + Eijkl(2, 2, 2, 2) * tanpsi * tanphi);
324 
325  if (1 - alpha * Eijkl(0, 2, 0, 2) / tau >= 0)
326  {
327  // returning to the "surface" of the cone
328  returned_stress(2, 2) = stress(2, 2) - alpha * Eijkl(2, 2, 2, 2) * norm[2];
329  returned_stress(0, 2) = returned_stress(2, 0) =
330  stress(0, 2) - alpha * 2.0 * Eijkl(0, 2, 0, 2) * norm[0];
331  returned_stress(1, 2) = returned_stress(2, 1) =
332  stress(1, 2) - alpha * 2.0 * Eijkl(1, 2, 1, 2) * norm[1];
333  }
334  else
335  {
336  // returning to the "tip" of the cone
337  returned_stress(2, 2) = cohesion(intnl) / tanphi;
338  returned_stress(0, 2) = returned_stress(2, 0) = returned_stress(1, 2) = returned_stress(2, 1) =
339  0.0;
340  }
341  returned_stress(0, 0) =
342  stress(0, 0) - Eijkl(0, 0, 2, 2) * (stress(2, 2) - returned_stress(2, 2)) / Eijkl(2, 2, 2, 2);
343  returned_stress(1, 1) =
344  stress(1, 1) - Eijkl(1, 1, 2, 2) * (stress(2, 2) - returned_stress(2, 2)) / Eijkl(2, 2, 2, 2);
345 
346  act[0] = true;
347 }

◆ cohesion()

Real TensorMechanicsPlasticWeakPlaneShear::cohesion ( const Real  internal_param) const
protectedvirtual

cohesion as a function of internal parameter

Definition at line 185 of file TensorMechanicsPlasticWeakPlaneShear.C.

186 {
187  return _cohesion.value(internal_param);
188 }

Referenced by activeConstraints(), TensorMechanicsPlasticWeakPlaneShear(), and yieldFunction().

◆ consistentTangentOperator()

RankFourTensor TensorMechanicsPlasticModel::consistentTangentOperator ( const RankTwoTensor trial_stress,
Real  intnl_old,
const RankTwoTensor stress,
Real  intnl,
const RankFourTensor E_ijkl,
const std::vector< Real > &  cumulative_pm 
) const
virtualinherited

Calculates a custom consistent tangent operator.

You may choose to over-ride this in your derived TensorMechanicsPlasticXXXX class.

(Note, if you over-ride returnMap, you will probably want to override consistentTangentOpertor too, otherwise it will default to E_ijkl.)

Parameters
stress_oldtrial stress before returning
intnl_oldinternal parameter before returning
stresscurrent returned stress state
intnlinternal parameter
E_ijklelasticity tensor
cumulative_pmthe cumulative plastic multipliers
Returns
the consistent tangent operator: E_ijkl if not over-ridden

Reimplemented in TensorMechanicsPlasticTensileMulti, TensorMechanicsPlasticDruckerPragerHyperbolic, TensorMechanicsPlasticMeanCapTC, and TensorMechanicsPlasticJ2.

Definition at line 254 of file TensorMechanicsPlasticModel.C.

261 {
262  return E_ijkl;
263 }

Referenced by TensorMechanicsPlasticJ2::consistentTangentOperator(), TensorMechanicsPlasticDruckerPragerHyperbolic::consistentTangentOperator(), TensorMechanicsPlasticMeanCapTC::consistentTangentOperator(), and TensorMechanicsPlasticTensileMulti::consistentTangentOperator().

◆ d2smooth()

Real TensorMechanicsPlasticWeakPlaneShear::d2smooth ( const RankTwoTensor stress) const
protectedvirtual

returns the d^2a/dstress(2,2)^2 - see doco for _tip_scheme

Definition at line 256 of file TensorMechanicsPlasticWeakPlaneShear.C.

257 {
258  Real d2smoother2 = 0;
259  if (_tip_scheme == "cap")
260  {
261  Real x = stress(2, 2) - _cap_start;
262  Real p = 0;
263  Real dp_dx = 0;
264  Real d2p_dx2 = 0;
265  if (x > 0)
266  {
267  const Real exp_cap_rate_x = std::exp(-_cap_rate * x);
268  p = x * (1.0 - exp_cap_rate_x);
269  dp_dx = (1.0 - exp_cap_rate_x) + x * _cap_rate * exp_cap_rate_x;
270  d2p_dx2 = 2.0 * _cap_rate * exp_cap_rate_x - x * Utility::pow<2>(_cap_rate) * exp_cap_rate_x;
271  }
272  d2smoother2 += 2.0 * Utility::pow<2>(dp_dx) + 2.0 * p * d2p_dx2;
273  }
274  return d2smoother2;
275 }

Referenced by dflowPotential_dstress().

◆ dcohesion()

Real TensorMechanicsPlasticWeakPlaneShear::dcohesion ( const Real  internal_param) const
protectedvirtual

d(cohesion)/d(internal_param)

Definition at line 191 of file TensorMechanicsPlasticWeakPlaneShear.C.

192 {
193  return _cohesion.derivative(internal_param);
194 }

Referenced by dyieldFunction_dintnl().

◆ df_dsig()

RankTwoTensor TensorMechanicsPlasticWeakPlaneShear::df_dsig ( const RankTwoTensor stress,
Real  _tan_phi_or_psi 
) const
protected

Function that's used in dyieldFunction_dstress and flowPotential.

Definition at line 94 of file TensorMechanicsPlasticWeakPlaneShear.C.

96 {
97  RankTwoTensor deriv; // the constructor zeroes this
98 
99  Real tau = std::sqrt(Utility::pow<2>((stress(0, 2) + stress(2, 0)) / 2.0) +
100  Utility::pow<2>((stress(1, 2) + stress(2, 1)) / 2.0) + smooth(stress));
101  // note that i explicitly symmeterise in preparation for Cosserat
102  if (tau == 0.0)
103  {
104  // the derivative is not defined here, but i want to set it nonzero
105  // because otherwise the return direction might be too crazy
106  deriv(0, 2) = deriv(2, 0) = 0.5;
107  deriv(1, 2) = deriv(2, 1) = 0.5;
108  }
109  else
110  {
111  deriv(0, 2) = deriv(2, 0) = 0.25 * (stress(0, 2) + stress(2, 0)) / tau;
112  deriv(1, 2) = deriv(2, 1) = 0.25 * (stress(1, 2) + stress(2, 1)) / tau;
113  deriv(2, 2) = 0.5 * dsmooth(stress) / tau;
114  }
115  deriv(2, 2) += _tan_phi_or_psi;
116  return deriv;
117 }

Referenced by dyieldFunction_dstress(), and flowPotential().

◆ dflowPotential_dintnl()

RankTwoTensor TensorMechanicsPlasticWeakPlaneShear::dflowPotential_dintnl ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

The derivative of the flow potential with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
Returns
dr_dintnl(i, j) = dr(i, j)/dintnl

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 176 of file TensorMechanicsPlasticWeakPlaneShear.C.

178 {
179  RankTwoTensor dr_dintnl;
180  dr_dintnl(2, 2) = dtan_psi(intnl);
181  return dr_dintnl;
182 }

◆ dflowPotential_dintnlV()

void TensorMechanicsPlasticModel::dflowPotential_dintnlV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankTwoTensor > &  dr_dintnl 
) const
virtualinherited

The derivative of the flow potential with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
[out]dr_dintnldr_dintnl[alpha](i, j) = dr[alpha](i, j)/dintnl

Reimplemented in TensorMechanicsPlasticMohrCoulombMulti, and TensorMechanicsPlasticTensileMulti.

Definition at line 139 of file TensorMechanicsPlasticModel.C.

142 {
143  return dr_dintnl.assign(1, dflowPotential_dintnl(stress, intnl));
144 }

◆ dflowPotential_dstress()

RankFourTensor TensorMechanicsPlasticWeakPlaneShear::dflowPotential_dstress ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

The derivative of the flow potential with respect to stress.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
Returns
dr_dstress(i, j, k, l) = dr(i, j)/dstress(k, l)

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 140 of file TensorMechanicsPlasticWeakPlaneShear.C.

142 {
143  RankFourTensor dr_dstress;
144  Real tau = std::sqrt(Utility::pow<2>((stress(0, 2) + stress(2, 0)) / 2.0) +
145  Utility::pow<2>((stress(1, 2) + stress(2, 1)) / 2.0) + smooth(stress));
146  if (tau == 0.0)
147  return dr_dstress;
148 
149  // note that i explicitly symmeterise
150  RankTwoTensor dtau;
151  dtau(0, 2) = dtau(2, 0) = 0.25 * (stress(0, 2) + stress(2, 0)) / tau;
152  dtau(1, 2) = dtau(2, 1) = 0.25 * (stress(1, 2) + stress(2, 1)) / tau;
153  dtau(2, 2) = 0.5 * dsmooth(stress) / tau;
154 
155  for (unsigned i = 0; i < 3; ++i)
156  for (unsigned j = 0; j < 3; ++j)
157  for (unsigned k = 0; k < 3; ++k)
158  for (unsigned l = 0; l < 3; ++l)
159  dr_dstress(i, j, k, l) = -dtau(i, j) * dtau(k, l) / tau;
160 
161  // note that i explicitly symmeterise
162  dr_dstress(0, 2, 0, 2) += 0.25 / tau;
163  dr_dstress(0, 2, 2, 0) += 0.25 / tau;
164  dr_dstress(2, 0, 0, 2) += 0.25 / tau;
165  dr_dstress(2, 0, 2, 0) += 0.25 / tau;
166  dr_dstress(1, 2, 1, 2) += 0.25 / tau;
167  dr_dstress(1, 2, 2, 1) += 0.25 / tau;
168  dr_dstress(2, 1, 1, 2) += 0.25 / tau;
169  dr_dstress(2, 1, 2, 1) += 0.25 / tau;
170  dr_dstress(2, 2, 2, 2) += 0.5 * d2smooth(stress) / tau;
171 
172  return dr_dstress;
173 }

◆ dflowPotential_dstressV()

void TensorMechanicsPlasticModel::dflowPotential_dstressV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankFourTensor > &  dr_dstress 
) const
virtualinherited

The derivative of the flow potential with respect to stress.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
[out]dr_dstressdr_dstress[alpha](i, j, k, l) = dr[alpha](i, j)/dstress(k, l)

Reimplemented in TensorMechanicsPlasticMohrCoulombMulti, and TensorMechanicsPlasticTensileMulti.

Definition at line 125 of file TensorMechanicsPlasticModel.C.

128 {
129  return dr_dstress.assign(1, dflowPotential_dstress(stress, intnl));
130 }

◆ dhardPotential_dintnl()

Real TensorMechanicsPlasticModel::dhardPotential_dintnl ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The derivative of the hardening potential with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlinternal parameter
Returns
the derivative

Reimplemented in TensorMechanicsPlasticMeanCapTC.

Definition at line 174 of file TensorMechanicsPlasticModel.C.

176 {
177  return 0.0;
178 }

Referenced by TensorMechanicsPlasticModel::dhardPotential_dintnlV().

◆ dhardPotential_dintnlV()

void TensorMechanicsPlasticModel::dhardPotential_dintnlV ( const RankTwoTensor stress,
Real  intnl,
std::vector< Real > &  dh_dintnl 
) const
virtualinherited

The derivative of the hardening potential with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlinternal parameter
[out]dh_dintnldh_dintnl[alpha] = dh[alpha]/dintnl

Definition at line 180 of file TensorMechanicsPlasticModel.C.

183 {
184  dh_dintnl.resize(numberSurfaces(), dhardPotential_dintnl(stress, intnl));
185 }

◆ dhardPotential_dstress()

RankTwoTensor TensorMechanicsPlasticModel::dhardPotential_dstress ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The derivative of the hardening potential with respect to stress.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlinternal parameter
Returns
dh_dstress(i, j) = dh/dstress(i, j)

Reimplemented in TensorMechanicsPlasticMeanCapTC.

Definition at line 160 of file TensorMechanicsPlasticModel.C.

162 {
163  return RankTwoTensor();
164 }

Referenced by TensorMechanicsPlasticModel::dhardPotential_dstressV().

◆ dhardPotential_dstressV()

void TensorMechanicsPlasticModel::dhardPotential_dstressV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankTwoTensor > &  dh_dstress 
) const
virtualinherited

The derivative of the hardening potential with respect to stress.

Parameters
stressthe stress at which to calculate the hardening potentials
intnlinternal parameter
[out]dh_dstressdh_dstress[alpha](i, j) = dh[alpha]/dstress(i, j)

Definition at line 166 of file TensorMechanicsPlasticModel.C.

169 {
170  dh_dstress.assign(numberSurfaces(), dhardPotential_dstress(stress, intnl));
171 }

◆ dsmooth()

Real TensorMechanicsPlasticWeakPlaneShear::dsmooth ( const RankTwoTensor stress) const
protectedvirtual

returns the da/dstress(2,2) - see doco for _tip_scheme

Definition at line 236 of file TensorMechanicsPlasticWeakPlaneShear.C.

237 {
238  Real dsmoother2 = 0;
239  if (_tip_scheme == "cap")
240  {
241  Real x = stress(2, 2) - _cap_start;
242  Real p = 0;
243  Real dp_dx = 0;
244  if (x > 0)
245  {
246  const Real exp_cap_rate_x = std::exp(-_cap_rate * x);
247  p = x * (1 - exp_cap_rate_x);
248  dp_dx = (1 - exp_cap_rate_x) + x * _cap_rate * exp_cap_rate_x;
249  }
250  dsmoother2 += 2 * p * dp_dx;
251  }
252  return dsmoother2;
253 }

Referenced by df_dsig(), and dflowPotential_dstress().

◆ dtan_phi()

Real TensorMechanicsPlasticWeakPlaneShear::dtan_phi ( const Real  internal_param) const
protectedvirtual

d(tan_phi)/d(internal_param);

Definition at line 203 of file TensorMechanicsPlasticWeakPlaneShear.C.

204 {
205  return _tan_phi.derivative(internal_param);
206 }

Referenced by dyieldFunction_dintnl().

◆ dtan_psi()

Real TensorMechanicsPlasticWeakPlaneShear::dtan_psi ( const Real  internal_param) const
protectedvirtual

d(tan_psi)/d(internal_param);

Definition at line 215 of file TensorMechanicsPlasticWeakPlaneShear.C.

216 {
217  return _tan_psi.derivative(internal_param);
218 }

Referenced by dflowPotential_dintnl().

◆ dyieldFunction_dintnl()

Real TensorMechanicsPlasticWeakPlaneShear::dyieldFunction_dintnl ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

The derivative of yield function with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
Returns
the derivative

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 127 of file TensorMechanicsPlasticWeakPlaneShear.C.

129 {
130  return stress(2, 2) * dtan_phi(intnl) - dcohesion(intnl);
131 }

◆ dyieldFunction_dintnlV()

void TensorMechanicsPlasticModel::dyieldFunction_dintnlV ( const RankTwoTensor stress,
Real  intnl,
std::vector< Real > &  df_dintnl 
) const
virtualinherited

The derivative of yield functions with respect to the internal parameter.

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
[out]df_dintnldf_dintnl[alpha] = df[alpha]/dintnl

Reimplemented in TensorMechanicsPlasticMohrCoulombMulti, and TensorMechanicsPlasticTensileMulti.

Definition at line 98 of file TensorMechanicsPlasticModel.C.

101 {
102  return df_dintnl.assign(1, dyieldFunction_dintnl(stress, intnl));
103 }

◆ dyieldFunction_dstress()

RankTwoTensor TensorMechanicsPlasticWeakPlaneShear::dyieldFunction_dstress ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

The derivative of yield function with respect to stress.

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
Returns
df_dstress(i, j) = dyieldFunction/dstress(i, j)

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 120 of file TensorMechanicsPlasticWeakPlaneShear.C.

122 {
123  return df_dsig(stress, tan_phi(intnl));
124 }

◆ dyieldFunction_dstressV()

void TensorMechanicsPlasticModel::dyieldFunction_dstressV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankTwoTensor > &  df_dstress 
) const
virtualinherited

The derivative of yield functions with respect to stress.

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
[out]df_dstressdf_dstress[alpha](i, j) = dyieldFunction[alpha]/dstress(i, j)

Reimplemented in TensorMechanicsPlasticMohrCoulombMulti, and TensorMechanicsPlasticTensileMulti.

Definition at line 84 of file TensorMechanicsPlasticModel.C.

87 {
88  df_dstress.assign(1, dyieldFunction_dstress(stress, intnl));
89 }

◆ execute()

void TensorMechanicsPlasticModel::execute ( )
inherited

Definition at line 47 of file TensorMechanicsPlasticModel.C.

48 {
49 }

◆ finalize()

void TensorMechanicsPlasticModel::finalize ( )
inherited

Definition at line 52 of file TensorMechanicsPlasticModel.C.

53 {
54 }

◆ flowPotential()

RankTwoTensor TensorMechanicsPlasticWeakPlaneShear::flowPotential ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

The flow potential.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
Returns
the flow potential

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 134 of file TensorMechanicsPlasticWeakPlaneShear.C.

135 {
136  return df_dsig(stress, tan_psi(intnl));
137 }

◆ flowPotentialV()

void TensorMechanicsPlasticModel::flowPotentialV ( const RankTwoTensor stress,
Real  intnl,
std::vector< RankTwoTensor > &  r 
) const
virtualinherited

The flow potentials.

Parameters
stressthe stress at which to calculate the flow potential
intnlinternal parameter
[out]rr[alpha] is the flow potential for the "alpha" yield function

Reimplemented in TensorMechanicsPlasticMohrCoulombMulti, and TensorMechanicsPlasticTensileMulti.

Definition at line 111 of file TensorMechanicsPlasticModel.C.

114 {
115  return r.assign(1, flowPotential(stress, intnl));
116 }

◆ hardPotential()

Real TensorMechanicsPlasticModel::hardPotential ( const RankTwoTensor stress,
Real  intnl 
) const
protectedvirtualinherited

The hardening potential.

Parameters
stressthe stress at which to calculate the hardening potential
intnlinternal parameter
Returns
the hardening potential

Reimplemented in TensorMechanicsPlasticMeanCapTC.

Definition at line 147 of file TensorMechanicsPlasticModel.C.

148 {
149  return -1.0;
150 }

Referenced by TensorMechanicsPlasticModel::hardPotentialV().

◆ hardPotentialV()

void TensorMechanicsPlasticModel::hardPotentialV ( const RankTwoTensor stress,
Real  intnl,
std::vector< Real > &  h 
) const
virtualinherited

The hardening potential.

Parameters
stressthe stress at which to calculate the hardening potential
intnlinternal parameter
[out]hh[alpha] is the hardening potential for the "alpha" yield function

Definition at line 152 of file TensorMechanicsPlasticModel.C.

155 {
156  h.assign(numberSurfaces(), hardPotential(stress, intnl));
157 }

◆ initialize()

void TensorMechanicsPlasticModel::initialize ( )
inherited

Definition at line 42 of file TensorMechanicsPlasticModel.C.

43 {
44 }

◆ KuhnTuckerSingleSurface()

bool TensorMechanicsPlasticModel::KuhnTuckerSingleSurface ( Real  yf,
Real  dpm,
Real  dpm_tol 
) const
inherited

Returns true if the Kuhn-Tucker conditions for the single surface are satisfied.

Parameters
yfYield function value
dpmplastic multiplier
dpm_toltolerance on plastic multiplier: viz dpm>-dpm_tol means "dpm is non-negative"

Definition at line 248 of file TensorMechanicsPlasticModel.C.

249 {
250  return (dpm == 0 && yf <= _f_tol) || (dpm > -dpm_tol && yf <= _f_tol && yf >= -_f_tol);
251 }

Referenced by TensorMechanicsPlasticMohrCoulombMulti::KuhnTuckerOK(), TensorMechanicsPlasticTensileMulti::KuhnTuckerOK(), and TensorMechanicsPlasticModel::returnMap().

◆ modelName()

std::string TensorMechanicsPlasticWeakPlaneShear::modelName ( ) const
overridevirtual

Implements TensorMechanicsPlasticModel.

Definition at line 350 of file TensorMechanicsPlasticWeakPlaneShear.C.

351 {
352  return "WeakPlaneShear";
353 }

◆ numberSurfaces()

unsigned TensorMechanicsPlasticModel::numberSurfaces ( ) const
virtualinherited

◆ returnMap()

bool TensorMechanicsPlasticModel::returnMap ( 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
virtualinherited

Performs a custom return-map.

You may choose to over-ride this in your derived TensorMechanicsPlasticXXXX class, and you may implement the return-map algorithm in any way that suits you. Eg, using a Newton-Raphson approach, or a radial-return, etc. This may also be used as a quick way of ascertaining whether (trial_stress, intnl_old) is in fact admissible.

For over-riding this function, please note the following.

(1) Denoting the return value of the function by "successful_return", the only possible output values should be: (A) trial_stress_inadmissible=false, successful_return=true. That is, (trial_stress, intnl_old) is in fact admissible (in the elastic domain). (B) trial_stress_inadmissible=true, successful_return=false. That is (trial_stress, intnl_old) is inadmissible (outside the yield surface), and you didn't return to the yield surface. (C) trial_stress_inadmissible=true, successful_return=true. That is (trial_stress, intnl_old) is inadmissible (outside the yield surface), but you did return to the yield surface. The default implementation only handles case (A) and (B): it does not attempt to do a return-map algorithm.

(2) you must correctly signal "successful_return" using the return value of this function. Don't assume the calling function will do Kuhn-Tucker checking and so forth!

(3) In cases (A) and (B) you needn't set returned_stress, returned_intnl, delta_dp, or dpm. This is for computational efficiency.

(4) In cases (A) and (B), you MUST place the yield function values at (trial_stress, intnl_old) into yf so the calling function can use this information optimally. You will have already calculated these yield function values, which can be quite expensive, and it's not very optimal for the calling function to have to re-calculate them.

(5) In case (C), you need to set: returned_stress (the returned value of stress) returned_intnl (the returned value of the internal variable) delta_dp (the change in plastic strain) dpm (the plastic multipliers needed to bring about the return) yf (yield function values at the returned configuration)

(Note, if you over-ride returnMap, you will probably want to override consistentTangentOpertor too, otherwise it will default to E_ijkl.)

Parameters
trial_stressThe trial stress
intnl_oldValue of the internal parameter
E_ijklElasticity tensor
ep_plastic_toleranceTolerance defined by the user for the plastic strain
[out]returned_stressIn case (C): lies on the yield surface after returning and produces the correct plastic strain (normality condition). Otherwise: not defined
[out]returned_intnlIn case (C): the value of the internal parameter after returning. Otherwise: not defined
[out]dpmIn case (C): the plastic multipliers needed to bring about the return. Otherwise: not defined
[out]delta_dpIn case (C): The change in plastic strain induced by the return process. Otherwise: not defined
[out]yfIn case (C): the yield function at (returned_stress, returned_intnl). Otherwise: the yield function at (trial_stress, intnl_old)
[out]trial_stress_inadmissibleShould be set to false if the trial_stress is admissible, and true if the trial_stress is inadmissible. This can be used by the calling prorgram
Returns
true if a successful return (or a return-map not needed), false if the trial_stress is inadmissible but the return process failed

Reimplemented in TensorMechanicsPlasticTensileMulti, TensorMechanicsPlasticMohrCoulombMulti, TensorMechanicsPlasticDruckerPragerHyperbolic, TensorMechanicsPlasticMeanCapTC, and TensorMechanicsPlasticJ2.

Definition at line 221 of file TensorMechanicsPlasticModel.C.

231 {
232  trial_stress_inadmissible = false;
233  yieldFunctionV(trial_stress, intnl_old, yf);
234 
235  for (unsigned sf = 0; sf < numberSurfaces(); ++sf)
236  if (yf[sf] > _f_tol)
237  trial_stress_inadmissible = true;
238 
239  // example of checking Kuhn-Tucker
240  std::vector<Real> dpm(numberSurfaces(), 0);
241  for (unsigned sf = 0; sf < numberSurfaces(); ++sf)
242  if (!KuhnTuckerSingleSurface(yf[sf], dpm[sf], 0))
243  return false;
244  return true;
245 }

Referenced by TensorMechanicsPlasticJ2::returnMap(), TensorMechanicsPlasticDruckerPragerHyperbolic::returnMap(), TensorMechanicsPlasticMeanCapTC::returnMap(), TensorMechanicsPlasticMohrCoulombMulti::returnMap(), and TensorMechanicsPlasticTensileMulti::returnMap().

◆ smooth()

Real TensorMechanicsPlasticWeakPlaneShear::smooth ( const RankTwoTensor stress) const
protectedvirtual

returns the 'a' parameter - see doco for _tip_scheme

Definition at line 221 of file TensorMechanicsPlasticWeakPlaneShear.C.

222 {
223  Real smoother2 = _small_smoother2;
224  if (_tip_scheme == "cap")
225  {
226  Real x = stress(2, 2) - _cap_start;
227  Real p = 0;
228  if (x > 0)
229  p = x * (1 - std::exp(-_cap_rate * x));
230  smoother2 += Utility::pow<2>(p);
231  }
232  return smoother2;
233 }

Referenced by df_dsig(), dflowPotential_dstress(), and yieldFunction().

◆ tan_phi()

Real TensorMechanicsPlasticWeakPlaneShear::tan_phi ( const Real  internal_param) const
protectedvirtual

tan_phi as a function of internal parameter

Definition at line 197 of file TensorMechanicsPlasticWeakPlaneShear.C.

198 {
199  return _tan_phi.value(internal_param);
200 }

Referenced by activeConstraints(), dyieldFunction_dstress(), TensorMechanicsPlasticWeakPlaneShear(), and yieldFunction().

◆ tan_psi()

Real TensorMechanicsPlasticWeakPlaneShear::tan_psi ( const Real  internal_param) const
protectedvirtual

tan_psi as a function of internal parameter

Definition at line 209 of file TensorMechanicsPlasticWeakPlaneShear.C.

210 {
211  return _tan_psi.value(internal_param);
212 }

Referenced by activeConstraints(), flowPotential(), and TensorMechanicsPlasticWeakPlaneShear().

◆ useCustomCTO()

bool TensorMechanicsPlasticModel::useCustomCTO ( ) const
virtualinherited

Returns false. You will want to override this in your derived class if you write a custom consistent tangent operator function.

Reimplemented in TensorMechanicsPlasticTensileMulti, TensorMechanicsPlasticMeanCapTC, TensorMechanicsPlasticDruckerPragerHyperbolic, and TensorMechanicsPlasticJ2.

Definition at line 215 of file TensorMechanicsPlasticModel.C.

216 {
217  return false;
218 }

◆ useCustomReturnMap()

bool TensorMechanicsPlasticModel::useCustomReturnMap ( ) const
virtualinherited

Returns false. You will want to override this in your derived class if you write a custom returnMap function.

Reimplemented in TensorMechanicsPlasticMohrCoulombMulti, TensorMechanicsPlasticTensileMulti, TensorMechanicsPlasticMeanCapTC, TensorMechanicsPlasticDruckerPragerHyperbolic, and TensorMechanicsPlasticJ2.

Definition at line 209 of file TensorMechanicsPlasticModel.C.

210 {
211  return false;
212 }

◆ validParams()

InputParameters TensorMechanicsPlasticWeakPlaneShear::validParams ( )
static

Definition at line 19 of file TensorMechanicsPlasticWeakPlaneShear.C.

20 {
21  InputParameters params = TensorMechanicsPlasticModel::validParams();
22  params.addRequiredParam<UserObjectName>(
23  "cohesion",
24  "A TensorMechanicsHardening UserObject that defines hardening of the cohesion. "
25  "Physically the cohesion should not be negative.");
26  params.addRequiredParam<UserObjectName>("tan_friction_angle",
27  "A TensorMechanicsHardening UserObject that defines "
28  "hardening of tan(friction angle). Physically the "
29  "friction angle should be between 0 and 90deg.");
30  params.addRequiredParam<UserObjectName>(
31  "tan_dilation_angle",
32  "A TensorMechanicsHardening UserObject that defines hardening of the "
33  "tan(dilation angle). Usually the dilation angle is not greater than "
34  "the friction angle, and it is between 0 and 90deg.");
35  MooseEnum tip_scheme("hyperbolic cap", "hyperbolic");
36  params.addParam<MooseEnum>(
37  "tip_scheme", tip_scheme, "Scheme by which the cone's tip will be smoothed.");
38  params.addRequiredRangeCheckedParam<Real>(
39  "smoother",
40  "smoother>=0",
41  "For the 'hyperbolic' tip_scheme, the cone vertex at shear-stress "
42  "= 0 will be smoothed by the given amount. For the 'cap' "
43  "tip_scheme, additional smoothing will occur. Typical value is "
44  "0.1*cohesion");
45  params.addParam<Real>(
46  "cap_start",
47  0.0,
48  "For the 'cap' tip_scheme, smoothing is performed in the stress_zz > cap_start region");
49  params.addRangeCheckedParam<Real>("cap_rate",
50  0.0,
51  "cap_rate>=0",
52  "For the 'cap' tip_scheme, this controls how quickly the cap "
53  "degenerates to a hemisphere: small values mean a slow "
54  "degeneration to a hemisphere (and zero means the 'cap' will "
55  "be totally inactive). Typical value is 1/cohesion");
56  params.addClassDescription("Non-associative finite-strain weak-plane shear perfect plasticity. "
57  "Here cohesion = 1, tan(phi) = 1 = tan(psi)");
58 
59  return params;
60 }

◆ yieldFunction()

Real TensorMechanicsPlasticWeakPlaneShear::yieldFunction ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

The following functions are what you should override when building single-plasticity models.

The yield function

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
Returns
the yield function

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 85 of file TensorMechanicsPlasticWeakPlaneShear.C.

86 {
87  // note that i explicitly symmeterise in preparation for Cosserat
88  return std::sqrt(Utility::pow<2>((stress(0, 2) + stress(2, 0)) / 2.0) +
89  Utility::pow<2>((stress(1, 2) + stress(2, 1)) / 2.0) + smooth(stress)) +
90  stress(2, 2) * tan_phi(intnl) - cohesion(intnl);
91 }

◆ yieldFunctionV()

void TensorMechanicsPlasticModel::yieldFunctionV ( const RankTwoTensor stress,
Real  intnl,
std::vector< Real > &  f 
) const
virtualinherited

Calculates the yield functions.

Note that for single-surface plasticity you don't want to override this - override the private yieldFunction below

Parameters
stressthe stress at which to calculate the yield function
intnlinternal parameter
[out]fthe yield functions

Reimplemented in TensorMechanicsPlasticMohrCoulombMulti, and TensorMechanicsPlasticTensileMulti.

Definition at line 69 of file TensorMechanicsPlasticModel.C.

72 {
73  f.assign(1, yieldFunction(stress, intnl));
74 }

Referenced by TensorMechanicsPlasticModel::returnMap().

Member Data Documentation

◆ _cap_rate

Real TensorMechanicsPlasticWeakPlaneShear::_cap_rate
protected

dictates how quickly the 'cap' degenerates to a hemisphere - see doco for _tip_scheme

Definition at line 79 of file TensorMechanicsPlasticWeakPlaneShear.h.

Referenced by d2smooth(), dsmooth(), and smooth().

◆ _cap_start

Real TensorMechanicsPlasticWeakPlaneShear::_cap_start
protected

smoothing parameter dictating when the 'cap' will start - see doco for _tip_scheme

Definition at line 76 of file TensorMechanicsPlasticWeakPlaneShear.h.

Referenced by d2smooth(), dsmooth(), and smooth().

◆ _cohesion

const TensorMechanicsHardeningModel& TensorMechanicsPlasticWeakPlaneShear::_cohesion
protected

Hardening model for cohesion.

Definition at line 42 of file TensorMechanicsPlasticWeakPlaneShear.h.

Referenced by cohesion(), and dcohesion().

◆ _f_tol

const Real TensorMechanicsPlasticModel::_f_tol
inherited

◆ _ic_tol

const Real TensorMechanicsPlasticModel::_ic_tol
inherited

Tolerance on internal constraint.

Definition at line 178 of file TensorMechanicsPlasticModel.h.

◆ _small_smoother2

Real TensorMechanicsPlasticWeakPlaneShear::_small_smoother2
protected

smoothing parameter for the cone's tip - see doco for _tip_scheme

Definition at line 73 of file TensorMechanicsPlasticWeakPlaneShear.h.

Referenced by smooth().

◆ _tan_phi

const TensorMechanicsHardeningModel& TensorMechanicsPlasticWeakPlaneShear::_tan_phi
protected

Hardening model for tan(phi)

Definition at line 45 of file TensorMechanicsPlasticWeakPlaneShear.h.

Referenced by dtan_phi(), and tan_phi().

◆ _tan_psi

const TensorMechanicsHardeningModel& TensorMechanicsPlasticWeakPlaneShear::_tan_psi
protected

Hardening model for tan(psi)

Definition at line 48 of file TensorMechanicsPlasticWeakPlaneShear.h.

Referenced by dtan_psi(), and tan_psi().

◆ _tip_scheme

MooseEnum TensorMechanicsPlasticWeakPlaneShear::_tip_scheme
protected

The yield function is modified to f = sqrt(s_xz^2 + s_yz^2 + a) + s_zz*_tan_phi - _cohesion where "a" depends on the tip_scheme.

Currently _tip_scheme is 'hyperbolic', where a = _small_smoother2 'cap' where a = _small_smoother2 + (p(stress(2,2) - _cap_start))^2 with the function p(x)=x(1-exp(-_cap_rate*x)) for x>0, and p=0 otherwise

Definition at line 70 of file TensorMechanicsPlasticWeakPlaneShear.h.

Referenced by d2smooth(), dsmooth(), and smooth().


The documentation for this class was generated from the following files:
TensorMechanicsPlasticModel::numberSurfaces
virtual unsigned int numberSurfaces() const
The number of yield surfaces for this plasticity model.
Definition: TensorMechanicsPlasticModel.C:57
TensorMechanicsPlasticWeakPlaneShear::tan_psi
virtual Real tan_psi(const Real internal_param) const
tan_psi as a function of internal parameter
Definition: TensorMechanicsPlasticWeakPlaneShear.C:209
TensorMechanicsPlasticWeakPlaneShear::_small_smoother2
Real _small_smoother2
smoothing parameter for the cone's tip - see doco for _tip_scheme
Definition: TensorMechanicsPlasticWeakPlaneShear.h:73
TensorMechanicsPlasticModel::dflowPotential_dintnl
virtual RankTwoTensor dflowPotential_dintnl(const RankTwoTensor &stress, Real intnl) const
The derivative of the flow potential with respect to the internal parameter.
Definition: TensorMechanicsPlasticModel.C:133
TensorMechanicsPlasticModel::validParams
static InputParameters validParams()
Definition: TensorMechanicsPlasticModel.C:18
TensorMechanicsPlasticModel::dyieldFunction_dstress
virtual RankTwoTensor dyieldFunction_dstress(const RankTwoTensor &stress, Real intnl) const
The derivative of yield function with respect to stress.
Definition: TensorMechanicsPlasticModel.C:77
TensorMechanicsPlasticModel::dhardPotential_dintnl
virtual Real dhardPotential_dintnl(const RankTwoTensor &stress, Real intnl) const
The derivative of the hardening potential with respect to the internal parameter.
Definition: TensorMechanicsPlasticModel.C:174
TensorMechanicsPlasticWeakPlaneShear::_cap_start
Real _cap_start
smoothing parameter dictating when the 'cap' will start - see doco for _tip_scheme
Definition: TensorMechanicsPlasticWeakPlaneShear.h:76
TensorMechanicsPlasticWeakPlaneShear::cohesion
virtual Real cohesion(const Real internal_param) const
cohesion as a function of internal parameter
Definition: TensorMechanicsPlasticWeakPlaneShear.C:185
TensorMechanicsPlasticWeakPlaneShear::_tip_scheme
MooseEnum _tip_scheme
The yield function is modified to f = sqrt(s_xz^2 + s_yz^2 + a) + s_zz*_tan_phi - _cohesion where "a"...
Definition: TensorMechanicsPlasticWeakPlaneShear.h:70
TensorMechanicsPlasticWeakPlaneShear::dtan_phi
virtual Real dtan_phi(const Real internal_param) const
d(tan_phi)/d(internal_param);
Definition: TensorMechanicsPlasticWeakPlaneShear.C:203
TensorMechanicsPlasticWeakPlaneShear::dtan_psi
virtual Real dtan_psi(const Real internal_param) const
d(tan_psi)/d(internal_param);
Definition: TensorMechanicsPlasticWeakPlaneShear.C:215
TensorMechanicsPlasticModel::KuhnTuckerSingleSurface
bool KuhnTuckerSingleSurface(Real yf, Real dpm, Real dpm_tol) const
Returns true if the Kuhn-Tucker conditions for the single surface are satisfied.
Definition: TensorMechanicsPlasticModel.C:248
TensorMechanicsPlasticWeakPlaneShear::d2smooth
virtual Real d2smooth(const RankTwoTensor &stress) const
returns the d^2a/dstress(2,2)^2 - see doco for _tip_scheme
Definition: TensorMechanicsPlasticWeakPlaneShear.C:256
TensorMechanicsPlasticWeakPlaneShear::_tan_phi
const TensorMechanicsHardeningModel & _tan_phi
Hardening model for tan(phi)
Definition: TensorMechanicsPlasticWeakPlaneShear.h:45
TensorMechanicsPlasticWeakPlaneShear::_cap_rate
Real _cap_rate
dictates how quickly the 'cap' degenerates to a hemisphere - see doco for _tip_scheme
Definition: TensorMechanicsPlasticWeakPlaneShear.h:79
TensorMechanicsPlasticModel::_f_tol
const Real _f_tol
Tolerance on yield function.
Definition: TensorMechanicsPlasticModel.h:175
TensorMechanicsPlasticModel::TensorMechanicsPlasticModel
TensorMechanicsPlasticModel(const InputParameters &parameters)
Definition: TensorMechanicsPlasticModel.C:34
TensorMechanicsHardeningModel::derivative
virtual Real derivative(Real intnl) const
Definition: TensorMechanicsHardeningModel.C:47
TensorMechanicsPlasticWeakPlaneShear::_cohesion
const TensorMechanicsHardeningModel & _cohesion
Hardening model for cohesion.
Definition: TensorMechanicsPlasticWeakPlaneShear.h:42
TensorMechanicsPlasticModel::hardPotential
virtual Real hardPotential(const RankTwoTensor &stress, Real intnl) const
The hardening potential.
Definition: TensorMechanicsPlasticModel.C:147
TensorMechanicsPlasticModel::yieldFunctionV
virtual void yieldFunctionV(const RankTwoTensor &stress, Real intnl, std::vector< Real > &f) const
Calculates the yield functions.
Definition: TensorMechanicsPlasticModel.C:69
TensorMechanicsHardeningModel::value
virtual Real value(Real intnl) const
Definition: TensorMechanicsHardeningModel.C:45
RankTwoTensor
RankTwoTensorTempl< Real > RankTwoTensor
Definition: ACGrGrElasticDrivingForce.h:17
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
TensorMechanicsPlasticWeakPlaneShear::smooth
virtual Real smooth(const RankTwoTensor &stress) const
returns the 'a' parameter - see doco for _tip_scheme
Definition: TensorMechanicsPlasticWeakPlaneShear.C:221
TensorMechanicsPlasticModel::dflowPotential_dstress
virtual RankFourTensor dflowPotential_dstress(const RankTwoTensor &stress, Real intnl) const
The derivative of the flow potential with respect to stress.
Definition: TensorMechanicsPlasticModel.C:119
TensorMechanicsPlasticModel::yieldFunction
virtual Real yieldFunction(const RankTwoTensor &stress, Real intnl) const
The following functions are what you should override when building single-plasticity models.
Definition: TensorMechanicsPlasticModel.C:63
TensorMechanicsPlasticModel::flowPotential
virtual RankTwoTensor flowPotential(const RankTwoTensor &stress, Real intnl) const
The flow potential.
Definition: TensorMechanicsPlasticModel.C:106
TensorMechanicsPlasticWeakPlaneShear::dcohesion
virtual Real dcohesion(const Real internal_param) const
d(cohesion)/d(internal_param)
Definition: TensorMechanicsPlasticWeakPlaneShear.C:191
TensorMechanicsPlasticWeakPlaneShear::tan_phi
virtual Real tan_phi(const Real internal_param) const
tan_phi as a function of internal parameter
Definition: TensorMechanicsPlasticWeakPlaneShear.C:197
RankTwoTensorTempl< Real >
TensorMechanicsPlasticWeakPlaneShear::_tan_psi
const TensorMechanicsHardeningModel & _tan_psi
Hardening model for tan(psi)
Definition: TensorMechanicsPlasticWeakPlaneShear.h:48
TensorMechanicsPlasticModel::dhardPotential_dstress
virtual RankTwoTensor dhardPotential_dstress(const RankTwoTensor &stress, Real intnl) const
The derivative of the hardening potential with respect to stress.
Definition: TensorMechanicsPlasticModel.C:160
TensorMechanicsPlasticWeakPlaneShear::df_dsig
RankTwoTensor df_dsig(const RankTwoTensor &stress, Real _tan_phi_or_psi) const
Function that's used in dyieldFunction_dstress and flowPotential.
Definition: TensorMechanicsPlasticWeakPlaneShear.C:94
TensorMechanicsPlasticWeakPlaneShear::dsmooth
virtual Real dsmooth(const RankTwoTensor &stress) const
returns the da/dstress(2,2) - see doco for _tip_scheme
Definition: TensorMechanicsPlasticWeakPlaneShear.C:236
TensorMechanicsPlasticModel::dyieldFunction_dintnl
virtual Real dyieldFunction_dintnl(const RankTwoTensor &stress, Real intnl) const
The derivative of yield function with respect to the internal parameter.
Definition: TensorMechanicsPlasticModel.C:92