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

Aux Kernel class to postprocess and output the strain and stress components and equivalents for peridynamic models excluding correspondence material models. More...

#include <NodalRankTwoPD.h>

Inheritance diagram for NodalRankTwoPD:
[legend]

Public Member Functions

 NodalRankTwoPD (const InputParameters &parameters)
 

Protected Member Functions

Real computeValue () override
 
virtual RankTwoTensor computeNodalTotalStrain ()
 Function to compute the total strain tensor at each discrete material point. More...
 
virtual RankTwoTensor computeNodalMechanicalStrain ()
 Function to compute the elastic strain tensor at each discrete material point. More...
 
virtual RankTwoTensor computeNodalStress ()
 Function to compute the stress tensor at each discrete materials point. More...
 

Protected Attributes

bool _has_temp
 
MooseVariable * _temp_var
 coupled temperature variable More...
 
MooseVariable * _bond_status_var
 bond_status variable More...
 
const bool _plane_stress
 plane stress problem or not More...
 
const Real _temp_ref
 reference temperature More...
 
std::string _rank_two_tensor
 name of rank two tensor to be processed: total_strain, mechanical_strain or stress More...
 
std::string _output_type
 output type: component or scalar More...
 
MooseEnum _scalar_type
 specific scalar type to be output More...
 
PeridynamicsMesh_pdmesh
 Reference to peridynamic mesh object. More...
 
const unsigned int _dim
 Problem dimension. More...
 
bool _scalar_out_of_plane_strain_coupled
 variable for generalized plane strain cases More...
 
const VariableValue & _scalar_out_of_plane_strain
 
const Real _youngs_modulus
 material constants More...
 
const Real _poissons_ratio
 
const Real _alpha
 
const unsigned int _i
 component index More...
 
const unsigned int _j
 
const Point _point1
 Points for direction dependent scalar output. More...
 
const Point _point2
 
std::vector< MooseVariable * > _disp_var
 @ More...
 
RankFourTensor _Cijkl
 elasticity tensor More...
 

Detailed Description

Aux Kernel class to postprocess and output the strain and stress components and equivalents for peridynamic models excluding correspondence material models.

Definition at line 26 of file NodalRankTwoPD.h.

Constructor & Destructor Documentation

◆ NodalRankTwoPD()

NodalRankTwoPD::NodalRankTwoPD ( const InputParameters &  parameters)

Definition at line 58 of file NodalRankTwoPD.C.

59  : AuxKernelBasePD(parameters),
60  _has_temp(isCoupled("temperature")),
61  _temp_var(_has_temp ? getVar("temperature", 0) : nullptr),
62  _bond_status_var(&_subproblem.getStandardVariable(_tid, "bond_status")),
63  _scalar_out_of_plane_strain_coupled(isCoupledScalar("scalar_out_of_plane_strain")),
65  ? coupledScalarValue("scalar_out_of_plane_strain")
66  : _zero),
67  _plane_stress(getParam<bool>("plane_stress")),
68  _youngs_modulus(getParam<Real>("youngs_modulus")),
69  _poissons_ratio(getParam<Real>("poissons_ratio")),
70  _alpha(getParam<Real>("thermal_expansion_coeff")),
71  _temp_ref(getParam<Real>("stress_free_temperature")),
72  _rank_two_tensor(getParam<std::string>("rank_two_tensor")),
73  _output_type(getParam<std::string>("output_type")),
74  _scalar_type(getParam<MooseEnum>("scalar_type")),
75  _i(getParam<unsigned int>("index_i")),
76  _j(getParam<unsigned int>("index_j")),
77  _point1(parameters.get<Point>("point1")),
78  _point2(parameters.get<Point>("point2"))
79 {
80  if (!_var.isNodal())
81  mooseError("NodalRankTwoPD operates on nodal variable!");
82 
83  if (coupledComponents("displacements") != _dim)
84  mooseError("Number of displacements should be the same as problem dimension!");
85  for (unsigned int i = 0; i < coupledComponents("displacements"); ++i)
86  _disp_var.push_back(getVar("displacements", i));
87 
88  if (_rank_two_tensor == "stress" && !isParamValid("youngs_modulus") &&
89  !isParamValid("poissons_ratio"))
90  mooseError("Both Young's modulus and Poisson's ratio must be provided for stress calculation");
91 
92  std::vector<Real> iso_const(2);
93  iso_const[0] =
94  _youngs_modulus * _poissons_ratio / ((1.0 + _poissons_ratio) * (1.0 - 2.0 * _poissons_ratio));
95  iso_const[1] = _youngs_modulus / (2.0 * (1.0 + _poissons_ratio));
96 
97  // fill elasticity tensor
98  _Cijkl.fillFromInputVector(iso_const, RankFourTensor::symmetric_isotropic);
99 
100  if (_has_temp && !isParamValid("thermal_expansion_coeff"))
101  mooseError("thermal_expansion_coeff is required when temperature is coupled");
102 
103  if (_output_type == "component" && !(isParamValid("index_i") || isParamValid("index_j")))
104  mooseError("The output_type is 'component', but no 'index_i' and 'index_j' are provided!");
105 
106  if (_output_type == "scalar" && !isParamValid("scalar_type"))
107  mooseError("The output_type is 'scalar', but no 'scalar_type' is provided!");
108 }

