SphereWallCollision.cpp 69.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
//======================================================================================================================
//
//  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 SphereWallCollision.cpp
//! \ingroup lbm_mesapd_coupling
//! \author Christoph Rettinger <christoph.rettinger@fau.de>
//
//======================================================================================================================

#include "blockforest/Initialization.h"
#include "blockforest/communication/UniformBufferedScheme.h"

#include "boundary/all.h"

#include "core/waLBerlaBuildInfo.h"
#include "core/DataTypes.h"
#include "core/Environment.h"
#include "core/SharedFunctor.h"
#include "core/debug/Debug.h"
#include "core/debug/TestSubsystem.h"
#include "core/math/all.h"
#include "core/timing/RemainingTimeLogger.h"
#include "core/logging/all.h"
#include "core/mpi/Broadcast.h"
#include "core/mpi/Reduce.h"

#include "domain_decomposition/SharedSweep.h"

#include "field/AddToStorage.h"
#include "field/adaptors/AdaptorCreators.h"
#include "field/communication/PackInfo.h"
#include "field/interpolators/FieldInterpolatorCreators.h"
#include "field/interpolators/NearestNeighborFieldInterpolator.h"

#include "lbm/boundary/all.h"
#include "lbm/communication/PdfFieldPackInfo.h"
#include "lbm/field/AddToStorage.h"
#include "lbm/field/Adaptors.h"
#include "lbm/field/PdfField.h"
#include "lbm/lattice_model/D3Q19.h"
#include "lbm/sweeps/CellwiseSweep.h"
#include "lbm/sweeps/SweepWrappers.h"

#include "lbm_mesapd_coupling/mapping/ParticleMapping.h"
#include "lbm_mesapd_coupling/momentum_exchange_method/MovingParticleMapping.h"
#include "lbm_mesapd_coupling/momentum_exchange_method/boundary/SimpleBB.h"
#include "lbm_mesapd_coupling/momentum_exchange_method/boundary/CurvedLinear.h"
#include "lbm_mesapd_coupling/momentum_exchange_method/reconstruction/Reconstructor.h"
#include "lbm_mesapd_coupling/momentum_exchange_method/reconstruction/ExtrapolationDirectionFinder.h"
#include "lbm_mesapd_coupling/momentum_exchange_method/reconstruction/PdfReconstructionManager.h"
#include "lbm_mesapd_coupling/utility/AddForceOnParticlesKernel.h"
#include "lbm_mesapd_coupling/utility/ParticleSelector.h"
#include "lbm_mesapd_coupling/DataTypes.h"
#include "lbm_mesapd_coupling/utility/AverageHydrodynamicForceTorqueKernel.h"
#include "lbm_mesapd_coupling/utility/AddHydrodynamicInteractionKernel.h"
#include "lbm_mesapd_coupling/utility/ResetHydrodynamicForceTorqueKernel.h"
#include "lbm_mesapd_coupling/utility/LubricationCorrectionKernel.h"
#include "lbm_mesapd_coupling/utility/OmegaBulkAdaption.h"
71
#include "lbm_mesapd_coupling/utility/InitializeHydrodynamicForceTorqueForAveragingKernel.h"
72
73
74
75
76
77
78
79
80
81
82

#include "mesa_pd/collision_detection/AnalyticContactDetection.h"
#include "mesa_pd/data/ParticleAccessorWithShape.h"
#include "mesa_pd/data/ParticleStorage.h"
#include "mesa_pd/data/ShapeStorage.h"
#include "mesa_pd/data/DataTypes.h"
#include "mesa_pd/data/shape/HalfSpace.h"
#include "mesa_pd/data/shape/Sphere.h"
#include "mesa_pd/domain/BlockForestDomain.h"
#include "mesa_pd/domain/BlockForestDataHandling.h"
#include "mesa_pd/kernel/DoubleCast.h"
83
#include "mesa_pd/kernel/ExplicitEuler.h"
84
85
86
#include "mesa_pd/kernel/LinearSpringDashpot.h"
#include "mesa_pd/kernel/NonLinearSpringDashpot.h"
#include "mesa_pd/kernel/ParticleSelector.h"
87
#include "mesa_pd/kernel/VelocityVerlet.h"
88
#include "mesa_pd/mpi/ClearNextNeighborSync.h"
89
90
91
92
93
#include "mesa_pd/mpi/SyncNextNeighbors.h"
#include "mesa_pd/mpi/ReduceProperty.h"
#include "mesa_pd/mpi/ReduceContactHistory.h"
#include "mesa_pd/mpi/ContactFilter.h"
#include "mesa_pd/mpi/notifications/ForceTorqueNotification.h"
94
#include "mesa_pd/mpi/notifications/HydrodynamicForceTorqueNotification.h"
95
96
97
98
99
100
101
102
103
104
#include "mesa_pd/vtk/ParticleVtkOutput.h"

#include "timeloop/SweepTimeloop.h"

#include "vtk/all.h"
#include "field/vtk/all.h"
#include "lbm/vtk/all.h"

#include <functional>

105
106
107
#ifdef WALBERLA_BUILD_WITH_CODEGEN
#include "GeneratedLBM.h"
#endif
108
109
110
111
112
113
114
115
116
117
118

namespace sphere_wall_collision
{

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

using namespace walberla;
using walberla::uint_t;

119
#ifdef WALBERLA_BUILD_WITH_CODEGEN
120
using LatticeModel_T = lbm::GeneratedLBM;
121
122
123
#else
using LatticeModel_T = lbm::D3Q19< lbm::collision_model::D3Q19MRT>;
#endif
124
125
126
127
128
129
130
131
132
133
134

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

using flag_t = walberla::uint8_t;
using FlagField_T = FlagField<flag_t>;

using ScalarField_T = GhostLayerField< real_t, 1>;

const uint_t FieldGhostLayers = 1;

135
136
137
138
#define USE_TRTLIKE
//#define USE_D3Q27TRTLIKE
//#define USE_CUMULANT

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
///////////
// FLAGS //
///////////

const FlagUID Fluid_Flag( "fluid" );
const FlagUID NoSlip_Flag( "no slip" );
const FlagUID MO_Flag( "moving obstacle" );
const FlagUID FormerMO_Flag( "former moving obstacle" );
const FlagUID SimplePressure_Flag("simple pressure");

/////////////////////////////////////
// BOUNDARY HANDLING CUSTOMIZATION //
/////////////////////////////////////
template <typename ParticleAccessor_T>
class MyBoundaryHandling
{
public:

   using NoSlip_T = lbm::NoSlip< LatticeModel_T, flag_t >;
   using MO_T = lbm_mesapd_coupling::CurvedLinear< LatticeModel_T, FlagField_T, ParticleAccessor_T >;
   using SimplePressure_T = lbm::SimplePressure< LatticeModel_T, flag_t >;
   using Type = BoundaryHandling< FlagField_T, Stencil_T, NoSlip_T, MO_T, SimplePressure_T >;

   MyBoundaryHandling( const BlockDataID & flagFieldID, const BlockDataID & pdfFieldID,
                       const BlockDataID & particleFieldID, const shared_ptr<ParticleAccessor_T>& ac,
                       bool applyOutflowBCAtTop) :
         flagFieldID_( flagFieldID ), pdfFieldID_( pdfFieldID ), particleFieldID_( particleFieldID ), ac_( ac ), applyOutflowBCAtTop_(applyOutflowBCAtTop) {}

   Type * operator()( IBlock* const block, const StructuredBlockStorage* const storage ) const
   {
      WALBERLA_ASSERT_NOT_NULLPTR( block );
      WALBERLA_ASSERT_NOT_NULLPTR( storage );

      auto * flagField     = block->getData< FlagField_T >( flagFieldID_ );
      auto *  pdfField     = block->getData< PdfField_T > ( pdfFieldID_ );
      auto * particleField = block->getData< lbm_mesapd_coupling::ParticleField_T > ( particleFieldID_ );

      const auto fluid = flagField->flagExists( Fluid_Flag ) ? flagField->getFlag( Fluid_Flag ) : flagField->registerFlag( Fluid_Flag );

      Type * handling = new Type( "moving obstacle boundary handling", flagField, fluid,
                                  NoSlip_T( "NoSlip", NoSlip_Flag, pdfField ),
                                  MO_T( "MO", MO_Flag, pdfField, flagField, particleField, ac_, fluid, *storage, *block ),
                                  SimplePressure_T( "SimplePressure", SimplePressure_Flag, pdfField, real_t(1) ) );

      if(applyOutflowBCAtTop_)
      {
         const auto simplePressure = flagField->getFlag( SimplePressure_Flag );

         CellInterval domainBB = storage->getDomainCellBB();
         domainBB.xMin() -= cell_idx_c( FieldGhostLayers );
         domainBB.xMax() += cell_idx_c( FieldGhostLayers );

         domainBB.yMin() -= cell_idx_c( FieldGhostLayers );
         domainBB.yMax() += cell_idx_c( FieldGhostLayers );

         domainBB.zMin() -= cell_idx_c( FieldGhostLayers );
         domainBB.zMax() += cell_idx_c( FieldGhostLayers );

         // TOP
         CellInterval top( domainBB.xMin(), domainBB.yMin(), domainBB.zMax(), domainBB.xMax(), domainBB.yMax(), domainBB.zMax() );
         storage->transformGlobalToBlockLocalCellInterval( top, *block );
         handling->forceBoundary( simplePressure, top );
      }

      handling->fillWithDomain( FieldGhostLayers );

      return handling;
   }

private:

   const BlockDataID flagFieldID_;
   const BlockDataID pdfFieldID_;
   const BlockDataID particleFieldID_;

   shared_ptr<ParticleAccessor_T> ac_;

