diff --git a/src/core/mpi/BufferSystem.cpp b/src/core/mpi/BufferSystem.cpp
index 551fce3e60dae2355ff448ab71e380fc9ca77437..cfe037d8fbf8a56282a9a9798630e1321920ee3a 100644
--- a/src/core/mpi/BufferSystem.cpp
+++ b/src/core/mpi/BufferSystem.cpp
@@ -53,6 +53,9 @@ void BufferSystem::iterator::operator++()
    currentRecvBuffer_ = bufferSystem_.waitForNext( currentSenderRank_ );
    if ( ! currentRecvBuffer_ ) {
       WALBERLA_ASSERT_EQUAL( currentSenderRank_, -1 );
+   } else
+   {
+      bufferSystem_.bytesReceived_ += currentRecvBuffer_->size() * sizeof(RecvBuffer::ElementType);
    }
 }
 
@@ -303,7 +306,10 @@ void BufferSystem::sendAll()
       if ( ! iter->second.alreadySent )
       {
          if ( iter->second.buffer.size() > 0 )
+         {
+            bytesSent_ += iter->second.buffer.size() * sizeof(SendBuffer::ElementType);
             currentComm_->send( iter->first, iter->second.buffer );
+         }
 
          iter->second.alreadySent = true;
       }
@@ -331,7 +337,10 @@ void BufferSystem::send( MPIRank rank )
    WALBERLA_ASSERT( ! iter->second.alreadySent ); // this buffer has already been sent
 
    if ( iter->second.buffer.size() > 0 )
+   {
+      bytesSent_ += iter->second.buffer.size() * sizeof(SendBuffer::ElementType);
       currentComm_->send( rank, iter->second.buffer );
+   }
 
    iter->second.alreadySent = true;
 }
@@ -362,6 +371,9 @@ void BufferSystem::startCommunication()
 
    currentComm_->scheduleReceives( recvInfos_ );
    communicationRunning_ = true;
+
+   bytesSent_     = 0;
+   bytesReceived_ = 0;
 }
 
 
diff --git a/src/core/mpi/BufferSystem.h b/src/core/mpi/BufferSystem.h
index a79322d1ab326518156a7114f70cd024d2a9b041..03ea8bcbe058513fbb04222df6ec089895a64f3f 100644
--- a/src/core/mpi/BufferSystem.h
+++ b/src/core/mpi/BufferSystem.h
@@ -193,6 +193,9 @@ public:
    //@}
    //*******************************************************************************************************************
 
+   int64_t getBytesSent() const { return bytesSent_; }
+   int64_t getBytesReceived() const { return bytesReceived_; }
+
 
    //* Rank Ranges     *************************************************************************************************
    /*! \name Rank Ranges  */
@@ -240,6 +243,9 @@ protected:
    //stores tags of running communications in debug mode to ensure that
    //each concurrently running communication uses different tags
    static std::set<int> activeTags_;
+
+   int64_t bytesSent_     = 0; ///< number of bytes sent during last communication
+   int64_t bytesReceived_ = 0; ///< number of bytes received during last communication
 };
 
 
diff --git a/tests/core/mpi/BufferSystemTest.cpp b/tests/core/mpi/BufferSystemTest.cpp
index a57d24123ed8269d30e71112f4780ce72737b8e6..d84d12da34e0ee84143d5ae1d7f59442f2bec2a6 100644
--- a/tests/core/mpi/BufferSystemTest.cpp
+++ b/tests/core/mpi/BufferSystemTest.cpp
@@ -110,6 +110,9 @@ void symmetricCommunication()
 
       WALBERLA_CHECK_EQUAL( receivedVal, it.rank() );
    }
+
+   WALBERLA_CHECK_EQUAL( bs.getBytesSent(), (MSG_SIZE * sizeof(int) + MSG_SIZE * mpi::BUFFER_DEBUG_OVERHEAD) * 2 );
+   WALBERLA_CHECK_EQUAL( bs.getBytesReceived(), (MSG_SIZE * sizeof(int) + MSG_SIZE * mpi::BUFFER_DEBUG_OVERHEAD) * 2 );
 }
 
 /**
@@ -175,6 +178,9 @@ void asymmetricCommunication()
          WALBERLA_CHECK( it.buffer().isEmpty() );
       }
    }
+
+   WALBERLA_CHECK_EQUAL( bs.getBytesSent(), int64_c(sizeof(int) + mpi::BUFFER_DEBUG_OVERHEAD) * int64_c(rank + rank) );
+   WALBERLA_CHECK_EQUAL( bs.getBytesReceived(), int64_c(sizeof(int) + mpi::BUFFER_DEBUG_OVERHEAD) * int64_c(leftNeighbor + rightNeighbor) );
 }