aboutsummaryrefslogtreecommitdiff
path: root/src/Write1D.c
blob: 764211ac8dfa5fe2c3da0dfa55236f316088485f (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
/*@@
   @file      Write1D.c
   @date      Fri May 21 1999
   @author    Thomas Radke, Gerd Lanfermann
   @desc 
   Output and trigger routines for 1D output into HDF5 files.
   @enddesc 
   @hendhistory
 @@*/


#include "cctk.h"
#include "cctk_Parameters.h"
#include "StoreNamedData.h"
#include "CactusBase/IOUtil/src/ioGH.h"
#include "ioHDF5GH.h"

#include <stdio.h>
#include <stdlib.h>
#include <assert.h>


/* local function prototypes */
void IOHDF5i_DumpParameters (cGH *GH, hid_t group);
void IOHDF5i_DumpGHExtensions (cGH *GH, hid_t group);
char *IOHDF5_Get_filename_ND (cGH *GH, int index,  ioHDF5Geo_t *slab_geo,
			      const char *name, int *isNewFile);
int NumDirection(int sdim, int vdim);
int SetDirection(int vdim, int sdim, int ci, int *direction);

/*@@
  @routine Write1D
  @author  Paul Walker 
  @date    Feb 1997
  @calledby IOHDF5_Output1DVarAs, IOHDF5_TriggerOutput1D
  @var     GH
  @vdesc   Pointer to CCTK GH
  @vtype   cGH
  @vio     in
  @vcomment
  @endvar
  @var     index
  @vdesc   index of variable to output
  @vtype   int
  @vio     in
  @vcomment
  @endvar
  @var     alias
  @vdesc   alias name of variable to output
  @vtype   const char *
  @vio     in
  @vcomment
  @endvar
@@*/
int IOHDF5_Write1D (cGH *GH, int index, const char *alias)
{
  DECLARE_CCTK_PARAMETERS
  int timelevel;
  ioGH *ioUtilGH; 
  ioHDF5GH *h5GH;
  ioHDF5Geo_t geo;

  hid_t fid, plist;
  int isNewFile;
  char *filename;

  int SDIM;
  int si,max_slabs;

  /* Get the handle for IO extensions */
  ioUtilGH = (ioGH *) GH->extensions [CCTK_GHExtensionHandle ("IO")];
  h5GH     = (ioHDF5GH *) GH->extensions [CCTK_GHExtensionHandle ("IOHDF5")];

  /* The dimension of the slab - HARDCODED */
  SDIM = 1;

  /* Get the slab geometry for 1D slabs */
  geo      = h5GH->out_geo[index][SDIM];
  /* Set geo.sdim to announce 1D output */
  geo.sdim = SDIM;
  /* get the dimension of the variable */
  geo.vdim = CCTK_GroupDimFromVarI(index);

  if (verbose)
    CCTK_VInfo (CCTK_THORNSTRING, "HDF5 %dD-output of variable '%s', "
                "downsampling parameters (%d, %d, %d)", SDIM, alias, 
                geo.downs[0],geo.downs[1], geo.downs[2]);

  /* get the current timelevel */
  timelevel = CCTK_NumTimeLevelsFromVarI (index) - 1;
  if (timelevel > 0)
    timelevel--;

  /* Maximal number of slabs (dimension sdim) in 
     given volume (dimension vdim) */
  if ((max_slabs=NumDirection(SDIM, geo.vdim))<0) 
  {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
		"Cannot get valid number of slabs (%d) for given volume (%d).\n",
	   SDIM, geo.vdim);
    return(-1);
  }
  
  /* Loop over all possible choices */
  for (si=0;si<max_slabs;si++) {

    /* Get the next set of slab directions */
    if (SetDirection(geo.vdim, SDIM, si, geo.direction)<0) {
         CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
	     "Cannot set direction for slab (#%d) in given volume (%d).\n",
	     si, geo.vdim);
	 continue;
    }

    /* get the filename for output */
    filename = IOHDF5_Get_filename_ND (GH, index, &geo,
				       alias, &isNewFile);

    /* open the output file on IO procs no parallel HDF5  */
    if (CCTK_MyProc (GH) != ioUtilGH->ioproc)
      fid = -1;
    else
    {
      
      if (verbose)
	CCTK_VInfo (CCTK_THORNSTRING, "%s HDF5 output file '%s'",
		    isNewFile ? "Opening" : "Appending", filename);
      
      CACTUS_IOHDF5_ERROR (plist = H5Pcreate (H5P_FILE_ACCESS));
      
      if (isNewFile)
	CACTUS_IOHDF5_ERROR (fid = H5Fcreate (filename, H5F_ACC_TRUNC,
					      H5P_DEFAULT, plist));
      else
	CACTUS_IOHDF5_ERROR (fid = H5Fopen (filename, H5F_ACC_RDWR, plist));
      
      H5Pclose (plist);
    }


    /* output global attributes */
    if (isNewFile && fid >= 0) {
      hid_t group;
      
      /* all GH extension variables and parameter stuff which is not
	 specific to any dataset goes into dedicated groups */
      if (out1D_parameters)
      {
	CACTUS_IOHDF5_ERROR (group = H5Gcreate (fid, GLOBAL_PARAMETERS_GROUP, 0));
	IOHDF5i_DumpParameters (GH, group);
	CACTUS_IOHDF5_ERROR (H5Gclose (group));
      }
      
      if (verbose)
	CCTK_INFO ("Dumping GH extensions ...");
      
      CACTUS_IOHDF5_ERROR (group = H5Gcreate (fid, GHEXTENSIONS_GROUP, 0));
      IOHDF5i_DumpGHExtensions (GH, group);
      CACTUS_IOHDF5_ERROR (H5Gclose (group));
    }
    
    IOHDF5_DumpVar(GH, index, timelevel, &geo, fid);
    
    /* close the file */
    if (fid >= 0)
    {
      if (verbose)
	CCTK_INFO ("Closing HDF5 output file from this iteration");
      CACTUS_IOHDF5_ERROR (H5Fclose (fid));
    }
    
    /* free the filename if it was not registered for re-opening
       in IOHDF5_Get1D_filename () */
    if (out1D_septimefiles)
      free (filename);
  }
  return(1);
}