CouetteFlow.cpp 59.7 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
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
386
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
450
451
452
453
454
455
456
457
458
459
460
461
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
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
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
701
702
703
704
705
706
707
708
709
710
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
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
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
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
//======================================================================================================================
//
//  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 CouetteFlow.cpp
//! \author Florian Schornbaum <florian.schornbaum@fau.de>
//
//======================================================================================================================

#include "blockforest/AABBRefinementSelection.h"
#include "blockforest/communication/NonUniformBufferedScheme.h"
#include "blockforest/SetupBlockForest.h"
#include "blockforest/StructuredBlockForest.h"
#include "blockforest/loadbalancing/StaticCurve.h"

#include "boundary/BoundaryHandling.h"

#include "core/Abort.h"
#include "core/DataTypes.h"
#include "core/SharedFunctor.h"
#include "core/cell/CellInterval.h"
#include "core/config/Config.h"
#include "core/debug/Debug.h"
#include "core/debug/TestSubsystem.h"
#include "core/logging/Logging.h"
#include "core/math/Sample.h"
#include "core/math/Vector3.h"
#include "core/mpi/Environment.h"
#include "core/mpi/Gatherv.h"
#include "core/mpi/MPIManager.h"
#include "core/mpi/RecvBuffer.h"
#include "core/mpi/Reduce.h"
#include "core/mpi/SendBuffer.h"
#include "core/timing/RemainingTimeLogger.h"
#include "core/timing/TimingPool.h"

#include "field/AccuracyEvaluation.h"
#include "field/AccuracyEvaluationLinePlot.h"
#include "field/AddToStorage.h"
#include "field/CellCounter.h"
#include "field/FlagField.h"
#include "field/FlagUID.h"
#include "field/StabilityChecker.h"
#include "field/VolumetricFlowRateEvaluation.h"
#include "field/adaptors/AdaptorCreators.h"
#include "field/iterators/FieldIterator.h"
#include "field/vtk/FlagFieldCellFilter.h"
#include "field/vtk/VTKWriter.h"

#include "lbm/BlockForestEvaluation.h"
#include "lbm/MassEvaluation.h"
#include "lbm/PerformanceEvaluation.h"
#include "lbm/boundary/NoSlip.h"
#include "lbm/boundary/SimpleUBB.h"
#include "lbm/field/Adaptors.h"
#include "lbm/field/AddToStorage.h"
#include "lbm/field/PdfField.h"
#include "lbm/lattice_model/CollisionModel.h"
#include "lbm/lattice_model/D3Q15.h"
#include "lbm/lattice_model/D3Q19.h"
#include "lbm/lattice_model/D3Q27.h"
#include "lbm/refinement/PdfFieldSyncPackInfo.h"
#include "lbm/refinement/TimeStep.h"
#include "lbm/sweeps/CellwiseSweep.h"
#include "lbm/sweeps/SplitPureSweep.h"
#include "lbm/sweeps/SplitSweep.h"
#include "lbm/vtk/Density.h"
#include "lbm/vtk/NonEquilibrium.h"
#include "lbm/vtk/Velocity.h"

#include "postprocessing/sqlite/SQLite.h"

#include "stencil/D3Q15.h"
#include "stencil/D3Q19.h"
#include "stencil/D3Q27.h"

#include "timeloop/SweepTimeloop.h"

#include "vtk/BlockCellDataWriter.h"
#include "vtk/Initialization.h"
#include "vtk/VTKOutput.h"

#include <boost/bind.hpp>
#include <boost/mpl/or.hpp>
#include <boost/type_traits/is_same.hpp>

#include <algorithm>
#include <cmath>
#include <cstdlib>
#include <cstring>
#include <iostream>
#include <utility>
#include <vector>



namespace couette_flow {

///////////
// USING //
///////////

using namespace walberla;
using walberla::uint_t;
using walberla::real_t;

//////////////
// TYPEDEFS //
//////////////

typedef lbm::D3Q15< lbm::collision_model::SRT,      false > D3Q15_SRT_INCOMP;
typedef lbm::D3Q15< lbm::collision_model::SRT,      true  > D3Q15_SRT_COMP;
typedef lbm::D3Q15< lbm::collision_model::TRT,      false > D3Q15_TRT_INCOMP;
typedef lbm::D3Q15< lbm::collision_model::TRT,      true  > D3Q15_TRT_COMP;

typedef lbm::D3Q19< lbm::collision_model::SRT,      false > D3Q19_SRT_INCOMP;
typedef lbm::D3Q19< lbm::collision_model::SRT,      true  > D3Q19_SRT_COMP;
typedef lbm::D3Q19< lbm::collision_model::TRT,      false > D3Q19_TRT_INCOMP;
typedef lbm::D3Q19< lbm::collision_model::TRT,      true  > D3Q19_TRT_COMP;
typedef lbm::D3Q19< lbm::collision_model::D3Q19MRT, false > D3Q19_MRT_INCOMP;

typedef lbm::D3Q27< lbm::collision_model::SRT,      false > D3Q27_SRT_INCOMP;
typedef lbm::D3Q27< lbm::collision_model::SRT,      true  > D3Q27_SRT_COMP;
typedef lbm::D3Q27< lbm::collision_model::TRT,      false > D3Q27_TRT_INCOMP;
typedef lbm::D3Q27< lbm::collision_model::TRT,      true  > D3Q27_TRT_COMP;

template< typename LatticeModel_T >
struct Types
{
   typedef typename LatticeModel_T::Stencil Stencil_T;
   typedef lbm::PdfField< LatticeModel_T >  PdfField_T;
};

typedef walberla::uint16_t  flag_t;
typedef FlagField< flag_t > FlagField_T;

const uint_t FieldGhostLayers  = uint_t(4);

///////////
// FLAGS //
///////////

const FlagUID  Fluid_Flag( "fluid" );
const FlagUID    UBB_Flag( "velocity bounce back" );
const FlagUID NoSlip_Flag( "no slip" );

/////////////////////
// OUTPUT HELPERS  //
/////////////////////

template< typename LatticeModel_T, class Enable = void >
struct StencilString;

template< typename LatticeModel_T >
struct StencilString< LatticeModel_T, typename boost::enable_if_c< boost::is_same< typename LatticeModel_T::Stencil, stencil::D3Q15 >::value >::type >
{
   static const char * str() { return "D3Q15"; }
};

template< typename LatticeModel_T >
struct StencilString< LatticeModel_T, typename boost::enable_if_c< boost::is_same< typename LatticeModel_T::Stencil, stencil::D3Q19 >::value >::type >
{
   static const char * str() { return "D3Q19"; }
};

template< typename LatticeModel_T >
struct StencilString< LatticeModel_T, typename boost::enable_if_c< boost::is_same< typename LatticeModel_T::Stencil, stencil::D3Q27 >::value >::type >
{
   static const char * str() { return "D3Q27"; }
};


template< typename LatticeModel_T, class Enable = void >
struct CollisionModelString;

template< typename LatticeModel_T >
struct CollisionModelString< LatticeModel_T, typename boost::enable_if_c< boost::is_same< typename LatticeModel_T::CollisionModel::tag,
                                                                                          lbm::collision_model::SRT_tag >::value >::type >
{
   static const char * str() { return "SRT"; }
};

template< typename LatticeModel_T >
struct CollisionModelString< LatticeModel_T, typename boost::enable_if_c< boost::is_same< typename LatticeModel_T::CollisionModel::tag,
                                                                                          lbm::collision_model::TRT_tag >::value >::type >
{
   static const char * str() { return "TRT"; }
};

template< typename LatticeModel_T >
struct CollisionModelString< LatticeModel_T, typename boost::enable_if_c< boost::is_same< typename LatticeModel_T::CollisionModel::tag,
                                                                                          lbm::collision_model::MRT_tag >::value >::type >
{
   static const char * str() { return "MRT"; }
};

//////////////////////
// Parameter Struct //
//////////////////////

struct Setup
{
   uint_t xBlocks;
   uint_t yBlocks;
   uint_t zBlocks;

