https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
WeakPlaneStressNOSPD Class Reference

Kernel class for weak plane stress formulation based on Form I of the horizon-stabilized peridynamic correspondence model. More...

#include <WeakPlaneStressNOSPD.h>

Inheritance diagram for WeakPlaneStressNOSPD:
[legend]

Public Types

typedef DerivativeMaterialPropertyNameInterface::SymbolName SymbolName
 

Public Member Functions

 WeakPlaneStressNOSPD (const InputParameters &parameters)
 
virtual void computeLocalResidual () override
 
virtual void computeLocalJacobian () override
 
virtual void computeLocalOffDiagJacobian (unsigned int jvar_num, unsigned int coupled_component) override
 Function to compute local contribution to the off-diagonal Jacobian at the current nodes. More...
 
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. More...
 
virtual void computeOffDiagJacobian (unsigned int jvar) override
 
virtual void initialSetup () override
 
virtual void prepare () override
 
const GenericMaterialProperty< U, is_ad > & getDefaultMaterialProperty (const std::string &name)
 
const GenericMaterialProperty< U, is_ad > & getDefaultMaterialPropertyByName (const std::string &name)
 
void validateDerivativeMaterialPropertyBase (const std::string &base)
 
const MaterialPropertyName derivativePropertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName derivativePropertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName derivativePropertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName derivativePropertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
GenericMaterialProperty< U, is_ad > & declarePropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, unsigned int v2, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, unsigned int v1, unsigned int v2=libMesh::invalid_uint, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, const SymbolName &c1, unsigned int v2, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivative (const std::string &base, unsigned int v1, unsigned int v2=libMesh::invalid_uint, unsigned int v3=libMesh::invalid_uint)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< VariableName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const std::vector< SymbolName > &c)
 
const GenericMaterialProperty< U, is_ad > & getMaterialPropertyDerivativeByName (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2="", const SymbolName &c3="")
 
void validateCoupling (const MaterialPropertyName &base, const std::vector< VariableName > &c, bool validate_aux=true)
 
void validateCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateCoupling (const MaterialPropertyName &base, const std::vector< VariableName > &c, bool validate_aux=true)
 
void validateCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateNonlinearCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
void validateNonlinearCoupling (const MaterialPropertyName &base, const VariableName &c1="", const VariableName &c2="", const VariableName &c3="")
 
const MaterialPropertyName propertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName propertyName (const MaterialPropertyName &base, const std::vector< SymbolName > &c) const
 
const MaterialPropertyName propertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName propertyNameFirst (const MaterialPropertyName &base, const SymbolName &c1) const
 
const MaterialPropertyName propertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName propertyNameSecond (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2) const
 
const MaterialPropertyName propertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 
const MaterialPropertyName propertyNameThird (const MaterialPropertyName &base, const SymbolName &c1, const SymbolName &c2, const SymbolName &c3) const
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

virtual RankTwoTensor computeDSDU (unsigned int component, unsigned int nd)
 Function to compute derivative of stress with respect to displacements for small strain problems. More...
 

Protected Attributes

std::vector< MooseVariable * > _disp_var
 displacement variables More...
 
unsigned int _ndisp
 number of displacement components More...
 
const std::vector< RealGradient > * _orientation
 Vector of bond in current configuration. More...
 
std::vector< dof_id_type_ivardofs
 Current variable dof numbers for nodes i and j. More...
 
std::vector< Real_weights
 weights used for the current element to obtain the nodal stress More...
 
RealGradient _current_vec
 Vector of bond in current configuration. More...
 
RealGradient _current_unit_vec
 Unit vector of bond in current configuration. More...
 
const MaterialProperty< Real > & _multi
 Material point based material properties. More...
 
const MaterialProperty< RankTwoTensor > & _stress
 
const MaterialProperty< RankTwoTensor > & _shape2
 
const MaterialProperty< RankTwoTensor > & _dgrad
 
