globalPoints.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "globalPoints.H"
30 #include "processorPolyPatch.H"
31 #include "cyclicPolyPatch.H"
32 #include "polyMesh.H"
33 #include "mapDistribute.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(globalPoints, 0);
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 Foam::label Foam::globalPoints::countPatchPoints
46 (
47  const polyBoundaryMesh& patches
48 )
49 {
50  label nTotPoints = 0;
51 
52  for (const polyPatch& pp : patches)
53  {
54  if (pp.coupled())
55  {
56  nTotPoints += pp.nPoints();
57  }
58  }
59  return nTotPoints;
60 }
61 
62 
63 Foam::label Foam::globalPoints::findSamePoint
64 (
65  const labelPairList& allInfo,
66  const labelPair& info
67 ) const
68 {
69  const label proci = globalTransforms_.processor(info);
70  const label index = globalTransforms_.index(info);
71 
72  forAll(allInfo, i)
73  {
74  if
75  (
76  globalTransforms_.processor(allInfo[i]) == proci
77  && globalTransforms_.index(allInfo[i]) == index
78  )
79  {
80  return i;
81  }
82  }
83  return -1;
84 }
85 
86 
87 Foam::labelPairList Foam::globalPoints::addSendTransform
88 (
89  const label patchi,
90  const labelPairList& info
91 ) const
92 {
93  scalar tol = refCast<const coupledPolyPatch>
94  (
95  mesh_.boundaryMesh()[patchi]
96  ).matchTolerance();
97 
98  labelPairList sendInfo(info.size());
99 
100  forAll(info, i)
101  {
102  //Pout<< " adding send transform to" << nl
103  // << " proc:" << globalTransforms_.processor(info[i])
104  // << nl
105  // << " index:" << globalTransforms_.index(info[i]) << nl
106  // << " trafo:"
107  // << globalTransforms_.decodeTransformIndex
108  // (globalTransforms_.transformIndex(info[i]))
109  // << endl;
110 
111  sendInfo[i] = globalTransforms_.encode
112  (
113  globalTransforms_.processor(info[i]),
114  globalTransforms_.index(info[i]),
115  globalTransforms_.addToTransformIndex
116  (
117  globalTransforms_.transformIndex(info[i]),
118  patchi,
119  true, // patchi is sending side
120  tol // tolerance for comparison
121  )
122  );
123  }
124  return sendInfo;
125 }
126 
127 
128 void Foam::globalPoints::addToSend
129 (
130  const polyPatch& pp,
131  const label patchPointi,
132  const labelPairList& knownInfo,
133 
134  DynamicList<label>& patchFaces,
135  DynamicList<label>& indexInFace,
136  DynamicList<labelPairList>& allInfo
137 ) const
138 {
139  // Collect all topological information about a point on a patch. (this
140  // information is the patch faces using the point and the relative position
141  // of the point in the face)
142 
143  const label meshPointi = pp.meshPoints()[patchPointi];
144 
145  // Add all faces using the point so we are sure we find it on the
146  // other side.
147  const labelList& pFaces = pp.pointFaces()[patchPointi];
148 
149  for (const label patchFacei : pFaces)
150  {
151  const face& f = pp[patchFacei];
152 
153  patchFaces.append(patchFacei);
154  indexInFace.append(f.find(meshPointi));
155 
156  // Add patch transformation
157  allInfo.append(addSendTransform(pp.index(), knownInfo));
158  }
159 }
160 
161 
162 bool Foam::globalPoints::mergeInfo
163 (
164  const labelPairList& nbrInfo,
165  const label localPointi,
166  labelPairList& myInfo
167 ) const
168 {
169  // Add nbrInfo to myInfo. Return true if anything changed. nbrInfo is for a
170  // point a list of all the global points using it
171 
172  bool anyChanged = false;
173 
174  // Extend to make space for the nbrInfo (trimmed later)
175  labelPairList newInfo(myInfo);
176  label newI = newInfo.size();
177  newInfo.setSize(newI + nbrInfo.size());
178 
179  forAll(nbrInfo, i)
180  {
181  // Check if already have information about nbr point. There are two
182  // possibilities:
183  // - information found about same point but different transform.
184  // Combine transforms
185  // - information not found.
186 
187  label index = findSamePoint(myInfo, nbrInfo[i]);
188 
189  if (index == -1)
190  {
191  // New point
192  newInfo[newI++] = nbrInfo[i];
193  anyChanged = true;
194  }
195  else
196  {
197  // Same point. So we already have a connection between localPointi
198  // and the nbrIndex. Two situations:
199  // - same transform
200  // - one transform takes two steps, the other just a single.
201  if (myInfo[index] == nbrInfo[i])
202  {
203  // Everything same (so also transform). Nothing changed.
204  }
205  else
206  {
207  label myTransform = globalTransforms_.transformIndex
208  (
209  myInfo[index]
210  );
211  label nbrTransform = globalTransforms_.transformIndex
212  (
213  nbrInfo[i]
214  );
215 
216  // Different transform. See which is 'simplest'.
217  label minTransform = globalTransforms_.minimumTransformIndex
218  (
219  myTransform,
220  nbrTransform
221  );
222 
223  if (minTransform != myTransform)
224  {
225  // Use nbr info.
226  newInfo[index] = nbrInfo[i];
227  anyChanged = true;
228  }
229  }
230  }
231  }
232 
233  newInfo.setSize(newI);
234  myInfo.transfer(newInfo);
235 
236  return anyChanged;
237 }
238 
239 
240 Foam::label Foam::globalPoints::meshToLocalPoint
241 (
242  const Map<label>& meshToPatchPoint, // from mesh point to local numbering
243  const label meshPointi
244 )
245 {
246  return
247  (
248  meshToPatchPoint.size() == 0
249  ? meshPointi
250  : meshToPatchPoint[meshPointi]
251  );
252 }
253 
254 
255 Foam::label Foam::globalPoints::localToMeshPoint
256 (
257  const labelList& patchToMeshPoint,
258  const label localPointi
259 )
260 {
261  return
262  (
263  patchToMeshPoint.size() == 0
264  ? localPointi
265  : patchToMeshPoint[localPointi]
266  );
267 }
268 
269 
270 bool Foam::globalPoints::mergeInfo
271 (
272  const labelPairList& nbrInfo,
273  const label localPointi
274 )
275 {
276  // Updates database of current information on meshpoints with nbrInfo. Uses
277  // mergeInfo above. Returns true if data kept for meshPointi changed.
278 
279  label infoChanged = false;
280 
281  // Get the index into the procPoints list.
282  const auto iter = meshToProcPoint_.cfind(localPointi);
283 
284  if (iter.found())
285  {
286  if (mergeInfo(nbrInfo, localPointi, procPoints_[iter.val()]))
287  {
288  infoChanged = true;
289  }
290  }
291  else
292  {
293  // Construct local index for point
294  labelPairList knownInfo
295  (
296  1,
297  globalTransforms_.encode
298  (
300  localPointi,
301  globalTransforms_.nullTransformIndex()
302  )
303  );
304 
305  if (mergeInfo(nbrInfo, localPointi, knownInfo))
306  {
307  // Update addressing from into procPoints
308  meshToProcPoint_.insert(localPointi, procPoints_.size());
309  // Insert into list of equivalences.
310  procPoints_.append(knownInfo);
311 
312  infoChanged = true;
313  }
314  }
315  return infoChanged;
316 }
317 
318 
319 bool Foam::globalPoints::storeInitialInfo
320 (
321  const labelPairList& nbrInfo,
322  const label localPointi
323 )
324 {
325  // Updates database of current information on meshpoints with nbrInfo. Uses
326  // mergeInfo above. Returns true if data kept for meshPointi changed.
327 
328  label infoChanged = false;
329 
330  // Get the index into the procPoints list.
331  const auto iter = meshToProcPoint_.find(localPointi);
332 
333  if (iter.found())
334  {
335  if (mergeInfo(nbrInfo, localPointi, procPoints_[iter.val()]))
336  {
337  infoChanged = true;
338  }
339  }
340  else
341  {
342  // Update addressing into procPoints
343  meshToProcPoint_.insert(localPointi, procPoints_.size());
344  // Insert into list of equivalences.
345  procPoints_.append(nbrInfo);
346 
347  infoChanged = true;
348  }
349  return infoChanged;
350 }
351 
352 
353 void Foam::globalPoints::printProcPoint
354 (
355  const labelList& patchToMeshPoint,
356  const labelPair& pointInfo
357 ) const
358 {
359  label proci = globalTransforms_.processor(pointInfo);
360  label index = globalTransforms_.index(pointInfo);
361  label trafoI = globalTransforms_.transformIndex(pointInfo);
362 
363  Pout<< " proc:" << proci;
364  Pout<< " localpoint:";
365  Pout<< index;
366  Pout<< " through transform:"
367  << trafoI << " bits:"
368  << globalTransforms_.decodeTransformIndex(trafoI);
369 
370  if (proci == Pstream::myProcNo())
371  {
372  label meshPointi = localToMeshPoint(patchToMeshPoint, index);
373  Pout<< " at:" << mesh_.points()[meshPointi];
374  }
375 }
376 
377 
378 void Foam::globalPoints::printProcPoints
379 (
380  const labelList& patchToMeshPoint,
381  const labelPairList& pointInfo
382 ) const
383 {
384  forAll(pointInfo, i)
385  {
386  printProcPoint(patchToMeshPoint, pointInfo[i]);
387  Pout<< endl;
388  }
389 }
390 
391 
392 void Foam::globalPoints::initOwnPoints
393 (
394  const Map<label>& meshToPatchPoint,
395  const bool allPoints,
396  labelHashSet& changedPoints
397 )
398 {
399  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
400 
401  forAll(patches, patchi)
402  {
403  const polyPatch& pp = patches[patchi];
404 
405  if (pp.coupled())
406  {
407  const labelList& meshPoints = pp.meshPoints();
408 
409  if (allPoints)
410  {
411  // All points on patch
412  forAll(meshPoints, patchPointi)
413  {
414  label meshPointi = meshPoints[patchPointi];
415  label localPointi = meshToLocalPoint
416  (
417  meshToPatchPoint,
418  meshPointi
419  );
420 
421  labelPairList knownInfo
422  (
423  1,
424  globalTransforms_.encode
425  (
427  localPointi,
428  globalTransforms_.nullTransformIndex()
429  )
430  );
431 
432  //Pout<< "For point "<< pp.points()[meshPointi]
433  // << " inserting info " << knownInfo
434  // << endl;
435 
436  // Update changedpoints info.
437  if (storeInitialInfo(knownInfo, localPointi))
438  {
439  changedPoints.insert(localPointi);
440  }
441  }
442  }
443  else
444  {
445  // Boundary points only
446  const labelList& boundaryPoints = pp.boundaryPoints();
447 
448  forAll(boundaryPoints, i)
449  {
450  label meshPointi = meshPoints[boundaryPoints[i]];
451  label localPointi = meshToLocalPoint
452  (
453  meshToPatchPoint,
454  meshPointi
455  );
456 
457  labelPairList knownInfo
458  (
459  1,
460  globalTransforms_.encode
461  (
463  localPointi,
464  globalTransforms_.nullTransformIndex()
465  )
466  );
467 
468  if (storeInitialInfo(knownInfo, localPointi))
469  {
470  changedPoints.insert(localPointi);
471  }
472  }
473  }
474  }
475  }
476 }
477 
478 
479 void Foam::globalPoints::sendPatchPoints
480 (
481  const bool mergeSeparated,
482  const Map<label>& meshToPatchPoint,
483  PstreamBuffers& pBufs,
484  const labelHashSet& changedPoints
485 ) const
486 {
487  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
488  const labelPairList& patchInfo = globalTransforms_.patchTransformSign();
489 
490  forAll(patches, patchi)
491  {
492  const polyPatch& pp = patches[patchi];
493 
494  // mergeSeparated=true : send from all processor patches
495  // =false: send from ones without transform
496 
497  if
498  (
499  (Pstream::parRun() && isA<processorPolyPatch>(pp))
500  && (mergeSeparated || patchInfo[patchi].first() == -1)
501  )
502  {
503  const processorPolyPatch& procPatch =
504  refCast<const processorPolyPatch>(pp);
505 
506  // Information to send:
507  // patch face
508  DynamicList<label> patchFaces(pp.nPoints());
509  // index in patch face
510  DynamicList<label> indexInFace(pp.nPoints());
511  // all information I currently hold about this patchPoint
512  DynamicList<labelPairList> allInfo(pp.nPoints());
513 
514 
515  // Now collect information on all points mentioned in
516  // changedPoints. Note that these points only should occur on
517  // processorPatches (or rather this is a limitation!).
518 
519  const labelList& meshPoints = pp.meshPoints();
520 
521  forAll(meshPoints, patchPointi)
522  {
523  label meshPointi = meshPoints[patchPointi];
524  label localPointi = meshToLocalPoint
525  (
526  meshToPatchPoint,
527  meshPointi
528  );
529 
530  if (changedPoints.found(localPointi))
531  {
532  label index = meshToProcPoint_[localPointi];
533 
534  const labelPairList& knownInfo = procPoints_[index];
535 
536  // Add my information about localPointi to the
537  // send buffers. Encode the transformation
538  addToSend
539  (
540  pp,
541  patchPointi,
542  knownInfo,
543 
544  patchFaces,
545  indexInFace,
546  allInfo
547  );
548  }
549  }
550 
551  // Send to neighbour
552  if (debug)
553  {
554  Pout<< " Sending from " << pp.name() << " to "
555  << procPatch.neighbProcNo() << " point information:"
556  << patchFaces.size() << endl;
557  }
558 
559  UOPstream toNeighbour(procPatch.neighbProcNo(), pBufs);
560  toNeighbour << patchFaces << indexInFace << allInfo;
561  }
562  }
563 }
564 
565 
566 void Foam::globalPoints::receivePatchPoints
567 (
568  const bool mergeSeparated,
569  const Map<label>& meshToPatchPoint,
570  const labelList& patchToMeshPoint,
571  PstreamBuffers& pBufs,
572  labelHashSet& changedPoints
573 )
574 {
575  // Receive all my neighbours' information and merge with mine.
576  // After finishing will have updated
577  // - procPoints_ : all neighbour information merged in.
578  // - meshToProcPoint_
579  // - changedPoints: all points for which something changed.
580 
581  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
582  const labelPairList& patchInfo = globalTransforms_.patchTransformSign();
583 
584  // Reset changed points
585  changedPoints.clear();
586 
587  forAll(patches, patchi)
588  {
589  const polyPatch& pp = patches[patchi];
590 
591  if
592  (
593  (Pstream::parRun() && isA<processorPolyPatch>(pp))
594  && (mergeSeparated || patchInfo[patchi].first() == -1)
595  )
596  {
597  const processorPolyPatch& procPatch =
598  refCast<const processorPolyPatch>(pp);
599 
600  labelList patchFaces;
601  labelList indexInFace;
602  List<labelPairList> nbrInfo;
603 
604  {
605  UIPstream fromNeighbour(procPatch.neighbProcNo(), pBufs);
606  fromNeighbour >> patchFaces >> indexInFace >> nbrInfo;
607  }
608 
609  if (debug)
610  {
611  Pout<< " On " << pp.name()
612  << " Received from "
613  << procPatch.neighbProcNo() << " point information:"
614  << patchFaces.size() << endl;
615  }
616 
617  forAll(patchFaces, i)
618  {
619  const face& f = pp[patchFaces[i]];
620 
621  // Get index in this face from index on face on other side.
622  label index = (f.size() - indexInFace[i]) % f.size();
623 
624  // Get the meshpoint on my side
625  label meshPointi = f[index];
626 
627  label localPointi = meshToLocalPoint
628  (
629  meshToPatchPoint,
630  meshPointi
631  );
632 
633  if (mergeInfo(nbrInfo[i], localPointi))
634  {
635  changedPoints.insert(localPointi);
636  }
637  }
638  }
639  else if
640  (
641  (
642  isA<cyclicPolyPatch>(pp)
643  && refCast<const cyclicPolyPatch>(pp).owner()
644  )
645  && (mergeSeparated || patchInfo[patchi].first() == -1)
646  )
647  {
648  // Handle cyclics: send lower half to upper half and vice versa.
649  // Or since they both are in memory just do it point by point.
650 
651  const cyclicPolyPatch& cycPatch =
652  refCast<const cyclicPolyPatch>(pp);
653 
654  //Pout<< "Patch:" << patchi << " name:" << pp.name() << endl;
655 
656  const labelList& meshPoints = pp.meshPoints();
657  const labelList coupledMeshPoints(reverseMeshPoints(cycPatch));
658 
659  forAll(meshPoints, i)
660  {
661  label meshPointA = meshPoints[i];
662  label meshPointB = coupledMeshPoints[i];
663 
664  if (meshPointA != meshPointB)
665  {
666  //Pout<< "Connection between point " << meshPointA
667  // << " at " << mesh_.points()[meshPointA]
668  // << " and " << meshPointB
669  // << " at " << mesh_.points()[meshPointB] << endl;
670 
671  label localA = meshToLocalPoint
672  (
673  meshToPatchPoint,
674  meshPointA
675  );
676  label localB = meshToLocalPoint
677  (
678  meshToPatchPoint,
679  meshPointB
680  );
681 
682 
683  // Do we have information on pointA?
684  const auto procPointA = meshToProcPoint_.cfind(localA);
685 
686  if (procPointA.found())
687  {
688  const labelPairList infoA = addSendTransform
689  (
690  cycPatch.index(),
691  procPoints_[procPointA()]
692  );
693 
694  if (mergeInfo(infoA, localB))
695  {
696  changedPoints.insert(localB);
697  }
698  }
699 
700  // Same for info on pointB
701  const auto procPointB = meshToProcPoint_.cfind(localB);
702 
703  if (procPointB.found())
704  {
705  const labelPairList infoB = addSendTransform
706  (
707  cycPatch.neighbPatchID(),
708  procPoints_[procPointB()]
709  );
710 
711  if (mergeInfo(infoB, localA))
712  {
713  changedPoints.insert(localA);
714  }
715  }
716  }
717  }
718  }
719  }
720 }
721 
722 
723 void Foam::globalPoints::remove
724 (
725  const labelList& patchToMeshPoint,
726  const Map<label>& directNeighbours
727 )
728 {
729  // Remove entries which are handled by normal face-face communication. I.e.
730  // those points where the equivalence list is only me and my (face)neighbour
731 
732  // Save old ones.
733  Map<label> oldMeshToProcPoint(std::move(meshToProcPoint_));
734  meshToProcPoint_.resize(oldMeshToProcPoint.size());
735  DynamicList<labelPairList> oldProcPoints(std::move(procPoints_));
736  procPoints_.setCapacity(oldProcPoints.size());
737 
738  // Go through all equivalences
739  forAllConstIters(oldMeshToProcPoint, iter)
740  {
741  const label localPointi = iter.key();
742  const labelPairList& pointInfo = oldProcPoints[iter.val()];
743 
744  if (pointInfo.size() == 2)
745  {
746  // I will be in this equivalence list.
747  // Check whether my direct (=face) neighbour
748  // is in it. This would be an ordinary connection and can be
749  // handled by normal face-face connectivity.
750 
751  label proc0 = globalTransforms_.processor(pointInfo[0]);
752  label proc1 = globalTransforms_.processor(pointInfo[1]);
753 
754  if
755  (
756  (
757  proc0 == Pstream::myProcNo()
758  && directNeighbours.found
759  (
760  globalTransforms_.index(pointInfo[0])
761  )
762  )
763  || (
764  proc1 == Pstream::myProcNo()
765  && directNeighbours.found
766  (
767  globalTransforms_.index(pointInfo[1])
768  )
769  )
770  )
771  {
772  // Normal faceNeighbours
773  if (proc0 == Pstream::myProcNo())
774  {
775  //Pout<< "Removing direct neighbour:"
776  // << mesh_.points()
777  // [globalTransforms_.index(pointInfo[0])]
778  // << endl;
779  }
780  else if (proc1 == Pstream::myProcNo())
781  {
782  //Pout<< "Removing direct neighbour:"
783  // << mesh_.points()
784  // [globalTransforms_.index(pointInfo[1])]
785  // << endl;
786  }
787  }
788  else
789  {
790  // This condition will be very rare: points are used by
791  // two processors which are not face-face connected.
792  // e.g.
793  // +------+------+
794  // | wall | B |
795  // +------+------+
796  // | A | wall |
797  // +------+------+
798  // Processor A and B share a point. Note that this only will
799  // be found if the two domains are face connected at all
800  // (not shown in the picture)
801 
802  meshToProcPoint_.insert(localPointi, procPoints_.size());
803  procPoints_.append(pointInfo);
804  }
805  }
806  else if (pointInfo.size() == 1)
807  {
808  // This happens for 'wedge' like cyclics where the two halves
809  // come together in the same point so share the same meshPoint.
810  // So this meshPoint will have info of size one only.
811  if
812  (
813  globalTransforms_.processor(pointInfo[0])
814  != Pstream::myProcNo()
815  || !directNeighbours.found
816  (
817  globalTransforms_.index(pointInfo[0])
818  )
819  )
820  {
821  meshToProcPoint_.insert(localPointi, procPoints_.size());
822  procPoints_.append(pointInfo);
823  }
824  }
825  else
826  {
827  meshToProcPoint_.insert(localPointi, procPoints_.size());
828  procPoints_.append(pointInfo);
829  }
830  }
831 
832  procPoints_.shrink();
833  meshToProcPoint_.resize(2*procPoints_.size());
834 }
835 
836 
837 Foam::labelList Foam::globalPoints::reverseMeshPoints
838 (
839  const cyclicPolyPatch& pp
840 )
841 {
842  const cyclicPolyPatch& nbrPatch = pp.neighbPatch();
843 
844  faceList masterFaces(nbrPatch.size());
845 
846  forAll(nbrPatch, facei)
847  {
848  masterFaces[facei] = nbrPatch[facei].reverseFace();
849  }
850 
851  return primitiveFacePatch
852  (
853  masterFaces,
854  nbrPatch.points()
855  ).meshPoints();
856 }
857 
858 
859 void Foam::globalPoints::calculateSharedPoints
860 (
861  const Map<label>& meshToPatchPoint, // from mesh point to local numbering
862  const labelList& patchToMeshPoint, // from local numbering to mesh point
863  const bool keepAllPoints,
864  const bool mergeSeparated
865 )
866 {
867  if (debug)
868  {
869  Pout<< "globalPoints::calculateSharedPoints(..) : "
870  << "doing processor to processor communication to get sharedPoints"
871  << endl
872  << " keepAllPoints :" << keepAllPoints << endl
873  << " mergeSeparated:" << mergeSeparated << endl
874  << endl;
875  }
876 
877 
878  labelHashSet changedPoints(2*nPatchPoints_);
879 
880  // Initialize procPoints with my patch points. Keep track of points
881  // inserted (in changedPoints)
882  // There are two possible forms of this:
883  // - initialize with all patch points (allPoints = true). This causes all
884  // patch points to be exchanged so a lot of information gets stored and
885  // transferred. This all gets filtered out later when removing the
886  // equivalence lists of size 2.
887  // - initialize with boundary points of patches only (allPoints = false).
888  // This should work for all decompositions except extreme ones where a
889  // shared point is not on the boundary of any processor patches using it.
890  // This would happen if a domain was pinched such that two patches share
891  // a point or edge.
892  initOwnPoints(meshToPatchPoint, true, changedPoints);
893 
894  // Do one exchange iteration to get neighbour points.
895  {
896  // Note: to use 'scheduled' would have to intersperse send and receive.
897  // So for now just use nonBlocking. Also globalPoints itself gets
898  // constructed by mesh.globalData().patchSchedule() so creates a loop.
899  PstreamBuffers pBufs
900  (
901  (
905  )
906  );
907  sendPatchPoints
908  (
909  mergeSeparated,
910  meshToPatchPoint,
911  pBufs,
912  changedPoints
913  );
914  pBufs.finishedSends();
915  receivePatchPoints
916  (
917  mergeSeparated,
918  meshToPatchPoint,
919  patchToMeshPoint,
920  pBufs,
921  changedPoints
922  );
923  }
924 
925  // Save neighbours reachable through face-face communication.
926  Map<label> neighbourList;
927  if (!keepAllPoints)
928  {
929  neighbourList = meshToProcPoint_;
930  }
931 
932  // Exchange until nothing changes on all processors.
933  bool changed = false;
934 
935  do
936  {
937  PstreamBuffers pBufs
938  (
939  (
943  )
944  );
945  sendPatchPoints
946  (
947  mergeSeparated,
948  meshToPatchPoint,
949  pBufs,
950  changedPoints
951  );
952  pBufs.finishedSends();
953  receivePatchPoints
954  (
955  mergeSeparated,
956  meshToPatchPoint,
957  patchToMeshPoint,
958  pBufs,
959  changedPoints
960  );
961 
962  changed = changedPoints.size() > 0;
963  reduce(changed, orOp<bool>());
964 
965  } while (changed);
966 
967 
968  //Pout<< "**ALL** connected points:" << endl;
969  //forAllConstIters(meshToProcPoint_, iter)
970  //{
971  // label localI = iter.key();
972  // const labelPairList& pointInfo = procPoints_[iter.val()];
973  // Pout<< "pointi:" << localI << " index:" << iter.val()
974  // << " coord:"
975  // << mesh_.points()[localToMeshPoint(patchToMeshPoint, localI)]
976  // << endl;
977  // printProcPoints(patchToMeshPoint, pointInfo);
978  // Pout<< endl;
979  //}
980 
981 
982  // Remove direct neighbours from point equivalences.
983  if (!keepAllPoints)
984  {
985  remove(patchToMeshPoint, neighbourList);
986  }
987 
988 
989  // Sort procPoints in incremental order. This will make
990  // the master the first element on all processors.
991  // Note: why not sort in decreasing order? Give more work to higher
992  // processors.
993  forAllConstIters(meshToProcPoint_, iter)
994  {
995  labelPairList& pointInfo = procPoints_[iter.val()];
996  sort(pointInfo, globalIndexAndTransform::less(globalTransforms_));
997  }
998 
999 
1000  // We now have - in procPoints_ - a list of points which are shared between
1001  // multiple processors. Filter into non-transformed and transformed
1002  // connections.
1003 
1004  pointPoints_.setSize(globalIndices_.localSize());
1005  List<labelPairList> transformedPoints(globalIndices_.localSize());
1006 
1007  forAllConstIters(meshToProcPoint_, iter)
1008  {
1009  const labelPairList& pointInfo = procPoints_[iter.val()];
1010 
1011  if (pointInfo.size() >= 2)
1012  {
1013  // Since sorted master point is the first element
1014  const labelPair& masterInfo = pointInfo[0];
1015 
1016  if
1017  (
1018  (
1019  globalTransforms_.processor(masterInfo)
1020  == Pstream::myProcNo()
1021  )
1022  && (globalTransforms_.index(masterInfo) == iter.key())
1023  )
1024  {
1025  labelList& pPoints = pointPoints_[iter.key()];
1026  pPoints.setSize(pointInfo.size()-1);
1027 
1028  labelPairList& trafoPPoints = transformedPoints[iter.key()];
1029  trafoPPoints.setSize(pointInfo.size()-1);
1030 
1031  label nonTransformI = 0;
1032  label transformI = 0;
1033 
1034  for (label i = 1; i < pointInfo.size(); i++)
1035  {
1036  const labelPair& info = pointInfo[i];
1037  label proci = globalTransforms_.processor(info);
1038  label index = globalTransforms_.index(info);
1039  label transform = globalTransforms_.transformIndex
1040  (
1041  info
1042  );
1043 
1044  if (transform == globalTransforms_.nullTransformIndex())
1045  {
1046  pPoints[nonTransformI++] = globalIndices_.toGlobal
1047  (
1048  proci,
1049  index
1050  );
1051  }
1052  else
1053  {
1054  trafoPPoints[transformI++] = info;
1055  }
1056  }
1057 
1058  pPoints.setSize(nonTransformI);
1059  trafoPPoints.setSize(transformI);
1060  }
1061  }
1062  }
1063 
1064 
1065  List<Map<label>> compactMap;
1066  map_.reset
1067  (
1068  new mapDistribute
1069  (
1070  globalIndices_,
1071  pointPoints_,
1072 
1073  globalTransforms_,
1074  transformedPoints,
1075  transformedPointPoints_,
1076 
1077  compactMap
1078  )
1079  );
1080 
1081  if (debug)
1082  {
1083  Pout<< "globalPoints::calculateSharedPoints(..) : "
1084  << "Finished global points" << endl;
1085  }
1086 }
1087 
1088 
1089 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1090 
1091 Foam::globalPoints::globalPoints
1093  const polyMesh& mesh,
1094  const bool keepAllPoints,
1095  const bool mergeSeparated
1096 )
1097 :
1098  mesh_(mesh),
1099  globalIndices_(mesh_.nPoints()),
1100  globalTransforms_(mesh),
1101  nPatchPoints_(countPatchPoints(mesh.boundaryMesh())),
1102  procPoints_(nPatchPoints_),
1103  meshToProcPoint_(nPatchPoints_)
1104 {
1105  // Empty patch maps to signal storing mesh point labels
1106  Map<label> meshToPatchPoint(0);
1107  labelList patchToMeshPoint(0);
1108 
1109  calculateSharedPoints
1110  (
1111  meshToPatchPoint,
1112  patchToMeshPoint,
1113  keepAllPoints,
1114  mergeSeparated
1115  );
1116 }
1117 
1118 
1119 Foam::globalPoints::globalPoints
1121  const polyMesh& mesh,
1122  const indirectPrimitivePatch& coupledPatch,
1123  const bool keepAllPoints,
1124  const bool mergeSeparated
1125 )
1126 :
1127  mesh_(mesh),
1128  globalIndices_(coupledPatch.nPoints()),
1129  globalTransforms_(mesh),
1130  nPatchPoints_(coupledPatch.nPoints()),
1131  procPoints_(nPatchPoints_),
1132  meshToProcPoint_(nPatchPoints_)
1133 {
1134  calculateSharedPoints
1135  (
1136  coupledPatch.meshPointMap(),
1137  coupledPatch.meshPoints(),
1138  keepAllPoints,
1139  mergeSeparated
1140  );
1141 }
1142 
1143 
1144 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
cyclicPolyPatch.H
Foam::primitiveFacePatch
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
Definition: primitiveFacePatch.H:51
Foam::less
static bool less(const vector &x, const vector &y)
To compare normals.
Definition: meshRefinementRefine.C:57
Foam::labelPairList
List< labelPair > labelPairList
List of labelPairs.
Definition: labelPair.H:64
Foam::Map< label >
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::UPstream::defaultCommsType
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:281
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::transform
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
polyMesh.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::labelPair
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:54
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::OSstream::name
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::sort
void sort(UList< T > &a)
Definition: UList.C:261
Foam::PrimitivePatch::nPoints
label nPoints() const
Number of points supporting patch faces.
Definition: PrimitivePatch.H:316
processorPolyPatch.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::UPstream::commsTypes::nonBlocking
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
mapDistribute.H
Foam::UPstream::commsTypes::scheduled
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::PrimitivePatch::meshPointMap
const Map< label > & meshPointMap() const
Mesh point map.
Definition: PrimitivePatch.C:343
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
Foam::PrimitivePatch::meshPoints
const labelList & meshPoints() const
Return labelList of mesh points in patch.
Definition: PrimitivePatch.C:330
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::HEX);labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells].reset(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Foam::dimensionSet::reset
void reset(const dimensionSet &ds)
Copy assign the exponents from the dimensionSet.
Definition: dimensionSet.C:149
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatch.H:79
globalPoints.H