aboutsummaryrefslogtreecommitdiff
path: root/src/Write.c
blob: 8c416c971e1852da9f98a3bc87af17486c542b32 (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
/*@@
   @file      Write.c
   @date      Thu 18 April 2002
   @author    Thomas Radke
   @desc
              Output two-dimensional slices in Jpeg image format.
   @enddesc
   @version   $Id$
 @@*/


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

#include "cctk.h"
#include "cctk_Parameters.h"
#include "CactusBase/IOUtil/src/ioutil_AdvertisedFiles.h"
#include "ioJpegGH.h"

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


/********************************************************************
 ********************    Internal Routines   ************************
 ********************************************************************/
static void WriteData (const cGH *GH, int vindex, const char *alias, int dim,
                       int dir, CCTK_REAL min, CCTK_REAL max,
                       const CCTK_INT hsize[2], const void *hdata);
static void *RefineData (CCTK_INT input_dims[2], const void *input_data);


/*@@
   @routine   IOJpeg_Write
   @date      Thu 18 April 2002
   @author    Thomas Radke
   @desc
              Writes the slices of a variable into separate Jpeg files.
   @enddesc
   @calls     Hyperslab_GlobalMappingByIndex
              Hyperslab_FreeMapping
              Hyperslab_GetList
              OpenFile
              WriteData

   @var       GH
   @vdesc     pointer to CCTK GH
   @vtype     const cGH *
   @vio       in
   @endvar
   @var       vindex
   @vdesc     index of variable to output
   @vtype     CCTK_INT
   @vio       in
   @endvar
   @var       alias
   @vdesc     alias name of variable to output
   @vtype     const char *
   @vio       in
   @endvar

   @returntype int
   @returndesc
                0 for success, or<BR>
               -1 if variable has no storage assigned<BR>
               -2 if output file couldn't be opened<BR>
               -3 if hyperslab for coordinates and/or variable couldn't be
                  extracted
   @endreturndesc
@@*/
int IOJpeg_Write (const cGH *GH, CCTK_INT vindex, const char *alias)
{
  const ioJpegGH *myGH;
  int i, mapping, total_hsize;
  int dir, dir_i, dir_j, maxdir, myproc, groupindex;
  cGroup gdata;
  char *fullname;
  int extent_int[3];
  CCTK_INT origin[3], extent[2], direction[6], hsize[2];
  void *hdata, *outdata;
  CCTK_REAL min, max;
  const CCTK_INT htype = CCTK_VARIABLE_REAL;
  DECLARE_CCTK_PARAMETERS


  /* get the variable name and group information */
  fullname = CCTK_FullName (vindex);
  groupindex = CCTK_GroupIndexFromVarI (vindex);
  CCTK_GroupData (groupindex, &gdata);

  /* check if variable has storage assigned */
  if (! CCTK_QueryGroupStorageI (GH, groupindex))
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "No IOJpeg output for '%s' (no storage)", fullname);
    free (fullname);
    return (-1);
  }

  /* get the handle for IOJpeg extensions */
  myproc = CCTK_MyProc (GH);
  myGH = (const ioJpegGH *) CCTK_GHExtension (GH, "IOJpeg");

  /* set the minimum/maximum for the colormap */
  if (CCTK_Equals (colormap, "custom"))
  {
    min = colormap_min;
    max = colormap_max;
  }
  else
  {
    i = CCTK_ReductionHandle ("minimum");
    CCTK_Reduce (GH, 0, i, 1, CCTK_VARIABLE_REAL, &min, 1, vindex);
    i = CCTK_ReductionHandle ("maximum");
    CCTK_Reduce (GH, 0, i, 1, CCTK_VARIABLE_REAL, &max, 1, vindex);
  }

  /* get the number of slices to output */
  /* in general: maxdir = gdata.dim * (gdata.dim - 1) / 2; */
  maxdir = gdata.dim == 2 ? 1 : 3;

  /* get the extents of the variable */
  CCTK_GroupgshVI (GH, 3, extent_int, vindex);

  /* now do the actual I/O looping over all directions */
  for (dir = 0; dir < maxdir; dir++)
  {
    /* get the directions to span the hyperslab */
    if (dir == 0)
    {
      dir_i = 0; dir_j = 1;   /* xy */
    }
    else if (dir == 1)
    {
      dir_i = 0; dir_j = 2;   /* xz */
    }
    else
    {
      dir_i = 1; dir_j = 2;   /* yz */
    }

    /* set the extent vector */
    extent[0] = extent_int[dir_i];
    extent[1] = extent_int[dir_j];

    /* set the origin using the slice center from IOUtil */
    memset (origin, 0, sizeof (origin));
    if (gdata.grouptype == CCTK_GF)
    {
      origin[maxdir-dir-1] = myGH->sp2xyz[gdata.dim - 1][dir];
    }

    /* set the direction vector */
    memset (direction, 0, sizeof (direction));
    direction[dir_i] = direction[gdata.dim + dir_j] = 1;

    mapping = Hyperslab_GlobalMappingByIndex (GH, vindex, 2,
                                              direction, origin, extent,
                                              NULL, -1, NULL, hsize);
    if (mapping < 0)
    {
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Failed to define hyperslab mapping for variable '%s'",
                  fullname);
      continue;
    }
    total_hsize = hsize[0] * hsize[1] * CCTK_VarTypeSize (gdata.vartype);
    if (total_hsize <= 0)
    {
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Selected hyperslab has zero size for variable '%s' "
                  "direction %d", fullname, dir);
      Hyperslab_FreeMapping (mapping);
      continue;
    }

    /* allocate hyperslab buffer */
    hdata = myproc == 0 ? malloc (total_hsize) : NULL;

    /* get the hyperslab */
    i = Hyperslab_Get (GH, mapping, 0, vindex, 0, htype, hdata);
    if (i)
    {
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Failed to extract hyperslab for variable '%s'", fullname);
    }

    /* release the mapping structure */
    Hyperslab_FreeMapping (mapping);

    /* and dump the data to file */
    if (myproc == 0)
    {
      if (i == 0)
      {
        /* refine the 2D slice if requested */
        outdata = refinement_factor > 1 ? RefineData (hsize, hdata) : hdata;
        if (! outdata)
        {
          CCTK_WARN (1, "Failed to refine 2D array");
          outdata = hdata;
        }
        WriteData (GH, vindex, alias, gdata.dim, dir, min, max, hsize, outdata);
        if (outdata != hdata)
        {
          free (outdata);
        }
      }

      /* clean up */
      free (hdata);
    }
  } /* end of looping through xyz directions */

  free (fullname);

  return (0);
}