const MaterialProperty< RankTwoTensor > & _ddgraddu
 
const MaterialProperty< RankTwoTensor > & _ddgraddv
 
const MaterialProperty< RankTwoTensor > & _ddgraddw
 
const MaterialProperty< RankFourTensor > & _Jacobian_mult
 
const std::vector< MaterialPropertyName > _eigenstrain_names
 
std::vector< const MaterialProperty< RankTwoTensor > * > _deigenstrain_dT
 
const bool _temp_coupled
 Temperature variable. More...
 
MooseVariable_temp_var
 
const bool _out_of_plane_strain_coupled
 Parameters for out-of-plane strain in weak plane stress formulation. More...
 
MooseVariable_out_of_plane_strain_var
 

Detailed Description

Kernel class for weak plane stress formulation based on Form I of the horizon-stabilized peridynamic correspondence model.

Definition at line 18 of file WeakPlaneStressNOSPD.h.

Constructor & Destructor Documentation

◆ WeakPlaneStressNOSPD()

WeakPlaneStressNOSPD::WeakPlaneStressNOSPD ( const InputParameters parameters)

Definition at line 27 of file WeakPlaneStressNOSPD.C.

28  : MechanicsBaseNOSPD(parameters)
29 {
30 }
MechanicsBaseNOSPD(const InputParameters &parameters)

Member Function Documentation

◆ computeDSDU()

RankTwoTensor MechanicsBaseNOSPD::computeDSDU ( unsigned int  component,
unsigned int  nd 
)
protectedvirtualinherited

Function to compute derivative of stress with respect to displacements for small strain problems.

Parameters
componentThe index of displacement component
ndThe local index of element node (either 1 or 2 for Edge2 element)
Returns
The calculated derivative

Reimplemented in MechanicsFiniteStrainBaseNOSPD.

Definition at line 47 of file MechanicsBaseNOSPD.C.

Referenced by GeneralizedPlaneStrainOffDiagNOSPD::computeDispFullOffDiagJacobianScalar(), GeneralizedPlaneStrainOffDiagNOSPD::computeDispPartialOffDiagJacobianScalar(), ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeLocalJacobian(), computeLocalOffDiagJacobian(), ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), and HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian().

48 {
49  // compute the derivative of stress w.r.t the solution components for small strain
50  RankTwoTensor dSdU;
51  if (component == 0)
52  dSdU = _Jacobian_mult[nd] * 0.5 * (_ddgraddu[nd].transpose() + _ddgraddu[nd]);
53  else if (component == 1)
54  dSdU = _Jacobian_mult[nd] * 0.5 * (_ddgraddv[nd].transpose() + _ddgraddv[nd]);
55  else if (component == 2)
56  dSdU = _Jacobian_mult[nd] * 0.5 * (_ddgraddw[nd].transpose() + _ddgraddw[nd]);
57 
58  return dSdU;
59 }
const MaterialProperty< RankTwoTensor > & _ddgraddu
const MaterialProperty< RankFourTensor > & _Jacobian_mult
static const std::string component
Definition: NS.h:153
const MaterialProperty< RankTwoTensor > & _ddgraddw
const MaterialProperty< RankTwoTensor > & _ddgraddv

◆ computeLocalJacobian()

void WeakPlaneStressNOSPD::computeLocalJacobian ( )
overridevirtual

Definition at line 40 of file WeakPlaneStressNOSPD.C.

41 {
42  for (unsigned int i = 0; i < _nnodes; ++i)
43  for (unsigned int j = 0; j < _nnodes; ++j)
44  _local_ke(i, j) += (i == j ? 1 : 0) * _Jacobian_mult[i](2, 2, 2, 2) * _weights[i] *
45  _node_vol[i] * _bond_status;
46 }
const MaterialProperty< RankFourTensor > & _Jacobian_mult
std::vector< Real > _weights
weights used for the current element to obtain the nodal stress
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")

◆ computeLocalOffDiagJacobian()

