11 #ifndef BEMBEL_SRC_LAPLACE_SINGLELAYERPOTENTIALGRADIENT_HPP_
12 #define BEMBEL_SRC_LAPLACE_SINGLELAYERPOTENTIALGRADIENT_HPP_
17 template <
typename LinOp>
18 class LaplaceSingleLayerPotentialGradient;
20 template <
typename LinOp>
22 typedef Eigen::VectorXd::Scalar Scalar;
23 static constexpr
int OutputSpaceDimension = 3;
29 template <
typename LinOp>
31 :
public PotentialBase<LaplaceSingleLayerPotentialGradient<LinOp>, LinOp> {
42 const Eigen::Vector3d &point,
45 auto s = p.segment<2>(0);
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);
56 auto x_kappa = x_f_dx.cross(x_f_dy).norm();
62 auto cauchy_value = fun_ev.evaluate(element, p);
65 auto integrand = kernel * cauchy_value * x_kappa * ws;
74 const Eigen::Vector3d &y)
const {
78 return -c / r3 / 4. / BEMBEL_PI;
The ElementTreeNode corresponds to an element in the element tree.
Eigen::VectorXd evaluateKernelGrad(const Eigen::Vector3d &x, const Eigen::Vector3d &y) const
Fundamental solution of Laplace problem.
Eigen::Matrix< double, 12, 1 > SurfacePoint
typedef of SurfacePoint
Routines for the evalutation of pointwise errors.
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.
Base case for specifying the return type of the potential.
struct containing specifications on the functional has to be specialized or derived for any particula...