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

Rate-independent non-associative Drucker Prager with hardening/softening. More...

#include <TensorMechanicsPlasticDruckerPrager.h>

Inheritance diagram for TensorMechanicsPlasticDruckerPrager:
[legend]

Public Types

enum  FrictionDilation { friction = 0, dilation = 1 }
 bbb (friction) and bbb_flow (dilation) are computed using the same function, onlyB, and this parameter tells that function whether to compute bbb or bbb_flow More...
 

Public Member Functions

 TensorMechanicsPlasticDruckerPrager (const InputParameters &parameters)
 
virtual std::string modelName () const override
 
void bothAB (Real intnl, Real &aaa, Real &bbb) const
 Calculates aaa and bbb as a function of the internal parameter intnl. More...
 
void dbothAB (Real intnl, Real &daaa, Real &dbbb) const
 Calculates d(aaa)/d(intnl) and d(bbb)/d(intnl) as a function of the internal parameter intnl. More...
 
void onlyB (Real intnl, int fd, Real &bbb) const
 Calculate bbb or bbb_flow. More...
 
void donlyB (Real intnl, int fd, Real &dbbb) const
 Calculate d(bbb)/d(intnl) or d(bbb_flow)/d(intnl) More...
 
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 void activeConstraints (const std::vector< Real > &f, const RankTwoTensor &stress, Real intnl, const RankFourTensor &Eijkl, std::vector< bool > &act, RankTwoTensor &returned_stress) const
 The active yield surfaces, given a vector of yield functions. 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

virtual Real yieldFunction (const RankTwoTensor &stress, Real intnl) const override
 The following functions are what you should override when building single-plasticity models. More...
 
virtual RankTwoTensor dyieldFunction_dstress (const RankTwoTensor &stress, Real intnl) const override
 The derivative of yield function with respect to stress. More...
 
virtual Real dyieldFunction_dintnl (const RankTwoTensor &stress, Real intnl) const override
 The derivative of yield function with respect to the internal parameter. More...
 
virtual RankTwoTensor flowPotential (const RankTwoTensor &stress, Real intnl) const override
 The flow potential. More...
 
virtual RankFourTensor dflowPotential_dstress (const RankTwoTensor &stress, Real intnl) const override
 The derivative of the flow potential with respect to stress. More...
 
virtual RankTwoTensor dflowPotential_dintnl (const RankTwoTensor &stress, Real intnl) const override
 The derivative of the flow potential with respect to the internal parameter. More...
 
virtual RankTwoTensor df_dsig (const RankTwoTensor &stress, Real bbb) const
 Function that's used in dyieldFunction_dstress and flowPotential. 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_mc_cohesion
 Hardening model for cohesion. More...
 
const TensorMechanicsHardeningModel_mc_phi
 Hardening model for tan(phi) More...
 
const TensorMechanicsHardeningModel_mc_psi
 Hardening model for tan(psi) More...
 
const MooseEnum _mc_interpolation_scheme
 The parameters aaa and bbb are chosen to closely match the Mohr-Coulomb yield surface. More...
 
const bool _zero_cohesion_hardening
 True if there is no hardening of cohesion. More...
 
const bool _zero_phi_hardening
 True if there is no hardening of friction angle. More...
 
const bool _zero_psi_hardening
 True if there is no hardening of dilation angle. More...
 

Private Member Functions

void initializeAandB (Real intnl, Real &aaa, Real &bbb) const
 Returns the Drucker-Prager parameters A nice reference on the different relationships between Drucker-Prager and Mohr-Coulomb is H Jiang and X Yongli, A note on the Mohr-Coulomb and Drucker-Prager strength criteria, Mechanics Research Communications 38 (2011) 309-314. More...
 
void initializeB (Real intnl, int fd, Real &bbb) const
 Returns the Drucker-Prager parameters A nice reference on the different relationships between Drucker-Prager and Mohr-Coulomb is H Jiang and X Yongli, A note on the Mohr-Coulomb and Drucker-Prager strength criteria, Mechanics Research Communications 38 (2011) 309-314. More...
 

Private Attributes

Real _aaa
 
Real _bbb
 
Real _bbb_flow
 

Detailed Description

Rate-independent non-associative Drucker Prager with hardening/softening.

The cone's tip is not smoothed. f = sqrt(J2) + bbb*Tr(stress) - aaa with aaa = "cohesion" (a non-negative number) and bbb = tan(friction angle) (a positive number)

Definition at line 27 of file TensorMechanicsPlasticDruckerPrager.h.

Member Enumeration Documentation

◆ FrictionDilation

bbb (friction) and bbb_flow (dilation) are computed using the same function, onlyB, and this parameter tells that function whether to compute bbb or bbb_flow

Enumerator
friction 
dilation 

Definition at line 46 of file TensorMechanicsPlasticDruckerPrager.h.

47  {
48  friction = 0,
49  dilation = 1
50  };

Constructor & Destructor Documentation

◆ TensorMechanicsPlasticDruckerPrager()

TensorMechanicsPlasticDruckerPrager::TensorMechanicsPlasticDruckerPrager ( const InputParameters &  parameters)

Definition at line 54 of file TensorMechanicsPlasticDruckerPrager.C.

56  : TensorMechanicsPlasticModel(parameters),
57  _mc_cohesion(getUserObject<TensorMechanicsHardeningModel>("mc_cohesion")),
58  _mc_phi(getUserObject<TensorMechanicsHardeningModel>("mc_friction_angle")),
59  _mc_psi(getUserObject<TensorMechanicsHardeningModel>("mc_dilation_angle")),
60  _mc_interpolation_scheme(getParam<MooseEnum>("mc_interpolation_scheme")),
61  _zero_cohesion_hardening(_mc_cohesion.modelName().compare("Constant") == 0),
62  _zero_phi_hardening(_mc_phi.modelName().compare("Constant") == 0),
63  _zero_psi_hardening(_mc_psi.modelName().compare("Constant") == 0)
64 {
65  if (_mc_phi.value(0.0) < 0.0 || _mc_psi.value(0.0) < 0.0 ||
66  _mc_phi.value(0.0) > libMesh::pi / 2.0 || _mc_psi.value(0.0) > libMesh::pi / 2.0)
67  mooseError("TensorMechanicsPlasticDruckerPrager: MC friction and dilation angles must lie in "
68  "[0, Pi/2]");
69  if (_mc_phi.value(0) < _mc_psi.value(0.0))
70  mooseError("TensorMechanicsPlasticDruckerPrager: MC friction angle must not be less than MC "
71  "dilation angle");
72  if (_mc_cohesion.value(0.0) < 0)
73  mooseError("TensorMechanicsPlasticDruckerPrager: MC cohesion should not be negative");
74 
75  initializeAandB(0.0, _aaa, _bbb);
77 }

Member Function Documentation

◆ activeConstraints()

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

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 in TensorMechanicsPlasticMohrCoulombMulti, TensorMechanicsPlasticTensileMulti, TensorMechanicsPlasticMeanCapTC, TensorMechanicsPlasticWeakPlaneShear, and TensorMechanicsPlasticWeakPlaneTensile.

Definition at line 188 of file TensorMechanicsPlasticModel.C.

194 {
195  mooseAssert(f.size() == numberSurfaces(),
196  "f incorrectly sized at " << f.size() << " in activeConstraints");
197  act.resize(numberSurfaces());
198  for (unsigned surface = 0; surface < numberSurfaces(); ++surface)
199  act[surface] = (f[surface] > _f_tol);
200 }

◆ bothAB()

void TensorMechanicsPlasticDruckerPrager::bothAB ( Real  intnl,
Real &  aaa,
Real &  bbb 
) const

◆ 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().

◆ dbothAB()

void TensorMechanicsPlasticDruckerPrager::dbothAB ( Real  intnl,
Real &  daaa,
Real &  dbbb 
) const

Calculates d(aaa)/d(intnl) and d(bbb)/d(intnl) as a function of the internal parameter intnl.

Definition at line 219 of file TensorMechanicsPlasticDruckerPrager.C.

