www.mooseframework.org
Public Member Functions | Static Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | Private 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
 

Static Public Member Functions

static InputParameters validParams ()
 

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< const Function * > _initial_stress_fcn
 initial stress components More...
 
MaterialProperty< RankFourTensor > & _Jacobian_mult
 derivative of stress w.r.t. strain (_dstress_dstrain) More...
 

Private Attributes

const InputParameters & _gc_params
 Parameters of the object with this interface. More...
 
FEProblemBase *const _gc_feproblem
 Reference to the FEProblemBase class. More...
 
BlockRestrictable *const _gc_block_restrict
 Access block restrictions of the object with this interface. 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 26 of file ComputeLinearElasticPFFractureStress.h.

Member Enumeration Documentation

◆ Decomposition_type

Decomposittion type.

Enumerator
strain_spectral 
strain_vol_dev 
stress_spectral 
none 

Definition at line 61 of file ComputeLinearElasticPFFractureStress.h.

62  {
63  strain_spectral,
64  strain_vol_dev,
65  stress_spectral,
66  none

Constructor & Destructor Documentation

◆ ComputeLinearElasticPFFractureStress()

ComputeLinearElasticPFFractureStress::ComputeLinearElasticPFFractureStress ( const InputParameters &  parameters)

Definition at line 31 of file ComputeLinearElasticPFFractureStress.C.

33  : ComputePFFractureStressBase(parameters),
34  GuaranteeConsumer(this),
35  _decomposition_type(getParam<MooseEnum>("decomposition_type").getEnum<Decomposition_type>())
36 {
37 }

Member Function Documentation

◆ computeQpProperties()

void ComputeStressBase::computeQpProperties ( )
overrideprotectedvirtualinherited

Definition at line 50 of file ComputeStressBase.C.

51 {
53 
54  // Add in extra stress
55  _stress[_qp] += _extra_stress[_qp];
56 }

◆ 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 194 of file ComputeLinearElasticPFFractureStress.C.

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 }

◆ 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 51 of file ComputeLinearElasticPFFractureStress.C.

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 }

Referenced by computeQpStress().

◆ 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 151 of file ComputeLinearElasticPFFractureStress.C.

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 }

Referenced by computeQpStress().

◆ 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 119 of file ComputeLinearElasticPFFractureStress.C.

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 }

Referenced by computeQpStress().

◆ hasGuaranteedMaterialProperty()

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

Definition at line 28 of file GuaranteeConsumer.C.

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<MaterialBase>> & 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 }

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

◆ initialSetup()

void ComputeLinearElasticPFFractureStress::initialSetup ( )
override

Definition at line 40 of file ComputeLinearElasticPFFractureStress.C.

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 }

◆ initQpStatefulProperties()

void ComputePFFractureStressBase::initQpStatefulProperties ( )
overrideprotectedvirtualinherited

Reimplemented from ComputeStressBase.

Definition at line 61 of file ComputePFFractureStressBase.C.

62 {
63  _H[_qp] = 0.0;
64 }

◆ validParams()

InputParameters ComputeLinearElasticPFFractureStress::validParams ( )
static

Definition at line 18 of file ComputeLinearElasticPFFractureStress.C.

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 }

Member Data Documentation

◆ _barrier

const MaterialProperty<Real>& ComputePFFractureStressBase::_barrier
protectedinherited

material property for fracture energy barrier

Definition at line 56 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 45 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 38 of file ComputePFFractureStressBase.h.

◆ _D

const MaterialProperty<Real>& ComputePFFractureStressBase::_D
protectedinherited

Material property for energetic degradation function.

Definition at line 74 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 80 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 65 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 71 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 77 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 62 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 68 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 59 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 35 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 33 of file ComputePFFractureStressBase.h.

Referenced by initialSetup().

◆ _extra_stress

const MaterialProperty<RankTwoTensor>& ComputeStressBase::_extra_stress
protectedinherited

Extra stress tensor.

Definition at line 55 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 44 of file ComputePFFractureStressBase.h.

◆ _gc_block_restrict

