decompositionMethod.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
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 "BitOps.H"
38 #include "FaceCellWave.H"
39 
40 // Compatibility (MAY-2014)
45 
46 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
47 
48 namespace Foam
49 {
50  defineTypeNameAndDebug(decompositionMethod, 0);
51  defineRunTimeSelectionTable(decompositionMethod, dictionary);
52 
53 } // End namespace Foam
54 
55 
56 // * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
57 
58 namespace Foam
59 {
60 
61 // Find named coefficents dictionary, or use default "coeffs"
62 static inline const dictionary* cfindCoeffsDict
63 (
64  const dictionary& dict,
65  const word& coeffsName,
66  const bool allowDefault
67 )
68 {
69  const dictionary* dictptr = dict.findDict(coeffsName);
70  if (!dictptr && allowDefault)
71  {
72  dictptr = dict.findDict("coeffs");
73  }
74  return dictptr;
75 }
76 
77 } // End namespace Foam
78 
79 
80 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
81 
83 (
84  const dictionary& decompDict,
85  const word& regionName
86 )
87 {
88  const label nDomainsGlobal = decompDict.get<label>("numberOfSubdomains");
89 
90  if (!regionName.empty())
91  {
92  const dictionary& regionDict =
93  optionalRegionDict(decompDict, regionName);
94 
95  label nDomainsRegion;
96  if (regionDict.readIfPresent("numberOfSubdomains", nDomainsRegion))
97  {
98  if (nDomainsRegion >= 1 && nDomainsRegion <= nDomainsGlobal)
99  {
100  return nDomainsRegion;
101  }
102 
104  << "ignoring out of range numberOfSubdomains "
105  << nDomainsRegion << " for region " << regionName
106  << nl << nl
107  << endl;
108  }
109  }
110 
111  return nDomainsGlobal;
112 }
113 
114 
116 (
117  const dictionary& decompDict,
118  const word& regionName
119 )
120 {
121  const dictionary* dictptr = nullptr;
122  if
123  (
124  !regionName.empty()
125  && (dictptr = decompDict.findDict("regions")) != nullptr
126  )
127  {
128  dictptr = dictptr->findDict(regionName);
129  }
130  return (dictptr ? *dictptr : dictionary::null);
131 }
132 
133 
134 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
135 
136 bool Foam::decompositionMethod::constraintCompat(const word& modelType) const
137 {
138  bool usable = decompDict_.found(modelType);
139  if (!usable)
140  {
141  return false;
142  }
143 
144  for (const auto& item : constraints_)
145  {
146  if (modelType == item.type())
147  {
148  usable = false;
149  break;
150  }
151  }
152 
153  if (usable)
154  {
155  Warning
156  << nl << " Using '" << modelType
157  << "' constraint specification." << nl;
158  }
159  else
160  {
161  Warning
162  << nl << " Ignoring '" << modelType
163  << "' constraint specification - was already specified." << nl;
164  }
165 
166  // The syntax changed MAY-2014
167  error::warnAboutAge("constraint keyword", 1406);
168 
169  return usable;
170 }
171 
172 
173 void Foam::decompositionMethod::readConstraints()
174 {
175  constraints_.clear();
176 
177  const dictionary* dictptr = decompDict_.findDict("constraints");
178 
179  if (dictptr)
180  {
181  for (const entry& dEntry : *dictptr)
182  {
183  if (!dEntry.isDict()) // safety
184  {
185  // Ignore or warn
186  continue;
187  }
188 
189  const dictionary& dict = dEntry.dict();
190 
191  if (dict.getOrDefault("enabled", true))
192  {
193  constraints_.append(decompositionConstraint::New(dict));
194  }
195  }
196  }
197 
198  // Backwards compatibility (MAY-2014)
199  if (constraintCompat("preserveBaffles"))
200  {
201  constraints_.append
202  (
203  new decompositionConstraints::preserveBaffles()
204  );
205  }
206 
207  if (constraintCompat("preservePatches"))
208  {
209  constraints_.append
210  (
211  new decompositionConstraints::preservePatches
212  (
213  decompDict_.get<wordRes>("preservePatches")
214  )
215  );
216  }
217 
218  if (constraintCompat("preserveFaceZones"))
219  {
220  constraints_.append
221  (
222  new decompositionConstraints::preserveFaceZones
223  (
224  decompDict_.get<wordRes>("preserveFaceZones")
225  )
226  );
227  }
228 
229  if (constraintCompat("singleProcessorFaceSets"))
230  {
231  constraints_.append
232  (
233  new decompositionConstraints::singleProcessorFaceSets
234  (
235  decompDict_.lookup("singleProcessorFaceSets")
236  )
237  );
238  }
239 }
240 
241 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
242 
244 (
245  const dictionary& dict,
246  const word& coeffsName,
247  int select
248 )
249 {
250  const bool allowDefault = !(select & selectionType::EXACT);
251 
252  const dictionary* dictptr =
253  cfindCoeffsDict(dict, coeffsName, allowDefault);
254 
255  if (dictptr)
256  {
257  return *dictptr;
258  }
259 
260  // Not found
261  if (select & selectionType::MANDATORY)
262  {
264  << "'" << coeffsName << "' dictionary not found in dictionary "
265  << dict.name() << endl
266  << abort(FatalIOError);
267  }
268 
269  if (select & selectionType::NULL_DICT)
270  {
271  return dictionary::null;
272  }
273 
274  return dict; // Return original dictionary
275 }
276 
277 
279 (
280  const word& coeffsName,
281  int select
282 ) const
283 {
284  const bool allowDefault = !(select & selectionType::EXACT);
285 
286  const dictionary* dictptr = nullptr;
287 
288  if (!decompRegionDict_.empty())
289  {
290  // Region-specific dictionary
291  dictptr = cfindCoeffsDict(decompRegionDict_, coeffsName, allowDefault);
292  }
293  if (!dictptr)
294  {
295  // General
296  dictptr = cfindCoeffsDict(decompDict_, coeffsName, allowDefault);
297  }
298 
299  if (dictptr)
300  {
301  return *dictptr;
302  }
303 
304  // Not found
305  if (select & selectionType::MANDATORY)
306  {
308  << "'" << coeffsName << "' dictionary not found in dictionary "
309  << decompDict_.name() << endl
310  << abort(FatalIOError);
311  }
312 
313  if (select & selectionType::NULL_DICT)
314  {
315  return dictionary::null;
316  }
317 
318  return decompDict_; // Return general dictionary
319 }
320 
321 
322 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
323 
324 Foam::decompositionMethod::decompositionMethod
325 (
326  const dictionary& decompDict,
327  const word& regionName
328 )
329 :
330  decompDict_(decompDict),
331  decompRegionDict_
332  (
333  optionalRegionDict(decompDict_, regionName)
334  ),
335  nDomains_(nDomains(decompDict, regionName))
336 {
337  readConstraints();
338 }
339 
340 
341 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
342 
344 (
345  const dictionary& decompDict,
346  const word& regionName
347 )
348 {
349  word methodType(decompDict.get<word>("method"));
350 
351  const dictionary& regionDict = optionalRegionDict(decompDict, regionName);
352  regionDict.readIfPresent("method", methodType);
353 
354  auto* ctorPtr = dictionaryConstructorTable(methodType);
355 
356  if (!ctorPtr)
357  {
359  (
360  decompDict,
361  "decompositionMethod",
362  methodType,
363  *dictionaryConstructorTablePtr_
364  ) << exit(FatalIOError);
365  }
366 
367  // verbose
368  {
369  Info<< "Decomposition method " << methodType
370  << " [" << (nDomains(decompDict, regionName)) << ']';
371 
372  if (!regionName.empty())
373  {
374  Info<< " (region " << regionName << ')';
375  }
376  Info<< endl;
377  }
378 
379  return autoPtr<decompositionMethod>(ctorPtr(decompDict, regionName));
380 }
381 
382 
384 (
385  const polyMesh& mesh,
386  const pointField& points
387 ) const
388 {
389  scalarField weights(points.size(), scalar(1));
390 
391  return decompose(mesh, points, weights);
392 }
393 
394 
396 (
397  const polyMesh& mesh,
398  const labelList& fineToCoarse,
399  const pointField& coarsePoints,
400  const scalarField& coarseWeights
401 ) const
402 {
403  CompactListList<label> coarseCellCells;
404  calcCellCells
405  (
406  mesh,
407  fineToCoarse,
408  coarsePoints.size(),
409  true, // use global cell labels
410  coarseCellCells
411  );
412 
413  // Decompose based on agglomerated points
414  labelList coarseDistribution
415  (
416  decompose
417  (
418  coarseCellCells(),
419  coarsePoints,
420  coarseWeights
421  )
422  );
423 
424  // Rework back into decomposition for original mesh_
425  labelList fineDistribution(fineToCoarse.size());
426 
427  forAll(fineDistribution, i)
428  {
429  fineDistribution[i] = coarseDistribution[fineToCoarse[i]];
430  }
431 
432  return fineDistribution;
433 }
434 
435 
437 (
438  const polyMesh& mesh,
439  const labelList& fineToCoarse,
440  const pointField& coarsePoints
441 ) const
442 {
443  scalarField weights(coarsePoints.size(), scalar(1));
444 
445  return decompose
446  (
447  mesh,
448  fineToCoarse,
449  coarsePoints,
450  weights
451  );
452 }
453 
454 
456 (
457  const labelListList& globalCellCells,
458  const pointField& cc
459 ) const
460 {
461  scalarField weights(cc.size(), scalar(1));
462 
463  return decompose(globalCellCells, cc, weights);
464 }
465 
466 
468 (
469  const polyMesh& mesh,
470  const labelList& agglom,
471  const label nLocalCoarse,
472  const bool parallel,
473  CompactListList<label>& cellCells
474 )
475 {
476  const labelList& faceOwner = mesh.faceOwner();
477  const labelList& faceNeighbour = mesh.faceNeighbour();
479 
480 
481  // Create global cell numbers
482  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
483 
484  globalIndex globalAgglom
485  (
486  nLocalCoarse,
489  parallel
490  );
491 
492 
493  // Get agglomerate owner on other side of coupled faces
494  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
495 
496  labelList globalNeighbour(mesh.nBoundaryFaces());
497 
498  for (const polyPatch& pp : patches)
499  {
500  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
501  {
502  label facei = pp.start();
503  label bFacei = pp.start() - mesh.nInternalFaces();
504 
505  forAll(pp, i)
506  {
507  globalNeighbour[bFacei] = globalAgglom.toGlobal
508  (
509  agglom[faceOwner[facei]]
510  );
511 
512  ++facei;
513  ++bFacei;
514  }
515  }
516  }
517 
518  // Get the cell on the other side of coupled patches
519  syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
520 
521 
522  // Count number of faces (internal + coupled)
523  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
524 
525  // Number of faces per coarse cell
526  labelList nFacesPerCell(nLocalCoarse, Zero);
527 
528  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
529  {
530  const label own = agglom[faceOwner[facei]];
531  const label nei = agglom[faceNeighbour[facei]];
532 
533  nFacesPerCell[own]++;
534  nFacesPerCell[nei]++;
535  }
536 
537  for (const polyPatch& pp : patches)
538  {
539  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
540  {
541  label facei = pp.start();
542  label bFacei = pp.start()-mesh.nInternalFaces();
543 
544  forAll(pp, i)
545  {
546  const label own = agglom[faceOwner[facei]];
547  const label globalNei = globalNeighbour[bFacei];
548 
549  if
550  (
551  !globalAgglom.isLocal(globalNei)
552  || globalAgglom.toLocal(globalNei) != own
553  )
554  {
555  nFacesPerCell[own]++;
556  }
557 
558  ++facei;
559  ++bFacei;
560  }
561  }
562  }
563 
564 
565  // Fill in offset and data
566  // ~~~~~~~~~~~~~~~~~~~~~~~
567 
568  cellCells.setSize(nFacesPerCell);
569 
570  nFacesPerCell = 0;
571 
572  labelList& m = cellCells.m();
573  const labelList& offsets = cellCells.offsets();
574 
575  // For internal faces is just offsetted owner and neighbour
576  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
577  {
578  const label own = agglom[faceOwner[facei]];
579  const label nei = agglom[faceNeighbour[facei]];
580 
581  m[offsets[own] + nFacesPerCell[own]++] = globalAgglom.toGlobal(nei);
582  m[offsets[nei] + nFacesPerCell[nei]++] = globalAgglom.toGlobal(own);
583  }
584 
585  // For boundary faces is offsetted coupled neighbour
586  for (const polyPatch& pp : patches)
587  {
588  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
589  {
590  label facei = pp.start();
591  label bFacei = pp.start()-mesh.nInternalFaces();
592 
593  forAll(pp, i)
594  {
595  const label own = agglom[faceOwner[facei]];
596  const label globalNei = globalNeighbour[bFacei];
597 
598  if
599  (
600  !globalAgglom.isLocal(globalNei)
601  || globalAgglom.toLocal(globalNei) != own
602  )
603  {
604  m[offsets[own] + nFacesPerCell[own]++] = globalNei;
605  }
606 
607  ++facei;
608  ++bFacei;
609  }
610  }
611  }
612 
613 
614  // Check for duplicates connections between cells
615  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
616  // Done as postprocessing step since we now have cellCells.
617 
618  if (cellCells.size() == 0)
619  {
620  return;
621  }
622 
623  label newIndex = 0;
624  labelHashSet nbrCells;
625 
626  label startIndex = cellCells.offsets()[0];
627 
628  forAll(cellCells, celli)
629  {
630  nbrCells.clear();
631  nbrCells.insert(globalAgglom.toGlobal(celli));
632 
633  const label endIndex = cellCells.offsets()[celli+1];
634 
635  for (label i = startIndex; i < endIndex; ++i)
636  {
637  if (nbrCells.insert(cellCells.m()[i]))
638  {
639  cellCells.m()[newIndex++] = cellCells.m()[i];
640  }
641  }
642  startIndex = endIndex;
643  cellCells.offsets()[celli+1] = newIndex;
644  }
645 
646  cellCells.m().setSize(newIndex);
647 
648  //forAll(cellCells, celli)
649  //{
650  // Pout<< "Original: Coarse cell " << celli << endl;
651  // forAll(mesh.cellCells()[celli], i)
652  // {
653  // Pout<< " nbr:" << mesh.cellCells()[celli][i] << endl;
654  // }
655  // Pout<< "Compacted: Coarse cell " << celli << endl;
656  // const labelUList cCells = cellCells[celli];
657  // forAll(cCells, i)
658  // {
659  // Pout<< " nbr:" << cCells[i] << endl;
660  // }
661  //}
662 }
663 
664 
666 (
667  const polyMesh& mesh,
668  const labelList& agglom,
669  const label nLocalCoarse,
670  const bool parallel,
671  CompactListList<label>& cellCells,
672  CompactListList<scalar>& cellCellWeights
673 )
674 {
675  const labelList& faceOwner = mesh.faceOwner();
676  const labelList& faceNeighbour = mesh.faceNeighbour();
678 
679 
680  // Create global cell numbers
681  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
682 
683  globalIndex globalAgglom
684  (
685  nLocalCoarse,
688  parallel
689  );
690 
691 
692  // Get agglomerate owner on other side of coupled faces
693  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
694 
695  labelList globalNeighbour(mesh.nBoundaryFaces());
696 
697  for (const polyPatch& pp : patches)
698  {
699  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
700  {
701  label facei = pp.start();
702  label bFacei = pp.start() - mesh.nInternalFaces();
703 
704  forAll(pp, i)
705  {
706  globalNeighbour[bFacei] = globalAgglom.toGlobal
707  (
708  agglom[faceOwner[facei]]
709  );
710 
711  ++facei;
712  ++bFacei;
713  }
714  }
715  }
716 
717  // Get the cell on the other side of coupled patches
718  syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
719 
720 
721  // Count number of faces (internal + coupled)
722  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
723 
724  // Number of faces per coarse cell
725  labelList nFacesPerCell(nLocalCoarse, Zero);
726 
727  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
728  {
729  const label own = agglom[faceOwner[facei]];
730  const label nei = agglom[faceNeighbour[facei]];
731 
732  nFacesPerCell[own]++;
733  nFacesPerCell[nei]++;
734  }
735 
736  for (const polyPatch& pp : patches)
737  {
738  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
739  {
740  label facei = pp.start();
741  label bFacei = pp.start() - mesh.nInternalFaces();
742 
743  forAll(pp, i)
744  {
745  const label own = agglom[faceOwner[facei]];
746  const label globalNei = globalNeighbour[bFacei];
747 
748  if
749  (
750  !globalAgglom.isLocal(globalNei)
751  || globalAgglom.toLocal(globalNei) != own
752  )
753  {
754  nFacesPerCell[own]++;
755  }
756 
757  ++facei;
758  ++bFacei;
759  }
760  }
761  }
762 
763 
764  // Fill in offset and data
765  // ~~~~~~~~~~~~~~~~~~~~~~~
766 
767  cellCells.setSize(nFacesPerCell);
768  cellCellWeights.setSize(nFacesPerCell);
769 
770  nFacesPerCell = 0;
771 
772  labelList& m = cellCells.m();
773  scalarList& w = cellCellWeights.m();
774  const labelList& offsets = cellCells.offsets();
775 
776  // For internal faces is just offsetted owner and neighbour
777  for (label facei = 0; facei < mesh.nInternalFaces(); ++facei)
778  {
779  const label own = agglom[faceOwner[facei]];
780  const label nei = agglom[faceNeighbour[facei]];
781 
782  const label ownIndex = offsets[own] + nFacesPerCell[own]++;
783  const label neiIndex = offsets[nei] + nFacesPerCell[nei]++;
784 
785  m[ownIndex] = globalAgglom.toGlobal(nei);
786  w[ownIndex] = mag(mesh.faceAreas()[facei]);
787  m[neiIndex] = globalAgglom.toGlobal(own);
788  w[ownIndex] = mag(mesh.faceAreas()[facei]);
789  }
790 
791  // For boundary faces is offsetted coupled neighbour
792  for (const polyPatch& pp : patches)
793  {
794  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
795  {
796  label facei = pp.start();
797  label bFacei = pp.start()-mesh.nInternalFaces();
798 
799  forAll(pp, i)
800  {
801  const label own = agglom[faceOwner[facei]];
802  const label globalNei = globalNeighbour[bFacei];
803 
804  if
805  (
806  !globalAgglom.isLocal(globalNei)
807  || globalAgglom.toLocal(globalNei) != own
808  )
809  {
810  const label ownIndex = offsets[own] + nFacesPerCell[own]++;
811  m[ownIndex] = globalNei;
812  w[ownIndex] = mag(mesh.faceAreas()[facei]);
813  }
814 
815  ++facei;
816  ++bFacei;
817  }
818  }
819  }
820 
821 
822  // Check for duplicates connections between cells
823  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
824  // Done as postprocessing step since we now have cellCells.
825 
826  if (cellCells.size() == 0)
827  {
828  return;
829  }
830 
831  label newIndex = 0;
832  labelHashSet nbrCells;
833 
834  label startIndex = cellCells.offsets()[0];
835 
836  forAll(cellCells, celli)
837  {
838  nbrCells.clear();
839  nbrCells.insert(globalAgglom.toGlobal(celli));
840 
841  const label endIndex = cellCells.offsets()[celli+1];
842 
843  for (label i = startIndex; i < endIndex; ++i)
844  {
845  if (nbrCells.insert(cellCells.m()[i]))
846  {
847  cellCells.m()[newIndex] = cellCells.m()[i];
848  cellCellWeights.m()[newIndex] = cellCellWeights.m()[i];
849  newIndex++;
850  }
851  }
852  startIndex = endIndex;
853  cellCells.offsets()[celli+1] = newIndex;
854  cellCellWeights.offsets()[celli+1] = newIndex;
855  }
856 
857  cellCells.m().setSize(newIndex);
858  cellCellWeights.m().setSize(newIndex);
859 }
860 
861 
862 // NOTE:
863 // - alternative calcCellCells that handled explicitConnections was
864 // deactivated (2014 or earlier) and finally removed APR-2018.
865 
867 (
868  const polyMesh& mesh,
869  const scalarField& cellWeights,
870 
871  //- Whether owner and neighbour should be on same processor
872  // (takes priority over explicitConnections)
873  const boolList& blockedFace,
874 
875  //- Whether whole sets of faces (and point neighbours) need to be kept
876  // on single processor
877  const PtrList<labelList>& specifiedProcessorFaces,
878  const labelList& specifiedProcessor,
879 
880  //- Additional connections between boundary faces
881  const List<labelPair>& explicitConnections
882 ) const
883 {
884  // Any weights specified?
885  const bool hasWeights = returnReduce(!cellWeights.empty(), orOp<bool>());
886 
887  if (hasWeights && cellWeights.size() != mesh.nCells())
888  {
890  << "Number of weights " << cellWeights.size()
891  << " differs from number of cells " << mesh.nCells()
892  << exit(FatalError);
893  }
894 
895  // Any faces not blocked?
896  const bool hasUnblocked =
898  (
899  (!blockedFace.empty() && !BitOps::all(blockedFace)),
900  orOp<bool>()
901  );
902 
903 
904  // Any non-mesh connections?
905  const label nConnections = returnReduce
906  (
907  explicitConnections.size(),
908  sumOp<label>()
909  );
910 
911 
912  // Any processor sets?
913  label nProcSets = 0;
914  for (const labelList& procset : specifiedProcessorFaces)
915  {
916  nProcSets += procset.size();
917  }
918  reduce(nProcSets, sumOp<label>());
919 
920 
921  // Either do decomposition on cell centres or on agglomeration
922 
923  if (!hasUnblocked && !nConnections && !nProcSets)
924  {
925  // No constraints, possibly weights
926 
927  return
928  (
929  hasWeights
930  ? decompose(mesh, mesh.cellCentres(), cellWeights)
931  : decompose(mesh, mesh.cellCentres())
932  );
933  }
934 
935 
936  // The harder work.
937  // When we have processor sets, connections, or blocked faces.
938 
939 
940  // Determine local regions, separated by blockedFaces
941  regionSplit localRegion(mesh, blockedFace, explicitConnections, false);
942 
943  if (debug)
944  {
945  // Only need to count unblocked faces for debugging
946  const label nUnblocked =
947  (
948  hasUnblocked
949  ? returnReduce
950  (
951  label(BitOps::count(blockedFace, false)),
952  sumOp<label>()
953  )
954  : 0
955  );
956 
957  Info<< "Constrained decomposition:" << nl
958  << " faces with same owner and neighbour processor : "
959  << nUnblocked << nl
960  << " baffle faces with same owner processor : "
961  << nConnections << nl
962  << " faces all on same processor : "
963  << nProcSets << nl
964  << " split into " << localRegion.nLocalRegions()
965  << " regions."
966  << endl;
967  }
968 
969 
970  // Gather region weights and determine region cell centres
971  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
972 
973  // For the region centre, just take the first cell in the region.
974  // If we average the region centre instead, cyclics could cause
975  // the average domain centre to be outside of domain.
976 
977  scalarField regionWeights(localRegion.nLocalRegions(), Zero);
978 
979  pointField regionCentres(localRegion.nLocalRegions(), point::max);
980 
981  if (hasWeights)
982  {
983  forAll(localRegion, celli)
984  {
985  const label regioni = localRegion[celli];
986 
987  regionWeights[regioni] += cellWeights[celli];
988 
989  if (regionCentres[regioni] == point::max)
990  {
991  regionCentres[regioni] = mesh.cellCentres()[celli];
992  }
993  }
994  }
995  else
996  {
997  forAll(localRegion, celli)
998  {
999  const label regioni = localRegion[celli];
1000 
1001  regionWeights[regioni] += 1.0;
1002 
1003  if (regionCentres[regioni] == point::max)
1004  {
1005  regionCentres[regioni] = mesh.cellCentres()[celli];
1006  }
1007  }
1008  }
1009 
1010  // Do decomposition on agglomeration
1011  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
1012 
1013  labelList finalDecomp =
1014  decompose
1015  (
1016  mesh,
1017  localRegion,
1018  regionCentres,
1019  regionWeights
1020  );
1021 
1022 
1023  // Apply explicitConnections since decompose did not know about them
1024  for (const labelPair& baffle : explicitConnections)
1025  {
1026  const label f0 = baffle.first();
1027  const label f1 = baffle.second();
1028 
1029  if (!blockedFace[f0] && !blockedFace[f1])
1030  {
1031  // Note: what if internal faces and owner and neighbour on
1032  // different processor?
1033  // So for now just push owner side proc
1034 
1035  const label proci = finalDecomp[mesh.faceOwner()[f0]];
1036 
1037  finalDecomp[mesh.faceOwner()[f1]] = proci;
1038  if (mesh.isInternalFace(f1))
1039  {
1040  finalDecomp[mesh.faceNeighbour()[f1]] = proci;
1041  }
1042  }
1043  else if (blockedFace[f0] != blockedFace[f1])
1044  {
1046  << "On explicit connection between faces " << f0
1047  << " and " << f1
1048  << " the two blockedFace status are not equal : "
1049  << blockedFace[f0] << " and " << blockedFace[f1]
1050  << exit(FatalError);
1051  }
1052  }
1053 
1054 
1055  // blockedFaces corresponding to processor faces need to be handled
1056  // separately since not handled by local regionSplit. We need to
1057  // walk now across coupled faces and make sure to move a whole
1058  // global region across
1059 
1060  // This additionally consolidates/compacts the regions numbers globally,
1061  // since that was skipped in the previous regionSplit.
1062  if (Pstream::parRun())
1063  {
1064  // Re-do regionSplit
1065 
1066  // Field on cells and faces.
1067  List<minData> cellData(mesh.nCells());
1068  List<minData> faceData(mesh.nFaces());
1069 
1070  // Take over blockedFaces by seeding a negative number
1071  // (so is always less than the decomposition)
1072  label nUnblocked = 0;
1073  forAll(blockedFace, facei)
1074  {
1075  if (blockedFace[facei])
1076  {
1077  faceData[facei] = minData(-123);
1078  }
1079  else
1080  {
1081  ++nUnblocked;
1082  }
1083  }
1084 
1085  // Seed unblocked faces with destination processor
1086  labelList seedFaces(nUnblocked);
1087  List<minData> seedData(nUnblocked);
1088  nUnblocked = 0;
1089 
1090  forAll(blockedFace, facei)
1091  {
1092  if (!blockedFace[facei])
1093  {
1094  const label own = mesh.faceOwner()[facei];
1095  seedFaces[nUnblocked] = facei;
1096  seedData[nUnblocked] = minData(finalDecomp[own]);
1097  nUnblocked++;
1098  }
1099  }
1100 
1101 
1102  // Propagate information inwards
1103  FaceCellWave<minData> deltaCalc
1104  (
1105  mesh,
1106  seedFaces,
1107  seedData,
1108  faceData,
1109  cellData,
1111  );
1112 
1113  // And extract
1114  forAll(finalDecomp, celli)
1115  {
1116  if (cellData[celli].valid(deltaCalc.data()))
1117  {
1118  finalDecomp[celli] = cellData[celli].data();
1119  }
1120  }
1121  }
1122 
1123 
1124  // For specifiedProcessorFaces rework the cellToProc to enforce
1125  // all on one processor since we can't guarantee that the input
1126  // to regionSplit was a single region.
1127  // E.g. faceSet 'a' with the cells split into two regions
1128  // by a notch formed by two walls
1129  //
1130  // \ /
1131  // \ /
1132  // ---a----+-----a-----
1133  //
1134  //
1135  // Note that reworking the cellToProc might make the decomposition
1136  // unbalanced.
1137  forAll(specifiedProcessorFaces, seti)
1138  {
1139  const labelList& set = specifiedProcessorFaces[seti];
1140 
1141  label proci = specifiedProcessor[seti];
1142  if (proci == -1)
1143  {
1144  // If no processor specified - use the one from the 0th element
1145  if (set.size())
1146  {
1147  proci = finalDecomp[mesh.faceOwner()[set[0]]];
1148  }
1149  else
1150  {
1151  // Zero-sized processor (e.g. from redistributePar)
1152  proci = 0;
1153  }
1154  }
1155 
1156  for (const label facei : set)
1157  {
1158  const face& f = mesh.faces()[facei];
1159  for (const label pointi : f)
1160  {
1161  const labelList& pFaces = mesh.pointFaces()[pointi];
1162  for (const label pFacei : pFaces)
1163  {
1164  finalDecomp[mesh.faceOwner()[pFacei]] = proci;
1165  if (mesh.isInternalFace(pFacei))
1166  {
1167  finalDecomp[mesh.faceNeighbour()[pFacei]] = proci;
1168  }
1169  }
1170  }
1171  }
1172  }
1173 
1174 
1175  if (debug && Pstream::parRun())
1176  {
1177  labelList nbrDecomp;
1178  syncTools::swapBoundaryCellList(mesh, finalDecomp, nbrDecomp);
1179 
1181  for (const polyPatch& pp : patches)
1182  {
1183  if (pp.coupled())
1184  {
1185  forAll(pp, i)
1186  {
1187  const label facei = pp.start()+i;
1188  const label own = mesh.faceOwner()[facei];
1189  const label bFacei = facei-mesh.nInternalFaces();
1190 
1191  if (!blockedFace[facei])
1192  {
1193  const label ownProc = finalDecomp[own];
1194  const label nbrProc = nbrDecomp[bFacei];
1195 
1196  if (ownProc != nbrProc)
1197  {
1199  << "patch:" << pp.name()
1200  << " face:" << facei
1201  << " at:" << mesh.faceCentres()[facei]
1202  << " ownProc:" << ownProc
1203  << " nbrProc:" << nbrProc
1204  << exit(FatalError);
1205  }
1206  }
1207  }
1208  }
1209  }
1210  }
1211 
1212  return finalDecomp;
1213 }
1214 
1215 
1218  const polyMesh& mesh,
1219  boolList& blockedFace,
1220  PtrList<labelList>& specifiedProcessorFaces,
1221  labelList& specifiedProcessor,
1222  List<labelPair>& explicitConnections
1223 ) const
1224 {
1225  blockedFace.setSize(mesh.nFaces());
1226  blockedFace = true;
1227 
1228  specifiedProcessorFaces.clear();
1229  explicitConnections.clear();
1230 
1231  for (const decompositionConstraint& decompConstraint : constraints_)
1232  {
1233  decompConstraint.add
1234  (
1235  mesh,
1236  blockedFace,
1237  specifiedProcessorFaces,
1238  specifiedProcessor,
1239  explicitConnections
1240  );
1241  }
1242 }
1243 
1244 
1247  const polyMesh& mesh,
1248  const boolList& blockedFace,
1249  const PtrList<labelList>& specifiedProcessorFaces,
1250  const labelList& specifiedProcessor,
1251  const List<labelPair>& explicitConnections,
1252  labelList& decomposition
1253 ) const
1254 {
1255  for (const decompositionConstraint& decompConstraint : constraints_)
1256  {
1257  decompConstraint.apply
1258  (
1259  mesh,
1260  blockedFace,
1261  specifiedProcessorFaces,
1262  specifiedProcessor,
1263  explicitConnections,
1264  decomposition
1265  );
1266  }
1267 }
1268 
1269 
1272  const polyMesh& mesh,
1273  const scalarField& cellWeights
1274 ) const
1275 {
1276  // Collect all constraints
1277 
1278  boolList blockedFace;
1279  PtrList<labelList> specifiedProcessorFaces;
1280  labelList specifiedProcessor;
1281  List<labelPair> explicitConnections;
1282  setConstraints
1283  (
1284  mesh,
1285  blockedFace,
1286  specifiedProcessorFaces,
1287  specifiedProcessor,
1288  explicitConnections
1289  );
1290 
1291 
1292  // Construct decomposition method and either do decomposition on
1293  // cell centres or on agglomeration
1294 
1295  labelList finalDecomp = decompose
1296  (
1297  mesh,
1298  cellWeights, // optional weights
1299  blockedFace, // any cells to be combined
1300  specifiedProcessorFaces,// any whole cluster of cells to be kept
1301  specifiedProcessor,
1302  explicitConnections // baffles
1303  );
1304 
1305 
1306  // Give any constraint the option of modifying the decomposition
1307 
1308  applyConstraints
1309  (
1310  mesh,
1311  blockedFace,
1312  specifiedProcessorFaces,
1313  specifiedProcessor,
1314  explicitConnections,
1315  finalDecomp
1316  );
1317 
1318  return finalDecomp;
1319 }
1320 
1321 
1322 // * * * * * * * * * * * * * * * Stub Functions * * * * * * * * * * * * * * //
1323 
1326  const pointField& points,
1327  const scalarField& pointWeights
1328 ) const
1329 {
1331  return labelList();
1332 }
1333 
1334 
1337  const pointField& points
1338 ) const
1339 {
1341  return labelList();
1342 }
1343 
1344 
1345 // ************************************************************************* //
Foam::CompactListList::offsets
const List< label > & offsets() const
Return the offset table (= size()+1)
Definition: CompactListListI.H:142
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::dictionary::findDict
dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX)
Find and return a sub-dictionary pointer if present.
Definition: dictionaryI.H:127
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
BitOps.H
preserveFaceZonesConstraint.H
Foam::CompactListList::m
const List< T > & m() const
Return the packed matrix of data.
Definition: CompactListListI.H:156
Foam::decompositionMethod::calcCellCells
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.
Definition: decompositionMethod.C:468
Foam::decompositionConstraint
Abstract class for handling decomposition constraints.
Definition: decompositionConstraint.H:58
Foam::decompositionMethod::decompDict_
const dictionary & decompDict_
Top-level decomposition dictionary (eg, decomposeParDict)
Definition: decompositionMethod.H:89
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::minData
For use with FaceCellWave. Transports minimum passive data.
Definition: minData.H:60
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::decompositionMethod::optionalRegionDict
static const dictionary & optionalRegionDict(const dictionary &decompDict, const word &regionName)
Definition: decompositionMethod.C:116
Tuple2.H
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::CompactListList
A packed storage unstructured matrix of objects of type <T> using an offset table for access.
Definition: CompactListList.H:63
globalIndex.H
Foam::Warning
messageStream Warning
Foam::primitiveMesh::pointFaces
const labelListList & pointFaces() const
Definition: primitiveMeshPointFaces.C:34
localPointRegion.H
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &faceValues)
Swap coupled boundary face values. Uses eqOp.
Definition: syncTools.H:445
Foam::cfindCoeffsDict
static const dictionary * cfindCoeffsDict(const dictionary &dict, const word &coeffsName, const bool allowDefault)
Definition: decompositionMethod.C:63
Foam::globalIndex::isLocal
bool isLocal(const label i) const
Is on local processor.
Definition: globalIndexI.H:224
Foam::FatalIOError
IOerror FatalIOError
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
decompositionMethod.H
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
Foam::decompositionMethod::setConstraints
void setConstraints(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections) const
Helper: extract constraints:
Definition: decompositionMethod.C:1217
Foam::dictionary::name
const fileName & name() const noexcept
The dictionary name.
Definition: dictionaryI.H:48
Foam::HashSet< label, Hash< label > >
syncTools.H
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::BitOps::all
bool all(const UList< bool > &bools)
True if all entries are 'true' or if the set is empty.
Definition: BitOps.H:84
Foam::sumOp
Definition: ops.H:213
Foam::decompositionConstraint::New
static autoPtr< decompositionConstraint > New(const dictionary &constraintDict)
Return a reference to the selected decompositionConstraint.
Definition: decompositionConstraint.C:155
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::dictionary::null
static const dictionary null
An empty dictionary, which is also the parent for all dictionaries.
Definition: dictionary.H:392
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:517
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
FatalIOErrorInLookup
#define FatalIOErrorInLookup(ios, lookupTag, lookupName, lookupTable)
Report an error message using Foam::FatalIOError.
Definition: error.H:478
Foam::Field< vector >
Foam::CompactListList::setSize
void setSize(const label mRows)
Reset size of CompactListList.
Definition: CompactListList.C:108
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
faceSet.H
regionSplit.H
Foam::primitiveMesh::nBoundaryFaces
label nBoundaryFaces() const noexcept
Number of boundary faces (== nFaces - nInternalFaces)
Definition: primitiveMeshI.H:84
singleProcessorFaceSetsConstraint.H
Foam::regionSplit
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:140
Foam::PtrList
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition: List.H:59
preservePatchesConstraint.H
Foam::syncTools::swapBoundaryCellList
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
Definition: syncToolsTemplates.C:1392
dict
dictionary dict
Definition: searchingEngine.H:14
minData.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
preserveBafflesConstraint.H
Foam::PtrList::clear
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:97
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::FaceCellWave
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:78
Foam::globalIndex::toLocal
label toLocal(const label i) const
From global to local on current processor.
Definition: globalIndexI.H:305
Foam::autoPtr< Foam::decompositionMethod >
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const noexcept
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:103
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::decompositionMethod::decompose
virtual labelList decompose(const pointField &points, const scalarField &pointWeights) const
Return the wanted processor number for every coordinate.
Definition: decompositionMethod.C:1325
Foam::UPstream::msgType
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:540
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Pair< label >
Foam::VectorSpace< Vector< Cmpt >, Cmpt, 3 >::max
static const Vector< Cmpt > max
Definition: VectorSpace.H:117
f
labelList f(nPoints)
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::UPstream::worldComm
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:293
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< label >
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:77
FaceCellWave.H
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::List::clear
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:116
Foam::decompositionMethod::New
static autoPtr< decompositionMethod > New(const dictionary &decompDict, const word &regionName="")
Definition: decompositionMethod.C:344
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::CompactListList::size
label size() const noexcept
The primary size (the number of rows)
Definition: CompactListListI.H:127
Foam::globalMeshData::nTotalCells
label nTotalCells() const noexcept
Return total number of cells in decomposed mesh.
Definition: globalMeshData.H:371
Foam::error::warnAboutAge
static bool warnAboutAge(const int version) noexcept
Test if an age warning should be emitted.
Definition: error.C:55
Foam::decompositionMethod::nDomains
label nDomains() const noexcept
Number of domains.
Definition: decompositionMethod.H:209
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::polyMesh::globalData
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1295
Foam::decompositionMethod::applyConstraints
void applyConstraints(const polyMesh &mesh, const boolList &blockedFace, const PtrList< labelList > &specifiedProcessorFaces, const labelList &specifiedProcessor, const List< labelPair > &explicitConnections, labelList &finalDecomp) const
Helper: apply constraints to a decomposition.
Definition: decompositionMethod.C:1246
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:148
Foam::orOp
Definition: ops.H:234
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::regionSplit::nLocalRegions
label nLocalRegions() const
Return local number of regions.
Definition: regionSplit.H:288
Foam::decompositionMethod::findCoeffsDict
static const dictionary & findCoeffsDict(const dictionary &dict, const word &coeffsName, int select=selectionType::DEFAULT)
Definition: decompositionMethod.C:244
pFaces
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=cellModel::ref(cellModel::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].reset(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::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::decompositionMethod::constraints_
PtrList< decompositionConstraint > constraints_
Optional constraints.
Definition: decompositionMethod.H:98
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:240
Foam::primitiveMesh::faceAreas
const vectorField & faceAreas() const
Definition: primitiveMeshFaceCentresAndAreas.C:89