| Line | Branch | Exec | Source |
|---|---|---|---|
| 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_FUNCTIONEVALUATOR_HPP_ | ||
| 12 | #define BEMBEL_SRC_ANSATZSPACE_FUNCTIONEVALUATOR_HPP_ | ||
| 13 | |||
| 14 | namespace Bembel { | ||
| 15 | /** | ||
| 16 | * \ingroup AnsatzSpace | ||
| 17 | * \brief The FunctionEvaluator provides means to evaluate coefficient vectors | ||
| 18 | * as functions on the geometry. | ||
| 19 | */ | ||
| 20 | template <typename Derived> | ||
| 21 | class FunctionEvaluator { | ||
| 22 | public: | ||
| 23 | ////////////////////////////////////////////////////////////////////////////// | ||
| 24 | // constructors | ||
| 25 | ////////////////////////////////////////////////////////////////////////////// | ||
| 26 |
2/4✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 6 times.
✗ Branch 6 not taken.
|
6 | FunctionEvaluator() {} |
| 27 | FunctionEvaluator(const FunctionEvaluator &other) = default; | ||
| 28 | FunctionEvaluator(FunctionEvaluator &&other) = default; | ||
| 29 | 6 | FunctionEvaluator &operator=(FunctionEvaluator other) { | |
| 30 |
1/2✓ Branch 2 taken 6 times.
✗ Branch 3 not taken.
|
6 | ansatz_space_ = other.ansatz_space_; |
| 31 | 6 | fun_ = other.fun_; | |
| 32 | 6 | polynomial_degree_plus_one_squared_ = | |
| 33 | 6 | other.polynomial_degree_plus_one_squared_; | |
| 34 | 6 | return *this; | |
| 35 | } | ||
| 36 |
2/4✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 10 times.
✗ Branch 6 not taken.
|
10 | explicit FunctionEvaluator(const AnsatzSpace<Derived> &ansatz_space) { |
| 37 |
1/2✓ Branch 1 taken 10 times.
✗ Branch 2 not taken.
|
10 | init_FunctionEvaluator(ansatz_space); |
| 38 | 10 | return; | |
| 39 | ✗ | } | |
| 40 | 256 | FunctionEvaluator( | |
| 41 | const AnsatzSpace<Derived> &ansatz_space, | ||
| 42 | const Eigen::Matrix<typename LinearOperatorTraits<Derived>::Scalar, | ||
| 43 |
2/4✓ Branch 2 taken 128 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 128 times.
✗ Branch 6 not taken.
|
256 | Eigen::Dynamic, 1> &fun) { |
| 44 |
1/2✓ Branch 1 taken 128 times.
✗ Branch 2 not taken.
|
256 | init_FunctionEvaluator(ansatz_space, fun); |
| 45 | 256 | return; | |
| 46 | ✗ | } | |
| 47 | ////////////////////////////////////////////////////////////////////////////// | ||
| 48 | // init_Ansatzspace | ||
| 49 | ////////////////////////////////////////////////////////////////////////////// | ||
| 50 | 10 | void init_FunctionEvaluator(const AnsatzSpace<Derived> &ansatz_space) { | |
| 51 |
1/2✓ Branch 2 taken 10 times.
✗ Branch 3 not taken.
|
10 | ansatz_space_ = ansatz_space; |
| 52 | 10 | auto polynomial_degree = ansatz_space_.get_polynomial_degree(); | |
| 53 | 10 | polynomial_degree_plus_one_squared_ = | |
| 54 | 10 | (polynomial_degree + 1) * (polynomial_degree + 1); | |
| 55 | 20 | reordering_vector_ = ansatz_space_.get_superspace() | |
| 56 | 10 | .get_mesh() | |
| 57 | 10 | .get_element_tree() | |
| 58 | .computeReorderingVector(); | ||
| 59 | 10 | return; | |
| 60 | } | ||
| 61 | 256 | void init_FunctionEvaluator( | |
| 62 | const AnsatzSpace<Derived> &ansatz_space, | ||
| 63 | const Eigen::Matrix<typename LinearOperatorTraits<Derived>::Scalar, | ||
| 64 | Eigen::Dynamic, 1> &fun) { | ||
| 65 |
1/2✓ Branch 2 taken 128 times.
✗ Branch 3 not taken.
|
256 | ansatz_space_ = ansatz_space; |
| 66 |
1/2✓ Branch 2 taken 128 times.
✗ Branch 3 not taken.
|
256 | set_function(fun); |
| 67 | 256 | auto polynomial_degree = ansatz_space_.get_polynomial_degree(); | |
| 68 | 256 | polynomial_degree_plus_one_squared_ = | |
| 69 | 256 | (polynomial_degree + 1) * (polynomial_degree + 1); | |
| 70 | 512 | reordering_vector_ = ansatz_space_.get_superspace() | |
| 71 | 256 | .get_mesh() | |
| 72 | 256 | .get_element_tree() | |
| 73 | .computeReorderingVector(); | ||
| 74 | 256 | return; | |
| 75 | } | ||
| 76 | ////////////////////////////////////////////////////////////////////////////// | ||
| 77 | // evaluators | ||
| 78 | ////////////////////////////////////////////////////////////////////////////// | ||
| 79 | Eigen::Matrix< | ||
| 80 | typename LinearOperatorTraits<Derived>::Scalar, | ||
| 81 | getFunctionSpaceOutputDimension<LinearOperatorTraits<Derived>::Form>(), 1> | ||
| 82 | 7680 | evaluateOnPatch(int patch, const Eigen::Vector2d &ref_point) const { | |
| 83 | 7680 | const int elements_per_direction = | |
| 84 |
1/2✓ Branch 1 taken 3840 times.
✗ Branch 2 not taken.
|
7680 | (1 << ansatz_space_.get_refinement_level()); |
| 85 | 7680 | const int elements_per_patch = | |
| 86 | elements_per_direction * elements_per_direction; | ||
| 87 | 7680 | const double h = 1. / ((double)(elements_per_direction)); | |
| 88 |
1/2✓ Branch 1 taken 3840 times.
✗ Branch 2 not taken.
|
7680 | const int x_idx_ = std::floor(ref_point(0) / h); |
| 89 |
1/2✓ Branch 1 taken 3840 times.
✗ Branch 2 not taken.
|
7680 | const int y_idx_ = std::floor(ref_point(1) / h); |
| 90 | 7680 | const int x_idx = std::min(std::max(x_idx_, 0), elements_per_direction - 1); | |
| 91 | 7680 | const int y_idx = std::min(std::max(y_idx_, 0), elements_per_direction - 1); | |
| 92 | 7680 | const int tp_idx = | |
| 93 | 7680 | x_idx + elements_per_direction * y_idx + patch * elements_per_patch; | |
| 94 | 7680 | const int et_idx = reordering_vector_[tp_idx]; | |
| 95 | const ElementTree &et = | ||
| 96 | 7680 | ansatz_space_.get_superspace().get_mesh().get_element_tree(); | |
| 97 | 7680 | auto it = et.cpbegin(); | |
| 98 |
1/2✓ Branch 1 taken 3840 times.
✗ Branch 2 not taken.
|
7680 | std::advance(it, et_idx); |
| 99 | 7680 | const ElementTreeNode &element = *it; | |
| 100 | |||
| 101 |
1/2✓ Branch 1 taken 3840 times.
✗ Branch 2 not taken.
|
7680 | SurfacePoint sp; |
| 102 |
2/4✓ Branch 4 taken 3840 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 3840 times.
✗ Branch 8 not taken.
|
7680 | ansatz_space_.get_superspace().get_geometry()[patch].updateSurfacePoint( |
| 103 | &sp, ref_point, 1, element.mapToReferenceElement(ref_point)); | ||
| 104 |
1/2✓ Branch 1 taken 3840 times.
✗ Branch 2 not taken.
|
15360 | return evaluate(element, sp); |
| 105 | } | ||
| 106 | Eigen::Matrix< | ||
| 107 | typename LinearOperatorTraits<Derived>::Scalar, | ||
| 108 | getFunctionSpaceOutputDimension<LinearOperatorTraits<Derived>::Form>(), 1> | ||
| 109 | 1002182 | evaluate(const ElementTreeNode &element, const SurfacePoint &p) const { | |
| 110 | return eval_.eval( | ||
| 111 | 1002182 | ansatz_space_.get_superspace(), polynomial_degree_plus_one_squared_, | |
| 112 | element, p, | ||
| 113 | 3006546 | fun_.block(polynomial_degree_plus_one_squared_ * element.id_, 0, | |
| 114 | 1002182 | polynomial_degree_plus_one_squared_, | |
| 115 | getFunctionSpaceVectorDimension< | ||
| 116 |
2/4✓ Branch 2 taken 998342 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 998342 times.
✗ Branch 7 not taken.
|
2004364 | LinearOperatorTraits<Derived>::Form>())); |
| 117 | } | ||
| 118 | 91800 | typename LinearOperatorTraits<Derived>::Scalar evaluateDiv( | |
| 119 | const ElementTreeNode &element, const SurfacePoint &p) const { | ||
| 120 |
2/4✓ Branch 1 taken 91800 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 91800 times.
✗ Branch 5 not taken.
|
275400 | return eval_.evalDiv( |
| 121 | 91800 | ansatz_space_.get_superspace(), polynomial_degree_plus_one_squared_, | |
| 122 | element, p, | ||
| 123 | 183600 | fun_.block(polynomial_degree_plus_one_squared_ * element.id_, 0, | |
| 124 | 91800 | polynomial_degree_plus_one_squared_, | |
| 125 | getFunctionSpaceVectorDimension< | ||
| 126 | 183600 | LinearOperatorTraits<Derived>::Form>())); | |
| 127 | } | ||
| 128 | ////////////////////////////////////////////////////////////////////////////// | ||
| 129 | // setters | ||
| 130 | ////////////////////////////////////////////////////////////////////////////// | ||
| 131 | 266 | void set_function( | |
| 132 | Eigen::Matrix<typename LinearOperatorTraits<Derived>::Scalar, | ||
| 133 | Eigen::Dynamic, 1> | ||
| 134 | fun) { | ||
| 135 | 266 | const auto vec_dim = | |
| 136 | getFunctionSpaceVectorDimension<LinearOperatorTraits<Derived>::Form>(); | ||
| 137 |
2/4✓ Branch 2 taken 138 times.
✗ Branch 3 not taken.
✓ Branch 5 taken 138 times.
✗ Branch 6 not taken.
|
266 | auto longfun = (ansatz_space_.get_transformation_matrix() * fun).eval(); |
| 138 |
1/2✓ Branch 1 taken 138 times.
✗ Branch 2 not taken.
|
266 | fun_ = |
| 139 |
5/8✓ Branch 1 taken 138 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 66 times.
✓ Branch 5 taken 72 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 66 times.
✓ Branch 8 taken 72 times.
✗ Branch 9 not taken.
|
396 | Eigen::Map<Eigen::Matrix<typename LinearOperatorTraits<Derived>::Scalar, |
| 140 | Eigen::Dynamic, vec_dim>>( | ||
| 141 | 130 | longfun.data(), longfun.rows() / vec_dim, vec_dim); | |
| 142 | 266 | } | |
| 143 | ////////////////////////////////////////////////////////////////////////////// | ||
| 144 | // private member variables | ||
| 145 | ////////////////////////////////////////////////////////////////////////////// | ||
| 146 | private: | ||
| 147 | std::vector<int> reordering_vector_; | ||
| 148 | AnsatzSpace<Derived> ansatz_space_; | ||
| 149 | Eigen::Matrix< | ||
| 150 | typename LinearOperatorTraits<Derived>::Scalar, Eigen::Dynamic, | ||
| 151 | getFunctionSpaceVectorDimension<LinearOperatorTraits<Derived>::Form>()> | ||
| 152 | fun_; | ||
| 153 | int polynomial_degree_plus_one_squared_; | ||
| 154 | FunctionEvaluatorEval<typename LinearOperatorTraits<Derived>::Scalar, | ||
| 155 | LinearOperatorTraits<Derived>::Form, Derived> | ||
| 156 | eval_; | ||
| 157 | }; | ||
| 158 | |||
| 159 | } // namespace Bembel | ||
| 160 | #endif // BEMBEL_SRC_ANSATZSPACE_FUNCTIONEVALUATOR_HPP_ | ||
| 161 |