fvMeshDistributorsLoadBalancer.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) 2021-2022 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 
27 #include "decompositionMethod.H"
28 #include "cpuLoad.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace fvMeshDistributors
36 {
39  (
42  fvMesh
43  );
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
50 void Foam::fvMeshDistributors::loadBalancer::readDict()
51 {
53 
54  const dictionary& distributorDict(dict());
55 
56  multiConstraint_ =
57  distributorDict.lookupOrDefault<Switch>("multiConstraint", true);
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
62 
64 :
65  distributor(mesh)
66 {
67  readDict();
68 }
69 
70 
71 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
72 
74 {}
75 
76 
77 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78 
80 {
81  const fvMesh& mesh = this->mesh();
82 
83  bool redistributed = false;
84 
85  if
86  (
87  Pstream::nProcs() > 1
88  && mesh.time().timeIndex() > 1
89  && timeIndex_ != mesh.time().timeIndex()
90  )
91  {
92  timeIndex_ = mesh.time().timeIndex();
93 
94  const scalar timeStepCpuTime = cpuTime_.cpuTimeIncrement();
95 
96  // CPU loads per cell
97  HashTable<cpuLoad*> cpuLoads(this->mesh().lookupClass<cpuLoad>());
98 
99  if (!cpuLoads.size())
100  {
102  << "No CPU loads have been allocated"
103  << exit(FatalError);
104  }
105 
106  if (mesh.time().timeIndex() % redistributionInterval_ == 0)
107  {
108  timeIndex_ = mesh.time().timeIndex();
109 
110  scalar sumCpuLoad = 0;
111 
112  forAllConstIter(HashTable<cpuLoad*>, cpuLoads, iter)
113  {
114  sumCpuLoad += sum(iter()->field());
115  }
116 
117  const scalar cellCFDCpuTime = returnReduce
118  (
119  (timeStepCpuTime - sumCpuLoad)/mesh.nCells(),
120  minOp<scalar>()
121  );
122 
123  // Total CPU time for this processor
124  const scalar processorCpuTime =
125  mesh.nCells()*cellCFDCpuTime + sumCpuLoad;
126 
127  // Average processor CPU time
128  const scalar averageProcessorCpuTime =
129  returnReduce(processorCpuTime, sumOp<scalar>())
130  /Pstream::nProcs();
131 
132  Pout<< "imbalance "
133  << " " << sumCpuLoad
134  << " " << mesh.nCells()*cellCFDCpuTime
135  << " " << processorCpuTime
136  << " " << averageProcessorCpuTime << endl;
137 
138  const scalar imbalance = returnReduce
139  (
140  mag(1 - processorCpuTime/averageProcessorCpuTime),
141  maxOp<scalar>()
142  );
143 
144  scalarField weights;
145 
146  if (multiConstraint_)
147  {
148  const int nWeights = cpuLoads.size() + 1;
149 
150  weights.setSize(nWeights*mesh.nCells());
151 
152  for (label i=0; i<mesh.nCells(); i++)
153  {
154  weights[nWeights*i] = cellCFDCpuTime;
155  }
156 
157  label loadi = 1;
158  forAllConstIter(HashTable<cpuLoad*>, cpuLoads, iter)
159  {
160  const scalarField& cpuLoadField = iter()->field();
161 
162  forAll(cpuLoadField, i)
163  {
164  weights[nWeights*i + loadi] = cpuLoadField[i];
165  }
166 
167  loadi++;
168  }
169  }
170  else
171  {
172  weights.setSize(mesh.nCells(), cellCFDCpuTime);
173 
174  forAllConstIter(HashTable<cpuLoad*>, cpuLoads, iter)
175  {
176  weights += iter()->field();
177  }
178  }
179 
180  if (imbalance > maxImbalance_)
181  {
182  Info<< "Redistributing mesh with imbalance "
183  << imbalance << endl;
184 
185  // Create new decomposition distribution
186  const labelList distribution
187  (
188  distributor_->decompose(mesh, weights)
189  );
190 
191  distribute(distribution);
192 
193  redistributed = true;
194  }
195  }
196 
197  forAllIter(HashTable<cpuLoad*>, cpuLoads, iter)
198  {
199  iter()->checkOut();
200  }
201  }
202 
203  return redistributed;
204 }
205 
206 
207 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
An STL-conforming hash table.
Definition: HashTable.H:127
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
void setSize(const label)
Reset size of List.
Definition: List.C:281
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
Accumulating histogram of values. Specified bin resolution automatic generation of bins.
Definition: distribution.H:61
Abstract base class for fvMesh movers.
const dictionary & dict() const
Return the dynamicMeshDict/distributor sub-dict.
Dynamic mesh redistribution using the distributor specified in decomposeParDict.
void readDict()
Read the projection parameters from dictionary.
Dynamic mesh redistribution using the distributor specified in decomposeParDict.
loadBalancer(fvMesh &mesh)
Construct from fvMesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:406
label nCells() const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
addToRunTimeSelectionTable(fvMeshDistributor, none, fvMesh)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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:251
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError