www.mooseframework.org
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
ComputeLinearElasticPFFractureStress Class Reference

Phase-field fracture This class computes the stress and energy contribution for the small strain Linear Elastic formulation of phase field fracture. More...

#include <ComputeLinearElasticPFFractureStress.h>

Inheritance diagram for ComputeLinearElasticPFFractureStress:
[legend]

Public Member Functions

 ComputeLinearElasticPFFractureStress (const InputParameters &parameters)
 
void initialSetup () override
 

Protected Types

enum  Decomposition_type { Decomposition_type::strain_spectral, Decomposition_type::strain_vol_dev, Decomposition_type::stress_spectral, Decomposition_type::none }
 Decomposittion type. More...
 

Protected Member Functions

virtual void computeQpStress () override
 Compute the stress and store it in the _stress material property for the current quadrature point. More...
 
void computeStrainSpectral (Real &F_pos, Real &F_neg)
 Method to split elastic energy based on strain spectral decomposition. More...
 
void computeStrainVolDev (Real &F_pos, Real &F_neg)
 Method to split elastic energy based on strain volumetric/deviatoric decomposition. More...
 
void computeStressSpectral (Real &F_pos, Real &F_neg)
 Method to split elastic energy based on stress spectral decomposition. More...
 
virtual void initQpStatefulProperties () override
 
virtual void computeQpProperties () override
 
bool hasGuaranteedMaterialProperty (const MaterialPropertyName &prop, Guarantee guarantee)
 

Protected Attributes

enum ComputeLinearElasticPFFractureStress::Decomposition_type _decomposition_type
 
const std::string _elasticity_tensor_name
 Name of the elasticity tensor material property. More...
 
const MaterialProperty< RankFourTensor > & _elasticity_tensor
 Elasticity tensor material property. More...
 
const VariableValue & _c
 Coupled order parameter defining the crack. More...
 
const MaterialProperty< Real > & _l
 Material property defining crack width, declared elsewhere. More...
 
const MaterialProperty< Real > & _gc
 Material property defining gc parameter, declared elsewhere. More...
 
bool _use_current_hist
 Use current value of history variable. More...
 
MaterialProperty< Real > & _H
 History variable that prevents crack healing, declared in this material. More...
 
const MaterialProperty< Real > & _H_old
 Old value of history variable. More...
 
const MaterialProperty< Real > & _barrier
 material property for fracture energy barrier More...
 
MaterialProperty< Real > & _E
 Material property for elastic energy. More...
 
MaterialProperty< Real > & _dEdc
 Derivative of elastic energy w.r.t damage variable. More...
 
MaterialProperty< Real > & _d2Ed2c
 Second-order derivative of elastic energy w.r.t damage variable. More...
 
MaterialProperty< RankTwoTensor > & _dstress_dc
 Derivative of stress w.r.t damage variable. More...
 
MaterialProperty< RankTwoTensor > & _d2Fdcdstrain
 Second-order derivative of elastic energy w.r.t damage variable and strain. More...
 
const MaterialProperty< Real > & _D
 Material property for energetic degradation function. More...
 
const MaterialProperty< Real > & _dDdc
 Derivative of degradation function w.r.t damage variable. More...
 
const MaterialProperty< Real > & _d2Dd2c
 Second-order derivative of degradation w.r.t damage variable. More...
 
const std::string _base_name
 Base name prepended to all material property names to allow for multi-material systems. More...
 
const MaterialProperty< RankTwoTensor > & _mechanical_strain
 Mechanical strain material property. More...
 
MaterialProperty< RankTwoTensor > & _stress
 Stress material property. More...
 
MaterialProperty< RankTwoTensor > & _elastic_strain
 Elastic strain material property. More...
 
const MaterialProperty< RankTwoTensor > & _extra_stress
 Extra stress tensor. More...
 
std::vector< Function * > _initial_stress_fcn
 initial stress components More...
 
MaterialProperty< RankFourTensor > & _Jacobian_mult
 derivative of stress w.r.t. strain (_dstress_dstrain) More...
 

Detailed Description

Phase-field fracture This class computes the stress and energy contribution for the small strain Linear Elastic formulation of phase field fracture.

Definition at line 27 of file ComputeLinearElasticPFFractureStress.h.

Member Enumeration Documentation

◆ Decomposition_type

Decomposittion type.

