HCSITSKernels.cpp 21.1 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
//======================================================================================================================
//
//  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   HCSITSKernels.cpp
//! \author Tobias Leemann <tobias.leemann@fau.de>
//
//======================================================================================================================


/** Test Collision Detection and Insertion of contacts into the contact storage */

#include <mesa_pd/data/DataTypes.h>
#include <mesa_pd/data/ContactStorage.h>
#include <mesa_pd/data/ParticleStorage.h>
#include <mesa_pd/data/ShapeStorage.h>
#include <mesa_pd/kernel/DetectAndStoreContacts.h>
#include <mesa_pd/domain/InfiniteDomain.h>
#include <mesa_pd/data/ContactAccessor.h>
#include <mesa_pd/kernel/ParticleSelector.h>
//HCSITS Kernels

#include "mesa_pd/kernel/InitContactsForHCSITS.h"
#include "mesa_pd/kernel/InitParticlesForHCSITS.h"
#include "mesa_pd/kernel/HCSITSRelaxationStep.h"
#include "mesa_pd/kernel/IntegrateParticlesHCSITS.h"

// Communication kernels
#include "mesa_pd/mpi/ReduceProperty.h"
#include "mesa_pd/mpi/BroadcastProperty.h"

#include "mesa_pd/mpi/notifications/VelocityUpdateNotification.h"
#include "mesa_pd/mpi/notifications/VelocityCorrectionNotification.h"

#include <mesa_pd/data/ParticleAccessor.h>
#include <core/Environment.h>
#include <core/logging/Logging.h>

#include <iostream>

namespace walberla {
Sebastian Eibl's avatar
Sebastian Eibl committed
53
namespace mesa_pd {
54

Sebastian Eibl's avatar
Sebastian Eibl committed
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
class ParticleAccessorWithShape : public data::ParticleAccessor
{
public:
   ParticleAccessorWithShape(std::shared_ptr<data::ParticleStorage>& ps, std::shared_ptr<data::ShapeStorage>& ss)
      : ParticleAccessor(ps)
      , ss_(ss)
   {}

   const real_t& getMass(const size_t p_idx) const {return ss_->shapes[ps_->getShapeIDRef(p_idx)]->getMass();}
   const real_t& getInvMass(const size_t p_idx) const {return ss_->shapes[ps_->getShapeIDRef(p_idx)]->getInvMass();}

   const Mat3& getInertia(const size_t p_idx) const {return ss_->shapes[ps_->getShapeIDRef(p_idx)]->getInertiaBF();}
   const Mat3& getInvInertia(const size_t p_idx) const {return ss_->shapes[ps_->getShapeIDRef(p_idx)]->getInvInertiaBF();}

   data::BaseShape* getShape(const size_t p_idx) const {return ss_->shapes[ps_->getShapeIDRef(p_idx)].get();}
private:
   std::shared_ptr<data::ShapeStorage> ss_;
};

template<typename PStorage, typename CStorage, typename PAccessor, typename CAccessor>
class TestHCSITSKernel {
public:
   TestHCSITSKernel(PStorage &ps_, CStorage &cs_, PAccessor &pa_, CAccessor &ca_) : ps(ps_), cs(cs_), pa(pa_), ca(ca_),
      erp(real_t(1.0)), model(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticFrictionlessContact), contactThreshold(0), globalAcc(0) {}

   void operator()(real_t dt){
      // Perform Collision detection (call kernel, that stores contacts into cs)
      kernel::DetectAndStoreContacts detectAndStore(cs);
      cs.clear();

      domain::InfiniteDomain domain;
      collision_detection::AnalyticContactDetection acd;
      acd.getContactThreshold() = contactThreshold;
      ps.forEachParticlePairHalf(false, kernel::ExcludeInfiniteInfinite(), pa, detectAndStore, pa, domain, acd);

      // Create Kernels
      kernel::InitContactsForHCSITS initContacts(1);
      initContacts.setFriction(0,0,real_t(0.2));
      initContacts.setErp(real_t(erp));

      kernel::InitParticlesForHCSITS initParticles;
      initParticles.setGlobalAcceleration(globalAcc);

      kernel::HCSITSRelaxationStep relaxationStep;
      relaxationStep.setRelaxationModel(model);
      relaxationStep.setCor(real_t(0.6)); // Only effective for PGSM

      kernel::IntegrateParticlesHCSITS integration;

      mesa_pd::mpi::ReduceProperty reductionKernel;
      mesa_pd::mpi::BroadcastProperty broadcastKernel;

      // Run the HCSITS loop
      cs.forEachContact(false, kernel::SelectAll(), ca, initContacts, ca, pa);
      ps.forEachParticle(false, kernel::SelectAll(), pa, initParticles, pa, dt);

      VelocityUpdateNotification::Parameters::relaxationParam = real_t(1.0);
      reductionKernel.operator()<VelocityCorrectionNotification>(ps);
      broadcastKernel.operator()<VelocityUpdateNotification>(ps);

      VelocityUpdateNotification::Parameters::relaxationParam = real_t(0.8);
      for(int i = 0; i < 10; i++){
         cs.forEachContact(false, kernel::SelectAll(), ca, relaxationStep, ca, pa, dt);
118
119
120
         reductionKernel.operator()<VelocityCorrectionNotification>(ps);
         broadcastKernel.operator()<VelocityUpdateNotification>(ps);
      }
Sebastian Eibl's avatar
Sebastian Eibl committed
121
122
      ps.forEachParticle(false, kernel::SelectAll(), pa, integration, pa, dt);
   }
123

Sebastian Eibl's avatar
Sebastian Eibl committed
124
125
126
127
128
private:
   PStorage  &ps;
   CStorage  &cs;
   PAccessor &pa;
   CAccessor &ca;
129

Sebastian Eibl's avatar
Sebastian Eibl committed
130
131
132
133
134
135
public:
   real_t erp;
   kernel::HCSITSRelaxationStep::RelaxationModel model;
   real_t contactThreshold;
   Vec3 globalAcc;
};
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


void normalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel model)
{
   //init data structures
   auto ps = std::make_shared<data::ParticleStorage>(100);
   auto cs = std::make_shared<data::ContactStorage>(100);
   auto ss = std::make_shared<data::ShapeStorage>();
   ParticleAccessorWithShape paccessor(ps, ss);
   data::ContactAccessor caccessor(cs);
   auto density = real_t(7.874);
   auto radius = real_t(1.1);

   //Geometries for sphere and half space.
   auto smallSphere = ss->create<data::Sphere>( radius );
   auto halfSpace = ss->create<data::HalfSpace>(Vec3(0,0,1));

   ss->shapes[halfSpace]->updateMassAndInertia( density );
   ss->shapes[smallSphere]->updateMassAndInertia( density );

   // Create four slightly overlapping spheres in a row (located at x=0,2)
   auto p = ps->create();
   p->getPositionRef()        = Vec3(real_t(0), real_t(0), real_t(0));
   p->getShapeIDRef()         = smallSphere;
   p->getOwnerRef()           = walberla::mpi::MPIManager::instance()->rank();
   p->getLinearVelocityRef()  = Vec3(real_t(5), real_t(5), real_t(6));
   p->getTypeRef()            = 0;
   //WALBERLA_LOG_INFO(paccessor.ParticleAccessorWithShape::getInvMass(0));


   auto p2 = ps->create();
   p2->getPositionRef()          = Vec3(real_t(5), real_t(5), real_t(5));
   p2->getShapeIDRef()           = halfSpace;
   p2->getOwnerRef()             = walberla::mpi::MPIManager::instance()->rank();
   p2->getTypeRef()               = 0;
   data::particle_flags::set(p2->getFlagsRef(), data::particle_flags::FIXED);
   data::particle_flags::set(p2->getFlagsRef(), data::particle_flags::GLOBAL);

   TestHCSITSKernel<data::ParticleStorage, data::ContactStorage, ParticleAccessorWithShape, data::ContactAccessor> testHCSITS(*ps, *cs, paccessor, caccessor);
   testHCSITS.model = model;

   WALBERLA_LOG_INFO(paccessor.getInvMass(0))
Sebastian Eibl's avatar
Sebastian Eibl committed
178
         WALBERLA_LOG_INFO(paccessor.getInvMass(1))
179

Sebastian Eibl's avatar
Sebastian Eibl committed
180
181
182
         // plane at 5,5,5
         // radius 1.1
         p->setPosition(  Vec3(5,5,6) );
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
   p->setLinearVelocity( Vec3(0,0,0) );
   testHCSITS( real_c( real_t(1.0) ) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,real_t(6.1)) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,real_t(0.1)) );
   WALBERLA_LOG_INFO(p->getPosition());
   WALBERLA_LOG_INFO(p->getLinearVelocity());

   p->setPosition(  Vec3(5,5,6) );
   p->setLinearVelocity( Vec3(0,0,0) );
   testHCSITS.erp = real_t(0.5) ;
   testHCSITS( real_t(1.0) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,real_t(6.05)) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,real_t(0.05)) );
   WALBERLA_LOG_INFO(p->getPosition());
   WALBERLA_LOG_INFO(p->getLinearVelocity());

   p->setPosition(  Vec3(5,5,6) );
   p->setLinearVelocity( Vec3(0,0,0) );
   testHCSITS.erp = real_t(0.0) ;
   testHCSITS( real_t(1.0) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,6) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,0) );
   WALBERLA_LOG_INFO(p->getPosition());
   WALBERLA_LOG_INFO(p->getLinearVelocity());

   p->setPosition(  Vec3(5,5,6) );
   p->setLinearVelocity( Vec3(0,0,-1) );
   testHCSITS.erp = real_t(1.0) ;
   testHCSITS( real_t(1.0) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,real_t(6.1)) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,real_t(0.1)) );
   WALBERLA_LOG_INFO(p->getPosition());
   WALBERLA_LOG_INFO(p->getLinearVelocity());

   p->setPosition(  Vec3(5,5,6) );
   p->setLinearVelocity( Vec3(0,0,-1) );
   testHCSITS.erp = real_t(0.5) ;
   testHCSITS( real_t(1.0) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,real_t(6.05)) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,real_t(0.05)) );
   WALBERLA_LOG_INFO(p->getPosition());
   WALBERLA_LOG_INFO(p->getLinearVelocity());

   p->setPosition(  Vec3(5,5,6) );
   p->setLinearVelocity( Vec3(0,0,-1) );
   testHCSITS.erp = real_t(0.0) ;
   testHCSITS( real_t(1.0) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,6) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,0) );
   WALBERLA_LOG_INFO(p->getPosition());
   WALBERLA_LOG_INFO(p->getLinearVelocity());

   // No collision
   p->setPosition(  Vec3(5,5,real_t(6.2)) );
   p->setLinearVelocity( Vec3(0,0,real_t(-0.2)) );
   testHCSITS.erp = real_t(1.0) ;
   testHCSITS( real_t(1.0) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,real_t(6.0)) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,real_t(-0.2)) );

   testHCSITS.contactThreshold = real_t(1.0);
   p->setPosition(  Vec3(5,5,real_t(6.2)) );
   p->setLinearVelocity( Vec3(0,0,real_t(-0.2)) );
   testHCSITS.erp = real_t(1.0) ;
   testHCSITS( real_t(1.0) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,real_t(6.1)) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,real_t(-0.1)) );

   p->setPosition(  Vec3(5,5,real_t(6.1)) );
   p->setLinearVelocity( Vec3(0,0,real_t(-0.1)) );
   testHCSITS.erp = real_t(1.0) ;
   testHCSITS( real_t(1.0) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getPosition() , Vec3(5,5,real_t(6.1)) );
   WALBERLA_CHECK_FLOAT_EQUAL( p->getLinearVelocity(), Vec3(0,0,real_t(0)) );
   WALBERLA_LOG_INFO(p->getPosition());
   WALBERLA_LOG_INFO(p->getLinearVelocity());
}


