aboutsummaryrefslogtreecommitdiff
path: root/src/RadiationBoundary.c
diff options
context:
space:
mode:
authorrideout <rideout@6a38eb6e-646e-4a02-a296-d141613ad6c4>2002-08-29 08:32:52 +0000
committerrideout <rideout@6a38eb6e-646e-4a02-a296-d141613ad6c4>2002-08-29 08:32:52 +0000
commit83fa18f1a0f12565d2ac322d8b324d5a5589f086 (patch)
treef91c8f37ed09b3315a6bbf25b7660d10e200661f /src/RadiationBoundary.c
parent99a903e3f8a6370b753c35d5f177f810322f795d (diff)
Check for valid coordinate index before using it to get a data pointer.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/Boundary/trunk@181 6a38eb6e-646e-4a02-a296-d141613ad6c4
Diffstat (limited to 'src/RadiationBoundary.c')
-rw-r--r--src/RadiationBoundary.c14
1 files changed, 9 insertions, 5 deletions
diff --git a/src/RadiationBoundary.c b/src/RadiationBoundary.c
index 43c59a2..8cd220d 100644
--- a/src/RadiationBoundary.c
+++ b/src/RadiationBoundary.c
@@ -1359,7 +1359,7 @@ static int ApplyBndRadiative (const cGH *GH,
int num_vars)
{
DECLARE_CCTK_PARAMETERS
- int i, gdim;
+ int i, gdim, indx;
int var_to, var_from;
int timelvl_to, timelvl_from;
SymmetryGHex *sGHex;
@@ -1436,14 +1436,16 @@ static int ApplyBndRadiative (const cGH *GH,
sprintf (coord_system_name, "cart%dd", gdim);
for (i = 0; i < gdim; i++)
{
+
/* Radiative boundaries need the underlying Cartesian coordinates */
- coords[i] = GH->data[CCTK_CoordIndex (i + 1, NULL, coord_system_name)][0];
- if (coords[i] < 0)
+ indx = CCTK_CoordIndex (i + 1, NULL, coord_system_name);
+ if (indx < 0)
{
CCTK_VWarn(1,__LINE__,__FILE__,CCTK_THORNSTRING,
"Coordinate for system %s not found",coord_system_name);
return (-6);
}
+ coords[i] = GH->data[indx][0];
/* According to the Cactus spec, the true delta_space values for a
grid are calculated as follows: */
@@ -1455,13 +1457,15 @@ static int ApplyBndRadiative (const cGH *GH,
}
sprintf (coord_system_name, "spher%dd", gdim);
- coords[MAXDIM] = GH->data[CCTK_CoordIndex (-1, "r", coord_system_name)][0];
- if (coords[MAXDIM] < 0)
+ indx = CCTK_CoordIndex (-1, "r", coord_system_name);
+ if (indx < 0)
{
CCTK_VWarn(1,__LINE__,__FILE__,CCTK_THORNSTRING,
"Coordinate for system %s not found",coord_system_name);
return (-6);
}
+ coords[MAXDIM] = GH->data[indx][0];
+
/* see if we have a symmetry array */
sGHex = (SymmetryGHex *) CCTK_GHExtension (GH, "Symmetry");