www.mooseframework.org
GeneralizedPlaneStrainUserObjectNOSPD.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 "RankTwoTensor.h"
12 #include "RankFourTensor.h"
13 #include "Function.h"
14 
16 
17 template <>
18 InputParameters
20 {
21  InputParameters params = validParams<GeneralizedPlaneStrainUserObjectBasePD>();
22  params.addClassDescription("Class for calculating the scalar residual and diagonal Jacobian "
23  "entry of generalized plane strain in SNOSPD formulation");
24 
25  return params;
26 }
27 
29  const InputParameters & parameters)
31  _stress(getMaterialProperty<RankTwoTensor>("stress"))
32 {
33 }
34 
35 void
37 {
38  // dof_id_type for node i and j
39  dof_id_type node_i = _current_elem->node_id(0);
40  dof_id_type node_j = _current_elem->node_id(1);
41 
42  // coordinates for node i and j
43  Point coord_i = *_pdmesh.nodePtr(node_i);
44  Point coord_j = *_pdmesh.nodePtr(node_j);
45 
46  // nodal area for node i and j
47  Real nv_i = _pdmesh.getPDNodeVolume(node_i);
48  Real nv_j = _pdmesh.getPDNodeVolume(node_j);
49 
50  // sum of volumes of material points used in bond-associated deformation gradient calculation
51  dof_id_type id_j_in_i = _pdmesh.getNeighborIndex(node_i, node_j);
52  dof_id_type id_i_in_j = _pdmesh.getNeighborIndex(node_j, node_i);
53 
54  Real dg_vol_frac_i = _pdmesh.getDefGradVolFraction(node_i, id_j_in_i);
55  Real dg_vol_frac_j = _pdmesh.getDefGradVolFraction(node_j, id_i_in_j);
56 
57  // residual
58  _residual += (_stress[0](2, 2) - _pressure.value(_t, coord_i) * _factor) * nv_i * dg_vol_frac_i;
59  _residual += (_stress[1](2, 2) - _pressure.value(_t, coord_j) * _factor) * nv_j * dg_vol_frac_j;
60 
61  // diagonal jacobian
62  _jacobian +=
63  _Cijkl[0](2, 2, 2, 2) * nv_i * dg_vol_frac_i + _Cijkl[1](2, 2, 2, 2) * nv_j * dg_vol_frac_j;
64 }
GeneralizedPlaneStrainUserObjectBasePD::_pressure
const Function & _pressure
Applied out-of-plane force parameters.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:52
PeridynamicsMesh::getPDNodeVolume
Real getPDNodeVolume(dof_id_type node_id)
Function to return nodal volume for node node_id.
Definition: PeridynamicsMesh.C:435
GeneralizedPlaneStrainUserObjectBasePD::_factor
const Real _factor
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:53
GeneralizedPlaneStrainUserObjectNOSPD::execute
virtual void execute() override
Definition: GeneralizedPlaneStrainUserObjectNOSPD.C:36
validParams< GeneralizedPlaneStrainUserObjectBasePD >
InputParameters validParams< GeneralizedPlaneStrainUserObjectBasePD >()
Definition: GeneralizedPlaneStrainUserObjectBasePD.C:15
GeneralizedPlaneStrainUserObjectBasePD::_jacobian
Real _jacobian
Jacobian parameter.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:60
GeneralizedPlaneStrainUserObjectNOSPD.h
GeneralizedPlaneStrainUserObjectBasePD
Base userObject class to compute the residual and diagonal Jacobian components for scalar out-of-plan...
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:23
GeneralizedPlaneStrainUserObjectNOSPD::_stress
const MaterialProperty< RankTwoTensor > & _stress
Materials property stress.
Definition: GeneralizedPlaneStrainUserObjectNOSPD.h:33
GeneralizedPlaneStrainUserObjectNOSPD::GeneralizedPlaneStrainUserObjectNOSPD
GeneralizedPlaneStrainUserObjectNOSPD(const InputParameters &parameters)
Definition: GeneralizedPlaneStrainUserObjectNOSPD.C:28
GeneralizedPlaneStrainUserObjectNOSPD
UserObject class to compute the residual and diagonal Jacobian components for scalar out-of-plane str...
Definition: GeneralizedPlaneStrainUserObjectNOSPD.h:24
registerMooseObject
registerMooseObject("PeridynamicsApp", GeneralizedPlaneStrainUserObjectNOSPD)
GeneralizedPlaneStrainUserObjectBasePD::_Cijkl
const MaterialProperty< RankFourTensor > & _Cijkl
Elasticity tensor.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:49
GeneralizedPlaneStrainUserObjectBasePD::_residual
Real _residual
Residual parameter.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:57
PeridynamicsMesh::getNeighborIndex
dof_id_type getNeighborIndex(dof_id_type node_i, dof_id_type node_j)
Function to return the local neighbor index of node_j from node_i's neighbor list.
Definition: PeridynamicsMesh.C:368
RankTwoTensorTempl< Real >
ElementUserObjectBasePD::_pdmesh
PeridynamicsMesh & _pdmesh
Reference to Peridynamic mesh.
Definition: ElementUserObjectBasePD.h:39
PeridynamicsMesh::getDefGradVolFraction
Real getDefGradVolFraction(dof_id_type node_id, dof_id_type neighbor_id)
Function to return summation of volumes of bond-associated neighbors used in the deformation gradient...
Definition: PeridynamicsMesh.C:453
validParams< GeneralizedPlaneStrainUserObjectNOSPD >
InputParameters validParams< GeneralizedPlaneStrainUserObjectNOSPD >()
Definition: GeneralizedPlaneStrainUserObjectNOSPD.C:19