aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetIOHDF5/src/iohdf5.cc
blob: 2e676081c30602e82889c86f5f638541e93b8f17 (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
#include <assert.h>
#include <limits.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <sys/stat.h>
#include <sys/types.h>

#include <algorithm>
#include <fstream>
#include <sstream>
#include <vector>

#include <hdf5.h>

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

extern "C" {
  static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetIOHDF5/src/iohdf5.cc,v 1.5 2004/03/08 22:50:41 cott Exp $";
  CCTK_FILEVERSION(Carpet_CarpetIOHDF5_iohdf5_cc);
}

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

#include "bbox.hh"
#include "data.hh"
#include "gdata.hh"
#include "ggf.hh"
#include "vect.hh"

#include "carpet.hh"

#include "iohdf5.hh"



namespace CarpetIOHDF5 {
  
  using namespace std;
  using namespace Carpet;
  
  
  
  // Variable definitions
  int GHExtension;
  int IOMethod;
  vector<bool> do_truncate;     // [var]
  vector<vector<vector<int> > > last_output; // [ml][rl][var]
  
  
  int CarpetIOHDF5Startup ()
  {
    int ierr;
    
    CCTK_RegisterBanner ("AMR 3D HDF5 I/O provided by CarpetIOHDF5");
    
    GHExtension = CCTK_RegisterGHExtension ("CarpetIOHDF5");
    CCTK_RegisterGHExtensionSetupGH (GHExtension, SetupGH);
    
    IOMethod = CCTK_RegisterIOMethod ("CarpetIOHDF5");
    CCTK_RegisterIOMethodOutputGH (IOMethod, OutputGH);
    CCTK_RegisterIOMethodOutputVarAs (IOMethod, OutputVarAs);
    CCTK_RegisterIOMethodTimeToOutput (IOMethod, TimeToOutput);
    CCTK_RegisterIOMethodTriggerOutput (IOMethod, TriggerOutput);
    
#if 0
    ierr = IOUtil_RegisterRecover ("CarpetIOHDF5", Recover);
    assert (! ierr);
#endif
    
    return 0;
  }
  
  
  
  void* SetupGH (tFleshConfig* const fc,
		 const int convLevel, cGH* const cctkGH)
  {
    DECLARE_CCTK_PARAMETERS;
    
    // Truncate all files if this is not a restart
    do_truncate.resize(CCTK_NumVars(), true);
    
    // No iterations have yet been output
    last_output.resize(mglevels);
    for (int ml=0; ml<mglevels; ++ml) {
      last_output.at(ml).resize(maxreflevels);
      for (int rl=0; rl<maxreflevels; ++rl) {
        last_output.at(ml).at(rl).resize(CCTK_NumVars(), INT_MIN);
      }
    }
    
    // We register only once, ergo we get only one handle.  We store
    // that statically, so there is no need to pass anything to
    // Cactus.
    return 0;
  }
  
  
  
  int OutputGH (const cGH* const cctkGH) {
    for (int vindex=0; vindex<CCTK_NumVars(); ++vindex) {
      if (TimeToOutput(cctkGH, vindex)) {
	TriggerOutput(cctkGH, vindex);
      }
    }
    return 0;
  }
  
  
  
  int OutputVarAs (const cGH* const cctkGH, const char* const varname,
		   const char* const alias) {
    DECLARE_CCTK_ARGUMENTS;
    DECLARE_CCTK_PARAMETERS;
    
    herr_t herr;
    
    const int n = CCTK_VarIndex(varname);
    assert (n>=0 && n<CCTK_NumVars());
    const int group = CCTK_GroupIndexFromVarI (n);
    assert (group>=0 && group<(int)Carpet::arrdata.size());
    const int n0 = CCTK_FirstVarIndexI(group);
    assert (n0>=0 && n0<CCTK_NumVars());
    const int var = n - n0;
    assert (var>=0 && var<CCTK_NumVars());
    const int tl = 0;
    
    // Check for storage
    if (! CCTK_QueryGroupStorageI(cctkGH, group)) {
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
		  "Cannot output variable \"%s\" because it has no storage",
		  varname);
      return 0;
    }
    
    const int grouptype = CCTK_GroupTypeI(group);
    switch (grouptype) {
    case CCTK_SCALAR:
    case CCTK_ARRAY:
      assert (do_global_mode);
      break;
    case CCTK_GF:
      /* do nothing */
      break;
    default:
      assert (0);
    }
    const int rl = grouptype==CCTK_GF ? reflevel : 0;
    
    // Get grid hierarchy extentsion from IOUtil
    const ioGH * const iogh = (const ioGH *)CCTK_GHExtension (cctkGH, "IO");
    assert (iogh);
    
    // Create the output directory
    const char* const myoutdir = GetStringParameter("out3D_dir", out_dir);
    if (CCTK_MyProc(cctkGH)==0) {
      CCTK_CreateDirectory (0755, myoutdir);
    }
    
    // Invent a file name
    ostringstream filenamebuf;
    filenamebuf << myoutdir << "/" << alias << out3D_extension;
    string filenamestr = filenamebuf.str();
    const char * const filename = filenamestr.c_str();
    
    hid_t writer = -1;
    
    // Write the file only on the root processor
    if (CCTK_MyProc(cctkGH)==0) {
      
      // If this is the first time, then create and truncate the file
      if (do_truncate.at(n)) {
	struct stat fileinfo;
	if (! iogh->recovered || stat(filename, &fileinfo)!=0) {
	  writer = H5Fcreate (filename, H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
          assert (writer>=0);
          herr = H5Fclose (writer);
          assert (!herr);
	  writer = -1;
	}
      }
      
      // Open the file 
      writer = H5Fopen (filename, H5F_ACC_RDWR, H5P_DEFAULT);
      assert (writer>=0);
      
//       const int gpdim = CCTK_GroupDimI(group);
      
//       // Set coordinate information
//       double origin[dim], delta[dim], timestep;
//       for (int d=0; d<dim; ++d) {
// 	origin[d] = cctkGH->cctk_origin_space[d];
// 	delta[d] = cctkGH->cctk_delta_space[d];
//       }
//       timestep = cctkGH->cctk_delta_time;
//       amrwriter->setTopLevelParameters
// 	(gpdim, origin, delta, timestep, maxreflevels);
      
//       // Set refinement information
//       int interlevel_timerefinement;
//       int interlevel_spacerefinement[dim];
//       int initial_gridplacementrefinement[dim];
//       interlevel_timerefinement = reffact;
//       for (int d=0; d<dim; ++d) {
// 	interlevel_spacerefinement[d] = reffact;
// 	initial_gridplacementrefinement[d] = 1;
//       }
//       amrwriter->setRefinement
// 	(interlevel_timerefinement, interlevel_spacerefinement,
// 	 initial_gridplacementrefinement);
      
//       // Set level
//       amrwriter->setLevel (rl);
      
//       // Set current time
//       amrwriter->setTime (cctk_iteration);
    }
    
    // Traverse all components
    BEGIN_MAP_LOOP(cctkGH, grouptype) {
      BEGIN_COMPONENT_LOOP(cctkGH, grouptype) {
        
        const ggf<dim>* ff = 0;
        
        ff = (ggf<dim>*)arrdata.at(group).at(Carpet::map).data.at(var);
        
        const gdata<dim>* const data = (*ff) (tl, rl, component, mglevel);
        
        // Make temporary copy on processor 0
        const ibbox ext = data->extent();
//         vect<int,dim> lo = ext.lower();
//         vect<int,dim> hi = ext.upper();
//         vect<int,dim> str = ext.stride();
        
        gdata<dim>* const tmp = data->make_typed (n);
        tmp->allocate (ext, 0);
        for (comm_state<dim> state; !state.done(); state.step()) {
          tmp->copy_from (state, data, ext);
        }
        
        // Write data
        if (CCTK_MyProc(cctkGH)==0) {
          
          hsize_t shape[dim];
          for (int d=0; d<dim; ++d) {
            shape[dim-1-d] = (ext.shape() / ext.stride())[d];
          }
          const hid_t dataspace = H5Screate_simple (dim, shape, NULL);
          assert (dataspace>=0);
          
          // Select datatype
          assert (true
                  || (CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL8
                      && sizeof(CCTK_REAL8) == sizeof(double))
                  || (CCTK_VarTypeI(n) == CCTK_VARIABLE_REAL
                      && sizeof(CCTK_REAL) == sizeof(double)));
          // TODO: Set datatype correctly
          const hid_t datatype = H5T_NATIVE_DOUBLE;
          
          ostringstream datasetnamebuf;
          datasetnamebuf << varname
                         << " it=" << cctk_iteration
                         << " ml=" << mglevel
                         << " rl=" << rl
                         << " m=" << Carpet::map
                         << " c=" << component;
          string datasetnamestr = datasetnamebuf.str();
          const char * const datasetname = datasetnamestr.c_str();
          const hid_t dataset = H5Dcreate (writer, datasetname, datatype, dataspace, H5P_DEFAULT);
          assert (dataset>=0);
          
          const void * const data = (void*)tmp->storage();
          herr = H5Dwrite (dataset, datatype, H5S_ALL, H5S_ALL, H5P_DEFAULT, data);
          assert (!herr);
          
          // Write FlexIO attributes
          WriteAttribute (dataset, "level", rl);
          {
            CCTK_REAL origin[dim], delta[dim];
            CCTK_REAL min_ext[dim], max_ext[dim];
            for (int d=0; d<dim; ++d) {
              origin[d] = CCTK_ORIGIN_SPACE(d) + cctk_lbnd[d] * delta[d];
              delta[d] = CCTK_DELTA_SPACE(d);
              min_ext[d] = origin[d];
              max_ext[d] = origin[d] + cctk_lsh[d] * delta[d];
            }
            WriteAttribute (dataset, "origin", origin, dim);
            WriteAttribute (dataset, "delta", delta, dim);
            WriteAttribute (dataset, "min_ext", min_ext, dim);
            WriteAttribute (dataset, "max_ext", max_ext, dim);
          }
          WriteAttribute (dataset, "time", cctk_time);
          WriteAttribute (dataset, "timestep", cctk_iteration);
          WriteAttribute (dataset, "level_timestep", cctk_iteration / reflevelfact);
          WriteAttribute (dataset, "persistence", maxreflevelfact / reflevelfact);
          {
            int time_refinement;
            int spatial_refinement[dim];
            int grid_placement_refinement[dim];
            time_refinement = reflevelfact;
            for (int d=0; d<dim; ++d) {
              spatial_refinement[d] = reflevelfact;
              grid_placement_refinement[d] = reflevelfact;
            }
            WriteAttribute (dataset, "time_refinement", time_refinement);
            WriteAttribute (dataset, "spatial_refinement", spatial_refinement, dim);
            WriteAttribute (dataset, "grid_placement_refinement", grid_placement_refinement, dim);
          }
          {
            int iorigin[dim];
            for (int d=0; d<dim; ++d) {
              iorigin[d] = (ext.lower() / ext.stride())[d];
            }
            WriteAttribute (dataset, "iorigin", iorigin, dim);
          }
          
          // Write some additional attributes
          
          // Legacy arguments
          {
            char * fullname = CCTK_FullName(n);
            assert (fullname);
            WriteAttribute (dataset, "name", fullname);
            free (fullname);
          }
          
          // Group arguments
          WriteAttribute (dataset, "group_version", 1);
          {
            char * fullname = CCTK_FullName(n);
            assert (fullname);
            WriteAttribute (dataset, "group_fullname", fullname);
            free (fullname);
          }
          WriteAttribute (dataset, "group_varname", CCTK_VarName(n));
          {
            char * groupname = CCTK_GroupName(group);
            assert (groupname);
            WriteAttribute (dataset, "group_groupname", groupname);
            free (groupname);
          }
          switch (grouptype) {
          case CCTK_GF:
            WriteAttribute (dataset, "group_grouptype", "CCTK_GF");
            break;
          case CCTK_ARRAY:
            WriteAttribute (dataset, "group_grouptype", "CCTK_ARRAY");
            break;
          case CCTK_SCALAR:
            WriteAttribute (dataset, "group_grouptype", "CCTK_SCALAR");
            break;
          default:
            assert (0);
          }
          WriteAttribute (dataset, "group_dim", CCTK_GroupDimI(group));
          WriteAttribute (dataset, "group_timelevel", tl);
          WriteAttribute (dataset, "group_numtimelevels", CCTK_NumTimeLevelsI(group));
          
          // Cactus arguments
          WriteAttribute (dataset, "cctk_version", 1);
          WriteAttribute (dataset, "cctk_dim", cctk_dim);
          WriteAttribute (dataset, "cctk_iteration", cctk_iteration);
// TODO: disable temporarily
//           WriteAttribute (dataset, "cctk_nmaps", cctk_nmaps);
//           WriteAttribute (dataset, "cctk_map", cctk_map);
          WriteAttribute (dataset, "cctk_gsh", cctk_gsh, dim);
          WriteAttribute (dataset, "cctk_lsh", cctk_lsh, dim);
          WriteAttribute (dataset, "cctk_lbnd", cctk_lbnd, dim);
          WriteAttribute (dataset, "cctk_delta_time", cctk_delta_time);
          WriteAttribute (dataset, "cctk_delta_space", cctk_delta_space, dim);
          WriteAttribute (dataset, "cctk_origin_space", cctk_origin_space, dim);
          WriteAttribute (dataset, "cctk_bbox", cctk_bbox, 2*dim);
          WriteAttribute (dataset, "cctk_levfac", cctk_levfac, dim);
          WriteAttribute (dataset, "cctk_levoff", cctk_levoff, dim);
          WriteAttribute (dataset, "cctk_levoffdenom", cctk_levoffdenom, dim);
          WriteAttribute (dataset, "cctk_timefac", cctk_timefac);
          WriteAttribute (dataset, "cctk_convlevel", cctk_convlevel);
          WriteAttribute (dataset, "cctk_convfac", cctk_convfac);
          WriteAttribute (dataset, "cctk_nghostzones", cctk_nghostzones, dim);
          WriteAttribute (dataset, "cctk_time", cctk_time);
          
          // Carpet arguments
          WriteAttribute (dataset, "carpet_version", 1);
          WriteAttribute (dataset, "carpet_dim", dim);
          WriteAttribute (dataset, "carpet_basemglevel", basemglevel);
          WriteAttribute (dataset, "carpet_mglevel", mglevel);
          WriteAttribute (dataset, "carpet_mglevels", mglevels);
          WriteAttribute (dataset, "carpet_mgface", mgfact);
          WriteAttribute (dataset, "carpet_reflevel", reflevel);
          WriteAttribute (dataset, "carpet_reflevels", reflevels);
          WriteAttribute (dataset, "carpet_reffact", reffact);
          WriteAttribute (dataset, "carpet_map", Carpet::map);
          WriteAttribute (dataset, "carpet_maps", maps);
          WriteAttribute (dataset, "carpet_component", component);
          WriteAttribute (dataset, "carpet_components", vhh.at(Carpet::map)->components(reflevel));
          
          herr = H5Dclose (dataset);
          assert (!herr);
          
          herr = H5Sclose (dataspace);
          assert (!herr);
          
        } // if on root processor
        
        // Delete temporary copy
        delete tmp;
        
      } END_COMPONENT_LOOP;
    } END_MAP_LOOP;
    
    // Close the file
    if (CCTK_MyProc(cctkGH)==0) {
      herr = H5Fclose (writer);
      assert (!herr);
      writer = -1;
    }
    
    // Don't truncate again
    do_truncate.at(n) = false;
    
    return 0;
  }
  
  
  
  int TimeToOutput (const cGH* const cctkGH, const int vindex) {
    DECLARE_CCTK_ARGUMENTS;
    DECLARE_CCTK_PARAMETERS;
    
    assert (vindex>=0 && vindex<CCTK_NumVars());

    const int grouptype = CCTK_GroupTypeFromVarI(vindex);
    switch (grouptype) {
    case CCTK_SCALAR:
    case CCTK_ARRAY:
      if (! do_global_mode) return 0;
      break;
    case CCTK_GF:
      /* do nothing */
      break;
    default:
      assert (0);
    }
    
    const int myoutevery = GetIntParameter("out3D_every", out_every);
    
    if (myoutevery < 0) {
      // Nothing should be output at all
      return 0;
    }
    
    if (cctk_iteration % myoutevery != 0) {
      // Nothing should be output during this iteration
      return 0;
    }
    
    if (! CheckForVariable(cctkGH, GetStringParameter("out3D_vars",""), vindex)) {
      // This variable should not be output
      return 0;
    }
    
    if (last_output.at(mglevel).at(reflevel).at(vindex) == cctk_iteration) {
      // Has already been output during this iteration
      char* varname = CCTK_FullName(vindex);
      CCTK_VWarn (5, __LINE__, __FILE__, CCTK_THORNSTRING,
		  "Skipping output for variable \"%s\", because this variable "
		  "has already been output during the current iteration -- "
		  "probably via a trigger during the analysis stage",
		  varname);
      free (varname);
      return 0;
    }
    
    assert (last_output.at(mglevel).at(reflevel).at(vindex) < cctk_iteration);
    
    // Should be output during this iteration
    return 1;
  }
  
  
  
  int TriggerOutput (const cGH* const cctkGH, const int vindex) {
    assert (vindex>=0 && vindex<CCTK_NumVars());
    
    char* varname = CCTK_FullName(vindex);
    const int retval = OutputVarAs (cctkGH, varname, CCTK_VarName(vindex));
    free (varname);
    
    last_output.at(mglevel).at(reflevel).at(vindex) = cctkGH->cctk_iteration;
    
    return retval;
  }
  
  
  

  int InputGH (const cGH* const cctkGH) {
    int retval = 0;
    for (int vindex=0; vindex<CCTK_NumVars(); ++vindex) {
      if (CheckForVariable(cctkGH, GetStringParameter("in3D_vars",""), vindex)) {
	char* varname = CCTK_FullName(vindex);
	retval = InputVarAs (cctkGH, varname, CCTK_VarName(vindex));
	free (varname);
	if (retval != 0) return retval;
      }
    }
    return retval;
  }
  
  
  
  int InputVarAs (const cGH* const cctkGH, const char* const varname, 
		  const char* const alias) {
    DECLARE_CCTK_PARAMETERS;
    
    const int n = CCTK_VarIndex(varname);
    assert (n>=0 && n<CCTK_NumVars());
    const int group = CCTK_GroupIndexFromVarI (n);
    assert (group>=0 && group<(int)Carpet::arrdata.size());
    const int n0 = CCTK_FirstVarIndexI(group);
    assert (n0>=0 && n0<CCTK_NumVars());
    const int var = n - n0;
    assert (var>=0 && var<CCTK_NumVars());
    const int tl = 0;

    herr_t herr;
    int have_dataset = 0;
    int want_dataset = 0;
    int did_read_something = 0;
    int ndatasets = 0;
    hid_t dataset = 0;

    char datasetname[1024];

    CCTK_REAL *h5data = NULL;

    // Check for storage
    if (! CCTK_QueryGroupStorageI(cctkGH, group)) {
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Cannot input variable \"%s\" because it has no storage",
                  varname);
      return 0;
    }
    
    const int grouptype = CCTK_GroupTypeI(group);
    const int rl = grouptype==CCTK_GF ? reflevel : 0;
    
    // Find the input directory
    const char* myindir = GetStringParameter("in3D_dir", ".");
    
    // Invent a file name
    ostringstream filenamebuf;
    filenamebuf << myindir << "/" << alias << in3D_extension;
    string filenamestr = filenamebuf.str();
    const char * const filename = filenamestr.c_str();
    
    hid_t reader = -1;
    
    const int gpdim = CCTK_GroupDimI(group);
    
    // Read the file only on the root processor
    if (CCTK_MyProc(cctkGH)==0) {
      
      // Open the file 
      if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Opening file \"%s\"", filename);
      reader = H5Fopen (filename, H5F_ACC_RDONLY, H5P_DEFAULT);
      if (reader<0) {
        CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "Could not open file \"%s\" for reading", filename);
      }
      assert (reader>=0);
      // get the number of datasets in the file 
      ndatasets=GetnDatasets(reader);
    }
    
    vector<ibset> regions_read(Carpet::maps);
    
    // Broadcast number of datasets
    MPI_Bcast (&ndatasets, 1, MPI_INT, 0, dist::comm);
    assert (ndatasets>=0);
      
    
    for (int datasetid=0; datasetid<ndatasets; ++datasetid) {
      if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Handling dataset #%d", datasetid);

        
      // Read data
      if (CCTK_MyProc(cctkGH)==0) {
	GetDatasetName(reader,datasetid,datasetname);
	cout << datasetname << "\n";
  
	dataset = H5Dopen (reader, datasetname);
	assert(dataset);
      }
    
      
      int amr_level;
      int amr_origin[dim];
      int amr_dims[dim];
	
      if (CCTK_MyProc(cctkGH)==0) {
	  
       // Read data
       char * name;
       //       cout << "reading name" << "\n";
       ReadAttribute (dataset, "name", name);
       // cout << "done reading name" << "\n";
       if (verbose) {
	 if (name) {
		CCTK_VInfo (CCTK_THORNSTRING, "Dataset name is \"%s\"", name);
	 }
       }
       want_dataset = name && CCTK_EQUALS(name, varname);
       free (name);
     
       if(want_dataset) {
	 
	 // get dataset dimensions
	 const hid_t dataspace = H5Dget_space(dataset);
	 assert (dataspace>=0);
	 hsize_t rank=H5Sget_simple_extent_ndims(dataspace);
	 hsize_t shape[rank];
	 int rank2 = H5Sget_simple_extent_dims (dataspace, shape, NULL);
	 herr = H5Sclose(dataspace);
	 assert(!herr);
	 assert (rank2 == rank);
	 assert (gpdim == rank);
	    
	 int datalength=1;
	 
	 for(int i=0;i<rank;i++) {
	   datalength=datalength*shape[i];
	 }
	 const hid_t datatype = H5T_NATIVE_DOUBLE;
	 
	 //	 cout << "datalength: " << datalength << " rank: " << rank << "\n";
	 // cout << shape[0] << " " << shape[1] << " " << shape[2] << "\n";
	 
	 h5data = (CCTK_REAL*) malloc(sizeof(double)*datalength);
	 herr = H5Dread(dataset,datatype,H5S_ALL, H5S_ALL, H5P_DEFAULT,(void*)h5data);
	 assert(!herr);
	 
	 ReadAttribute(dataset,"level",amr_level);
	 //cout << amr_level << "," << gpdim << "\n";
	 ReadAttribute(dataset,"iorigin",amr_origin,dim);
	 //    cout << amr_origin[0] << "\n";
	 
	 herr = H5Dclose(dataset);
	 assert(!herr);

	 for (int d=0; d<gpdim; ++d) {
	   amr_dims[d] = shape[d];
	 }
  
       } // want_dataset
     } // MyProc == 0
	
     MPI_Bcast (&want_dataset, 1, MPI_INT, 0, dist::comm);
     MPI_Bcast (&amr_level, 1, MPI_INT, 0, dist::comm);
     MPI_Bcast (amr_origin, dim, MPI_INT, 0, dist::comm);
     MPI_Bcast (amr_dims, dim, MPI_INT, 0, dist::comm);
	
     if (want_dataset && amr_level == reflevel) {
       did_read_something = true;
	  
       // Traverse all components on all levels
       BEGIN_MAP_LOOP(cctkGH, grouptype) {
	 BEGIN_COMPONENT_LOOP(cctkGH, grouptype) {
            
	   ggf<dim>* ff = 0;
	      
	   assert (var < (int)arrdata.at(group).at(Carpet::map).data.size());
	   ff = (ggf<dim>*)arrdata.at(group).at(Carpet::map).data.at(var);
            
	   gdata<dim>* const data = (*ff) (tl, rl, component, mglevel);
	   
	   // Create temporary data storage on processor 0
	   const vect<int,dim> str
	     = vect<int,dim>(maxreflevelfact/reflevelfact);
	   const vect<int,dim> lb = vect<int,dim>::ref(amr_origin) * str;
	   const vect<int,dim> ub
	     = lb + (vect<int,dim>::ref(amr_dims) - 1) * str;
	   const bbox<int,dim> ext(lb,ub,str);
	   
	   gdata<dim>* const tmp = data->make_typed (n);
	   
	   if (CCTK_MyProc(cctkGH)==0) {
	     tmp->allocate (ext, 0, h5data);
	   } else {
	     tmp->allocate (ext, 0);
	   }
            
	   // Initialise with what is found in the file -- this does
	   // not guarantee that everything is initialised.
	   const bbox<int,dim> overlap = tmp->extent() & data->extent();
	   regions_read.at(Carpet::map) |= overlap;
	   
	   // Copy into grid function
	   for (comm_state<dim> state; !state.done(); state.step()) {
	     data->copy_from (state, tmp, overlap);
	   }
	   
	   // Delete temporary copy
	   delete tmp;
	      
	 } END_COMPONENT_LOOP;
       } END_MAP_LOOP;
       
     } // if want_dataset && level == reflevel
       
     if (CCTK_MyProc(cctkGH)==0) {
       free (h5data);
     }
  
    } // loop over datasets
      
    // Close the file
    if (CCTK_MyProc(cctkGH)==0) {
      if (verbose) CCTK_VInfo (CCTK_THORNSTRING, "Closing file");
      herr = H5Fclose(reader);
      //	  cout << "blah! " << reader << "\n";
      // cout << "closing file " << herr << "\n";
      assert(!herr);
      reader=-1;
    }
  
    // Was everything initialised?
    if (did_read_something) {
      for (int m=0; m<Carpet::maps; ++m) {
	dh<dim>& thedd = *arrdata.at(group).at(m).dd;
	ibset all_exterior;
	for (size_t c=0; c<thedd.boxes.at(rl).size(); ++c) {
	  all_exterior |= thedd.boxes.at(rl).at(c).at(mglevel).exterior;
	    }
	if (regions_read.at(m) != all_exterior) {
	  CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
		      "Variable \"%s\" could not be initialised from file -- the file may be missing data",
		      varname);
	}
      }
    } // if did_read_something
    
    return did_read_something ? 0 : -1;