   bool applyOutflowBCAtTop_;
};
//*******************************************************************************************************************

//*******************************************************************************************************************
/*!\brief Evaluating the position and velocity of the sphere
 *
 */
//*******************************************************************************************************************
template< typename ParticleAccessor_T>
class SpherePropertyLogger
{
public:
   SpherePropertyLogger( const shared_ptr< ParticleAccessor_T > & ac, walberla::id_t sphereUid,
                         const std::string & fileName, bool fileIO,
                         real_t diameter, real_t uIn, uint_t numberOfSubSteps) :
      ac_( ac ), sphereUid_( sphereUid ), fileName_( fileName ), fileIO_(fileIO),
      diameter_( diameter ), uIn_( uIn ), numberOfSubSteps_(numberOfSubSteps),
      gap_( real_t(0) ), settlingVelocity_( real_t(0) )
   {
      if ( fileIO_ )
      {
         WALBERLA_ROOT_SECTION()
         {
            std::ofstream file;
            file.open( fileName_.c_str() );
            file << "#\t t/tref\t gapSize\t gapSize/D\t velZ\t velZ/uIn\n";
            file.close();
         }
      }
   }

   void operator()(const uint_t timeStep, const uint_t subStep )
   {
      real_t curTimestep = real_c(timeStep) + real_c(subStep) / real_c(numberOfSubSteps_);

      Vector3<real_t> pos(real_t(0));
      Vector3<real_t> transVel(real_t(0));

      size_t idx = ac_->uidToIdx(sphereUid_);
      if( idx != std::numeric_limits<size_t>::max())
      {
         if(!mesa_pd::data::particle_flags::isSet( ac_->getFlags(idx), mesa_pd::data::particle_flags::GHOST))
         {
            pos = ac_->getPosition(idx);
            transVel = ac_->getLinearVelocity(idx);
         }
      }

      WALBERLA_MPI_SECTION()
      {
         mpi::allReduceInplace( pos[0], mpi::SUM );
         mpi::allReduceInplace( pos[1], mpi::SUM );
         mpi::allReduceInplace( pos[2], mpi::SUM );
         mpi::allReduceInplace( transVel[2], mpi::SUM );
      }

      position_ = pos;
      gap_ = pos[2] - real_t(0.5) * diameter_;
      settlingVelocity_ = transVel[2];

      WALBERLA_ROOT_SECTION()
      {
         // coarsen the output to the result file
         if(subStep == 0)
         {
            gapLogging_.push_back(gap_);
            velocityLogging_.push_back(settlingVelocity_);
         }
      }

      if( fileIO_ )
         writeToFile( curTimestep, gap_, settlingVelocity_);
   }

   real_t getGapSize() const
   {
      return gap_;
   }

   real_t getSettlingVelocity() const
   {
      return settlingVelocity_;
   }

   Vector3<real_t> getPosition() const
   {
      return position_;
   }

306
   void writeResult(const std::string & filename, uint_t tImpact)
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
   {
      WALBERLA_ROOT_SECTION()
      {
         WALBERLA_CHECK_EQUAL(gapLogging_.size(), velocityLogging_.size());
         std::ofstream file;
         file.open( filename.c_str());

         file << "#\t t\t t/tref\t gapSize\t gapSize/D\t velZ\t velZ/uIn\n";
         real_t tref = diameter_ / uIn_;

         for( uint_t i = uint_t(0); i < gapLogging_.size(); ++i)
         {
            int timestep = int(i) - int(tImpact);
            file << timestep << "\t" << real_c(timestep) / tref<< "\t"
                 << "\t" << gapLogging_[i] << "\t" << gapLogging_[i] / diameter_
                 << "\t" << velocityLogging_[i] << "\t" << velocityLogging_[i] / uIn_
                 << "\n";
         }

         file.close();
      }
   }


private:
   void writeToFile( real_t timestep, real_t gap, real_t velocity )
   {
      WALBERLA_ROOT_SECTION()
      {
         std::ofstream file;
         file.open( fileName_.c_str(), std::ofstream::app );

         real_t tref = diameter_ / uIn_;

         file << timestep << "\t" << timestep / tref<< "\t"
              << "\t" << gap << "\t" << gap / diameter_
              << "\t" << velocity << "\t" << velocity / uIn_
              << "\n";
         file.close();
      }
   }


   shared_ptr< ParticleAccessor_T > ac_;
   const walberla::id_t sphereUid_;
   std::string fileName_;
   bool fileIO_;
   real_t diameter_, uIn_;
   uint_t numberOfSubSteps_;

   real_t gap_, settlingVelocity_;
   Vector3<real_t> position_;

   std::vector<real_t> gapLogging_, velocityLogging_;

};


void createPlaneSetup(const shared_ptr<mesa_pd::data::ParticleStorage> & ps, const shared_ptr<mesa_pd::data::ShapeStorage> & ss, const math::AABB & simulationDomain, bool applyOutflowBCAtTop )
{
   // create bounding planes
   mesa_pd::data::Particle&& p0 = *ps->create(true);
   p0.setPosition(simulationDomain.minCorner());
370
   p0.setInteractionRadius(std::numeric_limits<real_t>::infinity());
371
372
373
374
375
376
377
378
379
380
381
   p0.setShapeID(ss->create<mesa_pd::data::HalfSpace>( Vector3<real_t>(0,0,1) ));
   p0.setOwner(mpi::MPIManager::instance()->rank());
   p0.setType(0);
   mesa_pd::data::particle_flags::set(p0.getFlagsRef(), mesa_pd::data::particle_flags::INFINITE);
   mesa_pd::data::particle_flags::set(p0.getFlagsRef(), mesa_pd::data::particle_flags::FIXED);

   if(!applyOutflowBCAtTop)
   {
      //only create top plane when no outflow BC should be set there
      mesa_pd::data::Particle&& p1 = *ps->create(true);
      p1.setPosition(simulationDomain.maxCorner());
382
      p1.setInteractionRadius(std::numeric_limits<real_t>::infinity());
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
      p1.setShapeID(ss->create<mesa_pd::data::HalfSpace>( Vector3<real_t>(0,0,-1) ));
      p1.setOwner(mpi::MPIManager::instance()->rank());
      p1.setType(0);
      mesa_pd::data::particle_flags::set(p1.getFlagsRef(), mesa_pd::data::particle_flags::INFINITE);
      mesa_pd::data::particle_flags::set(p1.getFlagsRef(), mesa_pd::data::particle_flags::FIXED);
   }

}

struct NonCollidingSphereSelector
{
   explicit NonCollidingSphereSelector( real_t radius ) : radius_(radius) {}

   template< typename ParticleAccessor_T >
   bool inline operator()(const size_t particleIdx, const ParticleAccessor_T & ac) const
   {
      static_assert(std::is_base_of<mesa_pd::data::IAccessor, ParticleAccessor_T>::value, "Provide a valid accessor as template");
      return !mesa_pd::data::particle_flags::isSet( ac.getFlags(particleIdx), mesa_pd::data::particle_flags::FIXED) &&
             !mesa_pd::data::particle_flags::isSet( ac.getFlags(particleIdx), mesa_pd::data::particle_flags::GLOBAL) &&
             ac.getPosition(particleIdx)[2] > radius_;
   }
private:
   real_t radius_;
};


class ForcesOnSphereLogger
{
public:
   ForcesOnSphereLogger( const std::string & fileName, bool fileIO, real_t weightForce, uint_t numberOfSubSteps ) :
         fileName_( fileName ), fileIO_( fileIO ), weightForce_( weightForce ), numberOfSubSteps_( numberOfSubSteps )
   {
      if ( fileIO_ )
      {
         WALBERLA_ROOT_SECTION()
         {
            std::ofstream file;
            file.open( fileName_.c_str() );
            file << "#t\t BuoyAndGravF\t HydInteractionF\t LubricationF\t CollisionForce\n";
            file.close();
         }
      }
   }

   void operator()(const uint_t timeStep, const uint_t subStep, real_t hydForce, real_t lubForce, real_t collisionForce )
   {
      real_t curTimestep = real_c(timeStep) + real_c(subStep) / real_c(numberOfSubSteps_);

      if ( fileIO_ )
      {
         WALBERLA_ROOT_SECTION()
         {
            std::ofstream file;
            file.open( fileName_.c_str(), std::ofstream::app );

            file << std::setprecision(10) << curTimestep << "\t" << weightForce_ << "\t" << hydForce << "\t" << lubForce << "\t" << collisionForce << "\n";
            file.close();
         }
      }

   }


private:

   std::string fileName_;
   bool fileIO_;

