aboutsummaryrefslogtreecommitdiff
path: root/src/RecoverGH.c
blob: 5a674a0e00515d65b0547c7087d3de154390f09a (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
 /*@@
   @file      RecoverGH.c
   @date      Tue Oct 10 2000
   @author    Thomas Radke
   @desc
              Contains the routines to recover from a streamed HDF5 checkpoint.

              Currently can recover from:
                (1) One file containing recombined data
                (2) Multiple unrecombined files, where the current
                    number of processors and IO processors
                    match those used to write the files.
   @enddesc
   @version   $Id$
 @@*/


#include <stdlib.h>

#include "cctk.h"
#include "cctk_Parameters.h"
#include "CactusPUGH/PUGH/src/include/pugh.h"
#include "CactusBase/IOUtil/src/ioutil_CheckpointRecovery.h"
#include "ioStreamedHDF5GH.h"

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


/********************************************************************
 ********************    External Routines   ************************
 ********************************************************************/
int IOStreamedHDF5_RecoverParameters (void);


/********************************************************************
 ********************    Internal Routines   ************************
 ********************************************************************/
static int OpenFile (cGH *GH, const char *basefilename, int called_from,
                     fileinfo_t *fileinfo);


 /*@@
   @routine    IOStreamedHDF5_Recover
   @date       Tue Oct 10 2000
   @author     Thomas Radke
   @desc
               Recovers a GH from a streamed HDF5 checkpoint file.
               This routine is registered with IOUtil
               as IOStreamedHDF5's recovery routine.
   @enddesc
   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      cGH *
   @vio        in
   @endvar
   @var        basefilename
   @vdesc      the basefilename containing the 'host:port' of the sending process
   @vtype      const char *
   @vio        in
   @endvar
   @var        called_from
   @vdesc      flag indicating where this routine was called from
               (either RECOVER or filereader)
   @vtype      int
   @vio        in
   @endvar

   @calls      OpenFile
               IOHDF5Util_RecoverParameters
               IOHDF5Util_RecoverVariables
               IOHDF5Util_RecoverGHextensions
               CCTK_TimerStartI
               CCTK_TimerStopI
               IOUtil_PrintTimings

   @returntype int
   @returndesc
               -1 if there is no valid HDF5 file,
               or the returncode of
                 @seeroutine IOHDF5Util_RecoverParameters or
                 @seeroutine IOHDF5Util_RecoverVariables or
                 @seeroutine IOHDF5Util_RecoverGHextensions
   @endreturndesc
@@*/
int IOStreamedHDF5_Recover (cGH *GH, const char *basefilename, int called_from)
{
  int result;
  ioStreamedHDF5GH *myGH;
  static fileinfo_t fileinfo;  /* this is static because info is passed from
                                  CP_RECOVERY_PARAMETERS to CP_RECOVERY_DATA */
  DECLARE_CCTK_PARAMETERS


  /* start the recovery timer if we were called at CCTK_RECOVER */
  myGH = CCTK_GHExtension (GH, "IOStreamedHDF5");
  if (myGH && myGH->print_timing_info)
  {
    CCTK_TimerStartI (myGH->timers[RECOVERY_TIMER]);
  }

  /* open the file if it wasn't already opened at CCTK_RECOVER_PARAMETERS */
  /* FIXME Gab ... asymmetric levfac */
  if (called_from == CP_RECOVER_PARAMETERS ||
      called_from == FILEREADER_DATA ||
      (GH && (GH->cctk_levfac[0] > 1 || GH->cctk_convlevel > 0)))
  {
    if (OpenFile (GH, basefilename, called_from, &fileinfo) < 0)
    {
      return (-1);
    }
  }
  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 (! fileinfo.is_HDF5_file)
    {
      return (-1);
    }
  }

  /* if called at CCTK_RECOVER_PARAMETERS
     just do this and return (keeping the file open) */
  if (called_from == CP_RECOVER_PARAMETERS)
  {
    return (IOHDF5Util_RecoverParameters (&fileinfo));
  }

  /* Recover variables */
  if (CCTK_Equals (verbose, "full"))
  {
    CCTK_VInfo (CCTK_THORNSTRING, "Recovering %schunked data with ioproc %d, "
                "ioproc_every %d", fileinfo.unchunked ? "un" : "",
                fileinfo.ioproc, fileinfo.ioproc_every);
  }
  result = IOHDF5Util_RecoverVariables (GH, &fileinfo);

  /* Recover GH extensions */
  if (result == 0 && called_from == CP_RECOVER_DATA)
  {
    if (CCTK_Equals (verbose, "full"))
    {
      CCTK_INFO ("Recovering GH extensions");
    }
    result = IOHDF5Util_RecoverGHextensions (GH, &fileinfo);
  }

  /* Close the file */
  if (CCTK_MyProc (GH) == fileinfo.ioproc)
  {
    if (CCTK_Equals (verbose, "full"))
    {
      if (called_from == CP_RECOVER_DATA)
      {
        CCTK_VInfo (CCTK_THORNSTRING, "Closing checkpoint file '%s' after "
                                      "successful recovery", fileinfo.filename);
      }
      else
      {
        CCTK_VInfo (CCTK_THORNSTRING, "Closing data file '%s'",
                    fileinfo.filename);
      }
    }
    HDF5_ERROR (H5Fclose (fileinfo.file));
  }

  /* free the allocated filename */
  if (fileinfo.filename)
  {
    free (fileinfo.filename);
  }

  /* stop recovery timer and print timing info */
  if (called_from == CP_RECOVER_DATA && myGH->print_timing_info)
  {
    const char *timer_description = "Time to recover:";


    CCTK_TimerStopI (myGH->timers[RECOVERY_TIMER]);
    IOUtil_PrintTimings ("Timing information for recovery in IOStreamedHDF5:",
                         1, &myGH->timers[RECOVERY_TIMER], &timer_description);
  }

  /* print an info message */
  if (called_from == CP_RECOVER_DATA)
  {
    CCTK_VInfo (CCTK_THORNSTRING,
                "Restarting simulation at iteration %d (physical time %g)",
                GH->cctk_iteration, (double) GH->cctk_time);
  }

  return (result);
}


 /*@@
   @routine    IOStreamedHDF5_RecoverParameters
   @date       Tue Oct 10 2000
   @author     Thomas Radke
   @desc
               This routine is scheduled at CCTK_RECOVER_PARAMETERS.
               It recovers the parameters from an HDF5 checkpoint file.

               Note that it cannot be registered with IOUtil to be scheduled
               from there (as done with the IOStreamedHDF5_Recover routine)
               because the registration mechanism isn't activated yet
               at CCTK_RECOVER_PARAMETERS.
               Instead we call the generic parameter recovery routine
               from IOUtil here, and just pass the necessary callback function
               and its arguments.

               Note also that this routine doesn't get passed any parameters,
               not even a GH, because this doesn't exist yet at the time it is
               being called.
   @enddesc

   @calls      IOUtil_RecoverParameters

   @returntype int
   @returndesc
               -1 if an invalid reocvery mode was set
               or the returncode of
                 @seeroutine IOUtil_RecoverParameters
   @endreturndesc
@@*/
int IOStreamedHDF5_RecoverParameters (void)
{
  int retval;
  DECLARE_CCTK_PARAMETERS


  if (CCTK_Equals (recover, "auto") || CCTK_Equals (recover, "autoprobe"))
  {
    CCTK_WARN (2, "'IO::recover = \"auto|autoprobe\"' selected which doesn't "
                  "make sense for streamed HDF5 checkpoint files");
    retval = -1;
  }
  else
  {
    retval = IOUtil_RecoverParameters (IOStreamedHDF5_Recover, NULL, NULL);
  }

  return (retval);
}