/**Check hard contact constraints on two overlapping, colliding spheres
 * Works only for the solvers that really archieve seperation after a single
 * timestep. Use SphereSeperationTest to check for seperation after multiple
 * timesteps.
 * @param model The collision model to use.
 * */
void SphereSphereTest(kernel::HCSITSRelaxationStep::RelaxationModel model){

Sebastian Eibl's avatar
Sebastian Eibl committed
270
271
272
273
274
275
276
277
   //init data structures
   auto ps = std::make_shared<data::ParticleStorage>(100);
   auto cs = std::make_shared<data::ContactStorage>(100);
   auto ss = std::make_shared<data::ShapeStorage>();
   ParticleAccessorWithShape paccessor(ps, ss);
   data::ContactAccessor caccessor(cs);
   auto density = real_t(7.874);
   auto radius = real_t(1.1);
278

Sebastian Eibl's avatar
Sebastian Eibl committed
279
280
   auto smallSphere = ss->create<data::Sphere>( radius );
   ss->shapes[smallSphere]->updateMassAndInertia( density );
281

Sebastian Eibl's avatar
Sebastian Eibl committed
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
   auto dt = real_t(1);

   // Create two slightly overlapping spheres in a row (located at x=0,2)
   auto p = ps->create();
   p->getPositionRef()          = Vec3(real_t(0), real_t(0), real_t(0));
   p->getShapeIDRef()           = smallSphere;
   p->getOwnerRef()             = walberla::mpi::MPIManager::instance()->rank();
   p->getLinearVelocityRef()    = Vec3(real_t(1), real_t(0), real_t(0));
   p->getTypeRef()              = 0;
   auto p2 = ps->create();
   p2->getPositionRef()          = Vec3(real_t(2), real_t(0), real_t(0));
   p2->getShapeIDRef()           = smallSphere;
   p2->getOwnerRef()             = walberla::mpi::MPIManager::instance()->rank();
   p2->getLinearVelocityRef() = Vec3(real_t(-1), real_t(0), real_t(0));
   p2->getTypeRef()              = 0;
   TestHCSITSKernel<data::ParticleStorage, data::ContactStorage, ParticleAccessorWithShape, data::ContactAccessor> testHCSITS(*ps, *cs, paccessor, caccessor);
   testHCSITS.model = model;
   testHCSITS(dt);

   WALBERLA_CHECK_FLOAT_EQUAL(p->getPosition(), Vec3(real_t(-0.1),0,0));
   WALBERLA_CHECK_FLOAT_EQUAL(p->getLinearVelocity(), Vec3(real_t(-0.1),0,0));
   WALBERLA_CHECK_FLOAT_EQUAL(p->getAngularVelocity(), Vec3(0,0,0));
   WALBERLA_CHECK_FLOAT_EQUAL(p2->getPosition(), Vec3(real_t(2.1),0,0));
   WALBERLA_CHECK_FLOAT_EQUAL(p2->getLinearVelocity(), Vec3(real_t(0.1),0,0))
         WALBERLA_CHECK_FLOAT_EQUAL(p2->getAngularVelocity(), Vec3(0,0,0));

   WALBERLA_LOG_INFO(p->getPosition());
   WALBERLA_LOG_INFO(p->getLinearVelocity());
   WALBERLA_LOG_INFO(p2->getPosition());
   WALBERLA_LOG_INFO(p2->getLinearVelocity());
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
}

/**
 * Create two overlapping spheres with opposing velocities, that are in contact.
 * Assert that the collision is resolved after a certain (10) number of timesteps.
 * @param model The collision model to use.
 */
void SphereSeperationTest(kernel::HCSITSRelaxationStep::RelaxationModel model){

   //init data structures
   auto ps = std::make_shared<data::ParticleStorage>(100);
   auto cs = std::make_shared<data::ContactStorage>(100);
   auto ss = std::make_shared<data::ShapeStorage>();
   ParticleAccessorWithShape paccessor(ps, ss);
   data::ContactAccessor caccessor(cs);
   auto density = real_t(7.874);
   auto radius = real_t(1.1);

   auto smallSphere = ss->create<data::Sphere>( radius );

   ss->shapes[smallSphere]->updateMassAndInertia( density );
   auto dt = real_t(0.2);

   // Create two slightly overlapping spheres in a row (located at x=0,2)
   auto p = ps->create();
   p->getPositionRef()          = Vec3(real_t(0), real_t(0), real_t(0));
   p->getShapeIDRef()           = smallSphere;
   p->getOwnerRef()             = walberla::mpi::MPIManager::instance()->rank();
   p->getLinearVelocityRef()    = Vec3(real_t(1), real_t(0), real_t(0));
   p->getTypeRef()              = 0;
   auto p2 = ps->create();
   p2->getPositionRef()          = Vec3(real_t(2.0), real_t(0), real_t(0));
   p2->getShapeIDRef()           = smallSphere;
   p2->getOwnerRef()             = walberla::mpi::MPIManager::instance()->rank();
   p2->getLinearVelocityRef()    = Vec3(real_t(-1), real_t(0), real_t(0));
   p2->getTypeRef()              = 0;
   TestHCSITSKernel<data::ParticleStorage, data::ContactStorage, ParticleAccessorWithShape, data::ContactAccessor> testHCSITS(*ps, *cs, paccessor, caccessor);

   int solveCount = 0;
   testHCSITS.model = model;

   // Number of allowed iterations
   int maxIter = 10;
   while(p2->getPosition()[0]-p->getPosition()[0] < 2.2){
      testHCSITS(dt);
      WALBERLA_LOG_INFO(p->getPosition());
      WALBERLA_LOG_INFO(p->getLinearVelocity());
      WALBERLA_LOG_INFO(p2->getPosition());
      WALBERLA_LOG_INFO(p2->getLinearVelocity());
      solveCount ++;
      if(solveCount==maxIter){
         WALBERLA_CHECK(false, "Seperation did not occur after " << maxIter << " Iterations performed.");
      }
   }
   WALBERLA_LOG_INFO("Seperation achieved after " << solveCount << " iterations.");
}

/**
 * Create a sphere that slides on a plane with only linear velocity. It is
 * put into a spin by the frictional reactions, until it rolls with no slip.
 * Assert that the sphere rolls with all slip resolved and at the physically correct
 * speed after a certain number of timesteps. Use only with solvers that obey the coulomb friction law.
 * @param model The collision model to use.
 */
void SlidingSphereFrictionalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel model){

   //init data structures
   auto ps = std::make_shared<data::ParticleStorage>(100);
   auto cs = std::make_shared<data::ContactStorage>(100);
   auto ss = std::make_shared<data::ShapeStorage>();
   ParticleAccessorWithShape paccessor(ps, ss);
   data::ContactAccessor caccessor(cs);
   auto density = real_t(1);
   auto radius = real_t(1);

   auto smallSphere = ss->create<data::Sphere>( radius );
   auto halfSpace = ss->create<data::HalfSpace>(Vec3(0,0,1));
   ss->shapes[smallSphere]->updateMassAndInertia( density );
   ss->shapes[halfSpace]->updateMassAndInertia( density );
   auto dt = real_t(0.002);

   // Create a spheres (located at x=0, height = 1)
   auto p = ps->create();
   p->getPositionRef()          = Vec3(real_t(0), real_t(0), real_t(1));
   p->getShapeIDRef()           = smallSphere;
   p->getOwnerRef()             = walberla::mpi::MPIManager::instance()->rank();
   p->getLinearVelocityRef()    = Vec3(real_t(5), real_t(0), real_t(0));
   p->getTypeRef()              = 0;

   auto p2 = ps->create();
   p2->getPositionRef()          = Vec3(real_t(0), real_t(0), real_t(0));
   p2->getShapeIDRef()           = halfSpace;
   p2->getOwnerRef()             = walberla::mpi::MPIManager::instance()->rank();
   p2->getTypeRef()               = 0;
   data::particle_flags::set(p2->getFlagsRef(), data::particle_flags::FIXED);
   data::particle_flags::set(p2->getFlagsRef(), data::particle_flags::GLOBAL);

   TestHCSITSKernel<data::ParticleStorage, data::ContactStorage, ParticleAccessorWithShape, data::ContactAccessor> testHCSITS(*ps, *cs, paccessor, caccessor);
   testHCSITS.model = model;
   testHCSITS.globalAcc = Vec3(0,0,-10);
   int solveCount = 0;


   // Number of allowed iterations
   int maxIter = 500;
   while(!walberla::floatIsEqual(p->getAngularVelocity()[1],p->getLinearVelocity()[0], real_t(0.002))){
      testHCSITS(dt);
      if(solveCount % 50 == 0){
         WALBERLA_LOG_INFO(p->getAngularVelocity());
         WALBERLA_LOG_INFO(p->getLinearVelocity());
      }
      solveCount ++;
      if(solveCount==maxIter){
         WALBERLA_CHECK(false, "End of slip did not occur after " << maxIter << " Iterations performed.");
      }
   }
   WALBERLA_LOG_INFO("Rolling with no slip achieved after " << solveCount << " iterations.");

   // Check if the value obtained values equal the physically correct values
   // (which can be determined by newtons equation to be 25/7).
   WALBERLA_CHECK_FLOAT_EQUAL_EPSILON(real_t(25.0/7.0), p->getAngularVelocity()[1], real_t(0.01), "Angular velocity is not physically correct");
   WALBERLA_CHECK_FLOAT_EQUAL_EPSILON(real_t(25.0/7.0), p->getLinearVelocity()[0], real_t(0.01), "Linear velocity is not physically correct");
}

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

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

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


   WALBERLA_LOG_INFO_ON_ROOT("InelasticFrictionlessContact");
   normalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticFrictionlessContact);
   SphereSphereTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticFrictionlessContact);

   WALBERLA_LOG_INFO_ON_ROOT("ApproximateInelasticCoulombContactByDecoupling");
   normalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::ApproximateInelasticCoulombContactByDecoupling);
   SphereSphereTest(kernel::HCSITSRelaxationStep::RelaxationModel::ApproximateInelasticCoulombContactByDecoupling);
   SlidingSphereFrictionalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::ApproximateInelasticCoulombContactByDecoupling);

   WALBERLA_LOG_INFO_ON_ROOT("InelasticCoulombContactByDecoupling");
   normalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticCoulombContactByDecoupling);
   SphereSphereTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticCoulombContactByDecoupling);
   SlidingSphereFrictionalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticCoulombContactByDecoupling);

   WALBERLA_LOG_INFO_ON_ROOT("InelasticGeneralizedMaximumDissipationContact");
   normalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticGeneralizedMaximumDissipationContact);
   SphereSphereTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticGeneralizedMaximumDissipationContact);
   SlidingSphereFrictionalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticGeneralizedMaximumDissipationContact);

   WALBERLA_LOG_INFO_ON_ROOT("InelasticCoulombContactByOrthogonalProjections");
   SphereSeperationTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticCoulombContactByOrthogonalProjections);
   SlidingSphereFrictionalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticCoulombContactByOrthogonalProjections);

   WALBERLA_LOG_INFO_ON_ROOT("ApproximateInelasticCoulombContactByOrthogonalProjections");
   SphereSeperationTest(kernel::HCSITSRelaxationStep::RelaxationModel::ApproximateInelasticCoulombContactByOrthogonalProjections);
   SlidingSphereFrictionalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::ApproximateInelasticCoulombContactByOrthogonalProjections);

   WALBERLA_LOG_INFO_ON_ROOT("InelasticProjectedGaussSeidel");
   SphereSeperationTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticProjectedGaussSeidel);
   SlidingSphereFrictionalReactionTest(kernel::HCSITSRelaxationStep::RelaxationModel::InelasticProjectedGaussSeidel);
   return EXIT_SUCCESS;
}

Sebastian Eibl's avatar
Sebastian Eibl committed
479
} //namespace mesa_pd
480
481
482
483
} //namespace walberla

int main( int argc, char ** argv )
{
Sebastian Eibl's avatar
Sebastian Eibl committed
484
   return walberla::mesa_pd::main(argc, argv);
485
}