aboutsummaryrefslogtreecommitdiff
path: root/src/DumpVar.c
blob: f75a0a72978c2ff6c6b997963c4ee53029259a03 (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
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
/*@@
   @file      DumpVar.c
   @desc 
   Do the actual writing of a 3D grid function, for output or for checkpointing
   @enddesc 
   @history
   @hauthor Gabrielle Allen @hdate Oct 17 1998
   @hdesc Split into subroutines and tidied.
   @hauthor Gabrielle Allen @hdate 19 Oct 1998      
   @hdesc Changed names ready for thorn_IO
   @hauthor Thomas Radke @hdate 16 Mar 1999
   @hdesc Converted to Cactus 4.0
          introduced separate downsample values for each dimension
   @hauthor Thomas Radke @hdate 17 Mar 1999
   @hdesc Changed naming: IEEEIO -> IOFlexIO
   @hauthor Thomas Radke @hdate 17 Mar 1999
   @hdesc fixed errors for MPI
   @hauthor Thomas Radke @hdate 12 Apr 1999
   @hdesc Store full variable name in IOFlexIO_AddCommonAttributes()
   @hauthor Thomas Radke @hdate 16 Apr 1999
   @hdesc Added groupname, grouptype, ntimelevels, and current timelevel
          as common attributes
   @hendhistory
   @version $Id$
 @@*/

static char *rcsid = "$Id$";

#include <stdio.h>
#include <stdlib.h>
#ifdef SGI
#include <time.h>
#endif

#include "cctk.h"
#include "cctk_Parameters.h"
#include "CactusPUGH/PUGH/src/include/pugh.h"
#include "CactusBase/IOUtil/src/ioGH.h"
#include "ioFlexGH.h"


#define IOTAGBASE 20000 /* This may break on more than 2000 processors */


static char *char_time_date = NULL;

/* local function prototypes */
static void IOFlexIO_DumpScalar (cGH *GH, int index, int timelevel, IOFile iof,
               int ioflex_type);
static void IOFlexIO_DumpGF (cGH *GH, int index, int timelevel, IOFile iof,
               int ioflex_type, int mpi_type, int element_size);
static void IOFlexIO_DumpArray (cGH *GH, int index, int timelevel, IOFile iof,
               int ioflex_type, int mpi_type, int element_size);
static void IOFlexIO_getDumpData (cGH *GH, int index, int timelevel,                         void **outme, int *free_outme, CCTK_INT4 bnd[9], int element_size);
static void IOFlexIO_eachProcDump (cGH *GH, int index, int timelevel,
               void *outme, CCTK_INT4 bnd [9], IOFile iof, int ioflex_type);
static void IOFlexIO_procDump (IOFile iof, cGH *GH, int index, void *outme,
               CCTK_INT4 bnd [9], int ioflex_type);
static void IOFlexIO_AddChunkAttributes (cGH *GH, int index, CCTK_INT4 bnd [9],
               IOFile iof);
static void IOFlexIO_AddCommonAttributes (cGH *GH, int index, int timelevel,
               CCTK_INT4 gsz [3], IOFile iof);



/*@@
  @routine    IOFlexIO_DumpVar
  @date       16 Apr 1999
  @author     Thomas Radke
  @desc 
              Generic dump routine, just calls the appropriate dump routine
  @enddesc 
  @calls
  @calledby   IOFlexIO_DumpGH IOFlexIO_Write3D

  @var        GH
  @vdesc      Pointer to CCTK grid hierarchy
  @vtype      cGH
  @vio        in
  @endvar
  @var        index
  @vdesc      global index of the variable to be dumped
  @vtype      int
  @vio        in
  @endvar
  @var        timelevel
  @vdesc      the timelevel to store
  @vtype      int
  @vio        in
  @endvar
  @var        iof
  @vdesc      the IEEEIO file to dump to
  @vtype      IOFile
  @vio        in
  @endvar
  @history 
  @endhistory 

@@*/
void IOFlexIO_DumpVar (cGH *GH, int index, int timelevel, IOFile iof)
{
  ioGH *ioUtilGH;
  int variable_type;
  int ioflex_type, mpi_type = 0, element_size;


  /* Get the handle for IOUtil extensions */
  ioUtilGH = (ioGH *) GH->extensions [CCTK_GHExtensionHandle ("IO")];

  variable_type = CCTK_VarTypeI (index);

  switch (variable_type) {
    case CCTK_VARIABLE_REAL:
      ioflex_type = ioUtilGH->out_single ? FLEXIO_REAL4 : FLEXIO_REAL;
      element_size = ioUtilGH->out_single ?
                     sizeof (CCTK_REAL4) : sizeof (CCTK_REAL);
#ifdef MPI
      mpi_type = ioUtilGH->out_single ? PUGH_MPI_REAL4 : PUGH_MPI_REAL;
#endif
      break;

    case CCTK_VARIABLE_INT:
      ioflex_type = FLEXIO_INT;
      element_size = sizeof (CCTK_INT);
#ifdef MPI
      mpi_type = PUGH_MPI_INT;
#endif
      break;

    case CCTK_VARIABLE_CHAR:
      ioflex_type = FLEXIO_CHAR;
      element_size = sizeof (CCTK_CHAR);
#ifdef MPI
      mpi_type = PUGH_MPI_CHAR;
#endif
      break;

#if 0
/* FIXME: Don't know how to support COMPLEX types too !! */
    case CCTK_VARIABLE_COMPLEX:
      ioflex_type = ioUtilGH->out_single ? FLEXIO_REAL4 : FLEXIO_REAL;
      element_size = ioUtilGH->out_single ?
                     2 * sizeof (CCTK_REAL4) : 2 * sizeof (CCTK_REAL);
#ifdef MPI
      mpi_type = ioUtilGH->out_single ? PUGH_MPI_REAL4 : PUGH_MPI_REAL;
#endif
      break;
#endif

    default:
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Unsupported variable type %d", variable_type);
      return;
  }

  switch (CCTK_GroupTypeFromVarI (index)) {
    case GROUP_SCALAR:
      IOFlexIO_DumpScalar (GH, index, timelevel, iof, ioflex_type);
      break;

    case GROUP_ARRAY:
      IOFlexIO_DumpArray (GH, index, timelevel, iof, ioflex_type, mpi_type,
                          element_size);
      break;

    case GROUP_GF:
      IOFlexIO_DumpGF (GH, index, timelevel, iof, ioflex_type, mpi_type,
                       element_size);
      break;

    default:
      CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Invalid group type %d", CCTK_GroupTypeFromVarI (index));
      return;
  }
}


/************************** local routines ******************************/

/*@@
  @routine    IOFlexIO_DumpScalar
  @date       16 Apr 1999
  @author     Thomas Radke
  @desc 
              Dumps a SCALAR type variable into a IEEEIO file
  @enddesc 
  @calls
  @calledby   IOFlexIO_DumpVar

  @var        GH
  @vdesc      Pointer to CCTK grid hierarchy
  @vtype      cGH
  @vio        in
  @endvar
  @var        index
  @vdesc      global index of the variable to be dumped
  @vtype      int
  @vio        in
  @endvar
  @var        timelevel
  @vdesc      the timelevel to store
  @vtype      int
  @vio        in
  @endvar
  @var        iof
  @vdesc      the IEEEIO file to dump to
  @vtype      IOFile
  @vio        in
  @endvar
  @var        ioflex_type
  @vdesc      the IOFlexIO datatype to store
  @vtype      int
  @vio        in
  @endvar
  @history 
  @endhistory 

@@*/
static void IOFlexIO_DumpScalar (cGH *GH, int index, int timelevel, IOFile iof,
                                 int ioflex_type)
{
  ioGH *ioUtilGH;
  CCTK_INT gsz [] = {0, 0, 0};    /* not needed here ? */
  const int dim [] = {1};         /* size of scalar */


  /* Get the handle for IOUtil extensions */
  ioUtilGH = (ioGH *) GH->extensions [CCTK_GHExtensionHandle ("IO")];

  if (CCTK_MyProc (GH) != ioUtilGH->ioproc)
    return;

  /* first dump the data then add the attributes */
  CACTUS_IEEEIO_ERROR (IOwrite (iof, ioflex_type, 1, dim,
                       CCTK_VarDataPtrI (GH, timelevel, index)));
  IOFlexIO_AddCommonAttributes (GH, index, timelevel, gsz, iof);
}


/*@@
  @routine    IOFlexIO_DumpArray
  @date       16 Apr 1999
  @author     Thomas Radke
  @desc 
              Dumps a SCALAR type variable into a IEEEIO file
  @enddesc 
  @calls
  @calledby   IOFlexIO_DumpVar

  @var        GH
  @vdesc      Pointer to CCTK grid hierarchy
  @vtype      cGH
  @vio        in
  @endvar
  @var        index
  @vdesc      global index of the variable to be dumped
  @vtype      int
  @vio        in
  @endvar
  @var        timelevel
  @vdesc      the timelevel to store
  @vtype      int
  @vio        in
  @endvar
  @var        iof
  @vdesc      the IEEEIO file to dump to
  @vtype      IOFile
  @vio        in
  @endvar
  @var        ioflex_type
  @vdesc      the IOFlexIO datatype to store
  @vtype      int
  @vio        in
  @endvar
  @var        mpi_type
  @vdesc      the MPI datatype to communicate
  @vtype      int
  @vio        in
  @endvar
  @var        element_size
  @vdesc      the size of an element of variable
  @vtype      int
  @vio        in
  @endvar
  @history 
  @endhistory 

@@*/
static void IOFlexIO_DumpArray (cGH *GH, int index, int timelevel, IOFile iof,
                                int hdf5io_type, int mpi_type, int element_size)
{
  CCTK_WARN (1, "Dumping of arrays not yet supported");
}


/*@@
  @routine    IOFlexIO_DumpGF
  @date       July 1998
  @author     Paul Walker
  @desc 
  
  @enddesc 
  @calls     
  @calledby   IOFlexIO_DumpVar

  @var        GH
  @vdesc      Pointer to CCTK grid hierarchy
  @vtype      cGH
  @vio        in
  @endvar
  @var        index
  @vdesc      global index of the variable to be dumped
  @vtype      int
  @vio        in
  @endvar
  @var        timelevel
  @vdesc      the timelevel to store
  @vtype      int
  @vio        in
  @endvar
  @var        iof
  @vdesc      the IEEEIO file to dump to
  @vtype      IOFile
  @vio        in
  @endvar
  @var        ioflex_type
  @vdesc      the IOFlexIO datatype to store
  @vtype      int
  @vio        in
  @endvar
  @var        mpi_type
  @vdesc      the MPI datatype to communicate
  @vtype      int
  @vio        in
  @endvar
  @var        element_size
  @vdesc      the size of an element of variable
  @vtype      int
  @vio        in
  @endvar

  @history 
  @hdate Wed Sep  2 10:14:17 1998 @hauthor Tom Goodale
  @hdesc Merged in Szu-Wen's Panda changes 
  @hdate 26 Sep 1998
  @hauthor Gabrielle Allen
  @hdesc Moved into separate file, tidied up, changed parameter names.
         Removed zoom_active flag
  @hauthor    Thomas Radke
  @hdesc      New parameter timelevel specifies the timelevel to dump
  @endhistory 

@@*/
static void IOFlexIO_DumpGF (cGH *GH, int index, int timelevel, IOFile iof,
                             int ioflex_type, int mpi_type, int element_size)
{
  DECLARE_CCTK_PARAMETERS
  ioGH *ioUtilGH;
  pGH *pughGH;
  int myproc;
  int nprocs;
  CCTK_INT4 bnd [9]; /* Lower bounds, size and global size of data we send */
  void *outme;          /* pointer to the data to be dumped */
  int free_outme;       /* indicates whether data buffer needs freeing */
#ifdef SGI
  /* Get the date if we can */
  time_t t = time (NULL);
  char_time_date = asctime (localtime (&t));
#endif
#ifdef MPI
  int i;
  MPI_Status ms;
#endif


  /* Get the handles for pugh and IO extensions */
  ioUtilGH = (ioGH *) GH->extensions [CCTK_GHExtensionHandle ("IO")];
  pughGH = (pGH *) GH->extensions [CCTK_GHExtensionHandle ("PUGH")];

  myproc = CCTK_MyProc (GH);
  nprocs = CCTK_nProcs (GH);

  /* Get the pointer to the data we want to output (includes downsampling) */
  IOFlexIO_getDumpData (GH, index, timelevel, &outme, &free_outme, bnd,
                        element_size);

  if (ioUtilGH->ioproc_every == 1) {
  
    /* 
     * OUTPUT ON EVERY PROCESSOR
     */

    IOFlexIO_eachProcDump (GH, index, timelevel, outme, bnd, iof, ioflex_type);

  }

#ifdef MPI
  else if (myproc != ioUtilGH->ioproc) {

    int outgoing;

    outgoing = bnd [3] * bnd [4] * bnd [5];
      
    /* 
     * OUTPUT ON A SUBSET OF PROCESSORS
     */

    /* Send the geometry and data from the processors which aren't 
     * outputing to the ones that are 
     */

    /* GEOMETRY SEND */
    CACTUS_MPI_ERROR (MPI_Send (bnd, 6, PUGH_MPI_INT4, ioUtilGH->ioproc,
                      2*myproc+IOTAGBASE+1, pughGH->PUGH_COMM_WORLD));
    /* DATA SEND */
    CACTUS_MPI_ERROR (MPI_Send (outme, outgoing, mpi_type, ioUtilGH->ioproc,
                      2*myproc+IOTAGBASE, pughGH->PUGH_COMM_WORLD));
    if (verbose) {
      printf ("Sent %d data points\n", outgoing);
      fflush (stdout);
    }

  } else if (myproc == ioUtilGH->ioproc) {

    int lsz [3];

    /* First calculate the local size and reserve the chunk */
    if (ioUtilGH->ioproc_every >= nprocs) {

      /* One output file ... the local chunk size is the global size */
      lsz [0] = bnd [6]; lsz [1] = bnd [7]; lsz [2] = bnd [8];
    }

    /* Dump data held on output processor */

    if (ioUtilGH->unchunked)
      IOreserveChunk (iof, ioflex_type, 3, lsz);

    /* first the data then the attributes */
    IOFlexIO_procDump (iof, GH, index, outme, bnd, ioflex_type);
    /* delay adding attributes for unchunked files
       until all chunks were written (otherwise attributes get lost !) */
    if (! ioUtilGH->unchunked)
      IOFlexIO_AddCommonAttributes (GH, index, timelevel, &bnd [6], iof);

    /* Dump data from all other processors */

    for (i = myproc+1; i < myproc+ioUtilGH->ioproc_every && i < nprocs; i++) {

      /* Allocate temp to the right size (based on GH->ub [i] et..) */
      void *tmpd;
      int incoming;

      /* Receive bounds information */
      CACTUS_MPI_ERROR (MPI_Recv (bnd, 6, PUGH_MPI_INT4, i, 2*i+IOTAGBASE+1,
                        pughGH->PUGH_COMM_WORLD, &ms));

      incoming = bnd [3] * bnd [4] * bnd [5];

      tmpd = malloc (incoming * element_size);
      CACTUS_MPI_ERROR (MPI_Recv (tmpd, incoming, mpi_type, i,
                        2*i+IOTAGBASE, pughGH->PUGH_COMM_WORLD, &ms));

      IOFlexIO_procDump (iof, GH, index, tmpd, bnd, ioflex_type);
      free (tmpd);
      
    } /* End loop over processors */

    /* now add the attributes for unchunked files */
    if (ioUtilGH->unchunked)
      IOFlexIO_AddCommonAttributes (GH, index, timelevel, &bnd [6], iof);

  } /* End myproc = 0 */

  CACTUS_MPI_ERROR (MPI_Barrier (pughGH->PUGH_COMM_WORLD));
#endif /* MPI */

  if (free_outme)
    free (outme);

  USE_CCTK_PARAMETERS
}


/*@@
   @routine    IOFlexIO_AddCommonAttributes
   @date       July 1998
   @author     Paul Walker
   @desc 
   Add "Common" attributes, these are the GF name, the current date,
   simulation time, origin, bounding box, gridspacings (both downsampled
   and evolution), global grid size, number of processors and iteration
   number. Note that the datestamp should be turned of if you are byte 
   comparing two output files.
   @enddesc 
   @calls     
   @calledby   
   @history 
   @hdate Wed Sep  2 10:15:22 1998 @hauthor Tom Goodale
   @hdesc Abstracted WRITE_ATTRIBUTE to merge in Szu-Wen's Panda changes. 
   @hdate Apr 16 1999 @hauthor Thomas Radke
   @hdesc Added attributes groupname, grouptype, ntimelevels,
          and current timelevel to be stored
   @hdate May 02 1999 @hauthor Thomas Radke
   @hdesc Made chunked attribute nioprocs common so that it can be found by
          the recombiner even for unchunked datasets (eg. SCALARs)
   @hdate May 05 1999 @hauthor Thomas Radke
   @hdesc Added "unchunked" attribute to distinguish between chunked/unchunked
          output files
   @endhistory 
@@*/
static void IOFlexIO_AddCommonAttributes (cGH *GH, int index, int timelevel,
                                          CCTK_INT4 gsz [3], IOFile iof)
{
  DECLARE_CCTK_PARAMETERS
  CCTK_REAL d3_to_IO [6]; /* buffer for writing doubles to IEEEIO */
  CCTK_INT4 i_to_IO;      /* buffer for writing an int to IEEEIO */
  char *name, *gname;
  ioGH *ioUtilGH;


  /* Get the handle for IO extensions */
  ioUtilGH = (ioGH *) GH->extensions [CCTK_GHExtensionHandle ("IO")];

  name = CCTK_FullName (index);
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "name", FLEXIO_CHAR,
                       strlen (name) + 1, name));
  free (name);

  gname = CCTK_GroupNameFromVarI (index);
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "groupname", FLEXIO_CHAR,
                       strlen (gname) + 1, gname));
  free (gname);

  i_to_IO = CCTK_GroupTypeFromVarI (index);
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "grouptype", FLEXIO_INT4,
                       1, &i_to_IO));

  i_to_IO = CCTK_NumTimeLevelsFromVarI (index);
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "ntimelevels", FLEXIO_INT4,
                       1, &i_to_IO));

  i_to_IO = timelevel;
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "timelevel", FLEXIO_INT4,
                       1, &i_to_IO));

  if (char_time_date && out3D_datestamp)
    CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "date", FLEXIO_CHAR,
                         strlen (char_time_date) + 1, char_time_date));

  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "time", FLEXIO_REAL, 1,&GH->cctk_time));

  d3_to_IO [0] = CCTK_CoordOrigin ("x");
  d3_to_IO [1] = CCTK_CoordOrigin ("y");
  d3_to_IO [2] = CCTK_CoordOrigin ("z");
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "origin", FLEXIO_REAL,3,d3_to_IO));
  CCTK_CoordRange (GH, &d3_to_IO [0], &d3_to_IO [3], "x");
  CCTK_CoordRange (GH, &d3_to_IO [1], &d3_to_IO [4], "y");
  CCTK_CoordRange (GH, &d3_to_IO [2], &d3_to_IO [5], "z");
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "min_ext",FLEXIO_REAL,3,d3_to_IO));
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "max_ext",FLEXIO_REAL,3,d3_to_IO + 3));


  d3_to_IO [0] = GH->cctk_delta_space [0] * ioUtilGH->downsample_x;
  d3_to_IO [1] = GH->cctk_delta_space [1] * ioUtilGH->downsample_y;
  d3_to_IO [2] = GH->cctk_delta_space [2] * ioUtilGH->downsample_z;
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "delta", FLEXIO_REAL, 3,d3_to_IO));

  if (ioUtilGH->downsample_x > 1 ||
      ioUtilGH->downsample_y > 1 ||
      ioUtilGH->downsample_z > 1) {
    d3_to_IO [0] = GH->cctk_delta_space [0];
    d3_to_IO [1] = GH->cctk_delta_space [1];
    d3_to_IO [2] = GH->cctk_delta_space [2];
    CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "evolution_delta", FLEXIO_REAL,
                         3, d3_to_IO));
  }

  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "global_size", FLEXIO_INT4,
                       3, gsz));

  i_to_IO = CCTK_nProcs (GH);
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "nprocs", FLEXIO_INT4,
                       1, &i_to_IO));

  i_to_IO = ioUtilGH->ioproc_every;
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "ioproc_every", FLEXIO_INT4,
                       1, &i_to_IO));

  i_to_IO = ioUtilGH->unchunked;
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "unchunked", FLEXIO_INT4,
                       1, &i_to_IO));

  i_to_IO = GH->cctk_iteration;
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "iteration", FLEXIO_INT4,
                       1, &i_to_IO));

  USE_CCTK_PARAMETERS
}


