MESA_PD_LoadBalancing.cpp 27.5 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
//======================================================================================================================
//
//  This file is part of waLBerla. waLBerla is free software: you can
//  redistribute it and/or modify it under the terms of the GNU General Public
//  License as published by the Free Software Foundation, either version 3 of
//  the License, or (at your option) any later version.
//
//  waLBerla is distributed in the hope that it will be useful, but WITHOUT
//  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
//  FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
//  for more details.
//
//  You should have received a copy of the GNU General Public License along
//  with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>.
//
//! \file   MESA_PD_LoadBalancing.cpp
//! \author Sebastian Eibl <sebastian.eibl@fau.de>
//
//======================================================================================================================

#include "Accessor.h"
#include "check.h"
#include "Contact.h"
#include "CreateParticles.h"
#include "NodeTimings.h"
#include "Parameters.h"
#include "SelectProperty.h"
#include "sortParticleStorage.h"
#include "SQLProperties.h"

#include <mesa_pd/vtk/ParticleVtkOutput.h>

#include <mesa_pd/collision_detection/AnalyticContactDetection.h>
#include <mesa_pd/data/LinkedCells.h>
#include <mesa_pd/data/ParticleAccessor.h>
#include <mesa_pd/data/ParticleStorage.h>
#include <mesa_pd/data/ShapeStorage.h>
#include <mesa_pd/data/STLOverloads.h>
#include <mesa_pd/domain/BlockForestDataHandling.h>
#include <mesa_pd/domain/BlockForestDomain.h>
#include <mesa_pd/domain/InfoCollection.h>
#include <mesa_pd/kernel/AssocToBlock.h>
#include <mesa_pd/kernel/DoubleCast.h>
#include <mesa_pd/kernel/ExplicitEulerWithShape.h>
#include <mesa_pd/kernel/InsertParticleIntoLinkedCells.h>
#include <mesa_pd/kernel/ParticleSelector.h>
#include <mesa_pd/kernel/SpringDashpot.h>
#include <mesa_pd/mpi/ContactFilter.h>
#include <mesa_pd/mpi/ReduceProperty.h>
#include <mesa_pd/mpi/SyncNextNeighbors.h>
#include <mesa_pd/mpi/SyncNextNeighborsBlockForest.h>

#include <mesa_pd/mpi/notifications/ForceTorqueNotification.h>

#include <blockforest/BlockForest.h>
#include <blockforest/Initialization.h>
#include <blockforest/loadbalancing/DynamicCurve.h>
#include <blockforest/loadbalancing/DynamicParMetis.h>
#include <blockforest/loadbalancing/PODPhantomData.h>
#include <core/Abort.h>
#include <core/Environment.h>
#include <core/math/Random.h>
#include <core/mpi/Reduce.h>
#include <core/mpi/MPITextFile.h>
#include <core/grid_generator/SCIterator.h>
#include <core/logging/Logging.h>
#include <core/OpenMP.h>
#include <core/timing/Timer.h>
#include <core/waLBerlaBuildInfo.h>
#include <pe/amr/level_determination/MinMaxLevelDetermination.h>
#include <pe/amr/weight_assignment/MetisAssignmentFunctor.h>
#include <pe/amr/weight_assignment/WeightAssignmentFunctor.h>
Sebastian Eibl's avatar
Sebastian Eibl committed
73
#include <sqlite/SQLite.h>
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
#include <vtk/VTKOutput.h>

#include <functional>
#include <memory>
#include <string>
#include <type_traits>