   uint_t xCells;
   uint_t yCells;
   uint_t zCells;

   real_t Re;

   real_t viscosity_L; // on the coarsest grid
   real_t maxVelocity_L;
   real_t meanVelocity_L;
   real_t flowRate_L; // volumetric flow rate - in lattice units of the coarsest grid
};






/////////////////////
// BLOCK STRUCTURE //
/////////////////////

class BorderRefinementSelection
{
public:

   BorderRefinementSelection( const Setup & setup, const uint_t level, const real_t bufferDistance ) :
      setup_( setup ), level_( level ), bufferDistance_( bufferDistance ) {}

   void operator()( SetupBlockForest & forest )
   {
      for( auto block = forest.begin(); block != forest.end(); ++block )
      {
         if( block->getLevel() < level_ && !channelContains( forest, *block ) )
            block->setMarker( true );
      }
   }

private:

   bool channelContains( SetupBlockForest & forest, const SetupBlock & block )
   {
      AABB domain = forest.getDomain();
      domain.setAxisBounds( uint_t(1), domain.yMin() + bufferDistance_, domain.yMax() - bufferDistance_ );

      return domain.contains( block.getAABB() ) && !(forest.atDomainYMinBorder( block )) && !(forest.atDomainYMaxBorder( block ));
   }

   const Setup & setup_;
   uint_t level_;
   real_t bufferDistance_;

}; // class BorderRefinementSelection



static void workloadAndMemoryAssignment( SetupBlockForest& forest, const memory_t memoryPerBlock )
{
   for( auto block = forest.begin(); block != forest.end(); ++block )
   {
      block->setWorkload( numeric_cast< workload_t >( uint_t(1) << block->getLevel() ) );
      block->setMemory( memoryPerBlock );
   }
}



static shared_ptr< SetupBlockForest > createSetupBlockForest( const blockforest::RefinementSelectionFunctions & refinementSelectionFunctions,
                                                              const Setup & setup,
                                                              uint_t numberOfProcesses, const uint_t blocksPerProcess,
                                                              const memory_t memoryPerCell, const memory_t processMemoryLimit,
                                                              const bool outputSetupForest )
{
   shared_ptr< SetupBlockForest > forest = make_shared< SetupBlockForest >();

   const memory_t memoryPerBlock = numeric_cast< memory_t >( ( setup.xCells + uint_t(2) * FieldGhostLayers ) *
                                                             ( setup.yCells + uint_t(2) * FieldGhostLayers ) *
                                                             ( setup.zCells + uint_t(2) * FieldGhostLayers ) ) * memoryPerCell;

   forest->addRefinementSelectionFunction( refinementSelectionFunctions );
   forest->addWorkloadMemorySUIDAssignmentFunction( boost::bind( workloadAndMemoryAssignment, _1, memoryPerBlock ) );

   forest->init( AABB( real_c(0), real_c(0), real_c(0), real_c( setup.xBlocks * setup.xCells ),
                                                        real_c( setup.yBlocks * setup.yCells ),
                                                        real_c( setup.zBlocks * setup.zCells ) ),
                 setup.xBlocks, setup.yBlocks, setup.zBlocks, true, false, true );

   MPIManager::instance()->useWorldComm();

   if( blocksPerProcess != 0 )
      numberOfProcesses = uint_c( std::ceil( real_c( forest->getNumberOfBlocks() ) / real_c( blocksPerProcess ) ) );

   forest->balanceLoad( blockforest::StaticLevelwiseCurveBalance(true), numberOfProcesses, real_t(0), processMemoryLimit, true );

   if( outputSetupForest ) 
   {
      forest->writeVTKOutput( "domain_decomposition" );
      forest->writeCSV( "process_distribution" );
   }

   WALBERLA_LOG_INFO_ON_ROOT( "SetupBlockForest created successfully:\n" << *forest );

   return forest;
}



shared_ptr< blockforest::StructuredBlockForest > createStructuredBlockForest( const Config::BlockHandle & configBlock,
                                                                              const blockforest::RefinementSelectionFunctions & refinementSelectionFunctions,
                                                                              const Setup & setup,
                                                                              const memory_t memoryPerCell, const memory_t processMemoryLimit )
{
   if( configBlock.isDefined( "sbffile" ) )
   {
      std::string sbffile = configBlock.getParameter< std::string >( "sbffile" );

      WALBERLA_LOG_INFO_ON_ROOT( "Creating the block structure: loading from file \'" << sbffile << "\' ..." );

      MPIManager::instance()->useWorldComm();

      auto bf = shared_ptr< BlockForest >( new BlockForest( uint_c( MPIManager::instance()->rank() ), sbffile.c_str(), true, false ) );

      auto sbf = shared_ptr< StructuredBlockForest >( new StructuredBlockForest( bf, setup.xCells, setup.yCells, setup.zCells ) );
      sbf->createCellBoundingBoxes();

      return sbf;
   }

   WALBERLA_LOG_INFO_ON_ROOT( "Creating the block structure ..." );

   shared_ptr< SetupBlockForest > sforest = createSetupBlockForest( refinementSelectionFunctions, setup,
                                                                    uint_c( MPIManager::instance()->numProcesses() ), uint_t(0),
                                                                    memoryPerCell, processMemoryLimit,
                                                                    configBlock.getParameter< bool >( "outputSetupForest", false ) );

   auto bf = shared_ptr< blockforest::BlockForest >( new blockforest::BlockForest( uint_c( MPIManager::instance()->rank() ), *sforest, false ) );

   auto sbf = shared_ptr< blockforest::StructuredBlockForest >( new blockforest::StructuredBlockForest( bf, setup.xCells, setup.yCells, setup.zCells ) );
   sbf->createCellBoundingBoxes();
   
   return sbf;
}






///////////////////////
// BOUNDARY HANDLING //
///////////////////////

template< typename LatticeModel_T >
class MyBoundaryHandling
{
public:

