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