www.mooseframework.org
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...
 

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 27 of file TensorMechanicsPlasticMeanCapTC.h.

Constructor & Destructor Documentation

TensorMechanicsPlasticMeanCapTC::TensorMechanicsPlasticMeanCapTC ( const InputParameters &  parameters)

Definition at line 42 of file TensorMechanicsPlasticMeanCapTC.C.

43  : TensorMechanicsPlasticModel(parameters),
44  _max_iters(getParam<unsigned>("max_iterations")),
45  _use_custom_returnMap(getParam<bool>("use_custom_returnMap")),
46  _use_custom_cto(getParam<bool>("use_custom_cto")),
47  _strength(getUserObject<TensorMechanicsHardeningModel>("tensile_strength")),
48  _c_strength(getUserObject<TensorMechanicsHardeningModel>("compressive_strength"))
49 {
50  // cannot check the following for all values of the internal parameter, but this will catch most
51  // errors
52  if (_strength.value(0) <= _c_strength.value(0))
53  mooseError("MeanCapTC: tensile strength (which is usually positive) must not be less than "
54  "compressive strength (which is usually negative)");
55 }
TensorMechanicsPlasticModel(const InputParameters &parameters)
const TensorMechanicsHardeningModel & _c_strength
the compressive strength
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.
const unsigned _max_iters
max iters for custom return map loop
const bool _use_custom_cto
Whether to use the custom consistent tangent operator algorithm.
virtual Real value(Real intnl) const
const TensorMechanicsHardeningModel & _strength
the tensile strength

Member Function Documentation

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 240 of file TensorMechanicsPlasticMeanCapTC.C.

246 {
247  act.assign(1, false);
248 
249  if (f[0] <= _f_tol)
250  {
251  returned_stress = stress;
252  return;
253  }
254 
255  const Real tr = stress.trace();
256  const Real t_str = tensile_strength(intnl);
257  Real str;
258  Real dirn;
259  if (tr >= t_str)
260  {
261  str = t_str;
262  dirn = 1.0;
263  }
264  else
265  {
266  str = compressive_strength(intnl);
267  dirn = -1.0;
268  }
269 
270  RankTwoTensor n; // flow direction
271  for (unsigned i = 0; i < 3; ++i)
272  for (unsigned j = 0; j < 3; ++j)
273  for (unsigned k = 0; k < 3; ++k)
274  n(i, j) += dirn * Eijkl(i, j, k, k);
275 
276  // returned_stress = stress - gamma*n
277  // and taking the trace of this and using
278  // Tr(returned_stress) = str, gives
279  // gamma = (Tr(stress) - str)/Tr(n)
280  Real gamma = (stress.trace() - str) / n.trace();
281 
282  for (unsigned i = 0; i < 3; ++i)
283  for (unsigned j = 0; j < 3; ++j)
284  returned_stress(i, j) = stress(i, j) - gamma * n(i, j);
285 
286  act[0] = true;
287 }
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
const Real _f_tol
Tolerance on yield function.
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 228 of file TensorMechanicsPlasticMeanCapTC.C.

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

229 {
230  return _c_strength.value(internal_param);
231 }
const TensorMechanicsHardeningModel & _c_strength
the compressive strength
virtual Real value(Real intnl) const
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 392 of file TensorMechanicsPlasticMeanCapTC.C.

399 {
400  if (!_use_custom_cto)
402  trial_stress, intnl_old, stress, intnl, E_ijkl, cumulative_pm);
403 
404  Real df_dq;
405  Real alpha;
406  if (trial_stress.trace() >= tensile_strength(intnl_old))
407  {
408  df_dq = -dtensile_strength(intnl);
409  alpha = 1.0;
410  }
411  else
412  {
413  df_dq = dcompressive_strength(intnl);
414  alpha = -1.0;
415  }
416 
417  RankTwoTensor elas;
418  for (unsigned int i = 0; i < 3; ++i)
419  for (unsigned int j = 0; j < 3; ++j)
420  for (unsigned int k = 0; k < 3; ++k)
421  elas(i, j) += E_ijkl(i, j, k, k);
422 
423  const Real hw = -df_dq + alpha * elas.trace();
424 
425  return E_ijkl - alpha / hw * elas.outerProduct(elas);
426 }
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 ...
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.
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
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 ...
const bool _use_custom_cto
Whether to use the custom consistent tangent operator algorithm.
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 234 of file TensorMechanicsPlasticMeanCapTC.C.

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

