aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authoreschnett <eschnett@2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8>2013-01-15 02:38:41 +0000
committereschnett <eschnett@2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8>2013-01-15 02:38:41 +0000
commit4508a50e32aea320593740dea0f382d9273490d9 (patch)
treeed22201c510ed98588579da1ab43652a44eefbba
parent423826345276afcab96e60ec3b17133f8526fd3b (diff)
API Update: Support array padding (cctk_ash)
Add new field ash to structures holding lsh. Use ash when calculating array indices. NOTE: This changes the API; all callers must be updated. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/Slab/trunk@89 2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8
-rw-r--r--src/slab.cc76
-rw-r--r--src/slab.h4
-rw-r--r--src/slab.inc8
3 files changed, 58 insertions, 30 deletions
diff --git a/src/slab.cc b/src/slab.cc
index 0194397..f234a86 100644
--- a/src/slab.cc
+++ b/src/slab.cc
@@ -30,14 +30,14 @@
#include <cstring>
#include <vector>
-#include "cctk.h"
-#include "cctk_Arguments.h"
-#include "cctk_Parameters.h"
-#include "cctk_DefineThorn.h"
-#include "util_ErrorCodes.h"
-#include "util_Table.h"
+#include <cctk.h>
+#include <cctk_Arguments.h>
+#include <cctk_Parameters.h>
+#include <cctk_DefineThorn.h>
+#include <util_ErrorCodes.h>
+#include <util_Table.h>
-#include "loopcontrol.h"
+#include <loopcontrol.h>
#ifdef CCTK_MPI
# include <mpi.h>
@@ -543,6 +543,7 @@ mpi_type (int const cactustype)
struct bbox {
int off; // offset
int len; // length
+ int alen; // allocated length
int str; // stride
};
@@ -586,6 +587,7 @@ bbox_check (bbox const * restrict const bbox)
{
assert (bbox);
assert (bbox->len >= 0);
+ assert (bbox->alen >= bbox->len);
assert (bbox->str > 0);
}
@@ -598,6 +600,7 @@ global2bbox (slabinfo const * restrict const slab,
assert (slab->gsh >= 0);
bbox->off = 0;
bbox->len = slab->gsh;
+ bbox->alen = slab->gsh;
bbox->str = 1;
bbox_check (bbox);
}
@@ -610,9 +613,11 @@ local2bbox (slabinfo const * restrict const slab,
assert (bbox);
assert (slab->lbnd >= 0);
assert (slab->lsh >= 0);
+ assert (slab->ash >= slab->lsh);
assert (slab->lbnd + slab->lsh <= slab->gsh);
bbox->off = slab->lbnd;
bbox->len = slab->lsh;
+ bbox->alen = slab->ash;
bbox->str = 1;
bbox_check (bbox);
}
@@ -627,6 +632,7 @@ active2bbox (slabinfo const * restrict const slab,
assert (useghosts == 0 or useghosts == 1);
assert (slab->lbnd >= 0);
assert (slab->lsh >= 0);
+ assert (slab->ash >= slab->lsh);
assert (slab->lbnd + slab->lsh <= slab->gsh);
assert (slab->lbbox == 0 or slab->lbbox == 1);
assert (slab->ubbox == 0 or slab->ubbox == 1);
@@ -635,6 +641,7 @@ active2bbox (slabinfo const * restrict const slab,
int const nughostzones = slab->ubbox or useghosts ? 0 : slab->nghostzones;
bbox->off = slab->lbnd + nlghostzones;
bbox->len = slab->lsh - nlghostzones - nughostzones;
+ bbox->alen = slab->ash;
bbox->str = 1;
bbox_check (bbox);
}
@@ -647,6 +654,7 @@ slab2bbox (slabinfo const * restrict const slab,
assert (bbox);
bbox->off = slab->off;
bbox->len = slab->len;
+ bbox->alen = slab->len;
bbox->str = slab->str;
bbox_check (bbox);
}
@@ -716,6 +724,7 @@ bbox_xform (bbox * restrict const ydst,
assert ((ysrc_last - xsrc->off) % ysrc->str == 0);
assert (ydst_last - ydst->off + ydst->str >= 0);
ydst->len = (ydst_last - ydst->off + ydst->str) / ydst->str;
+ ydst->alen = xdst->alen;
bbox_check (ydst);
}
@@ -728,6 +737,7 @@ print_slabinfo (FILE * const out,
{
fprintf (out, " gsh: %d\n", slabinfo->gsh);
fprintf (out, " lbnd: %d, lsh: %d\n", slabinfo->lbnd, slabinfo->lsh);
+ fprintf (out, " ash: %d\n", slabinfo->ash);
fprintf (out, " lbbox: %d, ubbox: %d, nghostzones: %d\n",
slabinfo->lbbox, slabinfo->ubbox, slabinfo->nghostzones);
fprintf (out, " off: %d, str: %d, len: %d\n",
@@ -774,6 +784,10 @@ copy_data (const vector<xfer> &info,
int const srclenj = info[1].src.local.len;
int const srclenk = info[2].src.local.len;
+ int const srcaleni = info[0].src.local.alen;
+ int const srcalenj = info[1].src.local.alen;
+ int const srcalenk = info[2].src.local.alen;
+
int const srcdetailoffi = srcdetail[n*SLAB_MAXDIM+0].off;
int const srcdetailoffj = srcdetail[n*SLAB_MAXDIM+1].off;
int const srcdetailoffk = srcdetail[n*SLAB_MAXDIM+2].off;
@@ -815,7 +829,7 @@ copy_data (const vector<xfer> &info,
# pragma omp parallel
CCTK_LOOP3(Slab_copy_in_noxpose, i,j,k,
0,0,0, srcdetailleni,srcdetaillenj,srcdetaillenk,
- srcleni,srclenj,srclenk)
+ srcaleni,srcalenj,srcalenk)
{
int const srcindi = srcdetailoffi + i - srcoffi;
int const srcindj = srcdetailoffj + j - srcoffj;
@@ -824,7 +838,7 @@ copy_data (const vector<xfer> &info,
ifcheck assert (srcindj>=0 and srcindj<srclenj);
ifcheck assert (srcindk>=0 and srcindk<srclenk);
size_t const srcind =
- srcindi + srcleni * (srcindj + srclenj * srcindk);
+ srcindi + srcaleni * (srcindj + srcalenj * srcindk);
size_t const bufind =
i + dstdetailleni * (j + dstdetaillenj * k);
srcdataptr[bufind] = srcptr[srcind];
@@ -839,7 +853,7 @@ copy_data (const vector<xfer> &info,
// Interchange i and j loops
CCTK_LOOP3(Slab_copy_in_xposexy, j,i,k,
0,0,0, srcdetaillenj,srcdetailleni,srcdetaillenk,
- srclenj,srcleni,srclenk)
+ srcalenj,srcaleni,srcalenk)
{
int const srcindi = srcdetailoffi + i - srcoffi;
int const srcindj = srcdetailoffj + j - srcoffj;
@@ -848,7 +862,7 @@ copy_data (const vector<xfer> &info,
ifcheck assert (srcindj>=0 and srcindj<srclenj);
ifcheck assert (srcindk>=0 and srcindk<srclenk);
size_t const srcind =
- srcindi + srcleni * (srcindj + srclenj * srcindk);
+ srcindi + srcaleni * (srcindj + srcalenj * srcindk);
size_t const bufind =
j + dstdetailleni * (i + dstdetaillenj * k);
srcdataptr[bufind] = srcptr[srcind];
@@ -862,7 +876,7 @@ copy_data (const vector<xfer> &info,
# pragma omp parallel
CCTK_LOOP3(Slab_copy_in_xposegeneral, i,j,k,
0,0,0, srcdetailleni,srcdetaillenj,srcdetaillenk,
- srcleni,srclenj,srclenk)
+ srcaleni,srcalenj,srcalenk)
{
int ipos[SLAB_MAXDIM];
ipos[0] = i;
@@ -875,7 +889,9 @@ copy_data (const vector<xfer> &info,
ifcheck assert (srcindj>=0 and srcindj<srclenj);
ifcheck assert (srcindk>=0 and srcindk<srclenk);
size_t const srcind =
- srcindi + srcleni * (srcindj + srclenj * srcindk);
+ srcindi + srcaleni * (srcindj + srcalenj * srcindk);
+ // TODO: rewrite this, moving the index onto dstdetaillen and
+ // outside the loop; then collapse all loop specialisations
size_t const bufind =
(ipos[xpose_x] + dstdetailleni *
(ipos[xpose_y] + dstdetaillenj * ipos[xpose_z]));
@@ -913,6 +929,10 @@ copy_data_back (const vector<xfer> &info,
int const dstlenj = info[1].dst.local.len;
int const dstlenk = info[2].dst.local.len;
+ int const dstaleni = info[0].dst.local.alen;
+ int const dstalenj = info[1].dst.local.alen;
+ int const dstalenk = info[2].dst.local.alen;
+
int const dstdetailoffi = dstdetail[n*SLAB_MAXDIM+0].off;
int const dstdetailoffj = dstdetail[n*SLAB_MAXDIM+1].off;
int const dstdetailoffk = dstdetail[n*SLAB_MAXDIM+2].off;
@@ -937,7 +957,7 @@ copy_data_back (const vector<xfer> &info,
# pragma omp parallel
CCTK_LOOP3(Slab_copy_back_noflip, i,j,k,
0,0,0, dstdetailleni,dstdetaillenj,dstdetaillenk,
- dstleni,dstlenj,dstlenk)
+ dstaleni,dstalenj,dstalenk)
{
int const dstindi = dstdetailoffi + i - dstoffi;
int const dstindj = dstdetailoffj + j - dstoffj;
@@ -946,7 +966,7 @@ copy_data_back (const vector<xfer> &info,
ifcheck assert (dstindj>=0 and dstindj<dstlenj);
ifcheck assert (dstindk>=0 and dstindk<dstlenk);
size_t const dstind =
- dstindi + dstleni * (dstindj + dstlenj * dstindk);
+ dstindi + dstaleni * (dstindj + dstalenj * dstindk);
size_t const bufind =
i + dstdetailleni * (j + dstdetaillenj * k);
dstptr[dstind] = dstdataptr[bufind];
@@ -960,7 +980,7 @@ copy_data_back (const vector<xfer> &info,
# pragma omp parallel
CCTK_LOOP3(Slab_copy_back_flipx, i,j,k,
0,0,0, dstdetailleni,dstdetaillenj,dstdetaillenk,
- dstleni,dstlenj,dstlenk)
+ dstaleni,dstalenj,dstalenk)
{
int const dstindi = dstdetailoffi + (dstdetailleni - 1 - i) - dstoffi;
int const dstindj = dstdetailoffj + j - dstoffj;
@@ -969,7 +989,7 @@ copy_data_back (const vector<xfer> &info,
ifcheck assert (dstindj>=0 and dstindj<dstlenj);
ifcheck assert (dstindk>=0 and dstindk<dstlenk);
size_t const dstind =
- dstindi + dstleni * (dstindj + dstlenj * dstindk);
+ dstindi + dstaleni * (dstindj + dstalenj * dstindk);
size_t const bufind =
i + dstdetailleni * (j + dstdetaillenj * k);
dstptr[dstind] = dstdataptr[bufind];
@@ -983,7 +1003,7 @@ copy_data_back (const vector<xfer> &info,
# pragma omp parallel
CCTK_LOOP3(Slab_copy_back_flipy, i,j,k,
0,0,0, dstdetailleni,dstdetaillenj,dstdetaillenk,
- dstleni,dstlenj,dstlenk)
+ dstaleni,dstalenj,dstalenk)
{
int const dstindi = dstdetailoffi + i - dstoffi;
int const dstindj = dstdetailoffj + (dstdetaillenj - 1 - j) - dstoffj;
@@ -992,7 +1012,7 @@ copy_data_back (const vector<xfer> &info,
ifcheck assert (dstindj>=0 and dstindj<dstlenj);
ifcheck assert (dstindk>=0 and dstindk<dstlenk);
size_t const dstind =
- dstindi + dstleni * (dstindj + dstlenj * dstindk);
+ dstindi + dstaleni * (dstindj + dstalenj * dstindk);
size_t const bufind =
i + dstdetailleni * (j + dstdetaillenj * k);
dstptr[dstind] = dstdataptr[bufind];
@@ -1006,7 +1026,7 @@ copy_data_back (const vector<xfer> &info,
# pragma omp parallel
CCTK_LOOP3(Slab_copy_back_flipxy, i,j,k,
0,0,0, dstdetailleni,dstdetaillenj,dstdetaillenk,
- dstleni,dstlenj,dstlenk)
+ dstaleni,dstalenj,dstalenk)
{
int const dstindi = dstdetailoffi + (dstdetailleni - 1 - i) - dstoffi;
int const dstindj = dstdetailoffj + (dstdetaillenj - 1 - j) - dstoffj;
@@ -1015,7 +1035,7 @@ copy_data_back (const vector<xfer> &info,
ifcheck assert (dstindj>=0 and dstindj<dstlenj);
ifcheck assert (dstindk>=0 and dstindk<dstlenk);
size_t const dstind =
- dstindi + dstleni * (dstindj + dstlenj * dstindk);
+ dstindi + dstaleni * (dstindj + dstalenj * dstindk);
size_t const bufind =
i + dstdetailleni * (j + dstdetaillenj * k);
dstptr[dstind] = dstdataptr[bufind];
@@ -1029,8 +1049,10 @@ copy_data_back (const vector<xfer> &info,
# pragma omp parallel
CCTK_LOOP3(Slab_copy_back_flipgeneral, i,j,k,
0,0,0, dstdetailleni,dstdetaillenj,dstdetaillenk,
- dstleni,dstlenj,dstlenk)
+ dstaleni,dstalenj,dstalenk)
{
+ // TODO: rewrite this, moving the flipping onto dstlen and
+ // outside the loop; then collapse all loop specialisations
int const dstindi =
dstdetailoffi + (flip_x ? dstdetailleni - 1 - i : i) - dstoffi;
int const dstindj =
@@ -1041,7 +1063,7 @@ copy_data_back (const vector<xfer> &info,
ifcheck assert (dstindj>=0 and dstindj<dstlenj);
ifcheck assert (dstindk>=0 and dstindk<dstlenk);
size_t const dstind =
- dstindi + dstleni * (dstindj + dstlenj * dstindk);
+ dstindi + dstaleni * (dstindj + dstalenj * dstindk);
size_t const bufind =
i + dstdetailleni * (j + dstdetaillenj * k);
dstptr[dstind] = dstdataptr[bufind];
@@ -1670,7 +1692,7 @@ Slab_MultiTransfer_Apply
size_t bufind = 0;
for (int 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;
+ srcind = srcind * info[d].src.local.alen + srcipos[d] - info[d].src.local.off;
bufind = bufind * srcdetail[n*SLAB_MAXDIM+c].len + bufipos[d];
}
assert (srcind < srclentot);
@@ -1834,7 +1856,7 @@ Slab_MultiTransfer_Apply
size_t dstind = 0;
for (int 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;
+ dstind = dstind * info[d].dst.local.alen + dstipos[d] - info[d].dst.local.off;
}
ifcheck assert (bufind < (size_t)dstcount[n]);
ifcheck assert (dstind < dstlentot);
@@ -1943,6 +1965,7 @@ CCTK_FNAME(Slab_Transfer) (int * restrict const ierr,
int const * restrict const src_gsh,
int const * restrict const src_lbnd,
int const * restrict const src_lsh,
+ int const * restrict const src_ash,
int const * restrict const src_lbbox,
int const * restrict const src_ubbox,
int const * restrict const src_nghostzones,
@@ -1952,6 +1975,7 @@ CCTK_FNAME(Slab_Transfer) (int * restrict const ierr,
int const * restrict const dst_gsh,
int const * restrict const dst_lbnd,
int const * restrict const dst_lsh,
+ int const * restrict const dst_ash,
int const * restrict const dst_lbbox,
int const * restrict const dst_ubbox,
int const * restrict const dst_nghostzones,
@@ -1972,6 +1996,7 @@ CCTK_FNAME(Slab_Transfer) (int * restrict const ierr,
xferinfo[d].src.gsh = src_gsh[d];
xferinfo[d].src.lbnd = src_lbnd[d];
xferinfo[d].src.lsh = src_lsh[d];
+ xferinfo[d].src.ash = src_ash[d];
xferinfo[d].src.lbbox = src_lbbox[d];
xferinfo[d].src.ubbox = src_ubbox[d];
xferinfo[d].src.nghostzones = src_nghostzones[d];
@@ -1982,6 +2007,7 @@ CCTK_FNAME(Slab_Transfer) (int * restrict const ierr,
xferinfo[d].dst.gsh = dst_gsh[d];
xferinfo[d].dst.lbnd = dst_lbnd[d];
xferinfo[d].dst.lsh = dst_lsh[d];
+ xferinfo[d].dst.ash = dst_ash[d];
xferinfo[d].dst.lbbox = dst_lbbox[d];
xferinfo[d].dst.ubbox = dst_ubbox[d];
xferinfo[d].dst.nghostzones = dst_nghostzones[d];
diff --git a/src/slab.h b/src/slab.h
index 9718c42..ee252e0 100644
--- a/src/slab.h
+++ b/src/slab.h
@@ -8,7 +8,7 @@ extern "C"
#include <stdio.h>
-#include "cctk.h"
+#include <cctk.h>
/*
* Slab_Transfer copies a slab from one array into a slab of another
@@ -58,6 +58,7 @@ extern "C"
* gsh >= 0
* lbnd >= 0
* lsh >= 0
+ * ash >= lsh
* lbnd + lsh <= gsh
* lbbox and ubbox must be booleans, i.e., either 0 or 1
* nghostzones >= 0
@@ -80,6 +81,7 @@ extern "C"
struct slabinfo {
int gsh;
int lbnd, lsh;
+ int ash;
int lbbox, ubbox, nghostzones;
int off, str, len;
};
diff --git a/src/slab.inc b/src/slab.inc
index 791cfd8..60fdea4 100644
--- a/src/slab.inc
+++ b/src/slab.inc
@@ -2,10 +2,10 @@
interface
subroutine Slab_Transfer (ierr, cctkGH, dim, &
- src_gsh, src_lbnd, src_lsh, &
+ src_gsh, src_lbnd, src_lsh, src_ash, &
src_lbbox, src_ubbox, src_nghostzones, &
src_off, src_str, src_len, &
- dst_gsh, dst_lbnd, dst_lsh, &
+ dst_gsh, dst_lbnd, dst_lsh, dst_ash, &
dst_lbbox, dst_ubbox, dst_nghostzones, &
dst_off, dst_str, dst_len, &
xpose, flip, &
@@ -16,10 +16,10 @@ interface
integer ierr
CCTK_POINTER cctkGH
integer dim
- integer src_gsh(dim), src_lbnd(dim), src_lsh(dim)
+ integer src_gsh(dim), src_lbnd(dim), src_lsh(dim), src_ash(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_gsh(dim), dst_lbnd(dim), dst_lsh(dim), dst_ash(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)