Skip to content
Snippets Groups Projects
BufferSystem.cpp 15.8 KiB
Newer Older
1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463
//======================================================================================================================
//
//  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 BufferSystem.cpp
//! \ingroup core
//! \author Martin Bauer <martin.bauer@fau.de>
//
//======================================================================================================================

#include "BufferSystem.h"
#include "core/logging/Logging.h"
#include "core/mpi/MPIManager.h"


namespace walberla {
namespace mpi {


#ifndef NDEBUG
std::set<int> BufferSystem::activeTags_;
#endif

//======================================================================================================================
//
//  Iterator
//
//======================================================================================================================



BufferSystem::iterator::iterator( BufferSystem & bufferSystem, bool begin )
    : bufferSystem_( bufferSystem), currentRecvBuffer_( NULL ), currentSenderRank_( -1 )
{
   if ( begin ) // init iterator
      ++(*this);
}


void BufferSystem::iterator::operator++()
{
   currentRecvBuffer_ = bufferSystem_.waitForNext( currentSenderRank_ );
   if ( ! currentRecvBuffer_ ) {
      WALBERLA_ASSERT_EQUAL( currentSenderRank_, -1 );
   }
}

bool BufferSystem::iterator::operator==( const BufferSystem::iterator & other )
{
   // only equality checks with end iterators are allowed
   WALBERLA_ASSERT( other.currentSenderRank_ == -1 || currentSenderRank_ == -1 );

   return ( currentSenderRank_ == other.currentSenderRank_ );
}

bool BufferSystem::iterator::operator!=( const BufferSystem::iterator & other )
{
   // only equality checks with end iterators are allowed
   WALBERLA_ASSERT( other.currentSenderRank_ == -1 || currentSenderRank_ == -1 );

   return ( currentSenderRank_ != other.currentSenderRank_ );
}




//======================================================================================================================
//
//  Constructors
//
//======================================================================================================================


BufferSystem::BufferSystem( const MPI_Comm & communicator, int tag )
   : knownSizeComm_  ( communicator, tag ),
     unknownSizeComm_( communicator, tag ),
     noMPIComm_( communicator, tag ),
     currentComm_    ( NULL ),
     sizeChangesEverytime_( true ),
     communicationRunning_( false )
{
}

//======================================================================================================================
//
//  Receive Information Setup
//
//======================================================================================================================


//**********************************************************************************************************************
/*! Sets receiver information, when message sizes are unknown
*
* \param ranksToRecvFrom  set of all ranks where messages are received
* \param changingSize     true if the message size is different in each communication step.
*                         If false the message size is exchanged once and is expected to be constant
*                         If true the message size is exchanged before each communication step.
*                         The behavior can be changed later one using setReceiverInfo() or sizeHasChanged().
*/
//**********************************************************************************************************************
void BufferSystem::setReceiverInfo( const std::set<MPIRank> & ranksToRecvFrom, bool changingSize )
{
   WALBERLA_ASSERT( ! communicationRunning_ );

   recvInfos_.clear();
   for ( auto it = ranksToRecvFrom.begin(); it != ranksToRecvFrom.end(); ++it )
   {
      const MPIRank sender = *it;
      recvInfos_[ sender ].size = INVALID_SIZE;
   }

   sizeChangesEverytime_ = changingSize;
   setCommunicationType( false ); // no size information on first run -> UnknownSizeCommunication
}



//**********************************************************************************************************************
/*! Sets receiver information, when message sizes are known
*
* \param ranksToRecvFrom  Map containing all ranks, where messages are received from, as keys
*                         and the message sizes as values.
*                         The message sizes are expected to be constant for all communication step until
*                         behavior is changed with setReceiverInfo*() or sizeHasChanged()
*/
//**********************************************************************************************************************
void BufferSystem::setReceiverInfo( const std::map<MPIRank,MPISize> & ranksToRecvFrom )
{
   WALBERLA_ASSERT( ! communicationRunning_ );

   recvInfos_.clear();
   for ( auto it = ranksToRecvFrom.begin(); it != ranksToRecvFrom.end(); ++it )
   {
      const MPIRank sender       = it->first;
      const MPISize senderSize   = it->second;
      WALBERLA_ASSERT_GREATER( senderSize, 0 );
      recvInfos_[ sender ].size   = senderSize;
   }

   sizeChangesEverytime_ = false;
   setCommunicationType( true );
}



//**********************************************************************************************************************
/*! Sets receiver information, using SendBuffers (symmetric communication)
*
* Gives the BufferSystem the information that messages are received from the same processes that we
* send to (i.e. from all ranks where SendBuffers were already filled )
* sendBuffer() has to be called before, and corresponding SendBuffers have to be filled.
*
*
* \param useSizeFromSendBuffers  If true, the sizes are expected to be known and equal to the size
*                                of the SendBuffers. SendBuffers with zero size are ignored.
*                                If false, all SendBuffers (also with zero size) are registered as ranks where
*                                messages are received from. The size is unknown, and communicated before.
*
* \param changingSize            if true the size is communicated before every communication step.
*                                if useSizeFromSendBuffer==true and changingSize==true, the size is not
*                                communicated in the first step but in all following steps.
*/
//**********************************************************************************************************************
void BufferSystem::setReceiverInfoFromSendBufferState( bool useSizeFromSendBuffers, bool changingSize )
{
   WALBERLA_ASSERT( ! communicationRunning_ );

   recvInfos_.clear();
   for ( auto it = sendInfos_.begin(); it != sendInfos_.end(); ++it )
   {
      const MPIRank sender = it->first;
      const SendBuffer & buffer = it->second.buffer;

      if ( buffer.size() == 0 && useSizeFromSendBuffers )
         continue;

      recvInfos_[ sender ].size  = useSizeFromSendBuffers ? int_c( buffer.size() )  : INVALID_SIZE;
   }

   sizeChangesEverytime_ = changingSize;

   setCommunicationType( useSizeFromSendBuffers );
}



//**********************************************************************************************************************
/*! Notifies that BufferSystem that message sizes have changed ( and optionally are changing in all following steps)
*
* Useful when setReceiverInfo was set such that BufferSystem assumes constant message sizes for all steps.
* Can only be called if no communication is currently running.
*
* \param alwaysChangingSize  if true the message sizes is communicated in all following steps, if false
*                            only in the next step.
*/
//**********************************************************************************************************************
void BufferSystem::sizeHasChanged( bool alwaysChangingSize )
{
   WALBERLA_ASSERT( ! communicationRunning_ );

   sizeChangesEverytime_ = alwaysChangingSize;
   setCommunicationType( false );
}



//======================================================================================================================
//
//  Step 1: Schedule Receives and ISends
//
//======================================================================================================================



//**********************************************************************************************************************
/*! Returns an existing SendBuffer, or creates a new one (only if !isCommunicationRunning() )
*
* \param rank  the rank where the buffer should be sent to
*/
//**********************************************************************************************************************
SendBuffer & BufferSystem::sendBuffer( MPIRank rank )
{
   return sendInfos_[rank].buffer;
}





//**********************************************************************************************************************
/*! Sends filled ( nonzero length) SendBuffers to their destination ranks
*
* If some of the SendBuffers have been sent manually before using send(int rank) they are skipped,
* only the remaining buffers are sent.
*
* If communication was not started before, it is started with this function.
*/
//**********************************************************************************************************************
void BufferSystem::sendAll()
{
   WALBERLA_ASSERT_NOT_NULLPTR( currentComm_ ); // call setReceiverInfo first!

   if ( !communicationRunning_ )
      startCommunication();

   for( auto iter = sendInfos_.begin(); iter != sendInfos_.end(); ++iter )
   {
      if ( ! iter->second.alreadySent )
      {
         if ( iter->second.buffer.size() > 0 )
            currentComm_->send( iter->first, iter->second.buffer );

         iter->second.alreadySent = true;
      }
   }
}


//**********************************************************************************************************************
/*! Sends a single SendBuffer to its destination rank.
*
* If SendBuffer is empty no message is sent.
*
* If communication was not started before, it is started with this function.
*/
//**********************************************************************************************************************
void BufferSystem::send( MPIRank rank )
{
   WALBERLA_ASSERT_NOT_NULLPTR( currentComm_ ); // call setReceiverInfo first!

   if ( !communicationRunning_ )
      startCommunication();

   auto iter = sendInfos_.find(rank);
   WALBERLA_ASSERT( iter != sendInfos_.end() );   // no send buffer was created for this rank
   WALBERLA_ASSERT( ! iter->second.alreadySent ); // this buffer has already been sent

   if ( iter->second.buffer.size() > 0 )
      currentComm_->send( rank, iter->second.buffer );

   iter->second.alreadySent = true;
}



//======================================================================================================================
//
//  Private Helper Functions
//
//======================================================================================================================


//**********************************************************************************************************************
/*! Starts communication
*
* - schedules receives and reserves space for MPI_Request vectors in the currentComm_ member
*/
//**********************************************************************************************************************
void BufferSystem::startCommunication()
{
#ifndef NDEBUG
   const auto tag = currentComm_->getTag();
   WALBERLA_ASSERT_EQUAL(activeTags_.find(tag), activeTags_.end(),
                         "Another communication with the same MPI tag is currently in progress.");
  activeTags_.insert(tag);
#endif

   WALBERLA_ASSERT( ! communicationRunning_ );

   currentComm_->scheduleReceives( recvInfos_ );
   communicationRunning_ = true;
}



//**********************************************************************************************************************
/*! Cleanup after communication
*
* - wait for sends to complete
* - clear buffers
* - manage sizeChangesEverytime
*/
//**********************************************************************************************************************
void BufferSystem::endCommunication()
{
   WALBERLA_ASSERT( communicationRunning_ );
   currentComm_->waitForSends();

   // Clear send buffers
   for( auto iter = sendInfos_.begin(); iter != sendInfos_.end(); ++iter )
   {
      iter->second.alreadySent = false;
      iter->second.buffer.clear();
   }

   // Clear receive buffers
   for( auto iter = recvInfos_.begin(); iter != recvInfos_.end(); ++iter )  {
      iter->second.buffer.clear();
   }


   if( !sizeChangesEverytime_ )
      setCommunicationType( true );

   communicationRunning_ = false;

#ifndef NDEBUG
   activeTags_.erase( activeTags_.find( currentComm_->getTag() ) );
#endif
}



//**********************************************************************************************************************
/*! Helper function for iterator
*
*  See documentation of AbstractCommunication::waitForNext()
*/
//**********************************************************************************************************************
RecvBuffer * BufferSystem::waitForNext( MPIRank & fromRank )
{
   WALBERLA_ASSERT( communicationRunning_ );

   fromRank = currentComm_->waitForNextReceive( recvInfos_ );

   if( fromRank >= 0 )
      return & ( recvInfos_[fromRank].buffer );
   else
   {
      endCommunication();
      return NULL;
   }

}



//**********************************************************************************************************************
/*! Sets the communication type to known size, unknown size or NoMPI comm
*/
//**********************************************************************************************************************
void BufferSystem::setCommunicationType( const bool knownSize )
{
   WALBERLA_NON_MPI_SECTION()
   {
      currentComm_ = &noMPIComm_;
   }

   WALBERLA_MPI_SECTION()
   {
      if( knownSize )
         currentComm_ = &knownSizeComm_;
      else
         currentComm_ = &unknownSizeComm_;
   }
}


//======================================================================================================================
//
//  Rank Ranges
//
//======================================================================================================================

// using boost::counting_range didn't work on all supported compilers
// so the range is created explicitly

BufferSystem::RankRange BufferSystem::noRanks()
{
   return RankRange ( RankCountIter( 0 ),
                      RankCountIter( 0 ) );
}

BufferSystem::RankRange BufferSystem::allRanks()
{
   int numProcesses = MPIManager::instance()->numProcesses();
   return RankRange ( RankCountIter( 0            ),
                      RankCountIter( numProcesses ) );
}

BufferSystem::RankRange BufferSystem::allRanksButRoot()
{
   int numProcesses = MPIManager::instance()->numProcesses();
   return RankRange ( RankCountIter( 1            ),
                      RankCountIter( numProcesses ) );
}

BufferSystem::RankRange BufferSystem::onlyRank( int includedRank )
{
   WALBERLA_ASSERT_LESS( includedRank, MPIManager::instance()->numProcesses() );
   return RankRange ( RankCountIter( includedRank   ),
                      RankCountIter( includedRank+1 ) );
}


BufferSystem::RankRange BufferSystem::onlyRoot()
{
   return RankRange ( RankCountIter( 0 ),
                      RankCountIter( 1 ) );
}



} // namespace mpi
} // namespace walberla