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