aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetIOHDF5/src/Input.cc
blob: 6495eb276d67e74e743c2a4c3fc64cfef1521da0 (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
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
#include <cassert>
#include <cstring>
#include <sstream>
#include <vector>

#include "util_Table.h"
#include "cctk.h"
#include "cctk_Parameters.h"

#include "CarpetIOHDF5.hh"
#include "CactusBase/IOUtil/src/ioGH.h"
#include "CactusBase/IOUtil/src/ioutil_CheckpointRecovery.h"

#include "defs.hh"


namespace CarpetIOHDF5
{

using namespace std;
using namespace Carpet;

// structure describing a patch as a single dataset of an HDF5 file to read from
typedef struct {
  string patchname;
  int vindex;
  int map;
  int mglevel;
  int reflevel;
  int timestep;
  int timelevel;
  int component;
  int rank;
  ivect iorigin;
  vector<hsize_t> shape;
} patch_t;

// structure describing the contents of an HDF5 file to read from
typedef struct {
  char *filename;
  hid_t file;
  list<patch_t> patches;
} file_t;

// structure describing a set of HDF5 files
typedef struct {
  string setname;
  string basefilename;
  int first_ioproc;
  vector<file_t> files;            // [nioprocs]

  int nioprocs;
  int num_mglevels;
  int num_reflevels;
  int cctk_iteration;
  int main_loop_index;
  CCTK_REAL global_time;
  CCTK_REAL delta_time;
  vector<CCTK_REAL> mgleveltimes;  // [num_mglevels*num_reflevels]

  vector<grid_structure_t> grid_structure; // [maps]
} fileset_t;

// list of checkpoint/filereader files
static list<fileset_t> filesets;

// number of reflevels in the checkpoint
static int num_reflevels = -1;


static list<fileset_t>::iterator OpenFileSet (const cGH* const cctkGH,
                                              const string setname,
                                              const char *basefilename,
                                              int called_from);
static void ReadMetadata (fileset_t& fileset, hid_t file);
static herr_t BrowseDatasets (hid_t group, const char *objectname, void *arg);
static int ReadVar (const cGH* const cctkGH,
                    hid_t file,
                    list<patch_t>::const_iterator patch,
                    vector<ibset> &bboxes_read,
                    bool in_recovery);


//////////////////////////////////////////////////////////////////////////////
// Register with the Cactus Recovery Interface
//////////////////////////////////////////////////////////////////////////////
int CarpetIOHDF5_RecoverParameters ()
{
  return IOUtil_RecoverParameters (Recover, ".h5", "HDF5");
}

//////////////////////////////////////////////////////////////////////////////
// Recover the grid structure
//////////////////////////////////////////////////////////////////////////////
void CarpetIOHDF5_RecoverGridStructure (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  
  fileset_t & fileset = * filesets.begin();
  
  // Abort with an error if there is no grid structure in the
  // checkpoint file
  assert (fileset.grid_structure.size() == maps);
  
  for (int m = 0; m < maps; ++ m) {
    grid_structure_t const & grid_structure = fileset.grid_structure.at(m);
    
    int const rls = grid_structure.bbss.size();
    assert (grid_structure.obss.size() == rls);
    
    vector <vector <ibbox> >  bbss = grid_structure.bbss;
    vector <vector <bbvect> > obss = grid_structure.obss;
    vector <vector <int> >    pss (rls);
    
    for (int rl = 0; rl < rls; ++ rl) {
      
      vector <ibbox>  & bbs = bbss.at(rl);
      vector <bbvect> & obs = obss.at(rl);
      vector <int>    & ps  = pss.at(rl);
      
      // Make multiprocessor aware
      Carpet::SplitRegions (cctkGH, bbs, obs, ps);
      
    } // for rl
    
    // Make multigrid aware
    vector <vector <vector <ibbox> > > bbsss;
    Carpet::MakeMultigridBoxes (cctkGH, bbss, obss, bbsss);
    
    // Regrid
    Recompose (cctkGH, m, bbsss, obss, pss, false);
    
  } // for m
  
  PostRecompose ();
}

//////////////////////////////////////////////////////////////////////////////
// Overwrite the "CarpetRegrid::refinement_levels"
// with the number of levels given in the checkpoint file
//
// Note that this has to be done after parameter recovery in order to have
// any effect of steering "CarpetRegrid::refinement_levels".
//////////////////////////////////////////////////////////////////////////////
int CarpetIOHDF5_SetNumRefinementLevels ()
{
  DECLARE_CCTK_PARAMETERS;

  if (num_reflevels > 0) {
    if (not CCTK_Equals (verbose, "none")) {
      char *buffer = CCTK_ParameterValString ("refinement_levels",
                                              "CarpetRegrid");
      assert (buffer);
      CCTK_VInfo (CCTK_THORNSTRING, "Using %i reflevels from checkpoint file. "
                  "Ignoring value '%s' in parameter file.",
                  num_reflevels, buffer);
      free (buffer);
    }

    char buffer[32];
    snprintf (buffer, sizeof (buffer), "%d", num_reflevels);
    int const retval = CCTK_ParameterSet ("refinement_levels", "CarpetRegrid",
                                          buffer);
    assert (retval == 0);
  }
  
  return 0;
}


//////////////////////////////////////////////////////////////////////////////
// close all open checkpoint/filereader files after recovering grid variables
//////////////////////////////////////////////////////////////////////////////
void CarpetIOHDF5_CloseFiles (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  int error_count = 0;


  for (list<fileset_t>::const_iterator set = filesets.begin();
       set != filesets.end(); set++) {
    for (unsigned int i = 0; i < set->files.size(); i++) {
      if (set->files[i].file >= 0) {

        if (CCTK_Equals (verbose, "full")) {
          CCTK_VInfo (CCTK_THORNSTRING, "closing file '%s' after recovery",
                      set->files[i].filename);
        }

        free (set->files[i].filename);
        HDF5_ERROR (H5Fclose (set->files[i].file));
      }
    }
  }

  filesets.clear();
}


//////////////////////////////////////////////////////////////////////////////
// Top-level recovery/filereader routine
// Opens the checkpoint/data file on the first time through,
// recovers parameters and variables
//////////////////////////////////////////////////////////////////////////////
int Recover (cGH* cctkGH, const char *basefilename, int called_from)
{
  int error_count = 0;
  DECLARE_CCTK_PARAMETERS;


  assert (called_from == CP_RECOVER_PARAMETERS or
          called_from == CP_RECOVER_DATA or
          called_from == FILEREADER_DATA);
  const bool in_recovery = called_from == CP_RECOVER_PARAMETERS or
                           called_from == CP_RECOVER_DATA;

  const char* const dir = in_recovery ? recover_dir : filereader_ID_dir;
  string setname (dir);
  if (called_from == FILEREADER_DATA) {
    setname.append (basefilename);
  }

  list<fileset_t>::iterator fileset = filesets.begin();
  while (fileset != filesets.end()) {
    if (fileset->setname == setname) {
      break;
    }
    fileset++;
  }
  if (fileset == filesets.end()) {
    assert (called_from == CP_RECOVER_PARAMETERS or
            called_from == FILEREADER_DATA);

    fileset = OpenFileSet (cctkGH, setname, basefilename, called_from);
    if (fileset == filesets.end()) {
      return (-1);
    }

    if (called_from == CP_RECOVER_PARAMETERS) {
      // return here if only parameters should be recovered
      // (come back later with called_from == CP_RECOVER_DATA)
      return (1);
    }
  }

  // set global Cactus/Carpet variables
  if (in_recovery) {

    if (use_grid_structure_from_checkpoint) {
      // recover the grid structure only once
      static bool is_first = true;
      if (is_first) {
        is_first = false;
        CarpetIOHDF5_RecoverGridStructure (cctkGH);
      }
    }

    global_time = fileset->global_time;
    delta_time = fileset->delta_time;
    CCTK_SetMainLoopIndex (fileset->main_loop_index);

    cctkGH->cctk_iteration = fileset->cctk_iteration;
    int const idx = mglevel*fileset->num_reflevels + reflevel;
    cctkGH->cctk_time = fileset->mgleveltimes.at(idx);
  }

  if (not CCTK_Equals (verbose, "none")) {
    CCTK_VInfo (CCTK_THORNSTRING,
                "reading grid variables on mglevel %d reflevel %d",
                mglevel, reflevel);
  }

  // create a bbox set for each active timelevel of all variables
  // to mark how much has been read already
  const int numvars = CCTK_NumVars ();
  vector<vector<bool> > read_completely(numvars);
  vector<vector<vector<ibset> > > bboxes_read (numvars);
  for (unsigned int vindex = 0; vindex < bboxes_read.size(); vindex++) {
    const int timelevels = CCTK_ActiveTimeLevelsVI (cctkGH, vindex);
    read_completely[vindex].resize (timelevels, false);
    bboxes_read[vindex].resize (timelevels);
    for (int tl = 0; tl < timelevels; tl++) {
      bboxes_read[vindex][tl].resize (Carpet::maps);
    }
  }

  // query for groups which have the 'CHECKPOINT = "no"' option set
  // Such groups are not checked against being read completely.
  vector<bool> not_checkpointed(numvars);
  if (in_recovery) {
    for (unsigned int vindex = 0; vindex < not_checkpointed.size(); vindex++) {
      int gindex = CCTK_GroupIndexFromVarI (vindex);
      int tagstable = CCTK_GroupTagsTableI (gindex);
      int const len = Util_TableGetString (tagstable, 0, NULL, "checkpoint");
      if (len > 0) {
        char* value = new char[len + 1];
        Util_TableGetString (tagstable, len + 1, value, "checkpoint");
        if (len == sizeof ("no") - 1 and CCTK_Equals (value, "no")) {
          not_checkpointed[vindex] = true;
        }
        delete[] value;
      }
    }
  }

  // create a bbox set for each group to list how much needs to be read
  const int numgroups = CCTK_NumGroups ();
  vector<vector<ibset> > group_bboxes (numgroups);
  for (unsigned int gindex = 0; gindex < group_bboxes.size(); gindex++) {
    group_bboxes[gindex].resize (Carpet::maps);
    const int grouptype = CCTK_GroupTypeI (gindex);
    BEGIN_MAP_LOOP (cctkGH, grouptype) {
      struct arrdesc& data = arrdata.at(gindex).at(Carpet::map);

      BEGIN_LOCAL_COMPONENT_LOOP (cctkGH, grouptype) {
        if (grouptype == CCTK_GF or (mglevel == 0 and reflevel == 0)) {
          group_bboxes[gindex][Carpet::map] |=
            data.dd->boxes.at(mglevel).at(reflevel).at(component).exterior;
        }
      } END_LOCAL_COMPONENT_LOOP;
    } END_MAP_LOOP;
  }

  const ioGH* ioUtilGH = (const ioGH*) CCTK_GHExtension (cctkGH, "IO");

  // loop over all input files of this set
  for (unsigned int i = 0; i < fileset->files.size(); i++) {

    const int file_idx = (i + fileset->first_ioproc) % fileset->nioprocs;
    file_t& file = fileset->files[file_idx];

    // open the file (if it hasn't been already) and read its metadata
    if (file.file < 0) {
      file.filename =
        IOUtil_AssembleFilename (NULL, fileset->basefilename.c_str(),
                                 "", ".h5", called_from, file_idx, 0);
      assert (file.filename);
      HDF5_ERROR (file.file = H5Fopen (file.filename, H5F_ACC_RDONLY,
                                       H5P_DEFAULT));

      if (CCTK_Equals (verbose, "full")) {
        CCTK_VInfo (CCTK_THORNSTRING, "opening %s file '%s'",
                    in_recovery ? "checkpoint" : "input", file.filename);
      }

      // browse through all datasets contained in this file
      HDF5_ERROR (H5Giterate (file.file, "/", NULL, BrowseDatasets, &file));
    }
    assert (file.patches.size() > 0);

    // some optimisation for the case when all processors recover
    // from a single chunked checkpoint file:
    // reorder the list so that processor-local components come first
    if (fileset->nioprocs == 1) {
      list<patch_t>::iterator patch = file.patches.begin();
      while (patch != file.patches.end()) {
        if (patch->component == dist::rank()) {
          file.patches.push_front (*patch);
          patch = file.patches.erase (patch);
        } else {
          patch++;
        }
      }
    }

    // now loop over all patches contained in this file
    for (list<patch_t>::const_iterator patch = file.patches.begin();
         patch != file.patches.end();
         patch++) {

      // only recover grid variables for the current mglevel/reflevel
      if (patch->mglevel != mglevel or patch->reflevel != reflevel) {
        continue;
      }

      // When called by the filereader, IOUtil will set the 'do_inVars' vector
      // if the parameter 'IO::filereader_ID_vars' has been set.
      // In this parameter the user can specify individual variables to be read
      // as well as (optionally) a specific timestep.
      //
      // If 'IO::filereader_ID_vars' is left empty, 'do_inVars' will be NULL
      // which implicitely selects all variables for input.
      // If 'do_inVars[vindex]' is 0, variable 'vindex' was not selected.
      // If 'do_inVars[vindex]' is positive, a specific timestep of 'vindex'
      // was selected.
      // If 'do_inVars[vindex]' is negative, 'vindex' was selected with no
      // specific timestep; all timesteps will be read, the last one wins.
      if (ioUtilGH->do_inVars and
          ioUtilGH->do_inVars[patch->vindex] >= 0 and
          ioUtilGH->do_inVars[patch->vindex] != patch->timestep + 1) {
        continue;
      }

      // check the timelevel
      if ((unsigned int) patch->timelevel >=
          read_completely.at(patch->vindex).size()) {
        CCTK_VWarn (3, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "Ignoring dataset '%s' (invalid timelevel %d)",
                    patch->patchname.c_str(), patch->timelevel);
        continue;
      }

      // actually read the patch
      if (not read_completely.at(patch->vindex).at(patch->timelevel)) {
        ReadVar (cctkGH, file.file, patch,
                 bboxes_read.at(patch->vindex).at(patch->timelevel),
                 in_recovery);

        // update the read_completely entry
        const int gindex = CCTK_GroupIndexFromVarI (patch->vindex);
        read_completely.at(patch->vindex).at(patch->timelevel) =
          bboxes_read.at(patch->vindex).at(patch->timelevel) ==
          group_bboxes.at(gindex);
      }
    }

    // check if all variables have been read completely already
    bool all_done = true;
    for (unsigned int vindex = 0; vindex < read_completely.size(); vindex++) {

      // skip all variables which aren't expected to be recovered
      if (not_checkpointed[vindex] or
          (CCTK_GroupTypeFromVarI (vindex) != CCTK_GF and reflevel > 0)) {
        continue;
      }

      for (unsigned int tl = 0; tl < read_completely[vindex].size(); tl++) {
        all_done &= read_completely[vindex][tl];
      }
    }
    if (all_done) {
      break;
    }
  }

  // check that all variables have been read completely on this mglevel/reflevel
  int num_incomplete = 0;
  for (unsigned int vindex = 0; vindex < read_completely.size(); vindex++) {

    if (CCTK_GroupTypeFromVarI (vindex) != CCTK_GF and reflevel > 0) {
      continue;
    }

    for (unsigned int tl = 0; tl < read_completely[vindex].size(); tl++) {
      if (called_from == FILEREADER_DATA and not
          (ioUtilGH->do_inVars and ioUtilGH->do_inVars[vindex])) {
        continue;
      }

      if (not read_completely[vindex][tl]) {
        // check if the variable has been read partially
        size_t size = 0;
        for (int map = 0; map < Carpet::maps; map++) {
          size += bboxes_read[vindex][tl][map].size();
        }
        char* fullname = CCTK_FullName (vindex);
        if (size == 0) {
          if (not_checkpointed[vindex]) {
            CCTK_VWarn (4, __LINE__, __FILE__, CCTK_THORNSTRING,
                        "variable '%s' timelevel %d has not been read "
                        "(variable has option tag \"CHECKPOINT = 'no'\")",
                        fullname, tl);
          } else {
            CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                        "variable '%s' timelevel %d has not been read",
                        fullname, tl);
          }
        } else {
          CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                      "variable '%s' timelevel %d has been read only partially",
                      fullname, tl);
          num_incomplete++;
        }
        free (fullname);
      }
    }
  }
  if (num_incomplete) {
    CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                "%d variables on mglevel %d reflevel %d have been read "
                "only partially", num_incomplete, mglevel, reflevel);
  }

  if (in_recovery and not CCTK_Equals (verbose, "none")) {
    CCTK_VInfo (CCTK_THORNSTRING,
                "restarting simulation on mglevel %d reflevel %d "
                "at iteration %d (physical time %g)", mglevel, reflevel,
                cctkGH->cctk_iteration, (double) cctkGH->cctk_time);
  }

  return (0);
}


//////////////////////////////////////////////////////////////////////////////
// Open a set of checkpoint/filereader files
//////////////////////////////////////////////////////////////////////////////
static list<fileset_t>::iterator OpenFileSet (const cGH* const cctkGH,
                                              const string setname,
                                              const char *basefilename,
                                              int called_from)
{
  file_t file;
  fileset_t fileset;
  int error_count = 0;
  DECLARE_CCTK_PARAMETERS;

  // first try to open a chunked file written on this processor
  // (note that dist::rank() cannot be called yet during RECOVER_PARAMETERS)
  fileset.first_ioproc = CCTK_MyProc (cctkGH);
  file.filename = IOUtil_AssembleFilename (NULL, basefilename, "", ".h5",
                                           called_from, fileset.first_ioproc, 0);

  // try to open the file (prevent HDF5 error messages if it fails)
  H5E_BEGIN_TRY {
    file.file = H5Fopen (file.filename, H5F_ACC_RDONLY, H5P_DEFAULT);
  } H5E_END_TRY;

  // if that failed, try a chunked file written on processor 0
  // (which always is an I/O proc)
  if (file.file < 0) {
    free (file.filename);
    fileset.first_ioproc = 0;
    file.filename = IOUtil_AssembleFilename (NULL, basefilename, "", ".h5",
                                             called_from, fileset.first_ioproc, 0);
    H5E_BEGIN_TRY {
      file.file = H5Fopen (file.filename, H5F_ACC_RDONLY, H5P_DEFAULT);
    } H5E_END_TRY;
  }

  // if that still failed, try an unchunked file
  // (which is always written on processor 0)
  if (file.file < 0) {
    free (file.filename);
    file.filename = IOUtil_AssembleFilename (NULL, basefilename, "", ".h5",
                                             called_from, fileset.first_ioproc, 1);
    H5E_BEGIN_TRY {
      file.file = H5Fopen (file.filename, H5F_ACC_RDONLY, H5P_DEFAULT);
    } H5E_END_TRY;
  }

  // return if no valid checkpoint could be found
  if (file.file < 0) {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "No valid HDF5 file with basename '%s' found", basefilename);
    free (file.filename);
    return (filesets.end());
  }

  if (CCTK_Equals (verbose, "full")) {
    CCTK_VInfo (CCTK_THORNSTRING, "opening %s file '%s'",
                called_from == CP_RECOVER_PARAMETERS ?
                "checkpoint" : "input", file.filename);
  }

  // read all the metadata information
  ReadMetadata (fileset, file.file);

  // browse through all datasets contained in this file
  HDF5_ERROR (H5Giterate (file.file, "/", NULL, BrowseDatasets, &file));
  assert (file.patches.size() > 0);

  // recover parameters
  if (called_from == CP_RECOVER_PARAMETERS) {
    if (not CCTK_Equals (verbose, "none")) {
      CCTK_VInfo (CCTK_THORNSTRING, "Recovering parameters from checkpoint "
                  "file '%s'", file.filename);
    }
    hid_t dataset;
    HDF5_ERROR (dataset = H5Dopen (file.file,
                                   METADATA_GROUP "/" ALL_PARAMETERS));
    char *parameters = new char[H5Dget_storage_size (dataset) + 1];
    HDF5_ERROR (H5Dread (dataset, H5T_NATIVE_CHAR, H5S_ALL, H5S_ALL,
                         H5P_DEFAULT, parameters));
    HDF5_ERROR (H5Dclose (dataset));
    IOUtil_SetAllParameters (parameters);
    delete[] parameters;

    num_reflevels = fileset.num_reflevels;
  }

  // allocate and initialise the list of input files for this set
  fileset.files.resize (fileset.nioprocs);
  for (int i = 0; i < fileset.nioprocs; i++) {
    fileset.files[i].file = -1;
  }
  fileset.files[fileset.first_ioproc] = file;
  fileset.setname = setname;
  fileset.basefilename = basefilename;

  // add the new fileset to the list of sets and return an iterator to it
  filesets.push_front (fileset);
  return (filesets.begin());
}


//////////////////////////////////////////////////////////////////////////////
// Read the metadata for a set of input files
//////////////////////////////////////////////////////////////////////////////
static void ReadMetadata (fileset_t& fileset, hid_t file)
{
  int error_count = 0;
  DECLARE_CCTK_PARAMETERS;

  fileset.nioprocs = 1;
  hid_t metadata, attr, dataspace;
  H5E_BEGIN_TRY {
    metadata = H5Gopen (file, METADATA_GROUP);
  } H5E_END_TRY;
  if (metadata < 0) {
    // no metadata at all - this must be old-fashioned data output file
    return;
  }

  // in old times (before June 2005) the metadata attributes were attached
  // to the parameters dataset which was stored in checkpoint files only
  H5E_BEGIN_TRY {
    attr = H5Aopen_name (metadata, "nioprocs");
  } H5E_END_TRY;
  const bool is_old_fashioned_file = attr < 0;
  if (is_old_fashioned_file) {
    HDF5_ERROR (H5Gclose (metadata));
    HDF5_ERROR (metadata = H5Dopen (file,
                                    METADATA_GROUP "/" ALL_PARAMETERS));
  } else {
    HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &fileset.nioprocs));
    HDF5_ERROR (H5Aclose (attr));
  }
  HDF5_ERROR (attr = H5Aopen_name (metadata, "carpet_reflevels"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &fileset.num_reflevels));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (attr = H5Aopen_name (metadata, "numberofmgtimes"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &fileset.num_mglevels));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (attr = H5Aopen_name (metadata, "GH$iteration"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &fileset.cctk_iteration));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (attr = H5Aopen_name (metadata, "main loop index"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &fileset.main_loop_index));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (attr = H5Aopen_name (metadata, "carpet_global_time"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_DOUBLE, &fileset.global_time));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (attr = H5Aopen_name (metadata, "carpet_delta_time"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_DOUBLE, &fileset.delta_time));
  HDF5_ERROR (H5Aclose (attr));

  // Read grid structure if it is present
  hid_t dataset;
  H5E_BEGIN_TRY {
    dataset = H5Dopen (metadata, GRID_STRUCTURE);
  } H5E_END_TRY;
  if (dataset >= 0) {
    vector<char> gs_cstr (H5Dget_storage_size (dataset) + 1);
    HDF5_ERROR (H5Dread (dataset, H5T_NATIVE_CHAR, H5S_ALL, H5S_ALL,
                         H5P_DEFAULT, &gs_cstr.front()));
    HDF5_ERROR (H5Dclose (dataset));
    istringstream gs_buf (&gs_cstr.front());
    gs_buf >> fileset.grid_structure;
  }

  fileset.mgleveltimes.resize (fileset.num_mglevels * fileset.num_reflevels);
  for (int i = 0; i < fileset.num_mglevels; i++) {
    char buffer[32];

    snprintf (buffer, sizeof (buffer), "mgleveltimes %d", i);
    HDF5_ERROR (attr = H5Aopen_name (metadata, buffer));
    HDF5_ERROR (dataspace = H5Aget_space (attr));
    assert (H5Sget_simple_extent_npoints (dataspace) == fileset.num_reflevels);
    HDF5_ERROR (H5Sclose (dataspace));
    HDF5_ERROR (H5Aread (attr, H5T_NATIVE_DOUBLE,
                         &fileset.mgleveltimes[i * fileset.num_reflevels]));
    HDF5_ERROR (H5Aclose (attr));
  }

  if (is_old_fashioned_file) {
    HDF5_ERROR (H5Dclose (metadata));
  } else {
    HDF5_ERROR (H5Gclose (metadata));
  }
}


//////////////////////////////////////////////////////////////////////////////
// Browses through all the datasets in a checkpoint/data file
// and stores this information in a list of patch_t objects
//////////////////////////////////////////////////////////////////////////////
static herr_t BrowseDatasets (hid_t group, const char *objectname, void *arg)
{
  int error_count = 0;
  file_t *file = (file_t *) arg;
  patch_t patch;
  hid_t dataset, dataspace, attr, attrtype;
  H5G_stat_t object_info;


  // skip everything in the metadata group
  if (strcmp (objectname, METADATA_GROUP) == 0) {
    return (0);
  }
  // we are interested only in datasets
  HDF5_ERROR (H5Gget_objinfo (group, objectname, 0, &object_info));
  if (object_info.type != H5G_DATASET) {
    return (0);
  }

  HDF5_ERROR (dataset = H5Dopen (group, objectname));
  HDF5_ERROR (dataspace = H5Dget_space (dataset));
  patch.rank = H5Sget_simple_extent_ndims (dataspace);
  assert (patch.rank > 0 and patch.rank <= ::dim);
  patch.shape.resize (patch.rank);
  HDF5_ERROR (H5Sget_simple_extent_dims (dataspace, &patch.shape[0], NULL));
  HDF5_ERROR (H5Sclose (dataspace));
  HDF5_ERROR (attr = H5Aopen_name (dataset, "name"));
  HDF5_ERROR (attrtype = H5Aget_type (attr));
  size_t length = H5Tget_size (attrtype);
  assert (length > 0);
  vector<char> varname(length + 1);
  HDF5_ERROR (H5Aread (attr, attrtype, &varname[0]));
  HDF5_ERROR (H5Tclose (attrtype));
  HDF5_ERROR (H5Aclose (attr));
  patch.vindex = CCTK_VarIndex (&varname[0]);
  HDF5_ERROR (attr = H5Aopen_name (dataset, "carpet_mglevel"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &patch.mglevel));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (attr = H5Aopen_name (dataset, "level"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &patch.reflevel));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (attr = H5Aopen_name (dataset, "timestep"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &patch.timestep));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (attr = H5Aopen_name (dataset, "group_timelevel"));
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &patch.timelevel));
  HDF5_ERROR (H5Aclose (attr));
  // in old days (before February 2005)
  // timelevels were stored as negative numbers
  patch.timelevel = abs (patch.timelevel);
  HDF5_ERROR (attr = H5Aopen_name (dataset, "iorigin"));
  HDF5_ERROR (dataspace = H5Aget_space (attr));
  assert (H5Sget_simple_extent_npoints (dataspace) == patch.rank);
  HDF5_ERROR (H5Sclose (dataspace));
  patch.iorigin = 0;
  HDF5_ERROR (H5Aread (attr, H5T_NATIVE_INT, &patch.iorigin[0]));
  HDF5_ERROR (H5Aclose (attr));
  HDF5_ERROR (H5Dclose (dataset));

  // try to obtain the map and component number from the patch's name
  // (cleaner way would be to store attributes with the dataset)
  patch.map = 0;
  const char* map_string = strstr (objectname, " m=");
  if (map_string) {
    sscanf (map_string, " m=%d", &patch.map);
  }
  patch.component = -1;
  const char* component_string = strstr (objectname, " c=");
  if (component_string) {
    sscanf (component_string, " c=%d", &patch.component);
  }

  // add this patch to our list
  if (patch.vindex >= 0 and patch.vindex < CCTK_NumVars ()) {
    patch.patchname = objectname;
    file->patches.push_back (patch);
  } else {
    CCTK_VWarn (3, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Ignoring dataset '%s' (invalid variable name)", objectname);
  }

  return (0);
}


//////////////////////////////////////////////////////////////////////////////
// Reads a single grid variable from a checkpoint/data file
//////////////////////////////////////////////////////////////////////////////
static int ReadVar (const cGH* const cctkGH,
                    hid_t file,
                    list<patch_t>::const_iterator patch,
                    vector<ibset> &bboxes_read,
                    bool in_recovery)
{
  int error_count = 0;
  DECLARE_CCTK_PARAMETERS;

  const int gindex = CCTK_GroupIndexFromVarI (patch->vindex);
  assert (gindex >= 0 and (unsigned int) gindex < Carpet::arrdata.size());
  cGroup group;
  CCTK_GroupData (gindex, &group);


  // grid arrays and scalars are read once on the coarsest reflevel
  assert (group.grouptype == CCTK_GF or (mglevel == 0 and reflevel == 0));

  // Check for storage
  if (not CCTK_QueryGroupStorageI(cctkGH, gindex)) {
    char *fullname = CCTK_FullName (patch->vindex);
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Cannot input variable '%s' (no storage)", fullname);
    free (fullname);
    return 0;
  }

  // filereader reads the current timelevel
  int timelevel = in_recovery ? patch->timelevel : 0;

  // get dataset dimensions
  assert (group.dim == (group.grouptype == CCTK_SCALAR ?
                        patch->rank-1 : patch->rank));
  ivect shape(1);
  int elems = 1;
  for (int i = 0; i < patch->rank; i++) {
    shape[i] = patch->shape[patch->rank-i-1];
    elems   *= shape[i];
  }
  const hid_t datatype = CCTKtoHDF5_Datatype (cctkGH, group.vartype, 0);

  const ivect stride = group.grouptype == CCTK_GF ?
                       maxspacereflevelfact/spacereflevelfact : 1;
  ivect lower = patch->iorigin * stride;
  ivect upper = lower + (shape - 1) * stride;

  // Traverse all local components on all maps
  hid_t filespace = -1, dataset = -1;
  BEGIN_MAP_LOOP (cctkGH, group.grouptype) {

    // skip this dataset if it belongs to another map
    if (group.grouptype == CCTK_GF and patch->map != Carpet::map) {
      continue;
    }

    struct arrdesc& data = arrdata.at(gindex).at(Carpet::map);

    BEGIN_LOCAL_COMPONENT_LOOP (cctkGH, group.grouptype) {

      // reset bounds for DISTRIB=CONST variables
      if (group.disttype == CCTK_DISTRIB_CONSTANT) {
        assert (group.grouptype==CCTK_SCALAR or group.grouptype==CCTK_ARRAY);
        const int newlower = lower[patch->rank-1] +
          (upper[patch->rank-1]-lower[patch->rank-1]+1)*
          (data.hh->processors().at(reflevel).at(component));
        const int newupper = upper[patch->rank-1] +
          (upper[patch->rank-1]-lower[patch->rank-1]+1)*
          (data.hh->processors().at(reflevel).at(component));
        lower[patch->rank-1] = newlower;
        upper[patch->rank-1] = newupper;
      }
      const ibbox filebox(lower, upper, stride);

      ibbox& interior_membox =
        data.dd->boxes.at(mglevel).at(reflevel).at(component).interior;

      // skip this dataset if it doesn't overlap with this component's interior
      const ibbox interior_overlap = interior_membox & filebox;
      if (interior_overlap.empty()) {
        continue;
      }

      if (CCTK_Equals (verbose, "full")) {
        char *fullname = CCTK_FullName (patch->vindex);
        CCTK_VInfo (CCTK_THORNSTRING, "  reading '%s' from dataset '%s'",
                    fullname, patch->patchname.c_str());
        free (fullname);
      }

      // get the overlap with this component's exterior
      ibbox& membox =
        data.dd->boxes.at(mglevel).at(reflevel).at(component).exterior;
      const ibbox overlap = membox & filebox;
      bboxes_read.at(Carpet::map) |= overlap;

      // calculate hyperslab selection parameters
#if (H5_VERS_MAJOR == 1 && H5_VERS_MINOR >= 6 && H5_VERS_RELEASE >= 4)
      hsize_t memorigin[dim], fileorigin[dim];
#else
      hssize_t memorigin[dim], fileorigin[dim];
#endif
      hsize_t memdims[dim], count[dim];
      for (int i = 0; i < patch->rank; i++) {
        memdims[patch->rank-i-1] =
          (membox.upper() - membox.lower())[i] / stride[i] + 1;
        count[patch->rank-i-1] =
          (overlap.upper() - overlap.lower())[i] / stride[i] + 1;
        memorigin[patch->rank-i-1] =
          (overlap.lower() - membox.lower())[i] / stride[i];
        fileorigin[patch->rank-i-1] =
          group.disttype == CCTK_DISTRIB_CONSTANT ?
          0 : (overlap.lower() - lower)[i] / stride[i];
      }

      // open the dataset on the first time through
      if (dataset < 0) {
        HDF5_ERROR (dataset = H5Dopen (file, patch->patchname.c_str()));
        HDF5_ERROR (filespace = H5Dget_space (dataset));
      }

      // read the hyperslab
      hid_t memspace;
      HDF5_ERROR (memspace  = H5Screate_simple (patch->rank, memdims, NULL));
      HDF5_ERROR (H5Sselect_hyperslab (filespace, H5S_SELECT_SET, fileorigin,
                                       NULL, count, NULL));
      HDF5_ERROR (H5Sselect_hyperslab (memspace, H5S_SELECT_SET, memorigin,
                                       NULL, count, NULL));
      HDF5_ERROR (H5Dread (dataset, datatype, memspace, filespace, H5P_DEFAULT,
                           cctkGH->data[patch->vindex][timelevel]));
      HDF5_ERROR (H5Sclose (memspace));

    } END_LOCAL_COMPONENT_LOOP;

    if (in_recovery) {
      data.tt->set_time (reflevel, mglevel,
                         ((cctkGH->cctk_time - cctk_initial_time)
                          / (delta_time * mglevelfact)) );
    }

  } END_MAP_LOOP;

  if (dataset >= 0) {
    HDF5_ERROR (H5Dclose (dataset));
    HDF5_ERROR (H5Sclose (filespace));
  }

  return 1;
}

} // namespace CarpetIOHDF5