12 #ifndef BEMBEL_SRC_AUGMENTEDEFIE_EFIE_HPP_
13 #define BEMBEL_SRC_AUGMENTEDEFIE_EFIE_HPP_
18 template <
typename LinOp>
21 template <
typename LinOp>
23 typedef Eigen::VectorXcd::Scalar Scalar;
24 static constexpr
int OutputSpaceDimension = 3;
30 template <
typename LinOp>
38 std::complex<double>>::Scalar,
42 const Eigen::Vector3d &point,
45 auto s = p.segment<2>(0);
51 auto x_f = p.segment<3>(3);
54 auto h = element.
get_h();
61 auto scalar_part = fun_ev.evaluate(element, p);
62 auto divergence_part = fun_ev.evaluateDiv(element, p);
66 std::complex<double> omega =
67 wavenumber_ / std::sqrt(Constants::mu0 * Constants::eps0);
68 auto integrand = -(std::complex<double>(0., 1.) * omega * Constants::mu0 *
69 kernel * scalar_part -
70 kernel_gradient * divergence_part / Constants::eps0 /
71 std::complex<double>(0., 1.) / omega) *
81 const Eigen::Vector3d &y)
const {
82 auto r = (x - y).norm();
83 return std::exp(-std::complex<double>(0., 1.) * wavenumber_ * r) / 4. /
90 const Eigen::Vector3d &y)
const {
94 auto i = std::complex<double>(0., 1.);
95 return (std::exp(-i * wavenumber_ * r) * (-1. - i * wavenumber_ * r) / 4. /
102 void set_wavenumber(std::complex<double> wavenumber) {
103 wavenumber_ = wavenumber;
108 std::complex<double> get_wavenumber() {
return wavenumber_; }
111 std::complex<double> wavenumber_;
Eigen::VectorXcd evaluateKernelGrad(const Eigen::Vector3d &x, const Eigen::Vector3d &y) const
Gradient of fundamental solution of Helmholtz problem.
std::complex< double > evaluateKernel(const Eigen::Vector3d &x, const Eigen::Vector3d &y) const
Fundamental solution of Helmholtz problem.
The ElementTreeNode corresponds to an element in the element tree.
double get_h() const
getter
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...