10 #ifndef PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
11 #define PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
15 #include "Xpetra_Map.hpp"
17 #ifdef HAVE_XPETRA_EPETRA
18 #include "Epetra_Util.h"
21 #ifdef HAVE_XPETRA_TPETRA
22 #include "Tpetra_Import_Util2.hpp"
34 template <
class Scalar,
37 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
38 class CrsMatrixUtils {
39 #undef XPETRA_CRSMATRIXUTILS_SHORT
45 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
46 const Teuchos::ArrayView<Scalar>& CRS_vals,
49 #if defined(HAVE_XPETRA_EPETRA)
50 throw(
Xpetra::Exceptions::RuntimeError(
"Xpetra::CrsMatrixUtils::sortCrsEntries only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)"));
51 #endif // HAVE_XPETRA_EPETRA
53 #ifdef HAVE_XPETRA_TPETRA
54 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
55 #endif // HAVE_XPETRA_TPETRA
64 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
65 const Teuchos::ArrayView<Scalar>& CRS_vals,
68 #if defined(HAVE_XPETRA_EPETRA)
69 throw(
Xpetra::Exceptions::RuntimeError(
"Xpetra::CrsMatrixUtils::sortAndMergeCrsEntries only available for GO=int or GO=long long with EpetraNode (Serial or OpenMP depending on configuration)"));
70 #endif // HAVE_XPETRA_EPETRA
72 #ifdef HAVE_XPETRA_TPETRA
73 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
74 #endif // HAVE_XPETRA_TPETRA
82 #ifdef HAVE_XPETRA_EPETRA
90 #undef XPETRA_CRSMATRIXUTILS_SHORT
96 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
97 const Teuchos::ArrayView<Scalar>& CRS_vals,
100 #if defined(HAVE_XPETRA_EPETRA)
101 int rv = Epetra_Util::SortCrsEntries(Teuchos::as<int>(CRS_rowptr.size() - 1),
102 CRS_rowptr.getRawPtr(),
103 CRS_colind.getRawPtr(),
104 CRS_vals.getRawPtr());
108 #endif // HAVE_XPETRA_EPETRA
110 #ifdef HAVE_XPETRA_TPETRA
111 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
112 #endif // HAVE_XPETRA_TPETRA
121 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
122 const Teuchos::ArrayView<Scalar>& CRS_vals,
125 #if defined(HAVE_XPETRA_EPETRA)
126 int rv = Epetra_Util::SortAndMergeCrsEntries(Teuchos::as<int>(CRS_rowptr.size() - 1),
127 CRS_rowptr.getRawPtr(),
128 CRS_colind.getRawPtr(),
129 CRS_vals.getRawPtr());
133 #endif // HAVE_XPETRA_EPETRA
135 #ifdef HAVE_XPETRA_TPETRA
136 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
137 #endif // HAVE_XPETRA_TPETRA
152 #undef XPETRA_CRSMATRIXUTILS_SHORT
158 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
159 const Teuchos::ArrayView<Scalar>& CRS_vals,
162 #if defined(HAVE_XPETRA_EPETRA)
163 int rv = Epetra_Util::SortCrsEntries(Teuchos::as<int>(CRS_rowptr.size() - 1),
164 CRS_rowptr.getRawPtr(),
165 CRS_colind.getRawPtr(),
166 CRS_vals.getRawPtr());
170 #endif // HAVE_XPETRA_EPETRA
172 #ifdef HAVE_XPETRA_TPETRA
173 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
174 #endif // HAVE_XPETRA_TPETRA
183 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
184 const Teuchos::ArrayView<Scalar>& CRS_vals,
187 #if defined(HAVE_XPETRA_EPETRA)
188 int rv = Epetra_Util::SortAndMergeCrsEntries(Teuchos::as<int>(CRS_rowptr.size() - 1),
189 CRS_rowptr.getRawPtr(),
190 CRS_colind.getRawPtr(),
191 CRS_vals.getRawPtr());
195 #endif // HAVE_XPETRA_EPETRA
197 #ifdef HAVE_XPETRA_TPETRA
198 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
199 #endif // HAVE_XPETRA_TPETRA
206 #endif // HAVE_XPETRA_EPETRA for Epetra scpecialization
210 #define XPETRA_CRSMATRIXUTILS_SHORT
212 #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.
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode
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.