https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ComputePlaneSmallStrainNOSPD Class Reference

Material class for 2D correspondence material model for small strain: plane strain, generalized plane strain, weak plane stress. More...

#include <ComputePlaneSmallStrainNOSPD.h>

Inheritance diagram for ComputePlaneSmallStrainNOSPD:
[legend]

Public Types

typedef DerivativeMaterialPropertyNameInterface::SymbolName SymbolName
 

Public Member Functions

 ComputePlaneSmallStrainNOSPD (const InputParameters &parameters)
 
virtual void initQpStatefulProperties () override
 
const GenericMaterialProperty< U, is_ad > & getDefaultMaterialProperty (const std::string &name)
 
const GenericMaterialProperty< U, is_ad > & getDefaultMaterialPropertyByName (const std::string &name)
 
void validateDerivativeMaterialPropertyBase (const std::string &base)
 
const MaterialPropertyName derivativePropertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName derivativePropertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName derivativePropertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName derivativePropertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, unsigned int v2, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, unsigned int v1, unsigned int v2=libMesh::invalid_uint, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, unsigned int v2, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, unsigned int v1, unsigned int v2=libMesh::invalid_uint, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
void validateCoupling (const MaterialPropertyName &base, const std::vector< VariableName > &c, bool validate_aux=true)
 
void validateCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateCoupling (const MaterialPropertyName &base, const std::vector< VariableName > &c, bool validate_aux=true)
 
void validateCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateNonlinearCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateNonlinearCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
const MaterialPropertyName propertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName propertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName propertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName propertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName propertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName propertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName propertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 
const MaterialPropertyName propertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

virtual void computeQpTotalStrain () override
 Function to compute the total strain tensor for small strain case. More...
 
Real computeOutOfPlaneStrain ()
 Function to compute out-of-plane component of strain tensor for generalized plane strain and weak plane stress. More...
 
virtual void computeQpStrain () override
 Function to compute strain tensors. More...
 
virtual void computeProperties () override
 
virtual void computeBondStretch () override
 
virtual void computeQpDeformationGradient ()
 Function to compute deformation gradient for peridynamic correspondence model. More...
 
virtual void computeConventionalQpDeformationGradient ()
 Function to compute conventional nonlocal deformation gradient. More...
 
virtual void computeBondHorizonQpDeformationGradient ()
 Function to compute bond-associated horizon based deformation gradient. More...
 

Protected Attributes

const MooseEnum _stabilization
 Option of stabilization scheme for correspondence material model: FORCE, BOND_HORIZON_I or BOND_HORIZON_II. More...
 
const bool _plane_strain
 Plane strain problem or not, this is only used for mechanical stretch calculation. More...
 
const MaterialProperty< RankFourTensor > & _Cijkl
 Material properties to fetch. More...
 
std::vector< MaterialPropertyName > _eigenstrain_names
 
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
 
MaterialProperty< RankTwoTensor > & _shape2
 Material properties to store. More...
 
MaterialProperty< RankTwoTensor > & _deformation_gradient
 
MaterialProperty< RankTwoTensor > & _ddgraddu
 
MaterialProperty< RankTwoTensor > & _ddgraddv
 
MaterialProperty< RankTwoTensor > & _ddgraddw
 
MaterialProperty< RankTwoTensor > & _total_strain
 
MaterialProperty< RankTwoTensor > & _mechanical_strain
 
MaterialProperty< Real > & _multi
 
MaterialProperty< Real > & _sf_coeff
 fictitious force coefficient for force stabilized model More...
 

Private Attributes

const bool _scalar_out_of_plane_strain_coupled
 Scalar out-of-plane strain for generalized plane strain. More...
 
const VariableValue_scalar_out_of_plane_strain
 
const bool _out_of_plane_strain_coupled
 Out-of-plane strain for weak plane stress. More...
 
const VariableValue_out_of_plane_strain
 

Detailed Description

Material class for 2D correspondence material model for small strain: plane strain, generalized plane strain, weak plane stress.

Definition at line 19 of file ComputePlaneSmallStrainNOSPD.h.

Constructor & Destructor Documentation

◆ ComputePlaneSmallStrainNOSPD()

