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

ComputeCrackTipEnrichmentSmallStrain calculates the sum of standard strain and enrichement strain. More...

#include <ComputeCrackTipEnrichmentSmallStrain.h>

Inheritance diagram for ComputeCrackTipEnrichmentSmallStrain:
[legend]

Public Member Functions

 ComputeCrackTipEnrichmentSmallStrain (const InputParameters &parameters)
 
virtual ~ComputeCrackTipEnrichmentSmallStrain ()
 
void initialSetup () override
 
virtual unsigned int crackTipEnrichementFunctionAtPoint (const Point &point, std::vector< Real > &B)
 calculate the enrichment function values at point More...
 
virtual unsigned int crackTipEnrichementFunctionDerivativeAtPoint (const Point &point, std::vector< RealVectorValue > &dB)
 calculate the enrichment function derivatives at point More...
 
void rotateFromCrackFrontCoordsToGlobal (const RealVectorValue &vector, RealVectorValue &rotated_vector, const unsigned int point_index)
 rotate a vector from crack front coordinate to global cooridate More...
 

Static Public Member Functions

static InputParameters validParams ()
 

Protected Member Functions

virtual void computeProperties () override
 
virtual void computeQpProperties () override
 
virtual void initQpStatefulProperties () override
 
virtual void displacementIntegrityCheck ()
 

Protected Attributes

std::vector< Real > _enrich_disp
 enrichment displacement More...
 
std::vector< RealVectorValue > _grad_enrich_disp
 gradient of enrichment displacement More...
 
std::vector< std::vector< MooseVariableFEBase * > > _enrich_variable
 enrichment displacement variables More...
 
const VariablePhiValue & _phi
 the current shape functions More...
 
const VariablePhiGradient & _grad_phi
 gradient of the shape function More...
 
unsigned int _ndisp
 Coupled displacement variables. More...
 
std::vector< const VariableValue * > _disp
 
std::vector< const VariableGradient * > _grad_disp
 
const std::string _base_name
 
MaterialProperty< RankTwoTensor > & _mechanical_strain
 
MaterialProperty< RankTwoTensor > & _total_strain
 
std::vector< MaterialPropertyName > _eigenstrain_names
 
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
 
const MaterialProperty< RankTwoTensor > * _global_strain
 
const bool _volumetric_locking_correction
 
const Real & _current_elem_volume
 

Private Attributes

std::vector< Real > _B
 enrichment function value More...
 
std::vector< RealVectorValue > _dBX
 derivatives of enrichment function respect to global cooridnate More...
 
std::vector< RealVectorValue > _dBx
 derivatives of enrichment function respect to crack front cooridnate More...
 
std::vector< std::vector< Real > > _BI
 enrichment function at node I More...
 
const std::vector< std::vector< Real > > * _fe_phi
 shape function More...
 
const std::vector< std::vector< RealGradient > > * _fe_dphi
 gradient of shape function More...
 
NonlinearSystem * _nl
 
const NumericVector< Number > * _sln
 
const CrackFrontDefinition_crack_front_definition
 
Real _r
 
Real _theta
 

Detailed Description

ComputeCrackTipEnrichmentSmallStrain calculates the sum of standard strain and enrichement strain.

Definition at line 29 of file ComputeCrackTipEnrichmentSmallStrain.h.

Constructor & Destructor Documentation

◆ ComputeCrackTipEnrichmentSmallStrain()

ComputeCrackTipEnrichmentSmallStrain::ComputeCrackTipEnrichmentSmallStrain ( const InputParameters &  parameters)

Definition at line 33 of file ComputeCrackTipEnrichmentSmallStrain.C.

35  : ComputeStrainBase(parameters),
36  EnrichmentFunctionCalculation(&getUserObject<CrackFrontDefinition>("crack_front_definition")),
37  _enrich_disp(3),
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 }

◆ ~ComputeCrackTipEnrichmentSmallStrain()