//	CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,"stop!");
  }

  
  /** returns the number of recovered variables */
  int Recover (cGH* const cctkGH, const char *basefilename,
               const int called_from)
  {
    assert (cctkGH);
    assert (basefilename);
    assert (called_from == CP_INITIAL_DATA
            || called_from == CP_EVOLUTION_DATA
            || called_from == CP_RECOVER_PARAMETERS
            || called_from == CP_RECOVER_DATA
            || called_from == FILEREADER_DATA);
    
    // the other modes are not supported yet
    assert (called_from == FILEREADER_DATA);
    
    const ioGH * const iogh = (const ioGH *) CCTK_GHExtension (cctkGH, "IO");
    assert (iogh);
    
    int num_vars_read = 0;
    assert (iogh->do_inVars);
    for (int n=0; n<CCTK_NumVars(); ++n) {
      if (iogh->do_inVars[n]) {
        char * const fullname = CCTK_FullName(n);
        assert (fullname);
        const int ierr = InputVarAs (cctkGH, fullname, basefilename);
        if (! ierr) {
          ++ num_vars_read;
        }
        free (fullname);
      }
    }
    
    return num_vars_read;
  }

  
  
  int CarpetIOHDF5ReadData (CCTK_ARGUMENTS)
  {
    DECLARE_CCTK_ARGUMENTS;
    return InputGH(cctkGH);
  }

  
  
} // namespace CarpetIOHDF5