test_resource_copy.c 6.78 KB
Newer Older
1
2
3
4
5
6
7
8
9
#if defined (HAVE_CONFIG_H)
#include "config.h"
#endif

#include <stdio.h>

#ifdef USE_MPI
#include <mpi.h>
#include "cdi.h"
10
#include "cdipio.h"
11
#include "pio_util.h"
12
13
14
#include "resource_handle.h"
#include "resource_unpack.h"
#include "pio_serialize.h"
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

extern void   arrayDestroy          ( void );

enum {
  IOMode           = PIO_NONE,
  nProcsIO         = 1,
  nNamespaces      = 2,
  DOUBLE_PRECISION = 8,
  nlon             = 12,
  nlat             = 6,
  nlev             = 5,
  ntsteps          = 3 };

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[nlev] = {101300, 92500, 85000, 50000, 20000};

int defineGrid ()
{
  int gridID = CDI_UNDEFID;
  int mask_vec[nlon*nlat];
  const int * mp = &mask_vec[0];
  double area_vec[nlon*nlat];
  const double * ap = &area_vec[0];
  int i;

  gridID = gridCreate(GRID_LONLAT, nlon*nlat);
  gridDefXsize(gridID, nlon);
  gridDefYsize(gridID, nlat);
  gridDefXvals(gridID, lons);
  gridDefYvals(gridID, lats);
  gridDefNvertex ( gridID, 1 );
  gridDefXbounds ( gridID, lons );
  gridDefYbounds ( gridID, lats );
  for ( i = 0; i < nlon*nlat; i++ )
    mask_vec[i] = i % 2 ;
  gridDefMaskGME ( gridID, mp );
  for ( i = 0; i < nlon*nlat; i++ )
    mask_vec[i] = 1;
  gridDefMask ( gridID, mp );
  gridDefXname ( gridID, "myXname" );
  gridDefXlongname ( gridID, "myXlongname" );
  gridDefXunits ( gridID, "myXunits" );
  gridDefYname ( gridID, "myYname" );
  gridDefYlongname ( gridID, "myYlongname" );
  gridDefYunits ( gridID, "myYunits" );
  gridDefPrec ( gridID, DOUBLE_PRECISION );
  gridDefXpole ( gridID, 90.0 );
  gridDefYpole ( gridID, 180.0 );
  gridDefAngle ( gridID, 360.0 );
  gridDefTrunc ( gridID, 1 );
  gridDefGMEnd ( gridID, 2 );
  gridDefGMEni ( gridID, 3 );
  gridDefGMEni2 ( gridID, 4 );
  gridDefGMEni3 ( gridID, 5 );
  gridDefNumber ( gridID, 6 );
  gridDefPosition ( gridID, 7 );
  gridDefReference ( gridID, "myReference" );
/* gridDefLCC ( gridID, double originLon, double originLat,  */
/*         double lonParY, double lat1, double lat2, double xinc, double yinc, int projflag, int scanflag); */
/* gridDefLcc2 ( gridID, double earth_radius, double lon_0,  */
/*          double lat_0, double lat_1,double lat_2);*/
/* gridDefLaea ( gridID, double earth_radius, double lon_0, double lat_0); */
  for ( i = 0; i < nlon*nlat; i++ )
    area_vec[i] = 0.1 * i;
  gridDefArea ( gridID, ap );
  for ( i = 0; i < nlon*nlat; i++ )
    mask_vec[i] = i;
  gridDefRowlon ( gridID, nlon*nlat, mp );
  gridDefComplexPacking ( gridID, 1 );

  return gridID;
}

int defineZaxis ()
{
  int zaxisID = CDI_UNDEFID;
  double vct[3] = { 3.0, 3.3, 3.6 };

  zaxisID = zaxisCreate(ZAXIS_PRESSURE, nlev);
  zaxisDefLevels(zaxisID, levs);
  zaxisDefLevel ( zaxisID, 2, 8507.3 );
  zaxisDefName ( zaxisID, "myName" );
  zaxisDefLongname ( zaxisID, "myLongname" );
  zaxisDefUnits ( zaxisID, "myUnits" );
  zaxisDefPrec ( zaxisID, DOUBLE_PRECISION );
  zaxisDefLtype ( zaxisID, 1 );
  zaxisDefVct ( zaxisID, 3, vct );
  zaxisDefLbounds ( zaxisID, &levs[0] );
  zaxisDefUbounds ( zaxisID, &levs[0] );
  zaxisDefWeights ( zaxisID, &levs[0] );

  return zaxisID;
}

int defineTaxis ()
{
  int taxisID = CDI_UNDEFID;

  taxisID = taxisCreate(TAXIS_ABSOLUTE);

  taxisDefType  ( taxisID, 0 );
  taxisDefVdate ( taxisID, 1 );
  taxisDefVtime ( taxisID, 2 );
  taxisDefRdate ( taxisID, 3 );
  taxisDefRtime ( taxisID, 4 );
  taxisDefVdateBounds ( taxisID, 5, 6 );
  taxisDefVtimeBounds ( taxisID, 7, 8 );
  taxisDefCalendar ( taxisID, 1 );
  taxisDefTunit ( taxisID, 1 );
  taxisDefNumavg ( taxisID, 1 );

  return taxisID;
}

void defineStream ( int streamID, int vlistID )
{
  streamDefByteorder ( streamID, 1 );
  streamDefCompType  ( streamID, 2 );
  streamDefCompLevel ( streamID, 3 );
  streamDefVlist(streamID, vlistID);
}

int defineVlist ( int gridID, int zaxisID, int taxisID )
{
  int vlistID = CDI_UNDEFID;
  int zaxisID2 = zaxisCreate(ZAXIS_SURFACE, 1);
  int varID1, varID2;

  vlistID = vlistCreate();
  varID1 = vlistDefVar(vlistID, gridID, zaxisID, TIME_VARIABLE);
  varID2 = vlistDefVar(vlistID, gridID, zaxisID2, TIME_VARIABLE);
  vlistDefVarName(vlistID, varID1, "varname1");
  {
    int globfac[] = { 23, 42 };
    vlistDefAttInt(vlistID, varID1, "seer's globule factors", DATATYPE_INT16,
                   2, globfac);
  }
  vlistDefVarName(vlistID, varID2, "varname2");
  vlistDefAttTxt(vlistID, varID2, "txt demo", 6, "banana");
  vlistDefTaxis(vlistID, taxisID);
  return vlistID;
}

int defineInstitute ()
{
  int instID = CDI_UNDEFID;

  instID = institutDef( 0, 0,"MYINSTITUTE", "myInstitute");

  return instID;
}

int defineModel ( int instID )
{
  int modelID = CDI_UNDEFID;

  modelID = modelDef ( instID, 0, "myModel");

  return modelID;
}

177
int modelRun ( MPI_Comm comm )
178
179
180
181
182
183
184
{
  int gridID, zaxisID, taxisID, instID, modelID, vlistID, streamID;

  char * recvBuffer, * sendBuffer;
  int bufferSize, differ;

  pioNamespaceSetActive ( 0 );
185
186
  if (IOMode != PIO_NONE)
    serializeSetMPI();
187
188
189
190
191
192
193
194
195
196
197

  gridID  = defineGrid      ();
  zaxisID = defineZaxis     ();
  taxisID = defineTaxis     ();
  instID  = defineInstitute ();
  modelID = defineModel     ( instID );
  vlistID = defineVlist     ( gridID, zaxisID, taxisID);
  streamID = streamOpenWrite("example.grb", FILETYPE_GRB);
  if ( streamID < 0 ) xabort ( "Could not open file" );
  defineStream ( streamID, vlistID );

198
  reshPackBufferCreate ( &sendBuffer, &bufferSize, &comm );
Thomas Jahns's avatar
Thomas Jahns committed
199
  recvBuffer = xmalloc(bufferSize);
200
201
202
  xmpi(MPI_Sendrecv(sendBuffer, bufferSize, MPI_PACKED, 0, 0,
                    recvBuffer, bufferSize, MPI_PACKED, 0, 0,
                    MPI_COMM_SELF, MPI_STATUS_IGNORE));
203
  pioNamespaceSetActive ( 1 );
204
  reshUnpackResources(recvBuffer, bufferSize, &comm);
205
206
207
208
209
210
211
  free ( recvBuffer );
  reshPackBufferDestroy ( &sendBuffer );

  differ = reshListCompare ( 0, 1 );

  pioNamespaceSetActive ( 0 );
  streamClose(streamID);
212
  return differ;
213
214
215
216
217
218
}

#endif

int main (int argc, char *argv[])
{
219
  int exitCode = 77;
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
#ifdef USE_MPI
  int sizeGlob, pioNamespace;
  MPI_Comm commGlob, commModel;

  MPI_Init(&argc, &argv);
  xmpi ( MPI_Comm_dup ( MPI_COMM_WORLD, &commGlob ));
  xmpi ( MPI_Comm_set_errhandler ( commGlob, MPI_ERRORS_RETURN ));
  xmpi ( MPI_Comm_size ( commGlob, &sizeGlob ));

  if ( sizeGlob != 1 )
      xabort ( "test transition of resource array only with 1 PE." );

  if ( nProcsIO != 1 )
    xabort ( "bad distribution of tasks on PEs" );

235
236
  commModel = pioInit(commGlob, nProcsIO, IOMode, &pioNamespace, 1.0f,
                      cdiPioNoPostCommSetup);
237
238
  pioNamespaceSetActive(pioNamespace);

239
  exitCode = modelRun(commModel);
240

241
  xmpi(MPI_Finalize());
242
243
244
245
#else
  printf ( "Use MPI for this testprogram.\n" );
#endif

246
  return exitCode;
247
248
249
250
251
252
253
254
255
256
257
}

/*
 * Local Variables:
 * c-file-style: "Java"
 * c-basic-offset: 2
 * indent-tabs-mode: nil
 * show-trailing-whitespace: t
 * require-trailing-newline: t
 * End:
 */