virtual ComputeCrackTipEnrichmentSmallStrain::~ComputeCrackTipEnrichmentSmallStrain ( )
inlinevirtual

Definition at line 34 of file ComputeCrackTipEnrichmentSmallStrain.h.

34 {}

Member Function Documentation

◆ computeProperties()

void ComputeCrackTipEnrichmentSmallStrain::computeProperties ( )
overrideprotectedvirtual

Definition at line 123 of file ComputeCrackTipEnrichmentSmallStrain.C.

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 }

◆ computeQpProperties()

void ComputeCrackTipEnrichmentSmallStrain::computeQpProperties ( )
overrideprotectedvirtual

Definition at line 74 of file ComputeCrackTipEnrichmentSmallStrain.C.

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 }

◆ crackTipEnrichementFunctionAtPoint()

unsigned int EnrichmentFunctionCalculation::crackTipEnrichementFunctionAtPoint ( const Point &  point,
std::vector< Real > &  B 
)
virtualinherited

calculate the enrichment function values at point

Returns
the closest crack front index

Definition at line 19 of file EnrichmentFunctionCalculation.C.

21 {
22  unsigned int crack_front_point_index =
24 
25  if (MooseUtils::absoluteFuzzyEqual(_r, 0.0))
26  mooseError("EnrichmentFunctionCalculation: the distance between a point and the crack "
27  "tip/front is zero.");
28 
29  Real st = std::sin(_theta);
30  Real st2 = std::sin(_theta / 2.0);
31  Real ct2 = std::cos(_theta / 2.0);
32  Real sr = std::sqrt(_r);
33 
34  B[0] = sr * st2;
35  B[1] = sr * ct2;
36  B[2] = sr * st2 * st;
37  B[3] = sr * ct2 * st;
38 
39  return crack_front_point_index;
40 }

Referenced by computeProperties(), CrackTipEnrichmentStressDivergenceTensors::computeQpJacobian(), CrackTipEnrichmentStressDivergenceTensors::computeQpOffDiagJacobian(), computeQpProperties(), and CrackTipEnrichmentStressDivergenceTensors::computeQpResidual().

◆ crackTipEnrichementFunctionDerivativeAtPoint()

unsigned int EnrichmentFunctionCalculation::crackTipEnrichementFunctionDerivativeAtPoint ( const Point &  point,
std::vector< RealVectorValue > &  dB 
)
virtualinherited

calculate the enrichment function derivatives at point

Returns
the closest crack front index

Definition at line 43 of file EnrichmentFunctionCalculation.C.

45 {
46  unsigned int crack_front_point_index =
48 
49  if (MooseUtils::absoluteFuzzyEqual(_r, 0.0))
50  mooseError("EnrichmentFunctionCalculation: the distance between a point and the crack "
51  "tip/front is zero.");
52 
53  Real st = std::sin(_theta);
54  Real ct = std::cos(_theta);
55  Real st2 = std::sin(_theta / 2.0);
56  Real ct2 = std::cos(_theta / 2.0);
57  Real st15 = std::sin(1.5 * _theta);
58  Real ct15 = std::cos(1.5 * _theta);
59  Real sr = std::sqrt(_r);
60 
61  dB[0](0) = -0.5 / sr * st2;
62  dB[0](1) = 0.5 / sr * ct2;
63  dB[0](2) = 0.0;
64  dB[1](0) = 0.5 / sr * ct2;
65  dB[1](1) = 0.5 / sr * st2;
66  dB[1](2) = 0.0;
67  dB[2](0) = -0.5 / sr * st15 * st;
68  dB[2](1) = 0.5 / sr * (st2 + st15 * ct);
69  dB[2](2) = 0.0;
70  dB[3](0) = -0.5 / sr * ct15 * st;
71  dB[3](1) = 0.5 / sr * (ct2 + ct15 * ct);
72  dB[3](2) = 0.0;
73 
74  return crack_front_point_index;
75 }

