~ubuntu-branches/ubuntu/trusty/silo-llnl/trusty

« back to all changes in this revision

Viewing changes to tests/pmpio_silo_test_mesh.c

  • Committer: Bazaar Package Importer
  • Author(s): Alastair McKinstry
  • Date: 2011-01-02 00:03:01 UTC
  • Revision ID: james.westby@ubuntu.com-20110102000301-9s2hfsjrkguz6h4r
Tags: upstream-4.8
ImportĀ upstreamĀ versionĀ 4.8

Show diffs side-by-side

added added

removed removed

Lines of Context:
 
1
/*
 
2
Copyright (c) 1994 - 2010, Lawrence Livermore National Security, LLC.
 
3
LLNL-CODE-425250.
 
4
All rights reserved.
 
5
 
 
6
This file is part of Silo. For details, see silo.llnl.gov.
 
7
 
 
8
Redistribution and use in source and binary forms, with or without
 
9
modification, are permitted provided that the following conditions
 
10
are met:
 
11
 
 
12
   * Redistributions of source code must retain the above copyright
 
13
     notice, this list of conditions and the disclaimer below.
 
14
   * Redistributions in binary form must reproduce the above copyright
 
15
     notice, this list of conditions and the disclaimer (as noted
 
16
     below) in the documentation and/or other materials provided with
 
17
     the distribution.
 
18
   * Neither the name of the LLNS/LLNL nor the names of its
 
19
     contributors may be used to endorse or promote products derived
 
20
     from this software without specific prior written permission.
 
21
 
 
22
THIS SOFTWARE  IS PROVIDED BY  THE COPYRIGHT HOLDERS  AND CONTRIBUTORS
 
23
"AS  IS" AND  ANY EXPRESS  OR IMPLIED  WARRANTIES, INCLUDING,  BUT NOT
 
24
LIMITED TO, THE IMPLIED  WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
 
25
A  PARTICULAR  PURPOSE ARE  DISCLAIMED.  IN  NO  EVENT SHALL  LAWRENCE
 
26
LIVERMORE  NATIONAL SECURITY, LLC,  THE U.S.  DEPARTMENT OF  ENERGY OR
 
27
CONTRIBUTORS BE LIABLE FOR  ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
 
28
EXEMPLARY, OR  CONSEQUENTIAL DAMAGES  (INCLUDING, BUT NOT  LIMITED TO,
 
29
PROCUREMENT OF  SUBSTITUTE GOODS  OR SERVICES; LOSS  OF USE,  DATA, OR
 
30
PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
 
31
LIABILITY, WHETHER  IN CONTRACT, STRICT LIABILITY,  OR TORT (INCLUDING
 
32
NEGLIGENCE OR  OTHERWISE) ARISING IN  ANY WAY OUT  OF THE USE  OF THIS
 
33
SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
 
34
 
 
35
This work was produced at Lawrence Livermore National Laboratory under
 
36
Contract  No.   DE-AC52-07NA27344 with  the  DOE.  Neither the  United
 
37
States Government  nor Lawrence  Livermore National Security,  LLC nor
 
38
any of  their employees,  makes any warranty,  express or  implied, or
 
39
assumes   any   liability   or   responsibility  for   the   accuracy,
 
40
completeness, or usefulness of any information, apparatus, product, or
 
41
process  disclosed, or  represents  that its  use  would not  infringe
 
42
privately-owned   rights.  Any  reference   herein  to   any  specific
 
43
commercial products,  process, or  services by trade  name, trademark,
 
44
manufacturer or otherwise does not necessarily constitute or imply its
 
45
endorsement,  recommendation,   or  favoring  by   the  United  States
 
46
Government or Lawrence Livermore National Security, LLC. The views and
 
47
opinions  of authors  expressed  herein do  not  necessarily state  or
 
48
reflect those  of the United  States Government or  Lawrence Livermore
 
49
National  Security, LLC,  and shall  not  be used  for advertising  or
 
50
product endorsement purposes.
 
51
*/
 
52
#include <mpi.h>
 
53
#include <silo.h>
 
54
#include <pmpio.h>
 
55
#include <string.h>
 
56
#include <stdlib.h>
 
57
 
 
58
void WriteMultiXXXObjects(DBfile *siloFile, PMPIO_baton_t *bat, int size,
 
59
    const char *file_ext);
 
60
 
 
61
 
 
62
/*-----------------------------------------------------------------------------
 
63
 * Purpose:     Impliment the create callback to initialize pmpio
 
64
 *              Will create the silo file and the 'first' directory (namespace)
 
65
 *              in it. The driver type (DB_PDB or DB_HDF5) is passed as user
 
66
 *              data; a void pointer to the driver determined in main.
 
67
 *-----------------------------------------------------------------------------
 
68
 */
 