/********************************************************************
 ********************    Internal Routines   ************************
 ********************************************************************/
/*@@
   @routine    RefineData
   @date       Thu 4 November 2004
   @author     Thomas Radke
   @desc
               Refines a given 2D input array by a certain factor
   @enddesc

   @returntype void *
   @returndesc
               a pointer to the allocated refined output array (must be freed
               by the caller), or NULL in case of an error
   @endreturndesc
 @@*/
static void *RefineData (CCTK_INT input_dims[2], const void *input_data)
{
  const CCTK_REAL coord_origin[2] = {0.0, 0.0};
  CCTK_REAL coord_delta[2];
  const CCTK_INT array_type_codes[2] = {CCTK_VARIABLE_REAL, CCTK_VARIABLE_REAL};
  CCTK_REAL *interp_coords[2], *refined_data;
  int i, j, interp_handle, N_interp_points;
  CCTK_INT output_dims[2];
  DECLARE_CCTK_PARAMETERS


  interp_handle = CCTK_InterpHandle ("uniform cartesian");
  if (interp_handle < 0)
  {
    CCTK_WARN (1, "Couldn't get handle for interpolation operator 'uniform "
               "cartesian'. Did you forget to activate a thorn providing "
               "CCTK_InterpLocalUniform() ?");
    return (NULL);
  }

  coord_delta[0] = coord_delta[1] = refinement_factor;

  output_dims[0] = refinement_factor * (input_dims[0] - 1);
  output_dims[1] = refinement_factor * (input_dims[1] - 1);
  N_interp_points = output_dims[0] * output_dims[1];

  interp_coords[0] = malloc (2 * N_interp_points * sizeof (CCTK_REAL));
  interp_coords[1] = interp_coords[0];
  for (i = 0; i < output_dims[0]; i++)
  {
    for (j = 0; j < output_dims[1]; j++)
    {
      *interp_coords[0]++ = i;
      *interp_coords[1]++ = j;
    }
  }
  interp_coords[0] -= N_interp_points;
  interp_coords[1] -= N_interp_points;

  refined_data = malloc (N_interp_points * sizeof (CCTK_REAL));
  if (CCTK_InterpLocalUniform (2, interp_handle, -1,
                               coord_origin, coord_delta, N_interp_points,
                               CCTK_VARIABLE_REAL,
                               (const void *const *) interp_coords,
                               1, input_dims, array_type_codes,
                               &input_data, 1, array_type_codes,
                               (void *const *) &refined_data) != 0)
  {
    CCTK_WARN (1, "Failed to interpolate 2D array");
    free (refined_data);
    refined_data = NULL;
  }
  else
  {
    input_dims[0] = output_dims[0];
    input_dims[1] = output_dims[1];
  }

  free (interp_coords[0]);

  return (refined_data);
}