235 {
236  return _c_strength.derivative(internal_param);
237 }
const TensorMechanicsHardeningModel & _c_strength
the compressive strength
virtual Real derivative(Real intnl) const
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 103 of file TensorMechanicsPlasticMeanCapTC.C.

Referenced by dyieldFunction_dstress(), and flowPotential().

104 {
105  const Real tr = stress.trace();
106  const Real t_str = tensile_strength(intnl);
107  if (tr >= t_str)
108  return stress.dtrace();
109 
110  const Real c_str = compressive_strength(intnl);
111  if (tr <= c_str)
112  return -stress.dtrace();
113 
114  return -std::cos(libMesh::pi * (tr - c_str) / (t_str - c_str)) * stress.dtrace();
115 }
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 141 of file TensorMechanicsPlasticMeanCapTC.C.

143 {
144  const Real tr = stress.trace();
145  const Real t_str = tensile_strength(intnl);
146  if (tr >= t_str)
147  return RankTwoTensor();
148 
149  const Real c_str = compressive_strength(intnl);
150  if (tr <= c_str)
151  return RankTwoTensor();
152 
153  const Real dt = dtensile_strength(intnl);
154  const Real dc = dcompressive_strength(intnl);
155  return std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * stress.dtrace() * libMesh::pi /
156  Utility::pow<2>(t_str - c_str) * ((tr - t_str) * dc - (tr - c_str) * dt);
157 }
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 ...
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
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 ...
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 134 of file TensorMechanicsPlasticModel.C.

137 {
138  return dr_dintnl.assign(1, dflowPotential_dintnl(stress, intnl));
139 }
virtual RankTwoTensor dflowPotential_dintnl(const RankTwoTensor &stress, Real intnl) const
The derivative of the flow potential with respect to the internal parameter.
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 124 of file TensorMechanicsPlasticMeanCapTC.C.

126 {
127  const Real tr = stress.trace();
128  const Real t_str = tensile_strength(intnl);
129  if (tr >= t_str)
130  return RankFourTensor();
131 
132  const Real c_str = compressive_strength(intnl);
133  if (tr <= c_str)
134  return RankFourTensor();
135 
136  return libMesh::pi / (t_str - c_str) * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) *
137  stress.dtrace().outerProduct(stress.dtrace());
138 }
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 120 of file TensorMechanicsPlasticModel.C.

123 {
124  return dr_dstress.assign(1, dflowPotential_dstress(stress, intnl));
125 }
virtual RankFourTensor dflowPotential_dstress(const RankTwoTensor &stress, Real intnl) const
The derivative of the flow potential with respect to stress.
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 197 of file TensorMechanicsPlasticMeanCapTC.C.

199 {
200  const Real tr = stress.trace();
201  const Real t_str = tensile_strength(intnl);
202  if (tr >= t_str)
203  return 0.0;
204 
205  const Real c_str = compressive_strength(intnl);
206  if (tr <= c_str)
207  return 0.0;
208 
209  const Real dt = dtensile_strength(intnl);
210  const Real dc = dcompressive_strength(intnl);
211  return -std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * libMesh::pi /
212  Utility::pow<2>(t_str - c_str) * ((tr - t_str) * dc - (tr - c_str) * dt);
213 }
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 ...
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
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 ...
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 175 of file TensorMechanicsPlasticModel.C.

178 {
179  dh_dintnl.resize(numberSurfaces(), dhardPotential_dintnl(stress, intnl));
180 }
virtual unsigned int numberSurfaces() const
The number of yield surfaces for this plasticity model.
virtual Real dhardPotential_dintnl(const RankTwoTensor &stress, Real intnl) const
The derivative of the hardening potential with respect to the internal parameter. ...
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 180 of file TensorMechanicsPlasticMeanCapTC.C.