/********************************************************************
 ********************    Internal Routines   ************************
 ********************************************************************/
 /*@@
   @routine    OpenFile
   @date       Tue Oct 10 2000
   @author     Thomas Radke
   @desc
               Open a streamed HDF5 file given by its filename
               (which has format 'host:port') and checks whether
               we can recover from that file.
               The file information is broadcasted by the IO processor(s)
               to all other processors so that everyone knows how to proceed.
   @enddesc

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      cGH *
   @vio        in
   @endvar
   @var        basefilename
   @vdesc      basefilename of the HDF5 file to recover from
               For streamed files this should be of format 'host:port'.
   @vtype      int
   @vio        in
   @endvar
   @var        called_from
   @vdesc      flag indicating where this routine was called from
   @vtype      int
   @vio        in
   @endvar
   @var        fileinfo
   @vdesc      pointer to structure describing the file
   @vtype      fileinfo_t *
   @vio        out
   @endvar

   @returntype int
   @returndesc
               0 for success, -1 if file could not be opened
   @endreturndesc
@@*/
static int OpenFile (cGH *GH, const char *basefilename, int called_from,
                     fileinfo_t *fileinfo)
{
  int nprocs, myproc;
  hid_t fapl, group, version_attr;
#ifdef CCTK_MPI
  MPI_Comm comm;
  CCTK_INT4 info[4];
#endif
  DECLARE_CCTK_PARAMETERS


  /* suppress compiler warnings about unused variables */
  (void) (called_from + 0);

#ifdef CCTK_MPI
  /* Get the communicator for broadcasting the info structure */
  /* NOTE: When recovering parameters thorn PUGH is not yet initialized
           so that we have to use MPI_COMM_WORLD in this case */
  comm = CCTK_GHExtensionHandle ("PUGH") < 0 ?
           MPI_COMM_WORLD : PUGH_pGH (GH)->PUGH_COMM_WORLD;
#endif

  /* identify myself */
  nprocs = CCTK_nProcs (GH);
  myproc = CCTK_MyProc (GH);

  /* use the basename (which should be of format 'host:port') as the filename */
  fileinfo->filename = strdup (basefilename);

  if (myproc == 0)
  {
    if (CCTK_Equals (verbose, "full"))
    {
      CCTK_VInfo (CCTK_THORNSTRING, "Opening streamed HDF5 file from '%s'",
                  fileinfo->filename);
    }

    /* set up file access property list and select Stream VFD */
    HDF5_ERROR (fapl = H5Pcreate (H5P_FILE_ACCESS));
    HDF5_ERROR (H5Pset_fapl_stream (fapl, NULL));

    fileinfo->file = H5Fopen (fileinfo->filename, H5F_ACC_RDONLY, fapl);
    HDF5_ERROR (H5Pclose (fapl));
    if (fileinfo->file < 0)
    {
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Cannot open streamed file HDF5 file from '%s'",
                  fileinfo->filename);
    }
    else if ((group = H5Gopen (fileinfo->file, GLOBAL_ATTRIBUTES_GROUP)) < 0)
    {
      CCTK_WARN (1, "Can't find global attributes group. "
                    "Is this really a Cactus HDF5 datafile ?");
    }
    else
    {
      /* Determine how the data was written */
      READ_ATTRIBUTE (group, "nprocs", H5T_NATIVE_INT, &fileinfo->nprocs);
      READ_ATTRIBUTE (group, "unchunked", H5T_NATIVE_INT, &fileinfo->unchunked);
      READ_ATTRIBUTE (group, "ioproc_every", H5T_NATIVE_INT,
                      &fileinfo->ioproc_every);

      /* check if there exists a version attribute
         For this we temporarily turn off automatic error printing. */
      H5E_BEGIN_TRY
      {
        version_attr = H5Aopen_name (group, "Cactus version");
      } H5E_END_TRY;

      fileinfo->has_version = version_attr >= 0;
      if (version_attr >= 0)
      {
        HDF5_ERROR (H5Aclose (version_attr));
      }

      HDF5_ERROR (H5Gclose (group));

      /* If we recover from chunked file(s) the number of
       * writing processors must match the number of reading
       * processors, and the total number of processors must match.
       */

      if ((fileinfo->ioproc_every == nprocs && nprocs > 1) ||
          fileinfo->unchunked)
      {
        if (CCTK_Equals (verbose, "full"))
        {
          CCTK_VInfo (CCTK_THORNSTRING, "Recovering from one %schunked file",
                                        fileinfo->unchunked ? "un":"");
        }
        fileinfo->ioproc_every = nprocs;
        fileinfo->is_HDF5_file = 1;
      }
      else
      {
        if (fileinfo->nprocs != nprocs)
        {
          CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                      "Must restart on %d processors with chunked files "
                      "or recombine them", fileinfo->nprocs);
        }
        else
        {
          if (CCTK_Equals (verbose, "full"))
          {
            CCTK_VInfo (CCTK_THORNSTRING, "Recovering from %d chunked files",
                        nprocs / fileinfo->ioproc_every +
                        (nprocs % fileinfo->ioproc_every ? 1 : 0));
          }
          fileinfo->is_HDF5_file = 1;
        }
      }
    }
  }

  if (myproc == 0 && ! fileinfo->is_HDF5_file)
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "No valid HDF5 file '%s' found", fileinfo->filename);
  }

#ifdef CCTK_MPI
  /* Broadcast the file information to all processors
     Need to convert everything into CCTK_INTs which can be communicated. */
  info[0] = fileinfo->is_HDF5_file;
  info[1] = fileinfo->unchunked;
  info[2] = fileinfo->ioproc_every;
  info[3] = fileinfo->has_version;
  CACTUS_MPI_ERROR (MPI_Bcast (info, 4, PUGH_MPI_INT4, 0, comm));
  fileinfo->is_HDF5_file = info[0];
  fileinfo->unchunked = info[1];
  fileinfo->ioproc_every = info[2];
  fileinfo->has_version = info[3];
#endif

  /* return 0 for success otherwise negative */
  return (fileinfo->is_HDF5_file ? 0 : -1);
}