Bimetric 3+1 toolkit for spherical symmetry
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines
MPIWorld Class Reference

Implements the message passing interface (MPI) for high-performance computing (HPC). More...

#include <mpiWorld.h>

Inheritance diagram for MPIWorld:

Public Member Functions

 MPIWorld ()
 Initializes the MPI environment and finds our rank in the MPI world.
void cleanup ()
 Wait all outstanding MPI requests to finish, then finalize MPI.
 ~MPIWorld ()
void quit (int rc)
int size () const
int rank () const
bool isFirstInRank () const
bool isLastInRank () const
void defineGhostChunk (Int nGFs, Int nTotal, Int nGhost)
 Defines a MPI datatype for the ghost chunk, which consists of nGFs blocks where each block size is nGhost Reals, repeating every nTotal strides.
bool exchangeBoundaries (Real *leftGhost, Real *leftEdge, Real *rightEdge, Real *rightGhost, bool receiveEdges=true)
 Exchange the grid regions at boundaries with the other MPI ranks.
bool waitExchange ()
 Wait the other ranks to have received our edges.
void abortExchange (Real *dummyData)
 Nonblocking signal to all the ranks that we are quitting the integration.

Private Attributes

int worldSize
 The size of the MPI world.
int myRank
 Our rank in the MPI world.
int leftRank
 The rank of the left neighbor.
int rightRank
 The rank of the right neighbor.
MPI_Request waitLeft
 Outstanding MPI_Isend request to the left rank.
MPI_Request waitRight
 Outstanding MPI_Isend request to the right rank.
MPI_Datatype ghostChunk
 Describes a ghost-size block at each slice.

Detailed Description

Implements the message passing interface (MPI) for high-performance computing (HPC).

In the MPI environment, the calculations on the grid are split among the ranks so the edges are exchanged between each two neighbors in the following way:

        .---- left -----.-------.---- right ----.
        | ghost | edge  |       | edge  | ghost |
        +-------+-----------------------+-------+
   ...  |       |     grid chunk i      |       |
        +-------+-----------------------+-------+------------+-------+
                                |       |   grid chunk i+1   |       |  ...
                                +-------+--------------------+-------+
                                :     left      :
                                : ghost | edge  :
                                `-------v-------´
                         exchanged between two neighbors              

The overlapping edges of one rank become other's ghosts (of the same size). The ghost of the left-most rank is determined from the inner boundary conditions and the ghost of the right-most rank is determined from the outer boundary conditions.

When compiling, use mpic++ with -D_USEMPI. To run, use mpiexec or mpirun.

Warning:
The maximal slicing is not compliant with MPI since the boundary value problem for the slicing differential equation requires access to the whole grid!
Todo:
Implement the collective I/O with MPI_File_open and MPI_File_write.

Definition at line 38 of file mpiWorld.h.


Constructor & Destructor Documentation

MPIWorld::MPIWorld ( ) [inline]

Initializes the MPI environment and finds our rank in the MPI world.

Definition at line 52 of file mpiWorld.h.

    {
        // Initialize the MPI environment
        //
        MPI_Init( NULL, NULL );

        // Find out the rank and the size of the MPI
        //
        MPI_Comm_rank( MPI_COMM_WORLD, &myRank );
        MPI_Comm_size( MPI_COMM_WORLD, &worldSize );

        leftRank  = isFirstInRank() ? -1 : myRank - 1;
        rightRank = isLastInRank()  ? -1 : myRank + 1;

        waitLeft = waitRight = MPI_REQUEST_NULL;
    }
MPIWorld::~MPIWorld ( ) [inline]

Definition at line 83 of file mpiWorld.h.

    {
        cleanup ();
    }

Member Function Documentation

void MPIWorld::abortExchange ( Real dummyData) [inline]

Nonblocking signal to all the ranks that we are quitting the integration.

Definition at line 229 of file mpiWorld.h.

    {
        if( leftRank >= 0 ) {
            MPI_Request request;
            MPI_Isend( dummyData, 1, ghostChunk, leftRank, 1, MPI_COMM_WORLD, &request );
            MPI_Request_free( &request );
        }
        if( rightRank >= 0 ) {
            MPI_Request request;
            MPI_Isend( dummyData, 1, ghostChunk, rightRank, 1, MPI_COMM_WORLD, &request );
            MPI_Request_free( &request );
        }
    }
void MPIWorld::cleanup ( ) [inline]

Wait all outstanding MPI requests to finish, then finalize MPI.

Definition at line 71 of file mpiWorld.h.

    {
        if( waitLeft != MPI_REQUEST_NULL ) {
            MPI_Cancel( &waitLeft );
        }
        if( waitRight != MPI_REQUEST_NULL ) {
            MPI_Cancel( &waitRight );
        }
        MPI_Barrier( MPI_COMM_WORLD );
        MPI_Finalize ();
    }
void MPIWorld::defineGhostChunk ( Int  nGFs,
Int  nTotal,
Int  nGhost 
) [inline]

Defines a MPI datatype for the ghost chunk, which consists of nGFs blocks where each block size is nGhost Reals, repeating every nTotal strides.

We have: grid index (gf,m,n) = (m) * nGFs * nTotal + (gf) * nTotal + (n), hence:

Note that nGFs does not need to be GFCNT (we can transfer a truncated list of grid functions).

  MPI_Type_vector( count, blocklen, stride, oldtype, newtype )

     count    = nGFs,         number of blocks
     blocklen = nGhost,       number of elements in each block
     stride   = nTotal,       number of elements between start of each block
     oldtype  = MPI_DOUBLE,   old datatype (handle)
     newtype  = &ghostChunk,  new datatype (handle)
   |<- blocklen ->|     |<- blocklen ->|    ...    |<- blocklen ->|
   :                    :                                         :
   |<----- stride ----->|                                         :
   :                                                              :
    `------------------------------.-----------------------------´
                                 count
   Total size = count * blocklen                             

