aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authortradke <tradke@eff87b29-5268-4891-90a3-a07138403961>2002-04-19 08:19:25 +0000
committertradke <tradke@eff87b29-5268-4891-90a3-a07138403961>2002-04-19 08:19:25 +0000
commit0034155495b8a4e58b6f6f32becae974bceaa48c (patch)
treec84cb891011208d097dc45a0d5ab0a4e3dc22e97
parent2d91d4f970e1f31365da061077384eee8a839b28 (diff)
Bugfix for multiple processor runs. The reduction wasn't called as a global
operation so the simulation got stuck. git-svn-id: http://svn.cactuscode.org/arrangements/CactusIO/IOJpeg/trunk@60 eff87b29-5268-4891-90a3-a07138403961
-rw-r--r--src/Write.c56
1 files changed, 29 insertions, 27 deletions
diff --git a/src/Write.c b/src/Write.c
index 64d1928..f4a15db 100644
--- a/src/Write.c
+++ b/src/Write.c
@@ -27,7 +27,8 @@ CCTK_FILEVERSION(CactusIO_IOJpeg_Write_c)
******************** Internal Routines ************************
********************************************************************/
static void WriteData (const cGH *GH, int vindex, const char *alias, int dim,
- int dir, const CCTK_INT hsize[2], const void *hdata);
+ int dir, CCTK_REAL min, CCTK_REAL max,
+ const CCTK_INT hsize[2], const void *hdata);
/*@@
@@ -79,7 +80,14 @@ int IOJpeg_Write (const cGH *GH, CCTK_INT vindex, const char *alias)
CCTK_INT mapping;
CCTK_INT origin[3], extent[3], direction[6], hsize[2];
void *hdata;
+ CCTK_REAL min, max;
const CCTK_INT htype = CCTK_VARIABLE_REAL;
+ /*** FIXME: can CCTK_Reduce() have a 'const cGH *' parameter ?? ***/
+ union
+ {
+ const cGH *const_ptr;
+ cGH *non_const_ptr;
+ } GH_fake_const;
DECLARE_CCTK_PARAMETERS
@@ -101,6 +109,23 @@ int IOJpeg_Write (const cGH *GH, CCTK_INT vindex, const char *alias)
myproc = CCTK_MyProc (GH);
myGH = (ioJpegGH *) CCTK_GHExtension (GH, "IOJpeg");
+ /* set the minimum/maximum for the colormap */
+ if (CCTK_Equals (colormap, "custom"))
+ {
+ min = colormap_min;
+ max = colormap_max;
+ }
+ else
+ {
+ GH_fake_const.const_ptr = GH;
+ i = CCTK_ReductionHandle ("minimum");
+ CCTK_Reduce (GH_fake_const.non_const_ptr, 0, i, 1, CCTK_VARIABLE_REAL,
+ &min, 1, vindex);
+ i = CCTK_ReductionHandle ("maximum");
+ CCTK_Reduce (GH_fake_const.non_const_ptr, 0, i, 1, CCTK_VARIABLE_REAL,
+ &max, 1, vindex);
+ }
+
/* get the number of slices to output */
/* in general: maxdir = gdata.dim * (gdata.dim - 1) / 2; */
maxdir = gdata.dim == 2 ? 1 : 3;
@@ -186,7 +211,7 @@ int IOJpeg_Write (const cGH *GH, CCTK_INT vindex, const char *alias)
{
if (i == 1)
{
- WriteData (GH, vindex, alias, gdata.dim, dir, hsize, hdata);
+ WriteData (GH, vindex, alias, gdata.dim, dir, min, max, hsize, hdata);
}
/* clean up */
@@ -212,44 +237,21 @@ int IOJpeg_Write (const cGH *GH, CCTK_INT vindex, const char *alias)
@enddesc
@@*/
static void WriteData (const cGH *GH, int vindex, const char *alias, int dim,
- int dir, const CCTK_INT hsize[2], const void *hdata)
+ int dir, CCTK_REAL min, CCTK_REAL max,
+ const CCTK_INT hsize[2], const void *hdata)
{
FILE *file;
- CCTK_REAL min, max;
- int reduction_handle;
unsigned char *dataout;
ioJpegGH *myGH;
ioAdvertisedFileDesc advertised_file;
char *filename, *fullname;
char slicename[30];
const char *extensions[] = {"xy", "xz", "yz"};
- /*** FIXME: can CCTK_Reduce() have a 'const cGH *' parameter ?? ***/
- union
- {
- const cGH *const_ptr;
- cGH *non_const_ptr;
- } GH_fake_const;
DECLARE_CCTK_PARAMETERS
- GH_fake_const.const_ptr = GH;
myGH = (ioJpegGH *) CCTK_GHExtension (GH, "IOJpeg");
- if (CCTK_Equals (colormap, "custom"))
- {
- min = colormap_min;
- max = colormap_max;
- }
- else
- {
- reduction_handle = CCTK_ReductionHandle ("minimum");
- CCTK_Reduce (GH_fake_const.non_const_ptr, 0, reduction_handle, 1,
- CCTK_VARIABLE_REAL, &min, 1, vindex);
- reduction_handle = CCTK_ReductionHandle ("maximum");
- CCTK_Reduce (GH_fake_const.non_const_ptr, 0, reduction_handle, 1,
- CCTK_VARIABLE_REAL, &max, 1, vindex);
- }
-
/* allocate the RGB image buffer */
dataout = (unsigned char *) malloc (3 * hsize[0] * hsize[1]);