aboutsummaryrefslogtreecommitdiff
diff options
context:
space:
mode:
authorlanfer <lanfer@c78560ca-4b45-4335-b268-5f3340f3cb52>2000-04-18 14:13:47 +0000
committerlanfer <lanfer@c78560ca-4b45-4335-b268-5f3340f3cb52>2000-04-18 14:13:47 +0000
commit06e6bac2be0e5d99c061016a637ff8fad5cbf162 (patch)
tree67cec58146be19f959d25ade204fe926b27efc4f
parent3b2df5e16b6a1412db8f8435551337e412d5769b (diff)
symmetry now works with no_origin
git-svn-id: http://svn.cactuscode.org/arrangements/CactusBase/CartGrid3D/trunk@102 c78560ca-4b45-4335-b268-5f3340f3cb52
-rw-r--r--src/Symmetry.c104
-rw-r--r--src/SymmetryWrappers.c49
2 files changed, 137 insertions, 16 deletions
diff --git a/src/Symmetry.c b/src/Symmetry.c
index 0b0209c..23e2f88 100644
--- a/src/Symmetry.c
+++ b/src/Symmetry.c
@@ -1,4 +1,13 @@
+ /*@@
+ @file Symmetry.c
+ @date Tue Apr 18 14:14:16 2000
+ @author Gerd Lanfermann
+ @desc
+ Routines to apply the 1/2/3D Symmetries for
+ all symmetry domains (octant/bitant/quadrant).
+ @enddesc
+ @@*/
#include <stdio.h>
#include <assert.h>
@@ -9,12 +18,40 @@
/*#define SYM_DEBUG*/
-/* index convention:
+/*@@
+ @routine CartApplySym3Di
+ @date Tue Apr 18 14:17:23 2000
+ @author Gerd Lanfermann
+ @desc Apply Symmetry BC to 3D variables
+
+ Variables passed through:
+ cGH *GH pointer to cGH
+ int gdim dimension of the variable
+ int *doSym flags whether to apply a symmetries on a given face
+ size 2*dim, here we only check for lower faces:0,2,4
+ int *cntstag value used when the gridpoints are staggered
+ around the origin
+ int *lssh size of the domain,
+ int *ghostz size of the ghostzone
+ int *sym symmetry values
+ CCTK_REAL *var pointer to variable
+
+ index convention:
i ~ x ~ 0
j ~ y ~ 1
- k ~ z ~ 2 */
+ k ~ z ~ 2
-int CartApplySym3Di(cGH *GH, int gdim, int *doSym,
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+
+int CartApplySym3Di(cGH *GH, int gdim, int *doSym, int *cntstag,
int *lssh, int *ghostz, int *sym, CCTK_REAL *var)
{
@@ -28,6 +65,7 @@ int CartApplySym3Di(cGH *GH, int gdim, int *doSym,
printf(" lssh: %d %d %d sym: %d %d %d \n",
lssh[0],lssh[1],lssh[2], sym[0], sym[2], sym[4] );
printf(" ghostz %d %d %d \n",ghostz[0],ghostz[1],ghostz[2]);
+ printf(" cntstag: %d %d %d\n",cntstag[0],cntstag[1],cntstag[2]);
#endif
if (doSym[0] == 1)
@@ -39,7 +77,7 @@ int CartApplySym3Di(cGH *GH, int gdim, int *doSym,
for(i=0; i < ghostz[0]; i++)
{
var[CCTK_GFINDEX3D(GH,i,j,k)] =
- sym[0]*var[CCTK_GFINDEX3D(GH,2*ghostz[0]-1-i,j,k)];
+ sym[0]*var[CCTK_GFINDEX3D(GH,2*ghostz[0]-cntstag[0]-i,j,k)];
}
}
}
@@ -53,7 +91,7 @@ int CartApplySym3Di(cGH *GH, int gdim, int *doSym,
for(j=0; j < ghostz[1]; j++)
{
var[CCTK_GFINDEX3D(GH,i,j,k)] =
- sym[2]*var[CCTK_GFINDEX3D(GH,i,2*ghostz[1]-1-j,k)];
+ sym[2]*var[CCTK_GFINDEX3D(GH,i,2*ghostz[1]-cntstag[1]-j,k)];
}
}
}
@@ -67,15 +105,37 @@ int CartApplySym3Di(cGH *GH, int gdim, int *doSym,
for(k=0; k < ghostz[2]; k++)
{
var[CCTK_GFINDEX3D(GH,i,j,k)] =
- sym[4]*var[CCTK_GFINDEX3D(GH,i,j,2*ghostz[2]-1-k)];
+ sym[4]*var[CCTK_GFINDEX3D(GH,i,j,2*ghostz[2]-cntstag[2]-k)];
}
}
}
}
return(0);
}
+
+
+/*@@
+ @routine CartApplySym2Di
+ @date Tue Apr 18 14:17:23 2000
+ @author Gerd Lanfermann
+ @desc Apply Symmetry BC to 2D variables
+
+
+ index convention:
+ i ~ x ~ 0
+ j ~ y ~ 1
+ k ~ z ~ 2
-int CartApplySym2Di(cGH *GH, int gdim, int *doSym,
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+int CartApplySym2Di(cGH *GH, int gdim, int *doSym, int *cntstag,
int *lssh, int *ghostz, int *sym, CCTK_REAL *var)
{
int i,j;
@@ -87,7 +147,7 @@ int CartApplySym2Di(cGH *GH, int gdim, int *doSym,
for(i=0; i < ghostz[0]; i++)
{
var[CCTK_GFINDEX2D(GH,i,j)] =
- sym[0]*var[CCTK_GFINDEX2D(GH,2*ghostz[0]-1-i,j)];
+ sym[0]*var[CCTK_GFINDEX2D(GH,2*ghostz[0]-cntstag[0]-i,j)];
}
}
}
@@ -99,7 +159,7 @@ int CartApplySym2Di(cGH *GH, int gdim, int *doSym,
for(j=0; j < ghostz[1]; j++)
{
var[CCTK_GFINDEX2D(GH,i,j)] =
- sym[2]*var[CCTK_GFINDEX2D(GH,i,2*ghostz[1]-1-j)];
+ sym[2]*var[CCTK_GFINDEX2D(GH,i,2*ghostz[1]-cntstag[1]-j)];
}
}
}
@@ -107,7 +167,29 @@ int CartApplySym2Di(cGH *GH, int gdim, int *doSym,
return(0);
}
-int CartApplySym1Di(cGH *GH, int gdim, int *doSym,
+
+/*@@
+ @routine CartApplySym1Di
+ @date Tue Apr 18 14:17:23 2000
+ @author Gerd Lanfermann
+ @desc Apply Symmetry BC to 1D variables
+
+
+ index convention:
+ i ~ x ~ 0
+ j ~ y ~ 1
+ k ~ z ~ 2
+
+ @enddesc
+ @calls
+ @calledby
+ @history
+
+ @endhistory
+
+@@*/
+
+int CartApplySym1Di(cGH *GH, int gdim, int *doSym, int *cntstag,
int *lssh, int *ghostz, int *sym, CCTK_REAL *var)
{
int i;
@@ -117,7 +199,7 @@ int CartApplySym1Di(cGH *GH, int gdim, int *doSym,
for(i=0; i < ghostz[0]; i++)
{
var[CCTK_GFINDEX1D(GH,i)] =
- sym[0]*var[CCTK_GFINDEX1D(GH,2*ghostz[0]-1-i)];
+ sym[0]*var[CCTK_GFINDEX1D(GH,2*ghostz[0]-cntstag[0]-i)];
}
}
diff --git a/src/SymmetryWrappers.c b/src/SymmetryWrappers.c
index 92247cd..18e4ad0 100644
--- a/src/SymmetryWrappers.c
+++ b/src/SymmetryWrappers.c
@@ -13,6 +13,7 @@
int CartApplySym3Di(cGH *GH,
int gdim,
int *doSym,
+ int *cntstag,
int *lssh,
int *ghostz,
int *sym,
@@ -21,6 +22,7 @@ int CartApplySym3Di(cGH *GH,
int CartApplySym2Di(cGH *GH,
int gdim,
int *doSym,
+ int *cntstag,
int *lssh,
int *ghostz,
int *sym,
@@ -29,6 +31,7 @@ int CartApplySym2Di(cGH *GH,
int CartApplySym1Di(cGH *GH,
int gdim,
int *doSym,
+ int *cntstag,
int *lssh,
int *ghostz,
int *sym,
@@ -45,7 +48,7 @@ int CartSymGI(cGH *GH, int gi)
int idim, gdim;
int berr=-1,ierr=-1;
int time;
- int *doSym, *dstag, *lssh;
+ int *doSym, *dstag, *lssh, *cntstag;
SymmetryGHex *sGHex;
/* Get out if we are sure no symmetries should be applied */
@@ -68,9 +71,24 @@ int CartSymGI(cGH *GH, int gi)
doSym = (int *)malloc((2*gdim)*sizeof(int));
dstag = (int *)malloc(gdim*sizeof(int));
lssh = (int *)malloc(gdim*sizeof(int));
-
+ cntstag= (int *)malloc(gdim*sizeof(int));
+
/* get the directional staggering of the group */
ierr = CCTK_GroupStaggerDirArrayGI(dstag, gdim, gi);
+
+ /* Set value to one if grid is staggered around the center */
+ if (no_origin)
+ {
+ cntstag[0]=1;
+ cntstag[1]=1;
+ cntstag[2]=1;
+ }
+ else
+ {
+ cntstag[0]=0;
+ cntstag[1]=0;
+ cntstag[2]=0;
+ }
/* Use next time level, if present */
time = CCTK_NumTimeLevelsFromVarI(first_vi) - 1;
@@ -126,7 +144,8 @@ int CartSymGI(cGH *GH, int gi)
{
case 1: berr = CartApplySym1Di(GH,
gdim,
- doSym,
+ doSym,
+ cntstag,
lssh,
GH->cctk_nghostzones,
sGHex->GFSym[vi],
@@ -134,13 +153,15 @@ int CartSymGI(cGH *GH, int gi)
case 2: berr = CartApplySym2Di(GH,
gdim,
doSym,
+ cntstag,
lssh,
GH->cctk_nghostzones,
sGHex->GFSym[vi],
GH->data[vi][time]); break;
case 3: berr = CartApplySym3Di(GH,
gdim,
- doSym,
+ doSym,
+ cntstag,
lssh,
GH->cctk_nghostzones,
sGHex->GFSym[vi],
@@ -207,7 +228,7 @@ int CartSymVI(cGH *GH, int vi)
int idim, gdim;
int berr=-1, ierr=-1;
int time;
- int *doSym, *dstag, *lssh;
+ int *doSym, *dstag, *lssh, *cntstag;
SymmetryGHex *sGHex;
/* Get out if we are sure no symmetries should be applied */
@@ -231,9 +252,23 @@ int CartSymVI(cGH *GH, int vi)
doSym = (int *)malloc((2*gdim)*sizeof(int));
dstag = (int *)malloc(gdim*sizeof(int));
lssh = (int *)malloc(gdim*sizeof(int));
+ cntstag= (int *)malloc(gdim*sizeof(int));
/* get the directional staggering of the group */
ierr = CCTK_GroupStaggerDirArrayGI(dstag, gdim, gi);
+
+ if (no_origin)
+ {
+ cntstag[0]=1;
+ cntstag[1]=1;
+ cntstag[2]=1;
+ }
+ else
+ {
+ cntstag[0]=0;
+ cntstag[1]=0;
+ cntstag[2]=0;
+ }
/* Use next time level, if present */
time = CCTK_NumTimeLevelsFromVarI(vi) - 1;
@@ -289,6 +324,7 @@ int CartSymVI(cGH *GH, int vi)
case 1: berr = CartApplySym1Di(GH,
gdim,
doSym,
+ cntstag,
lssh,
GH->cctk_nghostzones,
sGHex->GFSym[vi],
@@ -296,6 +332,7 @@ int CartSymVI(cGH *GH, int vi)
case 2: berr = CartApplySym2Di(GH,
gdim,
doSym,
+ cntstag,
lssh,
GH->cctk_nghostzones,
sGHex->GFSym[vi],
@@ -303,6 +340,7 @@ int CartSymVI(cGH *GH, int vi)
case 3: berr = CartApplySym3Di(GH,
gdim,
doSym,
+ cntstag,
lssh,
GH->cctk_nghostzones,
sGHex->GFSym[vi],
@@ -312,6 +350,7 @@ int CartSymVI(cGH *GH, int vi)
free(dstag);
free(doSym);
+ free(cntstag);
return(berr);
}