   typedef lbm::NoSlip< LatticeModel_T, flag_t >    NoSlip_T;
   typedef lbm::SimpleUBB< LatticeModel_T, flag_t > UBB_T;

   typedef boost::tuples::tuple< NoSlip_T, UBB_T > BoundaryConditions_T;

   typedef BoundaryHandling< FlagField_T, typename Types<LatticeModel_T>::Stencil_T, BoundaryConditions_T > BoundaryHandling_T;



   MyBoundaryHandling( const BlockDataID & flagFieldId, const BlockDataID & pdfFieldId, const real_t topVelocity ) :
      flagFieldId_( flagFieldId ), pdfFieldId_( pdfFieldId ), topVelocity_( topVelocity ) {}

   BoundaryHandling_T * operator()( IBlock* const block ) const;

private:

   const BlockDataID flagFieldId_;
   const BlockDataID  pdfFieldId_;

   const real_t topVelocity_;

}; // class MyBoundaryHandling

template< typename LatticeModel_T >
typename MyBoundaryHandling<LatticeModel_T>::BoundaryHandling_T *
MyBoundaryHandling<LatticeModel_T>::operator()( IBlock * const block ) const
{
   typedef typename Types<LatticeModel_T>::PdfField_T PdfField_T;

   FlagField_T * flagField = block->getData< FlagField_T >( flagFieldId_ );
   PdfField_T *   pdfField = block->getData< PdfField_T > (  pdfFieldId_ );

   const flag_t fluid = flagField->registerFlag( Fluid_Flag );

   return new BoundaryHandling_T( "boundary handling", flagField, fluid,
         boost::tuples::make_tuple( NoSlip_T( "no slip", NoSlip_Flag, pdfField ),
                                       UBB_T( "velocity bounce back", UBB_Flag, pdfField, topVelocity_, real_t(0), real_t(0) ) ) );
}



template< typename LatticeModel_T, typename BoundaryHandling_T >
void setFlags( shared_ptr< StructuredBlockForest > & blocks, const BlockDataID & boundaryHandlingId )
{
   for( auto block = blocks->begin(); block != blocks->end(); ++block )
   {
      BoundaryHandling_T * boundaryHandling = block->template getData< BoundaryHandling_T >( boundaryHandlingId );
      
      const uint_t level = blocks->getLevel(*block);
      CellInterval domainBB = blocks->getDomainCellBB( level );
      blocks->transformGlobalToBlockLocalCellInterval( domainBB, *block );
         
      domainBB.xMin() -= cell_idx_c(FieldGhostLayers);
      domainBB.xMax() += cell_idx_c(FieldGhostLayers);
      domainBB.yMin() -= 1;
      domainBB.yMax() += 1;
      domainBB.zMin() -= cell_idx_c(FieldGhostLayers);
      domainBB.zMax() += cell_idx_c(FieldGhostLayers);

      // no slip BOTTOM
      CellInterval south( domainBB.xMin(), domainBB.yMin(), domainBB.zMin(), domainBB.xMax(), domainBB.yMin(), domainBB.zMax() );
      boundaryHandling->forceBoundary( NoSlip_Flag, south );

      // velocity TOP
      CellInterval north( domainBB.xMin(), domainBB.yMax(), domainBB.zMin(), domainBB.xMax(), domainBB.yMax(), domainBB.zMax() );
      boundaryHandling->forceBoundary( UBB_Flag, north );
      
      boundaryHandling->fillWithDomain( domainBB );
   }
}






/////////
// VTK //
/////////

template< typename LatticeModel_T, typename OutputType = float >
class ErrorVTKWriter : public vtk::BlockCellDataWriter< OutputType, 3 >
{
public:

   typedef lbm::PdfField< LatticeModel_T > PdfField_T;

   ErrorVTKWriter( const ConstBlockDataID & pdfFieldId, const std::string & id, const Setup & setup ) :
      vtk::BlockCellDataWriter< OutputType, 3 >( id ), bdid_( pdfFieldId ), pdf_( NULL ), setup_( setup ) {}

protected:

   void configure() { WALBERLA_ASSERT_NOT_NULLPTR( this->block_ ); pdf_ = this->block_->template getData< PdfField_T >( bdid_ ); }

