aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8>2005-05-10 08:15:24 +0000
committerschnetter <schnetter@2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8>2005-05-10 08:15:24 +0000
commit9186e24e77ab6b825e56d21191ec36bc9ea8d0e7 (patch)
tree77dd5f9e0f92939e76e32eeaafb1ada69abe6202
parent833f438bef6f5124db58796299ba99be3d6ba848 (diff)
Introduce new function Slab_MultiTransfer that transfers slabs of
several grid functions at once. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/Slab/trunk@40 2e825fa2-fb71-486d-8b7f-a5ff3f0f6cb8
-rw-r--r--src/slab.c168
-rw-r--r--src/slab.h10
2 files changed, 132 insertions, 46 deletions
diff --git a/src/slab.c b/src/slab.c
index f95f750..f4cea36 100644
--- a/src/slab.c
+++ b/src/slab.c
@@ -19,7 +19,7 @@
#undef CHECK
/* Omit all self-checks? (Overrides CHECK) */
-#undef NDEBUG
+#define NDEBUG
/* Byte value for poison checks: use 255 for nan, or e.g. 113 for a
large value */
@@ -294,13 +294,19 @@ get_mpi_comm (const cGH * restrict const cctkGH)
{
#ifdef CCTK_MPI
# if defined CARPET_CARPET || defined CARPETORIG_CARPET
- if (CCTK_IsThornActive ("Carpet")) {
- return CarpetMPIComm ();
+ {
+ static int Carpet_active = -1;
+ if (Carpet_active == -1) Carpet_active = CCTK_IsThornActive ("Carpet");
+ assert (Carpet_active >= 0);
+ if (Carpet_active) return CarpetMPIComm ();
}
# endif
# if defined CACTUSPUGH_PUGH
- if (CCTK_IsThornActive ("PUGH")) {
- return PUGH_pGH(cctkGH)->PUGH_COMM_WORLD;
+ {
+ static int PUGH_active = -1;
+ if (PUGH_active == -1) PUGH_active = CCTK_IsThornActive ("PUGH");
+ assert (PUGH_active >= 0);
+ if (PUGH_active) return PUGH_pGH(cctkGH)->PUGH_COMM_WORLD;
}
# endif
return MPI_COMM_WORLD;
@@ -612,14 +618,15 @@ static void bbox_xform (struct bbox * restrict const ydst,
-int Slab_Transfer (cGH const * restrict const cctkGH,
- int const dim,
- struct xferinfo const * restrict const xferinfo,
- int const options,
- int const srctype,
- void const * const srcptr,
- int const dsttype,
- void * const dstptr)
+int Slab_MultiTransfer (cGH const * const cctkGH,
+ int const dim,
+ struct xferinfo const * const xferinfo,
+ int const options,
+ int const nvars,
+ int const * const srctypes,
+ void const * const * const srcptrs,
+ int const * const dsttypes,
+ void * const * const dstptrs)
{
struct info * restrict info;
size_t srclentot, dstlentot;
@@ -628,8 +635,10 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
struct bbox * restrict srcdetail;
struct bbox * restrict dstdetail;
+ int * restrict srcelems;
int * restrict srccount;
int * restrict srcoffset;
+ int * restrict dstelems;
int * restrict dstcount;
int * restrict dstoffset;
@@ -642,6 +651,7 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
int size, rank;
MPI_Datatype srcdatatype, dstdatatype;
+ int var;
int i, j, k;
int n;
int d;
@@ -650,10 +660,15 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
assert (cctkGH);
assert (dim >= 0);
assert (xferinfo);
- assert (srctype >= 0);
-/* assert (srcptr); */
- assert (dsttype >= 0);
-/* assert (dstptr); */
+ assert (nvars >= 0);
+ assert (nvars==0 || srctypes);
+ for (var=0; var<nvars; ++var) assert (srctypes[var] >= 0);
+ assert (nvars==0 || srcptrs);
+/* for (var=0; var<nvars; ++var) assert (srcptrs[var]); */
+ assert (nvars==0 || dsttypes);
+ for (var=0; var<nvars; ++var) assert (dsttypes[var] >= 0);
+ assert (nvars==0 || dstptrs);
+/* for (var=0; var<nvars; ++var) assert (dstptrs[var]); */
CCTK_TimerStartI (timer_init);
@@ -790,6 +805,13 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
++ count;
}
+ assert (nvars >= 1);
+ int const srctype = srctypes[0];
+ int const dsttype = dsttypes[0];
+ for (var=0; var<nvars; ++var) {
+ assert (srctypes[var] == srctype);
+ assert (dsttypes[var] == dsttype);
+ }
assert (srctype == dsttype);
srctypesize = CCTK_VarTypeSize (srctype);
assert (srctypesize > 0);
@@ -887,29 +909,32 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
}
}
+ srcelems = malloc (size * sizeof *srcelems);
+ assert (srcelems);
srccount = malloc (size * sizeof *srccount);
assert (srccount);
srcoffset = malloc ((size + 1) * sizeof *srcoffset);
assert (srcoffset);
srcoffset[0] = 0;
for (n = 0; n < size; ++n) {
- srccount[n] = 1;
+ srcelems[n] = 1;
for (d=0; d<SLAB_MAXDIM; ++d) {
- srccount[n] *= srcdetail[n*SLAB_MAXDIM+d].len;
+ srcelems[n] *= srcdetail[n*SLAB_MAXDIM+d].len;
}
+ srccount[n] = nvars * srcelems[n];
ifdebug printf
("srccnt n=%d offset=%d count=%d\n", n, srcoffset[n], srccount[n]);
srcoffset[n+1] = srcoffset[n] + srccount[n];
}
srcdata = malloc (srcoffset[size] * srctypesize);
- assert (srcoffset[size] == 0 || srcdata);
+ assert (nvars==0 || srcdata);
ifcheck {
if (srctype == CCTK_VARIABLE_REAL) {
CCTK_REAL * restrict const srcdataptr = srcdata;
CCTK_REAL marker;
memset (&marker, POISON_VALUE, sizeof marker);
for (i = 0; i < srcoffset[size]; ++i) {
- memcpy (&srcdataptr[i], &marker, sizeof marker);
+ memcpy (&srcdataptr[i], &marker, sizeof marker);
}
}
}
@@ -956,29 +981,32 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
}
}
+ dstelems = malloc (size * sizeof *dstelems);
+ assert (dstelems);
dstcount = malloc (size * sizeof *dstcount);
assert (dstcount);
dstoffset = malloc ((size + 1) * sizeof *dstoffset);
assert (dstoffset);
dstoffset[0] = 0;
for (n = 0; n < size; ++n) {
- dstcount[n] = 1;
+ dstelems[n] = 1;
for (d=0; d<SLAB_MAXDIM; ++d) {
- dstcount[n] *= dstdetail[n*SLAB_MAXDIM+d].len;
+ dstelems[n] *= dstdetail[n*SLAB_MAXDIM+d].len;
}
+ dstcount[n] = nvars * dstelems[n];
ifdebug printf
("dstcnt n=%d offset=%d count=%d\n", n, dstoffset[n], dstcount[n]);
dstoffset[n+1] = dstoffset[n] + dstcount[n];
}
dstdata = malloc (dstoffset[size] * dsttypesize);
- assert (dstoffset[size] == 0 || dstdata);
+ assert (nvars==0 || dstdata);
ifcheck {
if (dsttype == CCTK_VARIABLE_REAL) {
CCTK_REAL * restrict const dstdataptr = dstdata;
CCTK_REAL marker;
memset (&marker, POISON_VALUE, sizeof marker);
for (i = 0; i < dstoffset[size]; ++i) {
- memcpy (&dstdataptr[i], &marker, sizeof marker);
+ memcpy (&dstdataptr[i], &marker, sizeof marker);
}
}
}
@@ -1047,7 +1075,9 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
ifcheck assert (srcindi>=0 && srcindi<srcleni);
ifcheck assert (srcindj>=0 && srcindj<srclenj);
ifcheck assert (srcindk>=0 && srcindk<srclenk);
- ((CCTK_REAL*)srcdata)[srcoffset[n] + bufind] = ((const CCTK_REAL*)srcptr)[srcind];
+ for (var=0; var<nvars; ++var) {
+ ((CCTK_REAL*)srcdata)[srcoffset[n] + var * srcelems[n] + bufind] = ((const CCTK_REAL*)srcptrs[var])[srcind];
+ }
}
}
}
@@ -1087,7 +1117,9 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
ifcheck assert (srcindi>=0 && srcindi<srcleni);
ifcheck assert (srcindj>=0 && srcindj<srclenj);
ifcheck assert (srcindk>=0 && srcindk<srclenk);
- ((CCTK_REAL*)srcdata)[srcoffset[n] + bufind] = ((const CCTK_REAL*)srcptr)[srcind];
+ for (var=0; var<nvars; ++var) {
+ ((CCTK_REAL*)srcdata)[srcoffset[n] + var * srcelems[n] + bufind] = ((const CCTK_REAL*)srcptrs[var])[srcind];
+ }
}
}
}
@@ -1126,8 +1158,10 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
}
assert (srcind < srclentot);
assert (bufind < (size_t)srccount[n]);
- ((CCTK_REAL*)srcdata)[srcoffset[n] + bufind]
- = ((const CCTK_REAL*)srcptr)[srcind];
+ for (var=0; var<nvars; ++var) {
+ ((CCTK_REAL*)srcdata)[srcoffset[n] + var * srcelems[n] + bufind]
+ = ((const CCTK_REAL*)srcptrs[var])[srcind];
+ }
}
}
}
@@ -1168,9 +1202,11 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
assert (bufind < (size_t)srccount[n]);
/* ((CCTK_REAL*)srcdata)[srcoffset[n] + bufind] */
/* = ((const CCTK_REAL*)srcptr)[srcind]; */
- memcpy ((char *) srcdata + srctypesize * (srcoffset[n] + bufind),
- (const char *) srcptr + srctypesize * srcind,
- srctypesize);
+ for (var=0; var<nvars; ++var) {
+ memcpy ((char *) srcdata + srctypesize * (srcoffset[n] + var * srcelems[n] + bufind),
+ (const char *) srcptrs[var] + srctypesize * srcind,
+ srctypesize);
+ }
}
}
}
@@ -1200,11 +1236,13 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
ifcheck {
if (dsttype == CCTK_VARIABLE_REAL) {
- const CCTK_REAL * restrict const dstdataptr = dstdata;
- CCTK_REAL marker;
- memset (&marker, POISON_VALUE, sizeof marker);
- for (i = 0; i < dstoffset[size]; ++i) {
- assert (memcmp(&dstdataptr[i], &marker, sizeof marker) != 0);
+ for (var=0; var<nvars; ++var) {
+ const CCTK_REAL * restrict const dstdataptr = dstdata;
+ CCTK_REAL marker;
+ memset (&marker, POISON_VALUE, sizeof marker);
+ for (i = 0; i < dstoffset[size]; ++i) {
+ assert (memcmp(&dstdataptr[i], &marker, sizeof marker) != 0);
+ }
}
}
}
@@ -1248,7 +1286,10 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
ifcheck assert (dstindi>=0 && dstindi<dstleni);
ifcheck assert (dstindj>=0 && dstindj<dstlenj);
ifcheck assert (dstindk>=0 && dstindk<dstlenk);
- ((CCTK_REAL*)dstptr)[dstind] = ((const CCTK_REAL*)dstdata)[dstoffset[n] + bufind];
+ for (var=0; var<nvars; ++var) {
+ ((CCTK_REAL*)dstptrs[var])[dstind]
+ = ((const CCTK_REAL*)dstdata)[dstoffset[n] + var * dstelems[n] + bufind];
+ }
}
}
}
@@ -1285,7 +1326,10 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
ifcheck assert (dstindi>=0 && dstindi<dstleni);
ifcheck assert (dstindj>=0 && dstindj<dstlenj);
ifcheck assert (dstindk>=0 && dstindk<dstlenk);
- ((CCTK_REAL*)dstptr)[dstind] = ((const CCTK_REAL*)dstdata)[dstoffset[n] + bufind];
+ for (var=0; var<nvars; ++var) {
+ ((CCTK_REAL*)dstptrs[var])[dstind]
+ = ((const CCTK_REAL*)dstdata)[dstoffset[n] + var * dstelems[n] + bufind];
+ }
}
}
}
@@ -1322,7 +1366,10 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
ifcheck assert (dstindi>=0 && dstindi<dstleni);
ifcheck assert (dstindj>=0 && dstindj<dstlenj);
ifcheck assert (dstindk>=0 && dstindk<dstlenk);
- ((CCTK_REAL*)dstptr)[dstind] = ((const CCTK_REAL*)dstdata)[dstoffset[n] + bufind];
+ for (var=0; var<nvars; ++var) {
+ ((CCTK_REAL*)dstptrs[var])[dstind]
+ = ((const CCTK_REAL*)dstdata)[dstoffset[n] + var * dstelems[n] + bufind];
+ }
}
}
}
@@ -1359,7 +1406,10 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
ifcheck assert (dstindi>=0 && dstindi<dstleni);
ifcheck assert (dstindj>=0 && dstindj<dstlenj);
ifcheck assert (dstindk>=0 && dstindk<dstlenk);
- ((CCTK_REAL*)dstptr)[dstind] = ((const CCTK_REAL*)dstdata)[dstoffset[n] + bufind];
+ for (var=0; var<nvars; ++var) {
+ ((CCTK_REAL*)dstptrs[var])[dstind]
+ = ((const CCTK_REAL*)dstdata)[dstoffset[n] + var * dstelems[n] + bufind];
+ }
}
}
}
@@ -1400,8 +1450,10 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
}
ifcheck assert (bufind < (size_t)dstcount[n]);
ifcheck assert (dstind < dstlentot);
- ((CCTK_REAL*)dstptr)[dstind]
- = ((const CCTK_REAL*)dstdata)[dstoffset[n] + bufind];
+ for (var=0; var<nvars; ++var) {
+ ((CCTK_REAL*)dstptrs[var])[dstind]
+ = ((const CCTK_REAL*)dstdata)[dstoffset[n] + var * dstelems[n] + bufind];
+ }
}
}
}
@@ -1444,9 +1496,11 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
ifcheck assert (dstind < dstlentot);
/* ((CCTK_REAL*)dstptr)[dstind] */
/* = ((const CCTK_REAL*)dstdata)[dstoffset[n] + bufind]; */
- memcpy ((char *) dstptr + dsttypesize * dstind,
- (const char *) dstdata + dsttypesize * (dstoffset[n] + bufind),
- dsttypesize);
+ for (var=0; var<nvars; ++var) {
+ memcpy ((char *) dstptrs[var] + dsttypesize * dstind,
+ (const char *) dstdata + dsttypesize * (dstoffset[n] + var * dstelems[n] + bufind),
+ dsttypesize);
+ }
}
}
}
@@ -1459,11 +1513,13 @@ int Slab_Transfer (cGH const * restrict const cctkGH,
free (dstdata);
+ free (dstelems);
free (dstcount);
free (dstoffset);
free (dstdetail);
free (srcdata);
+ free (srcelems);
free (srccount);
free (srcoffset);
free (srcdetail);
@@ -1549,3 +1605,23 @@ CCTK_FNAME(Slab_Transfer) (int * restrict const ierr,
free (xferinfo);
}
+
+
+
+int Slab_Transfer (cGH const * restrict const cctkGH,
+ int const dim,
+ struct xferinfo const * restrict const xferinfo,
+ int const options,
+ int const srctype,
+ void const * const srcptr,
+ int const dsttype,
+ void * const dstptr)
+{
+ int const nvars = 1;
+ int const srctypes[] = { srctype };
+ void const * const srcptrs[] = { srcptr };
+ int const dsttypes[] = { dsttype };
+ void * const dstptrs[] = { dstptr };
+ return Slab_MultiTransfer (cctkGH, dim, xferinfo, options,
+ nvars, srctypes, srcptrs, dsttypes, dstptrs);
+}
diff --git a/src/slab.h b/src/slab.h
index 0ac241b..a9e8b87 100644
--- a/src/slab.h
+++ b/src/slab.h
@@ -99,6 +99,16 @@ int Slab_Transfer (cGH const * const cctkGH,
int const dsttype,
void * const dstptr);
+int Slab_MultiTransfer (cGH const * const cctkGH,
+ int const dim,
+ struct xferinfo const * const xferinfo,
+ int const options,
+ int const nvars,
+ int const * const srctypes,
+ void const * const * const srcptrs,
+ int const * const dsttypes,
+ void * const * const dstptrs);
+
#ifdef __cplusplus
}
#endif