syncToolsTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "syncTools.H"
27 #include "polyMesh.H"
28 #include "processorPolyPatch.H"
29 #include "cyclicPolyPatch.H"
30 #include "globalMeshData.H"
31 #include "contiguous.H"
32 #include "transform.H"
33 #include "SubField.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
37 template<class T, class CombineOp>
38 void Foam::syncTools::combine
39 (
40  Map<T>& pointValues,
41  const CombineOp& cop,
42  const label index,
43  const T& val
44 )
45 {
46  typename Map<T>::iterator iter = pointValues.find(index);
47 
48  if (iter != pointValues.end())
49  {
50  cop(iter(), val);
51  }
52  else
53  {
54  pointValues.insert(index, val);
55  }
56 }
57 
58 
59 template<class T, class CombineOp>
60 void Foam::syncTools::combine
61 (
62  EdgeMap<T>& edgeValues,
63  const CombineOp& cop,
64  const edge& index,
65  const T& val
66 )
67 {
68  typename EdgeMap<T>::iterator iter = edgeValues.find(index);
69 
70  if (iter != edgeValues.end())
71  {
72  cop(iter(), val);
73  }
74  else
75  {
76  edgeValues.insert(index, val);
77  }
78 }
79 
80 
81 template<class T, class CombineOp, class TransformOp>
83 (
84  const polyMesh& mesh,
85  Map<T>& pointValues, // from mesh point label to value
86  const CombineOp& cop,
87  const TransformOp& top
88 )
89 {
90  const polyBoundaryMesh& patches = mesh.boundaryMesh();
91 
92  // Synchronise multiple shared points.
93  const globalMeshData& pd = mesh.globalData();
94 
95  // Values on shared points. Keyed on global shared index.
96  Map<T> sharedPointValues(0);
97 
98  if (pd.nGlobalPoints() > 0)
99  {
100  // meshPoint per local index
101  const labelList& sharedPtLabels = pd.sharedPointLabels();
102  // global shared index per local index
103  const labelList& sharedPtAddr = pd.sharedPointAddr();
104 
105  sharedPointValues.resize(sharedPtAddr.size());
106 
107  // Fill my entries in the shared points
108  forAll(sharedPtLabels, i)
109  {
110  label meshPointi = sharedPtLabels[i];
111 
112  typename Map<T>::const_iterator fnd =
113  pointValues.find(meshPointi);
114 
115  if (fnd != pointValues.end())
116  {
117  combine
118  (
119  sharedPointValues,
120  cop,
121  sharedPtAddr[i], // index
122  fnd() // value
123  );
124  }
125  }
126  }
127 
128 
129  if (Pstream::parRun())
130  {
132 
133  // Send
134 
135  forAll(patches, patchi)
136  {
137  if
138  (
139  isA<processorPolyPatch>(patches[patchi])
140  && patches[patchi].nPoints() > 0
141  )
142  {
143  const processorPolyPatch& procPatch =
144  refCast<const processorPolyPatch>(patches[patchi]);
145 
146  // Get data per patchPoint in neighbouring point numbers.
147 
148  const labelList& meshPts = procPatch.meshPoints();
149  const labelList& nbrPts = procPatch.nbrPoints();
150 
151  // Extract local values. Create map from nbrPoint to value.
152  // Note: how small initial size?
153  Map<T> patchInfo(meshPts.size() / 20);
154 
155  forAll(meshPts, i)
156  {
157  typename Map<T>::const_iterator iter =
158  pointValues.find(meshPts[i]);
159 
160  if (iter != pointValues.end())
161  {
162  patchInfo.insert(nbrPts[i], iter());
163  }
164  }
165 
166  UOPstream toNeighb(procPatch.neighbProcNo(), pBufs);
167  toNeighb << patchInfo;
168  }
169  }
170 
171  pBufs.finishedSends();
172 
173  // Receive and combine.
174 
175  forAll(patches, patchi)
176  {
177  if
178  (
179  isA<processorPolyPatch>(patches[patchi])
180  && patches[patchi].nPoints() > 0
181  )
182  {
183  const processorPolyPatch& procPatch =
184  refCast<const processorPolyPatch>(patches[patchi]);
185 
186  UIPstream fromNb(procPatch.neighbProcNo(), pBufs);
187  Map<T> nbrPatchInfo(fromNb);
188 
189  // Transform
190  top(procPatch, nbrPatchInfo);
191 
192  const labelList& meshPts = procPatch.meshPoints();
193 
194  // Only update those values which come from neighbour
195  forAllConstIter(typename Map<T>, nbrPatchInfo, nbrIter)
196  {
197  combine
198  (
199  pointValues,
200  cop,
201  meshPts[nbrIter.key()],
202  nbrIter()
203  );
204  }
205  }
206  }
207  }
208 
209  // Do the cyclics.
210  forAll(patches, patchi)
211  {
212  if (isA<cyclicPolyPatch>(patches[patchi]))
213  {
214  const cyclicPolyPatch& cycPatch =
215  refCast<const cyclicPolyPatch>(patches[patchi]);
216 
217  if (cycPatch.owner())
218  {
219  // Owner does all.
220 
221  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
222  const edgeList& coupledPoints = cycPatch.coupledPoints();
223  const labelList& meshPtsA = cycPatch.meshPoints();
224  const labelList& meshPtsB = nbrPatch.meshPoints();
225 
226  // Extract local values. Create map from coupled-edge to value.
227  Map<T> half0Values(meshPtsA.size() / 20);
228  Map<T> half1Values(half0Values.size());
229 
230  forAll(coupledPoints, i)
231  {
232  const edge& e = coupledPoints[i];
233 
234  typename Map<T>::const_iterator point0Fnd =
235  pointValues.find(meshPtsA[e[0]]);
236 
237  if (point0Fnd != pointValues.end())
238  {
239  half0Values.insert(i, point0Fnd());
240  }
241 
242  typename Map<T>::const_iterator point1Fnd =
243  pointValues.find(meshPtsB[e[1]]);
244 
245  if (point1Fnd != pointValues.end())
246  {
247  half1Values.insert(i, point1Fnd());
248  }
249  }
250 
251  // Transform to receiving side
252  top(cycPatch, half1Values);
253  top(nbrPatch, half0Values);
254 
255  forAll(coupledPoints, i)
256  {
257  const edge& e = coupledPoints[i];
258 
259  typename Map<T>::const_iterator half0Fnd =
260  half0Values.find(i);
261 
262  if (half0Fnd != half0Values.end())
263  {
264  combine
265  (
266  pointValues,
267  cop,
268  meshPtsB[e[1]],
269  half0Fnd()
270  );
271  }
272 
273  typename Map<T>::const_iterator half1Fnd =
274  half1Values.find(i);
275 
276  if (half1Fnd != half1Values.end())
277  {
278  combine
279  (
280  pointValues,
281  cop,
282  meshPtsA[e[0]],
283  half1Fnd()
284  );
285  }
286  }
287  }
288  }
289  }
290 
291  // Synchronise multiple shared points.
292  if (pd.nGlobalPoints() > 0)
293  {
294  // meshPoint per local index
295  const labelList& sharedPtLabels = pd.sharedPointLabels();
296  // global shared index per local index
297  const labelList& sharedPtAddr = pd.sharedPointAddr();
298 
299  // Reduce on master.
300 
301  if (Pstream::parRun())
302  {
303  if (Pstream::master())
304  {
305  // Receive the edges using shared points from the slave.
306  for
307  (
308  int slave=Pstream::firstSlave();
309  slave<=Pstream::lastSlave();
310  slave++
311  )
312  {
313  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
314  Map<T> nbrValues(fromSlave);
315 
316  // Merge neighbouring values with my values
317  forAllConstIter(typename Map<T>, nbrValues, iter)
318  {
319  combine
320  (
321  sharedPointValues,
322  cop,
323  iter.key(), // edge
324  iter() // value
325  );
326  }
327  }
328 
329  // Send back
330  for
331  (
332  int slave=Pstream::firstSlave();
333  slave<=Pstream::lastSlave();
334  slave++
335  )
336  {
338  toSlave << sharedPointValues;
339  }
340  }
341  else
342  {
343  // Slave: send to master
344  {
345  OPstream toMaster
346  (
349  );
350  toMaster << sharedPointValues;
351  }
352  // Receive merged values
353  {
354  IPstream fromMaster
355  (
358  );
359  fromMaster >> sharedPointValues;
360  }
361  }
362  }
363 
364 
365  // Merge sharedPointValues (keyed on sharedPointAddr) into
366  // pointValues (keyed on mesh points).
367 
368  // Map from global shared index to meshpoint
369  Map<label> sharedToMeshPoint(2*sharedPtAddr.size());
370  forAll(sharedPtAddr, i)
371  {
372  sharedToMeshPoint.insert(sharedPtAddr[i], sharedPtLabels[i]);
373  }
374 
375  forAllConstIter(Map<label>, sharedToMeshPoint, iter)
376  {
377  // Do I have a value for my shared point
378  typename Map<T>::const_iterator sharedFnd =
379  sharedPointValues.find(iter.key());
380 
381  if (sharedFnd != sharedPointValues.end())
382  {
383  pointValues.set(iter(), sharedFnd());
384  }
385  }
386  }
387 }
388 
389 
390 template<class T, class CombineOp, class TransformOp>
392 (
393  const polyMesh& mesh,
394  EdgeMap<T>& edgeValues,
395  const CombineOp& cop,
396  const TransformOp& top
397 )
398 {
399  const polyBoundaryMesh& patches = mesh.boundaryMesh();
400 
401 
402  // Do synchronisation without constructing globalEdge addressing
403  // (since this constructs mesh edge addressing)
404 
405 
406  // Swap proc patch info
407  // ~~~~~~~~~~~~~~~~~~~~
408 
409  if (Pstream::parRun())
410  {
412 
413  // Send
414 
415  forAll(patches, patchi)
416  {
417  if
418  (
419  isA<processorPolyPatch>(patches[patchi])
420  && patches[patchi].nEdges() > 0
421  )
422  {
423  const processorPolyPatch& procPatch =
424  refCast<const processorPolyPatch>(patches[patchi]);
425 
426 
427  // Get data per patch edge in neighbouring edge.
428 
429  const edgeList& edges = procPatch.edges();
430  const labelList& meshPts = procPatch.meshPoints();
431  const labelList& nbrPts = procPatch.nbrPoints();
432 
433  EdgeMap<T> patchInfo(edges.size() / 20);
434 
435  forAll(edges, i)
436  {
437  const edge& e = edges[i];
438  const edge meshEdge(meshPts[e[0]], meshPts[e[1]]);
439 
440  typename EdgeMap<T>::const_iterator iter =
441  edgeValues.find(meshEdge);
442 
443  if (iter != edgeValues.end())
444  {
445  const edge nbrEdge(nbrPts[e[0]], nbrPts[e[1]]);
446  patchInfo.insert(nbrEdge, iter());
447  }
448  }
449 
450  UOPstream toNeighb(procPatch.neighbProcNo(), pBufs);
451  toNeighb << patchInfo;
452  }
453  }
454 
455  pBufs.finishedSends();
456 
457  // Receive and combine.
458 
459  forAll(patches, patchi)
460  {
461  if
462  (
463  isA<processorPolyPatch>(patches[patchi])
464  && patches[patchi].nEdges() > 0
465  )
466  {
467  const processorPolyPatch& procPatch =
468  refCast<const processorPolyPatch>(patches[patchi]);
469 
470  EdgeMap<T> nbrPatchInfo;
471  {
472  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
473  fromNbr >> nbrPatchInfo;
474  }
475 
476  // Apply transform to convert to this side properties.
477  top(procPatch, nbrPatchInfo);
478 
479 
480  // Only update those values which come from neighbour
481  const labelList& meshPts = procPatch.meshPoints();
482 
483  forAllConstIter(typename EdgeMap<T>, nbrPatchInfo, nbrIter)
484  {
485  const edge& e = nbrIter.key();
486  const edge meshEdge(meshPts[e[0]], meshPts[e[1]]);
487 
488  combine
489  (
490  edgeValues,
491  cop,
492  meshEdge, // edge
493  nbrIter() // value
494  );
495  }
496  }
497  }
498  }
499 
500 
501  // Swap cyclic info
502  // ~~~~~~~~~~~~~~~~
503 
504  forAll(patches, patchi)
505  {
506  if (isA<cyclicPolyPatch>(patches[patchi]))
507  {
508  const cyclicPolyPatch& cycPatch =
509  refCast<const cyclicPolyPatch>(patches[patchi]);
510 
511  if (cycPatch.owner())
512  {
513  // Owner does all.
514 
515  const edgeList& coupledEdges = cycPatch.coupledEdges();
516  const labelList& meshPtsA = cycPatch.meshPoints();
517  const edgeList& edgesA = cycPatch.edges();
518  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
519  const labelList& meshPtsB = nbrPatch.meshPoints();
520  const edgeList& edgesB = nbrPatch.edges();
521 
522  // Extract local values. Create map from edge to value.
523  Map<T> half0Values(edgesA.size() / 20);
524  Map<T> half1Values(half0Values.size());
525 
526  forAll(coupledEdges, i)
527  {
528  const edge& twoEdges = coupledEdges[i];
529 
530  {
531  const edge& e0 = edgesA[twoEdges[0]];
532  const edge meshEdge0(meshPtsA[e0[0]], meshPtsA[e0[1]]);
533 
534  typename EdgeMap<T>::const_iterator iter =
535  edgeValues.find(meshEdge0);
536 
537  if (iter != edgeValues.end())
538  {
539  half0Values.insert(i, iter());
540  }
541  }
542  {
543  const edge& e1 = edgesB[twoEdges[1]];
544  const edge meshEdge1(meshPtsB[e1[0]], meshPtsB[e1[1]]);
545 
546  typename EdgeMap<T>::const_iterator iter =
547  edgeValues.find(meshEdge1);
548 
549  if (iter != edgeValues.end())
550  {
551  half1Values.insert(i, iter());
552  }
553  }
554  }
555 
556  // Transform to this side
557  top(cycPatch, half1Values);
558  top(nbrPatch, half0Values);
559 
560 
561  // Extract and combine information
562 
563  forAll(coupledEdges, i)
564  {
565  const edge& twoEdges = coupledEdges[i];
566 
567  typename Map<T>::const_iterator half1Fnd =
568  half1Values.find(i);
569 
570  if (half1Fnd != half1Values.end())
571  {
572  const edge& e0 = edgesA[twoEdges[0]];
573  const edge meshEdge0(meshPtsA[e0[0]], meshPtsA[e0[1]]);
574 
575  combine
576  (
577  edgeValues,
578  cop,
579  meshEdge0, // edge
580  half1Fnd() // value
581  );
582  }
583 
584  typename Map<T>::const_iterator half0Fnd =
585  half0Values.find(i);
586  if (half0Fnd != half0Values.end())
587  {
588  const edge& e1 = edgesB[twoEdges[1]];
589  const edge meshEdge1(meshPtsB[e1[0]], meshPtsB[e1[1]]);
590 
591  combine
592  (
593  edgeValues,
594  cop,
595  meshEdge1, // edge
596  half0Fnd() // value
597  );
598  }
599  }
600  }
601  }
602  }
603 
604  // Synchronise multiple shared points.
605  // Problem is that we don't want to construct shared edges so basically
606  // we do here like globalMeshData but then using sparse edge representation
607  // (EdgeMap instead of mesh.edges())
608 
609  const globalMeshData& pd = mesh.globalData();
610  const labelList& sharedPtAddr = pd.sharedPointAddr();
611  const labelList& sharedPtLabels = pd.sharedPointLabels();
612 
613  // 1. Create map from meshPoint to globalShared index.
614  Map<label> meshToShared(2*sharedPtLabels.size());
615  forAll(sharedPtLabels, i)
616  {
617  meshToShared.insert(sharedPtLabels[i], sharedPtAddr[i]);
618  }
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  label v0 = f[fp];
637  label v1 = f[f.fcIndex(fp)];
638 
639  Map<label>::const_iterator v0Fnd = meshToShared.find(v0);
640 
641  if (v0Fnd != meshToShared.end())
642  {
643  Map<label>::const_iterator v1Fnd = meshToShared.find(v1);
644 
645  if (v1Fnd != meshToShared.end())
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  typename EdgeMap<T>::const_iterator edgeFnd =
656  edgeValues.find(meshEdge);
657 
658  if (edgeFnd != edgeValues.end())
659  {
660  // edge exists in edgeValues. See if already in map
661  // (since on same processor, e.g. cyclic)
662  combine
663  (
664  sharedEdgeValues,
665  cop,
666  sharedEdge, // edge
667  edgeFnd() // value
668  );
669  }
670  }
671  }
672  }
673  }
674 
675 
676  // Now sharedEdgeValues will contain per potential sharedEdge the value.
677  // (potential since an edge having two shared points is not necessary a
678  // shared edge).
679  // Reduce this on the master.
680 
681  if (Pstream::parRun())
682  {
683  if (Pstream::master())
684  {
685  // Receive the edges using shared points from the slave.
686  for
687  (
688  int slave=Pstream::firstSlave();
689  slave<=Pstream::lastSlave();
690  slave++
691  )
692  {
693  IPstream fromSlave(Pstream::commsTypes::scheduled, slave);
694  EdgeMap<T> nbrValues(fromSlave);
695 
696  // Merge neighbouring values with my values
697  forAllConstIter(typename EdgeMap<T>, nbrValues, iter)
698  {
699  combine
700  (
701  sharedEdgeValues,
702  cop,
703  iter.key(), // edge
704  iter() // value
705  );
706  }
707  }
708 
709  // Send back
710  for
711  (
712  int slave=Pstream::firstSlave();
713  slave<=Pstream::lastSlave();
714  slave++
715  )
716  {
717 
719  toSlave << sharedEdgeValues;
720  }
721  }
722  else
723  {
724  // Send to master
725  {
726  OPstream toMaster
727  (
730  );
731  toMaster << sharedEdgeValues;
732  }
733  // Receive merged values
734  {
735  IPstream fromMaster
736  (
739  );
740  fromMaster >> sharedEdgeValues;
741  }
742  }
743  }
744 
745 
746  // Merge sharedEdgeValues (keyed on sharedPointAddr) into edgeValues
747  // (keyed on mesh points).
748 
749  // Loop over all my shared edges.
750  forAllConstIter(typename EdgeMap<edge>, potentialSharedEdge, iter)
751  {
752  const edge& sharedEdge = iter.key();
753  const edge& meshEdge = iter();
754 
755  // Do I have a value for the shared edge?
756  typename EdgeMap<T>::const_iterator sharedFnd =
757  sharedEdgeValues.find(sharedEdge);
758 
759  if (sharedFnd != sharedEdgeValues.end())
760  {
761  combine
762  (
763  edgeValues,
764  cop,
765  meshEdge, // edge
766  sharedFnd() // value
767  );
768  }
769  }
770 }
771 
772 
773 //template<class T, class CombineOp, class TransformOp>
774 //void Foam::syncTools::syncPointList
775 //(
776 // const polyMesh& mesh,
777 // List<T>& pointValues,
778 // const CombineOp& cop,
779 // const T& nullValue,
780 // const TransformOp& top
781 //)
782 //{
783 // if (pointValues.size() != mesh.nPoints())
784 // {
785 // FatalErrorInFunction
786 // << "Number of values " << pointValues.size()
787 // << " is not equal to the number of points in the mesh "
788 // << mesh.nPoints() << abort(FatalError);
789 // }
790 //
791 // const polyBoundaryMesh& patches = mesh.boundaryMesh();
792 //
793 // // Synchronise multiple shared points.
794 // const globalMeshData& pd = mesh.globalData();
795 //
796 // // Values on shared points.
797 // Field<T> sharedPts(0);
798 // if (pd.nGlobalPoints() > 0)
799 // {
800 // // Values on shared points.
801 // sharedPts.setSize(pd.nGlobalPoints(), nullValue);
802 //
803 // forAll(pd.sharedPointLabels(), i)
804 // {
805 // label meshPointi = pd.sharedPointLabels()[i];
806 // // Fill my entries in the shared points
807 // sharedPts[pd.sharedPointAddr()[i]] = pointValues[meshPointi];
808 // }
809 // }
810 //
811 // if (Pstream::parRun())
812 // {
813 // PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking);
814 //
815 // // Send
816 //
817 // forAll(patches, patchi)
818 // {
819 // if
820 // (
821 // isA<processorPolyPatch>(patches[patchi])
822 // && patches[patchi].nPoints() > 0
823 // )
824 // {
825 // const processorPolyPatch& procPatch =
826 // refCast<const processorPolyPatch>(patches[patchi]);
827 //
828 // // Get data per patchPoint in neighbouring point numbers.
829 // Field<T> patchInfo(procPatch.nPoints());
830 //
831 // const labelList& meshPts = procPatch.meshPoints();
832 // const labelList& nbrPts = procPatch.nbrPoints();
833 //
834 // forAll(nbrPts, pointi)
835 // {
836 // label nbrPointi = nbrPts[pointi];
837 // patchInfo[nbrPointi] = pointValues[meshPts[pointi]];
838 // }
839 //
840 // UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
841 // toNbr << patchInfo;
842 // }
843 // }
844 //
845 // pBufs.finishedSends();
846 //
847 // // Receive and combine.
848 //
849 // forAll(patches, patchi)
850 // {
851 // if
852 // (
853 // isA<processorPolyPatch>(patches[patchi])
854 // && patches[patchi].nPoints() > 0
855 // )
856 // {
857 // const processorPolyPatch& procPatch =
858 // refCast<const processorPolyPatch>(patches[patchi]);
859 //
860 // Field<T> nbrPatchInfo(procPatch.nPoints());
861 // {
862 // UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
863 // fromNbr >> nbrPatchInfo;
864 // }
865 //
866 // // Transform to this side
867 // top(procPatch, nbrPatchInfo);
868 //
869 // const labelList& meshPts = procPatch.meshPoints();
870 //
871 // forAll(meshPts, pointi)
872 // {
873 // label meshPointi = meshPts[pointi];
874 // cop(pointValues[meshPointi], nbrPatchInfo[pointi]);
875 // }
876 // }
877 // }
878 // }
879 //
880 // // Do the cyclics.
881 // forAll(patches, patchi)
882 // {
883 // if (isA<cyclicPolyPatch>(patches[patchi]))
884 // {
885 // const cyclicPolyPatch& cycPatch =
886 // refCast<const cyclicPolyPatch>(patches[patchi]);
887 //
888 // if (cycPatch.owner())
889 // {
890 // // Owner does all.
891 //
892 // const edgeList& coupledPoints = cycPatch.coupledPoints();
893 // const labelList& meshPts = cycPatch.meshPoints();
894 // const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
895 // const labelList& nbrMeshPoints = nbrPatch.meshPoints();
896 //
897 // Field<T> half0Values(coupledPoints.size());
898 // Field<T> half1Values(coupledPoints.size());
899 //
900 // forAll(coupledPoints, i)
901 // {
902 // const edge& e = coupledPoints[i];
903 // half0Values[i] = pointValues[meshPts[e[0]]];
904 // half1Values[i] = pointValues[nbrMeshPoints[e[1]]];
905 // }
906 //
907 // // SubField<T> slice0(half0Values, half0Values.size());
908 // // SubField<T> slice1(half1Values, half1Values.size());
909 // // top(cycPatch, reinterpret_cast<Field<T>&>(slice1));
910 // // top(nbrPatch, reinterpret_cast<Field<T>&>(slice0));
911 //
912 // top(cycPatch, half1Values);
913 // top(nbrPatch, half0Values);
914 //
915 // forAll(coupledPoints, i)
916 // {
917 // const edge& e = coupledPoints[i];
918 // cop(pointValues[meshPts[e[0]]], half1Values[i]);
919 // cop(pointValues[nbrMeshPoints[e[1]]], half0Values[i]);
920 // }
921 // }
922 // }
923 // }
924 //
925 // // Synchronise multiple shared points.
926 // const globalMeshData& pd = mesh.globalData();
927 //
928 // if (pd.nGlobalPoints() > 0)
929 // {
930 // // Combine on master.
931 // Pstream::listCombineGather(sharedPts, cop);
932 // Pstream::listCombineScatter(sharedPts);
933 //
934 // // Now we will all have the same information. Merge it back with
935 // // my local information.
936 // forAll(pd.sharedPointLabels(), i)
937 // {
938 // label meshPointi = pd.sharedPointLabels()[i];
939 // pointValues[meshPointi] = sharedPts[pd.sharedPointAddr()[i]];
940 // }
941 // }
942 //}
943 
944 
945 //template<class T, class CombineOp, class TransformOp>
946 //void Foam::syncTools::syncPointList
947 //(
948 // const polyMesh& mesh,
949 // const labelList& meshPoints,
950 // List<T>& pointValues,
951 // const CombineOp& cop,
952 // const T& nullValue,
953 // const TransformOp& top
954 //)
955 //{
956 // if (pointValues.size() != meshPoints.size())
957 // {
958 // FatalErrorInFunction
959 // << "Number of values " << pointValues.size()
960 // << " is not equal to the number of points "
961 // << meshPoints.size() << abort(FatalError);
962 // }
963 //
964 // if (!hasCouples(mesh.boundaryMesh()))
965 // {
966 // return;
967 // }
968 //
969 // Field<T> meshValues(mesh.nPoints(), nullValue);
970 //
971 // forAll(meshPoints, i)
972 // {
973 // meshValues[meshPoints[i]] = pointValues[i];
974 // }
975 //
976 // syncTools::syncPointList
977 // (
978 // mesh,
979 // meshValues,
980 // cop, // combine op
981 // nullValue, // null value
982 // top // position or field
983 // );
984 //
985 // forAll(meshPoints, i)
986 // {
987 // pointValues[i] = meshValues[meshPoints[i]];
988 // }
989 //}
990 
991 template<class T, class CombineOp, class TransformOp>
993 (
994  const polyMesh& mesh,
995  List<T>& pointValues,
996  const CombineOp& cop,
997  const T& nullValue,
998  const TransformOp& top
999 )
1000 {
1001  if (pointValues.size() != mesh.nPoints())
1002  {
1004  << "Number of values " << pointValues.size()
1005  << " is not equal to the number of points in the mesh "
1006  << mesh.nPoints() << abort(FatalError);
1007  }
1008 
1009  mesh.globalData().syncPointData(pointValues, cop, top);
1010 }
1011 
1012 
1013 //template<class CombineOp>
1014 //void Foam::syncTools::syncPointPositions
1015 //(
1016 // const polyMesh& mesh,
1017 // List<point>& pointValues,
1018 // const CombineOp& cop,
1019 // const point& nullValue
1020 //)
1021 //{
1022 // if (pointValues.size() != mesh.nPoints())
1023 // {
1024 // FatalErrorInFunction
1025 // << "Number of values " << pointValues.size()
1026 // << " is not equal to the number of points in the mesh "
1027 // << mesh.nPoints() << abort(FatalError);
1028 // }
1029 //
1030 // mesh.globalData().syncPointData(pointValues, cop, true);
1031 //}
1032 
1033 
1034 template<class T, class CombineOp, class TransformOp>
1037  const polyMesh& mesh,
1038  const labelList& meshPoints,
1039  List<T>& pointValues,
1040  const CombineOp& cop,
1041  const T& nullValue,
1042  const TransformOp& top
1043 )
1044 {
1045  if (pointValues.size() != meshPoints.size())
1046  {
1048  << "Number of values " << pointValues.size()
1049  << " is not equal to the number of meshPoints "
1050  << meshPoints.size() << abort(FatalError);
1051  }
1052  const globalMeshData& gd = mesh.globalData();
1053  const indirectPrimitivePatch& cpp = gd.coupledPatch();
1054  const Map<label>& mpm = cpp.meshPointMap();
1055 
1056  List<T> cppFld(cpp.nPoints(), nullValue);
1057 
1058  forAll(meshPoints, i)
1059  {
1060  label pointi = meshPoints[i];
1061  Map<label>::const_iterator iter = mpm.find(pointi);
1062  if (iter != mpm.end())
1063  {
1064  cppFld[iter()] = pointValues[i];
1065  }
1066  }
1067 
1069  (
1070  cppFld,
1071  gd.globalPointSlaves(),
1073  gd.globalPointSlavesMap(),
1074  gd.globalTransforms(),
1075  cop,
1076  top
1077  );
1078 
1079  forAll(meshPoints, i)
1080  {
1081  label pointi = meshPoints[i];
1082  Map<label>::const_iterator iter = mpm.find(pointi);
1083  if (iter != mpm.end())
1084  {
1085  pointValues[i] = cppFld[iter()];
1086  }
1087  }
1088 }
1089 
1090 
1091 //template<class CombineOp>
1092 //void Foam::syncTools::syncPointPositions
1093 //(
1094 // const polyMesh& mesh,
1095 // const labelList& meshPoints,
1096 // List<point>& pointValues,
1097 // const CombineOp& cop,
1098 // const point& nullValue
1099 //)
1100 //{
1101 // if (pointValues.size() != meshPoints.size())
1102 // {
1103 // FatalErrorInFunction
1104 // << "Number of values " << pointValues.size()
1105 // << " is not equal to the number of meshPoints "
1106 // << meshPoints.size() << abort(FatalError);
1107 // }
1108 // const globalMeshData& gd = mesh.globalData();
1109 // const indirectPrimitivePatch& cpp = gd.coupledPatch();
1110 // const Map<label>& mpm = cpp.meshPointMap();
1111 //
1112 // List<point> cppFld(cpp.nPoints(), nullValue);
1113 //
1114 // forAll(meshPoints, i)
1115 // {
1116 // label pointi = meshPoints[i];
1117 // Map<label>::const_iterator iter = mpm.find(pointi);
1118 // if (iter != mpm.end())
1119 // {
1120 // cppFld[iter()] = pointValues[i];
1121 // }
1122 // }
1123 //
1124 // globalMeshData::syncData
1125 // (
1126 // cppFld,
1127 // gd.globalPointSlaves(),
1128 // gd.globalPointTransformedSlaves(),
1129 // gd.globalPointSlavesMap(),
1130 // gd.globalTransforms(),
1131 // cop,
1132 // true, // position?
1133 // mapDistribute::transform() // not used
1134 // );
1135 //
1136 // forAll(meshPoints, i)
1137 // {
1138 // label pointi = meshPoints[i];
1139 // Map<label>::const_iterator iter = mpm.find(pointi);
1140 // if (iter != mpm.end())
1141 // {
1142 // pointValues[i] = cppFld[iter()];
1143 // }
1144 // }
1145 //}
1146 
1147 
1148 template<class T, class CombineOp, class TransformOp>
1151  const polyMesh& mesh,
1152  List<T>& edgeValues,
1153  const CombineOp& cop,
1154  const T& nullValue,
1155  const TransformOp& top
1156 )
1157 {
1158  if (edgeValues.size() != mesh.nEdges())
1159  {
1161  << "Number of values " << edgeValues.size()
1162  << " is not equal to the number of edges in the mesh "
1163  << mesh.nEdges() << abort(FatalError);
1164  }
1165 
1166  const globalMeshData& gd = mesh.globalData();
1167  const labelList& meshEdges = gd.coupledPatchMeshEdges();
1168  const globalIndexAndTransform& git = gd.globalTransforms();
1169  const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
1170 
1171  List<T> cppFld(UIndirectList<T>(edgeValues, meshEdges));
1172 
1174  (
1175  cppFld,
1176  gd.globalEdgeSlaves(),
1178  edgeMap,
1179  git,
1180  cop,
1181  top
1182  );
1183 
1184  // Extract back onto mesh
1185  forAll(meshEdges, i)
1186  {
1187  edgeValues[meshEdges[i]] = cppFld[i];
1188  }
1189 }
1190 
1191 
1192 //template<class CombineOp>
1193 //void Foam::syncTools::syncEdgePositions
1194 //(
1195 // const polyMesh& mesh,
1196 // List<point>& edgeValues,
1197 // const CombineOp& cop,
1198 // const point& nullValue
1199 //)
1200 //{
1201 // if (edgeValues.size() != mesh.nEdges())
1202 // {
1203 // FatalErrorInFunction
1204 // << "Number of values " << edgeValues.size()
1205 // << " is not equal to the number of edges in the mesh "
1206 // << mesh.nEdges() << abort(FatalError);
1207 // }
1208 //
1209 // const globalMeshData& gd = mesh.globalData();
1210 // const labelList& meshEdges = gd.coupledPatchMeshEdges();
1211 // const globalIndexAndTransform& git = gd.globalTransforms();
1212 // const mapDistribute& map = gd.globalEdgeSlavesMap();
1213 //
1214 // List<point> cppFld(UIndirectList<point>(edgeValues, meshEdges));
1215 //
1216 // globalMeshData::syncData
1217 // (
1218 // cppFld,
1219 // gd.globalEdgeSlaves(),
1220 // gd.globalEdgeTransformedSlaves(),
1221 // map,
1222 // git,
1223 // cop,
1224 // true, // position?
1225 // mapDistribute::transform() // not used
1226 // );
1227 //
1228 // // Extract back onto mesh
1229 // forAll(meshEdges, i)
1230 // {
1231 // edgeValues[meshEdges[i]] = cppFld[i];
1232 // }
1233 //}
1234 
1235 
1236 template<class T, class CombineOp, class TransformOp>
1239  const polyMesh& mesh,
1240  const labelList& meshEdges,
1241  List<T>& edgeValues,
1242  const CombineOp& cop,
1243  const T& nullValue,
1244  const TransformOp& top
1245 )
1246 {
1247  if (edgeValues.size() != meshEdges.size())
1248  {
1250  << "Number of values " << edgeValues.size()
1251  << " is not equal to the number of meshEdges "
1252  << meshEdges.size() << abort(FatalError);
1253  }
1254  const globalMeshData& gd = mesh.globalData();
1255  const indirectPrimitivePatch& cpp = gd.coupledPatch();
1256  const Map<label>& mpm = gd.coupledPatchMeshEdgeMap();
1257 
1258  List<T> cppFld(cpp.nEdges(), nullValue);
1259 
1260  forAll(meshEdges, i)
1261  {
1262  label edgeI = meshEdges[i];
1263  Map<label>::const_iterator iter = mpm.find(edgeI);
1264  if (iter != mpm.end())
1265  {
1266  cppFld[iter()] = edgeValues[i];
1267  }
1268  }
1269 
1271  (
1272  cppFld,
1273  gd.globalEdgeSlaves(),
1275  gd.globalEdgeSlavesMap(),
1276  gd.globalTransforms(),
1277  cop,
1278  top
1279  );
1280 
1281  forAll(meshEdges, i)
1282  {
1283  label edgeI = meshEdges[i];
1284  Map<label>::const_iterator iter = mpm.find(edgeI);
1285  if (iter != mpm.end())
1286  {
1287  edgeValues[i] = cppFld[iter()];
1288  }
1289  }
1290 }
1291 
1292 template<class T, class CombineOp, class TransformOp>
1295  const polyMesh& mesh,
1296  UList<T>& faceValues,
1297  const CombineOp& cop,
1298  const TransformOp& top,
1299  const bool parRun
1300 )
1301 {
1302  const label nBFaces = mesh.nFaces() - mesh.nInternalFaces();
1303 
1304  if (faceValues.size() != nBFaces)
1305  {
1307  << "Number of values " << faceValues.size()
1308  << " is not equal to the number of boundary faces in the mesh "
1309  << nBFaces << abort(FatalError);
1310  }
1311 
1312  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1313 
1314  if (parRun)
1315  {
1317 
1318  // Send
1319 
1320  forAll(patches, patchi)
1321  {
1322  if
1323  (
1324  isA<processorPolyPatch>(patches[patchi])
1325  && patches[patchi].size() > 0
1326  )
1327  {
1328  const processorPolyPatch& procPatch =
1329  refCast<const processorPolyPatch>(patches[patchi]);
1330 
1331  label patchStart = procPatch.start()-mesh.nInternalFaces();
1332 
1333  UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
1334  toNbr << SubField<T>(faceValues, procPatch.size(), patchStart);
1335  }
1336  }
1337 
1338 
1339  pBufs.finishedSends();
1340 
1341 
1342  // Receive and combine.
1343 
1344  forAll(patches, patchi)
1345  {
1346  if
1347  (
1348  isA<processorPolyPatch>(patches[patchi])
1349  && patches[patchi].size() > 0
1350  )
1351  {
1352  const processorPolyPatch& procPatch =
1353  refCast<const processorPolyPatch>(patches[patchi]);
1354 
1355  Field<T> nbrPatchInfo(procPatch.size());
1356 
1357  UIPstream fromNeighb(procPatch.neighbProcNo(), pBufs);
1358  fromNeighb >> nbrPatchInfo;
1359 
1360  top(procPatch, nbrPatchInfo);
1361 
1362  label bFacei = procPatch.start()-mesh.nInternalFaces();
1363 
1364  forAll(nbrPatchInfo, i)
1365  {
1366  cop(faceValues[bFacei++], nbrPatchInfo[i]);
1367  }
1368  }
1369  }
1370  }
1371 
1372  // Do the cyclics.
1373  forAll(patches, patchi)
1374  {
1375  if (isA<cyclicPolyPatch>(patches[patchi]))
1376  {
1377  const cyclicPolyPatch& cycPatch =
1378  refCast<const cyclicPolyPatch>(patches[patchi]);
1379 
1380  if (cycPatch.owner())
1381  {
1382  // Owner does all.
1383  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
1384  label ownStart = cycPatch.start()-mesh.nInternalFaces();
1385  label nbrStart = nbrPatch.start()-mesh.nInternalFaces();
1386 
1387  label sz = cycPatch.size();
1388 
1389  // Transform (copy of) data on both sides
1390  Field<T> ownVals(SubField<T>(faceValues, sz, ownStart));
1391  top(nbrPatch, ownVals);
1392 
1393  Field<T> nbrVals(SubField<T>(faceValues, sz, nbrStart));
1394  top(cycPatch, nbrVals);
1395 
1396  label i0 = ownStart;
1397  forAll(nbrVals, i)
1398  {
1399  cop(faceValues[i0++], nbrVals[i]);
1400  }
1401 
1402  label i1 = nbrStart;
1403  forAll(ownVals, i)
1404  {
1405  cop(faceValues[i1++], ownVals[i]);
1406  }
1407  }
1408  }
1409  }
1410 }
1411 
1412 
1413 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1414 
1415 template<unsigned nBits, class CombineOp>
1418  const polyMesh& mesh,
1419  PackedList<nBits>& faceValues,
1420  const CombineOp& cop,
1421  const bool parRun
1422 )
1423 {
1424  if (faceValues.size() != mesh.nFaces())
1425  {
1427  << "Number of values " << faceValues.size()
1428  << " is not equal to the number of faces in the mesh "
1429  << mesh.nFaces() << abort(FatalError);
1430  }
1431 
1432  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1433 
1434  if (parRun)
1435  {
1437 
1438  // Send
1439 
1440  forAll(patches, patchi)
1441  {
1442  if
1443  (
1444  isA<processorPolyPatch>(patches[patchi])
1445  && patches[patchi].size() > 0
1446  )
1447  {
1448  const processorPolyPatch& procPatch =
1449  refCast<const processorPolyPatch>(patches[patchi]);
1450 
1451  List<unsigned int> patchInfo(procPatch.size());
1452  forAll(procPatch, i)
1453  {
1454  patchInfo[i] = faceValues[procPatch.start()+i];
1455  }
1456 
1457  UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
1458  toNbr << patchInfo;
1459  }
1460  }
1461 
1462 
1463  pBufs.finishedSends();
1464 
1465  // Receive and combine.
1466 
1467  forAll(patches, patchi)
1468  {
1469  if
1470  (
1471  isA<processorPolyPatch>(patches[patchi])
1472  && patches[patchi].size() > 0
1473  )
1474  {
1475  const processorPolyPatch& procPatch =
1476  refCast<const processorPolyPatch>(patches[patchi]);
1477 
1478  List<unsigned int> patchInfo(procPatch.size());
1479  {
1480  UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
1481  fromNbr >> patchInfo;
1482  }
1483 
1484  // Combine (bitwise)
1485  forAll(procPatch, i)
1486  {
1487  unsigned int patchVal = patchInfo[i];
1488  label meshFacei = procPatch.start()+i;
1489  unsigned int faceVal = faceValues[meshFacei];
1490  cop(faceVal, patchVal);
1491  faceValues[meshFacei] = faceVal;
1492  }
1493  }
1494  }
1495  }
1496 
1497  // Do the cyclics.
1498  forAll(patches, patchi)
1499  {
1500  if (isA<cyclicPolyPatch>(patches[patchi]))
1501  {
1502  const cyclicPolyPatch& cycPatch =
1503  refCast<const cyclicPolyPatch>(patches[patchi]);
1504 
1505  if (cycPatch.owner())
1506  {
1507  // Owner does all.
1508  const cyclicPolyPatch& nbrPatch = cycPatch.nbrPatch();
1509 
1510  for (label i = 0; i < cycPatch.size(); i++)
1511  {
1512  label meshFace0 = cycPatch.start()+i;
1513  unsigned int val0 = faceValues[meshFace0];
1514  label meshFace1 = nbrPatch.start()+i;
1515  unsigned int val1 = faceValues[meshFace1];
1516 
1517  unsigned int t = val0;
1518  cop(t, val1);
1519  faceValues[meshFace0] = t;
1520 
1521  cop(val1, val0);
1522  faceValues[meshFace1] = val1;
1523  }
1524  }
1525  }
1526  }
1527 }
1528 
1529 
1530 template<class T>
1533  const polyMesh& mesh,
1534  const UList<T>& cellData,
1535  List<T>& neighbourCellData
1536 )
1537 {
1538  if (cellData.size() != mesh.nCells())
1539  {
1541  << "Number of cell values " << cellData.size()
1542  << " is not equal to the number of cells in the mesh "
1543  << mesh.nCells() << abort(FatalError);
1544  }
1545 
1546  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1547 
1548  label nBnd = mesh.nFaces()-mesh.nInternalFaces();
1549 
1550  neighbourCellData.setSize(nBnd);
1551 
1552  forAll(patches, patchi)
1553  {
1554  const polyPatch& pp = patches[patchi];
1555  const labelUList& faceCells = pp.faceCells();
1556  forAll(faceCells, i)
1557  {
1558  label bFacei = pp.start()+i-mesh.nInternalFaces();
1559  neighbourCellData[bFacei] = cellData[faceCells[i]];
1560  }
1561  }
1562  syncTools::swapBoundaryFaceList(mesh, neighbourCellData);
1563 }
1564 
1565 
1566 template<unsigned nBits>
1569  const polyMesh& mesh,
1570  PackedList<nBits>& faceValues
1571 )
1572 {
1573  syncFaceList(mesh, faceValues, eqOp<unsigned int>());
1574 }
1575 
1576 
1577 template<unsigned nBits, class CombineOp>
1580  const polyMesh& mesh,
1581  PackedList<nBits>& pointValues,
1582  const CombineOp& cop,
1583  const unsigned int nullValue
1584 )
1585 {
1586  if (pointValues.size() != mesh.nPoints())
1587  {
1589  << "Number of values " << pointValues.size()
1590  << " is not equal to the number of points in the mesh "
1591  << mesh.nPoints() << abort(FatalError);
1592  }
1593 
1594  const globalMeshData& gd = mesh.globalData();
1595  const labelList& meshPoints = gd.coupledPatch().meshPoints();
1596 
1598  forAll(meshPoints, i)
1599  {
1600  cppFld[i] = pointValues[meshPoints[i]];
1601  }
1602 
1604  (
1605  cppFld,
1606  gd.globalPointSlaves(),
1608  gd.globalPointSlavesMap(),
1609  cop
1610  );
1611 
1612  // Extract back to mesh
1613  forAll(meshPoints, i)
1614  {
1615  pointValues[meshPoints[i]] = cppFld[i];
1616  }
1617 }
1618 
1619 
1620 template<unsigned nBits, class CombineOp>
1623  const polyMesh& mesh,
1624  PackedList<nBits>& edgeValues,
1625  const CombineOp& cop,
1626  const unsigned int nullValue
1627 )
1628 {
1629  if (edgeValues.size() != mesh.nEdges())
1630  {
1632  << "Number of values " << edgeValues.size()
1633  << " is not equal to the number of edges in the mesh "
1634  << mesh.nEdges() << abort(FatalError);
1635  }
1636 
1637  const globalMeshData& gd = mesh.globalData();
1638  const labelList& meshEdges = gd.coupledPatchMeshEdges();
1639 
1641  forAll(meshEdges, i)
1642  {
1643  cppFld[i] = edgeValues[meshEdges[i]];
1644  }
1645 
1647  (
1648  cppFld,
1649  gd.globalEdgeSlaves(),
1651  gd.globalEdgeSlavesMap(),
1652  cop
1653  );
1654 
1655  // Extract back to mesh
1656  forAll(meshEdges, i)
1657  {
1658  edgeValues[meshEdges[i]] = cppFld[i];
1659  }
1660 }
1661 
1662 
1663 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
static void syncBoundaryFaceList(const polyMesh &, UList< T > &, const CombineOp &cop, const TransformOp &top, const bool parRun=Pstream::parRun())
Synchronise values on boundary faces only.
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
const Map< label > & coupledPatchMeshEdgeMap() const
Return map from mesh edges to coupledPatch edges.
label nPoints() const
Return number of points supporting patch faces.
const labelListList & globalPointSlaves() const
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:59
static int masterNo()
Process index of the master.
Definition: UPstream.H:417
const labelList & nbrPoints() const
Return neighbour point labels. WIP.
const labelListList & globalPointTransformedSlaves() const
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
An STL-conforming const_iterator.
Definition: HashTable.H:481
const edgeList & coupledPoints() const
Return connected points (from patch local to neighbour patch local)
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
label nInternalFaces() const
static int firstSlave()
Process index of first slave.
Definition: UPstream.H:446
const mapDistribute & globalEdgeSlavesMap() const
HashTable< T, label, Hash< label > >::iterator iterator
Definition: Map.H:56
label nFaces() const
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
label nCells() const
Template function to specify if the data of a type are contiguous.
const mapDistribute & globalPointSlavesMap() const
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
patches[0]
const labelListList & globalEdgeSlaves() const
const cyclicPolyPatch & nbrPatch() const
const labelList & coupledPatchMeshEdges() const
Return map from coupledPatch edges to mesh edges.
Pre-declare related SubField type.
Definition: Field.H:60
const edgeList & coupledEdges() const
Return connected edges (from patch local to neighbour patch local).
const labelList & meshPoints() const
Return labelList of mesh points in patch. They are constructed.
Neighbour processor patch.
Input inter-processor communications stream.
Definition: IPstream.H:50
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:53
A list of faces which address into the list of points.
virtual bool owner() const
Does this side own the patch ?
label nGlobalPoints() const
Return number of globally shared points.
const Map< label > & meshPointMap() const
Mesh point map. Given the global point index find its.
const globalIndexAndTransform & globalTransforms() const
Global transforms numbering.
Definition: ops.H:70
A dynamically allocatable list of packed unsigned integers.
Definition: PackedList.H:117
const labelUList & faceCells() const
Return face-cell addressing.
Definition: polyPatch.C:334
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
static void syncEdgeList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh edges.
label nPoints
Cyclic plane patch.
3D tensor transformation operations.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
static void swapFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled face values.
Definition: syncTools.H:463
static void syncPointMap(const polyMesh &, Map< T > &pointValues, const CombineOp &cop, const TransformOp &top)
Synchronise values on selected points.
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1156
int neighbProcNo() const
Return neighbour processor number.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
Foam::polyBoundaryMesh.
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:54
label nEdges() const
Return number of edges in patch.
Map from edge (expressed as its endpoints) to value.
Definition: EdgeMap.H:47
labelList f(nPoints)
Output inter-processor communications stream.
Definition: OPstream.H:50
Buffers for inter-processor communications streams (UOPstream, UIPstream).
static void syncEdgeMap(const polyMesh &, EdgeMap< T > &edgeValues, const CombineOp &cop, const TransformOp &top)
Synchronise values on selected edges.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const volScalarField & T
label nEdges() const
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
static void syncData(List< Type > &pointData, const labelListList &slaves, const labelListList &transformedSlaves, const mapDistribute &slavesMap, const globalIndexAndTransform &, const CombineOp &cop, const TransformOp &top)
Helper: synchronise data with transforms.
label constructSize() const
Constructed data size.
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:381
void setSize(const label)
Reset size of List.
Definition: List.C:281
bool set(const label &, const T &newElmt)
Assign a new hashedEntry, overwriting existing entries.
Definition: HashTableI.H:91
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label patchi
Class containing processor-to-processor mapping information.
A List with indirect addressing.
Definition: fvMatrix.H:106
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
static void syncPointList(const polyMesh &, List< T > &, const CombineOp &cop, const T &nullValue, const TransformOp &top)
Synchronise values on all mesh points.
label nPoints() const
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:432
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
label size() const
Number of entries.
Definition: PackedListI.H:711
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
label size() const
Return the number of elements in the UList.
Definition: UListI.H:311
friend class iterator
Declare friendship with the iterator.
Definition: HashTable.H:194
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points.
Determination and storage of the possible independent transforms introduced by coupledPolyPatches, as well as all of the possible permutations of these transforms generated by the presence of multiple coupledPolyPatches, i.e. more than one cyclic boundary. Note that any given point can be on maximum 3 transforms only (and these transforms have to be perpendicular)
static int lastSlave(const label communicator=0)
Process index of last slave.
Definition: UPstream.H:452
A HashTable to objects of type <T> with a label key.
Definition: Map.H:49
const labelListList & globalEdgeTransformedSlaves() const
void syncPointData(List< Type > &pointData, const CombineOp &cop, const TransformOp &top) const
Helper to synchronise coupled patch point data.