www.mooseframework.org
ComputeMultipleInelasticStress.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 
13 
14 #include "StressUpdateBase.h"
15 
17 class DamageBase;
18 
19 template <>
21 
39 {
40 public:
41  static InputParameters validParams();
42 
43  ComputeMultipleInelasticStress(const InputParameters & parameters);
44 
45  virtual void initialSetup() override;
46 
47 protected:
48  virtual void initQpStatefulProperties() override;
49 
50  virtual void computeQpStress() override;
51 
57 
64  virtual void finiteStrainRotation(const bool force_elasticity_rotation = false);
65 
78  virtual void updateQpState(RankTwoTensor & elastic_strain_increment,
79  RankTwoTensor & combined_inelastic_strain_increment);
80 
91  virtual void updateQpStateSingleModel(unsigned model_number,
92  RankTwoTensor & elastic_strain_increment,
93  RankTwoTensor & combined_inelastic_strain_increment);
94 
100  virtual void computeQpJacobianMult();
101 
116  virtual void computeAdmissibleState(unsigned model_number,
117  RankTwoTensor & elastic_strain_increment,
118  RankTwoTensor & inelastic_strain_increment,
119  RankFourTensor & consistent_tangent_operator);
120 
122  const unsigned int _max_iterations;
127 
130 
132  const MaterialProperty<RankTwoTensor> & _elastic_strain_old;
133  const MaterialProperty<RankTwoTensor> & _strain_increment;
135 
137  MaterialProperty<RankTwoTensor> & _inelastic_strain;
138 
140  const MaterialProperty<RankTwoTensor> & _inelastic_strain_old;
141 
144 
146  const unsigned _num_models;
147 
149  std::vector<bool> _tangent_computation_flag;
150 
153 
155  const std::vector<Real> _inelastic_weights;
156 
158  std::vector<RankFourTensor> _consistent_tangent_operator;
159 
161  const bool _cycle_models;
162 
163  MaterialProperty<Real> & _matl_timestep_limit;
164 
169 
177  std::vector<StressUpdateBase *> _models;
178 
181 
184 
187 
189 };
ComputeMultipleInelasticStress
ComputeMultipleInelasticStress computes the stress, the consistent tangent operator (or an approximat...
Definition: ComputeMultipleInelasticStress.h:38
ComputeMultipleInelasticStress::_inelastic_weights
const std::vector< Real > _inelastic_weights
_inelastic_strain = sum_i (_inelastic_weights_i * inelastic_strain_from_model_i)
Definition: ComputeMultipleInelasticStress.h:155
ComputeMultipleInelasticStress::_cycle_models
const bool _cycle_models
whether to cycle through the models, using only one model per timestep
Definition: ComputeMultipleInelasticStress.h:161
ComputeMultipleInelasticStress::_all_models_isotropic
bool _all_models_isotropic
are all inelastic models inherently isotropic? (not the case for e.g. weak plane plasticity models)
Definition: ComputeMultipleInelasticStress.h:183
ComputeFiniteStrainElasticStress
ComputeFiniteStrainElasticStress computes the stress following elasticity theory for finite strains.
Definition: ComputeFiniteStrainElasticStress.h:24
ComputeMultipleInelasticStress::validParams
static InputParameters validParams()
Definition: ComputeMultipleInelasticStress.C:21
ComputeMultipleInelasticStress::initialSetup
virtual void initialSetup() override
Definition: ComputeMultipleInelasticStress.C:122
ComputeMultipleInelasticStress::_num_models
const unsigned _num_models
number of plastic models
Definition: ComputeMultipleInelasticStress.h:146
validParams< ComputeMultipleInelasticStress >
InputParameters validParams< ComputeMultipleInelasticStress >()
ComputeMultipleInelasticStress::computeQpStress
virtual void computeQpStress() override
Compute the stress and store it in the _stress material property for the current quadrature point.
Definition: ComputeMultipleInelasticStress.C:187
ComputeMultipleInelasticStress::updateQpStateSingleModel
virtual void updateQpStateSingleModel(unsigned model_number, RankTwoTensor &elastic_strain_increment, RankTwoTensor &combined_inelastic_strain_increment)
An optimised version of updateQpState that gets used when the number of plastic models is unity,...
Definition: ComputeMultipleInelasticStress.C:415
ComputeMultipleInelasticStress::_inelastic_strain_old
const MaterialProperty< RankTwoTensor > & _inelastic_strain_old
old value of inelastic strain
Definition: ComputeMultipleInelasticStress.h:140
ComputeMultipleInelasticStress::_is_elasticity_tensor_guaranteed_isotropic
bool _is_elasticity_tensor_guaranteed_isotropic
is the elasticity tensor guaranteed to be isotropic?
Definition: ComputeMultipleInelasticStress.h:180
ComputeMultipleInelasticStress::_tangent_computation_flag
std::vector< bool > _tangent_computation_flag
Flags to compute tangent during updateState call.
Definition: ComputeMultipleInelasticStress.h:149
ComputeMultipleInelasticStress::_models
std::vector< StressUpdateBase * > _models
The user supplied list of inelastic models to use in the simulation.
Definition: ComputeMultipleInelasticStress.h:177
ComputeMultipleInelasticStress::_perform_finite_strain_rotations
const bool _perform_finite_strain_rotations
after updateQpState, rotate the stress, elastic_strain, inelastic_strain and Jacobian_mult using _rot...
Definition: ComputeMultipleInelasticStress.h:129
ComputeMultipleInelasticStress::TangentOperatorEnum::nonlinear
ComputeMultipleInelasticStress::TangentOperatorEnum::elastic
ComputeMultipleInelasticStress::_elastic_strain_old
const MaterialProperty< RankTwoTensor > & _elastic_strain_old
Strain tensors.
Definition: ComputeMultipleInelasticStress.h:132
ComputeMultipleInelasticStress::_undamaged_stress_old
RankTwoTensor _undamaged_stress_old
Definition: ComputeMultipleInelasticStress.h:188
ComputeMultipleInelasticStress::_internal_solve_full_iteration_history
const bool _internal_solve_full_iteration_history
Definition: ComputeMultipleInelasticStress.h:125
ComputeMultipleInelasticStress::_consistent_tangent_operator
std::vector< RankFourTensor > _consistent_tangent_operator
the consistent tangent operators computed by each plastic model
Definition: ComputeMultipleInelasticStress.h:158
ComputeMultipleInelasticStress::_inelastic_strain
MaterialProperty< RankTwoTensor > & _inelastic_strain
The sum of the inelastic strains that come from the plastic models.
Definition: ComputeMultipleInelasticStress.h:137
ComputeMultipleInelasticStress::_identity_symmetric_four
const RankFourTensor _identity_symmetric_four
Rank four symmetric identity tensor.
Definition: ComputeMultipleInelasticStress.h:168
ComputeMultipleInelasticStress::_strain_increment
const MaterialProperty< RankTwoTensor > & _strain_increment
Definition: ComputeMultipleInelasticStress.h:133
DamageBase
DamageBase is a base class for damage models, which modify the stress tensor computed by another mode...
Definition: DamageBase.h:28
ComputeMultipleInelasticStress::_damage_model
DamageBase * _damage_model
Pointer to the damage model.
Definition: ComputeMultipleInelasticStress.h:186
ComputeMultipleInelasticStress::_absolute_tolerance
const Real _absolute_tolerance
Definition: ComputeMultipleInelasticStress.h:124
ComputeMultipleInelasticStress::finiteStrainRotation
virtual void finiteStrainRotation(const bool force_elasticity_rotation=false)
Rotate _elastic_strain, _stress, _inelastic_strain, and _Jacobian_mult to the new configuration.
Definition: ComputeMultipleInelasticStress.C:255
ComputeMultipleInelasticStress::computeQpStressIntermediateConfiguration
virtual void computeQpStressIntermediateConfiguration()
Compute the stress for the current QP, but do not rotate tensors from the intermediate configuration ...
Definition: ComputeMultipleInelasticStress.C:215
ComputeMultipleInelasticStress::updateQpState
virtual void updateQpState(RankTwoTensor &elastic_strain_increment, RankTwoTensor &combined_inelastic_strain_increment)
Given the _strain_increment[_qp], iterate over all of the user-specified recompute materials in order...
Definition: ComputeMultipleInelasticStress.C:271
ComputeMultipleInelasticStress::initQpStatefulProperties
virtual void initQpStatefulProperties() override
Definition: ComputeMultipleInelasticStress.C:115
ComputeMultipleInelasticStress::computeAdmissibleState
virtual void computeAdmissibleState(unsigned model_number, RankTwoTensor &elastic_strain_increment, RankTwoTensor &inelastic_strain_increment, RankFourTensor &consistent_tangent_operator)
Given a trial stress (_stress[_qp]) and a strain increment (elastic_strain_increment) let the model_n...
Definition: ComputeMultipleInelasticStress.C:466
ComputeMultipleInelasticStress::_matl_timestep_limit
MaterialProperty< Real > & _matl_timestep_limit
Definition: ComputeMultipleInelasticStress.h:163
RankFourTensorTempl< Real >
ComputeMultipleInelasticStress::computeQpJacobianMult
virtual void computeQpJacobianMult()
Using _elasticity_tensor[_qp] and the consistent tangent operators, _consistent_tangent_operator[....
Definition: ComputeMultipleInelasticStress.C:393
ComputeMultipleInelasticStress::_tangent_calculation_method
TangentCalculationMethod _tangent_calculation_method
Calculation method for the tangent modulus.
Definition: ComputeMultipleInelasticStress.h:152
TangentCalculationMethod
TangentCalculationMethod
TangentCalculationMethod is an enum that determines the calculation method for the tangent operator.
Definition: StressUpdateBase.h:30
RankTwoTensorTempl< Real >
StressUpdateBase.h
ComputeMultipleInelasticStress::_tangent_operator_type
enum ComputeMultipleInelasticStress::TangentOperatorEnum _tangent_operator_type
ComputeMultipleInelasticStress::_max_iterations
const unsigned int _max_iterations
Input parameters associated with the recompute iteration to return the stress state to the yield surf...
Definition: ComputeMultipleInelasticStress.h:122
ComputeFiniteStrainElasticStress.h
ComputeMultipleInelasticStress::TangentOperatorEnum
TangentOperatorEnum
what sort of Tangent operator to calculate
Definition: ComputeMultipleInelasticStress.h:143
ComputeMultipleInelasticStress::ComputeMultipleInelasticStress
ComputeMultipleInelasticStress(const InputParameters &parameters)
Definition: ComputeMultipleInelasticStress.C:77
ComputeMultipleInelasticStress::_relative_tolerance
const Real _relative_tolerance
Definition: ComputeMultipleInelasticStress.h:123