aboutsummaryrefslogtreecommitdiff
path: root/CarpetDev/CarpetAdaptiveRegrid/src/CAR.cc
blob: 2e4388921206b8c31f9dd3860dbce114e5903406 (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
#include <cassert>
#include <cmath>
#include <algorithm>
#include <iostream>
#include <sstream>
#include <string>
#include <vector>
#include <stack>

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

#include "gh.hh"
#include "vect.hh"

#include "carpet.hh"
#include "CAR.hh"

extern "C" {
  static const char* rcsid = "$Header:$";
  CCTK_FILEVERSION(Carpet_CarpetAdaptiveregrid_regrid_cc);
}

//
// For the moment this is going to live as one large file with one
// large routine. However, it should be broken up into pieces later
// when possible.
//

namespace CarpetAdaptiveRegrid {

  //
  // The following "local_" variables contain copies of the list of
  // boxes and the outer boundaries. However, these should be
  // PROCESSOR INDEPENDENT; that is, these are never used in calls to
  // SplitRegions. All box computation is done with respect to these
  // local variables, so should always be independent of the number of
  // processors. The results are split across processors and put into
  // the "standard" variables before being passed back to Carpet.
  //

  static gh::mexts local_bbsss;
  static gh::rbnds local_obss;
  
  //
  // Keep track of the last iteration on which we were called. This
  // means that we can return if we've been called on this timestep
  // (because we want to ensure proper nesting we regrid all levels
  // finer than the one on which we are called. Therefore if we're
  // called on the same iteration, as long as we're called in
  // coarse->fine order, we only need to regrid once.) 
  //
  // FIXME: Note that this may well cause problems with checkpoint /
  // restart; should this be moved to a grid scalar?
  //

  static int last_iteration = -1;  
  
  using namespace std;
  using namespace Carpet;

  //
  // The following Fortran helper routines are defined in
  // CAR_utils.F90.
  //

  extern "C" {
    void CCTK_FCALL CCTK_FNAME(copy_mask)
      (const int& snx, const int& sny, const int& snz,
       const int* smask, const int sbbox[3][3],
       const int& dnx, const int& dny, const int& dnz,
       int* dmask, const int dbbox[3][3]);
    void CCTK_FCALL CCTK_FNAME(check_box)
      (const int& nx, const int& ny, const int& nz,
       const int* mask,
       int* sum_x, int* sum_y, int* sum_z,
       int* sig_x, int* sig_y, int* sig_z,
       const int bbox[3][3],
       int newbbox1[3][3], int newbbox2[3][3],
       const int& min_width, const CCTK_REAL& min_fraction,
       int& didit);
  }

  //
  // The following helper routines translate between real coordinates
  // and integer (Carpet bbox) coordinates on a given refinement
  // level. Basically copied from CarpetRegrid, "pos" <-> "position"
  // and "int" <-> "integer". Actually defined at the bottom of this
  // file.
  //

  ivect pos2int (const cGH* const cctkGH, const gh& hh,
                 const rvect & rpos, const int rl);
  rvect int2pos (const cGH* const cctkGH, const gh& hh,
                 const ivect & ipos, const int rl);
  
  //
  // The real routine that does all the work
  //

  CCTK_INT CarpetAdaptiveRegrid_Regrid (CCTK_POINTER_TO_CONST const cctkGH_,
                                        CCTK_POINTER const bbsss_,
                                        CCTK_POINTER const obss_,
                                        CCTK_POINTER const pss_,
                                        CCTK_INT force)
  {
    DECLARE_CCTK_PARAMETERS;
    
    const cGH * const cctkGH = (const cGH *) cctkGH_;

    //
    // The following variables are the important ones that Carpet uses
    // to set up the grids. They contain the bounding boxes (bb), the
    // outer boundary specifications (ob) and the processor numbers
    // (p). I believe the "s" stands for "set"; so the "obss" is a set
    // of sets of outer boundary specifications. However, each "s" is
    // implemented as a vector, not a set. For the ob and p the
    // indices on the vectors are the refinement levels and the maps;
    // for the bounding boxes the order is multigrid levels,
    // refinement levels, maps.
    //

    gh::mexts  & bbsss = * (gh::mexts  *) bbsss_;
    gh::rbnds  & obss  = * (gh::rbnds  *) obss_;
    gh::rprocs & pss   = * (gh::rprocs *) pss_;
    
    gh const & hh = *vhh.at(Carpet::map);
    
    assert (is_singlemap_mode());

    if (local_bbsss.empty()) { // It's the first call
      //
      // We set up the "local_" variables using just the base grid
      // which we get from looking at the base extent; of course, the
      // outer boundary set up is trivial.
      //
      // Note that this will need improving (in fact the whole
      // regridding mechanism need minor changes) in order to do mesh
      // refinement with multiple maps.
      //
      const ibbox& baseext = 
        vdd.at(Carpet::map)->bases.at(mglevel).at(reflevel).exterior;
      vector<ibbox> tmp_bbs;
      tmp_bbs.push_back (baseext);
      vector<bbvect> tmp_obs;
      tmp_obs.push_back (bbvect(true));
      vector<vector<ibbox> > tmp_bbss(1);
      vector<vector<bbvect> > tmp_obss(1);
      tmp_bbss.at(0) = tmp_bbs;
      tmp_obss.at(0) = tmp_obs;
      MakeMultigridBoxes(cctkGH, tmp_bbss, tmp_obss, local_bbsss);
      local_obss = tmp_obss;
      last_iteration = cctkGH->cctk_iteration;
      //
      // Having set up the base grid we then set any finer grids
      // according to the coordinate parameters, as if we were the
      // standard CarpetRegrid.
      //
      int do_recompose = 
        ManualCoordinateList (cctkGH, hh, bbsss, obss, pss, 
                              local_bbsss, local_obss);

      if (verbose) {
        ostringstream buf;
        buf << "Done with manual coordinate list. Total list is:"
            << endl << local_bbsss;
        CCTK_INFO(buf.str().c_str());
      }      

      return do_recompose;
    }

    // FIXME: We should check that the local reflevel "agrees"
    // with what is passed in.
    
    // In force mode (force == true) we do not check the
    // CarpetAdaptiveregrid parameters

    if (!force) {

      assert (regrid_every == -1 || regrid_every == 0
	      || regrid_every % maxmglevelfact == 0);
    
      // Return if no regridding is desired
      if (regrid_every == -1) return 0;
      
      // Return if we want to regrid during initial data only, and this
      // is not the time for initial data
      if (regrid_every == 0 && cctkGH->cctk_iteration != 0) return 0;

      // Return if we want to regrid regularly, but not at this time
      if (regrid_every > 0 && cctkGH->cctk_iteration != 0
	  && (cctkGH->cctk_iteration-1) % regrid_every != 0)
      {
	return 0;
      }
      
    }

    //
    // Even if force is set, there are times that I'm not going to do
    // any regridding, so be told.
    //

    if (reflevel == maxreflevels - 1) return 0;

    // Return if we want to regrid regularly, but not at this time
    if (regrid_every > 0 && cctkGH->cctk_iteration != 0
        && (cctkGH->cctk_iteration-1) % regrid_every != 0)
    {
      return 0;
    }

    // Return if we have already been called on this iteration
    if (cctkGH->cctk_iteration == last_iteration) {
      return 0;
    }
    else {
      last_iteration = cctkGH->cctk_iteration;
    }

    //
    // Set up a few local variables for later use. In particular we
    // need to keep track of the level and map on which we were
    // called. 
    //

    int do_recompose;
    do_recompose = 1;

    int sum_handle = CCTK_ReductionArrayHandle("sum");

    int called_on_ml = mglevel;
    int called_on_rl = reflevel;
    int called_on_grouptype = mc_grouptype;
    int called_on_map = carpetGH.map;
    
    int finest_current_rl = local_bbsss.at(0).size();
    finest_current_rl = min(finest_current_rl, (maxreflevels - 1));

    //
    // Loop over all levels finer than this one.
    //

    leave_singlemap_mode(const_cast<cGH *> (cctkGH));
    leave_level_mode(const_cast<cGH *> (cctkGH));
    for (int rl = called_on_rl; rl < finest_current_rl; ++rl) {
      enter_level_mode(const_cast<cGH *> (cctkGH), rl);
      enter_singlemap_mode(const_cast<cGH *> (cctkGH), called_on_map, called_on_grouptype);

      if (verbose) {
        ostringstream buf;
        buf << "Entering level " << rl << " (of " << finest_current_rl 
            << "), map " << called_on_map;
        CCTK_INFO(buf.str().c_str());
      }
      
      vector<ibbox> bbs = local_bbsss.at(mglevel).at(reflevel);
      
      stack<ibbox> final;
      
      vector<vector<ibbox> > bbss = bbsss.at(0);
      vector<vector<ibbox> > local_bbss = local_bbsss.at(0);
      
      ivect const reffact = spacereffacts.at(rl+1) / spacereffacts.at(rl);

      bool did_regrid = false;
            
      //
      // For later use get the physical domain specification in real
      // coordinates.
      //

      rvect physical_min, physical_max;
      rvect interior_min, interior_max;
      rvect exterior_min, exterior_max;
      rvect base_spacing;
      int ierr = GetDomainSpecification
        (dim, &physical_min[0], &physical_max[0],
         &interior_min[0], &interior_max[0],
         &exterior_min[0], &exterior_max[0], &base_spacing[0]);
      assert (!ierr);

      //
      // Loop over all boxes on this level.
      //
      
      for ( vector<ibbox>::const_iterator bbi = bbs.begin(); 
            bbi != bbs.end();
            ++bbi) 
      {
        
        ivect low = bbi->lower();
        ivect upp = bbi->upper();
        
        // low and upp now define the starting bbox.
        
        ibbox bb(low, upp, bbs.at(0).stride());
        
        if (verbose) {
          ostringstream buf;
          buf << "Found the local size of the box: " << endl << bb;
          CCTK_INFO(buf.str().c_str());
        }
        
        vector<int> local_mask(prod(bb.shape()/bb.stride()), 0),
          mask(prod(bb.shape()/bb.stride()), 0);
        
        if (veryverbose) {
          ostringstream buf;
          buf << "Allocated mask size: " << bb.shape()/bb.stride() 
              << " (points: " << prod(bb.shape()/bb.stride()) << ")";
          CCTK_INFO(buf.str().c_str());
        }

        //        
        // Setup the mask.
        // That is, for any errors that we can locally see that exceed
        // the threshold, set the mask to 1.
        //

        const ibbox& baseext = 
          vdd.at(Carpet::map)->bases.at(mglevel).at(reflevel).exterior;
        ivect imin = (bb.lower() - baseext.lower())/bb.stride(), 
          imax = (bb.upper() - baseext.lower())/bb.stride();
        
        BEGIN_LOCAL_COMPONENT_LOOP(cctkGH, CCTK_GF) 
        {
          const CCTK_REAL *error_var_ptr = 
            static_cast<const CCTK_REAL*>(CCTK_VarDataPtr(cctkGH, 
                                                          0, error_var));
          const CCTK_REAL *x_var_ptr = 
            static_cast<const CCTK_REAL*>(CCTK_VarDataPtr(cctkGH, 
                                                          0, "Grid::x"));
          const CCTK_REAL *y_var_ptr = 
            static_cast<const CCTK_REAL*>(CCTK_VarDataPtr(cctkGH, 
                                                          0, "Grid::y"));
          const CCTK_REAL *z_var_ptr = 
            static_cast<const CCTK_REAL*>(CCTK_VarDataPtr(cctkGH, 
                                                          0, "Grid::z"));
          

          // These can actually be negative if the parent shrinks.
          // Of course, the final grid should still be properly
          // nested...
          //          assert(all(imin >= 0));
          //          assert(all(imax >= 0));
          // FIXME: Why should the following assert be true?
          //      assert(all(imax < ivect::ref(cctkGH->cctk_lsh)));
          assert(all(imin <= imax));
          
          for (int k = 0; k < cctkGH->cctk_lsh[2]; ++k) {
            for (int j = 0; j < cctkGH->cctk_lsh[1]; ++j) {
              for (int i = 0; i < cctkGH->cctk_lsh[0]; ++i) {
                int index = CCTK_GFINDEX3D(cctkGH, i, j, k);
                CCTK_REAL local_error = abs(error_var_ptr[index]);
                if (local_error > max_error) {
                  int ii = i + cctkGH->cctk_lbnd[0] - imin[0];
                  int jj = j + cctkGH->cctk_lbnd[1] - imin[1];
                  int kk = k + cctkGH->cctk_lbnd[2] - imin[2];
                  // Check that this point actually intersects with 
                  // this box (if this component was actually a
                  // different grid on the same processor, it need not)
                  if ( (ii >= 0) and (jj >= 0) and (kk >= 0) and 
                       (ii <= imax[0] - imin[0]) and
                       (jj <= imax[1] - imin[1]) and
                       (kk <= imax[2] - imin[2]) )
                  {
                    assert (ii >= 0);
                    assert (jj >= 0);
                    assert (kk >= 0);
                    assert (ii <= imax[0] - imin[0]);
                    assert (jj <= imax[1] - imin[1]);
                    assert (kk <= imax[2] - imin[2]);
                    int mindex = ii + 
                      (imax[0] - imin[0] + 1)*
                      (jj + (imax[1] - imin[1] + 1) * kk);
                    local_mask[mindex] = 1;
                    if (veryverbose) {
                      CCTK_VInfo(CCTK_THORNSTRING, "In error at point"
                                 "\n(%g,%g,%g) [%d,%d,%d] [[%d,%d,%d]]",
                                 (double) x_var_ptr[index],
                                 (double) y_var_ptr[index],
                                 (double) z_var_ptr[index],
                                 ii, jj, kk,
                                 i, j, k);
                    }
                  }
                }
              }
            }
          }
        } END_LOCAL_COMPONENT_LOOP;
        
        //
        // Check the error on child level, if such a level exists
        // Also only worry if there's a grandchild level.
        // This should fix the "orphaned grandchild" problem
        //

        if ((int)local_bbss.size() > reflevel+2) {
                
          int currentml = mglevel;
          int currentrl = reflevel;
          int currentmap = carpetGH.map;

          leave_singlemap_mode(const_cast<cGH *> (cctkGH));
          leave_level_mode(const_cast<cGH *> (cctkGH));
      
          enter_level_mode(const_cast<cGH *> (cctkGH), currentrl + 1);
          enter_singlemap_mode(const_cast<cGH *> (cctkGH), currentmap, CCTK_GF);

          const ibbox& child_baseext = 
            vdd.at(Carpet::map)->bases.at(mglevel).at(reflevel).exterior;
          ivect child_levoff = child_baseext.lower()/(bb.stride()/reffact);
                          
          if (verbose) {
            ostringstream buf;
            buf << "Checking for errors on child level " 
                << reflevel << " map " << currentmap;
            CCTK_INFO(buf.str().c_str());
          }
      
          BEGIN_LOCAL_COMPONENT_LOOP(cctkGH, CCTK_GF) {
            const CCTK_REAL *error_var_ptr = 
              static_cast<const CCTK_REAL*>(CCTK_VarDataPtr(cctkGH, 
                                                            0, error_var));
            const CCTK_REAL *x_var_ptr = 
              static_cast<const CCTK_REAL*>(CCTK_VarDataPtr(cctkGH, 
                                                            0, "Grid::x"));
            const CCTK_REAL *y_var_ptr = 
              static_cast<const CCTK_REAL*>(CCTK_VarDataPtr(cctkGH, 
                                                            0, "Grid::y"));
            const CCTK_REAL *z_var_ptr = 
              static_cast<const CCTK_REAL*>(CCTK_VarDataPtr(cctkGH, 
                                                            0, "Grid::z"));
            
            //            assert(all(imin >= 0));
            //            assert(all(imax >= 0));
            // FIXME: Why should the following assert be true?
            //      assert(all(imax < ivect::ref(cctkGH->cctk_lsh)));
            assert(all(imin <= imax));

            if (verbose) {
              ostringstream buf;
              buf << "Checking for errors on child level " 
                  << reflevel << " map " << currentmap
                  << " component " << component << endl
                  << "Mask: " << imin << imax << endl
                  << "Component: " << ivect::ref(cctkGH->cctk_lbnd)
                  << ivect::ref(cctkGH->cctk_lsh) << endl
                  << "level offset: " << child_levoff;
              CCTK_INFO(buf.str().c_str());
            }
            
            for (int k = 0; k < cctkGH->cctk_lsh[2]; ++k) {
              for (int j = 0; j < cctkGH->cctk_lsh[1]; ++j) {
                for (int i = 0; i < cctkGH->cctk_lsh[0]; ++i) {
                  int index = CCTK_GFINDEX3D(cctkGH, i, j, k);
                  CCTK_REAL local_error = abs(error_var_ptr[index]);
                  if (local_error > max_error) {
                    // Correct for the change in level !
                    int ii = (i + cctkGH->cctk_lbnd[0] 
                                   - cctkGH->cctk_nghostzones[0]
                                   + child_levoff[0]) / reffact[0] - 
                      imin[0];
                    int jj = (j + cctkGH->cctk_lbnd[1] 
                                   - cctkGH->cctk_nghostzones[1]
                                   + child_levoff[1]) / reffact[1] - 
                      imin[1];
                    int kk = (k + cctkGH->cctk_lbnd[2] 
                                   - cctkGH->cctk_nghostzones[2]
                                   + child_levoff[1]) / reffact[2] - 
                      imin[2];
                    // Check that this point actually intersects with 
                    // this box (if this component was actually a
                    // different grid on the same processor, it need not)
                    if ( (ii >= 0) and (jj >= 0) and (kk >= 0) and 
                         (ii <= imax[0] - imin[0]) and
                         (jj <= imax[1] - imin[1]) and
                         (kk <= imax[2] - imin[2]) )
                    {
                      assert (ii >= 0);
                      assert (jj >= 0);
                      assert (kk >= 0);
                      assert (ii <= imax[0] - imin[0]);
                      assert (jj <= imax[1] - imin[1]);
                      assert (kk <= imax[2] - imin[2]);
                      int mindex = ii + 
                        (imax[0] - imin[0] + 1)*
                        (jj + (imax[1] - imin[1] + 1) * kk);
                      local_mask[mindex] = 1;
                      if (veryverbose) {
                        CCTK_VInfo(CCTK_THORNSTRING, "In error at point"
                                   "\n(%g,%g,%g) [%d,%d,%d] [[%d,%d,%d]]",
                                   (double) x_var_ptr[index],
                                   (double) y_var_ptr[index],
                                   (double) z_var_ptr[index],
                                   ii, jj, kk,
                                   i, j, k);
                      }
                    }
                  }
                }
              }
            }
          } END_LOCAL_COMPONENT_LOOP;
        
          leave_singlemap_mode(const_cast<cGH *> (cctkGH));
          leave_level_mode(const_cast<cGH *> (cctkGH));
          
          enter_level_mode(const_cast<cGH *> (cctkGH), currentrl);
          enter_singlemap_mode(const_cast<cGH *> (cctkGH), currentmap, CCTK_GF);

        }

        //        
        // Reduce errors (MPI sum).
        // This sets up the mask globally, although the mask is now
        // non-zero where in error instead of just being 1. It's still
        // 0 at points not in error.
        //

        int ierr = 
          CCTK_ReduceLocArrayToArray1D (cctkGH,
                                        -1,
                                        sum_handle,
                                        &local_mask.front(),
                                        &mask.front(),
                                        local_mask.size(),
                                        CCTK_VARIABLE_INT );
        
        if (ierr) {
          ostringstream buf;
          buf << "The reduction on level " << reflevel << " failed "
              << "with error " << ierr;
          CCTK_WARN(0, buf.str().c_str());
        }

        //
        // Set errors to 0/1; so, where the mask indicates that a
        // point is in error, set it precisely to 1.
        //

        for (int k = 0; k < imax[2] - imin[2] + 1; k++) {
          for (int j = 0; j < imax[1] - imin[1] + 1; j++) {
            for (int i = 0; i < imax[0] - imin[0] + 1; i++) {
              int index =  i + 
                (imax[0] - imin[0] + 1)*(j + (imax[1] - imin[1] + 1) * k);
              if (mask[index]) {
                mask[index] = 1;
              }
            }
          }
        }
        
        //
        // Pad the errors: stage 1 - buffer points marked as 2.
        //

        for (int k = 0; k < imax[2] - imin[2] + 1; k++) {
          for (int j = 0; j < imax[1] - imin[1] + 1; j++) {
            for (int i = 0; i < imax[0] - imin[0] + 1; i++) {
              int index =  i + 
                (imax[0] - imin[0] + 1)*(j + (imax[1] - imin[1] + 1) * k);
              if (mask[index] == 1) {
                for (int kk = max((int)(k - pad), 0);
                     kk < min((int)(k + pad + 1), imax[2] - imin[2] + 1);
                       ++kk)
                {
                  for (int jj = max((int)(j - pad), 0); 
                       jj < min((int)(j + pad + 1), imax[1] - imin[1] + 1);
                       ++jj)
                  {
                    for (int ii = max((int)(i - pad), 0); 
                         ii < min((int)(i + pad + 1), imax[0] - imin[0] + 1);
                         ++ii)
                    {
                      int mindex = ii + 
                        (imax[0] - imin[0] + 1)*
                        (jj + (imax[1] - imin[1] + 1) * kk);
                      if (!mask[mindex]) mask[mindex] = 2;
                    }
                  }
                }
              }
            }
          }
        }
        //
        // stage 2: all buffer points marked truly in error.
        // Also mark if there are any errors.
        //
        bool should_regrid = false;
        for (int k = 0; k < imax[2] - imin[2] + 1; k++) {
          if (maskpicture) {
            cout << endl << "k = " << k << endl;
          }
          for (int j = 0; j < imax[1] - imin[1] + 1; j++) {
            if (maskpicture) {
              cout << endl;
            }
            for (int i = 0; i < imax[0] - imin[0] + 1; i++) {
              int index =  i + 
                (imax[0]-imin[0] + 1)*(j + (imax[1] - imin[1] + 1) * k);
              if (mask[index] > 1) mask[index] = 1;
              if ((veryverbose) and (mask[index])) {
                CCTK_VInfo(CCTK_THORNSTRING, "Mask set at point"
                           "\n[%d,%d,%d]",
                           i, j, k);
              } 
              if (maskpicture) {
                if (mask[index]) {
                  cout << "x";
                }
                else {
                  cout << " ";
                }
              }
              should_regrid |= (mask[index]);
              did_regrid |= should_regrid;
            }
          }
        }    
        
        if (verbose) {
          ostringstream buf;
          buf << "Finished looking for errors on level " 
              << reflevel << endl << "should_regrid " << should_regrid
              << " did_regrid " << did_regrid;
          CCTK_INFO(buf.str().c_str());
        }
        
        //
        // For this box on this level we now have the marked
        // points. If there are any errors then we should actually
        // create the new boxes. Starting from the actual bbox, placed
        // in the "todo" stack, we iterate over the "todo" stack
        // creating new boxes which are either accepted (and hence
        // placed on the "final" stack) or still need work done to
        // them (hence placed on the "todo" stack).
        //

        if (should_regrid) {
          
          stack<ibbox> todo;
          
          todo.push(bb);
          
          // Define vector of masks (contains error mask)
          
          stack<vector<int> > masklist;
          
          masklist.push(mask);

          //          
          // Loop over all entries in todo:
          //   Setup appropriate 1d array memory
          //   Call fortran routine
          //   If return is:
          //     zero: add bbox to final
          //     one:  add new bbox to todo and assoc mask to masklist
          //     two:  add both new bboxs to todo and assoc masks to
          //           masklist 
          //

          while (!todo.empty())
          {
            
            ibbox bb = todo.top(); todo.pop();
            vector<int> mask = masklist.top(); masklist.pop();
            
            int nx = bb.shape()[0]/bb.stride()[0];
            int ny = bb.shape()[1]/bb.stride()[1];
            int nz = bb.shape()[2]/bb.stride()[2];
            
            if (verbose) {
              ostringstream buf;
              buf << "todo loop. Box: " << endl << bb;
              CCTK_INFO(buf.str().c_str());
            }
            
            vector<int> sum_x(nx, 0);
            vector<int> sig_x(nx, 0);
            vector<int> sum_y(ny, 0);
            vector<int> sig_y(ny, 0);
            vector<int> sum_z(nz, 0);
            vector<int> sig_z(nz, 0);
            
            int fbbox[3][3], fbbox1[3][3], fbbox2[3][3];
            
            for (int d = 0; d < 3; ++d) {
              fbbox[0][d] = bb.lower()[d];
              fbbox[1][d] = bb.upper()[d];
              fbbox[2][d] = bb.stride()[d];
            }
              
            int didit;

            //
            // This is the actual Fortran routine that does the work
            // of the Berger-Rigoutsos algorithm.
            //
              
            CCTK_FNAME(check_box)(nx, ny, nz,
                                  &mask.front(),
                                  &sum_x.front(), &sum_y.front(),
                                  &sum_z.front(),
                                  &sig_x.front(), &sig_y.front(),
                                  &sig_z.front(),
                                  fbbox,
                                  fbbox1, fbbox2,
                                  min_width, min_fraction,
                                  didit);
            
            if (didit == 0) { // Box was accepted
              
              final.push(bb);          
              
              if (verbose) {
                ostringstream buf;
                buf << "todo loop. Box pushed to final: " 
                    << endl << bb;
                CCTK_INFO(buf.str().c_str());
              }
            }
            else if (didit == 1) { // Box was replaced by a new single
                                   // box
              ibbox newbbox1(ivect::ref(&fbbox1[0][0]),
                             ivect::ref(&fbbox1[1][0]),
                             ivect::ref(&fbbox1[2][0]));
              todo.push(newbbox1);
              
              int dnx = newbbox1.shape()[0]/newbbox1.stride()[0];
              int dny = newbbox1.shape()[1]/newbbox1.stride()[1];
              int dnz = newbbox1.shape()[2]/newbbox1.stride()[2];
              
              vector<int>  
                newmask1(prod(newbbox1.shape()/newbbox1.stride()), 0);
              
              CCTK_FNAME(copy_mask)(nx, ny, nz,
                                    &mask.front(), fbbox,
                                    dnx, dny, dnz,
                                    &newmask1.front(), fbbox1);
              masklist.push(newmask1);
              
              if (verbose) {
                ostringstream buf;
                buf << "todo loop. New (single) box created: " 
                    << endl << newbbox1;
                CCTK_INFO(buf.str().c_str());
              }
            }
            else if (didit == 2) { // Box was replaced with two boxes
              
              ibbox newbbox1(ivect::ref(&fbbox1[0][0]),
                             ivect::ref(&fbbox1[1][0]),
                             ivect::ref(&fbbox1[2][0]));
              todo.push(newbbox1);
              ibbox newbbox2(ivect::ref(&fbbox2[0][0]),
                             ivect::ref(&fbbox2[1][0]),
                             ivect::ref(&fbbox2[2][0]));
              todo.push(newbbox2);
              
              int dnx = newbbox1.shape()[0]/newbbox1.stride()[0];
              int dny = newbbox1.shape()[1]/newbbox1.stride()[1];
              int dnz = newbbox1.shape()[2]/newbbox1.stride()[2];
              
              vector<int>  
                newmask1(prod(newbbox1.shape()/newbbox1.stride()), 0);
              
              CCTK_FNAME(copy_mask)(nx, ny, nz,
                                    &mask.front(), fbbox,
                                    dnx, dny, dnz,
                                    &newmask1.front(), fbbox1);
              masklist.push(newmask1);
              
              dnx = newbbox2.shape()[0]/newbbox2.stride()[0];
              dny = newbbox2.shape()[1]/newbbox2.stride()[1];
              dnz = newbbox2.shape()[2]/newbbox2.stride()[2];
              
              vector<int>  
                newmask2(prod(newbbox2.shape()/newbbox2.stride()), 0);
              
              CCTK_FNAME(copy_mask)(nx, ny, nz,
                                    &mask.front(), fbbox,
                                    dnx, dny, dnz,
                                    &newmask2.front(), fbbox2);
              masklist.push(newmask2);
              
              if (verbose) {
                ostringstream buf;
                buf << "todo loop. New (double) box created. Box 1: " 
                    << endl << newbbox1
                    << "                                     Box 2: "
                    << endl << newbbox2;
                CCTK_INFO(buf.str().c_str());
              }
            }
            else {
              CCTK_WARN(0, "The fortran routine must be confused.");
            }
            
          } // loop over todo vector (boxes needing to be done).
        } // should regrid.
      } // Loop over boxes on the parent grid.
      
      if (did_regrid) { // If we actually did something, reconvert the
                        // boxes to correct Carpet style, plus correct
                        // the boundaries.
        // Fixup the stride
        vector<ibbox> newbbs;
        vector<bbvect> obs;
        while (! final.empty()) {
          ibbox bb = final.top(); final.pop();
          
          if (veryverbose) {
            ostringstream buf;
            buf << "Looping over the final list. Box is:"
                << endl << bb;
            CCTK_INFO(buf.str().c_str());
          }
          
          ivect ilo = bb.lower();
          ivect ihi = bb.upper();
          rvect lo = int2pos(cctkGH, hh, ilo, reflevel);
          rvect hi = int2pos(cctkGH, hh, ihi, reflevel);
          rvect str = base_spacing * 
            ipow((CCTK_REAL)mgfact, basemglevel) /
            rvect (spacereffacts.at(reflevel));
          rbbox newbbcoord(lo, hi, str);
          
          if (veryverbose) {
            ostringstream buf;
            buf << "Dealing with boundaries. Coord box is:"
                << endl << newbbcoord;
            CCTK_INFO(buf.str().c_str());
          }
          
          // Set the correct ob here.
          
          bbvect ob(false);
          for (int d=0; d<dim; ++d) {
            assert (mglevel==0);
            
            // Find the size of the physical domain
            
            rvect const spacing = base_spacing * 
              ipow((CCTK_REAL)mgfact, basemglevel) /
              rvect (ipow(reffact, reflevel+1));
            ierr = ConvertFromPhysicalBoundary
              (dim, &physical_min[0], &physical_max[0],
               &interior_min[0], &interior_max[0],
               &exterior_min[0], &exterior_max[0], &spacing[0]);
            assert (!ierr);
            
            // If need be clip the domain
            
            rvect lo = newbbcoord.lower();
            if (newbbcoord.lower()[d] < physical_min[d]) {
              lo[d] = exterior_min[d];
            }
            rvect up = newbbcoord.upper();
            if (newbbcoord.upper()[d] > physical_max[d]) {
              up[d] = exterior_max[d];
            }
            rvect str = newbbcoord.stride();
            
            // Set the ob if outside the physical domain
            
            ob[d][0] = 
              abs(lo[d] - exterior_min[d]) < 1.0e-6 * spacing[d];
            ob[d][1] = 
              abs(up[d] - exterior_max[d]) < 1.0e-6 * spacing[d];
            
            if (veryverbose) {
              ostringstream buf;
              buf << "Done clipping domain:"
                  << endl << lo << endl << up << endl << str;
              CCTK_INFO(buf.str().c_str());
            } 
            
            // Check that the striding is correct.
            
            CCTK_REAL remainder = fmod((up[d] - lo[d]), str[d])/str[d];
            
            if ( abs(remainder) > 1.e-6 ) {
              if (ob[d][0]) {
                up[d] += str[d] * (1 - remainder);
              }
              else if (ob[d][1]) {
                lo[d] -= str[d] * remainder;
              }
            }
            
            if (veryverbose) {
              ostringstream buf;
              buf << "Corrected coords for striding:"
                  << endl << lo << endl << up << endl << str;
              CCTK_INFO(buf.str().c_str());
            } 
            
            newbbcoord = rbbox(lo, up, str);
          }
          if (verbose) {
            ostringstream buf;
            buf << "Done dealing with boundaries. Coord box is:"
                << endl << newbbcoord << endl
                << "obox is:" << endl << ob;
            CCTK_INFO(buf.str().c_str());
          }

          //          
          // Convert back to integer coordinates
          // We have to do this on the fine grid to ensure that
          // it is correct for an outer boundary with odd numbers
          // of ghost zones where the bbox does not align with the
          // parent. 
          //

          ilo = pos2int(cctkGH, hh, newbbcoord.lower(), reflevel+1);
          ihi = pos2int(cctkGH, hh, newbbcoord.upper(), reflevel+1);
          ivect istr = bb.stride() / reffact;
          
          // Check that the width is sufficient
          // This can only be too small if the domain was clipped
          for (int d=0; d < dim; ++d) {
            if (ihi[d] - ilo[d] < min_width * istr[d]) {
              if (ob[d][0]) {
                if (ob[d][1]) {
                  CCTK_WARN(0, "The domain is too small?!");
                }
                ihi[d] = ilo[d] + min_width * istr[d];
              }
              else if (ob[d][1]) {
                if (ob[d][0]) {
                  CCTK_WARN(0, "The domain is too small?!");
                }
                ilo[d] = ihi[d] - min_width * istr[d];
              }
              else {
                ostringstream buf;
                buf << "The grid is unclipped and too small?" << endl 
                    << ilo << endl << ihi << endl << istr << endl << d;
                CCTK_WARN(0, buf.str().c_str());
              }
            }
          }
          
          if (veryverbose) {
            ostringstream buf;
            buf << "Corrected integer coords for min_width:"
                << endl << ilo << endl << ihi << endl << istr;
            CCTK_INFO(buf.str().c_str());
          } 
          
          ibbox newbb(ilo, ihi, istr);          
          
          if (verbose) {
            ostringstream buf;
            buf << "After dealing with boundaries. Final box is:"
                << endl << newbb;
            CCTK_INFO(buf.str().c_str());
          }
          
          newbbs.push_back (newbb);
          obs.push_back(ob);
        }
        
        
        // FIXME: check if the newbbs is really different
        // from the current bbs
        //        if not, set do_recompose = 0
        bbs = newbbs;
        
        // Set local bbss

        //
        // This mess ensures that the local bbsss is correctly set for
        // multigrid levels, then splits over regions (remember that
        // the local_bbsss is processor independent), then goes
        // through the whole multigrid procedure with the processor
        // dependent bbsss.
        //
        
        if ((int)bbss.size() < reflevel+2) {
          if (verbose) {
            CCTK_INFO("Adding new refinement level");
          }
          local_bbss.resize(reflevel+2);
          bbss.resize(reflevel+2);
          local_obss.resize(reflevel+2);
          obss.resize(reflevel+2);
          pss.resize(reflevel+2);
        }
        local_bbss.at(reflevel+1) = bbs;
        local_obss.at(reflevel+1) = obs;
        MakeMultigridBoxes (cctkGH, local_bbss, local_obss, local_bbsss);
        
        // make multiprocessor aware
        gh::cprocs ps;
        SplitRegions (cctkGH, bbs, obs, ps);    
        
        bbss.at(reflevel+1) = bbs;
        obss.at(reflevel+1) = obs;
        pss.at(reflevel+1) = ps;
        
      } // did_regrid?
      else
      {
        if ((int)local_bbss.size() > reflevel+1) {
          if (verbose) {
            CCTK_INFO("Removing refinement level");
          }
        }
        local_bbss.resize(reflevel+1);
        bbss.resize(reflevel+1);
        local_obss.resize(reflevel+1);
        obss.resize(reflevel+1);
        // Set local bbsss
        MakeMultigridBoxes (cctkGH, local_bbss, local_obss, local_bbsss);
        
        pss.resize(reflevel+1);
        
        do_recompose = 1;
      }
        
      // make multigrid aware
      MakeMultigridBoxes (cctkGH, bbss, obss, bbsss);
      
      leave_singlemap_mode(const_cast<cGH *> (cctkGH));
      leave_level_mode(const_cast<cGH *> (cctkGH));
      
    } 

    enter_level_mode(const_cast<cGH *> (cctkGH), called_on_rl);
    enter_singlemap_mode(const_cast<cGH *> (cctkGH), called_on_map, called_on_grouptype);

    if (verbose) {
      ostringstream buf;
      buf << "Done with it all. Iteration " << cctkGH->cctk_iteration
          << " level " << reflevel << endl << "Total list is:"
          << endl << local_bbsss;
      CCTK_INFO(buf.str().c_str());
    }      
    
    return do_recompose;
  }
  
  ivect pos2int (const cGH* const cctkGH, const gh& hh,
                 const rvect & rpos, const int rl)
  {
    rvect global_lower, global_upper;
    for (int d=0; d<dim; ++d) {
      const int ierr = CCTK_CoordRange
	(cctkGH, &global_lower[d], &global_upper[d], d+1, 0, "cart3d");
      if (ierr<0) {
	global_lower[d] = 0;
	global_upper[d] = 1;
      }
    }
    const ivect global_extent (hh.baseextent.upper() - hh.baseextent.lower());
    
    const rvect scale  = rvect(global_extent) / (global_upper - global_lower);
    const ivect levfac = hh.reffacts.at(rl);
    assert (all (hh.baseextent.stride() % levfac == 0));
    const ivect istride = hh.baseextent.stride() / levfac;
    
    const ivect ipos
      = (ivect(floor((rpos - global_lower) * scale / rvect(istride)
                     + (CCTK_REAL) 0.5))
         * istride);
    
    return ipos;
  }
  
  rvect int2pos (const cGH* const cctkGH, const gh& hh,
                 const ivect & ipos, const int rl)
  {
    rvect global_lower, global_upper;
    for (int d=0; d<dim; ++d) {
      const int ierr = CCTK_CoordRange
	(cctkGH, &global_lower[d], &global_upper[d], d+1, 0, "cart3d");
      
      if (ierr<0) {
	global_lower[d] = 0;
	global_upper[d] = 1;
      }
    }
    const ivect global_extent (hh.baseextent.upper() - hh.baseextent.lower());
    
    const rvect scale  = rvect(global_extent) / (global_upper - global_lower);
    const ivect levfac = hh.reffacts.at(rl);
    assert (all (hh.baseextent.stride() % levfac == 0));
    const ivect istride = hh.baseextent.stride() / levfac;    
    
    const rvect rpos
      = rvect(ipos) / scale + global_lower;
    
    return rpos;
  }
  
} // namespace CarpetAdaptiveRegrid