aboutsummaryrefslogtreecommitdiff
path: root/src/apply.c
blob: 114c1f771a7df234fadf3a214bff901c028f61f1 (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
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
#include <assert.h>
#include <stdlib.h>
#include <string.h>

#include "cctk.h"
#include "cctk_Arguments.h"
#include "cctk_Parameters.h"

#include "util_ErrorCodes.h"
#include "util_Table.h"

#include "reflection.h"



static const char rcsid[] = "$Header$";
CCTK_FILEVERSION(AEIDevelopment_ReflectionSymmetry_apply_c);



#define COPY(VARTYPE)                                                   \
  static void                                                           \
  copy_##VARTYPE (VARTYPE const * restrict const srcvar,                \
                  VARTYPE * restrict const dstvar,                      \
                  int const ni, int const nj, int const nk,             \
                  int const imin, int const jmin, int const kmin,       \
                  int const imax, int const jmax, int const kmax,       \
                  int const ioff, int const joff, int const koff,       \
                  int const idir, int const jdir, int const kdir,       \
                  int const parity)                                     \
  {                                                                     \
    int i, j, k;                                                        \
                                                                        \
    for (k=kmin; k<kmax; ++k)                                           \
    {                                                                   \
      for (j=jmin; j<jmax; ++j)                                         \
      {                                                                 \
        for (i=imin; j<imax; ++i)                                       \
        {                                                               \
          int const dstind = i + ni * (j + nj * k);                     \
          int const ii = ioff + idir * i;                               \
          int const jj = joff + jdir * j;                               \
          int const kk = koff + kdir * k;                               \
          int const srcind = ii + ni * (jj + nj * kk);                  \
          REAL(dstvar[dstind] RE = parity * srcvar[srcind] RE;)         \
            IMAG(dstvar[dstind] IM = parity * srcvar[srcind] IM;)       \
            }                                                           \
      }                                                                 \
    }                                                                   \
  }

#define REAL(x) x
#define IMAG(x) /* nothing */
#define RE /* nothing */
#define IM /* nothing */

#ifdef CCTK_INT1
COPY(CCTK_INT1)
#endif

#ifdef CCTK_INT2
     COPY(CCTK_INT2)
#endif

#ifdef CCTK_INT4
     COPY(CCTK_INT4)
#endif

#ifdef CCTK_INT8
     COPY(CCTK_INT8)
#endif

#ifdef CCTK_REAL4
     COPY(CCTK_REAL4)
#endif

#ifdef CCTK_REAL8
     COPY(CCTK_REAL8)
#endif

#ifdef CCTK_REAL16
     COPY(CCTK_REAL16)
#endif

#undef REAL
#undef IMAG
#undef RE
#undef IM

#define REAL(x) x
#define IMAG(x) x
#define RE .Re
#define IM .Im

#ifdef CCTK_REAL4
     COPY(CCTK_COMPLEX8)
#endif

#ifdef CCTK_REAL8
     COPY(CCTK_COMPLEX16)
#endif

#ifdef CCTK_REAL16
     COPY(CCTK_COMPLEX32)
#endif

#undef REAL
#undef IMAG
#undef RE
#undef IM

#undef COPY



static int
BndReflectVI (cGH const * restrict const cctkGH,
              int const * restrict const stencil,
              int const vi)
{
  DECLARE_CCTK_PARAMETERS;
  
  int gi;
  cGroup group;
  cGroupDynamicData data;
  int firstvar, numvars;
  char * restrict fullname;
  char * restrict groupname;
  
  void * restrict varptr;
  
  int table;
  char tensortypealias[1000];
  enum tensortype { UNKNOWN, SCALAR, VECTOR, TENSOR };
  enum tensortype ttype;
  int tcomponent;
  
  int do_reflection[6];
  int do_stagger[6];
  
  int dir, face;
  
  int lsh[3], imin[3], imax[3], ioff[3], idir[3];
  
  int parity;
  
  int d;
  
  int ierr;
  
  
  
  /* Check arguments */
  if (! cctkGH)
  {
    CCTK_WARN (0, "Argument cctkGH is NULL");
  }
  if (! stencil)
  {
    CCTK_WARN (0, "Argument stencil is NULL");
  }
  if (vi < 0 || vi >= CCTK_NumVars())
  {
    CCTK_WARN (0, "Illegal variable index");
  }
  
  if (verbose)
  {
    fullname = CCTK_FullName (vi);
    if (! fullname)
    {
      CCTK_WARN (0, "Internal error in CCTK_FullName");
    }
    CCTK_VInfo (CCTK_THORNSTRING,
                "Applying reflection boundary conditions to \"%s\"",
                fullname);
    free (fullname);
  }
  
  
  
  /* Get and check group information */
  gi = CCTK_GroupIndexFromVarI (vi);
  if (gi < 0 || gi > CCTK_NumGroups())
  {
    CCTK_WARN (0, "Internal error in CCTK_GroupIndexFromVarI");
  }
  
  ierr = CCTK_GroupData (gi, &group);
  assert (!ierr);
  assert (group.grouptype == CCTK_GF);
  assert (group.disttype == CCTK_DISTRIB_DEFAULT);
  assert (group.stagtype == 0);
  
  firstvar = CCTK_FirstVarIndexI (gi);
  assert (firstvar>=0 && firstvar<CCTK_NumVars());
  numvars = CCTK_NumVarsInGroupI (gi);
  assert (numvars>=0);
  
  ierr = CCTK_GroupDynamicData (cctkGH, gi, &data);
  assert (!ierr);
  
  table = CCTK_GroupTagsTableI(gi);
  assert (table>=0);
  
  varptr = CCTK_VarDataPtrI (cctkGH, 0, vi);
  assert (varptr);
  
  
  
  /* Get and check tensor type information */
  ierr = Util_TableGetString
    (table, sizeof tensortypealias, tensortypealias, "tensortypealias");
  if (ierr == UTIL_ERROR_TABLE_NO_SUCH_KEY) {
    groupname = CCTK_GroupName(gi);
    assert (groupname);
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Tensor type alias not declared for group \"%s\" -- assuming a scalar",
                groupname);
    free (groupname);
    strcpy (tensortypealias, "scalar");
  } else if (ierr<0) {
    groupname = CCTK_GroupName(gi);
    assert (groupname);
    CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Error in tensor type alias declaration for group \"%s\"",
                groupname);
    free (groupname);
  }
  
  ttype = UNKNOWN;
  tcomponent = 0;
  if (CCTK_EQUALS (tensortypealias, "scalar"))
  {
    /* scalar */
    if (numvars != 1) {
      groupname = CCTK_GroupName(gi);
      assert (groupname);
      CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "Group \"%s\" has the tensor type alias \"scalar\", but contains more than 1 element",
                  groupname);
      free (groupname);
    }
    ttype = SCALAR;
    tcomponent = 0;
  }
  else if (CCTK_EQUALS (tensortypealias, "u"))
  {
    /* vector */
    assert (numvars == 3);
    ttype = VECTOR;
    tcomponent = vi - firstvar;
  } else if (CCTK_EQUALS (tensortypealias, "dd_sym")) {
    /* symmetric tensor */
    assert (numvars == 6);
    ttype = TENSOR;
    tcomponent = vi - firstvar;
  } else {
    char * groupname = CCTK_GroupName(gi);
    assert (groupname);
    CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Illegal tensor type alias for group \"%s\"",
                groupname);
    free (groupname);
  }
  
  switch (ttype)
  {
  case SCALAR:
    assert (tcomponent>=0 && tcomponent<1);
    break;
  case VECTOR:
    assert (tcomponent>=0 && tcomponent<3);
    break;
  case TENSOR:
    assert (tcomponent>=0 && tcomponent<6);
    break;
  default:
    assert (0);
  }
  
  
  
  /* Reflection symmetry information */
  do_reflection[0] = reflection_x;
  do_reflection[1] = reflection_upper_x;
  do_reflection[2] = reflection_y;
  do_reflection[3] = reflection_upper_y;
  do_reflection[4] = reflection_z;
  do_reflection[5] = reflection_upper_z;
  
  do_stagger[0] = avoid_origin_x;
  do_stagger[1] = avoid_origin_upper_x;
  do_stagger[2] = avoid_origin_y;
  do_stagger[3] = avoid_origin_upper_y;
  do_stagger[4] = avoid_origin_z;
  do_stagger[5] = avoid_origin_upper_z;
  
  
  
  /* Loop over all directions and faces */
  for (dir=0; dir<3; ++dir)
  {
    for (face=0; face<2; ++face)
    {
      /* If there is a reflection symmetry on that face */
      if (do_reflection[2*dir+face])
      {
        /* If we have the outer boundary of that face */
        if (cctkGH->cctk_bbox[2*dir+face])
        {
          
          /* Find parity */
          switch (ttype)
          {
          case SCALAR:
            parity = +1;
            break;
          case VECTOR:
            parity = dir == tcomponent ? -1 : +1;
            break;
          case TENSOR:
            switch (tcomponent)
            {
            case 0: parity = +1; break;
            case 1: parity = (dir == 0 || dir == 1) ? -1 : +1; break;
            case 2: parity = (dir == 0 || dir == 2) ? -1 : +1; break;
            case 3: parity = +1; break;
            case 4: parity = (dir == 1 || dir == 2) ? -1 : +1; break;
            case 5: parity = +1; break;
            default: assert (0);
            }
            break;
          default:
            assert (0);
          }
          
          /* Find region extent */
          for (d=0; d<3; ++d)
          {
            lsh[d] = cctkGH->cctk_lsh[d];
            imin[d] = 0;
            imax[d] = cctkGH->cctk_lsh[d];
            ioff[d] = 0;
          }
          if (face == 0)
          {
            imax[d] = cctkGH->cctk_nghostzones[d];
            ioff[d] = (+ cctkGH->cctk_nghostzones[d] 
                       + (do_stagger[2*d+face] ? 0 : 1));
          }
          else
          {
            imin[d] = cctkGH->cctk_lsh[d] - cctkGH->cctk_nghostzones[d];
            ioff[d] = (- cctkGH->cctk_nghostzones[d]
                       - (do_stagger[2*d+face] ? 0 : 1));
          }
          
          /* Copy region */
          switch (group.vartype)
          {
            
#ifdef CCTK_INT1
          case CCTK_VARIABLE_INT1:
#ifdef CCTK_INTER_PRECISION_1
          case CCTK_VARIABLE_INT:
#endif
            copy_CCTK_INT1 (varptr, varptr,
                            lsh[0], lsh[1], lsh[2],
                            imin[0], imin[1], imin[2],
                            imax[0], imax[1], imax[2],
                            ioff[0], ioff[1], ioff[2],
                            idir[0], idir[1], idir[2],
                            parity);
            break;
#endif
            
#ifdef CCTK_INT2
          case CCTK_VARIABLE_INT2:
#ifdef CCTK_INTER_PRECISION_2
          case CCTK_VARIABLE_INT:
#endif
            copy_CCTK_INT2 (varptr, varptr,
                            lsh[0], lsh[1], lsh[2],
                            imin[0], imin[1], imin[2],
                            imax[0], imax[1], imax[2],
                            ioff[0], ioff[1], ioff[2],
                            idir[0], idir[1], idir[2],
                            parity);
            break;
#endif
            
#ifdef CCTK_INT4
          case CCTK_VARIABLE_INT4:
#ifdef CCTK_INTER_PRECISION_4
          case CCTK_VARIABLE_INT:
#endif
            copy_CCTK_INT4 (varptr, varptr,
                            lsh[0], lsh[1], lsh[2],
                            imin[0], imin[1], imin[2],
                            imax[0], imax[1], imax[2],
                            ioff[0], ioff[1], ioff[2],
                            idir[0], idir[1], idir[2],
                            parity);
            break;
#endif
            
#ifdef CCTK_INT8
          case CCTK_VARIABLE_INT8:
#ifdef CCTK_INTER_PRECISION_8
          case CCTK_VARIABLE_INT:
#endif
            copy_CCTK_INT8 (varptr, varptr,
                            lsh[0], lsh[1], lsh[2],
                            imin[0], imin[1], imin[2],
                            imax[0], imax[1], imax[2],
                            ioff[0], ioff[1], ioff[2],
                            idir[0], idir[1], idir[2],
                            parity);
            break;
#endif
            
#ifdef CCTK_REAL4
          case CCTK_VARIABLE_REAL4:
#ifdef CCTK_REAL_PRECISION_4
          case CCTK_VARIABLE_REAL:
#endif
            copy_CCTK_REAL4 (varptr, varptr,
                             lsh[0], lsh[1], lsh[2],
                             imin[0], imin[1], imin[2],
                             imax[0], imax[1], imax[2],
                             ioff[0], ioff[1], ioff[2],
                             idir[0], idir[1], idir[2],
                             parity);
            break;
#endif
            
#ifdef CCTK_REAL8
          case CCTK_VARIABLE_REAL8:
#ifdef CCTK_REAL_PRECISION_8
          case CCTK_VARIABLE_REAL:
#endif
            copy_CCTK_REAL8 (varptr, varptr,
                             lsh[0], lsh[1], lsh[2],
                             imin[0], imin[1], imin[2],
                             imax[0], imax[1], imax[2],
                             ioff[0], ioff[1], ioff[2],
                             idir[0], idir[1], idir[2],
                             parity);
            break;
#endif
            
#ifdef CCTK_REAL16
          case CCTK_VARIABLE_REAL16:
#ifdef CCTK_REAL_PRECISION_16
          case CCTK_VARIABLE_REAL:
#endif
            copy_CCTK_REAL16 (varptr, varptr,
                              lsh[0], lsh[1], lsh[2],
                              imin[0], imin[1], imin[2],
                              imax[0], imax[1], imax[2],
                              ioff[0], ioff[1], ioff[2],
                              idir[0], idir[1], idir[2],
                              parity);
            break;
#endif
            
#ifdef CCTK_REAL4
          case CCTK_VARIABLE_COMPLEX8:
#ifdef CCTK_COMPLEX_PRECISION_8
          case CCTK_VARIABLE_COMPLEX:
#endif
            copy_CCTK_COMPLEX8 (varptr, varptr,
                                lsh[0], lsh[1], lsh[2],
                                imin[0], imin[1], imin[2],
                                imax[0], imax[1], imax[2],
                                ioff[0], ioff[1], ioff[2],
                                idir[0], idir[1], idir[2],
                                parity);
            break;
#endif
            
#ifdef CCTK_REAL8
          case CCTK_VARIABLE_COMPLEX16:
#ifdef CCTK_COMPLEX_PRECISION_16
          case CCTK_VARIABLE_COMPLEX:
#endif
            copy_CCTK_COMPLEX16 (varptr, varptr,
                                 lsh[0], lsh[1], lsh[2],
                                 imin[0], imin[1], imin[2],
                                 imax[0], imax[1], imax[2],
                                 ioff[0], ioff[1], ioff[2],
                                 idir[0], idir[1], idir[2],
                                 parity);
            break;
#endif
            
#ifdef CCTK_REAL16
          case CCTK_VARIABLE_COMPLEX32:
#ifdef CCTK_COMPLEX_PRECISION_32
          case CCTK_VARIABLE_COMPLEX:
#endif
            copy_CCTK_COMPLEX32 (varptr, varptr,
                                 lsh[0], lsh[1], lsh[2],
                                 imin[0], imin[1], imin[2],
                                 imax[0], imax[1], imax[2],
                                 ioff[0], ioff[1], ioff[2],
                                 idir[0], idir[1], idir[2],
                                 parity);
            break;
#endif
            
          default:
            CCTK_WARN (0, "Unsupported variable type");
          }
          
        } /* if cctk_bbox */
      } /* if do_reflection */
    } /* for face */
  } /* for dir */
  
  /* Success */
  return 0;
}