ComputePlaneSmallStrainNOSPD::ComputePlaneSmallStrainNOSPD ( const InputParameters parameters)

Definition at line 30 of file ComputePlaneSmallStrainNOSPD.C.

31  : ComputeSmallStrainNOSPD(parameters),
32  _scalar_out_of_plane_strain_coupled(isCoupledScalar("scalar_out_of_plane_strain")),
34  ? coupledScalarValue("scalar_out_of_plane_strain")
35  : _zero),
36  _out_of_plane_strain_coupled(isCoupled("out_of_plane_strain")),
37  _out_of_plane_strain(_out_of_plane_strain_coupled ? coupledValue("out_of_plane_strain") : _zero)
38 {
39 }
const VariableValue & _scalar_out_of_plane_strain
const bool _scalar_out_of_plane_strain_coupled
Scalar out-of-plane strain for generalized plane strain.
ComputeSmallStrainNOSPD(const InputParameters &parameters)
const bool _out_of_plane_strain_coupled
Out-of-plane strain for weak plane stress.

Member Function Documentation

◆ computeBondHorizonQpDeformationGradient()

void ComputeStrainBaseNOSPD::computeBondHorizonQpDeformationGradient ( )
protectedvirtualinherited

Function to compute bond-associated horizon based deformation gradient.

Definition at line 163 of file ComputeStrainBaseNOSPD.C.

Referenced by ComputeStrainBaseNOSPD::computeQpDeformationGradient().

164 {
165  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(_qp));
166  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(_qp));
167 
168  dof_id_type nb_index =
169  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - _qp)) -
170  neighbors.begin();
171  std::vector<dof_id_type> dg_neighbors =
172  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(_qp), nb_index);
173  Real dg_sub_vol = _pdmesh.getHorizonSubsetVolume(_current_elem->node_id(_qp), nb_index);
174  Real dg_sub_vol_sum = _pdmesh.getHorizonSubsetVolumeSum(_current_elem->node_id(_qp));
175 
176  // calculate the shape tensor and prepare the deformation gradient tensor
177  Real vol_nb, weight_nb;
178  RealGradient origin_vec_nb, current_vec_nb;
179 
180  for (unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
181  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
182  {
183  vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
184  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
185  _pdmesh.getNodeCoord(_current_elem->node_id(_qp));
186 
187  for (unsigned int k = 0; k < _dim; ++k)
188  current_vec_nb(k) =
189  origin_vec_nb(k) +
190  _disp_var[k]->getNodalValue(*_pdmesh.nodePtr(neighbors[dg_neighbors[nb]])) -
191  _disp_var[k]->getNodalValue(*_current_elem->node_ptr(_qp));
192 
193  weight_nb = _horizon_radius[_qp] / origin_vec_nb.norm();
194  for (unsigned int k = 0; k < _dim; ++k)
195  {
196  for (unsigned int l = 0; l < _dim; ++l)
197  {
198  _shape2[_qp](k, l) += weight_nb * origin_vec_nb(k) * origin_vec_nb(l) * vol_nb;
199  _deformation_gradient[_qp](k, l) +=
200  weight_nb * current_vec_nb(k) * origin_vec_nb(l) * vol_nb;
201  }
202  // calculate derivatives of deformation_gradient w.r.t displacements of node i
203  _ddgraddu[_qp](0, k) += -weight_nb * origin_vec_nb(k) * vol_nb;
204  _ddgraddv[_qp](1, k) += -weight_nb * origin_vec_nb(k) * vol_nb;
205  if (_dim == 3)
206  _ddgraddw[_qp](2, k) += -weight_nb * origin_vec_nb(k) * vol_nb;
207  }
208  }
209  // finalize the deformation gradient and its derivatives
210  if (_shape2[_qp].det() == 0.)
211  computeConventionalQpDeformationGradient(); // this will affect the corresponding jacobian
212  // calculation
213  else
214  {
215  _deformation_gradient[_qp] *= _shape2[_qp].inverse();
216  _ddgraddu[_qp] *= _shape2[_qp].inverse();
217  _ddgraddv[_qp] *= _shape2[_qp].inverse();
218  _ddgraddw[_qp] *= _shape2[_qp].inverse();
219  }
220 
221  // force state multiplier
222  if (_stabilization == "BOND_HORIZON_I")
223  _multi[_qp] = _horizon_radius[_qp] / _origin_vec.norm() * _node_vol[0] * _node_vol[1] *
224  dg_sub_vol / _horizon_vol[_qp];
225  else if (_stabilization == "BOND_HORIZON_II")
226  _multi[_qp] = _node_vol[_qp] * dg_sub_vol / dg_sub_vol_sum;
227 }
MaterialProperty< RankTwoTensor > & _shape2
Material properties to store.
auto norm() const -> decltype(std::norm(Real()))
MaterialProperty< RankTwoTensor > & _deformation_gradient
const MooseEnum _stabilization
Option of stabilization scheme for correspondence material model: FORCE, BOND_HORIZON_I or BOND_HORIZ...
MaterialProperty< RankTwoTensor > & _ddgraddv
MaterialProperty< RankTwoTensor > & _ddgraddw
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void computeConventionalQpDeformationGradient()
Function to compute conventional nonlocal deformation gradient.
MaterialProperty< Real > & _multi
static const std::string k
Definition: NS.h:130
MaterialProperty< RankTwoTensor > & _ddgraddu
uint8_t dof_id_type