220 {
222  {
223  daaa = 0;
224  dbbb = 0;
225  return;
226  }
227 
228  const Real C = _mc_cohesion.value(intnl);
229  const Real dC = _mc_cohesion.derivative(intnl);
230  const Real cosphi = std::cos(_mc_phi.value(intnl));
231  const Real dcosphi = -std::sin(_mc_phi.value(intnl)) * _mc_phi.derivative(intnl);
232  const Real sinphi = std::sin(_mc_phi.value(intnl));
233  const Real dsinphi = std::cos(_mc_phi.value(intnl)) * _mc_phi.derivative(intnl);
234  switch (_mc_interpolation_scheme)
235  {
236  case 0: // outer_tip
237  daaa = 2.0 * std::sqrt(3.0) *
238  (dC * cosphi / (3.0 - sinphi) + C * dcosphi / (3.0 - sinphi) +
239  C * cosphi * dsinphi / Utility::pow<2>(3.0 - sinphi));
240  dbbb = 2.0 / std::sqrt(3.0) *
241  (dsinphi / (3.0 - sinphi) + sinphi * dsinphi / Utility::pow<2>(3.0 - sinphi));
242  break;
243  case 1: // inner_tip
244  daaa = 2.0 * std::sqrt(3.0) *
245  (dC * cosphi / (3.0 + sinphi) + C * dcosphi / (3.0 + sinphi) -
246  C * cosphi * dsinphi / Utility::pow<2>(3.0 + sinphi));
247  dbbb = 2.0 / std::sqrt(3.0) *
248  (dsinphi / (3.0 + sinphi) - sinphi * dsinphi / Utility::pow<2>(3.0 + sinphi));
249  break;
250  case 2: // lode_zero
251  daaa = dC * cosphi + C * dcosphi;
252  dbbb = dsinphi / 3.0;
253  break;
254  case 3: // inner_edge
255  daaa = 3.0 * dC * cosphi / std::sqrt(9.0 + 3.0 * Utility::pow<2>(sinphi)) +
256  3.0 * C * dcosphi / std::sqrt(9.0 + 3.0 * Utility::pow<2>(sinphi)) -
257  3.0 * C * cosphi * 3.0 * sinphi * dsinphi /
258  std::pow(9.0 + 3.0 * Utility::pow<2>(sinphi), 1.5);
259  dbbb = dsinphi / std::sqrt(9.0 + 3.0 * Utility::pow<2>(sinphi)) -
260  3.0 * sinphi * sinphi * dsinphi / std::pow(9.0 + 3.0 * Utility::pow<2>(sinphi), 1.5);
261  break;
262  case 4: // native
263  daaa = dC;
264  dbbb = dsinphi / cosphi - sinphi * dcosphi / Utility::pow<2>(cosphi);
265  break;
266  }
267 }

Referenced by CappedDruckerPragerStressUpdate::computeAllQ(), dyieldFunction_dintnl(), and TensorMechanicsPlasticDruckerPragerHyperbolic::returnMap().

◆ df_dsig()

RankTwoTensor TensorMechanicsPlasticDruckerPrager::df_dsig ( const RankTwoTensor stress,
Real  bbb 
) const
protectedvirtual

Function that's used in dyieldFunction_dstress and flowPotential.

Reimplemented in TensorMechanicsPlasticDruckerPragerHyperbolic.

Definition at line 89 of file TensorMechanicsPlasticDruckerPrager.C.

90 {
91  return 0.5 * stress.dsecondInvariant() / std::sqrt(stress.secondInvariant()) +
92  stress.dtrace() * bbb;
93 }

Referenced by dyieldFunction_dstress(), and flowPotential().

◆ dflowPotential_dintnl()

RankTwoTensor TensorMechanicsPlasticDruckerPrager::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 134 of file TensorMechanicsPlasticDruckerPrager.C.

136 {
137  Real dbbb;
138  donlyB(intnl, dilation, dbbb);
139  return stress.dtrace() * dbbb;
140 }

◆ 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 TensorMechanicsPlasticDruckerPrager::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.

Reimplemented in TensorMechanicsPlasticDruckerPragerHyperbolic.

Definition at line 123 of file TensorMechanicsPlasticDruckerPrager.C.

125 {
126  RankFourTensor dr_dstress;
127  dr_dstress = 0.5 * stress.d2secondInvariant() / std::sqrt(stress.secondInvariant());
128  dr_dstress += -0.5 * 0.5 * stress.dsecondInvariant().outerProduct(stress.dsecondInvariant()) /
129  std::pow(stress.secondInvariant(), 1.5);
130  return dr_dstress;
131 }

