NonuniformGridCodegenPerformance.cpp 30 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
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
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
124
125
126
127
128
129
130
131
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
180
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
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
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
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
//======================================================================================================================
//
//  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 NonuniformGridCodegenPerformance.cpp
//! \author Frederik Hennig <frederik.hennig@fau.de>
//
//======================================================================================================================


#include "blockforest/all.h"

#include "core/all.h"

#include "domain_decomposition/all.h"

#include "lbm/all.h"

#include "field/all.h"

#include "geometry/all.h"

#include "timeloop/all.h"

#include "vtk/all.h"

#include "python_coupling/CreateConfig.h"
#include "python_coupling/PythonCallback.h"
#include "python_coupling/DictWrapper.h"

#include "RefinementCodegen.h"


namespace walberla{

using LatticeModel_T         = lbm::GeneratedLatticeModel;
using Stencil_T              = LatticeModel_T::Stencil;
using CommunicationStencil_T = LatticeModel_T::CommunicationStencil;
using PdfField_T             = lbm::PdfField< LatticeModel_T >;

using FlagField_T            = field::FlagField< uint8_t >;
using ScalarField_T          = GhostLayerField< real_t, 1 >;
using VectorField_T          = GhostLayerField< real_t, Stencil_T::D >;

using VelocityEvalFunc       = std::function< Vector3< real_t > ( Vector3< real_t > ) >;
using VelocityCallback       = std::function< Vector3<real_t> (const Cell &, const shared_ptr< StructuredBlockForest >&, IBlock&)>;

using stencil::Direction;
using blockforest::communication::NonUniformBufferedScheme;



//////////////////////////////////////////////////////////////////
///    Workaround for In-Place Streaming in Boundaries         ///
//////////////////////////////////////////////////////////////////

/*
 * Wrapper class around the generated boundary, to extract the timestep from
 * the lattice model.
 */
template< typename LM_T, typename Boundary_T, bool inplace = false >
struct BoundarySweep
{
 public:
   BoundarySweep(Boundary_T&, BlockDataID) {};

   void operator()(IBlock*) = 0;
};

template< typename LM_T, typename Boundary_T >
struct BoundarySweep< LM_T, Boundary_T, false >
{
 public:
   BoundarySweep(Boundary_T& boundary, BlockDataID pdfFieldId) : boundary_(boundary), pdfFieldId_(pdfFieldId)
   {}

   void operator()(IBlock* block)
   {
      boundary_.run(block, uint8_t(0));
   }

 private:
   Boundary_T& boundary_;
   BlockDataID pdfFieldId_;
};

template< typename LM_T, typename Boundary_T >
struct BoundarySweep< LM_T, Boundary_T, true >
{
 public:
   BoundarySweep(Boundary_T& boundary, BlockDataID pdfFieldId) : boundary_(boundary), pdfFieldId_(pdfFieldId)
   {}

   void operator()(IBlock* block)
   {
      auto field       = block->getData< lbm::PdfField< LM_T > >(pdfFieldId_);
      uint8_t timestep = field->latticeModel().getTimestep();
      boundary_.run(block, timestep);
   }

 private:
   Boundary_T& boundary_;
   BlockDataID pdfFieldId_;
};


/////////////////////////////////////////
///         Grid and Flow Setup       ///
/////////////////////////////////////////

using RefinementSelectionFunctor = SetupBlockForest::RefinementSelectionFunction;

class SetupBase {
 public:
   virtual RefinementSelectionFunctor refinementSelector() = 0;
   virtual void setupBoundaryFlagField(StructuredBlockForest & /*sbfs*/, const BlockDataID /*flagFieldID*/) {};


   virtual real_t shearRelaxationRate() const  = 0;
   virtual Vector3< real_t > acceleration() const = 0;
   virtual void getBoundarySweeps(std::shared_ptr< StructuredBlockForest > & /*sbfs*/,
                                 const BlockDataID /*pdfFieldID*/,
                                 const BlockDataID /*flagFieldID*/,
                                 const FlagUID /*fluidFlagUID*/,
                                 std::vector< std::function< void( Block * ) > > & /*sweeps*/) {};
};

class AABBRefinement
{
 private:
   AABB refinementRegion_;
   uint_t refinementDepth_;

   uint_t targetLevel(const AABB & blockAABB){
      if(refinementRegion_.contains(blockAABB.center())){
         return refinementDepth_;
      }else{
         return 0;
      }
   }

