www.mooseframework.org
TensileStressUpdate.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 
14 
16 
17 template <>
18 InputParameters validParams<TensileStressUpdate>();
19 
25 {
26 public:
27  static InputParameters validParams();
28 
29  TensileStressUpdate(const InputParameters & parameters);
30 
34  bool requiresIsotropicTensor() override { return true; }
35 
36 protected:
39 
41  const bool _perfect_guess;
42 
45 
47  virtual Real tensile_strength(const Real internal_param) const;
48 
50  virtual Real dtensile_strength(const Real internal_param) const;
51 
52  void computeStressParams(const RankTwoTensor & stress,
53  std::vector<Real> & stress_params) const override;
54 
55  std::vector<RankTwoTensor> dstress_param_dstress(const RankTwoTensor & stress) const override;
56 
57  std::vector<RankFourTensor> d2stress_param_dstress(const RankTwoTensor & stress) const override;
58 
59  virtual void setStressAfterReturnV(const RankTwoTensor & stress_trial,
60  const std::vector<Real> & stress_params,
61  Real gaE,
62  const std::vector<Real> & intnl,
63  const yieldAndFlow & smoothed_q,
64  const RankFourTensor & Eijkl,
65  RankTwoTensor & stress) const override;
66 
67  virtual void preReturnMapV(const std::vector<Real> & trial_stress_params,
68  const RankTwoTensor & stress_trial,
69  const std::vector<Real> & intnl_old,
70  const std::vector<Real> & yf,
71  const RankFourTensor & Eijkl) override;
72 
73  void setEffectiveElasticity(const RankFourTensor & Eijkl) override;
74 
75  void yieldFunctionValuesV(const std::vector<Real> & stress_params,
76  const std::vector<Real> & intnl,
77  std::vector<Real> & yf) const override;
78 
79  void computeAllQV(const std::vector<Real> & stress_params,
80  const std::vector<Real> & intnl,
81  std::vector<yieldAndFlow> & all_q) const override;
82 
83  void initializeVarsV(const std::vector<Real> & trial_stress_params,
84  const std::vector<Real> & intnl_old,
85  std::vector<Real> & stress_params,
86  Real & gaE,
87  std::vector<Real> & intnl) const override;
88 
89  void setIntnlValuesV(const std::vector<Real> & trial_stress_params,
90  const std::vector<Real> & current_stress_params,
91  const std::vector<Real> & intnl_old,
92  std::vector<Real> & intnl) const override;
93 
94  void setIntnlDerivativesV(const std::vector<Real> & trial_stress_params,
95  const std::vector<Real> & current_stress_params,
96  const std::vector<Real> & intnl,
97  std::vector<std::vector<Real>> & dintnl) const override;
98 
99  virtual void consistentTangentOperatorV(const RankTwoTensor & stress_trial,
100  const std::vector<Real> & trial_stress_params,
101  const RankTwoTensor & stress,
102  const std::vector<Real> & stress_params,
103  Real gaE,
104  const yieldAndFlow & smoothed_q,
105  const RankFourTensor & Eijkl,
106  bool compute_full_tangent_operator,
107  const std::vector<std::vector<Real>> & dvar_dtrial,
108  RankFourTensor & cto) override;
109 };
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
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
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
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
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
validParams< TensileStressUpdate >
InputParameters validParams< TensileStressUpdate >()
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
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
TensileStressUpdate::requiresIsotropicTensor
bool requiresIsotropicTensor() override
Does the model require the elasticity tensor to be isotropic?
Definition: TensileStressUpdate.h:34
RankFourTensorTempl< Real >
MultiParameterPlasticityStressUpdate.h
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
TensorMechanicsHardeningModel.h
TensorMechanicsHardeningModel
Hardening Model base class.
Definition: TensorMechanicsHardeningModel.h:27