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