 public:
   AABBRefinement(const AABB & refinementRegion, uint_t depth) : refinementRegion_(refinementRegion), refinementDepth_(depth) {}

   void operator() ( SetupBlockForest & forest ){
      std::vector< SetupBlock* > blocks;
      forest.getBlocks(blocks);

      for (auto b : blocks)
      {
         AABB aabb = b->getAABB();
         if (b->getLevel() < targetLevel(aabb)) b->setMarker(true);
      }
   }
};

class RefineEverywhere
{
 private:
   uint_t refinementDepth_;

   uint_t targetLevel(const AABB & /*blockAABB*/){
      return refinementDepth_;
   }

 public:
   RefineEverywhere(uint_t depth) : refinementDepth_(depth) {}

   void operator() ( SetupBlockForest & forest ){
      std::vector< SetupBlock* > blocks;
      forest.getBlocks(blocks);

      for (auto b : blocks)
      {
         AABB aabb = b->getAABB();
         if (b->getLevel() < targetLevel(aabb)) b->setMarker(true);
      }
   }
};

/////////////////////////////////////////
///         Lid Driven Cavity         ///
/////////////////////////////////////////

class LDCRefinement {
 private:
   const uint_t refinementDepth_;

 public:
   LDCRefinement(const uint_t depth) : refinementDepth_(depth) {};

   void operator() ( SetupBlockForest & forest ){
      std::vector< SetupBlock* > blocks;
      forest.getBlocks(blocks);

      for (auto b : blocks)
      {
         if(forest.atDomainZMaxBorder(*b) && ( forest.atDomainXMinBorder(*b) || forest.atDomainXMaxBorder(*b) )){
            if(b->getLevel() < refinementDepth_){
               b->setMarker(true);
            }
         }
      }
   }
};

class LDC : public SetupBase {
 private:
   const real_t lidVelocity_LU_;

   const real_t omegaShear_;

   const std::string refinementProfile_;
   const uint_t refinementDepth_;

   const FlagUID noSlipFlagUID_;
   const FlagUID ubbFlagUID_;

   std::shared_ptr< lbm::GeneratedNoSlip > noSlipObject;
   std::shared_ptr< lbm::GeneratedUBB > ubbObject;

 public:
   LDC(const real_t lidVelocity_LU, const real_t omegaShear,
       const std::string profile, const uint_t depth)
      : lidVelocity_LU_(lidVelocity_LU), omegaShear_(omegaShear),
        refinementProfile_(profile), refinementDepth_(depth),
        noSlipFlagUID_("NoSlip"), ubbFlagUID_("Lid")
   {};

   real_t shearRelaxationRate() const override {
      return omegaShear_;
   }

   Vector3< real_t > acceleration() const override { return Vector3< real_t >(0.0); }

   RefinementSelectionFunctor refinementSelector() override {
      if(refinementProfile_ == "refineEverywhere"){
         return RefineEverywhere(refinementDepth_);
      }
      if(refinementProfile_ == "LDC"){
         return LDCRefinement(refinementDepth_);
      }
      WALBERLA_ABORT("Invalid refinement profile.");
   }

   void setupBoundaryFlagField(StructuredBlockForest & sbfs, const BlockDataID flagFieldID) override {
      for(auto bIt = sbfs.begin(); bIt != sbfs.end(); ++bIt)
      {
         Block& b = dynamic_cast< Block& >(*bIt);
         uint_t level = b.getLevel();
         auto flagField = b.getData< FlagField_T >(flagFieldID);
         uint8_t noslipFlag = flagField->registerFlag(noSlipFlagUID_);
         uint8_t ubbFlag = flagField->registerFlag(ubbFlagUID_);
         for(auto cIt = flagField->beginWithGhostLayerXYZ(2); cIt != flagField->end(); ++cIt){
            Cell localCell = cIt.cell();
            Cell globalCell(localCell);
            sbfs.transformBlockLocalToGlobalCell(globalCell, b);
            Vector3< real_t > cellCenter(sbfs.getCellCenter(globalCell, level));
            if( globalCell.z() >= cell_idx_c(sbfs.getNumberOfZCells(level) )){
               flagField->addFlag(localCell, ubbFlag);
            }
            else if( globalCell.z() < 0
                   || globalCell.x() < 0
                   || globalCell.x() >= cell_idx_c(sbfs.getNumberOfXCells(level))){
                  flagField->addFlag(localCell, noslipFlag);
            }
         }
      }
   }

