ScaLAPACK 2.1  2.1
ScaLAPACK: Scalable Linear Algebra PACKage
pdgemrdrv.c
Go to the documentation of this file.
1 #include "redist.h"
2 /* $Id: pdgemrdrv.c,v 1.1.1.1 2000/02/15 18:04:11 susan Exp $
3  *
4  * pdgemrdrv.c :
5  *
6  *
7  * PURPOSE:
8  *
9  * this driver is testing the PDGEMR2D routine. It calls it to obtain a new
10  * scattered block data decomposition of a distributed DOUBLE PRECISION
11  * (block scattered) matrix. Then it calls PDGEMR2D for the inverse
12  * redistribution and checks the results with the initial data.
13  *
14  * Data are going from a Block Scattered nbrow0 x nbcol0 decomposition on the
15  * processor grid p0 x q0, to data distributed in a BS nbrow1 x nbcol1 on the
16  * processor grid p1 x q1, then back to the BS nbrow0 x nbcol0 decomposition
17  * on the processor grid p0 x q0.
18  *
19  * See pdgemr.c file for detailed info on the PDGEMR2D function.
20  *
21  *
22  * The testing parameters are read from the file GEMR2D.dat, see the file in the
23  * distribution to have an example.
24  *
25  * created by Bernard Tourancheau in April 1994.
26  *
27  * modifications : see sccs history
28  *
29  * ===================================
30  *
31  *
32  * NOTE :
33  *
34  * - the matrix elements are DOUBLE PRECISION
35  *
36  * - memory requirements : this procedure requires approximately 3 times the
37  * memory space of the initial data block in grid 0 (initial block, copy for
38  * test and second redistribution result) and 1 time the memory space of the
39  * result data block in grid 1. with the element size = sizeof(double)
40  * bytes,
41  *
42  *
43  * - use the procedures of the files:
44  *
45  * pdgemr.o pdgemr2.o pdgemraux.o
46  *
47  *
48  * ======================================
49  *
50  * WARNING ASSUMPTIONS :
51  *
52  *
53  * ========================================
54  *
55  *
56  * Planned changes:
57  *
58  *
59  *
60  * ========================================= */
61 #define static2 static
62 #if defined(Add_) || defined(f77IsF2C)
63 #define fortran_mr2d pdgemr2do_
64 #define fortran_mr2dnew pdgemr2d_
65 #elif defined(UpCase)
66 #define fortran_mr2dnew PDGEMR2D
67 #define fortran_mr2d PDGEMR2DO
68 #define dcopy_ DCOPY
69 #define dlacpy_ DLACPY
70 #else
71 #define fortran_mr2d pdgemr2do
72 #define fortran_mr2dnew pdgemr2d
73 #define dcopy_ dcopy
74 #define dlacpy_ dlacpy
75 #endif
76 #define Clacpy Cdgelacpy
77 void Clacpy();
78 typedef struct {
79  int desctype;
80  int ctxt;
81  int m;
82  int n;
83  int nbrow;
84  int nbcol;
85  int sprow;
86  int spcol;
87  int lda;
88 } MDESC;
89 #define BLOCK_CYCLIC_2D 1
90 typedef struct {
91  int lstart;
92  int len;
93 } IDESC;
94 #define SHIFT(row,sprow,nbrow) ((row)-(sprow)+ ((row) >= (sprow) ? 0 : (nbrow)))
95 #define max(A,B) ((A)>(B)?(A):(B))
96 #define min(A,B) ((A)>(B)?(B):(A))
97 #define DIVUP(a,b) ( ((a)-1) /(b)+1)
98 #define ROUNDUP(a,b) (DIVUP(a,b)*(b))
99 #ifdef MALLOCDEBUG
100 #define malloc mymalloc
101 #define free myfree
102 #define realloc myrealloc
103 #endif
104 /* Cblacs */
105 extern void Cblacs_pcoord();
106 extern int Cblacs_pnum();
107 extern void Csetpvmtids();
108 extern void Cblacs_get();
109 extern void Cblacs_pinfo();
110 extern void Cblacs_gridinfo();
111 extern void Cblacs_gridinit();
112 extern void Cblacs_exit();
113 extern void Cblacs_gridexit();
114 extern void Cblacs_setup();
115 extern void Cigebs2d();
116 extern void Cigebr2d();
117 extern void Cigesd2d();
118 extern void Cigerv2d();
119 extern void Cigsum2d();
120 extern void Cigamn2d();
121 extern void Cigamx2d();
122 extern void Cdgesd2d();
123 extern void Cdgerv2d();
124 /* lapack */
125 void dlacpy_();
126 /* aux fonctions */
127 extern int localindice();
128 extern void *mr2d_malloc();
129 extern int ppcm();
130 extern int localsize();
131 extern int memoryblocksize();
132 extern int changeorigin();
133 extern void paramcheck();
134 /* tools and others function */
135 #define scanD0 dgescanD0
136 #define dispmat dgedispmat
137 #define setmemory dgesetmemory
138 #define freememory dgefreememory
139 #define scan_intervals dgescan_intervals
140 extern void scanD0();
141 extern void dispmat();
142 extern void setmemory();
143 extern void freememory();
144 extern int scan_intervals();
145 extern void Cpdgemr2do();
146 extern void Cpdgemr2d();
147 /* some defines for Cpdgemr2do */
148 #define SENDBUFF 0
149 #define RECVBUFF 1
150 #define SIZEBUFF 2
151 #if 0
152 #define DEBUG
153 #endif
154 #ifndef DEBUG
155 #define NDEBUG
156 #endif
157 #include <stdio.h>
158 #include <stdlib.h>
159 #include <string.h>
160 #include <ctype.h>
161 #include <assert.h>
162 /* initblock: intialize the local part of a matrix with random data (well,
163  * not very random) */
164 static2 void
165 initblock(block, m, n)
166  double *block;
167  int m, n;
168 {
169  double *pdata;
170  int i;
171  pdata = block;
172  for (i = 0; i < m * n; i++, pdata++) {
173  (*pdata) = i;
174  };
175 }
176 /* getparam:read from a file a list of integer parameters, the end of the
177  * parameters to read is given by a NULL at the end of the args list */
178 #ifdef __STDC__
179 #include <stdarg.h>
180 static void
181 getparam(FILE * f,...)
182 {
183 #else
184 #include <varargs.h>
185 static void
186 getparam(va_alist)
187 va_dcl
188 {
189  FILE *f;
190 #endif
191  va_list ap;
192  int i;
193  static int nbline;
194  char *ptr, *next;
195  int *var;
196  static char buffer[200];
197 #ifdef __STDC__
198  va_start(ap, f);
199 #else
200  va_start(ap);
201  f = va_arg(ap, FILE *);
202 #endif
203  do {
204  next = fgets(buffer, 200, f);
205  if (next == NULL) {
206  fprintf(stderr, "bad configuration driver file:after line %d\n", nbline);
207  exit(1);
208  }
209  nbline += 1;
210  } while (buffer[0] == '#');
211  ptr = buffer;
212  var = va_arg(ap, int *);
213  while (var != NULL) {
214  *var = strtol(ptr, &next, 10);
215  if (ptr == next) {
216  fprintf(stderr, "bad configuration driver file:error line %d\n", nbline);
217  exit(1);
218  }
219  ptr = next;
220  var = va_arg(ap, int *);
221  }
222  va_end(ap);
223 }
224 void
225 initforpvm(argc, argv)
226  int argc;
227  char *argv[];
228 {
229  int pnum, nproc;
230  Cblacs_pinfo(&pnum, &nproc);
231  if (nproc < 1) { /* we are with PVM */
232  if (pnum == 0) {
233  if (argc < 2) {
234  fprintf(stderr, "usage with PVM:xdgemr nbproc\n\
235 \t where nbproc is the number of nodes to initialize\n");
236  exit(1);
237  }
238  nproc = atoi(argv[1]);
239  }
240  Cblacs_setup(&pnum, &nproc);
241  }
242 }
243 int
244 main(argc, argv)
245  int argc;
246  char *argv[];
247 {
248  /* We initialize the data-block on the current processor, then redistribute
249  * it, and perform the inverse redistribution to compare the local memory
250  * with the initial one. */
251  /* Data file */
252  FILE *fp;
253  int nbre, nbremax;
254  /* Data distribution 0 parameters */
255  int p0, /* # of rows in the processor grid */
256  q0; /* # of columns in the processor grid */
257  /* Data distribution 1 parameters */
258  int p1, q1;
259  /* # of parameter to be read on the keyboard */
260 #define nbparameter 24
261  /* General variables */
262  int blocksize0;
263  int mypnum, nprocs;
264  int parameters[nbparameter], nberrors;
265  int i;
266  int ia, ja, ib, jb, m, n;
267  int gcontext, context0, context1;
268  int myprow1, myprow0, mypcol0, mypcol1;
269  int dummy;
270  MDESC ma, mb;
271  double *ptrmyblock, *ptrsavemyblock, *ptrmyblockcopy, *ptrmyblockvide;
272 #ifdef UsingMpiBlacs
273  MPI_Init(&argc, &argv);
274 #endif
275  setvbuf(stdout, NULL, _IOLBF, 0);
276  setvbuf(stderr, NULL, _IOLBF, 0);
277 #ifdef T3D
278  free(malloc(14000000));
279 #endif
280  initforpvm(argc, argv);
281  /* Read physical parameters */
282  Cblacs_pinfo(&mypnum, &nprocs);
283  /* initialize BLACS for the parameter communication */
284  Cblacs_get(0, 0, &gcontext);
285  Cblacs_gridinit(&gcontext, "R", nprocs, 1);
286  Cblacs_gridinfo(gcontext, &dummy, &dummy, &mypnum, &dummy);
287  if (mypnum == 0) {
288  if ((fp = fopen("GEMR2D.dat", "r")) == NULL) {
289  fprintf(stderr, "Can't open GEMR2D.dat\n");
290  exit(1);
291  };
292  printf("\n// DGEMR2D TESTER for DOUBLE PRECISION //\n");
293  getparam(fp, &nbre, NULL);
294  printf("////////// %d tests \n\n", nbre);
295  parameters[0] = nbre;
296  Cigebs2d(gcontext, "All", "H", 1, 1, parameters, 1);
297  } else {
298  Cigebr2d(gcontext, "All", "H", 1, 1, parameters, 1, 0, 0);
299  nbre = parameters[0];
300  };
301  if (mypnum == 0) {
302  printf("\n m n m0 n0 sr0 sc0 i0 j0 p0 q0 nbr0 nbc0 \
303 m1 n1 sr1 sc1 i1 j1 p1 q1 nbr1 nbc1\n\n");
304  };
305  /****** TEST LOOP *****/
306  /* Here we are in grip 1xnprocs */
307  nbremax = nbre;
308 #ifdef DEBUG
309  fprintf(stderr, "bonjour,je suis le noeud %d\n", mypnum);
310 #endif
311  while (nbre-- != 0) { /* Loop on the serie of tests */
312  /* All the processors read the parameters so we have to be in a 1xnprocs
313  * grid at each iteration */
314  /* Read processors grid and matrices parameters */
315  if (mypnum == 0) {
316  int u, d;
317  getparam(fp,
318  &m, &n,
319  &ma.m, &ma.n, &ma.sprow, &ma.spcol,
320  &ia, &ja, &p0, &q0, &ma.nbrow, &ma.nbcol,
321  &mb.m, &mb.n, &mb.sprow, &mb.spcol,
322  &ib, &jb, &p1, &q1, &mb.nbrow, &mb.nbcol,
323  NULL);
324  printf("\t\t************* TEST # %d **********\n",
325  nbremax - nbre);
326  printf(" %3d %3d %3d %3d %3d %3d %3d %3d \
327 %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d %3d",
328  m, n,
329  ma.m, ma.n, ma.sprow, ma.spcol,
330  ia, ja, p0, q0, ma.nbrow, ma.nbcol,
331  mb.m, mb.n, mb.sprow, mb.spcol,
332  ib, jb, p1, q1, mb.nbrow, mb.nbcol);
333  printf("\n");
334  if (p0 * q0 > nprocs || p1 * q1 > nprocs) {
335  fprintf(stderr, "not enough nodes:%d processors required\n",
336  max(p0 * q0, p1 * q1));
337  exit(1);
338  }
339  parameters[0] = p0;
340  parameters[1] = q0;
341  parameters[2] = ma.nbrow;
342  parameters[3] = ma.nbcol;
343  parameters[4] = p1;
344  parameters[5] = q1;
345  parameters[6] = mb.nbrow;
346  parameters[7] = mb.nbcol;
347  parameters[8] = ma.m;
348  parameters[9] = ma.n;
349  parameters[10] = ma.sprow;
350  parameters[11] = ma.spcol;
351  parameters[12] = mb.sprow;
352  parameters[13] = mb.spcol;
353  parameters[14] = ia;
354  parameters[15] = ja;
355  parameters[16] = ib;
356  parameters[17] = jb;
357  parameters[18] = m;
358  parameters[19] = n;
359  parameters[20] = mb.m;
360  parameters[21] = mb.n;
361  Cigebs2d(gcontext, "All", "H", 1, nbparameter, parameters, 1);
362  } else {
363  Cigebr2d(gcontext, "All", "H", 1, nbparameter, parameters, 1, 0, 0);
364  p0 = parameters[0];
365  q0 = parameters[1];
366  ma.nbrow = parameters[2];
367  ma.nbcol = parameters[3];
368  p1 = parameters[4];
369  q1 = parameters[5];
370  mb.nbrow = parameters[6];
371  mb.nbcol = parameters[7];
372  ma.m = parameters[8];
373  ma.n = parameters[9];
374  ma.sprow = parameters[10];
375  ma.spcol = parameters[11];
376  mb.sprow = parameters[12];
377  mb.spcol = parameters[13];
378  ia = parameters[14];
379  ja = parameters[15];
380  ib = parameters[16];
381  jb = parameters[17];
382  m = parameters[18];
383  n = parameters[19];
384  mb.m = parameters[20];
385  mb.n = parameters[21];
388  };
389  Cblacs_get(0, 0, &context0);
390  Cblacs_gridinit(&context0, "R", p0, q0);
391  Cblacs_get(0, 0, &context1);
392  Cblacs_gridinit(&context1, "R", p1, q1);
393  Cblacs_gridinfo(context0, &dummy, &dummy, &myprow0, &mypcol0);
394  if (myprow0 >= p0 || mypcol0 >= q0)
395  myprow0 = mypcol0 = -1;
396  Cblacs_gridinfo(context1, &dummy, &dummy, &myprow1, &mypcol1);
397  if (myprow1 >= p1 || mypcol1 >= q1)
398  myprow1 = mypcol1 = -1;
399  assert((myprow0 < p0 && mypcol0 < q0) || (myprow0 == -1 && mypcol0 == -1));
400  assert((myprow1 < p1 && mypcol1 < q1) || (myprow1 == -1 && mypcol1 == -1));
401  ma.ctxt = context0;
402  mb.ctxt = context1;
403  /* From here, we are not assuming that only the processors working in the
404  * redistribution are calling xxMR2D, but the ones not concerned will do
405  * nothing. */
406  /* We compute the exact size of the local memory block for the memory
407  * allocations */
408  if (myprow0 >= 0 && mypcol0 >= 0) {
409  blocksize0 = memoryblocksize(&ma);
410  ma.lda = localsize(SHIFT(myprow0, ma.sprow, p0), p0, ma.nbrow, ma.m);
411  setmemory(&ptrmyblock, blocksize0);
412  initblock(ptrmyblock, 1, blocksize0);
413  setmemory(&ptrmyblockcopy, blocksize0);
414  memcpy((char *) ptrmyblockcopy, (char *) ptrmyblock,
415  blocksize0 * sizeof(double));
416  setmemory(&ptrmyblockvide, blocksize0);
417  for (i = 0; i < blocksize0; i++)
418  ptrmyblockvide[i] = -1;
419  }; /* if (mypnum < p0 * q0) */
420  if (myprow1 >= 0 && mypcol1 >= 0) {
421  setmemory(&ptrsavemyblock, memoryblocksize(&mb));
422  mb.lda = localsize(SHIFT(myprow1, mb.sprow, p1), p1, mb.nbrow, mb.m);
423  }; /* if (mypnum < p1 * q1) */
424  /* Redistribute the matrix from grid 0 to grid 1 (memory location
425  * ptrmyblock to ptrsavemyblock) */
426  Cpdgemr2d(m, n,
427  ptrmyblock, ia, ja, &ma,
428  ptrsavemyblock, ib, jb, &mb, gcontext);
429  /* Perform the inverse redistribution of the matrix from grid 1 to grid 0
430  * (memory location ptrsavemyblock to ptrmyblockvide) */
431  Cpdgemr2d(m, n,
432  ptrsavemyblock, ib, jb, &mb,
433  ptrmyblockvide, ia, ja, &ma, gcontext);
434  /* Check the differences */
435  nberrors = 0;
436  if (myprow0 >= 0 && mypcol0 >= 0) {
437  /* only for the processors that do have data at the begining */
438  for (i = 0; i < blocksize0; i++) {
439  int li, lj, gi, gj;
440  int in;
441  in = 1;
442  li = i % ma.lda;
443  lj = i / ma.lda;
444  gi = (li / ma.nbrow) * p0 * ma.nbrow +
445  SHIFT(myprow0, ma.sprow, p0) * ma.nbrow + li % ma.nbrow;
446  gj = (lj / ma.nbcol) * q0 * ma.nbcol +
447  SHIFT(mypcol0, ma.spcol, q0) * ma.nbcol + lj % ma.nbcol;
448  assert(gi < ma.m && gj < ma.n);
449  gi -= (ia - 1);
450  gj -= (ja - 1);
451  if (gi < 0 || gj < 0 || gi >= m || gj >= n)
452  in = 0;
453  if (!in) {
454  ptrmyblockcopy[i] = -1;
455  }
456  if (ptrmyblockvide[i] != ptrmyblockcopy[i]) {
457  nberrors++;
458  printf("Proc %d : Error element number %d, value = %f , initvalue =%f \n"
459  ,mypnum, i,
460  ptrmyblockvide[i], ptrmyblockcopy[i]);
461  };
462  };
463  if (nberrors > 0) {
464  printf("Processor %d, has tested %d DOUBLE PRECISION elements,\
465 Number of redistribution errors = %d \n",
466  mypnum, blocksize0, nberrors);
467  }
468  }
469  /* Look at the errors on all the processors at this point. */
470  Cigsum2d(gcontext, "All", "H", 1, 1, &nberrors, 1, 0, 0);
471  if (mypnum == 0)
472  if (nberrors)
473  printf(" => Total number of redistribution errors = %d \n",
474  nberrors);
475  else
476  printf("TEST PASSED OK\n");
477  /* release memory for the next iteration */
478  if (myprow0 >= 0 && mypcol0 >= 0) {
479  freememory((char *) ptrmyblock);
480  freememory((char *) ptrmyblockvide);
481  freememory((char *) ptrmyblockcopy);
482  }; /* if (mypnum < p0 * q0) */
483  /* release memory for the next iteration */
484  if (myprow1 >= 0 && mypcol1 >= 0) {
485  freememory((char *) ptrsavemyblock);
486  };
487  if (myprow0 >= 0)
488  Cblacs_gridexit(context0);
489  if (myprow1 >= 0)
490  Cblacs_gridexit(context1);
491  }; /* while nbre != 0 */
492  if (mypnum == 0) {
493  fclose(fp);
494  };
495  Cblacs_exit(0);
496  return 0;
497 }/* main */
Cigamx2d
void Cigamx2d()
Cblacs_pcoord
void Cblacs_pcoord()
Cblacs_pinfo
void Cblacs_pinfo()
localsize
int localsize()
Cblacs_get
void Cblacs_get()
main
int main(int argc, argv)
Definition: pdgemrdrv.c:244
MDESC::ctxt
int ctxt
Definition: pcgemr.c:165
Cpdgemr2d
void Cpdgemr2d()
BLOCK_CYCLIC_2D
#define BLOCK_CYCLIC_2D
Definition: pdgemrdrv.c:89
dispmat
#define dispmat
Definition: pdgemrdrv.c:136
MDESC::sprow
int sprow
Definition: pcgemr.c:170
MDESC::n
int n
Definition: pcgemr.c:167
MDESC::nbcol
int nbcol
Definition: pcgemr.c:169
MDESC
Definition: pcgemr.c:163
dlacpy_
#define dlacpy_
Definition: pdgemrdrv.c:74
Cdgesd2d
void Cdgesd2d()
Cblacs_gridinfo
void Cblacs_gridinfo()
paramcheck
void paramcheck()
mr2d_malloc
void * mr2d_malloc()
Cblacs_setup
void Cblacs_setup()
redist.h
MDESC::desctype
int desctype
Definition: pcgemr.c:164
initforpvm
void initforpvm(int argc, argv)
Definition: pdgemrdrv.c:225
setmemory
#define setmemory
Definition: pdgemrdrv.c:137
IDESC
Definition: pcgemr.c:175
Cigamn2d
void Cigamn2d()
Cigerv2d
void Cigerv2d()
Cblacs_exit
void Cblacs_exit()
Cigsum2d
void Cigsum2d()
initblock
static2 void initblock(double *block, int m, int n)
Definition: pdgemrdrv.c:165
Csetpvmtids
void Csetpvmtids()
scanD0
#define scanD0
Definition: pdgemrdrv.c:135
MDESC::spcol
int spcol
Definition: pcgemr.c:171
Cblacs_pnum
int Cblacs_pnum()
max
#define max(A, B)
Definition: pdgemrdrv.c:95
Cigesd2d
void Cigesd2d()
scan_intervals
#define scan_intervals
Definition: pdgemrdrv.c:139
Cpdgemr2do
void Cpdgemr2do()
nbparameter
#define nbparameter
Cigebs2d
void Cigebs2d()
changeorigin
int changeorigin()
Clacpy
#define Clacpy
Definition: pdgemrdrv.c:76
freememory
#define freememory
Definition: pdgemrdrv.c:138
SHIFT
#define SHIFT(row, sprow, nbrow)
Definition: pdgemrdrv.c:94
MDESC::lda
int lda
Definition: pcgemr.c:172
Cblacs_gridexit
void Cblacs_gridexit()
MDESC::nbrow
int nbrow
Definition: pcgemr.c:168
ppcm
int ppcm()
Cdgerv2d
void Cdgerv2d()
localindice
int localindice()
MDESC::m
int m
Definition: pcgemr.c:166
static2
#define static2
Definition: pdgemrdrv.c:61
memoryblocksize
int memoryblocksize()
Cigebr2d
void Cigebr2d()
Cblacs_gridinit
void Cblacs_gridinit()