aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetMask/src/mask_surface.cc
blob: e1a02079fac5a77785e1d2864e1550abc75f9c4f (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
#include <algorithm>
#include <cassert>
#include <cmath>
#include <vector>

#include <cctk.h>
#include <cctk_Arguments.h>
#include <cctk_Parameters.h>

#include "mask_surface.hh"



namespace CarpetMask {
  
  using namespace std;
  
  
  
  // Number of excluded regions which are defined in param.ccl
  int const num_excluded = 10;
  
  
  
  /**
   * Set the weight in the excluded regions to zero.
   */
  
  void
  CarpetSurfaceSetup (CCTK_ARGUMENTS)
  {
    DECLARE_CCTK_ARGUMENTS;
    DECLARE_CCTK_PARAMETERS;
    
    // Some state verbose output
    static vector <int> last_output;
    if (last_output.empty()) {
      last_output.resize (num_excluded, -1);
    }
    
    CCTK_REAL * const weight =
      static_cast <CCTK_REAL *>
      (CCTK_VarDataPtr (cctkGH, 0, "CarpetReduce::weight"));
    
    if (not weight) {
      CCTK_WARN (CCTK_WARN_ABORT,
                 "CarpetReduce is not active, or CarpetReduce::mask does not have storage");
    }
    
    for (int n = 0; n < num_excluded; ++ n) {
      
      int const sn = excluded_surface[n];
      if (sn >= 0) {
        
        assert (sn < nsurfaces);
        
        if (sf_active[sn]) {
          
          CCTK_REAL const shrink_factor = excluded_surface_factor[n];
          
          CCTK_REAL const x0 = sf_origin_x[sn];
          CCTK_REAL const y0 = sf_origin_y[sn];
          CCTK_REAL const z0 = sf_origin_z[sn];
          
          CCTK_REAL const theta0 = sf_origin_theta[sn];
          CCTK_REAL const phi0   = sf_origin_phi  [sn];
          CCTK_REAL const dtheta = sf_delta_theta[sn];
          CCTK_REAL const dphi   = sf_delta_phi  [sn];
          
          int const ntheta = sf_ntheta[sn];
          int const nphi   = sf_nphi  [sn];
          
          if (verbose) {
            if (cctk_iteration > last_output.at(n)) {
              last_output.at(n) = cctk_iteration;
              CCTK_VInfo (CCTK_THORNSTRING,
                          "Setting mask from surface %d at [%g,%g,%g], average radius %g",
                          sn,
                          double(x0), double(y0), double(z0),
                          double(sf_mean_radius[sn]));
            }
          }
          
          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);
                
                CCTK_REAL const dx = x[ind] - x0;
                CCTK_REAL const dy = y[ind] - y0;
                CCTK_REAL const dz = z[ind] - z0;
                
                CCTK_REAL const rho =
                  sqrt (pow (dx, 2) + pow (dy, 2) + pow (dz, 2));
                if (rho < 1.0e-12) {
                  // Always excise the surface origin
                  weight[ind] = 0.0;
                } else {
                  CCTK_REAL theta =
                    acos (min (CCTK_REAL (+1.0),
                               max (CCTK_REAL (-1.0), dz / rho)));
                  if (symmetric_z[sn]) {
                    if (theta > M_PI/2.0) {
                      theta = M_PI - theta;
                    }
                  }
                  
                  assert (not isnan (theta));
                  assert (theta >= 0);
                  assert (theta <= M_PI);
                  CCTK_REAL phi =
                    fmod (atan2 (dy, dx) + CCTK_REAL (2 * M_PI),
                          CCTK_REAL (2 * M_PI));
                  if (symmetric_x[sn] or symmetric_y[sn]) {
                    if (symmetric_x[sn] and symmetric_y[sn]) {
                      if (phi > M_PI / 2.0) {
                        phi = M_PI - phi;
                      }
                    } else {
                      if (phi > M_PI) {
                        phi = 2 * M_PI - phi;
                      }
                    }
                  }
                  assert (not isnan (phi));
                  assert (phi >= 0);
                  assert (phi < 2 * M_PI);
                  int const a = floor ((theta - theta0) / dtheta + 0.5);
                  assert (a >= 0);
                  assert (a < ntheta);
                  int const b = floor ((phi   - phi0  ) / dphi   + 0.5);
                  assert (b >= 0);
                  assert (b < nphi);
                  
                  assert (a >= 0 and a < ntheta);
                  assert (b >= 0 and b < nphi  );
                  
                  CCTK_REAL const dr =
                    sf_radius[a + maxntheta * (b + maxnphi * sn)];
                  if (rho <= dr * shrink_factor) {
                    weight[ind] = 0.0;
                  }
                }
                
              }
            }
          }
          
        } else {
          
          if (verbose) {
            if (cctk_iteration > last_output.at(n)) {
              last_output.at(n) = cctk_iteration;
              CCTK_VInfo (CCTK_THORNSTRING, "Surface %d is not active", sn);
            }
          }
          
        } // if surface is not active
        
      } // if sn >= 0
    }   // for n
    
  }
  
} // namespace CarpetMask