namespace walberla {
namespace mesa_pd {

int main( int argc, char ** argv )
{
   using namespace walberla::timing;

   Environment env(argc, argv);
   auto mpiManager = walberla::mpi::MPIManager::instance();
   mpiManager->useWorldComm();

92
93
   WALBERLA_LOG_DEVEL_ON_ROOT("MESA_PD_LoadBalancing" );

Sebastian Eibl's avatar
Sebastian Eibl committed
94
95
96
   //   logging::Logging::instance()->setStreamLogLevel(logging::Logging::INFO);
   //   logging::Logging::instance()->includeLoggingToFile("LoadBalancing");
   //   logging::Logging::instance()->setFileLogLevel(logging::Logging::DETAIL);
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123

   WALBERLA_LOG_INFO_ON_ROOT( "config file: " << argv[1] );
   WALBERLA_LOG_INFO_ON_ROOT( "waLBerla Revision: " << WALBERLA_GIT_SHA1 );

   math::seedRandomGenerator( static_cast<unsigned int>(1337 * mpiManager->worldRank()) );

   std::map< std::string, walberla::int64_t > integerProperties;
   std::map< std::string, double >            realProperties;
   std::map< std::string, std::string >       stringProperties;

   WALBERLA_LOG_INFO_ON_ROOT("*** READING CONFIG FILE ***");
   auto cfg = env.config();
   if (cfg == nullptr) WALBERLA_ABORT("No config specified!");
   const Config::BlockHandle mainConf  = cfg->getBlock( "GranularGas" );

   Parameters params;
   loadFromConfig(params, mainConf);

   WALBERLA_LOG_INFO_ON_ROOT("*** BLOCKFOREST ***");
   // create forest
   auto forest = blockforest::createBlockForestFromConfig( mainConf );
   if (!forest)
   {
      WALBERLA_LOG_INFO_ON_ROOT( "No BlockForest created ... exiting!");
      return EXIT_SUCCESS;
   }

124
125
126
127
   forest->recalculateBlockLevelsInRefresh( params.recalculateBlockLevelsInRefresh );
   forest->alwaysRebalanceInRefresh( params.alwaysRebalanceInRefresh );
   forest->reevaluateMinTargetLevelsAfterForcedRefinement( params.reevaluateMinTargetLevelsAfterForcedRefinement );
   forest->allowRefreshChangingDepth( params.allowRefreshChangingDepth );
128

129
130
131
   forest->allowMultipleRefreshCycles( params.allowMultipleRefreshCycles );
   forest->checkForEarlyOutInRefresh( params.checkForEarlyOutInRefresh );
   forest->checkForLateOutInRefresh( params.checkForLateOutInRefresh );
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179

   auto ic = make_shared<pe::InfoCollection>();

   //   pe::amr::MinMaxLevelDetermination regrid(ic, regridMin, regridMax);
   //   forest->setRefreshMinTargetLevelDeterminationFunction( regrid );

   bool bRebalance = true;
   if (params.LBAlgorithm == "None")
   {
      bRebalance = false;
   } else if (params.LBAlgorithm == "Morton")
   {
      forest->setRefreshPhantomBlockDataAssignmentFunction( pe::amr::WeightAssignmentFunctor( ic, params.baseWeight ) );
      forest->setRefreshPhantomBlockDataPackFunction( pe::amr::WeightAssignmentFunctor::PhantomBlockWeightPackUnpackFunctor() );
      forest->setRefreshPhantomBlockDataUnpackFunction( pe::amr::WeightAssignmentFunctor::PhantomBlockWeightPackUnpackFunctor() );

      auto prepFunc = blockforest::DynamicCurveBalance< pe::amr::WeightAssignmentFunctor::PhantomBlockWeight >( false, true, false );
      prepFunc.setMaxBlocksPerProcess( params.maxBlocksPerProcess );
      forest->setRefreshPhantomBlockMigrationPreparationFunction( prepFunc );
   } else if (params.LBAlgorithm == "Hilbert")
   {
      forest->setRefreshPhantomBlockDataAssignmentFunction( pe::amr::WeightAssignmentFunctor( ic, params.baseWeight ) );
      forest->setRefreshPhantomBlockDataPackFunction( pe::amr::WeightAssignmentFunctor::PhantomBlockWeightPackUnpackFunctor() );
      forest->setRefreshPhantomBlockDataUnpackFunction( pe::amr::WeightAssignmentFunctor::PhantomBlockWeightPackUnpackFunctor() );

      auto prepFunc = blockforest::DynamicCurveBalance< pe::amr::WeightAssignmentFunctor::PhantomBlockWeight >( true, true, false );
      prepFunc.setMaxBlocksPerProcess( params.maxBlocksPerProcess );
      forest->setRefreshPhantomBlockMigrationPreparationFunction( prepFunc );
   } else if (params.LBAlgorithm == "Metis")
   {
      auto assFunc = pe::amr::MetisAssignmentFunctor( ic, params.baseWeight );
      forest->setRefreshPhantomBlockDataAssignmentFunction( assFunc );
      forest->setRefreshPhantomBlockDataPackFunction( pe::amr::MetisAssignmentFunctor::PhantomBlockWeightPackUnpackFunctor() );
      forest->setRefreshPhantomBlockDataUnpackFunction( pe::amr::MetisAssignmentFunctor::PhantomBlockWeightPackUnpackFunctor() );

      auto alg     = blockforest::DynamicParMetis::stringToAlgorithm(    params.metisAlgorithm );
      auto vWeight = blockforest::DynamicParMetis::stringToWeightsToUse( params.metisWeightsToUse );
      auto eWeight = blockforest::DynamicParMetis::stringToEdgeSource(   params.metisEdgeSource );

      auto prepFunc = blockforest::DynamicParMetis( alg, vWeight, eWeight );
      prepFunc.setipc2redist(params.metisipc2redist);
      addParMetisPropertiesToSQL(prepFunc, integerProperties, realProperties, stringProperties);
      forest->setRefreshPhantomBlockMigrationPreparationFunction( prepFunc );
   } else if (params.LBAlgorithm == "Diffusive")
   {
      forest->setRefreshPhantomBlockDataAssignmentFunction( pe::amr::WeightAssignmentFunctor( ic, params.baseWeight ) );
      forest->setRefreshPhantomBlockDataPackFunction( pe::amr::WeightAssignmentFunctor::PhantomBlockWeightPackUnpackFunctor() );
      forest->setRefreshPhantomBlockDataUnpackFunction( pe::amr::WeightAssignmentFunctor::PhantomBlockWeightPackUnpackFunctor() );
180
      auto prepFunc = blockforest::DynamicDiffusionBalance< pe::amr::WeightAssignmentFunctor::PhantomBlockWeight >( 1, 1, false );
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
      //configure(cfg, prepFunc);
      //addDynamicDiffusivePropertiesToSQL(prepFunc, integerProperties, realProperties, stringProperties);
      forest->setRefreshPhantomBlockMigrationPreparationFunction(prepFunc);
   } else
   {
      WALBERLA_ABORT("Unknown LBAlgorithm: " << params.LBAlgorithm);
   }

   auto domain = std::make_shared<domain::BlockForestDomain>(forest);

   WALBERLA_LOG_INFO_ON_ROOT("*** SETUP - START ***");

   //init data structures
   auto ps = std::make_shared<data::ParticleStorage>(100);
   auto ss = std::make_shared<data::ShapeStorage>();
   ParticleAccessorWithShape accessor(ps, ss);
   auto lc = std::make_shared<data::LinkedCells>(domain->getUnionOfLocalAABBs().getExtended(params.spacing), params.spacing );
   forest->addBlockData(domain::createBlockForestDataHandling(ps), "Storage");

   auto center = forest->getDomain().center();
   auto  smallSphere = ss->create<data::Sphere>( params.radius );
   ss->shapes[smallSphere]->updateMassAndInertia(real_t(2707));
   for (auto& iBlk : *forest)
   {
      for (auto pt : grid_generator::SCGrid(iBlk.getAABB(),
206
                                            Vector3<real_t>(params.spacing) * real_c(0.5) + params.shift,
207
208
209
                                            params.spacing))
      {
         WALBERLA_CHECK(iBlk.getAABB().contains(pt));
Sebastian Eibl's avatar
Sebastian Eibl committed
210
         auto tmp = dot( (pt - center), params.normal );
211
212
213
214
215
216
217
218
219
220
221
222
223
         if (tmp < 0)
         {
            createSphere(*ps, pt, params.radius, smallSphere);
         }
      }
   }
   int64_t numParticles = int64_c(ps->size());
   walberla::mpi::reduceInplace(numParticles, walberla::mpi::SUM);
   WALBERLA_LOG_INFO_ON_ROOT("#particles created: " << numParticles);

   WALBERLA_LOG_INFO_ON_ROOT("*** SETUP - END ***");

   WALBERLA_LOG_INFO_ON_ROOT("*** VTK ***");
Sebastian Eibl's avatar
Sebastian Eibl committed
224
   auto vtkDomainOutput = walberla::vtk::createVTKOutput_DomainDecomposition( forest, "domain_decomposition", 1, params.vtk_out, "simulation_step" );
225
   auto vtkOutput       = make_shared<mesa_pd::vtk::ParticleVtkOutput>(ps) ;
Sebastian Eibl's avatar
Sebastian Eibl committed
226
   auto vtkWriter       = walberla::vtk::createVTKOutput_PointData(vtkOutput, "Bodies", 1, params.vtk_out, "simulation_step", false, false);
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
   vtkOutput->addOutput<SelectRank>("rank");
   vtkOutput->addOutput<data::SelectParticleOwner>("owner");
   vtkDomainOutput->write();

   WALBERLA_LOG_INFO_ON_ROOT("*** SIMULATION - START ***");
   // Init kernels
   kernel::ExplicitEulerWithShape        explicitEulerWithShape( params.dt );
   kernel::InsertParticleIntoLinkedCells ipilc;
   kernel::SpringDashpot                 dem(1);
   dem.setStiffness(0, 0, real_t(0));
   dem.setDampingN (0, 0, real_t(0));
   dem.setDampingT (0, 0, real_t(0));
   dem.setFriction (0, 0, real_t(0));
   collision_detection::AnalyticContactDetection              acd;
   kernel::AssocToBlock                  assoc(forest);
   kernel::DoubleCast                    double_cast;
   mpi::ContactFilter                    contact_filter;
   mpi::ReduceProperty                   RP;
   mpi::SyncNextNeighborsBlockForest     SNN;

   // initial sync
   ps->forEachParticle(true, kernel::SelectLocal(), accessor, assoc, accessor);
   SNN(*ps, forest, domain);
   sortParticleStorage(*ps, params.sorting, lc->domain_, uint_c(lc->numCellsPerDim_[0]));

Sebastian Eibl's avatar
Sebastian Eibl committed
252
253
254
255
256
   WcTimer      timerImbalanced;
   WcTimer      timerLoadBalancing;
   WcTimer      timerBalanced;
   WcTimingPool tpImbalanced;
   WcTimingPool tpBalanced;
257

Sebastian Eibl's avatar
Sebastian Eibl committed
258
259
260
261
262
263
264
265
   auto    SNNBytesSent     = SNN.getBytesSent();
   auto    SNNBytesReceived = SNN.getBytesReceived();
   auto    SNNSends         = SNN.getNumberOfSends();
   auto    SNNReceives      = SNN.getNumberOfReceives();
   auto    RPBytesSent      = RP.getBytesSent();
   auto    RPBytesReceived  = RP.getBytesReceived();
   auto    RPSends          = RP.getNumberOfSends();
   auto    RPReceives       = RP.getNumberOfReceives();
266
267
268
   int64_t imbalancedContactsChecked  = 0;
   int64_t imbalancedContactsDetected = 0;
   int64_t imbalancedContactsTreated  = 0;
Sebastian Eibl's avatar
Sebastian Eibl committed
269

270
   WALBERLA_MPI_BARRIER();
Sebastian Eibl's avatar
Sebastian Eibl committed
271
272
273
   WALBERLA_LOG_DEVEL_ON_ROOT("running imbalanced simulation");
   timerImbalanced.start();
   for (int64_t i=0; i < params.simulationSteps; ++i)
274
   {
Sebastian Eibl's avatar
Sebastian Eibl committed
275
276
277
278
279
280
281
282
      //WALBERLA_LOG_DEVEL_ON_ROOT("timestep: " << i << " / " << params.simulationSteps );
      //         if (i % params.visSpacing == 0)
      //         {
      //            vtkWriter->write();
      //         }

      tpImbalanced["AssocToBlock"].start();
      ps->forEachParticle(true, kernel::SelectLocal(), accessor, assoc, accessor);
283
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
Sebastian Eibl's avatar
Sebastian Eibl committed
284
285
286
287
288
289
290
291
292
      tpImbalanced["AssocToBlock"].end();

      tpImbalanced["GenerateLinkedCells"].start();
      lc->clear();
      ps->forEachParticle(true, kernel::SelectAll(), accessor, ipilc, accessor, *lc);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpImbalanced["GenerateLinkedCells"].end();

      tpImbalanced["DEM"].start();
293
294
295
      imbalancedContactsChecked  = 0;
      imbalancedContactsDetected = 0;
      imbalancedContactsTreated  = 0;
Sebastian Eibl's avatar
Sebastian Eibl committed
296
297
298
299
      lc->forEachParticlePairHalf(true,
                                  kernel::SelectAll(),
                                  accessor,
                                  [&](const size_t idx1, const size_t idx2, auto& ac)
300
      {
301
         ++imbalancedContactsChecked;
Sebastian Eibl's avatar
Sebastian Eibl committed
302
         if (double_cast(idx1, idx2, ac, acd, ac ))
303
         {
304
            ++imbalancedContactsDetected;
Sebastian Eibl's avatar
Sebastian Eibl committed
305
            if (contact_filter(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), *domain))
306
            {
307
               ++imbalancedContactsTreated;
Sebastian Eibl's avatar
Sebastian Eibl committed
308
               dem(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), acd.getContactNormal(), acd.getPenetrationDepth());
309
            }
Sebastian Eibl's avatar
Sebastian Eibl committed
310
311
312
313
314
         }
      },
      accessor );
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpImbalanced["DEM"].end();
315

Sebastian Eibl's avatar
Sebastian Eibl committed
316
317
318
319
      tpImbalanced["ReduceForce"].start();
      RP.operator()<ForceTorqueNotification>(*ps);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpImbalanced["ReduceForce"].end();
320

Sebastian Eibl's avatar
Sebastian Eibl committed
321
322
323
324
325
326
327
328
329
330
331
332
      tpImbalanced["Euler"].start();
      //ps->forEachParticle(false, [&](const size_t idx){WALBERLA_CHECK_EQUAL(ps->getForce(idx), Vec3(0,0,0), *(*ps)[idx] << "\n" << idx);});
      ps->forEachParticle(true, kernel::SelectLocal(), accessor, explicitEulerWithShape, accessor);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpImbalanced["Euler"].end();

