aboutsummaryrefslogtreecommitdiff
path: root/src/puncture_tracker.c
blob: ac79fe665ad694e3716b24f99f708a8cd8b3b97d (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
#include <assert.h>
#include <math.h>
#include <stdio.h>

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

#include "util_Table.h"



static int const max_num_tracked = 10;



void
PunctureTracker_Init (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  if (verbose) {
    CCTK_INFO ("Initializing PunctureTracker");
  }
  
  for (int n = 0; n < max_num_tracked; ++ n) {
    if (track[n]) {
      pt_loc_t[n] = cctk_time;
      pt_loc_x[n] = initial_x[n];
      pt_loc_y[n] = initial_y[n];
      pt_loc_z[n] = initial_z[n];
      pt_vel_t[n] = cctk_time;
      pt_vel_x[n] = 0.0;
      pt_vel_y[n] = 0.0;
      pt_vel_z[n] = 0.0;
    } else {
      // Initialise to some sensible but unimportant values
      pt_loc_t[n] = 0.0;
      pt_loc_x[n] = 0.0;
      pt_loc_y[n] = 0.0;
      pt_loc_z[n] = 0.0;
      pt_vel_t[n] = 0.0;
      pt_vel_x[n] = 0.0;
      pt_vel_y[n] = 0.0;
      pt_vel_z[n] = 0.0;
    }
    pt_loc_t_p[n] = 0.0;
    pt_loc_x_p[n] = 0.0;
    pt_loc_y_p[n] = 0.0;
    pt_loc_z_p[n] = 0.0;
  }
}



void
PunctureTracker_Track (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  // Do not track while setting up initial data;
  // time interpolation may fail
  
  if (cctk_iteration == 0) {
    return;
  }
  
  // Some output
  
  if (verbose) {
    CCTK_INFO ("Tracking punctures...");
  }
  
  if (verbose) {
    for (int n = 0; n < max_num_tracked; ++ n) {
      if (track[n]) {
        CCTK_VInfo (CCTK_THORNSTRING,
                    "Puncture #%d is at (%g,%g,%g)",
                    n,
                    (double)pt_loc_x[n],
                    (double)pt_loc_y[n],
                    (double)pt_loc_z[n]);
      }
    }
  }
  
  // Manual time level cycling
  
  for (int n = 0; n < max_num_tracked; ++ n) {
    if (track[n]) {
      pt_loc_t_p[n] = pt_loc_t[n];
      pt_loc_x_p[n] = pt_loc_x[n];
      pt_loc_y_p[n] = pt_loc_y[n];
      pt_loc_z_p[n] = pt_loc_z[n];
      
      pt_loc_t[n] = cctk_time;
      pt_vel_t[n] = cctk_time;
    }
  }
  
  // Interpolate
  
  // Dimensions
  int const dim = 3;
  
  // Interpolation operator
  int const operator_handle =
    CCTK_InterpHandle ("Lagrange polynomial interpolation");
  if (operator_handle < 0) {
    CCTK_WARN (CCTK_WARN_ALERT, "Can't get interpolation handle");
    goto label_return;
  }
  
  // Interpolation parameter table
  int const order = 4;
  int const param_table_handle = Util_TableCreateFromString ("order=4");
  if (param_table_handle < 0) {
    CCTK_WARN (CCTK_WARN_ALERT, "Can't create parameter table");
    goto label_return;
  }
  
  {
    
    // Interpolation coordinate system
    int const coordsys_handle = CCTK_CoordSystemHandle ("cart3d");
    if (coordsys_handle < 0) {
      CCTK_WARN (CCTK_WARN_ALERT, "Can't get coordinate system handle");
      goto label_free_param_table;
    }
    
    // Only processor 0 interpolates
    int const num_points = CCTK_MyProc(cctkGH) == 0 ? max_num_tracked : 0;
    
    // Interpolation coordinates
    assert (dim == 3);
    CCTK_POINTER_TO_CONST interp_coords[3];
    interp_coords[0] = pt_loc_x_p;
    interp_coords[1] = pt_loc_y_p;
    interp_coords[2] = pt_loc_z_p;
    
    // Number of interpolation variables
    int const num_vars = 3;
    
    // Interpolated variables
    assert (num_vars == 3);
    int input_array_indices[3];
    input_array_indices[0] = CCTK_VarIndex ("ADMBase::betax");
    input_array_indices[1] = CCTK_VarIndex ("ADMBase::betay");
    input_array_indices[2] = CCTK_VarIndex ("ADMBase::betaz");
    
    // Interpolation result types
    assert (num_vars == 3);
    CCTK_INT output_array_type_codes[3];
    output_array_type_codes[0] = CCTK_VARIABLE_REAL;
    output_array_type_codes[1] = CCTK_VARIABLE_REAL;
    output_array_type_codes[2] = CCTK_VARIABLE_REAL;
    
    // Interpolation result
    CCTK_REAL pt_betax[max_num_tracked];
    CCTK_REAL pt_betay[max_num_tracked];
    CCTK_REAL pt_betaz[max_num_tracked];
    
    assert (num_vars == 3);
    CCTK_POINTER output_arrays[3];
    output_arrays[0] = pt_betax;
    output_arrays[1] = pt_betay;
    output_arrays[2] = pt_betaz;
    
    // Interpolate
    int ierr;
    if (CCTK_IsFunctionAliased ("InterpGridArrays")) {
      // TODO: use correct array types
      // (CCTK_POINTER[] vs. CCTK_REAL[])
      ierr = InterpGridArrays
        (cctkGH, dim,
         order,
         num_points,
         interp_coords,
         num_vars, input_array_indices,
         num_vars, output_arrays);
    } else {
      ierr = CCTK_InterpGridArrays
        (cctkGH, dim,
         operator_handle, param_table_handle, coordsys_handle,
         num_points,
         CCTK_VARIABLE_REAL,
         interp_coords,
         num_vars, input_array_indices,
         num_vars, output_array_type_codes, output_arrays);
    }
    if (ierr < 0) {
      CCTK_WARN (CCTK_WARN_ALERT, "Interpolation error");
      goto label_free_param_table;
    }
    
    if (CCTK_MyProc(cctkGH) == 0) {
      
      // Some more output
      
      if (verbose && CCTK_MyProc(cctkGH) == 0) {
        for (int n = 0; n < max_num_tracked; ++ n) {
          if (track[n]) {
          CCTK_VInfo (CCTK_THORNSTRING,
                      "Shift at puncture #%d is at (%g,%g,%g)",
                      n,
                      (double)pt_betax[n],
                      (double)pt_betay[n],
                      (double)pt_betaz[n]);
          }
        }
      }

      // Check for NaNs and large shift components
      if (CCTK_MyProc(cctkGH) == 0) {
        for (int n = 0; n < max_num_tracked; ++ n) {
          if (track[n]) {
            CCTK_REAL norm = sqrt(pow((double) pt_betax[n],2) + 
                                  pow((double) pt_betay[n],2) + 
                                  pow((double) pt_betaz[n],2));

            if (!isfinite(norm) || norm > shift_limit) {
              CCTK_VWarn(0, __LINE__, __FILE__, CCTK_THORNSTRING,
                         "Shift at puncture #%d is (%g,%g,%g).  This likely indicates an error in the simulation.",
                         n,
                         (double) pt_betax[n],
                         (double) pt_betay[n],
                         (double) pt_betaz[n]);
            }
          }
        }
      }
      
      // Time evolution
      
      for (int n = 0; n < max_num_tracked; ++ n) {
        if (track[n]) {
          CCTK_REAL const dt = pt_loc_t[n] - pt_loc_t_p[n];
          // First order time integrator
          // Michael Koppitz says this works...
          // if it doesn't, we can make it second order accurate
          pt_loc_x[n] = pt_loc_x_p[n] + dt * (- pt_betax[n]);
          pt_loc_y[n] = pt_loc_y_p[n] + dt * (- pt_betay[n]);
          pt_loc_z[n] = pt_loc_z_p[n] + dt * (- pt_betaz[n]);
          pt_vel_x[n] = - pt_betax[n];
          pt_vel_y[n] = - pt_betay[n];
          pt_vel_z[n] = - pt_betaz[n];
        }
      }
      
    }
    
    // Broadcast result
    
    CCTK_REAL loc_local[6*max_num_tracked]; /* 3 components for location, 3 components for velocity */
    if (CCTK_MyProc(cctkGH) == 0) {
      for (int n = 0; n < max_num_tracked; ++ n) {
        loc_local[                  n] = pt_loc_x[n];
        loc_local[  max_num_tracked+n] = pt_loc_y[n];
        loc_local[2*max_num_tracked+n] = pt_loc_z[n];
        loc_local[3*max_num_tracked+n] = pt_vel_x[n];
        loc_local[4*max_num_tracked+n] = pt_vel_y[n];
        loc_local[5*max_num_tracked+n] = pt_vel_z[n];
      }
    } else {
      for (int n = 0; n < max_num_tracked; ++ n) {
        loc_local[                  n] = 0.0;
        loc_local[  max_num_tracked+n] = 0.0;
        loc_local[2*max_num_tracked+n] = 0.0;
        loc_local[3*max_num_tracked+n] = 0.0;
        loc_local[4*max_num_tracked+n] = 0.0;
        loc_local[5*max_num_tracked+n] = 0.0;
      }
    }
    
    CCTK_REAL loc_global[6*max_num_tracked]; /* 3 components for location, 3 components for velocity */
    
    int const handle_sum = CCTK_ReductionHandle ("sum");
    if (handle_sum < 0) {
      CCTK_WARN (CCTK_WARN_ALERT, "Can't get redunction handle");
      goto label_free_param_table;
    }
    
    int const ierr2 = CCTK_ReduceLocArrayToArray1D
      (cctkGH, -1, handle_sum,
       loc_local, loc_global, 6*max_num_tracked, CCTK_VARIABLE_REAL);
    if (ierr2 < 0) {
      CCTK_WARN (CCTK_WARN_ALERT, "Reduction error");
      goto label_free_param_table;
    }
    
    for (int n = 0; n < max_num_tracked; ++ n) {
      pt_loc_x[n] = loc_global[                  n];
      pt_loc_y[n] = loc_global[  max_num_tracked+n];
      pt_loc_z[n] = loc_global[2*max_num_tracked+n];
      pt_vel_x[n] = loc_global[3*max_num_tracked+n];
      pt_vel_y[n] = loc_global[4*max_num_tracked+n];
      pt_vel_z[n] = loc_global[5*max_num_tracked+n];
    }
    
  }
  
  // Done
  
  // Poor man's exception handling
 label_free_param_table:
  Util_TableDestroy (param_table_handle);
  
 label_return:
  return;
}



void
PunctureTracker_SetPositions (CCTK_ARGUMENTS)
{
  DECLARE_CCTK_ARGUMENTS;
  DECLARE_CCTK_PARAMETERS;
  
  CCTK_REAL dist;

  for (int n = 0; n < max_num_tracked; ++ n) {
    if (track[n]) {
      // store puncture location in spherical surface
      if (which_surface_to_store_info[n] != -1) {
	int sn = which_surface_to_store_info[n];
	
	sf_centroid_x[sn] = pt_loc_x[n];
	sf_centroid_y[sn] = pt_loc_y[n];
	sf_centroid_z[sn] = pt_loc_z[n];
	
	sf_active[sn] = 1;
	sf_valid[sn] = 1;

        if (verbose) {
          CCTK_VInfo (CCTK_THORNSTRING,
                      "Setting spherical surface %d centroid from puncture #%d to (%g,%g,%g)",
                      sn,
                      n,
                      (double)pt_loc_x[n],
                      (double)pt_loc_y[n],
                      (double)pt_loc_z[n]);
        }
      }
    }
  }

  if ( modify_puncture[0]>=0 && modify_puncture[0]<max_num_tracked &&
       modify_puncture[1]>=0 && modify_puncture[1]<max_num_tracked &&
       modify_puncture[0]!=modify_puncture[1] ) {

    if (track[modify_puncture[0]] && track[modify_puncture[1]]) {
 
      dist = sqrt ( pow ( pt_loc_x[modify_puncture[0]] 
                          - pt_loc_x[modify_puncture[1]], 2 )
                  + pow ( pt_loc_y[modify_puncture[0]] 
                          - pt_loc_y[modify_puncture[1]], 2 )
                  + pow ( pt_loc_z[modify_puncture[0]] 
                          - pt_loc_z[modify_puncture[1]], 2 ) );
 
      if ( dist < modify_distance ) {
 
        if ( new_reflevel_number[0] > -1 ) {
          if (verbose) {
            CCTK_VInfo (CCTK_THORNSTRING,
                      "Setting the number of refinement levels to %d for refinement region #%d",
                       new_reflevel_number[0],modify_puncture[0]);
          }
          num_levels[modify_puncture[0]] = new_reflevel_number[0];
        }

        if ( new_reflevel_number[1] > -1 ) {
          if (verbose) {
            CCTK_VInfo (CCTK_THORNSTRING,
                      "Setting the number of refinement levels to %d for refinement region #%d",
                       new_reflevel_number[1],modify_puncture[1]);
          }
          num_levels[modify_puncture[1]] = new_reflevel_number[1];
        }

      }
      
    }

  }
 
}