Enumerator
strain_spectral 
strain_vol_dev 
stress_spectral 
none 

Definition at line 60 of file ComputeLinearElasticPFFractureStress.h.

61  {
62  strain_spectral,
63  strain_vol_dev,
64  stress_spectral,
65  none
enum ComputeLinearElasticPFFractureStress::Decomposition_type _decomposition_type

Constructor & Destructor Documentation

◆ ComputeLinearElasticPFFractureStress()

ComputeLinearElasticPFFractureStress::ComputeLinearElasticPFFractureStress ( const InputParameters &  parameters)

Definition at line 30 of file ComputeLinearElasticPFFractureStress.C.

32  : ComputePFFractureStressBase(parameters),
33  GuaranteeConsumer(this),
34  _decomposition_type(getParam<MooseEnum>("decomposition_type").getEnum<Decomposition_type>())
35 {
36 }
GuaranteeConsumer(MooseObject *moose_object)
ComputePFFractureStressBase(const InputParameters &parameters)
enum ComputeLinearElasticPFFractureStress::Decomposition_type _decomposition_type

Member Function Documentation

◆ computeQpProperties()

void ComputeStressBase::computeQpProperties ( )
overrideprotectedvirtualinherited

Definition at line 49 of file ComputeStressBase.C.

50 {
52 
53  // Add in extra stress
54  _stress[_qp] += _extra_stress[_qp];
55 }
virtual void computeQpStress()=0
Compute the stress and store it in the _stress material property for the current quadrature point...
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
const MaterialProperty< RankTwoTensor > & _extra_stress
Extra stress tensor.

◆ computeQpStress()

void ComputeLinearElasticPFFractureStress::computeQpStress ( )
overrideprotectedvirtual

Compute the stress and store it in the _stress material property for the current quadrature point.

Implements ComputeStressBase.

Definition at line 193 of file ComputeLinearElasticPFFractureStress.C.

194 {
195  Real F_pos, F_neg;
196 
197  switch (_decomposition_type)
198  {
200  computeStrainSpectral(F_pos, F_neg);
201  break;
203  computeStrainVolDev(F_pos, F_neg);
204  break;
206  computeStressSpectral(F_pos, F_neg);
207  break;
208  default:
209  {
210  _stress[_qp] = _D[_qp] * _elasticity_tensor[_qp] * _mechanical_strain[_qp];
211  F_pos = (_stress[_qp]).doubleContraction(_mechanical_strain[_qp]) / 2.0;
212  F_neg = 0.0;
213  if (_use_current_hist)
214  _d2Fdcdstrain[_qp] = _stress[_qp] * _dDdc[_qp];
215 
216  _dstress_dc[_qp] = _stress[_qp] * _dDdc[_qp];
217  _Jacobian_mult[_qp] = _D[_qp] * _elasticity_tensor[_qp];
218  }
219  }
220 
221  // // Assign history variable
222  if (F_pos > _H_old[_qp])
223  _H[_qp] = F_pos;
224  else
225  _H[_qp] = _H_old[_qp];
226 
227  Real hist_variable = _H_old[_qp];
228  if (_use_current_hist)
229  hist_variable = _H[_qp];
230 
231  if (hist_variable < _barrier[_qp])
232  hist_variable = _barrier[_qp];
233 
234  // Elastic free energy density
235  _E[_qp] = hist_variable * _D[_qp] + F_neg;
236  _dEdc[_qp] = hist_variable * _dDdc[_qp];
237  _d2Ed2c[_qp] = hist_variable * _d2Dd2c[_qp];
238 }
MaterialProperty< Real > & _H
History variable that prevents crack healing, declared in this material.
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
const MaterialProperty< Real > & _d2Dd2c
Second-order derivative of degradation w.r.t damage variable.
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
MaterialProperty< Real > & _E
Material property for elastic energy.
MaterialProperty< Real > & _d2Ed2c
Second-order derivative of elastic energy w.r.t damage variable.
MaterialProperty< RankTwoTensor > & _dstress_dc
Derivative of stress w.r.t damage variable.
MaterialProperty< RankTwoTensor > & _d2Fdcdstrain
Second-order derivative of elastic energy w.r.t damage variable and strain.
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
void computeStrainSpectral(Real &F_pos, Real &F_neg)
Method to split elastic energy based on strain spectral decomposition.
const MaterialProperty< Real > & _barrier
material property for fracture energy barrier
void computeStressSpectral(Real &F_pos, Real &F_neg)
Method to split elastic energy based on stress spectral decomposition.
enum ComputeLinearElasticPFFractureStress::Decomposition_type _decomposition_type
const MaterialProperty< Real > & _dDdc
Derivative of degradation function w.r.t damage variable.
MaterialProperty< Real > & _dEdc
Derivative of elastic energy w.r.t damage variable.
bool _use_current_hist
Use current value of history variable.
const MaterialProperty< RankTwoTensor > & _mechanical_strain
Mechanical strain material property.
const MaterialProperty< Real > & _H_old
Old value of history variable.
const MaterialProperty< Real > & _D
Material property for energetic degradation function.
void computeStrainVolDev(Real &F_pos, Real &F_neg)
Method to split elastic energy based on strain volumetric/deviatoric decomposition.

◆ computeStrainSpectral()

void ComputeLinearElasticPFFractureStress::computeStrainSpectral ( Real &  F_pos,
Real &  F_neg 
)
protected

Method to split elastic energy based on strain spectral decomposition.

Parameters
F_postensile part of total elastic energy
F_negcompressive part of total elastic energy

Definition at line 50 of file ComputeLinearElasticPFFractureStress.C.

Referenced by computeQpStress().

51 {
52  // Isotropic elasticity is assumed and should be enforced
53  const Real lambda = _elasticity_tensor[_qp](0, 0, 1, 1);
54  const Real mu = _elasticity_tensor[_qp](0, 1, 0, 1);
55 
56  // Compute eigenvectors and eigenvalues of mechanical strain and projection tensor
57  RankTwoTensor eigvec;
58  std::vector<Real> eigval(LIBMESH_DIM);
59  RankFourTensor Ppos =
60  _mechanical_strain[_qp].positiveProjectionEigenDecomposition(eigval, eigvec);
61  RankFourTensor I4sym(RankFourTensor::initIdentitySymmetricFour);
62 
63  // Calculate tensors of outerproduct of eigen vectors
64  std::vector<RankTwoTensor> etens(LIBMESH_DIM);
65 
66  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
67  etens[i].vectorOuterProduct(eigvec.column(i), eigvec.column(i));
68 
69  // Separate out positive and negative eigen values
70  std::vector<Real> epos(LIBMESH_DIM), eneg(LIBMESH_DIM);
71  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
72  {
73  epos[i] = (std::abs(eigval[i]) + eigval[i]) / 2.0;
74  eneg[i] = -(std::abs(eigval[i]) - eigval[i]) / 2.0;
75  }
76 
77  // Seprate positive and negative sums of all eigenvalues
78  Real etr = 0.0;
79  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
80  etr += eigval[i];
81 
82  const Real etrpos = (std::abs(etr) + etr) / 2.0;
83  const Real etrneg = -(std::abs(etr) - etr) / 2.0;
84 
85  // Calculate the tensile (postive) and compressive (negative) parts of stress
86  RankTwoTensor stress0pos, stress0neg;
87  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
88  {
89  stress0pos += etens[i] * (lambda * etrpos + 2.0 * mu * epos[i]);
90  stress0neg += etens[i] * (lambda * etrneg + 2.0 * mu * eneg[i]);
91  }
92 
93  // sum squares of epos and eneg
94  Real pval(0.0), nval(0.0);
95  for (unsigned int i = 0; i < LIBMESH_DIM; ++i)
96  {
97  pval += epos[i] * epos[i];
98  nval += eneg[i] * eneg[i];
99  }
100 
101  _stress[_qp] = stress0pos * _D[_qp] + stress0neg;
102 
103  // Energy with positive principal strains
104  F_pos = lambda * etrpos * etrpos / 2.0 + mu * pval;
105  F_neg = -lambda * etrneg * etrneg / 2.0 + mu * nval;
106 
107  // 2nd derivative wrt c and strain = 0.0 if we used the previous step's history varible
108  if (_use_current_hist)
109  _d2Fdcdstrain[_qp] = stress0pos * _dDdc[_qp];
110 
111  // Used in StressDivergencePFFracTensors off-diagonal Jacobian
112  _dstress_dc[_qp] = stress0pos * _dDdc[_qp];
113 
114  _Jacobian_mult[_qp] = (I4sym - (1 - _D[_qp]) * Ppos) * _elasticity_tensor[_qp];
115 }
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
MaterialProperty< RankTwoTensor > & _dstress_dc
Derivative of stress w.r.t damage variable.
MaterialProperty< RankTwoTensor > & _d2Fdcdstrain
Second-order derivative of elastic energy w.r.t damage variable and strain.
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
const MaterialProperty< Real > & _dDdc
Derivative of degradation function w.r.t damage variable.
bool _use_current_hist
Use current value of history variable.
const MaterialProperty< RankTwoTensor > & _mechanical_strain
Mechanical strain material property.
const MaterialProperty< Real > & _D
Material property for energetic degradation function.

◆ computeStrainVolDev()

void ComputeLinearElasticPFFractureStress::computeStrainVolDev ( Real &  F_pos,
Real &  F_neg 
)
protected

Method to split elastic energy based on strain volumetric/deviatoric decomposition.

Parameters
F_postensile part of total elastic energy
F_negcompressive part of total elastic energy

Definition at line 150 of file ComputeLinearElasticPFFractureStress.C.

Referenced by computeQpStress().

151 {
152  // Isotropic elasticity is assumed and should be enforced
153  const Real lambda = _elasticity_tensor[_qp](0, 0, 1, 1);
154  const Real mu = _elasticity_tensor[_qp](0, 1, 0, 1);
155  const Real k = lambda + 2.0 * mu / LIBMESH_DIM;
156 
157  RankTwoTensor I2(RankTwoTensor::initIdentity);
158  RankFourTensor I2I2 = I2.outerProduct(I2);
159 
160  RankFourTensor Jacobian_pos, Jacobian_neg;
161  RankTwoTensor strain0vol, strain0dev;
162  RankTwoTensor stress0pos, stress0neg;
163  Real strain0tr, strain0tr_neg, strain0tr_pos;
164 
165  strain0dev = _mechanical_strain[_qp].deviatoric();
166  strain0vol = _mechanical_strain[_qp] - strain0dev;
167  strain0tr = _mechanical_strain[_qp].trace();
168  strain0tr_neg = std::min(strain0tr, 0.0);
169  strain0tr_pos = strain0tr - strain0tr_neg;
170  stress0neg = k * strain0tr_neg * I2;
171  stress0pos = _elasticity_tensor[_qp] * _mechanical_strain[_qp] - stress0neg;
172  // Energy with positive principal strains
173  RankTwoTensor strain0dev2 = strain0dev * strain0dev;
174  F_pos = 0.5 * k * strain0tr_pos * strain0tr_pos + mu * strain0dev2.trace();
175  F_neg = 0.5 * k * strain0tr_neg * strain0tr_neg;
176 
177  _stress[_qp] = stress0pos * _D[_qp] + stress0neg;
178 
179  // 2nd derivative wrt c and strain = 0.0 if we used the previous step's history varible
180  if (_use_current_hist)
181  _d2Fdcdstrain[_qp] = stress0pos * _dDdc[_qp];
182 
183  // Used in StressDivergencePFFracTensors off-diagonal Jacobian
184  _dstress_dc[_qp] = stress0pos * _dDdc[_qp];
185 
186  if (strain0tr < 0)
187  Jacobian_neg = k * I2I2;
188  Jacobian_pos = _elasticity_tensor[_qp] - Jacobian_neg;
189  _Jacobian_mult[_qp] = _D[_qp] * Jacobian_pos + Jacobian_neg;
190 }
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
MaterialProperty< RankTwoTensor > & _dstress_dc
Derivative of stress w.r.t damage variable.
MaterialProperty< RankTwoTensor > & _d2Fdcdstrain
Second-order derivative of elastic energy w.r.t damage variable and strain.
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
const MaterialProperty< Real > & _dDdc
Derivative of degradation function w.r.t damage variable.
bool _use_current_hist
Use current value of history variable.
const MaterialProperty< RankTwoTensor > & _mechanical_strain
Mechanical strain material property.
const MaterialProperty< Real > & _D
Material property for energetic degradation function.

◆ computeStressSpectral()

void ComputeLinearElasticPFFractureStress::computeStressSpectral ( Real &  F_pos,
Real &  F_neg 
)
protected

Method to split elastic energy based on stress spectral decomposition.

Parameters
F_postensile part of total elastic energy
F_negcompressive part of total elastic energy

Definition at line 118 of file ComputeLinearElasticPFFractureStress.C.

Referenced by computeQpStress().

119 {
120  // Compute Uncracked stress
122 
123  // Create the positive and negative projection tensors
124  RankFourTensor I4sym(RankFourTensor::initIdentitySymmetricFour);
125  std::vector<Real> eigval;
126  RankTwoTensor eigvec;
127  RankFourTensor Ppos = stress.positiveProjectionEigenDecomposition(eigval, eigvec);
128 
129  // Project the positive and negative stresses
130  RankTwoTensor stress0pos = Ppos * stress;
131  RankTwoTensor stress0neg = stress - stress0pos;
132 
133  // Compute the positive and negative elastic energies
134  F_pos = (stress0pos).doubleContraction(_mechanical_strain[_qp]) / 2.0;
135  F_neg = (stress0neg).doubleContraction(_mechanical_strain[_qp]) / 2.0;
136 
137  _stress[_qp] = stress0pos * _D[_qp] + stress0neg;
138 
139  // 2nd derivative wrt c and strain = 0.0 if we used the previous step's history varible
140  if (_use_current_hist)
141  _d2Fdcdstrain[_qp] = stress0pos * _dDdc[_qp];
142 
143  // Used in StressDivergencePFFracTensors off-diagonal Jacobian
144  _dstress_dc[_qp] = stress0pos * _dDdc[_qp];
145 
146  _Jacobian_mult[_qp] = (I4sym - (1 - _D[_qp]) * Ppos) * _elasticity_tensor[_qp];
147 }
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
const MaterialProperty< RankFourTensor > & _elasticity_tensor
Elasticity tensor material property.
MaterialProperty< RankTwoTensor > & _dstress_dc
Derivative of stress w.r.t damage variable.
MaterialProperty< RankTwoTensor > & _d2Fdcdstrain
Second-order derivative of elastic energy w.r.t damage variable and strain.
MaterialProperty< RankTwoTensor > & _stress
Stress material property.
const MaterialProperty< Real > & _dDdc
Derivative of degradation function w.r.t damage variable.
bool _use_current_hist
Use current value of history variable.
const MaterialProperty< RankTwoTensor > & _mechanical_strain
Mechanical strain material property.
const MaterialProperty< Real > & _D
Material property for energetic degradation function.

◆ hasGuaranteedMaterialProperty()

bool GuaranteeConsumer::hasGuaranteedMaterialProperty ( const MaterialPropertyName &  prop,
Guarantee  guarantee 
)
protectedinherited

Definition at line 28 of file GuaranteeConsumer.C.

Referenced by ComputeFiniteStrainElasticStress::initialSetup(), ComputeSmearedCrackingStress::initialSetup(), initialSetup(), and ComputeMultipleInelasticStress::initialSetup().

30 {
31  if (!_gc_feproblem->startedInitialSetup())
32  mooseError("hasGuaranteedMaterialProperty() needs to be called in initialSetup()");
33 
34  // Reference to MaterialWarehouse for testing and retrieving block ids
35  const auto & warehouse = _gc_feproblem->getMaterialWarehouse();
36 
37  // Complete set of ids that this object is active
38  const auto & ids = (_gc_block_restrict && _gc_block_restrict->blockRestricted())
39  ? _gc_block_restrict->blockIDs()
40  : _gc_feproblem->mesh().meshSubdomains();
41 
42  // Loop over each id for this object
43  for (const auto & id : ids)
44  {
45  // If block materials exist, look if any issue the required guarantee
46  if (warehouse.hasActiveBlockObjects(id))
47  {
48  const std::vector<std::shared_ptr<Material>> & mats = warehouse.getActiveBlockObjects(id);
49  for (const auto & mat : mats)
50  {
51  const auto & mat_props = mat->getSuppliedItems();
52  if (mat_props.count(prop_name))
53  {
54  auto guarantee_mat = dynamic_cast<GuaranteeProvider *>(mat.get());
55  if (guarantee_mat && !guarantee_mat->hasGuarantee(prop_name, guarantee))
56  {
57  // we found at least one material on the set of block we operate on
58  // that does _not_ provide the requested guarantee
59  return false;
60  }
61  }
62  }
63  }
64  }
65 
66  return true;
67 }
Add-on class that provides the functionality to issue guarantees for declared material properties...
BlockRestrictable *const _gc_block_restrict
Access block restrictions of the object with this interface.
FEProblemBase *const _gc_feproblem
Reference to the FEProblemBase class.

◆ initialSetup()

void ComputeLinearElasticPFFractureStress::initialSetup ( )
override

Definition at line 39 of file ComputeLinearElasticPFFractureStress.C.

40 {
44  mooseError("Decomposition approach of strain_vol_dev and strain_spectral can only be used with "
45  "isotropic elasticity tensor materials, use stress_spectral for anistropic "
46  "elasticity tensor materials");
47 }
enum ComputeLinearElasticPFFractureStress::Decomposition_type _decomposition_type
const std::string _elasticity_tensor_name
Name of the elasticity tensor material property.
bool hasGuaranteedMaterialProperty(const MaterialPropertyName &prop, Guarantee guarantee)

◆ initQpStatefulProperties()

void ComputePFFractureStressBase::initQpStatefulProperties ( )
overrideprotectedvirtualinherited

Reimplemented from ComputeStressBase.

Definition at line 60 of file ComputePFFractureStressBase.C.

61 {
62  _H[_qp] = 0.0;
63 }
MaterialProperty< Real > & _H
History variable that prevents crack healing, declared in this material.

Member Data Documentation

◆ _barrier

const MaterialProperty<Real>& ComputePFFractureStressBase::_barrier
protectedinherited

material property for fracture energy barrier

Definition at line 55 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress().

◆ _base_name

const std::string ComputeStressBase::_base_name
protectedinherited

Base name prepended to all material property names to allow for multi-material systems.

Definition at line 44 of file ComputeStressBase.h.

Referenced by ComputeLinearElasticStress::initialSetup(), and ComputeCosseratLinearElasticStress::initialSetup().

◆ _c

const VariableValue& ComputePFFractureStressBase::_c
protectedinherited

Coupled order parameter defining the crack.

Definition at line 37 of file ComputePFFractureStressBase.h.

◆ _D

const MaterialProperty<Real>& ComputePFFractureStressBase::_D
protectedinherited

Material property for energetic degradation function.

Definition at line 73 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress(), computeStrainSpectral(), computeStrainVolDev(), and computeStressSpectral().

◆ _d2Dd2c

const MaterialProperty<Real>& ComputePFFractureStressBase::_d2Dd2c
protectedinherited

Second-order derivative of degradation w.r.t damage variable.

Definition at line 79 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress().

◆ _d2Ed2c

MaterialProperty<Real>& ComputePFFractureStressBase::_d2Ed2c
protectedinherited

Second-order derivative of elastic energy w.r.t damage variable.

Definition at line 64 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress().

◆ _d2Fdcdstrain

MaterialProperty<RankTwoTensor>& ComputePFFractureStressBase::_d2Fdcdstrain
protectedinherited

Second-order derivative of elastic energy w.r.t damage variable and strain.

Definition at line 70 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress(), computeStrainSpectral(), computeStrainVolDev(), and computeStressSpectral().

◆ _dDdc

const MaterialProperty<Real>& ComputePFFractureStressBase::_dDdc
protectedinherited

Derivative of degradation function w.r.t damage variable.

Definition at line 76 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress(), computeStrainSpectral(), computeStrainVolDev(), and computeStressSpectral().

◆ _decomposition_type

enum ComputeLinearElasticPFFractureStress::Decomposition_type ComputeLinearElasticPFFractureStress::_decomposition_type
protected

Referenced by computeQpStress(), and initialSetup().

◆ _dEdc

MaterialProperty<Real>& ComputePFFractureStressBase::_dEdc
protectedinherited

Derivative of elastic energy w.r.t damage variable.

Definition at line 61 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress().

◆ _dstress_dc

MaterialProperty<RankTwoTensor>& ComputePFFractureStressBase::_dstress_dc
protectedinherited

Derivative of stress w.r.t damage variable.

Definition at line 67 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress(), computeStrainSpectral(), computeStrainVolDev(), and computeStressSpectral().

◆ _E

MaterialProperty<Real>& ComputePFFractureStressBase::_E
protectedinherited

Material property for elastic energy.

Definition at line 58 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress().

◆ _elastic_strain

MaterialProperty<RankTwoTensor>& ComputeStressBase::_elastic_strain
protectedinherited

◆ _elasticity_tensor

const MaterialProperty<RankFourTensor>& ComputePFFractureStressBase::_elasticity_tensor
protectedinherited

Elasticity tensor material property.

Definition at line 34 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress(), computeStrainSpectral(), computeStrainVolDev(), and computeStressSpectral().

◆ _elasticity_tensor_name

const std::string ComputePFFractureStressBase::_elasticity_tensor_name
protectedinherited

Name of the elasticity tensor material property.

Definition at line 32 of file ComputePFFractureStressBase.h.

Referenced by initialSetup().

◆ _extra_stress

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_extra_stress
protectedinherited

Extra stress tensor.

Definition at line 54 of file ComputeStressBase.h.

Referenced by ComputeStressBase::computeQpProperties().

◆ _gc

const MaterialProperty<Real>& ComputePFFractureStressBase::_gc
protectedinherited

Material property defining gc parameter, declared elsewhere.

Definition at line 43 of file ComputePFFractureStressBase.h.

◆ _H

MaterialProperty<Real>& ComputePFFractureStressBase::_H
protectedinherited

History variable that prevents crack healing, declared in this material.

Definition at line 49 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress(), and ComputePFFractureStressBase::initQpStatefulProperties().

◆ _H_old

const MaterialProperty<Real>& ComputePFFractureStressBase::_H_old
protectedinherited

Old value of history variable.

Definition at line 52 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress().

◆ _initial_stress_fcn

std::vector<Function *> ComputeStressBase::_initial_stress_fcn
protectedinherited

initial stress components

Definition at line 57 of file ComputeStressBase.h.

◆ _Jacobian_mult

MaterialProperty<RankFourTensor>& ComputeStressBase::_Jacobian_mult
protectedinherited

◆ _l

const MaterialProperty<Real>& ComputePFFractureStressBase::_l
protectedinherited

Material property defining crack width, declared elsewhere.

Definition at line 40 of file ComputePFFractureStressBase.h.

◆ _mechanical_strain

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_mechanical_strain
protectedinherited

◆ _stress

MaterialProperty<RankTwoTensor>& ComputeStressBase::_stress
protectedinherited

Stress material property.

Definition at line 49 of file ComputeStressBase.h.

Referenced by ComputeMultipleInelasticCosseratStress::computeAdmissibleState(), ComputeMultipleInelasticStress::computeAdmissibleState(), ComputeStressBase::computeQpProperties(), ComputeStrainIncrementBasedStress::computeQpStress(), ComputeLinearElasticStress::computeQpStress(), ComputeFiniteStrainElasticStress::computeQpStress(), ComputeCosseratLinearElasticStress::computeQpStress(), ComputeDamageStress::computeQpStress(), ComputeSmearedCrackingStress::computeQpStress(), FiniteStrainPlasticMaterial::computeQpStress(), computeQpStress(), ComputeMultiPlasticityStress::computeQpStress(), ComputeLinearViscoelasticStress::computeQpStress(), ComputeMultipleInelasticStress::computeQpStressIntermediateConfiguration(), computeStrainSpectral(), computeStrainVolDev(), computeStressSpectral(), ComputeMultipleInelasticStress::finiteStrainRotation(), ComputeStressBase::initQpStatefulProperties(), FiniteStrainCrystalPlasticity::initQpStatefulProperties(), FiniteStrainUObasedCP::initQpStatefulProperties(), FiniteStrainHyperElasticViscoPlastic::initQpStatefulProperties(), ComputeMultiPlasticityStress::postReturnMap(), FiniteStrainUObasedCP::postSolveQp(), FiniteStrainHyperElasticViscoPlastic::postSolveQp(), FiniteStrainCrystalPlasticity::postSolveQp(), ComputeSmearedCrackingStress::updateCrackingStateAndStress(), ComputeMultipleInelasticStress::updateQpState(), and ComputeMultipleInelasticStress::updateQpStateSingleModel().

◆ _use_current_hist

bool ComputePFFractureStressBase::_use_current_hist
protectedinherited

Use current value of history variable.

Definition at line 46 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress(), computeStrainSpectral(), computeStrainVolDev(), and computeStressSpectral().


The documentation for this class was generated from the following files: