/*@@ @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 #include #ifdef SGI #include #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 */ void IOFlexIO_DumpScalar (cGH *GH, int index, int timelevel, IOFile iof, int ioflex_type); void IOFlexIO_DumpGF (cGH *GH, int index, int timelevel, IOFile iof, int ioflex_type, int mpi_type, int element_size); void IOFlexIO_DumpArray (cGH *GH, int index, int timelevel, IOFile iof, int ioflex_type, int mpi_type, int element_size); void IOFlexIO_getDumpData (cGH *GH, int index, int timelevel, void **outme, int *free_outme, CCTK_INT4 bnd [9], int element_size); void IOFlexIO_eachProcDump (cGH *GH, int index, int timelevel, void *outme, CCTK_INT4 bnd [9], IOFile iof, int ioflex_type); void IOFlexIO_procDump (IOFile iof, cGH *GH, int index, void *outme, CCTK_INT4 bnd [9], int ioflex_type); void IOFlexIO_AddChunkAttributes (cGH *GH, int index, CCTK_INT4 bnd [9], IOFile iof); 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) { char msg [80]; 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; case CCTK_VARIABLE_COMPLEX: CCTK_WARN (1, "Don't know how to dump complex datatype in IOFlexIO_DumpScalar"); return; default: sprintf (msg, "Unknown variable type %d in IOFlexIO_DumpVar", variable_type); CCTK_WARN (1, msg); 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: sprintf (msg, "Invalid group type %d in IOFlexIO_DumpVar", CCTK_GroupTypeFromVarI (index)); CCTK_WARN (1, msg); return; } } /*@@ @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 @@*/ 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 ? */ int dim [] = {1}; /* size of dimensions for 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 @@*/ 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 @@*/ 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); /* If I'm not vertex centered, bail (for now) */ if (pughGH->stagger != PUGH_VERTEXCTR) { CCTK_WARN (1, "3D Output for staggered variables not available"); return; } /* 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 @@*/ 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 @@*/ 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 @@*/ 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_VERTEXCTR][i][0]) { start [i] += downsample; bnd [i] ++; } end [i] = ((GH->cctk_lbnd [i] + pughGH->ownership [PUGH_VERTEXCTR][i][1] - 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; 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 @@*/ 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 @@*/ 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); } }