aboutsummaryrefslogtreecommitdiff
path: root/src/OutputScalar.c
blob: a744164d6796cac06c06030d125cc47d30e788dd (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
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
 /*@@
   @file      Output.c
   @date      Mon 21 September
   @author    Gabrielle Allen
   @desc
              Functions to deal with scalar output of grid variables
   @enddesc
 @@*/

/* #define IOBASIC_DEBUG 1 */

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

#include "cctk.h"
#include "cctk_Parameters.h"
#include "util_String.h"
#include "iobasicGH.h"

static const char *rcsid = "$Header$";

CCTK_FILEVERSION(CactusBase_IOBasic_OutputScalar_c)


/********************************************************************
 ********************    Internal Routines   ************************
 ********************************************************************/
static void CheckSteerableParameters (iobasicGH *myGH);


/*@@
   @routine    IOBasic_ScalarOutputGH
   @date       Sat March 6 1999
   @author     Gabrielle Allen
   @desc
               Loops over all variables and outputs them if necessary
   @enddesc
   @calls      IOBasic_TimeForScalarOutput
               IOBasic_WriteScalar

   @var        GH
   @vdesc      Pointer to CCTK GH
   @vtype      const cGH *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               the number of variables which were output at this iteration
               (or 0 if it wasn't time to output yet)
   @endreturndesc
@@*/
int IOBasic_ScalarOutputGH (const cGH *GH)
{
  int vindex, retval;
  iobasicGH *myGH;


  myGH = CCTK_GHExtension (GH, "IOBasic");
  CheckSteerableParameters (myGH);

  /* return if no output is required */
  if (myGH->outScalar_every <= 0)
  {
    return (0);
  }

  retval = 0;

  /* loop over all variables */
  for (vindex = CCTK_NumVars () - 1; vindex >= 0; vindex--)
  {
    if (IOBasic_TimeForScalarOutput (GH, vindex) &&
        IOBasic_WriteScalar (GH, vindex, CCTK_VarName (vindex)) == 0)
    {
      /* register variable as having output this iteration */
      myGH->outScalar_last[vindex] = GH->cctk_iteration;
      retval++;
    }
  }

  return (retval);
}


/*@@
   @routine    IOBasic_OutputVarAs
   @date       Sun 8 Mar 2003
   @author     Thomas Radke
   @desc
               Unconditional output of a variable using the "Scalar" I/O method
   @enddesc
   @calls      IOBasic_WriteScalar

   @var        GH
   @vdesc      pointer to CCTK GH
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        fullname
   @vdesc      complete name of variable to output, maybe appended by an
               options string
   @vtype      const char *
   @vio        in
   @endvar
   @var        alias
   @vdesc      variable's alias name used to build the output filename
   @vtype      const char *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine IOBasic_WriteScalar, or<BR>
               -1 if variable cannot be output by "Scalar"
   @endreturndesc
@@*/
int IOBasic_OutputVarAs (const cGH *GH, const char *fullname, const char *alias)
{
  int vindex, old_num_reductions, retval;
  char *varname, *tmp;
  iobasic_parseinfo_t info;
  iobasic_reduction_t *reduction, *old_reductions, *next;
  iobasicGH *myGH;
  DECLARE_CCTK_PARAMETERS


  /* get the variable index from 'fullname' */
  varname = strdup (fullname);
  tmp = strchr (fullname, '[');
  if (tmp)
  {
    varname[tmp - fullname] = 0;
  }

  vindex = CCTK_VarIndex (varname);
  free (varname);
  if (vindex < -1)
  {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "IOBasic_OutputVarAs: invalid variable name in 'fullname' "
                "argument '%s'", fullname);
    return (-1);
  }

  /* save the old scalar_reductions[vindex] info ... */
  myGH = CCTK_GHExtension (GH, "IOBasic");
  old_num_reductions = myGH->scalar_reductions[vindex].num_reductions;
  old_reductions = myGH->scalar_reductions[vindex].reductions;

  /* ... and create a temporary one from 'fullname' */
  retval = 0;
  info.reduction_list = myGH->scalar_reductions;
  info.reductions_string = outScalar_reductions;
  if (CCTK_TraverseString (fullname, IOBasic_AssignReductionList, &info,
                           CCTK_VAR) < 0)
  {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "IOBasic_OutputVarAs: failed to parse 'fullname' argument '%s'",
                fullname);
    retval = -1;
  }
  else if (myGH->scalar_reductions[vindex].num_reductions > 0)
  {
    /* do the actual output */
    retval = IOBasic_WriteScalar (GH, vindex, alias);

    /* free temporary scalar_reductions[vindex] info and restore the old one */
    reduction = myGH->scalar_reductions[vindex].reductions;
    while (reduction)
    {
      next = reduction->next;
      free (reduction->name);
      free (reduction);
      reduction = next;
    }
  }
  myGH->scalar_reductions[vindex].num_reductions = old_num_reductions;
  myGH->scalar_reductions[vindex].reductions = old_reductions;

  return (retval);
}


 /*@@
   @routine    IOBasic_TimeForScalarOutput
   @date       Sat March 6 1999
   @author     Gabrielle Allen
   @desc
               Decides if it is time to do Scalar output.
   @enddesc

   @var        GH
   @vdesc      Pointer to CCTK GH
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        vindex
   @vdesc      index of variable
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               1 if output should take place at this iteration, or<BR>
               0 if not
   @endreturndesc
@@*/
int IOBasic_TimeForScalarOutput (const cGH *GH, int vindex)
{
  int retval;
  char *fullname;
  iobasicGH *myGH;


  /* Get the GH extension for IOBasic */
  myGH = CCTK_GHExtension (GH, "IOBasic");

  CheckSteerableParameters (myGH);

  /* check if any output was requested */
  if (myGH->outScalar_every <= 0 || GH->cctk_iteration % myGH->outScalar_every
      || myGH->scalar_reductions[vindex].num_reductions == 0)
  {
    retval = 0;
  }
  else
  {
    /* Check if variable wasn't already output this iteration */
    retval = myGH->outScalar_last[vindex] != GH->cctk_iteration;
    if (! retval)
    {
      fullname = CCTK_FullName (vindex);
      CCTK_VWarn (5, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Already done scalar output for '%s' in current iteration "
                  "(probably via triggers)", fullname);
      free (fullname);
    }
  }

  return (retval);
}


