1#ifndef TPETRA_DETAILS_MAKECOLMAP_DEF_HPP
2#define TPETRA_DETAILS_MAKECOLMAP_DEF_HPP
24#include "KokkosCompat_View.hpp"
25#include "Kokkos_Atomic.hpp"
26#include "Kokkos_Core_fwd.hpp"
27#include "Kokkos_Pair.hpp"
28#include "Teuchos_Assert.hpp"
29#include "Tpetra_CrsGraph.hpp"
31#include "Teuchos_Array.hpp"
32#include "Kokkos_Bitset.hpp"
33#include "Kokkos_UnorderedMap.hpp"
34#include "Kokkos_Sort.hpp"
39template <
class keys_view_type>
40struct StableSortComparator {
43 StableSortComparator(keys_view_type keys_)
46 KOKKOS_INLINE_FUNCTION
47 bool operator()(
size_t i,
size_t j)
const {
48 if (keys(i) == keys(j)) {
51 return keys(i) < keys(j);
56template <
class ExecutionSpace,
57 class keys_view_type_in,
58 class keys_view_type_out,
59 class values_view_type_in,
60 class values_view_type_out>
61void sort_by_key(
const ExecutionSpace& exec,
62 const keys_view_type_in& keys_in,
63 const values_view_type_in& values_in,
64 const keys_view_type_out& keys_out,
65 const values_view_type_out& values_out) {
66 auto numItems = keys_in.extent(0);
67 Kokkos::View<size_t*, typename keys_view_type_out::memory_space> indices(Kokkos::ViewAllocateWithoutInitializing(
""), numItems);
69 Kokkos::RangePolicy<ExecutionSpace>(0, numItems), KOKKOS_LAMBDA(
const size_t k) { indices(k) = k; });
70 StableSortComparator comparison(keys_in);
71 Kokkos::sort(indices, comparison);
74 Kokkos::RangePolicy<ExecutionSpace>(0, numItems), KOKKOS_LAMBDA(
const size_t k) {
75 auto idx = indices(k);
76 values_out(k) = values_in(idx);
77 keys_out(k) = keys_in(idx);
81template <
class LO,
class GO,
class NT>
83 Teuchos::Array<int>& remotePIDs,
85 size_t numLocalColGIDs,
86 size_t numRemoteColGIDs,
87 const Kokkos::View<GO*, typename NT::memory_space> RemoteGIDs,
88 const Kokkos::View<bool*, typename NT::memory_space>& GIDisLocal,
89 std::ostream* errStrm) {
92 using execution_space =
typename NT::execution_space;
93 using memory_space =
typename NT::memory_space;
96 const char prefix[] =
"Tpetra::Details::makeColMapImpl: ";
127 if (domMap->getComm()->getSize() == 1) {
128 if (numRemoteColGIDs != 0) {
130 if (errStrm != NULL) {
131 *errStrm << prefix <<
"The domain Map only has one process, but "
132 << numRemoteColGIDs <<
" column "
133 << (numRemoteColGIDs != 1 ?
"indices are" :
"index is")
134 <<
" not in the domain Map. Either these indices are "
135 "invalid or the domain Map is invalid. Remember that nonsquare "
136 "matrices, or matrices where the row and range Maps differ, "
137 "require calling the version of fillComplete that takes the "
138 "domain and range Maps as input."
142 if (numLocalColGIDs == domMap->getLocalNumElements()) {
152 Kokkos::View<GO*, memory_space> myColumns(Kokkos::ViewAllocateWithoutInitializing(
"myColumns"), numLocalColGIDs + numRemoteColGIDs);
154 auto LocalColGIDs = Kokkos::subview(myColumns, Kokkos::make_pair((
size_t)0, numLocalColGIDs));
155 auto remoteColGIDs = Kokkos::subview(myColumns, Kokkos::make_pair(numLocalColGIDs, numLocalColGIDs + numRemoteColGIDs));
158 Kokkos::deep_copy(remoteColGIDs, RemoteGIDs);
163 if (
static_cast<size_t>(remotePIDs.size()) != numRemoteColGIDs) {
164 remotePIDs.resize(numRemoteColGIDs);
169 auto remoteColGIDs_h = Kokkos::create_mirror_view(remoteColGIDs);
170 Kokkos::deep_copy(remoteColGIDs_h, remoteColGIDs);
171 auto remoteColGIDs_av = Kokkos::Compat::getArrayView(remoteColGIDs_h);
173 domMap->getRemoteIndexList(remoteColGIDs_av, remotePIDs());
183 if (errStrm != NULL) {
184 *errStrm << prefix <<
"Some column indices are not in the domain Map."
185 "Either these column indices are invalid or the domain Map is "
186 "invalid. Likely cause: For a nonsquare matrix, you must give the "
187 "domain and range Maps as input to fillComplete."
208 Kokkos::View<int*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> remotePIDs_h(remotePIDs.data(), remotePIDs.size());
210 auto remotePIDs_d = Kokkos::create_mirror_view_and_copy(execution_space(), remotePIDs_h);
211 Kokkos::View<GO*, memory_space> remoteColGIDsSorted(Kokkos::ViewAllocateWithoutInitializing(
"remoteColGIDsSorted"), numRemoteColGIDs);
212 Kokkos::View<int*, memory_space> remotePIDsSorted(Kokkos::ViewAllocateWithoutInitializing(
"remotePIDsSorted"), numRemoteColGIDs);
213 sort_by_key(execution_space(), remotePIDs_d, remoteColGIDs, remotePIDsSorted, remoteColGIDsSorted);
214 Kokkos::deep_copy(remoteColGIDs, remoteColGIDsSorted);
215 Kokkos::deep_copy(remotePIDs_h, remotePIDsSorted);
228 const size_t numDomainElts = domMap->getLocalNumElements();
229 if (numLocalColGIDs == numDomainElts) {
233 if (domMap->isContiguous()) {
238 auto domMapMinGlobalIndex = domMap->getMinGlobalIndex();
239 Kokkos::parallel_for(
240 Kokkos::RangePolicy<execution_space>(0, numLocalColGIDs), KOKKOS_LAMBDA(
const size_t k) { LocalColGIDs(k) = domMapMinGlobalIndex + k; });
242 Kokkos::deep_copy(LocalColGIDs, domMap->getMyGlobalIndicesDevice());
247 size_t numLocalCount = 0;
248 if (domMap->isContiguous()) {
253 auto domMapMinGlobalIndex = domMap->getMinGlobalIndex();
254 Kokkos::parallel_scan(
255 Kokkos::RangePolicy<execution_space>(0, numDomainElts), KOKKOS_LAMBDA(
const size_t i,
size_t& offset,
const bool is_final) {
258 LocalColGIDs(offset) = domMapMinGlobalIndex + i;
262 auto domainElts = domMap->getMyGlobalIndicesDevice();
263 Kokkos::parallel_scan(
264 Kokkos::RangePolicy<execution_space>(0, numDomainElts), KOKKOS_LAMBDA(
const size_t i,
size_t& offset,
const bool is_final) {
267 LocalColGIDs(offset) = domainElts(i);
271 if (numLocalCount != numLocalColGIDs) {
272 if (errStrm != NULL) {
273 *errStrm << prefix <<
"numLocalCount = " << numLocalCount
274 <<
" != numLocalColGIDs = " << numLocalColGIDs
275 <<
". This should never happen. "
276 "Please report this bug to the Tpetra developers."
333 Tpetra::Details::OrdinalTraits<global_size_t>::invalid();
338 const GO indexBase = domMap->getIndexBase();
339 colMap = rcp(
new map_type(INV, myColumns, indexBase, domMap->getComm()));
343template <
class ValueTypeView,
class ValuesIdxType>
345 using value_type =
typename ValueTypeView::non_const_value_type;
348 void op(ValueTypeView values, ValuesIdxType values_idx,
349 const value_type v)
const {
350 Kokkos::atomic_min(values.data() + values_idx, v);
354template <
class LO,
class GO,
class NT>
356 Teuchos::Array<int>& remotePIDs,
363 using execution_space =
typename NT::execution_space;
364 using memory_space =
typename NT::memory_space;
366 const char prefix[] =
"Tpetra::Details::makeColMap: ";
377 colMap = Teuchos::null;
381 Kokkos::View<GO*, memory_space>
myColumns;
382 if (
graph.isLocallyIndexed()) {
383 colMap =
graph.getColMap();
386 if (colMap.is_null() || !
graph.hasColMap()) {
389 *
errStrm <<
prefix <<
"The graph is locally indexed on the calling "
390 "process, but has no column Map (either getColMap() returns null, "
391 "or hasColMap() returns false)."
403 if (colMap->isContiguous()) {
405 const LO
numCurGids =
static_cast<LO
>(colMap->getLocalNumElements());
406 myColumns = Kokkos::View<GO*, memory_space>(Kokkos::ViewAllocateWithoutInitializing(
"myColumns"),
numCurGids);
408 Kokkos::parallel_for(
411 auto curGids =
graph.getColMap()->getMyGlobalIndicesDevice();
414 myColumns = Kokkos::View<GO*, memory_space>(Kokkos::ViewAllocateWithoutInitializing(
"myColumns"),
numCurGids);
418 }
else if (
graph.isGloballyIndexed()) {
442 const LO
LINV = Tpetra::Details::OrdinalTraits<LO>::invalid();
443 size_t numLocalColGIDs = 0;
444 size_t numRemoteColGIDs = 0;
448 Kokkos::View<bool*, memory_space>
GIDisLocal(
"GIDisLocal",
domMap->getLocalNumElements());
450 if (!
graph.getRowMap().is_null()) {
452 const LO
lclNumRows = rowMap.getLocalNumElements();
458 auto numRowEntries = Kokkos::create_mirror_view_and_copy(execution_space(),
graph.k_numRowEntries_);
460 using unordered_map_type = Kokkos::UnorderedMap<GO, size_t, typename NT::device_type>;
464 if (
domMap->getComm()->getSize() == 1) {
467 Kokkos::parallel_reduce(
488 Kokkos::parallel_reduce(
502 } }, numLocalColGIDs);
512 RemoteGIDs = Kokkos::View<GO*, memory_space>(Kokkos::ViewAllocateWithoutInitializing(
"RemoteGIDs"), numRemoteColGIDs);
514 auto RemoteGIDsOffsets = Kokkos::View<LO*, memory_space>(Kokkos::ViewAllocateWithoutInitializing(
"RemoteGIDsOffsets"), numRemoteColGIDs);
516 Kokkos::parallel_scan(
520 RemoteGIDs(offset) = RemoteGIDMap.key_at(i);
521 RemoteGIDsOffsets(offset) = RemoteGIDMap.value_at(i);
529 Kokkos::parallel_scan(
533 RemoteGIDs(offset) = RemoteGIDMap.key_at(i);
547 numLocalColGIDs, numRemoteColGIDs,
561 Tpetra::Details::OrdinalTraits<global_size_t>::invalid();
566 const GO indexBase = domMap->getIndexBase();
567 colMap = rcp(
new map_type(INV, myColumns, indexBase, domMap->getComm()));
571template <
typename GOView,
typename bitset_t>
572struct GatherPresentEntries {
573 using GO =
typename GOView::non_const_value_type;
575 GatherPresentEntries(GO minGID_,
const GOView& gids_,
const bitset_t& present_)
578 , present(present_) {}
580 KOKKOS_INLINE_FUNCTION
void operator()(
const GO i)
const {
581 present.set(gids(i) - minGID);
589template <
typename GO,
typename mem_space>
590struct MinMaxReduceFunctor {
591 using MinMaxValue =
typename Kokkos::MinMax<GO>::value_type;
592 using GOView = Kokkos::View<GO*, mem_space>;
594 MinMaxReduceFunctor(
const GOView& gids_)
597 KOKKOS_INLINE_FUNCTION
void operator()(
const GO i, MinMaxValue& lminmax)
const {
599 if (gid < lminmax.min_val)
600 lminmax.min_val = gid;
601 if (gid > lminmax.max_val)
602 lminmax.max_val = gid;
608template <
class LO,
class GO,
class NT>
611 Kokkos::View<GO*, typename NT::memory_space> gids,
613 using Kokkos::RangePolicy;
614 using Teuchos::Array;
615 using device_t =
typename NT::device_type;
616 using exec_space =
typename device_t::execution_space;
617 using memory_space =
typename device_t::memory_space;
623 using bitset_t = Kokkos::Bitset<typename exec_space::memory_space>;
624 using const_bitset_t = Kokkos::ConstBitset<typename exec_space::memory_space>;
625 using GOView = Kokkos::View<GO*, memory_space>;
627 GO minGID = Teuchos::OrdinalTraits<GO>::max();
629 using MinMaxValue =
typename Kokkos::MinMax<GO>::value_type;
642 Kokkos::View<bool*, memory_space>
GIDisLocal(
"GIDisLocal",
domMap->getLocalNumElements());
645 GO numLocalColGIDs = 0;
646 const auto LINV = ::Tpetra::Details::OrdinalTraits<LO>::invalid();
647 Kokkos::parallel_reduce(
650 auto gid = i + minGID;
651 auto lid = localDomMap.getLocalElement(gid);
652 bool isRemote = (lid == LINV);
654 GIDisLocal(lid) = true;
657 } }, numLocalColGIDs);
660 auto RemoteGIDs = Kokkos::View<GO*, memory_space>(
"RemoteGIDs", numRemoteColGIDs);
662 Kokkos::parallel_scan(
666 bool isRemote = (localDomMap.getLocalElement(gid) == LINV);
669 RemoteGIDs(lcount) = gid;
676 Array<int> remotePIDs;
678 return makeColMapImpl<LO, GO, NT>(
682 static_cast<size_t>(numLocalColGIDs),
683 static_cast<size_t>(numRemoteColGIDs),
696#define TPETRA_DETAILS_MAKECOLMAP_INSTANT(LO, GO, NT) \
697 namespace Details { \
699 makeColMap(Teuchos::RCP<const Tpetra::Map<LO, GO, NT>>&, \
700 Teuchos::Array<int>&, \
701 const Teuchos::RCP<const Tpetra::Map<LO, GO, NT>>&, \
702 const CrsGraph<LO, GO, NT>&, \
706 makeColMap(Teuchos::RCP<const Tpetra::Map<LO, GO, NT>>&, \
707 const Teuchos::RCP<const Tpetra::Map<LO, GO, NT>>&, \
708 Kokkos::View<GO*, typename NT::memory_space>, \
Stand-alone utility functions and macros.
Struct that holds views of the contents of a CrsMatrix.
Nonmember function that computes a residual Computes R = B - A * X.
int makeColMap(Teuchos::RCP< const Tpetra::Map< LO, GO, NT > > &colMap, Teuchos::Array< int > &remotePIDs, const Teuchos::RCP< const Tpetra::Map< LO, GO, NT > > &domMap, const CrsGraph< LO, GO, NT > &graph, const bool sortEachProcsGids=true, std::ostream *errStrm=NULL)
Make the graph's column Map.
LookupStatus
Return status of Map remote index lookup (getRemoteIndexList()).
size_t global_size_t
Global size_t object.