/*@@
  @routine IOFlexIO_AddChunkAttributes
  @author  Paul Walker 
  @date    Feb 1997
  @desc
  This routine adds chunk attributes to a data set. That is,
  is adds attributes to the chunks of data passed to the
  io processors. 
  @enddesc
  @history 
  @hauthor Gabrielle Allen
  @hdate Jul 10 1998
  @hdesc Added the name of the grid function 
  @hdate Wed Sep  2 10:08:31 1998 @hauthor Tom Goodale
  @hdesc Abstracted WRITE_ATTRIBUTE to merge in Szu-Wen's Panda calls.
  @endhistory 
@@*/
static void IOFlexIO_AddChunkAttributes (cGH *GH, int index, CCTK_INT4 bnd [9],
                                         IOFile iof)
{
  char *name;
  CCTK_INT4 i_to_IO;


  /* there is nothing to do for a serial run */
  if (CCTK_nProcs (GH) == 1)
    return;

  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "chunk_origin", FLEXIO_INT4,
                       3, &bnd [0]));
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "subchunk_lb", FLEXIO_INT4,
                       3, &bnd [0]));
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "global_size", FLEXIO_INT4,
                       3, &bnd [6]));

  i_to_IO = GH->cctk_iteration;
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "chunk_dataset", FLEXIO_INT4,
                       1, &i_to_IO));

  name = CCTK_FullName (index);
  CACTUS_IEEEIO_ERROR (IOwriteAttribute (iof, "name", FLEXIO_CHAR,
                       strlen (name)+1, name));
  free (name);

}


