aboutsummaryrefslogtreecommitdiff
path: root/src/Write.c
blob: bdb78fc8bac13ee55689ea9fa7969b6324915bb0 (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
/*@@
   @routine    Write.c
   @date       18th September 1999
   @author     Gabrielle Allen
   @desc 
   Writes scalar grid variable data.
   @enddesc 
   @calls    
   @calledby   
   @history
   @hauthor  @hdate 
   @hdesc 
   @hendhistory
@@*/

#include <stdio.h>

#include "cctk.h"
#include "cctk_parameters.h"
#include "iobasicGH.h"

void IOBasic_Write (cGH *GH, int index, const char *alias)
{
  DECLARE_CCTK_PARAMETERS
  char *openmode;
  FILE *file;
  CCTK_REAL *data_real;
  CCTK_INT *data_int;


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

    fullname = CCTK_FullName (index);
    msg = (char *) malloc (200 * sizeof (char) + strlen (fullname));
    sprintf (msg, "No scalar output for '%s' (no storage)", fullname);
    CCTK_WARN (2, msg);
    free (fullname);
    free (msg);
    return;
  }

  /* Open the file (we write only on proc0) */
  if (CCTK_MyProc (GH) == 0) 
  {
    char fname[256];
    iobasicGH *myGH;

    myGH = (iobasicGH *) GH->extensions [CCTK_GHExtensionHandle ("IOBasic")];

    /* see if output files for this alias name were already created */
    if (GetNamedData (myGH->filenameListScalar, alias) == NULL) {
      /* just store a non-NULL pointer in database */
      StoreNamedData (&myGH->filenameListScalar, alias, (void *) 1);
      openmode = "w";
    } else
      openmode = "a";

    sprintf (fname, "%s/%s.tl", myGH->outdirScalar, alias);
    file = fopen(fname,openmode);
    
    if (*openmode == 'w') {
      char title_start_char;

      if (CCTK_Equals(outScalar_style,"gnuplot")) 
        title_start_char = '#';
      else {
        if (! CCTK_Equals(outScalar_style,"xgraph"))
          CCTK_WARN(3,"Don't understand outScalar_style ... using xgraph");
        title_start_char = 34;
      }

      fprintf (file,"%c%s v time\n",title_start_char,alias);
    }

    switch (CCTK_VarTypeI(index)) {
    case CCTK_VARIABLE_REAL:
      data_real  = ((CCTK_REAL ***) GH->data) [index][0]; 
      fprintf(file,"%f %25.13f\n",GH->cctk_time,*data_real);
      break;
    case CCTK_VARIABLE_INT:
      data_int  = ((CCTK_INT ***) GH->data) [index][0]; 
      fprintf(file,"%f %d\n",GH->cctk_time,*data_int);
      break;
    default:
      CCTK_WARN (3, "Unsupported data type");
      return;
    }
    
    fclose(file);

  }

  USE_CCTK_PARAMETERS

}