www.mooseframework.org
Public Member Functions | Protected Attributes | List of all members
GeneralizedPlaneStrainUserObjectNOSPD Class Reference

UserObject class to compute the residual and diagonal Jacobian components for scalar out-of-plane strain variable of generalized plane strain model based on self-stablized non-ordinary state-based peridynamic model. More...

#include <GeneralizedPlaneStrainUserObjectNOSPD.h>

Inheritance diagram for GeneralizedPlaneStrainUserObjectNOSPD:
[legend]

Public Member Functions

 GeneralizedPlaneStrainUserObjectNOSPD (const InputParameters &parameters)
 
virtual void execute () override
 
virtual void initialize () override
 
virtual void threadJoin (const UserObject &uo) override
 
virtual void finalize () override
 
Real returnResidual () const
 Function to return the computed residual. More...
 
Real returnJacobian () const
 Function to return the computed diagonal Jacobian. More...
 

Protected Attributes

const MaterialProperty< RankTwoTensor > & _stress
 Materials property stress. More...
 
const MooseEnum _strain
 Option of strain formulation: SMALL or FINITE. More...
 
const MaterialProperty< RankFourTensor > & _Cijkl
 Elasticity tensor. More...
 
Real _residual
 Residual parameter. More...
 
Real _jacobian
 Jacobian parameter. More...
 
MooseVariable * _bond_status_var
 Bond status aux variable. More...
 
AuxiliarySystem & _aux
 Reference to auxiliary system. More...
 
NumericVector< Number > & _aux_sln
 Solution vector for aux variables. More...
 
PeridynamicsMesh_pdmesh
 Reference to Peridynamic mesh. More...
 
const Function & _pressure
 Applied out-of-plane force parameters. More...
 
const Real _factor
 

Detailed Description

UserObject class to compute the residual and diagonal Jacobian components for scalar out-of-plane strain variable of generalized plane strain model based on self-stablized non-ordinary state-based peridynamic model.

Definition at line 24 of file GeneralizedPlaneStrainUserObjectNOSPD.h.

Constructor & Destructor Documentation

◆ GeneralizedPlaneStrainUserObjectNOSPD()

GeneralizedPlaneStrainUserObjectNOSPD::GeneralizedPlaneStrainUserObjectNOSPD ( const InputParameters &  parameters)

Definition at line 28 of file GeneralizedPlaneStrainUserObjectNOSPD.C.

31  _stress(getMaterialProperty<RankTwoTensor>("stress"))
32 {
33 }

Member Function Documentation

◆ execute()

void GeneralizedPlaneStrainUserObjectNOSPD::execute ( )
overridevirtual

Definition at line 36 of file GeneralizedPlaneStrainUserObjectNOSPD.C.

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 }

◆ finalize()

void GeneralizedPlaneStrainUserObjectBasePD::finalize ( )
overridevirtualinherited

Definition at line 60 of file GeneralizedPlaneStrainUserObjectBasePD.C.

61 {
62  gatherSum(_residual);
63  gatherSum(_jacobian);
64 }

◆ initialize()

void GeneralizedPlaneStrainUserObjectBasePD::initialize ( )
overridevirtualinherited

Definition at line 44 of file GeneralizedPlaneStrainUserObjectBasePD.C.

45 {
46  _residual = 0;
47  _jacobian = 0;
48 }

◆ returnJacobian()

Real GeneralizedPlaneStrainUserObjectBasePD::returnJacobian ( ) const
inherited

Function to return the computed diagonal Jacobian.

Returns
The computed Jacobian

Definition at line 73 of file GeneralizedPlaneStrainUserObjectBasePD.C.

74 {
75  return _jacobian;
76 }

Referenced by GeneralizedPlaneStrainPD::computeJacobian().

◆ returnResidual()

Real GeneralizedPlaneStrainUserObjectBasePD::returnResidual ( ) const
inherited

Function to return the computed residual.

Returns
The computed residual

Definition at line 67 of file GeneralizedPlaneStrainUserObjectBasePD.C.

68 {
69  return _residual;
70 }

Referenced by GeneralizedPlaneStrainPD::computeResidual().

◆ threadJoin()

void GeneralizedPlaneStrainUserObjectBasePD::threadJoin ( const UserObject &  uo)
overridevirtualinherited

Definition at line 51 of file GeneralizedPlaneStrainUserObjectBasePD.C.

52 {
54  static_cast<const GeneralizedPlaneStrainUserObjectBasePD &>(uo);
55  _residual += gpsuo._residual;
56  _jacobian += gpsuo._jacobian;
57 }

Member Data Documentation

◆ _aux

AuxiliarySystem& ElementUserObjectBasePD::_aux
protectedinherited

Reference to auxiliary system.

Definition at line 33 of file ElementUserObjectBasePD.h.

Referenced by NodalAuxVariableUserObjectBasePD::execute(), and NodalAuxVariableUserObjectBasePD::initialize().

◆ _aux_sln

NumericVector<Number>& ElementUserObjectBasePD::_aux_sln
protectedinherited

◆ _bond_status_var

MooseVariable* ElementUserObjectBasePD::_bond_status_var
protectedinherited

◆ _Cijkl

const MaterialProperty<RankFourTensor>& GeneralizedPlaneStrainUserObjectBasePD::_Cijkl
protectedinherited

Elasticity tensor.

Definition at line 49 of file GeneralizedPlaneStrainUserObjectBasePD.h.

Referenced by execute(), and GeneralizedPlaneStrainUserObjectOSPD::execute().

◆ _factor

const Real GeneralizedPlaneStrainUserObjectBasePD::_factor
protectedinherited

◆ _jacobian

Real GeneralizedPlaneStrainUserObjectBasePD::_jacobian
protectedinherited

◆ _pdmesh

PeridynamicsMesh& ElementUserObjectBasePD::_pdmesh
protectedinherited

◆ _pressure

const Function& GeneralizedPlaneStrainUserObjectBasePD::_pressure
protectedinherited

Applied out-of-plane force parameters.

Definition at line 52 of file GeneralizedPlaneStrainUserObjectBasePD.h.

Referenced by execute(), and GeneralizedPlaneStrainUserObjectOSPD::execute().

◆ _residual

Real GeneralizedPlaneStrainUserObjectBasePD::_residual
protectedinherited

◆ _strain

const MooseEnum GeneralizedPlaneStrainUserObjectBasePD::_strain
protectedinherited

Option of strain formulation: SMALL or FINITE.

Definition at line 46 of file GeneralizedPlaneStrainUserObjectBasePD.h.

◆ _stress

const MaterialProperty<RankTwoTensor>& GeneralizedPlaneStrainUserObjectNOSPD::_stress
protected

Materials property stress.

Definition at line 33 of file GeneralizedPlaneStrainUserObjectNOSPD.h.

Referenced by execute().


The documentation for this class was generated from the following files:
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
GeneralizedPlaneStrainUserObjectBasePD::_jacobian
Real _jacobian
Jacobian parameter.
Definition: GeneralizedPlaneStrainUserObjectBasePD.h:60
GeneralizedPlaneStrainUserObjectBasePD::GeneralizedPlaneStrainUserObjectBasePD
GeneralizedPlaneStrainUserObjectBasePD(const InputParameters &parameters)
Definition: GeneralizedPlaneStrainUserObjectBasePD.C:33
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
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
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