aboutsummaryrefslogtreecommitdiff
path: root/Carpet/CarpetInterp2/src/interp2.cc
blob: 109819f08390e0c5799c555b4ded198da1be0c55 (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
#include <cctk.h>

#include "fasterp.hh"



namespace CarpetInterp2 {
  
  extern "C"
  CCTK_INT
  CarpetInterp2_InterpGridArrays (CCTK_POINTER_TO_CONST const cctkGH_,
                                  CCTK_INT const N_dims,
                                  CCTK_INT const order,
                                  CCTK_INT const N_interp_points,
                                  CCTK_POINTER_TO_CONST const interp_coords_,
                                  CCTK_INT const N_input_arrays,
                                  CCTK_INT const * const input_array_indices,
                                  CCTK_INT const N_output_arrays,
                                  CCTK_POINTER const output_arrays_)
  {
    // Check input values and convert types
    cGH const * const cctkGH = static_cast<cGH const *> (cctkGH_);
    assert (cctkGH);
    
    assert (N_dims == dim);
    
    assert (N_interp_points >= 0);
    
    CCTK_REAL const * const * const interp_coords =
      static_cast<CCTK_REAL const * const *> (interp_coords_);
    assert (interp_coords);
    for (int d=0; d<dim; ++d) {
      assert (N_interp_points==0 or interp_coords[d]);
    }
    
    assert (N_input_arrays >= 0);
    
    assert (input_array_indices);
    for (int n=0; n<N_input_arrays; ++n) {
      assert (input_array_indices[n]>=0 and
              input_array_indices[n]<CCTK_NumVars());
    }
    
    assert (N_output_arrays >= 0);
    assert (N_output_arrays == N_input_arrays);
    
    CCTK_REAL * const * const output_arrays =
      static_cast<CCTK_REAL * const *> (output_arrays_);
    assert (output_arrays);
    for (int n=0; n<N_output_arrays; ++n) {
      assert (N_output_arrays==0 or output_arrays[n]);
    }
    
    // Set up interpolation
    fasterp_glocs_t locations (N_interp_points);
    for (int d=0; d<dim; ++d) {
      for (int i=0; i<N_interp_points; ++i) {
        locations.coords[d].AT(i) = interp_coords[d][i];
      }
    }
    fasterp_setup_t const setup (cctkGH, locations, order);
    
    // Interpolate
    vector<int> varinds (N_input_arrays);
    for (int n=0; n<N_input_arrays; ++n) {
      varinds.AT(n) = input_array_indices[n];
    }
    vector<CCTK_REAL *> values (N_input_arrays);
    for (int n=0; n<N_input_arrays; ++n) {
      values.AT(n) = output_arrays[n];
    }
    setup.interpolate (cctkGH, varinds, values);
    
    // Done
    return 0;
  }
  
  extern "C"
  CCTK_INT
  CarpetInterp2_MultiPatchInterpGridArrays (CCTK_POINTER_TO_CONST const cctkGH_,
                                            CCTK_INT const N_dims,
                                            CCTK_INT const order,
                                            CCTK_INT const N_interp_points,
                                            CCTK_INT const * const interp_maps,
                                            CCTK_POINTER_TO_CONST const interp_coords_,
                                            CCTK_INT const N_input_arrays,
                                            CCTK_INT const * const input_array_indices,
                                            CCTK_INT const N_output_arrays,
                                            CCTK_POINTER const output_arrays_)
  {
    // Check input values and convert types
    cGH const * const cctkGH = static_cast<cGH const *> (cctkGH_);
    assert (cctkGH);
    
    assert (N_dims == dim);
    
    assert (N_interp_points >= 0);
    
    assert (interp_maps);
    CCTK_REAL const * const * const interp_coords =
      static_cast<CCTK_REAL const * const *> (interp_coords_);
    assert (interp_coords);
    for (int d=0; d<dim; ++d) {
      assert (N_interp_points==0 or interp_coords[d]);
    }
    
    assert (N_input_arrays >= 0);
    
    assert (input_array_indices);
    for (int n=0; n<N_input_arrays; ++n) {
      assert (input_array_indices[n]>=0 and
              input_array_indices[n]<CCTK_NumVars());
    }
    
    assert (N_output_arrays >= 0);
    assert (N_output_arrays == N_input_arrays);
    
    CCTK_REAL * const * const output_arrays =
      static_cast<CCTK_REAL * const *> (output_arrays_);
    assert (output_arrays);
    for (int n=0; n<N_output_arrays; ++n) {
      assert (N_output_arrays==0 or output_arrays[n]);
    }
    
    // Set up interpolation
    fasterp_llocs_t locations (N_interp_points);
    for (int i=0; i<N_interp_points; ++i) {
      locations.maps.AT(i) = interp_maps[i];
    }
    for (int d=0; d<dim; ++d) {
      for (int i=0; i<N_interp_points; ++i) {
        locations.coords[d].AT(i) = interp_coords[d][i];
      }
    }
    fasterp_setup_t const setup (cctkGH, locations, order);
    
    // Interpolate
    vector<int> varinds (N_input_arrays);
    for (int n=0; n<N_input_arrays; ++n) {
      varinds.AT(n) = input_array_indices[n];
    }
    vector<CCTK_REAL *> values (N_input_arrays);
    for (int n=0; n<N_input_arrays; ++n) {
      values.AT(n) = output_arrays[n];
    }
    setup.interpolate (cctkGH, varinds, values);
    
    // Done
    return 0;
  }
  
  
  
  extern "C"
  CCTK_POINTER
  CarpetInterp2_Interp2GridArraysSetup
  (CCTK_POINTER_TO_CONST const cctkGH_,
   CCTK_INT const N_dims,
   CCTK_INT const order,
   CCTK_INT const N_interp_points,
   CCTK_POINTER_TO_CONST const interp_coords_)
  {
    // Check input values and convert types
    cGH const * const cctkGH = static_cast<cGH const *> (cctkGH_);
    assert (cctkGH);
    
    assert (N_dims == dim);
    
    assert (N_interp_points >= 0);
    
    CCTK_REAL const * const * const interp_coords =
      static_cast<CCTK_REAL const * const *> (interp_coords_);
    assert (interp_coords);
    for (int d=0; d<dim; ++d) {
      assert (N_interp_points==0 or interp_coords[d]);
    }
    
    // Set up interpolation
    fasterp_glocs_t locations (N_interp_points);
    for (int d=0; d<dim; ++d) {
      for (int i=0; i<N_interp_points; ++i) {
        locations.coords[d].AT(i) = interp_coords[d][i];
      }
    }
    fasterp_setup_t * const setup =
      new fasterp_setup_t (cctkGH, locations, order);
    
    // Done
    return setup;
  }
  
  extern "C"
  CCTK_INT
  CarpetInterp2_Interp2GridArrays
  (CCTK_POINTER_TO_CONST const cctkGH_,
   CCTK_POINTER_TO_CONST const setup_,
   CCTK_INT const N_input_arrays,
   CCTK_INT const * const input_array_indices,
   CCTK_INT const N_output_arrays,
   CCTK_POINTER const output_arrays_)
  {
    // Check input values and convert types
    cGH const * const cctkGH = static_cast<cGH const *> (cctkGH_);
    assert (cctkGH);
    
    fasterp_setup_t const * const setup =
      static_cast<fasterp_setup_t const *> (setup_);
    assert (setup);
    
    assert (N_input_arrays >= 0);
    
    assert (input_array_indices);
    for (int n=0; n<N_input_arrays; ++n) {
      assert (input_array_indices[n]>=0 and
              input_array_indices[n]<CCTK_NumVars());
    }
    
    assert (N_output_arrays >= 0);
    assert (N_output_arrays == N_input_arrays);
    
    CCTK_REAL * const * const output_arrays =
      static_cast<CCTK_REAL * const *> (output_arrays_);
    assert (output_arrays);
    for (int n=0; n<N_output_arrays; ++n) {
      assert (N_output_arrays==0 or output_arrays[n]);
    }
    
    // Interpolate
    vector<int> varinds (N_input_arrays);
    for (int n=0; n<N_input_arrays; ++n) {
      varinds.AT(n) = input_array_indices[n];
    }
    vector<CCTK_REAL *> values (N_input_arrays);
    for (int n=0; n<N_input_arrays; ++n) {
      values.AT(n) = output_arrays[n];
    }
    setup->interpolate (cctkGH, varinds, values);
    
    // Done
    return 0;
  }
  
  extern "C"
  CCTK_INT
  CarpetInterp2_Interp2GridArraysFree
  (CCTK_POINTER_TO_CONST const cctkGH_,
   CCTK_POINTER          const setup_)
  {
    // Check input values and convert types
    cGH const * const cctkGH = static_cast<cGH const *> (cctkGH_);
    assert (cctkGH);
    
    fasterp_setup_t const * const setup =
      static_cast<fasterp_setup_t const *> (setup_);
    assert (setup);
    
    delete setup;
    
    // Done
    return 0;
  }
  
  
  
} // namespace CarpetInterp2