Member Function Documentation

◆ computeNodalMechanicalStrain()

RankTwoTensor NodalRankTwoPD::computeNodalMechanicalStrain ( )
protectedvirtual

Function to compute the elastic strain tensor at each discrete material point.

Returns
The calculated elastic strain tensor

Definition at line 207 of file NodalRankTwoPD.C.

208 {
209  RankTwoTensor total_strain = computeNodalTotalStrain();
210  RankTwoTensor delta(RankTwoTensor::initIdentity);
211  RankTwoTensor thermal_strain;
212 
213  if (_has_temp)
214  thermal_strain = _alpha * (_temp_var->getNodalValue(*_current_node) - _temp_ref) * delta;
215 
216  RankTwoTensor mechanical_strain = total_strain - thermal_strain;
217 
218  return mechanical_strain;
219 }

Referenced by computeNodalStress(), and computeValue().

◆ computeNodalStress()

RankTwoTensor NodalRankTwoPD::computeNodalStress ( )
protectedvirtual

Function to compute the stress tensor at each discrete materials point.

Returns
The calculated stress tensor

Definition at line 222 of file NodalRankTwoPD.C.

223 {
225 }

Referenced by computeValue().

◆ computeNodalTotalStrain()

RankTwoTensor NodalRankTwoPD::computeNodalTotalStrain ( )
protectedvirtual

Function to compute the total strain tensor at each discrete material point.

Returns
The calculated total strain tensor

Definition at line 137 of file NodalRankTwoPD.C.

