aboutsummaryrefslogtreecommitdiff
path: root/src/RadiationBoundary.c
blob: 479c8bf534d642fa6e26b67e0d09cc17c0dc1350 (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
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
 /*@@
   @file      RadiationBoundary.c
   @date      Mon Mar 15 15:09:00 1999
   @author    Miguel Alcubierre, Gabrielle Allen, Gerd Lanfermann
   @desc
              Routines for applying radiation boundary conditions

              The radiative boundary condition that is implemented is:

                f  =  f0  +  u(r - v*t) / r  +  h(r + v*t) / r

              That is, I assume outgoing radial waves with a 1/r
              fall off, and the correct asymptotic value f0, plus
              I include the possibility of incoming waves as well
              (these incoming waves should be modeled somehow).

              The condition above leads to the differential equation:

                (x / r) d f  +  v d f  + v x (f - f0) / r^2  =  v x H / r^2
                  i      t         i        i                      i

              where x_i is the normal direction to the given boundaries,
              and H = 2 dh(s)/ds.

              So at a given boundary I only worry about derivatives in
              the normal direction.  Notice that u(r-v*t) has dissapeared,
              but we still do not know the value of H.

              To get H what I do is the following:  I evaluate the
              expression one point in from the boundary and solve for H
              there.  We now need a way of extrapolation H to the boundary.
              For this I assume that H falls off as a power law:

                H = k/r**n   =>  d H  =  - n H/r
                                  i

              The value of n is is defined by the parameter "radpower".
              If this parameter is negative, H is forced to be zero (this
              corresponds to pure outgoing waves and is the default).

              The behaviour I have observed is the following:  Using H=0
              is very stable, but has a very bad initial transient. Taking
              n to be 0 or positive improves the initial behaviour considerably,
              but introduces a drift that can kill the evolution at very late
              times.  Empirically, the best value I have found is n=2, for
              which the initial behaviour is very nice, and the late time drift
              is quite small.

              Another problem with this condition is that it does not
              use the physical characteristic speed, but rather it assumes
              a wave speed of v, so the boundaries should be out in
              the region where the characteristic speed is constant.
              Notice that this speed does not have to be 1.  For gauge
              quantities {alpha,phi,trK} we can have a different asymptotic
              speed, which is why the value of v is passed as a parameter.
   @enddesc
   @history
   @hdate     unknown
   @hauthor   Gerd Lanfermann
   @hdesc     Ported to Cactus 4.0
   @hdate     Fri 6 Apr 2001
   @hauthor   Thomas Radke
   @hdesc     BC routines generalized for applying to arbitrary CCTK data types
   @endhistory
   @version   $Id$
 @@*/

#include <stdio.h>
#include <stdlib.h>
#include <string.h>

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

#include "Boundary.h"

/* #define DEBUG */

/* the rcs ID and its dummy function to use it */
static const char *rcsid = "$Header$";
CCTK_FILEVERSION(CactusBase_Boundary_RadiationBoundary_c);

static int ApplyBndRadiative (const cGH *GH,
                              int stencil_dir,
                              const CCTK_INT *stencil_alldirs,
                              int dir,
                              CCTK_REAL var0,
                              CCTK_REAL speed,
                              CCTK_INT first_var_to,
                              CCTK_INT first_var_from,
                              int num_vars);
static int OldApplyBndRadiative (const cGH *GH,
                                 int stencil_dir,
                                 const CCTK_INT *stencil_alldirs,
                                 int dir,
                                 CCTK_REAL var0,
                                 CCTK_REAL speed,
                                 int first_var_to,
                                 int first_var_from,
                                 int num_vars);

/********************************************************************
 ********************    External Routines   ************************
 ********************************************************************/

/*@@
   @routine    BndRadiative
   @date       6 Nov 2002
   @author     David Rideout
   @desc
               Top level function which is registered as handling
               the Radiative boundary condition
   @enddesc
   @calls      ApplyBndRadiative

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        num_vars
   @vdesc      number of variables passed in through var_indices[]
   @vtype      CCTK_INT
   @vio        in
   @endvar
   @var        var_indices
   @vdesc      array of variable indicies to which to apply this boundary
               condition
   @vtype      CCTK_INT *
   @vio        in
   @endvar
   @var        faces
   @vdesc      array of set of faces to which to apply the bc
   @vtype      CCTK_INT
   @vio        in
   @endvar
   @var        widths
   @vdesc      array of boundary widths for each variable
   @vtype      CCTK_INT
   @vio        in
   @endvar
   @var        table_handles
   @vdesc      array of table handles which hold extra arguments
   @vtype      CCTK_INT
   @vio        in
   @endvar
   @returntype CCTK_INT
   @returndesc
               return code of @seeroutine ApplyBndRadiative
               -21 error reading boundary width array from table
               -22 wrong size boundary width array in table
   @endreturndesc
@@*/

CCTK_INT BndRadiative(const cGH *GH, CCTK_INT num_vars, CCTK_INT *vars,
                      CCTK_INT *faces, CCTK_INT *widths, CCTK_INT *tables)
{
  int i, j, k, gi, gdim, max_gdim, err, retval;
  CCTK_INT value_type, value_size;
  char *prev_time_level_name;

  /* variables to pass to ApplyBndRadiative */
  CCTK_INT *width_alldirs;  /* width of boundary in all directions */
  int dir;             /* direction in which to apply bc */
  CCTK_REAL limit, speed;
  CCTK_INT prev_time_level; /* variable index which holds the previous time level */

#ifdef DEBUG
  printf("BndRadiative() got passed: GH=%p, num_vars=%d:\n", (const void *) GH,
         num_vars);
  printf("var index  var name  table handle\n");
  for (i=0; i<num_vars; ++i)
  {
    printf("%d  %12s  %d\n", vars[i], CCTK_VarName(vars[i]), tables[i]);
  }
  printf("end of table\n");

  /*  CCTK_WARN(0, "stopping code"); */
#endif

  /* Initialize variables */
  retval = 0; width_alldirs = NULL; max_gdim = 0;

  /* loop through variables, j at a time */
  for (i=0; i<num_vars; i+=j) {
    /* find other adjacent vars which are selected for identical bcs */
    j=1;
    /* Since GFs are allowed to have different staggering, the best we
       can do is find variables of the same group which are selected
       for identical bcs.  If all GFs had the same staggering then we
       could groups many GFs together. */
    gi = CCTK_GroupIndexFromVarI (vars[i]);
#ifdef DEBUG
    printf("starting increment computation with group %d:\n", gi);
    printf("this group had %d members\n", CCTK_NumVarsInGroupI(gi));
#endif
    while (i+j<num_vars && vars[i+j]==vars[i]+j &&
           CCTK_GroupIndexFromVarI(vars[i+j])==gi && tables[i+j]==tables[i]
           && faces[i+j]==faces[i] && widths[i+j]==widths[i])
    {
      ++j;
    }

    /* Check to see if faces specification is valid */
    if (faces[i] != CCTK_ALL_FACES)
    {
      CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
                 "Faces specification %d for Radiative boundary conditions on "
                 "%s is not implemented yet.  "
                 "Applying Radiative bcs to all (external) faces.", faces[i],
                 CCTK_VarName(vars[i]));
    }
    dir = 0; /* apply bc to all faces */

    /* Set up default arguments for ApplyBndRadiative */
    /* Defaults for remainder of arguments */
    limit = 0.;
    prev_time_level = vars[i];
    speed = 1.;

    /* Look on table for possible non-default arguments
     * (If any of these table look-ups fail, the value will be unchanged
     * from its default value)
     */
    /* Asymptotic value of function at infinity */
    err = Util_TableGetReal(tables[i], &limit, "LIMIT");
    if (err == UTIL_ERROR_BAD_HANDLE)
    {
      CCTK_VWarn(5, __LINE__, __FILE__, CCTK_THORNSTRING,
                 "Invalid table handle passed for Radiative boundary "
                 "conditions for %s.  Using all default values.",
                 CCTK_VarName(vars[i]));
    } else
    {
      /* Previous time level (for GVs which don't use Cactus time levels)
         (to be deprecated) */
      if (Util_TableQueryValueInfo(tables[i], &value_type, &value_size,
                                   "PREVIOUS_TIME_LEVEL"))
      {
        /* Key for previous time level is set */
        if (value_type==CCTK_VARIABLE_STRING)
        {
          prev_time_level_name = malloc(value_size*sizeof(char));
          Util_TableGetString(tables[i], value_size, prev_time_level_name,
                              "PREVIOUS_TIME_LEVEL");
          prev_time_level = (CCTK_INT) CCTK_VarIndex(prev_time_level_name);
          free(prev_time_level_name);
        } else if (value_type==CCTK_VARIABLE_INT)
        {
          Util_TableGetInt(tables[i], &prev_time_level, "PREVIOUS_TIME_LEVEL");
        } else
        {
          CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
                     "Invalid data type for key \"PREVIOUS_TIME_LEVEL\" "
                     "Please use CCTK_STRING for the variable name, "
                     "or CCTK_INT for the variable index.");
        }
      }

      /* Wave speed */
      Util_TableGetReal(tables[i], &speed, "SPEED");
    }

    /* Determine boundary width on all faces */
    /* allocate memory for buffer */
    gdim = CCTK_GroupDimI(gi);
    if (gdim > max_gdim)
    {
      width_alldirs = realloc(width_alldirs, 2*gdim*sizeof(CCTK_INT));
      max_gdim = gdim;
    }

    /* fill it with values, either from table or the boundary_width
       parameter */
    if (widths[i]<0)
    {
      err = Util_TableGetIntArray(tables[i], gdim, width_alldirs,
                                  "BOUNDARY_WIDTH");
      if (err<0)
      {
        CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
                   "Error %d when reading boundary width array from table "
                   "for %s", err, CCTK_VarName(vars[i]));
        return -21;
      } else if (err!=2*gdim)
      {
        CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
                   "Boundary width array for %s has %d elements, but %d "
                   "expected", CCTK_VarName(vars[i]), err, 2*gdim);
        return -22;
      }
    } else
    {
      for (k=0; k<2*gdim; ++k)
      {
        width_alldirs[k] = widths[i];
      }
    }

    /* Apply the boundary condition */
    if ((retval = ApplyBndRadiative(GH, 0, width_alldirs, dir, limit, speed,
                                    vars[i], prev_time_level, j)) < 0)
    {
      CCTK_VWarn(1, __LINE__, __FILE__, CCTK_THORNSTRING,
                 "ApplyBndRadiative() returned %d", retval);
    }
  }
#ifdef DEBUG
  printf("BndRadiative(): returning %d\n",retval);
#endif
  free(width_alldirs);

  return retval;
}

/* prototypes for external C routines are declared in header Boundary.h
   here only follow the fortran wrapper prototypes */
void CCTK_FCALL CCTK_FNAME (BndRadiativeDirGI)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil_size,
                            const int *dir,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            const int *gi_to,
                            const int *gi_from);
void CCTK_FCALL CCTK_FNAME (BndRadiativeGI)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            const int *gi_to,
                            const int *gi_from);
void CCTK_FCALL CCTK_FNAME (BndRadiativeDirGN)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil_size,
                            const int *dir,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            TWO_FORTSTRING_ARG);
void CCTK_FCALL CCTK_FNAME (BndRadiativeGN)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            TWO_FORTSTRING_ARG);
void CCTK_FCALL CCTK_FNAME (BndRadiativeDirVI)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil_size,
                            const int *dir,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            const int *vi_to,
                            const int *vi_from);
void CCTK_FCALL CCTK_FNAME (BndRadiativeVI)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            const int *vi_to,
                            const int *vi_from);
void CCTK_FCALL CCTK_FNAME (BndRadiativeDirVN)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil_size,
                            const int *dir,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            TWO_FORTSTRING_ARG);
void CCTK_FCALL CCTK_FNAME (BndRadiativeVN)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            TWO_FORTSTRING_ARG);


/********************************************************************
 ********************    Internal Routines   ************************
 ********************************************************************/
