Xpetra  Version of the Day
Xpetra_CrsMatrixWrap_def.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_DEF_HPP
50 #define XPETRA_CRSMATRIXWRAP_DEF_HPP
51 
52 #include <Kokkos_DefaultNode.hpp>
53 
55 #include <Teuchos_Hashtable.hpp>
56 
57 #include "Xpetra_ConfigDefs.hpp"
58 #include "Xpetra_Exceptions.hpp"
59 
60 #include "Xpetra_MultiVector.hpp"
61 #include "Xpetra_CrsGraph.hpp"
62 #include "Xpetra_CrsMatrix.hpp"
64 
65 #include "Xpetra_Matrix.hpp"
66 
68 
69 namespace Xpetra {
70  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
72  : finalDefaultView_ (false)
73  {
76  }
77 
78  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
80  size_t maxNumEntriesPerRow)
81  : finalDefaultView_ (false)
82  {
83  matrixData_ = CrsMatrixFactory::Build (rowMap, maxNumEntriesPerRow);
85  }
86 
87  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
89  const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc)
90  : finalDefaultView_ (false)
91  {
92  matrixData_ = CrsMatrixFactory::Build(rowMap, NumEntriesPerRowToAlloc);
94  }
95 
96  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
98  : finalDefaultView_(false)
99  {
100  // Set matrix data
101  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, maxNumEntriesPerRow);
102 
103  // Default view
105  }
106 
107  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
109  : finalDefaultView_(false)
110  {
111  // Set matrix data
112  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, NumEntriesPerRowToAlloc);
113 
114  // Default view
116  }
117 
118 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
119 #ifdef HAVE_XPETRA_TPETRA
120  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
121  CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>::CrsMatrixWrap(const RCP<const Map> &rowMap, const RCP<const Map>& colMap, const local_matrix_type& lclMatrix, const Teuchos::RCP<Teuchos::ParameterList>& params)
122  : finalDefaultView_(false)
123  {
124  // Set matrix data
125  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, lclMatrix, params);
126 
127  // Default view
129  }
130 
131  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
132  CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>::CrsMatrixWrap(const local_matrix_type& lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map>& colMap,
133  const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap,
135  : finalDefaultView_(false)
136  {
137  // Set matrix data
138  matrixData_ = CrsMatrixFactory::Build(lclMatrix, rowMap, colMap, domainMap, rangeMap, params);
139 
140  // Default view
142  }
143 #else
144 #ifdef __GNUC__
145 #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."
146 #endif
147 #endif
148 #endif
149 
150  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
152  : finalDefaultView_(matrix->isFillComplete())
153  {
154  // Set matrix data
155  matrixData_ = matrix;
156 
157  // Default view
159  }
160 
161  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
163  : finalDefaultView_(false)
164  {
165  // Set matrix data
166  matrixData_ = CrsMatrixFactory::Build(graph, paramList);
167 
168  // Default view
170  }
171 
172  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
174 
175  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
177  matrixData_->insertGlobalValues(globalRow, cols, vals);
178  }
179 
180  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
182  matrixData_->insertLocalValues(localRow, cols, vals);
183  }
184 
185  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
187  const ArrayView<const GlobalOrdinal> &cols,
188  const ArrayView<const Scalar> &vals) { matrixData_->replaceGlobalValues(globalRow, cols, vals); }
189 
190  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
192  const ArrayView<const LocalOrdinal> &cols,
193  const ArrayView<const Scalar> &vals) { matrixData_->replaceLocalValues(localRow, cols, vals); }
194 
195  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
196  void CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setAllToScalar(const Scalar &alpha) { matrixData_->setAllToScalar(alpha); }
197 
198  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
200  matrixData_->scale(alpha);
201  }
202 
203  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
205  matrixData_->resumeFill(params);
206  }
207 
208  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
210  matrixData_->fillComplete(domainMap, rangeMap, params);
211 
212  // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
213  updateDefaultView();
214  }
215 
216  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
218  matrixData_->fillComplete(params);
219 
220  // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
221  updateDefaultView();
222  }
223 
224  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
226  return matrixData_->getGlobalNumRows();
227  }
228 
229  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
231  return matrixData_->getGlobalNumCols();
232  }
233 
234  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
236  return matrixData_->getNodeNumRows();
237  }
238 
239  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
241  return matrixData_->getGlobalNumEntries();
242  }
243 
244  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246  return matrixData_->getNodeNumEntries();
247  }
248 
249  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
251  return matrixData_->getNumEntriesInLocalRow(localRow);
252  }
253 
254  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
256  return matrixData_->getNumEntriesInGlobalRow(globalRow);
257  }
258 
259  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
261  return matrixData_->getGlobalMaxNumRowEntries();
262  }
263 
264  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
266  return matrixData_->getNodeMaxNumRowEntries();
267  }
268 
269  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
271  return matrixData_->isLocallyIndexed();
272  }
273 
274  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
276  return matrixData_->isGloballyIndexed();
277  }
278 
279  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281  return matrixData_->isFillComplete();
282  }
283 
284  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
286  const ArrayView<LocalOrdinal> &Indices,
287  const ArrayView<Scalar> &Values,
288  size_t &NumEntries
289  ) const {
290  matrixData_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
291  }
292 
293  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
295  matrixData_->getGlobalRowView(GlobalRow, indices, values);
296  }
297 
298  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300  matrixData_->getLocalRowView(LocalRow, indices, values);
301  }
302 
303  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
305  matrixData_->getLocalDiagCopy(diag);
306  }
307 
308  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
310  matrixData_->getLocalDiagOffsets(offsets);
311  }
312 
313  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
315  matrixData_->getLocalDiagCopy(diag,offsets);
316  }
317 
318  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
320  return matrixData_->getFrobeniusNorm();
321  }
322 
323  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
325  matrixData_->leftScale(x);
326  }
327 
328  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
330  matrixData_->rightScale(x);
331  }
332 
333  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
335  return matrixData_->haveGlobalConstants();
336  }
337 
338  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
341  Teuchos::ETransp mode,
342  Scalar alpha,
343  Scalar beta) const {
344 
345  matrixData_->apply(X,Y,mode,alpha,beta);
346  }
347 
348  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
351  Teuchos::ETransp mode,
352  Scalar alpha,
353  Scalar beta,
354  bool sumInterfaceValues,
355  const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> >& regionInterfaceImporter,
356  const Teuchos::ArrayRCP<LocalOrdinal>& regionInterfaceLIDs
357  ) const{
358  matrixData_->apply(X,Y,mode,alpha,beta,sumInterfaceValues,regionInterfaceImporter,regionInterfaceLIDs);
359  }
360 
361  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
363  return matrixData_->getDomainMap();
364  }
365 
366  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
368  return matrixData_->getRangeMap();
369  }
370 
371  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
373 
374  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
376  TEUCHOS_TEST_FOR_EXCEPTION(Matrix::operatorViewTable_.containsKey(viewLabel) == false, Xpetra::Exceptions::RuntimeError, "Xpetra::Matrix.GetColMap(): view '" + viewLabel + "' does not exist.");
377  updateDefaultView(); // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated.
378  return Matrix::operatorViewTable_.get(viewLabel)->GetColMap();
379  }
380 
381  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
383  matrixData_->removeEmptyProcessesInPlace(newMap);
384  this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetRowMap(matrixData_->getRowMap());
385  this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetColMap(matrixData_->getColMap());
386  }
387 
388  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
390  return matrixData_->getMap();
391  }
392 
393  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
396  const CrsMatrixWrap & sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
397  matrixData_->doImport(*sourceWrp.getCrsMatrix(), importer, CM);
398  }
399 
400  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
403  const CrsMatrixWrap & destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
404  matrixData_->doExport(*destWrp.getCrsMatrix(), importer, CM);
405  }
406 
407  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
410  const CrsMatrixWrap & sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
411  matrixData_->doImport(*sourceWrp.getCrsMatrix(), exporter, CM);
412  }
413 
414  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
417  const CrsMatrixWrap & destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
418  matrixData_->doExport(*destWrp.getCrsMatrix(), exporter, CM);
419  }
420 
421  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
423  return "Xpetra::CrsMatrixWrap";
424  }
425 
426  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
428  // Teuchos::EVerbosityLevel vl = verbLevel;
429  // if (vl == VERB_DEFAULT) vl = VERB_LOW;
430  // RCP<const Comm<int> > comm = this->getComm();
431  // const int myImageID = comm->getRank(),
432  // numImages = comm->getSize();
433 
434  // if (myImageID == 0) out << this->description() << std::endl;
435 
436  matrixData_->describe(out,verbLevel);
437 
438  // Teuchos::OSTab tab(out);
439  }
440 
441  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
444  matrixData_->setObjectLabel(objectLabel);
445  }
446 
447 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
448 #ifdef HAVE_XPETRA_TPETRA
449 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
450  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
453  return getLocalMatrixDevice();
454  }
455 #endif
456  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
458  CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalMatrixHost () const {
459  return matrixData_->getLocalMatrixHost();
460  }
461  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
463  CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>::getLocalMatrixDevice () const {
464  return matrixData_->getLocalMatrixDevice();
465  }
466 #else
467 #ifdef __GNUC__
468 #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."
469 #endif
470 #endif
471 #endif
472 
473  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
475 
476 
477  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
479 
480 
481  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
483 
484 // Default view is created after fillComplete()
485  // Because ColMap might not be available before fillComplete().
486  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
488 
489  // Create default view
490  this->defaultViewLabel_ = "point";
491  this->CreateView(this->GetDefaultViewLabel(), matrixData_->getRowMap(), matrixData_->getColMap());
492 
493  // Set current view
494  this->currentViewLabel_ = this->GetDefaultViewLabel();
495  }
496 
497  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
499  if ((finalDefaultView_ == false) && matrixData_->isFillComplete() ) {
500  // Update default view with the colMap
501  Matrix::operatorViewTable_.get(Matrix::GetDefaultViewLabel())->SetColMap(matrixData_->getColMap());
502  finalDefaultView_ = true;
503  }
504  }
505 
506 
507  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
512  matrixData_->residual(X,B,R);
513  }
514 
515 
516 } //namespace Xpetra
517 
518 #endif //ifndef XPETRA_CRSMATRIXWRAP_DEF_HPP
virtual void setObjectLabel(const std::string &objectLabel)
static RCP< CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap)
Constructor for empty matrix (intended use is an import/export target - can't insert entries directly...
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...
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.
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.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
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.
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.
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
Exception throws to report errors in the internal logical of the program.
Xpetra-specific matrix class.
const viewLabel_t & GetDefaultViewLabel() const
const viewLabel_t & GetCurrentViewLabel() const
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Xpetra namespace
std::string viewLabel_t
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.