GRAVI Pipeline Reference Manual 1.9.2
Loading...
Searching...
No Matches
gravi_astrometry.c
Go to the documentation of this file.
1/*
2 * This file is part of the GRAVI Pipeline
3 * Copyright (C) 2002,2003 European Southern Observatory
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
18 */
19
29#include "gravi_astrometry.h"
30#include "gravi_dfs.h"
31#include "gravi_pfits.h"
32#include "gravi_utils.h"
33
34#include <assert.h>
35#include <ctype.h>
36#include <math.h>
37
38#include <gsl/gsl_matrix.h>
39#include <gsl/gsl_vector.h>
40#include <gsl/gsl_blas.h>
41#include <gsl/gsl_complex.h>
42#include <gsl/gsl_complex_math.h>
43#include <gsl/gsl_statistics.h>
44#include <gsl/gsl_multimin.h>
45
46#define PI 3.14159265359
47#define TWOPI 6.28318530718
48#define MAS_TO_RAD (2.0 * PI / (3600 * 360 * 1000))
49
50/*-----------------------------------------------------------------------------
51 Private prototypes
52 -----------------------------------------------------------------------------*/
53
57 double sobj_x, sobj_y;
58 double dit;
59 double mjd;
60
63 cpl_boolean swap;
64
65 gsl_vector *wave;
66 gsl_vector *u, *v;
67 gsl_matrix *ucoord, *vcoord;
68
69 gsl_matrix_complex *visdata;
70 gsl_matrix_complex *viserr;
71 gsl_matrix_complex *visdata_ft;
72
73 int nflag;
74 gsl_matrix_int *flag;
75 // gsl_matrix_int **flag_cov;
76
77 gsl_matrix *phase_ref;
78 gsl_matrix *opd_disp;
79 gsl_matrix *phase_met_telfc;
80 gsl_matrix_complex *vis_ref;
81
82 gsl_vector **vis_ref_cov;
83 gsl_vector_complex **vis_ref_pcov;
84
85 gsl_matrix *amp_ref_astro;
86 gsl_matrix_complex *phase_ref_astro;
87};
88
89typedef enum _find_closest_mode {
94
95static gsl_vector *load_vector(cpl_table *table, const char *name) CPL_ATTR_ALLOC;
96static gsl_matrix *load_matrix(cpl_table *table, const char *name) CPL_ATTR_ALLOC;
97static gsl_matrix_complex *load_matrix_complex(cpl_table *table, const char *name) CPL_ATTR_ALLOC;
98
99static cpl_error_code gravi_astrometry_scale_visibilities(astro_data *self, double factor);
100static cpl_error_code gravi_astrometry_mul_visibilities(astro_data *self, const gsl_vector *factor);
101static cpl_error_code gravi_astrometry_add_phase(astro_data *self, const gsl_matrix *phase);
102static cpl_error_code gravi_astrometry_recentre_phase(astro_data *self, double xvalue, double yvalue);
103static cpl_size gravi_astrometry_find_closest_mjd(astro_data *self, astro_data **others, cpl_size n_other, find_closest_mode mode);
104
105static gsl_vector *average_vector_over_dits(gsl_vector *v, int ndit, int nchannel) CPL_ATTR_ALLOC;
106static gsl_matrix *average_matrix_over_dits(gsl_matrix *m, int ndit, int nchannel, int nwave, gsl_matrix_int *flag) CPL_ATTR_ALLOC;
107static gsl_matrix_complex *average_matrix_complex_over_dits(gsl_matrix_complex *m, int ndit, int nchannel, int nwave, gsl_matrix_int *flag) CPL_ATTR_ALLOC;
108
109static astro_data *gravi_astrometry_average_over_dits(astro_data *self) CPL_ATTR_ALLOC;
110
113 int n1, n2;
115
116static gsl_matrix_complex *gravi_astrometry_calculate_visref_swap(astro_data **data, cpl_size ndata, double ra, double dec) CPL_ATTR_ALLOC;
117static gsl_matrix_complex *gravi_astrometry_create_swap_reference(astro_data **swap_data, cpl_size nswap) CPL_ATTR_ALLOC;
118static double gravi_astrometry_calculate_chi2(const gsl_vector *X, void *params);
119static cpl_error_code gravi_astrometry_minimise_chi2_grid(gravi_astrometry_model_params *params, const gsl_vector *ra_grid, const gsl_vector *dec_grid, gsl_vector *ra_dec_out, gsl_matrix *chi2_map);
120static cpl_error_code gravi_astrometry_minimise_chi2_descent(gravi_astrometry_model_params *params, double ra_guess, double dec_guess, gsl_vector *Xsolve);
121
122/*-----------------------------------------------------------------------------
123 Function code
124 -----------------------------------------------------------------------------*/
125
126static int _gsl_vector_int_sum(const gsl_vector_int *a)
127{
128 const size_t N = a->size;
129 const size_t stride = a->stride;
130 int sum = 0.0;
131 size_t i;
132
133 for (i = 0; i < N; i++)
134 sum += a->data[i * stride];
135 return sum;
136}
137
141static gsl_vector *load_vector(cpl_table *table, const char *name)
142{
143 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
144 cpl_ensure(name != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
145
146 cpl_size depth = cpl_table_get_column_depth(table, name);
147 if (depth == - 1) {
148 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "Field %s not found", name);
149 return NULL;
150 } else if (depth > 0) {
151 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE, "Field %s not of scalar type", name);
152 return NULL;
153 }
154
155 cpl_size nelem = cpl_table_get_nrow(table);
156 double * cpl_vec = cpl_table_get_data_double(table, name);
157
158 gsl_vector *gsl_vec = gsl_vector_alloc(nelem);
159 memcpy(gsl_vec->data, cpl_vec, nelem * sizeof(double));
160
161 return gsl_vec;
162}
163
167static gsl_matrix *load_matrix(cpl_table *table, const char *name)
168{
169 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
170 cpl_ensure(name != NULL, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
171
172 cpl_size depth = cpl_table_get_column_depth(table, name);
173 if (depth == -1) {
174 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "Field %s not found", name);
175 return NULL;
176 } else if (depth == 0) {
177 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE, "Field %s not of 2D array type", name);
178 return NULL;
179 }
180
181 cpl_size nrow = cpl_table_get_nrow(table);
182 cpl_size ncol = cpl_table_get_column_dimension(table, name, 0);
183 const cpl_array ** cpl_matrix = cpl_table_get_data_array_const(table, name);
184
185 gsl_matrix *gsl_matrix = gsl_matrix_alloc(nrow, ncol);
186 for (int i = 0; i < nrow; i++) {
187 gsl_vector_const_view row_view = gsl_vector_const_view_array(
188 cpl_array_get_data_double_const(cpl_matrix[i]), ncol);
189 gsl_matrix_set_row(gsl_matrix, i, &row_view.vector);
190 }
191
192 return gsl_matrix;
193}
194
198static gsl_matrix_complex *load_matrix_complex(cpl_table *table, const char *name)
199{
200 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
201 cpl_ensure(name != NULL, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
202
203 cpl_size depth = cpl_table_get_column_depth(table, name);
204 if (depth == - 1) {
205 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "Field %s not found", name);
206 return NULL;
207 } else if (depth == 0) {
208 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE, "Field %s not of 2D array type", name);
209 return NULL;
210 }
211
212 cpl_size nrow = cpl_table_get_nrow(table);
213 cpl_size ncol = cpl_table_get_column_dimension(table, name, 0);
214 const cpl_array ** cpl_matrix = cpl_table_get_data_array_const(table, name);
215
216 gsl_matrix_complex *gsl_matrix = gsl_matrix_complex_alloc(nrow, ncol);
217 for (int i = 0; i < nrow; i++) {
218 for (int j = 0; j < ncol; j++) {
219 double complex cpl_matrix_val = cpl_array_get_complex(cpl_matrix[i], j, NULL);
220 gsl_complex gsl_matrix_val = gsl_complex_rect(creal(cpl_matrix_val), cimag(cpl_matrix_val));
221 gsl_matrix_complex_set(gsl_matrix, i, j, gsl_matrix_val);
222 }
223 }
224
225 return gsl_matrix;
226}
227
234static cpl_error_code gravi_astrometry_scale_visibilities(astro_data *self, double factor)
235{
236 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
237
238 gsl_complex complex_factor = gsl_complex_rect(factor, 0);
239 gsl_matrix_complex_scale(self->visdata, complex_factor);
240 gsl_matrix_complex_scale(self->visdata_ft, complex_factor);
241 gsl_matrix_complex_scale(self->vis_ref, complex_factor);
242
243 double factor2 = factor * factor;
244 for (int i = 0; i < self->ndit * self->nchannel; i++) {
245 gsl_vector_scale(self->vis_ref_cov[i], factor2);
246 gsl_vector_complex_scale(self->vis_ref_pcov[i], gsl_complex_rect(factor2, 0));
247 }
248
249 return CPL_ERROR_NONE;
250}
251
258static cpl_error_code gravi_astrometry_mul_visibilities(astro_data *self, const gsl_vector *factor)
259{
260 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
261 cpl_ensure_code(factor, CPL_ERROR_NULL_INPUT);
262 cpl_ensure_code((size_t) (self->ndit * self->nchannel) == factor->size, CPL_ERROR_INCOMPATIBLE_INPUT);
263
264 gsl_matrix_complex *complex_factor = gsl_matrix_complex_alloc(self->ndit * self->nchannel, self->nwave);
265 gsl_matrix_complex *complex_factor_ft = gsl_matrix_complex_alloc(self->ndit * self->nchannel, self->nwave_ft);
266 for (int i = 0; i < self->ndit * self->nchannel; i++) {
267 gsl_complex cf = gsl_complex_rect(gsl_vector_get(factor, i), 0);
268 for (int j = 0; j < self->nwave; j++)
269 gsl_matrix_complex_set(complex_factor, i, j, cf);
270 for (int j = 0; j < self->nwave_ft; j++)
271 gsl_matrix_complex_set(complex_factor_ft, i, j, cf);
272 }
273
274 gsl_matrix_complex_mul_elements(self->visdata, complex_factor);
275 gsl_matrix_complex_mul_elements(self->vis_ref, complex_factor);
276 gsl_matrix_complex_mul_elements(self->visdata_ft, complex_factor_ft);
277
278 for (int i = 0; i < self->ndit * self->nchannel; i++) {
279 double factor2 = gsl_vector_get(factor, i);
280 factor2 *= factor2;
281 gsl_vector_scale(self->vis_ref_cov[i], factor2);
282 gsl_vector_complex_scale(self->vis_ref_pcov[i], gsl_complex_rect(factor2, 0));
283 }
284
285 FREE(gsl_matrix_complex_free, complex_factor);
286 FREE(gsl_matrix_complex_free, complex_factor_ft);
287
288 return CPL_ERROR_NONE;
289}
290
294static cpl_error_code gravi_astrometry_add_phase(astro_data *self, const gsl_matrix *phase)
295{
296 gsl_matrix_complex *phi = gsl_matrix_complex_alloc(self->ndit * self->nchannel, self->nwave);
297 for (int i = 0; i < self->ndit * self->nchannel; i++) {
298 for (int j = 0; j < self->nwave; j++) {
299 gsl_complex phi_val = gsl_complex_polar(1.0, gsl_matrix_get(phase, i, j));
300 gsl_matrix_complex_set(phi, i, j, phi_val);
301 }
302 }
303
304 gsl_matrix_complex_mul_elements(self->vis_ref, phi);
305
306 gsl_vector_complex *cov = gsl_vector_complex_alloc(self->nwave);
307 gsl_vector_complex *pcov = gsl_vector_complex_alloc(self->nwave);
308 gsl_vector_complex *tmp = gsl_vector_complex_alloc(self->nwave);
309
310 for (int i = 0; i < self->ndit * self->nchannel; i++) {
311 gsl_vector_complex_set_zero(cov);
312 gsl_vector_complex_set_zero(pcov);
313 gsl_vector_complex_set_zero(tmp);
314
315 // Conjugate vector of phase
316 gsl_vector_complex_const_view phi_row = gsl_matrix_complex_const_row(phi, i);
317 gsl_vector_complex_memcpy(tmp, &phi_row.vector);
318 gsl_vector_view phi_imag = gsl_vector_complex_imag(tmp);
319 gsl_vector_scale(&phi_imag.vector, -1);
320
321 // Need complex vector of C
322 gsl_vector_view cov_real = gsl_vector_complex_real(cov);
323 gsl_vector_memcpy(&cov_real.vector, self->vis_ref_cov[i]);
324
325 // Adjust Cov
326 gsl_vector_complex_mul(tmp, cov); // D = conj(phi) * C
327 gsl_vector_complex_mul(tmp, &phi_row.vector); // D * phi
328 cov_real = gsl_vector_complex_real(tmp);
329 gsl_vector_memcpy(self->vis_ref_cov[i], &cov_real.vector);
330
331 // Adjust Pcov
332 gsl_vector_complex_memcpy(pcov, self->vis_ref_pcov[i]);
333 gsl_vector_complex_mul(pcov, &phi_row.vector);
334 gsl_vector_complex_mul(pcov, &phi_row.vector);
335 gsl_vector_complex_memcpy(self->vis_ref_pcov[i], pcov);
336 }
337
338 FREE(gsl_vector_complex_free, cov);
339 FREE(gsl_vector_complex_free, pcov);
340 FREE(gsl_vector_complex_free, tmp);
341 FREE(gsl_matrix_complex_free, phi);
342
343 return CPL_ERROR_NONE;
344}
345
349static cpl_error_code gravi_astrometry_recentre_phase(astro_data *self, double xvalue, double yvalue)
350{
351 gsl_matrix *this_u = gsl_matrix_alloc(self->ndit * self->nchannel, self->nwave);
352 gsl_matrix_memcpy(this_u, self->ucoord);
353 gsl_matrix_scale(this_u, xvalue);
354
355 gsl_matrix *this_v = gsl_matrix_alloc(self->ndit * self->nchannel, self->nwave);
356 gsl_matrix_memcpy(this_v, self->vcoord);
357 gsl_matrix_scale(this_v, yvalue);
358
359 gsl_matrix_add(this_u, this_v);
360 gsl_matrix_scale(this_u, TWOPI * MAS_TO_RAD);
361
362 gravi_astrometry_add_phase(self, this_u);
363
364 FREE(gsl_matrix_free, this_u);
365 FREE(gsl_matrix_free, this_v);
366
367 return CPL_ERROR_NONE;
368}
369
373cpl_error_code gravi_astrometry_filter_ftflux(astro_data *self, double threshold)
374{
375 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
376
377 int npoints = 0;
378 for (int i = 0; i < self->ndit * self->nchannel; i++) {
379 gsl_vector_complex_const_view visft_row = gsl_matrix_complex_const_row(self->visdata_ft, i);
380 double row_abs_sum = 0.0;
381 for (int j = 0; j < self->nwave_ft; j++)
382 row_abs_sum += gsl_complex_abs(gsl_vector_complex_get(&visft_row.vector, j));
383 row_abs_sum /= self->nwave_ft;
384
385 if (row_abs_sum < threshold) {
386 npoints += 233;
387 for (int k = 0; k < self->nwave; k++) {
388 gsl_matrix_int_set(self->flag, i, k, 1);
389 gsl_matrix_complex_set(self->visdata, i, k, GSL_COMPLEX_ZERO);
390 // for (int k = 0; k < self->nwave; k++) {
391 // gsl_matrix_int_set(self->flag_cov[i], j, k, 1);
392 // gsl_matrix_int_set(self->flag_cov[i], k, j, 1);
393 // }
394 }
395 }
396 // if (row_abs_sum < threshold) {
397 // npoints += 233;
398 // int dit = i / self->nchannel;
399 // // int cha = i % self->nchannel;
400 // // printf("%d %d %d\n", i, dit, cha);
401 // for (int j = 0; j < self->nchannel; j++) {
402 // int z = j + dit * self->nchannel;
403 // for (int k = 0; k < self->nwave; k++) {
404 // gsl_matrix_int_set(self->flag, z, k, 1);
405 // gsl_matrix_complex_set(self->visdata, z, k, GSL_COMPLEX_ZERO);
406 // // for (int k = 0; k < self->nwave; k++) {
407 // // gsl_matrix_int_set(self->flag_cov[i], j, k, 1);
408 // // gsl_matrix_int_set(self->flag_cov[i], k, j, 1);
409 // // }
410 // }
411 // }
412 // }
413 }
414 cpl_msg_debug(cpl_func, "Flagging %d points below FT threshold %.3e", npoints, threshold);
415 self->nflag += npoints;
416
417 return CPL_ERROR_NONE;
418}
419
424{
425 gsl_vector *ftflux_scale = gsl_vector_alloc(self->ndit * self->nchannel);
426 for (int i = 0; i < self->ndit * self->nchannel; i++) {
427 double scale = 0;
428 for (int j = 0; j < self->nwave_ft; j++)
429 scale += gsl_complex_abs(gsl_matrix_complex_get(self->visdata_ft, i, j));
430 gsl_vector_set(ftflux_scale, i, self->nwave / scale);
431 }
432
433 gravi_astrometry_mul_visibilities(self, ftflux_scale);
434 FREE(gsl_vector_free, ftflux_scale);
435 return CPL_ERROR_NONE;
436}
437
442{
443 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
444
445 cpl_propertylist *hdr = gravi_data_get_header(data);
446
447 astro_data *self = cpl_malloc(sizeof(astro_data));
448
449 const char *filename = cpl_propertylist_get_string(hdr, "PIPEFILE");
450 self->filename = cpl_strdup(filename);
451
452 const char *insname = cpl_propertylist_get_string(hdr, "TELESCOP");
453 self->insname = cpl_strdup(insname);
454
455 self->nchannel = 6;
456 if (cpl_propertylist_has(hdr, "ESO TPL NDIT OBJECT"))
457 self->ndit = cpl_propertylist_get_int(hdr, "ESO TPL NDIT OBJECT");
458 else
459 self->ndit = cpl_propertylist_get_int(hdr, "ESO TPL NDIT SKY");
460 self->dit = cpl_propertylist_get_double(hdr, "ESO DET2 SEQ1 DIT");
461
462 self->swap = (!strcmp(cpl_propertylist_get_string(hdr, "ESO INS SOBJ SWAP"), "YES"));
463 cpl_msg_debug(cpl_func, "SWAP keyword is: %s", cpl_propertylist_get_string(hdr, "ESO INS SOBJ SWAP"));
464 // self->swap_astrometry_guess = NULL;
465 // self->swap_astrometry_fit = NULL;
466
467 self->sobj_x = cpl_propertylist_get_double(hdr, "ESO INS SOBJ X");
468 self->sobj_y = cpl_propertylist_get_double(hdr, "ESO INS SOBJ Y");
469 cpl_msg_debug(cpl_func, "Fiber position is: (%f, %f)", self->sobj_x, self->sobj_y);
470
471 self->mjd = cpl_propertylist_get_double(hdr, "MJD-OBS");
472
473 /* Assigned by gravi_astrometry_create_phase_reference */
474 self->amp_ref_astro = NULL;
475 self->phase_ref_astro = NULL;
476
477 /* Get size of wavelength axis */
478 cpl_propertylist *wave_plist = gravi_data_get_oi_wave_plist(data, GRAVI_SC, 0, 1);
479 self->nwave = cpl_propertylist_get_int(wave_plist, "NWAVE");
480
481 wave_plist = gravi_data_get_oi_wave_plist(data, GRAVI_FT, 0, 1);
482 self->nwave_ft = cpl_propertylist_get_int(wave_plist, "NWAVE");
483
484 /* Get wavelength and cast to double */
485 cpl_table *wave_table = gravi_data_get_oi_wave(data, GRAVI_SC, 0, 1);
486 cpl_table_cast_column(wave_table, "EFF_WAVE", "EFF_WAVE", CPL_TYPE_DOUBLE);
487 self->wave = load_vector(wave_table, "EFF_WAVE");
488
489 /* Load visibility data from VIS_OI */
490 cpl_table *oivis_table = gravi_data_get_oi_vis(data, GRAVI_SC, 0, 1);
491
492 self->u = load_vector(oivis_table, "UCOORD");
493 self->v = load_vector(oivis_table, "VCOORD");
494
495 /* uv divided by wavelength is more generally useful */
496 self->ucoord = gsl_matrix_alloc(self->ndit * self->nchannel, self->nwave);
497 self->vcoord = gsl_matrix_alloc(self->ndit * self->nchannel, self->nwave);
498 for (int j = 0; j < self->nwave; j++) {
499 gsl_vector_view col_view;
500 col_view = gsl_matrix_column(self->ucoord, j);
501 gsl_vector_memcpy(&col_view.vector, self->u);
502 col_view = gsl_matrix_column(self->vcoord, j);
503 gsl_vector_memcpy(&col_view.vector, self->v);
504 }
505
506 for (int i = 0; i < self->ndit * self->nchannel; i++) {
507 gsl_vector_view row_view;
508 row_view = gsl_matrix_row(self->ucoord, i);
509 gsl_vector_div(&row_view.vector, self->wave);
510 row_view = gsl_matrix_row(self->vcoord, i);
511 gsl_vector_div(&row_view.vector, self->wave);
512 }
513
514 self->visdata = load_matrix_complex(oivis_table, "VISDATA");
515 self->viserr = load_matrix_complex(oivis_table, "VISERR");
516 self->visdata_ft = load_matrix_complex(oivis_table, "VISDATA_FT");
517
518 self->phase_ref = load_matrix(oivis_table, "PHASE_REF");
519 self->opd_disp = load_matrix(oivis_table, "OPD_DISP");
520 self->phase_met_telfc = load_matrix(oivis_table, "PHASE_MET_TELFC");
521
522 /* Load FLAG */
523 const cpl_array ** cpl_flag = cpl_table_get_data_array_const(oivis_table, "FLAG");
524 const int* cpl_rej = NULL;
525 if (cpl_table_has_column(oivis_table, "REJECTION_FLAG"))
526 cpl_rej = cpl_table_get_data_int_const(oivis_table, "REJECTION_FLAG");
527
528 self->nflag = 0;
529 self->flag = gsl_matrix_int_alloc(self->ndit * self->nchannel, self->nwave);
530 for (int i = 0; i < self->ndit * self->nchannel; i++) {
531 // int rej_val = ((cpl_rej ? cpl_rej[i] : 0) & 19) > 0;
532 int rej_val = (cpl_rej ? cpl_rej[i] : 0) > 0;
533 for (int j = 0; j < self->nwave; j++) {
534 int flag_val = cpl_array_get_int(cpl_flag[i], j, NULL);
535 flag_val |= rej_val;
536 gsl_matrix_int_set(self->flag, i, j, flag_val);
537 self->nflag += flag_val;
538 }
539 }
540
541 /* Clean flagged visData */
542 int ncleaned = 0;
543 for (int i = 0; i < self->ndit * self->nchannel; i++) {
544 for (int j = 0; j < self->nwave; j++) {
545 if (gsl_matrix_int_get(self->flag, i, j)) {
546 ++ncleaned;
547 gsl_matrix_complex_set(self->visdata, i, j, GSL_COMPLEX_ZERO);
548 }
549 }
550 }
551 cpl_msg_debug(cpl_func, "Cleaned %d flagged values in VISDATA", ncleaned);
552
553 /* Construct flags for covariance matrices (unused) */
554 // self->flag_cov = cpl_malloc(self->ndit * self->nchannel * sizeof(gsl_matrix_int*));
555 // for (int i = 0; i < self->ndit * self->nchannel; i++) {
556 // self->flag_cov[i] = gsl_matrix_int_alloc(self->nwave, self->nwave);
557 // for (int j = 0; j < self->nwave; j++) {
558 // if (gsl_matrix_int_get(self->flag, i, j)) {
559 // for (int k = 0; k < self->nwave; k++) {
560 // gsl_matrix_int_set(self->flag_cov[i], j, k, 1);
561 // gsl_matrix_int_set(self->flag_cov[i], k, j, 1);
562 // }
563 // }
564 // }
565 // }
566
567 /* Compute VIS_REF (visData_phasedFT per manual §10.26) */
568 self->vis_ref = gsl_matrix_complex_alloc(self->ndit * self->nchannel, self->nwave);
569 for (int i = 0; i < self->ndit * self->nchannel; i++) {
570 for (int j = 0; j < self->nwave; j++) {
571 gsl_complex vis_ref_val = gsl_matrix_complex_get(self->visdata, i, j);
572 gsl_complex pha = gsl_complex_polar(1.0, gsl_matrix_get(self->phase_ref, i, j));
573 vis_ref_val = gsl_complex_mul(vis_ref_val, pha);
574 gsl_matrix_complex_set(self->vis_ref, i, j, vis_ref_val);
575 }
576 }
577
578 /* Construct covariance and psuedo-covariance matrices for VIS_REF */
579 /* These are diagonal, so we only store a vector of the diagonal terms */
580 self->vis_ref_cov = cpl_malloc(self->ndit * self->nchannel * sizeof(gsl_vector*));
581 self->vis_ref_pcov = cpl_malloc(self->ndit * self->nchannel * sizeof(gsl_vector_complex*));
582 gsl_vector *vreal = gsl_vector_alloc(self->nwave);
583 gsl_vector *vimag = gsl_vector_alloc(self->nwave);
584 for (int i = 0; i < self->ndit * self->nchannel; i++) {
585 // Diagonal of |err|^2 for this dit and channel
586 gsl_vector_complex_const_view err_vals = gsl_matrix_complex_const_row(self->viserr, i);
587 gsl_vector_const_view real_vals = gsl_vector_complex_const_real(&err_vals.vector);
588 gsl_vector_const_view imag_vals = gsl_vector_complex_const_imag(&err_vals.vector);
589 // Real part squared
590 gsl_vector_memcpy(vreal, &real_vals.vector);
591 gsl_vector_mul(vreal, vreal);
592 // Imag part squared
593 gsl_vector_memcpy(vimag, &imag_vals.vector);
594 gsl_vector_mul(vimag, vimag);
595
596 // Form Cov
597 gsl_vector *cov_diag = gsl_vector_alloc(self->nwave);
598 gsl_vector_memcpy(cov_diag, vreal);
599 gsl_vector_add(cov_diag, vimag);
600 self->vis_ref_cov[i] = cov_diag;
601
602 // Form PCov
603 gsl_vector_complex *pcov_diag = gsl_vector_complex_calloc(self->nwave);
604 gsl_vector_view pcov_diag_real = gsl_vector_complex_real(pcov_diag);
605 gsl_vector_memcpy(&pcov_diag_real.vector, vreal);
606 gsl_vector_sub(&pcov_diag_real.vector, vimag);
607
608 for (int j = 0; j < self->nwave; j++) {
609 // pha = pcov_diag * exp(2i*phase_ref)
610 gsl_complex pcov_val = gsl_complex_polar(1.0, 2 * gsl_matrix_get(self->phase_ref, i, j));
611 pcov_val = gsl_complex_mul(gsl_vector_complex_get(pcov_diag, j), pcov_val);
612 gsl_vector_complex_set(pcov_diag, j, pcov_val);
613 }
614 self->vis_ref_pcov[i] = pcov_diag;
615 }
616 FREE(gsl_vector_free, vreal);
617 FREE(gsl_vector_free, vimag);
618
619 /* Normalise visibilities by the DIT */
620 gravi_astrometry_scale_visibilities(self, 1.0 / self->dit);
621
622 /* Compute metrology and dispersion phase corrections */
623 gsl_matrix *wave_opd_disp = gsl_matrix_alloc(self->ndit * self->nchannel, self->nwave);
624 gsl_matrix_memcpy(wave_opd_disp, self->opd_disp);
625
626 gsl_vector *wavenumber = gsl_vector_alloc(self->nwave);
627 gsl_vector_set_all(wavenumber, TWOPI);
628 gsl_vector_div(wavenumber, self->wave);
629
630 /* 2*pi/wave * opdDisp */
631 for (int i = 0; i < self->ndit * self->nchannel; i++) {
632 gsl_vector_view wave_opd_dist_row = gsl_matrix_row(wave_opd_disp, i);
633 gsl_vector_mul(&wave_opd_dist_row.vector, wavenumber);
634 }
635
636 /* Total correction should be subtracted from phase, so negate first */
637 gsl_matrix_add(wave_opd_disp, self->phase_met_telfc);
638 gsl_matrix_scale(wave_opd_disp, -1);
639 gravi_astrometry_add_phase(self, wave_opd_disp);
640
641 FREE(gsl_vector_free, wavenumber);
642 FREE(gsl_matrix_free, wave_opd_disp);
643
644 return self;
645}
646
647cpl_error_code gravi_astrometry_dump(astro_data *self, FILE *handle)
648{
649 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
650 cpl_ensure_code(handle, CPL_ERROR_NULL_INPUT);
651
652 fprintf(handle, "%s\n", self->filename);
653 fprintf(handle, "%s\n", self->insname);
654
655 fprintf(handle, "%d %d %d\n", self->ndit, self->nchannel, self->nwave);
656
657 gsl_matrix_fprintf(handle, self->ucoord, "%f");
658 gsl_matrix_fprintf(handle, self->vcoord, "%f");
659 gsl_matrix_int_fprintf(handle, self->flag, "%d");
660 gsl_matrix_complex_fprintf(handle, self->visdata, "%f");
661 gsl_matrix_complex_fprintf(handle, self->vis_ref, "%f");
662
663 return CPL_ERROR_NONE;
664}
665
667{
668 FREE(cpl_free, self->filename);
669 FREE(cpl_free, self->insname);
670 FREE(gsl_vector_free, self->wave);
671 FREE(gsl_vector_free, self->u);
672 FREE(gsl_vector_free, self->v);
673
674 FREE(gsl_matrix_complex_free, self->visdata);
675 FREE(gsl_matrix_complex_free, self->visdata_ft);
676 FREE(gsl_matrix_complex_free, self->viserr);
677 FREE(gsl_matrix_free, self->phase_ref);
678 FREE(gsl_matrix_free, self->phase_met_telfc);
679 FREE(gsl_matrix_free, self->opd_disp);
680 FREE(gsl_matrix_complex_free, self->vis_ref);
681 FREE(gsl_matrix_int_free, self->flag);
682
683 // FREELOOP(gsl_matrix_int_free, self->flag_cov, self->ndit * self->nchannel);
684 FREELOOP(gsl_vector_free, self->vis_ref_cov, self->ndit * self->nchannel);
685 FREELOOP(gsl_vector_complex_free, self->vis_ref_pcov, self->ndit * self->nchannel);
686
687 // FREE(gsl_vector_free, self->swap_astrometry_guess);
688 // FREE(gsl_vector_free, self->swap_astrometry_fit);
689 FREE(gsl_matrix_free, self->amp_ref_astro);
690 FREE(gsl_matrix_complex_free, self->phase_ref_astro);
691
692 cpl_free(self);
693}
694
696{
697 cpl_ensure(self, CPL_ERROR_NULL_INPUT, -1);
698
699 double abs_val = 0;
700 for (int i = 0; i < self->ndit * self->nchannel; i++) {
701 for (int j = 0; j < self->nwave_ft; j++) {
702 abs_val += gsl_complex_abs(gsl_matrix_complex_get(self->visdata_ft, i, j));
703 }
704 }
705 return abs_val / (self->ndit * self->nchannel * self->nwave_ft);
706}
707
719cpl_size gravi_astrometry_find_closest_mjd(astro_data *self, astro_data **others, cpl_size n_other, find_closest_mode mode)
720{
721 cpl_ensure(self, CPL_ERROR_NULL_INPUT, -1);
722 cpl_ensure(others, CPL_ERROR_NULL_INPUT, -1);
723
724 cpl_boolean any_before = FALSE, any_after = FALSE;
725 for (int i = 0; i < n_other; i++) {
726 double delta = self->mjd - others[i]->mjd;
727 if (delta >= 0)
728 any_before = CPL_TRUE;
729 if (delta <= 0)
730 any_after = CPL_TRUE;
731 }
732
733 cpl_size imin = -1;
734 double min_delta = INFINITY;
735
736 for (int i = 0; i < n_other; i++) {
737 double delta = self->mjd - others[i]->mjd;
738 switch (mode) {
740 if (any_before && delta >= 0 && delta < min_delta) {
741 imin = i;
742 min_delta = delta;
743 }
744 break;
745case FIND_MODE_AFTER:
746 if (any_after && delta <= 0 && -delta < min_delta) {
747 imin = i;
748 min_delta = -delta;
749 }
750 break;
751default:
752 if (delta <= 0 && fabs(delta) < min_delta) {
753 imin = i;
754 min_delta = delta;
755 }
756 break;
757 }
758 }
759 return imin;
760}
761
772cpl_error_code gravi_astrometry_create_phase_reference(astro_data *self, astro_data **phase_refs, cpl_size nphase, astro_data **swaps, cpl_size nswap, cpl_parameterlist *parlist)
773{
774 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
775 cpl_ensure_code(phase_refs, CPL_ERROR_NULL_INPUT);
776 if (nswap > 0)
777 cpl_ensure_code(swaps, CPL_ERROR_NULL_INPUT);
778 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
779
780 gsl_matrix *amp_ref = NULL, *amp_ref_tmp = NULL;
781 gsl_matrix_complex *vis_ref = NULL, *vis_ref_tmp = NULL;
782
783 /* Parameters */
784 char *calib_strategy = cpl_strdup(cpl_parameter_get_string(
785 cpl_parameterlist_find(parlist, "gravity.astrometry.calib-strategy")));
786 /* just in case: tolerate mixed/lower case arguments */
787 for(int i = 0; calib_strategy[i]; i++)
788 calib_strategy[i] = toupper(calib_strategy[i]);
789
790 CPLCHECK_INT("Could not get parameters");
791
792 cpl_msg_debug(cpl_func, "Calibrating amplitude reference with '%s' strategy", calib_strategy);
793
794 if (!strcmp(calib_strategy, "NONE")) {
795 amp_ref = gsl_matrix_alloc(self->nchannel, self->nwave);
796 vis_ref = gsl_matrix_complex_alloc(self->nchannel, self->nwave);
797 gsl_matrix_set_all(amp_ref, 1.0);
798 gsl_matrix_complex_set_all(vis_ref, GSL_COMPLEX_ONE);
799 } else if (!strcmp(calib_strategy, "SELF")) {
800 // abs first, then mean over dits
801 amp_ref_tmp = gsl_matrix_alloc(self->ndit * self->nchannel, self->nwave);
802 for (int i = 0; i < self->ndit * self->nchannel; i++) {
803 for (int j = 0; j < self->nwave; j++) {
804 gsl_complex vis_ref_val = gsl_matrix_complex_get(self->vis_ref, i, j);
805 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
806 }
807 }
808 amp_ref = average_matrix_over_dits(amp_ref_tmp, self->ndit, self->nchannel, self->nwave, self->flag);
809 FREE(gsl_matrix_free, amp_ref_tmp);
810
811 vis_ref = gsl_matrix_complex_alloc(self->nchannel, self->nwave);
812 for (int i = 0; i < self->nchannel; i++) {
813 for (int j = 0; j < self->nwave; j++) {
814 double amp_ref_val = gsl_matrix_get(amp_ref, i, j);
815 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_rect(amp_ref_val, 0));
816 }
817 }
818 } else {
819 int nphase_used;
820 gsl_vector_int *phase_indices;
821
822 if (!strcmp(calib_strategy, "ALL")) {
823 nphase_used = nphase;
824 phase_indices = gsl_vector_int_alloc(nphase_used);
825 for (int n = 0; n < nphase_used; n++)
826 gsl_vector_int_set(phase_indices, n, n);
827 } else if (!strcmp(calib_strategy, "NEAREST")) {
828 nphase_used = 2;
829 phase_indices = gsl_vector_int_alloc(nphase_used);
830 cpl_size idx_before = gravi_astrometry_find_closest_mjd(self, phase_refs, nphase, FIND_MODE_BEFORE);
831 cpl_size idx_after = gravi_astrometry_find_closest_mjd(self, phase_refs, nphase, FIND_MODE_AFTER);
832 // cpl_msg_debug(cpl_func, "File %s", self->filename);
833 // cpl_msg_debug(cpl_func, "Preceding on-star file %s", phase_refs[idx_before]->filename);
834 // cpl_msg_debug(cpl_func, "Succeeding on-star file %s", phase_refs[idx_after]->filename);
835 if (idx_before == -1 || idx_after == -1) {
836 cpl_error_set_message(cpl_func, CPL_ERROR_ACCESS_OUT_OF_RANGE, "Cannot find time-bracketing phase references");
837 return CPL_ERROR_ACCESS_OUT_OF_RANGE;
838 }
839 gsl_vector_int_set(phase_indices, 0, idx_before);
840 gsl_vector_int_set(phase_indices, 1, idx_after);
841 } else {
842 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT, "Unknown calibration strategy %s", calib_strategy);
843 return CPL_ERROR_ILLEGAL_INPUT;
844 }
845
846 amp_ref = gsl_matrix_calloc(self->nchannel, self->nwave);
847 vis_ref = gsl_matrix_complex_calloc(self->nchannel, self->nwave);
848 amp_ref_tmp = gsl_matrix_alloc(self->nchannel, self->nwave);
849
850 for (int n = 0; n < nphase_used; n++) {
851 cpl_size idx = gsl_vector_int_get(phase_indices, n);
852 astro_data *soi = phase_refs[idx];
853 cpl_msg_debug(cpl_func, "Take idx %lld, filename %s", idx, soi->filename);
854
855 // Average over DITs to get object of shape (n_channel, n_wave)
856 vis_ref_tmp = average_matrix_complex_over_dits(soi->vis_ref, soi->ndit, soi->nchannel, soi->nwave, soi->flag);
857
858 for (int i = 0; i < soi->nchannel; i++) {
859 for (int j = 0; j < soi->nwave; j++) {
860 gsl_complex vis_ref_val = gsl_matrix_complex_get(vis_ref_tmp, i, j);
861 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
862 }
863 }
864
865 gsl_matrix_add(amp_ref, amp_ref_tmp);
866 gsl_matrix_complex_add(vis_ref, vis_ref_tmp);
867 FREE(gsl_matrix_complex_free, vis_ref_tmp);
868 }
869
870 gsl_matrix_scale(amp_ref, 1.0 / nphase_used);
871 gsl_matrix_complex_scale(vis_ref, gsl_complex_rect(1.0 / nphase_used, 0.0));
872
873 for (int i = 0; i < self->nchannel; i++) {
874 for (int j = 0; j < self->nwave; j++) {
875 double theta = gsl_complex_arg(gsl_matrix_complex_get(vis_ref, i, j));
876 double r = gsl_matrix_get(amp_ref, i, j);
877 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_polar(r, theta));
878 }
879 }
880 FREE(gsl_matrix_free, amp_ref_tmp);
881 FREE(gsl_vector_int_free, phase_indices);
882 }
883
884 self->amp_ref_astro = amp_ref;
885 self->phase_ref_astro = vis_ref;
886
887 /* extract phase reference from swaps */
888 if (nswap > 0) {
889 cpl_msg_debug(cpl_func, "Extract phase reference from swaps");
890 /* Should already exist, but will be modified and/or overwritten */
891 amp_ref = self->amp_ref_astro;
892 vis_ref = self->phase_ref_astro;
893
894 /* Take the amplitude reference from the swaps */
895 if (!strcmp(calib_strategy, "SWAP")) {
896 cpl_size n_swapon = 0;
897
898 gsl_matrix_set_zero(amp_ref);
899 gsl_matrix_complex_set_zero(vis_ref);
900
901 for (int n = 0; n < nswap; n++) {
902 astro_data *soi = swaps[n];
903
904 if (soi->swap) {
905 n_swapon++;
906 // abs first, then mean over dits
907 amp_ref_tmp = gsl_matrix_alloc(soi->ndit * soi->nchannel, soi->nwave);
908 for (int i = 0; i < soi->ndit * soi->nchannel; i++) {
909 for (int j = 0; j < soi->nwave; j++) {
910 gsl_complex vis_ref_val = gsl_matrix_complex_get(soi->vis_ref, i, j);
911 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
912 }
913 }
914 gsl_matrix *amp_ref_avg = average_matrix_over_dits(amp_ref_tmp, soi->ndit, soi->nchannel, soi->nwave, soi->flag);
915 gsl_matrix_add(amp_ref, amp_ref_avg);
916
917 FREE(gsl_matrix_free, amp_ref_tmp);
918 FREE(gsl_matrix_free, amp_ref_avg);
919 }
920 }
921
922 gsl_matrix_scale(amp_ref, 1.0 / n_swapon);
923 for (int i = 0; i < self->nchannel; i++) {
924 for (int j = 0; j < self->nwave; j++) {
925 /* Phase ref will be added later. Divide by 2 is because flux is 2x higher in swap */
926 double amp_ref_val = gsl_matrix_get(amp_ref, i, j);
927 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_rect(0.5 * amp_ref_val, 0));
928 }
929 }
930 }
931
932 /* Finally calculate and add the phase ref */
933 gsl_matrix_complex *swap_ref = swaps[0]->phase_ref_astro;
934 for (int i = 0; i < self->nchannel; i++) {
935 for (int j = 0; j < self->nwave; j++) {
936 gsl_complex phi = gsl_matrix_complex_get(swap_ref, i, j);
937 double argphi = gsl_complex_arg(phi);
938
939 gsl_complex vis_ref_val = gsl_matrix_complex_get(vis_ref, i, j);
940 double vis_ref_abs = gsl_complex_abs(vis_ref_val);
941 /* Factor 2 because of the beamsplitter for on-star observations*/
942 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_polar(2 * vis_ref_abs, argphi));
943 }
944 }
945 }
946
947 cpl_free(calib_strategy);
948 return CPL_ERROR_NONE;
949}
950
954// cpl_table *gravi_astrometry_get_phase_reference(astro_data *self)
955// {
956// cpl_table *table = cpl_table_new(self->nchannel);
957// cpl_table_new_column_array(table, "ASTRO_VISREF", CPL_TYPE_DOUBLE_COMPLEX, self->nwave);
958// cpl_table_new_column_array(table, "ASTRO_AMPREF", CPL_TYPE_DOUBLE, self->nwave);
959// cpl_table_new_column_array(table, "ASTRO_COV", CPL_TYPE_DOUBLE, self->ndit * self->nwave);
960// cpl_table_new_column_array(table, "ASTRO_PCOV", CPL_TYPE_DOUBLE_COMPLEX, self->ndit * self->nwave);
961
962// cpl_array *tmp_arr = cpl_array_new(self->nwave, CPL_TYPE_DOUBLE_COMPLEX);
963// for (int i = 0; i < self->nchannel; i++) {
964// for (int j = 0; j < self->nwave; j++) {
965// gsl_complex val = gsl_matrix_complex_get(self->phase_ref_astro, i, j);
966// cpl_array_set_complex(tmp_arr, j, GSL_REAL(val) + I * GSL_IMAG(val));
967// }
968// cpl_table_set_array(table, "ASTRO_VISREF", i, tmp_arr);
969// }
970// FREE(cpl_array_delete, tmp_arr);
971
972// tmp_arr = cpl_array_new(self->nwave, CPL_TYPE_DOUBLE);
973// for (int i = 0; i < self->nchannel; i++) {
974// for (int j = 0; j < self->nwave; j++) {
975// double val = gsl_matrix_get(self->amp_ref_astro, i, j);
976// cpl_array_set(tmp_arr, j, val);
977// }
978// cpl_table_set_array(table, "ASTRO_AMPREF", i, tmp_arr);
979// }
980// FREE(cpl_array_delete, tmp_arr);
981
982// // vis_ref_cov is array of size ndit * nchannel, of vectors of length nwave
983// tmp_arr = cpl_array_new(self->ndit * self->nwave, CPL_TYPE_DOUBLE);
984// for (int i = 0; i < self->nchannel; i++) {
985// for (int j = 0; j < self->ndit; j++) {
986// for (int k = 0; k < self->nwave; k++) {
987// double val = gsl_vector_get(self->vis_ref_cov[i + j * self->nchannel], k);
988// cpl_array_set(tmp_arr, k + j * self->nwave, val);
989// }
990// }
991// cpl_table_set_array(table, "ASTRO_COV", i, tmp_arr);
992// }
993// CPLCHECK_NUL("Could not add ASTRO_COV");
994// FREE(cpl_array_delete, tmp_arr);
995
996// // vis_ref_pcov is array of size ndit * nchannel, of complex vectors of length nwave
997// tmp_arr = cpl_array_new(self->ndit * self->nwave, CPL_TYPE_DOUBLE_COMPLEX);
998// for (int i = 0; i < self->nchannel; i++) {
999// for (int j = 0; j < self->ndit; j++) {
1000// for (int k = 0; k < self->nwave; k++) {
1001// gsl_complex val = gsl_vector_complex_get(self->vis_ref_pcov[i + j * self->nchannel], k);
1002// cpl_array_set_complex(tmp_arr, k + j * self->nwave, GSL_REAL(val) + I * GSL_IMAG(val));
1003// }
1004// }
1005// cpl_table_set_array(table, "ASTRO_PCOV", i, tmp_arr);
1006// }
1007// CPLCHECK_NUL("Could not add ASTRO_PCOV");
1008// FREE(cpl_array_delete, tmp_arr);
1009
1010// return table;
1011// }
1012
1021{
1022 cpl_table *table = cpl_table_new(self->nwave);
1023 cpl_table_new_column_array(table, "ASTRO_VISREF", CPL_TYPE_DOUBLE_COMPLEX, self->nchannel);
1024 cpl_table_new_column_array(table, "ASTRO_AMPREF", CPL_TYPE_DOUBLE, self->nchannel);
1025 cpl_table_new_column_array(table, "ASTRO_COV", CPL_TYPE_DOUBLE, self->ndit * self->nchannel);
1026 cpl_table_new_column_array(table, "ASTRO_PCOV", CPL_TYPE_DOUBLE_COMPLEX, self->ndit * self->nchannel);
1027
1028 cpl_array *tmp_arr = cpl_array_new(self->nchannel, CPL_TYPE_DOUBLE_COMPLEX);
1029 for (int j = 0; j < self->nwave; j++) {
1030 for (int i = 0; i < self->nchannel; i++) {
1031 gsl_complex val = gsl_matrix_complex_get(self->phase_ref_astro, i, j);
1032 cpl_array_set_complex(tmp_arr, i, GSL_REAL(val) + I * GSL_IMAG(val));
1033 }
1034 cpl_table_set_array(table, "ASTRO_VISREF", j, tmp_arr);
1035 }
1036 CPLCHECK_NUL("Could not add ASTRO_VISREF");
1037 FREE(cpl_array_delete, tmp_arr);
1038
1039 tmp_arr = cpl_array_new(self->nchannel, CPL_TYPE_DOUBLE);
1040 for (int j = 0; j < self->nwave; j++) {
1041 for (int i = 0; i < self->nchannel; i++) {
1042 double val = gsl_matrix_get(self->amp_ref_astro, i, j);
1043 cpl_array_set(tmp_arr, i, val);
1044 }
1045 cpl_table_set_array(table, "ASTRO_AMPREF", j, tmp_arr);
1046 }
1047 CPLCHECK_NUL("Could not add ASTRO_AMPREF");
1048 FREE(cpl_array_delete, tmp_arr);
1049
1050 // vis_ref_cov is array of size ndit * nchannel, of vectors of length nwave
1051 tmp_arr = cpl_array_new(self->ndit * self->nchannel, CPL_TYPE_DOUBLE);
1052 for (int j = 0; j < self->nwave; j++) {
1053 for (int i = 0; i < self->ndit * self->nchannel; i++) {
1054 double val = gsl_vector_get(self->vis_ref_cov[i], j);
1055 cpl_array_set(tmp_arr, i, val);
1056 }
1057 cpl_table_set_array(table, "ASTRO_COV", j, tmp_arr);
1058 }
1059 CPLCHECK_NUL("Could not add ASTRO_COV");
1060 FREE(cpl_array_delete, tmp_arr);
1061
1062 // vis_ref_pcov is array of size ndit * nchannel, of vectors of length nwave
1063 tmp_arr = cpl_array_new(self->ndit * self->nchannel, CPL_TYPE_DOUBLE_COMPLEX);
1064 for (int j = 0; j < self->nwave; j++) {
1065 for (int i = 0; i < self->ndit * self->nchannel; i++) {
1066 gsl_complex val = gsl_vector_complex_get(self->vis_ref_pcov[i], j);
1067 cpl_array_set_complex(tmp_arr, i, GSL_REAL(val) + I * GSL_IMAG(val));
1068 }
1069 cpl_table_set_array(table, "ASTRO_PCOV", j, tmp_arr);
1070 }
1071 CPLCHECK_NUL("Could not add ASTRO_PCOV");
1072 FREE(cpl_array_delete, tmp_arr);
1073
1074 return table;
1075}
1076
1077static gsl_vector *average_vector_over_dits(gsl_vector *v, int ndit, int nchannel)
1078{
1079 gsl_vector *vavg = gsl_vector_alloc(nchannel);
1080
1081 for (int i = 0; i < nchannel; i++) {
1082 double avg = gsl_stats_mean(v->data + i, nchannel, ndit);
1083 gsl_vector_set(vavg, i, avg);
1084 }
1085
1086 return vavg;
1087}
1088
1089static gsl_matrix *average_matrix_over_dits(gsl_matrix *m, int ndit, int nchannel, int nwave, gsl_matrix_int *flag)
1090{
1091 gsl_matrix *mavg = gsl_matrix_alloc(nchannel, nwave);
1092
1093 for (int i = 0; i < nchannel; i++) {
1094 for (int j = 0; j < nwave; j++) {
1095 cpl_size start = i * nwave + j;
1096 gsl_vector_const_view dit_view = gsl_vector_const_view_array_with_stride(
1097 m->data + start, nwave * nchannel, ndit);
1098 gsl_vector_int_const_view flag_view = flag ?
1099 gsl_vector_int_const_view_array_with_stride(flag->data + start, nwave * nchannel, ndit) :
1100 (gsl_vector_int_const_view){};
1101
1102 double avg = 0.0;
1103 cpl_size Nvalid = 0;
1104 for (int k = 0; k < ndit; k++) {
1105 if (!flag || gsl_vector_int_get(&flag_view.vector, k) == 0) {
1106 avg += gsl_vector_get(&dit_view.vector, k);
1107 Nvalid++;
1108 }
1109 }
1110 if (Nvalid > 0) {
1111 avg /= Nvalid;
1112 gsl_matrix_set(mavg, i, j, avg);
1113 } else {
1114 gsl_matrix_set(mavg, i, j, 0.0);
1115 }
1116 }
1117 }
1118
1119 return mavg;
1120}
1121
1122static gsl_matrix_complex *average_matrix_complex_over_dits(gsl_matrix_complex *m, int ndit, int nchannel, int nwave, gsl_matrix_int *flag)
1123{
1124 gsl_matrix_complex *mavg = gsl_matrix_complex_alloc(nchannel, nwave);
1125
1126 for (int i = 0; i < nchannel; i++) {
1127 for (int j = 0; j < nwave; j++) {
1128 cpl_size start = i * nwave + j;
1129 gsl_vector_complex_const_view dit_view = gsl_vector_complex_const_view_array_with_stride(
1130 m->data + 2 * start, nwave * nchannel, ndit); /* factor 2 for complex */
1131 gsl_vector_int_const_view flag_view = flag ?
1132 gsl_vector_int_const_view_array_with_stride(flag->data + start, nwave * nchannel, ndit) :
1133 (gsl_vector_int_const_view){};
1134
1135 gsl_complex avg = GSL_COMPLEX_ZERO;
1136 cpl_size Nvalid = 0;
1137 for (int k = 0; k < ndit; k++) {
1138 if (!flag || gsl_vector_int_get(&flag_view.vector, k) == 0) {
1139 avg = gsl_complex_add(avg, gsl_vector_complex_get(&dit_view.vector, k));
1140 Nvalid++;
1141 }
1142 }
1143 if (Nvalid > 0) {
1144 avg = gsl_complex_div_real(avg, Nvalid);
1145 gsl_matrix_complex_set(mavg, i, j, avg);
1146 } else {
1147 gsl_matrix_complex_set(mavg, i, j, GSL_COMPLEX_ZERO);
1148 }
1149 }
1150 }
1151
1152 return mavg;
1153}
1154
1159{
1160 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
1161 cpl_msg_debug(cpl_func, "Averaging file %s", self->filename);
1162
1163 astro_data *flat = cpl_malloc(sizeof(astro_data));
1164
1165 flat->filename = cpl_strdup(self->filename);
1166 flat->insname = cpl_strdup(self->insname);
1167
1168 flat->nchannel = self->nchannel;
1169 flat->ndit = 1;
1170 flat->dit = self->dit;
1171 flat->sobj_x = self->sobj_x;
1172 flat->sobj_y = self->sobj_y;
1173 flat->mjd = self->mjd;
1174 flat->swap = self->swap;
1175 flat->nwave = self->nwave;
1176 flat->nwave_ft = self->nwave_ft;
1177
1178 // flat->swap_astrometry_guess = NULL;
1179 // flat->swap_astrometry_fit = NULL;
1180
1181 flat->wave = gsl_vector_alloc(flat->nwave);
1182 gsl_vector_memcpy(flat->wave, self->wave);
1183
1184 /* Calculate mean in SC frame */
1185 gravi_astrometry_recentre_phase(self, self->sobj_x, self->sobj_y);
1186
1187 flat->u = average_vector_over_dits(self->u, self->ndit, self->nchannel);
1188 flat->v = average_vector_over_dits(self->u, self->ndit, self->nchannel);
1189 flat->ucoord = average_matrix_over_dits(self->ucoord, self->ndit, self->nchannel, self->nwave, self->flag);
1190 flat->vcoord = average_matrix_over_dits(self->vcoord, self->ndit, self->nchannel, self->nwave, self->flag);
1191
1192 flat->visdata = average_matrix_complex_over_dits(self->visdata, self->ndit, self->nchannel, self->nwave, self->flag);
1193 flat->vis_ref = average_matrix_complex_over_dits(self->vis_ref, self->ndit, self->nchannel, self->nwave, self->flag);
1194
1195 /* cov, pcov are array of shape (ndit * nchannel,) of vectors with shape (nwave,)
1196 sum over dits, divide by ndit**2 */
1197 flat->vis_ref_cov = cpl_malloc(flat->nchannel * sizeof(gsl_vector *));
1198 flat->vis_ref_pcov = cpl_malloc(flat->nchannel * sizeof(gsl_vector_complex *));
1199
1200 for (int i = 0; i < self->nchannel; i++) {
1201 flat->vis_ref_cov[i] = gsl_vector_calloc(flat->nwave);
1202 flat->vis_ref_pcov[i] = gsl_vector_complex_calloc(flat->nwave);
1203
1204 for (int j = 0; j < self->ndit; j++) {
1205 cpl_size idx = j * self->nchannel + i;
1206 gsl_vector_add(flat->vis_ref_cov[i], self->vis_ref_cov[idx]);
1207 gsl_vector_complex_add(flat->vis_ref_pcov[i], self->vis_ref_pcov[idx]);
1208 }
1209 gsl_vector_scale(flat->vis_ref_cov[i], 1.0 / (self->ndit * self->ndit));
1210 gsl_vector_complex_scale(flat->vis_ref_pcov[i], gsl_complex_rect(1.0 / (self->ndit * self->ndit), 0));
1211 }
1212
1213 flat->nflag = 0;
1214 flat->flag = gsl_matrix_int_calloc(flat->nchannel, flat->nwave);
1215 for (int i = 0; i < self->nchannel; i++) {
1216 for (int j = 0; j < self->nwave; j++) {
1217 cpl_size start = i * self->nwave + j;
1218 gsl_vector_int_const_view flag_view = gsl_vector_int_const_view_array_with_stride(
1219 self->flag->data + start, self->nwave * self->nchannel, self->ndit);
1220
1221 if (_gsl_vector_int_sum(&flag_view.vector) == self->ndit) {
1222 gsl_matrix_int_set(flat->flag, i, j, 1);
1223 flat->nflag++;
1224 }
1225 }
1226 }
1227
1228 /* Shift flattened back to FT frame */
1229 gravi_astrometry_recentre_phase(flat, -flat->sobj_x, -flat->sobj_y);
1230
1231 /* Shift original back to FT frame */
1232 gravi_astrometry_recentre_phase(self, -self->sobj_x, -self->sobj_y);
1233
1234 return flat;
1235}
1236
1237static gsl_matrix_complex *gravi_astrometry_calculate_visref_swap(astro_data **data, cpl_size ndata, double ra, double dec) {
1238 gsl_matrix_complex *visref = gsl_matrix_complex_calloc(data[0]->nchannel, data[0]->nwave);
1239 gsl_matrix_complex *visref_ditsum = gsl_matrix_complex_alloc(data[0]->nchannel, data[0]->nwave);
1240
1241 gsl_matrix_int *ngood = gsl_matrix_int_calloc(data[0]->nchannel, data[0]->nwave);
1242 gsl_matrix_int *ngood_ditsum = gsl_matrix_int_alloc(data[0]->nchannel, data[0]->nwave);
1243
1244 for (int n = 0; n < ndata; n++) {
1245 astro_data *self = data[n];
1246
1247 for (int j = 0; j < self->nchannel; j++) {
1248 for (int k = 0; k < self->nwave; k++) {
1249 int ngood_jk = 0;
1250 gsl_complex visref_jk = GSL_COMPLEX_ZERO;
1251
1252 for (int i = 0; i < self->ndit; i++) {
1253 cpl_size idx = i * self->nchannel + j;
1254 double uv = gsl_matrix_get(self->ucoord, idx, k) * ra + gsl_matrix_get(self->vcoord, idx, k) * dec;
1255 double phase = uv * TWOPI * MAS_TO_RAD;
1256
1257 gsl_complex phi = gsl_complex_polar(1.0, phase);
1258 if (self->swap)
1259 phi = gsl_complex_conjugate(phi);
1260 gsl_complex visref_ijk = gsl_complex_mul(phi, gsl_matrix_complex_get(self->vis_ref, idx, k));
1261
1262 int flag = gsl_matrix_int_get(self->flag, idx, k);
1263 if (flag == 0) {
1264 ngood_jk++;
1265 visref_jk = gsl_complex_add(visref_jk, visref_ijk);
1266 }
1267 }
1268 gsl_matrix_complex_set(visref_ditsum, j, k, visref_jk);
1269 gsl_matrix_int_set(ngood_ditsum, j, k, ngood_jk);
1270 }
1271 }
1272 gsl_matrix_complex_add(visref, visref_ditsum);
1273 gsl_matrix_int_add(ngood, ngood_ditsum);
1274 }
1275
1276 for (int j = 0; j < data[0]->nchannel; j++) {
1277 for (int k = 0; k < data[0]->nwave; k++) {
1278 int ngood_jk = gsl_matrix_int_get(ngood, j, k);
1279 if (ngood_jk > 0) {
1280 gsl_complex visref_jk = gsl_matrix_complex_get(visref, j, k);
1281 visref_jk = gsl_complex_div_real(visref_jk, ngood_jk);
1282 gsl_matrix_complex_set(visref, j, k, visref_jk);
1283 }
1284 }
1285 }
1286
1287 FREE(gsl_matrix_complex_free, visref_ditsum);
1288 return visref;
1289}
1290
1297static double gravi_astrometry_calculate_chi2(const gsl_vector *X, void *params)
1298{
1300 astro_data** s1 = p->group1;
1301 astro_data** s2 = p->group2;
1302
1303 cpl_size nchannel = s1[0]->nchannel;
1304 cpl_size nwave = s1[0]->nwave;
1305
1306 double ra = gsl_vector_get(X, 0);
1307 double dec = gsl_vector_get(X, 1);
1308
1309 gsl_matrix_complex *visref_s1 = gravi_astrometry_calculate_visref_swap(s1, p->n1, ra, dec);
1310 gsl_matrix_complex *visref_s2 = gravi_astrometry_calculate_visref_swap(s2, p->n2, ra, dec);
1311
1312 double chi2 = 0.0;
1313 for (int i = 0; i < nchannel; i++) {
1314 double chi2_channel = 0;
1315 for (int k = 0; k < nwave; k++) {
1316 gsl_complex visref_swap = gsl_complex_mul(
1317 gsl_matrix_complex_get(visref_s1, i, k),
1318 gsl_complex_conjugate(gsl_matrix_complex_get(visref_s2, i, k))
1319 );
1320 visref_swap = gsl_complex_sqrt(visref_swap);
1321 chi2_channel += GSL_IMAG(visref_swap) * GSL_IMAG(visref_swap);
1322 }
1323 chi2 += chi2_channel;
1324 }
1325
1326 FREE(gsl_matrix_complex_free, visref_s1);
1327 FREE(gsl_matrix_complex_free, visref_s2);
1328 return chi2;
1329}
1330
1332 gravi_astrometry_model_params *params, const gsl_vector *ra_grid, const gsl_vector *dec_grid, gsl_vector *ra_dec_out, gsl_matrix *chi2_map)
1333{
1334 cpl_ensure_code(ra_grid, CPL_ERROR_NULL_INPUT);
1335 cpl_ensure_code(dec_grid, CPL_ERROR_NULL_INPUT);
1336 cpl_ensure_code(ra_dec_out, CPL_ERROR_NULL_INPUT);
1337
1338 cpl_size n_ra = ra_grid->size;
1339 cpl_size n_dec = dec_grid->size;
1340 double ra_dec[2] = {0};
1341 gsl_vector_const_view ra_dec_view = gsl_vector_const_view_array(ra_dec, 2);
1342 double best_ra, best_dec, best_chi2 = INFINITY;
1343
1344 for (int i = 0; i < n_ra; i++) {
1345 ra_dec[0] = gsl_vector_get(ra_grid, i);
1346 cpl_msg_debug(cpl_func, "calculating chi2 for ra=%f", ra_dec[0]);
1347 for (int j = 0; j < n_dec; j++) {
1348 ra_dec[1] = gsl_vector_get(dec_grid, j);
1349
1350 double chi2 = gravi_astrometry_calculate_chi2(
1351 &ra_dec_view.vector, params
1352 );
1353
1354 if (chi2_map)
1355 gsl_matrix_set(chi2_map, i, j, chi2);
1356
1357 if (chi2 < best_chi2) {
1358 best_ra = ra_dec[0];
1359 best_dec = ra_dec[1];
1360 best_chi2 = chi2;
1361 }
1362 }
1363 }
1364
1365 cpl_msg_debug(cpl_func, "Best chi2 is %f at [%f, %f]",
1366 best_chi2, best_ra, best_dec);
1367
1368 gsl_vector_set(ra_dec_out, 0, best_ra);
1369 gsl_vector_set(ra_dec_out, 1, best_dec);
1370 return CPL_ERROR_NONE;
1371}
1372
1373static cpl_error_code gravi_astrometry_minimise_chi2_descent(gravi_astrometry_model_params *params, double ra_guess, double dec_guess, gsl_vector *Xsolve)
1374{
1375 /* The minimisation has 2 parameters: ra, dec */
1376 const int dim = 2;
1377
1378 const int MAX_ITERATIONS = 1000;
1379 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex2rand;
1380 gsl_multimin_fminimizer *mini = gsl_multimin_fminimizer_alloc(T, dim);
1381 int status;
1382 double size;
1383
1384 gsl_multimin_function func = {
1386 .n = dim,
1387 .params = params
1388 };
1389
1390 double initial_guess[] = {ra_guess, dec_guess};
1391 gsl_vector_view guess_view = gsl_vector_view_array(initial_guess, dim);
1392 double step_size[] = {0.1, 0.1};
1393 gsl_vector_view step_view = gsl_vector_view_array(step_size, dim);
1394 gsl_multimin_fminimizer_set(mini, &func, &guess_view.vector, &step_view.vector);
1395
1396 int iterations = 0;
1397 do {
1398 iterations++;
1399 status = gsl_multimin_fminimizer_iterate(mini);
1400
1401 /* Error taking optimisation step */
1402 if (status)
1403 break;
1404
1405 size = gsl_multimin_fminimizer_size(mini);
1406 status = gsl_multimin_test_size(size, 1e-3);
1407 } while (status == GSL_CONTINUE && iterations < MAX_ITERATIONS);
1408
1409 if (status == GSL_SUCCESS)
1410 gsl_vector_memcpy(Xsolve, mini->x);
1411
1412 if (status == GSL_CONTINUE)
1413 cpl_msg_warning(cpl_func, "model-fitting did not converge");
1414 else if (status != GSL_SUCCESS)
1415 cpl_msg_error(cpl_func, "model-fitting failed with error code %d (%s)", status, gsl_strerror(status));
1416 else
1417 cpl_msg_debug(cpl_func, "converged with best-fit:\nra\tdec\tchi2\n%.3f\t%.3f\t%.4f", mini->x->data[0], mini->x->data[1], mini->fval);
1418
1419 FREE(gsl_multimin_fminimizer_free, mini);
1420 return status ? CPL_ERROR_ILLEGAL_OUTPUT : CPL_ERROR_NONE;
1421}
1422
1423cpl_error_code gravi_astrometry_reduce_swaps(astro_data **swap_data, cpl_size nswap, cpl_parameterlist *parlist)
1424{
1425 cpl_ensure_code(swap_data, CPL_ERROR_NULL_INPUT);
1426 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
1427
1428 astro_data **group1 = NULL, **group2 = NULL;
1429 cpl_size n1 = 0, n2 = 0;
1430
1431 gsl_vector *ra_grid = NULL, *dec_grid = NULL;
1432
1433 gsl_vector *ra_dec_fit = NULL;
1434 gsl_matrix *chi2_map = NULL;
1435
1436 /* Parameters */
1437 cpl_boolean use_swap_fiber_pos = cpl_parameter_get_bool(
1438 cpl_parameterlist_find(parlist, "gravity.astrometry.use-swap-fiber-pos"));
1439
1440 cpl_boolean go_fast = cpl_parameter_get_bool(
1441 cpl_parameterlist_find(parlist, "gravity.astrometry.average-over-dits"));
1442
1443 double ra_lim = cpl_parameter_get_double(
1444 cpl_parameterlist_find(parlist, "gravity.astrometry.ra-lim-swap"));
1445 int n_ra = cpl_parameter_get_int(
1446 cpl_parameterlist_find(parlist, "gravity.astrometry.nra-swap"));
1447
1448 double dec_lim = cpl_parameter_get_double(
1449 cpl_parameterlist_find(parlist, "gravity.astrometry.ra-lim-swap"));
1450 int n_dec = cpl_parameter_get_int(
1451 cpl_parameterlist_find(parlist, "gravity.astrometry.ndec-swap"));
1452
1453 double zoom = cpl_parameter_get_double(
1454 cpl_parameterlist_find(parlist, "gravity.astrometry.zoom-factor"));
1455
1456 CPLCHECK_CLEAN("Could not get parameters");
1457
1458 if (use_swap_fiber_pos) {
1459 cpl_msg_warning(cpl_func, "No astrometric solution is computed for swaps. Default to fiber position");
1460 for (int i = 0; i < nswap; i++) {
1461 astro_data *swap = swap_data[i];
1462 // swap->swap_astrometry_guess = gsl_vector_alloc(2);
1463 // gsl_vector_set(swap->swap_astrometry_guess, 0, swap->sobj_x);
1464 // gsl_vector_set(swap->swap_astrometry_guess, 1, swap->sobj_y);
1465 // swap->swap_astrometry_fit = gsl_vector_alloc(2);
1466 // gsl_vector_memcpy(swap->swap_astrometry_fit, swap->swap_astrometry_guess);
1467 swap->swap_astrometry_guess[0] = swap->swap_astrometry_fit[0] = swap->sobj_x;
1468 swap->swap_astrometry_guess[1] = swap->swap_astrometry_fit[1] = swap->sobj_y;
1469 }
1470 return CPL_ERROR_NONE;
1471 }
1472
1473 /* Divide the swaps into on/off-axis groups */
1474 for (int i = 0; i < nswap; i++) {
1475 if (swap_data[i]->swap)
1476 n1++;
1477 else
1478 n2++;
1479 }
1480
1481 if (n1 == 0)
1482 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT, "No files with SWAP=YES");
1483 else if (n2 == 0)
1484 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT, "No files with SWAP=NO");
1485 CPLCHECK_CLEAN("Could not partition swaps");
1486
1487 group1 = cpl_malloc(n1 * sizeof(astro_data*));
1488 group2 = cpl_malloc(n2 * sizeof(astro_data*));
1489 cpl_size i1 = 0, i2 = 0;
1490
1491 for (int i = 0; i < nswap; i++) {
1492 if (swap_data[i]->swap) {
1493 group1[i1++] = go_fast ? gravi_astrometry_average_over_dits(swap_data[i]) : swap_data[i];
1494 } else {
1495 group2[i2++] = go_fast ? gravi_astrometry_average_over_dits(swap_data[i]) : swap_data[i];
1496 }
1497 }
1498
1500 .group1 = group1, .group2 = group2, .n1 = n1, .n2 = n2
1501 };
1502
1503 /* Determine RA/dec grid */
1504 int field;
1505 if (strstr(swap_data[0]->insname, "U1234"))
1506 field = 60;
1507 else
1508 field = 240;
1509
1510 if (ra_lim < 0)
1511 ra_lim = field / 2;
1512 if (dec_lim < 0)
1513 dec_lim = field / 2;
1514
1515 double avg_sx = 0.0, avg_sy = 0.0;
1516 for (int i = 0; i < nswap; i++) {
1517 avg_sx += swap_data[i]->sobj_x * (swap_data[i]->swap ? -1 : 1);
1518 avg_sy += swap_data[i]->sobj_y * (swap_data[i]->swap ? -1 : 1);
1519 }
1520 avg_sx /= nswap;
1521 avg_sy /= nswap;
1522 cpl_msg_debug(cpl_func, "Centre: [%f, %f] mas", avg_sx, avg_sy);
1523
1524 double ra_fit, dec_fit;
1525 double ra_min = avg_sx - ra_lim, ra_max = avg_sx + ra_lim;
1526 double dec_min = avg_sy - dec_lim, dec_max = avg_sy + dec_lim;
1527 double d_ra = (ra_max - ra_min) / (n_ra - 1);
1528 double d_dec = (dec_max - dec_min) / (n_dec - 1);
1529
1530 ra_grid = gsl_vector_alloc(n_ra);
1531 dec_grid = gsl_vector_alloc(n_dec);
1532 for (int i = 0; i < n_ra; i++)
1533 gsl_vector_set(ra_grid, i, ra_min + i * d_ra);
1534 for (int j = 0; j < n_dec; j++)
1535 gsl_vector_set(dec_grid, j, dec_min + j * d_dec);
1536
1537 ra_dec_fit = gsl_vector_alloc(2);
1538 chi2_map = gsl_matrix_alloc(n_ra, n_dec);
1539
1540 if (zoom > 1.0) {
1541 cpl_msg_debug(cpl_func, "RA grid: [%f, %f] with %d points", ra_min, ra_max, n_ra);
1542 cpl_msg_debug(cpl_func, "dec grid: [%f, %f] with %d points", dec_min, dec_max, n_dec);
1543
1544 gravi_astrometry_minimise_chi2_grid(&params, ra_grid, dec_grid, ra_dec_fit, chi2_map);
1545 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1546 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1547 cpl_msg_debug(cpl_func, "big chi2 map minimsed for ra=%f, dec=%f", ra_fit, dec_fit);
1548
1549 // FILE *f = fopen("chi2_map_initial.dat", "w");
1550 // gsl_matrix_fprintf(f, chi2_map, "%g");
1551 // fclose(f);
1552
1553 ra_lim = (ra_max - ra_min) / (4 * zoom);
1554 dec_lim = (dec_max - dec_min) / (4 * zoom);
1555
1556 ra_min = ra_fit - ra_lim;
1557 ra_max = ra_fit + ra_lim;
1558 dec_min = dec_fit - dec_lim;
1559 dec_max = dec_fit + dec_lim;
1560
1561 d_ra = (ra_max - ra_min) / (n_ra - 1);
1562 d_dec = (dec_max - dec_min) / (n_dec - 1);
1563
1564 for (int i = 0; i < n_ra; i++)
1565 gsl_vector_set(ra_grid, i, ra_min + i * d_ra);
1566 for (int j = 0; j < n_dec; j++)
1567 gsl_vector_set(dec_grid, j, dec_min + j * d_dec);
1568 }
1569
1570 cpl_msg_debug(cpl_func, "RA grid: [%f, %f] with %d points", ra_min, ra_max, n_ra);
1571 cpl_msg_debug(cpl_func, "dec grid: [%f, %f] with %d points", dec_min, dec_max, n_dec);
1572
1573 gravi_astrometry_minimise_chi2_grid(&params, ra_grid, dec_grid, ra_dec_fit, chi2_map);
1574 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1575 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1576 cpl_msg_debug(cpl_func, "zoomed chi2 map minimised for ra=%f, dec=%f", ra_fit, dec_fit);
1577
1578 // FILE *f = fopen("chi2_map.dat", "w");
1579 // gsl_matrix_fprintf(f, chi2_map, "%g");
1580 // fclose(f);
1581
1582 for (int i = 0; i < nswap; i++) {
1583 // swap_data[i]->swap_astrometry_guess = gsl_vector_alloc(2);
1584 // gsl_vector_memcpy(swap_data[i]->swap_astrometry_guess, ra_dec_fit);
1585 swap_data[i]->swap_astrometry_guess[0] = ra_fit;
1586 swap_data[i]->swap_astrometry_guess[1] = dec_fit;
1587 }
1588
1589 gravi_astrometry_minimise_chi2_descent(&params, ra_fit, dec_fit, ra_dec_fit);
1590 CPLCHECK_CLEAN("Could not minimise swap astrometry");
1591
1592 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1593 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1594 cpl_msg_debug(cpl_func, "final astrometric solution after gradient descent is ra=%f, dec=%f", ra_fit, dec_fit);
1595
1596 for (int i = 0; i < nswap; i++) {
1597 // swap_data[i]->swap_astrometry_fit = gsl_vector_alloc(2);
1598 // gsl_vector_memcpy(swap_data[i]->swap_astrometry_fit, ra_dec_fit);
1599 swap_data[i]->swap_astrometry_fit[0] = ra_fit;
1600 swap_data[i]->swap_astrometry_fit[1] = dec_fit;
1601 }
1602
1603 /* Compute the swap phaseref, which is common to all the swap data */
1604 cpl_msg_debug(cpl_func, "Create the swap phase reference");
1605 gsl_matrix_complex *swap_ref = gravi_astrometry_create_swap_reference(swap_data, nswap);
1606 CPLCHECK_MSG("Could not extract swap phase reference");
1607 swap_data[0]->phase_ref_astro = swap_ref;
1608 for (int i = 1; i < nswap; i++)
1609 swap_data[i]->phase_ref_astro = gsl_matrix_complex_alloc_from_matrix(swap_ref, 0, 0, swap_ref->size1, swap_ref->size2);
1610
1611cleanup:
1612 if (go_fast) {
1613 FREELOOP(gravi_astrometry_delete, group1, n1);
1614 FREELOOP(gravi_astrometry_delete, group2, n2);
1615 } else {
1616 FREE(cpl_free, group1);
1617 FREE(cpl_free, group2);
1618 }
1619
1620 FREE(gsl_matrix_free, chi2_map);
1621 FREE(gsl_vector_free, ra_grid);
1622 FREE(gsl_vector_free, dec_grid);
1623 FREE(gsl_vector_free, ra_dec_fit);
1624
1625 return CPL_ERROR_NONE;
1626}
1627
1634static gsl_matrix_complex *gravi_astrometry_create_swap_reference(astro_data **swap_data, cpl_size nswap)
1635{
1636 cpl_ensure(swap_data, CPL_ERROR_NULL_INPUT, NULL);
1637
1638 int nchannel = swap_data[0]->nchannel;
1639 int nwave = swap_data[0]->nwave;
1640
1641 gsl_matrix_complex *phase_ref_s1 = gsl_matrix_complex_calloc(nchannel, nwave);
1642 gsl_matrix_complex *phase_ref_s2 = gsl_matrix_complex_calloc(nchannel, nwave);
1643 gsl_matrix_complex *phase_ref_avg = gsl_matrix_complex_calloc(nchannel, nwave);
1644
1645 for (int n = 0; n < nswap; n++) {
1646 astro_data *soi = swap_data[n];
1647 cpl_ensure(
1648 soi->nchannel == nchannel &&
1649 soi->nwave == nwave,
1650 CPL_ERROR_INCOMPATIBLE_INPUT,
1651 NULL
1652 );
1653
1654 /* First shift swap visibilities to zero OPD using the fitted ra and dec */
1655 // double swap_ra = gsl_vector_get(soi->swap_astrometry_fit, 0) * (soi->swap ? -1 : 1);
1656 // double swap_dec = gsl_vector_get(soi->swap_astrometry_fit, 1) * (soi->swap ? -1 : 1);
1657 double swap_ra = soi->swap_astrometry_fit[0] * (soi->swap ? -1 : 1);
1658 double swap_dec =soi->swap_astrometry_fit[1] * (soi->swap ? -1 : 1);
1659 gravi_astrometry_recentre_phase(soi, swap_ra, swap_dec);
1660
1661 /* Now take mean of vis and extract phase, for each position of the swap */
1662 gsl_matrix_complex *vis_ref_avg = average_matrix_complex_over_dits(soi->vis_ref, soi->ndit, soi->nchannel, soi->nwave, NULL);
1663 if (soi->swap)
1664 gsl_matrix_complex_add(phase_ref_s1, vis_ref_avg);
1665 else
1666 gsl_matrix_complex_add(phase_ref_s2, vis_ref_avg);
1667 FREE(gsl_matrix_complex_free, vis_ref_avg);
1668
1669 }
1670
1671 gsl_matrix_complex_memcpy(phase_ref_avg, phase_ref_s1);
1672 gsl_matrix_complex_add(phase_ref_avg, phase_ref_s2);
1673 gsl_matrix_complex_scale(phase_ref_avg, gsl_complex_rect(0.5, 0));
1674
1675 FREE(gsl_matrix_complex_free, phase_ref_s1);
1676 FREE(gsl_matrix_complex_free, phase_ref_s2);
1677
1678 return phase_ref_avg;
1679}
static cpl_error_code gravi_astrometry_mul_visibilities(astro_data *self, const gsl_vector *factor)
Scale VISDATA, VISREF, errors by vector factor, broadcasting over wavelength axis.
void gravi_astrometry_delete(astro_data *self)
cpl_error_code gravi_astrometry_create_phase_reference(astro_data *self, astro_data **phase_refs, cpl_size nphase, astro_data **swaps, cpl_size nswap, cpl_parameterlist *parlist)
Compute the final astrometric phase reference.
cpl_table * gravi_astrometry_get_phase_reference(astro_data *self)
double gravi_astrometry_get_mean_ftflux(astro_data *self)
static gsl_matrix * load_matrix(cpl_table *table, const char *name) CPL_ATTR_ALLOC
Load data from table into GSL matrix. Each table row is copied into the corresponding matrix row.
static gsl_vector * load_vector(cpl_table *table, const char *name) CPL_ATTR_ALLOC
Load data from table into GSL vector.
static gsl_matrix_complex * load_matrix_complex(cpl_table *table, const char *name) CPL_ATTR_ALLOC
Load data from table into GSL matrix. Each table row is copied into the corresponding matrix row.
static cpl_error_code gravi_astrometry_minimise_chi2_descent(gravi_astrometry_model_params *params, double ra_guess, double dec_guess, gsl_vector *Xsolve)
#define TWOPI
static gsl_matrix_complex * gravi_astrometry_calculate_visref_swap(astro_data **data, cpl_size ndata, double ra, double dec) CPL_ATTR_ALLOC
static int _gsl_vector_int_sum(const gsl_vector_int *a)
#define MAS_TO_RAD
cpl_error_code gravi_astrometry_dump(astro_data *self, FILE *handle)
_find_closest_mode
@ FIND_MODE_AFTER
@ FIND_MODE_BEFORE
@ FIND_MODE_NONE
struct _gravi_astrometry_model_params_ gravi_astrometry_model_params
static astro_data * gravi_astrometry_average_over_dits(astro_data *self) CPL_ATTR_ALLOC
Average input astro_data over DITs and return new flattened astro_data.
cpl_error_code gravi_astrometry_filter_ftflux(astro_data *self, double threshold)
Filter based on FT flux threshold and normalise.
static cpl_error_code gravi_astrometry_minimise_chi2_grid(gravi_astrometry_model_params *params, const gsl_vector *ra_grid, const gsl_vector *dec_grid, gsl_vector *ra_dec_out, gsl_matrix *chi2_map)
static double gravi_astrometry_calculate_chi2(const gsl_vector *X, void *params)
Evaluate chi^2 statistic.
static cpl_error_code gravi_astrometry_recentre_phase(astro_data *self, double xvalue, double yvalue)
Recentre visRef and visRefError on the given position (in mas) by shifting the phases.
cpl_error_code gravi_astrometry_normalise_to_ft(astro_data *self)
Normalise visibilities to average FT flux.
static gsl_matrix_complex * average_matrix_complex_over_dits(gsl_matrix_complex *m, int ndit, int nchannel, int nwave, gsl_matrix_int *flag) CPL_ATTR_ALLOC
static gsl_matrix * average_matrix_over_dits(gsl_matrix *m, int ndit, int nchannel, int nwave, gsl_matrix_int *flag) CPL_ATTR_ALLOC
static cpl_size gravi_astrometry_find_closest_mjd(astro_data *self, astro_data **others, cpl_size n_other, find_closest_mode mode)
Given observation and list of reference observations, find closest in time from the list.
enum _find_closest_mode find_closest_mode
static gsl_matrix_complex * gravi_astrometry_create_swap_reference(astro_data **swap_data, cpl_size nswap) CPL_ATTR_ALLOC
Centre swap data on zero OPD and extract average phase.
astro_data * gravi_astrometry_load(gravi_data *data)
Load data for astrometry from a gravi_data.
static gsl_vector * average_vector_over_dits(gsl_vector *v, int ndit, int nchannel) CPL_ATTR_ALLOC
cpl_error_code gravi_astrometry_reduce_swaps(astro_data **swap_data, cpl_size nswap, cpl_parameterlist *parlist)
static cpl_error_code gravi_astrometry_add_phase(astro_data *self, const gsl_matrix *phase)
Add phase to visRef, updating errors accordingly.
static cpl_error_code gravi_astrometry_scale_visibilities(astro_data *self, double factor)
Scale VISDATA, VISREF, errors by scalar factor.
typedefCPL_BEGIN_DECLS struct _gravi_data_ gravi_data
Definition: gravi_data.h:39
#define gravi_data_get_header(data)
Definition: gravi_data.h:75
#define gravi_data_get_oi_wave(data, type, pol, npol)
Definition: gravi_data.h:45
#define gravi_data_get_oi_vis(data, type, pol, npol)
Definition: gravi_data.h:46
#define gravi_data_get_oi_wave_plist(data, type, pol, npol)
Definition: gravi_data.h:69
cpl_msg_debug(cpl_func, "Spectra has <50 pixels -> don't flat")
#define GRAVI_SC
Definition: gravi_pfits.h:165
#define GRAVI_FT
Definition: gravi_pfits.h:166
#define CPLCHECK_INT(msg)
Definition: gravi_utils.h:51
#define CPLCHECK_CLEAN(msg)
Definition: gravi_utils.h:54
#define FREE(function, variable)
Definition: gravi_utils.h:69
#define CPLCHECK_NUL(msg)
Definition: gravi_utils.h:48
#define CPLCHECK_MSG(msg)
Definition: gravi_utils.h:45
#define FREELOOP(function, variable, n)
Definition: gravi_utils.h:72
gsl_matrix_complex * phase_ref_astro
gsl_vector ** vis_ref_cov
gsl_matrix * vcoord
gsl_matrix_complex * visdata
gsl_matrix * ucoord
gsl_matrix * phase_ref
gsl_vector_complex ** vis_ref_pcov
gsl_matrix * amp_ref_astro
double swap_astrometry_fit[2]
double swap_astrometry_guess[2]
gsl_matrix_int * flag
gsl_matrix * opd_disp
gsl_matrix_complex * vis_ref
gsl_matrix_complex * viserr
gsl_matrix_complex * visdata_ft
gsl_matrix * phase_met_telfc