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

Rate-independent associative mean-cap tensile AND compressive failure with hardening/softening of the tensile and compressive strength. More...

#include <TensorMechanicsPlasticMeanCapTC.h>

Inheritance diagram for TensorMechanicsPlasticMeanCapTC:
[legend]

Public Member Functions

 TensorMechanicsPlasticMeanCapTC (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 bool useCustomReturnMap () const override
 Returns false. You will want to override this in your derived class if you write a custom returnMap function. More...
 
virtual bool useCustomCTO () const override
 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 override
 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 override
 Calculates a custom consistent tangent operator. 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...
 
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...
 
Real hardPotential (const RankTwoTensor &stress, Real intnl) const override
 The hardening potential. More...
 
virtual RankTwoTensor dhardPotential_dstress (const RankTwoTensor &stress, Real intnl) const override
 The derivative of the hardening potential with respect to stress. More...
 
virtual Real dhardPotential_dintnl (const RankTwoTensor &stress, Real intnl) const override
 The derivative of the hardening potential with respect to the internal parameter. More...
 
RankTwoTensor df_dsig (const RankTwoTensor &stress, Real intnl) const
 Derivative of the yield function with respect to stress. More...
 
virtual Real tensile_strength (const Real internal_param) const
 tensile strength as a function of residual value, rate, and internal_param More...
 
virtual Real dtensile_strength (const Real internal_param) const
 d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param More...
 
virtual Real compressive_strength (const Real internal_param) const
 compressive strength as a function of residual value, rate, and internal_param More...
 
virtual Real dcompressive_strength (const Real internal_param) const
 d(compressive strength)/d(internal_param) as a function of residual value, rate, and internal_param More...
 

Protected Attributes

const unsigned _max_iters
 max iters for custom return map loop More...
 
const bool _use_custom_returnMap
 Whether to use the custom return-map algorithm. More...
 
const bool _use_custom_cto
 Whether to use the custom consistent tangent operator algorithm. More...
 
const TensorMechanicsHardeningModel_strength
 the tensile strength More...
 
const TensorMechanicsHardeningModel_c_strength
 the compressive strength More...
 

Detailed Description

Rate-independent associative mean-cap tensile AND compressive failure with hardening/softening of the tensile and compressive strength.

The key point here is that the internal parameter is equal to the volumetric plastic strain. This means that upon tensile failure, the compressive strength can soften (using a TensorMechanicsHardening object) which physically means that a subsequent compressive stress can easily squash the material

Definition at line 29 of file TensorMechanicsPlasticMeanCapTC.h.

Constructor & Destructor Documentation

◆ TensorMechanicsPlasticMeanCapTC()

TensorMechanicsPlasticMeanCapTC::TensorMechanicsPlasticMeanCapTC ( const InputParameters &  parameters)

Definition at line 49 of file TensorMechanicsPlasticMeanCapTC.C.

50  : TensorMechanicsPlasticModel(parameters),
51  _max_iters(getParam<unsigned>("max_iterations")),
52  _use_custom_returnMap(getParam<bool>("use_custom_returnMap")),
53  _use_custom_cto(getParam<bool>("use_custom_cto")),
54  _strength(getUserObject<TensorMechanicsHardeningModel>("tensile_strength")),
55  _c_strength(getUserObject<TensorMechanicsHardeningModel>("compressive_strength"))
56 {
57  // cannot check the following for all values of the internal parameter, but this will catch most
58  // errors
59  if (_strength.value(0) <= _c_strength.value(0))
60  mooseError("MeanCapTC: tensile strength (which is usually positive) must not be less than "
61  "compressive strength (which is usually negative)");
62 }

Member Function Documentation

◆ activeConstraints()

void TensorMechanicsPlasticMeanCapTC::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 247 of file TensorMechanicsPlasticMeanCapTC.C.

253 {
254  act.assign(1, false);
255 
256  if (f[0] <= _f_tol)
257  {
258  returned_stress = stress;
259  return;
260  }
261 
262  const Real tr = stress.trace();
263  const Real t_str = tensile_strength(intnl);
264  Real str;
265  Real dirn;
266  if (tr >= t_str)
267  {
268  str = t_str;
269  dirn = 1.0;
270  }
271  else
272  {
273  str = compressive_strength(intnl);
274  dirn = -1.0;
275  }
276 
277  RankTwoTensor n; // flow direction
278  for (unsigned i = 0; i < 3; ++i)
279  for (unsigned j = 0; j < 3; ++j)
280  for (unsigned k = 0; k < 3; ++k)
281  n(i, j) += dirn * Eijkl(i, j, k, k);
282 
283  // returned_stress = stress - gamma*n
284  // and taking the trace of this and using
285  // Tr(returned_stress) = str, gives
286  // gamma = (Tr(stress) - str)/Tr(n)
287  Real gamma = (stress.trace() - str) / n.trace();
288 
289  for (unsigned i = 0; i < 3; ++i)
290  for (unsigned j = 0; j < 3; ++j)
291  returned_stress(i, j) = stress(i, j) - gamma * n(i, j);
292 
293  act[0] = true;
294 }

◆ compressive_strength()

Real TensorMechanicsPlasticMeanCapTC::compressive_strength ( const Real  internal_param) const
protectedvirtual

compressive strength as a function of residual value, rate, and internal_param

Definition at line 235 of file TensorMechanicsPlasticMeanCapTC.C.

236 {
237  return _c_strength.value(internal_param);
238 }

Referenced by activeConstraints(), df_dsig(), dflowPotential_dintnl(), dflowPotential_dstress(), dhardPotential_dintnl(), dhardPotential_dstress(), dyieldFunction_dintnl(), hardPotential(), returnMap(), and yieldFunction().

◆ consistentTangentOperator()

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

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 from TensorMechanicsPlasticModel.

Definition at line 399 of file TensorMechanicsPlasticMeanCapTC.C.

406 {
407  if (!_use_custom_cto)
409  trial_stress, intnl_old, stress, intnl, E_ijkl, cumulative_pm);
410 
411  Real df_dq;
412  Real alpha;
413  if (trial_stress.trace() >= tensile_strength(intnl_old))
414  {
415  df_dq = -dtensile_strength(intnl);
416  alpha = 1.0;
417  }
418  else
419  {
420  df_dq = dcompressive_strength(intnl);
421  alpha = -1.0;
422  }
423 
424  RankTwoTensor elas;
425  for (unsigned int i = 0; i < 3; ++i)
426  for (unsigned int j = 0; j < 3; ++j)
427  for (unsigned int k = 0; k < 3; ++k)
428  elas(i, j) += E_ijkl(i, j, k, k);
429 
430  const Real hw = -df_dq + alpha * elas.trace();
431 
432  return E_ijkl - alpha / hw * elas.outerProduct(elas);
433 }