   OutputType evaluate( const cell_idx_t x, const cell_idx_t y, const cell_idx_t z, const cell_idx_t f )
   {
      WALBERLA_ASSERT_NOT_NULLPTR( pdf_ );

      const auto & domain = this->blockStorage_->getDomain();
      const auto center = this->blockStorage_->getBlockLocalCellCenter( *(this->block_), Cell(x,y,z) );

      const real_t refVelocity_x = setup_.maxVelocity_L * ( center[1] - domain.yMin() ) / domain.ySize();

      const auto velocity = pdf_->getVelocity(x,y,z);

      // error is always a relative error:
      // 1st error component -> base = reference velocity
      // 2nd error component -> base = maximum velocity
      // 3rd error component -> base = mean velocity

      if( f == cell_idx_t(0) )
         return numeric_cast< OutputType >( std::fabs( ( refVelocity_x - velocity[0] ) / refVelocity_x ) );
      else if( f == cell_idx_t(1) )
         return numeric_cast< OutputType >( std::fabs( ( refVelocity_x - velocity[0] ) / setup_.maxVelocity_L ) );
      return numeric_cast< OutputType >( std::fabs( ( refVelocity_x - velocity[0] ) / setup_.meanVelocity_L ) );
   }

   const ConstBlockDataID bdid_;
   const PdfField_T * pdf_;
   
   const Setup & setup_;

}; // class ErrorVTKWriter



template< typename LatticeModel_T >
class MyVTKOutput {

public:

   MyVTKOutput( const ConstBlockDataID & pdfField, const ConstBlockDataID & flagField,
                vtk::VTKOutput::BeforeFunction pdfGhostLayerSync, const Setup & setup ) :
      setup_( setup ), pdfField_( pdfField ), flagField_( flagField ), pdfGhostLayerSync_( pdfGhostLayerSync ) {}

   void operator()( std::vector< shared_ptr<vtk::BlockCellDataWriterInterface> > & writers,
                    std::map< std::string, vtk::VTKOutput::CellFilter > &          filters,
                    std::map< std::string, vtk::VTKOutput::BeforeFunction > &      beforeFunctions );

private:

   const Setup & setup_;

   const ConstBlockDataID pdfField_;
   const ConstBlockDataID flagField_;

