aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8>2003-02-28 16:05:44 +0000
committerschnetter <schnetter@2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8>2003-02-28 16:05:44 +0000
commit342ae4c153d733e9ffb3ce655437e3f17b91c129 (patch)
tree5bf90093e832cb74ac6c9fd9889f7061b2266408
parenta49d0b6ae9a948563470998c99e3ea30587efc61 (diff)
Completed the work for the previous commit.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/Slab/trunk@15 2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8
-rw-r--r--src/slab.c182
-rw-r--r--src/slab.inc18
2 files changed, 107 insertions, 93 deletions
diff --git a/src/slab.c b/src/slab.c
index 304adf4..026b002 100644
--- a/src/slab.c
+++ b/src/slab.c
@@ -620,25 +620,25 @@ int Slab_Transfer (cGH * restrict const cctkGH,
}
{
- int iflag[dim];
- for (d=0; d<dim; ++d) {
+ int iflag[SLAB_MAXDIM];
+ for (d=0; d<SLAB_MAXDIM; ++d) {
iflag[d] = 0;
}
- for (d=0; d<dim; ++d) {
+ for (d=0; d<SLAB_MAXDIM; ++d) {
assert (! iflag[info[d].xpose]);
iflag[info[d].xpose] = 1;
}
- for (d=0; d<dim; ++d) {
+ for (d=0; d<SLAB_MAXDIM; ++d) {
assert (iflag[d]);
}
- for (d=0; d<dim; ++d) {
+ for (d=0; d<SLAB_MAXDIM; ++d) {
assert (info[info[d].xpose].src.slab.len == info[d].dst.slab.len);
}
}
srclentot = 1;
dstlentot = 1;
- for (d=0; d<dim; ++d) {
+ for (d=0; d<SLAB_MAXDIM; ++d) {
srclentot *= info[d].src.local.len;
dstlentot *= info[d].dst.local.len;
}
@@ -667,72 +667,85 @@ int Slab_Transfer (cGH * restrict const cctkGH,
- allinfo = malloc (size * dim * sizeof *allinfo);
+ allinfo = malloc (size * SLAB_MAXDIM * sizeof *allinfo);
assert (allinfo);
{
int const info_nints = sizeof(struct info) / sizeof(int);
ifdebug fflush (stdout);
MPI_Allgather
- (info, dim * info_nints, MPI_INT,
- allinfo, dim * info_nints, MPI_INT, comm);
+ (info, SLAB_MAXDIM * info_nints, MPI_INT,
+ allinfo, SLAB_MAXDIM * info_nints, MPI_INT, comm);
}
for (n = 0; n < size; ++n) {
- for (d=0; d<dim; ++d) {
- assert (allinfo[n*dim+d].src.global.off == info[d].src.global.off);
- assert (allinfo[n*dim+d].src.global.len == info[d].src.global.len);
- assert (allinfo[n*dim+d].src.global.str == info[d].src.global.str);
- assert (allinfo[n*dim+d].dst.global.off == info[d].dst.global.off);
- assert (allinfo[n*dim+d].dst.global.len == info[d].dst.global.len);
- assert (allinfo[n*dim+d].dst.global.str == info[d].dst.global.str);
- assert (allinfo[n*dim+d].src.local.str == info[d].src.local.str);
- assert (allinfo[n*dim+d].dst.local.str == info[d].dst.local.str);
- assert (allinfo[n*dim+d].src.active.str == info[d].src.active.str);
- assert (allinfo[n*dim+d].dst.active.str == info[d].dst.active.str);
- assert (allinfo[n*dim+d].src.slab.str == info[d].src.slab.str);
- assert (allinfo[n*dim+d].dst.slab.str == info[d].dst.slab.str);
- assert (allinfo[n*dim+d].xpose == info[d].xpose);
- assert (allinfo[n*dim+d].flip == info[d].flip);
+ for (d=0; d<SLAB_MAXDIM; ++d) {
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].src.global.off == info[d].src.global.off);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].src.global.len == info[d].src.global.len);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].src.global.str == info[d].src.global.str);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].dst.global.off == info[d].dst.global.off);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].dst.global.len == info[d].dst.global.len);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].dst.global.str == info[d].dst.global.str);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].src.local.str == info[d].src.local.str);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].dst.local.str == info[d].dst.local.str);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].src.active.str == info[d].src.active.str);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].dst.active.str == info[d].dst.active.str);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].src.slab.str == info[d].src.slab.str);
+ assert
+ (allinfo[n*SLAB_MAXDIM+d].dst.slab.str == info[d].dst.slab.str);
+ assert (allinfo[n*SLAB_MAXDIM+d].xpose == info[d].xpose);
+ assert (allinfo[n*SLAB_MAXDIM+d].flip == info[d].flip);
}
}
- srcdetail = malloc (size * dim * sizeof *srcdetail);
+ srcdetail = malloc (size * SLAB_MAXDIM * sizeof *srcdetail);
assert (srcdetail);
for (n = 0; n < size; ++n) {
ifdebug printf ("srcdetail n=%d:\n", n);
- for (d=0; d<dim; ++d) {
- srcdetail[n*dim+d] = allinfo[n*dim+d].src.slab;
+ for (d=0; d<SLAB_MAXDIM; ++d) {
+ srcdetail[n*SLAB_MAXDIM+d] = allinfo[n*SLAB_MAXDIM+d].src.slab;
ifdebug printf (" src.slab d=%d ", d);
- ifdebug bbox_print (&srcdetail[n*dim+d]);
+ ifdebug bbox_print (&srcdetail[n*SLAB_MAXDIM+d]);
ifdebug printf ("\n");
- bbox_clip (&srcdetail[n*dim+d], &info[d].src.active);
+ bbox_clip (&srcdetail[n*SLAB_MAXDIM+d], &info[d].src.active);
ifdebug printf (" clipped with src.active d=%d ", d);
- ifdebug bbox_print (&srcdetail[n*dim+d]);
+ ifdebug bbox_print (&srcdetail[n*SLAB_MAXDIM+d]);
ifdebug printf ("\n");
}
- for (d=0; d<dim; ++d) {
+ for (d=0; d<SLAB_MAXDIM; ++d) {
struct bbox whereto;
struct bbox wherefrom;
- whereto = allinfo[n*dim+d].dst.slab;
+ whereto = allinfo[n*SLAB_MAXDIM+d].dst.slab;
ifdebug printf (" dst.slab d=%d ", info[d].xpose);
ifdebug bbox_print (&whereto);
ifdebug printf ("\n");
- bbox_clip (&whereto, &allinfo[n*dim+d].dst.active);
+ bbox_clip (&whereto, &allinfo[n*SLAB_MAXDIM+d].dst.active);
ifdebug printf (" whereto d=%d ", info[d].xpose);
ifdebug bbox_print (&whereto);
ifdebug printf ("\n");
bbox_xform
(&wherefrom, &whereto,
- &allinfo[n*dim+info[d].xpose].src.slab, &allinfo[n*dim+d].dst.slab,
+ &allinfo[n*SLAB_MAXDIM+info[d].xpose].src.slab,
+ &allinfo[n*SLAB_MAXDIM+d].dst.slab,
info[d].flip);
ifdebug printf (" wherefrom d=%d ", info[d].xpose);
ifdebug bbox_print (&wherefrom);
ifdebug printf ("\n");
- bbox_clip (&srcdetail[n*dim+info[d].xpose], &wherefrom);
+ bbox_clip (&srcdetail[n*SLAB_MAXDIM+info[d].xpose], &wherefrom);
ifdebug printf (" clipped with wherefrom d=%d ", info[d].xpose);
- ifdebug bbox_print (&srcdetail[n*dim+info[d].xpose]);
+ ifdebug bbox_print (&srcdetail[n*SLAB_MAXDIM+info[d].xpose]);
ifdebug printf ("\n");
}
}
@@ -744,8 +757,8 @@ int Slab_Transfer (cGH * restrict const cctkGH,
srcoffset[0] = 0;
for (n = 0; n < size; ++n) {
srccount[n] = 1;
- for (d=0; d<dim; ++d) {
- srccount[n] *= srcdetail[n*dim+d].len;
+ for (d=0; d<SLAB_MAXDIM; ++d) {
+ srccount[n] *= srcdetail[n*SLAB_MAXDIM+d].len;
}
ifdebug printf
("srccnt n=%d offset=%d count=%d\n", n, srcoffset[n], srccount[n]);
@@ -764,39 +777,40 @@ int Slab_Transfer (cGH * restrict const cctkGH,
}
}
- dstdetail = malloc (size * dim * sizeof *dstdetail);
+ dstdetail = malloc (size * SLAB_MAXDIM * sizeof *dstdetail);
assert (dstdetail);
for (n = 0; n < size; ++n) {
ifdebug printf ("dstdetail n=%d:\n", n);
- for (d=0; d<dim; ++d) {
+ for (d=0; d<SLAB_MAXDIM; ++d) {
struct bbox wherefrom;
struct bbox whereto;
- dstdetail[n*dim+d] = allinfo[n*dim+d].dst.slab;
+ dstdetail[n*SLAB_MAXDIM+d] = allinfo[n*SLAB_MAXDIM+d].dst.slab;
ifdebug printf (" dst.slab d=%d ", d);
- ifdebug bbox_print (&dstdetail[n*dim+d]);
+ ifdebug bbox_print (&dstdetail[n*SLAB_MAXDIM+d]);
ifdebug printf ("\n");
- bbox_clip (&dstdetail[n*dim+d], &info[d].dst.active);
+ bbox_clip (&dstdetail[n*SLAB_MAXDIM+d], &info[d].dst.active);
ifdebug printf (" clipped with dst.active d=%d ", d);
- ifdebug bbox_print (&dstdetail[n*dim+d]);
+ ifdebug bbox_print (&dstdetail[n*SLAB_MAXDIM+d]);
ifdebug printf ("\n");
- wherefrom = allinfo[n*dim+info[d].xpose].src.slab;
+ wherefrom = allinfo[n*SLAB_MAXDIM+info[d].xpose].src.slab;
ifdebug printf (" src.slab d=%d ", d);
- ifdebug bbox_print (&dstdetail[n*dim+d]);
+ ifdebug bbox_print (&dstdetail[n*SLAB_MAXDIM+d]);
ifdebug printf ("\n");
- bbox_clip (&wherefrom, &allinfo[n*dim+info[d].xpose].src.active);
+ bbox_clip (&wherefrom, &allinfo[n*SLAB_MAXDIM+info[d].xpose].src.active);
ifdebug printf (" wherefrom d=%d ", d);
- ifdebug bbox_print (&dstdetail[n*dim+d]);
+ ifdebug bbox_print (&dstdetail[n*SLAB_MAXDIM+d]);
ifdebug printf ("\n");
bbox_xform
(&whereto, &wherefrom,
- &allinfo[n*dim+d].dst.slab, &allinfo[n*dim+info[d].xpose].src.slab,
+ &allinfo[n*SLAB_MAXDIM+d].dst.slab,
+ &allinfo[n*SLAB_MAXDIM+info[d].xpose].src.slab,
info[d].flip);
ifdebug printf (" whereto d=%d ", d);
- ifdebug bbox_print (&dstdetail[n*dim+d]);
+ ifdebug bbox_print (&dstdetail[n*SLAB_MAXDIM+d]);
ifdebug printf ("\n");
- bbox_clip (&dstdetail[n*dim+d], &whereto);
+ bbox_clip (&dstdetail[n*SLAB_MAXDIM+d], &whereto);
ifdebug printf (" clipped with whereto d=%d ", d);
- ifdebug bbox_print (&dstdetail[n*dim+d]);
+ ifdebug bbox_print (&dstdetail[n*SLAB_MAXDIM+d]);
ifdebug printf ("\n");
}
}
@@ -808,8 +822,8 @@ int Slab_Transfer (cGH * restrict const cctkGH,
dstoffset[0] = 0;
for (n = 0; n < size; ++n) {
dstcount[n] = 1;
- for (d=0; d<dim; ++d) {
- dstcount[n] *= dstdetail[n*dim+d].len;
+ for (d=0; d<SLAB_MAXDIM; ++d) {
+ dstcount[n] *= dstdetail[n*SLAB_MAXDIM+d].len;
}
ifdebug printf
("dstcnt n=%d offset=%d count=%d\n", n, dstoffset[n], dstcount[n]);
@@ -849,38 +863,38 @@ int Slab_Transfer (cGH * restrict const cctkGH,
for (n = 0; n < size; ++n) {
- assert (dim == 3);
- for (k = 0; k < srcdetail[n*dim+info[2].xpose].len; ++k) {
- for (j = 0; j < srcdetail[n*dim+info[1].xpose].len; ++j) {
- for (i = 0; i < srcdetail[n*dim+info[0].xpose].len; ++i) {
- int ipos[3];
- int srcipos[3];
- int bufipos[3];
+ assert (SLAB_MAXDIM == 3);
+ for (k = 0; k < srcdetail[n*SLAB_MAXDIM+info[2].xpose].len; ++k) {
+ for (j = 0; j < srcdetail[n*SLAB_MAXDIM+info[1].xpose].len; ++j) {
+ for (i = 0; i < srcdetail[n*SLAB_MAXDIM+info[0].xpose].len; ++i) {
+ int ipos[SLAB_MAXDIM];
+ int srcipos[SLAB_MAXDIM];
+ int bufipos[SLAB_MAXDIM];
size_t srcind;
size_t bufind;
ipos[0] = i;
ipos[1] = j;
ipos[2] = k;
- for (d=0; d<dim; ++d) {
+ for (d=0; d<SLAB_MAXDIM; ++d) {
int const c = info[d].xpose;
- srcipos[c] = srcdetail[n*dim+c].off + ipos[d] * srcdetail[n*dim+c].str;
+ srcipos[c] = srcdetail[n*SLAB_MAXDIM+c].off + ipos[d] * srcdetail[n*SLAB_MAXDIM+c].str;
assert (srcipos[c] >= info[c].src.local.off
&& srcipos[c] < info[c].src.local.off + info[c].src.local.len);
- if (! (srcipos[c] >= allinfo[n*dim+c].src.slab.off
- && srcipos[c] <= allinfo[n*dim+c].src.slab.off + (allinfo[n*dim+c].src.slab.len - 1) * allinfo[n*dim+c].src.slab.str)) {
+ if (! (srcipos[c] >= allinfo[n*SLAB_MAXDIM+c].src.slab.off
+ && srcipos[c] <= allinfo[n*SLAB_MAXDIM+c].src.slab.off + (allinfo[n*SLAB_MAXDIM+c].src.slab.len - 1) * allinfo[n*SLAB_MAXDIM+c].src.slab.str)) {
}
- assert (srcipos[c] >= allinfo[n*dim+c].src.slab.off
- && srcipos[c] <= allinfo[n*dim+c].src.slab.off + (allinfo[n*dim+c].src.slab.len - 1) * allinfo[n*dim+c].src.slab.str);
- assert ((srcipos[c] - allinfo[n*dim+c].src.slab.off) % allinfo[n*dim+c].src.slab.str == 0);
+ assert (srcipos[c] >= allinfo[n*SLAB_MAXDIM+c].src.slab.off
+ && srcipos[c] <= allinfo[n*SLAB_MAXDIM+c].src.slab.off + (allinfo[n*SLAB_MAXDIM+c].src.slab.len - 1) * allinfo[n*SLAB_MAXDIM+c].src.slab.str);
+ assert ((srcipos[c] - allinfo[n*SLAB_MAXDIM+c].src.slab.off) % allinfo[n*SLAB_MAXDIM+c].src.slab.str == 0);
bufipos[d] = ipos[d];
- assert (bufipos[d] >= 0 && bufipos[d] < srcdetail[n*dim+c].len);
+ assert (bufipos[d] >= 0 && bufipos[d] < srcdetail[n*SLAB_MAXDIM+c].len);
}
srcind = 0;
bufind = 0;
- for (d=dim-1; d>=0; --d) {
+ for (d=SLAB_MAXDIM-1; d>=0; --d) {
int const c = info[d].xpose;
srcind = srcind * info[d].src.local.len + srcipos[d] - info[d].src.local.off;
- bufind = bufind * srcdetail[n*dim+c].len + bufipos[d];
+ bufind = bufind * srcdetail[n*SLAB_MAXDIM+c].len + bufipos[d];
}
assert (srcind < srclentot);
assert (bufind < (size_t)srccount[n]);
@@ -922,26 +936,26 @@ int Slab_Transfer (cGH * restrict const cctkGH,
}
for (n = 0; n < size; ++n) {
- assert (dim == 3);
- for (k = 0; k < dstdetail[n*dim+2].len; ++k) {
- for (j = 0; j < dstdetail[n*dim+1].len; ++j) {
- for (i = 0; i < dstdetail[n*dim+0].len; ++i) {
- int ipos[3];
- int bufipos[3];
- int dstipos[3];
+ assert (SLAB_MAXDIM == 3);
+ for (k = 0; k < dstdetail[n*SLAB_MAXDIM+2].len; ++k) {
+ for (j = 0; j < dstdetail[n*SLAB_MAXDIM+1].len; ++j) {
+ for (i = 0; i < dstdetail[n*SLAB_MAXDIM+0].len; ++i) {
+ int ipos[SLAB_MAXDIM];
+ int bufipos[SLAB_MAXDIM];
+ int dstipos[SLAB_MAXDIM];
size_t bufind;
size_t dstind;
ipos[0] = i;
ipos[1] = j;
ipos[2] = k;
- for (d=0; d<dim; ++d) {
+ for (d=0; d<SLAB_MAXDIM; ++d) {
if (! info[d].flip) {
bufipos[d] = ipos[d];
} else {
- bufipos[d] = dstdetail[n*dim+d].len - 1 - ipos[d];
+ bufipos[d] = dstdetail[n*SLAB_MAXDIM+d].len - 1 - ipos[d];
}
- assert (bufipos[d] >= 0 && bufipos[d] < dstdetail[n*dim+d].len);
- dstipos[d] = dstdetail[n*dim+d].off + ipos[d] * info[d].dst.slab.str;
+ assert (bufipos[d] >= 0 && bufipos[d] < dstdetail[n*SLAB_MAXDIM+d].len);
+ dstipos[d] = dstdetail[n*SLAB_MAXDIM+d].off + ipos[d] * info[d].dst.slab.str;
assert (dstipos[d] >= info[d].dst.local.off
&& dstipos[d] < info[d].dst.local.off + info[d].dst.local.len);
assert (dstipos[d] >= info[d].dst.slab.off
@@ -950,8 +964,8 @@ int Slab_Transfer (cGH * restrict const cctkGH,
}
bufind = 0;
dstind = 0;
- for (d=dim-1; d>=0; --d) {
- bufind = bufind * dstdetail[n*dim+d].len + bufipos[d];
+ for (d=SLAB_MAXDIM-1; d>=0; --d) {
+ bufind = bufind * dstdetail[n*SLAB_MAXDIM+d].len + bufipos[d];
dstind = dstind * info[d].dst.local.len + dstipos[d] - info[d].dst.local.off;
}
assert (bufind < (size_t)dstcount[n]);
diff --git a/src/slab.inc b/src/slab.inc
index a0f9bdc..395985a 100644
--- a/src/slab.inc
+++ b/src/slab.inc
@@ -17,17 +17,17 @@ interface
integer ierr
CCTK_POINTER cctkGH
integer dim
- integer src_gsh, src_lbnd, src_lsh
- integer src_lbbox, src_ubbox, src_nghostzones
- integer src_off, src_str, src_len
- integer dst_gsh, dst_lbnd, dst_lsh
- integer dst_lbbox, dst_ubbox, dst_nghostzones
- integer dst_off, dst_str, dst_len
- integer xpose, flip
+ integer src_gsh(dim), src_lbnd(dim), src_lsh(dim)
+ integer src_lbbox(dim), src_ubbox(dim), src_nghostzones(dim)
+ integer src_off(dim), src_str(dim), src_len(dim)
+ integer dst_gsh(dim), dst_lbnd(dim), dst_lsh(dim)
+ integer dst_lbbox(dim), dst_ubbox(dim), dst_nghostzones(dim)
+ integer dst_off(dim), dst_str(dim), dst_len(dim)
+ integer xpose(dim), flip(dim)
integer options
integer srctype
- CCTK_POINTER srcptr
+ CCTK_REAL srcptr(*)
integer dsttype
- CCTK_POINTER dstptr
+ CCTK_REAL dstptr(*)
end subroutine Slab_Transfer
end interface