www.mooseframework.org
ComputeLinearElasticPFFractureStress.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 "MathUtils.h"
12 
14 
16 
17 InputParameters
19 {
20  InputParameters params = ComputePFFractureStressBase::validParams();
21  params.addClassDescription("Computes the stress and free energy derivatives for the phase field "
22  "fracture model, with small strain");
23  MooseEnum Decomposition("strain_spectral strain_vol_dev stress_spectral none", "none");
24  params.addParam<MooseEnum>("decomposition_type",
25  Decomposition,
26  "Decomposition approaches. Choices are: " +
27  Decomposition.getRawNames());
28  return params;
29 }
30 
32  const InputParameters & parameters)
33  : ComputePFFractureStressBase(parameters),
34  GuaranteeConsumer(this),
35  _decomposition_type(getParam<MooseEnum>("decomposition_type").getEnum<Decomposition_type>())
36 {
37 }
38 
39 void
41 {
45  mooseError("Decomposition approach of strain_vol_dev and strain_spectral can only be used with "
46  "isotropic elasticity tensor materials, use stress_spectral for anistropic "
47  "elasticity tensor materials");
48 }
49 
50 void
52 {
53  // Isotropic elasticity is assumed and should be enforced
54  const Real lambda = _elasticity_tensor[_qp](0, 0, 1, 1);
55  const Real mu = _elasticity_tensor[_qp](0, 1, 0, 1);
56 
57  // Compute eigenvectors and eigenvalues of mechanical strain and projection tensor
58  RankTwoTensor eigvec;
59  std::vector<Real> eigval(LIBMESH_DIM);
60  RankFourTensor Ppos =
61  _mechanical_strain[_qp].positiveProjectionEigenDecomposition(eigval, eigvec);
62  RankFourTensor I4sym(RankFourTensor::initIdentitySymmetricFour);
63 
64  // Calculate tensors of outerproduct of eigen vectors
65  std::vector<RankTwoTensor> etens(LIBMESH_DIM);
66 
67  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
68  etens[i].vectorOuterProduct(eigvec.column(i), eigvec.column(i));
69 
70  // Separate out positive and negative eigen values
71  std::vector<Real> epos(LIBMESH_DIM), eneg(LIBMESH_DIM);
72  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
73  {
74  epos[i] = (std::abs(eigval[i]) + eigval[i]) / 2.0;
75  eneg[i] = -(std::abs(eigval[i]) - eigval[i]) / 2.0;
76  }
77 
78  // Seprate positive and negative sums of all eigenvalues
79  Real etr = 0.0;
80  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
81  etr += eigval[i];
82 
83  const Real etrpos = (std::abs(etr) + etr) / 2.0;
84  const Real etrneg = -(std::abs(etr) - etr) / 2.0;
85 
86  // Calculate the tensile (postive) and compressive (negative) parts of stress
87  RankTwoTensor stress0pos, stress0neg;
88  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
89  {
90  stress0pos += etens[i] * (lambda * etrpos + 2.0 * mu * epos[i]);
91  stress0neg += etens[i] * (lambda * etrneg + 2.0 * mu * eneg[i]);
92  }
93 
94  // sum squares of epos and eneg
95  Real pval(0.0), nval(0.0);
96  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
97  {
98  pval += epos[i] * epos[i];
99  nval += eneg[i] * eneg[i];
100  }
101 
102  _stress[_qp] = stress0pos * _D[_qp] + stress0neg;
103 
104  // Energy with positive principal strains
105  F_pos = lambda * etrpos * etrpos / 2.0 + mu * pval;
106  F_neg = -lambda * etrneg * etrneg / 2.0 + mu * nval;
107 
108  // 2nd derivative wrt c and strain = 0.0 if we used the previous step's history varible
109  if (_use_current_hist)
110  _d2Fdcdstrain[_qp] = stress0pos * _dDdc[_qp];
111 
112  // Used in StressDivergencePFFracTensors off-diagonal Jacobian
113  _dstress_dc[_qp] = stress0pos * _dDdc[_qp];
114 
115  _Jacobian_mult[_qp] = (I4sym - (1 - _D[_qp]) * Ppos) * _elasticity_tensor[_qp];
116 }
117 
118 void
120 {
121  // Compute Uncracked stress
123 
124  // Create the positive and negative projection tensors
125  RankFourTensor I4sym(RankFourTensor::initIdentitySymmetricFour);
126  std::vector<Real> eigval;
127  RankTwoTensor eigvec;
128  RankFourTensor Ppos = stress.positiveProjectionEigenDecomposition(eigval, eigvec);
129 
130  // Project the positive and negative stresses
131  RankTwoTensor stress0pos = Ppos * stress;
132  RankTwoTensor stress0neg = stress - stress0pos;
133 
134  // Compute the positive and negative elastic energies
135  F_pos = (stress0pos).doubleContraction(_mechanical_strain[_qp]) / 2.0;
136  F_neg = (stress0neg).doubleContraction(_mechanical_strain[_qp]) / 2.0;
137 
138  _stress[_qp] = stress0pos * _D[_qp] + stress0neg;
139 
140  // 2nd derivative wrt c and strain = 0.0 if we used the previous step's history varible
141  if (_use_current_hist)
142  _d2Fdcdstrain[_qp] = stress0pos * _dDdc[_qp];
143 
144  // Used in StressDivergencePFFracTensors off-diagonal Jacobian
145  _dstress_dc[_qp] = stress0pos * _dDdc[_qp];
146 
147  _Jacobian_mult[_qp] = (I4sym - (1 - _D[_qp]) * Ppos) * _elasticity_tensor[_qp];
148 }
149 
150 void
152 {
153  // Isotropic elasticity is assumed and should be enforced
154  const Real lambda = _elasticity_tensor[_qp](0, 0, 1, 1);
155  const Real mu = _elasticity_tensor[_qp](0, 1, 0, 1);
156  const Real k = lambda + 2.0 * mu / LIBMESH_DIM;
157 
158  RankTwoTensor I2(RankTwoTensor::initIdentity);
159  RankFourTensor I2I2 = I2.outerProduct(I2);
160 
161  RankFourTensor Jacobian_pos, Jacobian_neg;
162  RankTwoTensor strain0vol, strain0dev;
163  RankTwoTensor stress0pos, stress0neg;
164  Real strain0tr, strain0tr_neg, strain0tr_pos;
165 
166  strain0dev = _mechanical_strain[_qp].deviatoric();
167  strain0vol = _mechanical_strain[_qp] - strain0dev;
168  strain0tr = _mechanical_strain[_qp].trace();
169  strain0tr_neg = std::min(strain0tr, 0.0);
170  strain0tr_pos = strain0tr - strain0tr_neg;
171  stress0neg = k * strain0tr_neg * I2;
172  stress0pos = _elasticity_tensor[_qp] * _mechanical_strain[_qp] - stress0neg;
173  // Energy with positive principal strains
174  RankTwoTensor strain0dev2 = strain0dev * strain0dev;
175  F_pos = 0.5 * k * strain0tr_pos * strain0tr_pos + mu * strain0dev2.trace();
176  F_neg = 0.5 * k * strain0tr_neg * strain0tr_neg;
177 
178  _stress[_qp] = stress0pos * _D[_qp] + stress0neg;
179 
180  // 2nd derivative wrt c and strain = 0.0 if we used the previous step's history varible
181  if (_use_current_hist)
182  _d2Fdcdstrain[_qp] = stress0pos * _dDdc[_qp];
183 
184  // Used in StressDivergencePFFracTensors off-diagonal Jacobian
185  _dstress_dc[_qp] = stress0pos * _dDdc[_qp];
186 
187  if (strain0tr < 0)
188  Jacobian_neg = k * I2I2;
189  Jacobian_pos = _elasticity_tensor[_qp] - Jacobian_neg;
190  _Jacobian_mult[_qp] = _D[_qp] * Jacobian_pos + Jacobian_neg;
191 }
192 
193 void
195 {
196  Real F_pos, F_neg;
197 
198  switch (_decomposition_type)
199  {
201  computeStrainSpectral(F_pos, F_neg);
202  break;
204  computeStrainVolDev(F_pos, F_neg);
205  break;
207  computeStressSpectral(F_pos, F_neg);
208  break;
209  default:
210  {
211  _stress[_qp] = _D[_qp] * _elasticity_tensor[_qp] * _mechanical_strain[_qp];
212  F_pos = (_stress[_qp]).doubleContraction(_mechanical_strain[_qp]) / 2.0;
213  F_neg = 0.0;
214  if (_use_current_hist)
215  _d2Fdcdstrain[_qp] = _stress[_qp] * _dDdc[_qp];
216 
217  _dstress_dc[_qp] = _stress[_qp] * _dDdc[_qp];
218  _Jacobian_mult[_qp] = _D[_qp] * _elasticity_tensor[_qp];
219  }
220  }
221 
222  // // Assign history variable
223  if (F_pos > _H_old[_qp])
224  _H[_qp] = F_pos;
225  else
226  _H[_qp] = _H_old[_qp];
227 
228  Real hist_variable = _H_old[_qp];
229  if (_use_current_hist)
230  hist_variable = _H[_qp];
231 
232  if (hist_variable < _barrier[_qp])
233  hist_variable = _barrier[_qp];
234 
235  // Elastic free energy density
236  _E[_qp] = hist_variable * _D[_qp] + F_neg;
237  _dEdc[_qp] = hist_variable * _dDdc[_qp];
238  _d2Ed2c[_qp] = hist_variable * _d2Dd2c[_qp];
239 }
defineLegacyParams
defineLegacyParams(ComputeLinearElasticPFFractureStress)
ComputePFFractureStressBase
ComputePFFractureStressBase is the base class for stress in phase field fracture model.
Definition: ComputePFFractureStressBase.h:22
ComputePFFractureStressBase::_dstress_dc
MaterialProperty< RankTwoTensor > & _dstress_dc
Derivative of stress w.r.t damage variable.
Definition: ComputePFFractureStressBase.h:68
ComputePFFractureStressBase::validParams
static InputParameters validParams()
Definition: ComputePFFractureStressBase.C:15
ComputeStressBase::_stress
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
Definition: ComputeStressBase.h:50
ComputePFFractureStressBase::_d2Dd2c
const MaterialProperty< Real > & _d2Dd2c
Second-order derivative of degradation w.r.t damage variable.
Definition: ComputePFFractureStressBase.h:80
ComputeLinearElasticPFFractureStress::computeQpStress
virtual void computeQpStress() override
Compute the stress and store it in the _stress material property for the current quadrature point.
Definition: ComputeLinearElasticPFFractureStress.C:194
ComputeStressBase::_Jacobian_mult
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
Definition: ComputeStressBase.h:61
ComputeLinearElasticPFFractureStress::initialSetup
void initialSetup() override
Definition: ComputeLinearElasticPFFractureStress.C:40
ComputeLinearElasticPFFractureStress::Decomposition_type::strain_vol_dev
ComputePFFractureStressBase::_barrier
const MaterialProperty< Real > & _barrier
material property for fracture energy barrier
Definition: ComputePFFractureStressBase.h:56
ComputeLinearElasticPFFractureStress
Phase-field fracture This class computes the stress and energy contribution for the small strain Line...
Definition: ComputeLinearElasticPFFractureStress.h:26
ComputePFFractureStressBase::_d2Fdcdstrain
MaterialProperty< RankTwoTensor > & _d2Fdcdstrain
Second-order derivative of elastic energy w.r.t damage variable and strain.
Definition: ComputePFFractureStressBase.h:71
ComputePFFractureStressBase::_H
MaterialProperty< Real > & _H
History variable that prevents crack healing, declared in this material.
Definition: ComputePFFractureStressBase.h:50
ComputeLinearElasticPFFractureStress::computeStrainSpectral
void computeStrainSpectral(Real &F_pos, Real &F_neg)
Method to split elastic energy based on strain spectral decomposition.
Definition: ComputeLinearElasticPFFractureStress.C:51
ComputePFFractureStressBase::_dEdc
MaterialProperty< Real > & _dEdc
Derivative of elastic energy w.r.t damage variable.
Definition: ComputePFFractureStressBase.h:62
ComputePFFractureStressBase::_elasticity_tensor
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
Definition: ComputePFFractureStressBase.h:35
ComputePFFractureStressBase::_E
MaterialProperty< Real > & _E
Material property for elastic energy.
Definition: ComputePFFractureStressBase.h:59
ComputePFFractureStressBase::_d2Ed2c
MaterialProperty< Real > & _d2Ed2c
Second-order derivative of elastic energy w.r.t damage variable.
Definition: ComputePFFractureStressBase.h:65
ComputeLinearElasticPFFractureStress::computeStrainVolDev
void computeStrainVolDev(Real &F_pos, Real &F_neg)
Method to split elastic energy based on strain volumetric/deviatoric decomposition.
Definition: ComputeLinearElasticPFFractureStress.C:151
ComputeLinearElasticPFFractureStress::Decomposition_type::strain_spectral
ComputeLinearElasticPFFractureStress::Decomposition_type::stress_spectral
ComputeLinearElasticPFFractureStress.h
ComputeLinearElasticPFFractureStress::ComputeLinearElasticPFFractureStress
ComputeLinearElasticPFFractureStress(const InputParameters &parameters)
Definition: ComputeLinearElasticPFFractureStress.C:31
ComputePFFractureStressBase::_H_old
const MaterialProperty< Real > & _H_old
Old value of history variable.
Definition: ComputePFFractureStressBase.h:53
GuaranteeConsumer
Add-on class that provides the functionality to check if guarantees for material properties are provi...
Definition: GuaranteeConsumer.h:25
ComputePFFractureStressBase::_use_current_hist
bool _use_current_hist
Use current value of history variable.
Definition: ComputePFFractureStressBase.h:47
ComputePFFractureStressBase::_elasticity_tensor_name
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
Definition: ComputePFFractureStressBase.h:33
GuaranteeConsumer::hasGuaranteedMaterialProperty
bool hasGuaranteedMaterialProperty(const MaterialPropertyName &prop, Guarantee guarantee)
Definition: GuaranteeConsumer.C:28
registerMooseObject
registerMooseObject("TensorMechanicsApp", ComputeLinearElasticPFFractureStress)
RankFourTensorTempl< Real >
ComputeLinearElasticPFFractureStress::validParams
static InputParameters validParams()
Definition: ComputeLinearElasticPFFractureStress.C:18
ComputePFFractureStressBase::_D
const MaterialProperty< Real > & _D
Material property for energetic degradation function.
Definition: ComputePFFractureStressBase.h:74
ComputeLinearElasticPFFractureStress::Decomposition_type
Decomposition_type
Decomposittion type.
Definition: ComputeLinearElasticPFFractureStress.h:61
ComputeStressBase::_mechanical_strain
const MaterialProperty< RankTwoTensor > & _mechanical_strain
Mechanical strain material property.
Definition: ComputeStressBase.h:48
RankTwoTensorTempl< Real >
ComputePFFractureStressBase::_dDdc
const MaterialProperty< Real > & _dDdc
Derivative of degradation function w.r.t damage variable.
Definition: ComputePFFractureStressBase.h:77
ComputeLinearElasticPFFractureStress::_decomposition_type
enum ComputeLinearElasticPFFractureStress::Decomposition_type _decomposition_type
ComputeLinearElasticPFFractureStress::computeStressSpectral
void computeStressSpectral(Real &F_pos, Real &F_neg)
Method to split elastic energy based on stress spectral decomposition.
Definition: ComputeLinearElasticPFFractureStress.C:119
Guarantee::ISOTROPIC