◆ computeBondStretch()

void ComputeStrainBaseNOSPD::computeBondStretch ( )
overrideprotectedvirtualinherited

Definition at line 241 of file ComputeStrainBaseNOSPD.C.

Referenced by ComputeStrainBaseNOSPD::computeProperties().

242 {
243  _total_stretch[0] = _current_len / _origin_vec.norm() - 1.0;
244  _total_stretch[1] = _total_stretch[0];
245  _mechanical_stretch[0] = _total_stretch[0];
246 
247  Real factor = 1.0;
248  if (_plane_strain)
250 
251  for (auto es : _eigenstrains)
252  _mechanical_stretch[0] -= 0.5 * factor * ((*es)[0](2, 2) + (*es)[1](2, 2));
253 
254  _mechanical_stretch[1] = _mechanical_stretch[0];
255 }
const bool _plane_strain
Plane strain problem or not, this is only used for mechanical stretch calculation.
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
T getIsotropicPoissonsRatio(const RankFourTensorTempl< T > &elasticity_tensor)
Get the Poisson&#39;s modulus for an isotropic elasticity tensor param elasticity_tensor the tensor (must...
const MaterialProperty< RankFourTensor > & _Cijkl
Material properties to fetch.

◆ computeConventionalQpDeformationGradient()

void ComputeStrainBaseNOSPD::computeConventionalQpDeformationGradient ( )
protectedvirtualinherited

Function to compute conventional nonlocal deformation gradient.

Definition at line 113 of file ComputeStrainBaseNOSPD.C.

Referenced by ComputeStrainBaseNOSPD::computeBondHorizonQpDeformationGradient(), and ComputeStrainBaseNOSPD::computeQpDeformationGradient().