   vtk::VTKOutput::BeforeFunction pdfGhostLayerSync_;

}; // class MyVTKOutput

template< typename LatticeModel_T >
void MyVTKOutput<LatticeModel_T>::operator()( std::vector< shared_ptr<vtk::BlockCellDataWriterInterface> > & writers,
                                              std::map< std::string, vtk::VTKOutput::CellFilter > &          filters,
                                              std::map< std::string, vtk::VTKOutput::BeforeFunction > &      beforeFunctions )
{
   // block data writers

   writers.push_back( make_shared< lbm::VelocityVTKWriter<LatticeModel_T>            >( pdfField_, "VelocityFromPDF" ) );
   writers.push_back( make_shared< lbm::VelocityMagnitudeVTKWriter<LatticeModel_T>   >( pdfField_, "VelocityMagnitudeFromPDF" ) );
   writers.push_back( make_shared< lbm::DensityVTKWriter<LatticeModel_T>             >( pdfField_, "DensityFromPDF" ) );
   writers.push_back( make_shared< lbm::NonEqulibriumVTKWriter<LatticeModel_T>       >( pdfField_, "NonEquPart" ) );
   writers.push_back( make_shared< field::VTKWriter< lbm::PdfField<LatticeModel_T> > >( pdfField_, "PDF" ) );
   writers.push_back( make_shared< ErrorVTKWriter<LatticeModel_T>                    >( pdfField_, "Error", setup_ ) );

   writers.push_back( make_shared< field::VTKWriter< FlagField_T > >( flagField_, "FlagField" ) );

   // cell filters

   field::FlagFieldCellFilter<FlagField_T> fluidFilter( flagField_ );
   fluidFilter.addFlag( Fluid_Flag );
   filters[ "FluidFilter" ] = fluidFilter;

   field::FlagFieldCellFilter<FlagField_T> obstacleFilter( flagField_ );
   obstacleFilter.addFlag( NoSlip_Flag );
   obstacleFilter.addFlag(    UBB_Flag );
   filters[ "ObstacleFilter" ] = obstacleFilter;

   // before functions

   beforeFunctions[ "PDFGhostLayerSync" ] = pdfGhostLayerSync_;
}






////////////////
// Evaluation //
////////////////

real_t exactFlowRate( const real_t flowRate )
{
   return flowRate;
}

Vector3< real_t > exactVelocity( const Vector3< real_t > & p, const math::AABB & domain, const real_t maxLatticeVelocity )
{
   return Vector3< real_t >( maxLatticeVelocity * ( p[1] - domain.yMin() ) / domain.ySize(), real_t(0), real_t(0) );
}






////////////////////
// THE SIMULATION //
////////////////////

template< typename LatticeModel_T, typename Sweep_T >
void addRefinementTimeStep( SweepTimeloop & timeloop, shared_ptr< blockforest::StructuredBlockForest > & blocks,
                            const BlockDataID & pdfFieldId, const BlockDataID & boundaryHandlingId,
                            const shared_ptr<WcTimingPool> & timingPool, const shared_ptr<WcTimingPool> & levelwiseTimingPool,
                            const bool syncComm, const bool fullComm, const bool linearExplosion,
                            shared_ptr< Sweep_T > & sweep, const std::string & info )
{
   typedef typename MyBoundaryHandling< LatticeModel_T >::BoundaryHandling_T BH_T;

   auto ts = lbm::refinement::makeTimeStep< LatticeModel_T, BH_T >( blocks, sweep, pdfFieldId, boundaryHandlingId );
   ts->asynchronousCommunication( !syncComm );
   ts->optimizeCommunication( !fullComm );
   ts->performLinearExplosion( linearExplosion );
   ts->enableTiming( timingPool, levelwiseTimingPool );

   timeloop.addFuncBeforeTimeStep( makeSharedFunctor(ts), info );
}

template< typename LatticeModel_T, class Enable = void >
struct AddRefinementTimeStep
{
   static void add( SweepTimeloop & timeloop, shared_ptr< blockforest::StructuredBlockForest > & blocks,
                    const BlockDataID & pdfFieldId, const BlockDataID & flagFieldId, const BlockDataID & boundaryHandlingId,
                    const shared_ptr<WcTimingPool> & timingPool, const shared_ptr<WcTimingPool> & levelwiseTimingPool,
                    const bool split, const bool pure, const bool syncComm, const bool fullComm, const bool linearExplosion )
   {
      if( split )
      {
         if( pure )
         {
            typedef lbm::SplitPureSweep< LatticeModel_T > Sweep_T;
            auto mySweep = make_shared< Sweep_T >( pdfFieldId );

            addRefinementTimeStep< LatticeModel_T, Sweep_T >( timeloop, blocks, pdfFieldId, boundaryHandlingId, timingPool, levelwiseTimingPool,
                                                              syncComm, fullComm, linearExplosion, mySweep, "LBM refinement time step (split pure LB sweep)" );
         }
         else
         {
            typedef lbm::SplitSweep< LatticeModel_T, FlagField_T > Sweep_T;
            auto mySweep = make_shared< Sweep_T >( pdfFieldId, flagFieldId, Fluid_Flag );

            addRefinementTimeStep< LatticeModel_T, Sweep_T >( timeloop, blocks, pdfFieldId, boundaryHandlingId, timingPool, levelwiseTimingPool,
                                                              syncComm, fullComm, linearExplosion, mySweep, "LBM refinement time step (split LB sweep)" );
         }
      }
      else
      {
         auto mySweep = lbm::makeCellwiseSweep< LatticeModel_T, FlagField_T >( pdfFieldId, flagFieldId, Fluid_Flag );

         addRefinementTimeStep< LatticeModel_T >( timeloop, blocks, pdfFieldId, boundaryHandlingId, timingPool, levelwiseTimingPool,
                                                  syncComm, fullComm, linearExplosion, mySweep, "LBM refinement time step (default LB sweep)" );
      }
   }
};

template< typename LatticeModel_T  >
struct AddRefinementTimeStep< LatticeModel_T, typename boost::enable_if< boost::mpl::or_< boost::is_same< typename LatticeModel_T::CollisionModel::tag,
                                                                                                          lbm::collision_model::MRT_tag >,
                                                                                          boost::is_same< typename LatticeModel_T::Stencil, stencil::D3Q15 >,
                                                                                          boost::is_same< typename LatticeModel_T::Stencil, stencil::D3Q27 > > >::type >
{
   static void add( SweepTimeloop & timeloop, shared_ptr< blockforest::StructuredBlockForest > & blocks,
                    const BlockDataID & pdfFieldId, const BlockDataID & flagFieldId, const BlockDataID & boundaryHandlingId,
                    const shared_ptr<WcTimingPool> & timingPool, const shared_ptr<WcTimingPool> & levelwiseTimingPool,
                    const bool /*split*/, const bool /*pure*/, const bool syncComm, const bool fullComm, const bool linearExplosion )
   {
      auto mySweep = lbm::makeCellwiseSweep< LatticeModel_T, FlagField_T >( pdfFieldId, flagFieldId, Fluid_Flag );

      addRefinementTimeStep< LatticeModel_T >( timeloop, blocks, pdfFieldId, boundaryHandlingId, timingPool, levelwiseTimingPool,
                                               syncComm, fullComm, linearExplosion, mySweep, "LBM refinement time step (default LB sweep)" );
   }
};



template< typename LatticeModel_T >
void run( const shared_ptr< Config > & config, const LatticeModel_T & latticeModel,
          const bool split, const bool pure, const bool fzyx, const bool syncComm, const bool fullComm, const bool linearExplosion,
          const blockforest::RefinementSelectionFunctions & refinementSelectionFunctions, Setup & setup,
          const memory_t memoryPerCell, const memory_t processMemoryLimit )
{
   Config::BlockHandle configBlock = config->getBlock( "CouetteFlow" );

   setup.viscosity_L    = latticeModel.collisionModel().viscosity( uint_t(0) );
   setup.meanVelocity_L = ( setup.Re * setup.viscosity_L ) / real_c( setup.yBlocks * setup.yCells );
   setup.maxVelocity_L  = real_t(2) * setup.meanVelocity_L;
   setup.flowRate_L     = ( setup.maxVelocity_L * real_c( setup.yBlocks * setup.yCells ) * real_c( setup.zBlocks * setup.zCells ) ) / real_t(2);

   // creating the block structure

   auto blocks = createStructuredBlockForest( configBlock, refinementSelectionFunctions, setup, memoryPerCell, processMemoryLimit );

   // add pdf field to blocks

   const real_t initVelocity = ( configBlock.getParameter< bool >( "initWithMeanVelocity", false ) ) ? setup.meanVelocity_L : real_t(0);

   BlockDataID pdfFieldId = fzyx ? lbm::addPdfFieldToStorage( blocks, "pdf field (fzyx)", latticeModel,
                                                              Vector3< real_t >( initVelocity, real_c(0), real_c(0) ), real_t(1),
                                                              FieldGhostLayers, field::fzyx ) :
                                   lbm::addPdfFieldToStorage( blocks, "pdf field (zyxf)", latticeModel,
                                                              Vector3< real_t >( initVelocity, real_c(0), real_c(0) ), real_t(1),
                                                              FieldGhostLayers, field::zyxf );

   typedef typename lbm::Adaptor<LatticeModel_T>::VelocityVector  VelocityAdaptor_T;
   typedef typename lbm::Adaptor<LatticeModel_T>::Density          DensityAdaptor_T;
   BlockDataID velocityAdaptorId = field::addFieldAdaptor< VelocityAdaptor_T >( blocks, pdfFieldId, "velocity adaptor" );
   BlockDataID  densityAdaptorId = field::addFieldAdaptor<  DensityAdaptor_T >( blocks, pdfFieldId, "density adaptor" );

   // add flag field to blocks

   BlockDataID flagFieldId = field::addFlagFieldToStorage< FlagField_T >( blocks, "flag field", FieldGhostLayers );

   // add LB boundary handling to blocks

   BlockDataID boundaryHandlingId = blocks->template addBlockData< typename MyBoundaryHandling< LatticeModel_T >::BoundaryHandling_T >(
            MyBoundaryHandling< LatticeModel_T >( flagFieldId, pdfFieldId, setup.maxVelocity_L ), "boundary handling" );

   setFlags< LatticeModel_T, typename MyBoundaryHandling< LatticeModel_T >::BoundaryHandling_T >( blocks, boundaryHandlingId );

   // creating the time loop

   const uint_t outerTimeSteps = configBlock.getParameter< uint_t >( "outerTimeSteps", uint_c(10) );
   const uint_t innerTimeSteps = configBlock.getParameter< uint_t >( "innerTimeSteps", uint_c(10) );

   SweepTimeloop timeloop( blocks->getBlockStorage(), ( outerTimeSteps * innerTimeSteps ) + uint_t(1) );

   // VTK

   blockforest::communication::NonUniformBufferedScheme< stencil::D3Q27 > pdfGhostLayerSync( blocks );
   pdfGhostLayerSync.addPackInfo( make_shared< lbm::refinement::PdfFieldSyncPackInfo< LatticeModel_T > >( pdfFieldId ) );

   MyVTKOutput< LatticeModel_T > myVTKOutput( pdfFieldId, flagFieldId, pdfGhostLayerSync, setup );

   std::map< std::string, vtk::SelectableOutputFunction > vtkOutputFunctions;
   vtk::initializeVTKOutput( vtkOutputFunctions, myVTKOutput, blocks, config );   
   
   const bool vtkBeforeTimeStep = configBlock.getParameter< bool >( "vtkBeforeTimeStep", true );

   if( vtkBeforeTimeStep )
   {
      for( auto output = vtkOutputFunctions.begin(); output != vtkOutputFunctions.end(); ++output )
         timeloop.addFuncBeforeTimeStep( output->second.outputFunction, std::string("VTK: ") + output->first,
                                         output->second.requiredGlobalStates, output->second.incompatibleGlobalStates );
   }

   // add 'refinement' LB time step to time loop

   shared_ptr<WcTimingPool> refinementTimeStepTiming = make_shared<WcTimingPool>();
   shared_ptr<WcTimingPool> refinementTimeStepLevelwiseTiming = make_shared<WcTimingPool>();

   AddRefinementTimeStep< LatticeModel_T >::add( timeloop, blocks, pdfFieldId, flagFieldId, boundaryHandlingId, refinementTimeStepTiming,
                                                 refinementTimeStepLevelwiseTiming, split, pure, syncComm, fullComm, linearExplosion );

   // evaluation

   const auto exactSolutionFunction = boost::bind( exactVelocity, _1, blocks->getDomain(), setup.maxVelocity_L );

   auto volumetricFlowRate = field::makeVolumetricFlowRateEvaluation< VelocityAdaptor_T, FlagField_T >( configBlock, blocks, velocityAdaptorId,
                                                                                                        flagFieldId, Fluid_Flag,
                                                                                                        boost::bind( exactFlowRate, setup.flowRate_L ),
                                                                                                        exactSolutionFunction );
   volumetricFlowRate->setNormalizationFactor( real_t(1) / setup.maxVelocity_L );
   volumetricFlowRate->setDomainNormalization( Vector3<real_t>( real_t(1) ) );

   timeloop.addFuncBeforeTimeStep( makeSharedFunctor( volumetricFlowRate ), "volumetric flow rate evaluation" );

   auto accuracyEvaluation = field::makeAccuracyEvaluation< VelocityAdaptor_T >( configBlock, blocks, velocityAdaptorId, exactSolutionFunction );
   accuracyEvaluation->setNormalizationFactor( real_t(1) / setup.maxVelocity_L );

   timeloop.addFuncBeforeTimeStep( makeSharedFunctor( accuracyEvaluation ), "accuracy evaluation" );

   auto linePlot = field::makeAccuracyEvaluationLinePlot< VelocityAdaptor_T >( configBlock, blocks, velocityAdaptorId, exactSolutionFunction );
   linePlot->setNormalizationFactor( real_t(1) / setup.maxVelocity_L );

   timeloop.addFuncBeforeTimeStep( makeSharedFunctor( field::makeAccuracyEvaluationLinePlotter( configBlock, linePlot ) ), "accuracy evaluation (line plot)" );

   timeloop.addFuncBeforeTimeStep( makeSharedFunctor( lbm::makeMassEvaluation< DensityAdaptor_T >( configBlock, blocks, uint_t(0), densityAdaptorId ) ), "mass evaluation" );

   // stability check (non-finite values in the PDF field?)

   timeloop.addFuncAfterTimeStep( makeSharedFunctor( field::makeStabilityChecker< lbm::PdfField< LatticeModel_T >, FlagField_T >(
                                                        configBlock, blocks, pdfFieldId, flagFieldId, Fluid_Flag ) ),
                                  "LBM stability check" );

   // VTK

   if( !vtkBeforeTimeStep )
   {
      for( auto output = vtkOutputFunctions.begin(); output != vtkOutputFunctions.end(); ++output )
         timeloop.addFuncAfterTimeStep( output->second.outputFunction, std::string("VTK: ") + output->first,
                                        output->second.requiredGlobalStates, output->second.incompatibleGlobalStates );
   }

   // remaining time logger

   const double remainingTimeLoggerFrequency = configBlock.getParameter< double >( "remainingTimeLoggerFrequency", 3.0 );
   timeloop.addFuncAfterTimeStep( timing::RemainingTimeLogger( timeloop.getNrOfTimeSteps(), remainingTimeLoggerFrequency ), "Remaining time logger" );

   // logging right before the simulation starts

   lbm::BlockForestEvaluation< FlagField_T > blockForest( blocks, flagFieldId, Fluid_Flag );
   blockForest.logInfoOnRoot();

   field::CellCounter< FlagField_T > fluidCells( blocks, flagFieldId, Fluid_Flag );
   fluidCells();

   WALBERLA_LOG_INFO_ON_ROOT( "Benchmark run data:"
                              "\n- simulation parameters:"
                              "\n   + collision model:  " << CollisionModelString< LatticeModel_T >::str() <<
                              "\n   + stencil:          " << StencilString< LatticeModel_T >::str() <<
                              "\n   + compressible:     " << ( LatticeModel_T::compressible ? "yes" : "no" ) <<
                              "\n   + split kernel:     " << ( split ? "yes" : "no" ) <<
                              "\n   + pure kernel:      " << ( pure ? "yes (collision is also performed within obstacle cells)" : "no" ) <<
                              "\n   + data layout:      " << ( fzyx ? "fzyx (structure of arrays [SoA])" : "zyxf (array of structures [AoS])" ) <<
                              "\n   + communication:    " << ( fullComm ? ( syncComm ? "synchronous, full synchronization" : "asynchronous, full synchronization" ) :
                                                                          ( syncComm ? "synchronous, block neighborhood and direction-aware optimizations" : "asynchronous, block neighborhood and direction-aware optimizations" ) ) <<
                              "\n   + linear explosion: " << ( linearExplosion ? "yes" : "no" ) <<
                              "\n- simulation properties:"
                              "\n   + fluid cells:      " << fluidCells.numberOfCells() << " (in total on all levels)" <<
                              "\n   + Reynolds number:  " << setup.Re <<
                              "\n   + kin. viscosity:   " << setup.viscosity_L << " (on the coarsest grid)" <<
                              "\n   + maximum velocity: " << setup.maxVelocity_L <<
                              "\n   + mean velocity:    " << setup.meanVelocity_L <<
                              "\n   + flow rate:        " << setup.flowRate_L << " (in lattice units of the coarsest grid)" <<
                              "\n   + #time steps:      " << timeloop.getNrOfTimeSteps() << " (on the coarsest grid)" );

   // run the simulation

   lbm::PerformanceEvaluation< FlagField_T > performance( blocks, flagFieldId, Fluid_Flag );

   for( uint_t outerRun = 0; outerRun < outerTimeSteps; ++outerRun )
   {
      WcTimingPool timeloopTiming;

      WALBERLA_MPI_WORLD_BARRIER();
      WcTimer timer;
      timer.start();

      for( uint_t innerRun = 0; innerRun < innerTimeSteps; ++innerRun )
         timeloop.singleStep( timeloopTiming );

      timer.end();

      double time = timer.max();
      mpi::reduceInplace( time, mpi::MAX );

      const auto reducedTimeloopTiming = timeloopTiming.getReduced();
      const auto reducedRTSTiming  = refinementTimeStepTiming->getReduced();
      const auto reducedRTSLTiming = refinementTimeStepLevelwiseTiming->getReduced();
      refinementTimeStepTiming->clear();
      refinementTimeStepLevelwiseTiming->clear();

      WALBERLA_LOG_RESULT_ON_ROOT( "Time loop timing:\n" << *reducedTimeloopTiming );
      WALBERLA_LOG_RESULT_ON_ROOT( "Refinement time step timing:\n" << *reducedRTSTiming );
      WALBERLA_LOG_RESULT_ON_ROOT( "Refinement time step timing (one timer per level):\n" << *reducedRTSLTiming );

      performance.logResultOnRoot( innerTimeSteps, time );

      WALBERLA_ROOT_SECTION()
      {
         // logging in SQL database

         if( configBlock.getParameter< bool >( "logToSqlDB", true ) )
         {
            const std::string sqlFile = configBlock.getParameter< std::string >( "sqlFile", "performance.sqlite" );

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

            performance.getResultsForSQLOnRoot( integerProperties, realProperties, stringProperties, innerTimeSteps, time );
            blockForest.getResultsForSQLOnRoot( integerProperties, realProperties, stringProperties );

            stringProperties[ "collisionModel" ]    = CollisionModelString< LatticeModel_T >::str();
            stringProperties[ "stencil" ]           = StencilString< LatticeModel_T >::str();
            stringProperties[ "compressible" ]      = ( LatticeModel_T::compressible ? "yes" : "no" );
            stringProperties[ "splitKernel" ]       = ( split ? "yes" : "no" );
            stringProperties[ "pureKernel" ]        = ( pure ? "yes" : "no" );
            stringProperties[ "dataLayout" ]        = ( fzyx ? "fzyx" : "zyxf" );
            stringProperties[ "syncCommunication" ] = ( syncComm ? "yes" : "no" );
            stringProperties[ "fullCommunication" ] = ( fullComm ? "yes" : "no" );
            stringProperties[ "linearExplosion" ]   = ( linearExplosion ? "yes" : "no" );

            realProperties[ "Re" ]             = double_c( setup.Re );
            realProperties[ "viscosity_L" ]    = double_c( setup.viscosity_L );
            realProperties[ "maxVelocity_L" ]  = double_c( setup.maxVelocity_L );
            realProperties[ "meanVelocity_L" ] = double_c( setup.meanVelocity_L );
            realProperties[ "flowRate_L" ]     = double_c( setup.flowRate_L );
            
            integerProperties[ "simulationTimeSteps" ] = int_c( timeloop.getNrOfTimeSteps() );
            
            realProperties[ "simulationProgress" ] = double_c( ( outerRun + uint_t(1) ) * innerTimeSteps ) / double_c( outerTimeSteps * innerTimeSteps );

            auto runId = postprocessing::storeRunInSqliteDB( sqlFile, integerProperties, stringProperties, realProperties );
            postprocessing::storeTimingPoolInSqliteDB( sqlFile, runId, *reducedTimeloopTiming, "Timeloop" );
            postprocessing::storeTimingPoolInSqliteDB( sqlFile, runId, *reducedRTSTiming, "RefinementTimeStep" );
            postprocessing::storeTimingPoolInSqliteDB( sqlFile, runId, *reducedRTSLTiming, "RefinementTimeStepLevelwise" );
         }
      }
   }

   // Do one more step so that we see the final output
   timeloop.singleStep();

   // logging once again at the end of the simulation, identical to logging at the beginning :-)

   blockForest.logInfoOnRoot();

   WALBERLA_LOG_INFO_ON_ROOT( "Benchmark run data:"
                              "\n- simulation parameters:"
                              "\n   + collision model:  " << CollisionModelString< LatticeModel_T >::str() <<
                              "\n   + stencil:          " << StencilString< LatticeModel_T >::str() <<
                              "\n   + compressible:     " << ( LatticeModel_T::compressible ? "yes" : "no" ) <<
                              "\n   + split kernel:     " << ( split ? "yes" : "no" ) <<
                              "\n   + pure kernel:      " << ( pure ? "yes (collision is also performed within obstacle cells)" : "no" ) <<
                              "\n   + data layout:      " << ( fzyx ? "fzyx (structure of arrays [SoA])" : "zyxf (array of structures [AoS])" ) <<
                              "\n   + communication:    " << ( fullComm ? ( syncComm ? "synchronous, full synchronization" : "asynchronous, full synchronization" ) :
                                                                          ( syncComm ? "synchronous, block neighborhood and direction-aware optimizations" : "asynchronous, block neighborhood and direction-aware optimizations" ) ) <<
                              "\n   + linear explosion: " << ( linearExplosion ? "yes" : "no" ) <<
                              "\n- simulation properties:"
                              "\n   + fluid cells:      " << fluidCells.numberOfCells() << " (in total on all levels)" <<
                              "\n   + Reynolds number:  " << setup.Re <<
                              "\n   + kin. viscosity:   " << setup.viscosity_L << " (on the coarsest grid)" <<
                              "\n   + maximum velocity: " << setup.maxVelocity_L <<
                              "\n   + mean velocity:    " << setup.meanVelocity_L <<
                              "\n   + flow rate:        " << setup.flowRate_L << " (in lattice units of the coarsest grid)" <<
                              "\n   + #time steps:      " << timeloop.getNrOfTimeSteps() << " (on the coarsest grid)" );

   auto accuracyEvaluationLinePlotBlock = configBlock.getBlock("AccuracyEvaluationLinePlot");
   if( accuracyEvaluationLinePlotBlock && accuracyEvaluationLinePlotBlock.isDefined("filename") )
      (*linePlot)( accuracyEvaluationLinePlotBlock.template getParameter<std::string>("filename") );

   WALBERLA_ROOT_SECTION()
   {
      if( configBlock.getParameter< bool >( "check", false ) )
      {
         if( configBlock.isDefined( "checkFlowRateError" ) )
         {
            const real_t checkFlowRateError = configBlock.getParameter< real_t >( "checkFlowRateError" );
            WALBERLA_CHECK_LESS( std::abs( ( volumetricFlowRate->flowRate() - volumetricFlowRate->solution() ) / volumetricFlowRate->solution() ),
                                 checkFlowRateError );
         }
         if( configBlock.isDefined( "checkErrorL1" ) )
         {
            const real_t checkErrorL1 = configBlock.getParameter< real_t >( "checkErrorL1" );
            WALBERLA_CHECK_LESS( accuracyEvaluation->L1(), checkErrorL1 );
         }
         if( configBlock.isDefined( "checkErrorL2" ) )
         {
            const real_t checkErrorL2 = configBlock.getParameter< real_t >( "checkErrorL2" );
            WALBERLA_CHECK_LESS( accuracyEvaluation->L2(), checkErrorL2 );
         }
         if( configBlock.isDefined( "checkErrorLmax" ) )
         {
            const real_t checkErrorLmax = configBlock.getParameter< real_t >( "checkErrorLmax" );
            WALBERLA_CHECK_LESS( accuracyEvaluation->Lmax(), checkErrorLmax );
         }
      }
   }
}



//////////
// MAIN //
//////////

enum CM { CMSRT, CMTRT, CMMRT };
enum LM { LMD3Q15, LMD3Q19, LMD3Q27 };

int main( int argc, char **argv )
{
   mpi::Environment env( argc, argv );
   
   if( argc < 2 )
   {
      WALBERLA_ROOT_SECTION()
      {
         std::cout << "Usage: " << argv[0] << " path-to-configuration-file [--trt | --mrt] [--d3q15 | --d3q27] [--comp] [--split [--pure]] [--fzyx] [--sync-comm] [--full-comm] [--linear-exp]\n"
                      "\n"
                      "By default, SRT is selected as collision model, an asynchronous communication scheme with block neighborhood and\n"
                      "direction-aware optimizations is chosen, and an incompressible, basic D3Q19 LB kernel is executed on a PDF field with\n"
                      "layout 'zyxf' (= array of structures [AoS]).\n"
                      "\n"
                      "Optional arguments:\n"
                      " --trt:        collision model = TRT\n"
                      " --mrt:        collision model = MRT\n"
                      " --d3q15:      A D3Q15 model is used.\n"
                      " --d3q27:      A D3Q27 model is used.\n"
                      " --comp:       LB kernel is switched from incompressible to compressible\n"
                      " --split:      LB kernel split by PDF direction\n"
                      "               Should always be combined with --fzyx.\n"
                      " --pure:       LB kernel is executed in every cell (including obstacle/boundary cells)\n"
                      "               Only available in combination with --split.\n"
                      " --fzyx:       data layout switched to 'fzyx' (structure of arrays [SoA])\n"
                      " --sync-comm:  A synchronous communication scheme is used instead of an asynchronous scheme\n"
                      "               which is used by default.\n"
                      " --full-comm:  A full synchronization of neighboring blocks is performed instead of using a communication\n"
                      "               that uses block neighborhood and direction-aware optimizations.\n"
                      " --linear-exp: When communicating from coarse to fine grids, a linear interpolation scheme is used\n"
                      "               instead of a uniform distribution of a coarse cell to eight fine cells.\n"
                      "\n"
                      "Please note: Depending on the underlying hardware and the configuration of the benchmark (more precisely: the number of cells\n"
                      "             in each block), the best performance may be achieved with split, pure kernels combined with a structure of arrays\n"
                      "             ('fzyx') data layout!" << std::endl;
      }
      return EXIT_SUCCESS;
   }   

   logging::Logging::printHeaderOnStream();
   //WALBERLA_ROOT_SECTION() { logging::Logging::instance()->setLogLevel( logging::Logging::PROGRESS ); }

#ifdef _OPENMP
   if( std::getenv( "OMP_NUM_THREADS" ) == NULL )
      WALBERLA_ABORT( "If you are using a version of the program that was compiled with OpenMP you have to "
                      "specify the environment variable \'OMP_NUM_THREADS\' accordingly!" );
#endif
For faster browsing, not all history is shown. View entire blame