www.mooseframework.org
RateDepSmearCrackModel.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 "ConstitutiveModel.h"
13 #include "SymmElasticityTensor.h"
14 
15 #include "petscsys.h"
16 #include "petscblaslapack.h"
17 
18 #if PETSC_VERSION_LESS_THAN(3, 5, 0)
19 extern "C" void FORTRAN_CALL(dgetri)(...); // matrix inversion routine from LAPACK
20 #endif
21 
23 
24 template <>
25 InputParameters validParams<RateDepSmearCrackModel>();
26 
33 {
34 public:
35  RateDepSmearCrackModel(const InputParameters & parameters);
36 
37  virtual ~RateDepSmearCrackModel();
38 
39 protected:
40  virtual void computeStress(const Elem & current_elem,
41  const SymmElasticityTensor & elasticity_tensor,
42  const SymmTensor & stress_old,
43  SymmTensor & strain_increment,
44  SymmTensor & stress_new);
45 
46  virtual void initQpStatefulProperties();
47 
48  virtual void initVariables();
49 
57  virtual void solve();
61  virtual void postSolveVariables();
66  virtual void postSolveStress();
67 
72  virtual void updateVariables();
73 
77  bool getConvergeVar();
78 
83  virtual void calcResidual();
84 
89  virtual void calcStateIncr();
90 
94  virtual void calcJacobian();
95 
96  int matrixInversion(std::vector<Real> & A, int n) const;
97 
99  unsigned int _nstate;
100  Real _exponent;
101  unsigned int _maxiter;
102  Real _tol;
103  Real _zero_tol;
107 
108  MaterialProperty<std::vector<Real>> & _intvar;
109  const MaterialProperty<std::vector<Real>> & _intvar_old;
110 
111  MaterialProperty<SymmTensor> & _stress_undamaged;
112  const MaterialProperty<SymmTensor> & _stress_undamaged_old;
113 
114  std::vector<Real> _intvar_incr;
115  std::vector<Real> _intvar_tmp, _intvar_old_tmp;
116  std::vector<Real> _resid;
117  std::vector<Real> _jac;
118  std::vector<Real> _dvar;
119 
123  bool _nconv;
124  bool _err_tol;
125 
126 private:
127 };
128 
dgetri
void FORTRAN_CALL() dgetri(...)
SymmElasticityTensor.h
RateDepSmearCrackModel::_intvar_incr_tol
Real _intvar_incr_tol
Tolerance for zero.
Definition: RateDepSmearCrackModel.h:104
RateDepSmearCrackModel::postSolveStress
virtual void postSolveStress()
This function updates the stress after solve.
Definition: RateDepSmearCrackModel.C:194
RateDepSmearCrackModel::_stress_old
SymmTensor _stress_old
Definition: RateDepSmearCrackModel.h:121
RateDepSmearCrackModel::_stress_new
SymmTensor _stress_new
Definition: RateDepSmearCrackModel.h:121
RateDepSmearCrackModel::_elasticity
SymmElasticityTensor _elasticity
Definition: RateDepSmearCrackModel.h:120
ConstitutiveModel.h
RateDepSmearCrackModel::initVariables
virtual void initVariables()
Definition: RateDepSmearCrackModel.C:108
RateDepSmearCrackModel
RateDepSmearCrackModel is the base class for rate dependent continuum damage model.
Definition: RateDepSmearCrackModel.h:32
ConstitutiveModel
Definition: ConstitutiveModel.h:22
RateDepSmearCrackModel::_intvar_tmp
std::vector< Real > _intvar_tmp
Definition: RateDepSmearCrackModel.h:115
RateDepSmearCrackModel::_maxiter
unsigned int _maxiter
Definition: RateDepSmearCrackModel.h:101
RateDepSmearCrackModel::postSolveVariables
virtual void postSolveVariables()
This function updates the internal variables after solve.
Definition: RateDepSmearCrackModel.C:187
RateDepSmearCrackModel::_zero_tol
Real _zero_tol
Relative tolerance factor for convergence of the Newton Raphson solve.
Definition: RateDepSmearCrackModel.h:103
RateDepSmearCrackModel::_dstrain
SymmTensor _dstrain
Definition: RateDepSmearCrackModel.h:121
RateDepSmearCrackModel::_stress_undamaged_old
const MaterialProperty< SymmTensor > & _stress_undamaged_old
Definition: RateDepSmearCrackModel.h:112
RateDepSmearCrackModel::_err_tol
bool _err_tol
Convergence flag.
Definition: RateDepSmearCrackModel.h:124
RateDepSmearCrackModel::_intvar_incr
std::vector< Real > _intvar_incr
Definition: RateDepSmearCrackModel.h:114
RateDepSmearCrackModel::RateDepSmearCrackModel
RateDepSmearCrackModel(const InputParameters &parameters)
Definition: RateDepSmearCrackModel.C:41
RateDepSmearCrackModel::_jac
std::vector< Real > _jac
Definition: RateDepSmearCrackModel.h:117
RateDepSmearCrackModel::updateVariables
virtual void updateVariables()
This function updates variables during solve a(i+1) = a(i) + da(i+1)
Definition: RateDepSmearCrackModel.C:145
RateDepSmearCrackModel::_ref_damage_rate
Real _ref_damage_rate
Definition: RateDepSmearCrackModel.h:98
SymmElasticityTensor
This class defines a basic set of capabilities any elasticity tensor should have.
Definition: SymmElasticityTensor.h:55
RateDepSmearCrackModel::_stress_undamaged
MaterialProperty< SymmTensor > & _stress_undamaged
Definition: RateDepSmearCrackModel.h:111
RateDepSmearCrackModel::computeStress
virtual void computeStress(const Elem &current_elem, const SymmElasticityTensor &elasticity_tensor, const SymmTensor &stress_old, SymmTensor &strain_increment, SymmTensor &stress_new)
Definition: RateDepSmearCrackModel.C:76
RateDepSmearCrackModel::_rndm_scale_var
Real _rndm_scale_var
Flag to specify scaling parameter to generate random stress.
Definition: RateDepSmearCrackModel.h:106
validParams< RateDepSmearCrackModel >
InputParameters validParams< RateDepSmearCrackModel >()
Definition: RateDepSmearCrackModel.C:16
RateDepSmearCrackModel::~RateDepSmearCrackModel
virtual ~RateDepSmearCrackModel()
Definition: RateDepSmearCrackModel.C:246
RateDepSmearCrackModel::_nstate
unsigned int _nstate
reference damage rate
Definition: RateDepSmearCrackModel.h:99
RateDepSmearCrackModel::solve
virtual void solve()
This function solves the state variables.
Definition: RateDepSmearCrackModel.C:122
RateDepSmearCrackModel::_dvar
std::vector< Real > _dvar
Definition: RateDepSmearCrackModel.h:118
RateDepSmearCrackModel::_resid
std::vector< Real > _resid
Definition: RateDepSmearCrackModel.h:116
RateDepSmearCrackModel::_stress0
SymmTensor _stress0
Definition: RateDepSmearCrackModel.h:122
RateDepSmearCrackModel::_nconv
bool _nconv
Definition: RateDepSmearCrackModel.h:123
RateDepSmearCrackModel::initQpStatefulProperties
virtual void initQpStatefulProperties()
Definition: RateDepSmearCrackModel.C:69
RateDepSmearCrackModel::calcJacobian
virtual void calcJacobian()
This function calculates the Jacobian.
Definition: RateDepSmearCrackModel.C:216
RateDepSmearCrackModel::_intvar_old_tmp
std::vector< Real > _intvar_old_tmp
Definition: RateDepSmearCrackModel.h:115
RateDepSmearCrackModel::_intvar
MaterialProperty< std::vector< Real > > & _intvar
Variable value.
Definition: RateDepSmearCrackModel.h:108
RateDepSmearCrackModel::calcResidual
virtual void calcResidual()
This function calculates the residual as r = v - v_old - dv.
Definition: RateDepSmearCrackModel.C:199
RateDepSmearCrackModel::_exponent
Real _exponent
Number of state variables.
Definition: RateDepSmearCrackModel.h:100
RateDepSmearCrackModel::matrixInversion
int matrixInversion(std::vector< Real > &A, int n) const
Definition: RateDepSmearCrackModel.C:221
RateDepSmearCrackModel::_dstress0
SymmTensor _dstress0
Definition: RateDepSmearCrackModel.h:122
SymmTensor
Definition: SymmTensor.h:21
RateDepSmearCrackModel::_input_rndm_scale_var
bool _input_rndm_scale_var
Allowable relative increment size of state variables (dv)
Definition: RateDepSmearCrackModel.h:105
RateDepSmearCrackModel::getConvergeVar
bool getConvergeVar()
This function returns true if convergence is not achieved.
Definition: RateDepSmearCrackModel.C:163
RateDepSmearCrackModel::_tol
Real _tol
Maximum number of Newton Raphson iteration.
Definition: RateDepSmearCrackModel.h:102
RateDepSmearCrackModel::_intvar_old
const MaterialProperty< std::vector< Real > > & _intvar_old
Definition: RateDepSmearCrackModel.h:109
RateDepSmearCrackModel::calcStateIncr
virtual void calcStateIncr()
This function calculated thw increment of the state variables (dv) used to form the residual.
Definition: RateDepSmearCrackModel.C:211