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-2022 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  bool infoChanged = false;
280 
281  // Get the index into the procPoints list.
282  const auto iter = meshToProcPoint_.cfind(localPointi);
283 
284  if (iter.good())
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  bool infoChanged = false;
329 
330  // Get the index into the procPoints list.
331  const auto iter = meshToProcPoint_.find(localPointi);
332 
333  if (iter.good())
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  // Reduce communication by only sending non-zero data,
491  // but with multiply-connected processor/processor
492  // (eg, processorCyclic) also need to send zero information
493  // to keep things synchronised
494 
495  // Reset buffers, initialize for registerSend() bookkeeping
496  pBufs.clear();
497  pBufs.initRegisterSend();
498 
499 
500  // Information to send:
501 
502  // The patch face
503  DynamicList<label> patchFaces;
504 
505  // Index in patch face
506  DynamicList<label> indexInFace;
507 
508  // All information I currently hold about the patchPoint
509  DynamicList<labelPairList> allInfo;
510 
511 
512  forAll(patches, patchi)
513  {
514  const polyPatch& pp = patches[patchi];
515 
516  // mergeSeparated=true : send from all processor patches
517  // =false: send from ones without transform
518 
519  if
520  (
521  (UPstream::parRun() && isA<processorPolyPatch>(pp))
522  && (mergeSeparated || patchInfo[patchi].first() == -1)
523  )
524  {
525  const auto& procPatch = refCast<const processorPolyPatch>(pp);
526  const label nbrProci = procPatch.neighbProcNo();
527 
528  patchFaces.clear();
529  patchFaces.reserve(pp.nPoints());
530 
531  indexInFace.clear();
532  indexInFace.reserve(pp.nPoints());
533 
534  allInfo.clear();
535  allInfo.reserve(pp.nPoints());
536 
537  // Now collect information on all points mentioned in
538  // changedPoints. Note that these points only should occur on
539  // processorPatches (or rather this is a limitation!).
540 
541  const labelList& meshPoints = pp.meshPoints();
542 
543  forAll(meshPoints, patchPointi)
544  {
545  label meshPointi = meshPoints[patchPointi];
546  label localPointi = meshToLocalPoint
547  (
548  meshToPatchPoint,
549  meshPointi
550  );
551 
552  if (changedPoints.found(localPointi))
553  {
554  label index = meshToProcPoint_[localPointi];
555 
556  const labelPairList& knownInfo = procPoints_[index];
557 
558  // Add my information about localPointi to the
559  // send buffers. Encode the transformation
560  addToSend
561  (
562  pp,
563  patchPointi,
564  knownInfo,
565 
566  patchFaces,
567  indexInFace,
568  allInfo
569  );
570  }
571  }
572 
573 
574  // Send to neighbour
575  {
576  UOPstream toNbr(nbrProci, pBufs);
577  toNbr << patchFaces << indexInFace << allInfo;
578 
579  // Record if send is required (data are non-zero)
580  pBufs.registerSend(nbrProci, !patchFaces.empty());
581 
582  if (debug)
583  {
584  Pout<< "Sending from " << pp.name() << " to proc:"
585  << nbrProci << " point information:"
586  << patchFaces.size() << endl;
587  }
588  }
589  }
590  }
591 
592  // Discard unnecessary (unregistered) sends
593  pBufs.clearUnregistered();
594 }
595 
596 
597 void Foam::globalPoints::receivePatchPoints
598 (
599  const bool mergeSeparated,
600  const Map<label>& meshToPatchPoint,
601  const labelList& patchToMeshPoint,
602  PstreamBuffers& pBufs,
603  labelHashSet& changedPoints
604 )
605 {
606  // Receive all my neighbours' information and merge with mine.
607  // After finishing will have updated
608  // - procPoints_ : all neighbour information merged in.
609  // - meshToProcPoint_
610  // - changedPoints: all points for which something changed.
611 
612  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
613  const labelPairList& patchInfo = globalTransforms_.patchTransformSign();
614 
615  // Reset changed points
616  changedPoints.clear();
617 
618  forAll(patches, patchi)
619  {
620  const polyPatch& pp = patches[patchi];
621 
622  if
623  (
624  (UPstream::parRun() && isA<processorPolyPatch>(pp))
625  && (mergeSeparated || patchInfo[patchi].first() == -1)
626  )
627  {
628  const auto& procPatch = refCast<const processorPolyPatch>(pp);
629  const label nbrProci = procPatch.neighbProcNo();
630 
631  if (!pBufs.recvDataCount(nbrProci))
632  {
633  // Nothing to receive
634  continue;
635  }
636 
637  labelList patchFaces;
638  labelList indexInFace;
639  List<labelPairList> nbrInfo;
640 
641  {
642  UIPstream fromNbr(nbrProci, pBufs);
643  fromNbr >> patchFaces >> indexInFace >> nbrInfo;
644  }
645 
646  if (debug)
647  {
648  Pout<< " On " << pp.name()
649  << " Received from "
650  << nbrProci << " point information:"
651  << patchFaces.size() << endl;
652  }
653 
654  forAll(patchFaces, i)
655  {
656  const face& f = pp[patchFaces[i]];
657 
658  // Get index in this face from index on face on other side.
659  label index = (f.size() - indexInFace[i]) % f.size();
660 
661  // Get the meshpoint on my side
662  label meshPointi = f[index];
663 
664  label localPointi = meshToLocalPoint
665  (
666  meshToPatchPoint,
667  meshPointi
668  );
669 
670  if (mergeInfo(nbrInfo[i], localPointi))
671  {
672  changedPoints.insert(localPointi);
673  }
674  }
675  }
676  else if
677  (
678  (
679  isA<cyclicPolyPatch>(pp)
680  && refCast<const cyclicPolyPatch>(pp).owner()
681  )
682  && (mergeSeparated || patchInfo[patchi].first() == -1)
683  )
684  {
685  // Handle cyclics: send lower half to upper half and vice versa.
686  // Or since they both are in memory just do it point by point.
687 
688  const cyclicPolyPatch& cycPatch =
689  refCast<const cyclicPolyPatch>(pp);
690 
691  //Pout<< "Patch:" << patchi << " name:" << pp.name() << endl;
692 
693  const labelList& meshPoints = pp.meshPoints();
694  const labelList coupledMeshPoints(reverseMeshPoints(cycPatch));
695 
696  forAll(meshPoints, i)
697  {
698  label meshPointA = meshPoints[i];
699  label meshPointB = coupledMeshPoints[i];
700 
701  if (meshPointA != meshPointB)
702  {
703  //Pout<< "Connection between point " << meshPointA
704  // << " at " << mesh_.points()[meshPointA]
705  // << " and " << meshPointB
706  // << " at " << mesh_.points()[meshPointB] << endl;
707 
708  label localA = meshToLocalPoint
709  (
710  meshToPatchPoint,
711  meshPointA
712  );
713  label localB = meshToLocalPoint
714  (
715  meshToPatchPoint,
716  meshPointB
717  );
718 
719 
720  // Do we have information on pointA?
721  const auto procPointA = meshToProcPoint_.cfind(localA);
722 
723  if (procPointA.good())
724  {
725  const labelPairList infoA = addSendTransform
726  (
727  cycPatch.index(),
728  procPoints_[procPointA()]
729  );
730 
731  if (mergeInfo(infoA, localB))
732  {
733  changedPoints.insert(localB);
734  }
735  }
736 
737  // Same for info on pointB
738  const auto procPointB = meshToProcPoint_.cfind(localB);
739 
740  if (procPointB.good())
741  {
742  const labelPairList infoB = addSendTransform
743  (
744  cycPatch.neighbPatchID(),
745  procPoints_[procPointB()]
746  );
747 
748  if (mergeInfo(infoB, localA))
749  {
750  changedPoints.insert(localA);
751  }
752  }
753  }
754  }
755  }
756  }
757 }
758 
759 
760 void Foam::globalPoints::remove
761 (
762  const labelList& patchToMeshPoint,
763  const Map<label>& directNeighbours
764 )
765 {
766  // Remove entries which are handled by normal face-face communication. I.e.
767  // those points where the equivalence list is only me and my (face)neighbour
768 
769  // Save old ones.
770  Map<label> oldMeshToProcPoint(std::move(meshToProcPoint_));
771  meshToProcPoint_.reserve(oldMeshToProcPoint.size());
772  DynamicList<labelPairList> oldProcPoints(std::move(procPoints_));
773  procPoints_.setCapacity(oldProcPoints.size());
774 
775  // Go through all equivalences
776  forAllConstIters(oldMeshToProcPoint, iter)
777  {
778  const label localPointi = iter.key();
779  const labelPairList& pointInfo = oldProcPoints[iter.val()];
780 
781  if (pointInfo.size() == 2)
782  {
783  // I will be in this equivalence list.
784  // Check whether my direct (=face) neighbour
785  // is in it. This would be an ordinary connection and can be
786  // handled by normal face-face connectivity.
787 
788  label proc0 = globalTransforms_.processor(pointInfo[0]);
789  label proc1 = globalTransforms_.processor(pointInfo[1]);
790 
791  if
792  (
793  (
794  proc0 == Pstream::myProcNo()
795  && directNeighbours.found
796  (
797  globalTransforms_.index(pointInfo[0])
798  )
799  )
800  || (
801  proc1 == Pstream::myProcNo()
802  && directNeighbours.found
803  (
804  globalTransforms_.index(pointInfo[1])
805  )
806  )
807  )
808  {
809  // Normal faceNeighbours
810  if (proc0 == Pstream::myProcNo())
811  {
812  //Pout<< "Removing direct neighbour:"
813  // << mesh_.points()
814  // [globalTransforms_.index(pointInfo[0])]
815  // << endl;
816  }
817  else if (proc1 == Pstream::myProcNo())
818  {
819  //Pout<< "Removing direct neighbour:"
820  // << mesh_.points()
821  // [globalTransforms_.index(pointInfo[1])]
822  // << endl;
823  }
824  }
825  else
826  {
827  // This condition will be very rare: points are used by
828  // two processors which are not face-face connected.
829  // e.g.
830  // +------+------+
831  // | wall | B |
832  // +------+------+
833  // | A | wall |
834  // +------+------+
835  // Processor A and B share a point. Note that this only will
836  // be found if the two domains are face connected at all
837  // (not shown in the picture)
838 
839  meshToProcPoint_.insert(localPointi, procPoints_.size());
840  procPoints_.append(pointInfo);
841  }
842  }
843  else if (pointInfo.size() == 1)
844  {
845  // This happens for 'wedge' like cyclics where the two halves
846  // come together in the same point so share the same meshPoint.
847  // So this meshPoint will have info of size one only.
848  if
849  (
850  globalTransforms_.processor(pointInfo[0])
851  != Pstream::myProcNo()
852  || !directNeighbours.found
853  (
854  globalTransforms_.index(pointInfo[0])
855  )
856  )
857  {
858  meshToProcPoint_.insert(localPointi, procPoints_.size());
859  procPoints_.append(pointInfo);
860  }
861  }
862  else
863  {
864  meshToProcPoint_.insert(localPointi, procPoints_.size());
865  procPoints_.append(pointInfo);
866  }
867  }
868 
869  procPoints_.shrink();
870  meshToProcPoint_.reserve(procPoints_.size());
871 }
872 
873 
874 Foam::labelList Foam::globalPoints::reverseMeshPoints
875 (
876  const cyclicPolyPatch& pp
877 )
878 {
879  const cyclicPolyPatch& nbrPatch = pp.neighbPatch();
880 
881  faceList masterFaces(nbrPatch.size());
882 
883  forAll(nbrPatch, facei)
884  {
885  masterFaces[facei] = nbrPatch[facei].reverseFace();
886  }
887 
888  return primitiveFacePatch
889  (
890  masterFaces,
891  nbrPatch.points()
892  ).meshPoints();
893 }
894 
895 
896 void Foam::globalPoints::calculateSharedPoints
897 (
898  const Map<label>& meshToPatchPoint, // from mesh point to local numbering
899  const labelList& patchToMeshPoint, // from local numbering to mesh point
900  const bool keepAllPoints,
901  const bool mergeSeparated
902 )
903 {
904  if (debug)
905  {
906  Pout<< "globalPoints::calculateSharedPoints(..) : "
907  << "doing processor to processor communication to get sharedPoints"
908  << endl
909  << " keepAllPoints :" << keepAllPoints << endl
910  << " mergeSeparated:" << mergeSeparated << endl
911  << endl;
912  }
913 
914 
915  labelHashSet changedPoints(2*nPatchPoints_);
916 
917  // Initialize procPoints with my patch points. Keep track of points
918  // inserted (in changedPoints)
919  // There are two possible forms of this:
920  // - initialize with all patch points (allPoints = true). This causes all
921  // patch points to be exchanged so a lot of information gets stored and
922  // transferred. This all gets filtered out later when removing the
923  // equivalence lists of size 2.
924  // - initialize with boundary points of patches only (allPoints = false).
925  // This should work for all decompositions except extreme ones where a
926  // shared point is not on the boundary of any processor patches using it.
927  // This would happen if a domain was pinched such that two patches share
928  // a point or edge.
929  initOwnPoints(meshToPatchPoint, true, changedPoints);
930 
931  // Note: to use 'scheduled' would have to intersperse send and receive.
932  // So for now just use nonBlocking. Also globalPoints itself gets
933  // constructed by mesh.globalData().patchSchedule() so creates a loop.
934  PstreamBuffers pBufs
935  (
936  (
940  )
941  );
942 
943  // Don't clear storage on persistent buffer
944  pBufs.allowClearRecv(false);
945 
946  // Do one exchange iteration to get neighbour points.
947  {
948  sendPatchPoints
949  (
950  mergeSeparated,
951  meshToPatchPoint,
952  pBufs,
953  changedPoints
954  );
955 
956  pBufs.finishedSends();
957 
958  receivePatchPoints
959  (
960  mergeSeparated,
961  meshToPatchPoint,
962  patchToMeshPoint,
963  pBufs,
964  changedPoints
965  );
966  }
967 
968  // Save neighbours reachable through face-face communication.
969  Map<label> neighbourList;
970  if (!keepAllPoints)
971  {
972  neighbourList = meshToProcPoint_;
973  }
974 
975  // Exchange until nothing changes on any processors.
976 
977  do
978  {
979  sendPatchPoints
980  (
981  mergeSeparated,
982  meshToPatchPoint,
983  pBufs,
984  changedPoints
985  );
986 
987  pBufs.finishedSends();
988 
989  receivePatchPoints
990  (
991  mergeSeparated,
992  meshToPatchPoint,
993  patchToMeshPoint,
994  pBufs,
995  changedPoints
996  );
997 
998  } while (returnReduceOr(changedPoints.size()));
999 
1000 
1001  //Pout<< "**ALL** connected points:" << endl;
1002  //forAllConstIters(meshToProcPoint_, iter)
1003  //{
1004  // label localI = iter.key();
1005  // const labelPairList& pointInfo = procPoints_[iter.val()];
1006  // Pout<< "pointi:" << localI << " index:" << iter.val()
1007  // << " coord:"
1008  // << mesh_.points()[localToMeshPoint(patchToMeshPoint, localI)]
1009  // << endl;
1010  // printProcPoints(patchToMeshPoint, pointInfo);
1011  // Pout<< endl;
1012  //}
1013 
1014 
1015  // Remove direct neighbours from point equivalences.
1016  if (!keepAllPoints)
1017  {
1018  remove(patchToMeshPoint, neighbourList);
1019  }
1020 
1021 
1022  // Sort procPoints in incremental order. This will make
1023  // the master the first element on all processors.
1024  // Note: why not sort in decreasing order? Give more work to higher
1025  // processors.
1026  forAllConstIters(meshToProcPoint_, iter)
1027  {
1028  labelPairList& pointInfo = procPoints_[iter.val()];
1029  sort(pointInfo, globalIndexAndTransform::less(globalTransforms_));
1030  }
1031 
1032 
1033  // We now have - in procPoints_ - a list of points which are shared between
1034  // multiple processors. Filter into non-transformed and transformed
1035  // connections.
1036 
1037  pointPoints_.setSize(globalIndices_.localSize());
1038  List<labelPairList> transformedPoints(globalIndices_.localSize());
1039 
1040  forAllConstIters(meshToProcPoint_, iter)
1041  {
1042  const labelPairList& pointInfo = procPoints_[iter.val()];
1043 
1044  if (pointInfo.size() >= 2)
1045  {
1046  // Since sorted master point is the first element
1047  const labelPair& masterInfo = pointInfo[0];
1048 
1049  if
1050  (
1051  (
1052  globalTransforms_.processor(masterInfo)
1053  == Pstream::myProcNo()
1054  )
1055  && (globalTransforms_.index(masterInfo) == iter.key())
1056  )
1057  {
1058  labelList& pPoints = pointPoints_[iter.key()];
1059  pPoints.setSize(pointInfo.size()-1);
1060 
1061  labelPairList& trafoPPoints = transformedPoints[iter.key()];
1062  trafoPPoints.setSize(pointInfo.size()-1);
1063 
1064  label nonTransformI = 0;
1065  label transformI = 0;
1066 
1067  for (label i = 1; i < pointInfo.size(); i++)
1068  {
1069  const labelPair& info = pointInfo[i];
1070  label proci = globalTransforms_.processor(info);
1071  label index = globalTransforms_.index(info);
1072  label transform = globalTransforms_.transformIndex
1073  (
1074  info
1075  );
1076 
1077  if (transform == globalTransforms_.nullTransformIndex())
1078  {
1079  pPoints[nonTransformI++] = globalIndices_.toGlobal
1080  (
1081  proci,
1082  index
1083  );
1084  }
1085  else
1086  {
1087  trafoPPoints[transformI++] = info;
1088  }
1089  }
1090 
1091  pPoints.setSize(nonTransformI);
1092  trafoPPoints.setSize(transformI);
1093  }
1094  }
1095  }
1096 
1097 
1098  List<Map<label>> compactMap;
1099  map_.reset
1100  (
1101  new mapDistribute
1102  (
1103  globalIndices_,
1104  pointPoints_,
1105 
1106  globalTransforms_,
1107  transformedPoints,
1108  transformedPointPoints_,
1109 
1110  compactMap
1111  )
1112  );
1113 
1114  if (debug)
1115  {
1116  Pout<< "globalPoints::calculateSharedPoints(..) : "
1117  << "Finished global points" << endl;
1118  }
1119 }
1120 
1121 
1122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1123 
1124 Foam::globalPoints::globalPoints
1125 (
1126  const polyMesh& mesh,
1127  const bool keepAllPoints,
1128  const bool mergeSeparated
1129 )
1130 :
1131  mesh_(mesh),
1132  globalIndices_(mesh_.nPoints()),
1133  globalTransforms_(mesh),
1134  nPatchPoints_(countPatchPoints(mesh.boundaryMesh())),
1135  procPoints_(nPatchPoints_),
1136  meshToProcPoint_(nPatchPoints_)
1137 {
1138  // Empty patch maps to signal storing mesh point labels
1139  Map<label> meshToPatchPoint(0);
1140  labelList patchToMeshPoint(0);
1141 
1142  calculateSharedPoints
1143  (
1144  meshToPatchPoint,
1145  patchToMeshPoint,
1146  keepAllPoints,
1147  mergeSeparated
1148  );
1149 }
1150 
1151 
1152 Foam::globalPoints::globalPoints
1153 (
1154  const polyMesh& mesh,
1155  const indirectPrimitivePatch& coupledPatch,
1156  const bool keepAllPoints,
1157  const bool mergeSeparated
1158 )
1159 :
1160  mesh_(mesh),
1161  globalIndices_(coupledPatch.nPoints()),
1162  globalTransforms_(mesh),
1163  nPatchPoints_(coupledPatch.nPoints()),
1164  procPoints_(nPatchPoints_),
1165  meshToProcPoint_(nPatchPoints_)
1166 {
1167  calculateSharedPoints
1168  (
1169  coupledPatch.meshPointMap(),
1170  coupledPatch.meshPoints(),
1171  keepAllPoints,
1172  mergeSeparated
1173  );
1174 }
1175 
1176 
1177 // ************************************************************************* //
label nPoints() const
Number of points supporting patch faces.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
const labelList & boundaryPoints() const
Return list of boundary points, address into LOCAL point list.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
List< labelPair > labelPairList
List of labelPair.
Definition: labelPair.H:33
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
static bool less(const vector &x, const vector &y)
To compare normals.
const labelList & meshPoints() const
Return labelList of mesh points in patch.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
Definition: HashSet.H:85
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
A list of faces which address into the list of points.
const Map< label > & meshPointMap() const
Mesh point map.
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
"scheduled" : (MPI_Send, MPI_Recv)
label nPoints
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
label find(const T &val) const
Find index of the first occurrence of the value.
Definition: UList.C:173
void reset(const dimensionSet &ds)
Copy assign the exponents from the dimensionSet.
Definition: dimensionSet.C:142
const labelListList & pointFaces() const
Return point-face addressing.
void clear()
Clear the patch list and all demand-driven data.
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
PrimitivePatch< List< face >, const pointField & > primitiveFacePatch
A PrimitivePatch with List storage for the faces, const reference for the point field.
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
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:385
const polyBoundaryMesh & patches
"nonBlocking" : (MPI_Isend, MPI_Irecv)
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28