Bembel
HypersingularOperator.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 #ifndef BEMBEL_SRC_LAPLACE_HYPERSINGULAROPERATOR_HPP_
12 #define BEMBEL_SRC_LAPLACE_HYPERSINGULAROPERATOR_HPP_
13 
14 namespace Bembel {
15 // forward declaration of class LaplaceHypersingularOperator in order to define
16 // traits
17 class LaplaceHypersingularOperator;
18 
19 template <>
21  typedef Eigen::VectorXd EigenType;
22  typedef Eigen::VectorXd::Scalar Scalar;
23  enum {
24  OperatorOrder = 1,
25  Form = DifferentialForm::Continuous,
26  NumberOfFMMComponents = 2
27  };
28 };
29 
34  : public LinearOperatorBase<LaplaceHypersingularOperator> {
35  // implementation of the kernel evaluation, which may be based on the
36  // information available from the superSpace
37  public:
39  template <class T>
40  void evaluateIntegrand_impl(
41  const T &super_space, const SurfacePoint &p1, const SurfacePoint &p2,
42  Eigen::Matrix<
44  Eigen::Dynamic, Eigen::Dynamic> *intval) const {
45  auto polynomial_degree = super_space.get_polynomial_degree();
46  auto polynomial_degree_plus_one_squared =
47  (polynomial_degree + 1) * (polynomial_degree + 1);
48 
49  // get evaluation points on unit square
50  auto s = p1.segment<2>(0);
51  auto t = p2.segment<2>(0);
52 
53  // get quadrature weights
54  auto ws = p1(2);
55  auto wt = p2(2);
56 
57  // get points on geometry and tangential derivatives
58  auto x_f = p1.segment<3>(3);
59  auto x_f_dx = p1.segment<3>(6);
60  auto x_f_dy = p1.segment<3>(9);
61  auto y_f = p2.segment<3>(3);
62  auto y_f_dx = p2.segment<3>(6);
63  auto y_f_dy = p2.segment<3>(9);
64 
65  // compute surface measures from tangential derivatives
66  auto x_kappa = x_f_dx.cross(x_f_dy).norm();
67  auto y_kappa = y_f_dx.cross(y_f_dy).norm();
68 
69  // compute h
70  auto h = 1. / (1 << super_space.get_refinement_level()); // h = 1 ./ (2^M)
71 
72  // integrand without basis functions
73  auto integrand =
74  evaluateKernel(x_f, y_f) * x_kappa * y_kappa * ws * wt / h / h;
75 
76  // multiply basis functions with integrand and add to intval, this is an
77  // efficient implementation of
78  super_space.addScaledSurfaceCurlInteraction(intval, integrand, p1, p2);
79 
80  return;
81  }
82 
83  Eigen::Matrix<double, 2, 2> evaluateFMMInterpolation_impl(
84  const SurfacePoint &p1, const SurfacePoint &p2) const {
85  // get evaluation points on unit square
86  auto s = p1.segment<2>(0);
87  auto t = p2.segment<2>(0);
88 
89  // get points on geometry and tangential derivatives
90  auto x_f = p1.segment<3>(3);
91  auto x_f_dx = p1.segment<3>(6);
92  auto x_f_dy = p1.segment<3>(9);
93  auto y_f = p2.segment<3>(3);
94  auto y_f_dx = p2.segment<3>(6);
95  auto y_f_dy = p2.segment<3>(9);
96 
97  // evaluate kernel
98  auto kernel = evaluateKernel(x_f, y_f);
99 
100  // interpolation
101  Eigen::Matrix<double, 2, 2> intval;
102  intval.setZero();
103  intval(0, 0) = kernel * x_f_dy.dot(y_f_dy);
104  intval(0, 1) = -kernel * x_f_dy.dot(y_f_dx);
105  intval(1, 0) = -kernel * x_f_dx.dot(y_f_dy);
106  intval(1, 1) = kernel * x_f_dx.dot(y_f_dx);
107 
108  return intval;
109  }
110 
114  double evaluateKernel(const Eigen::Vector3d &x,
115  const Eigen::Vector3d &y) const {
116  return 1. / 4. / BEMBEL_PI / (x - y).norm();
117  }
118 };
119 
125 template <typename InterpolationPoints>
126 struct H2Multipole::Moment2D<InterpolationPoints,
128  static std::vector<Eigen::MatrixXd> compute2DMoment(
129  const SuperSpace<LaplaceHypersingularOperator> &super_space,
130  const int cluster_level, const int cluster_refinements,
131  const int number_of_points) {
132  Eigen::MatrixXd moment_dx = moment2DComputer<
135  super_space, cluster_level, cluster_refinements, number_of_points);
136  Eigen::MatrixXd moment_dy = moment2DComputer<
139  super_space, cluster_level, cluster_refinements, number_of_points);
140 
141  Eigen::MatrixXd moment(moment_dx.rows() + moment_dy.rows(),
142  moment_dx.cols());
143  moment << moment_dx, moment_dy;
144 
145  std::vector<Eigen::MatrixXd> vector_of_moments;
146  vector_of_moments.push_back(moment);
147 
148  return vector_of_moments;
149  }
150 };
151 
152 } // namespace Bembel
153 #endif // BEMBEL_SRC_LAPLACE_HYPERSINGULAROPERATOR_HPP_
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
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