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>
15#include "Tpetra_Details_residual.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>
37 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(graph), params))) {}
38
39template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
42
43template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
49 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
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>
67 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
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>
83 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
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>
103 typedef Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> MyTpetraCrsMatrix;
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
115#ifdef HAVE_XPETRA_TPETRA
116template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
119 const local_matrix_type &lclMatrix,
121 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(toTpetra(rowMap), toTpetra(colMap), lclMatrix, params))) {}
122
123template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
125 const local_matrix_type &lclMatrix,
131 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), params))) {}
132
133template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
135 const local_matrix_type &lclMatrix,
143 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(lclMatrix, toTpetra(rowMap), toTpetra(colMap), toTpetra(domainMap), toTpetra(rangeMap), toTpetra(importer), toTpetra(exporter), params))) {}
144#endif
145
146template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148
149template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
151 XPETRA_MONITOR("TpetraCrsMatrix::insertGlobalValues");
152 mtx_->insertGlobalValues(globalRow, cols, vals);
153}
154
155template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
157 XPETRA_MONITOR("TpetraCrsMatrix::insertLocalValues");
158 mtx_->insertLocalValues(localRow, cols, vals);
159}
160
161template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
163 XPETRA_MONITOR("TpetraCrsMatrix::replaceGlobalValues");
164 mtx_->replaceGlobalValues(globalRow, cols, vals);
165}
166
167template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
170 const ArrayView<const Scalar> &vals) {
171 XPETRA_MONITOR("TpetraCrsMatrix::replaceLocalValues");
172 typedef typename ArrayView<const LocalOrdinal>::size_type size_type;
173 const LocalOrdinal numValid =
174 mtx_->replaceLocalValues(localRow, cols, vals);
176 static_cast<size_type>(numValid) != cols.size(), std::runtime_error,
177 "replaceLocalValues returned " << numValid << " != cols.size() = " << cols.size() << ".");
178}
179
180template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
182 XPETRA_MONITOR("TpetraCrsMatrix::setAllToScalar");
183 mtx_->setAllToScalar(alpha);
184}
185
186template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
188 XPETRA_MONITOR("TpetraCrsMatrix::scale");
189 mtx_->scale(alpha);
190}
191
192template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
194 XPETRA_MONITOR("TpetraCrsMatrix::allocateAllValues");
195 rowptr.resize(getLocalNumRows() + 1);
196 colind.resize(numNonZeros);
197 values.resize(numNonZeros);
198}
199
200template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
202 XPETRA_MONITOR("TpetraCrsMatrix::setAllValues");
203 mtx_->setAllValues(rowptr, colind, values);
204}
205
206template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
208 XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
209 // TODO: Change getAllValues interface to return Kokkos views,
210 // TODO: or add getLocalRowPtrsHost, getLocalIndicesHost, etc., interfaces
211 // TODO: and use them in MueLu
212 rowptr = Kokkos::Compat::persistingView(mtx_->getLocalRowPtrsHost());
213 colind = Kokkos::Compat::persistingView(mtx_->getLocalIndicesHost());
214 values = Teuchos::arcp_reinterpret_cast<const Scalar>(
215 Kokkos::Compat::persistingView(mtx_->getLocalValuesHost(
216 Tpetra::Access::ReadOnly)));
217}
218
219template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
221 XPETRA_MONITOR("TpetraCrsMatrix::getAllValues");
222 // TODO: Change getAllValues interface to return Kokkos view,
223 // TODO: or add getLocalValuesDevice() interfaces
224 values = Teuchos::arcp_reinterpret_cast<Scalar>(
225 Kokkos::Compat::persistingView(mtx_->getLocalValuesDevice(
226 Tpetra::Access::ReadWrite)));
227}
228
229template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
231
232template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
234 XPETRA_MONITOR("TpetraCrsMatrix::resumeFill");
235 mtx_->resumeFill(params);
236}
237
238template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
240 XPETRA_MONITOR("TpetraCrsMatrix::fillComplete");
241 mtx_->fillComplete(toTpetra(domainMap), toTpetra(rangeMap), params);
242}
243
244template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246 XPETRA_MONITOR("TpetraCrsMatrix::fillComplete");
247 mtx_->fillComplete(params);
248}
249
250template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
252 XPETRA_MONITOR("TpetraCrsMatrix::replaceDomainMapAndImporter");
253 XPETRA_DYNAMIC_CAST(const TpetraImportClass, *newImporter, tImporter, "Xpetra::TpetraCrsMatrix::replaceDomainMapAndImporter only accepts Xpetra::TpetraImport.");
254 RCP<const Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> myImport = tImporter.getTpetra_Import();
255 mtx_->replaceDomainMapAndImporter(toTpetra(newDomainMap), myImport);
256}
257
258template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
260 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> &rangeMap,
261 const RCP<const Import<LocalOrdinal, GlobalOrdinal, Node>> &importer,
262 const RCP<const Export<LocalOrdinal, GlobalOrdinal, Node>> &exporter,
263 const RCP<ParameterList> &params) {
264 XPETRA_MONITOR("TpetraCrsMatrix::expertStaticFillComplete");
267
268 if (importer != Teuchos::null) {
269 XPETRA_DYNAMIC_CAST(const TpetraImportClass, *importer, tImporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraImport.");
270 myImport = tImporter.getTpetra_Import();
271 }
272 if (exporter != Teuchos::null) {
273 XPETRA_DYNAMIC_CAST(const TpetraExportClass, *exporter, tExporter, "Xpetra::TpetraCrsMatrix::expertStaticFillComplete only accepts Xpetra::TpetraExport.");
274 myExport = tExporter.getTpetra_Export();
275 }
276
277 mtx_->expertStaticFillComplete(toTpetra(domainMap), toTpetra(rangeMap), myImport, myExport, params);
278}
279
280template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
285
286template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
291
292template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
297
298template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
300 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumRows");
301 return mtx_->getGlobalNumRows();
302}
303
304template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
306 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumCols");
307 return mtx_->getGlobalNumCols();
308}
309
310template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
312 XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumRows");
313 return mtx_->getLocalNumRows();
314}
315
316template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
318 XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumCols");
319 return mtx_->getLocalNumCols();
320}
321
322template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
324 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalNumEntries");
325 return mtx_->getGlobalNumEntries();
326}
327
328template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
330 XPETRA_MONITOR("TpetraCrsMatrix::getLocalNumEntries");
331 return mtx_->getLocalNumEntries();
332}
333
334template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
336 XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInLocalRow");
337 return mtx_->getNumEntriesInLocalRow(localRow);
338}
339
340template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
342 XPETRA_MONITOR("TpetraCrsMatrix::getNumEntriesInGlobalRow");
343 return mtx_->getNumEntriesInGlobalRow(globalRow);
344}
345
346template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
348 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalMaxNumRowEntries");
349 return mtx_->getGlobalMaxNumRowEntries();
350}
351
352template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
354 XPETRA_MONITOR("TpetraCrsMatrix::getLocalMaxNumRowEntries");
355 return mtx_->getLocalMaxNumRowEntries();
356}
357
358template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
360 XPETRA_MONITOR("TpetraCrsMatrix::isLocallyIndexed");
361 return mtx_->isLocallyIndexed();
362}
363
364template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
366 XPETRA_MONITOR("TpetraCrsMatrix::isGloballyIndexed");
367 return mtx_->isGloballyIndexed();
368}
369
370template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
372 XPETRA_MONITOR("TpetraCrsMatrix::isFillComplete");
373 return mtx_->isFillComplete();
374}
375
376template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
378 XPETRA_MONITOR("TpetraCrsMatrix::isFillActive");
379 return mtx_->isFillActive();
380}
381
382template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
384 XPETRA_MONITOR("TpetraCrsMatrix::getFrobeniusNorm");
385 return mtx_->getFrobeniusNorm();
386}
387
388template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
390 XPETRA_MONITOR("TpetraCrsMatrix::supportsRowViews");
391 return mtx_->supportsRowViews();
392}
393
394template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
395void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalRowCopy(LocalOrdinal LocalRow, const ArrayView<LocalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const {
396 XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowCopy");
397 typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::nonconst_local_inds_host_view_type indices("indices", Indices.size());
398 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type values("values", Values.size());
399
400 mtx_->getLocalRowCopy(LocalRow, indices, values, NumEntries);
401 for (size_t i = 0; i < NumEntries; ++i) {
402 Indices[i] = indices(i);
403 Values[i] = values(i);
404 }
405}
406
407template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
408void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getGlobalRowCopy(GlobalOrdinal GlobalRow, const ArrayView<GlobalOrdinal> &Indices, const ArrayView<Scalar> &Values, size_t &NumEntries) const {
409 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowCopy");
410 typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::nonconst_global_inds_host_view_type indices("indices", Indices.size());
411 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::nonconst_values_host_view_type values("values", Values.size());
412
413 mtx_->getGlobalRowCopy(GlobalRow, indices, values, NumEntries);
414 for (size_t i = 0; i < NumEntries; ++i) {
415 Indices[i] = indices(i);
416 Values[i] = values(i);
417 }
418}
419
420template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
422 XPETRA_MONITOR("TpetraCrsMatrix::getGlobalRowView");
423 typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::global_inds_host_view_type indices;
424 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type values;
425
426 mtx_->getGlobalRowView(GlobalRow, indices, values);
427 Indices = ArrayView<const GlobalOrdinal>(indices.data(), indices.extent(0));
428 Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
429}
430
431template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
433 XPETRA_MONITOR("TpetraCrsMatrix::getLocalRowView");
434 typename Tpetra::CrsGraph<LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type indices;
435 typename Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::values_host_view_type values;
436
437 mtx_->getLocalRowView(LocalRow, indices, values);
438 Indices = ArrayView<const LocalOrdinal>(indices.data(), indices.extent(0));
439 Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
440}
441
442template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
444 XPETRA_MONITOR("TpetraCrsMatrix::apply");
445 mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
446}
447
448template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
449void 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 {
450 XPETRA_MONITOR("TpetraCrsMatrix::apply(region)");
451 RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>> regionInterfaceMap = regionInterfaceImporter->getTargetMap();
452 mtx_->localApply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
453 if (sumInterfaceValues) {
454 // preform communication to propagate local interface
455 // values to all the processor that share interfaces.
457 matvecInterfaceTmp->doImport(Y, *regionInterfaceImporter, INSERT);
458
459 // sum all contributions to interface values
460 // on all ranks
462 ArrayRCP<Scalar> interfaceData = matvecInterfaceTmp->getDataNonConst(0);
463 for (LocalOrdinal interfaceIdx = 0; interfaceIdx < static_cast<LocalOrdinal>(interfaceData.size()); ++interfaceIdx) {
464 YData[regionInterfaceLIDs[interfaceIdx]] += interfaceData[interfaceIdx];
465 }
466 }
467}
468
469template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
471 XPETRA_MONITOR("TpetraCrsMatrix::getDomainMap");
472 return toXpetra(mtx_->getDomainMap());
473}
474
475template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
480
481template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
483 XPETRA_MONITOR("TpetraCrsMatrix::description");
484 return mtx_->description();
485}
486
487template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
489 XPETRA_MONITOR("TpetraCrsMatrix::describe");
490 mtx_->describe(out, verbLevel);
491}
492
493template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
495 XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel");
497 mtx_->setObjectLabel(objectLabel);
498}
499
500template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
502 : mtx_(Teuchos::rcp(new Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(*(matrix.mtx_), Teuchos::Copy))) {}
503
504template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
506 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
507 XPETRA_DYNAMIC_CAST(TpetraVectorClass, diag, tDiag, "Xpetra::TpetraCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
508 mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
509}
510
511template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
513 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagOffsets");
514 mtx_->getLocalDiagOffsets(offsets);
515}
516
517template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
519 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
520 mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
521}
522
523template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
524void TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::getLocalDiagCopy(Vector &diag, const Kokkos::View<const size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets) const {
525 XPETRA_MONITOR("TpetraCrsMatrix::getLocalDiagCopy");
526 mtx_->getLocalDiagCopy(*(toTpetra(diag)), offsets);
527}
528
529template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
531 XPETRA_MONITOR("TpetraCrsMatrix::replaceDiag");
532 Tpetra::replaceDiagonalCrsMatrix(*mtx_, *(toTpetra(diag)));
533}
534
535template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
537 XPETRA_MONITOR("TpetraCrsMatrix::leftScale");
538 mtx_->leftScale(*(toTpetra(x)));
539}
540
541template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
543 XPETRA_MONITOR("TpetraCrsMatrix::rightScale");
544 mtx_->rightScale(*(toTpetra(x)));
545}
546
547template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
552
553template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
556 XPETRA_MONITOR("TpetraCrsMatrix::doImport");
557
558 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
560 // mtx_->doImport(toTpetraCrsMatrix(source), *tImporter.getTpetra_Import(), toTpetra(CM));
561 mtx_->doImport(*v, toTpetra(importer), toTpetra(CM));
562}
563
565template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
568 XPETRA_MONITOR("TpetraCrsMatrix::doExport");
569
570 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
572 mtx_->doExport(*v, toTpetra(importer), toTpetra(CM));
573}
574
575template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
578 XPETRA_MONITOR("TpetraCrsMatrix::doImport");
579
580 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, source, tSource, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
582 mtx_->doImport(*v, toTpetra(exporter), toTpetra(CM));
583}
584
585template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
588 XPETRA_MONITOR("TpetraCrsMatrix::doExport");
589
590 XPETRA_DYNAMIC_CAST(const TpetraCrsMatrixClass, dest, tDest, "Xpetra::TpetraCrsMatrix::doImport only accept Xpetra::TpetraCrsMatrix as input arguments."); // TODO: remove and use toTpetra()
592 mtx_->doExport(*v, toTpetra(exporter), toTpetra(CM));
593}
594
595template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
597 XPETRA_MONITOR("TpetraCrsMatrix::removeEmptyProcessesInPlace");
598 mtx_->removeEmptyProcessesInPlace(toTpetra(newMap));
599}
600
601template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
603 return !mtx_.is_null();
604}
605
606template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
607TpetraCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::TpetraCrsMatrix(const Teuchos::RCP<Tpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>> &mtx)
608 : mtx_(mtx) {}
609
610template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
612
614template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
616
618template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
620 const MultiVector &B,
621 MultiVector &R) const {
622 Tpetra::Details::residual(*mtx_, toTpetra(X), toTpetra(B), toTpetra(R));
623}
624
627// End of TpetrCrsMatrix class definition //
630
631} // namespace Xpetra
632
633#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
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.
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.
Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
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)
size_t global_size_t
Global size_t object.
RCP< const CrsGraph< int, GlobalOrdinal, Node > > toXpetra(const Epetra_CrsGraph &g)
RCP< const Tpetra::CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > toTpetra(const RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > &graph)
CombineMode
Xpetra::Combine Mode enumerable type.