www.mooseframework.org
HyperElasticPhaseFieldIsoDamage.C
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 
11 #include "libmesh/utility.h"
12 
14 
16 
17 InputParameters
19 {
20  InputParameters params = FiniteStrainHyperElasticViscoPlastic::validParams();
21  params.addParam<bool>("numerical_stiffness", false, "Flag for numerical stiffness");
22  params.addParam<Real>("damage_stiffness", 1e-8, "Avoid zero after complete damage");
23  params.addParam<Real>("zero_tol", 1e-12, "Tolerance for numerical zero");
24  params.addParam<Real>(
25  "zero_perturb", 1e-8, "Perturbation value when strain value less than numerical zero");
26  params.addParam<Real>("perturbation_scale_factor", 1e-5, "Perturbation scale factor");
27  params.addRequiredCoupledVar("c", "Damage variable");
28  params.addParam<bool>(
29  "use_current_history_variable", false, "Use the current value of the history variable.");
30  params.addParam<MaterialPropertyName>(
31  "F_name", "E_el", "Name of material property storing the elastic energy");
32  params.addClassDescription(
33  "Computes damaged stress and energy in the intermediate configuration assuming isotropy");
34 
35  return params;
36 }
37 
40  _num_stiffness(getParam<bool>("numerical_stiffness")),
41  _kdamage(getParam<Real>("damage_stiffness")),
42  _use_current_hist(getParam<bool>("use_current_history_variable")),
43  _l(getMaterialProperty<Real>("l")),
44  _gc(getMaterialProperty<Real>("gc_prop")),
45  _zero_tol(getParam<Real>("zero_tol")),
46  _zero_pert(getParam<Real>("zero_perturb")),
47  _pert_val(getParam<Real>("perturbation_scale_factor")),
48  _c(coupledValue("c")),
49  _save_state(false),
50  _dstress_dc(
51  declarePropertyDerivative<RankTwoTensor>(_base_name + "stress", getVar("c", 0)->name())),
52  _etens(LIBMESH_DIM),
53  _F(declareProperty<Real>(getParam<MaterialPropertyName>("F_name"))),
54  _dFdc(declarePropertyDerivative<Real>(getParam<MaterialPropertyName>("F_name"),
55  getVar("c", 0)->name())),
56  _d2Fdc2(declarePropertyDerivative<Real>(
57  getParam<MaterialPropertyName>("F_name"), getVar("c", 0)->name(), getVar("c", 0)->name())),
58  _d2Fdcdstrain(declareProperty<RankTwoTensor>("d2Fdcdstrain")),
59  _hist(declareProperty<Real>("hist")),
60  _hist_old(getMaterialPropertyOld<Real>("hist"))
61 {
62 }
63 
64 void
66 {
68 
69  _save_state = true;
71  _pk2[_qp] = _pk2_tmp;
72 
73  _save_state = false;
74  if (_num_stiffness)
76 
77  if (_num_stiffness)
79 
80  _dce_dfe.zero();
81  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
82  for (unsigned int j = 0; j < LIBMESH_DIM; ++j)
83  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
84  {
85  _dce_dfe(i, j, k, i) = _dce_dfe(i, j, k, i) + _fe(k, j);
86  _dce_dfe(i, j, k, j) = _dce_dfe(i, j, k, j) + _fe(k, i);
87  }
88 
90 }
91 
92 void
94 {
95  Real lambda = _elasticity_tensor[_qp](0, 0, 1, 1);
96  Real mu = _elasticity_tensor[_qp](0, 1, 0, 1);
97 
98  Real c = _c[_qp];
99  Real xfac = Utility::pow<2>(1.0 - c) + _kdamage;
100 
101  std::vector<Real> eigval;
102  RankTwoTensor evec;
103  _ee.symmetricEigenvaluesEigenvectors(eigval, evec);
104 
105  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
106  _etens[i].vectorOuterProduct(evec.column(i), evec.column(i));
107 
108  Real etr = 0.0;
109  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
110  etr += eigval[i];
111 
112  Real etrpos = (std::abs(etr) + etr) / 2.0;
113  Real etrneg = (std::abs(etr) - etr) / 2.0;
114 
115  RankTwoTensor pk2pos, pk2neg;
116 
117  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
118  {
119  pk2pos += _etens[i] * (lambda * etrpos + 2.0 * mu * (std::abs(eigval[i]) + eigval[i]) / 2.0);
120  pk2neg += _etens[i] * (lambda * etrneg + 2.0 * mu * (std::abs(eigval[i]) - eigval[i]) / 2.0);
121  }
122 
123  _pk2_tmp = pk2pos * xfac - pk2neg;
124 
125  if (_save_state)
126  {
127  std::vector<Real> epos(LIBMESH_DIM), eneg(LIBMESH_DIM);
128  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
129  {
130  epos[i] = (std::abs(eigval[i]) + eigval[i]) / 2.0;
131  eneg[i] = (std::abs(eigval[i]) - eigval[i]) / 2.0;
132  }
133 
134  // sum squares of epos and eneg
135  Real pval(0.0), nval(0.0);
136  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
137  {
138  pval += epos[i] * epos[i];
139  nval += eneg[i] * eneg[i];
140  }
141 
142  // Energy with positive principal strains
143  const Real G0_pos = lambda * etrpos * etrpos / 2.0 + mu * pval;
144  const Real G0_neg = lambda * etrneg * etrneg / 2.0 + mu * nval;
145 
146  // Assign history variable and derivative
147  if (G0_pos > _hist_old[_qp])
148  _hist[_qp] = G0_pos;
149  else
150  _hist[_qp] = _hist_old[_qp];
151 
152  Real hist_variable = _hist_old[_qp];
153  if (_use_current_hist)
154  hist_variable = _hist[_qp];
155 
156  // Elastic free energy density
157  _F[_qp] = hist_variable * xfac - G0_neg + _gc[_qp] / (2 * _l[_qp]) * c * c;
158 
159  // derivative of elastic free energy density wrt c
160  _dFdc[_qp] = -hist_variable * 2.0 * (1.0 - c) * (1 - _kdamage) + _gc[_qp] / _l[_qp] * c;
161 
162  // 2nd derivative of elastic free energy density wrt c
163  _d2Fdc2[_qp] = hist_variable * 2.0 * (1 - _kdamage) + _gc[_qp] / _l[_qp];
164 
165  _dG0_dee = pk2pos;
166 
167  _dpk2_dc = -pk2pos * 2.0 * (1.0 - c);
168  }
169 }
170 
171 void
173 {
174  RankTwoTensor ee_tmp;
175 
176  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
177  for (unsigned int j = i; j < LIBMESH_DIM; ++j)
178  {
179  Real ee_pert = _zero_pert;
180  if (std::abs(_ee(i, j)) > _zero_tol)
181  ee_pert = _pert_val * std::abs(_ee(i, j));
182 
183  ee_tmp = _ee;
184  _ee(i, j) += ee_pert;
185 
187 
188  for (unsigned int k = 0; k < LIBMESH_DIM; ++k)
189  for (unsigned int l = 0; l < LIBMESH_DIM; ++l)
190  {
191  _dpk2_dee(k, l, i, j) = (_pk2_tmp(k, l) - _pk2[_qp](k, l)) / ee_pert;
192  _dpk2_dee(k, l, j, i) = (_pk2_tmp(k, l) - _pk2[_qp](k, l)) / ee_pert;
193  }
194  _ee = ee_tmp;
195  }
196 }
197 
198 void
200 {
202 
203  RankTwoTensor dG0_dce = _dee_dce.innerProductTranspose(_dG0_dee);
204  RankTwoTensor dG0_dfe = _dce_dfe.innerProductTranspose(dG0_dce);
205  RankTwoTensor dG0_df = _dfe_df.innerProductTranspose(dG0_dfe);
206 
207  // 2nd derivative wrt c and strain = 0.0 if we used the previous step's history varible
208  if (_use_current_hist)
209  _d2Fdcdstrain[_qp] =
210  -_df_dstretch_inc.innerProductTranspose(dG0_df) * 2.0 * (1.0 - _c[_qp]) * (1 - _kdamage);
211 
212  _dstress_dc[_qp] = _fe.mixedProductIkJl(_fe) * _dpk2_dc;
213 }
HyperElasticPhaseFieldIsoDamage::_num_stiffness
bool _num_stiffness
Flag to compute numerical stiffness.
Definition: HyperElasticPhaseFieldIsoDamage.h:47
HyperElasticPhaseFieldIsoDamage::_hist
MaterialProperty< Real > & _hist
History variable that prevents crack healing, declared in this material.
Definition: HyperElasticPhaseFieldIsoDamage.h:80
HyperElasticPhaseFieldIsoDamage::validParams
static InputParameters validParams()
Definition: HyperElasticPhaseFieldIsoDamage.C:18
FiniteStrainHyperElasticViscoPlastic::_dfe_df
RankFourTensor _dfe_df
Definition: FiniteStrainHyperElasticViscoPlastic.h:221
FiniteStrainHyperElasticViscoPlastic::computeElasticStrain
virtual void computeElasticStrain()
Computes elastic Lagrangian strain.
Definition: FiniteStrainHyperElasticViscoPlastic.C:482
FiniteStrainHyperElasticViscoPlastic::_df_dstretch_inc
RankFourTensor _df_dstretch_inc
Definition: FiniteStrainHyperElasticViscoPlastic.h:222
FiniteStrainHyperElasticViscoPlastic::computeQpJacobian
virtual void computeQpJacobian()
This function computes the Jacobian.
Definition: FiniteStrainHyperElasticViscoPlastic.C:549
HyperElasticPhaseFieldIsoDamage::computeQpJacobian
virtual void computeQpJacobian()
This function computes tensors used to construct diagonal and off-diagonal Jacobian.
Definition: HyperElasticPhaseFieldIsoDamage.C:199
HyperElasticPhaseFieldIsoDamage::_dG0_dee
RankTwoTensor _dG0_dee
Definition: HyperElasticPhaseFieldIsoDamage.h:69
registerMooseObject
registerMooseObject("TensorMechanicsApp", HyperElasticPhaseFieldIsoDamage)
HyperElasticPhaseFieldIsoDamage::_pk2_tmp
RankTwoTensor _pk2_tmp
Definition: HyperElasticPhaseFieldIsoDamage.h:69
FiniteStrainHyperElasticViscoPlastic::_dpk2_dce
RankFourTensor _dpk2_dce
Definition: FiniteStrainHyperElasticViscoPlastic.h:220
FiniteStrainHyperElasticViscoPlastic::validParams
static InputParameters validParams()
Definition: FiniteStrainHyperElasticViscoPlastic.C:18
HyperElasticPhaseFieldIsoDamage::computePK2StressAndDerivative
virtual void computePK2StressAndDerivative()
This function computes PK2 stress.
Definition: HyperElasticPhaseFieldIsoDamage.C:65
HyperElasticPhaseFieldIsoDamage::_hist_old
const MaterialProperty< Real > & _hist_old
Old value of history variable.
Definition: HyperElasticPhaseFieldIsoDamage.h:83
FiniteStrainHyperElasticViscoPlastic::_elasticity_tensor
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
Definition: FiniteStrainHyperElasticViscoPlastic.h:204
HyperElasticPhaseFieldIsoDamage::_use_current_hist
bool _use_current_hist
Use current value of history variable.
Definition: HyperElasticPhaseFieldIsoDamage.h:51
HyperElasticPhaseFieldIsoDamage::HyperElasticPhaseFieldIsoDamage
HyperElasticPhaseFieldIsoDamage(const InputParameters &parameters)
Definition: HyperElasticPhaseFieldIsoDamage.C:38
HyperElasticPhaseFieldIsoDamage::_save_state
bool _save_state
Flag to save couple material properties.
Definition: HyperElasticPhaseFieldIsoDamage.h:65
HyperElasticPhaseFieldIsoDamage::_dstress_dc
MaterialProperty< RankTwoTensor > & _dstress_dc
Definition: HyperElasticPhaseFieldIsoDamage.h:67
FiniteStrainHyperElasticViscoPlastic::_dee_dce
RankFourTensor _dee_dce
Definition: FiniteStrainHyperElasticViscoPlastic.h:221
defineLegacyParams
defineLegacyParams(HyperElasticPhaseFieldIsoDamage)
HyperElasticPhaseFieldIsoDamage::_zero_pert
Real _zero_pert
Perturbation value for near zero or zero strain components.
Definition: HyperElasticPhaseFieldIsoDamage.h:59
HyperElasticPhaseFieldIsoDamage::_d2Fdcdstrain
MaterialProperty< RankTwoTensor > & _d2Fdcdstrain
Definition: HyperElasticPhaseFieldIsoDamage.h:77
HyperElasticPhaseFieldIsoDamage::computeNumStiffness
virtual void computeNumStiffness()
This function computes numerical stiffness.
Definition: HyperElasticPhaseFieldIsoDamage.C:172
FiniteStrainHyperElasticViscoPlastic::_dce_dfe
RankFourTensor _dce_dfe
Definition: FiniteStrainHyperElasticViscoPlastic.h:221
name
const std::string name
Definition: Setup.h:21
HyperElasticPhaseFieldIsoDamage::_dpk2_dc
RankTwoTensor _dpk2_dc
Definition: HyperElasticPhaseFieldIsoDamage.h:69
FiniteStrainHyperElasticViscoPlastic::_dpk2_dfe
RankFourTensor _dpk2_dfe
Definition: FiniteStrainHyperElasticViscoPlastic.h:220
HyperElasticPhaseFieldIsoDamage::_dFdc
MaterialProperty< Real > & _dFdc
Definition: HyperElasticPhaseFieldIsoDamage.h:75
HyperElasticPhaseFieldIsoDamage::_pert_val
Real _pert_val
Perturbation value for strain components.
Definition: HyperElasticPhaseFieldIsoDamage.h:61
HyperElasticPhaseFieldIsoDamage::_etens
std::vector< RankTwoTensor > _etens
Definition: HyperElasticPhaseFieldIsoDamage.h:71
FiniteStrainHyperElasticViscoPlastic::_ee
RankTwoTensor _ee
Definition: FiniteStrainHyperElasticViscoPlastic.h:218
HyperElasticPhaseFieldIsoDamage::_l
const MaterialProperty< Real > & _l
Material property defining crack width, declared elsewhere.
Definition: HyperElasticPhaseFieldIsoDamage.h:53
FiniteStrainHyperElasticViscoPlastic::_fe
RankTwoTensor _fe
Definition: FiniteStrainHyperElasticViscoPlastic.h:218
HyperElasticPhaseFieldIsoDamage::_gc
const MaterialProperty< Real > & _gc
Material property defining gc parameter, declared elsewhere.
Definition: HyperElasticPhaseFieldIsoDamage.h:55
HyperElasticPhaseFieldIsoDamage::_F
MaterialProperty< Real > & _F
Elastic energy and derivatives, declared in this material.
Definition: HyperElasticPhaseFieldIsoDamage.h:74
RankTwoTensorTempl< Real >
HyperElasticPhaseFieldIsoDamage::_dpk2_dee
RankFourTensor _dpk2_dee
Definition: HyperElasticPhaseFieldIsoDamage.h:70
HyperElasticPhaseFieldIsoDamage::_c
const VariableValue & _c
Compupled damage variable.
Definition: HyperElasticPhaseFieldIsoDamage.h:63
FiniteStrainHyperElasticViscoPlastic::_pk2
MaterialProperty< RankTwoTensor > & _pk2
Definition: FiniteStrainHyperElasticViscoPlastic.h:196
HyperElasticPhaseFieldIsoDamage
This class solves visco plastic model based on isotropically damaged stress The damage parameter is o...
Definition: HyperElasticPhaseFieldIsoDamage.h:25
HyperElasticPhaseFieldIsoDamage::computeDamageStress
virtual void computeDamageStress()
This function computes PK2 stress modified to account for damage Computes numerical stiffness if flag...
Definition: HyperElasticPhaseFieldIsoDamage.C:93
FiniteStrainHyperElasticViscoPlastic
This class solves the viscoplastic flow rate equations in the total form Involves 4 different types o...
Definition: FiniteStrainHyperElasticViscoPlastic.h:33
HyperElasticPhaseFieldIsoDamage::_zero_tol
Real _zero_tol
Used in numerical stiffness calculation to check near zero values.
Definition: HyperElasticPhaseFieldIsoDamage.h:57
HyperElasticPhaseFieldIsoDamage::_kdamage
Real _kdamage
Small stiffness of completely damaged material point.
Definition: HyperElasticPhaseFieldIsoDamage.h:49
HyperElasticPhaseFieldIsoDamage.h
HyperElasticPhaseFieldIsoDamage::_d2Fdc2
MaterialProperty< Real > & _d2Fdc2
Definition: HyperElasticPhaseFieldIsoDamage.h:76