lightning user document heading  
NCAR
Last update: 02/11/2005

Lightning user doc contents

SPMD parallel job example using MPI-standard message passing

This example uses MPI SPMD (Message Passing Interface-standard Single Processor Multiple Data structures) to sum a set of 1,000,000 exponentials. The Message Passing Interface is a standardized programming paradigm for parallel computers designed to allow separate processes to efficiently communicate with each other via messages.

Note that this example starts two identical processes (process 1 and process 0) that have different execution paths (tasks). The execution path in process 1 performs the calculation and sends that sum to process 0. Process 0 waits until it receives the sum from process 1, then prints the sum. The code for both processes is identical, and the execution path is determined by testing the value of the process id. This is noted by the comments in the code.

This example consists of four parts:

  • A script for the LSF batch subsystem that submits the job to lightning
  • Fortran code that runs the example
  • C code that runs the same example
  • C++ code that runs the same example

To run this example, you have two choices:

  • In the /usr/local/examples/lsf/batch/ directory on lightning:
    cp mpi.* $PWD
    Submit the example codes to the LSF batch subsystem by entering:
    bsub < mpi.lsf

  • Copy the codes on this page and paste them into your own files:

    1. Copy the LSF batch job script below and paste it into a file named mpi.lsf in your working directory on lightning.
    2. Copy the Fortran code below and paste it into a file named mpi.f in your working directory on lightning.
    3. Copy the C code below and paste it into a file named mpi.c in your working directory on lightning.
    4. Copy the C++ code below and paste it into a file named mpi.cc in your working directory on lightning.
    5. Submit the example codes to the LSF batch subsystem by entering:
      bsub < mpi.lsf

Studying this example will help you prepare your own SPMD jobs for submittal to lightning via LSF.

LSF batch job script to submit the SPMD job

#!/bin/ksh
#
# LSF batch script to run an SPMD/MPI code
#
#BSUB -a mpich_gm                       # select the mpich-gm elim
#BSUB -x                                # exlusive use of node (not_shared)
#BSUB -n 2                              # number of total tasks
#BSUB -R "span[ptile=1]"                # run 1 task per node
#BSUB -o mpilsf.out                     # output filename (%J to add job id)
#BSUB -e mpilsf.err                     # error filename
#BSUB -J mpilsf.test                    # job name
#BSUB -q regular                        # queue

# Run this executable in SPMD mode:
# Fortran example
mpif90 -Mextend -o mpif mpi.f
mpirun.lsf ./mpif
rm mpif

# C example
mpicc -o mpic mpi.c
mpirun.lsf ./mpic
rm mpic

# C++ example
mpicxx --no_auto_instantiation -o mpicc mpi.cc
mpirun.lsf ./mpicc
rm .mpicc

SPMD example code in Fortran

      program main
      implicit none
      include 'mpif.h'

      character (len=MPI_MAX_PROCESSOR_NAME) nodename
      integer name_len, ierr

! Establish the number of the process being run with the variable "rank"
      integer i
      integer rank,error,tag,length,status(MPI_STATUS_SIZE)
      integer hz, clock0, clock1, t
      real(kind=8)::  sum, buf(2), elapsed

      call mpi_init(error)
! Assign the process number to the current process with "rank"
      call mpi_comm_rank(MPI_COMM_WORLD,rank,error)
      length=2
      tag=999

! Test for the process being run, then calculate sum (process 1) or
! print sum (process 0)
      call system_clock(count_rate = hz)
      call system_clock(count = clock0)
      if (rank .eq. 1) then
          sum=0.0
          do i=1,1000000
              sum=sum+exp(.00000001*i)
          end do
          call system_clock(count = clock1)
          elapsed = real(clock1 - clock0) / hz
          buf(1)=sum
          buf(2)=elapsed/1000.0
          call mpi_send(buf,length,MPI_REAL8,0,tag,MPI_COMM_WORLD,error)
          call MPI_GET_PROCESSOR_NAME(nodename, name_len, ierr)
          print *,'f90 MPI task ', rank,' data sent from node ', nodename(1:name_len)
      elseif (rank .eq. 0) then
          call MPI_GET_PROCESSOR_NAME(nodename, name_len, ierr)
          print *,'f90 MPI task ', rank,' data received on node ', nodename(1:name_len)
          call  mpi_recv(buf,length,MPI_REAL8,1,tag,MPI_COMM_WORLD,status,error)
          print 10, buf(1), buf(2)
10        format(' f90 MPI Results: Sum = ',1pe12.6,' Loop time = ',  0pf12.8)
      end if

      call mpi_finalize(error)
      stop
      end

SPMD example code in C

#include <stdlib.h>
#include <stdio.h>
#include <math.h>
#include <mpi.h>

/* Establish the number of the process being run with the variable  "rank" */
main(int argc,char **argv)
{
        MPI_Status status;
        int i;
        int rank,error,tag=42,length=2;
        double buf[2], sum=0.0, elapsed, rtc();

        int name_len;
        char nodename[MPI_MAX_PROCESSOR_NAME];

        error = MPI_Init(&argc,&argv);
/* Assign the process number to the current process with "rank" */
        error=MPI_Comm_rank(MPI_COMM_WORLD,&rank);

/* Test for the process being run, then calculate sum (process 1) or */
/* print sum (process 0) */
        if(rank==1)
        {
                elapsed=rtc();
                for(i=1; i<1000000; i++)
                {
                        sum += exp( .00000001 * (double)i );
                }
                elapsed=rtc()-elapsed;
                buf[0]=sum;
                buf[1]=elapsed;
                error = MPI_Send(&buf, length, MPI_DOUBLE,0,tag,MPI_COMM_WORLD);
                MPI_Get_processor_name(nodename, &name_len);
                printf("   c MPI task %4d data sent from node %s\n",  rank, nodename);
        }

        else
        {
                error=MPI_Recv(&buf,length,MPI_DOUBLE,1,tag,MPI_COMM_WORLD,&status);
                MPI_Get_processor_name(nodename, &name_len);
                printf("   c MPI task %4d data received on node  %s\n",  rank, nodename);
                printf( "   c MPI Results: Sum = %11e Loop time = %8f  \n",buf[0],buf[1] );

        }
        error=MPI_Finalize();
        exit (0);
}
double rtc()
{
static struct timeval time;
gettimeofday(&time,NULL);
return ( (double)(time.tv_sec)+(double)(time.tv_usec)/1000000.0 );
}

SPMD example code in C++

#include <iostream>
#include <string>
#include <math.h>
#include <mpi.h>
#include <sys/time.h>
using namespace std;

//---------------Class  Defs-----------------------------------------------
class exp_sum {
public:
        double elapsed;
        double sum();
private:
        double summer;
        struct timeval time;
        double rtc();
};
double exp_sum::sum()
{
        int i;
        summer=0.0;
        elapsed=rtc();
        for(i=1; i<1000000; i++)
        {
                summer += exp( .00000001 * (double)i );
        }
        elapsed=rtc()-elapsed;
        return summer;
};
double exp_sum:: rtc()
{
        gettimeofday(&time,NULL);
        return ( (double)(time.tv_sec*1000000+time.tv_usec)/1000000 );
};
//----------------------------------------------------------------------- ---



//--------------Main  Program--------------------------------------------------
int main(int argc,char **argv)
{
        exp_sum total;
        MPI_Status status;
        int rank,error,tag=99,length=2;
        double buf[2];


        int name_len;
        char nodename[MPI_MAX_PROCESSOR_NAME];

        error = MPI_Init(&argc,&argv);
        error=MPI_Comm_rank(MPI_COMM_WORLD,&rank);
        if(rank==1)
        {
                buf[0]=total.sum();
                buf[1]=total.elapsed;

                error = MPI_Send(&buf, length,  MPI_DOUBLE,0,tag,MPI_COMM_WORLD);
                MPI_Get_processor_name(nodename, &name_len);
                cout << " c++ MPI task " << rank << " data sent from  node " << \
                          nodename << endl;
        }
        else
        {
                 error=MPI_Recv(&buf,length,MPI_DOUBLE,1,tag,MPI_COMM_WORLD,&status);
                MPI_Get_processor_name(nodename, &name_len);
                cout << " c++ MPI task " << rank << " data received on  node " << \
                          nodename  << endl;
                cout << " c++ mpi Results: Sum = " << buf[0] << " Loop  time = " << buf[1] << '\n';
        }
        error=MPI_Finalize();

        return 0;
}

Next page | IBM Linux cluster systems fundamentals - Table of contents

If you have questions about this document, please contact SCD Customer Support. You can also reach us by telephone 24 hours a day, seven days a week at 303-497-1278. Additional contact methods: consult1@ucar.edu and during business hours in NCAR Mesa Lab Suite 39.

© Copyright 2005. University Corporation for Atmospheric Research (UCAR). All Rights Reserved.

Address of this page: http://www.scd.ucar.edu/docs/lightning/examples/spmd.jsp