aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorallen <allen@b589c3ab-70e8-4b4d-a09f-cba2dd200880>2000-07-19 12:20:56 +0000
committerallen <allen@b589c3ab-70e8-4b4d-a09f-cba2dd200880>2000-07-19 12:20:56 +0000
commit5ef5219e049640848f134969fec13edf850606c7 (patch)
tree34d6f00ec5b9ae1600c7cb25eecb54a815b35faf /src
parent4708b9c442dcf8919a614866d599e8a0a239ebcf (diff)
Only output scalar's if their values have been successfully
calculated git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/IOBasic/trunk@63 b589c3ab-70e8-4b4d-a09f-cba2dd200880
Diffstat (limited to 'src')
-rw-r--r--src/OutputInfo.c55
-rw-r--r--src/Write.c5
-rw-r--r--src/WriteGF.c10
-rw-r--r--src/WriteInfo.c33
-rw-r--r--src/iobasicGH.h2
5 files changed, 75 insertions, 30 deletions
diff --git a/src/OutputInfo.c b/src/OutputInfo.c
index 70bb8a1..b0d8c0e 100644
--- a/src/OutputInfo.c
+++ b/src/OutputInfo.c
@@ -20,7 +20,8 @@
#include "iobasicGH.h"
/* extern prototypes */
-CCTK_REAL IOBasic_WriteInfo (cGH *GH, int index, const char *operator,
+int IOBasic_WriteInfo (cGH *GH, CCTK_REAL *val,int index,
+ const char *operator,
const char *alias);
@@ -58,10 +59,12 @@ int IOBasic_OutputInfoGH (cGH *GH)
DECLARE_CCTK_PARAMETERS
int i;
+ int ierr1=-1;
+ int ierr2=-1;
const char *name;
iobasicGH *myGH;
const cParamData *paramdata;
-
+ CCTK_REAL val;
/* Get the GH extensions for IOBasic */
myGH = (iobasicGH *) GH->extensions [CCTK_GHExtensionHandle ("IOBasic")];
@@ -69,7 +72,9 @@ int IOBasic_OutputInfoGH (cGH *GH)
/* How often to output */
myGH->outInfo_every = out_every > 0 ? out_every : -1;
if (outInfo_every > 0)
+ {
myGH->outInfo_every = outInfo_every;
+ }
/* Return if no output is required */
if (myGH->outInfo_every < 1 ||
@@ -168,29 +173,46 @@ int IOBasic_OutputInfoGH (cGH *GH)
#endif
/* Make the IO call */
- myGH->infovals[i][0] = IOBasic_WriteInfo (GH, i, "minimum", name);
- myGH->infovals[i][1] = IOBasic_WriteInfo (GH, i, "maximum", name);
+ ierr1 = IOBasic_WriteInfo (GH, &val, i, "minimum", name);
+ myGH->infovals[i][0] = val;
+ ierr2 = IOBasic_WriteInfo (GH, &val, i, "maximum", name);
+ myGH->infovals[i][1] = val;
/* Register GF as having info output this iteration */
myGH->outInfo_last [i] = GH->cctk_iteration;
}
- /* finally print out the stuff */
- if (USE_DECIMAL_NOTATION (myGH->infovals [i][0]))
+ if (ierr1 == 0)
{
- printf ("%12.8f |", myGH->infovals [i][0]);
+ /* finally print out the stuff */
+ if (USE_DECIMAL_NOTATION (myGH->infovals [i][0]))
+ {
+ printf ("%12.8f |", myGH->infovals [i][0]);
+ }
+ else
+ {
+ printf ("%10.6e |", myGH->infovals [i][0]);
+ }
}
else
{
- printf ("%10.6e |", myGH->infovals [i][0]);
+ printf(" ----------- |");
}
- if (USE_DECIMAL_NOTATION (myGH->infovals [i][1]))
+
+ if (ierr2 == 0)
{
- printf ("%12.8f |", myGH->infovals [i][1]);
+ if (USE_DECIMAL_NOTATION (myGH->infovals [i][1]))
+ {
+ printf ("%12.8f |", myGH->infovals [i][1]);
+ }
+ else
+ {
+ printf ("%10.6e |", myGH->infovals [i][1]);
+ }
}
else
{
- printf ("%10.6e |", myGH->infovals [i][1]);
+ printf(" ----------- |");
}
} /* end of loop over all variables */
@@ -322,7 +344,8 @@ int IOBasic_TriggerOutputInfo (cGH *GH, int variable)
{
const char *var;
iobasicGH *myGH;
-
+ int ierr1,ierr2;
+ CCTK_REAL val;
/* Get the GH extension for IOBasic */
myGH = (iobasicGH *) GH->extensions [CCTK_GHExtensionHandle ("IOBasic")];
@@ -336,9 +359,11 @@ int IOBasic_TriggerOutputInfo (cGH *GH, int variable)
printf(" Variable = -%s-\n",var);
#endif
- myGH->infovals[variable][0]=IOBasic_WriteInfo (GH, variable, "minimum",var);
- myGH->infovals[variable][1]=IOBasic_WriteInfo (GH, variable, "maximum",var);
-
+ ierr1=IOBasic_WriteInfo (GH, &val, variable, "minimum",var);
+ myGH->infovals[variable][0] = val;
+ ierr2=IOBasic_WriteInfo (GH, &val, variable, "maximum",var);
+ myGH->infovals[variable][1] = val;
+
/* Register GF as having Info output this iteration */
myGH->outInfo_last [variable] = GH->cctk_iteration;
diff --git a/src/Write.c b/src/Write.c
index 65be41d..1faa4a0 100644
--- a/src/Write.c
+++ b/src/Write.c
@@ -40,7 +40,9 @@ void IOBasic_Write (cGH *GH, int index, const char *alias)
/* output is done by processor 0 only */
if (CCTK_MyProc (GH) != 0)
+ {
return;
+ }
/* first, check if variable has storage assigned */
if (! CCTK_QueryGroupStorageI (GH, CCTK_GroupIndexFromVarI (index)))
@@ -49,7 +51,8 @@ void IOBasic_Write (cGH *GH, int index, const char *alias)
fullname = CCTK_FullName (index);
CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
- "No scalar output for '%s' (no storage)", fullname);
+ "IOBasic_Write: No scalar output for '%s' (no storage)",
+ fullname);
free (fullname);
return;
}
diff --git a/src/WriteGF.c b/src/WriteGF.c
index 6f02c63..aa3a1d9 100644
--- a/src/WriteGF.c
+++ b/src/WriteGF.c
@@ -28,6 +28,7 @@
void IOBasic_WriteGF (cGH *GH, int index, const char *alias)
{
DECLARE_CCTK_PARAMETERS
+ int ierr;
int i;
ioGH *ioUtilGH;
iobasicGH *myGH;
@@ -92,7 +93,8 @@ void IOBasic_WriteGF (cGH *GH, int index, const char *alias)
/* get the reduction handle from the reduction operator */
reduction_handle = CCTK_ReductionHandle (reductions [i].operator);
- if (reduction_handle < 0) {
+ if (reduction_handle < 0)
+ {
CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
"IOBasic_WriteGF: Invalid reduction operator '%s'",
reductions [i].operator);
@@ -100,11 +102,11 @@ void IOBasic_WriteGF (cGH *GH, int index, const char *alias)
}
/* do the reduction (all processors) */
- CCTK_Reduce (GH, 0, reduction_handle, 1, CCTK_VARIABLE_REAL,
- &reduction_value, 1, index);
+ ierr = CCTK_Reduce (GH, 0, reduction_handle, 1, CCTK_VARIABLE_REAL,
+ &reduction_value, 1, index);
/* dump the reduction value to file by processor 0 */
- if (CCTK_MyProc (GH) == 0)
+ if (ierr == 0 && CCTK_MyProc (GH) == 0)
{
/* build the filename */
diff --git a/src/WriteInfo.c b/src/WriteInfo.c
index ca06ae9..504b7a9 100644
--- a/src/WriteInfo.c
+++ b/src/WriteInfo.c
@@ -19,34 +19,49 @@
#include "cctk.h"
#include "CactusBase/IOUtil/src/ioGH.h"
-CCTK_REAL IOBasic_WriteInfo (cGH *GH, int index, const char *operator, const char *alias)
+int IOBasic_WriteInfo (cGH *GH,
+ CCTK_REAL *val,
+ int index,
+ const char *operator,
+ const char *alias)
{
+ int ierr;
int reduce_handle;
- CCTK_REAL retval = 0;
-
+ int retval = 0;
/* first, check if variable has storage assigned */
- if (! CCTK_QueryGroupStorageI (GH, CCTK_GroupIndexFromVarI (index))) {
+ if (! CCTK_QueryGroupStorageI (GH, CCTK_GroupIndexFromVarI (index)))
+ {
char *fullname;
fullname = CCTK_FullName (index);
CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
- "No info output for '%s' (no storage)", fullname);
+ "IOBasic_WriteInfo: No info output for '%s' (no storage)",
+ fullname);
free (fullname);
- return (0);
+ return -3;
}
reduce_handle = CCTK_ReductionHandle(operator);
if (reduce_handle > -1)
{
- CCTK_Reduce(GH,0,reduce_handle,1,CCTK_VARIABLE_REAL,&retval,1,index);
+ ierr = CCTK_Reduce(GH,0,reduce_handle,1,
+ CCTK_VARIABLE_REAL,val,1,index);
+ if (ierr < 0)
+ {
+ CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "IOBasic_WriteInfo: Internal error in reduction '%s'",
+ operator);
+ retval = -2;
+ }
}
else
{
CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
- "Reduction operator '%s' doesn't exist or failed internally", operator);
- retval = -27;
+ "IOBasic_WriteInfo: Reduction operator '%s' not found",
+ operator);
+ retval = -1;
}
return retval;
diff --git a/src/iobasicGH.h b/src/iobasicGH.h
index 1c21c61..3f41cf5 100644
--- a/src/iobasicGH.h
+++ b/src/iobasicGH.h
@@ -40,4 +40,4 @@ typedef struct IOBASICGH {
/* function prototypes */
void IOBasic_Write (cGH *GH, int index, const char *alias);
void IOBasic_WriteGF (cGH *GH, int index, const char *alias);
-CCTK_REAL IOBasic_WriteInfo (cGH *GH, int index, const char *operator, const char *alias);
+int IOBasic_WriteInfo (cGH *GH, CCTK_REAL *val,int index, const char *operator, const char *alias);