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