aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authortradke <tradke@b32723a9-ab3a-4a60-88e2-2e5d99d7c17a>2004-05-09 15:50:27 +0000
committertradke <tradke@b32723a9-ab3a-4a60-88e2-2e5d99d7c17a>2004-05-09 15:50:27 +0000
commit43ab51690c928f935a44da4f6460a55fe4fc3bb1 (patch)
treef439649d854c6162150baf019dc7def53bb60d08
parentdd0c56d782f2442d7be5b2195cce0395660c8fb5 (diff)
Call aliased function Coord_GetDefaultSystem() to get the default coordinate
system for setting up 1D line and 2D slice centers. This closes PR IO/1681: "IOUtil still uses the old CCTK_Coord functions rather than the CoordBase functions". git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/IOUtil/trunk@204 b32723a9-ab3a-4a60-88e2-2e5d99d7c17a
-rw-r--r--interface.ccl10
-rw-r--r--src/Utils.c204
2 files changed, 121 insertions, 93 deletions
diff --git a/interface.ccl b/interface.ccl
index b2f16ad..40f0d99 100644
--- a/interface.ccl
+++ b/interface.ccl
@@ -1,3 +1,13 @@
# Interface definitions for thorn IOUtil
+# $Header$
implements: IO
+
+
+# aliased functions required from Coordinate base thorn
+
+CCTK_INT FUNCTION Coord_GetDefaultSystem \
+ (CCTK_POINTER_TO_CONST IN GH, \
+ CCTK_INT IN systemdim)
+
+REQUIRES FUNCTION Coord_GetDefaultSystem
diff --git a/src/Utils.c b/src/Utils.c
index d2bae38..af74e7a 100644
--- a/src/Utils.c
+++ b/src/Utils.c
@@ -260,8 +260,7 @@ void IOUtil_ParseOutputFrequency (const char *method_name,
of that line on the grid.
@enddesc
- @calls CCTK_CoordSystemHandle
- CCTK_CoordRange
+ @calls Coord_GetDefaultSystem
@var GH
@vdesc Pointer to CCTK GH
@@ -301,64 +300,74 @@ int IOUtil_1DLines (const cGH *GH,
CCTK_REAL *const *const origin_phys,
int *const *slice_center)
{
- int dim, dir, handle, no_origin_index;
- char coord_system_name[20];
- CCTK_REAL *lower_range, *upper_range;
-
-
- /* get the appropriate coordinate system name */
- sprintf (coord_system_name, "cart%dd", num_dims);
- handle = CCTK_CoordSystemHandle (coord_system_name);
+ int dim, dir, coord_system_handle, have_coords, have_origin_index, len;
+ char *coord_system_name;
+ CCTK_INT *coord_handles;
+ CCTK_REAL *lower, *upper;
+
+
+ /* get the default coordinate system associated with this grid dimension */
+ coord_system_name = NULL;
+ coord_system_handle = Coord_GetDefaultSystem (GH, num_dims);
+ coord_handles = malloc (num_dims * sizeof (CCTK_INT));
+ have_coords = coord_system_handle >= 0 &&
+ Util_TableGetIntArray (coord_system_handle, num_dims,
+ coord_handles, "COORDINATES") >= 0;
+ if (have_coords)
+ {
+ len = Util_TableGetString (coord_system_handle, 0, NULL, "NAME");
+ have_coords = len > 0;
+ if (have_coords)
+ {
+ coord_system_name = malloc (len);
+ Util_TableGetString (coord_system_handle, len, coord_system_name, "NAME");
+ have_coords = Util_StrCmpi (coord_system_name, "cart");
+ }
+ }
- /* if no coordinate system is available then origin_index[] must be given */
- no_origin_index = origin_index == NULL;
- if (handle < 0 && origin_index)
+ /* check that origin_index[] is valid */
+ have_origin_index = origin_index != NULL;
+ for (dir = 0; dir < num_dims && have_origin_index; dir++)
{
- /* check that origin_index[] is valid */
- for (dir = 0; dir < num_dims; dir++)
+ for (dim = 0; dim < num_dims && have_origin_index; dim++)
{
- for (dim = 0; dim < num_dims; dim++)
- {
- if (origin_index[dir][dim] < 0)
- {
- no_origin_index = 1;
- break;
- }
- }
+ have_origin_index = origin_index[dir][dim] >= 0;
}
}
- if (handle < 0 && no_origin_index)
+ /* if no coordinate info is available then origin_index[] must be given */
+ if (! have_coords && ! have_origin_index)
{
CCTK_VWarn (3, __LINE__, __FILE__, CCTK_THORNSTRING,
- "IOUtil_1DLines: Cartesian coordinate system '%s' not found "
- "and no slice center index coordinates set - slice center will "
- "not be set up for output of 1D lines from %dD variables",
- coord_system_name, num_dims);
+ "IOUtil_1DLines: Found no default Cartesian coordinate system "
+ "associated with grid variables of dimension %d, and no slice "
+ "center index coordinates were given either - slice center "
+ "will not be set up for output of 1D lines from %dD variables",
+ num_dims, num_dims);
for (dir = 0; dir < num_dims; dir++)
{
memset (slice_center[dir], 0, num_dims * sizeof (int));
}
+ free (coord_handles);
+ free (coord_system_name);
+
return (-1);
}
- lower_range = upper_range = NULL;
- if (handle >= 0)
+ /* get the ranges in every direction */
+ lower = calloc (2 * num_dims, sizeof (CCTK_REAL));
+ upper = lower + num_dims;
+ if (! have_origin_index)
{
- /* allocate arrays for ranges */
- lower_range = calloc (2 * num_dims, sizeof (CCTK_REAL));
- upper_range = lower_range + num_dims;
-
- /* get the ranges in every direction */
for (dir = 0; dir < num_dims; dir++)
{
- if (CCTK_CoordRange (GH, &lower_range[dir], &upper_range[dir],
- dir + 1, NULL, coord_system_name) < 0)
+ if (Util_TableGetReal (coord_handles[dir], &lower[dir], "COMPMIN") < 0 ||
+ Util_TableGetReal (coord_handles[dir], &upper[dir], "COMPMAX") < 0)
{
- CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
+ CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
"IOUtil_1DLines: Could not get ranges for %c-direction "
- "of coordinate system '%s'",
+ "of associated default coordinate system '%s'",
'x' + dir, coord_system_name);
}
}
@@ -375,13 +384,13 @@ int IOUtil_1DLines (const cGH *GH,
/* line always starts at the first point */
slice_center[dir][dim] = 0;
}
- else if (origin_index && origin_index[dir][dim] >= 0)
+ else if (have_origin_index)
{
/* FIXME: check upper index bounds also ?? */
slice_center[dir][dim] = origin_index[dir][dim];
}
- else if (lower_range[dim] > origin_phys[dir][dim] ||
- upper_range[dim] < origin_phys[dir][dim])
+ else if (lower[dim] > origin_phys[dir][dim] ||
+ upper[dim] < origin_phys[dir][dim])
{
CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
"IOUtil_1DLines: %c-coordinate (%f) for slice center of "
@@ -389,7 +398,7 @@ int IOUtil_1DLines (const cGH *GH,
"grid coordinates range (%f, %f)",
'x' + dim, (double) origin_phys[dir][dim],
'x' + dir, num_dims,
- (double) lower_range[dim], (double) upper_range[dim]);
+ (double) lower[dim], (double) upper[dim]);
CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
"IOUtil_1DLines: slice center will default to %c-index 0",
'x' + dir);
@@ -399,7 +408,7 @@ int IOUtil_1DLines (const cGH *GH,
{
/* Find index for first point above the chosen coordinate */
slice_center[dir][dim] =
- ceil ((origin_phys[dir][dim] - lower_range[dim]) /
+ ceil ((origin_phys[dir][dim] - lower[dim]) /
GH->cctk_delta_space[dim] - 1e-6);
#ifdef DEBUG_IOUTIL
@@ -410,11 +419,9 @@ int IOUtil_1DLines (const cGH *GH,
}
}
- /* free allocated resources */
- if (lower_range)
- {
- free (lower_range);
- }
+ free (lower);
+ free (coord_handles);
+ free (coord_system_name);
return (0);
}
@@ -429,7 +436,7 @@ int IOUtil_1DLines (const cGH *GH,
for output on a multidimensional regular Cartesian unigrid.
@enddesc
- @calls CCTK_CoordRange
+ @calls Coord_GetDefaultSystem
@var GH
@vdesc Pointer to CCTK GH
@@ -469,55 +476,68 @@ int IOUtil_2DPlanes (const cGH *GH,
const CCTK_REAL *origin_phys,
int *slice_center)
{
- int handle, dir;
- char coord_system_name[20];
- CCTK_REAL *lower_range, *upper_range;
-
-
- /* get the appropriate coordinate system name */
- sprintf (coord_system_name, "cart%dd", num_dims);
- handle = CCTK_CoordSystemHandle (coord_system_name);
+ int have_coords, len, coord_system_handle, dir;
+ char *coord_system_name;
+ CCTK_INT *coord_handles;
+ CCTK_REAL *lower, *upper;
+
+
+ /* get the default coordinate system associated with this grid dimension */
+ coord_system_name = NULL;
+ coord_system_handle = Coord_GetDefaultSystem (GH, num_dims);
+ coord_handles = malloc (num_dims * sizeof (CCTK_INT));
+ have_coords = coord_system_handle >= 0 &&
+ Util_TableGetIntArray (coord_system_handle, num_dims,
+ coord_handles, "COORDINATES") >= 0;
+ if (have_coords)
+ {
+ len = Util_TableGetString (coord_system_handle, 0, NULL, "NAME");
+ have_coords = len > 0;
+ if (have_coords)
+ {
+ coord_system_name = malloc (len);
+ Util_TableGetString (coord_system_handle, len, coord_system_name, "NAME");
+ have_coords = Util_StrCmpi (coord_system_name, "cart");
+ }
+ }
- /* if no coordinate system is available then origin_index[] must be given */
- if (handle < 0 && origin_index)
+ for (dir = 0; dir < num_dims && origin_index; dir++)
{
- /* check that origin_index[] is valid */
- for (dir = 0; dir < num_dims; dir++)
+ if (origin_index[dir] < 0)
{
- if (origin_index[dir] < 0)
- {
- origin_index = NULL;
- break;
- }
+ origin_index = NULL;
}
}
- if (handle < 0 && ! origin_index)
+ /* if no coordinate system is available then origin_index[] must be given */
+ if (! have_coords && ! origin_index)
{
CCTK_VWarn (3, __LINE__, __FILE__, CCTK_THORNSTRING,
- "IOUtil_2DPlanes: Cartesian coordinate system '%s' not found "
- "and no slice center index coordinates set - slice center will "
- "not be set up for output of 2D planes from %dD variables",
- coord_system_name, num_dims);
+ "IOUtil_2DPlanes: Found no default Cartesian coordinate system "
+ "associated with grid variables of dimension %d, and no slice "
+ "center index coordinates were given either - slice center "
+ "will not be set up for output of 2D planes from %dD variables",
+ num_dims, num_dims);
+
+ free (coord_handles);
+ free (coord_system_name);
+
return (-1);
}
- lower_range = upper_range = NULL;
- if (handle >= 0)
+ /* get the ranges in every direction */
+ lower = calloc (2 * num_dims, sizeof (CCTK_REAL));
+ upper = lower + num_dims;
+ if (! origin_index)
{
- /* allocate arrays for ranges */
- lower_range = calloc (2 * num_dims, sizeof (CCTK_REAL));
- upper_range = lower_range + num_dims;
-
- /* get the ranges in every direction */
for (dir = 0; dir < num_dims; dir++)
{
- if (CCTK_CoordRange (GH, &lower_range[dir], &upper_range[dir],
- num_dims - dir, NULL, coord_system_name) < 0)
+ if (Util_TableGetReal (coord_handles[dir], &lower[dir], "COMPMIN") < 0 ||
+ Util_TableGetReal (coord_handles[dir], &upper[dir], "COMPMAX") < 0)
{
- CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
+ CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
"IOUtil_2DPlanes: Could not get ranges for %c-direction "
- "of coordinate system '%s'",
+ "of associated coordinate system '%s'",
'x' + (num_dims - dir - 1), coord_system_name);
}
}
@@ -527,18 +547,18 @@ int IOUtil_2DPlanes (const cGH *GH,
according to origin_index[] or origin_phys[] */
for (dir = 0; dir < num_dims; dir++)
{
- if (origin_index && origin_index[dir] >= 0)
+ if (origin_index)
{
slice_center[dir] = origin_index[dir];
}
- else if (lower_range[dir] > origin_phys[dir] ||
- upper_range[dir] < origin_phys[dir])
+ else if (lower[dir] > origin_phys[dir] ||
+ upper[dir] < origin_phys[dir])
{
CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
"IOUtil_2DPlanes: %c-coordinate for slice center of 2D "
"planes (%f) is out of grid coordinates range (%f, %f)",
'x' + (num_dims - dir - 1), (double) origin_phys[dir],
- (double) lower_range[dir], (double) upper_range[dir]);
+ (double) lower[dir], (double) upper[dir]);
CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
"IOUtil_2DPlanes: slice center will default to %c-index 0",
'x' + (num_dims - dir - 1));
@@ -547,7 +567,7 @@ int IOUtil_2DPlanes (const cGH *GH,
else
{
/* Find index for first point above the chosen coordinate */
- slice_center[dir] = ceil ((origin_phys[dir] - lower_range[dir]) /
+ slice_center[dir] = ceil ((origin_phys[dir] - lower[dir]) /
GH->cctk_delta_space[(num_dims - dir - 1)] -
1e-6);
@@ -558,11 +578,9 @@ int IOUtil_2DPlanes (const cGH *GH,
}
}
- /* free allocated resources */
- if (lower_range)
- {
- free (lower_range);
- }
+ free (lower);
+ free (coord_handles);
+ free (coord_system_name);
return (0);
}