void WeakPlaneStressNOSPD::computeLocalOffDiagJacobian ( unsigned int  ,
unsigned int   
)
overridevirtual

Function to compute local contribution to the off-diagonal Jacobian at the current nodes.

Parameters
coupled_componentThe coupled variable number

Reimplemented from MechanicsBasePD.

Definition at line 49 of file WeakPlaneStressNOSPD.C.

51 {
52  _local_ke.zero();
53  if (_temp_coupled && jvar_num == _temp_var->number()) // for coupled temperature variable
54  {
55  std::vector<RankTwoTensor> dSdT(_nnodes);
56  for (unsigned int nd = 0; nd < _nnodes; nd++)
57  for (unsigned int es = 0; es < _deigenstrain_dT.size(); es++)
58  dSdT[nd] = -_Jacobian_mult[nd] * (*_deigenstrain_dT[es])[nd];
59 
60  for (unsigned int i = 0; i < _nnodes; ++i)
61  for (unsigned int j = 0; j < _nnodes; ++j)
62  _local_ke(i, j) +=
63  (i == j ? 1 : 0) * dSdT[i](2, 2) * _weights[i] * _node_vol[i] * _bond_status;
64  }
65  else // for coupled displacement variables
66  {
67  // dSidUi and dSjdUj are considered as local off-diagonal Jacobian
68  // contributions to their neighbors are considered as nonlocal off-diagonal
69  for (unsigned int i = 0; i < _nnodes; ++i)
70  for (unsigned int j = 0; j < _nnodes; ++j)
71  _local_ke(i, j) += (i == j ? 1 : 0) * computeDSDU(coupled_component, j)(2, 2) *
72  _weights[j] * _node_vol[j] * _bond_status;
73  }
74 }
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...
const MaterialProperty< RankFourTensor > & _Jacobian_mult
unsigned int number() const
std::vector< const MaterialProperty< RankTwoTensor > * > _deigenstrain_dT
std::vector< Real > _weights
weights used for the current element to obtain the nodal stress
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
MooseVariable * _temp_var

◆ computeLocalResidual()

void WeakPlaneStressNOSPD::computeLocalResidual ( )
overridevirtual

Definition at line 33 of file WeakPlaneStressNOSPD.C.

34 {
35  _local_re(0) = _stress[0](2, 2) * _weights[0] * _node_vol[0] * _bond_status;
36  _local_re(1) = _stress[1](2, 2) * _weights[1] * _node_vol[1] * _bond_status;
37 }
const MaterialProperty< RankTwoTensor > & _stress
std::vector< Real > _weights
weights used for the current element to obtain the nodal stress

◆ computeOffDiagJacobian()

void MechanicsBasePD::computeOffDiagJacobian ( unsigned int  jvar)
overridevirtualinherited

Definition at line 73 of file MechanicsBasePD.C.

74 {
75  prepare();
76 
77  if (jvar_num == _var.number())
78  computeJacobian();
79  else
80  {
81  unsigned int coupled_component = 0;
82  bool active = false;
83 
84  for (unsigned int i = 0; i < _dim; ++i)
85  if (jvar_num == _disp_var[i]->number())
86  {
87  coupled_component = i;
88  active = true;
89  }
90 
91  if (_temp_coupled && jvar_num == _temp_var->number())
92  active = true;
93 
95  active = true;
96 
97  if (active)
98  {
99  prepareMatrixTag(_assembly, _var.number(), jvar_num);
100  computeLocalOffDiagJacobian(jvar_num, coupled_component);
101  accumulateTaggedLocalMatrix();
102 
103  if (_use_full_jacobian)
104  computePDNonlocalOffDiagJacobian(jvar_num, coupled_component);
105  }
106  }
107 }
const bool _temp_coupled
Temperature variable.
unsigned int number() const
virtual void computeLocalOffDiagJacobian(unsigned int, unsigned int)
Function to compute local contribution to the off-diagonal Jacobian at the current nodes...
const bool _out_of_plane_strain_coupled
Parameters for out-of-plane strain in weak plane stress formulation.
std::vector< MooseVariable * > _disp_var
displacement variables
virtual void prepare() override
MooseVariable * _out_of_plane_strain_var
virtual void computePDNonlocalOffDiagJacobian(unsigned int, unsigned int)
Function to compute nonlocal contribution to the off-diagonal Jacobian at the current nodes...
MooseVariable * _temp_var

