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-2019 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 "contiguous.H"
35 #include "transformList.H"
36 #include "SubField.H"
37 
38 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
39 
40 template<class T, class CombineOp>
41 void Foam::syncTools::combine
42 (
43  Map<T>& pointValues,
44  const CombineOp& cop,
45  const label index,
46  const T& val
47 )
48 {
49  auto iter = pointValues.find(index);
50 
51  if (iter.found())
52  {
53  cop(*iter, val);
54  }
55  else
56  {
57  pointValues.insert(index, val);
58  }
59 }
60 
61 
62 template<class T, class CombineOp>
63 void Foam::syncTools::combine
64 (
65  EdgeMap<T>& edgeValues,
66  const CombineOp& cop,
67  const edge& index,
68  const T& val
69 )
70 {
71  auto iter = edgeValues.find(index);
72 
73  if (iter.found())
74  {
75  cop(*iter, val);
76  }
77  else
78  {
79  edgeValues.insert(index, val);
80  }
81 }
82 
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(0);
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.resize(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.found())
117  {
118  combine
119  (
120  sharedPointValues,
121  cop,
122  sharedPtAddr[i], // index
123  *fnd // value
124  );
125  }
126  }
127  }
128 
129 
130  if (Pstream::parRun())
131  {
133 
134  // Send
135 
136  for (const polyPatch& pp : patches)
137  {
138  if (isA<processorPolyPatch>(pp) && pp.nPoints() > 0)
139  {
140  const processorPolyPatch& procPatch =
141  refCast<const processorPolyPatch>(pp);
142 
143  // Get data per patchPoint in neighbouring point numbers.
144 
145  const labelList& meshPts = procPatch.meshPoints();
146  const labelList& nbrPts = procPatch.neighbPoints();
147 
148  // Extract local values. Create map from nbrPoint to value.
149  // Note: how small initial size?
150  Map<T> patchInfo(meshPts.size() / 20);
151 
152  forAll(meshPts, i)
153  {
154  const auto iter = pointValues.cfind(meshPts[i]);
155 
156  if (iter.found())
157  {
158  patchInfo.insert(nbrPts[i], *iter);
159  }
160  }
161 
162  UOPstream toNeighb(procPatch.neighbProcNo(), pBufs);
163  toNeighb << patchInfo;
164  }
165  }
166 
167  pBufs.finishedSends();
168 
169  // Receive and combine.
170 
171  for (const polyPatch& pp : patches)
172  {
173  if (isA<processorPolyPatch>(pp) && pp.nPoints() > 0)
174  {
175  const processorPolyPatch& procPatch =
176  refCast<const processorPolyPatch>(pp);
177 
178  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
179  Map<T> nbrPatchInfo(fromNbr);
180 
181  // Transform
182  top(procPatch, nbrPatchInfo);
183 
184  const labelList& meshPts = procPatch.meshPoints();
185 
186  // Only update those values which come from neighbour
187  forAllConstIters(nbrPatchInfo, nbrIter)
188  {
189  combine
190  (
191  pointValues,
192  cop,
193  meshPts[nbrIter.key()],
194  nbrIter.val()
195  );
196  }
197  }
198  }
199  }
200 
201  // Do the cyclics.
202  for (const polyPatch& pp : patches)
203  {
204  if (isA<cyclicPolyPatch>(pp))
205  {
206  const cyclicPolyPatch& cycPatch =
207  refCast<const cyclicPolyPatch>(pp);
208 
209  if (cycPatch.owner())
210  {
211  // Owner does all.
212 
213  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
214  const edgeList& coupledPoints = cycPatch.coupledPoints();
215  const labelList& meshPtsA = cycPatch.meshPoints();
216  const labelList& meshPtsB = nbrPatch.meshPoints();
217 
218  // Extract local values. Create map from coupled-edge to value.
219  Map<T> half0Values(meshPtsA.size() / 20);
220  Map<T> half1Values(half0Values.size());
221 
222  forAll(coupledPoints, i)
223  {
224  const edge& e = coupledPoints[i];
225 
226  const auto point0Fnd = pointValues.cfind(meshPtsA[e[0]]);
227 
228  if (point0Fnd.found())
229  {
230  half0Values.insert(i, *point0Fnd);
231  }
232 
233  const auto point1Fnd = pointValues.cfind(meshPtsB[e[1]]);
234 
235  if (point1Fnd.found())
236  {
237  half1Values.insert(i, *point1Fnd);
238  }
239  }
240 
241  // Transform to receiving side
242  top(cycPatch, half1Values);
243  top(nbrPatch, half0Values);
244 
245  forAll(coupledPoints, i)
246  {
247  const edge& e = coupledPoints[i];
248 
249  const auto half0Fnd = half0Values.cfind(i);
250 
251  if (half0Fnd.found())
252  {
253  combine
254  (
255  pointValues,
256  cop,
257  meshPtsB[e[1]],
258  *half0Fnd
259  );
260  }
261 
262  const auto half1Fnd = half1Values.cfind(i);
263 
264  if (half1Fnd.found())
265  {
266  combine
267  (
268  pointValues,
269  cop,
270  meshPtsA[e[0]],
271  *half1Fnd
272  );
273  }
274  }
275  }
276  }
277  }
278 
279  // Synchronize multiple shared points.
280  if (pd.nGlobalPoints() > 0)
281  {
282  // meshPoint per local index
283  const labelList& sharedPtLabels = pd.sharedPointLabels();
284  // global shared index per local index
285  const labelList& sharedPtAddr = pd.sharedPointAddr();
286 
287  // Reduce on master.
288 
289  if (Pstream::parRun())
290  {
291  if (Pstream::master())
292  {
293  // Receive the edges using shared points from the slave.
294  for
295  (
296  int slave=Pstream::firstSlave();
297  slave<=Pstream::lastSlave();
298  slave++
299  )
300  {
301  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
302  Map<T> nbrValues(fromSlave);
303 
304  // Merge neighbouring values with my values
305  forAllConstIters(nbrValues, iter)
306  {
307  combine
308  (
309  sharedPointValues,
310  cop,
311  iter.key(), // edge
312  iter.val() // value
313  );
314  }
315  }
316 
317  // Send back
318  for
319  (
320  int slave=Pstream::firstSlave();
321  slave<=Pstream::lastSlave();
322  slave++
323  )
324  {
326  toSlave << sharedPointValues;
327  }
328  }
329  else
330  {
331  // Slave: send to master
332  {
333  OPstream toMaster
334  (
337  );
338  toMaster << sharedPointValues;
339  }
340  // Receive merged values
341  {
342  IPstream fromMaster
343  (
346  );
347  fromMaster >> sharedPointValues;
348  }
349  }
350  }
351 
352 
353  // Merge sharedPointValues (keyed on sharedPointAddr) into
354  // pointValues (keyed on mesh points).
355 
356  // Map from global shared index to meshpoint
357  Map<label> sharedToMeshPoint(2*sharedPtAddr.size());
358  forAll(sharedPtAddr, i)
359  {
360  sharedToMeshPoint.insert(sharedPtAddr[i], sharedPtLabels[i]);
361  }
362 
363  forAllConstIters(sharedToMeshPoint, iter)
364  {
365  // Do I have a value for my shared point
366  const auto sharedFnd = sharedPointValues.cfind(iter.key());
367 
368  if (sharedFnd.found())
369  {
370  pointValues.set(*iter, *sharedFnd);
371  }
372  }
373  }
374 }
375 
376 
377 template<class T, class CombineOp, class TransformOp>
379 (
380  const polyMesh& mesh,
381  EdgeMap<T>& edgeValues,
382  const CombineOp& cop,
383  const TransformOp& top
384 )
385 {
387 
388 
389  // Do synchronisation without constructing globalEdge addressing
390  // (since this constructs mesh edge addressing)
391 
392 
393  // Swap proc patch info
394  // ~~~~~~~~~~~~~~~~~~~~
395 
396  if (Pstream::parRun())
397  {
399 
400  // Send
401 
402  for (const polyPatch& pp : patches)
403  {
404  if (isA<processorPolyPatch>(pp) && pp.nEdges() > 0)
405  {
406  const processorPolyPatch& procPatch =
407  refCast<const processorPolyPatch>(pp);
408 
409 
410  // Get data per patch edge in neighbouring edge.
411 
412  const edgeList& edges = procPatch.edges();
413  const labelList& meshPts = procPatch.meshPoints();
414  const labelList& nbrPts = procPatch.neighbPoints();
415 
416  EdgeMap<T> patchInfo(edges.size() / 20);
417 
418  for (const edge& e : edges)
419  {
420  const edge meshEdge(meshPts[e[0]], meshPts[e[1]]);
421 
422  const auto iter = edgeValues.cfind(meshEdge);
423 
424  if (iter.found())
425  {
426  const edge nbrEdge(nbrPts[e[0]], nbrPts[e[1]]);
427  patchInfo.insert(nbrEdge, *iter);
428  }
429  }
430 
431  UOPstream toNeighb(procPatch.neighbProcNo(), pBufs);
432  toNeighb << patchInfo;
433  }
434  }
435 
436  pBufs.finishedSends();
437 
438  // Receive and combine.
439 
440  for (const polyPatch& pp : patches)
441  {
442  if (isA<processorPolyPatch>(pp) && pp.nEdges() > 0)
443  {
444  const processorPolyPatch& procPatch =
445  refCast<const processorPolyPatch>(pp);
446 
447  EdgeMap<T> nbrPatchInfo;
448  {
449  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
450  fromNbr >> nbrPatchInfo;
451  }
452 
453  // Apply transform to convert to this side properties.
454  top(procPatch, nbrPatchInfo);
455 
456 
457  // Only update those values which come from neighbour
458  const labelList& meshPts = procPatch.meshPoints();
459 
460  forAllConstIters(nbrPatchInfo, nbrIter)
461  {
462  const edge& e = nbrIter.key();
463  const edge meshEdge(meshPts[e[0]], meshPts[e[1]]);
464 
465  combine
466  (
467  edgeValues,
468  cop,
469  meshEdge, // edge
470  nbrIter.val() // value
471  );
472  }
473  }
474  }
475  }
476 
477 
478  // Swap cyclic info
479  // ~~~~~~~~~~~~~~~~
480 
481  for (const polyPatch& pp : patches)
482  {
483  if (isA<cyclicPolyPatch>(pp))
484  {
485  const cyclicPolyPatch& cycPatch =
486  refCast<const cyclicPolyPatch>(pp);
487 
488  if (cycPatch.owner())
489  {
490  // Owner does all.
491 
492  const edgeList& coupledEdges = cycPatch.coupledEdges();
493  const labelList& meshPtsA = cycPatch.meshPoints();
494  const edgeList& edgesA = cycPatch.edges();
495  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
496  const labelList& meshPtsB = nbrPatch.meshPoints();
497  const edgeList& edgesB = nbrPatch.edges();
498 
499  // Extract local values. Create map from edge to value.
500  Map<T> half0Values(edgesA.size() / 20);
501  Map<T> half1Values(half0Values.size());
502 
503  forAll(coupledEdges, i)
504  {
505  const edge& twoEdges = coupledEdges[i];
506 
507  {
508  const edge& e0 = edgesA[twoEdges[0]];
509  const edge meshEdge0(meshPtsA[e0[0]], meshPtsA[e0[1]]);
510 
511  const auto iter = edgeValues.cfind(meshEdge0);
512 
513  if (iter.found())
514  {
515  half0Values.insert(i, *iter);
516  }
517  }
518  {
519  const edge& e1 = edgesB[twoEdges[1]];
520  const edge meshEdge1(meshPtsB[e1[0]], meshPtsB[e1[1]]);
521 
522  const auto iter = edgeValues.cfind(meshEdge1);
523 
524  if (iter.found())
525  {
526  half1Values.insert(i, *iter);
527  }
528  }
529  }
530 
531  // Transform to this side
532  top(cycPatch, half1Values);
533  top(nbrPatch, half0Values);
534 
535 
536  // Extract and combine information
537 
538  forAll(coupledEdges, i)
539  {
540  const edge& twoEdges = coupledEdges[i];
541 
542  const auto half1Fnd = half1Values.cfind(i);
543 
544  if (half1Fnd.found())
545  {
546  const edge& e0 = edgesA[twoEdges[0]];
547  const edge meshEdge0(meshPtsA[e0[0]], meshPtsA[e0[1]]);
548 
549  combine
550  (
551  edgeValues,
552  cop,
553  meshEdge0, // edge
554  *half1Fnd // value
555  );
556  }
557 
558  const auto half0Fnd = half0Values.cfind(i);
559 
560  if (half0Fnd.found())
561  {
562  const edge& e1 = edgesB[twoEdges[1]];
563  const edge meshEdge1(meshPtsB[e1[0]], meshPtsB[e1[1]]);
564 
565  combine
566  (
567  edgeValues,
568  cop,
569  meshEdge1, // edge
570  *half0Fnd // value
571  );
572  }
573  }
574  }
575  }
576  }
577 
578  // Synchronize multiple shared points.
579  // Problem is that we don't want to construct shared edges so basically
580  // we do here like globalMeshData but then using sparse edge representation
581  // (EdgeMap instead of mesh.edges())
582 
583  const globalMeshData& pd = mesh.globalData();
584  const labelList& sharedPtAddr = pd.sharedPointAddr();
585  const labelList& sharedPtLabels = pd.sharedPointLabels();
586 
587  // 1. Create map from meshPoint to globalShared index.
588  Map<label> meshToShared(2*sharedPtLabels.size());
589  forAll(sharedPtLabels, i)
590  {
591  meshToShared.insert(sharedPtLabels[i], sharedPtAddr[i]);
592  }
593 
594  // Values on shared points. From two sharedPtAddr indices to a value.
595  EdgeMap<T> sharedEdgeValues(meshToShared.size());
596 
597  // From shared edge to mesh edge. Used for merging later on.
598  EdgeMap<edge> potentialSharedEdge(meshToShared.size());
599 
600  // 2. Find any edges using two global shared points. These will always be
601  // on the outside of the mesh. (though might not be on coupled patch
602  // if is single edge and not on coupled face)
603  // Store value (if any) on sharedEdgeValues
604  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
605  {
606  const face& f = mesh.faces()[facei];
607 
608  forAll(f, fp)
609  {
610  const label v0 = f[fp];
611  const label v1 = f[f.fcIndex(fp)];
612 
613  const auto v0Fnd = meshToShared.cfind(v0);
614 
615  if (v0Fnd.found())
616  {
617  const auto v1Fnd = meshToShared.cfind(v1);
618 
619  if (v1Fnd.found())
620  {
621  const edge meshEdge(v0, v1);
622 
623  // edge in shared point labels
624  const edge sharedEdge(*v0Fnd, *v1Fnd);
625 
626  // Store mesh edge as a potential shared edge.
627  potentialSharedEdge.insert(sharedEdge, meshEdge);
628 
629  const auto edgeFnd = edgeValues.cfind(meshEdge);
630 
631  if (edgeFnd.found())
632  {
633  // edge exists in edgeValues. See if already in map
634  // (since on same processor, e.g. cyclic)
635  combine
636  (
637  sharedEdgeValues,
638  cop,
639  sharedEdge, // edge
640  *edgeFnd // value
641  );
642  }
643  }
644  }
645  }
646  }
647 
648 
649  // Now sharedEdgeValues will contain per potential sharedEdge the value.
650  // (potential since an edge having two shared points is not necessary a
651  // shared edge).
652  // Reduce this on the master.
653 
654  if (Pstream::parRun())
655  {
656  if (Pstream::master())
657  {
658  // Receive the edges using shared points from the slave.
659  for
660  (
661  int slave=Pstream::firstSlave();
662  slave<=Pstream::lastSlave();
663  slave++
664  )
665  {
666  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
667  EdgeMap<T> nbrValues(fromSlave);
668 
669  // Merge neighbouring values with my values
670  forAllConstIters(nbrValues, iter)
671  {
672  combine
673  (
674  sharedEdgeValues,
675  cop,
676  iter.key(), // edge
677  iter.val() // value
678  );
679  }
680  }
681 
682  // Send back
683  for
684  (
685  int slave=Pstream::firstSlave();
686  slave<=Pstream::lastSlave();
687  slave++
688  )
689  {
690 
692  toSlave << sharedEdgeValues;
693  }
694  }
695  else
696  {
697  // Send to master
698  {
699  OPstream toMaster
700  (
703  );
704  toMaster << sharedEdgeValues;
705  }
706  // Receive merged values
707  {
708  IPstream fromMaster
709  (
712  );
713  fromMaster >> sharedEdgeValues;
714  }
715  }
716  }
717 
718 
719  // Merge sharedEdgeValues (keyed on sharedPointAddr) into edgeValues
720  // (keyed on mesh points).
721 
722  // Loop over all my shared edges.
723  forAllConstIters(potentialSharedEdge, iter)
724  {
725  const edge& sharedEdge = iter.key();
726  const edge& meshEdge = iter.val();
727 
728  // Do I have a value for the shared edge?
729  const auto sharedFnd = sharedEdgeValues.cfind(sharedEdge);
730 
731  if (sharedFnd.found())
732  {
733  combine
734  (
735  edgeValues,
736  cop,
737  meshEdge, // edge
738  *sharedFnd // value
739  );
740  }
741  }
742 }
743 
744 
745 template<class T, class CombineOp, class TransformOp>
747 (
748  const polyMesh& mesh,
749  List<T>& pointValues,
750  const CombineOp& cop,
751  const T& nullValue,
752  const TransformOp& top
753 )
754 {
755  if (pointValues.size() != mesh.nPoints())
756  {
758  << "Number of values " << pointValues.size()
759  << " is not equal to the number of points in the mesh "
760  << mesh.nPoints() << abort(FatalError);
761  }
762 
763  mesh.globalData().syncPointData(pointValues, cop, top);
764 }
765 
766 
767 template<class T, class CombineOp, class TransformOp>
769 (
770  const polyMesh& mesh,
771  const labelUList& meshPoints,
772  List<T>& pointValues,
773  const CombineOp& cop,
774  const T& nullValue,
775  const TransformOp& top
776 )
777 {
778  if (pointValues.size() != meshPoints.size())
779  {
781  << "Number of values " << pointValues.size()
782  << " is not equal to the number of meshPoints "
783  << meshPoints.size() << abort(FatalError);
784  }
785  const globalMeshData& gd = mesh.globalData();
786  const indirectPrimitivePatch& cpp = gd.coupledPatch();
787  const Map<label>& mpm = cpp.meshPointMap();
788 
789  List<T> cppFld(cpp.nPoints(), nullValue);
790 
791  forAll(meshPoints, i)
792  {
793  const auto iter = mpm.cfind(meshPoints[i]);
794 
795  if (iter.found())
796  {
797  cppFld[*iter] = pointValues[i];
798  }
799  }
800 
802  (
803  cppFld,
804  gd.globalPointSlaves(),
807  gd.globalTransforms(),
808  cop,
809  top
810  );
811 
812  forAll(meshPoints, i)
813  {
814  const auto iter = mpm.cfind(meshPoints[i]);
815 
816  if (iter.found())
817  {
818  pointValues[i] = cppFld[*iter];
819  }
820  }
821 }
822 
823 
824 template<class T, class CombineOp, class TransformOp>
826 (
827  const polyMesh& mesh,
828  List<T>& edgeValues,
829  const CombineOp& cop,
830  const T& nullValue,
831  const TransformOp& top
832 )
833 {
834  if (edgeValues.size() != mesh.nEdges())
835  {
837  << "Number of values " << edgeValues.size()
838  << " is not equal to the number of edges in the mesh "
839  << mesh.nEdges() << abort(FatalError);
840  }
841 
842  const globalMeshData& gd = mesh.globalData();
843  const labelList& meshEdges = gd.coupledPatchMeshEdges();
844  const globalIndexAndTransform& git = gd.globalTransforms();
845  const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
846 
847  List<T> cppFld(UIndirectList<T>(edgeValues, meshEdges));
848 
850  (
851  cppFld,
852  gd.globalEdgeSlaves(),
854  edgeMap,
855  git,
856  cop,
857  top
858  );
859 
860  // Extract back onto mesh
861  forAll(meshEdges, i)
862  {
863  edgeValues[meshEdges[i]] = cppFld[i];
864  }
865 }
866 
867 
868 template<class T, class CombineOp, class TransformOp>
870 (
871  const polyMesh& mesh,
872  const labelList& meshEdges,
873  List<T>& edgeValues,
874  const CombineOp& cop,
875  const T& nullValue,
876  const TransformOp& top
877 )
878 {
879  if (edgeValues.size() != meshEdges.size())
880  {
882  << "Number of values " << edgeValues.size()
883  << " is not equal to the number of meshEdges "
884  << meshEdges.size() << abort(FatalError);
885  }
886  const globalMeshData& gd = mesh.globalData();
887  const indirectPrimitivePatch& cpp = gd.coupledPatch();
888  const Map<label>& mpm = gd.coupledPatchMeshEdgeMap();
889 
890  List<T> cppFld(cpp.nEdges(), nullValue);
891 
892  forAll(meshEdges, i)
893  {
894  const auto iter = mpm.cfind(meshEdges[i]);
895 
896  if (iter.found())
897  {
898  cppFld[*iter] = edgeValues[i];
899  }
900  }
901 
903  (
904  cppFld,
905  gd.globalEdgeSlaves(),
907  gd.globalEdgeSlavesMap(),
908  gd.globalTransforms(),
909  cop,
910  top
911  );
912 
913  forAll(meshEdges, i)
914  {
915  const auto iter = mpm.cfind(meshEdges[i]);
916 
917  if (iter.found())
918  {
919  edgeValues[i] = cppFld[*iter];
920  }
921  }
922 }
923 
924 
925 template<class T, class CombineOp, class TransformOp>
927 (
928  const polyMesh& mesh,
929  UList<T>& faceValues,
930  const CombineOp& cop,
931  const TransformOp& top,
932  const bool parRun
933 )
934 {
935  // Offset (global to local) for start of boundaries
936  const label boundaryOffset = mesh.nInternalFaces();
937 
938  if (faceValues.size() != mesh.nBoundaryFaces())
939  {
941  << "Number of values " << faceValues.size()
942  << " is not equal to the number of boundary faces in the mesh "
943  << mesh.nBoundaryFaces() << nl
944  << abort(FatalError);
945  }
946 
948 
949  if (parRun)
950  {
952 
953  // Send
954 
955  for (const polyPatch& pp : patches)
956  {
957  if (isA<processorPolyPatch>(pp) && pp.size() > 0)
958  {
959  const processorPolyPatch& procPatch =
960  refCast<const processorPolyPatch>(pp);
961 
962  const label patchStart = procPatch.start()-boundaryOffset;
963 
964  // Send slice of values on the patch
965  UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
966  toNbr<< SubList<T>(faceValues, procPatch.size(), patchStart);
967  }
968  }
969 
970 
971  pBufs.finishedSends();
972 
973 
974  // Receive and combine.
975 
976  for (const polyPatch& pp : patches)
977  {
978  if (isA<processorPolyPatch>(pp) && pp.size() > 0)
979  {
980  const processorPolyPatch& procPatch =
981  refCast<const processorPolyPatch>(pp);
982 
983  Field<T> nbrVals(procPatch.size());
984 
985  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
986  fromNbr >> nbrVals;
987 
988  top(procPatch, nbrVals);
989 
990  label bFacei = procPatch.start()-boundaryOffset;
991 
992  for (T& nbrVal : nbrVals)
993  {
994  cop(faceValues[bFacei++], nbrVal);
995  }
996  }
997  }
998  }
999 
1000  // Do the cyclics.
1001  for (const polyPatch& pp : patches)
1002  {
1003  if (isA<cyclicPolyPatch>(pp))
1004  {
1005  const cyclicPolyPatch& cycPatch =
1006  refCast<const cyclicPolyPatch>(pp);
1007 
1008  if (cycPatch.owner())
1009  {
1010  // Owner does all.
1011  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
1012 
1013  const label patchSize = cycPatch.size();
1014  const label ownStart = cycPatch.start()-boundaryOffset;
1015  const label nbrStart = nbrPatch.start()-boundaryOffset;
1016 
1017  // Transform (copy of) data on both sides
1018  Field<T> ownVals(SubList<T>(faceValues, patchSize, ownStart));
1019  top(nbrPatch, ownVals);
1020 
1021  Field<T> nbrVals(SubList<T>(faceValues, patchSize, nbrStart));
1022  top(cycPatch, nbrVals);
1023 
1024  label bFacei = ownStart;
1025  for (T& nbrVal : nbrVals)
1026  {
1027  cop(faceValues[bFacei++], nbrVal);
1028  }
1029 
1030  bFacei = nbrStart;
1031  for (T& ownVal : ownVals)
1032  {
1033  cop(faceValues[bFacei++], ownVal);
1034  }
1035  }
1036  }
1037  }
1038 }
1039 
1040 
1041 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1042 
1043 template<unsigned Width, class CombineOp>
1046  const polyMesh& mesh,
1047  const bool isBoundaryOnly,
1048  PackedList<Width>& faceValues,
1049  const CombineOp& cop,
1050  const bool parRun
1051 )
1052 {
1053  // Offset (global to local) for start of boundaries
1054  const label boundaryOffset = (isBoundaryOnly ? mesh.nInternalFaces() : 0);
1055 
1056  if
1057  (
1058  faceValues.size()
1059  != (isBoundaryOnly ? mesh.nBoundaryFaces() : mesh.nFaces())
1060  )
1061  {
1063  << "Number of values " << faceValues.size()
1064  << " is not equal to the number of "
1065  << (isBoundaryOnly ? "boundary" : "mesh") << " faces "
1066  << (isBoundaryOnly ? mesh.nBoundaryFaces() : mesh.nFaces()) << nl
1067  << abort(FatalError);
1068  }
1069 
1071 
1072  if (parRun)
1073  {
1075 
1076  // Send
1077 
1078  for (const polyPatch& pp : patches)
1079  {
1080  if (isA<processorPolyPatch>(pp) && pp.size())
1081  {
1082  const processorPolyPatch& procPatch =
1083  refCast<const processorPolyPatch>(pp);
1084 
1085  const labelRange range
1086  (
1087  procPatch.start()-boundaryOffset,
1088  procPatch.size()
1089  );
1090 
1091  // Send slice of values on the patch
1092  UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
1093  toNbr<< PackedList<Width>(faceValues, range);
1094  }
1095  }
1096 
1097  pBufs.finishedSends();
1098 
1099 
1100  // Receive and combine.
1101 
1102  for (const polyPatch& pp : patches)
1103  {
1104  if (isA<processorPolyPatch>(pp) && pp.size())
1105  {
1106  const processorPolyPatch& procPatch =
1107  refCast<const processorPolyPatch>(pp);
1108 
1109  const label patchSize = procPatch.size();
1110 
1111  // Recv slice of values on the patch
1112  PackedList<Width> recvInfo(patchSize);
1113  {
1114  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
1115  fromNbr >> recvInfo;
1116  }
1117 
1118  // Combine (bitwise)
1119  label bFacei = procPatch.start()-boundaryOffset;
1120  for (label i = 0; i < patchSize; ++i)
1121  {
1122  unsigned int recvVal = recvInfo[i];
1123  unsigned int faceVal = faceValues[bFacei];
1124 
1125  cop(faceVal, recvVal);
1126  faceValues.set(bFacei, faceVal);
1127 
1128  ++bFacei;
1129  }
1130  }
1131  }
1132  }
1133 
1134 
1135  // Do the cyclics.
1136  for (const polyPatch& pp : patches)
1137  {
1138  if (isA<cyclicPolyPatch>(pp))
1139  {
1140  const cyclicPolyPatch& cycPatch =
1141  refCast<const cyclicPolyPatch>(pp);
1142 
1143  if (cycPatch.owner())
1144  {
1145  // Owner does all.
1146  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
1147 
1148  const label patchSize = cycPatch.size();
1149 
1150  label face0 = cycPatch.start()-boundaryOffset;
1151  label face1 = nbrPatch.start()-boundaryOffset;
1152  for (label i = 0; i < patchSize; ++i)
1153  {
1154  unsigned int val0 = faceValues[face0];
1155  unsigned int val1 = faceValues[face1];
1156 
1157  unsigned int t = val0;
1158  cop(t, val1);
1159  faceValues[face0] = t;
1160 
1161  cop(val1, val0);
1162  faceValues[face1] = val1;
1163 
1164  ++face0;
1165  ++face1;
1166  }
1167  }
1168  }
1169  }
1170 }
1171 
1172 
1173 template<class T>
1176  const polyMesh& mesh,
1177  const UList<T>& cellData,
1178  List<T>& neighbourCellData
1179 )
1180 {
1181  if (cellData.size() != mesh.nCells())
1182  {
1184  << "Number of cell values " << cellData.size()
1185  << " is not equal to the number of cells in the mesh "
1186  << mesh.nCells() << abort(FatalError);
1187  }
1188 
1190 
1191  neighbourCellData.resize(mesh.nBoundaryFaces());
1192 
1193  for (const polyPatch& pp : patches)
1194  {
1195  label bFacei = pp.start()-mesh.nInternalFaces();
1196 
1197  const labelUList& faceCells = pp.faceCells();
1198 
1199  for (const label celli : faceCells)
1200  {
1201  neighbourCellData[bFacei] = cellData[celli];
1202  ++bFacei;
1203  }
1204  }
1205  syncTools::swapBoundaryFaceList(mesh, neighbourCellData);
1206 }
1207 
1208 
1209 template<unsigned Width, class CombineOp>
1212  const polyMesh& mesh,
1213  PackedList<Width>& faceValues,
1214  const CombineOp& cop,
1215  const bool parRun
1216 )
1217 {
1218  syncFaceList(mesh, false, faceValues, cop, parRun);
1219 }
1220 
1221 
1222 template<unsigned Width, class CombineOp>
1225  const polyMesh& mesh,
1226  PackedList<Width>& faceValues,
1227  const CombineOp& cop,
1228  const bool parRun
1229 )
1230 {
1231  syncFaceList(mesh, true, faceValues, cop, parRun);
1232 }
1233 
1234 
1235 template<unsigned Width>
1238  const polyMesh& mesh,
1239  PackedList<Width>& faceValues
1240 )
1241 {
1242  syncFaceList(mesh, faceValues, eqOp<unsigned int>());
1243 }
1244 
1245 
1246 template<unsigned Width>
1249  const polyMesh& mesh,
1250  PackedList<Width>& faceValues
1251 )
1252 {
1253  syncBoundaryFaceList(mesh, faceValues, eqOp<unsigned int>());
1254 }
1255 
1256 
1257 template<unsigned Width, class CombineOp>
1260  const polyMesh& mesh,
1261  PackedList<Width>& pointValues,
1262  const CombineOp& cop,
1263  const unsigned int nullValue
1264 )
1265 {
1266  if (pointValues.size() != mesh.nPoints())
1267  {
1269  << "Number of values " << pointValues.size()
1270  << " is not equal to the number of points in the mesh "
1271  << mesh.nPoints() << abort(FatalError);
1272  }
1273 
1274  const globalMeshData& gd = mesh.globalData();
1275  const labelList& meshPoints = gd.coupledPatch().meshPoints();
1276 
1278  forAll(meshPoints, i)
1279  {
1280  cppFld[i] = pointValues[meshPoints[i]];
1281  }
1282 
1284  (
1285  cppFld,
1286  gd.globalPointSlaves(),
1288  gd.globalPointSlavesMap(),
1289  cop
1290  );
1291 
1292  // Extract back to mesh
1293  forAll(meshPoints, i)
1294  {
1295  pointValues[meshPoints[i]] = cppFld[i];
1296  }
1297 }
1298 
1299 
1300 template<unsigned Width, class CombineOp>
1303  const polyMesh& mesh,
1304  PackedList<Width>& edgeValues,
1305  const CombineOp& cop,
1306  const unsigned int nullValue
1307 )
1308 {
1309  if (edgeValues.size() != mesh.nEdges())
1310  {
1312  << "Number of values " << edgeValues.size()
1313  << " is not equal to the number of edges in the mesh "
1314  << mesh.nEdges() << abort(FatalError);
1315  }
1316 
1317  const globalMeshData& gd = mesh.globalData();
1318  const labelList& meshEdges = gd.coupledPatchMeshEdges();
1319 
1321  forAll(meshEdges, i)
1322  {
1323  cppFld[i] = edgeValues[meshEdges[i]];
1324  }
1325 
1327  (
1328  cppFld,
1329  gd.globalEdgeSlaves(),
1331  gd.globalEdgeSlavesMap(),
1332  cop
1333  );
1334 
1335  // Extract back to mesh
1336  forAll(meshEdges, i)
1337  {
1338  edgeValues[meshEdges[i]] = cppFld[i];
1339  }
1340 }
1341 
1342 
1343 // ************************************************************************* //
Foam::globalMeshData::globalPointTransformedSlaves
const labelListList & globalPointTransformedSlaves() const
Definition: globalMeshData.C:2200
Foam::val
label ListType::const_reference val
Definition: ListOps.H:407
Foam::UPstream::commsTypes::nonBlocking
Foam::globalMeshData::globalPointSlaves
const labelListList & globalPointSlaves() const
Definition: globalMeshData.C:2190
Foam::UPstream::masterNo
static constexpr int masterNo() noexcept
Process index of the master.
Definition: UPstream.H:432
Foam::UOPstream
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:57
SubField.H
Foam::globalMeshData::coupledPatchMeshEdgeMap
const Map< label > & coupledPatchMeshEdgeMap() const
Return map from mesh edges to coupledPatch edges.
Definition: globalMeshData.C:2147
Foam::syncTools::syncPointMap
static void syncPointMap(const polyMesh &mesh, Map< T > &pointValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected points.
Definition: syncToolsTemplates.C:86
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatch.C:238
Foam::globalMeshData::globalEdgeSlavesMap
const mapDistribute & globalEdgeSlavesMap() const
Definition: globalMeshData.C:2265
cyclicPolyPatch.H
Foam::cyclicPolyPatch
Cyclic plane patch.
Definition: cyclicPolyPatch.H:65
Foam::OPstream
Output inter-processor communications stream.
Definition: OPstream.H:52
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::cyclicPolyPatch::neighbPatch
const cyclicPolyPatch & neighbPatch() const
Definition: cyclicPolyPatch.H:330
globalMeshData.H
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::PrimitivePatch::nEdges
label nEdges() const
Return number of edges in patch.
Definition: PrimitivePatch.H:317
Foam::syncTools::syncBoundaryFaceList
static void syncBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronize values on boundary faces only.
Definition: syncToolsTemplates.C:927
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:88
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::primitiveMesh::nEdges
label nEdges() const
Number of mesh edges.
Definition: primitiveMeshI.H:67
Foam::Map
A HashTable to objects of type <T> with a label key.
Definition: HashTableFwd.H:46
Foam::processorPolyPatch::neighbProcNo
int neighbProcNo() const
Return neighbour processor number.
Definition: processorPolyPatch.H:274
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:439
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:435
Foam::globalMeshData::globalPointSlavesMap
const mapDistribute & globalPointSlavesMap() const
Definition: globalMeshData.C:2211
polyMesh.H
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::syncTools::syncPointList
static void syncPointList(const polyMesh &mesh, List< T > &pointValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh points.
Definition: syncToolsTemplates.C:747
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::cyclicPolyPatch::coupledPoints
const edgeList & coupledPoints() const
Return connected points (from patch local to neighbour patch local)
Definition: cyclicPolyPatch.C:1008
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:390
Foam::globalMeshData::nGlobalPoints
label nGlobalPoints() const
Return number of globally shared points.
Definition: globalMeshData.C:2006
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:62
Foam::Field< T >
Foam::UPstream::lastSlave
static int lastSlave(const label communicator=0)
Process index of last slave.
Definition: UPstream.H:467
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::ListListOps::combine
AccessType combine(const UList< T > &lists, AccessOp aop=accessOp< T >())
Combines sub-lists into a single list.
Definition: ListListOps.C:69
Foam::globalMeshData::globalTransforms
const globalIndexAndTransform & globalTransforms() const
Global transforms numbering.
Definition: globalMeshData.C:2180
Foam::globalMeshData::coupledPatchMeshEdges
const labelList & coupledPatchMeshEdges() const
Return map from coupledPatch edges to mesh edges.
Definition: globalMeshData.C:2127
Foam::syncTools::syncEdgeList
static void syncEdgeList(const polyMesh &mesh, List< T > &edgeValues, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronize values on all mesh edges.
Definition: syncToolsTemplates.C:826
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:163
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::labelRange
A range or interval of labels defined by a start and a size.
Definition: labelRange.H:58
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:58
Foam::globalMeshData::globalEdgeSlaves
const labelListList & globalEdgeSlaves() const
Definition: globalMeshData.C:2234
Foam::globalMeshData::syncData
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.
Definition: globalMeshDataTemplates.C:37
Foam::List::resize
void resize(const label newSize)
Adjust allocated size of list.
Definition: ListI.H:139
Foam::syncTools::swapBoundaryCellList
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
Definition: syncToolsTemplates.C:1175
Foam::PstreamBuffers::finishedSends
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
Definition: PstreamBuffers.C:80
Foam::UPstream::commsTypes::scheduled
Foam::cyclicPolyPatch::coupledEdges
const edgeList & coupledEdges() const
Return connected edges (from patch local to neighbour patch local).
Definition: cyclicPolyPatch.C:1089
Foam::PrimitivePatch::nPoints
label nPoints() const
Return number of points supporting patch faces.
Definition: PrimitivePatch.H:311
Foam::cyclicPolyPatch::owner
virtual bool owner() const
Does this side own the patch ?
Definition: cyclicPolyPatch.H:320
Foam::FatalError
error FatalError
processorPolyPatch.H
Foam::primitiveMesh::nBoundaryFaces
label nBoundaryFaces() const
Number of boundary faces (== nFaces - nInternalFaces)
Definition: primitiveMeshI.H:84
Foam::eqOp
Definition: ops.H:71
Foam::globalMeshData
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
Definition: globalMeshData.H:106
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:311
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::UPstream::firstSlave
static constexpr int firstSlave() noexcept
Process index of first slave.
Definition: UPstream.H:461
Foam::mapDistributeBase::constructSize
label constructSize() const
Constructed data size.
Definition: mapDistributeBase.H:270
Foam::EdgeMap
Map from edge (expressed as its endpoints) to value. For easier forward declaration it is currently i...
Definition: EdgeMap.H:51
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1063
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::processorPolyPatch::neighbPoints
const labelList & neighbPoints() const
Return neighbour point labels. WIP.
Definition: processorPolyPatch.C:529
Foam::HashTable< T, edge, Hash< edge > >::cfind
const_iterator cfind(const edge &key) const
Find and return an const_iterator set at the hashed entry.
Definition: HashTableI.H:141
range
scalar range
Definition: LISASMDCalcMethod1.H:12
Foam::nl
constexpr char nl
Definition: Ostream.H:372
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::globalMeshData::globalEdgeTransformedSlaves
const labelListList & globalEdgeTransformedSlaves() const
Definition: globalMeshData.C:2244
f
labelList f(nPoints)
Foam::syncTools::swapFaceList
static void swapFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled face values. Uses eqOp.
Definition: syncTools.H:472
contiguous.H
Foam::List< label >
Foam::globalMeshData::coupledPatch
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
Definition: globalMeshData.C:2066
Foam::PackedList
A dynamic list of packed unsigned integers, with the number of bits per item specified by the <Width>...
Definition: PackedList.H:108
Foam::syncTools::syncEdgeMap
static void syncEdgeMap(const polyMesh &mesh, EdgeMap< T > &edgeValues, const CombineOp &cop, const TransformOp &top)
Synchronize values on selected edges.
Definition: syncToolsTemplates.C:379
Foam::PackedList::size
label size() const noexcept
Number of entries.
Definition: PackedListI.H:366
Foam::UList< label >
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::primitiveMesh::nPoints
label nPoints() const
Number of mesh points.
Definition: primitiveMeshI.H:37
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::UList::size
void size(const label n) noexcept
Override size to be inconsistent with allocated storage.
Definition: UListI.H:360
Foam::globalMeshData::sharedPointLabels
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
Definition: globalMeshData.C:2016
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
Foam::PackedList::set
bool set(const label i, unsigned int val=~0u)
Set value at index i, default value set is the max_value.
Definition: PackedListI.H:608
Foam::globalMeshData::sharedPointAddr
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points.
Definition: globalMeshData.C:2026
Foam::IPstream
Input inter-processor communications stream.
Definition: IPstream.H:52
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
Foam::UIPstream
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:56
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1241
Foam::PrimitivePatch::meshPointMap
const Map< label > & meshPointMap() const
Mesh point map.
Definition: PrimitivePatch.C:438
transformList.H
Spatial transformation functions for list of values and primitive fields.
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:418
Foam::globalIndexAndTransform
Determination and storage of the possible independent transforms introduced by coupledPolyPatches,...
Definition: globalIndexAndTransform.H:64
Foam::faceCells
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:56
Foam::globalMeshData::syncPointData
void syncPointData(List< Type > &pointData, const CombineOp &cop, const TransformOp &top) const
Helper to synchronise coupled patch point data.
Definition: globalMeshDataTemplates.C:159
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:90