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>
40#define PI 3.14159265359
41#define TWOPI 6.28318530718
42#define MAS_TO_RAD (2.0 * PI / (3600 * 360 * 1000))
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;
100static gsl_matrix *
average_matrix_over_dits(gsl_matrix *m,
int ndit,
int nchannel,
int nwave, gsl_matrix_int *flag) CPL_ATTR_ALLOC;
122 const size_t N = a->size;
123 const size_t stride = a->stride;
127 for (i = 0; i < N; i++)
128 sum += a->data[i * stride];
137 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
138 cpl_ensure(name != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
140 cpl_size depth = cpl_table_get_column_depth(table, name);
142 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
144 }
else if (depth > 0) {
145 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of scalar type", name);
149 cpl_size nelem = cpl_table_get_nrow(table);
150 double * cpl_vec = cpl_table_get_data_double(table, name);
152 gsl_vector *gsl_vec = gsl_vector_alloc(nelem);
153 memcpy(gsl_vec->data, cpl_vec, nelem *
sizeof(
double));
163 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
164 cpl_ensure(name != NULL, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
166 cpl_size depth = cpl_table_get_column_depth(table, name);
168 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
170 }
else if (depth == 0) {
171 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of 2D array type", name);
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);
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);
194 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
195 cpl_ensure(name != NULL, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
197 cpl_size depth = cpl_table_get_column_depth(table, name);
199 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
201 }
else if (depth == 0) {
202 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of 2D array type", name);
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);
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);
230 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
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);
237 double factor2 = factor * factor;
240 gsl_vector_complex_scale(self->
vis_ref_pcov[i], gsl_complex_rect(factor2, 0));
243 return CPL_ERROR_NONE;
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);
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);
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);
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);
273 double factor2 = gsl_vector_get(factor, i);
276 gsl_vector_complex_scale(self->
vis_ref_pcov[i], gsl_complex_rect(factor2, 0));
279 FREE(gsl_matrix_complex_free, complex_factor);
280 FREE(gsl_matrix_complex_free, complex_factor_ft);
282 return CPL_ERROR_NONE;
290 gsl_matrix_complex *phi = gsl_matrix_complex_alloc(self->
ndit * self->
nchannel, self->
nwave);
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);
298 gsl_matrix_complex_mul_elements(self->
vis_ref, phi);
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);
305 gsl_vector_complex_set_zero(cov);
306 gsl_vector_complex_set_zero(pcov);
307 gsl_vector_complex_set_zero(tmp);
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);
316 gsl_vector_view cov_real = gsl_vector_complex_real(cov);
317 gsl_vector_memcpy(&cov_real.vector, self->
vis_ref_cov[i]);
320 gsl_vector_complex_mul(tmp, cov);
321 gsl_vector_complex_mul(tmp, &phi_row.vector);
322 cov_real = gsl_vector_complex_real(tmp);
323 gsl_vector_memcpy(self->
vis_ref_cov[i], &cov_real.vector);
327 gsl_vector_complex_mul(pcov, &phi_row.vector);
328 gsl_vector_complex_mul(pcov, &phi_row.vector);
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);
337 return CPL_ERROR_NONE;
346 gsl_matrix_memcpy(this_u, self->
ucoord);
347 gsl_matrix_scale(this_u, xvalue);
350 gsl_matrix_memcpy(this_v, self->
vcoord);
351 gsl_matrix_scale(this_v, yvalue);
353 gsl_matrix_add(this_u, this_v);
358 FREE(gsl_matrix_free, this_u);
359 FREE(gsl_matrix_free, this_v);
361 return CPL_ERROR_NONE;
369 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
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));
379 if (row_abs_sum < threshold) {
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);
408 cpl_msg_debug(cpl_func,
"Flagging %d points below FT threshold %.3e", npoints, threshold);
409 self->
nflag += npoints;
411 return CPL_ERROR_NONE;
419 gsl_vector *ftflux_scale = gsl_vector_alloc(self->
ndit * self->
nchannel);
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);
428 FREE(gsl_vector_free, ftflux_scale);
429 return CPL_ERROR_NONE;
437 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
443 const char *filename = cpl_propertylist_get_string(hdr,
"PIPEFILE");
444 self->
filename = cpl_strdup(filename);
446 const char *insname = cpl_propertylist_get_string(hdr,
"TELESCOP");
447 self->
insname = cpl_strdup(insname);
450 if (cpl_propertylist_has(hdr,
"ESO TPL NDIT OBJECT"))
451 self->
ndit = cpl_propertylist_get_int(hdr,
"ESO TPL NDIT OBJECT");
453 self->
ndit = cpl_propertylist_get_int(hdr,
"ESO TPL NDIT SKY");
454 self->
dit = cpl_propertylist_get_double(hdr,
"ESO DET2 SEQ1 DIT");
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"));
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");
465 self->
mjd = cpl_propertylist_get_double(hdr,
"MJD-OBS");
473 self->
nwave = cpl_propertylist_get_int(wave_plist,
"NWAVE");
476 self->
nwave_ft = cpl_propertylist_get_int(wave_plist,
"NWAVE");
480 cpl_table_cast_column(wave_table,
"EFF_WAVE",
"EFF_WAVE", CPL_TYPE_DOUBLE);
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);
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);
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");
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);
530 gsl_matrix_int_set(self->
flag, i, j, flag_val);
531 self->
nflag += flag_val;
538 for (
int j = 0; j < self->
nwave; j++) {
539 if (gsl_matrix_int_get(self->
flag, i, j)) {
541 gsl_matrix_complex_set(self->
visdata, i, j, GSL_COMPLEX_ZERO);
545 cpl_msg_debug(cpl_func,
"Cleaned %d flagged values in VISDATA", ncleaned);
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);
576 gsl_vector *vreal = gsl_vector_alloc(self->
nwave);
577 gsl_vector *vimag = gsl_vector_alloc(self->
nwave);
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);
584 gsl_vector_memcpy(vreal, &real_vals.vector);
585 gsl_vector_mul(vreal, vreal);
587 gsl_vector_memcpy(vimag, &imag_vals.vector);
588 gsl_vector_mul(vimag, vimag);
591 gsl_vector *cov_diag = gsl_vector_alloc(self->
nwave);
592 gsl_vector_memcpy(cov_diag, vreal);
593 gsl_vector_add(cov_diag, vimag);
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);
602 for (
int j = 0; j < self->
nwave; j++) {
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);
610 FREE(gsl_vector_free, vreal);
611 FREE(gsl_vector_free, vimag);
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);
620 gsl_vector *wavenumber = gsl_vector_alloc(self->
nwave);
621 gsl_vector_set_all(wavenumber,
TWOPI);
622 gsl_vector_div(wavenumber, self->
wave);
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);
632 gsl_matrix_scale(wave_opd_disp, -1);
635 FREE(gsl_vector_free, wavenumber);
636 FREE(gsl_matrix_free, wave_opd_disp);
643 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
644 cpl_ensure_code(handle, CPL_ERROR_NULL_INPUT);
646 fprintf(handle,
"%s\n", self->
filename);
647 fprintf(handle,
"%s\n", self->
insname);
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");
657 return CPL_ERROR_NONE;
665 FREE(gsl_vector_free, self->
u);
666 FREE(gsl_vector_free, self->
v);
675 FREE(gsl_matrix_int_free, self->
flag);
691 cpl_ensure(self, CPL_ERROR_NULL_INPUT, -1);
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));
715 cpl_ensure(self, CPL_ERROR_NULL_INPUT, -1);
716 cpl_ensure(others, CPL_ERROR_NULL_INPUT, -1);
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;
722 any_before = CPL_TRUE;
724 any_after = CPL_TRUE;
728 double min_delta = INFINITY;
730 for (
int i = 0; i < n_other; i++) {
731 double delta = self->
mjd - others[i]->
mjd;
734 if (any_before && delta >= 0 && delta < min_delta) {
740 if (any_after && delta <= 0 && -delta < min_delta) {
746 if (delta <= 0 && fabs(delta) < min_delta) {
768 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
769 cpl_ensure_code(phase_refs, CPL_ERROR_NULL_INPUT);
771 cpl_ensure_code(swaps, CPL_ERROR_NULL_INPUT);
772 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
774 gsl_matrix *amp_ref = NULL, *amp_ref_tmp = NULL;
775 gsl_matrix_complex *vis_ref = NULL, *vis_ref_tmp = NULL;
778 const char *calib_strategy = cpl_parameter_get_string(
779 cpl_parameterlist_find(parlist,
"gravity.astrometry.calib-strategy"));
782 cpl_msg_debug(cpl_func,
"Calibrating amplitude reference with '%s' strategy", calib_strategy);
784 if (!strcmp(calib_strategy,
"none")) {
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")) {
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));
799 FREE(gsl_matrix_free, amp_ref_tmp);
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));
810 gsl_vector_int *phase_indices;
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")) {
819 phase_indices = gsl_vector_int_alloc(nphase_used);
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;
829 gsl_vector_int_set(phase_indices, 0, idx_before);
830 gsl_vector_int_set(phase_indices, 1, idx_after);
832 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"Unknown calibration strategy %s", calib_strategy);
833 return CPL_ERROR_ILLEGAL_INPUT;
837 vis_ref = gsl_matrix_complex_calloc(self->
nchannel, self->
nwave);
840 for (
int n = 0; n < nphase_used; n++) {
841 cpl_size idx = gsl_vector_int_get(phase_indices, n);
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));
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);
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));
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));
870 FREE(gsl_matrix_free, amp_ref_tmp);
871 FREE(gsl_vector_int_free, phase_indices);
879 cpl_msg_debug(cpl_func,
"Extract phase reference from swaps");
885 if (!strcmp(calib_strategy,
"swap")) {
886 cpl_size n_swapon = 0;
888 gsl_matrix_set_zero(amp_ref);
889 gsl_matrix_complex_set_zero(vis_ref);
891 for (
int n = 0; n < nswap; n++) {
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));
905 gsl_matrix_add(amp_ref, amp_ref_avg);
907 FREE(gsl_matrix_free, amp_ref_tmp);
908 FREE(gsl_matrix_free, amp_ref_avg);
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++) {
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));
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);
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);
932 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_polar(2 * vis_ref_abs, argphi));
937 return CPL_ERROR_NONE;
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);
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));
1023 cpl_table_set_array(table,
"ASTRO_VISREF", j, tmp_arr);
1026 FREE(cpl_array_delete, tmp_arr);
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++) {
1032 cpl_array_set(tmp_arr, i, val);
1034 cpl_table_set_array(table,
"ASTRO_AMPREF", j, tmp_arr);
1037 FREE(cpl_array_delete, tmp_arr);
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);
1046 cpl_table_set_array(table,
"ASTRO_COV", j, tmp_arr);
1049 FREE(cpl_array_delete, tmp_arr);
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));
1058 cpl_table_set_array(table,
"ASTRO_PCOV", j, tmp_arr);
1061 FREE(cpl_array_delete, tmp_arr);
1068 gsl_vector *vavg = gsl_vector_alloc(nchannel);
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);
1080 gsl_matrix *mavg = gsl_matrix_alloc(nchannel, nwave);
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){};
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);
1101 gsl_matrix_set(mavg, i, j, avg);
1103 gsl_matrix_set(mavg, i, j, 0.0);
1113 gsl_matrix_complex *mavg = gsl_matrix_complex_alloc(nchannel, nwave);
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);
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){};
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));
1133 avg = gsl_complex_div_real(avg, Nvalid);
1134 gsl_matrix_complex_set(mavg, i, j, avg);
1136 gsl_matrix_complex_set(mavg, i, j, GSL_COMPLEX_ZERO);
1149 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
1170 flat->
wave = gsl_vector_alloc(flat->
nwave);
1171 gsl_vector_memcpy(flat->
wave, self->
wave);
1189 for (
int i = 0; i < self->
nchannel; i++) {
1193 for (
int j = 0; j < self->
ndit; j++) {
1194 cpl_size idx = j * self->
nchannel + i;
1199 gsl_vector_complex_scale(flat->
vis_ref_pcov[i], gsl_complex_rect(1.0 / (self->
ndit * self->
ndit), 0));
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(
1211 gsl_matrix_int_set(flat->
flag, i, j, 1);
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);
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);
1233 for (
int n = 0; n < ndata; n++) {
1236 for (
int j = 0; j < self->
nchannel; j++) {
1237 for (
int k = 0; k < self->
nwave; k++) {
1239 gsl_complex visref_jk = GSL_COMPLEX_ZERO;
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;
1246 gsl_complex phi = gsl_complex_polar(1.0, phase);
1248 phi = gsl_complex_conjugate(phi);
1249 gsl_complex visref_ijk = gsl_complex_mul(phi, gsl_matrix_complex_get(self->
vis_ref, idx, k));
1251 int flag = gsl_matrix_int_get(self->
flag, idx, k);
1254 visref_jk = gsl_complex_add(visref_jk, visref_ijk);
1257 gsl_matrix_complex_set(visref_ditsum, j, k, visref_jk);
1258 gsl_matrix_int_set(ngood_ditsum, j, k, ngood_jk);
1261 gsl_matrix_complex_add(visref, visref_ditsum);
1262 gsl_matrix_int_add(ngood, ngood_ditsum);
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);
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);
1276 FREE(gsl_matrix_complex_free, visref_ditsum);
1292 cpl_size nchannel = s1[0]->
nchannel;
1293 cpl_size nwave = s1[0]->
nwave;
1295 double ra = gsl_vector_get(X, 0);
1296 double dec = gsl_vector_get(X, 1);
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))
1309 visref_swap = gsl_complex_sqrt(visref_swap);
1310 chi2_channel += GSL_IMAG(visref_swap) * GSL_IMAG(visref_swap);
1312 chi2 += chi2_channel;
1315 FREE(gsl_matrix_complex_free, visref_s1);
1316 FREE(gsl_matrix_complex_free, visref_s2);
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);
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;
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);
1340 &ra_dec_view.vector, params
1344 gsl_matrix_set(chi2_map, i, j, chi2);
1346 if (chi2 < best_chi2) {
1347 best_ra = ra_dec[0];
1348 best_dec = ra_dec[1];
1355 best_chi2, best_ra, best_dec);
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;
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);
1373 gsl_multimin_function func = {
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);
1388 status = gsl_multimin_fminimizer_iterate(mini);
1394 size = gsl_multimin_fminimizer_size(mini);
1395 status = gsl_multimin_test_size(size, 1e-3);
1396 }
while (status == GSL_CONTINUE && iterations < MAX_ITERATIONS);
1398 if (status == GSL_SUCCESS)
1399 gsl_vector_memcpy(Xsolve, mini->x);
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));
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);
1408 FREE(gsl_multimin_fminimizer_free, mini);
1409 return status ? CPL_ERROR_ILLEGAL_OUTPUT : CPL_ERROR_NONE;
1414 cpl_ensure_code(swap_data, CPL_ERROR_NULL_INPUT);
1415 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
1418 cpl_size n1 = 0, n2 = 0;
1420 gsl_vector *ra_grid = NULL, *dec_grid = NULL;
1422 gsl_vector *ra_dec_fit = NULL;
1423 gsl_matrix *chi2_map = NULL;
1426 cpl_boolean use_swap_fiber_pos = cpl_parameter_get_bool(
1427 cpl_parameterlist_find(parlist,
"gravity.astrometry.use-swap-fiber-pos"));
1429 cpl_boolean go_fast = cpl_parameter_get_bool(
1430 cpl_parameterlist_find(parlist,
"gravity.astrometry.average-over-dits"));
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"));
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"));
1442 double zoom = cpl_parameter_get_double(
1443 cpl_parameterlist_find(parlist,
"gravity.astrometry.zoom-factor"));
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++) {
1459 return CPL_ERROR_NONE;
1463 for (
int i = 0; i < nswap; i++) {
1464 if (swap_data[i]->swap)
1471 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"No files with SWAP=YES");
1473 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"No files with SWAP=NO");
1476 group1 = cpl_malloc(n1 *
sizeof(
astro_data*));
1477 group2 = cpl_malloc(n2 *
sizeof(
astro_data*));
1478 cpl_size i1 = 0, i2 = 0;
1480 for (
int i = 0; i < nswap; i++) {
1481 if (swap_data[i]->swap) {
1489 .
group1 = group1, .group2 = group2, .n1 = n1, .n2 = n2
1494 if (strstr(swap_data[0]->insname,
"U1234"))
1502 dec_lim = field / 2;
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);
1511 cpl_msg_debug(cpl_func,
"Centre: [%f, %f] mas", avg_sx, avg_sy);
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);
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);
1526 ra_dec_fit = gsl_vector_alloc(2);
1527 chi2_map = gsl_matrix_alloc(n_ra, n_dec);
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);
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);
1542 ra_lim = (ra_max - ra_min) / (4 * zoom);
1543 dec_lim = (dec_max - dec_min) / (4 * zoom);
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;
1550 d_ra = (ra_max - ra_min) / (n_ra - 1);
1551 d_dec = (dec_max - dec_min) / (n_dec - 1);
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);
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);
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);
1571 for (
int i = 0; i < nswap; i++) {
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);
1585 for (
int i = 0; i < nswap; i++) {
1595 CPLCHECK_MSG(
"Could not extract swap phase reference");
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);
1605 FREE(cpl_free, group1);
1606 FREE(cpl_free, group2);
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);
1614 return CPL_ERROR_NONE;
1625 cpl_ensure(swap_data, CPL_ERROR_NULL_INPUT, NULL);
1627 int nchannel = swap_data[0]->
nchannel;
1628 int nwave = swap_data[0]->
nwave;
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);
1634 for (
int n = 0; n < nswap; n++) {
1638 soi->
nwave == nwave,
1639 CPL_ERROR_INCOMPATIBLE_INPUT,
1653 gsl_matrix_complex_add(phase_ref_s1, vis_ref_avg);
1655 gsl_matrix_complex_add(phase_ref_s2, vis_ref_avg);
1656 FREE(gsl_matrix_complex_free, vis_ref_avg);
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));
1664 FREE(gsl_matrix_complex_free, phase_ref_s1);
1665 FREE(gsl_matrix_complex_free, phase_ref_s2);
1667 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