Referenced by CrackTipEnrichmentStressDivergenceTensors::computeQpJacobian(), CrackTipEnrichmentStressDivergenceTensors::computeQpOffDiagJacobian(), computeQpProperties(), and CrackTipEnrichmentStressDivergenceTensors::computeQpResidual().

◆ displacementIntegrityCheck()

void ComputeStrainBase::displacementIntegrityCheck ( )
protectedvirtualinherited

Reimplemented in Compute2DFiniteStrain, Compute2DIncrementalStrain, and Compute2DSmallStrain.

Definition at line 94 of file ComputeStrainBase.C.

95 {
96  // Checking for consistency between mesh size and length of the provided displacements vector
97  if (_ndisp != _mesh.dimension())
98  paramError(
99  "displacements",
100  "The number of variables supplied in 'displacements' must match the mesh dimension.");
101 }

Referenced by ComputeStrainBase::initialSetup().

◆ initialSetup()

void ComputeStrainBase::initialSetup ( )
overrideinherited

Definition at line 75 of file ComputeStrainBase.C.

76 {
78  // fetch coupled variables and gradients (as stateful properties if necessary)
79  for (unsigned int i = 0; i < _ndisp; ++i)
80  {
81  _disp[i] = &coupledValue("displacements", i);
82  _grad_disp[i] = &coupledGradient("displacements", i);
83  }
84 
85  // set unused dimensions to zero
86  for (unsigned i = _ndisp; i < 3; ++i)
87  {
88  _disp[i] = &_zero;
89  _grad_disp[i] = &_grad_zero;
90  }
91 }

Referenced by ComputeIncrementalStrainBase::initialSetup(), and ComputeAxisymmetric1DSmallStrain::initialSetup().

◆ initQpStatefulProperties()

void ComputeStrainBase::initQpStatefulProperties ( )
overrideprotectedvirtualinherited

Reimplemented in ComputeIncrementalStrainBase, and ComputeCosseratIncrementalSmallStrain.

Definition at line 104 of file ComputeStrainBase.C.

105 {
106  _mechanical_strain[_qp].zero();
107  _total_strain[_qp].zero();
108 }

◆ rotateFromCrackFrontCoordsToGlobal()

void EnrichmentFunctionCalculation::rotateFromCrackFrontCoordsToGlobal ( const RealVectorValue &  vector,
RealVectorValue &  rotated_vector,
const unsigned int  point_index 
)
inherited

rotate a vector from crack front coordinate to global cooridate

Parameters
rotated_vectorrotated vector

Definition at line 78 of file EnrichmentFunctionCalculation.C.

81 {
82  rotated_vector = _crack_front_definition->rotateFromCrackFrontCoordsToGlobal(vector, point_index);
83 }

Referenced by CrackTipEnrichmentStressDivergenceTensors::computeQpJacobian(), CrackTipEnrichmentStressDivergenceTensors::computeQpOffDiagJacobian(), computeQpProperties(), and CrackTipEnrichmentStressDivergenceTensors::computeQpResidual().

◆ validParams()

InputParameters ComputeStrainBase::validParams ( )
staticinherited

Definition at line 17 of file ComputeStrainBase.C.

18 {
19  InputParameters params = Material::validParams();
20  params.addRequiredCoupledVar(
21  "displacements",
22  "The displacements appropriate for the simulation geometry and coordinate system");
23  params.addParam<std::string>("base_name",
24  "Optional parameter that allows the user to define "
25  "multiple mechanics material systems on the same "
26  "block, i.e. for multiple phases");
27  params.addParam<bool>(
28  "volumetric_locking_correction", false, "Flag to correct volumetric locking");
29  params.addParam<std::vector<MaterialPropertyName>>(
30  "eigenstrain_names", "List of eigenstrains to be applied in this strain calculation");
31  params.addParam<MaterialPropertyName>("global_strain",
32  "Optional material property holding a global strain "
33  "tensor applied to the mesh as a whole");
34  params.suppressParameter<bool>("use_displaced_mesh");
35  return params;
36 }