114 {
115  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(_qp));
116  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(_qp));
117 
118  // calculate the shape tensor and prepare the deformation gradient tensor
119  Real vol_nb, weight_nb;
120  RealGradient origin_vec_nb, current_vec_nb;
121 
122  for (unsigned int nb = 0; nb < neighbors.size(); ++nb)
123  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[nb])) > 0.5)
124  {
125  vol_nb = _pdmesh.getNodeVolume(neighbors[nb]);
126  origin_vec_nb =
127  _pdmesh.getNodeCoord(neighbors[nb]) - _pdmesh.getNodeCoord(_current_elem->node_id(_qp));
128 
129  for (unsigned int k = 0; k < _dim; ++k)
130  current_vec_nb(k) = origin_vec_nb(k) +
131  _disp_var[k]->getNodalValue(*_pdmesh.nodePtr(neighbors[nb])) -
132  _disp_var[k]->getNodalValue(*_current_elem->node_ptr(_qp));
133 
134  weight_nb = _horizon_radius[_qp] / origin_vec_nb.norm();
135  for (unsigned int k = 0; k < _dim; ++k)
136  {
137  for (unsigned int l = 0; l < _dim; ++l)
138  {
139  _shape2[_qp](k, l) += weight_nb * origin_vec_nb(k) * origin_vec_nb(l) * vol_nb;
140  _deformation_gradient[_qp](k, l) +=
141  weight_nb * current_vec_nb(k) * origin_vec_nb(l) * vol_nb;
142  }
143  // calculate derivatives of deformation_gradient w.r.t displacements of node nb
144  _ddgraddu[_qp](0, k) += -weight_nb * origin_vec_nb(k) * vol_nb;
145  _ddgraddv[_qp](1, k) += -weight_nb * origin_vec_nb(k) * vol_nb;
146  if (_dim == 3)
147  _ddgraddw[_qp](2, k) += -weight_nb * origin_vec_nb(k) * vol_nb;
148  }
149  }
150 
151  // finalize the deformation gradient tensor
152  if (_shape2[_qp].det() == 0.)
153  mooseError("Singular shape tensor is detected in ComputeStrainBaseNOSPD! Use "
154  "SingularShapeTensorEliminatorUserObjectPD to avoid singular shape tensor!");
155 
156  _deformation_gradient[_qp] *= _shape2[_qp].inverse();
157  _ddgraddu[_qp] *= _shape2[_qp].inverse();
158  _ddgraddv[_qp] *= _shape2[_qp].inverse();
159  _ddgraddw[_qp] *= _shape2[_qp].inverse();
160 }
MaterialProperty< RankTwoTensor > & _shape2
Material properties to store.
auto norm() const -> decltype(std::norm(Real()))
MaterialProperty< RankTwoTensor > & _deformation_gradient
void mooseError(Args &&... args)
MaterialProperty< RankTwoTensor > & _ddgraddv
MaterialProperty< RankTwoTensor > & _ddgraddw
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string k
Definition: NS.h:130
MaterialProperty< RankTwoTensor > & _ddgraddu

◆ computeOutOfPlaneStrain()

Real ComputePlaneSmallStrainNOSPD::computeOutOfPlaneStrain ( )
protected

Function to compute out-of-plane component of strain tensor for generalized plane strain and weak plane stress.

Returns
The value of out-of-plane strain

Definition at line 51 of file ComputePlaneSmallStrainNOSPD.C.

Referenced by computeQpTotalStrain().

52 {
55  else
56  return _out_of_plane_strain[_qp];
57 }
const VariableValue & _scalar_out_of_plane_strain
const bool _scalar_out_of_plane_strain_coupled
Scalar out-of-plane strain for generalized plane strain.

◆ computeProperties()

void ComputeStrainBaseNOSPD::computeProperties ( )
overrideprotectedvirtualinherited

Definition at line 230 of file ComputeStrainBaseNOSPD.C.

231 {
232  setupMeshRelatedData(); // function from base class
233  computeBondCurrentLength(); // current length of a bond from base class
234  computeBondStretch(); // stretch of a bond
235 
236  for (_qp = 0; _qp < _nnodes; ++_qp)
237  computeQpStrain();
238 }
virtual void computeBondStretch() override
virtual void computeQpStrain()=0
Function to compute strain tensors.

◆ computeQpDeformationGradient()

void ComputeStrainBaseNOSPD::computeQpDeformationGradient ( )
protectedvirtualinherited

Function to compute deformation gradient for peridynamic correspondence model.

Definition at line 75 of file ComputeStrainBaseNOSPD.C.

Referenced by ComputeSmallStrainNOSPD::computeQpStrain(), and ComputeFiniteStrainNOSPD::computeQpStrain().