/*@@
   @routine    BndRadiativeDirGI
   @date       Sat Jan 20 2001
   @author     Gabrielle Allen
   @desc
               Aply radiative BC's by group index in given direction
   @enddesc
   @calls      ApplyBndRadiative

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        stencil_size
   @vdesc      stencil size in this direction
   @vtype      int
   @vio        in
   @endvar
   @var        dir
   @vdesc      direction to apply BC
   @vtype      int
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        gi_to
   @vdesc      index of group to apply BC to
   @vtype      int
   @vio        in
   @endvar
   @var        gi_from
   @vdesc      index of group to apply BC from
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine ApplyBndRadiative <BR>
               -1 if invalid group indices are given
   @endreturndesc
@@*/
int BndRadiativeDirGI (const cGH *GH,
                       int stencil_size,
                       int dir,
                       CCTK_REAL var0,
                       CCTK_REAL speed,
                       int gi_to,
                       int gi_from)
{
  int first_vi_to, first_vi_from, retval;


  first_vi_to   = CCTK_FirstVarIndexI (gi_to);
  first_vi_from = CCTK_FirstVarIndexI (gi_from);
  if (first_vi_to >= 0 && first_vi_from >= 0)
  {
    retval = ApplyBndRadiative (GH, stencil_size, NULL, dir, var0, speed,
                                first_vi_to, first_vi_from,
                                CCTK_NumVarsInGroupI (gi_to));
  }
  else
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Invalid group indices %d and/or %d in BndRadiativeDirGI",
                gi_to, gi_from);
    retval = -1;
  }

  return (retval);
}

void CCTK_FCALL CCTK_FNAME (BndRadiativeDirGI)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil_size,
                            const int *dir,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            const int *gi_to,
                            const int *gi_from)
{
  *ierr = BndRadiativeDirGI (*GH, *stencil_size, *dir, *var0, *speed,
                             *gi_to, *gi_from);
}


/*@@
   @routine    BndRadiativeGI
   @date       Tue Jul 18 18:04:07 2000
   @author     Gerd Lanfermann
   @desc
               Aply radiative BC's by group index
   @enddesc
   @calls      ApplyBndRadiative

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        stencil
   @vdesc      stencil width
   @vtype      int [ dimension of group ]
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        gi_to
   @vdesc      index of group to apply BC to
   @vtype      int
   @vio        in
   @endvar
   @var        gi_from
   @vdesc      index of group to apply BC from
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine ApplyBndRadiative
               -1 if invalid group indices are given
   @endreturndesc
@@*/
int BndRadiativeGI (const cGH *GH,
                    const int *stencil,
                    CCTK_REAL var0,
                    CCTK_REAL speed,
                    int gi_to,
                    int gi_from)
{
  int first_vi_to, first_vi_from, retval;


  first_vi_to   = CCTK_FirstVarIndexI (gi_to);
  first_vi_from = CCTK_FirstVarIndexI (gi_from);
  if (first_vi_to >= 0 && first_vi_from >= 0)
  {
    retval = OldApplyBndRadiative (GH, -1, (const CCTK_INT *) stencil, 0, var0,
                                   speed, first_vi_to, first_vi_from,
                                   CCTK_NumVarsInGroupI (gi_to));
  }
  else
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Invalid group indices %d and/or %d in BndRadiativeGI",
                gi_to, gi_from);
    retval = -1;
  }

  return (retval);
}

void CCTK_FCALL CCTK_FNAME (BndRadiativeGI)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            const int *gi_to,
                            const int *gi_from)
{
  *ierr = BndRadiativeGI (*GH, stencil, *var0, *speed, *gi_to, *gi_from);
}


/* ===================================================================== */

/*@@
   @routine    BndRadiativeDirGN
   @date       Sat Jan 20 2001
   @author     Gabrielle Allen
   @desc
               Apply radiative BC's by group name in given direction
   @enddesc
   @calls      BndRadiativeDirGI

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        stencil_size
   @vdesc      stencil size in this direction
   @vtype      int
   @vio        in
   @endvar
   @var        dir
   @vdesc      direction to apply BC
   @vtype      int
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        gname_to
   @vdesc      name of group to apply BC to
   @vtype      const char *
   @vio        in
   @endvar
   @var        gname_from
   @vdesc      name of group to apply BC from
   @vtype      const char *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine BndRadiativeDirGI
               -1 if invalid group names are given
   @endreturndesc
@@*/
int BndRadiativeDirGN (const cGH *GH,
                       int stencil_size,
                       int dir,
                       CCTK_REAL var0,
                       CCTK_REAL speed,
                       const char *gname_to,
                       const char *gname_from)
{
  int gi_to, gi_from, retval;


  gi_to   = CCTK_GroupIndex (gname_to);
  gi_from = CCTK_GroupIndex (gname_from);
  if (gi_to >= 0 && gi_from >= 0)
  {
    retval = BndRadiativeDirGI (GH, stencil_size, dir, var0, speed,
                                gi_to, gi_from);
  }
  else
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Invalid group names '%s' and/or '%s' in BndRadiativeDirGN",
                gname_to, gname_from);
    retval = -1;
  }

  return (retval);
}

void CCTK_FCALL CCTK_FNAME (BndRadiativeDirGN)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil_size,
                            const int *dir,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            TWO_FORTSTRING_ARG)
{
  TWO_FORTSTRINGS_CREATE (gname_to, gname_from)
  *ierr = BndRadiativeDirGN (*GH, *stencil_size, *dir, *var0, *speed,
                             gname_to, gname_from);
  free (gname_to);
  free (gname_from);
}


/*@@
   @routine    BndRadiativeGN
   @date       Tue Jul 18 18:04:07 2000
   @author     Gerd Lanfermann
   @desc
               Aply radiative BC's by group name
   @enddesc
   @calls      BndRadiativeGI

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        stencil
   @vdesc      stencil width
   @vtype      int [ dimension of group ]
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        gname_to
   @vdesc      name of group to apply BC to
   @vtype      const char *
   @vio        in
   @endvar
   @var        gname_from
   @vdesc      name of group to apply BC from
   @vtype      const char *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine BndRadiativeGI <BR>
               -1 if invalid group names are given
   @endreturndesc
@@*/
int BndRadiativeGN (const cGH *GH,
                    const int *stencil,
                    CCTK_REAL var0,
                    CCTK_REAL speed,
                    const char *gname_to,
                    const char *gname_from)
{
  int gi_to, gi_from, retval;


  gi_to   = CCTK_GroupIndex (gname_to);
  gi_from = CCTK_GroupIndex (gname_from);
  if (gi_to >= 0 && gi_from >= 0)
  {
    retval = BndRadiativeGI (GH, stencil, var0, speed, gi_to, gi_from);
  }
  else
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Invalid group names '%s' and/or '%s' in BndRadiativeGN",
                gname_to, gname_from);
    retval = -1;
  }

  return (retval);
}

