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-2024 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "decompositionMethod.H"
27 #include "Time.H"
28 #include "globalIndex.H"
29 #include "syncTools.H"
30 #include "Tuple2.H"
31 #include "faceSet.H"
32 #include "regionSplit.H"
33 #include "localPointRegion.H"
34 #include "minData.H"
35 #include "FaceCellWave.H"
36 
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
49 }
50 
51 
52 // * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * //
53 
55 (
56  const pointField& points,
57  const scalarField& pointWeights
58 ) const
59 {
60  const label localnWeights =
61  points.size() ? pointWeights.size()/points.size() : 0;
62 
63  const label nWeights = returnReduce(localnWeights, maxOp<label>());
64 
65  if (localnWeights && localnWeights != nWeights)
66  {
68  << "Number of weights on this processor " << localnWeights
69  << " does not equal the maximum number of weights " << nWeights
70  << exit(FatalError);
71  }
72 
73  return nWeights;
74 }
75 
76 
78 (
79  const pointField& points,
80  const scalarField& pointWeights
81 ) const
82 {
83  const label nWeights = this->nWeights(points, pointWeights);
84 
85  if (nWeights > 1)
86  {
88  << "decompositionMethod " << type()
89  << " does not support multiple constraints"
90  << exit(FatalError);
91  }
92 
93  return nWeights;
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
98 
100 (
101  const dictionary& decompositionDict
102 )
103 :
104  decompositionDict_(decompositionDict),
105  nProcessors_
106  (
107  decompositionDict.lookup<label>("numberOfSubdomains")
108  )
109 {
110  // Read any constraints
111  wordList constraintTypes_;
112 
113  if (decompositionDict_.found("constraints"))
114  {
115  const dictionary& constraintsList = decompositionDict_.subDict
116  (
117  "constraints"
118  );
119  forAllConstIter(dictionary, constraintsList, iter)
120  {
121  const dictionary& dict = iter().dict();
122 
123  constraintTypes_.append(dict.lookup("type"));
124 
125  constraints_.append
126  (
128  (
129  dict,
130  constraintTypes_.last()
131  )
132  );
133  }
134  }
135 }
136 
137 
138 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
139 
142 (
143  const dictionary& decompositionDict
144 )
145 {
146  const word methodType
147  (
148  decompositionDict.lookupBackwardsCompatible<word>
149  (
150  {"decomposer", "method"}
151  )
152  );
153 
154  Info<< "Selecting decomposer " << methodType << endl;
155 
156  libs.open
157  (
158  decompositionDict,
159  "libs",
160  decomposerConstructorTablePtr_
161  );
162 
163  decomposerConstructorTable::iterator cstrIter =
164  decomposerConstructorTablePtr_->find(methodType);
165 
166  if (cstrIter == decomposerConstructorTablePtr_->end())
167  {
169  << "Unknown decomposer "
170  << methodType << nl << nl
171  << "Valid decomposers are : " << endl
172  << decomposerConstructorTablePtr_->sortedToc()
173  << exit(FatalError);
174  }
175 
176  return autoPtr<decompositionMethod>(cstrIter()(decompositionDict));
177 }
178 
179 
182 (
183  const dictionary& distributionDict
184 )
185 {
186  const word methodType
187  (
188  distributionDict.lookupBackwardsCompatible<word>
189  (
190  {"distributor", "method"}
191  )
192  );
193 
194  Info<< "Selecting distributor " << methodType << endl;
195 
196  libs.open
197  (
198  distributionDict,
199  "libs",
200  distributorConstructorTablePtr_
201  );
202 
203  distributorConstructorTable::iterator cstrIter =
204  distributorConstructorTablePtr_->find(methodType);
205 
206  if (cstrIter == distributorConstructorTablePtr_->end())
207  {
209  << "Unknown distributor "
210  << methodType << nl << nl
211  << "Valid distributors are : " << endl
212  << distributorConstructorTablePtr_->sortedToc()
213  << exit(FatalError);
214  }
215 
216  return autoPtr<decompositionMethod>(cstrIter()(distributionDict));
217 }
218 
219 
221 {
222  return IOdictionary
223  (
224  IOobject
225  (
226  "decomposeParDict",
227  time.system(),
228  time,
231  false
232  )
233  );
234 }
235 
236 
238 (
239  const polyMesh& mesh,
240  const pointField& points
241 )
242 {
243  scalarField weights(points.size(), 1.0);
244 
245  return decompose(mesh, points, weights);
246 }
247 
248 
250 (
251  const polyMesh& mesh,
252  const labelList& fineToCoarse,
253  const pointField& coarsePoints,
254  const scalarField& coarseWeights
255 )
256 {
257  CompactListList<label> coarseCellCells;
258  calcCellCells
259  (
260  mesh,
261  fineToCoarse,
262  coarsePoints.size(),
263  true, // use global cell labels
264  coarseCellCells
265  );
266 
267  // Decompose based on agglomerated points
268  labelList coarseDistribution
269  (
270  decompose
271  (
272  coarseCellCells.list(),
273  coarsePoints,
274  coarseWeights
275  )
276  );
277 
278  // Rework back into decomposition for original mesh_
279  labelList fineDistribution(fineToCoarse.size());
280 
281  forAll(fineDistribution, i)
282  {
283  fineDistribution[i] = coarseDistribution[fineToCoarse[i]];
284  }
285 
286  return fineDistribution;
287 }
288 
289 
291 (
292  const polyMesh& mesh,
293  const labelList& fineToCoarse,
294  const pointField& coarsePoints
295 )
296 {
297  scalarField cellWeights(coarsePoints.size(), 1.0);
298 
299  return decompose
300  (
301  mesh,
302  fineToCoarse,
303  coarsePoints,
304  cellWeights
305  );
306 }
307 
308 
310 (
311  const labelListList& globalCellCells,
312  const pointField& cellCentres
313 )
314 {
315  scalarField cellWeights(cellCentres.size(), 1.0);
316 
317  return decompose(globalCellCells, cellCentres, cellWeights);
318 }
319 
320 
322 (
323  const scalarField& weights,
324  label& nWeights,
325  const bool distributed
326 )
327 {
328  labelList intWeights;
329 
330  if (nWeights > 0)
331  {
332  // Calculate the scalar -> integer scaling factor
333  const scalar sumWeights
334  (
335  distributed
336  ? gSum(weights)
337  : sum(weights)
338  );
339  const scalar scale = labelMax/(2*sumWeights);
340 
341  // Convert weights to integer
342  intWeights.setSize(weights.size());
343  forAll(intWeights, i)
344  {
345  intWeights[i] = ceil(scale*weights[i]);
346  }
347 
348  /*
349  // Alternatively calculate a separate scaling factor
350  // for each weight, it is not clear from the parMETIS or Scotch manuals
351  // which method should be used.
352  //
353  // Calculate the scalar -> integer scaling factors
354  scalarList sumWeights(nWeights, 0.0);
355 
356  forAll(weights, i)
357  {
358  sumWeights[i % nWeights] += weights[i];
359  }
360 
361  if (distributed)
362  {
363  reduce(sumWeights, ListOp<sumOp<scalar>>());
364  }
365 
366  scalarList scale(nWeights, 0.0);
367  forAll(scale, i)
368  {
369  if (sumWeights[i] > small)
370  {
371  scale[i] = labelMax/(2*sumWeights[i]);
372  }
373  }
374 
375  // Convert weights to integer
376  intWeights.setSize(weights.size());
377  forAll(intWeights, i)
378  {
379  intWeights[i] = ceil(scale[i % nWeights]*weights[i]);
380  }
381  */
382 
383  // Calculate the sum over all processors of each weight
384  labelList sumIntWeights(nWeights, 0);
385  forAll(weights, i)
386  {
387  sumIntWeights[i % nWeights] += intWeights[i];
388  }
389 
390  if (distributed)
391  {
392  reduce(sumIntWeights, ListOp<sumOp<label>>());
393  }
394 
395  // Check that the sum of each weight is non-zero
396  boolList nonZeroWeights(nWeights, false);
397  label nNonZeroWeights = 0;
398  forAll(sumIntWeights, i)
399  {
400  if (sumIntWeights[i] != 0)
401  {
402  nonZeroWeights[i] = true;
403  nNonZeroWeights++;
404  }
405  }
406 
407  // If there are zero weights remove them from the weights list
408  if (nNonZeroWeights != nWeights)
409  {
410  label j = 0;
411  forAll(intWeights, i)
412  {
413  if (nonZeroWeights[i % nWeights])
414  {
415  intWeights[j++] = intWeights[i];
416  }
417  }
418 
419  // Resize the weights list
420  intWeights.setSize(nNonZeroWeights*(intWeights.size()/nWeights));
421 
422  // Reset the number of weight to the number of non-zero weights
423  nWeights = nNonZeroWeights;
424  }
425  }
426 
427  return intWeights;
428 }
429 
430 
432 (
433  const polyMesh& mesh,
434  const labelList& agglom,
435  const label nLocalCoarse,
436  const bool parallel,
437  CompactListList<label>& cellCells
438 )
439 {
440  const labelList& faceOwner = mesh.faceOwner();
441  const labelList& faceNeighbour = mesh.faceNeighbour();
442  const polyBoundaryMesh& patches = mesh.boundaryMesh();
443 
444 
445  // Create global cell numbers
446  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
447 
448  globalIndex globalAgglom
449  (
450  nLocalCoarse,
453  parallel
454  );
455 
456 
457  // Get agglomerate owner on other side of coupled faces
458  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
459 
460  labelList globalNeighbour(mesh.nFaces() - mesh.nInternalFaces());
461 
463  {
464  const polyPatch& pp = patches[patchi];
465 
466  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
467  {
468  label facei = pp.start();
469  label bFacei = pp.start() - mesh.nInternalFaces();
470 
471  forAll(pp, i)
472  {
473  globalNeighbour[bFacei] = globalAgglom.toGlobal
474  (
475  agglom[faceOwner[facei]]
476  );
477 
478  bFacei++;
479  facei++;
480  }
481  }
482  }
483 
484  // Get the cell on the other side of coupled patches
485  syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
486 
487 
488  // Count number of faces (internal + coupled)
489  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
490 
491  // Number of faces per coarse cell
492  labelList nFacesPerCell(nLocalCoarse, 0);
493 
494  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
495  {
496  label own = agglom[faceOwner[facei]];
497  label nei = agglom[faceNeighbour[facei]];
498 
499  nFacesPerCell[own]++;
500  nFacesPerCell[nei]++;
501  }
502 
504  {
505  const polyPatch& pp = patches[patchi];
506 
507  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
508  {
509  label facei = pp.start();
510  label bFacei = pp.start() - mesh.nInternalFaces();
511 
512  forAll(pp, i)
513  {
514  label own = agglom[faceOwner[facei]];
515 
516  label globalNei = globalNeighbour[bFacei];
517  if
518  (
519  !globalAgglom.isLocal(globalNei)
520  || globalAgglom.toLocal(globalNei) != own
521  )
522  {
523  nFacesPerCell[own]++;
524  }
525 
526  facei++;
527  bFacei++;
528  }
529  }
530  }
531 
532 
533  // Fill in offset and data
534  // ~~~~~~~~~~~~~~~~~~~~~~~
535 
536  cellCells.setSize(nFacesPerCell);
537 
538  nFacesPerCell = 0;
539 
540  labelUList& m = cellCells.m();
541  const labelList& offsets = cellCells.offsets();
542 
543  // For internal faces is just offsetted owner and neighbour
544  for (label facei = 0; facei < mesh.nInternalFaces(); facei++)
545  {
546  label own = agglom[faceOwner[facei]];
547  label nei = agglom[faceNeighbour[facei]];
548 
549  m[offsets[own] + nFacesPerCell[own]++] = globalAgglom.toGlobal(nei);
550  m[offsets[nei] + nFacesPerCell[nei]++] = globalAgglom.toGlobal(own);
551  }
552 
553  // For boundary faces is offsetted coupled neighbour
555  {
556  const polyPatch& pp = patches[patchi];
557 
558  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
559  {
560  label facei = pp.start();
561  label bFacei = pp.start() - mesh.nInternalFaces();
562 
563  forAll(pp, i)
564  {
565  label own = agglom[faceOwner[facei]];
566 
567  label globalNei = globalNeighbour[bFacei];
568 
569  if
570  (
571  !globalAgglom.isLocal(globalNei)
572  || globalAgglom.toLocal(globalNei) != own
573  )
574  {
575  m[offsets[own] + nFacesPerCell[own]++] = globalNei;
576  }
577 
578  facei++;
579  bFacei++;
580  }
581  }
582  }
583 
584 
585  // Check for duplicates connections between cells
586  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
587  // Done as postprocessing step since we now have cellCells.
588  label newIndex = 0;
589  labelHashSet nbrCells;
590 
591 
592  if (cellCells.size() == 0)
593  {
594  return;
595  }
596 
597  label startIndex = cellCells.offsets()[0];
598 
599  forAll(cellCells, celli)
600  {
601  nbrCells.clear();
602  nbrCells.insert(globalAgglom.toGlobal(celli));
603 
604  label endIndex = cellCells.offsets()[celli+1];
605 
606  for (label i = startIndex; i < endIndex; i++)
607  {
608  if (nbrCells.insert(cellCells.m()[i]))
609  {
610  cellCells.m()[newIndex++] = cellCells.m()[i];
611  }
612  }
613  startIndex = endIndex;
614  cellCells.offsets()[celli+1] = newIndex;
615  }
616 
617  cellCells.setSize(cellCells.size(), newIndex);
618 }
619 
620 
622 (
623  const polyMesh& mesh,
624  const labelList& agglom,
625  const label nLocalCoarse,
626  const bool parallel,
627  CompactListList<label>& cellCells,
628  CompactListList<scalar>& cellCellWeights
629 )
630 {
631  const labelList& faceOwner = mesh.faceOwner();
632  const labelList& faceNeighbour = mesh.faceNeighbour();
633  const polyBoundaryMesh& patches = mesh.boundaryMesh();
634 
635 
636  // Create global cell numbers
637  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
638 
639  globalIndex globalAgglom
640  (
641  nLocalCoarse,
644  parallel
645  );
646 
647 
648  // Get agglomerate owner on other side of coupled faces
649  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
650 
651  labelList globalNeighbour(mesh.nFaces() - mesh.nInternalFaces());
652 
654  {
655  const polyPatch& pp = patches[patchi];
656 
657  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
658  {
659  label faceI = pp.start();
660  label bFaceI = pp.start() - mesh.nInternalFaces();
661 
662  forAll(pp, i)
663  {
664  globalNeighbour[bFaceI] = globalAgglom.toGlobal
665  (
666  agglom[faceOwner[faceI]]
667  );
668 
669  bFaceI++;
670  faceI++;
671  }
672  }
673  }
674 
675  // Get the cell on the other side of coupled patches
676  syncTools::swapBoundaryFaceList(mesh, globalNeighbour);
677 
678 
679  // Count number of faces (internal + coupled)
680  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
681 
682  // Number of faces per coarse cell
683  labelList nFacesPerCell(nLocalCoarse, 0);
684 
685  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
686  {
687  label own = agglom[faceOwner[faceI]];
688  label nei = agglom[faceNeighbour[faceI]];
689 
690  nFacesPerCell[own]++;
691  nFacesPerCell[nei]++;
692  }
693 
695  {
696  const polyPatch& pp = patches[patchi];
697 
698  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
699  {
700  label faceI = pp.start();
701  label bFaceI = pp.start()-mesh.nInternalFaces();
702 
703  forAll(pp, i)
704  {
705  label own = agglom[faceOwner[faceI]];
706 
707  label globalNei = globalNeighbour[bFaceI];
708  if
709  (
710  !globalAgglom.isLocal(globalNei)
711  || globalAgglom.toLocal(globalNei) != own
712  )
713  {
714  nFacesPerCell[own]++;
715  }
716 
717  faceI++;
718  bFaceI++;
719  }
720  }
721  }
722 
723 
724  // Fill in offset and data
725  // ~~~~~~~~~~~~~~~~~~~~~~~
726 
727  cellCells.setSize(nFacesPerCell);
728  cellCellWeights.setSize(nFacesPerCell);
729 
730  nFacesPerCell = 0;
731 
732  labelUList& m = cellCells.m();
733  scalarUList& w = cellCellWeights.m();
734  const labelList& offsets = cellCells.offsets();
735 
736  // For internal faces is just offsetted owner and neighbour
737  for (label faceI = 0; faceI < mesh.nInternalFaces(); faceI++)
738  {
739  label own = agglom[faceOwner[faceI]];
740  label nei = agglom[faceNeighbour[faceI]];
741 
742  label ownIndex = offsets[own] + nFacesPerCell[own]++;
743  label neiIndex = offsets[nei] + nFacesPerCell[nei]++;
744 
745  m[ownIndex] = globalAgglom.toGlobal(nei);
746  w[ownIndex] = mesh.magFaceAreas()[faceI];
747  m[neiIndex] = globalAgglom.toGlobal(own);
748  w[ownIndex] = mesh.magFaceAreas()[faceI];
749  }
750 
751  // For boundary faces is offsetted coupled neighbour
753  {
754  const polyPatch& pp = patches[patchi];
755 
756  if (pp.coupled() && (parallel || !isA<processorPolyPatch>(pp)))
757  {
758  label faceI = pp.start();
759  label bFaceI = pp.start()-mesh.nInternalFaces();
760 
761  forAll(pp, i)
762  {
763  label own = agglom[faceOwner[faceI]];
764 
765  label globalNei = globalNeighbour[bFaceI];
766 
767  if
768  (
769  !globalAgglom.isLocal(globalNei)
770  || globalAgglom.toLocal(globalNei) != own
771  )
772  {
773  label ownIndex = offsets[own] + nFacesPerCell[own]++;
774  m[ownIndex] = globalNei;
775  w[ownIndex] = mesh.magFaceAreas()[faceI];
776  }
777 
778  faceI++;
779  bFaceI++;
780  }
781  }
782  }
783 
784 
785  // Check for duplicates connections between cells
786  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
787  // Done as postprocessing step since we now have cellCells.
788  label newIndex = 0;
789  labelHashSet nbrCells;
790 
791 
792  if (cellCells.size() == 0)
793  {
794  return;
795  }
796 
797  label startIndex = cellCells.offsets()[0];
798 
799  forAll(cellCells, cellI)
800  {
801  nbrCells.clear();
802  nbrCells.insert(globalAgglom.toGlobal(cellI));
803 
804  label endIndex = cellCells.offsets()[cellI+1];
805 
806  for (label i = startIndex; i < endIndex; i++)
807  {
808  if (nbrCells.insert(cellCells.m()[i]))
809  {
810  cellCells.m()[newIndex] = cellCells.m()[i];
811  cellCellWeights.m()[newIndex] = cellCellWeights.m()[i];
812  newIndex++;
813  }
814  }
815  startIndex = endIndex;
816  cellCells.offsets()[cellI+1] = newIndex;
817  cellCellWeights.offsets()[cellI+1] = newIndex;
818  }
819 
820  cellCells.setSize(cellCells.size(), newIndex);
821  cellCellWeights.setSize(cellCells.size(), newIndex);
822 }
823 
824 
826 (
827  const polyMesh& mesh,
828  const scalarField& cellWeights,
829  const boolList& blockedFace,
830  const PtrList<labelList>& specifiedProcessorFaces,
831  const labelList& specifiedProcessor,
832  const List<labelPair>& explicitConnections
833 )
834 {
835  // Any weights specified?
836  const bool hasWeights =
837  returnReduce(cellWeights.size(), sumOp<label>()) > 0;
838 
839  // Any processor sets?
840  label nProcSets = 0;
841  forAll(specifiedProcessorFaces, setI)
842  {
843  nProcSets += specifiedProcessorFaces[setI].size();
844  }
845  reduce(nProcSets, sumOp<label>());
846 
847  // Any non-mesh connections?
848  label nConnections = returnReduce
849  (
850  explicitConnections.size(),
851  sumOp<label>()
852  );
853 
854  // Any faces not blocked?
855  label nUnblocked = 0;
856  forAll(blockedFace, facei)
857  {
858  if (!blockedFace[facei])
859  {
860  nUnblocked++;
861  }
862  }
863  reduce(nUnblocked, sumOp<label>());
864 
865 
866  // Either do decomposition on cell centres or on agglomeration
867 
868  labelList finalDecomp;
869 
870 
871  if (nProcSets+nConnections+nUnblocked == 0)
872  {
873  // No constraints, possibly weights
874 
875  if (hasWeights)
876  {
877  finalDecomp = decompose
878  (
879  mesh,
880  mesh.cellCentres(),
881  cellWeights
882  );
883  }
884  else
885  {
886  finalDecomp = decompose(mesh, mesh.cellCentres());
887  }
888  }
889  else
890  {
891  if (debug)
892  {
893  Info<< "Constrained decomposition:" << endl
894  << " faces with same owner and neighbour processor : "
895  << nUnblocked << endl
896  << " baffle faces with same owner processor : "
897  << nConnections << endl
898  << " faces all on same processor : "
899  << nProcSets << endl << endl;
900  }
901 
902  // Determine local regions, separated by blockedFaces
903  regionSplit localRegion(mesh, blockedFace, explicitConnections, false);
904 
905 
906  if (debug)
907  {
908  Info<< "Constrained decomposition:" << endl
909  << " split into " << localRegion.nLocalRegions()
910  << " regions."
911  << endl;
912  }
913 
914  // Determine region cell centres
915  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
916 
917  // This just takes the first cell in the region. Otherwise the problem
918  // is with cyclics - if we'd average the region centre might be
919  // somewhere in the middle of the domain which might not be anywhere
920  // near any of the cells.
921 
922  pointField regionCentres(localRegion.nLocalRegions(), point::max);
923 
924  forAll(localRegion, celli)
925  {
926  label regionI = localRegion[celli];
927 
928  if (regionCentres[regionI] == point::max)
929  {
930  regionCentres[regionI] = mesh.cellCentres()[celli];
931  }
932  }
933 
934  // Do decomposition on agglomeration
935  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
936 
937  scalarField regionWeights(localRegion.nLocalRegions(), 0);
938 
939  if (hasWeights)
940  {
941  forAll(localRegion, celli)
942  {
943  const label regionI = localRegion[celli];
944 
945  regionWeights[regionI] += cellWeights[celli];
946  }
947  }
948  else
949  {
950  forAll(localRegion, celli)
951  {
952  label regionI = localRegion[celli];
953 
954  regionWeights[regionI] += 1.0;
955  }
956  }
957 
958  finalDecomp = decompose
959  (
960  mesh,
961  localRegion,
962  regionCentres,
963  regionWeights
964  );
965 
966 
967 
968  // Implement the explicitConnections since above decompose
969  // does not know about them
970  forAll(explicitConnections, i)
971  {
972  const labelPair& baffle = explicitConnections[i];
973  label f0 = baffle.first();
974  label f1 = baffle.second();
975 
976  if (!blockedFace[f0] && !blockedFace[f1])
977  {
978  // Note: what if internal faces and owner and neighbour on
979  // different processor? So for now just push owner side
980  // proc
981 
982  const label proci = finalDecomp[mesh.faceOwner()[f0]];
983 
984  finalDecomp[mesh.faceOwner()[f1]] = proci;
985  if (mesh.isInternalFace(f1))
986  {
987  finalDecomp[mesh.faceNeighbour()[f1]] = proci;
988  }
989  }
990  else if (blockedFace[f0] != blockedFace[f1])
991  {
993  << "On explicit connection between faces " << f0
994  << " and " << f1
995  << " the two blockedFace status are not equal : "
996  << blockedFace[f0] << " and " << blockedFace[f1]
997  << exit(FatalError);
998  }
999  }
1000 
1001 
1002  // blockedFaces corresponding to processor faces need to be handled
1003  // separately since not handled by local regionSplit. We need to
1004  // walk now across coupled faces and make sure to move a whole
1005  // global region across
1006  if (Pstream::parRun())
1007  {
1008  // Re-do regionSplit
1009 
1010  // Field on cells and faces.
1011  List<minData> cellData(mesh.nCells());
1012  List<minData> faceData(mesh.nFaces());
1013 
1014  // Take over blockedFaces by seeding a negative number
1015  // (so is always less than the decomposition)
1016  label nUnblocked = 0;
1017  forAll(blockedFace, facei)
1018  {
1019  if (blockedFace[facei])
1020  {
1021  faceData[facei] = minData(-123);
1022  }
1023  else
1024  {
1025  nUnblocked++;
1026  }
1027  }
1028 
1029  // Seed unblocked faces with destination processor
1030  labelList seedFaces(nUnblocked);
1031  List<minData> seedData(nUnblocked);
1032  nUnblocked = 0;
1033 
1034  forAll(blockedFace, facei)
1035  {
1036  if (!blockedFace[facei])
1037  {
1038  label own = mesh.faceOwner()[facei];
1039  seedFaces[nUnblocked] = facei;
1040  seedData[nUnblocked] = minData(finalDecomp[own]);
1041  nUnblocked++;
1042  }
1043  }
1044 
1045 
1046  // Propagate information inwards
1047  FaceCellWave<minData> deltaCalc
1048  (
1049  mesh,
1050  seedFaces,
1051  seedData,
1052  faceData,
1053  cellData,
1054  mesh.globalData().nTotalCells()+1
1055  );
1056 
1057  // And extract
1058  forAll(finalDecomp, celli)
1059  {
1060  if (cellData[celli].valid(deltaCalc.data()))
1061  {
1062  finalDecomp[celli] = cellData[celli].data();
1063  }
1064  }
1065  }
1066 
1067 
1068  // For specifiedProcessorFaces rework the cellToProc to enforce
1069  // all on one processor since we can't guarantee that the input
1070  // to regionSplit was a single region.
1071  // E.g. faceSet 'a' with the cells split into two regions
1072  // by a notch formed by two walls
1073  //
1074  // \ /
1075  // \ /
1076  // ---a----+-----a-----
1077  //
1078  //
1079  // Note that reworking the cellToProc might make the decomposition
1080  // unbalanced.
1081  forAll(specifiedProcessorFaces, setI)
1082  {
1083  const labelList& set = specifiedProcessorFaces[setI];
1084 
1085  label proci = specifiedProcessor[setI];
1086  if (proci == -1)
1087  {
1088  // If no processor specified use the one from the
1089  // 0th element
1090  proci = finalDecomp[mesh.faceOwner()[set[0]]];
1091  }
1092 
1093  forAll(set, fI)
1094  {
1095  const face& f = mesh.faces()[set[fI]];
1096  forAll(f, fp)
1097  {
1098  const labelList& pFaces = mesh.pointFaces()[f[fp]];
1099  forAll(pFaces, i)
1100  {
1101  label facei = pFaces[i];
1102 
1103  finalDecomp[mesh.faceOwner()[facei]] = proci;
1104  if (mesh.isInternalFace(facei))
1105  {
1106  finalDecomp[mesh.faceNeighbour()[facei]] = proci;
1107  }
1108  }
1109  }
1110  }
1111  }
1112 
1113 
1114  if (debug && Pstream::parRun())
1115  {
1116  labelList nbrDecomp;
1117  syncTools::swapBoundaryCellList(mesh, finalDecomp, nbrDecomp);
1118 
1119  const polyBoundaryMesh& patches = mesh.boundaryMesh();
1121  {
1122  const polyPatch& pp = patches[patchi];
1123  if (pp.coupled())
1124  {
1125  forAll(pp, i)
1126  {
1127  label facei = pp.start()+i;
1128  label own = mesh.faceOwner()[facei];
1129  label bFacei = facei-mesh.nInternalFaces();
1130 
1131  if (!blockedFace[facei])
1132  {
1133  label ownProc = finalDecomp[own];
1134  label nbrProc = nbrDecomp[bFacei];
1135  if (ownProc != nbrProc)
1136  {
1138  << "patch:" << pp.name()
1139  << " face:" << facei
1140  << " at:" << mesh.faceCentres()[facei]
1141  << " ownProc:" << ownProc
1142  << " nbrProc:" << nbrProc
1143  << exit(FatalError);
1144  }
1145  }
1146  }
1147  }
1148  }
1149  }
1150  }
1151 
1152  return finalDecomp;
1153 }
1154 
1155 
1157 (
1158  const polyMesh& mesh,
1159  boolList& blockedFace,
1160  PtrList<labelList>& specifiedProcessorFaces,
1161  labelList& specifiedProcessor,
1162  List<labelPair>& explicitConnections
1163 )
1164 {
1165  blockedFace.setSize(mesh.nFaces());
1166  blockedFace = true;
1167 
1168  specifiedProcessorFaces.clear();
1169  explicitConnections.clear();
1170 
1171  forAll(constraints_, constraintI)
1172  {
1173  constraints_[constraintI].add
1174  (
1175  mesh,
1176  blockedFace,
1177  specifiedProcessorFaces,
1178  specifiedProcessor,
1179  explicitConnections
1180  );
1181  }
1182 }
1183 
1184 
1186 (
1187  const polyMesh& mesh,
1188  const boolList& blockedFace,
1189  const PtrList<labelList>& specifiedProcessorFaces,
1190  const labelList& specifiedProcessor,
1191  const List<labelPair>& explicitConnections,
1192  labelList& decomposition
1193 )
1194 {
1195  forAll(constraints_, constraintI)
1196  {
1197  constraints_[constraintI].apply
1198  (
1199  mesh,
1200  blockedFace,
1201  specifiedProcessorFaces,
1202  specifiedProcessor,
1203  explicitConnections,
1204  decomposition
1205  );
1206  }
1207 }
1208 
1209 
1211 (
1212  const polyMesh& mesh,
1213  const scalarField& cellWeights
1214 )
1215 {
1216  // Collect all constraints
1217 
1218  boolList blockedFace;
1219  PtrList<labelList> specifiedProcessorFaces;
1220  labelList specifiedProcessor;
1221  List<labelPair> explicitConnections;
1222  setConstraints
1223  (
1224  mesh,
1225  blockedFace,
1226  specifiedProcessorFaces,
1227  specifiedProcessor,
1228  explicitConnections
1229  );
1230 
1231 
1232  // Construct decomposition method and either do decomposition on
1233  // cell centres or on agglomeration
1234 
1235  labelList finalDecomp = decompose
1236  (
1237  mesh,
1238  cellWeights, // optional weights
1239  blockedFace, // any cells to be combined
1240  specifiedProcessorFaces,// any whole cluster of cells to be kept
1241  specifiedProcessor,
1242  explicitConnections // baffles
1243  );
1244 
1245 
1246  // Give any constraint the option of modifying the decomposition
1247 
1248  applyConstraints
1249  (
1250  mesh,
1251  blockedFace,
1252  specifiedProcessorFaces,
1253  specifiedProcessor,
1254  explicitConnections,
1255  finalDecomp
1256  );
1257 
1258  return finalDecomp;
1259 }
1260 
1261 
1262 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
void setSize(const label mRows)
Reset size of CompactListList.
Wave propagation of information through grid. Every iteration information goes through one layer of c...
Definition: FaceCellWave.H:79
const TrackingData & data() const
Additional data to be passed into container.
Definition: FaceCellWave.H:308
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
void clear()
Clear all entries from table.
Definition: HashTable.C:542
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
void setSize(const label)
Reset size of List.
Definition: List.C:281
const Type & second() const
Return second.
Definition: Pair.H:110
const Type & first() const
Return first.
Definition: Pair.H:98
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
void clear()
Clear the PtrList, i.e. set size to zero deleting all the.
Definition: PtrList.C:198
static const word & system()
Return system name.
Definition: TimePaths.H:112
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
List< Container > list() const
Convert to List<Container>
const UList< label > & offsets() const
Return the offset table (= size()+1)
const UList< T > & m() const
Return the packed matrix of data.
label size() const
Return the primary size, i.e. the number of rows.
T * data()
Return a pointer to the first data element,.
Definition: UListI.H:149
T & last()
Return the last element of the list.
Definition: UListI.H:128
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:278
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:476
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static const Form max
Definition: VectorSpace.H:120
static autoPtr< decompositionConstraint > New(const dictionary &constraintsDict, const word &type)
Return a reference to the selected decompositionConstraint.
Abstract base class for decomposition.
decompositionMethod(const dictionary &decompositionDict)
Construct given the decomposition dictionary.
static IOdictionary decomposeParDict(const Time &time)
Read and return the decomposeParDict.
PtrList< decompositionConstraint > constraints_
Optional constraints.
label nWeights(const pointField &points, const scalarField &pointWeights) const
Return the number of weights per point.
virtual labelList decompose(const pointField &points, const scalarField &pointWeights)
Return for every coordinate the wanted processor number.
void setConstraints(const polyMesh &mesh, boolList &blockedFace, PtrList< labelList > &specifiedProcessorFaces, labelList &specifiedProcessor, List< labelPair > &explicitConnections)
Helper: extract constraints:
label checkWeights(const pointField &points, const scalarField &pointWeights) const
Check the weights against the points.
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.
static autoPtr< decompositionMethod > NewDistributor(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
static labelList scaleWeights(const scalarField &weights, label &nWeights, const bool distributed=true)
Convert the given scalar weights to labels.
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.
static autoPtr< decompositionMethod > NewDecomposer(const dictionary &decompositionDict)
Return a reference to the selected decomposition method.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
ITstream & lookupBackwardsCompatible(const wordList &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream, trying a list of keywords.
Definition: dictionary.C:721
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
bool open(const fileName &libName, const bool verbose=true)
Open the named library, optionally with warnings if problems occur.
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:64
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
label toLocal(const label i) const
From global to local on current processor.
Definition: globalIndexI.H:117
bool isLocal(const label i) const
Is on local processor.
Definition: globalIndexI.H:95
label nTotalCells() const
Return total number of cells in decomposed mesh.
For use with FaceCellWave. Transports minimum passive data.
Definition: minData.H:64
const word & name() const
Return name.
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1326
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1339
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1515
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:404
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1345
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:290
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:280
label nCells() const
const vectorField & faceCentres() const
const vectorField & cellCentres() const
label nInternalFaces() const
label nFaces() const
const scalarField & magFaceAreas() const
const labelListList & pointFaces() const
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
This class separates the mesh into distinct unconnected regions, each of which is then given a label ...
Definition: regionSplit.H:118
label nLocalRegions() const
Return local number of regions.
Definition: regionSplit.H:202
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:436
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const pointField & points
const fvPatchList & patches
bool valid(const PtrList< ModelType > &l)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dlLibraryTable libs
Table of loaded dynamic libraries.
Type gSum(const FieldField< Field, Type > &f)
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
defineTypeNameAndDebug(combustionModel, 0)
error FatalError
static const label labelMax
Definition: label.H:62
static const char nl
Definition: Ostream.H:266
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
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, &mergedCyclicPolyPatch::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:229
dictionary dict
Operator to apply a binary operation to a pair of lists.
Definition: ListOps.H:288