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