www.mooseframework.org
GeneralizedPlaneStrainUserObjectOSPD.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 "RankFourTensor.h"
12 #include "Function.h"
13 
15 
16 template <>
17 InputParameters
19 {
20  InputParameters params = validParams<GeneralizedPlaneStrainUserObjectBasePD>();
21  params.addClassDescription("Class for calculating the scalar residual and diagonal Jacobian "
22  "entry of generalized plane strain in OSPD formulation");
23 
24  params.addCoupledVar("out_of_plane_stress_variable",
25  "Auxiliary variable name for out-of-plane stress in GPS simulation");
26 
27  return params;
28 }
29 
31  const InputParameters & parameters)
33  _out_of_plane_stress_var(getVar("out_of_plane_stress_variable", 0))
34 {
35 }
36 
37 void
39 {
40  // dof_id_type for node i and j
41  dof_id_type node_i = _current_elem->node_id(0);
42  dof_id_type node_j = _current_elem->node_id(1);
43 
44  // coordinates for node i and j
45  Point coord_i = *_pdmesh.nodePtr(node_i);
46  Point coord_j = *_pdmesh.nodePtr(node_j);
47 
48  // nodal area for node i and j
49  Real nv_i = _pdmesh.getPDNodeVolume(node_i);
50  Real nv_j = _pdmesh.getPDNodeVolume(node_j);
51 
52  // number of neighbors for node i and j, used to avoid repeated accounting nodal stress in
53  // element-wise loop
54  unsigned int nn_i = _pdmesh.getNeighbors(node_i).size();
55  unsigned int nn_j = _pdmesh.getNeighbors(node_j).size();
56 
57  // residual
58  _residual += (_out_of_plane_stress_var->getNodalValue(*_current_elem->node_ptr(0)) -
59  _pressure.value(_t, coord_i) * _factor) *
60  nv_i / nn_i;
61  _residual += (_out_of_plane_stress_var->getNodalValue(*_current_elem->node_ptr(1)) -
62  _pressure.value(_t, coord_j) * _factor) *
63  nv_j / nn_j;
64 
65  // diagonal jacobian
66  _jacobian += _Cijkl[0](2, 2, 2, 2) * nv_i / nn_i + _Cijkl[0](2, 2, 2, 2) * nv_j / nn_j;
67 }
GeneralizedPlaneStrainUserObjectBasePD::_pressure
const Function & _pressure
Applied out-of-plane force parameters.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:52
GeneralizedPlaneStrainUserObjectOSPD::execute
virtual void execute() override
Definition: GeneralizedPlaneStrainUserObjectOSPD.C:38
PeridynamicsMesh::getPDNodeVolume
Real getPDNodeVolume(dof_id_type node_id)
Function to return nodal volume for node node_id.
Definition: PeridynamicsMesh.C:435
GeneralizedPlaneStrainUserObjectOSPD.h
GeneralizedPlaneStrainUserObjectBasePD::_factor
const Real _factor
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:53
GeneralizedPlaneStrainUserObjectOSPD::GeneralizedPlaneStrainUserObjectOSPD
GeneralizedPlaneStrainUserObjectOSPD(const InputParameters &parameters)
Definition: GeneralizedPlaneStrainUserObjectOSPD.C:30
validParams< GeneralizedPlaneStrainUserObjectBasePD >
InputParameters validParams< GeneralizedPlaneStrainUserObjectBasePD >()
Definition: GeneralizedPlaneStrainUserObjectBasePD.C:15
validParams< GeneralizedPlaneStrainUserObjectOSPD >
InputParameters validParams< GeneralizedPlaneStrainUserObjectOSPD >()
Definition: GeneralizedPlaneStrainUserObjectOSPD.C:18
GeneralizedPlaneStrainUserObjectBasePD::_jacobian
Real _jacobian
Jacobian parameter.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:60
GeneralizedPlaneStrainUserObjectBasePD
Base userObject class to compute the residual and diagonal Jacobian components for scalar out-of-plan...
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:23
GeneralizedPlaneStrainUserObjectBasePD::_Cijkl
const MaterialProperty< RankFourTensor > & _Cijkl
Elasticity tensor.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:49
PeridynamicsMesh::getNeighbors
std::vector< dof_id_type > getNeighbors(dof_id_type node_id)
Function to return neighbor nodes indices for node node_id.
Definition: PeridynamicsMesh.C:362
registerMooseObject
registerMooseObject("PeridynamicsApp", GeneralizedPlaneStrainUserObjectOSPD)
GeneralizedPlaneStrainUserObjectBasePD::_residual
Real _residual
Residual parameter.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:57
GeneralizedPlaneStrainUserObjectOSPD::_out_of_plane_stress_var
MooseVariable * _out_of_plane_stress_var
Variable for out-of-plane stress component.
Definition: GeneralizedPlaneStrainUserObjectOSPD.h:33
GeneralizedPlaneStrainUserObjectOSPD
UserObject class to compute the residual and diagonal Jacobian components for scalar out-of-plane str...
Definition: GeneralizedPlaneStrainUserObjectOSPD.h:24
ElementUserObjectBasePD::_pdmesh
PeridynamicsMesh & _pdmesh
Reference to Peridynamic mesh.
Definition: ElementUserObjectBasePD.h:39