aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetIOHDF5/src/CarpetIOHDF5.hh
blob: 1c1d7b52bda9424a99b15194482132d06909b791 (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
#ifndef CARPETIOHDF5_HH
#define CARPETIOHDF5_HH

#define H5_USE_16_API 1
#include <hdf5.h>

#include "cctk_Arguments.h"
#include "CactusBase/IOUtil/src/ioutil_Utils.h"
#include "carpet.hh"





// some macros for HDF5 group names
#define METADATA_GROUP "Parameters and Global Attributes"
#define ALL_PARAMETERS "All Parameters"
#define GRID_STRUCTURE "Grid Structure"

// atomic HDF5 datatypes for the generic CCTK datatypes
// (the one for CCTK_COMPLEX is created at startup as a compound HDF5 datatype)
#define HDF5_CHAR   H5T_NATIVE_CHAR

#ifdef  CCTK_REAL_PRECISION_16
#define HDF5_REAL   H5T_NATIVE_LDOUBLE
#elif   CCTK_REAL_PRECISION_8
#define HDF5_REAL   H5T_NATIVE_DOUBLE
#elif   CCTK_REAL_PRECISION_4
#define HDF5_REAL   H5T_NATIVE_FLOAT
#endif

#ifdef  CCTK_INTEGER_PRECISION_8
#define HDF5_INT    H5T_NATIVE_LLONG
#elif   CCTK_INTEGER_PRECISION_4
#define HDF5_INT    H5T_NATIVE_INT
#elif   CCTK_INTEGER_PRECISION_2
#define HDF5_INT    H5T_NATIVE_SHORT
#elif   CCTK_INTEGER_PRECISION_1
#define HDF5_INT    H5T_NATIVE_CHAR
#endif


// check return code of HDF5 call and print a warning in case of an error
#define HDF5_ERROR(fn_call)                                                   \
        do {                                                                  \
          int _error_code = fn_call;                                          \
                                                                              \
                                                                              \
          if (_error_code < 0)                                                \
          {                                                                   \
            CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,              \
                        "HDF5 call '%s' returned error code %d",              \
                        #fn_call, _error_code);                               \
            error_count++;                                                    \
          }                                                                   \
        } while (0)


// CarpetIOHDF5 GH extension structure
typedef struct
{
  // default number of times to output
  int out_every_default;

  // list of variables to output
  char *out_vars;

  // stop on I/O parameter parsing errors ?
  int stop_on_parse_errors;

  // I/O request description list (for all variables)
  vector<ioRequest*> requests;

  // directory in which to output
  char *out_dir;

  // ring buffer for list of successfully created cp files
  int    checkpoint_keep;
  int    cp_filename_index;
  char **cp_filename_list;

  // list of recovery files to remove
  int    recovery_num_filenames;
  char **recovery_filename_list;

  // iteration number of the last checkpoint
  int last_checkpoint_iteration;

  // hdf5 datatype for complex variables; to be set at run time
  hid_t HDF5_COMPLEX, HDF5_COMPLEX8, HDF5_COMPLEX16, HDF5_COMPLEX32;
  
} CarpetIOHDF5GH;


namespace CarpetIOHDF5
{
  // callback routine registered for recovery/filereader
  int Recover (cGH* cctkGH, const char *basefilename, int called_from);

  // worker routines to write a single variable
  int WriteVarUnchunked (const cGH* const cctkGH,
                         hid_t file,
                         CCTK_REAL & io_bytes,
                         const ioRequest* const request,
                         bool called_from_checkpoint);
  int WriteVarChunkedSequential (const cGH* const cctkGH,
                                 hid_t file,
                                 CCTK_REAL & io_bytes,
                                 const ioRequest* const request,
                                 bool called_from_checkpoint);
  int WriteVarChunkedParallel (const cGH* const cctkGH,
                               hid_t file,
                               CCTK_REAL & io_bytes,
                               const ioRequest* const request,
                               bool called_from_checkpoint);

  // returns an HDF5 datatype corresponding to the given CCTK datatype
  hid_t CCTKtoHDF5_Datatype (const cGH* const cctkGH,
                             int cctk_type,
                             bool single_precision);

  // scheduled routines (must be declared as C according to schedule.ccl)
  extern "C" {

    int CarpetIOHDF5_Startup (void);
    void CarpetIOHDF5_Init (CCTK_ARGUMENTS);
    void CarpetIOHDF5_InitCheckpointingIntervals (CCTK_ARGUMENTS);
    int CarpetIOHDF5_SetNumRefinementLevels (void);
    void CarpetIOHDF5_CloseFiles (CCTK_ARGUMENTS);
    void CarpetIOHDF5_InitialDataCheckpoint (CCTK_ARGUMENTS);
    void CarpetIOHDF5_EvolutionCheckpoint (CCTK_ARGUMENTS);
    void CarpetIOHDF5_TerminationCheckpoint (CCTK_ARGUMENTS);
    int CarpetIOHDF5_RecoverParameters (void);
    void CarpetIOHDF5_RecoverGridStructure (CCTK_ARGUMENTS);

  } // extern "C"

} // namespace CarpetIOHDF5

#endif // !defined(CARPETIOHDF5_HH)