www.mooseframework.org
ComputeSmallStrainMaterialBaseOSPD.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 
12 template <>
13 InputParameters
15 {
16  InputParameters params = validParams<ParametricMaterialBasePD>();
17  params.addClassDescription("Base class for ordinary state-based peridynamic mechanics models");
18 
19  return params;
20 }
21 
23  const InputParameters & parameters)
24  : ParametricMaterialBasePD(parameters),
25  _bond_force_i_j(declareProperty<Real>("bond_force_i_j")),
26  _bond_dfdU_i_j(declareProperty<Real>("bond_dfdU_i_j")),
27  _bond_dfdT_i_j(declareProperty<Real>("bond_dfdT_i_j")),
28  _bond_dfdE_i_j(declareProperty<Real>("bond_dfdE_i_j")),
29  _d(2)
30 {
31 }
32 
33 void
35 {
37  {
38  _bond_force_ij[_qp] = 2.0 * _b *
39  (_mechanical_stretch[_qp] +
41  _alpha * (0.5 * (_temp[0] + _temp[1]) - _temp_ref))) /
43  _bond_dfdT_ij[_qp] = -2.0 * _b * (1.0 + _poissons_ratio) * 0.5 * _alpha / _origin_length *
44  _node_vol[0] * _node_vol[1];
45 
46  _bond_force_i_j[_qp] =
47  2.0 * _a * _d[_qp] * _d[_qp] *
48  (_mechanical_stretch[_qp] + _alpha * (0.5 * (_temp[0] + _temp[1]) - _temp_ref) -
49  _alpha * (_temp[_qp] - _temp_ref) +
51  _node_vol[0] * _node_vol[1];
52  _bond_dfdT_i_j[_qp] = -2.0 * _a * _d[_qp] * _d[_qp] * (1.0 + _poissons_ratio) * _alpha *
53  _node_vol[0] * _node_vol[1];
54  }
55  else
56  {
57  _bond_force_ij[_qp] =
58  2.0 * _b * _mechanical_stretch[_qp] / _origin_length * _node_vol[0] * _node_vol[1];
59  _bond_dfdT_ij[_qp] = -2.0 * _b * 0.5 * _alpha / _origin_length * _node_vol[0] * _node_vol[1];
60 
61  _bond_force_i_j[_qp] =
62  2.0 * _a * _d[_qp] * _d[_qp] *
63  (_mechanical_stretch[_qp] + _alpha * (0.5 * (_temp[0] + _temp[1]) - _temp_ref) -
64  _alpha * (_temp[_qp] - _temp_ref)) *
65  _node_vol[0] * _node_vol[1];
66  _bond_dfdT_i_j[_qp] = -2.0 * _a * _d[_qp] * _d[_qp] * _alpha * _node_vol[0] * _node_vol[1];
67  }
68 
70  _bond_dfdU_i_j[_qp] = 2.0 * _a * _d[_qp] * _d[_qp] / _origin_length * _node_vol[0] * _node_vol[1];
71 
73  _bond_dfdE_i_j[_qp] =
74  2.0 * _a * _d[_qp] * _d[_qp] * _poissons_ratio * _node_vol[0] * _node_vol[1];
75 }
PeridynamicsMaterialBase::_origin_length
Real _origin_length
Definition: PeridynamicsMaterialBase.h:43
MechanicsMaterialBasePD::_mechanical_stretch
MaterialProperty< Real > & _mechanical_stretch
Definition: MechanicsMaterialBasePD.h:48
ParametricMaterialBasePD::_bond_dfdE_ij
MaterialProperty< Real > & _bond_dfdE_ij
Definition: ParametricMaterialBasePD.h:70
ComputeSmallStrainMaterialBaseOSPD::_b
Real _b
Definition: ComputeSmallStrainMaterialBaseOSPD.h:39
ComputeSmallStrainMaterialBaseOSPD::_bond_dfdU_i_j
MaterialProperty< Real > & _bond_dfdU_i_j
Definition: ComputeSmallStrainMaterialBaseOSPD.h:32
ComputeSmallStrainMaterialBaseOSPD::_bond_force_i_j
MaterialProperty< Real > & _bond_force_i_j
Material properties to store.
Definition: ComputeSmallStrainMaterialBaseOSPD.h:31
ParametricMaterialBasePD::_scalar_out_of_plane_strain
const VariableValue & _scalar_out_of_plane_strain
Definition: ParametricMaterialBasePD.h:53
ParametricMaterialBasePD::_temp_ref
const Real _temp_ref
Reference temperature.
Definition: ParametricMaterialBasePD.h:60
ParametricMaterialBasePD::_poissons_ratio
Real _poissons_ratio
Definition: ParametricMaterialBasePD.h:79
ComputeSmallStrainMaterialBaseOSPD::_bond_dfdE_i_j
MaterialProperty< Real > & _bond_dfdE_i_j
Definition: ComputeSmallStrainMaterialBaseOSPD.h:34
ParametricMaterialBasePD::_bond_dfdT_ij
MaterialProperty< Real > & _bond_dfdT_ij
Definition: ParametricMaterialBasePD.h:69
ParametricMaterialBasePD::_scalar_out_of_plane_strain_coupled
const bool _scalar_out_of_plane_strain_coupled
Scalar out-of-plane component of strain tensor for generalized plane strain.
Definition: ParametricMaterialBasePD.h:52
ComputeSmallStrainMaterialBaseOSPD::computeBondForce
virtual void computeBondForce() override
Function to compute force of a bond.
Definition: ComputeSmallStrainMaterialBaseOSPD.C:34
validParams< ParametricMaterialBasePD >
InputParameters validParams< ParametricMaterialBasePD >()
Definition: ParametricMaterialBasePD.C:17
ComputeSmallStrainMaterialBaseOSPD::_a
Real _a
Model parameters.
Definition: ComputeSmallStrainMaterialBaseOSPD.h:38
validParams< ComputeSmallStrainMaterialBaseOSPD >
InputParameters validParams< ComputeSmallStrainMaterialBaseOSPD >()
Definition: ComputeSmallStrainMaterialBaseOSPD.C:14
PeridynamicsMaterialBase::_node_vol
std::vector< Real > _node_vol
Definition: PeridynamicsMaterialBase.h:39
ParametricMaterialBasePD::_bond_force_ij
MaterialProperty< Real > & _bond_force_ij
Material properties to store.
Definition: ParametricMaterialBasePD.h:67
ParametricMaterialBasePD
Base material class for bond-based and ordinary state-based peridynamic models, i....
Definition: ParametricMaterialBasePD.h:24
ComputeSmallStrainMaterialBaseOSPD::_bond_dfdT_i_j
MaterialProperty< Real > & _bond_dfdT_i_j
Definition: ComputeSmallStrainMaterialBaseOSPD.h:33
ParametricMaterialBasePD::_bond_dfdU_ij
MaterialProperty< Real > & _bond_dfdU_ij
Definition: ParametricMaterialBasePD.h:68
ComputeSmallStrainMaterialBaseOSPD::_d
std::vector< Real > _d
Definition: ComputeSmallStrainMaterialBaseOSPD.h:40
ComputeSmallStrainMaterialBaseOSPD.h
ParametricMaterialBasePD::_alpha
Real _alpha
Definition: ParametricMaterialBasePD.h:64
ComputeSmallStrainMaterialBaseOSPD::ComputeSmallStrainMaterialBaseOSPD
ComputeSmallStrainMaterialBaseOSPD(const InputParameters &parameters)
Definition: ComputeSmallStrainMaterialBaseOSPD.C:22
ParametricMaterialBasePD::_temp
std::vector< Real > _temp
Temperature variable.
Definition: ParametricMaterialBasePD.h:57