Bembel
VectorPotential.hpp
1 // This file is part of Bembel, the higher order C++ boundary element library.
2 //
3 // Copyright (C) 2024 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_AUGMENTEDEFIE_VECTORPOTENTIAL_HPP_
13 #define BEMBEL_SRC_AUGMENTEDEFIE_VECTORPOTENTIAL_HPP_
14 
15 namespace Bembel {
16 // forward declaration of class VectorPotential in order to define
17 // traits
18 template <typename LinOp>
19 class VectorPotential;
20 
21 template <typename LinOp>
23  typedef Eigen::VectorXcd::Scalar Scalar;
24  static constexpr int OutputSpaceDimension = 3;
25 };
26 
30 template <typename LinOp>
31 class VectorPotential : public PotentialBase<VectorPotential<LinOp>, LinOp> {
32  // implementation of the kernel evaluation, which may be based on the
33  // information available from the superSpace
34  public:
35  VectorPotential() {}
36  Eigen::Matrix<typename PotentialReturnScalar<
38  std::complex<double>>::Scalar,
39  3, 1>
40  evaluateIntegrand_impl(const FunctionEvaluator<LinOp> &fun_ev,
41  const ElementTreeNode &element,
42  const Eigen::Vector3d &point,
43  const SurfacePoint &p) const {
44  // get evaluation points on unit square
45  auto s = p.segment<2>(0);
46 
47  // get quadrature weights
48  auto ws = p(2);
49 
50  // get points on geometry and tangential derivatives
51  auto x_f = p.segment<3>(3);
52 
53  // evaluate kernel
54  auto kernel = evaluateKernel(point, x_f);
55 
56  // assemble Galerkin solution
57  auto scalar_part = fun_ev.evaluate(element, p);
58 
59  // integrand without basis functions, note that the surface measure
60  // disappears for the divergence
61  auto integrand = kernel * scalar_part * ws;
62 
63  return integrand;
64  }
65 
69  std::complex<double> evaluateKernel(const Eigen::Vector3d &x,
70  const Eigen::Vector3d &y) const {
71  auto r = (x - y).norm();
72  return std::exp(-std::complex<double>(0., 1.) * wavenumber_ * r) / 4. /
73  BEMBEL_PI / r;
74  }
76  // setters
78  void set_wavenumber(std::complex<double> wavenumber) {
79  wavenumber_ = wavenumber;
80  }
82  // getters
84  std::complex<double> get_wavenumber() { return wavenumber_; }
85 
86  private:
87  std::complex<double> wavenumber_;
88 };
89 
90 } // namespace Bembel
91 #endif // BEMBEL_SRC_AUGMENTEDEFIE_VECTORPOTENTIAL_HPP_
The ElementTreeNode corresponds to an element in the element tree.
std::complex< double > evaluateKernel(const Eigen::Vector3d &x, const Eigen::Vector3d &y) const
Fundamental solution of Helmholtz problem.
Eigen::Matrix< double, 12, 1 > SurfacePoint
typedef of SurfacePoint
Routines for the evalutation of pointwise errors.
Definition: AnsatzSpace.hpp:14
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