root/ompi/tools/mpisync/mpigclock.c

/* [<][>][^][v][top][bottom][index][help] */

DEFINITIONS

This source file includes following definitions.
  1. mpigclock_sync_linear
  2. mpigclock_sync_log
  3. mpigclock_measure_offset_adaptive

   1 /*
   2  * Copyright (c) 2010-2011, Siberian State University of Telecommunications
   3  *                          and Information Sciences. All rights reserved.
   4  * Copyright (c) 2010-2011, A.V. Rzhanov Institute of Semiconductor Physics SB RAS.
   5  *                          All rights reserved.
   6  *
   7  * mpigclock.c: MPI clock synchronization.
   8  * http://mpiperf.cpct.sibsutis.ru/index.php
   9  *
  10  * Copyright (C) 2011 Mikhail Kurnosov <mkurnosov@gmail.com>
  11  *
  12  * Redistribution and use in source and binary forms, with or without
  13  * modification, are permitted provided that the following conditions are met:
  14  * Redistributions of source code must retain the above copyright
  15  * notice, this list of conditions and the following disclaimer.
  16  *
  17  * Redistributions in binary form must reproduce the above copyright
  18  * notice, this list of conditions and the following disclaimer in the
  19  * documentation and/or other materials provided with the distribution.
  20  *
  21  * Neither the name of the the copyright holders nor the
  22  * names of its contributors may be used to endorse or promote products
  23  * derived from this software without specific prior written permission.
  24  *
  25  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
  26  * ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
  27  * WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
  28  * DISCLAIMED. IN NO EVENT SHALL <COPYRIGHT HOLDER> BE LIABLE FOR ANY
  29  * DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
  30  * (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
  31  * LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
  32  * ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
  33  * (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
  34  * SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
  35  *
  36  * $COPYRIGHT$
  37  *
  38  * Additional copyrights may follow
  39  *
  40  * $HEADER$
  41  */
  42 
  43 #include <stdio.h>
  44 #include <stdlib.h>
  45 
  46 #include <mpi.h>
  47 
  48 #include "mpigclock.h"
  49 #include "hpctimer.h"
  50 
  51 #define INVALIDTIME -1.0
  52 #define MPIGCLOCK_RTTMIN_NOTCHANGED_MAX 100
  53 #define MPIGCLOCK_MSGTAG 128
  54 
  55 static double mpigclock_measure_offset_adaptive(MPI_Comm comm, int root, int peer, double *min_rtt, double root_offset);
  56 
  57 
  58 /*
  59  * mpigclock_sync_linear: Clock synchronization algorithm with O(n) steps.
  60  */
  61 double mpigclock_sync_linear(MPI_Comm comm, int root, double *rtt)
  62 {
  63     int peer, rank, commsize;
  64     double ret = 0;
  65 
  66     MPI_Comm_rank(comm, &rank);
  67     MPI_Comm_size(comm, &commsize);
  68 
  69     for (peer = 1; peer < commsize; peer++) {
  70         MPI_Barrier(comm);
  71         if (rank == root || rank == peer) {
  72             ret = mpigclock_measure_offset_adaptive(comm, root, peer, rtt, 0.0);
  73         }
  74     }
  75     return ret;
  76 }
  77 
  78 /*
  79  * mpigclock_sync_log: Clock synchronization algorithm with O(logn) steps.
  80  * rtt argumrnt does not have a meaning
  81  */
  82 double mpigclock_sync_log(MPI_Comm comm, int root, double *rtt)
  83 {
  84     int peer, rank, commsize;
  85     double root_offset = 0;
  86     double ret = 0;
  87 
  88     MPI_Comm_rank(comm, &rank);
  89     MPI_Comm_size(comm, &commsize);
  90 
  91     /* I am a peer */
  92     if (rank != root) {
  93         root = (rank - 1) / 2;
  94         MPI_Recv(&root_offset, 1, MPI_DOUBLE, root, MPIGCLOCK_MSGTAG, comm, MPI_STATUS_IGNORE);
  95         ret = mpigclock_measure_offset_adaptive(comm, root, rank, rtt, root_offset);
  96     }
  97 
  98     root_offset = ret;
  99 
 100     /* I am a root */
 101     *rtt = 0;
 102     peer = 2 * rank + 1;
 103     if (peer < commsize) {
 104         MPI_Send(&root_offset, 1, MPI_DOUBLE, peer, MPIGCLOCK_MSGTAG, comm);
 105         mpigclock_measure_offset_adaptive(comm, rank, peer, rtt, root_offset);
 106     }
 107 
 108     *rtt = 0;
 109     peer = 2 * rank + 2;
 110     if (peer < commsize) {
 111         MPI_Send(&root_offset, 1, MPI_DOUBLE, peer, MPIGCLOCK_MSGTAG, comm);
 112         mpigclock_measure_offset_adaptive(comm, rank, peer, rtt, root_offset);
 113     }
 114 
 115     return ret;
 116 }
 117 
 118 /* mpigclock_measure_offset_adaptive: Measures clock's offset of peer. */
 119 static double mpigclock_measure_offset_adaptive(MPI_Comm comm, int root, int peer, double *min_rtt, double root_offset)
 120 {
 121     int rank, commsize, rttmin_notchanged = 0;
 122     double starttime, stoptime, peertime, rtt, rttmin = 1E12,
 123            invalidtime = INVALIDTIME, offset;
 124 
 125     MPI_Comm_rank(comm, &rank);
 126     MPI_Comm_size(comm, &commsize);
 127 
 128     offset = 0.0;
 129     for (;;) {
 130         if (rank != root) {
 131             /* Peer process */
 132             starttime = hpctimer_wtime();
 133             MPI_Send(&starttime, 1, MPI_DOUBLE, root, MPIGCLOCK_MSGTAG, comm);
 134             MPI_Recv(&peertime, 1, MPI_DOUBLE, root, MPIGCLOCK_MSGTAG, comm,
 135                      MPI_STATUS_IGNORE);
 136             stoptime = hpctimer_wtime();
 137             rtt = stoptime - starttime;
 138 
 139             if (rtt < rttmin) {
 140                 rttmin = rtt;
 141                 rttmin_notchanged = 0;
 142                 offset =  peertime - rtt / 2.0 - starttime;
 143             } else {
 144                 if (++rttmin_notchanged == MPIGCLOCK_RTTMIN_NOTCHANGED_MAX) {
 145                     MPI_Send(&invalidtime, 1, MPI_DOUBLE, root, MPIGCLOCK_MSGTAG,
 146                              comm);
 147                     break;
 148                 }
 149             }
 150         } else {
 151             /* Root process */
 152             MPI_Recv(&starttime, 1, MPI_DOUBLE, peer, MPIGCLOCK_MSGTAG, comm,
 153                      MPI_STATUS_IGNORE);
 154             peertime = hpctimer_wtime() + root_offset;
 155             if (starttime < 0.0) {
 156                 break;
 157             }
 158             MPI_Send(&peertime, 1, MPI_DOUBLE, peer, MPIGCLOCK_MSGTAG, comm);
 159         }
 160     } /* for */
 161 
 162     if( rank != root ){
 163         *min_rtt = rttmin;
 164     } else {
 165         rtt = 0.0;
 166     }
 167     return offset;
 168 }

/* [<][>][^][v][top][bottom][index][help] */