BlockRestrictable* const GuaranteeConsumer::_gc_block_restrict
privateinherited

Access block restrictions of the object with this interface.

Definition at line 41 of file GuaranteeConsumer.h.

Referenced by GuaranteeConsumer::hasGuaranteedMaterialProperty().

◆ _gc_feproblem

FEProblemBase* const GuaranteeConsumer::_gc_feproblem
privateinherited

Reference to the FEProblemBase class.

Definition at line 38 of file GuaranteeConsumer.h.

Referenced by GuaranteeConsumer::hasGuaranteedMaterialProperty().

◆ _gc_params

const InputParameters& GuaranteeConsumer::_gc_params
privateinherited

Parameters of the object with this interface.

Definition at line 35 of file GuaranteeConsumer.h.

◆ _H

MaterialProperty<Real>& ComputePFFractureStressBase::_H
protectedinherited

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

Definition at line 50 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 53 of file ComputePFFractureStressBase.h.

Referenced by computeQpStress().

◆ _initial_stress_fcn

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

initial stress components

Definition at line 58 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 41 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 50 of file ComputeStressBase.h.

Referenced by ComputeMultipleInelasticCosseratStress::computeAdmissibleState(), ComputeMultipleInelasticStress::computeAdmissibleState(), ComputeStressBase::computeQpProperties(), ComputeStrainIncrementBasedStress::computeQpStress(), ComputeLinearElasticStress::computeQpStress(), ComputeDamageStress::computeQpStress(), ComputeFiniteStrainElasticStress::computeQpStress(), ComputeCosseratLinearElasticStress::computeQpStress(), ComputeSmearedCrackingStress::computeQpStress(), computeQpStress(), FiniteStrainPlasticMaterial::computeQpStress(), ComputeMultiPlasticityStress::computeQpStress(), ComputeLinearViscoelasticStress::computeQpStress(), ComputeMultipleInelasticStress::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 47 of file ComputePFFractureStressBase.h.

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


The documentation for this class was generated from the following files:
ComputePFFractureStressBase::ComputePFFractureStressBase
ComputePFFractureStressBase(const InputParameters &parameters)
Definition: ComputePFFractureStressBase.C:34
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
ComputeStressBase::_extra_stress
const MaterialProperty< RankTwoTensor > & _extra_stress
Extra stress tensor.
Definition: ComputeStressBase.h:55
ComputeStressBase::_Jacobian_mult
MaterialProperty< RankFourTensor > & _Jacobian_mult
derivative of stress w.r.t. strain (_dstress_dstrain)
Definition: ComputeStressBase.h:61
GuaranteeConsumer::GuaranteeConsumer
GuaranteeConsumer(MooseObject *moose_object)
Definition: GuaranteeConsumer.C:20
ComputeLinearElasticPFFractureStress::Decomposition_type::strain_vol_dev
ComputePFFractureStressBase::_barrier
const MaterialProperty< Real > & _barrier
material property for fracture energy barrier
Definition: ComputePFFractureStressBase.h:56
GuaranteeConsumer::_gc_block_restrict
BlockRestrictable *const _gc_block_restrict
Access block restrictions of the object with this interface.
Definition: GuaranteeConsumer.h:41
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
ComputeStressBase::computeQpStress
virtual void computeQpStress()=0
Compute the stress and store it in the _stress material property for the current quadrature point.
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
GuaranteeConsumer::_gc_feproblem
FEProblemBase *const _gc_feproblem
Reference to the FEProblemBase class.
Definition: GuaranteeConsumer.h:38
ComputeLinearElasticPFFractureStress::Decomposition_type::strain_spectral
ComputeLinearElasticPFFractureStress::Decomposition_type::stress_spectral
ComputePFFractureStressBase::_H_old
const MaterialProperty< Real > & _H_old
Old value of history variable.
Definition: ComputePFFractureStressBase.h:53
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
RankFourTensorTempl
Definition: ACGrGrElasticDrivingForce.h:20
ComputePFFractureStressBase::_D
const MaterialProperty< Real > & _D
Material property for energetic degradation function.
Definition: ComputePFFractureStressBase.h:74
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