| Line | Branch | Exec | Source | 
|---|---|---|---|
| 1 | // This file is part of Bembel, the higher order C++ boundary element library. | ||
| 2 | // | ||
| 3 | // Copyright (C) 2022 see <http://www.bembel.eu> | ||
| 4 | // | ||
| 5 | // It was written as part of a cooperation of J. Doelz, H. Harbrecht, S. Kurz, | ||
| 6 | // M. Multerer, S. Schoeps, and F. Wolf at Technische Universitaet Darmstadt, | ||
| 7 | // Universitaet Basel, and Universita della Svizzera italiana, Lugano. This | ||
| 8 | // source code is subject to the GNU General Public License version 3 and | ||
| 9 | // provided WITHOUT ANY WARRANTY, see <http://www.bembel.eu> for further | ||
| 10 | // information. | ||
| 11 | // | ||
| 12 | #ifndef BEMBEL_SRC_SPLINE_SHAPEFUNCTIONS_HPP_ | ||
| 13 | #define BEMBEL_SRC_SPLINE_SHAPEFUNCTIONS_HPP_ | ||
| 14 | |||
| 15 | namespace Bembel { | ||
| 16 | namespace Basis { | ||
| 17 | |||
| 18 | using funptr_doubleOut_doubleptrDoubleIn = double (*)(double*, double); | ||
| 19 | using funptr_voidOut_doubleptrDoubleIn = void (*)(double*, double); | ||
| 20 | |||
| 21 | /** | ||
| 22 | * \ingroup Spline | ||
| 23 | * \brief These routines implement a template recursion that allows to choose a | ||
| 24 | *compile time instantiation of a basis-evaluation routine with a runtime p. To | ||
| 25 | *replace the underlying basis, only these routines should be changed. | ||
| 26 | **/ | ||
| 27 | template <int P> | ||
| 28 | class PSpecificShapeFunctionHandler { | ||
| 29 | public: | ||
| 30 | 15300 | inline static double evalCoef(int p, double* ar, double x) { | |
| 31 | 2/2✓ Branch 0 taken 476 times. ✓ Branch 1 taken 7174 times. | 15300 | return p == P ? Bembel::Basis::EvalBernstein<double, P>(ar, x) | 
| 32 | 15300 | : PSpecificShapeFunctionHandler<P - 1>::evalCoef(p, ar, x); | |
| 33 | } | ||
| 34 | 4620 | inline static double evalDerCoef(int p, double* ar, double x) { | |
| 35 | 2/2✓ Branch 0 taken 220 times. ✓ Branch 1 taken 2090 times. | 4620 | return p == P ? Bembel::Basis::EvalBernsteinDer<double, P>(ar, x) | 
| 36 | 4620 | : PSpecificShapeFunctionHandler<P - 1>::evalDerCoef(p, ar, x); | |
| 37 | } | ||
| 38 | 7949353304 | inline static void evalBasis(int p, double* ar, double x) { | |
| 39 | 2/2✓ Branch 0 taken 243962956 times. ✓ Branch 1 taken 3730713696 times. | 7949353304 | return p == P ? Bembel::Basis::EvalBernsteinBasis<double, P>(ar, x) | 
| 40 | 7949353304 | : PSpecificShapeFunctionHandler<P - 1>::evalBasis(p, ar, x); | |
| 41 | } | ||
| 42 | 7875393386 | inline static void evalDerBasis(int p, double* ar, double x) { | |
| 43 | return p == P | ||
| 44 | 2/2✓ Branch 0 taken 210090490 times. ✓ Branch 1 taken 3727606203 times. | 7875393386 | ? Bembel::Basis::EvalBernsteinDerBasis<double, P>(ar, x) | 
| 45 | 7875393386 | : PSpecificShapeFunctionHandler<P - 1>::evalDerBasis(p, ar, x); | |
| 46 | } | ||
| 47 | inline static constexpr funptr_doubleOut_doubleptrDoubleIn ptrEvalCoef( | ||
| 48 | int p) { | ||
| 49 | return p == P ? &Bembel::Basis::EvalBernstein<double, P> | ||
| 50 | : PSpecificShapeFunctionHandler<P - 1>::ptrEvalCoef(p); | ||
| 51 | } | ||
| 52 | inline static constexpr funptr_doubleOut_doubleptrDoubleIn ptrEvalDerCoef( | ||
| 53 | int p) { | ||
| 54 | return p == P ? &Bembel::Basis::EvalBernsteinDer<double, P> | ||
| 55 | : PSpecificShapeFunctionHandler<P - 1>::ptrEvalDerCoef(p); | ||
| 56 | } | ||
| 57 | inline static constexpr funptr_voidOut_doubleptrDoubleIn ptrEvalBasis(int p) { | ||
| 58 | return p == P ? &Bembel::Basis::EvalBernsteinBasis<double, P> | ||
| 59 | : PSpecificShapeFunctionHandler<P - 1>::ptrEvalBasis(p); | ||
| 60 | } | ||
| 61 | inline static constexpr funptr_voidOut_doubleptrDoubleIn ptrEvalDerBasis( | ||
| 62 | int p) { | ||
| 63 | return p == P ? &Bembel::Basis::EvalBernsteinDerBasis<double, P> | ||
| 64 | : PSpecificShapeFunctionHandler<P - 1>::ptrEvalDerBasis(p); | ||
| 65 | } | ||
| 66 | inline static constexpr bool checkP(int p) { | ||
| 67 | static_assert(P > 0, "Polynomial degree must be larger than zero"); | ||
| 68 | return p <= Constants::MaxP; | ||
| 69 | } | ||
| 70 | }; | ||
| 71 | |||
| 72 | template <> | ||
| 73 | class PSpecificShapeFunctionHandler<0> { | ||
| 74 | public: | ||
| 75 | 11 | inline static double evalCoef(int p, double* ar, double x) { | |
| 76 | 11 | return Bembel::Basis::EvalBernstein<double, 0>(ar, x); | |
| 77 | } | ||
| 78 | ✗ | inline static double evalDerCoef(int p, double* ar, double x) { | |
| 79 | ✗ | return Bembel::Basis::EvalBernsteinDer<double, 0>(ar, x); | |
| 80 | } | ||
| 81 | 11330670 | inline static void evalBasis(int p, double* ar, double x) { | |
| 82 | 11330670 | return Bembel::Basis::EvalBernsteinBasis<double, 0>(ar, x); | |
| 83 | } | ||
| 84 | ✗ | inline static void evalDerBasis(int p, double* ar, double x) { | |
| 85 | ✗ | return Bembel::Basis::EvalBernsteinDerBasis<double, 0>(ar, x); | |
| 86 | } | ||
| 87 | inline static constexpr funptr_doubleOut_doubleptrDoubleIn ptrEvalCoef( | ||
| 88 | int p) { | ||
| 89 | return &Bembel::Basis::EvalBernstein<double, 0>; | ||
| 90 | } | ||
| 91 | inline static constexpr funptr_doubleOut_doubleptrDoubleIn ptrEvalDerCoef( | ||
| 92 | int p) { | ||
| 93 | return &Bembel::Basis::EvalBernsteinDer<double, 0>; | ||
| 94 | } | ||
| 95 | inline static constexpr funptr_voidOut_doubleptrDoubleIn ptrEvalBasis(int p) { | ||
| 96 | return &Bembel::Basis::EvalBernsteinBasis<double, 0>; | ||
| 97 | } | ||
| 98 | inline static constexpr funptr_voidOut_doubleptrDoubleIn ptrEvalDerBasis( | ||
| 99 | int p) { | ||
| 100 | return &Bembel::Basis::EvalBernsteinDerBasis<double, 0>; | ||
| 101 | } | ||
| 102 | inline static constexpr bool checkP(int p) { return Constants::MaxP >= 0; } | ||
| 103 | }; | ||
| 104 | |||
| 105 | using ShapeFunctionHandler = PSpecificShapeFunctionHandler<Constants::MaxP>; | ||
| 106 | |||
| 107 | } // namespace Basis | ||
| 108 | } // namespace Bembel | ||
| 109 | #endif // BEMBEL_SRC_SPLINE_SHAPEFUNCTIONS_HPP_ | ||
| 110 |