HyperHDG
shape_function.hxx
Go to the documentation of this file.
1 #pragma once // Ensure that file is included only once in a single compilation.
2 
4 
5 namespace TPP
6 {
7 /*!*************************************************************************************************
8  * \brief Struct that handles different types of evaluation of shape functions.
9  *
10  * \tparam shape_t Type of shape functions.
11  *
12  * \authors Andreas Rupp, Heidelberg University, 2021.
13  **************************************************************************************************/
14 template <typename shape_t>
16 {
17  /*!***********************************************************************************************
18  * \brief Make type of shape functions accessable to everyone.
19  ************************************************************************************************/
20  typedef shape_t shape_fun_t;
21  /*!***********************************************************************************************
22  * \brief Number of shape functions that span the local polynomial space.
23  ************************************************************************************************/
24  static constexpr unsigned int n_fun() { return shape_t::n_fun(); }
25  /*!***********************************************************************************************
26  * \brief Dimension of abscissas of the polynomials.
27  ************************************************************************************************/
28  static constexpr unsigned int dim() { return shape_t::dim(); }
29  /*!***********************************************************************************************
30  * \brief Maximum degree of all polynomials.
31  ************************************************************************************************/
32  static constexpr unsigned int degree() { return shape_t::degree(); }
33  /*!***********************************************************************************************
34  * \brief Evaluate value of shape function.
35  ************************************************************************************************/
36  template <typename return_t, typename point_t>
37  static constexpr return_t fct_val(const unsigned int index, const point_t& point)
38  {
39  return shape_t::template fct_val<return_t>(index, point);
40  }
41  /*!***********************************************************************************************
42  * \brief Evaluate value of partial derivative of shape function.
43  ************************************************************************************************/
44  template <typename return_t, typename point_t>
45  static constexpr return_t der_val(const unsigned int index,
46  const point_t& point,
47  const unsigned int der_dim)
48  {
49  return shape_t::template der_val<return_t>(index, point, der_dim);
50  }
51 
52  /*!***********************************************************************************************
53  * \brief Evaluate value of gradient of shape function.
54  ************************************************************************************************/
55  template <typename return_t, typename point_t>
56  static constexpr return_t grad_val(const unsigned int index, const point_t& point)
57  {
58  static_assert(return_t::size() == point_t::size(), "Gradient and point have same dimension.");
59  static_assert(point_t::size() == dim(), "Gradient and dimension need to be equal.");
60 
61  return_t value;
62  for (unsigned int der_dim = 0; der_dim < dim(); ++der_dim)
63  value[der_dim] = der_val<return_t::value_type>(index, point, der_dim);
64  return value;
65  }
66  /*!***********************************************************************************************
67  * \brief Evaluate value of divergence of shape function
68  ************************************************************************************************/
69  template <typename return_t, typename point_t>
70  static constexpr return_t div_val(const unsigned int index, const point_t& point)
71  {
72  static_assert(point_t::size() == dim(), "Gradient and dimension need to be equal.");
73 
74  return_t value = 0.;
75  for (unsigned int der_dim = 0; der_dim < dim(); ++der_dim)
76  value += der_val<return_t>(index, point, der_dim);
77  return value;
78  }
79 
80  /*!***********************************************************************************************
81  * \brief Evaluate value linear combination of shape functions.
82  ************************************************************************************************/
83  template <typename return_t, typename coeffs_t, typename point_t>
84  static constexpr return_t lin_comb_fct_val(const coeffs_t& coeffs, const point_t& point)
85  {
86  static_assert(point_t::size() == dim(), "Point needs to have correct dimension.");
87  static_assert(coeffs_t::size() == n_fun(), "Coeffs size must coincide with poly space!");
88 
89  return_t value = 0.;
90  for (unsigned int index = 0; index < coeffs.size(); ++index)
91  value += coeffs[index] * fct_val<return_t>(index, point);
92  return value;
93  }
94  /*!***********************************************************************************************
95  * \brief Evaluate value of partial derivative of linear combination of shape functions.
96  ************************************************************************************************/
97  template <typename return_t, typename coeffs_t, typename point_t>
98  static constexpr return_t lin_comb_der_val(const coeffs_t& coeffs,
99  const point_t& point,
100  const unsigned int der_dim)
101  {
102  static_assert(point_t::size() == dim(), "Point needs to have correct dimension.");
103  static_assert(coeffs_t::size() == n_fun(), "Coeffs size must coincide with poly space!");
104 
105  return_t value = 0.;
106  for (unsigned int index = 0; index < coeffs.size(); ++index)
107  value += coeffs[index] * der_val<return_t>(index, point, der_dim);
108  return value;
109  }
110 }; // end of struct ShapeFunction
111 
112 } // end of namespace TPP
check_push_test.index
index
Definition: check_push_test.py:10
TPP::ShapeFunction::n_fun
static constexpr unsigned int n_fun()
Number of shape functions that span the local polynomial space.
Definition: shape_function.hxx:24
tensorial.hxx
TPP::ShapeFunction::dim
static constexpr unsigned int dim()
Dimension of abscissas of the polynomials.
Definition: shape_function.hxx:28
TPP::ShapeFunction::lin_comb_fct_val
static constexpr return_t lin_comb_fct_val(const coeffs_t &coeffs, const point_t &point)
Evaluate value linear combination of shape functions.
Definition: shape_function.hxx:84
TPP::ShapeFunction::lin_comb_der_val
static constexpr return_t lin_comb_der_val(const coeffs_t &coeffs, const point_t &point, const unsigned int der_dim)
Evaluate value of partial derivative of linear combination of shape functions.
Definition: shape_function.hxx:98
TPP::ShapeFunction::degree
static constexpr unsigned int degree()
Maximum degree of all polynomials.
Definition: shape_function.hxx:32
TPP
Contains the functionalities to evaluate and integrate tensor product type polynomials.
Definition: compile_time_tricks.hxx:7
TPP::ShapeFunction
Struct that handles different types of evaluation of shape functions.
Definition: shape_function.hxx:15
TPP::ShapeFunction::der_val
static constexpr return_t der_val(const unsigned int index, const point_t &point, const unsigned int der_dim)
Evaluate value of partial derivative of shape function.
Definition: shape_function.hxx:45
TPP::ShapeFunction::fct_val
static constexpr return_t fct_val(const unsigned int index, const point_t &point)
Evaluate value of shape function.
Definition: shape_function.hxx:37
TPP::ShapeFunction::shape_fun_t
shape_t shape_fun_t
Make type of shape functions accessable to everyone.
Definition: shape_function.hxx:20
TPP::ShapeFunction::div_val
static constexpr return_t div_val(const unsigned int index, const point_t &point)
Evaluate value of divergence of shape function.
Definition: shape_function.hxx:70
TPP::ShapeFunction::grad_val
static constexpr return_t grad_val(const unsigned int index, const point_t &point)
Evaluate value of gradient of shape function.
Definition: shape_function.hxx:56