/*@@
  @routine IOFlexIO_getDumpData
  @author  Paul Walker 
  @date    Feb 1997
  @desc
  Bounds and data to be output, takes into account downsampling
  @enddesc
  @history 
  @hauthor Gabrielle Allen @hdate Oct 5 1998
  @hdesc Made into subroutine
  @endhistory 
  @var     GH
  @vdesc   Identifies grid hierachy
  @vtype   pGH *
  @vio     in
  @vcomment 
  @endvar 
  @var     GF
  @vdesc   Identifies grid function
  @vtype   pGF *
  @vio     in
  @vcomment 
  @endvar 
  @var     outme
  @vdesc   data segment to output
  @vtype   void **
  @vio     out
  @vcomment This is just GF->data if there is no downsampling
  @endvar 
  @var     free_outme
  @vdesc   Specifies whether or not to free the storage associated with outme
  @vtype   int *
  @vio     out
  @vcomment 1: free storage; 0: don't free storage
  @endvar 
  @var     bnd
  @vdesc   bounds, local size and global shape of the output
  @vtype   CCTK_INT4 [9]
  @vio     out
  @vcomment bnd [0..2] lower bounds and bnd[3],bnd[4],bnd[5] 
            bnd [3..5] local size of data to send
            bnd [6..8] global shape
  @endvar 
  @var     element_size
  @vdesc   the size of an element of variable
  @vtype   int
  @vio     in
  @endvar
@@*/
static void IOFlexIO_getDumpData (cGH *GH, int index, int timelevel,
              void **outme, int *free_outme, CCTK_INT4 bnd[9], int element_size)
{
  DECLARE_CCTK_PARAMETERS
  int i;
  ioGH *ioUtilGH;
  pGH *pughGH;
  CCTK_REAL4 *single_ptr;
  CCTK_REAL *real_ptr;
  CCTK_CHAR *char_ptr;
  CCTK_INT *int_ptr;
  void *data = CCTK_VarDataPtrI (GH, timelevel, index);


  /* to make the compiler happy */
  single_ptr = NULL;
  real_ptr = NULL;
  char_ptr = NULL;
  int_ptr = NULL;

  ioUtilGH = (ioGH *) GH->extensions [CCTK_GHExtensionHandle ("IO")];
  pughGH = (pGH *) GH->extensions [CCTK_GHExtensionHandle ("PUGH")];

  if (ioUtilGH->downsample_x == 1 &&
      ioUtilGH->downsample_y == 1 &&
      ioUtilGH->downsample_z == 1) {

    if (ioUtilGH->out_single) {
      single_ptr = (CCTK_REAL4 *) malloc (pughGH->npoints*sizeof (CCTK_REAL4));

      for (i = 0; i < pughGH->npoints; i++)
        single_ptr [i] = (CCTK_REAL4) ((CCTK_REAL *) data) [i];

      *outme = single_ptr;
      *free_outme = 1;
    } else {
      *outme = data;
      *free_outme = 0;
    }

    for (i = 0; i < 3; i++) {
      bnd [i] = GH->cctk_lbnd[i];     /* the bounds */
      bnd [i+3] = GH->cctk_lsh[i];    /* the sizes */
      bnd [i+6] = GH->cctk_gsh[i];    /* the global space */
    }

  } else {

    int start [3], end [3];
    int i, j, k, l;

    /* Downsampling code ... */
    bnd [6] = GH->cctk_gsh[0] / ioUtilGH->downsample_x;
    if (GH->cctk_gsh[0] % ioUtilGH->downsample_x)
      bnd [6]++;
    bnd [7] = GH->cctk_gsh[1] / ioUtilGH->downsample_y;
    if (GH->cctk_gsh[1] % ioUtilGH->downsample_y)
      bnd [7]++;
    bnd [8] = GH->cctk_gsh[2] / ioUtilGH->downsample_z;
    if (GH->cctk_gsh[2] % ioUtilGH->downsample_z)
      bnd [8]++;

    if (verbose) 
      printf ("Downsampled sizes (%d, %d, %d) -> (%d, %d, %d)\n",
              GH->cctk_gsh[0], GH->cctk_gsh[1], GH->cctk_gsh[2],
              (int) bnd [6], (int) bnd [7], (int) bnd [8]);
    
    /* Now figure out the local downsampling */
    /* The local starts are the lb modded into the downsample */
    for (i = 0; i < 3; i++) {
      int downsample;

      if (i == 0)
        downsample = ioUtilGH->downsample_x;
      else if (i == 1)
        downsample = ioUtilGH->downsample_y;
      else
        downsample = ioUtilGH->downsample_z;

      bnd [i]   = GH->cctk_lbnd[i] / downsample;
      start [i] = bnd [i] * downsample;
      if (start [i] <
          GH->cctk_lbnd[i] + pughGH->ownership [PUGH_NO_STAGGER][0][i]) {
        start [i] += downsample;
        bnd [i] ++;
      }
      end [i]   = ((GH->cctk_lbnd [i] +
                   pughGH->ownership [PUGH_NO_STAGGER][1][i] - 1) / downsample)
                   * downsample;
      bnd [i+3] = (end [i] - start [i]) / downsample + 1;
    }

    if (verbose) {
      printf ("Downsample ranges (%d, %d, %d) -> (%d, %d, %d)\n",
        start [0], start [1], start [2],
        end [0], end [1], end [2]);
      printf ("Local size/bound  (%d, %d, %d)  (%d, %d, %d)\n",
        (int) bnd [3], (int) bnd [4], (int) bnd [5],
        (int) bnd [0], (int) bnd [1], (int) bnd [2]);
    }

    /* compute local ranges */
    for (i = 0; i < 3; i++) {
      start [i] -= GH->cctk_lbnd [i];
      end [i] -= GH->cctk_lbnd [i];
    }

    *outme = malloc (bnd [3] * bnd [4] * bnd [5] * element_size);
    *free_outme = 1;

    /* I hate it to repeat the loops for each case label
       but that way produces much more efficient code */
    l = 0;
    switch (CCTK_VarTypeI (index)) {
      case CCTK_VARIABLE_CHAR:
        char_ptr = (CCTK_CHAR *) *outme;
        for (k = start [2]; k <= end [2]; k += ioUtilGH->downsample_z)
          for (j = start [1]; j <= end [1]; j += ioUtilGH->downsample_y)
            for (i = start [0]; i <= end [0]; i += ioUtilGH->downsample_x)
              char_ptr [l++] = ((CCTK_CHAR *) data) [DATINDEX (pughGH, i, j, k)];
        break;

      case CCTK_VARIABLE_INT:
        int_ptr = (CCTK_INT *) *outme;
        for (k = start [2]; k <= end [2]; k += ioUtilGH->downsample_z)
          for (j = start [1]; j <= end [1]; j += ioUtilGH->downsample_y)
            for (i = start [0]; i <= end [0]; i += ioUtilGH->downsample_x)
              int_ptr [l++] = ((CCTK_INT *) data) [DATINDEX (pughGH, i, j, k)];
        break;

      case CCTK_VARIABLE_REAL:
        if (ioUtilGH->out_single)
          single_ptr = (CCTK_REAL4 *) *outme;
        else
          real_ptr = (CCTK_REAL *) *outme;
        for (k = start [2]; k <= end [2]; k += ioUtilGH->downsample_z)
          for (j = start [1]; j <= end [1]; j += ioUtilGH->downsample_y)
            for (i = start [0]; i <= end [0]; i += ioUtilGH->downsample_x)
              if (ioUtilGH->out_single)
                single_ptr [l++] = (CCTK_REAL4)
                      (((CCTK_REAL *) data) [DATINDEX (pughGH, i, j, k)]);
              else
                real_ptr [l++] =
                      ((CCTK_REAL *) data) [DATINDEX (pughGH, i, j, k)];
        break;

#if 0
/* FIXME: Don't know how to support COMPLEX type too !! */
      case CCTK_VARIABLE_COMPLEX:
        if (ioUtilGH->out_single)
          single_ptr = (CCTK_REAL4 *) *outme;
        else
          cmplx_ptr = (CCTK_COMPLEX *) *outme;
        for (k = start [2]; k <= end [2]; k += ioUtilGH->downsample_z)
          for (j = start [1]; j <= end [1]; j += ioUtilGH->downsample_y)
            for (i = start [0]; i <= end [0]; i += ioUtilGH->downsample_x)
              if (ioUtilGH->out_single) {
                single_ptr [l++] = (CCTK_REAL4)
                       (((CCTK_REAL *) data) [0 + DATINDEX (pughGH, i, j, k)]);
                single_ptr [l++] = (CCTK_REAL4)
                       (((CCTK_REAL *) data) [1 + DATINDEX (pughGH, i, j, k)]);
              } else
                cmplx_ptr [l++] =
                       ((CCTK_COMPLEX *) data) [DATINDEX (pughGH, i, j, k)];
        break;
#endif

      default:
        CCTK_WARN (1, "Unsupported variable type in IOFlexIO_getDumpData");
        return;
    }
  }

  if (verbose) {
    printf ("Global size: %d %d %d\n",
            (int) bnd [6], (int) bnd [7], (int) bnd [8]);
    printf ("Lower bound: %d %d %d\n",
            (int) bnd [0], (int) bnd [1], (int) bnd [2]);
    printf ("Chunk size : %d %d %d\n",
            (int) bnd [3], (int) bnd [4], (int) bnd [5]);
  }

  USE_CCTK_PARAMETERS
}