69
void *CreateSiloFile(const char *fname, const char *nsname, void *userData)
 
70
{
 
71
    int driver = *((int*) userData);
 
72
    DBfile *siloFile = DBCreate(fname, DB_CLOBBER, DB_LOCAL, "pmpio testing", driver);
 
73
    if (siloFile)
 
74
    {
 
75
        DBMkDir(siloFile, nsname);
 
76
        DBSetDir(siloFile, nsname);
 
77
    }
 
78
    return (void *) siloFile;
 
79
}
 
80
 
 
81
/*-----------------------------------------------------------------------------
 
82
 * Purpose:     Impliment the open callback to initialize pmpio
 
83
 *              Will open the silo file and, for write, create the new
 
84
 *              directory or, for read, just cd into the right directory.
 
85
 *-----------------------------------------------------------------------------
 
86
 */
 
87
void *OpenSiloFile(const char *fname, const char *nsname, PMPIO_iomode_t ioMode,
 
88
    void *userData)
 
89
{
 
90
    DBfile *siloFile = DBOpen(fname, DB_UNKNOWN,
 
91
        ioMode == PMPIO_WRITE ? DB_APPEND : DB_READ);
 
92
    if (siloFile)
 
93
    {
 
94
        if (ioMode == PMPIO_WRITE)
 
95
            DBMkDir(siloFile, nsname);
 
96
        DBSetDir(siloFile, nsname);
 
97
    }
 
98
    return (void *) siloFile;
 
99
}
 
100
 
 
101
/*-----------------------------------------------------------------------------
 
102
 * Purpose:     Impliment the close callback for pmpio
 
103
 *-----------------------------------------------------------------------------
 
104
 */
 
105
void CloseSiloFile(void *file, void *userData)
 
106
{
 
107
    DBfile *siloFile = (DBfile *) file;
 
108
    if (siloFile)
 
109
        DBClose(siloFile);
 
110
}
 
111
 
 
112
/*-----------------------------------------------------------------------------
 
113
 * Audience:    Public
 
114
 * Chapter:     Main 
 
115
 * Purpose:     Demonstrate use of PMPIO 
 
116
 * Description:
 
117
 * This simple program demonstrates the use of PMPIO to write a set of silo
 
118
 * files for a very simple quad mesh. The mesh will be N x N where N is the
 
119
 * number of processors. Each processor will write a single 1 x N strip of
 
120
 * zones (2 nodes and 1 zone wide). The second command line argument 
 
121
 * indicates the total number of files that will be produced. Note that this
 
122
 * is totally independent of the number of processors. The resulting silo
 
123
 * file can be visualized, in parallel, in VisIt by opening the 'root' file
 
124
 * named "silo_000". Note that PMPIO's role is merely to coordinate and
 
125
 * manage access to the silo file(s) that get created.
 
126
 *
 
127
 * By default, this example will use Silo's PDB driver. However, if you pass
 
128
 * "DB_HDF5" anywhere on the command line, it will use the HDF5 driver. Any
 
129
 * integer appearing on the command line is taken to be the total number of
 
130
 * files.
 
131
 *
 
132
 * An example of how you would invoke this example is...
 
133
 *
 
134
 *     mpirun -np 17 pmpio_silo_test_mesh 3 DB_HDF5     
 
135
 *
 
136
 * which would run on 17 processors, creating a 17x17 mesh but writing it to
 
137
 * 3 files using the HDF5 driver.
 
138
 *-----------------------------------------------------------------------------
 
139
 */
 
140
int main(int argc, char **argv)
 
141
{
 
142
    int size, rank;
 
143
    int numGroups = -1;
 
144
    int driver = DB_PDB;
 
145
    DBfile *siloFile;
 
146
    char *file_ext = "pdb";
 
147
    char fileName[256], nsName[256];
 
148
    int i, j, dims[2], ndims = 2;
 
149
    char *coordnames[2], *varnames[2];
 
150
    float *x, *y, *coords[2], *vars[2];
 
151
    float *vx, *vy;
 
152
    float *temp;
 
153
    PMPIO_baton_t *bat; 
 
154
 
 
155
    /* specify the number of Silo files to create and driver */
 
156
    for (i = 1; i < argc; i++)
 
157
    {
 
158
        if (!strcmp(argv[i], "DB_HDF5"))
 
159
        {
 
160
            driver = DB_HDF5;
 
161
            file_ext = "h5";
 
162
        }
 
163
        else if (strtol(argv[i], 0, 10) > 0)
 
164
        {
 
165
            numGroups = strtol(argv[i], 0, 10);
 
166
        }
 
167
        else if (argv[i][0] != '\0')
 
168
        {
 
169
            fprintf(stderr, "%s: ignored argument `%s'\n", argv[0], argv[i]);
 
170
        }
 
171
    }
 
172
 
 
173
    /* standard MPI initialization stuff */
 
174
    MPI_Init(&argc, &argv);
 
175
    MPI_Comm_size(MPI_COMM_WORLD, &size);
 
176
    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
 
177
    if (rank == 0)
 
178
    {
 
179
        if (numGroups < 0)
 
180
            numGroups = 3; 
 
181
    }
 
182
 
 
183
    /* Ensure all procs agree on numGroups, driver and file_ext */
 
184
    MPI_Bcast(&numGroups, 1, MPI_INT, 0, MPI_COMM_WORLD);
 
185
    MPI_Bcast(&driver, 1, MPI_INT, 0, MPI_COMM_WORLD);
 
186
    if (driver == DB_HDF5)
 
187
        file_ext = "h5";
 
188
 
 
189
    /* Initialize PMPIO, pass a pointer to the driver type as the
 
190
       user data. */
 
191
    bat = PMPIO_Init(numGroups, PMPIO_WRITE, MPI_COMM_WORLD, 1,
 
192
        CreateSiloFile, OpenSiloFile, CloseSiloFile, &driver);
 
193
 
 
194
    /* Construct names for the silo files and the directories in them */
 
195
    sprintf(fileName, "silo_%03d.%s", PMPIO_GroupRank(bat, rank), file_ext);
 
196
    sprintf(nsName, "domain_%03d", rank);
 
197
 
 
198
    /* Wait for write access to the file. All processors call this.
 
199
     * Some processors (the first in each group) return immediately
 
200
     * with write access to the file. Other processors wind up waiting
 
201
     * until they are given control by the preceeding processor in 
 
202
     * the group when that processor calls "HandOffBaton" */
 
203
    siloFile = (DBfile *) PMPIO_WaitForBaton(bat, fileName, nsName);
 
204
 
 
205
    /* Do some work to construct the mesh data for this processor. */
 
206
    dims[0] = 2;
 
207
    dims[1] = size;
 
208
    coordnames[0] = "xcoords";
 
209
    coordnames[1] = "ycoords";
 
210
    x = (float *) malloc(2 * sizeof(float));
 
211
    x[0] = (float) rank;
 
212
    x[1] = (float) rank+1;
 
213
    y = (float *) malloc(size * sizeof(float));
 
214
    for (i = 0; i < size; i++)
 
215
        y[i] = (float) i;
 
216
    coords[0] = x;
 
217
    coords[1] = y;
 
218
 
 
219
    /* Do some work to create some simple variable data for this processor. */
 
220
    vx = (float *) malloc(2 * size * sizeof(float));
 
221
    vy = (float *) malloc(2 * size * sizeof(float));
 
222
    for (j = 0; j < size; j++)
 
223
    {
 
224
        for (i = 0; i < 2; i++)
 
225
        {
 
226
            vx[j*2+i] = ((float) (size-(rank+i)))/((float) (2*size));
 
227
            vy[j*2+i] = ((float) ((size/2) - (rank+i))) / ((float) size);
 
228
        }
 
229
    }
 
230
    temp = (float *) malloc(size * sizeof(float));
 
231
    for (i = 0; i < size; i++)
 
232
        temp[i] = (float) i * rank;
 
233
 
 
234
    /*
 
235
     *
 
236
     * BEGIN THIS PROCESSOR'S LOCAL WORK ON A SILO FILE
 
237
     *
 
238
     */
 
239
 
 
240
    /* This processor's local work on the file */
 
241
    DBPutQuadmesh(siloFile, "qmesh", coordnames, coords, dims, ndims,
 
242
        DB_FLOAT, DB_COLLINEAR, 0);
 
243
    free(x);
 
244
    free(y);
 
245
 
 
246
    /* Output velocity on the mesh */
 
247
    vars[0] = vx;
 
248
    vars[1] = vy;
 
249
    varnames[0] = "vx";
 
250
    varnames[1] = "vy";
 
251
    DBPutQuadvar(siloFile, "velocity", "qmesh", 2, varnames, vars,
 
252
        dims, ndims, NULL, 0, DB_FLOAT, DB_NODECENT, 0);
 
253
    free(vx);
 
254
    free(vy);
 
255
 
 
256
    /* Output temp on the mesh */
 
257
    dims[0] = 1;
 
258
    dims[1] = size-1;
 
259
    DBPutQuadvar1(siloFile, "temp", "qmesh", temp, dims, ndims,
 
260
        NULL, 0, DB_FLOAT, DB_ZONECENT, 0);
 
261
    free(temp);
 
262
 
 
263
    /*
 
264
     *
 
265
     * END THIS PROCESSORS LOCAL WORK ON THE FILE
 
266
     *
 
267
     */
 
268
 
 
269
    /* If this is the 'root' processor, also write Silo's multi-XXX objects */
 
270
    if (rank == 0)
 
271
        WriteMultiXXXObjects(siloFile, bat, size, file_ext);
 
272
 
 
273
    /* Hand off the baton to the next processor. This winds up closing
 
274
     * the file so that the next processor that opens it can be assured
 
275
     * of getting a consistent and up to date view of the file's contents. */
 
276
    PMPIO_HandOffBaton(bat, siloFile);
 
277
 
 
278
    /* We're done using PMPIO, so finish it off */
 
279
    PMPIO_Finish(bat);
 
280
 
 
281
    /* Standard MPI finalization */
 
282
    MPI_Finalize();
 
283
 
 
284
    return 0;
 
285
}
 
286
 
 
287
 
 
288
void WriteMultiXXXObjects(DBfile *siloFile, PMPIO_baton_t *bat, int size,
 
289
    const char *file_ext)
 
290
{
 
291
    int i;
 
292
    char **meshBlockNames = (char **) malloc(size * sizeof(char*));
 
293
    char **tempBlockNames = (char **) malloc(size * sizeof(char*));
 
294
    char **velBlockNames = (char **) malloc(size * sizeof(char*));
 
295
    int *blockTypes = (int *) malloc(size * sizeof(int));
 
296
    int *varTypes = (int *) malloc(size * sizeof(int));
 
297
 
 
298
    /* Go to root directory in the silo file */
 
299
    DBSetDir(siloFile, "/");
 
300
 
 
301
    /* Construct the lists of individual object names */
 
302
    for (i = 0; i < size; i++)
 
303
    {
 
304
        int groupRank = PMPIO_GroupRank(bat, i);
 
305
        meshBlockNames[i] = (char *) malloc(1024);
 
306
        velBlockNames[i] = (char *) malloc(1024);
 
307
        tempBlockNames[i] = (char *) malloc(1024);
 
308
        if (groupRank == 0)
 
309
        {
 
310
            /* this mesh block is in the file 'root' owns */
 
311
            sprintf(meshBlockNames[i], "/domain_%03d/qmesh", i);
 
312
            sprintf(velBlockNames[i], "/domain_%03d/velocity", i);
 
313
            sprintf(tempBlockNames[i], "/domain_%03d/temp", i);
 
314
        }
 
315
        else
 
316
        {
 
317
            /* this mesh block is another file */ 
 
318
            sprintf(meshBlockNames[i], "silo_%03d.%s:/domain_%03d/qmesh",
 
319
                groupRank, file_ext, i);
 
320
            sprintf(velBlockNames[i], "silo_%03d.%s:/domain_%03d/velocity",
 
321
                groupRank, file_ext, i);
 
322
            sprintf(tempBlockNames[i], "silo_%03d.%s:/domain_%03d/temp",
 
323
                groupRank, file_ext, i);
 
324
        }
 
325
        blockTypes[i] = DB_QUADMESH;
 
326
        varTypes[i] = DB_QUADVAR;
 
327
    }
 
328
 
 
329
    /* Write the multi-block objects */
 
330
    DBPutMultimesh(siloFile, "multi_qmesh", size, meshBlockNames, blockTypes, 0);
 
331
    DBPutMultivar(siloFile, "multi_velocity", size, velBlockNames, varTypes, 0);
 
332
    DBPutMultivar(siloFile, "multi_temp", size, tempBlockNames, varTypes, 0);
 
333
 
 
334
    /* Clean up */
 
335
    for (i = 0; i < size; i++)
 
336
    {
 
337
        free(meshBlockNames[i]);
 
338
        free(velBlockNames[i]);
 
339
        free(tempBlockNames[i]);
 
340
    }
 
341
    free(meshBlockNames);
 
342
    free(velBlockNames);
 
343
    free(tempBlockNames);
 
344
    free(blockTypes);
 
345
    free(varTypes);
 
346
}