38#include <gsl/gsl_matrix.h>
39#include <gsl/gsl_vector.h>
40#include <gsl/gsl_blas.h>
41#include <gsl/gsl_complex.h>
42#include <gsl/gsl_complex_math.h>
43#include <gsl/gsl_statistics.h>
44#include <gsl/gsl_multimin.h>
46#define PI 3.14159265359
47#define TWOPI 6.28318530718
48#define MAS_TO_RAD (2.0 * PI / (3600 * 360 * 1000))
95static gsl_vector *
load_vector(cpl_table *table,
const char *name) CPL_ATTR_ALLOC;
96static gsl_matrix *
load_matrix(cpl_table *table,
const char *name) CPL_ATTR_ALLOC;
97static gsl_matrix_complex *
load_matrix_complex(cpl_table *table,
const char *name) CPL_ATTR_ALLOC;
106static gsl_matrix *
average_matrix_over_dits(gsl_matrix *m,
int ndit,
int nchannel,
int nwave, gsl_matrix_int *flag) CPL_ATTR_ALLOC;
128 const size_t N = a->size;
129 const size_t stride = a->stride;
133 for (i = 0; i < N; i++)
134 sum += a->data[i * stride];
143 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
144 cpl_ensure(name != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
146 cpl_size depth = cpl_table_get_column_depth(table, name);
148 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
150 }
else if (depth > 0) {
151 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of scalar type", name);
155 cpl_size nelem = cpl_table_get_nrow(table);
156 double * cpl_vec = cpl_table_get_data_double(table, name);
158 gsl_vector *gsl_vec = gsl_vector_alloc(nelem);
159 memcpy(gsl_vec->data, cpl_vec, nelem *
sizeof(
double));
169 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
170 cpl_ensure(name != NULL, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
172 cpl_size depth = cpl_table_get_column_depth(table, name);
174 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
176 }
else if (depth == 0) {
177 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of 2D array type", name);
181 cpl_size nrow = cpl_table_get_nrow(table);
182 cpl_size ncol = cpl_table_get_column_dimension(table, name, 0);
183 const cpl_array ** cpl_matrix = cpl_table_get_data_array_const(table, name);
185 gsl_matrix *gsl_matrix = gsl_matrix_alloc(nrow, ncol);
186 for (
int i = 0; i < nrow; i++) {
187 gsl_vector_const_view row_view = gsl_vector_const_view_array(
188 cpl_array_get_data_double_const(cpl_matrix[i]), ncol);
189 gsl_matrix_set_row(gsl_matrix, i, &row_view.vector);
200 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
201 cpl_ensure(name != NULL, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
203 cpl_size depth = cpl_table_get_column_depth(table, name);
205 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
207 }
else if (depth == 0) {
208 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of 2D array type", name);
212 cpl_size nrow = cpl_table_get_nrow(table);
213 cpl_size ncol = cpl_table_get_column_dimension(table, name, 0);
214 const cpl_array ** cpl_matrix = cpl_table_get_data_array_const(table, name);
216 gsl_matrix_complex *gsl_matrix = gsl_matrix_complex_alloc(nrow, ncol);
217 for (
int i = 0; i < nrow; i++) {
218 for (
int j = 0; j < ncol; j++) {
219 double complex cpl_matrix_val = cpl_array_get_complex(cpl_matrix[i], j, NULL);
220 gsl_complex gsl_matrix_val = gsl_complex_rect(creal(cpl_matrix_val), cimag(cpl_matrix_val));
221 gsl_matrix_complex_set(gsl_matrix, i, j, gsl_matrix_val);
236 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
238 gsl_complex complex_factor = gsl_complex_rect(factor, 0);
239 gsl_matrix_complex_scale(self->
visdata, complex_factor);
240 gsl_matrix_complex_scale(self->
visdata_ft, complex_factor);
241 gsl_matrix_complex_scale(self->
vis_ref, complex_factor);
243 double factor2 = factor * factor;
246 gsl_vector_complex_scale(self->
vis_ref_pcov[i], gsl_complex_rect(factor2, 0));
249 return CPL_ERROR_NONE;
260 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
261 cpl_ensure_code(factor, CPL_ERROR_NULL_INPUT);
262 cpl_ensure_code((
size_t) (self->
ndit * self->
nchannel) == factor->size, CPL_ERROR_INCOMPATIBLE_INPUT);
264 gsl_matrix_complex *complex_factor = gsl_matrix_complex_alloc(self->
ndit * self->
nchannel, self->
nwave);
265 gsl_matrix_complex *complex_factor_ft = gsl_matrix_complex_alloc(self->
ndit * self->
nchannel, self->
nwave_ft);
267 gsl_complex cf = gsl_complex_rect(gsl_vector_get(factor, i), 0);
268 for (
int j = 0; j < self->
nwave; j++)
269 gsl_matrix_complex_set(complex_factor, i, j, cf);
270 for (
int j = 0; j < self->
nwave_ft; j++)
271 gsl_matrix_complex_set(complex_factor_ft, i, j, cf);
274 gsl_matrix_complex_mul_elements(self->
visdata, complex_factor);
275 gsl_matrix_complex_mul_elements(self->
vis_ref, complex_factor);
276 gsl_matrix_complex_mul_elements(self->
visdata_ft, complex_factor_ft);
279 double factor2 = gsl_vector_get(factor, i);
282 gsl_vector_complex_scale(self->
vis_ref_pcov[i], gsl_complex_rect(factor2, 0));
285 FREE(gsl_matrix_complex_free, complex_factor);
286 FREE(gsl_matrix_complex_free, complex_factor_ft);
288 return CPL_ERROR_NONE;
296 gsl_matrix_complex *phi = gsl_matrix_complex_alloc(self->
ndit * self->
nchannel, self->
nwave);
298 for (
int j = 0; j < self->
nwave; j++) {
299 gsl_complex phi_val = gsl_complex_polar(1.0, gsl_matrix_get(phase, i, j));
300 gsl_matrix_complex_set(phi, i, j, phi_val);
304 gsl_matrix_complex_mul_elements(self->
vis_ref, phi);
306 gsl_vector_complex *cov = gsl_vector_complex_alloc(self->
nwave);
307 gsl_vector_complex *pcov = gsl_vector_complex_alloc(self->
nwave);
308 gsl_vector_complex *tmp = gsl_vector_complex_alloc(self->
nwave);
311 gsl_vector_complex_set_zero(cov);
312 gsl_vector_complex_set_zero(pcov);
313 gsl_vector_complex_set_zero(tmp);
316 gsl_vector_complex_const_view phi_row = gsl_matrix_complex_const_row(phi, i);
317 gsl_vector_complex_memcpy(tmp, &phi_row.vector);
318 gsl_vector_view phi_imag = gsl_vector_complex_imag(tmp);
319 gsl_vector_scale(&phi_imag.vector, -1);
322 gsl_vector_view cov_real = gsl_vector_complex_real(cov);
323 gsl_vector_memcpy(&cov_real.vector, self->
vis_ref_cov[i]);
326 gsl_vector_complex_mul(tmp, cov);
327 gsl_vector_complex_mul(tmp, &phi_row.vector);
328 cov_real = gsl_vector_complex_real(tmp);
329 gsl_vector_memcpy(self->
vis_ref_cov[i], &cov_real.vector);
333 gsl_vector_complex_mul(pcov, &phi_row.vector);
334 gsl_vector_complex_mul(pcov, &phi_row.vector);
338 FREE(gsl_vector_complex_free, cov);
339 FREE(gsl_vector_complex_free, pcov);
340 FREE(gsl_vector_complex_free, tmp);
341 FREE(gsl_matrix_complex_free, phi);
343 return CPL_ERROR_NONE;
352 gsl_matrix_memcpy(this_u, self->
ucoord);
353 gsl_matrix_scale(this_u, xvalue);
356 gsl_matrix_memcpy(this_v, self->
vcoord);
357 gsl_matrix_scale(this_v, yvalue);
359 gsl_matrix_add(this_u, this_v);
364 FREE(gsl_matrix_free, this_u);
365 FREE(gsl_matrix_free, this_v);
367 return CPL_ERROR_NONE;
375 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
379 gsl_vector_complex_const_view visft_row = gsl_matrix_complex_const_row(self->
visdata_ft, i);
380 double row_abs_sum = 0.0;
381 for (
int j = 0; j < self->
nwave_ft; j++)
382 row_abs_sum += gsl_complex_abs(gsl_vector_complex_get(&visft_row.vector, j));
385 if (row_abs_sum < threshold) {
387 for (
int k = 0; k < self->
nwave; k++) {
388 gsl_matrix_int_set(self->
flag, i, k, 1);
389 gsl_matrix_complex_set(self->
visdata, i, k, GSL_COMPLEX_ZERO);
414 cpl_msg_debug(cpl_func,
"Flagging %d points below FT threshold %.3e", npoints, threshold);
415 self->
nflag += npoints;
417 return CPL_ERROR_NONE;
425 gsl_vector *ftflux_scale = gsl_vector_alloc(self->
ndit * self->
nchannel);
428 for (
int j = 0; j < self->
nwave_ft; j++)
429 scale += gsl_complex_abs(gsl_matrix_complex_get(self->
visdata_ft, i, j));
430 gsl_vector_set(ftflux_scale, i, self->
nwave / scale);
434 FREE(gsl_vector_free, ftflux_scale);
435 return CPL_ERROR_NONE;
443 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
449 const char *filename = cpl_propertylist_get_string(hdr,
"PIPEFILE");
450 self->
filename = cpl_strdup(filename);
452 const char *insname = cpl_propertylist_get_string(hdr,
"TELESCOP");
453 self->
insname = cpl_strdup(insname);
456 if (cpl_propertylist_has(hdr,
"ESO TPL NDIT OBJECT"))
457 self->
ndit = cpl_propertylist_get_int(hdr,
"ESO TPL NDIT OBJECT");
459 self->
ndit = cpl_propertylist_get_int(hdr,
"ESO TPL NDIT SKY");
460 self->
dit = cpl_propertylist_get_double(hdr,
"ESO DET2 SEQ1 DIT");
462 self->
swap = (!strcmp(cpl_propertylist_get_string(hdr,
"ESO INS SOBJ SWAP"),
"YES"));
463 cpl_msg_debug(cpl_func,
"SWAP keyword is: %s", cpl_propertylist_get_string(hdr,
"ESO INS SOBJ SWAP"));
467 self->
sobj_x = cpl_propertylist_get_double(hdr,
"ESO INS SOBJ X");
468 self->
sobj_y = cpl_propertylist_get_double(hdr,
"ESO INS SOBJ Y");
471 self->
mjd = cpl_propertylist_get_double(hdr,
"MJD-OBS");
479 self->
nwave = cpl_propertylist_get_int(wave_plist,
"NWAVE");
482 self->
nwave_ft = cpl_propertylist_get_int(wave_plist,
"NWAVE");
486 cpl_table_cast_column(wave_table,
"EFF_WAVE",
"EFF_WAVE", CPL_TYPE_DOUBLE);
498 for (
int j = 0; j < self->
nwave; j++) {
499 gsl_vector_view col_view;
500 col_view = gsl_matrix_column(self->
ucoord, j);
501 gsl_vector_memcpy(&col_view.vector, self->
u);
502 col_view = gsl_matrix_column(self->
vcoord, j);
503 gsl_vector_memcpy(&col_view.vector, self->
v);
507 gsl_vector_view row_view;
508 row_view = gsl_matrix_row(self->
ucoord, i);
509 gsl_vector_div(&row_view.vector, self->
wave);
510 row_view = gsl_matrix_row(self->
vcoord, i);
511 gsl_vector_div(&row_view.vector, self->
wave);
523 const cpl_array ** cpl_flag = cpl_table_get_data_array_const(oivis_table,
"FLAG");
524 const int* cpl_rej = NULL;
525 if (cpl_table_has_column(oivis_table,
"REJECTION_FLAG"))
526 cpl_rej = cpl_table_get_data_int_const(oivis_table,
"REJECTION_FLAG");
532 int rej_val = (cpl_rej ? cpl_rej[i] : 0) > 0;
533 for (
int j = 0; j < self->
nwave; j++) {
534 int flag_val = cpl_array_get_int(cpl_flag[i], j, NULL);
536 gsl_matrix_int_set(self->
flag, i, j, flag_val);
537 self->
nflag += flag_val;
544 for (
int j = 0; j < self->
nwave; j++) {
545 if (gsl_matrix_int_get(self->
flag, i, j)) {
547 gsl_matrix_complex_set(self->
visdata, i, j, GSL_COMPLEX_ZERO);
551 cpl_msg_debug(cpl_func,
"Cleaned %d flagged values in VISDATA", ncleaned);
570 for (
int j = 0; j < self->
nwave; j++) {
571 gsl_complex vis_ref_val = gsl_matrix_complex_get(self->
visdata, i, j);
572 gsl_complex pha = gsl_complex_polar(1.0, gsl_matrix_get(self->
phase_ref, i, j));
573 vis_ref_val = gsl_complex_mul(vis_ref_val, pha);
574 gsl_matrix_complex_set(self->
vis_ref, i, j, vis_ref_val);
582 gsl_vector *vreal = gsl_vector_alloc(self->
nwave);
583 gsl_vector *vimag = gsl_vector_alloc(self->
nwave);
586 gsl_vector_complex_const_view err_vals = gsl_matrix_complex_const_row(self->
viserr, i);
587 gsl_vector_const_view real_vals = gsl_vector_complex_const_real(&err_vals.vector);
588 gsl_vector_const_view imag_vals = gsl_vector_complex_const_imag(&err_vals.vector);
590 gsl_vector_memcpy(vreal, &real_vals.vector);
591 gsl_vector_mul(vreal, vreal);
593 gsl_vector_memcpy(vimag, &imag_vals.vector);
594 gsl_vector_mul(vimag, vimag);
597 gsl_vector *cov_diag = gsl_vector_alloc(self->
nwave);
598 gsl_vector_memcpy(cov_diag, vreal);
599 gsl_vector_add(cov_diag, vimag);
603 gsl_vector_complex *pcov_diag = gsl_vector_complex_calloc(self->
nwave);
604 gsl_vector_view pcov_diag_real = gsl_vector_complex_real(pcov_diag);
605 gsl_vector_memcpy(&pcov_diag_real.vector, vreal);
606 gsl_vector_sub(&pcov_diag_real.vector, vimag);
608 for (
int j = 0; j < self->
nwave; j++) {
610 gsl_complex pcov_val = gsl_complex_polar(1.0, 2 * gsl_matrix_get(self->
phase_ref, i, j));
611 pcov_val = gsl_complex_mul(gsl_vector_complex_get(pcov_diag, j), pcov_val);
612 gsl_vector_complex_set(pcov_diag, j, pcov_val);
616 FREE(gsl_vector_free, vreal);
617 FREE(gsl_vector_free, vimag);
623 gsl_matrix *wave_opd_disp = gsl_matrix_alloc(self->
ndit * self->
nchannel, self->
nwave);
624 gsl_matrix_memcpy(wave_opd_disp, self->
opd_disp);
626 gsl_vector *wavenumber = gsl_vector_alloc(self->
nwave);
627 gsl_vector_set_all(wavenumber,
TWOPI);
628 gsl_vector_div(wavenumber, self->
wave);
632 gsl_vector_view wave_opd_dist_row = gsl_matrix_row(wave_opd_disp, i);
633 gsl_vector_mul(&wave_opd_dist_row.vector, wavenumber);
638 gsl_matrix_scale(wave_opd_disp, -1);
641 FREE(gsl_vector_free, wavenumber);
642 FREE(gsl_matrix_free, wave_opd_disp);
649 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
650 cpl_ensure_code(handle, CPL_ERROR_NULL_INPUT);
652 fprintf(handle,
"%s\n", self->
filename);
653 fprintf(handle,
"%s\n", self->
insname);
657 gsl_matrix_fprintf(handle, self->
ucoord,
"%f");
658 gsl_matrix_fprintf(handle, self->
vcoord,
"%f");
659 gsl_matrix_int_fprintf(handle, self->
flag,
"%d");
660 gsl_matrix_complex_fprintf(handle, self->
visdata,
"%f");
661 gsl_matrix_complex_fprintf(handle, self->
vis_ref,
"%f");
663 return CPL_ERROR_NONE;
671 FREE(gsl_vector_free, self->
u);
672 FREE(gsl_vector_free, self->
v);
681 FREE(gsl_matrix_int_free, self->
flag);
697 cpl_ensure(self, CPL_ERROR_NULL_INPUT, -1);
701 for (
int j = 0; j < self->
nwave_ft; j++) {
702 abs_val += gsl_complex_abs(gsl_matrix_complex_get(self->
visdata_ft, i, j));
721 cpl_ensure(self, CPL_ERROR_NULL_INPUT, -1);
722 cpl_ensure(others, CPL_ERROR_NULL_INPUT, -1);
724 cpl_boolean any_before = FALSE, any_after = FALSE;
725 for (
int i = 0; i < n_other; i++) {
726 double delta = self->
mjd - others[i]->
mjd;
728 any_before = CPL_TRUE;
730 any_after = CPL_TRUE;
734 double min_delta = INFINITY;
736 for (
int i = 0; i < n_other; i++) {
737 double delta = self->
mjd - others[i]->
mjd;
740 if (any_before && delta >= 0 && delta < min_delta) {
746 if (any_after && delta <= 0 && -delta < min_delta) {
752 if (delta <= 0 && fabs(delta) < min_delta) {
774 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
775 cpl_ensure_code(phase_refs, CPL_ERROR_NULL_INPUT);
777 cpl_ensure_code(swaps, CPL_ERROR_NULL_INPUT);
778 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
780 gsl_matrix *amp_ref = NULL, *amp_ref_tmp = NULL;
781 gsl_matrix_complex *vis_ref = NULL, *vis_ref_tmp = NULL;
784 char *calib_strategy = cpl_strdup(cpl_parameter_get_string(
785 cpl_parameterlist_find(parlist,
"gravity.astrometry.calib-strategy")));
787 for(
int i = 0; calib_strategy[i]; i++)
788 calib_strategy[i] = toupper(calib_strategy[i]);
792 cpl_msg_debug(cpl_func,
"Calibrating amplitude reference with '%s' strategy", calib_strategy);
794 if (!strcmp(calib_strategy,
"NONE")) {
796 vis_ref = gsl_matrix_complex_alloc(self->
nchannel, self->
nwave);
797 gsl_matrix_set_all(amp_ref, 1.0);
798 gsl_matrix_complex_set_all(vis_ref, GSL_COMPLEX_ONE);
799 }
else if (!strcmp(calib_strategy,
"SELF")) {
803 for (
int j = 0; j < self->
nwave; j++) {
804 gsl_complex vis_ref_val = gsl_matrix_complex_get(self->
vis_ref, i, j);
805 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
809 FREE(gsl_matrix_free, amp_ref_tmp);
811 vis_ref = gsl_matrix_complex_alloc(self->
nchannel, self->
nwave);
812 for (
int i = 0; i < self->
nchannel; i++) {
813 for (
int j = 0; j < self->
nwave; j++) {
814 double amp_ref_val = gsl_matrix_get(amp_ref, i, j);
815 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_rect(amp_ref_val, 0));
820 gsl_vector_int *phase_indices;
822 if (!strcmp(calib_strategy,
"ALL")) {
823 nphase_used = nphase;
824 phase_indices = gsl_vector_int_alloc(nphase_used);
825 for (
int n = 0; n < nphase_used; n++)
826 gsl_vector_int_set(phase_indices, n, n);
827 }
else if (!strcmp(calib_strategy,
"NEAREST")) {
829 phase_indices = gsl_vector_int_alloc(nphase_used);
835 if (idx_before == -1 || idx_after == -1) {
836 cpl_error_set_message(cpl_func, CPL_ERROR_ACCESS_OUT_OF_RANGE,
"Cannot find time-bracketing phase references");
837 return CPL_ERROR_ACCESS_OUT_OF_RANGE;
839 gsl_vector_int_set(phase_indices, 0, idx_before);
840 gsl_vector_int_set(phase_indices, 1, idx_after);
842 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"Unknown calibration strategy %s", calib_strategy);
843 return CPL_ERROR_ILLEGAL_INPUT;
847 vis_ref = gsl_matrix_complex_calloc(self->
nchannel, self->
nwave);
850 for (
int n = 0; n < nphase_used; n++) {
851 cpl_size idx = gsl_vector_int_get(phase_indices, n);
858 for (
int i = 0; i < soi->
nchannel; i++) {
859 for (
int j = 0; j < soi->
nwave; j++) {
860 gsl_complex vis_ref_val = gsl_matrix_complex_get(vis_ref_tmp, i, j);
861 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
865 gsl_matrix_add(amp_ref, amp_ref_tmp);
866 gsl_matrix_complex_add(vis_ref, vis_ref_tmp);
867 FREE(gsl_matrix_complex_free, vis_ref_tmp);
870 gsl_matrix_scale(amp_ref, 1.0 / nphase_used);
871 gsl_matrix_complex_scale(vis_ref, gsl_complex_rect(1.0 / nphase_used, 0.0));
873 for (
int i = 0; i < self->
nchannel; i++) {
874 for (
int j = 0; j < self->
nwave; j++) {
875 double theta = gsl_complex_arg(gsl_matrix_complex_get(vis_ref, i, j));
876 double r = gsl_matrix_get(amp_ref, i, j);
877 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_polar(r, theta));
880 FREE(gsl_matrix_free, amp_ref_tmp);
881 FREE(gsl_vector_int_free, phase_indices);
889 cpl_msg_debug(cpl_func,
"Extract phase reference from swaps");
895 if (!strcmp(calib_strategy,
"SWAP")) {
896 cpl_size n_swapon = 0;
898 gsl_matrix_set_zero(amp_ref);
899 gsl_matrix_complex_set_zero(vis_ref);
901 for (
int n = 0; n < nswap; n++) {
909 for (
int j = 0; j < soi->
nwave; j++) {
910 gsl_complex vis_ref_val = gsl_matrix_complex_get(soi->
vis_ref, i, j);
911 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
915 gsl_matrix_add(amp_ref, amp_ref_avg);
917 FREE(gsl_matrix_free, amp_ref_tmp);
918 FREE(gsl_matrix_free, amp_ref_avg);
922 gsl_matrix_scale(amp_ref, 1.0 / n_swapon);
923 for (
int i = 0; i < self->
nchannel; i++) {
924 for (
int j = 0; j < self->
nwave; j++) {
926 double amp_ref_val = gsl_matrix_get(amp_ref, i, j);
927 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_rect(0.5 * amp_ref_val, 0));
934 for (
int i = 0; i < self->
nchannel; i++) {
935 for (
int j = 0; j < self->
nwave; j++) {
936 gsl_complex phi = gsl_matrix_complex_get(swap_ref, i, j);
937 double argphi = gsl_complex_arg(phi);
939 gsl_complex vis_ref_val = gsl_matrix_complex_get(vis_ref, i, j);
940 double vis_ref_abs = gsl_complex_abs(vis_ref_val);
942 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_polar(2 * vis_ref_abs, argphi));
947 cpl_free(calib_strategy);
948 return CPL_ERROR_NONE;
1022 cpl_table *table = cpl_table_new(self->
nwave);
1023 cpl_table_new_column_array(table,
"ASTRO_VISREF", CPL_TYPE_DOUBLE_COMPLEX, self->
nchannel);
1024 cpl_table_new_column_array(table,
"ASTRO_AMPREF", CPL_TYPE_DOUBLE, self->
nchannel);
1025 cpl_table_new_column_array(table,
"ASTRO_COV", CPL_TYPE_DOUBLE, self->
ndit * self->
nchannel);
1026 cpl_table_new_column_array(table,
"ASTRO_PCOV", CPL_TYPE_DOUBLE_COMPLEX, self->
ndit * self->
nchannel);
1028 cpl_array *tmp_arr = cpl_array_new(self->
nchannel, CPL_TYPE_DOUBLE_COMPLEX);
1029 for (
int j = 0; j < self->
nwave; j++) {
1030 for (
int i = 0; i < self->
nchannel; i++) {
1031 gsl_complex val = gsl_matrix_complex_get(self->
phase_ref_astro, i, j);
1032 cpl_array_set_complex(tmp_arr, i, GSL_REAL(val) + I * GSL_IMAG(val));
1034 cpl_table_set_array(table,
"ASTRO_VISREF", j, tmp_arr);
1037 FREE(cpl_array_delete, tmp_arr);
1039 tmp_arr = cpl_array_new(self->
nchannel, CPL_TYPE_DOUBLE);
1040 for (
int j = 0; j < self->
nwave; j++) {
1041 for (
int i = 0; i < self->
nchannel; i++) {
1043 cpl_array_set(tmp_arr, i, val);
1045 cpl_table_set_array(table,
"ASTRO_AMPREF", j, tmp_arr);
1048 FREE(cpl_array_delete, tmp_arr);
1051 tmp_arr = cpl_array_new(self->
ndit * self->
nchannel, CPL_TYPE_DOUBLE);
1052 for (
int j = 0; j < self->
nwave; j++) {
1053 for (
int i = 0; i < self->
ndit * self->
nchannel; i++) {
1054 double val = gsl_vector_get(self->
vis_ref_cov[i], j);
1055 cpl_array_set(tmp_arr, i, val);
1057 cpl_table_set_array(table,
"ASTRO_COV", j, tmp_arr);
1060 FREE(cpl_array_delete, tmp_arr);
1063 tmp_arr = cpl_array_new(self->
ndit * self->
nchannel, CPL_TYPE_DOUBLE_COMPLEX);
1064 for (
int j = 0; j < self->
nwave; j++) {
1065 for (
int i = 0; i < self->
ndit * self->
nchannel; i++) {
1066 gsl_complex val = gsl_vector_complex_get(self->
vis_ref_pcov[i], j);
1067 cpl_array_set_complex(tmp_arr, i, GSL_REAL(val) + I * GSL_IMAG(val));
1069 cpl_table_set_array(table,
"ASTRO_PCOV", j, tmp_arr);
1072 FREE(cpl_array_delete, tmp_arr);
1079 gsl_vector *vavg = gsl_vector_alloc(nchannel);
1081 for (
int i = 0; i < nchannel; i++) {
1082 double avg = gsl_stats_mean(v->data + i, nchannel, ndit);
1083 gsl_vector_set(vavg, i, avg);
1091 gsl_matrix *mavg = gsl_matrix_alloc(nchannel, nwave);
1093 for (
int i = 0; i < nchannel; i++) {
1094 for (
int j = 0; j < nwave; j++) {
1095 cpl_size start = i * nwave + j;
1096 gsl_vector_const_view dit_view = gsl_vector_const_view_array_with_stride(
1097 m->data + start, nwave * nchannel, ndit);
1098 gsl_vector_int_const_view flag_view = flag ?
1099 gsl_vector_int_const_view_array_with_stride(flag->data + start, nwave * nchannel, ndit) :
1100 (gsl_vector_int_const_view){};
1103 cpl_size Nvalid = 0;
1104 for (
int k = 0; k < ndit; k++) {
1105 if (!flag || gsl_vector_int_get(&flag_view.vector, k) == 0) {
1106 avg += gsl_vector_get(&dit_view.vector, k);
1112 gsl_matrix_set(mavg, i, j, avg);
1114 gsl_matrix_set(mavg, i, j, 0.0);
1124 gsl_matrix_complex *mavg = gsl_matrix_complex_alloc(nchannel, nwave);
1126 for (
int i = 0; i < nchannel; i++) {
1127 for (
int j = 0; j < nwave; j++) {
1128 cpl_size start = i * nwave + j;
1129 gsl_vector_complex_const_view dit_view = gsl_vector_complex_const_view_array_with_stride(
1130 m->data + 2 * start, nwave * nchannel, ndit);
1131 gsl_vector_int_const_view flag_view = flag ?
1132 gsl_vector_int_const_view_array_with_stride(flag->data + start, nwave * nchannel, ndit) :
1133 (gsl_vector_int_const_view){};
1135 gsl_complex avg = GSL_COMPLEX_ZERO;
1136 cpl_size Nvalid = 0;
1137 for (
int k = 0; k < ndit; k++) {
1138 if (!flag || gsl_vector_int_get(&flag_view.vector, k) == 0) {
1139 avg = gsl_complex_add(avg, gsl_vector_complex_get(&dit_view.vector, k));
1144 avg = gsl_complex_div_real(avg, Nvalid);
1145 gsl_matrix_complex_set(mavg, i, j, avg);
1147 gsl_matrix_complex_set(mavg, i, j, GSL_COMPLEX_ZERO);
1160 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
1181 flat->
wave = gsl_vector_alloc(flat->
nwave);
1182 gsl_vector_memcpy(flat->
wave, self->
wave);
1200 for (
int i = 0; i < self->
nchannel; i++) {
1204 for (
int j = 0; j < self->
ndit; j++) {
1205 cpl_size idx = j * self->
nchannel + i;
1210 gsl_vector_complex_scale(flat->
vis_ref_pcov[i], gsl_complex_rect(1.0 / (self->
ndit * self->
ndit), 0));
1215 for (
int i = 0; i < self->
nchannel; i++) {
1216 for (
int j = 0; j < self->
nwave; j++) {
1217 cpl_size start = i * self->
nwave + j;
1218 gsl_vector_int_const_view flag_view = gsl_vector_int_const_view_array_with_stride(
1222 gsl_matrix_int_set(flat->
flag, i, j, 1);
1238 gsl_matrix_complex *visref = gsl_matrix_complex_calloc(data[0]->nchannel, data[0]->nwave);
1239 gsl_matrix_complex *visref_ditsum = gsl_matrix_complex_alloc(data[0]->nchannel, data[0]->nwave);
1241 gsl_matrix_int *ngood = gsl_matrix_int_calloc(data[0]->nchannel, data[0]->nwave);
1242 gsl_matrix_int *ngood_ditsum = gsl_matrix_int_alloc(data[0]->nchannel, data[0]->nwave);
1244 for (
int n = 0; n < ndata; n++) {
1247 for (
int j = 0; j < self->
nchannel; j++) {
1248 for (
int k = 0; k < self->
nwave; k++) {
1250 gsl_complex visref_jk = GSL_COMPLEX_ZERO;
1252 for (
int i = 0; i < self->
ndit; i++) {
1253 cpl_size idx = i * self->
nchannel + j;
1254 double uv = gsl_matrix_get(self->
ucoord, idx, k) * ra + gsl_matrix_get(self->
vcoord, idx, k) * dec;
1257 gsl_complex phi = gsl_complex_polar(1.0, phase);
1259 phi = gsl_complex_conjugate(phi);
1260 gsl_complex visref_ijk = gsl_complex_mul(phi, gsl_matrix_complex_get(self->
vis_ref, idx, k));
1262 int flag = gsl_matrix_int_get(self->
flag, idx, k);
1265 visref_jk = gsl_complex_add(visref_jk, visref_ijk);
1268 gsl_matrix_complex_set(visref_ditsum, j, k, visref_jk);
1269 gsl_matrix_int_set(ngood_ditsum, j, k, ngood_jk);
1272 gsl_matrix_complex_add(visref, visref_ditsum);
1273 gsl_matrix_int_add(ngood, ngood_ditsum);
1276 for (
int j = 0; j < data[0]->
nchannel; j++) {
1277 for (
int k = 0; k < data[0]->
nwave; k++) {
1278 int ngood_jk = gsl_matrix_int_get(ngood, j, k);
1280 gsl_complex visref_jk = gsl_matrix_complex_get(visref, j, k);
1281 visref_jk = gsl_complex_div_real(visref_jk, ngood_jk);
1282 gsl_matrix_complex_set(visref, j, k, visref_jk);
1287 FREE(gsl_matrix_complex_free, visref_ditsum);
1303 cpl_size nchannel = s1[0]->
nchannel;
1304 cpl_size nwave = s1[0]->
nwave;
1306 double ra = gsl_vector_get(X, 0);
1307 double dec = gsl_vector_get(X, 1);
1313 for (
int i = 0; i < nchannel; i++) {
1314 double chi2_channel = 0;
1315 for (
int k = 0; k < nwave; k++) {
1316 gsl_complex visref_swap = gsl_complex_mul(
1317 gsl_matrix_complex_get(visref_s1, i, k),
1318 gsl_complex_conjugate(gsl_matrix_complex_get(visref_s2, i, k))
1320 visref_swap = gsl_complex_sqrt(visref_swap);
1321 chi2_channel += GSL_IMAG(visref_swap) * GSL_IMAG(visref_swap);
1323 chi2 += chi2_channel;
1326 FREE(gsl_matrix_complex_free, visref_s1);
1327 FREE(gsl_matrix_complex_free, visref_s2);
1334 cpl_ensure_code(ra_grid, CPL_ERROR_NULL_INPUT);
1335 cpl_ensure_code(dec_grid, CPL_ERROR_NULL_INPUT);
1336 cpl_ensure_code(ra_dec_out, CPL_ERROR_NULL_INPUT);
1338 cpl_size n_ra = ra_grid->size;
1339 cpl_size n_dec = dec_grid->size;
1340 double ra_dec[2] = {0};
1341 gsl_vector_const_view ra_dec_view = gsl_vector_const_view_array(ra_dec, 2);
1342 double best_ra, best_dec, best_chi2 = INFINITY;
1344 for (
int i = 0; i < n_ra; i++) {
1345 ra_dec[0] = gsl_vector_get(ra_grid, i);
1346 cpl_msg_debug(cpl_func,
"calculating chi2 for ra=%f", ra_dec[0]);
1347 for (
int j = 0; j < n_dec; j++) {
1348 ra_dec[1] = gsl_vector_get(dec_grid, j);
1351 &ra_dec_view.vector, params
1355 gsl_matrix_set(chi2_map, i, j, chi2);
1357 if (chi2 < best_chi2) {
1358 best_ra = ra_dec[0];
1359 best_dec = ra_dec[1];
1366 best_chi2, best_ra, best_dec);
1368 gsl_vector_set(ra_dec_out, 0, best_ra);
1369 gsl_vector_set(ra_dec_out, 1, best_dec);
1370 return CPL_ERROR_NONE;
1378 const int MAX_ITERATIONS = 1000;
1379 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex2rand;
1380 gsl_multimin_fminimizer *mini = gsl_multimin_fminimizer_alloc(T, dim);
1384 gsl_multimin_function func = {
1390 double initial_guess[] = {ra_guess, dec_guess};
1391 gsl_vector_view guess_view = gsl_vector_view_array(initial_guess, dim);
1392 double step_size[] = {0.1, 0.1};
1393 gsl_vector_view step_view = gsl_vector_view_array(step_size, dim);
1394 gsl_multimin_fminimizer_set(mini, &func, &guess_view.vector, &step_view.vector);
1399 status = gsl_multimin_fminimizer_iterate(mini);
1405 size = gsl_multimin_fminimizer_size(mini);
1406 status = gsl_multimin_test_size(size, 1e-3);
1407 }
while (status == GSL_CONTINUE && iterations < MAX_ITERATIONS);
1409 if (status == GSL_SUCCESS)
1410 gsl_vector_memcpy(Xsolve, mini->x);
1412 if (status == GSL_CONTINUE)
1413 cpl_msg_warning(cpl_func,
"model-fitting did not converge");
1414 else if (status != GSL_SUCCESS)
1415 cpl_msg_error(cpl_func,
"model-fitting failed with error code %d (%s)", status, gsl_strerror(status));
1417 cpl_msg_debug(cpl_func,
"converged with best-fit:\nra\tdec\tchi2\n%.3f\t%.3f\t%.4f", mini->x->data[0], mini->x->data[1], mini->fval);
1419 FREE(gsl_multimin_fminimizer_free, mini);
1420 return status ? CPL_ERROR_ILLEGAL_OUTPUT : CPL_ERROR_NONE;
1425 cpl_ensure_code(swap_data, CPL_ERROR_NULL_INPUT);
1426 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
1429 cpl_size n1 = 0, n2 = 0;
1431 gsl_vector *ra_grid = NULL, *dec_grid = NULL;
1433 gsl_vector *ra_dec_fit = NULL;
1434 gsl_matrix *chi2_map = NULL;
1437 cpl_boolean use_swap_fiber_pos = cpl_parameter_get_bool(
1438 cpl_parameterlist_find(parlist,
"gravity.astrometry.use-swap-fiber-pos"));
1440 cpl_boolean go_fast = cpl_parameter_get_bool(
1441 cpl_parameterlist_find(parlist,
"gravity.astrometry.average-over-dits"));
1443 double ra_lim = cpl_parameter_get_double(
1444 cpl_parameterlist_find(parlist,
"gravity.astrometry.ra-lim-swap"));
1445 int n_ra = cpl_parameter_get_int(
1446 cpl_parameterlist_find(parlist,
"gravity.astrometry.nra-swap"));
1448 double dec_lim = cpl_parameter_get_double(
1449 cpl_parameterlist_find(parlist,
"gravity.astrometry.ra-lim-swap"));
1450 int n_dec = cpl_parameter_get_int(
1451 cpl_parameterlist_find(parlist,
"gravity.astrometry.ndec-swap"));
1453 double zoom = cpl_parameter_get_double(
1454 cpl_parameterlist_find(parlist,
"gravity.astrometry.zoom-factor"));
1458 if (use_swap_fiber_pos) {
1459 cpl_msg_warning(cpl_func,
"No astrometric solution is computed for swaps. Default to fiber position");
1460 for (
int i = 0; i < nswap; i++) {
1470 return CPL_ERROR_NONE;
1474 for (
int i = 0; i < nswap; i++) {
1475 if (swap_data[i]->swap)
1482 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"No files with SWAP=YES");
1484 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"No files with SWAP=NO");
1487 group1 = cpl_malloc(n1 *
sizeof(
astro_data*));
1488 group2 = cpl_malloc(n2 *
sizeof(
astro_data*));
1489 cpl_size i1 = 0, i2 = 0;
1491 for (
int i = 0; i < nswap; i++) {
1492 if (swap_data[i]->swap) {
1500 .
group1 = group1, .group2 = group2, .n1 = n1, .n2 = n2
1505 if (strstr(swap_data[0]->insname,
"U1234"))
1513 dec_lim = field / 2;
1515 double avg_sx = 0.0, avg_sy = 0.0;
1516 for (
int i = 0; i < nswap; i++) {
1517 avg_sx += swap_data[i]->
sobj_x * (swap_data[i]->
swap ? -1 : 1);
1518 avg_sy += swap_data[i]->
sobj_y * (swap_data[i]->
swap ? -1 : 1);
1522 cpl_msg_debug(cpl_func,
"Centre: [%f, %f] mas", avg_sx, avg_sy);
1524 double ra_fit, dec_fit;
1525 double ra_min = avg_sx - ra_lim, ra_max = avg_sx + ra_lim;
1526 double dec_min = avg_sy - dec_lim, dec_max = avg_sy + dec_lim;
1527 double d_ra = (ra_max - ra_min) / (n_ra - 1);
1528 double d_dec = (dec_max - dec_min) / (n_dec - 1);
1530 ra_grid = gsl_vector_alloc(n_ra);
1531 dec_grid = gsl_vector_alloc(n_dec);
1532 for (
int i = 0; i < n_ra; i++)
1533 gsl_vector_set(ra_grid, i, ra_min + i * d_ra);
1534 for (
int j = 0; j < n_dec; j++)
1535 gsl_vector_set(dec_grid, j, dec_min + j * d_dec);
1537 ra_dec_fit = gsl_vector_alloc(2);
1538 chi2_map = gsl_matrix_alloc(n_ra, n_dec);
1541 cpl_msg_debug(cpl_func,
"RA grid: [%f, %f] with %d points", ra_min, ra_max, n_ra);
1542 cpl_msg_debug(cpl_func,
"dec grid: [%f, %f] with %d points", dec_min, dec_max, n_dec);
1545 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1546 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1547 cpl_msg_debug(cpl_func,
"big chi2 map minimsed for ra=%f, dec=%f", ra_fit, dec_fit);
1553 ra_lim = (ra_max - ra_min) / (4 * zoom);
1554 dec_lim = (dec_max - dec_min) / (4 * zoom);
1556 ra_min = ra_fit - ra_lim;
1557 ra_max = ra_fit + ra_lim;
1558 dec_min = dec_fit - dec_lim;
1559 dec_max = dec_fit + dec_lim;
1561 d_ra = (ra_max - ra_min) / (n_ra - 1);
1562 d_dec = (dec_max - dec_min) / (n_dec - 1);
1564 for (
int i = 0; i < n_ra; i++)
1565 gsl_vector_set(ra_grid, i, ra_min + i * d_ra);
1566 for (
int j = 0; j < n_dec; j++)
1567 gsl_vector_set(dec_grid, j, dec_min + j * d_dec);
1570 cpl_msg_debug(cpl_func,
"RA grid: [%f, %f] with %d points", ra_min, ra_max, n_ra);
1571 cpl_msg_debug(cpl_func,
"dec grid: [%f, %f] with %d points", dec_min, dec_max, n_dec);
1574 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1575 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1576 cpl_msg_debug(cpl_func,
"zoomed chi2 map minimised for ra=%f, dec=%f", ra_fit, dec_fit);
1582 for (
int i = 0; i < nswap; i++) {
1592 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1593 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1594 cpl_msg_debug(cpl_func,
"final astrometric solution after gradient descent is ra=%f, dec=%f", ra_fit, dec_fit);
1596 for (
int i = 0; i < nswap; i++) {
1606 CPLCHECK_MSG(
"Could not extract swap phase reference");
1608 for (
int i = 1; i < nswap; i++)
1609 swap_data[i]->phase_ref_astro = gsl_matrix_complex_alloc_from_matrix(swap_ref, 0, 0, swap_ref->size1, swap_ref->size2);
1616 FREE(cpl_free, group1);
1617 FREE(cpl_free, group2);
1620 FREE(gsl_matrix_free, chi2_map);
1621 FREE(gsl_vector_free, ra_grid);
1622 FREE(gsl_vector_free, dec_grid);
1623 FREE(gsl_vector_free, ra_dec_fit);
1625 return CPL_ERROR_NONE;
1636 cpl_ensure(swap_data, CPL_ERROR_NULL_INPUT, NULL);
1638 int nchannel = swap_data[0]->
nchannel;
1639 int nwave = swap_data[0]->
nwave;
1641 gsl_matrix_complex *phase_ref_s1 = gsl_matrix_complex_calloc(nchannel, nwave);
1642 gsl_matrix_complex *phase_ref_s2 = gsl_matrix_complex_calloc(nchannel, nwave);
1643 gsl_matrix_complex *phase_ref_avg = gsl_matrix_complex_calloc(nchannel, nwave);
1645 for (
int n = 0; n < nswap; n++) {
1649 soi->
nwave == nwave,
1650 CPL_ERROR_INCOMPATIBLE_INPUT,
1664 gsl_matrix_complex_add(phase_ref_s1, vis_ref_avg);
1666 gsl_matrix_complex_add(phase_ref_s2, vis_ref_avg);
1667 FREE(gsl_matrix_complex_free, vis_ref_avg);
1671 gsl_matrix_complex_memcpy(phase_ref_avg, phase_ref_s1);
1672 gsl_matrix_complex_add(phase_ref_avg, phase_ref_s2);
1673 gsl_matrix_complex_scale(phase_ref_avg, gsl_complex_rect(0.5, 0));
1675 FREE(gsl_matrix_complex_free, phase_ref_s1);
1676 FREE(gsl_matrix_complex_free, phase_ref_s2);
1678 return phase_ref_avg;
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)
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)
cpl_error_code gravi_astrometry_dump(astro_data *self, FILE *handle)
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
#define gravi_data_get_header(data)
#define gravi_data_get_oi_wave(data, type, pol, npol)
#define gravi_data_get_oi_vis(data, type, pol, npol)
#define gravi_data_get_oi_wave_plist(data, type, pol, npol)
cpl_msg_debug(cpl_func, "Spectra has <50 pixels -> don't flat")
#define CPLCHECK_INT(msg)
#define CPLCHECK_CLEAN(msg)
#define FREE(function, variable)
#define CPLCHECK_NUL(msg)
#define CPLCHECK_MSG(msg)
#define FREELOOP(function, variable, n)
gsl_matrix_complex * phase_ref_astro
gsl_vector ** vis_ref_cov
gsl_matrix_complex * visdata
gsl_vector_complex ** vis_ref_pcov
gsl_matrix * amp_ref_astro
double swap_astrometry_fit[2]
double swap_astrometry_guess[2]
gsl_matrix_complex * vis_ref
gsl_matrix_complex * viserr
gsl_matrix_complex * visdata_ft
gsl_matrix * phase_met_telfc