GRAVI Pipeline Reference Manual  1.2.3
gravi_wave.c
1 /* $Id: gravi_calib.c,v 1.10 2012/03/23 15:10:40 nazouaoui Exp $
2  *
3  * This file is part of the GRAVI Pipeline
4  * Copyright (C) 2002,2003 European Southern Observatory
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program; if not, write to the Free Software
18  * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19  */
20 
43 /*-----------------------------------------------------------------------------
44  Includes
45  -----------------------------------------------------------------------------*/
46 
47 #ifdef HAVE_CONFIG_H
48 #include <config.h>
49 #endif
50 
51 #include <cpl.h>
52 #include "cpl_wlcalib.h"
53 #include "irplib_wavecal.h"
54 #include <string.h>
55 #include <stdio.h>
56 #include <time.h>
57 #include <math.h>
58 #include <complex.h>
59 
60 #include "gravi_data.h"
61 #include "gravi_dfs.h"
62 #include "gravi_pfits.h"
63 #include "gravi_cpl.h"
64 
65 #include "gravi_utils.h"
66 #include "gravi_ellipse.h"
67 #include "gravi_signal.h"
68 
69 #include "gravi_preproc.h"
70 #include "gravi_calib.h"
71 #include "gravi_wave.h"
72 #include "gravi_metrology.h"
73 
74 /*----------------------------------------------------------------------------
75  DEBUG
76  -----------------------------------------------------------------------------*/
77 
78 #define PLOT_WAVE_PHASE_VS_OPD 0
79 #define WAVE_TO_PLOT 12
80 
81 /*-----------------------------------------------------------------------------
82  Private prototypes
83  -----------------------------------------------------------------------------*/
84 
85 cpl_table * gravi_opds_compute_guess (cpl_table * spectrumsc_table,
86  cpl_table * opdft_table,
87  cpl_table * met_table,
88  double dit_sc,
89  double lbd_met);
90 
91 cpl_table * gravi_opds_calibration (cpl_table * spectrum_table,
92  cpl_table * detector_table,
93  cpl_table * opdguess_table);
94 
95 cpl_error_code gravi_opds_correct_closures (cpl_table * opd_sc,
96  const char *name);
97 
98 cpl_vector * gravi_opds_fit_opdmet (cpl_table *opd_ft, double lbd_met);
99 
100 cpl_table * gravi_wave_fibre (cpl_table * spectrum_table,
101  cpl_table * detector_table,
102  cpl_table * opd_table);
103 
104 cpl_table * gravi_wave_fit_2d (cpl_table * wavefibre_table,
105  cpl_table * detector_table,
106  cpl_size fullstartx,
107  int spatial_order,
108  int spectral_order,
109  double * rms_residuals);
110 
111 cpl_table * gravi_wave_fit_individual (cpl_table * wave_individual_table,
112  cpl_table * weight_individual_table,
113  cpl_table * wave_fitted_table,
114  cpl_table * opd_table,
115  cpl_table * spectrum_table,
116  cpl_table * detector_table,
117  cpl_size fullstartx,
118  double n0, double n1, double n2,
119  double * rms_residuals);
120 
121 cpl_error_code gravi_wave_correct_dispersion (cpl_table * wave_fibre,
122  double n0, double n1, double n2);
123 
124 cpl_imagelist * gravi_wave_test_image (cpl_table * wavedata_table,
125  cpl_table * wavefibre_table,
126  cpl_table * profile_table,
127  cpl_table * detector_table);
128 
129 
130 
131 
132 /*-----------------------------------------------------------------------------
133  Functions code
134  -----------------------------------------------------------------------------*/
135 
136 /*----------------------------------------------------------------------------*/
147 /*----------------------------------------------------------------------------*/
148 
149 cpl_table * gravi_compute_argon_wave (cpl_table * spectrum_table)
150 {
151  gravi_msg_function_start(1);
152  cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
153 
154  /* Get the SC spectrums */
155  cpl_size nwave = cpl_table_get_column_depth (spectrum_table,"DATA1");
156  cpl_size nregion = gravi_spectrum_get_nregion (spectrum_table);
157  CPLCHECK_NUL ("Cannot get data");
158 
159  cpl_msg_info (cpl_func, "spectrum_table = %lld",cpl_table_get_column_depth (spectrum_table, "DATA1"));
160 
161  /* Prepare the outputs */
162  cpl_table * wave_data_sc = cpl_table_new (1);
163  CPLCHECK_NUL ("Cannot create output");
164 
165  /* Init the list of lines*/
166  cpl_size nlines = 14;
167  double plines_str[] = {60,
168  160,
169  90,
170  1800,
171  50,
172  30,
173  30,
174  1500,
175  130,
176  510,
177  240,
178  169,
179  64,
180  73
181  };
182 
183  double plines_lbd[] = {1.982291,
184  1.997118,
185  2.032256,
186  2.062186,
187  2.065277,
188  2.073922,
189  2.081672,
190  2.099184,
191  2.133871,
192  2.154009,
193  2.208321,
194  2.313952,
195  2.385154,
196  2.397306
197  };
198 
199  cpl_vector * lines_lbd = cpl_vector_wrap (nlines,plines_lbd);
200  cpl_vector * lines_str = cpl_vector_wrap (nlines,plines_str);
201  cpl_bivector * lines = cpl_bivector_wrap_vectors (lines_lbd, lines_str);
202 
203  /* Init the slit and lamp model */
204  cpl_wlcalib_slitmodel * model = cpl_wlcalib_slitmodel_new ();
205  cpl_wlcalib_slitmodel_set_wslit (model, 0.1);
206  cpl_wlcalib_slitmodel_set_wfwhm (model, 1.25);
207  cpl_wlcalib_slitmodel_set_threshold (model, 5.0);
208  cpl_wlcalib_slitmodel_set_catalog (model, lines);
209 
210  /* Starting point of dispersion */
211  double values_hr[] = {1.97048,0.2228e-3,2.7728e-08,-4.72e-12};
212  double values_mr[] = {1.9697,0.2179e-2,2.5e-07,-5e-10};
213  double * values = nwave > 1000 ? values_hr : values_mr;
214 
215  /* Fill the dispersion polynomial */
216  cpl_polynomial * dispersion0 = cpl_polynomial_new (1);
217  for (cpl_size order = 0 ; order < 4 ; order ++) {
218  cpl_polynomial_set_coeff (dispersion0, &order, values[order]);
219  }
220 
221  /* QC parameters */
222  double minwave = -1e10;
223  double maxwave = +1e10;
224 
225  cpl_vector * spectrum = cpl_vector_new (nwave);
226  cpl_vector * input_spectrum = cpl_vector_new (nwave);
227  cpl_vector * model_spectrum = cpl_vector_new (nwave);
228 
229  /* Prepare output table */
230  cpl_table * fit_table = cpl_table_new (nregion);
231  gravi_table_new_column_array (fit_table, "DATA", "adu", CPL_TYPE_DOUBLE, nwave);
232  gravi_table_new_column_array (fit_table, "DATA_MODEL", "adu", CPL_TYPE_DOUBLE, nwave);
233  gravi_table_new_column_array (fit_table, "DATA_WAVE", "adu", CPL_TYPE_DOUBLE, nwave);
234 
235  /* The starting point will be the fit of previous region */
236  cpl_polynomial * dispersion = cpl_polynomial_duplicate (dispersion0);
237 
238  /* Loop on regions */
239  for (cpl_size reg = 0; reg < nregion ; reg ++ ) {
240  cpl_msg_info (cpl_func,"Fit region %lld over %lld", reg+1, nregion);
241 
242  /* Copy the spectrum of this region into a vector */
243  const cpl_array * data = cpl_table_get_array (spectrum_table, GRAVI_DATA[reg], 0);
244  for (cpl_size wave=0; wave < nwave; wave ++) {
245  cpl_vector_set (spectrum, wave, log (CPL_MAX (cpl_array_get (data, wave, NULL), 1.0)));
246  cpl_vector_set (input_spectrum, wave, cpl_array_get (data, wave, NULL));
247  }
248 
249  CPLCHECK_NUL ("Cannot prepare data");
250 
251  /* Parameters: Explore a large hsize [pixel]
252  * only for the first region, since we keep
253  * the starting point */
254  int maxdeg = 3, nmaxima = 0, linelim = 99;
255  int maxite = 100 * maxdeg, maxfail = 10, maxcont = 10;
256  int hsize = (reg == 0 ? 40 : 5);
257  double pixtol = 1e-6, pixstep = 0.5, pxc = 0.0;
258 
259  /* Use the log of spectrum to put
260  * less weight on the intensity */
261  irplib_polynomial_find_1d_from_correlation_all (dispersion, maxdeg, spectrum,
262  nmaxima, linelim,
263  (irplib_base_spectrum_model*)model,
264  &irplib_vector_fill_logline_spectrum,
265  pixtol, pixstep,
266  hsize, maxite, maxfail, maxcont,
267  CPL_FALSE,&pxc);
268  CPLCHECK_NUL ("Cannot fit polynomial");
269 
270  /* Fill in the final fit */
271  irplib_vector_fill_line_spectrum (model_spectrum, dispersion, (void *)model);
272  cpl_polynomial_dump (dispersion, stdout);
273 
274  /* Evaluate the polynomial for the output */
275  cpl_array * wave_data = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
276  for (cpl_size pix = 0; pix < nwave ; pix++)
277  cpl_array_set (wave_data, pix, cpl_polynomial_eval_1d (dispersion, (double)pix, NULL) * 1e-6);
278 
279  /* Check the limits for the QC */
280  minwave = CPL_MAX (cpl_array_get_min (wave_data), minwave);
281  maxwave = CPL_MIN (cpl_array_get_max (wave_data), maxwave);
282 
283  /* Fill the WAVE_DATA_SC */
284  cpl_table_new_column_array (wave_data_sc, GRAVI_DATA[reg], CPL_TYPE_DOUBLE, nwave);
285  cpl_table_set_array (wave_data_sc, GRAVI_DATA[reg], 0, wave_data);
286  FREE (cpl_array_delete, wave_data);
287 
288  /* Save the test tables */
289  cpl_array * tmp_array;
290  tmp_array = cpl_array_wrap_double (cpl_vector_get_data (input_spectrum), nwave);
291  cpl_table_set_array (fit_table, "DATA", reg, tmp_array);
292  cpl_array_unwrap (tmp_array);
293  tmp_array = cpl_array_wrap_double (cpl_vector_get_data (model_spectrum), nwave);
294  cpl_table_set_array (fit_table, "DATA_MODEL", reg, tmp_array);
295  cpl_array_unwrap (tmp_array);
296  tmp_array = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
297  for (cpl_size wave = 0; wave < nwave; wave ++)
298  cpl_array_set (tmp_array, wave, cpl_polynomial_eval_1d (dispersion, (double)wave, NULL) * 1e-6);
299  cpl_table_set_array (fit_table, "DATA_WAVE", reg, tmp_array);
300  cpl_array_delete (tmp_array);
301 
302  CPLCHECK_NUL ("Cannot set data");
303  } /* End loop on regions */
304 
305  // cpl_table_save (fit_table, NULL, NULL, "fit_wave.fits", CPL_IO_CREATE);
306  FREE (cpl_table_delete, fit_table);
307 
308  /* Delete */
309  FREE (cpl_vector_delete, spectrum);
310  FREE (cpl_vector_delete, input_spectrum);
311  FREE (cpl_vector_delete, model_spectrum);
312  FREE (cpl_polynomial_delete, dispersion0);
313  FREE (cpl_polynomial_delete, dispersion);
314 
315  /* This memory desalocation does not work: */
316  // FREE (cpl_wlcalib_slitmodel_delete, model);
317 
318  gravi_msg_function_exit(1);
319  return wave_data_sc;
320 }
321 
322 /*----------------------------------------------------------------------------*/
352 /*----------------------------------------------------------------------------*/
353 
354 cpl_error_code gravi_opds_correct_closures (cpl_table * phase_table,
355  const char *name)
356 {
357  gravi_msg_function_start(1);
358  cpl_ensure_code (phase_table, CPL_ERROR_NULL_INPUT);
359  cpl_ensure_code (name, CPL_ERROR_NULL_INPUT);
360 
361  int nbase = 6, nclo = 4;
362  cpl_size nrow = cpl_table_get_nrow (phase_table) / nbase;
363 
364  /* Get the data */
365  double * phase = cpl_table_get_data_double (phase_table, name);
366  CPLCHECK_MSG ("Cannot get data");
367 
368  /* Model and right-hand-side for the lineary system
369  * (unfilled matrix are 0.0) */
370  cpl_matrix * rhs_matrix = cpl_matrix_new (nrow * nclo, 1);
371  cpl_matrix * model_matrix = cpl_matrix_new (nrow * nclo, nbase-1 + nclo);
372 
373  for (cpl_size row = 0; row < nrow; row++) {
374  for (int clo = 0; clo < nclo; clo++) {
375  int b0 = GRAVI_CLO_BASE[clo][0];
376  int b1 = GRAVI_CLO_BASE[clo][1];
377  int b2 = GRAVI_CLO_BASE[clo][2];
378 
379  /* Fill the rhs with measured closure */
380  double closure = phase[row*nbase+b0] +
381  phase[row*nbase+b1] -
382  phase[row*nbase+b2];
383  cpl_matrix_set (rhs_matrix, row*nclo+clo, 0, closure);
384  CPLCHECK_MSG ("Cannot fill rhs_matrix");
385 
386  /* Fill the k (unfilled are 0.0) */
387  cpl_matrix_set (model_matrix, row*nclo+clo, clo, 1.0);
388  CPLCHECK_MSG ("Cannot fill k in model_matrix");
389 
390  /* Fill the f in of the model_matrix */
391  if (b0 < nbase-1) cpl_matrix_set (model_matrix, row*nclo+clo, nclo+b0,
392  +phase[row*nbase+b0]);
393  if (b1 < nbase-1) cpl_matrix_set (model_matrix, row*nclo+clo, nclo+b1,
394  +phase[row*nbase+b1]);
395  if (b2 < nbase-1) cpl_matrix_set (model_matrix, row*nclo+clo, nclo+b2,
396  -phase[row*nbase+b2]);
397  CPLCHECK_MSG ("Cannot fill phase in model_matrix");
398  }
399  } /* End loop on clo and rows */
400 
401  /* Solve the system */
402  cpl_matrix * res_matrix = cpl_matrix_solve_normal (model_matrix, rhs_matrix);
403  CPLCHECK_MSG ("Cannot solve system");
404 
405  /* Compute residuals */
406  cpl_matrix * residual_matrix = cpl_matrix_product_create (model_matrix, res_matrix);
407  cpl_matrix_subtract (residual_matrix, rhs_matrix);
408  double rms_fit = cpl_matrix_get_stdev (residual_matrix);
409 
410  /* Correct baseline with (1-f). Last baseline is kept unchanged */
411  for (int base = 0; base < nbase - 1; base ++) {
412  double f = cpl_matrix_get (res_matrix, nclo + base, 0);
413  cpl_msg_info (cpl_func,"correction factor f = 1 %+.20f", -1*f);
414  for (cpl_size row = 0; row < nrow; row++) phase[row*nbase+base] *= 1 - f;
415  }
416 
417  /* Dump the residual with their units */
418  const char * unit = cpl_table_get_column_unit (phase_table, name);
419  cpl_msg_info (cpl_func, "residual stdev = %.20g [%s]", rms_fit, unit);
420 
421  /* Delete matrix */
422  FREE (cpl_matrix_delete, residual_matrix);
423  FREE (cpl_matrix_delete, res_matrix);
424  FREE (cpl_matrix_delete, model_matrix);
425  FREE (cpl_matrix_delete, rhs_matrix);
426 
427  gravi_msg_function_exit(1);
428  return CPL_ERROR_NONE;
429 }
430 
431 /*----------------------------------------------------------------------------*/
457 /*----------------------------------------------------------------------------*/
458 
459 cpl_vector * gravi_opds_fit_opdmet (cpl_table * ft_table, double lbd_met)
460 {
461  gravi_msg_function_start(1);
462  cpl_ensure (ft_table, CPL_ERROR_NULL_INPUT, NULL);
463 
464  int nbase = 6;
465 
466  /* Get the number of acquisitions */
467  cpl_size nrow = cpl_table_get_nrow (ft_table) / nbase;
468 
469  /* Get the pointer to data */
470  double * opd_sc = cpl_table_get_data_double (ft_table, "OPD_SC");
471  double * opd_ft = cpl_table_get_data_double (ft_table, "OPD");
472  double * phase_met = cpl_table_get_data_double (ft_table, "PHASE_MET_FC");
473  CPLCHECK_NUL ("Cannot get data");
474 
475  /* Number of valid rows */
476  cpl_size nrow_valid = 0;
477  for (cpl_size row = 0; row < nrow; row++) if (opd_sc[row*nbase] != 0) nrow_valid++;
478  cpl_msg_info (cpl_func,"nrow_valid = %lld", nrow_valid);
479 
480  cpl_ensure (nrow_valid > 100, CPL_ERROR_ILLEGAL_INPUT, NULL);
481 
482  /* Model and right-hand-side for the lineary system
483  * (unfilled matrix are 0.0) */
484  cpl_matrix * rhs_matrix = cpl_matrix_new (nrow_valid*nbase, 1);
485  cpl_matrix * model_matrix = cpl_matrix_new (nrow_valid*nbase, nbase + 2);
486 
487  for (int base = 0; base < nbase; base++) {
488  cpl_size row_valid = 0;
489 
490  for (cpl_size row=0; row<nrow; row++) {
491  if (opd_sc[row*nbase] == 0) continue;
492 
493  int idv = row_valid * nbase + base;
494  int id = row * nbase + base;
495 
496  /* Fill the OPD metrology */
497  cpl_matrix_set (rhs_matrix, idv, 0, phase_met[id] * lbd_met / CPL_MATH_2PI);
498  CPLCHECK_NUL ("Cannot set OPD");
499 
500  /* Fill the model b and c */
501  cpl_matrix_set (model_matrix, idv, 0, 1*opd_sc[id]);
502  cpl_matrix_set (model_matrix, idv, 1, -1*opd_ft[id]);
503  CPLCHECK_NUL ("Cannot set SC or FT");
504 
505  /* Fill the model Aij (unfilled matrix are 0.0) */
506  cpl_matrix_set (model_matrix, idv, 2 + base, 1.0);
507  CPLCHECK_NUL ("Cannot set the zero points");
508 
509  row_valid++;
510  } /* End loop rows */
511  } /* End loop on base*/
512 
513  /* Solve the system */
514  cpl_matrix * res_matrix = cpl_matrix_solve_normal (model_matrix, rhs_matrix);
515 
516  /* Compute residuals */
517  cpl_matrix * residual_matrix = cpl_matrix_product_create (model_matrix, res_matrix);
518  cpl_matrix_subtract (residual_matrix, rhs_matrix);
519  double rms_fit = cpl_matrix_get_stdev (residual_matrix);
520 
521  /* Verbose */
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);
525 
526  /* Set the Results */
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);
531 
532  /* Delete matrix */
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);
537 
538  gravi_msg_function_exit(1);
539  return opd_coeff;
540 }
541 
542 /*----------------------------------------------------------------------------*/
564 /*----------------------------------------------------------------------------*/
565 
566 cpl_table * gravi_opds_compute_guess (cpl_table * spectrumsc_table,
567  cpl_table * ft_table,
568  cpl_table * vismet_table,
569  double dit_sc,
570  double lbd_met)
571 {
572  gravi_msg_function_start(1);
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);
577 
578  int nbase = 6, ntel = 4;
579 
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");
584  CPLCHECK_NUL ("Cannot get data");
585 
586  /* Pointer on MET data */
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");
589  CPLCHECK_NUL ("Cannot get MET data");
590 
591  /* Pointer on OPD data */
592  int * time_FT = cpl_table_get_data_int (ft_table, "TIME");
593  double * opd_FT = cpl_table_get_data_double (ft_table, "OPD");
594  CPLCHECK_NUL ("Cannot get FT data");
595 
596  /* Create table */
597  cpl_table * guess_table = cpl_table_new (nrow * nbase);
598  gravi_table_new_column (guess_table, "OPD", "m", CPL_TYPE_DOUBLE);
599 
600  /* Loop on base */
601  for (int base = 0; base < nbase; base++) {
602  int tel0 = GRAVI_BASE_TEL[base][0];
603  int tel1 = GRAVI_BASE_TEL[base][1];
604 
605  /* Loop on the SC frames */
606  cpl_size row_met = 0, row_ft = 0;
607  for (cpl_size row = 0 ; row < nrow ; row++) {
608 
609  /*
610  * Average the MET
611  */
612  int counter_met = 0;
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];
617  counter_met ++;
618  }
619 
620  /* If not enough data to synchronize */
621  if (row_met > nrow_met - 2) {
622  cpl_msg_warning (cpl_func,"Not enough data to synchronise the MET with SC");
623  break;
624  }
625  row_met ++;
626  }
627  /* End sum the MET over the current SC DIT */
628 
629  opd_met = opd_met * lbd_met / CPL_MATH_2PI / counter_met; // [m]
630 
631  /*
632  * Average the FT
633  */
634  int counter_ft = 0;
635  double opd_ft = 0.0;
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];
639  counter_ft ++;
640  }
641 
642  /* If not enough data to synchronize */
643  if (row_ft > nrow_ft - 2) {
644  cpl_msg_warning (cpl_func,"Not enough data to synchronise the FT with SC");
645  break;
646  }
647  row_ft ++;
648  }
649  /* End sum the FT over the current SC DIT */
650 
651  opd_ft = opd_ft / counter_ft;
652 
653  /* Set the total guess OPD as OPD_FT - OPD_MET */
654  cpl_table_set (guess_table, "OPD", row*nbase+base, opd_ft - opd_met);
655  CPLCHECK_NUL ("Cannot compute opd guess");
656  } /* End loop on row */
657  } /* End loop on base */
658 
659  /*
660  * Remove the mean of each base
661  */
662  cpl_msg_info (cpl_func, "Remove the mean from opdguess table");
663 
664  for (int base = 0; base < 6; base++) {
665  double mean = gravi_table_get_column_mean (guess_table, "OPD", base, nbase);
666  gravi_table_add_scalar (guess_table, "OPD", base, nbase, -1*mean);
667  }
668 
669  gravi_msg_function_exit(1);
670  return guess_table;
671 }
672 
673 /*----------------------------------------------------------------------------*/
696 /*----------------------------------------------------------------------------*/
697 
698 cpl_table * gravi_opds_calibration (cpl_table * spectrum_table,
699  cpl_table * detector_table,
700  cpl_table * guess_table)
701 {
702  int nbase = 6;
703 
704  /* Verbose */
705  gravi_msg_function_start(1);
706  cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
707  cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
708 
709  /* Verbose the phase algorithm selected */
710  if (USE_LINEAR_ELLIPSE) {
711  cpl_msg_info (cpl_func, "Solve linearly the ellipse X,Y");
712  }
713  else {
714  cpl_msg_info (cpl_func, "Solve non-linearly the ellipse X,Y");
715  }
716 
717  /* Get the size of the vectors */
718  cpl_size nrow = cpl_table_get_nrow (spectrum_table);
719  cpl_size nwave = gravi_spectrum_get_nwave (spectrum_table);
720  cpl_size npol = gravi_spectrum_get_npol (spectrum_table);
721  CPLCHECK_NUL ("Cannot get size");
722 
723  /* Create a fake OI_WAVE table, to average the different
724  * spectral channels [m] */
725 
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);
733  }
734  }
735 
736  /* Create column in phase table */
737  cpl_table * output_table = cpl_table_new (nbase * nrow);
738  gravi_table_new_column (output_table,"OPD", "m", CPL_TYPE_DOUBLE);
739  gravi_table_new_column (output_table,"TIME", "us", CPL_TYPE_INT);
740 
741  /* Set the time */
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);
746  }
747 
748  /* Loop on base */
749  for (int base = 0; base < 6; base ++) {
750 
751  /* Check if a guess exists */
752  cpl_vector * opd_guess = NULL;
753  if (guess_table) {
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);
758  }
759  }
760 
761  /* Recover the mean OPD modulation (averaved over pol
762  * and channels) using the ellipses. In [m]. However since oiwave
763  * is not known... this is to a scaling factor correction. */
764  cpl_vector * mean_opd;
765  mean_opd = gravi_ellipse_meanopd_create (spectrum_table, detector_table,
766  oiwave_tables, opd_guess, base);
767  CPLCHECK_NUL ("Cannot compute opd");
768 
769  /* Save the mean opd [m] */
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);
773  }
774  CPLCHECK_NUL ("Cannot set phase");
775 
776  FREE (cpl_vector_delete, mean_opd);
777  FREE (cpl_vector_delete, opd_guess);
778 
779  } /* End loop on base */
780 
781  FREELOOP (cpl_table_delete, oiwave_tables, npol);
782 
783  /* Verbose */
784  gravi_msg_function_exit(1);
785  return output_table;
786 }
787 
788 /*----------------------------------------------------------------------------*/
820 /*----------------------------------------------------------------------------*/
821 
822 cpl_error_code gravi_wave_compute_opds (gravi_data * spectrum_data,
823  cpl_table * met_table,
824  enum gravi_detector_type det_type)
825 {
826  gravi_msg_function_start(1);
827  cpl_ensure_code (spectrum_data, CPL_ERROR_NULL_INPUT);
828  cpl_ensure_code (met_table, CPL_ERROR_NULL_INPUT);
829 
830  /* Get DIT */
831  cpl_propertylist * spectrum_header = gravi_data_get_header (spectrum_data);
832  double dit_sc = gravi_pfits_get_dit_sc (spectrum_header)*1e6; // [us]
833  CPLCHECK_MSG ("Cannot get DIT");
834 
835  /* Get the input */
836  cpl_table * spectrumft_table = gravi_data_get_spectrum_data (spectrum_data, GRAVI_FT);
837  cpl_table * detectorft_table = gravi_data_get_imaging_detector (spectrum_data, GRAVI_FT);
838  cpl_table * spectrumsc_table = gravi_data_get_spectrum_data (spectrum_data, GRAVI_SC);
839  cpl_table * detectorsc_table = gravi_data_get_imaging_detector (spectrum_data, GRAVI_SC);
840  CPLCHECK_MSG ("Cannot get data");
841 
842 
843  /*
844  * Reduce the raw METROLOGY into VIS_MET
845  */
846  cpl_msg_info (cpl_func, "Compute the phase of MET_FC");
847 
848  cpl_table * vismet_table = gravi_metrology_create (met_table, spectrum_header);
849  gravi_metrology_drs (met_table, vismet_table, spectrum_header);
850 
851  /* Compute the mean LBD_MET for this file */
852  double lbd_met = gravi_pfits_get_met_wavelength_mean (spectrum_header, met_table);
853 
854  /*
855  * Compute the phase of SC and FT. For the SC, we use a guess
856  * of the opd modulation, computed from FT and MET, to unwrap.
857  */
858  cpl_msg_info (cpl_func, "Compute OPD of FT from ellipse");
859 
860  cpl_table * ft_table;
861  ft_table = gravi_opds_calibration (spectrumft_table, detectorft_table, NULL);
862  gravi_opds_correct_closures (ft_table, "OPD");
863 
864  cpl_msg_info (cpl_func, "Compute OPD of SC from ellipse");
865 
866  cpl_table * guess_table;
867  guess_table = gravi_opds_compute_guess (spectrumsc_table, ft_table, vismet_table, dit_sc, lbd_met);
868 
869  cpl_table * sc_table;
870  sc_table = gravi_opds_calibration (spectrumsc_table, detectorsc_table, guess_table);
871  gravi_opds_correct_closures (sc_table, "OPD");
872  FREE (cpl_table_delete, guess_table);
873 
874  CPLCHECK_MSG ("Cannot calibrate phase");
875 
876  /*
877  * Interpolate MET and SC at the sampling frequency of FT
878  * Results are saved in the ft_table table
879  */
880  cpl_msg_info (cpl_func, "Fit MET = a.SC - b.FT + c to get absolute modulation");
881 
882  gravi_vis_create_met_ft (ft_table, vismet_table);
883  gravi_vis_create_opdsc_ft (ft_table, sc_table, dit_sc);
884 
885  CPLCHECK_MSG ("Cannot resample SC or MET at the FT frequency");
886 
887  /*
888  * Compute the scaling coefficients of OPDs by fitting:
889  * OPD_MET_ijt = a.OPD_SC_ijt - b.OPD_FT_ijt + c_ij
890  */
891  cpl_vector * coeff_vector = gravi_opds_fit_opdmet (ft_table, lbd_met);
892 
893  CPLCHECK_MSG ("Cannot fit opdmet");
894 
895  /* Set the CHI2 of the fit in the QC parameters */
896  cpl_propertylist_update_double (spectrum_header, "ESO QC PHASE_CALIBRATION_CHI2",
897  cpl_vector_get (coeff_vector, 2));
898  cpl_propertylist_set_comment (spectrum_header, "ESO QC PHASE_CALIBRATION_CHI2",
899  "chi2 of a.SC-b.FT+c=MET");
900 
901  /* Add the OPD COEFF in QC parameters */
902  for (int type_data = 0; type_data < 2; type_data++) {
903  double tmp = cpl_vector_get (coeff_vector, type_data);
904  cpl_propertylist_update_float (spectrum_header, OPD_COEFF_SIGN(type_data), tmp);
905  cpl_propertylist_set_comment (spectrum_header, OPD_COEFF_SIGN(type_data), "wavelength correction");
906  }
907 
908  /* Set a warning */
909  if (cpl_vector_get (coeff_vector, 2) > 1e-7) {
910  gravi_pfits_add_check (spectrum_header,"Residu of fit MET=a.SC-b.FT+c are high");
911 
912  cpl_msg_info (cpl_func, "*************************************************");
913  cpl_msg_warning (cpl_func, "**** !!! residuals of the fit too high !!! ****");
914  cpl_msg_warning (cpl_func, "**** SC and RMN may be desynchronized ****");
915  cpl_msg_warning (cpl_func, "**** (or out of the envelope in LOW) ****");
916  cpl_msg_info (cpl_func, "*************************************************");
917  }
918 
919  CPLCHECK_MSG ("Cannot set QC parameter");
920 
921  /*
922  * Correct opd from overall scaling by fitting MET
923  */
924  double coeff_sc = cpl_vector_get (coeff_vector, GRAVI_SC);
925  cpl_table_multiply_scalar (sc_table, "OPD", coeff_sc);
926 
927  double coeff_ft = cpl_vector_get (coeff_vector, GRAVI_FT);
928  cpl_table_multiply_scalar (ft_table, "OPD", coeff_ft);
929 
930  FREE (cpl_vector_delete, coeff_vector);
931  CPLCHECK_MSG ("Cannot correct OPDs from scaling coefficients");
932 
933  /*
934  * Fill the output
935  */
936  if ((det_type == GRAVI_DET_SC || det_type == GRAVI_DET_ALL))
937  {
938  gravi_data_add_table (spectrum_data, NULL, "OPD_SC", sc_table);
939  CPLCHECK_MSG ("Cannot set OPD_SC table");
940  }
941  else
942  FREE (cpl_table_delete, sc_table);
943 
944 
945  if ((det_type == GRAVI_DET_FT || det_type == GRAVI_DET_ALL))
946  {
947  gravi_data_add_table (spectrum_data, NULL, "OPD_FT", ft_table);
948  CPLCHECK_MSG ("Cannot set OPD_FT table");
949  }
950  else
951  FREE (cpl_table_delete, ft_table);
952 
953  gravi_data_add_table (spectrum_data, NULL, GRAVI_OI_VIS_MET_EXT, vismet_table);
954  CPLCHECK_MSG ("Cannot set OI_VIS_MET table");
955 
956  gravi_msg_function_exit(1);
957  return CPL_ERROR_NONE;
958 }
959 
960 /*----------------------------------------------------------------------------*/
984 /*----------------------------------------------------------------------------*/
985 
986 cpl_table * gravi_wave_fibre (cpl_table * spectrum_table,
987  cpl_table * detector_table,
988  cpl_table * opd_table)
989 {
990  gravi_msg_function_start(1);
991  cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
992  cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
993  cpl_ensure (opd_table, CPL_ERROR_NULL_INPUT, NULL);
994 
995  int nbase = 6;
996  char name[100];
997 
998  /* Create the output table */
999  cpl_table * wave_fibre = cpl_table_new (1);
1000 
1001  /* Get the number of wavelength, region, polarisation... */
1002  cpl_size nwave = cpl_table_get_column_depth (spectrum_table, "DATA1");
1003  cpl_size n_region = cpl_table_get_nrow (detector_table);
1004  cpl_size nrow = cpl_table_get_nrow (spectrum_table);
1005 
1006  int npol = (n_region > 24 ? 2 : 1);
1007 
1008  /*
1009  * Calibration of each polarization and base
1010  */
1011 
1012  for (int pol = 0; pol < npol; pol++) {
1013  for (int base = 0; base < nbase; base ++) {
1014  cpl_msg_info (cpl_func, "Compute wave fibre for pol %i over %i, base %i over %i",
1015  pol+1, npol, base+1, nbase);
1016 
1017  /* Get the index of the ABCD. */
1018  int iA = gravi_get_region (detector_table, base, 'A', pol);
1019  int iB = gravi_get_region (detector_table, base, 'B', pol);
1020  int iC = gravi_get_region (detector_table, base, 'C', pol);
1021  int iD = gravi_get_region (detector_table, base, 'D', pol);
1022  if (iA<0 || iB<0 || iC<0 || iD<0){
1023  cpl_msg_warning (cpl_func, "Don't found the A, B, C or D !!!");
1024  }
1025 
1026  /* Sign of this baseline */
1027  double phi_sign = gravi_region_get_base_sign (detector_table, base);
1028 
1029  /* Get the OPD into various flarous */
1030  cpl_matrix * opd_matrix = cpl_matrix_new (1, nrow);
1031  cpl_vector * opd_vector = cpl_vector_new (nrow);
1032  for (cpl_size row = 0; row < nrow; row ++ ) {
1033  double value = cpl_table_get (opd_table, "OPD", row*nbase+base, NULL);
1034  cpl_matrix_set (opd_matrix, 0, row, value);
1035  cpl_vector_set (opd_vector, row, value);
1036  }
1037  CPLCHECK_NUL ("Cannot extract the OPD");
1038 
1039  /* Create array to fill the wavelenghts */
1040  cpl_array * wavelength = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1041  cpl_array * wavechi2 = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1042  cpl_array_fill_window (wavelength, 0, nwave, 0.0);
1043  cpl_array_fill_window (wavechi2, 0, nwave, 1e10);
1044 
1045  /*
1046  * Loop on spectral channel
1047  */
1048  for (cpl_size wave = 0; wave < nwave; wave++) {
1049 
1050  cpl_vector * vector_T = NULL, * vector_X, * vector_Y;
1051 
1052  /* Define the vector_X = C - A */
1053  vector_X = gravi_table_get_vector (spectrum_table, wave, GRAVI_DATA[iC]);
1054  vector_T = gravi_table_get_vector (spectrum_table, wave, GRAVI_DATA[iA]);
1055  cpl_vector_subtract (vector_X, vector_T);
1056  FREE (cpl_vector_delete, vector_T);
1057 
1058  /* Define the vector_Y = D - B */
1059  vector_Y = gravi_table_get_vector (spectrum_table, wave, GRAVI_DATA[iD]);
1060  vector_T = gravi_table_get_vector (spectrum_table, wave, GRAVI_DATA[iB]);
1061  cpl_vector_subtract (vector_Y, vector_T);
1062  FREE (cpl_vector_delete, vector_T);
1063 
1064  CPLCHECK_NUL ("Cannot compute vector_X or Y");
1065 
1066  /* Compute envelope from OPD for this channel */
1067  cpl_vector * envelope_vector = gravi_compute_envelope (opd_vector, wave, nwave);
1068 
1069  /* Compute the phase from the ellipse */
1070  cpl_vector * phase;
1071  phase = gravi_ellipse_phase_create (vector_X, vector_Y,
1072  envelope_vector);
1073  FREE (cpl_vector_delete, vector_X);
1074  FREE (cpl_vector_delete, vector_Y);
1075  FREE (cpl_vector_delete, envelope_vector);
1076 
1077  /* If the computation of the ellipse fails, we continue with next wave */
1078  if (phase == NULL) {
1079  cpl_msg_warning (cpl_func, "Cannot compute wave for %lld and base %d", wave, base);
1080  continue;
1081  }
1082 
1083  /* Multiply by the sign */
1084  cpl_vector_multiply_scalar (phase, phi_sign);
1085 
1086  /* Unwrap phase from the OPD */
1087  double lbd_channel = 1.95e-6 + (2.46e-6 - 1.95e-6) / nwave * wave;
1088  gravi_vector_unwrap_with_guess (phase, opd_vector, CPL_MATH_2PI / lbd_channel);
1089 
1090  /* Fit the slope of the phase versus OPD gives the
1091  * wavelength of the spectral element */
1092  const cpl_size mindeg = 0, maxdeg = 1;
1093  cpl_polynomial * fit_slope = cpl_polynomial_new (1);
1094  cpl_polynomial_fit (fit_slope, opd_matrix, NULL, phase, NULL, CPL_FALSE, &mindeg, &maxdeg);
1095 
1096  /* Compute residuals */
1097  cpl_vector * residuals = cpl_vector_new (nwave);
1098  double rechisq;
1099  cpl_vector_fill_polynomial_fit_residual (residuals, phase, NULL, fit_slope, opd_matrix, &rechisq);
1100 
1101  if(PLOT_WAVE_PHASE_VS_OPD && (wave == WAVE_TO_PLOT))
1102  {
1103  char gnuplot_str[200];
1104  sprintf (gnuplot_str, "set title 'Wavelength (base %d)'; set xlabel 'Phase [rad]'; set ylabel 'OPD (m)';", base);
1105  cpl_plot_vector (gnuplot_str, NULL, NULL, opd_vector);
1106  sprintf (gnuplot_str, "set title 'Wavelength residuals (base %d)'; set xlabel 'Phase [rad]'; set ylabel 'OPD (m)';", base);
1107  cpl_plot_vector (gnuplot_str, NULL, NULL, residuals);
1108  CPLCHECK_NUL ("Cannot plot OPD versus phase");
1109  }
1110 
1111  CPLCHECK_NUL ("Cannot fit the OPD and phase");
1112 
1113  /* Get the slope */
1114  const cpl_size pow_slope = 1;
1115  double slope = cpl_polynomial_get_coeff (fit_slope, &pow_slope);
1116 
1117  /* Check slope sign. Should be positive */
1118  if (slope < 0.0 && wave == 0) {
1119  cpl_msg_warning (cpl_func, "Negative wavelength!! "
1120  "Report to DRS team");
1121  }
1122 
1123  /* Set the wavelength */
1124  cpl_array_set (wavechi2, wave, sqrt(rechisq));
1125  cpl_array_set (wavelength, wave, CPL_MATH_2PI / fabs(slope));
1126 
1127  /* Free memory */
1128  cpl_vector_delete (phase);
1129  cpl_vector_delete (residuals);
1130  cpl_polynomial_delete (fit_slope);
1131 
1132  } /* End loop on wave */
1133 
1134  cpl_matrix_delete (opd_matrix);
1135  cpl_vector_delete (opd_vector);
1136 
1137  /* Set the wavelength array in the output table */
1138  sprintf (name, "BASE_%s_%s", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol,npol));
1139  cpl_table_new_column_array (wave_fibre, name, CPL_TYPE_DOUBLE, nwave);
1140  cpl_table_set_column_unit (wave_fibre, name, "m");
1141  cpl_table_set_array (wave_fibre, name, 0, wavelength);
1142  cpl_array_delete (wavelength);
1143 
1144  /* Set the chi2 array in the output table */
1145  sprintf (name, "BASE_%s_%s_CHI2", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol,npol));
1146  cpl_table_new_column_array (wave_fibre, name, CPL_TYPE_DOUBLE, nwave);
1147  cpl_table_set_array (wave_fibre, name, 0, wavechi2);
1148  cpl_array_delete (wavechi2);
1149 
1150  } /* End loop on base*/
1151  } /* End loop on polar */
1152 
1153  gravi_msg_function_exit(1);
1154  return wave_fibre;
1155 }
1156 
1157 /*----------------------------------------------------------------------------*/
1175 /*----------------------------------------------------------------------------*/
1176 
1177 cpl_table * gravi_wave_fit_2d (cpl_table * wavefibre_table,
1178  cpl_table * detector_table,
1179  cpl_size fullstartx,
1180  int spatial_order,
1181  int spectral_order,
1182  double * rms_residuals)
1183 {
1184  gravi_msg_function_start(1);
1185  cpl_ensure (wavefibre_table, CPL_ERROR_NULL_INPUT, NULL);
1186  cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
1187 
1188  int nbase = 6;
1189  char name[100];
1190  *rms_residuals = 0;
1191 
1192  /* Get numbers */
1193  cpl_size n_region = cpl_table_get_nrow (detector_table);
1194  int npol = n_region > 24 ? 2 : 1;
1195  cpl_size nwave = cpl_table_get_column_depth (wavefibre_table, npol > 1 ? "BASE_12_S" : "BASE_12_C");
1196  CPLCHECK_NUL ("Cannot get data");
1197 
1198  /* Odd index, for SC only */
1199  cpl_vector * odd_index = cpl_vector_new (nwave);
1200  for (int i = fullstartx; i < fullstartx + nwave; i++) {
1201  if (nwave > GRAVI_LBD_FTSC) cpl_vector_set (odd_index, i - fullstartx, ((i/64)%2 == 0) ? 0 : 1);
1202  else cpl_vector_set (odd_index, i - fullstartx, 0);
1203  }
1204 
1205  CPLCHECK_NUL ("Cannot buid odd_index");
1206 
1207  /* Save the 2D coeficient of each polar */
1208  cpl_polynomial ** coef_poly = cpl_calloc (npol, sizeof (cpl_polynomial*));
1209 
1210  /* Loop on polarisation */
1211  for (int pol = 0; pol < npol; pol++) {
1212 
1213  /* Prepare the 2D coordinates
1214  * and values to fit */
1215  cpl_vector * coord_X = cpl_vector_new (nbase * nwave);
1216  cpl_vector * coord_Y = cpl_vector_new (nbase * nwave);
1217 
1218  cpl_vector * all_wavelength = cpl_vector_new (nbase * nwave);
1219  cpl_vector * all_wavechi2 = cpl_vector_new (nbase * nwave);
1220  cpl_vector * all_valid = cpl_vector_new (nbase * nwave);
1221 
1222  /*
1223  * Loop on base and wave to get all
1224  * wavelenght and coordinates
1225  */
1226  for (int base = 0; base < nbase; base ++) {
1227 
1228  /* Mean position of this baseline */
1229  int iA = gravi_get_region (detector_table, base, 'A', pol);
1230  int iB = gravi_get_region (detector_table, base, 'B', pol);
1231  int iC = gravi_get_region (detector_table, base, 'C', pol);
1232  int iD = gravi_get_region (detector_table, base, 'D', pol);
1233 
1234  /* WAVE_FIBRE data */
1235  sprintf (name, "BASE_%s_%s", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol,npol));
1236  cpl_array * wavelength = cpl_table_get_data_array (wavefibre_table, name)[0];
1237 
1238  sprintf (name, "BASE_%s_%s_CHI2", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol,npol));
1239  cpl_array * wavechi2 = cpl_table_get_data_array (wavefibre_table, name)[0];
1240 
1241  /* Loop on wave */
1242  for (cpl_size wave = 0; wave < nwave; wave++) {
1243  cpl_size pos = base * nwave + wave;
1244 
1245  /* Get the values */
1246  int nv = 0;
1247  double wave_value = cpl_array_get (wavelength, wave, &nv);
1248  double chi2_value = cpl_array_get (wavechi2, wave, &nv);
1249 
1250  /* The FT accept all channel for the 2D fit */
1251  cpl_vector_set (all_valid, pos, 1);
1252 
1253  if (nwave > 100) {
1254  /* SC HIGH and MEDIUM */
1255  if ((chi2_value > M_PI_4 ||
1256  wave_value < 2.01e-6 ||
1257  wave_value > 2.43e-6))
1258  cpl_vector_set (all_valid, pos, 0);
1259  }
1260  else if (nwave > GRAVI_LBD_FTSC) {
1261  /* SC LOW */
1262  if ((chi2_value > M_PI_4 ||
1263  //wave_value < 2.01e-6 ||
1264  wave_value < 1.99e-6 ||
1265  wave_value > 2.5e-6 ||
1266  wave == 0 || wave == nwave-1))
1267  cpl_vector_set (all_valid, pos, 0);
1268  }
1269 
1270  /* Set the wavelength */
1271  cpl_vector_set (all_wavelength, pos, wave_value);
1272  cpl_vector_set (all_wavechi2, pos, chi2_value);
1273 
1274  /* Set the X position as the mean of the 4 regions */
1275  cpl_vector_set (coord_X, pos, (double)(iA + iB + iC + iD) / 4.);
1276 
1277  /* Set the Y position. Add a shift of 0.15 pixels
1278  * on odd output for SC detector */
1279  cpl_vector_set (coord_Y, pos, wave + cpl_vector_get (odd_index, wave)*0.15);
1280 
1281  } /* End loop on wave */
1282  } /* End loop on base */
1283 
1284  CPLCHECK_NUL ("Error get all wavelength");
1285 
1286  /*
1287  * Reformat the valid point in vector and matrix
1288  */
1289  cpl_size nvalid = cpl_vector_get_sum (all_valid);
1290 
1291  cpl_msg_info (cpl_func, "Remove %lld/%lld bad wave",
1292  nbase*nwave - nvalid, nbase * nwave);
1293 
1294  cpl_vector * vector = cpl_vector_new (nvalid);
1295  cpl_matrix * matrix = cpl_matrix_new (2, nvalid);
1296 
1297  for (cpl_size c = 0, i = 0 ; i < nwave * nbase; i ++) {
1298  if (!cpl_vector_get (all_valid, i)) continue;
1299  cpl_vector_set (vector, c, cpl_vector_get (all_wavelength, i));
1300  cpl_matrix_set (matrix, 0, c, cpl_vector_get (coord_X, i));
1301  cpl_matrix_set (matrix, 1, c, cpl_vector_get (coord_Y, i));
1302  c++;
1303  }
1304 
1305  CPLCHECK_NUL ("Error cropping all wavelength");
1306 
1307  /*
1308  * Perform a 2D fit with a polynomial model
1309  * between the position and the wavelength = F(y, x)
1310  */
1311  cpl_size deg2d[2] = {2, 3};
1312  if ( (nwave < 20) && (nwave > 8) ) {deg2d[0] = 2; deg2d[1] = 7;} // FIXME Useless line ?
1313  deg2d[0] = spatial_order;
1314  deg2d[1] = spectral_order;
1315 
1316  cpl_msg_info (cpl_func, "Fit a 2d polynomial {%lli..%lli} to the wavelengths map", deg2d[0], deg2d[1]);
1317 
1318  cpl_polynomial * fit2d = cpl_polynomial_new (2);
1319  cpl_polynomial_fit (fit2d, matrix, NULL, vector, NULL, CPL_TRUE, NULL, deg2d);
1320  coef_poly[pol] = fit2d;
1321 
1322  CPLCHECK_NUL ("Cannot fit 2D");
1323 
1324  /*
1325  * Compute residuals
1326  */
1327  double rechisq = 0.0;
1328  cpl_vector * residuals = cpl_vector_new (nvalid);
1329  cpl_vector_fill_polynomial_fit_residual (residuals, vector, NULL, fit2d, matrix, &rechisq);
1330  *rms_residuals += cpl_vector_get_stdev(residuals)/npol;
1331  FREE (cpl_vector_delete, residuals);
1332  CPLCHECK_NUL ("Cannot compute residuals");
1333 
1334 
1335  FREE (cpl_matrix_delete, matrix);
1336  FREE (cpl_vector_delete, vector);
1337  FREE (cpl_vector_delete, all_wavelength);
1338  FREE (cpl_vector_delete, all_wavechi2);
1339  FREE (cpl_vector_delete, all_valid);
1340  FREE (cpl_vector_delete, coord_X);
1341  FREE (cpl_vector_delete, coord_Y);
1342 
1343  } /* End loop on polarisation */
1344 
1345 
1346  /*
1347  * Create and fill the interpolated WAVE_DATA table
1348  */
1349 
1350  cpl_table * wavedata_table = cpl_table_new (1);
1351  cpl_vector * pos = cpl_vector_new (2);
1352  cpl_array * value = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1353 
1354  for (cpl_size region = 0 ; region < n_region; region ++) {
1355 
1356  int pol = gravi_region_get_pol (detector_table, region);
1357  /* Loop on wave to evaluate the 2D polynome */
1358  for (cpl_size wave = 0; wave < nwave; wave ++) {
1359 
1360  /* Evaluate */
1361  cpl_vector_set (pos, 0, region);
1362  cpl_vector_set (pos, 1, wave + cpl_vector_get (odd_index, wave)*0.15);
1363 
1364  double result = cpl_polynomial_eval (coef_poly[pol], pos);
1365  cpl_array_set (value, wave, result);
1366 
1367  }
1368  /* ensure cresent wavelength */
1369  double previous_wave = cpl_array_get(value, nwave/2, NULL);
1370  for (cpl_size wave_loop = nwave/2 ; wave_loop >= 0 ; wave_loop --){
1371  if (previous_wave < cpl_array_get(value, wave_loop, NULL))
1372  cpl_array_set(value, wave_loop, previous_wave);
1373  else previous_wave = cpl_array_get(value, wave_loop, NULL);
1374  }
1375 
1376  previous_wave = cpl_array_get(value, nwave/2, NULL);
1377  for (cpl_size wave_loop = nwave/2 ; wave_loop < nwave ; wave_loop ++){
1378  if (previous_wave > cpl_array_get(value, wave_loop, NULL))
1379  cpl_array_set(value, wave_loop, previous_wave);
1380  else previous_wave = cpl_array_get(value, wave_loop, NULL);
1381  }
1382 
1383  /* Add column */
1384  char * data_x = GRAVI_DATA[region];
1385  cpl_table_new_column_array (wavedata_table, data_x, CPL_TYPE_DOUBLE, nwave);
1386  cpl_table_set_array (wavedata_table, data_x, 0, value);
1387 
1388  } /* End loop on regions */
1389 
1390  /* Delete allocations */
1391  FREE (cpl_vector_delete, pos);
1392  FREE (cpl_array_delete, value);
1393  FREELOOP (cpl_polynomial_delete, coef_poly, npol);
1394  FREE (cpl_vector_delete, odd_index);
1395 
1396  gravi_msg_function_exit(1);
1397  return wavedata_table;
1398 }
1399 
1400 
1401 
1402 /*----------------------------------------------------------------------------*/
1403 /*
1404  * @brief TBD
1405  */
1406 /*----------------------------------------------------------------------------*/
1407 
1408 cpl_table * gravi_wave_fit_individual (cpl_table * wave_individual_table,
1409  cpl_table * weight_individual_table,
1410  cpl_table * wave_fitted_table,
1411  cpl_table * opd_table,
1412  cpl_table * spectrum_table,
1413  cpl_table * detector_table,
1414  cpl_size fullstartx,
1415  double n0, double n1, double n2,
1416  double * rms_residuals)
1417 {
1418 
1419  gravi_msg_function_start(1);
1420 
1421  cpl_ensure (wave_individual_table, CPL_ERROR_NULL_INPUT, NULL);
1422  cpl_ensure (weight_individual_table, CPL_ERROR_NULL_INPUT, NULL);
1423  cpl_ensure (wave_fitted_table, CPL_ERROR_NULL_INPUT, NULL);
1424  cpl_ensure (opd_table, CPL_ERROR_NULL_INPUT, NULL);
1425  cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
1426  cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
1427 
1428  /* Get the number of wavelength, region, polarisation... */
1429  cpl_size nwave = cpl_table_get_column_depth (spectrum_table, "DATA1");
1430  cpl_size n_region = cpl_table_get_nrow (detector_table);
1431  cpl_size nrow = cpl_table_get_nrow (spectrum_table);
1432  int npol = (n_region > 24 ? 2 : 1);
1433  int nbase = 6;
1434  cpl_size nwave_ref=3000;
1435  if (nwave<10) nwave_ref=600;
1436 
1437  CPLCHECK_NUL ("Cannot buid odd_index");
1438 
1439  /* Get OPD Table */
1440 
1441 
1442 
1443  /* create temporary variables */
1444 
1445  cpl_array * wave_individual_array = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1446  cpl_array * weight_individual_array = cpl_array_new (nwave, CPL_TYPE_DOUBLE);
1447  cpl_matrix * data_flux_matrix = cpl_matrix_new (nrow, nwave);
1448  cpl_matrix * vis_to_flux_matrix = cpl_matrix_new (nrow, 3);
1449  cpl_matrix * signal_matrix = cpl_matrix_new (nwave, nwave_ref);
1450  cpl_matrix * residual_matrix = cpl_matrix_new (nwave, nwave_ref);
1451  cpl_array * wave_reference_array = cpl_array_new (nwave_ref, CPL_TYPE_DOUBLE);
1452 
1453  /*initialise arrays */
1454 
1455  cpl_matrix_fill_column(vis_to_flux_matrix,1,2);
1456  for (cpl_size wave_ref = 0; wave_ref < nwave_ref; wave_ref+=1)
1457  {
1458  // make matrix todo
1459  double wave_value=1.95e-6+wave_ref*0.6e-6/((double) nwave_ref);
1460  cpl_array_set_double(wave_reference_array,wave_ref,wave_value);
1461  }
1462 
1463  CPLCHECK_NUL ("Cannot initialize arrays for wavelength fit");
1464 
1465 
1466  for (cpl_size region = 0 ; region < n_region; region ++) {
1467 
1468  int base=gravi_region_get_base (detector_table, region);
1469  char * data_x = GRAVI_DATA[region];
1470 
1471  cpl_msg_info_overwritable (cpl_func, "Least square fitting of wavelength for region %s", data_x);
1472  // get data_flux_matrix
1473 
1474  for (cpl_size row = 0; row < nrow; row ++ ) {
1475  cpl_array * flux_array= cpl_table_get_data_array(spectrum_table,data_x)[row];
1476  for (cpl_size wave = 0; wave < nwave; wave ++) {
1477  cpl_matrix_set (data_flux_matrix, row, wave, cpl_array_get(flux_array,wave, NULL));
1478  }
1479  }
1480 
1481 
1482  for (cpl_size wave_ref = 0; wave_ref < nwave_ref; wave_ref+=1)
1483  {
1484  // make matrix
1485  double wave_value=cpl_array_get(wave_reference_array,wave_ref,NULL);
1486 
1487  for (cpl_size row = 0; row < nrow; row ++ ) {
1488  double opd = cpl_table_get (opd_table, "OPD", row*nbase+base, NULL);
1489  double coherence_loss=1;
1490  if (fabs(opd) > 1e-9)
1491  {
1492  if (nwave <30)
1493  {
1494  coherence_loss=sin(opd*19500)/(opd*19500);
1495  }
1496  }
1497  cpl_matrix_set(vis_to_flux_matrix,row,0,cos(opd*6.28318530718/wave_value)*coherence_loss);
1498  cpl_matrix_set(vis_to_flux_matrix,row,1,sin(opd*6.28318530718/wave_value)*coherence_loss);
1499  }
1500 
1501  cpl_matrix * coef_vis = cpl_matrix_solve_normal(vis_to_flux_matrix,data_flux_matrix); // coef_vis is 3xnwave
1502  cpl_matrix * data_flux_fit = cpl_matrix_product_create(vis_to_flux_matrix,coef_vis);
1503  cpl_matrix * residuals_fit = cpl_matrix_duplicate(data_flux_fit);
1504  cpl_matrix_subtract (residuals_fit,data_flux_matrix); //residuals_fit is nrowxnwave
1505 
1506  for (cpl_size wave = 0; wave < nwave; wave ++ ) {
1507  cpl_matrix * temp_matrix = cpl_matrix_extract_column (data_flux_fit, wave);
1508  cpl_matrix_set(signal_matrix,wave,wave_ref,cpl_matrix_get_stdev(temp_matrix));
1509  FREE (cpl_matrix_delete, temp_matrix);
1510  cpl_matrix * temp_matrix2 = cpl_matrix_extract_column (residuals_fit, wave);
1511  cpl_matrix_set(residual_matrix,wave,wave_ref,cpl_matrix_get_stdev(temp_matrix2));
1512  FREE (cpl_matrix_delete, temp_matrix2);
1513  }
1514 
1515  FREE (cpl_matrix_delete, coef_vis);
1516  FREE (cpl_matrix_delete, data_flux_fit);
1517  FREE (cpl_matrix_delete, residuals_fit);
1518 
1519  CPLCHECK_NUL ("Cannot do Matrix inversion to calculate optimum wavelength");
1520 
1521  }
1522 
1523  // get minimum chi2 and amplitude signal for each wave
1524  cpl_size wave_ref=1;
1525  cpl_size discarded=1;
1526  for (cpl_size wave = 0; wave < nwave; wave ++ ) {
1527 
1528  cpl_matrix * chi2_extract=cpl_matrix_extract_row(residual_matrix,wave);
1529 
1530  cpl_matrix_get_minpos(chi2_extract,&discarded, &wave_ref );
1531 
1532  double wave_value = cpl_array_get(wave_reference_array, wave_ref, NULL );
1533  double weight_value = cpl_matrix_get(signal_matrix, wave , wave_ref )/(0.1+cpl_matrix_get(residual_matrix, wave, wave_ref ));
1534 
1535  cpl_array_set_double (wave_individual_array, wave, wave_value);
1536  cpl_array_set_double (weight_individual_array, wave, weight_value);
1537 
1538  FREE (cpl_matrix_delete, chi2_extract);
1539 
1540  }
1541 
1542  /* Add column */
1543  cpl_table_new_column_array (wave_individual_table, data_x, CPL_TYPE_DOUBLE, nwave);
1544  cpl_table_set_array (wave_individual_table, data_x, 0, wave_individual_array);
1545  cpl_table_new_column_array (weight_individual_table, data_x, CPL_TYPE_DOUBLE, nwave);
1546  cpl_table_set_array (weight_individual_table, data_x, 0, weight_individual_array);
1547  cpl_table_new_column_array (wave_fitted_table, data_x, CPL_TYPE_DOUBLE, nwave);
1548  cpl_table_set_array (wave_fitted_table, data_x, 0, wave_individual_array);
1549  }
1550 
1551  CPLCHECK_NUL ("Cannot get individual wavelength for each pixel");
1552 
1553  FREE (cpl_array_delete ,wave_individual_array);
1554  FREE (cpl_array_delete ,weight_individual_array);
1555  FREE (cpl_array_delete ,wave_reference_array);
1556  FREE (cpl_matrix_delete ,data_flux_matrix);
1557  FREE (cpl_matrix_delete ,vis_to_flux_matrix);
1558  FREE (cpl_matrix_delete ,signal_matrix);
1559  FREE (cpl_matrix_delete ,residual_matrix);
1560 
1561  cpl_msg_info (cpl_func, "Now fitting polynomials on wavelength channels");
1562 
1563  cpl_matrix * coef_to_wave = cpl_matrix_new (n_region / npol,5);
1564  cpl_matrix * coef_to_wave_weight = cpl_matrix_new (n_region / npol,n_region / npol);
1565  cpl_matrix * wavelength = cpl_matrix_new(n_region / npol,1);
1566 
1567  // set coordinates
1568  for (cpl_size region = 0 ; region < n_region/ npol; region ++)
1569  {
1570  double mean_region = region - (n_region/npol-1)*0.5;
1571  cpl_matrix_set (coef_to_wave, region, 0, 1);
1572  cpl_matrix_set (coef_to_wave, region, 1, mean_region);
1573  cpl_matrix_set (coef_to_wave, region, 2, mean_region*mean_region);
1574  cpl_matrix_set (coef_to_wave, region, 3, mean_region*mean_region*mean_region);
1575  cpl_matrix_set (coef_to_wave, region, 4, mean_region*mean_region*mean_region*mean_region);
1576  }
1577 
1578  for (int pol = 0; pol < npol; pol++) {
1579 
1580  cpl_msg_info (cpl_func, "Looping for polyfit now, with pol: %i",(int) pol);
1581 
1582  for (cpl_size wave = 0; wave < nwave; wave ++) {
1583 
1584  // get the data for a common wave row
1585  for (cpl_size region = 0 ; region < n_region/ npol; region ++) {
1586 
1587  // Get the pointers to the table arrays
1588  char * data_x = GRAVI_DATA[region*npol+pol];
1589 
1590 
1591  const cpl_array * wave_array = cpl_table_get_array (wave_individual_table, data_x, 0);
1592  const cpl_array * weight_array = cpl_table_get_array (weight_individual_table, data_x, 0);
1593 
1594 
1595  cpl_matrix_set (wavelength, region, 0, cpl_array_get(wave_array,wave,NULL));
1596  double weight_value=cpl_array_get(weight_array,wave,NULL);
1597  cpl_matrix_set (coef_to_wave_weight, region, region, weight_value*weight_value);
1598 
1599  }
1600 
1601 
1602  cpl_matrix * coef_to_wave2 = cpl_matrix_product_create(coef_to_wave_weight,coef_to_wave);
1603  cpl_matrix * wavelength2 = cpl_matrix_product_create(coef_to_wave_weight,wavelength);
1604 
1605  // Fit a second order polynomial
1606  cpl_matrix * coeff = cpl_matrix_solve_normal(coef_to_wave2, wavelength2); // 5 x 1
1607  cpl_matrix * wavelength_fitted = cpl_matrix_product_create(coef_to_wave, coeff); // n_region / npol x 1
1608 
1609  CPLCHECK_NUL ("Cannot do Matrix inversion to calculate optimum polynomial for wavelength");
1610 
1611  // Write result to new wave table
1612  for (cpl_size region = 0 ; region < n_region/ npol; region ++) {
1613 
1614  char * data_x = GRAVI_DATA[region*npol+pol];
1615  //cpl_msg_info (cpl_func, "Writing: region %i",region);
1616  cpl_array * wave_array = cpl_table_get_data_array (wave_fitted_table, data_x)[0];
1617 
1618  cpl_array_set_double(wave_array,wave,cpl_matrix_get(wavelength_fitted,region,0));
1619 
1620  }
1621 
1622  FREE (cpl_matrix_delete ,coef_to_wave2);
1623  FREE (cpl_matrix_delete ,wavelength2);
1624  FREE (cpl_matrix_delete ,coeff);
1625  FREE (cpl_matrix_delete ,wavelength_fitted);
1626  }
1627 
1628  }
1629  FREE (cpl_matrix_delete ,coef_to_wave);
1630  FREE (cpl_matrix_delete ,coef_to_wave_weight);
1631  FREE (cpl_matrix_delete ,wavelength);
1632 
1633  CPLCHECK_NUL ("Cannot fit individual wavelength with 3rd order polynomial");
1634 
1635  cpl_msg_info (cpl_func, "Correcting for wavelength error");
1636 
1637 
1638  /* Correct the computed wavelength from the dispersion */
1639  for (cpl_size region = 0 ; region < n_region; region ++)
1640  {
1641  char * data_x = GRAVI_DATA[region];
1642  cpl_array * wavelength = cpl_table_get_data_array (wave_fitted_table, data_x)[0];
1643  cpl_size nwave = cpl_array_get_size (wavelength);
1644  for (cpl_size wave = 0 ; wave < nwave ; wave ++ ) {
1645 
1646  double result = cpl_array_get (wavelength, wave, NULL);
1647  double d_met = (result - LAMBDA_MET) / LAMBDA_MET;
1648  cpl_array_set (wavelength, wave, result * (n0 + n1*d_met + n2*d_met*d_met));
1649 
1650  }
1651  for (cpl_size wave = nwave/2 ; wave < nwave-1 ; wave ++ ) {
1652  double result = cpl_array_get (wavelength, wave, NULL);
1653  double result2 = cpl_array_get (wavelength, wave+1, NULL);
1654  if (result2<result+2e-10) {
1655  result2=result+2e-10;
1656  cpl_array_set (wavelength, wave+1, result2);
1657  }
1658  }
1659  for (cpl_size wave = nwave/2 ; wave > 0 ; wave -- ) {
1660  double result = cpl_array_get (wavelength, wave, NULL);
1661  double result2 = cpl_array_get (wavelength, wave-1, NULL);
1662  if (result2>result-2e-10) {
1663  result2=result-2e-10;
1664  cpl_array_set (wavelength, wave-1, result2);
1665  }
1666  }
1667  }
1668 
1669  CPLCHECK_NUL ("Error in correcting for dispersion");
1670 
1671  gravi_msg_function_exit(1);
1672  return wave_fitted_table;
1673 }
1674 
1675 
1676 
1677 
1678 /*----------------------------------------------------------------------------*/
1685 /*----------------------------------------------------------------------------*/
1686 
1687 cpl_error_code gravi_wave_correct_dispersion (cpl_table * wave_fibre,
1688  double n0, double n1, double n2)
1689 {
1690  gravi_msg_function_start(1);
1691  cpl_ensure_code (wave_fibre, CPL_ERROR_NULL_INPUT);
1692 
1693  int nbase = 6;
1694  char name[100];
1695 
1696  /* Get the number of polarisation */
1697  int npol = cpl_table_has_column (wave_fibre, "BASE_12_P") ? 2 : 1;
1698 
1699  /* Loop on columns in the table */
1700  for (int pol = 0; pol < npol; pol ++) {
1701  for (int base = 0; base < nbase; base ++) {
1702 
1703  /* Get data of this region */
1704  sprintf (name, "BASE_%s_%s", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol, npol));
1705  cpl_array * wavelength = cpl_table_get_data_array (wave_fibre, name)[0];
1706  CPLCHECK_MSG ("Cannot get data");
1707 
1708  /* Loop on wave */
1709  cpl_size nwave = cpl_array_get_size (wavelength);
1710  for (cpl_size wave = 0 ; wave < nwave ; wave ++ ) {
1711 
1712  /* Correct the computed wavelength from the dispersion */
1713  double result = cpl_array_get (wavelength, wave, NULL);
1714  double d_met = (result - LAMBDA_MET) / LAMBDA_MET;
1715  cpl_array_set (wavelength, wave, result * (n0 + n1*d_met + n2*d_met*d_met));
1716 
1717  }
1718  } /* End loop on base */
1719  } /* End loop on pol */
1720 
1721  gravi_msg_function_exit(1);
1722  return CPL_ERROR_NONE;
1723 }
1724 
1725 /*----------------------------------------------------------------------------*/
1731 /*----------------------------------------------------------------------------*/
1732 
1733 cpl_error_code gravi_wave_correct_color (gravi_data * vis_data)
1734 {
1735  gravi_msg_function_start(1);
1736  cpl_ensure_code (vis_data, CPL_ERROR_NULL_INPUT);
1737 
1738  cpl_propertylist * primary_header = gravi_data_get_header (vis_data);
1739  cpl_propertylist * oiwave_header = NULL;
1740  cpl_table * oiwave_table = NULL;
1741 
1742  /* For each type of data SC / FT */
1743  int ntype_data = 2;
1744  for (int type_data = 0; type_data < ntype_data ; type_data ++) {
1745 
1746  /* Loop on polarisation */
1747  int npol = gravi_pfits_get_pola_num (primary_header, type_data);
1748  for (int pol = 0 ; pol<npol ; pol++) {
1749  oiwave_table = cpl_table_duplicate ( gravi_data_get_oi_wave ( vis_data, type_data, pol, npol ) );
1750  oiwave_header = cpl_propertylist_duplicate ( gravi_data_get_oi_wave_plist ( vis_data, type_data, pol, npol ) );
1751 
1752  /* here you can do what you want on this duplicated oi_wave_table
1753  * to get OI_FLUX table :
1754  * cpl_table * oiflux_table = gravi_data_get_oi_flux(vis_data, type_data, pol, npol)
1755  * */
1756 
1757 
1758  /* save the new extension with name OI_WAVELENGTH_CORR */
1759  gravi_data_add_table(vis_data, oiwave_header, "OI_WAVELENGTH_CORR", oiwave_table);
1760 
1761  CPLCHECK_MSG("Cannot apply color wave correction");
1762  }
1763  /* End loop on polarisation */
1764  }
1765  /* End loop on data_type */
1766 
1767  gravi_msg_function_exit(1);
1768  return CPL_ERROR_NONE;
1769 }
1770 
1771 
1772 /*----------------------------------------------------------------------------*/
1784 /*----------------------------------------------------------------------------*/
1785 
1786 cpl_imagelist * gravi_wave_test_image (cpl_table * wavedata_table,
1787  cpl_table * wavefibre_table,
1788  cpl_table * profile_table,
1789  cpl_table * detector_table)
1790 {
1791  gravi_msg_function_start(1);
1792  cpl_ensure (wavedata_table, CPL_ERROR_NULL_INPUT, NULL);
1793  cpl_ensure (wavefibre_table, CPL_ERROR_NULL_INPUT, NULL);
1794  cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
1795  cpl_ensure (profile_table, CPL_ERROR_NULL_INPUT, NULL);
1796 
1797  int nv = 0;
1798  char name[100];
1799 
1800  cpl_size n_region = cpl_table_get_nrow (detector_table);
1801  int npol = (n_region > 24 ? 2 : 1);
1802 
1803  CPLCHECK_NUL ("Cannot get data");
1804 
1805  /* Init the output images */
1806  cpl_size sizex = cpl_table_get_column_dimension (profile_table, "DATA1", 0);
1807  cpl_size sizey = cpl_table_get_column_dimension (profile_table, "DATA1", 1);
1808 
1809  cpl_image * profilesum_image = cpl_image_new (sizex, sizey, CPL_TYPE_DOUBLE);
1810  cpl_image_fill_window (profilesum_image, 1, 1, sizex, sizey, 0.0);
1811 
1812  cpl_image * wave_image = cpl_image_new (sizex, sizey, CPL_TYPE_DOUBLE);
1813  cpl_image_fill_window (wave_image, 1, 1, sizex, sizey, 0.0);
1814 
1815  cpl_image * realwave_image = cpl_image_new (sizex, sizey, CPL_TYPE_DOUBLE);
1816  cpl_image_fill_window (realwave_image, 1, 1, sizex, sizey, 0.0);
1817 
1818  CPLCHECK_NUL ("Cannot prepare output image");
1819 
1820  /* Loop on regions */
1821  for (cpl_size region = 0 ; region < n_region; region ++) {
1822 
1823  /* Load the profile of this region as an image */
1824  cpl_imagelist * profile_imglist = gravi_imagelist_wrap_column (profile_table, GRAVI_DATA[region]);
1825  cpl_image * profile_image = cpl_imagelist_get (profile_imglist, 0);
1826 
1827  CPLCHECK_NUL ("Cannot get data");
1828 
1829  /* Sum all profils of all regions */
1830  cpl_image_add (profilesum_image, profile_image);
1831 
1832  /*
1833  * Define an image of each region with its WAVE_DATA
1834  */
1835  const cpl_array * wavelength;
1836  wavelength = cpl_table_get_array (wavedata_table, GRAVI_DATA[region], 0);
1837  CPLCHECK_NUL ("Cannot get data");
1838 
1839  for (cpl_size x = 0; x < sizex; x ++){
1840  for (cpl_size y = 0; y < sizey; y ++){
1841  if (cpl_image_get (profile_image, x+1, y+1, &nv) > 0.01)
1842  cpl_image_set (wave_image, x+1, y+1,
1843  cpl_array_get (wavelength, x, NULL));
1844  }
1845  }
1846  CPLCHECK_NUL ("Cannot make image of wave_map");
1847 
1848  /*
1849  * Define an image of each region with its WAVE_FIBER
1850  */
1851  int base = gravi_region_get_base (detector_table, region);
1852  int pol = gravi_region_get_pol (detector_table, region);
1853  sprintf (name, "BASE_%s_%s", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol, npol));
1854  wavelength = cpl_table_get_array (wavefibre_table, name, 0);
1855  CPLCHECK_NUL ("Cannot get data");
1856 
1857  for (cpl_size x = 0; x < sizex; x ++){
1858  for (cpl_size y = 0; y < sizey; y ++){
1859  if (cpl_image_get (profile_image, x+1, y+1, &nv) > 0.01)
1860  cpl_image_set (realwave_image, x+1, y+1,
1861  cpl_array_get (wavelength, x, NULL));
1862  }
1863  }
1864 
1865  gravi_imagelist_unwrap_images (profile_imglist);
1866  CPLCHECK_NUL ("Cannot make image of wave_fibre");
1867  } /* End loop on regions */
1868 
1869  /* Create an imagelist */
1870  cpl_imagelist * testwave_imglist = cpl_imagelist_new ();
1871  cpl_imagelist_set (testwave_imglist, wave_image, 0);
1872  cpl_imagelist_set (testwave_imglist, profilesum_image, 1);
1873  cpl_imagelist_set (testwave_imglist, realwave_image, 2);
1874 
1875  gravi_msg_function_exit(1);
1876  return testwave_imglist;
1877 }
1878 
1879 /*----------------------------------------------------------------------------*/
1889 /*----------------------------------------------------------------------------*/
1890 
1891 cpl_error_code gravi_wave_qc (gravi_data * wave_map, gravi_data * profile_map)
1892 {
1893  gravi_msg_function_start(1);
1894  cpl_ensure_code (wave_map, CPL_ERROR_NULL_INPUT);
1895 
1896  char name[100];
1897 
1898  cpl_propertylist * wave_header = gravi_data_get_header (wave_map);
1899 
1900  /* Loop on extensions (thus SC/FT) */
1901  for (int type_data = 0; type_data < 2; type_data ++ ) {
1902 
1903  /* Get WAVE_FIBRE and IMAGING_DETECTOR tables */
1904  cpl_table * detector_table = gravi_data_get_imaging_detector (wave_map, type_data);
1905  cpl_table * wave_data = gravi_data_get_wave_data (wave_map, type_data);
1906  cpl_size n_region = cpl_table_get_nrow (detector_table);
1907 
1908  /* Get the full STARTX */
1909  cpl_propertylist * plist = gravi_data_get_wave_data_plist (wave_map, type_data);
1910  int fullstartx = gravi_pfits_get_fullstartx (plist);
1911 
1912  CPLCHECK_MSG ("Cannot get data");
1913 
1914  /*
1915  * Compute the min and max wave over all regions
1916  */
1917  double minwave = -1e10;
1918  double maxwave = 1e10;
1919 
1920  for (int region = 0; region < n_region; region++) {
1921  const cpl_array * wavelength = cpl_table_get_array (wave_data, GRAVI_DATA[region], 0);
1922 
1923  minwave = CPL_MAX (minwave, cpl_array_get_min (wavelength));
1924  maxwave = CPL_MIN (maxwave, cpl_array_get_max (wavelength));
1925  } /* End loop on regions */
1926 
1927  cpl_msg_info (cpl_func,"%s = %g [m]", QC_MINWAVE(type_data), minwave);
1928  cpl_msg_info (cpl_func,"%s = %g [m]", QC_MAXWAVE(type_data), maxwave);
1929  cpl_propertylist_update_double (wave_header, QC_MINWAVE(type_data), minwave);
1930  cpl_propertylist_update_double (wave_header, QC_MAXWAVE(type_data), maxwave);
1931 
1932  CPLCHECK_MSG ("Cannot compute minwave or maxwave");
1933 
1934 
1935  /*
1936  * Compute the pixel position of QC specified argon wavelength.
1937  * Having two lines allow to check the position and dispersion
1938  */
1939  int qc_reg[3] = {0,23,47};
1940  double qc_wave[2] = {2.099184e-6, 2.313952e-6};
1941 
1942  /* Loop on regions */
1943  for (int reg = 0; reg < (n_region > 24 ? 3 : 2); reg++) {
1944  cpl_size region = qc_reg[reg];
1945 
1946  const cpl_array * wavelength = cpl_table_get_array (wave_data, GRAVI_DATA[region], 0);
1947  cpl_size nwave = cpl_array_get_size (wavelength);
1948  const double * wave_tab = cpl_array_get_data_double_const (wavelength);
1949  CPLCHECK_MSG ("Cannot get wavelength data");
1950 
1951  /* Loop on the two argon lines */
1952  for (int iqc = 0 ; iqc < 2 ; iqc++) {
1953  cpl_size l2 = 0;
1954  if ( wave_tab[0] < wave_tab[nwave-1]) {while (wave_tab[l2] < qc_wave[iqc]) l2 ++;}
1955  else {while (wave_tab[l2] > qc_wave[iqc]) l2 ++;}
1956 
1957  if (l2-1 < 0 || l2 > nwave-1) {
1958  cpl_msg_error (cpl_func, "Cannot find the QC position for lbd=%g", qc_wave[iqc]);
1959  continue;
1960  }
1961 
1962  /* Position on full detector */
1963  double qc_pos = 0.0;
1964  qc_pos = fullstartx + (l2-1) + (qc_wave[iqc] - wave_tab[l2-1]) / (wave_tab[l2] - wave_tab[l2-1]);
1965 
1966  sprintf (name, "ESO QC REFWAVE%i", iqc+1);
1967  cpl_propertylist_update_double (wave_header, name, qc_wave[iqc]);
1968  cpl_propertylist_set_comment (wave_header, name, "[m] value of ref wave");
1969 
1970  sprintf (name, "ESO QC REFPOS%i %s%lli", iqc+1, GRAVI_TYPE(type_data),region+1);
1971  cpl_propertylist_update_double (wave_header, name, qc_pos);
1972  cpl_propertylist_set_comment (wave_header, name, "[pix] position of ref wave");
1973 
1974  cpl_msg_info (cpl_func, "%s = %f [pix] for %e [m]", name, qc_pos, qc_wave[iqc]);
1975 
1976  CPLCHECK_MSG ("Cannot set QC");
1977  }
1978  /* End loop on the 2 argon lines */
1979  } /* End loop on 2 or 3 regions */
1980 
1981  } /* End loop on SC / FT */
1982 
1983 
1984  /*
1985  * Create the test image for SC (only used for debug)
1986  */
1987  cpl_table * wavedata_table = gravi_data_get_wave_data (wave_map, GRAVI_SC);
1988  cpl_table * wavefibre_table = gravi_data_get_wave_fibre (wave_map, GRAVI_SC);
1989  cpl_table * profile_table = gravi_data_get_table (profile_map, GRAVI_PROFILE_DATA_EXT);
1990  cpl_table * detector_table = gravi_data_get_imaging_detector (wave_map, GRAVI_SC);
1991 
1992  cpl_imagelist * testwave_imglist;
1993  testwave_imglist = gravi_wave_test_image (wavedata_table, wavefibre_table,
1994  profile_table, detector_table);
1995  CPLCHECK_MSG ("Cannot compute TEST_WAVE");
1996 
1997  /* Set TEST_WAVE in output */
1998  cpl_propertylist * plist = gravi_data_get_plist (profile_map, GRAVI_PROFILE_DATA_EXT);
1999  gravi_data_add_cube (wave_map, cpl_propertylist_duplicate (plist),
2000  "TEST_WAVE", testwave_imglist);
2001 
2002  CPLCHECK_MSG ("Cannot set data");
2003 
2004  gravi_msg_function_exit(1);
2005  return CPL_ERROR_NONE;
2006 }
2007 
2008 /*----------------------------------------------------------------------------*/
2036 /*----------------------------------------------------------------------------*/
2037 
2038 cpl_error_code gravi_compute_wave (gravi_data * wave_map,
2039  gravi_data * spectrum_data,
2040  int type_data, const cpl_parameterlist * parlist)
2041 {
2042  gravi_msg_function_start(1);
2043  cpl_ensure_code (wave_map, CPL_ERROR_NULL_INPUT);
2044  cpl_ensure_code (spectrum_data, CPL_ERROR_NULL_INPUT);
2045  cpl_ensure_code (type_data == GRAVI_SC || type_data == GRAVI_FT,
2046  CPL_ERROR_ILLEGAL_INPUT);
2047 
2048  /* Get headers */
2049  cpl_propertylist * raw_header = gravi_data_get_header (spectrum_data);
2050  cpl_propertylist * wave_header = gravi_data_get_header (wave_map);
2051 
2052  /* Dump full header of wave data */
2053  if (type_data == GRAVI_FT) {
2054  cpl_propertylist_append (wave_header, raw_header);
2055  }
2056 
2057 
2058  /* Copy IMAGING_DETECTOR in output WAVE */
2059  const char * extname = (type_data == GRAVI_SC) ? GRAVI_IMAGING_DETECTOR_SC_EXT : GRAVI_IMAGING_DETECTOR_FT_EXT;
2060  gravi_data_copy_ext (wave_map, spectrum_data, extname);
2061 
2062  /*
2063  * Compute WAVE_FIBRE map
2064  */
2065  cpl_msg_info (cpl_func, "Compute WAVE_FIBRE for %s", GRAVI_TYPE(type_data));
2066 
2067  cpl_table * spectrum_table = gravi_data_get_spectrum_data (spectrum_data, type_data);
2068  cpl_table * detector_table = gravi_data_get_imaging_detector (spectrum_data, type_data);
2069  cpl_table * opd_table = gravi_data_get_table (spectrum_data, type_data==GRAVI_SC ? "OPD_SC" : "OPD_FT");
2070 
2071  cpl_table * wavefibre_table;
2072  wavefibre_table = gravi_wave_fibre (spectrum_table, detector_table, opd_table);
2073  CPLCHECK_MSG ("Cannot compute the WAVE_FIBRE");
2074 
2075  /*
2076  * Correct WAVE_FIBRE from the dispersion model
2077  */
2078  cpl_msg_info (cpl_func, "Correct dispersion in WAVE_FIBRE of %s", GRAVI_TYPE(type_data));
2079 
2080  /* Hardcoded values to correct for the fiber dispersion */
2081  double n0 = 1.0, n1 = -0.0165448, n2 = 0.00256002;
2082  cpl_msg_info (cpl_func,"Rescale wavelengths with dispersion (%g,%g,%g)",n0,n1,n2);
2083 
2084  cpl_propertylist_update_string (wave_header, "ESO QC WAVE_CORR", "lbd*(N0+N1*(lbd-lbd0)/lbd0+N2*(lbd-lbd0)^2/lbd0^2)");
2085  cpl_propertylist_update_double (wave_header, "ESO QC WAVE_CORR N0", n0);
2086  cpl_propertylist_update_double (wave_header, "ESO QC WAVE_CORR N1", n1);
2087  cpl_propertylist_update_double (wave_header, "ESO QC WAVE_CORR N2", n2);
2088  CPLCHECK_MSG ("Cannot set Keywords");
2089 
2090  gravi_wave_correct_dispersion (wavefibre_table, n0, n1, n2);
2091  CPLCHECK_MSG ("Cannot correct dispersion in wave");
2092 
2093  /*
2094  * Fit the 2d dispersion WAVE_FIBRE into WAVE_DATA
2095  */
2096  cpl_msg_info (cpl_func,"Fit WAVE_DATA map for %s", GRAVI_TYPE(type_data));
2097 
2098  /* Get the fullstartx */
2099  cpl_propertylist * spectrum_plist = gravi_data_get_spectrum_data_plist (spectrum_data, type_data);
2100  int fullstartx = gravi_pfits_get_fullstartx (spectrum_plist);
2101 
2102  /* Interpolate table 2D */
2103  cpl_table * wavedata_table;
2104  int spatial_order=2; // default spatial order
2105  int spectral_order=3; // default spectral order
2106  double rms_residuals;
2107  if (type_data == GRAVI_FT && gravi_param_get_bool(parlist, "gravity.calib.force-wave-ft-equal")) {
2108  spatial_order = 0;
2109  cpl_msg_info (cpl_func, "Option force-waveFT-equal applied");
2110  }
2111 // Keep default value
2112 // if (type_data == GRAVI_SC) {
2113 // spectral_order = gravi_param_get_int(parlist, "gravity.calib.wave-spectral-order");
2114 // cpl_msg_info (cpl_func, "Option set_spectral order to %d", spectral_order);
2115 // }
2116 
2117  wavedata_table = gravi_wave_fit_2d (wavefibre_table,
2118  detector_table,
2119  fullstartx, spatial_order, spectral_order, &rms_residuals);
2120  cpl_propertylist_update_double (wave_header, QC_RMS_RESIDUALS(type_data), rms_residuals);
2121 
2122  CPLCHECK_MSG ("Cannot fit 2d data");
2123 
2124  /*
2125  * New Wavelength interpolation made by sylvestre on January 30 2018
2126  */
2127 
2128  if (type_data == GRAVI_SC && !strcmp (gravi_pfits_get_spec_res(raw_header), "LOW"))
2129  {
2130  cpl_msg_info (cpl_func, "Additional Wavelength Fit");
2131  cpl_table * wave_individual_table = cpl_table_new (1);
2132  cpl_table * weight_individual_table = cpl_table_new (1);
2133  cpl_table * wave_fitted_table = cpl_table_new (1);
2134 
2135  gravi_wave_fit_individual (wave_individual_table,
2136  weight_individual_table,
2137  wave_fitted_table,
2138  opd_table,
2139  spectrum_table,
2140  detector_table,
2141  fullstartx,
2142  n0,n1,n2,
2143  &rms_residuals);
2144 
2145  cpl_msg_info (cpl_func,"Add tables in wave_map");
2146 
2147  gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2148  "WAVE_INDIV_SC", wave_individual_table);
2149  gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2150  "WAVE_WEIGHT_SC", weight_individual_table);
2151  gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2152  "WAVE_FITTED_SC", wavedata_table);
2153 
2154  gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2155  GRAVI_WAVE_FIBRE_EXT(type_data), wavefibre_table);
2156  gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2157  GRAVI_WAVE_DATA_EXT(type_data), wave_fitted_table);
2158  } else {
2159  /*
2160  * Add the WAVE_FIBRE and WAVE_DATA table in the wave_map
2161  */
2162  cpl_msg_info (cpl_func,"Add WAVE_FIBRE and WAVE_DATA in wave_map");
2163 
2164  gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2165  GRAVI_WAVE_FIBRE_EXT(type_data), wavefibre_table);
2166 
2167  gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2168  GRAVI_WAVE_DATA_EXT(type_data), wavedata_table);
2169 
2170  }
2171 
2172 
2173 
2174 
2175  CPLCHECK_MSG ("Cannot set data");
2176 
2177  gravi_msg_function_exit(1);
2178  return CPL_ERROR_NONE;
2179 }
2180 
2181 
cpl_error_code gravi_wave_qc(gravi_data *wave_map, gravi_data *profile_map)
Compute the QC parameters of the WAVE product.
Definition: gravi_wave.c:1891
int gravi_region_get_pol(cpl_table *imaging_detector, int region)
Return the polarisation id of a region.
Definition: gravi_utils.c:445
int gravi_region_get_base(cpl_table *imaging_detector, int region)
Return the base of a region.
Definition: gravi_utils.c:361
cpl_error_code gravi_wave_correct_color(gravi_data *vis_data)
Create a OI_WAVELENGTH_CORR table with color corrected wavelength.
Definition: gravi_wave.c:1733
cpl_error_code gravi_data_copy_ext(gravi_data *output, gravi_data *input, const char *name)
Copy extensions from one data to another.
Definition: gravi_data.c:1325
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_base_sign(cpl_table *imaging_detector, int base)
Return the sign of a base by looking at the PORT order.
Definition: gravi_utils.c:399
cpl_table * gravi_metrology_create(cpl_table *metrology_table, cpl_propertylist *header)
Create the VIS_MET table.
cpl_error_code gravi_pfits_add_check(cpl_propertylist *header, char *msg)
Add a QC.CHECK keyword to the header.
Definition: gravi_pfits.c:1398
cpl_table * gravi_wave_fit_2d(cpl_table *wavefibre_table, cpl_table *detector_table, 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...
Definition: gravi_wave.c:1177
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.
Definition: gravi_cpl.c:1870
cpl_vector * gravi_ellipse_phase_create(cpl_vector *vectCA, cpl_vector *vectDB, cpl_vector *envelope)
Compute the phase atan{X&#39;,Y&#39;}, unwraped from first sample.
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.
Definition: gravi_wave.c:986
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.
Definition: gravi_wave.c:1687
cpl_error_code gravi_wave_compute_opds(gravi_data *spectrum_data, cpl_table *met_table, enum gravi_detector_type det_type)
Recover the OPD modulation from a spectrum_data and vismet_table.
Definition: gravi_wave.c:822
cpl_table * gravi_compute_argon_wave(cpl_table *spectrum_table)
Compute a WAVE calibration from the ARGON data (SC only)
Definition: gravi_wave.c:149
cpl_table * gravi_data_get_table(gravi_data *self, const char *extname)
Return a pointer on a table extension by its EXTNAME.
Definition: gravi_data.c:1716
cpl_vector * gravi_compute_envelope(const cpl_vector *opd, int wave, int nwave)
Compute the envelope value.
Definition: gravi_utils.c:1081
int gravi_get_region(cpl_table *img_det, int base, char phase, int pol)
Find the region matching base, phase and pol.
Definition: gravi_utils.c:584
cpl_error_code gravi_compute_wave(gravi_data *wave_map, gravi_data *spectrum_data, int type_data, const cpl_parameterlist *parlist)
Create the WAVE calibration map.
Definition: gravi_wave.c:2038
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.
Definition: gravi_wave.c:1786
cpl_error_code gravi_imagelist_unwrap_images(cpl_imagelist *imglist)
Unwrap an imagelist an all its images.
Definition: gravi_cpl.c:1496
cpl_error_code gravi_opds_correct_closures(cpl_table *opd_sc, const char *name)
Correct the input table to ensure constant closure phase.
Definition: gravi_wave.c:354
cpl_error_code gravi_metrology_drs(cpl_table *metrology_table, cpl_table *vismet_table, cpl_propertylist *header)
Fill the VIS_MET table with the DRS algorithm.
cpl_propertylist * gravi_data_get_plist(gravi_data *self, const char *extname)
Get the propertylist from EXTNAME.
Definition: gravi_data.c:1684
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_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...
Definition: gravi_cpl.c:2284
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.
Definition: gravi_cpl.c:2381
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.
Definition: gravi_data.c:1972
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.
Definition: gravi_wave.c:566
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.
Definition: gravi_wave.c:698
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.
Definition: gravi_data.c:1909
cpl_error_code gravi_vis_create_met_ft(cpl_table *vis_FT, cpl_table *vis_MET)
Compute the resampled MET signal for each FT DIT per base.
cpl_imagelist * gravi_imagelist_wrap_column(cpl_table *table_data, const char *data_x)
Wrap a column array of a table into an imagelist.
Definition: gravi_cpl.c:1526
cpl_vector * gravi_opds_fit_opdmet(cpl_table *opd_ft, double lbd_met)
Compute the absolute scaling coefficients of SC and FT OPDs.
Definition: gravi_wave.c:459