138 {
139  std::vector<dof_id_type> neighbors = _pdmesh.getNeighbors(_current_node->id());
140  std::vector<dof_id_type> bonds = _pdmesh.getBonds(_current_node->id());
141  Real horizon = _pdmesh.getHorizon(_current_node->id());
142 
143  // calculate the shape tensor and prepare the deformation gradient tensor
144  RankTwoTensor shape, dgrad, delta(RankTwoTensor::initIdentity);
145  shape.zero();
146  dgrad.zero();
147  if (_dim == 2)
148  shape(2, 2) = dgrad(2, 2) = 1.0;
149 
150  RealGradient origin_vector(3), current_vector(3);
151  origin_vector = 0.0;
152  current_vector = 0.0;
153 
154  for (unsigned int j = 0; j < neighbors.size(); ++j)
155  {
156  Node * node_j = _pdmesh.nodePtr(neighbors[j]);
157  Real vol_j = _pdmesh.getPDNodeVolume(neighbors[j]);
158  for (unsigned int k = 0; k < _dim; ++k)
159  {
160  origin_vector(k) = (*node_j)(k) - (*_pdmesh.nodePtr(_current_node->id()))(k);
161  current_vector(k) = origin_vector(k) + _disp_var[k]->getNodalValue(*node_j) -
162  _disp_var[k]->getNodalValue(*_current_node);
163  }
164  Real origin_length = origin_vector.norm();
165  Real bond_status_ij = _bond_status_var->getElementalValue(_pdmesh.elemPtr(bonds[j]));
166  for (unsigned int k = 0; k < _dim; ++k)
167  for (unsigned int l = 0; l < _dim; ++l)
168  {
169  shape(k, l) +=
170  vol_j * horizon / origin_length * origin_vector(k) * origin_vector(l) * bond_status_ij;
171  dgrad(k, l) +=
172  vol_j * horizon / origin_length * current_vector(k) * origin_vector(l) * bond_status_ij;
173  }
174  }
175 
176  // finalize the deformation gradient tensor
177  dgrad *= shape.inverse();
178 
179  // the green-lagrange strain tensor
180  RankTwoTensor total_strain = 0.5 * (dgrad.transpose() * dgrad - delta);
181 
183  total_strain(2, 2) = _scalar_out_of_plane_strain[0];
184  else if (_plane_stress)
185  {
186  if (_has_temp)
187  {
188  Real mstrain00 =
189  (total_strain(0, 0) - _alpha * (_temp_var->getNodalValue(*_current_node) - _temp_ref));
190  Real mstrain11 =
191  (total_strain(1, 1) - _alpha * (_temp_var->getNodalValue(*_current_node) - _temp_ref));
192  Real mstrain22 =
193  -(_Cijkl(2, 2, 0, 0) * mstrain00 + _Cijkl(2, 2, 1, 1) * mstrain11) / _Cijkl(2, 2, 2, 2);
194  total_strain(2, 2) =
195  mstrain22 + _alpha * (_temp_var->getNodalValue(*_current_node) - _temp_ref);
196  }
197  else
198  total_strain(2, 2) =
199  -(_Cijkl(2, 2, 0, 0) * total_strain(0, 0) + _Cijkl(2, 2, 1, 1) * total_strain(1, 1)) /
200  _Cijkl(2, 2, 2, 2);
201  }
202 
203  return total_strain;
204 }

Referenced by computeNodalMechanicalStrain(), and computeValue().

◆ computeValue()

Real NodalRankTwoPD::computeValue ( )
overrideprotected

Definition at line 111 of file NodalRankTwoPD.C.

112 {
113  Real val = 0.0;
114  RankTwoTensor tensor;
115  Point p = Point(0, 0, 0);
116  if (_rank_two_tensor == "total_strain")
117  tensor = computeNodalTotalStrain();
118  else if (_rank_two_tensor == "mechanical_strain")
119  tensor = computeNodalMechanicalStrain();
120  else if (_rank_two_tensor == "stress")
121  tensor = computeNodalStress();
122  else
123  mooseError("NodalRankTwoPD Error: Pass valid rank two tensor: total_strain, "
124  "mechanical_strain or stress");
125 
126  if (_output_type == "component")
127  val = RankTwoScalarTools::component(tensor, _i, _j);
128  else if (_output_type == "scalar")
129  val = RankTwoScalarTools::getQuantity(tensor, _scalar_type, _point1, _point2, _q_point[_qp], p);
130  else
131  mooseError("NodalRankTwoPD Error: Pass valid output type - component or scalar");
132 
133  return val;
134 }

Member Data Documentation

◆ _alpha

const Real NodalRankTwoPD::_alpha
protected

Definition at line 71 of file NodalRankTwoPD.h.

Referenced by computeNodalMechanicalStrain(), and computeNodalTotalStrain().

◆ _bond_status_var

MooseVariable* NodalRankTwoPD::_bond_status_var
protected

bond_status variable

Definition at line 58 of file NodalRankTwoPD.h.

Referenced by computeNodalTotalStrain().

◆ _Cijkl

RankFourTensor NodalRankTwoPD::_Cijkl
protected

elasticity tensor

Definition at line 100 of file NodalRankTwoPD.h.

Referenced by computeNodalStress(), computeNodalTotalStrain(), and NodalRankTwoPD().

◆ _dim

const unsigned int AuxKernelBasePD::_dim
protectedinherited

Problem dimension.

Definition at line 33 of file AuxKernelBasePD.h.

Referenced by computeNodalTotalStrain(), and NodalRankTwoPD().

◆ _disp_var

std::vector<MooseVariable *> NodalRankTwoPD::_disp_var
protected

@

displacement variables

Definition at line 97 of file NodalRankTwoPD.h.

Referenced by computeNodalTotalStrain(), and NodalRankTwoPD().

◆ _has_temp

