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

#include <carpet.hh>

#include "boundary.hh"
#include "property.hh"

#include <typeinfo>

// Consistency properties for the grid structure



namespace CarpetRegrid2 {
  
  using namespace std;
  using namespace Carpet;
  
  
  
  // Each property consists of a test, which returns true or false
  // depending on whether the property is satisfied, and an action
  // that enforces the property.
  
  bool property::
  test (gh const& hh, dh const& dd,
        level_boundary const& bnd,
        vector<ibset> const& regions, int const rl)
  {
    assert (rl>=0 and rl<int(regions.size()));
    return test_impl (hh, dd, bnd, regions, rl);
  }
  
  void property::
  enforce (gh const& hh, dh const& dd,
           level_boundary const& bnd,
           vector<ibset>& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    enforce_impl (hh, dd, bnd, regions, rl);
    if (not test(hh, dd, bnd, regions, rl)) {
      cout << "Property " << typeid(*this).name() << "\n";
      CCTK_WARN (CCTK_WARN_ABORT,
                 "Property does not hold after being enforced");
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  // Ensure that this grid contains the next finer grid
  //////////////////////////////////////////////////////////////////////////////
  
  ibset proper_nesting::
  enlarged_fine_grid (gh const& hh, dh const& dd,
                      level_boundary const& bnd,
                      vector<ibset> const& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (rl+1 < int(regions.size()));
    
    // The minimum amount of space required between the boundaries of
    // this and the next finer grid. We need a certain amount of space
    // on the coarse and a certain amount on the fine grid.
    i2vect const fdistance = dd.ghost_widths.at(rl+1);
    i2vect const cdistance =
      i2vect(min_distance + dd.prolongation_stencil_size(rl));
    
    ibset enlarged;
    
    // Loop over all bboxes that make up the next finer level
    for (ibset::const_iterator ibb = regions.at(rl+1).begin();
         ibb != regions.at(rl+1).end();
         ++ ibb)
    {
      ibbox const& fbb = *ibb;
      
      // Find out which faces are on a boundary
      bvect const lower_is_outer = fbb.lower() <= bnd.level_physical_ilower;
      bvect const upper_is_outer = fbb.upper() >= bnd.level_physical_iupper;
      b2vect const ob (lower_is_outer, upper_is_outer);
      
      ibbox const domext = hh.baseextent(0,rl);
      
      // Enlarge the bbox, first on the fine grid, then transfer it to
      // the coarse grid, then enlarge it again
      ibbox const ebb = fbb.expand (i2vect(not ob) * fdistance);
      ibbox const cbb = ebb.expanded_for (domext);
      ibbox const ecbb = cbb.expand (i2vect(not ob) * cdistance);
      
      // Add it
      enlarged |= ecbb;
    }
    
    return enlarged;
  }
  
  bool proper_nesting::
  test_impl (gh const& hh, const dh& dd,
             level_boundary const& bnd,
             vector<ibset> const& regions, int const rl)
  {
    // This should not be tested because it has to be applied
    // unconditionally and only once
    return true;
  }
  
  void proper_nesting::
  enforce_impl (gh const& hh, dh const& dd,
                level_boundary const& bnd,
                vector<ibset>& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (not ensure_proper_nesting) return;
    if (rl == int(regions.size()) - 1) return;
    
    if (veryverbose) {
      cout << "Refinement level " << rl << ": ensuring proper nesting...\n";
    }
    
    // Enlarge the level
    regions.AT(rl) |= enlarged_fine_grid (hh, dd, bnd, regions, rl);
    
    if (veryverbose) {
      cout << "   New regions are " << regions.at(rl) << "\n";
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  // Add buffer zones (do this only once)
  //////////////////////////////////////////////////////////////////////////////
  
  ibset add_buffers::
  buffered_regions (gh const& hh, dh const& dd,
                    level_boundary const& bnd,
                    vector<ibset> const& regions, int const rl)
  {
    ibset buffered;
    for (ibset::const_iterator
           ibb = regions.at(rl).begin(); ibb != regions.at(rl).end(); ++ ibb)
    {
      ibbox const& bb = *ibb;
      buffered |= bb.expand (dd.buffer_widths.at(rl));
    }
    return buffered;
  }
  
  bool add_buffers::
  test_impl (gh const& hh, dh const& dd,
             level_boundary const& bnd,
             vector<ibset> const& regions, int const rl)
  {
    // This should not be tested because it has to be applied
    // unconditionally and only once
    return true;
  }
  
  void add_buffers::
  enforce_impl (gh const& hh, dh const& dd,
                level_boundary const& bnd,
                vector<ibset>& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (veryverbose) {
      cout << "Refinement level " << rl << ": adding buffer zones...\n";
    }
    
    regions.at(rl) = buffered_regions (hh, dd, bnd, regions, rl);
    
    if (veryverbose) {
      cout << "   New regions are " << regions.at(rl) << "\n";
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  // Combine all regions into a single region, if this is worthwhile
  //////////////////////////////////////////////////////////////////////////////
  
  ibbox combine_regions::
  combined_regions (gh const& hh, dh const& dd,
                    level_boundary const& bnd,
                    vector<ibset> const& regions, int const rl)
  {
    ibbox single;
    for (ibset::const_iterator
           ibb = regions.at(rl).begin(); ibb != regions.at(rl).end(); ++ ibb)
    {
      ibbox const& bb = *ibb;
      single = single.expanded_containing (bb);
    }
    return single;
  }
  
  bool combine_regions::
  test_impl (gh const& hh, dh const& dd,
             level_boundary const& bnd,
             vector<ibset> const& regions, int const rl)
  {
    // This should not be tested because it has to be applied
    // unconditionally and only once
    return true;
  }
  
  void combine_regions::
  enforce_impl (gh const& hh, dh const& dd,
                level_boundary const& bnd,
                vector<ibset>& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (veryverbose) {
      cout << "Refinement level " << rl << ": combining regions...\n";
    }
    
    ibbox const combined = combined_regions (hh, dd, bnd, regions, rl);
    
    CCTK_REAL const regions_size =
      static_cast <CCTK_REAL> (regions.at(rl).size());
    CCTK_REAL const combined_size =
      static_cast <CCTK_REAL> (combined.size());
    
    // Would a single bbox be efficient enough?
    // TODO: Check this also for pairs of regions
    if (min_fraction * combined_size <= regions_size) {
      regions.at(rl) = combined;
    }
    
    if (veryverbose) {
      cout << "   New regions are " << regions.at(rl) << "\n";
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  // Align the boxes with the next coarser grid
  //////////////////////////////////////////////////////////////////////////////
  
  ibset snap_coarse::
  snapped_regions (gh const& hh, dh const& dd, level_boundary const& bnd,
                   vector<ibset> const& regions, int const rl)
  {
    assert (rl>0);
    
    ibbox const& base  = hh.baseextent(0,rl);
    ibbox const& cbase = hh.baseextent(0,rl-1);
    assert (all (cbase.stride() % base.stride() == 0));
    ivect const reffact = cbase.stride() / base.stride();
    i2vect const& buffers = dd.buffer_widths.at(rl);
    
    ibset snapped;
    
    for (ibset::const_iterator
           ibb = regions.at(rl).begin(); ibb != regions.at(rl).end(); ++ ibb)
    {
      ibbox const& bb = *ibb;
      
      // We want to align the current level (without its buffer zones)
      // with the next coarser level. Conceptually, we therefore
      // subtract the buffer zones, align, and add the buffer zones
      // again.
      
      // In detail, we first expand by reffact-1-N points, then expand
      // to the next coarser grid, then expand back to the current
      // grid, and expand by N points again. This sequence is correct
      // for both vertex and cell centred grids, and N is determined
      // by the number of buffer zones.
      
      // N is the number of buffer zones modulo the refinement factor.
      // We cannot shrink the domain (since we cannot shrink
      // bboxsets). For alignment, only operations modulo the
      // refinement factor are relevant.
      
      snapped |= bb.
        expand(reffact-1 - buffers % reffact).
        contracted_for(cbase).
        expanded_for(base).
        expand(buffers % reffact);
    }
    
    return snapped;
  }
  
  bool snap_coarse::
  test_impl (gh const& hh, dh const& dd,
             level_boundary const& bnd,
             vector<ibset> const& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (not snap_to_coarse) return true;
    
    ibset const snapped = snapped_regions (hh, dd, bnd, regions, rl);
    return regions.AT(rl) == snapped;
  }
  
  void snap_coarse::
  enforce_impl (gh const& hh, dh const& dd,
                level_boundary const& bnd,
                vector<ibset>& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (snap_to_coarse);
    
    if (veryverbose) {
      cout << "Refinement level " << rl << ": aligning regions with next coarser grid...\n";
    }
    
    regions.AT(rl) = snapped_regions (hh, dd, bnd, regions, rl);
    
    if (veryverbose) {
      cout << "   New regions are " << regions.at(rl) << "\n";
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  // Make the boxes rotating-90 symmetric
  //////////////////////////////////////////////////////////////////////////////
  
  ibset rotsym90::
  symmetrised_regions (gh const& hh, dh const& dd, level_boundary const& bnd,
                       vector<ibset> const& regions, int const rl)
  {
    ibset symmetrised;
    for (ibset::const_iterator
           ibb = regions.at(rl).begin(); ibb != regions.at(rl).end(); ++ ibb)
    {
      ibbox const& bb = *ibb;
      
      bvect const lower_is_outside_lower =
        bb.lower() - bnd.min_bnd_dist_away[0] * bb.stride() <=
        bnd.level_physical_ilower;
      
      // Treat both x and y directions
      for (int dir=0; dir<=1; ++dir) {
        if (lower_is_outside_lower[dir]) {
          ivect const ilo = bb.lower();
          ivect const iup = bb.upper();
          ivect const istr = bb.stride();
          
          // Origin
          rvect const axis (bnd.physical_lower[0],
                            bnd.physical_lower[1],
                            bnd.physical_lower[2]); // z component is unused
          ivect const iaxis0 = rpos2ipos (axis, bnd.origin, bnd.scale, hh, rl);
          assert (all (iaxis0 % istr == 0));
          ivect const iaxis1 = rpos2ipos1 (axis, bnd.origin, bnd.scale, hh, rl);
          assert (all (iaxis1 % istr == 0));
          ivect const offset = iaxis1 - iaxis0;
          assert (all (offset % istr == 0));
          assert (all (offset >= 0 and offset < 2*istr));
          assert (all ((iaxis0 + iaxis1 - offset) % (2*istr) == 0));
          ivect const iaxis = (iaxis0 + iaxis1 - offset) / 2;
          // negated (reflected) domain boundaries
          ivect const neg_ilo = (2*iaxis+offset) - ilo;
          ivect const neg_iup = (2*iaxis+offset) - iup;
          // offset to add when permuting directions
          ivect const permute01 (-iaxis[0]+iaxis[1], -iaxis[1]+iaxis[0], 0);
          
          // Rotate 90 degrees about z axis
          ivect new_ilo, new_iup;
          if (dir==0) {
            // rotate clockwise
            new_ilo = ivect (ilo[1], neg_iup[0], ilo[2]) + permute01;
            new_iup = ivect (iup[1], neg_ilo[0], iup[2]) + permute01;
          }
          if (dir==1) {
            // rotate counterclockwise
            new_ilo = ivect (neg_iup[1], ilo[0],  ilo[2]) + permute01;
            new_iup = ivect (neg_ilo[1],  iup[0],  iup[2]) + permute01;
          }
          ivect const new_istr (istr);
          
          ibbox const new_bb (new_ilo, new_iup, new_istr);
          // Will be clipped later
          // assert (new_bb.is_contained_in (baseextent));
          
          symmetrised |= new_bb;
        }
      }
    }
    
    return symmetrised;
  }
  
  bool rotsym90::
  test_impl (gh const& hh, dh const& dd,
             level_boundary const& bnd,
             vector<ibset> const& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (not symmetry_rotating90) return true;
    
    ibset const symmetrised = symmetrised_regions (hh, dd, bnd, regions, rl);
    return regions.AT(rl) == symmetrised;
  }
  
  void rotsym90::
  enforce_impl (gh const& hh, dh const& dd,
                level_boundary const& bnd,
                vector<ibset>& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (symmetry_rotating90);
    
    if (veryverbose) {
      cout << "Refinement level " << rl << ": making regions rotating-90 symmetric...\n";
    }
    
    regions.AT(rl) = symmetrised_regions (hh, dd, bnd, regions, rl);
    
    if (veryverbose) {
      cout << "   New regions are " << regions.at(rl) << "\n";
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  // Make the boxes rotating-180 symmetric
  //////////////////////////////////////////////////////////////////////////////
  
  ibset rotsym180::
  symmetrised_regions (gh const& hh, dh const& dd,
                       level_boundary const& bnd,
                       vector<ibset> const& regions, int const rl)
  {
    ibbox const& baseextent = hh.baseextent(0,rl);
    
    ibset symmetrised;
    for (ibset::const_iterator
           ibb = regions.at(rl).begin(); ibb != regions.at(rl).end(); ++ ibb)
    {
      ibbox const& bb = *ibb;
      
      bvect const lower_is_outside_lower =
        bb.lower() - bnd.min_bnd_dist_away[0] * bb.stride() <=
        bnd.level_physical_ilower;
      
      // Treat x direction
      int const dir = 0;
      if (lower_is_outside_lower[dir]) {
        ivect const ilo = bb.lower();
        ivect const iup = bb.upper();
        ivect const istr = bb.stride();
        assert (istr[0] == istr[1]);
        
        // Origin
        assert (hh.refcent == vertex_centered or all (istr % 2 == 0));
        rvect const axis (bnd.physical_lower[0],
                          (bnd.physical_lower[1] + bnd.physical_upper[1]) / 2,
                          bnd.physical_lower[2]); // z component is unused
        ivect const iaxis0 = rpos2ipos (axis, bnd.origin, bnd.scale, hh, rl);
        assert (all ((iaxis0 - baseextent.lower()) % istr == 0));
        ivect const iaxis1 = rpos2ipos1 (axis, bnd.origin, bnd.scale, hh, rl);
        assert (all ((iaxis1 - baseextent.lower()) % istr == 0));
        ivect const offset = iaxis1 - iaxis0;
        assert (all (offset % istr == 0));
        if (hh.refcent == vertex_centered) {
          assert (all (offset >= 0 and offset < 2*istr));
          assert (all ((iaxis0 + iaxis1 - offset) % (2*istr) == 0));
        } else {
          // The offset may be negative because both boundaries are
          // shifted inwards by 1/2 grid spacing, and therefore iaxis0
          // < iaxis1 + istr
          assert (all (offset >= -istr and offset < istr));
          assert (all ((iaxis0 + iaxis1 - offset) % (2*istr) == istr));
          assert (all (istr % 2 == 0));
        }
        ivect const iaxis = (iaxis0 + iaxis1 - offset) / 2;
        ivect const neg_ilo = (2*iaxis+offset) - ilo;
        ivect const neg_iup = (2*iaxis+offset) - iup;
        
        // Rotate 180 degrees about z axis
        ivect const new_ilo (neg_iup[0], neg_iup[1], ilo[2]);
        ivect const new_iup (neg_ilo[0], neg_ilo[1], iup[2]);
        ivect const new_istr (istr);
        
        ibbox const new_bb (new_ilo, new_iup, new_istr);
        // Will be clipped later
        // assert (new_bb.is_contained_in (baseextent));
        
        symmetrised |= new_bb;
      }
    }
    
    return symmetrised;
  }
  
  bool rotsym180::
  test_impl (gh const& hh, dh const& dd,
             level_boundary const& bnd,
             vector<ibset> const& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (not symmetry_rotating180) return true;
    
    ibset const symmetrised = symmetrised_regions (hh, dd, bnd, regions, rl);
    return regions.AT(rl) == symmetrised;
  }
  
  void rotsym180::
  enforce_impl (gh const& hh, dh const& dd,
                level_boundary const& bnd,
                vector<ibset>& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (symmetry_rotating180);
    
    if (veryverbose) {
      cout << "Refinement level " << rl << ": making regions rotating-180 symmetric...\n";
    }
    
    regions.AT(rl) = symmetrised_regions (hh, dd, bnd, regions, rl);
    
    if (veryverbose) {
      cout << "   New regions are " << regions.at(rl) << "\n";
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  // Clip at the outer boundary
  //////////////////////////////////////////////////////////////////////////////
  
  ibset boundary_clip::
  clipped_regions (gh const& hh, dh const& dd,
                   level_boundary const& bnd,
                   vector<ibset> const& regions, int const rl)
  {
    ibbox const& baseextent = hh.baseextent(0,rl);
    
    ibset clipped;
    for (ibset::const_iterator
           ibb = regions.at(rl).begin(); ibb != regions.at(rl).end(); ++ ibb)
    {
      ibbox const& bb = *ibb;
      
      // Clip boxes that extend outside the boundary. Enlarge boxes
      // that are inside but too close to the outer boundary.
      bvect const lower_is_outside_lower =
        bb.lower() - bnd.min_bnd_dist_away[0] * bb.stride() <=
        bnd.level_physical_ilower;
      // Remove bboxes that are completely outside.
      bvect const upper_is_outside_lower =
        bb.upper() < bnd.level_physical_ilower;
      // Enlarge bboxes that extend not far enough inwards.
      bvect const upper_is_almost_outside_lower =
        bb.upper() <
        bnd.level_physical_ilower + bnd.min_bnd_dist_incl[0] * bb.stride();
      
      // Ditto for the upper boundary.
      bvect const upper_is_outside_upper =
        bb.upper() + bnd.min_bnd_dist_away[1] * bb.stride() >=
        bnd.level_physical_iupper;
      bvect const lower_is_outside_upper =
        bb.lower() > bnd.level_physical_iupper;
      bvect const lower_is_almost_outside_upper =
        bb.lower() >
        bnd.level_physical_iupper - bnd.min_bnd_dist_incl[1] * bb.stride();
      
      assert (not any (lower_is_almost_outside_upper and
                       lower_is_outside_lower));
      assert (not any (upper_is_almost_outside_lower and
                       upper_is_outside_upper));
      
      if (any (upper_is_outside_lower or lower_is_outside_upper)) {
        // The box is completely outside. Ignore it.
        continue;
      }
          
      if (any ((lower_is_outside_lower and
                bnd.boundary_staggering_mismatch[0]) or
               (upper_is_outside_upper and
                bnd.boundary_staggering_mismatch[1])))
      {
        ostringstream msg;
        msg << "Level " << rl << " of the refinement hierarchy has inconsistent bountary staggering."
            << "  The refined region extends up to the boundary, but the staggering of the boundary is different from the staggering of the mesh refinement."
            << "  lower_is_outside_lower=" << lower_is_outside_lower
            << "  upper_is_outside_upper=" << upper_is_outside_upper
            << "  boundary_staggering_mismatch=" << bnd.boundary_staggering_mismatch
            << "  level_physical_ilower=" << bnd.level_physical_ilower
            << "  level_physical_iupper=" << bnd.level_physical_iupper
            << "  baseextent=" << baseextent;
        CCTK_WARN (CCTK_WARN_ABORT, msg.str().c_str());
      }
      
      ibbox const clipped_bb
        (either (lower_is_outside_lower,
                 bnd.level_exterior_ilower,
                 either (lower_is_almost_outside_upper,
                         (bnd.level_physical_iupper -
                          bnd.min_bnd_dist_incl[1] * bb.stride()),
                         bb.lower())),
         either (upper_is_outside_upper,
                 bnd.level_exterior_iupper,
                 either (upper_is_almost_outside_lower,
                         (bnd.level_physical_ilower +
                          bnd.min_bnd_dist_incl[0] * bb.stride()),
                         bb.upper())),
         bb.stride());
      if (not clipped_bb.is_contained_in (baseextent)) {
        ostringstream msg;
        msg << "Level " << rl << " of the refinement hierarchy is not contained in the simulation domain."
            << "  (There may be too many ghost or buffer zones.)"
            << "  One bbox is " << clipped_bb << "."
            << "  lower_is_outside_lower=" << lower_is_outside_lower
            << "  upper_is_outside_upper=" << upper_is_outside_upper
            << "  lower_is_almost_outside_upper=" << lower_is_almost_outside_upper
            << "  upper_is_almost_outside_lower=" << upper_is_almost_outside_lower
            << "  level_exterior_ilower=" << bnd.level_exterior_ilower
            << "  level_exterior_iupper=" << bnd.level_exterior_iupper
            << "  level_physical_ilower=" << bnd.level_physical_ilower
            << "  level_physical_iupper=" << bnd.level_physical_iupper
            << "  baseextent=" << baseextent;
        CCTK_WARN (CCTK_WARN_ABORT, msg.str().c_str());
      }
      assert (clipped_bb.is_contained_in (baseextent));
      
      clipped |= clipped_bb;
    }
    
    return clipped;
  }
  
  bool boundary_clip::
  test_impl (gh const& hh, dh const& dd,
             level_boundary const& bnd,
             vector<ibset> const& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    ibset const clipped = clipped_regions (hh, dd, bnd, regions, rl);
    return regions.AT(rl) == clipped;
  }

  void boundary_clip::
  enforce_impl (gh const& hh, dh const& dd,
                level_boundary const& bnd,
                vector<ibset>& regions, int const rl)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (veryverbose) {
      cout << "Refinement level " << rl << ": clipping at outer boundary...\n";
    }
    
    regions.AT(rl) = clipped_regions (hh, dd, bnd, regions, rl);
    
    if (veryverbose) {
      cout << "   New regions are " << regions.at(rl) << "\n";
    }
  }
  
  
  
  //////////////////////////////////////////////////////////////////////////////
  // Ensure that this grid is contained in the domain
  //////////////////////////////////////////////////////////////////////////////
  
  bool in_domain::
  test_impl (gh const& hh, dh const& dd,
             level_boundary const& bnd,
             vector<ibset> const& regions, int const rl)
  {
    return regions.at(rl) <= hh.baseextent(0,rl);
  }
  
  void in_domain::
  enforce_impl (gh const& hh, dh const& dd,
                level_boundary const& bnd,
                vector<ibset>& regions, int const rl)
  {
    // There is nothing we can do here, since we can't enlarge the
    // domain
    CCTK_WARN (CCTK_WARN_ABORT, "internal error");
  }
  
  
  
} // namespace CarpetRegrid2