Bembel
SingleLayerPotential.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_HELMHOLTZ_SINGLELAYERPOTENTIAL_HPP_
12 #define BEMBEL_SRC_HELMHOLTZ_SINGLELAYERPOTENTIAL_HPP_
13 
14 namespace Bembel {
15 // forward declaration of class HelmholtzSingleLayerPotential in order to define
16 // traits
17 template <typename LinOp>
18 class HelmholtzSingleLayerPotential;
19 
23 template <typename LinOp>
25  typedef Eigen::VectorXcd::Scalar Scalar;
26  static constexpr int OutputSpaceDimension = 1;
27 };
28 
34 template <typename LinOp>
36  : public PotentialBase<HelmholtzSingleLayerPotential<LinOp>, LinOp> {
37  // implementation of the kernel evaluation, which may be based on the
38  // information available from the superSpace
39  public:
41  Eigen::Matrix<typename PotentialReturnScalar<
43  std::complex<double>>::Scalar,
44  1, 1>
45  evaluateIntegrand_impl(const FunctionEvaluator<LinOp> &fun_ev,
46  const ElementTreeNode &element,
47  const Eigen::Vector3d &point,
48  const SurfacePoint &p) const {
49  // get evaluation points on unit square
50  auto s = p.segment<2>(0);
51 
52  // get quadrature weights
53  auto ws = p(2);
54 
55  // get points on geometry and tangential derivatives
56  auto x_f = p.segment<3>(3);
57  auto x_f_dx = p.segment<3>(6);
58  auto x_f_dy = p.segment<3>(9);
59 
60  // compute surface measures from tangential derivatives
61  auto x_kappa = x_f_dx.cross(x_f_dy).norm();
62 
63  // evaluate kernel
64  auto kernel = evaluateKernel(point, x_f);
65 
66  // assemble Galerkin solution
67  auto cauchy_value = fun_ev.evaluate(element, p);
68 
69  // integrand without basis functions
70  auto integrand = kernel * cauchy_value * x_kappa * ws;
71 
72  return integrand;
73  }
74 
78  std::complex<double> evaluateKernel(const Eigen::Vector3d &x,
79  const Eigen::Vector3d &y) const {
80  auto r = (x - y).norm();
81  return std::exp(-std::complex<double>(0., 1.) * wavenumber_ * r) / 4. /
82  BEMBEL_PI / r;
83  }
85  // setters
87  void set_wavenumber(std::complex<double> wavenumber) {
88  wavenumber_ = wavenumber;
89  }
91  // getters
93  std::complex<double> get_wavenumber() { return wavenumber_; }
94 
95  private:
96  std::complex<double> wavenumber_;
97 };
98 
99 } // namespace Bembel
100 #endif // BEMBEL_SRC_HELMHOLTZ_SINGLELAYERPOTENTIAL_HPP_
The ElementTreeNode corresponds to an element in the element tree.
This class implements the specification of the integration for the single layer potential for Helmhol...
std::complex< double > evaluateKernel(const Eigen::Vector3d &x, const Eigen::Vector3d &y) const
Fundamental solution of Laplace 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