aboutsummaryrefslogtreecommitdiff
path: root/src/WriteInfo.c
diff options
context:
space:
mode:
authortradke <tradke@b589c3ab-70e8-4b4d-a09f-cba2dd200880>2001-08-01 11:38:17 +0000
committertradke <tradke@b589c3ab-70e8-4b4d-a09f-cba2dd200880>2001-08-01 11:38:17 +0000
commit1ec7c99cb05a7a5dde7ed81fbc185be658822bda (patch)
tree40ea9282884da6a01d6cbae0eea39cc0968abb8f /src/WriteInfo.c
parentad7050ddb6009847f4b345dbb6c32e1e02e54466 (diff)
For info output, you can now specify the reduction values to print by setting
IOBasic::outInfo_reductions or add an option string to the variable name in IOBasic::outInfo_vars. Scalar variables will be output by their value. The option string syntax is IOBasic::outInfo_vars = "var1[reductions=<minimum maximum>] var2" and might still change slightly in the future. git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/IOBasic/trunk@96 b589c3ab-70e8-4b4d-a09f-cba2dd200880
Diffstat (limited to 'src/WriteInfo.c')
-rw-r--r--src/WriteInfo.c185
1 files changed, 141 insertions, 44 deletions
diff --git a/src/WriteInfo.c b/src/WriteInfo.c
index 4c22f9e..3345161 100644
--- a/src/WriteInfo.c
+++ b/src/WriteInfo.c
@@ -1,71 +1,168 @@
/*@@
- @routine WriteInfo
- @date June 31 1999
- @author Gabrielle Allen, Paul Walker
- @desc
- Dumps the Info data.
- @enddesc
- @calls
- @calledby
- @history
- @hauthor Thomas Radke @hdate 17 Mar 1999
- @hdesc included "cctk_Comm.h" to overload CCTK_MyProc() properly
- @hendhistory
+ @file WriteInfo.c
+ @date June 31 1999
+ @author Gabrielle Allen, Paul Walker, Thomas Radke
+ @desc
+ Gets the data for IOBasic's info output.
+ @enddesc
+ @version $Id: /cactusdevcvs/CactusBase/IOBasic/src/OutputInfo.c,v 1.24 2001
@@*/
-#include <stdio.h>
#include <stdlib.h>
#include "cctk.h"
+#include "iobasicGH.h"
+/* the rcs ID and its dummy function to use it */
static const char *rcsid = "$Header$";
-
CCTK_FILEVERSION(CactusBase_IOBasic_WriteInfo_c)
-int IOBasic_WriteInfo (cGH *GH,
- CCTK_REAL *val,
- int index,
- const char *operator,
- const char *alias)
+
+ /*@@
+ @routine IOBasic_WriteInfo
+ @date Tue 31 July 2001
+ @author Thomas Radke
+ @desc
+ Gets the data values to output for a given CCTK variable.
+ For CCTK_SCALAR variables, their value is taken,
+ for CCTK_GF and CCTK_ARRAY variables the appropriate
+ reduction operations are performed and their results taken.
+ @enddesc
+
+ @calls CCTK_VarDataPtrI
+ CCTK_Reduce
+@@*/
+void IOBasic_WriteInfo (cGH *GH, int vindex)
{
- int ierr;
- int reduce_handle;
- int retval = 0;
+ int vtype;
+ char *fullname;
+ iobasicGH *myGH;
+ void *ptr;
+ iobasic_reduction_t *reduction;
+
+
+ myGH = (iobasicGH *) CCTK_GHExtension (GH, "IOBasic");
+ reduction = myGH->info_reductions[vindex].reductions;
/* first, check if variable has storage assigned */
- if (! CCTK_QueryGroupStorageI (GH, CCTK_GroupIndexFromVarI (index)))
+ if (! CCTK_QueryGroupStorageI (GH, CCTK_GroupIndexFromVarI (vindex)))
{
- char *fullname;
-
- fullname = CCTK_FullName (index);
+ fullname = CCTK_FullName (vindex);
CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
- "IOBasic_WriteInfo: No info output for '%s' (no storage)",
+ "IOBasic_WriteInfo: No info output for '%s' (no storage)",
fullname);
free (fullname);
- return -3;
+ /* invalidate data buffer for this variable */
+ while (reduction)
+ {
+ reduction->is_valid = 0;
+ reduction = reduction->next;
+ }
+ return;
}
- reduce_handle = CCTK_ReductionHandle(operator);
-
- if (reduce_handle > -1)
+ /* CCTK scalars are printed by their value, cast into a C double */
+ if (CCTK_GroupTypeFromVarI (vindex) == CCTK_SCALAR)
{
- ierr = CCTK_Reduce(GH, 0, reduce_handle, 1, CCTK_VARIABLE_REAL, val,
- 1, index);
- if (ierr < 0)
+ reduction->is_valid = 1;
+
+ /* get the variable's data type and data pointer */
+ ptr = CCTK_VarDataPtrI (GH, 0, vindex);
+ vtype = CCTK_VarTypeI (vindex);
+ if (vtype == CCTK_VARIABLE_CHAR)
+ {
+ reduction->value = (double) *(CCTK_CHAR *) ptr;
+ }
+ else if (vtype == CCTK_VARIABLE_INT)
+ {
+ reduction->value = (double) *(CCTK_INT *) ptr;
+ }
+ else if (vtype == CCTK_VARIABLE_REAL)
+ {
+ reduction->value = (double) *(CCTK_REAL *) ptr;
+ }
+ else if (vtype == CCTK_VARIABLE_COMPLEX)
+ {
+ reduction->value = (double) ((CCTK_REAL *) ptr)[0];
+ reduction->next->value = (double) ((CCTK_REAL *) ptr)[1];
+ reduction->next->is_valid = 1;
+ }
+#ifdef CCTK_INT2
+ else if (vtype == CCTK_VARIABLE_INT2)
+ {
+ reduction->value = (double) *(CCTK_INT2 *) ptr;
+ }
+#endif
+#ifdef CCTK_INT4
+ else if (vtype == CCTK_VARIABLE_INT4)
{
- CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
- "IOBasic_WriteInfo: Internal error in reduction '%s'",
- operator);
- retval = -2;
+ reduction->value = (double) *(CCTK_INT4 *) ptr;
+ }
+#endif
+#ifdef CCTK_INT8
+ else if (vtype == CCTK_VARIABLE_INT8)
+ {
+ reduction->value = (double) *(CCTK_INT8 *) ptr;
+ }
+#endif
+#ifdef CCTK_REAL4
+ else if (vtype == CCTK_VARIABLE_REAL4)
+ {
+ reduction->value = (double) *(CCTK_REAL4 *) ptr;
+ }
+ else if (vtype == CCTK_VARIABLE_COMPLEX8)
+ {
+ reduction->value = (double) ((CCTK_REAL4 *) ptr)[0];
+ reduction->next->value = (double) ((CCTK_REAL4 *) ptr)[1];
+ reduction->next->is_valid = 1;
+ }
+#endif
+#ifdef CCTK_REAL8
+ else if (vtype == CCTK_VARIABLE_REAL8)
+ {
+ reduction->value = (double) *(CCTK_REAL8 *) ptr;
+ }
+ else if (vtype == CCTK_VARIABLE_COMPLEX16)
+ {
+ reduction->value = (double) ((CCTK_REAL8 *) ptr)[0];
+ reduction->next->value = (double) ((CCTK_REAL8 *) ptr)[1];
+ reduction->next->is_valid = 1;
+ }
+#endif
+#ifdef CCTK_REAL16
+ else if (vtype == CCTK_VARIABLE_REAL16)
+ {
+ reduction->value = (double) *(CCTK_REAL16 *) ptr;
+ }
+ else if (vtype == CCTK_VARIABLE_COMPLEX32)
+ {
+ reduction->value = (double) ((CCTK_REAL16 *) ptr)[0];
+ reduction->next->value = (double) ((CCTK_REAL16 *) ptr)[1];
+ reduction->next->is_valid = 1;
+ }
+#endif
+ else
+ {
+ CCTK_WARN (3, "IOBasic_WriteInfo: Unsupported data type");
+ reduction->is_valid = 0;
}
}
else
{
- CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
- "IOBasic_WriteInfo: Reduction operator '%s' not found",
- operator);
- retval = -1;
- }
+ /* for CCTK_GF and CCTK_ARRAY variables: loop over all reductions */
+ while (reduction)
+ {
+ reduction->is_valid = CCTK_Reduce (GH, 0, reduction->handle, 1,
+ CCTK_VARIABLE_REAL,
+ &reduction->value, 1, vindex) == 0;
+ if (! reduction->is_valid)
+ {
+ CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "IOBasic_WriteInfo: Internal error in reduction '%s'",
+ reduction->name);
+ }
- return retval;
+ reduction = reduction->next;
+ }
+ }
}