www.mooseframework.org
MultiParameterPlasticityStressUpdate.h
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 #pragma once
11 
12 #include "StressUpdateBase.h"
13 
14 #include <array>
15 
95 {
96 public:
98 
100  unsigned num_sp,
101  unsigned num_yf,
102  unsigned num_intnl);
103 
104 protected:
105  virtual void initQpStatefulProperties() override;
107  virtual void updateState(RankTwoTensor & strain_increment,
108  RankTwoTensor & inelastic_strain_increment,
109  const RankTwoTensor & rotation_increment,
110  RankTwoTensor & stress_new,
111  const RankTwoTensor & stress_old,
113  const RankTwoTensor & elastic_strain_old,
114  bool compute_full_tangent_operator,
115  RankFourTensor & tangent_operator) override;
116 
117  virtual void propagateQpStatefulProperties() override;
118 
120  {
122  }
123 
125  constexpr static unsigned _tensor_dimensionality = 3;
126 
128  const unsigned _num_sp;
129 
131  const std::vector<Real> _definitely_ok_sp;
132 
134  std::vector<std::vector<Real>> _Eij;
135 
138 
140  std::vector<std::vector<Real>> _Cij;
141 
143  const unsigned _num_yf;
144 
146  const unsigned _num_intnl;
147 
149  const unsigned _max_nr_its;
150 
153 
156 
159 
161  const Real _f_tol;
162 
164  const Real _f_tol2;
165 
172 
174  bool _step_one;
175 
178 
181 
184 
187 
190 
193 
196 
199 
202 
205 
211  {
212  Real f; // yield function value
213  std::vector<Real> df; // df/d(stress_param[i])
214  std::vector<Real> df_di; // df/d(intnl_variable[a])
215  std::vector<Real> dg; // d(flow)/d(stress_param[i])
216  std::vector<std::vector<Real>> d2g; // d^2(flow)/d(sp[i])/d(sp[j])
217  std::vector<std::vector<Real>> d2g_di; // d^2(flow)/d(sp[i])/d(intnl[a])
218 
220 
221  yieldAndFlow(unsigned num_var, unsigned num_intnl)
222  : f(0.0),
223  df(num_var),
224  df_di(num_intnl),
225  dg(num_var),
226  d2g(num_var, std::vector<Real>(num_var, 0.0)),
227  d2g_di(num_var, std::vector<Real>(num_intnl, 0.0))
228  {
229  }
230 
231  // this may be involved in the smoothing of a group of yield functions
232  bool operator<(const yieldAndFlow & fd) const { return f < fd.f; }
233  };
234 
242  Real yieldF(const std::vector<Real> & stress_params, const std::vector<Real> & intnl) const;
243 
249  Real yieldF(const std::vector<Real> & yfs) const;
250 
262  Real ismoother(Real f_diff) const;
263 
267  Real smoother(Real f_diff) const;
268 
272  Real dsmoother(Real f_diff) const;
273 
281  yieldAndFlow smoothAllQuantities(const std::vector<Real> & stress_params,
282  const std::vector<Real> & intnl) const;
283 
323  int lineSearch(Real & res2,
324  std::vector<Real> & stress_params,
325  Real & gaE,
326  const std::vector<Real> & trial_stress_params,
327  yieldAndFlow & smoothed_q,
328  const std::vector<Real> & intnl_ok,
329  std::vector<Real> & intnl,
330  std::vector<Real> & rhs,
331  Real & linesearch_needed) const;
332 
346  int nrStep(const yieldAndFlow & smoothed_q,
347  const std::vector<Real> & trial_stress_params,
348  const std::vector<Real> & stress_params,
349  const std::vector<Real> & intnl,
350  Real gaE,
351  std::vector<Real> & rhs) const;
352 
358  Real calculateRes2(const std::vector<Real> & rhs) const;
359 
376  void calculateRHS(const std::vector<Real> & trial_stress_params,
377  const std::vector<Real> & stress_params,
378  Real gaE,
379  const yieldAndFlow & smoothed_q,
380  std::vector<Real> & rhs) const;
381 
395  void dnRHSdVar(const yieldAndFlow & smoothed_q,
396  const std::vector<std::vector<Real>> & dintnl,
397  const std::vector<Real> & stress_params,
398  Real gaE,
399  std::vector<double> & jac) const;
400 
405  virtual void errorHandler(const std::string & message) const;
406 
415  virtual void yieldFunctionValuesV(const std::vector<Real> & stress_params,
416  const std::vector<Real> & intnl,
417  std::vector<Real> & yf) const = 0;
418 
431  virtual void computeAllQV(const std::vector<Real> & stress_params,
432  const std::vector<Real> & intnl,
433  std::vector<yieldAndFlow> & all_q) const = 0;
434 
446  virtual void preReturnMapV(const std::vector<Real> & trial_stress_params,
447  const RankTwoTensor & stress_trial,
448  const std::vector<Real> & intnl_old,
449  const std::vector<Real> & yf,
450  const RankFourTensor & Eijkl);
451 
467  virtual void initializeVarsV(const std::vector<Real> & trial_stress_params,
468  const std::vector<Real> & intnl_old,
469  std::vector<Real> & stress_params,
470  Real & gaE,
471  std::vector<Real> & intnl) const;
472 
483  virtual void setIntnlValuesV(const std::vector<Real> & trial_stress_params,
484  const std::vector<Real> & current_stress_params,
485  const std::vector<Real> & intnl_old,
486  std::vector<Real> & intnl) const = 0;
487 
498  virtual void setIntnlDerivativesV(const std::vector<Real> & trial_stress_params,
499  const std::vector<Real> & current_stress_params,
500  const std::vector<Real> & intnl,
501  std::vector<std::vector<Real>> & dintnl) const = 0;
502 
509  virtual void computeStressParams(const RankTwoTensor & stress,
510  std::vector<Real> & stress_params) const = 0;
511 
519  virtual void initializeReturnProcess();
520 
527  virtual void finalizeReturnProcess(const RankTwoTensor & rotation_increment);
528 
546  virtual void setStressAfterReturnV(const RankTwoTensor & stress_trial,
547  const std::vector<Real> & stress_params,
548  Real gaE,
549  const std::vector<Real> & intnl,
550  const yieldAndFlow & smoothed_q,
551  const RankFourTensor & Eijkl,
552  RankTwoTensor & stress) const = 0;
553 
569  virtual void
571  Real gaE,
572  const yieldAndFlow & smoothed_q,
573  const RankFourTensor & elasticity_tensor,
574  const RankTwoTensor & returned_stress,
575  RankTwoTensor & inelastic_strain_increment) const;
576 
598  virtual void consistentTangentOperatorV(const RankTwoTensor & stress_trial,
599  const std::vector<Real> & trial_stress_params,
600  const RankTwoTensor & stress,
601  const std::vector<Real> & stress_params,
602  Real gaE,
603  const yieldAndFlow & smoothed_q,
604  const RankFourTensor & Eijkl,
605  bool compute_full_tangent_operator,
606  const std::vector<std::vector<Real>> & dvar_dtrial,
607  RankFourTensor & cto);
608 
614  virtual std::vector<RankTwoTensor> dstress_param_dstress(const RankTwoTensor & stress) const = 0;
615 
621  virtual std::vector<RankFourTensor>
622  d2stress_param_dstress(const RankTwoTensor & stress) const = 0;
623 
625  virtual void setEffectiveElasticity(const RankFourTensor & Eijkl) = 0;
626 
648  void dVardTrial(bool elastic_only,
649  const std::vector<Real> & trial_stress_params,
650  const std::vector<Real> & stress_params,
651  Real gaE,
652  const std::vector<Real> & intnl,
653  const yieldAndFlow & smoothed_q,
654  Real step_size,
655  bool compute_full_tangent_operator,
656  std::vector<std::vector<Real>> & dvar_dtrial) const;
657 
665  bool precisionLoss(const std::vector<Real> & solution,
666  const std::vector<Real> & stress_params,
667  Real gaE) const;
668 
669 private:
677  std::vector<Real> _trial_sp;
678 
685 
694  std::vector<Real> _rhs;
695 
699  std::vector<std::vector<Real>> _dvar_dtrial;
700 
709  std::vector<Real> _ok_sp;
710 
714  std::vector<Real> _ok_intnl;
715 
723  std::vector<Real> _del_stress_params;
724 
728  std::vector<Real> _current_sp;
729 
733  std::vector<Real> _current_intnl;
734 
735 private:
740  {
741  cos,
742  poly1,
743  poly2,
744  poly3
746 };
const Real _smoothing_tol2
Square of the smoothing tolerance.
MaterialProperty< std::vector< Real > > & _yf
yield functions
std::vector< Real > _trial_sp
"Trial" value of stress_params that initializes the return-map process This is derived from stress = ...
void dnRHSdVar(const yieldAndFlow &smoothed_q, const std::vector< std::vector< Real >> &dintnl, const std::vector< Real > &stress_params, Real gaE, std::vector< double > &jac) const
Derivative of -RHS with respect to the stress_params and gaE, placed into an array ready for solving ...
MultiParameterPlasticityStressUpdate(const InputParameters &parameters, unsigned num_sp, unsigned num_yf, unsigned num_intnl)
virtual void setEffectiveElasticity(const RankFourTensor &Eijkl)=0
Sets _Eij and _En and _Cij.
const bool _warn_about_precision_loss
Output a warning message if precision loss is encountered during the return-map process.
virtual void errorHandler(const std::string &message) const
Performs any necessary cleaning-up, then throw MooseException(message)
virtual 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
Sets (stress_params, intnl) at "good guesses" of the solution to the Return-Map algorithm.
int lineSearch(Real &res2, std::vector< Real > &stress_params, Real &gaE, const std::vector< Real > &trial_stress_params, yieldAndFlow &smoothed_q, const std::vector< Real > &intnl_ok, std::vector< Real > &intnl, std::vector< Real > &rhs, Real &linesearch_needed) const
Performs a line-search to find stress_params and gaE Upon entry:
const unsigned _num_intnl
Number of internal parameters.
Struct designed to hold info about a single yield function and its derivatives, as well as the flow d...
std::vector< std::vector< Real > > _Cij
_Cij[i, j] * _Eij[j, k] = 1 iff j == k
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)
Calculates the consistent tangent operator.
virtual 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 =0
Sets the internal parameters based on the trial values of stress_params, their current values...
std::vector< Real > _ok_intnl
The state (ok_sp, ok_intnl) is known to be admissible.
const Real _f_tol2
Square of the yield-function tolerance.
const Real _min_step_size
In order to help the Newton-Raphson procedure, the applied strain increment may be applied in sub-inc...
Real smoother(Real f_diff) const
Derivative of ismoother.
int nrStep(const yieldAndFlow &smoothed_q, const std::vector< Real > &trial_stress_params, const std::vector< Real > &stress_params, const std::vector< Real > &intnl, Real gaE, std::vector< Real > &rhs) const
Performs a Newton-Raphson step to attempt to zero rhs Upon return, rhs will contain the solution...
void calculateRHS(const std::vector< Real > &trial_stress_params, const std::vector< Real > &stress_params, Real gaE, const yieldAndFlow &smoothed_q, std::vector< Real > &rhs) const
Calculates the RHS in the following 0 = rhs[0] = S[0] - S[0]^trial + ga * E[0, j] * dg/dS[j] 0 = rhs[...
virtual void updateState(RankTwoTensor &strain_increment, RankTwoTensor &inelastic_strain_increment, const RankTwoTensor &rotation_increment, RankTwoTensor &stress_new, const RankTwoTensor &stress_old, const RankFourTensor &elasticity_tensor, const RankTwoTensor &elastic_strain_old, bool compute_full_tangent_operator, RankFourTensor &tangent_operator) override
MaterialProperty< std::vector< Real > > & _intnl
internal parameters
virtual void updateState(GR2 &strain_increment, GR2 &inelastic_strain_increment, const GR2 &rotation_increment, GR2 &stress_new, const RankTwoTensor &stress_old, const GR4 &elasticity_tensor, const RankTwoTensor &elastic_strain_old, bool compute_full_tangent_operator=false, RankFourTensor &tangent_operator=StressUpdateBaseTempl< is_ad >::_identityTensor)
Given a strain increment that results in a trial stress, perform some procedure (such as an iterative...
yieldAndFlow smoothAllQuantities(const std::vector< Real > &stress_params, const std::vector< Real > &intnl) const
Calculates all yield functions and derivatives, and then performs the smoothing scheme.
Real elasticity_tensor(unsigned int i, unsigned int j, unsigned int k, unsigned int l)
std::vector< std::vector< Real > > _Eij
E[i, j] in the system of equations to be solved.
const Real _f_tol
The yield-function tolerance.
const unsigned _num_yf
Number of yield functions.
virtual void computeAllQV(const std::vector< Real > &stress_params, const std::vector< Real > &intnl, std::vector< yieldAndFlow > &all_q) const =0
Completely fills all_q with correct values.
virtual void setInelasticStrainIncrementAfterReturn(const RankTwoTensor &stress_trial, Real gaE, const yieldAndFlow &smoothed_q, const RankFourTensor &elasticity_tensor, const RankTwoTensor &returned_stress, RankTwoTensor &inelastic_strain_increment) const
Sets inelastic strain increment from the returned configuration This is called after the return-map p...
const std::vector< Real > _definitely_ok_sp
An admissible value of stress_params at the initial time.
const MaterialProperty< std::vector< Real > > & _intnl_old
old values of internal parameters
RankTwoTensor _stress_trial
"Trial" value of stress that is set at the beginning of the return-map process.
virtual void propagateQpStatefulProperties() override
If updateState is not called during a timestep, this will be.
virtual void computeStressParams(const RankTwoTensor &stress, std::vector< Real > &stress_params) const =0
Computes stress_params, given stress.
const MaterialProperty< Real > & _max_iter_used_old
Old value of maximum number of Newton-Raphson iterations used in the return-map during the course of ...
virtual TangentCalculationMethod getTangentCalculationMethod() override
Real calculateRes2(const std::vector< Real > &rhs) const
Calculates the residual-squared for the Newton-Raphson + line-search.
virtual 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 =0
Sets the derivatives of internal parameters, based on the trial values of stress_params, their current values, and the current values of the internal parameters.
virtual void yieldFunctionValuesV(const std::vector< Real > &stress_params, const std::vector< Real > &intnl, std::vector< Real > &yf) const =0
Computes the values of the yield functions, given stress_params and intnl parameters.
StressUpdateBase is a material that is not called by MOOSE because of the compute=false flag set in t...
Real dsmoother(Real f_diff) const
Derivative of smoother.
enum MultiParameterPlasticityStressUpdate::SmootherFunctionType _smoother_function_type
MaterialProperty< Real > & _iter
Number of Newton-Raphson iterations used in the return-map.
MaterialProperty< RankTwoTensor > & _plastic_strain
plastic strain
virtual std::vector< RankFourTensor > d2stress_param_dstress(const RankTwoTensor &stress) const =0
d2(stress_param[i])/d(stress)/d(stress) at given stress
const Real _smoothing_tol
Smoothing tolerance: edges of the yield surface get smoothed by this amount.
std::vector< Real > _current_intnl
The current values of the internal params during the Newton-Raphson.
virtual void finalizeReturnProcess(const RankTwoTensor &rotation_increment)
Derived classes may use this to perform calculations after the return-map process has completed succe...
const bool _perform_finite_strain_rotations
Whether to perform finite-strain rotations.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string message
MaterialProperty< Real > & _linesearch_needed
Whether a line-search was needed in the latest Newton-Raphson process (1 if true, 0 otherwise) ...
MultiParameterPlasticityStressUpdate performs the return-map algorithm and associated stress updates ...
std::vector< std::vector< Real > > _dvar_dtrial
d({stress_param[i], gaE})/d(trial_stress_param[j])
std::vector< Real > _current_sp
The current values of the stress params during the Newton-Raphson.
bool _step_one
handles case of initial_stress that is inadmissible being supplied
void dVardTrial(bool elastic_only, const std::vector< Real > &trial_stress_params, const std::vector< Real > &stress_params, Real gaE, const std::vector< Real > &intnl, const yieldAndFlow &smoothed_q, Real step_size, bool compute_full_tangent_operator, std::vector< std::vector< Real >> &dvar_dtrial) const
Calculates derivatives of the stress_params and gaE with repect to the trial values of the stress_par...
const InputParameters & parameters() const
virtual std::vector< RankTwoTensor > dstress_param_dstress(const RankTwoTensor &stress) const =0
d(stress_param[i])/d(stress) at given stress
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 =0
Sets stress from the admissible parameters.
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)
Derived classes may employ this function to record stuff or do other computations prior to the return...
MaterialProperty< Real > & _max_iter_used
Maximum number of Newton-Raphson iterations used in the return-map during the course of the entire si...
bool precisionLoss(const std::vector< Real > &solution, const std::vector< Real > &stress_params, Real gaE) const
Check whether precision loss has occurred.
std::vector< Real > _del_stress_params
_del_stress_params = trial_stress_params - ok_sp This is fixed at the beginning of the return-map pro...
TangentCalculationMethod
TangentCalculationMethod is an enum that determines the calculation method for the tangent operator...
const unsigned _max_nr_its
Maximum number of Newton-Raphson iterations allowed in the return-map process.
virtual void initializeReturnProcess()
Derived classes may use this to perform calculations before any return-map process is performed...
const MaterialProperty< RankTwoTensor > & _plastic_strain_old
Old value of plastic strain.
static constexpr unsigned _tensor_dimensionality
Internal dimensionality of tensors (currently this is 3 throughout tensor_mechanics) ...
Real yieldF(const std::vector< Real > &stress_params, const std::vector< Real > &intnl) const
Computes the smoothed yield function.
std::vector< Real > _rhs
0 = rhs[0] = S[0] - S[0]^trial + ga * E[0, i] * dg/dS[i] 0 = rhs[1] = S[1] - S[1]^trial + ga * E[1...
std::vector< Real > _ok_sp
The state (ok_sp, ok_intnl) is known to be admissible, so ok_sp are stress_params that are "OK"...
const unsigned _num_sp
Number of stress parameters.
Real ismoother(Real f_diff) const
Smooths yield functions.