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 35 of file TemplateIsoparametric.h.

#include <TemplateIsoparametric.h>

Public Types

using ShapeFunctionType = ShapeFunctionType_
 
using ShapeMatrices = typename ShapeMatrixTypes_::ShapeMatrices
 Coordinate mapping matrices type. More...
 
using MeshElementType = typename ShapeFunctionType_::MeshElement
 Type of the underlying geometrical element. More...
 
using NaturalCoordsMappingType = NaturalCoordinatesMapping< MeshElementType, ShapeFunctionType, ShapeMatrices >
 

Public Member Functions

 TemplateIsoparametric ()
 
 TemplateIsoparametric (const MeshElementType &e)
 
const MeshElementTypegetMeshElement () const
 return current mesh element More...
 
void setMeshElement (const MeshElementType &e)
 Sets the mesh element. More...
 
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. More...
 

Private Member Functions

void computeIntegralMeasure (bool is_axially_symmetric, ShapeMatrices &shape) const
 

Private Attributes

const MeshElementType_ele
 

Member Typedef Documentation

◆ MeshElementType

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

Type of the underlying geometrical element.

Definition at line 44 of file TemplateIsoparametric.h.

◆ NaturalCoordsMappingType

template<class ShapeFunctionType_ , class ShapeMatrixTypes_ >
using NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::NaturalCoordsMappingType = NaturalCoordinatesMapping<MeshElementType, ShapeFunctionType, ShapeMatrices>

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

Definition at line 48 of file TemplateIsoparametric.h.

◆ ShapeFunctionType

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

Definition at line 38 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 41 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 56 of file TemplateIsoparametric.h.

56 : _ele(nullptr) {}

◆ TemplateIsoparametric() [2/2]

template<class ShapeFunctionType_ , class ShapeMatrixTypes_ >
NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::TemplateIsoparametric ( const MeshElementType e)
inlineexplicit

Construct this object for the given mesh element.

Parameters
eMesh element object

Definition at line 63 of file TemplateIsoparametric.h.

63 : _ele(&e) {}

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 142 of file TemplateIsoparametric.h.

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  }
double interpolateZerothCoordinate(typename ShapeMatrices::ShapeType const &N) const
static const double r

References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::interpolateZerothCoordinate(), and MathLib::r.

Referenced by NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::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 79 of file TemplateIsoparametric.h.

82  {
84  global_dim);
85  computeIntegralMeasure(is_axially_symmetric, shape);
86  }
void computeIntegralMeasure(bool is_axially_symmetric, ShapeMatrices &shape) const
static void computeShapeMatrices(const T_MESH_ELEMENT &ele, const double *natural_pt, T_SHAPE_MATRICES &shapemat, const unsigned global_dim)

References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::_ele, NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::computeIntegralMeasure(), and NumLib::NaturalCoordinatesMapping< T_MESH_ELEMENT, T_SHAPE_FUNC, T_SHAPE_MATRICES >::computeShapeMatrices().

◆ 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 99 of file TemplateIsoparametric.h.

102  {
103  NaturalCoordsMappingType::template computeShapeMatrices<
104  T_SHAPE_MATRIX_TYPE>(*_ele, natural_pt, shape, global_dim);
105  computeIntegralMeasure(is_axially_symmetric, shape);
106  }
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)

References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::_ele, NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::computeIntegralMeasure(), and NumLib::computeShapeMatrices().

◆ getMeshElement()

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

return current mesh element

Definition at line 66 of file TemplateIsoparametric.h.

66 { return _ele; }

References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::_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 124 of file TemplateIsoparametric.h.

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  }

References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::_ele, and NumLib::N.

◆ 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 110 of file TemplateIsoparametric.h.

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  }

References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::_ele, NumLib::N, and MathLib::r.

Referenced by NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::computeIntegralMeasure().

◆ setMeshElement()

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

Sets the mesh element.

Definition at line 69 of file TemplateIsoparametric.h.

69 { this->_ele = &e; }

References NumLib::TemplateIsoparametric< ShapeFunctionType_, ShapeMatrixTypes_ >::_ele.

Member Data Documentation

◆ _ele


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