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 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237
|
#if defined(HAVE_CONFIG_H)
#include "config.h"
#endif
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#ifdef USE_MPI
#include <unistd.h>
#include <mpi.h>
#include <yaxt.h>
#include "pio_util.h"
#else
typedef int MPI_Comm;
#define MPI_COMM_NULL 0
#endif
#include "cdi.h"
#include "error.h"
#ifdef USE_MPI
#include "cdipio.h"
static int uniform_partition_start(int set_interval[2], int nparts, int part_idx);
#endif
static void
modelRun(MPI_Comm commModel)
{
enum
{
filetype = CDI_FILETYPE_GRB,
ntfiles = 2,
ntsteps = 3,
nVars = 5,
nlon = 12,
nlat = 6,
maxlev = 5
};
static int nlev[nVars] = { 1, 1, 5, 5, 2 };
static const char *name = "example";
int gridID, zaxisID[nVars], taxisID;
int vlistID, varID[nVars], streamID, tsID, tfID = 0;
int i, j, nmiss = 0;
double lons[nlon] = { 0, 30, 60, 90, 120, 150, 180, 210, 240, 270, 300, 330 };
double lats[nlat] = { -75, -45, -15, 15, 45, 75 };
double levs[maxlev] = { 101300, 92500, 85000, 50000, 20000 };
double var[nlon * nlat * maxlev];
int vdate = 19850101, vtime = 120000;
char filename[1024];
size_t varSize[nVars];
#if USE_MPI
int rank, comm_size;
struct var1DDeco
{
int chunkSize, start;
Xt_idxlist partDesc;
} varDeco[nVars];
#endif
#ifndef USE_MPI
(void) commModel;
#endif
gridID = gridCreate(GRID_LONLAT, nlon * nlat);
gridDefXsize(gridID, nlon);
gridDefYsize(gridID, nlat);
gridDefXvals(gridID, lons);
gridDefYvals(gridID, lats);
for (i = 0; i < nVars; i++)
{
zaxisID[i] = zaxisCreate(ZAXIS_PRESSURE, nlev[i]);
zaxisDefLevels(zaxisID[i], levs);
varSize[i] = nlon * nlat * (size_t) nlev[i];
}
vlistID = vlistCreate();
#if USE_MPI
xmpi(MPI_Comm_rank(commModel, &rank));
xmpi(MPI_Comm_size(commModel, &comm_size));
#endif
for (i = 0; i < nVars; i++)
{
varID[i] = vlistDefVar(vlistID, gridID, zaxisID[i], TIME_VARIABLE);
#ifdef USE_MPI
{
int start = uniform_partition_start((int[2]){ 0, (int) varSize[i] - 1 }, comm_size, rank),
chunkSize = uniform_partition_start((int[2]){ 0, (int) varSize[i] - 1 }, comm_size, rank + 1) - start;
fprintf(stderr, "%d: start=%d, chunkSize = %d\n", rank, start, chunkSize);
Xt_idxlist idxlist = xt_idxstripes_new(&(struct Xt_stripe){ .start = start, .nstrides = chunkSize, .stride = 1 }, 1);
varDeco[i] = (struct var1DDeco){ .start = start, .chunkSize = chunkSize, .partDesc = idxlist };
}
#endif
}
taxisID = taxisCreate(TAXIS_ABSOLUTE);
vlistDefTaxis(vlistID, taxisID);
snprintf(&filename[0], sizeof(filename), "%s_%d.grb", name, tfID);
streamID = streamOpenWrite(filename, filetype);
xassert(streamID >= 0);
streamDefVlist(streamID, vlistID);
#ifdef USE_MPI
pioEndDef();
#endif
for (tfID = 0; tfID < ntfiles; tfID++)
{
/* if ( tfID > 0 ) */
{
streamClose(streamID);
snprintf(&filename[0], sizeof(filename), "%s_%d.grb", name, tfID);
streamID = streamOpenWrite(filename, filetype);
xassert(streamID >= 0);
streamDefVlist(streamID, vlistID);
}
for (tsID = 0; tsID < ntsteps; tsID++)
{
taxisDefVdate(taxisID, vdate);
taxisDefVtime(taxisID, vtime);
streamDefTimestep(streamID, tsID);
for (i = 0; i < nVars; i++)
{
#ifdef USE_MPI
int chunk = varDeco[i].chunkSize;
#else
int chunk = (int) varSize[i];
#endif
for (j = 0; j < chunk; ++j) var[j] = 2.2;
#ifdef USE_MPI
streamWriteVarPart(streamID, varID[i], var, nmiss, varDeco[i].partDesc);
#else
streamWriteVar(streamID, varID[i], var, nmiss);
#endif
}
#ifdef USE_MPI
pioWriteTimestep();
#endif
}
}
streamClose(streamID);
vlistDestroy(vlistID);
taxisDestroy(taxisID);
for (i = 0; i < nVars; i++) zaxisDestroy(zaxisID[i]);
gridDestroy(gridID);
#ifdef USE_MPI
for (int varIdx = 0; varIdx < nVars; ++varIdx) xt_idxlist_delete(varDeco[varIdx].partDesc);
MPI_Barrier(commModel);
#endif
}
int
main(int argc, char *argv[])
{
MPI_Comm commModel = MPI_COMM_NULL;
#ifdef USE_MPI
MPI_Comm commGlob;
int sizeGlob, rankGlob, pioNamespace;
int IOMode = PIO_WRITER;
int nProcsIO = 2;
xmpi(MPI_Init(&argc, &argv));
commGlob = MPI_COMM_WORLD;
xt_initialize(commGlob);
xmpi(MPI_Comm_set_errhandler(commGlob, MPI_ERRORS_RETURN));
xmpi(MPI_Comm_size(commGlob, &sizeGlob));
xmpi(MPI_Comm_rank(commGlob, &rankGlob));
{
int opt;
while ((opt = getopt(argc, argv, "p:w:")) != -1) switch (opt)
{
case 'p':
IOMode = cdiPioStr2IOMode(optarg);
if (IOMode < 0)
{
fprintf(stderr, "Unsupported PIO mode requested: %s\n", optarg);
exit(EXIT_FAILURE);
}
break;
case 'w':
{
long temp = strtol(optarg, NULL, 0);
if (temp < 0 || temp > INT_MAX / 2)
{
fprintf(stderr, "Unsupported number of I/O servers: %ld\n", temp);
exit(EXIT_FAILURE);
}
nProcsIO = (int) temp;
break;
}
}
}
commModel = pioInit(commGlob, nProcsIO, IOMode, &pioNamespace, 1.0f, cdiPioNoPostCommSetup);
if (commModel != MPI_COMM_NULL)
{
namespaceSetActive(pioNamespace);
#else
(void) argc;
(void) argv;
#endif
modelRun(commModel);
#ifdef USE_MPI
}
pioFinalize();
xt_finalize();
MPI_Finalize();
#endif
return 0;
}
#ifdef USE_MPI
static int
uniform_partition_start(int set_interval[2], int nparts, int part_idx)
{
int part_offset
= (int) ((((long long) set_interval[1] - (long long) set_interval[0] + 1LL) * (long long) part_idx) / (long long) nparts);
int start = set_interval[0] + part_offset;
return start;
}
#endif
/*
* Local Variables:
* c-file-style: "Java"
* c-basic-offset: 2
* indent-tabs-mode: nil
* show-trailing-whitespace: t
* require-trailing-newline: t
* End:
*/
|