File: comm_init.c

package info (click to toggle)
ga 5.9.2-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 18,472 kB
  • sloc: ansic: 192,963; fortran: 53,761; f90: 11,218; cpp: 5,784; makefile: 2,248; sh: 1,945; python: 1,734; perl: 534; csh: 134; asm: 106
file content (124 lines) | stat: -rw-r--r-- 2,765 bytes parent folder | download | duplicates (3)
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

#if HAVE_CONFIG_H
#   include "config.h"
#endif

#include <stdio.h>
#include <stdlib.h>
#include <math.h>

#include "mpi.h"
#include "ga.h"

#include "mp3.h"
#include "macdecls.h"

#define DIM 256

void comm_test(MPI_Comm comm)
{
  double *buf;
  int i, j, ii, jj, ld;
  int lo[2], hi[2];
  int dims[2];
  int two = 2;
  int g_a;
  int me, nprocs, prem;
  int nelems, ok;
  int rank;
  MPI_Comm_rank(comm,&rank);
  if (GA_Initialize_comm(comm)) {
    me = GA_Nodeid();
    nprocs = GA_Nnodes();
    prem = (me+1)%nprocs;
    if (me == 0) {
      printf("\nRunning comm_test on %d processors\n\n",nprocs);
    }

    g_a = NGA_Create_handle();
    dims[0] = DIM;
    dims[1] = DIM;
    NGA_Set_data(g_a,two,dims,C_DBL);
    NGA_Allocate(g_a);
    GA_Zero(g_a);

    /* Find out what section of global array I'm writing to */
    NGA_Distribution(g_a,prem,lo,hi);
    nelems = (hi[0]-lo[0]+1)*(hi[1]-lo[1]+1);
    ld = hi[1]-lo[1]+1;
    buf = (double*)malloc(nelems*sizeof(double));
    for (i=lo[0]; i<=hi[0]; i++) {
      ii = i-lo[0];
      for (j=lo[1]; j<=hi[1]; j++) {
        jj = j-lo[1];
        buf[ii*ld+jj] = (double)(i*DIM+j);
      }
    }
    NGA_Put(g_a,lo,hi,buf,&ld);
    GA_Sync();
    free(buf);
    /* Copy data from a different section of global array and check values */
    prem = me-1;
    if (prem < 0) prem = nprocs-1;
    NGA_Distribution(g_a,prem,lo,hi);
    nelems = (hi[0]-lo[0]+1)*(hi[1]-lo[1]+1);
    ld = hi[1]-lo[1]+1;
    buf = (double*)malloc(nelems*sizeof(double));
    for (i=0; i<nelems; i++) buf[i] = 0.0;
    NGA_Get(g_a,lo,hi,buf,&ld);

    ok = 1;
    for (i=lo[0]; i<=hi[0]; i++) {
      ii = i-lo[0];
      for (j=lo[1]; j<=hi[1]; j++) {
        jj = j-lo[1];
        if (buf[ii*ld+jj] != (double)(i*DIM+j)) {
          if (ok == 1) {
            printf("p[%d] Mismatch for (%d,%d) expected: %f actual: %f\n",
                me,i,j,buf[ii*ld+jj],(double)(i*DIM+j));
            ok = 0;
          }
        }
      }
    }
    free(buf);
    if (me == 0 && ok) {
      printf("Test passed\n");
    } else if (!ok) {
      printf("Test failed\n");
    }
    GA_Terminate();
  }
}

int main(int argc, char **argv)
{
  int rank, size;
  int color;
  MPI_Comm group;

  MP_INIT(argc, argv);

  MPI_Comm_size(MPI_COMM_WORLD, &size);
  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
  if (rank == 0) {
    printf("MPI COMM WORLD has %d processors\n",size);
  }
  color = 1;
  if (rank < size/2) color = 0;
  MPI_Comm_split(MPI_COMM_WORLD,color,rank,&group);
  if (rank == 0) {
    printf("\nRun test on first group\n");
  }
  if (color == 0) {
    comm_test(group);
  }
  MPI_Barrier(MPI_COMM_WORLD);
  if (rank == 0) {
    printf("\nRun test on second group\n");
  }
  if (color == 1) {
    comm_test(group);
  }
  MPI_Finalize();
}