www.mooseframework.org
TensileStressUpdate.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 
10 #include "TensileStressUpdate.h"
11 #include "libmesh/utility.h"
12 
13 registerMooseObject("TensorMechanicsApp", TensileStressUpdate);
14 
16 
17 InputParameters
19 {
20  InputParameters params = MultiParameterPlasticityStressUpdate::validParams();
21  params.addRequiredParam<UserObjectName>(
22  "tensile_strength",
23  "A TensorMechanicsHardening UserObject that defines hardening of the tensile strength");
24  params.addParam<bool>("perfect_guess",
25  true,
26  "Provide a guess to the Newton-Raphson procedure "
27  "that is the result from perfect plasticity. With "
28  "severe hardening/softening this may be "
29  "suboptimal.");
30  params.addClassDescription(
31  "Associative, smoothed, tensile (Rankine) plasticity with hardening/softening");
32  return params;
33 }
34 
35 TensileStressUpdate::TensileStressUpdate(const InputParameters & parameters)
36  : MultiParameterPlasticityStressUpdate(parameters, 3, 3, 1),
37  _strength(getUserObject<TensorMechanicsHardeningModel>("tensile_strength")),
38  _perfect_guess(getParam<bool>("perfect_guess")),
39  _eigvecs(RankTwoTensor())
40 {
41 }
42 
43 void
45  std::vector<Real> & stress_params) const
46 {
47  // stress_params[0] = smallest eigenvalue, stress_params[2] = largest eigenvalue
48  stress.symmetricEigenvalues(stress_params);
49 }
50 
51 std::vector<RankTwoTensor>
53 {
54  std::vector<Real> sp;
55  std::vector<RankTwoTensor> dsp;
56  stress.dsymmetricEigenvalues(sp, dsp);
57  return dsp;
58 }
59 
60 std::vector<RankFourTensor>
62 {
63  std::vector<RankFourTensor> d2;
64  stress.d2symmetricEigenvalues(d2);
65  return d2;
66 }
67 
68 void
69 TensileStressUpdate::preReturnMapV(const std::vector<Real> & /*trial_stress_params*/,
70  const RankTwoTensor & stress_trial,
71  const std::vector<Real> & /*intnl_old*/,
72  const std::vector<Real> & /*yf*/,
73  const RankFourTensor & /*Eijkl*/)
74 {
75  std::vector<Real> eigvals;
76  stress_trial.symmetricEigenvaluesEigenvectors(eigvals, _eigvecs);
77 }
78 
79 void
81  const std::vector<Real> & stress_params,
82  Real /*gaE*/,
83  const std::vector<Real> & /*intnl*/,
84  const yieldAndFlow & /*smoothed_q*/,
85  const RankFourTensor & /*Eijkl*/,
86  RankTwoTensor & stress) const
87 {
88  // form the diagonal stress
89  stress = RankTwoTensor(stress_params[0], stress_params[1], stress_params[2], 0.0, 0.0, 0.0);
90  // rotate to the original frame
91  stress = _eigvecs * stress * (_eigvecs.transpose());
92 }
93 
94 void
95 TensileStressUpdate::yieldFunctionValuesV(const std::vector<Real> & stress_params,
96  const std::vector<Real> & intnl,
97  std::vector<Real> & yf) const
98 {
99  const Real ts = tensile_strength(intnl[0]);
100  // The smoothing strategy means that the last yield function
101  // gives the smallest value when returning to a line/point where,
102  // without smoothing, the yield functions are equal.
103  // Therefore, i use the smallest eigenvalue in this yield function
104  yf[0] = stress_params[2] - ts; // use largest eigenvalue
105  yf[1] = stress_params[1] - ts;
106  yf[2] = stress_params[0] - ts; // use smallest eigenvalue
107 }
108 
109 void
110 TensileStressUpdate::computeAllQV(const std::vector<Real> & stress_params,
111  const std::vector<Real> & intnl,
112  std::vector<yieldAndFlow> & all_q) const
113 {
114  const Real ts = tensile_strength(intnl[0]);
115  const Real dts = dtensile_strength(intnl[0]);
116 
117  // yield functions. See comment in yieldFunctionValuesV
118  all_q[0].f = stress_params[2] - ts;
119  all_q[1].f = stress_params[1] - ts;
120  all_q[2].f = stress_params[0] - ts;
121 
122  // d(yield function)/d(stress_params)
123  for (unsigned yf = 0; yf < _num_yf; ++yf)
124  for (unsigned a = 0; a < _num_sp; ++a)
125  all_q[yf].df[a] = 0.0;
126  all_q[0].df[2] = 1.0;
127  all_q[1].df[1] = 1.0;
128  all_q[2].df[0] = 1.0;
129 
130  // d(yield function)/d(intnl)
131  for (unsigned yf = 0; yf < _num_yf; ++yf)
132  all_q[yf].df_di[0] = -dts;
133 
134  // the flow potential is just the yield function
135  // d(flow potential)/d(stress_params)
136  for (unsigned yf = 0; yf < _num_yf; ++yf)
137  for (unsigned a = 0; a < _num_sp; ++a)
138  all_q[yf].dg[a] = all_q[yf].df[a];
139 
140  // d(flow potential)/d(stress_params)/d(intnl)
141  for (unsigned yf = 0; yf < _num_yf; ++yf)
142  for (unsigned a = 0; a < _num_sp; ++a)
143  all_q[yf].d2g_di[a][0] = 0.0;
144 
145  // d(flow potential)/d(stress_params)/d(stress_params)
146  for (unsigned yf = 0; yf < _num_yf; ++yf)
147  for (unsigned a = 0; a < _num_sp; ++a)
148  for (unsigned b = 0; b < _num_sp; ++b)
149  all_q[yf].d2g[a][b] = 0.0;
150 }
151 
152 void
154 {
155  // Eijkl is required to be isotropic, so we can use the
156  // frame where stress is diagonal
157  for (unsigned a = 0; a < _num_sp; ++a)
158  for (unsigned b = 0; b < _num_sp; ++b)
159  _Eij[a][b] = Eijkl(a, a, b, b);
160  _En = _Eij[2][2];
161  const Real denom = _Eij[0][0] * (_Eij[0][0] + _Eij[0][1]) - 2 * Utility::pow<2>(_Eij[0][1]);
162  for (unsigned a = 0; a < _num_sp; ++a)
163  {
164  _Cij[a][a] = (_Eij[0][0] + _Eij[0][1]) / denom;
165  for (unsigned b = 0; b < a; ++b)
166  _Cij[a][b] = _Cij[b][a] = -_Eij[0][1] / denom;
167  }
168 }
169 
170 void
171 TensileStressUpdate::initializeVarsV(const std::vector<Real> & trial_stress_params,
172  const std::vector<Real> & intnl_old,
173  std::vector<Real> & stress_params,
174  Real & gaE,
175  std::vector<Real> & intnl) const
176 {
177  if (!_perfect_guess)
178  {
179  for (unsigned i = 0; i < _num_sp; ++i)
180  stress_params[i] = trial_stress_params[i];
181  gaE = 0.0;
182  }
183  else
184  {
185  const Real ts = tensile_strength(intnl_old[0]);
186  stress_params[2] = ts; // largest eigenvalue
187  stress_params[1] = std::min(stress_params[1], ts);
188  stress_params[0] = std::min(stress_params[0], ts);
189  gaE = trial_stress_params[2] - stress_params[2];
190  }
191  setIntnlValuesV(trial_stress_params, stress_params, intnl_old, intnl);
192 }
193 
194 void
195 TensileStressUpdate::setIntnlValuesV(const std::vector<Real> & trial_stress_params,
196  const std::vector<Real> & current_stress_params,
197  const std::vector<Real> & intnl_old,
198  std::vector<Real> & intnl) const
199 {
200  intnl[0] = intnl_old[0] + (trial_stress_params[2] - current_stress_params[2]) / _Eij[2][2];
201 }
202 
203 void
204 TensileStressUpdate::setIntnlDerivativesV(const std::vector<Real> & /*trial_stress_params*/,
205  const std::vector<Real> & /*current_stress_params*/,
206  const std::vector<Real> & /*intnl*/,
207  std::vector<std::vector<Real>> & dintnl) const
208 {
209  dintnl[0][0] = 0.0;
210  dintnl[0][1] = 0.0;
211  dintnl[0][2] = -1.0 / _Eij[2][2];
212 }
213 
214 Real
215 TensileStressUpdate::tensile_strength(const Real internal_param) const
216 {
217  return _strength.value(internal_param);
218 }
219 
220 Real
221 TensileStressUpdate::dtensile_strength(const Real internal_param) const
222 {
223  return _strength.derivative(internal_param);
224 }
225 
226 void
228  const std::vector<Real> & trial_stress_params,
229  const RankTwoTensor & /*stress*/,
230  const std::vector<Real> & stress_params,
231  Real /*gaE*/,
232  const yieldAndFlow & /*smoothed_q*/,
233  const RankFourTensor & elasticity_tensor,
234  bool compute_full_tangent_operator,
235  const std::vector<std::vector<Real>> & dvar_dtrial,
236  RankFourTensor & cto)
237 {
238  cto = elasticity_tensor;
239  if (!compute_full_tangent_operator)
240  return;
241 
242  // dvar_dtrial has been computed already, so
243  // d(stress)/d(trial_stress) = d(eigvecs * stress_params * eigvecs.transpose())/d(trial_stress)
244  // eigvecs is a rotation matrix, rot(i, j) = e_j(i) = i^th component of j^th eigenvector
245  // d(rot_ij)/d(stress_kl) = d(e_j(i))/d(stress_kl)
246  // = sum_a 0.5 * e_a(i) * (e_a(k)e_j(l) + e_a(l)e_j(k)) / (la_j - la_a)
247  // = sum_a 0.5 * rot(i,a) * (rot(k,a)rot(l,j) + rot(l,a)*rot(k,j)) / (la_j - la_a)
248  RankFourTensor drot_dstress;
249  for (unsigned i = 0; i < _tensor_dimensionality; ++i)
250  for (unsigned j = 0; j < _tensor_dimensionality; ++j)
251  for (unsigned k = 0; k < _tensor_dimensionality; ++k)
252  for (unsigned l = 0; l < _tensor_dimensionality; ++l)
253  for (unsigned a = 0; a < _num_sp; ++a)
254  {
255  if (trial_stress_params[a] == trial_stress_params[j])
256  continue;
257  drot_dstress(i, j, k, l) += 0.5 * _eigvecs(i, a) * (_eigvecs(k, a) * _eigvecs(l, j) +
258  _eigvecs(l, a) * _eigvecs(k, j)) /
259  (trial_stress_params[j] - trial_stress_params[a]);
260  }
261 
262  const RankTwoTensor eT = _eigvecs.transpose();
263 
264  RankFourTensor dstress_dtrial;
265  for (unsigned i = 0; i < _tensor_dimensionality; ++i)
266  for (unsigned j = 0; j < _tensor_dimensionality; ++j)
267  for (unsigned k = 0; k < _tensor_dimensionality; ++k)
268  for (unsigned l = 0; l < _tensor_dimensionality; ++l)
269  for (unsigned a = 0; a < _num_sp; ++a)
270  dstress_dtrial(i, j, k, l) +=
271  drot_dstress(i, a, k, l) * stress_params[a] * eT(a, j) +
272  _eigvecs(i, a) * stress_params[a] * drot_dstress(j, a, k, l);
273 
274  const std::vector<RankTwoTensor> dsp_trial = dstress_param_dstress(stress_trial);
275  for (unsigned i = 0; i < _tensor_dimensionality; ++i)
276  for (unsigned j = 0; j < _tensor_dimensionality; ++j)
277  for (unsigned k = 0; k < _tensor_dimensionality; ++k)
278  for (unsigned l = 0; l < _tensor_dimensionality; ++l)
279  for (unsigned a = 0; a < _num_sp; ++a)
280  for (unsigned b = 0; b < _num_sp; ++b)
281  dstress_dtrial(i, j, k, l) +=
282  _eigvecs(i, a) * dvar_dtrial[a][b] * dsp_trial[b](k, l) * eT(a, j);
283 
284  cto = dstress_dtrial * elasticity_tensor;
285 }
TensileStressUpdate::setStressAfterReturnV
virtual void setStressAfterReturnV(const RankTwoTensor &stress_trial, const std::vector< Real > &stress_params, Real gaE, const std::vector< Real > &intnl, const yieldAndFlow &smoothed_q, const RankFourTensor &Eijkl, RankTwoTensor &stress) const override
Sets stress from the admissible parameters.
Definition: TensileStressUpdate.C:80
MultiParameterPlasticityStressUpdate::yieldAndFlow
Struct designed to hold info about a single yield function and its derivatives, as well as the flow d...
Definition: MultiParameterPlasticityStressUpdate.h:214
TensileStressUpdate
TensileStressUpdate implements rate-independent associative tensile failure ("Rankine" plasticity) wi...
Definition: TensileStressUpdate.h:24
MultiParameterPlasticityStressUpdate::_Cij
std::vector< std::vector< Real > > _Cij
_Cij[i, j] * _Eij[j, k] = 1 iff j == k
Definition: MultiParameterPlasticityStressUpdate.h:144
TensileStressUpdate::computeStressParams
void computeStressParams(const RankTwoTensor &stress, std::vector< Real > &stress_params) const override
Computes stress_params, given stress.
Definition: TensileStressUpdate.C:44
TensileStressUpdate::_perfect_guess
const bool _perfect_guess
Whether to provide an estimate of the returned stress, based on perfect plasticity.
Definition: TensileStressUpdate.h:41
TensileStressUpdate::computeAllQV
void computeAllQV(const std::vector< Real > &stress_params, const std::vector< Real > &intnl, std::vector< yieldAndFlow > &all_q) const override
Completely fills all_q with correct values.
Definition: TensileStressUpdate.C:110
TensileStressUpdate::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: TensileStressUpdate.C:221
TensileStressUpdate::TensileStressUpdate
TensileStressUpdate(const InputParameters &parameters)
Definition: TensileStressUpdate.C:35
MultiParameterPlasticityStressUpdate::validParams
static InputParameters validParams()
Definition: MultiParameterPlasticityStressUpdate.C:23
MultiParameterPlasticityStressUpdate::_Eij
std::vector< std::vector< Real > > _Eij
E[i, j] in the system of equations to be solved.
Definition: MultiParameterPlasticityStressUpdate.h:138
TensileStressUpdate::d2stress_param_dstress
std::vector< RankFourTensor > d2stress_param_dstress(const RankTwoTensor &stress) const override
d2(stress_param[i])/d(stress)/d(stress) at given stress
Definition: TensileStressUpdate.C:61
TensileStressUpdate::dstress_param_dstress
std::vector< RankTwoTensor > dstress_param_dstress(const RankTwoTensor &stress) const override
d(stress_param[i])/d(stress) at given stress
Definition: TensileStressUpdate.C:52
TensileStressUpdate::validParams
static InputParameters validParams()
Definition: TensileStressUpdate.C:18
TensileStressUpdate::consistentTangentOperatorV
virtual void consistentTangentOperatorV(const RankTwoTensor &stress_trial, const std::vector< Real > &trial_stress_params, const RankTwoTensor &stress, const std::vector< Real > &stress_params, Real gaE, const yieldAndFlow &smoothed_q, const RankFourTensor &Eijkl, bool compute_full_tangent_operator, const std::vector< std::vector< Real >> &dvar_dtrial, RankFourTensor &cto) override
Calculates the consistent tangent operator.
Definition: TensileStressUpdate.C:227
defineLegacyParams
defineLegacyParams(TensileStressUpdate)
registerMooseObject
registerMooseObject("TensorMechanicsApp", TensileStressUpdate)
MultiParameterPlasticityStressUpdate::_num_sp
const unsigned _num_sp
Number of stress parameters.
Definition: MultiParameterPlasticityStressUpdate.h:132
TensileStressUpdate::preReturnMapV
virtual void preReturnMapV(const std::vector< Real > &trial_stress_params, const RankTwoTensor &stress_trial, const std::vector< Real > &intnl_old, const std::vector< Real > &yf, const RankFourTensor &Eijkl) override
Derived classes may employ this function to record stuff or do other computations prior to the return...
Definition: TensileStressUpdate.C:69
TensileStressUpdate::_strength
const TensorMechanicsHardeningModel & _strength
Hardening model for tensile strength.
Definition: TensileStressUpdate.h:38
MultiParameterPlasticityStressUpdate::_num_yf
const unsigned _num_yf
Number of yield functions.
Definition: MultiParameterPlasticityStressUpdate.h:147
TensileStressUpdate::tensile_strength
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
Definition: TensileStressUpdate.C:215
TensorMechanicsHardeningModel::derivative
virtual Real derivative(Real intnl) const
Definition: TensorMechanicsHardeningModel.C:47
TensileStressUpdate.h
TensileStressUpdate::initializeVarsV
void initializeVarsV(const std::vector< Real > &trial_stress_params, const std::vector< Real > &intnl_old, std::vector< Real > &stress_params, Real &gaE, std::vector< Real > &intnl) const override
Sets (stress_params, intnl) at "good guesses" of the solution to the Return-Map algorithm.
Definition: TensileStressUpdate.C:171
TensorMechanicsHardeningModel::value
virtual Real value(Real intnl) const
Definition: TensorMechanicsHardeningModel.C:45
RankTwoTensor
RankTwoTensorTempl< Real > RankTwoTensor
Definition: ACGrGrElasticDrivingForce.h:17
TensileStressUpdate::yieldFunctionValuesV
void yieldFunctionValuesV(const std::vector< Real > &stress_params, const std::vector< Real > &intnl, std::vector< Real > &yf) const override
Computes the values of the yield functions, given stress_params and intnl parameters.
Definition: TensileStressUpdate.C:95
MultiParameterPlasticityStressUpdate::_En
Real _En
normalising factor
Definition: MultiParameterPlasticityStressUpdate.h:141
RankFourTensorTempl< Real >
TensileStressUpdate::setEffectiveElasticity
void setEffectiveElasticity(const RankFourTensor &Eijkl) override
Sets _Eij and _En and _Cij.
Definition: TensileStressUpdate.C:153
RankTwoTensorTempl< Real >
TensileStressUpdate::_eigvecs
RankTwoTensor _eigvecs
Eigenvectors of the trial stress as a RankTwoTensor, in order to rotate the returned stress back to s...
Definition: TensileStressUpdate.h:44
TensileStressUpdate::setIntnlDerivativesV
void setIntnlDerivativesV(const std::vector< Real > &trial_stress_params, const std::vector< Real > &current_stress_params, const std::vector< Real > &intnl, std::vector< std::vector< Real >> &dintnl) const override
Sets the derivatives of internal parameters, based on the trial values of stress_params,...
Definition: TensileStressUpdate.C:204
TensileStressUpdate::setIntnlValuesV
void setIntnlValuesV(const std::vector< Real > &trial_stress_params, const std::vector< Real > &current_stress_params, const std::vector< Real > &intnl_old, std::vector< Real > &intnl) const override
Sets the internal parameters based on the trial values of stress_params, their current values,...
Definition: TensileStressUpdate.C:195
MultiParameterPlasticityStressUpdate
MultiParameterPlasticityStressUpdate performs the return-map algorithm and associated stress updates ...
Definition: MultiParameterPlasticityStressUpdate.h:99
MultiParameterPlasticityStressUpdate::_tensor_dimensionality
constexpr static unsigned _tensor_dimensionality
Internal dimensionality of tensors (currently this is 3 throughout tensor_mechanics)
Definition: MultiParameterPlasticityStressUpdate.h:129
TensorMechanicsHardeningModel
Hardening Model base class.
Definition: TensorMechanicsHardeningModel.h:27