46 #ifndef XPETRA_TPETRACRSMATRIX_DEF_HPP
47 #define XPETRA_TPETRACRSMATRIX_DEF_HPP
53 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
55 : mtx_ (Teuchos::
rcp (new Tpetra::
CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (
toTpetra(rowMap), maxNumEntriesPerRow,
toTpetra(pftype), params))) { }
57 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
59 : mtx_(Teuchos::
rcp(new Tpetra::
CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> (
toTpetra(rowMap), NumEntriesPerRowToAlloc(),
toTpetra(pftype), params))) { }
61 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
62 TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(
const Teuchos::RCP<
const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap,
const Teuchos::RCP<
const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap,
size_t maxNumEntriesPerRow,
ProfileType pftype,
const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
65 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
66 TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::TpetraCrsMatrix(
const Teuchos::RCP<
const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap,
const Teuchos::RCP<
const Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap,
const ArrayRCP< const size_t > &NumEntriesPerRowToAlloc,
ProfileType pftype,
const Teuchos::RCP< Teuchos::ParameterList > ¶ms)
67 : mtx_(Teuchos::
rcp(new Tpetra::
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(rowMap),
toTpetra(colMap), NumEntriesPerRowToAlloc(),
toTpetra(pftype), params))) { }
69 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
71 : mtx_(Teuchos::
rcp(new Tpetra::
CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(
toTpetra(graph), params))) { }
75 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
82 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
88 mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),
toTpetra(importer),myDomainMap,myRangeMap,params);
89 bool restrictComm=
false;
90 if(!params.
is_null()) restrictComm = params->
get(
"Restrict Communicator",restrictComm);
91 if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
95 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
102 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
108 mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),
toTpetra(exporter),myDomainMap,myRangeMap,params);
112 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
120 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
127 mtx_=Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),
toTpetra(RowImporter),
toTpetra(*DomainImporter),myDomainMap,myRangeMap,params);
128 bool restrictComm=
false;
129 if(!params.
is_null()) restrictComm = params->
get(
"Restrict Communicator",restrictComm);
130 if(restrictComm && mtx_->getRowMap().is_null()) mtx_=Teuchos::null;
133 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
141 typedef Tpetra::CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node> MyTpetraCrsMatrix;
148 mtx_=Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(),
toTpetra(RowExporter),
toTpetra(*DomainExporter),myDomainMap,myRangeMap,params);
154 #ifdef HAVE_XPETRA_KOKKOS_REFACTOR
155 #ifdef HAVE_XPETRA_TPETRA
156 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
159 const local_matrix_type& lclMatrix,
161 : mtx_(Teuchos::
rcp(new Tpetra::
CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(
toTpetra(rowMap),
toTpetra(colMap), lclMatrix, params))) { }
163 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
165 const local_matrix_type& lclMatrix,
166 const Teuchos::RCP<
const Map<LocalOrdinal,GlobalOrdinal,Node> >& rowMap,
167 const Teuchos::RCP<
const Map<LocalOrdinal,GlobalOrdinal,Node> >& colMap,
168 const Teuchos::RCP<
const Map<LocalOrdinal,GlobalOrdinal,Node> >& domainMap,
169 const Teuchos::RCP<
const Map<LocalOrdinal,GlobalOrdinal,Node> >& rangeMap,
171 : mtx_(Teuchos::
rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix,
toTpetra(rowMap),
toTpetra(colMap),
toTpetra(domainMap),
toTpetra(rangeMap), params))) { }
175 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
178 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
181 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
184 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
187 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
195 const LocalOrdinal numValid =
196 mtx_->replaceLocalValues (localRow, cols, vals);
198 static_cast<size_type> (numValid) != cols.
size (), std::runtime_error,
199 "replaceLocalValues returned " << numValid <<
" != cols.size() = " <<
200 cols.
size () <<
".");
203 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
206 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
209 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
213 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
215 {
XPETRA_MONITOR(
"TpetraCrsMatrix::setAllValues"); mtx_->setAllValues(rowptr,colind,values); }
217 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
219 {
XPETRA_MONITOR(
"TpetraCrsMatrix::getAllValues"); mtx_->getAllValues(rowptr,colind,values); }
221 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
223 {
return mtx_->haveGlobalConstants();}
225 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
228 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
229 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::fillComplete(
const RCP<
const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap,
const RCP<
const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap,
const RCP< ParameterList > ¶ms) {
XPETRA_MONITOR(
"TpetraCrsMatrix::fillComplete"); mtx_->fillComplete(
toTpetra(domainMap),
toTpetra(rangeMap), params); }
231 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
235 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
240 mtx_->replaceDomainMapAndImporter(
toTpetra(newDomainMap),myImport);
243 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
253 if(importer!=Teuchos::null) {
255 myImport = tImporter.getTpetra_Import();
257 if(exporter!=Teuchos::null) {
259 myExport = tExporter.getTpetra_Export();
262 mtx_->expertStaticFillComplete(
toTpetra(domainMap),
toTpetra(rangeMap),myImport,myExport,params);
265 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
268 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
271 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
274 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
277 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
280 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
283 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
286 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
289 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
292 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
295 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
298 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
301 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
304 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
307 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
310 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
313 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
316 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
319 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
322 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
325 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
328 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
331 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
332 void TpetraCrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>::apply(
const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X,
MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &Y,
Teuchos::ETransp mode, Scalar alpha, Scalar beta)
const {
XPETRA_MONITOR(
"TpetraCrsMatrix::apply"); mtx_->apply(
toTpetra(X),
toTpetra(Y), mode, alpha, beta); }
334 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
337 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
340 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
343 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
346 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
349 Teuchos::LabeledObject::setObjectLabel(objectLabel);
350 mtx_->setObjectLabel(objectLabel);
355 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
357 : mtx_(Teuchos::
rcp(new Tpetra::
CrsMatrix<Scalar,LocalOrdinal,GlobalOrdinal,Node>(*(matrix.mtx_),Teuchos::
Copy))) {}
359 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
363 mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
366 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
369 mtx_->getLocalDiagOffsets(offsets);
372 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
375 mtx_->getLocalDiagCopy(*(
toTpetra(diag)), offsets);
378 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
381 Tpetra::replaceDiagonalCrsMatrix(*mtx_, *(
toTpetra(diag)));
384 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
390 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
396 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
399 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
411 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
422 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
433 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
444 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
447 mtx_->removeEmptyProcessesInPlace(
toTpetra(newMap));
450 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
452 return ! mtx_.is_null ();
455 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
458 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
462 template <
class Scalar,
class LocalOrdinal,
class GlobalOrdinal,
class Node>
473 #endif // XPETRA_TPETRACRSMATRIX_DEF_HPP
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
std::string description() const
A simple one-line description of this object.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getDomainMap() const
Returns the Map associated with the domain of this operator. This will be null until fillComplete() i...
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
T & get(const std::string &name, T def_value)
void apply(const MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &X, 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.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void allocateAllValues(size_t numNonZeros, ArrayRCP< size_t > &rowptr, ArrayRCP< LocalOrdinal > &colind, ArrayRCP< Scalar > &values)
Allocates and returns ArrayRCPs of the Crs arrays — This is an Xpetra-only routine.
void resumeFill(const RCP< ParameterList > ¶ms=null)
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
size_t getNodeNumEntries() const
Returns the local number of entries in this matrix.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, ProfileType pftype=DynamicProfile, const Teuchos::RCP< Teuchos::ParameterList > ¶ms=Teuchos::null)
Constructor specifying fixed number of entries for each row.
void setObjectLabel(const std::string &objectLabel)
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
size_t getNodeNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix. ...
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const
Gets the 1D pointer arrays of the graph.
void replaceDomainMapAndImporter(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newDomainMap, Teuchos::RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &newImporter)
Replaces the current domainMap and importer with the user-specified objects.
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
void resize(const size_type n, const T &val=T())
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
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 expertStaticFillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< const Import< LocalOrdinal, GlobalOrdinal, Node > > &importer=Teuchos::null, const RCP< const Export< LocalOrdinal, GlobalOrdinal, Node > > &exporter=Teuchos::null, const RCP< ParameterList > ¶ms=Teuchos::null)
Expert static fill complete.
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 describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Print the object with some verbosity level to an FancyOStream object.
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_t global_size_t
Global size_t object.
virtual ~TpetraCrsMatrix()
Destructor.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRangeMap() const
Returns the Map associated with the range of this operator, which must be compatible with Y...
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
bool hasMatrix() const
Does this have an underlying matrix.
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
CombineMode
Xpetra::Combine Mode enumerable type.
#define XPETRA_MONITOR(funcName)
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
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.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
size_t getNodeNumRows() const
Returns the number of matrix rows owned on the calling node.
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > ¶ms=null)
Signal that data entry is complete, specifying domain and range maps.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView< GlobalOrdinal > &indices, const ArrayView< Scalar > &values, size_t &numEntries) const
Extract a list of entries in a specified global row of this matrix. Put into pre-allocated storage...
size_t getNodeMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
bool isFillActive() const
Returns true if the matrix is in edit mode.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.