Definition at line 136 of file mpiWorld.h.

    {
        MPI_Datatype data_type = sizeof(Real) == sizeof(double)
                               ? MPI_DOUBLE : MPI_LONG_DOUBLE;
        MPI_Type_vector( nGFs, nGhost, nTotal, data_type, &ghostChunk );
        MPI_Type_commit( &ghostChunk );
    }
bool MPIWorld::exchangeBoundaries ( Real leftGhost,
Real leftEdge,
Real rightEdge,
Real rightGhost,
bool  receiveEdges = true 
) [inline]

Exchange the grid regions at boundaries with the other MPI ranks.

  • Send our edges (they become other's ghosts)
  • Optionally receive other's edges so they become our ghosts
Parameters:
leftGhostWhere to put the right edge of the left neighbor
leftEdgePointer to the left edge of our data
rightEdgePointer to the right edge of our data
rightGhostWhere to put the left edge of the right neighbor
receiveEdgesWhether to receive the edges (or just send ghosts)

Definition at line 147 of file mpiWorld.h.

    {
        waitLeft = waitRight = MPI_REQUEST_NULL;

        /// - Send our edges (they become other's ghosts)
        ///
        if( leftRank >= 0 ) {
            if( MPI_SUCCESS != MPI_Isend( leftEdge, 1, ghostChunk, leftRank, 0,
                                          MPI_COMM_WORLD, &waitLeft ) ) {
                return false;
            }
        }
        if( rightRank >= 0 ) {
            if( MPI_SUCCESS != MPI_Isend( rightEdge, 1, ghostChunk, rightRank, 0,
                       MPI_COMM_WORLD, &waitRight ) ) {
                return false;
            }
        }

        /// - Optionally receive other's edges so they become our ghosts
        ///
        if( ! receiveEdges ) {
            return true;
        }

        if( leftRank >= 0 )
        {
            MPI_Status status;
            if( MPI_SUCCESS != MPI_Recv( leftGhost, 1, ghostChunk,
                                         leftRank, MPI_ANY_TAG,
                                         MPI_COMM_WORLD, &status )
                || status.MPI_TAG != 0 )
            {
                // std::cerr << "[MPI #" << rank() << " got abort from left]";
                abortExchange( leftGhost );
                return false;
            }
        }
        if( rightRank >= 0 )
        {
            MPI_Status status;
            if( MPI_SUCCESS != MPI_Recv( rightGhost, 1, ghostChunk,
                                         rightRank, MPI_ANY_TAG,
                                         MPI_COMM_WORLD, &status )
                ||status.MPI_TAG != 0 )
            {
                // std::cerr << "[MPI #" << rank() << " got abort from right]";
                abortExchange( leftGhost );
                return false;
            }
        }

        return true;
    }
bool MPIWorld::isFirstInRank ( ) const [inline]

Definition at line 102 of file mpiWorld.h.

                                {
        return myRank == 0;
    }
bool MPIWorld::isLastInRank ( ) const [inline]

Definition at line 106 of file mpiWorld.h.

                               {
        return myRank == worldSize - 1;
    }
void MPIWorld::quit ( int  rc) [inline]

Definition at line 88 of file mpiWorld.h.

    {
        cleanup ();
        std::exit( rc );
    }
int MPIWorld::rank ( ) const [inline]

Definition at line 98 of file mpiWorld.h.

                      {
        return myRank;
    }
int MPIWorld::size ( ) const [inline]

Definition at line 94 of file mpiWorld.h.

                      {
        return worldSize;
    }
bool MPIWorld::waitExchange ( ) [inline]

Wait the other ranks to have received our edges.

Definition at line 210 of file mpiWorld.h.

    {
        if( waitLeft != MPI_REQUEST_NULL ) {
            MPI_Status status;
            if ( MPI_SUCCESS != MPI_Wait( &waitLeft, &status ) ) {
                return false;
            }
        }
        if( waitRight != MPI_REQUEST_NULL ) {
            MPI_Status status;
            if( MPI_SUCCESS != MPI_Wait( &waitRight, &status ) ) {
                return false;
            }
        }
        return true;
    }

Field Documentation

MPI_Datatype MPIWorld::ghostChunk [private]

Describes a ghost-size block at each slice.

Definition at line 46 of file mpiWorld.h.

int MPIWorld::leftRank [private]

The rank of the left neighbor.

Definition at line 42 of file mpiWorld.h.

int MPIWorld::myRank [private]

Our rank in the MPI world.

Definition at line 41 of file mpiWorld.h.

int MPIWorld::rightRank [private]

The rank of the right neighbor.

Definition at line 43 of file mpiWorld.h.

MPI_Request MPIWorld::waitLeft [private]

Outstanding MPI_Isend request to the left rank.

Definition at line 44 of file mpiWorld.h.

MPI_Request MPIWorld::waitRight [private]

Outstanding MPI_Isend request to the right rank.

Definition at line 45 of file mpiWorld.h.

int MPIWorld::worldSize [private]

The size of the MPI world.

Definition at line 40 of file mpiWorld.h.


The documentation for this class was generated from the following file: