diff options
author | tradke <tradke@94b1c47f-dcfd-45ef-a468-0854c0e9e350> | 2002-04-17 08:40:21 +0000 |
---|---|---|
committer | tradke <tradke@94b1c47f-dcfd-45ef-a468-0854c0e9e350> | 2002-04-17 08:40:21 +0000 |
commit | bed4e30fff00ce199f97fb781999fa332798ce8a (patch) | |
tree | 9756a8b1a19667efab70617cb40531f6337b6cde | |
parent | 8c87ce629147b4040461ed28048fb15ca1f4a2c4 (diff) |
Bugfix for computing the hyperslab extent. This caused a few testsuites
(eg. BenchADM) to fail because of missing lines in the 1D output files.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/IOASCII/trunk@106 94b1c47f-dcfd-45ef-a468-0854c0e9e350
-rw-r--r-- | src/Write1D.c | 55 |
1 files changed, 24 insertions, 31 deletions
diff --git a/src/Write1D.c b/src/Write1D.c index b497d89..7905c1c 100644 --- a/src/Write1D.c +++ b/src/Write1D.c @@ -144,8 +144,8 @@ int IOASCII_Write1D (const cGH *GH, int vindex, const char *alias) FILE *file[8]; CCTK_REAL offset; CCTK_REAL coord_lower[3]; - CCTK_INT *origin, *direction, *extent; - CCTK_INT mapping, hsize; + CCTK_INT *origin, *direction; + CCTK_INT mapping, hsize, extent; CCTK_INT vindices[2]; void *hdata[2]; DECLARE_CCTK_PARAMETERS @@ -284,12 +284,22 @@ int IOASCII_Write1D (const cGH *GH, int vindex, const char *alias) have_coords ? coord_lower : NULL, num_files, file); } - origin = (CCTK_INT *) calloc (3*gdata.dim, sizeof (CCTK_INT)); - direction = origin + 1*gdata.dim; - extent = origin + 2*gdata.dim; - extent_int = (int *) malloc (gdata.dim * sizeof (int)); + origin = (CCTK_INT *) calloc (2*gdata.dim, sizeof (CCTK_INT)); + direction = origin + gdata.dim; + extent_int = (int *) malloc ((gdata.dim + 1) * sizeof (int)); + /* get the variable's extents, compute the extent for 3D-diagonals as the + minimum of grid points in each direction */ CCTK_GroupgshVI (GH, gdata.dim, extent_int, vindex); + if (gdata.dim == 3) + { + extent_int[3] = extent_int[0] < extent_int[1] ? + extent_int[0] : extent_int[1]; + if (extent_int[2] < extent_int[3]) + { + extent_int[3] = extent_int[2]; + } + } /* now do the actual I/O looping over all directions */ for (dir = 0; dir < 4; dir++) @@ -310,11 +320,8 @@ int IOASCII_Write1D (const cGH *GH, int vindex, const char *alias) direction[i] = (dir == i || dir == 3) ? 1 : 0; } - /* set the extent vector */ - for (i = 0; i < gdata.dim; i++) - { - extent[i] = extent_int[i]; - } + /* set the extent */ + extent = extent_int[dir]; /* set the origin of the line */ if (gdata.grouptype == CCTK_GF && dir < 3) @@ -322,15 +329,15 @@ int IOASCII_Write1D (const cGH *GH, int vindex, const char *alias) for (i = 0; i < gdata.dim; i++) { origin[i] = myGH->spxyz[gdata.dim-1][dir][i]; - extent[i] -= origin[i]; } + extent -= origin[dir]; #if 0 if (CCTK_CoordRangePhysIndex (GH, &lower, &upper, dir + 1, NULL, coord_system) >= 0) { -fprintf (stderr, "CCTK_CoordRangePhysIndex() resets origin/extent in dir %d from %d/%d to ", dir, origin[dir], extent[dir]); - origin[dir] = lower; extent[dir] = upper - lower + 1; -fprintf (stderr, "%d/%d\n", origin[dir], extent[dir]); +fprintf (stderr, "CCTK_CoordRangePhysIndex() resets origin/extent in dir %d from %d/%d to ", dir, origin[dir], extent); + origin[dir] = lower; extent = upper - lower + 1; +fprintf (stderr, "%d/%d\n", origin[dir], extent); } else { @@ -340,29 +347,15 @@ fprintf (stderr, "%d/%d\n", origin[dir], extent[dir]); } #endif } - else /* origin for CCTK_ARRAYS and diagonals is always (0, 0, 0) */ + else /* origin for CCTK_ARRAYS is always (0, 0, 0) */ { memset (origin, 0, gdata.dim * sizeof (CCTK_INT)); - - /* for diagonals: get the extent as the minimum of points - in all direction */ - if (dir == 3) - { - if (extent[0] > extent[1]) - { - extent[0] = extent[1]; - } - if (extent[0] > extent[2]) - { - extent[0] = extent[2]; - } - } } mapping = Hyperslab_DefineGlobalMappingByIndex (GH, vindex, 1, direction, origin, - extent, + &extent, NULL, /* downsample */ -1, /* table handle */ NULL /* conversion fn */, |