https://mooseframework.inl.gov
HorizonStabilizedFormISmallStrainMechanicsNOSPD.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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 
14 
17 {
19  params.addClassDescription(
20  "Class for calculating the residual and the Jacobian for Form I of the horizon-stabilized"
21  "peridynamic correspondence model under small strain assumptions");
22 
23  params.addRequiredParam<unsigned int>(
24  "component",
25  "An integer corresponding to the variable this kernel acts on (0 for x, 1 for y, 2 for z)");
26 
27  return params;
28 }
29 
31  const InputParameters & parameters)
32  : MechanicsBaseNOSPD(parameters), _component(getParam<unsigned int>("component"))
33 {
34 }
35 
36 void
38 {
39  // For small strain assumptions, stress measures, i.e., Cauchy stress (Sigma), the first
40  // Piola-Kirchhoff stress (P), and the second Piola-Kirchhoff stress (S) are approximately the
41  // same. Thus, the nodal force state tensors are calculated using the Cauchy stresses,
42  // i.e., T = Sigma * inv(Shape) * xi * multi.
43  // Cauchy stress is calculated as Sigma = C * E in the ComputeSmallStrainNOSPD material class.
44 
45  for (unsigned int nd = 0; nd < _nnodes; ++nd)
46  for (unsigned int i = 0; i < _nnodes; ++i)
47  _local_re(i) += (i == 0 ? -1 : 1) * _multi[nd] *
48  (_stress[nd] * _shape2[nd].inverse()).row(_component) * _origin_vec *
49  _bond_status;
50 }
51 
52 void
54 {
55  // excludes dTi/dUj and dTj/dUi contributions, which were considered as nonlocal contribution
56  for (unsigned int i = 0; i < _nnodes; ++i)
57  for (unsigned int j = 0; j < _nnodes; ++j)
58  _local_ke(i, j) += (i == 0 ? -1 : 1) * _multi[j] *
59  (computeDSDU(_component, j) * _shape2[j].inverse()).row(_component) *
60  _origin_vec * _bond_status;
61 }
62 
63 void
65 {
66  // includes dTi/dUj and dTj/dUi contributions
67  // excludes contributions to node i and j, which were considered as local contributions
68  for (unsigned int nd = 0; nd < _nnodes; ++nd)
69  {
70  // calculation of jacobian contribution to current_node's neighbors
71  std::vector<dof_id_type> ivardofs(_nnodes);
72  ivardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
73  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
74  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
75 
76  dof_id_type nb_index =
77  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
78  neighbors.begin();
79  std::vector<dof_id_type> dg_neighbors =
80  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
81 
82  Real vol_nb;
83  RealGradient origin_vec_nb;
84  RankTwoTensor dFdUk, dPxdUkx;
85 
86  for (unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
87  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
88  {
89  ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
90  ->dof_number(_sys.number(), _var.number(), 0);
91  vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
92 
93  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
94  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
95 
96  dFdUk.zero();
97  for (unsigned int i = 0; i < _dim; ++i)
98  dFdUk(_component, i) =
99  _horizon_radius[nd] / origin_vec_nb.norm() * origin_vec_nb(i) * vol_nb;
100 
101  dFdUk *= _shape2[nd].inverse();
102 
103  dPxdUkx = _Jacobian_mult[nd] * 0.5 * (dFdUk.transpose() + dFdUk);
104 
105  for (unsigned int i = 0; i < _nnodes; ++i)
106  for (unsigned int j = 0; j < _nnodes; ++j)
107  _local_ke(i, j) = (i == 0 ? -1 : 1) * (j == 0 ? 0 : 1) * _multi[nd] *
108  (dPxdUkx * _shape2[nd].inverse()).row(_component) * _origin_vec *
109  _bond_status;
110 
111  addJacobian(_assembly, _local_ke, _ivardofs, ivardofs, _var.scalingFactor());
112 
113  if (_has_diag_save_in)
114  {
115  unsigned int rows = _nnodes;
116  DenseVector<Real> diag(rows);
117  for (unsigned int i = 0; i < rows; ++i)
118  diag(i) = _local_ke(i, i);
119 
120  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
121  for (unsigned int i = 0; i < _diag_save_in.size(); ++i)
122  {
123  std::vector<dof_id_type> diag_save_in_dofs(2);
124  diag_save_in_dofs[0] = _current_elem->node_ptr(nd)->dof_number(
125  _diag_save_in[i]->sys().number(), _diag_save_in[i]->number(), 0);
126  diag_save_in_dofs[1] =
127  _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
128  ->dof_number(_diag_save_in[i]->sys().number(), _diag_save_in[i]->number(), 0);
129 
130  _diag_save_in[i]->sys().solution().add_vector(diag, diag_save_in_dofs);
131  }
132  }
133  }
134  }
135 }
136 
137 void
139  unsigned int jvar_num, unsigned int coupled_component)
140 {
141  _local_ke.zero();
142  if (_temp_coupled && jvar_num == _temp_var->number()) // temperature is coupled
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 (unsigned int i = 0; i < _nnodes; ++i)
150  for (unsigned int j = 0; j < _nnodes; ++j)
151  _local_ke(i, j) += (i == 0 ? -1 : 1) * _multi[j] *
152  (dSdT[j] * _shape2[j].inverse()).row(_component) * _origin_vec *
153  _bond_status;
154  }
155  else if (_out_of_plane_strain_coupled &&
156  jvar_num == _out_of_plane_strain_var
157  ->number()) // weak plane stress case, out_of_plane_strain is coupled
158  {
159  std::vector<RankTwoTensor> dSdE33(_nnodes);
160  for (unsigned int nd = 0; nd < _nnodes; ++nd)
161  for (unsigned int i = 0; i < 3; ++i)
162  for (unsigned int j = 0; j < 3; ++j)
163  dSdE33[nd](i, j) = _Jacobian_mult[nd](i, j, 2, 2);
164 
165  for (unsigned int i = 0; i < _nnodes; ++i)
166  for (unsigned int j = 0; j < _nnodes; ++j)
167  _local_ke(i, j) += (i == 0 ? -1 : 1) * _multi[j] *
168  (dSdE33[j] * _shape2[j].inverse()).row(_component) * _origin_vec *
169  _bond_status;
170  }
171  else // off-diagonal Jacobian with respect to other displacement variables
172  {
173  // ONLY consider the contributions to node i and j
174  // contributions to their neighbors are considered as nonlocal off-diagonal
175  for (unsigned int i = 0; i < _nnodes; ++i)
176  for (unsigned int j = 0; j < _nnodes; ++j)
177  _local_ke(i, j) +=
178  (i == 0 ? -1 : 1) * _multi[j] *
179  (computeDSDU(coupled_component, j) * _shape2[j].inverse()).row(_component) *
180  _origin_vec * _bond_status;
181  }
182 }
183 
184 void
186  unsigned int jvar_num, unsigned int coupled_component)
187 {
188  if (_temp_coupled && jvar_num == _temp_var->number())
189  {
190  // no nonlocal contribution from temperature
191  }
193  {
194  // no nonlocal contribution from out of plane strain
195  }
196  else
197  {
198  for (unsigned int nd = 0; nd < _nnodes; ++nd)
199  {
200  // calculation of jacobian contribution to current_node's neighbors
201  // NOT including the contribution to nodes i and j, which is considered as local off-diagonal
202  std::vector<dof_id_type> jvardofs(_nnodes);
203  jvardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), jvar_num, 0);
204  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
205  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
206 
207  dof_id_type nb_index =
208  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
209  neighbors.begin();
210  std::vector<dof_id_type> dg_neighbors =
211  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
212 
213  Real vol_nb;
214  RealGradient origin_vec_nb;
215  RankTwoTensor dFdUk, dPxdUky;
216 
217  for (unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
218  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
219  {
220  jvardofs[1] =
221  _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])->dof_number(_sys.number(), jvar_num, 0);
222  vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
223 
224  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
225  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
226 
227  dFdUk.zero();
228  for (unsigned int i = 0; i < _dim; ++i)
229  dFdUk(coupled_component, i) =
230  _horizon_radius[nd] / origin_vec_nb.norm() * origin_vec_nb(i) * vol_nb;
231 
232  dFdUk *= _shape2[nd].inverse();
233 
234  dPxdUky = _Jacobian_mult[nd] * 0.5 * (dFdUk.transpose() + dFdUk);
235 
236  for (unsigned int i = 0; i < _nnodes; ++i)
237  for (unsigned int j = 0; j < _nnodes; ++j)
238  _local_ke(i, j) = (i == 0 ? -1 : 1) * (j == 0 ? 0 : 1) * _multi[nd] *
239  (dPxdUky * _shape2[nd].inverse()).row(_component) * _origin_vec *
240  _bond_status;
241 
242  addJacobian(_assembly, _local_ke, _ivardofs, jvardofs, _var.scalingFactor());
243  }
244  }
245  }
246 }
const bool _temp_coupled
Temperature variable.
virtual RankTwoTensor computeDSDU(unsigned int component, unsigned int nd)
Function to compute derivative of stress with respect to displacements for small strain problems...
auto norm() const -> decltype(std::norm(Real()))
const MaterialProperty< RankFourTensor > & _Jacobian_mult
unsigned int number() const
const MaterialProperty< RankTwoTensor > & _stress
virtual 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...
const unsigned int _component
The index of displacement component.
void inverse(const std::vector< std::vector< Real >> &m, std::vector< std::vector< Real >> &m_inv)
const bool _out_of_plane_strain_coupled
Parameters for out-of-plane strain in weak plane stress formulation.
void addRequiredParam(const std::string &name, const std::string &doc_string)
std::vector< const MaterialProperty< RankTwoTensor > * > _deigenstrain_dT
Kernel class for Form I of the horizon-stabilized peridynamic correspondence model for small strain...
static InputParameters validParams()
virtual void computeLocalOffDiagJacobian(unsigned int jvar_num, unsigned int coupled_component) override
Function to compute local contribution to the off-diagonal Jacobian at the current nodes...
MooseVariable * _out_of_plane_strain_var
std::vector< dof_id_type > _ivardofs
Current variable dof numbers for nodes i and j.
const MaterialProperty< RankTwoTensor > & _shape2
const MaterialProperty< Real > & _multi
Material point based material properties.
RankTwoTensorTempl< Real > transpose() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
registerMooseObject("PeridynamicsApp", HorizonStabilizedFormISmallStrainMechanicsNOSPD)
void addClassDescription(const std::string &doc_string)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
MooseVariable * _temp_var
void ErrorVector unsigned int
Base kernel class for bond-associated correspondence material models.
uint8_t dof_id_type