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-2020 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.found())
50  {
51  cop(*iter, 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.found())
72  {
73  cop(*iter, 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 {
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.resize(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.found())
115  {
116  combine
117  (
118  sharedPointValues,
119  cop,
120  sharedPtAddr[i], // index
121  *fnd // value
122  );
123  }
124  }
125  }
126 
127 
128  if (Pstream::parRun())
129  {
131 
132  // Send
133 
134  for (const polyPatch& pp : patches)
135  {
136  if (isA<processorPolyPatch>(pp) && pp.nPoints() > 0)
137  {
138  const processorPolyPatch& procPatch =
139  refCast<const processorPolyPatch>(pp);
140 
141  // Get data per patchPoint in neighbouring point numbers.
142 
143  const labelList& meshPts = procPatch.meshPoints();
144  const labelList& nbrPts = procPatch.neighbPoints();
145 
146  // Extract local values. Create map from nbrPoint to value.
147  // Note: how small initial size?
148  Map<T> patchInfo(meshPts.size() / 20);
149 
150  forAll(meshPts, i)
151  {
152  const auto iter = pointValues.cfind(meshPts[i]);
153 
154  if (iter.found())
155  {
156  patchInfo.insert(nbrPts[i], *iter);
157  }
158  }
159 
160  UOPstream toNeighb(procPatch.neighbProcNo(), pBufs);
161  toNeighb << patchInfo;
162  }
163  }
164 
165  pBufs.finishedSends();
166 
167  // Receive and combine.
168 
169  for (const polyPatch& pp : patches)
170  {
171  if (isA<processorPolyPatch>(pp) && pp.nPoints() > 0)
172  {
173  const processorPolyPatch& procPatch =
174  refCast<const processorPolyPatch>(pp);
175 
176  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
177  Map<T> nbrPatchInfo(fromNbr);
178 
179  // Transform
180  top(procPatch, nbrPatchInfo);
181 
182  const labelList& meshPts = procPatch.meshPoints();
183 
184  // Only update those values which come from neighbour
185  forAllConstIters(nbrPatchInfo, nbrIter)
186  {
187  combine
188  (
189  pointValues,
190  cop,
191  meshPts[nbrIter.key()],
192  nbrIter.val()
193  );
194  }
195  }
196  }
197  }
198 
199  // Do the cyclics.
200  for (const polyPatch& pp : patches)
201  {
202  if (isA<cyclicPolyPatch>(pp))
203  {
204  const cyclicPolyPatch& cycPatch =
205  refCast<const cyclicPolyPatch>(pp);
206 
207  if (cycPatch.owner())
208  {
209  // Owner does all.
210 
211  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
212  const edgeList& coupledPoints = cycPatch.coupledPoints();
213  const labelList& meshPtsA = cycPatch.meshPoints();
214  const labelList& meshPtsB = nbrPatch.meshPoints();
215 
216  // Extract local values. Create map from coupled-edge to value.
217  Map<T> half0Values(meshPtsA.size() / 20);
218  Map<T> half1Values(half0Values.size());
219 
220  forAll(coupledPoints, i)
221  {
222  const edge& e = coupledPoints[i];
223 
224  const auto point0Fnd = pointValues.cfind(meshPtsA[e[0]]);
225 
226  if (point0Fnd.found())
227  {
228  half0Values.insert(i, *point0Fnd);
229  }
230 
231  const auto point1Fnd = pointValues.cfind(meshPtsB[e[1]]);
232 
233  if (point1Fnd.found())
234  {
235  half1Values.insert(i, *point1Fnd);
236  }
237  }
238 
239  // Transform to receiving side
240  top(cycPatch, half1Values);
241  top(nbrPatch, half0Values);
242 
243  forAll(coupledPoints, i)
244  {
245  const edge& e = coupledPoints[i];
246 
247  const auto half0Fnd = half0Values.cfind(i);
248 
249  if (half0Fnd.found())
250  {
251  combine
252  (
253  pointValues,
254  cop,
255  meshPtsB[e[1]],
256  *half0Fnd
257  );
258  }
259 
260  const auto half1Fnd = half1Values.cfind(i);
261 
262  if (half1Fnd.found())
263  {
264  combine
265  (
266  pointValues,
267  cop,
268  meshPtsA[e[0]],
269  *half1Fnd
270  );
271  }
272  }
273  }
274  }
275  }
276 
277  // Synchronize multiple shared points.
278  if (pd.nGlobalPoints() > 0)
279  {
280  // meshPoint per local index
281  const labelList& sharedPtLabels = pd.sharedPointLabels();
282  // global shared index per local index
283  const labelList& sharedPtAddr = pd.sharedPointAddr();
284 
285  // Reduce on master.
286 
287  if (Pstream::parRun())
288  {
289  if (Pstream::master())
290  {
291  // Receive the edges using shared points from the slave.
292  for (const int slave : Pstream::subProcs())
293  {
294  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
295  Map<T> nbrValues(fromSlave);
296 
297  // Merge neighbouring values with my values
298  forAllConstIters(nbrValues, iter)
299  {
300  combine
301  (
302  sharedPointValues,
303  cop,
304  iter.key(), // edge
305  iter.val() // value
306  );
307  }
308  }
309 
310  // Send back
311  for (const int slave : Pstream::subProcs())
312  {
314  toSlave << sharedPointValues;
315  }
316  }
317  else
318  {
319  // Slave: send to master
320  {
321  OPstream toMaster
322  (
325  );
326  toMaster << sharedPointValues;
327  }
328  // Receive merged values
329  {
330  IPstream fromMaster
331  (
334  );
335  fromMaster >> sharedPointValues;
336  }
337  }
338  }
339 
340 
341  // Merge sharedPointValues (keyed on sharedPointAddr) into
342  // pointValues (keyed on mesh points).
343 
344  // Map from global shared index to meshpoint
345  Map<label> sharedToMeshPoint(2*sharedPtAddr.size());
346  forAll(sharedPtAddr, i)
347  {
348  sharedToMeshPoint.insert(sharedPtAddr[i], sharedPtLabels[i]);
349  }
350 
351  forAllConstIters(sharedToMeshPoint, iter)
352  {
353  // Do I have a value for my shared point
354  const auto sharedFnd = sharedPointValues.cfind(iter.key());
355 
356  if (sharedFnd.found())
357  {
358  pointValues.set(*iter, *sharedFnd);
359  }
360  }
361  }
362 }
363 
364 
365 template<class T, class CombineOp, class TransformOp>
367 (
368  const polyMesh& mesh,
369  EdgeMap<T>& edgeValues,
370  const CombineOp& cop,
371  const TransformOp& top
372 )
373 {
375 
376 
377  // Do synchronisation without constructing globalEdge addressing
378  // (since this constructs mesh edge addressing)
379 
380 
381  // Swap proc patch info
382  // ~~~~~~~~~~~~~~~~~~~~
383 
384  if (Pstream::parRun())
385  {
387 
388  // Send
389 
390  for (const polyPatch& pp : patches)
391  {
392  if (isA<processorPolyPatch>(pp) && pp.nEdges() > 0)
393  {
394  const processorPolyPatch& procPatch =
395  refCast<const processorPolyPatch>(pp);
396 
397 
398  // Get data per patch edge in neighbouring edge.
399 
400  const edgeList& edges = procPatch.edges();
401  const labelList& meshPts = procPatch.meshPoints();
402  const labelList& nbrPts = procPatch.neighbPoints();
403 
404  EdgeMap<T> patchInfo(edges.size() / 20);
405 
406  for (const edge& e : edges)
407  {
408  const edge meshEdge(meshPts[e[0]], meshPts[e[1]]);
409 
410  const auto iter = edgeValues.cfind(meshEdge);
411 
412  if (iter.found())
413  {
414  const edge nbrEdge(nbrPts[e[0]], nbrPts[e[1]]);
415  patchInfo.insert(nbrEdge, *iter);
416  }
417  }
418 
419  UOPstream toNeighb(procPatch.neighbProcNo(), pBufs);
420  toNeighb << patchInfo;
421  }
422  }
423 
424  pBufs.finishedSends();
425 
426  // Receive and combine.
427 
428  for (const polyPatch& pp : patches)
429  {
430  if (isA<processorPolyPatch>(pp) && pp.nEdges() > 0)
431  {
432  const processorPolyPatch& procPatch =
433  refCast<const processorPolyPatch>(pp);
434 
435  EdgeMap<T> nbrPatchInfo;
436  {
437  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
438  fromNbr >> nbrPatchInfo;
439  }
440 
441  // Apply transform to convert to this side properties.
442  top(procPatch, nbrPatchInfo);
443 
444 
445  // Only update those values which come from neighbour
446  const labelList& meshPts = procPatch.meshPoints();
447 
448  forAllConstIters(nbrPatchInfo, nbrIter)
449  {
450  const edge& e = nbrIter.key();
451  const edge meshEdge(meshPts[e[0]], meshPts[e[1]]);
452 
453  combine
454  (
455  edgeValues,
456  cop,
457  meshEdge, // edge
458  nbrIter.val() // value
459  );
460  }
461  }
462  }
463  }
464 
465 
466  // Swap cyclic info
467  // ~~~~~~~~~~~~~~~~
468 
469  for (const polyPatch& pp : patches)
470  {
471  if (isA<cyclicPolyPatch>(pp))
472  {
473  const cyclicPolyPatch& cycPatch =
474  refCast<const cyclicPolyPatch>(pp);
475 
476  if (cycPatch.owner())
477  {
478  // Owner does all.
479 
480  const edgeList& coupledEdges = cycPatch.coupledEdges();
481  const labelList& meshPtsA = cycPatch.meshPoints();
482  const edgeList& edgesA = cycPatch.edges();
483  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
484  const labelList& meshPtsB = nbrPatch.meshPoints();
485  const edgeList& edgesB = nbrPatch.edges();
486 
487  // Extract local values. Create map from edge to value.
488  Map<T> half0Values(edgesA.size() / 20);
489  Map<T> half1Values(half0Values.size());
490 
491  forAll(coupledEdges, i)
492  {
493  const edge& twoEdges = coupledEdges[i];
494 
495  {
496  const edge& e0 = edgesA[twoEdges[0]];
497  const edge meshEdge0(meshPtsA[e0[0]], meshPtsA[e0[1]]);
498 
499  const auto iter = edgeValues.cfind(meshEdge0);
500 
501  if (iter.found())
502  {
503  half0Values.insert(i, *iter);
504  }
505  }
506  {
507  const edge& e1 = edgesB[twoEdges[1]];
508  const edge meshEdge1(meshPtsB[e1[0]], meshPtsB[e1[1]]);
509 
510  const auto iter = edgeValues.cfind(meshEdge1);
511 
512  if (iter.found())
513  {
514  half1Values.insert(i, *iter);
515  }
516  }
517  }
518 
519  // Transform to this side
520  top(cycPatch, half1Values);
521  top(nbrPatch, half0Values);
522 
523 
524  // Extract and combine information
525 
526  forAll(coupledEdges, i)
527  {
528  const edge& twoEdges = coupledEdges[i];
529 
530  const auto half1Fnd = half1Values.cfind(i);
531 
532  if (half1Fnd.found())
533  {
534  const edge& e0 = edgesA[twoEdges[0]];
535  const edge meshEdge0(meshPtsA[e0[0]], meshPtsA[e0[1]]);
536 
537  combine
538  (
539  edgeValues,
540  cop,
541  meshEdge0, // edge
542  *half1Fnd // value
543  );
544  }
545 
546  const auto half0Fnd = half0Values.cfind(i);
547 
548  if (half0Fnd.found())
549  {
550  const edge& e1 = edgesB[twoEdges[1]];
551  const edge meshEdge1(meshPtsB[e1[0]], meshPtsB[e1[1]]);
552 
553  combine
554  (
555  edgeValues,
556  cop,
557  meshEdge1, // edge
558  *half0Fnd // value
559  );
560  }
561  }
562  }
563  }
564  }
565 
566  // Synchronize multiple shared points.
567  // Problem is that we don't want to construct shared edges so basically
568  // we do here like globalMeshData but then using sparse edge representation
569  // (EdgeMap instead of mesh.edges())
570 
571  const globalMeshData& pd = mesh.globalData();
572  const labelList& sharedPtAddr = pd.sharedPointAddr();
573  const labelList& sharedPtLabels = pd.sharedPointLabels();
574 
575  // 1. Create map from meshPoint to globalShared index.
576  Map<label> meshToShared(2*sharedPtLabels.size());
577  forAll(sharedPtLabels, i)
578  {
579  meshToShared.insert(sharedPtLabels[i], sharedPtAddr[i]);
580  }
581 
582  // Values on shared points. From two sharedPtAddr indices to a value.
583  EdgeMap<T> sharedEdgeValues(meshToShared.size());
584 
585  // From shared edge to mesh edge. Used for merging later on.
586  EdgeMap<edge> potentialSharedEdge(meshToShared.size());
587 
588  // 2. Find any edges using two global shared points. These will always be
589  // on the outside of the mesh. (though might not be on coupled patch
590  // if is single edge and not on coupled face)
591  // Store value (if any) on sharedEdgeValues
592  for (label facei = mesh.nInternalFaces(); facei < mesh.nFaces(); ++facei)
593  {
594  const face& f = mesh.faces()[facei];
595 
596  forAll(f, fp)
597  {
598  const label v0 = f[fp];
599  const label v1 = f[f.fcIndex(fp)];
600 
601  const auto v0Fnd = meshToShared.cfind(v0);
602 
603  if (v0Fnd.found())
604  {
605  const auto v1Fnd = meshToShared.cfind(v1);
606 
607  if (v1Fnd.found())
608  {
609  const edge meshEdge(v0, v1);
610 
611  // edge in shared point labels
612  const edge sharedEdge(*v0Fnd, *v1Fnd);
613 
614  // Store mesh edge as a potential shared edge.
615  potentialSharedEdge.insert(sharedEdge, meshEdge);
616 
617  const auto edgeFnd = edgeValues.cfind(meshEdge);
618 
619  if (edgeFnd.found())
620  {
621  // edge exists in edgeValues. See if already in map
622  // (since on same processor, e.g. cyclic)
623  combine
624  (
625  sharedEdgeValues,
626  cop,
627  sharedEdge, // edge
628  *edgeFnd // value
629  );
630  }
631  }
632  }
633  }
634  }
635 
636 
637  // Now sharedEdgeValues will contain per potential sharedEdge the value.
638  // (potential since an edge having two shared points is not necessary a
639  // shared edge).
640  // Reduce this on the master.
641 
642  if (Pstream::parRun())
643  {
644  if (Pstream::master())
645  {
646  // Receive the edges using shared points from the slave.
647  for (const int slave : Pstream::subProcs())
648  {
649  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
650  EdgeMap<T> nbrValues(fromSlave);
651 
652  // Merge neighbouring values with my values
653  forAllConstIters(nbrValues, iter)
654  {
655  combine
656  (
657  sharedEdgeValues,
658  cop,
659  iter.key(), // edge
660  iter.val() // value
661  );
662  }
663  }
664 
665  // Send back
666  for (const int slave : Pstream::subProcs())
667  {
669  toSlave << sharedEdgeValues;
670  }
671  }
672  else
673  {
674  // Send to master
675  {
676  OPstream toMaster
677  (
680  );
681  toMaster << sharedEdgeValues;
682  }
683  // Receive merged values
684  {
685  IPstream fromMaster
686  (
689  );
690  fromMaster >> sharedEdgeValues;
691  }
692  }
693  }
694 
695 
696  // Merge sharedEdgeValues (keyed on sharedPointAddr) into edgeValues
697  // (keyed on mesh points).
698 
699  // Loop over all my shared edges.
700  forAllConstIters(potentialSharedEdge, iter)
701  {
702  const edge& sharedEdge = iter.key();
703  const edge& meshEdge = iter.val();
704 
705  // Do I have a value for the shared edge?
706  const auto sharedFnd = sharedEdgeValues.cfind(sharedEdge);
707 
708  if (sharedFnd.found())
709  {
710  combine
711  (
712  edgeValues,
713  cop,
714  meshEdge, // edge
715  *sharedFnd // value
716  );
717  }
718  }
719 }
720 
721 
722 template<class T, class CombineOp, class TransformOp>
724 (
725  const polyMesh& mesh,
726  List<T>& pointValues,
727  const CombineOp& cop,
728  const T& nullValue,
729  const TransformOp& top
730 )
731 {
732  if (pointValues.size() != mesh.nPoints())
733  {
735  << "Number of values " << pointValues.size()
736  << " is not equal to the number of points in the mesh "
737  << mesh.nPoints() << abort(FatalError);
738  }
739 
740  mesh.globalData().syncPointData(pointValues, cop, top);
741 }
742 
743 
744 template<class T, class CombineOp, class TransformOp>
746 (
747  const polyMesh& mesh,
748  const labelUList& meshPoints,
749  List<T>& pointValues,
750  const CombineOp& cop,
751  const T& nullValue,
752  const TransformOp& top
753 )
754 {
755  if (pointValues.size() != meshPoints.size())
756  {
758  << "Number of values " << pointValues.size()
759  << " is not equal to the number of meshPoints "
760  << meshPoints.size() << abort(FatalError);
761  }
762  const globalMeshData& gd = mesh.globalData();
763  const indirectPrimitivePatch& cpp = gd.coupledPatch();
764  const Map<label>& mpm = cpp.meshPointMap();
765 
766  List<T> cppFld(cpp.nPoints(), nullValue);
767 
768  forAll(meshPoints, i)
769  {
770  const auto iter = mpm.cfind(meshPoints[i]);
771 
772  if (iter.found())
773  {
774  cppFld[*iter] = pointValues[i];
775  }
776  }
777 
779  (
780  cppFld,
781  gd.globalPointSlaves(),
784  gd.globalTransforms(),
785  cop,
786  top
787  );
788 
789  forAll(meshPoints, i)
790  {
791  const auto iter = mpm.cfind(meshPoints[i]);
792 
793  if (iter.found())
794  {
795  pointValues[i] = cppFld[*iter];
796  }
797  }
798 }
799 
800 
801 template<class T, class CombineOp, class TransformOp>
803 (
804  const polyMesh& mesh,
805  List<T>& edgeValues,
806  const CombineOp& cop,
807  const T& nullValue,
808  const TransformOp& top
809 )
810 {
811  if (edgeValues.size() != mesh.nEdges())
812  {
814  << "Number of values " << edgeValues.size()
815  << " is not equal to the number of edges in the mesh "
816  << mesh.nEdges() << abort(FatalError);
817  }
818 
819  const globalMeshData& gd = mesh.globalData();
820  const labelList& meshEdges = gd.coupledPatchMeshEdges();
821  const globalIndexAndTransform& git = gd.globalTransforms();
822  const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
823 
824  List<T> cppFld(UIndirectList<T>(edgeValues, meshEdges));
825 
827  (
828  cppFld,
829  gd.globalEdgeSlaves(),
831  edgeMap,
832  git,
833  cop,
834  top
835  );
836 
837  // Extract back onto mesh
838  forAll(meshEdges, i)
839  {
840  edgeValues[meshEdges[i]] = cppFld[i];
841  }
842 }
843 
844 
845 template<class T, class CombineOp, class TransformOp>
847 (
848  const polyMesh& mesh,
849  const labelList& meshEdges,
850  List<T>& edgeValues,
851  const CombineOp& cop,
852  const T& nullValue,
853  const TransformOp& top
854 )
855 {
856  if (edgeValues.size() != meshEdges.size())
857  {
859  << "Number of values " << edgeValues.size()
860  << " is not equal to the number of meshEdges "
861  << meshEdges.size() << abort(FatalError);
862  }
863  const globalMeshData& gd = mesh.globalData();
864  const indirectPrimitivePatch& cpp = gd.coupledPatch();
865  const Map<label>& mpm = gd.coupledPatchMeshEdgeMap();
866 
867  List<T> cppFld(cpp.nEdges(), nullValue);
868 
869  forAll(meshEdges, i)
870  {
871  const auto iter = mpm.cfind(meshEdges[i]);
872 
873  if (iter.found())
874  {
875  cppFld[*iter] = edgeValues[i];
876  }
877  }
878 
880  (
881  cppFld,
882  gd.globalEdgeSlaves(),
884  gd.globalEdgeSlavesMap(),
885  gd.globalTransforms(),
886  cop,
887  top
888  );
889 
890  forAll(meshEdges, i)
891  {
892  const auto iter = mpm.cfind(meshEdges[i]);
893 
894  if (iter.found())
895  {
896  edgeValues[i] = cppFld[*iter];
897  }
898  }
899 }
900 
901 
902 template<class T, class CombineOp, class TransformOp>
904 (
905  const polyMesh& mesh,
906  UList<T>& faceValues,
907  const CombineOp& cop,
908  const TransformOp& top,
909  const bool parRun
910 )
911 {
912  // Offset (global to local) for start of boundaries
913  const label boundaryOffset = mesh.nInternalFaces();
914 
915  if (faceValues.size() != mesh.nBoundaryFaces())
916  {
918  << "Number of values " << faceValues.size()
919  << " is not equal to the number of boundary faces in the mesh "
920  << mesh.nBoundaryFaces() << nl
921  << abort(FatalError);
922  }
923 
925 
926  if (parRun)
927  {
929 
930  // Send
931 
932  for (const polyPatch& pp : patches)
933  {
934  if (isA<processorPolyPatch>(pp) && pp.size() > 0)
935  {
936  const processorPolyPatch& procPatch =
937  refCast<const processorPolyPatch>(pp);
938 
939  const label patchStart = procPatch.start()-boundaryOffset;
940 
941  // Send slice of values on the patch
942  UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
943  toNbr<< SubList<T>(faceValues, procPatch.size(), patchStart);
944  }
945  }
946 
947 
948  pBufs.finishedSends();
949 
950 
951  // Receive and combine.
952 
953  for (const polyPatch& pp : patches)
954  {
955  if (isA<processorPolyPatch>(pp) && pp.size() > 0)
956  {
957  const processorPolyPatch& procPatch =
958  refCast<const processorPolyPatch>(pp);
959 
960  Field<T> nbrVals(procPatch.size());
961 
962  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
963  fromNbr >> nbrVals;
964 
965  top(procPatch, nbrVals);
966 
967  label bFacei = procPatch.start()-boundaryOffset;
968 
969  for (T& nbrVal : nbrVals)
970  {
971  cop(faceValues[bFacei++], nbrVal);
972  }
973  }
974  }
975  }
976 
977  // Do the cyclics.
978  for (const polyPatch& pp : patches)
979  {
980  if (isA<cyclicPolyPatch>(pp))
981  {
982  const cyclicPolyPatch& cycPatch =
983  refCast<const cyclicPolyPatch>(pp);
984 
985  if (cycPatch.owner())
986  {
987  // Owner does all.
988  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
989 
990  const label patchSize = cycPatch.size();
991  const label ownStart = cycPatch.start()-boundaryOffset;
992  const label nbrStart = nbrPatch.start()-boundaryOffset;
993 
994  // Transform (copy of) data on both sides
995  Field<T> ownVals(SubList<T>(faceValues, patchSize, ownStart));
996  top(nbrPatch, ownVals);
997 
998  Field<T> nbrVals(SubList<T>(faceValues, patchSize, nbrStart));
999  top(cycPatch, nbrVals);
1000 
1001  label bFacei = ownStart;
1002  for (T& nbrVal : nbrVals)
1003  {
1004  cop(faceValues[bFacei++], nbrVal);
1005  }
1006 
1007  bFacei = nbrStart;
1008  for (T& ownVal : ownVals)
1009  {
1010  cop(faceValues[bFacei++], ownVal);
1011  }
1012  }
1013  }
1014  }
1015 }
1016 
1017 
1018 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1019 
1020 template<unsigned Width, class CombineOp>
1023  const polyMesh& mesh,
1024  const bool isBoundaryOnly,
1025  PackedList<Width>& faceValues,
1026  const CombineOp& cop,
1027  const bool parRun
1028 )
1029 {
1030  // Offset (global to local) for start of boundaries
1031  const label boundaryOffset = (isBoundaryOnly ? mesh.nInternalFaces() : 0);
1032 
1033  if
1034  (
1035  faceValues.size()
1036  != (isBoundaryOnly ? mesh.nBoundaryFaces() : mesh.nFaces())
1037  )
1038  {
1040  << "Number of values " << faceValues.size()
1041  << " is not equal to the number of "
1042  << (isBoundaryOnly ? "boundary" : "mesh") << " faces "
1043  << (isBoundaryOnly ? mesh.nBoundaryFaces() : mesh.nFaces()) << nl
1044  << abort(FatalError);
1045  }
1046 
1048 
1049  if (parRun)
1050  {
1052 
1053  // Send
1054 
1055  for (const polyPatch& pp : patches)
1056  {
1057  if (isA<processorPolyPatch>(pp) && pp.size())
1058  {
1059  const processorPolyPatch& procPatch =
1060  refCast<const processorPolyPatch>(pp);
1061 
1062  const labelRange range
1063  (
1064  procPatch.start()-boundaryOffset,
1065  procPatch.size()
1066  );
1067 
1068  // Send slice of values on the patch
1069  UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
1070  toNbr<< PackedList<Width>(faceValues, range);
1071  }
1072  }
1073 
1074  pBufs.finishedSends();
1075 
1076 
1077  // Receive and combine.
1078 
1079  for (const polyPatch& pp : patches)
1080  {
1081  if (isA<processorPolyPatch>(pp) && pp.size())
1082  {
1083  const processorPolyPatch& procPatch =
1084  refCast<const processorPolyPatch>(pp);
1085 
1086  const label patchSize = procPatch.size();
1087 
1088  // Recv slice of values on the patch
1089  PackedList<Width> recvInfo(patchSize);
1090  {
1091  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
1092  fromNbr >> recvInfo;
1093  }
1094 
1095  // Combine (bitwise)
1096  label bFacei = procPatch.start()-boundaryOffset;
1097  for (label i = 0; i < patchSize; ++i)
1098  {
1099  unsigned int recvVal = recvInfo[i];
1100  unsigned int faceVal = faceValues[bFacei];
1101 
1102  cop(faceVal, recvVal);
1103  faceValues.set(bFacei, faceVal);
1104 
1105  ++bFacei;
1106  }
1107  }
1108  }
1109  }
1110 
1111 
1112  // Do the cyclics.
1113  for (const polyPatch& pp : patches)
1114  {
1115  if (isA<cyclicPolyPatch>(pp))
1116  {
1117  const cyclicPolyPatch& cycPatch =
1118  refCast<const cyclicPolyPatch>(pp);
1119 
1120  if (cycPatch.owner())
1121  {
1122  // Owner does all.
1123  const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
1124 
1125  const label patchSize = cycPatch.size();
1126 
1127  label face0 = cycPatch.start()-boundaryOffset;
1128  label face1 = nbrPatch.start()-boundaryOffset;
1129  for (label i = 0; i < patchSize; ++i)
1130  {
1131  unsigned int val0 = faceValues[face0];
1132  unsigned int val1 = faceValues[face1];
1133 
1134  unsigned int t = val0;
1135  cop(t, val1);
1136  faceValues[face0] = t;
1137 
1138  cop(val1, val0);
1139  faceValues[face1] = val1;
1140 
1141  ++face0;
1142  ++face1;
1143  }
1144  }
1145  }
1146  }
1147 }
1148 
1149 
1150 template<class T>
1153  const polyMesh& mesh,
1154  const UList<T>& cellData,
1155  List<T>& neighbourCellData
1156 )
1157 {
1158  if (cellData.size() != mesh.nCells())
1159  {
1161  << "Number of cell values " << cellData.size()
1162  << " is not equal to the number of cells in the mesh "
1163  << mesh.nCells() << abort(FatalError);
1164  }
1165 
1167 
1168  neighbourCellData.resize(mesh.nBoundaryFaces());
1169 
1170  for (const polyPatch& pp : patches)
1171  {
1172  label bFacei = pp.start()-mesh.nInternalFaces();
1173 
1174  const labelUList& faceCells = pp.faceCells();
1175 
1176  for (const label celli : faceCells)
1177  {
1178  neighbourCellData[bFacei] = cellData[celli];
1179  ++bFacei;
1180  }
1181  }
1182  syncTools::swapBoundaryFaceList(mesh, neighbourCellData);
1183 }
1184 
1185 
1186 template<unsigned Width, class CombineOp>
1189  const polyMesh& mesh,
1190  PackedList<Width>& faceValues,
1191  const CombineOp& cop,
1192  const bool parRun
1193 )
1194 {
1195  syncFaceList(mesh, false, faceValues, cop, parRun);
1196 }
1197 
1198 
1199 template<unsigned Width, class CombineOp>
1202  const polyMesh& mesh,
1203  PackedList<Width>& faceValues,
1204  const CombineOp& cop,
1205  const bool parRun
1206 )
1207 {
1208  syncFaceList(mesh, true, faceValues, cop, parRun);
1209 }
1210 
1211 
1212 template<unsigned Width>
1215  const polyMesh& mesh,
1216  PackedList<Width>& faceValues
1217 )
1218 {
1219  syncFaceList(mesh, faceValues, eqOp<unsigned int>());
1220 }
1221 
1222 
1223 template<unsigned Width>
1226  const polyMesh& mesh,
1227  PackedList<Width>& faceValues
1228 )
1229 {
1230  syncBoundaryFaceList(mesh, faceValues, eqOp<unsigned int>());
1231 }
1232 
1233 
1234 template<unsigned Width, class CombineOp>
1237  const polyMesh& mesh,
1238  PackedList<Width>& pointValues,
1239  const CombineOp& cop,
1240  const unsigned int nullValue
1241 )
1242 {
1243  if (pointValues.size() != mesh.nPoints())
1244  {
1246  << "Number of values " << pointValues.size()
1247  << " is not equal to the number of points in the mesh "
1248  << mesh.nPoints() << abort(FatalError);
1249  }
1250 
1251  const globalMeshData& gd = mesh.globalData();
1252  const labelList& meshPoints = gd.coupledPatch().meshPoints();
1253 
1255  forAll(meshPoints, i)
1256  {
1257  cppFld[i] = pointValues[meshPoints[i]];
1258  }
1259 
1261  (
1262  cppFld,
1263  gd.globalPointSlaves(),
1265  gd.globalPointSlavesMap(),
1266  cop
1267  );
1268 
1269  // Extract back to mesh
1270  forAll(meshPoints, i)
1271  {
1272  pointValues[meshPoints[i]] = cppFld[i];
1273  }
1274 }
1275 
1276 
1277 template<unsigned Width, class CombineOp>
1280  const polyMesh& mesh,
1281  PackedList<Width>& edgeValues,
1282  const CombineOp& cop,
1283  const unsigned int nullValue
1284 )
1285 {
1286  if (edgeValues.size() != mesh.nEdges())
1287  {
1289  << "Number of values " << edgeValues.size()
1290  << " is not equal to the number of edges in the mesh "
1291  << mesh.nEdges() << abort(FatalError);
1292  }
1293 
1294  const globalMeshData& gd = mesh.globalData();
1295  const labelList& meshEdges = gd.coupledPatchMeshEdges();
1296 
1298  forAll(meshEdges, i)
1299  {
1300  cppFld[i] = edgeValues[meshEdges[i]];
1301  }
1302 
1304  (
1305  cppFld,
1306  gd.globalEdgeSlaves(),
1308  gd.globalEdgeSlavesMap(),
1309  cop
1310  );
1311 
1312  // Extract back to mesh
1313  forAll(meshEdges, i)
1314  {
1315  edgeValues[meshEdges[i]] = cppFld[i];
1316  }
1317 }
1318 
1319 
1320 // ************************************************************************* //
Foam::globalMeshData::globalPointTransformedSlaves
const labelListList & globalPointTransformedSlaves() const
Definition: globalMeshData.C:2180
Foam::UPstream::commsTypes::nonBlocking
Foam::globalMeshData::globalPointSlaves
const labelListList & globalPointSlaves() const
Definition: globalMeshData.C:2170
Foam::UPstream::masterNo
static constexpr int masterNo() noexcept
Process index of the master (always 0)
Definition: UPstream.H:452
Foam::UOPstream
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:57
Foam::globalMeshData::coupledPatchMeshEdgeMap
const Map< label > & coupledPatchMeshEdgeMap() const
Return map from mesh edges to coupledPatch edges.
Definition: globalMeshData.C:2127
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:84
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:190
Foam::globalMeshData::globalEdgeSlavesMap
const mapDistribute & globalEdgeSlavesMap() const
Definition: globalMeshData.C:2245
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:322
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:904
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:87
Foam::UPstream::parRun
static bool & parRun()
Test if this a parallel run, or allow modify access.
Definition: UPstream.H:434
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: lumpedPointController.H:69
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:458
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:444
Foam::globalMeshData::globalPointSlavesMap
const mapDistribute & globalPointSlavesMap() const
Definition: globalMeshData.C:2191
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:724
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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::UPstream::subProcs
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:516
Foam::globalMeshData::nGlobalPoints
label nGlobalPoints() const
Return number of globally shared points.
Definition: globalMeshData.C:1986
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Field< T >
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
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:2160
Foam::globalMeshData::coupledPatchMeshEdges
const labelList & coupledPatchMeshEdges() const
Return map from coupledPatch edges to mesh edges.
Definition: globalMeshData.C:2107
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:803
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:55
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:58
Foam::globalMeshData::globalEdgeSlaves
const labelListList & globalEdgeSlaves() const
Definition: globalMeshData.C:2214
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:1152
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:316
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:107
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
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:313
Foam::mapDistributeBase::constructSize
label constructSize() const
Constructed data size.
Definition: mapDistributeBase.H:277
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:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
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:385
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::globalMeshData::globalEdgeTransformedSlaves
const labelListList & globalEdgeTransformedSlaves() const
Definition: globalMeshData.C:2224
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
Foam::List< label >
Foam::globalMeshData::coupledPatch
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
Definition: globalMeshData.C:2046
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:367
Foam::PackedList::size
label size() const noexcept
Number of entries.
Definition: PackedListI.H:377
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:1996
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
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:2006
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:72
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:1295
Foam::PrimitivePatch::meshPointMap
const Map< label > & meshPointMap() const
Mesh point map.
Definition: PrimitivePatch.C:323
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:310
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:85