11 #ifndef BEMBEL_SRC_ANSATZSPACE_FUNCTIONEVALUATOREVAL_HPP_
12 #define BEMBEL_SRC_ANSATZSPACE_FUNCTIONEVALUATOREVAL_HPP_
16 template <
typename Scalar,
unsigned int DF,
typename LinOp>
20 template <
typename Scalar,
typename LinOp>
23 getFunctionSpaceOutputDimension<DifferentialForm::Continuous>(),
26 const int polynomial_degree_plus_one_squared,
29 Scalar, Eigen::Dynamic,
30 getFunctionSpaceVectorDimension<DifferentialForm::Continuous>()>
32 auto s = p.segment<2>(0);
33 return coeff.transpose() * super_space.
basis(s) / element.
get_h();
38 template <
typename Scalar,
typename LinOp>
42 getFunctionSpaceOutputDimension<DifferentialForm::DivConforming>(), 1>
44 const int polynomial_degree_plus_one_squared,
47 Scalar, Eigen::Dynamic,
48 getFunctionSpaceVectorDimension<DifferentialForm::DivConforming>()>
50 auto s = p.segment<2>(0);
51 auto h = element.
get_h();
52 auto x_f_dx = p.segment<3>(6);
53 auto x_f_dy = p.segment<3>(9);
54 Eigen::Matrix<typename LinearOperatorTraits<LinOp>::Scalar, Eigen::Dynamic,
56 tangential_coefficients = coeff.transpose() * super_space.
basis(s);
57 return (x_f_dx * tangential_coefficients(0) +
58 x_f_dy * tangential_coefficients(1)) /
64 const int polynomial_degree_plus_one_squared,
67 Scalar, Eigen::Dynamic,
68 getFunctionSpaceVectorDimension<DifferentialForm::DivConforming>()>
70 auto s = p.segment<2>(0);
71 auto h = element.
get_h();
72 Eigen::Matrix<typename LinearOperatorTraits<LinOp>::Scalar, Eigen::Dynamic,
74 phiPhiVec_dx = super_space.
basisDx(s);
75 Eigen::Matrix<typename LinearOperatorTraits<LinOp>::Scalar, Eigen::Dynamic,
77 phiPhiVec_dy = super_space.
basisDy(s);
78 return (phiPhiVec_dx.dot(coeff.col(0)) + phiPhiVec_dy.dot(coeff.col(1))) /
84 template <
typename Scalar,
typename LinOp>
88 getFunctionSpaceOutputDimension<DifferentialForm::Discontinuous>(), 1>
90 const int polynomial_degree_plus_one_squared,
93 Scalar, Eigen::Dynamic,
94 getFunctionSpaceVectorDimension<DifferentialForm::Discontinuous>()>
96 auto s = p.segment<2>(0);
97 return coeff.transpose() * super_space.
basis(s) / element.
get_h();
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.
The superspace manages local polynomial bases on each element of the mesh and provides an interface t...
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > basisDx(const Eigen::Vector2d &s) const
Evaluate derivatives in x direction of local shape functions on the unit square at coordinate s.
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > basis(const Eigen::Vector2d &s) const
Evaluate local shape functions on the unit square at coordinate s.
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > basisDy(const Eigen::Vector2d &s) const
Evaluate derivatives in y direction of local shape functions on the unit square at coordinate s.