Tpetra parallel linear algebra Version of the Day
Loading...
Searching...
No Matches
Tpetra_Import_Util2.hpp
Go to the documentation of this file.
1// @HEADER
2// *****************************************************************************
3// Tpetra: Templated Linear Algebra Services Package
4//
5// Copyright 2008 NTESS and the Tpetra contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#ifndef TPETRA_IMPORT_UTIL2_HPP
11#define TPETRA_IMPORT_UTIL2_HPP
12
17
18#include "Tpetra_ConfigDefs.hpp"
19#include "Tpetra_Import.hpp"
20#include "Tpetra_HashTable.hpp"
21#include "Tpetra_Map.hpp"
22#include "Tpetra_Util.hpp"
23#include "Tpetra_Distributor.hpp"
25#include "Tpetra_Vector.hpp"
26#include "Kokkos_DualView.hpp"
27#include "KokkosSparse_SortCrs.hpp"
28#include <Teuchos_Array.hpp>
31#include <Kokkos_Core.hpp>
32#include <Kokkos_UnorderedMap.hpp>
33#include <Kokkos_Sort.hpp>
34#include <utility>
35#include <set>
36
38
39namespace {
40
41struct LocalRemoteCount {
42 size_t numLocalColGIDs;
43 size_t numRemoteColGIDs;
44
45 KOKKOS_INLINE_FUNCTION
46 LocalRemoteCount()
47 : numLocalColGIDs(0)
48 , numRemoteColGIDs(0) {}
49
50 KOKKOS_INLINE_FUNCTION
51 LocalRemoteCount& operator+=(const LocalRemoteCount& src) {
52 numLocalColGIDs += src.numLocalColGIDs;
53 numRemoteColGIDs += src.numRemoteColGIDs;
54 return *this;
55 }
56};
57
58} // namespace
59
60namespace Kokkos {
61
62template <>
63struct reduction_identity<LocalRemoteCount> {
64 KOKKOS_FORCEINLINE_FUNCTION static LocalRemoteCount sum() {
65 return LocalRemoteCount();
66 }
67};
68} // namespace Kokkos
69
70namespace Tpetra {
71namespace Import_Util {
72
75template <typename Scalar, typename Ordinal>
76void sortCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
77 const Teuchos::ArrayView<Ordinal>& CRS_colind,
78 const Teuchos::ArrayView<Scalar>& CRS_vals,
79 const ::KokkosSparse::SortAlgorithm option = ::KokkosSparse::SortAlgorithm::DEFAULT);
80
81template <typename Ordinal>
82void sortCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
83 const Teuchos::ArrayView<Ordinal>& CRS_colind,
84 const ::KokkosSparse::SortAlgorithm option = ::KokkosSparse::SortAlgorithm::DEFAULT);
85
86template <typename rowptr_array_type, typename colind_array_type, typename vals_array_type>
87void sortCrsEntries(const rowptr_array_type& CRS_rowptr,
88 const colind_array_type& CRS_colind,
89 const vals_array_type& CRS_vals,
90 const ::KokkosSparse::SortAlgorithm option = ::KokkosSparse::SortAlgorithm::DEFAULT);
91
92template <typename rowptr_array_type, typename colind_array_type>
93void sortCrsEntries(const rowptr_array_type& CRS_rowptr,
94 const colind_array_type& CRS_colind,
95 const ::KokkosSparse::SortAlgorithm option = ::KokkosSparse::SortAlgorithm::DEFAULT);
96
101template <typename Scalar, typename Ordinal>
102void sortAndMergeCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
103 const Teuchos::ArrayView<Ordinal>& CRS_colind,
104 const Teuchos::ArrayView<Scalar>& CRS_vals);
105
106template <typename Ordinal>
107void sortAndMergeCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
108 const Teuchos::ArrayView<Ordinal>& CRS_colind);
109
110template <class rowptr_view_type, class colind_view_type, class vals_view_type>
111void sortAndMergeCrsEntries(rowptr_view_type& CRS_rowptr,
112 colind_view_type& CRS_colind,
113 vals_view_type& CRS_vals,
114 const ::KokkosSparse::SortAlgorithm option);
115
131template <typename LocalOrdinal, typename GlobalOrdinal, typename Node>
133 const Teuchos::ArrayView<const size_t>& rowptr,
134 const Teuchos::ArrayView<LocalOrdinal>& colind_LID,
135 const Teuchos::ArrayView<GlobalOrdinal>& colind_GID,
136 const Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& domainMapRCP,
137 const Teuchos::ArrayView<const int>& owningPIDs,
138 Teuchos::Array<int>& remotePIDs,
139 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& colMap,
140 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
141
146template <typename LocalOrdinal, typename GlobalOrdinal, typename Node>
148 const Kokkos::View<const size_t*, typename Node::device_type> rowptr_view,
149 const Kokkos::View<LocalOrdinal*, typename Node::device_type> colind_LID_view,
150 const Kokkos::View<GlobalOrdinal*, typename Node::device_type> colind_GID_view,
151 const Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& domainMapRCP,
152 const Kokkos::View<const int*, typename Node::device_type> owningPIDs_view,
153 Teuchos::Array<int>& remotePIDs,
154 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& colMap,
155 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null);
156
170template <typename LocalOrdinal, typename GlobalOrdinal, typename Node>
171void getTwoTransferOwnershipVector(const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node>& transferThatDefinesOwnership,
172 bool useReverseModeForOwnership,
173 const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node>& transferForMigratingData,
174 bool useReverseModeForMigration,
176
177} // namespace Import_Util
178} // namespace Tpetra
179
180//
181// Implementations
182//
183
184namespace Tpetra {
185namespace Import_Util {
186
187template <typename PID, typename GlobalOrdinal>
188bool sort_PID_then_GID(const std::pair<PID, GlobalOrdinal>& a,
189 const std::pair<PID, GlobalOrdinal>& b) {
190 if (a.first != b.first)
191 return (a.first < b.first);
192 return (a.second < b.second);
193}
194
195template <typename PID,
196 typename GlobalOrdinal,
197 typename LocalOrdinal>
198bool sort_PID_then_pair_GID_LID(const std::pair<PID, std::pair<GlobalOrdinal, LocalOrdinal>>& a,
199 const std::pair<PID, std::pair<GlobalOrdinal, LocalOrdinal>>& b) {
200 if (a.first != b.first)
201 return a.first < b.first;
202 else
203 return (a.second.first < b.second.first);
204}
205
206template <typename Scalar,
207 typename LocalOrdinal,
208 typename GlobalOrdinal,
209 typename Node>
210void reverseNeighborDiscovery(const CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>& SourceMatrix,
211 const typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::row_ptrs_host_view_type& rowptr,
212 const typename CrsMatrix<Scalar, LocalOrdinal, GlobalOrdinal, Node>::local_inds_host_view_type& colind,
214 Teuchos::RCP<const Tpetra::Import<LocalOrdinal, GlobalOrdinal, Node>> MyImporter,
215 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> MyDomainMap,
216 Teuchos::ArrayRCP<int>& type3PIDs,
217 Teuchos::ArrayRCP<LocalOrdinal>& type3LIDs,
218 Teuchos::RCP<const Teuchos::Comm<int>>& rcomm) {
219#ifdef HAVE_TPETRACORE_MPI
220 using ::Tpetra::Details::Behavior;
221 typedef LocalOrdinal LO;
222 typedef GlobalOrdinal GO;
223 typedef std::pair<GO, GO> pidgidpair_t;
224 using Teuchos::RCP;
225 const std::string prefix{" Import_Util2::ReverseND:: "};
226 const std::string label("IU2::Neighbor");
227
228 // There can be no neighbor discovery if you don't have an importer
229 if (MyImporter.is_null()) return;
230
231 std::ostringstream errstr;
232 bool error = false;
233 auto const comm = MyDomainMap->getComm();
234
235 MPI_Comm rawComm = getRawMpiComm(*comm);
236 const int MyPID = rcomm->getRank();
237
238 // Things related to messages I am sending in forward mode (RowTransfer)
239 // *** Note: this will be incorrect for transferAndFillComplete if it is in reverse mode. FIXME cbl.
240 auto ExportPIDs = RowTransfer.getExportPIDs();
241 auto ExportLIDs = RowTransfer.getExportLIDs();
242 auto NumExportLIDs = RowTransfer.getNumExportIDs();
243
244 Distributor& Distor = MyImporter->getDistributor();
245 const size_t NumRecvs = Distor.getNumReceives();
246 const size_t NumSends = Distor.getNumSends();
247 auto RemoteLIDs = MyImporter->getRemoteLIDs();
248 auto const ProcsFrom = Distor.getProcsFrom();
249 auto const ProcsTo = Distor.getProcsTo();
250
251 auto LengthsFrom = Distor.getLengthsFrom();
252 auto MyColMap = SourceMatrix.getColMap();
253 const size_t numCols = MyColMap->getLocalNumElements();
254 RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> target = MyImporter->getTargetMap();
255
256 // Get the owning pids in a special way,
257 // s.t. ProcsFrom[RemotePIDs[i]] is the proc that owns RemoteLIDs[j]....
258 Teuchos::Array<int> RemotePIDOrder(numCols, -1);
259
260 // For each remote ID, record index into ProcsFrom, who owns it.
261 for (size_t i = 0, j = 0; i < NumRecvs; ++i) {
262 for (size_t k = 0; k < LengthsFrom[i]; ++k) {
263 const int pid = ProcsFrom[i];
264 if (pid != MyPID) {
265 RemotePIDOrder[RemoteLIDs[j]] = i;
266 }
267 j++;
268 }
269 }
270
271 // Step One: Start tacking the (GID,PID) pairs on the std sets
272 //
273 // For each index in ProcsFrom, we will insert into a set of (PID,
274 // GID) pairs, in order to build a list of such pairs for each of
275 // those processes. Since this is building a reverse, we will send
276 // to these processes.
277 Teuchos::Array<int> ReverseSendSizes(NumRecvs, 0);
278 // do this as C array to avoid Teuchos::Array value initialization of all reserved memory
279 Teuchos::Array<Teuchos::ArrayRCP<pidgidpair_t>> RSB(NumRecvs);
280
281 {
282 Tpetra::Details::ProfilingRegion set_all("Import_Util2::ReverseND::isMMallSetRSB");
283
284 // 25 Jul 2018: CBL
285 // todo:std::unordered_set (hash table),
286 // with an adequate prereservation ("bucket count").
287 // An onordered_set has to have a custom hasher for pid/gid pair
288 // However, when pidsets is copied to RSB, it will be in key
289 // order _not_ in pid,gid order. (unlike std::set).
290 // Impliment this with a reserve, and time BOTH building pidsets
291 // _and_ the sort after the receive. Even if unordered_set saves
292 // time, if it causes the sort to be longer, it's not a win.
293
294 Teuchos::Array<std::set<pidgidpair_t>> pidsets(NumRecvs);
295 {
296 Tpetra::Details::ProfilingRegion set_insert("Import_Util2::ReverseND::isMMallSetRSBinsert");
297 for (size_t i = 0; i < NumExportLIDs; i++) {
298 LO lid = ExportLIDs[i];
299 GO exp_pid = ExportPIDs[i];
300 for (auto j = rowptr[lid]; j < rowptr[lid + 1]; j++) {
301 int pid_order = RemotePIDOrder[colind[j]];
302 if (pid_order != -1) {
303 GO gid = MyColMap->getGlobalElement(colind[j]); // Epetra SM.GCID46 =>sm->graph-> {colmap(colind)}
304 auto tpair = pidgidpair_t(exp_pid, gid);
305 pidsets[pid_order].insert(pidsets[pid_order].end(), tpair);
306 }
307 }
308 }
309 }
310
311 {
312 Tpetra::Details::ProfilingRegion set_cpy("Import_Util2::ReverseND::isMMallSetRSBcpy");
313 int jj = 0;
314 for (auto&& ps : pidsets) {
315 auto s = ps.size();
316 RSB[jj] = Teuchos::arcp(new pidgidpair_t[s], 0, s, true);
317 std::copy(ps.begin(), ps.end(), RSB[jj]);
318 ReverseSendSizes[jj] = s;
319 ++jj;
320 }
321 }
322 } // end of set based packing.
323
324 Teuchos::Array<int> ReverseRecvSizes(NumSends, -1);
325 Teuchos::Array<MPI_Request> rawBreq(ProcsFrom.size() + ProcsTo.size(), MPI_REQUEST_NULL);
326 // 25 Jul 2018: MPI_TAG_UB is the largest tag value; could be < 32768.
327 const int mpi_tag_base_ = 3;
328
329 int mpireq_idx = 0;
330 for (int i = 0; i < ProcsTo.size(); ++i) {
331 int Rec_Tag = mpi_tag_base_ + ProcsTo[i];
332 int* thisrecv = (int*)(&ReverseRecvSizes[i]);
333 MPI_Request rawRequest = MPI_REQUEST_NULL;
334 MPI_Irecv(const_cast<int*>(thisrecv),
335 1,
336 MPI_INT,
337 ProcsTo[i],
338 Rec_Tag,
339 rawComm,
340 &rawRequest);
341 rawBreq[mpireq_idx++] = rawRequest;
342 }
343 for (int i = 0; i < ProcsFrom.size(); ++i) {
344 int Send_Tag = mpi_tag_base_ + MyPID;
345 int* mysend = (int*)(&ReverseSendSizes[i]);
346 MPI_Request rawRequest = MPI_REQUEST_NULL;
347 MPI_Isend(mysend,
348 1,
349 MPI_INT,
350 ProcsFrom[i],
351 Send_Tag,
352 rawComm,
353 &rawRequest);
354 rawBreq[mpireq_idx++] = rawRequest;
355 }
356 Teuchos::Array<MPI_Status> rawBstatus(rawBreq.size());
357#ifdef HAVE_TPETRA_DEBUG
358 const int err1 =
359#endif
360 MPI_Waitall(rawBreq.size(), rawBreq.getRawPtr(),
361 rawBstatus.getRawPtr());
362
363#ifdef HAVE_TPETRA_DEBUG
364 if (err1) {
365 errstr << MyPID << "sE1 reverseNeighborDiscovery Mpi_Waitall error on send ";
366 error = true;
367 std::cerr << errstr.str() << std::flush;
368 }
369#endif
370
371 int totalexportpairrecsize = 0;
372 for (size_t i = 0; i < NumSends; ++i) {
373 totalexportpairrecsize += ReverseRecvSizes[i];
374#ifdef HAVE_TPETRA_DEBUG
375 if (ReverseRecvSizes[i] < 0) {
376 errstr << MyPID << "E4 reverseNeighborDiscovery: got < 0 for receive size " << ReverseRecvSizes[i] << std::endl;
377 error = true;
378 }
379#endif
380 }
381 Teuchos::ArrayRCP<pidgidpair_t> AllReverseRecv = Teuchos::arcp(new pidgidpair_t[totalexportpairrecsize], 0, totalexportpairrecsize, true);
382 int offset = 0;
383 mpireq_idx = 0;
384 for (int i = 0; i < ProcsTo.size(); ++i) {
385 int recv_data_size = ReverseRecvSizes[i] * 2;
386 int recvData_MPI_Tag = mpi_tag_base_ * 2 + ProcsTo[i];
387 MPI_Request rawRequest = MPI_REQUEST_NULL;
388 GO* rec_bptr = (GO*)(&AllReverseRecv[offset]);
389 offset += ReverseRecvSizes[i];
390 MPI_Irecv(rec_bptr,
391 recv_data_size,
392 ::Tpetra::Details::MpiTypeTraits<GO>::getType(rec_bptr[0]),
393 ProcsTo[i],
394 recvData_MPI_Tag,
395 rawComm,
396 &rawRequest);
397 rawBreq[mpireq_idx++] = rawRequest;
398 }
399 for (int ii = 0; ii < ProcsFrom.size(); ++ii) {
400 GO* send_bptr = (GO*)(RSB[ii].getRawPtr());
401 MPI_Request rawSequest = MPI_REQUEST_NULL;
402 int send_data_size = ReverseSendSizes[ii] * 2; // 2 == count of pair
403 int sendData_MPI_Tag = mpi_tag_base_ * 2 + MyPID;
404 MPI_Isend(send_bptr,
405 send_data_size,
406 ::Tpetra::Details::MpiTypeTraits<GO>::getType(send_bptr[0]),
407 ProcsFrom[ii],
408 sendData_MPI_Tag,
409 rawComm,
410 &rawSequest);
411
412 rawBreq[mpireq_idx++] = rawSequest;
413 }
414#ifdef HAVE_TPETRA_DEBUG
415 const int err =
416#endif
417 MPI_Waitall(rawBreq.size(),
418 rawBreq.getRawPtr(),
419 rawBstatus.getRawPtr());
420#ifdef HAVE_TPETRA_DEBUG
421 if (err) {
422 errstr << MyPID << "E3.r reverseNeighborDiscovery Mpi_Waitall error on receive ";
423 error = true;
424 std::cerr << errstr.str() << std::flush;
425 }
426#endif
427 std::sort(AllReverseRecv.begin(), AllReverseRecv.end(), Tpetra::Import_Util::sort_PID_then_GID<GlobalOrdinal, GlobalOrdinal>);
428
429 auto newEndOfPairs = std::unique(AllReverseRecv.begin(), AllReverseRecv.end());
430 // don't resize to remove non-unique, just use the end-of-unique iterator
431 if (AllReverseRecv.begin() == newEndOfPairs) return;
432 int ARRsize = std::distance(AllReverseRecv.begin(), newEndOfPairs);
433 auto rPIDs = Teuchos::arcp(new int[ARRsize], 0, ARRsize, true);
434 auto rLIDs = Teuchos::arcp(new LocalOrdinal[ARRsize], 0, ARRsize, true);
435
436 int tsize = 0;
437 for (auto itr = AllReverseRecv.begin(); itr != newEndOfPairs; ++itr) {
438 if ((int)(itr->first) != MyPID) {
439 rPIDs[tsize] = (int)itr->first;
440 LocalOrdinal lid = MyDomainMap->getLocalElement(itr->second);
441 rLIDs[tsize] = lid;
442 tsize++;
443 }
444 }
445
446 type3PIDs = rPIDs.persistingView(0, tsize);
447 type3LIDs = rLIDs.persistingView(0, tsize);
448
449 if (error) {
450 std::cerr << errstr.str() << std::flush;
451 comm->barrier();
452 comm->barrier();
453 comm->barrier();
454 MPI_Abort(MPI_COMM_WORLD, -1);
455 }
456#endif
457}
458
459// Note: This should get merged with the other Tpetra sort routines eventually.
460template <typename Scalar, typename Ordinal>
461void sortCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
462 const Teuchos::ArrayView<Ordinal>& CRS_colind,
463 const Teuchos::ArrayView<Scalar>& CRS_vals,
464 const ::KokkosSparse::SortAlgorithm option) {
465 auto rowptr_k = Kokkos::View<size_t*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(CRS_rowptr.data(), CRS_rowptr.size());
466 auto colind_k = Kokkos::View<Ordinal*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(CRS_colind.data(), CRS_colind.size());
467 auto vals_k = Kokkos::View<Scalar*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(CRS_vals.data(), CRS_vals.size());
469}
470
471template <typename Ordinal>
472void sortCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
473 const Teuchos::ArrayView<Ordinal>& CRS_colind,
474 const ::KokkosSparse::SortAlgorithm option) {
475 auto rowptr_k = Kokkos::View<size_t*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(CRS_rowptr.data(), CRS_rowptr.size());
476 auto colind_k = Kokkos::View<Ordinal*, Kokkos::HostSpace, Kokkos::MemoryUnmanaged>(CRS_colind.data(), CRS_colind.size());
477 sortCrsEntries(rowptr_k, colind_k, option);
478}
479
480template <typename rowptr_array_type, typename colind_array_type, typename vals_array_type>
481void sortCrsEntries(const rowptr_array_type& CRS_rowptr,
482 const colind_array_type& CRS_colind,
483 const vals_array_type& CRS_vals,
484 const ::KokkosSparse::SortAlgorithm option) {
485 KokkosSparse::sort_crs_matrix(CRS_rowptr, CRS_colind, CRS_vals,
486 KokkosKernels::ArithTraits<typename colind_array_type::non_const_value_type>::max(),
487 option);
488}
489
490template <typename rowptr_array_type, typename colind_array_type>
491void sortCrsEntries(const rowptr_array_type& CRS_rowptr,
492 const colind_array_type& CRS_colind,
493 const ::KokkosSparse::SortAlgorithm option) {
494 KokkosSparse::sort_crs_graph(CRS_rowptr, CRS_colind,
495 KokkosKernels::ArithTraits<typename colind_array_type::non_const_value_type>::max(),
496 option);
497}
498
499// Note: This should get merged with the other Tpetra sort routines eventually.
500template <typename Scalar, typename Ordinal>
501void sortAndMergeCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
502 const Teuchos::ArrayView<Ordinal>& CRS_colind,
503 const Teuchos::ArrayView<Scalar>& CRS_vals) {
504 // For each row, sort column entries from smallest to largest,
505 // merging column ids that are identify by adding values. Use shell
506 // sort. Stable sort so it is fast if indices are already sorted.
507 // Code copied from Epetra_CrsMatrix::SortEntries()
508
509 if (CRS_rowptr.size() == 0) {
510 return; // no rows, so nothing to sort
511 }
512 const size_t NumRows = CRS_rowptr.size() - 1;
513 const size_t nnz = CRS_colind.size();
514 size_t new_curr = CRS_rowptr[0];
515 size_t old_curr = CRS_rowptr[0];
516
517 const bool permute_values_array = CRS_vals.size() > 0;
518
519 for (size_t i = 0; i < NumRows; i++) {
520 const size_t old_rowptr_i = CRS_rowptr[i];
522 if (old_rowptr_i >= nnz) continue;
523
524 size_t NumEntries = CRS_rowptr[i + 1] - old_rowptr_i;
525 Teuchos::ArrayRCP<Scalar> locValues;
527 locValues = Teuchos::arcp<Scalar>(&CRS_vals[old_rowptr_i], 0, NumEntries, false);
528 Teuchos::ArrayRCP<Ordinal> locIndices(&CRS_colind[old_rowptr_i], 0, NumEntries, false);
529
530 // Sort phase
532 Ordinal m = n / 2;
533
534 while (m > 0) {
535 Ordinal max = n - m;
536 for (Ordinal j = 0; j < max; j++) {
537 for (Ordinal k = j; k >= 0; k -= m) {
538 if (locIndices[k + m] >= locIndices[k])
539 break;
541 Scalar dtemp = locValues[k + m];
542 locValues[k + m] = locValues[k];
543 locValues[k] = dtemp;
544 }
546 locIndices[k + m] = locIndices[k];
547 locIndices[k] = itemp;
548 }
549 }
550 m = m / 2;
551 }
552
553 // Merge & shrink
554 for (size_t j = old_rowptr_i; j < CRS_rowptr[i + 1]; j++) {
555 if (j > old_rowptr_i && CRS_colind[j] == CRS_colind[new_curr - 1]) {
557 } else if (new_curr == j) {
558 new_curr++;
559 } else {
562 new_curr++;
563 }
564 }
566 }
567
569}
570
571template <typename Ordinal>
572void sortAndMergeCrsEntries(const Teuchos::ArrayView<size_t>& CRS_rowptr,
573 const Teuchos::ArrayView<Ordinal>& CRS_colind) {
574 Teuchos::ArrayView<Tpetra::Details::DefaultTypes::scalar_type> CRS_vals;
576}
577
578template <class rowptr_view_type, class colind_view_type, class vals_view_type>
579void sortAndMergeCrsEntries(rowptr_view_type& CRS_rowptr,
580 colind_view_type& CRS_colind,
581 vals_view_type& CRS_vals,
582 const ::KokkosSparse::SortAlgorithm option) {
583 using execution_space = typename vals_view_type::execution_space;
584
585 auto CRS_rowptr_in = CRS_rowptr;
586 auto CRS_colind_in = CRS_colind;
587 auto CRS_vals_in = CRS_vals;
588
589 KokkosSparse::sort_and_merge_matrix<execution_space, rowptr_view_type,
590 colind_view_type, vals_view_type>(CRS_rowptr_in, CRS_colind_in, CRS_vals_in,
591 CRS_rowptr, CRS_colind, CRS_vals,
592 KokkosKernels::ArithTraits<typename colind_view_type::non_const_value_type>::max(),
593 option);
594}
595
596template <typename LocalOrdinal, typename GlobalOrdinal, typename Node>
597void lowCommunicationMakeColMapAndReindexSerial(const Teuchos::ArrayView<const size_t>& rowptr,
598 const Teuchos::ArrayView<LocalOrdinal>& colind_LID,
599 const Teuchos::ArrayView<GlobalOrdinal>& colind_GID,
600 const Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& domainMapRCP,
601 const Teuchos::ArrayView<const int>& owningPIDs,
602 Teuchos::Array<int>& remotePIDs,
603 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& colMap,
604 const Teuchos::RCP<Teuchos::ParameterList>& params = Teuchos::null) {
605 using Teuchos::rcp;
606 typedef LocalOrdinal LO;
607 typedef GlobalOrdinal GO;
608 typedef Tpetra::global_size_t GST;
609 typedef Tpetra::Map<LO, GO, Node> map_type;
610 const char prefix[] = "lowCommunicationMakeColMapAndReindexSerial: ";
611
612 // The domainMap is an RCP because there is a shortcut for a
613 // (common) special case to return the columnMap = domainMap.
614 const map_type& domainMap = *domainMapRCP;
615
616 // Scan all column indices and sort into two groups:
617 // Local: those whose GID matches a GID of the domain map on this processor and
618 // Remote: All others.
619 const size_t numDomainElements = domainMap.getLocalNumElements();
620 Teuchos::Array<bool> LocalGIDs;
621 if (numDomainElements > 0) {
622 LocalGIDs.resize(numDomainElements, false); // Assume domain GIDs are not local
623 }
624
625 // In principle it is good to have RemoteGIDs and RemotGIDList be as
626 // long as the number of remote GIDs on this processor, but this
627 // would require two passes through the column IDs, so we make it
628 // the max of 100 and the number of block rows.
629 //
630 // FIXME (mfh 11 Feb 2015) Tpetra::Details::HashTable can hold at
631 // most INT_MAX entries, but it's possible to have more rows than
632 // that (if size_t is 64 bits and int is 32 bits).
633 const size_t numMyRows = rowptr.size() - 1;
634 const int hashsize = std::max(static_cast<int>(numMyRows), 100);
635
636 Tpetra::Details::HashTable<GO, LO> RemoteGIDs(hashsize);
637 Teuchos::Array<GO> RemoteGIDList;
638 RemoteGIDList.reserve(hashsize);
639 Teuchos::Array<int> PIDList;
640 PIDList.reserve(hashsize);
641
642 // Here we start using the *LocalOrdinal* colind_LID array. This is
643 // safe even if both columnIndices arrays are actually the same
644 // (because LocalOrdinal==GO). For *local* GID's set
645 // colind_LID with with their LID in the domainMap. For *remote*
646 // GIDs, we set colind_LID with (numDomainElements+NumRemoteColGIDs)
647 // before the increment of the remote count. These numberings will
648 // be separate because no local LID is greater than
649 // numDomainElements.
650
651 size_t NumLocalColGIDs = 0;
652 LO NumRemoteColGIDs = 0;
653 for (size_t i = 0; i < numMyRows; ++i) {
654 for (size_t j = rowptr[i]; j < rowptr[i + 1]; ++j) {
655 const GO GID = colind_GID[j];
656 // Check if GID matches a row GID
657 const LO LID = domainMap.getLocalElement(GID);
658 if (LID != -1) {
659 const bool alreadyFound = LocalGIDs[LID];
660 if (!alreadyFound) {
661 LocalGIDs[LID] = true; // There is a column in the graph associated with this domain map GID
662 NumLocalColGIDs++;
663 }
664 colind_LID[j] = LID;
665 } else {
666 const LO hash_value = RemoteGIDs.get(GID);
667 if (hash_value == -1) { // This means its a new remote GID
668 const int PID = owningPIDs[j];
669 TEUCHOS_TEST_FOR_EXCEPTION(
670 PID == -1, std::invalid_argument, prefix << "Cannot figure out if "
671 "PID is owned.");
672 colind_LID[j] = static_cast<LO>(numDomainElements + NumRemoteColGIDs);
673 RemoteGIDs.add(GID, NumRemoteColGIDs);
674 RemoteGIDList.push_back(GID);
675 PIDList.push_back(PID);
676 NumRemoteColGIDs++;
677 } else {
678 colind_LID[j] = static_cast<LO>(numDomainElements + hash_value);
679 }
680 }
681 }
682 }
683
684 // Possible short-circuit: If all domain map GIDs are present as
685 // column indices, then set ColMap=domainMap and quit.
686 if (domainMap.getComm()->getSize() == 1) {
687 // Sanity check: When there is only one process, there can be no
688 // remoteGIDs.
689 TEUCHOS_TEST_FOR_EXCEPTION(
690 NumRemoteColGIDs != 0, std::runtime_error, prefix << "There is only one "
691 "process in the domain Map's communicator, which means that there are no "
692 "\"remote\" indices. Nevertheless, some column indices are not in the "
693 "domain Map.");
694 if (static_cast<size_t>(NumLocalColGIDs) == numDomainElements) {
695 // In this case, we just use the domainMap's indices, which is,
696 // not coincidently, what we clobbered colind with up above
697 // anyway. No further reindexing is needed.
698 colMap = domainMapRCP;
699 return;
700 }
701 }
702
703 // Launch reduction to get global num entries in the column map
704 const LO numMyCols = NumLocalColGIDs + NumRemoteColGIDs;
705 GST numMyColsGST = static_cast<GST>(numMyCols);
706 GST numGlobalCols;
707 auto req = Details::iallreduce(numMyColsGST, numGlobalCols, Teuchos::REDUCE_SUM, *domainMap.getComm());
708
709 // Now build the array containing column GIDs
710 // Build back end, containing remote GIDs, first
711 Teuchos::Array<GO> ColIndices;
712 GO* RemoteColIndices = NULL;
713 if (numMyCols > 0) {
714 ColIndices.resize(numMyCols);
715 if (NumLocalColGIDs != static_cast<size_t>(numMyCols)) {
716 RemoteColIndices = &ColIndices[NumLocalColGIDs]; // Points to back half of ColIndices
717 }
718 }
719
720 for (LO i = 0; i < NumRemoteColGIDs; ++i) {
721 RemoteColIndices[i] = RemoteGIDList[i];
722 }
723
724 // Build permute array for *remote* reindexing.
725 Teuchos::Array<LO> RemotePermuteIDs(NumRemoteColGIDs);
726 for (LO i = 0; i < NumRemoteColGIDs; ++i) {
727 RemotePermuteIDs[i] = i;
728 }
729
730 // Sort External column indices so that all columns coming from a
731 // given remote processor are contiguous. This is a sort with two
732 // auxilary arrays: RemoteColIndices and RemotePermuteIDs.
733 Tpetra::sort3(PIDList.begin(), PIDList.end(),
734 ColIndices.begin() + NumLocalColGIDs,
735 RemotePermuteIDs.begin());
736
737 // Stash the RemotePIDs.
738 //
739 // Note: If Teuchos::Array had a shrink_to_fit like std::vector,
740 // we'd call it here.
741 remotePIDs = PIDList;
742
743 // Sort external column indices so that columns from a given remote
744 // processor are not only contiguous but also in ascending
745 // order. NOTE: I don't know if the number of externals associated
746 // with a given remote processor is known at this point ... so I
747 // count them here.
748
749 // NTS: Only sort the RemoteColIndices this time...
750 LO StartCurrent = 0, StartNext = 1;
751 while (StartNext < NumRemoteColGIDs) {
752 if (PIDList[StartNext] == PIDList[StartNext - 1]) {
753 StartNext++;
754 } else {
755 Tpetra::sort2(ColIndices.begin() + NumLocalColGIDs + StartCurrent,
756 ColIndices.begin() + NumLocalColGIDs + StartNext,
757 RemotePermuteIDs.begin() + StartCurrent);
758 StartCurrent = StartNext;
759 StartNext++;
760 }
761 }
762 Tpetra::sort2(ColIndices.begin() + NumLocalColGIDs + StartCurrent,
763 ColIndices.begin() + NumLocalColGIDs + StartNext,
764 RemotePermuteIDs.begin() + StartCurrent);
765
766 // Reverse the permutation to get the information we actually care about
767 Teuchos::Array<LO> ReverseRemotePermuteIDs(NumRemoteColGIDs);
768 for (LO i = 0; i < NumRemoteColGIDs; ++i) {
769 ReverseRemotePermuteIDs[RemotePermuteIDs[i]] = i;
770 }
771
772 // Build permute array for *local* reindexing.
773 bool use_local_permute = false;
774 Teuchos::Array<LO> LocalPermuteIDs(numDomainElements);
775
776 // Now fill front end. Two cases:
777 //
778 // (1) If the number of Local column GIDs is the same as the number
779 // of Local domain GIDs, we can simply read the domain GIDs into
780 // the front part of ColIndices, otherwise
781 //
782 // (2) We step through the GIDs of the domainMap, checking to see if
783 // each domain GID is a column GID. we want to do this to
784 // maintain a consistent ordering of GIDs between the columns
785 // and the domain.
786 Teuchos::ArrayView<const GO> domainGlobalElements = domainMap.getLocalElementList();
787 if (static_cast<size_t>(NumLocalColGIDs) == numDomainElements) {
788 if (NumLocalColGIDs > 0) {
789 // Load Global Indices into first numMyCols elements column GID list
790 std::copy(domainGlobalElements.begin(), domainGlobalElements.end(),
791 ColIndices.begin());
792 }
793 } else {
794 LO NumLocalAgain = 0;
795 use_local_permute = true;
796 for (size_t i = 0; i < numDomainElements; ++i) {
797 if (LocalGIDs[i]) {
798 LocalPermuteIDs[i] = NumLocalAgain;
799 ColIndices[NumLocalAgain++] = domainGlobalElements[i];
800 }
801 }
802 TEUCHOS_TEST_FOR_EXCEPTION(
803 static_cast<size_t>(NumLocalAgain) != NumLocalColGIDs,
804 std::runtime_error, prefix << "Local ID count test failed.");
805 }
806
807 // Make column Map
808 req->wait();
809 colMap = rcp(new map_type(numGlobalCols, ColIndices, domainMap.getIndexBase(),
810 domainMap.getComm(),
811 params));
812
813 // Low-cost reindex of the matrix
814 for (size_t i = 0; i < numMyRows; ++i) {
815 for (size_t j = rowptr[i]; j < rowptr[i + 1]; ++j) {
816 const LO ID = colind_LID[j];
817 if (static_cast<size_t>(ID) < numDomainElements) {
818 if (use_local_permute) {
819 colind_LID[j] = LocalPermuteIDs[colind_LID[j]];
820 }
821 // In the case where use_local_permute==false, we just copy
822 // the DomainMap's ordering, which it so happens is what we
823 // put in colind_LID to begin with.
824 } else {
825 colind_LID[j] = NumLocalColGIDs + ReverseRemotePermuteIDs[colind_LID[j] - numDomainElements];
826 }
827 }
828 }
829}
830
831template <typename LocalOrdinal, typename GlobalOrdinal, typename Node>
833 const Teuchos::ArrayView<const size_t>& rowptr,
834 const Teuchos::ArrayView<LocalOrdinal>& colind_LID,
835 const Teuchos::ArrayView<GlobalOrdinal>& colind_GID,
837 const Teuchos::ArrayView<const int>& owningPIDs,
838 Teuchos::Array<int>& remotePIDs,
839 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& colMap,
840 const Teuchos::RCP<Teuchos::ParameterList>& params) {
841 using DT = typename Node::device_type;
842 using execution_space = typename DT::execution_space;
843 execution_space exec;
844 DT device;
845
846 // Create device mirror and host mirror views from function parameters
847 // When we pass in views instead of Teuchos::ArrayViews, we can avoid copying views
848 auto colind_LID_view = Details::create_mirror_view_from_raw_host_array(device, colind_LID.getRawPtr(), colind_LID.size(), true, "colind_LID");
849 auto rowptr_view = Details::create_mirror_view_from_raw_host_array(device, rowptr.getRawPtr(), rowptr.size(), true, "rowptr");
850 auto colind_GID_view = Details::create_mirror_view_from_raw_host_array(device, colind_GID.getRawPtr(), colind_GID.size(), true, "colind_GID");
851 auto owningPIDs_view = Details::create_mirror_view_from_raw_host_array(device, owningPIDs.getRawPtr(), owningPIDs.size(), true, "owningPIDs");
852
853 typename decltype(colind_LID_view)::host_mirror_type colind_LID_host(colind_LID.getRawPtr(), colind_LID.size());
854
856
857 // For now, we copy back into colind_LID_host (which also overwrites the colind_LID Teuchos array)
858 // When colind_LID becomes a Kokkos View we can delete this
859 Kokkos::deep_copy(exec, colind_LID_host, colind_LID_view);
860}
861
862template <typename LocalOrdinal, typename GlobalOrdinal, typename Node>
864 const Kokkos::View<const size_t*, typename Node::device_type> rowptr_view,
865 const Kokkos::View<LocalOrdinal*, typename Node::device_type> colind_LID_view,
866 const Kokkos::View<GlobalOrdinal*, typename Node::device_type> colind_GID_view,
868 const Kokkos::View<const int*, typename Node::device_type> owningPIDs_view,
869 Teuchos::Array<int>& remotePIDs,
870 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>>& colMap,
871 const Teuchos::RCP<Teuchos::ParameterList>& params) {
872 using Teuchos::rcp;
873 typedef LocalOrdinal LO;
874 typedef GlobalOrdinal GO;
876 typedef Tpetra::Map<LO, GO, Node> map_type;
877 const char prefix[] = "lowCommunicationMakeColMapAndReindex: ";
878 Tpetra::Details::ProfilingRegion regOuter("lowCommunicationMakeColMapAndReindex");
879
880 typedef typename Node::device_type DT;
881 using execution_space = typename DT::execution_space;
882 using memory_space = typename DT::memory_space;
883 execution_space exec;
884 using team_policy = Kokkos::TeamPolicy<execution_space, Kokkos::Schedule<Kokkos::Dynamic>>;
885 typedef typename map_type::local_map_type local_map_type;
886
887 // The domainMap is an RCP because there is a shortcut for a
888 // (common) special case to return the columnMap = domainMap.
889 const map_type& domainMap = *domainMapRCP;
890
891 const size_t numMyRows = rowptr_view.size() - 1;
892 local_map_type domainMap_local = domainMap.getLocalMap();
893
894 const size_t numDomainElements = domainMap.getLocalNumElements();
895 Kokkos::View<bool*, DT> LocalGIDs_view("LocalGIDs", numDomainElements);
896
897 auto map_capacity = std::max(std::min(uint32_t(1.2 * colind_LID_view.size()),
899 uint32_t(5000));
900 Kokkos::UnorderedMap<GO, int, DT> GIDs_map(map_capacity);
901 size_t NumLocalColGIDs;
902 size_t NumRemoteColGIDs;
903 while (true) {
904 Tpetra::Details::ProfilingRegion reg("lowCommunicationMakeColMapAndReindex::fill_map");
905 const LO LINVALID = Teuchos::OrdinalTraits<LO>::invalid();
906
907 // Scan all column indices and sort into two groups:
908 // Local: those whose GID matches a GID of the domain map on this processor and
909 // Remote: All others.
910 // Kokkos::Parallel_reduce sums up NumLocalColGIDs, while we use the size of the Remote GIDs map to find NumRemoteColGIDs
911 LocalRemoteCount counts;
912 Kokkos::parallel_reduce(
913 team_policy(numMyRows, Kokkos::AUTO), KOKKOS_LAMBDA(const typename team_policy::member_type& member, LocalRemoteCount& update) {
914 const int i = member.league_rank();
915 LocalRemoteCount NumLocalColGIDs_temp;
916 size_t rowptr_start = rowptr_view[i];
917 size_t rowptr_end = rowptr_view[i + 1];
918 Kokkos::parallel_reduce(
919 Kokkos::TeamThreadRange(member, rowptr_start, rowptr_end), [&](const size_t j, LocalRemoteCount& innerUpdate) {
920 const GO GID = colind_GID_view[j];
921 // Check if GID matches a row GID in local domain map
922 const LO LID = domainMap_local.getLocalElement(GID);
923 if (LID != LINVALID) {
924 auto outcome = GIDs_map.insert(GID, -1);
925 // Fresh insert
926 if (outcome.success()) {
927 LocalGIDs_view[LID] = true;
928 innerUpdate.numLocalColGIDs++;
929 }
930 } else {
931 const int PID = owningPIDs_view[j];
932 auto outcome = GIDs_map.insert(GID, PID);
933 if (outcome.success()) {
934 innerUpdate.numRemoteColGIDs++;
935 if (PID == -1) {
936 Kokkos::abort("Cannot figure out if ID is owned.\n");
937 }
938 }
939 }
940 },
942 if (member.team_rank() == 0) update += NumLocalColGIDs_temp;
943 },
944 counts);
945 if (!GIDs_map.failed_insert()) {
946 NumLocalColGIDs = counts.numLocalColGIDs;
947 NumRemoteColGIDs = counts.numRemoteColGIDs;
948 break;
949 } else {
950 map_capacity *= 1.5;
951 GIDs_map = Kokkos::UnorderedMap<GO, int, DT>(map_capacity);
952 }
953 }
954
955 // Possible short-circuit: If all domain map GIDs are present as
956 // column indices, then set ColMap=domainMap and quit.
957 if (domainMap.getComm()->getSize() == 1) {
958 // Sanity check: When there is only one process, there can be no
959 // remoteGIDs.
961 NumRemoteColGIDs != 0, std::runtime_error, prefix << "There is only one "
962 "process in the domain Map's communicator, which means that there are no "
963 "\"remote\" indices. Nevertheless, some column indices are not in the "
964 "domain Map.");
965 if (static_cast<size_t>(NumLocalColGIDs) == numDomainElements) {
966 // In this case, we just use the domainMap's indices, which is,
967 // not coincidently, what we clobbered colind with up above
968 // anyway. No further reindexing is needed.
969 colMap = domainMapRCP;
970
971 // Fill out local colMap (which should only contain local GIDs)
972 auto localColMap = colMap->getLocalMap();
973 Kokkos::parallel_for(
974 Kokkos::RangePolicy<execution_space>(0, colind_GID_view.size()), KOKKOS_LAMBDA(const int i) {
975 colind_LID_view(i) = localColMap.getLocalElement(colind_GID_view(i));
976 });
977 return;
978 }
979 }
980
981 Kokkos::View<GO*, DT> RemoteColMapIndices_unsorted(Kokkos::ViewAllocateWithoutInitializing("RemoteColMapIndices_unsorted"), NumRemoteColGIDs);
982 Kokkos::View<int*, DT> RemotePIDs_unsorted(Kokkos::ViewAllocateWithoutInitializing("RemotePIDs_unsorted"), NumRemoteColGIDs);
983
984 // For each index in RemoteGIDs_map that contains a GID, use "update" to indicate the number of GIDs "before" this GID
985 // This maps each element in the RemoteGIDs hash table to an index in RemoteGIDList / PIDList without any overwriting or empty spaces between indices
986
987 {
988 Tpetra::Details::ProfilingRegion reg("lowCommunicationMakeColMapAndReindex::process_map");
989 size_t numEnteredRemotes;
990 Kokkos::parallel_scan(
991 Kokkos::RangePolicy<execution_space>(0, GIDs_map.capacity()), KOKKOS_LAMBDA(const int i, size_t& update, const bool final) {
992 if (GIDs_map.valid_at(i) && GIDs_map.value_at(i) != -1) {
993 if (final) {
994 RemoteColMapIndices_unsorted(update) = GIDs_map.key_at(i);
995 RemotePIDs_unsorted(update) = GIDs_map.value_at(i);
996 }
997 ++update;
998 }
999 },
1000 numEnteredRemotes);
1001 TEUCHOS_ASSERT(numEnteredRemotes == NumRemoteColGIDs);
1002 }
1003
1004 // Launch reduction to get global num entries in the column map
1005 const size_t numMyCols = NumLocalColGIDs + NumRemoteColGIDs;
1006 GST numMyColsGST = static_cast<GST>(numMyCols);
1007 GST numGlobalCols;
1008 auto req = Details::iallreduce(numMyColsGST, numGlobalCols, Teuchos::REDUCE_SUM, *domainMap.getComm());
1009
1010 // Now build the array containing column GIDs
1011 // Build back end, containing remote GIDs, first
1012 Kokkos::View<GO*, DT> ColMapIndices(Kokkos::ViewAllocateWithoutInitializing("ColMapIndices"), numMyCols);
1013
1014 // We don't need to load the back end of ColMapIndices or sort if there are no remote GIDs
1015 if (NumRemoteColGIDs > 0) {
1016 Tpetra::Details::ProfilingRegion reg("lowCommunicationMakeColMapAndReindex::sort");
1017 auto RemoteColMapIndices_sorted = Kokkos::subview(ColMapIndices, Kokkos::make_pair(NumLocalColGIDs, numMyCols));
1018
1019 // Create a permutation that sorts by first by PID and then by GID.
1020 Kokkos::View<size_t*, memory_space> index(Kokkos::ViewAllocateWithoutInitializing("index"), NumRemoteColGIDs);
1021 Kokkos::parallel_for(
1022 Kokkos::RangePolicy<execution_space>(0, NumRemoteColGIDs), KOKKOS_LAMBDA(const size_t i) {
1023 index(i) = i;
1024 });
1025 Kokkos::sort(
1026 exec, index, KOKKOS_LAMBDA(const size_t i, const size_t j) {
1027 auto pid_i = RemotePIDs_unsorted(i);
1028 auto pid_j = RemotePIDs_unsorted(j);
1029 if (pid_i == pid_j)
1030 return RemoteColMapIndices_unsorted(i) < RemoteColMapIndices_unsorted(j);
1031 else
1032 return pid_i < pid_j;
1033 });
1034 // Apply the permutation to remote GIDs and PIDs.
1035 Kokkos::View<int*, memory_space> RemotePIDs_sorted(Kokkos::ViewAllocateWithoutInitializing("RemotePIDs_new"), NumRemoteColGIDs);
1036 Kokkos::parallel_for(
1037 Kokkos::RangePolicy<execution_space>(0, NumRemoteColGIDs), KOKKOS_LAMBDA(const size_t i) {
1038 RemoteColMapIndices_sorted(i) = RemoteColMapIndices_unsorted(index(i));
1039 RemotePIDs_sorted(i) = RemotePIDs_unsorted(index(i));
1040 });
1041
1042 // Deep copy back from device to host
1043 // Stash the RemotePIDs. Once remotePIDs is changed to become a Kokkos view, we can remove this and copy directly.
1044 Teuchos::Array<int> PIDList(NumRemoteColGIDs);
1045 Kokkos::View<int*, Kokkos::HostSpace> RemotePIDs_sorted_host(PIDList.data(), NumRemoteColGIDs);
1046 Kokkos::deep_copy(exec, RemotePIDs_sorted_host, RemotePIDs_sorted);
1047 exec.fence();
1048
1049 remotePIDs = PIDList;
1050 }
1051
1052 // Build permute array for *local* reindexing.
1053
1054 // Now fill front end. Two cases:
1055 //
1056 // (1) If the number of Local column GIDs is the same as the number
1057 // of Local domain GIDs, we can simply read the domain GIDs into
1058 // the front part of ColIndices, otherwise
1059 //
1060 // (2) We step through the GIDs of the domainMap, checking to see if
1061 // each domain GID is a column GID. we want to do this to
1062 // maintain a consistent ordering of GIDs between the columns
1063 // and the domain.
1064 if (static_cast<size_t>(NumLocalColGIDs) == numDomainElements) {
1065 if (NumLocalColGIDs > 0) {
1066 // Load Global Indices into first numMyCols elements column GID list
1067 auto LocalColMapIndices = Kokkos::subview(ColMapIndices, Kokkos::make_pair((size_t)0, numDomainElements));
1068 Kokkos::deep_copy(LocalColMapIndices, domainMap.getMyGlobalIndicesDevice());
1069 }
1070 } else {
1071 // This part isn't actually tested in the unit tests
1072 LO NumLocalAgain = 0;
1073 Kokkos::parallel_scan(
1074 Kokkos::RangePolicy<execution_space>(0, numDomainElements), KOKKOS_LAMBDA(const int i, LO& update, const bool final) {
1075 if (LocalGIDs_view(i)) {
1076 if (final) {
1077 ColMapIndices(update) = domainMap_local.getGlobalElement(i);
1078 }
1079 update++;
1080 }
1081 },
1082 NumLocalAgain);
1083
1084 TEUCHOS_TEST_FOR_EXCEPTION(
1085 static_cast<size_t>(NumLocalAgain) != NumLocalColGIDs,
1086 std::runtime_error, prefix << "Local ID count test failed.");
1087 }
1088
1089 // Make column Map
1090 req->wait();
1091 colMap = rcp(new map_type(numGlobalCols, ColMapIndices, domainMap.getIndexBase(),
1092 domainMap.getComm(), params));
1093
1094 // Fill out colind_LID using local map
1095 auto localColMap = colMap->getLocalMap();
1096 Kokkos::parallel_for(
1097 Kokkos::RangePolicy<execution_space>(0, colind_GID_view.size()), KOKKOS_LAMBDA(const int i) {
1098 colind_LID_view(i) = localColMap.getLocalElement(colind_GID_view(i));
1099 });
1100}
1101
1102// Generates an list of owning PIDs based on two transfer (aka import/export objects)
1103// Let:
1104// OwningMap = useReverseModeForOwnership ? transferThatDefinesOwnership.getTargetMap() : transferThatDefinesOwnership.getSourceMap();
1105// MapAo = useReverseModeForOwnership ? transferThatDefinesOwnership.getSourceMap() : transferThatDefinesOwnership.getTargetMap();
1106// MapAm = useReverseModeForMigration ? transferThatDefinesMigration.getTargetMap() : transferThatDefinesMigration.getSourceMap();
1107// VectorMap = useReverseModeForMigration ? transferThatDefinesMigration.getSourceMap() : transferThatDefinesMigration.getTargetMap();
1108// Precondition:
1109// 1) MapAo.isSameAs(*MapAm) - map compatibility between transfers
1110// 2) VectorMap->isSameAs(*owningPIDs->getMap()) - map compabibility between transfer & vector
1111// 3) OwningMap->isOneToOne() - owning map is 1-to-1
1112// --- Precondition 3 is only checked in DEBUG mode ---
1113// Postcondition:
1114// owningPIDs[VectorMap->getLocalElement(GID i)] = j iff (OwningMap->isLocalElement(GID i) on rank j)
1115template <typename LocalOrdinal, typename GlobalOrdinal, typename Node>
1116void getTwoTransferOwnershipVector(const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node>& transferThatDefinesOwnership,
1118 const ::Tpetra::Details::Transfer<LocalOrdinal, GlobalOrdinal, Node>& transferThatDefinesMigration,
1123
1124 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> OwningMap = useReverseModeForOwnership ? transferThatDefinesOwnership.getTargetMap() : transferThatDefinesOwnership.getSourceMap();
1125 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> MapAo = useReverseModeForOwnership ? transferThatDefinesOwnership.getSourceMap() : transferThatDefinesOwnership.getTargetMap();
1126 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> MapAm = useReverseModeForMigration ? transferThatDefinesMigration.getTargetMap() : transferThatDefinesMigration.getSourceMap();
1127 Teuchos::RCP<const Tpetra::Map<LocalOrdinal, GlobalOrdinal, Node>> VectorMap = useReverseModeForMigration ? transferThatDefinesMigration.getSourceMap() : transferThatDefinesMigration.getTargetMap();
1128
1129 TEUCHOS_TEST_FOR_EXCEPTION(!MapAo->isSameAs(*MapAm), std::runtime_error, "Tpetra::Import_Util::getTwoTransferOwnershipVector map mismatch between transfers");
1130 TEUCHOS_TEST_FOR_EXCEPTION(!VectorMap->isSameAs(*owningPIDs.getMap()), std::runtime_error, "Tpetra::Import_Util::getTwoTransferOwnershipVector map mismatch transfer and vector");
1131#ifdef HAVE_TPETRA_DEBUG
1132 TEUCHOS_TEST_FOR_EXCEPTION(!OwningMap->isOneToOne(), std::runtime_error, "Tpetra::Import_Util::getTwoTransferOwnershipVector owner must be 1-to-1");
1133#endif
1134
1135 int rank = OwningMap->getComm()->getRank();
1136 // Generate "A" vector and fill it with owning information. We can read this from transferThatDefinesOwnership w/o communication
1137 // Note: Due to the 1-to-1 requirement, several of these options throw
1139 const import_type* ownAsImport = dynamic_cast<const import_type*>(&transferThatDefinesOwnership);
1140 const export_type* ownAsExport = dynamic_cast<const export_type*>(&transferThatDefinesOwnership);
1141
1142 Teuchos::ArrayRCP<int> pids = temp.getDataNonConst();
1143 Teuchos::ArrayView<int> v_pids = pids();
1145 TEUCHOS_TEST_FOR_EXCEPTION(1, std::runtime_error, "Tpetra::Import_Util::getTwoTransferOwnershipVector owner must be 1-to-1");
1147 getPids(*ownAsImport, v_pids, false);
1149 TEUCHOS_TEST_FOR_EXCEPTION(1, std::runtime_error, "Tpetra::Import_Util::getTwoTransferOwnershipVector this option not yet implemented");
1150 } else {
1151 TEUCHOS_TEST_FOR_EXCEPTION(1, std::runtime_error, "Tpetra::Import_Util::getTwoTransferOwnershipVector owner must be 1-to-1");
1152 }
1153
1154 const import_type* xferAsImport = dynamic_cast<const import_type*>(&transferThatDefinesMigration);
1155 const export_type* xferAsExport = dynamic_cast<const export_type*>(&transferThatDefinesMigration);
1156 TEUCHOS_TEST_FOR_EXCEPTION(!xferAsImport && !xferAsExport, std::runtime_error, "Tpetra::Import_Util::getTwoTransferOwnershipVector transfer undefined");
1157
1158 // Migrate from "A" vector to output vector
1159 owningPIDs.putScalar(rank);
1166 else
1168}
1169
1170} // namespace Import_Util
1171} // namespace Tpetra
1172
1173#endif // TPETRA_IMPORT_UTIL_HPP
Declaration of the Tpetra::CrsMatrix class.
Add specializations of Teuchos::Details::MpiTypeTraits for Kokkos::complex<float> and Kokkos::complex...
Declaration of Tpetra::Details::Profiling, a scope guard for Kokkos Profiling.
Functions that wrap Kokkos::create_mirror_view, in order to avoid deep copies when not necessary,...
void sortCrsEntries(const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< Ordinal > &CRS_colind, const Teuchos::ArrayView< Scalar > &CRS_vals, const ::KokkosSparse::SortAlgorithm option=::KokkosSparse::SortAlgorithm::DEFAULT)
Sort the entries of the (raw CSR) matrix by column index within each row.
void getTwoTransferOwnershipVector(const ::Tpetra::Details::Transfer< LocalOrdinal, GlobalOrdinal, Node > &transferThatDefinesOwnership, bool useReverseModeForOwnership, const ::Tpetra::Details::Transfer< LocalOrdinal, GlobalOrdinal, Node > &transferForMigratingData, bool useReverseModeForMigration, Tpetra::Vector< int, LocalOrdinal, GlobalOrdinal, Node > &owningPIDs)
Generates an list of owning PIDs based on two transfer (aka import/export objects) Let: OwningMap = u...
void sortAndMergeCrsEntries(const Teuchos::ArrayView< size_t > &CRS_rowptr, const Teuchos::ArrayView< Ordinal > &CRS_colind, const Teuchos::ArrayView< Scalar > &CRS_vals)
Sort and merge the entries of the (raw CSR) matrix by column index within each row.
void lowCommunicationMakeColMapAndReindex(const Teuchos::ArrayView< const size_t > &rowptr, const Teuchos::ArrayView< LocalOrdinal > &colind_LID, const Teuchos::ArrayView< GlobalOrdinal > &colind_GID, const Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &domainMapRCP, const Teuchos::ArrayView< const int > &owningPIDs, Teuchos::Array< int > &remotePIDs, Teuchos::RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &colMap, const Teuchos::RCP< Teuchos::ParameterList > &params=Teuchos::null)
lowCommunicationMakeColMapAndReindex
void getPids(const Tpetra::Import< LocalOrdinal, GlobalOrdinal, Node > &Importer, Teuchos::Array< int > &pids, bool use_minus_one_for_local)
Like getPidGidPairs, but just gets the PIDs, ordered by the column Map.
Stand-alone utility functions and macros.
Struct that holds views of the contents of a CrsMatrix.
Common base class of Import and Export.
Teuchos::ArrayView< const LO > getExportLIDs() const
List of entries in the source Map that will be sent to other processes.
size_t getNumExportIDs() const
Number of entries that must be sent by the calling process to other processes.
Teuchos::ArrayView< const int > getExportPIDs() const
List of processes to which entries will be sent.
Impl::CreateMirrorViewFromUnmanagedHostArray< ValueType, OutputDeviceType >::output_view_type create_mirror_view_from_raw_host_array(const OutputDeviceType &, ValueType *inPtr, const size_t inSize, const bool copy=true, const char label[]="")
Variant of Kokkos::create_mirror_view that takes a raw host 1-d array as input.
Namespace Tpetra contains the class and methods constituting the Tpetra library.
void sort2(const IT1 &first1, const IT1 &last1, const IT2 &first2, const bool stableSort=false)
Sort the first array, and apply the resulting permutation to the second array.
void sort3(const IT1 &first1, const IT1 &last1, const IT2 &first2, const IT3 &first3, const bool stableSort=false)
Sort the first array, and apply the same permutation to the second and third arrays.
@ REPLACE
Replace existing values with new values.