   void getBoundarySweeps(std::shared_ptr< StructuredBlockForest > & sbfs,
                          const BlockDataID pdfFieldID,
                          const BlockDataID flagFieldID,
                          const FlagUID fluidFlagUID,
                          std::vector< std::function< void( Block *) > > & sweeps) override {
      noSlipObject = std::make_shared< lbm::GeneratedNoSlip >(sbfs, pdfFieldID);
      noSlipObject->fillFromFlagField< FlagField_T >(sbfs, flagFieldID, noSlipFlagUID_, fluidFlagUID);
      BoundarySweep< LatticeModel_T, lbm::GeneratedNoSlip, LatticeModel_T::inplace > noSlipSweep(*noSlipObject, pdfFieldID);
      sweeps.emplace_back(noSlipSweep);

      VelocityCallback vCallback = [&](const Cell &, const shared_ptr< StructuredBlockForest >&, IBlock&){
        return Vector3< real_t > (lidVelocity_LU_, real_c(0), real_c(0));
      };

      ubbObject = std::make_shared< lbm::GeneratedUBB >(sbfs, pdfFieldID, vCallback);
      ubbObject->fillFromFlagField< FlagField_T >(sbfs, flagFieldID, ubbFlagUID_, fluidFlagUID);
      BoundarySweep< LatticeModel_T, lbm::GeneratedUBB, LatticeModel_T::inplace > ubbSweep(*ubbObject, pdfFieldID);
      sweeps.emplace_back(ubbSweep);
   }

};

Frederik Hennig's avatar
Frederik Hennig committed
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
/////////////////////////////////////////
///         Flow Around Cylinder      ///
/////////////////////////////////////////

class FlowAroundCylinder : public SetupBase
{
   friend class ConcentricCylinderWakeRefinement;

 private:
   Vector3< real_t > cylCenter_;
   const real_t cylRadius_;

   const real_t peakInflowVelocity_;
   const bool inflowParabolicProfile_;
   const bool walls_;

   const real_t omegaShear_;

   const std::string refinementProfile_;
   const uint_t refinementDepth_;

   Vector3< real_t > refinementRadii_;
   Vector3< real_t > refinementTailLengths_;

   const FlagUID noSlipFlagUID_;
   const FlagUID inflowFlagUID_;
   const FlagUID outflowFlagUID_;

   std::shared_ptr< lbm::GeneratedNoSlip > noSlipObject;
   std::shared_ptr< lbm::GeneratedUBB > inflowObject;
   std::shared_ptr< lbm::GeneratedOutflow > outflowObject;

 public:
   FlowAroundCylinder(const Vector3< real_t > cylCenter, const real_t cylRadius,
                      const real_t peakInflowVelocity, const bool inflowParabolicProfile, const bool walls,
                      const real_t omegaShear,
                      const std::string refinementProfile, const uint_t refinementDepth,
                      const Vector3< real_t > refinementRadii, const Vector3< real_t > refinementTailLengths)
   : cylRadius_(cylRadius),
     peakInflowVelocity_(peakInflowVelocity), inflowParabolicProfile_(inflowParabolicProfile), walls_(walls),
     omegaShear_(omegaShear),
     refinementProfile_(refinementProfile), refinementDepth_(refinementDepth),
     noSlipFlagUID_("NoSlip"), inflowFlagUID_("Inflow"), outflowFlagUID_("Outflow")
   {
      cylCenter_ = cylCenter;
      refinementRadii_ = refinementRadii;
      refinementTailLengths_ = refinementTailLengths;
   }

   real_t shearRelaxationRate() const override {
      return omegaShear_;
   }

   Vector3< real_t > acceleration() const override { return Vector3< real_t >(0.0); }

   RefinementSelectionFunctor refinementSelector() override;

