31template <
class ShapeFunctionType_,
class ShapeMatrixTypes_>
73 const unsigned global_dim,
74 bool is_axially_symmetric)
const
91 template <ShapeMatrixType T_SHAPE_MATRIX_TYPE>
93 const unsigned global_dim,
94 bool is_axially_symmetric)
const
97 T_SHAPE_MATRIX_TYPE>(*
_ele, natural_pt, shape, global_dim);
108 for (
int i = 0; i < rs.size(); ++i)
110 rs[i] = (*nodes[i])[0];
112 auto const r =
N.dot(rs);
122 std::array<double, 3> interpolated_coords;
123 for (
int d = 0; d < 3; ++d)
125 for (
int i = 0; i < rs.size(); ++i)
127 rs[i] = (*nodes[i])[d];
129 interpolated_coords[d] =
N.dot(rs);
131 return interpolated_coords;
138 if (!is_axially_symmetric)
140 shape.integralMeasure = 1.0;
151 shape.integralMeasure = 2. * std::numbers::pi * r;
159template <
typename ShapeFunction,
typename ShapeMatricesType>
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
ShapeFunctionType_ ShapeFunctionType
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)
const MeshLib::Element * _ele
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)