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
19#include "MeshLib/Node.h"
22
23namespace NumLib
24{
31template <class ShapeFunctionType_, class ShapeMatrixTypes_>
33{
34public:
35 using ShapeFunctionType = ShapeFunctionType_;
36
38 using ShapeMatrices = typename ShapeMatrixTypes_::ShapeMatrices;
39
44
50
56 explicit TemplateIsoparametric(const MeshLib::Element& e) : _ele(&e) {}
57
59 const MeshLib::Element* getMeshElement() const { return _ele; }
60
62 void setMeshElement(const MeshLib::Element& e) { this->_ele = &e; }
72 void computeShapeFunctions(const double* natural_pt, ShapeMatrices& shape,
73 const unsigned global_dim,
74 bool is_axially_symmetric) const
75 {
77 global_dim);
78 computeIntegralMeasure(is_axially_symmetric, shape);
79 }
80
91 template <ShapeMatrixType T_SHAPE_MATRIX_TYPE>
92 void computeShapeFunctions(const double* natural_pt, ShapeMatrices& shape,
93 const unsigned global_dim,
94 bool is_axially_symmetric) const
95 {
96 NaturalCoordsMappingType::template computeShapeMatrices<
97 T_SHAPE_MATRIX_TYPE>(*_ele, natural_pt, shape, global_dim);
98 computeIntegralMeasure(is_axially_symmetric, shape);
99 }
100
104 typename ShapeMatrices::ShapeType const& N) const
105 {
106 auto* nodes = _ele->getNodes();
107 typename ShapeMatrices::ShapeType rs(N.size());
108 for (int i = 0; i < rs.size(); ++i)
109 {
110 rs[i] = (*nodes[i])[0];
111 }
112 auto const r = N.dot(rs);
113 return r;
114 }
115
117 std::array<double, 3> interpolateCoordinates(
118 typename ShapeMatrices::ShapeType const& N) const
119 {
120 auto* nodes = _ele->getNodes();
121 typename ShapeMatrices::ShapeType rs(N.size());
122 std::array<double, 3> interpolated_coords;
123 for (int d = 0; d < 3; ++d)
124 {
125 for (int i = 0; i < rs.size(); ++i)
126 {
127 rs[i] = (*nodes[i])[d];
128 }
129 interpolated_coords[d] = N.dot(rs);
130 }
131 return interpolated_coords;
132 }
133
134private:
135 void computeIntegralMeasure(bool is_axially_symmetric,
136 ShapeMatrices& shape) const
137 {
138 if (!is_axially_symmetric)
139 {
140 shape.integralMeasure = 1.0;
141 return;
142 }
143
144 // Note: If an integration point is located at the rotation axis, r will
145 // be zero which might lead to problems with the assembled equation
146 // system.
147 // E.g., for triangle elements, if an integration point is
148 // located at edge of the unit triangle, it is possible that
149 // r becomes zero.
150 auto const r = interpolateZerothCoordinate(shape.N);
151 shape.integralMeasure = boost::math::constants::two_pi<double>() * r;
152 }
153
155};
156
159template <typename ShapeFunction, typename ShapeMatricesType>
168} // namespace NumLib
Definition of the Element class.
Definition of the Node class.
virtual Node *const * getNodes() const =0
Get array of element nodes.
Template class for isoparametric elements.
void computeShapeFunctions(const double *natural_pt, ShapeMatrices &shape, const unsigned global_dim, bool is_axially_symmetric) const
double interpolateZerothCoordinate(typename ShapeMatrices::ShapeType const &N) const
void setMeshElement(const MeshLib::Element &e)
Sets the mesh element.
void computeIntegralMeasure(bool is_axially_symmetric, ShapeMatrices &shape) const
std::array< double, 3 > interpolateCoordinates(typename ShapeMatrices::ShapeType const &N) const
Interpolates the coordinates of the element with the given shape matrix.
const MeshLib::Element * getMeshElement() const
return current mesh element
void computeShapeFunctions(const double *natural_pt, ShapeMatrices &shape, const unsigned global_dim, bool is_axially_symmetric) const
TemplateIsoparametric(const MeshLib::Element &e)
typename ShapeMatrixTypes_::ShapeMatrices ShapeMatrices
Coordinate mapping matrices type.
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)
NumLib::TemplateIsoparametric< ShapeFunction, ShapeMatricesType > createIsoparametricFiniteElement(MeshLib::Element const &e)
static void computeShapeMatrices(const MeshLib::Element &ele, const double *natural_pt, T_SHAPE_MATRICES &shapemat, const unsigned global_dim)