aboutsummaryrefslogtreecommitdiff
path: root/Carpet/Carpet/src/Recompose.cc
blob: 030ae6f76f062a2fccdff8855cacee118350494f (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
#include <cctk.h>
#include <cctk_Parameters.h>

#include <algorithm>
#include <cassert>
#include <cmath>
#include <cstdio>
#include <cstdlib>
#include <cstring>
#include <fstream>
#include <iomanip>
#include <iostream>
#include <limits>
#include <list>
#include <sstream>
#include <string>
#include <vector>

#include <sys/stat.h>
#include <sys/types.h>

#ifdef CCTK_MPI
#  include <mpi.h>
#else
#  include "nompi.h"
#endif

#include <loopcontrol.h>

#include <balance.hh>
#include <bbox.hh>
#include <bboxset.hh>
#include <defs.hh>
#include <dh.hh>
#include <gh.hh>
#include <region.hh>
#include <vect.hh>

#include <carpet.hh>
#include <modes.hh>
#include <variables.hh>
#include <Timers.hh>



namespace Carpet {
  
  using namespace std;
  
  
  
  // Helper routines for spliting regions automatically
  
  // The cost for a region, assuming a cost of 1 per interior point
  static rvect
  cost (region_t const & reg)
  {
    DECLARE_CCTK_PARAMETERS;
    static rvect costfactor;
    static bool initialised = false;
    if (not initialised) {
      costfactor = rvect(1.0);
      if (dim > 0) costfactor[0] = 1.0 / aspect_ratio_x;
      if (dim > 1) costfactor[1] = 1.0 / aspect_ratio_y;
      if (dim > 2) costfactor[2] = 1.0 / aspect_ratio_z;
    }
    if (reg.extent.empty()) return rvect(0);
    return rvect (reg.extent.shape() / reg.extent.stride()) * costfactor;
  }
  
  
  
  static void
  ClassifyPoints (cGH const * cctkGH, int rl);
  
  
  
  static void
  SplitRegionsMaps_Automatic_Recursively (bvect    const   & dims,
                                          int      const     firstproc,
                                          int      const     nprocs,
                                          region_t         & superreg,
                                          vector<region_t> & newregs);
  static void
  SplitRegions_AsSpecified (cGH const * cctkGH,
                            vector<region_t> & superregs,
                            vector<region_t> & regs);
  
  
  
  void
  CheckRegions (gh::mregs const & regsss)
  {
    char const * const where = "CheckRegions";
    static Timer timer (where);
    timer.start();
    
    // At least one multigrid level
    if (regsss.size() == 0) {
      CCTK_WARN (0, "I cannot set up a grid hierarchy with zero multigrid levels.");
    }
    assert (regsss.size() > 0);
    // At least one refinement level
    if (regsss.AT(0).size() == 0) {
      CCTK_WARN (0, "I cannot set up a grid hierarchy with zero refinement levels.");
    }
    assert (regsss.AT(0).size() > 0);
    // At most maxreflevels levels
    if ((int)regsss.AT(0).size() > maxreflevels) {
      CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
                  "I cannot set up a grid hierarchy with more than Carpet::max_refinement_levels refinement levels.  I found Carpet::max_refinement_levels=%d, while %d levels were requested.",
                  (int)maxreflevels, (int)regsss.AT(0).size());
    }
    assert ((int)regsss.AT(0).size() <= maxreflevels);
    for (int ml=0; ml<(int)regsss.size(); ++ml) {
      int num_regions = 0;
      for (int rl=0; rl<(int)regsss.AT(0).size(); ++rl) {
        // No empty levels
        // (but allow some empty maps)
        // assert (regsss.AT(ml).AT(rl).size() > 0);
        num_regions += regsss.AT(ml).AT(rl).size();
        for (int c=0; c<(int)regsss.AT(ml).AT(rl).size(); ++c) {
          // Check sizes
          // (but allow processors with zero grid points)
          // assert (all(regsss.AT(rl).AT(c).AT(ml).extent.lower() <=
          //             regsss.AT(rl).AT(c).AT(ml).extent.upper()));
          // Check strides
          ivect const str =
            (maxspacereflevelfact / spacereffacts.AT(rl) * ipow(mgfact, ml));
          assert (all(regsss.AT(ml).AT(rl).AT(c).extent.stride() % str == 0));
          // Check alignments
          assert (all(regsss.AT(ml).AT(rl).AT(c).extent.lower() % str == 0));
          assert (all(regsss.AT(ml).AT(rl).AT(c).extent.upper() % str == 0));
        }
      }
      // No empty levels
      assert (num_regions > 0);
    }
    
    timer.stop();
  }
  
  

  bool
  Regrid (cGH const * const cctkGH,
          bool const force_recompose,
          bool const do_init)
  {
    DECLARE_CCTK_PARAMETERS;
    
    char const * const where = "Regrid";
    static Timer timer (where);
    timer.start();
    
    Checkpoint ("Regridding level %d...", reflevel);
    
    assert (is_level_mode());
    
    bool const have_regrid     = CCTK_IsFunctionAliased ("Carpet_Regrid");
    bool const have_regridmaps = CCTK_IsFunctionAliased ("Carpet_RegridMaps");
    bool const use_regridmaps = regrid_in_level_mode and have_regridmaps;
    
    if (not use_regridmaps and not have_regrid) {
      static bool didtell = false;
      if (maxreflevels > 1 and not didtell) {
	CCTK_WARN (2, "The regridding routine Carpet_Regrid has not been provided.  There will be no regridding.  Maybe you forgot to activate a regridding thorn?");
	didtell = true;
      }
      timer.stop();
      return false;
    }
    
    if (regrid_in_level_mode and not have_regridmaps) {
      static bool didtell = false;
      if (maps > 1 and maxreflevels > 1 and not didtell) {
	CCTK_WARN (2, "The regridding routine Carpet_RegridMaps has not been provided.  Regridding will be performed in singlemap mode instead of level mode.");
	didtell = true;
      }
   }
    
    
    
    bool did_change = false;
    
    if (not use_regridmaps) {
      
      BEGIN_MAP_LOOP(cctkGH, CCTK_GF) {
        
        gh::rregs superregss = vhh.AT(map)->superregions;
        gh::mregs regsss;
        
        // Check whether to recompose
        CCTK_INT const do_recompose =
          Carpet_Regrid (cctkGH, & superregss, & regsss, force_recompose);
        assert (do_recompose >= 0);
        did_change = did_change or do_recompose;
        
        if (do_recompose) {
          RegridMap (cctkGH, map, superregss, regsss, do_init);
        }
        
      } END_MAP_LOOP;
      
    } else {
      
      vector<gh::rregs> superregsss (maps);
      vector<gh::mregs> regssss (maps);
      BEGIN_MAP_LOOP(cctkGH, CCTK_GF) {
        superregsss.AT(map) = vhh.AT(map)->superregions;
      } END_MAP_LOOP;
      
      // Check whether to recompose
      CCTK_INT const do_recompose =
        Carpet_RegridMaps (cctkGH, & superregsss, & regssss, force_recompose);
      assert (do_recompose >= 0);
      // TODO
#if 1 // #ifdef CARPET_DEBUG
      {
        int ival = do_recompose;
        MPI_Bcast (& ival, 1, MPI_INT, 0, dist::comm());
        assert (ival == do_recompose);
      }
#endif
      did_change = did_change or do_recompose;
      
      if (do_recompose) {
        BEGIN_MAP_LOOP(cctkGH, CCTK_GF) {
          gh::rregs const & superregss = superregsss.AT(map);
          gh::mregs const & regsss = regssss.AT(map);
          RegridMap (cctkGH, map, superregss, regsss, do_init);
        } END_MAP_LOOP;
      }
      
    } // if regrid in level mode
    
    
    
    if (did_change) {
      
      PostRegrid (cctkGH);
      
    } // if did change
    
    
      
    timer.stop();
    
    if (enable_no_storage) {
      CCTK_WARN (CCTK_WARN_ALERT,
                 "Carpet completed its internal setup, and would now normally go on to allocate memory.  Since the parameter Carpet::enable_no_storage has been set, Carpet will exit instead.");
      CCTK_Exit ((cGH*) cctkGH, 0);
    }
    
    return did_change;
  }
  
  
  
  void
  RegridMap (cGH const * const cctkGH,
             int const m,
             gh::rregs const & superregss,
             gh::mregs const & regsss,
             bool const do_init)
  {
    DECLARE_CCTK_PARAMETERS;
    
    char const * const where = "RegridMap";
    static Timer timer (where);
    timer.start();
    
    Waypoint ("Regridding map %d...", m);
    
    // TODO: keep levels fixed here
#if 0
    //
    // Keep this level fixed if it is not evolved
    //
    if (regrid_only_evolved_levels) {
      int type;
      bool const use_tapered_grids =
        * static_cast<CCTK_INT const *>
        (CCTK_ParameterGet ("use_tapered_grids", "Carpet", &type));
      assert (type == PARAMETER_BOOLEAN);
      
      int const do_every =
        use_tapered_grids ?
        maxtimereflevelfact / timereffacts.AT(max(0,rl-1)):
        maxtimereflevelfact / timereffacts.AT(      rl   );
      
      bool const regrid_this_level =
        (cctkGH->cctk_iteration - 1) % do_every == 0;
      
      if (not regrid_this_level) {
        // Set regions from current grid structure
        regions.AT(rl) = ...;
      }
    }
#endif
    
    // Check the regions
    CheckRegions (regsss);
    // TODO: check also that the current and all coarser levels did
    // not change
    
    // Regrid
    vhh.AT(m)->regrid (superregss, regsss, do_init);
    
    // Output grid structure to screen
    OutputSuperregions (cctkGH, m, * vhh.AT(m), superregss);
    OutputGrids (cctkGH, m, * vhh.AT(m), * vdd.AT(m));
    
    // Write grid structure to file
    OutputGridStructure (cctkGH, m, regsss);
    OutputGridCoordinates (cctkGH, m, regsss);
    
    Waypoint ("Done regridding map %d.", m);
    timer.stop();
  }
  
  
  
  void
  PostRegrid (cGH const * const cctkGH)
  {
    DECLARE_CCTK_PARAMETERS;
    
    char const * const where = "PostRegrid";
    static Timer timer (where);
    timer.start();
    
    // Calculate new number of levels
    int const oldreflevels = reflevels;
    reflevels = vhh.AT(0)->reflevels();
    for (int m=0; m<maps; ++m) {
      assert (vhh.AT(m)->reflevels() == reflevels);
    }
    
    // TODO
#if 1 // #ifdef CARPET_DEBUG
      {
        // All processes must use the same number of levels
        int ival = reflevels;
        MPI_Bcast (& ival, 1, MPI_INT, 0, dist::comm());
        assert (ival == reflevels);
      }
#endif
    
    // One cannot switch off the current level
    assert (reflevels > reflevel);
    
    // Set new number of active time levels
    for (int n=0; n<CCTK_NumGroups(); ++n) {
      int const grouptype = CCTK_GroupTypeI (n);
      if (grouptype == CCTK_GF) {
        for (int ml=0; ml<mglevels; ++ml) {
          groupdata.AT(n).activetimelevels.AT(ml).resize
            (reflevels, groupdata.AT(n).activetimelevels.AT(ml).AT(0));
        }
      }
    }
    
    // Set new level times
    // for (int ml=0; ml<mglevels; ++ml) {
    //   leveltimes.AT(ml).resize
    //     (reflevels, leveltimes.AT(ml).AT(oldreflevels-1));
    // }
    
    ++ regridding_epoch;
    // Mark all vanished levels as changed
    for (int rl=reflevels; rl<int(level_regridding_epochs.size()); ++rl) {
      ++ level_regridding_epochs.at(rl);
    }
    // Insert entries for new levels
    if (int(level_regridding_epochs.size()) < reflevels) {
      level_regridding_epochs.resize(reflevels, 0);
    }
    
    OutputGridStatistics (cctkGH);
    
    timer.stop();
  }
  
  
  
  bool
  Recompose (cGH const * const cctkGH,
             int const rl,
             bool const do_init)
  {
    char const * const where = "Recompose";
    static Timer timer (where);
    timer.start();
    
    bool did_recompose = false;
    
    for (int m=0; m<maps; ++m) {
      Waypoint ("Recomposing the grid hierarchy for map %d level %d...", m, rl);
      
      assert (rl>=0 and rl<vhh.AT(m)->reflevels());
      did_recompose |= vhh.AT(m)->recompose (rl, do_init);
      
      Waypoint ("Done recomposing the grid hierarchy for map %d level %d.",
                m, rl);
    }
      
    ClassifyPoints (cctkGH, rl);
    
    if (did_recompose) {
      ++ level_regridding_epochs.at(rl);
    }
    
    timer.stop();
    return did_recompose;
  }
  
  
  
  void
  RegridFree (cGH const * const cctkGH,
              bool const do_init)
  {
    char const * const where = "RegridFree";
    static Timer timer (where);
    timer.start();
    
    Checkpoint ("Freeing after regridding level %d...", reflevel);
    
    assert (is_level_mode());
    
    // Free old grid structure
    for (int m=0; m<maps; ++m) {
      vhh.AT(m)->regrid_free (do_init);
    }
    
    timer.stop();
  }
  
  
  
  void
  OutputSuperregions (cGH const * const cctkGH,
                      int const m,
                      gh const & hh,
                      gh::rregs const & superregss)
  {
    CCTK_INFO ("Grid structure (superregions, grid points):");
    for (int rl=0; rl<(int)superregss.size(); ++rl) {
      assert (rl < hh.reflevels());
      for (int c=0; c<(int)superregss.AT(rl).size(); ++c) {
        const ivect & levfact = spacereffacts.AT(rl);
        const ibbox & ext     = superregss.AT(rl).AT(c).extent;
        const ivect & lower   = ext.lower();
        const ivect & upper   = ext.upper();
        const ivect & stride  = ext.stride();
        assert (all(lower * levfact % maxspacereflevelfact == 0));
        assert (all(upper * levfact % maxspacereflevelfact == 0));
        cout << "   [" << rl << "][" << m << "][" << c << "]"
             << "   exterior: "
             << lower / stride
             << " : "
             << upper / stride
             << "   ("
             << (upper - lower) / stride + 1
             << " + PADDING) "
             << prod ((upper - lower) / stride + 1)
             << eol;
      }
    }
    
    CCTK_INFO ("Grid structure (superregions, coordinates):");
    for (int rl=0; rl<(int)superregss.size(); ++rl) {
      assert (rl < hh.reflevels());
      for (int c=0; c<(int)superregss.AT(rl).size(); ++c) {
        const rvect origin    = domainspecs.AT(m).exterior_min;
        const rvect delta     = (domainspecs.AT(m).exterior_max - domainspecs.AT(m).exterior_min) / rvect (domainspecs.AT(m).npoints - 1);
        const ibbox & ext     = superregss.AT(rl).AT(c).extent;
        const ivect & ilower  = ext.lower();
        const ivect & iupper  = ext.upper();
        const ivect & levfact = spacereffacts.AT(rl);
        const ibbox & base    = hh.baseextent(0,0);
        const ivect & bstride = base.stride();
        cout << "   [" << rl << "][" << m << "][" << c << "]"
             << "   exterior: "
             << origin + delta * rvect(ilower) / rvect(bstride)
             << " : "
             << origin + delta * rvect(iupper) / rvect(bstride)
             << " : "
             << delta / rvect(levfact) << eol;
      }
    }
    
    fflush (stdout);
  }
  
  
  
  void
  OutputGrids (cGH const * const cctkGH,
               int const m,
               gh const & hh,
               dh const & dd)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (verbose or veryverbose) {
      
      CCTK_INFO ("Grid structure (grid points):");
      for (int ml=0; ml<hh.mglevels(); ++ml) {
        for (int rl=0; rl<hh.reflevels(); ++rl) {
          for (int c=0; c<hh.components(rl); ++c) {
            const int convfact = ipow(mgfact, ml);
            const ivect levfact = spacereffacts.AT(rl);
            const ibbox   ext = hh.extent(ml,rl,c);
            const ivect & lower  = ext.lower();
            const ivect & upper  = ext.upper();
            const ivect & stride = ext.stride();
            assert (all(lower * levfact % maxspacereflevelfact == 0));
            assert (all(upper * levfact % maxspacereflevelfact == 0));
            assert (all(((upper - lower) * levfact / maxspacereflevelfact)
                        % convfact == 0));
            cout << "   [" << ml << "][" << rl << "][" << m << "][" << c << "]"
                 << "   exterior: "
                 << "proc "
                 << hh.processor(rl,c)
                 << "   "
                 << lower / stride
                 << " : "
                 << upper / stride
                 << "   ("
                 << (upper - lower) / stride + 1
                 << ") "
                 << prod ((upper - lower) / stride + 1)
                 << eol;
          }
        }
      }
      
      CCTK_INFO ("Grid structure (boundaries):");
      for (int rl=0; rl<hh.reflevels(); ++rl) {
        for (int c=0; c<hh.components(rl); ++c) {
          cout << "   [" << rl << "][" << m << "][" << c << "]"
               << "   bbox: "
               << hh.outer_boundaries(rl,c)
               << eol;
        }
      }
      
      CCTK_INFO ("Grid structure (coordinates):");
      for (int ml=0; ml<hh.mglevels(); ++ml) {
        for (int rl=0; rl<hh.reflevels(); ++rl) {
          for (int c=0; c<hh.components(rl); ++c) {
            const rvect origin = domainspecs.AT(m).exterior_min;
            const rvect delta  = (domainspecs.AT(m).exterior_max - domainspecs.AT(m).exterior_min) / rvect (domainspecs.AT(m).npoints - 1);
            const ibbox ext = hh.extent(ml,rl,c);
            const ivect & ilower = ext.lower();
            const ivect & iupper = ext.upper();
            const int convfact = ipow(mgfact, ml);
            const ivect levfact = spacereffacts.AT(rl);
            const ibbox & base    = hh.baseextent(ml,0);
            const ivect & bstride = base.stride();
            cout << "   [" << ml << "][" << rl << "][" << m << "][" << c << "]"
                 << "   exterior: "
                 << origin + delta * rvect(ilower) / rvect(bstride)
                 << " : "
                 << origin + delta * rvect(iupper) / rvect(bstride)
                 << " : "
                 << delta * rvect(convfact) / rvect(levfact) << eol;
          }
        }
      }
      
      CCTK_INFO ("Grid structure (coordinates, including ghosts):");
      for (int ml=0; ml<hh.mglevels(); ++ml) {
        for (int rl=0; rl<hh.reflevels(); ++rl) {
          for (int c=0; c<hh.components(rl); ++c) {
            const rvect origin = domainspecs.AT(m).exterior_min;
            const rvect delta  = (domainspecs.AT(m).exterior_max - domainspecs.AT(m).exterior_min) / rvect (domainspecs.AT(m).npoints - 1);
            const ivect lower = dd.light_boxes.AT(ml).AT(rl).AT(c).exterior.lower();
            const ivect upper = dd.light_boxes.AT(ml).AT(rl).AT(c).exterior.upper();
            const int convfact = ipow(mgfact, ml);
            const ivect levfact = spacereffacts.AT(rl);
            const ibbox & base    = hh.baseextent(ml,0);
            const ivect & bstride = base.stride();
            cout << "   [" << ml << "][" << rl << "][" << m << "][" << c << "]"
                 << "   exterior: "
                 << origin + delta * rvect(lower) / rvect(bstride)
                 << " : "
                 << origin + delta * rvect(upper) / rvect(bstride)
                 << " : "
                 << delta * rvect(convfact) / rvect(levfact) << eol;
          }
        }
      }
      
      CCTK_INFO ("Grid statistics:");
      const streamsize oldprecision = cout.precision();
      const ios_base::fmtflags oldflags = cout.flags();
      cout.setf (ios::fixed);
      for (int ml=0; ml<hh.mglevels(); ++ml) {
        CCTK_REAL coarsevolume = 0;
        for (int rl=0; rl<hh.reflevels(); ++rl) {
          
          const CCTK_REAL basevolume = hh.baseextents.AT(0).AT(0).size();
          CCTK_REAL countvolume = 0;
          CCTK_REAL totalvolume = 0;
          CCTK_REAL totalvolume2 = 0;
          
          for (int c=0; c<hh.components(rl); ++c) {
            const CCTK_REAL volume = hh.extent(ml,rl,c).size();
            ++ countvolume;
            totalvolume += volume;
            totalvolume2 += ipow(volume, 2);
          }
          
          const CCTK_REAL avgvolume = totalvolume / countvolume;
          const CCTK_REAL stddevvolume
            = sqrt (max (CCTK_REAL(0),
                         totalvolume2 / countvolume - ipow (avgvolume, 2)));
          
          for (int c=0; c<hh.components(rl); ++c) {
            const CCTK_REAL volume = hh.extent(ml,rl,c).size();
            cout << "   [" << ml << "][" << rl << "][" << m << "][" << c << "]"
                 << "   volume: " << setprecision(0) << volume
                 << "   of parent: " << setprecision(1) << 100 * volume / totalvolume << "%"
                 << "   of domain: " << setprecision(1) << 100 * volume / basevolume << "%"
                 << eol;
          }
          
          cout << "   [" << ml << "][" << rl << "][" << m << "]"
               << "   average volume: " << setprecision(0) << avgvolume
               << "   of parent: " << setprecision(1) << 100 * avgvolume / totalvolume << "%"
               << "   of domain: " << setprecision(1) << 100 * avgvolume / basevolume << "%"
               << eol;
          cout << "   [" << ml << "][" << rl << "][" << m << "]"
               << "   standard deviation: " << setprecision(0) << stddevvolume
               << "   of parent: " << setprecision(1) << 100 * stddevvolume / totalvolume << "%"
               << "   of domain: " << setprecision(1) << 100 * stddevvolume / basevolume << "%"
               << eol;
          
          // TODO: ghost points vs. interior points (and boundary
          // points)
          
          CCTK_REAL countquadrupole = 0;
          CCTK_REAL minquadrupole = 1;
          CCTK_REAL totalquadrupole = 0;
          CCTK_REAL totalquadrupole2 = 0;
          
          for (int c=0; c<hh.components(rl); ++c) {
            const ibbox ext = hh.extent(ml,rl,c);
            const ivect shape = ext.shape();
            const ivect stride = ext.stride();
            const CCTK_REAL minlength = minval (rvect (shape) / rvect (stride));
            const CCTK_REAL maxlength = maxval (rvect (shape) / rvect (stride));
            const CCTK_REAL quadrupole = minlength / maxlength;
            ++ countquadrupole;
            minquadrupole = min (minquadrupole, quadrupole);
            totalquadrupole += quadrupole;
            totalquadrupole2 += ipow (quadrupole, 2);
            cout << "   [" << ml << "][" << rl << "][" << m << "][" << c << "]"
                 << "   length ratio: " << setprecision(3) << quadrupole
                 << eol;
          }
          
          const CCTK_REAL avgquadrupole = totalquadrupole / countquadrupole;
          const CCTK_REAL stddevquadrupole
            = sqrt (max (CCTK_REAL(0),
                         (totalquadrupole2 / countquadrupole
                          - ipow (avgquadrupole, 2))));
          
          cout << "   [" << ml << "][" << rl << "][" << m << "]"
               << "   average length ratio: " << setprecision(3) << avgquadrupole
               << "   standard deviation: " << setprecision(3) << stddevquadrupole
               << eol;
          
          // TODO: processor distribution, average load, std deviation
          
          coarsevolume = totalvolume * prod (rvect (spacereflevelfact));
        }
      }
      cout.precision (oldprecision);
      cout.setf (oldflags);
      
      fflush (stdout);
    }
  }
  
  
  
  void
  OutputGridStructure (cGH const * const cctkGH,
                       int const m,
                       gh::mregs const & regsss)
  {
    DECLARE_CCTK_PARAMETERS;
    
    // Output only on the root processor
    if (CCTK_MyProc(cctkGH) != 0) return;
    
    // Output only if output is desired
    if (CCTK_EQUALS (grid_structure_filename, "")) return;
    
    // Create the output directory
    CCTK_CreateDirectory (0755, out_dir);
    
    ostringstream filenamebuf;
    filenamebuf << out_dir << "/" << grid_structure_filename;
    // we need a persistent temporary here
    string filenamestr = filenamebuf.str();
    const char * filename = filenamestr.c_str();
    
    ofstream file;
    
    static bool do_truncate = true;
    
    if (do_truncate) {
      do_truncate = false;
      struct stat fileinfo;
      if (IO_TruncateOutputFiles (cctkGH)
	  or stat(filename, &fileinfo)!=0) {
	file.open (filename, ios::out | ios::trunc);
	assert (file.good());
	file << "# grid structure" << eol
	     << "# format: map reflevel component mglevel   processor bounding-box is-outer-boundary" << eol;
	assert (file.good());
      }
    }
    if (not file.is_open()) {
      file.open (filename, ios::app);
      assert (file.good());
    }
    
    if (m == 0) {
      // Output a separator
      file << eol;
      file << "iteration " << cctkGH->cctk_iteration << eol;
      file << "maps " << maps << eol;
    }
    file << m << " mglevels " << regsss.size() << eol;
    for (int ml=0; ml<(int)regsss.size(); ++ml) {
      file << m << " " << ml << " reflevels " << regsss.AT(ml).size() << eol;
      for (int rl=0; rl<(int)regsss.AT(ml).size(); ++rl) {
        file << m << " " << ml << " " << rl << " components " << regsss.AT(ml).AT(rl).size() << eol;
        for (int c=0; c<(int)regsss.AT(ml).AT(rl).size(); ++c) {
          file << m << " " << ml << " " << rl << " " << c << "   "
               << regsss.AT(ml).AT(rl).AT(c).processor << " "
               << regsss.AT(ml).AT(rl).AT(c).extent << " "
               << regsss.AT(ml).AT(rl).AT(c).outer_boundaries << eol;
        }
      }
    }
    
    file.close();
    assert (file.good());
  }
  
  
  
  void
  OutputGridCoordinates (cGH const * const cctkGH,
                         int const m,
                         gh::mregs const & regsss)
  {
    DECLARE_CCTK_PARAMETERS;
    
    // Output only on the root processor
    if (CCTK_MyProc(cctkGH) != 0) return;
    
    // Output only if output is desired
    if (CCTK_EQUALS (grid_coordinates_filename, "")) return;
    
    // Create the output directory
    CCTK_CreateDirectory (0755, out_dir);
    
    ostringstream filenamebuf;
    filenamebuf << out_dir << "/" << grid_coordinates_filename;
    // we need a persistent temporary here
    string filenamestr = filenamebuf.str();
    const char * filename = filenamestr.c_str();
    
    ofstream file;
    
    static bool do_truncate = true;
    
    if (do_truncate) {
      do_truncate = false;
      struct stat fileinfo;
      if (IO_TruncateOutputFiles (cctkGH)
	  or stat(filename, &fileinfo)!=0) {
	file.open (filename, ios::out | ios::trunc);
	assert (file.good());
	file << "# grid coordinates" << eol
	     << "# format: map reflevel region mglevel   bounding-box" << eol;
	assert (file.good());
      }
    }
    if (not file.is_open()) {
      file.open (filename, ios::app);
      assert (file.good());
    }
    
    // Affine transformation between index space and coordinate space
    rvect const origin  = domainspecs.AT(m).exterior_min;
    rvect const spacing =
      (domainspecs.AT(m).exterior_max - domainspecs.AT(m).exterior_min) / 
        rvect (domainspecs.AT(m).npoints - 1);
    rvect const scale =
      rvect (vhh.AT(m)->baseextents.AT(0).AT(0).stride()) / spacing;
    
    
    
    file << "iteration " << cctkGH->cctk_iteration << eol;
    file << "maps " << maps << eol;
    file << m << " mglevels " << regsss.size() << eol;
    for (int ml=0; ml<(int)regsss.size(); ++ml) {
      file << m << " " << ml << " reflevels " << regsss.AT(ml).size() << eol;
      for (int rl=0; rl<(int)regsss.AT(ml).size(); ++rl) {
        ibset const extents (regsss.AT(ml).AT(rl), & region_t::extent);
        file << m << " " << ml << " " << rl << " regions " << extents.setsize() << eol;
        int c=0;
        for (ibset::const_iterator
               bi = extents.begin(); bi != extents.end(); ++ bi)
        {
#if 0
          ibbox const & ext       = * bi;
          ibbox const & baseext   = vhh.AT(m)->baseextents.AT(ml).AT(rl);
          ibbox const & coarseext = vhh.AT(m)->baseextents.AT(ml).AT(0 );
          
          // This is nice, but is wrong since CartGrid3D has not yet
          // initialised the coordinates
          ivect const cctk_levfac       = spacereffacts.AT(rl);
          ivect const cctk_levoff       = baseext.lower() - coarseext.lower();
          ivect const cctk_levoffdenom  = baseext.stride();
          
          ivect const cctk_lbnd =
            (ext.lower() - baseext.lower()) / ext.stride();
          ivect const cctk_ubnd =
            (ext.upper() - baseext.lower()) / ext.stride();
          
          rvect const cctk_origin_space =
            origin_space.AT(m).AT(ml);
          rvect const cctk_delta_space  =
            delta_space.AT(m) * rvect (ipow (mgfact, ml));
          
          rvect const CCTK_ORIGIN_SPACE =
            cctk_origin_space +
            cctk_delta_space / rvect (cctk_levfac) *
            rvect (cctk_levoff) / rvect (cctk_levoffdenom);
          rvect const CCTK_DELTA_SPACE =
            cctk_delta_space / rvect (cctk_levfac);
          
          rvect const rlower =
            CCTK_ORIGIN_SPACE + rvect (cctk_lbnd) * CCTK_DELTA_SPACE;
          rvect const rupper =
            CCTK_ORIGIN_SPACE + rvect (cctk_ubnd) * CCTK_DELTA_SPACE;
          rvect const rdelta =
            CCTK_DELTA_SPACE;
#endif
          
          ibbox const & ext = * bi;
          
          rvect const rlower = origin + rvect (ext.lower()) / scale;
          rvect const rupper = origin + rvect (ext.upper()) / scale;
          rvect const rdelta = rvect (ext.stride()) / scale;
          
          rbbox const rb (rlower, rupper, rdelta);
          file << m << " " << ml << " " << rl << " " << c << "   " << rb << eol;
          ++c;
        }
      }
    }
    file << eol;
    
    file.close();
    assert (file.good());
  }
  
  
  
  // TODO: this routine should go into CarpetRegrid (except maybe
  // SplitRegions_AlongZ for grid arrays)
  void
  SplitRegions (cGH const * const cctkGH,
                vector<region_t> & superregs,
                vector<region_t> & regs)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (regs.empty());
    if (CCTK_EQUALS (processor_topology, "along-z")) {
      SplitRegions_AlongZ (cctkGH, superregs, regs);
    } else if (CCTK_EQUALS (processor_topology, "along-dir")) {
      SplitRegions_AlongDir (cctkGH, superregs, regs, split_direction);
    } else if (CCTK_EQUALS (processor_topology, "automatic")) {
      SplitRegions_Automatic (cctkGH, superregs, regs);
    } else if (CCTK_EQUALS (processor_topology, "manual")) {
      SplitRegions_AsSpecified (cctkGH, superregs, regs);
    } else {
      assert (0);
    }
  }
  
  
  
  void
  OutputGridStatistics (cGH const * const cctkGH)
  {
    // Grid array statistics
    int num_arrays = 0;
    int num_gfs = 0;
    int size_gfs = 0;
    CCTK_REAL num_active_array_points = 0;
    CCTK_REAL num_total_array_points  = 0;
    CCTK_REAL size_total_array_points = 0;
    for (int g=0; g<CCTK_NumGroups(); ++g) {
      cGroup gdata;
      check (not CCTK_GroupData (g, &gdata));
      int const num_tl = CCTK_ActiveTimeLevelsGI (cctkGH, g);
      int const num_vars = gdata.numvars;
      int const size_vars = gdata.numvars * CCTK_VarTypeSize (gdata.vartype);
      switch (gdata.grouptype) {
      case CCTK_SCALAR:
      case CCTK_ARRAY: {
        num_arrays += num_tl * num_vars;
        gh const * const hh = arrdata.AT(g).AT(0).hh;
        dh const * const dd = arrdata.AT(g).AT(0).dd;
        for (int c=0; c<hh->components(0); ++c) {
          dh::light_dboxes const & b = dd->light_boxes.AT(0).AT(0).AT(c);
          num_active_array_points += num_tl * num_vars  * b.active_size;
          num_total_array_points  += num_tl * num_vars  * b.exterior_size;
          size_total_array_points += num_tl * size_vars * b.exterior_size;
        }
        break;
      }
      case CCTK_GF:
        num_gfs  += num_tl * num_vars;
        size_gfs += num_tl * size_vars;
        break;
      default:
        assert (0);
      }
    }
    
    // Grid function statistics
    int num_comps = 0;
    int num_steps = 0;
    CCTK_REAL num_active_mem_points = 0;
    CCTK_REAL num_owned_mem_points  = 0;
    CCTK_REAL num_total_mem_points  = 0;
    CCTK_REAL size_total_mem_points = 0;
    CCTK_REAL num_active_cpu_points = 0;
    CCTK_REAL num_owned_cpu_points  = 0;
    CCTK_REAL num_total_cpu_points  = 0;
    for (int m=0; m<maps; ++m) {
      gh const * const hh = vhh.AT(m);
      dh const * const dd = vdd.AT(m);
      for (int ml=0; ml<mglevels; ++ml) {
        for (int rl=0; rl<reflevels; ++rl) {
          int const trf = timereffacts.AT(rl);
          if (m==0 and ml==0) {
            num_steps += trf;
          }
          for (int c=0; c<hh->components(rl); ++c) {
            ++ num_comps;
            dh::light_dboxes const & b = dd->light_boxes.AT(ml).AT(rl).AT(c);
            num_active_mem_points += num_gfs  * b.active_size;
            num_owned_mem_points  += num_gfs  * b.owned_size;
            num_total_mem_points  += num_gfs  * b.exterior_size;
            size_total_mem_points += size_gfs * b.exterior_size;
            num_active_cpu_points += trf * b.active_size;
            num_owned_cpu_points  += trf * b.owned_size;
            num_total_cpu_points  += trf * b.exterior_size;
          }
        }
      }
    }
    num_active_cpu_points /= num_steps * delta_time;
    num_owned_cpu_points  /= num_steps * delta_time;
    num_total_cpu_points  /= num_steps * delta_time;
    
    // Load balance statistics
    CCTK_REAL const rzero = 0;
    CCTK_REAL const rmax  = numeric_limits<CCTK_REAL>::max();
    vector<CCTK_REAL> min_active (reflevels, rmax);
    vector<CCTK_REAL> max_active (reflevels, 0);
    vector<CCTK_REAL> avg_active (reflevels, 0);
    vector<CCTK_REAL> sdv_active (reflevels, 0);
    for (int rl=0; rl<reflevels; ++rl) {
      vector<CCTK_REAL> num_active_per_proc (dist::size(), 0);
      for (int m=0; m<maps; ++m) {
        gh const * const hh = vhh.AT(m);
        dh const * const dd = vdd.AT(m);
        for (int ml=0; ml<mglevels; ++ml) {
          for (int c=0; c<hh->components(rl); ++c) {
            int const p = hh->processor(rl,c);
            dh::light_dboxes const & b = dd->light_boxes.AT(ml).AT(rl).AT(c);
            num_active_per_proc.AT(p) += num_gfs * b.active_size;
          }
        }
      }
      for (int p=0; p<dist::size(); ++p) {
        min_active.AT(rl) = min (min_active.AT(rl), num_active_per_proc.AT(p));
        max_active.AT(rl) = max (max_active.AT(rl), num_active_per_proc.AT(p));
        avg_active.AT(rl) += num_active_per_proc.AT(p);
        sdv_active.AT(rl) += ipow (num_active_per_proc.AT(p), 2);
      }
      avg_active.AT(rl) /= dist::size();
      sdv_active.AT(rl) /= dist::size();
      sdv_active.AT(rl) =
        sqrt (max (rzero, sdv_active.AT(rl) - ipow (avg_active.AT(rl), 2)));
    } // for rl
    
    // Output
    CCTK_VInfo (CCTK_THORNSTRING,
                "Global grid structure statistics:");
    CCTK_VInfo (CCTK_THORNSTRING,
                "GF: rhs: %.0fk active, %.0fk owned (+%.0f%%), %.0fk total (+%.0f%%), %.3g steps/time",
                double (num_active_cpu_points / 1e+3),
                double (num_owned_cpu_points / 1e+3),
                double (num_active_cpu_points == 0
                        ? 0
                        : num_owned_cpu_points / num_active_cpu_points * 100 - 100),
                double (num_total_cpu_points / 1e+3),
                double (num_owned_cpu_points == 0
                        ? 0
                        : num_total_cpu_points / num_owned_cpu_points * 100 - 100),
                double (num_steps / delta_time));
    CCTK_VInfo (CCTK_THORNSTRING,
                "GF: vars: %d, pts: %.0fM active, %.0fM owned (+%.0f%%), %.0fM total (+%.0f%%), %.1f comp/proc",
                num_gfs,
                double (num_active_mem_points / 1e+6),
                double (num_owned_mem_points / 1e+6),
                double (num_active_mem_points == 0
                        ? 0
                        : num_owned_mem_points / num_active_mem_points * 100 - 100),
                double (num_total_mem_points / 1e+6),
                double (num_owned_mem_points == 0
                        ? 0
                        : num_total_mem_points / num_owned_mem_points * 100 - 100),
                double (num_comps / (reflevels * CCTK_nProcs (cctkGH))));
    CCTK_VInfo (CCTK_THORNSTRING,
                "GA: vars: %d, pts: %.0fM active, %.0fM total (+%.0f%%)",
                num_arrays,
                double (num_active_array_points / 1e+6),
                double (num_total_array_points / 1e+6),
                double (num_active_array_points == 0
                        ? 0
                        : num_total_array_points / num_active_array_points * 100 - 100));
    CCTK_VInfo (CCTK_THORNSTRING,
                "Total required memory: %.3f GByte (for GAs and currently active GFs)",
                double ((size_total_array_points + size_total_mem_points) / 1e+9));
    CCTK_VInfo (CCTK_THORNSTRING,
                "Load balance:  min     avg     max     sdv      max/avg-1");
    for (int rl=0; rl<reflevels; ++rl) {
      CCTK_VInfo (CCTK_THORNSTRING,
                  "Level %2d:    %4.0fM   %4.0fM   %4.0fM   %4.0fM active   %4.0f%%",
                  rl,
                  double (min_active.AT(rl) / 1e+6),
                  double (avg_active.AT(rl) / 1e+6),
                  double (max_active.AT(rl) / 1e+6),
                  double (sdv_active.AT(rl) / 1e+6),
                  double (100 * (max_active.AT(rl) / avg_active.AT(rl) - 1)));
    }
    
    // After this, we will begin to allocate memory for the grid
    // structure.  If we run out of memory, ensure that this output
    // still makes it to disk.
    fflush (stdout);
    fflush (stderr);
  }
  
  
  
  void
  SplitRegions_AlongZ (cGH const * const cctkGH,
                       vector<region_t> & superregs,
                       vector<region_t> & regs)
  {
    assert (regs.empty());
    SplitRegions_AlongDir (cctkGH, superregs, regs, dim - 1);
  }
  
  
  
  void
  SplitRegions_AlongDir (cGH const * const cctkGH,
                         vector<region_t> & superregs,
                         vector<region_t> & regs,
                         int const dir)
  {
    assert (regs.empty());
    
    // Something to do?
    if (superregs.size() == 0) {
      return;
    }
    
    const int nprocs = CCTK_nProcs (cctkGH);
    
    assert (superregs.size() == 1);
    regs = superregs;
    
    if (nprocs == 1) {
      regs.AT(0).processor = 0;
      pseudoregion_t const preg (regs.AT(0).extent, regs.AT(0).processor);
      assert (superregs.AT(0).processors == NULL);
      superregs.AT(0).processors = new ipfulltree (preg);
      return;
    }
    
    assert (dir>=0 and dir<dim);
    
    region_t const & reg0 = regs.AT(0);
    const ivect  rstr0  = reg0.extent.stride();
    const ivect  rlb0   = reg0.extent.lower();
    const ivect  rub0   = reg0.extent.upper() + rstr0;
    const b2vect obnd0 = reg0.outer_boundaries;
    
    regs.resize (nprocs);
    vector<int> bounds (nprocs+1);
    vector<ipfulltree *> subtrees (nprocs);
    for (int c=0; c<nprocs; ++c) {
      ivect cstr = rstr0;
      ivect clb  = rlb0;
      ivect cub  = rub0;
      const int glonpz = (rub0[dir] - rlb0[dir]) / cstr[dir];
      const int locnpz = (glonpz + nprocs - 1) / nprocs;
      const int zstep  = locnpz * cstr[dir];
      clb[dir] = rlb0[dir] + zstep *  c;
      cub[dir] = rlb0[dir] + zstep * (c+1);
      if (clb[dir] > rub0[dir]) clb[dir] = rub0[dir];
      if (cub[dir] > rub0[dir]) cub[dir] = rub0[dir];
      assert (clb[dir] <= cub[dir]);
      assert (cub[dir] <= rub0[dir]);
      region_t & reg = regs.AT(c);
      ibbox    & ext  = reg.extent;
      b2vect   & obnd = reg.outer_boundaries;
      int      & proc = reg.processor;
      ext = ibbox(clb, cub-cstr, cstr);
      obnd = obnd0;
      obnd[0][dir] &= clb[dir] == rlb0[dir];
      obnd[1][dir] &= cub[dir] == rub0[dir];
      proc = c;
      pseudoregion_t const preg (reg.extent, c);
      bounds.AT(c) = reg.extent.lower()[dir];
      subtrees.AT(c) = new ipfulltree (preg);
    }
    bounds.AT(nprocs) = rub0[dir] + rstr0[dir];
    
    assert (superregs.AT(0).processors == NULL);
    superregs.AT(0).processors = new ipfulltree (dir, bounds, subtrees);
    
    for (int c=0; c<(int)regs.size(); ++c) {
      assert (regs.AT(c).processor == c);
    }
  }
  
  
  
  void
  SplitRegions_Automatic (cGH const * const cctkGH,
                          vector<region_t> & superregs,
                          vector<region_t> & regs)
  {
    assert (regs.empty());
    vector<vector<region_t> > superregss (1, superregs);
    vector<vector<region_t> > regss (1);
    SplitRegionsMaps_Automatic (cctkGH, superregss, regss);
    assert (superregss.size() == 1);
    superregs = superregss.AT(0);
    assert (regss.size() == 1);
    regs = regss.AT(0);
  }
  
  
  
  static void
  SplitRegions_AsSpecified (cGH const * const cctkGH,
                            vector<region_t> & superregs,
                            vector<region_t> & regs)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (regs.empty());
    
    // Something to do?
    if (superregs.size() == 0) {
      return;
    }
    
    const int nprocs = CCTK_nProcs (cctkGH);
    
    assert (superregs.size() == 1);
    
    region_t const & reg0 = superregs.AT(0);
    const ivect rstr0  = reg0.extent.stride();
    const ivect rlb0   = reg0.extent.lower();
    const ivect rub0   = reg0.extent.upper() + rstr0;
    const b2vect obnd0 = reg0.outer_boundaries;
    
    const ivect nprocs_dir
      (processor_topology_3d_x,
       processor_topology_3d_y, 
       processor_topology_3d_z);
    assert (all (nprocs_dir > 0));
    if (prod (nprocs_dir) != nprocs) {
      CCTK_VWarn (0, __LINE__, __FILE__, CCTK_THORNSTRING,
		  "The specified processor topology [%d,%d,%d] requires %d processors, but there are %d processors",
                  nprocs_dir[0], nprocs_dir[1], nprocs_dir[2],
                  prod (nprocs_dir),
                  nprocs);
    }
    assert (prod (nprocs_dir) == nprocs);
    
    regs.resize (nprocs);
    const ivect cstr = rstr0;
    const ivect glonp = (rub0 - rlb0) / cstr;
