aboutsummaryrefslogtreecommitdiff
path: root/src/WriteGF.c
blob: 89cff1ee96922824a256d3c693e3e479e4245be1 (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
/*@@
   @routine    WriteGF
   @date       Tue Apr  1 16:45:35 1997
   @author     Paul Walker
   @desc 
   Dumps the scalar data.
   @enddesc 
   @calls    
   @calledby   
   @history
   @hauthor Thomas Radke @hdate 17 Mar 1999
   @hdesc included "cctk_Comm.h" to overload CCTK_GetMyProc() properly
   @hendhistory
@@*/

#include <stdio.h>
#include <stdlib.h>
#include <string.h>      /* strlen(3) */
#include <sys/types.h>
#include <sys/stat.h>    /* stat(2) */

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


void IOBasic_WriteGF (cGH *GH, int index, const char *alias)
{
  DECLARE_CCTK_PARAMETERS
  int ierr;
  int i;
  iobasicGH *myGH;
  FILE *file;
  char *openmode;
  char *filename;
  char title_start_char;
  int reduction_handle;
  char *format_str;
  struct stat fileinfo;
  
  CCTK_REAL reduction_value;
  /* add more reductions to this table if you want, telling it
       - the reduction operator
       - the label for the file header
       - the filename extension (suffixed by ".tl" for trace line) */
  const struct 
  {
    char *operator;
    char *label;
    char *extension;
  } reductions [] = {
    {"minimum", "min", "min"},
    {"maximum", "max", "max"},
    {"norm1", "norm1", "nm1"},
    {"norm2", "norm2", "nm2"},
  };

#define NUM_REDUCTIONS  (sizeof (reductions) / sizeof (reductions [0]))


  /* first, check if variable has storage assigned */
  if (! CCTK_QueryGroupStorageI (GH, CCTK_GroupIndexFromVarI (index)))
  {
    char *fullname;

    fullname = CCTK_FullName (index);
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "IOBasic_WriteGF: No 0D output for '%s' (no storage)", 
		fullname);
    free (fullname);
    return;
  }

  /* check if decimal or exponential format desired */
  format_str = CCTK_Equals (out_format, "f") ? "%f %25.13f\n" : "%e %25.13e\n";

  /* set the output style */
  title_start_char = CCTK_Equals (outScalar_style, "gnuplot") ? '#' : 34;

  /* get the GH extension handle for IOBasic */
  myGH = (iobasicGH *) GH->extensions [CCTK_GHExtensionHandle ("IOBasic")];

  /* the extra characters should be sufficient 
     for the longest filename extension */
  filename = (char *)
    malloc (strlen (myGH->outdirScalar) + strlen (alias) + 20);

  /* now loop over all IOBasic reduction operations */
  for (i = 0; i < NUM_REDUCTIONS; i++)
  {

    /* get the reduction handle from the reduction operator */
    reduction_handle = CCTK_ReductionHandle (reductions [i].operator);
    if (reduction_handle < 0) 
    {
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "IOBasic_WriteGF: Invalid reduction operator '%s'", 
		  reductions [i].operator);
      continue;
    }

    /* do the reduction (all processors) */
    ierr = CCTK_Reduce (GH, 0, reduction_handle, 1, CCTK_VARIABLE_REAL,
			&reduction_value, 1, index);

    /* dump the reduction value to file by processor 0 */
    if (ierr == 0 && CCTK_MyProc (GH) == 0)
    {

      /* build the filename */
      sprintf (filename, "%s/%s_%s.tl", myGH->outdirScalar, alias,
               reductions [i].extension);

      /* see if output files for this alias name were already created */
      if (GetNamedData (myGH->filenameListScalar, filename) == NULL)
      {
        /* if restart from recovery, all existing files are opened
           in append mode */
        if (IOUtil_RestartFromRecovery (GH))
        {
          openmode = stat (filename, &fileinfo) == 0 ? "a" : "w";
        }
        else
        {
          openmode = "w";
        }
      }
      else
      {
        openmode = "a";
      }

      file = fopen (filename, openmode);
      if (file == NULL)
      {
        CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "IOBasic_WriteGF: Could not open output file '%s'", 
		    filename);
        continue;
      }

      /* when creating the file, write the header,
         advertise the file, and save the filename in the database */
      if (*openmode == 'w')
      {
        ioAdvertisedFileDesc advertised_file;


        fprintf (file, "%c%s %s v time\n", title_start_char, alias,
                 reductions [i].label);

        /* advertise the file for downloading */
        advertised_file.slice       = reductions [i].extension;
        advertised_file.thorn       = CCTK_THORNSTRING;
        advertised_file.varname     = CCTK_FullName (index);
        advertised_file.description = "Reduction on Grid Functions";
        advertised_file.mimetype    = CCTK_Equals(outScalar_style,"gnuplot") ?
	  "application/gnuplot" : "application/x-graph";

        IOUtil_AdvertiseFile (GH, filename, &advertised_file);

        free (advertised_file.varname);

        /* just store a non-NULL pointer in database */
        StoreNamedData (&myGH->filenameListScalar, filename, (void *) 1);
      }

      /* write the data and close the file */
      fprintf (file, format_str, GH->cctk_time, reduction_value);
      fclose (file);
    }
  }

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