GenerateAnalyticContacts.cpp 6.61 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
//======================================================================================================================
//
//  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/>.
//
16
//! \file
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
//! \author Sebastian Eibl <sebastian.eibl@fau.de>
//
//======================================================================================================================

#include <mesa_pd/collision_detection/AnalyticContactDetection.h>
#include <mesa_pd/data/LinkedCells.h>
#include <mesa_pd/data/ParticleAccessor.h>
#include <mesa_pd/data/ParticleStorage.h>
#include <mesa_pd/data/ShapeStorage.h>
#include <mesa_pd/domain/BlockForestDomain.h>
#include <mesa_pd/kernel/DoubleCast.h>
#include <mesa_pd/kernel/InsertParticleIntoLinkedCells.h>
#include <mesa_pd/kernel/ParticleSelector.h>
#include <mesa_pd/mpi/ContactFilter.h>
#include <mesa_pd/mpi/notifications/ForceTorqueNotification.h>
#include <mesa_pd/mpi/ReduceProperty.h>
#include <mesa_pd/mpi/SyncNextNeighbors.h>

#include <blockforest/BlockForest.h>
#include <blockforest/Initialization.h>
#include <core/Environment.h>
#include <core/grid_generator/SCIterator.h>
#include <core/logging/Logging.h>
#include <core/mpi/Reduce.h>

#include <iostream>
#include <memory>

namespace walberla {
namespace mesa_pd {

class ParticleAccessorWithShape : public data::ParticleAccessor
{
public:
   ParticleAccessorWithShape(std::shared_ptr<data::ParticleStorage>& ps, std::shared_ptr<data::ShapeStorage>& ss)
52
53
         : ParticleAccessor(ps)
         , ss_(ss)
54
55
   {}

56
   const auto& getInvMass(const size_t p_idx) const {return ss_->shapes[ps_->getShapeID(p_idx)]->getInvMass();}
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

   const auto& getInvInertiaBF(const size_t p_idx) const {return ss_->shapes[ps_->getShapeID(p_idx)]->getInvInertiaBF();}

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

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

   //init domain partitioning
   auto forest = blockforest::createBlockForest( AABB(0,0,0,3,3,3), // simulation domain
                                                 Vector3<uint_t>(3,3,3), // blocks in each direction
                                                 Vector3<bool>(true, true, true) // periodicity
                                                 );
   domain::BlockForestDomain domain(forest);

   //init data structures
   auto ps = std::make_shared<data::ParticleStorage>(100);
   auto ss = std::make_shared<data::ShapeStorage>();
   ParticleAccessorWithShape accessor(ps, ss);
82
   data::LinkedCells      lc(math::AABB(-1,-1,-1,4,4,4), real_t(1.3));
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172

   //initialize particles
   const real_t radius  = real_t(0.6);
   const real_t spacing = real_t(1.0);
   auto smallSphere = ss->create<data::Sphere>( radius );
   ss->shapes[smallSphere]->updateMassAndInertia(real_t(2707));

   WALBERLA_CHECK_EQUAL(forest->size(), 1);
   const Block& blk = *static_cast<blockforest::Block*>(&*forest->begin());

   for (auto pt : grid_generator::SCGrid(blk.getAABB(), Vec3(spacing, spacing, spacing) * real_c(0.5), spacing))
   {
      data::Particle&& p          = *ps->create();
      p.getPositionRef()          = pt;
      p.getInteractionRadiusRef() = radius;
      p.getShapeIDRef()           = smallSphere;
      p.getOwnerRef()             = walberla::mpi::MPIManager::instance()->rank();
   }

   //init kernels
   kernel::InsertParticleIntoLinkedCells  insert_particle_into_linked_cells;
   mpi::ReduceProperty                    RP;
   mpi::SyncNextNeighbors                 SNN;


   SNN(*ps, domain);
   ps->forEachParticlePairHalf(false,
                               kernel::ExcludeInfiniteInfinite(),
                               accessor,
                               [&](const size_t idx1, const size_t idx2, auto& ac)
   {
      collision_detection::AnalyticContactDetection         acd;
      kernel::DoubleCast               double_cast;
      mpi::ContactFilter               contact_filter;
      if (double_cast(idx1, idx2, ac, acd, ac ))
      {
         if (contact_filter(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), domain))
         {
            ac.getForceRef(acd.getIdx1()) += acd.getContactNormal() + Vec3(1,0,0);
            ac.getForceRef(acd.getIdx2()) -= acd.getContactNormal() - Vec3(1,0,0);
         }
      }
   },
   accessor );

   RP.operator()<ForceTorqueNotification>(*ps);

   WALBERLA_CHECK_FLOAT_EQUAL(ps->getForce(0), Vec3(6,0,0));

   //TEST WITH LINKED CELLS
   ps->forEachParticle(false, kernel::SelectAll(), accessor, [](size_t idx, auto& ac) {ac.setForce(idx, Vec3(0,0,0));}, accessor);
   lc.clear();
   ps->forEachParticle(false, kernel::SelectAll(), accessor, insert_particle_into_linked_cells, accessor, lc);

   lc.forEachParticlePairHalf(false,
                              kernel::ExcludeInfiniteInfinite(),
                              accessor,
                              [&domain](const size_t idx1, const size_t idx2, auto& ac)
   {
      collision_detection::AnalyticContactDetection         acd;
      kernel::DoubleCast               double_cast;
      mpi::ContactFilter               contact_filter;
      if (double_cast(idx1, idx2, ac, acd, ac ))
      {
         if (contact_filter(acd.getIdx1(), acd.getIdx2(), ac, acd.getContactPoint(), domain))
         {
            ac.getForceRef(acd.getIdx1()) += acd.getContactNormal() + Vec3(1,0,0);
            ac.getForceRef(acd.getIdx2()) -= acd.getContactNormal() - Vec3(1,0,0);
         }
      }
   },
   accessor );

   RP.operator()<ForceTorqueNotification>(*ps);

   WALBERLA_CHECK_FLOAT_EQUAL(ps->getForce(0), Vec3(6,0,0));

   //ps->forEachParticle(false, [](data::ParticleStorage& ps, size_t idx) {WALBERLA_LOG_DEVEL_ON_ROOT(*ps[idx]);});
   //cs.forEachContact(false, [](data::ContactStorage& cs, size_t idx) {WALBERLA_LOG_DEVEL_ON_ROOT(*cs[idx]);});

   return EXIT_SUCCESS;
}

} //namespace mesa_pd
} //namespace walberla

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