Xpetra Version of the Day
Loading...
Searching...
No Matches
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 <Tpetra_KokkosCompat_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
69namespace Xpetra {
70 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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_TPETRA
119 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
121 : finalDefaultView_(false)
122 {
123 // Set matrix data
124 matrixData_ = CrsMatrixFactory::Build(rowMap, colMap, lclMatrix, params);
125
126 // Default view
128 }
129
130 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
132 const RCP<const Map>& domainMap, const RCP<const Map>& rangeMap,
134 : finalDefaultView_(false)
135 {
136 // Set matrix data
137 matrixData_ = CrsMatrixFactory::Build(lclMatrix, rowMap, colMap, domainMap, rangeMap, params);
138
139 // Default view
141 }
142#else
143#ifdef __GNUC__
144#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."
145#endif
146#endif
147
148 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
151 {
152 // Set matrix data
153 matrixData_ = matrix;
154
155 // Default view
157 }
158
159 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
161 : finalDefaultView_(false)
162 {
163 // Set matrix data
164 matrixData_ = CrsMatrixFactory::Build(graph, paramList);
165
166 // Default view
168 }
169
170
171 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
174
175 const RCP<ParameterList>& paramList)
176 : finalDefaultView_(false)
177 {
178 // Set matrix data
179 matrixData_ = CrsMatrixFactory::Build(graph, values, paramList);
180
181 // Default view
183 }
184
185
186 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
188
189 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
191 matrixData_->insertGlobalValues(globalRow, cols, vals);
192 }
193
194 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
196 matrixData_->insertLocalValues(localRow, cols, vals);
197 }
198
199 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
202 const ArrayView<const Scalar> &vals) { matrixData_->replaceGlobalValues(globalRow, cols, vals); }
203
204 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
207 const ArrayView<const Scalar> &vals) { matrixData_->replaceLocalValues(localRow, cols, vals); }
208
209 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
210 void CrsMatrixWrap<Scalar,LocalOrdinal,GlobalOrdinal,Node>::setAllToScalar(const Scalar &alpha) { matrixData_->setAllToScalar(alpha); }
211
212 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
214 matrixData_->scale(alpha);
215 }
216
217 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
221
222 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
224 matrixData_->fillComplete(domainMap, rangeMap, params);
225
226 // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
228 }
229
230 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
232 matrixData_->fillComplete(params);
233
234 // Update default view with the colMap because colMap can be <tt>null</tt> until fillComplete() is called.
236 }
237
238 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
242
243 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
247
248 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
250 return matrixData_->getLocalNumRows();
251 }
252
253 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
257
258 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
260 return matrixData_->getLocalNumEntries();
261 }
262
263 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
265 return matrixData_->getNumEntriesInLocalRow(localRow);
266 }
267
268 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
270 return matrixData_->getNumEntriesInGlobalRow(globalRow);
271 }
272
273 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
275 return matrixData_->getGlobalMaxNumRowEntries();
276 }
277
278 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
280 return matrixData_->getLocalMaxNumRowEntries();
281 }
282
283 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
287
288 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
292
293 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
297
298 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300 const ArrayView<LocalOrdinal> &Indices,
301 const ArrayView<Scalar> &Values,
302 size_t &NumEntries
303 ) const {
304 matrixData_->getLocalRowCopy(LocalRow, Indices, Values, NumEntries);
305 }
306
307 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
309 matrixData_->getGlobalRowView(GlobalRow, indices, values);
310 }
311
312 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
314 matrixData_->getLocalRowView(LocalRow, indices, values);
315 }
316
317 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
321
322 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
326
327 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
331
332 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
336
337 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
341
342 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
346
347 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
351
352 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
361
362 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
365 Teuchos::ETransp mode,
366 Scalar alpha,
367 Scalar beta,
368 bool sumInterfaceValues,
369 const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> >& regionInterfaceImporter,
370 const Teuchos::ArrayRCP<LocalOrdinal>& regionInterfaceLIDs
371 ) const{
372 matrixData_->apply(X,Y,mode,alpha,beta,sumInterfaceValues,regionInterfaceImporter,regionInterfaceLIDs);
373 }
374
375 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
379
380 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
384
385 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
387
388 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
390 TEUCHOS_TEST_FOR_EXCEPTION(Matrix::operatorViewTable_.containsKey(viewLabel) == false, Xpetra::Exceptions::RuntimeError, "Xpetra::Matrix.GetColMap(): view '" + viewLabel + "' does not exist.");
391 updateDefaultView(); // If CrsMatrix::fillComplete() have been used instead of CrsMatrixWrap::fillComplete(), the default view is updated.
392 return Matrix::operatorViewTable_.get(viewLabel)->GetColMap();
393 }
394
395 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
397 matrixData_->removeEmptyProcessesInPlace(newMap);
398 this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetRowMap(matrixData_->getRowMap());
399 this->operatorViewTable_.get(this->GetCurrentViewLabel())->SetColMap(matrixData_->getColMap());
400 }
401
402 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
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(), importer, 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(), importer, CM);
419 }
420
421 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
424 const CrsMatrixWrap & sourceWrp = dynamic_cast<const CrsMatrixWrap &>(source);
425 matrixData_->doImport(*sourceWrp.getCrsMatrix(), exporter, CM);
426 }
427
428 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
431 const CrsMatrixWrap & destWrp = dynamic_cast<const CrsMatrixWrap &>(dest);
432 matrixData_->doExport(*destWrp.getCrsMatrix(), exporter, CM);
433 }
434
435 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
437 return "Xpetra::CrsMatrixWrap";
438 }
439
440 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
442 // Teuchos::EVerbosityLevel vl = verbLevel;
443 // if (vl == VERB_DEFAULT) vl = VERB_LOW;
444 // RCP<const Comm<int> > comm = this->getComm();
445 // const int myImageID = comm->getRank(),
446 // numImages = comm->getSize();
447
448 // if (myImageID == 0) out << this->description() << std::endl;
449
450 matrixData_->describe(out,verbLevel);
451
452 // Teuchos::OSTab tab(out);
453 }
454
455 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
458 matrixData_->setObjectLabel(objectLabel);
459 }
460
461#ifdef HAVE_XPETRA_TPETRA
462 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
467 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
472#else
473#ifdef __GNUC__
474#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."
475#endif
476#endif
477
478 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
480
481
482 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
484
485
486 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
488
489// Default view is created after fillComplete()
490 // Because ColMap might not be available before fillComplete().
491 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
493
494 // Create default view
495 this->defaultViewLabel_ = "point";
496 this->CreateView(this->GetDefaultViewLabel(), matrixData_->getRowMap(), matrixData_->getColMap());
497
498 // Set current view
500 }
501
502 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
504 if ((finalDefaultView_ == false) && matrixData_->isFillComplete() ) {
505 // Update default view with the colMap
507 finalDefaultView_ = true;
508 }
509 }
510
511
512 // Expert only
513 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
515 // Clear the old view table
517 Matrix::operatorViewTable_ = dummy_table;
518
519 finalDefaultView_ = M->isFillComplete();
520 // Set matrix data
521 matrixData_ = M;
522
523
524 // Default view
526 }
527
528
529 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
531 return matrixData_->GetStorageBlockSize();
532 }
533
534 template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
541
542
543} //namespace Xpetra
544
545#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)
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale matrix using the given vector entries.
LocalOrdinal GetStorageBlockSize() const
Returns the block size of the storage mechanism, which is usually 1, except for Tpetra::BlockCrsMatri...
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.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
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.
virtual local_matrix_type::HostMirror getLocalMatrixHost() const
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
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.
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.
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.
virtual local_matrix_type getLocalMatrixDevice() const
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void replaceCrsMatrix(RCP< CrsMatrix > &M)
Expert only.
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.
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.
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.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
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.
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...
Exception throws to report errors in the internal logical of the program.
void CreateView(viewLabel_t viewLabel, const RCP< const Map > &rowMap, const RCP< const Map > &colMap)
Teuchos::Hashtable< viewLabel_t, RCP< MatrixView > > operatorViewTable_
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
basic_FancyOStream< char > FancyOStream
Xpetra namespace
std::string viewLabel_t
size_t global_size_t
Global size_t object.
CombineMode
Xpetra::Combine Mode enumerable type.