Statistics
| Branch: | Revision:

root / synthbench / skampi / synchronize.c

History | View | Annotate | Download (8.1 kB)

1
/* SKaMPI   MPI-Benchmark
2

3
Copyright 2003-2008 Werner Augustin, Lars Diesselberg
4
Lehrstuhl Informatik fuer Naturwissenschaftler und Ingenieure 
5
Fakultaet fuer Informatik
6
University of Karlsruhe
7

8
This program is free software; you can redistribute it and/or modify
9
it under the terms of version 2 of the GNU General Public License as
10
published by the Free Software Foundation
11

12
This program is distributed in the hope that it will be useful,
13
but WITHOUT ANY WARRANTY; without even the implied warranty of
14
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
15
GNU General Public License for more details.
16

17
You should have received a copy of the GNU General Public License
18
along with this program; if not, write to the Free Software
19
Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA */
20

    
21
#include <mpi.h>
22
#include <stdlib.h>
23
#include <stdio.h>
24
#include <math.h>
25
#include <unistd.h>
26
#include <assert.h>
27

    
28
#include "mpiversiontest.h"
29

    
30
#include "misc.h"
31
#include "private_misc.h"
32
#include "debug.h"
33
#include "synchronize.h"
34
#include "accounting.h"
35

    
36
double *ping_pong_min_time; /* ping_pong_min_time[i] is the minimum time of one ping_pong
37
                               between the current node and node i, negative value means 
38
                               time not yet determined */
39
double *tds;                /* tds[i] is the time difference between the 
40
                               current node and global node i */
