Xpetra  Version of the Day
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 <Kokkos_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 
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,
84  public Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
85 {
90 #ifdef HAVE_XPETRA_TPETRA
92 #endif
95 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
96 #ifdef HAVE_XPETRA_TPETRA
97  typedef typename CrsMatrix::local_matrix_type local_matrix_type;
98 #endif
99 #endif
100 
101 public:
103 
104 
106  CrsMatrixWrap (const RCP<const Map>& rowMap);
107 
109  CrsMatrixWrap (const RCP<const Map>& rowMap,
110  size_t maxNumEntriesPerRow);
111 
113  CrsMatrixWrap (const RCP<const Map>& rowMap,
114  const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc);
115 
117  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map>& colMap, size_t maxNumEntriesPerRow);
118 
120  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map>& colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc);
121 
122 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
123 #ifdef HAVE_XPETRA_TPETRA
125  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map>& colMap, const local_matrix_type& lclMatrix, const Teuchos::RCP<Teuchos::ParameterList>& params = null);
126 
128  CrsMatrixWrap(const local_matrix_type& lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map>& colMap,
129  const RCP<const Map>& domainMap = Teuchos::null, const RCP<const Map>& rangeMap = Teuchos::null,
130  const Teuchos::RCP<Teuchos::ParameterList>& params = null);
131 #else
132 #ifdef __GNUC__
133 #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."
134 #endif
135 #endif
136 #endif
137 
139 
140  CrsMatrixWrap(const RCP<const CrsGraph>& graph, const RCP<ParameterList>& paramList = Teuchos::null);
141 
143  virtual ~CrsMatrixWrap();
144 
146 
147 
149 
150 
152 
163  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals);
164 
166 
173  void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals);
174 
176 
181  void replaceGlobalValues(GlobalOrdinal globalRow,
182  const ArrayView<const GlobalOrdinal> &cols,
183  const ArrayView<const Scalar> &vals);
184 
186 
189  void replaceLocalValues(LocalOrdinal localRow,
190  const ArrayView<const LocalOrdinal> &cols,
191  const ArrayView<const Scalar> &vals);
192 
194  virtual void setAllToScalar(const Scalar &alpha);
195 
197  void scale(const Scalar &alpha);
198 
200 
202 
203 
212  void resumeFill(const RCP< ParameterList > &params=null);
213 
225  void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params = null);
226 
240  //TODO : Get ride of "Tpetra"::OptimizeOption
241  void fillComplete(const RCP<ParameterList> &params = null);
242 
244 
246 
249 
251 
254 
256  size_t getNodeNumRows() const;
257 
260 
262  size_t getNodeNumEntries() const;
263 
265 
266  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const;
267 
269 
270  size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const;
271 
273 
275  size_t getGlobalMaxNumRowEntries() const;
276 
278 
280  size_t getNodeMaxNumRowEntries() const;
281 
283  bool isLocallyIndexed() const;
284 
286  bool isGloballyIndexed() const;
287 
289  bool isFillComplete() const;
290 
292 
304  void getLocalRowCopy(LocalOrdinal LocalRow,
305  const ArrayView<LocalOrdinal> &Indices,
306  const ArrayView<Scalar> &Values,
307  size_t &NumEntries
308  ) const;
309 
311 
320  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const;
321 
323 
332  void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const;
333 
335 
338 
340  void getLocalDiagOffsets(Teuchos::ArrayRCP<size_t> &offsets) const;
341 
344 
347 
350 
353 
355  bool haveGlobalConstants() const;
356 
358 
360 
361 
363 
372  //TODO virtual=0 // TODO: Add default parameters ?
373 // void multiply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & X, MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y, Teuchos::ETransp trans, Scalar alpha, Scalar beta) const {
374 // matrixData_->multiply(X, Y, trans, alpha, beta);
375 // }
376 
378 
380 
381 
383 
389  Scalar alpha = ScalarTraits<Scalar>::one(),
390  Scalar beta = ScalarTraits<Scalar>::zero()) const;
391 
395  Teuchos::ETransp mode,
396  Scalar alpha,
397  Scalar beta,
398  bool sumInterfaceValues,
399  const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> >& regionInterfaceImporter,
400  const Teuchos::ArrayRCP<LocalOrdinal>& regionInterfaceLIDs) const;
401 
405 
409 
412  const RCP<const Map> & getColMap() const;
413 
415  const RCP<const Map> & getColMap(viewLabel_t viewLabel) const;
416 
418 
420 
422  //{@
423 
426 
428  void doImport(const Matrix &source,
430 
432  void doExport(const Matrix &dest,
434 
436  void doImport(const Matrix &source,
438 
440  void doExport(const Matrix &dest,
442 
443  // @}
444 
446 
447 
449  std::string description() const;
450 
453 
455 
456  void setObjectLabel( const std::string &objectLabel );
458 
459 
460 
461 
462 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
463 #ifdef HAVE_XPETRA_TPETRA
464 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
465  virtual local_matrix_type getLocalMatrix () const;
466 #endif
467 
468  virtual local_matrix_type getLocalMatrixDevice () const;
469  virtual typename local_matrix_type::HostMirror getLocalMatrixHost () const;
470 #else
471 #ifdef __GNUC__
472 #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."
473 #endif
474 #endif
475 #endif
476 
477  // JG: Added:
478 
479  bool hasCrsGraph() const;
480 
483 
485 
486 
491 
492 
494 private:
495 
496  // Default view is created after fillComplete()
497  // Because ColMap might not be available before fillComplete().
498  void CreateDefaultView();
499 
500 private:
501 
502  // 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.
503  // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated when getColMap() is called.
504  void updateDefaultView() const;
505  // The boolean finalDefaultView_ keep track of the status of the default view (= already updated or not)
506  // See also CrsMatrixWrap::updateDefaultView()
507  mutable bool finalDefaultView_;
508 
509  // The underlying matrix object
511 
512 }; // class CrsMatrixWrap
513 
514 } // namespace Xpetra
515 
516 #define XPETRA_CRSMATRIXWRAP_SHORT
517 #endif //XPETRA_CRSMATRIXWRAP_DECL_HPP
518 
519 //NOTE: if CrsMatrix and VbrMatrix share a common interface for fillComplete() etc, I can move some stuff in Xpetra_Matrix.hpp
520 //TODO: getUnderlyingMatrix() method
static const EVerbosityLevel verbLevel_default
Concrete implementation of Xpetra::Matrix.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale matrix using the given vector entries.
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.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified global row.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
Xpetra::Matrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > Matrix
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.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale matrix using the given vector entries.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
std::string description() const
Return a simple one-line description of this object.
Xpetra::CrsMatrixFactory< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrixFactory
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
Xpetra::TpetraCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > TpetraCrsMatrix
Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > Map
bool hasCrsGraph() const
Supports the getCrsGraph() call.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
void resumeFill(const RCP< ParameterList > &params=null)
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
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...
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
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.
global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
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...
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.
void setObjectLabel(const std::string &objectLabel)
const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
Xpetra::MatrixView< Scalar, LocalOrdinal, GlobalOrdinal, Node > MatrixView
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void doExport(const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
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.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
virtual ~CrsMatrixWrap()
Destructor.
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 insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
CrsMatrixWrap(const RCP< const Map > &rowMap)
Constructor for a dynamic profile matrix (Epetra only)
void doImport(const Matrix &source, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
RCP< CrsMatrix > getCrsMatrix() const
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.
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode.
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > CrsGraph
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
Xpetra-specific matrix class.
Xpetra namespace
std::string viewLabel_t
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.