www.mooseframework.org
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ComputeStrainBaseNOSPD Class Referenceabstract

Base material class for correspondence material model. More...

#include <ComputeStrainBaseNOSPD.h>

Inheritance diagram for ComputeStrainBaseNOSPD:
[legend]

Public Member Functions

 ComputeStrainBaseNOSPD (const InputParameters &parameters)
 
virtual void initQpStatefulProperties () override
 

Protected Member Functions

virtual void computeProperties () override
 
virtual void computeBondStretch () override
 
virtual void computeQpDeformationGradient ()
 Function to compute bond-associated deformation gradient. More...
 
virtual void computeQpStrain ()=0
 Function to compute strain tensors. More...
 

Protected Attributes

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
 

Detailed Description

Base material class for correspondence material model.

Definition at line 24 of file ComputeStrainBaseNOSPD.h.

Constructor & Destructor Documentation

◆ ComputeStrainBaseNOSPD()

ComputeStrainBaseNOSPD::ComputeStrainBaseNOSPD ( const InputParameters &  parameters)

Definition at line 34 of file ComputeStrainBaseNOSPD.C.

35  : DerivativeMaterialInterface<MechanicsMaterialBasePD>(parameters),
36  _plane_strain(getParam<bool>("plane_strain")),
37  _Cijkl(getMaterialProperty<RankFourTensor>("elasticity_tensor")),
38  _eigenstrain_names(getParam<std::vector<MaterialPropertyName>>("eigenstrain_names")),
40  _shape2(declareProperty<RankTwoTensor>("rank_two_shape_tensor")),
41  _deformation_gradient(declareProperty<RankTwoTensor>("deformation_gradient")),
42  _ddgraddu(declareProperty<RankTwoTensor>("ddeformation_gradient_du")),
43  _ddgraddv(declareProperty<RankTwoTensor>("ddeformation_gradient_dv")),
44  _ddgraddw(declareProperty<RankTwoTensor>("ddeformation_gradient_dw")),
45  _total_strain(declareProperty<RankTwoTensor>("total_strain")),
46  _mechanical_strain(declareProperty<RankTwoTensor>("mechanical_strain")),
47  _multi(declareProperty<Real>("multi"))
48 {
49  for (unsigned int i = 0; i < _eigenstrains.size(); ++i)
50  {
52  _eigenstrains[i] = &getMaterialProperty<RankTwoTensor>(_eigenstrain_names[i]);
53  }
54 }

Member Function Documentation

◆ computeBondStretch()

void ComputeStrainBaseNOSPD::computeBondStretch ( )
overrideprotectedvirtual

Definition at line 156 of file ComputeStrainBaseNOSPD.C.

157 {
158  _total_stretch[0] = _current_length / _origin_length - 1.0;
159  _total_stretch[1] = _total_stretch[0];
160  _mechanical_stretch[0] = _total_stretch[0];
161 
162  Real factor = 1.0;
163  if (_plane_strain)
165 
166  for (auto es : _eigenstrains)
167  _mechanical_stretch[0] -= 0.5 * factor * ((*es)[0](2, 2) + (*es)[1](2, 2));
168 
169  _mechanical_stretch[1] = _mechanical_stretch[0];
170 }

Referenced by computeProperties().

◆ computeProperties()

void ComputeStrainBaseNOSPD::computeProperties ( )
overrideprotectedvirtual

Definition at line 145 of file ComputeStrainBaseNOSPD.C.

146 {
147  setupMeshRelatedData(); // function from base class
148  computeBondCurrentLength(); // current length of a bond from base class
149  computeBondStretch(); // stretch of a bond
150 
151  for (_qp = 0; _qp < _nnodes; ++_qp)
152  computeQpStrain();
153 }

◆ computeQpDeformationGradient()

void ComputeStrainBaseNOSPD::computeQpDeformationGradient ( )
protectedvirtual

Function to compute bond-associated deformation gradient.

Reimplemented in ComputeForceStabilizedSmallStrainNOSPD.

Definition at line 71 of file ComputeStrainBaseNOSPD.C.

