Bembel
Potential.hpp
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 #ifndef BEMBEL_SRC_POTENTIAL_POTENTIAL_HPP_
12 #define BEMBEL_SRC_POTENTIAL_POTENTIAL_HPP_
13 
14 #include <Eigen/Dense>
15 
16 #include "../LinearOperator/DifferentialFormEnum.hpp"
17 #include "../LinearOperator/LinearOperatorTraits.hpp"
18 #include "../util/Macros.hpp"
19 
20 namespace Bembel {
27 template <typename Derived>
29  enum { YOU_DID_NOT_SPECIFY_POTENTIAL_TRAITS = 1 };
30 };
31 
35 template <typename S, typename T>
37  enum { RETURN_TYPE_ONLY_SPECIFIED_FOR_DOUBLE_OR_COMPLEX_DOUBLE = 1 };
38 };
43 template <>
44 struct PotentialReturnScalar<double, double> {
45  typedef double Scalar;
46 };
51 template <>
52 struct PotentialReturnScalar<std::complex<double>, double> {
53  typedef std::complex<double> Scalar;
54 };
59 template <>
60 struct PotentialReturnScalar<double, std::complex<double>> {
61  typedef std::complex<double> Scalar;
62 };
67 template <>
68 struct PotentialReturnScalar<std::complex<double>, std::complex<double>> {
69  typedef std::complex<double> Scalar;
70 };
71 
80 template <typename Derived, typename LinOp>
81 struct PotentialBase {
82  // Constructors
83  PotentialBase() {}
84 
85  // the user has to provide the implementation of this function, which
86  // tells
87  // is able to evaluate the integrand of the Galerkin formulation in a
88  // pair
89  // of quadrature points represented as a
90  // Surface point [xi; w; Chi(xi); dsChi(xi); dtChi(xi)]
91  template <typename T>
92  Eigen::Matrix<typename PotentialReturnScalar<
94  typename PotentialTraits<Derived>::Scalar>::Scalar,
96  evaluateIntegrand(
97  const AnsatzSpace<LinOp> &ansatz_space,
98  const Eigen::Matrix<
99  typename LinearOperatorTraits<LinOp>::Scalar, Eigen::Dynamic,
101  &coeff,
102  const Eigen::VectorXd &point, const SurfacePoint &p) const {
103  return static_cast<const Derived *>(this)->evaluatePotential_impl(
104  ansatz_space, coeff, point, p);
105  }
106  // pointer to the derived object
107  Derived &derived() { return *static_cast<Derived *>(this); }
108  // const pointer to the derived object
109  const Derived &derived() const { return *static_cast<const Derived *>(this); }
110 };
111 } // namespace Bembel
112 #endif // BEMBEL_SRC_POTENTIAL_POTENTIAL_HPP_
Eigen::Matrix< double, 12, 1 > SurfacePoint
typedef of SurfacePoint
Routines for the evalutation of pointwise errors.
Definition: AnsatzSpace.hpp:14
constexpr int getFunctionSpaceVectorDimension()
struct containing specifications on the linear operator has to be specialized or derived for any part...
functional base class. this serves as a common interface for existing functionals.
Definition: Potential.hpp:81
Base case for specifying the return type of the potential.
Definition: Potential.hpp:36
struct containing specifications on the functional has to be specialized or derived for any particula...
Definition: Potential.hpp:28