aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authortradke <tradke@6a38eb6e-646e-4a02-a296-d141613ad6c4>2003-10-06 15:50:03 +0000
committertradke <tradke@6a38eb6e-646e-4a02-a296-d141613ad6c4>2003-10-06 15:50:03 +0000
commit6d83ae21de4517c303319c91cfe67e79bc331e00 (patch)
tree84868c5da001ba1f92c5cb4507be99be578c85db
parent3268a17089b1c55495982f559d03eb5f36d2a0ea (diff)
Some optimizations: declare coordinates to be constant arrays (during the
boundary calculations), keep some pointers and indices in local variables. At least on the SX-5 this gives 50% speedup. git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/Boundary/trunk@248 6a38eb6e-646e-4a02-a296-d141613ad6c4
-rw-r--r--src/RadiationBoundary.c66
1 files changed, 34 insertions, 32 deletions
diff --git a/src/RadiationBoundary.c b/src/RadiationBoundary.c
index a4b15e2..e025c9e 100644
--- a/src/RadiationBoundary.c
+++ b/src/RadiationBoundary.c
@@ -1170,27 +1170,27 @@ void CCTK_FCALL CCTK_FNAME (BndRadiativeVN)
#define LOWER_RADIATIVE_BOUNDARY_3D(istart, jstart, kstart, dim, cctk_type) \
{ \
int _i, _j, _k; \
+ int _0 = 0*offset[dim], \
+ _1 = 1*offset[dim], \
+ _2 = 2*offset[dim]; \
\
\
for (_k = kstart-1; _k >= 0; _k--) \
{ \
for (_j = jstart-1; _j >= 0; _j--) \
{ \
- int _0 = 0*offset[dim], \
- _1 = 1*offset[dim], \
- _2 = 2*offset[dim]; \
int _idx = CCTK_GFINDEX3D (GH, istart-1, _j, _k); \
- CCTK_REAL *_r = coords[MAXDIM] + _idx, \
- *_xyz = coords[dim] + _idx; \
- cctk_type *_to = (cctk_type *) GH->data[var_to][timelvl_to] + _idx, \
- *_from = (cctk_type *) GH->data[var_from][timelvl_from] + _idx;\
+ const CCTK_REAL *_r = xyzr[MAXDIM] + _idx, \
+ *_xyz = xyzr[dim] + _idx; \
+ cctk_type *_to = (cctk_type *) to_ptr + _idx; \
+ const cctk_type *_from = (const cctk_type *) from_ptr + _idx; \
\
\
for (_i = istart-1; _i >= 0; _i--) \
{ \
CCTK_REAL _r0_inv = 1/_r[_0], _r1_inv = 1/_r[_1]; \
\
- if (radpower >= 0) \
+ if (radpower > 0) \
{ \
CCTK_REAL H; \
\
@@ -1253,27 +1253,27 @@ void CCTK_FCALL CCTK_FNAME (BndRadiativeVN)
#define UPPER_RADIATIVE_BOUNDARY_3D(istart, jstart, kstart, dim, cctk_type) \
{ \
int _i, _j, _k; \
+ int _0 = -0*offset[dim], \
+ _1 = -1*offset[dim], \
+ _2 = -2*offset[dim]; \
\
\
for (_k = kstart; _k < GH->cctk_lsh[2]; _k++) \
{ \
for (_j = jstart; _j < GH->cctk_lsh[1]; _j++) \
{ \
- int _0 = -0*offset[dim], \
- _1 = -1*offset[dim], \
- _2 = -2*offset[dim]; \
int _idx = CCTK_GFINDEX3D (GH, istart, _j, _k); \
- CCTK_REAL *_r = coords[MAXDIM] + _idx, \
- *_xyz = coords[dim] + _idx; \
- cctk_type *_to = (cctk_type *) GH->data[var_to][timelvl_to] + _idx, \
- *_from = (cctk_type *) GH->data[var_from][timelvl_from] + _idx;\
+ const CCTK_REAL *_r = xyzr[MAXDIM] + _idx, \
+ *_xyz = xyzr[dim] + _idx; \
+ cctk_type *_to = (cctk_type *) to_ptr + _idx; \
+ const cctk_type *_from = (const cctk_type *) from_ptr + _idx; \
\
\
for (_i = istart; _i < GH->cctk_lsh[0]; _i++) \
{ \
CCTK_REAL _r0_inv = 1/_r[_0], _r1_inv = 1/_r[_1]; \
\
- if (radpower >= 0) \
+ if (radpower > 0) \
{ \
CCTK_REAL H; \
\
@@ -1475,15 +1475,19 @@ static int ApplyBndRadiative (const cGH *GH,
CCTK_INT first_var_from,
int num_vars)
{
- DECLARE_CCTK_PARAMETERS;
int i, gdim, indx;
int var_to, var_from;
- int timelvl_to, timelvl_from;
+ int timelvl_from;
SymmetryGHex *sGHex;
char coord_system_name[10];
- CCTK_REAL dxyz[MAXDIM], rho[MAXDIM], *coords[MAXDIM+1];
+ CCTK_REAL dxyz[MAXDIM], rho[MAXDIM];
+ const CCTK_REAL *xyzr[MAXDIM+1];
int doBC[2*MAXDIM], widths[2*MAXDIM], offset[MAXDIM];
CCTK_REAL dtv, dtvh, dtvvar0, dtvvar0H;
+ void *to_ptr;
+ const void *from_ptr;
+ DECLARE_CCTK_PARAMETERS
+
/* check the direction parameter */
if (abs (dir) > MAXDIM)
@@ -1524,15 +1528,11 @@ static int ApplyBndRadiative (const cGH *GH,
}
/* Use next time level, if available */
- timelvl_to = 0;
+ timelvl_from = 0;
if (CCTK_MaxTimeLevelsVI(first_var_from) > 1)
{
timelvl_from = 1;
}
- else
- {
- timelvl_from = 0;
- }
/* Find Courant parameters. */
dtv = speed * GH->cctk_delta_time;
@@ -1543,7 +1543,6 @@ 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 */
indx = CCTK_CoordIndex (i + 1, NULL, coord_system_name);
if (indx < 0)
@@ -1552,7 +1551,7 @@ static int ApplyBndRadiative (const cGH *GH,
"Coordinate for system %s not found",coord_system_name);
return (-6);
}
- coords[i] = GH->data[indx][0];
+ xyzr[i] = GH->data[indx][0];
/* According to the Cactus spec, the true delta_space values for a
grid are calculated as follows: */
@@ -1563,26 +1562,29 @@ static int ApplyBndRadiative (const cGH *GH,
offset[i] = i == 0 ? 1 : offset[i-1] * GH->cctk_lsh[i-1];
}
- /* Append r grid variable to end of coords[] array */
+ /* Append r grid variable to end of xyzr[] array */
sprintf (coord_system_name, "spher%dd", gdim);
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);
+ CCTK_VWarn(1,__LINE__,__FILE__,CCTK_THORNSTRING,
+ "Coordinate for system %s not found",coord_system_name);
+ return (-6);
}
- coords[MAXDIM] = GH->data[indx][0];
+ xyzr[MAXDIM] = GH->data[indx][0];
/* see if we have a symmetry array */
- sGHex = (SymmetryGHex *) CCTK_GHExtension (GH, "Symmetry");
+ sGHex = CCTK_GHExtension (GH, "Symmetry");
/* now loop over all variables */
for (var_to = first_var_to, var_from = first_var_from;
var_to < first_var_to + num_vars;
var_to++, var_from++)
{
+ to_ptr = GH->data[var_to][0];
+ from_ptr = GH->data[var_from][timelvl_from];
+
/* Apply condition if:
+ boundary is not a symmetry boundary
(no symmetry or unset(=unsed))