◆ 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 }

◆ donlyB()

void TensorMechanicsPlasticDruckerPrager::donlyB ( Real  intnl,
int  fd,
Real &  dbbb 
) const

Calculate d(bbb)/d(intnl) or d(bbb_flow)/d(intnl)

Parameters
intnlthe internal parameter
fdif fd==friction then bbb is calculated. if fd==dilation then bbb_flow is calculated
bbbeither bbb or bbb_flow, depending on fd

Definition at line 177 of file TensorMechanicsPlasticDruckerPrager.C.

178 {
179  if (_zero_phi_hardening && (fd == friction))
180  {
181  dbbb = 0;
182  return;
183  }
184  if (_zero_psi_hardening && (fd == dilation))
185  {
186  dbbb = 0;
187  return;
188  }
189  const Real s = (fd == friction) ? std::sin(_mc_phi.value(intnl)) : std::sin(_mc_psi.value(intnl));
190  const Real ds = (fd == friction) ? std::cos(_mc_phi.value(intnl)) * _mc_phi.derivative(intnl)
191  : std::cos(_mc_psi.value(intnl)) * _mc_psi.derivative(intnl);
192  switch (_mc_interpolation_scheme)
193  {
194  case 0: // outer_tip
195  dbbb = 2.0 / std::sqrt(3.0) * (ds / (3.0 - s) + s * ds / Utility::pow<2>(3.0 - s));
196  break;
197  case 1: // inner_tip
198  dbbb = 2.0 / std::sqrt(3.0) * (ds / (3.0 + s) - s * ds / Utility::pow<2>(3.0 + s));
199  break;
200  case 2: // lode_zero
201  dbbb = ds / 3.0;
202  break;
203  case 3: // inner_edge
204  dbbb = ds / std::sqrt(9.0 + 3.0 * Utility::pow<2>(s)) -
205  3 * s * s * ds / std::pow(9.0 + 3.0 * Utility::pow<2>(s), 1.5);
206  break;
207  case 4: // native
208  const Real c =
209  (fd == friction) ? std::cos(_mc_phi.value(intnl)) : std::cos(_mc_psi.value(intnl));
210  const Real dc = (fd == friction)
211  ? -std::sin(_mc_phi.value(intnl)) * _mc_phi.derivative(intnl)
212  : -std::sin(_mc_psi.value(intnl)) * _mc_psi.derivative(intnl);
213  dbbb = ds / c - s * dc / Utility::pow<2>(c);
214  break;
215  }
216 }

Referenced by CappedDruckerPragerStressUpdate::computeAllQ(), TensorMechanicsPlasticDruckerPragerHyperbolic::consistentTangentOperator(), dflowPotential_dintnl(), TensorMechanicsPlasticDruckerPragerHyperbolic::returnMap(), and CappedDruckerPragerStressUpdate::setIntnlDerivatives().

◆ dyieldFunction_dintnl()

Real TensorMechanicsPlasticDruckerPrager::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 105 of file TensorMechanicsPlasticDruckerPrager.C.

107 {
108  Real daaa;
109  Real dbbb;
110  dbothAB(intnl, daaa, dbbb);
111  return stress.trace() * dbbb - daaa;
112 }

Referenced by TensorMechanicsPlasticDruckerPragerHyperbolic::consistentTangentOperator().

◆ 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 TensorMechanicsPlasticDruckerPrager::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 96 of file TensorMechanicsPlasticDruckerPrager.C.

98 {
99  Real bbb;
100  onlyB(intnl, friction, bbb);
101  return df_dsig(stress, bbb);
102 }

◆ 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 TensorMechanicsPlasticDruckerPrager::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 115 of file TensorMechanicsPlasticDruckerPrager.C.

116 {
117  Real bbb_flow;
118  onlyB(intnl, dilation, bbb_flow);
119  return df_dsig(stress, bbb_flow);
120 }

◆ 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 }

◆ initializeAandB()

void TensorMechanicsPlasticDruckerPrager::initializeAandB ( Real  intnl,
Real &  aaa,
Real &  bbb 
) const
private