◆ computePDNonlocalOffDiagJacobian()

void WeakPlaneStressNOSPD::computePDNonlocalOffDiagJacobian ( unsigned int  ,
unsigned int   
)
overridevirtual

Function to compute nonlocal contribution to the off-diagonal Jacobian at the current nodes.

Parameters
jvar_numThe number of the first coupled variable
coupled_componentThe component number of the second coupled variable

Reimplemented from MechanicsBasePD.

Definition at line 77 of file WeakPlaneStressNOSPD.C.

79 {
80  if (_temp_coupled && jvar_num == _temp_var->number())
81  {
82  // no nonlocal contribution from temperature
83  }
84  else
85  {
86  for (unsigned int nd = 0; nd < _nnodes; nd++)
87  {
88  // calculation of jacobian contribution to current_node's neighbors
89  // NOT including the contribution to nodes i and j, which is considered as local
90  // off-diagonal
91  std::vector<dof_id_type> jvardofs(_nnodes);
92  jvardofs[0] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), jvar_num, 0);
93  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_elem->node_id(nd));
94  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_elem->node_id(nd));
95 
96  dof_id_type nb_index =
97  std::find(neighbors.begin(), neighbors.end(), _current_elem->node_id(1 - nd)) -
98  neighbors.begin();
99  std::vector<dof_id_type> dg_neighbors =
100  _pdmesh.getBondDeformationGradientNeighbors(_current_elem->node_id(nd), nb_index);
101 
102  Real vol_nb, weight_nb;
103  RealGradient origin_vec_nb;
104  RankTwoTensor dFdUk, dPdUk;
105 
106  for (unsigned int nb = 0; nb < dg_neighbors.size(); nb++)
107  if (_bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[dg_neighbors[nb]])) > 0.5)
108  {
109  jvardofs[1] =
110  _pdmesh.nodePtr(neighbors[dg_neighbors[nb]])->dof_number(_sys.number(), jvar_num, 0);
111  vol_nb = _pdmesh.getNodeVolume(neighbors[dg_neighbors[nb]]);
112 
113  origin_vec_nb = _pdmesh.getNodeCoord(neighbors[dg_neighbors[nb]]) -
114  _pdmesh.getNodeCoord(_current_elem->node_id(nd));
115 
116  weight_nb = _horizon_radius[nd] / origin_vec_nb.norm();
117 
118  dFdUk.zero();
119  for (unsigned int i = 0; i < _dim; ++i)
120  dFdUk(coupled_component, i) = weight_nb * origin_vec_nb(i) * vol_nb;
121 
122  dFdUk *= _shape2[nd].inverse();
123 
124  dPdUk = _Jacobian_mult[nd] * 0.5 *
125  (dFdUk.transpose() * _dgrad[nd] + _dgrad[nd].transpose() * dFdUk);
126 
127  _local_ke.zero();
128  _local_ke(nd, 1) = dPdUk(2, 2) * _weights[nd] * _node_vol[nd] * _bond_status;
129 
130  addJacobian(_assembly, _local_ke, _ivardofs, jvardofs, _var.scalingFactor());
131  }
132  }
133  }
134 }
const bool _temp_coupled
Temperature variable.
auto norm() const -> decltype(std::norm(Real()))
const MaterialProperty< RankFourTensor > & _Jacobian_mult
unsigned int number() const
std::vector< Real > _weights
weights used for the current element to obtain the nodal stress
std::vector< dof_id_type > _ivardofs
Current variable dof numbers for nodes i and j.
const MaterialProperty< RankTwoTensor > & _shape2
RankTwoTensorTempl< Real > transpose() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
MooseVariable * _temp_var
uint8_t dof_id_type
const MaterialProperty< RankTwoTensor > & _dgrad