72 {
73  _shape2[_qp].zero();
74  _deformation_gradient[_qp].zero();
75  _ddgraddu[_qp].zero();
76  _ddgraddv[_qp].zero();
77  _ddgraddw[_qp].zero();
78  _multi[_qp] = 0.0;
79 
80  // for cases when current bond was broken, assign shape tensor and deformation gradient to unity
81  if (_bond_status_var->getElementalValue(_current_elem) < 0.5)
82  {
83  _shape2[_qp] = RankTwoTensor::initIdentity;
84  _deformation_gradient[_qp] = _shape2[_qp];
85  }
86  else
87  {
88  if (_dim == 2)
89  _shape2[_qp](2, 2) = _deformation_gradient[_qp](2, 2) = 1.0;
90 
91  const Node * cur_nd = _current_elem->node_ptr(_qp);
92  const Node * end_nd = _current_elem->node_ptr(1 - _qp); // two nodes for edge2 element
93  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(cur_nd->id());
94  std::vector<dof_id_type> bonds = _pdmesh.getBonds(cur_nd->id());
95 
96  unsigned int nb =
97  std::find(neighbors.begin(), neighbors.end(), end_nd->id()) - neighbors.begin();
98  std::vector<dof_id_type> dg_neighbors = _pdmesh.getDefGradNeighbors(cur_nd->id(), nb);
99 
100  // calculate the shape tensor and prepare the deformation gradient tensor
101  Real dgnodes_vsum = 0.0;
102  RealGradient ori_vec(_dim), cur_vec(_dim);
103 
104  for (unsigned int j = 0; j < dg_neighbors.size(); ++j)
105  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[j]])) > 0.5)
106  {
107  Node * node_j = _pdmesh.nodePtr(neighbors[dg_neighbors[j]]);
108  Real vol_j = _pdmesh.getPDNodeVolume(neighbors[dg_neighbors[j]]);
109  dgnodes_vsum += vol_j;
110  ori_vec = *node_j - *_pdmesh.nodePtr(cur_nd->id());
111 
112  for (unsigned int k = 0; k < _dim; ++k)
113  cur_vec(k) = ori_vec(k) + _disp_var[k]->getNodalValue(*node_j) -
114  _disp_var[k]->getNodalValue(*cur_nd);
115 
116  Real ori_len = ori_vec.norm();
117  for (unsigned int k = 0; k < _dim; ++k)
118  {
119  for (unsigned int l = 0; l < _dim; ++l)
120  {
121  _shape2[_qp](k, l) += _horiz_rad[_qp] / ori_len * ori_vec(k) * ori_vec(l) * vol_j;
122  _deformation_gradient[_qp](k, l) +=
123  _horiz_rad[_qp] / ori_len * cur_vec(k) * ori_vec(l) * vol_j;
124  }
125  // calculate derivatives of deformation_gradient w.r.t displacements of node i
126  _ddgraddu[_qp](0, k) += -_horiz_rad[_qp] / ori_len * ori_vec(k) * vol_j;
127  _ddgraddv[_qp](1, k) += -_horiz_rad[_qp] / ori_len * ori_vec(k) * vol_j;
128  if (_dim == 3)
129  _ddgraddw[_qp](2, k) += -_horiz_rad[_qp] / ori_len * ori_vec(k) * vol_j;
130  }
131  }
132  // finalize the deformation gradient and its derivatives
133  _deformation_gradient[_qp] *= _shape2[_qp].inverse();
134  _ddgraddu[_qp] *= _shape2[_qp].inverse();
135  _ddgraddv[_qp] *= _shape2[_qp].inverse();
136  _ddgraddw[_qp] *= _shape2[_qp].inverse();
137 
138  // force state multiplier
139  _multi[_qp] = _horiz_rad[_qp] / _origin_length * _node_vol[0] * _node_vol[1] * dgnodes_vsum /
140  _horiz_vol[_qp];
141  }
142 }

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

◆ computeQpStrain()

virtual void ComputeStrainBaseNOSPD::computeQpStrain ( )
protectedpure virtual

Function to compute strain tensors.

Implemented in ComputeFiniteStrainNOSPD, and ComputeSmallStrainNOSPD.

Referenced by computeProperties().

◆ initQpStatefulProperties()

void ComputeStrainBaseNOSPD::initQpStatefulProperties ( )
overridevirtual

Definition at line 57 of file ComputeStrainBaseNOSPD.C.

