www.mooseframework.org
ComputeMultiPlasticityStress.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 "ComputeStressBase.h"
14 
16 
17 template <>
19 
31 {
32 public:
33  static InputParameters validParams();
34 
35  ComputeMultiPlasticityStress(const InputParameters & parameters);
36 
37 protected:
38  virtual void computeQpStress();
39  virtual void initQpStatefulProperties();
40 
42  unsigned int _max_iter;
43 
46 
49 
52 
55  {
60 
62  Real _epp_tol;
63 
65  std::vector<Real> _dummy_pm;
66 
71  std::vector<Real> _cumulative_pm;
72 
73  /*
74  * Scheme by which constraints are deactivated.
75  * This enum is defined here for computational
76  * efficiency. If you add to this list you must
77  * also add to the MooseEnum in the .C file
78  */
80  {
89 
92 
94  RealVectorValue _n_input;
95 
97  RealTensorValue _rot;
98 
101 
103  const std::string _elasticity_tensor_name;
105  const MaterialProperty<RankFourTensor> & _elasticity_tensor;
106 
108  MaterialProperty<RankTwoTensor> & _plastic_strain;
109 
111  const MaterialProperty<RankTwoTensor> & _plastic_strain_old;
112 
114  MaterialProperty<std::vector<Real>> & _intnl;
115 
117  const MaterialProperty<std::vector<Real>> & _intnl_old;
118 
120  MaterialProperty<std::vector<Real>> & _yf;
121 
123  MaterialProperty<Real> & _iter;
124 
126  MaterialProperty<Real> & _linesearch_needed;
127 
129  MaterialProperty<Real> & _ld_encountered;
130 
132  MaterialProperty<Real> & _constraints_added;
133 
135  MaterialProperty<RealVectorValue> & _n;
136 
138  const MaterialProperty<RealVectorValue> & _n_old;
139 
141  const MaterialProperty<RankTwoTensor> & _strain_increment;
142 
144  const MaterialProperty<RankTwoTensor> & _total_strain_old;
145 
147  const MaterialProperty<RankTwoTensor> & _rotation_increment;
148 
150  const MaterialProperty<RankTwoTensor> & _stress_old;
151 
153  const MaterialProperty<RankTwoTensor> & _elastic_strain_old;
154 
156  bool _cosserat;
157 
159  const MaterialProperty<RankTwoTensor> * _curvature;
160 
162  const MaterialProperty<RankFourTensor> * _elastic_flexural_rigidity_tensor;
163 
165  MaterialProperty<RankTwoTensor> * _couple_stress;
166 
168  const MaterialProperty<RankTwoTensor> * _couple_stress_old;
169 
171  MaterialProperty<RankFourTensor> * _Jacobian_mult_couple;
172 
175 
178 
181 
184 
188  virtual bool reinstateLinearDependentConstraints(std::vector<bool> & deactivated_due_to_ld);
189 
193  virtual unsigned int numberActive(const std::vector<bool> & active);
194 
206  virtual Real residual2(const std::vector<Real> & pm,
207  const std::vector<Real> & f,
208  const RankTwoTensor & epp,
209  const std::vector<Real> & ic,
210  const std::vector<bool> & active,
211  const std::vector<bool> & deactivated_due_to_ld);
212 
249  virtual bool returnMap(const RankTwoTensor & stress_old,
250  RankTwoTensor & stress,
251  const std::vector<Real> & intnl_old,
252  std::vector<Real> & intnl,
253  const RankTwoTensor & plastic_strain_old,
254  RankTwoTensor & plastic_strain,
255  const RankFourTensor & E_ijkl,
256  const RankTwoTensor & strain_increment,
257  std::vector<Real> & f,
258  unsigned int & iter,
259  bool can_revert_to_dumb,
260  bool & linesearch_needed,
261  bool & ld_encountered,
262  bool & constraints_added,
263  bool final_step,
264  RankFourTensor & consistent_tangent_operator,
265  std::vector<Real> & cumulative_pm);
266 
299  virtual bool lineSearch(Real & nr_res2,
300  RankTwoTensor & stress,
301  const std::vector<Real> & intnl_old,
302  std::vector<Real> & intnl,
303  std::vector<Real> & pm,
304  const RankFourTensor & E_inv,
305  RankTwoTensor & delta_dp,
306  const RankTwoTensor & dstress,
307  const std::vector<Real> & dpm,
308  const std::vector<Real> & dintnl,
309  std::vector<Real> & f,
310  RankTwoTensor & epp,
311  std::vector<Real> & ic,
312  const std::vector<bool> & active,
313  const std::vector<bool> & deactivated_due_to_ld,
314  bool & linesearch_needed);
315 
342  virtual bool singleStep(Real & nr_res2,
343  RankTwoTensor & stress,
344  const std::vector<Real> & intnl_old,
345  std::vector<Real> & intnl,
346  std::vector<Real> & pm,
347  RankTwoTensor & delta_dp,
348  const RankFourTensor & E_inv,
349  std::vector<Real> & f,
350  RankTwoTensor & epp,
351  std::vector<Real> & ic,
352  std::vector<bool> & active,
353  DeactivationSchemeEnum deactivation_scheme,
354  bool & linesearch_needed,
355  bool & ld_encountered);
356 
364  virtual bool checkAdmissible(const RankTwoTensor & stress,
365  const std::vector<Real> & intnl,
366  std::vector<Real> & all_f);
367 
377  void buildDumbOrder(const RankTwoTensor & stress,
378  const std::vector<Real> & intnl,
379  std::vector<unsigned int> & dumb_order);
380 
392  virtual void incrementDumb(int & dumb_iteration,
393  const std::vector<unsigned int> & dumb_order,
394  std::vector<bool> & act);
395 
413  virtual bool checkKuhnTucker(const std::vector<Real> & f,
414  const std::vector<Real> & pm,
415  const std::vector<bool> & active);
416 
434  virtual void applyKuhnTucker(const std::vector<Real> & f,
435  const std::vector<Real> & pm,
436  std::vector<bool> & active);
437 
438  // gets called before any return-map
439  virtual void preReturnMap();
440 
441  // gets called after return-map
442  virtual void postReturnMap();
443 
446  {
449  };
450 
489  virtual bool quickStep(const RankTwoTensor & stress_old,
490  RankTwoTensor & stress,
491  const std::vector<Real> & intnl_old,
492  std::vector<Real> & intnl,
493  std::vector<Real> & pm,
494  std::vector<Real> & cumulative_pm,
495  const RankTwoTensor & plastic_strain_old,
496  RankTwoTensor & plastic_strain,
497  const RankFourTensor & E_ijkl,
498  const RankTwoTensor & strain_increment,
499  std::vector<Real> & yf,
500  unsigned int & iterations,
501  RankFourTensor & consistent_tangent_operator,
502  const quickStep_called_from_t called_from,
503  bool final_step);
504 
531  virtual bool plasticStep(const RankTwoTensor & stress_old,
532  RankTwoTensor & stress,
533  const std::vector<Real> & intnl_old,
534  std::vector<Real> & intnl,
535  const RankTwoTensor & plastic_strain_old,
536  RankTwoTensor & plastic_strain,
537  const RankFourTensor & E_ijkl,
538  const RankTwoTensor & strain_increment,
539  std::vector<Real> & yf,
540  unsigned int & iterations,
541  bool & linesearch_needed,
542  bool & ld_encountered,
543  bool & constraints_added,
544  RankFourTensor & consistent_tangent_operator);
545 
546  // bool checkAndModifyConstraints(bool nr_exit_condition, const RankTwoTensor & stress, const
547  // std::vector<Real> & intnl, const std::vector<Real> & pm, const std::vector<bool> &
548  // initial_act, bool can_revert_to_dumb, const RankTwoTensor & initial_stress, const
549  // std::vector<Real> & intnl_old, const std::vector<Real> & f, DeactivationSchemeEnum
550  // deact_scheme, std::vector<bool> & act, int & dumb_iteration, std::vector<unsigned int>
551  // dumb_order, bool & die);
552 
553  bool canChangeScheme(DeactivationSchemeEnum current_deactivation_scheme, bool can_revert_to_dumb);
554 
555  bool canIncrementDumb(int dumb_iteration);
556 
557  void changeScheme(const std::vector<bool> & initial_act,
558  bool can_revert_to_dumb,
559  const RankTwoTensor & initial_stress,
560  const std::vector<Real> & intnl_old,
561  DeactivationSchemeEnum & current_deactivation_scheme,
562  std::vector<bool> & act,
563  int & dumb_iteration,
564  std::vector<unsigned int> & dumb_order);
565 
566  bool canAddConstraints(const std::vector<bool> & act, const std::vector<Real> & all_f);
567 
568  unsigned int activeCombinationNumber(const std::vector<bool> & act);
569 
588  const std::vector<Real> & intnl,
589  const RankFourTensor & E_ijkl,
590  const std::vector<Real> & pm_this_step,
591  const std::vector<Real> & cumulative_pm);
592 
593 private:
594  RankTwoTensor rot(const RankTwoTensor & tens);
595 };
ComputeMultiPlasticityStress::_n
MaterialProperty< RealVectorValue > & _n
current value of transverse direction
Definition: ComputeMultiPlasticityStress.h:135
ComputeMultiPlasticityStress::safe
Definition: ComputeMultiPlasticityStress.h:82
ComputeMultiPlasticityStress::_linesearch_needed
MaterialProperty< Real > & _linesearch_needed
Whether a line-search was needed in the latest Newton-Raphson process (1 if true, 0 otherwise)
Definition: ComputeMultiPlasticityStress.h:126
ComputeMultiPlasticityStress::_dummy_pm
std::vector< Real > _dummy_pm
dummy "consistency parameters" (plastic multipliers) used in quickStep when called from computeQpStre...
Definition: ComputeMultiPlasticityStress.h:65
ComputeMultiPlasticityStress::elastic
Definition: ComputeMultiPlasticityStress.h:56
ComputeMultiPlasticityStress::_ld_encountered
MaterialProperty< Real > & _ld_encountered
Whether linear-dependence was encountered in the latest Newton-Raphson process (1 if true,...
Definition: ComputeMultiPlasticityStress.h:129
ComputeMultiPlasticityStress::changeScheme
void changeScheme(const std::vector< bool > &initial_act, bool can_revert_to_dumb, const RankTwoTensor &initial_stress, const std::vector< Real > &intnl_old, DeactivationSchemeEnum &current_deactivation_scheme, std::vector< bool > &act, int &dumb_iteration, std::vector< unsigned int > &dumb_order)
Definition: ComputeMultiPlasticityStress.C:1050
ComputeMultiPlasticityStress::checkKuhnTucker
virtual bool checkKuhnTucker(const std::vector< Real > &f, const std::vector< Real > &pm, const std::vector< bool > &active)
Checks Kuhn-Tucker conditions, and alters "active" if appropriate.
Definition: ComputeMultiPlasticityStress.C:1288
ComputeMultiPlasticityStress::lineSearch
virtual bool lineSearch(Real &nr_res2, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, std::vector< Real > &pm, const RankFourTensor &E_inv, RankTwoTensor &delta_dp, const RankTwoTensor &dstress, const std::vector< Real > &dpm, const std::vector< Real > &dintnl, std::vector< Real > &f, RankTwoTensor &epp, std::vector< Real > &ic, const std::vector< bool > &active, const std::vector< bool > &deactivated_due_to_ld, bool &linesearch_needed)
Performs a line search.
Definition: ComputeMultiPlasticityStress.C:1391
ComputeMultiPlasticityStress::rot
RankTwoTensor rot(const RankTwoTensor &tens)
Definition: ComputeMultiPlasticityStress.C:312
ComputeMultiPlasticityStress::canIncrementDumb
bool canIncrementDumb(int dumb_iteration)
Definition: ComputeMultiPlasticityStress.C:1567
ComputeMultiPlasticityStress
ComputeMultiPlasticityStress performs the return-map algorithm and associated stress updates for plas...
Definition: ComputeMultiPlasticityStress.h:30
ComputeMultiPlasticityStress::_tangent_operator_type
enum ComputeMultiPlasticityStress::TangentOperatorEnum _tangent_operator_type
ComputeMultiPlasticityStress::returnMap
virtual bool returnMap(const RankTwoTensor &stress_old, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, const RankTwoTensor &plastic_strain_old, RankTwoTensor &plastic_strain, const RankFourTensor &E_ijkl, const RankTwoTensor &strain_increment, std::vector< Real > &f, unsigned int &iter, bool can_revert_to_dumb, bool &linesearch_needed, bool &ld_encountered, bool &constraints_added, bool final_step, RankFourTensor &consistent_tangent_operator, std::vector< Real > &cumulative_pm)
Implements the return map.
Definition: ComputeMultiPlasticityStress.C:615
ComputeMultiPlasticityStress::incrementDumb
virtual void incrementDumb(int &dumb_iteration, const std::vector< unsigned int > &dumb_order, std::vector< bool > &act)
Increments "dumb_iteration" by 1, and sets "act" appropriately (act[alpha] = true iff alpha_th bit of...
Definition: ComputeMultiPlasticityStress.C:1555
ComputeMultiPlasticityStress::DeactivationSchemeEnum
DeactivationSchemeEnum
Definition: ComputeMultiPlasticityStress.h:79
ComputeMultiPlasticityStress::_n_input
RealVectorValue _n_input
the supplied transverse direction vector
Definition: ComputeMultiPlasticityStress.h:94
ComputeMultiPlasticityStress::_max_stepsize_for_dumb
Real _max_stepsize_for_dumb
"dumb" deactivation will only be used if the stepsize falls below this quantity
Definition: ComputeMultiPlasticityStress.h:48
ComputeMultiPlasticityStress::_couple_stress_old
const MaterialProperty< RankTwoTensor > * _couple_stress_old
the old value of Cosserat couple-stress
Definition: ComputeMultiPlasticityStress.h:168
ComputeMultiPlasticityStress::_epp_tol
Real _epp_tol
Tolerance on the plastic strain increment ("direction") constraint.
Definition: ComputeMultiPlasticityStress.h:62
ComputeMultiPlasticityStress::_intnl_old
const MaterialProperty< std::vector< Real > > & _intnl_old
old values of internal parameters
Definition: ComputeMultiPlasticityStress.h:117
ComputeMultiPlasticityStress::checkAdmissible
virtual bool checkAdmissible(const RankTwoTensor &stress, const std::vector< Real > &intnl, std::vector< Real > &all_f)
Checks whether the yield functions are in the admissible region.
Definition: ComputeMultiPlasticityStress.C:1271
ComputeMultiPlasticityStress::_min_stepsize
Real _min_stepsize
Minimum fraction of applied strain that may be applied during adaptive stepsizing.
Definition: ComputeMultiPlasticityStress.h:45
ComputeMultiPlasticityStress::ComputeMultiPlasticityStress
ComputeMultiPlasticityStress(const InputParameters &parameters)
Definition: ComputeMultiPlasticityStress.C:99
ComputeMultiPlasticityStress::residual2
virtual Real residual2(const std::vector< Real > &pm, const std::vector< Real > &f, const RankTwoTensor &epp, const std::vector< Real > &ic, const std::vector< bool > &active, const std::vector< bool > &deactivated_due_to_ld)
The residual-squared.
Definition: ComputeMultiPlasticityStress.C:1354
ComputeMultiPlasticityStress::dumb
Definition: ComputeMultiPlasticityStress.h:83
ComputeMultiPlasticityStress::_cumulative_pm
std::vector< Real > _cumulative_pm
the sum of the plastic multipliers over all the sub-steps.
Definition: ComputeMultiPlasticityStress.h:71
ComputeMultiPlasticityStress::_plastic_strain
MaterialProperty< RankTwoTensor > & _plastic_strain
plastic strain
Definition: ComputeMultiPlasticityStress.h:108
ComputeMultiPlasticityStress::computeQpStress
virtual void computeQpStress()
Compute the stress and store it in the _stress material property for the current quadrature point.
Definition: ComputeMultiPlasticityStress.C:221
ComputeMultiPlasticityStress::_my_curvature
RankTwoTensor _my_curvature
Curvature that can be rotated by this class, and split into multiple increments (ie,...
Definition: ComputeMultiPlasticityStress.h:183
ComputeMultiPlasticityStress::returnMap_function
Definition: ComputeMultiPlasticityStress.h:448
ComputeMultiPlasticityStress::_curvature
const MaterialProperty< RankTwoTensor > * _curvature
The Cosserat curvature strain.
Definition: ComputeMultiPlasticityStress.h:159
ComputeMultiPlasticityStress::consistentTangentOperator
RankFourTensor consistentTangentOperator(const RankTwoTensor &stress, const std::vector< Real > &intnl, const RankFourTensor &E_ijkl, const std::vector< Real > &pm_this_step, const std::vector< Real > &cumulative_pm)
Computes the consistent tangent operator (another name for the jacobian = d(stress_rate)/d(strain_rat...
Definition: ComputeMultiPlasticityStress.C:1585
ComputeMultiPlasticityStress::quickStep
virtual bool quickStep(const RankTwoTensor &stress_old, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, std::vector< Real > &pm, std::vector< Real > &cumulative_pm, const RankTwoTensor &plastic_strain_old, RankTwoTensor &plastic_strain, const RankFourTensor &E_ijkl, const RankTwoTensor &strain_increment, std::vector< Real > &yf, unsigned int &iterations, RankFourTensor &consistent_tangent_operator, const quickStep_called_from_t called_from, bool final_step)
Attempts to find an admissible (stress, intnl) by using the customized return-map algorithms defined ...
Definition: ComputeMultiPlasticityStress.C:371
ComputeMultiPlasticityStress::_deactivation_scheme
enum ComputeMultiPlasticityStress::DeactivationSchemeEnum _deactivation_scheme
ComputeStressBase
ComputeStressBase is the base class for stress tensors.
Definition: ComputeStressBase.h:26
ComputeMultiPlasticityStress::_elasticity_tensor
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
Definition: ComputeMultiPlasticityStress.h:105
ComputeMultiPlasticityStress::_plastic_strain_old
const MaterialProperty< RankTwoTensor > & _plastic_strain_old
Old value of plastic strain.
Definition: ComputeMultiPlasticityStress.h:111
ComputeMultiPlasticityStress::plasticStep
virtual bool plasticStep(const RankTwoTensor &stress_old, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, const RankTwoTensor &plastic_strain_old, RankTwoTensor &plastic_strain, const RankFourTensor &E_ijkl, const RankTwoTensor &strain_increment, std::vector< Real > &yf, unsigned int &iterations, bool &linesearch_needed, bool &ld_encountered, bool &constraints_added, RankFourTensor &consistent_tangent_operator)
performs a plastic step
Definition: ComputeMultiPlasticityStress.C:468
ComputeMultiPlasticityStress::_couple_stress
MaterialProperty< RankTwoTensor > * _couple_stress
the Cosserat couple-stress
Definition: ComputeMultiPlasticityStress.h:165
ComputeMultiPlasticityStress::_rot
RealTensorValue _rot
rotation matrix that takes _n to (0, 0, 1)
Definition: ComputeMultiPlasticityStress.h:97
ComputeMultiPlasticityStress::buildDumbOrder
void buildDumbOrder(const RankTwoTensor &stress, const std::vector< Real > &intnl, std::vector< unsigned int > &dumb_order)
Builds the order which "dumb" activation will take.
Definition: ComputeMultiPlasticityStress.C:1524
ComputeMultiPlasticityStress::_Jacobian_mult_couple
MaterialProperty< RankFourTensor > * _Jacobian_mult_couple
derivative of couple-stress w.r.t. curvature
Definition: ComputeMultiPlasticityStress.h:171
ComputeMultiPlasticityStress::numberActive
virtual unsigned int numberActive(const std::vector< bool > &active)
counts the number of active constraints
Definition: ComputeMultiPlasticityStress.C:1261
ComputeMultiPlasticityStress::_cosserat
bool _cosserat
whether Cosserat mechanics should be used
Definition: ComputeMultiPlasticityStress.h:156
ComputeMultiPlasticityStress::validParams
static InputParameters validParams()
Definition: ComputeMultiPlasticityStress.C:24
ComputeMultiPlasticityStress::_total_strain_old
const MaterialProperty< RankTwoTensor > & _total_strain_old
Old value of total strain (coming from ComputeIncrementalSmallStrain, for example)
Definition: ComputeMultiPlasticityStress.h:144
ComputeMultiPlasticityStress::_elastic_flexural_rigidity_tensor
const MaterialProperty< RankFourTensor > * _elastic_flexural_rigidity_tensor
The Cosserat elastic flexural rigidity tensor.
Definition: ComputeMultiPlasticityStress.h:162
ComputeMultiPlasticityStress::_perform_finite_strain_rotations
bool _perform_finite_strain_rotations
whether to perform the rotations necessary in finite-strain simulations
Definition: ComputeMultiPlasticityStress.h:100
validParams< ComputeMultiPlasticityStress >
InputParameters validParams< ComputeMultiPlasticityStress >()
ComputeMultiPlasticityStress::optimized_to_dumb
Definition: ComputeMultiPlasticityStress.h:87
ComputeMultiPlasticityStress::optimized_to_safe
Definition: ComputeMultiPlasticityStress.h:84
ComputeMultiPlasticityStress::reinstateLinearDependentConstraints
virtual bool reinstateLinearDependentConstraints(std::vector< bool > &deactivated_due_to_ld)
makes all deactivated_due_to_ld false, and if >0 of them were initially true, returns true
Definition: ComputeMultiPlasticityStress.C:1248
ComputeMultiPlasticityStress::quickStep_called_from_t
quickStep_called_from_t
The functions from which quickStep can be called.
Definition: ComputeMultiPlasticityStress.h:445
ComputeMultiPlasticityStress::_rotation_increment
const MaterialProperty< RankTwoTensor > & _rotation_increment
Rotation increment (coming from ComputeIncrementalSmallStrain, for example)
Definition: ComputeMultiPlasticityStress.h:147
ComputeMultiPlasticityStress::singleStep
virtual bool singleStep(Real &nr_res2, RankTwoTensor &stress, const std::vector< Real > &intnl_old, std::vector< Real > &intnl, std::vector< Real > &pm, RankTwoTensor &delta_dp, const RankFourTensor &E_inv, std::vector< Real > &f, RankTwoTensor &epp, std::vector< Real > &ic, std::vector< bool > &active, DeactivationSchemeEnum deactivation_scheme, bool &linesearch_needed, bool &ld_encountered)
Performs a single Newton-Raphson + linesearch step Constraints are deactivated and the step is re-don...
Definition: ComputeMultiPlasticityStress.C:1082
ComputeMultiPlasticityStress::_constraints_added
MaterialProperty< Real > & _constraints_added
Whether constraints were added in during the latest Newton-Raphson process (1 if true,...
Definition: ComputeMultiPlasticityStress.h:132
ComputeMultiPlasticityStress::preReturnMap
virtual void preReturnMap()
Definition: ComputeMultiPlasticityStress.C:320
ComputeMultiPlasticityStress::optimized
Definition: ComputeMultiPlasticityStress.h:81
ComputeMultiPlasticityStress::computeQpStress_function
Definition: ComputeMultiPlasticityStress.h:447
ComputeMultiPlasticityStress::safe_to_dumb
Definition: ComputeMultiPlasticityStress.h:85
MultiPlasticityDebugger.h
ComputeMultiPlasticityStress::canAddConstraints
bool canAddConstraints(const std::vector< bool > &act, const std::vector< Real > &all_f)
Definition: ComputeMultiPlasticityStress.C:1015
ComputeMultiPlasticityStress::_strain_increment
const MaterialProperty< RankTwoTensor > & _strain_increment
strain increment (coming from ComputeIncrementalSmallStrain, for example)
Definition: ComputeMultiPlasticityStress.h:141
MultiPlasticityDebugger
MultiPlasticityDebugger computes various finite-difference things to help developers remove bugs in t...
Definition: MultiPlasticityDebugger.h:24
ComputeMultiPlasticityStress::_n_old
const MaterialProperty< RealVectorValue > & _n_old
old value of transverse direction
Definition: ComputeMultiPlasticityStress.h:138
ComputeMultiPlasticityStress::_stress_old
const MaterialProperty< RankTwoTensor > & _stress_old
Old value of stress.
Definition: ComputeMultiPlasticityStress.h:150
ComputeMultiPlasticityStress::initQpStatefulProperties
virtual void initQpStatefulProperties()
Definition: ComputeMultiPlasticityStress.C:191
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
ComputeMultiPlasticityStress::_intnl
MaterialProperty< std::vector< Real > > & _intnl
internal parameters
Definition: ComputeMultiPlasticityStress.h:114
ComputeMultiPlasticityStress::_max_iter
unsigned int _max_iter
Maximum number of Newton-Raphson iterations allowed.
Definition: ComputeMultiPlasticityStress.h:42
ComputeStressBase.h
ComputeMultiPlasticityStress::canChangeScheme
bool canChangeScheme(DeactivationSchemeEnum current_deactivation_scheme, bool can_revert_to_dumb)
Definition: ComputeMultiPlasticityStress.C:1025
ComputeMultiPlasticityStress::activeCombinationNumber
unsigned int activeCombinationNumber(const std::vector< bool > &act)
Definition: ComputeMultiPlasticityStress.C:1574
ComputeMultiPlasticityStress::_ignore_failures
bool _ignore_failures
Even if the returnMap fails, return the best values found for stress and internal parameters.
Definition: ComputeMultiPlasticityStress.h:51
ComputeMultiPlasticityStress::_my_flexural_rigidity_tensor
RankFourTensor _my_flexural_rigidity_tensor
Flexual rigidity tensor that can be rotated by this class (ie, its not const)
Definition: ComputeMultiPlasticityStress.h:180
ComputeMultiPlasticityStress::TangentOperatorEnum
TangentOperatorEnum
The type of tangent operator to return. tangent operator = d(stress_rate)/d(strain_rate).
Definition: ComputeMultiPlasticityStress.h:54
ComputeMultiPlasticityStress::_my_elasticity_tensor
RankFourTensor _my_elasticity_tensor
Elasticity tensor that can be rotated by this class (ie, its not const)
Definition: ComputeMultiPlasticityStress.h:174
ComputeMultiPlasticityStress::_iter
MaterialProperty< Real > & _iter
Number of Newton-Raphson iterations used in the return-map.
Definition: ComputeMultiPlasticityStress.h:123
RankTwoTensorTempl< Real >
ComputeMultiPlasticityStress::linear
Definition: ComputeMultiPlasticityStress.h:57
ComputeMultiPlasticityStress::_my_strain_increment
RankTwoTensor _my_strain_increment
Strain increment that can be rotated by this class, and split into multiple increments (ie,...
Definition: ComputeMultiPlasticityStress.h:177
ComputeMultiPlasticityStress::_elasticity_tensor_name
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
Definition: ComputeMultiPlasticityStress.h:103
ComputeMultiPlasticityStress::applyKuhnTucker
virtual void applyKuhnTucker(const std::vector< Real > &f, const std::vector< Real > &pm, std::vector< bool > &active)
Checks Kuhn-Tucker conditions, and alters "active" if appropriate.
Definition: ComputeMultiPlasticityStress.C:1313
ComputeMultiPlasticityStress::_elastic_strain_old
const MaterialProperty< RankTwoTensor > & _elastic_strain_old
Old value of elastic strain.
Definition: ComputeMultiPlasticityStress.h:153
ComputeMultiPlasticityStress::nonlinear
Definition: ComputeMultiPlasticityStress.h:58
ComputeMultiPlasticityStress::postReturnMap
virtual void postReturnMap()
Definition: ComputeMultiPlasticityStress.C:339
ComputeMultiPlasticityStress::optimized_to_safe_to_dumb
Definition: ComputeMultiPlasticityStress.h:86
ComputeMultiPlasticityStress::_n_supplied
bool _n_supplied
User supplied the transverse direction vector.
Definition: ComputeMultiPlasticityStress.h:91
ComputeMultiPlasticityStress::_yf
MaterialProperty< std::vector< Real > > & _yf
yield functions
Definition: ComputeMultiPlasticityStress.h:120