41 #ifndef DTK_SPLINECOEFFICIENTMATRIX_IMPL_HPP 42 #define DTK_SPLINECOEFFICIENTMATRIX_IMPL_HPP 48 #include <Teuchos_Array.hpp> 56 template <
class Basis,
int DIM>
58 const Teuchos::RCP<
const Tpetra::Map<int, SupportId>> &operator_map,
59 const Teuchos::ArrayView<const double> &source_centers,
60 const Teuchos::ArrayView<const SupportId> &source_center_gids,
61 const Teuchos::ArrayView<const double> &dist_source_centers,
62 const Teuchos::ArrayView<const SupportId> &dist_source_center_gids,
65 DTK_CHECK( 0 == source_centers.size() % DIM );
66 DTK_CHECK( source_centers.size() / DIM == source_center_gids.size() );
67 DTK_CHECK( 0 == dist_source_centers.size() % DIM );
68 DTK_CHECK( dist_source_centers.size() / DIM ==
69 dist_source_center_gids.size() );
72 unsigned num_source_centers = source_center_gids.size();
76 Teuchos::RCP<Tpetra::MultiVector<double, int, SupportId>> P_vec =
77 Tpetra::createMultiVector<double, int, SupportId>( operator_map,
80 for (
unsigned i = 0; i < num_source_centers; ++i )
82 P_vec->replaceGlobalValue( source_center_gids[i], 0, 1.0 );
84 for (
int d = 0; d < DIM; ++d )
86 P_vec->replaceGlobalValue( source_center_gids[i], d + 1,
87 source_centers[di + d] );
94 Teuchos::ArrayRCP<SupportId> children_per_parent =
95 source_pairings.childrenPerParent();
96 SupportId max_entries_per_row = *std::max_element(
97 children_per_parent.begin(), children_per_parent.end() );
98 d_M = Teuchos::rcp(
new Tpetra::CrsMatrix<double, int, SupportId>(
99 operator_map, max_entries_per_row ) );
100 Teuchos::Array<SupportId> M_indices( max_entries_per_row );
101 Teuchos::Array<double> values( max_entries_per_row );
103 Teuchos::ArrayView<const unsigned> source_neighbors;
107 for (
unsigned i = 0; i < num_source_centers; ++i )
112 nsn = source_neighbors.size();
113 radius = source_pairings.parentSupportRadius( i );
116 for (
int j = 0; j < nsn; ++j )
118 dj = DIM * source_neighbors[j];
119 M_indices[j] = dist_source_center_gids[source_neighbors[j]];
122 &dist_source_centers[dj] );
124 values[j] = BP::evaluateValue( basis, radius, dist );
126 d_M->insertGlobalValues( source_center_gids[i], M_indices( 0, nsn ),
131 DTK_ENSURE( d_M->isFillComplete() );
140 #endif // end DTK_SPLINECOEFFICIENTMATRIX_IMPL_HPP Spline interpolation coefficient matrix.
SplineCoefficientMatrix(const Teuchos::RCP< const Tpetra::Map< int, SupportId >> &operator_map, const Teuchos::ArrayView< const double > &source_centers, const Teuchos::ArrayView< const SupportId > &source_center_gids, const Teuchos::ArrayView< const double > &dist_source_centers, const Teuchos::ArrayView< const SupportId > &dist_source_center_gids, const SplineInterpolationPairing< DIM > &source_pairings, const Basis &basis)
Constructor.
unsigned long int SupportId
Support id type.
Assertions and Design-by-Contract for error handling.
Eucliden distance function.
Vector apply implementation for polynomial matrices.
Teuchos::ArrayView< const unsigned > childCenterIds(const unsigned parent_id) const
Given a parent center local id get the ids of the child centers within the given radius.
Euclidean distance function.
Local child/parent center pairings.