void CCTK_FCALL CCTK_FNAME (BndRadiativeGN)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            TWO_FORTSTRING_ARG)
{
  TWO_FORTSTRINGS_CREATE (gname_to, gname_from)
  *ierr = BndRadiativeGN (*GH, stencil, *var0, *speed, gname_to, gname_from);
  free (gname_to);
  free (gname_from);
}


/* ===================================================================== */

/*@@
   @routine    BndRadiativeDirVI
   @date       Sat Jan 20 2001
   @author     Gabrielle Allen
   @desc
               Apply radiative BC's by variable index in given direction
   @enddesc
   @calls      ApplyBndRadiative

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        stencil_size
   @vdesc      stencil size in this direction
   @vtype      int
   @vio        in
   @endvar
   @var        dir
   @vdesc      direction to apply BC
   @vtype      int
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        vi_to
   @vdesc      index of variable to apply BC to
   @vtype      int
   @vio        in
   @endvar
   @var        vi_from
   @vdesc      index of variable to apply BC from
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine ApplyBndRadiative <BR>
               -1 if invalid variable indices are given
   @endreturndesc
@@*/
int BndRadiativeDirVI (const cGH *GH,
                       int stencil_size,
                       int dir,
                       CCTK_REAL var0,
                       CCTK_REAL speed,
                       int vi_to,
                       int vi_from)
{
  int retval, num_vars;


  num_vars = CCTK_NumVars ();
  if (vi_to >= 0 && vi_to < num_vars && vi_from >= 0 && vi_from < num_vars)
  {
    retval = ApplyBndRadiative (GH, stencil_size, NULL, dir, var0, speed,
                                vi_to, vi_from, 1);
  }
  else
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Invalid variable indices %d and/or %d in BndRadiativeDirVI",
                vi_to, vi_from);
    retval = -1;
  }

  return (retval);
}

void CCTK_FCALL CCTK_FNAME (BndRadiativeDirVI)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil_size,
                            const int *dir,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            const int *vi_to,
                            const int *vi_from)
{
  *ierr = BndRadiativeDirVI (*GH, *stencil_size, *dir, *var0, *speed,
                             *vi_to, *vi_from);
}


/*@@
   @routine    BndRadiativeVI
   @date       Tue Jul 18 18:04:07 2000
   @author     Gerd Lanfermann
   @desc
               Apply radiative BC's by variable index
   @enddesc
   @calls      ApplyBndRadiative

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        stencil
   @vdesc      stencil width
   @vtype      int [ dimension of variable ]
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        vi_to
   @vdesc      index of variable to apply BC to
   @vtype      int
   @vio        in
   @endvar
   @var        vi_from
   @vdesc      index of variable to apply BC from
   @vtype      int
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine ApplyBndRadiative <BR>
               -1 if invalid variable indices are given
   @endreturndesc
@@*/
int BndRadiativeVI (const cGH *GH,
                    const int *stencil,
                    CCTK_REAL var0,
                    CCTK_REAL speed,
                    int vi_to,
                    int vi_from)
{
  int retval, num_vars;


  num_vars = CCTK_NumVars ();
  if (vi_to >= 0 && vi_to < num_vars && vi_from >= 0 && vi_from < num_vars)
  {
    retval = OldApplyBndRadiative (GH, -1, (const CCTK_INT *) stencil, 0, var0, speed,
                                vi_to, vi_from, 1);
  }
  else
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Invalid variable indices %d and/or %d in BndRadiativeVI",
                vi_to, vi_from);
    retval = -1;
  }

  return (retval);
}

void CCTK_FCALL CCTK_FNAME (BndRadiativeVI)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            const int *vi_to,
                            const int *vi_from)
{
  *ierr = BndRadiativeVI (*GH, stencil, *var0, *speed, *vi_to, *vi_from);
}


/* ======================================================================= */

/*@@
   @routine    BndRadiativeDirVN
   @date       Sat Jan 20 2001
   @author     Gabrielle Allen
   @desc
               Apply radiative BC's by variable name in given direction
   @enddesc
   @calls      BndRadiativeDirVI

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        stencil_size
   @vdesc      stencil size in this direction
   @vtype      int
   @vio        in
   @endvar
   @var        dir
   @vdesc      direction to apply BC
   @vtype      int
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        vname_to
   @vdesc      name of variable to apply BC to
   @vtype      const char *
   @vio        in
   @endvar
   @var        vname_from
   @vdesc      name of variable to apply BC from
   @vtype      const char *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine BndRadiativeDirVI <BR>
               -1 if invalid variable names are given
   @endreturndesc
@@*/
int BndRadiativeDirVN (const cGH *GH,
                       int stencil_size,
                       int dir,
                       CCTK_REAL var0,
                       CCTK_REAL speed,
                       const char *vname_to,
                       const char *vname_from)
{
  int vi_to, vi_from, num_vars, retval;


  vi_to    = CCTK_VarIndex (vname_to);
  vi_from  = CCTK_VarIndex (vname_from);
  num_vars = CCTK_NumVars ();

  if (vi_to >= 0 && vi_to < num_vars && vi_from >= 0 && vi_from < num_vars)
  {
    retval = BndRadiativeDirVI (GH, stencil_size, dir, var0, speed,
                                vi_to, vi_from);
  }
  else
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Invalid variable names '%s' and/or '%s' in BndRadiativeDirVN",
                vname_to, vname_from);
    retval = -1;
  }

  return (retval);
}

