aboutsummaryrefslogtreecommitdiff
path: root/src
diff options
context:
space:
mode:
authorschnetter <schnetter@20f44201-0f4f-0410-9130-e5fc2714a787>2010-02-14 18:28:16 +0000
committerschnetter <schnetter@20f44201-0f4f-0410-9130-e5fc2714a787>2010-02-14 18:28:16 +0000
commitec93cf628d1c0d132c6e46c43202f234acd0c067 (patch)
tree999442aadb48d0251954340cc38083832059974d /src
parent53d0e813519bfa951e971c5b0a50e1dc3b45899e (diff)
Determine the symmetry boundary width from the registered boundary
width in CoordBase, not from the number of ghost zones. Abort if there are more than one components per MPI process; TAT/Slab cannot handle this. Begin to add support for CarpetSlab instead of TAT/Slab. Remove cvs $Header$ comments. git-svn-id: http://svn.cactuscode.org/arrangements/CactusNumerical/RotatingSymmetry180/trunk@48 20f44201-0f4f-0410-9130-e5fc2714a787
Diffstat (limited to 'src')
-rw-r--r--src/interpolate.c2
-rw-r--r--src/make.code.defn1
-rw-r--r--src/registersymmetry.c17
-rw-r--r--src/rotatingsymmetry180.c79
-rw-r--r--src/rotatingsymmetry180.h2
5 files changed, 80 insertions, 21 deletions
diff --git a/src/interpolate.c b/src/interpolate.c
index 717d352..ecf7477 100644
--- a/src/interpolate.c
+++ b/src/interpolate.c
@@ -1,5 +1,3 @@
-/* $Header$ */
-
#include <assert.h>
#include <math.h>
#include <stdio.h>
diff --git a/src/make.code.defn b/src/make.code.defn
index 03f55f4..bdd2ba6 100644
--- a/src/make.code.defn
+++ b/src/make.code.defn
@@ -1,5 +1,4 @@
# Main make.code.defn file for thorn RotatingSymmetry180
-# $Header$
# Source files in this directory
SRCS = interpolate.c rotatingsymmetry180.c registersymmetry.c
diff --git a/src/registersymmetry.c b/src/registersymmetry.c
index b326181..6ce17dc 100644
--- a/src/registersymmetry.c
+++ b/src/registersymmetry.c
@@ -1,5 +1,3 @@
-/* $Header$ */
-
#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"
@@ -13,19 +11,32 @@ void Rot180_RegisterSymmetry (CCTK_ARGUMENTS)
DECLARE_CCTK_ARGUMENTS;
DECLARE_CCTK_PARAMETERS;
+ CCTK_INT nboundaryzones[6];
+ CCTK_INT is_internal[6];
+ CCTK_INT is_staggered[6];
+ CCTK_INT shiftout[6];
+
int f;
CCTK_INT handle;
CCTK_INT faces[6];
CCTK_INT width[6];
CCTK_INT ierr;
+ /* Get the boundary specification */
+ ierr = GetBoundarySpecification
+ (6, nboundaryzones, is_internal, is_staggered, shiftout);
+ if (ierr < 0)
+ {
+ CCTK_WARN (0, "Could not get the boundary specification");
+ }
+
for (f=0; f<6; ++f) {
faces[f] = 0;
width[f] = 0;
}
faces[0] = 1;
- width[0] = cctk_nghostzones[0];
+ width[0] = nboundaryzones[0];
handle = SymmetryRegister ("rotating_symmetry_180");
if (handle < 0) {
diff --git a/src/rotatingsymmetry180.c b/src/rotatingsymmetry180.c
index 62b859d..ccafa3c 100644
--- a/src/rotatingsymmetry180.c
+++ b/src/rotatingsymmetry180.c
@@ -1,5 +1,3 @@
-/* $Header$ */
-
#include <assert.h>
#include <math.h>
#include <stdlib.h>
@@ -45,7 +43,6 @@ int BndRot180VI (cGH const * restrict const cctkGH,
int offset[3]; /* offset 0..1 due to avoid_origin */
struct xferinfo * restrict xferinfo;
- int options;
int var;
@@ -54,6 +51,7 @@ int BndRot180VI (cGH const * restrict const cctkGH,
int otherdir; /* the other direction of the rotation */
int q;
int d;
+ int icnt;
int ierr;
/* Check arguments */
@@ -451,16 +449,71 @@ int BndRot180VI (cGH const * restrict const cctkGH,
xferinfo[ dir].flip = 1;
xferinfo[otherdir].flip = 1;
- options = Util_TableCreateFromString ("useghosts=1");
- assert (options>=0);
-
- ierr = Slab_MultiTransfer
- (cctkGH, group.dim, xferinfo, options,
- nvars, vartypes, varptrs, vartypes, varptrs);
- assert (!ierr);
-
- ierr = Util_TableDestroy (options);
- assert (!ierr);
+ if (CCTK_EQUALS (hyperslabber, "TAT/Slab")) {
+
+ if (CCTK_IsFunctionAliased("GetLocalComponents")) {
+ int const num_local_components = GetLocalComponents(cctkGH);
+ if (num_local_components > 1) {
+ CCTK_WARN (CCTK_WARN_ABORT,
+ "TAT/Slab can only be used if there is a single local component per MPI process");
+ }
+ }
+
+ int options;
+ options = Util_TableCreateFromString ("useghosts=1");
+ assert (options>=0);
+
+ ierr = Slab_MultiTransfer
+ (cctkGH, group.dim, xferinfo, options,
+ nvars, vartypes, varptrs, vartypes, varptrs);
+ assert (!ierr);
+
+ ierr = Util_TableDestroy (options);
+ assert (!ierr);
+
+ } else if (CCTK_EQUALS (hyperslabber, "GetHyperslab")) {
+
+ CCTK_INT const direction[3*3] = {1, 0, 0, 0, 1, 0, 0, 0, 1};
+ CCTK_INT const origin[3] =
+ {xferinfo[0].src.off, xferinfo[1].src.off, xferinfo[2].src.off};
+ CCTK_INT const extent[3] =
+ {xferinfo[0].src.len, xferinfo[1].src.len, xferinfo[2].src.len};
+ CCTK_INT const downsample[3] =
+ {xferinfo[0].src.str, xferinfo[1].src.str, xferinfo[2].src.str};
+ int mapping;
+ CCTK_INT hsize[3], total_hsize;
+ CCTK_POINTER hdata[nvars];
+
+ assert (nvars > 0);
+ mapping = Hyperslab_GlobalMappingByIndex
+ (cctkGH, vis[0], 3, direction, origin, extent, downsample,
+ -1, NULL, &hsize);
+ assert (mapping>=0);
+
+ total_hsize = hsize[0] * hsize[1] * hsize[2];
+ assert (total_hsize >= 0);
+ for (var=0; var<nvars; ++var) {
+ hdata[var] = malloc(total_hsize * CCTK_VarTypeSize(CCTK_VarTypeI(var)));
+ assert (total_hsize==0 || hdata[var]);
+ }
+
+ icnt = Hyperslab_GetList
+ (cctkGH, mapping, nvars, NULL, vis, NULL, NULL, hdata, NULL);
+ assert (icnt == nvars);
+
+ // copy hyperslabs into grid functions
+ assert (0);
+
+ for (var=0; var<nvars; ++var) {
+ free (hdata[var]);
+ }
+
+ ierr = Hyperslab_FreeMapping (mapping);
+ assert (!ierr);
+
+ } else {
+ CCTK_WARN (CCTK_WARN_ABORT, "internal error");
+ }
if (poison_boundaries) {
/* check destination grid points for poison */
diff --git a/src/rotatingsymmetry180.h b/src/rotatingsymmetry180.h
index e33770a..d006812 100644
--- a/src/rotatingsymmetry180.h
+++ b/src/rotatingsymmetry180.h
@@ -1,5 +1,3 @@
-/* $Header$ */
-
#ifndef ROTATINGSYMMETRY180_H
#define ROTATINGSYMMETRY180_H