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

Material class for peridynamic correspondence model for finite strain. More...

#include <ComputeFiniteStrainNOSPD.h>

Inheritance diagram for ComputeFiniteStrainNOSPD:
[legend]

Public Types

typedef DerivativeMaterialPropertyNameInterface::SymbolName SymbolName
 

Public Member Functions

 ComputeFiniteStrainNOSPD (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 ()
 
static MooseEnum decompositionType ()
 

Protected Member Functions

virtual void computeQpStrain () override
 Function to compute strain tensors. More...
 
virtual void computeQpFhat ()
 
void computeQpStrainRotationIncrements (RankTwoTensor &e, RankTwoTensor &r)
 Function to compute strain and rotational increments. More...
 
void subtractEigenstrainIncrementFromStrain (RankTwoTensor &strain)
 Function to compute the mechanical strain tensor by subtracting thermal strain from the total strain. 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

std::vector< RankTwoTensor_Fhat
 'Incremental' deformation gradient More...
 
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...
 
MaterialProperty< RankTwoTensor > & _strain_rate
 Material properties to store. More...
 
MaterialProperty< RankTwoTensor > & _strain_increment
 
MaterialProperty< RankTwoTensor > & _rotation_increment
 
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
 Material properties to fetch. More...
 
const MaterialProperty< RankTwoTensor > & _mechanical_strain_old
 
const MaterialProperty< RankTwoTensor > & _total_strain_old
 
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains_old
 
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...
 
enum  DecompMethod { DecompMethod::TaylorExpansion, DecompMethod::EigenSolution }
 Method to decompose into rotation increment and strain increment. More...
 
const DecompMethod _decomposition_method
 

Detailed Description

Material class for peridynamic correspondence model for finite strain.

Definition at line 18 of file ComputeFiniteStrainNOSPD.h.

Member Enumeration Documentation

◆ DecompMethod

Method to decompose into rotation increment and strain increment.

Enumerator
TaylorExpansion 
EigenSolution 

Definition at line 60 of file ComputeFiniteStrainNOSPD.h.

61  {
62  TaylorExpansion,
63  EigenSolution
64  };

Constructor & Destructor Documentation

◆ ComputeFiniteStrainNOSPD()

ComputeFiniteStrainNOSPD::ComputeFiniteStrainNOSPD ( const InputParameters parameters)

Definition at line 37 of file ComputeFiniteStrainNOSPD.C.

38  : ComputeStrainBaseNOSPD(parameters),
39  _strain_rate(declareProperty<RankTwoTensor>("strain_rate")),
40  _strain_increment(declareProperty<RankTwoTensor>("strain_increment")),
41  _rotation_increment(declareProperty<RankTwoTensor>("rotation_increment")),
42  _deformation_gradient_old(getMaterialPropertyOld<RankTwoTensor>("deformation_gradient")),
43  _mechanical_strain_old(getMaterialPropertyOld<RankTwoTensor>("mechanical_strain")),
44  _total_strain_old(getMaterialPropertyOld<RankTwoTensor>("total_strain")),
46  _Fhat(2),
47  _decomposition_method(getParam<MooseEnum>("decomposition_method").getEnum<DecompMethod>())
48 {
49  for (unsigned int i = 0; i < _eigenstrains_old.size(); ++i)
50  _eigenstrains_old[i] = &getMaterialPropertyOld<RankTwoTensor>(_eigenstrain_names[i]);
51 }
const MaterialProperty< RankTwoTensor > & _total_strain_old
MaterialProperty< RankTwoTensor > & _strain_increment
const DecompMethod _decomposition_method
std::vector< RankTwoTensor > _Fhat
&#39;Incremental&#39; deformation gradient
ComputeStrainBaseNOSPD(const InputParameters &parameters)
const MaterialProperty< RankTwoTensor > & _mechanical_strain_old
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
Material properties to fetch.
std::vector< MaterialPropertyName > _eigenstrain_names
MaterialProperty< RankTwoTensor > & _strain_rate
Material properties to store.
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains_old
MaterialProperty< RankTwoTensor > & _rotation_increment

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

◆ 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 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

◆ computeQpFhat()

void ComputeFiniteStrainNOSPD::computeQpFhat ( )
protectedvirtual

Reimplemented in ComputePlaneFiniteStrainNOSPD.

Definition at line 97 of file ComputeFiniteStrainNOSPD.C.

Referenced by computeQpStrain().

98 {
99  // Incremental deformation gradient of current step w.r.t previous step:
100  // _Fhat = deformation_gradient * inv(deformation_gradient_old)
101  _Fhat[_qp] = _deformation_gradient[_qp] * _deformation_gradient_old[_qp].inverse();
102 }
MaterialProperty< RankTwoTensor > & _deformation_gradient
std::vector< RankTwoTensor > _Fhat
&#39;Incremental&#39; deformation gradient
const MaterialProperty< RankTwoTensor > & _deformation_gradient_old
Material properties to fetch.

◆ computeQpStrain()

void ComputeFiniteStrainNOSPD::computeQpStrain ( )
overrideprotectedvirtual

Function to compute strain tensors.

Implements ComputeStrainBaseNOSPD.

Definition at line 54 of file ComputeFiniteStrainNOSPD.C.

55 {
57 
58  computeQpFhat();
59 
60  RankTwoTensor total_strain_increment;
61 
62  // Two ways to calculate these increments: TaylorExpansion(default) or EigenSolution
63  computeQpStrainRotationIncrements(total_strain_increment, _rotation_increment[_qp]);
64 
65  _strain_increment[_qp] = total_strain_increment;
66 
67  // Remove the eigenstrain increment
69 
70  if (_dt > 0)
71  _strain_rate[_qp] = _strain_increment[_qp] / _dt;
72  else
73  _strain_rate[_qp].zero();
74 
75  // Update strain in intermediate configuration
77  _total_strain[_qp] = _total_strain_old[_qp] + total_strain_increment;
78 
79  // Rotate strain to current configuration
80  _mechanical_strain[_qp] =
81  _rotation_increment[_qp] * _mechanical_strain[_qp] * _rotation_increment[_qp].transpose();
82  _total_strain[_qp] =
83  _rotation_increment[_qp] * _total_strain[_qp] * _rotation_increment[_qp].transpose();
84 
85  // zero out all strain measures for broken bond
86  if (_bond_status_var->getElementalValue(_current_elem) < 0.5)
87  {
88  _strain_rate[_qp].zero();
89  _strain_increment[_qp].zero();
90  _rotation_increment[_qp].zero();
91  _mechanical_strain[_qp].zero();
92  _total_strain[_qp].zero();
93  }
94 }
void computeQpStrainRotationIncrements(RankTwoTensor &e, RankTwoTensor &r)
Function to compute strain and rotational increments.
virtual void computeQpDeformationGradient()
Function to compute deformation gradient for peridynamic correspondence model.
const MaterialProperty< RankTwoTensor > & _total_strain_old
MaterialProperty< RankTwoTensor > & _strain_increment
void subtractEigenstrainIncrementFromStrain(RankTwoTensor &strain)
Function to compute the mechanical strain tensor by subtracting thermal strain from the total strain...
const MaterialProperty< RankTwoTensor > & _mechanical_strain_old
MaterialProperty< RankTwoTensor > & _strain_rate
Material properties to store.
MaterialProperty< RankTwoTensor > & _rotation_increment
MaterialProperty< RankTwoTensor > & _total_strain
MaterialProperty< RankTwoTensor > & _mechanical_strain

◆ computeQpStrainRotationIncrements()

void ComputeFiniteStrainNOSPD::computeQpStrainRotationIncrements ( RankTwoTensor e,
RankTwoTensor r 
)
protected

Function to compute strain and rotational increments.

Definition at line 105 of file ComputeFiniteStrainNOSPD.C.

Referenced by computeQpStrain().

107 {
108  switch (_decomposition_method)
109  {
111  {
112  // inverse of _Fhat
113  RankTwoTensor invFhat(_Fhat[_qp].inverse());
114 
115  // A = I - _Fhat^-1
117  A -= invFhat;
118 
119  // Cinv - I = A A^T - A - A^T;
120  RankTwoTensor Cinv_I = A * A.transpose() - A - A.transpose();
121 
122  // strain rate D from Taylor expansion, Chat = (-1/2(Chat^-1 - I) + 1/4*(Chat^-1 - I)^2 + ...
123  total_strain_increment = -Cinv_I * 0.5 + Cinv_I * Cinv_I * 0.25;
124 
125  const Real a[3] = {invFhat(1, 2) - invFhat(2, 1),
126  invFhat(2, 0) - invFhat(0, 2),
127  invFhat(0, 1) - invFhat(1, 0)};
128 
129  Real q = (a[0] * a[0] + a[1] * a[1] + a[2] * a[2]) / 4.0;
130  Real trFhatinv_1 = invFhat.trace() - 1.0;
131  const Real p = trFhatinv_1 * trFhatinv_1 / 4.0;
132 
133  // cos theta_a
134  const Real C1_squared = p +
135  3.0 * Utility::pow<2>(p) * (1.0 - (p + q)) / Utility::pow<2>(p + q) -
136  2.0 * Utility::pow<3>(p) * (1.0 - (p + q)) / Utility::pow<3>(p + q);
137  if (C1_squared <= 0.0)
138  mooseException(
139  "Cannot take square root of a number less than or equal to zero in the calculation of "
140  "C1 for the Rashid approximation for the rotation tensor.");
141 
142  const Real C1 = std::sqrt(C1_squared);
143 
144  Real C2;
145  if (q > 0.01)
146  // (1-cos theta_a)/4q
147  C2 = (1.0 - C1) / (4.0 * q);
148  else
149  // alternate form for small q
150  C2 = 0.125 + q * 0.03125 * (Utility::pow<2>(p) - 12.0 * (p - 1.0)) / Utility::pow<2>(p) +
151  Utility::pow<2>(q) * (p - 2.0) * (Utility::pow<2>(p) - 10.0 * p + 32.0) /
152  Utility::pow<3>(p) +
153  Utility::pow<3>(q) *
154  (1104.0 - 992.0 * p + 376.0 * Utility::pow<2>(p) - 72.0 * Utility::pow<3>(p) +
155  5.0 * Utility::pow<4>(p)) /
156  (512.0 * Utility::pow<4>(p));
157 
158  const Real C3_test =
159  (p * q * (3.0 - q) + Utility::pow<3>(p) + Utility::pow<2>(q)) / Utility::pow<3>(p + q);
160 
161  if (C3_test <= 0.0)
162  mooseException(
163  "Cannot take square root of a number less than or equal to zero in the calculation of "
164  "C3_test for the Rashid approximation for the rotation tensor.");
165 
166  const Real C3 = 0.5 * std::sqrt(C3_test); // sin theta_a/(2 sqrt(q))
167 
168  // Calculate incremental rotation. Note that this value is the transpose of that from Rashid,
169  // 93, so we transpose it before storing
170  RankTwoTensor R_incr;
171  R_incr.addIa(C1);
172  for (unsigned int i = 0; i < 3; ++i)
173  for (unsigned int j = 0; j < 3; ++j)
174  R_incr(i, j) += C2 * a[i] * a[j];
175 
176  R_incr(0, 1) += C3 * a[2];
177  R_incr(0, 2) -= C3 * a[1];
178  R_incr(1, 0) -= C3 * a[2];
179  R_incr(1, 2) += C3 * a[0];
180  R_incr(2, 0) += C3 * a[1];
181  R_incr(2, 1) -= C3 * a[0];
182 
183  rotation_increment = R_incr.transpose();
184  break;
185  }
187  {
189  FactorizedRankTwoTensor Uhat = MathUtils::sqrt(Chat);
190  rotation_increment = _Fhat[_qp] * Uhat.inverse().get();
191  total_strain_increment = MathUtils::log(Uhat).get();
192  break;
193  }
194 
195  default:
196  mooseError("ComputeFiniteStrainNOSPD Error: Invalid decomposition type! Please specify : "
197  "TaylorExpansion or "
198  "EigenSolution.");
199  }
200 }
void mooseError(Args &&... args)
const DecompMethod _decomposition_method
void inverse(const std::vector< std::vector< Real >> &m, std::vector< std::vector< Real >> &m_inv)
std::vector< RankTwoTensor > _Fhat
&#39;Incremental&#39; deformation gradient
void addIa(const Real &a)
static RankTwoTensorTempl< Real > transposeTimes(const RankTwoTensorTempl< Real > &)
RankTwoTensorTempl< Real > transpose() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
FactorizedRankTwoTensorTempl< T > inverse() const

◆ decompositionType()

MooseEnum ComputeFiniteStrainNOSPD::decompositionType ( )
static

Definition at line 17 of file ComputeFiniteStrainNOSPD.C.

Referenced by validParams().

18 {
19  return MooseEnum("TaylorExpansion EigenSolution", "TaylorExpansion");
20 }

◆ 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

◆ subtractEigenstrainIncrementFromStrain()

void ComputeFiniteStrainNOSPD::subtractEigenstrainIncrementFromStrain ( RankTwoTensor strain)
protected

Function to compute the mechanical strain tensor by subtracting thermal strain from the total strain.

Definition at line 203 of file ComputeFiniteStrainNOSPD.C.

Referenced by computeQpStrain().

204 {
205  for (unsigned int i = 0; i < _eigenstrains.size(); ++i)
206  {
207  strain -= (*_eigenstrains[i])[_qp];
208  strain += (*_eigenstrains_old[i])[_qp];
209  }
210 }
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains_old

◆ validParams()

InputParameters ComputeFiniteStrainNOSPD::validParams ( )
static

Definition at line 23 of file ComputeFiniteStrainNOSPD.C.

Referenced by ComputePlaneFiniteStrainNOSPD::validParams().

24 {
26  params.addClassDescription(
27  "Class for computing nodal quantities for residual and jacobian calculation "
28  "for peridynamic correspondence models under finite strain assumptions");
29 
30  params.addParam<MooseEnum>("decomposition_method",
32  "Methods to calculate the strain and rotation increments");
33 
34  return params;
35 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
static InputParameters validParams()
void addClassDescription(const std::string &doc_string)
static MooseEnum decompositionType()

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

◆ _decomposition_method

const DecompMethod ComputeFiniteStrainNOSPD::_decomposition_method
private

Definition at line 65 of file ComputeFiniteStrainNOSPD.h.

Referenced by computeQpStrainRotationIncrements().

◆ _deformation_gradient

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_deformation_gradient
protectedinherited

◆ _deformation_gradient_old

const MaterialProperty<RankTwoTensor>& ComputeFiniteStrainNOSPD::_deformation_gradient_old
protected

Material properties to fetch.

Definition at line 49 of file ComputeFiniteStrainNOSPD.h.

Referenced by ComputePlaneFiniteStrainNOSPD::computeQpFhat(), and computeQpFhat().

◆ _eigenstrain_names

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

◆ _eigenstrains

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

◆ _eigenstrains_old

std::vector<const MaterialProperty<RankTwoTensor> *> ComputeFiniteStrainNOSPD::_eigenstrains_old
protected

◆ _Fhat

std::vector<RankTwoTensor> ComputeFiniteStrainNOSPD::_Fhat
protected

'Incremental' deformation gradient

Definition at line 56 of file ComputeFiniteStrainNOSPD.h.

Referenced by ComputePlaneFiniteStrainNOSPD::computeQpFhat(), computeQpFhat(), and computeQpStrainRotationIncrements().

◆ _mechanical_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_mechanical_strain
protectedinherited

◆ _mechanical_strain_old

const MaterialProperty<RankTwoTensor>& ComputeFiniteStrainNOSPD::_mechanical_strain_old
protected

Definition at line 50 of file ComputeFiniteStrainNOSPD.h.

Referenced by computeQpStrain().

◆ _multi

MaterialProperty<Real>& ComputeStrainBaseNOSPD::_multi
protectedinherited

◆ _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().

◆ _rotation_increment

MaterialProperty<RankTwoTensor>& ComputeFiniteStrainNOSPD::_rotation_increment
protected

Definition at line 45 of file ComputeFiniteStrainNOSPD.h.

Referenced by computeQpStrain().

◆ _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().

◆ _strain_increment

MaterialProperty<RankTwoTensor>& ComputeFiniteStrainNOSPD::_strain_increment
protected

Definition at line 44 of file ComputeFiniteStrainNOSPD.h.

Referenced by computeQpStrain().

◆ _strain_rate

MaterialProperty<RankTwoTensor>& ComputeFiniteStrainNOSPD::_strain_rate
protected

Material properties to store.

Definition at line 43 of file ComputeFiniteStrainNOSPD.h.

Referenced by computeQpStrain().

◆ _total_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_total_strain
protectedinherited

◆ _total_strain_old

const MaterialProperty<RankTwoTensor>& ComputeFiniteStrainNOSPD::_total_strain_old
protected

Definition at line 51 of file ComputeFiniteStrainNOSPD.h.

Referenced by computeQpStrain().


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