Loading...
Searching...
No Matches
syncToolsTemplates.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
27\*---------------------------------------------------------------------------*/
28
29#include "syncTools.H"
30#include "polyMesh.H"
31#include "processorPolyPatch.H"
32#include "cyclicPolyPatch.H"
33#include "globalMeshData.H"
34#include "transformList.H"
35
36// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37
38template<class T, class CombineOp>
39void Foam::syncTools::combine
40(
41 Map<T>& pointValues,
42 const CombineOp& cop,
43 const label index,
44 const T& val
45)
46{
47 auto iter = pointValues.find(index);
48
49 if (iter.good())
50 {
51 cop(iter.val(), val);
52 }
53 else
54 {
55 pointValues.insert(index, val);
56 }
57}
58
59
60template<class T, class CombineOp>
61void Foam::syncTools::combine
62(
63 EdgeMap<T>& edgeValues,
64 const CombineOp& cop,
65 const edge& index,
66 const T& val
67)
68{
69 auto iter = edgeValues.find(index);
70
71 if (iter.good())
72 {
73 cop(iter.val(), val);
74 }
75 else
76 {
77 edgeValues.insert(index, val);
78 }
79}
80
81
82// * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
83
84template<class T, class CombineOp, class TransformOp>
86(
87 const polyMesh& mesh,
88 Map<T>& pointValues, // from mesh point label to value
89 const CombineOp& cop,
90 const TransformOp& top
91)
92{
93 const polyBoundaryMesh& patches = mesh.boundaryMesh();
94
95 // Make sure we use unique message tag
96 const int oldTag = UPstream::incrMsgType();
97
98 // Synchronize multiple shared points.
99 const globalMeshData& pd = mesh.globalData();
100
101 // Values on shared points. Keyed on global shared index.
102 Map<T> sharedPointValues;
103
104 if (pd.nGlobalPoints() > 0)
105 {
106 // meshPoint per local index
107 const labelList& sharedPtLabels = pd.sharedPointLabels();
108
109 // global shared index per local index
110 const labelList& sharedPtAddr = pd.sharedPointAddr();
111
112 sharedPointValues.reserve(sharedPtAddr.size());
113
114 // Fill my entries in the shared points
115 forAll(sharedPtLabels, i)
116 {
117 const auto fnd = pointValues.cfind(sharedPtLabels[i]);
118
119 if (fnd.good())
120 {
121 combine
122 (
123 sharedPointValues,
124 cop,
125 sharedPtAddr[i], // index
126 fnd.val() // value
127 );
128 }
129 }
130 }
131
132
133 if (UPstream::parRun())
134 {
135 // Presize according to number of processor patches
136 // (global topology information may not yet be available...)
138 PstreamBuffers pBufs;
139
140 // Reduce communication by only sending non-zero data,
141 // but with multiply-connected processor/processor
142 // (eg, processorCyclic) also need to send zero information
143 // to keep things synchronised
144
145 // Initialize for registerSend() bookkeeping
146 pBufs.initRegisterSend();
147
148
149 // Sample and send.
150
151 for (const polyPatch& pp : patches)
152 {
153 const auto* ppp = isA<processorPolyPatch>(pp);
154
155 if (ppp && pp.nPoints())
156 {
157 const auto& procPatch = *ppp;
158 const label nbrProci = procPatch.neighbProcNo();
159
160 // Get data per patchPoint in neighbouring point numbers.
161 const labelList& meshPts = procPatch.meshPoints();
162 const labelList& nbrPts = procPatch.neighbPoints();
163
164 // Extract local values. Create map from nbrPoint to value.
165 // Note: how small initial size?
166 Map<T> patchInfo(meshPts.size() / 20);
167
168 forAll(meshPts, i)
169 {
170 const auto iter = pointValues.cfind(meshPts[i]);
171
172 if (iter.good())
173 {
174 patchInfo.insert(nbrPts[i], iter.val());
175 }
176 }
177
178 // Neighbour connectivity
179 neighbProcs.push_uniq(nbrProci);
180
181 // Send to neighbour
182 {
183 UOPstream toNbr(nbrProci, pBufs);
184 toNbr << patchInfo;
185
186 // Record if send is required (data are non-zero)
187 pBufs.registerSend(nbrProci, !patchInfo.empty());
188 }
189 }
190 }
191
192 // Limit exchange to involved procs.
193 // - automatically discards unnecessary (unregistered) sends
194 pBufs.finishedNeighbourSends(neighbProcs);
195
196
197 // Receive and combine.
198 for (const polyPatch& pp : patches)
199 {
200 const auto* ppp = isA<processorPolyPatch>(pp);
201
202 if (ppp && pp.nPoints())
203 {
204 const auto& procPatch = *ppp;
205 const label nbrProci = procPatch.neighbProcNo();
206
207 if (!pBufs.recvDataCount(nbrProci))
208 {
209 // Nothing to receive
210 continue;
211 }
212
213 Map<T> nbrPatchInfo;
214 {
215 UIPstream fromNbr(nbrProci, pBufs);
216 fromNbr >> nbrPatchInfo;
217 }
218
219 // Transform
220 top(procPatch, nbrPatchInfo);
221
222 const labelList& meshPts = procPatch.meshPoints();
223
224 // Only update those values which come from neighbour
225 forAllConstIters(nbrPatchInfo, nbrIter)
226 {
227 combine
228 (
229 pointValues,
230 cop,
231 meshPts[nbrIter.key()],
232 nbrIter.val()
233 );
234 }
235 }
236 }
237 }
238
239 // Do the cyclics.
240 for (const polyPatch& pp : patches)
241 {
243
244 if (cpp && cpp->owner())
245 {
246 // Owner does all.
247
248 const cyclicPolyPatch& cycPatch = *cpp;
249 const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
250
251 const edgeList& coupledPoints = cycPatch.coupledPoints();
252 const labelList& meshPtsA = cycPatch.meshPoints();
253 const labelList& meshPtsB = nbrPatch.meshPoints();
254
255 // Extract local values. Create map from coupled-edge to value.
256 Map<T> half0Values(meshPtsA.size() / 20);
257 Map<T> half1Values(half0Values.size());
258
259 forAll(coupledPoints, i)
260 {
261 const edge& e = coupledPoints[i];
262
263 const auto point0Fnd = pointValues.cfind(meshPtsA[e[0]]);
264
265 if (point0Fnd.good())
266 {
267 half0Values.insert(i, *point0Fnd);
268 }
269
270 const auto point1Fnd = pointValues.cfind(meshPtsB[e[1]]);
271
272 if (point1Fnd.good())
273 {
274 half1Values.insert(i, *point1Fnd);
275 }
276 }
277
278 // Transform to receiving side
279 top(cycPatch, half1Values);
280 top(nbrPatch, half0Values);
281
282 forAll(coupledPoints, i)
283 {
284 const edge& e = coupledPoints[i];
285
286 const auto half0Fnd = half0Values.cfind(i);
287
288 if (half0Fnd.good())
289 {
290 combine
291 (
292 pointValues,
293 cop,
294 meshPtsB[e[1]],
295 *half0Fnd
296 );
297 }
298
299 const auto half1Fnd = half1Values.cfind(i);
300
301 if (half1Fnd.good())
302 {
303 combine
304 (
305 pointValues,
306 cop,
307 meshPtsA[e[0]],
308 *half1Fnd
309 );
310 }
311 }
312 }
313 }
314
315 // Synchronize multiple shared points.
316 if (pd.nGlobalPoints() > 0)
317 {
318 // meshPoint per local index
319 const labelList& sharedPtLabels = pd.sharedPointLabels();
320 // global shared index per local index
321 const labelList& sharedPtAddr = pd.sharedPointAddr();
322
323 // Reduce on master.
324
325 if (UPstream::parRun())
326 {
327 if (UPstream::master())
328 {
329 // Receive the edges using shared points from other procs
330 for (const int proci : UPstream::subProcs())
331 {
332 Map<T> nbrValues;
333 IPstream::recv(nbrValues, proci);
334
335 // Merge neighbouring values with my values
336 forAllConstIters(nbrValues, iter)
337 {
338 combine
339 (
340 sharedPointValues,
341 cop,
342 iter.key(), // edge
343 iter.val() // value
344 );
345 }
346 }
347 }
348 else
349 {
350 // Send to master
351 OPstream::send(sharedPointValues, UPstream::masterNo());
352 }
353
354 // Broadcast: send merged values to all
355 Pstream::broadcast(sharedPointValues);
356 }
357
358 // Merge sharedPointValues (keyed on sharedPointAddr) into
359 // pointValues (keyed on mesh points).
360
361 // Map from global shared index to meshpoint
362 Map<label> sharedToMeshPoint(sharedPtAddr, sharedPtLabels);
363
364 forAllConstIters(sharedToMeshPoint, iter)
365 {
366 // Do I have a value for my shared point
367 const auto sharedFnd = sharedPointValues.cfind(iter.key());
368
369 if (sharedFnd.good())
370 {
371 pointValues.set(iter.val(), sharedFnd.val());
372 }
373 }
374 }
376 // Reset tag
377 UPstream::msgType(oldTag);
378}
379
380
381template<class T, class CombineOp, class TransformOp>
383(
384 const polyMesh& mesh,
385 EdgeMap<T>& edgeValues,
386 const CombineOp& cop,
387 const TransformOp& top
388)
389{
390 const polyBoundaryMesh& patches = mesh.boundaryMesh();
391
392
393 // Do synchronisation without constructing globalEdge addressing
394 // (since this constructs mesh edge addressing)
395
396 // Make sure we use unique message tag
397 const int oldTag = UPstream::incrMsgType();
398
399 // Swap proc patch info
400 // ~~~~~~~~~~~~~~~~~~~~
401
402 if (UPstream::parRun())
403 {
404 // Presize according to number of processor patches
405 // (global topology information may not yet be available...)
406 DynamicList<label> neighbProcs(patches.nProcessorPatches());
407 PstreamBuffers pBufs;
408
409 // Reduce communication by only sending non-zero data,
410 // but with multiply-connected processor/processor
411 // (eg, processorCyclic) also need to send zero information
412 // to keep things synchronised
413
414 // Initialize for registerSend() bookkeeping
415 pBufs.initRegisterSend();
416
417
418 // Sample and send.
419
420 for (const polyPatch& pp : patches)
421 {
422 const auto* ppp = isA<processorPolyPatch>(pp);
423
424 if (ppp && pp.nEdges())
425 {
426 const auto& procPatch = *ppp;
427 const label nbrProci = procPatch.neighbProcNo();
428
429 // Get data per patch edge in neighbouring edge.
430 const edgeList& edges = procPatch.edges();
431 const labelList& meshPts = procPatch.meshPoints();
432 const labelList& nbrPts = procPatch.neighbPoints();
433
434 EdgeMap<T> patchInfo(edges.size() / 20);
435
436 for (const edge& e : edges)
437 {
438 const edge meshEdge(meshPts, e);
439
440 const auto iter = edgeValues.cfind(meshEdge);
441
442 if (iter.good())
443 {
444 const edge nbrEdge(nbrPts, e);
445 patchInfo.insert(nbrEdge, iter.val());
446 }
447 }
448
449
450 // Neighbour connectivity
451 neighbProcs.push_uniq(nbrProci);
452
453 // Send to neighbour
454 {
455 UOPstream toNbr(nbrProci, pBufs);
456 toNbr << patchInfo;
457
458 // Record if send is required (data are non-zero)
459 pBufs.registerSend(nbrProci, !patchInfo.empty());
460 }
461 }
462 }
463
464 // Limit exchange to involved procs
465 // - automatically discards unnecessary (unregistered) sends
466 pBufs.finishedNeighbourSends(neighbProcs);
467
468
469 // Receive and combine.
470 for (const polyPatch& pp : patches)
471 {
472 const auto* ppp = isA<processorPolyPatch>(pp);
473
474 if (ppp && pp.nEdges())
475 {
476 const auto& procPatch = *ppp;
477 const label nbrProci = procPatch.neighbProcNo();
478
479 if (!pBufs.recvDataCount(nbrProci))
480 {
481 // Nothing to receive
482 continue;
483 }
484
485 EdgeMap<T> nbrPatchInfo;
486 {
487 UIPstream fromNbr(nbrProci, pBufs);
488 fromNbr >> nbrPatchInfo;
489 }
490
491 // Apply transform to convert to this side properties.
492 top(procPatch, nbrPatchInfo);
493
494
495 // Only update those values which come from neighbour
496 const labelList& meshPts = procPatch.meshPoints();
497
498 forAllConstIters(nbrPatchInfo, nbrIter)
499 {
500 const edge& e = nbrIter.key();
501 const edge meshEdge(meshPts, e);
502
503 combine
504 (
505 edgeValues,
506 cop,
507 meshEdge, // edge
508 nbrIter.val() // value
509 );
510 }
511 }
512 }
513 }
514
515
516 // Swap cyclic info
517 // ~~~~~~~~~~~~~~~~
518
519 for (const polyPatch& pp : patches)
520 {
522
523 if (cpp && cpp->owner())
524 {
525 // Owner does all.
526
527 const cyclicPolyPatch& cycPatch = *cpp;
528 const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
529
530 const edgeList& coupledEdges = cycPatch.coupledEdges();
531
532 const labelList& meshPtsA = cycPatch.meshPoints();
533 const edgeList& edgesA = cycPatch.edges();
534
535 const labelList& meshPtsB = nbrPatch.meshPoints();
536 const edgeList& edgesB = nbrPatch.edges();
537
538 // Extract local values. Create map from edge to value.
539 Map<T> half0Values(edgesA.size() / 20);
540 Map<T> half1Values(half0Values.size());
541
542 forAll(coupledEdges, edgei)
543 {
544 const edge& twoEdges = coupledEdges[edgei];
545
546 {
547 const edge& e0 = edgesA[twoEdges[0]];
548 const edge meshEdge0(meshPtsA, e0);
549
550 const auto iter = edgeValues.cfind(meshEdge0);
551
552 if (iter.good())
553 {
554 half0Values.insert(edgei, iter.val());
555 }
556 }
557 {
558 const edge& e1 = edgesB[twoEdges[1]];
559 const edge meshEdge1(meshPtsB, e1);
560
561 const auto iter = edgeValues.cfind(meshEdge1);
562
563 if (iter.good())
564 {
565 half1Values.insert(edgei, iter.val());
566 }
567 }
568 }
569
570 // Transform to this side
571 top(cycPatch, half1Values);
572 top(nbrPatch, half0Values);
573
574
575 // Extract and combine information
576
577 forAll(coupledEdges, edgei)
578 {
579 const edge& twoEdges = coupledEdges[edgei];
580
581 const auto half1Fnd = half1Values.cfind(edgei);
582
583 if (half1Fnd.good())
584 {
585 const edge& e0 = edgesA[twoEdges[0]];
586 const edge meshEdge0(meshPtsA, e0);
587
588 combine
589 (
590 edgeValues,
591 cop,
592 meshEdge0, // edge
593 *half1Fnd // value
594 );
595 }
596
597 const auto half0Fnd = half0Values.cfind(edgei);
598
599 if (half0Fnd.good())
600 {
601 const edge& e1 = edgesB[twoEdges[1]];
602 const edge meshEdge1(meshPtsB, e1);
603
604 combine
605 (
606 edgeValues,
607 cop,
608 meshEdge1, // edge
609 *half0Fnd // value
610 );
611 }
612 }
613 }
614 }
615
616 // Synchronize multiple shared points.
617 // Problem is that we don't want to construct shared edges so basically
618 // we do here like globalMeshData but then using sparse edge representation
619 // (EdgeMap instead of mesh.edges())
620
621 const globalMeshData& pd = mesh.globalData();
622 const labelList& sharedPtAddr = pd.sharedPointAddr();
623 const labelList& sharedPtLabels = pd.sharedPointLabels();
624
625 // 1. Create map from meshPoint to globalShared index.
626 Map<label> meshToShared(sharedPtLabels, sharedPtAddr);
627
628 // Values on shared points. From two sharedPtAddr indices to a value.
629 EdgeMap<T> sharedEdgeValues(meshToShared.size());
630
631 // From shared edge to mesh edge. Used for merging later on.
632 EdgeMap<edge> potentialSharedEdge(meshToShared.size());
633
634 // 2. Find any edges using two global shared points. These will always be
635 // on the outside of the mesh. (though might not be on coupled patch
636 // if is single edge and not on coupled face)
637 // Store value (if any) on sharedEdgeValues
638 for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
639 {
640 const face& f = mesh.faces()[facei];
641
642 forAll(f, fp)
643 {
644 const label v0 = f[fp];
645 const label v1 = f[f.fcIndex(fp)];
646
647 const auto v0Fnd = meshToShared.cfind(v0);
648
649 if (v0Fnd.good())
650 {
651 const auto v1Fnd = meshToShared.cfind(v1);
652
653 if (v1Fnd.good())
654 {
655 const edge meshEdge(v0, v1);
656
657 // edge in shared point labels
658 const edge sharedEdge(*v0Fnd, *v1Fnd);
659
660 // Store mesh edge as a potential shared edge.
661 potentialSharedEdge.insert(sharedEdge, meshEdge);
662
663 const auto edgeFnd = edgeValues.cfind(meshEdge);
664
665 if (edgeFnd.good())
666 {
667 // edge exists in edgeValues. See if already in map
668 // (since on same processor, e.g. cyclic)
669 combine
670 (
671 sharedEdgeValues,
672 cop,
673 sharedEdge, // edge
674 *edgeFnd // value
675 );
676 }
677 }
678 }
679 }
680 }
681
682
683 // Now sharedEdgeValues will contain per potential sharedEdge the value.
684 // (potential since an edge having two shared points is not necessary a
685 // shared edge).
686 // Reduce this on the master.
687
688 if (UPstream::parRun())
689 {
690 if (UPstream::master())
691 {
692 // Receive the edges using shared points from other procs
693 for (const int proci : UPstream::subProcs())
694 {
695 EdgeMap<T> nbrValues;
696 IPstream::recv(nbrValues, proci);
697
698 // Merge neighbouring values with my values
699 forAllConstIters(nbrValues, iter)
700 {
701 combine
702 (
703 sharedEdgeValues,
704 cop,
705 iter.key(), // edge
706 iter.val() // value
707 );
708 }
709 }
710 }
711 else
712 {
713 // Send to master
714 OPstream::send(sharedEdgeValues, UPstream::masterNo());
715 }
716
717 // Broadcast: send merged values to all
718 Pstream::broadcast(sharedEdgeValues);
719 }
720
721
722 // Merge sharedEdgeValues (keyed on sharedPointAddr) into edgeValues
723 // (keyed on mesh points).
724
725 // Loop over all my shared edges.
726 forAllConstIters(potentialSharedEdge, iter)
727 {
728 const edge& sharedEdge = iter.key();
729 const edge& meshEdge = iter.val();
730
731 // Do I have a value for the shared edge?
732 const auto sharedFnd = sharedEdgeValues.cfind(sharedEdge);
733
734 if (sharedFnd.good())
735 {
736 combine
737 (
738 edgeValues,
739 cop,
740 meshEdge, // edge
741 *sharedFnd // value
742 );
743 }
744 }
746 // Reset tag
747 UPstream::msgType(oldTag);
748}
749
750
751template<class T, class CombineOp, class TransformOp>
753(
754 const polyMesh& mesh,
755 List<T>& pointValues,
756 const CombineOp& cop,
757 const T& nullValue,
758 const TransformOp& top
759)
760{
761 if (pointValues.size() != mesh.nPoints())
762 {
764 << "Number of values " << pointValues.size()
765 << " != number of points " << mesh.nPoints() << nl
766 << abort(FatalError);
768
769 mesh.globalData().syncPointData(pointValues, cop, top);
770}
771
772
773template<class T, class CombineOp, class TransformOp>
775(
776 const polyMesh& mesh,
777 const labelUList& meshPoints,
778 List<T>& pointValues,
779 const CombineOp& cop,
780 const T& nullValue,
781 const TransformOp& top
782)
783{
784 if (pointValues.size() != meshPoints.size())
785 {
787 << "Number of values " << pointValues.size()
788 << " != number of meshPoints " << meshPoints.size() << nl
789 << abort(FatalError);
790 }
791 const globalMeshData& gd = mesh.globalData();
792 const indirectPrimitivePatch& cpp = gd.coupledPatch();
793 const Map<label>& mpm = cpp.meshPointMap();
794
795 List<T> cppFld(cpp.nPoints(), nullValue);
796
797 forAll(meshPoints, i)
798 {
799 const auto iter = mpm.cfind(meshPoints[i]);
800
801 if (iter.good())
802 {
803 cppFld[iter.val()] = pointValues[i];
804 }
805 }
806
808 (
809 cppFld,
810 gd.globalPointSlaves(),
811 gd.globalPointTransformedSlaves(),
812 gd.globalPointSlavesMap(),
813 gd.globalTransforms(),
814 cop,
815 top
816 );
817
818 forAll(meshPoints, i)
819 {
820 const auto iter = mpm.cfind(meshPoints[i]);
821
822 if (iter.good())
823 {
824 pointValues[i] = cppFld[iter.val()];
825 }
826 }
827}
828
829
830template<class T, class CombineOp, class TransformOp, class FlipOp>
832(
833 const polyMesh& mesh,
834 List<T>& edgeValues,
835 const CombineOp& cop,
836 const T& nullValue,
837 const TransformOp& top,
838 const FlipOp& fop
839)
840{
841 if (edgeValues.size() != mesh.nEdges())
842 {
844 << "Number of values " << edgeValues.size()
845 << " != number of edges " << mesh.nEdges() << nl
846 << abort(FatalError);
847 }
848
849 const edgeList& edges = mesh.edges();
850 const globalMeshData& gd = mesh.globalData();
851 const labelList& meshEdges = gd.coupledPatchMeshEdges();
852 const indirectPrimitivePatch& cpp = gd.coupledPatch();
853 const edgeList& cppEdges = cpp.edges();
854 const labelList& mp = cpp.meshPoints();
855 const globalIndexAndTransform& git = gd.globalTransforms();
856 const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
857 const bitSet& orientation = gd.globalEdgeOrientation();
858
859 List<T> cppFld(meshEdges.size());
860 forAll(meshEdges, i)
861 {
862 const edge& cppE = cppEdges[i];
863 const label meshEdgei = meshEdges[i];
864 const edge& meshE = edges[meshEdgei];
865
866 // 1. is cpp edge oriented as mesh edge
867 // 2. is cpp edge oriented same as master edge
868
869 const int dir = edge::compare(meshE, edge(mp, cppE));
870 if (dir == 0)
871 {
872 FatalErrorInFunction<< "Problem:"
873 << " mesh edge:" << meshE.line(mesh.points())
874 << " coupled edge:" << cppE.line(cpp.localPoints())
875 << exit(FatalError);
876 }
877
878 const bool sameOrientation = ((dir == 1) == orientation[i]);
879
880 if (sameOrientation)
881 {
882 cppFld[i] = edgeValues[meshEdgei];
883 }
884 else
885 {
886 cppFld[i] = fop(edgeValues[meshEdgei]);
887 }
888 }
889
890
892 (
893 cppFld,
894 gd.globalEdgeSlaves(),
895 gd.globalEdgeTransformedSlaves(),
896 edgeMap,
897 git,
898 cop,
899 top
900 );
901
902 // Extract back onto mesh
903 forAll(meshEdges, i)
904 {
905 const edge& cppE = cppEdges[i];
906 const label meshEdgei = meshEdges[i];
907 const edge& meshE = edges[meshEdgei];
908
909 // 1. is cpp edge oriented as mesh edge
910 // 2. is cpp edge oriented same as master edge
911
912 const int dir = edge::compare(meshE, edge(mp, cppE));
913 const bool sameOrientation = ((dir == 1) == orientation[i]);
914
915 if (sameOrientation)
916 {
917 edgeValues[meshEdges[i]] = cppFld[i];
918 }
919 else
920 {
921 edgeValues[meshEdges[i]] = fop(cppFld[i]);
922 }
923 }
924}
925
926
927template<class T, class CombineOp, class TransformOp, class FlipOp>
929(
930 const polyMesh& mesh,
931 const labelUList& meshEdges,
932 List<T>& edgeValues,
933 const CombineOp& cop,
934 const T& nullValue,
935 const TransformOp& top,
936 const FlipOp& fop
937)
938{
939 if (edgeValues.size() != meshEdges.size())
940 {
942 << "Number of values " << edgeValues.size()
943 << " != number of meshEdges " << meshEdges.size() << nl
944 << abort(FatalError);
945 }
946 const edgeList& edges = mesh.edges();
947 const globalMeshData& gd = mesh.globalData();
948 const indirectPrimitivePatch& cpp = gd.coupledPatch();
949 const edgeList& cppEdges = cpp.edges();
950 const labelList& mp = cpp.meshPoints();
951 const Map<label>& mpm = gd.coupledPatchMeshEdgeMap();
952 const bitSet& orientation = gd.globalEdgeOrientation();
953
954 List<T> cppFld(cpp.nEdges(), nullValue);
955
956 forAll(meshEdges, i)
957 {
958 const label meshEdgei = meshEdges[i];
959 const auto iter = mpm.cfind(meshEdgei);
960 if (iter.good())
961 {
962 const label cppEdgei = iter.val();
963 const edge& cppE = cppEdges[cppEdgei];
964 const edge& meshE = edges[meshEdgei];
965
966 // 1. is cpp edge oriented as mesh edge
967 // 2. is cpp edge oriented same as master edge
968
969 const int dir = edge::compare(meshE, edge(mp, cppE));
970 if (dir == 0)
971 {
972 FatalErrorInFunction<< "Problem:"
973 << " mesh edge:" << meshE.line(mesh.points())
974 << " coupled edge:" << cppE.line(cpp.localPoints())
975 << exit(FatalError);
976 }
977
978 const bool sameOrientation = ((dir == 1) == orientation[i]);
979
980 if (sameOrientation)
981 {
982 cppFld[cppEdgei] = edgeValues[i];
983 }
984 else
985 {
986 cppFld[cppEdgei] = fop(edgeValues[i]);
987 }
988 }
989 }
990
992 (
993 cppFld,
994 gd.globalEdgeSlaves(),
995 gd.globalEdgeTransformedSlaves(),
996 gd.globalEdgeSlavesMap(),
997 gd.globalTransforms(),
998 cop,
999 top
1000 );
1001
1002 forAll(meshEdges, i)
1003 {
1004 label meshEdgei = meshEdges[i];
1005 const auto iter = mpm.cfind(meshEdgei);
1006 if (iter.good())
1007 {
1008 label cppEdgei = iter.val();
1009 const edge& cppE = cppEdges[cppEdgei];
1010 const edge& meshE = edges[meshEdgei];
1011
1012 const int dir = edge::compare(meshE, edge(mp, cppE));
1013 const bool sameOrientation = ((dir == 1) == orientation[i]);
1014
1015 if (sameOrientation)
1016 {
1017 edgeValues[i] = cppFld[cppEdgei];
1018 }
1019 else
1020 {
1021 edgeValues[i] = fop(cppFld[cppEdgei]);
1023 }
1024 }
1025}
1026
1027
1028template<class T, class CombineOp, class TransformOp>
1030(
1031 const polyMesh& mesh,
1032 UList<T>& faceValues,
1033 const CombineOp& cop,
1034 const TransformOp& top,
1035 const bool parRun
1036)
1037{
1038 // Offset (global to local) for start of boundaries
1039 const label boundaryOffset = mesh.nInternalFaces();
1040
1041 if (faceValues.size() != mesh.nBoundaryFaces())
1042 {
1044 << "Number of values " << faceValues.size()
1045 << " != number of boundary faces " << mesh.nBoundaryFaces() << nl
1046 << abort(FatalError);
1047 }
1048
1050
1051 // Allocate unique tag for all comms
1052 const int oldTag = UPstream::incrMsgType();
1053
1054 if (parRun && UPstream::parRun())
1055 {
1056 // Avoid mesh.globalData() - possible race condition
1057
1058 if
1059 (
1062 )
1063 {
1064 const label startRequest = UPstream::nRequests();
1065
1066 // Receive buffer
1067 List<T> receivedValues(mesh.nBoundaryFaces());
1068
1069 // Set up reads
1070 for (const polyPatch& pp : patches)
1071 {
1072 const auto* ppp = isA<processorPolyPatch>(pp);
1073
1074 if (ppp && pp.size())
1075 {
1076 const auto& procPatch = *ppp;
1077
1079 (
1080 receivedValues,
1081 pp.size(),
1082 pp.start()-boundaryOffset
1083 );
1084
1086 (
1088 procPatch.neighbProcNo(),
1089 fld
1090 );
1091 }
1092 }
1093
1094 // Set up writes
1095 for (const polyPatch& pp : patches)
1096 {
1097 const auto* ppp = isA<processorPolyPatch>(pp);
1098
1099 if (ppp && pp.size())
1100 {
1101 const auto& procPatch = *ppp;
1102
1103 const SubList<T> fld
1104 (
1105 faceValues,
1106 pp.size(),
1107 pp.start()-boundaryOffset
1108 );
1109
1111 (
1113 procPatch.neighbProcNo(),
1114 fld
1115 );
1116 }
1117 }
1118
1119 // Wait for all comms to finish
1120 UPstream::waitRequests(startRequest);
1121
1122 // Combine with existing data
1123 for (const polyPatch& pp : patches)
1124 {
1125 const auto* ppp = isA<processorPolyPatch>(pp);
1126
1127 if (ppp && pp.size())
1128 {
1129 const auto& procPatch = *ppp;
1130
1131 SubList<T> recvFld
1132 (
1133 receivedValues,
1134 pp.size(),
1135 pp.start()-boundaryOffset
1136 );
1137
1138 top(procPatch, recvFld);
1139
1140 SubList<T> patchValues
1141 (
1142 faceValues,
1143 pp.size(),
1144 pp.start()-boundaryOffset
1145 );
1146
1147 forAll(patchValues, i)
1148 {
1149 cop(patchValues[i], recvFld[i]);
1150 }
1151 }
1152 }
1153 }
1154 else
1155 {
1156 DynamicList<label> neighbProcs;
1157 PstreamBuffers pBufs;
1158
1159 // Send
1160 for (const polyPatch& pp : patches)
1161 {
1162 const auto* ppp = isA<processorPolyPatch>(pp);
1163
1164 if (ppp && pp.size())
1165 {
1166 const auto& procPatch = *ppp;
1167 const label nbrProci = procPatch.neighbProcNo();
1168
1169 // Neighbour connectivity
1170 neighbProcs.push_uniq(nbrProci);
1171
1172 const SubList<T> fld
1173 (
1174 faceValues,
1175 pp.size(),
1176 pp.start()-boundaryOffset
1177 );
1178
1179 UOPstream toNbr(nbrProci, pBufs);
1180 toNbr << fld;
1181 }
1182 }
1183
1184 // Limit exchange to involved procs
1185 pBufs.finishedNeighbourSends(neighbProcs);
1186
1187
1188 // Receive and combine.
1189 for (const polyPatch& pp : patches)
1190 {
1191 const auto* ppp = isA<processorPolyPatch>(pp);
1192
1193 if (ppp && pp.size())
1194 {
1195 const auto& procPatch = *ppp;
1196 const label nbrProci = procPatch.neighbProcNo();
1197
1198 List<T> recvFld;
1199 {
1200 UIPstream fromNbr(nbrProci, pBufs);
1201 fromNbr >> recvFld;
1202 }
1203
1204 top(procPatch, recvFld);
1205
1206 SubList<T> patchValues
1207 (
1208 faceValues,
1209 pp.size(),
1210 pp.start()-boundaryOffset
1211 );
1212
1213 forAll(patchValues, i)
1214 {
1215 cop(patchValues[i], recvFld[i]);
1216 }
1217 }
1218 }
1219 }
1220 }
1221
1222 // Do the cyclics.
1223 for (const polyPatch& pp : patches)
1224 {
1226
1227 if (cpp && cpp->owner())
1228 {
1229 // Owner does all.
1230
1231 const cyclicPolyPatch& cycPatch = *cpp;
1232 const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
1233 const label patchSize = cycPatch.size();
1234
1235 SubList<T> ownPatchValues
1236 (
1237 faceValues,
1238 patchSize,
1239 cycPatch.start()-boundaryOffset
1240 );
1241
1242 SubList<T> nbrPatchValues
1243 (
1244 faceValues,
1245 patchSize,
1246 nbrPatch.start()-boundaryOffset
1247 );
1248
1249 // Transform (copy of) data on both sides
1250 List<T> ownVals(ownPatchValues);
1251 top(nbrPatch, ownVals);
1252
1253 List<T> nbrVals(nbrPatchValues);
1254 top(cycPatch, nbrVals);
1255
1256 forAll(ownPatchValues, i)
1257 {
1258 cop(ownPatchValues[i], nbrVals[i]);
1259 }
1260
1261 forAll(nbrPatchValues, i)
1262 {
1263 cop(nbrPatchValues[i], ownVals[i]);
1264 }
1265 }
1266 }
1267
1268 // Reset tag
1270}
1271
1272
1273// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1274
1275template<unsigned Width, class CombineOp>
1277(
1278 const polyMesh& mesh,
1279 const bool isBoundaryOnly,
1280 PackedList<Width>& faceValues,
1281 const CombineOp& cop,
1282 const bool parRun
1283)
1284{
1285 // Offset (global to local) for start of boundaries
1286 const label boundaryOffset = (isBoundaryOnly ? mesh.nInternalFaces() : 0);
1287
1288 // Check size
1289 if (faceValues.size() != (mesh.nFaces() - boundaryOffset))
1290 {
1292 << "Number of values " << faceValues.size()
1293 << " != number of "
1294 << (isBoundaryOnly ? "boundary" : "mesh") << " faces "
1295 << (mesh.nFaces() - boundaryOffset) << nl
1296 << abort(FatalError);
1297 }
1298
1299 const polyBoundaryMesh& patches = mesh.boundaryMesh();
1300
1301 // Allocate unique tag for all comms
1302 const int oldTag = UPstream::incrMsgType();
1303
1304 if (parRun && UPstream::parRun())
1305 {
1306 const label startRequest = UPstream::nRequests();
1307
1308 // Receive buffers
1309 PtrList<PackedList<Width>> recvBufs(patches.size());
1310
1311 // Set up reads
1312 for (const polyPatch& pp : patches)
1313 {
1314 const auto* ppp = isA<processorPolyPatch>(pp);
1315
1316 if (ppp && pp.size())
1317 {
1318 const auto& procPatch = *ppp;
1319 const label patchi = pp.index();
1320
1321 auto& recvbuf = recvBufs.emplace_set(patchi, pp.size());
1322
1324 (
1326 procPatch.neighbProcNo(),
1327 recvbuf.data_bytes(),
1328 recvbuf.size_bytes()
1329 );
1330 }
1331 }
1332
1333 // Send buffers
1335
1336 // Set up writes
1337 for (const polyPatch& pp : patches)
1338 {
1339 const auto* ppp = isA<processorPolyPatch>(pp);
1340
1341 if (ppp && pp.size())
1342 {
1343 const auto& procPatch = *ppp;
1344 const label patchi = pp.index();
1345
1346 const labelRange range(pp.start()-boundaryOffset, pp.size());
1347
1348 auto& sendbuf = sendBufs.emplace_set(patchi, faceValues, range);
1349
1351 (
1353 procPatch.neighbProcNo(),
1354 sendbuf.cdata_bytes(),
1355 sendbuf.size_bytes()
1356 );
1357 }
1358 }
1359
1360 // Wait for all comms to finish
1361 UPstream::waitRequests(startRequest);
1362
1363 // Combine with existing data
1364 for (const polyPatch& pp : patches)
1365 {
1366 const auto* ppp = isA<processorPolyPatch>(pp);
1367
1368 if (ppp && pp.size())
1369 {
1370 const label patchi = pp.index();
1371 const label patchSize = pp.size();
1372
1373 const auto& recvbuf = recvBufs[patchi];
1374
1375 // Combine (bitwise)
1376 label bFacei = pp.start()-boundaryOffset;
1377 for (label i = 0; i < patchSize; ++i)
1378 {
1379 unsigned int recvVal = recvbuf[i];
1380 unsigned int faceVal = faceValues[bFacei];
1381
1382 cop(faceVal, recvVal);
1383 faceValues.set(bFacei, faceVal);
1384
1385 ++bFacei;
1386 }
1387 }
1388 }
1389 }
1390
1391
1392 // Do the cyclics.
1393 for (const polyPatch& pp : patches)
1394 {
1396
1397 if (cpp && cpp->owner())
1398 {
1399 // Owner does all.
1400
1401 const cyclicPolyPatch& cycPatch = *cpp;
1402 const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
1403 const label patchSize = cycPatch.size();
1404
1405 label face0 = cycPatch.start()-boundaryOffset;
1406 label face1 = nbrPatch.start()-boundaryOffset;
1407 for (label i = 0; i < patchSize; ++i)
1408 {
1409 unsigned int val0 = faceValues[face0];
1410 unsigned int val1 = faceValues[face1];
1411
1412 unsigned int t = val0;
1413 cop(t, val1);
1414 faceValues[face0] = t;
1415
1416 cop(val1, val0);
1417 faceValues[face1] = val1;
1418
1419 ++face0;
1420 ++face1;
1421 }
1422 }
1423 }
1425 // Reset tag
1426 UPstream::msgType(oldTag);
1427}
1428
1429
1430template<class T>
1432(
1433 const polyMesh& mesh,
1434 const UList<T>& cellData,
1435 List<T>& neighbourCellData,
1436 const bool parRun
1437)
1438{
1439 if (cellData.size() != mesh.nCells())
1440 {
1442 << "Number of cell values " << cellData.size()
1443 << " != number of cells " << mesh.nCells() << nl
1444 << abort(FatalError);
1445 }
1446
1447 const polyBoundaryMesh& patches = mesh.boundaryMesh();
1448
1449 neighbourCellData.resize(mesh.nBoundaryFaces());
1450
1451 for (const polyPatch& pp : patches)
1452 {
1453 const auto& faceCells = pp.faceCells();
1454
1455 // ie, boundarySlice() = patchInternalList()
1456 SubList<T>
1457 (
1458 neighbourCellData,
1459 faceCells.size(),
1460 pp.offset()
1461 ) = UIndirectList<T>(cellData, faceCells);
1463
1464 syncTools::swapBoundaryFaceList(mesh, neighbourCellData, parRun);
1465}
1466
1467
1468template<unsigned Width, class CombineOp>
1470(
1471 const polyMesh& mesh,
1472 PackedList<Width>& faceValues,
1473 const CombineOp& cop,
1474 const bool parRun
1476{
1477 syncFaceList(mesh, false, faceValues, cop, parRun);
1478}
1479
1480
1481template<unsigned Width, class CombineOp>
1483(
1484 const polyMesh& mesh,
1485 PackedList<Width>& faceValues,
1486 const CombineOp& cop,
1487 const bool parRun
1489{
1490 syncFaceList(mesh, true, faceValues, cop, parRun);
1491}
1492
1493
1494template<unsigned Width>
1496(
1497 const polyMesh& mesh,
1498 PackedList<Width>& faceValues,
1499 const bool parRun
1500)
1501{
1502 syncFaceList
1503 (
1504 mesh,
1505 faceValues,
1507 parRun
1508 );
1509}
1510
1511
1512template<unsigned Width>
1514(
1515 const polyMesh& mesh,
1516 PackedList<Width>& faceValues,
1517 const bool parRun
1518)
1519{
1520 syncBoundaryFaceList
1521 (
1522 mesh,
1523 faceValues,
1525 parRun
1526 );
1527}
1528
1529
1530template<unsigned Width, class CombineOp>
1532(
1533 const polyMesh& mesh,
1534 PackedList<Width>& pointValues,
1535 const CombineOp& cop,
1536 const unsigned int nullValue
1537)
1538{
1539 if (pointValues.size() != mesh.nPoints())
1540 {
1542 << "Number of values " << pointValues.size()
1543 << " != number of points " << mesh.nPoints() << nl
1544 << abort(FatalError);
1545 }
1546
1547 const globalMeshData& gd = mesh.globalData();
1548 const labelList& meshPoints = gd.coupledPatch().meshPoints();
1549
1550 List<unsigned int> cppFld(gd.globalPointSlavesMap().constructSize());
1551 forAll(meshPoints, i)
1552 {
1553 cppFld[i] = pointValues[meshPoints[i]];
1554 }
1555
1557 (
1558 cppFld,
1559 gd.globalPointSlaves(),
1560 gd.globalPointTransformedSlaves(),
1561 gd.globalPointSlavesMap(),
1562 cop
1563 );
1564
1565 // Extract back to mesh
1566 forAll(meshPoints, i)
1568 pointValues[meshPoints[i]] = cppFld[i];
1569 }
1570}
1571
1572
1573template<unsigned Width, class CombineOp>
1575(
1576 const polyMesh& mesh,
1577 PackedList<Width>& edgeValues,
1578 const CombineOp& cop,
1579 const unsigned int nullValue
1580)
1581{
1582 if (edgeValues.size() != mesh.nEdges())
1583 {
1585 << "Number of values " << edgeValues.size()
1586 << " != number of edges " << mesh.nEdges() << nl
1587 << abort(FatalError);
1588 }
1589
1590 const globalMeshData& gd = mesh.globalData();
1591 const labelList& meshEdges = gd.coupledPatchMeshEdges();
1592
1593 List<unsigned int> cppFld(gd.globalEdgeSlavesMap().constructSize());
1594 forAll(meshEdges, i)
1595 {
1596 cppFld[i] = edgeValues[meshEdges[i]];
1597 }
1598
1600 (
1601 cppFld,
1602 gd.globalEdgeSlaves(),
1603 gd.globalEdgeTransformedSlaves(),
1604 gd.globalEdgeSlavesMap(),
1605 cop
1606 );
1607
1608 // Extract back to mesh
1609 forAll(meshEdges, i)
1610 {
1611 edgeValues[meshEdges[i]] = cppFld[i];
1612 }
1613}
1614
1615
1616// ************************************************************************* //
scalar range
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
label push_uniq(const T &val)
Append an element if not already in the list.
Map from edge (expressed as its endpoints) to value. Hashing (and ==) on an edge is symmetric.
Definition edgeHashes.H:59
bool set(const Key &key, const T &obj)
Copy assign a new entry, overwriting existing entries.
Definition HashTableI.H:174
const_iterator cfind(const Key &key) const
Find and return an const_iterator set at the hashed entry.
Definition HashTableI.H:113
void reserve(label numEntries)
Reserve space for at least the specified number of elements (not the number of buckets) and regenerat...
Definition HashTable.C:729
bool insert(const Key &key, const T &obj)
Copy insert a new entry, not overwriting existing entries.
Definition HashTableI.H:152
label size() const noexcept
The number of elements in table.
Definition HashTable.H:358
static void recv(Type &value, const int fromProcNo, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Receive and deserialize a value. Uses operator>> for de-serialization.
Definition IPstream.H:80
label size() const noexcept
The number of elements in the list.
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 resize(const label len)
Adjust allocated size of list.
Definition ListI.H:153
A HashTable to objects of type <T> with a label key.
Definition Map.H:54
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Definition PackedList.H:146
bool set(const label i, unsigned int val=~0u)
Set value at index i, default value set is the max_value.
label size() const noexcept
Number of entries.
Definition PackedList.H:392
label nEdges() const
Number of edges in patch.
label nPoints() const
Number of points supporting patch faces.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const Map< label > & meshPointMap() const
Mesh point map.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Buffers for inter-processor communications streams (UOPstream, UIPstream).
void registerSend(const label proci, const bool toggleOn=true)
Toggle an individual send buffer as 'registered'. The setting is sticky (does not turn off).
void initRegisterSend()
Initialise registerSend() bookkeeping by mark all send buffers as 'unregistered'.
label recvDataCount(const label proci) const
Number of unconsumed receive bytes for the specified processor. Must call finishedSends() or other fi...
void finishedNeighbourSends(const labelUList &neighProcs, const bool wait=true)
Mark the send phase as being finished, with communication being limited to a known subset of send/rec...
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
T & emplace_set(const label i, Args &&... args)
Construct and set a new element at given position, (discard old element at that location).
Definition PtrListI.H:191
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
Input inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UIPstream.H:313
static std::streamsize read(const UPstream::commsTypes commsType, const int fromProcNo, Type *buffer, std::streamsize count, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, UPstream::Request *req=nullptr)
Receive buffer contents (contiguous types) from given processor.
A List with indirect addressing. Like IndirectList but does not store addressing.
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
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
label fcIndex(const label i) const noexcept
The forward circular index. The next index in the list which returns to the first at the end of the l...
Definition UListI.H:99
Output inter-processor communications stream using MPI send/recv etc. - operating on external buffer.
Definition UOPstream.H:408
static bool write(const UPstream::commsTypes commsType, const int toProcNo, const Type *buffer, std::streamsize count, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, UPstream::Request *req=nullptr, const UPstream::sendModes sendMode=UPstream::sendModes::normal)
Write buffer contents (contiguous types only) to given processor.
bool send()
Send buffer contents now and not in destructor [advanced usage]. Returns true on success.
Definition OPstreams.C:84
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
Definition UPstream.H:1669
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition UPstream.H:1691
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition UPstream.H:1866
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
static void waitRequests()
Wait for all requests to finish.
Definition UPstream.H:2497
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
static commsTypes defaultCommsType
Default commsType.
Definition UPstream.H:1045
static bool & parRun() noexcept
Test if this a parallel run.
Definition UPstream.H:1681
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
Cyclic plane patch.
virtual bool owner() const
Does this side own the patch ?
const edgeList & coupledEdges() const
Return connected edges (from patch local to neighbour patch local).
const edgeList & coupledPoints() const
Return connected points (from patch local to neighbour patch local).
const cyclicPolyPatch & neighbPatch() const
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition edge.H:62
static int compare(const edge &a, const edge &b)
Compare edges.
Definition edgeI.H:26
linePointRef line(const UList< point > &pts) const
Return edge line.
Definition edgeI.H:463
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
Determination and storage of the possible independent transforms introduced by coupledPolyPatches,...
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points list.
const bitSet & globalEdgeOrientation() const
Is my edge same orientation as master edge.
const labelListList & globalEdgeTransformedSlaves() const
const Map< label > & coupledPatchMeshEdgeMap() const
Return map from mesh edges to coupledPatch edges.
label nGlobalPoints() const
Return number of globally shared points.
const mapDistribute & globalPointSlavesMap() const
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
const mapDistribute & globalEdgeSlavesMap() const
static void syncData(List< Type > &elems, const labelListList &slaves, const labelListList &transformedSlaves, const mapDistribute &slavesMap, const globalIndexAndTransform &, const CombineOp &cop, const TransformOp &top)
Helper: synchronise data with transforms.
const labelList & coupledPatchMeshEdges() const
Return map from coupledPatch edges to mesh edges.
const labelListList & globalPointSlaves() const
const globalIndexAndTransform & globalTransforms() const
Global transforms numbering.
const labelListList & globalPointTransformedSlaves() const
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
const labelListList & globalEdgeSlaves() const
A range or interval of labels defined by a start and a size.
Definition labelRange.H:66
label constructSize() const noexcept
Constructed data size.
Class containing processor-to-processor mapping information.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
label nProcessorPatches() const
The number of processorPolyPatch patches.
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
virtual const faceList & faces() const
Return raw faces.
Definition polyMesh.C:1088
const globalMeshData & globalData() const
Return parallel info (demand-driven).
Definition polyMesh.C:1296
virtual const pointField & points() const
Return raw points.
Definition polyMesh.C:1063
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label start() const noexcept
Return start label of this patch in the polyMesh face list.
Definition polyPatch.H:446
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces).
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
label nInternalFaces() const noexcept
Number of internal faces.
label nFaces() const noexcept
Number of mesh faces.
label nEdges() const
Number of mesh edges.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
static void syncEdgeMap(const polyMesh &mesh, EdgeMap< T > &edgeValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected edges.
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=UPstream::parRun())
Synchronize values on boundary faces only.
static void syncPointMap(const polyMesh &mesh, Map< T > &pointValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected points.
static void swapFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled face values. Uses eqOp.
Definition syncTools.H:567
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const bool parRun=UPstream::parRun())
Swap coupled boundary face values. Uses eqOp.
Definition syncTools.H:524
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition syncTools.H:465
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top, const FlipOp &fop)
Synchronize values on all mesh edges.
const volScalarField & T
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const dimensionedScalar mp
Proton mass.
List< edge > edgeList
List of edge.
Definition edgeList.H:32
List< label > labelList
A List of labels.
Definition List.H:62
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field.
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...
constexpr bool is_contiguous_v
The is_contiguous value of Type (after stripping of qualifiers).
Definition contiguous.H:77
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
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
labelList f(nPoints)
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Spatial transformation functions for list of values and primitive fields.