File: large_array.c

package info (click to toggle)
openmpi 5.0.8-4
  • links: PTS, VCS
  • area: main
  • in suites:
  • size: 201,684 kB
  • sloc: ansic: 613,078; makefile: 42,353; sh: 11,194; javascript: 9,244; f90: 7,052; java: 6,404; perl: 5,179; python: 1,859; lex: 740; fortran: 61; cpp: 20; tcl: 12
file content (170 lines) | stat: -rw-r--r-- 5,620 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
/*
 * Copyright (C) by Argonne National Laboratory
 *     See COPYRIGHT in top-level directory
 */

#include "mpi.h"
#include <stdio.h>
#include <string.h>
#include <stdlib.h>

/* Writes a 4-Gbyte distributed array, reads it back, and then deletes the
   file. Uses collective I/O. */
/* The file name is taken as a command-line argument. */
/* Run it only on a machine with sufficient memory and a file system
   on which ROMIO supports large files, i.e., any file system created after
   1999 */

/* This program will work only if the MPI implementation defines MPI_Aint
   as a 64-bit integer. */

static void handle_error(int errcode, const char *str)
{
    char msg[MPI_MAX_ERROR_STRING];
    int resultlen;
    MPI_Error_string(errcode, msg, &resultlen);
    fprintf(stderr, "%s: %s\n", str, msg);
    MPI_Abort(MPI_COMM_WORLD, 1);
}

#define MPI_CHECK(fn) { int errcode; errcode = (fn); if (errcode != MPI_SUCCESS) handle_error(errcode, #fn); }




int main(int argc, char **argv)
{
    MPI_Datatype newtype;
    int i, ndims, array_of_gsizes[3], array_of_distribs[3];
    int order, nprocs, len, flag, err;
    int array_of_dargs[3], array_of_psizes[3];
    int *readbuf, *writebuf, mynod;
    MPI_Count bufcount;
    char filename[1024];
    MPI_File fh;
    MPI_Status status;
    MPI_Aint size_with_aint;
    MPI_Offset size_with_offset;

    MPI_Init(&argc, &argv);
    MPI_Comm_rank(MPI_COMM_WORLD, &mynod);
    MPI_Comm_size(MPI_COMM_WORLD, &nprocs);

/* process 0 takes the file name as a command-line argument and
   broadcasts it to other processes */
    if (!mynod) {
        i = 1;
        while ((i < argc) && strcmp("-fname", *argv)) {
            i++;
            argv++;
        }
        if (i >= argc) {
            fprintf(stderr, "\n*#  Usage: large_array -fname filename\n\n");
            MPI_Abort(MPI_COMM_WORLD, 1);
        }
        argv++;
        len = strlen(*argv);
        strcpy(filename, *argv);
        MPI_Bcast(&len, 1, MPI_INT, 0, MPI_COMM_WORLD);
        MPI_Bcast(filename, len + 1, MPI_CHAR, 0, MPI_COMM_WORLD);
        fprintf(stderr,
                "This program creates a 4 Gbyte file. Don't run it if you don't have that much disk space!\n");
    } else {
        MPI_Bcast(&len, 1, MPI_INT, 0, MPI_COMM_WORLD);
        MPI_Bcast(filename, len + 1, MPI_CHAR, 0, MPI_COMM_WORLD);
    }

/* create the distributed array filetype */
    ndims = 3;
    order = MPI_ORDER_C;

    array_of_gsizes[0] = 1024;
    array_of_gsizes[1] = 1024;
    array_of_gsizes[2] = 4 * 1024 / sizeof(int);

    array_of_distribs[0] = MPI_DISTRIBUTE_BLOCK;
    array_of_distribs[1] = MPI_DISTRIBUTE_BLOCK;
    array_of_distribs[2] = MPI_DISTRIBUTE_BLOCK;

    array_of_dargs[0] = MPI_DISTRIBUTE_DFLT_DARG;
    array_of_dargs[1] = MPI_DISTRIBUTE_DFLT_DARG;
    array_of_dargs[2] = MPI_DISTRIBUTE_DFLT_DARG;

    for (i = 0; i < ndims; i++)
        array_of_psizes[i] = 0;
    MPI_Dims_create(nprocs, ndims, array_of_psizes);

/* check if MPI_Aint is large enough for size of global array.
   if not, complain. */

    size_with_aint = sizeof(int);
    for (i = 0; i < ndims; i++)
        size_with_aint *= array_of_gsizes[i];
    size_with_offset = sizeof(int);
    for (i = 0; i < ndims; i++)
        size_with_offset *= array_of_gsizes[i];
    if (size_with_aint != size_with_offset) {
        fprintf(stderr,
                "Can't use an array of this size unless the MPI implementation defines a 64-bit MPI_Aint\n");
        MPI_Abort(MPI_COMM_WORLD, 1);
    }

    MPI_Type_create_darray(nprocs, mynod, ndims, array_of_gsizes,
                           array_of_distribs, array_of_dargs,
                           array_of_psizes, order, MPI_INT, &newtype);
    MPI_Type_commit(&newtype);

/* initialize writebuf */

    MPI_Type_size_x(newtype, &bufcount);
    bufcount = bufcount / sizeof(int);
    writebuf = (int *) malloc(bufcount * sizeof(int));
    if (!writebuf)
        fprintf(stderr, "Process %d, not enough memory for writebuf\n", mynod);
    for (i = 0; i < bufcount; i++)
        writebuf[i] = mynod * 1024 + i;

    /* write the array to the file */
    MPI_CHECK(MPI_File_open(MPI_COMM_WORLD, filename,
                            MPI_MODE_CREATE | MPI_MODE_RDWR, MPI_INFO_NULL, &fh));
    MPI_CHECK(MPI_File_set_view(fh, 0, MPI_INT, newtype, "native", MPI_INFO_NULL));
    MPI_CHECK(MPI_File_write_all(fh, writebuf, bufcount, MPI_INT, &status));
    MPI_CHECK(MPI_File_close(&fh));

    free(writebuf);

    /* now read it back */
    readbuf = (int *) calloc(bufcount, sizeof(int));
    if (!readbuf)
        fprintf(stderr, "Process %d, not enough memory for readbuf\n", mynod);

    MPI_CHECK(MPI_File_open(MPI_COMM_WORLD, filename,
                            MPI_MODE_CREATE | MPI_MODE_RDWR, MPI_INFO_NULL, &fh));
    MPI_CHECK(MPI_File_set_view(fh, 0, MPI_INT, newtype, "native", MPI_INFO_NULL));
    MPI_CHECK(MPI_File_read_all(fh, readbuf, bufcount, MPI_INT, &status));
    MPI_CHECK(MPI_File_close(&fh));

    /* check the data read */
    flag = 0;
    for (i = 0; i < bufcount; i++)
        if (readbuf[i] != mynod * 1024 + i) {
            fprintf(stderr, "Process %d, readbuf=%d, writebuf=%d\n", mynod, readbuf[i],
                    mynod * 1024 + i);
            flag = 1;
        }
    if (!flag)
        fprintf(stderr, "Process %d: data read back is correct\n", mynod);

    MPI_Type_free(&newtype);
    free(readbuf);

    MPI_Barrier(MPI_COMM_WORLD);
    if (!mynod) {
        err = MPI_File_delete(filename, MPI_INFO_NULL);
        if (err == MPI_SUCCESS)
            fprintf(stderr, "file deleted\n");
    }

    MPI_Finalize();
    return 0;
}