Tpetra parallel linear algebra  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tpetra_Map_def.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ***********************************************************************
3 //
4 // Tpetra: Templated Linear Algebra Services Package
5 // Copyright (2008) Sandia Corporation
6 //
7 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
8 // the U.S. Government retains certain rights in this software.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
38 //
39 // ************************************************************************
40 // @HEADER
41 
46 
47 #ifndef TPETRA_MAP_DEF_HPP
48 #define TPETRA_MAP_DEF_HPP
49 
50 #include "Tpetra_Directory.hpp" // must include for implicit instantiation to work
51 #include "Tpetra_Details_FixedHashTable.hpp"
52 #include "Tpetra_Details_gathervPrint.hpp"
54 #include "Tpetra_Core.hpp"
55 #include "Tpetra_Util.hpp"
56 #include "Teuchos_as.hpp"
57 #include "Teuchos_TypeNameTraits.hpp"
58 #include "Teuchos_CommHelpers.hpp"
59 #include "Tpetra_Details_mpiIsInitialized.hpp"
60 #include "Tpetra_Details_extractMpiCommFromTeuchos.hpp" // teuchosCommIsAnMpiComm
62 #include <stdexcept>
63 #include <typeinfo>
64 
65 namespace Tpetra {
66 
67  template <class LocalOrdinal, class GlobalOrdinal, class Node>
69  Map () :
70  comm_ (new Teuchos::SerialComm<int> ()),
71  indexBase_ (0),
72  numGlobalElements_ (0),
73  numLocalElements_ (0),
74  minMyGID_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
75  maxMyGID_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
76  minAllGID_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
77  maxAllGID_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
78  firstContiguousGID_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
79  lastContiguousGID_ (Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ()),
80  uniform_ (false), // trivially
81  contiguous_ (false),
82  distributed_ (false), // no communicator yet
83  directory_ (new Directory<LocalOrdinal, GlobalOrdinal, Node> ())
84  {
86  }
87 
88 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
89  template <class LocalOrdinal, class GlobalOrdinal, class Node>
90  TPETRA_DEPRECATED
92  Map (const global_size_t numGlobalElements,
93  const global_ordinal_type indexBase,
94  const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
95  const LocalGlobal lOrG,
96  const Teuchos::RCP<Node> &/* node */) :
97  Map<LocalOrdinal,GlobalOrdinal,Node>::Map(numGlobalElements, indexBase,
98  comm, lOrG)
99  {}
100 #endif // TPETRA_ENABLE_DEPRECATED_CODE
101 
102  template <class LocalOrdinal, class GlobalOrdinal, class Node>
104  Map (const global_size_t numGlobalElements,
105  const global_ordinal_type indexBase,
106  const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
107  const LocalGlobal lOrG) :
108  comm_ (comm),
109  uniform_ (true),
110  directory_ (new Directory<LocalOrdinal, GlobalOrdinal, Node> ())
111  {
112  using Teuchos::as;
113  using Teuchos::broadcast;
114  using Teuchos::outArg;
115  using Teuchos::reduceAll;
116  using Teuchos::REDUCE_MIN;
117  using Teuchos::REDUCE_MAX;
118  using Teuchos::typeName;
119  using GO = global_ordinal_type;
120  using GST = global_size_t;
121  const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid ();
122 
124 
125 #ifdef HAVE_TPETRA_DEBUG
126  // In debug mode only, check whether numGlobalElements and
127  // indexBase are the same over all processes in the communicator.
128  {
129  GST proc0NumGlobalElements = numGlobalElements;
130  broadcast<int, GST> (*comm_, 0, outArg (proc0NumGlobalElements));
131  GST minNumGlobalElements = numGlobalElements;
132  GST maxNumGlobalElements = numGlobalElements;
133  reduceAll<int, GST> (*comm, REDUCE_MIN, numGlobalElements, outArg (minNumGlobalElements));
134  reduceAll<int, GST> (*comm, REDUCE_MAX, numGlobalElements, outArg (maxNumGlobalElements));
135  TEUCHOS_TEST_FOR_EXCEPTION(
136  minNumGlobalElements != maxNumGlobalElements || numGlobalElements != minNumGlobalElements,
137  std::invalid_argument,
138  "Tpetra::Map constructor: All processes must provide the same number "
139  "of global elements. Process 0 set numGlobalElements = "
140  << proc0NumGlobalElements << ". The calling process "
141  << comm->getRank () << " set numGlobalElements = " << numGlobalElements
142  << ". The min and max values over all processes are "
143  << minNumGlobalElements << " resp. " << maxNumGlobalElements << ".");
144 
145  GO proc0IndexBase = indexBase;
146  broadcast<int, GO> (*comm_, 0, outArg (proc0IndexBase));
147  GO minIndexBase = indexBase;
148  GO maxIndexBase = indexBase;
149  reduceAll<int, GO> (*comm, REDUCE_MIN, indexBase, outArg (minIndexBase));
150  reduceAll<int, GO> (*comm, REDUCE_MAX, indexBase, outArg (maxIndexBase));
151  TEUCHOS_TEST_FOR_EXCEPTION(
152  minIndexBase != maxIndexBase || indexBase != minIndexBase,
153  std::invalid_argument,
154  "Tpetra::Map constructor: "
155  "All processes must provide the same indexBase argument. "
156  "Process 0 set indexBase = " << proc0IndexBase << ". The calling "
157  "process " << comm->getRank () << " set indexBase = " << indexBase
158  << ". The min and max values over all processes are "
159  << minIndexBase << " resp. " << maxIndexBase << ".");
160  }
161 #endif // HAVE_TPETRA_DEBUG
162 
163  // Distribute the elements across the processes in the given
164  // communicator so that global IDs (GIDs) are
165  //
166  // - Nonoverlapping (only one process owns each GID)
167  // - Contiguous (the sequence of GIDs is nondecreasing, and no two
168  // adjacent GIDs differ by more than one)
169  // - As evenly distributed as possible (the numbers of GIDs on two
170  // different processes do not differ by more than one)
171 
172  // All processes have the same numGlobalElements, but we still
173  // need to check that it is valid. numGlobalElements must be
174  // positive and not the "invalid" value (GSTI).
175  //
176  // This comparison looks funny, but it avoids compiler warnings
177  // for comparing unsigned integers (numGlobalElements_in is a
178  // GST, which is unsigned) while still working if we
179  // later decide to make GST signed.
180  TEUCHOS_TEST_FOR_EXCEPTION(
181  (numGlobalElements < 1 && numGlobalElements != 0),
182  std::invalid_argument,
183  "Tpetra::Map constructor: numGlobalElements (= "
184  << numGlobalElements << ") must be nonnegative.");
185 
186  TEUCHOS_TEST_FOR_EXCEPTION(
187  numGlobalElements == GSTI, std::invalid_argument,
188  "Tpetra::Map constructor: You provided numGlobalElements = Teuchos::"
189  "OrdinalTraits<Tpetra::global_size_t>::invalid(). This version of the "
190  "constructor requires a valid value of numGlobalElements. You "
191  "probably mistook this constructor for the \"contiguous nonuniform\" "
192  "constructor, which can compute the global number of elements for you "
193  "if you set numGlobalElements to that value.");
194 
195  size_t numLocalElements = 0; // will set below
196  if (lOrG == GloballyDistributed) {
197  // Compute numLocalElements:
198  //
199  // If numGlobalElements == numProcs * B + remainder,
200  // then Proc r gets B+1 elements if r < remainder,
201  // and B elements if r >= remainder.
202  //
203  // This strategy is valid for any value of numGlobalElements and
204  // numProcs, including the following border cases:
205  // - numProcs == 1
206  // - numLocalElements < numProcs
207  //
208  // In the former case, remainder == 0 && numGlobalElements ==
209  // numLocalElements. In the latter case, remainder ==
210  // numGlobalElements && numLocalElements is either 0 or 1.
211  const GST numProcs = static_cast<GST> (comm_->getSize ());
212  const GST myRank = static_cast<GST> (comm_->getRank ());
213  const GST quotient = numGlobalElements / numProcs;
214  const GST remainder = numGlobalElements - quotient * numProcs;
215 
216  GO startIndex;
217  if (myRank < remainder) {
218  numLocalElements = static_cast<size_t> (1) + static_cast<size_t> (quotient);
219  // myRank was originally an int, so it should never overflow
220  // reasonable GO types.
221  startIndex = as<GO> (myRank) * as<GO> (numLocalElements);
222  } else {
223  numLocalElements = as<size_t> (quotient);
224  startIndex = as<GO> (myRank) * as<GO> (numLocalElements) +
225  as<GO> (remainder);
226  }
227 
228  minMyGID_ = indexBase + startIndex;
229  maxMyGID_ = indexBase + startIndex + numLocalElements - 1;
230  minAllGID_ = indexBase;
231  maxAllGID_ = indexBase + numGlobalElements - 1;
232  distributed_ = (numProcs > 1);
233  }
234  else { // lOrG == LocallyReplicated
235  numLocalElements = as<size_t> (numGlobalElements);
236  minMyGID_ = indexBase;
237  maxMyGID_ = indexBase + numGlobalElements - 1;
238  distributed_ = false;
239  }
240 
241  minAllGID_ = indexBase;
242  maxAllGID_ = indexBase + numGlobalElements - 1;
243  indexBase_ = indexBase;
244  numGlobalElements_ = numGlobalElements;
245  numLocalElements_ = numLocalElements;
246  firstContiguousGID_ = minMyGID_;
247  lastContiguousGID_ = maxMyGID_;
248  contiguous_ = true;
249 
250  // Create the Directory on demand in getRemoteIndexList().
251  //setupDirectory ();
252  }
253 
254 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
255  template <class LocalOrdinal, class GlobalOrdinal, class Node>
256  TPETRA_DEPRECATED
258  Map (const global_size_t numGlobalElements,
259  const size_t numLocalElements,
260  const global_ordinal_type indexBase,
261  const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
262  const Teuchos::RCP<Node> &/* node */) :
263  Map<LocalOrdinal,GlobalOrdinal,Node>::Map(numGlobalElements,
264  numLocalElements,
265  indexBase, comm)
266  {}
267 #endif // TPETRA_ENABLE_DEPRECATED_CODE
268 
269  template <class LocalOrdinal, class GlobalOrdinal, class Node>
271  Map (const global_size_t numGlobalElements,
272  const size_t numLocalElements,
273  const global_ordinal_type indexBase,
274  const Teuchos::RCP<const Teuchos::Comm<int> > &comm) :
275  comm_ (comm),
276  uniform_ (false),
277  directory_ (new Directory<LocalOrdinal, GlobalOrdinal, Node> ())
278  {
279  using Teuchos::as;
280  using Teuchos::broadcast;
281  using Teuchos::outArg;
282  using Teuchos::reduceAll;
283  using Teuchos::REDUCE_MIN;
284  using Teuchos::REDUCE_MAX;
285  using Teuchos::REDUCE_SUM;
286  using Teuchos::scan;
287  using GO = global_ordinal_type;
288  using GST = global_size_t;
289  const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid ();
290 
292 
293 #ifdef HAVE_TPETRA_DEBUG
294  // Global sum of numLocalElements over all processes.
295  // Keep this for later debug checks.
296  const GST debugGlobalSum =
297  initialNonuniformDebugCheck (numGlobalElements, numLocalElements,
298  indexBase, comm);
299 #endif // HAVE_TPETRA_DEBUG
300 
301  // Distribute the elements across the nodes so that they are
302  // - non-overlapping
303  // - contiguous
304 
305  // This differs from the first Map constructor (that only takes a
306  // global number of elements) in that the user has specified the
307  // number of local elements, so that the elements are not
308  // (necessarily) evenly distributed over the processes.
309 
310  // Compute my local offset. This is an inclusive scan, so to get
311  // the final offset, we subtract off the input.
312  GO scanResult = 0;
313  scan<int, GO> (*comm, REDUCE_SUM, numLocalElements, outArg (scanResult));
314  const GO myOffset = scanResult - numLocalElements;
315 
316  if (numGlobalElements != GSTI) {
317  numGlobalElements_ = numGlobalElements; // Use the user's value.
318  } else {
319  // Inclusive scan means that the last process has the final sum.
320  // Rather than doing a reduceAll to get the sum of
321  // numLocalElements, we can just have the last process broadcast
322  // its result. That saves us a round of log(numProcs) messages.
323  const int numProcs = comm->getSize ();
324  GST globalSum = scanResult;
325  if (numProcs > 1) {
326  broadcast (*comm, numProcs - 1, outArg (globalSum));
327  }
328  numGlobalElements_ = globalSum;
329 
330 #ifdef HAVE_TPETRA_DEBUG
331  // No need for an all-reduce here; both come from collectives.
332  TEUCHOS_TEST_FOR_EXCEPTION(
333  globalSum != debugGlobalSum, std::logic_error,
334  "Tpetra::Map constructor (contiguous nonuniform): "
335  "globalSum = " << globalSum << " != debugGlobalSum = " << debugGlobalSum
336  << ". Please report this bug to the Tpetra developers.");
337 #endif // HAVE_TPETRA_DEBUG
338  }
339  numLocalElements_ = numLocalElements;
340  indexBase_ = indexBase;
341  minAllGID_ = (numGlobalElements_ == 0) ?
342  std::numeric_limits<GO>::max () :
343  indexBase;
344  maxAllGID_ = (numGlobalElements_ == 0) ?
345  std::numeric_limits<GO>::lowest () :
346  indexBase + static_cast<GO> (numGlobalElements_) - static_cast<GO> (1);
347  minMyGID_ = (numLocalElements_ == 0) ?
348  std::numeric_limits<GO>::max () :
349  indexBase + static_cast<GO> (myOffset);
350  maxMyGID_ = (numLocalElements_ == 0) ?
351  std::numeric_limits<GO>::lowest () :
352  indexBase + myOffset + static_cast<GO> (numLocalElements) - static_cast<GO> (1);
353  firstContiguousGID_ = minMyGID_;
354  lastContiguousGID_ = maxMyGID_;
355  contiguous_ = true;
356  distributed_ = checkIsDist ();
357 
358  // Create the Directory on demand in getRemoteIndexList().
359  //setupDirectory ();
360  }
361 
362  template <class LocalOrdinal, class GlobalOrdinal, class Node>
365  initialNonuniformDebugCheck (const global_size_t numGlobalElements,
366  const size_t numLocalElements,
367  const GlobalOrdinal indexBase,
368  const Teuchos::RCP<const Teuchos::Comm<int> >& comm) const
369  {
370 #ifdef HAVE_TPETRA_DEBUG
371  using Teuchos::broadcast;
372  using Teuchos::outArg;
373  using Teuchos::ptr;
374  using Teuchos::REDUCE_MAX;
375  using Teuchos::REDUCE_MIN;
376  using Teuchos::REDUCE_SUM;
377  using Teuchos::reduceAll;
378  typedef GlobalOrdinal GO;
379  typedef global_size_t GST;
380  const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid ();
381 
382  // The user has specified the distribution of indices over the
383  // processes. The distribution is not necessarily contiguous or
384  // equally shared over the processes.
385  //
386  // We assume that the number of local elements can be stored in a
387  // size_t. The instance member numLocalElements_ is a size_t, so
388  // this variable and that should have the same type.
389 
390  GST debugGlobalSum = 0; // Will be global sum of numLocalElements
391  reduceAll<int, GST> (*comm, REDUCE_SUM, static_cast<GST> (numLocalElements),
392  outArg (debugGlobalSum));
393  // In debug mode only, check whether numGlobalElements and
394  // indexBase are the same over all processes in the communicator.
395  {
396  GST proc0NumGlobalElements = numGlobalElements;
397  broadcast<int, GST> (*comm_, 0, outArg (proc0NumGlobalElements));
398  GST minNumGlobalElements = numGlobalElements;
399  GST maxNumGlobalElements = numGlobalElements;
400  reduceAll<int, GST> (*comm, REDUCE_MIN, numGlobalElements,
401  outArg (minNumGlobalElements));
402  reduceAll<int, GST> (*comm, REDUCE_MAX, numGlobalElements,
403  outArg (maxNumGlobalElements));
404  TEUCHOS_TEST_FOR_EXCEPTION(
405  minNumGlobalElements != maxNumGlobalElements ||
406  numGlobalElements != minNumGlobalElements,
407  std::invalid_argument,
408  "Tpetra::Map constructor: All processes must provide the same number "
409  "of global elements. This is true even if that argument is Teuchos::"
410  "OrdinalTraits<global_size_t>::invalid() to signal that the Map should "
411  "compute the global number of elements. Process 0 set numGlobalElements"
412  " = " << proc0NumGlobalElements << ". The calling process "
413  << comm->getRank () << " set numGlobalElements = " << numGlobalElements
414  << ". The min and max values over all processes are "
415  << minNumGlobalElements << " resp. " << maxNumGlobalElements << ".");
416 
417  GO proc0IndexBase = indexBase;
418  broadcast<int, GO> (*comm_, 0, outArg (proc0IndexBase));
419  GO minIndexBase = indexBase;
420  GO maxIndexBase = indexBase;
421  reduceAll<int, GO> (*comm, REDUCE_MIN, indexBase, outArg (minIndexBase));
422  reduceAll<int, GO> (*comm, REDUCE_MAX, indexBase, outArg (maxIndexBase));
423  TEUCHOS_TEST_FOR_EXCEPTION(
424  minIndexBase != maxIndexBase || indexBase != minIndexBase,
425  std::invalid_argument,
426  "Tpetra::Map constructor: "
427  "All processes must provide the same indexBase argument. "
428  "Process 0 set indexBase = " << proc0IndexBase << ". The calling "
429  "process " << comm->getRank () << " set indexBase = " << indexBase
430  << ". The min and max values over all processes are "
431  << minIndexBase << " resp. " << maxIndexBase << ".");
432 
433  // Make sure that the sum of numLocalElements over all processes
434  // equals numGlobalElements.
435  TEUCHOS_TEST_FOR_EXCEPTION
436  (numGlobalElements != GSTI && debugGlobalSum != numGlobalElements,
437  std::invalid_argument, "Tpetra::Map constructor: The sum of each "
438  "process' number of indices over all processes, " << debugGlobalSum
439  << " != numGlobalElements = " << numGlobalElements << ". If you "
440  "would like this constructor to compute numGlobalElements for you, "
441  "you may set numGlobalElements = "
442  "Teuchos::OrdinalTraits<Tpetra::global_size_t>::invalid() on input. "
443  "Please note that this is NOT necessarily -1.");
444  }
445 
446  return debugGlobalSum;
447 #else
448  (void)numGlobalElements;
449  (void)numLocalElements;
450  (void)indexBase;
451  (void)comm;
452  return static_cast<global_size_t> (0);
453 #endif // HAVE_TPETRA_DEBUG
454  }
455 
456  template <class LocalOrdinal, class GlobalOrdinal, class Node>
457  void
458  Map<LocalOrdinal,GlobalOrdinal,Node>::
459  initWithNonownedHostIndexList (const global_size_t numGlobalElements,
460  const Kokkos::View<const GlobalOrdinal*,
461  Kokkos::LayoutLeft,
462  Kokkos::HostSpace,
463  Kokkos::MemoryUnmanaged>& entryList_host,
464  const GlobalOrdinal indexBase,
465  const Teuchos::RCP<const Teuchos::Comm<int> >& comm)
466  {
467  using Kokkos::LayoutLeft;
468  using Kokkos::subview;
469  using Kokkos::View;
470  using Teuchos::as;
471  using Teuchos::broadcast;
472  using Teuchos::outArg;
473  using Teuchos::ptr;
474  using Teuchos::REDUCE_MAX;
475  using Teuchos::REDUCE_MIN;
476  using Teuchos::REDUCE_SUM;
477  using Teuchos::reduceAll;
478  typedef LocalOrdinal LO;
479  typedef GlobalOrdinal GO;
480  typedef global_size_t GST;
481  const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid ();
482 
483  // Make sure that Kokkos has been initialized (Github Issue #513).
484  TEUCHOS_TEST_FOR_EXCEPTION
485  (! Kokkos::is_initialized (), std::runtime_error,
486  "Tpetra::Map constructor: The Kokkos execution space "
487  << Teuchos::TypeNameTraits<execution_space>::name ()
488  << " has not been initialized. "
489  "Please initialize it before creating a Map.")
490 
491  // The user has specified the distribution of indices over the
492  // processes, via the input array of global indices on each
493  // process. The distribution is not necessarily contiguous or
494  // equally shared over the processes.
495 
496  // The length of the input array on this process is the number of
497  // local indices to associate with this process, even though the
498  // input array contains global indices. We assume that the number
499  // of local indices on a process can be stored in a size_t;
500  // numLocalElements_ is a size_t, so this variable and that should
501  // have the same type.
502  const size_t numLocalElements = static_cast<size_t> (entryList_host.size ());
503 
504  initialNonuniformDebugCheck (numGlobalElements, numLocalElements,
505  indexBase, comm);
506 
507  // NOTE (mfh 20 Feb 2013, 10 Oct 2016) In some sense, this global
508  // reduction is redundant, since the directory Map will have to do
509  // the same thing. Thus, we could do the scan and broadcast for
510  // the directory Map here, and give the computed offsets to the
511  // directory Map's constructor. However, a reduction costs less
512  // than a scan and broadcast, so this still saves time if users of
513  // this Map don't ever need the Directory (i.e., if they never
514  // call getRemoteIndexList on this Map).
515  if (numGlobalElements != GSTI) {
516  numGlobalElements_ = numGlobalElements; // Use the user's value.
517  }
518  else { // The user wants us to compute the sum.
519  reduceAll<int, GST> (*comm, REDUCE_SUM,
520  static_cast<GST> (numLocalElements),
521  outArg (numGlobalElements_));
522  }
523 
524  // mfh 20 Feb 2013: We've never quite done the right thing for
525  // duplicate GIDs here. Duplicate GIDs have always been counted
526  // distinctly in numLocalElements_, and thus should get a
527  // different LID. However, we've always used std::map or a hash
528  // table for the GID -> LID lookup table, so distinct GIDs always
529  // map to the same LID. Furthermore, the order of the input GID
530  // list matters, so it's not desirable to sort for determining
531  // uniqueness.
532  //
533  // I've chosen for now to write this code as if the input GID list
534  // contains no duplicates. If this is not desired, we could use
535  // the lookup table itself to determine uniqueness: If we haven't
536  // seen the GID before, it gets a new LID and it's added to the
537  // LID -> GID and GID -> LID tables. If we have seen the GID
538  // before, it doesn't get added to either table. I would
539  // implement this, but it would cost more to do the double lookups
540  // in the table (one to check, and one to insert).
541  //
542  // More importantly, since we build the GID -> LID table in (a
543  // thread-) parallel (way), the order in which duplicate GIDs may
544  // get inserted is not defined. This would make the assignment of
545  // LID to GID nondeterministic.
546 
547  numLocalElements_ = numLocalElements;
548  indexBase_ = indexBase;
549 
550  minMyGID_ = indexBase_;
551  maxMyGID_ = indexBase_;
552 
553  // NOTE (mfh 27 May 2015): While finding the initial contiguous
554  // GID range requires looking at all the GIDs in the range,
555  // dismissing an interval of GIDs only requires looking at the
556  // first and last GIDs. Thus, we could do binary search backwards
557  // from the end in order to catch the common case of a contiguous
558  // interval followed by noncontiguous entries. On the other hand,
559  // we could just expose this case explicitly as yet another Map
560  // constructor, and avoid the trouble of detecting it.
561  if (numLocalElements_ > 0) {
562  // Find contiguous GID range, with the restriction that the
563  // beginning of the range starts with the first entry. While
564  // doing so, fill in the LID -> GID table.
565  View<GO*, LayoutLeft, device_type> lgMap ("lgMap", numLocalElements_);
566  auto lgMap_host = Kokkos::create_mirror_view (lgMap);
567 
568  // The input array entryList_host is already on host, so we
569  // don't need to take a host view of it.
570  // auto entryList_host = Kokkos::create_mirror_view (entryList);
571  // Kokkos::deep_copy (entryList_host, entryList);
572 
573  firstContiguousGID_ = entryList_host[0];
574  lastContiguousGID_ = firstContiguousGID_+1;
575 
576  // FIXME (mfh 23 Sep 2015) We need to copy the input GIDs
577  // anyway, so we have to look at them all. The logical way to
578  // find the first noncontiguous entry would thus be to "reduce,"
579  // where the local reduction result is whether entryList[i] + 1
580  // == entryList[i+1].
581 
582  lgMap_host[0] = firstContiguousGID_;
583  size_t i = 1;
584  for ( ; i < numLocalElements_; ++i) {
585  const GO curGid = entryList_host[i];
586  const LO curLid = as<LO> (i);
587 
588  if (lastContiguousGID_ != curGid) break;
589 
590  // Add the entry to the LID->GID table only after we know that
591  // the current GID is in the initial contiguous sequence, so
592  // that we don't repeat adding it in the first iteration of
593  // the loop below over the remaining noncontiguous GIDs.
594  lgMap_host[curLid] = curGid;
595  ++lastContiguousGID_;
596  }
597  --lastContiguousGID_;
598 
599  // [firstContiguousGID_, lastContigousGID_] is the initial
600  // sequence of contiguous GIDs. We can start the min and max
601  // GID using this range.
602  minMyGID_ = firstContiguousGID_;
603  maxMyGID_ = lastContiguousGID_;
604 
605  // Compute the GID -> LID lookup table, _not_ including the
606  // initial sequence of contiguous GIDs.
607  {
608  const std::pair<size_t, size_t> ncRange (i, entryList_host.extent (0));
609  auto nonContigGids_host = subview (entryList_host, ncRange);
610  TEUCHOS_TEST_FOR_EXCEPTION
611  (static_cast<size_t> (nonContigGids_host.extent (0)) !=
612  static_cast<size_t> (entryList_host.extent (0) - i),
613  std::logic_error, "Tpetra::Map noncontiguous constructor: "
614  "nonContigGids_host.extent(0) = "
615  << nonContigGids_host.extent (0)
616  << " != entryList_host.extent(0) - i = "
617  << (entryList_host.extent (0) - i) << " = "
618  << entryList_host.extent (0) << " - " << i
619  << ". Please report this bug to the Tpetra developers.");
620 
621  // FixedHashTable's constructor expects an owned device View,
622  // so we must deep-copy the subview of the input indices.
623  View<GO*, LayoutLeft, device_type>
624  nonContigGids ("nonContigGids", nonContigGids_host.size ());
625  Kokkos::deep_copy (nonContigGids, nonContigGids_host);
626 
627  glMap_ = global_to_local_table_type (nonContigGids,
628  firstContiguousGID_,
629  lastContiguousGID_,
630  static_cast<LO> (i));
631  }
632 
633  // FIXME (mfh 10 Oct 2016) When we construct the global-to-local
634  // table above, we have to look at all the (noncontiguous) input
635  // indices anyway. Thus, why not have the constructor compute
636  // and return the min and max?
637 
638  for ( ; i < numLocalElements_; ++i) {
639  const GO curGid = entryList_host[i];
640  const LO curLid = as<LO> (i);
641  lgMap_host[curLid] = curGid; // LID -> GID table
642 
643  // While iterating through entryList, we compute its
644  // (process-local) min and max elements.
645  if (curGid < minMyGID_) {
646  minMyGID_ = curGid;
647  }
648  if (curGid > maxMyGID_) {
649  maxMyGID_ = curGid;
650  }
651  }
652 
653  // We filled lgMap on host above; now sync back to device.
654  Kokkos::deep_copy (lgMap, lgMap_host);
655 
656  // "Commit" the local-to-global lookup table we filled in above.
657  lgMap_ = lgMap;
658  // We've already created this, so use it.
659  lgMapHost_ = lgMap_host;
660  }
661  else {
662  minMyGID_ = std::numeric_limits<GlobalOrdinal>::max();
663  maxMyGID_ = std::numeric_limits<GlobalOrdinal>::lowest();
664  // This insures tests for GIDs in the range
665  // [firstContiguousGID_, lastContiguousGID_] fail for processes
666  // with no local elements.
667  firstContiguousGID_ = indexBase_+1;
668  lastContiguousGID_ = indexBase_;
669  // glMap_ was default constructed, so it's already empty.
670  }
671 
672  // Compute the min and max of all processes' GIDs. If
673  // numLocalElements_ == 0 on this process, minMyGID_ and maxMyGID_
674  // are both indexBase_. This is wrong, but fixing it would
675  // require either a fancy sparse all-reduce, or a custom reduction
676  // operator that ignores invalid values ("invalid" means
677  // Tpetra::Details::OrdinalTraits<GO>::invalid()).
678  //
679  // Also, while we're at it, use the same all-reduce to figure out
680  // if the Map is distributed. "Distributed" means that there is
681  // at least one process with a number of local elements less than
682  // the number of global elements.
683  //
684  // We're computing the min and max of all processes' GIDs using a
685  // single MAX all-reduce, because min(x,y) = -max(-x,-y) (when x
686  // and y are signed). (This lets us combine the min and max into
687  // a single all-reduce.) If each process sets localDist=1 if its
688  // number of local elements is strictly less than the number of
689  // global elements, and localDist=0 otherwise, then a MAX
690  // all-reduce on localDist tells us if the Map is distributed (1
691  // if yes, 0 if no). Thus, we can append localDist onto the end
692  // of the data and get the global result from the all-reduce.
693  if (std::numeric_limits<GO>::is_signed) {
694  // Does my process NOT own all the elements?
695  const GO localDist =
696  (as<GST> (numLocalElements_) < numGlobalElements_) ? 1 : 0;
697 
698  GO minMaxInput[3];
699  minMaxInput[0] = -minMyGID_;
700  minMaxInput[1] = maxMyGID_;
701  minMaxInput[2] = localDist;
702 
703  GO minMaxOutput[3];
704  minMaxOutput[0] = 0;
705  minMaxOutput[1] = 0;
706  minMaxOutput[2] = 0;
707  reduceAll<int, GO> (*comm, REDUCE_MAX, 3, minMaxInput, minMaxOutput);
708  minAllGID_ = -minMaxOutput[0];
709  maxAllGID_ = minMaxOutput[1];
710  const GO globalDist = minMaxOutput[2];
711  distributed_ = (comm_->getSize () > 1 && globalDist == 1);
712  }
713  else { // unsigned; use two reductions
714  // This is always correct, no matter the signedness of GO.
715  reduceAll<int, GO> (*comm_, REDUCE_MIN, minMyGID_, outArg (minAllGID_));
716  reduceAll<int, GO> (*comm_, REDUCE_MAX, maxMyGID_, outArg (maxAllGID_));
717  distributed_ = checkIsDist ();
718  }
719 
720  contiguous_ = false; // "Contiguous" is conservative.
721 
722  TEUCHOS_TEST_FOR_EXCEPTION(
723  minAllGID_ < indexBase_,
724  std::invalid_argument,
725  "Tpetra::Map constructor (noncontiguous): "
726  "Minimum global ID = " << minAllGID_ << " over all process(es) is "
727  "less than the given indexBase = " << indexBase_ << ".");
728 
729  // Create the Directory on demand in getRemoteIndexList().
730  //setupDirectory ();
731  }
732 
733  template <class LocalOrdinal, class GlobalOrdinal, class Node>
735  Map (const global_size_t numGlobalElements,
736  const GlobalOrdinal indexList[],
737  const LocalOrdinal indexListSize,
738  const GlobalOrdinal indexBase,
739  const Teuchos::RCP<const Teuchos::Comm<int> >& comm) :
740  comm_ (comm),
741  uniform_ (false),
742  directory_ (new Directory<LocalOrdinal, GlobalOrdinal, Node> ())
743  {
745 
746  // Not quite sure if I trust all code to behave correctly if the
747  // pointer is nonnull but the array length is nonzero, so I'll
748  // make sure the raw pointer is null if the length is zero.
749  const GlobalOrdinal* const indsRaw = indexListSize == 0 ? NULL : indexList;
750  Kokkos::View<const GlobalOrdinal*,
751  Kokkos::LayoutLeft,
752  Kokkos::HostSpace,
753  Kokkos::MemoryUnmanaged> inds (indsRaw, indexListSize);
754  initWithNonownedHostIndexList (numGlobalElements, inds, indexBase, comm);
755  }
756 
757 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
758  template <class LocalOrdinal, class GlobalOrdinal, class Node>
760  Map (const global_size_t numGlobalElements,
761  const Teuchos::ArrayView<const GlobalOrdinal>& entryList,
762  const GlobalOrdinal indexBase,
763  const Teuchos::RCP<const Teuchos::Comm<int> >& comm,
764  const Teuchos::RCP<Node>& /* node */) :
765  Map<LocalOrdinal,GlobalOrdinal,Node>::Map(numGlobalElements, entryList,
766  indexBase, comm)
767  {}
768 #endif // TPETRA_ENABLE_DEPRECATED_CODE
769 
770  template <class LocalOrdinal, class GlobalOrdinal, class Node>
772  Map (const global_size_t numGlobalElements,
773  const Teuchos::ArrayView<const GlobalOrdinal>& entryList,
774  const GlobalOrdinal indexBase,
775  const Teuchos::RCP<const Teuchos::Comm<int> >& comm) :
776  comm_ (comm),
777  uniform_ (false),
778  directory_ (new Directory<LocalOrdinal, GlobalOrdinal, Node> ())
779  {
781 
782  const size_t numLclInds = static_cast<size_t> (entryList.size ());
783 
784  // Not quite sure if I trust both ArrayView and View to behave
785  // correctly if the pointer is nonnull but the array length is
786  // nonzero, so I'll make sure it's null if the length is zero.
787  const GlobalOrdinal* const indsRaw =
788  numLclInds == 0 ? NULL : entryList.getRawPtr ();
789  Kokkos::View<const GlobalOrdinal*,
790  Kokkos::LayoutLeft,
791  Kokkos::HostSpace,
792  Kokkos::MemoryUnmanaged> inds (indsRaw, numLclInds);
793  initWithNonownedHostIndexList (numGlobalElements, inds, indexBase, comm);
794  }
795 
796  template <class LocalOrdinal, class GlobalOrdinal, class Node>
798  Map (const global_size_t numGlobalElements,
799  const Kokkos::View<const GlobalOrdinal*, device_type>& entryList,
800  const GlobalOrdinal indexBase,
801  const Teuchos::RCP<const Teuchos::Comm<int> >& comm) :
802  comm_ (comm),
803  uniform_ (false),
804  directory_ (new Directory<LocalOrdinal, GlobalOrdinal, Node> ())
805  {
806  using Kokkos::LayoutLeft;
807  using Kokkos::subview;
808  using Kokkos::View;
809  using Teuchos::arcp;
810  using Teuchos::ArrayView;
811  using Teuchos::as;
812  using Teuchos::broadcast;
813  using Teuchos::outArg;
814  using Teuchos::ptr;
815  using Teuchos::REDUCE_MAX;
816  using Teuchos::REDUCE_MIN;
817  using Teuchos::REDUCE_SUM;
818  using Teuchos::reduceAll;
819  using Teuchos::typeName;
820  typedef LocalOrdinal LO;
821  typedef GlobalOrdinal GO;
822  typedef global_size_t GST;
823  const GST GSTI = Tpetra::Details::OrdinalTraits<GST>::invalid ();
824 
826 
827  // The user has specified the distribution of indices over the
828  // processes, via the input array of global indices on each
829  // process. The distribution is not necessarily contiguous or
830  // equally shared over the processes.
831 
832  // The length of the input array on this process is the number of
833  // local indices to associate with this process, even though the
834  // input array contains global indices. We assume that the number
835  // of local indices on a process can be stored in a size_t;
836  // numLocalElements_ is a size_t, so this variable and that should
837  // have the same type.
838  const size_t numLocalElements = static_cast<size_t> (entryList.size ());
839 
840  initialNonuniformDebugCheck (numGlobalElements, numLocalElements,
841  indexBase, comm);
842 
843  // NOTE (mfh 20 Feb 2013, 10 Oct 2016) In some sense, this global
844  // reduction is redundant, since the directory Map will have to do
845  // the same thing. Thus, we could do the scan and broadcast for
846  // the directory Map here, and give the computed offsets to the
847  // directory Map's constructor. However, a reduction costs less
848  // than a scan and broadcast, so this still saves time if users of
849  // this Map don't ever need the Directory (i.e., if they never
850  // call getRemoteIndexList on this Map).
851  if (numGlobalElements != GSTI) {
852  numGlobalElements_ = numGlobalElements; // Use the user's value.
853  }
854  else { // The user wants us to compute the sum.
855  reduceAll<int, GST> (*comm, REDUCE_SUM,
856  static_cast<GST> (numLocalElements),
857  outArg (numGlobalElements_));
858  }
859 
860  // mfh 20 Feb 2013: We've never quite done the right thing for
861  // duplicate GIDs here. Duplicate GIDs have always been counted
862  // distinctly in numLocalElements_, and thus should get a
863  // different LID. However, we've always used std::map or a hash
864  // table for the GID -> LID lookup table, so distinct GIDs always
865  // map to the same LID. Furthermore, the order of the input GID
866  // list matters, so it's not desirable to sort for determining
867  // uniqueness.
868  //
869  // I've chosen for now to write this code as if the input GID list
870  // contains no duplicates. If this is not desired, we could use
871  // the lookup table itself to determine uniqueness: If we haven't
872  // seen the GID before, it gets a new LID and it's added to the
873  // LID -> GID and GID -> LID tables. If we have seen the GID
874  // before, it doesn't get added to either table. I would
875  // implement this, but it would cost more to do the double lookups
876  // in the table (one to check, and one to insert).
877  //
878  // More importantly, since we build the GID -> LID table in (a
879  // thread-) parallel (way), the order in which duplicate GIDs may
880  // get inserted is not defined. This would make the assignment of
881  // LID to GID nondeterministic.
882 
883  numLocalElements_ = numLocalElements;
884  indexBase_ = indexBase;
885 
886  minMyGID_ = indexBase_;
887  maxMyGID_ = indexBase_;
888 
889  // NOTE (mfh 27 May 2015): While finding the initial contiguous
890  // GID range requires looking at all the GIDs in the range,
891  // dismissing an interval of GIDs only requires looking at the
892  // first and last GIDs. Thus, we could do binary search backwards
893  // from the end in order to catch the common case of a contiguous
894  // interval followed by noncontiguous entries. On the other hand,
895  // we could just expose this case explicitly as yet another Map
896  // constructor, and avoid the trouble of detecting it.
897  if (numLocalElements_ > 0) {
898  // Find contiguous GID range, with the restriction that the
899  // beginning of the range starts with the first entry. While
900  // doing so, fill in the LID -> GID table.
901  View<GO*, LayoutLeft, device_type> lgMap ("lgMap", numLocalElements_);
902  auto lgMap_host = Kokkos::create_mirror_view (lgMap);
903 
904  // Creating the mirror view is trivial, and the deep_copy is a
905  // no-op, if entryList is on host already.
906  auto entryList_host = Kokkos::create_mirror_view (entryList);
907  Kokkos::deep_copy (entryList_host, entryList);
908 
909  firstContiguousGID_ = entryList_host[0];
910  lastContiguousGID_ = firstContiguousGID_+1;
911 
912  // FIXME (mfh 23 Sep 2015) We need to copy the input GIDs
913  // anyway, so we have to look at them all. The logical way to
914  // find the first noncontiguous entry would thus be to "reduce,"
915  // where the local reduction result is whether entryList[i] + 1
916  // == entryList[i+1].
917 
918  lgMap_host[0] = firstContiguousGID_;
919  size_t i = 1;
920  for ( ; i < numLocalElements_; ++i) {
921  const GO curGid = entryList_host[i];
922  const LO curLid = as<LO> (i);
923 
924  if (lastContiguousGID_ != curGid) break;
925 
926  // Add the entry to the LID->GID table only after we know that
927  // the current GID is in the initial contiguous sequence, so
928  // that we don't repeat adding it in the first iteration of
929  // the loop below over the remaining noncontiguous GIDs.
930  lgMap_host[curLid] = curGid;
931  ++lastContiguousGID_;
932  }
933  --lastContiguousGID_;
934 
935  // [firstContiguousGID_, lastContigousGID_] is the initial
936  // sequence of contiguous GIDs. We can start the min and max
937  // GID using this range.
938  minMyGID_ = firstContiguousGID_;
939  maxMyGID_ = lastContiguousGID_;
940 
941  // Compute the GID -> LID lookup table, _not_ including the
942  // initial sequence of contiguous GIDs.
943  {
944  const std::pair<size_t, size_t> ncRange (i, entryList.extent (0));
945  auto nonContigGids = subview (entryList, ncRange);
946  TEUCHOS_TEST_FOR_EXCEPTION
947  (static_cast<size_t> (nonContigGids.extent (0)) !=
948  static_cast<size_t> (entryList.extent (0) - i),
949  std::logic_error, "Tpetra::Map noncontiguous constructor: "
950  "nonContigGids.extent(0) = "
951  << nonContigGids.extent (0)
952  << " != entryList.extent(0) - i = "
953  << (entryList.extent (0) - i) << " = "
954  << entryList.extent (0) << " - " << i
955  << ". Please report this bug to the Tpetra developers.");
956 
957  glMap_ = global_to_local_table_type (nonContigGids,
958  firstContiguousGID_,
959  lastContiguousGID_,
960  static_cast<LO> (i));
961  }
962 
963  // FIXME (mfh 10 Oct 2016) When we construct the global-to-local
964  // table above, we have to look at all the (noncontiguous) input
965  // indices anyway. Thus, why not have the constructor compute
966  // and return the min and max?
967 
968  for ( ; i < numLocalElements_; ++i) {
969  const GO curGid = entryList_host[i];
970  const LO curLid = static_cast<LO> (i);
971  lgMap_host[curLid] = curGid; // LID -> GID table
972 
973  // While iterating through entryList, we compute its
974  // (process-local) min and max elements.
975  if (curGid < minMyGID_) {
976  minMyGID_ = curGid;
977  }
978  if (curGid > maxMyGID_) {
979  maxMyGID_ = curGid;
980  }
981  }
982 
983  // We filled lgMap on host above; now sync back to device.
984  Kokkos::deep_copy (lgMap, lgMap_host);
985 
986  // "Commit" the local-to-global lookup table we filled in above.
987  lgMap_ = lgMap;
988  // We've already created this, so use it.
989  lgMapHost_ = lgMap_host;
990  }
991  else {
992  minMyGID_ = std::numeric_limits<GlobalOrdinal>::max();
993  maxMyGID_ = std::numeric_limits<GlobalOrdinal>::lowest();
994  // This insures tests for GIDs in the range
995  // [firstContiguousGID_, lastContiguousGID_] fail for processes
996  // with no local elements.
997  firstContiguousGID_ = indexBase_+1;
998  lastContiguousGID_ = indexBase_;
999  // glMap_ was default constructed, so it's already empty.
1000  }
1001 
1002  // Compute the min and max of all processes' GIDs. If
1003  // numLocalElements_ == 0 on this process, minMyGID_ and maxMyGID_
1004  // are both indexBase_. This is wrong, but fixing it would
1005  // require either a fancy sparse all-reduce, or a custom reduction
1006  // operator that ignores invalid values ("invalid" means
1007  // Tpetra::Details::OrdinalTraits<GO>::invalid()).
1008  //
1009  // Also, while we're at it, use the same all-reduce to figure out
1010  // if the Map is distributed. "Distributed" means that there is
1011  // at least one process with a number of local elements less than
1012  // the number of global elements.
1013  //
1014  // We're computing the min and max of all processes' GIDs using a
1015  // single MAX all-reduce, because min(x,y) = -max(-x,-y) (when x
1016  // and y are signed). (This lets us combine the min and max into
1017  // a single all-reduce.) If each process sets localDist=1 if its
1018  // number of local elements is strictly less than the number of
1019  // global elements, and localDist=0 otherwise, then a MAX
1020  // all-reduce on localDist tells us if the Map is distributed (1
1021  // if yes, 0 if no). Thus, we can append localDist onto the end
1022  // of the data and get the global result from the all-reduce.
1023  if (std::numeric_limits<GO>::is_signed) {
1024  // Does my process NOT own all the elements?
1025  const GO localDist =
1026  (as<GST> (numLocalElements_) < numGlobalElements_) ? 1 : 0;
1027 
1028  GO minMaxInput[3];
1029  minMaxInput[0] = -minMyGID_;
1030  minMaxInput[1] = maxMyGID_;
1031  minMaxInput[2] = localDist;
1032 
1033  GO minMaxOutput[3];
1034  minMaxOutput[0] = 0;
1035  minMaxOutput[1] = 0;
1036  minMaxOutput[2] = 0;
1037  reduceAll<int, GO> (*comm, REDUCE_MAX, 3, minMaxInput, minMaxOutput);
1038  minAllGID_ = -minMaxOutput[0];
1039  maxAllGID_ = minMaxOutput[1];
1040  const GO globalDist = minMaxOutput[2];
1041  distributed_ = (comm_->getSize () > 1 && globalDist == 1);
1042  }
1043  else { // unsigned; use two reductions
1044  // This is always correct, no matter the signedness of GO.
1045  reduceAll<int, GO> (*comm_, REDUCE_MIN, minMyGID_, outArg (minAllGID_));
1046  reduceAll<int, GO> (*comm_, REDUCE_MAX, maxMyGID_, outArg (maxAllGID_));
1047  distributed_ = checkIsDist ();
1048  }
1049 
1050  contiguous_ = false; // "Contiguous" is conservative.
1051 
1052  TEUCHOS_TEST_FOR_EXCEPTION(
1053  minAllGID_ < indexBase_,
1054  std::invalid_argument,
1055  "Tpetra::Map constructor (noncontiguous): "
1056  "Minimum global ID = " << minAllGID_ << " over all process(es) is "
1057  "less than the given indexBase = " << indexBase_ << ".");
1058 
1059  // Create the Directory on demand in getRemoteIndexList().
1060  //setupDirectory ();
1061  }
1062 
1063 
1064  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1066  {
1067  if (! Kokkos::is_initialized ()) {
1068  std::ostringstream os;
1069  os << "WARNING: Tpetra::Map destructor (~Map()) is being called after "
1070  "Kokkos::finalize() has been called. This is user error! There are "
1071  "two likely causes: " << std::endl <<
1072  " 1. You have a static Tpetra::Map (or RCP or shared_ptr of a Map)"
1073  << std::endl <<
1074  " 2. You declare and construct a Tpetra::Map (or RCP or shared_ptr "
1075  "of a Tpetra::Map) at the same scope in main() as Kokkos::finalize() "
1076  "or Tpetra::finalize()." << std::endl << std::endl <<
1077  "Don't do either of these! Please refer to GitHib Issue #2372."
1078  << std::endl;
1079  ::Tpetra::Details::printOnce (std::cerr, os.str (),
1080  this->getComm ().getRawPtr ());
1081  }
1082  else {
1083  using ::Tpetra::Details::mpiIsInitialized;
1084  using ::Tpetra::Details::mpiIsFinalized;
1085  using ::Tpetra::Details::teuchosCommIsAnMpiComm;
1086 
1087  Teuchos::RCP<const Teuchos::Comm<int> > comm = this->getComm ();
1088  if (! comm.is_null () && teuchosCommIsAnMpiComm (*comm) &&
1089  mpiIsInitialized () && mpiIsFinalized ()) {
1090  // Tpetra itself does not require MPI, even if building with
1091  // MPI. It is legal to create Tpetra objects that do not use
1092  // MPI, even in an MPI program. However, calling Tpetra stuff
1093  // after MPI_Finalize() has been called is a bad idea, since
1094  // some Tpetra defaults may use MPI if available.
1095  std::ostringstream os;
1096  os << "WARNING: Tpetra::Map destructor (~Map()) is being called after "
1097  "MPI_Finalize() has been called. This is user error! There are "
1098  "two likely causes: " << std::endl <<
1099  " 1. You have a static Tpetra::Map (or RCP or shared_ptr of a Map)"
1100  << std::endl <<
1101  " 2. You declare and construct a Tpetra::Map (or RCP or shared_ptr "
1102  "of a Tpetra::Map) at the same scope in main() as MPI_finalize() or "
1103  "Tpetra::finalize()." << std::endl << std::endl <<
1104  "Don't do either of these! Please refer to GitHib Issue #2372."
1105  << std::endl;
1106  ::Tpetra::Details::printOnce (std::cerr, os.str (), comm.getRawPtr ());
1107  }
1108  }
1109  // mfh 20 Mar 2018: We can't check Tpetra::isInitialized() yet,
1110  // because Tpetra does not yet require Tpetra::initialize /
1111  // Tpetra::finalize.
1112  }
1113 
1114 
1115  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1116  bool
1118  {
1119  TEUCHOS_TEST_FOR_EXCEPTION(
1120  getComm ().is_null (), std::logic_error, "Tpetra::Map::isOneToOne: "
1121  "getComm() returns null. Please report this bug to the Tpetra "
1122  "developers.");
1123 
1124  // This is a collective operation, if it hasn't been called before.
1125  setupDirectory ();
1126  return directory_->isOneToOne (*this);
1127  }
1128 
1129 
1130  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1131  LocalOrdinal
1133  getLocalElement (GlobalOrdinal globalIndex) const
1134  {
1135  if (isContiguous ()) {
1136  if (globalIndex < getMinGlobalIndex () ||
1137  globalIndex > getMaxGlobalIndex ()) {
1138  return Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid ();
1139  }
1140  return static_cast<LocalOrdinal> (globalIndex - getMinGlobalIndex ());
1141  }
1142  else if (globalIndex >= firstContiguousGID_ &&
1143  globalIndex <= lastContiguousGID_) {
1144  return static_cast<LocalOrdinal> (globalIndex - firstContiguousGID_);
1145  }
1146  else {
1147  // If the given global index is not in the table, this returns
1148  // the same value as OrdinalTraits<LocalOrdinal>::invalid().
1149  return glMap_.get (globalIndex);
1150  }
1151  }
1152 
1153  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1154  GlobalOrdinal
1156  getGlobalElement (LocalOrdinal localIndex) const
1157  {
1158  if (localIndex < getMinLocalIndex () || localIndex > getMaxLocalIndex ()) {
1159  return Tpetra::Details::OrdinalTraits<GlobalOrdinal>::invalid ();
1160  }
1161  if (isContiguous ()) {
1162  return getMinGlobalIndex () + localIndex;
1163  }
1164  else {
1165  // This is a host Kokkos::View access, with no RCP or ArrayRCP
1166  // involvement. As a result, it is thread safe.
1167  //
1168  // lgMapHost_ is a host pointer; this does NOT assume UVM.
1169  return lgMapHost_[localIndex];
1170  }
1171  }
1172 
1173  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1174  bool
1176  isNodeLocalElement (LocalOrdinal localIndex) const
1177  {
1178  if (localIndex < getMinLocalIndex () || localIndex > getMaxLocalIndex ()) {
1179  return false;
1180  } else {
1181  return true;
1182  }
1183  }
1184 
1185  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1186  bool
1188  isNodeGlobalElement (GlobalOrdinal globalIndex) const {
1189  return this->getLocalElement (globalIndex) !=
1190  Tpetra::Details::OrdinalTraits<LocalOrdinal>::invalid ();
1191  }
1192 
1193  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1195  return uniform_;
1196  }
1197 
1198  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1200  return contiguous_;
1201  }
1202 
1203 
1204  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1208  {
1209  return local_map_type (glMap_, lgMap_, getIndexBase (),
1210  getMinGlobalIndex (), getMaxGlobalIndex (),
1211  firstContiguousGID_, lastContiguousGID_,
1212  getNodeNumElements (), isContiguous ());
1213  }
1214 
1215  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1216  bool
1219  {
1220  using Teuchos::outArg;
1221  using Teuchos::REDUCE_MIN;
1222  using Teuchos::reduceAll;
1223  //
1224  // Tests that avoid the Boolean all-reduce below by using
1225  // globally consistent quantities.
1226  //
1227  if (this == &map) {
1228  // Pointer equality on one process always implies pointer
1229  // equality on all processes, since Map is immutable.
1230  return true;
1231  }
1232  else if (getComm ()->getSize () != map.getComm ()->getSize ()) {
1233  // The two communicators have different numbers of processes.
1234  // It's not correct to call isCompatible() in that case. This
1235  // may result in the all-reduce hanging below.
1236  return false;
1237  }
1238  else if (getGlobalNumElements () != map.getGlobalNumElements ()) {
1239  // Two Maps are definitely NOT compatible if they have different
1240  // global numbers of indices.
1241  return false;
1242  }
1243  else if (isContiguous () && isUniform () &&
1244  map.isContiguous () && map.isUniform ()) {
1245  // Contiguous uniform Maps with the same number of processes in
1246  // their communicators, and with the same global numbers of
1247  // indices, are always compatible.
1248  return true;
1249  }
1250  else if (! isContiguous () && ! map.isContiguous () &&
1251  lgMap_.extent (0) != 0 && map.lgMap_.extent (0) != 0 &&
1252  lgMap_.data () == map.lgMap_.data ()) {
1253  // Noncontiguous Maps whose global index lists are nonempty and
1254  // have the same pointer must be the same (and therefore
1255  // contiguous).
1256  //
1257  // Nonempty is important. For example, consider a communicator
1258  // with two processes, and two Maps that share this
1259  // communicator, with zero global indices on the first process,
1260  // and different nonzero numbers of global indices on the second
1261  // process. In that case, on the first process, the pointers
1262  // would both be NULL.
1263  return true;
1264  }
1265 
1266  TEUCHOS_TEST_FOR_EXCEPTION(
1267  getGlobalNumElements () != map.getGlobalNumElements (), std::logic_error,
1268  "Tpetra::Map::isCompatible: There's a bug in this method. We've already "
1269  "checked that this condition is true above, but it's false here. "
1270  "Please report this bug to the Tpetra developers.");
1271 
1272  // Do both Maps have the same number of indices on each process?
1273  const int locallyCompat =
1274  (getNodeNumElements () == map.getNodeNumElements ()) ? 1 : 0;
1275 
1276  int globallyCompat = 0;
1277  reduceAll<int, int> (*comm_, REDUCE_MIN, locallyCompat, outArg (globallyCompat));
1278  return (globallyCompat == 1);
1279  }
1280 
1281  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1282  bool
1285  {
1286  using Teuchos::ArrayView;
1287  typedef GlobalOrdinal GO;
1288  typedef typename ArrayView<const GO>::size_type size_type;
1289 
1290  // If both Maps are contiguous, we can compare their GID ranges
1291  // easily by looking at the min and max GID on this process.
1292  // Otherwise, we'll compare their GID lists. If only one Map is
1293  // contiguous, then we only have to call getNodeElementList() on
1294  // the noncontiguous Map. (It's best to avoid calling it on a
1295  // contiguous Map, since it results in unnecessary storage that
1296  // persists for the lifetime of the Map.)
1297 
1298  if (this == &map) {
1299  // Pointer equality on one process always implies pointer
1300  // equality on all processes, since Map is immutable.
1301  return true;
1302  }
1303  else if (getNodeNumElements () != map.getNodeNumElements ()) {
1304  return false;
1305  }
1306  else if (getMinGlobalIndex () != map.getMinGlobalIndex () ||
1307  getMaxGlobalIndex () != map.getMaxGlobalIndex ()) {
1308  return false;
1309  }
1310  else {
1311  if (isContiguous ()) {
1312  if (map.isContiguous ()) {
1313  return true; // min and max match, so the ranges match.
1314  }
1315  else { // *this is contiguous, but map is not contiguous
1316  TEUCHOS_TEST_FOR_EXCEPTION(
1317  ! this->isContiguous () || map.isContiguous (), std::logic_error,
1318  "Tpetra::Map::locallySameAs: BUG");
1319  ArrayView<const GO> rhsElts = map.getNodeElementList ();
1320  const GO minLhsGid = this->getMinGlobalIndex ();
1321  const size_type numRhsElts = rhsElts.size ();
1322  for (size_type k = 0; k < numRhsElts; ++k) {
1323  const GO curLhsGid = minLhsGid + static_cast<GO> (k);
1324  if (curLhsGid != rhsElts[k]) {
1325  return false; // stop on first mismatch
1326  }
1327  }
1328  return true;
1329  }
1330  }
1331  else if (map.isContiguous ()) { // *this is not contiguous, but map is
1332  TEUCHOS_TEST_FOR_EXCEPTION(
1333  this->isContiguous () || ! map.isContiguous (), std::logic_error,
1334  "Tpetra::Map::locallySameAs: BUG");
1335  ArrayView<const GO> lhsElts = this->getNodeElementList ();
1336  const GO minRhsGid = map.getMinGlobalIndex ();
1337  const size_type numLhsElts = lhsElts.size ();
1338  for (size_type k = 0; k < numLhsElts; ++k) {
1339  const GO curRhsGid = minRhsGid + static_cast<GO> (k);
1340  if (curRhsGid != lhsElts[k]) {
1341  return false; // stop on first mismatch
1342  }
1343  }
1344  return true;
1345  }
1346  else if (this->lgMap_.data () == map.lgMap_.data ()) {
1347  // Pointers to LID->GID "map" (actually just an array) are the
1348  // same, and the number of GIDs are the same.
1349  return this->getNodeNumElements () == map.getNodeNumElements ();
1350  }
1351  else { // we actually have to compare the GIDs
1352  if (this->getNodeNumElements () != map.getNodeNumElements ()) {
1353  return false; // We already checked above, but check just in case
1354  }
1355  else {
1356  ArrayView<const GO> lhsElts = getNodeElementList ();
1357  ArrayView<const GO> rhsElts = map.getNodeElementList ();
1358 
1359  // std::equal requires that the latter range is as large as
1360  // the former. We know the ranges have equal length, because
1361  // they have the same number of local entries.
1362  return std::equal (lhsElts.begin (), lhsElts.end (), rhsElts.begin ());
1363  }
1364  }
1365  }
1366  }
1367 
1368  template <class LocalOrdinal,class GlobalOrdinal, class Node>
1369  bool
1372  {
1373  if (this == &map)
1374  return true;
1375 
1376  // We are going to check if lmap1 is fitted into lmap2
1377  auto lmap1 = map.getLocalMap();
1378  auto lmap2 = this->getLocalMap();
1379 
1380  auto numLocalElements1 = lmap1.getNodeNumElements();
1381  auto numLocalElements2 = lmap2.getNodeNumElements();
1382 
1383  if (numLocalElements1 > numLocalElements2) {
1384  // There are more indices in the first map on this process than in second map.
1385  return false;
1386  }
1387 
1388  if (lmap1.isContiguous () && lmap2.isContiguous ()) {
1389  // When both Maps are contiguous, just check the interval inclusion.
1390  return ((lmap1.getMinGlobalIndex () == lmap2.getMinGlobalIndex ()) &&
1391  (lmap1.getMaxGlobalIndex () <= lmap2.getMaxGlobalIndex ()));
1392  }
1393 
1394  if (lmap1.getMinGlobalIndex () < lmap2.getMinGlobalIndex () ||
1395  lmap1.getMaxGlobalIndex () > lmap2.getMaxGlobalIndex ()) {
1396  // The second map does not include the first map bounds, and thus some of
1397  // the first map global indices are not in the second map.
1398  return false;
1399  }
1400 
1401  typedef Kokkos::RangePolicy<LocalOrdinal, typename Node::execution_space> range_type;
1402 
1403  // Check all elements.
1404  LocalOrdinal numDiff = 0;
1405  Kokkos::parallel_reduce("isLocallyFitted", range_type(0, numLocalElements1),
1406  KOKKOS_LAMBDA(const LocalOrdinal i, LocalOrdinal& diff) {
1407  diff += (lmap1.getGlobalElement(i) != lmap2.getGlobalElement(i));
1408  }, numDiff);
1409 
1410  return (numDiff == 0);
1411  }
1412 
1413  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1414  bool
1417  {
1418  using Teuchos::outArg;
1419  using Teuchos::REDUCE_MIN;
1420  using Teuchos::reduceAll;
1421  //
1422  // Tests that avoid the Boolean all-reduce below by using
1423  // globally consistent quantities.
1424  //
1425  if (this == &map) {
1426  // Pointer equality on one process always implies pointer
1427  // equality on all processes, since Map is immutable.
1428  return true;
1429  }
1430  else if (getComm ()->getSize () != map.getComm ()->getSize ()) {
1431  // The two communicators have different numbers of processes.
1432  // It's not correct to call isSameAs() in that case. This
1433  // may result in the all-reduce hanging below.
1434  return false;
1435  }
1436  else if (getGlobalNumElements () != map.getGlobalNumElements ()) {
1437  // Two Maps are definitely NOT the same if they have different
1438  // global numbers of indices.
1439  return false;
1440  }
1441  else if (getMinAllGlobalIndex () != map.getMinAllGlobalIndex () ||
1442  getMaxAllGlobalIndex () != map.getMaxAllGlobalIndex () ||
1443  getIndexBase () != map.getIndexBase ()) {
1444  // If the global min or max global index doesn't match, or if
1445  // the index base doesn't match, then the Maps aren't the same.
1446  return false;
1447  }
1448  else if (isDistributed () != map.isDistributed ()) {
1449  // One Map is distributed and the other is not, which means that
1450  // the Maps aren't the same.
1451  return false;
1452  }
1453  else if (isContiguous () && isUniform () &&
1454  map.isContiguous () && map.isUniform ()) {
1455  // Contiguous uniform Maps with the same number of processes in
1456  // their communicators, with the same global numbers of indices,
1457  // and with matching index bases and ranges, must be the same.
1458  return true;
1459  }
1460 
1461  // The two communicators must have the same number of processes,
1462  // with process ranks occurring in the same order. This uses
1463  // MPI_COMM_COMPARE. The MPI 3.1 standard (Section 6.4) says:
1464  // "Operations that access communicators are local and their
1465  // execution does not require interprocess communication."
1466  // However, just to be sure, I'll put this call after the above
1467  // tests that don't communicate.
1468  if (! ::Tpetra::Details::congruent (*comm_, * (map.getComm ()))) {
1469  return false;
1470  }
1471 
1472  // If we get this far, we need to check local properties and then
1473  // communicate local sameness across all processes.
1474  const int isSame_lcl = locallySameAs (map) ? 1 : 0;
1475 
1476  // Return true if and only if all processes report local sameness.
1477  int isSame_gbl = 0;
1478  reduceAll<int, int> (*comm_, REDUCE_MIN, isSame_lcl, outArg (isSame_gbl));
1479  return isSame_gbl == 1;
1480  }
1481 
1482  namespace { // (anonymous)
1483  template <class LO, class GO, class DT>
1484  class FillLgMap {
1485  public:
1486  FillLgMap (const Kokkos::View<GO*, DT>& lgMap,
1487  const GO startGid) :
1488  lgMap_ (lgMap), startGid_ (startGid)
1489  {
1490  Kokkos::RangePolicy<LO, typename DT::execution_space>
1491  range (static_cast<LO> (0), static_cast<LO> (lgMap.size ()));
1492  Kokkos::parallel_for (range, *this);
1493  }
1494 
1495  KOKKOS_INLINE_FUNCTION void operator () (const LO& lid) const {
1496  lgMap_(lid) = startGid_ + static_cast<GO> (lid);
1497  }
1498 
1499  private:
1500  const Kokkos::View<GO*, DT> lgMap_;
1501  const GO startGid_;
1502  };
1503 
1504  } // namespace (anonymous)
1505 
1506 
1507  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1508  typename Map<LocalOrdinal,GlobalOrdinal,Node>::global_indices_array_type
1510  {
1511  typedef LocalOrdinal LO;
1512  typedef GlobalOrdinal GO;
1513  typedef device_type DT;
1514 
1515  typedef decltype (lgMap_) const_lg_view_type;
1516  typedef typename const_lg_view_type::non_const_type lg_view_type;
1517 
1518  // If the local-to-global mapping doesn't exist yet, and if we
1519  // have local entries, then create and fill the local-to-global
1520  // mapping.
1521  const bool needToCreateLocalToGlobalMapping =
1522  lgMap_.extent (0) == 0 && numLocalElements_ > 0;
1523 
1524  if (needToCreateLocalToGlobalMapping) {
1525 #ifdef HAVE_TEUCHOS_DEBUG
1526  // The local-to-global mapping should have been set up already
1527  // for a noncontiguous map.
1528  TEUCHOS_TEST_FOR_EXCEPTION( ! isContiguous(), std::logic_error,
1529  "Tpetra::Map::getNodeElementList: The local-to-global mapping (lgMap_) "
1530  "should have been set up already for a noncontiguous Map. Please report"
1531  " this bug to the Tpetra team.");
1532 #endif // HAVE_TEUCHOS_DEBUG
1533 
1534  const LO numElts = static_cast<LO> (getNodeNumElements ());
1535 
1536  lg_view_type lgMap ("lgMap", numElts);
1537  FillLgMap<LO, GO, DT> fillIt (lgMap, minMyGID_);
1538 
1539  auto lgMapHost = Kokkos::create_mirror_view (lgMap);
1540  Kokkos::deep_copy (lgMapHost, lgMap);
1541 
1542  // "Commit" the local-to-global lookup table we filled in above.
1543  lgMap_ = lgMap;
1544  lgMapHost_ = lgMapHost;
1545 
1546  // lgMapHost_ may be a UVM View, so insert a fence to ensure
1547  // coherent host access. We only need to do this once, because
1548  // lgMapHost_ is immutable after initialization.
1549  execution_space::fence ();
1550  }
1551 
1552  return lgMap_;
1553  }
1554 
1555  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1556  Teuchos::ArrayView<const GlobalOrdinal>
1558  {
1559  typedef GlobalOrdinal GO; // convenient abbreviation
1560 
1561  // If the local-to-global mapping doesn't exist yet, and if we
1562  // have local entries, then create and fill the local-to-global
1563  // mapping.
1564  (void) this->getMyGlobalIndices ();
1565 
1566  // This does NOT assume UVM; lgMapHost_ is a host pointer.
1567  const GO* lgMapHostRawPtr = lgMapHost_.data ();
1568  // The third argument forces ArrayView not to try to track memory
1569  // in a debug build. We have to use it because the memory does
1570  // not belong to a Teuchos memory management class.
1571  return Teuchos::ArrayView<const GO> (lgMapHostRawPtr,
1572  lgMapHost_.extent (0),
1573  Teuchos::RCP_DISABLE_NODE_LOOKUP);
1574  }
1575 
1576  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1578  return distributed_;
1579  }
1580 
1581  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1583  using Teuchos::TypeNameTraits;
1584  std::ostringstream os;
1585 
1586  os << "Tpetra::Map: {"
1587  << "LocalOrdinalType: " << TypeNameTraits<LocalOrdinal>::name ()
1588  << ", GlobalOrdinalType: " << TypeNameTraits<GlobalOrdinal>::name ()
1589  << ", NodeType: " << TypeNameTraits<Node>::name ();
1590  if (this->getObjectLabel () != "") {
1591  os << ", Label: \"" << this->getObjectLabel () << "\"";
1592  }
1593  os << ", Global number of entries: " << getGlobalNumElements ()
1594  << ", Number of processes: " << getComm ()->getSize ()
1595  << ", Uniform: " << (isUniform () ? "true" : "false")
1596  << ", Contiguous: " << (isContiguous () ? "true" : "false")
1597  << ", Distributed: " << (isDistributed () ? "true" : "false")
1598  << "}";
1599  return os.str ();
1600  }
1601 
1606  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1607  std::string
1609  localDescribeToString (const Teuchos::EVerbosityLevel vl) const
1610  {
1611  typedef LocalOrdinal LO;
1612  using std::endl;
1613 
1614  // This preserves current behavior of Map.
1615  if (vl < Teuchos::VERB_HIGH) {
1616  return std::string ();
1617  }
1618  auto outStringP = Teuchos::rcp (new std::ostringstream ());
1619  Teuchos::RCP<Teuchos::FancyOStream> outp =
1620  Teuchos::getFancyOStream (outStringP);
1621  Teuchos::FancyOStream& out = *outp;
1622 
1623  auto comm = this->getComm ();
1624  const int myRank = comm->getRank ();
1625  const int numProcs = comm->getSize ();
1626  out << "Process " << myRank << " of " << numProcs << ":" << endl;
1627  Teuchos::OSTab tab1 (out);
1628 
1629  const LO numEnt = static_cast<LO> (this->getNodeNumElements ());
1630  out << "My number of entries: " << numEnt << endl
1631  << "My minimum global index: " << this->getMinGlobalIndex () << endl
1632  << "My maximum global index: " << this->getMaxGlobalIndex () << endl;
1633 
1634  if (vl == Teuchos::VERB_EXTREME) {
1635  out << "My global indices: [";
1636  const LO minLclInd = this->getMinLocalIndex ();
1637  for (LO k = 0; k < numEnt; ++k) {
1638  out << minLclInd + this->getGlobalElement (k);
1639  if (k + 1 < numEnt) {
1640  out << ", ";
1641  }
1642  }
1643  out << "]" << endl;
1644  }
1645 
1646  out.flush (); // make sure the ostringstream got everything
1647  return outStringP->str ();
1648  }
1649 
1650  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1651  void
1653  describe (Teuchos::FancyOStream &out,
1654  const Teuchos::EVerbosityLevel verbLevel) const
1655  {
1656  using Teuchos::TypeNameTraits;
1657  using Teuchos::VERB_DEFAULT;
1658  using Teuchos::VERB_NONE;
1659  using Teuchos::VERB_LOW;
1660  using Teuchos::VERB_HIGH;
1661  using std::endl;
1662  typedef LocalOrdinal LO;
1663  typedef GlobalOrdinal GO;
1664  const Teuchos::EVerbosityLevel vl =
1665  (verbLevel == VERB_DEFAULT) ? VERB_LOW : verbLevel;
1666 
1667  if (vl == VERB_NONE) {
1668  return; // don't print anything
1669  }
1670  // If this Map's Comm is null, then the Map does not participate
1671  // in collective operations with the other processes. In that
1672  // case, it is not even legal to call this method. The reasonable
1673  // thing to do in that case is nothing.
1674  auto comm = this->getComm ();
1675  if (comm.is_null ()) {
1676  return;
1677  }
1678  const int myRank = comm->getRank ();
1679  const int numProcs = comm->getSize ();
1680 
1681  // Only Process 0 should touch the output stream, but this method
1682  // in general may need to do communication. Thus, we may need to
1683  // preserve the current tab level across multiple "if (myRank ==
1684  // 0) { ... }" inner scopes. This is why we sometimes create
1685  // OSTab instances by pointer, instead of by value. We only need
1686  // to create them by pointer if the tab level must persist through
1687  // multiple inner scopes.
1688  Teuchos::RCP<Teuchos::OSTab> tab0, tab1;
1689 
1690  if (myRank == 0) {
1691  // At every verbosity level but VERB_NONE, Process 0 prints.
1692  // By convention, describe() always begins with a tab before
1693  // printing.
1694  tab0 = Teuchos::rcp (new Teuchos::OSTab (out));
1695  out << "\"Tpetra::Map\":" << endl;
1696  tab1 = Teuchos::rcp (new Teuchos::OSTab (out));
1697  {
1698  out << "Template parameters:" << endl;
1699  Teuchos::OSTab tab2 (out);
1700  out << "LocalOrdinal: " << TypeNameTraits<LO>::name () << endl
1701  << "GlobalOrdinal: " << TypeNameTraits<GO>::name () << endl
1702  << "Node: " << TypeNameTraits<Node>::name () << endl;
1703  }
1704  const std::string label = this->getObjectLabel ();
1705  if (label != "") {
1706  out << "Label: \"" << label << "\"" << endl;
1707  }
1708  out << "Global number of entries: " << getGlobalNumElements () << endl
1709  << "Minimum global index: " << getMinAllGlobalIndex () << endl
1710  << "Maximum global index: " << getMaxAllGlobalIndex () << endl
1711  << "Index base: " << getIndexBase () << endl
1712  << "Number of processes: " << numProcs << endl
1713  << "Uniform: " << (isUniform () ? "true" : "false") << endl
1714  << "Contiguous: " << (isContiguous () ? "true" : "false") << endl
1715  << "Distributed: " << (isDistributed () ? "true" : "false") << endl;
1716  }
1717 
1718  // This is collective over the Map's communicator.
1719  if (vl >= VERB_HIGH) { // VERB_HIGH or VERB_EXTREME
1720  const std::string lclStr = this->localDescribeToString (vl);
1721  Tpetra::Details::gathervPrint (out, lclStr, *comm);
1722  }
1723  }
1724 
1725  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1726  Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
1728  replaceCommWithSubset (const Teuchos::RCP<const Teuchos::Comm<int> >& newComm) const
1729  {
1730  using Teuchos::RCP;
1731  using Teuchos::rcp;
1732  typedef global_size_t GST;
1733  typedef LocalOrdinal LO;
1734  typedef GlobalOrdinal GO;
1735  typedef Map<LO, GO, Node> map_type;
1736 
1737  // mfh 26 Mar 2013: The lazy way to do this is simply to recreate
1738  // the Map by calling its ordinary public constructor, using the
1739  // original Map's data. This only involves O(1) all-reduces over
1740  // the new communicator, which in the common case only includes a
1741  // small number of processes.
1742 
1743  // Create the Map to return.
1744  if (newComm.is_null () || newComm->getSize () < 1) {
1745  return Teuchos::null; // my process does not participate in the new Map
1746  }
1747  else if (newComm->getSize () == 1) {
1748  // The case where the new communicator has only one process is
1749  // easy. We don't have to communicate to get all the
1750  // information we need. Use the default comm to create the new
1751  // Map, then fill in all the fields directly.
1752  RCP<map_type> newMap (new map_type ());
1753 
1754  newMap->comm_ = newComm;
1755  // mfh 07 Oct 2016: Preserve original behavior, even though the
1756  // original index base may no longer be the globally min global
1757  // index. See #616 for why this doesn't matter so much anymore.
1758  newMap->indexBase_ = this->indexBase_;
1759  newMap->numGlobalElements_ = this->numLocalElements_;
1760  newMap->numLocalElements_ = this->numLocalElements_;
1761  newMap->minMyGID_ = this->minMyGID_;
1762  newMap->maxMyGID_ = this->maxMyGID_;
1763  newMap->minAllGID_ = this->minMyGID_;
1764  newMap->maxAllGID_ = this->maxMyGID_;
1765  newMap->firstContiguousGID_ = this->firstContiguousGID_;
1766  newMap->lastContiguousGID_ = this->lastContiguousGID_;
1767  // Since the new communicator has only one process, neither
1768  // uniformity nor contiguity have changed.
1769  newMap->uniform_ = this->uniform_;
1770  newMap->contiguous_ = this->contiguous_;
1771  // The new communicator only has one process, so the new Map is
1772  // not distributed.
1773  newMap->distributed_ = false;
1774  newMap->lgMap_ = this->lgMap_;
1775  newMap->lgMapHost_ = this->lgMapHost_;
1776  newMap->glMap_ = this->glMap_;
1777  // It's OK not to initialize the new Map's Directory.
1778  // This is initialized lazily, on first call to getRemoteIndexList.
1779 
1780  return newMap;
1781  }
1782  else { // newComm->getSize() != 1
1783  // Even if the original Map is contiguous, the new Map might not
1784  // be, especially if the excluded processes have ranks != 0 or
1785  // newComm->getSize()-1. The common case for this method is to
1786  // exclude many (possibly even all but one) processes, so it
1787  // likely doesn't pay to do the global communication (over the
1788  // original communicator) to figure out whether we can optimize
1789  // the result Map. Thus, we just set up the result Map as
1790  // noncontiguous.
1791  //
1792  // TODO (mfh 07 Oct 2016) We don't actually need to reconstruct
1793  // the global-to-local table, etc. Optimize this code path to
1794  // avoid unnecessary local work.
1795 
1796  // Make Map (re)compute the global number of elements.
1797  const GST RECOMPUTE = Tpetra::Details::OrdinalTraits<GST>::invalid ();
1798  // TODO (mfh 07 Oct 2016) If we use any Map constructor, we have
1799  // to use the noncontiguous Map constructor, since the new Map
1800  // might not be contiguous. Even if the old Map was contiguous,
1801  // some process in the "middle" might have been excluded. If we
1802  // want to avoid local work, we either have to do the setup by
1803  // hand, or write a new Map constructor.
1804 #if 1
1805  // The disabled code here throws the following exception in
1806  // Map's replaceCommWithSubset test:
1807  //
1808  // Throw test that evaluated to true: static_cast<unsigned long long> (numKeys) > static_cast<unsigned long long> (::Kokkos::Details::ArithTraits<ValueType>::max ())
1809  // 10:
1810  // 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.
1811  // 10: Process 3: origComm->replaceCommWithSubset(subsetComm) threw an exception: /scratch/prj/Trilinos/Trilinos/packages/tpetra/core/src/Tpetra_Details_FixedHashTable_def.hpp:1044:
1812 
1813  auto lgMap = this->getMyGlobalIndices ();
1814  typedef typename std::decay<decltype (lgMap.extent (0)) >::type size_type;
1815  const size_type lclNumInds =
1816  static_cast<size_type> (this->getNodeNumElements ());
1817  using Teuchos::TypeNameTraits;
1818  TEUCHOS_TEST_FOR_EXCEPTION
1819  (lgMap.extent (0) != lclNumInds, std::logic_error,
1820  "Tpetra::Map::replaceCommWithSubset: Result of getMyGlobalIndices() "
1821  "has length " << lgMap.extent (0) << " (of type " <<
1822  TypeNameTraits<size_type>::name () << ") != this->getNodeNumElements()"
1823  " = " << this->getNodeNumElements () << ". The latter, upon being "
1824  "cast to size_type = " << TypeNameTraits<size_type>::name () << ", "
1825  "becomes " << lclNumInds << ". Please report this bug to the Tpetra "
1826  "developers.");
1827 #else
1828  Teuchos::ArrayView<const GO> lgMap = this->getNodeElementList ();
1829 #endif // 1
1830 
1831  const GO indexBase = this->getIndexBase ();
1832  return rcp (new map_type (RECOMPUTE, lgMap, indexBase, newComm));
1833  }
1834  }
1835 
1836  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1837  Teuchos::RCP<const Map<LocalOrdinal, GlobalOrdinal, Node> >
1840  {
1841  using Teuchos::Comm;
1842  using Teuchos::null;
1843  using Teuchos::outArg;
1844  using Teuchos::RCP;
1845  using Teuchos::rcp;
1846  using Teuchos::REDUCE_MIN;
1847  using Teuchos::reduceAll;
1848 
1849  // Create the new communicator. split() returns a valid
1850  // communicator on all processes. On processes where color == 0,
1851  // ignore the result. Passing key == 0 tells MPI to order the
1852  // processes in the new communicator by their rank in the old
1853  // communicator.
1854  const int color = (numLocalElements_ == 0) ? 0 : 1;
1855  // MPI_Comm_split must be called collectively over the original
1856  // communicator. We can't just call it on processes with color
1857  // one, even though we will ignore its result on processes with
1858  // color zero.
1859  RCP<const Comm<int> > newComm = comm_->split (color, 0);
1860  if (color == 0) {
1861  newComm = null;
1862  }
1863 
1864  // Create the Map to return.
1865  if (newComm.is_null ()) {
1866  return null; // my process does not participate in the new Map
1867  } else {
1868  // The default constructor that's useful for clone() above is
1869  // also useful here.
1870  RCP<Map> map = rcp (new Map ());
1871 
1872  map->comm_ = newComm;
1873  map->indexBase_ = indexBase_;
1874  map->numGlobalElements_ = numGlobalElements_;
1875  map->numLocalElements_ = numLocalElements_;
1876  map->minMyGID_ = minMyGID_;
1877  map->maxMyGID_ = maxMyGID_;
1878  map->minAllGID_ = minAllGID_;
1879  map->maxAllGID_ = maxAllGID_;
1880  map->firstContiguousGID_= firstContiguousGID_;
1881  map->lastContiguousGID_ = lastContiguousGID_;
1882 
1883  // Uniformity and contiguity have not changed. The directory
1884  // has changed, but we've taken care of that above.
1885  map->uniform_ = uniform_;
1886  map->contiguous_ = contiguous_;
1887 
1888  // If the original Map was NOT distributed, then the new Map
1889  // cannot be distributed.
1890  //
1891  // If the number of processes in the new communicator is 1, then
1892  // the new Map is not distributed.
1893  //
1894  // Otherwise, we have to check the new Map using an all-reduce
1895  // (over the new communicator). For example, the original Map
1896  // may have had some processes with zero elements, and all other
1897  // processes with the same number of elements as in the whole
1898  // Map. That Map is technically distributed, because of the
1899  // processes with zero elements. Removing those processes would
1900  // make the new Map locally replicated.
1901  if (! distributed_ || newComm->getSize () == 1) {
1902  map->distributed_ = false;
1903  } else {
1904  const int iOwnAllGids = (numLocalElements_ == numGlobalElements_) ? 1 : 0;
1905  int allProcsOwnAllGids = 0;
1906  reduceAll<int, int> (*newComm, REDUCE_MIN, iOwnAllGids, outArg (allProcsOwnAllGids));
1907  map->distributed_ = (allProcsOwnAllGids == 1) ? false : true;
1908  }
1909 
1910  map->lgMap_ = lgMap_;
1911  map->lgMapHost_ = lgMapHost_;
1912  map->glMap_ = glMap_;
1913 
1914  // Map's default constructor creates an uninitialized Directory.
1915  // The Directory will be initialized on demand in
1916  // getRemoteIndexList().
1917  //
1918  // FIXME (mfh 26 Mar 2013) It should be possible to "filter" the
1919  // directory more efficiently than just recreating it. If
1920  // directory recreation proves a bottleneck, we can always
1921  // revisit this. On the other hand, Directory creation is only
1922  // collective over the new, presumably much smaller
1923  // communicator, so it may not be worth the effort to optimize.
1924 
1925  return map;
1926  }
1927  }
1928 
1929  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1930  void
1932  {
1933  TEUCHOS_TEST_FOR_EXCEPTION(
1934  directory_.is_null (), std::logic_error, "Tpetra::Map::setupDirectory: "
1935  "The Directory is null. "
1936  "Please report this bug to the Tpetra developers.");
1937 
1938  // Only create the Directory if it hasn't been created yet.
1939  // This is a collective operation.
1940  if (! directory_->initialized ()) {
1941  directory_->initialize (*this);
1942  }
1943  }
1944 
1945  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1946  LookupStatus
1948  getRemoteIndexList (const Teuchos::ArrayView<const GlobalOrdinal>& GIDs,
1949  const Teuchos::ArrayView<int>& PIDs,
1950  const Teuchos::ArrayView<LocalOrdinal>& LIDs) const
1951  {
1952  using Tpetra::Details::OrdinalTraits;
1953  typedef Teuchos::ArrayView<int>::size_type size_type;
1954 
1955  // Empty Maps (i.e., containing no indices on any processes in the
1956  // Map's communicator) are perfectly valid. In that case, if the
1957  // input GID list is nonempty, we fill the output arrays with
1958  // invalid values, and return IDNotPresent to notify the caller.
1959  // It's perfectly valid to give getRemoteIndexList GIDs that the
1960  // Map doesn't own. SubmapImport test 2 needs this functionality.
1961  if (getGlobalNumElements () == 0) {
1962  if (GIDs.size () == 0) {
1963  return AllIDsPresent; // trivially
1964  } else {
1965  for (size_type k = 0; k < PIDs.size (); ++k) {
1966  PIDs[k] = OrdinalTraits<int>::invalid ();
1967  }
1968  for (size_type k = 0; k < LIDs.size (); ++k) {
1969  LIDs[k] = OrdinalTraits<LocalOrdinal>::invalid ();
1970  }
1971  return IDNotPresent;
1972  }
1973  }
1974 
1975  // getRemoteIndexList must be called collectively, and Directory
1976  // initialization is collective too, so it's OK to initialize the
1977  // Directory on demand.
1978  setupDirectory ();
1979  return directory_->getDirectoryEntries (*this, GIDs, PIDs, LIDs);
1980  }
1981 
1982  template <class LocalOrdinal, class GlobalOrdinal, class Node>
1983  LookupStatus
1985  getRemoteIndexList (const Teuchos::ArrayView<const GlobalOrdinal> & GIDs,
1986  const Teuchos::ArrayView<int> & PIDs) const
1987  {
1988  if (getGlobalNumElements () == 0) {
1989  if (GIDs.size () == 0) {
1990  return AllIDsPresent; // trivially
1991  } else {
1992  // The Map contains no indices, so all output PIDs are invalid.
1993  for (Teuchos::ArrayView<int>::size_type k = 0; k < PIDs.size (); ++k) {
1994  PIDs[k] = Tpetra::Details::OrdinalTraits<int>::invalid ();
1995  }
1996  return IDNotPresent;
1997  }
1998  }
1999 
2000  // getRemoteIndexList must be called collectively, and Directory
2001  // initialization is collective too, so it's OK to initialize the
2002  // Directory on demand.
2003  setupDirectory ();
2004  return directory_->getDirectoryEntries (*this, GIDs, PIDs);
2005  }
2006 
2007  template <class LocalOrdinal, class GlobalOrdinal, class Node>
2008  Teuchos::RCP<const Teuchos::Comm<int> >
2010  return comm_;
2011  }
2012 
2013 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2014  template <class LocalOrdinal, class GlobalOrdinal, class Node>
2015  TPETRA_DEPRECATED
2016  Teuchos::RCP<Node>
2018  // Node instances don't do anything any more, but sometimes it
2019  // helps for them to be nonnull.
2020  return Teuchos::rcp (new Node);
2021  }
2022 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2023 
2024  template <class LocalOrdinal,class GlobalOrdinal, class Node>
2025  bool Map<LocalOrdinal,GlobalOrdinal,Node>::checkIsDist() const {
2026  using Teuchos::as;
2027  using Teuchos::outArg;
2028  using Teuchos::REDUCE_MIN;
2029  using Teuchos::reduceAll;
2030 
2031  bool global = false;
2032  if (comm_->getSize () > 1) {
2033  // The communicator has more than one process, but that doesn't
2034  // necessarily mean the Map is distributed.
2035  int localRep = 0;
2036  if (numGlobalElements_ == as<global_size_t> (numLocalElements_)) {
2037  // The number of local elements on this process equals the
2038  // number of global elements.
2039  //
2040  // NOTE (mfh 22 Nov 2011) Does this still work if there were
2041  // duplicates in the global ID list on input (the third Map
2042  // constructor), so that the number of local elements (which
2043  // are not duplicated) on this process could be less than the
2044  // number of global elements, even if this process owns all
2045  // the elements?
2046  localRep = 1;
2047  }
2048  int allLocalRep;
2049  reduceAll<int, int> (*comm_, REDUCE_MIN, localRep, outArg (allLocalRep));
2050  if (allLocalRep != 1) {
2051  // At least one process does not own all the elements.
2052  // This makes the Map a distributed Map.
2053  global = true;
2054  }
2055  }
2056  // If the communicator has only one process, then the Map is not
2057  // distributed.
2058  return global;
2059  }
2060 
2061 } // namespace Tpetra
2062 
2063 template <class LocalOrdinal, class GlobalOrdinal>
2064 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal> >
2065 Tpetra::createLocalMap (const size_t numElements,
2066  const Teuchos::RCP<const Teuchos::Comm<int> >& comm)
2067 {
2068  typedef LocalOrdinal LO;
2069  typedef GlobalOrdinal GO;
2070  typedef typename ::Tpetra::Map<LO, GO>::node_type NT;
2071  return createLocalMapWithNode<LO, GO, NT> (numElements, comm);
2072 }
2073 
2074 template <class LocalOrdinal, class GlobalOrdinal>
2075 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal> >
2076 Tpetra::createUniformContigMap (const global_size_t numElements,
2077  const Teuchos::RCP<const Teuchos::Comm<int> >& comm)
2078 {
2079  typedef LocalOrdinal LO;
2080  typedef GlobalOrdinal GO;
2081  typedef typename ::Tpetra::Map<LO, GO>::node_type NT;
2082  return createUniformContigMapWithNode<LO, GO, NT> (numElements, comm);
2083 }
2084 
2085 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2086 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2087 TPETRA_DEPRECATED
2088 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2089 Tpetra::createUniformContigMapWithNode (const global_size_t numElements,
2090  const Teuchos::RCP<const Teuchos::Comm<int> >& comm,
2091  const Teuchos::RCP<Node>& /* node */
2092 )
2093 {
2094  return Tpetra::createUniformContigMapWithNode<LocalOrdinal,GlobalOrdinal,Node>
2095  (numElements, comm);
2096 }
2097 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2098 
2099 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2100 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2101 Tpetra::createUniformContigMapWithNode (const global_size_t numElements,
2102  const Teuchos::RCP<const Teuchos::Comm<int> >& comm
2103 )
2104 {
2105  using Teuchos::rcp;
2107  const GlobalOrdinal indexBase = static_cast<GlobalOrdinal> (0);
2108 
2109  return rcp (new map_type (numElements, indexBase, comm, GloballyDistributed));
2110 }
2111 
2112 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2113 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2114 TPETRA_DEPRECATED
2115 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2116 Tpetra::createLocalMapWithNode (const size_t numElements,
2117  const Teuchos::RCP<const Teuchos::Comm<int> >& comm,
2118  const Teuchos::RCP<Node>& /* node */
2119 )
2120 {
2121  return Tpetra::createLocalMapWithNode<LocalOrdinal,GlobalOrdinal,Node>
2122  (numElements, comm);
2123 }
2124 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2125 
2126 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2127 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2128 Tpetra::createLocalMapWithNode (const size_t numElements,
2129  const Teuchos::RCP<const Teuchos::Comm<int> >& comm
2130 )
2131 {
2132  using Tpetra::global_size_t;
2133  using Teuchos::rcp;
2135  const GlobalOrdinal indexBase = static_cast<GlobalOrdinal> (0);
2136  const global_size_t globalNumElts = static_cast<global_size_t> (numElements);
2137 
2138  return rcp (new map_type (globalNumElts, indexBase, comm, LocallyReplicated));
2139 }
2140 
2141 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2142 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2143 TPETRA_DEPRECATED
2144 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2146  const size_t localNumElements,
2147  const Teuchos::RCP<const Teuchos::Comm<int> >& comm,
2148  const Teuchos::RCP<Node>& /* node */
2149 )
2150 {
2151  return Tpetra::createContigMapWithNode<LocalOrdinal,GlobalOrdinal,Node>
2152  (numElements, localNumElements, comm);
2153 }
2154 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2155 
2156 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2157 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2159  const size_t localNumElements,
2160  const Teuchos::RCP<const Teuchos::Comm<int> >& comm
2161 )
2162 {
2163  using Teuchos::rcp;
2165  const GlobalOrdinal indexBase = static_cast<GlobalOrdinal> (0);
2166 
2167  return rcp (new map_type (numElements, localNumElements, indexBase, comm));
2168 }
2169 
2170 template <class LocalOrdinal, class GlobalOrdinal>
2171 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal> >
2173  const size_t localNumElements,
2174  const Teuchos::RCP<const Teuchos::Comm<int> >& comm)
2175 {
2176  typedef LocalOrdinal LO;
2177  typedef GlobalOrdinal GO;
2179 
2180  return Tpetra::createContigMapWithNode<LO, GO, NT> (numElements, localNumElements, comm);
2181 }
2182 
2183 
2184 template <class LocalOrdinal, class GlobalOrdinal>
2185 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal> >
2186 Tpetra::createNonContigMap(const Teuchos::ArrayView<const GlobalOrdinal>& elementList,
2187  const Teuchos::RCP<const Teuchos::Comm<int> >& comm)
2188 {
2189  typedef LocalOrdinal LO;
2190  typedef GlobalOrdinal GO;
2192 
2193  return Tpetra::createNonContigMapWithNode<LO, GO, NT> (elementList, comm);
2194 }
2195 
2196 
2197 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2198 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2199 TPETRA_DEPRECATED
2200 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2201 Tpetra::createNonContigMapWithNode (const Teuchos::ArrayView<const GlobalOrdinal>& elementList,
2202  const Teuchos::RCP<const Teuchos::Comm<int> >& comm,
2203  const Teuchos::RCP<Node>& /* node */
2204 )
2205 {
2206  return Tpetra::createNonContigMapWithNode<LocalOrdinal,GlobalOrdinal,Node>
2207  (elementList, comm);
2208 }
2209 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2210 
2211 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2212 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2213 Tpetra::createNonContigMapWithNode (const Teuchos::ArrayView<const GlobalOrdinal>& elementList,
2214  const Teuchos::RCP<const Teuchos::Comm<int> >& comm
2215 )
2216 {
2217  using Teuchos::rcp;
2219  using GST = Tpetra::global_size_t;
2220  const GST INV = Tpetra::Details::OrdinalTraits<GST>::invalid ();
2221  // FIXME (mfh 22 Jul 2016) This is what I found here, but maybe this
2222  // shouldn't be zero, given that the index base is supposed to equal
2223  // the globally min global index?
2224  const GlobalOrdinal indexBase = 0;
2225 
2226  return rcp (new map_type (INV, elementList, indexBase, comm));
2227 }
2228 
2229 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2230 template <class LocalOrdinal, class GlobalOrdinal, class Node>
2231 TPETRA_DEPRECATED
2232 Teuchos::RCP< const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2233 Tpetra::createWeightedContigMapWithNode (const int myWeight,
2234  const Tpetra::global_size_t numElements,
2235  const Teuchos::RCP<const Teuchos::Comm<int> >& comm,
2236  const Teuchos::RCP<Node>& /* node */
2237 )
2238 {
2239  Teuchos::RCP< Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> > map;
2240  int sumOfWeights, elemsLeft, localNumElements;
2241  const int numImages = comm->getSize();
2242  const int myImageID = comm->getRank();
2243  Teuchos::reduceAll<int>(*comm,Teuchos::REDUCE_SUM,myWeight,Teuchos::outArg(sumOfWeights));
2244  const double myShare = ((double)myWeight) / ((double)sumOfWeights);
2245  localNumElements = (int)std::floor( myShare * ((double)numElements) );
2246  // std::cout << "numElements: " << numElements << " myWeight: " << myWeight << " sumOfWeights: " << sumOfWeights << " myShare: " << myShare << std::endl;
2247  Teuchos::reduceAll<int>(*comm,Teuchos::REDUCE_SUM,localNumElements,Teuchos::outArg(elemsLeft));
2248  elemsLeft = numElements - elemsLeft;
2249  // std::cout << "(before) localNumElements: " << localNumElements << " elemsLeft: " << elemsLeft << std::endl;
2250  // i think this is true. just test it for now.
2251  TEUCHOS_TEST_FOR_EXCEPT(elemsLeft < -numImages || numImages < elemsLeft);
2252  if (elemsLeft < 0) {
2253  // last elemsLeft nodes lose an element
2254  if (myImageID >= numImages-elemsLeft) --localNumElements;
2255  }
2256  else if (elemsLeft > 0) {
2257  // first elemsLeft nodes gain an element
2258  if (myImageID < elemsLeft) ++localNumElements;
2259  }
2260  // std::cout << "(after) localNumElements: " << localNumElements << std::endl;
2261  return createContigMapWithNode<LocalOrdinal,GlobalOrdinal,Node>(numElements,localNumElements,comm);
2262 }
2263 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2264 
2265 
2266 template<class LO, class GO, class NT>
2267 Teuchos::RCP<const Tpetra::Map<LO, GO, NT> >
2268 Tpetra::createOneToOne (const Teuchos::RCP<const Tpetra::Map<LO, GO, NT> >& M)
2269 {
2270  using Teuchos::Array;
2271  using Teuchos::ArrayView;
2272  using Teuchos::as;
2273  using Teuchos::rcp;
2274  typedef Tpetra::Map<LO, GO, NT> map_type;
2275  typedef global_size_t GST;
2276  const GST GINV = Tpetra::Details::OrdinalTraits<GST>::invalid ();
2277  const int myRank = M->getComm ()->getRank ();
2278 
2279  // Bypasses for special cases where either M is known to be
2280  // one-to-one, or the one-to-one version of M is easy to compute.
2281  // This is why we take M as an RCP, not as a const reference -- so
2282  // that we can return M itself if it is 1-to-1.
2283  if (! M->isDistributed ()) {
2284  // For a locally replicated Map, we assume that users want to push
2285  // all the GIDs to Process 0.
2286 
2287  // mfh 05 Nov 2013: getGlobalNumElements() does indeed return what
2288  // you think it should return, in this special case of a locally
2289  // replicated contiguous Map.
2290  const GST numGlobalEntries = M->getGlobalNumElements ();
2291  if (M->isContiguous ()) {
2292  const size_t numLocalEntries =
2293  (myRank == 0) ? as<size_t> (numGlobalEntries) : static_cast<size_t> (0);
2294  return rcp (new map_type (numGlobalEntries, numLocalEntries,
2295  M->getIndexBase (), M->getComm ()));
2296  }
2297  else {
2298  ArrayView<const GO> myGids =
2299  (myRank == 0) ? M->getNodeElementList () : Teuchos::null;
2300  return rcp (new map_type (GINV, myGids (), M->getIndexBase (),
2301  M->getComm ()));
2302  }
2303  }
2304  else if (M->isContiguous ()) {
2305  // Contiguous, distributed Maps are one-to-one by construction.
2306  // (Locally replicated Maps can be contiguous.)
2307  return M;
2308  }
2309  else {
2311  const size_t numMyElems = M->getNodeNumElements ();
2312  ArrayView<const GO> myElems = M->getNodeElementList ();
2313  Array<int> owner_procs_vec (numMyElems);
2314 
2315  directory.getDirectoryEntries (*M, myElems, owner_procs_vec ());
2316 
2317  Array<GO> myOwned_vec (numMyElems);
2318  size_t numMyOwnedElems = 0;
2319  for (size_t i = 0; i < numMyElems; ++i) {
2320  const GO GID = myElems[i];
2321  const int owner = owner_procs_vec[i];
2322 
2323  if (myRank == owner) {
2324  myOwned_vec[numMyOwnedElems++] = GID;
2325  }
2326  }
2327  myOwned_vec.resize (numMyOwnedElems);
2328 
2329  return rcp (new map_type (GINV, myOwned_vec (), M->getIndexBase (),
2330  M->getComm ()));
2331  }
2332 }
2333 
2334 template<class LocalOrdinal, class GlobalOrdinal, class Node>
2335 Teuchos::RCP<const Tpetra::Map<LocalOrdinal,GlobalOrdinal,Node> >
2338 {
2339  using Teuchos::Array;
2340  using Teuchos::ArrayView;
2341  using Teuchos::rcp;
2342  typedef LocalOrdinal LO;
2343  typedef GlobalOrdinal GO;
2344  typedef Tpetra::Map<LO,GO,Node> map_type;
2345  int myID = M->getComm()->getRank();
2346 
2347  // FIXME (mfh 20 Feb 2013) We should have a bypass for contiguous
2348  // Maps (which are 1-to-1 by construction).
2349 
2350  //Based off Epetra's one to one.
2351 
2353  directory.initialize (*M, tie_break);
2354  size_t numMyElems = M->getNodeNumElements ();
2355  ArrayView<const GO> myElems = M->getNodeElementList ();
2356  Array<int> owner_procs_vec (numMyElems);
2357 
2358  directory.getDirectoryEntries (*M, myElems, owner_procs_vec ());
2359 
2360  Array<GO> myOwned_vec (numMyElems);
2361  size_t numMyOwnedElems = 0;
2362  for (size_t i = 0; i < numMyElems; ++i) {
2363  GO GID = myElems[i];
2364  int owner = owner_procs_vec[i];
2365 
2366  if (myID == owner) {
2367  myOwned_vec[numMyOwnedElems++] = GID;
2368  }
2369  }
2370  myOwned_vec.resize (numMyOwnedElems);
2371 
2372  // FIXME (mfh 08 May 2014) The above Directory should be perfectly
2373  // valid for the new Map. Why can't we reuse it?
2374  const global_size_t GINV =
2375  Tpetra::Details::OrdinalTraits<global_size_t>::invalid ();
2376  return rcp (new map_type (GINV, myOwned_vec (), M->getIndexBase (),
2377  M->getComm ()));
2378 }
2379 
2380 //
2381 // Explicit instantiation macro
2382 //
2383 // Must be expanded from within the Tpetra namespace!
2384 //
2385 
2387 #ifdef TPETRA_ENABLE_DEPRECATED_CODE
2388 
2389 #define TPETRA_MAP_INSTANT(LO,GO,NODE) \
2390  \
2391  template class Map< LO , GO , NODE >; \
2392  \
2393  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2394  createLocalMapWithNode<LO,GO,NODE> (const size_t numElements, \
2395  const Teuchos::RCP< const Teuchos::Comm< int > >& comm); \
2396  \
2397  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2398  createContigMapWithNode<LO,GO,NODE> (const global_size_t numElements, \
2399  const size_t localNumElements, \
2400  const Teuchos::RCP< const Teuchos::Comm< int > >& comm); \
2401  \
2402  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2403  createNonContigMapWithNode(const Teuchos::ArrayView<const GO> &elementList, \
2404  const Teuchos::RCP<const Teuchos::Comm<int> > &comm); \
2405  \
2406  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2407  createUniformContigMapWithNode<LO,GO,NODE> (const global_size_t numElements, \
2408  const Teuchos::RCP< const Teuchos::Comm< int > >& comm); \
2409  \
2410  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2411  createLocalMapWithNode<LO,GO,NODE> (const size_t numElements, \
2412  const Teuchos::RCP< const Teuchos::Comm< int > >& comm, \
2413  const Teuchos::RCP< NODE >& node); \
2414  \
2415  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2416  createContigMapWithNode<LO,GO,NODE> (const global_size_t numElements, \
2417  const size_t localNumElements, \
2418  const Teuchos::RCP< const Teuchos::Comm< int > >& comm, \
2419  const Teuchos::RCP< NODE > &node); \
2420  \
2421  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2422  createNonContigMapWithNode(const Teuchos::ArrayView<const GO> &elementList, \
2423  const Teuchos::RCP<const Teuchos::Comm<int> > &comm, \
2424  const Teuchos::RCP<NODE> &node); \
2425  \
2426  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2427  createUniformContigMapWithNode<LO,GO,NODE> (const global_size_t numElements, \
2428  const Teuchos::RCP< const Teuchos::Comm< int > >& comm, \
2429  const Teuchos::RCP< NODE > &node); \
2430  \
2431  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2432  createWeightedContigMapWithNode<LO,GO,NODE> (const int thisNodeWeight, \
2433  const global_size_t numElements, \
2434  const Teuchos::RCP< const Teuchos::Comm< int > >& comm, \
2435  const Teuchos::RCP< NODE >& node); \
2436  \
2437  template Teuchos::RCP<const Map<LO,GO,NODE> > \
2438  createOneToOne (const Teuchos::RCP<const Map<LO,GO,NODE> >& M); \
2439  \
2440  template Teuchos::RCP<const Map<LO,GO,NODE> > \
2441  createOneToOne (const Teuchos::RCP<const Map<LO,GO,NODE> >& M, \
2442  const Tpetra::Details::TieBreak<LO,GO>& tie_break); \
2443 
2444 #else // !TPETRA_ENABLE_DEPRECATED_CODE
2445 
2446 #define TPETRA_MAP_INSTANT(LO,GO,NODE) \
2447  \
2448  template class Map< LO , GO , NODE >; \
2449  \
2450  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2451  createLocalMapWithNode<LO,GO,NODE> (const size_t numElements, \
2452  const Teuchos::RCP< const Teuchos::Comm< int > >& comm); \
2453  \
2454  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2455  createContigMapWithNode<LO,GO,NODE> (const global_size_t numElements, \
2456  const size_t localNumElements, \
2457  const Teuchos::RCP< const Teuchos::Comm< int > >& comm); \
2458  \
2459  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2460  createNonContigMapWithNode(const Teuchos::ArrayView<const GO> &elementList, \
2461  const Teuchos::RCP<const Teuchos::Comm<int> > &comm); \
2462  \
2463  template Teuchos::RCP< const Map<LO,GO,NODE> > \
2464  createUniformContigMapWithNode<LO,GO,NODE> (const global_size_t numElements, \
2465  const Teuchos::RCP< const Teuchos::Comm< int > >& comm); \
2466  \
2467  template Teuchos::RCP<const Map<LO,GO,NODE> > \
2468  createOneToOne (const Teuchos::RCP<const Map<LO,GO,NODE> >& M); \
2469  \
2470  template Teuchos::RCP<const Map<LO,GO,NODE> > \
2471  createOneToOne (const Teuchos::RCP<const Map<LO,GO,NODE> >& M, \
2472  const Tpetra::Details::TieBreak<LO,GO>& tie_break); \
2473 
2474 #endif // TPETRA_ENABLE_DEPRECATED_CODE
2475 
2477 #define TPETRA_MAP_INSTANT_DEFAULTNODE(LO,GO) \
2478  template Teuchos::RCP< const Map<LO,GO> > \
2479  createLocalMap<LO,GO>( const size_t, const Teuchos::RCP< const Teuchos::Comm< int > > &); \
2480  \
2481  template Teuchos::RCP< const Map<LO,GO> > \
2482  createContigMap<LO,GO>( global_size_t, size_t, \
2483  const Teuchos::RCP< const Teuchos::Comm< int > > &); \
2484  \
2485  template Teuchos::RCP< const Map<LO,GO> > \
2486  createNonContigMap(const Teuchos::ArrayView<const GO> &, \
2487  const Teuchos::RCP<const Teuchos::Comm<int> > &); \
2488  \
2489  template Teuchos::RCP< const Map<LO,GO> > \
2490  createUniformContigMap<LO,GO>( const global_size_t, \
2491  const Teuchos::RCP< const Teuchos::Comm< int > > &); \
2492 
2493 #endif // TPETRA_MAP_DEF_HPP
Interface for breaking ties in ownership.
bool congruent(const Teuchos::Comm< int > &comm1, const Teuchos::Comm< int > &comm2)
Whether the two communicators are congruent.
Definition: Tpetra_Util.cpp:65
bool isLocallyFitted(const Map< LocalOrdinal, GlobalOrdinal, Node > &map) const
True if and only if map is locally fitted to this Map.
Declaration of Tpetra::Details::extractMpiCommFromTeuchos.
virtual ~Map()
Destructor (virtual for memory safety of derived classes).
GlobalOrdinal getMinGlobalIndex() const
The minimum global index owned by the calling process.
GlobalOrdinal global_ordinal_type
The type of global indices.
Declaration of Tpetra::Details::printOnce.
bool isNodeGlobalElement(GlobalOrdinal globalIndex) const
Whether the given global index is owned by this Map on the calling process.
global_indices_array_type getMyGlobalIndices() const
Return a view of the global indices owned by this process.
bool isContiguous() const
True if this Map is distributed contiguously, else false.
LookupStatus
Return status of Map remote index lookup (getRemoteIndexList()).
bool isCompatible(const Map< LocalOrdinal, GlobalOrdinal, Node > &map) const
True if and only if map is compatible with this Map.
bool isUniform() const
Whether the range of global indices is uniform.
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...
GlobalOrdinal getIndexBase() const
The index base for this Map.
bool isDistributed() const
Whether this Map is globally distributed or locally replicated.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > removeEmptyProcesses() const
Advanced methods.
Teuchos::RCP< const Teuchos::Comm< int > > getComm() const
Accessors for the Teuchos::Comm and Kokkos Node objects.
GlobalOrdinal getMinAllGlobalIndex() const
The minimum global index over all processes in the communicator.
bool locallySameAs(const Map< LocalOrdinal, GlobalOrdinal, node_type > &map) const
Is this Map locally the same as the input Map?
bool isNodeLocalElement(LocalOrdinal localIndex) const
Whether the given local index is valid for this Map on the calling process.
&quot;Local&quot; part of Map suitable for Kokkos kernels.
Declaration of Tpetra::Details::initializeKokkos.
GlobalOrdinal getMaxAllGlobalIndex() const
The maximum global index over all processes in the communicator.
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, in rank order.
void initializeKokkos()
Initialize Kokkos, using command-line arguments (if any) given to Teuchos::GlobalMPISession.
size_t global_size_t
Global size_t object.
void deep_copy(MultiVector< DS, DL, DG, DN > &dst, const MultiVector< SS, SL, SG, SN > &src)
Copy the contents of the MultiVector src into dst.
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, possibly nondefault Kokkos Node type.
Teuchos::RCP< const Map< LocalOrdinal, GlobalOrdinal, Node > > replaceCommWithSubset(const Teuchos::RCP< const Teuchos::Comm< int > > &newComm) const
Replace this Map&#39;s communicator with a subset communicator.
Functions for initializing and finalizing Tpetra.
GlobalOrdinal getGlobalElement(LocalOrdinal localIndex) const
The global index corresponding to the given local index.
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...
GlobalOrdinal getMaxGlobalIndex() const
The maximum global index owned by the calling process.
LookupStatus getRemoteIndexList(const Teuchos::ArrayView< const GlobalOrdinal > &GIDList, const Teuchos::ArrayView< int > &nodeIDList, const Teuchos::ArrayView< LocalOrdinal > &LIDList) const
Return the process ranks and corresponding local indices for the given global indices.
LookupStatus getDirectoryEntries(const map_type &map, const Teuchos::ArrayView< const GlobalOrdinal > &globalIDs, const Teuchos::ArrayView< int > &nodeIDs) const
Given a global ID list, return the list of their owning process IDs.
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...
LocalOrdinal getLocalElement(GlobalOrdinal globalIndex) const
The local index corresponding to the given global index.
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.
local_map_type getLocalMap() const
Get the local Map for Kokkos kernels.
Teuchos::ArrayView< const GlobalOrdinal > getNodeElementList() const
Return a NONOWNING view of the global indices owned by this process.
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...
size_t getNodeNumElements() const
The number of elements belonging to the calling process.
A parallel distribution of indices over processes.
bool isSameAs(const Map< LocalOrdinal, GlobalOrdinal, Node > &map) const
True if and only if map is identical to this Map.
Implement mapping from global ID to process ID and local ID.
Stand-alone utility functions and macros.
bool isOneToOne() const
Whether the Map is one to one.
void initialize(const map_type &map)
Initialize the Directory with its Map.
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.
typename Node::device_type device_type
This class&#39; Kokkos::Device specialization.
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.
LocalGlobal
Enum for local versus global allocation of Map entries.
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.
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 ...
Node node_type
Legacy typedef that will go away at some point.
global_size_t getGlobalNumElements() const
The number of elements in this Map.
std::string description() const
Implementation of Teuchos::Describable.
Map()
Default constructor (that does nothing).