76 {
77  _shape2[_qp].zero();
78  _deformation_gradient[_qp].zero();
79  _ddgraddu[_qp].zero();
80  _ddgraddv[_qp].zero();
81  _ddgraddw[_qp].zero();
82  _multi[_qp] = 0.0;
83  _sf_coeff[_qp] = 0.0;
84 
85  if (_dim == 2)
86  _shape2[_qp](2, 2) = _deformation_gradient[_qp](2, 2) = 1.0;
87 
88  if (_bond_status_var->getElementalValue(_current_elem) > 0.5)
89  {
90  if (_stabilization == "FORCE")
91  {
93 
95  _pdmesh.getNodeAverageSpacing(_current_elem->node_id(_qp)) *
96  _horizon_radius[_qp] / _origin_vec.norm();
97  _multi[_qp] = _horizon_radius[_qp] / _origin_vec.norm() * _node_vol[0] * _node_vol[1];
98  }
99  else if (_stabilization == "BOND_HORIZON_I" || _stabilization == "BOND_HORIZON_II")
101  else
102  paramError("stabilization",
103  "Unknown stabilization scheme for peridynamic correspondence model!");
104  }
105  else
106  {
107  _shape2[_qp].setToIdentity();
108  _deformation_gradient[_qp].setToIdentity();
109  }
110 }
MaterialProperty< RankTwoTensor > & _shape2
Material properties to store.
MaterialProperty< RankTwoTensor > & _deformation_gradient
const MooseEnum _stabilization
Option of stabilization scheme for correspondence material model: FORCE, BOND_HORIZON_I or BOND_HORIZ...
MaterialProperty< RankTwoTensor > & _ddgraddv
MaterialProperty< RankTwoTensor > & _ddgraddw
virtual void computeBondHorizonQpDeformationGradient()
Function to compute bond-associated horizon based deformation gradient.
MaterialProperty< Real > & _sf_coeff
fictitious force coefficient for force stabilized model
virtual void computeConventionalQpDeformationGradient()
Function to compute conventional nonlocal deformation gradient.
MaterialProperty< Real > & _multi
const MaterialProperty< RankFourTensor > & _Cijkl
Material properties to fetch.
T getIsotropicYoungsModulus(const RankFourTensorTempl< T > &elasticity_tensor)
Get the Young&#39;s modulus for an isotropic elasticity tensor param elasticity_tensor the tensor (must b...
MaterialProperty< RankTwoTensor > & _ddgraddu

◆ computeQpStrain()

void ComputeSmallStrainNOSPD::computeQpStrain ( )
overrideprotectedvirtualinherited

Function to compute strain tensors.

Implements ComputeStrainBaseNOSPD.

Definition at line 31 of file ComputeSmallStrainNOSPD.C.

32 {
34 
36 
38 
39  // Subtract Eigen strains
40  for (auto es : _eigenstrains)
41  _mechanical_strain[_qp] -= (*es)[_qp];
42 
43  // zero out all strain measures for broken bond
44  if (_bond_status_var->getElementalValue(_current_elem) < 0.5)
45  {
46  _mechanical_strain[_qp].zero();
47  _total_strain[_qp].zero();
48  }
49 }
virtual void computeQpDeformationGradient()
Function to compute deformation gradient for peridynamic correspondence model.
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
virtual void computeQpTotalStrain()
Function to compute the total strain tensor for small strain case.
MaterialProperty< RankTwoTensor > & _total_strain
MaterialProperty< RankTwoTensor > & _mechanical_strain

◆ computeQpTotalStrain()

void ComputePlaneSmallStrainNOSPD::computeQpTotalStrain ( )
overrideprotectedvirtual

Function to compute the total strain tensor for small strain case.

Reimplemented from ComputeSmallStrainNOSPD.

Definition at line 42 of file ComputePlaneSmallStrainNOSPD.C.

43 {
44  // the green-lagrange strain tensor
45  _total_strain[_qp] = 0.5 * (_deformation_gradient[_qp].transpose() + _deformation_gradient[_qp]) -
48 }
MaterialProperty< RankTwoTensor > & _deformation_gradient
Real computeOutOfPlaneStrain()
Function to compute out-of-plane component of strain tensor for generalized plane strain and weak pla...
MaterialProperty< RankTwoTensor > & _total_strain

◆ initQpStatefulProperties()

void ComputeStrainBaseNOSPD::initQpStatefulProperties ( )
overridevirtualinherited

Definition at line 62 of file ComputeStrainBaseNOSPD.C.

63 {
64  _mechanical_strain[_qp].zero();
65  _total_strain[_qp].zero();
66  _deformation_gradient[_qp].setToIdentity();
67 
68  if (_qrule->n_points() < 2) // Gauss_Lobatto: order = 2n-3 (n is number of qp)
69  mooseError(
70  "NOSPD models require Gauss_Lobatto rule and a minimum of 2 quadrature points, i.e., "
71  "order of FIRST");
72 }
MaterialProperty< RankTwoTensor > & _deformation_gradient
void mooseError(Args &&... args)
MaterialProperty< RankTwoTensor > & _total_strain
MaterialProperty< RankTwoTensor > & _mechanical_strain

