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