Line data Source code
1 0 : // Distributed under the MIT License. 2 : // See LICENSE.txt for details. 3 : 4 : #pragma once 5 : 6 : #include <array> 7 : #include <cstddef> 8 : 9 : #include "DataStructures/Matrix.hpp" 10 : #include "DataStructures/Tensor/TypeAliases.hpp" 11 : 12 : /// \cond 13 : class DataVector; 14 : template <size_t Dim> 15 : class Mesh; 16 : namespace PUP { 17 : class er; 18 : } // namespace PUP 19 : /// \endcond 20 : 21 : namespace intrp { 22 : /*! 23 : * \brief Interpolates by doing partial summation in each dimension using 24 : * one-dimensional interpolation 25 : * 26 : * \details The one-dimensional matrices used to do the interpolation depend 27 : * upon the Spectral::Basis used in each dimension: 28 : * - For a Chebyshev or Legendre basis, the matrices are given by 29 : * Spectral::fornberg_interpolation_matrix at the quadrature points of the 30 : * source_mesh. (These are equivalent to those returned by 31 : * Spectral::interpolation_matrix.) 32 : * - For a Fourier basis, the matrix is given by 33 : * Spectral::fourier_interpolation_matrix at the quadrature points of the 34 : * source_mesh 35 : * 36 : * For multidimensional bases such as the SphericalHarmonic or ZernikeB2, the 37 : * matrices used for interpolating cannot be applied per dimension but must 38 : * be handled specially. 39 : * 40 : */ 41 : template <size_t Dim> 42 1 : class Cardinal { 43 : public: 44 0 : Cardinal( 45 : const Mesh<Dim>& source_mesh, 46 : const tnsr::I<DataVector, Dim, Frame::ElementLogical>& target_points); 47 0 : Cardinal(const Mesh<Dim>& source_mesh, 48 : const tnsr::I<double, Dim, Frame::ElementLogical>& target_point); 49 : 50 0 : Cardinal(); 51 : 52 : /// Interpolates the function `f` provided on the `source_mesh` to the 53 : /// `target_points` with which the interpolator was constructed. 54 1 : DataVector interpolate(const DataVector& f) const; 55 : 56 : /// The one-dimensional interpolation matrices used to do the interpolation 57 1 : const std::array<Matrix, Dim>& interpolation_matrices() const; 58 : 59 : // NOLINTNEXTLINE(google-runtime-references) 60 0 : void pup(PUP::er& p); 61 : 62 : private: 63 : /// Precomputes `zernike_weights_`, which is all the work independent of 64 : /// `f_source`, to avoid redundant computations. This is only needed when 65 : /// the source mesh has a B2 basis 66 1 : void set_zernike_b2_weights(); 67 : 68 : /// General routine called by `interpolate()` for a mesh using B2 bases 69 1 : DataVector interpolate_zernike_b2(const DataVector& f_source) const; 70 : 71 : template <size_t LocalDim> 72 : // NOLINTNEXTILNE(readability-redundant-declaration) 73 0 : friend bool operator==(const Cardinal<LocalDim>& lhs, 74 : const Cardinal<LocalDim>& rhs); 75 : 76 0 : size_t n_target_points_ = 0; 77 0 : Mesh<Dim> source_mesh_{}; 78 0 : std::array<Matrix, Dim> interpolation_matrices_{}; 79 0 : bool using_spherical_harmonics_{false}; 80 0 : bool using_zernike_b2_{false}; 81 0 : Matrix zernike_weights_{}; 82 : }; 83 : 84 : template <size_t Dim> 85 0 : bool operator==(const Cardinal<Dim>& lhs, const Cardinal<Dim>& rhs); 86 : 87 : template <size_t Dim> 88 0 : bool operator!=(const Cardinal<Dim>& lhs, const Cardinal<Dim>& rhs); 89 : } // namespace intrp