   real_t weightForce_;
   uint_t numberOfSubSteps_;
};

template<typename ParticleAccessor_T>
real_t getForce(walberla::id_t uid, ParticleAccessor_T & ac)
{
   auto idx = ac.uidToIdx(uid);
   real_t force = real_t(0);
   if(idx != ac.getInvalidIdx())
   {
      force = ac.getForce(idx)[2];
   }
   WALBERLA_MPI_SECTION()
   {
      mpi::allReduceInplace( force, mpi::SUM );
   }
   return force;
}

template< typename DensityInterpolator_T>
void logDensityToFile(const shared_ptr<StructuredBlockStorage> & blocks, BlockDataID densityInterpolatorID, real_t surfaceDistance,
473
                      Vector3<real_t> spherePosition, real_t diameter, const std::string & fileName, uint_t timestep, real_t averageDensityInSystem )
474
475
476
477
{
   real_t evalRadius = real_t(0.5) * diameter + surfaceDistance;
   std::vector< Vector3<real_t> > evaluationPositions;
   evaluationPositions.push_back(spherePosition + Vector3<real_t>(real_t(0), real_t(0), evalRadius));
478
   evaluationPositions.push_back(spherePosition + Vector3<real_t>(evalRadius / real_c(sqrt(real_t(2))), real_t(0), real_c(evalRadius / sqrt(real_t(2)))));
479
   evaluationPositions.push_back(spherePosition + Vector3<real_t>(evalRadius, real_t(0), real_t(0)));
480
   evaluationPositions.push_back(spherePosition + Vector3<real_t>(evalRadius / real_c(sqrt(real_t(2))), real_t(0), -evalRadius / real_c(sqrt(real_t(2)))));
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
   evaluationPositions.push_back(spherePosition + Vector3<real_t>(real_t(0), real_t(0), -evalRadius));

   std::vector<real_t> densityAtPos(evaluationPositions.size(), real_t(0));

   for( auto & block: *blocks)
   {
      auto densityInterpolator = block.getData<DensityInterpolator_T>(densityInterpolatorID);
      for(auto i = uint_t(0); i < evaluationPositions.size(); ++i)
      {
         auto pos = evaluationPositions[i];
         if(block.getAABB().contains(pos))
         {
            densityInterpolator->get(pos, &densityAtPos[i]);
         }
      }
   }

   for(auto i = uint_t(0); i < densityAtPos.size(); ++i)
   {
      // reduce to root
      mpi::reduceInplace(densityAtPos[i], mpi::SUM);
   }

   WALBERLA_ROOT_SECTION()
   {
      std::ofstream file;
      file.open( fileName.c_str(), std::ofstream::app );

      file << std::setprecision(8) << timestep << "\t" << averageDensityInSystem;
      for(auto density: densityAtPos)
      {
         file << std::setprecision(8) << "\t" << density;
      }
      file << "\n";
      file.close();
   }
}

template< typename BoundaryHandling_T>
real_t getAverageDensityInSystem(const shared_ptr<StructuredBlockStorage> & blocks, BlockDataID pdfFieldID, BlockDataID boundaryHandlingID)
{
   real_t totalMass = real_t(0);
   uint_t count = uint_t(0);
   for( auto & block: *blocks)
   {
      auto pdfField = block.getData<PdfField_T>(pdfFieldID);
      auto boundaryHandling = block.getData<BoundaryHandling_T>(boundaryHandlingID);
      WALBERLA_FOR_ALL_CELLS_XYZ(pdfField,
                                 if(boundaryHandling->isDomain(x,y,z))
                                 {
                                    totalMass += pdfField->getDensity(x,y,z);
                                    ++count;
                                 }
      );
   }

   // reduce to root
   mpi::reduceInplace(totalMass, mpi::SUM);
   mpi::reduceInplace(count, mpi::SUM);

   return totalMass/real_c(count);
}


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

//*******************************************************************************************************************
550
/*!\brief PHYSICAL test case of a sphere settling and colliding with a wall submerged in a viscous fluid.
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
 *
 * The trajectory of the bouncing sphere is logged and compared to reference experiments.
 *
 * for experiments see: Gondret, Lance, Petit - "Bouncing motion of spherical particles in fluids" 2002
 * for simulations see e.g.: Biegert, Vowinckel, Meiburg - "A collision model for grain-resolving simulations of flows over
 * dense, mobile, polydisperse granular sediment beds" 2017
 */
//*******************************************************************************************************************


int main( int argc, char **argv )
{
   Environment env( argc, argv );

   if( !env.config() )
   {
      WALBERLA_ABORT("Usage: " << argv[0] << " path-to-configuration-file \n");
   }

   Config::BlockHandle simulationInput = env.config()->getBlock( "BouncingSphere" );

   // setup description from separate block
   const std::string setup = simulationInput.getParameter<std::string>("setup");
   Config::BlockHandle setupDescription = env.config()->getBlock( setup );

   const real_t Re = setupDescription.getParameter<real_t>("Re");
   const real_t densityFluid_SI = setupDescription.getParameter<real_t>("densityFluid_SI");
   const real_t dynamicViscosityFluid_SI = setupDescription.getParameter<real_t>("dynamicViscosityFluid_SI");
   const real_t densitySphere_SI = setupDescription.getParameter<real_t>("densitySphere_SI");
   const real_t diameter_SI = setupDescription.getParameter<real_t>("diameter_SI");
   const real_t gravitationalAcceleration_SI = setupDescription.getParameter<real_t>("gravitationalAcceleration_SI");
   const real_t restitutionCoeff = setupDescription.getParameter<real_t>("restitutionCoeff");
   const uint_t numberOfBouncesUntilTermination = setupDescription.getParameter<uint_t>("numberOfBouncesUntilTermination");
   const Vector3<real_t> domainSizeNonDim = setupDescription.getParameter<Vector3<real_t> >("domainSize");
   const Vector3<uint_t> numberOfBlocksPerDirection = setupDescription.getParameter<Vector3<uint_t> >("numberOfBlocksPerDirection");

   // numerical parameters

   // control
   const bool randomizeInitialSpherePosition = simulationInput.getParameter<bool>("randomizeInitialSpherePosition");
   const bool initializeSphereVelocity = simulationInput.getParameter<bool>("initializeSphereVelocity");
   const bool applyOutflowBCAtTop = simulationInput.getParameter<bool>("applyOutflowBCAtTop");
   bool artificiallyAccelerateSphere = simulationInput.getParameter<bool>("artificiallyAccelerateSphere");

   // LBM
   const real_t uIn = simulationInput.getParameter<real_t>("uIn");
   const real_t diameter = simulationInput.getParameter<real_t>("diameter");
   const real_t magicNumber = simulationInput.getParameter<real_t>("magicNumber");
   const real_t bulkViscRateFactor = simulationInput.getParameter<real_t>("bulkViscRateFactor");
   const bool useOmegaBulkAdaption = simulationInput.getParameter<bool>("useOmegaBulkAdaption");

   const uint_t numRPDSubCycles = simulationInput.getParameter<uint_t>("numRPDSubCycles");
   const bool useLubricationCorrection = simulationInput.getParameter<bool>("useLubricationCorrection");
   const real_t lubricationCutOffDistance = simulationInput.getParameter<real_t>("lubricationCutOffDistance");
   const real_t lubricationMinimalGapSizeNonDim = simulationInput.getParameter<real_t>("lubricationMinimalGapSizeNonDim");
   const bool useVelocityVerlet = simulationInput.getParameter<bool>("useVelocityVerlet");

   // Collision Response
   const real_t collisionTime = simulationInput.getParameter<real_t>("collisionTime");
   const real_t StCrit = simulationInput.getParameter<real_t>("StCrit");
   const bool useACTM = simulationInput.getParameter<bool>("useACTM");

   const bool averageForceTorqueOverTwoTimeSteps = simulationInput.getParameter<bool>("averageForceTorqueOverTwoTimeSteps");
   const bool conserveMomentum = simulationInput.getParameter<bool>("conserveMomentum");
   const bool disableFluidForceDuringContact = simulationInput.getParameter<bool>("disableFluidForceDuringContact");
   const std::string reconstructorType = simulationInput.getParameter<std::string>("reconstructorType");

   const bool fileIO = simulationInput.getParameter<bool>("fileIO");
   const uint_t vtkIOFreq = simulationInput.getParameter<uint_t>("vtkIOFreq");
   const std::string baseFolder = simulationInput.getParameter<std::string>("baseFolder");
   bool vtkOutputOnCollision = simulationInput.getParameter<bool>("vtkOutputOnCollision");

   bool writeCheckPointFile = simulationInput.getParameter<bool>("writeCheckPointFile", false);
   bool readFromCheckPointFile = simulationInput.getParameter<bool>("readFromCheckPointFile", false);

   WALBERLA_ROOT_SECTION()
   {
      if( fileIO )
      {
         // create base directory if it does not yet exist
         filesystem::path tpath( baseFolder );
         if( !filesystem::exists( tpath ) )
            filesystem::create_directory( tpath );
      }
   }


   //////////////////////////////////////
   // SIMULATION PROPERTIES in SI units//
   //////////////////////////////////////

   const real_t densityRatio = densitySphere_SI / densityFluid_SI;
   const real_t kinematicViscosityFluid_SI = dynamicViscosityFluid_SI / densityFluid_SI;
   // here, we do the analogy with experiments via the Galileo number which is a-priori well-defined in contrast to the Stokes number which is a result of the simulation
   // -> no, this is not a good idea as the sphere in the experiments has in some cases not yet reached the terminal settling velocity
   // thus, it is better to pose the reported resulting Re as an input parameter
   const real_t uIn_SI = Re * kinematicViscosityFluid_SI / diameter_SI;
   const real_t ug_SI = std::sqrt((densityRatio-real_t(1)) * gravitationalAcceleration_SI * diameter_SI);
   const real_t Ga = ug_SI * diameter_SI / kinematicViscosityFluid_SI;

   //////////////////////////
   // NUMERICAL PARAMETERS //
   //////////////////////////

   const real_t dx_SI = diameter_SI / diameter;
   const real_t sphereVolume = math::pi / real_t(6) * diameter * diameter * diameter;
   const real_t dt_SI = uIn / uIn_SI * dx_SI;
   const real_t viscosity =  kinematicViscosityFluid_SI * dt_SI / ( dx_SI * dx_SI );
   const real_t gravitationalAcceleration = gravitationalAcceleration_SI * dt_SI * dt_SI / dx_SI;
   const real_t omega = lbm::collision_model::omegaFromViscosity(viscosity);

   const real_t Re_p = diameter * uIn / viscosity;
   const real_t St = densityRatio * uIn * diameter / ( real_t(9) * viscosity );

   const real_t omegaBulk = lbm_mesapd_coupling::omegaBulkFromOmega(omega, bulkViscRateFactor);

   const real_t densityFluid = real_t(1);
   const real_t densitySphere = densityRatio;

   const real_t dx = real_t(1);

   const real_t responseTime = densityRatio * diameter * diameter / ( real_t(18) * viscosity );
   const real_t accelerationFactor = real_t(1) / (real_t(0.1) * responseTime);
   const real_t tref = diameter / uIn;

   const real_t particleMass = densityRatio * sphereVolume;
   const real_t Mij = particleMass; // * particleMass / ( real_t(2) * particleMass ); // Mij = M for sphere-wall collision

   const real_t uInCrit = real_t(9) * StCrit * viscosity / ( densityRatio * diameter);

   Vector3<uint_t> domainSize( uint_c(domainSizeNonDim[0] * diameter ), uint_c(domainSizeNonDim[1] * diameter ), uint_c(domainSizeNonDim[2] * diameter ) );

   real_t initialSpherePosition = real_c(domainSize[2]) - real_t(1.5) * diameter;

   if( randomizeInitialSpherePosition )
   {
      uint_t seed1 = uint_c(std::chrono::system_clock::now().time_since_epoch().count());
      mpi::broadcastObject(seed1); // root process chooses seed and broadcasts it
689
      std::mt19937 g1 (static_cast<unsigned int>(seed1));
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

      real_t initialPositionOffset = real_t(0.5)*math::realRandom<real_t>(real_t(-1), real_t(1), g1);
      initialSpherePosition += initialPositionOffset;
   }


   WALBERLA_LOG_INFO_ON_ROOT("Setup (in SI units):");
   WALBERLA_LOG_INFO_ON_ROOT(" - sphere: diameter = " << diameter_SI << ", densityRatio = " << densityRatio);
   WALBERLA_LOG_INFO_ON_ROOT(" - fluid: density = " << densityFluid_SI << ", dyn. visc = " << dynamicViscosityFluid_SI << ", kin. visc = " << kinematicViscosityFluid_SI );
   WALBERLA_LOG_INFO_ON_ROOT(" - ug = " << ug_SI );
   WALBERLA_LOG_INFO_ON_ROOT(" - Galileo number = " << Ga );
   WALBERLA_LOG_INFO_ON_ROOT(" - target Reynolds number = " << Re );
   WALBERLA_LOG_INFO_ON_ROOT(" - domain size = <" << real_c(domainSizeNonDim[0]) * diameter_SI << "," << real_c(domainSizeNonDim[1]) * diameter_SI  << ","<< real_c(domainSizeNonDim[2]) * diameter_SI  << ">");
   WALBERLA_LOG_INFO_ON_ROOT(" - dx = " << dx_SI);
   WALBERLA_LOG_INFO_ON_ROOT(" - dt = " << dt_SI);

   WALBERLA_LOG_INFO_ON_ROOT("Setup (in simulation, i.e. lattice, units):");
   WALBERLA_LOG_INFO_ON_ROOT(" - domain size = " << domainSize);
   if(applyOutflowBCAtTop) WALBERLA_LOG_INFO_ON_ROOT(" - outflow BC at top");
   WALBERLA_LOG_INFO_ON_ROOT(" - sphere: diameter = " << diameter << ", density = " << densitySphere );
   WALBERLA_LOG_INFO_ON_ROOT(" - initial sphere position = " << initialSpherePosition );
   WALBERLA_LOG_INFO_ON_ROOT(" - fluid: density = " << densityFluid << ", relaxation time (tau) = " << real_t(1) / omega << ", omega = " << omega << ", kin. visc = " << viscosity );
   WALBERLA_LOG_INFO_ON_ROOT(" - magic number = " << magicNumber);
   WALBERLA_LOG_INFO_ON_ROOT(" - omega bulk = " << omegaBulk << ", bulk visc = " << lbm_mesapd_coupling::bulkViscosityFromOmegaBulk(omegaBulk) << " ( bulk visc rate factor (conste) = " << bulkViscRateFactor << ")");
   if(useOmegaBulkAdaption) WALBERLA_LOG_INFO_ON_ROOT(" - using omega bulk adaption in vicinity of sphere");
   WALBERLA_LOG_INFO_ON_ROOT(" - gravitational acceleration = " << gravitationalAcceleration );
   WALBERLA_LOG_INFO_ON_ROOT(" - expected settling velocity = " << uIn );
   WALBERLA_LOG_INFO_ON_ROOT(" - target Reynolds number = " << Re_p );
   WALBERLA_LOG_INFO_ON_ROOT(" - target Stokes number = " << St );
   WALBERLA_LOG_INFO_ON_ROOT(" - tref = " << tref );
   WALBERLA_LOG_INFO_ON_ROOT(" - Stokes response time = " << responseTime );
   if( artificiallyAccelerateSphere )
   {
      WALBERLA_LOG_INFO_ON_ROOT(" - artificially accelerating sphere with factor " << accelerationFactor)
   }
   WALBERLA_LOG_INFO_ON_ROOT(" - integrator = " << (useVelocityVerlet ? "Velocity Verlet" : "Explicit Euler") );
   if( vtkIOFreq > 0 )
   {
      WALBERLA_LOG_INFO_ON_ROOT(" - writing vtk files to folder \"" << baseFolder << "\" with frequency " << vtkIOFreq);
   }
   if(vtkOutputOnCollision)
   {
      WALBERLA_LOG_INFO_ON_ROOT(" - writing vtk files during collision events");
   }

   WALBERLA_LOG_INFO_ON_ROOT("Collision Response properties:" );
   WALBERLA_LOG_INFO_ON_ROOT(" - collision time = " << collisionTime );
   if(useACTM)
   {
      WALBERLA_LOG_INFO_ON_ROOT(" - using nonlinear collision model with ACTM" );
   }
   else
   {
743
      WALBERLA_LOG_INFO_ON_ROOT(" - using linear collision model with fixed parameters" );
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
   }
   WALBERLA_LOG_INFO_ON_ROOT(" - coeff of restitution = " << restitutionCoeff );
   WALBERLA_LOG_INFO_ON_ROOT(" - number of RPD sub cycles = " << numRPDSubCycles );
   WALBERLA_LOG_INFO_ON_ROOT(" - lubrication correction = " << (useLubricationCorrection ? "yes" : "no") );
   if( useLubricationCorrection )
   {
      WALBERLA_LOG_INFO_ON_ROOT("  - lubrication correction cut off = " << lubricationCutOffDistance );
      WALBERLA_LOG_INFO_ON_ROOT("  - lubrication correction minimal gap size non dim = " << lubricationMinimalGapSizeNonDim  << " ( = " << lubricationMinimalGapSizeNonDim * diameter * real_t(0.5) << " cells )");
   }
   WALBERLA_LOG_INFO_ON_ROOT("Coupling properties:" );
   WALBERLA_LOG_INFO_ON_ROOT(" - disable hydrodynamic forces during contact = " << (disableFluidForceDuringContact ? "yes" : "no") );
   if(disableFluidForceDuringContact) WALBERLA_LOG_INFO_ON_ROOT("  - StCrit = " << StCrit << ", uInCrit = " << uInCrit );
   WALBERLA_LOG_INFO_ON_ROOT(" - average forces over two LBM time steps = " << (averageForceTorqueOverTwoTimeSteps ? "yes" : "no") );
   WALBERLA_LOG_INFO_ON_ROOT(" - conserve momentum = " << (conserveMomentum ? "yes" : "no") );
   WALBERLA_LOG_INFO_ON_ROOT(" - reconstructor = " << reconstructorType );

   std::string checkPointFileName = "checkPointingFile_" + setup
                                    + "_uIn" + std::to_string(uIn)
                                    + "_d" + std::to_string(diameter)
                                    + "_bulkViscFactor" + std::to_string(bulkViscRateFactor)
                                    + "_mn" + std::to_string(magicNumber);
   if(useOmegaBulkAdaption) checkPointFileName +="_omegaBulkAdaption";
   if(applyOutflowBCAtTop) checkPointFileName +="_outflowBC";
767
768
769
770
771
772
#if defined(USE_TRTLIKE)
   checkPointFileName += "_trtlike";
#elif defined(USE_CUMULANT)
   checkPointFileName += "_cumulant";
#elif defined(USE_D3Q27TRTLIKE)
   checkPointFileName += "_d3q27trtlike";
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
#endif

   WALBERLA_LOG_INFO_ON_ROOT("Checkpointing:");
   WALBERLA_LOG_INFO_ON_ROOT(" - read from checkpoint file = " << (readFromCheckPointFile ? "yes" : "no"));
   WALBERLA_LOG_INFO_ON_ROOT(" - write checkpoint file = " << (writeCheckPointFile ? "yes" : "no"));
   if( readFromCheckPointFile || writeCheckPointFile) WALBERLA_LOG_INFO_ON_ROOT(" - checkPointingFileName = " << checkPointFileName );

   if(readFromCheckPointFile && writeCheckPointFile )
   {
      // decide which option to choose
      if(filesystem::exists( checkPointFileName+"_lbm.txt"))
      {
         WALBERLA_LOG_INFO_ON_ROOT("Checkpoint file already exists! Will skip writing check point file and start from this check point!");
         writeCheckPointFile = false;
      } else
      {
         WALBERLA_LOG_INFO_ON_ROOT("Checkpoint file does not exists yet! Will skip reading check point file and just regularly start the simulation from the beginning!");
         readFromCheckPointFile = false;
      }
   }

   if(readFromCheckPointFile && artificiallyAccelerateSphere)
   {
      // accelerating after reading from a checkpoint file does not make sense, as current actual runtime is not known
      WALBERLA_LOG_INFO_ON_ROOT("NOTE: switching off artificial acceleration of sphere due to reading from check point file");
      artificiallyAccelerateSphere = false;
   }


   ///////////////////////////
   // BLOCK STRUCTURE SETUP //
   ///////////////////////////


   Vector3<uint_t> cellsPerBlockPerDirection( domainSize[0] / numberOfBlocksPerDirection[0],
                                              domainSize[1] / numberOfBlocksPerDirection[1],
                                              domainSize[2] / numberOfBlocksPerDirection[2] );
   for( uint_t i = 0; i < 3; ++i ) {
      WALBERLA_CHECK_EQUAL(cellsPerBlockPerDirection[i] * numberOfBlocksPerDirection[i], domainSize[i],
                           "Unmatching domain decomposition in direction " << i << "!");
   }

815
816
817
818
819
820
821
822
823
   shared_ptr< StructuredBlockForest > blocks;
   if( readFromCheckPointFile )
   {
      WALBERLA_LOG_INFO_ON_ROOT("Reading block forest from file!");
      blocks = blockforest::createUniformBlockGrid(checkPointFileName+"_forest.txt", cellsPerBlockPerDirection[0], cellsPerBlockPerDirection[1], cellsPerBlockPerDirection[2], false);
   }
   else
   {
      blocks = blockforest::createUniformBlockGrid( numberOfBlocksPerDirection[0], numberOfBlocksPerDirection[1], numberOfBlocksPerDirection[2],
824
825
826
827
                                                      cellsPerBlockPerDirection[0], cellsPerBlockPerDirection[1], cellsPerBlockPerDirection[2], dx,
                                                      0, false, false,
                                                      true, true, false, //periodicity
                                                      false );
828
829
830
831
832
833
834
835
836
   
      if(writeCheckPointFile)
      {
         WALBERLA_LOG_INFO_ON_ROOT("Writing block forest to file!");
         blocks->getBlockForest().saveToFile(checkPointFileName+"_forest.txt");
      }
   }
   
   WALBERLA_LOG_INFO_ON_ROOT("Domain partitioning:");
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
   WALBERLA_LOG_INFO_ON_ROOT(" - blocks per direction = " << numberOfBlocksPerDirection );
   WALBERLA_LOG_INFO_ON_ROOT(" - cells per block = " << cellsPerBlockPerDirection );

   //write domain decomposition to file
   if( vtkIOFreq > 0 )
   {
      vtk::writeDomainDecomposition( blocks, "initial_domain_decomposition", baseFolder );
   }


   //////////////////
   // RPD COUPLING //
   //////////////////

   auto rpdDomain = std::make_shared<mesa_pd::domain::BlockForestDomain>(blocks->getBlockForestPointer());

   //init data structures
   auto ps = walberla::make_shared<mesa_pd::data::ParticleStorage>(1);
855
856
857
   auto ss = walberla::make_shared<mesa_pd::data::ShapeStorage>();
   using ParticleAccessor_T = mesa_pd::data::ParticleAccessorWithShape;
   auto accessor = walberla::make_shared<ParticleAccessor_T >(ps, ss);
858
859
860
861
862
   BlockDataID particleStorageID;
   if(readFromCheckPointFile)
   {
      WALBERLA_LOG_INFO_ON_ROOT( "Initializing particles from checkpointing file!" );
      particleStorageID = blocks->loadBlockData( checkPointFileName + "_mesa.txt", mesa_pd::domain::createBlockForestDataHandling(ps), "Particle Storage" );
863
864
      mesa_pd::mpi::ClearNextNeighborSync CNNS;
      CNNS(*accessor);
865
866
   } else {
      particleStorageID = blocks->addBlockData(mesa_pd::domain::createBlockForestDataHandling(ps), "Particle Storage");
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
   }

   // bounding planes
   createPlaneSetup(ps,ss,blocks->getDomain(), applyOutflowBCAtTop);

   // create sphere and store Uid
   Vector3<real_t> initialPosition( real_t(0.5) * real_c(domainSize[0]), real_t(0.5) * real_c(domainSize[1]), initialSpherePosition );
   auto sphereShape = ss->create<mesa_pd::data::Sphere>( diameter * real_t(0.5) );
   ss->shapes[sphereShape]->updateMassAndInertia(densitySphere);

   walberla::id_t sphereUid = 0;
   if( readFromCheckPointFile )
   {
      for( auto pIt = ps->begin(); pIt != ps->end(); ++pIt )
      {
         // find sphere in loaded data structure and store uid for later reference
         if( pIt->getShapeID() == sphereShape )
         {
            sphereUid = pIt->getUid();
         }
      }
   } else
   {
      // create sphere
      if (rpdDomain->isContainedInProcessSubdomain( uint_c(mpi::MPIManager::instance()->rank()), initialPosition ))
      {
         mesa_pd::data::Particle&& p = *ps->create();
         p.setPosition(initialPosition);
         p.setInteractionRadius(diameter * real_t(0.5));
         p.setOwner(mpi::MPIManager::instance()->rank());
         p.setShapeID(sphereShape);
         p.setType(0);
         if( initializeSphereVelocity ) p.setLinearVelocity(Vector3<real_t>(real_t(0), real_t(0), -real_t(0.1) * uIn));
         sphereUid = p.getUid();
      }
   }
   mpi::allReduceInplace(sphereUid, mpi::SUM);

   if(sphereUid == 0) WALBERLA_ABORT("No sphere present - aborting!"); //something went wrong in the checkpointing probably

   ///////////////////////
   // ADD DATA TO BLOCKS //
   ////////////////////////

   // add omega bulk field
912
   BlockDataID omegaBulkFieldID = field::addToStorage<ScalarField_T>( blocks, "omega bulk field", omegaBulk, field::fzyx);
913
914
915
916

   // create the lattice model
   real_t lambda_e = lbm::collision_model::TRT::lambda_e( omega );
   real_t lambda_d = lbm::collision_model::TRT::lambda_d( omega, magicNumber );
917
#ifdef WALBERLA_BUILD_WITH_CODEGEN
918
919
920
921
#ifdef USE_CUMULANT
   WALBERLA_LOG_INFO_ON_ROOT("Using generated cumulant lattice model!");
   LatticeModel_T latticeModel = LatticeModel_T(omega);
#elif defined(USE_TRTLIKE) || defined(USE_D3Q27TRTLIKE)
922
   WALBERLA_LOG_INFO_ON_ROOT("Using generated TRT-like lattice model!");
923
   LatticeModel_T latticeModel = LatticeModel_T(omegaBulkFieldID, lambda_d, lambda_e);
924
#endif
925
#else
926
927
   WALBERLA_LOG_INFO_ON_ROOT("Using waLBerla built-in MRT lattice model and ignoring omega bulk field since not supported!");
   LatticeModel_T latticeModel = LatticeModel_T(lbm::collision_model::D3Q19MRT( omegaBulk, omegaBulk, lambda_d, lambda_e, lambda_e, lambda_d ));
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
#endif

   // add PDF field
   BlockDataID pdfFieldID;
   if( readFromCheckPointFile )
   {
      // add PDF field
      WALBERLA_LOG_INFO_ON_ROOT( "Initializing PDF Field from checkpointing file!" );
      shared_ptr< lbm::internal::PdfFieldHandling< LatticeModel_T > > dataHandling =
            make_shared< lbm::internal::PdfFieldHandling< LatticeModel_T > >( blocks, latticeModel, false,
                                                                              Vector3< real_t >( real_t(0) ), real_t(1),
                                                                              uint_t(1), field::fzyx );

      pdfFieldID = blocks->loadBlockData( checkPointFileName+"_lbm.txt", dataHandling, "pdf field" );

   } else {
      // add PDF field
      pdfFieldID = lbm::addPdfFieldToStorage< LatticeModel_T >( blocks, "pdf field", latticeModel,
                                                                Vector3< real_t >( real_t(0) ), real_t(1),
                                                                uint_t(1), field::fzyx );
   }


   // add flag field
   BlockDataID flagFieldID = field::addFlagFieldToStorage<FlagField_T>( blocks, "flag field" );

   // add particle field
   BlockDataID particleFieldID = field::addToStorage<lbm_mesapd_coupling::ParticleField_T>( blocks, "particle field", accessor->getInvalidUid(), field::fzyx, FieldGhostLayers );

   // add boundary handling
   using BoundaryHandling_T = MyBoundaryHandling<ParticleAccessor_T>::Type;
   BlockDataID boundaryHandlingID = blocks->addStructuredBlockData< BoundaryHandling_T >(MyBoundaryHandling<ParticleAccessor_T>( flagFieldID, pdfFieldID, particleFieldID, accessor, applyOutflowBCAtTop), "boundary handling" );

   // interpolation functionality
   using DensityAdaptor_T = typename lbm::Adaptor< LatticeModel_T >::Density;
   BlockDataID densityAdaptorID = field::addFieldAdaptor<  DensityAdaptor_T >( blocks, pdfFieldID, "density adaptor" );

   using DensityInterpolator_T = typename field::NearestNeighborFieldInterpolator<DensityAdaptor_T, FlagField_T>;
   BlockDataID densityInterpolatorID = field::addFieldInterpolator< DensityInterpolator_T, FlagField_T >( blocks, densityAdaptorID, flagFieldID, Fluid_Flag );

   // set up RPD functionality
   std::function<void(void)> syncCall = [&ps,&rpdDomain](){
      const real_t overlap = real_t( 1.5 );
      mesa_pd::mpi::SyncNextNeighbors syncNextNeighborFunc;
      syncNextNeighborFunc(*ps, *rpdDomain, overlap);
   };

   syncCall();

   real_t timeStepSizeRPD = real_t(1)/real_t(numRPDSubCycles);
978
979
980
   mesa_pd::kernel::ExplicitEuler explEulerIntegrator(timeStepSizeRPD);
   mesa_pd::kernel::VelocityVerletPreForceUpdate  vvIntegratorPreForce(timeStepSizeRPD);
   mesa_pd::kernel::VelocityVerletPostForceUpdate vvIntegratorPostForce(timeStepSizeRPD);
981
982
983

   // linear model
   mesa_pd::kernel::LinearSpringDashpot linearCollisionResponse(1);
984
   linearCollisionResponse.setStiffnessAndDamping(0,0,restitutionCoeff,collisionTime,real_t(0),Mij); // no response in tangential direction
985
986
987

   // nonlinear model for ACTM
   mesa_pd::kernel::NonLinearSpringDashpot nonLinearCollisionResponse(1, collisionTime);
988
   nonLinearCollisionResponse.setLnCORsqr(0,0,std::log(restitutionCoeff) * std::log(restitutionCoeff));
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
   nonLinearCollisionResponse.setMeff(0,0,Mij);

   mesa_pd::mpi::ReduceProperty reduceProperty;

   mesa_pd::mpi::ReduceContactHistory reduceAndSwapContactHistory;

   // set up coupling functionality
   Vector3<real_t> gravitationalForce( real_t(0), real_t(0), -(densitySphere - densityFluid) * gravitationalAcceleration * sphereVolume );
   lbm_mesapd_coupling::AddForceOnParticlesKernel addGravitationalForce(gravitationalForce);
   lbm_mesapd_coupling::ResetHydrodynamicForceTorqueKernel resetHydrodynamicForceTorque;
   lbm_mesapd_coupling::AverageHydrodynamicForceTorqueKernel averageHydrodynamicForceTorque;
   lbm_mesapd_coupling::LubricationCorrectionKernel lubricationCorrectionKernel(viscosity, [lubricationMinimalGapSizeNonDim](real_t r){ return lubricationMinimalGapSizeNonDim * r;}, lubricationCutOffDistance, real_t(0), real_t(0) ); // no tangential components needed
   lbm_mesapd_coupling::ParticleMappingKernel<BoundaryHandling_T> particleMappingKernel(blocks, boundaryHandlingID);
   lbm_mesapd_coupling::MovingParticleMappingKernel<BoundaryHandling_T> movingParticleMappingKernel(blocks, boundaryHandlingID, particleFieldID);

   ///////////////
   // TIME LOOP //
   ///////////////

   // map planes into the LBM simulation -> act as no-slip boundaries
   ps->forEachParticle(false, lbm_mesapd_coupling::GlobalParticlesSelector(), *accessor, particleMappingKernel, *accessor, NoSlip_Flag);

   // map particles into the LBM simulation
   ps->forEachParticle(false, lbm_mesapd_coupling::RegularParticlesSelector(), *accessor, movingParticleMappingKernel, *accessor, MO_Flag);

   // setup of the LBM communication for synchronizing the pdf field between neighboring blocks
   blockforest::communication::UniformBufferedScheme< Stencil_T > optimizedPDFCommunicationScheme( blocks );
   optimizedPDFCommunicationScheme.addPackInfo( make_shared< lbm::PdfFieldPackInfo< LatticeModel_T > >( pdfFieldID ) ); // optimized sync

   blockforest::communication::UniformBufferedScheme< Stencil_T > fullPDFCommunicationScheme( blocks );
   fullPDFCommunicationScheme.addPackInfo( make_shared< field::communication::PackInfo< PdfField_T > >( pdfFieldID ) ); // full sync

   // create the timeloop
   const uint_t timesteps = uint_c( real_t(3) * initialPosition[2] / uIn);

   SweepTimeloop timeloop( blocks->getBlockStorage(), timesteps );

   timeloop.addFuncBeforeTimeStep( RemainingTimeLogger( timeloop.getNrOfTimeSteps() ), "Remaining Time Logger" );

   // vtk output
   if( vtkIOFreq != uint_t(0) )
   {
      // sphere
      auto particleVtkOutput = make_shared<mesa_pd::vtk::ParticleVtkOutput>(ps);
      particleVtkOutput->addOutput<mesa_pd::data::SelectParticleOwner>("owner");
      particleVtkOutput->addOutput<mesa_pd::data::SelectParticleUid>("uid");
      particleVtkOutput->addOutput<mesa_pd::data::SelectParticleLinearVelocity>("velocity");
      particleVtkOutput->setParticleSelector( [sphereShape](const mesa_pd::data::ParticleStorage::iterator& pIt) {return pIt->getShapeID() == sphereShape;} ); //limit output to sphere
      auto particleVtkWriter = vtk::createVTKOutput_PointData(particleVtkOutput, "Particles", vtkIOFreq, baseFolder, "simulation_step");
      timeloop.addFuncBeforeTimeStep( vtk::writeFiles( particleVtkWriter ), "VTK (sphere data)" );

      // flag field (written only once in the first time step, ghost layers are also written)

      auto flagFieldVTK = vtk::createVTKOutput_BlockData( blocks, "flag_field", timesteps, FieldGhostLayers, false, baseFolder );
      flagFieldVTK->addCellDataWriter( make_shared< field::VTKWriter< FlagField_T > >( flagFieldID, "FlagField" ) );
      timeloop.addFuncBeforeTimeStep( vtk::writeFiles( flagFieldVTK ), "VTK (flag field data)" );


      // pdf field, as slice
      auto pdfFieldVTK = vtk::createVTKOutput_BlockData( blocks, "fluid_field", vtkIOFreq, 0, false, baseFolder );

      pdfFieldVTK->addBeforeFunction( fullPDFCommunicationScheme );

      AABB sliceAABB( real_t(0), real_c(domainSize[1])*real_t(0.5)-real_t(1), real_t(0),
                      real_c(domainSize[0]), real_c(domainSize[1])*real_t(0.5)+real_t(1), real_c(domainSize[2]) );
      vtk::AABBCellFilter aabbSliceFilter( sliceAABB );

      field::FlagFieldCellFilter< FlagField_T > fluidFilter( flagFieldID );
      fluidFilter.addFlag( Fluid_Flag );

      vtk::ChainedFilter combinedSliceFilter;
      combinedSliceFilter.addFilter( fluidFilter );
      combinedSliceFilter.addFilter( aabbSliceFilter );

      pdfFieldVTK->addCellInclusionFilter( combinedSliceFilter );

      pdfFieldVTK->addCellDataWriter( make_shared< lbm::VelocityVTKWriter< LatticeModel_T, float > >( pdfFieldID, "VelocityFromPDF" ) );
      pdfFieldVTK->addCellDataWriter( make_shared< lbm::DensityVTKWriter < LatticeModel_T, float > >( pdfFieldID, "DensityFromPDF" ) );

      timeloop.addFuncBeforeTimeStep( vtk::writeFiles( pdfFieldVTK ), "VTK (fluid field data)" );

      // omega bulk field
      //timeloop.addFuncBeforeTimeStep( field::createVTKOutput<ScalarField_T, float>( omegaBulkFieldID, *blocks, "omega_bulk_field", vtkIOFreq, uint_t(0), false, baseFolder ), "VTK (omega bulk field)" );
   }



   // update bulk omega in all cells to adapt to changed particle position
   if(useOmegaBulkAdaption)
   {
      using OmegaBulkAdapter_T = lbm_mesapd_coupling::OmegaBulkAdapter<ParticleAccessor_T, lbm_mesapd_coupling::RegularParticlesSelector>;
      real_t defaultOmegaBulk = lbm_mesapd_coupling::omegaBulkFromOmega(omega, real_t(1));
      real_t adaptionLayerSize = real_t(2);
      shared_ptr<OmegaBulkAdapter_T> omegaBulkAdapter = make_shared<OmegaBulkAdapter_T>(blocks, omegaBulkFieldID, accessor, defaultOmegaBulk, omegaBulk, adaptionLayerSize, lbm_mesapd_coupling::RegularParticlesSelector());

      // initial adaption
      for (auto blockIt = blocks->begin(); blockIt != blocks->end(); ++blockIt) {
         (*omegaBulkAdapter)(blockIt.get());
      }

      timeloop.add() << Sweep( makeSharedSweep(omegaBulkAdapter), "Omega Bulk Adapter");
   }


   // add LBM communication function and boundary handling sweep (does the hydro force calculations and the no-slip treatment)
   timeloop.add() << BeforeFunction( optimizedPDFCommunicationScheme, "LBM Communication" )
                  << Sweep( BoundaryHandling_T::getBlockSweep( boundaryHandlingID ), "Boundary Handling" );

   // stream + collide LBM step
1098
1099
1100
1101
1102
1103
1104
#ifdef WALBERLA_BUILD_WITH_CODEGEN
   auto lbmSweep = LatticeModel_T::Sweep( pdfFieldID );
   timeloop.add() << Sweep( lbmSweep, "LB sweep" );
#else
   auto lbmSweep = lbm::makeCellwiseSweep< LatticeModel_T, FlagField_T >( pdfFieldID, flagFieldID, Fluid_Flag );
   timeloop.add() << Sweep( makeSharedSweep( lbmSweep ), "cell-wise LB sweep" );
#endif
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167

   // this is carried out after the particle integration, it corrects the flag field and restores missing PDF information
   // then, the checkpointing file can be written, as otherwise some cells are invalid and can not be recovered
   SweepTimeloop timeloopAfterParticles( blocks->getBlockStorage(), timesteps );

   // sweep for updating the particle mapping into the LBM simulation
   timeloopAfterParticles.add() << Sweep( lbm_mesapd_coupling::makeMovingParticleMapping<PdfField_T,BoundaryHandling_T>(blocks, pdfFieldID,boundaryHandlingID, particleFieldID, accessor, MO_Flag, FormerMO_Flag, lbm_mesapd_coupling::RegularParticlesSelector(), conserveMomentum), "Particle Mapping" );

   // sweep for restoring PDFs in cells previously occupied by particles
   if( reconstructorType == "EAN" )
   {
      auto sphereNormalExtrapolationDirectionFinder = make_shared<lbm_mesapd_coupling::SphereNormalExtrapolationDirectionFinder>(blocks);
      auto equilibriumAndNonEquilibriumSphereNormalReconstructor = lbm_mesapd_coupling::makeEquilibriumAndNonEquilibriumReconstructor<BoundaryHandling_T>(blocks, boundaryHandlingID, sphereNormalExtrapolationDirectionFinder, uint_t(3), true);
      timeloopAfterParticles.add() << BeforeFunction( fullPDFCommunicationScheme, "LBM Communication" )
                                   << Sweep( makeSharedSweep(lbm_mesapd_coupling::makePdfReconstructionManager<PdfField_T,BoundaryHandling_T>(blocks, pdfFieldID, boundaryHandlingID, particleFieldID, accessor, FormerMO_Flag, Fluid_Flag, equilibriumAndNonEquilibriumSphereNormalReconstructor, conserveMomentum)), "PDF Restore" );
   }else if( reconstructorType == "Grad" )
   {
      bool recomputeTargetDensity = false;
      auto gradReconstructor = lbm_mesapd_coupling::makeGradsMomentApproximationReconstructor<BoundaryHandling_T>(blocks, boundaryHandlingID, omega, recomputeTargetDensity,true);

      timeloopAfterParticles.add() << BeforeFunction( fullPDFCommunicationScheme, "PDF Communication" )
                                   << Sweep( makeSharedSweep(lbm_mesapd_coupling::makePdfReconstructionManager<PdfField_T,BoundaryHandling_T>(blocks, pdfFieldID, boundaryHandlingID, particleFieldID, accessor, FormerMO_Flag, Fluid_Flag, gradReconstructor, conserveMomentum) ), "PDF Restore" );
   }else if( reconstructorType == "Eq" )
   {
      timeloopAfterParticles.add() << Sweep( makeSharedSweep(lbm_mesapd_coupling::makePdfReconstructionManager<PdfField_T,BoundaryHandling_T>(blocks, pdfFieldID, boundaryHandlingID, particleFieldID, accessor, FormerMO_Flag, Fluid_Flag, conserveMomentum)), "PDF Restore" );
   }else if( reconstructorType == "Ext")
   {
      auto sphereNormalExtrapolationDirectionFinder = make_shared<lbm_mesapd_coupling::SphereNormalExtrapolationDirectionFinder>(blocks);
      auto extrapolationReconstructor = lbm_mesapd_coupling::makeExtrapolationReconstructor<BoundaryHandling_T, lbm_mesapd_coupling::SphereNormalExtrapolationDirectionFinder, false>(blocks, boundaryHandlingID, sphereNormalExtrapolationDirectionFinder, uint_t(3), true);
      timeloopAfterParticles.add() << BeforeFunction( fullPDFCommunicationScheme, "LBM Communication" )
                                   << Sweep( makeSharedSweep(lbm_mesapd_coupling::makePdfReconstructionManager<PdfField_T,BoundaryHandling_T>(blocks, pdfFieldID, boundaryHandlingID, particleFieldID, accessor, FormerMO_Flag, Fluid_Flag, extrapolationReconstructor, conserveMomentum)), "PDF Restore" );
   } else
   {
      WALBERLA_ABORT("Reconstructor Type " << reconstructorType << " not implemented!");
   }


   // evaluation functionality
   std::string loggingFileName( baseFolder + "/LoggingSphereWallCollision.txt");
   if( fileIO  )
   {
      WALBERLA_LOG_INFO_ON_ROOT(" - writing logging output to file \"" << loggingFileName << "\"");
   }
   SpherePropertyLogger<ParticleAccessor_T> logger( accessor, sphereUid, loggingFileName, fileIO, diameter, uIn, numRPDSubCycles );

   std::string forceLoggingFileName( baseFolder + "/ForceLoggingSphereWallCollision.txt");
   if( fileIO  )
   {
      WALBERLA_LOG_INFO_ON_ROOT(" - writing force logging output to file \"" << forceLoggingFileName << "\"");
   }
   ForcesOnSphereLogger sphereForceLogger(forceLoggingFileName, fileIO, gravitationalForce[2], numRPDSubCycles);

   ////////////////////////
   // EXECUTE SIMULATION //
   ////////////////////////

   WcTimingPool timeloopTiming;

   // evaluation quantities
   uint_t numBounces = uint_t(0);
   uint_t tImpact = uint_t(0);
   std::vector<uint_t> impactTimes;

1168
1169
1170
   real_t curVel(real_t(0));
   real_t oldVel(real_t(0));
   real_t maxSettlingVel(real_t(0));
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
   std::vector<real_t> maxSettlingVelBetweenBounces;

   real_t minGapSize(real_t(0));
   real_t maxGapSize(real_t(0));
   std::vector<real_t> maxGapSizeBetweenBounces;

   real_t actualSt(real_t(0));
   real_t actualRe(real_t(0));

   WALBERLA_LOG_INFO_ON_ROOT("Running for maximum of " << timesteps << " timesteps!");
   
   const bool useOpenMP = false;

   // special vtk output
   auto pdfFieldVTKCollision = vtk::createVTKOutput_BlockData( blocks, "collision_fluid_field", 1, 1, false, baseFolder );

   blockforest::communication::UniformBufferedScheme< stencil::D3Q27 > pdfGhostLayerSync( blocks );
   pdfGhostLayerSync.addPackInfo( make_shared< field::communication::PackInfo< PdfField_T > >( pdfFieldID ) );
   pdfFieldVTKCollision->addBeforeFunction( pdfGhostLayerSync );

   AABB boxAABB( initialPosition[0] - diameter, initialPosition[1] - diameter, -real_t(1),
                 initialPosition[0] + diameter, initialPosition[1] + diameter, real_t(2) * diameter);
   vtk::AABBCellFilter aabbBoxFilter( boxAABB );

   field::FlagFieldCellFilter< FlagField_T > fluidFilter( flagFieldID );
   fluidFilter.addFlag( Fluid_Flag );

   vtk::ChainedFilter combinedSliceFilter;
   combinedSliceFilter.addFilter( fluidFilter );
   combinedSliceFilter.addFilter( aabbBoxFilter );

   pdfFieldVTKCollision->addCellInclusionFilter( combinedSliceFilter );

   pdfFieldVTKCollision->addCellDataWriter( make_shared< lbm::VelocityVTKWriter< LatticeModel_T, float > >( pdfFieldID, "VelocityFromPDF" ) );
   pdfFieldVTKCollision->addCellDataWriter( make_shared< lbm::DensityVTKWriter < LatticeModel_T, float > >( pdfFieldID, "DensityFromPDF" ) );


   // time loop
   for (uint_t i = 0; i < timesteps; ++i )
   {
      // perform a single simulation step -> this contains LBM and setting of the hydrodynamic interactions
      timeloop.singleStep( timeloopTiming );
1213
1214
         
      reduceProperty.operator()<mesa_pd::HydrodynamicForceTorqueNotification>(*ps);
1215

1216
      if( averageForceTorqueOverTwoTimeSteps )
1217
      {
1218
1219
1220
         if( i == 0 )
         {
            lbm_mesapd_coupling::InitializeHydrodynamicForceTorqueForAveragingKernel initializeHydrodynamicForceTorqueForAveragingKernel;
1221
            ps->forEachParticle(useOpenMP, mesa_pd::kernel::SelectLocal(), *accessor, initializeHydrodynamicForceTorqueForAveragingKernel, *accessor );
1222
         }
1223
         ps->forEachParticle(useOpenMP, mesa_pd::kernel::SelectLocal(), *accessor, averageHydrodynamicForceTorque, *accessor );
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
      }

      real_t hydForce(real_t(0));
      real_t lubForce(real_t(0));
      real_t collisionForce(real_t(0));

      for(auto subCycle = uint_t(0); subCycle < numRPDSubCycles; ++subCycle )
      {

         timeloopTiming["RPD"].start();

         if( useVelocityVerlet )
         {

            Vector3<real_t> oldForce;

            if(artificiallyAccelerateSphere)
            {
               // since the pre-force step of VV updates the position based on velocity and old force, we set oldForce to zero here for this step while artificially accelerating
               // to not perturb the step after which artificial acceleration is switched off (which requires valid oldForce values then) we store the old force and then re-apply it
               auto idx = accessor->uidToIdx(sphereUid);
               if(idx != accessor->getInvalidIdx())
               {
                  oldForce = accessor->getOldForce(idx);
                  accessor->setOldForce(idx, Vector3<real_t>(real_t(0)));
               }
            }

            ps->forEachParticle(useOpenMP, mesa_pd::kernel::SelectLocal(), *accessor, vvIntegratorPreForce, *accessor);
            syncCall();

            if(artificiallyAccelerateSphere)
            {
               // re-apply old force
               auto idx = accessor->uidToIdx(sphereUid);
               if(idx != accessor->getInvalidIdx())
               {
                  accessor->setOldForce(idx, oldForce);
               }
            }
         }

         // lubrication correction
         ps->forEachParticlePairHalf(useOpenMP, mesa_pd::kernel::ExcludeInfiniteInfinite(), *accessor,
                                     [&lubricationCorrectionKernel,&rpdDomain](const size_t idx1, const size_t idx2, auto& ac)
                                     {
                                        mesa_pd::collision_detection::AnalyticContactDetection acd;
                                        acd.getContactThreshold() = lubricationCorrectionKernel.getNormalCutOffDistance();
                                        mesa_pd::kernel::DoubleCast double_cast;
                                        mesa_pd::mpi::ContactFilter contact_filter;
                                        if (double_cast(idx1, idx2, ac, acd, ac ))
                                        {
                                           if (contact_filter(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), *rpdDomain))
                                           {
                                              double_cast(acd.getIdx1(), acd.getIdx2(), ac, lubricationCorrectionKernel, ac, acd.getContactNormal(), acd.getPenetrationDepth());
                                           }
                                        }
                                     },
                                     *accessor );

         lubForce = getForce(sphereUid,*accessor);

         // one could add linked cells here

         // collision response
         ps->forEachParticlePairHalf(useOpenMP, mesa_pd::kernel::ExcludeInfiniteInfinite(), *accessor,
                                     [&linearCollisionResponse, &nonLinearCollisionResponse, &rpdDomain, timeStepSizeRPD, useACTM]
                                     (const size_t idx1, const size_t idx2, auto& ac)
                                     {
                                        mesa_pd::collision_detection::AnalyticContactDetection acd;
                                        mesa_pd::kernel::DoubleCast double_cast;
                                        mesa_pd::mpi::ContactFilter contact_filter;
                                        if (double_cast(idx1, idx2, ac, acd, ac ))
                                        {
                                           if (contact_filter(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), *rpdDomain))
                                           {
                                              if(useACTM) nonLinearCollisionResponse(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), acd.getContactNormal(), acd.getPenetrationDepth(), timeStepSizeRPD);
                                              else        linearCollisionResponse(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), acd.getContactNormal(), acd.getPenetrationDepth(), timeStepSizeRPD);
                                           }
                                        }
                                     },
                                     *accessor );

         collisionForce = getForce(sphereUid,*accessor) - lubForce;

         reduceAndSwapContactHistory(*ps);

         // add hydrodynamic force
         if( disableFluidForceDuringContact )
         {
            lbm_mesapd_coupling::StokesNumberBasedSphereSelector sphereSelector(StCrit, densityFluid, densitySphere, viscosity);
            lbm_mesapd_coupling::AddHydrodynamicInteractionKernel addHydrodynamicInteraction;
            ps->forEachParticle(useOpenMP, sphereSelector, *accessor, addHydrodynamicInteraction, *accessor );
         }
         else
         {
            lbm_mesapd_coupling::AddHydrodynamicInteractionKernel addHydrodynamicInteraction;
1321
            ps->forEachParticle(useOpenMP, mesa_pd::kernel::SelectLocal(), *accessor, addHydrodynamicInteraction, *accessor );
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
         }

         hydForce = getForce(sphereUid,*accessor) - lubForce - collisionForce;

         ps->forEachParticle( useOpenMP, mesa_pd::kernel::SelectLocal(), *accessor, addGravitationalForce, *accessor );

         reduceProperty.operator()<mesa_pd::ForceTorqueNotification>(*ps);

         // integration
         if( useVelocityVerlet ) ps->forEachParticle( useOpenMP, mesa_pd::kernel::SelectLocal(), *accessor, vvIntegratorPostForce, *accessor);
         else ps->forEachParticle( useOpenMP, mesa_pd::kernel::SelectLocal(), *accessor, explEulerIntegrator, *accessor);

         syncCall();

         if( artificiallyAccelerateSphere)
         {
            // overwrite velocity of particle with prescribed one
            lbm_mesapd_coupling::RegularParticlesSelector sphereSelector;
1340
            real_t newSphereVel = uIn * (std::exp(-accelerationFactor * real_t(i) ) - real_t(1));
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
            ps->forEachParticle(useOpenMP, sphereSelector, *accessor, [newSphereVel](const size_t idx, ParticleAccessor_T& ac){ ac.setLinearVelocity(idx, Vector3<real_t>(real_t(0), real_t(0), newSphereVel));}, *accessor);
         }

         timeloopTiming["RPD"].end();

         // logging
         timeloopTiming["Logging"].start();
         logger(i, subCycle);
         sphereForceLogger(i, subCycle, hydForce, lubForce, collisionForce);
         timeloopTiming["Logging"].end();

      }

      ps->forEachParticle( useOpenMP, mesa_pd::kernel::SelectAll(), *accessor, resetHydrodynamicForceTorque, *accessor );


      // update particle mapping
      timeloopAfterParticles.singleStep(timeloopTiming);

      // check for termination
      oldVel = curVel;
      curVel = logger.getSettlingVelocity();

      maxSettlingVel = std::min(maxSettlingVel, curVel);
      minGapSize = std::min(minGapSize, logger.getGapSize());

      if(numBounces >= uint_t(1))
      {
         maxGapSize = std::max(maxGapSize, logger.getGapSize());
      }

      // detect a bounce
      if( oldVel < real_t(0) && curVel > real_t(0) && logger.getGapSize() < real_t(1))
      {

         ++numBounces;
         WALBERLA_LOG_INFO_ON_ROOT("Detected bounce number " << numBounces);
         WALBERLA_LOG_INFO_ON_ROOT("Max settling velocity " << maxSettlingVel << " -> St = "  << densityRatio * std::abs(maxSettlingVel) * diameter / ( real_t(9) * viscosity ) );

         if(numBounces == uint_t(1))
         {
            actualSt = densityRatio * std::abs(maxSettlingVel) * diameter / ( real_t(9) * viscosity );
            actualRe = std::abs(maxSettlingVel) * diameter / viscosity;
         }

         // reset and store quantities
         maxSettlingVelBetweenBounces.push_back(maxSettlingVel);
         maxSettlingVel = real_t(0);

         if(numBounces > uint_t(1))
         {
            maxGapSizeBetweenBounces.push_back(maxGapSize);
            maxGapSize = real_t(0);
         }

         if( numBounces >= numberOfBouncesUntilTermination )
         {
            break;
         }
      }

      // impact times are measured when the contact between sphere and wall is broken up again
      if( tImpact == uint_t(0) && numBounces == uint_t(1) && logger.getGapSize() > real_t(0) )
      {
         tImpact = i;
         WALBERLA_LOG_INFO_ON_ROOT("Detected impact time at time step " << tImpact);

         // switch off special vtk output after first bounce
         if(vtkOutputOnCollision) vtkOutputOnCollision = false;
      }

      if( numBounces > impactTimes.size() && logger.getGapSize() > real_t(0))
      {
         impactTimes.push_back(i);
      }

      // evaluate density around sphere
      if(fileIO)
      {
         std::string densityLoggingFileName(baseFolder + "/DensityLogging.txt");
         if(i == uint_t(0))
         {
            WALBERLA_ROOT_SECTION()
            {
               std::ofstream file;
               file.open( densityLoggingFileName.c_str());
               file << "# t \t avgDensity \t densitiesAtDifferentPositions \n";
               file.close();
            }
         }

         real_t surfaceDistance = real_t(2);
         real_t avgDensity = getAverageDensityInSystem<BoundaryHandling_T>(blocks, pdfFieldID, boundaryHandlingID);
         logDensityToFile<DensityInterpolator_T>(blocks, densityInterpolatorID, surfaceDistance, logger.getPosition(), diameter, densityLoggingFileName, i, avgDensity);
      }


      // check if sphere is close to bottom plane
      if( logger.getGapSize() < real_t(1) * diameter ) {

         // write a single checkpointing file before any collision relevant parameters take effect
         if( writeCheckPointFile )
         {
            WALBERLA_LOG_INFO_ON_ROOT("Writing checkpointing file in time step " << timeloop.getCurrentTimeStep());

            blocks->saveBlockData( checkPointFileName + "_lbm.txt", pdfFieldID );
            blocks->saveBlockData( checkPointFileName + "_mesa.txt", particleStorageID );

            writeCheckPointFile = false;
         }

         // switch off acceleration after writing of check point file
         artificiallyAccelerateSphere = false;

      }

      if(vtkOutputOnCollision && logger.getGapSize() < real_t(0))
      {
         pdfFieldVTKCollision->write();
      }

   }

   WALBERLA_LOG_INFO_ON_ROOT("Detected " << numBounces << " bounces, terminating simulation.");
   WALBERLA_LOG_INFO_ON_ROOT("Maximum settling velocities: " );
   for( auto vel : maxSettlingVelBetweenBounces)
   {
      WALBERLA_LOG_INFO_ON_ROOT(" - vel = " << vel << " -> St = "  << densityRatio * std::abs(vel) * diameter / ( real_t(9) * viscosity ) << ", Re = " << std::abs(vel) * diameter / viscosity );
   }
   WALBERLA_LOG_INFO_ON_ROOT("Maximum gap sizes between bounces:")
   for( auto gapSize : maxGapSizeBetweenBounces)
   {
      WALBERLA_LOG_INFO_ON_ROOT(" - gap size = " << gapSize << ", gap size / diameter = "  << gapSize / diameter );
   }


   std::string resultFile(baseFolder + "/ResultBouncingSphere.txt");
   WALBERLA_LOG_INFO_ON_ROOT("Writing logging file " << resultFile);
   logger.writeResult(resultFile, tImpact);

   std::string summaryFile(baseFolder + "/Summary.txt");
   WALBERLA_LOG_INFO_ON_ROOT("Writing summary file " << summaryFile);
   WALBERLA_ROOT_SECTION()
   {
      std::ofstream file;
      file.open( summaryFile.c_str());

      file << "waLBerla Revision = " << std::string(WALBERLA_GIT_SHA1).substr(0,8) << "\n";
      file << "\nInput parameters:\n";
      file << "Ga = " << Ga << "\n";
      file << "uIn = " << uIn << "\n";
      file << "LBM parameters:\n";
      file << " - magic number = " << magicNumber << "\n";
      file << " - bulk viscosity rate factor = " << bulkViscRateFactor << "\n";
      file << " - use omega bulk adaption = " << useOmegaBulkAdaption << " (layer size = 2)\n";
      file << "Collision parameters:\n";
      file << " - subCycles = " << numRPDSubCycles << "\n";
      file << " - collision time (Tc) = " << collisionTime << "\n";
      file << " - use ACTM = " << useACTM << "\n";
      file << "use lubrication correction = " << useLubricationCorrection << "\n";
      file << " - minimum gap size non dim = " << lubricationMinimalGapSizeNonDim << "\n";
      file << " - cut off distance = " << lubricationCutOffDistance << "\n";
      file << "switch off fluid interaction force during contact = " << disableFluidForceDuringContact << "\n";
      file << " - St_Crit = " << StCrit << "\n";
      file << "reconstructor type " << reconstructorType << "\n";
      file << "apply outflow BC at top = " << applyOutflowBCAtTop << "\n";
      file << "started from checkpoint file = " << readFromCheckPointFile << " ( " << checkPointFileName << ")\n";

      file << "\nOutput quantities:\n";
      file << "actual St = " << actualSt << "\n";
      file << "actual Re = " << actualRe << "\n";
      file << "impact times:\n";
      uint_t bounce = uint_t(0);
      for( auto impact : impactTimes)
      {
         file << " - " << bounce++ << ": impact time = " << impact << "\n";
      }
      file << "settling velocities:\n";
      bounce = uint_t(0);
      for( auto vel : maxSettlingVelBetweenBounces)
      {
        file << " - " << bounce++ << ": vel = " << vel << " -> St = "  << densityRatio * std::abs(vel) * diameter / ( real_t(9) * viscosity ) << ", Re = " << std::abs(vel) * diameter / viscosity << "\n";
      }
      file << "maximal overlap = " << std::abs(minGapSize) << " (" << std::abs(minGapSize)/diameter * real_t(100) << "%)\n";
      file << "maximum gap sizes:\n";
      bounce = uint_t(1);
      for( auto gapSize : maxGapSizeBetweenBounces)
      {
         file << " - " << bounce++ << ": gap size = " << gapSize << ", gap size / diameter = "  << gapSize / diameter << "\n";
      }

      file.close();
   }




   timeloopTiming.logResultOnRoot();

   return EXIT_SUCCESS;
}

} // namespace sphere_wall_collision

int main( int argc, char **argv ){
   sphere_wall_collision::main(argc, argv);
}