aboutsummaryrefslogtreecommitdiff
path: root/src/ParseVars.c
blob: 2c587032fe4e1a29211ba9eca763055f54e79b5e (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
 /*@@
   @file      GHExtension.c
   @date      Tue 9th Feb 1999
   @author    Gabrielle Allen
   @desc 
              IOUtil GH extension stuff.
   @enddesc 
   @version   $Id$
 @@*/

#include <stdlib.h>
#include <string.h>

#include "cctk.h"
#include "util_String.h"
#include "cctk_Parameters.h"
#include "cctk_GNU.h"
#include "CactusBase/IOUtil/src/ioGH.h"
#include "ioHDF5UtilGH.h"


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


typedef struct
{
  const cGH *GH;
  ioSlab **slablist;
} info_t;

/* prototypes of routines defined in this source file */
static void SetOutputVar (int vindex, const char *optstring, void *arg);
static ioSlab *DefaultIOHyperslab (const cGH *GH, int vindex);

/* prototypes of external routines for which no header files exist */
int CCTK_RegexMatch (const char *string,
                     const char *pattern,
                     const int nmatch,
                     regmatch_t *pmatch);


/* =============================================================== 
   utility routines used by other IO thorns
   ===============================================================*/

 /*@@
   @routine    IOUtil_ParseVarsForOutput
   @date       Sat March 6 1999
   @author     Gabrielle Allen
   @desc 
               Sets each flag in the do_output[] do_output to true 
               if var[i] could be found in the list of variable names.
               var_list my also contain group names of variables as well as the
               special keyword "all" which indicates that output is requested
               on all variables.
   @enddesc 

   @var        out_vars
   @vdesc      list of variables and/or group names
   @vtype      const char *
   @vio        in
   @endvar 
   @var        do_output
   @vdesc      do_output of flags indicating output was requested for var[i]
   @vtype      char []
   @vio        out
   @endvar 
@@*/
void IOHDF5Util_ParseVarsForOutput (const cGH *GH, const char *out_vars,
                                    ioSlab *slablist[])
{
  int i;
  info_t info;


  /* free current list of hyperslabs */
  for (i = CCTK_NumVars () - 1; i >= 0; i--)
  {
    if (slablist[i])
    {
      free (slablist[i]->vectors);
      free (slablist[i]);
      slablist[i] = NULL;
    }
  }

  /* generate new list of hyperslabs */
  info.GH = GH;
  info.slablist = slablist;
  CCTK_TraverseString (out_vars, SetOutputVar, &info, CCTK_GROUP_OR_VAR);
}


static void SetOutputVar (int vindex, const char *optstring, void *arg)
{
  info_t *info;
  regmatch_t gmatch[6], *dmatch;
  int i, j, bytes, matched;
  char *token, *separator, *fullname;
  char *substring, *parsestring, *regexstring;
  ioSlab *slab;
  DECLARE_CCTK_PARAMETERS


  info = (info_t *) arg;

  /* allocate a new I/O hyperslab structure and initialize it with defaults */
  slab = DefaultIOHyperslab (info->GH, vindex);

  if (! optstring)
  {
    info->slablist[vindex] = slab;
    return;
  }

  /* parse the hyperslab information */
  matched = CCTK_RegexMatch (optstring,
                             "\\{([0-9]*)\\}"            /* dimension  */
                             "\\{([0-9,()]+)*\\}"        /* direction  */
                             "\\{([0-9,]+)*\\}"          /* origin     */
                             "\\{([-0-9,]+)*\\}"         /* extent     */
                             "\\{([0-9,]+)*\\}",         /* downsample */
                             6, gmatch);
  if (matched <= 0)
  {
    fullname = CCTK_FullName (vindex);
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Couldn't parse hyperslab options '%s' for variable '%s'",
                optstring, fullname);
    free (fullname);
    free (slab);
    return;
  }

  /* SLAB DIMENSION */
  bytes = (int) (gmatch[1].rm_eo - gmatch[1].rm_so);
  if (gmatch[1].rm_so != -1 && bytes > 0)
  {
    slab->hdim = atoi (optstring + gmatch[1].rm_so);
    if (slab->hdim <= 0 || slab->hdim > slab->vdim)
    {
      fullname = CCTK_FullName (vindex);
      CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Invalid dimension given %d in hyperslab options '%s' "
                  "for variable '%s'", slab->hdim, optstring, fullname);
      free (fullname);
      free (slab);
      return;
    }
  }

  /* DIRECTION */
  bytes = (int) (gmatch[2].rm_eo - gmatch[2].rm_so);
  if (gmatch[2].rm_so != -1 && bytes > 0)
  {
    substring = strdup (optstring + gmatch[2].rm_so);
    substring[bytes] = 0;

    dmatch = (regmatch_t *) malloc ((slab->hdim + 1) *
                                    sizeof (regmatch_t));
    regexstring = (char *) malloc (slab->hdim * sizeof ("\\(([0-9,]+)\\)"));
    regexstring[0] = 0;
    for (i = 0; i < slab->hdim; i++)
    {
      strcat (regexstring, "\\(([0-9,]+)\\)");
    }
    matched = CCTK_RegexMatch (substring, regexstring, slab->hdim + 1,
                               dmatch);
    free (regexstring);
    if (matched <= 0)
    {
      fullname = CCTK_FullName (vindex);
      CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Couldn't parse direction vectors in hyperslab options '%s' "
                  "for variable '%s'.", optstring, fullname);
      free (fullname);
      free (dmatch);
      free (slab);
      return;
    }

    for (j = 0; j < slab->hdim; j++)
    {
      i = 0;
      bytes = (int) (dmatch[j + 1].rm_eo - dmatch[j + 1].rm_so);
      if (dmatch[j + 1].rm_so != -1 && bytes > 0)
      {
        parsestring = strdup (substring + dmatch[j + 1].rm_so);
        parsestring[bytes] = 0;

        token = parsestring;
        while ((separator = strchr (token, ',')) != NULL)
        {
          *separator = 0;
          slab->direction[j * slab->vdim + i] = atoi (token);
          if (++i >= slab->vdim)
          {
            break;
          }
          token = separator + 1;
        }
        slab->direction[j * slab->vdim + i++] = atoi (token);
        free (parsestring);
      }
      free (substring);
      if (i < slab->vdim)
      {
        fullname = CCTK_FullName (vindex);
        CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "Direction vectors are incomplete or missing in hyperslab "
                    "options '%s' for variable '%s'.", optstring, fullname);
        free (fullname);
        free (slab);
        return;
      }
    }
  }

  /* ORIGIN */
  bytes = (int) (gmatch[3].rm_eo - gmatch[3].rm_so);
  if (gmatch[3].rm_so != -1 && bytes > 0)
  {
    substring = strdup (optstring + gmatch[3].rm_so);
    substring[bytes] = 0;

    i = 0;
    token = substring;
    while ((separator = strchr (token, ',')) != NULL)
    {
      *separator = 0;
      slab->origin[i] = atoi (token);
      if (++i >= slab->vdim)
      {
        break;
      }
      token = separator + 1;
    }
    slab->origin[i++] = atoi (token);
    free (substring);

    if (i < slab->vdim)
    {
      fullname = CCTK_FullName (vindex);
      CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Origin vector is incomplete or missing in hyperslab "
                  "options '%s' for variable '%s'.", optstring, fullname);
      free (fullname);
      free (slab);
      return;
    }
  }

  /* LENGTH */
  bytes = (int) (gmatch[4].rm_eo - gmatch[4].rm_so);
  if(gmatch[4].rm_so != -1 && bytes > 0)
  {
    substring = strdup (optstring + gmatch[4].rm_so);
    substring[bytes] = 0;

    i = 0;
    token = substring;
    while ((separator = strchr (token, ',')) != NULL)
    {
      *separator = 0;
      slab->extent[i] = atoi (token);
      if (++i >= slab->hdim)
      {
        break;
      }
      token = separator + 1;
    }
    slab->extent[i++] = atoi (token);
    free (substring);

    if (i < slab->hdim)
    {
      fullname = CCTK_FullName (vindex);
      CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Length vector is incomplete or missing in hyperslab "
                  "options '%s' for variable '%s'.", optstring, fullname);
      free (fullname);
      free (slab);
      return;
    }
  }

  /* DOWNSAMPLING */
  bytes = (int) (gmatch[5].rm_eo - gmatch[5].rm_so);
  if(gmatch[5].rm_so != -1 && bytes > 0)
  {
    substring = strdup (optstring + gmatch[5].rm_so);
    substring[bytes] = 0;

    i = 0;
    token = substring;
    while ((separator = strchr (token, ',')) != NULL)
    {
      *separator = 0;
      slab->downsample[i] = atoi (token);
      if (++i >= slab->hdim)
      {
        break;
      }
      token = separator + 1;
    }
    slab->downsample[i++] = atoi (token);
    free (substring);

    if (i < slab->hdim)
    {
      fullname = CCTK_FullName (vindex);
      CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Downsampling vector is incomplete or missing in hyperslab "
                  "options '%s' for variable '%s'.", optstring, fullname);
      free (fullname);
      free (slab);
      return;
    }
  }

  /* assign the new hyperslab */
  info->slablist[vindex] = slab;
}


