Loading...
Searching...
No Matches
UPstream.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2015-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Class
28 Foam::UPstream
29
30Description
31 Inter-processor communications stream
32
33SourceFiles
34 UPstream.C
35 UPstream.txx
36 UPstreamCommsStruct.C
37 UPstreamReduceOffsets.H
38
39\*---------------------------------------------------------------------------*/
40
41#ifndef Foam_UPstream_H
42#define Foam_UPstream_H
43
44#include "wordList.H"
45#include "labelList.H"
46#include "DynamicList.H"
47#include "HashTable.H"
48#include "Map.H"
49#include "Enum.H"
50#include "ListOps.H"
51#include "OffsetRange.H"
52#include "OffsetRangeIO.H"
53
54// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55
56namespace Foam
58
59//- Implementation details for UPstream/Pstream/MPI etc.
60namespace PstreamDetail {}
61
62//- Interface handling for UPstream/Pstream/MPI etc.
63namespace PstreamUtils {}
64
65/*---------------------------------------------------------------------------*\
66 Class UPstream Declaration
67\*---------------------------------------------------------------------------*/
69class UPstream
70{
71public:
72
73 //- Int ranges are used for MPI ranks (processes)
76 //- Communications types
77 enum class commsTypes : char
78 {
79 buffered,
82 // Aliases
84 };
85
86 //- Enumerated names for the communication types
88
89 //- Different MPI-send modes (ignored for commsTypes::buffered)
90 enum class sendModes : char
91 {
93 sync
94 };
95
96 //- Mapping of some fundamental and aggregate types to MPI data types
97 enum class dataTypes : char
98 {
99 // NOTE: changes here require adjustment in
100 // PstreamGlobals, UPstreamTraits
101
102 // Fundamental Types [10]:
119 // User Types [6]:
126
127 // Internal markers [1]
128 invalid,
132 };
134 //- Mapping of some MPI op codes.
135 // Currently excluding min/max location until they are needed
136 enum class opCodes : char
138 // NOTE: changes here require adjustment in
139 // PstreamGlobals, UPstreamTraits
140
142
143 // Reduce or window operations [10]
145 op_max,
146 op_sum,
147 op_prod,
149 op_bool_or,
151 op_bit_and,
152 op_bit_or,
154
159 // Window-only operations [2]
163 // Internal markers [1]
168 };
169
170 //- Some bit masks corresponding to topology controls
171 // These selectively enable topology-aware handling
172 enum class topoControls : int
174 broadcast = 1,
175 reduce = 2,
176 gather = 4,
177 combine = 16,
179 gatherList = 64,
180 };
181
182
183 // Public Classes
184
185 //- Wrapper for OpenFOAM internal communicator index
186 class communicator; // Forward Declaration
188 //- Wrapper for MPI_Comm
189 class Communicator; // Forward Declaration
191 //- Wrapper for MPI_File routines
192 class File; // Forward Declaration
194 //- Wrapper for MPI_Request
195 class Request; // Forward Declaration
196
197 //- Wrapper for MPI_Win
198 class Window; // Forward Declaration
199
200 //- Structure for communicating between processors
201 class commsStruct
202 {
203 // Private Data
204
205 //- The procID of the processor \em directly above
206 int above_;
207
208 //- The procIDs of processors \em directly below
209 List<int> below_;
210
211 //- The procIDs of all processors below myProcNo,
212 //- not just directly below
213 List<int> allBelow_;
214
215 //- The procIDs of all processors not below myProcNo
216 //- (inverse of allBelow_ without myProcNo)
217 List<int> allNotBelow_;
218
219
220 public:
221
222 // Constructors
223
224 //- Default construct with above == -1
225 commsStruct() noexcept : above_(-1) {}
226
227 //- Move construct from components
229 (
230 const int above,
234 );
235
236 //- Copy construct from below, allBelow components
238 (
239 const int numProcs,
240 const int myProcID,
241 const int above,
242 const UList<int>& below,
243 const UList<int>& allBelow
244 );
245
246
247 // Member Functions
248
249 // Access
250
251 //- The procID of the processor \em directly above
252 int above() const noexcept { return above_; }
253
254 //- The procIDs of the processors \em directly below
255 const List<int>& below() const noexcept { return below_; }
256
257 //- The procIDs of \em all processors below
258 //- (so not just directly below)
259 const List<int>& allBelow() const noexcept { return allBelow_; }
260
261 //- The procIDs of all processors not below myProcNo.
262 //- The inverse set of allBelow without myProcNo.
263 const List<int>& allNotBelow() const noexcept
264 {
265 return allNotBelow_;
266 }
267
268 //- The number of processors addressed by the structure
269 int nProcs() const noexcept;
270
271
272 // Edit
273
274 //- Reset to default constructed state
275 void reset();
276
277 //- Reset to linear (flat) communication
278 void reset_linear
279 (
280 const int myProci,
281 const int numProcs
282 );
283
284 //- Reset to tree selection (linear when very small),
285 //- possibly with communicator-specific adjustments
286 void reset
287 (
288 const int myProci,
289 const int numProcs,
290 const int communicator
291 );
292
293
294 // Member / Friend Operators
296 bool operator==(const commsStruct&) const;
297 bool operator!=(const commsStruct&) const;
298
299 friend Ostream& operator<<(Ostream&, const commsStruct&);
300 };
301
302 //- Collection of communication structures
303 class commsStructList
304 {
305 // Private Data
307 //- The communicator index
308 int comm_;
309
310 //- Linear (flat) communication instead of tree communication
311 bool flat_;
313 //- The communication tree
314 List<commsStruct> tree_;
315
316 public:
317
318 // Constructors
319
320 //- Construct empty with invalid communicator
322 :
323 comm_(-1),
324 flat_(false)
325 {}
326
327 //- Construct empty with given communicator
328 explicit commsStructList(int comm, bool flat=false) noexcept
329 :
330 comm_(comm),
331 flat_(flat)
332 {}
333
334
335 // Static Functions
336
337 //- An empty structure. Used for placeholders etc.
338 static const commsStructList& null();
339
340
341 // Member Functions
342
343 //- True if communicator is non-negative (ie, was assigned)
344 bool good() const noexcept { return (comm_ >= 0); }
345
346 //- The communicator internal index
347 int comm() const noexcept { return comm_; }
348
349 //- Linear (flat) communication instead of tree communication
350 bool linear() const noexcept { return flat_; }
351
352 //- Change communication type (linear vs tree communication).
353 //- Resets existing structure.
354 void linear(bool on);
355
356 //- Clear the list
357 void clear() { return tree_.clear(); }
358
359 //- True if the list is empty
360 bool empty() const noexcept { return tree_.empty(); }
361
362 //- The number of entries
363 label size() const noexcept { return tree_.size(); }
364
365 //- Reset communicator index, fill tree with empty entries
366 void init(int communicator);
367
368 //- Reset communicator index, fill tree with empty entries.
369 //- Specify flat vs tree communication.
370 void init(int communicator, bool flat);
371
372 //- Reset communicator index, clear tree entries
373 void reset(int communicator);
374
375 //- Reset communicator index, clear tree entries.
376 //- Specify flat vs tree communication.
377 void reset(int communicator, bool flat);
378
379 //- Get existing or create (demand-driven) entry
380 const UPstream::commsStruct& get(int proci) const;
381
382 //- Get existing or create (demand-driven) entry
383 const UPstream::commsStruct& operator[](int proci) const
384 {
385 return get(proci);
386 }
387
388 //- Print un-directed graph in graphviz dot format
389 void printGraph(Ostream& os, int proci = 0) const;
390
391 //- Write the communication tree
392 Ostream& writeList(Ostream& os) const
393 {
394 return tree_.writeList(os);
395 }
396
398 // Member / Friend Operators
399
400 //- Write the communication tree
401 friend Ostream& operator<<
402 (
403 Ostream& os,
404 const commsStructList& rhs
405 )
406 {
407 return rhs.writeList(os);
408 }
409 };
410
411
412private:
413
414 // Private Data
415
416 //- Communications type of this stream
417 commsTypes commsType_;
418
419
420 // Private Static Data
421
422 //- By default this is not a parallel run
423 static bool parRun_;
424
425 //- Have support for threads?
426 static bool haveThreads_;
428 //- Initial MPI_Comm_dup(MPI_COMM_WORLD,...) disabled? (default: false)
429 static bool noInitialCommDup_;
430
431 //- Standard transfer message type
432 static int msgType_;
433
434 //- Index to the world-communicator as defined at startup
435 //- (after any multi-world definitions).
436 //- Is unaffected by any later changes to worldComm.
437 static int constWorldComm_;
439 //- Index to the inter-node communicator (between nodes),
440 //- defined based on constWorldComm_
441 static int commInterNode_;
442
443 //- Index to the intra-host communicator (within a node),
444 //- defined based on constWorldComm_
445 static int commLocalNode_;
446
447 //- The number of shared/host nodes in the (const) world communicator.
448 static int numNodes_;
449
450 //- Names of all worlds
451 static wordList allWorlds_;
452
453 //- Per processor the world index (into allWorlds_)
454 static labelList worldIDs_;
455
456
457 // Communicator specific data
458
459 //- My processor number
460 static DynamicList<int> myProcNo_;
461
462 //- List of process IDs
463 static DynamicList<List<int>> procIDs_;
464
465 //- Parent communicator
466 static DynamicList<label> parentComm_;
467
468 //- Free communicators
469 static DynamicList<label> freeComms_;
470
471 //- Linear communication schedule
472 static DynamicList<commsStructList> linearCommunication_;
473
474 //- Multi level communication schedule
475 static DynamicList<commsStructList> treeCommunication_;
476
477
478 // Private Member Functions
479
480 //- Set data for parallel running
481 static void setParRun(const int nProcs, const bool haveThreads);
482
483 //- Initialise entries for new communicator.
484 //
485 // Resets corresponding entry in myProcNo_, procIDs_,
486 // linearCommunication_, treeCommunication_
487 // \return the communicator index
488 static label getAvailableCommIndex(const label parentIndex);
489
490 //- Define inter-host/intra-host communicators (uses commConstWorld).
491 // Optionally specify a given number per node.
492 static bool setHostCommunicators(const int numPerNode = 0);
494 //- Define inter-host/intra-host communicators based on
495 //- shared-memory information. Uses comm-world.
496 static bool setSharedMemoryCommunicators();
497
498 //- Allocate MPI components of communicator with given index.
499 // This represents a "top-down" approach, creating a communicator
500 // based on the procIDs_ groupings.
501 //
502 // Modifies myProcNo_, reads and modifies procIDs_
503 static void allocateCommunicatorComponents
505 const label parentIndex,
506 const label index
507 );
508
509 //- Allocate MPI components as duplicate of the parent communicator
510 //
511 // Modifies myProcNo_, procIDs_
512 static void dupCommunicatorComponents
513 (
514 const label parentIndex,
515 const label index
516 );
517
518 //- Allocate MPI components for the given index by splitting
519 //- the parent communicator on the given \em colour.
520 // This represents a "bottom-up" approach, when the individual ranks
521 // only know which group they should belong to, but don't yet know
522 // which other ranks will be in their group.
523 //
524 // Modifies myProcNo_, procIDs_
525 static void splitCommunicatorComponents
526 (
527 const label parentIndex,
528 const label index,
529 const int colour,
531 const bool two_step = true
532 );
533
534 //- Free MPI components of communicator.
535 // Does not touch the first two communicators (SELF, WORLD)
536 static void freeCommunicatorComponents(const label index);
537
538
539 // Private implementation helpers
540
541 //- Test for communicator equality.
542 // True if they have the same index or address the same ranks
543 static bool sameProcs_impl(int comm1, int comm2)
544 {
545 return
546 (
547 (comm1 == comm2) ||
548 (
549 // With guard against bad index
550 (comm1 >= 0 && comm1 <= procIDs_.size())
551 && (comm2 >= 0 && comm2 <= procIDs_.size())
552 && (procIDs_[comm1] == procIDs_[comm2])
553 )
554 );
555 }
556
557 //- Test equality of communicator procs with the given list of ranks
558 template<typename IntType>
559 static bool sameProcs_impl
560 (
561 int communicator,
562 const UList<IntType>& procs
563 )
564 {
565 return
566 (
567 // With guard against bad index
568 (communicator >= 0 && communicator <= procIDs_.size())
569 && ListOps::equal(procIDs_[communicator], procs)
570 );
571 }
572
573 //- Test the equality of two lists of ranks
574 template<typename Type1, typename Type2>
575 static bool sameProcs_impl
576 (
577 const UList<Type1>& procs1,
578 const UList<Type2>& procs2
579 )
580 {
581 return ListOps::equal(procs1, procs2);
582 }
583
584
585protected:
586
587 // Protected Member Functions
588
589 // Static Functions
590
591 //- Broadcast buffer contents to all ranks (default: from rank=0).
592 //- The sizes must match on all processes.
593 // For \b non-parallel : do nothing.
594 //
595 // \note The method uses a \c void pointer and the required data type
596 // (as per MPI). This means it should almost never be called directly
597 // but always via a compile-time checked caller.
598 // \return True on success
599 static bool mpi_broadcast
600 (
601 void* buf,
602 std::streamsize count,
603 const UPstream::dataTypes dataTypeId,
604 const int communicator,
605 const int root = 0
606 );
607
608 //- In-place reduction of \c values with result on rank 0.
609 // Includes internal parallel guard and checks on data types, opcode.
610 //
611 // \note The method uses a \c void pointer and the required data type
612 // (as per MPI). This means it should almost never be called directly
613 // but always via a compile-time checked caller.
614 static void mpi_reduce
615 (
616 void* values,
617 int count,
618 const UPstream::dataTypes dataTypeId,
619 const UPstream::opCodes opCodeId,
620 const int communicator,
622 UPstream::Request* req = nullptr
623 );
624
625 //- In-place reduction of \c values with same result on all ranks.
626 // Includes internal parallel guard and checks on data types, opcode.
627 //
628 // \note The method uses a \c void pointer and the required data type
629 // (as per MPI). This means it should almost never be called directly
630 // but always via a compile-time checked caller.
631 static void mpi_allreduce
632 (
633 void* values,
634 int count,
635 const UPstream::dataTypes dataTypeId,
636 const UPstream::opCodes opCodeId,
637 const int communicator,
639 UPstream::Request* req = nullptr
640 );
641
642 //- In-place scan/exscan reduction of \c values.
643 // Includes internal parallel guard and checks on data types, opcode.
644 //
645 // \note The method uses a \c void pointer and the required data type
646 // (as per MPI). This means it should almost never be called directly
647 // but always via a compile-time checked caller.
648 static void mpi_scan_reduce
649 (
650 void* values,
651 int count,
652 const UPstream::dataTypes dataTypeId,
653 const UPstream::opCodes opCodeId,
654 const int communicator,
656 const bool exclusive
657 );
658
659 //- Send buffer contents of specified data type to given processor.
660 //
661 // \note The method uses a \c void pointer and the required data type
662 // (as per MPI). This means it should almost never be called directly
663 // but always via a compile-time checked caller.
664 // \return True on success (or Fatal)
665 static bool mpi_send
666 (
668 const void* buf,
669 std::streamsize count,
670 const UPstream::dataTypes dataTypeId,
671 const int toProcNo,
672 const int tag, // eg, UPstream::msgType()
673 const int communicator, // eg, UPstream::worldComm
675 UPstream::Request* req = nullptr,
678 );
679
680 //- Receive buffer contents of specified data type from given processor.
681 //
682 // \note The method uses a \c void pointer and the required data type
683 // (as per MPI). This means it should almost never be called directly
684 // but always via a compile-time checked caller.
685 // The commsType will be ignored if UPstream::Request is specified.
686 //
687 // \return number of elements read. May change in the future
688 static std::streamsize mpi_receive
689 (
691 void* buf,
692 std::streamsize count,
693 const UPstream::dataTypes dataTypeId,
694 const int fromProcNo,
695 const int tag, // eg, UPstream::msgType()
696 const int communicator, // eg, UPstream::worldComm
698 UPstream::Request* req = nullptr
699 );
700
701 //- Receive identically-sized (contiguous) data from all ranks,
702 //- placing the result on rank 0.
703 // Includes internal parallel guard.
704 // For non-parallel, does not copy any data.
705 // If needed, this must be done by the caller.
706 static void mpi_gather
707 (
709 const void* sendData,
712 void* recvData,
714 int count,
715 const UPstream::dataTypes dataTypeId,
716 const int communicator,
718 UPstream::Request* req = nullptr
719 );
720
721 //- Send identically-sized (contiguous) data from rank 0
722 //- to all other ranks.
723 // Includes internal parallel guard.
724 static void mpi_scatter
725 (
728 const void* sendData,
730 void* recvData,
732 int count,
733 const UPstream::dataTypes dataTypeId,
734 const int communicator,
736 UPstream::Request* req = nullptr
737 );
738
739 //- Gather/scatter identically-sized data
740 // Send data from proc slot, receive into all slots
741 static void mpi_allgather
742 (
744 void* allData,
746 int count,
747 const UPstream::dataTypes dataTypeId,
748 const int communicator,
750 UPstream::Request* req = nullptr
751 );
752
753 //- Receive variable length data from all ranks,
754 //- placing the result on rank 0.
755 //- (caution: known to scale poorly)
756 static void mpi_gatherv
757 (
758 const void* sendData,
759 int sendCount,
760 void* recvData,
761 const UList<int>& recvCounts,
762 const UList<int>& recvOffsets,
763 const UPstream::dataTypes dataTypeId,
764 const int communicator
765 );
766
767 //- Send variable length data from rank 0 to all ranks.
768 //- (caution: known to scale poorly)
769 static void mpi_scatterv
770 (
771 const void* sendData,
772 const UList<int>& sendCounts,
773 const UList<int>& sendOffsets,
774 void* recvData,
775 int recvCount,
776 const UPstream::dataTypes dataTypeId,
777 const int communicator
778 );
779
780public:
781
782 //- Declare name of the class and its debug switch
783 ClassName("UPstream");
784
785
786 // Static Data
787
788 //- Use of host/node topology-aware routines
789 // 0: disabled
790 // 1: split by hostname [default]
791 // 2: split by shared
792 // >=4: (debug) split with given number per node
793 static int nodeCommsControl_;
794
795 //- Minimum number of nodes before topology-aware routines are enabled
796 // <= 2 : always
797 // >= 3 : when there are more than N nodes
798 static int nodeCommsMin_;
799
800 //- Selection of topology-aware routines as a bitmask combination
801 //- of the topoControls enumerations
802 static int topologyControl_;
803
804 //- Test for selection of given topology-aware routine
805 static bool usingTopoControl(UPstream::topoControls ctrl) noexcept
806 {
807 return static_cast<bool>(topologyControl_ & int(ctrl));
808 }
809
810 //- Should compact transfer be used in which floats replace doubles
811 //- reducing the bandwidth requirement at the expense of some loss
812 //- in accuracy
813 static bool floatTransfer;
814
815 //- Number of processors to change from linear to tree communication
816 static int nProcsSimpleSum;
817
818 //- Number of processors to change to nonBlocking consensual
819 //- exchange (NBX). Ignored for zero or negative values.
820 static int nProcsNonblockingExchange;
821
822 //- Number of polling cycles in processor updates
823 static int nPollProcInterfaces;
824
825 //- Default commsType
827
828 //- Optional maximum message size (bytes)
829 static int maxCommsSize;
830
831 //- Tuning parameters for non-blocking exchange (NBX)
832 static int tuning_NBX_;
833
834 //- MPI buffer-size (bytes)
835 static const int mpiBufferSize;
836
837
838 // Standard Communicators
839
840 //- Communicator for all ranks.
841 //- May differ from commGlobal() if local worlds are in use
842 static label worldComm;
843
844 //- Debugging: warn for use of any communicator differing from warnComm
845 static label warnComm;
846
847 //- Communicator for all ranks, irrespective of any local worlds.
848 // This value \em never changes during a simulation.
849 static constexpr int commGlobal() noexcept { return 0; }
850
851 //- Communicator within the current rank only
852 // This value \em never changes during a simulation.
853 static constexpr int commSelf() noexcept { return 1; }
854
855 //- Communicator for all ranks (respecting any local worlds).
856 // This value \em never changes after startup. Unlike the commWorld()
857 // which can be temporarily overriden.
858 static int commConstWorld() noexcept { return constWorldComm_; }
859
860 //- Communicator for all ranks (respecting any local worlds)
861 static label commWorld() noexcept { return worldComm; }
862
863 //- Set world communicator. Negative values are a no-op.
864 // \returns old world communicator index
865 static label commWorld(const label communicator) noexcept
866 {
867 label old(worldComm);
868 if (communicator >= 0) worldComm = communicator;
869 return old;
870 }
871
872 //- Alter communicator debugging setting.
873 //- Warns for use of any communicator differing from specified.
874 //- Negative values disable.
875 // \returns the previous warn index
876 static label commWarn(const label communicator) noexcept
877 {
878 label old(warnComm);
879 warnComm = communicator;
880 return old;
881 }
882
883 //- Number of currently defined communicators
884 static label nComms() noexcept { return parentComm_.size(); }
885
886 //- Debugging: print the communication tree
887 static void printCommTree(int communicator, bool linear = false);
888
889
890 // Host Communicators
891
892 //- Communicator between nodes/hosts (respects any local worlds)
893 static int commInterNode() noexcept
894 {
895 return (parRun_ ? commInterNode_ : constWorldComm_);
896 }
897
898 //- Communicator within the node/host (respects any local worlds)
899 static int commLocalNode() noexcept
900 {
901 return (parRun_ ? commLocalNode_ : constWorldComm_);
902 }
903
904 //- Both inter-node and local-node communicators have been created
905 static bool hasNodeCommunicators() noexcept
906 {
907 return
908 (
909 (commInterNode_ > constWorldComm_)
910 && (commLocalNode_ > constWorldComm_)
911 );
912 }
913
914 //- True if node topology-aware routines have been enabled,
915 //- it is running in parallel, the starting point is the
916 //- world-communicator and it is not an odd corner case
917 //- (ie, all processes on one node, all processes on different nodes)
918 static bool usingNodeComms(const int communicator);
919
920
921 // Constructors
922
923 //- Construct for given communication type
924 explicit UPstream(const commsTypes commsType) noexcept
925 :
926 commsType_(commsType)
927 {}
928
929
930 // Member Functions
931
932 //- Create new communicator with sub-ranks on the parent communicator
933 static label newCommunicator
934 (
936 const label parent,
937
939 const labelRange& subRanks,
940
942 const bool withComponents = true
943 );
944
945 //- Create new communicator with sub-ranks on the parent communicator
946 static label newCommunicator
947 (
949 const label parent,
950
952 const labelUList& subRanks,
953
955 const bool withComponents = true
956 );
957
958 //- Duplicate the parent communicator
959 //
960 // Always calls dupCommunicatorComponents() internally
961 static label dupCommunicator
962 (
964 const label parent
965 );
966
967 //- Allocate a new communicator by splitting the parent communicator
968 //- on the given \em colour.
969 // Always calls splitCommunicatorComponents() internally
970 static label splitCommunicator
971 (
973 const label parent,
974
977 const int colour,
978
980 const bool two_step = true
981 );
983 //- Free a previously allocated communicator.
984 // Ignores placeholder (negative) communicators.
985 static void freeCommunicator
986 (
987 const label communicator,
988 const bool withComponents = true
989 );
990
991 //- Return physical processor number (i.e. processor number in
992 //- worldComm) given communicator and processor
993 static int baseProcNo(label comm, int procID);
994
995 //- Return processor number in communicator (given physical processor
996 //- number) (= reverse of baseProcNo)
997 static label procNo(const label comm, const int baseProcID);
998
999 //- Return processor number in communicator (given processor number
1000 //- and communicator)
1001 static label procNo
1002 (
1003 const label comm,
1004 const label currentComm,
1005 const int currentProcID
1006 );
1007
1008 //- Add the valid option this type of communications library
1009 //- adds/requires on the command line
1010 static void addValidParOptions(HashTable<string>& validParOptions);
1011
1012 //- Initialisation function called from main
1013 // Spawns sub-processes and initialises inter-communication
1014 static bool init(int& argc, char**& argv, const bool needsThread);
1015
1016 //- Special purpose initialisation function.
1017 // Performs a basic MPI_Init without any other setup.
1018 // Only used for applications that need MPI communication when
1019 // OpenFOAM is running in a non-parallel mode.
1020 // \note Behaves as a no-op if MPI has already been initialized.
1021 // Fatal if MPI has already been finalized.
1022 static bool initNull();
1023
1024 //- Impose a synchronisation barrier (optionally non-blocking)
1025 static void barrier
1026 (
1027 const int communicator,
1028 UPstream::Request* req = nullptr
1030
1031 //- Impose a point-to-point synchronisation barrier
1032 //- by sending a zero-byte \em "done" message to given rank.
1033 // A no-op for non-parallel
1034 static void send_done
1037 const int toProcNo,
1039 const int communicator,
1041 const int tag = UPstream::msgType()+1970
1042 );
1043
1044 //- Impose a point-to-point synchronisation barrier
1045 //- by receiving a zero-byte \em "done" message from given rank
1046 // A no-op for non-parallel
1047 // \returns the source rank (useful for ANY_SOURCE messages)
1048 // or -1 for non-parallel
1049 static int wait_done
1052 const int fromProcNo,
1054 const int communicator,
1056 const int tag = UPstream::msgType()+1970
1057 );
1058
1059 //- Probe for an incoming message.
1061 // \param commsType Non-blocking or not
1062 // \param fromProcNo The source rank (negative == ANY_SOURCE)
1063 // \param tag The source message tag
1064 // \param communicator The communicator index
1065 //
1066 // \returns source rank and message size (bytes)
1067 // and (-1, 0) on failure
1068 static std::pair<int,int64_t> probeMessage
1071 const int fromProcNo,
1072 const int tag = UPstream::msgType(),
1073 const int communicator = worldComm
1075
1076
1077 // Information
1078
1079 //- Report the node-communication settings
1080 static void printNodeCommsControl(Ostream& os);
1082 //- Report the topology routines settings
1083 static void printTopoControl(Ostream& os);
1084
1085
1086 // Requests (non-blocking comms).
1087 // Pending requests are usually handled as an internal (global) list,
1088 // since this simplifies the overall tracking and provides a convenient
1089 // wrapping to avoid exposing MPI-specific types, but can also handle
1090 // with a wrapped UPstream::Request as well.
1091
1092 //- Number of outstanding requests (on the internal list of requests)
1093 static label nRequests() noexcept;
1094
1095 //- Truncate outstanding requests to given length, which is
1096 //- expected to be in the range [0 to nRequests()].
1097 // A no-op for out-of-range values.
1098 static void resetRequests(const label n);
1099
1100 //- Transfer the (wrapped) MPI request to the internal global list
1101 //- and invalidate the parameter (ignores null requests)
1102 // A no-op for non-parallel
1103 static void addRequest(UPstream::Request& req);
1104
1105 //- Non-blocking comms: cancel and free outstanding request.
1106 //- Corresponds to MPI_Cancel() + MPI_Request_free()
1107 // A no-op if parRun() == false
1108 // if there are no pending requests,
1109 // or if the index is out-of-range (0 to nRequests)
1110 static void cancelRequest(const label i);
1111
1112 //- Non-blocking comms: cancel and free outstanding request.
1113 //- Corresponds to MPI_Cancel() + MPI_Request_free()
1114 // A no-op for non-parallel
1115 static void cancelRequest(UPstream::Request& req);
1116
1117 //- Non-blocking comms: cancel and free outstanding requests.
1118 //- Corresponds to MPI_Cancel() + MPI_Request_free()
1119 // A no-op if parRun() == false or list is empty
1120 static void cancelRequests(UList<UPstream::Request>& requests);
1121
1122 //- Non-blocking comms: cancel/free outstanding requests
1123 //- (from position onwards) and remove from internal list of requests.
1124 //- Corresponds to MPI_Cancel() + MPI_Request_free()
1125 // A no-op if parRun() == false,
1126 // if the position is out-of-range [0 to nRequests()],
1127 // or the internal list of requests is empty.
1128 //
1129 // \param pos starting position within the internal list of requests
1130 // \param len length of slice to remove (negative = until the end)
1131 static void removeRequests(label pos, label len = -1);
1133 //- Non-blocking comms: free outstanding request.
1134 //- Corresponds to MPI_Request_free()
1135 // A no-op if parRun() == false
1136 static void freeRequest(UPstream::Request& req);
1137
1138 //- Non-blocking comms: free outstanding requests.
1139 //- Corresponds to MPI_Request_free()
1140 // A no-op if parRun() == false or list is empty
1141 static void freeRequests(UList<UPstream::Request>& requests);
1142
1143 //- Wait until all requests (from position onwards) have finished.
1144 //- Corresponds to MPI_Waitall()
1145 // A no-op if parRun() == false,
1146 // if the position is out-of-range [0 to nRequests()],
1147 // or the internal list of requests is empty.
1148 //
1149 // If checking a trailing portion of the list, it will also trim
1150 // the list of outstanding requests as a side-effect.
1151 // This is a feature (not a bug) to conveniently manange the list.
1152 //
1153 // \param pos starting position within the internal list of requests
1154 // \param len length of slice to check (negative = until the end)
1155 static void waitRequests(label pos, label len = -1);
1156
1157 //- Wait until all requests have finished.
1158 //- Corresponds to MPI_Waitall()
1159 // A no-op if parRun() == false, or the list is empty.
1160 static void waitRequests(UList<UPstream::Request>& requests);
1162 //- Wait until any request (from position onwards) has finished.
1163 //- Corresponds to MPI_Waitany()
1164 // A no-op and returns false if parRun() == false,
1165 // if the position is out-of-range [0 to nRequests()],
1166 // or the internal list of requests is empty.
1167 //
1168 // \returns true if any pending request completed.
1169 // \returns false if all requests have already been handled.
1170 //
1171 // \param pos starting position within the internal list of requests
1172 // \param len length of slice to check (negative = until the end)
1173 static bool waitAnyRequest(label pos, label len = -1);
1174
1175 //- Wait until some requests (from position onwards) have finished.
1176 //- Corresponds to MPI_Waitsome()
1177 // A no-op and returns false if parRun() == false,
1178 // if the position is out-of-range [0 to nRequests],
1179 // or the internal list of requests is empty.
1180 //
1181 // \returns true if some pending requests completed.
1182 // \returns false if all requests have already been handled
1183 //
1184 // \param pos starting position within the internal list of requests
1185 // \param len length of slice to check (negative = until the end)
1186 // \param[out] indices the completed request indices relative to the
1187 // starting position. This is an optional parameter that can be
1188 // used to recover the indices or simply to avoid reallocations
1189 // when calling within a loop.
1190 static bool waitSomeRequests
1191 (
1192 label pos,
1193 label len = -1,
1194 DynamicList<int>* indices = nullptr
1195 );
1196
1197 //- Wait until some requests have finished.
1198 //- Corresponds to MPI_Waitsome()
1199 // A no-op and returns false if parRun() == false,
1200 // the list is empty,
1201 // or if all the requests have already been handled.
1202 //
1203 // \param requests the requests
1204 // \param[out] indices the completed request indices relative to the
1205 // starting position. This is an optional parameter that can be
1206 // used to recover the indices or simply to avoid reallocations
1207 // when calling within a loop.
1208 static bool waitSomeRequests
1209 (
1210 UList<UPstream::Request>& requests,
1211 DynamicList<int>* indices = nullptr
1212 );
1213
1214 //- Wait until any request has finished and return its index.
1215 //- Corresponds to MPI_Waitany()
1216 // Returns -1 if parRun() == false, or the list is empty,
1217 // or if all the requests have already been handled
1218 static int waitAnyRequest(UList<UPstream::Request>& requests);
1219
1220 //- Wait until request i has finished.
1221 //- Corresponds to MPI_Wait()
1222 // A no-op if parRun() == false,
1223 // if there are no pending requests,
1224 // or if the index is out-of-range (0 to nRequests)
1225 static void waitRequest(const label i);
1226
1227 //- Wait until specified request has finished.
1228 //- Corresponds to MPI_Wait()
1229 // A no-op if parRun() == false or for a null-request
1230 static void waitRequest(UPstream::Request& req);
1231
1232 //- Is request \p i active (!= MPI_REQUEST_NULL)?
1233 // False if there are no pending requests,
1234 // or if the index is out-of-range (0 to nRequests)
1235 static bool activeRequest(const label i);
1236
1237 //- Is request active (!= MPI_REQUEST_NULL)?
1238 static bool activeRequest(const UPstream::Request& req);
1239
1240 //- Non-blocking comms: has request i finished?
1241 //- Corresponds to MPI_Test()
1242 // A no-op and returns true if parRun() == false,
1243 // if there are no pending requests,
1244 // or if the index is out-of-range (0 to nRequests)
1245 static bool finishedRequest(const label i);
1246
1247 //- Non-blocking comms: has request finished?
1248 //- Corresponds to MPI_Test()
1249 // A no-op and returns true if parRun() == false
1250 // or for a null-request
1251 static bool finishedRequest(UPstream::Request& req);
1252
1253 //- Non-blocking comms: have all requests (from position onwards)
1254 //- finished?
1255 //- Corresponds to MPI_Testall()
1256 // A no-op and returns true if parRun() == false,
1257 // if there are no pending requests,
1258 // or if the index is out-of-range (0 to nRequests)
1259 // or the addressed range is empty etc.
1260 //
1261 // \param pos starting position within the internal list of requests
1262 // \param len length of slice to check (negative = until the end)
1263 static bool finishedRequests(label pos, label len = -1);
1264
1265 //- Non-blocking comms: have all requests finished?
1266 //- Corresponds to MPI_Testall()
1267 // A no-op and returns true if parRun() == false or list is empty
1268 static bool finishedRequests(UList<UPstream::Request>& requests);
1269
1270 //- Non-blocking comms: have both requests finished?
1271 //- Corresponds to pair of MPI_Test()
1272 // A no-op and returns true if parRun() == false,
1273 // if there are no pending requests,
1274 // or if the indices are out-of-range (0 to nRequests)
1275 // Each finished request parameter is set to -1 (ie, done).
1276 static bool finishedRequestPair(label& req0, label& req1);
1277
1278 //- Non-blocking comms: wait for both requests to finish.
1279 //- Corresponds to pair of MPI_Wait()
1280 // A no-op if parRun() == false,
1281 // if there are no pending requests,
1282 // or if the indices are out-of-range (0 to nRequests)
1283 // Each finished request parameter is set to -1 (ie, done).
1284 static void waitRequestPair(label& req0, label& req1);
1285
1286
1287 // General
1288
1289 //- Set as parallel run on/off.
1290 // \return the previous value
1291 static bool parRun(const bool on) noexcept
1292 {
1293 bool old(parRun_);
1294 parRun_ = on;
1295 return old;
1296 }
1297
1298 //- Test if this a parallel run
1299 // Modify access is deprecated
1300 static bool& parRun() noexcept { return parRun_; }
1301
1302 //- Have support for threads
1303 static bool haveThreads() noexcept { return haveThreads_; }
1304
1305 //- Relative rank for the master process - is always 0.
1306 static constexpr int masterNo() noexcept { return 0; }
1307
1308 //- Number of ranks in parallel run (for given communicator).
1309 //- It is 1 for serial run
1310 static label nProcs(const label communicator = worldComm)
1311 {
1312 return procIDs_[communicator].size();
1313 }
1314
1315 //- Rank of this process in the communicator (starting from masterNo()).
1316 //- Negative if the process is not a rank in the communicator.
1317 static int myProcNo(const label communicator = worldComm)
1318 {
1319 return myProcNo_[communicator];
1320 }
1321
1322 //- True if process corresponds to the master rank in the communicator
1323 static bool master(const label communicator = worldComm)
1324 {
1325 return myProcNo_[communicator] == masterNo();
1326 }
1327
1328 //- True if process corresponds to \b any rank (master or sub-rank)
1329 //- in the given communicator
1330 static bool is_rank(const label communicator = worldComm)
1331 {
1332 return myProcNo_[communicator] >= masterNo();
1333 }
1334
1335 //- True if process corresponds to a sub-rank in the given communicator
1336 static bool is_subrank(const label communicator = worldComm)
1337 {
1338 return myProcNo_[communicator] > masterNo();
1339 }
1340
1341 //- True if parallel algorithm or exchange is required.
1342 // This is when parRun() == true, the process corresponds to a rank
1343 // in the communicator and there is more than one rank in the
1344 // communicator
1345 static bool is_parallel(const label communicator = worldComm)
1346 {
1347 return
1348 (
1349 parRun_ && is_rank(communicator) && nProcs(communicator) > 1
1350 );
1351 }
1352
1353 //- The number of shared/host nodes in the (const) world communicator.
1354 static int numNodes() noexcept { return numNodes_; }
1355
1356 //- The parent communicator
1357 static label parent(int communicator)
1358 {
1359 return parentComm_(communicator);
1360 }
1361
1362 //- The list of ranks within a given communicator
1363 static List<int>& procID(int communicator)
1364 {
1365 return procIDs_[communicator];
1366 }
1367
1368
1369 //- Test for communicator equality.
1370 // True if they have the same index or address the same ranks
1371 static bool sameProcs(int communicator1, int communicator2)
1372 {
1373 return sameProcs_impl(communicator1, communicator2);
1374 }
1375
1376 //- Test equality of communicator procs with the given list of ranks.
1377 //- Includes a guard for the communicator index.
1378 template
1379 <
1380 typename T1,
1381 typename = std::void_t
1382 <std::enable_if_t<std::is_integral_v<T1>>>
1383 >
1384 static bool sameProcs(int communicator, const UList<T1>& procs)
1385 {
1386 return sameProcs_impl(communicator, procs);
1387 }
1388
1389 //- Test the equality of two lists of ranks
1390 template
1391 <
1392 typename T1,
1393 typename T2,
1394 typename = std::void_t
1395 <
1396 std::enable_if_t<std::is_integral_v<T1>>,
1397 std::enable_if_t<std::is_integral_v<T2>>
1398 >
1399 >
1400 static bool sameProcs(const UList<T1>& procs1, const UList<T2>& procs2)
1401 {
1402 return sameProcs_impl(procs1, procs2);
1403 }
1404
1405
1406 // Worlds
1407
1408 //- All worlds
1409 static const wordList& allWorlds() noexcept
1410 {
1411 return allWorlds_;
1412 }
1413
1414 //- The indices into allWorlds for all processes
1415 static const labelList& worldIDs() noexcept
1416 {
1417 return worldIDs_;
1418 }
1419
1420 //- My worldID
1421 static label myWorldID()
1422 {
1423 return worldIDs_[myProcNo(UPstream::commGlobal())];
1424 }
1425
1426 //- My world
1427 static const word& myWorld()
1428 {
1429 return allWorlds_[worldIDs_[myProcNo(UPstream::commGlobal())]];
1430 }
1431
1432
1433 // Rank addressing
1434
1435 //- Range of process indices for all processes
1436 static rangeType allProcs(const label communicator = worldComm)
1437 {
1438 // Proc 0 -> nProcs (int value)
1439 return rangeType(static_cast<int>(nProcs(communicator)));
1440 }
1441
1442 //- Range of process indices for sub-processes
1443 static rangeType subProcs(const label communicator = worldComm)
1444 {
1445 // Proc 1 -> nProcs (int value)
1446 return rangeType(1, static_cast<int>(nProcs(communicator)-1));
1447 }
1448
1449 //- Processor offsets corresponding to the inter-node communicator
1450 static const List<int>& interNode_offsets();
1451
1452 //- The range (start/size) of the commLocalNode ranks in terms of the
1453 //- (const) world communicator processors.
1454 static const rangeType& localNode_parentProcs();
1455
1456 //- Linear communication schedule (special case) for given communicator
1457 static const commsStructList& linearCommunication(int communicator);
1458
1459 //- Tree communication schedule (standard case) for given communicator
1460 static const commsStructList& treeCommunication(int communicator);
1461
1462 //- Communication schedule for all-to-master (proc 0) as
1463 //- linear/tree/none with switching based on UPstream::nProcsSimpleSum,
1464 //- the is_parallel() state and the optional \c linear parameter.
1465 static const commsStructList& whichCommunication
1466 (
1467 const int communicator,
1469 bool linear = false
1470 )
1471 {
1472 const label np
1473 (
1474 parRun_ && is_rank(communicator) // cf. is_parallel()
1475 ? nProcs(communicator)
1476 : label(0)
1477 );
1478
1479 return
1480 (
1481 np <= 1
1483 : (linear || np < UPstream::nProcsSimpleSum)
1484 ? linearCommunication(communicator)
1485 : treeCommunication(communicator)
1486 );
1487 }
1488
1489
1490 //- Message tag of standard messages
1491 static int& msgType() noexcept
1492 {
1493 return msgType_;
1494 }
1495
1496 //- Set the message tag for standard messages
1497 // \return the previous value
1498 static int msgType(int val) noexcept
1499 {
1500 int old(msgType_);
1501 msgType_ = val;
1502 return old;
1503 }
1504
1505 //- Increment the message tag for standard messages
1506 // \return the previous value
1507 static int incrMsgType(int val = 1) noexcept
1508 {
1509 int old(msgType_);
1510 msgType_ += val;
1511 return old;
1512 }
1513
1514 //- Get the communications type of the stream
1516 {
1517 return commsType_;
1518 }
1519
1520 //- Set the communications type of the stream
1521 // \return the previous value
1522 commsTypes commsType(const commsTypes ct) noexcept
1523 {
1524 commsTypes old(commsType_);
1525 commsType_ = ct;
1526 return old;
1527 }
1528
1529
1530 //- Shutdown (finalize) MPI as required.
1531 // Uses MPI_Abort instead of MPI_Finalize if errNo is non-zero
1532 static void shutdown(int errNo = 0);
1533
1534 //- Call MPI_Abort with no other checks or cleanup
1535 static void abort(int errNo = 1);
1536
1537 //- Shutdown (finalize) MPI as required and exit program with errNo.
1538 static void exit(int errNo = 1);
1539
1540
1541 #undef Pstream_CommonRoutines
1542 #define Pstream_CommonRoutines(Type) \
1543 \
1544 \
1545 \
1546 static void allToAll \
1547 ( \
1548 \
1549 const UList<Type>& sendData, \
1550 \
1551 UList<Type>& recvData, \
1552 const int communicator = UPstream::worldComm \
1553 ); \
1554 \
1555 \
1556 \
1557 \
1558 \
1559 \
1560 \
1561 \
1562 \
1563 static void allToAllConsensus \
1564 ( \
1565 \
1566 const UList<Type>& sendData, \
1567 \
1568 UList<Type>& recvData, \
1569 \
1570 const int tag, \
1571 const int communicator = UPstream::worldComm \
1572 ); \
1573 \
1574 \
1575 \
1576 \
1577 \
1578 \
1579 static void allToAllConsensus \
1580 ( \
1581 \
1582 const Map<Type>& sendData, \
1583 \
1584 Map<Type>& recvData, \
1585 \
1586 const int tag, \
1587 const int communicator = UPstream::worldComm \
1588 ); \
1589 \
1590 \
1591 \
1592 static Map<Type> allToAllConsensus \
1593 ( \
1594 \
1595 const Map<Type>& sendData, \
1596 \
1597 const int tag, \
1598 const int communicator = UPstream::worldComm \
1599 ) \
1600 { \
1601 Map<Type> recvData; \
1602 allToAllConsensus(sendData, recvData, tag, communicator); \
1603 return recvData; \
1604 }
1605
1606 Pstream_CommonRoutines(int32_t);
1607 Pstream_CommonRoutines(int64_t);
1608
1609 #undef Pstream_CommonRoutines
1610
1611
1612 // Low-level gather/scatter routines
1613
1614 //- Receive identically-sized (contiguous) data from all ranks
1615 template<class Type>
1616 static void mpiGather
1617 (
1619 const Type* sendData,
1622 Type* recvData,
1624 int count,
1625 const int communicator = UPstream::worldComm
1626 );
1627
1628 //- Send identically-sized (contiguous) data to all ranks
1629 template<class Type>
1630 static void mpiScatter
1631 (
1634 const Type* sendData,
1636 Type* recvData,
1638 int count,
1639 const int communicator = UPstream::worldComm
1640 );
1641
1642 //- Gather/scatter identically-sized data
1643 // Send data from proc slot, receive into all slots
1644 template<class Type>
1645 static void mpiAllGather
1646 (
1648 Type* allData,
1650 int count,
1651 const int communicator = UPstream::worldComm
1652 );
1653
1654 //- Receive variable length data from all ranks
1655 template<class Type>
1656 static void mpiGatherv
1657 (
1658 const Type* sendData,
1659 int sendCount,
1660 Type* recvData,
1661 const UList<int>& recvCounts,
1662 const UList<int>& recvOffsets,
1663 const int communicator = UPstream::worldComm
1664 );
1665
1666 //- Send variable length data to all ranks
1667 template<class Type>
1668 static void mpiScatterv
1670 const Type* sendData,
1671 const UList<int>& sendCounts,
1672 const UList<int>& sendOffsets,
1673 Type* recvData,
1674 int recvCount,
1676 );
1677
1678
1679 // Gather single, contiguous value(s)
1680
1681 //- Allgather individual values into list locations.
1682 // The returned list has size nProcs, identical on all ranks.
1683 template<class T>
1685 (
1686 const T& localValue,
1688 );
1689
1690 //- Gather individual values into list locations.
1691 // On master list length == nProcs, otherwise zero length.
1692 // \n
1693 // For \b non-parallel :
1694 // the returned list length is 1 with localValue.
1695 template<class T>
1698 const T& localValue,
1700 );
1701
1702 //- Scatter individual values from list locations.
1703 // On master input list length == nProcs, ignored on other procs.
1704 // \n
1705 // For \b non-parallel :
1706 // returns the first list element (or default initialized).
1707 template<class T>
1708 static T listScatterValues
1709 (
1710 const UList<T>& allValues,
1712 );
1713
1715 // Broadcast Functions
1716
1717 //- Broadcast buffer contents (contiguous types) to all ranks
1718 //- (default: from rank=0). The sizes must match on all processes.
1719 // For \b non-parallel : do nothing.
1720 // \return True on success
1721 template<class Type>
1722 inline static bool broadcast
1724 Type* buffer,
1725 std::streamsize count,
1726 const int communicator,
1727 const int root = UPstream::masterNo()
1728 );
1729
1730 //- Broadcast fixed-list content (contiguous types) to all ranks
1731 //- (default: from rank=0). The sizes must match on all processes.
1732 // For \b non-parallel : do nothing.
1733 // \return True on success
1734 template<class Type, unsigned N>
1735 inline static bool broadcast
1736 (
1737 FixedList<Type, N>& list,
1738 const int communicator,
1739 const int root = UPstream::masterNo()
1740 );
1741
1742
1743 // Reductions
1744
1745 //- MPI_Reduce (blocking) for known operators
1746 // For \b non-parallel : do nothing.
1747 template<class T>
1748 static void mpiReduce
1749 (
1750 T values[],
1751 int count,
1752 const UPstream::opCodes opCodeId,
1753 const int communicator
1755
1756 //- MPI_Reduce (blocking) for known operators
1757 // For \b non-parallel : do nothing.
1758 template<UPstream::opCodes opCode, class T>
1759 static void mpiReduce
1760 (
1761 T values[],
1762 int count,
1763 const int communicator
1764 );
1765
1766 //- MPI_Allreduce (blocking) for known operators
1767 // For \b non-parallel : do nothing.
1768 template<class T>
1769 static void mpiAllReduce
1770 (
1771 T values[],
1772 int count,
1773 const UPstream::opCodes opCodeId,
1774 const int communicator
1775 );
1776
1777 //- MPI_Allreduce (blocking) for known operators
1778 // For \b non-parallel : do nothing.
1779 template<UPstream::opCodes opCode, class T>
1780 static void mpiAllReduce
1781 (
1782 T values[],
1783 int count,
1784 const int communicator
1785 );
1786
1787 //- MPI_Iallreduce (non-blocking) for known operators
1788 // For \b non-parallel : do nothing.
1789 template<class T>
1790 static void mpiAllReduce
1791 (
1792 T values[],
1793 int count,
1794 const UPstream::opCodes opCodeId,
1795 const int communicator,
1797 );
1798
1799 //- MPI_Iallreduce (non-blocking) for known operators
1800 // For \b non-parallel : do nothing.
1801 template<UPstream::opCodes opCode, class T>
1802 static void mpiAllReduce
1803 (
1804 T values[],
1805 int count,
1806 const int communicator,
1808 );
1809
1810
1811 // Partial reductions
1812
1813 //- Inclusive/exclusive scan (in-place).
1814 // In exclusive mode, the degenerate value on rank=0 has no meaning
1815 // and will normally be treated like non-exclusive mode
1816 // (ie, original values).
1817 // However, for the \b \c opCodes::op_sum type we can provide a
1818 // sensible value and overwrite it with a zero-initialized value.
1819 //
1820 // \note For \b non-parallel : do nothing.
1821 template<Foam::UPstream::opCodes opCode, class T>
1822 static void mpiScan
1823 (
1824 T values[],
1825 int count,
1826 const int communicator,
1828 const bool exclusive = false
1829 );
1831 //- Inclusive/exclusive scan returning the result.
1832 //- In exclusive mode, the degenerate value on rank=0 has no meaning
1833 //- but is treated like non-exclusive mode (ie, original values)
1834 //
1835 // \note For \b non-parallel : return own value
1836 template<Foam::UPstream::opCodes opCode, class T>
1837 static T mpiScan
1839 const T& localValue,
1840 const int communicator,
1842 const bool exclusive = false
1843 );
1844
1845 #undef Pstream_ScanRoutines
1846 #define Pstream_ScanRoutines(OpName, OpCode) \
1847 \
1848 \
1849 \
1850 template<class T> \
1851 static void mpiScan_##OpName \
1852 ( \
1853 T values[], \
1854 int count, \
1855 const int communicator, \
1856 \
1857 const bool exclusive = false \
1858 ) \
1859 { \
1860 UPstream::mpiScan<OpCode>(values, count, communicator, exclusive);\
1861 } \
1862 \
1863 \
1864 \
1865 template<class T> \
1866 static void mpiExscan_##OpName \
1867 ( \
1868 T values[], \
1869 int count, \
1870 const int communicator \
1871 ) \
1872 { \
1873 UPstream::mpiScan<OpCode>(values, count, communicator, true); \
1874 } \
1875 \
1876 \
1877 \
1878 template<class T> \
1879 static T mpiScan_##OpName \
1880 ( \
1881 const T& value, \
1882 const int communicator, \
1883 \
1884 const bool exclusive = false \
1885 ) \
1886 { \
1887 return UPstream::mpiScan<OpCode>(value, communicator, exclusive); \
1888 } \
1889 \
1890 \
1891 \
1892 template<class T> \
1893 static T mpiExscan_##OpName \
1894 ( \
1895 const T& value, \
1896 const int communicator \
1897 ) \
1898 { \
1899 return UPstream::mpiScan<OpCode>(value, communicator, true); \
1900 }
1901
1905 #undef Pstream_ScanRoutines
1906
1907
1908 // Logical reductions
1909
1910 //- Logical (and) reduction (MPI_AllReduce)
1911 // For \b non-parallel : do nothing
1912 static void reduceAnd
1913 (
1914 bool& value,
1915 const int communicator = worldComm
1916 );
1917
1918 //- Logical (or) reduction (MPI_AllReduce)
1919 // For \b non-parallel : do nothing
1920 static void reduceOr
1921 (
1922 bool& value,
1923 const int communicator = worldComm
1924 );
1925
1926 //- Locate the first rank for which the \p condition is true,
1927 //- or -1 if no ranks satisfy the condition.
1928 static int find_first(bool condition, int communicator);
1929
1930 //- Locate the last rank for which the \p condition is true,
1931 //- or -1 if no ranks satisfy the condition.
1932 static int find_last(bool condition, int communicator);
1933
1934
1935 // Housekeeping
1937 //- \deprecated(2025-02) prefer newCommunicator
1938 FOAM_DEPRECATED_FOR(2025-02, "newCommunicator()")
1939 static label allocateCommunicator
1940 (
1941 const label parent,
1942 const labelRange& subRanks,
1943 const bool withComponents = true
1944 )
1945 {
1946 return newCommunicator(parent, subRanks, withComponents);
1947 }
1949 //- \deprecated(2025-02) prefer newCommunicator
1950 FOAM_DEPRECATED_FOR(2025-02, "newCommunicator()")
1951 static label allocateCommunicator
1952 (
1953 const label parent,
1954 const labelUList& subRanks,
1955 const bool withComponents = true
1956 )
1957 {
1958 return newCommunicator(parent, subRanks, withComponents);
1959 }
1960
1961 //- Communicator between nodes (respects any local worlds)
1962 FOAM_DEPRECATED_FOR(2025-02, "commInterNode()")
1963 static label commInterHost() noexcept { return commInterNode(); }
1964
1965 //- Communicator within the node (respects any local worlds)
1966 FOAM_DEPRECATED_FOR(2025-02, "commLocalNode()")
1967 static label commIntraHost() noexcept { return commLocalNode(); }
1969 //- Wait for all requests to finish.
1970 // \deprecated(2023-01) Probably not what you want.
1971 // Should normally be restricted to a particular starting request.
1972 FOAM_DEPRECATED_FOR(2023-01, "waitRequests(int) method")
1973 static void waitRequests() { waitRequests(0); }
1974
1975 //- \deprecated(2025-02) prefer mpiGatherv
1976 template<class Type>
1977 FOAM_DEPRECATED_FOR(2025-02, "mpiGatherv()")
1978 static void gather
1979 (
1980 const Type* send,
1981 int count,
1982 Type* recv,
1983 const UList<int>& counts,
1984 const UList<int>& offsets,
1985 const int comm = UPstream::worldComm
1986 )
1987 {
1988 UPstream::mpiGatherv(send, count, recv, counts, offsets, comm);
1989 }
1990
1991 //- \deprecated(2025-02) prefer mpiScatterv
1992 template<class Type>
1993 FOAM_DEPRECATED_FOR(2025-02, "mpiScatterv()")
1994 static void scatter
1996 const Type* send,
1997 const UList<int>& counts,
1998 const UList<int>& offsets,
1999 Type* recv,
2000 int count,
2001 const int comm = UPstream::worldComm
2002 )
2003 {
2004 UPstream::mpiScatterv(send, counts, offsets, recv, count, comm);
2005 }
2006};
2007
2008/*---------------------------------------------------------------------------*\
2009 Class UPstream::communicator Declaration
2010\*---------------------------------------------------------------------------*/
2011
2012//- Wrapper for internally indexed communicator label.
2013//- Always invokes UPstream::allocateCommunicatorComponents()
2014//- and UPstream::freeCommunicatorComponents()
2015class UPstream::communicator
2016{
2017 // Private Data
2018
2019 //- The communicator label
2020 label comm_;
2021
2022public:
2023
2024 // Generated Methods
2025
2026 //- No copy construct
2027 communicator(const communicator&) = delete;
2028
2029 //- No copy assignment
2030 void operator=(const communicator&) = delete;
2031
2032
2033 // Constructors
2034
2035 //- Default construct (a placeholder communicator)
2036 communicator() noexcept : comm_(-1) {}
2037
2038 //- Move construct, takes ownership
2039 communicator(communicator&& c) : comm_(c.comm_) { c.comm_ = -1; }
2040
2041 //- Allocate communicator for contiguous sub-ranks on given parent
2042 communicator
2043 (
2045 const label parentComm,
2047 const labelRange& subRanks
2048 )
2049 :
2050 comm_(UPstream::newCommunicator(parentComm, subRanks))
2051 {}
2052
2053 //- Allocate communicator for sub-ranks on given parent
2054 communicator
2055 (
2057 const label parentComm,
2058
2059 const labelUList& subRanks
2061 :
2062 comm_(UPstream::newCommunicator(parentComm, subRanks))
2063 {}
2064
2065
2066 // Destructor
2067
2068 //- Free allocated communicator
2069 ~communicator() { UPstream::freeCommunicator(comm_); }
2070
2072 // Factory Methods
2073
2074 //- Duplicate the given communicator
2075 static communicator duplicate(const label parentComm)
2076 {
2077 communicator c;
2078 c.comm_ = UPstream::dupCommunicator(parentComm);
2079 return c;
2080 }
2081
2082 //- Factory Method :
2083 //- Split the communicator on the given \em colour.
2084 static communicator split
2085 (
2087 const label parentComm,
2090 const int colour,
2092 const bool two_step = true
2093 )
2094 {
2095 communicator c;
2096 c.comm_ = UPstream::splitCommunicator(parentComm, colour, two_step);
2097 return c;
2098 }
2099
2100
2101 // Member Functions
2102
2103 //- True if communicator is non-negative (ie, was allocated)
2104 bool good() const noexcept { return (comm_ >= 0); }
2106 //- The communicator label
2107 label comm() const noexcept { return comm_; }
2108
2109 //- Return non-const reference to this
2110 communicator& constCast() const noexcept
2111 {
2112 return const_cast<communicator&>(*this);
2113 }
2114
2115 //- Release ownership of the communicator, return old value.
2116 // Leave further management to the caller
2117 label release() noexcept { label c(comm_); comm_ = -1; return c; }
2119 //- Free allocated communicator
2120 void reset() { UPstream::freeCommunicator(comm_); comm_ = -1; }
2121
2122 //- Allocate with contiguous sub-ranks of parent communicator
2123 void reset(label parent, const labelRange& subRanks)
2124 {
2126 comm_ = UPstream::newCommunicator(parent, subRanks);
2127 }
2128
2129 //- Allocate with sub-ranks of parent communicator
2130 void reset(label parent, const labelUList& subRanks)
2131 {
2133 comm_ = UPstream::newCommunicator(parent, subRanks);
2134 }
2135
2136 //- Take ownership, free managed communicator
2137 void reset(communicator&& c)
2138 {
2139 if (this == &c) return; // No self-assignment
2140 if (comm_ != c.comm_) UPstream::freeCommunicator(comm_);
2141 comm_ = c.comm_;
2142 c.comm_ = -1;
2143 }
2144
2145 //- Swap communicator labels
2146 void swap(communicator& c) { std::swap(comm_, c.comm_); }
2147
2148
2149 // Member Operators
2150
2151 //- Implicit cast to label - the same as comm()
2152 operator label() const noexcept { return comm_; }
2153
2154 //- Move assignment, takes ownership
2155 void operator=(communicator&& c) { reset(std::move(c)); }
2156
2157 //- Test for equality
2158 bool operator==(const communicator& c) const noexcept
2159 {
2160 return (comm_ == c.comm_);
2161 }
2162
2163 //- Test for inequality
2164 bool operator!=(const communicator& c) const noexcept
2165 {
2166 return (comm_ != c.comm_);
2167 }
2168};
2169
2170
2171/*---------------------------------------------------------------------------*\
2172 Class UPstream::Communicator Declaration
2173\*---------------------------------------------------------------------------*/
2174
2175//- An opaque wrapper for MPI_Comm with a vendor-independent
2176//- representation without any \c <mpi.h> header.
2177// The MPI standard states that MPI_Comm is always an opaque object.
2178// Generally it is either an integer (eg, mpich) or a pointer (eg, openmpi).
2179class UPstream::Communicator
2180{
2181public:
2182
2183 // Public Types
2184
2185 //- Storage for MPI_Comm (as integer or pointer)
2186 typedef std::intptr_t value_type;
2187
2188
2189private:
2190
2191 // Private Data
2192
2193 //- The MPI_Comm (as wrapped value)
2194 value_type value_;
2195
2196public:
2197
2198 // Generated Methods
2199
2200 //- Copy construct
2201 Communicator(const Communicator&) noexcept = default;
2202
2203 //- Move construct
2204 Communicator(Communicator&&) noexcept = default;
2205
2206 //- Copy assignment
2207 Communicator& operator=(const Communicator&) noexcept = default;
2208
2209 //- Move assignment
2210 Communicator& operator=(Communicator&&) noexcept = default;
2211
2212
2213 // Member Operators
2215 //- Test for equality
2216 bool operator==(const Communicator& rhs) const noexcept
2217 {
2218 return (value_ == rhs.value_);
2219 }
2220
2221 //- Test for inequality
2222 bool operator!=(const Communicator& rhs) const noexcept
2223 {
2224 return (value_ != rhs.value_);
2225 }
2226
2227
2228 // Constructors
2229
2230 //- Default construct as MPI_COMM_NULL
2232
2233 //- Construct from MPI_Comm (as pointer type)
2234 explicit Communicator(const void* p) noexcept
2235 :
2236 value_(reinterpret_cast<value_type>(p))
2237 {}
2238
2239 //- Construct from MPI_Comm (as integer type)
2240 explicit Communicator(value_type val) noexcept
2241 :
2242 value_(val)
2243 {}
2245
2246 // Factory Methods
2247
2248 //- Transcribe internally indexed communicator to wrapped value.
2249 //
2250 // Example,
2251 // \code
2252 // PstreamUtils::Cast::to_mpi
2253 // (
2254 // UPstream::Communicator::lookup(UPstream::worldComm)
2255 // )
2256 // \endcode
2258 (
2261 int communicator = -1
2262 );
2263
2264
2265 // Member Functions
2266
2267 // Access
2268
2269 //- Return raw value
2270 value_type value() const noexcept { return value_; }
2272 //- Return as pointer value
2273 const void* pointer() const noexcept
2274 {
2275 return reinterpret_cast<const void*>(value_);
2276 }
2277
2278
2279 // Query
2280
2281 //- True if not equal to MPI_COMM_NULL
2282 bool good() const noexcept;
2283
2284 //- Reset to default constructed value (MPI_COMM_NULL)
2285 void reset() noexcept;
2286
2287 //- The number of ranks associated with the communicator
2288 int size() const;
2289};
2290
2291
2292/*---------------------------------------------------------------------------*\
2293 Class UPstream::Request Declaration
2294\*---------------------------------------------------------------------------*/
2295
2296//- An opaque wrapper for MPI_Request with a vendor-independent
2297//- representation without any \c <mpi.h> header.
2298// The MPI standard states that MPI_Request is always an opaque object.
2299// Generally it is either an integer (eg, mpich) or a pointer (eg, openmpi).
2300class UPstream::Request
2301{
2302public:
2303
2304 // Public Types
2305
2306 //- Storage for MPI_Request (as integer or pointer)
2307 typedef std::intptr_t value_type;
2308
2309
2310private:
2311
2312 // Private Data
2313
2314 //- The MPI_Request (as wrapped value)
2315 value_type value_;
2316
2317public:
2318
2319 // Generated Methods
2320
2321 //- Copy construct
2322 Request(const Request&) noexcept = default;
2323
2324 //- Move construct
2325 Request(Request&&) noexcept = default;
2326
2327 //- Copy assignment
2328 Request& operator=(const Request&) noexcept = default;
2329
2330 //- Move assignment
2331 Request& operator=(Request&&) noexcept = default;
2332
2333
2334 // Member Operators
2335
2336 //- Test for equality
2337 bool operator==(const Request& rhs) const noexcept
2338 {
2339 return (value_ == rhs.value_);
2341
2342 //- Test for inequality
2343 bool operator!=(const Request& rhs) const noexcept
2344 {
2345 return (value_ != rhs.value_);
2346 }
2347
2348
2349 // Constructors
2350
2351 //- Default construct as MPI_REQUEST_NULL
2352 Request() noexcept;
2353
2354 //- Construct from MPI_Request (as pointer type)
2355 explicit Request(const void* p) noexcept
2356 :
2357 value_(reinterpret_cast<value_type>(p))
2358 {}
2359
2360 //- Construct from MPI_Request (as integer type)
2361 explicit Request(value_type val) noexcept
2362 :
2363 value_(val)
2364 {}
2365
2366
2367 // Member Functions
2368
2369 // Access
2370
2371 //- Return raw value
2372 value_type value() const noexcept { return value_; }
2373
2374 //- Return as pointer value
2375 const void* pointer() const noexcept
2376 {
2377 return reinterpret_cast<const void*>(value_);
2378 }
2379
2380
2381 // Basics
2382
2383 //- True if not equal to MPI_REQUEST_NULL
2384 bool good() const noexcept;
2385
2386 //- Reset to default constructed value (MPI_REQUEST_NULL)
2387 void reset() noexcept;
2388
2389 //- True if request is active (!= MPI_REQUEST_NULL)
2390 //- Same as good(). Same as calling UPstream::activeRequest()
2391 bool active() const noexcept { return good(); }
2392
2393 //- Same as calling UPstream::finishedRequest()
2394 // Uses MPI_Test()
2395 bool finished() { return UPstream::finishedRequest(*this); }
2396
2397 //- Same as calling UPstream::waitRequest().
2398 // Uses MPI_Wait()
2399 void wait() { UPstream::waitRequest(*this); }
2400
2401
2402 // Other
2403
2404 //- Same as calling UPstream::cancelRequest().
2405 // Uses MPI_Cancel(), MPI_Request_free()
2406 void cancel() { UPstream::cancelRequest(*this); }
2408 //- Same as calling UPstream::freeRequest().
2409 // Uses MPI_Request_free()
2410 void free() { UPstream::freeRequest(*this); }
2411};
2412
2413
2414// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
2415
2416Ostream& operator<<(Ostream&, const UPstream::commsStruct&);
2417
2418
2419// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2420
2421} // End namespace Foam
2422
2423// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2424
2425// Traits, nested classes etc
2426#include "UPstreamTraits.H"
2427#include "UPstreamWindow.H"
2428
2429// Locally include the following where required:
2430// - UPstreamFile.H
2431
2432// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2433
2434#ifdef NoRepository
2435 #include "UPstream.txx"
2436 #include "UPstreamReduceOffsets.H"
2437#endif
2438
2439// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
2440
2441#endif
2442
2443// ************************************************************************* //
Various functions to operate on Lists.
Implementation of IO operators for Foam::OffsetRange.
label n
A set of traits associated with UPstream communication.
#define Pstream_ScanRoutines(OpName, OpCode)
Definition UPstream.H:2349
#define Pstream_CommonRoutines(Type)
Definition UPstream.H:1995
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition Enum.H:57
A 1D vector of objects of type <T> with a fixed length <N>.
Definition FixedList.H:73
A HashTable similar to std::unordered_map.
Definition HashTable.H:124
An interval of (signed) integers defined by a start and a size.
Definition IntRange.H:69
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void clear()
Clear the list, i.e. set size to zero.
Definition ListI.H:133
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
An opaque wrapper for MPI_Comm with a vendor-independent representation without any <mpi....
Definition UPstream.H:2760
Communicator(const Communicator &) noexcept=default
Copy construct.
const void * pointer() const noexcept
Return as pointer value.
Definition UPstream.H:2882
value_type value() const noexcept
Return raw value.
Definition UPstream.H:2877
std::intptr_t value_type
Storage for MPI_Comm (as integer or pointer).
Definition UPstream.H:2768
bool operator!=(const Communicator &rhs) const noexcept
Test for inequality.
Definition UPstream.H:2818
Communicator(Communicator &&) noexcept=default
Move construct.
Communicator(value_type val) noexcept
Construct from MPI_Comm (as integer type).
Definition UPstream.H:2842
An opaque wrapper for MPI_File methods without any <mpi.h> header dependency.
An opaque wrapper for MPI_Request with a vendor-independent representation without any <mpi....
Definition UPstream.H:2919
void cancel()
Same as calling UPstream::cancelRequest().
Definition UPstream.H:3065
Request(const Request &) noexcept=default
Copy construct.
Request(Request &&) noexcept=default
Move construct.
const void * pointer() const noexcept
Return as pointer value.
Definition UPstream.H:3019
bool operator!=(const Request &rhs) const noexcept
Test for inequality.
Definition UPstream.H:2977
value_type value() const noexcept
Return raw value.
Definition UPstream.H:3014
bool good() const noexcept
True if not equal to MPI_REQUEST_NULL.
bool active() const noexcept
True if request is active (!= MPI_REQUEST_NULL) Same as good(). Same as calling UPstream::activeReque...
Definition UPstream.H:3041
std::intptr_t value_type
Storage for MPI_Request (as integer or pointer).
Definition UPstream.H:2927
void wait()
Same as calling UPstream::waitRequest().
Definition UPstream.H:3055
void free()
Same as calling UPstream::freeRequest().
Definition UPstream.H:3072
bool finished()
Same as calling UPstream::finishedRequest().
Definition UPstream.H:3048
Request(value_type val) noexcept
Construct from MPI_Request (as integer type).
Definition UPstream.H:3001
An opaque wrapper for MPI_Win with a vendor-independent representation and without any <mpi....
Collection of communication structures.
Definition UPstream.H:363
const UPstream::commsStruct & get(int proci) const
Get existing or create (demand-driven) entry.
bool linear() const noexcept
Linear (flat) communication instead of tree communication.
Definition UPstream.H:427
bool empty() const noexcept
True if the list is empty.
Definition UPstream.H:443
commsStructList() noexcept
Construct empty with invalid communicator.
Definition UPstream.H:388
bool good() const noexcept
True if communicator is non-negative (ie, was assigned).
Definition UPstream.H:417
const UPstream::commsStruct & operator[](int proci) const
Get existing or create (demand-driven) entry.
Definition UPstream.H:480
static const commsStructList & null()
An empty structure. Used for placeholders etc.
int comm() const noexcept
The communicator internal index.
Definition UPstream.H:422
Ostream & writeList(Ostream &os) const
Write the communication tree.
Definition UPstream.H:493
commsStructList(int comm, bool flat=false) noexcept
Construct empty with given communicator.
Definition UPstream.H:397
label size() const noexcept
The number of entries.
Definition UPstream.H:448
void clear()
Clear the list.
Definition UPstream.H:438
void printGraph(Ostream &os, int proci=0) const
Print un-directed graph in graphviz dot format.
Structure for communicating between processors.
Definition UPstream.H:229
int nProcs() const noexcept
The number of processors addressed by the structure.
void reset_linear(const int myProci, const int numProcs)
Reset to linear (flat) communication.
const List< int > & below() const noexcept
The procIDs of the processors directly below.
Definition UPstream.H:300
commsStruct() noexcept
Default construct with above == -1.
Definition UPstream.H:262
const List< int > & allNotBelow() const noexcept
The procIDs of all processors not below myProcNo. The inverse set of allBelow without myProcNo.
Definition UPstream.H:312
void reset()
Reset to default constructed state.
int above() const noexcept
The procID of the processor directly above.
Definition UPstream.H:295
const List< int > & allBelow() const noexcept
The procIDs of all processors below (so not just directly below).
Definition UPstream.H:306
Wrapper for internally indexed communicator label. Always invokes UPstream::allocateCommunicatorCompo...
Definition UPstream.H:2546
communicator(const label parentComm, const labelUList &subRanks)
Allocate communicator for sub-ranks on given parent.
Definition UPstream.H:2599
communicator(communicator &&c)
Move construct, takes ownership.
Definition UPstream.H:2579
void reset(communicator &&c)
Take ownership, free managed communicator.
Definition UPstream.H:2704
static communicator split(const label parentComm, const int colour, const bool two_step=true)
Factory Method : Split the communicator on the given colour.
Definition UPstream.H:2635
bool good() const noexcept
True if communicator is non-negative (ie, was allocated).
Definition UPstream.H:2656
communicator(const label parentComm, const labelRange &subRanks)
Allocate communicator for contiguous sub-ranks on given parent.
Definition UPstream.H:2585
void operator=(communicator &&c)
Move assignment, takes ownership.
Definition UPstream.H:2728
label release() noexcept
Release ownership of the communicator, return old value.
Definition UPstream.H:2676
void reset(label parent, const labelRange &subRanks)
Allocate with contiguous sub-ranks of parent communicator.
Definition UPstream.H:2686
communicator & constCast() const noexcept
Return non-const reference to this.
Definition UPstream.H:2666
bool operator!=(const communicator &c) const noexcept
Test for inequality.
Definition UPstream.H:2741
void swap(communicator &c)
Swap communicator labels.
Definition UPstream.H:2715
label comm() const noexcept
The communicator label.
Definition UPstream.H:2661
~communicator()
Free allocated communicator.
Definition UPstream.H:2615
communicator(const communicator &)=delete
No copy construct.
void reset()
Free allocated communicator.
Definition UPstream.H:2681
bool operator==(const communicator &c) const noexcept
Test for equality.
Definition UPstream.H:2733
static communicator duplicate(const label parentComm)
Duplicate the given communicator.
Definition UPstream.H:2623
communicator() noexcept
Default construct (a placeholder communicator).
Definition UPstream.H:2574
void operator=(const communicator &)=delete
No copy assignment.
void reset(label parent, const labelUList &subRanks)
Allocate with sub-ranks of parent communicator.
Definition UPstream.H:2695
Inter-processor communications stream.
Definition UPstream.H:69
static void freeRequest(UPstream::Request &req)
Non-blocking comms: free outstanding request. Corresponds to MPI_Request_free().
static void mpiAllReduce(T values[], int count, const int communicator)
MPI_Allreduce (blocking) for known operators.
static label commWorld() noexcept
Communicator for all ranks (respecting any local worlds).
Definition UPstream.H:1101
static constexpr int commSelf() noexcept
Communicator within the current rank only.
Definition UPstream.H:1088
static bool broadcast(FixedList< Type, N > &list, const int communicator, const int root=UPstream::masterNo())
Broadcast fixed-list content (contiguous types) to all ranks (default: from rank=0)....
static std::pair< int, int64_t > probeMessage(const UPstream::commsTypes commsType, const int fromProcNo, const int tag=UPstream::msgType(), const int communicator=worldComm)
Probe for an incoming message.
Definition UPstream.C:132
static void mpi_reduce(void *values, int count, const UPstream::dataTypes dataTypeId, const UPstream::opCodes opCodeId, const int communicator, UPstream::Request *req=nullptr)
In-place reduction of values with result on rank 0.
ClassName("UPstream")
Declare name of the class and its debug switch.
static bool waitSomeRequests(label pos, label len=-1, DynamicList< int > *indices=nullptr)
Wait until some requests (from position onwards) have finished. Corresponds to MPI_Waitsome().
static int wait_done(const int fromProcNo, const int communicator, const int tag=UPstream::msgType()+1970)
Impose a point-to-point synchronisation barrier by receiving a zero-byte "done" message from given ra...
Definition UPstream.C:120
static void printCommTree(int communicator, bool linear=false)
Debugging: print the communication tree.
Definition UPstream.C:737
static void mpiAllReduce(T values[], int count, const UPstream::opCodes opCodeId, const int communicator, UPstream::Request &req)
MPI_Iallreduce (non-blocking) for known operators.
static void mpiGatherv(const Type *sendData, int sendCount, Type *recvData, const UList< int > &recvCounts, const UList< int > &recvOffsets, const int communicator=UPstream::worldComm)
Receive variable length data from all ranks.
static bool finishedRequests(label pos, label len=-1)
Non-blocking comms: have all requests (from position onwards) finished? Corresponds to MPI_Testall().
static void mpiReduce(T values[], int count, const UPstream::opCodes opCodeId, const int communicator)
MPI_Reduce (blocking) for known operators.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static void mpi_scatter(const void *sendData, void *recvData, int count, const UPstream::dataTypes dataTypeId, const int communicator, UPstream::Request *req=nullptr)
Send identically-sized (contiguous) data from rank 0 to all other ranks.
static void addRequest(UPstream::Request &req)
Transfer the (wrapped) MPI request to the internal global list and invalidate the parameter (ignores ...
static const commsStructList & whichCommunication(const int communicator, bool linear=false)
Communication schedule for all-to-master (proc 0) as linear/tree/none with switching based on UPstrea...
Definition UPstream.H:1899
static List< int > & procID(int communicator)
The list of ranks within a given communicator.
Definition UPstream.H:1767
static bool init(int &argc, char **&argv, const bool needsThread)
Initialisation function called from main.
Definition UPstream.C:40
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
static label parent(int communicator)
The parent communicator.
Definition UPstream.H:1759
commsTypes
Communications types.
Definition UPstream.H:81
@ blocking
compatibility name for buffered
Definition UPstream.H:86
@ scheduled
"scheduled" (MPI standard) : (MPI_Send, MPI_Recv)
Definition UPstream.H:83
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
static label commInterHost() noexcept
Communicator between nodes (respects any local worlds).
Definition UPstream.H:2482
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static label warnComm
Debugging: warn for use of any communicator differing from warnComm.
Definition UPstream.H:1074
static void mpiReduce(T values[], int count, const int communicator)
MPI_Reduce (blocking) for known operators.
static bool initNull()
Special purpose initialisation function.
Definition UPstream.C:30
static void shutdown(int errNo=0)
Shutdown (finalize) MPI as required.
Definition UPstream.C:57
static const Enum< commsTypes > commsTypeNames
Enumerated names for the communication types.
Definition UPstream.H:92
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static bool floatTransfer
Should compact transfer be used in which floats replace doubles reducing the bandwidth requirement at...
Definition UPstream.H:1024
static void printTopoControl(Ostream &os)
Report the topology routines settings.
Definition UPstream.C:100
static bool mpi_broadcast(void *buf, std::streamsize count, const UPstream::dataTypes dataTypeId, const int communicator, const int root=0)
Broadcast buffer contents to all ranks (default: from rank=0). The sizes must match on all processes.
static void mpiScatterv(const Type *sendData, const UList< int > &sendCounts, const UList< int > &sendOffsets, Type *recvData, int recvCount, const int communicator=UPstream::worldComm)
Send variable length data to all ranks.
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition UPstream.H:1691
static int numNodes() noexcept
The number of shared/host nodes in the (const) world communicator.
Definition UPstream.H:1754
static void mpi_allgather(void *allData, int count, const UPstream::dataTypes dataTypeId, const int communicator, UPstream::Request *req=nullptr)
Gather/scatter identically-sized data.
static void mpiGather(const Type *sendData, Type *recvData, int count, const int communicator=UPstream::worldComm)
Receive identically-sized (contiguous) data from all ranks.
static bool broadcast(Type *buffer, std::streamsize count, const int communicator, const int root=UPstream::masterNo())
Broadcast buffer contents (contiguous types) to all ranks (default: from rank=0). The sizes must matc...
static void cancelRequest(const label i)
Non-blocking comms: cancel and free outstanding request. Corresponds to MPI_Cancel() + MPI_Request_fr...
static void mpiScatter(const Type *sendData, Type *recvData, int count, const int communicator=UPstream::worldComm)
Send identically-sized (contiguous) data to all ranks.
static int topologyControl_
Selection of topology-aware routines as a bitmask combination of the topoControls enumerations.
Definition UPstream.H:1009
static void send_done(const int toProcNo, const int communicator, const int tag=UPstream::msgType()+1970)
Impose a point-to-point synchronisation barrier by sending a zero-byte "done" message to given rank.
Definition UPstream.C:111
static std::streamsize mpi_receive(const UPstream::commsTypes commsType, void *buf, std::streamsize count, const UPstream::dataTypes dataTypeId, const int fromProcNo, const int tag, const int communicator, UPstream::Request *req=nullptr)
Receive buffer contents of specified data type from given processor.
static int nodeCommsMin_
Minimum number of nodes before topology-aware routines are enabled.
Definition UPstream.H:1003
static label dupCommunicator(const label parent)
Duplicate the parent communicator.
Definition UPstream.C:401
static List< T > listGatherValues(const T &localValue, const int communicator=UPstream::worldComm)
Gather individual values into list locations.
static void barrier(const int communicator, UPstream::Request *req=nullptr)
Impose a synchronisation barrier (optionally non-blocking).
Definition UPstream.C:106
opCodes
Mapping of some MPI op codes.
Definition UPstream.H:149
@ op_bool_and
Logical and.
Definition UPstream.H:160
@ Extra_end
(internal use) end marker [window types]
Definition UPstream.H:178
@ op_no_op
No-op (window only).
Definition UPstream.H:173
@ op_bit_xor
Bit-wise xor for (unsigned) integral types.
Definition UPstream.H:165
@ Extra_begin
(internal use) begin marker [window types]
Definition UPstream.H:168
@ OpCodes_end
(internal use) end marker [all types]
Definition UPstream.H:178
@ op_replace
Replace (window only).
Definition UPstream.H:172
@ op_bit_or
Bit-wise or for (unsigned) integral types.
Definition UPstream.H:164
@ op_bool_xor
Logical xor.
Definition UPstream.H:162
@ op_bit_and
Bit-wise and for (unsigned) integral types.
Definition UPstream.H:163
@ op_bool_or
Logical or.
Definition UPstream.H:161
static const int mpiBufferSize
MPI buffer-size (bytes).
Definition UPstream.H:1060
static const rangeType & localNode_parentProcs()
The range (start/size) of the commLocalNode ranks in terms of the (const) world communicator processo...
Definition UPstream.C:817
static void reduceAnd(bool &value, const int communicator=worldComm)
Logical (and) reduction (MPI_AllReduce).
static bool usingTopoControl(UPstream::topoControls ctrl) noexcept
Test for selection of given topology-aware routine.
Definition UPstream.H:1014
static int commLocalNode() noexcept
Communicator within the node/host (respects any local worlds).
Definition UPstream.H:1153
static void gather(const Type *send, int count, Type *recv, const UList< int > &counts, const UList< int > &offsets, const int comm=UPstream::worldComm)
Definition UPstream.H:2505
static void mpiAllReduce(T values[], int count, const UPstream::opCodes opCodeId, const int communicator)
MPI_Allreduce (blocking) for known operators.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static void reduceOr(bool &value, const int communicator=worldComm)
Logical (or) reduction (MPI_AllReduce).
static bool hasNodeCommunicators() noexcept
Both inter-node and local-node communicators have been created.
Definition UPstream.H:1161
static void resetRequests(const label n)
Truncate outstanding requests to given length, which is expected to be in the range [0 to nRequests()...
static bool is_subrank(const label communicator=worldComm)
True if process corresponds to a sub-rank in the given communicator.
Definition UPstream.H:1731
static const labelList & worldIDs() noexcept
The indices into allWorlds for all processes.
Definition UPstream.H:1830
static int nProcsNonblockingExchange
Number of processors to change to nonBlocking consensual exchange (NBX). Ignored for zero or negative...
Definition UPstream.H:1035
static label commWorld(const label communicator) noexcept
Set world communicator. Negative values are a no-op.
Definition UPstream.H:1108
static void exit(int errNo=1)
Shutdown (finalize) MPI as required and exit program with errNo.
Definition UPstream.C:61
static const commsStructList & linearCommunication(int communicator)
Linear communication schedule (special case) for given communicator.
Definition UPstream.C:699
static bool sameProcs(const UList< T1 > &procs1, const UList< T2 > &procs2)
Test the equality of two lists of ranks.
Definition UPstream.H:1811
sendModes
Different MPI-send modes (ignored for commsTypes::buffered).
Definition UPstream.H:98
@ sync
(MPI_Ssend, MPI_Issend)
Definition UPstream.H:100
@ normal
(MPI_Send, MPI_Isend)
Definition UPstream.H:99
static bool is_parallel(const label communicator=worldComm)
True if parallel algorithm or exchange is required.
Definition UPstream.H:1743
static void freeRequests(UList< UPstream::Request > &requests)
Non-blocking comms: free outstanding requests. Corresponds to MPI_Request_free().
static int tuning_NBX_
Tuning parameters for non-blocking exchange (NBX).
Definition UPstream.H:1055
commsTypes commsType(const commsTypes ct) noexcept
Set the communications type of the stream.
Definition UPstream.H:1968
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
commsTypes commsType() const noexcept
Get the communications type of the stream.
Definition UPstream.H:1958
static label allocateCommunicator(const label parent, const labelRange &subRanks, const bool withComponents=true)
Definition UPstream.H:2455
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
static void abort(int errNo=1)
Call MPI_Abort with no other checks or cleanup.
Definition UPstream.C:68
static label nComms() noexcept
Number of currently defined communicators.
Definition UPstream.H:1132
static bool finishedRequestPair(label &req0, label &req1)
Non-blocking comms: have both requests finished? Corresponds to pair of MPI_Test().
IntRange< int > rangeType
Int ranges are used for MPI ranks (processes).
Definition UPstream.H:75
static bool haveThreads() noexcept
Have support for threads.
Definition UPstream.H:1686
static label splitCommunicator(const label parent, const int colour, const bool two_step=true)
Allocate a new communicator by splitting the parent communicator on the given colour.
Definition UPstream.C:439
static label myWorldID()
My worldID.
Definition UPstream.H:1838
static void cancelRequests(UList< UPstream::Request > &requests)
Non-blocking comms: cancel and free outstanding requests. Corresponds to MPI_Cancel() + MPI_Request_f...
static void mpi_scatterv(const void *sendData, const UList< int > &sendCounts, const UList< int > &sendOffsets, void *recvData, int recvCount, const UPstream::dataTypes dataTypeId, const int communicator)
Send variable length data from rank 0 to all ranks. (caution: known to scale poorly).
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition UPstream.H:1857
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition UPstream.H:1866
static int find_last(bool condition, int communicator)
Locate the last rank for which the condition is true, or -1 if no ranks satisfy the condition.
Definition UPstream.C:1087
static int nProcsSimpleSum
Number of processors to change from linear to tree communication.
Definition UPstream.H:1029
static int commInterNode() noexcept
Communicator between nodes/hosts (respects any local worlds).
Definition UPstream.H:1145
static bool usingNodeComms(const int communicator)
True if node topology-aware routines have been enabled, it is running in parallel,...
Definition UPstream.C:751
static int nPollProcInterfaces
Number of polling cycles in processor updates.
Definition UPstream.H:1040
static bool is_rank(const label communicator=worldComm)
True if process corresponds to any rank (master or sub-rank) in the given communicator.
Definition UPstream.H:1723
static void mpi_allreduce(void *values, int count, const UPstream::dataTypes dataTypeId, const UPstream::opCodes opCodeId, const int communicator, UPstream::Request *req=nullptr)
In-place reduction of values with same result on all ranks.
static int nodeCommsControl_
Use of host/node topology-aware routines.
Definition UPstream.H:995
static const commsStructList & treeCommunication(int communicator)
Tree communication schedule (standard case) for given communicator.
Definition UPstream.C:718
static bool sameProcs(int communicator1, int communicator2)
Test for communicator equality.
Definition UPstream.H:1778
static void mpiScan(T values[], int count, const int communicator, const bool exclusive=false)
Inclusive/exclusive scan (in-place).
static void mpi_scan_reduce(void *values, int count, const UPstream::dataTypes dataTypeId, const UPstream::opCodes opCodeId, const int communicator, const bool exclusive)
In-place scan/exscan reduction of values.
static void mpi_gather(const void *sendData, void *recvData, int count, const UPstream::dataTypes dataTypeId, const int communicator, UPstream::Request *req=nullptr)
Receive identically-sized (contiguous) data from all ranks, placing the result on rank 0.
static int commConstWorld() noexcept
Communicator for all ranks (respecting any local worlds).
Definition UPstream.H:1096
static T listScatterValues(const UList< T > &allValues, const int communicator=UPstream::worldComm)
Scatter individual values from list locations.
static T mpiScan(const T &localValue, const int communicator, const bool exclusive=false)
Inclusive/exclusive scan returning the result. In exclusive mode, the degenerate value on rank=0 has ...
static void mpiAllGather(Type *allData, int count, const int communicator=UPstream::worldComm)
Gather/scatter identically-sized data.
static label commWarn(const label communicator) noexcept
Alter communicator debugging setting. Warns for use of any communicator differing from specified....
Definition UPstream.H:1122
static const wordList & allWorlds() noexcept
All worlds.
Definition UPstream.H:1822
static bool waitAnyRequest(label pos, label len=-1)
Wait until any request (from position onwards) has finished. Corresponds to MPI_Waitany().
static label commIntraHost() noexcept
Communicator within the node (respects any local worlds).
Definition UPstream.H:2488
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
static void addValidParOptions(HashTable< string > &validParOptions)
Add the valid option this type of communications library adds/requires on the command line.
Definition UPstream.C:26
static void removeRequests(label pos, label len=-1)
Non-blocking comms: cancel/free outstanding requests (from position onwards) and remove from internal...
static int find_first(bool condition, int communicator)
Locate the first rank for which the condition is true, or -1 if no ranks satisfy the condition.
Definition UPstream.C:1055
static const List< int > & interNode_offsets()
Processor offsets corresponding to the inter-node communicator.
Definition UPstream.C:776
static constexpr int commGlobal() noexcept
Communicator for all ranks, irrespective of any local worlds.
Definition UPstream.H:1081
static void waitRequestPair(label &req0, label &req1)
Non-blocking comms: wait for both requests to finish. Corresponds to pair of MPI_Wait().
static int msgType(int val) noexcept
Set the message tag for standard messages.
Definition UPstream.H:1936
static void waitRequests()
Wait for all requests to finish.
Definition UPstream.H:2497
static void waitRequests(label pos, label len=-1)
Wait until all requests (from position onwards) have finished. Corresponds to MPI_Waitall().
UPstream(const commsTypes commsType) noexcept
Construct for given communication type.
Definition UPstream.H:1184
static List< T > allGatherValues(const T &localValue, const int communicator=UPstream::worldComm)
Allgather individual values into list locations.
static label newCommunicator(const label parent, const labelRange &subRanks, const bool withComponents=true)
Create new communicator with sub-ranks on the parent communicator.
Definition UPstream.C:272
static int maxCommsSize
Optional maximum message size (bytes).
Definition UPstream.H:1050
static bool mpi_send(const UPstream::commsTypes commsType, const void *buf, std::streamsize count, const UPstream::dataTypes dataTypeId, const int toProcNo, const int tag, const int communicator, UPstream::Request *req=nullptr, const UPstream::sendModes sendMode=UPstream::sendModes::normal)
Send buffer contents of specified data type to given processor.
static label procNo(const label comm, const int baseProcID)
Return processor number in communicator (given physical processor number) (= reverse of baseProcNo).
Definition UPstream.C:670
static bool sameProcs(int communicator, const UList< T1 > &procs)
Test equality of communicator procs with the given list of ranks. Includes a guard for the communicat...
Definition UPstream.H:1793
static int baseProcNo(label comm, int procID)
Return physical processor number (i.e. processor number in worldComm) given communicator and processo...
Definition UPstream.C:657
topoControls
Some bit masks corresponding to topology controls.
Definition UPstream.H:188
@ mapGather
mapGather (reduction) [manual algorithm]
Definition UPstream.H:193
@ gatherList
gatherList [manual algorithm]
Definition UPstream.H:194
static bool activeRequest(const label i)
Is request i active (!= MPI_REQUEST_NULL)?
static void printNodeCommsControl(Ostream &os)
Report the node-communication settings.
Definition UPstream.C:56
static void scatter(const Type *send, const UList< int > &counts, const UList< int > &offsets, Type *recv, int count, const int comm=UPstream::worldComm)
Definition UPstream.H:2523
static void mpiAllReduce(T values[], int count, const int communicator, UPstream::Request &req)
MPI_Iallreduce (non-blocking) for known operators.
static void freeCommunicator(const label communicator, const bool withComponents=true)
Free a previously allocated communicator.
Definition UPstream.C:622
static commsTypes defaultCommsType
Default commsType.
Definition UPstream.H:1045
dataTypes
Mapping of some fundamental and aggregate types to MPI data types.
Definition UPstream.H:107
@ Basic_end
(internal use) end marker [basic types]
Definition UPstream.H:125
@ type_3double
3*double (eg, doubleVector)
Definition UPstream.H:130
@ DataTypes_end
(internal use) end marker [all types]
Definition UPstream.H:139
@ type_9float
9*float (eg, floatTensor)
Definition UPstream.H:133
@ User_begin
(internal use) begin marker [user types]
Definition UPstream.H:125
@ type_6float
6*float (eg, floatSymmTensor, complex vector)
Definition UPstream.H:131
@ type_6double
6*double (eg, doubleSymmTensor, complex vector)
Definition UPstream.H:132
@ type_3float
3*float (eg, floatVector)
Definition UPstream.H:129
@ type_byte
byte, char, unsigned char, ...
Definition UPstream.H:113
@ User_end
(internal use) end marker [user types]
Definition UPstream.H:139
@ Basic_begin
(internal use) begin marker [basic/all types]
Definition UPstream.H:112
@ type_9double
9*double (eg, doubleTensor)
Definition UPstream.H:134
static void mpi_gatherv(const void *sendData, int sendCount, void *recvData, const UList< int > &recvCounts, const UList< int > &recvOffsets, const UPstream::dataTypes dataTypeId, const int communicator)
Receive variable length data from all ranks, placing the result on rank 0. (caution: known to scale p...
static bool finishedRequest(const label i)
Non-blocking comms: has request i finished? Corresponds to MPI_Test().
static void waitRequest(const label i)
Wait until request i has finished. Corresponds to MPI_Wait().
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
static const word & myWorld()
My world.
Definition UPstream.H:1846
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
Central-differencing interpolation scheme class.
Definition linear.H:54
Lookup type of boundary radiation properties.
Definition lookup.H:60
A class for handling words, derived from Foam::string.
Definition word.H:66
#define ClassName(TypeNameString)
Add typeName information from argument TypeNameString to a class.
Definition className.H:74
volScalarField & p
limits reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL))
static bool split(const std::string &line, std::string &key, std::string &val)
Definition cpuInfo.C:32
OBJstream os(runTime.globalPath()/outputName)
surface1 clear()
bool equal(const UList< Type1 > &a, const UList< Type2 > &b)
Test for list equality with different but compatible data types. Eg, int32 and int64.
Implementation details for UPstream/Pstream/MPI etc.
Definition UPstream.H:57
Interface handling for UPstream/Pstream/MPI etc.
Definition UPstream.H:62
const dimensionedScalar c
Speed of light in a vacuum.
Namespace for OpenFOAM.
bool operator!=(const eddy &a, const eddy &b)
Definition eddy.H:297
dimensionedScalar pos(const dimensionedScalar &ds)
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< label > labelList
A List of labels.
Definition List.H:62
Ostream & operator<<(Ostream &, const boundaryPatch &p)
Write boundaryPatch as dictionary entries (without surrounding braces).
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
void rhs(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
const direction noexcept
Definition scalarImpl.H:265
UList< label > labelUList
A UList of labels.
Definition UList.H:75
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
#define FOAM_DEPRECATED_FOR(since, replacement)
Definition stdFoam.H:43