void CCTK_FCALL CCTK_FNAME (BndRadiativeDirVN)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil_size,
                            const int *dir,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            TWO_FORTSTRING_ARG)
{
  TWO_FORTSTRINGS_CREATE (vname_to, vname_from)
  *ierr = BndRadiativeDirVN (*GH, *stencil_size, *dir, *var0, *speed,
                             vname_to, vname_from);
  free (vname_to);
  free (vname_from);
}


/*@@
   @routine    BndRadiativeVN
   @date       Tue Jul 18 18:04:07 2000
   @author     Gerd Lanfermann
   @desc
               Apply radiative BC's by variable name
   @enddesc
   @calls      BndRadiativeVI

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        stencil
   @vdesc      stencil width
   @vtype      int [ dimension of variable ]
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        vname_to
   @vdesc      name of variable to apply BC to
   @vtype      const char *
   @vio        in
   @endvar
   @var        vname_from
   @vdesc      name of variable to apply BC from
   @vtype      const char *
   @vio        in
   @endvar

   @returntype int
   @returndesc
               return code of @seeroutine BndRadiativeVI <BR>
               -1 if invalid variable names are given
   @endreturndesc
@@*/
int BndRadiativeVN (const cGH *GH,
                    const int *stencil,
                    CCTK_REAL var0,
                    CCTK_REAL speed,
                    const char *vname_to,
                    const char *vname_from)
{
  int vi_to, vi_from, num_vars, retval;


  vi_to    = CCTK_VarIndex (vname_to);
  vi_from  = CCTK_VarIndex (vname_from);
  num_vars = CCTK_NumVars ();

  if (vi_to >= 0 && vi_to < num_vars && vi_from >= 0 && vi_from < num_vars)
  {
    retval = BndRadiativeVI (GH, stencil, var0, speed, vi_to, vi_from);
  }
  else
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Invalid variable names '%s' and/or '%s' in BndRadiativeVN",
                vname_to, vname_from);
    retval = -1;
  }

  return (retval);
}

void CCTK_FCALL CCTK_FNAME (BndRadiativeVN)
                           (int *ierr,
                            const cGH **GH,
                            const int *stencil,
                            const CCTK_REAL *var0,
                            const CCTK_REAL *speed,
                            TWO_FORTSTRING_ARG)
{
  TWO_FORTSTRINGS_CREATE (vname_to, vname_from)
  *ierr = BndRadiativeVN (*GH, stencil, *var0, *speed, vname_to, vname_from);
  free (vname_to);
  free (vname_from);
}


/********************************************************************
 *********************     Local Routines   *************************
 ********************************************************************/

/* shortcut for multiplying a variable with itself */
#define SQR(a) ((a) * (a))

/* the maximum dimension we can deal with */
#define MAXDIM 3


/*@@
   @routine    LOWER_RADIATIVE_BOUNDARY_3D
   @date       Mon 9 Apr 2001
   @author     Thomas Radke
   @desc
               Macro to apply radiative BC to a lower bound of a 3D variable
   @enddesc

   @var        istart, jstart, kstart
   @vdesc      start index for the x,y,z direction
   @vtype      int
   @vio        in
   @endvar
   @var        dim
   @vdesc      dimension to apply BC
   @vtype      int
   @vio        in
   @endvar
   @var        cctk_type
   @vdesc      CCTK datatypes of the source and target variable
   @vtype      <cctk_type>
   @vio        in
   @endvar
@@*/
#define LOWER_RADIATIVE_BOUNDARY_3D(istart, jstart, kstart, dim, cctk_type)   \
{                                                                             \
  int _i, _j, _k;                                                             \
  int _0 = 0*offset[dim],                                                     \
      _1 = 1*offset[dim],                                                     \
      _2 = 2*offset[dim];                                                     \
                                                                              \
                                                                              \
  for (_k = kstart-1; _k >= 0; _k--)                                          \
  {                                                                           \
    for (_j = jstart-1; _j >= 0; _j--)                                        \
    {                                                                         \
      int _idx = CCTK_GFINDEX3D (GH, istart-1, _j, _k);                       \
      const CCTK_REAL *_r    = xyzr[MAXDIM] + _idx,                           \
                      *_xyz  = xyzr[dim]    + _idx;                           \
      cctk_type *_to         = (cctk_type *) to_ptr         + _idx;           \
      const cctk_type *_from = (const cctk_type *) from_ptr + _idx;           \
                                                                              \
                                                                              \
      for (_i = istart-1; _i >= 0; _i--)                                      \
      {                                                                       \
        CCTK_REAL _r0_inv = 1/_r[_0], _r1_inv = 1/_r[_1];                     \
                                                                              \
        if (radpower > 0)                                                     \
        {                                                                     \
          CCTK_REAL H;                                                        \
                                                                              \
          H  = 0.25 * radpower * dxyz[dim]                                    \
               * (_xyz[_0]*SQR (_r0_inv) + _xyz[_1]*SQR (_r1_inv));           \
          H  = (1 + H) / (1 - H);                                             \
          H *= dtv * (0.25*(_to[_1] + _to[_2] + _from[_1] + _from[_2]) - var0)\
               + 0.5*(  _r[_1] * (_to[_1] - _from[_1])                        \
                      + _r[_2] * (_to[_2] - _from[_2]))                       \
               + 0.25*(_to[_2] - _to[_1] + _from[_2] - _from[_1]) * rho[dim]  \
                     *(  SQR (_r[_1]) / _xyz[_1]                              \
                       + SQR (_r[_2]) / _xyz[_2]);                            \
          dtvvar0H = dtvvar0 + H;                                             \
        }                                                                     \
                                                                              \
        _to[_0] = (cctk_type) (                                               \
                   (dtvvar0H * (  _xyz[_0] * SQR (_r0_inv)                    \
                                + _xyz[_1] * SQR (_r1_inv))                   \
                    - _to[_1]   * ( rho[dim] + _xyz[_1]*_r1_inv               \
                                                     *(1 + dtvh*_r1_inv))     \
                    + _from[_0] * ( rho[dim] + _xyz[_0]*_r0_inv               \
                                                     *(1 - dtvh*_r0_inv))     \
                    - _from[_1] * ( rho[dim] - _xyz[_1]*_r1_inv               \
                                                     *(1 - dtvh*_r1_inv))     \
                   )                                                          \
                   /              (-rho[dim] + _xyz[_0]*_r0_inv               \
                                                     *(1 + dtvh*_r0_inv))     \
                  );                                                          \
        _r--; _xyz--; _to--; _from--;                                         \
      }                                                                       \
    }                                                                         \
  }                                                                           \
}


