https://mooseframework.inl.gov
TensileStressUpdate.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 
14 
20 {
21 public:
23 
25 
29  bool requiresIsotropicTensor() override { return true; }
30 
31 protected:
34 
36  const bool _perfect_guess;
37 
40 
42  virtual Real tensile_strength(const Real internal_param) const;
43 
45  virtual Real dtensile_strength(const Real internal_param) const;
46 
47  void computeStressParams(const RankTwoTensor & stress,
48  std::vector<Real> & stress_params) const override;
49 
50  std::vector<RankTwoTensor> dstress_param_dstress(const RankTwoTensor & stress) const override;
51 
52  std::vector<RankFourTensor> d2stress_param_dstress(const RankTwoTensor & stress) const override;
53 
54  virtual void setStressAfterReturnV(const RankTwoTensor & stress_trial,
55  const std::vector<Real> & stress_params,
56  Real gaE,
57  const std::vector<Real> & intnl,
58  const yieldAndFlow & smoothed_q,
59  const RankFourTensor & Eijkl,
60  RankTwoTensor & stress) const override;
61 
62  virtual void preReturnMapV(const std::vector<Real> & trial_stress_params,
63  const RankTwoTensor & stress_trial,
64  const std::vector<Real> & intnl_old,
65  const std::vector<Real> & yf,
66  const RankFourTensor & Eijkl) override;
67 
68  void setEffectiveElasticity(const RankFourTensor & Eijkl) override;
69 
70  void yieldFunctionValuesV(const std::vector<Real> & stress_params,
71  const std::vector<Real> & intnl,
72  std::vector<Real> & yf) const override;
73 
74  void computeAllQV(const std::vector<Real> & stress_params,
75  const std::vector<Real> & intnl,
76  std::vector<yieldAndFlow> & all_q) const override;
77 
78  void initializeVarsV(const std::vector<Real> & trial_stress_params,
79  const std::vector<Real> & intnl_old,
80  std::vector<Real> & stress_params,
81  Real & gaE,
82  std::vector<Real> & intnl) const override;
83 
84  void setIntnlValuesV(const std::vector<Real> & trial_stress_params,
85  const std::vector<Real> & current_stress_params,
86  const std::vector<Real> & intnl_old,
87  std::vector<Real> & intnl) const override;
88 
89  void setIntnlDerivativesV(const std::vector<Real> & trial_stress_params,
90  const std::vector<Real> & current_stress_params,
91  const std::vector<Real> & intnl,
92  std::vector<std::vector<Real>> & dintnl) const override;
93 
94  virtual void consistentTangentOperatorV(const RankTwoTensor & stress_trial,
95  const std::vector<Real> & trial_stress_params,
96  const RankTwoTensor & stress,
97  const std::vector<Real> & stress_params,
98  Real gaE,
99  const yieldAndFlow & smoothed_q,
100  const RankFourTensor & Eijkl,
101  bool compute_full_tangent_operator,
102  const std::vector<std::vector<Real>> & dvar_dtrial,
103  RankFourTensor & cto) override;
104 };
TensileStressUpdate(const InputParameters &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) override
Derived classes may employ this function to record stuff or do other computations prior to the return...
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.
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 ...
Struct designed to hold info about a single yield function and its derivatives, as well as the flow d...
TensileStressUpdate implements rate-independent associative tensile failure ("Rankine" plasticity) wi...
void setEffectiveElasticity(const RankFourTensor &Eijkl) override
Sets _Eij and _En and _Cij.
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.
const SolidMechanicsHardeningModel & _strength
Hardening model for tensile strength.
RankTwoTensor _eigvecs
Eigenvectors of the trial stress as a RankTwoTensor, in order to rotate the returned stress back to s...
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, their current values, and the current values of the internal parameters.
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.
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.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< RankFourTensor > d2stress_param_dstress(const RankTwoTensor &stress) const override
d2(stress_param[i])/d(stress)/d(stress) at given stress
const bool _perfect_guess
Whether to provide an estimate of the returned stress, based on perfect plasticity.
MultiParameterPlasticityStressUpdate performs the return-map algorithm and associated stress updates ...
static InputParameters validParams()
std::vector< RankTwoTensor > dstress_param_dstress(const RankTwoTensor &stress) const override
d(stress_param[i])/d(stress) at given stress
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...
bool requiresIsotropicTensor() override
Does the model require the elasticity tensor to be isotropic?
const InputParameters & parameters() const
void computeStressParams(const RankTwoTensor &stress, std::vector< Real > &stress_params) const override
Computes stress_params, given stress.
virtual Real tensile_strength(const Real internal_param) const
tensile strength as a function of residual value, rate, and internal_param
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.