◆ initialSetup()

void MechanicsBasePD::initialSetup ( )
overridevirtualinherited

Definition at line 45 of file MechanicsBasePD.C.

46 {
47  _orientation = &_assembly.getFE(FEType(), 1)->get_dxyzdxi();
48 }
const std::vector< RealGradient > * _orientation
Vector of bond in current configuration.

◆ prepare()

void MechanicsBasePD::prepare ( )
overridevirtualinherited

Definition at line 51 of file MechanicsBasePD.C.

Referenced by MechanicsBasePD::computeOffDiagJacobian(), GeneralizedPlaneStrainOffDiagNOSPD::computeOffDiagJacobianScalar(), and GeneralizedPlaneStrainOffDiagOSPD::computeOffDiagJacobianScalar().

52 {
54 
55  _ivardofs.resize(_nnodes);
56  _weights.resize(_nnodes);
57  for (unsigned int nd = 0; nd < _nnodes; ++nd)
58  {
59  _ivardofs[nd] = _current_elem->node_ptr(nd)->dof_number(_sys.number(), _var.number(), 0);
60  _weights[nd] = _pdmesh.getNeighborWeight(
61  _current_elem->node_id(nd),
62  _pdmesh.getNeighborIndex(_current_elem->node_id(nd), _current_elem->node_id(1 - nd)));
63  }
64 
65  for (unsigned int i = 0; i < _dim; ++i)
66  _current_vec(i) = _origin_vec(i) + _disp_var[i]->getNodalValue(*_current_elem->node_ptr(1)) -
67  _disp_var[i]->getNodalValue(*_current_elem->node_ptr(0));
68 
70 }
auto norm() const -> decltype(std::norm(Real()))
RealGradient _current_unit_vec
Unit vector of bond in current configuration.
std::vector< Real > _weights
weights used for the current element to obtain the nodal stress
std::vector< MooseVariable * > _disp_var
displacement variables
RealGradient _current_vec
Vector of bond in current configuration.
std::vector< dof_id_type > _ivardofs
Current variable dof numbers for nodes i and j.
virtual void prepare()
Function to precalculate data which will be used in the derived classes.

◆ validParams()

InputParameters WeakPlaneStressNOSPD::validParams ( )
static

Definition at line 17 of file WeakPlaneStressNOSPD.C.

18 {
20  params.addClassDescription(
21  "Class for calculating the residual and the Jacobian for the peridynamic plane "
22  "stress model using weak formulation based on peridynamic correspondence models");
23 
24  return params;
25 }
static InputParameters validParams()
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _current_unit_vec

RealGradient MechanicsBasePD::_current_unit_vec
protectedinherited

◆ _current_vec

RealGradient MechanicsBasePD::_current_vec
protectedinherited

◆ _ddgraddu

const MaterialProperty<RankTwoTensor>& MechanicsBaseNOSPD::_ddgraddu
protectedinherited

◆ _ddgraddv

const MaterialProperty<RankTwoTensor>& MechanicsBaseNOSPD::_ddgraddv
protectedinherited

◆ _ddgraddw

const MaterialProperty<RankTwoTensor>& MechanicsBaseNOSPD::_ddgraddw
protectedinherited

◆ _deigenstrain_dT

std::vector<const MaterialProperty<RankTwoTensor> *> MechanicsBaseNOSPD::_deigenstrain_dT
protectedinherited

◆ _dgrad

const MaterialProperty<RankTwoTensor>& MechanicsBaseNOSPD::_dgrad
protectedinherited

◆ _disp_var

std::vector<MooseVariable *> MechanicsBasePD::_disp_var
protectedinherited

◆ _eigenstrain_names

const std::vector<MaterialPropertyName> MechanicsBaseNOSPD::_eigenstrain_names
protectedinherited