◆ dcompressive_strength()

Real TensorMechanicsPlasticMeanCapTC::dcompressive_strength ( const Real  internal_param) const
protectedvirtual

d(compressive strength)/d(internal_param) as a function of residual value, rate, and internal_param

Definition at line 241 of file TensorMechanicsPlasticMeanCapTC.C.

242 {
243  return _c_strength.derivative(internal_param);
244 }

Referenced by consistentTangentOperator(), dflowPotential_dintnl(), dhardPotential_dintnl(), dyieldFunction_dintnl(), and returnMap().

◆ df_dsig()

RankTwoTensor TensorMechanicsPlasticMeanCapTC::df_dsig ( const RankTwoTensor stress,
Real  intnl 
) const
protected

Derivative of the yield function with respect to stress.

This is also the flow potential.

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

Definition at line 110 of file TensorMechanicsPlasticMeanCapTC.C.

111 {
112  const Real tr = stress.trace();
113  const Real t_str = tensile_strength(intnl);
114  if (tr >= t_str)
115  return stress.dtrace();
116 
117  const Real c_str = compressive_strength(intnl);
118  if (tr <= c_str)
119  return -stress.dtrace();
120 
121  return -std::cos(libMesh::pi * (tr - c_str) / (t_str - c_str)) * stress.dtrace();
122 }

Referenced by dyieldFunction_dstress(), and flowPotential().

◆ dflowPotential_dintnl()

RankTwoTensor TensorMechanicsPlasticMeanCapTC::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 148 of file TensorMechanicsPlasticMeanCapTC.C.

150 {
151  const Real tr = stress.trace();
152  const Real t_str = tensile_strength(intnl);
153  if (tr >= t_str)
154  return RankTwoTensor();
155 
156  const Real c_str = compressive_strength(intnl);
157  if (tr <= c_str)
158  return RankTwoTensor();
159 
160  const Real dt = dtensile_strength(intnl);
161  const Real dc = dcompressive_strength(intnl);
162  return std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * stress.dtrace() * libMesh::pi /
163  Utility::pow<2>(t_str - c_str) * ((tr - t_str) * dc - (tr - c_str) * dt);
164 }

◆ 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 TensorMechanicsPlasticMeanCapTC::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 131 of file TensorMechanicsPlasticMeanCapTC.C.

133 {
134  const Real tr = stress.trace();
135  const Real t_str = tensile_strength(intnl);
136  if (tr >= t_str)
137  return RankFourTensor();
138 
139  const Real c_str = compressive_strength(intnl);
140  if (tr <= c_str)
141  return RankFourTensor();
142 
143  return libMesh::pi / (t_str - c_str) * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) *
144  stress.dtrace().outerProduct(stress.dtrace());
145 }

◆ 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 TensorMechanicsPlasticMeanCapTC::dhardPotential_dintnl ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

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 from TensorMechanicsPlasticModel.

Definition at line 204 of file TensorMechanicsPlasticMeanCapTC.C.

206 {
207  const Real tr = stress.trace();
208  const Real t_str = tensile_strength(intnl);
209  if (tr >= t_str)
210  return 0.0;
211 
212  const Real c_str = compressive_strength(intnl);
213  if (tr <= c_str)
214  return 0.0;
215 
216  const Real dt = dtensile_strength(intnl);
217  const Real dc = dcompressive_strength(intnl);
218  return -std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * libMesh::pi /
219  Utility::pow<2>(t_str - c_str) * ((tr - t_str) * dc - (tr - c_str) * dt);
220 }

◆ 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 TensorMechanicsPlasticMeanCapTC::dhardPotential_dstress ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

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 from TensorMechanicsPlasticModel.

Definition at line 187 of file TensorMechanicsPlasticMeanCapTC.C.

189 {
190  const Real tr = stress.trace();
191  const Real t_str = tensile_strength(intnl);
192  if (tr >= t_str)
193  return RankTwoTensor();
194 
195  const Real c_str = compressive_strength(intnl);
196  if (tr <= c_str)
197  return RankTwoTensor();
198 
199  return -std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * libMesh::pi / (t_str - c_str) *
200  stress.dtrace();
201 }

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

◆ dtensile_strength()

Real TensorMechanicsPlasticMeanCapTC::dtensile_strength ( const Real  internal_param) const
protectedvirtual

d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param

Definition at line 229 of file TensorMechanicsPlasticMeanCapTC.C.

230 {
231  return _strength.derivative(internal_param);
232 }

Referenced by consistentTangentOperator(), dflowPotential_dintnl(), dhardPotential_dintnl(), dyieldFunction_dintnl(), and returnMap().

◆ dyieldFunction_dintnl()

Real TensorMechanicsPlasticMeanCapTC::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 90 of file TensorMechanicsPlasticMeanCapTC.C.

92 {
93  const Real tr = stress.trace();
94  const Real t_str = tensile_strength(intnl);
95  if (tr >= t_str)
96  return -dtensile_strength(intnl);
97 
98  const Real c_str = compressive_strength(intnl);
99  if (tr <= c_str)
100  return dcompressive_strength(intnl);
101 
102  const Real dt = dtensile_strength(intnl);
103  const Real dc = dcompressive_strength(intnl);
104  return (dc - dt) / libMesh::pi * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) +
105  1.0 / (t_str - c_str) * std::cos(libMesh::pi * (tr - c_str) / (t_str - c_str)) *
106  ((tr - c_str) * dt - (tr - t_str) * dc);
107 }

◆ 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 TensorMechanicsPlasticMeanCapTC::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 83 of file TensorMechanicsPlasticMeanCapTC.C.

85 {
86  return df_dsig(stress, intnl);
87 }

◆ 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 TensorMechanicsPlasticMeanCapTC::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 125 of file TensorMechanicsPlasticMeanCapTC.C.

126 {
127  return df_dsig(stress, intnl);
128 }

◆ 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 TensorMechanicsPlasticMeanCapTC::hardPotential ( const RankTwoTensor stress,
Real  intnl 
) const
overrideprotectedvirtual

The hardening potential.

Note that it is -1 for stress.trace() > _strength, and +1 for stress.trace() < _c_strength. This implements the idea that tensile failure will cause a massive reduction in compressive strength

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

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 167 of file TensorMechanicsPlasticMeanCapTC.C.

