decompositionMethod.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 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 InClass
25  decompositionMethod
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "decompositionMethod.H"
30 #include "globalIndex.H"
31 #include "syncTools.H"
32 #include "Tuple2.H"
33 #include "faceSet.H"
34 #include "regionSplit.H"
35 #include "localPointRegion.H"
36 #include "minData.H"
37 #include "FaceCellWave.H"
38 
43 
44 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
45 
46 namespace Foam
47 {
48  defineTypeNameAndDebug(decompositionMethod, 0);
49  defineRunTimeSelectionTable(decompositionMethod, dictionary);
50 }
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 Foam::decompositionMethod::decompositionMethod
55 (
56  const dictionary& decompositionDict
57 )
58 :
59  decompositionDict_(decompositionDict),
60  nProcessors_
61  (
62  readLabel(decompositionDict.lookup("numberOfSubdomains"))
63  )
64 {
65  // Read any constraints
66  wordList constraintTypes_;
67  if (decompositionDict_.found("constraints"))
68  {
69  // PtrList<dictionary> constraintsList
70  //(
71  // decompositionDict_.lookup("constraints")
72  //);
73  // forAll(constraintsList, i)
74  //{
75  // const dictionary& dict = constraintsList[i];
76  const dictionary& constraintsList = decompositionDict_.subDict
77  (
78  "constraints"
79  );
80  forAllConstIter(dictionary, constraintsList, iter)
81  {
82  const dictionary& dict = iter().dict();
83 
84  constraintTypes_.append(dict.lookup("type"));
85 
86  constraints_.append
87  (
89  (
90  dict,
91  constraintTypes_.last()
92  )
93  );
94  }
95  }
96 
97  // Backwards compatibility
98  if
99  (
100  decompositionDict_.found("preserveBaffles")
101  && findIndex
102  (
103  constraintTypes_,
104  decompositionConstraints::preserveBafflesConstraint::typeName
105  ) == -1
106  )
107  {
108  constraints_.append
109  (
111  );
112  }
113 
114  if
115  (
116  decompositionDict_.found("preservePatches")
117  && findIndex
118  (
119  constraintTypes_,
120  decompositionConstraints::preservePatchesConstraint::typeName
121  ) == -1
122  )
123  {
124  const wordReList pNames(decompositionDict_.lookup("preservePatches"));
125 
126  constraints_.append
127  (
129  );
130  }
131 
132  if
133  (
134  decompositionDict_.found("preserveFaceZones")
135  && findIndex
136  (
137  constraintTypes_,
138  decompositionConstraints::preserveFaceZonesConstraint::typeName
139  ) == -1
140  )
141  {
142  const wordReList zNames(decompositionDict_.lookup("preserveFaceZones"));
143 
144  constraints_.append
145  (
147  );
148  }
149 
150  if
151  (
152  decompositionDict_.found("singleProcessorFaceSets")
153  && findIndex
154  (
155  constraintTypes_,
156  decompositionConstraints::preserveFaceZonesConstraint::typeName
157  ) == -1
158  )
159  {
160  const List<Tuple2<word, label>> zNameAndProcs
161  (
162  decompositionDict_.lookup("singleProcessorFaceSets")
163  );
164 
165  constraints_.append
166  (
168  (
169  zNameAndProcs
170  )
171  );
172  }
173 }
174 
175 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
176 
178 (
179  const dictionary& decompositionDict
180 )
181 {
182  word methodType(decompositionDict.lookup("method"));
183 
184  Info<< "Selecting decompositionMethod " << methodType << endl;
185 
186  dictionaryConstructorTable::iterator cstrIter =
187  dictionaryConstructorTablePtr_->find(methodType);
188 
189  if (cstrIter == dictionaryConstructorTablePtr_->end())
190  {
192  << "Unknown decompositionMethod "
193  << methodType << nl << nl
194  << "Valid decompositionMethods are : " << endl
195  << dictionaryConstructorTablePtr_->sortedToc()
196  << exit(FatalError);
197  }
198 
199  return autoPtr<decompositionMethod>(cstrIter()(decompositionDict));
200 }
201 
202 
204 (
205  const polyMesh& mesh,
206  const pointField& points
207 )
208 {
209  scalarField weights(points.size(), 1.0);
210 
211  return decompose(mesh, points, weights);
212 }
213 
214 
216 (
217  const polyMesh& mesh,
218  const labelList& fineToCoarse,
219  const pointField& coarsePoints,
220  const scalarField& coarseWeights
221 )
222 {
223  CompactListList<label> coarseCellCells;
224  calcCellCells
225  (
226  mesh,
227  fineToCoarse,
228  coarsePoints.size(),
229  true, // use global cell labels
230  coarseCellCells
231  );
232 
233  // Decompose based on agglomerated points
234  labelList coarseDistribution
235  (
236  decompose
237  (
238  coarseCellCells(),
239  coarsePoints,
240  coarseWeights
241  )
242  );
243 
244  // Rework back into decomposition for original mesh_
245  labelList fineDistribution(fineToCoarse.size());
246 
247  forAll(fineDistribution, i)
248  {
249  fineDistribution[i] = coarseDistribution[fineToCoarse[i]];
250  }
251 
252  return fineDistribution;
253 }
254 
255 
257 (
258  const polyMesh& mesh,
259  const labelList& fineToCoarse,
260  const pointField& coarsePoints
261 )
262 {
263  scalarField cWeights(coarsePoints.size(), 1.0);
264 
265  return decompose
266  (
267  mesh,
268  fineToCoarse,
269  coarsePoints,
270  cWeights
271  );
272 }
273 
274 
276 (
277  const labelListList& globalCellCells,
278  const pointField& cc
279 )
280 {
281  scalarField cWeights(cc.size(), 1.0);
282 
283  return decompose(globalCellCells, cc, cWeights);
284 }
285 
286 
288 (
289  const polyMesh& mesh,
290  const labelList& agglom,
291  const label nLocalCoarse,
292  const bool parallel,
293  CompactListList<label>& cellCells
294 )
295 {
296  const labelList& faceOwner = mesh.faceOwner();
297  const labelList& faceNeighbour = mesh.faceNeighbour();
298  const polyBoundaryMesh& patches = mesh.boundaryMesh();
299 
300 
301  // Create global cell numbers
302  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
303 
304  globalIndex globalAgglom
305  (
306  nLocalCoarse,
309  parallel
310  );
311 
312 
313  // Get agglomerate owner on other side of coupled faces
314  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
315 
316  labelList globalNeighbour(mesh.nFaces()-mesh.nInternalFaces());
317 
318  forAll(patches, patchi)
319  {
320  const polyPatch& pp = patches[patchi];
321 
322  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
323  {
324  label facei = pp.start();
325  label bFacei = pp.start() - mesh.nInternalFaces();
326 
327  forAll(pp, i)
328  {
329  globalNeighbour[bFacei] = globalAgglom.toGlobal
330  (
331  agglom[faceOwner[facei]]
332  );
333 
334  bFacei++;
335  facei++;
336  }
337  }
338  }
339 
340  // Get the cell on the other side of coupled patches
341  syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
342 
343 
344  // Count number of faces (internal + coupled)
345  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
346 
347  // Number of faces per coarse cell
348  labelList nFacesPerCell(nLocalCoarse, 0);
349 
350  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
351  {
352  label own = agglom[faceOwner[facei]];
353  label nei = agglom[faceNeighbour[facei]];
354 
355  nFacesPerCell[own]++;
356  nFacesPerCell[nei]++;
357  }
358 
359  forAll(patches, patchi)
360  {
361  const polyPatch& pp = patches[patchi];
362 
363  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
364  {
365  label facei = pp.start();
366  label bFacei = pp.start()-mesh.nInternalFaces();
367 
368  forAll(pp, i)
369  {
370  label own = agglom[faceOwner[facei]];
371 
372  label globalNei = globalNeighbour[bFacei];
373  if
374  (
375  !globalAgglom.isLocal(globalNei)
376  || globalAgglom.toLocal(globalNei) != own
377  )
378  {
379  nFacesPerCell[own]++;
380  }
381 
382  facei++;
383  bFacei++;
384  }
385  }
386  }
387 
388 
389  // Fill in offset and data
390  // ~~~~~~~~~~~~~~~~~~~~~~~
391 
392  cellCells.setSize(nFacesPerCell);
393 
394  nFacesPerCell = 0;
395 
396  labelList& m = cellCells.m();
397  const labelList& offsets = cellCells.offsets();
398 
399  // For internal faces is just offsetted owner and neighbour
400  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
401  {
402  label own = agglom[faceOwner[facei]];
403  label nei = agglom[faceNeighbour[facei]];
404 
405  m[offsets[own] + nFacesPerCell[own]++] = globalAgglom.toGlobal(nei);
406  m[offsets[nei] + nFacesPerCell[nei]++] = globalAgglom.toGlobal(own);
407  }
408 
409  // For boundary faces is offsetted coupled neighbour
410  forAll(patches, patchi)
411  {
412  const polyPatch& pp = patches[patchi];
413 
414  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
415  {
416  label facei = pp.start();
417  label bFacei = pp.start()-mesh.nInternalFaces();
418 
419  forAll(pp, i)
420  {
421  label own = agglom[faceOwner[facei]];
422 
423  label globalNei = globalNeighbour[bFacei];
424 
425  if
426  (
427  !globalAgglom.isLocal(globalNei)
428  || globalAgglom.toLocal(globalNei) != own
429  )
430  {
431  m[offsets[own] + nFacesPerCell[own]++] = globalNei;
432  }
433 
434  facei++;
435  bFacei++;
436  }
437  }
438  }
439 
440 
441  // Check for duplicates connections between cells
442  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
443  // Done as postprocessing step since we now have cellCells.
444  label newIndex = 0;
445  labelHashSet nbrCells;
446 
447 
448  if (cellCells.size() == 0)
449  {
450  return;
451  }
452 
453  label startIndex = cellCells.offsets()[0];
454 
455  forAll(cellCells, celli)
456  {
457  nbrCells.clear();
458  nbrCells.insert(globalAgglom.toGlobal(celli));
459 
460  label endIndex = cellCells.offsets()[celli+1];
461 
462  for (label i = startIndex; i < endIndex; i++)
463  {
464  if (nbrCells.insert(cellCells.m()[i]))
465  {
466  cellCells.m()[newIndex++] = cellCells.m()[i];
467  }
468  }
469  startIndex = endIndex;
470  cellCells.offsets()[celli+1] = newIndex;
471  }
472 
473  cellCells.m().setSize(newIndex);
474 
475  // forAll(cellCells, celli)
476  //{
477  // Pout<< "Original: Coarse cell " << celli << endl;
478  // forAll(mesh.cellCells()[celli], i)
479  // {
480  // Pout<< " nbr:" << mesh.cellCells()[celli][i] << endl;
481  // }
482  // Pout<< "Compacted: Coarse cell " << celli << endl;
483  // const labelUList cCells = cellCells[celli];
484  // forAll(cCells, i)
485  // {
486  // Pout<< " nbr:" << cCells[i] << endl;
487  // }
488  //}
489 }
490 
491 
493 (
494  const polyMesh& mesh,
495  const labelList& agglom,
496  const label nLocalCoarse,
497  const bool parallel,
498  CompactListList<label>& cellCells,
499  CompactListList<scalar>& cellCellWeights
500 )
501 {
502  const labelList& faceOwner = mesh.faceOwner();
503  const labelList& faceNeighbour = mesh.faceNeighbour();
504  const polyBoundaryMesh& patches = mesh.boundaryMesh();
505 
506 
507  // Create global cell numbers
508  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
509 
510  globalIndex globalAgglom
511  (
512  nLocalCoarse,
515  parallel
516  );
517 
518 
519  // Get agglomerate owner on other side of coupled faces
520  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
521 
522  labelList globalNeighbour(mesh.nFaces()-mesh.nInternalFaces());
523 
524  forAll(patches, patchi)
525  {
526  const polyPatch& pp = patches[patchi];
527 
528  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
529  {
530  label faceI = pp.start();
531  label bFaceI = pp.start() - mesh.nInternalFaces();
532 
533  forAll(pp, i)
534  {
535  globalNeighbour[bFaceI] = globalAgglom.toGlobal
536  (
537  agglom[faceOwner[faceI]]
538  );
539 
540  bFaceI++;
541  faceI++;
542  }
543  }
544  }
545 
546  // Get the cell on the other side of coupled patches
547  syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
548 
549 
550  // Count number of faces (internal + coupled)
551  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
552 
553  // Number of faces per coarse cell
554  labelList nFacesPerCell(nLocalCoarse, 0);
555 
556  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
557  {
558  label own = agglom[faceOwner[faceI]];
559  label nei = agglom[faceNeighbour[faceI]];
560 
561  nFacesPerCell[own]++;
562  nFacesPerCell[nei]++;
563  }
564 
565  forAll(patches, patchi)
566  {
567  const polyPatch& pp = patches[patchi];
568 
569  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
570  {
571  label faceI = pp.start();
572  label bFaceI = pp.start()-mesh.nInternalFaces();
573 
574  forAll(pp, i)
575  {
576  label own = agglom[faceOwner[faceI]];
577 
578  label globalNei = globalNeighbour[bFaceI];
579  if
580  (
581  !globalAgglom.isLocal(globalNei)
582  || globalAgglom.toLocal(globalNei) != own
583  )
584  {
585  nFacesPerCell[own]++;
586  }
587 
588  faceI++;
589  bFaceI++;
590  }
591  }
592  }
593 
594 
595  // Fill in offset and data
596  // ~~~~~~~~~~~~~~~~~~~~~~~
597 
598  cellCells.setSize(nFacesPerCell);
599  cellCellWeights.setSize(nFacesPerCell);
600 
601  nFacesPerCell = 0;
602 
603  labelList& m = cellCells.m();
604  scalarList& w = cellCellWeights.m();
605  const labelList& offsets = cellCells.offsets();
606 
607  // For internal faces is just offsetted owner and neighbour
608  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
609  {
610  label own = agglom[faceOwner[faceI]];
611  label nei = agglom[faceNeighbour[faceI]];
612 
613  label ownIndex = offsets[own] + nFacesPerCell[own]++;
614  label neiIndex = offsets[nei] + nFacesPerCell[nei]++;
615 
616  m[ownIndex] = globalAgglom.toGlobal(nei);
617  w[ownIndex] = mag(mesh.faceAreas()[faceI]);
618  m[neiIndex] = globalAgglom.toGlobal(own);
619  w[ownIndex] = mag(mesh.faceAreas()[faceI]);
620  }
621 
622  // For boundary faces is offsetted coupled neighbour
623  forAll(patches, patchi)
624  {
625  const polyPatch& pp = patches[patchi];
626 
627  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
628  {
629  label faceI = pp.start();
630  label bFaceI = pp.start()-mesh.nInternalFaces();
631 
632  forAll(pp, i)
633  {
634  label own = agglom[faceOwner[faceI]];
635 
636  label globalNei = globalNeighbour[bFaceI];
637 
638  if
639  (
640  !globalAgglom.isLocal(globalNei)
641  || globalAgglom.toLocal(globalNei) != own
642  )
643  {
644  label ownIndex = offsets[own] + nFacesPerCell[own]++;
645  m[ownIndex] = globalNei;
646  w[ownIndex] = mag(mesh.faceAreas()[faceI]);
647  }
648 
649  faceI++;
650  bFaceI++;
651  }
652  }
653  }
654 
655 
656  // Check for duplicates connections between cells
657  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
658  // Done as postprocessing step since we now have cellCells.
659  label newIndex = 0;
660  labelHashSet nbrCells;
661 
662 
663  if (cellCells.size() == 0)
664  {
665  return;
666  }
667 
668  label startIndex = cellCells.offsets()[0];
669 
670  forAll(cellCells, cellI)
671  {
672  nbrCells.clear();
673  nbrCells.insert(globalAgglom.toGlobal(cellI));
674 
675  label endIndex = cellCells.offsets()[cellI+1];
676 
677  for (label i = startIndex; i < endIndex; i++)
678  {
679  if (nbrCells.insert(cellCells.m()[i]))
680  {
681  cellCells.m()[newIndex] = cellCells.m()[i];
682  cellCellWeights.m()[newIndex] = cellCellWeights.m()[i];
683  newIndex++;
684  }
685  }
686  startIndex = endIndex;
687  cellCells.offsets()[cellI+1] = newIndex;
688  cellCellWeights.offsets()[cellI+1] = newIndex;
689  }
690 
691  cellCells.m().setSize(newIndex);
692  cellCellWeights.m().setSize(newIndex);
693 }
694 
695 
696 //void Foam::decompositionMethod::calcCellCells
697 //(
698 // const polyMesh& mesh,
699 // const boolList& blockedFace,
700 // const List<labelPair>& explicitConnections,
701 // const labelList& agglom,
702 // const label nLocalCoarse,
703 // const bool parallel,
704 // CompactListList<label>& cellCells
705 //)
706 //{
707 // const labelList& faceOwner = mesh.faceOwner();
708 // const labelList& faceNeighbour = mesh.faceNeighbour();
709 // const polyBoundaryMesh& patches = mesh.boundaryMesh();
710 //
711 //
712 // // Create global cell numbers
713 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~
714 //
715 // globalIndex globalAgglom
716 // (
717 // nLocalCoarse,
718 // Pstream::msgType(),
719 // Pstream::worldComm,
720 // parallel
721 // );
722 //
723 //
724 // // Get agglomerate owner on other side of coupled faces
725 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
726 //
727 // labelList globalNeighbour(mesh.nFaces()-mesh.nInternalFaces());
728 //
729 // forAll(patches, patchi)
730 // {
731 // const polyPatch& pp = patches[patchi];
732 //
733 // if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
734 // {
735 // label facei = pp.start();
736 // label bFacei = pp.start() - mesh.nInternalFaces();
737 //
738 // forAll(pp, i)
739 // {
740 // globalNeighbour[bFacei] = globalAgglom.toGlobal
741 // (
742 // agglom[faceOwner[facei]]
743 // );
744 //
745 // bFacei++;
746 // facei++;
747 // }
748 // }
749 // }
750 //
751 // // Get the cell on the other side of coupled patches
752 // syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
753 //
754 //
755 // // Count number of faces (internal + coupled)
756 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
757 //
758 // // Number of faces per coarse cell
759 // labelList nFacesPerCell(nLocalCoarse, 0);
760 //
761 // // 1. Internal faces
762 // for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
763 // {
764 // if (!blockedFace[facei])
765 // {
766 // label own = agglom[faceOwner[facei]];
767 // label nei = agglom[faceNeighbour[facei]];
768 //
769 // nFacesPerCell[own]++;
770 // nFacesPerCell[nei]++;
771 // }
772 // }
773 //
774 // // 2. Coupled faces
775 // forAll(patches, patchi)
776 // {
777 // const polyPatch& pp = patches[patchi];
778 //
779 // if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
780 // {
781 // label facei = pp.start();
782 // label bFacei = pp.start()-mesh.nInternalFaces();
783 //
784 // forAll(pp, i)
785 // {
786 // if (!blockedFace[facei])
787 // {
788 // label own = agglom[faceOwner[facei]];
789 //
790 // label globalNei = globalNeighbour[bFacei];
791 // if
792 // (
793 // !globalAgglom.isLocal(globalNei)
794 // || globalAgglom.toLocal(globalNei) != own
795 // )
796 // {
797 // nFacesPerCell[own]++;
798 // }
799 //
800 // facei++;
801 // bFacei++;
802 // }
803 // }
804 // }
805 // }
806 //
807 // // 3. Explicit connections between non-coupled boundary faces
808 // forAll(explicitConnections, i)
809 // {
810 // const labelPair& baffle = explicitConnections[i];
811 // label f0 = baffle.first();
812 // label f1 = baffle.second();
813 //
814 // if (!blockedFace[f0] && blockedFace[f1])
815 // {
816 // label f0Own = agglom[faceOwner[f0]];
817 // label f1Own = agglom[faceOwner[f1]];
818 //
819 // // Always count the connection between the two owner sides
820 // if (f0Own != f1Own)
821 // {
822 // nFacesPerCell[f0Own]++;
823 // nFacesPerCell[f1Own]++;
824 // }
825 //
826 // // Add any neighbour side connections
827 // if (mesh.isInternalFace(f0))
828 // {
829 // label f0Nei = agglom[faceNeighbour[f0]];
830 //
831 // if (mesh.isInternalFace(f1))
832 // {
833 // // Internal faces
834 // label f1Nei = agglom[faceNeighbour[f1]];
835 //
836 // if (f0Own != f1Nei)
837 // {
838 // nFacesPerCell[f0Own]++;
839 // nFacesPerCell[f1Nei]++;
840 // }
841 // if (f0Nei != f1Own)
842 // {
843 // nFacesPerCell[f0Nei]++;
844 // nFacesPerCell[f1Own]++;
845 // }
846 // if (f0Nei != f1Nei)
847 // {
848 // nFacesPerCell[f0Nei]++;
849 // nFacesPerCell[f1Nei]++;
850 // }
851 // }
852 // else
853 // {
854 // // f1 boundary face
855 // if (f0Nei != f1Own)
856 // {
857 // nFacesPerCell[f0Nei]++;
858 // nFacesPerCell[f1Own]++;
859 // }
860 // }
861 // }
862 // else
863 // {
864 // if (mesh.isInternalFace(f1))
865 // {
866 // label f1Nei = agglom[faceNeighbour[f1]];
867 // if (f0Own != f1Nei)
868 // {
869 // nFacesPerCell[f0Own]++;
870 // nFacesPerCell[f1Nei]++;
871 // }
872 // }
873 // }
874 // }
875 // }
876 //
877 //
878 // // Fill in offset and data
879 // // ~~~~~~~~~~~~~~~~~~~~~~~
880 //
881 // cellCells.setSize(nFacesPerCell);
882 //
883 // nFacesPerCell = 0;
884 //
885 // labelList& m = cellCells.m();
886 // const labelList& offsets = cellCells.offsets();
887 //
888 // // 1. For internal faces is just offsetted owner and neighbour
889 // for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
890 // {
891 // if (!blockedFace[facei])
892 // {
893 // label own = agglom[faceOwner[facei]];
894 // label nei = agglom[faceNeighbour[facei]];
895 //
896 // m[offsets[own] + nFacesPerCell[own]++] =
897 // globalAgglom.toGlobal(nei);
898 // m[offsets[nei] + nFacesPerCell[nei]++] =
899 // globalAgglom.toGlobal(own);
900 // }
901 // }
902 //
903 // // 2. For boundary faces is offsetted coupled neighbour
904 // forAll(patches, patchi)
905 // {
906 // const polyPatch& pp = patches[patchi];
907 //
908 // if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
909 // {
910 // label facei = pp.start();
911 // label bFacei = pp.start()-mesh.nInternalFaces();
912 //
913 // forAll(pp, i)
914 // {
915 // if (!blockedFace[facei])
916 // {
917 // label own = agglom[faceOwner[facei]];
918 //
919 // label globalNei = globalNeighbour[bFacei];
920 //
921 // if
922 // (
923 // !globalAgglom.isLocal(globalNei)
924 // || globalAgglom.toLocal(globalNei) != own
925 // )
926 // {
927 // m[offsets[own] + nFacesPerCell[own]++] = globalNei;
928 // }
929 //
930 // facei++;
931 // bFacei++;
932 // }
933 // }
934 // }
935 // }
936 //
937 // // 3. Explicit connections between non-coupled boundary faces
938 // forAll(explicitConnections, i)
939 // {
940 // const labelPair& baffle = explicitConnections[i];
941 // label f0 = baffle.first();
942 // label f1 = baffle.second();
943 //
944 // if (!blockedFace[f0] && blockedFace[f1])
945 // {
946 // label f0Own = agglom[faceOwner[f0]];
947 // label f1Own = agglom[faceOwner[f1]];
948 //
949 // // Always count the connection between the two owner sides
950 // if (f0Own != f1Own)
951 // {
952 // m[offsets[f0Own] + nFacesPerCell[f0Own]++] =
953 // globalAgglom.toGlobal(f1Own);
954 // m[offsets[f1Own] + nFacesPerCell[f1Own]++] =
955 // globalAgglom.toGlobal(f0Own);
956 // }
957 //
958 // // Add any neighbour side connections
959 // if (mesh.isInternalFace(f0))
960 // {
961 // label f0Nei = agglom[faceNeighbour[f0]];
962 //
963 // if (mesh.isInternalFace(f1))
964 // {
965 // // Internal faces
966 // label f1Nei = agglom[faceNeighbour[f1]];
967 //
968 // if (f0Own != f1Nei)
969 // {
970 // m[offsets[f0Own] + nFacesPerCell[f0Own]++] =
971 // globalAgglom.toGlobal(f1Nei);
972 // m[offsets[f1Nei] + nFacesPerCell[f1Nei]++] =
973 // globalAgglom.toGlobal(f1Nei);
974 // }
975 // if (f0Nei != f1Own)
976 // {
977 // m[offsets[f0Nei] + nFacesPerCell[f0Nei]++] =
978 // globalAgglom.toGlobal(f1Own);
979 // m[offsets[f1Own] + nFacesPerCell[f1Own]++] =
980 // globalAgglom.toGlobal(f0Nei);
981 // }
982 // if (f0Nei != f1Nei)
983 // {
984 // m[offsets[f0Nei] + nFacesPerCell[f0Nei]++] =
985 // globalAgglom.toGlobal(f1Nei);
986 // m[offsets[f1Nei] + nFacesPerCell[f1Nei]++] =
987 // globalAgglom.toGlobal(f0Nei);
988 // }
989 // }
990 // else
991 // {
992 // // f1 boundary face
993 // if (f0Nei != f1Own)
994 // {
995 // m[offsets[f0Nei] + nFacesPerCell[f0Nei]++] =
996 // globalAgglom.toGlobal(f1Own);
997 // m[offsets[f1Own] + nFacesPerCell[f1Own]++] =
998 // globalAgglom.toGlobal(f0Nei);
999 // }
1000 // }
1001 // }
1002 // else
1003 // {
1004 // if (mesh.isInternalFace(f1))
1005 // {
1006 // label f1Nei = agglom[faceNeighbour[f1]];
1007 // if (f0Own != f1Nei)
1008 // {
1009 // m[offsets[f0Own] + nFacesPerCell[f0Own]++] =
1010 // globalAgglom.toGlobal(f1Nei);
1011 // m[offsets[f1Nei] + nFacesPerCell[f1Nei]++] =
1012 // globalAgglom.toGlobal(f0Own);
1013 // }
1014 // }
1015 // }
1016 // }
1017 // }
1018 //
1019 //
1020 // // Check for duplicates connections between cells
1021 // // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1022 // // Done as postprocessing step since we now have cellCells.
1023 // label newIndex = 0;
1024 // labelHashSet nbrCells;
1025 //
1026 //
1027 // if (cellCells.size() == 0)
1028 // {
1029 // return;
1030 // }
1031 //
1032 // label startIndex = cellCells.offsets()[0];
1033 //
1034 // forAll(cellCells, celli)
1035 // {
1036 // nbrCells.clear();
1037 // nbrCells.insert(globalAgglom.toGlobal(celli));
1038 //
1039 // label endIndex = cellCells.offsets()[celli+1];
1040 //
1041 // for (label i = startIndex; i < endIndex; i++)
1042 // {
1043 // if (nbrCells.insert(cellCells.m()[i]))
1044 // {
1045 // cellCells.m()[newIndex++] = cellCells.m()[i];
1046 // }
1047 // }
1048 // startIndex = endIndex;
1049 // cellCells.offsets()[celli+1] = newIndex;
1050 // }
1051 //
1052 // cellCells.m().setSize(newIndex);
1053 //
1054 // // forAll(cellCells, celli)
1055 // //{
1056 // // Pout<< "Original: Coarse cell " << celli << endl;
1057 // // forAll(mesh.cellCells()[celli], i)
1058 // // {
1059 // // Pout<< " nbr:" << mesh.cellCells()[celli][i] << endl;
1060 // // }
1061 // // Pout<< "Compacted: Coarse cell " << celli << endl;
1062 // // const labelUList cCells = cellCells[celli];
1063 // // forAll(cCells, i)
1064 // // {
1065 // // Pout<< " nbr:" << cCells[i] << endl;
1066 // // }
1067 // //}
1068 //}
1069 
1070 
1073  const polyMesh& mesh,
1074  const scalarField& cellWeights,
1075 
1076  //- Whether owner and neighbour should be on same processor
1077  // (takes priority over explicitConnections)
1078  const boolList& blockedFace,
1079 
1080  //- Whether whole sets of faces (and point neighbours) need to be kept
1081  // on single processor
1082  const PtrList<labelList>& specifiedProcessorFaces,
1083  const labelList& specifiedProcessor,
1084 
1085  //- Additional connections between boundary faces
1086  const List<labelPair>& explicitConnections
1087 )
1088 {
1089  // Any weights specified?
1090  label nWeights = returnReduce(cellWeights.size(), sumOp<label>());
1091 
1092  if (nWeights > 0 && cellWeights.size() != mesh.nCells())
1093  {
1095  << "Number of weights " << cellWeights.size()
1096  << " differs from number of cells " << mesh.nCells()
1097  << exit(FatalError);
1098  }
1099 
1100 
1101  // Any processor sets?
1102  label nProcSets = 0;
1103  forAll(specifiedProcessorFaces, setI)
1104  {
1105  nProcSets += specifiedProcessorFaces[setI].size();
1106  }
1107  reduce(nProcSets, sumOp<label>());
1108 
1109  // Any non-mesh connections?
1110  label nConnections = returnReduce
1111  (
1112  explicitConnections.size(),
1113  sumOp<label>()
1114  );
1115 
1116  // Any faces not blocked?
1117  label nUnblocked = 0;
1118  forAll(blockedFace, facei)
1119  {
1120  if (!blockedFace[facei])
1121  {
1122  nUnblocked++;
1123  }
1124  }
1125  reduce(nUnblocked, sumOp<label>());
1126 
1127 
1128 
1129  // Either do decomposition on cell centres or on agglomeration
1130 
1131  labelList finalDecomp;
1132 
1133 
1134  if (nProcSets+nConnections+nUnblocked == 0)
1135  {
1136  // No constraints, possibly weights
1137 
1138  if (nWeights > 0)
1139  {
1140  finalDecomp = decompose
1141  (
1142  mesh,
1143  mesh.cellCentres(),
1144  cellWeights
1145  );
1146  }
1147  else
1148  {
1149  finalDecomp = decompose(mesh, mesh.cellCentres());
1150  }
1151  }
1152  else
1153  {
1154  if (debug)
1155  {
1156  Info<< "Constrained decomposition:" << endl
1157  << " faces with same owner and neighbour processor : "
1158  << nUnblocked << endl
1159  << " baffle faces with same owner processor : "
1160  << nConnections << endl
1161  << " faces all on same processor : "
1162  << nProcSets << endl << endl;
1163  }
1164 
1165  // Determine local regions, separated by blockedFaces
1166  regionSplit localRegion(mesh, blockedFace, explicitConnections, false);
1167 
1168 
1169  if (debug)
1170  {
1171  Info<< "Constrained decomposition:" << endl
1172  << " split into " << localRegion.nLocalRegions()
1173  << " regions."
1174  << endl;
1175  }
1176 
1177  // Determine region cell centres
1178  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1179 
1180  // This just takes the first cell in the region. Otherwise the problem
1181  // is with cyclics - if we'd average the region centre might be
1182  // somewhere in the middle of the domain which might not be anywhere
1183  // near any of the cells.
1184 
1185  pointField regionCentres(localRegion.nLocalRegions(), point::max);
1186 
1187  forAll(localRegion, celli)
1188  {
1189  label regionI = localRegion[celli];
1190 
1191  if (regionCentres[regionI] == point::max)
1192  {
1193  regionCentres[regionI] = mesh.cellCentres()[celli];
1194  }
1195  }
1196 
1197  // Do decomposition on agglomeration
1198  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1199 
1200  scalarField regionWeights(localRegion.nLocalRegions(), 0);
1201 
1202  if (nWeights > 0)
1203  {
1204  forAll(localRegion, celli)
1205  {
1206  label regionI = localRegion[celli];
1207 
1208  regionWeights[regionI] += cellWeights[celli];
1209  }
1210  }
1211  else
1212  {
1213  forAll(localRegion, celli)
1214  {
1215  label regionI = localRegion[celli];
1216 
1217  regionWeights[regionI] += 1.0;
1218  }
1219  }
1220 
1221  finalDecomp = decompose
1222  (
1223  mesh,
1224  localRegion,
1225  regionCentres,
1226  regionWeights
1227  );
1228 
1229 
1230 
1231  // Implement the explicitConnections since above decompose
1232  // does not know about them
1233  forAll(explicitConnections, i)
1234  {
1235  const labelPair& baffle = explicitConnections[i];
1236  label f0 = baffle.first();
1237  label f1 = baffle.second();
1238 
1239  if (!blockedFace[f0] && !blockedFace[f1])
1240  {
1241  // Note: what if internal faces and owner and neighbour on
1242  // different processor? So for now just push owner side
1243  // proc
1244 
1245  const label proci = finalDecomp[mesh.faceOwner()[f0]];
1246 
1247  finalDecomp[mesh.faceOwner()[f1]] = proci;
1248  if (mesh.isInternalFace(f1))
1249  {
1250  finalDecomp[mesh.faceNeighbour()[f1]] = proci;
1251  }
1252  }
1253  else if (blockedFace[f0] != blockedFace[f1])
1254  {
1256  << "On explicit connection between faces " << f0
1257  << " and " << f1
1258  << " the two blockedFace status are not equal : "
1259  << blockedFace[f0] << " and " << blockedFace[f1]
1260  << exit(FatalError);
1261  }
1262  }
1263 
1264 
1265  // blockedFaces corresponding to processor faces need to be handled
1266  // separately since not handled by local regionSplit. We need to
1267  // walk now across coupled faces and make sure to move a whole
1268  // global region across
1269  if (Pstream::parRun())
1270  {
1271  // Re-do regionSplit
1272 
1273  // Field on cells and faces.
1274  List<minData> cellData(mesh.nCells());
1275  List<minData> faceData(mesh.nFaces());
1276 
1277  // Take over blockedFaces by seeding a negative number
1278  // (so is always less than the decomposition)
1279  label nUnblocked = 0;
1280  forAll(blockedFace, facei)
1281  {
1282  if (blockedFace[facei])
1283  {
1284  faceData[facei] = minData(-123);
1285  }
1286  else
1287  {
1288  nUnblocked++;
1289  }
1290  }
1291 
1292  // Seed unblocked faces with destination processor
1293  labelList seedFaces(nUnblocked);
1294  List<minData> seedData(nUnblocked);
1295  nUnblocked = 0;
1296 
1297  forAll(blockedFace, facei)
1298  {
1299  if (!blockedFace[facei])
1300  {
1301  label own = mesh.faceOwner()[facei];
1302  seedFaces[nUnblocked] = facei;
1303  seedData[nUnblocked] = minData(finalDecomp[own]);
1304  nUnblocked++;
1305  }
1306  }
1307 
1308 
1309  // Propagate information inwards
1310  FaceCellWave<minData> deltaCalc
1311  (
1312  mesh,
1313  seedFaces,
1314  seedData,
1315  faceData,
1316  cellData,
1317  mesh.globalData().nTotalCells()+1
1318  );
1319 
1320  // And extract
1321  forAll(finalDecomp, celli)
1322  {
1323  if (cellData[celli].valid(deltaCalc.data()))
1324  {
1325  finalDecomp[celli] = cellData[celli].data();
1326  }
1327  }
1328  }
1329 
1330 
1331  // For specifiedProcessorFaces rework the cellToProc to enforce
1332  // all on one processor since we can't guarantee that the input
1333  // to regionSplit was a single region.
1334  // E.g. faceSet 'a' with the cells split into two regions
1335  // by a notch formed by two walls
1336  //
1337  // \ /
1338  // \ /
1339  // ---a----+-----a-----
1340  //
1341  //
1342  // Note that reworking the cellToProc might make the decomposition
1343  // unbalanced.
1344  forAll(specifiedProcessorFaces, setI)
1345  {
1346  const labelList& set = specifiedProcessorFaces[setI];
1347 
1348  label proci = specifiedProcessor[setI];
1349  if (proci == -1)
1350  {
1351  // If no processor specified use the one from the
1352  // 0th element
1353  proci = finalDecomp[mesh.faceOwner()[set[0]]];
1354  }
1355 
1356  forAll(set, fI)
1357  {
1358  const face& f = mesh.faces()[set[fI]];
1359  forAll(f, fp)
1360  {
1361  const labelList& pFaces = mesh.pointFaces()[f[fp]];
1362  forAll(pFaces, i)
1363  {
1364  label facei = pFaces[i];
1365 
1366  finalDecomp[mesh.faceOwner()[facei]] = proci;
1367  if (mesh.isInternalFace(facei))
1368  {
1369  finalDecomp[mesh.faceNeighbour()[facei]] = proci;
1370  }
1371  }
1372  }
1373  }
1374  }
1375 
1376 
1377  if (debug && Pstream::parRun())
1378  {
1379  labelList nbrDecomp;
1380  syncTools::swapBoundaryCellList(mesh, finalDecomp, nbrDecomp);
1381 
1382  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1383  forAll(patches, patchi)
1384  {
1385  const polyPatch& pp = patches[patchi];
1386  if (pp.coupled())
1387  {
1388  forAll(pp, i)
1389  {
1390  label facei = pp.start()+i;
1391  label own = mesh.faceOwner()[facei];
1392  label bFacei = facei-mesh.nInternalFaces();
1393 
1394  if (!blockedFace[facei])
1395  {
1396  label ownProc = finalDecomp[own];
1397  label nbrProc = nbrDecomp[bFacei];
1398  if (ownProc != nbrProc)
1399  {
1401  << "patch:" << pp.name()
1402  << " face:" << facei
1403  << " at:" << mesh.faceCentres()[facei]
1404  << " ownProc:" << ownProc
1405  << " nbrProc:" << nbrProc
1406  << exit(FatalError);
1407  }
1408  }
1409  }
1410  }
1411  }
1412  }
1413  }
1414 
1415  return finalDecomp;
1416 }
1417 
1418 
1421  const polyMesh& mesh,
1422  boolList& blockedFace,
1423  PtrList<labelList>& specifiedProcessorFaces,
1424  labelList& specifiedProcessor,
1425  List<labelPair>& explicitConnections
1426 )
1427 {
1428  blockedFace.setSize(mesh.nFaces());
1429  blockedFace = true;
1430 
1431  specifiedProcessorFaces.clear();
1432  explicitConnections.clear();
1433 
1434  forAll(constraints_, constraintI)
1435  {
1436  constraints_[constraintI].add
1437  (
1438  mesh,
1439  blockedFace,
1440  specifiedProcessorFaces,
1441  specifiedProcessor,
1442  explicitConnections
1443  );
1444  }
1445 }
1446 
1447 
1450  const polyMesh& mesh,
1451  const boolList& blockedFace,
1452  const PtrList<labelList>& specifiedProcessorFaces,
1453  const labelList& specifiedProcessor,
1454  const List<labelPair>& explicitConnections,
1455  labelList& decomposition
1456 )
1457 {
1458  forAll(constraints_, constraintI)
1459  {
1460  constraints_[constraintI].apply
1461  (
1462  mesh,
1463  blockedFace,
1464  specifiedProcessorFaces,
1465  specifiedProcessor,
1466  explicitConnections,
1467  decomposition
1468  );
1469  }
1470 }
1471 
1472 
1475  const polyMesh& mesh,
1476  const scalarField& cellWeights
1477 )
1478 {
1479  // Collect all constraints
1480 
1481  boolList blockedFace;
1482  PtrList<labelList> specifiedProcessorFaces;
1483  labelList specifiedProcessor;
1484  List<labelPair> explicitConnections;
1485  setConstraints
1486  (
1487  mesh,
1488  blockedFace,
1489  specifiedProcessorFaces,
1490  specifiedProcessor,
1491  explicitConnections
1492  );
1493 
1494 
1495  // Construct decomposition method and either do decomposition on
1496  // cell centres or on agglomeration
1497 
1498  labelList finalDecomp = decompose
1499  (
1500  mesh,
1501  cellWeights, // optional weights
1502  blockedFace, // any cells to be combined
1503  specifiedProcessorFaces,// any whole cluster of cells to be kept
1504  specifiedProcessor,
1505  explicitConnections // baffles
1506  );
1507 
1508 
1509  // Give any constraint the option of modifying the decomposition
1510 
1511  applyConstraints
1512  (
1513  mesh,
1514  blockedFace,
1515  specifiedProcessorFaces,
1516  specifiedProcessor,
1517  explicitConnections,
1518  finalDecomp
1519  );
1520 
1521  return finalDecomp;
1522 }
1523 
1524 
1525 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
bool isLocal(const label i) const
Is on local processor.
Definition: globalIndexI.H:95
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Return for every coordinate the wanted processor number.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:119
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const word & name() const
Return name.
label nLocalRegions() const
Return local number of regions.
Definition: regionSplit.H:195
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
label nInternalFaces() const
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1047
label nFaces() const
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
const List< T > & m() const
Return the packed matrix of data.
label nCells() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:76
label nTotalCells() const
Return total number of cells in decomposed mesh.
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:474
Constraint to keep/move owner and neighbour of faceZone onto same processor.
patches[0]
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:278
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
scalar f1
Definition: createFields.H:28
static void calcCellCells(const polyMesh &mesh, const labelList &agglom, const label nLocalCoarse, const bool global, CompactListList< label > &cellCells)
Helper: determine (local or global) cellCells from mesh.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Detects baffles and keeps owner and neighbour on same processor.
An ordered pair of two objects of type <T> with first() and second() elements.
Definition: contiguous.H:49
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:310
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:124
A class for handling words, derived from string.
Definition: word.H:59
void clear()
Clear all entries from table.
Definition: HashTable.C:468
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:184
void setConstraints(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections)
Helper: extract constraints:
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1041
label size() const
Return the primary size, i.e. the number of rows.
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1174
label toLocal(const label i) const
From global to local on current processor.
Definition: globalIndexI.H:117
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1028
const vectorField & cellCentres() const
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
label readLabel(Istream &is)
Definition: label.H:64
static autoPtr< decompositionMethod > New(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
Foam::polyBoundaryMesh.
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
A packed storage unstructured matrix of objects of type <T> using an offset table for access...
void setSize(const label mRows)
Reset size of CompactListList.
static const char nl
Definition: Ostream.H:265
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
static autoPtr< decompositionConstraint > New(const dictionary &constraintsDict, const word &type)
Return a reference to the selected decompositionConstraint.
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const Type & second() const
Return second.
Definition: Pair.H:99
Constraint to keep owner and neighbour of (cyclic) patch on same processor.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const vectorField & faceCentres() const
void setSize(const label)
Reset size of List.
Definition: List.C:281
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
T * data()
Return a pointer to the first data element,.
Definition: UListI.H:149
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:397
void applyConstraints(const polyMesh &mesh, const boolList &blockedFace, const PtrList< labelList > &specifiedProcessorFaces, const labelList &specifiedProcessor, const List< labelPair > &explicitConnections, labelList &finalDecomp)
Helper: apply constraints to a decomposition. This gives.
label patchi
const List< label > & offsets() const
Return the offset table (= size()+1)
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:63
const vectorField & faceAreas() const
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:174
messageStream Info
const labelListList & pointFaces() const
dimensioned< scalar > mag(const dimensioned< Type > &)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Constraint to keep all cells connected to face or point of faceSet on a single processor.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
T & last()
Return the last element of the list.
Definition: UListI.H:128
For use with FaceCellWave. Transports minimum passive data.
Definition: minData.H:61
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
const Type & first() const
Return first.
Definition: Pair.H:87
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576