aboutsummaryrefslogtreecommitdiff
path: root/src/WriteInfo.c
diff options
context:
space:
mode:
authoryye00 <yye00@b589c3ab-70e8-4b4d-a09f-cba2dd200880>2005-11-01 18:50:38 +0000
committeryye00 <yye00@b589c3ab-70e8-4b4d-a09f-cba2dd200880>2005-11-01 18:50:38 +0000
commitffa02dbbfb29511a6d60d2823fa36775fb6027d2 (patch)
tree7552784df92d6d1ddc21754e696c94f23fe28701 /src/WriteInfo.c
parentd675eee199ebdf7f889ab6780c56376d313eaeb0 (diff)
commits to change IOBasic to use the new reduction interface. This breaks the testsuite for testcomplex because of a new reduction scheme for complex numbers. WaveToy1DF77 also fails with inconsistency, failing at 3,7 processors with varying files and differences between them. Currently investigating these failurse in detail.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/IOBasic/trunk@169 b589c3ab-70e8-4b4d-a09f-cba2dd200880
Diffstat (limited to 'src/WriteInfo.c')
-rw-r--r--src/WriteInfo.c20
1 files changed, 14 insertions, 6 deletions
diff --git a/src/WriteInfo.c b/src/WriteInfo.c
index 163a3c7..cb73a24 100644
--- a/src/WriteInfo.c
+++ b/src/WriteInfo.c
@@ -30,7 +30,7 @@ CCTK_FILEVERSION(CactusBase_IOBasic_WriteInfo_c)
@enddesc
@calls CCTK_VarDataPtrI
- CCTK_Reduce
+ CCTK_ReduceGridArrays
@var GH
@vdesc Pointer to CCTK grid hierarchy
@@ -61,8 +61,16 @@ int IOBasic_WriteInfo (const cGH *GH, int vindex)
iobasicGH *myGH;
void *ptr;
iobasic_reduction_t *reduction;
-
-
+
+ int input_array[1];
+ CCTK_INT output_array_type_codes[1];
+ CCTK_REAL output_value;
+ void * reduction_value[1];
+
+ input_array[0] = vindex;
+ output_array_type_codes[0] = CCTK_VARIABLE_REAL;
+ reduction_value[0] = &output_value;
+
myGH = CCTK_GHExtension (GH, "IOBasic");
reduction = myGH->info_reductions[vindex].reductions;
@@ -169,9 +177,9 @@ int IOBasic_WriteInfo (const cGH *GH, int vindex)
/* for CCTK_GF and CCTK_ARRAY variables: loop over all reductions */
while (reduction)
{
- reduction->is_valid = CCTK_Reduce (GH, -1, reduction->handle, 1,
- CCTK_VARIABLE_REAL,
- &reduction->value, 1, vindex) == 0;
+ reduction->is_valid = CCTK_ReduceGridArrays (GH, -1, reduction->handle, -1, 1, input_array , 1, output_array_type_codes, reduction_value) == 0;
+ reduction->value = output_value;
+
if (! reduction->is_valid)
{
CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,