Returns the Drucker-Prager parameters A nice reference on the different relationships between Drucker-Prager and Mohr-Coulomb is H Jiang and X Yongli, A note on the Mohr-Coulomb and Drucker-Prager strength criteria, Mechanics Research Communications 38 (2011) 309-314.

This function uses Table1 in that reference, with the following translations aaa = k bbb = -alpha/3

Parameters
intnlThe internal parameter
[out]aaaThe Drucker-Prager aaa quantity
[out]bbbThe Drucker-Prager bbb quantity

Definition at line 270 of file TensorMechanicsPlasticDruckerPrager.C.

271 {
272  const Real C = _mc_cohesion.value(intnl);
273  const Real cosphi = std::cos(_mc_phi.value(intnl));
274  const Real sinphi = std::sin(_mc_phi.value(intnl));
275  switch (_mc_interpolation_scheme)
276  {
277  case 0: // outer_tip
278  aaa = 2.0 * std::sqrt(3.0) * C * cosphi / (3.0 - sinphi);
279  bbb = 2.0 * sinphi / std::sqrt(3.0) / (3.0 - sinphi);
280  break;
281  case 1: // inner_tip
282  aaa = 2.0 * std::sqrt(3.0) * C * cosphi / (3.0 + sinphi);
283  bbb = 2.0 * sinphi / std::sqrt(3.0) / (3.0 + sinphi);
284  break;
285  case 2: // lode_zero
286  aaa = C * cosphi;
287  bbb = sinphi / 3.0;
288  break;
289  case 3: // inner_edge
290  aaa = 3.0 * C * cosphi / std::sqrt(9.0 + 3.0 * Utility::pow<2>(sinphi));
291  bbb = sinphi / std::sqrt(9.0 + 3.0 * Utility::pow<2>(sinphi));
292  break;
293  case 4: // native
294  aaa = C;
295  bbb = sinphi / cosphi;
296  break;
297  }
298 }

Referenced by bothAB(), and TensorMechanicsPlasticDruckerPrager().

◆ initializeB()

void TensorMechanicsPlasticDruckerPrager::initializeB ( Real  intnl,
int  fd,
Real &  bbb 
) const
private

Returns the Drucker-Prager parameters A nice reference on the different relationships between Drucker-Prager and Mohr-Coulomb is H Jiang and X Yongli, A note on the Mohr-Coulomb and Drucker-Prager strength criteria, Mechanics Research Communications 38 (2011) 309-314.

This function uses Table1 in that reference, with the following translations aaa = k bbb = -alpha/3

Parameters
intnlThe internal parameter
fdIf fd == frction then the friction angle is used to set bbb, otherwise the dilation angle is used
[out]bbbThe Drucker-Prager bbb quantity

Definition at line 301 of file TensorMechanicsPlasticDruckerPrager.C.

302 {
303  const Real s = (fd == friction) ? std::sin(_mc_phi.value(intnl)) : std::sin(_mc_psi.value(intnl));
304  switch (_mc_interpolation_scheme)
305  {
306  case 0: // outer_tip
307  bbb = 2.0 * s / std::sqrt(3.0) / (3.0 - s);
308  break;
309  case 1: // inner_tip
310  bbb = 2.0 * s / std::sqrt(3.0) / (3.0 + s);
311  break;
312  case 2: // lode_zero
313  bbb = s / 3.0;
314  break;
315  case 3: // inner_edge
316  bbb = s / std::sqrt(9.0 + 3.0 * Utility::pow<2>(s));
317  break;
318  case 4: // native
319  const Real c =
320  (fd == friction) ? std::cos(_mc_phi.value(intnl)) : std::cos(_mc_psi.value(intnl));
321  bbb = s / c;
322  break;
323  }
324 }

Referenced by onlyB(), and TensorMechanicsPlasticDruckerPrager().

◆ 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 TensorMechanicsPlasticDruckerPrager::modelName ( ) const
overridevirtual

Implements TensorMechanicsPlasticModel.

Reimplemented in TensorMechanicsPlasticDruckerPragerHyperbolic.

Definition at line 143 of file TensorMechanicsPlasticDruckerPrager.C.

144 {
145  return "DruckerPrager";
146 }

◆ numberSurfaces()

unsigned TensorMechanicsPlasticModel::numberSurfaces ( ) const
virtualinherited

◆ onlyB()

