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