Xpetra  Version of the Day
Xpetra_CrsMatrixUtils.hpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Xpetra: A linear algebra interface package
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact
39 // Jonathan Hu (jhu@sandia.gov)
40 // Ray Tuminaro (rstumin@sandia.gov)
41 // Chris Siefert (csiefer@sandia.gov)
42 // Luc Berger-Vergoat (lberge@sandia.gov)
43 //
44 // ***********************************************************************
45 //
46 // @HEADER
47 #ifndef PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
48 #define PACKAGES_XPETRA_CRSMATRIX_UTILS_HPP_
49 
50 #include "Xpetra_ConfigDefs.hpp"
51 #include "Xpetra_Exceptions.hpp"
52 #include "Xpetra_Map.hpp" // definition of UnderlyingLib
53 
54 #ifdef HAVE_XPETRA_EPETRA
55 #include "Epetra_Util.h"
56 #endif
57 
58 #ifdef HAVE_XPETRA_TPETRA
59 #include "Tpetra_Import_Util2.hpp"
60 #endif
61 
62 namespace Xpetra {
63 
71  template <class Scalar,
72  class LocalOrdinal,
73  class GlobalOrdinal,
76 #undef XPETRA_CRSMATRIXUTILS_SHORT
77 
78  public:
79 
82  static void sortCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
83  const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
84  const Teuchos::ArrayView<Scalar>&CRS_vals,
85  const UnderlyingLib lib) {
86 
87  if(lib == Xpetra::UseEpetra) {
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
91  } else if(lib == Xpetra::UseTpetra) {
92 #ifdef HAVE_XPETRA_TPETRA
93  Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
94 #endif // HAVE_XPETRA_TPETRA
95  }
96 
97  return;
98  }
99 
102  static void sortAndMergeCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
103  const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
104  const Teuchos::ArrayView<Scalar>&CRS_vals,
105  const UnderlyingLib lib) {
106 
107  if(lib == Xpetra::UseEpetra) {
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
111  } else if(lib == Xpetra::UseTpetra) {
112 #ifdef HAVE_XPETRA_TPETRA
113  Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
114 #endif // HAVE_XPETRA_TPETRA
115  }
116 
117  return;
118  }
119 
120  }; // end class CrsMatrixUtils
121 
122 #ifdef HAVE_XPETRA_EPETRA
123 // Specialization for double, int, int, EpetraNode
124  template <>
125  class CrsMatrixUtils<double,int,int,EpetraNode> {
126  typedef double Scalar;
127  typedef int LocalOrdinal;
128  typedef int GlobalOrdinal;
129  typedef EpetraNode Node;
130 #undef XPETRA_CRSMATRIXUTILS_SHORT
131 
132  public:
133 
136  static void sortCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
137  const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
138  const Teuchos::ArrayView<Scalar>&CRS_vals,
139  const UnderlyingLib lib) {
140 
141  if(lib == Xpetra::UseEpetra) {
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());
147  if (rv != 0) {
148  throw Exceptions::RuntimeError("Epetra_Util::SortCrsEntries() return value of "
149  + Teuchos::toString(rv));
150  }
151 #endif // HAVE_XPETRA_EPETRA
152  } else if(lib == Xpetra::UseTpetra) {
153 #ifdef HAVE_XPETRA_TPETRA
154  Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
155 #endif // HAVE_XPETRA_TPETRA
156  }
157 
158  return;
159  }
160 
163  static void sortAndMergeCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
164  const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
165  const Teuchos::ArrayView<Scalar>&CRS_vals,
166  const UnderlyingLib lib) {
167 
168  if(lib == Xpetra::UseEpetra) {
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());
174  if (rv != 0) {
175  throw Exceptions::RuntimeError("Epetra_Util::SortCrsEntries() return value of "
176  + Teuchos::toString(rv));
177  }
178 #endif // HAVE_XPETRA_EPETRA
179  } else if(lib == Xpetra::UseTpetra) {
180 #ifdef HAVE_XPETRA_TPETRA
181  Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
182 #endif // HAVE_XPETRA_TPETRA
183  }
184 
185  return;
186  }
187 
188  }; // end class CrsMatrixUtils
189 
190 
191 // Specialization for double, int, long long, EpetraNode
192  template <>
193  class CrsMatrixUtils<double,int,long long,EpetraNode> {
194  typedef double Scalar;
195  typedef int LocalOrdinal;
196  typedef long long GlobalOrdinal;
197  typedef EpetraNode Node;
198 #undef XPETRA_CRSMATRIXUTILS_SHORT
199 
200  public:
201 
204  static void sortCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
205  const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
206  const Teuchos::ArrayView<Scalar>&CRS_vals,
207  const UnderlyingLib lib) {
208 
209  if(lib == Xpetra::UseEpetra) {
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());
215  if (rv != 0) {
216  throw Exceptions::RuntimeError("Epetra_Util::SortCrsEntries() return value of "
217  + Teuchos::toString(rv));
218  }
219 #endif // HAVE_XPETRA_EPETRA
220  } else if(lib == Xpetra::UseTpetra) {
221 #ifdef HAVE_XPETRA_TPETRA
222  Tpetra::Import_Util::sortCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
223 #endif // HAVE_XPETRA_TPETRA
224  }
225 
226  return;
227  }
228 
231  static void sortAndMergeCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
232  const Teuchos::ArrayView<LocalOrdinal>& CRS_colind,
233  const Teuchos::ArrayView<Scalar>&CRS_vals,
234  const UnderlyingLib lib) {
235 
236  if(lib == Xpetra::UseEpetra) {
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());
242  if (rv != 0) {
243  throw Exceptions::RuntimeError("Epetra_Util::SortCrsEntries() return value of "
244  + Teuchos::toString(rv));
245  }
246 #endif // HAVE_XPETRA_EPETRA
247  } else if(lib == Xpetra::UseTpetra) {
248 #ifdef HAVE_XPETRA_TPETRA
249  Tpetra::Import_Util::sortAndMergeCrsEntries(CRS_rowptr, CRS_colind, CRS_vals);
250 #endif // HAVE_XPETRA_TPETRA
251  }
252 
253  return;
254  }
255 
256  }; // end class CrsMatrixUtils
257 #endif // HAVE_XPETRA_EPETRA for Epetra scpecialization
258 
259 } // end namespace Xpetra
260 
261 #define XPETRA_CRSMATRIXUTILS_SHORT
262 
263 #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)
T * getRawPtr() const
size_type size() 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)
Xpetra namespace