www.mooseframework.org
ComputeCrackTipEnrichmentSmallStrain.h
Go to the documentation of this file.
1 /****************************************************************/
2 /* MOOSE - Multiphysics Object Oriented Simulation Environment */
3 /* */
4 /* All contents are licensed under LGPL V2.1 */
5 /* See LICENSE for full restrictions */
6 /****************************************************************/
7 #ifndef COMPUTECRACKTIPENRICHMENTSMALLSTRAIN_H
8 #define COMPUTECRACKTIPENRICHMENTSMALLSTRAIN_H
9 
10 #include "Material.h"
11 #include "RankTwoTensor.h"
12 #include "RankFourTensor.h"
13 #include "RotationTensor.h"
14 #include "ComputeStrainBase.h"
15 #include "Assembly.h"
16 #include "CrackFrontDefinition.h"
18 
24 {
25 public:
26  ComputeCrackTipEnrichmentSmallStrain(const InputParameters & parameters);
28 
29 protected:
30  virtual void computeProperties() override;
31 
32  virtual void computeQpProperties() override;
33 
35  std::vector<Real> _enrich_disp;
36 
38  std::vector<RealVectorValue> _grad_enrich_disp;
39 
41  std::vector<std::vector<MooseVariable *>> _enrich_variable;
42 
44  const VariablePhiValue & _phi;
45 
47  const VariablePhiGradient & _grad_phi;
48 
49 private:
51  std::vector<Real> _B;
53  std::vector<RealVectorValue> _dBX;
55  std::vector<RealVectorValue> _dBx;
57  std::vector<std::vector<Real>> _BI;
59  const std::vector<std::vector<Real>> * _fe_phi;
61  const std::vector<std::vector<RealGradient>> * _fe_dphi;
62  NonlinearSystem * _nl;
63  const NumericVector<Number> * _sln;
64 };
65 
66 #endif // COMPUTECRACKTIPENRICHMENTSMALLSTRAIN_H
const VariablePhiValue & _phi
the current shape functions
std::vector< Real > _B
enrichment function value
std::vector< RealVectorValue > _grad_enrich_disp
gradient of enrichment displacement
std::vector< RealVectorValue > _dBX
derivatives of enrichment function respect to global cooridnate
const VariablePhiGradient & _grad_phi
gradient of the shape function
const std::vector< std::vector< Real > > * _fe_phi
shape function
std::vector< std::vector< Real > > _BI
enrichment function at node I
std::vector< Real > _enrich_disp
enrichment displacement
Perform calculation of enrichment function values and derivatives.
std::vector< RealVectorValue > _dBx
derivatives of enrichment function respect to crack front cooridnate
ComputeCrackTipEnrichmentSmallStrain(const InputParameters &parameters)
ComputeStrainBase is the base class for strain tensors.
ComputeCrackTipEnrichmentSmallStrain calculates the sum of standard strain and enrichement strain...
const std::vector< std::vector< RealGradient > > * _fe_dphi
gradient of shape function
std::vector< std::vector< MooseVariable * > > _enrich_variable
enrichment displacement variables