Bembel
SingleLayerOperator.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 //
12 #ifndef BEMBEL_SRC_HELMHOLTZ_SINGLELAYEROPERATOR_HPP_
13 #define BEMBEL_SRC_HELMHOLTZ_SINGLELAYEROPERATOR_HPP_
14 
15 namespace Bembel {
16 // forward declaration of class HelmholtzSingleLayerOperator in order to define
17 // traits
18 class HelmholtzSingleLayerOperator;
22 template <>
24  typedef Eigen::VectorXcd EigenType;
25  typedef Eigen::VectorXcd::Scalar Scalar;
26  enum {
27  OperatorOrder = -1,
28  Form = DifferentialForm::Discontinuous,
29  NumberOfFMMComponents = 1
30  };
31 };
32 
39  : public LinearOperatorBase<HelmholtzSingleLayerOperator> {
40  // implementation of the kernel evaluation, which may be based on the
41  // information available from the superSpace
42  public:
44  template <class T>
45  void evaluateIntegrand_impl(
46  const T &super_space, const SurfacePoint &p1, const SurfacePoint &p2,
47  Eigen::Matrix<
49  Eigen::Dynamic, Eigen::Dynamic> *intval) const {
50  auto polynomial_degree = super_space.get_polynomial_degree();
51  auto polynomial_degree_plus_one_squared =
52  (polynomial_degree + 1) * (polynomial_degree + 1);
53 
54  // get evaluation points on unit square
55  auto s = p1.segment<2>(0);
56  auto t = p2.segment<2>(0);
57 
58  // get quadrature weights
59  auto ws = p1(2);
60  auto wt = p2(2);
61 
62  // get points on geometry and tangential derivatives
63  auto x_f = p1.segment<3>(3);
64  auto x_f_dx = p1.segment<3>(6);
65  auto x_f_dy = p1.segment<3>(9);
66  auto y_f = p2.segment<3>(3);
67  auto y_f_dx = p2.segment<3>(6);
68  auto y_f_dy = p2.segment<3>(9);
69 
70  // compute surface measures from tangential derivatives
71  auto x_kappa = x_f_dx.cross(x_f_dy).norm();
72  auto y_kappa = y_f_dx.cross(y_f_dy).norm();
73 
74  // integrand without basis functions
75  auto integrand = evaluateKernel(x_f, y_f) * x_kappa * y_kappa * ws * wt;
76 
77  // multiply basis functions with integrand and add to intval, this is an
78  // efficient implementation of
79  // (*intval) += super_space.BasisInteraction(s, t) * evaluateKernel(x_f,
80  // y_f)
81  // * x_kappa * y_kappa * ws * wt;
82  super_space.addScaledBasisInteraction(intval, integrand, s, t);
83 
84  return;
85  }
86 
87  Eigen::Matrix<std::complex<double>, 1, 1> evaluateFMMInterpolation_impl(
88  const SurfacePoint &p1, const SurfacePoint &p2) const {
89  // get evaluation points on unit square
90  auto s = p1.segment<2>(0);
91  auto t = p2.segment<2>(0);
92 
93  // get points on geometry and tangential derivatives
94  auto x_f = p1.segment<3>(3);
95  auto x_f_dx = p1.segment<3>(6);
96  auto x_f_dy = p1.segment<3>(9);
97  auto y_f = p2.segment<3>(3);
98  auto y_f_dx = p2.segment<3>(6);
99  auto y_f_dy = p2.segment<3>(9);
100 
101  // compute surface measures from tangential derivatives
102  auto x_kappa = x_f_dx.cross(x_f_dy).norm();
103  auto y_kappa = y_f_dx.cross(y_f_dy).norm();
104 
105  // interpolation
106  Eigen::Matrix<std::complex<double>, 1, 1> intval;
107  intval(0) = evaluateKernel(x_f, y_f) * x_kappa * y_kappa;
108 
109  return intval;
110  }
111 
115  std::complex<double> evaluateKernel(const Eigen::Vector3d &x,
116  const Eigen::Vector3d &y) const {
117  auto r = (x - y).norm();
118  return std::exp(-std::complex<double>(0., 1.) * wavenumber_ * r) / 4. /
119  BEMBEL_PI / r;
120  }
122  // setters
124  void set_wavenumber(std::complex<double> wavenumber) {
125  wavenumber_ = wavenumber;
126  }
128  // getters
130  std::complex<double> get_wavenumber() { return wavenumber_; }
131 
132  private:
133  std::complex<double> wavenumber_;
134 };
135 
136 } // namespace Bembel
137 #endif // BEMBEL_SRC_HELMHOLTZ_SINGLELAYEROPERATOR_HPP_
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 Helmholtz problem.
Eigen::Matrix< double, 12, 1 > SurfacePoint
typedef of SurfacePoint
Routines for the evalutation of pointwise errors.
Definition: AnsatzSpace.hpp:14
linear operator base class. this serves as a common interface for existing linear operators.
struct containing specifications on the linear operator has to be specialized or derived for any part...