aboutsummaryrefslogtreecommitdiff
path: root/src/SetBoundAxisym.c
blob: 37096e5675ae7b27752e48404da9fc2751d23aa6 (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
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
/* set_bound_axisym.c, BB */
/*@@
   @file      SetBoundAxisym.c
   @date      Sat Nov  6 16:35:46 MET 1999
   @author    Sai Iyer
   @desc 
              An implementation of Steve Brandt's idea for doing
              axisymmetry with 3d stencils.
              Cactus 4 version of Bernd Bruegmann's set_bound_axisym_GT.
   @enddesc
 @@*/

static char *rcsid = "$Id$"

#include "cctk.h"
#include "cctk_arguments.h"
#include "cctk_parameters.h"

#include "cartoon2D.h"

#include "math.h"

/* Number of components */
int tensor_type_n[N_TENSOR_TYPES] = {1, 3, 6};

/* set boundaries of a grid tensor assuming axisymmetry 
   - handles lower boundary in x
   - does not touch other boundaries
   - coordinates and rotation coefficients are independent of z and should
     be precomputed
   - this is also true for the constants in interpolator, but this may 
     be too messy
   - minimizes conceptual warpage, not computationally optimized
*/
/* uses rotation matrix and its inverse as linear transformation on
   arbitrary tensor indices -- I consider this a good compromise
   between doing index loops versus using explicit formulas in
   cos(phi) etc with messy signs 
*/
int SetBoundAxiSym(cGH *GH, int vi, int tensortype) {

{
  int i, j, k, di, dj, dk, ijk;
  int jj, n, s;
  Double x, y, rho;
  Double rxx, rxy, ryx, ryy;
  Double sxx, sxy, syx, syy;
  Double f[6], fx, fy, fz, fxx, fxy, fxz, fyy, fyz, fzz;
  Double *t, *tx, *ty, *tz, *txx, *txy, *txz, *tyy, *tyz, *tzz;

  if (0) printf("sba(%s)\n", GT->gf[0]->name);

  /* stencil width: in x-direction use ghost zone width,
                    in y-direction use ny/2, see below   */
  s = GH->stencil_width;

  /* Cactus loop in C:  grid points with y = 0 */
  /* compare thorn_BAM_Elliptic */
  /* strides used in stencils, should be provided by cactus */
  di = 1;
  dj = GH->lnx;
  dk = GH->lnx*GH->lny;

  /* y = 0 */
  j = GH->ny/2;

  /* z-direction: include lower and upper boundary */
  for (k = 0; k < GH->lnz; k++)

  /* y-direction: as many zombies as the evolution stencil needs */
  for (jj = 1, dj = jj*GH->lnx; jj <= GH->ny/2; jj++, dj = jj*GH->lnx) 

  /* x-direction: zombie for x < 0, including upper boundary for extrapol */
  for (i = -s; i < GH->lnx; i++) {

    /* index into linear memory */
    ijk = DI(GH,i,j,k);
   
    /* what a neat way to hack in the zombie for x < 0, y = 0 */
    if (i < 0) {i += s; ijk += s; dj = 0;}


    /* compute coordinates (could also use Cactus grid function) */
    x = GH->lx0 + GH->dx0 * i;
    y = (dj) ? GH->dy0 * jj : 0;
    rho = sqrt(x*x + y*y);

    /* compute rotation matrix
       note that this also works for x <= 0 (at lower boundary in x)
       note that rho is nonzero by definition if cactus checks dy ... 
    */
    rxx = x/rho;
    ryx = y/rho;
    rxy = -ryx;
    ryy = rxx;

    /* inverse rotation matrix, assuming detr = 1 */
    sxx = ryy;
    syx = -ryx;
    sxy = -rxy;
    syy = rxx;

    /* interpolate grid functions */
    for (n = 0; n < GT->n; n++)
      f[n] = axisym_interpolate(GH, GT->gf[n], i, di, ijk, rho); 

    /* rotate grid tensor by matrix multiplication */
    if (GT->tensor_type == TENSOR_TYPE_SCALAR) {
      t = GT->gf[0]->data;
      t[ijk+dj] = f[0];
      t[ijk-dj] = f[0];
    }
    else if (GT->tensor_type == TENSOR_TYPE_U) {
      tx = GT->gf[0]->data;
      ty = GT->gf[1]->data;
      tz = GT->gf[2]->data;
      fx = f[0];
      fy = f[1];
      fz = f[2];

      tx[ijk+dj] = rxx * fx + rxy * fy;
      ty[ijk+dj] = ryx * fx + ryy * fy;
      tz[ijk+dj] = fz;

      tx[ijk-dj] = sxx * fx + sxy * fy;
      ty[ijk-dj] = syx * fx + syy * fy;
      tz[ijk-dj] = fz;
    }
    else if (GT->tensor_type == TENSOR_TYPE_DDSYM) {
      txx = GT->gf[0]->data;
      txy = GT->gf[1]->data;
      txz = GT->gf[2]->data;
      tyy = GT->gf[3]->data;
      tyz = GT->gf[4]->data;
      tzz = GT->gf[5]->data;
      fxx = f[0];
      fxy = f[1];
      fxz = f[2];
      fyy = f[3];
      fyz = f[4];
      fzz = f[5];
      
      txx[ijk+dj] = fxx*sxx*sxx + 2*fxy*sxx*syx + fyy*syx*syx;
      tyy[ijk+dj] = fxx*sxy*sxy + 2*fxy*sxy*syy + fyy*syy*syy;
      txy[ijk+dj] = fxx*sxx*sxy + fxy*(sxy*syx+sxx*syy) + fyy*syx*syy;
      txz[ijk+dj] = fxz*sxx + fyz*syx;
      tyz[ijk+dj] = fxz*sxy + fyz*syy;
      tzz[ijk+dj] = fzz;

      txx[ijk-dj] = fxx*rxx*rxx + 2*fxy*rxx*ryx + fyy*ryx*ryx;
      tyy[ijk-dj] = fxx*rxy*rxy + 2*fxy*rxy*ryy + fyy*ryy*ryy;
      txy[ijk-dj] = fxx*rxx*rxy + fxy*(rxy*ryx+rxx*ryy) + fyy*ryx*ryy;
      txz[ijk-dj] = fxz*rxx + fyz*ryx;
      tyz[ijk-dj] = fxz*rxy + fyz*ryy;
      tzz[ijk-dj] = fzz;
    } 
    else {
      /* oops */
    }

    /* what a neat way to hack out the zombies for x < 0, y = 0 */
    if (dj == 0) {i -= s; ijk -= s; dj = jj*GH->lnx;}
  }

  /* syncs needed after interpolation (actually only for x direction) */
  for (n = 0; n < GT->n; n++)
    SyncGF(GH, GT->gf[n]); 
}





/* interpolation on x-axis */
Double axisym_interpolate(pGH *GH, pGF *GF, int i, int di, int ijk, Double x)
{
  Double *f = GF->data;
  Double c, d, x0;

  /* find i such that x(i) < x <= x(i+1)
     for rotation on entry always x > x(i), but sometimes also x > x(i+1) */
  while (x > GH->lx0 + GH->dx0 * (i+1)) {i++; ijk++;}

  /* the original demo */
  if (0) {
    c = (f[ijk+di] - f[ijk])/GH->dx0;
    d = f[ijk];
    x0 =  GH->lx0 + GH->dx0 * i; 
    return c*(x-x0) + d;
  }

  /* first attempt to interface to JT's interpolator */
  if (1) {
    double interpolate_local(int order, double x0, double dx, 
			     double y[], double x);
    double y[6], r;
    int n, offset;
    int order = cartoon_order;

    /* offset of stencil, note that rotation leads to x close to x0 
       for large x */
    offset = order/2;

    /* shift stencil at boundaries */
    /* note that for simplicity we don't distinguish between true boundaries
       and decomposition boundaries: the sync fixes that! */
    if (i-offset < 0) offset = i;
    if (i-offset+order >= GH->lnx) offset = i+order-GH->lnx+1;

    /* fill in info */
    /* fills in old data near axis, but as long as stencil at axis is 
       centered, no interpolation happens anyway */
    x0 = GH->lx0 + GH->dx0 * (i-offset);
    for (n = 0; n <= order; n++) {
      y[n] = f[ijk-offset+n];
      if (0) printf("y[%d]=%6.4f ", n, y[n]); 
    }

    /* call interpolator */
    r = interpolate_local(order, x0, GH->dx0, y, x);

    if (0) printf("i %d  x %7.4f  x0 %7.4f  r %7.4f\n", i, x, x0, r); 
    return r;
  }

  return 0;
}