www.mooseframework.org
ComputeCrackedStress.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 
10 #include "ComputeCrackedStress.h"
11 
12 registerMooseObject("TensorMechanicsApp", ComputeCrackedStress);
13 
15 
16 InputParameters
18 {
19  InputParameters params = Material::validParams();
20  params.addClassDescription("Computes energy and modifies the stress for phase field fracture");
21  params.addRequiredCoupledVar("c", "Order parameter for damage");
22  params.addParam<Real>("kdamage", 1e-9, "Stiffness of damaged matrix");
23  params.addParam<bool>("finite_strain_model", false, "The model is using finite strain");
24  params.addParam<bool>(
25  "use_current_history_variable", false, "Use the current value of the history variable.");
26  params.addParam<MaterialPropertyName>(
27  "F_name", "E_el", "Name of material property storing the elastic energy");
28  params.addParam<MaterialPropertyName>(
29  "kappa_name",
30  "kappa_op",
31  "Name of material property being created to store the interfacial parameter kappa");
32  params.addParam<MaterialPropertyName>(
33  "mobility_name", "L", "Name of material property being created to store the mobility L");
34  params.addParam<std::string>("base_name", "The base name used to save the cracked stress");
35  params.addRequiredParam<std::string>("uncracked_base_name",
36  "The base name used to calculate the original stress");
37  return params;
38 }
39 
40 ComputeCrackedStress::ComputeCrackedStress(const InputParameters & parameters)
41  : DerivativeMaterialInterface<Material>(parameters),
42  _base_name(isParamValid("base_name") ? getParam<std::string>("base_name") + "_" : ""),
43  _uncracked_base_name(getParam<std::string>("uncracked_base_name") + "_"),
44  _finite_strain_model(getParam<bool>("finite_strain_model")),
45  _use_current_hist(getParam<bool>("use_current_history_variable")),
46  _strain(
47  _finite_strain_model
48  ? getMaterialPropertyByName<RankTwoTensor>(_uncracked_base_name + "elastic_strain")
49  : getMaterialPropertyByName<RankTwoTensor>(_uncracked_base_name + "mechanical_strain")),
50  _uncracked_stress(getMaterialPropertyByName<RankTwoTensor>(_uncracked_base_name + "stress")),
51  _uncracked_Jacobian_mult(
52  getMaterialPropertyByName<RankFourTensor>(_uncracked_base_name + "Jacobian_mult")),
53  _c(coupledValue("c")),
54  _gc_prop(getMaterialProperty<Real>("gc_prop")),
55  _l(getMaterialProperty<Real>("l")),
56  _visco(getMaterialProperty<Real>("visco")),
57  _kdamage(getParam<Real>("kdamage")),
58  _stress(declareProperty<RankTwoTensor>(_base_name + "stress")),
59  _F(declareProperty<Real>(getParam<MaterialPropertyName>("F_name"))),
60  _dFdc(declarePropertyDerivative<Real>(getParam<MaterialPropertyName>("F_name"),
61  getVar("c", 0)->name())),
62  _d2Fdc2(declarePropertyDerivative<Real>(
63  getParam<MaterialPropertyName>("F_name"), getVar("c", 0)->name(), getVar("c", 0)->name())),
64  _d2Fdcdstrain(declareProperty<RankTwoTensor>("d2Fdcdstrain")),
65  _dstress_dc(declarePropertyDerivative<RankTwoTensor>("stress", getVar("c", 0)->name())),
66  _hist(declareProperty<Real>("hist")),
67  _hist_old(getMaterialPropertyOld<Real>("hist")),
68  _Jacobian_mult(declareProperty<RankFourTensor>(_base_name + "Jacobian_mult")),
69  _kappa(declareProperty<Real>(getParam<MaterialPropertyName>("kappa_name"))),
70  _L(declareProperty<Real>(getParam<MaterialPropertyName>("mobility_name")))
71 {
72 }
73 
74 void
76 {
77  _stress[_qp].zero();
78  _hist[_qp] = 0.0;
79 }
80 
81 void
83 {
84  const Real c = _c[_qp];
85 
86  // Zero out values when c > 1
87  Real cfactor = 1.0;
88  if (c > 1.0)
89  cfactor = 0.0;
90 
91  // Create the positive and negative projection tensors
92  RankFourTensor I4sym(RankFourTensor::initIdentitySymmetricFour);
93  std::vector<Real> eigval;
94  RankTwoTensor eigvec;
95  RankFourTensor Ppos = _uncracked_stress[_qp].positiveProjectionEigenDecomposition(eigval, eigvec);
96  RankFourTensor Pneg = I4sym - Ppos;
97 
98  // Project the positive and negative stresses
99  RankTwoTensor stress0pos = Ppos * _uncracked_stress[_qp];
100  RankTwoTensor stress0neg = Pneg * _uncracked_stress[_qp];
101 
102  // Compute the positive and negative elastic energies
103  Real G0_pos = (stress0pos).doubleContraction(_strain[_qp]) / 2.0;
104  Real G0_neg = (stress0neg).doubleContraction(_strain[_qp]) / 2.0;
105 
106  // Update the history variable
107  if (G0_pos > _hist_old[_qp])
108  _hist[_qp] = G0_pos;
109  else
110  _hist[_qp] = _hist_old[_qp];
111 
112  Real hist_variable = _hist_old[_qp];
113  if (_use_current_hist)
114  hist_variable = _hist[_qp];
115 
116  // Compute degredation function and derivatives
117  Real h = cfactor * (1.0 - c) * (1.0 - c) * (1.0 - _kdamage) + _kdamage;
118  Real dhdc = -2.0 * cfactor * (1.0 - c) * (1.0 - _kdamage);
119  Real d2hdc2 = 2.0 * cfactor * (1.0 - _kdamage);
120 
121  // Compute stress and its derivatives
122  _stress[_qp] = (Ppos * h + Pneg) * _uncracked_stress[_qp];
123  _dstress_dc[_qp] = stress0pos * dhdc;
124  _Jacobian_mult[_qp] = (Ppos * h + Pneg) * _uncracked_Jacobian_mult[_qp];
125 
126  // Compute energy and its derivatives
127  _F[_qp] = hist_variable * h - G0_neg + _gc_prop[_qp] * c * c / (2 * _l[_qp]);
128  _dFdc[_qp] = hist_variable * dhdc + _gc_prop[_qp] * c / _l[_qp];
129  _d2Fdc2[_qp] = hist_variable * d2hdc2 + _gc_prop[_qp] / _l[_qp];
130 
131  // 2nd derivative wrt c and strain = 0.0 if we used the previous step's history varible
132  if (_use_current_hist)
133  _d2Fdcdstrain[_qp] = stress0pos * dhdc;
134 
135  // Assign L and kappa
136  _kappa[_qp] = _gc_prop[_qp] * _l[_qp];
137  _L[_qp] = 1.0 / (_gc_prop[_qp] * _visco[_qp]);
138 }
ComputeCrackedStress::computeQpProperties
virtual void computeQpProperties()
Definition: ComputeCrackedStress.C:82
registerMooseObject
registerMooseObject("TensorMechanicsApp", ComputeCrackedStress)
ComputeCrackedStress::_d2Fdc2
MaterialProperty< Real > & _d2Fdc2
Definition: ComputeCrackedStress.h:78
ComputeCrackedStress::_uncracked_Jacobian_mult
const MaterialProperty< RankFourTensor > & _uncracked_Jacobian_mult
Uncracked Jacobian_mult calculated by another material.
Definition: ComputeCrackedStress.h:56
ComputeCrackedStress::_c
const VariableValue & _c
Variable defining the phase field damage parameter.
Definition: ComputeCrackedStress.h:59
ComputeCrackedStress::_F
MaterialProperty< Real > & _F
Definition: ComputeCrackedStress.h:76
ComputeCrackedStress::_Jacobian_mult
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
Definition: ComputeCrackedStress.h:87
ComputeCrackedStress::_strain
const MaterialProperty< RankTwoTensor > & _strain
Mechanical_strain if finite_strain_model = false, otherwise elastic_strain.
Definition: ComputeCrackedStress.h:50
ComputeCrackedStress::_hist
MaterialProperty< Real > & _hist
history variable storing the maximum positive deformation energy
Definition: ComputeCrackedStress.h:83
ComputeCrackedStress::_l
const MaterialProperty< Real > & _l
Characteristic length, controls damage zone thickness.
Definition: ComputeCrackedStress.h:65
ComputeCrackedStress::_dstress_dc
MaterialProperty< RankTwoTensor > & _dstress_dc
Definition: ComputeCrackedStress.h:80
defineLegacyParams
defineLegacyParams(ComputeCrackedStress)
ComputeCrackedStress::_d2Fdcdstrain
MaterialProperty< RankTwoTensor > & _d2Fdcdstrain
Definition: ComputeCrackedStress.h:79
ComputeCrackedStress::_kappa
MaterialProperty< Real > & _kappa
Property where the value for kappa will be defined.
Definition: ComputeCrackedStress.h:90
ComputeCrackedStress::initQpStatefulProperties
virtual void initQpStatefulProperties()
Definition: ComputeCrackedStress.C:75
ComputeCrackedStress::ComputeCrackedStress
ComputeCrackedStress(const InputParameters &parameters)
Definition: ComputeCrackedStress.C:40
validParams
InputParameters validParams()
ComputeCrackedStress::_stress
MaterialProperty< RankTwoTensor > & _stress
Stress being computed by this kernel.
Definition: ComputeCrackedStress.h:74
ComputeCrackedStress::_visco
const MaterialProperty< Real > & _visco
Definition: ComputeCrackedStress.h:68
ComputeCrackedStress.h
name
const std::string name
Definition: Setup.h:21
ComputeCrackedStress::_dFdc
MaterialProperty< Real > & _dFdc
Definition: ComputeCrackedStress.h:77
ComputeCrackedStress::_use_current_hist
bool _use_current_hist
Use current value of history variable.
Definition: ComputeCrackedStress.h:47
ComputeCrackedStress::_kdamage
Real _kdamage
Small number to avoid non-positive definiteness at or near complete damage.
Definition: ComputeCrackedStress.h:71
RankFourTensorTempl< Real >
ComputeCrackedStress
Computes energy and modifies the stress for phase field fracture.
Definition: ComputeCrackedStress.h:27
ComputeCrackedStress::_L
MaterialProperty< Real > & _L
Property where the value for L will be defined.
Definition: ComputeCrackedStress.h:93
ComputeCrackedStress::_gc_prop
const MaterialProperty< Real > & _gc_prop
Critical energy release rate for fracture.
Definition: ComputeCrackedStress.h:62
RankTwoTensorTempl
Definition: ACGrGrElasticDrivingForce.h:17
ComputeCrackedStress::validParams
static InputParameters validParams()
Definition: ComputeCrackedStress.C:17
ComputeCrackedStress::_hist_old
const MaterialProperty< Real > & _hist_old
Definition: ComputeCrackedStress.h:84
ComputeCrackedStress::_uncracked_stress
const MaterialProperty< RankTwoTensor > & _uncracked_stress
Uncracked stress calculated by another material.
Definition: ComputeCrackedStress.h:53