/*@@
  @routine IOFlexIO_eachProcDump
  @author  Paul Walker 
  @date    Feb 1997
  @desc
  Dump data from each processor
  @enddesc
  @history 
  @hauthor Gabrielle Allen @hdate Oct 5 1998
  @hdesc Made into subroutine
  @endhistory 
@@*/
static void IOFlexIO_eachProcDump (cGH *GH, int index, int timelevel,
                void *outme, CCTK_INT4 bnd [9], IOFile iof, int ioflex_type) 
{
  int chunk_dims [3];


  /* So set up the local shape. */
  chunk_dims [0] = bnd [3];
  chunk_dims [1] = bnd [4];
  chunk_dims [2] = bnd [5];

  /* Dump the data */
  CACTUS_IEEEIO_ERROR (IOwrite (iof, ioflex_type, 3, chunk_dims, outme));
  
  /* Add attributes for global space */
  IOFlexIO_AddCommonAttributes (GH, index, timelevel, &bnd [6], iof);
  
  /* Add chunk attributes */
  IOFlexIO_AddChunkAttributes (GH, index, bnd, iof);
}


/*@@
  @routine IOFlexIO_procDump
  @author  Paul Walker 
  @date    Feb 1997
  @desc
  Dump data 
  @enddesc
  @history 
  @hauthor Gabrielle Allen @hdate Oct 5 1998
  @hdesc Made into subroutine
  @endhistory 
@@*/
static void IOFlexIO_procDump (IOFile iof, cGH *GH, int index, void *outme,
                               CCTK_INT4 bnd [9], int ioflex_type)
{
  int chunk_dims [3], chunk_origin [3];      /* Chunk dim and origin */
  ioGH *ioUtilGH;


  ioUtilGH = (ioGH *) GH->extensions [CCTK_GHExtensionHandle ("IO")];

  chunk_origin [0] = bnd [0];
  chunk_origin [1] = bnd [1];
  chunk_origin [2] = bnd [2];
  chunk_dims [0] = bnd [3];
  chunk_dims [1] = bnd [4];
  chunk_dims [2] = bnd [5];

  if (ioUtilGH->unchunked) {
    
    /* Unchunked data */
    CACTUS_IEEEIO_ERROR (IOwriteChunk (iof, chunk_dims, chunk_origin, outme));
    
  } else {
    
    /* Chunked data */
    CACTUS_IEEEIO_ERROR (IOwrite (iof, ioflex_type, 3, chunk_dims, outme));

    /* Write chunk attributes */
    IOFlexIO_AddChunkAttributes (GH, index, bnd, iof);

  }
}