Xpetra Version of the Day
Loading...
Searching...
No Matches
Xpetra_TpetraBlockCrsMatrix_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_TPETRABLOCKCRSMATRIX_DEF_HPP
11#define XPETRA_TPETRABLOCKCRSMATRIX_DEF_HPP
12
14#include "Xpetra_TpetraCrsGraph.hpp"
15
16namespace Xpetra {
17
19template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
22 size_t maxNumEntriesPerRow,
24 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
25}
26
28template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
31 const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc,
33 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
34}
35
37template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
41 size_t maxNumEntriesPerRow,
43 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
44}
45
47template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
51 const ArrayRCP<const size_t> &NumEntriesPerRowToAlloc,
53 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
54}
55
57template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
61// : mtx_(Teuchos::rcp(new Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >(toTpetra(graph), params)))
62// * there is no Tpetra::BlockCrsMatrix(graph, params) c'tor. We throw anyways here so no need to set mtx_.
63{
64 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
65}
66
68template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
71 const LocalOrdinal blockSize)
72 : mtx_(Teuchos::rcp(new Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(*toTpetra(graph), blockSize))) {}
73
75template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
78 const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &pointDomainMap,
79 const Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &pointRangeMap,
80 const LocalOrdinal blockSize)
81 : mtx_(Teuchos::rcp(new Tpetra::BlockCrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>(*toTpetra(graph), *toTpetra(pointDomainMap), *toTpetra(pointRangeMap), blockSize))) {}
82
84template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
91 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
92}
93
95template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
102 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
103}
104
106template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
110 const Teuchos::RCP<const Import<LocalOrdinal, GlobalOrdinal, Node> > DomainImporter,
114 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
115}
116
118template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
122 const Teuchos::RCP<const Export<LocalOrdinal, GlobalOrdinal, Node> > DomainExporter,
126 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
127}
128
130template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
133
135
137template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
139 insertGlobalValues(GlobalOrdinal globalRow,
141 const ArrayView<const Scalar> &vals) {
142 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
143}
144
146template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
148 insertLocalValues(LocalOrdinal localRow,
150 const ArrayView<const Scalar> &vals) {
151 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
152}
153
155template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
157 replaceGlobalValues(GlobalOrdinal globalRow,
159 const ArrayView<const Scalar> &vals) {
160 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in" + std::string(__FILE__) + ":" + std::to_string(__LINE__));
161}
162
164template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
166 replaceLocalValues(LocalOrdinal localRow, const ArrayView<const LocalOrdinal> &cols, const ArrayView<const Scalar> &vals) {
167 XPETRA_MONITOR("TpetraBlockCrsMatrix::replaceLocalValues");
168 mtx_->replaceLocalValues(localRow, cols.getRawPtr(), vals.getRawPtr(), cols.size());
169}
170
172template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
174 setAllToScalar(const Scalar &alpha) {
175 XPETRA_MONITOR("TpetraBlockCrsMatrix::setAllToScalar");
176 mtx_->setAllToScalar(alpha);
177}
178
180template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
182 scale(const Scalar &alpha) {
183 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
184}
185
187//** \warning This is an expert-only routine and should not be called from user code. (not implemented)
188template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
190 allocateAllValues(size_t numNonZeros, ArrayRCP<size_t> &rowptr, ArrayRCP<LocalOrdinal> &colind, ArrayRCP<Scalar> &values) {
191 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
192}
193
195template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
197 setAllValues(const ArrayRCP<size_t> &rowptr, const ArrayRCP<LocalOrdinal> &colind, const ArrayRCP<Scalar> &values) {
198 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
199}
200
202template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
206 ArrayRCP<const Scalar> &values) const {
207 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
208}
209
211template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
214 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
215}
216
218
219// Transformational Methods
221
222template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
227
228template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
235
236template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
241
242template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
246 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
247}
248
249template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
252 const RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> > &rangeMap,
253 const RCP<const Import<LocalOrdinal, GlobalOrdinal, Node> > &importer,
254 const RCP<const Export<LocalOrdinal, GlobalOrdinal, Node> > &exporter,
255 const RCP<ParameterList> &params) {
256 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
257}
258
260
262
263template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
266 getRowMap() const {
267 XPETRA_MONITOR("TpetraBlockCrsMatrix::getRowMap");
268 return toXpetra(mtx_->getRowMap());
269}
270
271template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
274 getColMap() const {
275 XPETRA_MONITOR("TpetraBlockCrsMatrix::getColMap");
276 return toXpetra(mtx_->getColMap());
277}
278
279template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
282 getCrsGraph() const {
283 XPETRA_MONITOR("TpetraBlockCrsMatrix::getCrsGraph");
286 RCP<G_t> t_graph = Teuchos::rcp_const_cast<G_t>(Teuchos::rcpFromRef(mtx_->getCrsGraph()));
287 RCP<const G_x> x_graph = rcp(new G_x(t_graph));
288 return x_graph;
289}
290
291template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
294 getGlobalNumRows() const {
295 XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumRows");
296 return mtx_->getGlobalNumRows();
297}
298
299template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
302 getGlobalNumCols() const {
303 XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumCols");
304 return mtx_->getGlobalNumCols();
305}
306
307template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
308size_t
310 getLocalNumRows() const {
311 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalNumRows");
312 return mtx_->getLocalNumRows();
313}
314
315template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
316size_t
318 getLocalNumCols() const {
319 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalNumCols");
320 return mtx_->getLocalNumCols();
321}
322
323template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
326 getGlobalNumEntries() const {
327 XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalNumEntries");
328 return mtx_->getGlobalNumEntries();
329}
330
331template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
332size_t
334 getLocalNumEntries() const {
335 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalNumEntries");
336 return mtx_->getLocalNumEntries();
337}
338
339template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
340size_t
342 getNumEntriesInLocalRow(LocalOrdinal localRow) const {
343 XPETRA_MONITOR("TpetraBlockCrsMatrix::getNumEntriesInLocalRow");
344 return mtx_->getNumEntriesInLocalRow(localRow);
345}
346
347template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
348size_t
350 getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const {
351 XPETRA_MONITOR("TpetraBlockCrsMatrix::getNumEntriesInGlobalRow");
352 return mtx_->getNumEntriesInGlobalRow(globalRow);
353}
354
355template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
357 XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalMaxNumRowEntries");
358 return mtx_->getGlobalMaxNumRowEntries();
359}
360
361template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
363 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalMaxNumRowEntries");
364 return mtx_->getLocalMaxNumRowEntries();
365}
366
367template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
369 XPETRA_MONITOR("TpetraBlockCrsMatrix::isLocallyIndexed");
370 return mtx_->isLocallyIndexed();
371}
372
373template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
375 XPETRA_MONITOR("TpetraBlockCrsMatrix::isGloballyIndexed");
376 return mtx_->isGloballyIndexed();
377}
378
379template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
381 XPETRA_MONITOR("TpetraBlockCrsMatrix::isFillComplete");
382 return mtx_->isFillComplete();
383}
384
385template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
387 XPETRA_MONITOR("TpetraBlockCrsMatrix::isFillActive");
388 return false;
389}
390
391template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
393 XPETRA_MONITOR("TpetraBlockCrsMatrix::getFrobeniusNorm");
394 return mtx_->getFrobeniusNorm();
395}
396
397template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
399 XPETRA_MONITOR("TpetraBlockCrsMatrix::supportsRowViews");
400 return mtx_->supportsRowViews();
401}
402
403template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
405 getLocalRowCopy(LocalOrdinal LocalRow,
406 const ArrayView<LocalOrdinal> &Indices,
407 const ArrayView<Scalar> &Values,
408 size_t &NumEntries) const {
409 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalRowCopy");
412
413 mtx_->getLocalRowCopy(LocalRow, indices, values, NumEntries);
414 for (size_t i = 0; i < NumEntries; ++i) {
415 Indices[i] = indices(i);
416 Values[i] = values(i);
417 }
418}
419template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
421 getLocalRowView(LocalOrdinal LocalRow, ArrayView<const LocalOrdinal> &Indices,
422 ArrayView<const Scalar> &Values) const {
423 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalRowView");
426
427 mtx_->getLocalRowView(LocalRow, indices, values);
428 Indices = ArrayView<const LocalOrdinal>(indices.data(), indices.extent(0));
429 Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
430}
431
432template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
434 getGlobalRowView(GlobalOrdinal GlobalRow,
436 ArrayView<const Scalar> &Values) const {
437 XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalRowView");
440
441 mtx_->getGlobalRowView(GlobalRow, indices, values);
442 Indices = ArrayView<const GlobalOrdinal>(indices.data(), indices.extent(0));
443 Values = ArrayView<const Scalar>(reinterpret_cast<const Scalar *>(values.data()), values.extent(0));
444}
445
446template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
448 getGlobalRowCopy(GlobalOrdinal GlobalRow,
449 const ArrayView<GlobalOrdinal> &Indices,
450 const ArrayView<Scalar> &Values,
451 size_t &NumEntries) const {
452 XPETRA_MONITOR("TpetraBlockCrsMatrix::getGlobalRowCopy");
455
456 mtx_->getGlobalRowCopy(GlobalRow, indices, values, NumEntries);
457 for (size_t i = 0; i < NumEntries; ++i) {
458 Indices[i] = indices(i);
459 Values[i] = values(i);
460 }
461}
462
463template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
466
467template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
471 Teuchos::ETransp mode,
472 Scalar alpha,
473 Scalar beta) const {
474 XPETRA_MONITOR("TpetraBlockCrsMatrix::apply");
475 mtx_->apply(toTpetra(X), toTpetra(Y), mode, alpha, beta);
476}
477
478template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
482 Teuchos::ETransp mode,
483 Scalar alpha,
484 Scalar beta,
485 bool sumInterfaceValues,
486 const RCP<Import<LocalOrdinal, GlobalOrdinal, Node> > &regionInterfaceImporter,
487 const Teuchos::ArrayRCP<LocalOrdinal> &regionInterfaceLIDs) const {}
488
489template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
492 getDomainMap() const {
493 XPETRA_MONITOR("TpetraBlockCrsMatrix::getDomainMap");
494 return toXpetra(mtx_->getDomainMap());
495}
496
497template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
500 getRangeMap() const {
501 XPETRA_MONITOR("TpetraBlockCrsMatrix::getRangeMap");
502 return toXpetra(mtx_->getRangeMap());
503}
504
505template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
506std::string
508 description() const {
509 XPETRA_MONITOR("TpetraBlockCrsMatrix::description");
510 return mtx_->description();
511}
512
513template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
516 const Teuchos::EVerbosityLevel verbLevel) const {
517 XPETRA_MONITOR("TpetraBlockCrsMatrix::describe");
518 mtx_->describe(out, verbLevel);
519}
520
521template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
523 setObjectLabel(const std::string &objectLabel) {
524 XPETRA_MONITOR("TpetraCrsMatrix::setObjectLabel");
526 mtx_->setObjectLabel(objectLabel);
527}
528
529template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
532 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalDiagCopy");
534 diag,
535 tDiag,
536 "Xpetra::TpetraBlockCrsMatrix.getLocalDiagCopy() only accept Xpetra::TpetraVector as input arguments.");
537 mtx_->getLocalDiagCopy(*tDiag.getTpetra_Vector());
538}
539
541template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
544 const Teuchos::ArrayView<const size_t> &offsets) const {
545 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
546}
547
549template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
552 const Kokkos::View<const size_t *, typename Node::device_type, Kokkos::MemoryUnmanaged> &offsets) const {
553 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
554}
555
556template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
559 XPETRA_MONITOR("TpetraBlockCrsMatrix::getLocalDiagOffsets");
560
561 const size_t lclNumRows = mtx_->getGraph()->getLocalNumRows();
562 if (static_cast<size_t>(offsets.size()) < lclNumRows) {
563 offsets.resize(lclNumRows);
564 }
565
566 // The input ArrayRCP must always be a host pointer. Thus, if
567 // device_type::memory_space is Kokkos::HostSpace, it's OK for us
568 // to write to that allocation directly as a Kokkos::View.
569 typedef typename Node::device_type device_type;
570 typedef typename device_type::memory_space memory_space;
571 if (std::is_same<memory_space, Kokkos::HostSpace>::value) {
572 // It is always syntactically correct to assign a raw host
573 // pointer to a device View, so this code will compile correctly
574 // even if this branch never runs.
575 typedef Kokkos::View<size_t *, device_type, Kokkos::MemoryUnmanaged> output_type;
576 output_type offsetsOut(offsets.getRawPtr(), offsets.size());
577 mtx_->getLocalDiagOffsets(offsetsOut);
578 } else {
579 Kokkos::View<size_t *, device_type> offsetsTmp("diagOffsets", offsets.size());
580 mtx_->getLocalDiagOffsets(offsetsTmp);
581 typedef Kokkos::View<size_t *, Kokkos::HostSpace, Kokkos::MemoryUnmanaged> output_type;
582 output_type offsetsOut(offsets.getRawPtr(), offsets.size());
583 Kokkos::deep_copy(offsetsOut, offsetsTmp);
584 }
585}
586
587template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
590 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix::replaceDiag: function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
591}
592
593template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
596 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
597}
598
599template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
602 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
603}
604
605template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
608 getMap() const {
609 XPETRA_MONITOR("TpetraBlockCrsMatrix::getMap");
610 return rcp(new TpetraMap<LocalOrdinal, GlobalOrdinal, Node>(mtx_->getMap()));
611}
612
614template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
618 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
619}
620
622template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
626 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
627}
628
630template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
634 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
635}
636
638template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
642 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
643}
644
645template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
648 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix function not implemented in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
649}
650
651template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
653 hasMatrix() const {
654 return !mtx_.is_null();
655}
656
657template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
661
662template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
668
669// TODO: remove
670template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
676
677// was: typedef typename Xpetra::CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_type local_matrix_type;
678// using local_matrix_type = typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_matrix_type;
679
680template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
683 getLocalMatrixDevice() const {
684 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix does not support getLocalMatrix due to missing Kokkos::CrsMatrix in Tpetra's experimental implementation in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
685
686#ifndef __NVCC__
688#endif // __NVCC__
689
691}
692
693template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
696 getLocalMatrixHost() const {
697 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix does not support getLocalMatrix due to missing Kokkos::CrsMatrix in Tpetra's experimental implementation in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
698
699#ifndef __NVCC__
700 typename local_matrix_type::host_mirror_type ret;
701#endif // __NVCC__
702
704}
705
706template <class Scalar, class LocalOrdinal, class GlobalOrdinal, class Node>
708 setAllValues(const typename local_matrix_type::row_map_type &ptr,
709 const typename local_matrix_type::StaticCrsGraphType::entries_type::non_const_type &ind,
710 const typename local_matrix_type::values_type &val) {
711 throw std::runtime_error("Xpetra::TpetraBlockCrsMatrix does not support setAllValues due to missing Kokkos::CrsMatrix in Tpetra's experimental implementation in " + std::string(__FILE__) + ":" + std::to_string(__LINE__));
712}
713
714} // namespace Xpetra
715
716#endif // XPETRA_TPETRABLOCKCRSMATRIX_DEF_HPP
#define XPETRA_MONITOR(funcName)
#define XPETRA_DYNAMIC_CAST(type, obj, newObj, exceptionMsg)
size_type size() const
T * getRawPtr() const
void resize(const size_type n, const T &val=T())
size_type size() const
T * getRawPtr() const
virtual void setObjectLabel(const std::string &objectLabel)
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
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...
void insertGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using global IDs (not implemented)
TpetraBlockCrsMatrix(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 (not implemented)
void setObjectLabel(const std::string &objectLabel)
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getRowMap() const
Returns the Map that describes the row distribution in this matrix.
void setAllToScalar(const Scalar &alpha)
Set all matrix entries equal to scalarThis.
std::string description() const
A simple one-line description of this object.
virtual bool haveGlobalConstants() const
Returns true if globalConstants have been computed; false otherwise.
size_t getLocalNumCols() const
Returns the number of columns connected to the locally owned rows of this matrix.
typename Xpetra::CrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node >::local_matrix_type local_matrix_type
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 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...
size_t getGlobalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on all nodes.
local_matrix_type::host_mirror_type getLocalMatrixHost() const
bool supportsRowViews() const
Returns true if getLocalRowView() and getGlobalRowView() are valid for this class.
bool isLocallyIndexed() const
If matrix indices are in the local range, this function returns true. Otherwise, this function return...
global_size_t getGlobalNumRows() const
Number of global elements in the row map of this matrix.
size_t getNumEntriesInLocalRow(LocalOrdinal localRow) const
Returns the current number of entries on this node in the specified local row.
void leftScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Left scale operator with given vector values.
void getAllValues(ArrayRCP< const size_t > &rowptr, ArrayRCP< const LocalOrdinal > &colind, ArrayRCP< const Scalar > &values) const
Gets the 1D pointer arrays of the graph (not implemented)
size_t getLocalNumRows() const
Returns the number of matrix rows owned on the calling node.
void getLocalDiagCopy(Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag) const
Get a copy of the diagonal entries owned by this node, with local row idices.
bool hasMatrix() const
Does this have an underlying matrix.
void getLocalDiagOffsets(Teuchos::ArrayRCP< size_t > &offsets) const
Get offsets of the diagonal entries in the matrix.
void setAllValues(const ArrayRCP< size_t > &rowptr, const ArrayRCP< LocalOrdinal > &colind, const ArrayRCP< Scalar > &values)
Sets the 1D pointer arrays of the graph (not impelmented)
void replaceDiag(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &diag)
Replace the diagonal entries of the matrix.
global_size_t getGlobalNumEntries() const
Returns the global number of entries in this matrix.
bool isGloballyIndexed() const
If matrix indices are in the global range, this function returns true. Otherwise, this function retur...
void rightScale(const Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > &x)
Right scale operator with given vector values.
RCP< const Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrix() const
Get the underlying Tpetra matrix.
void scale(const Scalar &alpha)
Scale the current values of a matrix, this = alpha*this (not implemented)
const RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getColMap() const
Returns the Map that describes the column distribution in this matrix.
void insertLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Insert matrix entries, using local IDs (not implemented)
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.
RCP< const CrsGraph< LocalOrdinal, GlobalOrdinal, Node > > getCrsGraph() const
Returns the CrsGraph associated with this 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...
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....
size_t getLocalNumEntries() const
Returns the local number of entries in this 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.
ScalarTraits< Scalar >::magnitudeType getFrobeniusNorm() const
Returns the Frobenius norm of the matrix.
void resumeFill(const RCP< ParameterList > &params=null)
void doImport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &source, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Import.
global_size_t getGlobalNumCols() const
Number of global columns in the matrix.
size_t getNumEntriesInGlobalRow(GlobalOrdinal globalRow) const
Returns the current number of entries in the (locally owned) global row.
bool isFillComplete() const
Returns true if the matrix is in compute mode, i.e. if fillComplete() has been called.
void replaceLocalValues(LocalOrdinal localRow, const ArrayView< const LocalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using local 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.
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 doExport(const DistObject< char, LocalOrdinal, GlobalOrdinal, Node > &dest, const Import< LocalOrdinal, GlobalOrdinal, Node > &importer, CombineMode CM)
Export.
void replaceGlobalValues(GlobalOrdinal globalRow, const ArrayView< const GlobalOrdinal > &cols, const ArrayView< const Scalar > &vals)
Replace matrix entries, using global IDs (not implemented)
bool isFillActive() const
Returns true if the matrix is in edit mode.
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 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.
size_t getLocalMaxNumRowEntries() const
Returns the maximum number of entries across all rows/columns on this node.
RCP< Tpetra::BlockCrsMatrix< Scalar, LocalOrdinal, GlobalOrdinal, Node > > getTpetra_BlockCrsMatrixNonConst() const
Get the underlying Tpetra matrix.
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 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.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > getMap() const
Implements DistObject interface.
void removeEmptyProcessesInPlace(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &newMap)
#define TEUCHOS_UNREACHABLE_RETURN(dummyReturnVal)
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.