~ubuntu-branches/ubuntu/vivid/adios/vivid-proposed

« back to all changes in this revision

Viewing changes to examples/C/schema/structured2d_noxml.c

  • Committer: Package Import Robot
  • Author(s): Alastair McKinstry
  • Date: 2014-06-16 23:06:38 UTC
  • mfrom: (15.1.8 sid)
  • Revision ID: package-import@ubuntu.com-20140616230638-cxryhot6b8ge32l6
Tags: 1.7.0-1
* New upstream release.
* Add adios.pc pkgconfig file. adios_config now uses this.

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
/* 
 
2
 * ADIOS is freely available under the terms of the BSD license described
 
3
 * in the COPYING file in the top level directory of this source distribution.
 
4
 *
 
5
 * Copyright (c) 2008 - 2009.  UT-BATTELLE, LLC. All rights reserved.
 
6
 */
 
7
 
 
8
/* ADIOS C Example: write variables along with an rectilinear mesh. 
 
9
 */
 
10
#include <stdio.h>
 
11
#include <stdlib.h>
 
12
#include <string.h>
 
13
#include <unistd.h>
 
14
#include <fcntl.h>
 
15
#include <errno.h>
 
16
#include "mpi.h"
 
17
#include "public/adios.h"
 
18
#include "public/adios_types.h"
 
19
 
 
20
//will work with 12 cores, which are arranged by npx=4, npy=3 (4x3)
 
21
char   npx_str[256];       // # of procs in x dim (string value)
 
22
char   npy_str[256];       // # of procs in y dim (string value)
 
23
int    npx;                // # of procs in x direction
 
24
int    npy;                // # of procs in y direction
 
25
int    nproc;               // # of total procs
 
26
 
 
27
void printUsage(char *prgname)
 
28
{
 
29
    printf("Usage: mpirun -np <N> %s <nx> <ny>\n"
 
30
           "    <nx> <ny>  2D decomposition values in each dimension of an 2D array\n"
 
31
           "         The product of these number must be equal the number of processes\n"
 
32
           "         e.g. for N=12 you may use  4 3\n"
 
33
        ,prgname);
 
34
}
 
35
 
 
36
 
 
37
int processArgs(int argc, char ** argv)
 
38
{
 
39
    if (argc < 3) {
 
40
        printUsage (argv[0]);
 
41
        return 1;
 
42
    }
 
43
 
 
44
    strncpy(npx_str, argv[1], sizeof(npx_str));
 
45
    strncpy(npy_str, argv[2], sizeof(npy_str));
 
46
 
 
47
    npx = atoi(npx_str);
 
48
    npy = atoi(npy_str);
 
49
 
 
50
    if (npx*npy != nproc) {
 
51
        printf ("ERROR: Product of decomposition numbers in X and Y dimension %d != number of processes %d\n", npx*npy, nproc);
 
52
        printUsage(argv[0]);
 
53
        return 1;
 
54
    }
 
55
 
 
56
    return 0;
 
57
}
 
58
 
 
59
 
 
60
int main (int argc, char ** argv) 
 
