aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetIOHDF5/src/iohdf5chckpt_recover.cc
blob: 035d05f4cf2805197443fca4b6bdb45e6cdecffc (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
#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 <string>
#include <vector>

#include <hdf5.h>

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

extern "C" {
  static const char* rcsid = "$Header:$";
  CCTK_FILEVERSION(Carpet_CarpetIOHDF5_iohdf5chckpt_recover_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"
#include "iohdf5GH.h"

/* some macros for HDF5 group names */
#define PARAMETERS_GLOBAL_ATTRIBUTES_GROUP "Parameters and Global Attributes"
#define ALL_PARAMETERS "All Parameters"



namespace CarpetIOHDF5 {
  
  using namespace std;
  using namespace Carpet;

  // linked list for reading in the checkpoint file

  list<string> datasetnamelist;
  
  int Checkpoint (const cGH* const cctkGH, int called_from);
  int DumpParametersGHExtentions (const cGH *cctkGH, int all, hid_t writer);  
  
  int RecoverParameters (hid_t reader);
  int RecoverGHextensions (cGH* cctkGH, hid_t reader);
  int RecoverVariables (cGH* cctkGH, hid_t reader);

  void CarpetIOHDF5InitialDataCheckpoint( const cGH* const cgh){
    
    DECLARE_CCTK_PARAMETERS;
    
    if (checkpoint && checkpoint_ID)
      {
	CCTK_INFO ("---------------------------------------------------------");
	CCTK_INFO ("Dumping initial data checkpoint");
	CCTK_INFO ("---------------------------------------------------------");
	
	Checkpoint (cgh, CP_INITIAL_DATA);
      }
  } // CarpetIOHDF5_InitialDataCheckpoint


  void CarpetIOHDF5EvolutionCheckpoint( const cGH* const cgh){
    
    DECLARE_CCTK_PARAMETERS;

    // Test checkpoint_every, and adjust it. This is necessary until 
    // recovery is more flexible.
    int every_full = pow(2.0,maxreflevels-1);
    if (checkpoint_every>0 && (checkpoint_every % every_full) != 0) {
      int every = (checkpoint_every / every_full + 1) * every_full;
      ostringstream outevery;
      outevery << every;
      CCTK_ParameterSet("checkpoint_every","IOUtil",
			outevery.str().c_str());
      CCTK_VInfo (CCTK_THORNSTRING,"I have adjusted your checkpoint_every to %i.",every);
      
    }

      if (checkpoint &&
      ((checkpoint_every > 0 && cgh->cctk_iteration % checkpoint_every == 0) ||
       checkpoint_next))
      {
	CCTK_INFO ("---------------------------------------------------------");
	CCTK_VInfo (CCTK_THORNSTRING, "Dumping periodic checkpoint at "
		    "iteration %d", cgh->cctk_iteration);
	CCTK_INFO ("---------------------------------------------------------");

	Checkpoint (cgh, CP_EVOLUTION_DATA);

       	if (checkpoint_next)
	{
	  CCTK_ParameterSet ("checkpoint_next", CCTK_THORNSTRING, "no");
	}
      }
  } // CarpetIOHDF5_EvolutionCheckpoint


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

  int CarpetIOHDF5_Recover (cGH* cctkGH, const char *basefilename, int called_from) {
    
    DECLARE_CCTK_PARAMETERS;

    int result=0;
    int myproc=-1;

    CarpetIOHDF5GH *myGH;

    
    static hid_t reader=0; //this thing absolutely needs to be static!!!

    myGH = NULL;

    myproc = CCTK_MyProc (cctkGH);
    
  
    if (called_from == CP_RECOVER_PARAMETERS) {
      // Okay, let's see what we can do about the parameters
    
       // Invent a file name
      ostringstream filenamebuf;

      //      if(CCTK_nProcs(cctkGH) == 1)
      filenamebuf << recover_dir << "/" << basefilename << ".h5";
      //else
      //	filenamebuf << recover_dir << "/" << basefilename << ".file_0.h5";

      string filenamestr = filenamebuf.str();
      const char * const filename = filenamestr.c_str();

      if (myproc == 0) {
	// First, open the file
	if (verbose) 
	  CCTK_VInfo(CCTK_THORNSTRING, "Opening Checkpoint file %s for recovery",filename);
	reader = H5Fopen (filename, H5F_ACC_RDONLY, H5P_DEFAULT);
	if (reader<0) {
	  CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
		    "Could not recover from \"%s\"", filename);
	}
      } // myproc == 0
    }
    else {
	/* This is the case for CP_RECOVER_DATA.
	   CCTK_RECOVER_PARAMETERS must have been called before
	   and set up the file info structure. */
	if (myproc == 0) {
	  assert(reader>=0);
	}
    }
  
    if (called_from == CP_RECOVER_PARAMETERS)
      {
	return (RecoverParameters (reader));
      }
  
    if (called_from == CP_RECOVER_DATA) {
      CCTK_INT4 numberofmgtimes=0;

      CCTK_VInfo(CCTK_THORNSTRING,"Starting to recover data on reflevel %d!!!",reflevel);

      if (myproc == 0) {
	
	// Use refinement levels parameter from checkpointing file?
	if (use_reflevels_from_checkpoint && reflevel==0) {
	  
	  herr_t herr;
	  
	  hid_t group = H5Gopen (reader, PARAMETERS_GLOBAL_ATTRIBUTES_GROUP);
	  assert(group >= 0);
	  hid_t dataset = H5Dopen (group, ALL_PARAMETERS);
	  assert(dataset>= 0);
	  
	  int reflevels_chkpt;
	
	  ReadAttribute (dataset, "carpet_reflevels", reflevels_chkpt);
      
	  herr = H5Dclose(dataset);
	  assert(!herr);
	  herr = H5Gclose(group);
	  assert(!herr);
      
	  ostringstream reflevels_str;
	  reflevels_str << reflevels_chkpt;
	  CCTK_ParameterSet ("refinement_levels","CarpetRegrid",reflevels_str.str().c_str());
	  	 
	  CCTK_VInfo (CCTK_THORNSTRING, "Using %i reflevels read from checkpoint file. Ignoring value in parameter file.", reflevels_chkpt);
	  
	}     
	
	/* we need all the times on the individual levels */
	// these are a bit messy to extract

	// Actually, we do have to do this only once

	  hid_t group = H5Gopen (reader, PARAMETERS_GLOBAL_ATTRIBUTES_GROUP);
	  assert(group>=0);
	  hid_t dataset = H5Dopen (group, ALL_PARAMETERS);
	  assert(dataset >= 0);
	  hid_t attr = H5Aopen_name (dataset, "numberofmgtimes");
	  assert(attr >= 0);
	  hid_t atype = H5Aget_type (attr);
	  if(H5Tequal(atype, H5T_NATIVE_INT)) {
	    herr_t herr = H5Aread(attr, atype, &numberofmgtimes);
	    assert(!herr);
	    herr = H5Aclose(attr);
	    assert(numberofmgtimes==mglevels);
	    char buffer[100];
	    for(int lcv=0;lcv<numberofmgtimes;lcv++) {
	      sprintf(buffer,"mgleveltimes %d",lcv);
	      attr = H5Aopen_name(dataset, buffer);
	      assert (attr>=0);
	      atype = H5Aget_type (attr);
	      assert (atype>=0);
	      herr = H5Aread (attr, atype, &leveltimes.at(lcv).at(0));
	      assert(!herr);
	      herr = H5Aclose(attr);
	      assert(!herr);
	    }  
	    herr = H5Dclose(dataset);
	    assert(!herr);
	    herr = H5Gclose(group);
	    assert(!herr);
	  } else {
	    CCTK_WARN(0,"BAD BAD BAD! Can't read leveltimes!!");
	  }

      } // myproc == 0  	


      // communicate the time on all the mglevels and reflevels

      int mpierr = MPI_Bcast (&numberofmgtimes, 1, CARPET_MPI_INT4, 0,MPI_COMM_WORLD);
      assert(!mpierr);


      for(int i=0;i<numberofmgtimes;i++) {
	mpierr = MPI_Bcast (&(leveltimes.at(i).at(0)), reflevels, CARPET_MPI_REAL, 0, MPI_COMM_WORLD);
	assert(!mpierr);
      }

      if (verbose) {
	cout << "leveltimes: " << leveltimes << endl;
      }

      cctkGH->cctk_time = leveltimes.at(mglevel).at(reflevel);

      result += RecoverGHextensions(cctkGH,reader);

	  
      if (verbose) {
	cout << "reflevel: " << reflevel << endl;
      }

      result += RecoverVariables (cctkGH,reader);


      CCTK_VInfo (CCTK_THORNSTRING,
		  "Restarting simulation at iteration %d (physical time %g)",
		  cctkGH->cctk_iteration, (double) cctkGH->cctk_time);

    } // called_from == CP_RECOVER_DATA
  
    if (reflevel==maxreflevels) {	
      if(myproc == 0) {
	H5Fclose(reader);
      }
    }

    return (result);
  } // CarpetIOHDF5_Recover

  
  int RecoverVariables (cGH* cctkGH, hid_t reader) {

    DECLARE_CCTK_PARAMETERS;

    int retval = 0;
    int myproc = CCTK_MyProc (cctkGH);
    char * name;

    int ndatasets;

    int varindex; 

    double datasettime;
    double leveltime;
    static double totaltime;

    hid_t dataset;
    herr_t herr;

    list<string> refleveldatasetnamelist;

    if (reflevel==0) {
      totaltime = 0;
    }

    leveltime = MPI_Wtime();


    if(myproc==0) {
      ndatasets=GetnDatasets(reader);
      assert (ndatasets>=0);
    }

    // Broadcast number of datasets
    MPI_Bcast (&ndatasets, 1, MPI_INT, 0, dist::comm);

    assert ((ndatasets)>=0);
    
    //if (verbose && reflevel==0) cout << "ndatasets: " << ndatasets << endl;
    if ( reflevel == 0) {
      cout << "ndatasets: " << ndatasets << endl;
    }

    if (reflevel==0) {
      for (int currdataset=0;currdataset<ndatasets+1;currdataset++){
	char datasetname[256];
	if (myproc==0) {
	  GetDatasetName(reader,currdataset,datasetname);
	  datasetnamelist.push_back(datasetname);
	} //myproc = 0
	else {
	  datasetnamelist.push_back("blah");
	}
      }
    } 

    cout << "I have " << datasetnamelist.size() << endl;

    double comparetime = MPI_Wtime();
      
    if(myproc==0) {
      //      for(currdataset=datasetnamelist.begin();
      //	  currdataset!=datasetnamelist.end();
      //	  ++currdataset) {

      list<string>::iterator currdataset;
      currdataset=datasetnamelist.begin();
       while(currdataset!=datasetnamelist.end()) {
	char tempstr[10];
	sprintf(tempstr,"rl=%d m=",reflevel);
	if ( (*currdataset).find(tempstr) < (*currdataset).size() ) {
	  refleveldatasetnamelist.push_back((*currdataset).c_str());
	  currdataset = datasetnamelist.erase(currdataset);
	} else {
	  ++currdataset;
	} // if ...
      } // while ...
    } // if(myproc==0)



   
    static long reflevelnamenum; // static keeps this thing local

    if(myproc==0) {
	reflevelnamenum=refleveldatasetnamelist.size();
      }

    MPI_Bcast (&reflevelnamenum, 1, MPI_LONG, 0, dist::comm);
    assert ((reflevelnamenum)>=0);
    
    // fill bogus namelist for non-IO cpus
    if(myproc !=0) {
      for(long i=0;i < reflevelnamenum;i++) {
	refleveldatasetnamelist.push_back("blah");
      }
    }

    comparetime = MPI_Wtime() - comparetime;
    cout << "Time for string comparison: " << comparetime << endl;
    cout << "I have for this reflevel " << refleveldatasetnamelist.size() << endl;

    list<string>::iterator currdataset;

      currdataset=refleveldatasetnamelist.begin();
 
      while(currdataset!=refleveldatasetnamelist.end()) {

	//	cout << "name: " << (*currdataset).c_str() << endl;

      if (myproc==0) {
	dataset = H5Dopen (reader, (*currdataset).c_str());
	assert(dataset);
	// Read data
	ReadAttribute (dataset, "name", name);
	varindex = CCTK_VarIndex(name);
      }

      MPI_Bcast (&varindex, 1, MPI_INT, 0, dist::comm);

      name = CCTK_FullName(varindex);

      if (verbose) {
	cout << name << "  rl: " << reflevel << endl;
      }
      vector<ibset> regions_read(Carpet::maps);

      assert (varindex>=0 && varindex<CCTK_NumVars());
      const int group = CCTK_GroupIndexFromVarI (varindex);
      const int grouptype = CCTK_GroupTypeI(group);

      int did_read_something = ReadVar(cctkGH,name,dataset,regions_read,1);

      MPI_Bcast (&did_read_something, 1, MPI_INT, 0, dist::comm);

      if (did_read_something) {
      	currdataset = refleveldatasetnamelist.erase(currdataset);
      } else {
	++currdataset;
      }

      if(myproc==0) {
	herr = H5Dclose(dataset);
	assert(!herr);
      }
      free(name);

    }  // for (currdataset ... )
    leveltime = MPI_Wtime() - leveltime;
    totaltime = totaltime + leveltime;

    cout << "Timers: leveltime: " << leveltime << " totaltime: " << totaltime << endl; 

    return retval;
  }
   
  
  int RecoverGHextensions (cGH *cctkGH, hid_t reader)
  {
    const int myproc = CCTK_MyProc(cctkGH);
    CCTK_INT4 int4Buffer[3];
    CCTK_REAL realBuffer;
    CCTK_REAL realBuffer2;
    CCTK_INT4 intbuffer;

    int mpierr = 0;

    if (myproc==0)
      {
		
	// First open group and dataset
	hid_t group = H5Gopen (reader, PARAMETERS_GLOBAL_ATTRIBUTES_GROUP);
	assert(group>=0);
	hid_t dataset = H5Dopen (group, ALL_PARAMETERS);
	assert(dataset >= 0);

	ReadAttribute(dataset,"GH$iteration",int4Buffer[0]);
	ReadAttribute(dataset,"main loop index",int4Buffer[1]);
	ReadAttribute(dataset,"carpet_global_time",realBuffer);
	//	ReadAttribute(dataset,"carpet_reflevels",int4Buffer[2]);
	ReadAttribute(dataset,"carpet_delta_time",realBuffer2);

	herr_t herr = H5Dclose(dataset);
	assert(!herr);
	herr = H5Gclose(group);
	assert(!herr);
	
      }
  /* Broadcast the GH extensions to all processors */
  /* NOTE: We have to use MPI_COMM_WORLD here
     because PUGH_COMM_WORLD is not yet set up at parameter recovery time.
     We also assume that PUGH_MPI_INT4 is a compile-time defined datatype. */

    mpierr = MPI_Bcast (int4Buffer, 3, CARPET_MPI_INT4, 0,MPI_COMM_WORLD);
    assert(!mpierr);
    mpierr = MPI_Bcast (int4Buffer, 3, CARPET_MPI_INT4, 0,MPI_COMM_WORLD);
    assert(!mpierr);
    mpierr = MPI_Bcast (&realBuffer, 1, CARPET_MPI_REAL,0,MPI_COMM_WORLD);
    assert(!mpierr);
    mpierr = MPI_Bcast (&realBuffer2, 1, CARPET_MPI_REAL,0,MPI_COMM_WORLD);
    assert(!mpierr);

    global_time = (CCTK_REAL) realBuffer;
    delta_time = (CCTK_REAL) realBuffer2;
    //    reflevels = (int) int4Buffer[2];
    cctkGH->cctk_iteration = (int) int4Buffer[0];
    CCTK_SetMainLoopIndex ((int) int4Buffer[1]);


    return (0);

  } // RecoverGHExtensions

  int RecoverParameters(hid_t reader){

    DECLARE_CCTK_PARAMETERS; 
    
    int myproc, retval;
    char *parameters;
    CCTK_INT4 parameterSize;

    hid_t group,dataset;
    herr_t herr;

    int mpierr;

    myproc = CCTK_MyProc (NULL);

    if (myproc == 0){
      CCTK_VInfo (CCTK_THORNSTRING, "Recovering parameters from checkpoint ");
      
      group = H5Gopen (reader, PARAMETERS_GLOBAL_ATTRIBUTES_GROUP);
      assert(group >= 0);
      dataset = H5Dopen (group, ALL_PARAMETERS);
      assert(dataset>= 0);

      parameterSize = H5Dget_storage_size (dataset);
      parameters = (char *) malloc (parameterSize);
      herr = H5Dread (dataset, H5T_NATIVE_CHAR, H5S_ALL, H5S_ALL,
                             H5P_DEFAULT, parameters);
      assert(!herr);
      herr = H5Dclose(dataset);
      assert(!herr);
      herr = H5Gclose(group);
      assert(!herr);

      if(verbose) {
	CCTK_VInfo (CCTK_THORNSTRING, "\n%s\n",parameters);
      }
      
      CCTK_VInfo(CCTK_THORNSTRING, "Successfully recovered parameters!");
    } // myproc == 0

    /* Broadcast the parameter buffer size to all processors */
    /* NOTE: We have to use MPI_COMM_WORLD here
       because CARPET_COMM_WORLD is not yet set up at parameter recovery time.
       We also assume that CARPET_MPI_INT4 is a compile-time defined datatype. */
    mpierr = MPI_Bcast (&parameterSize, 1, CARPET_MPI_INT4, 0,
			MPI_COMM_WORLD);
    assert(!mpierr);

    if (parameterSize > 0) {
      if (myproc) {
	parameters = (char*) malloc (parameterSize + 1);
      }
	
      mpierr = MPI_Bcast (parameters, parameterSize + 1, CARPET_MPI_CHAR,
			  0,MPI_COMM_WORLD);
      assert(!mpierr);
      
      IOUtil_SetAllParameters (parameters);
    
      free (parameters);
    }
  
    /* return positive value for success otherwise negative */
    retval = (parameterSize > 0 ? 1 : -1);

    return (retval);

  } // RecoverParameters



  int Checkpoint (const cGH* const cctkGH, int called_from)
  {
    DECLARE_CCTK_ARGUMENTS;
    DECLARE_CCTK_PARAMETERS;

    herr_t herr;
    
    int retval = 0;

    const ioGH *ioUtilGH;
    
    ioRequest *request;

    CarpetIOHDF5GH *myGH;
    myGH = (CarpetIOHDF5GH *) CCTK_GHExtension (cctkGH, "CarpetIOHDF5");

    // check if CarpetIOHDF5 was registered as I/O method
    if (myGH == NULL) {
      CCTK_WARN (0, "No CarpetIOHDF5 I/O methods registered");
      return -1;
    }
    
    int const myproc = CCTK_MyProc (cctkGH);
    
    
    
    // Invent a file name
    ostringstream cp_tempname_buf;
    ostringstream cp_filename_buf;
    
    switch (called_from)
    {
    case CP_INITIAL_DATA:
      cp_tempname_buf << checkpoint_dir << "/tmp_" << checkpoint_ID_file << ".it_" << cctkGH->cctk_iteration << ".h5";
      cp_filename_buf << checkpoint_dir << "/" << checkpoint_ID_file << ".it_" << cctkGH->cctk_iteration << ".h5";
      break;
      
    case CP_EVOLUTION_DATA:
      cp_tempname_buf << checkpoint_dir << "/tmp_" << checkpoint_file << ".it_" << cctkGH->cctk_iteration << ".h5";
      cp_filename_buf << checkpoint_dir << "/" << checkpoint_file << ".it_" << cctkGH->cctk_iteration << ".h5";
      break;
      
    case CP_RECOVER_DATA:
    case CP_RECOVER_PARAMETERS:
    case FILEREADER_DATA:
      cp_tempname_buf << recover_dir << "/tmp_" << recover_file << ".h5";
      cp_filename_buf << recover_dir << "/" << recover_file << ".h5";
      break;
      
    default:
      CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "IOUtil_PrepareFilename: Unknown calling mode %d",
                  called_from);
      break;
    }
    
    string const cp_tempname_str = cp_tempname_buf.str();
    string const cp_filename_str = cp_filename_buf.str();
    char const * const cp_tempname = cp_tempname_str.c_str();
    char const * const cp_filename = cp_filename_str.c_str();
    
    
    
    hid_t writer = -1;

    if (myproc == 0) {

      if (verbose) {
	CCTK_VInfo (CCTK_THORNSTRING, "Creating temporary checkpoint file '%s'", cp_tempname);
      }

      CCTK_VInfo (CCTK_THORNSTRING, "Creating temporary checkpoint file '%s'", cp_tempname);
      
      CCTK_VInfo (CCTK_THORNSTRING, "Verbose = %d", verbose);

      writer = H5Fcreate (cp_tempname, H5F_ACC_TRUNC, H5P_DEFAULT, H5P_DEFAULT);
      assert (writer>=0);
      herr = H5Fclose (writer);
      assert (!herr);
      writer = -1;

      // Now open the file

      writer = H5Fopen (cp_tempname, H5F_ACC_RDWR, H5P_DEFAULT);
      assert (writer>=0);

      // Dump all parameters and GHExtentions
      retval += DumpParametersGHExtentions (cctkGH, 1, writer);
      assert(!retval);

    } // myproc == 0 


    // now dump the grid variables on all mglevels, reflevels, maps and components
    BEGIN_MGLEVEL_LOOP(cctkGH) {

      BEGIN_REFLEVEL_LOOP(cctkGH) {

	if (verbose)
	  {
	    CCTK_INFO ("Dumping Grid Variables ...");
	  }
	for (int group = CCTK_NumGroups () - 1; group >= 0; group--)
	  {
	    /* only dump groups which have storage assigned */

	    if (CCTK_QueryGroupStorageI (cctkGH, group) <= 0)
	      {
		continue;
	      }

	    const int grouptype = CCTK_GroupTypeI(group);

	    /* scalars and grid arrays only have 1 reflevel: */
	    if ( (grouptype != CCTK_GF) && (reflevel != 0) ) {
	      continue;
	    }
	    /* now check if there is any memory allocated
	       for GFs and GAs. GSs should always have 
	       memory allocated and there is at this point
	       no CCTK function to check this :/
	    */

	    if ( (grouptype == CCTK_GF) || (grouptype == CCTK_ARRAY)){
	      const int gpdim = CCTK_GroupDimI(group);
	      int gtotalsize=1;
	      vector<int> tlsh(gpdim);
	      assert(!CCTK_GrouplshGI(cctkGH,gpdim,&tlsh[0],group));
	      for(int i=0;i<gpdim;i++) {
		gtotalsize=tlsh[i];		  
	      }
	      if(gtotalsize == 0){
		if (verbose) CCTK_VInfo(CCTK_THORNSTRING, 
		    "Group %s is zero-sized. No checkpoint info written",CCTK_GroupName(group));
		continue;
	      }
	    }
	    
	    /* get the number of allocated timelevels */
            cGroup gdata;
	    CCTK_GroupData (group, &gdata);
	    gdata.numtimelevels = 0;
	    gdata.numtimelevels = CCTK_GroupStorageIncrease (cctkGH, 1, &group,
							     &gdata.numtimelevels,NULL);
	    
	    int first_vindex = CCTK_FirstVarIndexI (group);

	    /* get the default I/O request for this group */
	    request = IOUtil_DefaultIORequest (cctkGH, first_vindex, 1);
	    
	    /* disable checking for old data objects, disable datatype conversion
	       and downsampling */
	    request->check_exist = 0;
	    request->hdatatype = gdata.vartype;
	    for (request->hdim = 0; request->hdim < request->vdim; request->hdim++)
	      {
		request->downsample[request->hdim] = 1;
	      }
	    
	    /* loop over all variables in this group */
	    for (request->vindex = first_vindex;
		 request->vindex < first_vindex + gdata.numvars;
		 request->vindex++)
	      {
		/* loop over all timelevels of this variable */
		for (request->timelevel = 0;
		     request->timelevel < gdata.numtimelevels;
		     request->timelevel++)
		  {
		    if (verbose)
		      {
			char *fullname = CCTK_FullName (request->vindex);
			if (fullname != NULL) {
			    CCTK_VInfo (CCTK_THORNSTRING, "  %s (timelevel %d)",
					fullname, request->timelevel);
			    free (fullname);
			}
			else {
			  CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
				     "Invalid variable with varindex %d", request->vindex);
			}
		      }
		    // write the var
		    
		    if (grouptype == CCTK_ARRAY || grouptype == CCTK_GF || grouptype == CCTK_SCALAR)
		      {
			char* fullname = CCTK_FullName (request->vindex);
			if (verbose)
			  CCTK_VInfo (CCTK_THORNSTRING,"%s: reflevel: %d map: %d component: %d grouptype: %d ",
				      fullname,reflevel,Carpet::map,component,grouptype);
			free(fullname);
		    retval += WriteVar(cctkGH,writer,request,1);
		      }
		    else
		      {
			char *fullname = CCTK_FullName (request->vindex);
			CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
				    "Invalid group type %d for variable '%s'", grouptype, fullname);
			free(fullname);
			retval = -1;
		      }
		    
		  }
	      } /* end of loop over all variables */
	  } /* end of loop over all groups */
      } END_REFLEVEL_LOOP;
      
    } END_MGLEVEL_LOOP;
  

    if (myproc==0) {
      // Close the file
      herr = H5Fclose(writer);
      assert(!herr);
    }
    
    if (retval == 0) {
      if (myproc==0) {
	if (rename (cp_tempname, cp_filename))
	  {
	    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
			"Could not rename temporary checkpoint file '%s' to '%s'",
			cp_tempname, cp_filename);
	    retval = -1;
	  }
	else {
	  if (myGH->cp_filename_list[myGH->cp_filename_index]) {
	    if (checkpoint_keep > 0) {
	      remove (myGH->cp_filename_list[myGH->cp_filename_index]);
	    }
	    free (myGH->cp_filename_list[myGH->cp_filename_index]);
	  }
	  myGH->cp_filename_list[myGH->cp_filename_index] = strdup (cp_filename);
	  myGH->cp_filename_index = (myGH->cp_filename_index+1) % abs (checkpoint_keep);
	} // else
      } // myproc == 0
    } // retval == 0

    return retval;

  } // Checkpoint

  
  int DumpParametersGHExtentions (const cGH *cctkGH, int all, hid_t writer)
  {
    // large parts of this routine were taken from
    // Thomas Radke's IOHDF5Util. Thanks Thomas!

    DECLARE_CCTK_PARAMETERS;

    char *parameters;
    hid_t group, dataspace, dataset;
    hsize_t size;
    herr_t herr;
    
    CCTK_INT4 itmp;
    CCTK_REAL dtmp;
    const char *version;
    ioGH *ioUtilGH;
  
    if (verbose) {
      CCTK_INFO ("Dumping Parameters and GH Extentions...");
    }
  
    /* get the parameter string buffer */
    parameters = IOUtil_GetAllParameters (cctkGH, all);
  
    if (parameters)
    {
      size = strlen (parameters) + 1;
      group = H5Gcreate (writer, PARAMETERS_GLOBAL_ATTRIBUTES_GROUP, 0);
      assert(group>=0);
      dataspace = H5Screate_simple (1, &size, NULL);
      assert(dataspace>=0);
      dataset = H5Dcreate (group, ALL_PARAMETERS, H5T_NATIVE_UCHAR,
                                       dataspace, H5P_DEFAULT);
      assert(dataset>=0);
      herr = H5Dwrite (dataset, H5T_NATIVE_CHAR, H5S_ALL, H5S_ALL,
                            H5P_DEFAULT, parameters);
      assert(!herr);

      // now dump the GH Extentions

      /* get the handle for IOUtil extensions */
      ioUtilGH = (ioGH *) CCTK_GHExtension (cctkGH, "IO");
      
      itmp = CCTK_MainLoopIndex ();
      WriteAttribute(dataset,"main loop index",itmp);
      
      itmp = cctkGH->cctk_iteration;
      WriteAttribute(dataset,"GH$iteration",itmp);
      
      itmp = ioUtilGH->ioproc_every;
      WriteAttribute(dataset,"GH$ioproc_every",itmp);
      
      itmp = CCTK_nProcs (cctkGH);
      WriteAttribute(dataset,"GH$nprocs",itmp);
      
      dtmp = cctkGH->cctk_time;
      WriteAttribute(dataset,"GH$time", dtmp);

      dtmp = global_time;
      WriteAttribute(dataset,"carpet_global_time", dtmp);

      itmp = reflevels;
      WriteAttribute(dataset,"carpet_reflevels", itmp);

      dtmp = delta_time;
      WriteAttribute(dataset,"carpet_delta_time", dtmp);

      version = CCTK_FullVersion();
      WriteAttribute(dataset,"Cactus version", version);

      /* finally, we need all the times on the individual refinement levels */

      const int numberofmgtimes=mglevels;
      WriteAttribute(dataset,"numberofmgtimes",numberofmgtimes);
      for(int i=0;i < numberofmgtimes;i++) {
	char buffer[100];
	sprintf(buffer,"mgleveltimes %d",i);
	WriteAttribute(dataset,buffer,(double *) &leveltimes.at(i).at(0), reflevels);
      }

      herr = H5Dclose (dataset);
      assert(!herr);
      herr = H5Sclose (dataspace);
      assert(!herr);
      herr = H5Gclose (group);
      assert(!herr);
  
      free (parameters);
    }

    return 0;
  } // DumpParametersGHExtentions
  

    
  
} // namespace CarpetIOHDF5