52#include "cpl_wlcalib.h"
53#include "irplib_wavecal.h"
78#define PLOT_WAVE_PHASE_VS_OPD 0
86 cpl_table * opdft_table,
87 cpl_table * met_table,
92 cpl_table * detector_table,
93 cpl_table * opdguess_table);
101 cpl_table * detector_table,
102 cpl_table * opd_table);
105 cpl_table * detector_table,
110 double * rms_residuals);
113 cpl_table * weight_individual_table,
114 cpl_table * wave_fitted_table,
115 cpl_table * opd_table,
116 cpl_table * spectrum_table,
117 cpl_table * detector_table,
118 double n0,
double n1,
double n2);
121 double n0,
double n1,
double n2);
124 cpl_table * wavefibre_table,
125 cpl_table * profile_table,
126 cpl_table * detector_table);
151 cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
154 cpl_size nwave = cpl_table_get_column_depth (spectrum_table,
"DATA1");
158 cpl_msg_info (cpl_func,
"spectrum_table = %lld",cpl_table_get_column_depth (spectrum_table,
"DATA1"));
165 cpl_size nlines = 14;
166 double plines_str[] = {60,
182 double plines_lbd[] = {1.982291,
198 cpl_vector * lines_lbd = cpl_vector_wrap (nlines,plines_lbd);
199 cpl_vector * lines_str = cpl_vector_wrap (nlines,plines_str);
200 cpl_bivector * lines = cpl_bivector_wrap_vectors (lines_lbd, lines_str);
203 cpl_wlcalib_slitmodel * model = cpl_wlcalib_slitmodel_new ();
204 cpl_wlcalib_slitmodel_set_wslit (model, 0.1);
205 cpl_wlcalib_slitmodel_set_wfwhm (model, 1.25);
206 cpl_wlcalib_slitmodel_set_threshold (model, 5.0);
207 cpl_wlcalib_slitmodel_set_catalog (model, lines);
210 double values_hr[] = {1.97048,0.2228e-3,2.7728e-08,-4.72e-12};
211 double values_mr[] = {1.9697,0.2179e-2,2.5e-07,-5e-10};
212 double * values = nwave > 1000 ? values_hr : values_mr;
215 cpl_polynomial * dispersion0 = cpl_polynomial_new (1);
216 for (cpl_size order = 0 ; order < 4 ; order ++) {
217 cpl_polynomial_set_coeff (dispersion0, &order, values[order]);
221 double minwave = -1e10;
222 double maxwave = +1e10;
224 cpl_vector * spectrum = cpl_vector_new (nwave);
225 cpl_vector * input_spectrum = cpl_vector_new (nwave);
226 cpl_vector * model_spectrum = cpl_vector_new (nwave);
229 cpl_table * fit_table = cpl_table_new (nregion);
235 cpl_polynomial * dispersion = cpl_polynomial_duplicate (dispersion0);
238 for (cpl_size reg = 0; reg < nregion ; reg ++ ) {
239 cpl_msg_info (cpl_func,
"Fit region %lld over %lld", reg+1, nregion);
242 const cpl_array * data = cpl_table_get_array (spectrum_table,
GRAVI_DATA[reg], 0);
243 for (cpl_size wave=0; wave < nwave; wave ++) {
244 cpl_vector_set (spectrum, wave, log (CPL_MAX (cpl_array_get (data, wave, NULL), 1.0)));
245 cpl_vector_set (input_spectrum, wave, cpl_array_get (data, wave, NULL));
253 int maxdeg = 3, nmaxima = 0, linelim = 99;
254 int maxite = 100 * maxdeg, maxfail = 10, maxcont = 10;
255 int hsize = (reg == 0 ? 40 : 5);
256 double pixtol = 1e-6, pixstep = 0.5, pxc = 0.0;
260 irplib_polynomial_find_1d_from_correlation_all (dispersion, maxdeg, spectrum,
262 (irplib_base_spectrum_model*)model,
263 &irplib_vector_fill_logline_spectrum,
265 hsize, maxite, maxfail, maxcont,
270 irplib_vector_fill_line_spectrum (model_spectrum, dispersion, (
void *)model);
271 cpl_polynomial_dump (dispersion, stdout);
274 cpl_array *
wave_data = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
275 for (cpl_size pix = 0; pix < nwave ; pix++)
276 cpl_array_set (
wave_data, pix, cpl_polynomial_eval_1d (dispersion, (
double)pix, NULL) * 1e-6);
279 minwave = CPL_MAX (cpl_array_get_min (
wave_data), minwave);
280 maxwave = CPL_MIN (cpl_array_get_max (
wave_data), maxwave);
288 cpl_array * tmp_array;
289 tmp_array = cpl_array_wrap_double (cpl_vector_get_data (input_spectrum), nwave);
290 cpl_table_set_array (fit_table,
"DATA", reg, tmp_array);
291 cpl_array_unwrap (tmp_array);
292 tmp_array = cpl_array_wrap_double (cpl_vector_get_data (model_spectrum), nwave);
293 cpl_table_set_array (fit_table,
"DATA_MODEL", reg, tmp_array);
294 cpl_array_unwrap (tmp_array);
295 tmp_array = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
296 for (cpl_size wave = 0; wave < nwave; wave ++)
297 cpl_array_set (tmp_array, wave, cpl_polynomial_eval_1d (dispersion, (
double)wave, NULL) * 1e-6);
298 cpl_table_set_array (fit_table,
"DATA_WAVE", reg, tmp_array);
299 cpl_array_delete (tmp_array);
305 FREE (cpl_table_delete, fit_table);
308 FREE (cpl_vector_delete, spectrum);
309 FREE (cpl_vector_delete, input_spectrum);
310 FREE (cpl_vector_delete, model_spectrum);
311 FREE (cpl_polynomial_delete, dispersion0);
312 FREE (cpl_polynomial_delete, dispersion);
357 cpl_ensure_code (phase_table, CPL_ERROR_NULL_INPUT);
358 cpl_ensure_code (name, CPL_ERROR_NULL_INPUT);
360 int nbase = 6, nclo = 4;
361 cpl_size nrow = cpl_table_get_nrow (phase_table) / nbase;
364 double * phase = cpl_table_get_data_double (phase_table, name);
369 cpl_matrix * rhs_matrix = cpl_matrix_new (nrow * nclo, 1);
370 cpl_matrix * model_matrix = cpl_matrix_new (nrow * nclo, nbase-1 + nclo);
372 for (cpl_size row = 0; row < nrow; row++) {
373 for (
int clo = 0; clo < nclo; clo++) {
379 double closure = phase[row*nbase+b0] +
380 phase[row*nbase+b1] -
382 cpl_matrix_set (rhs_matrix, row*nclo+clo, 0, closure);
386 cpl_matrix_set (model_matrix, row*nclo+clo, clo, 1.0);
390 if (b0 < nbase-1) cpl_matrix_set (model_matrix, row*nclo+clo, nclo+b0,
391 +phase[row*nbase+b0]);
392 if (b1 < nbase-1) cpl_matrix_set (model_matrix, row*nclo+clo, nclo+b1,
393 +phase[row*nbase+b1]);
394 if (b2 < nbase-1) cpl_matrix_set (model_matrix, row*nclo+clo, nclo+b2,
395 -phase[row*nbase+b2]);
401 cpl_matrix * res_matrix = cpl_matrix_solve_normal (model_matrix, rhs_matrix);
405 cpl_matrix * residual_matrix = cpl_matrix_product_create (model_matrix, res_matrix);
406 cpl_matrix_subtract (residual_matrix, rhs_matrix);
407 double rms_fit = cpl_matrix_get_stdev (residual_matrix);
410 for (
int base = 0; base < nbase - 1; base ++) {
411 double f = cpl_matrix_get (res_matrix, nclo + base, 0);
412 cpl_msg_info (cpl_func,
"correction factor f = 1 %+.20f", -1*f);
413 for (cpl_size row = 0; row < nrow; row++) phase[row*nbase+base] *= 1 - f;
417 const char * unit = cpl_table_get_column_unit (phase_table, name);
418 cpl_msg_info (cpl_func,
"residual stdev = %.20g [%s]", rms_fit, unit);
421 FREE (cpl_matrix_delete, residual_matrix);
422 FREE (cpl_matrix_delete, res_matrix);
423 FREE (cpl_matrix_delete, model_matrix);
424 FREE (cpl_matrix_delete, rhs_matrix);
427 return CPL_ERROR_NONE;
461 cpl_ensure (ft_table, CPL_ERROR_NULL_INPUT, NULL);
466 cpl_size nrow = cpl_table_get_nrow (ft_table) / nbase;
469 double * opd_sc = cpl_table_get_data_double (ft_table,
"OPD_SC");
470 double * opd_ft = cpl_table_get_data_double (ft_table,
"OPD");
471 double * phase_met = cpl_table_get_data_double (ft_table,
"PHASE_MET_FC");
475 cpl_size nrow_valid = 0;
476 for (cpl_size row = 0; row < nrow; row++)
if (opd_sc[row*nbase] != 0) nrow_valid++;
477 cpl_msg_info (cpl_func,
"nrow_valid = %lld", nrow_valid);
479 cpl_ensure (nrow_valid > 100, CPL_ERROR_ILLEGAL_INPUT, NULL);
483 cpl_matrix * rhs_matrix = cpl_matrix_new (nrow_valid*nbase, 1);
484 cpl_matrix * model_matrix = cpl_matrix_new (nrow_valid*nbase, nbase + 2);
486 for (
int base = 0; base < nbase; base++) {
487 cpl_size row_valid = 0;
489 for (cpl_size row=0; row<nrow; row++) {
490 if (opd_sc[row*nbase] == 0)
continue;
492 int idv = row_valid * nbase + base;
493 int id = row * nbase + base;
496 cpl_matrix_set (rhs_matrix, idv, 0, phase_met[
id] * lbd_met / CPL_MATH_2PI);
500 cpl_matrix_set (model_matrix, idv, 0, 1*opd_sc[
id]);
501 cpl_matrix_set (model_matrix, idv, 1, -1*opd_ft[
id]);
505 cpl_matrix_set (model_matrix, idv, 2 + base, 1.0);
513 cpl_matrix * res_matrix = cpl_matrix_solve_normal (model_matrix, rhs_matrix);
516 cpl_matrix * residual_matrix = cpl_matrix_product_create (model_matrix, res_matrix);
517 cpl_matrix_subtract (residual_matrix, rhs_matrix);
518 double rms_fit = cpl_matrix_get_stdev (residual_matrix);
522 cpl_msg_info (cpl_func,
"coeff SC = %.20g ", cpl_matrix_get (res_matrix, 0, 0));
523 cpl_msg_info (cpl_func,
"coeff FT = %.20g ", cpl_matrix_get (res_matrix, 1, 0));
524 cpl_msg_info (cpl_func,
"residual stdev = %.20g [m]", rms_fit);
527 cpl_vector * opd_coeff = cpl_vector_new(3);
528 cpl_vector_set (opd_coeff,
GRAVI_SC, cpl_matrix_get (res_matrix, 0, 0));
529 cpl_vector_set (opd_coeff,
GRAVI_FT, cpl_matrix_get (res_matrix, 1, 0));
530 cpl_vector_set (opd_coeff, 2, rms_fit);
533 FREE (cpl_matrix_delete, residual_matrix);
534 FREE (cpl_matrix_delete, res_matrix);
535 FREE (cpl_matrix_delete, model_matrix);
536 FREE (cpl_matrix_delete, rhs_matrix);
567 cpl_table * ft_table,
568 cpl_table * vismet_table,
573 cpl_ensure (spectrumsc_table, CPL_ERROR_NULL_INPUT, NULL);
574 cpl_ensure (ft_table, CPL_ERROR_NULL_INPUT, NULL);
575 cpl_ensure (vismet_table, CPL_ERROR_NULL_INPUT, NULL);
576 cpl_ensure (dit_sc>0, CPL_ERROR_ILLEGAL_INPUT, NULL);
578 int nbase = 6,
ntel = 4;
580 cpl_size nrow = cpl_table_get_nrow (spectrumsc_table);
581 cpl_size nrow_met = cpl_table_get_nrow (vismet_table) /
ntel;
582 cpl_size nrow_ft = cpl_table_get_nrow (ft_table) / nbase;
583 int * time_SC = cpl_table_get_data_int (spectrumsc_table,
"TIME");
587 int * time_MET = cpl_table_get_data_int (vismet_table,
"TIME");
588 double * phase_MET = cpl_table_get_data_double (vismet_table,
"PHASE_FC_DRS");
592 int * time_FT = cpl_table_get_data_int (ft_table,
"TIME");
593 double * opd_FT = cpl_table_get_data_double (ft_table,
"OPD");
597 cpl_table * guess_table = cpl_table_new (nrow * nbase);
601 for (
int base = 0; base < nbase; base++) {
606 cpl_size row_met = 0, row_ft = 0;
607 for (cpl_size row = 0 ; row < nrow ; row++) {
613 double opd_met = 0.0;
614 while (time_MET[row_met*
ntel] < (time_SC[row] + dit_sc/2.)) {
615 if ((time_MET[row_met*
ntel] > (time_SC[row] - dit_sc/2.)) && (row_met < nrow_met)) {
616 opd_met += phase_MET[row_met*
ntel+tel1] - phase_MET[row_met*
ntel+tel0];
621 if (row_met > nrow_met - 2) {
622 cpl_msg_warning (cpl_func,
"Not enough data to synchronise the MET with SC");
629 opd_met = opd_met * lbd_met / CPL_MATH_2PI / counter_met;
636 while (time_FT[row_ft*nbase+base] < (time_SC[row] + dit_sc/2.)) {
637 if ((time_FT[row_ft*nbase+base] > (time_SC[row] - dit_sc/2.)) && (row_ft < nrow_ft)) {
638 opd_ft += opd_FT[row_ft*nbase+base];
643 if (row_ft > nrow_ft - 2) {
644 cpl_msg_warning (cpl_func,
"Not enough data to synchronise the FT with SC");
651 opd_ft = opd_ft / counter_ft;
654 cpl_table_set (guess_table,
"OPD", row*nbase+base, opd_ft - opd_met);
662 cpl_msg_info (cpl_func,
"Remove the mean from opdguess table");
664 for (
int base = 0; base < 6; base++) {
699 cpl_table * detector_table,
700 cpl_table * guess_table)
706 cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
707 cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
711 cpl_msg_info (cpl_func,
"Solve linearly the ellipse X,Y");
714 cpl_msg_info (cpl_func,
"Solve non-linearly the ellipse X,Y");
718 cpl_size nrow = cpl_table_get_nrow (spectrum_table);
726 cpl_table ** oiwave_tables = cpl_calloc (npol,
sizeof (cpl_table *));
727 for (
int pol = 0; pol < npol; pol++) {
728 oiwave_tables[pol] = cpl_table_new (nwave);
729 cpl_table_new_column (oiwave_tables[pol],
"EFF_WAVE", CPL_TYPE_DOUBLE);
730 for (cpl_size wave = 0; wave < nwave; wave++) {
731 double value = 1.97e-6 + (2.48e-6 - 1.97e-6) / nwave * wave;
732 cpl_table_set (oiwave_tables[pol],
"EFF_WAVE", wave, value);
737 cpl_table * output_table = cpl_table_new (nbase * nrow);
742 for (cpl_size row = 0; row < nrow; row++) {
743 double value = cpl_table_get (spectrum_table,
"TIME", row, NULL);
744 for (
int base = 0; base < nbase; base++)
745 cpl_table_set (output_table,
"TIME", row*nbase+base, value);
749 for (
int base = 0; base < 6; base ++) {
752 cpl_vector * opd_guess = NULL;
754 opd_guess = cpl_vector_new (nrow);
755 for (cpl_size row = 0; row < nrow; row++) {
756 double value = cpl_table_get (guess_table,
"OPD", row*nbase+base, NULL);
757 cpl_vector_set (opd_guess, row, value);
764 cpl_vector * mean_opd;
766 oiwave_tables, opd_guess, base);
770 for (cpl_size row = 0; row < nrow; row++) {
771 double value = cpl_vector_get (mean_opd, row);
772 cpl_table_set (output_table,
"OPD", row*nbase+base, value);
776 FREE (cpl_vector_delete, mean_opd);
777 FREE (cpl_vector_delete, opd_guess);
781 FREELOOP (cpl_table_delete, oiwave_tables, npol);
823 cpl_table * met_table,
824 const cpl_parameterlist * parlist)
827 cpl_ensure_code (spectrum_data, CPL_ERROR_NULL_INPUT);
828 cpl_ensure_code (met_table, CPL_ERROR_NULL_INPUT);
860 cpl_msg_info (cpl_func,
"Compute OPD of FT from ellipse");
862 cpl_table * ft_table;
866 cpl_msg_info (cpl_func,
"Compute OPD of SC from ellipse");
868 cpl_table * guess_table;
871 cpl_table * sc_table;
874 FREE (cpl_table_delete, guess_table);
882 cpl_msg_info (cpl_func,
"Fit MET = a.SC - b.FT + c to get absolute modulation");
887 CPLCHECK_MSG (
"Cannot resample SC or MET at the FT frequency");
899 cpl_vector_get (coeff_vector, 2));
900 cpl_propertylist_set_comment (spectrum_header,
QC_PHASECHI2,
901 "chi2 of a.SC-b.FT+c=MET");
904 for (
int type_data = 0; type_data < 2; type_data++) {
905 double tmp = cpl_vector_get (coeff_vector, type_data);
906 cpl_propertylist_update_float (spectrum_header,
OPD_COEFF_SIGN(type_data), tmp);
907 cpl_propertylist_set_comment (spectrum_header,
OPD_COEFF_SIGN(type_data),
"wavelength correction");
911 if (cpl_vector_get (coeff_vector, 2) > 1e-7) {
915 if (cpl_vector_get (coeff_vector, 2) > 1.2e-7) {
916 cpl_msg_warning (cpl_func,
"*************************************************");
917 cpl_msg_warning (cpl_func,
"**** !!! residuals of the fit too high !!! ****");
918 cpl_msg_warning (cpl_func,
"**** Residuals are:%7.0f nm ****",cpl_vector_get (coeff_vector, 2)*1e9);
919 cpl_msg_warning (cpl_func,
"**** SC and RMN may be desynchronized ****");
920 cpl_msg_warning (cpl_func,
"**** (or out of the envelope in LOW) ****");
921 cpl_msg_warning (cpl_func,
"*************************************************");
929 double coeff_sc = cpl_vector_get (coeff_vector,
GRAVI_SC);
930 cpl_table_multiply_scalar (sc_table,
"OPD", coeff_sc);
932 double coeff_ft = cpl_vector_get (coeff_vector,
GRAVI_FT);
933 cpl_table_multiply_scalar (ft_table,
"OPD", coeff_ft);
935 FREE (cpl_vector_delete, coeff_vector);
936 CPLCHECK_MSG (
"Cannot correct OPDs from scaling coefficients");
962 return CPL_ERROR_NONE;
992 cpl_table * detector_table,
993 cpl_table * opd_table)
996 cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
997 cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
998 cpl_ensure (opd_table, CPL_ERROR_NULL_INPUT, NULL);
1004 cpl_table * wave_fibre = cpl_table_new (1);
1007 cpl_size nwave = cpl_table_get_column_depth (spectrum_table,
"DATA1");
1008 cpl_size n_region = cpl_table_get_nrow (detector_table);
1009 cpl_size nrow = cpl_table_get_nrow (spectrum_table);
1011 int npol = (n_region > 24 ? 2 : 1);
1017 for (
int pol = 0; pol < npol; pol++) {
1018 for (
int base = 0; base < nbase; base ++) {
1019 cpl_msg_info (cpl_func,
"Compute wave fibre for pol %i over %i, base %i over %i",
1020 pol+1, npol, base+1, nbase);
1027 if (iA<0 || iB<0 || iC<0 || iD<0){
1028 cpl_msg_warning (cpl_func,
"Don't found the A, B, C or D !!!");
1035 cpl_matrix * opd_matrix = cpl_matrix_new (1, nrow);
1036 cpl_vector * opd_vector = cpl_vector_new (nrow);
1037 for (cpl_size row = 0; row < nrow; row ++ ) {
1038 double value = cpl_table_get (opd_table,
"OPD", row*nbase+base, NULL);
1039 cpl_matrix_set (opd_matrix, 0, row, value);
1040 cpl_vector_set (opd_vector, row, value);
1045 cpl_array * wavelength = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1046 cpl_array * wavechi2 = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1047 cpl_array_fill_window (wavelength, 0, nwave, 0.0);
1048 cpl_array_fill_window (wavechi2, 0, nwave, 1e10);
1053 for (cpl_size wave = 0; wave < nwave; wave++) {
1055 cpl_vector * vector_T = NULL, * vector_X, * vector_Y;
1060 cpl_vector_subtract (vector_X, vector_T);
1061 FREE (cpl_vector_delete, vector_T);
1066 cpl_vector_subtract (vector_Y, vector_T);
1067 FREE (cpl_vector_delete, vector_T);
1078 FREE (cpl_vector_delete, vector_X);
1079 FREE (cpl_vector_delete, vector_Y);
1080 FREE (cpl_vector_delete, envelope_vector);
1083 if (phase == NULL) {
1084 cpl_msg_warning (cpl_func,
"Cannot compute wave for %lld and base %d", wave, base);
1089 cpl_vector_multiply_scalar (phase, phi_sign);
1092 double lbd_channel = 1.95e-6 + (2.46e-6 - 1.95e-6) / nwave * wave;
1097 const cpl_size mindeg = 0, maxdeg = 1;
1098 cpl_polynomial * fit_slope = cpl_polynomial_new (1);
1099 cpl_polynomial_fit (fit_slope, opd_matrix, NULL, phase, NULL, CPL_FALSE, &mindeg, &maxdeg);
1102 cpl_vector * residuals = cpl_vector_new (nwave);
1104 cpl_vector_fill_polynomial_fit_residual (residuals, phase, NULL, fit_slope, opd_matrix, &rechisq);
1108 char gnuplot_str[200];
1109 sprintf (gnuplot_str,
"set title 'Wavelength (base %d)'; set xlabel 'Phase [rad]'; set ylabel 'OPD (m)';", base);
1110 cpl_plot_vector (gnuplot_str, NULL, NULL, opd_vector);
1111 sprintf (gnuplot_str,
"set title 'Wavelength residuals (base %d)'; set xlabel 'Phase [rad]'; set ylabel 'OPD (m)';", base);
1112 cpl_plot_vector (gnuplot_str, NULL, NULL, residuals);
1119 const cpl_size pow_slope = 1;
1120 double slope = cpl_polynomial_get_coeff (fit_slope, &pow_slope);
1123 if (slope < 0.0 && wave == 0) {
1124 cpl_msg_warning (cpl_func,
"Negative wavelength!! "
1125 "Report to DRS team");
1129 cpl_array_set (wavechi2, wave, sqrt(rechisq));
1130 cpl_array_set (wavelength, wave, CPL_MATH_2PI / fabs(slope));
1133 cpl_vector_delete (phase);
1134 cpl_vector_delete (residuals);
1135 cpl_polynomial_delete (fit_slope);
1139 cpl_matrix_delete (opd_matrix);
1140 cpl_vector_delete (opd_vector);
1144 cpl_table_new_column_array (wave_fibre, name, CPL_TYPE_DOUBLE, nwave);
1145 cpl_table_set_column_unit (wave_fibre, name,
"m");
1146 cpl_table_set_array (wave_fibre, name, 0, wavelength);
1147 cpl_array_delete (wavelength);
1151 cpl_table_new_column_array (wave_fibre, name, CPL_TYPE_DOUBLE, nwave);
1152 cpl_table_set_array (wave_fibre, name, 0, wavechi2);
1153 cpl_array_delete (wavechi2);
1183 cpl_table * detector_table,
1185 cpl_size fullstartx,
1188 double * rms_residuals)
1191 cpl_ensure (wavefibre_table, CPL_ERROR_NULL_INPUT, NULL);
1192 cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
1199 cpl_size n_region = cpl_table_get_nrow (detector_table);
1200 int npol = n_region > 24 ? 2 : 1;
1201 cpl_size nwave = cpl_table_get_column_depth (wavefibre_table, npol > 1 ?
"BASE_12_S" :
"BASE_12_C");
1208 cpl_vector * odd_index = cpl_vector_new (nwave);
1209 for (
int i = fullstartx; i < fullstartx + nwave; i++) {
1210 if (nwave >
GRAVI_LBD_FTSC) cpl_vector_set (odd_index, i - fullstartx, ((i/64)%2 == 0) ? 0 : 1);
1211 else cpl_vector_set (odd_index, i - fullstartx, 0);
1217 cpl_polynomial ** coef_poly = cpl_calloc (npol,
sizeof (cpl_polynomial*));
1220 for (
int pol = 0; pol < npol; pol++) {
1224 cpl_vector * coord_X = cpl_vector_new (nbase * nwave);
1225 cpl_vector * coord_Y = cpl_vector_new (nbase * nwave);
1227 cpl_vector * all_wavelength = cpl_vector_new (nbase * nwave);
1228 cpl_vector * all_wavechi2 = cpl_vector_new (nbase * nwave);
1229 cpl_vector * all_valid = cpl_vector_new (nbase * nwave);
1235 for (
int base = 0; base < nbase; base ++) {
1245 cpl_array * wavelength = cpl_table_get_data_array (wavefibre_table, name)[0];
1248 cpl_array * wavechi2 = cpl_table_get_data_array (wavefibre_table, name)[0];
1251 for (cpl_size wave = 0; wave < nwave; wave++) {
1252 cpl_size pos = base * nwave + wave;
1256 double wave_value = cpl_array_get (wavelength, wave, &nv);
1257 double chi2_value = cpl_array_get (wavechi2, wave, &nv);
1260 cpl_vector_set (all_valid, pos, 1);
1264 if ((chi2_value > M_PI_4 ||
1267 cpl_vector_set (all_valid, pos, 0);
1269 else if (nwave > 100) {
1271 if ((chi2_value > M_PI_4 ||
1274 cpl_vector_set (all_valid, pos, 0);
1278 if ((chi2_value > M_PI_4 ||
1282 wave == 0 || wave == nwave-1))
1283 cpl_vector_set (all_valid, pos, 0);
1287 if ((chi2_value > M_PI_4 ||
1289 wave_value < 1.99e-6 ||
1290 wave_value > 2.5e-6))
1291 cpl_vector_set (all_valid, pos, 0);
1295 cpl_vector_set (all_wavelength, pos, wave_value);
1296 cpl_vector_set (all_wavechi2, pos, chi2_value);
1299 cpl_vector_set (coord_X, pos, (
double)(iA + iB + iC + iD) / 4.);
1303 cpl_vector_set (coord_Y, pos, wave + cpl_vector_get (odd_index, wave)*0.15);
1313 cpl_size nvalid = cpl_vector_get_sum (all_valid);
1316 nbase*nwave - nvalid, nbase * nwave);
1318 cpl_vector * vector = cpl_vector_new (nvalid);
1319 cpl_matrix * matrix = cpl_matrix_new (2, nvalid);
1321 for (cpl_size c = 0, i = 0 ; i < nwave * nbase; i ++) {
1322 if (!cpl_vector_get (all_valid, i))
continue;
1323 cpl_vector_set (vector, c, cpl_vector_get (all_wavelength, i));
1324 cpl_matrix_set (matrix, 0, c, cpl_vector_get (coord_X, i));
1325 cpl_matrix_set (matrix, 1, c, cpl_vector_get (coord_Y, i));
1335 cpl_size deg2d[2] = {2, 3};
1336 if ( (nwave < 20) && (nwave > 8) ) {deg2d[0] = 2; deg2d[1] = 7;}
1337 deg2d[0] = spatial_order;
1338 deg2d[1] = spectral_order;
1340 cpl_msg_info (cpl_func,
"Fit a 2d polynomial {%lli..%lli} to the wavelengths map", deg2d[0], deg2d[1]);
1342 cpl_polynomial * fit2d = cpl_polynomial_new (2);
1343 cpl_polynomial_fit (fit2d, matrix, NULL, vector, NULL, CPL_TRUE, NULL, deg2d);
1344 coef_poly[pol] = fit2d;
1351 double rechisq = 0.0;
1352 cpl_vector * residuals = cpl_vector_new (nvalid);
1353 cpl_vector_fill_polynomial_fit_residual (residuals, vector, NULL, fit2d, matrix, &rechisq);
1354 *rms_residuals += cpl_vector_get_stdev(residuals)/npol;
1355 FREE (cpl_vector_delete, residuals);
1359 FREE (cpl_matrix_delete, matrix);
1360 FREE (cpl_vector_delete, vector);
1361 FREE (cpl_vector_delete, all_wavelength);
1362 FREE (cpl_vector_delete, all_wavechi2);
1363 FREE (cpl_vector_delete, all_valid);
1364 FREE (cpl_vector_delete, coord_X);
1365 FREE (cpl_vector_delete, coord_Y);
1374 cpl_table * wavedata_table = cpl_table_new (1);
1375 cpl_vector * pos = cpl_vector_new (2);
1376 cpl_array * value = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1378 for (cpl_size region = 0 ; region < n_region; region ++) {
1382 for (cpl_size wave = 0; wave < nwave; wave ++) {
1385 cpl_vector_set (pos, 0, region);
1386 cpl_vector_set (pos, 1, wave + cpl_vector_get (odd_index, wave)*0.15);
1388 double result = cpl_polynomial_eval (coef_poly[pol], pos);
1389 cpl_array_set (value, wave, result);
1393 double previous_wave = cpl_array_get(value, nwave/2, NULL);
1394 for (cpl_size wave_loop = nwave/2 ; wave_loop >= 0 ; wave_loop --){
1395 if (previous_wave < cpl_array_get(value, wave_loop, NULL))
1396 cpl_array_set(value, wave_loop, previous_wave);
1397 else previous_wave = cpl_array_get(value, wave_loop, NULL);
1400 previous_wave = cpl_array_get(value, nwave/2, NULL);
1401 for (cpl_size wave_loop = nwave/2 ; wave_loop < nwave ; wave_loop ++){
1402 if (previous_wave > cpl_array_get(value, wave_loop, NULL))
1403 cpl_array_set(value, wave_loop, previous_wave);
1404 else previous_wave = cpl_array_get(value, wave_loop, NULL);
1409 cpl_table_new_column_array (wavedata_table, data_x, CPL_TYPE_DOUBLE, nwave);
1410 cpl_table_set_array (wavedata_table, data_x, 0, value);
1415 FREE (cpl_vector_delete, pos);
1416 FREE (cpl_array_delete, value);
1417 FREELOOP (cpl_polynomial_delete, coef_poly, npol);
1418 FREE (cpl_vector_delete, odd_index);
1421 return wavedata_table;
1433 cpl_table * weight_individual_table,
1434 cpl_table * wave_fitted_table,
1435 cpl_table * opd_table,
1436 cpl_table * spectrum_table,
1437 cpl_table * detector_table,
1438 double n0,
double n1,
double n2)
1443 cpl_ensure (wave_individual_table, CPL_ERROR_NULL_INPUT, NULL);
1444 cpl_ensure (weight_individual_table, CPL_ERROR_NULL_INPUT, NULL);
1445 cpl_ensure (wave_fitted_table, CPL_ERROR_NULL_INPUT, NULL);
1446 cpl_ensure (opd_table, CPL_ERROR_NULL_INPUT, NULL);
1447 cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
1448 cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
1451 cpl_size nwave = cpl_table_get_column_depth (spectrum_table,
"DATA1");
1452 cpl_size n_region = cpl_table_get_nrow (detector_table);
1453 cpl_size nrow = cpl_table_get_nrow (spectrum_table);
1454 int npol = (n_region > 24 ? 2 : 1);
1456 cpl_size nwave_ref=3000;
1457 if (nwave<10) nwave_ref=600;
1467 cpl_array * wave_individual_array = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1468 cpl_array * weight_individual_array = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1469 cpl_matrix * data_flux_matrix = cpl_matrix_new (nrow, nwave);
1470 cpl_matrix * vis_to_flux_matrix = cpl_matrix_new (nrow, 3);
1471 cpl_matrix * signal_matrix = cpl_matrix_new (nwave, nwave_ref);
1472 cpl_matrix * residual_matrix = cpl_matrix_new (nwave, nwave_ref);
1473 cpl_array * wave_reference_array = cpl_array_new (nwave_ref, CPL_TYPE_DOUBLE);
1477 cpl_matrix_fill_column(vis_to_flux_matrix,1,2);
1478 for (cpl_size wave_ref = 0; wave_ref < nwave_ref; wave_ref+=1)
1481 double wave_value=1.95e-6+wave_ref*0.6e-6/((double) nwave_ref);
1482 cpl_array_set_double(wave_reference_array,wave_ref,wave_value);
1485 CPLCHECK_NUL (
"Cannot initialize arrays for wavelength fit");
1488 for (cpl_size region = 0 ; region < n_region; region ++) {
1493 cpl_msg_info_overwritable (cpl_func,
"Least square fitting of wavelength for region %s", data_x);
1496 for (cpl_size row = 0; row < nrow; row ++ ) {
1497 cpl_array * flux_array= cpl_table_get_data_array(spectrum_table,data_x)[row];
1498 for (cpl_size wave = 0; wave < nwave; wave ++) {
1499 cpl_matrix_set (data_flux_matrix, row, wave, cpl_array_get(flux_array,wave, NULL));
1504 for (cpl_size wave_ref = 0; wave_ref < nwave_ref; wave_ref+=1)
1507 double wave_value=cpl_array_get(wave_reference_array,wave_ref,NULL);
1509 for (cpl_size row = 0; row < nrow; row ++ ) {
1510 double opd = cpl_table_get (opd_table,
"OPD", row*nbase+base, NULL);
1511 double coherence_loss=1;
1512 if (fabs(opd) > 1e-9)
1516 coherence_loss=sin(opd*19500)/(opd*19500);
1519 cpl_matrix_set(vis_to_flux_matrix,row,0,cos(opd*6.28318530718/wave_value)*coherence_loss);
1520 cpl_matrix_set(vis_to_flux_matrix,row,1,sin(opd*6.28318530718/wave_value)*coherence_loss);
1523 cpl_matrix * coef_vis = cpl_matrix_solve_normal(vis_to_flux_matrix,data_flux_matrix);
1524 cpl_matrix * data_flux_fit = cpl_matrix_product_create(vis_to_flux_matrix,coef_vis);
1525 cpl_matrix * residuals_fit = cpl_matrix_duplicate(data_flux_fit);
1526 cpl_matrix_subtract (residuals_fit,data_flux_matrix);
1528 for (cpl_size wave = 0; wave < nwave; wave ++ ) {
1529 cpl_matrix * temp_matrix = cpl_matrix_extract_column (data_flux_fit, wave);
1530 cpl_matrix_set(signal_matrix,wave,wave_ref,cpl_matrix_get_stdev(temp_matrix));
1531 FREE (cpl_matrix_delete, temp_matrix);
1532 cpl_matrix * temp_matrix2 = cpl_matrix_extract_column (residuals_fit, wave);
1533 cpl_matrix_set(residual_matrix,wave,wave_ref,cpl_matrix_get_stdev(temp_matrix2));
1534 FREE (cpl_matrix_delete, temp_matrix2);
1537 FREE (cpl_matrix_delete, coef_vis);
1538 FREE (cpl_matrix_delete, data_flux_fit);
1539 FREE (cpl_matrix_delete, residuals_fit);
1541 CPLCHECK_NUL (
"Cannot do Matrix inversion to calculate optimum wavelength");
1546 cpl_size wave_ref=1;
1547 cpl_size discarded=1;
1548 for (cpl_size wave = 0; wave < nwave; wave ++ ) {
1550 cpl_matrix * chi2_extract=cpl_matrix_extract_row(residual_matrix,wave);
1552 cpl_matrix_get_minpos(chi2_extract,&discarded, &wave_ref );
1554 double wave_value = cpl_array_get(wave_reference_array, wave_ref, NULL );
1555 double weight_value = cpl_matrix_get(signal_matrix, wave , wave_ref )/(0.1+cpl_matrix_get(residual_matrix, wave, wave_ref ));
1557 cpl_array_set_double (wave_individual_array, wave, wave_value);
1558 cpl_array_set_double (weight_individual_array, wave, weight_value);
1560 FREE (cpl_matrix_delete, chi2_extract);
1565 cpl_table_new_column_array (wave_individual_table, data_x, CPL_TYPE_DOUBLE, nwave);
1566 cpl_table_set_array (wave_individual_table, data_x, 0, wave_individual_array);
1567 cpl_table_new_column_array (weight_individual_table, data_x, CPL_TYPE_DOUBLE, nwave);
1568 cpl_table_set_array (weight_individual_table, data_x, 0, weight_individual_array);
1569 cpl_table_new_column_array (wave_fitted_table, data_x, CPL_TYPE_DOUBLE, nwave);
1570 cpl_table_set_array (wave_fitted_table, data_x, 0, wave_individual_array);
1573 CPLCHECK_NUL (
"Cannot get individual wavelength for each pixel");
1575 FREE (cpl_array_delete ,wave_individual_array);
1576 FREE (cpl_array_delete ,weight_individual_array);
1577 FREE (cpl_array_delete ,wave_reference_array);
1578 FREE (cpl_matrix_delete ,data_flux_matrix);
1579 FREE (cpl_matrix_delete ,vis_to_flux_matrix);
1580 FREE (cpl_matrix_delete ,signal_matrix);
1581 FREE (cpl_matrix_delete ,residual_matrix);
1583 cpl_msg_info (cpl_func,
"Now fitting polynomials on wavelength channels");
1585 cpl_matrix * coef_to_wave = cpl_matrix_new (n_region / npol,5);
1586 cpl_matrix * coef_to_wave_weight = cpl_matrix_new (n_region / npol,n_region / npol);
1587 cpl_matrix * wavelength = cpl_matrix_new(n_region / npol,1);
1590 for (cpl_size region = 0 ; region < n_region/ npol; region ++)
1592 double mean_region = region - (n_region/npol-1)*0.5;
1593 cpl_matrix_set (coef_to_wave, region, 0, 1);
1594 cpl_matrix_set (coef_to_wave, region, 1, mean_region);
1595 cpl_matrix_set (coef_to_wave, region, 2, mean_region*mean_region);
1596 cpl_matrix_set (coef_to_wave, region, 3, mean_region*mean_region*mean_region);
1597 cpl_matrix_set (coef_to_wave, region, 4, mean_region*mean_region*mean_region*mean_region);
1600 for (
int pol = 0; pol < npol; pol++) {
1602 cpl_msg_info (cpl_func,
"Looping for polyfit now, with pol: %i",(
int) pol);
1604 for (cpl_size wave = 0; wave < nwave; wave ++) {
1607 for (cpl_size region = 0 ; region < n_region/ npol; region ++) {
1613 const cpl_array * wave_array = cpl_table_get_array (wave_individual_table, data_x, 0);
1614 const cpl_array * weight_array = cpl_table_get_array (weight_individual_table, data_x, 0);
1617 cpl_matrix_set (wavelength, region, 0, cpl_array_get(wave_array,wave,NULL));
1618 double weight_value=cpl_array_get(weight_array,wave,NULL);
1619 cpl_matrix_set (coef_to_wave_weight, region, region, weight_value*weight_value);
1624 cpl_matrix * coef_to_wave2 = cpl_matrix_product_create(coef_to_wave_weight,coef_to_wave);
1625 cpl_matrix * wavelength2 = cpl_matrix_product_create(coef_to_wave_weight,wavelength);
1628 cpl_matrix * coeff = cpl_matrix_solve_normal(coef_to_wave2, wavelength2);
1629 cpl_matrix * wavelength_fitted = cpl_matrix_product_create(coef_to_wave, coeff);
1631 CPLCHECK_NUL (
"Cannot do Matrix inversion to calculate optimum polynomial for wavelength");
1634 for (cpl_size region = 0 ; region < n_region/ npol; region ++) {
1638 cpl_array * wave_array = cpl_table_get_data_array (wave_fitted_table, data_x)[0];
1640 cpl_array_set_double(wave_array,wave,cpl_matrix_get(wavelength_fitted,region,0));
1644 FREE (cpl_matrix_delete ,coef_to_wave2);
1645 FREE (cpl_matrix_delete ,wavelength2);
1646 FREE (cpl_matrix_delete ,coeff);
1647 FREE (cpl_matrix_delete ,wavelength_fitted);
1651 FREE (cpl_matrix_delete ,coef_to_wave);
1652 FREE (cpl_matrix_delete ,coef_to_wave_weight);
1653 FREE (cpl_matrix_delete ,wavelength);
1655 CPLCHECK_NUL (
"Cannot fit individual wavelength with 3rd order polynomial");
1657 cpl_msg_info (cpl_func,
"Correcting for wavelength error");
1661 for (cpl_size region = 0 ; region < n_region; region ++)
1664 cpl_array * wavelength_fitted = cpl_table_get_data_array (wave_fitted_table, data_x)[0];
1665 cpl_size nwave_fitted = cpl_array_get_size (wavelength_fitted);
1666 for (cpl_size wave = 0 ; wave < nwave_fitted ; wave ++ ) {
1668 double result = cpl_array_get (wavelength_fitted, wave, NULL);
1670 cpl_array_set (wavelength_fitted, wave, result * (n0 + n1*d_met + n2*d_met*d_met));
1673 for (cpl_size wave = nwave_fitted/2 ; wave < nwave_fitted-1 ; wave ++ ) {
1674 double result = cpl_array_get (wavelength_fitted, wave, NULL);
1675 double result2 = cpl_array_get (wavelength_fitted, wave+1, NULL);
1676 if (result2<result+2e-10) {
1677 result2=result+2e-10;
1678 cpl_array_set (wavelength_fitted, wave+1, result2);
1681 for (cpl_size wave = nwave_fitted/2 ; wave > 0 ; wave -- ) {
1682 double result = cpl_array_get (wavelength_fitted, wave, NULL);
1683 double result2 = cpl_array_get (wavelength_fitted, wave-1, NULL);
1684 if (result2>result-2e-10) {
1685 result2=result-2e-10;
1686 cpl_array_set (wavelength_fitted, wave-1, result2);
1694 return wave_fitted_table;
1710 double n0,
double n1,
double n2)
1713 cpl_ensure_code (wave_fibre, CPL_ERROR_NULL_INPUT);
1719 int npol = cpl_table_has_column (wave_fibre,
"BASE_12_P") ? 2 : 1;
1722 for (
int pol = 0; pol < npol; pol ++) {
1723 for (
int base = 0; base < nbase; base ++) {
1727 cpl_array * wavelength = cpl_table_get_data_array (wave_fibre, name)[0];
1731 cpl_size nwave = cpl_array_get_size (wavelength);
1732 for (cpl_size wave = 0 ; wave < nwave ; wave ++ ) {
1735 double result = cpl_array_get (wavelength, wave, NULL);
1737 cpl_array_set (wavelength, wave, result * (n0 + n1*d_met + n2*d_met*d_met));
1744 return CPL_ERROR_NONE;
1758 cpl_ensure_code (vis_data, CPL_ERROR_NULL_INPUT);
1761 cpl_propertylist * oiwave_header = NULL;
1762 cpl_table * oiwave_table = NULL;
1766 for (
int type_data = 0; type_data < ntype_data ; type_data ++) {
1770 for (
int pol = 0 ; pol<npol ; pol++) {
1790 return CPL_ERROR_NONE;
1809 cpl_table * wavefibre_table,
1810 cpl_table * profile_table,
1811 cpl_table * detector_table)
1814 cpl_ensure (wavedata_table, CPL_ERROR_NULL_INPUT, NULL);
1815 cpl_ensure (wavefibre_table, CPL_ERROR_NULL_INPUT, NULL);
1816 cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
1817 cpl_ensure (profile_table, CPL_ERROR_NULL_INPUT, NULL);
1822 cpl_size n_region = cpl_table_get_nrow (detector_table);
1823 int npol = (n_region > 24 ? 2 : 1);
1828 cpl_size sizex = cpl_table_get_column_dimension (profile_table,
"DATA1", 0);
1829 cpl_size sizey = cpl_table_get_column_dimension (profile_table,
"DATA1", 1);
1831 cpl_image * profilesum_image = cpl_image_new (sizex, sizey, CPL_TYPE_DOUBLE);
1832 cpl_image_fill_window (profilesum_image, 1, 1, sizex, sizey, 0.0);
1834 cpl_image * wave_image = cpl_image_new (sizex, sizey, CPL_TYPE_DOUBLE);
1835 cpl_image_fill_window (wave_image, 1, 1, sizex, sizey, 0.0);
1837 cpl_image * realwave_image = cpl_image_new (sizex, sizey, CPL_TYPE_DOUBLE);
1838 cpl_image_fill_window (realwave_image, 1, 1, sizex, sizey, 0.0);
1843 for (cpl_size region = 0 ; region < n_region; region ++) {
1847 cpl_image * profile_image = cpl_imagelist_get (profile_imglist, 0);
1852 cpl_image_add (profilesum_image, profile_image);
1857 const cpl_array * wavelength;
1858 wavelength = cpl_table_get_array (wavedata_table,
GRAVI_DATA[region], 0);
1861 for (cpl_size x = 0; x < sizex; x ++){
1862 for (cpl_size y = 0; y < sizey; y ++){
1863 if (cpl_image_get (profile_image, x+1, y+1, &nv) > 0.01)
1864 cpl_image_set (wave_image, x+1, y+1,
1865 cpl_array_get (wavelength, x, NULL));
1876 wavelength = cpl_table_get_array (wavefibre_table, name, 0);
1879 for (cpl_size x = 0; x < sizex; x ++){
1880 for (cpl_size y = 0; y < sizey; y ++){
1881 if (cpl_image_get (profile_image, x+1, y+1, &nv) > 0.01)
1882 cpl_image_set (realwave_image, x+1, y+1,
1883 cpl_array_get (wavelength, x, NULL));
1892 cpl_imagelist * testwave_imglist = cpl_imagelist_new ();
1893 cpl_imagelist_set (testwave_imglist, wave_image, 0);
1894 cpl_imagelist_set (testwave_imglist, profilesum_image, 1);
1895 cpl_imagelist_set (testwave_imglist, realwave_image, 2);
1898 return testwave_imglist;
1916 cpl_ensure_code (wave_map, CPL_ERROR_NULL_INPUT);
1923 for (
int type_data = 0; type_data < 2; type_data ++ ) {
1928 cpl_size n_region = cpl_table_get_nrow (detector_table);
1939 double minwave = -1e10;
1940 double maxwave = 1e10;
1942 for (
int region = 0; region < n_region; region++) {
1945 minwave = CPL_MAX (minwave, cpl_array_get_min (wavelength));
1946 maxwave = CPL_MIN (maxwave, cpl_array_get_max (wavelength));
1963 int qc_reg[3] = {0,23,47};
1964 double qc_wave[2] = {2.099184e-6, 2.313952e-6};
1967 for (
int reg = 0; reg < (n_region > 24 ? 3 : 2); reg++) {
1968 cpl_size region = qc_reg[reg];
1971 cpl_size nwave = cpl_array_get_size (wavelength);
1972 const double * wave_tab = cpl_array_get_data_double_const (wavelength);
1976 for (
int iqc = 0 ; iqc < 2 ; iqc++) {
1978 if ( wave_tab[0] < wave_tab[nwave-1]) {
while (wave_tab[l2] < qc_wave[iqc]) l2 ++;}
1979 else {
while (wave_tab[l2] > qc_wave[iqc]) l2 ++;}
1981 if (l2-1 < 0 || l2 > nwave-1) {
1982 cpl_msg_error (cpl_func,
"Cannot find the QC position for lbd=%g", qc_wave[iqc]);
1987 double qc_pos = 0.0;
1988 qc_pos = fullstartx + (l2-1) + (qc_wave[iqc] - wave_tab[l2-1]) / (wave_tab[l2] - wave_tab[l2-1]);
1990 sprintf (name,
"ESO QC REFWAVE%i", iqc+1);
1992 cpl_propertylist_set_comment (wave_header, name,
"[m] value of ref wave");
1994 sprintf (name,
"ESO QC REFPOS%i %s%lli", iqc+1,
GRAVI_TYPE(type_data),region+1);
1996 cpl_propertylist_set_comment (wave_header, name,
"[pix] position of ref wave");
1998 cpl_msg_info (cpl_func,
"%s = %f [pix] for %e [m]", name, qc_pos, qc_wave[iqc]);
2014 cpl_imagelist * testwave_imglist;
2016 profile_table, detector_table);
2022 "TEST_WAVE", testwave_imglist);
2027 return CPL_ERROR_NONE;
2062 int type_data,
const cpl_parameterlist * parlist,
2066 cpl_ensure_code (wave_map, CPL_ERROR_NULL_INPUT);
2067 cpl_ensure_code (spectrum_data, CPL_ERROR_NULL_INPUT);
2069 CPL_ERROR_ILLEGAL_INPUT);
2077 cpl_propertylist_append (wave_header, raw_header);
2094 cpl_table * wavefibre_table;
2095 wavefibre_table =
gravi_wave_fibre (spectrum_table, detector_table, opd_table);
2104 double n0 = 1.0, n1 = -0.0165448, n2 = 0.00256002;
2105 cpl_msg_info (cpl_func,
"Rescale wavelengths with dispersion (%g,%g,%g)",n0,n1,n2);
2107 cpl_propertylist_update_string (wave_header,
"ESO QC WAVE_CORR",
"lbd*(N0+N1*(lbd-lbd0)/lbd0+N2*(lbd-lbd0)^2/lbd0^2)");
2126 cpl_table * wavedata_table;
2127 int spatial_order=2;
2128 int spectral_order=3;
2129 double rms_residuals;
2132 cpl_msg_info (cpl_func,
"Option force-waveFT-equal applied");
2143 fullstartx, spatial_order, spectral_order, &rms_residuals);
2148 double rms_fit=cpl_propertylist_get_double (raw_header,
QC_PHASECHI2);
2150 cpl_propertylist_set_comment (wave_header,
QC_CHI2WAVE(type_data),
"[nm]rms a.SC-b.FT+c=MET");
2161 cpl_table * wave_individual_table = cpl_table_new (1);
2162 cpl_table * weight_individual_table = cpl_table_new (1);
2163 cpl_table * wave_fitted_table = cpl_table_new (1);
2166 weight_individual_table,
2176 "WAVE_INDIV_SC", wave_individual_table);
2178 "WAVE_WEIGHT_SC", weight_individual_table);
2180 "WAVE_FITTED_SC", wavedata_table);
2190 cpl_msg_info (cpl_func,
"Add WAVE_FIBRE and WAVE_DATA in wave_map");
2206 return CPL_ERROR_NONE;
#define GRAVI_PRIMARY_HDR_EXT
typedefCPL_BEGIN_DECLS struct _gravi_data_ gravi_data
#define gravi_data_get_spectrum_data(data, type)
#define gravi_data_get_wave_data_plist(data, type)
#define gravi_data_get_wave_data(data, type)
#define gravi_data_get_spectrum_data_plist(data, type)
#define gravi_data_get_header(data)
#define gravi_data_get_oi_wave(data, type, pol, npol)
#define gravi_data_get_imaging_detector(data, type)
#define gravi_data_get_wave_fibre(data, type)
#define gravi_data_get_oi_wave_plist(data, type, pol, npol)
#define USE_LINEAR_ELLIPSE
cpl_msg_info(cpl_func, "Compute WAVE_SCAN for %s", GRAVI_TYPE(type_data))
cpl_propertylist_update_double(header, "ESO QC MINWAVE SC", cpl_propertylist_get_double(plist, "ESO QC MINWAVE SC"))
gravi_vis_create_met_ft(phase_ft, vis_met)
#define GRAVI_OI_VIS_MET_EXT
#define GRAVI_IMAGING_DETECTOR_FT_EXT
#define QC_MINWAVE_UM(type)
#define GRAVI_PROFILE_DATA_EXT
#define GRAVI_WAVE_FIBRE_EXT(type)
#define QC_MAXWAVE_UM(type)
#define OPD_COEFF_SIGN(type)
#define QC_CHI2WAVE(type)
#define GRAVI_WAVE_DATA_EXT(type)
#define QC_RMS_RESIDUALS_UM(type)
#define GRAVI_IMAGING_DETECTOR_SC_EXT
#define QC_RMS_RESIDUALS(type)
#define GRAVI_POLAR(pol, npol)
#define gravi_spectrum_get_npol(table)
#define gravi_msg_function_exit(flag)
#define FREE(function, variable)
#define CPLCHECK_NUL(msg)
#define gravi_msg_function_start(flag)
#define CPLCHECK_MSG(msg)
#define FREELOOP(function, variable, n)
double gravi_table_get_column_mean(cpl_table *table, const char *name, int base, int nbase)
cpl_error_code gravi_table_new_column(cpl_table *table, const char *name, const char *unit, cpl_type type)
cpl_error_code gravi_vector_unwrap_with_guess(cpl_vector *vector, cpl_vector *ref, double ref_to_phase)
Unwrap a phase vector following a guess vector. The difference is actually unwrap and shall thus be s...
cpl_imagelist * gravi_imagelist_wrap_column(cpl_table *table_data, const char *data_x)
Wrap a column array of a table into an imagelist.
cpl_error_code gravi_imagelist_unwrap_images(cpl_imagelist *imglist)
Unwrap an imagelist an all its images.
cpl_vector * gravi_table_get_vector(cpl_table *spectrum_data, cpl_size index, const char *regname)
Create a vector from the row index of the column regname.
cpl_error_code gravi_table_add_scalar(cpl_table *table, const char *name, int base, int nbase, double value)
Multply scalar or array column by scalar.
cpl_error_code gravi_table_new_column_array(cpl_table *table, const char *name, const char *unit, cpl_type type, cpl_size size)
cpl_error_code gravi_data_add_cube(gravi_data *self, cpl_propertylist *plist, const char *extname, cpl_imagelist *imglist)
Add an IMAGE (imagelist) extension in gravi_data.
cpl_propertylist * gravi_data_get_plist(gravi_data *self, const char *extname)
Get the propertylist from EXTNAME.
cpl_error_code gravi_data_add_table(gravi_data *self, cpl_propertylist *plist, const char *extname, cpl_table *table)
Add a BINTABLE extension in gravi_data.
cpl_table * gravi_data_get_table(gravi_data *self, const char *extname)
Return a pointer on a table extension by its EXTNAME.
cpl_error_code gravi_data_copy_ext(gravi_data *output, gravi_data *input, const char *name)
Copy extensions from one data to another.
int gravi_param_get_bool(const cpl_parameterlist *parlist, const char *name)
cpl_vector * gravi_ellipse_meanopd_create(cpl_table *spectrum_table, cpl_table *detector_table, cpl_table **oiwave_tables, cpl_vector *guess_vector, int base)
Compute the OPD modulation of a baseline from spectrum.
cpl_vector * gravi_ellipse_phase_create(cpl_vector *vectCA, cpl_vector *vectDB, cpl_vector *envelope)
Compute the phase atan{X',Y'}, unwraped from first sample.
cpl_table * gravi_metrology_create(cpl_table *metrology_table, cpl_propertylist *header)
Create the VIS_MET table.
cpl_error_code gravi_metrology_drs(cpl_table *metrology_table, cpl_table *vismet_table, cpl_propertylist *header, const cpl_parameterlist *parlist)
Fill the VIS_MET table with the DRS algorithm.
cpl_error_code gravi_metrology_tac(cpl_table *metrology_table, cpl_table *vismet_table, cpl_propertylist *header)
Compute the metrology signal from TAC algorithm.
int gravi_pfits_get_pola_num(const cpl_propertylist *plist, int type_data)
int gravi_pfits_get_fullstartx(const cpl_propertylist *plist)
cpl_error_code gravi_pfits_add_check(cpl_propertylist *header, const char *msg)
Add a QC.CHECK keyword to the header.
const char * gravi_pfits_get_spec_res(const cpl_propertylist *plist)
double gravi_pfits_get_dit_sc(const cpl_propertylist *plist)
double gravi_pfits_get_met_wavelength_mean(const cpl_propertylist *plist, cpl_table *met_table)
double gravi_pfits_get_double_default(const cpl_propertylist *plist, const char *name, double def)
cpl_error_code gravi_vis_create_opdsc_ft(cpl_table *vis_FT, cpl_table *vis_SC, double dit_sc)
Compute the resampled SC signal for each FT DIT per base.
int gravi_region_get_pol(cpl_table *imaging_detector, int region)
Return the polarisation id of a region.
int gravi_region_get_base(cpl_table *imaging_detector, int region)
Return the base of a region.
char GRAVI_BASE_NAME[6][3]
cpl_size gravi_spectrum_get_nwave(const cpl_table *table)
cpl_size gravi_spectrum_get_nregion(const cpl_table *table)
cpl_vector * gravi_compute_envelope(const cpl_vector *opd, int wave, int nwave)
Compute the envelope value.
int gravi_get_region(cpl_table *img_det, int base, char phase, int pol)
Find the region matching base, phase and pol.
int gravi_region_get_base_sign(cpl_table *imaging_detector, int base)
Return the sign of a base by looking at the PORT order.
cpl_table * gravi_opds_calibration(cpl_table *spectrum_table, cpl_table *detector_table, cpl_table *opdguess_table)
Compute the mean opd of each baseline from spectrum.
cpl_vector * gravi_opds_fit_opdmet(cpl_table *opd_ft, double lbd_met)
Compute the absolute scaling coefficients of SC and FT OPDs.
cpl_table * gravi_wave_fibre(cpl_table *spectrum_table, cpl_table *detector_table, cpl_table *opd_table)
Compute the wavelength of each channel for each 6 baselines.
cpl_table * gravi_wave_fit_2d(cpl_table *wavefibre_table, cpl_table *detector_table, gravi_data *wave_param, cpl_size fullstartx, int spatial_order, int spectral_order, double *rms_residuals)
Compute the WAVE_DATA table (1 wavelength per region) from the WAVE_FIBRE table (1 wavelength per bas...
#define PLOT_WAVE_PHASE_VS_OPD
cpl_error_code gravi_wave_compute_opds(gravi_data *spectrum_data, cpl_table *met_table, const cpl_parameterlist *parlist)
Recover the OPD modulation from a spectrum_data and vismet_table.
cpl_error_code gravi_opds_correct_closures(cpl_table *opd_sc, const char *name)
Correct the input table to ensure constant closure phase.
cpl_error_code gravi_compute_wave(gravi_data *wave_map, gravi_data *spectrum_data, int type_data, const cpl_parameterlist *parlist, gravi_data *wave_param)
Create the WAVE calibration map.
cpl_table * gravi_compute_argon_wave(cpl_table *spectrum_table)
Compute a WAVE calibration from the ARGON data (SC only)
cpl_table * gravi_wave_fit_individual(cpl_table *wave_individual_table, cpl_table *weight_individual_table, cpl_table *wave_fitted_table, cpl_table *opd_table, cpl_table *spectrum_table, cpl_table *detector_table, double n0, double n1, double n2)
cpl_error_code gravi_wave_qc(gravi_data *wave_map, gravi_data *profile_map)
Compute the QC parameters of the WAVE product.
cpl_error_code gravi_wave_correct_dispersion(cpl_table *wave_fibre, double n0, double n1, double n2)
Correct the WAVE_FIBRE table from a harcoded dispersion model.
cpl_error_code gravi_wave_correct_color(gravi_data *vis_data)
Create a OI_WAVELENGTH_CORR table with color corrected wavelength.
cpl_imagelist * gravi_wave_test_image(cpl_table *wavedata_table, cpl_table *wavefibre_table, cpl_table *profile_table, cpl_table *detector_table)
Compute the (useless) TEST_WAVE table from the WAVE_FIBRE and the PROFILE maps.
cpl_table * gravi_opds_compute_guess(cpl_table *spectrumsc_table, cpl_table *opdft_table, cpl_table *met_table, double dit_sc, double lbd_met)
Compute a guess of the OPD modulation of SC from FT and MET.