boundaryLayerCells.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | cfMesh: A library for mesh generation
4  \\ / O peration |
5  \\ / A nd | Author: Franjo Juretic (franjo.juretic@c-fields.com)
6  \\/ M anipulation | Copyright (C) Creative Fields, Ltd.
7 -------------------------------------------------------------------------------
8 License
9  This file is part of cfMesh.
10 
11  cfMesh is free software; you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by the
13  Free Software Foundation; either version 3 of the License, or (at your
14  option) any later version.
15 
16  cfMesh 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 cfMesh. If not, see <http://www.gnu.org/licenses/>.
23 
24 Description
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "boundaryLayers.H"
29 #include "meshSurfaceEngine.H"
30 #include "helperFunctions.H"
31 #include "helperFunctionsPar.H"
32 #include "demandDrivenData.H"
33 #include "VRWGraphList.H"
34 
35 #include "labelledPair.H"
36 #include "HashSet.H"
37 
38 #include <map>
39 
40 //#define DEBUGLayer
41 
42 # ifdef DEBUGLayer
43 #include "polyMeshGenAddressing.H"
44 # endif
45 
46 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
54 {
55  Info << "Starting creating layer cells" << endl;
56 
57  const meshSurfaceEngine& mse = surfaceEngine();
58  const faceList::subList& bFaces = mse.boundaryFaces();
59  const edgeList& edges = mse.edges();
60  const VRWGraph& faceEdges = mse.faceEdges();
61  const VRWGraph& edgeFaces = mse.edgeFaces();
62  const labelList& boundaryFacePatches = mse.boundaryFacePatches();
63  const labelList& faceOwners = mse.faceOwners();
64  const labelList& bp = mse.bp();
65  const VRWGraph& pointFaces = mse.pointFaces();
66 
68  const VRWGraph& pointPatches = mPart.pointPatches();
69 
70  //- mark patches which will be extruded into layer cells
71  boolList treatPatches(mesh_.boundaries().size(), false);
72  forAll(patchLabels, patchI)
73  {
74  const label pLabel = patchLabels[patchI];
75  forAll(treatPatchesWithPatch_[pLabel], i)
76  treatPatches[treatPatchesWithPatch_[pLabel][i]] = true;
77  }
78 
79  //- create new faces at parallel boundaries
80  const Map<label>* otherProcPatchPtr(NULL);
81  const Map<label>* otherFaceProcPtr(NULL);
82 
83  if( Pstream::parRun() )
84  {
85  createNewFacesParallel(treatPatches);
86 
87  otherProcPatchPtr = &mse.otherEdgeFacePatch();
88  otherFaceProcPtr = &mse.otherEdgeFaceAtProc();
89  }
90 
91  //- create lists for new boundary faces
92  VRWGraph newBoundaryFaces;
93  labelLongList newBoundaryOwners;
94  labelLongList newBoundaryPatches;
95 
96  //- create storage for new cells
97  VRWGraphList cellsToAdd;
98 
99  //- create layer cells and store boundary faces
100  const label nOldCells = mesh_.cells().size();
101  forAll(bFaces, bfI)
102  {
103  if( treatPatches[boundaryFacePatches[bfI]] )
104  {
105  const face& f = bFaces[bfI];
106  const label pKey = patchKey_[boundaryFacePatches[bfI]];
107 
108  DynList<DynList<label> > cellFaces;
109 
110  DynList<label> newF;
111 
112  //- store the current boundary face
113  newF.clear();
114  newF.append(f[0]);
115  for(label pI=f.size()-1;pI>0;--pI)
116  newF.append(f[pI]);
117  cellFaces.append(newF);
118 
119  //- create parallel face
120  forAll(f, pI)
121  newF[pI] = findNewNodeLabel(f[pI], pKey);
122 
123  cellFaces.append(newF);
124 
125  newBoundaryFaces.appendList(newF);
126  newBoundaryOwners.append(cellsToAdd.size() + nOldCells);
127  newBoundaryPatches.append(boundaryFacePatches[bfI]);
128 
129  //- create quad faces
130  newF.setSize(4);
131  forAll(f, pI)
132  {
133  newF[0] = f[pI];
134  newF[1] = f.nextLabel(pI);
135  newF[2] = findNewNodeLabel(newF[1], pKey);
136  newF[3] = findNewNodeLabel(f[pI], pKey);
137 
138  cellFaces.append(newF);
139 
140  //- check if the face is at the boundary
141  //- of the treated partitions
142  const label edgeI = faceEdges(bfI, pI);
143  if( edgeFaces.sizeOfRow(edgeI) == 2 )
144  {
145  label neiFace = edgeFaces(edgeI, 0);
146  if( neiFace == bfI )
147  neiFace = edgeFaces(edgeI, 1);
148 
149  if( !treatPatches[boundaryFacePatches[neiFace]] )
150  {
151  newBoundaryFaces.appendList(newF);
152  newBoundaryOwners.append(cellsToAdd.size() + nOldCells);
153  newBoundaryPatches.append(boundaryFacePatches[neiFace]);
154  }
155  }
156  else if( edgeFaces.sizeOfRow(edgeI) == 1 )
157  {
158  const Map<label>& otherProcPatch = *otherProcPatchPtr;
159  if( !treatPatches[otherProcPatch[edgeI]] )
160  {
161  //- face is a new boundary face
162  newBoundaryFaces.appendList(newF);
163  newBoundaryOwners.append(cellsToAdd.size() + nOldCells);
164  newBoundaryPatches.append(otherProcPatch[edgeI]);
165  }
166  }
167  }
168 
169  # ifdef DEBUGLayer
170  Info << "Adding cell " << cellFaces << endl;
171  # endif
172 
173  cellsToAdd.appendGraph(cellFaces);
174  }
175  else
176  {
177  # ifdef DEBUGLayer
178  Info << "Storing original boundary face "
179  << bfI << " into patch " << boundaryFacePatches[bfI] << endl;
180  # endif
181 
182  newBoundaryFaces.appendList(bFaces[bfI]);
183  newBoundaryOwners.append(faceOwners[bfI]);
184  newBoundaryPatches.append(boundaryFacePatches[bfI]);
185  }
186  }
187 
188  //- data for parallel execution
189  boolList procPoint;
190  LongList<DynList<label, 4> > pointProcFaces;
191  LongList<labelPair> faceAtPatches;
192  if( Pstream::parRun() )
193  {
194  procPoint.setSize(nPoints_);
195  procPoint = false;
196 
197  const Map<label>& globalToLocal = mse.globalToLocalBndPointAddressing();
198  const labelList& bPoints = mse.boundaryPoints();
199 
200  for
201  (
202  Map<label>::const_iterator iter=globalToLocal.begin();
203  iter!=globalToLocal.end();
204  ++iter
205  )
206  {
207  const label bpI = iter();
208  procPoint[bPoints[bpI]] = true;
209  }
210  }
211 
212  //- create cells at edges
213  forAll(edgeFaces, edgeI)
214  {
215  //- do not consider edges with no faces attached to it
216  if( edgeFaces.sizeOfRow(edgeI) == 0 )
217  continue;
218 
219  //- cells are generated at the processor with the lowest label
220  if(
221  (edgeFaces.sizeOfRow(edgeI) == 1) &&
222  (otherFaceProcPtr->operator[](edgeI) < Pstream::myProcNo())
223  )
224  continue;
225 
226  //- check if the edge is a feature edge
227  const label patchI = boundaryFacePatches[edgeFaces(edgeI, 0)];
228 
229  label patchJ;
230  if( otherProcPatchPtr && otherProcPatchPtr->found(edgeI) )
231  {
232  patchJ = otherProcPatchPtr->operator[](edgeI);
233  }
234  else
235  {
236  patchJ = boundaryFacePatches[edgeFaces(edgeI, 1)];
237  }
238 
239  if( patchI == patchJ )
240  continue;
241 
242  //- check if the faces attached to the edge have different keys
243  const label pKeyI = patchKey_[patchI];
244  const label pKeyJ = patchKey_[patchJ];
245 
246  if( pKeyI < 0 || pKeyJ < 0 )
247  {
248  continue;
250  (
251  "void boundaryLayers::createLayerCells(const labelList&)"
252  ) << "Patch key is negative at concave edge" << abort(FatalError);
253  }
254 
255  if( pKeyI == pKeyJ )
256  continue;
257 
258  const edge& e = edges[edgeI];
259  if( otherVrts_.find(e.start()) == otherVrts_.end() )
260  continue;
261  if( otherVrts_.find(e.end()) == otherVrts_.end() )
262  continue;
263 
264  //- generate faces of the bnd layer cell
265  FixedList<FixedList<label, 4>, 6> cellFaces;
266  createNewCellFromEdge(e, pKeyI, pKeyJ, cellFaces);
267 
268  //- store boundary faces
269  newBoundaryFaces.appendList(cellFaces[1]);
270  newBoundaryOwners.append(nOldCells+cellsToAdd.size());
271  newBoundaryPatches.append(patchJ);
272 
273  newBoundaryFaces.appendList(cellFaces[3]);
274  newBoundaryOwners.append(nOldCells+cellsToAdd.size());
275  newBoundaryPatches.append(patchI);
276 
277  //- check if face 5 is a boundary face or at an inter-processor boundary
278  const label bps = bp[e.start()];
279  label unusedPatch(-1);
280  forAllRow(pointPatches, bps, i)
281  {
282  const label ptchI = pointPatches(bps, i);
283 
284  if( ptchI == patchI )
285  continue;
286  if( ptchI == patchJ )
287  continue;
288  if( unusedPatch != -1 )
289  {
290  unusedPatch = -1;
291  break;
292  }
293 
294  unusedPatch = ptchI;
295  }
296 
297  if( unusedPatch != -1 && treatedPatch_[unusedPatch] )
298  {
299  //- add a face in the empty patch in case of 2D layer generation
300  newBoundaryFaces.appendList(cellFaces[5]);
301  newBoundaryOwners.append(nOldCells+cellsToAdd.size());
302  newBoundaryPatches.append(unusedPatch);
303  }
304  else if( Pstream::parRun() && procPoint[e.start()] )
305  {
306  //- add a face at inter-pocessor boundary
307  pointProcFaces.append(cellFaces[5]);
308  faceAtPatches.append(labelPair(patchI, patchJ));
309  }
310 
311  //- check if face 4 is a boundary face or at an inter-processor boundary
312  const label bpe = bp[e.end()];
313  unusedPatch = -1;
314  forAllRow(pointPatches, bpe, i)
315  {
316  const label ptchI = pointPatches(bpe, i);
317 
318  if( ptchI == patchI )
319  continue;
320  if( ptchI == patchJ )
321  continue;
322  if( unusedPatch != -1 )
323  {
324  unusedPatch = -1;
325  break;
326  }
327 
328  unusedPatch = ptchI;
329  }
330 
331  if( unusedPatch != -1 && treatedPatch_[unusedPatch] )
332  {
333  //- add a face in the empty patch in case of 2D layer generation
334  newBoundaryFaces.appendList(cellFaces[4]);
335  newBoundaryOwners.append(nOldCells+cellsToAdd.size());
336  newBoundaryPatches.append(unusedPatch);
337  }
338  else if( Pstream::parRun() && procPoint[e.end()] )
339  {
340  //- add a face at inter-pocessor boundary
341  pointProcFaces.append(cellFaces[4]);
342  faceAtPatches.append(labelPair(patchI, patchJ));
343  }
344 
345  # ifdef DEBUGLayer
346  Info << "Adding new cell at edge " << cellFaces << endl;
347  # endif
348 
349  //- append cell to the queue
350  cellsToAdd.appendGraph(cellFaces);
351  }
352 
353  //- create cells for corner nodes
354  typedef std::map<std::pair<label, label>, label> mPairToLabelType;
355  typedef std::map<label, mPairToLabelType> mPointsType;
356  typedef std::map<label, DynList<label, 3> > ppType;
357 
358  ppType nodePatches;
359  labelHashSet parPoint;
360 
361  if( Pstream::parRun() )
362  {
363  const labelList& bPoints = mse.boundaryPoints();
364  const VRWGraph& pProcs = mse.bpAtProcs();
365  const labelList& globalPointLabel = mse.globalBoundaryPointLabel();
366  const Map<label>& globalToLocal = mse.globalToLocalBndPointAddressing();
367 
368  std::map<label, labelLongList> facesToSend;
369 
370  typedef std::map<label, DynList<DynList<label, 8>, 8> > ppfType;
371 
372  ppfType parPointFaces;
373  ppType parPointPatches;
374 
375  forAllConstIter(mPointsType, otherVrts_, iter)
376  {
377  //- skip points on feature edges
378  if( iter->second.size() == 2 )
379  continue;
380 
381  const label bpI = bp[iter->first];
382 
383  if( pProcs.sizeOfRow(bpI) != 0 )
384  {
385  parPoint.insert(iter->first);
386 
387  //- point is at a parallel boundary
388  label pMin = pProcs(bpI, 0);
389  forAllRow(pProcs, bpI, i)
390  {
391  const label prI = pProcs(bpI, i);
392 
393  if( facesToSend.find(prI) == facesToSend.end() )
394  facesToSend.insert
395  (
396  std::make_pair(prI, labelLongList())
397  );
398 
399  if( prI < pMin )
400  pMin = prI;
401  }
402 
403  if( Pstream::myProcNo() == pMin )
404  {
405  DynList<label, 3>& pPatches = parPointPatches[bpI];
406  pPatches.setSize(pointFaces.sizeOfRow(bpI));
407 
408  DynList<DynList<label, 8>, 8>& pFaces = parPointFaces[bpI];
409  pFaces.setSize(pPatches.size());
410 
411  forAllRow(pointFaces, bpI, pfI)
412  {
413  const label bfI = pointFaces(bpI, pfI);
414  const face& bf = bFaces[bfI];
415 
416  pPatches[pfI] = boundaryFacePatches[bfI];
417 
418  DynList<label, 8>& bfCopy = pFaces[pfI];
419  bfCopy.setSize(bf.size());
420  forAll(bf, pI)
421  bfCopy[pI] = globalPointLabel[bp[bf[pI]]];
422  }
423 
424  continue;
425  }
426 
427  labelLongList& stp = facesToSend[pMin];
428 
429  //- send the data to the processor with the lowest label
430  //- data is flatenned as follows
431  //- 1. the number of faces and global point label
432  //- 2. number of points in the face
433  //- 3. patch label
434  //- 4. global labels of face points
435  stp.append(globalPointLabel[bpI]);
436  stp.append(pointFaces.sizeOfRow(bpI));
437  forAllRow(pointFaces, bpI, pfI)
438  {
439  const label bfI = pointFaces(bpI, pfI);
440  const face& bf = bFaces[bfI];
441 
442  stp.append(bf.size());
443  stp.append(boundaryFacePatches[bfI]);
444  forAll(bf, pI)
445  stp.append(globalPointLabel[bp[bf[pI]]]);
446  }
447  }
448  }
449 
450  //- exchange data with other processors
451  labelLongList receivedData;
452  help::exchangeMap(facesToSend, receivedData);
453 
454  label counter(0);
455  while( counter < receivedData.size() )
456  {
457  const label bpI = globalToLocal[receivedData[counter++]];
458  const label nFaces = receivedData[counter++];
459  for(label fI=0;fI<nFaces;++fI)
460  {
461  DynList<label, 8> f(receivedData[counter++]);
462  parPointPatches[bpI].append(receivedData[counter++]);
463  forAll(f, pI)
464  f[pI] = receivedData[counter++];
465  parPointFaces[bpI].append(f);
466  }
467  }
468 
469  //- sort faces sharing corners at the parallel boundaries
470  forAllIter(ppfType, parPointFaces, iter)
471  {
472  DynList<DynList<label, 8>, 8>& pFaces = iter->second;
473  DynList<label, 3>& fPatches = parPointPatches[iter->first];
474  const label gpI = globalPointLabel[iter->first];
475 
476  for(label i=0;i<pFaces.size();++i)
477  {
478  const DynList<label, 8>& bf = pFaces[i];
479  const label pos = bf.containsAtPosition(gpI);
480  const edge e(bf[pos], bf[bf.fcIndex(pos)]);
481 
482  for(label j=i+1;j<pFaces.size();++j)
483  {
484  const DynList<label, 8>& obf = pFaces[j];
485  if( obf.contains(e.start()) && obf.contains(e.end()) )
486  {
488  add = pFaces[i+1];
489  pFaces[i+1] = pFaces[j];
490  pFaces[j] = add;
491 
492  const label pAdd = fPatches[i+1];
493  fPatches[i+1] = fPatches[j];
494  fPatches[j] = pAdd;
495  break;
496  }
497  }
498  }
499 
500  DynList<label, 3> patchIDs;
501  forAll(fPatches, fpI)
502  patchIDs.appendIfNotIn(fPatches[fpI]);
503 
504  nodePatches.insert(std::make_pair(bPoints[iter->first], patchIDs));
505  }
506  }
507 
508  //- sort out point which are not at inter-processor boundaries
509  forAllConstIter(mPointsType, otherVrts_, iter)
510  {
511  if( iter->second.size() == 2 )
512  continue;
513 
514  if( parPoint.found(iter->first) )
515  continue;
516 
517  const label bpI = bp[iter->first];
518 
519  //- ensure correct orientation
520  DynList<label> pFaces(pointFaces.sizeOfRow(bpI));
521  forAll(pFaces, fI)
522  pFaces[fI] = pointFaces(bpI, fI);
523 
524  for(label i=0;i<pFaces.size();++i)
525  {
526  const face& bf = bFaces[pFaces[i]];
527  const edge e = bf.faceEdge(bf.which(iter->first));
528 
529  for(label j=i+1;j<pFaces.size();++j)
530  {
531  const face& obf = bFaces[pFaces[j]];
532  if(
533  (obf.which(e.start()) >= 0) &&
534  (obf.which(e.end()) >= 0)
535  )
536  {
537  const label add = pFaces[i+1];
538  pFaces[i+1] = pFaces[j];
539  pFaces[j] = add;
540  break;
541  }
542  }
543  }
544 
545  DynList<label, 3> patchIDs;
546  forAll(pFaces, patchI)
547  {
548  patchIDs.appendIfNotIn(boundaryFacePatches[pFaces[patchI]]);
549  }
550 
551  nodePatches.insert(std::make_pair(iter->first, patchIDs));
552  }
553 
554  //- create layer cells for corner nodes
555  forAllIter(ppType, nodePatches, iter)
556  {
557  const DynList<label, 3>& patchIDs = iter->second;
558  DynList<label, 3> pKeys;
559  forAll(patchIDs, patchI)
560  {
561  const label pKey = patchKey_[patchIDs[patchI]];
562 
563  if( pKey < 0 )
564  continue;
565 
566  pKeys.appendIfNotIn(pKey);
567  }
568 
569  if( pKeys.size() != 3 )
570  continue;
571 
572  # ifdef DEBUGLayer
573  Pout << "Creating corner cell at point " << iter->first << endl;
574  # endif
575 
576  FixedList<FixedList<label, 4>, 6> cellFaces;
577  createNewCellFromNode(iter->first, pKeys, cellFaces);
578 
579  //- store boundary faces
580  newBoundaryFaces.appendList(cellFaces[1]);
581  newBoundaryOwners.append(nOldCells+cellsToAdd.size());
582  newBoundaryPatches.append(patchIDs[0]);
583 
584  newBoundaryFaces.appendList(cellFaces[3]);
585  newBoundaryOwners.append(nOldCells+cellsToAdd.size());
586  newBoundaryPatches.append(patchIDs[1]);
587 
588  newBoundaryFaces.appendList(cellFaces[5]);
589  newBoundaryOwners.append(nOldCells+cellsToAdd.size());
590  newBoundaryPatches.append(patchIDs[2]);
591 
592  if( Pstream::parRun() )
593  {
594  if( procPoint[iter->first] )
595  {
596  pointProcFaces.append(cellFaces[0]);
597  faceAtPatches.append(labelPair(patchIDs[1], patchIDs[2]));
598 
599  pointProcFaces.append(cellFaces[2]);
600  faceAtPatches.append(labelPair(patchIDs[0], patchIDs[2]));
601 
602  pointProcFaces.append(cellFaces[4]);
603  faceAtPatches.append(labelPair(patchIDs[0], patchIDs[1]));
604  }
605  }
606 
607  # ifdef DEBUGLayer
608  Info << "Adding corner cell " << cellFaces << endl;
609  # endif
610 
611  //- append cell to the queue
612  cellsToAdd.appendGraph(cellFaces);
613  }
614 
615  if( Pstream::parRun() )
616  {
617  //- create faces at parallel boundaries created from
618  //- points at parallel boundaries
620  (
621  pointProcFaces,
622  faceAtPatches
623  );
624  }
625 
626  //- create mesh modifier
627  polyMeshGenModifier meshModifier(mesh_);
628 
629  meshModifier.addCells(cellsToAdd);
630 
631  cellsToAdd.clear();
632  meshModifier.reorderBoundaryFaces();
633  meshModifier.replaceBoundary
634  (
635  patchNames_,
636  newBoundaryFaces,
637  newBoundaryOwners,
638  newBoundaryPatches
639  );
640 
641  PtrList<boundaryPatch>& boundaries = meshModifier.boundariesAccess();
642  forAll(boundaries, patchI)
643  boundaries[patchI].patchType() = patchTypes_[patchI];
644 
645  //- delete meshSurfaceEngine
646  this->clearOut();
647 
648  Info << "Finished creating layer cells" << endl;
649 }
650 
652 (
653  const LongList<DynList<label, 4> >& faceCandidates,
654  const LongList<labelPair>& candidatePatches
655 )
656 {
657  const meshSurfaceEngine& mse = this->surfaceEngine();
658  const labelList& bPoints = mse.boundaryPoints();
659  const labelList& bp = mse.bp();
660  const VRWGraph& bpAtProcs = mse.bpAtProcs();
661  const labelList& globalPointLabel = mse.globalBoundaryPointLabel();
662  const Map<label>& globalToLocal = mse.globalToLocalBndPointAddressing();
663 
664  labelList otherFaceProc(faceCandidates.size(), -1);
665  //- some faces may appear more than once
666  //- such faces are ordinary internal faces
667  VRWGraph pointFaceCandidates(nPoints_);
668  forAll(faceCandidates, fI)
669  {
670  forAll(faceCandidates[fI], pI)
671  pointFaceCandidates.append(faceCandidates[fI][pI], fI);
672  }
673 
674  boolList duplicateFace(faceCandidates.size(), false);
675  List<labelledPair> pointOfOrigin(faceCandidates.size());
676  std::map<labelledPair, label> pointOfOriginToFaceLabel;
677  forAll(faceCandidates, fI)
678  {
679  const DynList<label, 4>& f = faceCandidates[fI];
680 
681  const label pointI = f[0];
682 
683  const labelledPair lp
684  (
685  globalPointLabel[bp[pointI]],
687  (
688  patchKey_[candidatePatches[fI][0]],
689  patchKey_[candidatePatches[fI][1]]
690  )
691  );
692 
693  if(
694  pointOfOriginToFaceLabel.find(lp) != pointOfOriginToFaceLabel.end()
695  )
696  {
697  duplicateFace[fI] = true;
698  pointOfOrigin[fI] = lp;
699  duplicateFace[pointOfOriginToFaceLabel[lp]] = true;
700  continue;
701  }
702 
703  pointOfOrigin[fI] = lp;
704 
705  pointOfOriginToFaceLabel.insert(std::make_pair(lp, fI));
706  }
707 
708  //- find the processor patch for each processor boundary face
709  //- the key of the algorithm is the point from which the face was created
710  //- by sending the point label and the associated patches, it will be
711  //- possible to find the other processor containing that face
712  std::map<label, LongList<labelledPair> > exchangeData;
713  const DynList<label>& neiProcs = mse.bpNeiProcs();
714  forAll(neiProcs, procI)
715  {
716  const label neiProcI = neiProcs[procI];
717 
718  if( neiProcI == Pstream::myProcNo() )
719  continue;
720 
721  if( exchangeData.find(neiProcI) == exchangeData.end() )
722  exchangeData.insert
723  (
724  std::make_pair(neiProcI, LongList<labelledPair>())
725  );
726  }
727 
728  forAll(faceCandidates, fI)
729  {
730  if( duplicateFace[fI] )
731  continue;
732 
733  const label bpI = bp[faceCandidates[fI][0]];
734 
735  forAllRow(bpAtProcs, bpI, procI)
736  {
737  const label neiProcNo = bpAtProcs(bpI, procI);
738  if( neiProcNo == Pstream::myProcNo() )
739  continue;
740 
741  LongList<labelledPair>& dataToSend = exchangeData[neiProcNo];
742  dataToSend.append(pointOfOrigin[fI]);
743  }
744  }
745 
746  //- exchange the data with other processors
747  std::map<label, List<labelledPair> > receivedMap;
748  help::exchangeMap(exchangeData, receivedMap);
749  exchangeData.clear();
750 
751  for
752  (
753  std::map<label, List<labelledPair> >::const_iterator
754  iter=receivedMap.begin();
755  iter!=receivedMap.end();
756  ++iter
757  )
758  {
759  const List<labelledPair>& receivedData = iter->second;
760 
761  forAll(receivedData, i)
762  {
763  const labelledPair& lpp = receivedData[i];
764  const label gpI = lpp.pairLabel();
765  const label pointI = bPoints[globalToLocal[gpI]];
766  const labelPair& lp = lpp.pair();
767 
768  forAllRow(pointFaceCandidates, pointI, i)
769  {
770  const label fI = pointFaceCandidates(pointI, i);
771  const DynList<label, 4>& f = faceCandidates[fI];
772 
773  const labelPair pk
774  (
775  patchKey_[candidatePatches[fI][0]],
776  patchKey_[candidatePatches[fI][1]]
777  );
778 
779  const labelPair rpk
780  (
781  patchKey_[candidatePatches[fI][1]],
782  patchKey_[candidatePatches[fI][0]]
783  );
784 
785  if(
786  (f[0] == pointI) && ((pk == lp) || (rpk == lp))
787  )
788  {
789  //- found the processor containing other face
790  otherFaceProc[pointOfOriginToFaceLabel[lpp]] = iter->first;
791  }
792  }
793  }
794  }
795  receivedMap.clear();
796 
797  //- sort the points in ascending order
798  //- this ensures the correct order of faces at the processor boundaries
799  sort(pointOfOrigin);
800 
801  Map<label> otherProcToProcPatch;
802  forAll(mesh_.procBoundaries(), patchI)
803  {
804  const processorBoundaryPatch& wp = mesh_.procBoundaries()[patchI];
805  otherProcToProcPatch.insert(wp.neiProcNo(), patchI);
806  }
807 
808  //- store processor faces
809  VRWGraph newProcFaces;
810  labelLongList newProc;
811 
812  forAll(pointOfOrigin, i)
813  {
814  const label fI = pointOfOriginToFaceLabel[pointOfOrigin[i]];
815 
816  if( duplicateFace[fI] || (otherFaceProc[fI] == -1) )
817  continue;
818 
819  if( !otherProcToProcPatch.found(otherFaceProc[fI]) )
820  {
821  otherProcToProcPatch.insert
822  (
823  otherFaceProc[fI],
824  polyMeshGenModifier(mesh_).addProcessorPatch
825  (
826  otherFaceProc[fI]
827  )
828  );
829  }
830 
831  newProcFaces.appendList(faceCandidates[fI]);
832  newProc.append(otherProcToProcPatch[otherFaceProc[fI]]);
833  }
834 
835  polyMeshGenModifier(mesh_).addProcessorFaces(newProcFaces, newProc);
836 }
837 
838 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
839 
840 } // End namespace Foam
841 
842 // ************************************************************************* //
Foam::meshSurfaceEngine::bpAtProcs
const VRWGraph & bpAtProcs() const
processors which contain the vertex
Definition: meshSurfaceEngineI.H:451
Foam::meshSurfaceEngine::bp
const labelList & bp() const
Definition: meshSurfaceEngineI.H:64
Foam::labelledPair::pair
const labelPair & pair() const
return the pair
Definition: labelledPair.H:91
Foam::LongList::append
void append(const T &e)
Append an element at the end of the list.
Definition: LongListI.H:265
Foam::boundaryLayers::patchNames_
wordList patchNames_
patch names
Definition: boundaryLayers.H:82
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::boundaryLayers::surfaceEngine
const meshSurfaceEngine & surfaceEngine() const
Return const reference to meshSurfaceEngine.
Definition: boundaryLayers.C:54
Foam::VRWGraphList::size
label size() const
Returns the number of graphs.
Definition: VRWGraphListI.H:96
helperFunctionsPar.H
Foam::boundaryLayers::nPoints_
label nPoints_
number of vertices in the mesh
Definition: boundaryLayers.H:105
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::meshSurfaceEngine::globalBoundaryPointLabel
const labelList & globalBoundaryPointLabel() const
global boundary point label
Definition: meshSurfaceEngineI.H:410
Foam::boundaryLayers::mesh_
polyMeshGen & mesh_
Reference to the mesh.
Definition: boundaryLayers.H:63
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::boundaryLayers::patchTypes_
wordList patchTypes_
patch types
Definition: boundaryLayers.H:85
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::face::faceEdge
edge faceEdge(const label n) const
Return n-th face edge.
Definition: faceI.H:110
Foam::polyMeshGenModifier::addProcessorFaces
void addProcessorFaces(const VRWGraph &procFaces, const labelLongList &facePatches)
add additional faces into processor patches
Definition: polyMeshGenModifierAddProcessorFaces.C:41
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::meshSurfaceEngine::bpNeiProcs
const DynList< label > & bpNeiProcs() const
communication matrix for sending point data
Definition: meshSurfaceEngineI.H:470
Foam::meshSurfaceEngine::pointFaces
const VRWGraph & pointFaces() const
Definition: meshSurfaceEngineI.H:162
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:58
pMin
dimensionedScalar pMin("pMin", dimPressure, fluid)
Foam::meshSurfaceEngine::boundaryFacePatches
const labelList & boundaryFacePatches() const
patch label for each boundary face
Definition: meshSurfaceEngineI.H:123
Foam::DynList::fcIndex
label fcIndex(const label index, const label offset=1) const
return forward and reverse circular indices
Definition: DynListI.H:486
Foam::meshSurfaceEngine::faceEdges
const VRWGraph & faceEdges() const
Definition: meshSurfaceEngineI.H:353
Foam::Map< label >
Foam::polyMeshGenModifier::replaceBoundary
void replaceBoundary(const wordList &patchNames, const VRWGraph &boundaryFaces, const labelLongList &faceOwners, const labelLongList &facePatches)
replace the boundary with new boundary faces
Definition: polyMeshGenModifierReplaceBoundary.C:42
Foam::boundaryLayers::otherVrts_
std::map< label, std::map< std::pair< label, label >, label > > otherVrts_
Definition: boundaryLayers.H:99
Foam::labelledPair
Definition: labelledPair.H:49
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::boundaryLayers::createNewFacesFromPointsParallel
void createNewFacesFromPointsParallel(const LongList< DynList< label, 4 > > &faceCandidates, const LongList< labelPair > &candidatePatches)
Definition: boundaryLayerCells.C:652
Foam::LongList::size
label size() const
Size of the active part of the list.
Definition: LongListI.H:203
Foam::boundaryLayers::createNewCellFromNode
void createNewCellFromNode(const label pointI, const DynList< label, 3 > &pKeys, FixedList< FixedList< label, 4 >, 6 > &cellFaces) const
creating hex cells near corners
Definition: boundaryLayersI.H:231
Foam::VRWGraphList
Definition: VRWGraphList.H:51
Foam::HashSet< label, Hash< label > >
Foam::DynList::contains
bool contains(const T &e) const
check if the element is in the list (takes linear time)
Definition: DynListI.H:339
Foam::meshSurfacePartitioner::pointPatches
const VRWGraph & pointPatches() const
Definition: meshSurfacePartitioner.H:123
Foam::VRWGraphList::appendGraph
void appendGraph(const GraphType &l)
Append a graph at the end of the graphList.
Definition: VRWGraphListI.H:123
Foam::face::which
label which(const label globalIndex) const
Navigation through face vertices.
Definition: face.C:630
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< SMALL) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &oldCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:235
Foam::polyMeshGenModifier::boundariesAccess
PtrList< boundaryPatch > & boundariesAccess()
access to boundary data
Definition: polyMeshGenModifier.H:131
Foam::LongList< label >
labelledPair.H
Foam::meshSurfaceEngine::boundaryFaces
const faceList::subList & boundaryFaces() const
Definition: meshSurfaceEngineI.H:103
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::polyMeshGenModifier::reorderBoundaryFaces
void reorderBoundaryFaces()
Definition: polyMeshGenModifierReorderBoundaryFaces.C:44
Foam::Info
messageStream Info
Foam::polyMeshGenCells::cells
const cellListPMG & cells() const
access to cells
Definition: polyMeshGenCellsI.H:39
Foam::meshSurfaceEngine::otherEdgeFaceAtProc
const Map< label > & otherEdgeFaceAtProc() const
Definition: meshSurfaceEngineI.H:568
Foam::DynList::appendIfNotIn
void appendIfNotIn(const T &e)
Definition: DynListI.H:324
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
forAllRow
#define forAllRow(graph, rowI, index)
Definition: VRWGraph.H:277
HashSet.H
Foam::VRWGraph::sizeOfRow
label sizeOfRow(const label rowI) const
Returns the number of elements in the given row.
Definition: VRWGraphI.H:127
Foam::boundaryLayers::patchKey_
labelList patchKey_
a key assigned to each patch. It is needed to search in otherVrts_
Definition: boundaryLayers.H:102
Foam::FatalError
error FatalError
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::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::boundaryLayers::createLayerCells
void createLayerCells(const labelList &patchLabels)
Definition: boundaryLayerCells.C:53
Foam::meshSurfaceEngine::edgeFaces
const VRWGraph & edgeFaces() const
Definition: meshSurfaceEngineI.H:334
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
meshSurfaceEngine.H
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::processorBoundaryPatch
Definition: processorBoundaryPatch.H:47
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::DynList
Definition: DynList.H:53
Foam::meshSurfacePartitioner
Definition: meshSurfacePartitioner.H:52
Foam::meshSurfaceEngine::edges
const edgeList & edges() const
Definition: meshSurfaceEngineI.H:296
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::boundaryLayers::clearOut
void clearOut()
delete meshSurfaceEngine
Definition: boundaryLayers.H:225
Foam::polyMeshGenModifier
Definition: polyMeshGenModifier.H:52
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::boundaryLayers::createNewCellFromEdge
void createNewCellFromEdge(const edge &e, const label pKeyI, const label pKeyJ, FixedList< FixedList< label, 4 >, 6 > &cellFaces) const
creating hex cells near feature edges
Definition: boundaryLayersI.H:81
Foam::DynList::containsAtPosition
label containsAtPosition(const T &e) const
Definition: DynListI.H:356
Foam::processorBoundaryPatch::neiProcNo
label neiProcNo() const
return the neighbour processor
Definition: processorBoundaryPatch.H:118
Foam::VRWGraphList::clear
void clear()
Clear the graph.
Definition: VRWGraphListI.H:115
Foam::labelledPair::pairLabel
label pairLabel() const
return pair label
Definition: labelledPair.H:85
Foam::Pair< label >
Foam::boundaryLayers::treatPatchesWithPatch_
List< DynList< label > > treatPatchesWithPatch_
extrude patches with patch
Definition: boundaryLayers.H:92
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
helperFunctions.H
f
labelList f(nPoints)
Foam::VRWGraph::appendList
void appendList(const ListType &l)
Append a list as a row at the end of the graph.
Definition: VRWGraphI.H:286
Foam::boundaryLayers::treatedPatch_
boolList treatedPatch_
Definition: boundaryLayers.H:89
Foam::polyMeshGenFaces::boundaries
const PtrList< boundaryPatch > & boundaries() const
ordinary boundaries
Definition: polyMeshGenFacesI.H:111
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::DynList::setSize
void setSize(const label)
Reset size of List.
Definition: DynListI.H:263
Foam::FixedList
A 1D vector of objects of type <T> with a fixed size <Size>.
Definition: FixedList.H:53
Foam::boundaryLayers::createNewFacesParallel
void createNewFacesParallel(const boolList &treatPatches)
Definition: boundaryLayersFacesAndCells.C:180
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: List.C:379
Foam::help::exchangeMap
void exchangeMap(const std::map< label, ListType > &m, LongList< T > &data, const Pstream::commsTypes commsType)
Definition: helperFunctionsPar.C:129
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::VRWGraph::append
void append(const label rowI, const label)
Append an element to the given row.
Definition: VRWGraphI.H:303
Foam::polyMeshGenModifier::addCells
void addCells(const LongList< faceList > &cellFaces)
add cells (vertices must be added)
Definition: polyMeshGenModifierAddCells.C:37
Foam::meshSurfaceEngine::boundaryPoints
const labelList & boundaryPoints() const
Definition: meshSurfaceEngineI.H:84
FatalErrorIn
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:313
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
Foam::DynList::size
label size() const
Definition: DynListI.H:235
Foam::meshSurfaceEngine::globalToLocalBndPointAddressing
const Map< label > & globalToLocalBndPointAddressing() const
global point label to local label. Only for processors points
Definition: meshSurfaceEngineI.H:431
Foam::sort
void sort(UList< T > &)
Definition: UList.C:107
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::meshSurfaceEngine::otherEdgeFacePatch
const Map< label > & otherEdgeFacePatch() const
Definition: meshSurfaceEngineI.H:588
Foam::VRWGraph
Definition: VRWGraph.H:101
VRWGraphList.H
Foam::boundaryLayers::surfacePartitioner
const meshSurfacePartitioner & surfacePartitioner() const
return const reference to meshSurfacePartitioner
Definition: boundaryLayers.C:62
Foam::labelLongList
LongList< label > labelLongList
Definition: labelLongList.H:46
polyMeshGenAddressing.H
Foam::labelPair
Pair< label > labelPair
Label pair.
Definition: labelPair.H:48
boundaryLayers.H
Foam::DynList::clear
void clear()
Clear the list, i.e. set next free to zero.
Definition: DynListI.H:279
Foam::meshSurfaceEngine
Definition: meshSurfaceEngine.H:54
Foam::cellListPMG::size
label size() const
return the number of used elements
Definition: cellListPMGI.H:56
Foam::boundaryLayers::findNewNodeLabel
label findNewNodeLabel(const label pointI, const label pKey) const
helper function finding a new face label for multiply extruded nodes
Definition: boundaryLayersI.H:38
Foam::DynList::append
void append(const T &e)
Append an element at the end of the list.
Definition: DynListI.H:304
Foam::meshSurfaceEngine::faceOwners
const labelList & faceOwners() const
Definition: meshSurfaceEngineI.H:143
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190