www.mooseframework.org
ComputeCrackTipEnrichmentSmallStrain.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 "MooseMesh.h"
12 #include "libmesh/fe_interface.h"
13 #include "libmesh/string_to_enum.h"
14 
15 #include "libmesh/quadrature_gauss.h"
16 
18 
19 template <>
20 InputParameters
22 {
23  InputParameters params = validParams<ComputeStrainBase>();
24  params.addClassDescription(
25  "Computes the crack tip enrichment at a point within a small strain formulation.");
26  params.addRequiredParam<std::vector<NonlinearVariableName>>("enrichment_displacements",
27  "The enrichment displacement");
28  params.addRequiredParam<UserObjectName>("crack_front_definition",
29  "The CrackFrontDefinition user object name");
30  return params;
31 }
32 
34  const InputParameters & parameters)
35  : ComputeStrainBase(parameters),
36  EnrichmentFunctionCalculation(&getUserObject<CrackFrontDefinition>("crack_front_definition")),
37  _enrich_disp(3),
38  _grad_enrich_disp(3),
39  _enrich_variable(4),
40  _phi(_assembly.phi()),
41  _grad_phi(_assembly.gradPhi()),
42  _B(4),
43  _dBX(4),
44  _dBx(4)
45 {
46  for (unsigned int i = 0; i < _enrich_variable.size(); ++i)
47  _enrich_variable[i].resize(_ndisp);
48 
49  const std::vector<NonlinearVariableName> & nl_vnames =
50  getParam<std::vector<NonlinearVariableName>>("enrichment_displacements");
51 
52  if (_ndisp == 2 && nl_vnames.size() != 8)
53  mooseError("The number of enrichment displacements should be total 8 for 2D.");
54  else if (_ndisp == 3 && nl_vnames.size() != 12)
55  mooseError("The number of enrichment displacements should be total 12 for 3D.");
56 
57  NonlinearSystem & nl = _fe_problem.getNonlinearSystem();
58  _nl = &(_fe_problem.getNonlinearSystem());
59 
60  for (unsigned int j = 0; j < _ndisp; ++j)
61  for (unsigned int i = 0; i < 4; ++i)
62  _enrich_variable[i][j] = &(nl.getVariable(0, nl_vnames[j * 4 + i]));
63 
64  if (_ndisp == 2)
65  _BI.resize(4); // QUAD4
66  else if (_ndisp == 3)
67  _BI.resize(8); // HEX8
68 
69  for (unsigned int i = 0; i < _BI.size(); ++i)
70  _BI[i].resize(4);
71 }
72 
73 void
75 {
76  crackTipEnrichementFunctionAtPoint(_q_point[_qp], _B);
77  unsigned int crack_front_point_index =
79 
80  for (unsigned int i = 0; i < 4; ++i)
81  rotateFromCrackFrontCoordsToGlobal(_dBx[i], _dBX[i], crack_front_point_index);
82 
83  _sln = _nl->currentSolution();
84 
85  for (unsigned int m = 0; m < _ndisp; ++m)
86  {
87  _enrich_disp[m] = 0.0;
88  _grad_enrich_disp[m].zero();
89  for (unsigned int i = 0; i < _current_elem->n_nodes(); ++i)
90  {
91  const Node * node_i = _current_elem->node_ptr(i);
92  for (unsigned int j = 0; j < 4; ++j)
93  {
94  dof_id_type dof = node_i->dof_number(_nl->number(), _enrich_variable[j][m]->number(), 0);
95  Real soln = (*_sln)(dof);
96  _enrich_disp[m] += (*_fe_phi)[i][_qp] * (_B[j] - _BI[i][j]) * soln;
97  RealVectorValue grad_B(_dBX[j]);
98  _grad_enrich_disp[m] +=
99  ((*_fe_dphi)[i][_qp] * (_B[j] - _BI[i][j]) + (*_fe_phi)[i][_qp] * grad_B) * soln;
100  }
101  }
102  }
103 
104  RankTwoTensor grad_tensor_enrich(
106 
107  RankTwoTensor enrich_strain = (grad_tensor_enrich + grad_tensor_enrich.transpose()) / 2.0;
108 
109  RankTwoTensor grad_tensor((*_grad_disp[0])[_qp], (*_grad_disp[1])[_qp], (*_grad_disp[2])[_qp]);
110 
111  _total_strain[_qp] = (grad_tensor + grad_tensor.transpose()) / 2.0;
112 
113  _total_strain[_qp] += enrich_strain;
114 
115  _mechanical_strain[_qp] = _total_strain[_qp];
116 
117  // Remove the Eigen strain
118  for (auto es : _eigenstrains)
119  _mechanical_strain[_qp] -= (*es)[_qp];
120 }
121 
122 void
124 {
125  FEType fe_type(Utility::string_to_enum<Order>("first"),
126  Utility::string_to_enum<FEFamily>("lagrange"));
127  const unsigned int dim = _current_elem->dim();
128  std::unique_ptr<FEBase> fe(FEBase::build(dim, fe_type));
129  fe->attach_quadrature_rule(const_cast<QBase *>(_qrule));
130  _fe_phi = &(fe->get_phi());
131  _fe_dphi = &(fe->get_dphi());
132 
133  if (isBoundaryMaterial())
134  fe->reinit(_current_elem, _current_side);
135  else
136  fe->reinit(_current_elem);
137 
138  for (unsigned int i = 0; i < _BI.size(); ++i)
139  crackTipEnrichementFunctionAtPoint(*(_current_elem->node_ptr(i)), _BI[i]);
140 
141  ComputeStrainBase::computeProperties();
142 }
EnrichmentFunctionCalculation::crackTipEnrichementFunctionDerivativeAtPoint
virtual unsigned int crackTipEnrichementFunctionDerivativeAtPoint(const Point &point, std::vector< RealVectorValue > &dB)
calculate the enrichment function derivatives at point
Definition: EnrichmentFunctionCalculation.C:43
ComputeCrackTipEnrichmentSmallStrain::_nl
NonlinearSystem * _nl
Definition: ComputeCrackTipEnrichmentSmallStrain.h:69
registerMooseObject
registerMooseObject("XFEMApp", ComputeCrackTipEnrichmentSmallStrain)
ComputeCrackTipEnrichmentSmallStrain::_BI
std::vector< std::vector< Real > > _BI
enrichment function at node I
Definition: ComputeCrackTipEnrichmentSmallStrain.h:64
ComputeCrackTipEnrichmentSmallStrain::_sln
const NumericVector< Number > * _sln
Definition: ComputeCrackTipEnrichmentSmallStrain.h:70
CrackFrontDefinition
Works on top of NodalNormalsPreprocessor.
Definition: CrackFrontDefinition.h:36
validParams< ComputeCrackTipEnrichmentSmallStrain >
InputParameters validParams< ComputeCrackTipEnrichmentSmallStrain >()
Definition: ComputeCrackTipEnrichmentSmallStrain.C:21
ComputeStrainBase
ComputeStrainBase is the base class for strain tensors.
Definition: ComputeStrainBase.h:26
ComputeCrackTipEnrichmentSmallStrain::_enrich_variable
std::vector< std::vector< MooseVariableFEBase * > > _enrich_variable
enrichment displacement variables
Definition: ComputeCrackTipEnrichmentSmallStrain.h:48
ComputeCrackTipEnrichmentSmallStrain::computeProperties
virtual void computeProperties() override
Definition: ComputeCrackTipEnrichmentSmallStrain.C:123
ComputeCrackTipEnrichmentSmallStrain::_grad_enrich_disp
std::vector< RealVectorValue > _grad_enrich_disp
gradient of enrichment displacement
Definition: ComputeCrackTipEnrichmentSmallStrain.h:45
ComputeCrackTipEnrichmentSmallStrain::_fe_dphi
const std::vector< std::vector< RealGradient > > * _fe_dphi
gradient of shape function
Definition: ComputeCrackTipEnrichmentSmallStrain.h:68
ComputeStrainBase::_ndisp
unsigned int _ndisp
Coupled displacement variables.
Definition: ComputeStrainBase.h:40
ComputeCrackTipEnrichmentSmallStrain::_B
std::vector< Real > _B
enrichment function value
Definition: ComputeCrackTipEnrichmentSmallStrain.h:58
ComputeStrainBase::_mechanical_strain
MaterialProperty< RankTwoTensor > & _mechanical_strain
Definition: ComputeStrainBase.h:46
ComputeCrackTipEnrichmentSmallStrain::ComputeCrackTipEnrichmentSmallStrain
ComputeCrackTipEnrichmentSmallStrain(const InputParameters &parameters)
Definition: ComputeCrackTipEnrichmentSmallStrain.C:33
ComputeCrackTipEnrichmentSmallStrain::_fe_phi
const std::vector< std::vector< Real > > * _fe_phi
shape function
Definition: ComputeCrackTipEnrichmentSmallStrain.h:66
ComputeCrackTipEnrichmentSmallStrain::_dBx
std::vector< RealVectorValue > _dBx
derivatives of enrichment function respect to crack front cooridnate
Definition: ComputeCrackTipEnrichmentSmallStrain.h:62
ComputeCrackTipEnrichmentSmallStrain
ComputeCrackTipEnrichmentSmallStrain calculates the sum of standard strain and enrichement strain.
Definition: ComputeCrackTipEnrichmentSmallStrain.h:29
ComputeStrainBase::_eigenstrains
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
Definition: ComputeStrainBase.h:51
EnrichmentFunctionCalculation::crackTipEnrichementFunctionAtPoint
virtual unsigned int crackTipEnrichementFunctionAtPoint(const Point &point, std::vector< Real > &B)
calculate the enrichment function values at point
Definition: EnrichmentFunctionCalculation.C:19
validParams< ComputeStrainBase >
InputParameters validParams< ComputeStrainBase >()
RankTwoTensorTempl< Real >
ComputeStrainBase::_total_strain
MaterialProperty< RankTwoTensor > & _total_strain
Definition: ComputeStrainBase.h:48
EnrichmentFunctionCalculation::rotateFromCrackFrontCoordsToGlobal
void rotateFromCrackFrontCoordsToGlobal(const RealVectorValue &vector, RealVectorValue &rotated_vector, const unsigned int point_index)
rotate a vector from crack front coordinate to global cooridate
Definition: EnrichmentFunctionCalculation.C:78
ComputeStrainBase::_grad_disp
std::vector< const VariableGradient * > _grad_disp
Definition: ComputeStrainBase.h:42
ComputeCrackTipEnrichmentSmallStrain::computeQpProperties
virtual void computeQpProperties() override
Definition: ComputeCrackTipEnrichmentSmallStrain.C:74
ComputeCrackTipEnrichmentSmallStrain::_dBX
std::vector< RealVectorValue > _dBX
derivatives of enrichment function respect to global cooridnate
Definition: ComputeCrackTipEnrichmentSmallStrain.h:60
ComputeCrackTipEnrichmentSmallStrain::_enrich_disp
std::vector< Real > _enrich_disp
enrichment displacement
Definition: ComputeCrackTipEnrichmentSmallStrain.h:42
EnrichmentFunctionCalculation
Perform calculation of enrichment function values and derivatives.
Definition: EnrichmentFunctionCalculation.h:22
ComputeCrackTipEnrichmentSmallStrain.h