Bembel
FunctionEvaluatorEval.hpp
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 #ifndef BEMBEL_SRC_ANSATZSPACE_FUNCTIONEVALUATOREVAL_HPP_
12 #define BEMBEL_SRC_ANSATZSPACE_FUNCTIONEVALUATOREVAL_HPP_
13 
14 namespace Bembel {
15 
16 template <typename Scalar, unsigned int DF, typename LinOp>
18 
19 // continuous
20 template <typename Scalar, typename LinOp>
21 struct FunctionEvaluatorEval<Scalar, DifferentialForm::Continuous, LinOp> {
22  Eigen::Matrix<Scalar,
23  getFunctionSpaceOutputDimension<DifferentialForm::Continuous>(),
24  1>
25  eval(const SuperSpace<LinOp> &super_space,
26  const int polynomial_degree_plus_one_squared,
27  const ElementTreeNode &element, const SurfacePoint &p,
28  const Eigen::Matrix<
29  Scalar, Eigen::Dynamic,
30  getFunctionSpaceVectorDimension<DifferentialForm::Continuous>()>
31  &coeff) const {
32  auto s = p.segment<2>(0);
33  return coeff.transpose() * super_space.basis(s) / element.get_h();
34  }
35 };
36 
37 // div-conforming
38 template <typename Scalar, typename LinOp>
39 struct FunctionEvaluatorEval<Scalar, DifferentialForm::DivConforming, LinOp> {
40  Eigen::Matrix<
41  Scalar,
42  getFunctionSpaceOutputDimension<DifferentialForm::DivConforming>(), 1>
43  eval(const SuperSpace<LinOp> &super_space,
44  const int polynomial_degree_plus_one_squared,
45  const ElementTreeNode &element, const SurfacePoint &p,
46  const Eigen::Matrix<
47  Scalar, Eigen::Dynamic,
48  getFunctionSpaceVectorDimension<DifferentialForm::DivConforming>()>
49  &coeff) const {
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,
55  1>
56  tangential_coefficients = coeff.transpose() * super_space.basis(s);
57  return (x_f_dx * tangential_coefficients(0) +
58  x_f_dy * tangential_coefficients(1)) /
59  h;
60  }
61 
62  Scalar evalDiv(
63  const SuperSpace<LinOp> &super_space,
64  const int polynomial_degree_plus_one_squared,
65  const ElementTreeNode &element, const SurfacePoint &p,
66  const Eigen::Matrix<
67  Scalar, Eigen::Dynamic,
68  getFunctionSpaceVectorDimension<DifferentialForm::DivConforming>()>
69  &coeff) const {
70  auto s = p.segment<2>(0);
71  auto h = element.get_h();
72  Eigen::Matrix<typename LinearOperatorTraits<LinOp>::Scalar, Eigen::Dynamic,
73  1>
74  phiPhiVec_dx = super_space.basisDx(s);
75  Eigen::Matrix<typename LinearOperatorTraits<LinOp>::Scalar, Eigen::Dynamic,
76  1>
77  phiPhiVec_dy = super_space.basisDy(s);
78  return (phiPhiVec_dx.dot(coeff.col(0)) + phiPhiVec_dy.dot(coeff.col(1))) /
79  h / h;
80  }
81 };
82 
83 // discontinuous
84 template <typename Scalar, typename LinOp>
85 struct FunctionEvaluatorEval<Scalar, DifferentialForm::Discontinuous, LinOp> {
86  Eigen::Matrix<
87  Scalar,
88  getFunctionSpaceOutputDimension<DifferentialForm::Discontinuous>(), 1>
89  eval(const SuperSpace<LinOp> &super_space,
90  const int polynomial_degree_plus_one_squared,
91  const ElementTreeNode &element, const SurfacePoint &p,
92  const Eigen::Matrix<
93  Scalar, Eigen::Dynamic,
94  getFunctionSpaceVectorDimension<DifferentialForm::Discontinuous>()>
95  &coeff) const {
96  auto s = p.segment<2>(0);
97  return coeff.transpose() * super_space.basis(s) / element.get_h();
98  }
99 };
100 } // namespace Bembel
101 #endif // BEMBEL_SRC_ANSATZSPACE_FUNCTIONEVALUATOREVAL_HPP_
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.
Definition: AnsatzSpace.hpp:14
Provides information about the discrete space required for the discretisation of a specific operator.
The superspace manages local polynomial bases on each element of the mesh and provides an interface t...
Definition: SuperSpace.hpp:22
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.
Definition: SuperSpace.hpp:330
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > basis(const Eigen::Vector2d &s) const
Evaluate local shape functions on the unit square at coordinate s.
Definition: SuperSpace.hpp:309
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.
Definition: SuperSpace.hpp:351