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