Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
HorizonStabilizedFormIISmallStrainMechanicsNOSPD.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 II "
21  "of the horizon-stabilized 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] * _horizon_radius[nd] / _origin_vec.norm() *
48  (_stress[nd] * _shape2[nd].inverse()).row(_component) * _origin_vec *
49  _node_vol[1 - nd] * _bond_status;
50 }
51 
52 void
54 {
55  for (unsigned int nd = 0; nd < _nnodes; ++nd)
56  {
57  // calculation of residual contribution to current_node's neighbors
58  std::vector<dof_id_type> ivardofs(_nnodes);
59  ivardofs[0] = _ivardofs[nd];
60  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
61  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
62 
63  dof_id_type nb_index =
64  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
65  neighbors.begin();
66  std::vector<dof_id_type> dg_neighbors =
67  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
68 
69  RealGradient origin_vec_nb;
70  Real node_vol_nb;
71 
72  for (unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
73  if (neighbors[dg_neighbors[nb]] != _current_elem->node_id(1 - nd) &&
74  _bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
75  {
76  ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
77  ->dof_number(_sys.number(), _var.number(), 0);
78  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
79  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
80  node_vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
81 
82  for (unsigned int i = 0; i < _nnodes; ++i)
83  _local_re(i) = (i == 0 ? -1 : 1) * _multi[nd] * _horizon_radius[nd] /
84  origin_vec_nb.norm() *
85  (_stress[nd] * _shape2[nd].inverse()).row(_component) * origin_vec_nb *
86  node_vol_nb * _bond_status;
87 
88  // cache the residual contribution
89  addResiduals(_assembly, _local_re, ivardofs, _var.scalingFactor());
90 
91  if (_has_save_in)
92  {
93  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
94  for (unsigned int i = 0; i < _save_in.size(); ++i)
95  {
96  std::vector<dof_id_type> save_in_dofs(_nnodes);
97  save_in_dofs[0] = _current_elem->node_ptr(nd)->dof_number(
98  _save_in[i]->sys().number(), _save_in[i]->number(), 0);
99  save_in_dofs[1] =
100  _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
101  ->dof_number(_save_in[i]->sys().number(), _save_in[i]->number(), 0);
102 
103  _save_in[i]->sys().solution().add_vector(_local_re, save_in_dofs);
104  }
105  }
106  }
107  }
108 }
109 
110 void
112 {
113  for (unsigned int nd = 0; nd < _nnodes; ++nd)
114  {
115  std::vector<dof_id_type> ivardofs(_nnodes);
116  ivardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
117  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
118  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
119 
120  dof_id_type nb_index =
121  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
122  neighbors.begin();
123  std::vector<dof_id_type> dg_neighbors =
124  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
125 
126  RealGradient origin_vec_nb;
127  Real node_vol_nb;
128 
129  for (unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
130  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
131  {
132  ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
133  ->dof_number(_sys.number(), _var.number(), 0);
134  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
135  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
136  node_vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
137 
138  for (unsigned int i = 0; i < _nnodes; ++i)
139  for (unsigned int j = 0; j < _nnodes; ++j)
140  _local_ke(i, j) =
141  (i == 0 ? -1 : 1) * (j == 0 ? 1 : 0) * _multi[nd] * _horizon_radius[nd] /
142  origin_vec_nb.norm() *
143  (computeDSDU(_component, nd) * _shape2[nd].inverse()).row(_component) *
144  origin_vec_nb * node_vol_nb * _bond_status;
145 
146  addJacobian(_assembly, _local_ke, ivardofs, ivardofs, _var.scalingFactor());
147  }
148  }
149  _local_ke.zero();
150 }
151 
152 void
154 {
155  // includes dTi/dUj and dTj/dUi contributions
156  for (unsigned int nd = 0; nd < _nnodes; ++nd)
157  {
158  // calculation of jacobian contribution to current_node's neighbors
159  std::vector<dof_id_type> ivardofs(_nnodes), jvardofs(_nnodes);
160  ivardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
161  jvardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
162  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
163  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
164 
165  dof_id_type nb_index =
166  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
167  neighbors.begin();
168  std::vector<dof_id_type> dg_neighbors =
169  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
170 
171  RealGradient origin_vec_nb1;
172  Real node_vol_nb1;
173 
174  for (unsigned int nb1 = 0; nb1 < dg_neighbors.size(); ++nb1)
175  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb1]])) > 0.5)
176  {
177  ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb1]])
178  ->dof_number(_sys.number(), _var.number(), 0);
179  origin_vec_nb1 = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb1]]) -
180  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
181  node_vol_nb1 = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb1]]);
182 
183  Real vol_nb2;
184  RealGradient origin_vec_nb2;
185  RankTwoTensor dFdUk, dPxdUkx;
186 
187  for (unsigned int nb2 = 0; nb2 < dg_neighbors.size(); ++nb2)
188  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb2]])) > 0.5)
189  {
190  jvardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb2]])
191  ->dof_number(_sys.number(), _var.number(), 0);
192  vol_nb2 = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb2]]);
193 
194  origin_vec_nb2 = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb2]]) -
195  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
196 
197  dFdUk.zero();
198  for (unsigned int i = 0; i < _dim; ++i)
199  dFdUk(_component, i) =
200  _horizon_radius[nd] / origin_vec_nb2.norm() * origin_vec_nb2(i) * vol_nb2;
201 
202  dFdUk *= _shape2[nd].inverse();
203  dPxdUkx = _Jacobian_mult[nd] * 0.5 * (dFdUk.transpose() + dFdUk);
204 
205  for (unsigned int i = 0; i < _nnodes; ++i)
206  for (unsigned int j = 0; j < _nnodes; ++j)
207  _local_ke(i, j) = (i == 0 ? -1 : 1) * (j == 0 ? 0 : 1) * _multi[nd] *
208  _horizon_radius[nd] / origin_vec_nb1.norm() *
209  (dPxdUkx * _shape2[nd].inverse()).row(_component) *
210  origin_vec_nb1 * node_vol_nb1 * _bond_status;
211 
212  addJacobian(_assembly, _local_ke, ivardofs, jvardofs, _var.scalingFactor());
213 
214  if (_has_diag_save_in)
215  {
216  unsigned int rows = _nnodes;
217  DenseVector<Real> diag(rows);
218  for (unsigned int i = 0; i < rows; ++i)
219  diag(i) = _local_ke(i, i);
220 
221  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
222  for (unsigned int i = 0; i < _diag_save_in.size(); ++i)
223  {
224  std::vector<dof_id_type> diag_save_in_dofs(2);
225  diag_save_in_dofs[0] = _current_elem->node_ptr(nd)->dof_number(
226  _diag_save_in[i]->sys().number(), _diag_save_in[i]->number(), 0);
227  diag_save_in_dofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb2]])
228  ->dof_number(_diag_save_in[i]->sys().number(),
229  _diag_save_in[i]->number(),
230  0);
231 
232  _diag_save_in[i]->sys().solution().add_vector(diag, diag_save_in_dofs);
233  }
234  }
235  }
236  }
237  }
238 }
239 
240 void
242  unsigned int jvar_num, unsigned int coupled_component)
243 {
244  if (_temp_coupled && jvar_num == _temp_var->number()) // temperature is coupled
245  {
246  std::vector<RankTwoTensor> dSdT(_nnodes);
247  for (unsigned int nd = 0; nd < _nnodes; ++nd)
248  for (unsigned int es = 0; es < _deigenstrain_dT.size(); ++es)
249  dSdT[nd] = -_Jacobian_mult[nd] * (*_deigenstrain_dT[es])[nd];
250 
251  for (unsigned int nd = 0; nd < _nnodes; ++nd)
252  {
253  std::vector<dof_id_type> ivardofs(_nnodes), jvardofs(_nnodes);
254  ivardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
255  jvardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), jvar_num, 0);
256  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
257  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
258 
259  dof_id_type nb_index =
260  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
261  neighbors.begin();
262  std::vector<dof_id_type> dg_neighbors =
263  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
264 
265  RealGradient origin_vec_nb;
266  Real node_vol_nb;
267 
268  for (unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
269  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
270  {
271  ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
272  ->dof_number(_sys.number(), _var.number(), 0);
273  jvardofs[1] =
274  _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])->dof_number(_sys.number(), jvar_num, 0);
275  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
276  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
277  node_vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
278 
279  for (unsigned int i = 0; i < _nnodes; ++i)
280  for (unsigned int j = 0; j < _nnodes; ++j)
281  _local_ke(i, j) = (i == 0 ? -1 : 1) * (j == 0 ? 1 : 0) * _multi[nd] *
282  _horizon_radius[nd] / origin_vec_nb.norm() *
283  (dSdT[nd] * _shape2[nd].inverse()).row(_component) * origin_vec_nb *
284  node_vol_nb * _bond_status;
285 
286  addJacobian(_assembly, _local_ke, ivardofs, jvardofs, _var.scalingFactor());
287  }
288  }
289  _local_ke.zero();
290  }
291  else if (_out_of_plane_strain_coupled &&
292  jvar_num == _out_of_plane_strain_var
293  ->number()) // weak plane stress case, out_of_plane_strain is coupled
294  {
295  std::vector<RankTwoTensor> dSdE33(_nnodes);
296  for (unsigned int nd = 0; nd < _nnodes; ++nd)
297  for (unsigned int i = 0; i < 3; ++i)
298  for (unsigned int j = 0; j < 3; ++j)
299  dSdE33[nd](i, j) = _Jacobian_mult[nd](i, j, 2, 2);
300 
301  for (unsigned int nd = 0; nd < _nnodes; ++nd)
302  {
303  std::vector<dof_id_type> ivardofs(_nnodes), jvardofs(_nnodes);
304  ivardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
305  jvardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), jvar_num, 0);
306  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
307  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
308 
309  dof_id_type nb_index =
310  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
311  neighbors.begin();
312  std::vector<dof_id_type> dg_neighbors =
313  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
314 
315  RealGradient origin_vec_nb;
316  Real node_vol_nb;
317 
318  for (unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
319  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
320  {
321  ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
322  ->dof_number(_sys.number(), _var.number(), 0);
323  jvardofs[1] =
324  _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])->dof_number(_sys.number(), jvar_num, 0);
325  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
326  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
327  node_vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
328 
329  for (unsigned int i = 0; i < _nnodes; ++i)
330  for (unsigned int j = 0; j < _nnodes; ++j)
331  _local_ke(i, j) = (i == 0 ? -1 : 1) * (j == 0 ? 1 : 0) * _multi[nd] *
332  _horizon_radius[nd] / origin_vec_nb.norm() *
333  (dSdE33[nd] * _shape2[nd].inverse()).row(_component) *
334  origin_vec_nb * node_vol_nb * _bond_status;
335 
336  addJacobian(_assembly, _local_ke, ivardofs, jvardofs, _var.scalingFactor());
337  }
338  }
339  _local_ke.zero();
340  }
341  else // off-diagonal Jacobian with respect to other displacement variables
342  {
343  for (unsigned int nd = 0; nd < _nnodes; ++nd)
344  {
345  std::vector<dof_id_type> ivardofs(_nnodes), jvardofs(_nnodes);
346  ivardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
347  jvardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), jvar_num, 0);
348  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
349  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
350 
351  dof_id_type nb_index =
352  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
353  neighbors.begin();
354  std::vector<dof_id_type> dg_neighbors =
355  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
356 
357  RealGradient origin_vec_nb;
358  Real node_vol_nb;
359 
360  for (unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
361  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
362  {
363  ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
364  ->dof_number(_sys.number(), _var.number(), 0);
365  jvardofs[1] =
366  _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])->dof_number(_sys.number(), jvar_num, 0);
367  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
368  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
369  node_vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
370 
371  for (unsigned int i = 0; i < _nnodes; ++i)
372  for (unsigned int j = 0; j < _nnodes; ++j)
373  _local_ke(i, j) =
374  (i == 0 ? -1 : 1) * (j == 0 ? 1 : 0) * _multi[nd] * _horizon_radius[nd] /
375  origin_vec_nb.norm() *
376  (computeDSDU(coupled_component, nd) * _shape2[nd].inverse()).row(_component) *
377  origin_vec_nb * node_vol_nb * _bond_status;
378 
379  addJacobian(_assembly, _local_ke, ivardofs, jvardofs, _var.scalingFactor());
380  }
381  }
382  _local_ke.zero();
383  }
384 }
385 
386 void
388  unsigned int jvar_num, unsigned int coupled_component)
389 {
390  if (_temp_coupled && jvar_num == _temp_var->number())
391  {
392  // no nonlocal contribution from temperature
393  }
395  {
396  // no nonlocal contribution from out of plane strain
397  }
398  else
399  {
400  for (unsigned int nd = 0; nd < _nnodes; ++nd)
401  {
402  // calculation of jacobian contribution to current_node's neighbors
403  // NOT including the contribution to nodes i and j, which is considered as local off-diagonal
404  std::vector<dof_id_type> ivardofs(_nnodes), jvardofs(_nnodes);
405  ivardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
406  jvardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), jvar_num, 0);
407  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
408  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
409 
410  dof_id_type nb_index =
411  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
412  neighbors.begin();
413  std::vector<dof_id_type> dg_neighbors =
414  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
415 
416  RealGradient origin_vec_nb1;
417  Real node_vol_nb1;
418 
419  for (unsigned int nb1 = 0; nb1 < dg_neighbors.size(); ++nb1)
420  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb1]])) > 0.5)
421  {
422  ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb1]])
423  ->dof_number(_sys.number(), _var.number(), 0);
424  origin_vec_nb1 = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb1]]) -
425  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
426  node_vol_nb1 = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb1]]);
427 
428  Real vol_nb2;
429  RealGradient origin_vec_nb2;
430  RankTwoTensor dFdUk, dPxdUky;
431 
432  for (unsigned int nb2 = 0; nb2 < dg_neighbors.size(); ++nb2)
433  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb2]])) >
434  0.5)
435  {
436  jvardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb2]])
437  ->dof_number(_sys.number(), jvar_num, 0);
438  vol_nb2 = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb2]]);
439 
440  origin_vec_nb2 = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb2]]) -
441  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
442 
443  dFdUk.zero();
444  for (unsigned int i = 0; i < _dim; ++i)
445  dFdUk(coupled_component, i) =
446  _horizon_radius[nd] / origin_vec_nb2.norm() * origin_vec_nb2(i) * vol_nb2;
447 
448  dFdUk *= _shape2[nd].inverse();
449  dPxdUky = _Jacobian_mult[nd] * 0.5 * (dFdUk.transpose() + dFdUk);
450 
451  for (unsigned int i = 0; i < _nnodes; ++i)
452  for (unsigned int j = 0; j < _nnodes; ++j)
453  _local_ke(i, j) = (i == 0 ? -1 : 1) * (j == 0 ? 0 : 1) * _multi[nd] *
454  _horizon_radius[nd] / origin_vec_nb1.norm() *
455  (dPxdUky * _shape2[nd].inverse()).row(_component) *
456  origin_vec_nb1 * node_vol_nb1 * _bond_status;
457 
458  addJacobian(_assembly, _local_ke, ivardofs, jvardofs, _var.scalingFactor());
459  }
460  }
461  }
462  }
463 }
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 unsigned int _component
The index of displacement component.
Kernel class for Form II of the horizon-associated peridynamic correspondence material model for smal...
const MaterialProperty< RankTwoTensor > & _stress
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
static InputParameters validParams()
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
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...
registerMooseObject("PeridynamicsApp", HorizonStabilizedFormIISmallStrainMechanicsNOSPD)
void addClassDescription(const std::string &doc_string)
virtual void computeLocalOffDiagJacobian(unsigned int, unsigned int coupled_component) override
Function to compute local contribution to the off-diagonal Jacobian at the current nodes...
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