www.mooseframework.org
TensorMechanicsPlasticMeanCapTC.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "RankFourTensor.h"
12 #include "libmesh/utility.h"
13 
15 
17 
18 InputParameters
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 }
48 
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 }
63 
64 Real
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 }
81 
84  Real intnl) const
85 {
86  return df_dsig(stress, intnl);
87 }
88 
89 Real
91  Real intnl) const
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 }
108 
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 }
123 
126 {
127  return df_dsig(stress, intnl);
128 }
129 
132  Real intnl) const
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 }
146 
149  Real intnl) const
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 }
165 
166 Real
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 }
185 
188  Real intnl) const
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 }
202 
203 Real
205  Real intnl) const
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 }
221 
222 Real
224 {
225  return _strength.value(internal_param);
226 }
227 
228 Real
230 {
231  return _strength.derivative(internal_param);
232 }
233 
234 Real
236 {
237  return _c_strength.value(internal_param);
238 }
239 
240 Real
242 {
243  return _c_strength.derivative(internal_param);
244 }
245 
246 void
248  const RankTwoTensor & stress,
249  Real intnl,
250  const RankFourTensor & Eijkl,
251  std::vector<bool> & act,
252  RankTwoTensor & returned_stress) const
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 }
295 
296 bool
298  Real intnl_old,
299  const RankFourTensor & E_ijkl,
300  Real ep_plastic_tolerance,
301  RankTwoTensor & returned_stress,
302  Real & returned_intnl,
303  std::vector<Real> & dpm,
304  RankTwoTensor & delta_dp,
305  std::vector<Real> & yf,
306  bool & trial_stress_inadmissible) const
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 }
397 
400  const RankTwoTensor & trial_stress,
401  Real intnl_old,
402  const RankTwoTensor & stress,
403  Real intnl,
404  const RankFourTensor & E_ijkl,
405  const std::vector<Real> & cumulative_pm) const
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 }
434 
435 bool
437 {
438  return _use_custom_returnMap;
439 }
440 
441 bool
443 {
444  return _use_custom_cto;
445 }
446 
447 std::string
449 {
450  return "MeanCapTC";
451 }
TensorMechanicsPlasticMeanCapTC::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 override
Performs a custom return-map.
Definition: TensorMechanicsPlasticMeanCapTC.C:297
TensorMechanicsPlasticMeanCapTC::TensorMechanicsPlasticMeanCapTC
TensorMechanicsPlasticMeanCapTC(const InputParameters &parameters)
Definition: TensorMechanicsPlasticMeanCapTC.C:49
TensorMechanicsPlasticModel::validParams
static InputParameters validParams()
Definition: TensorMechanicsPlasticModel.C:18
TensorMechanicsPlasticMeanCapTC.h
TensorMechanicsPlasticMeanCapTC::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 override
Calculates a custom consistent tangent operator.
Definition: TensorMechanicsPlasticMeanCapTC.C:399
TensorMechanicsPlasticMeanCapTC::_use_custom_returnMap
const bool _use_custom_returnMap
Whether to use the custom return-map algorithm.
Definition: TensorMechanicsPlasticMeanCapTC.h:73
TensorMechanicsPlasticMeanCapTC
Rate-independent associative mean-cap tensile AND compressive failure with hardening/softening of the...
Definition: TensorMechanicsPlasticMeanCapTC.h:29
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::dyieldFunction_dstress
RankTwoTensor dyieldFunction_dstress(const RankTwoTensor &stress, Real intnl) const override
The derivative of yield function with respect to stress.
Definition: TensorMechanicsPlasticMeanCapTC.C:83
TensorMechanicsPlasticMeanCapTC::hardPotential
Real hardPotential(const RankTwoTensor &stress, Real intnl) const override
The hardening potential.
Definition: TensorMechanicsPlasticMeanCapTC.C:167
registerMooseObject
registerMooseObject("TensorMechanicsApp", TensorMechanicsPlasticMeanCapTC)
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
TensorMechanicsPlasticMeanCapTC::useCustomCTO
virtual bool useCustomCTO() const override
Returns false. You will want to override this in your derived class if you write a custom consistent ...
Definition: TensorMechanicsPlasticMeanCapTC.C:442
TensorMechanicsPlasticMeanCapTC::dhardPotential_dstress
virtual RankTwoTensor dhardPotential_dstress(const RankTwoTensor &stress, Real intnl) const override
The derivative of the hardening potential with respect to stress.
Definition: TensorMechanicsPlasticMeanCapTC.C:187
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
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::modelName
virtual std::string modelName() const override
Definition: TensorMechanicsPlasticMeanCapTC.C:448
TensorMechanicsPlasticMeanCapTC::dhardPotential_dintnl
virtual Real dhardPotential_dintnl(const RankTwoTensor &stress, Real intnl) const override
The derivative of the hardening potential with respect to the internal parameter.
Definition: TensorMechanicsPlasticMeanCapTC.C:204
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::dflowPotential_dintnl
RankTwoTensor dflowPotential_dintnl(const RankTwoTensor &stress, Real intnl) const override
The derivative of the flow potential with respect to the internal parameter.
Definition: TensorMechanicsPlasticMeanCapTC.C:148
defineLegacyParams
defineLegacyParams(TensorMechanicsPlasticMeanCapTC)
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
TensorMechanicsPlasticMeanCapTC::validParams
static InputParameters validParams()
Definition: TensorMechanicsPlasticMeanCapTC.C:19
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
Plastic Model base class The virtual functions written below must be over-ridden in derived classes t...
Definition: TensorMechanicsPlasticModel.h:42
TensorMechanicsPlasticMeanCapTC::flowPotential
RankTwoTensor flowPotential(const RankTwoTensor &stress, Real intnl) const override
The flow potential.
Definition: TensorMechanicsPlasticMeanCapTC.C:125
TensorMechanicsPlasticMeanCapTC::activeConstraints
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.
Definition: TensorMechanicsPlasticMeanCapTC.C:247
TensorMechanicsPlasticMeanCapTC::dflowPotential_dstress
RankFourTensor dflowPotential_dstress(const RankTwoTensor &stress, Real intnl) const override
The derivative of the flow potential with respect to stress.
Definition: TensorMechanicsPlasticMeanCapTC.C:131
TensorMechanicsPlasticMeanCapTC::dyieldFunction_dintnl
Real dyieldFunction_dintnl(const RankTwoTensor &stress, Real intnl) const override
The derivative of yield function with respect to the internal parameter.
Definition: TensorMechanicsPlasticMeanCapTC.C:90
TensorMechanicsPlasticMeanCapTC::useCustomReturnMap
virtual bool useCustomReturnMap() const override
Returns false. You will want to override this in your derived class if you write a custom returnMap f...
Definition: TensorMechanicsPlasticMeanCapTC.C:436
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
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
TensorMechanicsHardeningModel
Hardening Model base class.
Definition: TensorMechanicsHardeningModel.h:27