182 {
183  const Real tr = stress.trace();
184  const Real t_str = tensile_strength(intnl);
185  if (tr >= t_str)
186  return RankTwoTensor();
187 
188  const Real c_str = compressive_strength(intnl);
189  if (tr <= c_str)
190  return RankTwoTensor();
191 
192  return -std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) * libMesh::pi / (t_str - c_str) *
193  stress.dtrace();
194 }
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 161 of file TensorMechanicsPlasticModel.C.

164 {
165  dh_dstress.assign(numberSurfaces(), dhardPotential_dstress(stress, intnl));
166 }
virtual unsigned int numberSurfaces() const
The number of yield surfaces for this plasticity model.
virtual RankTwoTensor dhardPotential_dstress(const RankTwoTensor &stress, Real intnl) const
The derivative of the hardening potential with respect to stress.
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 222 of file TensorMechanicsPlasticMeanCapTC.C.

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

223 {
224  return _strength.derivative(internal_param);
225 }
virtual Real derivative(Real intnl) const
const TensorMechanicsHardeningModel & _strength
the tensile strength
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 83 of file TensorMechanicsPlasticMeanCapTC.C.

85 {
86  const Real tr = stress.trace();
87  const Real t_str = tensile_strength(intnl);
88  if (tr >= t_str)
89  return -dtensile_strength(intnl);
90 
91  const Real c_str = compressive_strength(intnl);
92  if (tr <= c_str)
93  return dcompressive_strength(intnl);
94 
95  const Real dt = dtensile_strength(intnl);
96  const Real dc = dcompressive_strength(intnl);
97  return (dc - dt) / libMesh::pi * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str)) +
98  1.0 / (t_str - c_str) * std::cos(libMesh::pi * (tr - c_str) / (t_str - c_str)) *
99  ((tr - c_str) * dt - (tr - t_str) * dc);
100 }
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 ...
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
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 ...
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 93 of file TensorMechanicsPlasticModel.C.

96 {
97  return df_dintnl.assign(1, dyieldFunction_dintnl(stress, intnl));
98 }
virtual Real dyieldFunction_dintnl(const RankTwoTensor &stress, Real intnl) const
The derivative of yield function with respect to the internal parameter.
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 76 of file TensorMechanicsPlasticMeanCapTC.C.

78 {
79  return df_dsig(stress, intnl);
80 }
RankTwoTensor df_dsig(const RankTwoTensor &stress, Real intnl) const
Derivative of the yield function with respect to stress.
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 79 of file TensorMechanicsPlasticModel.C.

82 {
83  df_dstress.assign(1, dyieldFunction_dstress(stress, intnl));
84 }
virtual RankTwoTensor dyieldFunction_dstress(const RankTwoTensor &stress, Real intnl) const
The derivative of yield function with respect to stress.
void TensorMechanicsPlasticModel::execute ( )
inherited

Definition at line 42 of file TensorMechanicsPlasticModel.C.

43 {
44 }
void TensorMechanicsPlasticModel::finalize ( )
inherited

Definition at line 47 of file TensorMechanicsPlasticModel.C.

48 {
49 }
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 118 of file TensorMechanicsPlasticMeanCapTC.C.

119 {
120  return df_dsig(stress, intnl);
121 }
RankTwoTensor df_dsig(const RankTwoTensor &stress, Real intnl) const
Derivative of the yield function with respect to stress.
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 106 of file TensorMechanicsPlasticModel.C.

109 {
110  return r.assign(1, flowPotential(stress, intnl));
111 }
virtual RankTwoTensor flowPotential(const RankTwoTensor &stress, Real intnl) const
The flow potential.
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 160 of file TensorMechanicsPlasticMeanCapTC.C.

161 {
162  // This is the key for this whole class!
163  const Real tr = stress.trace();
164  const Real t_str = tensile_strength(intnl);
165 
166  if (tr >= t_str)
167  return -1.0; // this will serve to *increase* the internal parameter (so internal parameter will
168  // be a measure of volumetric strain)
169 
170  const Real c_str = compressive_strength(intnl);
171  if (tr <= c_str)
172  return 1.0; // this will serve to *decrease* the internal parameter (so internal parameter will
173  // be a measure of volumetric strain)
174 
175  return std::cos(libMesh::pi * (tr - c_str) /
176  (t_str - c_str)); // this interpolates C2 smoothly between 1 and -1
177 }
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 147 of file TensorMechanicsPlasticModel.C.

