Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Map_def.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
14
15#ifndef TPETRA_MAP_DEF_HPP
16#define TPETRA_MAP_DEF_HPP
17
18#include <memory>
19#include <sstream>
20#include <stdexcept>
21#include <typeinfo>
22
23#include "Teuchos_as.hpp"
24#include "Teuchos_TypeNameTraits.hpp"
25#include "Teuchos_CommHelpers.hpp"
26
27#include "Kokkos_Sort.hpp"
28
29#include "Tpetra_Directory.hpp" // must include for implicit instantiation to work
32#include "Tpetra_Details_FixedHashTable.hpp"
35#include "Tpetra_Core.hpp"
36#include "Tpetra_Map_decl.hpp"
37#include "Tpetra_Util.hpp"
38#include "Tpetra_Details_mpiIsInitialized.hpp"
39#include "Tpetra_Details_extractMpiCommFromTeuchos.hpp" // teuchosCommIsAnMpiComm
42
43namespace Tpetra {
44namespace Impl {
45
46inline void checkMapInputArray(const char ctorName[],
47 const void* indexList,
48 const size_t indexListSize,
49 const Teuchos::Comm<int>* const comm) {
51
52 const bool debug = Behavior::debug("Map");
53 if (debug) {
54 using std::endl;
55 using Teuchos::outArg;
56 using Teuchos::REDUCE_MIN;
57 using Teuchos::reduceAll;
58
59 const int myRank = comm == nullptr ? 0 : comm->getRank();
60 const bool verbose = Behavior::verbose("Map");
61 std::ostringstream lclErrStrm;
62 int lclSuccess = 1;
63
64 if (indexListSize != 0 && indexList == nullptr) {
65 lclSuccess = 0;
66 if (verbose) {
67 lclErrStrm << "Proc " << myRank << ": indexList is null, "
68 "but indexListSize="
69 << indexListSize << " != 0." << endl;
70 }
71 }
72 int gblSuccess = 0; // output argument
73 reduceAll(*comm, REDUCE_MIN, lclSuccess, outArg(gblSuccess));
74 if (gblSuccess != 1) {
75 std::ostringstream gblErrStrm;
76 gblErrStrm << "Tpetra::Map constructor " << ctorName << " detected a problem with the input array "
77 "(raw array, Teuchos::ArrayView, or Kokkos::View) "
78 "of global indices."
79 << endl;
80 if (verbose) {
81 using ::Tpetra::Details::gathervPrint;
82 gathervPrint(gblErrStrm, lclErrStrm.str(), *comm);
83 }
84 TEUCHOS_TEST_FOR_EXCEPTION(true, std::invalid_argument, gblErrStrm.str());
85 }
86 }
87}
88
89template <class LocalOrdinal, class GlobalOrdinal, class ViewType>
90void computeConstantsOnDevice(const ViewType& entryList, GlobalOrdinal& minMyGID, GlobalOrdinal& maxMyGID, GlobalOrdinal& firstContiguousGID, GlobalOrdinal& lastContiguousGID_val, LocalOrdinal& lastContiguousGID_loc) {
91 using LO = LocalOrdinal;
92 using GO = GlobalOrdinal;
93 using exec_space = typename ViewType::device_type::execution_space;
94 using range_policy = Kokkos::RangePolicy<exec_space, Kokkos::IndexType<LO>>;
95 const LO numLocalElements = entryList.extent(0);
96
97 // We're going to use the minloc backwards because we need to have it sort on the "location" and have the "value" along for the
98 // ride, rather than the other way around
99 typedef typename Kokkos::MinLoc<LO, GO>::value_type minloc_type;
100 minloc_type myMinLoc;
101
102 // Find the initial sequence of parallel gids
103 // To find the lastContiguousGID_, we find the first guy where entryList[i] - entryList[0] != i-0. That's the first non-contiguous guy.
104 // We want the one *before* that guy.
105 Kokkos::parallel_reduce(
106 range_policy(0, numLocalElements), KOKKOS_LAMBDA(const LO& i, GO& l_myMin, GO& l_myMax, GO& l_firstCont, minloc_type& l_lastCont) {
107 GO entry_0 = entryList[0];
108 GO entry_i = entryList[i];
109
110 // Easy stuff
111 l_myMin = (l_myMin < entry_i) ? l_myMin : entry_i;
112 l_myMax = (l_myMax > entry_i) ? l_myMax : entry_i;
113 l_firstCont = entry_0;
114
115 if (entry_i - entry_0 != i && l_lastCont.val >= i) {
116 // We're non-contiguous, so the guy before us could be the last contiguous guy
117 l_lastCont.val = i - 1;
118 l_lastCont.loc = entryList[i - 1];
119 } else if (i == numLocalElements - 1 && i < l_lastCont.val) {
120 // If we're last, we always think we're the last contiguous guy, unless someone non-contiguous is already here
121 l_lastCont.val = i;
122 l_lastCont.loc = entry_i;
123 }
124 },
125 Kokkos::Min<GO>(minMyGID), Kokkos::Max<GO>(maxMyGID), Kokkos::Min<GO>(firstContiguousGID), Kokkos::MinLoc<LO, GO>(myMinLoc));
126
127 // This switch is intentional, since we're using MinLoc backwards
128 lastContiguousGID_val = myMinLoc.loc;
129 lastContiguousGID_loc = myMinLoc.val;
130}
131
132} // namespace Impl
133
134template <class LocalOrdinal, class GlobalOrdinal, class Node>
136 Map()
137 : comm_(new Teuchos::SerialComm<int>())
138 , indexBase_(0)
139 , numGlobalElements_(0)
140 , numLocalElements_(0)
145 , firstContiguousGID_(Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid())
146 , lastContiguousGID_(Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid())
147 , uniform_(false)
148 , // trivially
149 contiguous_(false)
150 , distributed_(false)
151 , haveGlobalConstants_(true)
152 , // no communicator yet
153 directory_(new Directory<LocalOrdinal, GlobalOrdinal, Node>()) {
156}
157
158template <class LocalOrdinal, class GlobalOrdinal, class Node>
162 const Teuchos::RCP<const Teuchos::Comm<int>>& comm,
163 const LocalGlobal lOrG)
164 : comm_(comm)
165 , uniform_(true)
166 , directory_(new Directory<LocalOrdinal, GlobalOrdinal, Node>()) {
167 using std::endl;
168 using Teuchos::as;
169 using Teuchos::broadcast;
170 using Teuchos::outArg;
171 using Teuchos::REDUCE_MAX;
172 using Teuchos::REDUCE_MIN;
173 using Teuchos::reduceAll;
174 using Teuchos::typeName;
175 using GO = global_ordinal_type;
176 using GST = global_size_t;
177 const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid();
178 const char funcName[] = "Map(gblNumInds,indexBase,comm,LG)";
179 const char exPfx[] =
180 "Tpetra::Map::Map(gblNumInds,indexBase,comm,LG): ";
181
182 const bool debug = Details::Behavior::debug("Map");
183 const bool verbose = Details::Behavior::verbose("Map");
184 std::unique_ptr<std::string> prefix;
185 if (verbose) {
187 comm_.getRawPtr(), "Map", funcName);
188 std::ostringstream os;
189 os << *prefix << "Start" << endl;
190 std::cerr << os.str();
191 }
194
195 // In debug mode only, check whether numGlobalElements and
196 // indexBase are the same over all processes in the communicator.
197 if (debug) {
208 std::invalid_argument, exPfx << "All processes must "
209 "provide the same number of global elements. Process 0 set "
210 "numGlobalElements="
211 << proc0NumGlobalElements << ". The "
212 "calling process "
213 << comm->getRank() << " set "
214 "numGlobalElements="
215 << numGlobalElements << ". The min "
216 "and max values over all processes are "
217 << minNumGlobalElements << " resp. " << maxNumGlobalElements << ".");
218
226 std::invalid_argument, exPfx << "All processes must "
227 "provide the same indexBase argument. Process 0 set "
228 "indexBase="
229 << proc0IndexBase << ". The calling process " << comm->getRank() << " set indexBase=" << indexBase << ". The min and max values over all processes are " << minIndexBase << " resp. " << maxIndexBase << ".");
230 }
231
232 // Distribute the elements across the processes in the given
233 // communicator so that global IDs (GIDs) are
234 //
235 // - Nonoverlapping (only one process owns each GID)
236 // - Contiguous (the sequence of GIDs is nondecreasing, and no two
237 // adjacent GIDs differ by more than one)
238 // - As evenly distributed as possible (the numbers of GIDs on two
239 // different processes do not differ by more than one)
240
241 // All processes have the same numGlobalElements, but we still
242 // need to check that it is valid. numGlobalElements must be
243 // positive and not the "invalid" value (GSTI).
244 //
245 // This comparison looks funny, but it avoids compiler warnings
246 // for comparing unsigned integers (numGlobalElements_in is a
247 // GST, which is unsigned) while still working if we
248 // later decide to make GST signed.
251 std::invalid_argument, exPfx << "numGlobalElements (= " << numGlobalElements << ") must be nonnegative.");
252
253 TEUCHOS_TEST_FOR_EXCEPTION(numGlobalElements == GSTI, std::invalid_argument, exPfx << "You provided numGlobalElements = Teuchos::OrdinalTraits<"
254 "Tpetra::global_size_t>::invalid(). This version of the "
255 "constructor requires a valid value of numGlobalElements. "
256 "You probably mistook this constructor for the \"contiguous "
257 "nonuniform\" constructor, which can compute the global "
258 "number of elements for you if you set numGlobalElements to "
259 "Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid().");
260
261 size_t numLocalElements = 0; // will set below
262 if (lOrG == GloballyDistributed) {
263 // Compute numLocalElements:
264 //
265 // If numGlobalElements == numProcs * B + remainder,
266 // then Proc r gets B+1 elements if r < remainder,
267 // and B elements if r >= remainder.
268 //
269 // This strategy is valid for any value of numGlobalElements and
270 // numProcs, including the following border cases:
271 // - numProcs == 1
272 // - numLocalElements < numProcs
273 //
274 // In the former case, remainder == 0 && numGlobalElements ==
275 // numLocalElements. In the latter case, remainder ==
276 // numGlobalElements && numLocalElements is either 0 or 1.
277 const GST numProcs = static_cast<GST>(comm_->getSize());
278 const GST myRank = static_cast<GST>(comm_->getRank());
281
282 GO startIndex;
283 if (myRank < remainder) {
284 numLocalElements = static_cast<size_t>(1) + static_cast<size_t>(quotient);
285 // myRank was originally an int, so it should never overflow
286 // reasonable GO types.
288 } else {
292 }
293
294 minMyGID_ = indexBase + startIndex;
295 maxMyGID_ = indexBase + startIndex + numLocalElements - 1;
296 minAllGID_ = indexBase;
297 maxAllGID_ = indexBase + numGlobalElements - 1;
298 distributed_ = (numProcs > 1);
299 } else { // lOrG == LocallyReplicated
301 minMyGID_ = indexBase;
302 maxMyGID_ = indexBase + numGlobalElements - 1;
303 distributed_ = false;
304 }
305
306 minAllGID_ = indexBase;
307 maxAllGID_ = indexBase + numGlobalElements - 1;
308 indexBase_ = indexBase;
309 numGlobalElements_ = numGlobalElements;
310 numLocalElements_ = numLocalElements;
311 firstContiguousGID_ = minMyGID_;
312 lastContiguousGID_ = maxMyGID_;
313 contiguous_ = true;
314 haveGlobalConstants_ = true;
315
316 // Create the Directory on demand in getRemoteIndexList().
317 // setupDirectory ();
318
319 if (verbose) {
320 std::ostringstream os;
321 os << *prefix << "Done" << endl;
322 std::cerr << os.str();
323 }
324}
325
326template <class LocalOrdinal, class GlobalOrdinal, class Node>
329 const size_t numLocalElements,
331 const Teuchos::RCP<const Teuchos::Comm<int>>& comm)
332 : comm_(comm)
333 , uniform_(false)
334 , directory_(new Directory<LocalOrdinal, GlobalOrdinal, Node>()) {
335 using std::endl;
336 using Teuchos::as;
337 using Teuchos::broadcast;
338 using Teuchos::outArg;
339 using Teuchos::REDUCE_MAX;
340 using Teuchos::REDUCE_MIN;
341 using Teuchos::REDUCE_SUM;
342 using Teuchos::reduceAll;
343 using Teuchos::scan;
344 using GO = global_ordinal_type;
345 using GST = global_size_t;
346 const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid();
347 const char funcName[] =
348 "Map(gblNumInds,lclNumInds,indexBase,comm)";
349 const char exPfx[] =
350 "Tpetra::Map::Map(gblNumInds,lclNumInds,indexBase,comm): ";
351 const char suffix[] =
352 ". Please report this bug to the Tpetra developers.";
353
354 const bool debug = Details::Behavior::debug("Map");
355 const bool verbose = Details::Behavior::verbose("Map");
356 std::unique_ptr<std::string> prefix;
357 if (verbose) {
359 comm_.getRawPtr(), "Map", funcName);
360 std::ostringstream os;
361 os << *prefix << "Start" << endl;
362 std::cerr << os.str();
363 }
366
367 // Global sum of numLocalElements over all processes.
368 // Keep this for later debug checks.
370 if (debug) {
371 debugGlobalSum = initialNonuniformDebugCheck(exPfx,
373 }
374
375 // Distribute the elements across the nodes so that they are
376 // - non-overlapping
377 // - contiguous
378
379 // This differs from the first Map constructor (that only takes a
380 // global number of elements) in that the user has specified the
381 // number of local elements, so that the elements are not
382 // (necessarily) evenly distributed over the processes.
383
384 // Compute my local offset. This is an inclusive scan, so to get
385 // the final offset, we subtract off the input.
386 GO scanResult = 0;
389
390 if (numGlobalElements != GSTI) {
391 numGlobalElements_ = numGlobalElements; // Use the user's value.
392 } else {
393 // Inclusive scan means that the last process has the final sum.
394 // Rather than doing a reduceAll to get the sum of
395 // numLocalElements, we can just have the last process broadcast
396 // its result. That saves us a round of log(numProcs) messages.
397 const int numProcs = comm->getSize();
399 if (numProcs > 1) {
400 broadcast(*comm, numProcs - 1, outArg(globalSum));
401 }
402 numGlobalElements_ = globalSum;
403
404 if (debug) {
405 // No need for an all-reduce here; both come from collectives.
406 TEUCHOS_TEST_FOR_EXCEPTION(globalSum != debugGlobalSum, std::logic_error, exPfx << "globalSum = " << globalSum << " != debugGlobalSum = " << debugGlobalSum << suffix);
407 }
408 }
409 numLocalElements_ = numLocalElements;
410 indexBase_ = indexBase;
411 minAllGID_ = (numGlobalElements_ == 0) ? std::numeric_limits<GO>::max() : indexBase;
412 maxAllGID_ = (numGlobalElements_ == 0) ? std::numeric_limits<GO>::lowest() : indexBase + GO(numGlobalElements_) - GO(1);
413 minMyGID_ = (numLocalElements_ == 0) ? std::numeric_limits<GO>::max() : indexBase + GO(myOffset);
414 maxMyGID_ = (numLocalElements_ == 0) ? std::numeric_limits<GO>::lowest() : indexBase + myOffset + GO(numLocalElements) - GO(1);
415 firstContiguousGID_ = minMyGID_;
416 lastContiguousGID_ = maxMyGID_;
417 contiguous_ = true;
418 distributed_ = (comm->getSize() > 1);
419 haveGlobalConstants_ = true;
420
421 // Create the Directory on demand in getRemoteIndexList().
422 // setupDirectory ();
423
424 if (verbose) {
425 std::ostringstream os;
426 os << *prefix << "Done" << endl;
427 std::cerr << os.str();
428 }
429}
430
431template <class LocalOrdinal, class GlobalOrdinal, class Node>
433Map<LocalOrdinal, GlobalOrdinal, Node>::
434 initialNonuniformDebugCheck(
435 const char errorMessagePrefix[],
436 const global_size_t numGlobalElements,
437 const size_t numLocalElements,
438 const global_ordinal_type indexBase,
439 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) const {
440 const bool debug = Details::Behavior::debug("Map");
441 if (!debug) {
442 return global_size_t(0);
443 }
444
445 using Teuchos::broadcast;
446 using Teuchos::outArg;
447 using Teuchos::ptr;
448 using Teuchos::REDUCE_MAX;
449 using Teuchos::REDUCE_MIN;
450 using Teuchos::REDUCE_SUM;
451 using Teuchos::reduceAll;
452 using GO = global_ordinal_type;
453 using GST = global_size_t;
454 const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid();
455
456 // The user has specified the distribution of indices over the
457 // processes. The distribution is not necessarily contiguous or
458 // equally shared over the processes.
459 //
460 // We assume that the number of local elements can be stored in a
461 // size_t. The instance member numLocalElements_ is a size_t, so
462 // this variable and that should have the same type.
463
464 GST debugGlobalSum = 0; // Will be global sum of numLocalElements
467 // In debug mode only, check whether numGlobalElements and
468 // indexBase are the same over all processes in the communicator.
469 {
480 std::invalid_argument, errorMessagePrefix << "All processes "
481 "must provide the same number of global elements, even if "
482 "that argument is "
483 "Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid() "
484 "(which signals that the Map should compute the global "
485 "number of elements). Process 0 set numGlobalElements"
486 "="
487 << proc0NumGlobalElements << ". The calling process " << comm->getRank() << " set numGlobalElements=" << numGlobalElements << ". The min and max values over all "
488 "processes are "
489 << minNumGlobalElements << " resp. " << maxNumGlobalElements << ".");
490
498 std::invalid_argument, errorMessagePrefix << "All processes must provide the same indexBase argument. "
499 "Process 0 set indexBase = "
500 << proc0IndexBase << ". The "
501 "calling process "
502 << comm->getRank() << " set indexBase=" << indexBase << ". The min and max values over all "
503 "processes are "
504 << minIndexBase << " resp. " << maxIndexBase << ".");
505
506 // Make sure that the sum of numLocalElements over all processes
507 // equals numGlobalElements.
510 std::invalid_argument,
511 errorMessagePrefix << "The sum of each process' number of "
512 "indices over all processes, "
513 << debugGlobalSum << ", != "
514 << "numGlobalElements=" << numGlobalElements << ". If you "
515 "would like this constructor to compute numGlobalElements "
516 "for you, you may set numGlobalElements="
517 "Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid() "
518 "on input. Please note that this is NOT necessarily -1.");
519 }
520 return debugGlobalSum;
521}
522
523template <class LocalOrdinal, class GlobalOrdinal, class Node>
524void Map<LocalOrdinal, GlobalOrdinal, Node>::
525 initWithNonownedHostIndexList(
526 const char errorMessagePrefix[],
527 const global_size_t numGlobalElements,
528 const Kokkos::View<const global_ordinal_type*,
529 Kokkos::LayoutLeft,
530 Kokkos::HostSpace,
531 Kokkos::MemoryUnmanaged>& entryList_host,
532 const global_ordinal_type indexBase,
533 const Teuchos::RCP<const Teuchos::Comm<int>>& comm,
534 const Teuchos::RCP<Teuchos::ParameterList>& params) {
535 Tpetra::Details::ProfilingRegion pr("Map::initWithNonownedHostIndexList()");
536
537 using Kokkos::LayoutLeft;
538 using Kokkos::subview;
539 using Kokkos::View;
540 using Kokkos::view_alloc;
541 using Kokkos::WithoutInitializing;
542 using Teuchos::as;
543 using Teuchos::broadcast;
544 using Teuchos::outArg;
545 using Teuchos::ptr;
546 using Teuchos::REDUCE_MAX;
547 using Teuchos::REDUCE_MIN;
548 using Teuchos::REDUCE_SUM;
549 using Teuchos::reduceAll;
550 using LO = local_ordinal_type;
551 using GO = global_ordinal_type;
553 const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid();
554 // Make sure that Kokkos has been initialized (Github Issue #513).
555 TEUCHOS_TEST_FOR_EXCEPTION(!Kokkos::is_initialized(), std::runtime_error,
556 errorMessagePrefix << "The Kokkos execution space "
557 << Teuchos::TypeNameTraits<execution_space>::name()
558 << " has not been initialized. "
559 "Please initialize it before creating a Map.")
560
561 // The user has specified the distribution of indices over the
562 // processes, via the input array of global indices on each
563 // process. The distribution is not necessarily contiguous or
564 // equally shared over the processes.
565
566 // The length of the input array on this process is the number of
567 // local indices to associate with this process, even though the
568 // input array contains global indices. We assume that the number
569 // of local indices on a process can be stored in a size_t;
570 // numLocalElements_ is a size_t, so this variable and that should
571 // have the same type.
573
574 initialNonuniformDebugCheck(errorMessagePrefix, numGlobalElements,
576
577 // NOTE (mfh 20 Feb 2013, 10 Oct 2016) In some sense, this global
578 // reduction is redundant, since the directory Map will have to do
579 // the same thing. Thus, we could do the scan and broadcast for
580 // the directory Map here, and give the computed offsets to the
581 // directory Map's constructor. However, a reduction costs less
582 // than a scan and broadcast, so this still saves time if users of
583 // this Map don't ever need the Directory (i.e., if they never
584 // call getRemoteIndexList on this Map).
585 std::shared_ptr<Details::CommRequest> req;
588 numGlobalElements_ = numGlobalElements; // Use the user's value.
589 } else { // The user wants us to compute the sum.
591 req = Details::iallreduce(numLocalElementsGST,
592 numGlobalElements_, REDUCE_SUM, *comm);
593 }
594
595 // mfh 20 Feb 2013: We've never quite done the right thing for
596 // duplicate GIDs here. Duplicate GIDs have always been counted
597 // distinctly in numLocalElements_, and thus should get a
598 // different LID. However, we've always used std::map or a hash
599 // table for the GID -> LID lookup table, so distinct GIDs always
600 // map to the same LID. Furthermore, the order of the input GID
601 // list matters, so it's not desirable to sort for determining
602 // uniqueness.
603 //
604 // I've chosen for now to write this code as if the input GID list
605 // contains no duplicates. If this is not desired, we could use
606 // the lookup table itself to determine uniqueness: If we haven't
607 // seen the GID before, it gets a new LID and it's added to the
608 // LID -> GID and GID -> LID tables. If we have seen the GID
609 // before, it doesn't get added to either table. I would
610 // implement this, but it would cost more to do the double lookups
611 // in the table (one to check, and one to insert).
612 //
613 // More importantly, since we build the GID -> LID table in (a
614 // thread-) parallel (way), the order in which duplicate GIDs may
615 // get inserted is not defined. This would make the assignment of
616 // LID to GID nondeterministic.
617
618 numLocalElements_ = numLocalElements;
619 indexBase_ = indexBase;
620
621 minMyGID_ = indexBase_;
622 maxMyGID_ = indexBase_;
623
624 // NOTE (mfh 27 May 2015): While finding the initial contiguous
625 // GID range requires looking at all the GIDs in the range,
626 // dismissing an interval of GIDs only requires looking at the
627 // first and last GIDs. Thus, we could do binary search backwards
628 // from the end in order to catch the common case of a contiguous
629 // interval followed by noncontiguous entries. On the other hand,
630 // we could just expose this case explicitly as yet another Map
631 // constructor, and avoid the trouble of detecting it.
632 if (numLocalElements_ > 0) {
633 Tpetra::Details::ProfilingRegion prLcl("Map::initWithNonownedHostIndexList::local");
634 // Find contiguous GID range, with the restriction that the
635 // beginning of the range starts with the first entry. While
636 // doing so, fill in the LID -> GID table.
637 typename decltype(lgMap_)::non_const_type lgMap(view_alloc("lgMap", WithoutInitializing), numLocalElements_);
638 auto lgMap_host =
639 Kokkos::create_mirror_view(Kokkos::HostSpace(), lgMap);
640
641 // The input array entryList_host is already on host, so we
642 // don't need to take a host view of it.
643 // auto entryList_host =
644 // Kokkos::create_mirror_view (Kokkos::HostSpace (), entryList);
645 // Kokkos::deep_copy (entryList_host, entryList);
646
647 firstContiguousGID_ = entryList_host[0];
648 lastContiguousGID_ = firstContiguousGID_ + 1;
649
650 // FIXME (mfh 23 Sep 2015) We need to copy the input GIDs
651 // anyway, so we have to look at them all. The logical way to
652 // find the first noncontiguous entry would thus be to "reduce,"
653 // where the local reduction result is whether entryList[i] + 1
654 // == entryList[i+1].
655
656 lgMap_host[0] = firstContiguousGID_;
657 size_t i = 1;
658 for (; i < numLocalElements_; ++i) {
659 const GO curGid = entryList_host[i];
660 const LO curLid = as<LO>(i);
661
662 if (lastContiguousGID_ != curGid) break;
663
664 // Add the entry to the LID->GID table only after we know that
665 // the current GID is in the initial contiguous sequence, so
666 // that we don't repeat adding it in the first iteration of
667 // the loop below over the remaining noncontiguous GIDs.
669 ++lastContiguousGID_;
670 }
671 --lastContiguousGID_;
672 // NOTE: i is the first non-contiguous index.
673
674 // [firstContiguousGID_, lastContigousGID_] is the initial
675 // sequence of contiguous GIDs. We can start the min and max
676 // GID using this range.
677 minMyGID_ = firstContiguousGID_;
678 maxMyGID_ = lastContiguousGID_;
679
680 // Compute the GID -> LID lookup table, _not_ including the
681 // initial sequence of contiguous GIDs.
683 {
684 const std::pair<size_t, size_t> ncRange(i, entryList_host.extent(0));
686 TEUCHOS_TEST_FOR_EXCEPTION(static_cast<size_t>(nonContigGids_host.extent(0)) !=
687 static_cast<size_t>(entryList_host.extent(0) - i),
688 std::logic_error,
689 "Tpetra::Map noncontiguous constructor: "
690 "nonContigGids_host.extent(0) = "
691 << nonContigGids_host.extent(0)
692 << " != entryList_host.extent(0) - i = "
693 << (entryList_host.extent(0) - i) << " = "
694 << entryList_host.extent(0) << " - " << i
695 << ". Please report this bug to the Tpetra developers.");
696
697 // FixedHashTable's constructor expects an owned device View,
698 // so we must deep-copy the subview of the input indices.
701 nonContigGids_host.size());
702
703 // DEEP_COPY REVIEW - HOST-TO-DEVICE
704 Kokkos::deep_copy(execution_space(), nonContigGids, nonContigGids_host);
705 Kokkos::fence("Map::initWithNonownedHostIndexList"); // for UVM issues below - which will be refatored soon so FixedHashTable can build as pure CudaSpace - then I think remove this fence
706
709 // Make host version - when memory spaces match these just do trivial assignment
710 glMapHost_ = global_to_local_table_host_type(glMap_);
711 }
712
713 // FIXME (mfh 10 Oct 2016) When we construct the global-to-local
714 // table above, we have to look at all the (noncontiguous) input
715 // indices anyway. Thus, why not have the constructor compute
716 // and return the min and max?
717
718 for (; i < numLocalElements_; ++i) {
719 const GO curGid = entryList_host[i];
720 const LO curLid = static_cast<LO>(i);
721 lgMap_host[curLid] = curGid; // LID -> GID table
722
723 // While iterating through entryList, we compute its
724 // (process-local) min and max elements.
725 if (curGid < minMyGID_) {
726 minMyGID_ = curGid;
727 }
728 if (curGid > maxMyGID_) {
729 maxMyGID_ = curGid;
730 }
731 }
732
733 // We filled lgMap on host above; now sync back to device.
734 // DEEP_COPY REVIEW - HOST-TO-DEVICE
735 Kokkos::deep_copy(execution_space(), lgMap, lgMap_host);
736
737 // "Commit" the local-to-global lookup table we filled in above.
738 lgMap_ = lgMap;
739 // We've already created this, so use it.
740 lgMapHost_ = lgMap_host;
741
742 } else {
743 minMyGID_ = std::numeric_limits<GlobalOrdinal>::max();
744 maxMyGID_ = std::numeric_limits<GlobalOrdinal>::lowest();
745 // This insures tests for GIDs in the range
746 // [firstContiguousGID_, lastContiguousGID_] fail for processes
747 // with no local elements.
748 firstContiguousGID_ = indexBase_ + 1;
749 lastContiguousGID_ = indexBase_;
750 // glMap_ was default constructed, so it's already empty.
751 }
752
753 contiguous_ = false; // "Contiguous" is conservative.
754
755 if (req) req->wait();
756
757 const bool callComputeGlobalConstants = params.get() == nullptr ||
758 params->get("compute global constants", true);
759
760 if (callComputeGlobalConstants)
761 computeGlobalConstants();
762 else
763 distributed_ = params->get("distributed", true) && (comm->getSize() > 1);
764
765 // Create the Directory on demand in getRemoteIndexList().
766 // setupDirectory ();
767}
768
769template <class LocalOrdinal, class GlobalOrdinal, class Node>
772 const GlobalOrdinal indexList[],
775 const Teuchos::RCP<const Teuchos::Comm<int>>& comm,
776 const Teuchos::RCP<Teuchos::ParameterList>& params)
777 : comm_(comm)
778 , uniform_(false)
779 , directory_(new Directory<LocalOrdinal, GlobalOrdinal, Node>()) {
780 using std::endl;
781 const char funcName[] =
782 "Map(gblNumInds,indexList,indexListSize,indexBase,comm)";
783
784 const bool verbose = Details::Behavior::verbose("Map");
785 std::unique_ptr<std::string> prefix;
786 if (verbose) {
788 comm_.getRawPtr(), "Map", funcName);
789 std::ostringstream os;
790 os << *prefix << "Start" << endl;
791 std::cerr << os.str();
792 }
796 Impl::checkMapInputArray("(GST, const GO[], LO, GO, comm)",
797 indexList, static_cast<size_t>(indexListSize),
798 comm.getRawPtr());
799 // Not quite sure if I trust all code to behave correctly if the
800 // pointer is nonnull but the array length is nonzero, so I'll
801 // make sure the raw pointer is null if the length is zero.
802 const GlobalOrdinal* const indsRaw = indexListSize == 0 ? NULL : indexList;
803 Kokkos::View<const GlobalOrdinal*,
804 Kokkos::LayoutLeft,
805 Kokkos::HostSpace,
806 Kokkos::MemoryUnmanaged>
808 initWithNonownedHostIndexList(funcName, numGlobalElements, inds,
809 indexBase, comm, params);
810 if (verbose) {
811 std::ostringstream os;
812 os << *prefix << "Done" << endl;
813 std::cerr << os.str();
814 }
815}
817template <class LocalOrdinal, class GlobalOrdinal, class Node>
820 const Teuchos::ArrayView<const GlobalOrdinal>& entryList,
822 const Teuchos::RCP<const Teuchos::Comm<int>>& comm,
823 const Teuchos::RCP<Teuchos::ParameterList>& params)
824 : comm_(comm)
825 , uniform_(false)
826 , directory_(new Directory<LocalOrdinal, GlobalOrdinal, Node>()) {
827 using std::endl;
828 const char* funcName = "Map(gblNumInds,entryList(Teuchos::ArrayView),indexBase,comm)";
829
830 const bool verbose = Details::Behavior::verbose("Map");
831 std::unique_ptr<std::string> prefix;
832 if (verbose) {
834 comm_.getRawPtr(), "Map", funcName);
835 std::ostringstream os;
836 os << *prefix << "Start" << endl;
837 std::cerr << os.str();
838 }
842 const size_t numLclInds = static_cast<size_t>(entryList.size());
843 Impl::checkMapInputArray("(GST, ArrayView, GO, comm)",
844 entryList.getRawPtr(), numLclInds,
845 comm.getRawPtr());
846 // Not quite sure if I trust both ArrayView and View to behave
847 // correctly if the pointer is nonnull but the array length is
848 // nonzero, so I'll make sure it's null if the length is zero.
849 const GlobalOrdinal* const indsRaw =
850 numLclInds == 0 ? NULL : entryList.getRawPtr();
851 Kokkos::View<const GlobalOrdinal*,
852 Kokkos::LayoutLeft,
853 Kokkos::HostSpace,
854 Kokkos::MemoryUnmanaged>
856 initWithNonownedHostIndexList(funcName, numGlobalElements, inds,
857 indexBase, comm, params);
858 if (verbose) {
859 std::ostringstream os;
860 os << *prefix << "Done" << endl;
861 std::cerr << os.str();
862 }
863}
864
865template <class LocalOrdinal, class GlobalOrdinal, class Node>
868 const Kokkos::View<const GlobalOrdinal*, device_type>& entryList,
870 const Teuchos::RCP<const Teuchos::Comm<int>>& comm,
871 const Teuchos::RCP<Teuchos::ParameterList>& params)
872 : comm_(comm)
873 , uniform_(false)
874 , directory_(new Directory<LocalOrdinal, GlobalOrdinal, Node>()) {
875 using Kokkos::LayoutLeft;
876 using Kokkos::subview;
877 using Kokkos::View;
878 using Kokkos::view_alloc;
879 using Kokkos::WithoutInitializing;
880 using std::endl;
881 using Teuchos::arcp;
882 using Teuchos::ArrayView;
883 using Teuchos::as;
884 using Teuchos::broadcast;
885 using Teuchos::outArg;
886 using Teuchos::ptr;
887 using Teuchos::REDUCE_MAX;
888 using Teuchos::REDUCE_MIN;
889 using Teuchos::REDUCE_SUM;
890 using Teuchos::reduceAll;
891 using Teuchos::typeName;
892 using LO = local_ordinal_type;
893 using GST = global_size_t;
894 const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid();
895 const char funcName[] =
896 "Map(gblNumInds,entryList(Kokkos::View),indexBase,comm)";
897
898 const bool verbose = Details::Behavior::verbose("Map");
899 std::unique_ptr<std::string> prefix;
900 if (verbose) {
902 comm_.getRawPtr(), "Map", funcName);
903 std::ostringstream os;
904 os << *prefix << "Start" << endl;
905 std::cerr << os.str();
906 }
910 Impl::checkMapInputArray("(GST, Kokkos::View, GO, comm)",
911 entryList.data(),
912 static_cast<size_t>(entryList.extent(0)),
913 comm.getRawPtr());
914
915 // The user has specified the distribution of indices over the
916 // processes, via the input array of global indices on each
917 // process. The distribution is not necessarily contiguous or
918 // equally shared over the processes.
919
920 // The length of the input array on this process is the number of
921 // local indices to associate with this process, even though the
922 // input array contains global indices. We assume that the number
923 // of local indices on a process can be stored in a size_t;
924 // numLocalElements_ is a size_t, so this variable and that should
925 // have the same type.
926 const size_t numLocalElements(entryList.size());
927
928 initialNonuniformDebugCheck(funcName, numGlobalElements,
930
931 // NOTE (mfh 20 Feb 2013, 10 Oct 2016) In some sense, this global
932 // reduction is redundant, since the directory Map will have to do
933 // the same thing. Thus, we could do the scan and broadcast for
934 // the directory Map here, and give the computed offsets to the
935 // directory Map's constructor. However, a reduction costs less
936 // than a scan and broadcast, so this still saves time if users of
937 // this Map don't ever need the Directory (i.e., if they never
938 // call getRemoteIndexList on this Map).
939 std::shared_ptr<Details::CommRequest> req;
941 if (numGlobalElements != GSTI) {
942 numGlobalElements_ = numGlobalElements; // Use the user's value.
943 } else { // The user wants us to compute the sum.
945 req = Details::iallreduce(numLocalElementsGST,
946 numGlobalElements_, REDUCE_SUM, *comm);
948
949 // mfh 20 Feb 2013: We've never quite done the right thing for
950 // duplicate GIDs here. Duplicate GIDs have always been counted
951 // distinctly in numLocalElements_, and thus should get a
952 // different LID. However, we've always used std::map or a hash
953 // table for the GID -> LID lookup table, so distinct GIDs always
954 // map to the same LID. Furthermore, the order of the input GID
955 // list matters, so it's not desirable to sort for determining
956 // uniqueness.
957 //
958 // I've chosen for now to write this code as if the input GID list
959 // contains no duplicates. If this is not desired, we could use
960 // the lookup table itself to determine uniqueness: If we haven't
961 // seen the GID before, it gets a new LID and it's added to the
962 // LID -> GID and GID -> LID tables. If we have seen the GID
963 // before, it doesn't get added to either table. I would
964 // implement this, but it would cost more to do the double lookups
965 // in the table (one to check, and one to insert).
966 //
967 // More importantly, since we build the GID -> LID table in (a
968 // thread-) parallel (way), the order in which duplicate GIDs may
969 // get inserted is not defined. This would make the assignment of
970 // LID to GID nondeterministic.
971
972 numLocalElements_ = numLocalElements;
973 indexBase_ = indexBase;
974
975 minMyGID_ = indexBase_;
976 maxMyGID_ = indexBase_;
977
978 // NOTE (mfh 27 May 2015): While finding the initial contiguous
979 // GID range requires looking at all the GIDs in the range,
980 // dismissing an interval of GIDs only requires looking at the
981 // first and last GIDs. Thus, we could do binary search backwards
982 // from the end in order to catch the common case of a contiguous
983 // interval followed by noncontiguous entries. On the other hand,
984 // we could just expose this case explicitly as yet another Map
985 // constructor, and avoid the trouble of detecting it.
986 if (numLocalElements_ > 0) {
987 // Find contiguous GID range, with the restriction that the
988 // beginning of the range starts with the first entry. While
989 // doing so, fill in the LID -> GID table.
990 typename decltype(lgMap_)::non_const_type lgMap(view_alloc("lgMap2", WithoutInitializing), numLocalElements_);
991
992 // Because you can't use lambdas in constructors on CUDA. Or using private/protected data.
993 // DEEP_COPY REVIEW - DEVICE-TO-DEVICE
994 Kokkos::deep_copy(typename device_type::execution_space(), lgMap, entryList);
996 Impl::computeConstantsOnDevice(entryList, minMyGID_, maxMyGID_, firstContiguousGID_, lastContiguousGID_, lastContiguousGID_loc);
998 auto nonContigGids = Kokkos::subview(entryList, std::pair<size_t, size_t>(firstNonContiguous_loc, entryList.extent(0)));
999
1000 // NOTE: We do not fill the glMapHost_ and lgMapHost_ views here. They will be filled lazily later
1003
1004 // "Commit" the local-to-global lookup table we filled in above.
1005 lgMap_ = lgMap;
1006
1007 } else {
1008 minMyGID_ = std::numeric_limits<GlobalOrdinal>::max();
1009 maxMyGID_ = std::numeric_limits<GlobalOrdinal>::lowest();
1010 // This insures tests for GIDs in the range
1011 // [firstContiguousGID_, lastContiguousGID_] fail for processes
1012 // with no local elements.
1013 firstContiguousGID_ = indexBase_ + 1;
1014 lastContiguousGID_ = indexBase_;
1015 // glMap_ was default constructed, so it's already empty.
1016 }
1017
1018 if (req) req->wait();
1019
1020 const bool callComputeGlobalConstants = params.get() == nullptr ||
1021 params->get("compute global constants", true);
1022
1023 if (callComputeGlobalConstants)
1025 else
1026 distributed_ = params->get("distributed", true) && (comm->getSize() > 1);
1027
1028 contiguous_ = false; // "Contiguous" is conservative.
1029
1030 // Create the Directory on demand in getRemoteIndexList().
1031 // setupDirectory ();
1032
1033 if (verbose) {
1034 std::ostringstream os;
1035 os << *prefix << "Done" << endl;
1036 std::cerr << os.str();
1037 }
1038}
1039
1040template <class LocalOrdinal, class GlobalOrdinal, class Node>
1042 using GO = global_ordinal_type;
1043 using GST = global_size_t;
1044
1045 if (haveGlobalConstants_)
1046 return;
1047
1048 Tpetra::Details::ProfilingRegion pr("Tpetra::Map::computeGlobalConstants");
1049
1050 // Compute the min and max of all processes' GIDs. If
1051 // numLocalElements_ == 0 on this process, minMyGID_ and maxMyGID_
1052 // are both indexBase_. This is wrong, but fixing it would
1053 // require either a fancy sparse all-reduce, or a custom reduction
1054 // operator that ignores invalid values ("invalid" means
1055 // Tpetra::Details::OrdinalTraits<GO>::invalid()).
1056 //
1057 // Also, while we're at it, use the same all-reduce to figure out
1058 // if the Map is distributed. "Distributed" means that there is
1059 // at least one process with a number of local elements less than
1060 // the number of global elements.
1061 //
1062 // We're computing the min and max of all processes' GIDs using a
1063 // single MAX all-reduce, because min(x,y) = -max(-x,-y) (when x
1064 // and y are signed). (This lets us combine the min and max into
1065 // a single all-reduce.) If each process sets localDist=1 if its
1066 // number of local elements is strictly less than the number of
1067 // global elements, and localDist=0 otherwise, then a MAX
1068 // all-reduce on localDist tells us if the Map is distributed (1
1069 // if yes, 0 if no). Thus, we can append localDist onto the end
1070 // of the data and get the global result from the all-reduce.
1071 Kokkos::View<GO*, Kokkos::HostSpace> minMaxInput(Kokkos::ViewAllocateWithoutInitializing("minMaxInput"), 3);
1072 Kokkos::View<GO*, Kokkos::HostSpace> minMaxOutput(Kokkos::ViewAllocateWithoutInitializing("minMaxOutput"), 3);
1073
1074 minMaxInput[0] = std::numeric_limits<GO>::max() - minMyGID_;
1075 minMaxInput[1] = maxMyGID_;
1076 minMaxInput[2] = std::numeric_limits<GO>::max() - static_cast<GO>(numLocalElements_);
1077
1078 Teuchos::reduceAll<int, GO>(*comm_, Teuchos::REDUCE_MAX, 3, minMaxInput.data(), minMaxOutput.data());
1079
1080 minAllGID_ = std::numeric_limits<GO>::max() - minMaxOutput[0];
1081 maxAllGID_ = minMaxOutput[1];
1082 const GO minNumLocalElements = std::numeric_limits<GO>::max() - minMaxOutput[2];
1083
1084 distributed_ = (comm_->getSize() > 1 && (static_cast<GST>(minNumLocalElements) < numGlobalElements_));
1085
1086 haveGlobalConstants_ = true;
1087
1089 minAllGID_ < indexBase_,
1090 std::invalid_argument,
1091 "Tpetra::Map constructor (noncontiguous): "
1092 "Minimum global ID = "
1093 << minAllGID_ << " over all process(es) is "
1094 "less than the given indexBase = "
1095 << indexBase_ << ".");
1096}
1097
1098template <class LocalOrdinal, class GlobalOrdinal, class Node>
1100 if (!Kokkos::is_initialized()) {
1101 std::ostringstream os;
1102 os << "WARNING: Tpetra::Map destructor (~Map()) is being called after "
1103 "Kokkos::finalize() has been called. This is user error! There are "
1104 "two likely causes: "
1105 << std::endl
1106 << " 1. You have a static Tpetra::Map (or RCP or shared_ptr of a Map)"
1107 << std::endl
1108 << " 2. You declare and construct a Tpetra::Map (or RCP or shared_ptr "
1109 "of a Tpetra::Map) at the same scope in main() as Kokkos::finalize() "
1110 "or Tpetra::finalize()."
1111 << std::endl
1112 << std::endl
1113 << "Don't do either of these! Please refer to GitHib Issue #2372."
1114 << std::endl;
1115 ::Tpetra::Details::printOnce(std::cerr, os.str(),
1116 this->getComm().getRawPtr());
1117 } else {
1118 using ::Tpetra::Details::mpiIsFinalized;
1119 using ::Tpetra::Details::mpiIsInitialized;
1120 using ::Tpetra::Details::teuchosCommIsAnMpiComm;
1121
1122 Teuchos::RCP<const Teuchos::Comm<int>> comm = this->getComm();
1123 if (!comm.is_null() && teuchosCommIsAnMpiComm(*comm) &&
1124 mpiIsInitialized() && mpiIsFinalized()) {
1125 // Tpetra itself does not require MPI, even if building with
1126 // MPI. It is legal to create Tpetra objects that do not use
1127 // MPI, even in an MPI program. However, calling Tpetra stuff
1128 // after MPI_Finalize() has been called is a bad idea, since
1129 // some Tpetra defaults may use MPI if available.
1130 std::ostringstream os;
1131 os << "WARNING: Tpetra::Map destructor (~Map()) is being called after "
1132 "MPI_Finalize() has been called. This is user error! There are "
1133 "two likely causes: "
1134 << std::endl
1135 << " 1. You have a static Tpetra::Map (or RCP or shared_ptr of a Map)"
1136 << std::endl
1137 << " 2. You declare and construct a Tpetra::Map (or RCP or shared_ptr "
1138 "of a Tpetra::Map) at the same scope in main() as MPI_finalize() or "
1139 "Tpetra::finalize()."
1140 << std::endl
1141 << std::endl
1142 << "Don't do either of these! Please refer to GitHib Issue #2372."
1143 << std::endl;
1144 ::Tpetra::Details::printOnce(std::cerr, os.str(), comm.getRawPtr());
1145 }
1146 }
1147 // mfh 20 Mar 2018: We can't check Tpetra::isInitialized() yet,
1148 // because Tpetra does not yet require Tpetra::initialize /
1149 // Tpetra::finalize.
1150}
1151
1152template <class LocalOrdinal, class GlobalOrdinal, class Node>
1155 getComm().is_null(), std::logic_error,
1156 "Tpetra::Map::isOneToOne: "
1157 "getComm() returns null. Please report this bug to the Tpetra "
1158 "developers.");
1159
1160 // This is a collective operation, if it hasn't been called before.
1161 setupDirectory();
1162 return directory_->isOneToOne(*this);
1163}
1164
1165template <class LocalOrdinal, class GlobalOrdinal, class Node>
1169 if (isContiguous()) {
1170 if (globalIndex < getMinGlobalIndex() ||
1171 globalIndex > getMaxGlobalIndex()) {
1172 return Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid();
1173 }
1174 return static_cast<LocalOrdinal>(globalIndex - getMinGlobalIndex());
1175 } else if (globalIndex >= firstContiguousGID_ &&
1176 globalIndex <= lastContiguousGID_) {
1177 return static_cast<LocalOrdinal>(globalIndex - firstContiguousGID_);
1178 } else {
1179 // If the given global index is not in the table, this returns
1180 // the same value as OrdinalTraits<LocalOrdinal>::invalid().
1181 // glMapHost_ is Host and does not assume UVM
1182 lazyPushToHost();
1183 return glMapHost_.get(globalIndex);
1184 }
1185}
1186
1187template <class LocalOrdinal, class GlobalOrdinal, class Node>
1191 if (localIndex < getMinLocalIndex() || localIndex > getMaxLocalIndex()) {
1192 return Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid();
1193 }
1194 if (isContiguous()) {
1195 return getMinGlobalIndex() + localIndex;
1196 } else {
1197 // This is a host Kokkos::View access, with no RCP or ArrayRCP
1198 // involvement. As a result, it is thread safe.
1199 //
1200 // lgMapHost_ is a host pointer; this does NOT assume UVM.
1201 lazyPushToHost();
1202 return lgMapHost_[localIndex];
1203 }
1204}
1205
1206template <class LocalOrdinal, class GlobalOrdinal, class Node>
1208 const local_ordinal_type localIndices[], size_t numEntries, global_ordinal_type globalIndices[]) const {
1209 auto const minGI = getMinGlobalIndex();
1210 auto const minLI = getMinLocalIndex();
1211 auto const maxLI = getMaxLocalIndex();
1212 if (isContiguous()) {
1213 for (size_t i = 0; i < numEntries; i++) {
1216 return true;
1217 }
1219 }
1220 } else {
1221 // This is a host Kokkos::View access, with no RCP or ArrayRCP
1222 // involvement. As a result, it is thread safe.
1223 //
1224 // lgMapHost_ is a host pointer; this does NOT assume UVM.
1225 lazyPushToHost();
1226 for (size_t i = 0; i < numEntries; i++) {
1227 auto lclInd = localIndices[i];
1228 if (lclInd < minLI || lclInd > maxLI) {
1229 return true;
1230 }
1231 globalIndices[i] = lgMapHost_[lclInd];
1232 }
1233 }
1234 return false;
1235}
1236
1237template <class LocalOrdinal, class GlobalOrdinal, class Node>
1240 if (localIndex < getMinLocalIndex() || localIndex > getMaxLocalIndex()) {
1241 return false;
1242 } else {
1243 return true;
1244 }
1245}
1246
1247template <class LocalOrdinal, class GlobalOrdinal, class Node>
1250 return this->getLocalElement(globalIndex) !=
1251 Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid();
1252}
1253
1254template <class LocalOrdinal, class GlobalOrdinal, class Node>
1256 return uniform_;
1257}
1258
1259template <class LocalOrdinal, class GlobalOrdinal, class Node>
1261 return contiguous_;
1262}
1263
1264template <class LocalOrdinal, class GlobalOrdinal, class Node>
1267 getLocalMap() const {
1268 return local_map_type(glMap_, lgMap_, getIndexBase(),
1269 getMinGlobalIndex(), getMaxGlobalIndex(),
1270 firstContiguousGID_, lastContiguousGID_,
1271 getLocalNumElements(), isContiguous());
1272}
1273
1274template <class LocalOrdinal, class GlobalOrdinal, class Node>
1277 using Teuchos::outArg;
1278 using Teuchos::REDUCE_MIN;
1279 using Teuchos::reduceAll;
1280 //
1281 // Tests that avoid the Boolean all-reduce below by using
1282 // globally consistent quantities.
1283 //
1284 if (this == &map) {
1285 // Pointer equality on one process always implies pointer
1286 // equality on all processes, since Map is immutable.
1287 return true;
1288 } else if (getComm()->getSize() != map.getComm()->getSize()) {
1289 // The two communicators have different numbers of processes.
1290 // It's not correct to call isCompatible() in that case. This
1291 // may result in the all-reduce hanging below.
1292 return false;
1293 } else if (getGlobalNumElements() != map.getGlobalNumElements()) {
1294 // Two Maps are definitely NOT compatible if they have different
1295 // global numbers of indices.
1296 return false;
1297 } else if (isContiguous() && isUniform() &&
1298 map.isContiguous() && map.isUniform()) {
1299 // Contiguous uniform Maps with the same number of processes in
1300 // their communicators, and with the same global numbers of
1301 // indices, are always compatible.
1302 return true;
1303 } else if (!isContiguous() && !map.isContiguous() &&
1304 lgMap_.extent(0) != 0 && map.lgMap_.extent(0) != 0 &&
1305 lgMap_.data() == map.lgMap_.data()) {
1306 // Noncontiguous Maps whose global index lists are nonempty and
1307 // have the same pointer must be the same (and therefore
1308 // contiguous).
1309 //
1310 // Nonempty is important. For example, consider a communicator
1311 // with two processes, and two Maps that share this
1312 // communicator, with zero global indices on the first process,
1313 // and different nonzero numbers of global indices on the second
1314 // process. In that case, on the first process, the pointers
1315 // would both be NULL.
1316 return true;
1317 }
1318
1320 getGlobalNumElements() != map.getGlobalNumElements(), std::logic_error,
1321 "Tpetra::Map::isCompatible: There's a bug in this method. We've already "
1322 "checked that this condition is true above, but it's false here. "
1323 "Please report this bug to the Tpetra developers.");
1324
1325 // Do both Maps have the same number of indices on each process?
1326 const int locallyCompat =
1327 (getLocalNumElements() == map.getLocalNumElements()) ? 1 : 0;
1328
1329 int globallyCompat = 0;
1331 return (globallyCompat == 1);
1332}
1333
1334template <class LocalOrdinal, class GlobalOrdinal, class Node>
1337 using Teuchos::ArrayView;
1338 using GO = global_ordinal_type;
1339 using size_type = typename ArrayView<const GO>::size_type;
1340
1341 // If both Maps are contiguous, we can compare their GID ranges
1342 // easily by looking at the min and max GID on this process.
1343 // Otherwise, we'll compare their GID lists. If only one Map is
1344 // contiguous, then we only have to call getLocalElementList() on
1345 // the noncontiguous Map. (It's best to avoid calling it on a
1346 // contiguous Map, since it results in unnecessary storage that
1347 // persists for the lifetime of the Map.)
1348
1349 if (this == &map) {
1350 // Pointer equality on one process always implies pointer
1351 // equality on all processes, since Map is immutable.
1352 return true;
1353 } else if (getLocalNumElements() != map.getLocalNumElements()) {
1354 return false;
1355 } else if (getMinGlobalIndex() != map.getMinGlobalIndex() ||
1356 getMaxGlobalIndex() != map.getMaxGlobalIndex()) {
1357 return false;
1358 } else {
1359 if (isContiguous()) {
1360 if (map.isContiguous()) {
1361 return true; // min and max match, so the ranges match.
1362 } else { // *this is contiguous, but map is not contiguous
1364 !this->isContiguous() || map.isContiguous(), std::logic_error,
1365 "Tpetra::Map::locallySameAs: BUG");
1366 ArrayView<const GO> rhsElts = map.getLocalElementList();
1367 const GO minLhsGid = this->getMinGlobalIndex();
1368 const size_type numRhsElts = rhsElts.size();
1369 for (size_type k = 0; k < numRhsElts; ++k) {
1370 const GO curLhsGid = minLhsGid + static_cast<GO>(k);
1371 if (curLhsGid != rhsElts[k]) {
1372 return false; // stop on first mismatch
1373 }
1374 }
1375 return true;
1376 }
1377 } else if (map.isContiguous()) { // *this is not contiguous, but map is
1379 this->isContiguous() || !map.isContiguous(), std::logic_error,
1380 "Tpetra::Map::locallySameAs: BUG");
1381 ArrayView<const GO> lhsElts = this->getLocalElementList();
1382 const GO minRhsGid = map.getMinGlobalIndex();
1383 const size_type numLhsElts = lhsElts.size();
1384 for (size_type k = 0; k < numLhsElts; ++k) {
1385 const GO curRhsGid = minRhsGid + static_cast<GO>(k);
1386 if (curRhsGid != lhsElts[k]) {
1387 return false; // stop on first mismatch
1388 }
1389 }
1390 return true;
1391 } else if (this->lgMap_.data() == map.lgMap_.data()) {
1392 // Pointers to LID->GID "map" (actually just an array) are the
1393 // same, and the number of GIDs are the same.
1394 return this->getLocalNumElements() == map.getLocalNumElements();
1395 } else { // we actually have to compare the GIDs
1396 if (this->getLocalNumElements() != map.getLocalNumElements()) {
1397 return false; // We already checked above, but check just in case
1398 } else {
1399 using range_type = Kokkos::RangePolicy<LocalOrdinal, typename node_type::execution_space>;
1400
1401 auto lhsLclMap = getLocalMap();
1402 auto rhsLclMap = map.getLocalMap();
1403
1405 Kokkos::parallel_reduce(
1406 "Tpetra::Map::locallySameAs",
1407 range_type(0, this->getLocalNumElements()),
1409 if (lhsLclMap.getGlobalElement(lid) != rhsLclMap.getGlobalElement(lid))
1410 ++numMismatches;
1411 },
1413
1414 return (numMismatchedElements == 0);
1415 }
1416 }
1417 }
1418}
1419
1420template <class LocalOrdinal, class GlobalOrdinal, class Node>
1423 if (this == &map)
1424 return true;
1425
1426 // We are going to check if lmap1 is fitted into lmap2:
1427 // Is lmap1 (map) a subset of lmap2 (this)?
1428 // And do the first lmap1.getLocalNumElements() global elements
1429 // of lmap1,lmap2 owned on each process exactly match?
1430 auto lmap1 = map.getLocalMap();
1431 auto lmap2 = this->getLocalMap();
1432
1433 auto numLocalElements1 = lmap1.getLocalNumElements();
1434 auto numLocalElements2 = lmap2.getLocalNumElements();
1435
1437 // There are more indices in the first map on this process than in second map.
1438 return false;
1439 }
1440
1441 if (lmap1.isContiguous() && lmap2.isContiguous()) {
1442 // When both Maps are contiguous, just check the interval inclusion.
1443 return ((lmap1.getMinGlobalIndex() == lmap2.getMinGlobalIndex()) &&
1444 (lmap1.getMaxGlobalIndex() <= lmap2.getMaxGlobalIndex()));
1445 }
1446
1447 if (lmap1.getMinGlobalIndex() < lmap2.getMinGlobalIndex() ||
1448 lmap1.getMaxGlobalIndex() > lmap2.getMaxGlobalIndex()) {
1449 // The second map does not include the first map bounds, and thus some of
1450 // the first map global indices are not in the second map.
1451 return false;
1452 }
1453
1454 using LO = local_ordinal_type;
1455 using range_type =
1456 Kokkos::RangePolicy<LO, typename node_type::execution_space>;
1457
1458 // Check all elements.
1459 LO numDiff = 0;
1460 Kokkos::parallel_reduce(
1461 "isLocallyFitted",
1462 range_type(0, numLocalElements1),
1463 KOKKOS_LAMBDA(const LO i, LO& diff) {
1464 diff += (lmap1.getGlobalElement(i) != lmap2.getGlobalElement(i));
1465 },
1466 numDiff);
1467
1468 return (numDiff == 0);
1469}
1470
1471template <class LocalOrdinal, class GlobalOrdinal, class Node>
1474 using Teuchos::outArg;
1475 using Teuchos::REDUCE_MIN;
1476 using Teuchos::reduceAll;
1477 //
1478 // Tests that avoid the Boolean all-reduce below by using
1479 // globally consistent quantities.
1480 //
1481 if (this == &map) {
1482 // Pointer equality on one process always implies pointer
1483 // equality on all processes, since Map is immutable.
1484 return true;
1485 } else if (getComm()->getSize() != map.getComm()->getSize()) {
1486 // The two communicators have different numbers of processes.
1487 // It's not correct to call isSameAs() in that case. This
1488 // may result in the all-reduce hanging below.
1489 return false;
1490 } else if (getGlobalNumElements() != map.getGlobalNumElements()) {
1491 // Two Maps are definitely NOT the same if they have different
1492 // global numbers of indices.
1493 return false;
1494 } else if (haveGlobalConstants() && map.haveGlobalConstants() && (getMinAllGlobalIndex() != map.getMinAllGlobalIndex() || getMaxAllGlobalIndex() != map.getMaxAllGlobalIndex() || getIndexBase() != map.getIndexBase())) {
1495 // If the global min or max global index doesn't match, or if
1496 // the index base doesn't match, then the Maps aren't the same.
1497 return false;
1498 } else if (haveGlobalConstants() && map.haveGlobalConstants() && (isDistributed() != map.isDistributed())) {
1499 // One Map is distributed and the other is not, which means that
1500 // the Maps aren't the same.
1501 return false;
1502 } else if (isContiguous() && isUniform() &&
1503 map.isContiguous() && map.isUniform()) {
1504 // Contiguous uniform Maps with the same number of processes in
1505 // their communicators, with the same global numbers of indices,
1506 // and with matching index bases and ranges, must be the same.
1507 return true;
1508 }
1509
1510 // The two communicators must have the same number of processes,
1511 // with process ranks occurring in the same order. This uses
1512 // MPI_COMM_COMPARE. The MPI 3.1 standard (Section 6.4) says:
1513 // "Operations that access communicators are local and their
1514 // execution does not require interprocess communication."
1515 // However, just to be sure, I'll put this call after the above
1516 // tests that don't communicate.
1517 if (!::Tpetra::Details::congruent(*comm_, *(map.getComm()))) {
1518 return false;
1519 }
1520
1521 // If we get this far, we need to check local properties and then
1522 // communicate local sameness across all processes.
1523 const int isSame_lcl = locallySameAs(map) ? 1 : 0;
1524
1525 // Return true if and only if all processes report local sameness.
1526 int isSame_gbl = 0;
1528 return isSame_gbl == 1;
1529}
1530
1531namespace { // (anonymous)
1532template <class LO, class GO, class DT>
1533class FillLgMap {
1534 public:
1535 FillLgMap(const Kokkos::View<GO*, DT>& lgMap,
1536 const GO startGid)
1537 : lgMap_(lgMap)
1538 , startGid_(startGid) {
1539 Kokkos::RangePolicy<LO, typename DT::execution_space>
1540 range(static_cast<LO>(0), static_cast<LO>(lgMap.size()));
1541 Kokkos::parallel_for(range, *this);
1542 }
1543
1544 KOKKOS_INLINE_FUNCTION void operator()(const LO& lid) const {
1545 lgMap_(lid) = startGid_ + static_cast<GO>(lid);
1546 }
1547
1548 private:
1549 const Kokkos::View<GO*, DT> lgMap_;
1550 const GO startGid_;
1551};
1552
1553} // namespace
1554
1555template <class LocalOrdinal, class GlobalOrdinal, class Node>
1556typename Map<LocalOrdinal, GlobalOrdinal, Node>::global_indices_array_type
1558 using std::endl;
1559 using LO = local_ordinal_type;
1560 using GO = global_ordinal_type;
1561 using const_lg_view_type = decltype(lgMap_);
1562 using lg_view_type = typename const_lg_view_type::non_const_type;
1563 const bool debug = Details::Behavior::debug("Map");
1564 const bool verbose = Details::Behavior::verbose("Map");
1565
1566 std::unique_ptr<std::string> prefix;
1567 if (verbose) {
1569 comm_.getRawPtr(), "Map", "getMyGlobalIndices");
1570 std::ostringstream os;
1571 os << *prefix << "Start" << endl;
1572 std::cerr << os.str();
1573 }
1574
1575 // If the local-to-global mapping doesn't exist yet, and if we
1576 // have local entries, then create and fill the local-to-global
1577 // mapping.
1579 lgMap_.extent(0) == 0 && numLocalElements_ > 0;
1580
1582 if (verbose) {
1583 std::ostringstream os;
1584 os << *prefix << "Need to create lgMap" << endl;
1585 std::cerr << os.str();
1586 }
1587 if (debug) {
1588 // The local-to-global mapping should have been set up already
1589 // for a noncontiguous map.
1590 TEUCHOS_TEST_FOR_EXCEPTION(!isContiguous(), std::logic_error,
1591 "Tpetra::Map::getMyGlobalIndices: The local-to-global "
1592 "mapping (lgMap_) should have been set up already for a "
1593 "noncontiguous Map. Please report this bug to the Tpetra "
1594 "developers.");
1595 }
1596 const LO numElts = static_cast<LO>(getLocalNumElements());
1597
1598 using Kokkos::view_alloc;
1599 using Kokkos::WithoutInitializing;
1600 lg_view_type lgMap("lgMap3", numElts);
1601 if (verbose) {
1602 std::ostringstream os;
1603 os << *prefix << "Fill lgMap" << endl;
1604 std::cerr << os.str();
1605 }
1607
1608 if (verbose) {
1609 std::ostringstream os;
1610 os << *prefix << "Copy lgMap to lgMapHost" << endl;
1611 std::cerr << os.str();
1612 }
1613
1614 auto lgMapHost = Kokkos::create_mirror_view(Kokkos::HostSpace(), lgMap);
1615 // DEEP_COPY REVIEW - DEVICE-TO-HOST
1617 Kokkos::deep_copy(exec_instance, lgMapHost, lgMap);
1618
1619 // There's a non-trivial chance we'll grab this on the host,
1620 // so let's make sure the copy finishes
1621 exec_instance.fence();
1622
1623 // "Commit" the local-to-global lookup table we filled in above.
1624 lgMap_ = lgMap;
1625 lgMapHost_ = lgMapHost;
1626 } else {
1627 lazyPushToHost();
1628 }
1629
1630 if (verbose) {
1631 std::ostringstream os;
1632 os << *prefix << "Done" << endl;
1633 std::cerr << os.str();
1634 }
1635 return lgMapHost_;
1636}
1637
1638template <class LocalOrdinal, class GlobalOrdinal, class Node>
1639typename Map<LocalOrdinal, GlobalOrdinal, Node>::global_indices_array_device_type
1641 using std::endl;
1642 using LO = local_ordinal_type;
1643 using GO = global_ordinal_type;
1644 using const_lg_view_type = decltype(lgMap_);
1645 using lg_view_type = typename const_lg_view_type::non_const_type;
1646 const bool debug = Details::Behavior::debug("Map");
1647 const bool verbose = Details::Behavior::verbose("Map");
1648
1649 std::unique_ptr<std::string> prefix;
1650 if (verbose) {
1652 comm_.getRawPtr(), "Map", "getMyGlobalIndicesDevice");
1653 std::ostringstream os;
1654 os << *prefix << "Start" << endl;
1655 std::cerr << os.str();
1656 }
1657
1658 // If the local-to-global mapping doesn't exist yet, and if we
1659 // have local entries, then create and fill the local-to-global
1660 // mapping.
1662 lgMap_.extent(0) == 0 && numLocalElements_ > 0;
1663
1665 if (verbose) {
1666 std::ostringstream os;
1667 os << *prefix << "Need to create lgMap" << endl;
1668 std::cerr << os.str();
1669 }
1670 if (debug) {
1671 // The local-to-global mapping should have been set up already
1672 // for a noncontiguous map.
1673 TEUCHOS_TEST_FOR_EXCEPTION(!isContiguous(), std::logic_error,
1674 "Tpetra::Map::getMyGlobalIndices: The local-to-global "
1675 "mapping (lgMap_) should have been set up already for a "
1676 "noncontiguous Map. Please report this bug to the Tpetra "
1677 "developers.");
1678 }
1679 const LO numElts = static_cast<LO>(getLocalNumElements());
1680
1681 using Kokkos::view_alloc;
1682 using Kokkos::WithoutInitializing;
1683 lg_view_type lgMap("lgMap4", numElts);
1684 if (verbose) {
1685 std::ostringstream os;
1686 os << *prefix << "Fill lgMap" << endl;
1687 std::cerr << os.str();
1688 }
1690
1691 // "Commit" the local-to-global lookup table we filled in above.
1692 lgMap_ = lgMap;
1693 }
1694
1695 if (verbose) {
1696 std::ostringstream os;
1697 os << *prefix << "Done" << endl;
1698 std::cerr << os.str();
1699 }
1700 return lgMap_;
1701}
1702
1703template <class LocalOrdinal, class GlobalOrdinal, class Node>
1704Teuchos::ArrayView<const GlobalOrdinal>
1706 using GO = global_ordinal_type;
1707
1708 // If the local-to-global mapping doesn't exist yet, and if we
1709 // have local entries, then create and fill the local-to-global
1710 // mapping.
1711 (void)this->getMyGlobalIndices();
1712
1713 // This does NOT assume UVM; lgMapHost_ is a host pointer.
1714 lazyPushToHost();
1715 const GO* lgMapHostRawPtr = lgMapHost_.data();
1716 // The third argument forces ArrayView not to try to track memory
1717 // in a debug build. We have to use it because the memory does
1718 // not belong to a Teuchos memory management class.
1719 return Teuchos::ArrayView<const GO>(
1721 lgMapHost_.extent(0),
1722 Teuchos::RCP_DISABLE_NODE_LOOKUP);
1723}
1724
1725template <class LocalOrdinal, class GlobalOrdinal, class Node>
1727 return distributed_;
1728}
1729
1730template <class LocalOrdinal, class GlobalOrdinal, class Node>
1732 using Teuchos::TypeNameTraits;
1733 std::ostringstream os;
1734
1735 os << "Tpetra::Map: {"
1736 << "LocalOrdinalType: " << TypeNameTraits<LocalOrdinal>::name()
1737 << ", GlobalOrdinalType: " << TypeNameTraits<GlobalOrdinal>::name()
1738 << ", NodeType: " << TypeNameTraits<Node>::name();
1739 if (this->getObjectLabel() != "") {
1740 os << ", Label: \"" << this->getObjectLabel() << "\"";
1741 }
1742 os << ", Global number of entries: " << getGlobalNumElements()
1743 << ", Number of processes: " << getComm()->getSize()
1744 << ", Uniform: " << (isUniform() ? "true" : "false")
1745 << ", Contiguous: " << (isContiguous() ? "true" : "false")
1746 << ", Distributed: " << (isDistributed() ? "true" : "false")
1747 << "}";
1748 return os.str();
1749}
1750
1755template <class LocalOrdinal, class GlobalOrdinal, class Node>
1756std::string
1758 localDescribeToString(const Teuchos::EVerbosityLevel vl) const {
1759 using LO = local_ordinal_type;
1760 using std::endl;
1761
1762 // This preserves current behavior of Map.
1763 if (vl < Teuchos::VERB_HIGH) {
1764 return std::string();
1765 }
1766 auto outStringP = Teuchos::rcp(new std::ostringstream());
1767 Teuchos::RCP<Teuchos::FancyOStream> outp =
1768 Teuchos::getFancyOStream(outStringP);
1769 Teuchos::FancyOStream& out = *outp;
1770
1771 auto comm = this->getComm();
1772 const int myRank = comm->getRank();
1773 const int numProcs = comm->getSize();
1774 out << "Process " << myRank << " of " << numProcs << ":" << endl;
1775 Teuchos::OSTab tab1(out);
1776
1777 const LO numEnt = static_cast<LO>(this->getLocalNumElements());
1778 out << "My number of entries: " << numEnt << endl
1779 << "My minimum global index: " << this->getMinGlobalIndex() << endl
1780 << "My maximum global index: " << this->getMaxGlobalIndex() << endl;
1781
1782 if (vl == Teuchos::VERB_EXTREME) {
1783 out << "My global indices: [";
1784 const LO minLclInd = this->getMinLocalIndex();
1785 for (LO k = 0; k < numEnt; ++k) {
1786 out << minLclInd + this->getGlobalElement(k);
1787 if (k + 1 < numEnt) {
1788 out << ", ";
1789 }
1790 }
1791 out << "]" << endl;
1792 }
1793
1794 out.flush(); // make sure the ostringstream got everything
1795 return outStringP->str();
1796}
1797
1798template <class LocalOrdinal, class GlobalOrdinal, class Node>
1800 describe(Teuchos::FancyOStream& out,
1801 const Teuchos::EVerbosityLevel verbLevel) const {
1802 using std::endl;
1803 using Teuchos::TypeNameTraits;
1804 using Teuchos::VERB_DEFAULT;
1805 using Teuchos::VERB_HIGH;
1806 using Teuchos::VERB_LOW;
1807 using Teuchos::VERB_NONE;
1808 using LO = local_ordinal_type;
1809 using GO = global_ordinal_type;
1810 const Teuchos::EVerbosityLevel vl =
1812
1813 if (vl == VERB_NONE) {
1814 return; // don't print anything
1815 }
1816 // If this Map's Comm is null, then the Map does not participate
1817 // in collective operations with the other processes. In that
1818 // case, it is not even legal to call this method. The reasonable
1819 // thing to do in that case is nothing.
1820 auto comm = this->getComm();
1821 if (comm.is_null()) {
1822 return;
1823 }
1824 const int myRank = comm->getRank();
1825 const int numProcs = comm->getSize();
1826
1827 // Only Process 0 should touch the output stream, but this method
1828 // in general may need to do communication. Thus, we may need to
1829 // preserve the current tab level across multiple "if (myRank ==
1830 // 0) { ... }" inner scopes. This is why we sometimes create
1831 // OSTab instances by pointer, instead of by value. We only need
1832 // to create them by pointer if the tab level must persist through
1833 // multiple inner scopes.
1834 Teuchos::RCP<Teuchos::OSTab> tab0, tab1;
1835
1836 if (myRank == 0) {
1837 // At every verbosity level but VERB_NONE, Process 0 prints.
1838 // By convention, describe() always begins with a tab before
1839 // printing.
1840 tab0 = Teuchos::rcp(new Teuchos::OSTab(out));
1841 out << "\"Tpetra::Map\":" << endl;
1842 tab1 = Teuchos::rcp(new Teuchos::OSTab(out));
1843 {
1844 out << "Template parameters:" << endl;
1845 Teuchos::OSTab tab2(out);
1846 out << "LocalOrdinal: " << TypeNameTraits<LO>::name() << endl
1847 << "GlobalOrdinal: " << TypeNameTraits<GO>::name() << endl
1848 << "Node: " << TypeNameTraits<Node>::name() << endl;
1849 }
1850 const std::string label = this->getObjectLabel();
1851 if (label != "") {
1852 out << "Label: \"" << label << "\"" << endl;
1853 }
1854 out << "Global number of entries: " << getGlobalNumElements() << endl
1855 << "Minimum global index: " << getMinAllGlobalIndex() << endl
1856 << "Maximum global index: " << getMaxAllGlobalIndex() << endl
1857 << "Index base: " << getIndexBase() << endl
1858 << "Number of processes: " << numProcs << endl
1859 << "Uniform: " << (isUniform() ? "true" : "false") << endl
1860 << "Contiguous: " << (isContiguous() ? "true" : "false") << endl
1861 << "Distributed: " << (isDistributed() ? "true" : "false") << endl;
1862 }
1863
1864 // This is collective over the Map's communicator.
1865 if (vl >= VERB_HIGH) { // VERB_HIGH or VERB_EXTREME
1866 const std::string lclStr = this->localDescribeToString(vl);
1868 }
1869}
1870
1871template <class LocalOrdinal, class GlobalOrdinal, class Node>
1872Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>>
1874 replaceCommWithSubset(const Teuchos::RCP<const Teuchos::Comm<int>>& newComm) const {
1875 using Teuchos::RCP;
1876 using Teuchos::rcp;
1877 using GST = global_size_t;
1878 using LO = local_ordinal_type;
1879 using GO = global_ordinal_type;
1880 using map_type = Map<LO, GO, Node>;
1881
1882 // mfh 26 Mar 2013: The lazy way to do this is simply to recreate
1883 // the Map by calling its ordinary public constructor, using the
1884 // original Map's data. This only involves O(1) all-reduces over
1885 // the new communicator, which in the common case only includes a
1886 // small number of processes.
1887
1888 // Create the Map to return.
1889 if (newComm.is_null() || newComm->getSize() < 1) {
1890 return Teuchos::null; // my process does not participate in the new Map
1891 } else if (newComm->getSize() == 1) {
1892 lazyPushToHost();
1893
1894 // The case where the new communicator has only one process is
1895 // easy. We don't have to communicate to get all the
1896 // information we need. Use the default comm to create the new
1897 // Map, then fill in all the fields directly.
1898 RCP<map_type> newMap(new map_type());
1899
1900 newMap->comm_ = newComm;
1901 // mfh 07 Oct 2016: Preserve original behavior, even though the
1902 // original index base may no longer be the globally min global
1903 // index. See #616 for why this doesn't matter so much anymore.
1904 newMap->indexBase_ = this->indexBase_;
1905 newMap->numGlobalElements_ = this->numLocalElements_;
1906 newMap->numLocalElements_ = this->numLocalElements_;
1907 newMap->minMyGID_ = this->minMyGID_;
1908 newMap->maxMyGID_ = this->maxMyGID_;
1909 newMap->minAllGID_ = this->minMyGID_;
1910 newMap->maxAllGID_ = this->maxMyGID_;
1911 newMap->firstContiguousGID_ = this->firstContiguousGID_;
1912 newMap->lastContiguousGID_ = this->lastContiguousGID_;
1913 newMap->haveGlobalConstants_ = this->haveGlobalConstants_;
1914 // Since the new communicator has only one process, neither
1915 // uniformity nor contiguity have changed.
1916 newMap->uniform_ = this->uniform_;
1917 newMap->contiguous_ = this->contiguous_;
1918 // The new communicator only has one process, so the new Map is
1919 // not distributed.
1920 newMap->distributed_ = false;
1921 newMap->lgMap_ = this->lgMap_;
1922 newMap->lgMapHost_ = this->lgMapHost_;
1923 newMap->glMap_ = this->glMap_;
1924 newMap->glMapHost_ = this->glMapHost_;
1925 // It's OK not to initialize the new Map's Directory.
1926 // This is initialized lazily, on first call to getRemoteIndexList.
1927
1928 return newMap;
1929 } else { // newComm->getSize() != 1
1930 // Even if the original Map is contiguous, the new Map might not
1931 // be, especially if the excluded processes have ranks != 0 or
1932 // newComm->getSize()-1. The common case for this method is to
1933 // exclude many (possibly even all but one) processes, so it
1934 // likely doesn't pay to do the global communication (over the
1935 // original communicator) to figure out whether we can optimize
1936 // the result Map. Thus, we just set up the result Map as
1937 // noncontiguous.
1938 //
1939 // TODO (mfh 07 Oct 2016) We don't actually need to reconstruct
1940 // the global-to-local table, etc. Optimize this code path to
1941 // avoid unnecessary local work.
1942
1943 // Make Map (re)compute the global number of elements.
1944 const GST RECOMPUTE = Tpetra::Details::OrdinalTraits<GST>::invalid();
1945 // TODO (mfh 07 Oct 2016) If we use any Map constructor, we have
1946 // to use the noncontiguous Map constructor, since the new Map
1947 // might not be contiguous. Even if the old Map was contiguous,
1948 // some process in the "middle" might have been excluded. If we
1949 // want to avoid local work, we either have to do the setup by
1950 // hand, or write a new Map constructor.
1951#if 1
1952 // The disabled code here throws the following exception in
1953 // Map's replaceCommWithSubset test:
1954 //
1955 // Throw test that evaluated to true: static_cast<unsigned long long> (numKeys) > static_cast<unsigned long long> (::Kokkos::ArithTraits<ValueType>::max ())
1956 // 10:
1957 // 10: Tpetra::Details::FixedHashTable: The number of keys -3 is greater than the maximum representable ValueType value 2147483647. This means that it is not possible to use this constructor.
1958 // 10: Process 3: origComm->replaceCommWithSubset(subsetComm) threw an exception: /scratch/prj/Trilinos/Trilinos/packages/tpetra/core/src/Tpetra_Details_FixedHashTable_def.hpp:1044:
1959
1960 auto lgMap = this->getMyGlobalIndices();
1961 using size_type =
1962 typename std::decay<decltype(lgMap.extent(0))>::type;
1963 const size_type lclNumInds =
1964 static_cast<size_type>(this->getLocalNumElements());
1965 using Teuchos::TypeNameTraits;
1966 TEUCHOS_TEST_FOR_EXCEPTION(lgMap.extent(0) != lclNumInds, std::logic_error,
1967 "Tpetra::Map::replaceCommWithSubset: Result of getMyGlobalIndices() "
1968 "has length "
1969 << lgMap.extent(0) << " (of type " << TypeNameTraits<size_type>::name() << ") != this->getLocalNumElements()"
1970 " = "
1971 << this->getLocalNumElements() << ". The latter, upon being "
1972 "cast to size_type = "
1974 "becomes "
1975 << lclNumInds << ". Please report this bug to the Tpetra "
1976 "developers.");
1977#else
1978 Teuchos::ArrayView<const GO> lgMap = this->getLocalElementList();
1979#endif // 1
1980
1981 const GO indexBase = this->getIndexBase();
1982 // map stores HostSpace of CudaSpace but constructor is still CudaUVMSpace
1983 auto lgMap_device = Kokkos::create_mirror_view_and_copy(device_type(), lgMap);
1984 return rcp(new map_type(RECOMPUTE, lgMap_device, indexBase, newComm));
1985 }
1986}
1987
1988template <class LocalOrdinal, class GlobalOrdinal, class Node>
1989Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node>>
1991 removeEmptyProcesses() const {
1992 using Teuchos::Comm;
1993 using Teuchos::null;
1994 using Teuchos::outArg;
1995 using Teuchos::RCP;
1996 using Teuchos::rcp;
1997 using Teuchos::REDUCE_MIN;
1998 using Teuchos::reduceAll;
1999
2000 // Create the new communicator. split() returns a valid
2001 // communicator on all processes. On processes where color == 0,
2002 // ignore the result. Passing key == 0 tells MPI to order the
2003 // processes in the new communicator by their rank in the old
2004 // communicator.
2005 const int color = (numLocalElements_ == 0) ? 0 : 1;
2006 // MPI_Comm_split must be called collectively over the original
2007 // communicator. We can't just call it on processes with color
2008 // one, even though we will ignore its result on processes with
2009 // color zero.
2010 RCP<const Comm<int>> newComm = comm_->split(color, 0);
2011 if (color == 0) {
2012 newComm = null;
2013 }
2014
2015 TEUCHOS_TEST_FOR_EXCEPTION(!haveGlobalConstants(), std::logic_error, "\"Map::removeEmptyProcesses\" called, but global constants have not been computed with \"Map::computeGlobalConstants\".");
2016
2017 // Create the Map to return.
2018 if (newComm.is_null()) {
2019 return null; // my process does not participate in the new Map
2020 } else {
2021 RCP<Map> map = rcp(new Map());
2022
2023 map->comm_ = newComm;
2024 map->indexBase_ = indexBase_;
2025 map->numGlobalElements_ = numGlobalElements_;
2026 map->numLocalElements_ = numLocalElements_;
2027 map->minMyGID_ = minMyGID_;
2028 map->maxMyGID_ = maxMyGID_;
2029 map->minAllGID_ = minAllGID_;
2030 map->maxAllGID_ = maxAllGID_;
2031 map->firstContiguousGID_ = firstContiguousGID_;
2032 map->lastContiguousGID_ = lastContiguousGID_;
2033 map->haveGlobalConstants_ = haveGlobalConstants_;
2034
2035 // Uniformity and contiguity have not changed. The directory
2036 // has changed, but we've taken care of that above.
2037 map->uniform_ = uniform_;
2038 map->contiguous_ = contiguous_;
2039
2040 // If the original Map was NOT distributed, then the new Map
2041 // cannot be distributed.
2042 //
2043 // If the number of processes in the new communicator is 1, then
2044 // the new Map is not distributed.
2045 //
2046 // Otherwise, we have to check the new Map using an all-reduce
2047 // (over the new communicator). For example, the original Map
2048 // may have had some processes with zero elements, and all other
2049 // processes with the same number of elements as in the whole
2050 // Map. That Map is technically distributed, because of the
2051 // processes with zero elements. Removing those processes would
2052 // make the new Map locally replicated.
2053 if (!distributed_ || newComm->getSize() == 1) {
2054 map->distributed_ = false;
2055 } else {
2056 const int iOwnAllGids = (numLocalElements_ == numGlobalElements_) ? 1 : 0;
2057 int allProcsOwnAllGids = 0;
2059 map->distributed_ = (allProcsOwnAllGids == 1) ? false : true;
2060 }
2061
2062 map->lgMap_ = lgMap_;
2063 map->lgMapHost_ = lgMapHost_;
2064 map->glMap_ = glMap_;
2065 map->glMapHost_ = glMapHost_;
2066
2067 // Map's default constructor creates an uninitialized Directory.
2068 // The Directory will be initialized on demand in
2069 // getRemoteIndexList().
2070 //
2071 // FIXME (mfh 26 Mar 2013) It should be possible to "filter" the
2072 // directory more efficiently than just recreating it. If
2073 // directory recreation proves a bottleneck, we can always
2074 // revisit this. On the other hand, Directory creation is only
2075 // collective over the new, presumably much smaller
2076 // communicator, so it may not be worth the effort to optimize.
2077
2078 return map;
2079 }
2080}
2081
2082template <class LocalOrdinal, class GlobalOrdinal, class Node>
2085 directory_.is_null(), std::logic_error,
2086 "Tpetra::Map::setupDirectory: "
2087 "The Directory is null. "
2088 "Please report this bug to the Tpetra developers.");
2089
2090 // Only create the Directory if it hasn't been created yet.
2091 // This is a collective operation.
2092 if (!directory_->initialized()) {
2093 directory_->initialize(*this);
2094 }
2095}
2096
2097template <class LocalOrdinal, class GlobalOrdinal, class Node>
2100 getRemoteIndexList(const Teuchos::ArrayView<const GlobalOrdinal>& GIDs,
2101 const Teuchos::ArrayView<int>& PIDs,
2102 const Teuchos::ArrayView<LocalOrdinal>& LIDs) const {
2104 using std::endl;
2105 using Tpetra::Details::OrdinalTraits;
2106 using size_type = Teuchos::ArrayView<int>::size_type;
2107
2108 const bool verbose = Details::Behavior::verbose("Map");
2110 std::unique_ptr<std::string> prefix;
2111 if (verbose) {
2112 prefix = Details::createPrefix(comm_.getRawPtr(),
2113 "Map", "getRemoteIndexList(GIDs,PIDs,LIDs)");
2114 std::ostringstream os;
2115 os << *prefix << "Start: ";
2116 verbosePrintArray(os, GIDs, "GIDs", maxNumToPrint);
2117 os << endl;
2118 std::cerr << os.str();
2119 }
2120
2121 // Empty Maps (i.e., containing no indices on any processes in the
2122 // Map's communicator) are perfectly valid. In that case, if the
2123 // input GID list is nonempty, we fill the output arrays with
2124 // invalid values, and return IDNotPresent to notify the caller.
2125 // It's perfectly valid to give getRemoteIndexList GIDs that the
2126 // Map doesn't own. SubmapImport test 2 needs this functionality.
2127 if (getGlobalNumElements() == 0) {
2128 if (GIDs.size() == 0) {
2129 if (verbose) {
2130 std::ostringstream os;
2131 os << *prefix << "Done; both Map & input are empty" << endl;
2132 std::cerr << os.str();
2133 }
2134 return AllIDsPresent; // trivially
2135 } else {
2136 if (verbose) {
2137 std::ostringstream os;
2138 os << *prefix << "Done: Map is empty on all processes, "
2139 "so all output PIDs & LIDs are invalid (-1)."
2140 << endl;
2141 std::cerr << os.str();
2142 }
2143 for (size_type k = 0; k < PIDs.size(); ++k) {
2145 }
2146 for (size_type k = 0; k < LIDs.size(); ++k) {
2148 }
2149 return IDNotPresent;
2150 }
2151 }
2152
2153 // getRemoteIndexList must be called collectively, and Directory
2154 // initialization is collective too, so it's OK to initialize the
2155 // Directory on demand.
2156
2157 if (verbose) {
2158 std::ostringstream os;
2159 os << *prefix << "Call setupDirectory" << endl;
2160 std::cerr << os.str();
2161 }
2162 setupDirectory();
2163 if (verbose) {
2164 std::ostringstream os;
2165 os << *prefix << "Call directory_->getDirectoryEntries" << endl;
2166 std::cerr << os.str();
2167 }
2169 directory_->getDirectoryEntries(*this, GIDs, PIDs, LIDs);
2170 if (verbose) {
2171 std::ostringstream os;
2172 os << *prefix << "Done; getDirectoryEntries returned "
2173 << (retVal == IDNotPresent ? "IDNotPresent" : "AllIDsPresent")
2174 << "; ";
2175 verbosePrintArray(os, PIDs, "PIDs", maxNumToPrint);
2176 os << ", ";
2177 verbosePrintArray(os, LIDs, "LIDs", maxNumToPrint);
2178 os << endl;
2179 std::cerr << os.str();
2180 }
2181 return retVal;
2182}
2183
2184template <class LocalOrdinal, class GlobalOrdinal, class Node>
2187 getRemoteIndexList(const Teuchos::ArrayView<const GlobalOrdinal>& GIDs,
2188 const Teuchos::ArrayView<int>& PIDs) const {
2190 using std::endl;
2191
2192 const bool verbose = Details::Behavior::verbose("Map");
2194 std::unique_ptr<std::string> prefix;
2195 if (verbose) {
2196 prefix = Details::createPrefix(comm_.getRawPtr(),
2197 "Map", "getRemoteIndexList(GIDs,PIDs)");
2198 std::ostringstream os;
2199 os << *prefix << "Start: ";
2200 verbosePrintArray(os, GIDs, "GIDs", maxNumToPrint);
2201 os << endl;
2202 std::cerr << os.str();
2203 }
2204
2205 if (getGlobalNumElements() == 0) {
2206 if (GIDs.size() == 0) {
2207 if (verbose) {
2208 std::ostringstream os;
2209 os << *prefix << "Done; both Map & input are empty" << endl;
2210 std::cerr << os.str();
2211 }
2212 return AllIDsPresent; // trivially
2213 } else {
2214 if (verbose) {
2215 std::ostringstream os;
2216 os << *prefix << "Done: Map is empty on all processes, "
2217 "so all output PIDs are invalid (-1)."
2218 << endl;
2219 std::cerr << os.str();
2220 }
2221 for (Teuchos::ArrayView<int>::size_type k = 0; k < PIDs.size(); ++k) {
2222 PIDs[k] = Tpetra::Details::OrdinalTraits<int>::invalid();
2223 }
2224 return IDNotPresent;
2225 }
2226 }
2227
2228 // getRemoteIndexList must be called collectively, and Directory
2229 // initialization is collective too, so it's OK to initialize the
2230 // Directory on demand.
2231
2232 if (verbose) {
2233 std::ostringstream os;
2234 os << *prefix << "Call setupDirectory" << endl;
2235 std::cerr << os.str();
2236 }
2237 setupDirectory();
2238 if (verbose) {
2239 std::ostringstream os;
2240 os << *prefix << "Call directory_->getDirectoryEntries" << endl;
2241 std::cerr << os.str();
2242 }
2244 directory_->getDirectoryEntries(*this, GIDs, PIDs);
2245 if (verbose) {
2246 std::ostringstream os;
2247 os << *prefix << "Done; getDirectoryEntries returned "
2248 << (retVal == IDNotPresent ? "IDNotPresent" : "AllIDsPresent")
2249 << "; ";
2250 verbosePrintArray(os, PIDs, "PIDs", maxNumToPrint);
2251 os << endl;
2252 std::cerr << os.str();
2253 }
2254 return retVal;
2255}
2256
2257template <class LocalOrdinal, class GlobalOrdinal, class Node>
2259 using exec_space = typename Node::device_type::execution_space;
2260 if (lgMap_.extent(0) != lgMapHost_.extent(0)) {
2261 Tpetra::Details::ProfilingRegion pr("Map::lazyPushToHost() - pushing data");
2262 // NOTE: We check lgMap_ and not glMap_, since the latter can
2263 // be somewhat error prone for contiguous maps
2264
2265 // create_mirror_view preserves const-ness. create_mirror does not
2266 auto lgMap_host = Kokkos::create_mirror(Kokkos::HostSpace(), lgMap_);
2267
2268 // Since this was computed on the default stream, we can copy on the stream and then fence
2269 // the stream
2270 Kokkos::deep_copy(exec_space(), lgMap_host, lgMap_);
2271 exec_space().fence();
2272 lgMapHost_ = lgMap_host;
2273
2274 // Make host version - when memory spaces match these just do trivial assignment
2275 glMapHost_ = global_to_local_table_host_type(glMap_);
2276 }
2277}
2278
2279template <class LocalOrdinal, class GlobalOrdinal, class Node>
2280Teuchos::RCP<const Teuchos::Comm<int>>
2282 return comm_;
2283}
2284
2285template <class LocalOrdinal, class GlobalOrdinal, class Node>
2287 checkIsDist() const {
2288 using std::endl;
2289 using Teuchos::as;
2290 using Teuchos::outArg;
2291 using Teuchos::REDUCE_MIN;
2292 using Teuchos::reduceAll;
2293
2294 const bool verbose = Details::Behavior::verbose("Map");
2295 std::unique_ptr<std::string> prefix;
2296 if (verbose) {
2298 comm_.getRawPtr(), "Map", "checkIsDist");
2299 std::ostringstream os;
2300 os << *prefix << "Start" << endl;
2301 std::cerr << os.str();
2302 }
2303
2304 bool global = false;
2305 if (comm_->getSize() > 1) {
2306 // The communicator has more than one process, but that doesn't
2307 // necessarily mean the Map is distributed.
2308 int localRep = 0;
2309 if (numGlobalElements_ == as<global_size_t>(numLocalElements_)) {
2310 // The number of local elements on this process equals the
2311 // number of global elements.
2312 //
2313 // NOTE (mfh 22 Nov 2011) Does this still work if there were
2314 // duplicates in the global ID list on input (the third Map
2315 // constructor), so that the number of local elements (which
2316 // are not duplicated) on this process could be less than the
2317 // number of global elements, even if this process owns all
2318 // the elements?
2319 localRep = 1;
2320 }
2321 int allLocalRep;
2322 reduceAll<int, int>(*comm_, REDUCE_MIN, localRep, outArg(allLocalRep));
2323 if (allLocalRep != 1) {
2324 // At least one process does not own all the elements.
2325 // This makes the Map a distributed Map.
2326 global = true;
2327 }
2328 }
2329 // If the communicator has only one process, then the Map is not
2330 // distributed.
2331
2332 if (verbose) {
2333 std::ostringstream os;
2334 os << *prefix << "Done; global=" << (global ? "true" : "false")
2335 << endl;
2336 std::cerr << os.str();
2337 }
2338 return global;
2339}
2340
2341} // namespace Tpetra
2342
2343template <class LocalOrdinal, class GlobalOrdinal>
2344Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal>>
2345Tpetra::createLocalMap(const size_t numElements,
2346 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) {
2347 typedef LocalOrdinal LO;
2348 typedef GlobalOrdinal GO;
2349 using NT = typename ::Tpetra::Map<LO, GO>::node_type;
2350 return createLocalMapWithNode<LO, GO, NT>(numElements, comm);
2351}
2352
2353template <class LocalOrdinal, class GlobalOrdinal>
2354Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal>>
2355Tpetra::createUniformContigMap(const global_size_t numElements,
2356 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) {
2357 typedef LocalOrdinal LO;
2358 typedef GlobalOrdinal GO;
2359 using NT = typename ::Tpetra::Map<LO, GO>::node_type;
2360 return createUniformContigMapWithNode<LO, GO, NT>(numElements, comm);
2361}
2362
2363template <class LocalOrdinal, class GlobalOrdinal, class Node>
2364Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
2365Tpetra::createUniformContigMapWithNode(const global_size_t numElements,
2366 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) {
2367 using Teuchos::rcp;
2369 const GlobalOrdinal indexBase = static_cast<GlobalOrdinal>(0);
2370
2371 return rcp(new map_type(numElements, indexBase, comm, GloballyDistributed));
2372}
2373
2374template <class LocalOrdinal, class GlobalOrdinal, class Node>
2375Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
2376Tpetra::createLocalMapWithNode(const size_t numElements,
2377 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) {
2378 using Teuchos::rcp;
2381 const GlobalOrdinal indexBase = 0;
2382 const global_size_t globalNumElts = static_cast<global_size_t>(numElements);
2383
2384 return rcp(new map_type(globalNumElts, indexBase, comm, LocallyReplicated));
2385}
2386
2387template <class LocalOrdinal, class GlobalOrdinal, class Node>
2388Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
2390 const size_t localNumElements,
2391 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) {
2392 using Teuchos::rcp;
2394 const GlobalOrdinal indexBase = 0;
2395
2396 return rcp(new map_type(numElements, localNumElements, indexBase, comm));
2397}
2398
2399template <class LocalOrdinal, class GlobalOrdinal>
2400Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal>>
2402 const size_t localNumElements,
2403 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) {
2404 typedef LocalOrdinal LO;
2405 typedef GlobalOrdinal GO;
2406 using NT = typename Tpetra::Map<LO, GO>::node_type;
2407
2408 return Tpetra::createContigMapWithNode<LO, GO, NT>(numElements, localNumElements, comm);
2409}
2410
2411template <class LocalOrdinal, class GlobalOrdinal>
2412Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal>>
2413Tpetra::createNonContigMap(const Teuchos::ArrayView<const GlobalOrdinal>& elementList,
2414 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) {
2415 typedef LocalOrdinal LO;
2416 typedef GlobalOrdinal GO;
2417 using NT = typename Tpetra::Map<LO, GO>::node_type;
2418
2419 return Tpetra::createNonContigMapWithNode<LO, GO, NT>(elementList, comm);
2420}
2421
2422template <class LocalOrdinal, class GlobalOrdinal, class Node>
2423Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
2424Tpetra::createNonContigMapWithNode(const Teuchos::ArrayView<const GlobalOrdinal>& elementList,
2425 const Teuchos::RCP<const Teuchos::Comm<int>>& comm) {
2426 using Teuchos::rcp;
2428 using GST = Tpetra::global_size_t;
2429 const GST INV = Tpetra::Details::OrdinalTraits<GST>::invalid();
2430 // FIXME (mfh 22 Jul 2016) This is what I found here, but maybe this
2431 // shouldn't be zero, given that the index base is supposed to equal
2432 // the globally min global index?
2433 const GlobalOrdinal indexBase = 0;
2434
2435 return rcp(new map_type(INV, elementList, indexBase, comm));
2436}
2437
2438template <class LO, class GO, class NT>
2439Teuchos::RCP<const Tpetra::Map<LO, GO, NT>>
2440Tpetra::createOneToOne(const Teuchos::RCP<const Tpetra::Map<LO, GO, NT>>& M) {
2441 using Details::verbosePrintArray;
2442 using std::cerr;
2443 using std::endl;
2444 using Teuchos::Array;
2445 using Teuchos::ArrayView;
2446 using Teuchos::as;
2447 using Teuchos::rcp;
2448 using map_type = Tpetra::Map<LO, GO, NT>;
2449 using GST = global_size_t;
2450
2451 const bool verbose = Details::Behavior::verbose("Map");
2452 std::unique_ptr<std::string> prefix;
2453 if (verbose) {
2454 auto comm = M.is_null() ? Teuchos::null : M->getComm();
2455 prefix = Details::createPrefix(
2456 comm.getRawPtr(), "createOneToOne(Map)");
2457 std::ostringstream os;
2458 os << *prefix << "Start" << endl;
2459 cerr << os.str();
2460 }
2461 const size_t maxNumToPrint = verbose ? Details::Behavior::verbosePrintCountThreshold() : size_t(0);
2462 const GST GINV = Tpetra::Details::OrdinalTraits<GST>::invalid();
2463 const int myRank = M->getComm()->getRank();
2464
2465 // Bypasses for special cases where either M is known to be
2466 // one-to-one, or the one-to-one version of M is easy to compute.
2467 // This is why we take M as an RCP, not as a const reference -- so
2468 // that we can return M itself if it is 1-to-1.
2469 if (!M->isDistributed()) {
2470 // For a locally replicated Map, we assume that users want to push
2471 // all the GIDs to Process 0.
2472
2473 // mfh 05 Nov 2013: getGlobalNumElements() does indeed return what
2474 // you think it should return, in this special case of a locally
2475 // replicated contiguous Map.
2476 const GST numGlobalEntries = M->getGlobalNumElements();
2477 if (M->isContiguous()) {
2478 const size_t numLocalEntries =
2479 (myRank == 0) ? as<size_t>(numGlobalEntries) : size_t(0);
2480 if (verbose) {
2481 std::ostringstream os;
2482 os << *prefix << "Input is locally replicated & contiguous; "
2483 "numLocalEntries="
2484 << numLocalEntries << endl;
2485 cerr << os.str();
2486 }
2487 auto retMap =
2488 rcp(new map_type(numGlobalEntries, numLocalEntries,
2489 M->getIndexBase(), M->getComm()));
2490 if (verbose) {
2491 std::ostringstream os;
2492 os << *prefix << "Done" << endl;
2493 cerr << os.str();
2494 }
2495 return retMap;
2496 } else {
2497 if (verbose) {
2498 std::ostringstream os;
2499 os << *prefix << "Input is locally replicated & noncontiguous"
2500 << endl;
2501 cerr << os.str();
2502 }
2503 ArrayView<const GO> myGids =
2504 (myRank == 0) ? M->getLocalElementList() : Teuchos::null;
2505 auto retMap =
2506 rcp(new map_type(GINV, myGids(), M->getIndexBase(),
2507 M->getComm()));
2508 if (verbose) {
2509 std::ostringstream os;
2510 os << *prefix << "Done" << endl;
2511 cerr << os.str();
2512 }
2513 return retMap;
2514 }
2515 } else if (M->isContiguous()) {
2516 if (verbose) {
2517 std::ostringstream os;
2518 os << *prefix << "Input is distributed & contiguous" << endl;
2519 cerr << os.str();
2520 }
2521 // Contiguous, distributed Maps are one-to-one by construction.
2522 // (Locally replicated Maps can be contiguous.)
2523 return M;
2524 } else {
2525 if (verbose) {
2526 std::ostringstream os;
2527 os << *prefix << "Input is distributed & noncontiguous" << endl;
2528 cerr << os.str();
2529 }
2531 const size_t numMyElems = M->getLocalNumElements();
2532 ArrayView<const GO> myElems = M->getLocalElementList();
2533 Array<int> owner_procs_vec(numMyElems);
2534
2535 if (verbose) {
2536 std::ostringstream os;
2537 os << *prefix << "Call Directory::getDirectoryEntries: ";
2538 verbosePrintArray(os, myElems, "GIDs", maxNumToPrint);
2539 os << endl;
2540 cerr << os.str();
2541 }
2542 directory.getDirectoryEntries(*M, myElems, owner_procs_vec());
2543 if (verbose) {
2544 std::ostringstream os;
2545 os << *prefix << "getDirectoryEntries result: ";
2546 verbosePrintArray(os, owner_procs_vec, "PIDs", maxNumToPrint);
2547 os << endl;
2548 cerr << os.str();
2549 }
2550
2551 Array<GO> myOwned_vec(numMyElems);
2552 size_t numMyOwnedElems = 0;
2553 for (size_t i = 0; i < numMyElems; ++i) {
2554 const GO GID = myElems[i];
2555 const int owner = owner_procs_vec[i];
2556
2557 if (myRank == owner) {
2558 myOwned_vec[numMyOwnedElems++] = GID;
2559 }
2560 }
2561 myOwned_vec.resize(numMyOwnedElems);
2562
2563 if (verbose) {
2564 std::ostringstream os;
2565 os << *prefix << "Create Map: ";
2566 verbosePrintArray(os, myOwned_vec, "GIDs", maxNumToPrint);
2567 os << endl;
2568 cerr << os.str();
2569 }
2570 auto retMap = rcp(new map_type(GINV, myOwned_vec(),
2571 M->getIndexBase(), M->getComm()));
2572 if (verbose) {
2573 std::ostringstream os;
2574 os << *prefix << "Done" << endl;
2575 cerr << os.str();
2576 }
2577 return retMap;
2578 }
2579}
2580
2581template <class LocalOrdinal, class GlobalOrdinal, class Node>
2582Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>
2585 using Details::Behavior;
2586 using Details::verbosePrintArray;
2587 using std::cerr;
2588 using std::endl;
2589 using Teuchos::Array;
2590 using Teuchos::ArrayView;
2591 using Teuchos::RCP;
2592 using Teuchos::rcp;
2593 using Teuchos::toString;
2594 using LO = LocalOrdinal;
2595 using GO = GlobalOrdinal;
2596 using map_type = Tpetra::Map<LO, GO, Node>;
2597
2598 const bool verbose = Behavior::verbose("Map");
2599 std::unique_ptr<std::string> prefix;
2600 if (verbose) {
2601 auto comm = M.is_null() ? Teuchos::null : M->getComm();
2602 prefix = Details::createPrefix(
2603 comm.getRawPtr(), "createOneToOne(Map,TieBreak)");
2604 std::ostringstream os;
2605 os << *prefix << "Start" << endl;
2606 cerr << os.str();
2607 }
2608 const size_t maxNumToPrint = verbose ? Behavior::verbosePrintCountThreshold() : size_t(0);
2609
2610 // FIXME (mfh 20 Feb 2013) We should have a bypass for contiguous
2611 // Maps (which are 1-to-1 by construction).
2612
2614 if (verbose) {
2615 std::ostringstream os;
2616 os << *prefix << "Initialize Directory" << endl;
2617 cerr << os.str();
2618 }
2619 directory.initialize(*M, tie_break);
2620 if (verbose) {
2621 std::ostringstream os;
2622 os << *prefix << "Done initializing Directory" << endl;
2623 cerr << os.str();
2624 }
2625 size_t numMyElems = M->getLocalNumElements();
2626 ArrayView<const GO> myElems = M->getLocalElementList();
2627 Array<int> owner_procs_vec(numMyElems);
2628 if (verbose) {
2629 std::ostringstream os;
2630 os << *prefix << "Call Directory::getDirectoryEntries: ";
2631 verbosePrintArray(os, myElems, "GIDs", maxNumToPrint);
2632 os << endl;
2633 cerr << os.str();
2634 }
2635 directory.getDirectoryEntries(*M, myElems, owner_procs_vec());
2636 if (verbose) {
2637 std::ostringstream os;
2638 os << *prefix << "getDirectoryEntries result: ";
2639 verbosePrintArray(os, owner_procs_vec, "PIDs", maxNumToPrint);
2640 os << endl;
2641 cerr << os.str();
2642 }
2643
2644 const int myRank = M->getComm()->getRank();
2645 Array<GO> myOwned_vec(numMyElems);
2646 size_t numMyOwnedElems = 0;
2647 for (size_t i = 0; i < numMyElems; ++i) {
2648 const GO GID = myElems[i];
2649 const int owner = owner_procs_vec[i];
2650 if (myRank == owner) {
2651 myOwned_vec[numMyOwnedElems++] = GID;
2652 }
2653 }
2654 myOwned_vec.resize(numMyOwnedElems);
2655
2656 // FIXME (mfh 08 May 2014) The above Directory should be perfectly
2657 // valid for the new Map. Why can't we reuse it?
2658 const global_size_t GINV =
2659 Tpetra::Details::OrdinalTraits<global_size_t>::invalid();
2660 if (verbose) {
2661 std::ostringstream os;
2662 os << *prefix << "Create Map: ";
2663 verbosePrintArray(os, myOwned_vec, "GIDs", maxNumToPrint);
2664 os << endl;
2665 cerr << os.str();
2666 }
2667 RCP<const map_type> retMap(new map_type(GINV, myOwned_vec(), M->getIndexBase(),
2668 M->getComm()));
2669 if (verbose) {
2670 std::ostringstream os;
2671 os << *prefix << "Done" << endl;
2672 cerr << os.str();
2673 }
2674 return retMap;
2675}
2676
2677namespace Tpetra::Details {
2678
2679template <class pids_view_type>
2680struct SortToFit {
2681 int myRank;
2682 pids_view_type pids;
2683
2684 SortToFit(int _myRank, pids_view_type _pids)
2685 : myRank(_myRank)
2686 , pids(_pids) {}
2687
2688 KOKKOS_FUNCTION
2689 bool operator()(size_t i, size_t j) const {
2690 if (pids(i) == myRank) {
2691 if (pids(j) == myRank)
2692 return i < j;
2693 else
2694 return true;
2695 } else {
2696 if (pids(j) == myRank)
2697 return false;
2698 else
2699 return i < j;
2700 }
2701 }
2702};
2703
2704} // namespace Tpetra::Details
2705
2706template <class LO, class GO, class NT>
2707Teuchos::RCP<const Tpetra::Map<LO, GO, NT>>
2709 Teuchos::RCP<const Tpetra::Map<LO, GO, NT>> constM = M;
2711
2712 const int locallyFitted = M->isLocallyFitted(*owned_node_map);
2713 int globallyLocallyFitted = 0;
2714 reduceAll(*M->getComm(), Teuchos::REDUCE_MIN, locallyFitted, Teuchos::outArg(globallyLocallyFitted));
2715
2716 if (globallyLocallyFitted == 0) {
2717 using pid_type = int;
2718
2719 const pid_type myRank = M->getComm()->getRank();
2720 const size_t numMyElems = M->getLocalNumElements();
2721
2722 Kokkos::View<pid_type*, typename NT::memory_space> pids("pids", numMyElems);
2723 {
2724 auto gids_vec = M->getLocalElementList();
2725 Teuchos::Array<pid_type> pids_vec(numMyElems);
2726 M->getRemoteIndexList(gids_vec, pids_vec);
2727 Kokkos::View<pid_type*, Kokkos::HostSpace, Kokkos::MemoryTraits<Kokkos::Unmanaged>> pids_h(pids_vec.data(), pids_vec.size());
2728 Kokkos::deep_copy(pids, pids_h);
2729 }
2730
2731 auto gids = M->getMyGlobalIndicesDevice();
2732
2733 auto policy = Kokkos::RangePolicy<size_t, typename NT::execution_space>(0, numMyElems);
2734
2735 Kokkos::View<size_t*, typename NT::memory_space> idx(Kokkos::ViewAllocateWithoutInitializing("idx"), numMyElems);
2736 Kokkos::parallel_for(
2737 policy, KOKKOS_LAMBDA(const size_t i) { idx(i) = i; });
2738
2739 Tpetra::Details::SortToFit cmp(myRank, pids);
2740 Kokkos::sort(typename NT::execution_space(), idx, cmp);
2741
2742 Kokkos::View<GO*, typename NT::memory_space> new_gids(Kokkos::ViewAllocateWithoutInitializing("new_gids"), numMyElems);
2743 Kokkos::parallel_for(
2744 policy, KOKKOS_LAMBDA(const size_t i) { new_gids(i) = gids(idx(i)); });
2745
2746 Teuchos::RCP<const Tpetra::Map<LO, GO, NT>> shared_node_map =
2747 Teuchos::rcp(new Tpetra::Map<LO, GO, NT>(M->getGlobalNumElements(),
2748 new_gids,
2749 M->getIndexBase(),
2750 M->getComm()));
2751 M.swap(shared_node_map);
2752 }
2753
2754 return owned_node_map;
2755}
2756
2757//
2758// Explicit instantiation macro
2759//
2760// Must be expanded from within the Tpetra namespace!
2761//
2762
2764
2765#define TPETRA_MAP_INSTANT(LO, GO, NODE) \
2766 \
2767 template class Map<LO, GO, NODE>; \
2768 \
2769 template Teuchos::RCP<const Map<LO, GO, NODE>> \
2770 createLocalMapWithNode<LO, GO, NODE>(const size_t numElements, \
2771 const Teuchos::RCP<const Teuchos::Comm<int>>& comm); \
2772 \
2773 template Teuchos::RCP<const Map<LO, GO, NODE>> \
2774 createContigMapWithNode<LO, GO, NODE>(const global_size_t numElements, \
2775 const size_t localNumElements, \
2776 const Teuchos::RCP<const Teuchos::Comm<int>>& comm); \
2777 \
2778 template Teuchos::RCP<const Map<LO, GO, NODE>> \
2779 createNonContigMapWithNode(const Teuchos::ArrayView<const GO>& elementList, \
2780 const Teuchos::RCP<const Teuchos::Comm<int>>& comm); \
2781 \
2782 template Teuchos::RCP<const Map<LO, GO, NODE>> \
2783 createUniformContigMapWithNode<LO, GO, NODE>(const global_size_t numElements, \
2784 const Teuchos::RCP<const Teuchos::Comm<int>>& comm); \
2785 \
2786 template Teuchos::RCP<const Map<LO, GO, NODE>> \
2787 createOneToOne(const Teuchos::RCP<const Map<LO, GO, NODE>>& M); \
2788 \
2789 template Teuchos::RCP<const Map<LO, GO, NODE>> \
2790 createOneToOne(const Teuchos::RCP<const Map<LO, GO, NODE>>& M, \
2791 const Tpetra::Details::TieBreak<LO, GO>& tie_break); \
2792 \
2793 template Teuchos::RCP<const Map<LO, GO, NODE>> \
2794 createOneToOneAndMakeOverlappingMapFitted(Teuchos::RCP<const Map<LO, GO, NODE>>& M);
2795
2797#define TPETRA_MAP_INSTANT_DEFAULTNODE(LO, GO) \
2798 template Teuchos::RCP<const Map<LO, GO>> \
2799 createLocalMap<LO, GO>(const size_t, const Teuchos::RCP<const Teuchos::Comm<int>>&); \
2800 \
2801 template Teuchos::RCP<const Map<LO, GO>> \
2802 createContigMap<LO, GO>(global_size_t, size_t, \
2803 const Teuchos::RCP<const Teuchos::Comm<int>>&); \
2804 \
2805 template Teuchos::RCP<const Map<LO, GO>> \
2806 createNonContigMap(const Teuchos::ArrayView<const GO>&, \
2807 const Teuchos::RCP<const Teuchos::Comm<int>>&); \
2808 \
2809 template Teuchos::RCP<const Map<LO, GO>> \
2810 createUniformContigMap<LO, GO>(const global_size_t, \
2811 const Teuchos::RCP<const Teuchos::Comm<int>>&);
2812
2813#endif // TPETRA_MAP_DEF_HPP
Functions for initializing and finalizing Tpetra.
Declaration of Tpetra::Details::Behavior, a class that describes Tpetra's behavior.
Declaration of Tpetra::Details::Profiling, a scope guard for Kokkos Profiling.
Declaration of Tpetra::Details::extractMpiCommFromTeuchos.
Declaration of a function that prints strings from each process.
Declaration of Tpetra::Details::iallreduce.
Declaration of Tpetra::Details::initializeKokkos.
Declaration of Tpetra::Details::printOnce.
Declaration of the Tpetra::Map class and related nonmember constructors.
Stand-alone utility functions and macros.
Struct that holds views of the contents of a CrsMatrix.
Description of Tpetra's behavior.
static void reject_unrecognized_env_vars()
Search the environment for TPETRA_ variables and reject unrecognized ones.
static bool debug()
Whether Tpetra is in debug mode.
static bool verbose()
Whether Tpetra is in verbose mode.
static size_t verbosePrintCountThreshold()
Number of entries below which arrays, lists, etc. will be printed in debug mode.
"Local" part of Map suitable for Kokkos kernels.
Interface for breaking ties in ownership.
Implement mapping from global ID to process ID and local ID.
A parallel distribution of indices over processes.
bool isDistributed() const
Whether this Map is globally distributed or locally replicated.
bool isOneToOne() const
Whether the Map is one to one.
std::string description() const
Implementation of Teuchos::Describable.
Teuchos::ArrayView< const global_ordinal_type > getLocalElementList() const
Return a NONOWNING view of the global indices owned by this process.
::Tpetra::Details::FixedHashTable< global_ordinal_type, local_ordinal_type, device_type > global_to_local_table_type
Type of a mapping from global IDs to local IDs.
global_ordinal_type getGlobalElement(local_ordinal_type localIndex) const
The global index corresponding to the given local index.
Node node_type
Legacy typedef that will go away at some point.
Kokkos::View< const global_ordinal_type *, Kokkos::LayoutLeft, device_type > lgMap_
A mapping from local IDs to global IDs.
Map()
Default constructor (that does nothing).
GlobalOrdinal global_ordinal_type
The type of global indices.
LookupStatus getRemoteIndexList(const Teuchos::ArrayView< const global_ordinal_type > &GIDList, const Teuchos::ArrayView< int > &nodeIDList, const Teuchos::ArrayView< local_ordinal_type > &LIDList) const
Return the process ranks and corresponding local indices for the given global indices.
bool isNodeLocalElement(local_ordinal_type localIndex) const
Whether the given local index is valid for this Map on the calling process.
bool isUniform() const
Whether the range of global indices is uniform.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Accessors for the Teuchos::Comm and Kokkos Node objects.
typename device_type::execution_space execution_space
The Kokkos execution space.
LO local_ordinal_type
The type of local indices.
Teuchos::RCP< const Map< local_ordinal_type, global_ordinal_type, Node > > removeEmptyProcesses() const
Advanced methods.
void lazyPushToHost() const
Push the device data to host, if needed.
bool isCompatible(const Map< local_ordinal_type, global_ordinal_type, Node > &map) const
True if and only if map is compatible with this Map.
bool locallySameAs(const Map< local_ordinal_type, global_ordinal_type, node_type > &map) const
Is this Map locally the same as the input Map?
bool isLocallyFitted(const Map< local_ordinal_type, global_ordinal_type, Node > &map) const
True if and only if map is locally fitted to this Map.
virtual ~Map()
Destructor (virtual for memory safety of derived classes).
global_indices_array_device_type getMyGlobalIndicesDevice() const
Return a view of the global indices owned by this process on the Map's device.
local_ordinal_type getLocalElement(global_ordinal_type globalIndex) const
The local index corresponding to the given global index.
Teuchos::RCP< const Map< local_ordinal_type, global_ordinal_type, Node > > replaceCommWithSubset(const Teuchos::RCP< const Teuchos::Comm< int > > &newComm) const
Replace this Map's communicator with a subset communicator.
bool isContiguous() const
True if this Map is distributed contiguously, else false.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel=Teuchos::Describable::verbLevel_default) const
Describe this object in a human-readable way to the given output stream.
void computeGlobalConstants()
Compute global constants for the map. This method needs to be called collectively over all processes ...
bool isNodeGlobalElement(global_ordinal_type globalIndex) const
Whether the given global index is owned by this Map on the calling process.
bool isSameAs(const Map< local_ordinal_type, global_ordinal_type, Node > &map) const
True if and only if map is identical to this Map.
global_to_local_table_type glMap_
A mapping from global IDs to local IDs.
local_map_type getLocalMap() const
Get the LocalMap for Kokkos-Kernels.
global_indices_array_type getMyGlobalIndices() const
Return a view of the global indices owned by this process.
typename Node::device_type device_type
This class' Kokkos::Device specialization.
Implementation details of Tpetra.
Nonmember function that computes a residual Computes R = B - A * X.
void verbosePrintArray(std::ostream &out, const ArrayType &x, const char name[], const size_t maxNumToPrint)
Print min(x.size(), maxNumToPrint) entries of x.
void printOnce(std::ostream &out, const std::string &s, const Teuchos::Comm< int > *comm)
Print on one process of the given communicator, or at least try to do so (if MPI is not initialized).
std::unique_ptr< std::string > createPrefix(const int myRank, const char prefix[])
Create string prefix for each line of verbose output.
bool congruent(const Teuchos::Comm< int > &comm1, const Teuchos::Comm< int > &comm2)
Whether the two communicators are congruent.
void initializeKokkos()
Initialize Kokkos, using command-line arguments (if any) given to Teuchos::GlobalMPISession.
void gathervPrint(std::ostream &out, const std::string &s, const Teuchos::Comm< int > &comm)
On Process 0 in the given communicator, print strings from each process in that communicator,...
Namespace Tpetra contains the class and methods constituting the Tpetra library.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createLocalMap(const size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a locally replicated Map with the default Kokkos Node.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createUniformContigMap(const global_size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Non-member constructor for a uniformly distributed, contiguous Map with the default Kokkos Node.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createNonContigMap(const Teuchos::ArrayView< const GlobalOrdinal > &elementList, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a non-contiguous Map using the default Kokkos::Device type.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createUniformContigMapWithNode(const global_size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Non-member constructor for a uniformly distributed, contiguous Map with a user-specified Kokkos Node.
Teuchos::RCP< const Tpetra::Map< LO, GO, NT > > createOneToOneAndMakeOverlappingMapFitted(Teuchos::RCP< const Tpetra::Map< LO, GO, NT > > &M)
Creates a one-to-one version of the given Map where each GID lives on only one process....
LookupStatus
Return status of Map remote index lookup (getRemoteIndexList()).
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal > > createContigMap(const global_size_t numElements, const size_t localNumElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Non-member constructor for a (potentially) non-uniformly distributed, contiguous Map using the defaul...
size_t global_size_t
Global size_t object.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createOneToOne(const Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > &M)
Nonmember constructor for a contiguous Map with user-defined weights and a user-specified,...
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createLocalMapWithNode(const size_t numElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a locally replicated Map with a specified Kokkos Node.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createNonContigMapWithNode(const Teuchos::ArrayView< const GlobalOrdinal > &elementList, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a noncontiguous Map with a user-specified, possibly nondefault Kokkos Node ...
LocalGlobal
Enum for local versus global allocation of Map entries.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > createContigMapWithNode(const global_size_t numElements, const size_t localNumElements, const Teuchos::RCP< const Teuchos::Comm< int > > &comm)
Nonmember constructor for a (potentially) nonuniformly distributed, contiguous Map for a user-specifi...