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_MAXWELL_SINGLELAYEROPERATOR_HPP_
13 #define BEMBEL_SRC_MAXWELL_SINGLELAYEROPERATOR_HPP_
14 
15 namespace Bembel {
16 // forward declaration of class MaxwellSingleLayerOperator in order to define
17 // traits
18 class MaxwellSingleLayerOperator;
22 template <>
24  typedef Eigen::VectorXcd EigenType;
25  typedef Eigen::VectorXcd::Scalar Scalar;
26  enum {
27  OperatorOrder = -1,
28  Form = DifferentialForm::DivConforming,
29  NumberOfFMMComponents = 2
30  };
31 };
32 
39  : public LinearOperatorBase<MaxwellSingleLayerOperator> {
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 h = 1. / (1 << super_space.get_refinement_level()); // h = 1 ./ (2^M)
72 
73  // integrand without basis functions
74  auto kernel_evaluation = evaluateKernel(x_f, y_f) * ws * wt;
75  auto integrand_vector = kernel_evaluation;
76  auto integrand_divergence = -kernel_evaluation / wavenumber2_ / h / h;
77 
78  // vector part: mulitply shape functions with integrand and add to buffer
79  super_space.addScaledVectorBasisInteraction(intval, integrand_vector, s, t,
80  x_f_dx, x_f_dy, y_f_dx, y_f_dy);
81 
82  // divergence part: multiply shape functions with integrand and add to
83  // buffer
84  super_space.addScaledVectorBasisDivergenceInteraction(
85  intval, integrand_divergence, s, t);
86 
87  return;
88  }
89 
90  Eigen::Matrix<std::complex<double>, 4, 4> evaluateFMMInterpolation_impl(
91  const SurfacePoint &p1, const SurfacePoint &p2) const {
92  // get evaluation points on unit square
93  auto s = p1.segment<2>(0);
94  auto t = p2.segment<2>(0);
95 
96  // get points on geometry and tangential derivatives
97  auto x_f = p1.segment<3>(3);
98  auto x_f_dx = p1.segment<3>(6);
99  auto x_f_dy = p1.segment<3>(9);
100  auto y_f = p2.segment<3>(3);
101  auto y_f_dx = p2.segment<3>(6);
102  auto y_f_dy = p2.segment<3>(9);
103 
104  // evaluate kernel
105  auto kernel = evaluateKernel(x_f, y_f);
106 
107  // interpolation
108  Eigen::Matrix<std::complex<double>, 4, 4> intval;
109  intval.setZero();
110  intval(0, 0) = kernel * x_f_dx.dot(y_f_dx);
111  intval(0, 2) = kernel * x_f_dx.dot(y_f_dy);
112  intval(2, 0) = kernel * x_f_dy.dot(y_f_dx);
113  intval(2, 2) = kernel * x_f_dy.dot(y_f_dy);
114  intval(1, 1) = -kernel / wavenumber2_;
115  intval(1, 3) = -kernel / wavenumber2_;
116  intval(3, 1) = -kernel / wavenumber2_;
117  intval(3, 3) = -kernel / wavenumber2_;
118 
119  return intval;
120  }
121 
125  std::complex<double> evaluateKernel(const Eigen::Vector3d &x,
126  const Eigen::Vector3d &y) const {
127  auto r = (x - y).norm();
128  return std::exp(-std::complex<double>(0., 1.) * wavenumber_ * r) / 4. /
129  BEMBEL_PI / r;
130  }
132  // setters
134  void set_wavenumber(std::complex<double> wavenumber) {
135  wavenumber_ = wavenumber;
136  wavenumber2_ = wavenumber_ * wavenumber_;
137  }
139  // getters
141  std::complex<double> get_wavenumber() { return wavenumber_; }
142 
143  private:
144  std::complex<double> wavenumber_;
145  std::complex<double> wavenumber2_;
146 };
147 
153 template <typename InterpolationPoints>
154 struct H2Multipole::Moment2D<InterpolationPoints, MaxwellSingleLayerOperator> {
155  static std::vector<Eigen::MatrixXd> compute2DMoment(
156  const SuperSpace<MaxwellSingleLayerOperator> &super_space,
157  const int cluster_level, const int cluster_refinements,
158  const int number_of_points) {
159  Eigen::MatrixXd moment = moment2DComputer<
162  super_space, cluster_level, cluster_refinements, number_of_points);
163  Eigen::MatrixXd moment_dx = moment2DComputer<
166  super_space, cluster_level, cluster_refinements, number_of_points);
167  Eigen::MatrixXd moment_dy = moment2DComputer<
170  super_space, cluster_level, cluster_refinements, number_of_points);
171 
172  Eigen::MatrixXd moment1(moment.rows() + moment_dx.rows(), moment.cols());
173  moment1 << moment, moment_dx;
174  Eigen::MatrixXd moment2(moment.rows() + moment_dy.rows(), moment.cols());
175  moment2 << moment, moment_dy;
176 
177  std::vector<Eigen::MatrixXd> vector_of_moments;
178  vector_of_moments.push_back(moment1);
179  vector_of_moments.push_back(moment2);
180 
181  return vector_of_moments;
182  }
183 };
184 
185 } // namespace Bembel
186 #endif // BEMBEL_SRC_MAXWELL_SINGLELAYEROPERATOR_HPP_
This class implements the specification of the integration for the Electric Field Integral Equation.
std::complex< double > evaluateKernel(const Eigen::Vector3d &x, const Eigen::Vector3d &y) const
Fundamental solution of Helmholtz/Maxwell problem.
Eigen::Matrix< double, 12, 1 > SurfacePoint
typedef of SurfacePoint
Eigen::MatrixXd moment2DComputer(const SuperSpace< LinOp > &super_space, const int cluster_level, const int cluster_refinements, const int number_of_points)
Computes a single 2D moment for the FMM by tensorisation of the 1D moments.
Routines for the evalutation of pointwise errors.
Definition: AnsatzSpace.hpp:14
Computes 1D moment for FMM using derivatives of the basis functions. All calculations ar performed on...
Computes 1D moment for FMM. All calculations ar performed on [0,1].
Computes all 2D moment for the FMM by tensorisation of the 1D moments. Specialice this for your linea...
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...
The superspace manages local polynomial bases on each element of the mesh and provides an interface t...
Definition: SuperSpace.hpp:22