58 {
59  _mechanical_strain[_qp].zero();
60  _total_strain[_qp].zero();
61  _deformation_gradient[_qp].zero();
62  _deformation_gradient[_qp].addIa(1.0);
63 
64  if (_qrule->n_points() < 2) // Gauss_Lobatto: order = 2n-3 (n is number of qp)
65  mooseError(
66  "NOSPD models require Gauss_Lobatto rule and a minimum of 2 quadrature points, i.e., "
67  "order of FIRST");
68 }

Member Data Documentation

◆ _Cijkl

const MaterialProperty<RankFourTensor>& ComputeStrainBaseNOSPD::_Cijkl
protected

Material properties to fetch.

Definition at line 48 of file ComputeStrainBaseNOSPD.h.

Referenced by computeBondStretch().

◆ _ddgraddu

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_ddgraddu
protected

◆ _ddgraddv

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_ddgraddv
protected

◆ _ddgraddw

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_ddgraddw
protected

◆ _deformation_gradient

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_deformation_gradient
protected

◆ _eigenstrain_names

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

◆ _eigenstrains

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

◆ _mechanical_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_mechanical_strain
protected

◆ _multi

MaterialProperty<Real>& ComputeStrainBaseNOSPD::_multi
protected

◆ _plane_strain

const bool ComputeStrainBaseNOSPD::_plane_strain
protected

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

Definition at line 45 of file ComputeStrainBaseNOSPD.h.

Referenced by computeBondStretch().

◆ _shape2

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_shape2
protected

◆ _total_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBaseNOSPD::_total_strain
protected

The documentation for this class was generated from the following files:
ComputeStrainBaseNOSPD::_mechanical_strain
MaterialProperty< RankTwoTensor > & _mechanical_strain
Definition: ComputeStrainBaseNOSPD.h:62
ComputeStrainBaseNOSPD::_ddgraddv
MaterialProperty< RankTwoTensor > & _ddgraddv
Definition: ComputeStrainBaseNOSPD.h:58
ComputeStrainBaseNOSPD::_shape2
MaterialProperty< RankTwoTensor > & _shape2
Material properties to store.
Definition: ComputeStrainBaseNOSPD.h:54
libMesh::RealGradient
VectorValue< Real > RealGradient
Definition: GrainForceAndTorqueInterface.h:17
ComputeStrainBaseNOSPD::_plane_strain
const bool _plane_strain
Plane strain problem or not, this is only used for mechanical stretch calculation.
Definition: ComputeStrainBaseNOSPD.h:45
ComputeStrainBaseNOSPD::_deformation_gradient
MaterialProperty< RankTwoTensor > & _deformation_gradient
Definition: ComputeStrainBaseNOSPD.h:55
ElasticityTensorTools::getIsotropicPoissonsRatio
T getIsotropicPoissonsRatio(const RankFourTensorTempl< T > &elasticity_tensor)
Get the Poisson's modulus for an isotropic elasticity tensor param elasticity_tensor the tensor (must...
Definition: ElasticityTensorTools.h:115
ComputeStrainBaseNOSPD::_Cijkl
const MaterialProperty< RankFourTensor > & _Cijkl
Material properties to fetch.
Definition: ComputeStrainBaseNOSPD.h:48
ComputeStrainBaseNOSPD::_eigenstrain_names
std::vector< MaterialPropertyName > _eigenstrain_names
Definition: ComputeStrainBaseNOSPD.h:49
ComputeStrainBaseNOSPD::_multi
MaterialProperty< Real > & _multi
Definition: ComputeStrainBaseNOSPD.h:64
ComputeStrainBaseNOSPD::_ddgraddw
MaterialProperty< RankTwoTensor > & _ddgraddw
Definition: ComputeStrainBaseNOSPD.h:59
ComputeStrainBaseNOSPD::_eigenstrains
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
Definition: ComputeStrainBaseNOSPD.h:50
ComputeStrainBaseNOSPD::computeQpStrain
virtual void computeQpStrain()=0
Function to compute strain tensors.
ComputeStrainBaseNOSPD::computeBondStretch
virtual void computeBondStretch() override
Definition: ComputeStrainBaseNOSPD.C:156
ComputeStrainBaseNOSPD::_ddgraddu
MaterialProperty< RankTwoTensor > & _ddgraddu
Definition: ComputeStrainBaseNOSPD.h:57
ComputeStrainBaseNOSPD::_total_strain
MaterialProperty< RankTwoTensor > & _total_strain
Definition: ComputeStrainBaseNOSPD.h:61