168 {
169  // This is the key for this whole class!
170  const Real tr = stress.trace();
171  const Real t_str = tensile_strength(intnl);
172 
173  if (tr >= t_str)
174  return -1.0; // this will serve to *increase* the internal parameter (so internal parameter will
175  // be a measure of volumetric strain)
176 
177  const Real c_str = compressive_strength(intnl);
178  if (tr <= c_str)
179  return 1.0; // this will serve to *decrease* the internal parameter (so internal parameter will
180  // be a measure of volumetric strain)
181 
182  return std::cos(libMesh::pi * (tr - c_str) /
183  (t_str - c_str)); // this interpolates C2 smoothly between 1 and -1
184 }

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

Implements TensorMechanicsPlasticModel.

Definition at line 448 of file TensorMechanicsPlasticMeanCapTC.C.

449 {
450  return "MeanCapTC";
451 }

◆ numberSurfaces()

unsigned TensorMechanicsPlasticModel::numberSurfaces ( ) const
virtualinherited

◆ returnMap()

bool TensorMechanicsPlasticMeanCapTC::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
overridevirtual

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 from TensorMechanicsPlasticModel.

Definition at line 297 of file TensorMechanicsPlasticMeanCapTC.C.

307 {
308  if (!(_use_custom_returnMap))
309  return TensorMechanicsPlasticModel::returnMap(trial_stress,
310  intnl_old,
311  E_ijkl,
312  ep_plastic_tolerance,
313  returned_stress,
314  returned_intnl,
315  dpm,
316  delta_dp,
317  yf,
318  trial_stress_inadmissible);
319 
320  yf.resize(1);
321 
322  Real yf_orig = yieldFunction(trial_stress, intnl_old);
323 
324  yf[0] = yf_orig;
325 
326  if (yf_orig < _f_tol)
327  {
328  // the trial_stress is admissible
329  trial_stress_inadmissible = false;
330  return true;
331  }
332 
333  trial_stress_inadmissible = true;
334 
335  // In the following we want to solve
336  // trial_stress - stress = E_ijkl * dpm * r ...... (1)
337  // and either
338  // stress.trace() = tensile_strength(intnl) ...... (2a)
339  // intnl = intnl_old + dpm ...... (3a)
340  // or
341  // stress.trace() = compressive_strength(intnl) ... (2b)
342  // intnl = intnl_old - dpm ...... (3b)
343  // The former (2a and 3a) are chosen if
344  // trial_stress.trace() > tensile_strength(intnl_old)
345  // while the latter (2b and 3b) are chosen if
346  // trial_stress.trace() < compressive_strength(intnl_old)
347  // The variables we want to solve for are stress, dpm
348  // and intnl. We do this using a Newton approach, starting
349  // with stress=trial_stress and intnl=intnl_old and dpm=0
350  const bool tensile_failure = (trial_stress.trace() >= tensile_strength(intnl_old));
351  const Real dirn = (tensile_failure ? 1.0 : -1.0);
352 
353  RankTwoTensor n; // flow direction, which is E_ijkl * r
354  for (unsigned i = 0; i < 3; ++i)
355  for (unsigned j = 0; j < 3; ++j)
356  for (unsigned k = 0; k < 3; ++k)
357  n(i, j) += dirn * E_ijkl(i, j, k, k);
358  const Real n_trace = n.trace();
359 
360  // Perform a Newton-Raphson to find dpm when
361  // residual = trial_stress.trace() - tensile_strength(intnl) - dpm * n.trace() [for
362  // tensile_failure=true]
363  // or
364  // residual = trial_stress.trace() - compressive_strength(intnl) - dpm * n.trace() [for
365  // tensile_failure=false]
366  Real trial_trace = trial_stress.trace();
367  Real residual;
368  Real jac;
369  dpm[0] = 0;
370  unsigned int iter = 0;
371  do
372  {
373  if (tensile_failure)
374  {
375  residual = trial_trace - tensile_strength(intnl_old + dpm[0]) - dpm[0] * n_trace;
376  jac = -dtensile_strength(intnl_old + dpm[0]) - n_trace;
377  }
378  else
379  {
380  residual = trial_trace - compressive_strength(intnl_old - dpm[0]) - dpm[0] * n_trace;
381  jac = -dcompressive_strength(intnl_old - dpm[0]) - n_trace;
382  }
383  dpm[0] += -residual / jac;
384  if (iter > _max_iters) // not converging
385  return false;
386  iter++;
387  } while (residual * residual > _f_tol * _f_tol);
388 
389  // set the returned values
390  yf[0] = 0;
391  returned_intnl = intnl_old + dirn * dpm[0];
392  returned_stress = trial_stress - dpm[0] * n;
393  delta_dp = dpm[0] * dirn * returned_stress.dtrace();
394 
395  return true;
396 }