/*@@
   @routine    UPPER_RADIATIVE_BOUNDARY_3D
   @date       Mon 9 Apr 2001
   @author     Thomas Radke
   @desc
               Macro to apply radiative BC to an upper bound of a 3D variable
   @enddesc

   @var        istart, jstart, kstart
   @vdesc      start index for the x,y,z direction
   @vtype      int
   @vio        in
   @endvar
   @var        dim
   @vdesc      dimension to apply BC
   @vtype      int
   @vio        in
   @endvar
   @var        cctk_type
   @vdesc      CCTK datatypes of the source and target variable
   @vtype      <cctk_type>
   @vio        in
   @endvar
@@*/
#define UPPER_RADIATIVE_BOUNDARY_3D(istart, jstart, kstart, dim, cctk_type)   \
{                                                                             \
  int _i, _j, _k;                                                             \
  int _0 = -0*offset[dim],                                                    \
      _1 = -1*offset[dim],                                                    \
      _2 = -2*offset[dim];                                                    \
                                                                              \
                                                                              \
  for (_k = kstart; _k < GH->cctk_lsh[2]; _k++)                               \
  {                                                                           \
    for (_j = jstart; _j < GH->cctk_lsh[1]; _j++)                             \
    {                                                                         \
      int _idx = CCTK_GFINDEX3D (GH, istart, _j, _k);                         \
      const CCTK_REAL *_r    = xyzr[MAXDIM] + _idx,                           \
                      *_xyz  = xyzr[dim]    + _idx;                           \
      cctk_type *_to         = (cctk_type *) to_ptr         + _idx;           \
      const cctk_type *_from = (const cctk_type *) from_ptr + _idx;           \
                                                                              \
                                                                              \
      for (_i = istart; _i < GH->cctk_lsh[0]; _i++)                           \
      {                                                                       \
        CCTK_REAL _r0_inv = 1/_r[_0], _r1_inv = 1/_r[_1];                     \
                                                                              \
        if (radpower > 0)                                                     \
        {                                                                     \
          CCTK_REAL H;                                                        \
                                                                              \
          H  = 0.25 * radpower * dxyz[dim]                                    \
               * (_xyz[_0]*SQR (_r0_inv) + _xyz[_1]*SQR (_r1_inv));           \
          H  = (1 - H) / (1 + H);                                             \
          H *= dtv * (0.25*(_to[_1] + _to[_2] + _from[_1] + _from[_2]) - var0)\
               + 0.5 *(  _r[_1] * (_to[_1] - _from[_1])                       \
                       + _r[_2] * (_to[_2] - _from[_2]))                      \
               + 0.25*(_to[_1] - _to[_2] + _from[_1] - _from[_2]) * rho[dim]  \
                     * (  SQR (_r[_1]) / _xyz[_1]                             \
                        + SQR (_r[_2]) / _xyz[_2]);                           \
          dtvvar0H = dtvvar0 + H;                                             \
        }                                                                     \
                                                                              \
        _to[_0] = (cctk_type) (                                               \
                   (dtvvar0H * (  _xyz[_0] * (SQR (_r0_inv))                  \
                                + _xyz[_1] * (SQR (_r1_inv)))                 \
                    + _to[_1]   * ( rho[dim] - _xyz[_1]*_r1_inv               \
                                                     *(1 + dtvh*_r1_inv))     \
                    + _from[_0] * (-rho[dim] + _xyz[_0]*_r0_inv               \
                                                     *(1 - dtvh*_r0_inv))     \
                    + _from[_1] * ( rho[dim] + _xyz[_1]*_r1_inv               \
                                                     *(1 - dtvh*_r1_inv))     \
                   )                                                          \
                   /              ( rho[dim] + _xyz[_0]*_r0_inv               \
                                                     *(1 + dtvh*_r0_inv))     \
                  );                                                          \
        _r++; _xyz++; _to++; _from++;                                         \
      }                                                                       \
    }                                                                         \
  }                                                                           \
}


/*@@
   @routine    RADIATIVE_BOUNDARY
   @date       Mon 9 Apr 2001
   @author     Thomas Radke
   @desc
               Macro to apply radiative BC to a variable
               Currently it is limited to 3D variables only.
   @enddesc
   @calls      LOWER_RADIATIVE_BOUNDARY_3D
               UPPER_RADIATIVE_BOUNDARY_3D

   @var        lsh
   @vdesc      local shape of the variable
   @vtype      int [ dim ]
   @vio        in
   @endvar
   @var        stencil
   @vdesc      stencils in every direction
   @vtype      int [ 2*dim ]
   @vio        in
   @endvar
   @var        cctk_type
   @vdesc      CCTK datatypes of the source and target variable
   @vtype      <cctk_type>
   @vio        in
   @endvar
@@*/
#define RADIATIVE_BOUNDARY(lsh, stencil, cctk_type)                           \
{                                                                             \
  /* check the dimensionality */                                              \
  if (gdim != 3)                                                              \
  {                                                                           \
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,                      \
                "ApplyBndRadiative: variable dimension of %d not supported",  \
                gdim);                                                        \
    return (-5);                                                              \
  }                                                                           \
                                                                              \
  /* Lower x-bound */                                                         \
  if (doBC[0])                                                                \
  {                                                                           \
    LOWER_RADIATIVE_BOUNDARY_3D (stencil[0], lsh[1], lsh[2], 0, cctk_type);   \
  }                                                                           \
                                                                              \
  /* Upper x-bound */                                                         \
  if (doBC[1])                                                                \
  {                                                                           \
    UPPER_RADIATIVE_BOUNDARY_3D (lsh[0] - stencil[1], 0, 0, 0, cctk_type);    \
  }                                                                           \
                                                                              \
  /* Lower y-bound */                                                         \
  if (doBC[2])                                                                \
  {                                                                           \
    LOWER_RADIATIVE_BOUNDARY_3D (lsh[0], stencil[2], lsh[2], 1, cctk_type);   \
  }                                                                           \
                                                                              \
  /* Upper y-bound */                                                         \
  if (doBC[3])                                                                \
  {                                                                           \
    UPPER_RADIATIVE_BOUNDARY_3D (0, lsh[1] - stencil[3], 0, 1, cctk_type);    \
  }                                                                           \
                                                                              \
  /* Lower z-bound */                                                         \
  if (doBC[4])                                                                \
  {                                                                           \
    LOWER_RADIATIVE_BOUNDARY_3D (lsh[0], lsh[1], stencil[4], 2, cctk_type);   \
  }                                                                           \
                                                                              \
  /* Upper z-bound */                                                         \
  if (doBC[5])                                                                \
  {                                                                           \
    UPPER_RADIATIVE_BOUNDARY_3D (0, 0, lsh[2] - stencil[5], 2, cctk_type);    \
  }                                                                           \
}


/*@@
   @routine    ApplyBndRadiative
   @date       Tue Jul 18 18:04:07 2000
   @author     Gerd Lanfermann
   @desc
               Apply radiation boundary conditions to a group of grid functions
               given by their indices
               This routine is called by the various BndRadiativeXXX wrappers.

               Although it is currently limited to handle 3D variables only
               it can easily be extended for other dimensions
               by adapting the appropriate macros.
   @enddesc

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        width_dir
   @vdesc      boundary width in direction dir
   @vtype      int
   @vio        in
   @endvar
   @var        in_widths
   @vdesc      boundary widths for all directions
   @vtype      const CCTK_INT [ dimension of variable(s) ]
   @vio        in
   @endvar
   @var        dir
   @vdesc      direction to copy boundaries (0 for copying all directions)
   @vtype      int
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        first_var_to
   @vdesc      index of first variable to copy boundaries to
   @vtype      CCTK_INT
   @vio        in
   @endvar
   @var        first_var_from
   @vdesc      index of first variable to copy boundaries from
   @vtype      CCTK_INT
   @vio        in
   @endvar
   @var        num_vars
   @vdesc      number of variables
   @vtype      int
   @vio        in
   @endvar
   @calls      CCTK_VarTypeI
               CCTK_GroupDimFromVarI
               RADIATIVE_BOUNDARY
   @history
   @hdate      Mon 9 Apr 2001
   @hauthor    Thomas Radke
   @hdesc      Merged separate routines for 1D, 2D, and 3D
               into a single generic routine
   @endhistory

   @returntype int
   @returndesc
                0 for success
               -1 if abs(direction) is greater than variables' dimension
               -2 if variable dimension is not supported
               -3 if NULL pointer passed as boundary width array
               -4 if variable type is not supported
               -5 if variable dimension is other than 3D
               -6 if a coordinate is not found
   @endreturndesc
@@*/
static int ApplyBndRadiative (const cGH *GH,
                              int width_dir,
                              const CCTK_INT *in_widths,
                              int dir,
                              CCTK_REAL var0,
                              CCTK_REAL speed,
                              CCTK_INT first_var_to,
                              CCTK_INT first_var_from,
                              int num_vars)
{
  int i, gdim, indx;
  int var_to, var_from;
  int timelvl_from;
  char coord_system_name[10];
  CCTK_REAL dxyz[MAXDIM], rho[MAXDIM];
  const CCTK_REAL *xyzr[MAXDIM+1];
  int doBC[2*MAXDIM], widths[2*MAXDIM], offset[MAXDIM];
  CCTK_INT symtable;
  CCTK_INT symbnd[2*MAXDIM];
  CCTK_INT is_physical[2*MAXDIM];
  CCTK_INT ierr;
  CCTK_REAL dtv, dtvh, dtvvar0, dtvvar0H;
  void *to_ptr;
  const void *from_ptr;
  DECLARE_CCTK_PARAMETERS


  /* check the direction parameter */
  if (abs (dir) > MAXDIM)
  {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "ApplyBndRadiative: direction %d is greater than maximum "
                "dimension %d", dir, MAXDIM);
    return (-1);
  }

  /* get the dimensionality */
  gdim = CCTK_GroupDimFromVarI (first_var_to);

  /* check the dimensionality */
  if (gdim > MAXDIM)
  {
    CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                "ApplyBndRadiative: variable dimension of %d not supported",
                gdim);
    return (-2);
  }

  /* set up boundary width array */
  if (dir)
  {
    widths[2*(abs(dir)-1)] = width_dir;
    widths[2*(abs(dir)-1)+1] = width_dir;
  }
  else if (in_widths)
  {
    memcpy (widths, in_widths, 2* gdim * sizeof (CCTK_INT));
  }
  else
  {
    CCTK_WARN (1, "ApplyBndRadiative: NULL pointer passed "
                  "for boundary width array");
    return (-3);
  }

  /* Use next time level, if available */
  timelvl_from = 0;
  if (CCTK_MaxTimeLevelsVI(first_var_from) > 1)
  {
    timelvl_from = 1;
  }

  /* Find Courant parameters. */
  dtv      = speed * GH->cctk_delta_time;
  dtvh     = 0.5 * dtv;
  dtvvar0  = dtv * var0;
  dtvvar0H = dtvvar0;

  sprintf (coord_system_name, "cart%dd", gdim);
  for (i = 0; i < gdim; i++)
  {
    /* Radiative boundaries need the underlying Cartesian coordinates */
    indx = CCTK_CoordIndex (i + 1, NULL, coord_system_name);
    if (indx < 0)
    {
      CCTK_VWarn(1,__LINE__,__FILE__,CCTK_THORNSTRING,
                "Coordinate for system %s not found",coord_system_name);
      return (-6);
    }
    xyzr[i] = GH->data[indx][0];

    /* According to the Cactus spec, the true delta_space values for a
       grid are calculated as follows: */
    dxyz[i] = GH->cctk_delta_space[i] / GH->cctk_levfac[i];

    rho[i] = dtv / dxyz[i];

    offset[i] = i == 0 ? 1 : offset[i-1] * GH->cctk_lsh[i-1];
  }

  /* Append r grid variable to end of xyzr[] array */
  sprintf (coord_system_name, "spher%dd", gdim);
  indx = CCTK_CoordIndex (-1, "r", coord_system_name);
  if (indx < 0)
  {
    CCTK_VWarn(1,__LINE__,__FILE__,CCTK_THORNSTRING,
              "Coordinate for system %s not found",coord_system_name);
    return (-6);
  }
  xyzr[MAXDIM] = GH->data[indx][0];


  /* see if we have a physical boundary */
  symtable = SymmetryTableHandleForGrid (GH);
  if (symtable < 0) CCTK_WARN (0, "internal error");
  ierr = Util_TableGetIntArray (symtable, 2 * gdim, symbnd, "symmetry_handle");
  if (ierr != 2 * gdim) CCTK_WARN (0, "internal error");
  for (i = 0; i < 2 * gdim; i++)
  {
    is_physical[i] = symbnd[i] < 0;
  }

  /* now loop over all variables */
  for (var_to = first_var_to, var_from = first_var_from;
       var_to < first_var_to + num_vars;
       var_to++, var_from++)
  {
    to_ptr   = GH->data[var_to][0];
    from_ptr = GH->data[var_from][timelvl_from];

    /* Apply condition if:
       + boundary is a physical boundary
       + boundary is an outer boundary
       + have enough grid points
    */
    for (i = 0; i < 2 * MAXDIM; i++)
    {
      doBC[i] = is_physical[i];
    }
    for (i = 0; i < MAXDIM; i++)
    {
      doBC[i*2]   &= GH->cctk_lsh[i] > widths[i*2] && GH->cctk_bbox[i*2];
      doBC[i*2+1] &= GH->cctk_lsh[i] > widths[i*2+1] && GH->cctk_bbox[i*2+1];
      if (dir != 0)
      {
        doBC[i*2]   &= (dir < 0 && (i + 1 == abs (dir)));
        doBC[i*2+1] &= (dir > 0 && (i + 1 == abs (dir)));
      }
    }

    switch (CCTK_VarTypeI (var_to))
    {
      case CCTK_VARIABLE_REAL:
        RADIATIVE_BOUNDARY (GH->cctk_lsh, widths, CCTK_REAL);
        break;

#ifdef CCTK_REAL4
      case CCTK_VARIABLE_REAL4:
        RADIATIVE_BOUNDARY (GH->cctk_lsh, widths, CCTK_REAL4);
        break;
#endif

#ifdef CCTK_REAL8
      case CCTK_VARIABLE_REAL8:
        RADIATIVE_BOUNDARY (GH->cctk_lsh, widths, CCTK_REAL8);
        break;
#endif

#ifdef CCTK_REAL16
      case CCTK_VARIABLE_REAL16:
        RADIATIVE_BOUNDARY (GH->cctk_lsh, widths, CCTK_REAL16);
        break;
#endif

      default:
        CCTK_VWarn (1, __LINE__, __FILE__, CCTK_THORNSTRING,
                    "Unsupported variable type %d for variable '%s'",
                    CCTK_VarTypeI (var_to), CCTK_VarName (var_to));
        return (-4);
    }
  }

  return (0);
}


