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 
17 
18 template <>
20 
100 {
101 public:
102  static InputParameters validParams();
103 
104  MultiParameterPlasticityStressUpdate(const InputParameters & parameters,
105  unsigned num_sp,
106  unsigned num_yf,
107  unsigned num_intnl);
108 
109 protected:
110  virtual void initQpStatefulProperties() override;
111  virtual void updateState(RankTwoTensor & strain_increment,
112  RankTwoTensor & inelastic_strain_increment,
113  const RankTwoTensor & rotation_increment,
114  RankTwoTensor & stress_new,
115  const RankTwoTensor & stress_old,
116  const RankFourTensor & elasticity_tensor,
117  const RankTwoTensor & elastic_strain_old,
118  bool compute_full_tangent_operator,
119  RankFourTensor & tangent_operator) override;
120 
121  virtual void propagateQpStatefulProperties() override;
122 
124  {
126  }
127 
129  constexpr static unsigned _tensor_dimensionality = 3;
130 
132  const unsigned _num_sp;
133 
135  const std::vector<Real> _definitely_ok_sp;
136 
138  std::vector<std::vector<Real>> _Eij;
139 
141  Real _En;
142 
144  std::vector<std::vector<Real>> _Cij;
145 
147  const unsigned _num_yf;
148 
150  const unsigned _num_intnl;
151 
153  const unsigned _max_nr_its;
154 
157 
159  const Real _smoothing_tol;
160 
162  const Real _smoothing_tol2;
163 
165  const Real _f_tol;
166 
168  const Real _f_tol2;
169 
175  const Real _min_step_size;
176 
178  bool _step_one;
179 
182 
184  MaterialProperty<RankTwoTensor> & _plastic_strain;
185 
187  const MaterialProperty<RankTwoTensor> & _plastic_strain_old;
188 
190  MaterialProperty<std::vector<Real>> & _intnl;
191 
193  const MaterialProperty<std::vector<Real>> & _intnl_old;
194 
196  MaterialProperty<std::vector<Real>> & _yf;
197 
199  MaterialProperty<Real> & _iter;
200 
202  MaterialProperty<Real> & _max_iter_used;
203 
205  const MaterialProperty<Real> & _max_iter_used_old;
206 
208  MaterialProperty<Real> & _linesearch_needed;
209 
215  {
216  Real f; // yield function value
217  std::vector<Real> df; // df/d(stress_param[i])
218  std::vector<Real> df_di; // df/d(intnl_variable[a])
219  std::vector<Real> dg; // d(flow)/d(stress_param[i])
220  std::vector<std::vector<Real>> d2g; // d^2(flow)/d(sp[i])/d(sp[j])
221  std::vector<std::vector<Real>> d2g_di; // d^2(flow)/d(sp[i])/d(intnl[a])
222 
224 
225  yieldAndFlow(unsigned num_var, unsigned num_intnl)
226  : f(0.0),
227  df(num_var),
228  df_di(num_intnl),
229  dg(num_var),
230  d2g(num_var, std::vector<Real>(num_var, 0.0)),
231  d2g_di(num_var, std::vector<Real>(num_intnl, 0.0))
232  {
233  }
234 
235  // this may be involved in the smoothing of a group of yield functions
236  bool operator<(const yieldAndFlow & fd) const { return f < fd.f; }
237  };
238 
246  Real yieldF(const std::vector<Real> & stress_params, const std::vector<Real> & intnl) const;
247 
253  Real yieldF(const std::vector<Real> & yfs) const;
254 
266  Real ismoother(Real f_diff) const;
267 
271  Real smoother(Real f_diff) const;
272 
276  Real dsmoother(Real f_diff) const;
277 
285  yieldAndFlow smoothAllQuantities(const std::vector<Real> & stress_params,
286  const std::vector<Real> & intnl) const;
287 
327  int lineSearch(Real & res2,
328  std::vector<Real> & stress_params,
329  Real & gaE,
330  const std::vector<Real> & trial_stress_params,
331  yieldAndFlow & smoothed_q,
332  const std::vector<Real> & intnl_ok,
333  std::vector<Real> & intnl,
334  std::vector<Real> & rhs,
335  Real & linesearch_needed) const;
336 
350  int nrStep(const yieldAndFlow & smoothed_q,
351  const std::vector<Real> & trial_stress_params,
352  const std::vector<Real> & stress_params,
353  const std::vector<Real> & intnl,
354  Real gaE,
355  std::vector<Real> & rhs) const;
356 
362  Real calculateRes2(const std::vector<Real> & rhs) const;
363 
380  void calculateRHS(const std::vector<Real> & trial_stress_params,
381  const std::vector<Real> & stress_params,
382  Real gaE,
383  const yieldAndFlow & smoothed_q,
384  std::vector<Real> & rhs) const;
385 
399  void dnRHSdVar(const yieldAndFlow & smoothed_q,
400  const std::vector<std::vector<Real>> & dintnl,
401  const std::vector<Real> & stress_params,
402  Real gaE,
403  std::vector<double> & jac) const;
404 
409  virtual void errorHandler(const std::string & message) const;
410 
419  virtual void yieldFunctionValuesV(const std::vector<Real> & stress_params,
420  const std::vector<Real> & intnl,
421  std::vector<Real> & yf) const = 0;
422 
435  virtual void computeAllQV(const std::vector<Real> & stress_params,
436  const std::vector<Real> & intnl,
437  std::vector<yieldAndFlow> & all_q) const = 0;
438 
450  virtual void preReturnMapV(const std::vector<Real> & trial_stress_params,
451  const RankTwoTensor & stress_trial,
452  const std::vector<Real> & intnl_old,
453  const std::vector<Real> & yf,
454  const RankFourTensor & Eijkl);
455 
471  virtual void initializeVarsV(const std::vector<Real> & trial_stress_params,
472  const std::vector<Real> & intnl_old,
473  std::vector<Real> & stress_params,
474  Real & gaE,
475  std::vector<Real> & intnl) const;
476 
487  virtual void setIntnlValuesV(const std::vector<Real> & trial_stress_params,
488  const std::vector<Real> & current_stress_params,
489  const std::vector<Real> & intnl_old,
490  std::vector<Real> & intnl) const = 0;
491 
502  virtual void setIntnlDerivativesV(const std::vector<Real> & trial_stress_params,
503  const std::vector<Real> & current_stress_params,
504  const std::vector<Real> & intnl,
505  std::vector<std::vector<Real>> & dintnl) const = 0;
506 
513  virtual void computeStressParams(const RankTwoTensor & stress,
514  std::vector<Real> & stress_params) const = 0;
515 
523  virtual void initializeReturnProcess();
524 
531  virtual void finalizeReturnProcess(const RankTwoTensor & rotation_increment);
532 
550  virtual void setStressAfterReturnV(const RankTwoTensor & stress_trial,
551  const std::vector<Real> & stress_params,
552  Real gaE,
553  const std::vector<Real> & intnl,
554  const yieldAndFlow & smoothed_q,
555  const RankFourTensor & Eijkl,
556  RankTwoTensor & stress) const = 0;
557 
573  virtual void
575  Real gaE,
576  const yieldAndFlow & smoothed_q,
577  const RankFourTensor & elasticity_tensor,
578  const RankTwoTensor & returned_stress,
579  RankTwoTensor & inelastic_strain_increment) const;
580 
602  virtual void consistentTangentOperatorV(const RankTwoTensor & stress_trial,
603  const std::vector<Real> & trial_stress_params,
604  const RankTwoTensor & stress,
605  const std::vector<Real> & stress_params,
606  Real gaE,
607  const yieldAndFlow & smoothed_q,
608  const RankFourTensor & Eijkl,
609  bool compute_full_tangent_operator,
610  const std::vector<std::vector<Real>> & dvar_dtrial,
611  RankFourTensor & cto);
612 
618  virtual std::vector<RankTwoTensor> dstress_param_dstress(const RankTwoTensor & stress) const = 0;
619 
625  virtual std::vector<RankFourTensor>
626  d2stress_param_dstress(const RankTwoTensor & stress) const = 0;
627 
629  virtual void setEffectiveElasticity(const RankFourTensor & Eijkl) = 0;
630 
652  void dVardTrial(bool elastic_only,
653  const std::vector<Real> & trial_stress_params,
654  const std::vector<Real> & stress_params,
655  Real gaE,
656  const std::vector<Real> & intnl,
657  const yieldAndFlow & smoothed_q,
658  Real step_size,
659  bool compute_full_tangent_operator,
660  std::vector<std::vector<Real>> & dvar_dtrial) const;
661 
669  bool precisionLoss(const std::vector<Real> & solution,
670  const std::vector<Real> & stress_params,
671  Real gaE) const;
672 
673 private:
681  std::vector<Real> _trial_sp;
682 
689 
698  std::vector<Real> _rhs;
699 
703  std::vector<std::vector<Real>> _dvar_dtrial;
704 
713  std::vector<Real> _ok_sp;
714 
718  std::vector<Real> _ok_intnl;
719 
727  std::vector<Real> _del_stress_params;
728 
732  std::vector<Real> _current_sp;
733 
737  std::vector<Real> _current_intnl;
738 
739 private:
744  {
745  cos,
746  poly1,
747  poly2,
748  poly3
750 };
MultiParameterPlasticityStressUpdate::computeAllQV
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.
MultiParameterPlasticityStressUpdate::_max_nr_its
const unsigned _max_nr_its
Maximum number of Newton-Raphson iterations allowed in the return-map process.
Definition: MultiParameterPlasticityStressUpdate.h:153
MultiParameterPlasticityStressUpdate::_warn_about_precision_loss
const bool _warn_about_precision_loss
Output a warning message if precision loss is encountered during the return-map process.
Definition: MultiParameterPlasticityStressUpdate.h:181
MultiParameterPlasticityStressUpdate::_stress_trial
RankTwoTensor _stress_trial
"Trial" value of stress that is set at the beginning of the return-map process.
Definition: MultiParameterPlasticityStressUpdate.h:688
MultiParameterPlasticityStressUpdate::ismoother
Real ismoother(Real f_diff) const
Smooths yield functions.
Definition: MultiParameterPlasticityStressUpdate.C:484
MultiParameterPlasticityStressUpdate::MultiParameterPlasticityStressUpdate
MultiParameterPlasticityStressUpdate(const InputParameters &parameters, unsigned num_sp, unsigned num_yf, unsigned num_intnl)
Definition: MultiParameterPlasticityStressUpdate.C:82
MultiParameterPlasticityStressUpdate::yieldAndFlow::d2g
std::vector< std::vector< Real > > d2g
Definition: MultiParameterPlasticityStressUpdate.h:220
MultiParameterPlasticityStressUpdate::_current_intnl
std::vector< Real > _current_intnl
The current values of the internal params during the Newton-Raphson.
Definition: MultiParameterPlasticityStressUpdate.h:737
MultiParameterPlasticityStressUpdate::SmootherFunctionType::poly1
MultiParameterPlasticityStressUpdate::_perform_finite_strain_rotations
const bool _perform_finite_strain_rotations
Whether to perform finite-strain rotations.
Definition: MultiParameterPlasticityStressUpdate.h:156
MultiParameterPlasticityStressUpdate::SmootherFunctionType::cos
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
MultiParameterPlasticityStressUpdate::_smoothing_tol
const Real _smoothing_tol
Smoothing tolerance: edges of the yield surface get smoothed by this amount.
Definition: MultiParameterPlasticityStressUpdate.h:159
MultiParameterPlasticityStressUpdate::yieldAndFlow::dg
std::vector< Real > dg
Definition: MultiParameterPlasticityStressUpdate.h:219
TangentCalculationMethod::FULL
MultiParameterPlasticityStressUpdate::finalizeReturnProcess
virtual void finalizeReturnProcess(const RankTwoTensor &rotation_increment)
Derived classes may use this to perform calculations after the return-map process has completed succe...
Definition: MultiParameterPlasticityStressUpdate.C:672
MultiParameterPlasticityStressUpdate::_smoother_function_type
enum MultiParameterPlasticityStressUpdate::SmootherFunctionType _smoother_function_type
MultiParameterPlasticityStressUpdate::_Cij
std::vector< std::vector< Real > > _Cij
_Cij[i, j] * _Eij[j, k] = 1 iff j == k
Definition: MultiParameterPlasticityStressUpdate.h:144
MultiParameterPlasticityStressUpdate::lineSearch
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:
Definition: MultiParameterPlasticityStressUpdate.C:553
MultiParameterPlasticityStressUpdate::SmootherFunctionType::poly3
MultiParameterPlasticityStressUpdate::_iter
MaterialProperty< Real > & _iter
Number of Newton-Raphson iterations used in the return-map.
Definition: MultiParameterPlasticityStressUpdate.h:199
MultiParameterPlasticityStressUpdate::_f_tol2
const Real _f_tol2
Square of the yield-function tolerance.
Definition: MultiParameterPlasticityStressUpdate.h:168
MultiParameterPlasticityStressUpdate::_min_step_size
const Real _min_step_size
In order to help the Newton-Raphson procedure, the applied strain increment may be applied in sub-inc...
Definition: MultiParameterPlasticityStressUpdate.h:175
MultiParameterPlasticityStressUpdate::setEffectiveElasticity
virtual void setEffectiveElasticity(const RankFourTensor &Eijkl)=0
Sets _Eij and _En and _Cij.
MultiParameterPlasticityStressUpdate::d2stress_param_dstress
virtual std::vector< RankFourTensor > d2stress_param_dstress(const RankTwoTensor &stress) const =0
d2(stress_param[i])/d(stress)/d(stress) at given stress
MultiParameterPlasticityStressUpdate::initQpStatefulProperties
virtual void initQpStatefulProperties() override
Definition: MultiParameterPlasticityStressUpdate.C:141
MultiParameterPlasticityStressUpdate::_ok_intnl
std::vector< Real > _ok_intnl
The state (ok_sp, ok_intnl) is known to be admissible.
Definition: MultiParameterPlasticityStressUpdate.h:718
MultiParameterPlasticityStressUpdate::updateState
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
Given a strain increment that results in a trial stress, perform some procedure (such as an iterative...
Definition: MultiParameterPlasticityStressUpdate.C:160
MultiParameterPlasticityStressUpdate::_linesearch_needed
MaterialProperty< Real > & _linesearch_needed
Whether a line-search was needed in the latest Newton-Raphson process (1 if true, 0 otherwise)
Definition: MultiParameterPlasticityStressUpdate.h:208
MultiParameterPlasticityStressUpdate::dstress_param_dstress
virtual std::vector< RankTwoTensor > dstress_param_dstress(const RankTwoTensor &stress) const =0
d(stress_param[i])/d(stress) at given stress
MultiParameterPlasticityStressUpdate::dVardTrial
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...
Definition: MultiParameterPlasticityStressUpdate.C:873
MultiParameterPlasticityStressUpdate::initializeVarsV
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.
Definition: MultiParameterPlasticityStressUpdate.C:713
MultiParameterPlasticityStressUpdate::precisionLoss
bool precisionLoss(const std::vector< Real > &solution, const std::vector< Real > &stress_params, Real gaE) const
Check whether precision loss has occurred.
Definition: MultiParameterPlasticityStressUpdate.C:983
MultiParameterPlasticityStressUpdate::validParams
static InputParameters validParams()
Definition: MultiParameterPlasticityStressUpdate.C:23
MultiParameterPlasticityStressUpdate::_plastic_strain_old
const MaterialProperty< RankTwoTensor > & _plastic_strain_old
Old value of plastic strain.
Definition: MultiParameterPlasticityStressUpdate.h:187
MultiParameterPlasticityStressUpdate::_Eij
std::vector< std::vector< Real > > _Eij
E[i, j] in the system of equations to be solved.
Definition: MultiParameterPlasticityStressUpdate.h:138
MultiParameterPlasticityStressUpdate::yieldF
Real yieldF(const std::vector< Real > &stress_params, const std::vector< Real > &intnl) const
Computes the smoothed yield function.
Definition: MultiParameterPlasticityStressUpdate.C:688
MultiParameterPlasticityStressUpdate::_del_stress_params
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...
Definition: MultiParameterPlasticityStressUpdate.h:727
MultiParameterPlasticityStressUpdate::_intnl_old
const MaterialProperty< std::vector< Real > > & _intnl_old
old values of internal parameters
Definition: MultiParameterPlasticityStressUpdate.h:193
MultiParameterPlasticityStressUpdate::computeStressParams
virtual void computeStressParams(const RankTwoTensor &stress, std::vector< Real > &stress_params) const =0
Computes stress_params, given stress.
MultiParameterPlasticityStressUpdate::initializeReturnProcess
virtual void initializeReturnProcess()
Derived classes may use this to perform calculations before any return-map process is performed,...
Definition: MultiParameterPlasticityStressUpdate.C:667
MultiParameterPlasticityStressUpdate::smoother
Real smoother(Real f_diff) const
Derivative of ismoother.
Definition: MultiParameterPlasticityStressUpdate.C:510
MultiParameterPlasticityStressUpdate::_rhs
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,...
Definition: MultiParameterPlasticityStressUpdate.h:698
MultiParameterPlasticityStressUpdate::_num_sp
const unsigned _num_sp
Number of stress parameters.
Definition: MultiParameterPlasticityStressUpdate.h:132
MultiParameterPlasticityStressUpdate::yieldFunctionValuesV
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.
MultiParameterPlasticityStressUpdate::_yf
MaterialProperty< std::vector< Real > > & _yf
yield functions
Definition: MultiParameterPlasticityStressUpdate.h:196
MultiParameterPlasticityStressUpdate::_num_yf
const unsigned _num_yf
Number of yield functions.
Definition: MultiParameterPlasticityStressUpdate.h:147
MultiParameterPlasticityStressUpdate::dnRHSdVar
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 ...
Definition: MultiParameterPlasticityStressUpdate.C:828
MultiParameterPlasticityStressUpdate::_plastic_strain
MaterialProperty< RankTwoTensor > & _plastic_strain
plastic strain
Definition: MultiParameterPlasticityStressUpdate.h:184
MultiParameterPlasticityStressUpdate::_definitely_ok_sp
const std::vector< Real > _definitely_ok_sp
An admissible value of stress_params at the initial time.
Definition: MultiParameterPlasticityStressUpdate.h:135
MultiParameterPlasticityStressUpdate::propagateQpStatefulProperties
virtual void propagateQpStatefulProperties() override
If updateState is not called during a timestep, this will be.
Definition: MultiParameterPlasticityStressUpdate.C:152
validParams< MultiParameterPlasticityStressUpdate >
InputParameters validParams< MultiParameterPlasticityStressUpdate >()
MultiParameterPlasticityStressUpdate::SmootherFunctionType
SmootherFunctionType
The type of smoother function.
Definition: MultiParameterPlasticityStressUpdate.h:743
MultiParameterPlasticityStressUpdate::_trial_sp
std::vector< Real > _trial_sp
"Trial" value of stress_params that initializes the return-map process This is derived from stress = ...
Definition: MultiParameterPlasticityStressUpdate.h:681
MultiParameterPlasticityStressUpdate::yieldAndFlow::yieldAndFlow
yieldAndFlow(unsigned num_var, unsigned num_intnl)
Definition: MultiParameterPlasticityStressUpdate.h:225
MultiParameterPlasticityStressUpdate::_num_intnl
const unsigned _num_intnl
Number of internal parameters.
Definition: MultiParameterPlasticityStressUpdate.h:150
MultiParameterPlasticityStressUpdate::_step_one
bool _step_one
handles case of initial_stress that is inadmissible being supplied
Definition: MultiParameterPlasticityStressUpdate.h:178
MultiParameterPlasticityStressUpdate::dsmoother
Real dsmoother(Real f_diff) const
Derivative of smoother.
Definition: MultiParameterPlasticityStressUpdate.C:533
MultiParameterPlasticityStressUpdate::calculateRes2
Real calculateRes2(const std::vector< Real > &rhs) const
Calculates the residual-squared for the Newton-Raphson + line-search.
Definition: MultiParameterPlasticityStressUpdate.C:802
MultiParameterPlasticityStressUpdate::setIntnlDerivativesV
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,...
MultiParameterPlasticityStressUpdate::yieldAndFlow::yieldAndFlow
yieldAndFlow()
Definition: MultiParameterPlasticityStressUpdate.h:223
MultiParameterPlasticityStressUpdate::errorHandler
virtual void errorHandler(const std::string &message) const
Performs any necessary cleaning-up, then throw MooseException(message)
Definition: MultiParameterPlasticityStressUpdate.C:661
MultiParameterPlasticityStressUpdate::_current_sp
std::vector< Real > _current_sp
The current values of the stress params during the Newton-Raphson.
Definition: MultiParameterPlasticityStressUpdate.h:732
MultiParameterPlasticityStressUpdate::_En
Real _En
normalising factor
Definition: MultiParameterPlasticityStressUpdate.h:141
MultiParameterPlasticityStressUpdate::yieldAndFlow::f
Real f
Definition: MultiParameterPlasticityStressUpdate.h:216
RankFourTensorTempl< Real >
MultiParameterPlasticityStressUpdate::yieldAndFlow::df
std::vector< Real > df
Definition: MultiParameterPlasticityStressUpdate.h:217
MultiParameterPlasticityStressUpdate::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)
Calculates the consistent tangent operator.
Definition: MultiParameterPlasticityStressUpdate.C:725
MultiParameterPlasticityStressUpdate::_max_iter_used
MaterialProperty< Real > & _max_iter_used
Maximum number of Newton-Raphson iterations used in the return-map during the course of the entire si...
Definition: MultiParameterPlasticityStressUpdate.h:202
StressUpdateBase
StressUpdateBase is a material that is not called by MOOSE because of the compute=false flag set in t...
Definition: StressUpdateBase.h:52
MultiParameterPlasticityStressUpdate::SmootherFunctionType::poly2
MultiParameterPlasticityStressUpdate::setIntnlValuesV
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,...
MultiParameterPlasticityStressUpdate::_dvar_dtrial
std::vector< std::vector< Real > > _dvar_dtrial
d({stress_param[i], gaE})/d(trial_stress_param[j])
Definition: MultiParameterPlasticityStressUpdate.h:703
MultiParameterPlasticityStressUpdate::smoothAllQuantities
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.
Definition: MultiParameterPlasticityStressUpdate.C:400
MultiParameterPlasticityStressUpdate::calculateRHS
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[...
Definition: MultiParameterPlasticityStressUpdate.C:811
TangentCalculationMethod
TangentCalculationMethod
TangentCalculationMethod is an enum that determines the calculation method for the tangent operator.
Definition: StressUpdateBase.h:30
RankTwoTensorTempl< Real >
MultiParameterPlasticityStressUpdate::yieldAndFlow::d2g_di
std::vector< std::vector< Real > > d2g_di
Definition: MultiParameterPlasticityStressUpdate.h:221
StressUpdateBase.h
MultiParameterPlasticityStressUpdate::_ok_sp
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".
Definition: MultiParameterPlasticityStressUpdate.h:713
MultiParameterPlasticityStressUpdate::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 =0
Sets stress from the admissible parameters.
MultiParameterPlasticityStressUpdate::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)
Derived classes may employ this function to record stuff or do other computations prior to the return...
Definition: MultiParameterPlasticityStressUpdate.C:678
MultiParameterPlasticityStressUpdate
MultiParameterPlasticityStressUpdate performs the return-map algorithm and associated stress updates ...
Definition: MultiParameterPlasticityStressUpdate.h:99
MultiParameterPlasticityStressUpdate::yieldAndFlow::df_di
std::vector< Real > df_di
Definition: MultiParameterPlasticityStressUpdate.h:218
MultiParameterPlasticityStressUpdate::_smoothing_tol2
const Real _smoothing_tol2
Square of the smoothing tolerance.
Definition: MultiParameterPlasticityStressUpdate.h:162
MultiParameterPlasticityStressUpdate::nrStep
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.
Definition: MultiParameterPlasticityStressUpdate.C:637
MultiParameterPlasticityStressUpdate::_tensor_dimensionality
constexpr static unsigned _tensor_dimensionality
Internal dimensionality of tensors (currently this is 3 throughout tensor_mechanics)
Definition: MultiParameterPlasticityStressUpdate.h:129
MultiParameterPlasticityStressUpdate::_intnl
MaterialProperty< std::vector< Real > > & _intnl
internal parameters
Definition: MultiParameterPlasticityStressUpdate.h:190
MultiParameterPlasticityStressUpdate::_max_iter_used_old
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 ...
Definition: MultiParameterPlasticityStressUpdate.h:205
MultiParameterPlasticityStressUpdate::getTangentCalculationMethod
virtual TangentCalculationMethod getTangentCalculationMethod() override
Definition: MultiParameterPlasticityStressUpdate.h:123
MultiParameterPlasticityStressUpdate::setInelasticStrainIncrementAfterReturn
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...
Definition: MultiParameterPlasticityStressUpdate.C:786
MultiParameterPlasticityStressUpdate::yieldAndFlow::operator<
bool operator<(const yieldAndFlow &fd) const
Definition: MultiParameterPlasticityStressUpdate.h:236
MultiParameterPlasticityStressUpdate::_f_tol
const Real _f_tol
The yield-function tolerance.
Definition: MultiParameterPlasticityStressUpdate.h:165