aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorschnetter <schnetter@20f44201-0f4f-0410-9130-e5fc2714a787>2005-05-10 08:17:28 +0000
committerschnetter <schnetter@20f44201-0f4f-0410-9130-e5fc2714a787>2005-05-10 08:17:28 +0000
commitcf42eef60ff3335c9d5915837ed270faf4466c54 (patch)
treec13274a88fa2886feca12ba7d5fd166c3b8207c5 /src
parent9c0c79442fcc16921165c2479144b30b5fef4e47 (diff)
Use the new function Slab_MultiTransfer to transfer all grid functions
at once. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/RotatingSymmetry180/trunk@19 20f44201-0f4f-0410-9130-e5fc2714a787
Diffstat (limited to 'src')
-rw-r--r--src/rotatingsymmetry180.c228
-rw-r--r--src/rotatingsymmetry180.h4
2 files changed, 134 insertions, 98 deletions
diff --git a/src/rotatingsymmetry180.c b/src/rotatingsymmetry180.c
index f25beb9..e0140d1 100644
--- a/src/rotatingsymmetry180.c
+++ b/src/rotatingsymmetry180.c
@@ -19,19 +19,19 @@
int BndRot180VI (cGH const * restrict const cctkGH,
- int const * restrict const stencil,
- int const vi)
+ int const nvars,
+ CCTK_INT const * restrict const vis)
{
DECLARE_CCTK_ARGUMENTS;
DECLARE_CCTK_PARAMETERS;
- int gi;
+ int * restrict gis;
cGroup group;
cGroupDynamicData data;
char * restrict fullname;
- void * restrict varptr;
+ void * restrict * restrict varptrs;
- int parities[3];
+ int (*paritiess)[3];
int global_bbox[6];
int global_lbnd[3], global_ubnd[3];
@@ -43,6 +43,9 @@ int BndRot180VI (cGH const * restrict const cctkGH,
int offset[3]; /* offset 0..1 due to avoid_origin */
struct xferinfo * restrict xferinfo;
+ int * restrict vartypes;
+
+ int var;
int alldirs[2];
int dir; /* direction of the symmetry face */
@@ -53,33 +56,45 @@ int BndRot180VI (cGH const * restrict const cctkGH,
/* Check arguments */
assert (cctkGH);
- assert (stencil);
- assert (vi>=0 && vi<CCTK_NumVars());
+ assert (nvars >= 0);
+ assert (nvars==0 || vis);
+ for (var=0; var<nvars; ++var) {
+ assert (vis[var]>=0 && vis[var]<CCTK_NumVars());
+ }
if (verbose) {
- fullname = CCTK_FullName(vi);
- assert (fullname);
- CCTK_VInfo (CCTK_THORNSTRING,
- "Applying 180 degree rotating symmetry boundary conditions to \"%s\"",
- fullname);
- free (fullname);
+ for (var=0; var<nvars; ++var) {
+ fullname = CCTK_FullName(vis[var]);
+ assert (fullname);
+ CCTK_VInfo (CCTK_THORNSTRING,
+ "Applying 180 degree rotating symmetry boundary conditions to \"%s\"",
+ fullname);
+ free (fullname);
+ }
}
/* Get and check group info */
- gi = CCTK_GroupIndexFromVarI (vi);
- assert (gi>=0 && gi<CCTK_NumGroups());
-
- ierr = CCTK_GroupData (gi, &group);
- assert (!ierr);
- assert (group.grouptype == CCTK_GF);
- assert (group.disttype == CCTK_DISTRIB_DEFAULT);
- assert (group.stagtype == 0);
-
- ierr = CCTK_GroupDynamicData (cctkGH, gi, &data);
- assert (!ierr);
-
- varptr = CCTK_VarDataPtrI (cctkGH, 0, vi);
- assert (varptr);
+ assert (nvars>0);
+ gis = malloc (nvars * sizeof *gis);
+ assert (nvars==0 || gis);
+ varptrs = malloc (nvars * sizeof *varptrs);
+ assert (nvars==0 || varptrs);
+ for (var=0; var<nvars; ++var) {
+ gis[var] = CCTK_GroupIndexFromVarI (vis[var]);
+ assert (gis[var]>=0 && gis[var]<CCTK_NumGroups());
+
+ ierr = CCTK_GroupData (gis[var], &group);
+ assert (!ierr);
+ assert (group.grouptype == CCTK_GF);
+ assert (group.disttype == CCTK_DISTRIB_DEFAULT);
+ assert (group.stagtype == 0);
+
+ ierr = CCTK_GroupDynamicData (cctkGH, gis[var], &data);
+ assert (!ierr);
+
+ varptrs[var] = CCTK_VarDataPtrI (cctkGH, 0, vis[var]);
+ assert (varptrs[var]);
+ }
for (d=0; d<group.dim; ++d) {
x0[d] = CCTK_ORIGIN_SPACE(d);
@@ -87,19 +102,23 @@ int BndRot180VI (cGH const * restrict const cctkGH,
}
/* find parity */
- {
+ paritiess = malloc (nvars * sizeof *paritiess);
+ assert (nvars==0 || paritiess);
+
+ for (var=0; var<nvars; ++var) {
+
char tensortypealias[100];
int firstvar, numvars;
int table;
int index;
- numvars = CCTK_NumVarsInGroupI(gi);
+ numvars = CCTK_NumVarsInGroupI(gis[var]);
assert (numvars>0);
- firstvar = CCTK_FirstVarIndexI(gi);
+ firstvar = CCTK_FirstVarIndexI(gis[var]);
assert (firstvar>=0);
- index = vi - firstvar;
+ index = vis[var] - firstvar;
assert (index>=0 && index<numvars);
- table = CCTK_GroupTagsTableI(gi);
+ table = CCTK_GroupTagsTableI(gis[var]);
assert (table>=0);
ierr = Util_TableGetString
@@ -108,7 +127,7 @@ int BndRot180VI (cGH const * restrict const cctkGH,
/* assume a scalar */
strcpy (tensortypealias, "scalar");
} else if (ierr<0) {
- char * groupname = CCTK_GroupName(gi);
+ char * groupname = CCTK_GroupName(gis[var]);
assert (groupname);
CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
"Error in tensor type alias declaration for group \"%s\"",
@@ -117,27 +136,27 @@ int BndRot180VI (cGH const * restrict const cctkGH,
}
if (CCTK_EQUALS (tensortypealias, "scalar")) {
- parities[0] = parities[1] = parities[2] = +1;
+ paritiess[var][0] = paritiess[var][1] = paritiess[var][2] = +1;
} else if (CCTK_EQUALS (tensortypealias, "u")
|| CCTK_EQUALS (tensortypealias, "d"))
{
assert (numvars == 3);
- parities[0] = parities[1] = parities[2] = +1;
- parities[index] = -1;
+ paritiess[var][0] = paritiess[var][1] = paritiess[var][2] = +1;
+ paritiess[var][index] = -1;
} else if (CCTK_EQUALS (tensortypealias, "dd_sym")) {
assert (numvars == 6);
- parities[0] = parities[1] = parities[2] = +1;
+ paritiess[var][0] = paritiess[var][1] = paritiess[var][2] = +1;
switch (index) {
case 0: break;
- case 1: parities[0] = parities[1] = -1; break;
- case 2: parities[0] = parities[2] = -1; break;
+ case 1: paritiess[var][0] = paritiess[var][1] = -1; break;
+ case 2: paritiess[var][0] = paritiess[var][2] = -1; break;
case 3: break;
- case 4: parities[1] = parities[2] = -1; break;
+ case 4: paritiess[var][1] = paritiess[var][2] = -1; break;
case 5: break;
default: assert(0);
}
} else {
- char * groupname = CCTK_GroupName(gi);
+ char * groupname = CCTK_GroupName(gis[var]);
assert (groupname);
CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
"Illegal tensor type alias for group \"%s\"",
@@ -145,11 +164,13 @@ int BndRot180VI (cGH const * restrict const cctkGH,
free (groupname);
}
- assert (abs(parities[0])==1 && abs(parities[1])==1
- && abs(parities[2])==1);
- }
+ assert (abs(paritiess[var][0])==1 && abs(paritiess[var][1])==1
+ && abs(paritiess[var][2])==1);
+
+ } /* for var */
{
+#if 0
int min_handle, max_handle;
CCTK_REAL local[6], global[6];
min_handle = CCTK_ReductionArrayHandle ("minimum");
@@ -171,6 +192,21 @@ int BndRot180VI (cGH const * restrict const cctkGH,
ierr = CCTK_ReduceLocArrayToArray1D
(cctkGH, -1, max_handle, local, global, 3, CCTK_VARIABLE_REAL);
for (d=0; d<3; ++d) global_ubnd[d] = (int)global[d];
+#else
+ int max_handle;
+ CCTK_INT local[12], global[12];
+ max_handle = CCTK_ReductionArrayHandle ("maximum");
+ if (max_handle<0) CCTK_WARN (0, "Could not obtain reduction handle");
+
+ for (d=0; d<6; ++d) local[ d] = cctkGH->cctk_bbox[d];
+ for (d=0; d<3; ++d) local[6+d] = -cctkGH->cctk_lbnd[d];
+ for (d=0; d<3; ++d) local[9+d] = cctkGH->cctk_ubnd[d];
+ ierr = CCTK_ReduceLocArrayToArray1D
+ (cctkGH, -1, max_handle, local, global, 12, CCTK_VARIABLE_INT);
+ for (d=0; d<6; ++d) global_bbox[d] = global[ d];
+ for (d=0; d<3; ++d) global_lbnd[d] = -global[6+d];
+ for (d=0; d<3; ++d) global_ubnd[d] = global[9+d];
+#endif
for (d=0; d<3; ++d) {
fake_bbox[2*d ] = data.lbnd[d] == global_lbnd[d];
@@ -205,7 +241,7 @@ int BndRot180VI (cGH const * restrict const cctkGH,
dir = 0;
otherdir = 1;
- assert (stencil[dir] >= 0);
+ assert (data.nghostzones[dir] >= 0);
if (global_bbox[2*dir]) {
@@ -215,14 +251,15 @@ int BndRot180VI (cGH const * restrict const cctkGH,
"xyz"[otherdir]);
}
- if (data.gsh[dir] < 2*stencil[dir] + offset[dir]) {
+ if (data.gsh[dir] < 2*data.nghostzones[dir] + offset[dir]) {
+ assert (nvars > 0);
CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
"The group \"%s\" has in the %c-direction only %d grid points. "
"This is not large enough for a 180 degree rotating symmetry boundary that is %d grid points wide. "
"The group needs to have at least %d grid points in that direction.",
- CCTK_GroupNameFromVarI(vi), "xyz"[dir], data.gsh[dir],
- stencil[dir],
- 2*stencil[dir] + offset);
+ CCTK_GroupNameFromVarI(vis[var]), "xyz"[dir], data.gsh[dir],
+ data.nghostzones[dir],
+ 2*data.nghostzones[dir] + offset);
}
/* Allocate slab transfer description */
@@ -255,47 +292,57 @@ int BndRot180VI (cGH const * restrict const cctkGH,
xferinfo[d].flip = 0;
}
- xferinfo[dir].src.off = stencil[dir] + offset[dir];
- xferinfo[dir].src.len = stencil[dir];
+ xferinfo[dir].src.off = data.nghostzones[dir] + offset[dir];
+ xferinfo[dir].src.len = data.nghostzones[dir];
xferinfo[dir].dst.off = 0;
- xferinfo[dir].dst.len = stencil[dir];
+ xferinfo[dir].dst.len = data.nghostzones[dir];
xferinfo[ dir].flip = 1;
xferinfo[otherdir].flip = 1;
- ierr = Slab_Transfer
+ vartypes = malloc (nvars * sizeof *vartypes);
+ assert (nvars==0 || vartypes);
+ for (var=0; var<nvars; ++var) {
+ vartypes[var] = group.vartype;
+ }
+
+ ierr = Slab_MultiTransfer
(cctkGH, group.dim, xferinfo, -1,
- group.vartype, varptr, group.vartype, varptr);
+ nvars, vartypes, varptrs, vartypes, varptrs);
assert (!ierr);
+ free (vartypes);
+
/* take parity into account */
if (cctkGH->cctk_bbox[2*dir]) {
- int parity;
-
- parity = parities[0] * parities[1];
- assert (abs(parity) == 1);
- if (parity == -1) {
- int imin[3], imax[3];
- int i, j, k;
- for (d=0; d<3; ++d) {
- imin[d] = xferinfo[d].dst.off;
- imax[d] = xferinfo[d].dst.off + xferinfo[d].dst.len;
- imin[d] -= cctk_lbnd[d];
- imax[d] -= cctk_lbnd[d];
- if (imin[d] < 0) imin[d] = 0;
- if (imax[d] >= cctk_lsh[d]) imax[d] = cctk_lsh[d];
- }
- assert (group.dim == 3);
- assert (group.vartype == CCTK_VARIABLE_REAL);
- for (k=imin[0]; k<imax[2]; ++k) {
- for (j=imin[1]; j<imax[1]; ++j) {
- for (i=imin[2]; i<imax[0]; ++i) {
- const int ind = CCTK_GFINDEX3D(cctkGH, i, j, k);
- ((CCTK_REAL *) varptr) [ind] *= -1;
+ for (var=0; var<nvars; ++var) {
+ int parity;
+
+ parity = paritiess[var][0] * paritiess[var][1];
+ assert (abs(parity) == 1);
+ if (parity == -1) {
+ int imin[3], imax[3];
+ int i, j, k;
+ for (d=0; d<3; ++d) {
+ imin[d] = xferinfo[d].dst.off;
+ imax[d] = xferinfo[d].dst.off + xferinfo[d].dst.len;
+ imin[d] -= cctk_lbnd[d];
+ imax[d] -= cctk_lbnd[d];
+ if (imin[d] < 0) imin[d] = 0;
+ if (imax[d] >= cctk_lsh[d]) imax[d] = cctk_lsh[d];
+ }
+ assert (group.dim == 3);
+ assert (group.vartype == CCTK_VARIABLE_REAL);
+ for (k=imin[0]; k<imax[2]; ++k) {
+ for (j=imin[1]; j<imax[1]; ++j) {
+ for (i=imin[2]; i<imax[0]; ++i) {
+ const int ind = CCTK_GFINDEX3D(cctkGH, i, j, k);
+ ((CCTK_REAL *) varptrs[var]) [ind] *= -1;
+ }
}
}
}
- }
+ } /* for var */
free (xferinfo);
@@ -303,6 +350,10 @@ int BndRot180VI (cGH const * restrict const cctkGH,
} /* if bbox */
+ free (gis);
+ free (varptrs);
+ free (paritiess);
+
return 0;
}
@@ -317,9 +368,6 @@ void Rot180_ApplyBC (CCTK_ARGUMENTS)
CCTK_INT * restrict faces;
CCTK_INT * restrict widths;
CCTK_INT * restrict tables;
- int vi;
- int dim;
- int * restrict stencil;
int i;
int ierr;
@@ -342,25 +390,13 @@ void Rot180_ApplyBC (CCTK_ARGUMENTS)
assert (ierr == nvars);
for (i=0; i<nvars; ++i) {
- vi = indices[i];
- assert (vi>=0 && vi<CCTK_NumVars());
-
+ assert (indices[i]>=0 && indices[i]<CCTK_NumVars());
assert (widths[i] >= 0);
-
- dim = CCTK_GroupDimFromVarI (vi);
- assert (dim>=0);
-
- stencil = malloc (dim * sizeof *stencil);
- assert (stencil);
- ierr = CCTK_GroupnghostzonesVI (cctkGH, dim, stencil, vi);
- assert (!ierr);
-
- ierr = BndRot180VI (cctkGH, stencil, vi);
- assert (!ierr);
-
- free (stencil);
}
+ ierr = BndRot180VI (cctkGH, nvars, indices);
+ assert (!ierr);
+
free (indices);
free (faces);
free (widths);
diff --git a/src/rotatingsymmetry180.h b/src/rotatingsymmetry180.h
index e560f20..9e01060 100644
--- a/src/rotatingsymmetry180.h
+++ b/src/rotatingsymmetry180.h
@@ -9,8 +9,8 @@
void Rot180_RegisterSymmetry (CCTK_ARGUMENTS);
int BndRot180VI (cGH const * restrict const cctkGH,
- int const * restrict const stencil,
- int const vi);
+ int const nvars,
+ CCTK_INT const * restrict const vis);
void Rot180_ApplyBC (CCTK_ARGUMENTS);
CCTK_INT