      tpImbalanced["SNN"].start();
      SNN(*ps, forest, domain);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpImbalanced["SNN"].end();
   }
   timerImbalanced.end();
333

334
   WALBERLA_MPI_BARRIER();
Sebastian Eibl's avatar
Sebastian Eibl committed
335
336
337
338
339
340
   timerLoadBalancing.start();
   if (bRebalance)
   {
      WALBERLA_LOG_DEVEL_ON_ROOT("running load balancing");
      domain::createWithNeighborhood( accessor, *forest, *ic );
      for (auto pIt = ps->begin(); pIt != ps->end(); )
341
      {
Sebastian Eibl's avatar
Sebastian Eibl committed
342
343
         using namespace walberla::mesa_pd::data::particle_flags;
         if (isSet(pIt->getFlags(), GHOST))
344
         {
Sebastian Eibl's avatar
Sebastian Eibl committed
345
            pIt = ps->erase(pIt);
346
347
         } else
         {
Sebastian Eibl's avatar
Sebastian Eibl committed
348
349
            pIt->getGhostOwnersRef().clear();
            ++pIt;
350
351
         }
      }
Sebastian Eibl's avatar
Sebastian Eibl committed
352
353
354
355
356
357
358
359
360
361
362
      forest->refresh();
      domain->refresh();
      lc = std::make_shared<data::LinkedCells>(domain->getUnionOfLocalAABBs().getExtended(params.spacing), params.spacing );
      ps->forEachParticle(true, kernel::SelectLocal(), accessor, assoc, accessor);
      SNN(*ps, forest, domain);
      sortParticleStorage(*ps, params.sorting, lc->domain_, uint_c(lc->numCellsPerDim_[0]));
   }
   timerLoadBalancing.end();