/*@@
   @routine    IOBasic_TriggerScalarOutput
   @date       Sat March 6 1999
   @author     Gabrielle Allen
   @desc
               Triggers the output of a variable using IOBasic's Scalar
               output method
   @enddesc

   @var        GH
   @vdesc      Pointer to CCTK GH
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        vindex
   @vdesc      index of variable to output
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine IOBasic_WriteScalar
   @endreturndesc
@@*/
int IOBasic_TriggerScalarOutput (const cGH *GH, int vindex)
{
  int retval;
  const char *name;
  iobasicGH *myGH;


  /* Get the GH extension for IOBasic */
  myGH = CCTK_GHExtension (GH, "IOBasic");

  name = CCTK_VarName (vindex);

#ifdef IOBASIC_DEBUG
  printf ("\nIn IOBasic_TriggerScalarOutput\n---------------------\n");
  printf ("  Index = %d\n", vindex);
  printf ("  Variable = -%s-\n", name);
#endif

  /* do the Scalar output */
  retval = IOBasic_WriteScalar (GH, vindex, name);
  if (retval == 0)
  {
    /* register variable as having Scalar output this iteration */
    myGH->outScalar_last[vindex] = GH->cctk_iteration;
  }

  return (retval);
}


/********************************************************************
 ********************    Internal Routines   ************************
 ********************************************************************/
 /*@@
   @routine    CheckSteerableParameters
   @date       Tue 31 Jul 2001
   @author     Thomas Radke
   @desc
               Re-evaluates 'IOBasic::outScalar_every' and/or 'IO::out_every'
               resp. to set myGH->outScalar_every to the frequency of scalar
               output. Re-evaluates 'IOBasic::outScalar_vars' and
               'IOBasic::outScalar_reductions'.
   @enddesc
   @calls      CCTK_ParameterQueryTimesSet

   @var        myGH
   @vdesc      Pointer to IOBasic's GH extension
   @vtype      iobasicGH *
   @vio        in
   @endvar
@@*/
static void CheckSteerableParameters (iobasicGH *myGH)
{
  int vindex, out_old, times_set, update_reductions_list, num_vars;
  iobasic_reduction_t *reduction, *next;
  static int outScalar_vars_lastset       = -1;
  static int outScalar_reductions_lastset = -1;
  iobasic_parseinfo_t info;
  char *msg, *oldmsg, *fullname;
  DECLARE_CCTK_PARAMETERS


  /* how often to output */
  out_old = myGH->outScalar_every;
  myGH->outScalar_every = outScalar_every >= 0 ? outScalar_every : out_every;

  /* report if frequency changed */
  if (myGH->outScalar_every != out_old && ! CCTK_Equals (verbose, "none"))
  {
    if (myGH->outScalar_every > 0)
    {
      CCTK_VInfo (CCTK_THORNSTRING, "Scalar: Periodic output every %d "
                  "iterations", myGH->outScalar_every);
    }
    else
    {
      CCTK_INFO ("Scalar: Periodic output turned off");
    }
  }

  /* return if there's nothing to do */
  if (myGH->outScalar_every <= 0)
  {
    return;
  }

  /* check if the 'outScalar_reductions' parameter if it was changed */
  times_set = CCTK_ParameterQueryTimesSet ("outScalar_reductions",
                                           CCTK_THORNSTRING);
  update_reductions_list = times_set != outScalar_reductions_lastset;
  outScalar_reductions_lastset = times_set;

  /* check if the 'outScalar_vars' parameter if it was changed */
  times_set = CCTK_ParameterQueryTimesSet ("outScalar_vars", CCTK_THORNSTRING);
  update_reductions_list |= times_set != outScalar_vars_lastset;
  outScalar_vars_lastset = times_set;

  if (update_reductions_list)
  {
    /* free old reduction lists ... */
    num_vars = CCTK_NumVars ();
    for (vindex = 0; vindex < num_vars; vindex++)
    {
      if (myGH->scalar_reductions[vindex].num_reductions > 0)
      {
        myGH->scalar_reductions[vindex].num_reductions = 0;
        reduction = myGH->scalar_reductions[vindex].reductions;
        while (reduction)
        {
          next = reduction->next;
          free (reduction->name);
          free (reduction);
          reduction = next;
        }
      }
    }

    /* ... and create new ones */
    info.reduction_list = myGH->scalar_reductions;
    info.reductions_string = outScalar_reductions;
    if (CCTK_TraverseString (outScalar_vars, IOBasic_AssignReductionList, &info,
                             CCTK_GROUP_OR_VAR) < 0)
    {
      CCTK_WARN (1, "Failed to parse 'IOBasic::outScalar_vars' parameter");
    }
    else if (! CCTK_Equals (verbose, "none"))
    {
      msg = NULL;
      for (vindex = 0; vindex < num_vars; vindex++)
      {
        if (myGH->scalar_reductions[vindex].num_reductions)
        {
          fullname = CCTK_FullName (vindex);
          if (! msg)
          {
            Util_asprintf (&msg, "Scalar: Periodic output requested for '%s'",
                           fullname);
          }
          else
          {
            oldmsg = msg;
            Util_asprintf (&msg, "%s, '%s'", oldmsg, fullname);
            free (oldmsg);
          }
          free (fullname);
        }
      }
      if (msg)
      {
        CCTK_INFO (msg);
        free (msg);
      }
    }
  }
}