/*@@
   @routine    OldApplyBndRadiative
   @date       5 May 2003
   @author     David Rideout
   @desc
               The new boundary API expects a 2d-element array for the
               boundary_widths (d=dimension of grid variable), while
               the old API expects a d-element array.  This function
               converts the old array to the new format.
   @enddesc

   @var        GH
   @vdesc      Pointer to CCTK grid hierarchy
   @vtype      const cGH *
   @vio        in
   @endvar
   @var        width_dir
   @vdesc      boundary width in direction dir
   @vtype      int
   @vio        in
   @endvar
   @var        stencil_alldirs
   @vdesc      boundary widths for all directions
   @vtype      const CCTK_INT [ dimension of variable(s) ]
   @vio        in
   @endvar
   @var        dir
   @vdesc      direction to copy boundaries (0 for copying all directions)
   @vtype      int
   @vio        in
   @endvar
   @var        var0
   @vdesc      asymptotic value of function at infinity
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        speed
   @vdesc      wave speed
   @vtype      CCTK_REAL
   @vio        in
   @endvar
   @var        first_var_to
   @vdesc      index of first variable to copy boundaries to
   @vtype      int
   @vio        in
   @endvar
   @var        first_var_from
   @vdesc      index of first variable to copy boundaries from
   @vtype      int
   @vio        in
   @endvar
   @var        num_vars
   @vdesc      number of variables
   @vtype      int
   @vio        in
   @endvar

   @calls      CCTK_GroupIndexFromVarI
               ApplyBndRadiative
   @returntype int
   @returndesc
               returncode from @seeroutine ApplyBndRadiative
   @endreturndesc
@@*/

static int OldApplyBndRadiative(const cGH *GH,
                                int width_dir,
                                const CCTK_INT *stencil_alldirs,
                                int dir,
                                CCTK_REAL var0,
                                CCTK_REAL speed,
                                int first_var_to,
                                int first_var_from,
                                int num_vars)
{
  int i, dim, retval;
  CCTK_INT *boundary_widths;
  static int warned;

  /* Convert stencil_alldirs to new format */
  dim = CCTK_GroupDimFromVarI(first_var_to);
  boundary_widths = malloc(2*dim*sizeof(CCTK_INT));
  for (i=0; i<2*dim; ++i)
  {
    boundary_widths[i] = stencil_alldirs[i/2];
  }

  /* Bug people for using the old interface */
  if (!warned)
  {
    CCTK_VWarn (2, __LINE__, __FILE__, CCTK_THORNSTRING,
                "Copied older d-element array of boundary widths into the "
                "newer 2d-element format.  Please use the new boundary "
                "interface to avoid this.");
    warned = 1;
  }

  /* Call ApplyBnd... with new boundary width array */
  retval = ApplyBndRadiative(GH, width_dir, boundary_widths, dir, var0, speed,
                             first_var_to, first_var_from, num_vars);

  free(boundary_widths);
  return retval;
}