Loading...
Searching...
No Matches
UPstream.C
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
27Note
28 Included by global/globals.C
29
30\*---------------------------------------------------------------------------*/
31
32#include "UPstream.H"
33#include "debug.H"
34#include "registerSwitch.H"
35#include "dictionary.H"
36#include "SHA1.H"
37#include "OSspecific.H" // for hostName()
38#include "IOstreams.H"
39
40// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41
42namespace Foam
43{
45}
46
47const Foam::Enum
48<
50>
52({
53 { commsTypes::buffered, "buffered" }, // "buffered"
54 { commsTypes::scheduled, "scheduled" },
55 { commsTypes::nonBlocking, "nonBlocking" }, // "immediate"
56 // compatibility names
57 { commsTypes::buffered, "blocking" },
58});
59
60
61// * * * * * * * * * * * * * Controls Information * * * * * * * * * * * * * //
62
64{
66 {
68 {
69 os << "on [";
70 }
71 else
72 {
73 os << "off [";
74 }
76 {
77 os << "min=" << UPstream::nodeCommsMin_ << ",";
78 }
79 os << "type=";
80
81 // 1: split by hostname [default]
82 // 2: split by shared
83 // >=4: (debug/manual) split with given number per node
85 {
87 }
88 else if (UPstream::nodeCommsControl_ == 2)
89 {
90 os << "shared";
91 }
92 else
93 {
94 os << "host";
95 }
96 os << "]";
97 }
98 else
99 {
100 os << "disabled";
101 }
102 os << " (" << UPstream::nProcs(UPstream::worldComm) << " ranks, "
103 << UPstream::numNodes() << " nodes)";
104}
105
106
108{
109 unsigned count = 0;
110
112 {
113 #undef PrintControl
114 #define PrintControl(Ctrl, Name) \
115 if (UPstream::usingTopoControl(topoControls::Ctrl)) \
116 { \
117 os << (count++ ? ' ' : '(') << Name; \
118 }
119
120 PrintControl(broadcast, "broadcast");
121 PrintControl(reduce, "reduce");
122 PrintControl(gather, "gather");
123 PrintControl(combine, "combine");
124 PrintControl(mapGather, "mapGather");
125 PrintControl(gatherList, "gatherList");
126
127 #undef PrintControl
128 }
129
130 if (count)
131 {
132 os << ')'; // End the list
133 }
134 else
135 {
136 os << "none";
137 }
138}
139
140
141// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
142
143void Foam::UPstream::setParRun(const int nProcs, const bool haveThreads)
144{
145 parRun_ = (nProcs > 0);
146 haveThreads_ = haveThreads;
147
148 label comm = -1;
149 labelRange singleProc(1);
150
151 // Redo communicators that were created during static initialisation.
152 // When parRun == true, redo with MPI components
153 // When parRun == false, just redo in case of future changes
154
155 if (!parRun_)
156 {
157 // These are already correct from the static initialisation,
158 // but just in case of future changes
159
160 // Using (world, self) ordering
161 freeCommunicator(UPstream::commSelf());
162 freeCommunicator(UPstream::commGlobal());
163
164 // 0: COMM_WORLD : commWorld() / commGlobal()
165 comm = newCommunicator(-1, singleProc, false);
166 if (comm != UPstream::commGlobal())
167 {
168 // Failed sanity check
170 << "problem : comm:" << comm
171 << " expected comm-global:" << UPstream::commGlobal()
173 }
174
175 // 1: COMM_SELF
176 comm = newCommunicator(-2, singleProc, false);
177 if (comm != UPstream::commSelf())
178 {
179 // Failed sanity check
181 << "problem : comm:" << comm
182 << " expected comm-self:" << UPstream::commSelf()
184 }
185
186 Pout.prefix().clear();
187 Perr.prefix().clear();
188 }
189 else
190 {
191 // Redo communicators that were created during static initialisation
192 // but this time with MPI components
193
194 // Using (world, self) ordering
195 freeCommunicator(UPstream::commSelf());
196 freeCommunicator(UPstream::commGlobal());
197
198 // 0: COMM_WORLD : commWorld() / commGlobal()
199 comm = newCommunicator(-1, labelRange(nProcs), true);
200 if (comm != UPstream::commGlobal())
201 {
202 // Failed sanity check
204 << "problem : comm:" << comm
205 << " expected comm-global:" << UPstream::commGlobal()
207 }
208
209 const int globalRanki = UPstream::myProcNo(UPstream::commGlobal());
210
211 // 1: COMM_SELF
212 // - Processor number wrt world communicator
213 singleProc.start() = globalRanki;
214 comm = newCommunicator(-2, singleProc, true);
215 if (comm != UPstream::commSelf())
216 {
217 // Failed sanity check
219 << "problem : comm:" << comm
220 << " expected comm-self:" << UPstream::commSelf()
222 }
223
224 Pout.prefix() = '[' + std::to_string(globalRanki) + "] ";
225 Perr.prefix() = Pout.prefix();
226 }
227
228 if (debug)
229 {
230 Perr<< "UPstream::setParRun :"
231 << " nProcs:" << nProcs
232 << " haveThreads:" << haveThreads
233 << endl;
234 }
235}
236
237
238Foam::label Foam::UPstream::getAvailableCommIndex(const label parentIndex)
239{
240 label index;
241 if (!freeComms_.empty())
242 {
243 // LIFO pop
244 index = freeComms_.back();
245 freeComms_.pop_back();
246
247 // Reset existing
248 myProcNo_[index] = -1;
249 parentComm_[index] = parentIndex;
250
251 procIDs_[index].clear();
252 // Sizing and filling are demand-driven
253 linearCommunication_[index].clear();
254 treeCommunication_[index].clear();
255 }
256 else
257 {
258 // Extend storage
259 index = parentComm_.size();
260
261 myProcNo_.push_back(-1);
262 parentComm_.push_back(parentIndex);
263
264 procIDs_.emplace_back();
265 // Sizing and filling are demand-driven
266 linearCommunication_.emplace_back(index);
267 treeCommunication_.emplace_back(index);
268 }
269
270 // Set the communication pattern
271 linearCommunication_[index].linear(true);
272 treeCommunication_[index].linear(false);
273
274 return index;
275}
276
277
279(
280 const label parentIndex,
281 const labelRange& subRanks,
282 const bool withComponents
283)
284{
285 const label index = getAvailableCommIndex(parentIndex);
286
287 if (debug)
288 {
289 Perr<< "Allocate communicator ["
290 << index << "] from [" << parentIndex
291 << "] ranks : " << subRanks << nl
292 << endl;
293 }
294
295 // Initially treat as master,
296 // overwritten by allocateCommunicatorComponents
297 myProcNo_[index] = UPstream::masterNo();
298 auto& procIds = procIDs_[index];
299
300 // The selected sub-ranks.
301 // - transcribe from label to int
302 // - already in monotonic order
303 if
304 (
305 (withComponents && UPstream::parRun())
306 ? (parentIndex < 0 || subRanks.contains(myProcNo_[parentIndex]))
307 : !subRanks.empty()
308 )
309 {
310 procIds.resize_nocopy(subRanks.size());
311 std::iota(procIds.begin(), procIds.end(), subRanks.start());
312 }
313 else
314 {
315 // Not involved
316 procIds.clear();
317 }
318
319 if (withComponents && UPstream::parRun())
320 {
321 allocateCommunicatorComponents(parentIndex, index);
322 }
323
324 return index;
325}
326
327
329(
330 const label parentIndex,
331 const labelUList& subRanks,
332 const bool withComponents
333)
334{
335 const label index = getAvailableCommIndex(parentIndex);
336
337 if (debug)
338 {
339 Perr<< "Allocate communicator ["
340 << index << "] from [" << parentIndex
341 << "] ranks : " << flatOutput(subRanks) << nl
342 << endl;
343 }
344
345 // Initially treat as master,
346 // overwritten by allocateCommunicatorComponents
347 myProcNo_[index] = UPstream::masterNo();
348 auto& procIds = procIDs_[index];
349
350 // The selected sub-ranks.
351 // - transcribe from label to int
352 // - sort into monotonic order (if needed)
353 if
354 (
355 (withComponents && UPstream::parRun())
356 ? (parentIndex < 0 || subRanks.contains(myProcNo_[parentIndex]))
357 : !subRanks.empty()
358 )
359 {
360 procIds.resize_nocopy(subRanks.size());
361
362 label count = 0;
363 bool monotonicOrder = true;
364 for (const auto ranki : subRanks)
365 {
366 if (ranki < 0)
367 {
368 continue;
369 }
370 // Could also flag/ignore out-of-range ranks
371 // (ranki >= numProcs)
372
373 if (monotonicOrder && count)
374 {
375 monotonicOrder = (procIds[count-1] < ranki);
376 }
377
378 procIds[count] = ranki;
379 ++count;
380 }
381
382 if (!monotonicOrder)
383 {
384 auto last = procIds.begin() + count;
385 std::sort(procIds.begin(), last);
386 last = std::unique(procIds.begin(), last);
387 count = label(last - procIds.begin());
388 }
389
390 procIds.resize(count);
391 }
392 else
393 {
394 // Not involved
395 procIds.clear();
396 }
397
398 if (withComponents && UPstream::parRun())
399 {
400 allocateCommunicatorComponents(parentIndex, index);
401 }
402
403 return index;
404}
405
406
408(
409 const label parentIndex
410)
411{
412 #ifdef FULLDEBUG
413 if (FOAM_UNLIKELY(parentIndex < 0))
414 {
415 // Failed sanity check
417 << "Attempted to duplicate an invalid communicator: "
418 << parentIndex
420 }
421 #endif
422
423 const label index = getAvailableCommIndex(parentIndex);
424
425 if (debug)
426 {
427 Perr<< "Duplicate communicator ["
428 << index << "] from [" << parentIndex << "]" << endl;
429 }
430
431 // Initially treat as unknown,
432 // overwritten by dupCommunicatorComponents
433 myProcNo_[index] = -1;
434 procIDs_[index].clear();
435
436 if (UPstream::parRun())
437 {
438 dupCommunicatorComponents(parentIndex, index);
439 }
440
441 return index;
442}
443
444
446(
447 const label parentIndex,
448 const int colour,
449 const bool two_step
450)
451{
452 #ifdef FULLDEBUG
453 if (FOAM_UNLIKELY(parentIndex < 0))
454 {
455 // Failed sanity check
457 << "Attempted to split an invalid communicator: "
458 << parentIndex
460 }
461 #endif
462
463 const label index = getAvailableCommIndex(parentIndex);
464
465 if (debug)
466 {
467 Perr<< "Split communicator ["
468 << index << "] from [" << parentIndex
469 << "] using colour=" << colour
470 << " (two_step=" << two_step << ")" << endl;
471 }
472
473 // Initially treat as unknown,
474 // overwritten by splitCommunicatorComponents
475 myProcNo_[index] = -1;
476 procIDs_[index].clear();
477
478 if (UPstream::parRun())
479 {
480 splitCommunicatorComponents(parentIndex, index, colour, two_step);
481 }
482
483 return index;
484}
485
486
487bool Foam::UPstream::setHostCommunicators(const int numPerNode)
488{
489 // Uses the world communicator (not global communicator)
490
491 // Skip if non-parallel
492 if (!UPstream::parRun())
493 {
494 numNodes_ = 1;
495 return false;
496 }
497
498 if (FOAM_UNLIKELY(commInterNode_ >= 0 || commLocalNode_ >= 0))
499 {
500 // Failed sanity check
502 << "Node communicator(s) already created!" << endl
504 return false;
505 }
506
507 commInterNode_ = getAvailableCommIndex(constWorldComm_);
508 commLocalNode_ = getAvailableCommIndex(constWorldComm_);
509
510 // Overwritten later
511 myProcNo_[commInterNode_] = UPstream::masterNo();
512 myProcNo_[commLocalNode_] = UPstream::masterNo();
513
514 // Sorted order, purely cosmetic
515 if (commLocalNode_ < commInterNode_)
516 {
517 std::swap(commLocalNode_, commInterNode_);
518 }
519
520 if (debug)
521 {
522 Perr<< "Allocating node communicators "
523 << commInterNode_ << ", " << commLocalNode_
524 << " on parent : " << constWorldComm_ << nl
525 << endl;
526 }
527
528 const int worldRank = UPstream::myProcNo(constWorldComm_);
529 const int worldSize = UPstream::nProcs(constWorldComm_);
530
531 if (numPerNode > 1)
532 {
533 // Manual splitting based on given number of ranks per node
534 const int myNodeId = (worldRank/numPerNode);
535
536 // Establish the topology
537 {
538 DynamicList<int> nodeGroup(numPerNode);
539 DynamicList<int> nodeLeader(1+worldSize/numPerNode);
540
541 for (int proci = 0; proci < worldSize; ++proci)
542 {
543 if (myNodeId == (proci/numPerNode))
544 {
545 nodeGroup.push_back(proci);
546 }
547
548 if ((proci % numPerNode) == 0)
549 {
550 // Local rank 0 is a node leader
551 nodeLeader.push_back(proci);
552 }
553 }
554
555 procIDs_[commInterNode_] = std::move(nodeLeader);
556 procIDs_[commLocalNode_] = std::move(nodeGroup);
557 }
558 }
559 else
560 {
561 // Determine inter-host/inter-host grouping based on the SHA1 of the
562 // hostnames. This allows a single initial Allgather to establish
563 // the overall topology. The alternative is to use MPI_Split_comm_type()
564 // on SHARED and then MPI_Comm_split() on the leader ranks.
565
566 // Could also add lowercase etc, but since hostName()
567 // will be consistent within the same node, there is no need.
568 const SHA1Digest myDigest(SHA1(hostName()).digest());
569
570 List<SHA1Digest> digests(worldSize);
571 digests[worldRank] = myDigest;
572
573 // The fixed-length digest allows use of MPI_Allgather.
575 (
576 digests.data_bytes(), // Send/Rev
577 SHA1Digest::size_bytes(), // Num send/recv data per rank
578 UPstream::constWorldComm_
579 );
580
581 // Establish the topology
582 {
583 DynamicList<int> nodeGroup(64);
584 DynamicList<int> nodeLeader(64);
585 DynamicList<SHA1Digest> uniqDigests(64);
586
587 for (int proci = 0; proci < worldSize; ++proci)
588 {
589 const auto& dig = digests[proci];
590
591 if (myDigest == dig)
592 {
593 nodeGroup.push_back(proci);
594 }
595
596 if (!uniqDigests.contains(dig))
597 {
598 // First appearance of host
599 uniqDigests.push_back(dig);
600 nodeLeader.push_back(proci);
601 }
602 }
603
604 procIDs_[commInterNode_] = std::move(nodeLeader);
605 procIDs_[commLocalNode_] = std::move(nodeGroup);
606 }
607 }
608
609
610 // Capture the size (number of nodes) before doing anything further
611 numNodes_ = procIDs_[commInterNode_].size();
612
613 // ~~~~~~~~~
614 // IMPORTANT
615 // ~~~~~~~~~
616 // Always retain knowledge of the inter-node leaders,
617 // even if this process is not on that communicator.
618 // This will help when constructing topology-aware communication.
619
620 // Allocate backend MPI components
621 allocateCommunicatorComponents(constWorldComm_, commInterNode_);
622 allocateCommunicatorComponents(constWorldComm_, commLocalNode_);
623
624 return true;
625}
626
627
629(
630 const label communicator,
631 const bool withComponents
632)
633{
634 // Filter out any placeholders
635 if (communicator < 0)
636 {
637 return;
638 }
639
640 if (debug)
641 {
642 Perr<< "Communicators : Freeing communicator " << communicator
643 << " parent: " << parentComm_[communicator]
644 << " myProcNo: " << myProcNo_[communicator]
645 << endl;
646 }
647
648 if (withComponents && parRun())
649 {
650 freeCommunicatorComponents(communicator);
651 }
652
653 myProcNo_[communicator] = -1;
654 parentComm_[communicator] = -1;
655 //procIDs_[communicator].clear();
656 linearCommunication_[communicator].clear();
657 treeCommunication_[communicator].clear();
658
659 // LIFO push
660 freeComms_.push_back(communicator);
661}
662
663
664int Foam::UPstream::baseProcNo(label comm, int procID)
665{
666 while (UPstream::parent(comm) >= 0 && procID >= 0)
667 {
668 const auto& parentRanks = UPstream::procID(comm);
669 procID = parentRanks[procID];
670 comm = parent(comm);
671 }
672
673 return procID;
674}
675
676
677Foam::label Foam::UPstream::procNo(const label comm, const int baseProcID)
678{
679 const auto& parentRanks = UPstream::procID(comm);
680 label parentComm = UPstream::parent(comm);
681
682 int procID = baseProcID;
683
684 if (parentComm >= 0)
685 {
686 procID = UPstream::procNo(parentComm, baseProcID);
687 }
688
689 return parentRanks.find(procID);
690}
691
692
693Foam::label Foam::UPstream::procNo
694(
695 const label comm,
696 const label currentComm,
697 const int currentProcID
698)
700 label physProcID = UPstream::baseProcNo(currentComm, currentProcID);
701 return UPstream::procNo(comm, physProcID);
702}
703
704
707{
708 // linear = true
709
710 if (linearCommunication_[communicator].empty())
711 {
712 linearCommunication_[communicator].init(communicator, true);
713 }
714 // Probably not needed
715 // else
716 // {
717 // linearCommunication_[communicator].linear(true);
718 // }
719
720 return linearCommunication_[communicator];
721}
722
723
725Foam::UPstream::treeCommunication(int communicator)
726{
727 // linear = false
728
729 if (treeCommunication_[communicator].empty())
730 {
731 treeCommunication_[communicator].init(communicator, false);
732 }
733 // Probably not needed
734 // else
735 // {
736 // treeCommunication_[communicator].linear(false);
737 // }
738
739 return treeCommunication_[communicator];
740}
741
742
744(
745 int communicator,
746 bool linear
747)
748{
749 const auto& comms = UPstream::whichCommunication(communicator, linear);
750
752 {
753 comms.printGraph(Info());
754 }
755}
756
757
758bool Foam::UPstream::usingNodeComms(const int communicator)
759{
760 // Starting point must be "real" world-communicator
761 // ("real" means without any local trickery with worldComm)
762 // Avoid corner cases:
763 // - everthing is on one node
764 // - everthing is on different nodes
765
766 return
767 (
768 parRun_ && (constWorldComm_ == communicator)
769 && (nodeCommsControl_ > 0)
770
771 // More than one node and above defined threshold
772 && (numNodes_ > 1) && (numNodes_ >= nodeCommsMin_)
773 // Some processes do share nodes
774 && (numNodes_ < procIDs_[constWorldComm_].size())
775
776 // Extra paranoid (guard against calling during startup)
777 && (commInterNode_ > constWorldComm_)
778 && (commLocalNode_ > constWorldComm_)
779 );
780}
781
782
784{
785 static std::unique_ptr<List<int>> singleton;
786
787 if (!singleton)
788 {
789 // Extra paranoid (guard against calling during startup)
790 if
791 (
792 (commInterNode_ <= constWorldComm_)
793 || (commInterNode_ >= procIDs_.size())
794 )
795 {
796 return List<int>::null();
797 }
798
799 singleton = std::make_unique<List<int>>();
800 auto& offsets = *singleton;
801
802 const auto& procs = procIDs_[commInterNode_];
803
804 // The procIDs_ are already the offsets, but missing the end offset
805 if (!procs.empty())
806 {
807 const auto count = procs.size();
808
809 offsets.resize(count+1);
810 std::copy_n
811 (
812 procs.begin(),
813 count,
814 offsets.begin()
815 );
816 offsets[count] = UPstream::nProcs(constWorldComm_);
818 }
819
820 return *singleton;
821}
822
823
825{
826 static UPstream::rangeType singleton;
827
828 if (singleton.empty())
829 {
830 // The inter-node offsets [in const world comm] also include a
831 // startup guard
832 const auto& offsets = UPstream::interNode_offsets();
833
834 const auto nodei =
836 (
837 offsets,
838 // My place within const world comm
839 UPstream::myProcNo(constWorldComm_)+1
840 );
841
842 if (nodei >= 0)
843 {
844 singleton.reset
845 (
846 offsets[nodei],
847 offsets[nodei+1] - offsets[nodei]
848 );
849 }
850 }
851
852 return singleton;
853}
854
855
856// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
857
858bool Foam::UPstream::parRun_(false);
859
860bool Foam::UPstream::haveThreads_(false);
861
862bool Foam::UPstream::noInitialCommDup_(false);
863
864int Foam::UPstream::msgType_(1);
865
866
867Foam::wordList Foam::UPstream::allWorlds_(Foam::one{}, "");
868Foam::labelList Foam::UPstream::worldIDs_(Foam::one{}, 0);
869
870
871Foam::DynamicList<int> Foam::UPstream::myProcNo_(16);
872Foam::DynamicList<Foam::List<int>> Foam::UPstream::procIDs_(16);
873
874Foam::DynamicList<Foam::label> Foam::UPstream::parentComm_(16);
875Foam::DynamicList<Foam::label> Foam::UPstream::freeComms_;
876
877Foam::DynamicList<Foam::UPstream::commsStructList>
878Foam::UPstream::linearCommunication_(16);
879
880Foam::DynamicList<Foam::UPstream::commsStructList>
881Foam::UPstream::treeCommunication_(16);
882
883
884int Foam::UPstream::constWorldComm_(0);
885int Foam::UPstream::commInterNode_(-1);
886int Foam::UPstream::commLocalNode_(-1);
887int Foam::UPstream::numNodes_(1);
889Foam::label Foam::UPstream::worldComm(0); // Initially same as constWorldComm_
890Foam::label Foam::UPstream::warnComm(-1);
891
892
893// Predefine world and self communicator slots.
894// These are overwritten in parallel mode (by UPstream::setParRun())
895const int nPredefinedComm = []()
896{
897 // 0: COMM_WORLD : commGlobal(), constWorldComm_, worldComm
899
900 // 1: COMM_SELF
902
903 return Foam::UPstream::nComms();
904}();
905
906
908(
909 Foam::debug::optimisationSwitch("nodeComms", 1)
910);
912(
913 "nodeComms",
914 int,
916);
917
919(
920 Foam::debug::optimisationSwitch("nodeComms.min", 0)
921);
923(
924 "nodeComms.min",
925 int,
927);
928
930(
931 Foam::debug::optimisationSwitch("topoControl", 0)
932);
934(
935 "topoControl",
936 int,
938);
939
941(
942 Foam::debug::optimisationSwitch("floatTransfer", 0)
943);
945(
946 "floatTransfer",
947 bool,
949);
950
952(
953 Foam::debug::optimisationSwitch("nProcsSimpleSum", 0)
954);
956(
957 "nProcsSimpleSum",
958 int,
960);
961
963(
965);
967(
968 "nbx.min",
969 int,
971);
972
974(
975 Foam::debug::optimisationSwitch("nbx.tuning", 0)
976);
978(
979 "nbx.tuning",
980 int,
983
984
986(
987 Foam::debug::optimisationSwitch("nPollProcInterfaces", 0)
988);
990(
991 "nPollProcInterfaces",
992 int,
994);
995
996
998(
999 commsTypeNames.get
1000 (
1001 "commsType",
1003 )
1004);
1005
1006
1008namespace Foam
1009{
1010 //- Registered reader for UPstream::defaultCommsType
1011 class addcommsTypeToOpt
1012 :
1014 {
1015 public:
1016
1017 addcommsTypeToOpt(const char* name)
1018 :
1019 ::Foam::simpleRegIOobject(Foam::debug::addOptimisationObject, name)
1020 {}
1021
1022 virtual ~addcommsTypeToOpt() = default;
1023
1024 virtual void readData(Foam::Istream& is)
1025 {
1026 UPstream::defaultCommsType =
1027 UPstream::commsTypeNames.read(is);
1028 }
1029
1030 virtual void writeData(Foam::Ostream& os) const
1031 {
1032 os << UPstream::commsTypeNames[UPstream::defaultCommsType];
1033 }
1034 };
1035
1036 addcommsTypeToOpt addcommsTypeToOpt_("commsType");
1037}
1041(
1042 Foam::debug::optimisationSwitch("maxCommsSize", 0)
1043);
1045(
1046 "maxCommsSize",
1047 int,
1049);
1050
1051
1053(
1054 Foam::debug::optimisationSwitch("mpiBufferSize", 0)
1056
1057
1058// * * * * * * * * * * * * * * * * Reduction * * * * * * * * * * * * * * * * //
1059
1060int Foam::UPstream::find_first(bool condition, int communicator)
1061{
1063 {
1064 // Ensure we use exact and known data types
1065 typedef int32_t IntType;
1066 constexpr auto dataType = UPstream::dataTypes::type_int32;
1067 constexpr auto opCodeId = UPstream::opCodes::op_min;
1068
1069 auto proci = static_cast<IntType>(UPstream::myProcNo(communicator));
1070 auto nproc = static_cast<IntType>(UPstream::nProcs(communicator));
1071
1072 if (!condition)
1073 {
1074 // The out-of-range value (min reduction)
1075 proci = nproc;
1076 }
1077
1078 UPstream::mpi_allreduce(&proci, 1, dataType, opCodeId, communicator);
1079
1080 // Found: the first rank (min). Not found: -1
1081 return (proci < nproc ? proci : -1);
1082 }
1083 else
1084 {
1085 // Not parallel (serial), so we are rank=0
1086 // return 0 (found) or -1 (not found)
1087 return (condition ? 0 : -1);
1088 }
1089}
1090
1091
1092int Foam::UPstream::find_last(bool condition, int communicator)
1093{
1094 if (UPstream::is_parallel(communicator))
1095 {
1096 // Ensure we use exact and known data types
1097 typedef int32_t IntType;
1098 constexpr auto dataType = UPstream::dataTypes::type_int32;
1099 constexpr auto opCodeId = UPstream::opCodes::op_max;
1100
1101 auto proci = static_cast<IntType>(UPstream::myProcNo(communicator));
1102
1103 if (!condition)
1104 {
1105 // The out-of-range value (max reduction)
1106 proci = -1;
1107 }
1108
1109 UPstream::mpi_allreduce(&proci, 1, dataType, opCodeId, communicator);
1110
1111 // Found: the last rank (max). Not found: -1
1112 return (proci >= 0 ? proci : -1);
1113 }
1114 else
1115 {
1116 // Not parallel (serial), so we are rank=0
1117 // return 0 (found) or -1 (not found)
1118 return (condition ? 0 : -1);
1119 }
1120}
1121
1122
1123// ************************************************************************* //
Useful combination of include files which define Sin, Sout and Serr and the use of IO streams general...
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
#define PrintControl(Ctrl, Name)
const int nPredefinedComm
Definition UPstream.C:888
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
void reset() noexcept
Reset to (0,0).
Definition IntRangeI.H:99
bool empty() const noexcept
True if range is empty (zero-sized).
Definition IntRange.H:198
IntType start() const noexcept
The (inclusive) lower value of the range.
Definition IntRange.H:218
IntType size() const noexcept
The size of the range.
Definition IntRange.H:208
bool contains(IntType value) const noexcept
True if the (global) value is within the range.
Definition IntRangeI.H:143
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
static const List< T > & null() noexcept
Return a null List (reference to a nullObject). Behaves like an empty List.
Definition List.H:138
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
The SHA1 message digest.
Definition SHA1Digest.H:58
static constexpr unsigned size_bytes() noexcept
The number of bytes in digest (20).
Definition SHA1Digest.H:176
Functions to compute SHA1 message digest according to the NIST specification FIPS-180-1.
Definition SHA1.H:57
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
bool contains(const T &val) const
True if the value is contained in the list.
Definition UListI.H:302
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Collection of communication structures.
Definition UPstream.H:363
Wrapper for internally indexed communicator label. Always invokes UPstream::allocateCommunicatorCompo...
Definition UPstream.H:2546
Inter-processor communications stream.
Definition UPstream.H:69
static constexpr int commSelf() noexcept
Communicator within the current rank only.
Definition UPstream.H:1088
static void printCommTree(int communicator, bool linear=false)
Debugging: print the communication tree.
Definition UPstream.C:737
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 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 label parent(int communicator)
The parent communicator.
Definition UPstream.H:1759
commsTypes
Communications types.
Definition UPstream.H:81
@ buffered
"buffered" : (MPI_Bsend, MPI_Recv)
Definition UPstream.H:82
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 const Enum< commsTypes > commsTypeNames
Enumerated names for the communication types.
Definition UPstream.H:92
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 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 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 int topologyControl_
Selection of topology-aware routines as a bitmask combination of the topoControls enumerations.
Definition UPstream.H:1009
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 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 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 bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static int nProcsNonblockingExchange
Number of processors to change to nonBlocking consensual exchange (NBX). Ignored for zero or negative...
Definition UPstream.H:1035
static const commsStructList & linearCommunication(int communicator)
Linear communication schedule (special case) for given communicator.
Definition UPstream.C:699
static bool is_parallel(const label communicator=worldComm)
True if parallel algorithm or exchange is required.
Definition UPstream.H:1743
static int tuning_NBX_
Tuning parameters for non-blocking exchange (NBX).
Definition UPstream.H:1055
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
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition UPstream.H:1069
static label nComms() noexcept
Number of currently defined communicators.
Definition UPstream.H:1132
IntRange< int > rangeType
Int ranges are used for MPI ranks (processes).
Definition UPstream.H:75
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 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 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 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 void mpiAllGather(Type *allData, int count, const int communicator=UPstream::worldComm)
Gather/scatter identically-sized data.
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 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 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 int baseProcNo(label comm, int procID)
Return physical processor number (i.e. processor number in worldComm) given communicator and processo...
Definition UPstream.C:657
static void printNodeCommsControl(Ostream &os)
Report the node-communication settings.
Definition UPstream.C:56
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
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
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
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
Abstract base class for registered object with I/O. Used in debug symbol registration.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition BitOps.H:73
Namespace for handling debugging switches.
Definition debug.C:45
dictionary & optimisationSwitches()
The OptimisationSwitches sub-dictionary in the central controlDict(s).
Definition debug.C:216
int optimisationSwitch(const char *name, const int deflt=0)
Lookup optimisation switch or add default value.
Definition debug.C:234
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
prefixOSstream Perr
OSstream wrapped stderr (std::cerr) with parallel prefix.
messageStream Info
Information stream (stdout output on master, null elsewhere).
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
string hostName()
Return the system's host name, as per hostname(1).
Definition POSIX.C:373
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
errorManip< error > abort(error &err)
Definition errorManip.H:139
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define registerOptSwitch(Name, Type, SwitchVar)
#define FOAM_UNLIKELY(cond)
Definition stdFoam.H:64
const bool writeData(pdfDictionary.get< bool >("writeData"))