◆ validParams()

InputParameters ComputePlaneSmallStrainNOSPD::validParams ( )
static

Definition at line 15 of file ComputePlaneSmallStrainNOSPD.C.

16 {
18  params.addClassDescription(
19  "Class for computing nodal quantities for residual and jacobian calculation "
20  "for peridynamic correspondence models under planar small strain assumptions");
21 
22  params.addCoupledVar("scalar_out_of_plane_strain",
23  "Scalar out-of-plane strain variable for generalized plane strain");
24  params.addCoupledVar("out_of_plane_strain",
25  "Nonlinear out-of-plane strain variable for plane stress condition");
26 
27  return params;
28 }
static InputParameters validParams()
void addCoupledVar(const std::string &name, const std::string &doc_string)
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _Cijkl

const MaterialProperty<RankFourTensor>& ComputeStrainBaseNOSPD::_Cijkl
protectedinherited

◆ _ddgraddu

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_ddgraddu
protectedinherited

◆ _ddgraddv

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_ddgraddv
protectedinherited

◆ _ddgraddw

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_ddgraddw
protectedinherited

◆ _deformation_gradient

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_deformation_gradient
protectedinherited

◆ _eigenstrain_names

std::vector<MaterialPropertyName> ComputeStrainBaseNOSPD::_eigenstrain_names
protectedinherited

◆ _eigenstrains

std::vector<const MaterialProperty<RankTwoTensor> *> ComputeStrainBaseNOSPD::_eigenstrains
protectedinherited

◆ _mechanical_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_mechanical_strain
protectedinherited

◆ _multi

MaterialProperty<Real>& ComputeStrainBaseNOSPD::_multi
protectedinherited

◆ _out_of_plane_strain

const VariableValue& ComputePlaneSmallStrainNOSPD::_out_of_plane_strain
private

Definition at line 44 of file ComputePlaneSmallStrainNOSPD.h.

Referenced by computeOutOfPlaneStrain().

◆ _out_of_plane_strain_coupled

const bool ComputePlaneSmallStrainNOSPD::_out_of_plane_strain_coupled
private

Out-of-plane strain for weak plane stress.

Definition at line 43 of file ComputePlaneSmallStrainNOSPD.h.

◆ _plane_strain

const bool ComputeStrainBaseNOSPD::_plane_strain
protectedinherited

Plane strain problem or not, this is only used for mechanical stretch calculation.

Definition at line 56 of file ComputeStrainBaseNOSPD.h.

Referenced by ComputeStrainBaseNOSPD::computeBondStretch().

◆ _scalar_out_of_plane_strain

const VariableValue& ComputePlaneSmallStrainNOSPD::_scalar_out_of_plane_strain
private

Definition at line 39 of file ComputePlaneSmallStrainNOSPD.h.

Referenced by computeOutOfPlaneStrain().

◆ _scalar_out_of_plane_strain_coupled

const bool ComputePlaneSmallStrainNOSPD::_scalar_out_of_plane_strain_coupled
private

Scalar out-of-plane strain for generalized plane strain.

Definition at line 38 of file ComputePlaneSmallStrainNOSPD.h.

Referenced by computeOutOfPlaneStrain().

◆ _sf_coeff

MaterialProperty<Real>& ComputeStrainBaseNOSPD::_sf_coeff
protectedinherited

fictitious force coefficient for force stabilized model

Definition at line 78 of file ComputeStrainBaseNOSPD.h.

Referenced by ComputeStrainBaseNOSPD::computeQpDeformationGradient().

◆ _shape2

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_shape2
protectedinherited

◆ _stabilization

const MooseEnum ComputeStrainBaseNOSPD::_stabilization
protectedinherited

Option of stabilization scheme for correspondence material model: FORCE, BOND_HORIZON_I or BOND_HORIZON_II.

Definition at line 53 of file ComputeStrainBaseNOSPD.h.

Referenced by ComputeStrainBaseNOSPD::computeBondHorizonQpDeformationGradient(), and ComputeStrainBaseNOSPD::computeQpDeformationGradient().

◆ _total_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_total_strain
protectedinherited

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