   void setupBoundaryFlagField(StructuredBlockForest & sbfs, const BlockDataID flagFieldID) override {
      for(auto bIt = sbfs.begin(); bIt != sbfs.end(); ++bIt)
      {
         Block& b = dynamic_cast< Block& >(*bIt);
         uint_t level = b.getLevel();

         auto flagField = b.getData< FlagField_T >(flagFieldID);
         uint8_t noslipFlag = flagField->registerFlag(noSlipFlagUID_);
         uint8_t inflowFlag = flagField->registerFlag(inflowFlagUID_);
         uint8_t outflowFlag = flagField->registerFlag(outflowFlagUID_);

         for(auto cIt = flagField->beginWithGhostLayerXYZ(2); cIt != flagField->end(); ++cIt){
            Cell localCell = cIt.cell();
            Cell globalCell(localCell);
            sbfs.transformBlockLocalToGlobalCell(globalCell, b);
            Vector3< real_t > cellCenter(sbfs.getCellCenter(globalCell, level));

            if(walls_ && ( globalCell.y() < 0 || globalCell.y() >= cell_idx_c(sbfs.getNumberOfYCells(level)) )){
               flagField->addFlag(localCell, noslipFlag);
            } else {
               // No wall -> Inflow or Outflow
               if(globalCell.x() < 0){
                  flagField->addFlag(localCell, inflowFlag);
               }
               else if(globalCell.x() >= cell_idx_c(sbfs.getNumberOfXCells(level))){
386
                  flagField->addFlag(localCell, outflowFlag);
Frederik Hennig's avatar
Frederik Hennig committed
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
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
444
445
446
447
448
449
               }
               else{
                  Vector3< real_t > distanceInPlane(cellCenter - cylCenter_);
                  distanceInPlane[2] = real_t(0.0);
                  if(distanceInPlane.length() <= cylRadius_){
                     flagField->addFlag(localCell, noslipFlag);
                  }
               }
            }
         }
      }
   }

   void getBoundarySweeps(std::shared_ptr< StructuredBlockForest > & sbfs,
                          const BlockDataID pdfFieldID,
                          const BlockDataID flagFieldID,
                          const FlagUID fluidFlagUID,
                          std::vector< std::function< void( Block *) > > & sweeps) override {
      // NoSlip Walls
      noSlipObject = std::make_shared< lbm::GeneratedNoSlip >(sbfs, pdfFieldID);
      noSlipObject->fillFromFlagField< FlagField_T >(sbfs, flagFieldID, noSlipFlagUID_, fluidFlagUID);
      BoundarySweep< LatticeModel_T, lbm::GeneratedNoSlip, LatticeModel_T::inplace > noSlipSweep(*noSlipObject, pdfFieldID);
      sweeps.emplace_back(noSlipSweep);

      real_t ySize = sbfs->getDomain().ySize();
      real_t R = ySize / real_t(2.0);

      // Inflow Velocity Profile

      VelocityCallback vCallback = [this, R](const Cell & cell, const shared_ptr< StructuredBlockForest >& blocks, IBlock& ib){
         if(this->walls_){
            Block & b = dynamic_cast< Block & >(ib);
            Cell globalCell(cell);
            blocks->transformBlockLocalToGlobalCell(globalCell, b);
            Vector3< real_t > cellCenter(blocks->getCellCenter(globalCell, b.getLevel()));
            real_t r = abs(cellCenter[1] - R);
            real_t velocityFactor = (R*R - r*r) / (R*R);
            return Vector3< real_t > (peakInflowVelocity_ * velocityFactor, real_c(0.0), real_c(0.0));
         } else {
            return Vector3< real_t > (this->peakInflowVelocity_, real_c(0.0), real_c(0.0));
         }
      };

      inflowObject = std::make_shared< lbm::GeneratedUBB >(sbfs, pdfFieldID, vCallback);
      inflowObject->fillFromFlagField< FlagField_T >(sbfs, flagFieldID, inflowFlagUID_, fluidFlagUID);
      BoundarySweep< LatticeModel_T, lbm::GeneratedUBB, LatticeModel_T::inplace > inflowSweep(*inflowObject, pdfFieldID);
      sweeps.emplace_back(inflowSweep);

      // Outflow
      outflowObject = std::make_shared< lbm::GeneratedOutflow >(sbfs, pdfFieldID);
      outflowObject->fillFromFlagField< FlagField_T >(sbfs, flagFieldID, outflowFlagUID_, fluidFlagUID);
      BoundarySweep< LatticeModel_T, lbm::GeneratedOutflow, LatticeModel_T::inplace > outflowSweep(*outflowObject, pdfFieldID);
      sweeps.emplace_back(outflowSweep);
   }
};


class ConcentricCylinderWakeRefinement {
 private:
   const FlowAroundCylinder & fs_;

