globalMeshData.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) 2015-2024 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 "globalMeshData.H"
30 #include "globalPoints.H"
31 #include "polyMesh.H"
32 #include "mapDistribute.H"
33 #include "labelIOList.H"
34 #include "mergePoints.H"
35 #include "processorPolyPatch.H"
36 #include "processorTopologyNew.H"
38 #include "ListOps.H"
39 #include "Pstream.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 defineTypeNameAndDebug(globalMeshData, 0);
46 
47 const scalar globalMeshData::matchTol_ = 1e-8;
48 
49 template<>
50 class minEqOp<labelPair>
51 {
52 public:
53  void operator()(labelPair& x, const labelPair& y) const
54  {
55  x[0] = min(x[0], y[0]);
56  x[1] = min(x[1], y[1]);
57  }
58 };
59 }
60 
61 
62 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
63 
64 void Foam::globalMeshData::initProcAddr()
65 {
66  processorPatchIndices_.resize_nocopy(mesh_.boundaryMesh().size());
67  processorPatchIndices_ = -1;
68 
69  processorPatchNeighbours_.resize_nocopy(mesh_.boundaryMesh().size());
70  processorPatchNeighbours_ = -1;
71 
72  // Construct processor patch indexing. processorPatchNeighbours_ only
73  // set if running in parallel!
74  processorPatches_.resize_nocopy(mesh_.boundaryMesh().size());
75 
76  label nNeighbours = 0;
77 
78  forAll(mesh_.boundaryMesh(), patchi)
79  {
80  if (isA<processorPolyPatch>(mesh_.boundaryMesh()[patchi]))
81  {
82  processorPatches_[nNeighbours] = patchi;
83  processorPatchIndices_[patchi] = nNeighbours++;
84  }
85  }
86  processorPatches_.resize(nNeighbours);
87 
88 
89  if (UPstream::parRun())
90  {
91  PstreamBuffers pBufs(mesh_.comm());
92 
93  // Send indices of my processor patches to my neighbours
94  for (const label patchi : processorPatches_)
95  {
96  UOPstream toNeighbour
97  (
98  refCast<const processorPolyPatch>
99  (
100  mesh_.boundaryMesh()[patchi]
101  ).neighbProcNo(),
102  pBufs
103  );
104 
105  toNeighbour << processorPatchIndices_[patchi];
106  }
107 
108  pBufs.finishedSends();
109 
110  for (const label patchi : processorPatches_)
111  {
112  UIPstream fromNeighbour
113  (
114  refCast<const processorPolyPatch>
115  (
116  mesh_.boundaryMesh()[patchi]
117  ).neighbProcNo(),
118  pBufs
119  );
120 
121  fromNeighbour >> processorPatchNeighbours_[patchi];
122  }
123  }
124 }
125 
126 
127 void Foam::globalMeshData::calcSharedPoints() const
128 {
129  if
130  (
131  nGlobalPoints_ != -1
132  || sharedPointLabelsPtr_
133  || sharedPointAddrPtr_
134  )
135  {
137  << "Shared point addressing already done" << abort(FatalError);
138  }
139 
140  // Calculate all shared points (exclude points that are only
141  // on two coupled patches). This does all the hard work.
142  const globalPoints parallelPoints(mesh_, false, true);
143 
144  // Count the number of master points
145  label nMaster = 0;
146  forAll(parallelPoints.pointPoints(), i)
147  {
148  const labelList& pPoints = parallelPoints.pointPoints()[i];
149  const labelList& transPPoints =
150  parallelPoints.transformedPointPoints()[i];
151 
152  if (pPoints.size()+transPPoints.size() > 0)
153  {
154  nMaster++;
155  }
156  }
157 
158  // Allocate global numbers
159  const globalIndex masterNumbering(nMaster);
160 
161  nGlobalPoints_ = masterNumbering.totalSize();
162 
163 
164  // Push master number to slaves
165  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
166  // 1. Fill master and slave slots
167  nMaster = 0;
168  labelList master(parallelPoints.map().constructSize(), -1);
169  forAll(parallelPoints.pointPoints(), i)
170  {
171  const labelList& pPoints = parallelPoints.pointPoints()[i];
172  const labelList& transPPoints =
173  parallelPoints.transformedPointPoints()[i];
174 
175  if (pPoints.size()+transPPoints.size() > 0)
176  {
177  master[i] = masterNumbering.toGlobal(nMaster);
178 
179  labelUIndList(master, pPoints) = master[i];
180  labelUIndList(master, transPPoints) = master[i];
181 
182  ++nMaster;
183  }
184  }
185 
186 
187  // 2. Push slave slots back to local storage on originating processor
188  // For all the four types of points:
189  // - local master : already set
190  // - local transformed slave point : the reverse transform at
191  // reverseDistribute will have copied it back to its originating local
192  // point
193  // - remote untransformed slave point : sent back to originating processor
194  // - remote transformed slave point : the reverse transform will
195  // copy it back into the remote slot which then gets sent back to
196  // originating processor
197 
198  parallelPoints.map().reverseDistribute
199  (
200  parallelPoints.map().constructSize(),
201  master
202  );
203 
204 
205  // Collect all points that are a master or refer to a master.
206  nMaster = 0;
207  forAll(parallelPoints.pointPoints(), i)
208  {
209  if (master[i] != -1)
210  {
211  nMaster++;
212  }
213  }
214 
215  sharedPointLabelsPtr_.reset(new labelList(nMaster));
216  labelList& sharedPointLabels = sharedPointLabelsPtr_();
217  sharedPointAddrPtr_.reset(new labelList(nMaster));
218  labelList& sharedPointAddr = sharedPointAddrPtr_();
219  nMaster = 0;
220 
221  forAll(parallelPoints.pointPoints(), i)
222  {
223  if (master[i] != -1)
224  {
225  // I am master or slave
226  sharedPointLabels[nMaster] = i;
227  sharedPointAddr[nMaster] = master[i];
228  nMaster++;
229  }
230  }
231 
232  if (debug)
233  {
234  Pout<< "globalMeshData : nGlobalPoints_:" << nGlobalPoints_ << nl
235  << "globalMeshData : sharedPointLabels_:"
236  << sharedPointLabelsPtr_().size() << nl
237  << "globalMeshData : sharedPointAddr_:"
238  << sharedPointAddrPtr_().size() << endl;
239  }
240 }
241 
242 
243 void Foam::globalMeshData::countSharedEdges
244 (
245  const EdgeMap<labelList>& procSharedEdges,
246  EdgeMap<label>& globalShared,
247  label& sharedEdgeI
248 )
249 {
250  // Count occurrences of procSharedEdges in global shared edges table.
251  forAllConstIters(procSharedEdges, iter)
252  {
253  const edge& e = iter.key();
254 
255  auto globalFnd = globalShared.find(e);
256 
257  if (globalFnd.good())
258  {
259  if (globalFnd() == -1)
260  {
261  // Second time occurrence of this edge.
262  // Assign proper edge label.
263  globalFnd() = sharedEdgeI++;
264  }
265  }
266  else
267  {
268  // First time occurrence of this edge. Check how many we are adding.
269  if (iter().size() == 1)
270  {
271  // Only one edge. Mark with special value.
272  globalShared.insert(e, -1);
273  }
274  else
275  {
276  // Edge used more than once (even by local shared edges alone)
277  // so allocate proper shared edge label.
278  globalShared.insert(e, sharedEdgeI++);
279  }
280  }
281  }
282 }
283 
284 
285 void Foam::globalMeshData::calcSharedEdges() const
286 {
287  // Shared edges are shared between multiple processors. By their nature both
288  // of their endpoints are shared points. (but not all edges using two shared
289  // points are shared edges! There might e.g. be an edge between two
290  // unrelated clusters of shared points)
291 
292  if
293  (
294  nGlobalEdges_ != -1
295  || sharedEdgeLabelsPtr_
296  || sharedEdgeAddrPtr_
297  )
298  {
300  << "Shared edge addressing already done" << abort(FatalError);
301  }
302 
303 
304  const labelList& sharedPtAddr = sharedPointAddr();
305  const labelList& sharedPtLabels = sharedPointLabels();
306 
307  // Since don't want to construct pointEdges for whole mesh create
308  // Map for all shared points.
309  Map<label> meshToShared(invertToMap(sharedPtLabels));
310 
311 
312  // Find edges using shared points. Store correspondence to local edge
313  // numbering. Note that multiple local edges can have the same shared
314  // points! (for cyclics or separated processor patches)
315  EdgeMap<labelList> localShared(2*sharedPtAddr.size());
316 
317  const edgeList& edges = mesh_.edges();
318 
319  forAll(edges, edgei)
320  {
321  const edge& e = edges[edgei];
322 
323  const auto e0Fnd = meshToShared.cfind(e[0]);
324 
325  if (e0Fnd.good())
326  {
327  const auto e1Fnd = meshToShared.cfind(e[1]);
328 
329  if (e1Fnd.good())
330  {
331  // Found edge which uses shared points. Probably shared.
332 
333  // Construct the edge in shared points (or rather global indices
334  // of the shared points)
335  edge sharedEdge
336  (
337  sharedPtAddr[e0Fnd.val()],
338  sharedPtAddr[e1Fnd.val()]
339  );
340 
341  // Add this edge to list of edge labels
342  localShared(sharedEdge).push_back(edgei);
343  }
344  }
345  }
346 
347 
348  // Now we have a table on every processors which gives its edges which use
349  // shared points. Send this all to the master and have it allocate
350  // global edge numbers for it. But only allocate a global edge number for
351  // edge if it is used more than once!
352  // Note that we are now sending the whole localShared to the master whereas
353  // we only need the local count (i.e. the number of times a global edge is
354  // used). But then this only gets done once so not too bothered about the
355  // extra global communication.
356 
357  EdgeMap<label> globalShared(2*nGlobalPoints());
358 
359  if (UPstream::master())
360  {
361  label sharedEdgeI = 0;
362 
363  // Merge my shared edges into the global list
364  if (debug)
365  {
366  Pout<< "globalMeshData::calcSharedEdges : Merging in from proc0 : "
367  << localShared.size() << endl;
368  }
369  countSharedEdges(localShared, globalShared, sharedEdgeI);
370 
371  // Receive data and insert
372  if (UPstream::parRun())
373  {
374  for (const int proci : UPstream::subProcs())
375  {
376  // Receive the edges using shared points from the slave.
377  EdgeMap<labelList> procSharedEdges;
378  IPstream::recv(procSharedEdges, proci);
379 
380  if (debug)
381  {
382  Pout<< "globalMeshData::calcSharedEdges : "
383  << "Merging in from proc"
384  << proci << " : " << procSharedEdges.size()
385  << endl;
386  }
387  countSharedEdges(procSharedEdges, globalShared, sharedEdgeI);
388  }
389  }
390 
391  // Now our globalShared should have some edges with -1 as edge label
392  // These were only used once so are not proper shared edges.
393  // Remove them.
394  {
395  EdgeMap<label> oldSharedEdges(std::move(globalShared));
396  globalShared.clear();
397 
398  forAllConstIters(oldSharedEdges, iter)
399  {
400  if (iter.val() != -1)
401  {
402  globalShared.insert(iter.key(), iter.val());
403  }
404  }
405 
406  if (debug)
407  {
408  Pout<< "globalMeshData::calcSharedEdges : Filtered "
409  << oldSharedEdges.size()
410  << " down to " << globalShared.size() << endl;
411  }
412  }
413  }
414  else
415  {
416  if (UPstream::parRun())
417  {
418  // send local edges to master
419  OPstream::send(localShared, UPstream::masterNo());
420  }
421  }
422 
423  // Broadcast: merged edges to all
424  Pstream::broadcast(globalShared); // == worldComm;
425 
426 
427  // Now use the global shared edges list (globalShared) to classify my local
428  // ones (localShared)
429 
430  nGlobalEdges_ = globalShared.size();
431 
432  DynamicList<label> dynSharedEdgeLabels(globalShared.size());
433  DynamicList<label> dynSharedEdgeAddr(globalShared.size());
434 
435  forAllConstIters(localShared, iter)
436  {
437  const edge& e = iter.key();
438 
439  const auto edgeFnd = globalShared.cfind(e);
440 
441  if (edgeFnd.good())
442  {
443  // My local edge is indeed a shared one. Go through all local edge
444  // labels with this point combination.
445  const labelList& edgeLabels = iter.val();
446 
447  for (const label edgei : edgeLabels)
448  {
449  // Store label of local mesh edge
450  dynSharedEdgeLabels.append(edgei);
451 
452  // Store label of shared edge
453  dynSharedEdgeAddr.append(edgeFnd());
454  }
455  }
456  }
457 
458 
459  sharedEdgeLabelsPtr_.reset
460  (
461  new labelList(std::move(dynSharedEdgeLabels))
462  );
463 
464  sharedEdgeAddrPtr_.reset
465  (
466  new labelList(std::move(dynSharedEdgeAddr))
467  );
468 
469  if (debug)
470  {
471  Pout<< "globalMeshData : nGlobalEdges_:" << nGlobalEdges_ << nl
472  << "globalMeshData : sharedEdgeLabels:"
473  << sharedEdgeLabelsPtr_().size() << nl
474  << "globalMeshData : sharedEdgeAddr:"
475  << sharedEdgeAddrPtr_().size() << endl;
476  }
477 }
478 
479 
480 void Foam::globalMeshData::calcGlobalPointSlaves() const
481 {
482  if (debug)
483  {
484  Pout<< "globalMeshData::calcGlobalPointSlaves() :"
485  << " calculating coupled master to slave point addressing."
486  << endl;
487  }
488 
489  // Calculate connected points for master points.
490  globalPoints globalData(mesh_, coupledPatch(), true, true);
491 
492  globalPointSlavesPtr_.reset
493  (
494  new labelListList
495  (
496  std::move(globalData.pointPoints())
497  )
498  );
499  globalPointTransformedSlavesPtr_.reset
500  (
501  new labelListList
502  (
503  std::move(globalData.transformedPointPoints())
504  )
505  );
506 
507  globalPointSlavesMapPtr_.reset
508  (
509  new mapDistribute
510  (
511  std::move(globalData.map())
512  )
513  );
514 }
515 
516 
517 void Foam::globalMeshData::calcPointConnectivity
518 (
519  List<labelPairList>& allPointConnectivity
520 ) const
521 {
522  const globalIndexAndTransform& transforms = globalTransforms();
523  const labelListList& slaves = globalPointSlaves();
524  const labelListList& transformedSlaves = globalPointTransformedSlaves();
525 
526 
527  // Create field with my local data
528  labelPairList myData(globalPointSlavesMap().constructSize());
529  forAll(slaves, pointi)
530  {
531  myData[pointi] = transforms.encode
532  (
534  pointi,
535  transforms.nullTransformIndex()
536  );
537  }
538  // Send to master
539  globalPointSlavesMap().distribute(myData);
540 
541 
542  // String of connected points with their transform
543  allPointConnectivity.setSize(globalPointSlavesMap().constructSize());
544  allPointConnectivity = labelPairList(0);
545 
546  // Pass1: do the master points since these also update local slaves
547  // (e.g. from local cyclics)
548  forAll(slaves, pointi)
549  {
550  // Reconstruct string of connected points
551  const labelList& pSlaves = slaves[pointi];
552  const labelList& pTransformSlaves = transformedSlaves[pointi];
553 
554  if (pSlaves.size()+pTransformSlaves.size())
555  {
556  labelPairList& pConnectivity = allPointConnectivity[pointi];
557 
558  pConnectivity.setSize(1+pSlaves.size()+pTransformSlaves.size());
559  label connI = 0;
560 
561  // Add myself
562  pConnectivity[connI++] = myData[pointi];
563  // Add untransformed points
564  forAll(pSlaves, i)
565  {
566  pConnectivity[connI++] = myData[pSlaves[i]];
567  }
568  // Add transformed points.
569  forAll(pTransformSlaves, i)
570  {
571  // Get transform from index
572  label transformI = globalPointSlavesMap().whichTransform
573  (
574  pTransformSlaves[i]
575  );
576  // Add transform to connectivity
577  const labelPair& n = myData[pTransformSlaves[i]];
578  label proci = transforms.processor(n);
579  label index = transforms.index(n);
580  pConnectivity[connI++] = transforms.encode
581  (
582  proci,
583  index,
584  transformI
585  );
586  }
587 
588  // Put back in slots
589  forAll(pSlaves, i)
590  {
591  allPointConnectivity[pSlaves[i]] = pConnectivity;
592  }
593  forAll(pTransformSlaves, i)
594  {
595  allPointConnectivity[pTransformSlaves[i]] = pConnectivity;
596  }
597  }
598  }
599 
600 
601  // Pass2: see if anything is still unset (should not be the case)
602  forAll(slaves, pointi)
603  {
604  labelPairList& pConnectivity = allPointConnectivity[pointi];
605 
606  if (pConnectivity.size() == 0)
607  {
608  pConnectivity.setSize(1, myData[pointi]);
609  }
610  }
611 
612 
613  globalPointSlavesMap().reverseDistribute
614  (
615  slaves.size(),
616  allPointConnectivity
617  );
618 }
619 
620 
621 void Foam::globalMeshData::calcGlobalPointEdges
622 (
623  labelListList& globalPointEdges,
624  List<labelPairList>& globalPointPoints
625 ) const
626 {
627  const edgeList& edges = coupledPatch().edges();
628  const labelListList& pointEdges = coupledPatch().pointEdges();
629  const globalIndex& globalEdgeNumbers = globalEdgeNumbering();
630  const labelListList& slaves = globalPointSlaves();
631  const labelListList& transformedSlaves = globalPointTransformedSlaves();
632  const globalIndexAndTransform& transforms = globalTransforms();
633 
634 
635  // Create local version
636  globalPointEdges.setSize(globalPointSlavesMap().constructSize());
637  globalPointPoints.setSize(globalPointSlavesMap().constructSize());
638  forAll(pointEdges, pointi)
639  {
640  const labelList& pEdges = pointEdges[pointi];
641  globalPointEdges[pointi] = globalEdgeNumbers.toGlobal(pEdges);
642 
643  labelPairList& globalPPoints = globalPointPoints[pointi];
644  globalPPoints.setSize(pEdges.size());
645  forAll(pEdges, i)
646  {
647  label otherPointi = edges[pEdges[i]].otherVertex(pointi);
648  globalPPoints[i] = transforms.encode
649  (
651  otherPointi,
652  transforms.nullTransformIndex()
653  );
654  }
655  }
656 
657  // Pull slave data to master. Dummy transform.
658  globalPointSlavesMap().distribute(globalPointEdges);
659  globalPointSlavesMap().distribute(globalPointPoints);
660  // Add all pointEdges
661  forAll(slaves, pointi)
662  {
663  const labelList& pSlaves = slaves[pointi];
664  const labelList& pTransformSlaves = transformedSlaves[pointi];
665 
666  label n = 0;
667  forAll(pSlaves, i)
668  {
669  n += globalPointEdges[pSlaves[i]].size();
670  }
671  forAll(pTransformSlaves, i)
672  {
673  n += globalPointEdges[pTransformSlaves[i]].size();
674  }
675 
676  // Add all the point edges of the slaves to those of the (master) point
677  {
678  labelList& globalPEdges = globalPointEdges[pointi];
679  label sz = globalPEdges.size();
680  globalPEdges.setSize(sz+n);
681  forAll(pSlaves, i)
682  {
683  const labelList& otherData = globalPointEdges[pSlaves[i]];
684  forAll(otherData, j)
685  {
686  globalPEdges[sz++] = otherData[j];
687  }
688  }
689  forAll(pTransformSlaves, i)
690  {
691  const labelList& otherData =
692  globalPointEdges[pTransformSlaves[i]];
693  forAll(otherData, j)
694  {
695  globalPEdges[sz++] = otherData[j];
696  }
697  }
698 
699  // Put back in slots
700  forAll(pSlaves, i)
701  {
702  globalPointEdges[pSlaves[i]] = globalPEdges;
703  }
704  forAll(pTransformSlaves, i)
705  {
706  globalPointEdges[pTransformSlaves[i]] = globalPEdges;
707  }
708  }
709 
710 
711  // Same for corresponding pointPoints
712  {
713  labelPairList& globalPPoints = globalPointPoints[pointi];
714  label sz = globalPPoints.size();
715  globalPPoints.setSize(sz + n);
716 
717  // Add untransformed points
718  forAll(pSlaves, i)
719  {
720  const labelPairList& otherData = globalPointPoints[pSlaves[i]];
721  forAll(otherData, j)
722  {
723  globalPPoints[sz++] = otherData[j];
724  }
725  }
726  // Add transformed points.
727  forAll(pTransformSlaves, i)
728  {
729  // Get transform from index
730  label transformI = globalPointSlavesMap().whichTransform
731  (
732  pTransformSlaves[i]
733  );
734 
735  const labelPairList& otherData =
736  globalPointPoints[pTransformSlaves[i]];
737  forAll(otherData, j)
738  {
739  // Add transform to connectivity
740  const labelPair& n = otherData[j];
741  label proci = transforms.processor(n);
742  label index = transforms.index(n);
743  globalPPoints[sz++] = transforms.encode
744  (
745  proci,
746  index,
747  transformI
748  );
749  }
750  }
751 
752  // Put back in slots
753  forAll(pSlaves, i)
754  {
755  globalPointPoints[pSlaves[i]] = globalPPoints;
756  }
757  forAll(pTransformSlaves, i)
758  {
759  globalPointPoints[pTransformSlaves[i]] = globalPPoints;
760  }
761  }
762  }
763  // Push back
764  globalPointSlavesMap().reverseDistribute
765  (
766  slaves.size(),
767  globalPointEdges
768  );
769  // Push back
770  globalPointSlavesMap().reverseDistribute
771  (
772  slaves.size(),
773  globalPointPoints
774  );
775 }
776 
777 
778 Foam::label Foam::globalMeshData::findTransform
779 (
780  const labelPairList& info,
781  const labelPair& remotePoint,
782  const label localPoint
783 ) const
784 {
785  const globalIndexAndTransform& transforms = globalTransforms();
786 
787  const label remoteProci = transforms.processor(remotePoint);
788  const label remoteIndex = transforms.index(remotePoint);
789 
790  label remoteTransformI = -1;
791  label localTransformI = -1;
792  forAll(info, i)
793  {
794  label proci = transforms.processor(info[i]);
795  label pointi = transforms.index(info[i]);
796  label transformI = transforms.transformIndex(info[i]);
797 
798  if (proci == Pstream::myProcNo() && pointi == localPoint)
799  {
800  localTransformI = transformI;
801  //Pout<< "For local :" << localPoint
802  // << " found transform:" << localTransformI
803  // << endl;
804  }
805  if (proci == remoteProci && pointi == remoteIndex)
806  {
807  remoteTransformI = transformI;
808  //Pout<< "For remote:" << remotePoint
809  // << " found transform:" << remoteTransformI
810  // << " at index:" << i
811  // << endl;
812  }
813  }
814 
815  if (remoteTransformI == -1 || localTransformI == -1)
816  {
818  << "Problem. Cannot find " << remotePoint
819  << " or " << localPoint << " "
820  << coupledPatch().localPoints()[localPoint]
821  << " in " << info
822  << endl
823  << "remoteTransformI:" << remoteTransformI << endl
824  << "localTransformI:" << localTransformI
825  << abort(FatalError);
826  }
827 
828  return transforms.subtractTransformIndex
829  (
830  remoteTransformI,
831  localTransformI
832  );
833 }
834 
835 
836 void Foam::globalMeshData::calcGlobalEdgeSlaves() const
837 {
838  if (debug)
839  {
840  Pout<< "globalMeshData::calcGlobalEdgeSlaves() :"
841  << " calculating coupled master to slave edge addressing." << endl;
842  }
843 
844  const edgeList& edges = coupledPatch().edges();
845  const globalIndex& globalEdgeNumbers = globalEdgeNumbering();
846  const globalIndexAndTransform& transforms = globalTransforms();
847 
848 
849  // The whole problem with deducting edge-connectivity from
850  // point-connectivity is that one of the endpoints might be
851  // a local master but the other endpoint might not. So we first
852  // need to make sure that all points know about connectivity and
853  // the transformations.
854 
855 
856  // 1. collect point connectivity - basically recreating globalPoints output.
857  // All points will now have a string of coupled points. The transforms are
858  // in respect to the master.
859  List<labelPairList> allPointConnectivity;
860  calcPointConnectivity(allPointConnectivity);
861 
862 
863  // 2. Get all pointEdges and pointPoints
864  // Coupled point to global coupled edges and corresponding endpoint.
865  labelListList globalPointEdges;
866  List<labelPairList> globalPointPoints;
867  calcGlobalPointEdges(globalPointEdges, globalPointPoints);
868 
869 
870  // 3. Now all points have
871  // - all the connected points with original transform
872  // - all the connected global edges
873 
874  // Now all we need to do is go through all the edges and check
875  // both endpoints. If there is a edge between the two which is
876  // produced by transforming both points in the same way it is a shared
877  // edge.
878 
879  // Collect strings of connected edges.
880  List<labelPairList> allEdgeConnectivity(edges.size());
881 
882  forAll(edges, edgeI)
883  {
884  const edge& e = edges[edgeI];
885  const labelList& pEdges0 = globalPointEdges[e[0]];
886  const labelPairList& pPoints0 = globalPointPoints[e[0]];
887  const labelList& pEdges1 = globalPointEdges[e[1]];
888  const labelPairList& pPoints1 = globalPointPoints[e[1]];
889 
890  // Most edges will be size 2
891  DynamicList<labelPair> eEdges(2);
892  // Append myself.
893  eEdges.append
894  (
895  transforms.encode
896  (
898  edgeI,
899  transforms.nullTransformIndex()
900  )
901  );
902 
903  forAll(pEdges0, i)
904  {
905  forAll(pEdges1, j)
906  {
907  if
908  (
909  pEdges0[i] == pEdges1[j]
910  && pEdges0[i] != globalEdgeNumbers.toGlobal(edgeI)
911  )
912  {
913  // Found a shared edge. Now check if the endpoints
914  // go through the same transformation.
915  // Local: e[0] remote:pPoints1[j]
916  // Local: e[1] remote:pPoints0[i]
917 
918 
919  // Find difference in transforms to go from point on remote
920  // edge (pPoints1[j]) to this point.
921 
922  label transform0 = findTransform
923  (
924  allPointConnectivity[e[0]],
925  pPoints1[j],
926  e[0]
927  );
928  label transform1 = findTransform
929  (
930  allPointConnectivity[e[1]],
931  pPoints0[i],
932  e[1]
933  );
934 
935  if (transform0 == transform1)
936  {
937  label proci = globalEdgeNumbers.whichProcID(pEdges0[i]);
938  eEdges.append
939  (
940  transforms.encode
941  (
942  proci,
943  globalEdgeNumbers.toLocal(proci, pEdges0[i]),
944  transform0
945  )
946  );
947  }
948  }
949  }
950  }
951 
952  allEdgeConnectivity[edgeI].transfer(eEdges);
953  Foam::sort
954  (
955  allEdgeConnectivity[edgeI],
957  );
958  }
959 
960  // We now have - in allEdgeConnectivity - a list of edges which are shared
961  // between multiple processors. Filter into non-transformed and transformed
962  // connections.
963 
964  globalEdgeSlavesPtr_.reset(new labelListList(edges.size()));
965  labelListList& globalEdgeSlaves = globalEdgeSlavesPtr_();
966  List<labelPairList> transformedEdges(edges.size());
967  forAll(allEdgeConnectivity, edgeI)
968  {
969  const labelPairList& edgeInfo = allEdgeConnectivity[edgeI];
970  if (edgeInfo.size() >= 2)
971  {
972  const labelPair& masterInfo = edgeInfo[0];
973 
974  // Check if master edge (= first element (since sorted)) is me.
975  if
976  (
977  (
978  transforms.processor(masterInfo)
979  == Pstream::myProcNo()
980  )
981  && (transforms.index(masterInfo) == edgeI)
982  )
983  {
984  // Sort into transformed and untransformed
985  labelList& eEdges = globalEdgeSlaves[edgeI];
986  eEdges.setSize(edgeInfo.size()-1);
987 
988  labelPairList& trafoEEdges = transformedEdges[edgeI];
989  trafoEEdges.setSize(edgeInfo.size()-1);
990 
991  label nonTransformI = 0;
992  label transformI = 0;
993 
994  for (label i = 1; i < edgeInfo.size(); i++)
995  {
996  const labelPair& info = edgeInfo[i];
997  label proci = transforms.processor(info);
998  label index = transforms.index(info);
999  label transform = transforms.transformIndex
1000  (
1001  info
1002  );
1003 
1004  if (transform == transforms.nullTransformIndex())
1005  {
1006  eEdges[nonTransformI++] = globalEdgeNumbers.toGlobal
1007  (
1008  proci,
1009  index
1010  );
1011  }
1012  else
1013  {
1014  trafoEEdges[transformI++] = info;
1015  }
1016  }
1017 
1018  eEdges.setSize(nonTransformI);
1019  trafoEEdges.setSize(transformI);
1020  }
1021  }
1022  }
1023 
1024 
1025  // Construct map
1026  globalEdgeTransformedSlavesPtr_.reset(new labelListList());
1027 
1028  List<Map<label>> compactMap(Pstream::nProcs());
1029  globalEdgeSlavesMapPtr_.reset
1030  (
1031  new mapDistribute
1032  (
1033  globalEdgeNumbers,
1034  globalEdgeSlaves,
1035 
1036  transforms,
1037  transformedEdges,
1038  globalEdgeTransformedSlavesPtr_(),
1039 
1040  compactMap
1041  )
1042  );
1043 
1044 
1045  if (debug)
1046  {
1047  Pout<< "globalMeshData::calcGlobalEdgeSlaves() :"
1048  << " coupled edges:" << edges.size()
1049  << " additional coupled edges:"
1050  << globalEdgeSlavesMapPtr_().constructSize() - edges.size()
1051  << endl;
1052  }
1053 }
1054 
1055 
1056 void Foam::globalMeshData::calcGlobalEdgeOrientation() const
1057 {
1058  if (debug)
1059  {
1060  Pout<< "globalMeshData::calcGlobalEdgeOrientation() :"
1061  << " calculating edge orientation w.r.t. master edge." << endl;
1062  }
1063 
1064  const globalIndex& globalPoints = globalPointNumbering();
1065 
1066  // 1. Determine master point
1067  labelList masterPoint;
1068  {
1069  const mapDistribute& map = globalPointSlavesMap();
1070 
1071  masterPoint.setSize(map.constructSize());
1072  masterPoint = labelMax;
1073 
1074  for (label pointi = 0; pointi < coupledPatch().nPoints(); pointi++)
1075  {
1076  masterPoint[pointi] = globalPoints.toGlobal(pointi);
1077  }
1078  syncData
1079  (
1080  masterPoint,
1081  globalPointSlaves(),
1082  globalPointTransformedSlaves(),
1083  map,
1084  minEqOp<label>()
1085  );
1086  }
1087 
1088  // Now all points should know who is master by comparing their global
1089  // pointID with the masterPointID. We now can use this information
1090  // to find the orientation of the master edge.
1091 
1092  {
1093  const mapDistribute& map = globalEdgeSlavesMap();
1094  const labelListList& slaves = globalEdgeSlaves();
1095  const labelListList& transformedSlaves = globalEdgeTransformedSlaves();
1096 
1097  // Distribute orientation of master edge (in masterPoint numbering)
1098  labelPairList masterEdgeVerts(map.constructSize());
1099  masterEdgeVerts = labelPair(labelMax, labelMax);
1100 
1101  for (label edgeI = 0; edgeI < coupledPatch().nEdges(); edgeI++)
1102  {
1103  if
1104  (
1105  (
1106  slaves[edgeI].size()
1107  + transformedSlaves[edgeI].size()
1108  )
1109  > 0
1110  )
1111  {
1112  // I am master. Fill in my masterPoint equivalent.
1113 
1114  const edge& e = coupledPatch().edges()[edgeI];
1115  masterEdgeVerts[edgeI] = labelPair
1116  (
1117  masterPoint[e[0]],
1118  masterPoint[e[1]]
1119  );
1120  }
1121  }
1122  syncData
1123  (
1124  masterEdgeVerts,
1125  slaves,
1126  transformedSlaves,
1127  map,
1128  minEqOp<labelPair>()
1129  );
1130 
1131  // Now check my edges on how they relate to the master's edgeVerts
1132  globalEdgeOrientationPtr_.reset
1133  (
1134  new bitSet(coupledPatch().nEdges())
1135  );
1136  bitSet& globalEdgeOrientation = globalEdgeOrientationPtr_();
1137 
1138  forAll(coupledPatch().edges(), edgeI)
1139  {
1140  // Test that edge is not single edge on cyclic baffle
1141  if (masterEdgeVerts[edgeI] != labelPair(labelMax, labelMax))
1142  {
1143  const edge& e = coupledPatch().edges()[edgeI];
1144  const labelPair masterE
1145  (
1146  masterPoint[e[0]],
1147  masterPoint[e[1]]
1148  );
1149 
1150  const int stat = labelPair::compare
1151  (
1152  masterE,
1153  masterEdgeVerts[edgeI]
1154  );
1155  if (stat == 0)
1156  {
1158  << "problem : my edge:" << e
1159  << " in master points:" << masterE
1160  << " v.s. masterEdgeVerts:" << masterEdgeVerts[edgeI]
1161  << exit(FatalError);
1162  }
1163  else
1164  {
1165  globalEdgeOrientation.set(edgeI, (stat == 1));
1166  }
1167  }
1168  else
1169  {
1170  globalEdgeOrientation.set(edgeI, true);
1171  }
1172  }
1173  }
1174 
1175  if (debug)
1176  {
1177  Pout<< "globalMeshData::calcGlobalEdgeOrientation() :"
1178  << " finished calculating edge orientation."
1179  << endl;
1180  }
1181 }
1182 
1183 
1184 void Foam::globalMeshData::calcPointBoundaryFaces
1185 (
1186  labelListList& pointBoundaryFaces
1187 ) const
1188 {
1189  const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
1190  const Map<label>& meshPointMap = coupledPatch().meshPointMap();
1191 
1192  // 1. Count
1193 
1194  labelList nPointFaces(coupledPatch().nPoints(), Zero);
1195 
1196  for (const polyPatch& pp : bMesh)
1197  {
1198  if (!pp.coupled())
1199  {
1200  for (const face& f : pp)
1201  {
1202  forAll(f, fp)
1203  {
1204  const auto iter = meshPointMap.cfind(f[fp]);
1205  if (iter.good())
1206  {
1207  nPointFaces[iter.val()]++;
1208  }
1209  }
1210  }
1211  }
1212  }
1213 
1214 
1215  // 2. Size
1216 
1217  pointBoundaryFaces.setSize(coupledPatch().nPoints());
1218  forAll(nPointFaces, pointi)
1219  {
1220  pointBoundaryFaces[pointi].setSize(nPointFaces[pointi]);
1221  }
1222  nPointFaces = 0;
1223 
1224 
1225  // 3. Fill
1226 
1227  forAll(bMesh, patchi)
1228  {
1229  const polyPatch& pp = bMesh[patchi];
1230 
1231  if (!pp.coupled())
1232  {
1233  forAll(pp, i)
1234  {
1235  const face& f = pp[i];
1236  forAll(f, fp)
1237  {
1238  const auto iter = meshPointMap.cfind(f[fp]);
1239 
1240  if (iter.good())
1241  {
1242  label bFacei =
1243  pp.start() + i - mesh_.nInternalFaces();
1244  pointBoundaryFaces[iter()][nPointFaces[iter()]++] =
1245  bFacei;
1246  }
1247  }
1248  }
1249  }
1250  }
1251 }
1252 
1253 
1254 void Foam::globalMeshData::calcGlobalPointBoundaryFaces() const
1255 {
1256  if (debug)
1257  {
1258  Pout<< "globalMeshData::calcGlobalPointBoundaryFaces() :"
1259  << " calculating coupled point to boundary face addressing."
1260  << endl;
1261  }
1262 
1263  const label myProci = UPstream::myProcNo();
1264 
1265  // Construct local point to (uncoupled)boundaryfaces.
1266  labelListList pointBoundaryFaces;
1267  calcPointBoundaryFaces(pointBoundaryFaces);
1268 
1269 
1270  // Global indices for boundary faces
1271  globalBoundaryFaceNumberingPtr_.reset
1272  (
1273  new globalIndex(mesh_.nBoundaryFaces())
1274  );
1275  const auto& globalIndices = *globalBoundaryFaceNumberingPtr_;
1276 
1277 
1278  // Convert local boundary faces to global numbering
1279  globalPointBoundaryFacesPtr_.reset
1280  (
1281  new labelListList(globalPointSlavesMap().constructSize())
1282  );
1283  auto& globalPointBoundaryFaces = *globalPointBoundaryFacesPtr_;
1284 
1285  forAll(pointBoundaryFaces, pointi)
1286  {
1287  globalPointBoundaryFaces[pointi] = globalIndices.toGlobal
1288  (
1289  myProci,
1290  pointBoundaryFaces[pointi]
1291  );
1292  }
1293 
1294 
1295  // Pull slave pointBoundaryFaces to master
1296  globalPointSlavesMap().distribute
1297  (
1298  globalPointBoundaryFaces,
1299  true // put data on transformed points into correct slots
1300  );
1301 
1302 
1303  // Merge slave labels into master globalPointBoundaryFaces.
1304  // Split into untransformed and transformed values.
1305  const labelListList& pointSlaves = globalPointSlaves();
1306  const labelListList& pointTransformSlaves =
1307  globalPointTransformedSlaves();
1308  const globalIndexAndTransform& transforms = globalTransforms();
1309 
1310 
1311  // Any faces coming in through transformation
1312  List<labelPairList> transformedFaces(pointSlaves.size());
1313 
1314 
1315  forAll(pointSlaves, pointi)
1316  {
1317  const labelList& slaves = pointSlaves[pointi];
1318  const labelList& transformedSlaves = pointTransformSlaves[pointi];
1319 
1320  if (slaves.size() > 0)
1321  {
1322  labelList& myBFaces = globalPointBoundaryFaces[pointi];
1323  label sz = myBFaces.size();
1324 
1325  // Count
1326  label n = 0;
1327  forAll(slaves, i)
1328  {
1329  n += globalPointBoundaryFaces[slaves[i]].size();
1330  }
1331  // Fill
1332  myBFaces.setSize(sz+n);
1333  n = sz;
1334  forAll(slaves, i)
1335  {
1336  const labelList& slaveBFaces =
1337  globalPointBoundaryFaces[slaves[i]];
1338 
1339  // Add all slaveBFaces. Note that need to check for
1340  // uniqueness only in case of cyclics.
1341 
1342  for (const label slave : slaveBFaces)
1343  {
1344  if (!SubList<label>(myBFaces, sz).found(slave))
1345  {
1346  myBFaces[n++] = slave;
1347  }
1348  }
1349  }
1350  myBFaces.setSize(n);
1351  }
1352 
1353 
1354  if (transformedSlaves.size() > 0)
1355  {
1356  const labelList& untrafoFaces = globalPointBoundaryFaces[pointi];
1357 
1358  labelPairList& myBFaces = transformedFaces[pointi];
1359  label sz = myBFaces.size();
1360 
1361  // Count
1362  label n = 0;
1363  forAll(transformedSlaves, i)
1364  {
1365  n += globalPointBoundaryFaces[transformedSlaves[i]].size();
1366  }
1367  // Fill
1368  myBFaces.setSize(sz+n);
1369  n = sz;
1370  forAll(transformedSlaves, i)
1371  {
1372  label transformI = globalPointSlavesMap().whichTransform
1373  (
1374  transformedSlaves[i]
1375  );
1376 
1377  const labelList& slaveBFaces =
1378  globalPointBoundaryFaces[transformedSlaves[i]];
1379 
1380  for (const label slave : slaveBFaces)
1381  {
1382  // Check that same face not already present untransformed
1383  if (!untrafoFaces.found(slave))
1384  {
1385  label proci = globalIndices.whichProcID(slave);
1386  label facei = globalIndices.toLocal(proci, slave);
1387 
1388  myBFaces[n++] = transforms.encode
1389  (
1390  proci,
1391  facei,
1392  transformI
1393  );
1394  }
1395  }
1396  }
1397  myBFaces.setSize(n);
1398  }
1399 
1400 
1401  if (slaves.size() + transformedSlaves.size() == 0)
1402  {
1403  globalPointBoundaryFaces[pointi].clear();
1404  }
1405  }
1406 
1407  // Construct a map to get the face data directly
1408  List<Map<label>> compactMap(Pstream::nProcs());
1409 
1410  globalPointTransformedBoundaryFacesPtr_.reset
1411  (
1412  new labelListList(transformedFaces.size())
1413  );
1414 
1415  globalPointBoundaryFacesMapPtr_.reset
1416  (
1417  new mapDistribute
1418  (
1419  globalIndices,
1420  globalPointBoundaryFaces,
1421 
1422  transforms,
1423  transformedFaces,
1424  globalPointTransformedBoundaryFacesPtr_(),
1425 
1426  compactMap
1427  )
1428  );
1429  globalPointBoundaryFaces.setSize(coupledPatch().nPoints());
1430  globalPointTransformedBoundaryFacesPtr_().setSize(coupledPatch().nPoints());
1431 
1432  if (debug)
1433  {
1434  Pout<< "globalMeshData::calcGlobalPointBoundaryFaces() :"
1435  << " coupled points:" << coupledPatch().nPoints()
1436  << " local boundary faces:" << globalIndices.localSize()
1437  << " additional coupled faces:"
1438  << globalPointBoundaryFacesMapPtr_().constructSize()
1439  - globalIndices.localSize()
1440  << endl;
1441  }
1442 }
1443 
1444 
1445 void Foam::globalMeshData::calcGlobalPointBoundaryCells() const
1446 {
1447  if (debug)
1448  {
1449  Pout<< "globalMeshData::calcGlobalPointBoundaryCells() :"
1450  << " calculating coupled point to boundary cell addressing."
1451  << endl;
1452  }
1453 
1454  const label myProci = UPstream::myProcNo();
1455 
1456  // Create map of boundary cells and point-cell addressing
1457  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1458 
1459  label bCelli = 0;
1460  Map<label> meshCellMap(4*coupledPatch().nPoints());
1461  DynamicList<label> cellMap(meshCellMap.size());
1462 
1463  // Create addressing for point to boundary cells (local)
1464  labelListList pointBoundaryCells(coupledPatch().nPoints());
1465 
1466  forAll(coupledPatch().meshPoints(), pointi)
1467  {
1468  label meshPointi = coupledPatch().meshPoints()[pointi];
1469  const labelList& pCells = mesh_.pointCells(meshPointi);
1470 
1471  labelList& bCells = pointBoundaryCells[pointi];
1472  bCells.setSize(pCells.size());
1473 
1474  forAll(pCells, i)
1475  {
1476  const label celli = pCells[i];
1477  const auto fnd = meshCellMap.cfind(celli);
1478 
1479  if (fnd.good())
1480  {
1481  bCells[i] = fnd();
1482  }
1483  else
1484  {
1485  meshCellMap.insert(celli, bCelli);
1486  cellMap.append(celli);
1487  bCells[i] = bCelli;
1488  bCelli++;
1489  }
1490  }
1491  }
1492 
1493 
1494  boundaryCellsPtr_.reset(new labelList(std::move(cellMap)));
1495  labelList& boundaryCells = boundaryCellsPtr_();
1496 
1497 
1498  // Convert point-cells to global (boundary)cell numbers
1499  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1500 
1501  globalBoundaryCellNumberingPtr_.reset
1502  (
1503  new globalIndex(boundaryCells.size())
1504  );
1505  const auto& globalIndices = *globalBoundaryCellNumberingPtr_;
1506 
1507 
1508  // Convert local boundary cells to global numbering
1509  globalPointBoundaryCellsPtr_.reset
1510  (
1511  new labelListList(globalPointSlavesMap().constructSize())
1512  );
1513  auto& globalPointBoundaryCells = *globalPointBoundaryCellsPtr_;
1514 
1515  forAll(pointBoundaryCells, pointi)
1516  {
1517  globalPointBoundaryCells[pointi] = globalIndices.toGlobal
1518  (
1519  myProci,
1520  pointBoundaryCells[pointi]
1521  );
1522  }
1523 
1524 
1525  // Pull slave pointBoundaryCells to master
1526  globalPointSlavesMap().distribute
1527  (
1528  globalPointBoundaryCells,
1529  true // put data on transformed points into correct slots
1530  );
1531 
1532 
1533  // Merge slave labels into master globalPointBoundaryCells
1534  const labelListList& pointSlaves = globalPointSlaves();
1535  const labelListList& pointTransformSlaves =
1536  globalPointTransformedSlaves();
1537  const globalIndexAndTransform& transforms = globalTransforms();
1538 
1539  List<labelPairList> transformedCells(pointSlaves.size());
1540 
1541 
1542  forAll(pointSlaves, pointi)
1543  {
1544  const labelList& slaves = pointSlaves[pointi];
1545  const labelList& transformedSlaves = pointTransformSlaves[pointi];
1546 
1547  if (slaves.size() > 0)
1548  {
1549  labelList& myBCells = globalPointBoundaryCells[pointi];
1550  label sz = myBCells.size();
1551 
1552  // Count
1553  label n = 0;
1554  forAll(slaves, i)
1555  {
1556  n += globalPointBoundaryCells[slaves[i]].size();
1557  }
1558  // Fill
1559  myBCells.setSize(sz+n);
1560  n = sz;
1561  forAll(slaves, i)
1562  {
1563  const labelList& slaveBCells =
1564  globalPointBoundaryCells[slaves[i]];
1565 
1566  // Add all slaveBCells. Note that need to check for
1567  // uniqueness only in case of cyclics.
1568 
1569  for (const label slave : slaveBCells)
1570  {
1571  if (!SubList<label>(myBCells, sz).found(slave))
1572  {
1573  myBCells[n++] = slave;
1574  }
1575  }
1576  }
1577  myBCells.setSize(n);
1578  }
1579 
1580 
1581  if (transformedSlaves.size() > 0)
1582  {
1583  const labelList& untrafoCells = globalPointBoundaryCells[pointi];
1584 
1585  labelPairList& myBCells = transformedCells[pointi];
1586  label sz = myBCells.size();
1587 
1588  // Count
1589  label n = 0;
1590  forAll(transformedSlaves, i)
1591  {
1592  n += globalPointBoundaryCells[transformedSlaves[i]].size();
1593  }
1594  // Fill
1595  myBCells.setSize(sz+n);
1596  n = sz;
1597  forAll(transformedSlaves, i)
1598  {
1599  label transformI = globalPointSlavesMap().whichTransform
1600  (
1601  transformedSlaves[i]
1602  );
1603 
1604  const labelList& slaveBCells =
1605  globalPointBoundaryCells[transformedSlaves[i]];
1606 
1607  for (const label slave : slaveBCells)
1608  {
1609  // Check that same cell not already present untransformed
1610  if (!untrafoCells.found(slave))
1611  {
1612  label proci = globalIndices.whichProcID(slave);
1613  label celli = globalIndices.toLocal(proci, slave);
1614  myBCells[n++] = transforms.encode
1615  (
1616  proci,
1617  celli,
1618  transformI
1619  );
1620  }
1621  }
1622  }
1623  myBCells.setSize(n);
1624  }
1625 
1626  if (slaves.size() + transformedSlaves.size() == 0)
1627  {
1628  globalPointBoundaryCells[pointi].clear();
1629  }
1630  }
1631 
1632  // Construct a map to get the cell data directly
1633  List<Map<label>> compactMap(Pstream::nProcs());
1634 
1635  globalPointTransformedBoundaryCellsPtr_.reset
1636  (
1637  new labelListList(transformedCells.size())
1638  );
1639 
1640  globalPointBoundaryCellsMapPtr_.reset
1641  (
1642  new mapDistribute
1643  (
1644  globalIndices,
1645  globalPointBoundaryCells,
1646 
1647  transforms,
1648  transformedCells,
1649  globalPointTransformedBoundaryCellsPtr_(),
1650 
1651  compactMap
1652  )
1653  );
1654  globalPointBoundaryCells.setSize(coupledPatch().nPoints());
1655  globalPointTransformedBoundaryCellsPtr_().setSize(coupledPatch().nPoints());
1656 
1657  if (debug)
1658  {
1659  Pout<< "globalMeshData::calcGlobalPointBoundaryCells() :"
1660  << " coupled points:" << coupledPatch().nPoints()
1661  << " local boundary cells:" << globalIndices.localSize()
1662  << " additional coupled cells:"
1663  << globalPointBoundaryCellsMapPtr_().constructSize()
1664  - globalIndices.localSize()
1665  << endl;
1666  }
1667 }
1668 
1669 
1670 void Foam::globalMeshData::calcGlobalCoPointSlaves() const
1671 {
1672  if (debug)
1673  {
1674  Pout<< "globalMeshData::calcGlobalCoPointSlaves() :"
1675  << " calculating coupled master to collocated"
1676  << " slave point addressing." << endl;
1677  }
1678 
1679  // Calculate connected points for master points.
1680  globalPoints globalData(mesh_, coupledPatch(), true, false);
1681 
1682  globalCoPointSlavesPtr_.reset
1683  (
1684  new labelListList
1685  (
1686  std::move(globalData.pointPoints())
1687  )
1688  );
1689  globalCoPointSlavesMapPtr_.reset
1690  (
1691  new mapDistribute
1692  (
1693  std::move(globalData.map())
1694  )
1695  );
1696 
1697  if (debug)
1698  {
1699  Pout<< "globalMeshData::calcGlobalCoPointSlaves() :"
1700  << " finished calculating coupled master to collocated"
1701  << " slave point addressing." << endl;
1702  }
1703 }
1704 
1705 
1706 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
1707 
1709 :
1710  mesh_(mesh),
1711  globalMeshPointAddr_(),
1712  globalMeshFaceAddr_(),
1713  globalMeshCellAddr_(),
1714  processorTopology_
1715  (
1716  processorTopology::New<processorPolyPatch>
1717  (
1718  mesh.boundaryMesh(),
1719  mesh_.comm()
1720  )
1721  ),
1722  processorPatches_(),
1723  processorPatchIndices_(),
1724  processorPatchNeighbours_(),
1725  nGlobalPoints_(-1),
1726  sharedPointLabelsPtr_(nullptr),
1727  sharedPointAddrPtr_(nullptr),
1728  sharedPointGlobalLabelsPtr_(nullptr),
1729  nGlobalEdges_(-1),
1730  sharedEdgeLabelsPtr_(nullptr),
1731  sharedEdgeAddrPtr_(nullptr)
1732 {
1734 }
1735 
1736 
1737 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1738 
1739 // A non-default destructor since we had incomplete types in the header
1741 {}
1742 
1743 
1745 {
1746  // Point
1747  nGlobalPoints_ = -1;
1748  sharedPointLabelsPtr_.clear();
1749  sharedPointAddrPtr_.clear();
1750  sharedPointGlobalLabelsPtr_.clear();
1751 
1752  // Edge
1753  nGlobalEdges_ = -1;
1754  sharedEdgeLabelsPtr_.clear();
1755  sharedEdgeAddrPtr_.clear();
1756 
1757  // Coupled patch
1758  coupledPatchPtr_.clear();
1759  coupledPatchMeshEdgesPtr_.clear();
1760  coupledPatchMeshEdgeMapPtr_.clear();
1761  globalTransformsPtr_.clear();
1762 
1763  // Point
1764  globalPointNumberingPtr_.clear();
1765  globalPointSlavesPtr_.clear();
1766  globalPointTransformedSlavesPtr_.clear();
1767  globalPointSlavesMapPtr_.clear();
1768 
1769  // Edge
1770  globalEdgeNumberingPtr_.clear();
1771  globalEdgeSlavesPtr_.clear();
1772  globalEdgeTransformedSlavesPtr_.clear();
1773  globalEdgeOrientationPtr_.clear();
1774  globalEdgeSlavesMapPtr_.clear();
1775 
1776  // Face
1777  globalBoundaryFaceNumberingPtr_.clear();
1778  globalPointBoundaryFacesPtr_.clear();
1779  globalPointTransformedBoundaryFacesPtr_.clear();
1780  globalPointBoundaryFacesMapPtr_.clear();
1781 
1782  // Cell
1783  boundaryCellsPtr_.clear();
1784  globalBoundaryCellNumberingPtr_.clear();
1785  globalPointBoundaryCellsPtr_.clear();
1786  globalPointTransformedBoundaryCellsPtr_.clear();
1787  globalPointBoundaryCellsMapPtr_.clear();
1788 
1789  // Other: collocated points
1790  globalCoPointSlavesPtr_.clear();
1791  globalCoPointSlavesMapPtr_.clear();
1792 }
1793 
1794 
1795 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1796 
1798 {
1799  if (!sharedPointGlobalLabelsPtr_)
1800  {
1801  sharedPointGlobalLabelsPtr_.reset
1802  (
1803  new labelList(sharedPointLabels().size())
1804  );
1805  labelList& sharedPointGlobalLabels = sharedPointGlobalLabelsPtr_();
1806 
1807  IOobject addrHeader
1808  (
1809  "pointProcAddressing",
1810  mesh_.facesInstance()/polyMesh::meshSubDir,
1811  mesh_,
1813  );
1814 
1815  if (addrHeader.typeHeaderOk<labelIOList>(true))
1816  {
1817  // There is a pointProcAddressing file so use it to get labels
1818  // on the original mesh
1819  Pout<< "globalMeshData::sharedPointGlobalLabels : "
1820  << "Reading pointProcAddressing" << endl;
1821 
1822  labelIOList pointProcAddressing(addrHeader);
1823 
1824  const labelList& pointLabels = sharedPointLabels();
1825 
1826  forAll(pointLabels, i)
1827  {
1828  // Get my mesh point
1829  label pointi = pointLabels[i];
1830 
1831  // Map to mesh point of original mesh
1832  sharedPointGlobalLabels[i] = pointProcAddressing[pointi];
1833  }
1834  }
1835  else
1836  {
1837  Pout<< "globalMeshData::sharedPointGlobalLabels :"
1838  << " Setting pointProcAddressing to -1" << endl;
1839 
1840  sharedPointGlobalLabels = -1;
1841  }
1842  }
1843 
1844  return *sharedPointGlobalLabelsPtr_;
1845 }
1846 
1847 
1849 {
1850  // Get all processors to send their shared points to master.
1851  // (not very efficient)
1852 
1853  pointField sharedPoints(nGlobalPoints());
1854  const labelList& pointAddr = sharedPointAddr();
1855  const labelList& pointLabels = sharedPointLabels();
1856 
1857  if (UPstream::master())
1858  {
1859  // Master:
1860  // insert my own data first
1861  forAll(pointLabels, i)
1862  {
1863  label sharedPointi = pointAddr[i];
1864 
1865  sharedPoints[sharedPointi] = mesh_.points()[pointLabels[i]];
1866  }
1867 
1868  // Receive data and insert
1869  for (const int proci : UPstream::subProcs())
1870  {
1871  IPstream fromProc(UPstream::commsTypes::scheduled, proci);
1872 
1873  labelList nbrSharedPointAddr;
1874  pointField nbrSharedPoints;
1875  fromProc >> nbrSharedPointAddr >> nbrSharedPoints;
1876 
1877  forAll(nbrSharedPointAddr, i)
1878  {
1879  label sharedPointi = nbrSharedPointAddr[i];
1880 
1881  sharedPoints[sharedPointi] = nbrSharedPoints[i];
1882  }
1883  }
1884  }
1885  else
1886  {
1887  if (UPstream::parRun())
1888  {
1889  // Send address and points
1890  OPstream toMaster
1891  (
1894  );
1895  toMaster
1896  << pointAddr
1897  << pointField(mesh_.points(), pointLabels);
1898  }
1899  }
1900 
1901  // Broadcast: sharedPoints to all
1902  Pstream::broadcast(sharedPoints); // == worldComm
1903 
1904 
1905  return sharedPoints;
1906 }
1907 
1908 
1910 {
1911  // Get coords of my shared points
1912  pointField sharedPoints(mesh_.points(), sharedPointLabels());
1913 
1914  // Append from all processors, globally consistent
1916 
1917  // Merge tolerance
1918  scalar tolDim = matchTol_ * mesh_.bounds().mag();
1919 
1920  labelList pointMap;
1922  (
1923  sharedPoints, // coordinates to merge
1924  tolDim, // tolerance
1925  false, // verbosity
1926  pointMap
1927  );
1928 
1929  return sharedPoints;
1930 }
1931 
1932 
1933 Foam::label Foam::globalMeshData::nGlobalPoints() const
1934 {
1935  if (nGlobalPoints_ == -1)
1936  {
1937  calcSharedPoints();
1938  }
1939  return nGlobalPoints_;
1940 }
1941 
1942 
1944 {
1945  if (!sharedPointLabelsPtr_)
1946  {
1947  calcSharedPoints();
1948  }
1949  return *sharedPointLabelsPtr_;
1950 }
1951 
1952 
1954 {
1955  if (!sharedPointAddrPtr_)
1956  {
1957  calcSharedPoints();
1958  }
1959  return *sharedPointAddrPtr_;
1960 }
1961 
1962 
1963 Foam::label Foam::globalMeshData::nGlobalEdges() const
1964 {
1965  if (nGlobalEdges_ == -1)
1966  {
1967  calcSharedEdges();
1968  }
1969  return nGlobalEdges_;
1970 }
1971 
1972 
1974 {
1975  if (!sharedEdgeLabelsPtr_)
1976  {
1977  calcSharedEdges();
1978  }
1979  return *sharedEdgeLabelsPtr_;
1980 }
1981 
1982 
1984 {
1985  if (!sharedEdgeAddrPtr_)
1986  {
1987  calcSharedEdges();
1988  }
1989  return *sharedEdgeAddrPtr_;
1990 }
1991 
1992 
1994 {
1995  if (!coupledPatchPtr_)
1996  {
1997  const polyBoundaryMesh& bMesh = mesh_.boundaryMesh();
1998 
1999  label nCoupled = 0;
2000 
2001  forAll(bMesh, patchi)
2002  {
2003  const polyPatch& pp = bMesh[patchi];
2004 
2005  if (pp.coupled())
2006  {
2007  nCoupled += pp.size();
2008  }
2009  }
2010  labelList coupledFaces(nCoupled);
2011  nCoupled = 0;
2012 
2013  forAll(bMesh, patchi)
2014  {
2015  const polyPatch& pp = bMesh[patchi];
2016 
2017  if (pp.coupled())
2018  {
2019  label facei = pp.start();
2020 
2021  forAll(pp, i)
2022  {
2023  coupledFaces[nCoupled++] = facei++;
2024  }
2025  }
2026  }
2027 
2028  coupledPatchPtr_.reset
2029  (
2031  (
2032  IndirectList<face>
2033  (
2034  mesh_.faces(),
2035  coupledFaces
2036  ),
2037  mesh_.points()
2038  )
2039  );
2040 
2041  if (debug)
2042  {
2043  Pout<< "globalMeshData::coupledPatch() :"
2044  << " constructed coupled faces patch:"
2045  << " faces:" << coupledPatchPtr_().size()
2046  << " points:" << coupledPatchPtr_().nPoints()
2047  << endl;
2048  }
2049  }
2050  return *coupledPatchPtr_;
2051 }
2052 
2053 
2055 {
2056  if (!coupledPatchMeshEdgesPtr_)
2057  {
2058  coupledPatchMeshEdgesPtr_.reset
2059  (
2060  new labelList
2061  (
2062  coupledPatch().meshEdges
2063  (
2064  mesh_.edges(),
2065  mesh_.pointEdges()
2066  )
2067  )
2068  );
2069  }
2070  return *coupledPatchMeshEdgesPtr_;
2071 }
2072 
2073 
2075 const
2076 {
2077  if (!coupledPatchMeshEdgeMapPtr_)
2078  {
2079  const labelList& me = coupledPatchMeshEdges();
2080 
2081  coupledPatchMeshEdgeMapPtr_.reset
2082  (
2084  );
2085  }
2086  return *coupledPatchMeshEdgeMapPtr_;
2087 }
2088 
2089 
2091 {
2092  if (!globalPointNumberingPtr_)
2093  {
2094  globalPointNumberingPtr_.reset
2095  (
2096  new globalIndex(coupledPatch().nPoints())
2097  );
2098  }
2099  return *globalPointNumberingPtr_;
2100 }
2101 
2102 
2105 {
2106  if (!globalTransformsPtr_)
2107  {
2108  globalTransformsPtr_.reset(new globalIndexAndTransform(mesh_));
2109  }
2110  return *globalTransformsPtr_;
2111 }
2112 
2113 
2115 {
2116  if (!globalPointSlavesPtr_)
2117  {
2118  calcGlobalPointSlaves();
2119  }
2120  return *globalPointSlavesPtr_;
2121 }
2122 
2123 
2125 const
2126 {
2127  if (!globalPointTransformedSlavesPtr_)
2128  {
2129  calcGlobalPointSlaves();
2130  }
2131  return *globalPointTransformedSlavesPtr_;
2132 }
2133 
2134 
2136 {
2137  if (!globalPointSlavesMapPtr_)
2138  {
2139  calcGlobalPointSlaves();
2140  }
2141  return *globalPointSlavesMapPtr_;
2142 }
2143 
2144 
2146 {
2147  if (!globalEdgeNumberingPtr_)
2148  {
2149  globalEdgeNumberingPtr_.reset
2150  (
2151  new globalIndex(coupledPatch().nEdges())
2152  );
2153  }
2154  return *globalEdgeNumberingPtr_;
2155 }
2156 
2157 
2159 {
2160  if (!globalEdgeSlavesPtr_)
2161  {
2162  calcGlobalEdgeSlaves();
2163  }
2164  return *globalEdgeSlavesPtr_;
2165 }
2166 
2167 
2169 const
2170 {
2171  if (!globalEdgeTransformedSlavesPtr_)
2172  {
2173  calcGlobalEdgeSlaves();
2174  }
2175  return *globalEdgeTransformedSlavesPtr_;
2176 }
2177 
2178 
2180 {
2181  if (!globalEdgeOrientationPtr_)
2182  {
2183  calcGlobalEdgeOrientation();
2184  }
2185  return *globalEdgeOrientationPtr_;
2186 }
2187 
2188 
2190 {
2191  if (!globalEdgeSlavesMapPtr_)
2192  {
2193  calcGlobalEdgeSlaves();
2194  }
2195  return *globalEdgeSlavesMapPtr_;
2196 }
2197 
2198 
2200 const
2201 {
2202  if (!globalBoundaryFaceNumberingPtr_)
2203  {
2204  calcGlobalPointBoundaryFaces();
2205  }
2206  return *globalBoundaryFaceNumberingPtr_;
2207 }
2208 
2209 
2211 const
2212 {
2213  if (!globalPointBoundaryFacesPtr_)
2214  {
2215  calcGlobalPointBoundaryFaces();
2216  }
2217  return *globalPointBoundaryFacesPtr_;
2218 }
2219 
2220 
2221 const Foam::labelListList&
2223 {
2224  if (!globalPointTransformedBoundaryFacesPtr_)
2225  {
2226  calcGlobalPointBoundaryFaces();
2227  }
2228  return *globalPointTransformedBoundaryFacesPtr_;
2229 }
2230 
2231 
2233 const
2234 {
2235  if (!globalPointBoundaryFacesMapPtr_)
2236  {
2237  calcGlobalPointBoundaryFaces();
2238  }
2239  return *globalPointBoundaryFacesMapPtr_;
2240 }
2241 
2242 
2244 {
2245  if (!boundaryCellsPtr_)
2246  {
2247  calcGlobalPointBoundaryCells();
2248  }
2249  return *boundaryCellsPtr_;
2250 }
2251 
2252 
2254 const
2255 {
2256  if (!globalBoundaryCellNumberingPtr_)
2257  {
2258  calcGlobalPointBoundaryCells();
2259  }
2260  return *globalBoundaryCellNumberingPtr_;
2261 }
2262 
2263 
2265 const
2266 {
2267  if (!globalPointBoundaryCellsPtr_)
2268  {
2269  calcGlobalPointBoundaryCells();
2270  }
2271  return *globalPointBoundaryCellsPtr_;
2272 }
2273 
2274 
2275 const Foam::labelListList&
2277 {
2278  if (!globalPointTransformedBoundaryCellsPtr_)
2279  {
2280  calcGlobalPointBoundaryCells();
2281  }
2282  return *globalPointTransformedBoundaryCellsPtr_;
2283 }
2284 
2285 
2287 const
2288 {
2289  if (!globalPointBoundaryCellsMapPtr_)
2290  {
2291  calcGlobalPointBoundaryCells();
2292  }
2293  return *globalPointBoundaryCellsMapPtr_;
2294 }
2295 
2296 
2298 {
2299  if (!globalCoPointSlavesPtr_)
2300  {
2301  calcGlobalCoPointSlaves();
2302  }
2303  return *globalCoPointSlavesPtr_;
2304 }
2305 
2306 
2308 {
2309  if (!globalCoPointSlavesMapPtr_)
2310  {
2311  calcGlobalCoPointSlaves();
2312  }
2313  return *globalCoPointSlavesMapPtr_;
2314 }
2315 
2316 
2318 (
2319  labelList& pointToGlobal,
2320  labelList& uniquePoints
2321 ) const
2322 {
2323  const indirectPrimitivePatch& cpp = coupledPatch();
2324  const globalIndex& globalCoupledPoints = globalPointNumbering();
2325  // Use collocated only
2326  const labelListList& pointSlaves = globalCoPointSlaves();
2327  const mapDistribute& pointSlavesMap = globalCoPointSlavesMap();
2328 
2329 
2330  // Points are either
2331  // - master with slaves
2332  // - slave with a master
2333  // - other (since e.g. non-collocated cyclics not connected)
2334 
2335  labelList masterGlobalPoint(cpp.nPoints(), -1);
2336  forAll(masterGlobalPoint, pointi)
2337  {
2338  const labelList& slavePoints = pointSlaves[pointi];
2339  if (slavePoints.size() > 0)
2340  {
2341  masterGlobalPoint[pointi] = globalCoupledPoints.toGlobal(pointi);
2342  }
2343  }
2344 
2345  // Sync by taking max
2346  syncData
2347  (
2348  masterGlobalPoint,
2349  pointSlaves,
2350  labelListList(0), // no transforms
2351  pointSlavesMap,
2352  maxEqOp<label>()
2353  );
2354 
2355 
2356  // 1. Count number of masters on my processor.
2357  label nMaster = 0;
2358  bitSet isMaster(mesh_.nPoints(), true);
2359  forAll(pointSlaves, pointi)
2360  {
2361  if (masterGlobalPoint[pointi] == -1)
2362  {
2363  // unconnected point (e.g. from separated cyclic)
2364  nMaster++;
2365  }
2366  else if
2367  (
2368  masterGlobalPoint[pointi]
2369  == globalCoupledPoints.toGlobal(pointi)
2370  )
2371  {
2372  // connected master
2373  nMaster++;
2374  }
2375  else
2376  {
2377  // connected slave point
2378  isMaster.unset(cpp.meshPoints()[pointi]);
2379  }
2380  }
2381 
2382  label myUniquePoints = mesh_.nPoints() - cpp.nPoints() + nMaster;
2383 
2384  //Pout<< "Points :" << nl
2385  // << " mesh : " << mesh_.nPoints() << nl
2386  // << " of which coupled : " << cpp.nPoints() << nl
2387  // << " of which master : " << nMaster << nl
2388  // << endl;
2389 
2390 
2391  // 2. Create global indexing for unique points.
2392  autoPtr<globalIndex> globalPointsPtr(new globalIndex(myUniquePoints));
2393 
2394 
2395  // 3. Assign global point numbers. Keep slaves unset.
2396  pointToGlobal.setSize(mesh_.nPoints());
2397  pointToGlobal = -1;
2398  uniquePoints.setSize(myUniquePoints);
2399  nMaster = 0;
2400 
2401  forAll(isMaster, meshPointi)
2402  {
2403  if (isMaster[meshPointi])
2404  {
2405  pointToGlobal[meshPointi] = globalPointsPtr().toGlobal(nMaster);
2406  uniquePoints[nMaster] = meshPointi;
2407  nMaster++;
2408  }
2409  }
2410 
2411 
2412  // 4. Push global index for coupled points to slaves.
2413  {
2414  labelList masterToGlobal(pointSlavesMap.constructSize(), -1);
2415 
2416  forAll(pointSlaves, pointi)
2417  {
2418  const labelList& slaves = pointSlaves[pointi];
2419 
2420  if (slaves.size() > 0)
2421  {
2422  // Duplicate master globalpoint into slave slots
2423  label meshPointi = cpp.meshPoints()[pointi];
2424  masterToGlobal[pointi] = pointToGlobal[meshPointi];
2425  forAll(slaves, i)
2426  {
2427  masterToGlobal[slaves[i]] = masterToGlobal[pointi];
2428  }
2429  }
2430  }
2431 
2432  // Send back
2433  pointSlavesMap.reverseDistribute(cpp.nPoints(), masterToGlobal);
2434 
2435  // On slave copy master index into overal map.
2436  forAll(pointSlaves, pointi)
2437  {
2438  label meshPointi = cpp.meshPoints()[pointi];
2439 
2440  if (!isMaster[meshPointi])
2441  {
2442  pointToGlobal[meshPointi] = masterToGlobal[pointi];
2443  }
2444  }
2445  }
2446 
2447  return globalPointsPtr;
2448 }
2449 
2450 
2452 (
2453  const labelUList& meshPoints,
2454  const Map<label>& /* unused: meshPointMap */,
2455  labelList& pointToGlobal,
2456  labelList& uniqueMeshPoints
2457 ) const
2458 {
2459  const indirectPrimitivePatch& cpp = coupledPatch();
2460  const labelListList& pointSlaves = globalCoPointSlaves();
2461  const mapDistribute& pointSlavesMap = globalCoPointSlavesMap();
2462 
2463 
2464  // The patch points come in two variants:
2465  // - not on a coupled patch so guaranteed unique
2466  // - on a coupled patch
2467  // If the point is on a coupled patch the problem is that the
2468  // master-slave structure (globalPointSlaves etc.) assigns one of the
2469  // coupled points to be the master but this master point is not
2470  // necessarily on the patch itself! (it might just be connected to the
2471  // patch point via coupled patches).
2472 
2473 
2474  // Determine mapping:
2475  // - from patch point to coupled point (or -1)
2476  // - from coupled point to global patch point
2477  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2478 
2479  const globalIndex globalPPoints(meshPoints.size());
2480 
2481  labelList patchToCoupled(meshPoints.size(), -1);
2482  labelList coupledToGlobalPatch(pointSlavesMap.constructSize(), -1);
2483  //label nCoupled = 0;
2484 
2485  // Note: loop over patch since usually smaller
2486  forAll(meshPoints, patchPointi)
2487  {
2488  label meshPointi = meshPoints[patchPointi];
2489 
2490  const auto iter = cpp.meshPointMap().cfind(meshPointi);
2491 
2492  if (iter.good())
2493  {
2494  patchToCoupled[patchPointi] = iter();
2495  coupledToGlobalPatch[iter()] = globalPPoints.toGlobal(patchPointi);
2496  //++nCoupled;
2497  }
2498  }
2499 
2500  //Pout<< "Patch:" << nl
2501  // << " points:" << meshPoints.size() << nl
2502  // << " of which on coupled patch:" << nCoupled << endl;
2503 
2504 
2505  // Determine master of connected points
2506  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2507  // Problem is that the coupled master might not be on the patch. So
2508  // work out the best patch-point master from all connected points.
2509  // - if the coupled master is on the patch it becomes the patch-point master
2510  // - else the slave with the lowest numbered patch point label
2511 
2512  // Get all data on master
2513  pointSlavesMap.distribute(coupledToGlobalPatch);
2514  forAll(pointSlaves, coupledPointi)
2515  {
2516  const labelList& slaves = pointSlaves[coupledPointi];
2517 
2518  if (slaves.size() > 0)
2519  {
2520  // I am master. What is the best candidate for patch-point master
2521  label masterI = labelMax;
2522  if (coupledToGlobalPatch[coupledPointi] != -1)
2523  {
2524  // I am master and on the coupled patch. Use me.
2525  masterI = coupledToGlobalPatch[coupledPointi];
2526  }
2527  else
2528  {
2529  // Get min of slaves as master.
2530  forAll(slaves, i)
2531  {
2532  label slavePp = coupledToGlobalPatch[slaves[i]];
2533  if (slavePp != -1 && slavePp < masterI)
2534  {
2535  masterI = slavePp;
2536  }
2537  }
2538  }
2539 
2540  if (masterI != labelMax)
2541  {
2542  // Push back
2543  coupledToGlobalPatch[coupledPointi] = masterI;
2544  forAll(slaves, i)
2545  {
2546  coupledToGlobalPatch[slaves[i]] = masterI;
2547  }
2548  }
2549  }
2550  }
2551  pointSlavesMap.reverseDistribute(cpp.nPoints(), coupledToGlobalPatch);
2552 
2553 
2554  // Generate compact numbering for master points
2555  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2556  // Now coupledToGlobalPatch is the globalIndex of the master point.
2557  // Now every processor can check whether they hold it and generate a
2558  // compact numbering.
2559 
2560  label nMasters = 0;
2561  forAll(meshPoints, patchPointi)
2562  {
2563  if (patchToCoupled[patchPointi] == -1)
2564  {
2565  nMasters++;
2566  }
2567  else
2568  {
2569  label coupledPointi = patchToCoupled[patchPointi];
2570  if
2571  (
2572  globalPPoints.toGlobal(patchPointi)
2573  == coupledToGlobalPatch[coupledPointi]
2574  )
2575  {
2576  // I am the master
2577  nMasters++;
2578  }
2579  }
2580  }
2581 
2582  autoPtr<globalIndex> globalPointsPtr(new globalIndex(nMasters));
2583 
2584  //Pout<< "Patch:" << nl
2585  // << " points:" << meshPoints.size() << nl
2586  // << " of which on coupled patch:" << nCoupled << nl
2587  // << " of which master:" << nMasters << endl;
2588 
2589 
2590 
2591  // Push back compact numbering for master point onto slaves
2592  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
2593 
2594  pointToGlobal.setSize(meshPoints.size());
2595  pointToGlobal = -1;
2596  uniqueMeshPoints.setSize(nMasters);
2597 
2598  // Sync master in global point numbering so all know the master point.
2599  // Initialise globalMaster to be -1 except at a globalMaster.
2600  labelList globalMaster(cpp.nPoints(), -1);
2601 
2602  nMasters = 0;
2603  forAll(meshPoints, patchPointi)
2604  {
2605  if (patchToCoupled[patchPointi] == -1)
2606  {
2607  uniqueMeshPoints[nMasters++] = meshPoints[patchPointi];
2608  }
2609  else
2610  {
2611  label coupledPointi = patchToCoupled[patchPointi];
2612  if
2613  (
2614  globalPPoints.toGlobal(patchPointi)
2615  == coupledToGlobalPatch[coupledPointi]
2616  )
2617  {
2618  globalMaster[coupledPointi] =
2619  globalPointsPtr().toGlobal(nMasters);
2620  uniqueMeshPoints[nMasters++] = meshPoints[patchPointi];
2621  }
2622  }
2623  }
2624 
2625 
2626  // Sync by taking max
2627  syncData
2628  (
2629  globalMaster,
2630  pointSlaves,
2631  labelListList(0), // no transforms
2632  pointSlavesMap,
2633  maxEqOp<label>()
2634  );
2635 
2636 
2637  // Now everyone has the master point in globalPointsPtr numbering. Fill
2638  // in the pointToGlobal map.
2639  nMasters = 0;
2640  forAll(meshPoints, patchPointi)
2641  {
2642  if (patchToCoupled[patchPointi] == -1)
2643  {
2644  pointToGlobal[patchPointi] = globalPointsPtr().toGlobal(nMasters++);
2645  }
2646  else
2647  {
2648  label coupledPointi = patchToCoupled[patchPointi];
2649  pointToGlobal[patchPointi] = globalMaster[coupledPointi];
2650 
2651  if
2652  (
2653  globalPPoints.toGlobal(patchPointi)
2654  == coupledToGlobalPatch[coupledPointi]
2655  )
2656  {
2657  nMasters++;
2658  }
2659  }
2660  }
2661 
2662  return globalPointsPtr;
2663 }
2664 
2665 
2666 void Foam::globalMeshData::movePoints(const pointField& newPoints)
2667 {
2668  // Topology does not change and we don't store any geometry so nothing
2669  // needs to be done.
2670  // Only global transformations might change but this is not really
2671  // supported.
2672 }
2673 
2674 
2676 {
2677  // Clear out old data
2678  clearOut();
2679 
2680  // Do processor patch addressing
2681  initProcAddr();
2682 
2683  scalar tolDim = matchTol_ * mesh_.bounds().mag();
2684 
2685  if (debug)
2686  {
2687  Pout<< "globalMeshData : merge dist:" << tolDim << endl;
2688  }
2689 
2690 
2691  const label comm = mesh_.comm();
2692  const label oldWarnComm = UPstream::commWarn(comm);
2693 
2694  if (UPstream::is_parallel(comm))
2695  {
2696  const label myProci = UPstream::myProcNo(comm);
2697  const label numProc = UPstream::nProcs(comm);
2698 
2699  // Gather all three sizes together
2700  labelList allSizes(3*numProc);
2701  {
2702  label* tup = allSizes.begin(3*myProci);
2703  tup[0] = mesh_.nPoints();
2704  tup[1] = mesh_.nFaces();
2705  tup[2] = mesh_.nCells();
2706  }
2707 
2708  UPstream::mpiAllGather(allSizes.data(), 3, comm);
2709 
2710  // Extract counts per mesh entity
2711  // TBD: check for label overflow?
2712 
2713  labelList counts(numProc);
2714  for (label proci = 0, idx = 0; proci < numProc; ++proci, idx += 3)
2715  {
2716  counts[proci] = allSizes[idx];
2717  }
2718  globalMeshPointAddr_.reset(counts);
2719 
2720  for (label proci = 0, idx = 1; proci < numProc; ++proci, idx += 3)
2721  {
2722  counts[proci] = allSizes[idx];
2723  }
2724  globalMeshFaceAddr_.reset(counts);
2725 
2726  for (label proci = 0, idx = 2; proci < numProc; ++proci, idx += 3)
2727  {
2728  counts[proci] = allSizes[idx];
2729  }
2730  globalMeshCellAddr_.reset(counts);
2731  }
2732  else
2733  {
2734  globalMeshPointAddr_.reset(globalIndex::gatherNone{}, mesh_.nPoints());
2735  globalMeshFaceAddr_.reset(globalIndex::gatherNone{}, mesh_.nFaces());
2736  globalMeshCellAddr_.reset(globalIndex::gatherNone{}, mesh_.nCells());
2737  }
2738 
2739  // Restore communicator settings
2740  UPstream::commWarn(oldWarnComm);
2741 
2742  if (debug)
2743  {
2744  Info<< "globalMeshData : Total points/faces/cells : ("
2745  << nTotalPoints() << ' '
2746  << nTotalFaces() << ' '
2747  << nTotalCells() << ')' << endl;
2748  }
2749 }
2750 
2751 
2752 // ************************************************************************* //
const labelList & sharedPointLabels() const
Return indices of local points that are globally shared.
const Map< label > & coupledPatchMeshEdgeMap() const
Return map from mesh edges to coupledPatch edges.
void reset(const label localSize, const label comm=UPstream::worldComm, const bool parallel=UPstream::parRun())
Reset from local size, using gather/broadcast with default/specified communicator if parallel...
Definition: globalIndex.C:346
label comm() const noexcept
The communicator used for parallel communication.
Definition: polyMesh.H:701
const labelListList & globalPointSlaves() const
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
const mapDistribute & globalCoPointSlavesMap() const
const labelListList & globalPointTransformedSlaves() const
pointField geometricSharedPoints() const
Like sharedPoints but keeps cyclic points separate. (does geometric merging; uses matchTol_*bb as mer...
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
const globalIndex & globalBoundaryCellNumbering() const
Numbering of boundary cells is according to boundaryCells()
labelList pointLabels(nPoints, -1)
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:153
const labelList & sharedPointGlobalLabels() const
Return shared point global labels. Tries to read &#39;pointProcAddressing&#39; and returns list or -1 if none...
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:521
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
const labelListList & globalPointTransformedBoundaryFaces() const
const mapDistribute & globalEdgeSlavesMap() const
static const Foam::scalar matchTol_
Geometric tolerance (fraction of bounding box)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void clearOut()
Remove all demand driven data.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:65
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:411
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void mpiAllGather(char *allData, int count, const label communicator=worldComm)
Gather/scatter identically-sized char data.
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1061
const mapDistribute & globalPointSlavesMap() const
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
const bitSet & globalEdgeOrientation() const
Is my edge same orientation as master edge.
const labelListList & globalPointBoundaryCells() const
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:168
List< labelPair > labelPairList
List of labelPair.
Definition: labelPair.H:33
autoPtr< globalIndex > mergePoints(labelList &pointToGlobal, labelList &uniquePoints) const
Helper for merging (collocated!) mesh point data.
const labelListList & globalEdgeSlaves() const
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:1086
const labelList & coupledPatchMeshEdges() const
Return map from coupledPatch edges to mesh edges.
PrimitivePatch< IndirectList< face >, const pointField & > indirectPrimitivePatch
A PrimitivePatch with an IndirectList for the faces, const reference for the point field...
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
const labelListList & globalPointBoundaryFaces() const
const dimensionedScalar me
Electron mass.
static bool less(const vector &x, const vector &y)
To compare normals.
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const labelList & boundaryCells() const
From boundary cell to mesh cell.
static label commWarn(const label communicator) noexcept
Alter communicator debugging setting. Warns for use of any communicator differing from specified...
Definition: UPstream.H:461
label nGlobalEdges() const
Return number of globally shared edges.
List helper to append y elements onto the end of x.
Definition: ListOps.H:712
scalar y
A list of faces which address into the list of points.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:61
const globalIndex & globalEdgeNumbering() const
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run. ...
Definition: UPstream.H:1077
label nGlobalPoints() const
Return number of globally shared points.
static void recv(Type &value, const int fromProcNo, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm, IOstreamOption::streamFormat fmt=IOstreamOption::BINARY)
Receive and deserialize a value. Uses operator>> for de-serialization.
Definition: IPstream.H:81
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const labelList & sharedEdgeAddr() const
Return addressing into the complete globally shared edge list.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
void setSize(const label n)
Alias for resize()
Definition: List.H:320
dynamicFvMesh & mesh
const globalIndexAndTransform & globalTransforms() const
Global transforms numbering.
"scheduled" (MPI standard) : (MPI_Send, MPI_Recv)
static int compare(const Pair< label > &a, const Pair< label > &b)
Compare Pairs.
Definition: PairI.H:24
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
label nPoints
void sort(UList< T > &list)
Sort the list.
Definition: UList.C:296
static void combineReduce(T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors...
static constexpr int masterNo() noexcept
Relative rank for the master process - is always 0.
Definition: UPstream.H:1071
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
const globalIndex & globalPointNumbering() const
Numbering of coupled points is according to coupledPatch.
static bool is_parallel(const label communicator=worldComm)
True if parallel algorithm or exchange is required.
Definition: UPstream.H:1123
label inplaceMergePoints(PointList &points, const scalar mergeTol, const bool verbose, labelList &pointToUnique)
Inplace merge points, preserving the original point order. All points closer/equal mergeTol are to be...
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void reset(const dimensionSet &ds)
Copy assign the exponents from the dimensionSet.
Definition: dimensionSet.C:142
int debug
Static debugging option.
Pair< label > labelPair
A pair of labels.
Definition: Pair.H:51
~globalMeshData()
Destructor.
defineTypeNameAndDebug(combustionModel, 0)
Geometric merging of points. See below.
labelList f(nPoints)
Define the processor-processor connection table by walking a list of patches and detecting the proces...
void movePoints(const pointField &newPoints)
Update for moving points.
const indirectPrimitivePatch & coupledPatch() const
Return patch of all coupled faces.
const labelListList & globalPointTransformedBoundaryCells() const
void operator()(T &x, const T &y) const
Definition: ops.H:76
Class containing processor-to-processor mapping information.
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:59
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
void updateMesh()
Change global mesh data given a topological change. Does a.
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
label n
static rangeType subProcs(const label communicator=worldComm)
Range of process indices for sub-processes.
Definition: UPstream.H:1197
List< label > labelList
A List of labels.
Definition: List.H:62
const mapDistribute & globalPointBoundaryFacesMap() const
dimensionSet transform(const dimensionSet &ds)
Return the argument; transformations do not change the dimensions.
Definition: dimensionSet.C:521
const labelList & sharedPointAddr() const
Return addressing into the complete globally shared points list.
bool send()
Send buffer contents now and not in destructor [advanced usage]. Returns true on success.
Definition: OPstreams.C:84
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
pointField sharedPoints() const
Collect coordinates of shared points on all processors. (does parallel communication!) ...
PrimitivePatch< List< face >, const pointField > bMesh
Holder of faceList and points. (v.s. e.g. primitivePatch which references points) ...
Definition: bMesh.H:39
Map< label > invertToMap(const labelUList &values)
Create inverse mapping, which is a lookup table into the given list.
Definition: ListOps.C:107
const mapDistribute & globalPointBoundaryCellsMap() const
bool found
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
const labelList & sharedEdgeLabels() const
Return indices of local edges that are globally shared.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Determination and storage of the possible independent transforms introduced by coupledPolyPatches, as well as all of the possible permutations of these transforms generated by the presence of multiple coupledPolyPatches, i.e. more than one cyclic boundary. Note that any given point can be on maximum 3 transforms only (and these transforms have to be perpendicular)
const labelListList & globalCoPointSlaves() const
A HashTable to objects of type <T> with a label key.
const labelListList & globalEdgeTransformedSlaves() const
globalMeshData(const globalMeshData &)=delete
No copy construct.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const globalIndex & globalBoundaryFaceNumbering() const
Numbering of boundary faces is face-mesh.nInternalFaces()