/*@@
   @routine    WriteData
   @date       Thu 18 April 2002
   @author     Thomas Radke
   @desc
               Writes the given hyperslab into a Jpeg output file.
   @enddesc
 @@*/
static void WriteData (const cGH *GH, int vindex, const char *alias, int dim,
                       int dir, CCTK_REAL min, CCTK_REAL max,
                       const CCTK_INT hsize[2], const void *hdata)
{
  FILE *file;
  unsigned char *dataout;
  const ioJpegGH *myGH;
  ioAdvertisedFileDesc advertised_file;
  char *filename, *fullname;
  char slicename[30];
  const char *extensions[] = {"xy", "xz", "yz"};
  DECLARE_CCTK_PARAMETERS


  myGH = (const ioJpegGH *) CCTK_GHExtension (GH, "IOJpeg");

  /* allocate the RGB image buffer */
  dataout = (unsigned char *) malloc (3 * hsize[0] * hsize[1]);

  AutoColorDataSlice (hsize[0], hsize[1], hdata, dataout, min, max,
                      colormap_bias, colormap_factor);

  /* open the file */
  if (dim == 2)
  {
    strcpy (slicename, "2D");
  }
  else
  {
    /* give the slice origin as range [1 .. n] */
    sprintf (slicename, "%s_[%d]", extensions[dir], myGH->sp2xyz[dim-1][dir]);
  }

  filename = (char *) malloc (strlen (myGH->out_dir) + strlen (alias) +
                              sizeof (slicename) + 20);

  if (CCTK_Equals (mode, "remove"))
  {
    sprintf (filename, "%s%s_%s.jpeg", myGH->out_dir, alias, slicename);
  }
  else
  {
    sprintf (filename, "%s%s_%s.it_%d.jpeg", myGH->out_dir, alias, slicename,
             GH->cctk_iteration);
  }

  file = fopen (filename, "w");
  if (file)
  {
    /* write the data */
    WriteJPEGToFileRGB (hsize[0], hsize[1], dataout, colormap_quality, file);

    /* advertise the file for downloading */
    if (CCTK_Equals (mode, "remove") && myGH->out_last[vindex] < 0)
    {
      fullname = CCTK_FullName (vindex);
      advertised_file.slice = slicename;
      advertised_file.thorn = CCTK_THORNSTRING;
      advertised_file.varname = fullname;
      advertised_file.description = "Jpegs of slices";
      advertised_file.mimetype = "image/jpeg";

      IOUtil_AdvertiseFile (GH, filename, &advertised_file);

      free (fullname);
    }

    /* close the file */
    fclose (file);
  }
  else
  {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Cannot open IOJpeg output file '%s'", filename);
  }

  /* clean up */
  free (dataout);
  free (filename);
}