All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Xpetra_CrsMatrixWrap.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_HPP
50 #define XPETRA_CRSMATRIXWRAP_HPP
51 
52 #include <Kokkos_DefaultNode.hpp>
53 
54 #include "Xpetra_ConfigDefs.hpp"
55 #include "Xpetra_Exceptions.hpp"
56 
57 #include "Xpetra_MultiVector.hpp"
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 = Matrix<>::scalar_type,
80  class LocalOrdinal = typename Matrix<Scalar>::local_ordinal_type,
81  class GlobalOrdinal =
82  typename Matrix<Scalar, LocalOrdinal>::global_ordinal_type,
83  class Node =
84  typename Matrix<Scalar, LocalOrdinal, GlobalOrdinal>::node_type>
85 class CrsMatrixWrap :
86  public Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>
87 {
92 #ifdef HAVE_XPETRA_TPETRA
94 #endif
97 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
98 #ifdef HAVE_XPETRA_TPETRA
99  typedef typename CrsMatrix::local_matrix_type local_matrix_type;
100 #endif
101 #endif
102 
103 public:
105 
106 
109  size_t maxNumEntriesPerRow,
111  : finalDefaultView_ (false)
112  {
113  matrixData_ = CrsMatrixFactory::Build (rowMap, maxNumEntriesPerRow, pftype);
115  }
116 
119  const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc,
121  : finalDefaultView_ (false)
122  {
123  matrixData_ = CrsMatrixFactory::Build(rowMap, NumEntriesPerRowToAlloc, pftype);
125  }
126 
128  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map>& colMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype = Xpetra::DynamicProfile)
129  : finalDefaultView_(false)
130  {
131  // Set matrix data
132  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, maxNumEntriesPerRow, pftype);
133 
134  // Default view
136  }
137 
139  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map>& colMap, const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc, Xpetra::ProfileType pftype = Xpetra::DynamicProfile)
140  : finalDefaultView_(false)
141  {
142  // Set matrix data
143  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, NumEntriesPerRowToAlloc, pftype);
144 
145  // Default view
147  }
148 
149 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
150 #ifdef HAVE_XPETRA_TPETRA
151  CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map>& colMap, const local_matrix_type& lclMatrix, const Teuchos::RCP<Teuchos::ParameterList>& params = null)
153  : finalDefaultView_(false)
154  {
155  // Set matrix data
156  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, lclMatrix, params);
157 
158  // Default view
160  }
162  CrsMatrixWrap(const local_matrix_type& lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map>& colMap,
163  const RCP<const Map>& domainMap = Teuchos::null, const RCP<const Map>& rangeMap = Teuchos::null,
164  const Teuchos::RCP<Teuchos::ParameterList>& params = null)
165  : finalDefaultView_(false)
166  {
167  // Set matrix data
168  matrixData_ = CrsMatrixFactory::Build(lclMatrix, rowMap, colMap, domainMap, rangeMap, params);
169 
170  // Default view
172  }
173 #else
174 #ifdef __GNUC__
175 #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."
176 #endif
177 #endif
178 #endif
179 
181  : finalDefaultView_(matrix->isFillComplete())
182  {
183  // Set matrix data
184  matrixData_ = matrix;
185 
186  // Default view
188  }
189 
190  CrsMatrixWrap(const RCP<const CrsGraph>& graph, const RCP<ParameterList>& paramList = Teuchos::null)
191  : finalDefaultView_(false)
192  {
193  // Set matrix data
194  matrixData_ = CrsMatrixFactory::Build(graph, paramList);
195 
196  // Default view
198  }
199 
201  virtual ~CrsMatrixWrap() {}
202 
204 
205 
207 
208 
210 
221  void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView<const GlobalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
222  matrixData_->insertGlobalValues(globalRow, cols, vals);
223  }
224 
226 
233  void insertLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
234  matrixData_->insertLocalValues(localRow, cols, vals);
235  }
236 
238 
243  void replaceGlobalValues(GlobalOrdinal globalRow,
244  const ArrayView<const GlobalOrdinal> &cols,
245  const ArrayView<const Scalar> &vals) { matrixData_->replaceGlobalValues(globalRow, cols, vals); }
246 
248 
251  void replaceLocalValues(LocalOrdinal localRow,
252  const ArrayView<const LocalOrdinal> &cols,
253  const ArrayView<const Scalar> &vals) { matrixData_->replaceLocalValues(localRow, cols, vals); }
254 
256  virtual void setAllToScalar(const Scalar &alpha) { matrixData_->setAllToScalar(alpha); }
257 
259  void scale(const Scalar &alpha) {
260  matrixData_->scale(alpha);
261  }
262 
264 
266 
267 
276  void resumeFill(const RCP< ParameterList > &params=null) {
277  matrixData_->resumeFill(params);
278  }
279 
291  void fillComplete(const RCP<const Map> &domainMap, const RCP<const Map> &rangeMap, const RCP<Teuchos::ParameterList> &params = null) {
292  matrixData_->fillComplete(domainMap, rangeMap, params);
293 
294  // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
296  }
297 
311  //TODO : Get ride of "Tpetra"::OptimizeOption
312  void fillComplete(const RCP<ParameterList> &params = null) {
313  matrixData_->fillComplete(params);
314 
315  // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
317  }
318 
320 
322 
325  return matrixData_->getGlobalNumRows();
326  }
327 
329 
332  return matrixData_->getGlobalNumCols();
333  }
334 
336  size_t getNodeNumRows() const {
337  return matrixData_->getNodeNumRows();
338  }
339 
342  return matrixData_->getGlobalNumEntries();
343  }
344 
346  size_t getNodeNumEntries() const {
347  return matrixData_->getNodeNumEntries();
348  }
349 
351 
352  size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const {
353  return matrixData_->getNumEntriesInLocalRow(localRow);
354  }
355 
357 
359  size_t getGlobalMaxNumRowEntries() const {
360  return matrixData_->getGlobalMaxNumRowEntries();
361  }
362 
364 
366  size_t getNodeMaxNumRowEntries() const {
367  return matrixData_->getNodeMaxNumRowEntries();
368  }
369 
371  bool isLocallyIndexed() const {
372  return matrixData_->isLocallyIndexed();
373  }
374 
376  bool isGloballyIndexed() const {
377  return matrixData_->isGloballyIndexed();
378  }
379 
381  bool isFillComplete() const {
382  return matrixData_->isFillComplete();
383  }
384 
386 
398  void getLocalRowCopy(LocalOrdinal LocalRow,
399  const ArrayView<LocalOrdinal> &Indices,
400  const ArrayView<Scalar> &Values,
401  size_t &NumEntries
402  ) const {
403  matrixData_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
404  }
405 
407 
416  void getGlobalRowView(GlobalOrdinal GlobalRow, ArrayView<const GlobalOrdinal> &indices, ArrayView<const Scalar> &values) const {
417  matrixData_->getGlobalRowView(GlobalRow, indices, values);
418  }
419 
421 
430  void getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &indices, ArrayView<const Scalar> &values) const {
431  matrixData_->getLocalRowView(LocalRow, indices, values);
432  }
433 
435 
438  matrixData_->getLocalDiagCopy(diag);
439  }
440 
443  matrixData_->getLocalDiagOffsets(offsets);
444  }
445 
448  matrixData_->getLocalDiagCopy(diag,offsets);
449  }
450 
453  return matrixData_->getFrobeniusNorm();
454  }
455 
458  matrixData_->leftScale(x);
459  }
460 
463  matrixData_->rightScale(x);
464  }
465 
467  bool haveGlobalConstants() const {
468  return matrixData_->haveGlobalConstants();
469  }
470 
472 
474 
475 
477 
486  //TODO virtual=0 // TODO: Add default parameters ?
487 // void multiply(const MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> & X, MultiVector<Scalar,LocalOrdinal,GlobalOrdinal,Node> &Y, Teuchos::ETransp trans, Scalar alpha, Scalar beta) const {
488 // matrixData_->multiply(X, Y, trans, alpha, beta);
489 // }
490 
492 
494 
495 
497 
503  Scalar alpha = ScalarTraits<Scalar>::one(),
504  Scalar beta = ScalarTraits<Scalar>::zero()) const {
505 
506  matrixData_->apply(X,Y,mode,alpha,beta);
507  }
508 
512  return matrixData_->getDomainMap();
513  }
514 
518  return matrixData_->getRangeMap();
519  }
520 
524 
526  const RCP<const Map> & getColMap(viewLabel_t viewLabel) const {
527  TEUCHOS_TEST_FOR_EXCEPTION(Matrix::operatorViewTable_.containsKey(viewLabel) == false, Xpetra::Exceptions::RuntimeError, "Xpetra::Matrix.GetColMap(): view '" + viewLabel + "' does not exist.");
528  updateDefaultView(); // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated.
529  return Matrix::operatorViewTable_.get(viewLabel)->GetColMap();
530  }
531 
533  matrixData_->removeEmptyProcessesInPlace(newMap);
534  this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetRowMap(matrixData_->getRowMap());
535  this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetColMap(matrixData_->getColMap());
536  }
537 
539 
541  //{@
542 
545  return matrixData_->getMap();
546  }
547 
549  void doImport(const Matrix &source,
551  const CrsMatrixWrap & sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
552  matrixData_->doImport(*sourceWrp.getCrsMatrix(), importer, CM);
553  }
554 
556  void doExport(const Matrix &dest,
558  const CrsMatrixWrap & destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
559  matrixData_->doExport(*destWrp.getCrsMatrix(), importer, CM);
560  }
561 
563  void doImport(const Matrix &source,
565  const CrsMatrixWrap & sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
566  matrixData_->doImport(*sourceWrp.getCrsMatrix(), exporter, CM);
567  }
568 
570  void doExport(const Matrix &dest,
572  const CrsMatrixWrap & destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
573  matrixData_->doExport(*destWrp.getCrsMatrix(), exporter, CM);
574  }
575 
576  // @}
577 
579 
580 
582  std::string description() const {
583  return "Xpetra::CrsMatrixWrap";
584  }
585 
588  // Teuchos::EVerbosityLevel vl = verbLevel;
589  // if (vl == VERB_DEFAULT) vl = VERB_LOW;
590  // RCP<const Comm<int> > comm = this->getComm();
591  // const int myImageID = comm->getRank(),
592  // numImages = comm->getSize();
593 
594  // if (myImageID == 0) out << this->description() << std::endl;
595 
596  matrixData_->describe(out,verbLevel);
597 
598  // Teuchos::OSTab tab(out);
599  }
600 
602 
603  void setObjectLabel( const std::string &objectLabel ) {
604  Teuchos::LabeledObject::setObjectLabel(objectLabel);
605  matrixData_->setObjectLabel(objectLabel);
606  }
608 
609 
610 
611 
612 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
613 #ifdef HAVE_XPETRA_TPETRA
614  local_matrix_type getLocalMatrix () const {
616  return matrixData_->getLocalMatrix();
617  }
618 #else
619 #ifdef __GNUC__
620 #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."
621 #endif
622 #endif
623 #endif
624 
625  // JG: Added:
626 
627  bool hasCrsGraph() const {return true;}
628 
630  RCP<const CrsGraph> getCrsGraph() const { return matrixData_->getCrsGraph(); }
631 
633 
635 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
636  template<class Node2>
637  RCP<Xpetra::Matrix<Scalar,LocalOrdinal,GlobalOrdinal,Node2> > XPETRA_DEPRECATED clone(const RCP<Node2> &node2) const {
638 #ifdef HAVE_XPETRA_TPETRA
641  if (tMatrix == Teuchos::null)
642  throw Xpetra::Exceptions::RuntimeError("clone() functionality is only available for Tpetra");
643 
645  // TODO: inherit strided maps/views ?
646 #else
647  return Teuchos::null;
648 #endif
649  }
650 #endif
651 
652 private:
653 
654  // Default view is created after fillComplete()
655  // Because ColMap might not be available before fillComplete().
657 
658  // Create default view
659  this->defaultViewLabel_ = "point";
660  this->CreateView(this->GetDefaultViewLabel(), matrixData_->getRowMap(), matrixData_->getColMap());
661 
662  // Set current view
663  this->currentViewLabel_ = this->GetDefaultViewLabel();
664  }
665 
666 private:
667 
668  // 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.
669  // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated when getColMap() is called.
670  void updateDefaultView() const {
671  if ((finalDefaultView_ == false) && matrixData_->isFillComplete() ) {
672  // Update default view with the colMap
674  finalDefaultView_ = true;
675  }
676  }
677  // The boolean finalDefaultView_ keep track of the status of the default view (= already updated or not)
678  // See also CrsMatrixWrap::updateDefaultView()
679  mutable bool finalDefaultView_;
680 
681  // The underlying matrix object
683 
684 }; // class CrsMatrixWrap
685 
686 } // namespace Xpetra
687 
688 #define XPETRA_CRSMATRIXWRAP_SHORT
689 #endif //XPETRA_CRSMATRIXWRAP_DECL_HPP
690 
691 //NOTE: if CrsMatrix and VbrMatrix share a common interface for fillComplete() etc, I can move some stuff in Xpetra_Matrix.hpp
692 //TODO: getUnderlyingMatrix() method
void setObjectLabel(const std::string &objectLabel)
RCP< Map< LocalOrdinal, GlobalOrdinal, Node2 > > clone(const Map< LocalOrdinal, GlobalOrdinal, Node1 > &map, const RCP< Node2 > &node2)
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale matrix using the given vector entries.
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
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
virtual ~CrsMatrixWrap()
Destructor.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
CrsMatrixWrap(const RCP< const Map > &rowMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying (possibly different) number of entries in each row.
RCP< CrsMatrix > getCrsMatrix() const
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
void doExport(const Matrix &dest, const Xpetra::Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Export (using an Importer).
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.
global_size_t getGlobalNumRows() const
Returns the number of global rows in this matrix.
const Teuchos::RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype=Xpetra::DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > &plist=Teuchos::null)
Constructor specifying fixed number of entries for each row.
bool hasCrsGraph() const
Supports the getCrsGraph() call.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
CrsMatrixWrap(RCP< CrsMatrix > matrix)
void getLocalDiagCopy(Xpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag, const Teuchos::ArrayView< const size_t > &offsets) const
Get a copy of the diagonal entries owned by this node, with local row indices, using row offsets...
Exception throws to report errors in the internal logical of the program.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Get Frobenius norm of the matrix.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
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.
CrsMatrixWrap(const RCP< const Map > &rowMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying fixed number of entries for each row.
RCP< const CrsGraph > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
void doExport(const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
const RCP< const Map > & getColMap() const
Returns the Map that describes the column distribution in this matrix. This might be null until fillC...
Xpetra::MatrixView< LocalOrdinal, GlobalOrdinal, Node > MatrixView
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...
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void doImport(const Matrix &source, const Xpetra::Export< LocalOrdinal, GlobalOrdinal, Node > &exporter, CombineMode CM)
Import (using an Exporter).
void doImport(const Matrix &source, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
CrsMatrixWrap(const RCP< const Map > &rowMap, const RCP< const Map > &colMap, size_t maxNumEntriesPerRow, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying fixed number of entries for each row and column map.
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.
viewLabel_t defaultViewLabel_
Xpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > CrsGraph
const viewLabel_t & GetCurrentViewLabel() const
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.
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
viewLabel_t currentViewLabel_
size_t getNodeNumEntries() const
Returns the local number of entries in 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.
static const EVerbosityLevel verbLevel_default
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
CrsMatrixWrap(const RCP< const Map > &rowMap, const RCP< const Map > &colMap, const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc, Xpetra::ProfileType pftype=Xpetra::DynamicProfile)
Constructor specifying fixed number of entries for each row and column map.
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.
void resumeFill(const RCP< ParameterList > &params=null)
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
CrsMatrixWrap(const RCP< const CrsGraph > &graph, const RCP< ParameterList > &paramList=Teuchos::null)
Concrete implementation of Xpetra::Matrix.
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.
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...
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > CrsMatrix
const viewLabel_t & GetDefaultViewLabel() const
void fillComplete(const RCP< ParameterList > &params=null)
Signal that data entry is complete.
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
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...
const RCP< const Map > & getColMap(viewLabel_t viewLabel) const
Returns the Map that describes the column distribution in this matrix.