20 "Class for calculating the residual and the Jacobian for Form I " 21 "of the horizon-stabilized peridynamic correspondence model under finite strain assumptions");
25 "An integer corresponding to the variable this kernel acts on (0 for x, 1 for y, 2 for z)");
46 for (
unsigned int nd = 0; nd < _nnodes; ++nd)
47 for (
unsigned int i = 0; i < _nnodes; ++i)
48 _local_re(i) += (i == 0 ? -1 : 1) *
_multi[nd] *
52 _origin_vec * _bond_status;
59 std::vector<RankTwoTensor> dPxdUx(_nnodes);
60 for (
unsigned int nd = 0; nd < _nnodes; ++nd)
65 for (
unsigned int i = 0; i < _nnodes; ++i)
66 for (
unsigned int j = 0;
j < _nnodes; ++
j)
67 _local_ke(i,
j) += (i == 0 ? -1 : 1) *
_multi[
j] *
76 for (
unsigned int nd = 0; nd < _nnodes; ++nd)
82 std::vector<dof_id_type> ivardofs(_nnodes);
83 ivardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
84 std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
85 std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
88 std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
90 std::vector<dof_id_type> dg_neighbors =
91 _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
97 for (
unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
98 if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
100 ivardofs[1] = _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
101 ->dof_number(_sys.number(), _var.number(), 0);
102 vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
104 origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
105 _pdmesh.getNodeCoord(_current_elem->node_id(nd));
108 for (
unsigned int i = 0; i < _dim; ++i)
110 _horizon_radius[nd] / origin_vec_nb.
norm() * origin_vec_nb(i) * vol_nb;
112 dFdUk *=
_shape2[nd].inverse();
116 for (
unsigned int i = 0; i < 3; ++i)
117 for (
unsigned int j = 0;
j < 3; ++
j)
118 dJdU += detF * invF(
j, i) * dFdUk(i,
j);
121 dSdU = dSdFhat * dFdUk *
_dgrad_old[nd].inverse();
125 for (
unsigned int i = 0; i < 3; ++i)
126 for (
unsigned int J = 0; J < 3; ++J)
127 for (
unsigned int k = 0;
k < 3; ++
k)
128 for (
unsigned int L = 0; L < 3; ++L)
129 dinvFTdU(i, J) += -invF(J,
k) * invF(L, i) * dFdUk(
k, L);
135 for (
unsigned int i = 0; i < _nnodes; ++i)
136 for (
unsigned int j = 0;
j < _nnodes; ++
j)
137 _local_ke(i,
j) = (i == 0 ? -1 : 1) * (
j == 0 ? 0 : 1) *
_multi[nd] *
141 addJacobian(_assembly, _local_ke,
_ivardofs, ivardofs, _var.scalingFactor());
143 if (_has_diag_save_in)
145 unsigned int rows = _nnodes;
146 DenseVector<Real> diag(rows);
147 for (
unsigned int i = 0; i < rows; ++i)
148 diag(i) = _local_ke(i, i);
150 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
151 for (
unsigned int i = 0; i < _diag_save_in.size(); ++i)
153 std::vector<dof_id_type> diag_save_in_dofs(2);
154 diag_save_in_dofs[0] = _current_elem->node_ptr(nd)->dof_number(
155 _diag_save_in[i]->sys().number(), _diag_save_in[i]->number(), 0);
156 diag_save_in_dofs[1] =
157 _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])
158 ->dof_number(_diag_save_in[i]->sys().number(), _diag_save_in[i]->number(), 0);
160 _diag_save_in[i]->sys().solution().add_vector(diag, diag_save_in_dofs);
169 unsigned int jvar_num,
unsigned int coupled_component)
174 std::vector<RankTwoTensor> dSdT(_nnodes);
175 for (
unsigned int nd = 0; nd < _nnodes; ++nd)
180 for (
unsigned int i = 0; i < _nnodes; ++i)
181 for (
unsigned int j = 0;
j < _nnodes; ++
j)
182 _local_ke(i,
j) += (i == 0 ? -1 : 1) *
_multi[
j] *
190 std::vector<RankTwoTensor> dSdE33(_nnodes);
191 for (
unsigned int nd = 0; nd < _nnodes; ++nd)
193 for (
unsigned int i = 0; i < 3; ++i)
194 for (
unsigned int j = 0;
j < 3; ++
j)
197 dSdE33[nd] =
_dgrad[nd].det() * dSdE33[nd] *
_dgrad[nd].inverse().transpose();
200 for (
unsigned int i = 0; i < _nnodes; ++i)
201 for (
unsigned int j = 0;
j < _nnodes; ++
j)
202 _local_ke(i,
j) += (i == 0 ? -1 : 1) *
_multi[
j] *
208 std::vector<RankTwoTensor> dPxdUy(_nnodes);
209 for (
unsigned int nd = 0; nd < _nnodes; ++nd)
215 for (
unsigned int i = 0; i < _nnodes; ++i)
216 for (
unsigned int j = 0;
j < _nnodes; ++
j)
217 _local_ke(i,
j) += (i == 0 ? -1 : 1) *
_multi[
j] *
225 unsigned int jvar_num,
unsigned int coupled_component)
237 for (
unsigned int nd = 0; nd < _nnodes; ++nd)
243 std::vector<dof_id_type> jvardofs(_nnodes);
244 jvardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), jvar_num, 0);
245 std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
246 std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
249 std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
251 std::vector<dof_id_type> dg_neighbors =
252 _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
258 for (
unsigned int nb = 0; nb < dg_neighbors.size(); ++nb)
259 if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
262 _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])->dof_number(_sys.number(), jvar_num, 0);
263 vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
265 origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
266 _pdmesh.getNodeCoord(_current_elem->node_id(nd));
269 for (
unsigned int i = 0; i < _dim; ++i)
270 dFdUk(coupled_component, i) =
271 _horizon_radius[nd] / origin_vec_nb.
norm() * origin_vec_nb(i) * vol_nb;
273 dFdUk *=
_shape2[nd].inverse();
277 for (
unsigned int i = 0; i < 3; ++i)
278 for (
unsigned int j = 0;
j < 3; ++
j)
279 dJdU += detF * invF(
j, i) * dFdUk(i,
j);
282 dSdU = dSdFhat * dFdUk *
_dgrad_old[nd].inverse();
286 for (
unsigned int i = 0; i < 3; ++i)
287 for (
unsigned int J = 0; J < 3; ++J)
288 for (
unsigned int k = 0;
k < 3; ++
k)
289 for (
unsigned int L = 0; L < 3; ++L)
290 dinvFTdU(i, J) += -invF(J,
k) * invF(L, i) * dFdUk(
k, L);
296 for (
unsigned int i = 0; i < _nnodes; ++i)
297 for (
unsigned int j = 0;
j < _nnodes; ++
j)
298 _local_ke(i,
j) = (i == 0 ? -1 : 1) * (
j == 0 ? 0 : 1) *
_multi[nd] *
302 addJacobian(_assembly, _local_ke,
_ivardofs, jvardofs, _var.scalingFactor());
RankFourTensor computeDSDFhat(unsigned int nd)
Function to compute derivative of stress with respect to derived deformation gradient.
Base kernel class for finite strain correspondence models.
const bool _temp_coupled
Temperature variable.
auto norm() const -> decltype(std::norm(Real()))
RankTwoTensor computeDinvFTDU(unsigned int component, unsigned int nd)
Function to compute derivative of deformation gradient inverse with respect to displacements.
const MaterialProperty< RankFourTensor > & _Jacobian_mult
unsigned int number() const
const MaterialProperty< RankTwoTensor > & _stress
const bool _out_of_plane_strain_coupled
Parameters for out-of-plane strain in weak plane stress formulation.
std::vector< const MaterialProperty< RankTwoTensor > * > _deigenstrain_dT
const MaterialProperty< RankTwoTensor > & _dgrad_old
Material point based material property.
virtual RankTwoTensor computeDSDU(unsigned int component, unsigned int nd) override
Function to compute derivative of stress with respect to displacements for small strain problems...
Real computeDJDU(unsigned int component, unsigned int nd)
Function to compute derivative of determinant of deformation gradient with respect to displacements...
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
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
MooseVariable * _temp_var
static const std::string k
void ErrorVector unsigned int
static InputParameters validParams()
const MaterialProperty< RankTwoTensor > & _dgrad