Definition at line 43 of file MechanicsBaseNOSPD.h.

Referenced by MechanicsBaseNOSPD::MechanicsBaseNOSPD().

◆ _ivardofs

std::vector<dof_id_type> MechanicsBasePD::_ivardofs
protectedinherited

◆ _Jacobian_mult

const MaterialProperty<RankFourTensor>& MechanicsBaseNOSPD::_Jacobian_mult
protectedinherited

◆ _multi

const MaterialProperty<Real>& MechanicsBaseNOSPD::_multi
protectedinherited

Material point based material properties.

Definition at line 35 of file MechanicsBaseNOSPD.h.

Referenced by GeneralizedPlaneStrainOffDiagNOSPD::computeDispFullOffDiagJacobianScalar(), GeneralizedPlaneStrainOffDiagNOSPD::computeDispPartialOffDiagJacobianScalar(), ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeLocalJacobian(), ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalResidual(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computeLocalResidual(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeLocalResidual(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeLocalResidual(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeLocalResidual(), ForceStabilizedSmallStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeNonlocalResidual(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeNonlocalResidual(), ForceStabilizedSmallStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian(), and HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian().

◆ _ndisp

unsigned int MechanicsBasePD::_ndisp
protectedinherited

number of displacement components

Definition at line 55 of file MechanicsBasePD.h.

Referenced by MechanicsBasePD::MechanicsBasePD().

◆ _orientation

const std::vector<RealGradient>* MechanicsBasePD::_orientation
protectedinherited

Vector of bond in current configuration.

Definition at line 63 of file MechanicsBasePD.h.

Referenced by MechanicsBasePD::initialSetup().

◆ _out_of_plane_strain_coupled

const bool MechanicsBasePD::_out_of_plane_strain_coupled
protectedinherited

◆ _out_of_plane_strain_var

MooseVariable* MechanicsBasePD::_out_of_plane_strain_var
protectedinherited

◆ _shape2

const MaterialProperty<RankTwoTensor>& MechanicsBaseNOSPD::_shape2
protectedinherited

Definition at line 37 of file MechanicsBaseNOSPD.h.

Referenced by GeneralizedPlaneStrainOffDiagNOSPD::computeDispFullOffDiagJacobianScalar(), GeneralizedPlaneStrainOffDiagNOSPD::computeDispPartialOffDiagJacobianScalar(), ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeLocalJacobian(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeLocalJacobian(), ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeLocalOffDiagJacobian(), ForceStabilizedSmallStrainMechanicsNOSPD::computeLocalResidual(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computeLocalResidual(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeLocalResidual(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeLocalResidual(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeLocalResidual(), ForceStabilizedSmallStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeNonlocalJacobian(), HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computeNonlocalResidual(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computeNonlocalResidual(), computePDNonlocalOffDiagJacobian(), ForceStabilizedSmallStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian(), HorizonStabilizedFormISmallStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian(), HorizonStabilizedFormIFiniteStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian(), HorizonStabilizedFormIIFiniteStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian(), and HorizonStabilizedFormIISmallStrainMechanicsNOSPD::computePDNonlocalOffDiagJacobian().

◆ _stress

const MaterialProperty<RankTwoTensor>& MechanicsBaseNOSPD::_stress
protectedinherited

◆ _temp_coupled

const bool MechanicsBasePD::_temp_coupled
protectedinherited

◆ _temp_var

MooseVariable* MechanicsBasePD::_temp_var
protectedinherited

◆ _weights

std::vector<Real> MechanicsBasePD::_weights
protectedinherited

weights used for the current element to obtain the nodal stress

Definition at line 69 of file MechanicsBasePD.h.

Referenced by computeLocalJacobian(), computeLocalOffDiagJacobian(), computeLocalResidual(), computePDNonlocalOffDiagJacobian(), and MechanicsBasePD::prepare().


The documentation for this class was generated from the following files: