aboutsummaryrefslogtreecommitdiff
path: root/src/Cartoon2DBC.c
blob: 9fcbfc4c1de668680b3426a4cef0fe2e41ca9a56 (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
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
/*@@
   @file      Cartoon2DBC.c
   @date      Thu Nov  4 13:35:00 MET 1999
   @author    Sai Iyer
   @desc 
              Apply Cartoon2D boundary conditions
              An implementation of Steve Brandt's idea for doing
              axisymmetry with 3d stencils.
              Cactus 4 version of Bernd Bruegmann's code.             
   @enddesc
 @@*/

#include <stdio.h>
#include <assert.h>
#include <stdlib.h>
#include <ctype.h>
#include <stdarg.h>
#include <string.h>
#include <math.h>

#include "cctk.h"
#include "cctk_Parameters.h"
#include "cctk_FortranString.h"

#include "Cartoon2D_tensors.h"

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

CCTK_FILEVERSION(Development_Cartoon2D_Cartoon2DBC_c)

CCTK_INT BndCartoon2DVN(cGH *GH, CCTK_INT tensortype, const char *var);
CCTK_INT BndCartoon2DVI(cGH *GH, CCTK_INT tensortype, CCTK_INT vi);
CCTK_REAL Cartoon2DInterp(cGH *GH, CCTK_REAL *f, CCTK_INT i, CCTK_INT di,
                             CCTK_INT ijk, CCTK_REAL x);
CCTK_REAL interpolate_local(CCTK_INT order, CCTK_REAL x0, CCTK_REAL dx, 
			 CCTK_REAL y[], CCTK_REAL x);
void CCTK_FCALL CCTK_FNAME(BndCartoon2DVI)(CCTK_INT *retval, cGH *GH, 
                CCTK_INT *tensortype, CCTK_INT *vi);
void CCTK_FCALL CCTK_FNAME(BndCartoon2DVN)(CCTK_INT *retval, cGH *GH,
                CCTK_INT *tensortype, ONE_FORTSTRING_ARG);


/* 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 
*/

CCTK_INT BndCartoon2DVI(cGH *GH, CCTK_INT tensortype, CCTK_INT vi)

{
  DECLARE_CCTK_PARAMETERS

  static int xindx=-1;

  CCTK_INT i, j, k, di, dj, ijk;
  CCTK_INT jj, n, s;
  CCTK_INT lnx, lny, lnz, ny;
  CCTK_INT n_tensortype;
  CCTK_INT timelevel;
  CCTK_REAL *xp;
  CCTK_REAL x, y, rho;
  CCTK_REAL lx0, dx0, dy0;
  CCTK_REAL rxx, rxy, ryx, ryy;
  CCTK_REAL sxx, sxy, syx, syy;
  CCTK_REAL f[6], fx, fy, fz, fxx, fxy, fxz, fyy, fyz, fzz;
  CCTK_REAL *t, *tx, *ty, *tz, *txx, *txy, *txz, *tyy, *tyz, *tzz;


  switch (tensortype) {
  case TENSORTYPE_SCALAR:
    n_tensortype = N_TENSORTYPE_SCALAR;
    break;
  case TENSORTYPE_U:
    n_tensortype = N_TENSORTYPE_U;
    break;
  case TENSORTYPE_DDSYM:
    n_tensortype = N_TENSORTYPE_DDSYM;
    break;
  default:
    CCTK_WARN(0,"Tensor type not recognized by Cartoon2D.");
  }
    
  s = GH->cctk_nghostzones[0];
  lnx = GH->cctk_lsh[0];
  lny = GH->cctk_lsh[1];
  lnz = GH->cctk_lsh[2];
  ny  = GH->cctk_gsh[1];

  xindx = (xindx<0)?CCTK_CoordIndex(-1,"x","cart3d"):xindx;
  xp    = GH->data[xindx][0];

/*  xp  = GH->data[CCTK_CoordIndex(-1,"x","cart3d")][0];*/
  lx0 = xp[0];
  dx0 = GH->cctk_delta_space[0]/GH->cctk_levfac[0];
  dy0 = GH->cctk_delta_space[1]/GH->cctk_levfac[1];

  timelevel = 0; 

  /* 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 = lnx;

  /* y = 0 */
  j = ny/2;

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

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

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

    /* index into linear memory */
    ijk = CCTK_GFINDEX3D(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 = lx0 + dx0 * i;
    y = (dj) ? 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 < n_tensortype; n++)
      f[n] = Cartoon2DInterp(GH, GH->data[vi+n][timelevel], i, di, ijk, rho); 

    /* rotate grid tensor by matrix multiplication */
    if (tensortype == TENSORTYPE_SCALAR) {
      t = GH->data[vi][timelevel];
      t[ijk+dj] = f[0];
      t[ijk-dj] = f[0];
    }
    else if (tensortype == TENSORTYPE_U) {
      tx = GH->data[vi][timelevel];
      ty = GH->data[vi+1][timelevel];
      tz = GH->data[vi+2][timelevel];
      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 (tensortype == TENSORTYPE_DDSYM) {
      txx = GH->data[vi][timelevel];
      txy = GH->data[vi+1][timelevel];
      txz = GH->data[vi+2][timelevel];
      tyy = GH->data[vi+3][timelevel];
      tyz = GH->data[vi+4][timelevel];
      tzz = GH->data[vi+5][timelevel];
      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 {
      return(-1);
    }

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

  /* syncs needed after interpolation (actually only for x direction) */
    CCTK_SyncGroup(GH, CCTK_GroupNameFromVarI(vi)); 

  return(0);
}

void CCTK_FCALL CCTK_FNAME(BndCartoon2DVI)(CCTK_INT *retval, cGH *GH, CCTK_INT *tensortype, CCTK_INT *vi) {
  *retval = BndCartoon2DVI(GH, *tensortype, *vi);
}

CCTK_INT BndCartoon2DVN(cGH *GH, CCTK_INT tensortype, const char *impvarname)
{
  DECLARE_CCTK_PARAMETERS

  CCTK_INT vi;

  if (verbose) printf("cartoon2D called for %s\n", impvarname);
  vi = CCTK_VarIndex(impvarname);

  return(BndCartoon2DVI(GH, tensortype, vi));
}

void CCTK_FCALL CCTK_FNAME(BndCartoon2DVN)(CCTK_INT *retval, cGH *GH, CCTK_INT *tensortype, ONE_FORTSTRING_ARG) {
  ONE_FORTSTRING_CREATE(impvarname)
  *retval = BndCartoon2DVN(GH, *tensortype, impvarname);
  free(impvarname);
}

CCTK_INT BndCartoon2DGN(cGH *GH, CCTK_INT tensortype, const char *group)
{
  DECLARE_CCTK_PARAMETERS

  CCTK_INT n_tensortype;

  if (verbose) printf("cartoon2D called for %s\n", group);

  switch (tensortype) {
  case TENSORTYPE_SCALAR:
    n_tensortype = N_TENSORTYPE_SCALAR;
    break;
  case TENSORTYPE_U:
    n_tensortype = N_TENSORTYPE_U;
    break;
  case TENSORTYPE_DDSYM:
    n_tensortype = N_TENSORTYPE_DDSYM;
    break;
  default:
    CCTK_WARN(0,"Tensor type not recognized by Cartoon2D.");
  }

  if(CCTK_NumVarsInGroup(group) != n_tensortype){
    char *msg;
    msg = malloc(1024*sizeof(char));
    sprintf(msg, "%s should have only %d component(s).\n", 
	    group, n_tensortype);
    CCTK_WARN(0,msg);
    free(msg);
  }
  
  return(BndCartoon2DVI(GH, tensortype, CCTK_FirstVarIndex(group)));
}
void CCTK_FCALL CCTK_FNAME(BndCartoon2DGN)(CCTK_INT *retval, cGH *GH, CCTK_INT *tensortype, ONE_FORTSTRING_ARG) {
  ONE_FORTSTRING_CREATE(group)
  *retval = BndCartoon2DGN(GH, *tensortype, group);
  free(group);
}

/* interpolation on x-axis */
CCTK_REAL Cartoon2DInterp(cGH *GH, CCTK_REAL *f, CCTK_INT i, CCTK_INT di, 
              CCTK_INT ijk, CCTK_REAL x)
{
  DECLARE_CCTK_PARAMETERS

  static int xindx=-1;

  CCTK_REAL x0;
  CCTK_REAL lx0, dx0;
  CCTK_REAL y[6], r;
  CCTK_REAL *xp;
  CCTK_INT lnx;
  CCTK_INT n, offset;

  lnx = GH->cctk_lsh[0];

  xindx = (xindx<0)?CCTK_CoordIndex(-1,"x","cart3d"):xindx;
  xp    = GH->data[xindx][0];


/*  xp = GH->data[CCTK_CoordIndex(-1,"x","cart3d")][0];*/
  lx0 = xp[0];
  dx0 = GH->cctk_delta_space[0]/GH->cctk_levfac[0];

  /* 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 > lx0 + dx0 * (i+1)) {i++; ijk++;}

  /* first attempt to interface to JT's interpolator */

  /* 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 >= lnx) offset = i+order-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 = lx0 + dx0 * (i-offset);
  for (n = 0; n <= order; n++) {
    y[n] = f[ijk-offset+n];
  }
  
  /* call interpolator */
  r = interpolate_local(order, x0, dx0, y, x);
  
  return r;
}