150 {
151  h.assign(numberSurfaces(), hardPotential(stress, intnl));
152 }
virtual Real hardPotential(const RankTwoTensor &stress, Real intnl) const
The hardening potential.
virtual unsigned int numberSurfaces() const
The number of yield surfaces for this plasticity model.
void TensorMechanicsPlasticModel::initialize ( )
inherited

Definition at line 37 of file TensorMechanicsPlasticModel.C.

38 {
39 }
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 243 of file TensorMechanicsPlasticModel.C.

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

244 {
245  return (dpm == 0 && yf <= _f_tol) || (dpm > -dpm_tol && yf <= _f_tol && yf >= -_f_tol);
246 }
const Real _f_tol
Tolerance on yield function.
std::string TensorMechanicsPlasticMeanCapTC::modelName ( ) const
overridevirtual

Implements TensorMechanicsPlasticModel.

Definition at line 441 of file TensorMechanicsPlasticMeanCapTC.C.

442 {
443  return "MeanCapTC";
444 }
unsigned TensorMechanicsPlasticModel::numberSurfaces ( ) const
virtualinherited
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 290 of file TensorMechanicsPlasticMeanCapTC.C.

300 {
301  if (!(_use_custom_returnMap))
302  return TensorMechanicsPlasticModel::returnMap(trial_stress,
303  intnl_old,
304  E_ijkl,
305  ep_plastic_tolerance,
306  returned_stress,
307  returned_intnl,
308  dpm,
309  delta_dp,
310  yf,
311  trial_stress_inadmissible);
312 
313  yf.resize(1);
314 
315  Real yf_orig = yieldFunction(trial_stress, intnl_old);
316 
317  yf[0] = yf_orig;
318 
319  if (yf_orig < _f_tol)
320  {
321  // the trial_stress is admissible
322  trial_stress_inadmissible = false;
323  return true;
324  }
325 
326  trial_stress_inadmissible = true;
327 
328  // In the following we want to solve
329  // trial_stress - stress = E_ijkl * dpm * r ...... (1)
330  // and either
331  // stress.trace() = tensile_strength(intnl) ...... (2a)
332  // intnl = intnl_old + dpm ...... (3a)
333  // or
334  // stress.trace() = compressive_strength(intnl) ... (2b)
335  // intnl = intnl_old - dpm ...... (3b)
336  // The former (2a and 3a) are chosen if
337  // trial_stress.trace() > tensile_strength(intnl_old)
338  // while the latter (2b and 3b) are chosen if
339  // trial_stress.trace() < compressive_strength(intnl_old)
340  // The variables we want to solve for are stress, dpm
341  // and intnl. We do this using a Newton approach, starting
342  // with stress=trial_stress and intnl=intnl_old and dpm=0
343  const bool tensile_failure = (trial_stress.trace() >= tensile_strength(intnl_old));
344  const Real dirn = (tensile_failure ? 1.0 : -1.0);
345 
346  RankTwoTensor n; // flow direction, which is E_ijkl * r
347  for (unsigned i = 0; i < 3; ++i)
348  for (unsigned j = 0; j < 3; ++j)
349  for (unsigned k = 0; k < 3; ++k)
350  n(i, j) += dirn * E_ijkl(i, j, k, k);
351  const Real n_trace = n.trace();
352 
353  // Perform a Newton-Raphson to find dpm when
354  // residual = trial_stress.trace() - tensile_strength(intnl) - dpm * n.trace() [for
355  // tensile_failure=true]
356  // or
357  // residual = trial_stress.trace() - compressive_strength(intnl) - dpm * n.trace() [for
358  // tensile_failure=false]
359  Real trial_trace = trial_stress.trace();
360  Real residual;
361  Real jac;
362  dpm[0] = 0;
363  unsigned int iter = 0;
364  do
365  {
366  if (tensile_failure)
367  {
368  residual = trial_trace - tensile_strength(intnl_old + dpm[0]) - dpm[0] * n_trace;
369  jac = -dtensile_strength(intnl_old + dpm[0]) - n_trace;
370  }
371  else
372  {
373  residual = trial_trace - compressive_strength(intnl_old - dpm[0]) - dpm[0] * n_trace;
374  jac = -dcompressive_strength(intnl_old - dpm[0]) - n_trace;
375  }
376  dpm[0] += -residual / jac;
377  if (iter > _max_iters) // not converging
378  return false;
379  iter++;
380  } while (residual * residual > _f_tol * _f_tol);
381 
382  // set the returned values
383  yf[0] = 0;
384  returned_intnl = intnl_old + dirn * dpm[0];
385  returned_stress = trial_stress - dpm[0] * n;
386  delta_dp = dpm[0] * dirn * returned_stress.dtrace();
387 
388  return true;
389 }
Real yieldFunction(const RankTwoTensor &stress, Real intnl) const override
The following functions are what you should override when building single-plasticity models...
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 ...
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
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.
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.
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 ...
const unsigned _max_iters
max iters for custom return map loop
const Real _f_tol
Tolerance on yield function.
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 216 of file TensorMechanicsPlasticMeanCapTC.C.

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

