aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetLib/src/dh.cc
blob: 7b82508396595498184b31b944fd9d089fd698f3 (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
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
2093
2094
2095
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108
2109
2110
2111
2112
2113
2114
2115
2116
2117
2118
2119
2120
2121
2122
2123
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138
2139
2140
2141
2142
2143
2144
2145
2146
2147
2148
2149
2150
2151
2152
2153
2154
2155
2156
2157
2158
2159
2160
2161
2162
2163
2164
2165
2166
2167
2168
2169
2170
2171
2172
2173
2174
2175
2176
2177
2178
2179
2180
2181
2182
2183
2184
2185
2186
2187
2188
2189
2190
2191
2192
2193
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
2206
2207
2208
2209
2210
2211
2212
2213
2214
2215
2216
2217
2218
2219
2220
2221
2222
2223
2224
2225
2226
2227
2228
2229
2230
2231
2232
2233
2234
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
2263
2264
2265
2266
2267
2268
2269
2270
2271
2272
2273
2274
2275
2276
2277
2278
2279
2280
2281
2282
2283
2284
2285
2286
2287
2288
2289
2290
2291
2292
2293
2294
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308
2309
2310
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324
2325
2326
2327
2328
2329
2330
2331
2332
2333
2334
2335
2336
2337
2338
2339
2340
2341
2342
2343
2344
2345
2346
2347
2348
2349
2350
2351
2352
2353
2354
2355
2356
2357
2358
2359
2360
2361
2362
2363
2364
2365
2366
2367
2368
2369
2370
2371
2372
2373
2374
2375
2376
2377
2378
2379
2380
2381
2382
2383
2384
2385
2386
2387
2388
2389
2390
2391
2392
2393
2394
2395
2396
2397
2398
2399
2400
2401
2402
2403
2404
2405
2406
2407
2408
2409
2410
2411
2412
2413
2414
2415
2416
2417
2418
2419
2420
2421
2422
2423
2424
2425
2426
2427
2428
2429
2430
2431
2432
2433
2434
2435
2436
2437
2438
2439
2440
2441
2442
2443
2444
2445
2446
2447
2448
2449
2450
2451
2452
2453
2454
2455
2456
2457
2458
2459
2460
2461
2462
2463
2464
2465
2466
2467
2468
2469
2470
2471
2472
2473
2474
2475
2476
2477
2478
2479
2480
2481
2482
2483
2484
2485
2486
2487
2488
2489
2490
2491
2492
2493
2494
2495
2496
2497
2498
2499
2500
2501
2502
2503
2504
2505
2506
2507
2508
2509
2510
2511
2512
2513
2514
2515
2516
2517
2518
2519
2520
2521
2522
2523
2524
2525
2526
2527
2528
2529
2530
2531
2532
2533
2534
2535
2536
2537
2538
2539
2540
2541
2542
2543
2544
2545
2546
2547
2548
2549
2550
2551
2552
2553
2554
2555
2556
2557
2558
2559
2560
2561
2562
2563
2564
2565
2566
2567
2568
2569
2570
2571
2572
2573
2574
2575
2576
2577
2578
2579
2580
2581
2582
2583
2584
2585
2586
2587
2588
2589
2590
2591
2592
2593
2594
2595
2596
2597
2598
2599
2600
2601
2602
2603
2604
2605
2606
2607
2608
2609
2610
2611
2612
2613
2614
2615
2616
2617
2618
2619
2620
2621
2622
2623
2624
2625
2626
2627
2628
2629
2630
2631
2632
2633
2634
2635
2636
2637
2638
2639
2640
2641
2642
2643
2644
2645
2646
2647
2648
2649
2650
2651
2652
2653
2654
2655
2656
2657
2658
2659
2660
2661
2662
2663
2664
2665
2666
2667
2668
2669
2670
2671
2672
2673
2674
2675
2676
2677
2678
2679
2680
2681
2682
2683
2684
2685
2686
2687
2688
2689
2690
2691
2692
2693
2694
2695
2696
2697
2698
2699
2700
2701
2702
2703
2704
2705
2706
2707
2708
2709
2710
2711
2712
2713
2714
2715
2716
2717
2718
2719
2720
2721
2722
2723
2724
2725
2726
2727
2728
2729
2730
2731
2732
2733
2734
2735
2736
2737
2738
2739
2740
2741
2742
2743
2744
2745
2746
2747
2748
2749
2750
2751
2752
2753
2754
2755
2756
2757
2758
2759
2760
2761
2762
2763
2764
2765
2766
2767
2768
2769
2770
2771
2772
2773
2774
2775
2776
2777
2778
2779
2780
2781
2782
2783
2784
2785
2786
2787
2788
2789
2790
2791
2792
2793
2794
2795
2796
2797
2798
2799
2800
2801
2802
2803
2804
2805
2806
2807
2808
2809
2810
2811
2812
2813
2814
2815
2816
2817
2818
2819
2820
2821
2822
2823
2824
2825
2826
2827
2828
2829
2830
2831
2832
2833
2834
2835
2836
2837
2838
2839
2840
2841
2842
2843
2844
2845
2846
2847
2848
2849
2850
2851
2852
2853
2854
2855
2856
2857
2858
2859
2860
2861
2862
2863
2864
2865
2866
2867
2868
2869
2870
2871
2872
2873
2874
2875
2876
2877
2878
#include <cctk.h>
#include <cctk_Parameters.h>

#include <cassert>
#include <cstddef>
#include <sstream>

#include <Timer.hh>

#include "mpi_string.hh"
#include "bbox.hh"
#include "bboxset.hh"
#include "defs.hh"
#include "dist.hh"
#include "ggf.hh"
#include "timestat.hh"
#include "vect.hh"

#include "dh.hh"

using namespace std;
using namespace CarpetLib;



set<dh*> dh::alldh;



// Indexing over neighbours

static int num_nbs ()
{
  return ipow(3,dim);
}

static ivect ind2nb (int ind)
{
  ivect nb;
  for (int d=0; d<dim; ++d) {
    nb[d] = ind % 3 - 1;
    ind /= 3;
  }
  assert (ind == 0);
  return nb;
}



vect<vect<dh::srpvect dh::fast_dboxes::*,2>,dim>
dh::fast_dboxes::
fast_ref_refl_sendrecv;

void
dh::fast_dboxes::
init_fast_ref_refl_sendrecv ()
{
  static bool initialised = false;
  if (initialised) return;
  initialised = true;
  fast_ref_refl_sendrecv[0][0] = & fast_dboxes::fast_ref_refl_sendrecv_0_0;
  fast_ref_refl_sendrecv[0][1] = & fast_dboxes::fast_ref_refl_sendrecv_0_1;
  fast_ref_refl_sendrecv[1][0] = & fast_dboxes::fast_ref_refl_sendrecv_1_0;
  fast_ref_refl_sendrecv[1][1] = & fast_dboxes::fast_ref_refl_sendrecv_1_1;
  fast_ref_refl_sendrecv[2][0] = & fast_dboxes::fast_ref_refl_sendrecv_2_0;
  fast_ref_refl_sendrecv[2][1] = & fast_dboxes::fast_ref_refl_sendrecv_2_1;
}

vect<vect<dh::srpvect dh::fast_dboxes::*,2>,dim>
dh::fast_dboxes::
fast_ref_refl_prol_sendrecv;

void
dh::fast_dboxes::
init_fast_ref_refl_prol_sendrecv ()
{
  static bool initialised = false;
  if (initialised) return;
  initialised = true;
  fast_ref_refl_prol_sendrecv[0][0] =
    & fast_dboxes::fast_ref_refl_prol_sendrecv_0_0;
  fast_ref_refl_prol_sendrecv[0][1] =
    & fast_dboxes::fast_ref_refl_prol_sendrecv_0_1;
  fast_ref_refl_prol_sendrecv[1][0] =
    & fast_dboxes::fast_ref_refl_prol_sendrecv_1_0;
  fast_ref_refl_prol_sendrecv[1][1] =
    & fast_dboxes::fast_ref_refl_prol_sendrecv_1_1;
  fast_ref_refl_prol_sendrecv[2][0] =
    & fast_dboxes::fast_ref_refl_prol_sendrecv_2_0;
  fast_ref_refl_prol_sendrecv[2][1] =
    & fast_dboxes::fast_ref_refl_prol_sendrecv_2_1;
}



// Constructors
dh::
dh (gh & h_,
    vector<i2vect> const & ghost_widths_,
    vector<i2vect> const & buffer_widths_,
    vector<i2vect> const & overlap_widths_,
    vector<int> const & prolongation_orders_space_)
  : h(h_),
    ghost_widths(ghost_widths_),
    buffer_widths(buffer_widths_),
    overlap_widths(overlap_widths_),
    prolongation_orders_space(prolongation_orders_space_)
{
  fast_dboxes::init_fast_ref_refl_sendrecv();
  fast_dboxes::init_fast_ref_refl_prol_sendrecv();
  size_t const maxreflevels = h.reffacts.size();
  assert (ghost_widths.size() >= maxreflevels);
  assert (buffer_widths.size() >= maxreflevels);
  assert (overlap_widths.size() >= maxreflevels);
  assert (prolongation_orders_space.size() >= maxreflevels);
  for (size_t rl=0; rl<maxreflevels; ++rl) {
    assert (all (all (ghost_widths.AT(rl) >= 0)));
    assert (all (all (buffer_widths.AT(rl) >= 0)));
    assert (all (all (overlap_widths.AT(rl) >= 0)));
    assert (prolongation_orders_space.AT(rl) >= 0);
  }
  
  alldh.insert(this);
  h.insert (this);
  CHECKPOINT;
  regrid (false);
  for (int rl = 0; rl < h.reflevels(); ++ rl) {
    recompose (rl, false);
  }
  regrid_free (false);
}



// Destructors
dh::
~dh ()
{
  CHECKPOINT;
  assert (gfs.empty());
  h.erase (this);
  alldh.erase (this);
}



// Helpers
int
dh::
prolongation_stencil_size (int const rl)
  const
{
  assert (prolongation_orders_space.AT(rl) >= 0);
  return prolongation_orders_space.AT(rl) / 2;
}



// Modifiers

// Calculate this quantity on this process? It does not need to be
// calculated if it won't be used later on.
inline
int
dh::this_proc (int const rl, int const c) const
{
  return h.processor (rl, c);
}

inline
bool
dh::on_this_proc (int const rl, int const c) const
{
  return this_proc (rl, c) == dist::rank();
}

inline
int
dh::this_oldproc (int const rl, int const c) const
{
  return h.old_processor (rl, c);
}

inline
bool
dh::on_this_oldproc (int const rl, int const c) const
{
  return this_oldproc (rl, c) == dist::rank();
}



bool there_was_an_error = false;

static
void
assert_error (char const * restrict const checkstring,
              char const * restrict const file, int const line,
              int const ml, int const rl,
              char const * restrict const message)
{
  CCTK_VWarn (CCTK_WARN_ALERT, __LINE__, __FILE__, CCTK_THORNSTRING,
              "\n%s:%d:\n   [ml=%d rl=%d] The following grid structure consistency check failed:\n   %s\n   %s",
              file, line, ml, rl, message, checkstring);
  there_was_an_error = true;
}

static
void
assert_error (char const * restrict const checkstring,
              char const * restrict const file, int const line,
              int const ml, int const rl, int const c,
              char const * restrict const message)
{
  CCTK_VWarn (CCTK_WARN_ALERT, __LINE__, __FILE__, CCTK_THORNSTRING,
              "\n%s:%d:\n   [ml=%d rl=%d c=%d] The following grid structure consistency check failed:\n   %s\n   %s",
              file, line, ml, rl, c, message, checkstring);
  there_was_an_error = true;
}

static
void
assert_error (char const * restrict const checkstring,
              char const * restrict const file, int const line,
              int const ml, int const rl, int const c, int const cc,
              char const * restrict const message)
{
  CCTK_VWarn (CCTK_WARN_ALERT, __LINE__, __FILE__, CCTK_THORNSTRING,
              "\n%s:%d:\n   [ml=%d rl=%d c=%d cc=%d] The following grid structure consistency check failed:\n   %s\n   %s",
              file, line, ml, rl, c, cc, message, checkstring);
  there_was_an_error = true;
}

#ifdef CARPET_OPTIMISE

// For highest efficiency, omit all self-checks
#define ASSERT_rl(check, message)
#define ASSERT_c(check, message)
#define ASSERT_cc(check, message)

#else

#define ASSERT_rl(check, message)                                       \
  do {                                                                  \
    if (not (check)) {                                                  \
      assert_error (#check, __FILE__, __LINE__, ml, rl, message);       \
    }                                                                   \
  } while (false)

#define ASSERT_c(check, message)                                        \
  do {                                                                  \
    if (not (check)) {                                                  \
      assert_error (#check, __FILE__, __LINE__, ml, rl, c, message);    \
    }                                                                   \
  } while (false)

#define ASSERT_cc(check, message)                                       \
  do {                                                                  \
    if (not (check)) {                                                  \
      assert_error (#check, __FILE__, __LINE__, ml, rl, c, cc, message); \
    }                                                                   \
  } while (false)

#endif



void
dh::
regrid (bool const do_init)
{
  DECLARE_CCTK_PARAMETERS;
    
  static Timers::Timer timer ("CarpetLib::dh::regrid");
  timer.start();
  
  CHECKPOINT;
  
  static Timer total ("CarpetLib::dh::regrid");
  total.start ();
  
  light_mboxes old_light_boxes;
  swap (light_boxes, old_light_boxes);
  
  full_mboxes full_boxes;
  
  fast_boxes.clear();
  local_boxes.clear();
  
  
  
  light_boxes.resize (h.mglevels());
  local_boxes.resize (h.mglevels());
  level_boxes.resize (h.mglevels());
  full_boxes.resize (h.mglevels());
  fast_boxes.resize (h.mglevels());
  for (int ml = 0; ml < h.mglevels(); ++ ml) {
    light_boxes.AT(ml).resize (h.reflevels());
    local_boxes.AT(ml).resize (h.reflevels());
    level_boxes.AT(ml).resize (h.reflevels());
    full_boxes.AT(ml).resize (h.reflevels());
    fast_boxes.AT(ml).resize (h.reflevels());
    for (int rl = 0; rl < h.reflevels(); ++ rl) {
      light_boxes.AT(ml).AT(rl).resize (h.components(rl));
      local_boxes.AT(ml).AT(rl).resize (h.local_components(rl));
      full_boxes.AT(ml).AT(rl).resize (h.components(rl));
      
      light_cboxes & light_level = light_boxes.AT(ml).AT(rl);
      local_cboxes & local_level = local_boxes.AT(ml).AT(rl);
      level_dboxes & level_level = level_boxes.AT(ml).AT(rl);
      full_cboxes & full_level = full_boxes.AT(ml).AT(rl);
      fast_dboxes & fast_level = fast_boxes.AT(ml).AT(rl);
      
      vector<fast_dboxes> fast_level_otherprocs (dist::size());
      
      
      
      i2vect const& ghost_width = ghost_widths.AT(rl);
      i2vect const& buffer_width = buffer_widths.AT(rl);
      i2vect const& overlap_width = overlap_widths.AT(rl);
      
      
      
      // Domain:
      
      static Timers::Timer timer_domain ("domain");
      timer_domain.start();
      
      ibbox const & domain_exterior = h.baseextent(ml,rl);
      // Variables may have size zero
      // ASSERT_rl (not domain_exterior.empty(),
      //            "The exterior of the domain must not be empty");
      
      i2vect const & boundary_width = h.boundary_width;
      ASSERT_rl (all (all (boundary_width >= 0)),
                 "The gh boundary widths must not be negative");
      
      ibbox const domain_active = domain_exterior.expand (- boundary_width);
      // Variables may have size zero
      // ASSERT_rl (not domain_active.empty(),
      //            "The active part of the domain must not be empty");
      ASSERT_rl (domain_active <= domain_exterior,
                 "The active part of the domain must be contained in the exterior part of the domain");
      
      ibset const domain_boundary = domain_exterior - domain_active;
      
      timer_domain.stop();
      
      
      
      static Timers::Timer timer_region ("region");
      timer_region.start();
      
      for (int c = 0; c < h.components(rl); ++ c) {
        
        full_dboxes & box = full_level.AT(c);
        
        
        
        // Interior:
        
        ibbox & intr = box.interior;
        intr = ibbox::poison();
        
        // The interior of the grid has the extent as specified by the
        // regridding thorn
        intr = h.extent (ml,rl,c);
        
        // (The interior must not be empty)
        // Variables may have size zero
        // ASSERT_c (not intr.empty(),
        //           "The interior must not be empty");
        
        // The interior must be contained in the domain
        ASSERT_c (intr <= h.baseextent(ml,rl),
                  "The interior must be contained in the domain");
        
        // All interiors must be disjunct
#ifdef CARPET_DEBUG
        for (int cc = 0; cc < c; ++ cc) {
          ASSERT_cc (not intr.intersects (full_level.AT(cc).interior),
                     "All interiors must be disjunct");
        }
#endif
        
        
        
        // Outer boundary faces:
        
        b2vect & is_outer_boundary = box.is_outer_boundary;
        
        // The outer boundary faces are where the interior extends up
        // to the outer boundary of the domain.  It is not possible to
        // check whether it extends past the active part of the
        // domain, since this would be wrong when the outer boundary
        // width is zero.
        is_outer_boundary[0] = intr.lower() == domain_exterior.lower(); 
        is_outer_boundary[1] = intr.upper() == domain_exterior.upper(); 
        
        
        
        // Exterior:
        
        ibbox & extr = box.exterior;
        extr = ibbox::poison();
        
        ASSERT_c (all (all (ghost_width >= 0)),
                  "The gh ghost widths must not be negative");
        extr = intr.expand (i2vect (not is_outer_boundary) * ghost_width);
        
        // (The exterior must not be empty)
        // Variables may have size zero
        // ASSERT_c (not extr.empty(),
        //           "The exterior must not be empty");
        
        // The exterior must be contained in the domain
        ASSERT_c (extr <= domain_exterior,
                  "The exterior must be contained in the domain");
        
        
        
        // Cactus ghost zones (which include outer boundaries):
        
        ibset & ghosts = box.ghosts;
        ghosts = ibset::poison();
        
        ghosts = extr - intr;
        
        // The ghosts must be contained in the domain.  Different from
        // the boundaries, the ghost can include part of the outer
        // boundary of the domain.
        ASSERT_c (ghosts <= domain_exterior,
                  "The ghost zones must be contained in the domain");
        
        
        
        // Communicated region:
        
        ibbox & comm = box.communicated;
        comm = ibbox::poison();
        
        comm = extr.expand (i2vect (is_outer_boundary) * (- boundary_width));
        
        // (The communicated region must not be empty)
        // Variables may have size zero
        // ASSERT_c (not comm.empty(),
        //           "The communicated region must not be empty");
        
        // The communicated region must be contained in the active
        // part of the domain
        ASSERT_c (comm <= domain_active,
                  "The communicated region must be contained in the active part of the domain");
        
        
        
        // Outer boundary:
        
        ibset & outer_boundaries = box.outer_boundaries;
        outer_boundaries = ibset::poison();
        
        outer_boundaries = extr - comm;
        
        // The outer boundary must be contained in the outer boundary
        // of the domain
        ASSERT_c (outer_boundaries <= domain_boundary,
                  "The outer boundary must be contained in the outer boundary of the domain");
        
        
        
        // Owned region:
        
        ibbox & owned = box.owned;
        owned = ibbox::poison();
        
        owned = intr.expand (i2vect (is_outer_boundary) * (- boundary_width));
        
        // (The owned region must not be empty)
        // Variables may have size zero
        // ASSERT_c (not owned.empty(),
        //           "The owned region must not be empty");
        
        // The owned region must be contained in the active part of
        // the domain
        ASSERT_c (owned <= domain_active,
                  "The owned region must be contained in the active part of the domain");
        
        // All owned regions must be disjunct
#ifdef CARPET_DEBUG
        for (int cc = 0; cc < c; ++ cc) {
          ASSERT_cc (not owned.intersects (full_level.AT(cc).owned),
                     "All owned regions must be disjunct");
        }
#endif
        
        
        
        // Boundary (Carpet ghost zones, which do not include outer
        // boundaries):
        
        ibset & boundaries = box.boundaries;
        boundaries = ibset::poison();
        
        boundaries = comm - owned;
        
        // The boundary must be contained in the active part of the
        // domain.  This prevents that a region is too close to the
        // outer boundary, so that it has ghost zones overlapping with
        // the outer boundary.
        ASSERT_c (boundaries <= domain_active,
                  "The boundary must be contained in the active part of the domain");
        
      } // for c
      
      timer_region.stop();
      
      
      
      // Conjunction of all buffer zones:
      
      static Timers::Timer timer_buffers ("buffers");
      timer_buffers.start();
      
      // Enlarge active part of domain
      i2vect const safedist = i2vect (0);
      ibbox const domain_enlarged = domain_active.expand (safedist);
      
      // All owned regions
      ibset const allowned (full_level, & full_dboxes::owned);
      ASSERT_rl (allowned <= domain_active,
                 "The owned regions must be contained in the active part of the domain");
      
      // All not-owned regions
      ibset const notowned = domain_enlarged - allowned;
      
      // All not-active points
      ibset const notactive = notowned.expand (buffer_width + overlap_width);
      
      // All not-active points, in stages
      int const num_substeps =
        any (any (ghost_width == 0)) ?
        0 :
        minval (minval (buffer_width / ghost_width));
      if (not all (all (buffer_width == num_substeps * ghost_width))) {
        ostringstream buf;
        buf << "The buffer width " << buffer_width << " is not a multiple of the ghost width " << ghost_width << " on level " << rl;
        CCTK_WARN (CCTK_WARN_COMPLAIN, buf.str().c_str());
      }
      
#ifdef CARPET_DEBUG
      vector<ibset> notactive_stepped (num_substeps+1);
      notactive_stepped.AT(0) = notowned;
      for (int substep = 1; substep <= num_substeps; ++ substep) {
        notactive_stepped.AT(substep) =
          notactive_stepped.AT(substep-1).expand (ghost_width);
      }
      ibset const notactive_overlaps =
        notactive_stepped.AT(num_substeps).expand (overlap_width);
      if (all (all (buffer_width == num_substeps * ghost_width))) {
        ASSERT_rl (notactive_overlaps == notactive,
                   "The stepped not-owned region including overlaps must be equal to the not-active region");
      }
#endif
      
      // All buffer zones
      //ibset const allbuffers = allowned & notowned.expand (buffer_width);
      ibset& allbuffers = level_level.buffers;
      allbuffers = allowned & notowned.expand (buffer_width);
      
      // All overlap zones
      ibset const alloverlaps = (allowned & notactive) - allbuffers;
      
      // All active points
      ibset& allactive = level_level.active;
      allactive = allowned - notactive;
      
#ifdef CARPET_DEBUG
      // All stepped buffer zones
      vector<ibset> allbuffers_stepped (num_substeps);
      ibset allbuffers_stepped_combined;
      for (int substep = 0; substep < num_substeps; ++ substep) {
        allbuffers_stepped.AT(substep) = 
          allowned &
          (notactive_stepped.AT(substep+1) - notactive_stepped.AT(substep));
        allbuffers_stepped_combined += allbuffers_stepped.AT(substep);
      }
      if (all (all (buffer_width == num_substeps * ghost_width))) {
        ASSERT_rl (allbuffers_stepped_combined == allbuffers,
                   "The stepped buffer zones must be equal to the buffer zones");
      }
#endif
      
      // Overlap zones and buffer zones must be in the active part of
      // the domain
      ASSERT_rl (allactive <= domain_active,
                 "The active region must be in the active part of the domain");
      ASSERT_rl (alloverlaps <= domain_active,
                 "The overlap zones must be in the active part of the domain");
      ASSERT_rl (allbuffers <= domain_active,
                 "The buffer zones must be in the active part of the domain");
      ASSERT_rl ((allactive & alloverlaps).empty(), 
                 "The active points and the overlap zones cannot overlap");
      ASSERT_rl ((allactive & allbuffers).empty(), 
                 "The active points and the buffer zones cannot overlap");
      ASSERT_rl ((alloverlaps & allbuffers).empty(), 
                 "The overlap zones and the buffer zones cannot overlap");
      ASSERT_rl (allactive + alloverlaps + allbuffers == allowned,
                 "The active points, the overlap points, and buffer points together must be exactly the owned region");
      
      
      
      for (int c = 0; c < h.components(rl); ++ c) {
        full_dboxes & box = full_level.AT(c);
        
        // Buffer zones:
        box.buffers = box.owned & allbuffers;
        
        // Overlap zones:
        box.overlaps = box.owned & alloverlaps;
        
        // Active region:
        box.active = box.owned & allactive;
        ASSERT_c (box.active == box.owned - (box.buffers + box.overlaps),
                  "The active region must equal the owned region minus the (buffer zones plus overlap zones)");
      } // for c
      
      for (int lc = 0; lc < h.local_components(rl); ++ lc) {
        int const c = h.get_component (rl, lc);
        local_dboxes & local_box = local_level.AT(lc);
        full_dboxes const& box = full_level.AT(c);
        
        local_box.buffers = box.buffers;
        local_box.overlaps = box.overlaps;
        
#if 0
        // Stepped buffer zones:
        local_box.buffers_stepped.resize (num_substeps);
        for (int substep = 0; substep < num_substeps; ++ substep) {
          local_box.buffers_stepped.AT(substep) =
            box.owned & allbuffers_stepped.AT(substep);
        }
#endif
        
        local_box.active = box.active;
      } // for lc
      
      
      
      // The conjunction of all buffer zones must equal allbuffers
      ibset const allbuffers1 (full_level, & full_dboxes::buffers);
      ASSERT_rl (allbuffers1 == allbuffers,
                 "Buffer zone consistency check");
      
      timer_buffers.stop();
      
      
      
      // Test constituency relations:
      
      static Timers::Timer timer_test ("test");
      timer_test.start();
      
      for (int c = 0; c < h.components(rl); ++ c) {
        full_dboxes const & box = full_level.AT(c);
        
        ASSERT_c ((box.active & box.buffers).empty(),
                  "Consistency check");
        ASSERT_c ((box.active & box.overlaps).empty(),
                  "Consistency check");
        ASSERT_c ((box.overlaps & box.buffers).empty(),
                  "Consistency check");
        ASSERT_c ((box.active | box.overlaps | box.buffers) == box.owned,
                  "Consistency check");
        
        ASSERT_c ((box.owned & box.boundaries).empty(),
                  "Consistency check");
        ASSERT_c ((box.owned | box.boundaries) == box.communicated,
                  "Consistency check");
        
        ASSERT_c ((box.communicated & box.outer_boundaries).empty(),
                  "Consistency check");
        ASSERT_c ((box.communicated | box.outer_boundaries) == box.exterior,
                  "Consistency check");
        
        ASSERT_c (box.boundaries <= box.ghosts,
                  "Consistency check");
        
        ASSERT_c ((box.interior & box.ghosts).empty(),
                  "Consistency check");
        ASSERT_c ((box.interior | box.ghosts) == box.exterior,
                  "Consistency check");
        
      } // for c
      
      timer_test.stop();
      
      
      
      // Communication schedule:
      
      static Timers::Timer timer_comm ("comm");
      timer_comm.start();
      
      static Timers::Timer timer_comm_mgrest ("mgrest");
      static Timers::Timer timer_comm_mgprol ("mgprol");
      static Timers::Timer timer_comm_refprol ("refprol");
      static Timers::Timer timer_comm_sync ("sync");
      static Timers::Timer timer_comm_refbndprol ("refbndprol");
      timer_comm_mgrest.instantiate ();
      timer_comm_mgprol.instantiate ();
      timer_comm_refprol.instantiate ();
      timer_comm_sync.instantiate ();
      timer_comm_refbndprol.instantiate ();
      
      for (int lc = 0; lc < h.local_components(rl); ++ lc) {
        int const c = h.get_component (rl, lc);
        
        full_dboxes & box = full_level.AT(c);
        
        
        
        // Multigrid restriction:
        
        timer_comm_mgrest.start();
        
        if (ml > 0) {
          int const oml = ml - 1;
          
          // Multigrid restriction must fill all active points
          
          full_dboxes const & obox = full_boxes.AT(oml).AT(rl).AT(c);
          
          ibset needrecv = box.active;
          
          ibset const contracted_oactive
            (obox.active.contracted_for (box.interior));
          ibset const ovlp = needrecv & contracted_oactive;
          
          for (ibset::const_iterator
                 ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
          {
            ibbox const & recv = * ri;
            ibbox const send = recv.expanded_for (obox.interior);
            ASSERT_c (send <= obox.exterior,
                      "Multigrid restriction: Send region must be contained in exterior");
            fast_level.fast_mg_rest_sendrecv.push_back
              (sendrecv_pseudoregion_t (send, c, recv, c));
          }
          
          needrecv -= ovlp;
          
          // All points must have been received
          ASSERT_c (needrecv.empty(),
                    "Multigrid restriction: All points must have been received");
          
        } // if ml > 0
        
        timer_comm_mgrest.stop();
        
        
        
        // Multigrid prolongation:
        
        timer_comm_mgprol.start();
        
        if (ml > 0) {
          int const oml = ml - 1;
          
          // Multigrid prolongation must fill all active points
          // (this could probably be relaxed)
          
          full_dboxes const & obox = full_boxes.AT(oml).AT(rl).AT(c);
          
          ibset oneedrecv = obox.active;
          
          i2vect const stencil_size = i2vect (prolongation_stencil_size(rl));
          
          ibset const expanded_active (box.active.expanded_for (obox.interior));
          ibset const ovlp = oneedrecv & expanded_active;
          
          for (ibset::const_iterator
                 ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
          {
            ibbox const & recv = * ri;
            ibbox const send =
              recv.expanded_for (box.interior).expand (stencil_size);
            ASSERT_c (send <= box.exterior,
                      "Multigrid prolongation: Send region must be contained in exterior");
            fast_level.fast_mg_prol_sendrecv.push_back
              (sendrecv_pseudoregion_t (send, c, recv, c));
          }
          
          oneedrecv -= ovlp;
          
          // All points must have been received
          ASSERT_c (oneedrecv.empty(),
                    "Multigrid prolongation: All points must have been received");
          
        } // if ml > 0
        
        timer_comm_mgprol.stop();
        
        
        
        // Refinement prolongation:
        
        timer_comm_refprol.start();
        
        if (rl > 0) {
          int const orl = rl - 1;
          
          // Refinement prolongation must fill all active points
          
          ibset needrecv = box.active + box.overlaps;
          
          i2vect const stencil_size = i2vect (prolongation_stencil_size(rl));
          
          ASSERT_c (all (h.reffacts.at(rl) % h.reffacts.at(orl) == 0),
                    "Refinement factors must be integer multiples of each other");
          i2vect const reffact =
            i2vect (h.reffacts.at(rl) / h.reffacts.at(orl));
          
          for (int cc = 0; cc < h.components(orl); ++ cc) {
            full_dboxes const & obox = full_boxes.AT(ml).AT(orl).AT(cc);
            
#if 0
            // This does not work for cell centering; if the domain is
            // 1 cell wide, the contraction disappears it
            ibset const expanded_oactive
              (obox.active.contracted_for (box.interior).expand (reffact));
#else
            ibset const expanded_oactive
              (obox.active.expanded_for (box.interior).expand
               (h.refcent == vertex_centered ? reffact : reffact-1));
#endif
            ibset const ovlp = needrecv & expanded_oactive;
            
            for (ibset::const_iterator
                   ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
            {
              ibbox const & recv = * ri;
              ibbox const send =
                recv.expanded_for (obox.interior).expand (stencil_size);
              ASSERT_c (send <= obox.exterior,
                        "Refinement prolongation: Send region must be contained in exterior");
              fast_level.fast_ref_prol_sendrecv.push_back
                (sendrecv_pseudoregion_t (send, cc, recv, c));
              if (not on_this_proc (orl, cc)) {
                fast_dboxes & fast_level_otherproc =
                  fast_level_otherprocs.AT(this_proc(orl, cc));
                fast_level_otherproc.fast_ref_prol_sendrecv.push_back
                  (sendrecv_pseudoregion_t (send, cc, recv, c));
              }
            }
            
            needrecv -= ovlp;
            
          } // for cc
          
          // All points must have been received
          if (not (needrecv.empty())) {
            cerr << "box.active=" << box.active << "\n"
                 << "needrecv=" << needrecv << "\n";
          }
          ASSERT_c (needrecv.empty(),
                    "Refinement prolongation: All points must have been received");
          
        } // if rl > 0
        
        timer_comm_refprol.stop();
        
        
        
        // Synchronisation:
        
        timer_comm_sync.start();
        
        {
          
          // Synchronisation should fill as many boundary points as
          // possible
          
#if 0
          // Outer boundaries are not synchronised, since they cannot
          // be filled by boundary prolongation either, and therefore
          // the user code must set them anyway.
          ibset needrecv = box.boundaries;
#else
          // Outer boundaries are synchronised for backward
          // compatibility.
          ibset needrecv = box.ghosts;
#endif
          
          ibset & sync = box.sync;
          
          for (int cc = 0; cc < h.components(rl); ++ cc) {
            full_dboxes const & obox = full_level.AT(cc);
            
#if 0
            ibset const ovlp = needrecv & obox.owned;
#else
            ibset const ovlp = needrecv & obox.interior;
#endif
            
            if (cc == c) {
              ASSERT_cc (ovlp.empty(),
                         "A region may not synchronise from itself");
            }
            
            for (ibset::const_iterator
                   ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
            {
              ibbox const & recv = * ri;
              ibbox const & send = recv;
              fast_level.fast_sync_sendrecv.push_back
                (sendrecv_pseudoregion_t (send, cc, recv, c));
              if (not on_this_proc (rl, cc)) {
                fast_dboxes & fast_level_otherproc =
                  fast_level_otherprocs.AT(this_proc(rl, cc));
                fast_level_otherproc.fast_sync_sendrecv.push_back
                  (sendrecv_pseudoregion_t (send, cc, recv, c));
              }
            }
            
            needrecv -= ovlp;
            sync += ovlp;
            
          } // for cc
          
        }
        
        timer_comm_sync.stop();
        
        
        
        // Boundary prolongation:
        
        timer_comm_refbndprol.start();
        
        if (rl > 0) {
          int const orl = rl - 1;
          
#if 0
          // Outer boundaries are not synchronised, since they cannot
          // be filled by boundary prolongation either, and therefore
          // the user code must set them anyway.
          ibset needrecv = box.boundaries;
#else
          // Outer boundaries are synchronised for backward
          // compatibility.
          ibset needrecv = box.ghosts;
#endif
          
          // Points which are synchronised need not be boundary
          // prolongated
          needrecv -= box.sync;
          
          // Outer boundary points cannot be boundary prolongated
          needrecv &= box.communicated;
          
          // Prolongation must fill what cannot be synchronised, and
          // also all buffer zones
          needrecv += box.buffers;
          
          ibset & bndref = box.bndref;
          
          i2vect const stencil_size = i2vect (prolongation_stencil_size(rl));
          
          ASSERT_c (all (h.reffacts.at(rl) % h.reffacts.at(orl) == 0),
                    "Refinement factors must be integer multiples of each other");
          i2vect const reffact =
            i2vect (h.reffacts.at(rl) / h.reffacts.at(orl));
          
          for (int cc = 0; cc < h.components(orl); ++ cc) {
            full_dboxes const & obox = full_boxes.AT(ml).AT(orl).AT(cc);
            
            // See "refinement prolongation"
            ibset const expanded_oactive
              (obox.active.expanded_for (box.interior).expand
               (h.refcent == vertex_centered ? reffact : reffact-1));
            ibset const ovlp = needrecv & expanded_oactive;
            
            for (ibset::const_iterator
                   ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
            {
              ibbox const & recv = * ri;
              ibbox const send =
                recv.expanded_for (obox.interior).expand (stencil_size);
              ASSERT_c (send <= obox.exterior,
                        "Boundary prolongation: Send region must be contained in exterior");
              fast_level.fast_ref_bnd_prol_sendrecv.push_back
                (sendrecv_pseudoregion_t (send, cc, recv, c));
              if (not on_this_proc (orl, cc)) {
                fast_dboxes & fast_level_otherproc =
                  fast_level_otherprocs.AT(this_proc(orl, cc));
                fast_level_otherproc.fast_ref_bnd_prol_sendrecv.push_back
                  (sendrecv_pseudoregion_t (send, cc, recv, c));
              }
            }
            
            needrecv -= ovlp;
            bndref += ovlp;
            
          } // for cc
          
          // All points must now have been received, either through
          // synchronisation or through boundary prolongation
          ASSERT_c (needrecv.empty(),
                    "Synchronisation and boundary prolongation: All points must have been received");
          
        } // if rl > 0
        
        timer_comm_refbndprol.stop();
        
      } // for lc
      
      
      
      // Refinement restriction:
      
      static Timers::Timer timer_comm_refrest ("refrest");
      timer_comm_refrest.start();
      
      if (rl > 0) {
        int const orl = rl - 1;
        full_cboxes const& full_olevel = full_boxes.AT(ml).AT(orl);
        fast_dboxes& fast_olevel = fast_boxes.AT(ml).AT(orl);
        ibbox const& odomext = h.baseextent(ml,orl);
        
        // Refinement restriction may fill all coarse interior points,
        // and must use all fine active points
        
        ibset allrestricted;
        switch (h.refcent) {
        case vertex_centered:
          allrestricted = allactive.contracted_for(odomext);
          break;
        case cell_centered: {
          ibset const& source = allactive;
          ibbox const& target = odomext;
          ibbox const all_source = allactive.container().expand(10,10);
          ibbox const all_target = all_source.contracted_for(target);
          ibset const tmp0 = source;
          ibset const tmp1 = tmp0.expanded_for(target);
          ibset const tmp2 = all_target - tmp1;
          ibset const tmp3 = tmp2.expand(1,1);
          ibset const tmp4 = all_target - tmp3;
          allrestricted = tmp4;
          break;
        }
        default:
          assert(0);
        }
        
        for (int olc = 0; olc < h.local_components(orl); ++ olc) {
          int const oc = h.get_component(orl, olc);
          full_dboxes const& obox = full_olevel.AT(oc);
          
          ibset needrecv = allrestricted & obox.owned;
          if(use_higher_order_restriction) {
            // NOTE: change in behaviour (affects only outer boundaries I think)!!!!
            // NOTE: b/c of this we need a low-level sync after the restrict
            needrecv = allrestricted & obox.interior;
          }
          // Cannot restrict into buffer zones
          assert ((allrestricted & obox.buffers).empty());
          
          for (int c = 0; c < h.components(rl); ++ c) {
            full_dboxes const& box = full_level.AT(c);
            
            // If we make a mistake expanding the domain of dependence here, it
            // _should_ be caught be by the expand()ed is_contained_in(srcbbox)
            // test in the actual operator. 
            int const shrink_by =
              use_higher_order_restriction and (h.refcent == cell_centered) ?
              restriction_order_space/2  : 0;
            ibbox const contracted_exterior = 
              box.exterior.expand(ivect(-shrink_by)).contracted_for(odomext);
            ibset const ovlp = needrecv & contracted_exterior;
            
            for (ibset::const_iterator
                   ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
            {
              ibbox const& recv = *ri;
              ibbox const send = recv.expanded_for(box.exterior).expand(ivect(shrink_by));
              ASSERT_c (send <= box.exterior,
                        "Refinement restriction: Send region must be contained in exterior");
              
              sendrecv_pseudoregion_t const preg (send, c, recv, oc);
              fast_olevel.fast_ref_rest_sendrecv.push_back(preg);
              if (not on_this_proc(rl, c)) {
                fast_dboxes& fast_level_otherproc =
                  fast_level_otherprocs.AT(this_proc(rl, c));
                fast_level_otherproc.fast_ref_rest_sendrecv.push_back(preg);
              }
            }
            
            needrecv -= ovlp;
            
          } // for c
          
            // All points must have been received
          ASSERT_rl (needrecv.empty(),
                     "Refinement restriction: All points must have been received");
          
        } // for olc
        
      } // if rl > 0
      
      timer_comm_refrest.stop();
      
      
      
      // Refinement refluxing:
      
      static Timers::Timer timer_comm_reflux ("reflux");
      timer_comm_reflux.start();
      
      // If there is no coarser level, do nothing
      if (rl > 0) {
        int const orl = rl - 1;
        local_cboxes & local_olevel = local_boxes.AT(ml).AT(orl);
        full_cboxes const& full_olevel = full_boxes.AT(ml).AT(orl);
        fast_dboxes & fast_olevel = fast_boxes.AT(ml).AT(orl);
        
        // This works only with cell centred grids
        if (h.refcent == cell_centered) {
          
          
          
          // Fine grids
          ibset const& fine_level = allactive;
          
          assert (all (h.reffacts.AT(rl) % h.reffacts.AT(orl) == 0));
          
          vect<vect<ibset,2>,dim> all_fine_boundary;
          
          for (int dir=0; dir<dim; ++dir) {
            // Unit vector
            ivect const idir = ivect::dir(dir);
            
            // Left and right faces
            ibset const fine_level_minus = fine_level.shift (-idir, 2);
            ibset const fine_level_plus  = fine_level.shift (+idir, 2);
            
            // Fine boundaries
            all_fine_boundary[dir][0] = fine_level_minus - fine_level_plus ;
            all_fine_boundary[dir][1] = fine_level_plus  - fine_level_minus;
          } // for dir
          
          // Coarse grids
          ibbox const& coarse_domain_exterior = h.baseextent(ml, orl);
          ibbox const& coarse_ext = coarse_domain_exterior;
          ibbox const& fine_domain_exterior = h.baseextent(ml, rl);
          ibbox const& fine_ext = fine_domain_exterior;
          
          // Coarse equivalent of fine boundary
          vect<vect<ibset,2>,dim> all_coarse_boundary;
          for (int dir=0; dir<3; ++dir) {
            // Unit vector
            ivect const idir = ivect::dir(dir);
            ibbox const coarse_faces = coarse_ext.shift(-idir,2);
            ibbox const fine_faces = fine_ext.shift(-idir,2);
            for (int face=0; face<2; ++face) {
              if (verbose) {
                cout << "REFREF rl=" << rl << " dir=" << dir << " face=" << face << "\n"
                     << "   all_fine_boundary=" << all_fine_boundary[dir][face] << "\n"
                     << "   coarse_ext=" << coarse_ext << "\n"
                     << "   coarse_faces=" << coarse_faces << "\n";
              }
              all_coarse_boundary[dir][face] =
                all_fine_boundary[dir][face].contracted_for (coarse_faces);
              if (verbose) {
                cout << "   all_coarse_boundary=" << all_coarse_boundary[dir][face] << "\n";
              }
              ASSERT_rl (all_coarse_boundary[dir][face].expanded_for(fine_faces) == all_fine_boundary[dir][face],
                         "Refluxing: Coarse grid boundaries must be consistent with fine grid boundaries");
            } // for face
          }   // for dir
          
          
          
          // Split onto components
          for (int lc = 0; lc < h.local_components(rl); ++ lc) {
            int const c = h.get_component(rl, lc);
            full_dboxes & box = full_level.AT(c);
            local_dboxes & local_box = local_level.AT(lc);
            
            for (int dir = 0; dir < dim; ++ dir) {
              // Unit vector
              ivect const idir = ivect::dir(dir);
              for (int face = 0; face < 2; ++ face) {
                
                // This is not really used (only for debugging)
                local_box.fine_boundary[dir][face] =
                  box.exterior.shift(-idir,2) & all_fine_boundary[dir][face];
                if (verbose) {
                  cout << "REFREF rl=" << rl << " lc=" << lc << " dir=" << dir << " face=" << face << "\n"
                       << "   local.fine_boundary=" << local_box.fine_boundary[dir][face] << "\n";
                }
                
#if 0                           // OFFSETS,SIZE
                ibset const& boxes =
                  local_box.fine_boundary[dir][face];
                vector<int>& offsets =
                  local_box.fine_boundary_offsets[dir][face];
                assert(offsets.empty());
                offsets.reserve(boxes.setsize());
                int offset = 0;
                for (ibset::const_iterator
                       bi=boxes.begin(); bi!=boxes.end(); ++bi)
                {
                  offsets.push_back(offset);
                  offset += (*bi).size();
                }
                local_box.fine_boundary_size[dir][face] = offset;
#endif
                
              } // for face
            }   // for dir
            
          } // for lc
          
#ifdef CARPET_DEBUG
          for (int dir = 0; dir < dim; ++ dir) {
            ivect const idir = ivect::dir(dir); // Unit vector
            for (int face = 0; face < 2; ++ face) {
              ibset all_fine_boundary_combined;
              for (int c = 0; c < h.components(rl); ++ c) {
                full_dboxes const & box = full_level.AT(c);
                all_fine_boundary_combined +=
                  box.exterior.shift(-idir,2) & all_fine_boundary[dir][face];
              } // for c
              ASSERT_rl (all_fine_boundary_combined == all_fine_boundary[dir][face],
                         "Refluxing: Union of all fine grid boundaries must be the total fine grid boundary ");
            }
          }
#endif
          
          for (int lc = 0; lc < h.local_components(orl); ++ lc) {
            int const c = h.get_component(orl, lc);
            full_dboxes const & obox = full_olevel.AT(c);
            local_dboxes & local_obox = local_olevel.AT(lc);
            
            for (int dir = 0; dir < dim; ++ dir) {
              // Unit vector
              ivect const idir = ivect::dir(dir);
              for (int face = 0; face < 2; ++ face) {
                
                // This is used for post-processing the fluxes
                // (calculating the difference between coarse and fine
                // fluxes, adjusting the state)
                local_obox.coarse_boundary[dir][face] =
                  obox.owned.shift(-idir,2) & all_coarse_boundary[dir][face];
                // Cannot reflux into buffer zones
                if (verbose) {
                  cout << "REFREF orl=" << orl << " lc=" << lc << " dir=" << dir << " face=" << face << "\n"
                       << "   local.coarse_boundary=" << local_obox.coarse_boundary[dir][face] << "\n";
                }
                assert ((obox.buffers.shift(-idir,2) & all_coarse_boundary[dir][face]).empty());
                
#if 0                           // OFFSETS,SIZE
                ibset const& boxes =
                  local_obox.coarse_boundary[dir][face];
                vector<int>& offsets =
                  local_obox.coarse_boundary_offsets[dir][face];
                assert(offsets.empty());
                offsets.reserve(boxes.setsize());
                int offset = 0;
                for (ibset::const_iterator
                       bi=boxes.begin(); bi!=boxes.end(); ++bi)
                {
                  offsets.push_back(offset);
                  offset += (*bi).size();
                }
                local_obox.coarse_boundary_size[dir][face] = offset;
#endif
                
              } // for face
            }   // for dir
            
          } // for lc
          
#ifdef CARPET_DEBUG
          for (int dir = 0; dir < dim; ++ dir) {
            ivect const idir = ivect::dir(dir); // Unit vector
            for (int face = 0; face < 2; ++ face) {
              ibset all_coarse_boundary_combined;
              for (int c = 0; c < h.components(orl); ++ c) {
                full_dboxes const & obox = full_olevel.AT(c);
                all_coarse_boundary_combined +=
                  obox.owned.shift(face ? +idir : -idir,2) & all_coarse_boundary[dir][face];
              } // for c
              ASSERT_rl (all_coarse_boundary_combined == all_coarse_boundary[dir][face],
                         "Refluxing: Union of all coarse grid boundaries must be the total coarse grid boundary ");
            }
          }
#endif
          
          
          
          // Set up communication schedule
          for (int olc = 0; olc < h.local_components(orl); ++ olc) {
            int const oc = h.get_component(orl, olc);
            local_dboxes & local_obox = local_olevel.AT(olc);
            if (verbose) {
              cout << "REF ref_refl ml=" << ml << " rl=" << rl << " olc=" << olc << " oc=" << oc << "\n";
            }
            
            for (int dir = 0; dir < dim; ++ dir) {
              // Unit vector
              ivect const idir = ivect::dir(dir);
              ibbox const coarse_faces = coarse_ext.shift(-idir,2);
              for (int face = 0; face < 2; ++ face) {
                
                srpvect fast_dboxes::* const fast_ref_refl_sendrecv =
                  fast_dboxes::fast_ref_refl_sendrecv[dir][face];
                srpvect fast_dboxes::* const fast_ref_refl_prol_sendrecv =
                  fast_dboxes::fast_ref_refl_prol_sendrecv[dir][face];
                
                // Refluxing must fill all coarse refluxing boundary
                // points, and may use all fine points
                
                ibset needrecv (local_obox.coarse_boundary[dir][face]);
                
                for (int c = 0; c < h.components(rl); ++ c) {
                  full_dboxes const & box = full_level.AT(c);
                  if (verbose) {
                    cout << "REF ref_refl ml=" << ml << " rl=" << rl << " olc=" << olc << " oc=" << oc << " c=" << c << " dir=" << dir << " face=" << face << "\n";
                  }
                  
                  ibbox const contracted_exterior =
                    box.exterior.shift(-idir, 2).contracted_for(coarse_faces);
                  if (verbose) {
                    cout << "   exterior=" << box.exterior.shift(-idir, 2) << "\n"
                         << "   contracted=" << contracted_exterior << "\n";
                  }
                  ibset const ovlp = needrecv & contracted_exterior;
                  if (verbose) {
                    cout << "   ovlp=" << ovlp << "\n";
                  }
                  
                  for (ibset::const_iterator
                         ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
                  {
                    ibbox const & recv = * ri;
                    ibbox const send =
                      recv.expanded_for (box.exterior.shift(-idir, 2));
                    ASSERT_c (send <= box.exterior.shift(-idir, 2),
                              "Refinement restriction: Send region must be contained in exterior");
                    ibbox const shifted_recv = recv.shift(idir, 2);
                    ibbox const shifted_send = send.shift(idir, 2);
                    sendrecv_pseudoregion_t const preg
                      (shifted_send, c, shifted_recv, oc);
                    if (verbose) {
                      cout << "REF ref_refl ml=" << ml << " rl=" << rl << " olc=" << olc << " c=" << c << " oc=" << oc << " dir=" << dir << " face=" << face << "\n"
                           << "   preg=" << preg << "\n";
                    }
                    (fast_olevel.*fast_ref_refl_sendrecv).push_back (preg);
                    if (not on_this_proc (rl, c)) {
                      fast_dboxes & fast_level_otherproc =
                        fast_level_otherprocs.AT(this_proc(rl, c));
                      (fast_level_otherproc.*fast_ref_refl_sendrecv).
                        push_back (preg);
                    }
                    
                    // reflux-prolongation
                    sendrecv_pseudoregion_t const preg_prol
                      (shifted_recv, oc, shifted_send, c);
                    if (verbose) {
                      cout << "REF ref_refL_prol ml=" << ml << " rl=" << rl << " olc=" << olc << " c=" << c << " oc=" << oc << " dir=" << dir << " face=" << face << "\n"
                           << "   preg=" << preg_prol << "\n";
                    }
                    (fast_level.*fast_ref_refl_prol_sendrecv).push_back
                      (preg_prol);
                    if (not on_this_proc (orl, oc)) {
                      fast_dboxes & fast_level_otherproc =
                        fast_level_otherprocs.AT(this_proc(orl, oc));
                      (fast_level_otherproc.*fast_ref_refl_prol_sendrecv).
                        push_back (preg_prol);
                    }
                  }
                  
                  needrecv -= ovlp;
                  
                } // for c
                
                // All points must have been received
                if (not needrecv.empty()) {
                  cout << "needrecv=" << needrecv << endl;
                }
                ASSERT_rl (needrecv.empty(),
                           "Refinement refluxing: All points must have been received");
                
              } // for face
            }   // for dir
            
          } // for olc
          
        } // if cell centered
        
      } // if rl > 0
      
      timer_comm_reflux.stop();
      
      timer_comm.stop();
      
      
      
      // Reduction mask:
      
      static Timers::Timer timer_mask ("mask");
      timer_mask.start();
      
      for (int lc=0; lc<h.local_components(rl); ++lc) {
        local_dboxes& local_box = local_level.AT(lc);
        local_box.prolongation_boundary.clear();
        local_box.restriction_boundary .clear();
        local_box.prolongation_boundary.resize(num_nbs());
        local_box.restriction_boundary .resize(num_nbs());
      }
      
      if (rl > 0) {
        int const orl = rl - 1;
        
        // Set the weight in the interior of the notactive and the
        // allactive regions to zero, and set the weight on the
        // boundary of the notactive and allactive regions to 1/2.
        //
        // For the prolongation region, the "boundary" is the first
        // layer outside of the region, and the "region" consists of
        // the inactive points. This corresponds to the outermost
        // layer of active grid points.
        //
        // For the restricted region, the "boundary" is the outermost
        // layer of grid points if this layer is aligned with the next
        // coarser (i.e. the current) grid; otherwise, the boundary is
        // empty. The "region" consists of the active points.
        
        // Note: We calculate the prolongation information for the
        // current level, and the restriction information for the next
        // coarser level. We do it this way because both calculations
        // start from the current level's set of active grid points.
        
        full_cboxes const& full_olevel = full_boxes.AT(ml).AT(orl);
        // Local boxes are not communicated or post-processed, and
        // thus can be modified even for coarser levels
        local_cboxes& local_olevel = local_boxes.AT(ml).AT(orl);
        
        if (verbose) {
          ostringstream buf;
          buf << "Setting prolongation region " << notactive << " on level " << rl;
          CCTK_INFO (buf.str().c_str());
        }
        if (verbose) {
          ostringstream buf;
          buf << "Setting restriction region " << allactive << " on level " << orl;
          CCTK_INFO (buf.str().c_str());
        }
        
        // ibset test_boxes;
        // ibset test_cfboxes;
        
        for (int neighbour=0; neighbour<num_nbs(); ++neighbour) {
          ivect const shift = ind2nb(neighbour);
          
          // In this loop, shift [1,1,1] denotes a convex corner of
          // the region which should be masked out, i.e. a region
          // where only a small part (1/8) of the region should be
          // masked out. Concave edges are treated implicitly
          // (sequentially), i.e. larger chunks are cut out multiple
          // times: a concave xy edge counts as both an x face and a y
          // face.
          
          // Prolongation boundary of this level: start with inactive
          // (prolongated from the next coarser level) points
          ibset boxes  = notactive;
          // Restriction boundary of the next coarser level: start
          // with this level's active (restricted to the next coarser)
          // points
          ibset fboxes = allactive;
          
          switch (h.refcent) {
          case vertex_centered:
            for (int d=0; d<dim; ++d) {
              ivect const dir = ivect::dir(d);
              fboxes = fboxes.shift(-dir) & fboxes & fboxes.shift(+dir);
            }
            for (int d=0; d<dim; ++d) {
              // Calculate the boundary in direction d
              ivect const dir = ivect::dir(d);
              switch (shift[d]) {
              case -1: {
                // left boundary
                boxes  = boxes .shift(-dir) - boxes;
                fboxes = fboxes.shift(-dir) - fboxes;
                break;
              }
              case 0: {
                // interior
                // do nothing
                break;
              }
              case +1: {
                // right boundary
                boxes  = boxes .shift(+dir) - boxes;
                fboxes = fboxes.shift(+dir) - fboxes;
                break;
              }
              default:
                assert (0);
              }
            }
            break;
          case cell_centered:
            // For cell-centred grids we assume that all cell
            // boundaries are aligned. This may not actually be the
            // case.
            // TODO: assert this.
            if (all(shift == 0)) {
              // do nothing
            } else {
              boxes  = ibset();
              fboxes = ibset();
            }
            break;
          default:
            assert (0);
          } // switch centering
          
          ibbox const& odomext = h.baseextent(ml, orl);
          ibset const cfboxes = fboxes.contracted_for(odomext);
          // test_boxes   |= boxes;
          // test_cfboxes |= cfboxes;
          
          if (verbose) {
            ostringstream buf;
            buf << "Setting boundary " << shift << ": prolongation region " << boxes;
            CCTK_INFO (buf.str().c_str());
          }
          if (verbose) {
            ostringstream buf;
            buf << "Setting boundary " << shift << ": restriction region " << fboxes;
            CCTK_INFO (buf.str().c_str());
          }
          
          // For the current level, the light boxes do not (yet) exist
          // here, so we use the full boxes
          for (int lc=0; lc<h.local_components(rl); ++lc) {
            int const c = h.get_component(rl, lc);
            full_dboxes const& full_box = full_level.AT(c);
            local_dboxes& local_box = local_level.AT(lc);
            local_box.prolongation_boundary.AT(neighbour) =
              boxes & full_box.exterior;
          }
          // For the next coarser level, the full boxes do not exist
          // (any more), so we use the light boxes
          for (int olc=0; olc<h.local_components(orl); ++olc) {
            int const oc = h.get_component(orl, olc);
            full_dboxes const& full_obox = full_olevel.AT(oc);
            local_dboxes& local_obox = local_olevel.AT(olc);
            local_obox.restriction_boundary.AT(neighbour) =
              cfboxes & full_obox.exterior;
          }
          
        } // for neighbour
        
      } // if not coarsest level
      
      timer_mask.stop();
      
      
      
      // Mask for unused points on coarser level (which do not
      // influence the future evolution provided regridding is done at
      // the right times):
      static Timers::Timer timer_overwrittenmask ("unusedpoints_mask");
      timer_overwrittenmask.start();
      
      // Declare this here to save it for later use. Contains all the boxes
      // which are active minus the boundary
      ibset all_refined;
      
      if (rl > 0) {
        int const orl = rl - 1;
        full_cboxes const& full_olevel = full_boxes.AT(ml).AT(orl);
        // Local boxes are not communicated or post-processed, and
        // thus can be modified even for coarser levels
        local_cboxes& local_olevel = local_boxes.AT(ml).AT(orl);
        
        ivect const reffact = h.reffacts.AT(rl) / h.reffacts.AT(orl);
        
        // This works only when the refinement factor is 2
        assert (all (reffact == 2));
        
        ibbox const& base = domain_exterior;
        ibbox const& obase = h.baseextent(ml,orl);
        
        // Calculate the union of all coarse regions
        ibset const allointr (full_olevel, & full_dboxes::interior);
        
        // Project to current level
        ivect const rf(reffact);
        ibset const parent (allointr.expanded_for(base));
        
        // Subtract the active region
        ibset const notrefined = parent - allactive;
        
        // Enlarge this set
        vect<vect<ibset,2>,dim> enlarged;
        for (int d=0; d<dim; ++d) {
          for (int f=0; f<2; ++f) {
            switch (h.refcent) {
            case vertex_centered: {
              ivect const dir = ivect::dir(d);
              enlarged[d][f] =
                ibset (notrefined.expand(f==1?dir:0, f==0?dir:0));
              break;
            }
            case cell_centered: {
              enlarged[d][f] = notrefined;
              // TODO: restriction boundaries are wrong (they are
              // empty, but should not be) with cell centring when
              // fine cell cut coarse cells
              bool const old_there_was_an_error = there_was_an_error;
              ASSERT_rl (notrefined.contracted_for(obase).expanded_for(base) ==
                         notrefined,
                         "Refinement mask: Fine grid boundaries must be aligned with coarse grid cells");
              // Do not abort because of this problem
              there_was_an_error = old_there_was_an_error;
              break;
            }
            default:
              assert (0);
            }
          }
        }
        
        // Intersect with the active region
        vect<vect<ibset,2>,dim> all_boundaries;
        for (int d=0; d<dim; ++d) {
          for (int f=0; f<2; ++f) {
            all_boundaries[d][f] = allactive & enlarged[d][f];
          }
        }
        
        // TODO: Ensure that the prolongation boundaries
        // all_boundaries are contained in the boundary prolongated
        // region
        
        // TODO: Ensure that the restriction boundaries and the
        // restricted region are contained in the restricted region
        
        // Subtract the boundaries from the refined region
        all_refined = allactive;
        for (int d=0; d<dim; ++d) {
          for (int f=0; f<2; ++f) {
            all_refined -= all_boundaries[d][f];
          }
        }
        
        for (int lc = 0; lc < h.local_components(orl); ++ lc) {
          int const c = h.get_component(orl, lc);
          full_dboxes const& obox = full_olevel.AT(c);
          // Local boxes are not communicated or post-processed, and
          // thus can be modified even for coarser levels
          local_dboxes & local_obox = local_olevel.AT(lc);
          
          // Set restriction information for next coarser level
          local_obox.restricted_region =
            all_refined.contracted_for(obox.exterior) & obox.owned;
        } // for lc
        
      } // if not coarsest level
      
      if (rl > 0) {
        int const orl = rl - 1;
        full_cboxes const& full_olevel = full_boxes.AT(ml).AT(orl);
        // Local boxes are not communicated or post-processed, and
        // thus can be modified even for coarser levels
        local_cboxes& local_olevel = local_boxes.AT(ml).AT(orl);

        // This works only when the refinement factor is 2
        ivect const reffact = h.reffacts.AT(rl) / h.reffacts.AT(orl);
        if (all (reffact == 2)) {
          // use the already computed 'all_refined' to get region from where
          // no information will be used later (overwritten)
          // First: get the region which will get restricted, on the coarse level
          ibset restricted_region = all_refined.contracted_for(h.baseextent(ml,orl));
          // This is too big - during MoL-substeps information within this
          // region will be used to update points outside -> need to
          // shrink it by some points
          // The way we shrink it is to invert it, expand that, and invert
          // again. To invert we define an enclosing box and subtract it from that.
          i2vect to_shrink = buffer_widths[orl] + ghost_widths[orl];
          ibbox enclosing = restricted_region.container().expand(ivect(1)+to_shrink);
          ibset unused_region = enclosing - (enclosing - restricted_region).expand(to_shrink);
          // Now we have the interesting region in 'unused_region' and need to store
          // the union of this and the local regions
          for (int lc = 0; lc < h.local_components(orl); ++ lc) {
            int const c = h.get_component(orl, lc);
            full_dboxes const& obox = full_olevel.AT(c);
            // Local boxes are not communicated or post-processed, and
            // thus can be modified even for coarser levels
            local_dboxes & local_obox = local_olevel.AT(lc);
            // Set unused information for next coarser level
            local_obox.unused_region = unused_region & obox.owned;
          } // for lc
        } // if reffact != 2
      } // if not coarsest level

      timer_overwrittenmask.stop();
      
      
      
      // Regridding schedule:
      
      fast_level.do_init = do_init;
      if (do_init) {
        
        static Timers::Timer timer_regrid ("regrid");
        timer_regrid.start();
        
        static Timers::Timer timer_regrid_sync ("sync");
        static Timers::Timer timer_regrid_prolongate ("prolongate");
        
        for (int lc = 0; lc < h.local_components(rl); ++ lc) {
          int const c = h.get_component (rl, lc);
          
          full_dboxes & box = full_level.AT(c);
          
          ibset needrecv = box.active + box.overlaps;
          
          
          
          // Regridding synchronisation:
          
          timer_regrid_sync.start();
          
          if (int (old_light_boxes.size()) > ml and
              int (old_light_boxes.AT(ml).size()) > rl)
          {
            
            int const oldcomponents = old_light_boxes.AT(ml).AT(rl).size();
            
            // Synchronisation copies from the same level of the old
            // grid structure.  It should fill as many active points
            // as possible.
            
            for (int cc = 0; cc < oldcomponents; ++ cc) {
              light_dboxes const & obox = old_light_boxes.AT(ml).AT(rl).AT(cc);
              
              ibset const ovlp = needrecv & obox.owned;
              
              for (ibset::const_iterator
                     ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
              {
                ibbox const & recv = * ri;
                ibbox const & send = recv;
                fast_level.fast_old2new_sync_sendrecv.push_back
                  (sendrecv_pseudoregion_t (send, cc, recv, c));
                if (not on_this_oldproc (rl, cc)) {
                  fast_dboxes & fast_level_otherproc =
                    fast_level_otherprocs.AT(this_oldproc(rl, cc));
                  fast_level_otherproc.fast_old2new_sync_sendrecv.push_back
                    (sendrecv_pseudoregion_t (send, cc, recv, c));
                }
              }
              
              needrecv -= ovlp;
              
            } // for cc
            
          } // if not old_light_boxes.empty
          
          timer_regrid_sync.stop();
          
          
          
          // Regridding prolongation:
          
          timer_regrid_prolongate.start();
          
          if (rl > 0) {
            int const orl = rl - 1;
            
            // Prolongation interpolates from the next coarser level
            // of the new grid structure.  It must fill what cannot be
            // synchronised.
            
            i2vect const stencil_size = i2vect (prolongation_stencil_size(rl));
            
            ASSERT_c (all (h.reffacts.at(rl) % h.reffacts.at(orl) == 0),
                      "Refinement factors must be integer multiples of each other");
            i2vect const reffact =
              i2vect (h.reffacts.at(rl) / h.reffacts.at(orl));
            
            for (int cc = 0; cc < h.components(orl); ++ cc) {
              full_dboxes const & obox = full_boxes.AT(ml).AT(orl).AT(cc);
              
              // See "refinement prolongation"
              ibset const expanded_oactive
                (obox.active.expanded_for (box.interior).expand
                 (h.refcent == vertex_centered ? reffact : reffact-1));
              ibset const ovlp = needrecv & expanded_oactive;
              
              for (ibset::const_iterator
                     ri = ovlp.begin(); ri != ovlp.end(); ++ ri)
              {
                ibbox const & recv = * ri;
                ibbox const send =
                  recv.expanded_for (obox.interior).expand (stencil_size);
                if (not (send <= obox.exterior)) {
                  cerr << "stencil_size=" << stencil_size << "\n"
                       << "obox=" << obox
                       << "recv=" << recv << "\n"
                       << "send=" << send << "\n";
                }
                ASSERT_c (send <= obox.exterior,
                          "Regridding prolongation: Send region must be contained in exterior");
                fast_level.fast_old2new_ref_prol_sendrecv.push_back
                  (sendrecv_pseudoregion_t (send, cc, recv, c));
                if (not on_this_proc (orl, cc)) {
                  fast_dboxes & fast_level_otherproc =
                    fast_level_otherprocs.AT(this_proc(orl, cc));
                  fast_level_otherproc.fast_old2new_ref_prol_sendrecv.
                    push_back (sendrecv_pseudoregion_t (send, cc, recv, c));
                }
              }
              
              needrecv -= ovlp;
              
            } // for cc
            
          } // if rl > 0
          
          if (int (old_light_boxes.size()) > ml and
              int (old_light_boxes.AT(ml).size()) > 0)
          {
            // All points must now have been received, either through
            // synchronisation or through prolongation
            ASSERT_c (needrecv.empty(),
                      "Regridding prolongation: All points must have been received");
          }
          
          timer_regrid_prolongate.stop();
          
        } // for lc
        
        timer_regrid.stop();
        
      } // if do_init
      
      
      
      for (int lc = 0; lc < h.local_components(rl); ++ lc) {
        int const c = h.get_component (rl, lc);
        
        light_level.AT(c).exterior = full_level.AT(c).exterior;
        light_level.AT(c).owned    = full_level.AT(c).owned;
        light_level.AT(c).interior = full_level.AT(c).interior;
        
        light_level.AT(c).exterior_size = full_level.AT(c).exterior.size();
        light_level.AT(c).owned_size    = full_level.AT(c).owned.size();
        light_level.AT(c).active_size   = full_level.AT(c).active.size();
        
      } // for lc
      
      
      
      // Broadcast grid structure and communication schedule
      
      {
        
        static Timers::Timer timer_bcast_boxes ("bcast_boxes");
        timer_bcast_boxes.start();
        
        int const count_send = h.local_components(rl);
        vector<light_dboxes> level_send (count_send);
        for (int lc = 0; lc < h.local_components(rl); ++ lc) {
          int const c = h.get_component (rl, lc);
          level_send.AT(lc) = light_level.AT(c);
        }
        vector<vector<light_dboxes> > const level_recv =
          allgatherv (dist::comm(), level_send);
        vector<int> count_recv (dist::size(), 0);
        for (int c = 0; c < h.components(rl); ++ c) {
          int const p = this_proc (rl, c);
          if (p != dist::rank()) {
            light_level.AT(c) = level_recv.AT(p).AT(count_recv.AT(p));
            ++ count_recv.AT(p);
          }
        }
        for (int p = 0; p < dist::size(); ++ p) {
          if (p != dist::rank()) {
            assert (count_recv.AT(p) == int(level_recv.AT(p).size()));
          }
        }
        
        timer_bcast_boxes.stop();
        
      }
      
      {
        
        static Timers::Timer timer_bcast_comm ("bcast_comm");
        timer_bcast_comm.start();
        
        static Timers::Timer timer_bcast_comm_ref_prol ("ref_prol");
        timer_bcast_comm_ref_prol.start();
        broadcast_schedule (fast_level_otherprocs, fast_level,
                            & fast_dboxes::fast_ref_prol_sendrecv);
        timer_bcast_comm_ref_prol.stop();
        
        static Timers::Timer timer_bcast_comm_sync ("sync");
        timer_bcast_comm_sync.start();
        broadcast_schedule (fast_level_otherprocs, fast_level,
                            & fast_dboxes::fast_sync_sendrecv);
        timer_bcast_comm_sync.stop();
        
        static Timers::Timer timer_bcast_comm_ref_bnd_prol ("ref_bnd_prol");
        timer_bcast_comm_ref_bnd_prol.start();
        broadcast_schedule (fast_level_otherprocs, fast_level,
                            & fast_dboxes::fast_ref_bnd_prol_sendrecv);
        timer_bcast_comm_ref_bnd_prol.stop();
        
        if (rl > 0) {
          int const orl = rl - 1;
          fast_dboxes & fast_olevel = fast_boxes.AT(ml).AT(orl);
          static Timers::Timer timer_bcast_comm_ref_rest ("ref_rest");
          timer_bcast_comm_ref_rest.start();
          broadcast_schedule (fast_level_otherprocs, fast_olevel,
                              & fast_dboxes::fast_ref_rest_sendrecv);
          timer_bcast_comm_ref_rest.stop();
        }
        
        if (rl > 0) {
          int const orl = rl - 1;
          fast_dboxes & fast_olevel = fast_boxes.AT(ml).AT(orl);
          static Timers::Timer timer_bcast_comm_ref_refl ("ref_refl");
          timer_bcast_comm_ref_refl.start();
          for (int dir = 0; dir < dim; ++ dir) {
            for (int face = 0; face < 2; ++ face) {
              srpvect fast_dboxes::* const fast_ref_refl_sendrecv =
                fast_dboxes::fast_ref_refl_sendrecv[dir][face];
              broadcast_schedule (fast_level_otherprocs, fast_olevel,
                                  fast_ref_refl_sendrecv);
            }
          }
          timer_bcast_comm_ref_refl.stop();
        }
        
        static Timers::Timer timer_bcast_comm_ref_refl_prol
          ("ref_refl_prol");
        timer_bcast_comm_ref_refl_prol.start();
        for (int dir = 0; dir < dim; ++ dir) {
          for (int face = 0; face < 2; ++ face) {
            srpvect fast_dboxes::* const fast_ref_refl_prol_sendrecv =
              fast_dboxes::fast_ref_refl_prol_sendrecv[dir][face];
            broadcast_schedule (fast_level_otherprocs, fast_level,
                                fast_ref_refl_prol_sendrecv);
          }
        }
        timer_bcast_comm_ref_refl_prol.stop();
        
        // TODO: Maybe broadcast old2new schedule only if do_init is
        // set
        static Timers::Timer timer_bcast_comm_old2new_sync ("old2new_sync");
        timer_bcast_comm_old2new_sync.start();
        broadcast_schedule (fast_level_otherprocs, fast_level,
                            & fast_dboxes::fast_old2new_sync_sendrecv);
        timer_bcast_comm_old2new_sync.stop();
        
        static Timers::Timer timer_bcast_comm_old2new_ref_prol
          ("old2new_ref_prol");
        timer_bcast_comm_old2new_ref_prol.start();
        broadcast_schedule (fast_level_otherprocs, fast_level,
                            & fast_dboxes::fast_old2new_ref_prol_sendrecv);
        timer_bcast_comm_old2new_ref_prol.stop();
        
        timer_bcast_comm.stop();
        
      }
      
      
      
      // Output:
      if (output_bboxes or there_was_an_error) {
        
        cout << eol;
        cout << "ml=" << ml << " rl=" << rl << eol;
        cout << "baseextent=" << h.baseextent(ml,rl) << eol;
        
        for (int c = 0; c < h.components(rl); ++ c) {
          cout << eol;
          cout << "ml=" << ml << " rl=" << rl << " c=" << c << eol;
          cout << "extent=" << h.extent(ml,rl,c) << eol;
          cout << "outer_boundaries=" << h.outer_boundaries(rl,c) << eol;
          cout << "processor=" << h.processor(rl,c) << eol;
        } // for c
        
        for (int c = 0; c < h.components(rl); ++ c) {
          full_dboxes const & box = full_boxes.AT(ml).AT(rl).AT(c);
          cout << eol;
          cout << "ml=" << ml << " rl=" << rl << " c=" << c << eol;
          cout << box;
        } // for c
        
        // light, local, and fast boxes are output later
        
      } // if output_bboxes
      
      
      
      // Free memory early to save space
      if (int (old_light_boxes.size()) > ml and
          int (old_light_boxes.AT(ml).size()) > rl)
      {
        old_light_boxes.AT(ml).AT(rl).clear();
      }
      
      if (ml > 0) {
        if (rl > 0) {
          full_boxes.AT(ml-1).AT(rl-1).clear();
        }
        if (rl == h.reflevels()-1) {
          full_boxes.AT(ml-1).AT(rl).clear();
        }
      }
      if (ml == h.mglevels()-1) {
        if (rl > 0) {
          full_boxes.AT(ml).AT(rl-1).clear();
        }
        if (rl == h.reflevels()-1) {
          full_boxes.AT(ml).AT(rl).clear();
        }
      }
      
    } // for rl
    
    if (ml > 0) {
      full_boxes.AT(ml-1).clear();
    }
    if (ml == h.mglevels()-1) {
      full_boxes.AT(ml).clear();
    }
    
  } // for ml
  
  
  
  // Output:
  if (output_bboxes or there_was_an_error) {
    
    for (int ml = 0; ml < h.mglevels(); ++ ml) {
      for (int rl = 0; rl < h.reflevels(); ++ rl) {
        
        cout << eol;
        cout << "ml=" << ml << " rl=" << rl << eol;
        cout << "baseextent=" << h.baseextent(ml,rl) << eol;
        
        for (int c = 0; c < h.components(rl); ++ c) {
          cout << eol;
          cout << "ml=" << ml << " rl=" << rl << " c=" << c << eol;
          cout << "extent=" << h.extent(ml,rl,c) << eol;
          cout << "outer_boundaries=" << h.outer_boundaries(rl,c) << eol;
          cout << "processor=" << h.processor(rl,c) << eol;
        } // for c
        
        // full boxes have already been output (and deallocated)
    
        for (int c = 0; c < h.components(rl); ++ c) {
          light_dboxes const & box = light_boxes.AT(ml).AT(rl).AT(c);
          cout << eol;
          cout << "ml=" << ml << " rl=" << rl << " c=" << c << eol;
          cout << box;
        } // for c
        
        for (int lc = 0; lc < h.local_components(rl); ++ lc) {
          int const c = h.get_component (rl, lc);
          local_dboxes const & box = local_boxes.AT(ml).AT(rl).AT(lc);
          cout << eol;
          cout << "ml=" << ml << " rl=" << rl << " lc=" << lc << " c=" << c << eol;
          cout << box;
        } // for lc
        
        {
          level_dboxes const & box = level_boxes.AT(ml).AT(rl);
          cout << eol;
          cout << "ml=" << ml << " rl=" << rl << eol;
          cout << box;
        }
        
        fast_dboxes const & fast_box = fast_boxes.AT(ml).AT(rl);
        cout << eol;
        cout << "ml=" << ml << " rl=" << rl << eol;
        cout << fast_box;
        
      } // for rl
    }   // for ml
    
    cout << eol;
    cout << "memoryof(gh)=" << memoryof(h) << eol;
    cout << "memoryof(dh)=" << memoryof(*this) << eol;
    cout << "memoryof(dh.light_boxes)=" << memoryof(light_boxes) << eol;
    cout << "memoryof(dh.local_boxes)=" << memoryof(local_boxes) << eol;
    cout << "memoryof(dh.level_boxes)=" << memoryof(level_boxes) << eol;
    cout << "memoryof(dh.fast_boxes)=" << memoryof(fast_boxes) << eol;
    int gfcount = 0;
    size_t gfmemory = 0;
    for (map<int,ggf*>::const_iterator gfi = gfs.begin(); gfi != gfs.end(); ++ gfi)
    {
      ++ gfcount;
      gfmemory += memoryof(*gfi->second);
    }
    cout << "#gfs=" << gfcount << eol;
    cout << "memoryof(gfs)=" << gfmemory << eol;
    
  } // if output_bboxes
  
  if (there_was_an_error) {
    CCTK_WARN (CCTK_WARN_ABORT,
               "The grid structure is inconsistent.  It is impossible to continue.");
  }
  
  
  
  total.stop (0);
  timer.stop();
}



void
dh::
broadcast_schedule (vector<fast_dboxes> & fast_level_otherprocs,
                    fast_dboxes & fast_level,
                    srpvect fast_dboxes::* const schedule_item)
{
  static Timers::Timer timer_bs1 ("CarpetLib::dh::bs1");
  timer_bs1.start();
  vector <srpvect> send (dist::size());
  for (int p=0; p<dist::size(); ++p) {
    swap (send.AT(p), fast_level_otherprocs.AT(p).*schedule_item);
  }
  timer_bs1.stop();
  
  static Timers::Timer timer_bs2 ("CarpetLib::dh::bs2");
  timer_bs2.start();
  srpvect const recv = alltoallv1 (dist::comm(), send);
  timer_bs2.stop();
  
  static Timers::Timer timer_bs3 ("CarpetLib::dh::bs3");
  timer_bs3.start();
  (fast_level.*schedule_item).insert
    ((fast_level.*schedule_item).end(), recv.begin(), recv.end());
  timer_bs3.stop();
}



void
dh::
regrid_free (bool const do_init)
{
  if (do_init) {
    for (int ml = 0; ml < h.mglevels(); ++ ml) {
      for (int rl = 0; rl < h.reflevels(); ++ rl) {
        fast_boxes.AT(ml).AT(rl).fast_old2new_sync_sendrecv.clear();
        fast_boxes.AT(ml).AT(rl).fast_old2new_ref_prol_sendrecv.clear();
      }
    }
  } else {
    for (int ml = 0; ml < h.mglevels(); ++ ml) {
      for (int rl = 0; rl < h.reflevels(); ++ rl) {
        assert (fast_boxes.AT(ml).AT(rl).fast_old2new_sync_sendrecv.empty());
        assert (fast_boxes.AT(ml).AT(rl).fast_old2new_ref_prol_sendrecv.empty());
      }
    }
  }
}



void
dh::
recompose (int const rl, bool const do_prolongate)
{
  DECLARE_CCTK_PARAMETERS;
  
  assert (rl>=0 and rl<h.reflevels());
  
  static Timers::Timer timer ("CarpetLib::dh::recompose");
  timer.start ();
  
  for (map<int,ggf*>::iterator f=gfs.begin(); f!=gfs.end(); ++f) {
    f->second->recompose_crop ();
  }
  
  if (combine_recompose) {
    // Recompose all grid functions of this refinement levels at once.
    // This may be faster, but requires more memory.
    for (map<int,ggf*>::iterator f=gfs.begin(); f!=gfs.end(); ++f) {
      f->second->recompose_allocate (rl);
    }
    // TODO: If this works, rename do_prolongate to do_init here, and
    // remove the do_prolongate parameter from ggf::recompose_fill
#if 0
    for (comm_state state; not state.done(); state.step()) {
      for (map<int,ggf*>::iterator f=gfs.begin(); f!=gfs.end(); ++f) {
        f->second->recompose_fill (state, rl, do_prolongate);
      }
    }
#endif
    if (do_prolongate) {
      for (comm_state state; not state.done(); state.step()) {
        for (map<int,ggf*>::iterator f=gfs.begin(); f!=gfs.end(); ++f) {
          f->second->recompose_fill (state, rl, true);
        }
      }
    }
    for (map<int,ggf*>::iterator f=gfs.begin(); f!=gfs.end(); ++f) {
      f->second->recompose_free_old (rl);
    }
  } else {
    // Recompose the grid functions sequentially.  This may be slower,
    // but requires less memory.  This is the default.
    for (map<int,ggf*>::iterator f=gfs.begin(); f!=gfs.end(); ++f) {
      f->second->recompose_allocate (rl);
#if 0
      for (comm_state state; not state.done(); state.step()) {
        f->second->recompose_fill (state, rl, do_prolongate);
      }
#endif
      if (do_prolongate) {
        for (comm_state state; not state.done(); state.step()) {
          f->second->recompose_fill (state, rl, true);
        }
      }
      f->second->recompose_free_old (rl);
    }
  }
  
  timer.stop ();
}



// Grid function management
void
dh::
insert (ggf * const f)
{
  CHECKPOINT;
  assert (f);
  assert (f->varindex >= 0);
  assert (f->varindex < CCTK_NumVars());
  const bool inserted =
    gfs.insert (pair<int,ggf*> (f->varindex, f)).second;
  assert (inserted);
}

void
dh::
erase (ggf * const f)
{
  CHECKPOINT;
  assert (f);
  assert (f->varindex >= 0);
  assert (f->varindex < CCTK_NumVars());
  const size_t erased = gfs.erase (f->varindex);
  assert (erased == 1);
}



// Equality

bool
dh::full_dboxes::
operator== (full_dboxes const & b) const
{
  return
    exterior         == b.exterior         and
    all(all(is_outer_boundary == b.is_outer_boundary)) and
    outer_boundaries == b.outer_boundaries and
    communicated     == b.communicated     and
    boundaries       == b.boundaries       and
    owned            == b.owned            and
    buffers          == b.buffers          and
    overlaps         == b.overlaps         and
    active           == b.active           and
    sync             == b.sync             and
    bndref           == b.bndref           and
    ghosts           == b.ghosts           and
    interior         == b.interior;
}



// MPI datatypes

MPI_Datatype
mpi_datatype (dh::light_dboxes const &)
{
  static bool initialised = false;
  static MPI_Datatype newtype;
  if (not initialised) {
    static dh::light_dboxes s;
#define ENTRY(type, name)                                               \
    {                                                                   \
      sizeof s.name / sizeof(type), /* count elements */                \
        (char*)&s.name - (char*)&s, /* offsetof doesn't work (why?) */  \
        dist::mpi_datatype<type>(), /* find MPI datatype */             \
        STRINGIFY(name),    /* field name */                            \
        STRINGIFY(type),    /* type name */                             \
        }
    dist::mpi_struct_descr_t const descr[] = {
      ENTRY(int, exterior),
      ENTRY(int, owned),
      ENTRY(int, interior),
      ENTRY(dh::light_dboxes::size_type, exterior_size),
      ENTRY(dh::light_dboxes::size_type, owned_size),
      ENTRY(dh::light_dboxes::size_type, active_size),
      {1, sizeof s, MPI_UB, "MPI_UB", "MPI_UB"}
    };
#undef ENTRY
    newtype =
      dist::create_mpi_datatype (sizeof descr / sizeof descr[0], descr,
                                 "dh::light::dboxes", sizeof s);
#if 0
    int type_size;
    MPI_Type_size (newtype, & type_size);
    assert (type_size <= sizeof s);
    MPI_Aint type_lb, type_ub;
    MPI_Type_lb (newtype, & type_lb);
    MPI_Type_ub (newtype, & type_ub);
    assert (type_ub - type_lb == sizeof s);
#endif
    initialised = true;
  }
  return newtype;
}

MPI_Datatype
mpi_datatype (dh::fast_dboxes const &)
{
  static bool initialised = false;
  static MPI_Datatype newtype;
  if (not initialised) {
    static dh::fast_dboxes s;
#define ENTRY(type, name)                                               \
    {                                                                   \
      sizeof s.name / sizeof(type), /* count elements */                \
        (char*)&s.name - (char*)&s, /* offsetof doesn't work (why?) */  \
        dist::mpi_datatype<type>(), /* find MPI datatype */             \
        STRINGIFY(name),    /* field name */                            \
        STRINGIFY(type),    /* type name */                             \
        }
    dist::mpi_struct_descr_t const descr[] = {
      ENTRY (dh::srpvect, fast_mg_rest_sendrecv),
      ENTRY (dh::srpvect, fast_mg_prol_sendrecv),
      ENTRY (dh::srpvect, fast_ref_prol_sendrecv),
      ENTRY (dh::srpvect, fast_ref_rest_sendrecv),
      ENTRY (dh::srpvect, fast_sync_sendrecv),
      ENTRY (dh::srpvect, fast_ref_bnd_prol_sendrecv),
      ENTRY (dh::srpvect, fast_old2new_sync_sendrecv),
      ENTRY (dh::srpvect, fast_old2new_ref_prol_sendrecv),
      ENTRY (dh::srpvect, fast_ref_refl_sendrecv_0_0),
      ENTRY (dh::srpvect, fast_ref_refl_sendrecv_0_1),
      ENTRY (dh::srpvect, fast_ref_refl_sendrecv_1_0),
      ENTRY (dh::srpvect, fast_ref_refl_sendrecv_1_1),
      ENTRY (dh::srpvect, fast_ref_refl_sendrecv_2_0),
      ENTRY (dh::srpvect, fast_ref_refl_sendrecv_2_1),
      ENTRY (dh::srpvect, fast_ref_refl_prol_sendrecv_0_0),
      ENTRY (dh::srpvect, fast_ref_refl_prol_sendrecv_0_1),
      ENTRY (dh::srpvect, fast_ref_refl_prol_sendrecv_1_0),
      ENTRY (dh::srpvect, fast_ref_refl_prol_sendrecv_1_1),
      ENTRY (dh::srpvect, fast_ref_refl_prol_sendrecv_2_0),
      ENTRY (dh::srpvect, fast_ref_refl_prol_sendrecv_2_1),
      {1, sizeof s, MPI_UB, "MPI_UB", "MPI_UB"}
    };
#undef ENTRY
    newtype =
      dist::create_mpi_datatype (sizeof descr / sizeof descr[0], descr,
                                 "dh::fast_dboxes", sizeof s);
    initialised = true;
  }
  return newtype;
}



// Memory usage

size_t
dh::
memory ()
  const
{
  return
    memoryof (ghost_widths) +
    memoryof (buffer_widths) +
    memoryof (overlap_widths) +
    memoryof (prolongation_orders_space) +
    memoryof (light_boxes) +
    memoryof (local_boxes) +
    memoryof (level_boxes) +
    memoryof (fast_boxes) +
    memoryof (gfs);
}

size_t
dh::
allmemory ()
{
  size_t mem = memoryof(alldh);
  for (set<dh*>::const_iterator dhi = alldh.begin(); dhi != alldh.end(); ++ dhi)
  {
    mem += memoryof(**dhi);
  }
  return mem;
}

size_t
dh::light_dboxes::
memory ()
  const
{
  return
    memoryof (exterior) +
    memoryof (owned) +
    memoryof (interior) +
    memoryof (exterior_size) +
    memoryof (owned_size) +
    memoryof (active_size);
}

size_t
dh::local_dboxes::
memory ()
  const
{
  return
    memoryof (buffers) +
    memoryof (overlaps) +
    memoryof (active) +
#if 0
    memoryof (buffers_stepped) +
#endif
#if 0
    memoryof (fine_active) +
#endif
    memoryof (prolongation_boundary) +
    memoryof (restriction_boundary) +
    memoryof (restricted_region) +
    memoryof (unused_region) +
    memoryof (coarse_boundary) +
    memoryof (fine_boundary);
#if 0                           // OFFSETS,SIZE
    memoryof (coarse_boundary_offsets) +
    memoryof (fine_boundary_offsets) +
    memoryof (coarse_boundary_size) +
    memoryof (fine_boundary_size);
#endif
}

size_t
dh::level_dboxes::
memory ()
  const
{
  return
    memoryof (active);
}

size_t
dh::full_dboxes::
memory ()
  const
{
  return
    memoryof (exterior) +
    memoryof (is_outer_boundary) +
    memoryof (outer_boundaries) +
    memoryof (communicated) +
    memoryof (boundaries) +
    memoryof (owned) +
    memoryof (buffers) +
    memoryof (overlaps) +
    memoryof (active) +
    memoryof (sync) +
    memoryof (bndref) +
    memoryof (ghosts) +
    memoryof (interior);
}

size_t
dh::fast_dboxes::
memory ()
  const
{
  return
    memoryof (fast_mg_rest_sendrecv) +
    memoryof (fast_mg_prol_sendrecv) +
    memoryof (fast_ref_prol_sendrecv) +
    memoryof (fast_ref_rest_sendrecv) +
    memoryof (fast_sync_sendrecv) +
    memoryof (fast_ref_bnd_prol_sendrecv) +
    memoryof (fast_ref_refl_sendrecv_0_0) +
    memoryof (fast_ref_refl_sendrecv_0_1) +
    memoryof (fast_ref_refl_sendrecv_1_0) +
    memoryof (fast_ref_refl_sendrecv_1_1) +
    memoryof (fast_ref_refl_sendrecv_2_0) +
    memoryof (fast_ref_refl_sendrecv_2_1) +
    memoryof (fast_ref_refl_prol_sendrecv_0_0) +
    memoryof (fast_ref_refl_prol_sendrecv_0_1) +
    memoryof (fast_ref_refl_prol_sendrecv_1_0) +
    memoryof (fast_ref_refl_prol_sendrecv_1_1) +
    memoryof (fast_ref_refl_prol_sendrecv_2_0) +
    memoryof (fast_ref_refl_prol_sendrecv_2_1) +
    memoryof (do_init) +
    memoryof (fast_old2new_sync_sendrecv) +
    memoryof (fast_old2new_ref_prol_sendrecv);
}



// Input

istream &
dh::light_dboxes::
input (istream & is)
{
  // Regions:
  try {
    skipws (is);
    consume (is, "dh::light_dboxes:{");
    skipws (is);
    consume (is, "exterior:");
    is >> exterior;
    exterior_size = exterior.size();
    skipws (is);
    consume (is, "owned:");
    is >> owned;
    owned_size = owned.size();
    skipws (is);
    consume (is, "interior:");
    is >> interior;
    skipws (is);
    consume (is, "active_size:");
    is >> active_size;
    skipws (is);
    consume (is, "}");
  } catch (input_error & err) {
    cout << "Input error while reading a dh::light_dboxes" << endl;
    throw err;
  }
  return is;
}

istream &
dh::local_dboxes::
input (istream & is)
{
  // Regions:
  try {
    skipws (is);
    consume (is, "dh::local_dboxes:{");
    skipws (is);
    consume (is, "buffers:");
    is >> buffers;
    skipws (is);
    consume (is, "overlaps:");
    is >> overlaps;
    skipws (is);
    consume (is, "active:");
    is >> active;
#if 0
    skipws (is);
    consume (is, "buffers_stepped:");
    is >> buffers_stepped;
#endif
#if 0
    skipws (is);
    consume (is, "fine_active:");
    is >> fine_active;
#endif
    skipws (is);
    consume (is, "prolongation_boundary:");
    is >> prolongation_boundary;
    skipws (is);
    consume (is, "restriction_boundary:");
    is >> restriction_boundary;
    skipws (is);
    consume (is, "restricted_region:");
    is >> restricted_region;
    skipws (is);
    consume (is, "unused_region:");
    is >> unused_region;
    skipws (is);
    consume (is, "coarse_boundary:");
    is >> coarse_boundary;
    skipws (is);
    consume (is, "fine_boundary:");
    is >> fine_boundary;
    skipws (is);
    // TODO: read boundary sizes and boundary offsets
    consume (is, "}");
  } catch (input_error & err) {
    cout << "Input error while reading a dh::local_dboxes" << endl;
    throw err;
  }
  return is;
}

istream &
dh::level_dboxes::
input (istream & is)
{
  // Regions:
  try {
    skipws (is);
    consume (is, "dh::level_dboxes:{");
    skipws (is);
    consume (is, "active:");
    is >> active;
    skipws (is);
    consume (is, "}");
  } catch (input_error & err) {
    cout << "Input error while reading a dh::level_dboxes" << endl;
    throw err;
  }
  return is;
}

istream &
dh::full_dboxes::
input (istream & is)
{
  // Regions:
  try {
    skipws (is);
    consume (is, "dh::full_dboxes:{");
    skipws (is);
    consume (is, "exterior:");
    is >> exterior;
    skipws (is);
    consume (is, "is_outer_boundary:");
    is >> is_outer_boundary;
    skipws (is);
    consume (is, "outer_boundaries:");
    is >> outer_boundaries;
    skipws (is);
    consume (is, "communicated:");
    is >> communicated;
    skipws (is);
    consume (is, "boundaries:");
    is >> boundaries;
    skipws (is);
    consume (is, "owned:");
    is >> owned;
    skipws (is);
    consume (is, "buffers:");
    is >> buffers;
    skipws (is);
    consume (is, "overlaps:");
    is >> overlaps;
    skipws (is);
    consume (is, "active:");
    is >> active;
    skipws (is);
    consume (is, "sync:");
    is >> sync;
    skipws (is);
    consume (is, "bndref:");
    is >> bndref;
    skipws (is);
    consume (is, "ghosts:");
    is >> ghosts;
    skipws (is);
    consume (is, "interior:");
    is >> interior;
    skipws (is);
    consume (is, "}");
  } catch (input_error & err) {
    cout << "Input error while reading a dh::full_dboxes" << endl;
    throw err;
  }
  return is;
}

istream &
dh::fast_dboxes::
input (istream & is)
{
  // Communication schedule:
  try {
    skipws (is);
    consume (is, "dh::fast_dboxes:{");
    skipws (is);
    consume (is, "fast_mg_rest_sendrecv:");
    is >> fast_mg_rest_sendrecv;
    skipws (is);
    consume (is, "fast_mg_prol_sendrecv:");
    is >> fast_mg_prol_sendrecv;
    skipws (is);
    consume (is, "fast_ref_prol_sendrecv:");
    is >> fast_ref_prol_sendrecv;
    skipws (is);
    consume (is, "fast_ref_rest_sendrecv:");
    is >> fast_ref_rest_sendrecv;
    skipws (is);
    consume (is, "fast_sync_sendrecv:");
    is >> fast_sync_sendrecv;
    skipws (is);
    consume (is, "fast_ref_bnd_prol_sendrecv:");
    is >> fast_ref_bnd_prol_sendrecv;
    skipws (is);
    consume (is, "fast_old2new_sync_sendrecv:");
    is >> fast_old2new_sync_sendrecv;
    skipws (is);
    consume (is, "fast_old2new_ref_prol_sendrecv:");
    is >> fast_old2new_ref_prol_sendrecv;
    skipws (is);
    consume (is, "}");
  } catch (input_error & err) {
    cout << "Input error while reading a dh::fast_dboxes" << endl;
    throw err;
  }
  return is;
}



// Output

ostream &
dh::
output (ostream & os)
  const
{
  os << "dh:"
     << "ghost_widths=" << ghost_widths << ","
     << "buffer_widths=" << buffer_widths << ","
     << "overlap_widths=" << overlap_widths << ","
     << "prolongation_orders_space=" << prolongation_orders_space << ","
     << "light_boxes=" << light_boxes << ","
     << "local_boxes=" << local_boxes << ","
     << "level_boxes=" << level_boxes << ","
     << "fast_boxes=" << fast_boxes << ","
     << "gfs={";
  {
    bool isfirst = true;
    for (map<int,ggf*>::const_iterator
           f = gfs.begin(); f != gfs.end(); ++ f, isfirst = false)
    {
      if (not isfirst) os << ",";
      os << *f;
    }
  }
  os << "}";
  return os;
}

ostream &
dh::light_dboxes::
output (ostream & os)
  const
{
  // Regions:
  os << "dh::light_dboxes:{" << eol
     << "   exterior: " << exterior << eol
     << "   owned: " << owned << eol
     << "   interior: " << interior << eol
     << "   active_size: " << active_size << eol
     << "}" << eol;
  return os;
}

ostream &
dh::local_dboxes::
output (ostream & os)
  const
{
  // Regions:
  os << "dh::local_dboxes:{" << eol
     << "   buffers: " << buffers << eol
     << "   overlaps: " << overlaps << eol
     << "   active: " << active << eol
#if 0
     << "   buffers_stepped: " << buffers_stepped << eol
#endif
#if 0
     << "   fine_active: " << fine_active << eol
#endif
     << "   prolongation_boundary: " << prolongation_boundary << eol
     << "   restriction_boundary: " << restriction_boundary << eol
     << "   restricted_region: " << restricted_region << eol
     << "   unused_region: " << unused_region << eol
     << "   coarse_boundary: " << coarse_boundary << eol
     << "   fine_boundary: " << fine_boundary << eol
#if 0                           // OFFSETS,SIZE
     << "   coarse_boundary_offsets: " << coarse_boundary_offsets << eol
     << "   fine_boundary_offsets: " << fine_boundary_offsets << eol
     << "   coarse_boundary_size: " << coarse_boundary_size << eol
     << "   fine_boundary_size: " << fine_boundary_size << eol
#endif
     << "}" << eol;
  return os;
}

ostream &
dh::level_dboxes::
output (ostream & os)
  const
{
  // Regions:
  os << "dh::level_dboxes:{" << eol
     << "   active: " << active << eol
     << "}" << eol;
  return os;
}

ostream &
dh::full_dboxes::
output (ostream & os)
  const
{
  // Regions:
  os << "dh::full_dboxes:{" << eol
     << "   exterior: " << exterior << eol
     << "   is_outer_boundary: " << is_outer_boundary << eol
     << "   outer_boundaries: " << outer_boundaries << eol
     << "   communicated: " << communicated << eol
     << "   boundaries: " << boundaries << eol
     << "   owned: " << owned << eol
     << "   buffers: " << buffers << eol
     << "   overlaps: " << overlaps << eol
     << "   active: " << active << eol
     << "   sync: " << sync << eol
     << "   bndref: " << bndref << eol
     << "   ghosts: " << ghosts << eol
     << "   interior: " << interior << eol
     << "}" << eol;
  return os;
}

ostream &
dh::fast_dboxes::
output (ostream & os)
  const
{
  // Communication schedule:
  os << "dh::fast_dboxes:{" << eol
     << "   fast_mg_rest_sendrecv: " << fast_mg_rest_sendrecv << eol
     << "   fast_mg_prol_sendrecv: " << fast_mg_prol_sendrecv << eol
     << "   fast_ref_prol_sendrecv: " << fast_ref_prol_sendrecv << eol
     << "   fast_ref_rest_sendrecv: " << fast_ref_rest_sendrecv << eol
     << "   fast_sync_sendrecv: " << fast_sync_sendrecv << eol
     << "   fast_ref_bnd_prol_sendrecv: " << fast_ref_bnd_prol_sendrecv << eol
     << "   fast_ref_refl_sendrecv_0_0: " << fast_ref_refl_sendrecv_0_0 << eol
     << "   fast_ref_refl_sendrecv_0_1: " << fast_ref_refl_sendrecv_0_1 << eol
     << "   fast_ref_refl_sendrecv_1_0: " << fast_ref_refl_sendrecv_1_0 << eol
     << "   fast_ref_refl_sendrecv_1_1: " << fast_ref_refl_sendrecv_1_1 << eol
     << "   fast_ref_refl_sendrecv_2_0: " << fast_ref_refl_sendrecv_2_0 << eol
     << "   fast_ref_refl_sendrecv_2_1: " << fast_ref_refl_sendrecv_2_1 << eol
     << "   fast_ref_refl_prol_sendrecv_0_0: " << fast_ref_refl_prol_sendrecv_0_0 << eol
     << "   fast_ref_refl_prol_sendrecv_0_1: " << fast_ref_refl_prol_sendrecv_0_1 << eol
     << "   fast_ref_refl_prol_sendrecv_1_0: " << fast_ref_refl_prol_sendrecv_1_0 << eol
     << "   fast_ref_refl_prol_sendrecv_1_1: " << fast_ref_refl_prol_sendrecv_1_1 << eol
     << "   fast_ref_refl_prol_sendrecv_2_0: " << fast_ref_refl_prol_sendrecv_2_0 << eol
     << "   fast_ref_refl_prol_sendrecv_2_1: " << fast_ref_refl_prol_sendrecv_2_1 << eol
     << "   do_init: " << do_init << eol
     << "   fast_old2new_sync_sendrecv: " << fast_old2new_sync_sendrecv << eol
     << "   fast_old2new_ref_prol_sendrecv: " << fast_old2new_ref_prol_sendrecv << eol
     << "}" << eol;
  return os;
}