void TensorMechanicsPlasticDruckerPrager::onlyB ( Real  intnl,
int  fd,
Real &  bbb 
) const

Calculate bbb or bbb_flow.

Parameters
intnlthe internal parameter
fdif fd==friction then bbb is calculated. if fd==dilation then bbb_flow is calculated
bbbeither bbb or bbb_flow, depending on fd

Definition at line 161 of file TensorMechanicsPlasticDruckerPrager.C.

162 {
163  if (_zero_phi_hardening && (fd == friction))
164  {
165  bbb = _bbb;
166  return;
167  }
168  if (_zero_psi_hardening && (fd == dilation))
169  {
170  bbb = _bbb_flow;
171  return;
172  }
173  initializeB(intnl, fd, bbb);
174 }

Referenced by CappedDruckerPragerStressUpdate::computeAllQ(), TensorMechanicsPlasticDruckerPragerHyperbolic::consistentTangentOperator(), dyieldFunction_dstress(), flowPotential(), CappedDruckerPragerStressUpdate::initializeVars(), TensorMechanicsPlasticDruckerPragerHyperbolic::returnMap(), CappedDruckerPragerStressUpdate::setIntnlDerivatives(), and CappedDruckerPragerStressUpdate::setIntnlValues().

◆ 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().

◆ 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 TensorMechanicsPlasticDruckerPrager::validParams ( )
static

Definition at line 19 of file TensorMechanicsPlasticDruckerPrager.C.

20 {
21  InputParameters params = TensorMechanicsPlasticModel::validParams();
22  MooseEnum mc_interpolation_scheme("outer_tip=0 inner_tip=1 lode_zero=2 inner_edge=3 native=4",
23  "lode_zero");
24  params.addParam<MooseEnum>(
25  "mc_interpolation_scheme",
26  mc_interpolation_scheme,
27  "Scheme by which the Drucker-Prager cohesion, friction angle and dilation angle are set from "
28  "the Mohr-Coulomb parameters mc_cohesion, mc_friction_angle and mc_dilation_angle. Consider "
29  "the DP and MC yield surfaces on the deviatoric (octahedral) plane. Outer_tip: the DP "
30  "circle touches the outer tips of the MC hex. Inner_tip: the DP circle touches the inner "
31  "tips of the MC hex. Lode_zero: the DP circle intersects the MC hex at lode angle=0. "
32  "Inner_edge: the DP circle is the largest circle that wholly fits inside the MC hex. "
33  "Native: The DP cohesion, friction angle and dilation angle are set equal to the mc_ "
34  "parameters entered.");
35  params.addRequiredParam<UserObjectName>(
36  "mc_cohesion",
37  "A TensorMechanicsHardening UserObject that defines hardening of the "
38  "Mohr-Coulomb cohesion. Physically this should not be negative.");
39  params.addRequiredParam<UserObjectName>(
40  "mc_friction_angle",
41  "A TensorMechanicsHardening UserObject that defines hardening of the "
42  "Mohr-Coulomb friction angle (in radians). Physically this should be "
43  "between 0 and Pi/2.");
44  params.addRequiredParam<UserObjectName>(
45  "mc_dilation_angle",
46  "A TensorMechanicsHardening UserObject that defines hardening of the "
47  "Mohr-Coulomb dilation angle (in radians). Usually the dilation angle "
48  "is not greater than the friction angle, and it is between 0 and Pi/2.");
49  params.addClassDescription(
50  "Non-associative Drucker Prager plasticity with no smoothing of the cone tip.");
51  return params;
52 }

Referenced by TensorMechanicsPlasticDruckerPragerHyperbolic::validParams().

◆ yieldFunction()

Real TensorMechanicsPlasticDruckerPrager::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.

Reimplemented in TensorMechanicsPlasticDruckerPragerHyperbolic.

Definition at line 80 of file TensorMechanicsPlasticDruckerPrager.C.

81 {
82  Real aaa;
83  Real bbb;
84  bothAB(intnl, aaa, bbb);
85  return std::sqrt(stress.secondInvariant()) + stress.trace() * bbb - aaa;
86 }

◆ 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

◆ _aaa

Real TensorMechanicsPlasticDruckerPrager::_aaa
private

◆ _bbb

