aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8>2003-05-23 23:47:00 +0000
committerschnetter <schnetter@2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8>2003-05-23 23:47:00 +0000
commitea02c23168e25979a327f2bb1401ea0af68dc97d (patch)
tree700f24845cb96faab9caf4b700b4cd5221e77a24
parent20b08143e9ba60a2a15d4def1c4b3804b5b59ff8 (diff)
Fix optimisations that didn't work on multiple processors.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/Slab/trunk@23 2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8
-rw-r--r--src/slab.c73
1 files changed, 53 insertions, 20 deletions
diff --git a/src/slab.c b/src/slab.c
index 7b3fed4..f68df9b 100644
--- a/src/slab.c
+++ b/src/slab.c
@@ -12,7 +12,21 @@
-#define NDEBUG
+/* Print debug information? */
+#undef DEBUG
+
+/* Perform expensive self-checks? */
+#undef CHECK
+
+/* Omit all self-checks? (Overrides CHECK) */
+#undef NDEBUG
+
+/* Byte value for poison checks: use 255 for nan, or e.g. 113 for a
+ large value */
+#define POISON_VALUE 254
+
+
+
#include <assert.h>
#include <stdio.h>
#include <stdlib.h>
@@ -35,18 +49,6 @@ CCTK_FILEVERSION(TAT_Slab_slab_c);
-/* Print debug information? */
-#undef DEBUG
-
-/* Perform expensive self-checks? */
-#undef CHECK
-
-/* Byte value for poison checks: use 255 for nan, or e.g. 113 for a
- large value */
-#define POISON_VALUE 254
-
-
-
#ifdef DEBUG
# define ifdebug
#else
@@ -615,6 +617,7 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
/* assert (dstptr); */
CCTK_TimerStartI (timer_init);
+
assert (dim <= SLAB_MAXDIM);
info = malloc (SLAB_MAXDIM * sizeof *info);
assert (info);
@@ -737,6 +740,17 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
MPI_Comm_size (comm, &size);
MPI_Comm_rank (comm, &rank);
+ ifcheck {
+ static int count = 424242;
+ int mincount, maxcount;
+ ifdebug fflush (stdout);
+ MPI_Allreduce (&count, &mincount, 1, MPI_INT, MPI_MIN, comm);
+ MPI_Allreduce (&count, &maxcount, 1, MPI_INT, MPI_MAX, comm);
+ assert (mincount == count);
+ assert (maxcount == count);
+ ++ count;
+ }
+
assert (srctype == dsttype);
srctypesize = CCTK_VarTypeSize (srctype);
assert (srctypesize > 0);
@@ -930,6 +944,8 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
}
}
+ assert (srccount[rank] == dstcount[rank]);
+
ifcheck {
int * restrict src2count;
int * restrict dst2count;
@@ -962,6 +978,10 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
&& srctype == CCTK_VARIABLE_REAL) {
/* Optimised version for a special case */
+ int const srcoffi = info[0].src.local.off;
+ int const srcoffj = info[1].src.local.off;
+ int const srcoffk = info[2].src.local.off;
+
int const srcleni = info[0].src.local.len;
int const srclenj = info[1].src.local.len;
int const srclenk = info[2].src.local.len;
@@ -974,12 +994,18 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
int const srcdetaillenj = srcdetail[n*SLAB_MAXDIM+1].len;
int const srcdetaillenk = srcdetail[n*SLAB_MAXDIM+2].len;
+if (n==0) assert (srcoffset[n]==0);
+if (n<size-1) assert (srcoffset[n+1]==srcoffset[n]+srcdetailleni*srcdetaillenj*srcdetaillenk);
+
for (k = 0; k < srcdetaillenk; ++k) {
for (j = 0; j < srcdetaillenj; ++j) {
for (i = 0; i < srcdetailleni; ++i) {
- int const srcindi = srcdetailoffi + i;
- int const srcindj = srcdetailoffj + j;
- int const srcindk = srcdetailoffk + k;
+ int const srcindi = srcdetailoffi + i - srcoffi;
+ int const srcindj = srcdetailoffj + j - srcoffj;
+ int const srcindk = srcdetailoffk + k - srcoffk;
+ ifcheck assert (srcindi>=0 && srcindi<srcleni);
+ ifcheck assert (srcindj>=0 && srcindj<srclenj);
+ ifcheck assert (srcindk>=0 && srcindk<srclenk);
size_t const srcind = srcindi + srcleni * (srcindj + srclenj * srcindk);
size_t const bufind = i + srcdetailleni * (j + srcdetaillenj * k);
((CCTK_REAL*)srcdata)[srcoffset[n] + bufind] = ((const CCTK_REAL*)srcptr)[srcind];
@@ -1080,9 +1106,13 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
&& dsttype == CCTK_VARIABLE_REAL) {
/* Optimised version for a special case */
+ int const dstoffi = info[0].dst.local.off;
+ int const dstoffj = info[1].dst.local.off;
+ int const dstoffk = info[2].dst.local.off;
+
int const dstleni = info[0].dst.local.len;
int const dstlenj = info[1].dst.local.len;
- int const dstlenk = info[2].src.local.len;
+ int const dstlenk = info[2].dst.local.len;
int const dstdetailoffi = dstdetail[n*SLAB_MAXDIM+0].off;
int const dstdetailoffj = dstdetail[n*SLAB_MAXDIM+1].off;
@@ -1096,9 +1126,12 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
for (j = 0; j < dstdetaillenj; ++j) {
for (i = 0; i < dstdetailleni; ++i) {
size_t const bufind = i + dstdetailleni * (j + dstdetaillenj * k);
- int const dstindi = dstdetailoffi + i;
- int const dstindj = dstdetailoffj + j;
- int const dstindk = dstdetailoffk + k;
+ int const dstindi = dstdetailoffi + i - dstoffi;
+ int const dstindj = dstdetailoffj + j - dstoffj;
+ int const dstindk = dstdetailoffk + k - dstoffk;
+ ifcheck assert (dstindi>=0 && dstindi<dstleni);
+ ifcheck assert (dstindj>=0 && dstindj<dstlenj);
+ ifcheck assert (dstindk>=0 && dstindk<dstlenk);
size_t const dstind = dstindi + dstleni * (dstindj + dstlenj * dstindk);
((CCTK_REAL*)dstptr)[dstind] = ((const CCTK_REAL*)dstdata)[dstoffset[n] + bufind];
}