11 #ifndef BEMBEL_SRC_LAPLACE_HYPERSINGULAROPERATOR_HPP_
12 #define BEMBEL_SRC_LAPLACE_HYPERSINGULAROPERATOR_HPP_
17 class LaplaceHypersingularOperator;
21 typedef Eigen::VectorXd EigenType;
22 typedef Eigen::VectorXd::Scalar Scalar;
25 Form = DifferentialForm::Continuous,
26 NumberOfFMMComponents = 2
40 void evaluateIntegrand_impl(
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);
50 auto s = p1.segment<2>(0);
51 auto t = p2.segment<2>(0);
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);
66 auto x_kappa = x_f_dx.cross(x_f_dy).norm();
67 auto y_kappa = y_f_dx.cross(y_f_dy).norm();
70 auto h = 1. / (1 << super_space.get_refinement_level());
78 super_space.addScaledSurfaceCurlInteraction(intval, integrand, p1, p2);
83 Eigen::Matrix<double, 2, 2> evaluateFMMInterpolation_impl(
86 auto s = p1.segment<2>(0);
87 auto t = p2.segment<2>(0);
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);
101 Eigen::Matrix<double, 2, 2> intval;
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);
115 const Eigen::Vector3d &y)
const {
116 return 1. / 4. / BEMBEL_PI / (x - y).norm();
125 template <
typename InterpolationPo
ints>
128 static std::vector<Eigen::MatrixXd> compute2DMoment(
130 const int cluster_level,
const int cluster_refinements,
131 const int number_of_points) {
135 super_space, cluster_level, cluster_refinements, number_of_points);
139 super_space, cluster_level, cluster_refinements, number_of_points);
141 Eigen::MatrixXd moment(moment_dx.rows() + moment_dy.rows(),
143 moment << moment_dx, moment_dy;
145 std::vector<Eigen::MatrixXd> vector_of_moments;
146 vector_of_moments.push_back(moment);
148 return vector_of_moments;
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.
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...