HyperHDG
tensorial.hxx
Go to the documentation of this file.
1 #pragma once // Ensure that file is included only once in a single compilation.
2 
3 #include <tpp/hypercube.hxx>
5 #include <tpp/tpp_assert.hxx>
6 
7 #include <array>
8 
9 namespace TPP
10 {
11 namespace ShapeType
12 {
13 /*!*************************************************************************************************
14  * \brief Struct that handles the evaluation of tensorial shape functions.
15  *
16  * \tparam shape_t The typename of the one dimensional shape functions.
17  * \tparam dimT Local dimension of the shape function's domain or the point.
18  *
19  * \authors Andreas Rupp, Heidelberg University, 2021.
20  **************************************************************************************************/
21 template <typename shape_t, unsigned int dimT>
22 struct Tensorial
23 {
24  /*!***********************************************************************************************
25  * \brief Make type of one-dimensional shape functions accessable to everyone.
26  ************************************************************************************************/
27  typedef shape_t shape_fun_1d;
28  /*!***********************************************************************************************
29  * \brief Number of shape functions that span the local polynomial space.
30  ************************************************************************************************/
31  static constexpr unsigned int n_fun() { return Hypercube<dimT>::pow(shape_fun_1d::n_fun()); }
32  /*!***********************************************************************************************
33  * \brief Dimension of abscissas of the polynomials.
34  ************************************************************************************************/
35  static constexpr unsigned int dim() { return dimT; }
36  /*!***********************************************************************************************
37  * \brief Maximum degree of all polynomials.
38  ************************************************************************************************/
39  static constexpr unsigned int degree() { return shape_t::degree(); }
40  /*!***********************************************************************************************
41  * \brief Evaluate value of one shape function.
42  *
43  * Evaluates value of the \c index one-dimensional shape function on the reference interval
44  * \f$[0,1]\f$ at abscissas \c point.
45  *
46  * \tparam return_t Floating type specification for return value.
47  * \tparam point_t Type of abscissa value.
48  * \param index Index of evaluated shape function.
49  * \param point Abscissa / function argument.
50  * \retval fct_value Evaluated value of shape function.
51  ************************************************************************************************/
52  template <typename return_t, typename point_t>
53  static constexpr return_t fct_val(__attribute__((unused)) const unsigned int index,
54  __attribute__((unused)) const point_t& point)
55  {
56  static_assert(point_t::size() == dim(), "Point needs to have correct dimension.");
57 
58  return_t value = 1.;
59  if constexpr (dimT == 0)
60  return (return_t)1.;
61  else
62  {
63  std::array<unsigned int, std::max(dimT, 1U)> index_dim =
64  Hypercube<dimT>::index_decompose(index, shape_fun_1d::n_fun());
65  for (unsigned int dim = 0; dim < dimT; ++dim)
66  value *= shape_fun_1d::template fct_val<return_t>(index_dim[dim], point[dim]);
67  }
68  return value;
69  }
70  /*!***********************************************************************************************
71  * \brief Evaluate value of one shape function.
72  *
73  * Evaluates value of the \c index one-dimensional shape function on the reference interval
74  * \f$[0,1]\f$ at abscissas \c x_val.
75  *
76  * \tparam return_t Floating type specification for return value.
77  * \tparam point_t Type of abscissa value.
78  * \param index Index of evaluated shape function.
79  * \param point Abscissa / function argument.
80  * \param der_dim Dimension with respect to which derivative is evaluated.
81  * \retval fct_value Evaluated value of shape function's derivative.
82  ************************************************************************************************/
83  template <typename return_t, typename point_t>
84  static constexpr return_t der_val(const unsigned int index,
85  const point_t& point,
86  const unsigned int der_dim)
87  {
88  static_assert(point_t::size() == dim(), "Point needs to have correct dimension.");
89  tpp_assert(der_dim < dimT, "The derivative needs to be with respect to a valid dimension.");
90 
91  return_t value = 1.;
92  std::array<unsigned int, std::max(dimT, 1U)> index_dim =
93  Hypercube<dimT>::index_decompose(index, shape_fun_1d::n_fun());
94  for (unsigned int dim = 0; dim < dimT; ++dim)
95  if (der_dim == dim)
96  value *= shape_fun_1d::template der_val<return_t>(index_dim[dim], point[dim]);
97  else
98  value *= shape_fun_1d::template fct_val<return_t>(index_dim[dim], point[dim]);
99 
100  return value;
101  }
102 }; // end of struct Tensorial
103 
104 } // end of namespace ShapeType
105 
106 } // end of namespace TPP
check_push_test.index
index
Definition: check_push_test.py:10
one_dimensional.hxx
tpp_assert
#define tpp_assert(Expr, Msg)
The assertion to be used within tpp — deactivate using -DNDEBUG compile flag.
Definition: tpp_assert.hxx:33
hypercube.hxx
TPP::ShapeType::Tensorial
Struct that handles the evaluation of tensorial shape functions.
Definition: tensorial.hxx:22
TPP::ShapeType::Tensorial::der_val
static constexpr return_t der_val(const unsigned int index, const point_t &point, const unsigned int der_dim)
Evaluate value of one shape function.
Definition: tensorial.hxx:84
tpp_assert.hxx
This file provides the function tpp_assert.
TPP
Contains the functionalities to evaluate and integrate tensor product type polynomials.
Definition: compile_time_tricks.hxx:7
TPP::ShapeType::Tensorial::shape_fun_1d
shape_t shape_fun_1d
Make type of one-dimensional shape functions accessable to everyone.
Definition: tensorial.hxx:27
TPP::ShapeType::Tensorial::n_fun
static constexpr unsigned int n_fun()
Number of shape functions that span the local polynomial space.
Definition: tensorial.hxx:31
TPP::ShapeType::Tensorial::dim
static constexpr unsigned int dim()
Dimension of abscissas of the polynomials.
Definition: tensorial.hxx:35
TPP::Hypercube::pow
static constexpr unsigned int pow(unsigned int base)
Return n to the power dim.
Definition: hypercube.hxx:29
TPP::ShapeType::Tensorial::fct_val
static constexpr return_t fct_val(__attribute__((unused)) const unsigned int index, __attribute__((unused)) const point_t &point)
Evaluate value of one shape function.
Definition: tensorial.hxx:53
TPP::Hypercube::index_decompose
static constexpr std::array< unsigned int, std::max(dimT, 1U)> index_decompose(unsigned int index, const unsigned int range)
Decompose global index of tensorial structure with respect to its dimensions.
Definition: hypercube.hxx:44
TPP::ShapeType::Tensorial::degree
static constexpr unsigned int degree()
Maximum degree of all polynomials.
Definition: tensorial.hxx:39