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