41

    
42

    
43
enum {
44
  Number_ping_pongs = 100,
45
  Minimum_ping_pongs = 8
46
};
47

    
48
bool mpi_wtime_is_global;
49

    
50
extern int wait_till(double time_stamp, double *last_time_stamp);
51
extern double should_wait_till(int counter, double interval, double offset);
52

    
53

    
54
void init_synchronization_module(void)
55
{
56
  int i, flag;
57
  int *mpi_wtime_is_global_ptr;
58

    
59
  tds = skampi_malloc_doubles(get_global_size());
60
  for(i = 0; i < get_global_size(); i++) tds[i] = 0.0;
61
  ping_pong_min_time = skampi_malloc_doubles(get_global_size());
62
  for( i = 0; i < get_global_size(); i++) ping_pong_min_time[i] = -1.0;
63

    
64
#if MPI_VERSION < 2
65
  MPI_Attr_get(MPI_COMM_WORLD, MPI_WTIME_IS_GLOBAL, &mpi_wtime_is_global_ptr, &flag);
66
#else
67
  MPI_Comm_get_attr(MPI_COMM_WORLD, MPI_WTIME_IS_GLOBAL, &mpi_wtime_is_global_ptr, &flag);
68
#endif
69
  if( flag == 0 ) 
70
    mpi_wtime_is_global = False;
71
  else
72
    mpi_wtime_is_global = *mpi_wtime_is_global_ptr;
73

    
74
  rdebug(DBG_TIMEDIFF, "mpi_wtime_is_global = %d\n", mpi_wtime_is_global);
75
}
76

    
77
void print_global_time_differences(void)
78
{
79
  int i, p, name_len, pid;
80
  char my_name[MPI_MAX_PROCESSOR_NAME];
81

    
82

    
83
  double *all_tds = NULL;
84
  char *names = NULL;
85
  int *pids = NULL;
86

    
87
  if( grootproc() ) {
88
    all_tds = skampi_malloc_doubles(get_global_size()*get_global_size());
89
    names = skampi_malloc_chars(get_global_size()*MPI_MAX_PROCESSOR_NAME);
90
    pids = skampi_malloc_ints(get_global_size());
91
  }
92

    
93
  MPI_Get_processor_name(my_name, &name_len);
94
  my_name[name_len] = '\0';
95
  pid = getpid();
96

    
97
  MPI_Gather(tds, get_global_size(), MPI_DOUBLE, 
98
             all_tds, get_global_size(), MPI_DOUBLE, 
99
             0, MPI_COMM_WORLD);
100
  MPI_Gather(my_name, MPI_MAX_PROCESSOR_NAME, MPI_CHAR,
101
             names, MPI_MAX_PROCESSOR_NAME, MPI_CHAR, 0, MPI_COMM_WORLD);
102
  MPI_Gather(&pid, 1, MPI_INT, pids, 1, MPI_INT, 0, MPI_COMM_WORLD);
103

    
104
  if( grootproc() ) {
105
    for( i = 0; i < get_global_size(); i++)
106
      printf("name[%d] = \"%s\" pid=%d\n", i, &(names[i*MPI_MAX_PROCESSOR_NAME]), pids[i]);
107

    
108
    for( p = 0; p < get_global_size(); p++) {
109
      printf("tds[%d -> ..] = [%9.6f", p, all_tds[p*get_global_size() + 0]);
110
      for( i = 1; i < get_global_size(); i++) 
111
        printf(", %9.6f", all_tds[p*get_global_size() + i]);
112
      printf("]\n");
113
    }
114

    
115
    free(all_tds);
116
    free(names);
117
    free(pids);
118
  }
119
}
120

    
121
/*----------------------------------------------------------------------------*/
122

    
123
static void ping_pong(int p1, int p2)
124
{
125
  int i, other_global_id;
126
  double my_time, my_last_time, other_time;
127
  double td_min, td_max;
128
  double invalid_time = -1.0;
129
  MPI_Status status;
130
  int pp_tag = 43;
131

    
132

    
133
  /* I had to unroll the main loop because I didn't find a portable way
134
     to define the initial td_min and td_max with INFINITY and NINFINITY */
135
  if( get_measurement_rank() == p1 ) {
136
    other_global_id = get_global_rank(p2);
137
    my_last_time = wtime();
138
    MPI_Send(&my_last_time, 1, MPI_DOUBLE, p2, pp_tag, get_measurement_comm());
139
    MPI_Recv(&other_time, 1, MPI_DOUBLE, p2, pp_tag, get_measurement_comm(), &status);
140
    my_time = wtime();
141
    td_min = other_time - my_time;
142
    td_max = other_time - my_last_time;
143
    MPI_Send(&my_time, 1, MPI_DOUBLE, p2, pp_tag, get_measurement_comm());
144
  } else {
145
    other_global_id = get_global_rank(p1);
146
    MPI_Recv(&other_time, 1, MPI_DOUBLE, p1, pp_tag, get_measurement_comm(), &status);
147
    my_last_time = wtime();
148
    td_min = other_time - my_last_time;
149
    MPI_Send(&my_last_time, 1, MPI_DOUBLE, p1, pp_tag, get_measurement_comm());
150
    MPI_Recv(&other_time, 1, MPI_DOUBLE, p1, pp_tag, get_measurement_comm(), &status);
151
    my_time = wtime();
152
    td_min = fmax2(td_min, other_time - my_time);
153
    td_max = other_time - my_last_time;
154
  }
155

    
156
  if( get_measurement_rank() == p1 ) {
157
    i = 1;
158
    while( 1 ) {
159
      MPI_Recv(&other_time, 1, MPI_DOUBLE, p2, pp_tag, get_measurement_comm(), &status);
160
      if( other_time < 0.0 ) break; 
161
      my_last_time = my_time; 
162
      my_time = wtime();
163
      td_min = fmax2(td_min, other_time - my_time);
164
      td_max = fmin2(td_max, other_time - my_last_time);
165
      if( ping_pong_min_time[other_global_id] >= 0.0  && 
166
          i >= Minimum_ping_pongs && 
167
          my_time - my_last_time < ping_pong_min_time[other_global_id]*1.10 ) {
168
        MPI_Send(&invalid_time, 1, MPI_DOUBLE, p2, pp_tag, get_measurement_comm());
169
        break;
170
      }
171
      i++;
172
      if( i == Number_ping_pongs ) {
173
        MPI_Send(&invalid_time, 1, MPI_DOUBLE, p2, pp_tag, get_measurement_comm());
174
        break;
175
      }
176
      MPI_Send(&my_time, 1, MPI_DOUBLE, p2, pp_tag, get_measurement_comm());
177

    
178
    }
179
  } else {
180
    i = 1;
181
    while( 1 ) {
182
      MPI_Send(&my_time, 1, MPI_DOUBLE, p1, pp_tag, get_measurement_comm());
183
      MPI_Recv(&other_time, 1, MPI_DOUBLE, p1, pp_tag, get_measurement_comm(), &status);
184
      if( other_time < 0.0 ) break; 
185
      my_last_time = my_time;
186
      my_time = wtime();
187
      td_min = fmax2(td_min, other_time - my_time);
188
      td_max = fmin2(td_max, other_time - my_last_time);
189
      if( ping_pong_min_time[other_global_id] >= 0.0 && 
190
          i >= Minimum_ping_pongs &&
191
          my_time - my_last_time < ping_pong_min_time[other_global_id]*1.10 ) {
192
        MPI_Send(&invalid_time, 1, MPI_DOUBLE, p1, pp_tag, get_measurement_comm());
193
        break;
194
      }
195
      i++;
196
    }
197
  }
198
  
199
  if( ping_pong_min_time[other_global_id] < 0.0) 
200
    ping_pong_min_time[other_global_id] = td_max-td_min;
201
  else 
202
    ping_pong_min_time[other_global_id] = 
203
      fmin2(ping_pong_min_time[other_global_id], td_max-td_min);
204

    
205
  tds[other_global_id] = (td_min+td_max)/2.0;
206
}
207

    
208
void determine_time_differences(void)
209
{
210
  int i;
211
  double *tmp_tds;
212
  double start_time, finish_time;
213

    
214
  if( lrootproc() ) 
215
    start_time = acc_start_sync();
216
  else
217
    start_time = MPI_Wtime();
218
  logging(DBG_SYNC, "determine_time_differences\n");
219
  for( i = 1; i < get_measurement_size(); i++) {
220
    MPI_Barrier(get_measurement_comm());
221
    if( get_measurement_rank() == 0 || get_measurement_rank() == i ) 
222
      ping_pong(0,i);
223
  }
224
  tmp_tds = skampi_malloc_doubles(get_measurement_size());
225
  if( lrootproc() ) {
226
    for( i = 1; i < get_measurement_size(); i++) 
227
            tmp_tds[i] = tds[get_global_rank(i)];
228
  }
229
  
230
  assert(get_measurement_size() - 1 >= 0 );
231
  MPI_Bcast(&(tmp_tds[1]), get_measurement_size()-1, MPI_DOUBLE, 0, get_measurement_comm());
232

    
233
  if( get_measurement_rank() != 0 ) {
234
    logging(DBG_TIMEDIFF, "tds[%d] = %9.1f us\n", 
235
            get_global_rank(0), tds[get_global_rank(0)]*1.0e6);
236
    for( i = 1; i < get_measurement_size(); i++) {
237
      tds[get_global_rank(i)] = tmp_tds[i] + tds[get_global_rank(0)];
238
      logging(DBG_TIMEDIFF, "tds[%d] = %9.1f us\n", 
239
              get_global_rank(i), tds[get_global_rank(i)]*1.0e6);
240
    }
241
  } 
242
  free(tmp_tds);
243
  MPI_Barrier(get_measurement_comm());
244
  if( lrootproc() )
245
    finish_time = acc_stop_sync();
246
  else
247
    finish_time = MPI_Wtime();
248
  rdebug(DBG_TIMEDIFF, "determine_time_differences() on %d processes: %f sec.\n", 
249
         get_measurement_size(), finish_time - start_time);
250
}