https://mooseframework.inl.gov
ComputeMultipleCrystalPlasticityStress.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 
13 
16 
17 #include "RankTwoTensor.h"
18 #include "RankFourTensor.h"
19 
31 {
32 public:
34 
36 
37  virtual void initialSetup() override;
38 
39 protected:
40  virtual void computeQpStress() override;
41 
48  virtual void updateStress(RankTwoTensor & cauchy_stress, RankFourTensor & jacobian_mult);
49 
55  virtual void initQpStatefulProperties() override;
56 
62 
67  void preSolveQp();
68 
73  void solveQp();
74 
78  void postSolveQp(RankTwoTensor & stress_new, RankFourTensor & jacobian_mult);
79 
84  void solveStateVariables();
85 
89  void solveStress();
90 
97  void calculateResidual();
98 
104  void calculateJacobian();
105 
107  void calcTangentModuli(RankFourTensor & jacobian_mult);
108  void elasticTangentModuli(RankFourTensor & jacobian_mult);
109  void elastoPlasticTangentModuli(RankFourTensor & jacobian_mult);
111 
113  bool lineSearchUpdate(const Real & rnorm_prev, const RankTwoTensor & dpk2);
114 
119 
121  const unsigned _num_models;
122 
124  std::vector<CrystalPlasticityStressUpdateBase *> _models;
125 
127  const unsigned _num_eigenstrains;
128 
130  std::vector<ComputeCrystalPlasticityEigenstrainBase *> _eigenstrains;
131 
133  const std::string _base_name;
134 
137 
142 
147 
149  unsigned int _maxiter;
151  unsigned int _maxiterg;
152 
154  const enum class TangentModuliType { EXACT, NONE } _tan_mod_type;
155 
157  unsigned int _max_substep_iter;
158 
161 
164 
167 
170 
173 
176 
181 
186 
191 
196 
199 
205 
211 
219 
222 
225 
230 
233 };
const MaterialProperty< RankTwoTensor > & _pk2_old
void solveStateVariables()
Solves the internal variables stress as a function of the slip specified by the constitutive model de...
virtual void computeQpStress() override
Compute the stress and store it in the _stress material property for the current quadrature point...
MaterialProperty< RankTwoTensor > & _total_lagrangian_strain
Lagrangian total strain measure for the entire crystal.
void postSolveQp(RankTwoTensor &stress_new, RankFourTensor &jacobian_mult)
Save the final stress and internal variable values after the iterative solve.
RankTwoTensor _temporary_deformation_gradient
Helper deformation gradient tensor variables used in iterative solve.
Real _line_search_tolerance
Line search bisection method tolerance.
void calculateEigenstrainDeformationGrad()
Calculates the deformation gradient due to eigenstrain.
bool _convergence_failed
Flag to check whether convergence is achieved or if substepping is needed.
void calculateResidual()
Calculate stress residual as the difference between the stored material property PK2 stress and the e...
unsigned int _maxiter
Maximum number of iterations for stress update.
ComputeMultipleCrystalPlasticityStress(const InputParameters &parameters)
virtual void updateStress(RankTwoTensor &cauchy_stress, RankFourTensor &jacobian_mult)
Updates the stress (PK2) at a quadrature point by calling constiutive relationship as defined in a ch...
MaterialProperty< RankTwoTensor > * _eigenstrain_deformation_gradient
Generalized eigenstrain deformation gradient RankTwoTensor for the crystal.
const MaterialProperty< RankTwoTensor > & _crysrot
Crystal rotation in the original, or reference, configuration as defined by Euler angle arguments in ...
bool lineSearchUpdate(const Real &rnorm_prev, const RankTwoTensor &dpk2)
performs the line search update
virtual void initQpStatefulProperties() override
initializes the stateful properties such as PK2 stress, resolved shear stress, plastic deformation gr...
enum ComputeMultipleCrystalPlasticityStress::LineSearchMethod _line_search_method
void calculateJacobian()
Calculates the jacobian as ${J} = {I} - {C} {d{E}^e}{d{F}^e} {d{F}^e}{d{F}^P^{-1}} {d{F}^P^{-1}}{d{PK...
const MaterialProperty< RankTwoTensor > & _deformation_gradient
Total deformation gradient RankTwoTensor for the crystal.
void elastoPlasticTangentModuli(RankFourTensor &jacobian_mult)
Real _min_line_search_step_size
Minimum line search step size.
MaterialProperty< RankTwoTensor > & _pk2
Second Piola-Kirchoff stress measure.
Real _abs_tol
Stress residual equation absolute tolerance.
const MaterialProperty< RankTwoTensor > * _eigenstrain_deformation_gradient_old
ComputeFiniteStrainElasticStress computes the stress following elasticity theory for finite strains...
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor as defined by a separate class.
unsigned int _line_search_max_iterations
Line search bisection method maximum iteration number.
MaterialProperty< RankTwoTensor > & _updated_rotation
Tracks the rotation of the crystal during deformation Note: this rotation tensor is not applied to th...
unsigned int _max_substep_iter
Maximum number of substep iterations.
MaterialProperty< RankTwoTensor > & _plastic_deformation_gradient
Plastic deformation gradient RankTwoTensor for the crystal.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
enum ComputeMultipleCrystalPlasticityStress::TangentModuliType _tan_mod_type
std::vector< CrystalPlasticityStressUpdateBase * > _models
The user supplied cyrstal plasticity consititutive models.
unsigned int _maxiterg
Maximum number of iterations for internal variable update.
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
const MaterialProperty< RankTwoTensor > & _plastic_deformation_gradient_old
Real _dfgrd_scale_factor
Scales the substepping increment to obtain deformation gradient at a substep iteration.
std::vector< ComputeCrystalPlasticityEigenstrainBase * > _eigenstrains
The user supplied cyrstal plasticity eigenstrains.
const InputParameters & parameters() const
void calculateResidualAndJacobian()
Calls the residual and jacobian functions used in the stress update algorithm.
void solveQp()
Solve the stress and internal state variables (e.g.
ComputeMultipleCrystalPlasticityStress (used together with CrystalPlasticityStressUpdateBase) uses th...
void solveStress()
solves for stress, updates plastic deformation gradient.
const bool _print_convergence_message
Flag to print to console warning messages on stress, constitutive model convergence.
void calcTangentModuli(RankFourTensor &jacobian_mult)
Calculates the tangent moduli for use as a preconditioner, using the elastic or elastic-plastic optio...
RankTwoTensor _delta_deformation_gradient
Used for substepping; Uniformly divides the increment in deformation gradient.
Real _rtol
Stress residual equation relative tolerance.
const std::string _base_name
optional parameter to define several mechanical systems on the same block, e.g. multiple phases ...
void preSolveQp()
Reset the PK2 stress and the inverse deformation gradient to old values and provide an interface for ...