All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 
71  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
73  size_t maxNumEntriesPerRow,
74  Xpetra::ProfileType pftype)
75  : finalDefaultView_ (false)
76  {
77  matrixData_ = CrsMatrixFactory::Build (rowMap, maxNumEntriesPerRow, pftype);
79  }
80 
81  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
83  const ArrayRCP<const size_t>& NumEntriesPerRowToAlloc,
84  ProfileType pftype)
85  : finalDefaultView_ (false)
86  {
87  matrixData_ = CrsMatrixFactory::Build(rowMap, NumEntriesPerRowToAlloc, pftype);
89  }
90 
91  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
93  : finalDefaultView_(false)
94  {
95  // Set matrix data
96  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, maxNumEntriesPerRow, pftype);
97 
98  // Default view
100  }
101 
102  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
104  : finalDefaultView_(false)
105  {
106  // Set matrix data
107  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, NumEntriesPerRowToAlloc, pftype);
108 
109  // Default view
111  }
112 
113 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
114 #ifdef HAVE_XPETRA_TPETRA
115  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
116  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)
117  : finalDefaultView_(false)
118  {
119  // Set matrix data
120  matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, lclMatrix, params);
121 
122  // Default view
124  }
125 
126  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
127  CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>::CrsMatrixWrap(const local_matrix_type& lclMatrix, const RCP<const Map> &rowMap, const RCP<const Map>& colMap,
128  const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap,
130  : finalDefaultView_(false)
131  {
132  // Set matrix data
133  matrixData_ = CrsMatrixFactory::Build(lclMatrix, rowMap, colMap, domainMap, rangeMap, params);
134 
135  // Default view
137  }
138 #else
139 #ifdef __GNUC__
140 #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."
141 #endif
142 #endif
143 #endif
144 
145  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
147  : finalDefaultView_(matrix->isFillComplete())
148  {
149  // Set matrix data
150  matrixData_ = matrix;
151 
152  // Default view
154  }
155 
156  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
158  : finalDefaultView_(false)
159  {
160  // Set matrix data
161  matrixData_ = CrsMatrixFactory::Build(graph, paramList);
162 
163  // Default view
165  }
166 
167  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
169 
170  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
172  matrixData_->insertGlobalValues(globalRow, cols, vals);
173  }
174 
175  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
177  matrixData_->insertLocalValues(localRow, cols, vals);
178  }
179 
180  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
182  const ArrayView<const GlobalOrdinal> &cols,
183  const ArrayView<const Scalar> &vals) { matrixData_->replaceGlobalValues(globalRow, cols, vals); }
184 
185  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
187  const ArrayView<const LocalOrdinal> &cols,
188  const ArrayView<const Scalar> &vals) { matrixData_->replaceLocalValues(localRow, cols, vals); }
189 
190  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
191  void CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setAllToScalar(const Scalar &alpha) { matrixData_->setAllToScalar(alpha); }
192 
193  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
195  matrixData_->scale(alpha);
196  }
197 
198  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
200  matrixData_->resumeFill(params);
201  }
202 
203  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
205  matrixData_->fillComplete(domainMap, rangeMap, params);
206 
207  // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
208  updateDefaultView();
209  }
210 
211  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
213  matrixData_->fillComplete(params);
214 
215  // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
216  updateDefaultView();
217  }
218 
219  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
221  return matrixData_->getGlobalNumRows();
222  }
223 
224  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
226  return matrixData_->getGlobalNumCols();
227  }
228 
229  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
231  return matrixData_->getNodeNumRows();
232  }
233 
234  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
236  return matrixData_->getGlobalNumEntries();
237  }
238 
239  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
241  return matrixData_->getNodeNumEntries();
242  }
243 
244  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246  return matrixData_->getNumEntriesInLocalRow(localRow);
247  }
248 
249  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
251  return matrixData_->getNumEntriesInGlobalRow(globalRow);
252  }
253 
254  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
256  return matrixData_->getGlobalMaxNumRowEntries();
257  }
258 
259  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
261  return matrixData_->getNodeMaxNumRowEntries();
262  }
263 
264  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
266  return matrixData_->isLocallyIndexed();
267  }
268 
269  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
271  return matrixData_->isGloballyIndexed();
272  }
273 
274  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
276  return matrixData_->isFillComplete();
277  }
278 
279  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
281  const ArrayView<LocalOrdinal> &Indices,
282  const ArrayView<Scalar> &Values,
283  size_t &NumEntries
284  ) const {
285  matrixData_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
286  }
287 
288  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
290  matrixData_->getGlobalRowView(GlobalRow, indices, values);
291  }
292 
293  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
295  matrixData_->getLocalRowView(LocalRow, indices, values);
296  }
297 
298  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300  matrixData_->getLocalDiagCopy(diag);
301  }
302 
303  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
305  matrixData_->getLocalDiagOffsets(offsets);
306  }
307 
308  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
310  matrixData_->getLocalDiagCopy(diag,offsets);
311  }
312 
313  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
315  return matrixData_->getFrobeniusNorm();
316  }
317 
318  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
320  matrixData_->leftScale(x);
321  }
322 
323  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
325  matrixData_->rightScale(x);
326  }
327 
328  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
330  return matrixData_->haveGlobalConstants();
331  }
332 
333  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
336  Teuchos::ETransp mode,
337  Scalar alpha,
338  Scalar beta) const {
339 
340  matrixData_->apply(X,Y,mode,alpha,beta);
341  }
342 
343  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
345  return matrixData_->getDomainMap();
346  }
347 
348  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
350  return matrixData_->getRangeMap();
351  }
352 
353  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
355 
356  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
358  TEUCHOS_TEST_FOR_EXCEPTION(Matrix::operatorViewTable_.containsKey(viewLabel) == false, Xpetra::Exceptions::RuntimeError, "Xpetra::Matrix.GetColMap(): view '" + viewLabel + "' does not exist.");
359  updateDefaultView(); // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated.
360  return Matrix::operatorViewTable_.get(viewLabel)->GetColMap();
361  }
362 
363  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
365  matrixData_->removeEmptyProcessesInPlace(newMap);
366  this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetRowMap(matrixData_->getRowMap());
367  this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetColMap(matrixData_->getColMap());
368  }
369 
370  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
372  return matrixData_->getMap();
373  }
374 
375  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
378  const CrsMatrixWrap & sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
379  matrixData_->doImport(*sourceWrp.getCrsMatrix(), importer, CM);
380  }
381 
382  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
385  const CrsMatrixWrap & destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
386  matrixData_->doExport(*destWrp.getCrsMatrix(), importer, CM);
387  }
388 
389  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
392  const CrsMatrixWrap & sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
393  matrixData_->doImport(*sourceWrp.getCrsMatrix(), exporter, CM);
394  }
395 
396  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
399  const CrsMatrixWrap & destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
400  matrixData_->doExport(*destWrp.getCrsMatrix(), exporter, CM);
401  }
402 
403  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
405  return "Xpetra::CrsMatrixWrap";
406  }
407 
408  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
410  // Teuchos::EVerbosityLevel vl = verbLevel;
411  // if (vl == VERB_DEFAULT) vl = VERB_LOW;
412  // RCP<const Comm<int> > comm = this->getComm();
413  // const int myImageID = comm->getRank(),
414  // numImages = comm->getSize();
415 
416  // if (myImageID == 0) out << this->description() << std::endl;
417 
418  matrixData_->describe(out,verbLevel);
419 
420  // Teuchos::OSTab tab(out);
421  }
422 
423  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
425  Teuchos::LabeledObject::setObjectLabel(objectLabel);
426  matrixData_->setObjectLabel(objectLabel);
427  }
428 
429 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
430 #ifdef HAVE_XPETRA_TPETRA
431  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
434  return matrixData_->getLocalMatrix();
435  }
436 #else
437 #ifdef __GNUC__
438 #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."
439 #endif
440 #endif
441 #endif
442 
443  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
445 
446 
447  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
449 
450 
451  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
453 
454 #ifdef XPETRA_ENABLE_DEPRECATED_CODE
455  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
456  template <class Node2>
458 #ifdef HAVE_XPETRA_TPETRA
460  Teuchos::rcp_dynamic_cast<const Xpetra::TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> >(matrixData_);
461  if (tMatrix == Teuchos::null)
462  throw Xpetra::Exceptions::RuntimeError("clone() functionality is only available for Tpetra");
463 
465  // TODO: inherit strided maps/views ?
466 #else
467  return Teuchos::null;
468 #endif
469  }
470 #endif
471 
472  // Default view is created after fillComplete()
473  // Because ColMap might not be available before fillComplete().
474  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
476 
477  // Create default view
478  this->defaultViewLabel_ = "point";
479  this->CreateView(this->GetDefaultViewLabel(), matrixData_->getRowMap(), matrixData_->getColMap());
480 
481  // Set current view
482  this->currentViewLabel_ = this->GetDefaultViewLabel();
483  }
484 
485  template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
487  if ((finalDefaultView_ == false) && matrixData_->isFillComplete() ) {
488  // Update default view with the colMap
489  Matrix::operatorViewTable_.get(Matrix::GetDefaultViewLabel())->SetColMap(matrixData_->getColMap());
490  finalDefaultView_ = true;
491  }
492  }
493 
494 } //namespace Xpetra
495 
496 #endif //ifndef XPETRA_CRSMATRIXWRAP_DEF_HPP
void setObjectLabel(const std::string &objectLabel)
RCP< CrsMatrix > getCrsMatrix() const
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.
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...
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.
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.
Exception throws to report errors in the internal logical of the program.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map > &newMap)
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.
void doExport(const Matrix &dest, const Xpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
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.
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.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
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)
bool isFillComplete() const
Returns true if fillComplete() has been called and the matrix is in compute mode. ...
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...
virtual void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalar.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
Concrete implementation of Xpetra::Matrix.
global_size_t getGlobalNumCols() const
Returns the number of global columns in the matrix.
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...
const viewLabel_t & GetDefaultViewLabel() 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.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the specified global row.
RCP< const Xpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
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 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...