meshToMeshParallelOps.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) 2012-2014 OpenFOAM Foundation
6  \\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
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 "meshToMesh.H"
27 #include "OFstream.H"
28 #include "Time.H"
29 #include "globalIndex.H"
30 #include "mergePoints.H"
31 #include "processorPolyPatch.H"
32 #include "SubField.H"
33 #include "AABBTree.H"
34 
35 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
36 
38 (
39  const polyMesh& src,
40  const polyMesh& tgt
41 ) const
42 {
43  label procI = 0;
44 
45  if (Pstream::parRun())
46  {
47  List<label> cellsPresentOnProc(Pstream::nProcs(), 0);
48  if ((src.nCells() > 0) || (tgt.nCells() > 0))
49  {
50  cellsPresentOnProc[Pstream::myProcNo()] = 1;
51  }
52  else
53  {
54  cellsPresentOnProc[Pstream::myProcNo()] = 0;
55  }
56 
57  Pstream::gatherList(cellsPresentOnProc);
58  Pstream::scatterList(cellsPresentOnProc);
59 
60  label nHaveCells = sum(cellsPresentOnProc);
61 
62  if (nHaveCells > 1)
63  {
64  procI = -1;
65  if (debug)
66  {
67  Info<< "meshToMesh::calcDistribution: "
68  << "Meshes split across multiple processors" << endl;
69  }
70  }
71  else if (nHaveCells == 1)
72  {
73  procI = findIndex(cellsPresentOnProc, 1);
74  if (debug)
75  {
76  Info<< "meshToMesh::calcDistribution: "
77  << "Meshes local to processor" << procI << endl;
78  }
79  }
80  }
81 
82  return procI;
83 }
84 
85 
87 (
88  const List<treeBoundBoxList>& procBb,
89  const boundBox& bb,
90  boolList& overlaps
91 ) const
92 {
93  overlaps = false;
94 
95  label nOverlaps = 0;
96 
97  forAll(procBb, procI)
98  {
99  const treeBoundBoxList& bbp = procBb[procI];
100 
101  forAll(bbp, bbI)
102  {
103  if (bbp[bbI].overlaps(bb))
104  {
105  overlaps[procI] = true;
106  nOverlaps++;
107  break;
108  }
109  }
110  }
111 
112  return nOverlaps;
113 }
114 
115 
117 (
118  const polyMesh& src,
119  const polyMesh& tgt
120 ) const
121 {
122  // get decomposition of cells on src mesh
123  List<treeBoundBoxList> procBb(Pstream::nProcs());
124 
125  if (src.nCells() > 0)
126  {
127  procBb[Pstream::myProcNo()] = AABBTree<labelList>
128  (
129  src.cellPoints(),
130  src.points(),
131  false
132  ).boundBoxes();
133  }
134  else
135  {
136  procBb[Pstream::myProcNo()] = treeBoundBoxList();
137  }
138 
139 
140  Pstream::gatherList(procBb);
141  Pstream::scatterList(procBb);
142 
143 
144  if (debug)
145  {
146  Info<< "Determining extent of src mesh per processor:" << nl
147  << "\tproc\tbb" << endl;
148  forAll(procBb, procI)
149  {
150  Info<< '\t' << procI << '\t' << procBb[procI] << endl;
151  }
152  }
153 
154 
155  // determine which cells of tgt mesh overlaps src mesh per proc
156  const cellList& cells = tgt.cells();
157  const faceList& faces = tgt.faces();
158  const pointField& points = tgt.points();
159 
160  labelListList sendMap;
161 
162  {
163  // per processor indices into all segments to send
164  List<DynamicList<label> > dynSendMap(Pstream::nProcs());
165  label iniSize = floor(tgt.nCells()/Pstream::nProcs());
166 
167  forAll(dynSendMap, procI)
168  {
169  dynSendMap[procI].setCapacity(iniSize);
170  }
171 
172  // work array - whether src processor bb overlaps the tgt cell bounds
173  boolList procBbOverlaps(Pstream::nProcs());
174  forAll(cells, cellI)
175  {
176  const cell& c = cells[cellI];
177 
178  // determine bounding box of tgt cell
179  boundBox cellBb(point::max, point::min);
180  forAll(c, faceI)
181  {
182  const face& f = faces[c[faceI]];
183  forAll(f, fp)
184  {
185  cellBb.min() = min(cellBb.min(), points[f[fp]]);
186  cellBb.max() = max(cellBb.max(), points[f[fp]]);
187  }
188  }
189 
190  // find the overlapping tgt cells on each src processor
191  (void)calcOverlappingProcs(procBb, cellBb, procBbOverlaps);
192 
193  forAll(procBbOverlaps, procI)
194  {
195  if (procBbOverlaps[procI])
196  {
197  dynSendMap[procI].append(cellI);
198  }
199  }
200  }
201 
202  // convert dynamicList to labelList
203  sendMap.setSize(Pstream::nProcs());
204  forAll(sendMap, procI)
205  {
206  sendMap[procI].transfer(dynSendMap[procI]);
207  }
208  }
209 
210  // debug printing
211  if (debug)
212  {
213  Pout<< "Of my " << cells.size() << " target cells I need to send to:"
214  << nl << "\tproc\tcells" << endl;
215  forAll(sendMap, procI)
216  {
217  Pout<< '\t' << procI << '\t' << sendMap[procI].size() << endl;
218  }
219  }
220 
221 
222  // send over how many tgt cells I need to receive from each processor
223  labelListList sendSizes(Pstream::nProcs());
224  sendSizes[Pstream::myProcNo()].setSize(Pstream::nProcs());
225  forAll(sendMap, procI)
226  {
227  sendSizes[Pstream::myProcNo()][procI] = sendMap[procI].size();
228  }
229  Pstream::gatherList(sendSizes);
230  Pstream::scatterList(sendSizes);
231 
232 
233  // determine order of receiving
234  labelListList constructMap(Pstream::nProcs());
235 
236  label segmentI = 0;
237  forAll(constructMap, procI)
238  {
239  // what I need to receive is what other processor is sending to me
240  label nRecv = sendSizes[procI][Pstream::myProcNo()];
241  constructMap[procI].setSize(nRecv);
242 
243  for (label i = 0; i < nRecv; i++)
244  {
245  constructMap[procI][i] = segmentI++;
246  }
247  }
248 
250  (
251  new mapDistribute
252  (
253  segmentI, // size after construction
254  sendMap.xfer(),
255  constructMap.xfer()
256  )
257  );
258 
259  return mapPtr;
260 }
261 
262 
264 (
265  const mapDistribute& map,
266  const polyMesh& tgtMesh,
267  const globalIndex& globalI,
269  List<label>& nInternalFaces,
270  List<faceList>& faces,
271  List<labelList>& faceOwner,
272  List<labelList>& faceNeighbour,
273  List<labelList>& cellIDs,
274  List<labelList>& nbrProcIDs,
275  List<labelList>& procLocalFaceIDs
276 ) const
277 {
278  PstreamBuffers pBufs(Pstream::nonBlocking);
279 
280  points.setSize(Pstream::nProcs());
281  nInternalFaces.setSize(Pstream::nProcs(), 0);
282  faces.setSize(Pstream::nProcs());
283  faceOwner.setSize(Pstream::nProcs());
284  faceNeighbour.setSize(Pstream::nProcs());
285  cellIDs.setSize(Pstream::nProcs());
286 
287  nbrProcIDs.setSize(Pstream::nProcs());;
288  procLocalFaceIDs.setSize(Pstream::nProcs());;
289 
290 
291  for (label domain = 0; domain < Pstream::nProcs(); domain++)
292  {
293  const labelList& sendElems = map.subMap()[domain];
294 
295  if (sendElems.size())
296  {
297  // reverse cell map
298  labelList reverseCellMap(tgtMesh.nCells(), -1);
299  forAll(sendElems, subCellI)
300  {
301  reverseCellMap[sendElems[subCellI]] = subCellI;
302  }
303 
304  DynamicList<face> subFaces(tgtMesh.nFaces());
305  DynamicList<label> subFaceOwner(tgtMesh.nFaces());
306  DynamicList<label> subFaceNeighbour(tgtMesh.nFaces());
307 
308  DynamicList<label> subNbrProcIDs(tgtMesh.nFaces());
309  DynamicList<label> subProcLocalFaceIDs(tgtMesh.nFaces());
310 
311  label nInternal = 0;
312 
313  // internal faces
314  forAll(tgtMesh.faceNeighbour(), faceI)
315  {
316  label own = tgtMesh.faceOwner()[faceI];
317  label nbr = tgtMesh.faceNeighbour()[faceI];
318  label subOwn = reverseCellMap[own];
319  label subNbr = reverseCellMap[nbr];
320 
321  if (subOwn != -1 && subNbr != -1)
322  {
323  nInternal++;
324 
325  if (subOwn < subNbr)
326  {
327  subFaces.append(tgtMesh.faces()[faceI]);
328  subFaceOwner.append(subOwn);
329  subFaceNeighbour.append(subNbr);
330  subNbrProcIDs.append(-1);
331  subProcLocalFaceIDs.append(-1);
332  }
333  else
334  {
335  subFaces.append(tgtMesh.faces()[faceI].reverseFace());
336  subFaceOwner.append(subNbr);
337  subFaceNeighbour.append(subOwn);
338  subNbrProcIDs.append(-1);
339  subProcLocalFaceIDs.append(-1);
340  }
341  }
342  }
343 
344  // boundary faces for new region
345  forAll(tgtMesh.faceNeighbour(), faceI)
346  {
347  label own = tgtMesh.faceOwner()[faceI];
348  label nbr = tgtMesh.faceNeighbour()[faceI];
349  label subOwn = reverseCellMap[own];
350  label subNbr = reverseCellMap[nbr];
351 
352  if (subOwn != -1 && subNbr == -1)
353  {
354  subFaces.append(tgtMesh.faces()[faceI]);
355  subFaceOwner.append(subOwn);
356  subFaceNeighbour.append(subNbr);
357  subNbrProcIDs.append(-1);
358  subProcLocalFaceIDs.append(-1);
359  }
360  else if (subOwn == -1 && subNbr != -1)
361  {
362  subFaces.append(tgtMesh.faces()[faceI].reverseFace());
363  subFaceOwner.append(subNbr);
364  subFaceNeighbour.append(subOwn);
365  subNbrProcIDs.append(-1);
366  subProcLocalFaceIDs.append(-1);
367  }
368  }
369 
370  // boundary faces of existing region
371  forAll(tgtMesh.boundaryMesh(), patchI)
372  {
373  const polyPatch& pp = tgtMesh.boundaryMesh()[patchI];
374 
375  label nbrProcI = -1;
376 
377  // store info for faces on processor patches
378  if (isA<processorPolyPatch>(pp))
379  {
380  const processorPolyPatch& ppp =
381  dynamic_cast<const processorPolyPatch&>(pp);
382 
383  nbrProcI = ppp.neighbProcNo();
384  }
385 
386  forAll(pp, i)
387  {
388  label faceI = pp.start() + i;
389  label own = tgtMesh.faceOwner()[faceI];
390 
391  if (reverseCellMap[own] != -1)
392  {
393  subFaces.append(tgtMesh.faces()[faceI]);
394  subFaceOwner.append(reverseCellMap[own]);
395  subFaceNeighbour.append(-1);
396  subNbrProcIDs.append(nbrProcI);
397  subProcLocalFaceIDs.append(i);
398  }
399  }
400  }
401 
402  // reverse point map
403  labelList reversePointMap(tgtMesh.nPoints(), -1);
404  DynamicList<point> subPoints(tgtMesh.nPoints());
405  forAll(subFaces, subFaceI)
406  {
407  face& f = subFaces[subFaceI];
408  forAll(f, fp)
409  {
410  label pointI = f[fp];
411  if (reversePointMap[pointI] == -1)
412  {
413  reversePointMap[pointI] = subPoints.size();
414  subPoints.append(tgtMesh.points()[pointI]);
415  }
416 
417  f[fp] = reversePointMap[pointI];
418  }
419  }
420 
421  // tgt cells into global numbering
422  labelList globalElems(sendElems.size());
423  forAll(sendElems, i)
424  {
425  if (debug > 1)
426  {
427  Pout<< "tgtProc:" << Pstream::myProcNo()
428  << " sending tgt cell " << sendElems[i]
429  << "[" << globalI.toGlobal(sendElems[i]) << "]"
430  << " to srcProc " << domain << endl;
431  }
432 
433  globalElems[i] = globalI.toGlobal(sendElems[i]);
434  }
435 
436  // pass data
437  if (domain == Pstream::myProcNo())
438  {
439  // allocate my own data
440  points[Pstream::myProcNo()] = subPoints;
441  nInternalFaces[Pstream::myProcNo()] = nInternal;
442  faces[Pstream::myProcNo()] = subFaces;
443  faceOwner[Pstream::myProcNo()] = subFaceOwner;
444  faceNeighbour[Pstream::myProcNo()] = subFaceNeighbour;
445  cellIDs[Pstream::myProcNo()] = globalElems;
446  nbrProcIDs[Pstream::myProcNo()] = subNbrProcIDs;
447  procLocalFaceIDs[Pstream::myProcNo()] = subProcLocalFaceIDs;
448  }
449  else
450  {
451  // send data to other processor domains
452  UOPstream toDomain(domain, pBufs);
453 
454  toDomain
455  << subPoints
456  << nInternal
457  << subFaces
458  << subFaceOwner
459  << subFaceNeighbour
460  << globalElems
461  << subNbrProcIDs
462  << subProcLocalFaceIDs;
463  }
464  }
465  }
466 
467  // Start receiving
468  pBufs.finishedSends();
469 
470  // Consume
471  for (label domain = 0; domain < Pstream::nProcs(); domain++)
472  {
473  const labelList& recvElems = map.constructMap()[domain];
474 
475  if (domain != Pstream::myProcNo() && recvElems.size())
476  {
477  UIPstream str(domain, pBufs);
478 
479  str >> points[domain]
480  >> nInternalFaces[domain]
481  >> faces[domain]
482  >> faceOwner[domain]
483  >> faceNeighbour[domain]
484  >> cellIDs[domain]
485  >> nbrProcIDs[domain]
486  >> procLocalFaceIDs[domain];
487  }
488 
489  if (debug)
490  {
491  Pout<< "Target mesh send sizes[" << domain << "]"
492  << ": points="<< points[domain].size()
493  << ", faces=" << faces[domain].size()
494  << ", nInternalFaces=" << nInternalFaces[domain]
495  << ", faceOwn=" << faceOwner[domain].size()
496  << ", faceNbr=" << faceNeighbour[domain].size()
497  << ", cellIDs=" << cellIDs[domain].size() << endl;
498  }
499  }
500 }
501 
502 
504 (
505  const mapDistribute& map,
506  const polyMesh& tgt,
507  const globalIndex& globalI,
508  pointField& tgtPoints,
509  faceList& tgtFaces,
510  labelList& tgtFaceOwners,
511  labelList& tgtFaceNeighbours,
512  labelList& tgtCellIDs
513 ) const
514 {
515  // Exchange per-processor data
517  List<label> allNInternalFaces;
518  List<faceList> allFaces;
519  List<labelList> allFaceOwners;
520  List<labelList> allFaceNeighbours;
521  List<labelList> allTgtCellIDs;
522 
523  // Per target mesh face the neighbouring proc and index in
524  // processor patch (all -1 for normal boundary face)
525  List<labelList> allNbrProcIDs;
526  List<labelList> allProcLocalFaceIDs;
527 
528  distributeCells
529  (
530  map,
531  tgt,
532  globalI,
533  allPoints,
534  allNInternalFaces,
535  allFaces,
536  allFaceOwners,
537  allFaceNeighbours,
538  allTgtCellIDs,
539  allNbrProcIDs,
540  allProcLocalFaceIDs
541  );
542 
543  // Convert lists into format that can be used to generate a valid polyMesh
544  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
545  //
546  // Points and cells are collected into single flat lists:
547  // - i.e. proc0, proc1 ... procN
548  //
549  // Faces need to be sorted after collection to that internal faces are
550  // contiguous, followed by all boundary faces
551  //
552  // Processor patch faces between included cells on neighbouring processors
553  // are converted into internal faces
554  //
555  // Face list structure:
556  // - Per processor:
557  // - internal faces
558  // - processor faces that have been converted into internal faces
559  // - Followed by all boundary faces
560  // - from 'normal' boundary faces
561  // - from singularly-sided processor patch faces
562 
563 
564  // Number of internal+coupled faces
565  labelList allNIntCoupledFaces(allNInternalFaces);
566 
567  // Starting offset for points
568  label nPoints = 0;
569  labelList pointOffset(Pstream::nProcs(), 0);
570  forAll(allPoints, procI)
571  {
572  pointOffset[procI] = nPoints;
573  nPoints += allPoints[procI].size();
574  }
575 
576  // Starting offset for cells
577  label nCells = 0;
578  labelList cellOffset(Pstream::nProcs(), 0);
579  forAll(allTgtCellIDs, procI)
580  {
581  cellOffset[procI] = nCells;
582  nCells += allTgtCellIDs[procI].size();
583  }
584 
585  // Count any coupled faces
586  typedef FixedList<label, 3> label3;
587  typedef HashTable<label, label3, label3::Hash<> > procCoupleInfo;
588  procCoupleInfo procFaceToGlobalCell;
589 
590  forAll(allNbrProcIDs, procI)
591  {
592  const labelList& nbrProcI = allNbrProcIDs[procI];
593  const labelList& localFaceI = allProcLocalFaceIDs[procI];
594 
595  forAll(nbrProcI, i)
596  {
597  if (nbrProcI[i] != -1 && localFaceI[i] != -1)
598  {
599  label3 key;
600  key[0] = min(procI, nbrProcI[i]);
601  key[1] = max(procI, nbrProcI[i]);
602  key[2] = localFaceI[i];
603 
604  procCoupleInfo::const_iterator fnd =
605  procFaceToGlobalCell.find(key);
606 
607  if (fnd == procFaceToGlobalCell.end())
608  {
609  procFaceToGlobalCell.insert(key, -1);
610  }
611  else
612  {
613  if (debug > 1)
614  {
615  Pout<< "Additional internal face between procs:"
616  << key[0] << " and " << key[1]
617  << " across local face " << key[2] << endl;
618  }
619 
620  allNIntCoupledFaces[procI]++;
621  }
622  }
623  }
624  }
625 
626 
627  // Starting offset for internal faces
628  label nIntFaces = 0;
629  label nFacesTotal = 0;
630  labelList internalFaceOffset(Pstream::nProcs(), 0);
631  forAll(allNIntCoupledFaces, procI)
632  {
633  label nCoupledFaces =
634  allNIntCoupledFaces[procI] - allNInternalFaces[procI];
635 
636  internalFaceOffset[procI] = nIntFaces;
637  nIntFaces += allNIntCoupledFaces[procI];
638  nFacesTotal += allFaceOwners[procI].size() - nCoupledFaces;
639  }
640 
641  tgtPoints.setSize(nPoints);
642  tgtFaces.setSize(nFacesTotal);
643  tgtFaceOwners.setSize(nFacesTotal);
644  tgtFaceNeighbours.setSize(nFacesTotal);
645  tgtCellIDs.setSize(nCells);
646 
647  // Insert points
648  forAll(allPoints, procI)
649  {
650  const pointField& pts = allPoints[procI];
651  SubList<point>(tgtPoints, pts.size(), pointOffset[procI]).assign(pts);
652  }
653 
654  // Insert cellIDs
655  forAll(allTgtCellIDs, procI)
656  {
657  const labelList& cellIDs = allTgtCellIDs[procI];
658  SubList<label>(tgtCellIDs, cellIDs.size(), cellOffset[procI]).assign
659  (
660  cellIDs
661  );
662  }
663 
664 
665  // Insert internal faces (from internal faces)
666  forAll(allFaces, procI)
667  {
668  const faceList& fcs = allFaces[procI];
669  const labelList& faceOs = allFaceOwners[procI];
670  const labelList& faceNs = allFaceNeighbours[procI];
671 
672  SubList<face> slice
673  (
674  tgtFaces,
675  allNInternalFaces[procI],
676  internalFaceOffset[procI]
677  );
678  slice.assign(SubList<face>(fcs, allNInternalFaces[procI]));
679  forAll(slice, i)
680  {
681  add(slice[i], pointOffset[procI]);
682  }
683 
684  SubField<label> ownSlice
685  (
686  tgtFaceOwners,
687  allNInternalFaces[procI],
688  internalFaceOffset[procI]
689  );
690  ownSlice.assign(SubField<label>(faceOs, allNInternalFaces[procI]));
691  add(ownSlice, cellOffset[procI]);
692 
693  SubField<label> nbrSlice
694  (
695  tgtFaceNeighbours,
696  allNInternalFaces[procI],
697  internalFaceOffset[procI]
698  );
699  nbrSlice.assign(SubField<label>(faceNs, allNInternalFaces[procI]));
700  add(nbrSlice, cellOffset[procI]);
701 
702  internalFaceOffset[procI] += allNInternalFaces[procI];
703  }
704 
705 
706  // Insert internal faces (from coupled face-pairs)
707  forAll(allNbrProcIDs, procI)
708  {
709  const labelList& nbrProcI = allNbrProcIDs[procI];
710  const labelList& localFaceI = allProcLocalFaceIDs[procI];
711  const labelList& faceOs = allFaceOwners[procI];
712  const faceList& fcs = allFaces[procI];
713 
714  forAll(nbrProcI, i)
715  {
716  if (nbrProcI[i] != -1 && localFaceI[i] != -1)
717  {
718  label3 key;
719  key[0] = min(procI, nbrProcI[i]);
720  key[1] = max(procI, nbrProcI[i]);
721  key[2] = localFaceI[i];
722 
723  procCoupleInfo::iterator fnd = procFaceToGlobalCell.find(key);
724 
725  if (fnd != procFaceToGlobalCell.end())
726  {
727  label tgtFaceI = fnd();
728  if (tgtFaceI == -1)
729  {
730  // on first visit store the new cell on this side
731  fnd() = cellOffset[procI] + faceOs[i];
732  }
733  else
734  {
735  // get owner and neighbour in new cell numbering
736  label newOwn = cellOffset[procI] + faceOs[i];
737  label newNbr = fnd();
738  label tgtFaceI = internalFaceOffset[procI]++;
739 
740  if (debug > 1)
741  {
742  Pout<< " proc " << procI
743  << "\tinserting face:" << tgtFaceI
744  << " connection between owner " << newOwn
745  << " and neighbour " << newNbr
746  << endl;
747  }
748 
749  if (newOwn < newNbr)
750  {
751  // we have correct orientation
752  tgtFaces[tgtFaceI] = fcs[i];
753  tgtFaceOwners[tgtFaceI] = newOwn;
754  tgtFaceNeighbours[tgtFaceI] = newNbr;
755  }
756  else
757  {
758  // reverse orientation
759  tgtFaces[tgtFaceI] = fcs[i].reverseFace();
760  tgtFaceOwners[tgtFaceI] = newNbr;
761  tgtFaceNeighbours[tgtFaceI] = newOwn;
762  }
763 
764  add(tgtFaces[tgtFaceI], pointOffset[procI]);
765 
766  // mark with unique value
767  fnd() = -2;
768  }
769  }
770  }
771  }
772  }
773 
774 
775  forAll(allNbrProcIDs, procI)
776  {
777  const labelList& nbrProcI = allNbrProcIDs[procI];
778  const labelList& localFaceI = allProcLocalFaceIDs[procI];
779  const labelList& faceOs = allFaceOwners[procI];
780  const labelList& faceNs = allFaceNeighbours[procI];
781  const faceList& fcs = allFaces[procI];
782 
783  forAll(nbrProcI, i)
784  {
785  // coupled boundary face
786  if (nbrProcI[i] != -1 && localFaceI[i] != -1)
787  {
788  label3 key;
789  key[0] = min(procI, nbrProcI[i]);
790  key[1] = max(procI, nbrProcI[i]);
791  key[2] = localFaceI[i];
792 
793  label tgtFaceI = procFaceToGlobalCell[key];
794 
795  if (tgtFaceI == -1)
796  {
798  << "Unvisited " << key
799  << abort(FatalError);
800  }
801  else if (tgtFaceI != -2)
802  {
803  label newOwn = cellOffset[procI] + faceOs[i];
804  label tgtFaceI = nIntFaces++;
805 
806  if (debug > 1)
807  {
808  Pout<< " proc " << procI
809  << "\tinserting boundary face:" << tgtFaceI
810  << " from coupled face " << key
811  << endl;
812  }
813 
814  tgtFaces[tgtFaceI] = fcs[i];
815  add(tgtFaces[tgtFaceI], pointOffset[procI]);
816 
817  tgtFaceOwners[tgtFaceI] = newOwn;
818  tgtFaceNeighbours[tgtFaceI] = -1;
819  }
820  }
821  // normal boundary face
822  else
823  {
824  label own = faceOs[i];
825  label nbr = faceNs[i];
826  if ((own != -1) && (nbr == -1))
827  {
828  label newOwn = cellOffset[procI] + faceOs[i];
829  label tgtFaceI = nIntFaces++;
830 
831  tgtFaces[tgtFaceI] = fcs[i];
832  add(tgtFaces[tgtFaceI], pointOffset[procI]);
833 
834  tgtFaceOwners[tgtFaceI] = newOwn;
835  tgtFaceNeighbours[tgtFaceI] = -1;
836  }
837  }
838  }
839  }
840 
841 
842  if (debug)
843  {
844  // only merging points in debug mode
845 
846  labelList oldToNew;
847  pointField newTgtPoints;
848  bool hasMerged = mergePoints
849  (
850  tgtPoints,
851  SMALL,
852  false,
853  oldToNew,
854  newTgtPoints
855  );
856 
857  if (hasMerged)
858  {
859  if (debug)
860  {
861  Pout<< "Merged from " << tgtPoints.size()
862  << " down to " << newTgtPoints.size() << " points" << endl;
863  }
864 
865  tgtPoints.transfer(newTgtPoints);
866  forAll(tgtFaces, i)
867  {
868  inplaceRenumber(oldToNew, tgtFaces[i]);
869  }
870  }
871  }
872 }
873 
874 
875 // ************************************************************************* //
Foam::mapDistributeBase::subMap
const labelListList & subMap() const
From subsetted data back to original data.
Definition: mapDistributeBase.H:256
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:979
Foam::UOPstream
Output inter-processor communications stream operating on external buffer.
Definition: UOPstream.H:54
SubField.H
Foam::boundBox::max
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::primitiveMesh::cellPoints
const labelListList & cellPoints() const
Definition: primitiveMeshCellPoints.C:31
Foam::meshToMesh::distributeAndMergeCells
void distributeAndMergeCells(const mapDistribute &map, const polyMesh &tgt, const globalIndex &globalI, pointField &tgtPoints, faceList &tgtFaces, labelList &tgtFaceOwners, labelList &tgtFaceNeighbours, labelList &tgtCellIDs) const
Collect pieces of tgt mesh from other procssors and restructure.
Definition: meshToMeshParallelOps.C:504
Foam::findIndex
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurence of given element and return index,.
Foam::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::List::xfer
Xfer< List< T > > xfer()
Transfer contents to the Xfer container.
Definition: ListI.H:90
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
globalIndex.H
Foam::meshToMesh::distributeCells
void distributeCells(const mapDistribute &map, const polyMesh &tgtMesh, const globalIndex &globalI, List< pointField > &points, List< label > &nInternalFaces, List< faceList > &faces, List< labelList > &faceOwner, List< labelList > &faceNeighbour, List< labelList > &cellIDs, List< labelList > &nbrProcIDs, List< labelList > &procLocalFaceIDs) const
Distribute mesh info from 'my' processor to others.
Definition: meshToMeshParallelOps.C:264
Foam::PstreamBuffers
Buffers for inter-processor communications streams (UOPstream, UIPstream).
Definition: PstreamBuffers.H:85
meshToMesh.H
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:136
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Foam::AABBTree::boundBoxes
const List< treeBoundBox > & boundBoxes() const
Return the bounding boxes making up the tree.
Definition: AABBTree.C:428
Foam::processorPolyPatch::neighbProcNo
int neighbProcNo() const
Return neigbour processor number.
Definition: processorPolyPatch.H:255
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::meshToMesh::calcDistribution
label calcDistribution(const polyMesh &src, const polyMesh &tgt) const
Determine whether the meshes are split across multiple pocessors.
Definition: meshToMeshParallelOps.C:38
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
OFstream.H
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::SubField
Pre-declare related SubField type.
Definition: Field.H:61
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
Definition: ListOpsTemplates.C:57
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::meshToMesh::calcProcMap
autoPtr< mapDistribute > calcProcMap(const polyMesh &src, const polyMesh &tgt) const
Calculate the mapping between processors.
Definition: meshToMeshParallelOps.C:117
Foam::primitiveMesh::nPoints
label nPoints() const
Definition: primitiveMeshI.H:35
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::List::append
void append(const T &)
Append an element at the end of the list.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::mapDistribute
Class containing processor-to-processor mapping information.
Definition: mapDistribute.H:152
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
Foam::UList::assign
void assign(const UList< T > &)
Assign elements to those from UList.
Definition: UList.C:37
AABBTree.H
Foam::processorPolyPatch
Neighbour processor patch.
Definition: processorPolyPatch.H:55
Foam::boundBox::min
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:54
Foam::meshToMesh::calcOverlappingProcs
label calcOverlappingProcs(const List< treeBoundBoxList > &procBb, const boundBox &bb, boolList &overlaps) const
Determine which processor bounding-boxes overlap.
Definition: meshToMeshParallelOps.C:87
Foam::PstreamBuffers::finishedSends
void finishedSends(const bool block=true)
Mark all sends as having been done. This will start receives.
Definition: PstreamBuffers.C:82
Foam::FatalError
error FatalError
processorPolyPatch.H
Foam::add
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
Definition: FieldFieldFunctions.C:870
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::HashTable
An STL-conforming hash table.
Definition: HashTable.H:61
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::autoPtr< Foam::mapDistribute >
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1004
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
Foam::AABBTree
Templated tree of axis-aligned bounding boxes (AABB)
Definition: AABBTree.H:57
Foam::List< label >
Foam::FixedList< label, 3 >
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
Foam::boundBox
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:55
mergePoints.H
Merge points. See below.
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::UIPstream
Input inter-processor communications stream operating on external buffer.
Definition: UIPstream.H:53
Foam::mapDistributeBase::constructMap
const labelListList & constructMap() const
From subsetted data to new reconstructed data.
Definition: mapDistributeBase.H:268
Foam::treeBoundBoxList
List< treeBoundBox > treeBoundBoxList
List of bounding boxes.
Definition: treeBoundBoxList.H:44
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
Foam::mergePoints
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.
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82