61
{
 
62
    MPI_Comm    comm = MPI_COMM_WORLD;
 
63
    int         rank;
 
64
    int         ndx, ndy;             // size of array per processor
 
65
    double      * data;
 
66
 
 
67
    double      *X;                   //X coordinate
 
68
    double      *Y;                   //Y coordinate
 
69
 
 
70
    // Offsets and sizes
 
71
    int         offs_x, offs_y;       //offset in x and y direction
 
72
    int         nx_local, ny_local;   //local address
 
73
    int         nx_global, ny_global; //global address
 
74
    int         posx, posy;           // position index in the array
 
75
    int         i,j;
 
76
 
 
77
    int64_t     m_adios_group;
 
78
    uint64_t    adios_groupsize, adios_totalsize;
 
79
    int64_t     adios_handle;
 
80
 
 
81
    MPI_Init (&argc, &argv);
 
82
    MPI_Comm_rank (comm, &rank);
 
83
    MPI_Comm_size (comm, &nproc);
 
84
 
 
85
    if (processArgs(argc, argv)) {
 
86
        return 1;
 
87
    }
 
88
 
 
89
    //will work with each core writing ndx = 65, ndy = 129, (65*4,129*3) global
 
90
    ndx = 65;
 
91
    ndy = 129;
 
92
 
 
93
    //2D array with block,block decomposition
 
94
    posx = rank%npx;           // 1st dim
 
95
    posy = rank/npx;           // 2nd dim
 
96
    offs_x = posx * ndx;
 
97
    offs_y = posy * ndy;
 
98
    nx_local = ndx;
 
99
    ny_local = ndy;
 
100
    nx_global = npx * ndx;
 
101
    ny_global = npy * ndy;
 
102
 
 
103
    data = malloc (ndx * ndy * sizeof(double));
 
104
    for( i = 0; i < ndx; i++ )
 
105
        for( j = 0; j < ndy; j++)
 
106
            data[i*ndy + j] = 1.0*rank;
 
107
 
 
108
    X = malloc (ndx * ndy * sizeof(double));
 
109
    for( i = 0; i < ndx; i++ )
 
110
        for( j = 0; j < ndy; j++)
 
111
            X[i*ndy + j] = offs_x + posy*ndx + i*ndx/ndx + (double)ndx*j/ndy;
 
112
 
 
113
    Y = malloc (ndx * ndy * sizeof(double));
 
114
    Y[0] = offs_y;
 
115
    for( i = 0; i < ndx; i++ )
 
116
        for( j = 0; j < ndy; j++)
 
117
            Y[i*ndy + j] = offs_y + ndy*j/ndy;
 
118
  
 
119
    char * schema_version = "1.1";
 
120
    char * dimemsions = "nx_global,ny_global";
 
121
 
 
122
        adios_init_noxml (comm);
 
123
    adios_allocate_buffer (ADIOS_BUFFER_ALLOC_NOW, 50);
 
124
 
 
125
    adios_declare_group (&m_adios_group, "structured2d", "", adios_flag_yes);
 
126
    adios_select_method (m_adios_group, "MPI", "", "");
 
127
 
 
128
    adios_define_var (m_adios_group, "nx_global"
 
129
                        ,"", adios_integer
 
130
                        ,0, 0, 0);
 
131
    adios_define_var (m_adios_group, "ny_global"
 
132
            ,"", adios_integer
 
133
            ,0, 0, 0);
 
134
    adios_define_var (m_adios_group, "nproc"
 
135
                ,"", adios_integer                
 
136
                ,0, 0, 0);
 
137
    adios_define_var (m_adios_group, "offs_x"
 
138
                ,"", adios_integer
 
139
                ,0, 0, 0);
 
140
    adios_define_var (m_adios_group, "offs_y"
 
141
                ,"", adios_integer
 
142
                ,0, 0, 0);
 
143
    adios_define_var (m_adios_group, "nx_local"
 
144
                ,"", adios_integer
 
145
                ,0, 0, 0);
 
146
    adios_define_var (m_adios_group, "ny_local"
 
147
                ,"", adios_integer
 
148
                ,0, 0, 0);
 
149
    adios_define_var (m_adios_group, "X"
 
150
                    ,"", adios_double
 
151
                    ,"nx_local,ny_local", "nx_global,ny_global", "offs_x,offs_y");
 
152
    adios_define_var (m_adios_group, "Y"
 
153
                    ,"", adios_double
 
154
                    ,"nx_local,ny_local", "nx_global,ny_global", "offs_x,offs_y");
 
155
    adios_define_var (m_adios_group, "data"
 
156
                    ,"", adios_double
 
157
                    ,"nx_local,ny_local", "nx_global,ny_global", "offs_x,offs_y");
 
158
 
 
159
    adios_define_schema_version (m_adios_group, schema_version);
 
160
    adios_define_mesh_structured (dimemsions, "X,Y", "2", m_adios_group, "structuredmesh");
 
161
    adios_define_mesh_timevarying ("no", m_adios_group, "structuredmesh");
 
162
    adios_define_var_mesh (m_adios_group, "data", "structuredmesh");
 
163
    adios_define_var_centering (m_adios_group, "data", "point");
 
164
 
 
165
    adios_open (&adios_handle, "structured2d", "structured2d_noxml.bp", "w", comm);
 
166
 
 
167
    adios_groupsize = 7*sizeof(int) \
 
168
    + 3*sizeof(double) * (nx_local*ny_local);
 
169
 
 
170
    adios_group_size (adios_handle, adios_groupsize, &adios_totalsize);
 
171
    adios_write (adios_handle, "nproc", &nproc);
 
172
    adios_write (adios_handle, "nx_global", &nx_global);
 
173
    adios_write (adios_handle, "ny_global", &ny_global);
 
174
    adios_write (adios_handle, "offs_x", &offs_x);
 
175
    adios_write (adios_handle, "offs_y", &offs_y);
 
176
    adios_write (adios_handle, "nx_local", &nx_local);
 
177
    adios_write (adios_handle, "ny_local", &ny_local);
 
178
    adios_write (adios_handle, "X", X);
 
179
    adios_write (adios_handle, "Y", Y);
 
180
    adios_write (adios_handle, "data", data);
 
181
 
 
182
    adios_close (adios_handle);
 
183
 
 
184
    MPI_Barrier (comm);
 
185
 
 
186
    free (data);
 
187
    free (X);
 
188
    free (Y);
 
189
 
 
190
        adios_finalize (rank);
 
191
 
 
192
        MPI_Finalize ();
 
193
        return 0;
 
194
}