www.mooseframework.org
ForceStabilizedSmallStrainMechanicsNOSPD.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 #include "PeridynamicsMesh.h"
12 #include "RankTwoTensor.h"
13 #include "RankFourTensor.h"
14 
16 
17 template <>
18 InputParameters
20 {
21  InputParameters params = validParams<MechanicsBaseNOSPD>();
22  params.addClassDescription("Class for calculating residual and Jacobian for Non-Ordinary "
23  "State-based PeriDynamic solid mechanics formulation using a "
24  "fictitious force method for stabilization.");
25 
26  params.addRequiredParam<unsigned int>(
27  "component",
28  "An integer corresponding to the variable this kernel acts on (0 for x, 1 for y, 2 for z)");
29 
30  return params;
31 }
32 
34  const InputParameters & parameters)
35  : MechanicsBaseNOSPD(parameters),
36  _sf_coeff(getMaterialProperty<Real>("stabilization_force_coeff")),
37  _component(getParam<unsigned int>("component"))
38 {
39 }
40 
41 void
43 {
44  Real sforce = 0.5 * (_sf_coeff[0] + _sf_coeff[1]) *
45  (_disp_var[_component]->getNodalValue(*_current_elem->node_ptr(1)) -
46  _disp_var[_component]->getNodalValue(*_current_elem->node_ptr(0)));
47 
48  // For small strain assumptions, stress measures, i.e., Cauchy stress (Sigma), the first
49  // Piola-Kirchhoff stress (P), and the second Piola-Kirchhoff stress (S) are approximately the
50  // same. Thus, the nodal force state tensors are calculated using the Cauchy stresses,
51  // i.e., T = Sigma * inv(Shape) * xi * multi.
52  // Cauchy stress is calculated as Sigma = C * E.
53 
54  std::vector<Real> nodal_force_comp(_nnodes);
55  for (unsigned int nd = 0; nd < _nnodes; ++nd)
56  nodal_force_comp[nd] = _multi[nd] * (_stress[nd] * _shape2[nd].inverse()).row(_component) *
57  (nd == 0 ? 1 : -1) * _origin_vec_ij;
58 
59  _local_re(0) = -(nodal_force_comp[0] - nodal_force_comp[1] + sforce) * _bond_status_ij;
60  _local_re(1) = -_local_re(0);
61 }
62 
63 void
65 {
66  // excludes dTi/dUj and dTj/dUi contribution which was considered as nonlocal contribution
67  for (_i = 0; _i < _test.size(); ++_i)
68  for (_j = 0; _j < _phi.size(); ++_j)
69  _local_ke(_i, _j) += (_i == 0 ? -1 : 1) * _multi[_j] *
70  (computeDSDU(_component, _j) * _shape2[_j].inverse()).row(_component) *
71  _origin_vec_ij * _bond_status_ij;
72 
73  // compute local stabilization force jacobian
74  Real diag = 0.5 * (_sf_coeff[0] + _sf_coeff[1]);
75  for (_i = 0; _i < _test.size(); ++_i)
76  for (_j = 0; _j < _phi.size(); ++_j)
77  _local_ke(_i, _j) += (_i == _j ? 1 : -1) * diag * _bond_status_ij;
78 }
79 
80 void
82 {
83  // includes dTi/dUj and dTj/dUi contributions
84  for (unsigned int cur_nd = 0; cur_nd < _nnodes; ++cur_nd)
85  {
86  // calculation of jacobian contribution to current_node's neighbors
87  std::vector<dof_id_type> dof(_nnodes);
88  dof[0] = _current_elem->node_ptr(cur_nd)->dof_number(_sys.number(), _var.number(), 0);
89  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(cur_nd));
90  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(cur_nd));
91  for (unsigned int k = 0; k < neighbors.size(); ++k)
92  {
93  Node * node_k = _pdmesh.nodePtr(neighbors[k]);
94  dof[1] = node_k->dof_number(_sys.number(), _var.number(), 0);
95  Real vol_k = _pdmesh.getPDNodeVolume(neighbors[k]);
96 
97  // obtain bond ik's origin vector
98  RealGradient origin_vec_ijk = *node_k - *_pdmesh.nodePtr(_current_elem->node_id(cur_nd));
99 
100  RankTwoTensor dFdUk;
101  dFdUk.zero();
102  for (unsigned int j = 0; j < _dim; ++j)
103  dFdUk(_component, j) =
104  vol_k * _horiz_rad[cur_nd] / origin_vec_ijk.norm() * origin_vec_ijk(j);
105 
106  dFdUk *= _shape2[cur_nd].inverse();
107  RankTwoTensor dSxdUkx = _Jacobian_mult[cur_nd] * 0.5 * (dFdUk.transpose() + dFdUk);
108 
109  // bond status for bond k
110  Real bond_status_ijk = _bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[k]));
111 
112  _local_ke.resize(_test.size(), _phi.size());
113  _local_ke.zero();
114  for (_i = 0; _i < _test.size(); ++_i)
115  for (_j = 0; _j < _phi.size(); ++_j)
116  _local_ke(_i, _j) = (_i == 0 ? -1 : 1) * (_j == 0 ? 0 : 1) * _multi[cur_nd] *
117  (dSxdUkx * _shape2[cur_nd].inverse()).row(_component) *
118  _origin_vec_ij * _bond_status_ij * bond_status_ijk;
119 
120  _assembly.cacheJacobianBlock(_local_ke, _ivardofs_ij, dof, _var.scalingFactor());
121 
122  if (_has_diag_save_in)
123  {
124  unsigned int rows = _test.size();
125  DenseVector<Real> diag(rows);
126  for (unsigned int i = 0; i < rows; ++i)
127  diag(i) = _local_ke(i, i);
128 
129  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
130  for (unsigned int i = 0; i < _diag_save_in.size(); ++i)
131  _diag_save_in[i]->sys().solution().add_vector(diag, _diag_save_in[i]->dofIndices());
132  }
133  }
134  }
135 }
136 
137 void
139  unsigned int coupled_component)
140 {
141  _local_ke.zero();
142  if (coupled_component == 3)
143  {
144  std::vector<RankTwoTensor> dSdT(_nnodes);
145  for (unsigned int nd = 0; nd < _nnodes; ++nd)
146  for (unsigned int es = 0; es < _deigenstrain_dT.size(); ++es)
147  dSdT[nd] = -(_Jacobian_mult[nd] * (*_deigenstrain_dT[es])[nd]);
148 
149  for (_i = 0; _i < _test.size(); ++_i)
150  for (_j = 0; _j < _phi.size(); ++_j)
151  _local_ke(_i, _j) += (_i == 0 ? -1 : 1) * _multi[_j] *
152  (dSdT[_j] * _shape2[_j].inverse()).row(_component) * _origin_vec_ij *
153  _bond_status_ij;
154  }
155  else
156  {
157  for (_i = 0; _i < _test.size(); ++_i)
158  for (_j = 0; _j < _phi.size(); ++_j)
159  _local_ke(_i, _j) +=
160  (_i == 0 ? -1 : 1) * _multi[_j] *
161  (computeDSDU(coupled_component, _j) * _shape2[_j].inverse()).row(_component) *
162  _origin_vec_ij * _bond_status_ij;
163  }
164 }
165 
166 void
168  unsigned int jvar_num, unsigned int coupled_component)
169 {
170  if (coupled_component != 3)
171  {
172  for (unsigned int cur_nd = 0; cur_nd < _nnodes; ++cur_nd)
173  {
174  // calculation of jacobian contribution to current_node's neighbors
175  std::vector<dof_id_type> jvardofs_ijk(_nnodes);
176  jvardofs_ijk[0] = _current_elem->node_ptr(cur_nd)->dof_number(_sys.number(), jvar_num, 0);
177  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(cur_nd));
178  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(cur_nd));
179  for (unsigned int k = 0; k < neighbors.size(); ++k)
180  {
181  Node * node_k = _pdmesh.nodePtr(neighbors[k]);
182  jvardofs_ijk[1] = node_k->dof_number(_sys.number(), jvar_num, 0);
183  Real vol_k = _pdmesh.getPDNodeVolume(neighbors[k]);
184 
185  // obtain bond k's origin vector
186  RealGradient origin_vec_ijk = *node_k - *_pdmesh.nodePtr(_current_elem->node_id(cur_nd));
187 
188  RankTwoTensor dFdUk;
189  dFdUk.zero();
190  for (unsigned int j = 0; j < _dim; ++j)
191  dFdUk(coupled_component, j) =
192  _horiz_rad[cur_nd] / origin_vec_ijk.norm() * origin_vec_ijk(j) * vol_k;
193 
194  dFdUk *= _shape2[cur_nd].inverse();
195  RankTwoTensor dSxdUky = _Jacobian_mult[cur_nd] * 0.5 * (dFdUk.transpose() + dFdUk);
196 
197  // bond status for bond k
198  Real bond_status_ijk = _bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[k]));
199 
200  _local_ke.zero();
201  for (_i = 0; _i < _test.size(); ++_i)
202  for (_j = 0; _j < _phi.size(); ++_j)
203  _local_ke(_i, _j) = (_i == 0 ? -1 : 1) * (_j == 0 ? 0 : 1) * _multi[cur_nd] *
204  (dSxdUky * _shape2[cur_nd].inverse()).row(_component) *
205  _origin_vec_ij * _bond_status_ij * bond_status_ijk;
206 
207  _assembly.cacheJacobianBlock(_local_ke, _ivardofs_ij, jvardofs_ijk, _var.scalingFactor());
208  }
209  }
210  }
211 }
MechanicsBaseNOSPD
Base kernel class for bond-associated correspondence material models.
Definition: MechanicsBaseNOSPD.h:22
ForceStabilizedSmallStrainMechanicsNOSPD::_sf_coeff
const MaterialProperty< Real > & _sf_coeff
Bond based material property for fictitious stabilization force.
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.h:38
MechanicsBaseNOSPD::_shape2
const MaterialProperty< RankTwoTensor > & _shape2
Definition: MechanicsBaseNOSPD.h:39
libMesh::RealGradient
VectorValue< Real > RealGradient
Definition: GrainForceAndTorqueInterface.h:17
MechanicsBaseNOSPD::_multi
const MaterialProperty< Real > & _multi
Material point based material properties.
Definition: MechanicsBaseNOSPD.h:37
MechanicsBaseNOSPD::_Jacobian_mult
const MaterialProperty< RankFourTensor > & _Jacobian_mult
Definition: MechanicsBaseNOSPD.h:44
ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian
void computeLocalOffDiagJacobian(unsigned int coupled_component) override
Function to compute local contribution to the off-diagonal Jacobian at the current nodes.
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.C:138
ForceStabilizedSmallStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian
void computePDNonlocalOffDiagJacobian(unsigned int jvar_num, unsigned int coupled_component) override
Function to compute nonlocal contribution to the off-diagonal Jacobian at the current nodes.
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.C:167
MechanicsBaseNOSPD::_stress
const MaterialProperty< RankTwoTensor > & _stress
Definition: MechanicsBaseNOSPD.h:38
ForceStabilizedSmallStrainMechanicsNOSPD::_component
const unsigned int _component
The index of displacement component.
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.h:41
ForceStabilizedSmallStrainMechanicsNOSPD.h
ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalResidual
virtual void computeLocalResidual() override
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.C:42
PeridynamicsMesh.h
MechanicsBaseNOSPD::computeDSDU
virtual RankTwoTensor computeDSDU(unsigned int component, unsigned int nd)
Function to compute derivative of stress with respect to displacements.
Definition: MechanicsBaseNOSPD.C:49
MechanicsBasePD::_ivardofs_ij
std::vector< dof_id_type > _ivardofs_ij
Current variable dof numbers for nodes i and j.
Definition: MechanicsBasePD.h:69
registerMooseObject
registerMooseObject("PeridynamicsApp", ForceStabilizedSmallStrainMechanicsNOSPD)
validParams< MechanicsBaseNOSPD >
InputParameters validParams< MechanicsBaseNOSPD >()
Definition: MechanicsBaseNOSPD.C:16
ForceStabilizedSmallStrainMechanicsNOSPD::ForceStabilizedSmallStrainMechanicsNOSPD
ForceStabilizedSmallStrainMechanicsNOSPD(const InputParameters &parameters)
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.C:33
RankTwoTensorTempl< Real >
MechanicsBaseNOSPD::_deigenstrain_dT
std::vector< const MaterialProperty< RankTwoTensor > * > _deigenstrain_dT
Definition: MechanicsBaseNOSPD.h:46
ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalJacobian
virtual void computeLocalJacobian() override
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.C:64
ForceStabilizedSmallStrainMechanicsNOSPD
Kernel class for fictitious force stabilized conventional correspondence material model for small str...
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.h:23
MechanicsBasePD::_disp_var
std::vector< MooseVariable * > _disp_var
displacement variables
Definition: MechanicsBasePD.h:50
ForceStabilizedSmallStrainMechanicsNOSPD::computeNonlocalJacobian
virtual void computeNonlocalJacobian() override
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.C:81
validParams< ForceStabilizedSmallStrainMechanicsNOSPD >
InputParameters validParams< ForceStabilizedSmallStrainMechanicsNOSPD >()
Definition: ForceStabilizedSmallStrainMechanicsNOSPD.C:19