   uint_t targetLevel(const SetupBlock * b){
      AABB blockAabb(b->getAABB());
450
451
452
453
454
455
456
457
458
459
460
461
      Vector3< real_t > aabbCenter(blockAabb.center());
      Vector3 axisReferencePoint(std::max(aabbCenter[0], fs_.cylCenter_[0]),
                                 fs_.cylCenter_[1],
                                 aabbCenter[2]);

      for(int k = int(fs_.refinementDepth_) - 1; k >= 0 ; --k)
      {
         real_t tailLength = fs_.refinementTailLengths_[uint_c(k)];
         real_t refRadius  = fs_.refinementRadii_[uint_c(k)];

         if(blockAabb.xMax() - fs_.cylCenter_[0] <= tailLength){
            if(blockAabb.distance(axisReferencePoint) < refRadius){
Frederik Hennig's avatar
Frederik Hennig committed
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
               return uint_c(k + 1);
            }
         }
      }

      return 0;
   }

 public:
   ConcentricCylinderWakeRefinement(const FlowAroundCylinder & flowSetup) : fs_(flowSetup) {};

   void operator() ( SetupBlockForest & forest ){
      std::vector< SetupBlock* > blocks;
      forest.getBlocks(blocks);

      for (const auto & b : blocks)
      {
         if(b->getLevel() < targetLevel(b)){
            b->setMarker(true);
         }
      }
   }
};

RefinementSelectionFunctor FlowAroundCylinder::refinementSelector() {
   if(refinementProfile_ == "refineEverywhere"){
      return RefineEverywhere(refinementDepth_);
   }
   if(refinementProfile_ == "ConcentricCylinderWake"){
      return ConcentricCylinderWakeRefinement(*this);
   }
   WALBERLA_ABORT("Invalid refinement profile.");
}

496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
/////////////////////////////////////////
///        Flow Setup Selection       ///
/////////////////////////////////////////

std::shared_ptr< SetupBase > getSetup(Config::BlockHandle & cblock){
   std::string setupName = cblock.getParameter<std::string>("setup");

   if(setupName ==  "LDC"){
      return std::make_shared< LDC >(
         cblock.getParameter< real_t >("lidVelocity_LU"),
         cblock.getParameter< real_t >("omegaShear"),
         cblock.getParameter< std::string >("refinementProfile"),
         cblock.getParameter< uint_t >("refinementDepth")
      );
   }

Frederik Hennig's avatar
Frederik Hennig committed
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
   if(setupName == "FlowAroundCylinder"){
      return std::make_shared< FlowAroundCylinder >(
         cblock.getParameter< Vector3< real_t > >("cylCenter"),
         cblock.getParameter< real_t >("cylRadius"),
         cblock.getParameter< real_t >("peakInflowVelocity"),
         cblock.getParameter< bool >("inflowParabolicProfile"),
         cblock.getParameter< bool >("walls"),
         cblock.getParameter< real_t >("omegaShear"),
         cblock.getParameter< std::string >("refinementProfile"),
         cblock.getParameter< uint_t >("refinementDepth"),
         cblock.getParameter< Vector3< real_t > >("refinementRadii", Vector3< real_t >(real_c(0.0))),
         cblock.getParameter< Vector3< real_t > >("refinementTailLengths", Vector3< real_t >(real_c(0.0)))
      );
   }

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
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
   WALBERLA_ABORT("No valid setup was specified.");
}

inline void setupUniformVelocityField(const std::shared_ptr< StructuredBlockForest > & sbfs,
                                   const BlockDataID pdfFieldId,
                                   const BlockDataID densityFieldId,
                                   const BlockDataID velocityFieldId,
                                   const Vector3< real_t > velocity){
   lbm::GeneratedMacValueKernels< LatticeModel_T >::Setter init(pdfFieldId, densityFieldId, velocityFieldId);
   for(auto bIt = sbfs->begin(); bIt != sbfs->end(); ++bIt){
      Block & b = dynamic_cast< Block & > (*bIt);
      auto densityField = b.getData< ScalarField_T >(densityFieldId);
      auto velocityField = b.getData< VectorField_T >(velocityFieldId);

      densityField->set(1.0);

      for(auto cIt = velocityField->beginWithGhostLayerXYZ(); cIt != velocityField->end(); ++cIt){
         cIt.getF(0) = velocity[0];
         cIt.getF(1) = velocity[1];
         cIt.getF(2) = velocity[2];
      }

      auto pdfField = b.getData< PdfField_T >(pdfFieldId);
      pdfField->setWithGhostLayer(0.0);
      init(bIt.get());
   }
}


/////////////////////////////////////////
///          Setup Block Forest       ///
/////////////////////////////////////////

static void createSetupBlockForest(SetupBlockForest & setupBfs,
                                   const Config::BlockHandle & domainSetup,
                                   SetupBase & setup)
{
   Vector3< real_t > domainSize = domainSetup.getParameter< Vector3< real_t > >("domainSize");
   Vector3< uint_t > rootBlocks = domainSetup.getParameter< Vector3< uint_t > >("rootBlocks");
   Vector3< bool > periodic     = domainSetup.getParameter< Vector3< bool > >("periodic");

   auto refSelection = setup.refinementSelector();
   setupBfs.addRefinementSelectionFunction(std::function< void(SetupBlockForest &) > (refSelection) );
   AABB domain(real_t(0.0), real_t(0.0), real_t(0.0), domainSize[0], domainSize[1], domainSize[2]);
   setupBfs.init(domain, rootBlocks[0], rootBlocks[1], rootBlocks[2], periodic[0], periodic[1], periodic[2]);
   setupBfs.balanceLoad(blockforest::StaticLevelwiseCurveBalance(true), uint_c(MPIManager::instance()->numProcesses()));
}


void run(std::shared_ptr< Config > & config) {
   auto params = config->getBlock("Parameters");
   if(!params) WALBERLA_ABORT("Simulation Parameters missing!");

   const std::string vtkBaseFolder = params.getParameter< std::string >("vtkBaseFolder", "vtk_out");

   const uint_t runs        = params.getParameter< uint_t >("runs", uint_c(500));
   const uint_t stepsPerRun = params.getParameter< uint_t >("stepsPerRun", uint_c(500));

   const uint_t vtkWriteFrequency = params.getParameter< uint_t >("vtkWriteFrequency", uint_c(0));
   const double remainingTimeLoggerFrequency =
      params.getParameter< double >("remainingTimeLoggerFrequency", 3.0); // in seconds

   // Setup Procedure -- make sure the desired geometry is correctly constructed
   auto flowSetup = getSetup(params);

   SetupBlockForest setupBfs;
   WALBERLA_LOG_INFO_ON_ROOT("Generating SetupBlockForest...");
   createSetupBlockForest(setupBfs, params, *flowSetup);

   // Create structured block forest
   Vector3< uint_t > cellsPerBlock = params.getParameter< Vector3< uint_t > >("cellsPerBlock");

   WALBERLA_LOG_INFO_ON_ROOT("Creating structured block forest...");
   auto bfs  = std::make_shared< BlockForest >(uint_c(MPIManager::instance()->worldRank()), setupBfs);
   auto sbfs = std::make_shared< StructuredBlockForest >(bfs, cellsPerBlock[0], cellsPerBlock[1], cellsPerBlock[2]);
   sbfs->createCellBoundingBoxes();

   vtk::writeDomainDecomposition(sbfs, "domainDecomposition", vtkBaseFolder);

   // Setup Boundary Handling
   WALBERLA_LOG_INFO_ON_ROOT("Setting up boundary conditions...");
   BlockDataID flagFieldID = field::addFlagFieldToStorage< FlagField_T >(sbfs, "Boundary Flag Field", uint_c(3));

   const FlagUID fluidFlagUID("Fluid");

   flowSetup->setupBoundaryFlagField(*sbfs, flagFieldID);
   geometry::setNonBoundaryCellsToDomain< FlagField_T >(*sbfs, flagFieldID, fluidFlagUID, 2);

   auto flagFieldWriter = field::createVTKOutput< FlagField_T, uint8_t >(flagFieldID, *sbfs, "boundaryFlagFieldOutput", 1, 0, false, vtkBaseFolder);
   flagFieldWriter();

   // Maybe exit early
   if(stepsPerRun == 0) {
      WALBERLA_LOG_INFO_ON_ROOT("Domain setup complete and written to VTK. Terminating, since stepsPerRun == 0.");
      return;
   }

   // Density and Velocity Fields
   WALBERLA_LOG_INFO_ON_ROOT("Creating velocity field...");

   BlockDataID densityFieldID = field::addToStorage< ScalarField_T >(sbfs, "rho", real_t(1), field::fzyx);
   BlockDataID velocityFieldID = field::addToStorage< VectorField_T >(sbfs, "u", real_t(0), field::fzyx);

   // Lattice Model and PDF field
   WALBERLA_LOG_INFO_ON_ROOT("Creating LB data structures...");

   Vector3< real_t > accelVec = flowSetup->acceleration();
   real_t omegaShear = flowSetup->shearRelaxationRate();

   LatticeModel_T latticeModel = LatticeModel_T(accelVec[0], accelVec[1], accelVec[2], omegaShear);
   BlockDataID pdfFieldID      = lbm::addPdfFieldToStorage< LatticeModel_T >(sbfs, "pdfs", latticeModel, 2, field::fzyx);

   // Boundary Handling
   std::vector< std::function< void( Block *) > > boundarySweeps;
   flowSetup->getBoundarySweeps(sbfs, pdfFieldID, flagFieldID, fluidFlagUID, boundarySweeps);
   std::function< void ( Block * ) > boundaryFunctor = [&](Block * b) {
      for(auto & func : boundarySweeps){
         func(b);
      }
   };

   // Setup Communication
   WALBERLA_LOG_INFO_ON_ROOT("Setting up communication...");
   auto comm = std::make_shared< NonUniformBufferedScheme< CommunicationStencil_T > >(sbfs);
   auto packInfo = lbm::setupNonuniformPdfCommunication< LatticeModel_T, lbm::GeneratedPackingKernels >(sbfs, pdfFieldID);
   comm->addPackInfo(packInfo);

   // Velocity Field Setup
   WALBERLA_LOG_INFO_ON_ROOT("Setting up velocity field...");
   setupUniformVelocityField(sbfs, pdfFieldID, densityFieldID, velocityFieldID, Vector3< real_t > (0.0));

   // Recursive Time Step
   lbm::GeneratedLatticeModel::Sweep sweep(pdfFieldID);
   lbm::BasicRecursiveTimeStep< LatticeModel_T, lbm::GeneratedPackingKernels, LatticeModel_T::Sweep > timestep(
      sbfs, pdfFieldID, sweep, boundaryFunctor, comm, packInfo);

   // Timeloop
   WALBERLA_LOG_INFO_ON_ROOT("Setting up Time Loop...");
   SweepTimeloop timeloop(sbfs, stepsPerRun);

   // VTK Output
   if( vtkWriteFrequency > 0){
      lbm::GeneratedMacValueKernels< LatticeModel_T >::Getter macValueGetter(pdfFieldID, densityFieldID, velocityFieldID);

      std::function< void() > getterSweep = [&](){
        for(auto bIt = sbfs->begin(); bIt != sbfs->end(); ++bIt){
           macValueGetter(bIt.get());
        }
      };

      field::FlagFieldCellFilter< FlagField_T > fluidFilter(flagFieldID, fluidFlagUID);
      auto densityWriter = std::make_shared< field::VTKWriter< ScalarField_T, float > >(densityFieldID, "density");
      auto velocityWriter = std::make_shared< field::VTKWriter< VectorField_T , float > >(velocityFieldID, "velocity");

      auto flowFieldOutput = vtk::createVTKOutput_BlockData(sbfs, "FlowFieldOutput", vtkWriteFrequency, 0, false, vtkBaseFolder);

      flowFieldOutput->addCellInclusionFilter(fluidFilter);
      flowFieldOutput->addCellDataWriter(densityWriter);
      flowFieldOutput->addCellDataWriter(velocityWriter);
      flowFieldOutput->addBeforeFunction(getterSweep);

      auto flowFieldWriter = vtk::writeFiles(flowFieldOutput);
      timeloop.addFuncBeforeTimeStep(flowFieldWriter);
   }

   timeloop.addFuncAfterTimeStep(timestep);

   // Aaaand... Go!

   uint_t totalCellsPerBlock(1);
   for(uint_t i = 0; i < 3; ++i){
      totalCellsPerBlock *= sbfs->getNumberOfCellsPerBlock(i);
   }

701
702
703
   const uint_t numProcesses = uint_c( MPIManager::instance()->numProcesses() );
   const uint_t ompNumThreads      = uint_c(omp_get_max_threads());

704
   const uint_t processLocalNumberOfBlocks = sbfs->getNumberOfBlocks();
705
706
   const uint_t globalNumberOfBlocks = setupBfs.getNumberOfBlocks();
   const real_t avgNumberOfBlocksPerProcess = real_c(globalNumberOfBlocks) / real_c(numProcesses);
707
708
709

   const uint_t processLocalNumberOfCells = processLocalNumberOfBlocks * totalCellsPerBlock;
   const uint_t globalNumberOfCells = globalNumberOfBlocks * totalCellsPerBlock;
710
   const real_t avgNumberOfCellsPerProcess = avgNumberOfBlocksPerProcess * real_c(totalCellsPerBlock);
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756

#define MLUPS (real_c(1e-6))

   for(uint_t run = 0; run < runs; ++run){
      timeloop.setCurrentTimeStepToZero();

      WALBERLA_MPI_WORLD_BARRIER();
      WcTimer simTimer;
      simTimer.start();
      timeloop.run();
      simTimer.end();

      double localTime = simTimer.last();
      double globalMaxTime = mpi::reduce(localTime, mpi::MAX);

      double localMlups = ( real_c( processLocalNumberOfCells * stepsPerRun ) / localTime ) * MLUPS;

      double localMlupsMin = mpi::reduce(localMlups, mpi::MIN);
      double localMlupsMax = mpi::reduce(localMlups, mpi::MAX);
      double localMlupsAvg = mpi::reduce(localMlups, mpi::SUM) / double( numProcesses );

      WALBERLA_ROOT_SECTION()
      {
         // Total MLUPs
         double globalMlups = ( real_c( globalNumberOfCells * stepsPerRun ) / globalMaxTime) * MLUPS;

         // Callback
         python_coupling::PythonCallback resultsCallback("results_callback");
         if( !resultsCallback.isCallable() ) WALBERLA_ABORT("Results callback not available.");

         /**
          * Need:
          *  - Stencil
          *  - Streaming Pattern
          *  - Collision Operator
          *  - Per-Process MLUPs MIN/MAX/AVG
          *  - Total MLUPs ( by localTime of slowest process )
          */

         resultsCallback.data().exposeValue("mpiNumProcesses", numProcesses);
         resultsCallback.data().exposeValue("ompNumThreads", ompNumThreads);

         resultsCallback.data().exposeValue("stencil", stencilName);
         resultsCallback.data().exposeValue("streamingPattern", streamingPattern);
         resultsCallback.data().exposeValue("collisionSetup", collisionSetup);

757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
         resultsCallback.data().exposeValue("globalNumberOfBlocks", globalNumberOfBlocks);
         resultsCallback.data().exposeValue("avgNumberOfBlocksPerProcess", avgNumberOfBlocksPerProcess);
         resultsCallback.data().exposeValue("globalNumberOfCells", globalNumberOfCells);
         resultsCallback.data().exposeValue("avgNumberOfCellsPerProcess", avgNumberOfCellsPerProcess);

         for(uint_t l = 0; l < setupBfs.getNumberOfLevels(); ++l){
            uint_t numBlocks = setupBfs.getNumberOfBlocks(l);
            std::string keyTotal("globalNumberOfLevel" + std::to_string(l) + "Blocks");
            resultsCallback.data().exposeValue(&keyTotal[0], numBlocks);

            real_t avgNumBlocks = real_c(numBlocks) / real_c(numProcesses);
            std::string keyAvg("avgNumberOfLevel" + std::to_string(l) + "BlocksPerProcess");
            resultsCallback.data().exposeValue(&keyAvg[0], avgNumBlocks);
         }

772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
         resultsCallback.data().exposeValue("globalMLUPs", globalMlups);
         resultsCallback.data().exposeValue("localMLUPsMin", localMlupsMin);
         resultsCallback.data().exposeValue("localMLUPsMax", localMlupsMax);
         resultsCallback.data().exposeValue("localMLUPsAvg", localMlupsAvg);

         resultsCallback();
      }
   }

}

} // namespace walberla

int main(int argc, char ** argv){
   walberla::Environment env(argc, argv);
   walberla::MPIManager::instance()->useWorldComm();

   for( auto cfg = walberla::python_coupling::configBegin(argc, argv);
        cfg != walberla::python_coupling::configEnd();
        ++cfg )
   {
      auto config = *cfg;
      walberla::run(config);
   }

   return EXIT_SUCCESS;
}