//     const ivect locnp = (glonp + nprocs_dir - 1) / nprocs_dir;
    const ivect locnp = glonp / nprocs_dir;
    const ivect rem = glonp % nprocs_dir;
    const ivect step = locnp * cstr;
    assert (dim==3);
    
    vector<int> boundsz(nprocs_dir[2]);
    vector<ipfulltree *> subtreesz(nprocs_dir[2]+1);
    for (int k=0; k<nprocs_dir[2]; ++k) {
      vector<int> boundsy(nprocs_dir[2]);
      vector<ipfulltree *> subtreesy(nprocs_dir[2]+1);
      for (int j=0; j<nprocs_dir[1]; ++j) {
        vector<int> boundsx(nprocs_dir[2]);
        vector<ipfulltree *> subtreesx(nprocs_dir[2]+1);
	for (int i=0; i<nprocs_dir[0]; ++i) {
          
	  const int c = i + nprocs_dir[0] * (j + nprocs_dir[1] * k);
	  const ivect ipos (i, j, k);
	  ivect clb = rlb0 + step *  ipos;
	  ivect cub = rlb0 + step * (ipos+1);
// 	  clb = min (clb, rub);
// 	  cub = min (cub, rub);
          for (int d=0; d<dim; ++d) {
            if (ipos[d]<rem[d]) {
              clb[d] += cstr[d] * ipos[d];
              cub[d] += cstr[d] * (ipos[d]+1);
            } else {
              clb[d] += cstr[d] * rem[d];
              cub[d] += cstr[d] * rem[d];
            }
          }
	  assert (all (clb >= 0));
	  assert (all (clb <= cub));
	  assert (all (cub <= rub0));
          assert (all (not (ipos==0) or clb==rlb0));
          assert (all (not (ipos==nprocs_dir-1) or cub==rub0));
          region_t & reg  = regs.AT(c);
          ibbox    & ext  = reg.extent;
          b2vect   & obnd = reg.outer_boundaries;
          int      & proc = reg.processor;
          ext = ibbox(clb, cub-cstr, cstr);
	  obnd = obnd0;
          obnd[0] &= clb == rlb0;
          obnd[1] &= cub == rub0;
          proc = c;
          
          pseudoregion_t preg (reg.extent, c);
          subtreesx.AT(i) = new ipfulltree (preg);
	}
        boundsx.AT(nprocs_dir[0]) = rub0[0] + rstr0[0];
        subtreesy.AT(j) = new ipfulltree (0, boundsx, subtreesx);
      }
      boundsy.AT(nprocs_dir[1]) = rub0[1] + rstr0[1];
      subtreesz.AT(k) = new ipfulltree (1, boundsy, subtreesy);
    }
    boundsz.AT(nprocs_dir[2]) = rub0[2] + rstr0[2];
    
    assert (superregs.AT(0).processors == NULL);
    superregs.AT(0).processors = new ipfulltree (2, boundsz, subtreesz);
    
    for (int c=0; c<(int)regs.size(); ++c) {
      assert (regs.AT(c).processor == c);
    }
  }
  
  
  
  // TODO: this routine should go into CarpetRegrid (except maybe
  // SplitRegions_AlongZ for grid arrays)
  void
  SplitRegionsMaps (cGH const * const cctkGH,
                    vector<vector<region_t> > & superregss,
                    vector<vector<region_t> > & regss)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (regss.size() == superregss.size());
    for (size_t m=0; m<regss.size(); ++m) {
      assert (regss.AT(m).empty());
    }
    if (CCTK_EQUALS (processor_topology, "along-z")) {
      assert (0);
//       SplitRegionsMaps_AlongZ (cctkGH, superregss, regss);
    } else if (CCTK_EQUALS (processor_topology, "along-dir")) {
      assert (0);
//       SplitRegionsMaps_AlongDir (cctkGH, superregss, regss, split_direction);
    } else if (CCTK_EQUALS (processor_topology, "automatic")) {
      SplitRegionsMaps_Automatic (cctkGH, superregss, regss);
    } else if (CCTK_EQUALS (processor_topology, "recursive")) {
      SplitRegionsMaps_Recursively (cctkGH, superregss, regss);
    } else if (CCTK_EQUALS (processor_topology, "balanced")) {
      SplitRegionsMaps_Balanced (cctkGH, superregss, regss);
    } else if (CCTK_EQUALS (processor_topology, "manual")) {
      assert (0);
//       SplitRegionsMaps_AsSpecified (cctkGH, superregss, regss);
    } else {
      assert (0);
    }
    
    for (size_t m=0; m<superregss.size(); ++m) {
      for (size_t r=0; r<superregss.AT(m).size(); ++r) {
        bool const good_superregs = superregss.AT(m).AT(r).check_region(true);
        if (not good_superregs) {
          cout << "superregs[" << m << "][" << r << "]:\n"
               << superregss.AT(m).AT(r) << "\n";
          cout << "all superregions:\n" << superregss << "\n";
          CCTK_VWarn(CCTK_WARN_ALERT, __LINE__, __FILE__, CCTK_THORNSTRING,
                     "Superregions failed self-check for map %d superregion %d",
                     int(m), int(r));
        }
      }
    }
    assert(regss.size() == superregss.size());
    for (size_t m=0; m<regss.size(); ++m) {
      for (size_t r=0; r<regss.AT(m).size(); ++r) {
        bool const good_regs = regss.AT(m).AT(r).check_region(false);
        if (not good_regs) {
          cout << "regs[" << m << "][" << r << "]:\n"
               << regss.AT(m).AT(r) << "\n";
          cout << "all superregions:\n" << superregss << "\n";
          cout << "all regions:\n" << regss << "\n";
          CCTK_VWarn(CCTK_WARN_ALERT, __LINE__, __FILE__, CCTK_THORNSTRING,
                     "Regions failed self-check for map %d superregion %d",
                     int(m), int(r));
        }
      }
    }
  }
  
  
  
  static void
  SplitRegionsMaps_Automatic_Recursively (bvect    const   & dims,
                                          int      const     firstproc,
                                          int      const     nprocs,
                                          region_t         & superreg,
                                          vector<region_t> & newregs)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (recompose_verbose) cout << "SRMAR enter" << endl;
    // Check preconditions
    assert (nprocs >= 1);
    
    // Remember old size
    size_t const oldsize = newregs.size();
    
    // Are we done?
    if (all(dims)) {
      if (recompose_verbose) cout << "SRMAR bottom" << endl;
      
      // Check precondition
      assert (nprocs == 1);
      
      // Return argument
      region_t newreg = superreg;
      newreg.processor = firstproc;
      newregs.push_back (newreg);
      pseudoregion_t const preg (newreg.extent, newreg.processor);
      superreg.processors = new ipfulltree (preg);
      
      // Check postcondition
      assert (newregs.size() == oldsize + nprocs);
      
      if (recompose_verbose) cout << "SRMAR exit" << endl;
      return;
    }
    
    // Special case
    if (superreg.extent.empty()) {
      if (recompose_verbose) cout << "SRMAR empty" << endl;
      
      // Choose a direction
      int mydim = -1;
      for (int d=0; d<dim; ++d) {
        if (not dims[d]) {
          mydim = d;
          break;
        }
      }
      assert (mydim>=0 and mydim<dim);
      
      // Create a new region
      region_t newreg (superreg);
      newreg.outer_boundaries = b2vect(false);
      if (recompose_verbose) cout << "SRMAR newreg " << newreg << endl;
      
#if 0
      // Store
      for (int p=0; p<nprocs; ++p) {
        newreg.processor = firstproc + p;
        newregs.push_back (newreg);
      }
      assert (superreg.processors == NULL);
      superreg.processors = new ipfulltree ();
#endif
      
      // Store
      vector<int> bounds (nprocs+1);
      vector<ipfulltree *> subtrees (nprocs);
      for (int p=0; p<nprocs; ++p) {
        newreg.processor = firstproc + p;
        newregs.push_back (newreg);
        bounds.AT(p) = newreg.extent.lower()[mydim];
        pseudoregion_t const preg (newreg.extent, newreg.processor);
        subtrees.AT(p) = new ipfulltree (preg);
      }
      bounds.AT(nprocs) = newreg.extent.lower()[mydim];
      assert (superreg.processors == NULL);
      superreg.processors = new ipfulltree (mydim, bounds, subtrees);
      
      // Check postcondition
      assert (newregs.size() == oldsize + nprocs);
      
      if (recompose_verbose) cout << "SRMAR exit" << endl;
      return;
    }
    
    // Calculate cost factors
    rvect rcost = cost (superreg);
    CCTK_REAL const rfact = pow (nprocs / prod(rcost), CCTK_REAL(1)/dim);
    rcost *= rfact;
    assert (abs (prod (rcost) - nprocs) <= 1.0e-6);
    if (recompose_verbose) cout << "SRMA shapes " << rcost << endl;
    
    // Choose a direction
    int mydim = -1;
    int nslices = -1;
    if (no_split_direction!=-1 and not dims[no_split_direction]) {
      mydim = no_split_direction;
      nslices = 1;
    } else {
      int alldims = 0;
      CCTK_REAL mycost = 0;
      CCTK_REAL totalcost = 1;
      // Prefer to split in the z direction
      for (int d=dim-1; d>=0; --d) {
        if (not dims[d]) {
          ++ alldims;
          CCTK_REAL const thiscost = rcost[d];
          if (thiscost >= 0.999999 * mycost) {
            mydim = d;
            mycost = thiscost;
          }
          totalcost *= thiscost;
        }
      }
      assert (mydim>=0 and mydim<dim);
      assert (mycost>=0);
      if (recompose_verbose) cout << "SRMAR mydim " << mydim << endl;
      if (recompose_verbose) cout << "SRMAR mycost " << mycost << endl;
      
      // Choose a number of slices for this direction
      CCTK_REAL const mycost1 =
        mycost * pow(nprocs / totalcost, CCTK_REAL(1) / alldims);
      nslices = min (nprocs, int (floor (mycost1 + CCTK_REAL(0.5))));
    }
    assert (nslices <= nprocs);
    if (recompose_verbose) cout << "SRMAR " << mydim << " nprocs " << nprocs << endl;
    if (recompose_verbose) cout << "SRMAR " << mydim << " nslices " << nslices << endl;
    
    // Mark this direction as done
    assert (not dims[mydim]);
    bvect const newdims = dims.replace(mydim, true);
    
    // Split the remaining processors
    vector<int> mynprocs(nslices);
    int const mynprocs_base = nprocs / nslices;
    int const mynprocs_left = nprocs - nslices * mynprocs_base;
    int sum_mynprocs = 0;
    for (int n=0; n<nslices; ++n) {
      mynprocs.AT(n) = mynprocs_base + int (n < mynprocs_left);
      sum_mynprocs += mynprocs.AT(n);
    }
    assert (sum_mynprocs == nprocs);
    if (recompose_verbose) cout << "SRMAR " << mydim << " mynprocs " << mynprocs << endl;
    
    // Split the region
    vector<int> mynpoints(nslices);
    int const npoints =
      (superreg.extent.shape() / superreg.extent.stride())[mydim];
    int const npoints_bnd_lo =
      superreg.outer_boundaries[0][mydim] ? granularity_boundary : 0;
    int const npoints_bnd_hi =
      superreg.outer_boundaries[1][mydim] ? granularity_boundary : 0;
    
    // Keep track of how many points and processors we have left to
    // distribute
    int npoints_left = npoints;
    int nprocs_left  = nprocs;
    for (int n=0; n<nslices; ++n) {
      assert (nprocs_left > 0);
      CCTK_REAL const npoints1 =
        CCTK_REAL(1) * npoints_left * mynprocs.AT(n) / nprocs_left;
      mynpoints.AT(n) =
        int (floor (npoints1 / granularity + CCTK_REAL(0.5))) * granularity;
      if (n == 0)
        mynpoints.AT(n) += npoints_bnd_lo;
      if (npoints_left - mynpoints.AT(n) <= npoints_bnd_hi)
        mynpoints.AT(n) += npoints_bnd_hi;
      mynpoints.AT(n) = min(npoints_left, mynpoints.AT(n));
      assert (mynpoints.AT(n) >= 0 and mynpoints.AT(n) <= npoints_left);
      npoints_left -= mynpoints.AT(n);
      nprocs_left -= mynprocs.AT(n);
      assert (npoints_left >= 0);
      assert (nprocs_left  >= 0);
    }
    assert (npoints_left == 0);
    assert (nprocs_left  == 0);
    if (recompose_verbose) cout << "SRMAR " << mydim << " mynpoints " << mynpoints << endl;
    
    // Create the regions and recurse
    if (recompose_verbose) cout << "SRMAR " << mydim << ": create bboxes" << endl;
    ivect lo = superreg.extent.lower();
    ivect up = superreg.extent.upper();
    ivect const str = superreg.extent.stride();
    int p = firstproc;
    vector<int> bounds (nslices+1);
    vector<ipfulltree *> subtrees (nslices);
    for (int n=0; n<nslices; ++n) {
      if (recompose_verbose) cout << "SRMAR " << mydim << " n " << n << endl;
      
      // Create a new region
      up[mydim] = lo[mydim] + (mynpoints.AT(n) - 1) * str[mydim];
      b2vect newob (superreg.outer_boundaries);
      if (n > 0) {
        newob[0][mydim] = false;
      }
      if (n < nslices-1) {
        up[mydim] = lo[mydim] + (mynpoints.AT(n) - 1) * str[mydim];
        newob[1][mydim] = false;
      }
      region_t newreg (superreg);
      newreg.extent = ibbox(lo, up, str);
      newreg.outer_boundaries = newob;
      if (recompose_verbose) cout << "SRMAR " << mydim << " newreg " << newreg << endl;
      
      // Recurse
      bounds.AT(n) = lo[mydim];
      SplitRegionsMaps_Automatic_Recursively
        (newdims, p, mynprocs.AT(n), newreg, newregs);
      if (recompose_verbose) cout << "SRMAR " << mydim << " newregs " << newregs << endl;
      
      assert (newreg.processors != NULL);
      subtrees.AT(n) = newreg.processors;
      newreg.processors = NULL;
      newreg.processor = p;
      
      // Next slice
      lo[mydim] = up[mydim] + str[mydim];
      p += mynprocs.AT(n);
    }
    assert (up[mydim] == superreg.extent.upper()[mydim]);
    assert (p == firstproc + nprocs);
    bounds.AT(nslices) = up[mydim] + str[mydim];
    
    // Create tree
    assert (superreg.processors == NULL);
    superreg.processors = new ipfulltree (mydim, bounds, subtrees);
    
    // Check postcondition
    assert (newregs.size() == oldsize + nprocs);
    
    if (recompose_verbose) cout << "SRMAR exit" << endl;
  }
  
  
  
  void
  SplitRegionsMaps_Automatic (cGH const * const cctkGH,
                              vector<vector<region_t> > & superregss,
                              vector<vector<region_t> > & regss)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (recompose_verbose) cout << "SRMA enter" << endl;
    
    int const nmaps = superregss.size();
    int minmap = 1000000000;
    for (int m=0; m<nmaps; ++m) {
      for (int r=0; r<int(superregss.AT(m).size()); ++r) {
        minmap = min (minmap, superregss.AT(m).AT(r).map);
      }
    }
    int nregs = 0;
    for (int m=0; m<nmaps; ++m) {
      nregs += superregss.AT(m).size();
    }
    if (recompose_verbose) cout << "SRMA nregs " << nregs << endl;
    
    // Something to do?
    if (nregs == 0) {
      return;
    }
    
    // Collect slices
    vector<region_t> superregs;
    {
      for (int m=0; m<nmaps; ++m) {
        combine_regions (superregss.AT(m), superregs);
      }
      nregs = superregs.size();
      
      // If the last region was removed, add a new empty region again.
      // A set of regions (corresponding to a refinement level or a
      // grid array) cannot be empty.
      if (nregs == 0) {
        assert (nmaps == 1);    // we should only be here for grid
                                // arrays
        region_t reg;
        reg.extent           = ibbox (ivect (0), ivect (-1), ivect (1));
        reg.outer_boundaries = b2vect (bvect (true), bvect (true));
        reg.map              = 0;
        superregs.push_back (reg);
        nregs = superregs.size();
      }
    }
    
    int const real_nprocs = CCTK_nProcs (cctkGH);
    if (recompose_verbose) cout << "SRMA real_nprocs " << real_nprocs << endl;
    
    // Deactivate some processors if there are too many
    int nprocs;
    if (min_points_per_proc == 0) {
      nprocs = real_nprocs;
    } else {
      CCTK_REAL mycost = 0;
      for (int r=0; r<nregs; ++r) {
        mycost += prod (cost (superregs.AT(r)));
      }
      int const goodnprocs = int (floor (mycost / min_points_per_proc));
      nprocs = max (1, min (real_nprocs, goodnprocs));
    }
    if (recompose_verbose) cout << "SRMA nprocs " << nprocs << endl;
    
    // ncomps: number of components per processor
    int const ncomps = (nregs + nprocs - 1) / nprocs;
    if (recompose_verbose) cout << "SRMA ncomps " << ncomps << endl;
    assert (ncomps > 0);
    
    int const newnregs = nprocs * ncomps;
    if (recompose_verbose) cout << "SRMA newnregs " << newnregs << endl;
    
    if (recompose_verbose) cout << "SRMA: distributing processors to regions" << endl;
    vector<CCTK_REAL> mycosts(nregs);
    for (int r=0; r<nregs; ++r) {
      mycosts.AT(r) = prod (cost (superregs.AT(r)));
    }
    int nregs_left = newnregs;
    vector<int> mynprocs(nregs);
    for (int r=0; r<nregs; ++r) {
      mynprocs.AT(r) = 1;
      -- nregs_left;
    }
    // TODO: split regions if necessary
    while (nregs_left > 0) {
      if (recompose_verbose) cout << "SRMA nregs_left " << nregs_left << endl;
      int maxr = -1;
      CCTK_REAL maxratio = -1;
      for (int r=0; r<nregs; ++r) {
        CCTK_REAL const ratio = mycosts.AT(r) / mynprocs.AT(r);
        if (ratio > maxratio) { maxr=r; maxratio=ratio; }
      }
      assert (maxr>=0 and maxr<nregs);
      ++ mynprocs.AT(maxr);
      if (recompose_verbose) cout << "SRMA maxr " << maxr << endl;
      if (recompose_verbose) cout << "SRMA mynprocs[maxr] " << mynprocs.AT(maxr) << endl;
      -- nregs_left;
    }
    assert (nregs_left == 0);
    int sum_nprocs = 0;
    for (int r=0; r<nregs; ++r) {
      sum_nprocs += mynprocs.AT(r);
    }
    assert (sum_nprocs == newnregs);
    if (recompose_verbose) cout << "SRMA mynprocs " << mynprocs << endl;
    
    if (recompose_verbose) cout << "SRMA: splitting work units" << endl;
    // TODO: rename newregs to regs
    // TODO: rename nregs to nsuperregs
    // TODO: rename newnregs to nregs
    vector<region_t> newregs;
    newregs.reserve (newnregs);
    for (int r=0, p=0; r<nregs; p+=mynprocs.AT(r), ++r) {
      if (recompose_verbose) cout << "SRMA superreg[" << r << "] " << superregs.AT(r) << endl;
      // bvect const dims = false;
      bvect const dims = not split_components;
      SplitRegionsMaps_Automatic_Recursively
        (dims, p, mynprocs.AT(r), superregs.AT(r), newregs);
    } // for r
    if (recompose_verbose) cout << "SRMA newregs " << newregs << endl;
    assert (int(newregs.size()) == newnregs);
    
    // Count components per map
    vector<int> myncomps(nmaps, 0);
#if 0
    vector<int> empty_comps(nmaps, 0);
#endif
    for (int r=0; r<newnregs; ++r) {
      int const m = newregs.AT(r).map - minmap;
      assert (m>=0 and m<nmaps);
#if 0
      if (not newregs.AT(r).extent.empty()) {
        // Ignore empty regions, which may come from empty grid arrays
        ++ myncomps.AT(m);
      } else {
        ++ empty_comps.AT(m);
      }
#endif
      ++ myncomps.AT(m);
    }
    vector<int> mynregs(nmaps, 0);
#if 0
    vector<int> empty_regs(nmaps, 0);
#endif
    for (int r=0; r<nregs; ++r) {
      int const m = superregs.AT(r).map - minmap;
      assert (m>=0 and m<nmaps);
      ++ mynregs.AT(m);
#if 0
      if (not superregs.AT(r).extent.empty()) {
        // Ignore empty regions, which may come from empty grid arrays
        ++ mynregs.AT(m);
      } else {
        ++ empty_regs.AT(m);
      }
#endif
    }
    
    // Convert virtual to real processors
    for (int r=0; r<newnregs; ++r) {
      newregs.AT(r).processor /= ncomps;
      assert (newregs.AT(r).processor >= 0 and
              newregs.AT(r).processor < nprocs);
    }
    {
      vector<int> tmpncomps(nmaps, 0);
      for (int r=0; r<nregs; ++r) {
        int const m = superregs.AT(r).map - minmap;
        ipfulltree * const regf = superregs.AT(r).processors;
        assert (regf != NULL);
        for (ipfulltree::iterator fti (* regf); not fti.done(); ++ fti) {
          pseudoregion_t & preg = (* fti).payload();
          preg.component = tmpncomps.AT(m)++;
        }
      }
      for (int m=0; m<nmaps; ++m) {
        if (not (tmpncomps.AT(m) == myncomps.AT(m))) {
          cout << "Recompose.cc" << endl
               << "superregss=" << superregss << endl
               << "regss=" << regss << endl
               << "nregs=" << nregs << endl
               << "newregs=" << newregs << endl
               << "newnregs=" << newnregs << endl
               << "m=" << m << endl
               << "tmpncomps.AT(m)=" << tmpncomps.AT(m) << endl
               << "myncomps.AT(m)=" << myncomps.AT(m) << endl;
          cout << "newregs:" << endl;
          for (int r=0; r<newnregs; ++r) {
            int const m = newregs.AT(r).map - minmap;
            assert (m>=0 and m<nmaps);
            cout << "r=" << r << endl
                 << "newregs=" << newregs.at(r) << endl
                 << "newregs.extent.size()=" << newregs.at(r).extent.size() << endl
                 << "newregs.extent.empty()=" << newregs.at(r).extent.empty() << endl;
          }
          cout << "*regf:" << endl;
          for (int r=0; r<nregs; ++r) {
            int const m = superregs.AT(r).map - minmap;
            ipfulltree * const regf = superregs.AT(r).processors;
            cout << "r=" << r << endl
                 << "m=" << m << endl;
            cout << "*regf=" << *regf << endl;
          }
          cout.flush();
        }
        assert (tmpncomps.AT(m) == myncomps.AT(m));
      }
    }
    // Count componets per process
    {
      vector<int> ncompsperproc(real_nprocs, 0);
      for (int r=0; r<newnregs; ++r) {
        int const p = newregs.AT(r).processor;
        ++ ncompsperproc.at(p);
      }
      if (same_number_of_components_on_each_process) {
        for (int p=0; p<real_nprocs; ++p) {
          assert (ncompsperproc.at(p) == ncomps);
        }
      }
    }
    
    // Distribute regions
    // Allocate regions
    assert ((int)regss.size() == nmaps);
    for (int m=0; m<nmaps; ++m) {
      assert (regss.AT(m).empty());
      regss.AT(m).reserve (myncomps.AT(m));
      superregss.AT(m).clear();
      superregss.AT(m).reserve (mynregs.AT(m));
    }
    // Assign regions
    for (int r=0; r<newnregs; ++r) {
      int const m = newregs.AT(r).map - minmap;
      assert (m>=0 and m<nmaps);
      regss.AT(m).push_back (newregs.AT(r));
    }
    for (int r=0; r<nregs; ++r) {
      int const m = superregs.AT(r).map - minmap;
      assert (m>=0 and m<nmaps);
      superregss.AT(m).push_back (superregs.AT(r));
    }
    // Output regions
    if (recompose_verbose) {
      cout << "SRMA superregss " << superregss << endl;
      cout << "SRMA regss " << regss << endl;
    }
    // Check sizes
    for (int m=0; m<nmaps; ++m) {
#if 0
      assert (int(regss.AT(m).size()) == myncomps.AT(m) + empty_comps.AT(m));
#endif
      assert (int(regss.AT(m).size()) == myncomps.AT(m));
      assert (int(superregss.AT(m).size()) == mynregs.AT(m));
    }
    
    if (recompose_verbose) cout << "SRMA exit" << endl;
  }
  
  
  
  void
  SplitRegionsMaps_Balanced (cGH const *const cctkGH,
                             vector<vector<region_t> >& superregss,
                             vector<vector<region_t> >& regss)
  {
    DECLARE_CCTK_PARAMETERS;
    
     // Find map numbers
    int const nmaps = superregss.size();
    vector<int> map_index(maps, -1);
    for (int m=0; m<nmaps; ++m) {
      assert(not superregss.AT(m).empty());
      if (not superregss.AT(m).empty()) {
        assert(map_index.AT(superregss.AT(m).AT(0).map) == -1);
        map_index.AT(superregss.AT(m).AT(0).map) = m;
      }
    }
    
    // Combine regions from all maps
    vector<region_t> regs;
    for (int m=0; m<nmaps; ++m) {
      combine_regions(superregss.AT(m), regs);
    }

    // Create superregion tree structure
    {
      int const ncomps = regs.size();
      for (int c=0; c<ncomps; ++c) {
        assert(regs.AT(c).processor == -1);
        assert(not regs.AT(c).processors);
        pseudoregion_t const preg(regs.AT(c).extent, -1);
        regs.AT(c).processors = new ipfulltree(preg);
      }
    }
    
    // Set up new superregions
    // Note: The tree structure is now reachable from both the
    // superregions and the regions
    for (int m=0; m<nmaps; ++m) {
      superregss.AT(m).clear();
    }
    {
      int const ncomps = regs.size();
      for (int c=0; c<ncomps; ++c) {
        int const m = map_index.AT(regs.AT(c).map);
        superregss.AT(m).push_back(regs.AT(c));
      }
    }
    
    // Split onto processes
    vector<vector<region_t> > split_regss;
    int const nprocs = CCTK_nProcs(cctkGH);
    int const nworkers = nprocs;
    balance(regs, split_regss, nworkers,
            maximum_imbalance, same_number_of_components_on_each_process);
    
    // Assign process numbers, and make the tree structure
    // inaccessible from the regions
    for (int p=0; p<nprocs; ++p) {
      int const ncomps = split_regss.AT(p).size();
      for (int c=0; c<ncomps; ++c) {
        // Set process number in tree structure
        assert(split_regss.AT(p).AT(c).processors);
        assert(split_regss.AT(p).AT(c).processors->is_leaf());
        assert(split_regss.AT(p).AT(c).processors->payload().component == -1);
        split_regss.AT(p).AT(c).processors->payload().component = p;
        // Set process number in region
        split_regss.AT(p).AT(c).processors = NULL;
        assert(split_regss.AT(p).AT(c).processor == -1);
        split_regss.AT(p).AT(c).processor = p;
      }
    }
    
    // Distribute by map, implicitly assigning component numbers
    for (int p=0; p<nprocs; ++p) {
      int const ncomps = split_regss.AT(p).size();
      for (int c=0; c<ncomps; ++c) {
        int const m = map_index.AT(split_regss.AT(p).AT(c).map);
        regss.AT(m).push_back(split_regss.AT(p).AT(c));
      }
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  
  
  
  static void
  MakeMultigridBoxes (cGH const * const cctkGH,
                      int const m,
                      ibbox    const & base,
                      region_t const & reg,
                      vector<region_t> & regs)
  {
    DECLARE_CCTK_PARAMETERS;
    
    regs.resize (mglevels, reg);
    if (mglevels > 1) {
      // boundary offsets
      jjvect nboundaryzones, is_internal, is_staggered, shiftout;
      if (domain_from_multipatch and
          CCTK_IsFunctionAliased ("MultiPatch_GetBoundarySpecification"))
      {
        check (not MultiPatch_GetBoundarySpecification
               (m, 2*dim, &nboundaryzones[0][0], &is_internal[0][0],
                &is_staggered[0][0], &shiftout[0][0]));
      } else {
        check (not GetBoundarySpecification
               (2*dim, &nboundaryzones[0][0], &is_internal[0][0],
                &is_staggered[0][0], &shiftout[0][0]));
      }
      // (distance in grid points between the exterior and the physical boundary)
      iivect offset;
      for (int d=0; d<dim; ++d) {
        for (int f=0; f<2; ++f) {
          assert (not is_staggered[d][f]);
          offset[d][f] = (+ (is_internal[d][f] ? 0 : nboundaryzones[d][f] - 1)
                          + shiftout[d][f]);
        }
      }
      vector<ibbox> bases(mglevels);
      bases.AT(0) = base;
      for (int ml=1; ml<mglevels; ++ml) {
        // next finer base
        ivect const fbaselo = bases.AT(ml-1).lower();
        ivect const fbasehi = bases.AT(ml-1).upper();
        ivect const fbasestr = bases.AT(ml-1).stride();
        // this base
        ivect const basestr = fbasestr * mgfact;
        ivect const baselo = fbaselo + (xpose(offset)[0] - ivect(mgfact) * xpose(offset)[0]) * fbasestr;
        ivect const basehi = fbasehi + (xpose(offset)[1] - ivect(mgfact) * xpose(offset)[1]) * fbasestr;
        ivect const baselo1 = baselo;
        ivect const basehi1 = baselo1 + (basehi - baselo1) / basestr * basestr;
        bases.AT(ml) = ibbox(baselo1, basehi1, basestr);
        // next finer grid
        ivect const flo = regs.AT(ml-1).extent.lower();
        ivect const fhi = regs.AT(ml-1).extent.upper();
        ivect const fstr = regs.AT(ml-1).extent.stride();
        // this grid
        ivect const str = fstr * mgfact;
        ivect const lo = flo + either (reg.outer_boundaries[0],   (xpose(offset)[0] - ivect(mgfact) * xpose(offset)[0]) * fstr, ivect(0));
        ivect const hi = fhi + either (reg.outer_boundaries[1], - (xpose(offset)[1] - ivect(mgfact) * xpose(offset)[1]) * fstr, ivect(0));
        ivect const lo1 = baselo1 + (lo - baselo1 + str - 1) / str * str;
        ivect const hi1 = lo1 + (hi - lo1) / str * str;
        regs.AT(ml).extent = ibbox(lo1, hi1, str);
      }
    }
  }
  
  void
  MakeMultigridBoxes (cGH const * const cctkGH,
                      int const m,
                      vector<vector<region_t> > const & regss,
                      vector<vector<vector<region_t> > > & regsss)
  {
    regsss.resize (mglevels);
    for (int ml=0; ml<mglevels; ++ml) {
      regsss.AT(ml).resize (regss.size());
      for (int rl=0; rl<(int)regss.size(); ++rl) {
        regsss.AT(ml).AT(rl).resize (regss.AT(rl).size());
      }
    }
    
    for (int rl=0; rl<(int)regss.size(); ++rl) {
      
      ibbox base;
      for (int c=0; c<(int)regss.AT(rl).size(); ++c) {
        base = base.expanded_containing(regss.AT(rl).AT(c).extent);
      }
      
      for (int c=0; c<(int)regss.AT(rl).size(); ++c) {
        vector<region_t> mg_regs;
        MakeMultigridBoxes (cctkGH, m, base, regss.AT(rl).AT(c), mg_regs);
        assert ((int)mg_regs.size() == mglevels);
        for (int ml=0; ml<mglevels; ++ml) {
          regsss.AT(ml).AT(rl).AT(c) = mg_regs.AT(ml);
        }
        
      } // for c
      
    } // for rl
  }
  
  void
  MakeMultigridBoxesMaps (cGH const * const cctkGH,
                          vector<vector<vector<region_t> > > const & regsss,
                          vector<vector<vector<vector<region_t> > > > & regssss)
  {
    for (int m = 0; m < maps; ++m) {
      MakeMultigridBoxes (cctkGH, m, regsss.AT(m), regssss.AT(m));
    } // for m
  }
  
  
  
  static
  void
  ClassifyPoints (cGH const * const cctkGH, int const rl)
  {
    // negative: needs to be set explicitly (e.g. boundary)
    // zero:     unused (e.g. ghost)
    // positive: needs to be evolved
    // -1:      boundary point (needs to be set explicitly)
    //  0:      unused (e.g. ghost point, or restriction target)
    //  n=1..N: evolved, used for integrator substeps i<=n
    //          (i=N..1, counting backwards; see MoL documentation) 
    //          i.e.: n=1: used between time steps (i.e., should be visualised)
    //                n>1: used only while time stepping (e.g. buffer zones)
    
    BEGIN_META_MODE(cctkGH) {
      BEGIN_MGLEVEL_LOOP (cctkGH) {
        ENTER_LEVEL_MODE (cctkGH, rl) {
          BEGIN_LOCAL_MAP_LOOP(cctkGH, CCTK_GF) {
            BEGIN_LOCAL_COMPONENT_LOOP(cctkGH, CCTK_GF) {
              DECLARE_CCTK_ARGUMENTS;
#pragma omp parallel
              CCTK_LOOP3_ALL(CarpetClassifyPoints, cctkGH, i,j,k) {
                int const ind = CCTK_GFINDEX3D (cctkGH, i, j, k);
                point_class[ind] = 1;
              } CCTK_ENDLOOP3_ALL(CarpetClassifyPoints);
            } END_LOCAL_COMPONENT_LOOP;
          } END_LOCAL_MAP_LOOP;
        } LEAVE_LEVEL_MODE;
      } END_MGLEVEL_LOOP;
    } END_META_MODE;
  }
  
} // namespace Carpet