Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_CrsMatrixUtils.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
11#define PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
12
13#include "Xpetra_ConfigDefs.hpp"
14#include "Xpetra_Exceptions.hpp"
15#include "Xpetra_Map.hpp" // definition of UnderlyingLib
16
17#ifdef HAVE_XPETRA_EPETRA
18#include "Epetra_Util.h"
19#endif
20
21#ifdef HAVE_XPETRA_TPETRA
22#include "Tpetra_Import_Util2.hpp"
23#endif
24
25namespace Xpetra {
26
34template <class Scalar,
35 class LocalOrdinal,
36 class GlobalOrdinal,
37 class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
39#undef XPETRA_CRSMATRIXUTILS_SHORT
40
41 public:
44 static void sortCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
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)"));
51#endif // HAVE_XPETRA_EPETRA
52 } else if (lib == Xpetra::UseTpetra) {
53#ifdef HAVE_XPETRA_TPETRA
54 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
55#endif // HAVE_XPETRA_TPETRA
56 }
57
58 return;
59 }
60
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)"));
70#endif // HAVE_XPETRA_EPETRA
71 } else if (lib == Xpetra::UseTpetra) {
72#ifdef HAVE_XPETRA_TPETRA
73 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
74#endif // HAVE_XPETRA_TPETRA
75 }
76
77 return;
78 }
79
80}; // end class CrsMatrixUtils
81
82#ifdef HAVE_XPETRA_EPETRA
83// Specialization for double, int, int, EpetraNode
84template <>
85class CrsMatrixUtils<double, int, int, EpetraNode> {
86 typedef double Scalar;
87 typedef int LocalOrdinal;
88 typedef int GlobalOrdinal;
90#undef XPETRA_CRSMATRIXUTILS_SHORT
91
92 public:
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)
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());
105 if (rv != 0) {
106 throw Exceptions::RuntimeError("Epetra_Util::SortCrsEntries() return value of " + Teuchos::toString(rv));
107 }
108#endif // HAVE_XPETRA_EPETRA
109 } else if (lib == Xpetra::UseTpetra) {
110#ifdef HAVE_XPETRA_TPETRA
111 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
112#endif // HAVE_XPETRA_TPETRA
113 }
114
115 return;
116 }
117
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)
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());
130 if (rv != 0) {
131 throw Exceptions::RuntimeError("Epetra_Util::SortCrsEntries() return value of " + Teuchos::toString(rv));
132 }
133#endif // HAVE_XPETRA_EPETRA
134 } else if (lib == Xpetra::UseTpetra) {
135#ifdef HAVE_XPETRA_TPETRA
136 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
137#endif // HAVE_XPETRA_TPETRA
138 }
139
140 return;
141 }
142
143}; // end class CrsMatrixUtils
144
145// Specialization for double, int, long long, EpetraNode
146template <>
147class CrsMatrixUtils<double, int, long long, EpetraNode> {
148 typedef double Scalar;
149 typedef int LocalOrdinal;
150 typedef long long GlobalOrdinal;
152#undef XPETRA_CRSMATRIXUTILS_SHORT
153
154 public:
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)
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());
167 if (rv != 0) {
168 throw Exceptions::RuntimeError("Epetra_Util::SortCrsEntries() return value of " + Teuchos::toString(rv));
169 }
170#endif // HAVE_XPETRA_EPETRA
171 } else if (lib == Xpetra::UseTpetra) {
172#ifdef HAVE_XPETRA_TPETRA
173 Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
174#endif // HAVE_XPETRA_TPETRA
175 }
176
177 return;
178 }
179
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)
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());
192 if (rv != 0) {
193 throw Exceptions::RuntimeError("Epetra_Util::SortCrsEntries() return value of " + Teuchos::toString(rv));
194 }
195#endif // HAVE_XPETRA_EPETRA
196 } else if (lib == Xpetra::UseTpetra) {
197#ifdef HAVE_XPETRA_TPETRA
198 Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
199#endif // HAVE_XPETRA_TPETRA
200 }
201
202 return;
203 }
204
205}; // end class CrsMatrixUtils
206#endif // HAVE_XPETRA_EPETRA for Epetra scpecialization
207
208} // end namespace Xpetra
209
210#define XPETRA_CRSMATRIXUTILS_SHORT
211
212#endif // PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
static int SortAndMergeCrsEntries(int NumRows, int *CRS_rowptr, int *CRS_colind, double *CRS_vals)
static int SortCrsEntries(int NumRows, const int *CRS_rowptr, int *CRS_colind, double *CRS_vals)
size_type size() const
T * getRawPtr() const
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.
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.
Xpetra utility class for CrsMatrix-related routines.
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.
Exception throws to report errors in the internal logical of the program.
std::string toString(const T &t)
Tpetra::KokkosCompat::KokkosSerialWrapperNode EpetraNode