aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetReduce/src/mask_carpet.cc
blob: d0861320f269646646ea01949228d1fcc750bf18 (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
#include <cassert>
#include <sstream>

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

#include "carpet.hh"

#include "mask_carpet.hh"

extern "C" {
  static const char* rcsid = "$Header: /home/eschnett/C/carpet/Carpet/Carpet/CarpetReduce/src/mask_carpet.cc,v 1.5 2004/08/05 11:15:49 schnetter Exp $";
  CCTK_FILEVERSION(Carpet_CarpetMask_Mask_cc);
}



namespace CarpetMask {
  
  using namespace std;
  using namespace Carpet;
  
  
  
  /**
   * Reduce the weight on the current and the next coarser level to
   * make things consistent.  Set the weight to 0 inside the
   * restriction region of the next coarser level, maybe to 1/2 near
   * the boundary of that region, and also to 1/2 near the
   * prolongation boundary of this level.
   */
  
  void
  CarpetMaskSetup (CCTK_ARGUMENTS)
  {
    DECLARE_CCTK_PARAMETERS;
    
    assert (reffact == 2);
    
    if (! is_singlemap_mode()) {
      CCTK_WARN (0, "This routine may only be called in singlemap mode");
    }
    
    if (reflevel > 0) {
      
      ivect const izero = ivect(0);
      ivect const ione  = ivect(1);
      
      gh<dim> const & hh = *vhh.at(Carpet::map);
      dh<dim> const & dd = *vdd.at(Carpet::map);
      
      ibbox const & base = hh.bases.at(reflevel).at(mglevel);
      
      
      
      // Calculate the union of all refined regions
      ibset refined;
      for (int c=0; c<hh.components(reflevel); ++c) {
        refined |= hh.extents.at(reflevel).at(c).at(mglevel);
      }
      refined.normalize();
      
      // Calculate the union of all coarse regions
      ibset parent;
      for (int c=0; c<hh.components(reflevel-1); ++c) {
        parent |= hh.extents.at(reflevel-1).at(c).at(mglevel).expanded_for(base);
      }
      parent.normalize();
      
      // Subtract the refined region
      ibset notrefined = parent - refined;
      notrefined.normalize();
      
      // Enlarge this set
      ibset enlarged[dim];
      for (int d=0; d<dim; ++d) {
        for (ibset::const_iterator bi = notrefined.begin();
             bi != notrefined.end();
             ++bi)
        {
          enlarged[d] |= (*bi).expand(ivect::dir(d), ivect::dir(d));
        }
        enlarged[d].normalize();
      }
      
      // Intersect with the original union
      ibset boundaries[dim];
      for (int d=0; d<dim; ++d) {
        boundaries[d] = refined & enlarged[d];
        boundaries[d].normalize();
      }
      
      // Subtract the boundaries from the refined region
      for (int d=0; d<dim; ++d) {
        refined -= boundaries[d];
      }
      refined.normalize();
      
      
      
      // Set prolongation boundaries of this level
      {
        BEGIN_LOCAL_COMPONENT_LOOP (cctkGH, CCTK_GF) {
          
          DECLARE_CCTK_ARGUMENTS;
          
          ibbox const & ext
            = dd.boxes.at(reflevel).at(component).at(mglevel).exterior;
          
          for (int d=0; d<dim; ++d) {
            for (ibset::const_iterator bi = boundaries[d].begin();
                 bi != boundaries[d].end();
                 ++bi)
            {
              
              ibbox const & box = (*bi) & ext;
              if (! box.empty()) {
                
                assert (all ((box.lower() - ext.lower()               ) >= 0));
                assert (all ((box.upper() - ext.lower() + ext.stride()) >= 0));
                assert (all ((box.lower() - ext.lower()               ) % ext.stride() == 0));
                assert (all ((box.upper() - ext.lower() + ext.stride()) % ext.stride() == 0));
                ivect const imin = (box.lower() - ext.lower()               ) / ext.stride();
                ivect const imax = (box.upper() - ext.lower() + ext.stride()) / ext.stride();
                assert (all (izero <= imin));
                assert (box.empty() || all (imin <= imax));
                assert (all (imax <= ivect::ref(cctk_lsh)));
                
                if (verbose) {
                  ostringstream buf;
                  buf << "Setting prolongation boundary on level " << reflevel << " direction " << d << " to weight 1/2: " << imin << ":" << imax-ione;
                  CCTK_INFO (buf.str().c_str());
                }
                
                // Set weight on the boundary to 1/2
                assert (dim == 3);
                for (int k=imin[2]; k<imax[2]; ++k) {
                  for (int j=imin[1]; j<imax[1]; ++j) {
                    for (int i=imin[0]; i<imax[0]; ++i) {
                      int const ind = CCTK_GFINDEX3D (cctkGH, i, j, k);
                      weight[ind] *= 0.5;
                    }
                  }
                }
                
              } // if box not empty
              
            } // for box
          } // for d
          
        } END_LOCAL_COMPONENT_LOOP;
      }
      
      
      
      // Set restriction region on next coarser level
      {
        int const oldreflevel = reflevel;
        int const oldmap = Carpet::map;
        leave_singlemap_mode (cctkGH);
        leave_level_mode (cctkGH);
        enter_level_mode (cctkGH, oldreflevel-1);
        enter_singlemap_mode (cctkGH, oldmap);
        
        BEGIN_LOCAL_COMPONENT_LOOP (cctkGH, CCTK_GF) {
          
          DECLARE_CCTK_ARGUMENTS;
          
          ibbox const & ext
            = dd.boxes.at(reflevel).at(component).at(mglevel).exterior;
          
          for (ibset::const_iterator bi = refined.begin();
               bi != refined.end();
               ++bi)
          {
            
            ibbox const & box = (*bi).contracted_for(ext) & ext;
            if (! box.empty()) {
              
              assert (all ((box.lower() - ext.lower()               ) >= 0));
              assert (all ((box.upper() - ext.lower() + ext.stride()) >= 0));
              assert (all ((box.lower() - ext.lower()               ) % ext.stride() == 0));
              assert (all ((box.upper() - ext.lower() + ext.stride()) % ext.stride() == 0));
              ivect const imin = (box.lower() - ext.lower()               ) / ext.stride();
              ivect const imax = (box.upper() - ext.lower() + ext.stride()) / ext.stride();
              assert (all (izero <= imin));
              assert (box.empty() || all (imin <= imax));
              assert (all (imax <= ivect::ref(cctk_lsh)));
              
              if (verbose) {
                ostringstream buf;
                buf << "Setting restricted region on level " << reflevel << " to weight 0: " << imin << ":" << imax-ione;
                CCTK_INFO (buf.str().c_str());
              }
              
              // Set weight in the restricted region to 0
              assert (dim == 3);
              for (int k=imin[2]; k<imax[2]; ++k) {
                for (int j=imin[1]; j<imax[1]; ++j) {
                  for (int i=imin[0]; i<imax[0]; ++i) {
                    int const ind = CCTK_GFINDEX3D (cctkGH, i, j, k);
                    weight[ind] = 0;
                  }
                }
              }
              
            } // if box not empty
              
          } // for box
          
          assert (dim == 3);
          vector<int> mask (cctk_lsh[0] * cctk_lsh[1] * cctk_lsh[2]);
          
          assert (dim == 3);
          for (int k=0; k<cctk_lsh[2]; ++k) {
            for (int j=0; j<cctk_lsh[1]; ++j) {
              for (int i=0; i<cctk_lsh[0]; ++i) {
                int const ind = CCTK_GFINDEX3D (cctkGH, i, j, k);
                mask[ind] = 0;
              }
            }
          }
          
          for (int d=0; d<dim; ++d) {
            for (ibset::const_iterator bi = boundaries[d].begin();
                 bi != boundaries[d].end();
                 ++bi)
            {
              
              ibbox const & box = (*bi).contracted_for(ext) & ext;
              if (! box.empty()) {
                
                assert (all ((box.lower() - ext.lower()               ) >= 0));
                assert (all ((box.upper() - ext.lower() + ext.stride()) >= 0));
                assert (all ((box.lower() - ext.lower()               ) % ext.stride() == 0));
                assert (all ((box.upper() - ext.lower() + ext.stride()) % ext.stride() == 0));
                ivect const imin = (box.lower() - ext.lower()               ) / ext.stride();
                ivect const imax = (box.upper() - ext.lower() + ext.stride()) / ext.stride();
                assert (all (izero <= imin));
                assert (box.empty() || all (imin <= imax));
                assert (all (imax <= ivect::ref(cctk_lsh)));
                
                if (verbose) {
                  ostringstream buf;
                  buf << "Setting restriction boundary on level " << reflevel << " direction " << d << " to weight 1/2: " << imin << ":" << imax-ione;
                  CCTK_INFO (buf.str().c_str());
                }
                
                // Set weight on the boundary to 1/2
                assert (dim == 3);
                for (int k=imin[2]; k<imax[2]; ++k) {
                  for (int j=imin[1]; j<imax[1]; ++j) {
                    for (int i=imin[0]; i<imax[0]; ++i) {
                      int const ind = CCTK_GFINDEX3D (cctkGH, i, j, k);
                      if (mask[ind] == 0) {
                        mask[ind] = 1;
                      }
                      mask[ind] *= 2;
                    }
                  }
                }
                
              } // if box not empty
              
            } // for box
          } // for d
          
          assert (dim == 3);
          for (int k=0; k<cctk_lsh[2]; ++k) {
            for (int j=0; j<cctk_lsh[1]; ++j) {
              for (int i=0; i<cctk_lsh[0]; ++i) {
                int const ind = CCTK_GFINDEX3D (cctkGH, i, j, k);
                if (mask[ind] > 0) {
                  weight[ind] *= 1.0 - 1.0 / mask[ind];
                }
              }
            }
          }
          
        } END_LOCAL_COMPONENT_LOOP;
        
        leave_singlemap_mode (cctkGH);
        leave_level_mode (cctkGH);
        enter_level_mode (cctkGH, oldreflevel);
        enter_singlemap_mode (cctkGH, oldmap);
      }
      
    } // if reflevel>0
  }
  
} // namespace CarpetMask