static ioSlab *DefaultIOHyperslab (const cGH *GH, int vindex)
{
  ioSlab *slab;
  int *extent_int;
  const ioGH *ioUtilGH;
  DECLARE_CCTK_PARAMETERS


  ioUtilGH = (const ioGH *) CCTK_GHExtension (GH, "IO");

  /* allocate a new I/O hyperslab structure */
  slab = (ioSlab *) malloc (sizeof (ioSlab));

  /* fill out the basics */
  slab->vindex = vindex;
  slab->timelevel = 0;
  slab->check_exist = ioUtilGH->recovered;

  /* get the hyperslab datatype (will be single-precision if requested) */
  slab->hdatatype = CCTK_VarTypeI (vindex);
  if (ioUtilGH->out_single)
  {
    if (slab->hdatatype == CCTK_VARIABLE_REAL)
    {
      slab->hdatatype = CCTK_VARIABLE_REAL4;
    }
    else if (slab->hdatatype == CCTK_VARIABLE_COMPLEX)
    {
      slab->hdatatype = CCTK_VARIABLE_COMPLEX8;
    }
#ifdef CCTK_INT2
    else if (slab->hdatatype == CCTK_VARIABLE_INT)
    {
      slab->hdatatype = CCTK_VARIABLE_INT2;
    }
#endif
  }

  /* get the variable's dimension and extents */
  slab->vdim = CCTK_GroupDimFromVarI (vindex);
  extent_int = (int *) malloc (slab->vdim * sizeof (int));
  CCTK_GroupgshVI (GH, slab->vdim, extent_int, vindex);

  /* allocate the arrays all in one go
     only initialize those which are mandatory for the Hyperslab API */
  slab->vectors = (CCTK_INT *) calloc ((slab->vdim + 6) * slab->vdim + 1,
                                       sizeof (CCTK_INT));
  slab->hoffset      = slab->vectors + 0*slab->vdim;
  slab->hsize        = slab->vectors + 1*slab->vdim;
  slab->hsize_chunk  = slab->vectors + 2*slab->vdim;
  slab->origin       = slab->vectors + 3*slab->vdim + 1;
  slab->extent       = slab->vectors + 4*slab->vdim + 1;
  slab->downsample   = slab->vectors + 5*slab->vdim + 1;
  slab->direction    = slab->vectors + 6*slab->vdim + 1;

  for (slab->hdim = 0; slab->hdim < slab->vdim; slab->hdim++)
  {
    slab->extent[slab->hdim] = extent_int[slab->hdim];
    slab->direction[slab->hdim * (slab->vdim + 1)] = 1;

    /* take the downsampling parameters from IOUtil */
    if (slab->hdim == 0)
    {
      slab->downsample[slab->hdim] = out3D_downsample_x;
    }
    else if (slab->hdim == 1)
    {
      slab->downsample[slab->hdim] = out3D_downsample_y;
    }
    else if (slab->hdim == 2)
    {
      slab->downsample[slab->hdim] = out3D_downsample_z;
    }
    else
    {
      slab->downsample[slab->hdim] = 1;
    }
  }

  /* clean up */
  free (extent_int);

  return (slab);
}