10 #ifndef TPETRA_SOLVERMAP_CRSMATRIX_DEF_HPP
11 #define TPETRA_SOLVERMAP_CRSMATRIX_DEF_HPP
30 template <
class Scalar,
35 : StructuralSameTypeTransform<
CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node> >()
36 , newColMap_(Teuchos::null)
37 , newGraph_(Teuchos::null) {
41 template <
class Scalar,
49 template <
class Scalar,
53 typename SolverMap_CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::NewType
57 assert(!origMatrix->isGloballyIndexed());
59 this->origObj_ = origMatrix;
64 Teuchos::RCP<map_t const> origDomainMap = origMatrix->getDomainMap();
65 Teuchos::RCP<map_t const> origColMap = origMatrix->getColMap();
67 typename map_t::local_map_type localOrigDomainMap(origDomainMap->getLocalMap());
68 typename map_t::local_map_type localOrigColMap(origColMap->getLocalMap());
70 if (origDomainMap->isLocallyFitted(*origColMap)) {
71 this->newObj_ = this->origObj_;
82 size_t const origDomainMap_localSize = origDomainMap->getLocalNumElements();
83 size_t newColMap_localSize(origDomainMap_localSize);
86 size_t newColMap_extraSize(0);
87 size_t const origColMap_localSize = origColMap->getLocalNumElements();
89 Kokkos::parallel_reduce(
90 "Tpetra::SolverMap_CrsMatrix::construct::newColMap_extraSize", Kokkos::RangePolicy<typename Node::device_type::execution_space, size_t>(0, origColMap_localSize), KOKKOS_LAMBDA(
size_t const i,
size_t& sizeToUpdate)->
void {
91 GlobalOrdinal
const globalColIndex(localOrigColMap.getGlobalElement(i));
92 if (localOrigDomainMap.getLocalElement(globalColIndex) == ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid()) {
98 newColMap_localSize += newColMap_extraSize;
101 Kokkos::View<GlobalOrdinal*, typename Node::device_type> newColMap_globalColIndices(
"", newColMap_localSize);
105 Kokkos::parallel_for(
106 "Tpetra::SolverMap_CrsMatrix::construct::copyDomainMapToNewColMap", Kokkos::RangePolicy<typename Node::device_type::execution_space, size_t>(0, origDomainMap_localSize), KOKKOS_LAMBDA(
size_t const i)->
void {
107 newColMap_globalColIndices(i) = localOrigDomainMap.getGlobalElement(i);
116 Kokkos::parallel_scan(
117 "Tpetra::SolverMap_CrsMatrix::construct::appendNewColMap", Kokkos::RangePolicy<typename Node::device_type::execution_space, size_t>(0, origColMap_localSize), KOKKOS_LAMBDA(
size_t const i,
size_t& jToUpdate,
bool const final)->
void {
118 GlobalOrdinal
const globalColIndex(localOrigColMap.getGlobalElement(i));
119 if (localOrigDomainMap.getLocalElement(globalColIndex) == ::Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid()) {
121 newColMap_globalColIndices(origDomainMap_localSize + jToUpdate) = globalColIndex;
131 Teuchos::RCP<map_t const> origRowMap = origMatrix->getRowMap();
132 Teuchos::RCP<Teuchos::Comm<int>
const> Comm = origRowMap->getComm();
133 size_t const newColMap_localNumCols = newColMap_globalColIndices.size();
134 size_t newColMap_globalNumCols(0);
135 Teuchos::reduceAll(*Comm, Teuchos::REDUCE_SUM, 1, &newColMap_localNumCols, &newColMap_globalNumCols);
137 newColMap_ = Teuchos::rcp<map_t>(
new map_t(newColMap_globalNumCols, newColMap_globalColIndices, origDomainMap->getIndexBase(), Comm));
142 cg_t
const* origGraph =
dynamic_cast<cg_t const*
>(origMatrix->getGraph().get());
143 newGraph_ = Teuchos::rcp<cg_t>(
new cg_t(*origGraph));
144 newGraph_->resumeFill();
145 newGraph_->reindexColumns(newColMap_);
146 newGraph_->fillComplete();
151 Teuchos::RCP<cm_t> newMatrix = Teuchos::rcp<cm_t>(
new cm_t(*origMatrix, newGraph_));
156 this->newObj_ = newMatrix;
159 return this->newObj_;
168 #define TPETRA_SOLVERMAPCRSMATRIX_INSTANT(SCALAR, LO, GO, NODE) \
169 template class SolverMap_CrsMatrix<SCALAR, LO, GO, NODE>;
173 #endif // TPETRA_SOLVERMAP_CRSMATRIX_DEF_HPP
Sparse matrix that presents a row-oriented interface that lets users read or modify entries...
Declaration of the Tpetra::SolverMap_CrsMatrix class.
NewType operator()(OriginalType const &origMatrix)
A distributed graph accessed by rows (adjacency lists) and stored sparsely.
A parallel distribution of indices over processes.