41 #ifndef DTK_SPLINEEVALUATIONMATRIX_HPP 42 #define DTK_SPLINEEVALUATIONMATRIX_HPP 48 #include <Teuchos_ArrayView.hpp> 49 #include <Teuchos_RCP.hpp> 51 #include <Tpetra_CrsMatrix.hpp> 52 #include <Tpetra_Map.hpp> 53 #include <Tpetra_MultiVector.hpp> 54 #include <Tpetra_Operator.hpp> 64 template <
class Basis,
int DIM>
75 const Teuchos::RCP<
const Tpetra::Map<int, SupportId>> &domain_map,
76 const Teuchos::RCP<
const Tpetra::Map<int, SupportId>> &range_map,
77 const Teuchos::ArrayView<const double> &target_centers,
78 const Teuchos::ArrayView<const SupportId> &target_center_gids,
79 const Teuchos::ArrayView<const double> &dist_source_centers,
80 const Teuchos::ArrayView<const SupportId> &dist_source_center_gids,
85 Teuchos::RCP<Tpetra::Operator<double, int, SupportId>> getN()
91 Teuchos::RCP<Tpetra::Operator<double, int, SupportId>> getQ()
98 Teuchos::RCP<Tpetra::CrsMatrix<double, int, SupportId>> d_N;
101 Teuchos::RCP<PolynomialMatrix> d_Q;
110 #endif // end DTK_SPLINEEVALUATIONMATRIX_HPP SplineEvaluationMatrix(const Teuchos::RCP< const Tpetra::Map< int, SupportId >> &domain_map, const Teuchos::RCP< const Tpetra::Map< int, SupportId >> &range_map, const Teuchos::ArrayView< const double > &target_centers, const Teuchos::ArrayView< const SupportId > &target_center_gids, const Teuchos::ArrayView< const double > &dist_source_centers, const Teuchos::ArrayView< const SupportId > &dist_source_center_gids, const SplineInterpolationPairing< DIM > &target_pairings, const Basis &basis)
Constructor.
Local source/parent center pairings.
Sparse spline transformation operator (the A matrix). A = N + Q.
RadialBasisPolicy< Basis > BP
Typedefs.
Policy class for spline interpolation basis functions.
Local child/parent center pairings.