Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TpetraCrsMatrix_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Xpetra: A linear algebra interface package
4//
5// Copyright 2012 NTESS and the Xpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef XPETRA_TPETRACRSMATRIX_DEF_HPP
11#define XPETRA_TPETRACRSMATRIX_DEF_HPP
12
13#include <Xpetra_MultiVectorFactory.hpp>
16
17namespace Xpetra {
18
19template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
21 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), maxNumEntriesPerRow, params))) {}
22
23template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
25 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), NumEntriesPerRowToAlloc(), params))) {}
26
27template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
29 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), maxNumEntriesPerRow, params))) {}
30
31template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
33 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), NumEntriesPerRowToAlloc(), params))) {}
34
35template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
38
39template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42
43template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
50 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
52
53 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
54 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
55 mtx_ = Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(importer), myDomainMap, myRangeMap, params);
56 bool restrictComm = false;
57 if (!params.is_null()) restrictComm = params->get("Restrict Communicator", restrictComm);
58 if (restrictComm && mtx_->getRowMap().is_null()) mtx_ = Teuchos::null;
59}
60
61template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
68 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
70
71 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
72 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
73 mtx_ = Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(exporter), myDomainMap, myRangeMap, params);
74}
75
76template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
84 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
86
87 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
88 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
89
90 mtx_ = Tpetra::importAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(RowImporter), toTpetra(*DomainImporter), myDomainMap, myRangeMap, params);
91 bool restrictComm = false;
92 if (!params.is_null()) restrictComm = params->get("Restrict Communicator", restrictComm);
93 if (restrictComm && mtx_->getRowMap().is_null()) mtx_ = Teuchos::null;
94}
95
96template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
104 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, *sourceMatrix, tSourceMatrix, "Xpetra::TpetraCrsMatrix constructor only accepts Xpetra::TpetraCrsMatrix as the input argument."); // TODO: remove and use toTpetra()
105 RCP<const Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> v = tSourceMatrix.getTpetra_CrsMatrix();
106
107 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myDomainMap = domainMap != Teuchos::null ? toTpetra(domainMap) : Teuchos::null;
108 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> myRangeMap = rangeMap != Teuchos::null ? toTpetra(rangeMap) : Teuchos::null;
109
110 mtx_ = Tpetra::exportAndFillCompleteCrsMatrix<MyTpetraCrsMatrix>(tSourceMatrix.getTpetra_CrsMatrix(), toTpetra(RowExporter), toTpetra(*DomainExporter), myDomainMap, myRangeMap, params);
111}
112
114
115template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
118 const local_matrix_type &lclMatrix,
120 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclMatrix, params))) {}
121
122template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
124 const local_matrix_type &lclMatrix,
130 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) {}
131
132template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
134 const local_matrix_type &lclMatrix,
142 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), toTpetra(importer), toTpetra(exporter), params))) {}
143
144template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
146
147template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
149 XPETRA_MONITOR("TpetraCrsMatrix::insertGlobalValues");
150 mtx_->insertGlobalValues(globalRow, cols, vals);
151}
152
153template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
155 XPETRA_MONITOR("TpetraCrsMatrix::insertLocalValues");
156 mtx_->insertLocalValues(localRow, cols, vals);
157}
158
159template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
161 XPETRA_MONITOR("TpetraCrsMatrix::replaceGlobalValues");
162 mtx_->replaceGlobalValues(globalRow, cols, vals);
163}
164
165template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
168 const ArrayView<const Scalar> &vals) {
169 XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
170 typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
171 const LocalOrdinal numValid =
172 mtx_->replaceLocalValues(localRow, cols, vals);
174 static_cast<size_type>(numValid) != cols.size(), std::runtime_error,
175 "replaceLocalValues returned " << numValid << " != cols.size() = " << cols.size() << ".");
176}
177
178template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
180 XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar");
181 mtx_->setAllToScalar(alpha);
182}
183
184template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
186 XPETRA_MONITOR("TpetraCrsMatrix::scale");
187 mtx_->scale(alpha);
188}
189
190template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
192 XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues");
193 rowptr.resize(getLocalNumRows() + 1);
194 colind.resize(numNonZeros);
195 values.resize(numNonZeros);
196}
197
198template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
200 XPETRA_MONITOR("TpetraCrsMatrix::setAllValues");
201 mtx_->setAllValues(rowptr, colind, values);
202}
203
204template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
206 XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
207 // TODO: Change getAllValues interface to return Kokkos views,
208 // TODO: or add getLocalRowPtrsHost, getLocalIndicesHost, etc., interfaces
209 // TODO: and use them in MueLu
210 rowptr = Kokkos::Compat::persistingView(mtx_->getLocalRowPtrsHost());
211 colind = Kokkos::Compat::persistingView(mtx_->getLocalIndicesHost());
212 values = Teuchos::arcp_reinterpret_cast<const Scalar>(
213 Kokkos::Compat::persistingView(mtx_->getLocalValuesHost(
214 Tpetra::Access::ReadOnly)));
215}
216
217template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
219 XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
220 // TODO: Change getAllValues interface to return Kokkos view,
221 // TODO: or add getLocalValuesDevice() interfaces
222 values = Teuchos::arcp_reinterpret_cast<Scalar>(
223 Kokkos::Compat::persistingView(mtx_->getLocalValuesDevice(
224 Tpetra::Access::ReadWrite)));
225}
226
227template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
229
230template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
232 XPETRA_MONITOR("TpetraCrsMatrix::resumeFill");
233 mtx_->resumeFill(params);
234}
235
236template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
238 XPETRA_MONITOR("TpetraCrsMatrix::fillComplete");
239 mtx_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params);
240}
241
242template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
244 XPETRA_MONITOR("TpetraCrsMatrix::fillComplete");
245 mtx_->fillComplete(params);
246}
247
248template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
250 XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
251 XPETRA_DYNAMIC_CAST(const TpetraImportClass, *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
252 RCP<const Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> myImport = tImporter.getTpetra_Import();
253 mtx_->replaceDomainMapAndImporter(toTpetra(newDomainMap), myImport);
254}
255
256template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
258 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
259 const RCP<const Import<LocalOrdinal, GlobalOrdinal, Node>> &importer,
260 const RCP<const Export<LocalOrdinal, GlobalOrdinal, Node>> &exporter,
261 const RCP<ParameterList> &params) {
262 XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
265
266 if (importer != Teuchos::null) {
267 XPETRA_DYNAMIC_CAST(const TpetraImportClass, *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
268 myImport = tImporter.getTpetra_Import();
269 }
270 if (exporter != Teuchos::null) {
271 XPETRA_DYNAMIC_CAST(const TpetraExportClass, *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
272 myExport = tExporter.getTpetra_Export();
273 }
274
275 mtx_->expertStaticFillComplete(toTpetra(domainMap), toTpetra(rangeMap), myImport, myExport, params);
276}
277
278template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
283
284template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
289
290template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
295
296template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
298 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows");
299 return mtx_->getGlobalNumRows();
300}
301
302template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
304 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols");
305 return mtx_->getGlobalNumCols();
306}
307
308template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
310 XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumRows");
311 return mtx_->getLocalNumRows();
312}
313
314template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
316 XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumCols");
317 return mtx_->getLocalNumCols();
318}
319
320template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
322 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries");
323 return mtx_->getGlobalNumEntries();
324}
325
326template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
328 XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumEntries");
329 return mtx_->getLocalNumEntries();
330}
331
332template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
334 XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow");
335 return mtx_->getNumEntriesInLocalRow(localRow);
336}
337
338template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
340 XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInGlobalRow");
341 return mtx_->getNumEntriesInGlobalRow(globalRow);
342}
343
344template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
346 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries");
347 return mtx_->getGlobalMaxNumRowEntries();
348}
349
350template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
352 XPETRA_MONITOR("TpetraCrsMatrix::getLocalMaxNumRowEntries");
353 return mtx_->getLocalMaxNumRowEntries();
354}
355
356template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
358 XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed");
359 return mtx_->isLocallyIndexed();
360}
361
362template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
364 XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed");
365 return mtx_->isGloballyIndexed();
366}
367
368template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
370 XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete");
371 return mtx_->isFillComplete();
372}
373
374template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
376 XPETRA_MONITOR("TpetraCrsMatrix::isFillActive");
377 return mtx_->isFillActive();
378}
379
380template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
382 XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm");
383 return mtx_->getFrobeniusNorm();
384}
385
386template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
388 XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews");
389 return mtx_->supportsRowViews();
390}
391
392template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
393void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView<LocalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const {
394 XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowCopy");
397
398 mtx_->getLocalRowCopy(LocalRow, indices, values, NumEntries);
399 for (size_t i = 0; i < NumEntries; ++i) {
400 Indices[i] = indices(i);
401 Values[i] = values(i);
402 }
403}
404
405template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
406void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView<GlobalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const {
407 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowCopy");
410
411 mtx_->getGlobalRowCopy(GlobalRow, indices, values, NumEntries);
412 for (size_t i = 0; i < NumEntries; ++i) {
413 Indices[i] = indices(i);
414 Values[i] = values(i);
415 }
416}
417
418template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
420 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView");
423
424 mtx_->getGlobalRowView(GlobalRow, indices, values);
425 Indices = ArrayView<const GlobalOrdinal>(indices.data(), indices.extent(0));
426 Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
427}
428
429template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
431 XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView");
434
435 mtx_->getLocalRowView(LocalRow, indices, values);
436 Indices = ArrayView<const LocalOrdinal>(indices.data(), indices.extent(0));
437 Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
438}
439
440template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
442 XPETRA_MONITOR("TpetraCrsMatrix::apply");
443 mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
444}
445
446template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
447void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode, Scalar alpha, Scalar beta, bool sumInterfaceValues, const RCP<Import<LocalOrdinal, GlobalOrdinal, Node>> &regionInterfaceImporter, const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const {
448 XPETRA_MONITOR("TpetraCrsMatrix::apply(region)");
449 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> regionInterfaceMap = regionInterfaceImporter->getTargetMap();
450 mtx_->localApply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
451 if (sumInterfaceValues) {
452 // preform communication to propagate local interface
453 // values to all the processor that share interfaces.
455 matvecInterfaceTmp->doImport(Y, *regionInterfaceImporter, INSERT);
456
457 // sum all contributions to interface values
458 // on all ranks
460 ArrayRCP<Scalar> interfaceData = matvecInterfaceTmp->getDataNonConst(0);
461 for (LocalOrdinal interfaceIdx = 0; interfaceIdx < static_cast<LocalOrdinal>(interfaceData.size()); ++interfaceIdx) {
462 YData[regionInterfaceLIDs[interfaceIdx]] += interfaceData[interfaceIdx];
463 }
464 }
465}
466
467template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
469 XPETRA_MONITOR("TpetraCrsMatrix::getDomainMap");
470 return toXpetra(mtx_->getDomainMap());
471}
472
473template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
478
479template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
481 XPETRA_MONITOR("TpetraCrsMatrix::description");
482 return mtx_->description();
483}
484
485template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
487 XPETRA_MONITOR("TpetraCrsMatrix::describe");
488 mtx_->describe(out, verbLevel);
489}
490
491template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
493 XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel");
495 mtx_->setObjectLabel(objectLabel);
496}
497
498template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
500 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(*(matrix.mtx_), Teuchos::Copy))) {}
501
502template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
504 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
505 XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
506 mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
507}
508
509template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
511 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
512 mtx_->getLocalDiagOffsets(offsets);
513}
514
515template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
517 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
518 mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
519}
520
521template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
522void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalDiagCopy(Vector &diag, const Kokkos::View<const size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets) const {
523 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
524 mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
525}
526
527template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
532
533template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
535 XPETRA_MONITOR("TpetraCrsMatrix::leftScale");
536 mtx_->leftScale(*(toTpetra(x)));
537}
538
539template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
541 XPETRA_MONITOR("TpetraCrsMatrix::rightScale");
542 mtx_->rightScale(*(toTpetra(x)));
543}
544
545template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
550
551template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
554 XPETRA_MONITOR("TpetraCrsMatrix::doImport");
555
556 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
558 // mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
559 mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
560}
561
563template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
566 XPETRA_MONITOR("TpetraCrsMatrix::doExport");
567
568 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
570 mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
571}
572
573template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
576 XPETRA_MONITOR("TpetraCrsMatrix::doImport");
577
578 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
580 mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
581}
582
583template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
586 XPETRA_MONITOR("TpetraCrsMatrix::doExport");
587
588 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
590 mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
591}
592
593template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
595 XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
596 mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
597}
598
599template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
601 return !mtx_.is_null();
602}
603
604template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
607
608template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
610
612template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
614
616template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
622
625// End of TpetrCrsMatrix class definition //
628
629} // namespace Xpetra
630
631#endif // XPETRA_TPETRACRSMATRIX_DEF_HPP
#define XPETRA_MONITOR(funcName)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_type size() const
void resize(const size_type n, const T &val=T())
size_type size() const
virtual void setObjectLabel(const std::string &objectLabel)
bool is_null() const
T * get() const
void getLocalRowView(LocalOrdinal LocalRow, local_inds_host_view_type &indices, values_host_view_type &values) const override
void getGlobalRowView(GlobalOrdinal GlobalRow, global_inds_host_view_type &indices, values_host_view_type &values) const override
static Teuchos::RCP< MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > Build(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &map, size_t NumVectors, bool zeroOut=true)
Constructor specifying the number of non-zeros for all rows.
virtual Teuchos::ArrayRCP< Scalar > getDataNonConst(size_t j)=0
View of the local values in a particular vector of this multivector.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs.
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
RCP< Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrixNonConst() const
Get the underlying Tpetra matrix.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local IDs.
TpetraCrsMatrix(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rowMap, size_t maxNumEntriesPerRow, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
Constructor specifying fixed number of entries for each row.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
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.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs.
void getLocalDiagCopy(Vector &diag) const
Get a copy of the diagonal entries owned by this node, with local row indices.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the (locally owned) global row.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this.
void rightScale(const Vector &x)
Right scale operator with given vector values.
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...
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
void fillComplete(const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMap, const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &rangeMap, const RCP< ParameterList > &params=null)
Signal that data entry is complete, specifying domain and range maps.
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
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 residual(const MultiVector &X, const MultiVector &B, MultiVector &R) const
Compute a residual R = B - (*this) * X.
void doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void resumeFill(const RCP< ParameterList > &params=null)
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs.
void leftScale(const Vector &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.
bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
bool hasMatrix() const
Does this have an underlying matrix.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
std::string description() const
A simple one-line description of this object.
RCP< const Tpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_CrsMatrix() const
Get the underlying Tpetra matrix.
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
bool isFillActive() const
Returns true if the matrix is in edit mode.
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.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
size_t getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in 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 doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
global_size_t getGlobalNumRows() const
Number of global elements in the row map of 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 getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
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.
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 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 > &params=Teuchos::null)
Expert static fill complete.
void replaceDiag(const Vector &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.
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
void setObjectLabel(const std::string &objectLabel)
void apply(const MultiVector &X, MultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, Scalar alpha=ScalarTraits< Scalar >::one(), Scalar beta=ScalarTraits< Scalar >::zero()) const
Computes the sparse matrix-multivector multiplication.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this matrix.
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....
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
size_t getLocalNumEntries() const
Returns the local number of entries in this matrix.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
void residual(const Operator< SC, LO, GO, NO > &Aop, const MultiVector< SC, LO, GO, NO > &X_in, const MultiVector< SC, LO, GO, NO > &B_in, MultiVector< SC, LO, GO, NO > &R_in)
LO replaceDiagonalCrsMatrix(::Tpetra::CrsMatrix< SC, LO, GO, NT > &matrix, const ::Tpetra::Vector< SC, LO, GO, NT > &newDiag)
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toXpetra(RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > graph)
size_t global_size_t
Global size_t object.
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
CombineMode
Xpetra::Combine Mode enumerable type.