aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetIOHDF5/src/Checkpoint.cc
blob: 9a8d1ab9e21173365c047f71136343fe06455ea0 (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
#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"

#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 "CarpetIOHDF5.hh"


namespace CarpetIOHDF5
{

using namespace std;
using namespace Carpet;

// when was the last checkpoint written ?
static int last_checkpoint_iteration = -1;


static int Checkpoint (const cGH* const cctkGH, int called_from);
static int DumpParametersGHExtentions (const cGH *cctkGH, int all, hid_t writer);


int CarpetIOHDF5_InitialDataCheckpoint (const cGH* const cctkGH)
{
  DECLARE_CCTK_PARAMETERS


  if (! CCTK_Equals (verbose, "none"))
  {
    CCTK_INFO ("---------------------------------------------------------");
    CCTK_INFO ("Dumping initial data checkpoint");
    CCTK_INFO ("---------------------------------------------------------");
  }
  int retval = Checkpoint (cctkGH, CP_INITIAL_DATA);

  return (retval);
}


int CarpetIOHDF5_EvolutionCheckpoint (const cGH* const cctkGH)
{
  int retval = 0;
  DECLARE_CCTK_PARAMETERS


  if (checkpoint &&
    ((checkpoint_every > 0 && cctkGH->cctk_iteration % checkpoint_every == 0) ||
     checkpoint_next))
  {
    if (! CCTK_Equals (verbose, "none"))
    {
      CCTK_INFO ("---------------------------------------------------------");
      CCTK_VInfo (CCTK_THORNSTRING, "Dumping periodic checkpoint at "
                  "iteration %d", cctkGH->cctk_iteration);
      CCTK_INFO ("---------------------------------------------------------");
    }

    retval = Checkpoint (cctkGH, CP_EVOLUTION_DATA);

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

  return (retval);
}


int CarpetIOHDF5_TerminationCheckpoint (const cGH *const GH)
{
  int retval = 0;
  DECLARE_CCTK_PARAMETERS


  if (checkpoint && checkpoint_on_terminate)
  {
    if (last_checkpoint_iteration < GH->cctk_iteration)
    {
      if (! CCTK_Equals (verbose, "none"))
      {
        CCTK_INFO ("---------------------------------------------------------");
        CCTK_VInfo (CCTK_THORNSTRING, "Dumping termination checkpoint at "
                    "iteration %d", GH->cctk_iteration);
        CCTK_INFO ("---------------------------------------------------------");
      }

      retval = Checkpoint (GH, CP_EVOLUTION_DATA);
    }
    else if (! CCTK_Equals (verbose, "none"))
    {
      CCTK_INFO ("---------------------------------------------------------");
      CCTK_VInfo (CCTK_THORNSTRING, "Termination checkpoint already dumped "
                  "as last evolution checkpoint at iteration %d",
                  last_checkpoint_iteration);
      CCTK_INFO ("---------------------------------------------------------");
    }
  }

  return (retval);
}


static int Checkpoint (const cGH* const cctkGH, int called_from)
{
  int retval = 0;
  DECLARE_CCTK_PARAMETERS


  CarpetIOHDF5GH *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);


  /* get the filenames for both the temporary and real checkpoint file */
  char *filename = IOUtil_AssembleFilename (cctkGH, NULL, "", ".h5",
                                            called_from, 0, 1);
  char *tempname = IOUtil_AssembleFilename (cctkGH, NULL, ".tmp", ".h5",
                                            called_from, 0, 1);

  hid_t writer = -1;

  if (myproc == 0)
  {
    if (CCTK_Equals (verbose, "full"))
    {
      CCTK_VInfo (CCTK_THORNSTRING, "Creating temporary checkpoint file '%s'",
                  tempname);
    }

    HDF5_ERROR (writer = H5Fcreate (tempname, H5F_ACC_TRUNC, H5P_DEFAULT,
                                    H5P_DEFAULT));

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

  }

  // now dump the grid variables on all mglevels, reflevels, maps and components
  BEGIN_MGLEVEL_LOOP (cctkGH)
  {
    BEGIN_REFLEVEL_LOOP (cctkGH)
    {
      if (CCTK_Equals (verbose, "full"))
      {
        CCTK_VInfo (CCTK_THORNSTRING,
                    "Dumping grid variables on mglevel %d reflevel %d ...",
                    mglevel, reflevel);
      }

      for (int group = CCTK_NumGroups () - 1; group >= 0; group--)
      {
        /* only dump groups which have storage assigned */
        if (CCTK_QueryGroupStorageI (cctkGH, group) <= 0)
        {
          continue;
        }

        cGroup gdata;
        CCTK_GroupData (group, &gdata);
        assert (gdata.grouptype == CCTK_ARRAY || gdata.grouptype == CCTK_GF ||
                gdata.grouptype == CCTK_SCALAR);

        // scalars and grid arrays only have one reflevel
        if (gdata.grouptype != CCTK_GF && reflevel > 0)
        {
          continue;
        }

        /* get the number of allocated timelevels */
        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 */
        ioRequest *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++)
        {
          char *fullname = CCTK_FullName (request->vindex);
          assert (fullname);

          /* loop over all timelevels of this variable */
          for (request->timelevel = 0;
               request->timelevel < gdata.numtimelevels;
               request->timelevel++)
          {
            if (CCTK_Equals (verbose, "full"))
            {
              CCTK_VInfo (CCTK_THORNSTRING, "  %s (timelevel %d)",
                          fullname, request->timelevel);
            }

            // write the var
            retval += WriteVar (cctkGH, writer, request, 1);
          }
          free (fullname);

        } /* end of loop over all variables */

        // free I/O request structure
        IOUtil_FreeIORequest (&request);

      } /* end of loop over all groups */
    } END_REFLEVEL_LOOP;

  } END_MGLEVEL_LOOP;


  // Close the file
  if (writer >= 0)
  {
    HDF5_ERROR (H5Fclose(writer));
  }

  if (retval == 0 && myproc == 0)
  {
    if (rename (tempname, filename))
    {
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Could not rename temporary checkpoint file '%s' to '%s'",
                  tempname, 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 (filename);
      myGH->cp_filename_index = (myGH->cp_filename_index+1) % abs (checkpoint_keep);
    }
  }

  // save the iteration number of this checkpoint
  last_checkpoint_iteration = cctkGH->cctk_iteration;

  // free allocated resources
  free (tempname);
  free (filename);

  return retval;

} // Checkpoint


static 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, dspace, dset;
  hsize_t size;

  int itmp;
  double dtmp;
  const char *version;
  const ioGH *ioUtilGH;

  if (CCTK_Equals (verbose, "full"))
  {
    CCTK_INFO ("Dumping Parameters and GH Extentions...");
  }

  /* get the parameter string buffer */
  parameters = IOUtil_GetAllParameters (cctkGH, all);
  if (parameters)
  {
    size = strlen (parameters) + 1;
    HDF5_ERROR (group = H5Gcreate (writer, METADATA_GROUP, 0));
    HDF5_ERROR (dspace = H5Screate_simple (1, &size, NULL));
    HDF5_ERROR (dset = H5Dcreate (group, ALL_PARAMETERS, H5T_NATIVE_UCHAR,
                                  dspace, H5P_DEFAULT));
    HDF5_ERROR (H5Dwrite (dset, H5T_NATIVE_CHAR, H5S_ALL, H5S_ALL,
                          H5P_DEFAULT, parameters));
    free (parameters);

    // now dump the GH Extentions

    /* get the handle for IOUtil extensions */
    ioUtilGH = (const ioGH *) CCTK_GHExtension (cctkGH, "IO");

    itmp = CCTK_MainLoopIndex ();
    WriteAttribute(dset,"main loop index",itmp);

    itmp = cctkGH->cctk_iteration;
    WriteAttribute(dset,"GH$iteration",itmp);

    itmp = ioUtilGH->ioproc_every;
    WriteAttribute(dset,"GH$ioproc_every",itmp);

    itmp = CCTK_nProcs (cctkGH);
    WriteAttribute(dset,"GH$nprocs",itmp);

    dtmp = cctkGH->cctk_time;
    WriteAttribute(dset,"GH$time", dtmp);

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

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

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

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

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

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

    HDF5_ERROR (H5Dclose (dset));
    HDF5_ERROR (H5Sclose (dspace));
    HDF5_ERROR (H5Gclose (group));
  }

  return 0;
} // DumpParametersGHExtentions

} // namespace CarpetIOHDF5