bool NodalRankTwoPD::_has_temp
protected

◆ _i

const unsigned int NodalRankTwoPD::_i
protected

component index

Definition at line 87 of file NodalRankTwoPD.h.

Referenced by computeValue().

◆ _j

const unsigned int NodalRankTwoPD::_j
protected

Definition at line 88 of file NodalRankTwoPD.h.

Referenced by computeValue().

◆ _output_type

std::string NodalRankTwoPD::_output_type
protected

output type: component or scalar

Definition at line 81 of file NodalRankTwoPD.h.

Referenced by computeValue(), and NodalRankTwoPD().

◆ _pdmesh

PeridynamicsMesh& AuxKernelBasePD::_pdmesh
protectedinherited

Reference to peridynamic mesh object.

Definition at line 30 of file AuxKernelBasePD.h.

Referenced by computeNodalTotalStrain(), BoundaryOffsetPD::computeValue(), and NodalVolumePD::computeValue().

◆ _plane_stress

const bool NodalRankTwoPD::_plane_stress
protected

plane stress problem or not

Definition at line 66 of file NodalRankTwoPD.h.

Referenced by computeNodalTotalStrain().

◆ _point1

const Point NodalRankTwoPD::_point1
protected

Points for direction dependent scalar output.

Definition at line 92 of file NodalRankTwoPD.h.

Referenced by computeValue().

◆ _point2

const Point NodalRankTwoPD::_point2
protected

Definition at line 93 of file NodalRankTwoPD.h.

Referenced by computeValue().

◆ _poissons_ratio

const Real NodalRankTwoPD::_poissons_ratio
protected

Definition at line 70 of file NodalRankTwoPD.h.

Referenced by NodalRankTwoPD().

◆ _rank_two_tensor

std::string NodalRankTwoPD::_rank_two_tensor
protected

name of rank two tensor to be processed: total_strain, mechanical_strain or stress

Definition at line 78 of file NodalRankTwoPD.h.

Referenced by computeValue(), and NodalRankTwoPD().

◆ _scalar_out_of_plane_strain

const VariableValue& NodalRankTwoPD::_scalar_out_of_plane_strain
protected

Definition at line 62 of file NodalRankTwoPD.h.

Referenced by computeNodalTotalStrain().

◆ _scalar_out_of_plane_strain_coupled

bool NodalRankTwoPD::_scalar_out_of_plane_strain_coupled
protected

variable for generalized plane strain cases

Definition at line 61 of file NodalRankTwoPD.h.

Referenced by computeNodalTotalStrain().

◆ _scalar_type

MooseEnum NodalRankTwoPD::_scalar_type
protected

specific scalar type to be output

Definition at line 84 of file NodalRankTwoPD.h.

Referenced by computeValue().

◆ _temp_ref

const Real NodalRankTwoPD::_temp_ref
protected

reference temperature

Definition at line 75 of file NodalRankTwoPD.h.

Referenced by computeNodalMechanicalStrain(), and computeNodalTotalStrain().

◆ _temp_var

MooseVariable* NodalRankTwoPD::_temp_var
protected

coupled temperature variable

Definition at line 55 of file NodalRankTwoPD.h.

Referenced by computeNodalMechanicalStrain(), and computeNodalTotalStrain().

◆ _youngs_modulus

const Real NodalRankTwoPD::_youngs_modulus
protected

material constants

Definition at line 69 of file NodalRankTwoPD.h.

Referenced by NodalRankTwoPD().