Real TensorMechanicsPlasticDruckerPrager::_bbb
private

◆ _bbb_flow

Real TensorMechanicsPlasticDruckerPrager::_bbb_flow
private

◆ _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.

◆ _mc_cohesion

const TensorMechanicsHardeningModel& TensorMechanicsPlasticDruckerPrager::_mc_cohesion
protected

Hardening model for cohesion.

Definition at line 70 of file TensorMechanicsPlasticDruckerPrager.h.

Referenced by dbothAB(), initializeAandB(), and TensorMechanicsPlasticDruckerPrager().

◆ _mc_interpolation_scheme

const MooseEnum TensorMechanicsPlasticDruckerPrager::_mc_interpolation_scheme
protected

The parameters aaa and bbb are chosen to closely match the Mohr-Coulomb yield surface.

Various matching schemes may be used and this parameter holds the user's choice.

Definition at line 98 of file TensorMechanicsPlasticDruckerPrager.h.

Referenced by dbothAB(), donlyB(), initializeAandB(), and initializeB().

◆ _mc_phi

const TensorMechanicsHardeningModel& TensorMechanicsPlasticDruckerPrager::_mc_phi
protected

Hardening model for tan(phi)

Definition at line 73 of file TensorMechanicsPlasticDruckerPrager.h.

Referenced by dbothAB(), donlyB(), initializeAandB(), initializeB(), and TensorMechanicsPlasticDruckerPrager().

◆ _mc_psi

const TensorMechanicsHardeningModel& TensorMechanicsPlasticDruckerPrager::_mc_psi
protected

Hardening model for tan(psi)

Definition at line 76 of file TensorMechanicsPlasticDruckerPrager.h.

Referenced by donlyB(), initializeB(), and TensorMechanicsPlasticDruckerPrager().

◆ _zero_cohesion_hardening

const bool TensorMechanicsPlasticDruckerPrager::_zero_cohesion_hardening
protected

True if there is no hardening of cohesion.

Definition at line 101 of file TensorMechanicsPlasticDruckerPrager.h.

Referenced by bothAB(), and dbothAB().

◆ _zero_phi_hardening

const bool TensorMechanicsPlasticDruckerPrager::_zero_phi_hardening
protected

True if there is no hardening of friction angle.

Definition at line 104 of file TensorMechanicsPlasticDruckerPrager.h.

Referenced by bothAB(), dbothAB(), donlyB(), and onlyB().

◆ _zero_psi_hardening

const bool TensorMechanicsPlasticDruckerPrager::_zero_psi_hardening
protected

True if there is no hardening of dilation angle.

Definition at line 107 of file TensorMechanicsPlasticDruckerPrager.h.

Referenced by donlyB(), and onlyB().


