Bembel
DoubleLayerPotential.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_DOUBLELAYERPOTENTIAL_HPP_
12 #define BEMBEL_SRC_LAPLACE_DOUBLELAYERPOTENTIAL_HPP_
13 
14 namespace Bembel {
15 // forward declaration of class LaplaceDoubleLayerPotential in order to define
16 // traits
17 template <typename LinOp>
18 class LaplaceDoubleLayerPotential;
19 
20 template <typename LinOp>
22  typedef Eigen::VectorXd::Scalar Scalar;
23  static constexpr int OutputSpaceDimension = 1;
24 };
25 
29 template <typename LinOp>
31  : public PotentialBase<LaplaceDoubleLayerPotential<LinOp>, LinOp> {
32  // implementation of the kernel evaluation, which may be based on the
33  // information available from the superSpace
34  public:
36  Eigen::Matrix<
37  typename PotentialReturnScalar<
38  typename LinearOperatorTraits<LinOp>::Scalar, double>::Scalar,
39  1, 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  auto x_f_dx = p.segment<3>(6);
53  auto x_f_dy = p.segment<3>(9);
54 
55  // compute unnormalized normal from tangential derivatives
56  auto x_n = x_f_dx.cross(x_f_dy);
57 
58  // assemble Galerkin solution
59  auto cauchy_value = fun_ev.evaluate(element, p);
60 
61  // integrand without basis functions
62  auto integrand =
63  evaluateKernelGrad(point, x_f).dot(x_n) * cauchy_value * ws;
64 
65  return integrand;
66  }
67 
71  Eigen::Vector3d evaluateKernelGrad(const Eigen::Vector3d &x,
72  const Eigen::Vector3d &y) const {
73  auto c = x - y;
74  auto r = c.norm();
75  auto r3 = r * r * r;
76  return c / 4. / BEMBEL_PI / r3;
77  }
78 };
79 
80 } // namespace Bembel
81 #endif // BEMBEL_SRC_LAPLACE_DOUBLELAYERPOTENTIAL_HPP_
The ElementTreeNode corresponds to an element in the element tree.
Eigen::Vector3d evaluateKernelGrad(const Eigen::Vector3d &x, const Eigen::Vector3d &y) const
Gradient of 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