37 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
39#undef XPETRA_CRSMATRIXUTILS_SHORT
45 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
46 const Teuchos::ArrayView<Scalar>& CRS_vals,
47 const UnderlyingLib lib) {
48 if (lib == Xpetra::UseEpetra) {
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)"));
52 }
else if (lib == Xpetra::UseTpetra) {
53#ifdef HAVE_XPETRA_TPETRA
54 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
64 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
65 const Teuchos::ArrayView<Scalar>& CRS_vals,
66 const UnderlyingLib lib) {
67 if (lib == Xpetra::UseEpetra) {
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)"));
71 }
else if (lib == Xpetra::UseTpetra) {
72#ifdef HAVE_XPETRA_TPETRA
73 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
82#ifdef HAVE_XPETRA_EPETRA
85class CrsMatrixUtils<double, int, int,
EpetraNode> {
90#undef XPETRA_CRSMATRIXUTILS_SHORT
95 static void sortCrsEntries(
const Teuchos::ArrayView<size_t>& CRS_rowptr,
96 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
97 const Teuchos::ArrayView<Scalar>& CRS_vals,
98 const UnderlyingLib lib) {
99 if (lib == Xpetra::UseEpetra) {
100#if defined(HAVE_XPETRA_EPETRA)
102 CRS_rowptr.getRawPtr(),
103 CRS_colind.getRawPtr(),
104 CRS_vals.getRawPtr());
106 throw Exceptions::RuntimeError(
"Epetra_Util::SortCrsEntries() return value of " + Teuchos::toString(rv));
109 }
else if (lib == Xpetra::UseTpetra) {
110#ifdef HAVE_XPETRA_TPETRA
111 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
121 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
122 const Teuchos::ArrayView<Scalar>& CRS_vals,
123 const UnderlyingLib lib) {
124 if (lib == Xpetra::UseEpetra) {
125#if defined(HAVE_XPETRA_EPETRA)
127 CRS_rowptr.getRawPtr(),
128 CRS_colind.getRawPtr(),
129 CRS_vals.getRawPtr());
131 throw Exceptions::RuntimeError(
"Epetra_Util::SortCrsEntries() return value of " + Teuchos::toString(rv));
134 }
else if (lib == Xpetra::UseTpetra) {
135#ifdef HAVE_XPETRA_TPETRA
136 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
147class CrsMatrixUtils<double, int, long long,
EpetraNode> {
152#undef XPETRA_CRSMATRIXUTILS_SHORT
157 static void sortCrsEntries(
const Teuchos::ArrayView<size_t>& CRS_rowptr,
158 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
159 const Teuchos::ArrayView<Scalar>& CRS_vals,
160 const UnderlyingLib lib) {
161 if (lib == Xpetra::UseEpetra) {
162#if defined(HAVE_XPETRA_EPETRA)
164 CRS_rowptr.getRawPtr(),
165 CRS_colind.getRawPtr(),
166 CRS_vals.getRawPtr());
168 throw Exceptions::RuntimeError(
"Epetra_Util::SortCrsEntries() return value of " + Teuchos::toString(rv));
171 }
else if (lib == Xpetra::UseTpetra) {
172#ifdef HAVE_XPETRA_TPETRA
173 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
183 const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
184 const Teuchos::ArrayView<Scalar>& CRS_vals,
185 const UnderlyingLib lib) {
186 if (lib == Xpetra::UseEpetra) {
187#if defined(HAVE_XPETRA_EPETRA)
189 CRS_rowptr.getRawPtr(),
190 CRS_colind.getRawPtr(),
191 CRS_vals.getRawPtr());
193 throw Exceptions::RuntimeError(
"Epetra_Util::SortCrsEntries() return value of " + Teuchos::toString(rv));
196 }
else if (lib == Xpetra::UseTpetra) {
197#ifdef HAVE_XPETRA_TPETRA
198 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);