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