aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorschnetter <schnetter@c3c03602-0f4f-0410-b3fa-d2c81c8a7dc5>2005-05-10 08:50:05 +0000
committerschnetter <schnetter@c3c03602-0f4f-0410-b3fa-d2c81c8a7dc5>2005-05-10 08:50:05 +0000
commit05e5b613bac22a6bc1e98054a2c16c98d536fce1 (patch)
tree34731a6261c80970d05859dd571d68a6fb8d445c
parent99d623dcc21908b93893011c8ebba19940cd9035 (diff)
Use new function Slab_MultiTransfer to transfer all variables at once.
git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/RotatingSymmetry90/trunk@20 c3c03602-0f4f-0410-b3fa-d2c81c8a7dc5
-rw-r--r--src/rotatingsymmetry90.c368
-rw-r--r--src/rotatingsymmetry90.h8
2 files changed, 208 insertions, 168 deletions
diff --git a/src/rotatingsymmetry90.c b/src/rotatingsymmetry90.c
index 8db4b93..e797cac 100644
--- a/src/rotatingsymmetry90.c
+++ b/src/rotatingsymmetry90.c
@@ -57,25 +57,25 @@ static int convert_index (int const step,
int BndRot90VI (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 firstvar;
char tensortypealias[100];
int numvars;
int index;
int srcvi;
- void const * restrict srcptr;
- int parity;
+ void const * restrict * restrict srcptrs;
+ int * restrict parities;
int global_bbox[6];
int global_lbnd[3], global_ubnd[3];
@@ -87,9 +87,12 @@ int BndRot90VI (cGH const * restrict const cctkGH,
int offset[3]; /* offset 0..1 due to avoid_origin */
struct xferinfo * restrict xferinfo;
+ int * restrict vartypes;
int have_global_bbox, have_local_bbox;
+ int var;
+
int alldirs[2];
int step;
int ndirs; /* 1 for face, 2 for edge */
@@ -101,87 +104,52 @@ int BndRot90VI (cGH const * restrict const cctkGH,
/* Check arguments */
assert (cctkGH);
- assert (stencil);
- assert (vi>=0 && vi<CCTK_NumVars());
+ 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 90 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 90 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);
dx[d] = CCTK_DELTA_SPACE(d);
}
- /* Find tensor type, source variable, and parity */
- {
- int table;
-
- numvars = CCTK_NumVarsInGroupI(gi);
- assert (numvars>0);
- firstvar = CCTK_FirstVarIndexI(gi);
- assert (firstvar>=0);
- index = vi - firstvar;
- assert (index>=0 && index<numvars);
- table = CCTK_GroupTagsTableI(gi);
- assert (table>=0);
-
- ierr = Util_TableGetString
- (table, sizeof tensortypealias, tensortypealias, "tensortypealias");
- if (ierr == UTIL_ERROR_TABLE_NO_SUCH_KEY) {
- /* assume a scalar */
- strcpy (tensortypealias, "scalar");
- } else if (ierr<0) {
- char * groupname = CCTK_GroupName(gi);
- assert (groupname);
- CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
- "Error in tensor type alias declaration for group \"%s\"",
- groupname);
- free (groupname);
- }
-
- if (CCTK_EQUALS (tensortypealias, "scalar")) {
- /* scalar */
- } else if (CCTK_EQUALS (tensortypealias, "u")
- || CCTK_EQUALS (tensortypealias, "d"))
- {
- /* vector */
- assert (numvars == 3);
- } else if (CCTK_EQUALS (tensortypealias, "dd_sym")) {
- /* symmetric tensor */
- assert (numvars == 6);
- } else {
- char * groupname = CCTK_GroupName(gi);
- assert (groupname);
- CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
- "Illegal tensor type alias for group \"%s\"",
- groupname);
- free (groupname);
- }
- }
-
{
+#if 0
int min_handle, max_handle;
CCTK_REAL local[6], global[6];
min_handle = CCTK_ReductionArrayHandle ("minimum");
@@ -203,6 +171,21 @@ int BndRot90VI (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];
@@ -232,6 +215,11 @@ int BndRot90VI (cGH const * restrict const cctkGH,
offset[d] = avoid_origin[d] ? 0 : 1;
}
+ parities = malloc (nvars * sizeof *parities);
+ assert (nvars==0 || parities);
+ srcptrs = malloc (nvars * sizeof *srcptrs);
+ assert (nvars==0 || srcptrs);
+
for (step=0; step<3; ++step) {
switch (step) {
case 0:
@@ -258,42 +246,93 @@ int BndRot90VI (cGH const * restrict const cctkGH,
assert (0);
}
- parity = +1;
- if (CCTK_EQUALS (tensortypealias, "scalar")) {
- /* scalar */
- srcvi = vi;
- /* do nothing */
- } else if (CCTK_EQUALS (tensortypealias, "u")
- || CCTK_EQUALS (tensortypealias, "d"))
- {
- /* vector */
- int srcindex;
- srcindex = convert_index (step, index, alldirs, &parity);
- srcvi = firstvar + srcindex;
- } else if (CCTK_EQUALS (tensortypealias, "dd_sym")) {
- /* symmetric tensor */
- int index1, index2;
- int srcindex1, srcindex2;
- int srcindex;
+ for (var=0; var<nvars; ++var) {
+
+ /* Find tensor type, source variable, and parity */
{
- int const expand1[6] = { 0,0,0,1,1,2 };
- int const expand2[6] = { 0,1,2,1,2,2 };
- index1 = expand1[index];
- index2 = expand2[index];
+ int table;
+
+ numvars = CCTK_NumVarsInGroupI(gis[var]);
+ assert (numvars>0);
+ firstvar = CCTK_FirstVarIndexI(gis[var]);
+ assert (firstvar>=0);
+ index = vis[var] - firstvar;
+ assert (index>=0 && index<numvars);
+ table = CCTK_GroupTagsTableI(gis[var]);
+ assert (table>=0);
+
+ ierr = Util_TableGetString
+ (table, sizeof tensortypealias, tensortypealias, "tensortypealias");
+ if (ierr == UTIL_ERROR_TABLE_NO_SUCH_KEY) {
+ /* assume a scalar */
+ strcpy (tensortypealias, "scalar");
+ } else if (ierr<0) {
+ char * groupname = CCTK_GroupName(gis[var]);
+ assert (groupname);
+ CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "Error in tensor type alias declaration for group \"%s\"",
+ groupname);
+ free (groupname);
+ }
+
+ if (CCTK_EQUALS (tensortypealias, "scalar")) {
+ /* scalar */
+ } else if (CCTK_EQUALS (tensortypealias, "u")
+ || CCTK_EQUALS (tensortypealias, "d"))
+ {
+ /* vector */
+ assert (numvars == 3);
+ } else if (CCTK_EQUALS (tensortypealias, "dd_sym")) {
+ /* symmetric tensor */
+ assert (numvars == 6);
+ } else {
+ char * groupname = CCTK_GroupName(gis[var]);
+ assert (groupname);
+ CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
+ "Illegal tensor type alias for group \"%s\"",
+ groupname);
+ free (groupname);
+ }
}
- srcindex1 = convert_index (step, index1, alldirs, &parity);
- srcindex2 = convert_index (step, index2, alldirs, &parity);
+
+ parities[var] = +1;
+ if (CCTK_EQUALS (tensortypealias, "scalar")) {
+ /* scalar */
+ srcvi = vis[var];
+ /* do nothing */
+ } else if (CCTK_EQUALS (tensortypealias, "u")
+ || CCTK_EQUALS (tensortypealias, "d"))
{
- int const compact[3][3] = { { 0,1,2 }, { 1,3,4 }, { 2,4,5 } };
- srcindex = compact[srcindex1][srcindex2];
+ /* vector */
+ int srcindex;
+ srcindex = convert_index (step, index, alldirs, &parities[var]);
+ srcvi = firstvar + srcindex;
+ } else if (CCTK_EQUALS (tensortypealias, "dd_sym")) {
+ /* symmetric tensor */
+ int index1, index2;
+ int srcindex1, srcindex2;
+ int srcindex;
+ {
+ int const expand1[6] = { 0,0,0,1,1,2 };
+ int const expand2[6] = { 0,1,2,1,2,2 };
+ index1 = expand1[index];
+ index2 = expand2[index];
+ }
+ srcindex1 = convert_index (step, index1, alldirs, &parities[var]);
+ srcindex2 = convert_index (step, index2, alldirs, &parities[var]);
+ {
+ int const compact[3][3] = { { 0,1,2 }, { 1,3,4 }, { 2,4,5 } };
+ srcindex = compact[srcindex1][srcindex2];
+ }
+ srcvi = firstvar + srcindex;
+ } else {
+ assert (0);
}
- srcvi = firstvar + srcindex;
- } else {
- assert (0);
- }
-
- srcptr = CCTK_VarDataPtrI (cctkGH, 0, srcvi);
- assert (srcptr);
+
+ srcptrs[var] = CCTK_VarDataPtrI (cctkGH, 0, srcvi);
+ assert (srcptrs[var]);
+
+ } /* for var */
have_global_bbox = 1;
for (q=0; q<ndirs; ++q) {
@@ -302,16 +341,17 @@ int BndRot90VI (cGH const * restrict const cctkGH,
if (have_global_bbox) {
for (q=0; q<ndirs; ++q) {
- assert (stencil[dir[q]] >= 0);
+ assert (data.nghostzones[dir[q]] >= 0);
- if (data.gsh[otherdir[q]] < stencil[dir[q]] + stencil[otherdir[q]] + offset[otherdir[q]]) {
+ if (data.gsh[otherdir[q]] < data.nghostzones[dir[q]] + data.nghostzones[otherdir[q]] + offset[otherdir[q]]) {
+ 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 90 degree rotating symmetry boundary that is %d grid points wide in the %c-direction. "
- "The group needs to have at least %d grid points in the %c-direction.",
- CCTK_GroupNameFromVarI(vi), "xyz"[dir[q]], data.gsh[dir[q]],
- stencil[dir[q]], "xyz"[otherdir[q]],
- stencil[dir[q]] + stencil[otherdir[q]] + offset[otherdir[q]], "xyz"[dir[q]]);
+ "The group needs to have at least %d grid points in the %c-direction.",
+ CCTK_GroupNameFromVarI(vis[0]), "xyz"[dir[q]], data.gsh[dir[q]],
+ data.nghostzones[dir[q]], "xyz"[otherdir[q]],
+ data.nghostzones[dir[q]] + data.nghostzones[otherdir[q]] + offset[otherdir[q]], "xyz"[dir[q]]);
}
}
@@ -349,10 +389,10 @@ int BndRot90VI (cGH const * restrict const cctkGH,
case 1:
/* face */
for (q=0; q<ndirs; ++q) {
- xferinfo[otherdir[q]].src.off = stencil[otherdir[q]] + offset[otherdir[q]];
- xferinfo[otherdir[q]].src.len = stencil[ dir[q]];
+ xferinfo[otherdir[q]].src.off = data.nghostzones[otherdir[q]] + offset[otherdir[q]];
+ xferinfo[otherdir[q]].src.len = data.nghostzones[ dir[q]];
xferinfo[ dir[q]].dst.off = 0;
- xferinfo[ dir[q]].dst.len = stencil[ dir[q]];
+ xferinfo[ dir[q]].dst.len = data.nghostzones[ dir[q]];
xferinfo[ dir[q]].xpose = otherdir[q];
xferinfo[otherdir[q]].xpose = dir[q];
@@ -364,10 +404,10 @@ int BndRot90VI (cGH const * restrict const cctkGH,
case 2:
/* edge */
for (q=0; q<ndirs; ++q) {
- xferinfo[otherdir[q]].src.off = stencil[otherdir[q]] + offset[otherdir[q]];
- xferinfo[otherdir[q]].src.len = stencil[ dir[q]];
+ xferinfo[otherdir[q]].src.off = data.nghostzones[otherdir[q]] + offset[otherdir[q]];
+ xferinfo[otherdir[q]].src.len = data.nghostzones[ dir[q]];
xferinfo[ dir[q]].dst.off = 0;
- xferinfo[ dir[q]].dst.len = stencil[ dir[q]];
+ xferinfo[ dir[q]].dst.len = data.nghostzones[ dir[q]];
xferinfo[dir[q]].flip = 1;
}
@@ -376,40 +416,50 @@ int BndRot90VI (cGH const * restrict const cctkGH,
assert (0);
}
- 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, srcptr, group.vartype, varptr);
+ nvars, vartypes, srcptrs, vartypes, varptrs);
assert (!ierr);
-
+
+ free (vartypes);
+
/* take parity into account */
have_local_bbox = 1;
for (q=0; q<ndirs; ++q) {
have_local_bbox = have_local_bbox && cctk_bbox[2*dir[q]];
}
if (have_local_bbox) {
- 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) {
+ assert (abs(parities[var]) == 1);
+ if (parities[var] == -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 */
} /* if bbox */
@@ -419,6 +469,11 @@ int BndRot90VI (cGH const * restrict const cctkGH,
} /* for step */
+ free (gis);
+ free (varptrs);
+ free (parities);
+ free (srcptrs);
+
return 0;
}
@@ -433,10 +488,7 @@ void Rot90_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 var;
int ierr;
assert (cctkGH);
@@ -457,26 +509,14 @@ void Rot90_ApplyBC (CCTK_ARGUMENTS)
(cctkGH, nvars, indices, faces, widths, tables, 0);
assert (ierr == nvars);
- for (i=0; i<nvars; ++i) {
- vi = indices[i];
- assert (vi>=0 && vi<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 = BndRot90VI (cctkGH, stencil, vi);
- assert (!ierr);
-
- free (stencil);
+ for (var=0; var<nvars; ++var) {
+ assert (indices[var]>=0 && indices[var]<CCTK_NumVars());
+ assert (widths[var] >= 0);
}
+ ierr = BndRot90VI (cctkGH, nvars, indices);
+ assert (!ierr);
+
free (indices);
free (faces);
free (widths);
diff --git a/src/rotatingsymmetry90.h b/src/rotatingsymmetry90.h
index f7b408d..aa4561b 100644
--- a/src/rotatingsymmetry90.h
+++ b/src/rotatingsymmetry90.h
@@ -9,10 +9,10 @@
void Rot90_RegisterSymmetry (CCTK_ARGUMENTS);
int BndRot90VI (cGH const * restrict const cctkGH,
- int const * restrict const stencil,
- int const vi);
-void
-Rot90_ApplyBC (CCTK_ARGUMENTS);
+ int const nvars,
+ CCTK_INT const * restrict const vis);
+
+void Rot90_ApplyBC (CCTK_ARGUMENTS);
CCTK_INT
Rot90_SymmetryInterpolate (CCTK_POINTER_TO_CONST const cctkGH_,