diff options
author | tradke <tradke@ff385933-4943-42dc-877b-ffc776028de6> | 2001-05-26 09:18:39 +0000 |
---|---|---|
committer | tradke <tradke@ff385933-4943-42dc-877b-ffc776028de6> | 2001-05-26 09:18:39 +0000 |
commit | 3f44a94319cdd9ee4901a3051c80b5817fba4811 (patch) | |
tree | 3f19e8f93d2a9d392a57fb2d7fe66a2cd3a665dc | |
parent | bf1a0a985be510aec7cf25b12f18aa83cc4e89ba (diff) |
Changed the name of the 'NaNChecker::check_max' parameter into
'NaNChecker::report_max' since this sound more appropriate.
Also changed its default from 0 to -1 to allow disabling any level 2 warnings.
Added output of physical coordinate location for grid functions in level 2
warnings if coordinates are available.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusUtils/NaNChecker/trunk@11 ff385933-4943-42dc-877b-ffc776028de6
-rw-r--r-- | README | 2 | ||||
-rw-r--r-- | param.ccl | 8 | ||||
-rw-r--r-- | src/NaNCheck.c | 98 |
3 files changed, 79 insertions, 29 deletions
@@ -25,7 +25,7 @@ Currently these actions can be to and how many (for grid array variables). For grid arrays it will also print level 2 warnings with the array index (in fortran order) for all NaN elements. You can limit the number of - such warnings by setting the 'NanChecker::check_max' parameter. + such warnings by setting the 'NanChecker::report_max' parameter. * also set the CCTK termination flag so that Cactus will stop the evolution loop and gracefully terminate at the next time possible (giving you the @@ -12,11 +12,11 @@ INT check_every "How often to check for NaN's" STEERABLE = ALWAYS 1:* :: "Every so many iterations" } 0 -INT check_max "How many NaN's to report for a single variable" STEERABLE = ALWAYS +INT report_max "How many NaN's to report for a single variable" STEERABLE = ALWAYS { - 0 :: "Report all (default)" - 1:* :: "Do not report more than check_max number of NaN's" -} 0 + -1 :: "Report all (default)" + 0:* :: "Do not report more than report_max number of NaN's" +} -1 STRING check_vars "Groups and/or variables to check for NaN's" STEERABLE = ALWAYS { diff --git a/src/NaNCheck.c b/src/NaNCheck.c index 79235cd..a688671 100644 --- a/src/NaNCheck.c +++ b/src/NaNCheck.c @@ -36,6 +36,7 @@ static void NaNCheck (int vindex, const char *optstring, void *arg); static void PrintWarning (const char *error_type, int linear_index, int fp_type, + const CCTK_REAL *const coords[], const char *fullname, const cGroupDynamicData *gdata); @@ -95,6 +96,8 @@ int NaNChecker_NaNCheck (cGH *GH) at the given processor-local linear index. The warning includes the variable's fullname along with the global index of the NaN element in fortran order. + If coordinates are available, the NaN's location on the grid + is also output. @enddesc @calls CCTK_VWarn @@ -103,13 +106,13 @@ int NaNChecker_NaNCheck (cGH *GH) @vtype const char * @vio in @endvar - @var fp_type - @vdesc indicates if variable of of real or complex type + @var linear_index + @vdesc processor-local linear index of the NaN/Inf @vtype int @vio in @endvar - @var linear_index - @vdesc processor-local linear index of the NaN/Inf + @var fp_type + @vdesc indicates if variable of of real or complex type @vtype int @vio in @endvar @@ -118,6 +121,11 @@ int NaNChecker_NaNCheck (cGH *GH) @vtype const char * @vio in @endvar + @var coords + @vdesc array of coordinates + @vtype const CCTK_REAL *const [ dimension of variable ] + @vio in + @endvar @var gdata @vdesc size information on the variable to compute the global index @vtype const cGroupDynamicData * @@ -127,11 +135,12 @@ int NaNChecker_NaNCheck (cGH *GH) static void PrintWarning (const char *error_type, int linear_index, int fp_type, + const CCTK_REAL *const coords[], const char *fullname, const cGroupDynamicData *gdata) { int i; - char *index_string; + char *index_string, *coord_string; const char *complex_part; @@ -153,21 +162,42 @@ static void PrintWarning (const char *error_type, } else { - /* assume max. 10 characters per index number (including separators) */ - index_string = (char *) malloc (10 * gdata->dim); + /* assume max. 10 characters per index number and 40 characters per + coordinate value (including separators) */ + index_string = (char *) malloc (5 * 10 * gdata->dim); + coord_string = index_string + 10 * gdata->dim; sprintf (index_string, "%d", (linear_index % gdata->lsh[0]) + gdata->lbnd[0] + 1); + if (coords) + { + sprintf (coord_string, "%5.3e", (double) coords[0][linear_index]); + } for (i = 1; i < gdata->dim; i++) { linear_index /= gdata->lsh[i - 1]; sprintf (index_string, "%s, %d", index_string, (linear_index % gdata->lsh[i]) + gdata->lbnd[i] + 1); + if (coords) + { + sprintf (coord_string, "%s, %5.3e", coord_string, + (double) coords[i][linear_index]); + } } - CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING, - "%s caught in %svariable '%s' at (%s)", - error_type, complex_part, fullname, index_string); + if (coords) + { + CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING, + "%s caught in %svariable '%s' at index (%s) with coordinates " + "(%s)", error_type, complex_part, fullname, index_string, + coord_string); + } + else + { + CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING, + "%s caught in %svariable '%s' at (%s)", + error_type, complex_part, fullname, index_string); + } free (index_string); } @@ -205,10 +235,10 @@ static void PrintWarning (const char *error_type, if (! finite ((double) _typed_data[_i])) \ { \ nans_found++; \ - if (check_max == 0 || nans_found <= check_max) \ + if (report_max < 0 || nans_found <= report_max) \ { \ PrintWarning (isnan ((double) _typed_data[_i]) ? "NaN" : "Inf", \ - _i, fp_type, fullname, &gdata); \ + _i, fp_type, coords, fullname, &gdata); \ } \ } \ } \ @@ -216,8 +246,6 @@ static void PrintWarning (const char *error_type, #else -#ifdef HAVE_ISNAN - #define CHECK_DATA(cctk_type) \ { \ int _i; \ @@ -230,20 +258,14 @@ static void PrintWarning (const char *error_type, if (isnan ((double) _typed_data[_i])) \ { \ nans_found++; \ - if (check_max == 0 || nans_found <= check_max) \ + if (report_max < 0 || nans_found <= report_max) \ { \ - PrintWarning ("NaN", _i, fp_type, fullname, &gdata); \ + PrintWarning ("NaN", _i, fp_type, coords, fullname, &gdata); \ } \ } \ } \ } -#else - -#error Unable to check for NaNs on this architecture yet - -#endif /* HAVE_ISNAN */ - #endif /* HAVE_FINITE */ @@ -282,9 +304,11 @@ static void NaNCheck (int vindex, const char *optstring, void *_GH) DECLARE_CCTK_PARAMETERS const cGH *GH; int i, fp_type, nans_found; - int vtype, gindex, nelems; + int vtype, gtype, gindex, nelems; char *fullname; const char *vtypename; + char coord_system_name[10]; + const CCTK_REAL **coords; cGroupDynamicData gdata; const void *data; @@ -322,12 +346,33 @@ static void NaNCheck (int vindex, const char *optstring, void *_GH) /* get the number of elements to check for this variable */ nelems = 1; gdata.dim = 0; - if (CCTK_GroupTypeI (vindex) != CCTK_SCALAR) + coords = NULL; + gtype = CCTK_GroupTypeI (gindex); + if (gtype != CCTK_SCALAR) { CCTK_GroupDynamicData (GH, gindex, &gdata); + if (gtype == CCTK_GF) + { + sprintf (coord_system_name, "cart%dd", gdata.dim); + if (CCTK_CoordSystemHandle (coord_system_name) >= 0) + { + coords = (const CCTK_REAL **) + malloc (gdata.dim * sizeof (CCTK_REAL *)); + } + } for (i = 0; i < gdata.dim; i++) { nelems *= gdata.lsh[i]; + if (coords) + { + coords[i] = (const CCTK_REAL *) CCTK_VarDataPtrI (GH, 0, + CCTK_CoordIndex (i + 1, NULL, coord_system_name)); + if (! coords[i]) + { + free (coords); + coords = NULL; + } + } } } @@ -390,6 +435,11 @@ static void NaNCheck (int vindex, const char *optstring, void *_GH) CCTK_Abort (NULL, 0); } } + + if (coords) + { + free (coords); + } } else { |