   WALBERLA_MPI_BARRIER();
   WALBERLA_LOG_DEVEL_ON_ROOT("running balanced simulation");
363
364
365
   int64_t balancedContactsChecked  = 0;
   int64_t balancedContactsDetected = 0;
   int64_t balancedContactsTreated  = 0;
Sebastian Eibl's avatar
Sebastian Eibl committed
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
   timerBalanced.start();
   for (int64_t i=0; i < params.simulationSteps; ++i)
   {
      //WALBERLA_LOG_DEVEL_ON_ROOT("timestep: " << i << " / " << params.simulationSteps );
      //         if (i % params.visSpacing == 0)
      //         {
      //            vtkWriter->write();
      //         }

      tpBalanced["AssocToBlock"].start();
      ps->forEachParticle(true, kernel::SelectLocal(), accessor, assoc, accessor);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpBalanced["AssocToBlock"].end();

      tpBalanced["GenerateLinkedCells"].start();
      lc->clear();
      ps->forEachParticle(true, kernel::SelectAll(), accessor, ipilc, accessor, *lc);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpBalanced["GenerateLinkedCells"].end();

      tpBalanced["DEM"].start();
387
388
389
      balancedContactsChecked  = 0;
      balancedContactsDetected = 0;
      balancedContactsTreated  = 0;
Sebastian Eibl's avatar
Sebastian Eibl committed
390
391
392
393
394
      lc->forEachParticlePairHalf(true,
                                  kernel::SelectAll(),
                                  accessor,
                                  [&](const size_t idx1, const size_t idx2, auto& ac)
      {
395
         ++balancedContactsChecked;
Sebastian Eibl's avatar
Sebastian Eibl committed
396
397
         if (double_cast(idx1, idx2, ac, acd, ac ))
         {
398
            ++balancedContactsDetected;
Sebastian Eibl's avatar
Sebastian Eibl committed
399
400
            if (contact_filter(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), *domain))
            {
401
               ++balancedContactsTreated;
Sebastian Eibl's avatar
Sebastian Eibl committed
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
               dem(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), acd.getContactNormal(), acd.getPenetrationDepth());
            }
         }
      },
      accessor );
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpBalanced["DEM"].end();

      tpBalanced["ReduceForce"].start();
      RP.operator()<ForceTorqueNotification>(*ps);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpBalanced["ReduceForce"].end();

      tpBalanced["Euler"].start();
      //ps->forEachParticle(false, [&](const size_t idx){WALBERLA_CHECK_EQUAL(ps->getForce(idx), Vec3(0,0,0), *(*ps)[idx] << "\n" << idx);});
      ps->forEachParticle(true, kernel::SelectLocal(), accessor, explicitEulerWithShape, accessor);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpBalanced["Euler"].end();

      tpBalanced["SNN"].start();
      SNN(*ps, forest, domain);
      if (params.bBarrier) WALBERLA_MPI_BARRIER();
      tpBalanced["SNN"].end();
   }
   timerBalanced.end();

   SNNBytesSent     = SNN.getBytesSent();
   SNNBytesReceived = SNN.getBytesReceived();
   SNNSends         = SNN.getNumberOfSends();
   SNNReceives      = SNN.getNumberOfReceives();
   RPBytesSent      = RP.getBytesSent();
   RPBytesReceived  = RP.getBytesReceived();
   RPSends          = RP.getNumberOfSends();
   RPReceives       = RP.getNumberOfReceives();
   walberla::mpi::reduceInplace(SNNBytesSent, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(SNNBytesReceived, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(SNNSends, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(SNNReceives, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(RPBytesSent, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(RPBytesReceived, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(RPSends, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(RPReceives, walberla::mpi::SUM);
444
445
446
   auto cC = walberla::mpi::reduce(balancedContactsChecked, walberla::mpi::SUM);
   auto cD = walberla::mpi::reduce(balancedContactsDetected, walberla::mpi::SUM);
   auto cT = walberla::mpi::reduce(balancedContactsTreated, walberla::mpi::SUM);
Sebastian Eibl's avatar
Sebastian Eibl committed
447
448
449
450
451
452
453
454
455
456
   WALBERLA_LOG_DEVEL_ON_ROOT( "SNN bytes communicated:   " << SNNBytesSent << " / " << SNNBytesReceived );
   WALBERLA_LOG_DEVEL_ON_ROOT( "SNN communication partners: " << SNNSends << " / " << SNNReceives );
   WALBERLA_LOG_DEVEL_ON_ROOT( "RP bytes communicated:  " << RPBytesSent << " / " << RPBytesReceived );
   WALBERLA_LOG_DEVEL_ON_ROOT( "RP communication partners: " << RPSends << " / " << RPReceives );
   WALBERLA_LOG_DEVEL_ON_ROOT( "contacts checked/detected/treated: " << cC << " / " << cD << " / " << cT );

   auto timerImbalancedReduced = walberla::timing::getReduced(timerImbalanced, REDUCE_TOTAL, 0);
   double PUpSImbalanced = 0.0;
   WALBERLA_ROOT_SECTION()
   {
457
      WALBERLA_LOG_INFO_ON_ROOT("IMBALANCED " << *timerImbalancedReduced);
Sebastian Eibl's avatar
Sebastian Eibl committed
458
459
460
461
462
463
464
465
      PUpSImbalanced = double_c(numParticles) * double_c(params.simulationSteps) / double_c(timerImbalancedReduced->max());
      WALBERLA_LOG_INFO_ON_ROOT("PUpS: " << PUpSImbalanced);
   }

   auto timerBalancedReduced = walberla::timing::getReduced(timerBalanced, REDUCE_TOTAL, 0);
   double PUpSBalanced = 0.0;
   WALBERLA_ROOT_SECTION()
   {
466
      WALBERLA_LOG_INFO_ON_ROOT("BALANCED " << *timerBalancedReduced);
Sebastian Eibl's avatar
Sebastian Eibl committed
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
      PUpSBalanced = double_c(numParticles) * double_c(params.simulationSteps) / double_c(timerBalancedReduced->max());
      WALBERLA_LOG_INFO_ON_ROOT("PUpS: " << PUpSBalanced);
   }

   auto timerLoadBalancingReduced = walberla::timing::getReduced(timerLoadBalancing, REDUCE_TOTAL, 0);

   auto tpImbalancedReduced = tpImbalanced.getReduced();
   WALBERLA_LOG_INFO_ON_ROOT(*tpImbalancedReduced);

   auto tpBalancedReduced = tpBalanced.getReduced();
   WALBERLA_LOG_INFO_ON_ROOT(*tpBalancedReduced);
   WALBERLA_LOG_INFO_ON_ROOT("*** SIMULATION - END ***");

   if (params.checkSimulation)
   {
      check(*ps, *forest, params.spacing);
   }
484

Sebastian Eibl's avatar
Sebastian Eibl committed
485
486
487
488
489
490
491
492
493
494
495
496
   WALBERLA_LOG_INFO_ON_ROOT("*** SQL OUTPUT - START ***");
   numParticles = 0;
   int64_t numGhostParticles = 0;
   ps->forEachParticle(false,
                       kernel::SelectAll(),
                       accessor,
                       [&numParticles, &numGhostParticles](const size_t idx, auto& ac)
   {
      if (data::particle_flags::isSet( ac.getFlagsRef(idx), data::particle_flags::GHOST))
      {
         ++numGhostParticles;
      } else
497
      {
Sebastian Eibl's avatar
Sebastian Eibl committed
498
         ++numParticles;
499
      }
Sebastian Eibl's avatar
Sebastian Eibl committed
500
501
   },
   accessor);
502
503
504
   auto minParticles = walberla::mpi::reduce(numParticles, walberla::mpi::MIN);
   auto maxParticles = walberla::mpi::reduce(numParticles, walberla::mpi::MAX);
   WALBERLA_LOG_DEVEL_ON_ROOT("particle ratio: " << minParticles << " / " << maxParticles);
Sebastian Eibl's avatar
Sebastian Eibl committed
505
506
   walberla::mpi::reduceInplace(numParticles, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(numGhostParticles, walberla::mpi::SUM);
507
508
509
510
511
512
   walberla::mpi::reduceInplace(imbalancedContactsChecked, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(imbalancedContactsDetected, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(imbalancedContactsTreated, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(balancedContactsChecked, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(balancedContactsDetected, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(balancedContactsTreated, walberla::mpi::SUM);
Sebastian Eibl's avatar
Sebastian Eibl committed
513
514
515
516
517
518
519
520
521
522
523
524
525
526
   double linkedCellsVolume = lc->domain_.volume();
   walberla::mpi::reduceInplace(linkedCellsVolume, walberla::mpi::SUM);
   size_t numLinkedCells = lc->cells_.size();
   walberla::mpi::reduceInplace(numLinkedCells, walberla::mpi::SUM);
   size_t local_aabbs         = domain->getNumLocalAABBs();
   size_t neighbor_subdomains = domain->getNumNeighborSubdomains();
   size_t neighbor_processes  = domain->getNumNeighborProcesses();
   walberla::mpi::reduceInplace(local_aabbs, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(neighbor_subdomains, walberla::mpi::SUM);
   walberla::mpi::reduceInplace(neighbor_processes, walberla::mpi::SUM);

   uint_t runId = uint_c(-1);
   WALBERLA_ROOT_SECTION()
   {
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
      stringProperties["walberla_git"]                  = WALBERLA_GIT_SHA1;
      stringProperties["tag"]                           = "mesa_pd";
      integerProperties["mpi_num_processes"]            = mpiManager->numProcesses();
      integerProperties["omp_max_threads"]              = omp_get_max_threads();
      realProperties["imbalanced_PUpS"]                 = double_c(PUpSImbalanced);
      realProperties["imbalanced_timer_min"]            = timerImbalancedReduced->min();
      realProperties["imbalanced_timer_max"]            = timerImbalancedReduced->max();
      realProperties["imbalanced_timer_average"]        = timerImbalancedReduced->average();
      realProperties["imbalanced_timer_total"]          = timerImbalancedReduced->total();
      realProperties["loadbalancing_timer_min"]         = timerLoadBalancingReduced->min();
      realProperties["loadbalancing_timer_max"]         = timerLoadBalancingReduced->max();
      realProperties["loadbalancing_timer_average"]     = timerLoadBalancingReduced->average();
      realProperties["loadbalancing_timer_total"]       = timerLoadBalancingReduced->total();
      realProperties["balanced_PUpS"]                   = double_c(PUpSBalanced);
      realProperties["balanced_timer_min"]              = timerBalancedReduced->min();
      realProperties["balanced_timer_max"]              = timerBalancedReduced->max();
      realProperties["balanced_timer_average"]          = timerBalancedReduced->average();
      realProperties["balanced_timer_total"]            = timerBalancedReduced->total();
      integerProperties["num_particles"]                = numParticles;
      integerProperties["num_ghost_particles"]          = numGhostParticles;
      integerProperties["minParticles"]                 = minParticles;
      integerProperties["maxParticles"]                 = maxParticles;
      integerProperties["imbalancedContactsChecked"]    = imbalancedContactsChecked;
      integerProperties["imbalancedContactsDetected"]   = imbalancedContactsDetected;
      integerProperties["imbalancedContactsTreated"]    = imbalancedContactsTreated;
      integerProperties["balancedContactsChecked"]      = balancedContactsChecked;
      integerProperties["balancedContactsDetected"]     = balancedContactsDetected;
      integerProperties["balancedContactsTreated"]      = balancedContactsTreated;
      integerProperties["local_aabbs"]                  = int64_c(local_aabbs);
      integerProperties["neighbor_subdomains"]          = int64_c(neighbor_subdomains);
      integerProperties["neighbor_processes"]           = int64_c(neighbor_processes);
      integerProperties["SNNBytesSent"]                 = SNNBytesSent;
      integerProperties["SNNBytesReceived"]             = SNNBytesReceived;
      integerProperties["SNNSends"]                     = SNNSends;
      integerProperties["SNNReceives"]                  = SNNReceives;
      integerProperties["RPBytesSent"]                  = RPBytesSent;
      integerProperties["RPBytesReceived"]              = RPBytesReceived;
      integerProperties["RPSends"]                      = RPSends;
      integerProperties["RPReceives"]                   = RPReceives;
      realProperties["linkedCellsVolume"]               = linkedCellsVolume;
      integerProperties["numLinkedCells"]               = int64_c(numLinkedCells);
Sebastian Eibl's avatar
Sebastian Eibl committed
568
569
570
571
572
573

      addBuildInfoToSQL( integerProperties, realProperties, stringProperties );
      saveToSQL(params, integerProperties, realProperties, stringProperties );
      addDomainPropertiesToSQL(*forest, integerProperties, realProperties, stringProperties);
      addSlurmPropertiesToSQL(integerProperties, realProperties, stringProperties);

Sebastian Eibl's avatar
Sebastian Eibl committed
574
575
576
      runId = sqlite::storeRunInSqliteDB( params.sqlFile, integerProperties, stringProperties, realProperties );
      sqlite::storeTimingPoolInSqliteDB( params.sqlFile, runId, *tpImbalancedReduced, "imbalanced" );
      sqlite::storeTimingPoolInSqliteDB( params.sqlFile, runId, *tpImbalancedReduced, "balanced" );
Sebastian Eibl's avatar
Sebastian Eibl committed
577
578
579
580
581
582
   }

   if (params.storeNodeTimings)
   {
      storeNodeTimings(runId, params.sqlFile, "NodeTimingImbalanced", tpImbalanced);
      storeNodeTimings(runId, params.sqlFile, "NodeTimingBalanced", tpBalanced);
583
   }
Sebastian Eibl's avatar
Sebastian Eibl committed
584
   WALBERLA_LOG_INFO_ON_ROOT("*** SQL OUTPUT - END ***");
585
586
587
588
589
590
591
592
593
594
595

   return EXIT_SUCCESS;
}

} // namespace mesa_pd
} // namespace walberla

int main( int argc, char* argv[] )
{
   return walberla::mesa_pd::main( argc, argv );
}