Referenced by ComputeCosseratSmallStrain::validParams(), ComputeSmallStrain::validParams(), and ComputeIncrementalStrainBase::validParams().

Member Data Documentation

◆ _B

std::vector<Real> ComputeCrackTipEnrichmentSmallStrain::_B
private

enrichment function value

Definition at line 58 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by computeQpProperties().

◆ _base_name

const std::string ComputeStrainBase::_base_name
protectedinherited

Definition at line 44 of file ComputeStrainBase.h.

Referenced by ComputeStrainBase::ComputeStrainBase().

◆ _BI

std::vector<std::vector<Real> > ComputeCrackTipEnrichmentSmallStrain::_BI
private

enrichment function at node I

Definition at line 64 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by ComputeCrackTipEnrichmentSmallStrain(), computeProperties(), and computeQpProperties().

◆ _crack_front_definition

const CrackFrontDefinition* EnrichmentFunctionCalculation::_crack_front_definition
privateinherited

◆ _current_elem_volume

const Real& ComputeStrainBase::_current_elem_volume
protectedinherited

◆ _dBX

std::vector<RealVectorValue> ComputeCrackTipEnrichmentSmallStrain::_dBX
private

derivatives of enrichment function respect to global cooridnate

Definition at line 60 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by computeQpProperties().

◆ _dBx

std::vector<RealVectorValue> ComputeCrackTipEnrichmentSmallStrain::_dBx
private

derivatives of enrichment function respect to crack front cooridnate

Definition at line 62 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by computeQpProperties().

◆ _disp

std::vector<const VariableValue *> ComputeStrainBase::_disp
protectedinherited

◆ _eigenstrain_names

std::vector<MaterialPropertyName> ComputeStrainBase::_eigenstrain_names
protectedinherited

◆ _eigenstrains

std::vector<const MaterialProperty<RankTwoTensor> *> ComputeStrainBase::_eigenstrains
protectedinherited

◆ _enrich_disp

std::vector<Real> ComputeCrackTipEnrichmentSmallStrain::_enrich_disp
protected

enrichment displacement

Definition at line 42 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by computeQpProperties().

◆ _enrich_variable

std::vector<std::vector<MooseVariableFEBase *> > ComputeCrackTipEnrichmentSmallStrain::_enrich_variable
protected

enrichment displacement variables

Definition at line 48 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by ComputeCrackTipEnrichmentSmallStrain(), and computeQpProperties().

◆ _fe_dphi

const std::vector<std::vector<RealGradient> >* ComputeCrackTipEnrichmentSmallStrain::_fe_dphi
private

gradient of shape function

Definition at line 68 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by computeProperties().

◆ _fe_phi

const std::vector<std::vector<Real> >* ComputeCrackTipEnrichmentSmallStrain::_fe_phi
private

shape function

Definition at line 66 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by computeProperties(), and computeQpProperties().

◆ _global_strain

const MaterialProperty<RankTwoTensor>* ComputeStrainBase::_global_strain
protectedinherited

◆ _grad_disp

std::vector<const VariableGradient *> ComputeStrainBase::_grad_disp
protectedinherited

◆ _grad_enrich_disp

std::vector<RealVectorValue> ComputeCrackTipEnrichmentSmallStrain::_grad_enrich_disp
protected

gradient of enrichment displacement

Definition at line 45 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by computeQpProperties().

◆ _grad_phi

const VariablePhiGradient& ComputeCrackTipEnrichmentSmallStrain::_grad_phi
protected

gradient of the shape function

Definition at line 54 of file ComputeCrackTipEnrichmentSmallStrain.h.

◆ _mechanical_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBase::_mechanical_strain
protectedinherited

◆ _ndisp

unsigned int ComputeStrainBase::_ndisp
protectedinherited