The documentation for this class was generated from the following files:
TensorMechanicsPlasticDruckerPrager::initializeB
void initializeB(Real intnl, int fd, Real &bbb) const
Returns the Drucker-Prager parameters A nice reference on the different relationships between Drucker...
Definition: TensorMechanicsPlasticDruckerPrager.C:301
TensorMechanicsPlasticModel::numberSurfaces
virtual unsigned int numberSurfaces() const
The number of yield surfaces for this plasticity model.
Definition: TensorMechanicsPlasticModel.C:57
TensorMechanicsPlasticDruckerPrager::_bbb_flow
Real _bbb_flow
Definition: TensorMechanicsPlasticDruckerPrager.h:115
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
TensorMechanicsPlasticDruckerPrager::_zero_psi_hardening
const bool _zero_psi_hardening
True if there is no hardening of dilation angle.
Definition: TensorMechanicsPlasticDruckerPrager.h:107
TensorMechanicsPlasticModel::validParams
static InputParameters validParams()
Definition: TensorMechanicsPlasticModel.C:18
TensorMechanicsPlasticDruckerPrager::onlyB
void onlyB(Real intnl, int fd, Real &bbb) const
Calculate bbb or bbb_flow.
Definition: TensorMechanicsPlasticDruckerPrager.C:161
pow
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
Definition: ExpressionBuilder.h:673
TensorMechanicsPlasticDruckerPrager::donlyB
void donlyB(Real intnl, int fd, Real &dbbb) const
Calculate d(bbb)/d(intnl) or d(bbb_flow)/d(intnl)
Definition: TensorMechanicsPlasticDruckerPrager.C:177
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
TensorMechanicsPlasticDruckerPrager::df_dsig
virtual RankTwoTensor df_dsig(const RankTwoTensor &stress, Real bbb) const
Function that's used in dyieldFunction_dstress and flowPotential.
Definition: TensorMechanicsPlasticDruckerPrager.C:89
TensorMechanicsPlasticDruckerPrager::dilation
Definition: TensorMechanicsPlasticDruckerPrager.h:49
TensorMechanicsPlasticDruckerPrager::_zero_phi_hardening
const bool _zero_phi_hardening
True if there is no hardening of friction angle.
Definition: TensorMechanicsPlasticDruckerPrager.h:104
TensorMechanicsPlasticDruckerPrager::_bbb
Real _bbb
Definition: TensorMechanicsPlasticDruckerPrager.h:114
TensorMechanicsPlasticDruckerPrager::_mc_phi
const TensorMechanicsHardeningModel & _mc_phi
Hardening model for tan(phi)
Definition: TensorMechanicsPlasticDruckerPrager.h:73
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
TensorMechanicsPlasticDruckerPrager::_aaa
Real _aaa
Definition: TensorMechanicsPlasticDruckerPrager.h:113
TensorMechanicsPlasticDruckerPrager::friction
Definition: TensorMechanicsPlasticDruckerPrager.h:48
TensorMechanicsPlasticModel::_f_tol
const Real _f_tol
Tolerance on yield function.
Definition: TensorMechanicsPlasticModel.h:175
TensorMechanicsHardeningModel::modelName
virtual std::string modelName() const =0
TensorMechanicsPlasticModel::TensorMechanicsPlasticModel
TensorMechanicsPlasticModel(const InputParameters &parameters)
Definition: TensorMechanicsPlasticModel.C:34
TensorMechanicsPlasticDruckerPrager::_zero_cohesion_hardening
const bool _zero_cohesion_hardening
True if there is no hardening of cohesion.
Definition: TensorMechanicsPlasticDruckerPrager.h:101
TensorMechanicsHardeningModel::derivative
virtual Real derivative(Real intnl) const
Definition: TensorMechanicsHardeningModel.C:47
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
TensorMechanicsPlasticDruckerPrager::_mc_psi
const TensorMechanicsHardeningModel & _mc_psi
Hardening model for tan(psi)
Definition: TensorMechanicsPlasticDruckerPrager.h:76
TensorMechanicsHardeningModel::value
virtual Real value(Real intnl) const
Definition: TensorMechanicsHardeningModel.C:45
RankTwoTensor
RankTwoTensorTempl< Real > RankTwoTensor
Definition: ACGrGrElasticDrivingForce.h:17
TensorMechanicsPlasticDruckerPrager::_mc_interpolation_scheme
const MooseEnum _mc_interpolation_scheme
The parameters aaa and bbb are chosen to closely match the Mohr-Coulomb yield surface.
Definition: TensorMechanicsPlasticDruckerPrager.h:98
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
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
TensorMechanicsPlasticDruckerPrager::_mc_cohesion
const TensorMechanicsHardeningModel & _mc_cohesion
Hardening model for cohesion.
Definition: TensorMechanicsPlasticDruckerPrager.h:70
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
TensorMechanicsPlasticDruckerPrager::dbothAB
void dbothAB(Real intnl, Real &daaa, Real &dbbb) const
Calculates d(aaa)/d(intnl) and d(bbb)/d(intnl) as a function of the internal parameter intnl.
Definition: TensorMechanicsPlasticDruckerPrager.C:219
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
TensorMechanicsPlasticDruckerPrager::initializeAandB
void initializeAandB(Real intnl, Real &aaa, Real &bbb) const
Returns the Drucker-Prager parameters A nice reference on the different relationships between Drucker...
Definition: TensorMechanicsPlasticDruckerPrager.C:270
TensorMechanicsPlasticDruckerPrager::bothAB
void bothAB(Real intnl, Real &aaa, Real &bbb) const
Calculates aaa and bbb as a function of the internal parameter intnl.
Definition: TensorMechanicsPlasticDruckerPrager.C:149
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