~chaffra/+junk/trilinos

« back to all changes in this revision

Viewing changes to packages/zoltan/example/simple/simpleBLOCK.c

  • Committer: Bazaar Package Importer
  • Author(s): Christophe Prud'homme, Christophe Prud'homme, Johannes Ring
  • Date: 2009-12-13 12:53:22 UTC
  • mfrom: (5.1.2 sid)
  • Revision ID: james.westby@ubuntu.com-20091213125322-in0nrdjc55deqsw9
Tags: 10.0.3.dfsg-1
[Christophe Prud'homme]
* New upstream release

[Johannes Ring]
* debian/patches/libname.patch: Add prefix 'libtrilinos_' to all
  libraries. 
* debian/patches/soname.patch: Add soversion to libraries.
* debian/watch: Update download URL.
* debian/control:
  - Remove python-numeric from Build-Depends (virtual package).
  - Remove automake and autotools from Build-Depends and add cmake to
    reflect switch to CMake.
  - Add python-support to Build-Depends.
* debian/rules: 
  - Cleanup and updates for switch to CMake.

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
/* Basic example of using Zoltan to compute a quick partitioning
 
2
** of a set of objects.
 
3
***************************************************************/
 
4
 
 
5
#include <mpi.h>
 
6
#include <stdio.h>
 
7
#include <stdlib.h>
 
8
#include <ctype.h>
 
9
#include "zoltan.h"
 
10
 
 
11
/* Name of file containing objects to be partitioned */
 
12
 
 
13
static char *fname="objects.txt";
 
14
 
 
15
/* Structure to hold object data */
 
16
 
 
17
typedef struct{
 
18
  int numGlobalObjects;
 
19
  int numMyObjects;
 
20
  int *myGlobalIDs;
 
21
} OBJECT_DATA;
 
22
 
 
23
static int get_number_of_objects(void *data, int *ierr);
 
24
 
 
25
static void get_object_list(void *data, int sizeGID, int sizeLID,
 
26
            ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
 
27
                  int wgt_dim, float *obj_wgts, int *ierr);
 
28
 
 
29
static int get_next_line(FILE *fp, char *buf, int bufsize);
 
30
 
 
31
static void input_file_error(int numProcs, int tag, int startProc);
 
32
 
 
33
static void showSimpleMeshPartitions(int myProc, int numIDs, int *GIDs, int *parts);
 
34
 
 
35
static void read_input_objects(int myRank, int numProcs, char *fname, OBJECT_DATA *myData);
 
36
 
 
37
int main(int argc, char *argv[])
 
38
{
 
39
  int rc, i;
 
40
  int myRank, numProcs;
 
41
  float ver;
 
42
  struct Zoltan_Struct *zz;
 
43
  int changes, numGidEntries, numLidEntries, numImport, numExport;
 
44
  ZOLTAN_ID_PTR importGlobalGids, importLocalGids;
 
45
  ZOLTAN_ID_PTR exportGlobalGids, exportLocalGids; 
 
46
  int *importProcs, *importToPart, *exportProcs, *exportToPart;
 
47
  int *parts = NULL;
 
48
 
 
49
  FILE *fp;
 
50
  OBJECT_DATA myData;
 
51
 
 
52
  /******************************************************************
 
53
  ** Initialize MPI and Zoltan
 
54
  ******************************************************************/
 
55
 
 
56
  MPI_Init(&argc, &argv);
 
57
  MPI_Comm_rank(MPI_COMM_WORLD, &myRank);
 
58
  MPI_Comm_size(MPI_COMM_WORLD, &numProcs);
 
59
 
 
60
  rc = Zoltan_Initialize(argc, argv, &ver);
 
61
 
 
62
  if (rc != ZOLTAN_OK){
 
63
    printf("Error initializing Zoltan\n");
 
64
    MPI_Finalize();
 
65
    exit(0);
 
66
  }
 
67
 
 
68
  /******************************************************************
 
69
  ** Read objects from input file and distribute them unevenly
 
70
  ******************************************************************/
 
71
 
 
72
  fp = fopen(fname, "r");
 
73
  if (!fp){
 
74
    if (myRank == 0) fprintf(stderr,"ERROR: Can not open %s\n",fname);
 
75
    MPI_Finalize();
 
76
    exit(1);
 
77
  }
 
78
  fclose(fp);
 
79
 
 
80
  read_input_objects(myRank, numProcs, fname, &myData);
 
81
 
 
82
  /******************************************************************
 
83
  ** Create a Zoltan library structure for this instance of load
 
84
  ** balancing.  Set the parameters and query functions.
 
85
  ******************************************************************/
 
86
 
 
87
  zz = Zoltan_Create(MPI_COMM_WORLD);
 
88
 
 
89
  /* General parameters */
 
90
 
 
91
  Zoltan_Set_Param(zz, "LB_METHOD", "BLOCK");  /* Zoltan method: "BLOCK" */
 
92
  Zoltan_Set_Param(zz, "NUM_GID_ENTRIES", "1"); /* global ID is 1 integer */
 
93
  Zoltan_Set_Param(zz, "NUM_LID_ENTRIES", "1"); /* local ID is 1 integer */
 
94
  Zoltan_Set_Param(zz, "OBJ_WEIGHT_DIM", "0"); /* we omit object weights */
 
95
 
 
96
  /* Query functions */
 
97
 
 
98
  Zoltan_Set_Num_Obj_Fn(zz, get_number_of_objects, &myData);
 
99
  Zoltan_Set_Obj_List_Fn(zz, get_object_list, &myData);
 
100
 
 
101
  /******************************************************************
 
102
  ** Call Zoltan to partition the objects.
 
103
  ******************************************************************/
 
104
 
 
105
  rc = Zoltan_LB_Partition(zz, /* input (all remaining fields are output) */
 
106
        &changes,        /* 1 if partitioning was changed, 0 otherwise */ 
 
107
        &numGidEntries,  /* Number of integers used for a global ID */
 
108
        &numLidEntries,  /* Number of integers used for a local ID */
 
109
        &numImport,      /* Number of objects to be sent to me */
 
110
        &importGlobalGids,  /* Global IDs of objects to be sent to me */
 
111
        &importLocalGids,   /* Local IDs of objects to be sent to me */
 
112
        &importProcs,    /* Process rank for source of each incoming object */
 
113
        &importToPart,   /* New partition for each incoming object */
 
114
        &numExport,      /* Number of objects I must send to other processes*/
 
115
        &exportGlobalGids,  /* Global IDs of the objects I must send */
 
116
        &exportLocalGids,   /* Local IDs of the objects I must send */
 
117
        &exportProcs,    /* Process to which I send each of the objects */
 
118
        &exportToPart);  /* Partition to which each object will belong */
 
119
 
 
120
  if (rc != ZOLTAN_OK){
 
121
    printf("Error in Zoltan library\n");
 
122
    MPI_Finalize();
 
123
    Zoltan_Destroy(&zz);
 
124
    exit(0);
 
125
  }
 
126
 
 
127
  /******************************************************************
 
128
  ** Visualize the object partitioning before and after calling Zoltan.
 
129
  **
 
130
  ** In this example, partition number equals process rank.
 
131
  ******************************************************************/
 
132
 
 
133
  parts = (int *)malloc(sizeof(int) * myData.numMyObjects);
 
134
 
 
135
  for (i=0; i < myData.numMyObjects; i++){
 
136
    parts[i] = myRank;
 
137
  }
 
138
 
 
139
  if (myRank== 0){
 
140
    printf("\nObject partition assignments before calling Zoltan\n");
 
141
  }
 
142
 
 
143
  showSimpleMeshPartitions(myRank, myData.numMyObjects, myData.myGlobalIDs, parts);
 
144
 
 
145
  for (i=0; i < numExport; i++){
 
146
    parts[exportLocalGids[i]] = exportToPart[i];
 
147
  }
 
148
 
 
149
  if (myRank == 0){
 
150
    printf("Object partition assignments after calling Zoltan\n");
 
151
  }
 
152
 
 
153
  showSimpleMeshPartitions(myRank, myData.numMyObjects, myData.myGlobalIDs, parts);
 
154
 
 
155
  /******************************************************************
 
156
  ** Free the arrays allocated by Zoltan_LB_Partition, and free
 
157
  ** the storage allocated for the Zoltan structure.
 
158
  ******************************************************************/
 
159
 
 
160
  Zoltan_LB_Free_Part(&importGlobalGids, &importLocalGids, 
 
161
                      &importProcs, &importToPart);
 
162
  Zoltan_LB_Free_Part(&exportGlobalGids, &exportLocalGids, 
 
163
                      &exportProcs, &exportToPart);
 
164
 
 
165
  Zoltan_Destroy(&zz);
 
166
 
 
167
  MPI_Finalize();
 
168
 
 
169
  return 0;
 
170
}
 
171
/* Application defined query functions */
 
172
 
 
173
static int get_number_of_objects(void *data, int *ierr)
 
174
{
 
175
  *ierr = ZOLTAN_OK;
 
176
  OBJECT_DATA *objects = (OBJECT_DATA *)data;
 
177
  return objects->numMyObjects;
 
178
}
 
179
 
 
180
static void get_object_list(void *data, int sizeGID, int sizeLID,
 
181
            ZOLTAN_ID_PTR globalID, ZOLTAN_ID_PTR localID,
 
182
                  int wgt_dim, float *obj_wgts, int *ierr)
 
183
{
 
184
int i;
 
185
  *ierr = ZOLTAN_OK;
 
186
  OBJECT_DATA *objects = (OBJECT_DATA *)data;
 
187
 
 
188
  /* In this example, return the IDs of our objects, but no weights.
 
189
   * Zoltan will assume equally weighted objects.
 
190
   */
 
191
 
 
192
  for (i=0; i<objects->numMyObjects; i++){
 
193
    globalID[i] = objects->myGlobalIDs[i];
 
194
    localID[i] = i;
 
195
  }
 
196
}
 
197
 
 
198
/* Function to find next line of information in input file */
 
199
 
 
200
static int get_next_line(FILE *fp, char *buf, int bufsize)
 
201
{
 
202
int i, cval, len;
 
203
char *c;
 
204
 
 
205
  while (1){
 
206
 
 
207
    c = fgets(buf, bufsize, fp);
 
208
 
 
209
    if (c == NULL)
 
210
      return 0;  /* end of file */
 
211
 
 
212
    len = strlen(c);
 
213
 
 
214
    for (i=0, c=buf; i < len; i++, c++){
 
215
      cval = (int)*c; 
 
216
      if (isspace(cval) == 0) break;
 
217
    }
 
218
    if (i == len) continue;   /* blank line */
 
219
    if (*c == '#') continue;  /* comment */
 
220
 
 
221
    if (c != buf){
 
222
      strcpy(buf, c);
 
223
    }
 
224
    break;
 
225
  }
 
226
 
 
227
  return strlen(buf);  /* number of characters */
 
228
}
 
229
 
 
230
/* Proc 0 notifies others of error and exits */
 
231
 
 
232
static void input_file_error(int numProcs, int tag, int startProc)
 
233
{
 
234
int i, val;
 
235
 
 
236
  val = -1;
 
237
 
 
238
  fprintf(stderr,"ERROR in input file.\n");
 
239
 
 
240
  for (i=startProc; i < numProcs; i++){
 
241
    /* these procs have posted receive for "tag" */
 
242
    MPI_Send(&val, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
 
243
  }
 
244
  for (i=1; i < startProc; i++){
 
245
    /* these procs are done */
 
246
    MPI_Send(&val, 1, MPI_INT, i, 0, MPI_COMM_WORLD);
 
247
  }
 
248
 
 
249
  MPI_Finalize();
 
250
  exit(1);
 
251
}
 
252
 
 
253
/* Draw the partition assignments of the objects */
 
254
 
 
255
void showSimpleMeshPartitions(int myProc, int numIDs, int *GIDs, int *parts)
 
256
{
 
257
int partAssign[25], allPartAssign[25];
 
258
int i, j, part;
 
259
 
 
260
  memset(partAssign, 0, sizeof(int) * 25);
 
261
 
 
262
  for (i=0; i < numIDs; i++){
 
263
    partAssign[GIDs[i]-1] = parts[i];
 
264
  }
 
265
 
 
266
  MPI_Reduce(partAssign, allPartAssign, 25, MPI_INT, MPI_MAX, 0, MPI_COMM_WORLD);
 
267
 
 
268
  if (myProc == 0){
 
269
 
 
270
    for (i=20; i >= 0; i-=5){
 
271
      for (j=0; j < 5; j++){
 
272
        part = allPartAssign[i + j];
 
273
        if (j < 4)
 
274
          printf("%d-----",part);
 
275
        else
 
276
          printf("%d\n",part);
 
277
      }
 
278
      if (i > 0)
 
279
        printf("|     |     |     |     |\n");
 
280
    }
 
281
    printf("\n");
 
282
  }
 
283
}
 
284
 
 
285
/* Proc 0 reads the objects in the input file and divides them across processes */
 
286
 
 
287
void read_input_objects(int myRank, int numProcs, char *fname, OBJECT_DATA *myData)
 
288
{
 
289
char *buf;
 
290
int bufsize = 512;
 
291
int num, nobj, remainingObj, ack=0;
 
292
int i, j;
 
293
int *gids;
 
294
FILE *fp;
 
295
MPI_Status status;
 
296
int obj_ack_tag = 5, obj_count_tag = 10, obj_id_tag = 15;
 
297
 
 
298
  if (myRank == 0){
 
299
 
 
300
    buf = (char *)malloc(sizeof(char) * bufsize);
 
301
    fp = fopen(fname, "r");
 
302
 
 
303
    num = get_next_line(fp, buf, bufsize);
 
304
    if (num == 0) input_file_error(numProcs, obj_count_tag, 1);
 
305
    num = sscanf(buf, "%d", &myData->numGlobalObjects);
 
306
    if (num != 1) input_file_error(numProcs, obj_count_tag, 1);
 
307
 
 
308
    if (numProcs > 1){
 
309
      nobj = myData->numGlobalObjects / 2;
 
310
      remainingObj = myData->numGlobalObjects - nobj;
 
311
    }
 
312
    else{
 
313
      nobj = myData->numGlobalObjects;
 
314
      remainingObj = 0;
 
315
    }
 
316
 
 
317
    myData->myGlobalIDs = (int *)malloc(sizeof(int) * nobj);
 
318
    myData->numMyObjects = nobj;
 
319
 
 
320
    for (i=0; i < nobj; i++){
 
321
 
 
322
      num = get_next_line(fp, buf, bufsize);
 
323
      if (num == 0) input_file_error(numProcs, obj_count_tag, 1);
 
324
      num = sscanf(buf, "%d", myData->myGlobalIDs + i);
 
325
      if (num != 1) input_file_error(numProcs, obj_count_tag, 1);
 
326
  
 
327
    }
 
328
 
 
329
    gids = (int *)malloc(sizeof(int) * (nobj + 1));
 
330
 
 
331
    for (i=1; i < numProcs; i++){
 
332
    
 
333
      if (remainingObj > 1){
 
334
        nobj = remainingObj / 2;
 
335
        remainingObj -= nobj;
 
336
      }
 
337
      else if (remainingObj == 1){
 
338
        nobj = 1;
 
339
        remainingObj = 0;
 
340
      }
 
341
      else{
 
342
        nobj = 0;
 
343
      }
 
344
 
 
345
      if ((i == numProcs - 1) && (remainingObj > 0))
 
346
        nobj += remainingObj;
 
347
 
 
348
      if (nobj > 0){
 
349
        for (j=0; j < nobj; j++){
 
350
          num = get_next_line(fp, buf, bufsize);
 
351
          if (num == 0) input_file_error(numProcs, obj_count_tag, i);
 
352
          num = sscanf(buf, "%d", gids + j);
 
353
          if (num != 1) input_file_error(numProcs, obj_count_tag, i);
 
354
        }
 
355
      }
 
356
 
 
357
      MPI_Send(&nobj, 1, MPI_INT, i, obj_count_tag, MPI_COMM_WORLD);
 
358
      MPI_Recv(&ack, 1, MPI_INT, i, obj_ack_tag, MPI_COMM_WORLD, &status);
 
359
 
 
360
      if (nobj > 0)
 
361
        MPI_Send(gids, nobj, MPI_INT, i, obj_id_tag, MPI_COMM_WORLD);
 
362
      
 
363
    }
 
364
 
 
365
    free(gids);
 
366
    fclose(fp);
 
367
    free(buf);
 
368
 
 
369
    /* signal all procs it is OK to go on */
 
370
    ack = 0;
 
371
    for (i=1; i < numProcs; i++){
 
372
      MPI_Send(&ack, 1, MPI_INT, i, 0, MPI_COMM_WORLD);
 
373
    }
 
374
  }
 
375
  else{
 
376
 
 
377
    MPI_Recv(&myData->numMyObjects, 1, MPI_INT, 0, obj_count_tag, MPI_COMM_WORLD, &status);
 
378
    ack = 0;
 
379
    if (myData->numMyObjects > 0){
 
380
      myData->myGlobalIDs = (int *)malloc(sizeof(int) * myData->numMyObjects);
 
381
      MPI_Send(&ack, 1, MPI_INT, 0, obj_ack_tag, MPI_COMM_WORLD);
 
382
      MPI_Recv(myData->myGlobalIDs, myData->numMyObjects, MPI_INT, 0, 
 
383
               obj_id_tag, MPI_COMM_WORLD, &status);
 
384
    }
 
385
    else if (myData->numMyObjects == 0){
 
386
      MPI_Send(&ack, 1, MPI_INT, 0, obj_ack_tag, MPI_COMM_WORLD);
 
387
    }
 
388
    else{
 
389
      MPI_Finalize();
 
390
      exit(1);
 
391
    }
 
392
 
 
393
    MPI_Recv(&ack, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, &status);
 
394
    if (ack < 0){
 
395
      MPI_Finalize();
 
396
      exit(1);
 
397
    }
 
398
  }
 
399
}