void
ReflectionSymmetry_Apply (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  
  int nvars;
  CCTK_INT * restrict indices;
  CCTK_INT * restrict faces;
  CCTK_INT * restrict widths;
  CCTK_INT * restrict tables;
  int vi;
  int dim;
  int * restrict stencil;
  int i;
  int istat;
  int ierr;
  
  if (! cctkGH)
  {
    CCTK_WARN (0, "Argument cctkGH is NULL");
  }
  
  nvars = Boundary_SelectedGVs (cctkGH, 0, NULL, NULL, NULL, NULL, NULL);
  if (nvars < 0)
  {
    CCTK_WARN (0, "Internal error in Boundary_SelectedGVs");
  }
  
  indices = malloc (nvars * sizeof *indices);
  if (! indices)
  {
    CCTK_WARN (0, "Out of memory");
  }
  faces = malloc (nvars * sizeof *faces);
  if (! faces)
  {
    CCTK_WARN (0, "Out of memory");
  }
  widths = malloc (nvars * sizeof *widths);
  if (! widths)
  {
    CCTK_WARN (0, "Out of memory");
  }
  tables = malloc (nvars * sizeof *tables);
  if (! tables)
  {
    CCTK_WARN (0, "Out of memory");
  }
  
  istat =  Boundary_SelectedGVs
    (cctkGH, nvars, indices, faces, widths, tables, 0);
  if (istat != nvars)
  {
    CCTK_WARN (0, "Internal error in Boundary_SelectedGVs");
  }
  
  for (i=0; i<nvars; ++i)
  {
    vi = indices[i];
    if (vi < 0 || vi >= CCTK_NumVars())
    {
      CCTK_WARN (0, "Illegal variable index");
    }
    
    if (widths[i] < 0)
    {
      CCTK_WARN (0, "Illegal boundary width");
    }
    
    dim = CCTK_GroupDimFromVarI (vi);
    if (dim < 0)
    {
      CCTK_WARN (0, "Illegal dimension");
    }
    
    stencil = malloc (dim * sizeof *stencil);
    if (! stencil)
    {
      CCTK_WARN (0, "Out of memory");
    }
    ierr = CCTK_GroupnghostzonesVI (cctkGH, dim, stencil, vi);
    if (ierr)
    {
      CCTK_WARN (0, "Internal error in CCTK_GroupnghostzonesVI");
    }
    
    ierr = BndReflectVI (cctkGH, stencil, vi);
    if (ierr)
    {
      CCTK_WARN (0, "Internal error in BndReflectVI");
    }
    
    free (stencil);
  }
  
  free (indices);
  free (faces);
  free (widths);
  free (tables);
}