aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@082bdb00-0f4f-0410-b49e-b1835e5f2039>2009-02-10 22:30:11 +0000
committerschnetter <schnetter@082bdb00-0f4f-0410-b49e-b1835e5f2039>2009-02-10 22:30:11 +0000
commit6e0b30c8e62c310435cd937906dd4abc7d2b91ad (patch)
tree1fae8cb6f1800b3830b8d0fa1aab79c6ecd427f1
parent6df1f900b3683003747ec42a4b4e1b0c7904f655 (diff)
Split COPY macro into COPY_PRE, COPY_LOOP, and COPY_POST to allow
adding OpenMP #pragma statements. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/ReflectionSymmetry/trunk@25 082bdb00-0f4f-0410-b49e-b1835e5f2039
-rw-r--r--src/apply.c73
1 files changed, 53 insertions, 20 deletions
diff --git a/src/apply.c b/src/apply.c
index d5223da..9375b56 100644
--- a/src/apply.c
+++ b/src/apply.c
@@ -18,7 +18,7 @@ CCTK_FILEVERSION(AEIDevelopment_ReflectionSymmetry_apply_c);
-#define COPY(VARTYPE) \
+#define COPY_PRE(VARTYPE) \
static void \
copy_##VARTYPE (VARTYPE const * restrict const srcvar, \
VARTYPE * restrict const dstvar, \
@@ -29,8 +29,6 @@ CCTK_FILEVERSION(AEIDevelopment_ReflectionSymmetry_apply_c);
int const idir, int const jdir, int const kdir, \
int const parity) \
{ \
- int i, j, k; \
- \
int const iioff = ioff + (1 - idir) * imin; \
int const jjoff = joff + (1 - jdir) * jmin; \
int const kkoff = koff + (1 - kdir) * kmin; \
@@ -48,11 +46,12 @@ CCTK_FILEVERSION(AEIDevelopment_ReflectionSymmetry_apply_c);
assert(kkmin>=0 && kkmax<=nk); \
assert(iimax>=-1 && iimin<ni); \
assert(jjmax>=-1 && jjmin<nj); \
- assert(kkmax>=-1 && kkmin<nk); \
- \
- for (k=kmin; k<kmax; ++k) { \
- for (j=jmin; j<jmax; ++j) { \
- for (i=imin; i<imax; ++i) { \
+ assert(kkmax>=-1 && kkmin<nk);
+
+#define COPY_LOOP(VARTYPE) \
+ for (int k=kmin; k<kmax; ++k) { \
+ for (int j=jmin; j<jmax; ++j) { \
+ for (int i=imin; i<imax; ++i) { \
int const dstind = i + ni * (j + nj * k); \
int const ii = iioff + idir * i; \
int const jj = jjoff + jdir * j; \
@@ -63,6 +62,8 @@ CCTK_FILEVERSION(AEIDevelopment_ReflectionSymmetry_apply_c);
} \
} \
} \
+
+#define COPY_POST(VARTYPE) \
}
#define REAL(x) x
@@ -71,31 +72,52 @@ CCTK_FILEVERSION(AEIDevelopment_ReflectionSymmetry_apply_c);
#define IM /* nothing */
#ifdef HAVE_CCTK_INT1
-COPY(CCTK_INT1)
+COPY_PRE(CCTK_INT1)
+#pragma omp parallel for
+COPY_LOOP(CCTK_INT1)
+COPY_POST(CCTK_INT1)
#endif
#ifdef HAVE_CCTK_INT2
-COPY(CCTK_INT2)
+COPY_PRE(CCTK_INT2)
+#pragma omp parallel for
+COPY_LOOP(CCTK_INT2)
+COPY_POST(CCTK_INT2)
#endif
#ifdef HAVE_CCTK_INT4
-COPY(CCTK_INT4)
+COPY_PRE(CCTK_INT4)
+#pragma omp parallel for
+COPY_LOOP(CCTK_INT4)
+COPY_POST(CCTK_INT4)
#endif
#ifdef HAVE_CCTK_INT8
-COPY(CCTK_INT8)
+COPY_PRE(CCTK_INT8)
+#pragma omp parallel for
+COPY_LOOP(CCTK_INT8)
+COPY_POST(CCTK_INT8)
#endif
#ifdef HAVE_CCTK_REAL4
-COPY(CCTK_REAL4)
+COPY_PRE(CCTK_REAL4)
+#pragma omp parallel for
+COPY_LOOP(CCTK_REAL4)
+COPY_POST(CCTK_REAL4)
#endif
#ifdef HAVE_CCTK_REAL8
-COPY(CCTK_REAL8)
+COPY_PRE(CCTK_REAL8)
+#pragma omp parallel for
+COPY_LOOP(CCTK_REAL8)
+COPY_POST(CCTK_REAL8)
#endif
#ifdef HAVE_CCTK_REAL16
-COPY(CCTK_REAL16)
+COPY_PRE(CCTK_REAL16)
+#pragma omp parallel for
+COPY_LOOP(CCTK_REAL16)
+COPY_POST(CCTK_REAL16)
#endif
#undef REAL
@@ -109,15 +131,24 @@ COPY(CCTK_REAL16)
#define IM .Im
#ifdef HAVE_CCTK_COMPLEX8
-COPY(CCTK_COMPLEX8)
+COPY_PRE(CCTK_COMPLEX8)
+#pragma omp parallel for
+COPY_LOOP(CCTK_COMPLEX8)
+COPY_POST(CCTK_COMPLEX8)
#endif
#ifdef HAVE_CCTK_COMPLEX16
-COPY(CCTK_COMPLEX16)
+COPY_PRE(CCTK_COMPLEX16)
+#pragma omp parallel for
+COPY_LOOP(CCTK_COMPLEX16)
+COPY_POST(CCTK_COMPLEX16)
#endif
#ifdef HAVE_CCTK_COMPLEX32
-COPY(CCTK_COMPLEX32)
+COPY_PRE(CCTK_COMPLEX32)
+#pragma omp parallel for
+COPY_LOOP(CCTK_COMPLEX32)
+COPY_POST(CCTK_COMPLEX32)
#endif
#undef REAL
@@ -125,7 +156,9 @@ COPY(CCTK_COMPLEX32)
#undef RE
#undef IM
-#undef COPY
+#undef COPY_PRE
+#undef COPY_LOOP
+#undef COPY_POST
@@ -647,7 +680,7 @@ CheckBoundaryParameters (cGH const * restrict const cctkGH,
ptr = CCTK_ParameterGet ("type", "CartGrid3D", & type);
assert (ptr != 0);
assert (type == PARAMETER_KEYWORD);
- coordtype = * (char const * *) ptr;
+ coordtype = * (char const * const *) ptr;
if (! CCTK_EQUALS (coordtype, "coordbase")) return;
/* Get the boundary specification */