Xpetra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Xpetra_CrsMatrixWrap_decl.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 // Andrey Prokopenko (aprokop@sandia.gov)
41 // Ray Tuminaro (rstumin@sandia.gov)
42 //
43 // ***********************************************************************
44 //
45 // @HEADER
46 
47 // WARNING: This code is experimental. Backwards compatibility should not be expected.
48 
49 #ifndef XPETRA_CRSMATRIXWRAP_DECL_HPP
50 #define XPETRA_CRSMATRIXWRAP_DECL_HPP
51 
52 #include <Tpetra_KokkosCompat_DefaultNode.hpp>
53 
54 #include "Xpetra_ConfigDefs.hpp"
55 #include "Xpetra_Exceptions.hpp"
56 
58 #include "Xpetra_CrsGraph.hpp"
59 #include "Xpetra_CrsMatrix.hpp"
61 
62 #include "Xpetra_Matrix.hpp"
63 
64 #include <Teuchos_SerialDenseMatrix.hpp>
65 #include <Teuchos_Hashtable.hpp>
66 
71 namespace Xpetra {
72 
73 typedef std::string viewLabel_t;
74 
79 template <class Scalar,
80  class LocalOrdinal,
81  class GlobalOrdinal,
82  class Node = Tpetra::KokkosClassic::DefaultNode::DefaultNodeType>
83 class CrsMatrixWrap : public Matrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> {
88 #ifdef HAVE_XPETRA_TPETRA
90 #endif
93 #ifdef HAVE_XPETRA_TPETRA
95 #endif
96 
97  public:
99 
100 
102  CrsMatrixWrap(const RCP<const Map> &rowMap);
103 
105  CrsMatrixWrap(const RCP<const Map> &rowMap,
106  size_t maxNumEntriesPerRow);
107 
109  CrsMatrixWrap(const RCP<const Map> &rowMap,
110  const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
111 
113  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, size_t maxNumEntriesPerRow);
114 
116  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
117 
118 #ifdef HAVE_XPETRA_TPETRA
119  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map> &colMap, const local_matrix_type &lclMatrix, const Teuchos::RCP<Teuchos::ParameterList> &params = null);
121 
123  CrsMatrixWrap(const local_matrix_type &lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map> &colMap,
124  const RCP<const Map> &domainMap = Teuchos::null, const RCP<const Map> &rangeMap = Teuchos::null,
125  const Teuchos::RCP<Teuchos::ParameterList> &params = null);
126 #else
127 #ifdef __GNUC__
128 #warning "Xpetra Kokkos interface for CrsMatrix is enabled (HAVE_XPETRA_KOKKOS_REFACTOR) but Tpetra is disabled. The Kokkos interface needs Tpetra to be enabled, too."
129 #endif
130 #endif
131 
132  CrsMatrixWrap(RCP<CrsMatrix> matrix);
133 
134  CrsMatrixWrap(const RCP<const CrsGraph> &graph, const RCP<ParameterList> &paramList = Teuchos::null);
135 
136  CrsMatrixWrap(const RCP<const CrsGraph> &graph,
138  const RCP<ParameterList> &paramList = Teuchos::null);
139 
141  virtual ~CrsMatrixWrap();
142 
144 
146 
147 
149 
160  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals);
161 
163 
170  void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals);
171 
173 
178  void replaceGlobalValues(GlobalOrdinal globalRow,
179  const ArrayView<const GlobalOrdinal> &cols,
180  const ArrayView<const Scalar> &vals);
181 
183 
186  void replaceLocalValues(LocalOrdinal localRow,
187  const ArrayView<const LocalOrdinal> &cols,
188  const ArrayView<const Scalar> &vals);
189 
191  virtual void setAllToScalar(const Scalar &alpha);
192 
194  void scale(const Scalar &alpha);
195 
197 
199 
200 
209  void resumeFill(const RCP<ParameterList> &params = null);
210 
222  void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params = null);
223 
237  // TODO : Get ride of "Tpetra"::OptimizeOption
238  void fillComplete(const RCP<ParameterList> &params = null);
239 
241 
243 
246 
248 
251 
253  size_t getLocalNumRows() const;
254 
257 
259  size_t getLocalNumEntries() const;
260 
262 
263  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
264 
266 
267  size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
268 
270 
272  size_t getGlobalMaxNumRowEntries() const;
273 
275 
277  size_t getLocalMaxNumRowEntries() const;
278 
280  bool isLocallyIndexed() const;
281 
283  bool isGloballyIndexed() const;
284 
286  bool isFillComplete() const;
287 
289 
301  void getLocalRowCopy(LocalOrdinal LocalRow,
302  const ArrayView<LocalOrdinal> &Indices,
303  const ArrayView<Scalar> &Values,
304  size_t &NumEntries) const;
305 
307 
316  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const;
317 
319 
328  void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const;
329 
331 
334 
336  void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const;
337 
339  void getLocalDiagCopy(Xpetra::Vector<Scalar, LocalOrdinal, GlobalOrdinal, Node> &diag, const Teuchos::ArrayView<const size_t> &offsets) const;
340 
342  typename ScalarTraits<Scalar>::magnitudeType getFrobeniusNorm() const;
343 
346 
349 
351  bool haveGlobalConstants() const;
352 
354 
356 
357 
359 
368  // TODO virtual=0 // TODO: Add default parameters ?
369  // void multiply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & X, MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y, Teuchos::ETransp trans, Scalar alpha, Scalar beta) const {
370  // matrixData_->multiply(X, Y, trans, alpha, beta);
371  // }
372 
374 
376 
377 
379 
384  Teuchos::ETransp mode = Teuchos::NO_TRANS,
385  Scalar alpha = ScalarTraits<Scalar>::one(),
386  Scalar beta = ScalarTraits<Scalar>::zero()) const;
387 
391  Teuchos::ETransp mode,
392  Scalar alpha,
393  Scalar beta,
394  bool sumInterfaceValues,
395  const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> > &regionInterfaceImporter,
396  const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const;
397 
400  const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getDomainMap() const;
401 
404  const RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getRangeMap() const;
405 
408  const RCP<const Map> &getColMap() const;
409 
411  const RCP<const Map> &getColMap(viewLabel_t viewLabel) const;
412 
413  void removeEmptyProcessesInPlace(const Teuchos::RCP<const Map> &newMap);
414 
416 
418  //{@
419 
421  const Teuchos::RCP<const Xpetra::Map<LocalOrdinal, GlobalOrdinal, Node> > getMap() const;
422 
424  void doImport(const Matrix &source,
426 
428  void doExport(const Matrix &dest,
430 
432  void doImport(const Matrix &source,
434 
436  void doExport(const Matrix &dest,
438 
439  // @}
440 
442 
443 
445  std::string description() const;
446 
448  void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel = Teuchos::Describable::verbLevel_default) const;
449 
451 
452  void setObjectLabel(const std::string &objectLabel);
454 
455 #ifdef HAVE_XPETRA_TPETRA
456  virtual local_matrix_type getLocalMatrixDevice() const;
457  virtual typename local_matrix_type::HostMirror getLocalMatrixHost() const;
458 #else
459 #ifdef __GNUC__
460 #warning "Xpetra Kokkos interface for CrsMatrix is enabled (HAVE_XPETRA_KOKKOS_REFACTOR) but Tpetra is disabled. The Kokkos interface needs Tpetra to be enabled, too."
461 #endif
462 #endif
463 
464  // JG: Added:
465 
466  bool hasCrsGraph() const;
467 
469  RCP<const CrsGraph> getCrsGraph() const;
470 
471  RCP<CrsMatrix> getCrsMatrix() const;
472 
474  LocalOrdinal GetStorageBlockSize() const;
475 
480 
482  void replaceCrsMatrix(RCP<CrsMatrix> &M);
483 
485  private:
486  // Default view is created after fillComplete()
487  // Because ColMap might not be available before fillComplete().
488  void CreateDefaultView();
489 
490  private:
491  // The colMap can be <tt>null</tt> until fillComplete() is called. The default view of the Matrix have to be updated when fillComplete() is called.
492  // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated when getColMap() is called.
493  void updateDefaultView() const;
494  // The boolean finalDefaultView_ keep track of the status of the default view (= already updated or not)
495  // See also CrsMatrixWrap::updateDefaultView()
496  mutable bool finalDefaultView_;
497 
498  // The underlying matrix object
499  RCP<CrsMatrix> matrixData_;
500 
501 }; // class CrsMatrixWrap
502 
503 } // namespace Xpetra
504 
505 #define XPETRA_CRSMATRIXWRAP_SHORT
506 #endif // XPETRA_CRSMATRIXWRAP_DECL_HPP
507 
508 // NOTE: if CrsMatrix and VbrMatrix share a common interface for fillComplete() etc, I can move some stuff in Xpetra_Matrix.hpp
509 // TODO: getUnderlyingMatrix() method
void setObjectLabel(const std::string &objectLabel)
RCP< CrsMatrix > getCrsMatrix() const
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale matrix using the given vector entries.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Xpetra::TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrix
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
virtual local_matrix_type getLocalMatrixDevice() const
global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
bool hasCrsGraph() const
Supports the getCrsGraph() call.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
LocalOrdinal GetStorageBlockSize() const
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
Xpetra::MatrixView< Scalar, LocalOrdinal, GlobalOrdinal, Node > MatrixView
void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView< const GlobalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of global indices in a specified row of the matrix.
void doExport(const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
void doImport(const Matrix &source, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
std::string description() const
Return a simple one-line description of this object.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > CrsGraph
CrsMatrix::local_matrix_type local_matrix_type
void getLocalDiagCopy(Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
CrsMatrixWrap(const RCP< const Map > &rowMap)
Constructor for a dynamic profile matrix (Epetra only)
void residual(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &B, MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &R) const
Compute a residual R = B - (*this) * X.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
std::string viewLabel_t
size_t global_size_t
Global size_t object.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void fillComplete(const RCP< const Map > &domainMap, const RCP< const Map > &rangeMap, const RCP< Teuchos::ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
virtual ~CrsMatrixWrap()
Destructor.
const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
void resumeFill(const RCP< ParameterList > &params=null)
void replaceCrsMatrix(RCP< CrsMatrix > &M)
Expert only.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
const RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
CombineMode
Xpetra::Combine Mode enumerable type.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale matrix using the given vector entries.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
virtual void apply(const Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, Xpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
KokkosSparse::CrsMatrix< impl_scalar_type, LocalOrdinal, execution_space, void, typename local_graph_type::size_type > local_matrix_type
The specialization of Kokkos::CrsMatrix that represents the part of the sparse matrix on each MPI pro...
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified global row.
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
void getLocalRowView(LocalOrdinal LocalRow, ArrayView< const LocalOrdinal > &indices, ArrayView< const Scalar > &values) const
Extract a const, non-persisting view of local indices in a specified row of the matrix.
Xpetra::CrsMatrixFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixFactory
virtual local_matrix_type::HostMirror getLocalMatrixHost() const
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
Xpetra-specific matrix class.
void getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView< LocalOrdinal > &Indices, const ArrayView< Scalar > &Values, size_t &NumEntries) const
Extract a list of entries in a specified local row of the matrix. Put into storage allocated by calli...