OGS
TemplateIsoparametric.h
Go to the documentation of this file.
1 
13 #pragma once
14 
15 #include <boost/math/constants/constants.hpp>
16 #include <cassert>
17 
20 
21 namespace MeshLib
22 {
23 class Element;
24 }
25 
26 namespace NumLib
27 {
34 template <class ShapeFunctionType_, class ShapeMatrixTypes_>
36 {
37 public:
38  using ShapeFunctionType = ShapeFunctionType_;
39 
41  using ShapeMatrices = typename ShapeMatrixTypes_::ShapeMatrices;
42 
44  using MeshElementType = typename ShapeFunctionType_::MeshElement;
45 
51 
56  TemplateIsoparametric() : _ele(nullptr) {}
57 
63  explicit TemplateIsoparametric(const MeshElementType& e) : _ele(&e) {}
64 
66  const MeshElementType* getMeshElement() const { return _ele; }
67 
69  void setMeshElement(const MeshElementType& e) { this->_ele = &e; }
79  void computeShapeFunctions(const double* natural_pt, ShapeMatrices& shape,
80  const unsigned global_dim,
81  bool is_axially_symmetric) const
82  {
84  global_dim);
85  computeIntegralMeasure(is_axially_symmetric, shape);
86  }
87 
98  template <ShapeMatrixType T_SHAPE_MATRIX_TYPE>
99  void computeShapeFunctions(const double* natural_pt, ShapeMatrices& shape,
100  const unsigned global_dim,
101  bool is_axially_symmetric) const
102  {
103  NaturalCoordsMappingType::template computeShapeMatrices<
104  T_SHAPE_MATRIX_TYPE>(*_ele, natural_pt, shape, global_dim);
105  computeIntegralMeasure(is_axially_symmetric, shape);
106  }
107 
111  typename ShapeMatrices::ShapeType const& N) const
112  {
113  auto* nodes = _ele->getNodes();
114  typename ShapeMatrices::ShapeType rs(N.size());
115  for (int i = 0; i < rs.size(); ++i)
116  {
117  rs[i] = (*nodes[i])[0];
118  }
119  auto const r = N.dot(rs);
120  return r;
121  }
122 
124  std::array<double, 3> interpolateCoordinates(
125  typename ShapeMatrices::ShapeType const& N) const
126  {
127  auto* nodes = _ele->getNodes();
128  typename ShapeMatrices::ShapeType rs(N.size());
129  std::array<double, 3> interpolated_coords;
130  for (int d = 0; d < 3; ++d)
131  {
132  for (int i = 0; i < rs.size(); ++i)
133  {
134  rs[i] = (*nodes[i])[d];
135  }
136  interpolated_coords[d] = N.dot(rs);
137  }
138  return interpolated_coords;
139  }
140 
141 private:
142  void computeIntegralMeasure(bool is_axially_symmetric,
143  ShapeMatrices& shape) const
144  {
145  if (!is_axially_symmetric)
146  {
147  shape.integralMeasure = 1.0;
148  return;
149  }
150 
151  // Note: If an integration point is located at the rotation axis, r will
152  // be zero which might lead to problems with the assembled equation
153  // system.
154  // E.g., for triangle elements, if an integration point is
155  // located at edge of the unit triangle, it is possible that
156  // r becomes zero.
157  auto const r = interpolateZerothCoordinate(shape.N);
158  shape.integralMeasure = boost::math::constants::two_pi<double>() * r;
159  }
160 
162 };
163 
166 template <typename ShapeFunction, typename ShapeMatricesType>
169 {
170  using FemType =
172 
173  return FemType{
174  *static_cast<const typename ShapeFunction::MeshElement*>(&e)};
175 }
176 } // namespace NumLib
Template class for isoparametric elements.
void computeShapeFunctions(const double *natural_pt, ShapeMatrices &shape, const unsigned global_dim, bool is_axially_symmetric) const
const MeshElementType * getMeshElement() const
return current mesh element
double interpolateZerothCoordinate(typename ShapeMatrices::ShapeType const &N) const
std::array< double, 3 > interpolateCoordinates(typename ShapeMatrices::ShapeType const &N) const
Interpolates the coordinates of the element with the given shape matrix.
void setMeshElement(const MeshElementType &e)
Sets the mesh element.
void computeIntegralMeasure(bool is_axially_symmetric, ShapeMatrices &shape) const
TemplateIsoparametric(const MeshElementType &e)
typename ShapeFunctionType_::MeshElement MeshElementType
Type of the underlying geometrical element.
void computeShapeFunctions(const double *natural_pt, ShapeMatrices &shape, const unsigned global_dim, bool is_axially_symmetric) const
typename ShapeMatrixTypes_::ShapeMatrices ShapeMatrices
Coordinate mapping matrices type.
static const double r
NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatricesType > createIsoparametricFiniteElement(MeshLib::Element const &e)
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
static void computeShapeMatrices(const T_MESH_ELEMENT &ele, const double *natural_pt, T_SHAPE_MATRICES &shapemat, const unsigned global_dim)