217 {
218  return _strength.value(internal_param);
219 }
virtual Real value(Real intnl) const
const TensorMechanicsHardeningModel & _strength
the tensile strength
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 435 of file TensorMechanicsPlasticMeanCapTC.C.

436 {
437  return _use_custom_cto;
438 }
const bool _use_custom_cto
Whether to use the custom consistent tangent operator algorithm.
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 429 of file TensorMechanicsPlasticMeanCapTC.C.

430 {
431  return _use_custom_returnMap;
432 }
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.
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 58 of file TensorMechanicsPlasticMeanCapTC.C.

Referenced by returnMap().

59 {
60  const Real tr = stress.trace();
61  const Real t_str = tensile_strength(intnl);
62  if (tr >= t_str)
63  return tr - t_str;
64 
65  const Real c_str = compressive_strength(intnl);
66  if (tr <= c_str)
67  return -(tr - c_str);
68  // the following is zero at tr = t_str, and at tr = c_str
69  // it also has derivative = 1 at tr = t_str, and derivative = -1 at tr = c_str
70  // it also has second derivative = 0, at these points.
71  // This makes the complete yield function C2 continuous.
72  return (c_str - t_str) / libMesh::pi * std::sin(libMesh::pi * (tr - c_str) / (t_str - c_str));
73 }
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
virtual Real compressive_strength(const Real internal_param) const
compressive strength as a function of residual value, rate, and internal_param
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 64 of file TensorMechanicsPlasticModel.C.

Referenced by TensorMechanicsPlasticModel::returnMap().

67 {
68  f.assign(1, yieldFunction(stress, intnl));
69 }
virtual Real yieldFunction(const RankTwoTensor &stress, Real intnl) const
The following functions are what you should override when building single-plasticity models...

Member Data Documentation

const TensorMechanicsHardeningModel& TensorMechanicsPlasticMeanCapTC::_c_strength
protected

the compressive strength

Definition at line 78 of file TensorMechanicsPlasticMeanCapTC.h.

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

const Real TensorMechanicsPlasticModel::_f_tol
inherited
const Real TensorMechanicsPlasticModel::_ic_tol
inherited

Tolerance on internal constraint.

Definition at line 174 of file TensorMechanicsPlasticModel.h.

const unsigned TensorMechanicsPlasticMeanCapTC::_max_iters
protected

max iters for custom return map loop

Definition at line 66 of file TensorMechanicsPlasticMeanCapTC.h.

Referenced by returnMap().

const TensorMechanicsHardeningModel& TensorMechanicsPlasticMeanCapTC::_strength
protected

the tensile strength

Definition at line 75 of file TensorMechanicsPlasticMeanCapTC.h.

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

const bool TensorMechanicsPlasticMeanCapTC::_use_custom_cto
protected

Whether to use the custom consistent tangent operator algorithm.

Definition at line 72 of file TensorMechanicsPlasticMeanCapTC.h.

Referenced by consistentTangentOperator(), and useCustomCTO().

const bool TensorMechanicsPlasticMeanCapTC::_use_custom_returnMap
protected

Whether to use the custom return-map algorithm.

Definition at line 69 of file TensorMechanicsPlasticMeanCapTC.h.

Referenced by returnMap(), and useCustomReturnMap().


The documentation for this class was generated from the following files: