aboutsummaryrefslogtreecommitdiff
path: root/src/DumpUtils.c
blob: d44873c402a1edd2f840910d5c61c93ed5aed83f (plain)
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
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
 /*@@
   @file      DumpUtils.c
   @date      Fri Oct 6 2000
   @author    Thomas Radke
   @desc 
              Utility routines for dumping data into HDF5 files.
   @enddesc 
   @version   $Id$
 @@*/


#include <stdlib.h>

#include "cctk.h"
#include "cctk_Version.h"
#include "cctk_Parameters.h"
#include "CactusBase/IOUtil/src/ioGH.h"
#include "CactusBase/IOUtil/src/ioutil_CheckpointRecovery.h"
#include "ioHDF5UtilGH.h"

/* the rcs ID and its dummy function to use it */
static const char *rcsid = "$Id$";
CCTK_FILEVERSION(BetaThorns_IOHDF5Util_DumpUtils_c)


/*@@
   @routine    IOHDF5Util_DumpGH
   @date       Tue Oct 10 2000
   @author     Thomas Radke
   @desc
               Dump all variables of the given GH along with
               the global attributes and parameters
               to the (already opened) HDF5 checkpoint file.
   @enddesc
   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        timers
   @vdesc      array of timers for timing the checkpointing
   @vtype      const int *
   @vio        in
   @endvar
   @var        file
   @vdesc      HDF5 file handle
   @vtype      hid_t
   @vio        in
   @endvar

   @returntype int
   @returndesc
               always 0 (for success)
   @endreturndesc
@@*/
int IOHDF5Util_DumpGH (const cGH *GH,
                       const int *timers,
                       hid_t file)
{
  DECLARE_CCTK_PARAMETERS
  int maxdim;
  int first_vindex, vindex, gindex;
  int timelevel, last_timelevel;
  cGroup gdata;
  int old_out_single;
  ioGH *ioUtilGH;
  ioHDF5Geo_t request;
  /*** FIXME: can CCTK_SyncGroup() have a 'const cGH *' parameter ?? ***/
  union
  {
    const cGH *const_ptr;
    cGH *non_const_ptr;
  } GH_fake_const;


  GH_fake_const.const_ptr = GH;

  /* Get the GH extension for IOUtil */
  ioUtilGH = (ioGH *) CCTK_GHExtension (GH, "IO");

  maxdim = CCTK_MaxDim ();
  request.origin     = (int *) calloc (4 + maxdim, maxdim * sizeof (int));
  request.length     = request.origin + 1 * maxdim;
  request.downsample = request.origin + 2 * maxdim;
  request.actlen     = request.origin + 3 * maxdim;
  request.direction  = request.origin + 4 * maxdim;

  for (vindex = 0; vindex < maxdim; vindex++)
  {
    request.length[vindex]                      = -1;
    request.downsample[vindex]                  = 1;
  }

  /* disable output in single precision */
  old_out_single = ioUtilGH->out_single;
  ioUtilGH->out_single = 0;

  /* sync all active groups */
  for (vindex = CCTK_NumGroups () - 1; vindex >= 0; vindex--)
  {
    if (CCTK_IsImplementationActive (
           CCTK_ImpFromVarI (CCTK_FirstVarIndexI (vindex))))
    {
      CCTK_SyncGroupI (GH_fake_const.non_const_ptr, vindex);
    }
  }

  /* start CP_PARAMETERS_TIMER timer */
  if (timers)
  {
    CCTK_TimerResetI (timers[CP_PARAMETERS_TIMER]);
    CCTK_TimerStartI (timers[CP_PARAMETERS_TIMER]);
  }

  /* Great; Now start dumping away! */
  if (file >= 0)
  {
    /* all GH extension variables and parameter stuff which is not
       specific to any dataset goes into dedicated groups */
    if (out3D_parameters)
    {
      IOHDF5Util_DumpParameters (GH, file);
    }

    IOHDF5Util_DumpGHExtensions (GH, file);
  }

  if (verbose)
  {
    CCTK_INFO ("Dumping variables ...");
  }

  /* stop CP_PARAMETERS_TIMER timer and start CP_VARIABLES_TIMER */
  if (timers)
  {
    CCTK_TimerStopI  (timers[CP_PARAMETERS_TIMER]);
    CCTK_TimerResetI (timers[CP_VARIABLES_TIMER]);
    CCTK_TimerStartI (timers[CP_VARIABLES_TIMER]);
  }

  /* ... now the variables, sorted by groups */
  for (gindex = CCTK_NumGroups () - 1; gindex >= 0; gindex--)
  {
    /* only dump groups which have storage assigned */
    if (CCTK_QueryGroupStorageI (GH, gindex) <= 0)
    {
      continue;
    }

    /* sync the group */
    CCTK_SyncGroupI (GH_fake_const.non_const_ptr, gindex);

    /* dump all timelevels except the oldest (for multi-level groups) */
    CCTK_GroupData (gindex, &gdata);
    last_timelevel = gdata.numtimelevels;
    if (last_timelevel > 1)
    {
      last_timelevel--;
    }
    request.sdim = request.vdim = gdata.dim;

    /* set the hyperslab directions (orthogonal to all axes) */
    memset (request.direction, 0, request.sdim * request.vdim);
    for (vindex = 0; vindex < request.sdim; vindex++)
    {
      request.direction[vindex * request.vdim + vindex] = 1;
    }

    /* loop over all variables in this group */
    first_vindex = CCTK_FirstVarIndexI (gindex);
    for (vindex = first_vindex; vindex < first_vindex + gdata.numvars; vindex++)
    {
      if (verbose && file >= 0)
      {
        CCTK_VInfo (CCTK_THORNSTRING, "  %s", CCTK_VarName (vindex));
      }

      /* loop over all timelevels of this variable */
      for (timelevel = 0; timelevel < last_timelevel; timelevel++)
      {
        IOHDF5Util_DumpVar (GH, vindex, timelevel, &request, file, 0);
      }

    } /* end of loop over all variables */
  } /* end of loop over all groups */

  /* stop CP_VARIABLES_TIMER timer */
  if (timers)
  {
    CCTK_TimerStopI (timers[CP_VARIABLES_TIMER]);
  }

  /* restore output precision flag */
  ioUtilGH->out_single = old_out_single;

  /* free temporary resources */
  free (request.origin);

  return (0);
}


/*@@
   @routine    IOHDF5Util_DumpCommonAttributes
   @date       May 21 1999
   @author     Thomas Radke
   @desc
               Add "Common" attributes to the given dataset, these are:
               <ul>
                 <li> the variable's groupname
                 <li> the grouptype
                 <li> number of timelevels
                 <li> simulation time
                 <li> origin
                 <li> bounding box
                 <li> gridspacings
                 <li> global grid size
               </ul>
   @enddesc
   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        vindex
   @vdesc      CCTK index of the variable to dump
   @vtype      int
   @vio        in
   @endvar
   @var        timelevel
   @vdesc      timelevel of the variable to dump
   @vtype      int
   @vio        in
   @endvar
   @var        global_shape
   @vdesc      global shape of the grid
   @vtype      CCTK_INT []
   @vio        in
   @endvar
   @var        request
   @vdesc      pointer to IO request structure describing the hyperslab
   @vtype      ioHDF5Geo_t *
   @vio        in
   @endvar
   @var        dataset
   @vdesc      the dataset handle where the attributes should be attached to
   @vtype      hid_t
   @vio        in
   @endvar
@@*/
void IOHDF5Util_DumpCommonAttributes (const cGH *GH,
                                      int vindex,
                                      int timelevel,
                                      const CCTK_INT *global_shape,
                                      const ioHDF5Geo_t *request,
                                      hid_t dataset)
{
  DECLARE_CCTK_PARAMETERS
  char *groupname;
  int dim, vdim;
  CCTK_INT attr_int;
  CCTK_REAL *attr_real;
  ioHDF5UtilGH *myGH;
  char coord_system_name[20];


  /* prevent compiler warning about unused parameter */
  timelevel = timelevel;

  /* Get the handle for IOHDF5Util extensions */
  myGH = (ioHDF5UtilGH *) CCTK_GHExtension (GH, "IOHDF5Util");

  /* attributes describing the variable */
  groupname = CCTK_GroupNameFromVarI (vindex);
  WRITE_ATTRIBUTE ("groupname", groupname, dataset,
                   myGH->scalar_dataspace, 0, myGH->IOHDF5_STRING);
  free (groupname);
  attr_int = CCTK_GroupTypeFromVarI (vindex);
  WRITE_ATTRIBUTE ("grouptype", &attr_int, dataset,
                   myGH->scalar_dataspace, 0, IOHDF5_INT);
  attr_int = CCTK_NumTimeLevelsFromVarI (vindex);
  WRITE_ATTRIBUTE ("ntimelevels", &attr_int, dataset,
                   myGH->scalar_dataspace, 0, IOHDF5_INT);
  WRITE_ATTRIBUTE ("global_size", global_shape, dataset,
                   myGH->array_dataspace, request->sdim, IOHDF5_INT);
  WRITE_ATTRIBUTE ("time", &GH->cctk_time, dataset,
                   myGH->scalar_dataspace, 0, IOHDF5_REAL);

  /* attributes describing the underlying grid
     These are only stored for CCTK_GF variables if they are associated
     with coordinates. */
  /* FIXME: This is hardcoded for cartesian coordinate systems.
            A better solution would be to be able to query the coordinate
            system which is associated with the variable. */
  vdim = CCTK_GroupDimFromVarI (vindex);
  sprintf (coord_system_name, "cart%dd", vdim);
  if (CCTK_GroupTypeFromVarI (vindex) == CCTK_GF &&
      CCTK_CoordSystemHandle (coord_system_name) >= 0)
  {
    attr_real = (CCTK_REAL *) malloc (2 * vdim * sizeof (CCTK_REAL));
    for (dim = 0; dim < vdim; dim++)
    {
      CCTK_CoordRange (GH, &attr_real[dim], &attr_real[dim + vdim], dim + 1,
                       NULL, coord_system_name);
    }

    WRITE_ATTRIBUTE ("origin", attr_real, dataset,
                     myGH->array_dataspace, vdim, IOHDF5_REAL);
    WRITE_ATTRIBUTE ("min_ext", attr_real, dataset,
                     myGH->array_dataspace, vdim, IOHDF5_REAL);
    WRITE_ATTRIBUTE ("max_ext", attr_real + vdim, dataset,
                     myGH->array_dataspace, vdim, IOHDF5_REAL);
    WRITE_ATTRIBUTE ("delta", GH->cctk_delta_space, dataset,
                     myGH->array_dataspace, vdim, IOHDF5_REAL);
    free (attr_real);
  }

#if 0
  /* attributes describing the hyperslab */
  /* FIXME: what attributes are really needed here ?? ***/
  if (request)
  {
    attr_real = (CCTK_REAL *) malloc (4 * request->sdim * sizeof (CCTK_REAL));
    for (dim = 0; dim < request->sdim; dim++)
    {
      attr_real[dim + 0*request->sdim] =
        request->origin[request->direction[dim]] *
        GH->cctk_delta_space[request->direction[dim]];
      attr_real[dim + 1*request->sdim] =
        request->origin[dim] * GH->cctk_delta_space[dim];
      attr_real[dim + 2*request->sdim] =
        (request->origin[dim] + request->actlen[dim]-1) *
        GH->cctk_delta_space[dim] * request->downsample[dim];
      attr_real[dim + 3*request->sdim] =
        GH->cctk_delta_space[request->direction[dim]] *
        request->downsample[request->direction[dim]];
    }
    WRITE_ATTRIBUTE ("origin_slab",  attr_real + 0*dim, dataset,
                     myGH->array_dataspace, dim, IOHDF5_REAL);
    WRITE_ATTRIBUTE ("min_ext_slab", attr_real + 1*dim, dataset,
                     myGH->array_dataspace, dim, IOHDF5_REAL);
    WRITE_ATTRIBUTE ("max_ext_slab", attr_real + 2*dim, dataset,
                     myGH->array_dataspace, dim, IOHDF5_REAL);
    WRITE_ATTRIBUTE ("delta_slab",   attr_real + 3*dim, dataset,
                     myGH->array_dataspace, dim, IOHDF5_REAL);
    free (attr_real);
  }
#endif
}


 /*@@
   @routine    IOHDF5Util_DumpParameters
   @date       Thu Jan 20 2000
   @author     Thomas Radke
   @desc
               Collects the parameters of all active implementations
               into a single string and writes it as an attribute
               attached to the CACTUS_PARAMETERS_GROUP group in the HDF5 file.
   @enddesc
   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        file
   @vdesc      the HDF5 file handle
   @vtype      hid_t
   @vio        in
   @endvar
@@*/
void IOHDF5Util_DumpParameters (const cGH *GH,
                                hid_t file)
{
  DECLARE_CCTK_PARAMETERS
  ioHDF5UtilGH *myGH;
  char *parameters;
  hid_t group;


  if (verbose)
  {
    CCTK_INFO ("Dumping GH extensions ...");
  }

  /* Get the parameter string buffer */
  parameters = IOUtil_GetAllParameters (GH);

  if (parameters)
  {
    /* Get the handle for IOHDF5Util extensions */
    myGH = (ioHDF5UtilGH *) CCTK_GHExtension (GH, "IOHDF5Util");

    IOHDF5_ERROR (group = H5Gcreate (file, CACTUS_PARAMETERS_GROUP, 0));
    WRITE_ATTRIBUTE (ALL_PARAMETERS, parameters, group,
                     myGH->scalar_dataspace, 0, myGH->IOHDF5_STRING);
    IOHDF5_ERROR (H5Gclose (group));

    free (parameters);
  }
}


 /*@@
   @routine    IOHDF5Util_DumpGHExtensions
   @date       Thu Jan 20 2000
   @author     Thomas Radke
   @desc
               Attaches the current values of important flesh and IO variables
               as attributes to the GLOBAL_ATTRIBUTES_GROUP group
               in the HDF5 file.
   @enddesc
   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        file
   @vdesc      the HDF5 file handle
   @vtype      hid_t
   @vio        in
   @endvar
@@*/
void IOHDF5Util_DumpGHExtensions (const cGH *GH,
                                  hid_t file)
{
  DECLARE_CCTK_PARAMETERS
  int value;
  hid_t group;
  char buffer[128];
  const char *version;
  ioGH *ioUtilGH;
  ioHDF5UtilGH *myGH;


  if (verbose)
  {
    CCTK_INFO ("Dumping GH extensions ...");
  }

  /* Get the handles for IOUtil and IOHDF5 extensions */
  ioUtilGH = (ioGH *) CCTK_GHExtension (GH, "IO");
  myGH = (ioHDF5UtilGH *) CCTK_GHExtension (GH, "IOHDF5Util");

  IOHDF5_ERROR (group = H5Gcreate (file, GLOBAL_ATTRIBUTES_GROUP, 0));

  value = CCTK_MainLoopIndex ();
  WRITE_ATTRIBUTE ("main_loop_index", &value, group,
                   myGH->scalar_dataspace, 0, H5T_NATIVE_INT);
  value = CCTK_nProcs (GH);
  WRITE_ATTRIBUTE ("nprocs", &value, group,
                   myGH->scalar_dataspace, 0, H5T_NATIVE_INT);
  WRITE_ATTRIBUTE ("ioproc_every", &ioUtilGH->ioproc_every, group,
                   myGH->scalar_dataspace, 0, H5T_NATIVE_INT);
  WRITE_ATTRIBUTE ("unchunked", &ioUtilGH->unchunked, group,
                   myGH->scalar_dataspace, 0, H5T_NATIVE_INT);
  WRITE_ATTRIBUTE ("cctk_time", &GH->cctk_time, group,
                   myGH->scalar_dataspace, 0, IOHDF5_REAL);
  WRITE_ATTRIBUTE ("cctk_iteration", &GH->cctk_iteration, group,
                   myGH->scalar_dataspace, 0, H5T_NATIVE_INT);
  version = CCTK_FullVersion ();
  WRITE_ATTRIBUTE ("Cactus version", version, group,
                   myGH->scalar_dataspace, 0, myGH->IOHDF5_STRING);

  /* add the parameter filename and the creation date
     as file identification attributes */
  if (CCTK_Equals (out_fileinfo, "parameter filename") ||
      CCTK_Equals (out_fileinfo, "all"))
  {
    buffer[0] = 0;
    CCTK_ParameterFilename (sizeof (buffer), buffer);
    WRITE_ATTRIBUTE ("parameter file", buffer, group,
                     myGH->scalar_dataspace, 0, myGH->IOHDF5_STRING);
  }
  if (CCTK_Equals (out_fileinfo, "creation date") ||
      CCTK_Equals (out_fileinfo, "all"))
  {
    buffer[0] = 0;
    Util_CurrentDate (sizeof (buffer), buffer);
    value = strlen (buffer) + 1;
    buffer[value-1] = ' ';
    Util_CurrentTime (sizeof (buffer) - value, buffer + value);
    WRITE_ATTRIBUTE ("creation date", buffer, group,
                     myGH->scalar_dataspace, 0, myGH->IOHDF5_STRING);
  }

  IOHDF5_ERROR (H5Gclose (group));
}


/*@@
   @routine    IOHDF5Util_DataType
   @author     Thomas Radke
   @date       Mon 11 June 2001
   @desc
               Returns the HDF5 datatype for a given CCTK datatype
   @enddesc

   @var        myGH
   @vdesc      Pointer to IOHDF5Util's GH extensions
   @vtype      ioHDF5UtilGH *
   @vio        in
   @endvar
   @var        cctk_type
   @vdesc      CCTK datatype
   @vtype      int
   @vio        in
   @endvar

   @returntype hid_t
   @returndesc
               the appropriate HDF5 datatype for success, or -1 otherwise
   @endreturndesc
@@*/
hid_t IOHDF5Util_DataType (const ioHDF5UtilGH *myGH, int cctk_type)
{
  hid_t retval;


  switch (cctk_type)
  {
    case CCTK_VARIABLE_CHAR:      retval = IOHDF5_CHAR; break;
    case CCTK_VARIABLE_INT:       retval = IOHDF5_INT; break;
    case CCTK_VARIABLE_REAL:      retval = IOHDF5_REAL; break;
    case CCTK_VARIABLE_COMPLEX:   retval = myGH->IOHDF5_COMPLEX; break;
#ifdef CCTK_INT2
    case CCTK_VARIABLE_INT2:      retval = IOHDF5_INT2; break;
#endif
#ifdef CCTK_INT4
    case CCTK_VARIABLE_INT4:      retval = IOHDF5_INT4; break;
#endif
#ifdef CCTK_INT8
    case CCTK_VARIABLE_INT8:      retval = IOHDF5_INT8; break;
#endif
#ifdef CCTK_REAL4
    case CCTK_VARIABLE_REAL4:     retval = IOHDF5_REAL4; break;
    case CCTK_VARIABLE_COMPLEX8:  retval = myGH->IOHDF5_COMPLEX8; break;
#endif
#ifdef CCTK_REAL8
    case CCTK_VARIABLE_REAL8:     retval = IOHDF5_REAL8; break;
    case CCTK_VARIABLE_COMPLEX16: retval = myGH->IOHDF5_COMPLEX16; break;
#endif
#ifdef CCTK_REAL16
    case CCTK_VARIABLE_REAL16:    retval = IOHDF5_REAL16; break;
    case CCTK_VARIABLE_COMPLEX32: retval = myGH->IOHDF5_COMPLEX32; break;
#endif

    default: CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                         "Unsupported CCTK variable datatype %d", cctk_type);
             retval = -1;
  }

  return (retval);
}