◆ tensile_strength()

Real TensorMechanicsPlasticMeanCapTC::tensile_strength ( const Real  internal_param) const
protectedvirtual

tensile strength as a function of residual value, rate, and internal_param

Definition at line 223 of file TensorMechanicsPlasticMeanCapTC.C.

224 {
225  return _strength.value(internal_param);
226 }

Referenced by activeConstraints(), consistentTangentOperator(), df_dsig(), dflowPotential_dintnl(), dflowPotential_dstress(), dhardPotential_dintnl(), dhardPotential_dstress(), dyieldFunction_dintnl(), hardPotential(), returnMap(), and yieldFunction().

◆ useCustomCTO()

bool TensorMechanicsPlasticMeanCapTC::useCustomCTO ( ) const
overridevirtual

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

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 442 of file TensorMechanicsPlasticMeanCapTC.C.

443 {
444  return _use_custom_cto;
445 }

◆ useCustomReturnMap()

bool TensorMechanicsPlasticMeanCapTC::useCustomReturnMap ( ) const
overridevirtual

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

Reimplemented from TensorMechanicsPlasticModel.

Definition at line 436 of file TensorMechanicsPlasticMeanCapTC.C.

437 {
438  return _use_custom_returnMap;
439 }

◆ validParams()

InputParameters TensorMechanicsPlasticMeanCapTC::validParams ( )
static

Definition at line 19 of file TensorMechanicsPlasticMeanCapTC.C.

20 {
21  InputParameters params = TensorMechanicsPlasticModel::validParams();
22  params.addRangeCheckedParam<unsigned>("max_iterations",
23  10,
24  "max_iterations>0",
25  "Maximum iterations for custom MeanCapTC return map");
26  params.addParam<bool>(
27  "use_custom_returnMap", true, "Whether to use the custom MeanCapTC returnMap algorithm.");
28  params.addParam<bool>("use_custom_cto",
29  true,
30  "Whether to use the custom consistent tangent operator computations.");
31  params.addRequiredParam<UserObjectName>("tensile_strength",
32  "A TensorMechanicsHardening UserObject that defines "
33  "hardening of the mean-cap tensile strength (it will "
34  "typically be positive). Yield function = trace(stress) "
35  "- tensile_strength for trace(stress)>tensile_strength.");
36  params.addRequiredParam<UserObjectName>(
37  "compressive_strength",
38  "A TensorMechanicsHardening UserObject that defines hardening of the "
39  "mean-cap compressive strength. This should always be less than "
40  "tensile_strength (it will typically be negative). Yield function = "
41  "- (trace(stress) - compressive_strength) for "
42  "trace(stress)<compressive_strength.");
43  params.addClassDescription(
44  "Associative mean-cap tensile and compressive plasticity with hardening/softening");
45 
46  return params;
47 }

◆ yieldFunction()

Real TensorMechanicsPlasticMeanCapTC::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 65 of file TensorMechanicsPlasticMeanCapTC.C.

66 {
67  const Real tr = stress.trace();
68  const Real t_str = tensile_strength(intnl);
69  if (tr >= t_str)
70  return tr - t_str;
71 
72  const Real c_str = compressive_strength(intnl);
73  if (tr <= c_str)
74  return -(tr - c_str);
75  // the following is zero at tr = t_str, and at tr = c_str
76  // it also has derivative = 1 at tr = t_str, and derivative = -1 at tr = c_str
77  // it also has second derivative = 0, at these points.
78  // This makes the complete yield function C2 continuous.
79  return (c_str - t_str) / libMesh::pi * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str));
80 }

Referenced by returnMap().

◆ 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

◆ _c_strength

const TensorMechanicsHardeningModel& TensorMechanicsPlasticMeanCapTC::_c_strength
protected

the compressive strength

Definition at line 82 of file TensorMechanicsPlasticMeanCapTC.h.

Referenced by compressive_strength(), dcompressive_strength(), and TensorMechanicsPlasticMeanCapTC().

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

◆ _max_iters

const unsigned TensorMechanicsPlasticMeanCapTC::_max_iters
protected

max iters for custom return map loop

Definition at line 70 of file TensorMechanicsPlasticMeanCapTC.h.

Referenced by returnMap().

