| Line | Branch | Exec | Source |
|---|---|---|---|
| 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_LAPLACE_SINGLELAYEROPERATOR_HPP_ | ||
| 13 | #define BEMBEL_SRC_LAPLACE_SINGLELAYEROPERATOR_HPP_ | ||
| 14 | |||
| 15 | namespace Bembel { | ||
| 16 | // forward declaration of class LaplaceSingleLayerOperator in order to define | ||
| 17 | // traits | ||
| 18 | class LaplaceSingleLayerOperator; | ||
| 19 | /** | ||
| 20 | * \brief Specification of the LinerOperatorTraits for Laplace. | ||
| 21 | */ | ||
| 22 | template <> | ||
| 23 | struct LinearOperatorTraits<LaplaceSingleLayerOperator> { | ||
| 24 | typedef Eigen::VectorXd EigenType; | ||
| 25 | typedef Eigen::VectorXd::Scalar Scalar; | ||
| 26 | enum { | ||
| 27 | OperatorOrder = -1, | ||
| 28 | Form = DifferentialForm::Discontinuous, | ||
| 29 | NumberOfFMMComponents = 1 | ||
| 30 | }; | ||
| 31 | }; | ||
| 32 | |||
| 33 | /** | ||
| 34 | * \ingroup Laplace | ||
| 35 | * \brief This class implements the specification of the integration for the | ||
| 36 | * single layer operator for Laplace. | ||
| 37 | */ | ||
| 38 | class LaplaceSingleLayerOperator | ||
| 39 | : public LinearOperatorBase<LaplaceSingleLayerOperator> { | ||
| 40 | // implementation of the kernel evaluation, which may be based on the | ||
| 41 | // information available from the superSpace | ||
| 42 | public: | ||
| 43 | 2 | LaplaceSingleLayerOperator() {} | |
| 44 | template <class T> | ||
| 45 | 795097 | void evaluateIntegrand_impl( | |
| 46 | const T &super_space, const SurfacePoint &p1, const SurfacePoint &p2, | ||
| 47 | Eigen::Matrix< | ||
| 48 | typename LinearOperatorTraits<LaplaceSingleLayerOperator>::Scalar, | ||
| 49 | Eigen::Dynamic, Eigen::Dynamic> *intval) const { | ||
| 50 | 795097 | auto polynomial_degree = super_space.get_polynomial_degree(); | |
| 51 | 795097 | auto polynomial_degree_plus_one_squared = | |
| 52 | 795097 | (polynomial_degree + 1) * (polynomial_degree + 1); | |
| 53 | |||
| 54 | // get evaluation points on unit square | ||
| 55 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto s = p1.segment<2>(0); |
| 56 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto t = p2.segment<2>(0); |
| 57 | |||
| 58 | // get quadrature weights | ||
| 59 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto ws = p1(2); |
| 60 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto wt = p2(2); |
| 61 | |||
| 62 | // get points on geometry and tangential derivatives | ||
| 63 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto x_f = p1.segment<3>(3); |
| 64 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto x_f_dx = p1.segment<3>(6); |
| 65 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto x_f_dy = p1.segment<3>(9); |
| 66 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto y_f = p2.segment<3>(3); |
| 67 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto y_f_dx = p2.segment<3>(6); |
| 68 |
1/2✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
|
795097 | auto y_f_dy = p2.segment<3>(9); |
| 69 | |||
| 70 | // compute surface measures from tangential derivatives | ||
| 71 |
2/4✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 795097 times.
✗ Branch 5 not taken.
|
795097 | auto x_kappa = x_f_dx.cross(x_f_dy).norm(); |
| 72 |
2/4✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 795097 times.
✗ Branch 5 not taken.
|
795097 | auto y_kappa = y_f_dx.cross(y_f_dy).norm(); |
| 73 | |||
| 74 | // integrand without basis functions | ||
| 75 |
3/6✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 795097 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 795097 times.
✗ Branch 8 not taken.
|
795097 | 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 |
3/6✓ Branch 1 taken 795097 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 795097 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 795097 times.
✗ Branch 8 not taken.
|
795097 | super_space.addScaledBasisInteraction(intval, integrand, s, t); |
| 83 | |||
| 84 | 1590194 | return; | |
| 85 | } | ||
| 86 | |||
| 87 | 1259712 | Eigen::Matrix<double, 1, 1> evaluateFMMInterpolation_impl( | |
| 88 | const SurfacePoint &p1, const SurfacePoint &p2) const { | ||
| 89 | // get evaluation points on unit square | ||
| 90 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | auto s = p1.segment<2>(0); |
| 91 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | auto t = p2.segment<2>(0); |
| 92 | |||
| 93 | // get points on geometry and tangential derivatives | ||
| 94 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | auto x_f = p1.segment<3>(3); |
| 95 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | auto x_f_dx = p1.segment<3>(6); |
| 96 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | auto x_f_dy = p1.segment<3>(9); |
| 97 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | auto y_f = p2.segment<3>(3); |
| 98 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | auto y_f_dx = p2.segment<3>(6); |
| 99 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | auto y_f_dy = p2.segment<3>(9); |
| 100 | |||
| 101 | // compute surface measures from tangential derivatives | ||
| 102 |
2/4✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1259712 times.
✗ Branch 5 not taken.
|
1259712 | auto x_kappa = x_f_dx.cross(x_f_dy).norm(); |
| 103 |
2/4✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1259712 times.
✗ Branch 5 not taken.
|
1259712 | auto y_kappa = y_f_dx.cross(y_f_dy).norm(); |
| 104 | |||
| 105 | // interpolation | ||
| 106 |
1/2✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
|
1259712 | Eigen::Matrix<double, 1, 1> intval; |
| 107 |
4/8✓ Branch 1 taken 1259712 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1259712 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1259712 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1259712 times.
✗ Branch 11 not taken.
|
1259712 | intval(0) = evaluateKernel(x_f, y_f) * x_kappa * y_kappa; |
| 108 | |||
| 109 | 2519424 | return intval; | |
| 110 | } | ||
| 111 | |||
| 112 | /** | ||
| 113 | * \brief Fundamental solution of Laplace problem | ||
| 114 | */ | ||
| 115 | 2054809 | double evaluateKernel(const Eigen::Vector3d &x, | |
| 116 | const Eigen::Vector3d &y) const { | ||
| 117 |
1/2✓ Branch 2 taken 2054809 times.
✗ Branch 3 not taken.
|
2054809 | return 1. / 4. / BEMBEL_PI / (x - y).norm(); |
| 118 | } | ||
| 119 | }; | ||
| 120 | |||
| 121 | } // namespace Bembel | ||
| 122 | #endif // BEMBEL_SRC_LAPLACE_SINGLELAYEROPERATOR_HPP_ | ||
| 123 |