| 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_HOMOGENISEDLAPLACE_SINGLELAYEROPERATOR_HPP_ | ||
| 13 | #define BEMBEL_SRC_HOMOGENISEDLAPLACE_SINGLELAYEROPERATOR_HPP_ | ||
| 14 | |||
| 15 | namespace Bembel { | ||
| 16 | // forward declaration of class HomogenisedLaplaceSingleLayerOperator | ||
| 17 | // in order to define traits | ||
| 18 | class HomogenisedLaplaceSingleLayerOperator; | ||
| 19 | |||
| 20 | /** | ||
| 21 | * \brief Specification of the LinerOperatorTraits for the Homogenised Laplace. | ||
| 22 | */ | ||
| 23 | template<> | ||
| 24 | struct LinearOperatorTraits<HomogenisedLaplaceSingleLayerOperator> { | ||
| 25 | typedef Eigen::VectorXd EigenType; | ||
| 26 | typedef Eigen::VectorXd::Scalar Scalar; | ||
| 27 | enum { | ||
| 28 | OperatorOrder = -1, | ||
| 29 | Form = DifferentialForm::Discontinuous, | ||
| 30 | NumberOfFMMComponents = 1 | ||
| 31 | }; | ||
| 32 | }; | ||
| 33 | |||
| 34 | /** | ||
| 35 | * \ingroup HomogenisedLaplace | ||
| 36 | * \brief This class implements the specification of the integration for the | ||
| 37 | * single layer operator for the homogenised Laplace. | ||
| 38 | */ | ||
| 39 | class HomogenisedLaplaceSingleLayerOperator : public LinearOperatorBase< | ||
| 40 | HomogenisedLaplaceSingleLayerOperator> { | ||
| 41 | // implementation of the kernel evaluation, which may be based on the | ||
| 42 | // information available from the superSpace | ||
| 43 | public: | ||
| 44 | /** | ||
| 45 | * \brief Constructs an object initialising the coefficients and the degree | ||
| 46 | * via the static variable precision. | ||
| 47 | */ | ||
| 48 | 2 | HomogenisedLaplaceSingleLayerOperator() { | |
| 49 | 2 | this->deg = getDegree(HomogenisedLaplaceSingleLayerOperator::precision); | |
| 50 |
1/2✓ Branch 1 taken 2 times.
✗ Branch 2 not taken.
|
4 | this->cs = getCoefficients( |
| 51 | 2 | HomogenisedLaplaceSingleLayerOperator::precision); | |
| 52 | 2 | } | |
| 53 | |||
| 54 | template<class T> | ||
| 55 | 787320 | void evaluateIntegrand_impl(const T &super_space, const SurfacePoint &p1, | |
| 56 | const SurfacePoint &p2, | ||
| 57 | Eigen::Matrix< | ||
| 58 | typename LinearOperatorTraits<HomogenisedLaplaceSingleLayerOperator | ||
| 59 | >::Scalar, Eigen::Dynamic, Eigen::Dynamic> *intval) const { | ||
| 60 | 787320 | auto polynomial_degree = super_space.get_polynomial_degree(); | |
| 61 | 787320 | auto polynomial_degree_plus_one_squared = (polynomial_degree + 1) | |
| 62 | 787320 | * (polynomial_degree + 1); | |
| 63 | |||
| 64 | // get evaluation points on unit square | ||
| 65 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto s = p1.segment < 2 > (0); |
| 66 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto t = p2.segment < 2 > (0); |
| 67 | |||
| 68 | // get quadrature weights | ||
| 69 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto ws = p1(2); |
| 70 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto wt = p2(2); |
| 71 | |||
| 72 | // get points on geometry and tangential derivatives | ||
| 73 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto x_f = p1.segment < 3 > (3); |
| 74 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto x_f_dx = p1.segment < 3 > (6); |
| 75 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto x_f_dy = p1.segment < 3 > (9); |
| 76 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto y_f = p2.segment < 3 > (3); |
| 77 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto y_f_dx = p2.segment < 3 > (6); |
| 78 |
1/2✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
|
787320 | auto y_f_dy = p2.segment < 3 > (9); |
| 79 | |||
| 80 | // compute surface measures from tangential derivatives | ||
| 81 |
2/4✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 787320 times.
✗ Branch 5 not taken.
|
787320 | auto x_kappa = x_f_dx.cross(x_f_dy).norm(); |
| 82 |
2/4✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 787320 times.
✗ Branch 5 not taken.
|
787320 | auto y_kappa = y_f_dx.cross(y_f_dy).norm(); |
| 83 | |||
| 84 | // integrand without basis functions | ||
| 85 |
3/6✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 787320 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 787320 times.
✗ Branch 8 not taken.
|
787320 | auto integrand = evaluateKernel(x_f, y_f) * x_kappa * y_kappa * ws * wt; |
| 86 | |||
| 87 | // multiply basis functions with integrand and add to intval, this is an | ||
| 88 | // efficient implementation of | ||
| 89 | // (*intval) += super_space.basisInteraction(s, t) | ||
| 90 | // * evaluateKernel(x_f, y_f) * x_kappa * y_kappa * ws * wt; | ||
| 91 |
3/6✓ Branch 1 taken 787320 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 787320 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 787320 times.
✗ Branch 8 not taken.
|
787320 | super_space.addScaledBasisInteraction(intval, integrand, s, t); |
| 92 | |||
| 93 | 1574640 | return; | |
| 94 | } | ||
| 95 | |||
| 96 | 1889568 | Eigen::Matrix<double, 1, 1> evaluateFMMInterpolation_impl( | |
| 97 | const SurfacePoint &p1, const SurfacePoint &p2) const { | ||
| 98 | // get evaluation points on unit square | ||
| 99 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | auto s = p1.segment < 2 > (0); |
| 100 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | auto t = p2.segment < 2 > (0); |
| 101 | |||
| 102 | // get points on geometry and tangential derivatives | ||
| 103 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | auto x_f = p1.segment < 3 > (3); |
| 104 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | auto x_f_dx = p1.segment < 3 > (6); |
| 105 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | auto x_f_dy = p1.segment < 3 > (9); |
| 106 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | auto y_f = p2.segment < 3 > (3); |
| 107 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | auto y_f_dx = p2.segment < 3 > (6); |
| 108 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | auto y_f_dy = p2.segment < 3 > (9); |
| 109 | |||
| 110 | // compute surface measures from tangential derivatives | ||
| 111 |
2/4✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1889568 times.
✗ Branch 5 not taken.
|
1889568 | auto x_kappa = x_f_dx.cross(x_f_dy).norm(); |
| 112 |
2/4✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1889568 times.
✗ Branch 5 not taken.
|
1889568 | auto y_kappa = y_f_dx.cross(y_f_dy).norm(); |
| 113 | |||
| 114 | // interpolation | ||
| 115 |
1/2✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
|
1889568 | Eigen::Matrix<double, 1, 1> intval; |
| 116 |
4/8✓ Branch 1 taken 1889568 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 1889568 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1889568 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 1889568 times.
✗ Branch 11 not taken.
|
1889568 | intval(0) = evaluateKernel(x_f, y_f) * x_kappa * y_kappa; |
| 117 | |||
| 118 | 3779136 | return intval; | |
| 119 | } | ||
| 120 | |||
| 121 | /** | ||
| 122 | * \brief Fundamental solution of the Homogenised Laplace problem | ||
| 123 | */ | ||
| 124 | 2676888 | double evaluateKernel(const Eigen::Vector3d &x, | |
| 125 | const Eigen::Vector3d &y) const { | ||
| 126 |
2/4✓ Branch 2 taken 2676888 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 2676888 times.
✗ Branch 6 not taken.
|
2676888 | return k_mod(x - y) |
| 127 |
4/8✓ Branch 1 taken 2676888 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2676888 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 2676888 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 2676888 times.
✗ Branch 11 not taken.
|
2676888 | + evaluate_solid_sphericals(x - y, this->cs, this->deg, false); |
| 128 | } | ||
| 129 | |||
| 130 | /** | ||
| 131 | * \brief sets the precision of the periodicity of the kernel | ||
| 132 | */ | ||
| 133 | 2 | static void setPrecision(double p) { | |
| 134 | 2 | HomogenisedLaplaceSingleLayerOperator::precision = p; | |
| 135 | 2 | } | |
| 136 | |||
| 137 | /** | ||
| 138 | * \brief returns the precision of the periodicity of the kernel | ||
| 139 | */ | ||
| 140 | 4 | static double getPrecision() { | |
| 141 | 4 | return HomogenisedLaplaceSingleLayerOperator::precision; | |
| 142 | } | ||
| 143 | |||
| 144 | private: | ||
| 145 | /** The degree of the spherical harmonics expansion */ | ||
| 146 | unsigned int deg; | ||
| 147 | /** The coefficients of the spherical harmonics expansion */ | ||
| 148 | Eigen::VectorXd cs; | ||
| 149 | /** The precision of the periodicity of the kernel */ | ||
| 150 | static double precision; | ||
| 151 | }; | ||
| 152 | |||
| 153 | double HomogenisedLaplaceSingleLayerOperator::precision = 0; | ||
| 154 | |||
| 155 | } // namespace Bembel | ||
| 156 | |||
| 157 | #endif // BEMBEL_SRC_HOMOGENISEDLAPLACE_SINGLELAYEROPERATOR_HPP_ | ||
| 158 |