SyncNextNeighbors.cpp 10.1 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
//! \author Sebastian Eibl <sebastian.eibl@fau.de>
//
//======================================================================================================================

//======================================================================================================================
//
//  THIS FILE IS GENERATED - PLEASE CHANGE THE TEMPLATE !!!
//
//======================================================================================================================

#include "SyncNextNeighbors.h"

Sebastian Eibl's avatar
Sebastian Eibl committed
29
30
#include <mesa_pd/mpi/RemoveAndNotify.h>

31
32
33
34
35
36
37
38
39
40
namespace walberla {
namespace mesa_pd {
namespace mpi {

void SyncNextNeighbors::operator()(data::ParticleStorage& ps,
                                   const domain::IDomain& domain,
                                   const real_t dx) const
{
   if (numProcesses_ == 1) return;

41
42
   walberla::mpi::BufferSystem bs( walberla::mpi::MPIManager::instance()->comm() );

43
44
45
46
47
48
49
50
51
   neighborRanks_ = domain.getNeighborProcesses();
   for( uint_t nbProcessRank : neighborRanks_ )
   {
      if (bs.sendBuffer(nbProcessRank).isEmpty())
      {
         // fill empty buffers with a dummy byte to force transmission
         bs.sendBuffer(nbProcessRank) << walberla::uint8_c(0);
      }
   }
52
   generateSynchronizationMessages(bs, ps, domain, dx);
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70

   // size of buffer is unknown and changes with each send
   bs.setReceiverInfoFromSendBufferState(false, true);
   bs.sendAll();

   // Receiving the updates for the remote rigid bodies from the connected processes
   WALBERLA_LOG_DETAIL( "Parsing of particle synchronization response starts..." );
   ParseMessage parseMessage;
   for( auto it = bs.begin(); it != bs.end(); ++it )
   {
      walberla::uint8_t tmp;
      it.buffer() >> tmp;
      while( !it.buffer().isEmpty() )
      {
         parseMessage(it.rank(), it.buffer(), ps, domain);
      }
   }
   WALBERLA_LOG_DETAIL( "Parsing of particle synchronization response ended." );
71
72
73
74
75

   bytesSent_ = bs.getBytesSent();
   bytesReceived_ = bs.getBytesReceived();
   numberOfSends_ = bs.getNumberOfSends();
   numberOfReceives_ = bs.getNumberOfReceives();
76
77
}

78
79
void SyncNextNeighbors::generateSynchronizationMessages(walberla::mpi::BufferSystem& bs,
                                                        data::ParticleStorage& ps,
80
81
82
83
84
85
86
87
88
89
                                                        const domain::IDomain& domain,
                                                        const real_t dx) const
{
   const uint_t ownRank = uint_c(rank_);

   WALBERLA_LOG_DETAIL( "Assembling of particle synchronization message starts..." );

   // position update
   for( auto pIt = ps.begin(); pIt != ps.end(); )
   {
90
91
      WALBERLA_ASSERT_GREATER(pIt->getInteractionRadius(), 0_r, "Did you forget to set the interaction radius?");

92
93
94
95
96
97
98
99
100
101
102
103
104
105
      //skip all ghost particles
      if (data::particle_flags::isSet( pIt->getFlags(), data::particle_flags::GHOST))
      {
         ++pIt;
         continue;
      }

      //skip all particles that do not communicate (create ghost particles) on other processes
      if (data::particle_flags::isSet( pIt->getFlags(), data::particle_flags::NON_COMMUNICATING))
      {
         ++pIt;
         continue;
      }

106
      if (domain.isContainedInLocalSubdomain(pIt->getPosition(), pIt->getInteractionRadius() + dx))
107
108
109
110
111
112
      {
         //no sync needed
         //just delete ghost particles if there are any

         for (const auto& ghostOwner : pIt->getGhostOwners() )
         {
Sebastian Eibl's avatar
Sebastian Eibl committed
113
            auto& buffer( bs.sendBuffer(static_cast<walberla::mpi::MPIRank>(ghostOwner)) );
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

            WALBERLA_LOG_DETAIL( "Sending removal notification for particle " << pIt->getUid() << " to process " << ghostOwner );

            packNotification(buffer, ParticleRemovalNotification( *pIt ));
         }

         pIt->getGhostOwnersRef().clear();

         ++pIt;
         continue;
      }

      //correct position to make sure particle is always inside the domain!
      //everything is decided by the master particle therefore ghost particles are not touched
      if (!data::particle_flags::isSet( pIt->getFlags(), data::particle_flags::FIXED) &&
          !data::particle_flags::isSet( pIt->getFlags(), data::particle_flags::GHOST))
      {
         domain.periodicallyMapToDomain( pIt->getPositionRef() );
      }

      // Note: At this point we know that the particle was locally owned before the position update.
      WALBERLA_CHECK_EQUAL(pIt->getOwner(), ownRank);

      WALBERLA_LOG_DETAIL( "Processing local particle " << pIt->getUid() );

      // Update nearest neighbor processes.
      for( uint_t nbProcessRank : neighborRanks_ )
      {
         if( domain.intersectsWithProcessSubdomain( nbProcessRank, pIt->getPosition(), pIt->getInteractionRadius() + dx ) )
         {
            auto ghostOwnerIt = std::find( pIt->getGhostOwners().begin(), pIt->getGhostOwners().end(), nbProcessRank );
            if( ghostOwnerIt != pIt->getGhostOwners().end() )
            {
               // already a ghost there -> update
               auto& buffer( bs.sendBuffer(nbProcessRank) );
               WALBERLA_LOG_DETAIL( "Sending update notification for particle " << pIt->getUid() << " to process " << (nbProcessRank) );
               packNotification(buffer, ParticleUpdateNotification( *pIt ));
            } else
            {
               // no ghost there -> create ghost
               auto& buffer( bs.sendBuffer(nbProcessRank) );
155
156
               WALBERLA_LOG_DETAIL( "Sending ghost copy notification for particle " << pIt->getUid() << " to process " << (nbProcessRank) );
               packNotification(buffer, ParticleGhostCopyNotification( *pIt ));
Sebastian Eibl's avatar
Sebastian Eibl committed
157
               pIt->getGhostOwnersRef().insert( int_c(nbProcessRank) );
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
            }
         }
         else
         {
            //no overlap with neighboring process -> delete if ghost is there
            auto ghostOwnerIt = std::find( pIt->getGhostOwners().begin(), pIt->getGhostOwners().end(), nbProcessRank );
            if( ghostOwnerIt != pIt->getGhostOwners().end() )
            {
               // In case the rigid particle no longer intersects the remote process nor interacts with it but is registered,
               // send removal notification.
               auto& buffer( bs.sendBuffer(nbProcessRank) );

               WALBERLA_LOG_DETAIL( "Sending removal notification for particle " << pIt->getUid() << " to process " << nbProcessRank );

               packNotification(buffer, ParticleRemovalNotification( *pIt ));

               pIt->getGhostOwnersRef().erase(ghostOwnerIt);
            }
         }
      }

      //particle has left subdomain?
      const auto ownerRank = domain.findContainingProcessRank( pIt->getPosition() );
      if( ownerRank != int_c(ownRank) )
      {
         WALBERLA_LOG_DETAIL( "Local particle " << pIt->getUid() << " is no longer on process " << ownRank << " but on process " << ownerRank );

         if( ownerRank < 0 ) {
            // No owner found: Outflow condition.
            WALBERLA_LOG_DETAIL( "Sending deletion notifications for particle " << pIt->getUid() << " due to outflow." );

            // Registered processes receive removal notification in the remove() routine.
            pIt = removeAndNotify( bs, ps, pIt );

            continue;
         }

         WALBERLA_LOG_DETAIL( "Sending migration notification for particle " << pIt->getUid() << " to process " << ownerRank << "." );
         //WALBERLA_LOG_DETAIL( "Process registration list before migration: " << pIt->getGhostOwners() );

         // Set new owner and transform to ghost particle
         pIt->setOwner(ownerRank);
         data::particle_flags::set( pIt->getFlagsRef(), data::particle_flags::GHOST );

         // currently position is mapped to periodically to global domain,
         // this might not be the correct position for a ghost particle
         domain.correctParticlePosition( pIt->getPositionRef() );

         // Correct registration list (exclude new owner and us - the old owner) and
         // notify registered processes (except for new owner) of (remote) migration since they possess a ghost particle.
         auto ownerIt = std::find( pIt->getGhostOwners().begin(), pIt->getGhostOwners().end(), ownerRank );
         WALBERLA_CHECK_UNEQUAL(ownerIt, pIt->getGhostOwners().end(), "New owner has to be former ghost owner!" );

         pIt->getGhostOwnersRef().erase( ownerIt );

         for( auto ghostRank : pIt->getGhostOwners() )
         {
Sebastian Eibl's avatar
Sebastian Eibl committed
215
            auto& buffer( bs.sendBuffer(static_cast<walberla::mpi::MPIRank>(ghostRank)) );
216
217
218
219
220
221
222

            WALBERLA_LOG_DETAIL( "Sending remote migration notification for particle " << pIt->getUid() <<
                                 " to process " << ghostRank );

            packNotification(buffer, ParticleRemoteMigrationNotification( *pIt, ownerRank ));
         }

Sebastian Eibl's avatar
Sebastian Eibl committed
223
         pIt->getGhostOwnersRef().insert( int_c(ownRank) );
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247

         // Send migration notification to new owner
         auto& buffer( bs.sendBuffer(ownerRank) );
         packNotification(buffer, ParticleMigrationNotification( *pIt ));

         pIt->getGhostOwnersRef().clear();

         continue;

      } else
      {
         // particle still is locally owned after position update.
         WALBERLA_LOG_DETAIL( "Owner of particle " << pIt->getUid() << " is still process " << pIt->getOwner() );
      }

      ++pIt;
   }

   WALBERLA_LOG_DETAIL( "Assembling of particle synchronization message ended." );
}

}  // namespace mpi
}  // namespace mesa_pd
}  // namespace walberla