The documentation for this class was generated from the following files:
NodalRankTwoPD::_poissons_ratio
const Real _poissons_ratio
Definition: NodalRankTwoPD.h:70
PeridynamicsMesh::getPDNodeVolume
Real getPDNodeVolume(dof_id_type node_id)
Function to return nodal volume for node node_id.
Definition: PeridynamicsMesh.C:435
NodalRankTwoPD::_bond_status_var
MooseVariable * _bond_status_var
bond_status variable
Definition: NodalRankTwoPD.h:58
AuxKernelBasePD::_pdmesh
PeridynamicsMesh & _pdmesh
Reference to peridynamic mesh object.
Definition: AuxKernelBasePD.h:30
NodalRankTwoPD::computeNodalTotalStrain
virtual RankTwoTensor computeNodalTotalStrain()
Function to compute the total strain tensor at each discrete material point.
Definition: NodalRankTwoPD.C:137
libMesh::RealGradient
VectorValue< Real > RealGradient
Definition: GrainForceAndTorqueInterface.h:17
NodalRankTwoPD::_has_temp
bool _has_temp
Definition: NodalRankTwoPD.h:52
NodalRankTwoPD::_j
const unsigned int _j
Definition: NodalRankTwoPD.h:88
AuxKernelBasePD::AuxKernelBasePD
AuxKernelBasePD(const InputParameters &parameters)
Definition: AuxKernelBasePD.C:23
RankTwoScalarTools::getQuantity
T getQuantity(const RankTwoTensorTempl< T > &tensor, const MooseEnum &scalar_type, const Point &point1, const Point &point2, const Point &curr_point, Point &direction)
Definition: RankTwoScalarTools.h:420
NodalRankTwoPD::computeNodalStress
virtual RankTwoTensor computeNodalStress()
Function to compute the stress tensor at each discrete materials point.
Definition: NodalRankTwoPD.C:222
NodalRankTwoPD::_i
const unsigned int _i
component index
Definition: NodalRankTwoPD.h:87
NodalRankTwoPD::_temp_var
MooseVariable * _temp_var
coupled temperature variable
Definition: NodalRankTwoPD.h:55
NodalRankTwoPD::_Cijkl
RankFourTensor _Cijkl
elasticity tensor
Definition: NodalRankTwoPD.h:100
NodalRankTwoPD::_disp_var
std::vector< MooseVariable * > _disp_var
@
Definition: NodalRankTwoPD.h:97
NodalRankTwoPD::_alpha
const Real _alpha
Definition: NodalRankTwoPD.h:71
RankTwoScalarTools::component
T component(const RankTwoTensorTempl< T > &r2tensor, unsigned int i, unsigned int j)
Definition: RankTwoScalarTools.h:31
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
PeridynamicsMesh::getHorizon
Real getHorizon(dof_id_type node_id)
Function to return horizon size.
Definition: PeridynamicsMesh.C:471
NodalRankTwoPD::_plane_stress
const bool _plane_stress
plane stress problem or not
Definition: NodalRankTwoPD.h:66
PeridynamicsMesh::getBonds
std::vector< dof_id_type > getBonds(dof_id_type node_id)
Function to return the bond number connected with node node_id.
Definition: PeridynamicsMesh.C:382
NodalRankTwoPD::_point1
const Point _point1
Points for direction dependent scalar output.
Definition: NodalRankTwoPD.h:92
NodalRankTwoPD::_youngs_modulus
const Real _youngs_modulus
material constants
Definition: NodalRankTwoPD.h:69
NodalRankTwoPD::_point2
const Point _point2
Definition: NodalRankTwoPD.h:93
NodalRankTwoPD::_scalar_type
MooseEnum _scalar_type
specific scalar type to be output
Definition: NodalRankTwoPD.h:84
AuxKernelBasePD::_dim
const unsigned int _dim
Problem dimension.
Definition: AuxKernelBasePD.h:33
NodalRankTwoPD::computeNodalMechanicalStrain
virtual RankTwoTensor computeNodalMechanicalStrain()
Function to compute the elastic strain tensor at each discrete material point.
Definition: NodalRankTwoPD.C:207
RankTwoTensorTempl< Real >
NodalRankTwoPD::_output_type
std::string _output_type
output type: component or scalar
Definition: NodalRankTwoPD.h:81
NodalRankTwoPD::_scalar_out_of_plane_strain
const VariableValue & _scalar_out_of_plane_strain
Definition: NodalRankTwoPD.h:62
NodalRankTwoPD::_scalar_out_of_plane_strain_coupled
bool _scalar_out_of_plane_strain_coupled
variable for generalized plane strain cases
Definition: NodalRankTwoPD.h:61
NodalRankTwoPD::_temp_ref
const Real _temp_ref
reference temperature
Definition: NodalRankTwoPD.h:75
NodalRankTwoPD::_rank_two_tensor
std::string _rank_two_tensor
name of rank two tensor to be processed: total_strain, mechanical_strain or stress
Definition: NodalRankTwoPD.h:78