38#include <gsl/gsl_blas.h>
39#include <gsl/gsl_complex.h>
40#include <gsl/gsl_complex_math.h>
41#include <gsl/gsl_statistics.h>
42#include <gsl/gsl_multimin.h>
44#define PI 3.14159265359
45#define TWOPI 6.28318530718
46#define MAS_TO_RAD (2.0 * PI / (3600 * 360 * 1000))
93static gsl_vector *
load_vector(cpl_table *table,
const char *name) CPL_ATTR_ALLOC;
94static gsl_matrix *
load_matrix(cpl_table *table,
const char *name) CPL_ATTR_ALLOC;
95static gsl_matrix_complex *
load_matrix_complex(cpl_table *table,
const char *name) CPL_ATTR_ALLOC;
104static gsl_matrix *
average_matrix_over_dits(gsl_matrix *m,
int ndit,
int nchannel,
int nwave, gsl_matrix_int *flag) CPL_ATTR_ALLOC;
126 const size_t N = a->size;
127 const size_t stride = a->stride;
131 for (i = 0; i < N; i++)
132 sum += a->data[i * stride];
141 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
142 cpl_ensure(name != NULL, CPL_ERROR_ILLEGAL_INPUT, NULL);
144 cpl_size depth = cpl_table_get_column_depth(table, name);
146 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
148 }
else if (depth > 0) {
149 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of scalar type", name);
153 cpl_size nelem = cpl_table_get_nrow(table);
154 double * cpl_vec = cpl_table_get_data_double(table, name);
156 gsl_vector *gsl_vec = gsl_vector_alloc(nelem);
157 memcpy(gsl_vec->data, cpl_vec, nelem *
sizeof(
double));
167 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
168 cpl_ensure(name != NULL, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
170 cpl_size depth = cpl_table_get_column_depth(table, name);
172 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
174 }
else if (depth == 0) {
175 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of 2D array type", name);
179 cpl_size nrow = cpl_table_get_nrow(table);
180 cpl_size ncol = cpl_table_get_column_dimension(table, name, 0);
181 const cpl_array ** cpl_matrix = cpl_table_get_data_array_const(table, name);
183 gsl_matrix *gsl_matrix = gsl_matrix_alloc(nrow, ncol);
184 for (
int i = 0; i < nrow; i++) {
185 gsl_vector_const_view row_view = gsl_vector_const_view_array(
186 cpl_array_get_data_double_const(cpl_matrix[i]), ncol);
187 gsl_matrix_set_row(gsl_matrix, i, &row_view.vector);
198 cpl_ensure(table != NULL, CPL_ERROR_NULL_INPUT, NULL);
199 cpl_ensure(name != NULL, CPL_ERROR_ACCESS_OUT_OF_RANGE, NULL);
201 cpl_size depth = cpl_table_get_column_depth(table, name);
203 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
"Field %s not found", name);
205 }
else if (depth == 0) {
206 cpl_error_set_message(cpl_func, CPL_ERROR_INVALID_TYPE,
"Field %s not of 2D array type", name);
210 cpl_size nrow = cpl_table_get_nrow(table);
211 cpl_size ncol = cpl_table_get_column_dimension(table, name, 0);
212 const cpl_array ** cpl_matrix = cpl_table_get_data_array_const(table, name);
214 gsl_matrix_complex *gsl_matrix = gsl_matrix_complex_alloc(nrow, ncol);
215 for (
int i = 0; i < nrow; i++) {
216 for (
int j = 0; j < ncol; j++) {
217 double complex cpl_matrix_val = cpl_array_get_complex(cpl_matrix[i], j, NULL);
218 gsl_complex gsl_matrix_val = gsl_complex_rect(creal(cpl_matrix_val), cimag(cpl_matrix_val));
219 gsl_matrix_complex_set(gsl_matrix, i, j, gsl_matrix_val);
234 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
236 gsl_complex complex_factor = gsl_complex_rect(factor, 0);
237 gsl_matrix_complex_scale(self->
visdata, complex_factor);
238 gsl_matrix_complex_scale(self->
visdata_ft, complex_factor);
239 gsl_matrix_complex_scale(self->
vis_ref, complex_factor);
241 double factor2 = factor * factor;
244 gsl_vector_complex_scale(self->
vis_ref_pcov[i], gsl_complex_rect(factor2, 0));
247 return CPL_ERROR_NONE;
258 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
259 cpl_ensure_code(factor, CPL_ERROR_NULL_INPUT);
260 cpl_ensure_code((
size_t) (self->
ndit * self->
nchannel) == factor->size, CPL_ERROR_INCOMPATIBLE_INPUT);
262 gsl_matrix_complex *complex_factor = gsl_matrix_complex_alloc(self->
ndit * self->
nchannel, self->
nwave);
263 gsl_matrix_complex *complex_factor_ft = gsl_matrix_complex_alloc(self->
ndit * self->
nchannel, self->
nwave_ft);
265 gsl_complex cf = gsl_complex_rect(gsl_vector_get(factor, i), 0);
266 for (
int j = 0; j < self->
nwave; j++)
267 gsl_matrix_complex_set(complex_factor, i, j, cf);
268 for (
int j = 0; j < self->
nwave_ft; j++)
269 gsl_matrix_complex_set(complex_factor_ft, i, j, cf);
272 gsl_matrix_complex_mul_elements(self->
visdata, complex_factor);
273 gsl_matrix_complex_mul_elements(self->
vis_ref, complex_factor);
274 gsl_matrix_complex_mul_elements(self->
visdata_ft, complex_factor_ft);
277 double factor2 = gsl_vector_get(factor, i);
280 gsl_vector_complex_scale(self->
vis_ref_pcov[i], gsl_complex_rect(factor2, 0));
283 FREE(gsl_matrix_complex_free, complex_factor);
284 FREE(gsl_matrix_complex_free, complex_factor_ft);
286 return CPL_ERROR_NONE;
294 gsl_matrix_complex *phi = gsl_matrix_complex_alloc(self->
ndit * self->
nchannel, self->
nwave);
296 for (
int j = 0; j < self->
nwave; j++) {
297 gsl_complex phi_val = gsl_complex_polar(1.0, gsl_matrix_get(phase, i, j));
298 gsl_matrix_complex_set(phi, i, j, phi_val);
302 gsl_matrix_complex_mul_elements(self->
vis_ref, phi);
304 gsl_vector_complex *cov = gsl_vector_complex_alloc(self->
nwave);
305 gsl_vector_complex *pcov = gsl_vector_complex_alloc(self->
nwave);
306 gsl_vector_complex *tmp = gsl_vector_complex_alloc(self->
nwave);
309 gsl_vector_complex_set_zero(cov);
310 gsl_vector_complex_set_zero(pcov);
311 gsl_vector_complex_set_zero(tmp);
314 gsl_vector_complex_const_view phi_row = gsl_matrix_complex_const_row(phi, i);
315 gsl_vector_complex_memcpy(tmp, &phi_row.vector);
316 gsl_vector_view phi_imag = gsl_vector_complex_imag(tmp);
317 gsl_vector_scale(&phi_imag.vector, -1);
320 gsl_vector_view cov_real = gsl_vector_complex_real(cov);
321 gsl_vector_memcpy(&cov_real.vector, self->
vis_ref_cov[i]);
324 gsl_vector_complex_mul(tmp, cov);
325 gsl_vector_complex_mul(tmp, &phi_row.vector);
326 cov_real = gsl_vector_complex_real(tmp);
327 gsl_vector_memcpy(self->
vis_ref_cov[i], &cov_real.vector);
331 gsl_vector_complex_mul(pcov, &phi_row.vector);
332 gsl_vector_complex_mul(pcov, &phi_row.vector);
336 FREE(gsl_vector_complex_free, cov);
337 FREE(gsl_vector_complex_free, pcov);
338 FREE(gsl_vector_complex_free, tmp);
339 FREE(gsl_matrix_complex_free, phi);
341 return CPL_ERROR_NONE;
350 gsl_matrix_memcpy(this_u, self->
ucoord);
351 gsl_matrix_scale(this_u, xvalue);
354 gsl_matrix_memcpy(this_v, self->
vcoord);
355 gsl_matrix_scale(this_v, yvalue);
357 gsl_matrix_add(this_u, this_v);
362 FREE(gsl_matrix_free, this_u);
363 FREE(gsl_matrix_free, this_v);
365 return CPL_ERROR_NONE;
373 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
377 gsl_vector_complex_const_view visft_row = gsl_matrix_complex_const_row(self->
visdata_ft, i);
378 double row_abs_sum = 0.0;
379 for (
int j = 0; j < self->
nwave_ft; j++)
380 row_abs_sum += gsl_complex_abs(gsl_vector_complex_get(&visft_row.vector, j));
383 if (row_abs_sum < threshold) {
385 for (
int k = 0; k < self->
nwave; k++) {
386 gsl_matrix_int_set(self->
flag, i, k, 1);
387 gsl_matrix_complex_set(self->
visdata, i, k, GSL_COMPLEX_ZERO);
412 cpl_msg_debug(cpl_func,
"Flagging %d points below FT threshold %.3e", npoints, threshold);
413 self->
nflag += npoints;
415 return CPL_ERROR_NONE;
423 gsl_vector *ftflux_scale = gsl_vector_alloc(self->
ndit * self->
nchannel);
426 for (
int j = 0; j < self->
nwave_ft; j++)
427 scale += gsl_complex_abs(gsl_matrix_complex_get(self->
visdata_ft, i, j));
428 gsl_vector_set(ftflux_scale, i, self->
nwave / scale);
432 FREE(gsl_vector_free, ftflux_scale);
433 return CPL_ERROR_NONE;
441 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, NULL);
447 const char *filename = cpl_propertylist_get_string(hdr,
"PIPEFILE");
448 self->
filename = cpl_strdup(filename);
450 const char *insname = cpl_propertylist_get_string(hdr,
"TELESCOP");
451 self->
insname = cpl_strdup(insname);
454 if (cpl_propertylist_has(hdr,
"ESO TPL NDIT OBJECT"))
455 self->
ndit = cpl_propertylist_get_int(hdr,
"ESO TPL NDIT OBJECT");
457 self->
ndit = cpl_propertylist_get_int(hdr,
"ESO TPL NDIT SKY");
458 self->
dit = cpl_propertylist_get_double(hdr,
"ESO DET2 SEQ1 DIT");
460 self->
swap = (!strcmp(cpl_propertylist_get_string(hdr,
"ESO INS SOBJ SWAP"),
"YES"));
461 cpl_msg_debug(cpl_func,
"SWAP keyword is: %s", cpl_propertylist_get_string(hdr,
"ESO INS SOBJ SWAP"));
465 self->
sobj_x = cpl_propertylist_get_double(hdr,
"ESO INS SOBJ X");
466 self->
sobj_y = cpl_propertylist_get_double(hdr,
"ESO INS SOBJ Y");
469 self->
mjd = cpl_propertylist_get_double(hdr,
"MJD-OBS");
477 self->
nwave = cpl_propertylist_get_int(wave_plist,
"NWAVE");
480 self->
nwave_ft = cpl_propertylist_get_int(wave_plist,
"NWAVE");
484 cpl_table_cast_column(wave_table,
"EFF_WAVE",
"EFF_WAVE", CPL_TYPE_DOUBLE);
496 for (
int j = 0; j < self->
nwave; j++) {
497 gsl_vector_view col_view;
498 col_view = gsl_matrix_column(self->
ucoord, j);
499 gsl_vector_memcpy(&col_view.vector, self->
u);
500 col_view = gsl_matrix_column(self->
vcoord, j);
501 gsl_vector_memcpy(&col_view.vector, self->
v);
505 gsl_vector_view row_view;
506 row_view = gsl_matrix_row(self->
ucoord, i);
507 gsl_vector_div(&row_view.vector, self->
wave);
508 row_view = gsl_matrix_row(self->
vcoord, i);
509 gsl_vector_div(&row_view.vector, self->
wave);
521 const cpl_array ** cpl_flag = cpl_table_get_data_array_const(oivis_table,
"FLAG");
522 const int* cpl_rej = NULL;
523 if (cpl_table_has_column(oivis_table,
"REJECTION_FLAG"))
524 cpl_rej = cpl_table_get_data_int_const(oivis_table,
"REJECTION_FLAG");
530 int rej_val = (cpl_rej ? cpl_rej[i] : 0) > 0;
531 for (
int j = 0; j < self->
nwave; j++) {
532 int flag_val = cpl_array_get_int(cpl_flag[i], j, NULL);
534 gsl_matrix_int_set(self->
flag, i, j, flag_val);
535 self->
nflag += flag_val;
542 for (
int j = 0; j < self->
nwave; j++) {
543 if (gsl_matrix_int_get(self->
flag, i, j)) {
545 gsl_matrix_complex_set(self->
visdata, i, j, GSL_COMPLEX_ZERO);
549 cpl_msg_debug(cpl_func,
"Cleaned %d flagged values in VISDATA", ncleaned);
568 for (
int j = 0; j < self->
nwave; j++) {
569 gsl_complex vis_ref_val = gsl_matrix_complex_get(self->
visdata, i, j);
570 gsl_complex pha = gsl_complex_polar(1.0, gsl_matrix_get(self->
phase_ref, i, j));
571 vis_ref_val = gsl_complex_mul(vis_ref_val, pha);
572 gsl_matrix_complex_set(self->
vis_ref, i, j, vis_ref_val);
580 gsl_vector *vreal = gsl_vector_alloc(self->
nwave);
581 gsl_vector *vimag = gsl_vector_alloc(self->
nwave);
584 gsl_vector_complex_const_view err_vals = gsl_matrix_complex_const_row(self->
viserr, i);
585 gsl_vector_const_view real_vals = gsl_vector_complex_const_real(&err_vals.vector);
586 gsl_vector_const_view imag_vals = gsl_vector_complex_const_imag(&err_vals.vector);
588 gsl_vector_memcpy(vreal, &real_vals.vector);
589 gsl_vector_mul(vreal, vreal);
591 gsl_vector_memcpy(vimag, &imag_vals.vector);
592 gsl_vector_mul(vimag, vimag);
595 gsl_vector *cov_diag = gsl_vector_alloc(self->
nwave);
596 gsl_vector_memcpy(cov_diag, vreal);
597 gsl_vector_add(cov_diag, vimag);
601 gsl_vector_complex *pcov_diag = gsl_vector_complex_calloc(self->
nwave);
602 gsl_vector_view pcov_diag_real = gsl_vector_complex_real(pcov_diag);
603 gsl_vector_memcpy(&pcov_diag_real.vector, vreal);
604 gsl_vector_sub(&pcov_diag_real.vector, vimag);
606 for (
int j = 0; j < self->
nwave; j++) {
608 gsl_complex pcov_val = gsl_complex_polar(1.0, 2 * gsl_matrix_get(self->
phase_ref, i, j));
609 pcov_val = gsl_complex_mul(gsl_vector_complex_get(pcov_diag, j), pcov_val);
610 gsl_vector_complex_set(pcov_diag, j, pcov_val);
614 FREE(gsl_vector_free, vreal);
615 FREE(gsl_vector_free, vimag);
621 gsl_matrix *wave_opd_disp = gsl_matrix_alloc(self->
ndit * self->
nchannel, self->
nwave);
622 gsl_matrix_memcpy(wave_opd_disp, self->
opd_disp);
624 gsl_vector *wavenumber = gsl_vector_alloc(self->
nwave);
625 gsl_vector_set_all(wavenumber,
TWOPI);
626 gsl_vector_div(wavenumber, self->
wave);
630 gsl_vector_view wave_opd_dist_row = gsl_matrix_row(wave_opd_disp, i);
631 gsl_vector_mul(&wave_opd_dist_row.vector, wavenumber);
636 gsl_matrix_scale(wave_opd_disp, -1);
639 FREE(gsl_vector_free, wavenumber);
640 FREE(gsl_matrix_free, wave_opd_disp);
647 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
648 cpl_ensure_code(handle, CPL_ERROR_NULL_INPUT);
650 fprintf(handle,
"%s\n", self->
filename);
651 fprintf(handle,
"%s\n", self->
insname);
655 gsl_matrix_fprintf(handle, self->
ucoord,
"%f");
656 gsl_matrix_fprintf(handle, self->
vcoord,
"%f");
657 gsl_matrix_int_fprintf(handle, self->
flag,
"%d");
658 gsl_matrix_complex_fprintf(handle, self->
visdata,
"%f");
659 gsl_matrix_complex_fprintf(handle, self->
vis_ref,
"%f");
661 return CPL_ERROR_NONE;
669 FREE(gsl_vector_free, self->
u);
670 FREE(gsl_vector_free, self->
v);
679 FREE(gsl_matrix_int_free, self->
flag);
695 cpl_ensure(self, CPL_ERROR_NULL_INPUT, -1);
699 for (
int j = 0; j < self->
nwave_ft; j++) {
700 abs_val += gsl_complex_abs(gsl_matrix_complex_get(self->
visdata_ft, i, j));
719 cpl_ensure(self, CPL_ERROR_NULL_INPUT, -1);
720 cpl_ensure(others, CPL_ERROR_NULL_INPUT, -1);
722 cpl_boolean any_before = FALSE, any_after = FALSE;
723 for (
int i = 0; i < n_other; i++) {
724 double delta = self->
mjd - others[i]->
mjd;
726 any_before = CPL_TRUE;
728 any_after = CPL_TRUE;
732 double min_delta = INFINITY;
734 for (
int i = 0; i < n_other; i++) {
735 double delta = self->
mjd - others[i]->
mjd;
738 if (any_before && delta >= 0 && delta < min_delta) {
744 if (any_after && delta <= 0 && -delta < min_delta) {
750 if (delta <= 0 && fabs(delta) < min_delta) {
772 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
773 cpl_ensure_code(phase_refs, CPL_ERROR_NULL_INPUT);
775 cpl_ensure_code(swaps, CPL_ERROR_NULL_INPUT);
776 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
778 gsl_matrix *amp_ref = NULL, *amp_ref_tmp = NULL;
779 gsl_matrix_complex *vis_ref = NULL, *vis_ref_tmp = NULL;
782 char *calib_strategy = cpl_strdup(cpl_parameter_get_string(
783 cpl_parameterlist_find(parlist,
"gravity.astrometry.calib-strategy")));
785 for(
int i = 0; calib_strategy[i]; i++)
786 calib_strategy[i] = toupper(calib_strategy[i]);
790 cpl_msg_debug(cpl_func,
"Calibrating amplitude reference with '%s' strategy", calib_strategy);
792 if (!strcmp(calib_strategy,
"NONE")) {
794 vis_ref = gsl_matrix_complex_alloc(self->
nchannel, self->
nwave);
795 gsl_matrix_set_all(amp_ref, 1.0);
796 gsl_matrix_complex_set_all(vis_ref, GSL_COMPLEX_ONE);
797 }
else if (!strcmp(calib_strategy,
"SELF")) {
801 for (
int j = 0; j < self->
nwave; j++) {
802 gsl_complex vis_ref_val = gsl_matrix_complex_get(self->
vis_ref, i, j);
803 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
807 FREE(gsl_matrix_free, amp_ref_tmp);
809 vis_ref = gsl_matrix_complex_alloc(self->
nchannel, self->
nwave);
810 for (
int i = 0; i < self->
nchannel; i++) {
811 for (
int j = 0; j < self->
nwave; j++) {
812 double amp_ref_val = gsl_matrix_get(amp_ref, i, j);
813 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_rect(amp_ref_val, 0));
818 gsl_vector_int *phase_indices;
820 if (!strcmp(calib_strategy,
"ALL")) {
821 nphase_used = nphase;
822 phase_indices = gsl_vector_int_alloc(nphase_used);
823 for (
int n = 0; n < nphase_used; n++)
824 gsl_vector_int_set(phase_indices, n, n);
825 }
else if (!strcmp(calib_strategy,
"NEAREST")) {
827 phase_indices = gsl_vector_int_alloc(nphase_used);
833 if (idx_before == -1 || idx_after == -1) {
834 cpl_error_set_message(cpl_func, CPL_ERROR_ACCESS_OUT_OF_RANGE,
"Cannot find time-bracketing phase references");
835 return CPL_ERROR_ACCESS_OUT_OF_RANGE;
837 gsl_vector_int_set(phase_indices, 0, idx_before);
838 gsl_vector_int_set(phase_indices, 1, idx_after);
840 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"Unknown calibration strategy %s", calib_strategy);
841 return CPL_ERROR_ILLEGAL_INPUT;
845 vis_ref = gsl_matrix_complex_calloc(self->
nchannel, self->
nwave);
848 for (
int n = 0; n < nphase_used; n++) {
849 cpl_size idx = gsl_vector_int_get(phase_indices, n);
856 for (
int i = 0; i < soi->
nchannel; i++) {
857 for (
int j = 0; j < soi->
nwave; j++) {
858 gsl_complex vis_ref_val = gsl_matrix_complex_get(vis_ref_tmp, i, j);
859 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
863 gsl_matrix_add(amp_ref, amp_ref_tmp);
864 gsl_matrix_complex_add(vis_ref, vis_ref_tmp);
865 FREE(gsl_matrix_complex_free, vis_ref_tmp);
868 gsl_matrix_scale(amp_ref, 1.0 / nphase_used);
869 gsl_matrix_complex_scale(vis_ref, gsl_complex_rect(1.0 / nphase_used, 0.0));
871 for (
int i = 0; i < self->
nchannel; i++) {
872 for (
int j = 0; j < self->
nwave; j++) {
873 double theta = gsl_complex_arg(gsl_matrix_complex_get(vis_ref, i, j));
874 double r = gsl_matrix_get(amp_ref, i, j);
875 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_polar(r, theta));
878 FREE(gsl_matrix_free, amp_ref_tmp);
879 FREE(gsl_vector_int_free, phase_indices);
887 cpl_msg_debug(cpl_func,
"Extract phase reference from swaps");
893 if (!strcmp(calib_strategy,
"SWAP")) {
894 cpl_size n_swapon = 0;
896 gsl_matrix_set_zero(amp_ref);
897 gsl_matrix_complex_set_zero(vis_ref);
899 for (
int n = 0; n < nswap; n++) {
907 for (
int j = 0; j < soi->
nwave; j++) {
908 gsl_complex vis_ref_val = gsl_matrix_complex_get(soi->
vis_ref, i, j);
909 gsl_matrix_set(amp_ref_tmp, i, j, gsl_complex_abs(vis_ref_val));
913 gsl_matrix_add(amp_ref, amp_ref_avg);
915 FREE(gsl_matrix_free, amp_ref_tmp);
916 FREE(gsl_matrix_free, amp_ref_avg);
920 gsl_matrix_scale(amp_ref, 1.0 / n_swapon);
921 for (
int i = 0; i < self->
nchannel; i++) {
922 for (
int j = 0; j < self->
nwave; j++) {
924 double amp_ref_val = gsl_matrix_get(amp_ref, i, j);
925 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_rect(0.5 * amp_ref_val, 0));
932 for (
int i = 0; i < self->
nchannel; i++) {
933 for (
int j = 0; j < self->
nwave; j++) {
934 gsl_complex phi = gsl_matrix_complex_get(swap_ref, i, j);
935 double argphi = gsl_complex_arg(phi);
937 gsl_complex vis_ref_val = gsl_matrix_complex_get(vis_ref, i, j);
938 double vis_ref_abs = gsl_complex_abs(vis_ref_val);
940 gsl_matrix_complex_set(vis_ref, i, j, gsl_complex_polar(2 * vis_ref_abs, argphi));
945 cpl_free(calib_strategy);
946 return CPL_ERROR_NONE;
1020 cpl_table *table = cpl_table_new(self->
nwave);
1021 cpl_table_new_column_array(table,
"ASTRO_VISREF", CPL_TYPE_DOUBLE_COMPLEX, self->
nchannel);
1022 cpl_table_new_column_array(table,
"ASTRO_AMPREF", CPL_TYPE_DOUBLE, self->
nchannel);
1023 cpl_table_new_column_array(table,
"ASTRO_COV", CPL_TYPE_DOUBLE, self->
ndit * self->
nchannel);
1024 cpl_table_new_column_array(table,
"ASTRO_PCOV", CPL_TYPE_DOUBLE_COMPLEX, self->
ndit * self->
nchannel);
1026 cpl_array *tmp_arr = cpl_array_new(self->
nchannel, CPL_TYPE_DOUBLE_COMPLEX);
1027 for (
int j = 0; j < self->
nwave; j++) {
1028 for (
int i = 0; i < self->
nchannel; i++) {
1029 gsl_complex val = gsl_matrix_complex_get(self->
phase_ref_astro, i, j);
1030 cpl_array_set_complex(tmp_arr, i, GSL_REAL(val) + I * GSL_IMAG(val));
1032 cpl_table_set_array(table,
"ASTRO_VISREF", j, tmp_arr);
1035 FREE(cpl_array_delete, tmp_arr);
1037 tmp_arr = cpl_array_new(self->
nchannel, CPL_TYPE_DOUBLE);
1038 for (
int j = 0; j < self->
nwave; j++) {
1039 for (
int i = 0; i < self->
nchannel; i++) {
1041 cpl_array_set(tmp_arr, i, val);
1043 cpl_table_set_array(table,
"ASTRO_AMPREF", j, tmp_arr);
1046 FREE(cpl_array_delete, tmp_arr);
1049 tmp_arr = cpl_array_new(self->
ndit * self->
nchannel, CPL_TYPE_DOUBLE);
1050 for (
int j = 0; j < self->
nwave; j++) {
1051 for (
int i = 0; i < self->
ndit * self->
nchannel; i++) {
1052 double val = gsl_vector_get(self->
vis_ref_cov[i], j);
1053 cpl_array_set(tmp_arr, i, val);
1055 cpl_table_set_array(table,
"ASTRO_COV", j, tmp_arr);
1058 FREE(cpl_array_delete, tmp_arr);
1061 tmp_arr = cpl_array_new(self->
ndit * self->
nchannel, CPL_TYPE_DOUBLE_COMPLEX);
1062 for (
int j = 0; j < self->
nwave; j++) {
1063 for (
int i = 0; i < self->
ndit * self->
nchannel; i++) {
1064 gsl_complex val = gsl_vector_complex_get(self->
vis_ref_pcov[i], j);
1065 cpl_array_set_complex(tmp_arr, i, GSL_REAL(val) + I * GSL_IMAG(val));
1067 cpl_table_set_array(table,
"ASTRO_PCOV", j, tmp_arr);
1070 FREE(cpl_array_delete, tmp_arr);
1077 gsl_vector *vavg = gsl_vector_alloc(nchannel);
1079 for (
int i = 0; i < nchannel; i++) {
1080 double avg = gsl_stats_mean(v->data + i, nchannel, ndit);
1081 gsl_vector_set(vavg, i, avg);
1089 gsl_matrix *mavg = gsl_matrix_alloc(nchannel, nwave);
1091 for (
int i = 0; i < nchannel; i++) {
1092 for (
int j = 0; j < nwave; j++) {
1093 cpl_size start = i * nwave + j;
1094 gsl_vector_const_view dit_view = gsl_vector_const_view_array_with_stride(
1095 m->data + start, nwave * nchannel, ndit);
1096 gsl_vector_int_const_view flag_view = flag ?
1097 gsl_vector_int_const_view_array_with_stride(flag->data + start, nwave * nchannel, ndit) :
1098 (gsl_vector_int_const_view){};
1101 cpl_size Nvalid = 0;
1102 for (
int k = 0; k < ndit; k++) {
1103 if (!flag || gsl_vector_int_get(&flag_view.vector, k) == 0) {
1104 avg += gsl_vector_get(&dit_view.vector, k);
1110 gsl_matrix_set(mavg, i, j, avg);
1112 gsl_matrix_set(mavg, i, j, 0.0);
1122 gsl_matrix_complex *mavg = gsl_matrix_complex_alloc(nchannel, nwave);
1124 for (
int i = 0; i < nchannel; i++) {
1125 for (
int j = 0; j < nwave; j++) {
1126 cpl_size start = i * nwave + j;
1127 gsl_vector_complex_const_view dit_view = gsl_vector_complex_const_view_array_with_stride(
1128 m->data + 2 * start, nwave * nchannel, ndit);
1129 gsl_vector_int_const_view flag_view = flag ?
1130 gsl_vector_int_const_view_array_with_stride(flag->data + start, nwave * nchannel, ndit) :
1131 (gsl_vector_int_const_view){};
1133 gsl_complex avg = GSL_COMPLEX_ZERO;
1134 cpl_size Nvalid = 0;
1135 for (
int k = 0; k < ndit; k++) {
1136 if (!flag || gsl_vector_int_get(&flag_view.vector, k) == 0) {
1137 avg = gsl_complex_add(avg, gsl_vector_complex_get(&dit_view.vector, k));
1142 avg = gsl_complex_div_real(avg, Nvalid);
1143 gsl_matrix_complex_set(mavg, i, j, avg);
1145 gsl_matrix_complex_set(mavg, i, j, GSL_COMPLEX_ZERO);
1158 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
1179 flat->
wave = gsl_vector_alloc(flat->
nwave);
1180 gsl_vector_memcpy(flat->
wave, self->
wave);
1198 for (
int i = 0; i < self->
nchannel; i++) {
1202 for (
int j = 0; j < self->
ndit; j++) {
1203 cpl_size idx = j * self->
nchannel + i;
1208 gsl_vector_complex_scale(flat->
vis_ref_pcov[i], gsl_complex_rect(1.0 / (self->
ndit * self->
ndit), 0));
1213 for (
int i = 0; i < self->
nchannel; i++) {
1214 for (
int j = 0; j < self->
nwave; j++) {
1215 cpl_size start = i * self->
nwave + j;
1216 gsl_vector_int_const_view flag_view = gsl_vector_int_const_view_array_with_stride(
1220 gsl_matrix_int_set(flat->
flag, i, j, 1);
1236 gsl_matrix_complex *visref = gsl_matrix_complex_calloc(data[0]->nchannel, data[0]->nwave);
1237 gsl_matrix_complex *visref_ditsum = gsl_matrix_complex_alloc(data[0]->nchannel, data[0]->nwave);
1239 gsl_matrix_int *ngood = gsl_matrix_int_calloc(data[0]->nchannel, data[0]->nwave);
1240 gsl_matrix_int *ngood_ditsum = gsl_matrix_int_alloc(data[0]->nchannel, data[0]->nwave);
1242 for (
int n = 0; n < ndata; n++) {
1245 for (
int j = 0; j < self->
nchannel; j++) {
1246 for (
int k = 0; k < self->
nwave; k++) {
1248 gsl_complex visref_jk = GSL_COMPLEX_ZERO;
1250 for (
int i = 0; i < self->
ndit; i++) {
1251 cpl_size idx = i * self->
nchannel + j;
1252 double uv = gsl_matrix_get(self->
ucoord, idx, k) * ra + gsl_matrix_get(self->
vcoord, idx, k) * dec;
1255 gsl_complex phi = gsl_complex_polar(1.0, phase);
1257 phi = gsl_complex_conjugate(phi);
1258 gsl_complex visref_ijk = gsl_complex_mul(phi, gsl_matrix_complex_get(self->
vis_ref, idx, k));
1260 int flag = gsl_matrix_int_get(self->
flag, idx, k);
1263 visref_jk = gsl_complex_add(visref_jk, visref_ijk);
1266 gsl_matrix_complex_set(visref_ditsum, j, k, visref_jk);
1267 gsl_matrix_int_set(ngood_ditsum, j, k, ngood_jk);
1270 gsl_matrix_complex_add(visref, visref_ditsum);
1271 gsl_matrix_int_add(ngood, ngood_ditsum);
1274 for (
int j = 0; j < data[0]->
nchannel; j++) {
1275 for (
int k = 0; k < data[0]->
nwave; k++) {
1276 int ngood_jk = gsl_matrix_int_get(ngood, j, k);
1278 gsl_complex visref_jk = gsl_matrix_complex_get(visref, j, k);
1279 visref_jk = gsl_complex_div_real(visref_jk, ngood_jk);
1280 gsl_matrix_complex_set(visref, j, k, visref_jk);
1285 FREE(gsl_matrix_complex_free, visref_ditsum);
1301 cpl_size nchannel = s1[0]->
nchannel;
1302 cpl_size nwave = s1[0]->
nwave;
1304 double ra = gsl_vector_get(X, 0);
1305 double dec = gsl_vector_get(X, 1);
1311 for (
int i = 0; i < nchannel; i++) {
1312 double chi2_channel = 0;
1313 for (
int k = 0; k < nwave; k++) {
1314 gsl_complex visref_swap = gsl_complex_mul(
1315 gsl_matrix_complex_get(visref_s1, i, k),
1316 gsl_complex_conjugate(gsl_matrix_complex_get(visref_s2, i, k))
1318 visref_swap = gsl_complex_sqrt(visref_swap);
1319 chi2_channel += GSL_IMAG(visref_swap) * GSL_IMAG(visref_swap);
1321 chi2 += chi2_channel;
1324 FREE(gsl_matrix_complex_free, visref_s1);
1325 FREE(gsl_matrix_complex_free, visref_s2);
1332 cpl_ensure_code(ra_grid, CPL_ERROR_NULL_INPUT);
1333 cpl_ensure_code(dec_grid, CPL_ERROR_NULL_INPUT);
1334 cpl_ensure_code(ra_dec_out, CPL_ERROR_NULL_INPUT);
1336 cpl_size n_ra = ra_grid->size;
1337 cpl_size n_dec = dec_grid->size;
1338 double ra_dec[2] = {0};
1339 gsl_vector_const_view ra_dec_view = gsl_vector_const_view_array(ra_dec, 2);
1340 double best_ra, best_dec, best_chi2 = INFINITY;
1342 for (
int i = 0; i < n_ra; i++) {
1343 ra_dec[0] = gsl_vector_get(ra_grid, i);
1344 cpl_msg_debug(cpl_func,
"calculating chi2 for ra=%f", ra_dec[0]);
1345 for (
int j = 0; j < n_dec; j++) {
1346 ra_dec[1] = gsl_vector_get(dec_grid, j);
1349 &ra_dec_view.vector, params
1353 gsl_matrix_set(chi2_map, i, j, chi2);
1355 if (chi2 < best_chi2) {
1356 best_ra = ra_dec[0];
1357 best_dec = ra_dec[1];
1364 best_chi2, best_ra, best_dec);
1366 gsl_vector_set(ra_dec_out, 0, best_ra);
1367 gsl_vector_set(ra_dec_out, 1, best_dec);
1368 return CPL_ERROR_NONE;
1376 const int MAX_ITERATIONS = 1000;
1377 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex2rand;
1378 gsl_multimin_fminimizer *mini = gsl_multimin_fminimizer_alloc(T, dim);
1382 gsl_multimin_function func = {
1388 double initial_guess[] = {ra_guess, dec_guess};
1389 gsl_vector_view guess_view = gsl_vector_view_array(initial_guess, dim);
1390 double step_size[] = {0.1, 0.1};
1391 gsl_vector_view step_view = gsl_vector_view_array(step_size, dim);
1392 gsl_multimin_fminimizer_set(mini, &func, &guess_view.vector, &step_view.vector);
1397 status = gsl_multimin_fminimizer_iterate(mini);
1403 size = gsl_multimin_fminimizer_size(mini);
1404 status = gsl_multimin_test_size(size, 1e-3);
1405 }
while (status == GSL_CONTINUE && iterations < MAX_ITERATIONS);
1407 if (status == GSL_SUCCESS)
1408 gsl_vector_memcpy(Xsolve, mini->x);
1410 if (status == GSL_CONTINUE)
1411 cpl_msg_warning(cpl_func,
"model-fitting did not converge");
1412 else if (status != GSL_SUCCESS)
1413 cpl_msg_error(cpl_func,
"model-fitting failed with error code %d (%s)", status, gsl_strerror(status));
1415 cpl_msg_debug(cpl_func,
"converged with best-fit:\nra\tdec\tchi2\n%.3f\t%.3f\t%.4f", mini->x->data[0], mini->x->data[1], mini->fval);
1417 FREE(gsl_multimin_fminimizer_free, mini);
1418 return status ? CPL_ERROR_ILLEGAL_OUTPUT : CPL_ERROR_NONE;
1423 cpl_ensure_code(swap_data, CPL_ERROR_NULL_INPUT);
1424 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
1427 cpl_size n1 = 0, n2 = 0;
1429 gsl_vector *ra_grid = NULL, *dec_grid = NULL;
1431 gsl_vector *ra_dec_fit = NULL;
1432 gsl_matrix *chi2_map = NULL;
1435 cpl_boolean use_swap_fiber_pos = cpl_parameter_get_bool(
1436 cpl_parameterlist_find(parlist,
"gravity.astrometry.use-swap-fiber-pos"));
1438 cpl_boolean go_fast = cpl_parameter_get_bool(
1439 cpl_parameterlist_find(parlist,
"gravity.astrometry.average-over-dits"));
1441 double ra_lim = cpl_parameter_get_double(
1442 cpl_parameterlist_find(parlist,
"gravity.astrometry.ra-lim-swap"));
1443 int n_ra = cpl_parameter_get_int(
1444 cpl_parameterlist_find(parlist,
"gravity.astrometry.nra-swap"));
1446 double dec_lim = cpl_parameter_get_double(
1447 cpl_parameterlist_find(parlist,
"gravity.astrometry.ra-lim-swap"));
1448 int n_dec = cpl_parameter_get_int(
1449 cpl_parameterlist_find(parlist,
"gravity.astrometry.ndec-swap"));
1451 double zoom = cpl_parameter_get_double(
1452 cpl_parameterlist_find(parlist,
"gravity.astrometry.zoom-factor"));
1456 if (use_swap_fiber_pos) {
1457 cpl_msg_warning(cpl_func,
"No astrometric solution is computed for swaps. Default to fiber position");
1458 for (
int i = 0; i < nswap; i++) {
1468 return CPL_ERROR_NONE;
1472 for (
int i = 0; i < nswap; i++) {
1473 if (swap_data[i]->swap)
1480 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"No files with SWAP=YES");
1482 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
"No files with SWAP=NO");
1485 group1 = cpl_malloc(n1 *
sizeof(
astro_data*));
1486 group2 = cpl_malloc(n2 *
sizeof(
astro_data*));
1487 cpl_size i1 = 0, i2 = 0;
1489 for (
int i = 0; i < nswap; i++) {
1490 if (swap_data[i]->swap) {
1498 .
group1 = group1, .group2 = group2, .n1 = n1, .n2 = n2
1503 if (strstr(swap_data[0]->insname,
"U1234"))
1511 dec_lim = field / 2;
1513 double avg_sx = 0.0, avg_sy = 0.0;
1514 for (
int i = 0; i < nswap; i++) {
1515 avg_sx += swap_data[i]->
sobj_x * (swap_data[i]->
swap ? -1 : 1);
1516 avg_sy += swap_data[i]->
sobj_y * (swap_data[i]->
swap ? -1 : 1);
1520 cpl_msg_debug(cpl_func,
"Centre: [%f, %f] mas", avg_sx, avg_sy);
1522 double ra_fit, dec_fit;
1523 double ra_min = avg_sx - ra_lim, ra_max = avg_sx + ra_lim;
1524 double dec_min = avg_sy - dec_lim, dec_max = avg_sy + dec_lim;
1525 double d_ra = (ra_max - ra_min) / (n_ra - 1);
1526 double d_dec = (dec_max - dec_min) / (n_dec - 1);
1528 ra_grid = gsl_vector_alloc(n_ra);
1529 dec_grid = gsl_vector_alloc(n_dec);
1530 for (
int i = 0; i < n_ra; i++)
1531 gsl_vector_set(ra_grid, i, ra_min + i * d_ra);
1532 for (
int j = 0; j < n_dec; j++)
1533 gsl_vector_set(dec_grid, j, dec_min + j * d_dec);
1535 ra_dec_fit = gsl_vector_alloc(2);
1536 chi2_map = gsl_matrix_alloc(n_ra, n_dec);
1539 cpl_msg_debug(cpl_func,
"RA grid: [%f, %f] with %d points", ra_min, ra_max, n_ra);
1540 cpl_msg_debug(cpl_func,
"dec grid: [%f, %f] with %d points", dec_min, dec_max, n_dec);
1543 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1544 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1545 cpl_msg_debug(cpl_func,
"big chi2 map minimsed for ra=%f, dec=%f", ra_fit, dec_fit);
1551 ra_lim = (ra_max - ra_min) / (4 * zoom);
1552 dec_lim = (dec_max - dec_min) / (4 * zoom);
1554 ra_min = ra_fit - ra_lim;
1555 ra_max = ra_fit + ra_lim;
1556 dec_min = dec_fit - dec_lim;
1557 dec_max = dec_fit + dec_lim;
1559 d_ra = (ra_max - ra_min) / (n_ra - 1);
1560 d_dec = (dec_max - dec_min) / (n_dec - 1);
1562 for (
int i = 0; i < n_ra; i++)
1563 gsl_vector_set(ra_grid, i, ra_min + i * d_ra);
1564 for (
int j = 0; j < n_dec; j++)
1565 gsl_vector_set(dec_grid, j, dec_min + j * d_dec);
1568 cpl_msg_debug(cpl_func,
"RA grid: [%f, %f] with %d points", ra_min, ra_max, n_ra);
1569 cpl_msg_debug(cpl_func,
"dec grid: [%f, %f] with %d points", dec_min, dec_max, n_dec);
1572 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1573 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1574 cpl_msg_debug(cpl_func,
"zoomed chi2 map minimised for ra=%f, dec=%f", ra_fit, dec_fit);
1580 for (
int i = 0; i < nswap; i++) {
1590 ra_fit = gsl_vector_get(ra_dec_fit, 0);
1591 dec_fit = gsl_vector_get(ra_dec_fit, 1);
1592 cpl_msg_debug(cpl_func,
"final astrometric solution after gradient descent is ra=%f, dec=%f", ra_fit, dec_fit);
1594 for (
int i = 0; i < nswap; i++) {
1604 CPLCHECK_MSG(
"Could not extract swap phase reference");
1606 for (
int i = 1; i < nswap; i++)
1607 swap_data[i]->phase_ref_astro = gsl_matrix_complex_alloc_from_matrix(swap_ref, 0, 0, swap_ref->size1, swap_ref->size2);
1614 FREE(cpl_free, group1);
1615 FREE(cpl_free, group2);
1618 FREE(gsl_matrix_free, chi2_map);
1619 FREE(gsl_vector_free, ra_grid);
1620 FREE(gsl_vector_free, dec_grid);
1621 FREE(gsl_vector_free, ra_dec_fit);
1623 return CPL_ERROR_NONE;
1634 cpl_ensure(swap_data, CPL_ERROR_NULL_INPUT, NULL);
1636 int nchannel = swap_data[0]->
nchannel;
1637 int nwave = swap_data[0]->
nwave;
1639 gsl_matrix_complex *phase_ref_s1 = gsl_matrix_complex_calloc(nchannel, nwave);
1640 gsl_matrix_complex *phase_ref_s2 = gsl_matrix_complex_calloc(nchannel, nwave);
1641 gsl_matrix_complex *phase_ref_avg = gsl_matrix_complex_calloc(nchannel, nwave);
1643 for (
int n = 0; n < nswap; n++) {
1647 soi->
nwave == nwave,
1648 CPL_ERROR_INCOMPATIBLE_INPUT,
1662 gsl_matrix_complex_add(phase_ref_s1, vis_ref_avg);
1664 gsl_matrix_complex_add(phase_ref_s2, vis_ref_avg);
1665 FREE(gsl_matrix_complex_free, vis_ref_avg);
1669 gsl_matrix_complex_memcpy(phase_ref_avg, phase_ref_s1);
1670 gsl_matrix_complex_add(phase_ref_avg, phase_ref_s2);
1671 gsl_matrix_complex_scale(phase_ref_avg, gsl_complex_rect(0.5, 0));
1673 FREE(gsl_matrix_complex_free, phase_ref_s1);
1674 FREE(gsl_matrix_complex_free, phase_ref_s2);
1676 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