1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
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.
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
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/>.
24 \*----------------------------------------------------------------------------*/
26 #include "syncTools.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"
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 // Combine val with existing value at index
39 template <class T, class CombineOp>
40 void Foam::syncTools::combine
48 typename Map<T>::iterator iter = pointValues.find(index);
50 if (iter != pointValues.end())
56 pointValues.insert(index, val);
61 // Combine val with existing value at (implicit index) e.
62 template <class T, class CombineOp>
63 void Foam::syncTools::combine
65 EdgeMap<T>& edgeValues,
71 typename EdgeMap<T>::iterator iter = edgeValues.find(index);
73 if (iter != edgeValues.end())
79 edgeValues.insert(index, val);
84 template <class T, class CombineOp, class TransformOp>
85 void Foam::syncTools::syncPointMap
88 Map<T>& pointValues, // from mesh point label to value
90 const TransformOp& top
93 const polyBoundaryMesh& patches = mesh.boundaryMesh();
95 // Synchronize multiple shared points.
96 const globalMeshData& pd = mesh.globalData();
98 // Values on shared points. Keyed on global shared index.
99 Map<T> sharedPointValues(0);
101 if (pd.nGlobalPoints() > 0)
103 // meshPoint per local index
104 const labelList& sharedPtLabels = pd.sharedPointLabels();
105 // global shared index per local index
106 const labelList& sharedPtAddr = pd.sharedPointAddr();
108 sharedPointValues.resize(sharedPtAddr.size());
110 // Fill my entries in the shared points
111 forAll(sharedPtLabels, i)
113 label meshPointI = sharedPtLabels[i];
115 typename Map<T>::const_iterator fnd =
116 pointValues.find(meshPointI);
118 if (fnd != pointValues.end())
124 sharedPtAddr[i], // index
132 if (Pstream::parRun())
134 PstreamBuffers pBufs(Pstream::nonBlocking);
138 forAll(patches, patchI)
142 isA<processorPolyPatch>(patches[patchI])
143 && patches[patchI].nPoints() > 0
146 const processorPolyPatch& procPatch =
147 refCast<const processorPolyPatch>(patches[patchI]);
149 // Get data per patchPoint in neighbouring point numbers.
151 const labelList& meshPts = procPatch.meshPoints();
152 const labelList& nbrPts = procPatch.neighbPoints();
154 // Extract local values. Create map from nbrPoint to value.
155 // Note: how small initial size?
156 Map<T> patchInfo(meshPts.size() / 20);
160 typename Map<T>::const_iterator iter =
161 pointValues.find(meshPts[i]);
163 if (iter != pointValues.end())
165 patchInfo.insert(nbrPts[i], iter());
169 UOPstream toNeighb(procPatch.neighbProcNo(), pBufs);
170 toNeighb << patchInfo;
174 pBufs.finishedSends();
176 // Receive and combine.
178 forAll(patches, patchI)
182 isA<processorPolyPatch>(patches[patchI])
183 && patches[patchI].nPoints() > 0
186 const processorPolyPatch& procPatch =
187 refCast<const processorPolyPatch>(patches[patchI]);
189 UIPstream fromNb(procPatch.neighbProcNo(), pBufs);
190 Map<T> nbrPatchInfo(fromNb);
193 top(procPatch, nbrPatchInfo);
195 const labelList& meshPts = procPatch.meshPoints();
197 // Only update those values which come from neighbour
198 forAllConstIter(typename Map<T>, nbrPatchInfo, nbrIter)
204 meshPts[nbrIter.key()],
213 forAll(patches, patchI)
215 if (isA<cyclicPolyPatch>(patches[patchI]))
217 const cyclicPolyPatch& cycPatch =
218 refCast<const cyclicPolyPatch>(patches[patchI]);
220 if (cycPatch.owner())
224 const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
225 const edgeList& coupledPoints = cycPatch.coupledPoints();
226 const labelList& meshPtsA = cycPatch.meshPoints();
227 const labelList& meshPtsB = nbrPatch.meshPoints();
229 // Extract local values. Create map from coupled-edge to value.
230 Map<T> half0Values(meshPtsA.size() / 20);
231 Map<T> half1Values(half0Values.size());
233 forAll(coupledPoints, i)
235 const edge& e = coupledPoints[i];
237 typename Map<T>::const_iterator point0Fnd =
238 pointValues.find(meshPtsA[e[0]]);
240 if (point0Fnd != pointValues.end())
242 half0Values.insert(i, point0Fnd());
245 typename Map<T>::const_iterator point1Fnd =
246 pointValues.find(meshPtsB[e[1]]);
248 if (point1Fnd != pointValues.end())
250 half1Values.insert(i, point1Fnd());
254 // Transform to receiving side
255 top(cycPatch, half1Values);
256 top(nbrPatch, half0Values);
258 forAll(coupledPoints, i)
260 const edge& e = coupledPoints[i];
262 typename Map<T>::const_iterator half0Fnd =
265 if (half0Fnd != half0Values.end())
276 typename Map<T>::const_iterator half1Fnd =
279 if (half1Fnd != half1Values.end())
294 // Synchronize multiple shared points.
295 if (pd.nGlobalPoints() > 0)
297 // meshPoint per local index
298 const labelList& sharedPtLabels = pd.sharedPointLabels();
299 // global shared index per local index
300 const labelList& sharedPtAddr = pd.sharedPointAddr();
304 if (Pstream::parRun())
306 if (Pstream::master())
308 // Receive the edges using shared points from the slave.
311 int slave=Pstream::firstSlave();
312 slave<=Pstream::lastSlave();
316 IPstream fromSlave(Pstream::scheduled, slave);
317 Map<T> nbrValues(fromSlave);
319 // Merge neighbouring values with my values
320 forAllConstIter(typename Map<T>, nbrValues, iter)
335 int slave=Pstream::firstSlave();
336 slave<=Pstream::lastSlave();
340 OPstream toSlave(Pstream::scheduled, slave);
341 toSlave << sharedPointValues;
346 // Slave: send to master
348 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
349 toMaster << sharedPointValues;
351 // Receive merged values
358 fromMaster >> sharedPointValues;
364 // Merge sharedPointValues (keyed on sharedPointAddr) into
365 // pointValues (keyed on mesh points).
367 // Map from global shared index to meshpoint
368 Map<label> sharedToMeshPoint(2*sharedPtAddr.size());
369 forAll(sharedPtAddr, i)
371 sharedToMeshPoint.insert(sharedPtAddr[i], sharedPtLabels[i]);
374 forAllConstIter(Map<label>, sharedToMeshPoint, iter)
376 // Do I have a value for my shared point
377 typename Map<T>::const_iterator sharedFnd =
378 sharedPointValues.find(iter.key());
380 if (sharedFnd != sharedPointValues.end())
382 pointValues.set(iter(), sharedFnd());
389 template <class T, class CombineOp, class TransformOp>
390 void Foam::syncTools::syncEdgeMap
392 const polyMesh& mesh,
393 EdgeMap<T>& edgeValues,
394 const CombineOp& cop,
395 const TransformOp& top
398 const polyBoundaryMesh& patches = mesh.boundaryMesh();
401 // Do synchronisation without constructing globalEdge addressing
402 // (since this constructs mesh edge addressing)
405 // Swap proc patch info
406 // ~~~~~~~~~~~~~~~~~~~~
408 if (Pstream::parRun())
410 PstreamBuffers pBufs(Pstream::nonBlocking);
414 forAll(patches, patchI)
418 isA<processorPolyPatch>(patches[patchI])
419 && patches[patchI].nEdges() > 0
422 const processorPolyPatch& procPatch =
423 refCast<const processorPolyPatch>(patches[patchI]);
426 // Get data per patch edge in neighbouring edge.
428 const edgeList& edges = procPatch.edges();
429 const labelList& meshPts = procPatch.meshPoints();
430 const labelList& nbrPts = procPatch.neighbPoints();
432 EdgeMap<T> patchInfo(edges.size() / 20);
436 const edge& e = edges[i];
437 const edge meshEdge(meshPts[e[0]], meshPts[e[1]]);
439 typename EdgeMap<T>::const_iterator iter =
440 edgeValues.find(meshEdge);
442 if (iter != edgeValues.end())
444 const edge nbrEdge(nbrPts[e[0]], nbrPts[e[1]]);
445 patchInfo.insert(nbrEdge, iter());
449 UOPstream toNeighb(procPatch.neighbProcNo(), pBufs);
450 toNeighb << patchInfo;
454 pBufs.finishedSends();
456 // Receive and combine.
458 forAll(patches, patchI)
462 isA<processorPolyPatch>(patches[patchI])
463 && patches[patchI].nEdges() > 0
466 const processorPolyPatch& procPatch =
467 refCast<const processorPolyPatch>(patches[patchI]);
469 EdgeMap<T> nbrPatchInfo;
471 UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
472 fromNbr >> nbrPatchInfo;
475 // Apply transform to convert to this side properties.
476 top(procPatch, nbrPatchInfo);
479 // Only update those values which come from neighbour
480 const labelList& meshPts = procPatch.meshPoints();
482 forAllConstIter(typename EdgeMap<T>, nbrPatchInfo, nbrIter)
484 const edge& e = nbrIter.key();
485 const edge meshEdge(meshPts[e[0]], meshPts[e[1]]);
503 forAll(patches, patchI)
505 if (isA<cyclicPolyPatch>(patches[patchI]))
507 const cyclicPolyPatch& cycPatch =
508 refCast<const cyclicPolyPatch>(patches[patchI]);
510 if (cycPatch.owner())
514 const edgeList& coupledEdges = cycPatch.coupledEdges();
515 const labelList& meshPtsA = cycPatch.meshPoints();
516 const edgeList& edgesA = cycPatch.edges();
517 const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
518 const labelList& meshPtsB = nbrPatch.meshPoints();
519 const edgeList& edgesB = nbrPatch.edges();
521 // Extract local values. Create map from edge to value.
522 Map<T> half0Values(edgesA.size() / 20);
523 Map<T> half1Values(half0Values.size());
525 forAll(coupledEdges, i)
527 const edge& twoEdges = coupledEdges[i];
530 const edge& e0 = edgesA[twoEdges[0]];
531 const edge meshEdge0(meshPtsA[e0[0]], meshPtsA[e0[1]]);
533 typename EdgeMap<T>::const_iterator iter =
534 edgeValues.find(meshEdge0);
536 if (iter != edgeValues.end())
538 half0Values.insert(i, iter());
542 const edge& e1 = edgesB[twoEdges[1]];
543 const edge meshEdge1(meshPtsB[e1[0]], meshPtsB[e1[1]]);
545 typename EdgeMap<T>::const_iterator iter =
546 edgeValues.find(meshEdge1);
548 if (iter != edgeValues.end())
550 half1Values.insert(i, iter());
555 // Transform to this side
556 top(cycPatch, half1Values);
557 top(nbrPatch, half0Values);
560 // Extract and combine information
562 forAll(coupledEdges, i)
564 const edge& twoEdges = coupledEdges[i];
566 typename Map<T>::const_iterator half1Fnd =
569 if (half1Fnd != half1Values.end())
571 const edge& e0 = edgesA[twoEdges[0]];
572 const edge meshEdge0(meshPtsA[e0[0]], meshPtsA[e0[1]]);
583 typename Map<T>::const_iterator half0Fnd =
585 if (half0Fnd != half0Values.end())
587 const edge& e1 = edgesB[twoEdges[1]];
588 const edge meshEdge1(meshPtsB[e1[0]], meshPtsB[e1[1]]);
603 // Synchronize multiple shared points.
604 // Problem is that we don't want to construct shared edges so basically
605 // we do here like globalMeshData but then using sparse edge representation
606 // (EdgeMap instead of mesh.edges())
608 const globalMeshData& pd = mesh.globalData();
609 const labelList& sharedPtAddr = pd.sharedPointAddr();
610 const labelList& sharedPtLabels = pd.sharedPointLabels();
612 // 1. Create map from meshPoint to globalShared index.
613 Map<label> meshToShared(2*sharedPtLabels.size());
614 forAll(sharedPtLabels, i)
616 meshToShared.insert(sharedPtLabels[i], sharedPtAddr[i]);
619 // Values on shared points. From two sharedPtAddr indices to a value.
620 EdgeMap<T> sharedEdgeValues(meshToShared.size());
622 // From shared edge to mesh edge. Used for merging later on.
623 EdgeMap<edge> potentialSharedEdge(meshToShared.size());
625 // 2. Find any edges using two global shared points. These will always be
626 // on the outside of the mesh. (though might not be on coupled patch
627 // if is single edge and not on coupled face)
628 // Store value (if any) on sharedEdgeValues
629 for (label faceI = mesh.nInternalFaces(); faceI < mesh.nFaces(); faceI++)
631 const face& f = mesh.faces()[faceI];
636 label v1 = f[f.fcIndex(fp)];
638 Map<label>::const_iterator v0Fnd = meshToShared.find(v0);
640 if (v0Fnd != meshToShared.end())
642 Map<label>::const_iterator v1Fnd = meshToShared.find(v1);
644 if (v1Fnd != meshToShared.end())
646 const edge meshEdge(v0, v1);
648 // edge in shared point labels
649 const edge sharedEdge(v0Fnd(), v1Fnd());
651 // Store mesh edge as a potential shared edge.
652 potentialSharedEdge.insert(sharedEdge, meshEdge);
654 typename EdgeMap<T>::const_iterator edgeFnd =
655 edgeValues.find(meshEdge);
657 if (edgeFnd != edgeValues.end())
659 // edge exists in edgeValues. See if already in map
660 // (since on same processor, e.g. cyclic)
675 // Now sharedEdgeValues will contain per potential sharedEdge the value.
676 // (potential since an edge having two shared points is not nessecary a
678 // Reduce this on the master.
680 if (Pstream::parRun())
682 if (Pstream::master())
684 // Receive the edges using shared points from the slave.
687 int slave=Pstream::firstSlave();
688 slave<=Pstream::lastSlave();
692 IPstream fromSlave(Pstream::scheduled, slave);
693 EdgeMap<T> nbrValues(fromSlave);
695 // Merge neighbouring values with my values
696 forAllConstIter(typename EdgeMap<T>, nbrValues, iter)
711 int slave=Pstream::firstSlave();
712 slave<=Pstream::lastSlave();
717 OPstream toSlave(Pstream::scheduled, slave);
718 toSlave << sharedEdgeValues;
725 OPstream toMaster(Pstream::scheduled, Pstream::masterNo());
726 toMaster << sharedEdgeValues;
728 // Receive merged values
730 IPstream fromMaster(Pstream::scheduled, Pstream::masterNo());
731 fromMaster >> sharedEdgeValues;
737 // Merge sharedEdgeValues (keyed on sharedPointAddr) into edgeValues
738 // (keyed on mesh points).
740 // Loop over all my shared edges.
741 forAllConstIter(typename EdgeMap<edge>, potentialSharedEdge, iter)
743 const edge& sharedEdge = iter.key();
744 const edge& meshEdge = iter();
746 // Do I have a value for the shared edge?
747 typename EdgeMap<T>::const_iterator sharedFnd =
748 sharedEdgeValues.find(sharedEdge);
750 if (sharedFnd != sharedEdgeValues.end())
764 //template <class T, class CombineOp, class TransformOp>
765 //void Foam::syncTools::syncPointList
767 // const polyMesh& mesh,
768 // List<T>& pointValues,
769 // const CombineOp& cop,
770 // const T& nullValue,
771 // const TransformOp& top
774 // if (pointValues.size() != mesh.nPoints())
778 // "syncTools<class T, class CombineOp>::syncPointList"
779 // "(const polyMesh&, List<T>&, const CombineOp&, const T&"
781 // ) << "Number of values " << pointValues.size()
782 // << " is not equal to the number of points in the mesh "
783 // << mesh.nPoints() << abort(FatalError);
786 // const polyBoundaryMesh& patches = mesh.boundaryMesh();
788 // // Synchronize multiple shared points.
789 // const globalMeshData& pd = mesh.globalData();
791 // // Values on shared points.
792 // Field<T> sharedPts(0);
793 // if (pd.nGlobalPoints() > 0)
795 // // Values on shared points.
796 // sharedPts.setSize(pd.nGlobalPoints(), nullValue);
798 // forAll(pd.sharedPointLabels(), i)
800 // label meshPointI = pd.sharedPointLabels()[i];
801 // // Fill my entries in the shared points
802 // sharedPts[pd.sharedPointAddr()[i]] = pointValues[meshPointI];
806 // if (Pstream::parRun())
808 // PstreamBuffers pBufs(Pstream::nonBlocking);
812 // forAll(patches, patchI)
816 // isA<processorPolyPatch>(patches[patchI])
817 // && patches[patchI].nPoints() > 0
820 // const processorPolyPatch& procPatch =
821 // refCast<const processorPolyPatch>(patches[patchI]);
823 // // Get data per patchPoint in neighbouring point numbers.
824 // Field<T> patchInfo(procPatch.nPoints());
826 // const labelList& meshPts = procPatch.meshPoints();
827 // const labelList& nbrPts = procPatch.neighbPoints();
829 // forAll(nbrPts, pointI)
831 // label nbrPointI = nbrPts[pointI];
832 // patchInfo[nbrPointI] = pointValues[meshPts[pointI]];
835 // UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
836 // toNbr << patchInfo;
840 // pBufs.finishedSends();
842 // // Receive and combine.
844 // forAll(patches, patchI)
848 // isA<processorPolyPatch>(patches[patchI])
849 // && patches[patchI].nPoints() > 0
852 // const processorPolyPatch& procPatch =
853 // refCast<const processorPolyPatch>(patches[patchI]);
855 // Field<T> nbrPatchInfo(procPatch.nPoints());
857 // UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
858 // fromNbr >> nbrPatchInfo;
861 // // Transform to this side
862 // top(procPatch, nbrPatchInfo);
864 // const labelList& meshPts = procPatch.meshPoints();
866 // forAll(meshPts, pointI)
868 // label meshPointI = meshPts[pointI];
869 // cop(pointValues[meshPointI], nbrPatchInfo[pointI]);
875 // // Do the cyclics.
876 // forAll(patches, patchI)
878 // if (isA<cyclicPolyPatch>(patches[patchI]))
880 // const cyclicPolyPatch& cycPatch =
881 // refCast<const cyclicPolyPatch>(patches[patchI]);
883 // if (cycPatch.owner())
885 // // Owner does all.
887 // const edgeList& coupledPoints = cycPatch.coupledPoints();
888 // const labelList& meshPts = cycPatch.meshPoints();
889 // const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
890 // const labelList& nbrMeshPoints = nbrPatch.meshPoints();
892 // Field<T> half0Values(coupledPoints.size());
893 // Field<T> half1Values(coupledPoints.size());
895 // forAll(coupledPoints, i)
897 // const edge& e = coupledPoints[i];
898 // half0Values[i] = pointValues[meshPts[e[0]]];
899 // half1Values[i] = pointValues[nbrMeshPoints[e[1]]];
902 // //SubField<T> slice0(half0Values, half0Values.size());
903 // //SubField<T> slice1(half1Values, half1Values.size());
904 // //top(cycPatch, reinterpret_cast<Field<T>&>(slice1));
905 // //top(nbrPatch, reinterpret_cast<Field<T>&>(slice0));
907 // top(cycPatch, half1Values);
908 // top(nbrPatch, half0Values);
910 // forAll(coupledPoints, i)
912 // const edge& e = coupledPoints[i];
913 // cop(pointValues[meshPts[e[0]]], half1Values[i]);
914 // cop(pointValues[nbrMeshPoints[e[1]]], half0Values[i]);
920 // // Synchronize multiple shared points.
921 // const globalMeshData& pd = mesh.globalData();
923 // if (pd.nGlobalPoints() > 0)
925 // // Combine on master.
926 // Pstream::listCombineGather(sharedPts, cop);
927 // Pstream::listCombineScatter(sharedPts);
929 // // Now we will all have the same information. Merge it back with
930 // // my local information.
931 // forAll(pd.sharedPointLabels(), i)
933 // label meshPointI = pd.sharedPointLabels()[i];
934 // pointValues[meshPointI] = sharedPts[pd.sharedPointAddr()[i]];
940 //template <class T, class CombineOp, class TransformOp>
941 //void Foam::syncTools::syncPointList
943 // const polyMesh& mesh,
944 // const labelList& meshPoints,
945 // List<T>& pointValues,
946 // const CombineOp& cop,
947 // const T& nullValue,
948 // const TransformOp& top
951 // if (pointValues.size() != meshPoints.size())
955 // "syncTools<class T, class CombineOp>::syncPointList"
956 // "(const polyMesh&, const labelList&, List<T>&, const CombineOp&"
957 // ", const T&, const bool)"
958 // ) << "Number of values " << pointValues.size()
959 // << " is not equal to the number of points "
960 // << meshPoints.size() << abort(FatalError);
963 // if (!hasCouples(mesh.boundaryMesh()))
968 // Field<T> meshValues(mesh.nPoints(), nullValue);
970 // forAll(meshPoints, i)
972 // meshValues[meshPoints[i]] = pointValues[i];
975 // syncTools::syncPointList
979 // cop, // combine op
980 // nullValue, // null value
981 // top // position or field
984 // forAll(meshPoints, i)
986 // pointValues[i] = meshValues[meshPoints[i]];
990 template <class T, class CombineOp, class TransformOp>
991 void Foam::syncTools::syncPointList
993 const polyMesh& mesh,
994 List<T>& pointValues,
995 const CombineOp& cop,
997 const TransformOp& top
1000 if (pointValues.size() != mesh.nPoints())
1004 "syncTools<class T, class CombineOp>::syncPointList"
1005 "(const polyMesh&, List<T>&, const CombineOp&, const T&"
1007 ) << "Number of values " << pointValues.size()
1008 << " is not equal to the number of points in the mesh "
1009 << mesh.nPoints() << abort(FatalError);
1012 mesh.globalData().syncPointData(pointValues, cop, top);
1016 //template <class CombineOp>
1017 //void Foam::syncTools::syncPointPositions
1019 // const polyMesh& mesh,
1020 // List<point>& pointValues,
1021 // const CombineOp& cop,
1022 // const point& nullValue
1025 // if (pointValues.size() != mesh.nPoints())
1029 // "syncTools<class CombineOp>::syncPointPositions"
1030 // "(const polyMesh&, List<point>&, const CombineOp&, const point&)"
1031 // ) << "Number of values " << pointValues.size()
1032 // << " is not equal to the number of points in the mesh "
1033 // << mesh.nPoints() << abort(FatalError);
1036 // mesh.globalData().syncPointData(pointValues, cop, true);
1040 template <class T, class CombineOp, class TransformOp>
1041 void Foam::syncTools::syncPointList
1043 const polyMesh& mesh,
1044 const labelList& meshPoints,
1045 List<T>& pointValues,
1046 const CombineOp& cop,
1048 const TransformOp& top
1051 if (pointValues.size() != meshPoints.size())
1055 "syncTools<class T, class CombineOp>::syncPointList"
1056 "(const polyMesh&, List<T>&, const CombineOp&, const T&)"
1057 ) << "Number of values " << pointValues.size()
1058 << " is not equal to the number of meshPoints "
1059 << meshPoints.size() << abort(FatalError);
1061 const globalMeshData& gd = mesh.globalData();
1062 const indirectPrimitivePatch& cpp = gd.coupledPatch();
1063 const Map<label>& mpm = cpp.meshPointMap();
1065 List<T> cppFld(cpp.nPoints(), nullValue);
1067 forAll(meshPoints, i)
1069 label pointI = meshPoints[i];
1070 Map<label>::const_iterator iter = mpm.find(pointI);
1071 if (iter != mpm.end())
1073 cppFld[iter()] = pointValues[i];
1077 globalMeshData::syncData
1080 gd.globalPointSlaves(),
1081 gd.globalPointTransformedSlaves(),
1082 gd.globalPointSlavesMap(),
1083 gd.globalTransforms(),
1088 forAll(meshPoints, i)
1090 label pointI = meshPoints[i];
1091 Map<label>::const_iterator iter = mpm.find(pointI);
1092 if (iter != mpm.end())
1094 pointValues[i] = cppFld[iter()];
1100 //template <class CombineOp>
1101 //void Foam::syncTools::syncPointPositions
1103 // const polyMesh& mesh,
1104 // const labelList& meshPoints,
1105 // List<point>& pointValues,
1106 // const CombineOp& cop,
1107 // const point& nullValue
1110 // if (pointValues.size() != meshPoints.size())
1114 // "syncTools<class CombineOp>::syncPointList"
1115 // "(const polyMesh&, List<point>&, const CombineOp&, const point&)"
1116 // ) << "Number of values " << pointValues.size()
1117 // << " is not equal to the number of meshPoints "
1118 // << meshPoints.size() << abort(FatalError);
1120 // const globalMeshData& gd = mesh.globalData();
1121 // const indirectPrimitivePatch& cpp = gd.coupledPatch();
1122 // const Map<label>& mpm = cpp.meshPointMap();
1124 // List<point> cppFld(cpp.nPoints(), nullValue);
1126 // forAll(meshPoints, i)
1128 // label pointI = meshPoints[i];
1129 // Map<label>::const_iterator iter = mpm.find(pointI);
1130 // if (iter != mpm.end())
1132 // cppFld[iter()] = pointValues[i];
1136 // globalMeshData::syncData
1139 // gd.globalPointSlaves(),
1140 // gd.globalPointTransformedSlaves(),
1141 // gd.globalPointSlavesMap(),
1142 // gd.globalTransforms(),
1144 // true, //position?
1145 // mapDistribute::transform() // not used
1148 // forAll(meshPoints, i)
1150 // label pointI = meshPoints[i];
1151 // Map<label>::const_iterator iter = mpm.find(pointI);
1152 // if (iter != mpm.end())
1154 // pointValues[i] = cppFld[iter()];
1160 template <class T, class CombineOp, class TransformOp>
1161 void Foam::syncTools::syncEdgeList
1163 const polyMesh& mesh,
1164 List<T>& edgeValues,
1165 const CombineOp& cop,
1167 const TransformOp& top
1170 if (edgeValues.size() != mesh.nEdges())
1174 "syncTools<class T, class CombineOp>::syncEdgeList"
1175 "(const polyMesh&, List<T>&, const CombineOp&, const T&)"
1176 ) << "Number of values " << edgeValues.size()
1177 << " is not equal to the number of edges in the mesh "
1178 << mesh.nEdges() << abort(FatalError);
1181 const globalMeshData& gd = mesh.globalData();
1182 const labelList& meshEdges = gd.coupledPatchMeshEdges();
1183 const globalIndexAndTransform& git = gd.globalTransforms();
1184 const mapDistribute& edgeMap = gd.globalEdgeSlavesMap();
1186 List<T> cppFld(UIndirectList<T>(edgeValues, meshEdges));
1188 globalMeshData::syncData
1191 gd.globalEdgeSlaves(),
1192 gd.globalEdgeTransformedSlaves(),
1199 // Extract back onto mesh
1200 forAll(meshEdges, i)
1202 edgeValues[meshEdges[i]] = cppFld[i];
1207 //template <class CombineOp>
1208 //void Foam::syncTools::syncEdgePositions
1210 // const polyMesh& mesh,
1211 // List<point>& edgeValues,
1212 // const CombineOp& cop,
1213 // const point& nullValue
1216 // if (edgeValues.size() != mesh.nEdges())
1220 // "syncTools<class CombineOp>::syncEdgePositions"
1221 // "(const polyMesh&, List<point>&, const CombineOp&, const point&)"
1222 // ) << "Number of values " << edgeValues.size()
1223 // << " is not equal to the number of edges in the mesh "
1224 // << mesh.nEdges() << abort(FatalError);
1227 // const globalMeshData& gd = mesh.globalData();
1228 // const labelList& meshEdges = gd.coupledPatchMeshEdges();
1229 // const globalIndexAndTransform& git = gd.globalTransforms();
1230 // const mapDistribute& map = gd.globalEdgeSlavesMap();
1232 // List<point> cppFld(UIndirectList<point>(edgeValues, meshEdges));
1234 // globalMeshData::syncData
1237 // gd.globalEdgeSlaves(),
1238 // gd.globalEdgeTransformedSlaves(),
1242 // true, //position?
1243 // mapDistribute::transform() // not used
1246 // // Extract back onto mesh
1247 // forAll(meshEdges, i)
1249 // edgeValues[meshEdges[i]] = cppFld[i];
1254 template <class T, class CombineOp, class TransformOp>
1255 void Foam::syncTools::syncBoundaryFaceList
1257 const polyMesh& mesh,
1258 UList<T>& faceValues,
1259 const CombineOp& cop,
1260 const TransformOp& top
1263 const label nBFaces = mesh.nFaces() - mesh.nInternalFaces();
1265 if (faceValues.size() != nBFaces)
1269 "syncTools<class T, class CombineOp>::syncBoundaryFaceList"
1270 "(const polyMesh&, UList<T>&, const CombineOp&"
1272 ) << "Number of values " << faceValues.size()
1273 << " is not equal to the number of boundary faces in the mesh "
1274 << nBFaces << abort(FatalError);
1277 const polyBoundaryMesh& patches = mesh.boundaryMesh();
1279 if (Pstream::parRun())
1281 PstreamBuffers pBufs(Pstream::nonBlocking);
1285 forAll(patches, patchI)
1289 isA<processorPolyPatch>(patches[patchI])
1290 && patches[patchI].size() > 0
1293 const processorPolyPatch& procPatch =
1294 refCast<const processorPolyPatch>(patches[patchI]);
1296 label patchStart = procPatch.start()-mesh.nInternalFaces();
1298 UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
1299 toNbr << SubField<T>(faceValues, procPatch.size(), patchStart);
1304 pBufs.finishedSends();
1307 // Receive and combine.
1309 forAll(patches, patchI)
1313 isA<processorPolyPatch>(patches[patchI])
1314 && patches[patchI].size() > 0
1317 const processorPolyPatch& procPatch =
1318 refCast<const processorPolyPatch>(patches[patchI]);
1320 Field<T> nbrPatchInfo(procPatch.size());
1322 UIPstream fromNeighb(procPatch.neighbProcNo(), pBufs);
1323 fromNeighb >> nbrPatchInfo;
1325 top(procPatch, nbrPatchInfo);
1327 label bFaceI = procPatch.start()-mesh.nInternalFaces();
1329 forAll(nbrPatchInfo, i)
1331 cop(faceValues[bFaceI++], nbrPatchInfo[i]);
1338 forAll(patches, patchI)
1340 if (isA<cyclicPolyPatch>(patches[patchI]))
1342 const cyclicPolyPatch& cycPatch =
1343 refCast<const cyclicPolyPatch>(patches[patchI]);
1345 if (cycPatch.owner())
1348 const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
1349 label ownStart = cycPatch.start()-mesh.nInternalFaces();
1350 label nbrStart = nbrPatch.start()-mesh.nInternalFaces();
1352 label sz = cycPatch.size();
1354 // Transform (copy of) data on both sides
1355 Field<T> ownVals(SubField<T>(faceValues, sz, ownStart));
1356 top(nbrPatch, ownVals);
1358 Field<T> nbrVals(SubField<T>(faceValues, sz, nbrStart));
1359 top(cycPatch, nbrVals);
1361 label i0 = ownStart;
1364 cop(faceValues[i0++], nbrVals[i]);
1367 label i1 = nbrStart;
1370 cop(faceValues[i1++], ownVals[i]);
1378 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1380 template <unsigned nBits, class CombineOp>
1381 void Foam::syncTools::syncFaceList
1383 const polyMesh& mesh,
1384 PackedList<nBits>& faceValues,
1385 const CombineOp& cop
1388 if (faceValues.size() != mesh.nFaces())
1392 "syncTools<unsigned nBits, class CombineOp>::syncFaceList"
1393 "(const polyMesh&, PackedList<nBits>&, const CombineOp&)"
1394 ) << "Number of values " << faceValues.size()
1395 << " is not equal to the number of faces in the mesh "
1396 << mesh.nFaces() << abort(FatalError);
1399 const polyBoundaryMesh& patches = mesh.boundaryMesh();
1401 if (Pstream::parRun())
1403 PstreamBuffers pBufs(Pstream::nonBlocking);
1407 forAll(patches, patchI)
1411 isA<processorPolyPatch>(patches[patchI])
1412 && patches[patchI].size() > 0
1415 const processorPolyPatch& procPatch =
1416 refCast<const processorPolyPatch>(patches[patchI]);
1418 List<unsigned int> patchInfo(procPatch.size());
1419 forAll(procPatch, i)
1421 patchInfo[i] = faceValues[procPatch.start()+i];
1424 UOPstream toNbr(procPatch.neighbProcNo(), pBufs);
1430 pBufs.finishedSends();
1432 // Receive and combine.
1434 forAll(patches, patchI)
1438 isA<processorPolyPatch>(patches[patchI])
1439 && patches[patchI].size() > 0
1442 const processorPolyPatch& procPatch =
1443 refCast<const processorPolyPatch>(patches[patchI]);
1445 List<unsigned int> patchInfo(procPatch.size());
1447 UIPstream fromNbr(procPatch.neighbProcNo(), pBufs);
1448 fromNbr >> patchInfo;
1451 // Combine (bitwise)
1452 forAll(procPatch, i)
1454 unsigned int patchVal = patchInfo[i];
1455 label meshFaceI = procPatch.start()+i;
1456 unsigned int faceVal = faceValues[meshFaceI];
1457 cop(faceVal, patchVal);
1458 faceValues[meshFaceI] = faceVal;
1465 forAll(patches, patchI)
1467 if (isA<cyclicPolyPatch>(patches[patchI]))
1469 const cyclicPolyPatch& cycPatch =
1470 refCast<const cyclicPolyPatch>(patches[patchI]);
1472 if (cycPatch.owner())
1475 const cyclicPolyPatch& nbrPatch = cycPatch.neighbPatch();
1477 for (label i = 0; i < cycPatch.size(); i++)
1479 label meshFace0 = cycPatch.start()+i;
1480 unsigned int val0 = faceValues[meshFace0];
1481 label meshFace1 = nbrPatch.start()+i;
1482 unsigned int val1 = faceValues[meshFace1];
1484 unsigned int t = val0;
1486 faceValues[meshFace0] = t;
1489 faceValues[meshFace1] = val1;
1497 template <unsigned nBits>
1498 void Foam::syncTools::swapFaceList
1500 const polyMesh& mesh,
1501 PackedList<nBits>& faceValues
1504 syncFaceList(mesh, faceValues, eqOp<unsigned int>());
1508 template <unsigned nBits, class CombineOp>
1509 void Foam::syncTools::syncPointList
1511 const polyMesh& mesh,
1512 PackedList<nBits>& pointValues,
1513 const CombineOp& cop,
1514 const unsigned int nullValue
1517 if (pointValues.size() != mesh.nPoints())
1521 "syncTools<unsigned nBits, class CombineOp>::syncPointList"
1522 "(const polyMesh&, PackedList<nBits>&, const CombineOp&"
1523 ", const unsigned int)"
1524 ) << "Number of values " << pointValues.size()
1525 << " is not equal to the number of points in the mesh "
1526 << mesh.nPoints() << abort(FatalError);
1529 const globalMeshData& gd = mesh.globalData();
1530 const labelList& meshPoints = gd.coupledPatch().meshPoints();
1532 List<unsigned int> cppFld(gd.globalPointSlavesMap().constructSize());
1533 forAll(meshPoints, i)
1535 cppFld[i] = pointValues[meshPoints[i]];
1538 globalMeshData::syncData
1541 gd.globalPointSlaves(),
1542 gd.globalPointTransformedSlaves(),
1543 gd.globalPointSlavesMap(),
1547 // Extract back to mesh
1548 forAll(meshPoints, i)
1550 pointValues[meshPoints[i]] = cppFld[i];
1555 template <unsigned nBits, class CombineOp>
1556 void Foam::syncTools::syncEdgeList
1558 const polyMesh& mesh,
1559 PackedList<nBits>& edgeValues,
1560 const CombineOp& cop,
1561 const unsigned int nullValue
1564 if (edgeValues.size() != mesh.nEdges())
1568 "syncTools<unsigned nBits, class CombineOp>::syncEdgeList"
1569 "(const polyMesh&, PackedList<nBits>&, const CombineOp&"
1570 ", const unsigned int)"
1571 ) << "Number of values " << edgeValues.size()
1572 << " is not equal to the number of edges in the mesh "
1573 << mesh.nEdges() << abort(FatalError);
1576 const globalMeshData& gd = mesh.globalData();
1577 const labelList& meshEdges = gd.coupledPatchMeshEdges();
1579 List<unsigned int> cppFld(gd.globalEdgeSlavesMap().constructSize());
1580 forAll(meshEdges, i)
1582 cppFld[i] = edgeValues[meshEdges[i]];
1585 globalMeshData::syncData
1588 gd.globalEdgeSlaves(),
1589 gd.globalEdgeTransformedSlaves(),
1590 gd.globalEdgeSlavesMap(),
1594 // Extract back to mesh
1595 forAll(meshEdges, i)
1597 edgeValues[meshEdges[i]] = cppFld[i];
1602 // ************************************************************************* //