deal.II version GIT relicensing-2167-g9622207b8f 2024-11-21 12:40:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
vector_data_exchange.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2020 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
16#include <deal.II/base/mpi.h>
17#include <deal.II/base/mpi.templates.h>
20#include <deal.II/base/timer.h>
21
23
24#include <boost/serialization/utility.hpp>
25
26#include <map>
27#include <vector>
28
29
31
32namespace internal
33{
34 namespace MatrixFreeFunctions
35 {
36 namespace VectorDataExchange
37 {
39 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner)
40 : partitioner(partitioner)
41 {}
42
43
44
45 unsigned int
47 {
48 return partitioner->locally_owned_size();
49 }
50
51
52
53 unsigned int
55 {
56 return partitioner->n_ghost_indices();
57 }
58
59
60
61 unsigned int
63 {
64 return partitioner->n_import_indices();
65 }
66
67
68
69 unsigned int
71 {
72 return 0;
73 }
74
75
76
79 {
80 return partitioner->size();
81 }
82
83
84
85 void
87 const unsigned int communication_channel,
88 const ArrayView<const double> &locally_owned_array,
89 const std::vector<ArrayView<const double>> &shared_arrays,
90 const ArrayView<double> &ghost_array,
91 const ArrayView<double> &temporary_storage,
92 std::vector<MPI_Request> &requests) const
93 {
94 (void)shared_arrays;
95#ifndef DEAL_II_WITH_MPI
96 (void)communication_channel;
97 (void)locally_owned_array;
98 (void)ghost_array;
99 (void)temporary_storage;
100 (void)requests;
101#else
102 partitioner->export_to_ghosted_array_start(communication_channel,
103 locally_owned_array,
104 temporary_storage,
105 ghost_array,
106 requests);
107#endif
108 }
109
110
111
112 void
114 const ArrayView<const double> &locally_owned_array,
115 const std::vector<ArrayView<const double>> &shared_arrays,
116 const ArrayView<double> &ghost_array,
117 std::vector<MPI_Request> &requests) const
118 {
119 (void)locally_owned_array;
120 (void)shared_arrays;
121#ifndef DEAL_II_WITH_MPI
122 (void)ghost_array;
123 (void)requests;
124#else
125 partitioner->export_to_ghosted_array_finish(ghost_array, requests);
126#endif
127 }
128
129
130
131 void
133 const VectorOperation::values vector_operation,
134 const unsigned int communication_channel,
135 const ArrayView<const double> &locally_owned_array,
136 const std::vector<ArrayView<const double>> &shared_arrays,
137 const ArrayView<double> &ghost_array,
138 const ArrayView<double> &temporary_storage,
139 std::vector<MPI_Request> &requests) const
140 {
141 (void)locally_owned_array;
142 (void)shared_arrays;
143#ifndef DEAL_II_WITH_MPI
144 (void)vector_operation;
145 (void)communication_channel;
146 (void)ghost_array;
147 (void)temporary_storage;
148 (void)requests;
149#else
150 partitioner->import_from_ghosted_array_start(vector_operation,
151 communication_channel,
152 ghost_array,
153 temporary_storage,
154 requests);
155#endif
156 }
157
158
159
160 void
162 const VectorOperation::values vector_operation,
163 const ArrayView<double> &locally_owned_storage,
164 const std::vector<ArrayView<const double>> &shared_arrays,
165 const ArrayView<double> &ghost_array,
166 const ArrayView<const double> &temporary_storage,
167 std::vector<MPI_Request> &requests) const
168 {
169 (void)shared_arrays;
170#ifndef DEAL_II_WITH_MPI
171 (void)vector_operation;
172 (void)locally_owned_storage;
173 (void)ghost_array;
174 (void)temporary_storage;
175 (void)requests;
176#else
177 partitioner->import_from_ghosted_array_finish(vector_operation,
178 temporary_storage,
179 locally_owned_storage,
180 ghost_array,
181 requests);
182#endif
183 }
184
185
186
187 void
189 const ArrayView<double> &ghost_array) const
190 {
191 reset_ghost_values_impl(ghost_array);
192 }
193
194
195
196 void
198 const unsigned int communication_channel,
199 const ArrayView<const float> &locally_owned_array,
200 const std::vector<ArrayView<const float>> &shared_arrays,
201 const ArrayView<float> &ghost_array,
202 const ArrayView<float> &temporary_storage,
203 std::vector<MPI_Request> &requests) const
204 {
205 (void)shared_arrays;
206#ifndef DEAL_II_WITH_MPI
207 (void)communication_channel;
208 (void)locally_owned_array;
209 (void)ghost_array;
210 (void)temporary_storage;
211 (void)requests;
212#else
213 partitioner->export_to_ghosted_array_start(communication_channel,
214 locally_owned_array,
215 temporary_storage,
216 ghost_array,
217 requests);
218#endif
219 }
220
221
222
223 void
225 const ArrayView<const float> &locally_owned_array,
226 const std::vector<ArrayView<const float>> &shared_arrays,
227 const ArrayView<float> &ghost_array,
228 std::vector<MPI_Request> &requests) const
229 {
230 (void)locally_owned_array;
231 (void)shared_arrays;
232#ifndef DEAL_II_WITH_MPI
233 (void)ghost_array;
234 (void)requests;
235#else
236 partitioner->export_to_ghosted_array_finish(ghost_array, requests);
237#endif
238 }
239
240
241
242 void
244 const VectorOperation::values vector_operation,
245 const unsigned int communication_channel,
246 const ArrayView<const float> &locally_owned_array,
247 const std::vector<ArrayView<const float>> &shared_arrays,
248 const ArrayView<float> &ghost_array,
249 const ArrayView<float> &temporary_storage,
250 std::vector<MPI_Request> &requests) const
251 {
252 (void)locally_owned_array;
253 (void)shared_arrays;
254#ifndef DEAL_II_WITH_MPI
255 (void)vector_operation;
256 (void)communication_channel;
257 (void)ghost_array;
258 (void)temporary_storage;
259 (void)requests;
260#else
261 partitioner->import_from_ghosted_array_start(vector_operation,
262 communication_channel,
263 ghost_array,
264 temporary_storage,
265 requests);
266#endif
267 }
268
269
270
271 void
273 const VectorOperation::values vector_operation,
274 const ArrayView<float> &locally_owned_storage,
275 const std::vector<ArrayView<const float>> &shared_arrays,
276 const ArrayView<float> &ghost_array,
277 const ArrayView<const float> &temporary_storage,
278 std::vector<MPI_Request> &requests) const
279 {
280 (void)shared_arrays;
281#ifndef DEAL_II_WITH_MPI
282 (void)vector_operation;
283 (void)locally_owned_storage;
284 (void)ghost_array;
285 (void)temporary_storage;
286 (void)requests;
287#else
288 partitioner->import_from_ghosted_array_finish(vector_operation,
289 temporary_storage,
290 locally_owned_storage,
291 ghost_array,
292 requests);
293#endif
294 }
295
296
297
298 void
300 const ArrayView<float> &ghost_array) const
301 {
302 reset_ghost_values_impl(ghost_array);
303 }
304
305
306
307 template <typename Number>
308 void
310 const ArrayView<Number> &ghost_array) const
311 {
312 for (const auto &my_ghosts :
313 partitioner->ghost_indices_within_larger_ghost_set())
314 for (unsigned int j = my_ghosts.first; j < my_ghosts.second; ++j)
315 ghost_array[j] = 0.;
316 }
317
318
319
320 namespace internal
321 {
322 std::pair<std::vector<unsigned int>,
323 std::vector<std::pair<unsigned int, unsigned int>>>
325 const std::vector<unsigned int> &sm_export_ptr,
326 const std::vector<unsigned int> &sm_export_indices)
327 {
328 std::vector<unsigned int> recv_ptr = {0};
329 std::vector<unsigned int> recv_indices;
330 std::vector<unsigned int> recv_len;
331
332 for (unsigned int i = 0; i + 1 < sm_export_ptr.size(); ++i)
333 {
334 if (sm_export_ptr[i] != sm_export_ptr[i + 1])
335 {
336 recv_indices.push_back(sm_export_indices[sm_export_ptr[i]]);
337 recv_len.push_back(1);
338
339 for (unsigned int j = sm_export_ptr[i] + 1;
340 j < sm_export_ptr[i + 1];
341 j++)
342 if (recv_indices.back() + recv_len.back() !=
343 sm_export_indices[j])
344 {
345 recv_indices.push_back(sm_export_indices[j]);
346 recv_len.push_back(1);
347 }
348 else
349 recv_len.back()++;
350 }
351 recv_ptr.push_back(recv_indices.size());
352 }
353
354 std::pair<std::vector<unsigned int>,
355 std::vector<std::pair<unsigned int, unsigned int>>>
356 result;
357
358 result.first = recv_ptr;
359
360 for (unsigned int i = 0; i < recv_indices.size(); ++i)
361 result.second.emplace_back(recv_indices[i], recv_len[i]);
362
363 return result;
364 }
365
366 } // namespace internal
367
368
369
371 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner,
372 const MPI_Comm communicator_sm)
373 : comm(partitioner->get_mpi_communicator())
374 , comm_sm(communicator_sm)
375 , n_local_elements(partitioner->locally_owned_range().n_elements())
376 , n_ghost_elements(partitioner->ghost_indices().n_elements())
377 , n_global_elements(partitioner->locally_owned_range().size())
378 {
379#ifndef DEAL_II_WITH_MPI
380 Assert(false, ExcNeedsMPI());
381#else
383 return; // nothing to do in serial case
384
385 const auto &is_locally_owned = partitioner->locally_owned_range();
386 const auto &is_locally_ghost = partitioner->ghost_indices();
387 const auto &ghost_indices_within_larger_ghost_set =
388 partitioner->ghost_indices_within_larger_ghost_set();
389
390 // temporal data structures
391 std::vector<unsigned int> n_ghost_indices_in_larger_set_by_remote_rank;
392
393 std::vector<std::array<unsigned int, 3>> ghost_targets_data;
394
395 std::vector<std::array<unsigned int, 3>> import_targets_data;
396
397 std::vector<unsigned int> sm_ghost_ranks;
398
399 std::vector<unsigned int> sm_import_ranks;
400
401 // temporary uncompressed data structures for ghost_indices_subset_data
402 std::vector<unsigned int> ghost_indices_subset_data_ptr = {0};
403 std::vector<unsigned int> ghost_indices_subset_data_indices;
404
405 // ... for import_indices_data
406 std::vector<unsigned int> import_indices_data_ptr = {0};
407 std::vector<unsigned int> import_indices_data_indices;
408
409 // ... for sm_export_data
410 std::vector<unsigned int> sm_export_data_ptr = {0};
411 std::vector<unsigned int> sm_export_data_indices;
412
413 // ... for sm_export_data_this
414 std::vector<unsigned int> sm_export_data_this_ptr = {0};
415 std::vector<unsigned int> sm_export_data_this_indices;
416
417 // ... for sm_import_data
418 std::vector<unsigned int> sm_import_data_ptr = {};
419 std::vector<unsigned int> sm_import_data_indices;
420
421 // ... for sm_import_data_this
422 std::vector<unsigned int> sm_import_data_this_ptr = {0};
423 std::vector<unsigned int> sm_import_data_this_indices;
424
425 // collect ranks of processes of shared-memory domain
426 const auto sm_ranks =
428
429 // determine owners of ghost indices and determine requesters
430 const auto [owning_ranks_of_ghosts, rank_to_global_indices] =
432 is_locally_ghost,
433 comm);
434
435 // decompress ghost_indices_within_larger_ghost_set for simpler
436 // data access during setup
437 std::vector<unsigned int> shifts_indices;
438 for (const auto &pair : ghost_indices_within_larger_ghost_set)
439 for (unsigned int k = pair.first; k < pair.second; ++k)
440 shifts_indices.push_back(k);
441
442 // process ghost indices
443 {
444 // collect ghost indices according to owning rank
445 std::map<unsigned int, std::vector<types::global_dof_index>>
446 rank_to_local_indices;
447
448 for (unsigned int i = 0; i < owning_ranks_of_ghosts.size(); ++i)
449 rank_to_local_indices[owning_ranks_of_ghosts[i]].push_back(i);
450
451 unsigned int compressed_offset = 0;
452
453 for (const auto &rank_and_local_indices : rank_to_local_indices)
454 {
455 const auto sm_ranks_ptr = std::find(sm_ranks.begin(),
456 sm_ranks.end(),
457 rank_and_local_indices.first);
458
459 if (sm_ranks_ptr == sm_ranks.end()) // remote process
460 {
461 ghost_targets_data.emplace_back(std::array<unsigned int, 3>{{
462 rank_and_local_indices.first, // rank
463 shifts_indices[compressed_offset], // offset
464 static_cast<unsigned int>(
465 rank_and_local_indices.second.size()) // length
466 }});
467
468 for (unsigned int i = 0;
469 i < rank_and_local_indices.second.size();
470 ++i)
471 ghost_indices_subset_data_indices.push_back(
472 shifts_indices[i + compressed_offset]);
473
474 ghost_indices_subset_data_ptr.push_back(
475 ghost_indices_subset_data_indices.size());
476
477 ghost_indices_subset_data.first.push_back(compressed_offset);
478
479 unsigned int i =
481
483 (shifts_indices[ghost_indices_subset_data.first[i] +
484 (ghost_targets_data[i][2] - 1)] -
485 shifts_indices[ghost_indices_subset_data.first[i]]) +
486 1);
487 }
488 else // shared process
489 {
490 sm_ghost_ranks.push_back(
491 std::distance(sm_ranks.begin(), sm_ranks_ptr));
492
493 sm_export_data_ptr.push_back(
494 sm_export_data_ptr.back() +
495 rank_and_local_indices.second.size());
496
497 for (unsigned int i = compressed_offset;
498 i <
499 rank_and_local_indices.second.size() + compressed_offset;
500 ++i)
501 sm_export_data_this_indices.push_back(
502 shifts_indices[i] + is_locally_owned.n_elements());
503
504 sm_export_data_this_ptr.push_back(
505 sm_export_data_this_indices.size());
506 }
507 compressed_offset += rank_and_local_indices.second.size();
508 }
509
510 sm_export_data_indices.resize(sm_export_data_ptr.back());
511 }
512
513 // process requesters
514 {
515 for (const auto &rank_and_global_indices : rank_to_global_indices)
516 {
517 const auto sm_ranks_ptr =
518 std::find(sm_ranks.begin(),
519 sm_ranks.end(),
520 rank_and_global_indices.first);
521
522 if (sm_ranks_ptr == sm_ranks.end()) // remote process
523 {
524 import_targets_data.emplace_back(std::array<unsigned int, 3>{{
525 rank_and_global_indices.first, // rank
526 static_cast<unsigned int>(
527 import_indices_data_indices.size()), // offset
528 static_cast<unsigned int>(
529 rank_and_global_indices.second.n_elements()) // length
530 }});
531
532 for (const auto i : rank_and_global_indices.second)
533 import_indices_data_indices.push_back(
534 is_locally_owned.index_within_set(i));
535
536 import_indices_data_ptr.push_back(
537 import_indices_data_indices.size());
538 }
539 else // shared process
540 {
541 sm_import_ranks.push_back(
542 std::distance(sm_ranks.begin(), sm_ranks_ptr));
543
544 for (const auto i : rank_and_global_indices.second)
545 sm_import_data_this_indices.push_back(
546 is_locally_owned.index_within_set(i));
547
548 sm_import_data_this_ptr.push_back(
549 sm_import_data_this_indices.size());
550 }
551 }
552
553 sm_import_data_ptr = sm_import_data_this_ptr;
554 sm_import_data_indices.resize(sm_import_data_this_ptr.back());
555 }
556
557 // send sm_export_data_this to sm-neighbor -> sm_import_data
558 {
559 std::vector<MPI_Request> requests(sm_ghost_ranks.size() +
560 sm_import_ranks.size());
561
562 for (unsigned int i = 0; i < sm_ghost_ranks.size(); ++i)
563 {
564 const int ierr = MPI_Isend(sm_export_data_this_indices.data() +
565 sm_export_data_this_ptr[i],
566 sm_export_data_this_ptr[i + 1] -
567 sm_export_data_this_ptr[i],
568 MPI_UNSIGNED,
570 4,
571 comm_sm,
572 requests.data() + i);
573 AssertThrowMPI(ierr);
574 }
575
576 for (unsigned int i = 0; i < sm_import_ranks.size(); ++i)
577 {
578 const int ierr =
579 MPI_Irecv(sm_import_data_indices.data() + sm_import_data_ptr[i],
580 sm_import_data_ptr[i + 1] - sm_import_data_ptr[i],
581 MPI_UNSIGNED,
583 4,
584 comm_sm,
585 requests.data() + sm_ghost_ranks.size() + i);
586 AssertThrowMPI(ierr);
587 }
588
589 const int ierr =
590 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
591 AssertThrowMPI(ierr);
592 }
593
594 // send sm_import_data_this to sm-neighbor -> sm_export_data_indices
595 {
596 std::vector<MPI_Request> requests(sm_import_ranks.size() +
597 sm_ghost_ranks.size());
598
599 for (unsigned int i = 0; i < sm_import_ranks.size(); ++i)
600 {
601 const int ierr = MPI_Isend(sm_import_data_this_indices.data() +
602 sm_import_data_this_ptr[i],
603 sm_import_data_this_ptr[i + 1] -
604 sm_import_data_this_ptr[i],
605 MPI_UNSIGNED,
607 2,
608 comm_sm,
609 requests.data() + i);
610 AssertThrowMPI(ierr);
611 }
612
613 for (unsigned int i = 0; i < sm_ghost_ranks.size(); ++i)
614 {
615 const int ierr =
616 MPI_Irecv(sm_export_data_indices.data() + sm_export_data_ptr[i],
617 sm_export_data_ptr[i + 1] - sm_export_data_ptr[i],
618 MPI_UNSIGNED,
620 2,
621 comm_sm,
622 requests.data() + sm_import_ranks.size() + i);
623 AssertThrowMPI(ierr);
624 }
625
626 const int ierr =
627 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
628 AssertThrowMPI(ierr);
629 }
630
631 // store data structures and, if needed, compress them
632 this->n_ghost_indices_in_larger_set_by_remote_rank =
634
637 ghost_indices_subset_data_ptr, ghost_indices_subset_data_indices);
638
639 this->ghost_targets_data = ghost_targets_data;
640
641 this->import_targets_data = import_targets_data;
642
643 this->import_indices_data =
644 internal::compress_to_contiguous_ranges(import_indices_data_ptr,
645 import_indices_data_indices);
646
647 this->sm_ghost_ranks = sm_ghost_ranks;
648
649 this->sm_export_data =
651 sm_export_data_indices);
652
653 this->sm_export_data_this =
654 internal::compress_to_contiguous_ranges(sm_export_data_this_ptr,
655 sm_export_data_this_indices);
656
657 this->sm_import_ranks = sm_import_ranks;
658
659 this->sm_import_data =
661 sm_import_data_indices);
662
663 this->sm_import_data_this =
664 internal::compress_to_contiguous_ranges(sm_import_data_this_ptr,
665 sm_import_data_this_indices);
666
667#endif
668 }
669
670
671
672 void
674 const unsigned int communication_channel,
675 const ArrayView<const double> &locally_owned_array,
676 const std::vector<ArrayView<const double>> &shared_arrays,
677 const ArrayView<double> &ghost_array,
678 const ArrayView<double> &temporary_storage,
679 std::vector<MPI_Request> &requests) const
680 {
681 export_to_ghosted_array_start_impl(communication_channel,
682 locally_owned_array,
683 shared_arrays,
684 ghost_array,
685 temporary_storage,
686 requests);
687 }
688
689
690
691 void
693 const ArrayView<const double> &locally_owned_array,
694 const std::vector<ArrayView<const double>> &shared_arrays,
695 const ArrayView<double> &ghost_array,
696 std::vector<MPI_Request> &requests) const
697 {
698 export_to_ghosted_array_finish_impl(locally_owned_array,
699 shared_arrays,
700 ghost_array,
701 requests);
702 }
703
704
705
706 void
708 const VectorOperation::values vector_operation,
709 const unsigned int communication_channel,
710 const ArrayView<const double> &locally_owned_array,
711 const std::vector<ArrayView<const double>> &shared_arrays,
712 const ArrayView<double> &ghost_array,
713 const ArrayView<double> &temporary_storage,
714 std::vector<MPI_Request> &requests) const
715 {
717 communication_channel,
718 locally_owned_array,
719 shared_arrays,
720 ghost_array,
721 temporary_storage,
722 requests);
723 }
724
725
726
727 void
729 const VectorOperation::values vector_operation,
730 const ArrayView<double> &locally_owned_storage,
731 const std::vector<ArrayView<const double>> &shared_arrays,
732 const ArrayView<double> &ghost_array,
733 const ArrayView<const double> &temporary_storage,
734 std::vector<MPI_Request> &requests) const
735 {
737 locally_owned_storage,
738 shared_arrays,
739 ghost_array,
740 temporary_storage,
741 requests);
742 }
743
744
745
746 void
748 const unsigned int communication_channel,
749 const ArrayView<const float> &locally_owned_array,
750 const std::vector<ArrayView<const float>> &shared_arrays,
751 const ArrayView<float> &ghost_array,
752 const ArrayView<float> &temporary_storage,
753 std::vector<MPI_Request> &requests) const
754 {
755 export_to_ghosted_array_start_impl(communication_channel,
756 locally_owned_array,
757 shared_arrays,
758 ghost_array,
759 temporary_storage,
760 requests);
761 }
762
763
764
765 void
767 const ArrayView<const float> &locally_owned_array,
768 const std::vector<ArrayView<const float>> &shared_arrays,
769 const ArrayView<float> &ghost_array,
770 std::vector<MPI_Request> &requests) const
771 {
772 export_to_ghosted_array_finish_impl(locally_owned_array,
773 shared_arrays,
774 ghost_array,
775 requests);
776 }
777
778
779
780 void
782 const VectorOperation::values vector_operation,
783 const unsigned int communication_channel,
784 const ArrayView<const float> &locally_owned_array,
785 const std::vector<ArrayView<const float>> &shared_arrays,
786 const ArrayView<float> &ghost_array,
787 const ArrayView<float> &temporary_storage,
788 std::vector<MPI_Request> &requests) const
789 {
791 communication_channel,
792 locally_owned_array,
793 shared_arrays,
794 ghost_array,
795 temporary_storage,
796 requests);
797 }
798
799
800
801 void
803 const VectorOperation::values vector_operation,
804 const ArrayView<float> &locally_owned_storage,
805 const std::vector<ArrayView<const float>> &shared_arrays,
806 const ArrayView<float> &ghost_array,
807 const ArrayView<const float> &temporary_storage,
808 std::vector<MPI_Request> &requests) const
809 {
811 locally_owned_storage,
812 shared_arrays,
813 ghost_array,
814 temporary_storage,
815 requests);
816 }
817
818
819
820 template <typename Number>
821 void
823 const unsigned int communication_channel,
824 const ArrayView<const Number> &data_this,
825 const std::vector<ArrayView<const Number>> &data_others,
826 const ArrayView<Number> &buffer,
827 const ArrayView<Number> &temporary_storage,
828 std::vector<MPI_Request> &requests) const
829 {
830#ifndef DEAL_II_WITH_MPI
831 Assert(false, ExcNeedsMPI());
832
833 (void)communication_channel;
834 (void)data_this;
835 (void)data_others;
836 (void)buffer;
837 (void)temporary_storage;
838 (void)requests;
839#else
840 (void)data_others;
841
842 requests.resize(sm_import_ranks.size() + sm_ghost_ranks.size() +
844
845 int dummy;
846 // receive a signal that relevant sm neighbors are ready
847 for (unsigned int i = 0; i < sm_ghost_ranks.size(); ++i)
848 {
849 const int ierr =
850 MPI_Irecv(&dummy,
851 0,
852 MPI_INT,
854 communication_channel + 0,
855 comm_sm,
856 requests.data() + sm_import_ranks.size() + i);
857 AssertThrowMPI(ierr);
858 }
859
860 // signal to all relevant sm neighbors that this process is ready
861 for (unsigned int i = 0; i < sm_import_ranks.size(); ++i)
862 {
863 const int ierr = MPI_Isend(&dummy,
864 0,
865 MPI_INT,
867 communication_channel + 0,
868 comm_sm,
869 requests.data() + i);
870 AssertThrowMPI(ierr);
871 }
872
873 // receive data from remote processes
874 for (unsigned int i = 0; i < ghost_targets_data.size(); ++i)
875 {
876 const unsigned int offset =
878 ghost_targets_data[i][2];
879
880 const int ierr = MPI_Irecv(
881 buffer.data() + ghost_targets_data[i][1] + offset,
882 ghost_targets_data[i][2],
883 Utilities::MPI::mpi_type_id_for_type<decltype(*buffer.data())>,
884 ghost_targets_data[i][0],
885 communication_channel + 1,
886 comm,
887 requests.data() + sm_import_ranks.size() + sm_ghost_ranks.size() +
888 i);
889 AssertThrowMPI(ierr);
890 }
891
892 // send data to remote processes
893 for (unsigned int i = 0, k = 0; i < import_targets_data.size(); ++i)
894 {
895 for (unsigned int j = import_indices_data.first[i];
896 j < import_indices_data.first[i + 1];
897 j++)
898 for (unsigned int l = 0; l < import_indices_data.second[j].second;
899 l++, k++)
900 temporary_storage[k] =
901 data_this[import_indices_data.second[j].first + l];
902
903 // send data away
904 const int ierr = MPI_Isend(
905 temporary_storage.data() + import_targets_data[i][1],
907 Utilities::MPI::mpi_type_id_for_type<decltype(*data_this.data())>,
909 communication_channel + 1,
910 comm,
911 requests.data() + sm_import_ranks.size() + sm_ghost_ranks.size() +
912 ghost_targets_data.size() + i);
913 AssertThrowMPI(ierr);
914 }
915#endif
916 }
917
918
919
920 template <typename Number>
921 void
923 const ArrayView<const Number> &data_this,
924 const std::vector<ArrayView<const Number>> &data_others,
925 const ArrayView<Number> &ghost_array,
926 std::vector<MPI_Request> &requests) const
927 {
928 (void)data_this;
929
930#ifndef DEAL_II_WITH_MPI
931 Assert(false, ExcNeedsMPI());
932
933 (void)data_others;
934 (void)ghost_array;
935 (void)requests;
936#else
937
938 AssertDimension(requests.size(),
939 sm_import_ranks.size() + sm_ghost_ranks.size() +
940 ghost_targets_data.size() +
941 import_targets_data.size());
942
943 const auto split =
944 [&](const unsigned int i) -> std::pair<unsigned int, unsigned int> {
946 (sm_ghost_ranks.size() + ghost_targets_data.size()));
947
948 if (i < sm_ghost_ranks.size())
949 return {0, i};
950 else
951 return {1, i - sm_ghost_ranks.size()};
952 };
953
954 for (unsigned int c = 0;
955 c < sm_ghost_ranks.size() + ghost_targets_data.size();
956 c++)
957 {
958 int i;
959 const int ierr =
960 MPI_Waitany(sm_ghost_ranks.size() + ghost_targets_data.size(),
961 requests.data() + sm_import_ranks.size(),
962 &i,
963 MPI_STATUS_IGNORE);
964 AssertThrowMPI(ierr);
965
966 const auto s = split(i);
967 i = s.second;
968
969 if (s.first == 0)
970 {
971 const Number *DEAL_II_RESTRICT data_others_ptr =
972 data_others[sm_ghost_ranks[i]].data();
973 Number *DEAL_II_RESTRICT data_this_ptr = ghost_array.data();
974
975 for (unsigned int lo = sm_export_data.first[i],
976 ko = sm_export_data_this.first[i],
977 li = 0,
978 ki = 0;
979 (lo < sm_export_data.first[i + 1]) &&
980 (ko < sm_export_data_this.first[i + 1]);)
981 {
982 for (; (li < sm_export_data.second[lo].second) &&
983 (ki < sm_export_data_this.second[ko].second);
984 ++li, ++ki)
985 data_this_ptr[sm_export_data_this.second[ko].first + ki -
987 data_others_ptr[sm_export_data.second[lo].first + li];
988
989 if (li == sm_export_data.second[lo].second)
990 {
991 lo++; // increment outer counter
992 li = 0; // reset inner counter
993 }
994
995 if (ki == sm_export_data_this.second[ko].second)
996 {
997 ko++; // increment outer counter
998 ki = 0; // reset inner counter
999 }
1000 }
1001 }
1002 else /*if(s.second == 1)*/
1003 {
1004 const unsigned int offset =
1006 ghost_targets_data[i][2];
1007
1008 for (unsigned int c = 0,
1009 ko = ghost_indices_subset_data.first[i],
1010 ki = 0;
1011 c < ghost_targets_data[i][2];
1012 ++c)
1013 {
1015 ghost_indices_subset_data.second.size());
1016
1017 const unsigned int idx_1 =
1018 ghost_indices_subset_data.second[ko].first + ki;
1019 const unsigned int idx_2 =
1020 ghost_targets_data[i][1] + c + offset;
1021
1022 AssertIndexRange(idx_1, ghost_array.size());
1023 AssertIndexRange(idx_2, ghost_array.size());
1024
1025 if (idx_1 == idx_2)
1026 {
1027 // noting to do
1028 }
1029 else if (idx_1 < idx_2)
1030 {
1031 ghost_array[idx_1] = ghost_array[idx_2];
1032 ghost_array[idx_2] = 0.0;
1033 }
1034 else
1035 {
1037 }
1038
1039 ++ki;
1040
1041 if (ki == ghost_indices_subset_data.second[ko].second)
1042 {
1043 ko++; // increment outer counter
1044 ki = 0; // reset inner counter
1045 }
1046 }
1047 }
1048 }
1049
1050 const int ierr =
1051 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
1052 AssertThrowMPI(ierr);
1053
1054#endif
1055 }
1056
1057
1058
1059 template <typename Number>
1060 void
1062 const VectorOperation::values operation,
1063 const unsigned int communication_channel,
1064 const ArrayView<const Number> &data_this,
1065 const std::vector<ArrayView<const Number>> &data_others,
1066 const ArrayView<Number> &buffer,
1067 const ArrayView<Number> &temporary_storage,
1068 std::vector<MPI_Request> &requests) const
1069 {
1070 (void)data_this;
1071
1072#ifndef DEAL_II_WITH_MPI
1073 Assert(false, ExcNeedsMPI());
1074
1075 (void)operation;
1076 (void)communication_channel;
1077 (void)data_others;
1078 (void)buffer;
1079 (void)temporary_storage;
1080 (void)requests;
1081#else
1082 // return;
1083
1084 (void)data_others;
1085 (void)operation;
1086
1088
1089 requests.resize(sm_ghost_ranks.size() + sm_import_ranks.size() +
1090 ghost_targets_data.size() + import_targets_data.size());
1091
1092 int dummy;
1093 for (unsigned int i = 0; i < sm_ghost_ranks.size(); ++i)
1094 {
1095 const int ierr = MPI_Isend(&dummy,
1096 0,
1097 MPI_INT,
1098 sm_ghost_ranks[i],
1099 communication_channel + 1,
1100 comm_sm,
1101 requests.data() + i);
1102 AssertThrowMPI(ierr);
1103 }
1104
1105 for (unsigned int i = 0; i < sm_import_ranks.size(); ++i)
1106 {
1107 const int ierr =
1108 MPI_Irecv(&dummy,
1109 0,
1110 MPI_INT,
1111 sm_import_ranks[i],
1112 communication_channel + 1,
1113 comm_sm,
1114 requests.data() + sm_ghost_ranks.size() + i);
1115 AssertThrowMPI(ierr);
1116 }
1117
1118 for (unsigned int i = 0; i < ghost_targets_data.size(); ++i)
1119 {
1120 for (unsigned int c = 0,
1121 ko = ghost_indices_subset_data.first[i],
1122 ki = 0;
1123 c < ghost_targets_data[i][2];
1124 ++c)
1125 {
1127
1128 const unsigned int idx_1 =
1129 ghost_indices_subset_data.second[ko].first + ki;
1130 const unsigned int idx_2 = ghost_targets_data[i][1] + c;
1131
1132 AssertIndexRange(idx_1, buffer.size());
1133 AssertIndexRange(idx_2, buffer.size());
1134
1135 if (idx_1 == idx_2)
1136 {
1137 // nothing to do
1138 }
1139 else if (idx_2 < idx_1)
1140 {
1141 buffer[idx_2] = buffer[idx_1];
1142 buffer[idx_1] = 0.0;
1143 }
1144 else
1145 {
1147 }
1148
1149 if (++ki == ghost_indices_subset_data.second[ko].second)
1150 {
1151 ko++; // increment outer counter
1152 ki = 0; // reset inner counter
1153 }
1154 }
1155
1156 const int ierr = MPI_Isend(
1157 buffer.data() + ghost_targets_data[i][1],
1158 ghost_targets_data[i][2],
1159 Utilities::MPI::mpi_type_id_for_type<decltype(*buffer.data())>,
1160 ghost_targets_data[i][0],
1161 communication_channel + 0,
1162 comm,
1163 requests.data() + sm_ghost_ranks.size() + sm_import_ranks.size() +
1164 i);
1165 AssertThrowMPI(ierr);
1166 }
1167
1168 for (unsigned int i = 0; i < import_targets_data.size(); ++i)
1169 {
1170 const int ierr =
1171 MPI_Irecv(temporary_storage.data() + import_targets_data[i][1],
1172 import_targets_data[i][2],
1174 decltype(*temporary_storage.data())>,
1175 import_targets_data[i][0],
1176 communication_channel + 0,
1177 comm,
1178 requests.data() + sm_ghost_ranks.size() +
1179 sm_import_ranks.size() + ghost_targets_data.size() +
1180 i);
1181 AssertThrowMPI(ierr);
1182 }
1183#endif
1184 }
1185
1186
1187
1188 template <typename Number>
1189 void
1191 const VectorOperation::values operation,
1192 const ArrayView<Number> &data_this,
1193 const std::vector<ArrayView<const Number>> &data_others,
1194 const ArrayView<Number> &buffer,
1195 const ArrayView<const Number> &temporary_storage,
1196 std::vector<MPI_Request> &requests) const
1197 {
1198#ifndef DEAL_II_WITH_MPI
1199 Assert(false, ExcNeedsMPI());
1200
1201 (void)operation;
1202 (void)data_this;
1203 (void)data_others;
1204 (void)buffer;
1205 (void)temporary_storage;
1206 (void)requests;
1207#else
1208
1209 (void)operation;
1210
1212
1213 AssertDimension(requests.size(),
1214 sm_ghost_ranks.size() + sm_import_ranks.size() +
1215 ghost_targets_data.size() +
1216 import_targets_data.size());
1217
1218 const auto split =
1219 [&](const unsigned int i) -> std::pair<unsigned int, unsigned int> {
1221 (sm_import_ranks.size() + ghost_targets_data.size() +
1222 import_targets_data.size()));
1223
1224 if (i < sm_import_ranks.size())
1225 return {0, i};
1226 else if (i < (sm_import_ranks.size() + ghost_targets_data.size()))
1227 return {2, i - sm_import_ranks.size()};
1228 else
1229 return {1, i - sm_import_ranks.size() - ghost_targets_data.size()};
1230 };
1231
1232 for (unsigned int c = 0;
1233 c < sm_import_ranks.size() + import_targets_data.size() +
1234 ghost_targets_data.size();
1235 c++)
1236 {
1237 int i;
1238 const int ierr =
1239 MPI_Waitany(sm_import_ranks.size() + import_targets_data.size() +
1240 ghost_targets_data.size(),
1241 requests.data() + sm_ghost_ranks.size(),
1242 &i,
1243 MPI_STATUS_IGNORE);
1244 AssertThrowMPI(ierr);
1245
1246 const auto &s = split(i);
1247 i = s.second;
1248
1249 if (s.first == 0)
1250 {
1251 Number *DEAL_II_RESTRICT data_others_ptr =
1252 const_cast<Number *>(data_others[sm_import_ranks[i]].data());
1253 Number *DEAL_II_RESTRICT data_this_ptr = data_this.data();
1254
1255 for (unsigned int lo = sm_import_data_this.first[i],
1256 ko = sm_import_data.first[i],
1257 li = 0,
1258 ki = 0;
1259 (lo < sm_import_data_this.first[i + 1]) &&
1260 (ko < sm_import_data.first[i + 1]);)
1261 {
1262 for (; (li < sm_import_data_this.second[lo].second) &&
1263 (ki < sm_import_data.second[ko].second);
1264 ++li, ++ki)
1265 {
1266 data_this_ptr[sm_import_data_this.second[lo].first +
1267 li] +=
1268 data_others_ptr[sm_import_data.second[ko].first + ki];
1269 data_others_ptr[sm_import_data.second[ko].first + ki] =
1270 0.0;
1271 }
1272
1273 if (li == sm_import_data_this.second[lo].second)
1274 {
1275 lo++; // increment outer counter
1276 li = 0; // reset inner counter
1277 }
1278 if (ki == sm_import_data.second[ko].second)
1279 {
1280 ko++; // increment outer counter
1281 ki = 0; // reset inner counter
1282 }
1283 }
1284 }
1285 else if (s.first == 1)
1286 {
1287 for (unsigned int j = import_indices_data.first[i],
1288 k = import_targets_data[i][1];
1289 j < import_indices_data.first[i + 1];
1290 j++)
1291 for (unsigned int l = 0;
1292 l < import_indices_data.second[j].second;
1293 l++)
1294 data_this[import_indices_data.second[j].first + l] +=
1295 temporary_storage[k++];
1296 }
1297 else /*if (s.first == 2)*/
1298 {
1299 std::memset(buffer.data() + ghost_targets_data[i][1],
1300 0.0,
1301 (ghost_targets_data[i][2]) * sizeof(Number));
1302 }
1303 }
1304
1305 const int ierr =
1306 MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
1307 AssertThrowMPI(ierr);
1308#endif
1309 }
1310
1311
1312
1313 unsigned int
1315 {
1316 return n_local_elements;
1317 }
1318
1319
1320
1321 unsigned int
1323 {
1324 return n_ghost_elements;
1325 }
1326
1327
1328
1329 unsigned int
1331 {
1332 if (import_targets_data.empty())
1333 return 0;
1334 return import_targets_data.back()[1] + import_targets_data.back()[2];
1335 }
1336
1337
1338
1339 unsigned int
1341 {
1342 return sm_import_ranks.size() + sm_ghost_ranks.size(); // TODO
1343 }
1344
1345
1346
1349 {
1350 return n_global_elements;
1351 }
1352
1353
1354
1355 MPI_Comm
1357 {
1358 return this->comm_sm;
1359 }
1360
1361
1362
1363 void
1365 {
1366 reset_ghost_values_impl(ghost_array);
1367 }
1368
1369
1370
1371 void
1373 {
1374 reset_ghost_values_impl(ghost_array);
1375 }
1376
1377
1378
1379 template <typename Number>
1380 void
1382 {
1383 // reset ghost values coming from shared-memory neighbors
1384 // TODO: only needed if values are buffered
1385 for (const auto &i : sm_export_data_this.second)
1386 std::memset(ghost_array.data() + (i.first - n_local_elements),
1387 0,
1388 sizeof(Number) * i.second);
1389
1390 // reset ghost values coming from remote neighbors
1391 for (const auto &i : ghost_indices_subset_data.second)
1392 std::memset(ghost_array.data() + i.first,
1393 0,
1394 sizeof(Number) * i.second);
1395 }
1396
1397
1398
1399 } // namespace VectorDataExchange
1400 } // namespace MatrixFreeFunctions
1401} // namespace internal
1402
1403
value_type * data() const noexcept
Definition array_view.h:661
std::size_t size() const
Definition array_view.h:684
void export_to_ghosted_array_start(const unsigned int communication_channel, const ArrayView< const double > &locally_owned_array, const std::vector< ArrayView< const double > > &shared_arrays, const ArrayView< double > &ghost_array, const ArrayView< double > &temporary_storage, std::vector< MPI_Request > &requests) const override
void reset_ghost_values(const ArrayView< double > &ghost_array) const override
std::vector< std::array< unsigned int, 3 > > ghost_targets_data
std::pair< std::vector< unsigned int >, std::vector< std::pair< unsigned int, unsigned int > > > sm_export_data_this
void export_to_ghosted_array_finish(const ArrayView< const double > &locally_owned_array, const std::vector< ArrayView< const double > > &shared_arrays, const ArrayView< double > &ghost_array, std::vector< MPI_Request > &requests) const override
void import_from_ghosted_array_start(const VectorOperation::values vector_operation, const unsigned int communication_channel, const ArrayView< const double > &locally_owned_array, const std::vector< ArrayView< const double > > &shared_arrays, const ArrayView< double > &ghost_array, const ArrayView< double > &temporary_storage, std::vector< MPI_Request > &requests) const override
void export_to_ghosted_array_start_impl(const unsigned int communication_channel, const ArrayView< const Number > &locally_owned_array, const std::vector< ArrayView< const Number > > &shared_arrays, const ArrayView< Number > &ghost_array, const ArrayView< Number > &temporary_storage, std::vector< MPI_Request > &requests) const
void import_from_ghosted_array_finish_impl(const VectorOperation::values vector_operation, const ArrayView< Number > &locally_owned_storage, const std::vector< ArrayView< const Number > > &shared_arrays, const ArrayView< Number > &ghost_array, const ArrayView< const Number > &temporary_storage, std::vector< MPI_Request > &requests) const
std::pair< std::vector< unsigned int >, std::vector< std::pair< unsigned int, unsigned int > > > sm_export_data
std::vector< std::array< unsigned int, 3 > > import_targets_data
std::pair< std::vector< unsigned int >, std::vector< std::pair< unsigned int, unsigned int > > > ghost_indices_subset_data
void import_from_ghosted_array_finish(const VectorOperation::values vector_operation, const ArrayView< double > &locally_owned_storage, const std::vector< ArrayView< const double > > &shared_arrays, const ArrayView< double > &ghost_array, const ArrayView< const double > &temporary_storage, std::vector< MPI_Request > &requests) const override
void export_to_ghosted_array_finish_impl(const ArrayView< const Number > &locally_owned_array, const std::vector< ArrayView< const Number > > &shared_arrays, const ArrayView< Number > &ghost_array, std::vector< MPI_Request > &requests) const
std::pair< std::vector< unsigned int >, std::vector< std::pair< unsigned int, unsigned int > > > sm_import_data
void import_from_ghosted_array_start_impl(const VectorOperation::values vector_operation, const unsigned int communication_channel, const ArrayView< const Number > &locally_owned_array, const std::vector< ArrayView< const Number > > &shared_arrays, const ArrayView< Number > &ghost_array, const ArrayView< Number > &temporary_storage, std::vector< MPI_Request > &requests) const
Full(const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner, const MPI_Comm communicator_sm)
void reset_ghost_values_impl(const ArrayView< Number > &ghost_array) const
std::pair< std::vector< unsigned int >, std::vector< std::pair< unsigned int, unsigned int > > > sm_import_data_this
virtual types::global_dof_index size() const override
std::pair< std::vector< unsigned int >, std::vector< std::pair< unsigned int, unsigned int > > > import_indices_data
const std::shared_ptr< const Utilities::MPI::Partitioner > partitioner
void export_to_ghosted_array_finish(const ArrayView< const double > &locally_owned_array, const std::vector< ArrayView< const double > > &shared_arrays, const ArrayView< double > &ghost_array, std::vector< MPI_Request > &requests) const override
PartitionerWrapper(const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner)
void export_to_ghosted_array_start(const unsigned int communication_channel, const ArrayView< const double > &locally_owned_array, const std::vector< ArrayView< const double > > &shared_arrays, const ArrayView< double > &ghost_array, const ArrayView< double > &temporary_storage, std::vector< MPI_Request > &requests) const override
void reset_ghost_values_impl(const ArrayView< Number > &ghost_array) const
void import_from_ghosted_array_finish(const VectorOperation::values vector_operation, const ArrayView< double > &locally_owned_storage, const std::vector< ArrayView< const double > > &shared_arrays, const ArrayView< double > &ghost_array, const ArrayView< const double > &temporary_storage, std::vector< MPI_Request > &requests) const override
void reset_ghost_values(const ArrayView< double > &ghost_array) const override
void import_from_ghosted_array_start(const VectorOperation::values vector_operation, const unsigned int communication_channel, const ArrayView< const double > &locally_owned_array, const std::vector< ArrayView< const double > > &shared_arrays, const ArrayView< double > &ghost_array, const ArrayView< double > &temporary_storage, std::vector< MPI_Request > &requests) const override
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_RESTRICT
Definition config.h:110
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
static ::ExceptionBase & ExcNotImplemented()
static ::ExceptionBase & ExcNeedsMPI()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
#define DEAL_II_NOT_IMPLEMENTED()
std::size_t size
Definition mpi.cc:734
const MPI_Comm comm
Definition mpi.cc:913
std::pair< std::vector< unsigned int >, std::map< unsigned int, IndexSet > > compute_index_owner_and_requesters(const IndexSet &owned_indices, const IndexSet &indices_to_look_up, const MPI_Comm &comm)
Definition mpi.cc:1860
bool job_supports_mpi()
Definition mpi.cc:681
const MPI_Datatype mpi_type_id_for_type
Definition mpi.h:1626
std::vector< unsigned int > mpi_processes_within_communicator(const MPI_Comm comm_large, const MPI_Comm comm_small)
Definition mpi.cc:123
std::pair< std::vector< unsigned int >, std::vector< std::pair< unsigned int, unsigned int > > > compress_to_contiguous_ranges(const std::vector< unsigned int > &sm_export_ptr, const std::vector< unsigned int > &sm_export_indices)