◆ _strength

const TensorMechanicsHardeningModel& TensorMechanicsPlasticMeanCapTC::_strength
protected

the tensile strength

Definition at line 79 of file TensorMechanicsPlasticMeanCapTC.h.

Referenced by dtensile_strength(), tensile_strength(), and TensorMechanicsPlasticMeanCapTC().

◆ _use_custom_cto

const bool TensorMechanicsPlasticMeanCapTC::_use_custom_cto
protected

Whether to use the custom consistent tangent operator algorithm.

Definition at line 76 of file TensorMechanicsPlasticMeanCapTC.h.

Referenced by consistentTangentOperator(), and useCustomCTO().

◆ _use_custom_returnMap

const bool TensorMechanicsPlasticMeanCapTC::_use_custom_returnMap
protected

Whether to use the custom return-map algorithm.

Definition at line 73 of file TensorMechanicsPlasticMeanCapTC.h.

Referenced by returnMap(), and useCustomReturnMap().


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
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
TensorMechanicsPlasticMeanCapTC::_use_custom_returnMap
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.
Definition: TensorMechanicsPlasticMeanCapTC.h:73
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
TensorMechanicsPlasticMeanCapTC::_strength
const TensorMechanicsHardeningModel & _strength
the tensile strength
Definition: TensorMechanicsPlasticMeanCapTC.h:79
TensorMechanicsPlasticMeanCapTC::_c_strength
const TensorMechanicsHardeningModel & _c_strength
the compressive strength
Definition: TensorMechanicsPlasticMeanCapTC.h:82
TensorMechanicsPlasticMeanCapTC::_max_iters
const unsigned _max_iters
max iters for custom return map loop
Definition: TensorMechanicsPlasticMeanCapTC.h:70
TensorMechanicsPlasticModel::_f_tol
const Real _f_tol
Tolerance on yield function.
Definition: TensorMechanicsPlasticModel.h:175
TensorMechanicsPlasticMeanCapTC::df_dsig
RankTwoTensor df_dsig(const RankTwoTensor &stress, Real intnl) const
Derivative of the yield function with respect to stress.
Definition: TensorMechanicsPlasticMeanCapTC.C:110
TensorMechanicsPlasticModel::TensorMechanicsPlasticModel
TensorMechanicsPlasticModel(const InputParameters &parameters)
Definition: TensorMechanicsPlasticModel.C:34
TensorMechanicsPlasticMeanCapTC::compressive_strength
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
Definition: TensorMechanicsPlasticMeanCapTC.C:235
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
TensorMechanicsPlasticMeanCapTC::_use_custom_cto
const bool _use_custom_cto
Whether to use the custom consistent tangent operator algorithm.
Definition: TensorMechanicsPlasticMeanCapTC.h:76
TensorMechanicsPlasticModel::returnMap
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.
Definition: TensorMechanicsPlasticModel.C:221
TensorMechanicsHardeningModel::value
virtual Real value(Real intnl) const
Definition: TensorMechanicsHardeningModel.C:45
RankTwoTensor
RankTwoTensorTempl< Real > RankTwoTensor
Definition: ACGrGrElasticDrivingForce.h:17
TensorMechanicsPlasticMeanCapTC::yieldFunction
Real yieldFunction(const RankTwoTensor &stress, Real intnl) const override
The following functions are what you should override when building single-plasticity models.
Definition: TensorMechanicsPlasticMeanCapTC.C:65
TensorMechanicsPlasticModel::consistentTangentOperator
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.
Definition: TensorMechanicsPlasticModel.C:254
TensorMechanicsPlasticMeanCapTC::dtensile_strength
virtual Real dtensile_strength(const Real internal_param) const
d(tensile strength)/d(internal_param) as a function of residual value, rate, and internal_param
Definition: TensorMechanicsPlasticMeanCapTC.C:229
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
RankTwoTensorTempl< Real >
TensorMechanicsPlasticMeanCapTC::tensile_strength
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
Definition: TensorMechanicsPlasticMeanCapTC.C:223
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
RankFourTensor
RankFourTensorTempl< Real > RankFourTensor
Definition: ACGrGrElasticDrivingForce.h:20
TensorMechanicsPlasticMeanCapTC::dcompressive_strength
virtual Real dcompressive_strength(const Real internal_param) const
d(compressive strength)/d(internal_param) as a function of residual value, rate, and internal_param
Definition: TensorMechanicsPlasticMeanCapTC.C:241
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