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