47 #ifndef PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
48 #define PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
52 #include "Xpetra_Map.hpp"
54 #ifdef HAVE_XPETRA_EPETRA
55 #include "Epetra_Util.h"
58 #ifdef HAVE_XPETRA_TPETRA
59 #include "Tpetra_Import_Util2.hpp"
71 template <
class Scalar,
75 class CrsMatrixUtils {
76 #undef XPETRA_CRSMATRIXUTILS_SHORT
83 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
84 const Teuchos::ArrayView<Scalar>&CRS_vals,
88 #if defined(HAVE_XPETRA_EPETRA)
89 throw(
Xpetra::Exceptions::RuntimeError(
"Xpetra::CrsMatrixUtils::sortCrsEntries only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)"));
90 #endif // HAVE_XPETRA_EPETRA
92 #ifdef HAVE_XPETRA_TPETRA
93 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
94 #endif // HAVE_XPETRA_TPETRA
103 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
104 const Teuchos::ArrayView<Scalar>&CRS_vals,
108 #if defined(HAVE_XPETRA_EPETRA)
109 throw(
Xpetra::Exceptions::RuntimeError(
"Xpetra::CrsMatrixUtils::sortAndMergeCrsEntries only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)"));
110 #endif // HAVE_XPETRA_EPETRA
112 #ifdef HAVE_XPETRA_TPETRA
113 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
114 #endif // HAVE_XPETRA_TPETRA
122 #ifdef HAVE_XPETRA_EPETRA
130 #undef XPETRA_CRSMATRIXUTILS_SHORT
137 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
138 const Teuchos::ArrayView<Scalar>&CRS_vals,
142 #if defined(HAVE_XPETRA_EPETRA)
143 int rv = Epetra_Util::SortCrsEntries(Teuchos::as<int>(CRS_rowptr.size()-1),
144 CRS_rowptr.getRawPtr(),
145 CRS_colind.getRawPtr(),
146 CRS_vals.getRawPtr());
151 #endif // HAVE_XPETRA_EPETRA
153 #ifdef HAVE_XPETRA_TPETRA
154 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
155 #endif // HAVE_XPETRA_TPETRA
164 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
165 const Teuchos::ArrayView<Scalar>&CRS_vals,
169 #if defined(HAVE_XPETRA_EPETRA)
170 int rv = Epetra_Util::SortAndMergeCrsEntries(Teuchos::as<int>(CRS_rowptr.size()-1),
171 CRS_rowptr.getRawPtr(),
172 CRS_colind.getRawPtr(),
173 CRS_vals.getRawPtr());
178 #endif // HAVE_XPETRA_EPETRA
180 #ifdef HAVE_XPETRA_TPETRA
181 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
182 #endif // HAVE_XPETRA_TPETRA
198 #undef XPETRA_CRSMATRIXUTILS_SHORT
205 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
206 const Teuchos::ArrayView<Scalar>&CRS_vals,
210 #if defined(HAVE_XPETRA_EPETRA)
211 int rv = Epetra_Util::SortCrsEntries(Teuchos::as<int>(CRS_rowptr.size()-1),
212 CRS_rowptr.getRawPtr(),
213 CRS_colind.getRawPtr(),
214 CRS_vals.getRawPtr());
219 #endif // HAVE_XPETRA_EPETRA
221 #ifdef HAVE_XPETRA_TPETRA
222 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
223 #endif // HAVE_XPETRA_TPETRA
232 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
233 const Teuchos::ArrayView<Scalar>&CRS_vals,
237 #if defined(HAVE_XPETRA_EPETRA)
238 int rv = Epetra_Util::SortAndMergeCrsEntries(Teuchos::as<int>(CRS_rowptr.size()-1),
239 CRS_rowptr.getRawPtr(),
240 CRS_colind.getRawPtr(),
241 CRS_vals.getRawPtr());
246 #endif // HAVE_XPETRA_EPETRA
248 #ifdef HAVE_XPETRA_TPETRA
249 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
250 #endif // HAVE_XPETRA_TPETRA
257 #endif // HAVE_XPETRA_EPETRA for Epetra scpecialization
261 #define XPETRA_CRSMATRIXUTILS_SHORT
263 #endif // PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
static void sortAndMergeCrsEntries(const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< LocalOrdinal > &CRS_colind, const Teuchos::ArrayView< Scalar > &CRS_vals, const UnderlyingLib lib)
Sort and merge the entries of the (raw CSR) matrix by column index within each row.
std::string toString(Xpetra::UnderlyingLib lib)
Convert a Xpetra::UnderlyingLib to a std::string.
Xpetra utility class for CrsMatrix-related routines.
Exception throws to report errors in the internal logical of the program.
static void sortCrsEntries(const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< LocalOrdinal > &CRS_colind, const Teuchos::ArrayView< Scalar > &CRS_vals, const UnderlyingLib lib)
Sort the entries of the (raw CSR) matrix by column index within each row.
static void sortCrsEntries(const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< LocalOrdinal > &CRS_colind, const Teuchos::ArrayView< Scalar > &CRS_vals, const UnderlyingLib lib)
Sort the entries of the (raw CSR) matrix by column index within each row.
static void sortAndMergeCrsEntries(const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< LocalOrdinal > &CRS_colind, const Teuchos::ArrayView< Scalar > &CRS_vals, const UnderlyingLib lib)
Sort and merge the entries of the (raw CSR) matrix by column index within each row.
static void sortCrsEntries(const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< LocalOrdinal > &CRS_colind, const Teuchos::ArrayView< Scalar > &CRS_vals, const UnderlyingLib lib)
Sort the entries of the (raw CSR) matrix by column index within each row.
static void sortAndMergeCrsEntries(const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< LocalOrdinal > &CRS_colind, const Teuchos::ArrayView< Scalar > &CRS_vals, const UnderlyingLib lib)
Sort and merge the entries of the (raw CSR) matrix by column index within each row.