◆ _nl

NonlinearSystem* ComputeCrackTipEnrichmentSmallStrain::_nl
private

◆ _phi

const VariablePhiValue& ComputeCrackTipEnrichmentSmallStrain::_phi
protected

the current shape functions

Definition at line 51 of file ComputeCrackTipEnrichmentSmallStrain.h.

◆ _r

Real EnrichmentFunctionCalculation::_r
privateinherited

◆ _sln

const NumericVector<Number>* ComputeCrackTipEnrichmentSmallStrain::_sln
private

Definition at line 70 of file ComputeCrackTipEnrichmentSmallStrain.h.

Referenced by computeQpProperties().

◆ _theta

Real EnrichmentFunctionCalculation::_theta
privateinherited

◆ _total_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBase::_total_strain
protectedinherited

◆ _volumetric_locking_correction

const bool ComputeStrainBase::_volumetric_locking_correction
protectedinherited

The documentation for this class was generated from the following files:
ComputeCrackTipEnrichmentSmallStrain::_phi
const VariablePhiValue & _phi
the current shape functions
Definition: ComputeCrackTipEnrichmentSmallStrain.h:51
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
ComputeStrainBase::ComputeStrainBase
ComputeStrainBase(const InputParameters &parameters)
Definition: ComputeStrainBase.C:38
ComputeStrainBase::displacementIntegrityCheck
virtual void displacementIntegrityCheck()
Definition: ComputeStrainBase.C:94
ComputeStrainBase::_disp
std::vector< const VariableValue * > _disp
Definition: ComputeStrainBase.h:41
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
EnrichmentFunctionCalculation::_theta
Real _theta
Definition: EnrichmentFunctionCalculation.h:50
ComputeCrackTipEnrichmentSmallStrain::_enrich_variable
std::vector< std::vector< MooseVariableFEBase * > > _enrich_variable
enrichment displacement variables
Definition: ComputeCrackTipEnrichmentSmallStrain.h:48
EnrichmentFunctionCalculation::_r
Real _r
Definition: EnrichmentFunctionCalculation.h:49
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
ComputeCrackTipEnrichmentSmallStrain::_grad_phi
const VariablePhiGradient & _grad_phi
gradient of the shape function
Definition: ComputeCrackTipEnrichmentSmallStrain.h:54
ComputeStrainBase::_mechanical_strain
MaterialProperty< RankTwoTensor > & _mechanical_strain
Definition: ComputeStrainBase.h:46
CrackFrontDefinition::rotateFromCrackFrontCoordsToGlobal
RealVectorValue rotateFromCrackFrontCoordsToGlobal(const RealVectorValue vector, const unsigned int point_index) const
rotate a vector from crack front cartesian coordinate to global cartesian coordinate
Definition: CrackFrontDefinition.C:1133
CrackFrontDefinition::calculateRThetaToCrackFront
void calculateRThetaToCrackFront(const Point qp, const unsigned int point_index, Real &r, Real &theta) const
calculate r and theta in the crack front polar coordinates
Definition: CrackFrontDefinition.C:1150
validParams
InputParameters validParams()
EnrichmentFunctionCalculation::_crack_front_definition
const CrackFrontDefinition * _crack_front_definition
Definition: EnrichmentFunctionCalculation.h:48
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
ComputeStrainBase::_eigenstrains
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
Definition: ComputeStrainBase.h:51
EnrichmentFunctionCalculation::EnrichmentFunctionCalculation
EnrichmentFunctionCalculation(const CrackFrontDefinition *crack_front_definition)
Definition: EnrichmentFunctionCalculation.C:12
EnrichmentFunctionCalculation::crackTipEnrichementFunctionAtPoint
virtual unsigned int crackTipEnrichementFunctionAtPoint(const Point &point, std::vector< Real > &B)
calculate the enrichment function values at point
Definition: EnrichmentFunctionCalculation.C:19
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::_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