File: testGraph_Bcast.c

package info (click to toggle)
spooles 2.2-9
  • links: PTS, VCS
  • area: main
  • in suites: wheezy
  • size: 19,012 kB
  • sloc: ansic: 146,834; csh: 3,615; makefile: 2,040; perl: 74
file content (280 lines) | stat: -rw-r--r-- 8,499 bytes parent folder | download | duplicates (7)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
/*  testGraph_Bcast.c  */

#include "../spoolesMPI.h"
#include "../../Drand.h"
#include "../../timings.h"

/*--------------------------------------------------------------------*/
int
main ( int argc, char *argv[] )
/*
   --------------------------------------------------------------------
   this program tests the Graph_MPI_Bcast() method

   (1) process root generates a random Graph object
       and computes its checksum
   (2) process root broadcasts the Graph object to the other processors
   (3) each process computes the checksum of its Graph object
   (4) the checksums are compared on root

   created -- 98sep10, cca
   --------------------------------------------------------------------
*/
{
char         *buffer ;
double       chksum, t1, t2 ;
double       *sums ;
Drand        drand ;
int          iproc, length, loc, msglvl, myid, nitem, nproc, 
             nvtx, root, seed, size, type, v ;
int          *list ;
FILE         *msgFile ;
Graph        *graph ;
/*
   ---------------------------------------------------------------
   find out the identity of this process and the number of process
   ---------------------------------------------------------------
*/
MPI_Init(&argc, &argv) ;
MPI_Comm_rank(MPI_COMM_WORLD, &myid) ;
MPI_Comm_size(MPI_COMM_WORLD, &nproc) ;
fprintf(stdout, "\n process %d of %d, argc = %d", myid, nproc, argc) ;
fflush(stdout) ;
if ( argc != 8 ) {
   fprintf(stdout, 
           "\n\n usage : %s msglvl msgFile type nvtx nitem root seed "
           "\n    msglvl      -- message level"
           "\n    msgFile     -- message file"
           "\n    type        -- type of graph"
           "\n    nvtx        -- # of vertices"
           "\n    nitem       -- # of items used to generate graph"
           "\n    root        -- root processor for broadcast"
           "\n    seed        -- random number seed"
           "\n", argv[0]) ;
   return(0) ;
}
msglvl = atoi(argv[1]) ;
if ( strcmp(argv[2], "stdout") == 0 ) {
   msgFile = stdout ;
} else {
   length = strlen(argv[2]) + 1 + 4 ;
   buffer = CVinit(length, '\0') ;
   sprintf(buffer, "%s.%d", argv[2], myid) ;
   if ( (msgFile = fopen(buffer, "w")) == NULL ) {
      fprintf(stderr, "\n fatal error in %s"
              "\n unable to open file %s\n",
              argv[0], argv[2]) ;
      return(-1) ;
   }
   CVfree(buffer) ;
}
type  = atoi(argv[3]) ;
nvtx  = atoi(argv[4]) ;
nitem = atoi(argv[5]) ;
root  = atoi(argv[6]) ;
seed  = atoi(argv[7]) ;
fprintf(msgFile, 
        "\n %s "
        "\n msglvl  -- %d" 
        "\n msgFile -- %s" 
        "\n type    -- %d" 
        "\n nvtx    -- %d" 
        "\n nitem   -- %d" 
        "\n root    -- %d" 
        "\n seed    -- %d" 
        "\n",
        argv[0], msglvl, argv[2], type, nvtx, nitem, root, seed) ;
fflush(msgFile) ;
/*
   -----------------------
   set up the Graph object
   -----------------------
*/
MARKTIME(t1) ;
graph = Graph_new() ;
if ( myid == root ) {
   InpMtx   *inpmtx ;
   int      nedges, totewght, totvwght, v ;
   int      *adj, *vwghts ;
   IVL      *adjIVL, *ewghtIVL ;
/*
   -----------------------
   generate a random graph
   -----------------------
*/
   inpmtx = InpMtx_new() ;
   InpMtx_init(inpmtx, INPMTX_BY_ROWS, INPMTX_INDICES_ONLY, nitem, 0) ;
   Drand_setDefaultFields(&drand) ;
   Drand_setSeed(&drand, seed) ;
   Drand_setUniform(&drand, 0, nvtx) ;
   Drand_fillIvector(&drand, nitem, InpMtx_ivec1(inpmtx)) ;
   Drand_fillIvector(&drand, nitem, InpMtx_ivec2(inpmtx)) ;
   InpMtx_setNent(inpmtx, nitem) ;
   InpMtx_changeStorageMode(inpmtx, INPMTX_BY_VECTORS) ;
   if ( msglvl > 2 ) {
      fprintf(msgFile, "\n\n inpmtx mtx filled with raw entries") ;
      InpMtx_writeForHumanEye(inpmtx, msgFile) ;
      fflush(msgFile) ;
   }
   adjIVL = InpMtx_fullAdjacency(inpmtx) ;
   if ( msglvl > 2 ) {
      fprintf(msgFile, "\n\n full adjacency structure") ;
      IVL_writeForHumanEye(adjIVL, msgFile) ;
      fflush(msgFile) ;
   }
   nedges = adjIVL->tsize ;
   if ( type == 1 || type == 3 ) {
      Drand_setUniform(&drand, 1, 10) ;
      vwghts = IVinit(nvtx, 0) ;
      Drand_fillIvector(&drand, nvtx, vwghts) ;
      totvwght = IVsum(nvtx, vwghts) ;
      if ( msglvl > 2 ) {
         fprintf(msgFile, "\n\n vertex weights") ;
         IVfprintf(msgFile, nvtx, vwghts) ;
         fflush(msgFile) ;
      }
   } else {
      vwghts = NULL ;
      totvwght = nvtx ;
   }
   if ( msglvl > 2 ) {
      fprintf(msgFile, "\n\n totvwght %d", totvwght) ;
      fflush(msgFile) ;
   }
   if ( type == 2 || type == 3 ) {
      ewghtIVL = IVL_new() ;
      IVL_init1(ewghtIVL, IVL_CHUNKED, nvtx) ;
      Drand_setUniform(&drand, 1, 100) ;
      totewght = 0 ;
      for ( v = 0 ; v < nvtx ; v++ ) {
         IVL_listAndSize(adjIVL, v, &size, &adj) ;
         IVL_setList(ewghtIVL, v, size, NULL) ;
         IVL_listAndSize(ewghtIVL, v, &size, &adj) ;
         Drand_fillIvector(&drand, size, adj) ;
         totewght += IVsum(size, adj) ;
      }
      if ( msglvl > 2 ) {
         fprintf(msgFile, "\n\n ewghtIVL") ;
         IVL_writeForHumanEye(ewghtIVL, msgFile) ;
         fflush(msgFile) ;
      }
   } else {
      ewghtIVL = NULL ;
      totewght = nedges ;
   }
   if ( msglvl > 2 ) {
      fprintf(msgFile, "\n\n totewght %d", totewght) ;
      fflush(msgFile) ;
   }
   Graph_init2(graph, type, nvtx, 0, nedges, totvwght, totewght,
               adjIVL, vwghts, ewghtIVL) ;
   InpMtx_free(inpmtx) ;
}
MARKTIME(t2) ;
fprintf(msgFile, 
        "\n CPU %8.3f : initialize the Graph object", t2 - t1) ;
fflush(msgFile) ;
if ( msglvl > 2 ) {
   Graph_writeForHumanEye(graph, msgFile) ;
} else {
   Graph_writeStats(graph, msgFile) ;
}
fflush(msgFile) ;
if ( myid == root ) {
/*
   ----------------------------------------
   compute the checksum of the Graph object
   ----------------------------------------
*/
   chksum = graph->type + graph->nvtx + graph->nvbnd 
          + graph->nedges + graph->totvwght + graph->totewght ;
   for ( v = 0 ; v < nvtx ; v++ ) {
      IVL_listAndSize(graph->adjIVL, v, &size, &list) ;
      chksum += 1 + v + size + IVsum(size, list) ;
   }
   if ( graph->vwghts != NULL ) {
      chksum += IVsum(nvtx, graph->vwghts) ;
   }
   if ( graph->ewghtIVL != NULL ) {
      for ( v = 0 ; v < nvtx ; v++ ) {
         IVL_listAndSize(graph->ewghtIVL, v, &size, &list) ;
         chksum += 1 + v + size + IVsum(size, list) ;
      }
   }
   fprintf(msgFile, "\n\n local chksum = %12.4e", chksum) ;
   fflush(msgFile) ;
}
/*
   --------------------------
   broadcast the Graph object
   --------------------------
*/
MARKTIME(t1) ;
graph = Graph_MPI_Bcast(graph, root, msglvl, msgFile, MPI_COMM_WORLD) ;
MARKTIME(t2) ;
fprintf(msgFile, "\n CPU %8.3f : broadcast the Graph object", t2 - t1) ;
if ( msglvl > 2 ) {
   Graph_writeForHumanEye(graph, msgFile) ;
} else {
   Graph_writeStats(graph, msgFile) ;
}
/*
   ----------------------------------------
   compute the checksum of the Graph object
   ----------------------------------------
*/
chksum = graph->type + graph->nvtx + graph->nvbnd 
       + graph->nedges + graph->totvwght + graph->totewght ;
for ( v = 0 ; v < nvtx ; v++ ) {
   IVL_listAndSize(graph->adjIVL, v, &size, &list) ;
   chksum += 1 + v + size + IVsum(size, list) ;
}
if ( graph->vwghts != NULL ) {
   chksum += IVsum(nvtx, graph->vwghts) ;
}
if ( graph->ewghtIVL != NULL ) {
   for ( v = 0 ; v < nvtx ; v++ ) {
      IVL_listAndSize(graph->ewghtIVL, v, &size, &list) ;
      chksum += 1 + v + size + IVsum(size, list) ;
   }
}
fprintf(msgFile, "\n\n local chksum = %12.4e", chksum) ;
fflush(msgFile) ;
/*
   ---------------------------------------
   gather the checksums from the processes
   ---------------------------------------
*/
sums = DVinit(nproc, 0.0) ;
MPI_Gather((void *) &chksum, 1, MPI_DOUBLE, 
           (void *) sums, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD) ;
if ( myid == 0 ) {
   fprintf(msgFile, "\n\n sums") ;
   DVfprintf(msgFile, nproc, sums) ;
   for ( iproc = 0 ; iproc < nproc ; iproc++ ) {
      sums[iproc] -= chksum ;
   }
   fprintf(msgFile, "\n\n errors") ;
   DVfprintf(msgFile, nproc, sums) ;
   fprintf(msgFile, "\n\n maxerror = %12.4e", DVmax(nproc, sums, &loc));
}
/*
   ----------------
   free the objects
   ----------------
*/
DVfree(sums) ;
Graph_free(graph) ;
/*
   ------------------------
   exit the MPI environment
   ------------------------
*/
MPI_Finalize() ;

fprintf(msgFile, "\n") ;
fclose(msgFile) ;

return(0) ; }

/*--------------------------------------------------------------------*/