aboutsummaryrefslogtreecommitdiff
path: root/src/SymmetryCondition.c
blob: c6c8ce6e6fcfca0d6c489ec7bab9874f6d1cb2f0 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
 /*@@
   @file      SymmetryCondition.c
   @date      Thu Oct  7 16:45:19 1999
   @author    Tom Goodale
   @desc 
   C version of Gerd's symmetry stuff.
   @enddesc 
 @@*/

#include "cctk.h"

static const char *rcsid = "$Header$";

CCTK_FILEVERSION(CactusBase_CartGrid3D_SymmetryCondition_c)


/*@@
   @routine    SymmetryCondition
   @date       Mon Mar 15 15:51:57 1999
   @author     Gerd Lanfermann
   @desc 
               Routine performs the symmetry boundary operations.
   @enddesc 
   @calls     
   @calledby   
   @history 
   @hdate Thu Oct  7 16:47:35 1999 @hauthor Tom Goodale
   @hdesc Converted to C 
   @endhistory 
@@*/


#define GFINDEX3D(sh,i,j,k) ((i) + sh[0]*((j)+sh[1]*(k)))

void SymmetryCondition(int nxyz[],CCTK_REAL var[], int nghostzones,int sym[], int doSym[])
{
  int i,j,k;

  int sw;

  /*      
   *      Apply symmetry if 
   *       * the grid chunk has a physical boundary (bbox)
   *       * its size in a direction is bigger than one (sh)
   *       * we actually want a symmetry (sx.ne.ESYM_UNSET)
   */
  if (doSym[1] == 1 || doSym[3] == 1 || doSym[5] == 1)
  {
    CCTKi_NotYetImplemented("Right hand side boundary conditions");
  }

  if (doSym[0] == 1)
  {
    for(k=0; k < nxyz[2]; k++)
    {
      for(j=0; j < nxyz[1]; j++)
      {
        for(sw=0; sw < nghostzones; sw++)
        {
          var[GFINDEX3D(nxyz,sw,j,k)] = 
	    sym[0]*var[GFINDEX3D(nxyz,2*nghostzones-1-sw,j,k)];
        }
      }
    }
  }

  if (doSym[2] == 1)
  {
    for(k=0; k < nxyz[2]; k++)
    {
      for(sw=0; sw < nghostzones; sw++)
      {
        for(i=0; i < nxyz[0]; i++)
        {
          var[GFINDEX3D(nxyz,i,sw,k)] = 
	    sym[2]*var[GFINDEX3D(nxyz,i,2*nghostzones-1-sw,k)];
        }
      }
    }
  }

  if (doSym[4] == 1)
  {
    for(sw=0; sw < nghostzones; sw++)
    {
      for(j=0; j < nxyz[1]; j++)
      {
        for(i=0; i < nxyz[0]; i++)
        {
          var[GFINDEX3D(nxyz,i,j,sw)] = 
	    sym[4]*var[GFINDEX3D(nxyz,i,j,2*nghostzones-1-sw)];
        }
      }
    }
  }

  return;
}

void CCTK_FNAME(SymmetryCondition)(int nxyz[],CCTK_REAL var[], int *nghostzones,int sym[], int doSym[])
{
  SymmetryCondition(nxyz, var, *nghostzones, sym, doSym);
}