aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetRegrid/src/manualcoordinates.cc
blob: e6d2d2ae06b6372129c3b267fdb24e930c16b276 (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
#include <cassert>
#include <vector>

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

#include "gh.hh"

#include "carpet.hh"
#include "regrid.hh"



namespace CarpetRegrid {
  
  using namespace std;
  using namespace Carpet;
  
  
  
  int ManualCoordinates (cGH const * const cctkGH,
                         gh const & hh,
                         gh::mregs & regsss)
  {
    DECLARE_CCTK_PARAMETERS;
    
    if (refinement_levels > 4) {
      CCTK_WARN (0, "Cannot currently specify manual refinement regions for more than 4 refinement levels");
    }
    assert (refinement_levels >= 1 and refinement_levels <= 4);
    
    // do nothing if the levels already exist
    if (reflevel == refinement_levels) return 0;
    
    assert (regsss.size() >= 1);
    vector<vector<region_t> > regss = regsss.at(0);
    
    regss.resize (refinement_levels);
    
    vector<rvect> lower(3), upper(3);
    lower.at(0) = rvect (l1xmin, l1ymin, l1zmin);
    upper.at(0) = rvect (l1xmax, l1ymax, l1zmax);
    lower.at(1) = rvect (l2xmin, l2ymin, l2zmin);
    upper.at(1) = rvect (l2xmax, l2ymax, l2zmax);
    lower.at(2) = rvect (l3xmin, l3ymin, l3zmin);
    upper.at(2) = rvect (l3xmax, l3ymax, l3zmax);
    
    assert (! smart_outer_boundaries);
    
    for (size_t rl=1; rl<regss.size(); ++rl) {
      
      b2vect const ob (false);
      region_t reg;
      reg.map = Carpet::map;
      reg.outer_boundaries = ob;
      
      vector<region_t> regs;
      ManualCoordinates_OneLevel
        (cctkGH, hh, rl, refinement_levels,
         lower.at(rl-1), upper.at(rl-1), reg, regs);
      
      // make multiprocessor aware
      SplitRegions (cctkGH, regs);
      
      regss.at(rl) = regs;
      
    } // for rl
    
    // make multigrid aware
    MakeMultigridBoxes (cctkGH, Carpet::map, regss, regsss);
    
    return 1;
  }
  
  
  
  void ManualCoordinates_OneLevel (const cGH * const cctkGH,
                                   const gh & hh,
                                   const int rl,
                                   const int numrl,
                                   const rvect lower,
                                   const rvect upper,
                                   const region_t & reg,
                                   vector<region_t> & regs)
  {
    if (rl >= numrl) return;
    
    jvect const ilower = pos2int (cctkGH, hh, lower, rl);
    jvect const iupper = pos2int (cctkGH, hh, upper, rl);
    
    ManualGridpoints_OneLevel
      (cctkGH, hh, rl, numrl, ilower, iupper, reg, regs);
  }
  
  
  
  ivect delta2int (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.baseextents.at(0).at(0).upper() -
                               hh.baseextents.at(0).at(0).lower());
    
    const rvect scale  = rvect(global_extent) / (global_upper - global_lower);
    const ivect levfac = hh.reffacts.at(rl);
    assert (all (hh.baseextents.at(0).at(0).stride() % levfac == 0));
    const ivect istride = hh.baseextents.at(0).at(0).stride() / levfac;
    
    const ivect ipos
      = ivect(floor(rpos * scale / rvect(istride) + (CCTK_REAL) 0.5)) * istride;
    
    const rvect apos = rpos * scale;
    assert (all(abs(apos - rvect(ipos)) < rvect(istride) * (CCTK_REAL) 0.01));
    
    return ipos;
  }
  
  
  
  ivect pos2int (const cGH* const cctkGH, const gh& hh,
                 const rvect & rpos, const int rl)
  {
    rvect global_lower, global_upper;
#if 0
    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;
      }
    }
#endif
    assert (Carpet::map >= 0);
    global_lower = domainspecs.at(Carpet::map).exterior_min;
    global_upper = domainspecs.at(Carpet::map).exterior_max;
    const ivect global_extent (hh.baseextents.at(0).at(0).upper() -
                               hh.baseextents.at(0).at(0).lower());
    
    const rvect scale  = rvect(global_extent) / (global_upper - global_lower);
    const ivect levfac = hh.reffacts.at(rl);
    assert (all (hh.baseextents.at(0).at(0).stride() % levfac == 0));
    const ivect istride = hh.baseextents.at(0).at(0).stride() / levfac;
    
    const ivect ipos
      = (ivect(floor((rpos - global_lower) * scale / rvect(istride)
                     + (CCTK_REAL) 0.5))
         * istride);
    
    return ipos;
  }
  
} // namespace CarpetRegrid