OGS
NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ > Class Template Reference

Detailed Description

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
class NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >

Template class for isoparametric elements.

Template Parameters
ShapeFunctionType_The shape function type.
ShapeMatrixTypes_An aggregate of shape matrix types.

Definition at line 32 of file TemplateIsoparametric.h.

#include <TemplateIsoparametric.h>

Collaboration diagram for NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >:
[legend]

Public Types

using ShapeFunctionType = ShapeFunctionType_
using ShapeMatrices = typename ShapeMatrixTypes_::ShapeMatrices
 Coordinate mapping matrices type.
using NaturalCoordsMappingType

Public Member Functions

 TemplateIsoparametric ()
 TemplateIsoparametric (const MeshLib::Element &e)
const MeshLib::ElementgetMeshElement () const
 return current mesh element
void setMeshElement (const MeshLib::Element &e)
 Sets the mesh element.
void computeShapeFunctions (const double *natural_pt, ShapeMatrices &shape, const unsigned global_dim, bool is_axially_symmetric) const
template<ShapeMatrixType T_SHAPE_MATRIX_TYPE>
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
std::array< double, 3 > interpolateCoordinates (typename ShapeMatrices::ShapeType const &N) const
 Interpolates the coordinates of the element with the given shape matrix.

Private Member Functions

void computeIntegralMeasure (bool is_axially_symmetric, ShapeMatrices &shape) const

Private Attributes

const MeshLib::Element_ele

Member Typedef Documentation

◆ NaturalCoordsMappingType

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
using NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::NaturalCoordsMappingType
Initial value:

Natural coordinates mapping tools specialization for specific MeshElement and ShapeFunction types.

Definition at line 42 of file TemplateIsoparametric.h.

◆ ShapeFunctionType

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
using NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::ShapeFunctionType = ShapeFunctionType_

Definition at line 35 of file TemplateIsoparametric.h.

◆ ShapeMatrices

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
using NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::ShapeMatrices = typename ShapeMatrixTypes_::ShapeMatrices

Coordinate mapping matrices type.

Definition at line 38 of file TemplateIsoparametric.h.

Constructor & Destructor Documentation

◆ TemplateIsoparametric() [1/2]

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::TemplateIsoparametric ( )
inline

Constructor without specifying a mesh element. setMeshElement() must be called afterwards.

Definition at line 49 of file TemplateIsoparametric.h.

49: _ele(nullptr) {}

◆ TemplateIsoparametric() [2/2]

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::TemplateIsoparametric ( const MeshLib::Element & e)
inlineexplicit

Construct this object for the given mesh element.

Parameters
eMesh element object

Definition at line 56 of file TemplateIsoparametric.h.

56: _ele(&e) {}
Template class for isoparametric elements.

Member Function Documentation

◆ computeIntegralMeasure()

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
void NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::computeIntegralMeasure ( bool is_axially_symmetric,
ShapeMatrices & shape ) const
inlineprivate

Definition at line 135 of file TemplateIsoparametric.h.

137 {
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 = 2. * std::numbers::pi * r;
152 }
double interpolateZerothCoordinate(typename ShapeMatrices::ShapeType const &N) const

Referenced by NumLib::TemplateIsoparametric< ShapePoint1, T_SHAPE_MATRIX_POLICY< ShapePoint1 > >::computeShapeFunctions().

◆ computeShapeFunctions() [1/2]

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
void NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::computeShapeFunctions ( const double * natural_pt,
ShapeMatrices & shape,
const unsigned global_dim,
bool is_axially_symmetric ) const
inline

compute shape functions

Parameters
natural_ptposition in natural coordinates
shapeevaluated shape function matrices
global_dimglobal dimension
is_axially_symmetricif true, the integral measure for cylinder coordinates is used, and 1 otherwise.

Definition at line 72 of file TemplateIsoparametric.h.

75 {
79 }
void computeIntegralMeasure(bool is_axially_symmetric, ShapeMatrices &shape) const
static void computeShapeMatrices(const MeshLib::Element &ele, const double *natural_pt, ShapeMatrices &shapemat, const unsigned global_dim)

◆ computeShapeFunctions() [2/2]

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
template<ShapeMatrixType T_SHAPE_MATRIX_TYPE>
void NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::computeShapeFunctions ( const double * natural_pt,
ShapeMatrices & shape,
const unsigned global_dim,
bool is_axially_symmetric ) const
inline

compute shape functions

Template Parameters
T_SHAPE_MATRIX_TYPEshape matrix types to be calculated
Parameters
natural_ptposition in natural coordinates
shapeevaluated shape function matrices
global_dimglobal dimension
is_axially_symmetricif true, the integral measure for cylinder coordinates is used, and 1 otherwise.

Definition at line 92 of file TemplateIsoparametric.h.

◆ getMeshElement()

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
const MeshLib::Element * NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::getMeshElement ( ) const
inline

return current mesh element

Definition at line 59 of file TemplateIsoparametric.h.

59{ return _ele; }

◆ interpolateCoordinates()

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
std::array< double, 3 > NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::interpolateCoordinates ( typename ShapeMatrices::ShapeType const & N) const
inline

Interpolates the coordinates of the element with the given shape matrix.

Definition at line 117 of file TemplateIsoparametric.h.

119 {
120 auto* nodes = _ele->getNodes();
121 typename ShapeMatrices::ShapeType rs(N.size());
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 }

Referenced by ProcessLib::BoundaryConditionAndSourceTerm::Python::BcAndStLocalAssemblerImpl< BcOrStData, ShapeFunction, LowerOrderShapeFunction, GlobalDim >::interpolateCoords().

◆ interpolateZerothCoordinate()

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
double NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::interpolateZerothCoordinate ( typename ShapeMatrices::ShapeType const & N) const
inline

Interpolates the x coordinate of the element with the given shape matrix.

Definition at line 103 of file TemplateIsoparametric.h.

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 }

Referenced by NumLib::TemplateIsoparametric< ShapePoint1, T_SHAPE_MATRIX_POLICY< ShapePoint1 > >::computeIntegralMeasure().

◆ setMeshElement()

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
void NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::setMeshElement ( const MeshLib::Element & e)
inline

Sets the mesh element.

Definition at line 62 of file TemplateIsoparametric.h.

62{ this->_ele = &e; }

Member Data Documentation

◆ _ele

template<class ShapeFunctionType_, class ShapeMatrixTypes_>
const MeshLib::Element* NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::_ele
private

Definition at line 154 of file TemplateIsoparametric.h.


The documentation for this class was generated from the following file: