GRAVI Pipeline Reference Manual 1.8.0
Loading...
Searching...
No Matches
gravi_wave.c
Go to the documentation of this file.
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 5
80
81/*-----------------------------------------------------------------------------
82 Private prototypes
83 -----------------------------------------------------------------------------*/
84
85cpl_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
91cpl_table * gravi_opds_calibration (cpl_table * spectrum_table,
92 cpl_table * detector_table,
93 cpl_table * opdguess_table);
94
95cpl_error_code gravi_opds_correct_closures (cpl_table * opd_sc,
96 const char *name);
97
98cpl_vector * gravi_opds_fit_opdmet (cpl_table *opd_ft, double lbd_met);
99
100cpl_table * gravi_wave_fibre (cpl_table * spectrum_table,
101 cpl_table * detector_table,
102 cpl_table * opd_table);
103
104cpl_table * gravi_wave_fit_2d (cpl_table * wavefibre_table,
105 cpl_table * detector_table,
106 gravi_data * wave_param,
107 cpl_size fullstartx,
108 int spatial_order,
109 int spectral_order,
110 double * rms_residuals);
111
112cpl_table * gravi_wave_fit_individual (cpl_table * wave_individual_table,
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);
119
120cpl_error_code gravi_wave_correct_dispersion (cpl_table * wave_fibre,
121 double n0, double n1, double n2);
122
123cpl_imagelist * gravi_wave_test_image (cpl_table * wavedata_table,
124 cpl_table * wavefibre_table,
125 cpl_table * profile_table,
126 cpl_table * detector_table);
127
128
129
130
131/*-----------------------------------------------------------------------------
132 Functions code
133 -----------------------------------------------------------------------------*/
134
135/*----------------------------------------------------------------------------*/
146/*----------------------------------------------------------------------------*/
147
148cpl_table * gravi_compute_argon_wave (cpl_table * spectrum_table)
149{
151 cpl_ensure (spectrum_table, CPL_ERROR_NULL_INPUT, NULL);
152
153 /* Get the SC spectrums */
154 cpl_size nwave = cpl_table_get_column_depth (spectrum_table,"DATA1");
155 cpl_size nregion = gravi_spectrum_get_nregion (spectrum_table);
156 CPLCHECK_NUL ("Cannot get data");
157
158 cpl_msg_info (cpl_func, "spectrum_table = %lld",cpl_table_get_column_depth (spectrum_table, "DATA1"));
159
160 /* Prepare the outputs */
161 cpl_table * wave_data_sc = cpl_table_new (1);
162 CPLCHECK_NUL ("Cannot create output");
163
164 /* Init the list of lines*/
165 cpl_size nlines = 14;
166 double plines_str[] = {60,
167 160,
168 90,
169 1800,
170 50,
171 30,
172 30,
173 1500,
174 130,
175 510,
176 240,
177 169,
178 64,
179 73
180 };
181
182 double plines_lbd[] = {1.982291,
183 1.997118,
184 2.032256,
185 2.062186,
186 2.065277,
187 2.073922,
188 2.081672,
189 2.099184,
190 2.133871,
191 2.154009,
192 2.208321,
193 2.313952,
194 2.385154,
195 2.397306
196 };
197
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);
201
202 /* Init the slit and lamp model */
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);
208
209 /* Starting point of dispersion */
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;
213
214 /* Fill the dispersion polynomial */
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]);
218 }
219
220 /* QC parameters */
221 double minwave = -1e10;
222 double maxwave = +1e10;
223
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);
227
228 /* Prepare output table */
229 cpl_table * fit_table = cpl_table_new (nregion);
230 gravi_table_new_column_array (fit_table, "DATA", "adu", CPL_TYPE_DOUBLE, nwave);
231 gravi_table_new_column_array (fit_table, "DATA_MODEL", "adu", CPL_TYPE_DOUBLE, nwave);
232 gravi_table_new_column_array (fit_table, "DATA_WAVE", "adu", CPL_TYPE_DOUBLE, nwave);
233
234 /* The starting point will be the fit of previous region */
235 cpl_polynomial * dispersion = cpl_polynomial_duplicate (dispersion0);
236
237 /* Loop on regions */
238 for (cpl_size reg = 0; reg < nregion ; reg ++ ) {
239 cpl_msg_info (cpl_func,"Fit region %lld over %lld", reg+1, nregion);
240
241 /* Copy the spectrum of this region into a vector */
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));
246 }
247
248 CPLCHECK_NUL ("Cannot prepare data");
249
250 /* Parameters: Explore a large hsize [pixel]
251 * only for the first region, since we keep
252 * the starting point */
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;
257
258 /* Use the log of spectrum to put
259 * less weight on the intensity */
260 irplib_polynomial_find_1d_from_correlation_all (dispersion, maxdeg, spectrum,
261 nmaxima, linelim,
262 (irplib_base_spectrum_model*)model,
263 &irplib_vector_fill_logline_spectrum,
264 pixtol, pixstep,
265 hsize, maxite, maxfail, maxcont,
266 CPL_FALSE,&pxc);
267 CPLCHECK_NUL ("Cannot fit polynomial");
268
269 /* Fill in the final fit */
270 irplib_vector_fill_line_spectrum (model_spectrum, dispersion, (void *)model);
271 cpl_polynomial_dump (dispersion, stdout);
272
273 /* Evaluate the polynomial for the output */
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);
277
278 /* Check the limits for the QC */
279 minwave = CPL_MAX (cpl_array_get_min (wave_data), minwave);
280 maxwave = CPL_MIN (cpl_array_get_max (wave_data), maxwave);
281
282 /* Fill the WAVE_DATA_SC */
283 cpl_table_new_column_array (wave_data_sc, GRAVI_DATA[reg], CPL_TYPE_DOUBLE, nwave);
284 cpl_table_set_array (wave_data_sc, GRAVI_DATA[reg], 0, wave_data);
285 FREE (cpl_array_delete, wave_data);
286
287 /* Save the test tables */
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);
300
301 CPLCHECK_NUL ("Cannot set data");
302 } /* End loop on regions */
303
304 // cpl_table_save (fit_table, NULL, NULL, "fit_wave.fits", CPL_IO_CREATE);
305 FREE (cpl_table_delete, fit_table);
306
307 /* Delete */
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);
313
314 /* This memory desalocation does not work: */
315 // FREE (cpl_wlcalib_slitmodel_delete, model);
316
318 return wave_data_sc;
319}
320
321/*----------------------------------------------------------------------------*/
351/*----------------------------------------------------------------------------*/
352
353cpl_error_code gravi_opds_correct_closures (cpl_table * phase_table,
354 const char *name)
355{
357 cpl_ensure_code (phase_table, CPL_ERROR_NULL_INPUT);
358 cpl_ensure_code (name, CPL_ERROR_NULL_INPUT);
359
360 int nbase = 6, nclo = 4;
361 cpl_size nrow = cpl_table_get_nrow (phase_table) / nbase;
362
363 /* Get the data */
364 double * phase = cpl_table_get_data_double (phase_table, name);
365 CPLCHECK_MSG ("Cannot get data");
366
367 /* Model and right-hand-side for the lineary system
368 * (unfilled matrix are 0.0) */
369 cpl_matrix * rhs_matrix = cpl_matrix_new (nrow * nclo, 1);
370 cpl_matrix * model_matrix = cpl_matrix_new (nrow * nclo, nbase-1 + nclo);
371
372 for (cpl_size row = 0; row < nrow; row++) {
373 for (int clo = 0; clo < nclo; clo++) {
374 int b0 = GRAVI_CLO_BASE[clo][0];
375 int b1 = GRAVI_CLO_BASE[clo][1];
376 int b2 = GRAVI_CLO_BASE[clo][2];
377
378 /* Fill the rhs with measured closure */
379 double closure = phase[row*nbase+b0] +
380 phase[row*nbase+b1] -
381 phase[row*nbase+b2];
382 cpl_matrix_set (rhs_matrix, row*nclo+clo, 0, closure);
383 CPLCHECK_MSG ("Cannot fill rhs_matrix");
384
385 /* Fill the k (unfilled are 0.0) */
386 cpl_matrix_set (model_matrix, row*nclo+clo, clo, 1.0);
387 CPLCHECK_MSG ("Cannot fill k in model_matrix");
388
389 /* Fill the f in of the model_matrix */
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]);
396 CPLCHECK_MSG ("Cannot fill phase in model_matrix");
397 }
398 } /* End loop on clo and rows */
399
400 /* Solve the system */
401 cpl_matrix * res_matrix = cpl_matrix_solve_normal (model_matrix, rhs_matrix);
402 CPLCHECK_MSG ("Cannot solve system");
403
404 /* Compute residuals */
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);
408
409 /* Correct baseline with (1-f). Last baseline is kept unchanged */
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;
414 }
415
416 /* Dump the residual with their units */
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);
419
420 /* Delete matrix */
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);
425
427 return CPL_ERROR_NONE;
428}
429
430/*----------------------------------------------------------------------------*/
456/*----------------------------------------------------------------------------*/
457
458cpl_vector * gravi_opds_fit_opdmet (cpl_table * ft_table, double lbd_met)
459{
461 cpl_ensure (ft_table, CPL_ERROR_NULL_INPUT, NULL);
462
463 int nbase = 6;
464
465 /* Get the number of acquisitions */
466 cpl_size nrow = cpl_table_get_nrow (ft_table) / nbase;
467
468 /* Get the pointer to data */
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");
472 CPLCHECK_NUL ("Cannot get data");
473
474 /* Number of valid rows */
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);
478
479 cpl_ensure (nrow_valid > 100, CPL_ERROR_ILLEGAL_INPUT, NULL);
480
481 /* Model and right-hand-side for the lineary system
482 * (unfilled matrix are 0.0) */
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);
485
486 for (int base = 0; base < nbase; base++) {
487 cpl_size row_valid = 0;
488
489 for (cpl_size row=0; row<nrow; row++) {
490 if (opd_sc[row*nbase] == 0) continue;
491
492 int idv = row_valid * nbase + base;
493 int id = row * nbase + base;
494
495 /* Fill the OPD metrology */
496 cpl_matrix_set (rhs_matrix, idv, 0, phase_met[id] * lbd_met / CPL_MATH_2PI);
497 CPLCHECK_NUL ("Cannot set OPD");
498
499 /* Fill the model b and c */
500 cpl_matrix_set (model_matrix, idv, 0, 1*opd_sc[id]);
501 cpl_matrix_set (model_matrix, idv, 1, -1*opd_ft[id]);
502 CPLCHECK_NUL ("Cannot set SC or FT");
503
504 /* Fill the model Aij (unfilled matrix are 0.0) */
505 cpl_matrix_set (model_matrix, idv, 2 + base, 1.0);
506 CPLCHECK_NUL ("Cannot set the zero points");
507
508 row_valid++;
509 } /* End loop rows */
510 } /* End loop on base*/
511
512 /* Solve the system */
513 cpl_matrix * res_matrix = cpl_matrix_solve_normal (model_matrix, rhs_matrix);
514
515 /* Compute residuals */
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);
519 //cpl_plot_vector(NULL, NULL, NULL, cpl_vector_wrap(nrow_valid*nbase,cpl_matrix_get_data(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
539 return opd_coeff;
540}
541
542/*----------------------------------------------------------------------------*/
564/*----------------------------------------------------------------------------*/
565
566cpl_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{
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_DRS");
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
670 return guess_table;
671}
672
673/*----------------------------------------------------------------------------*/
696/*----------------------------------------------------------------------------*/
697
698cpl_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 */
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 */
785 return output_table;
786}
787
788/*----------------------------------------------------------------------------*/
820/*----------------------------------------------------------------------------*/
821
822cpl_error_code gravi_wave_compute_opds (gravi_data * spectrum_data,
823 cpl_table * met_table,
824 const cpl_parameterlist * parlist)
825{
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, parlist);
850 /* also apply the TAC reduction */
851 gravi_metrology_tac(met_table, vismet_table, spectrum_header);
852
853 /* Compute the mean LBD_MET for this file */
854 double lbd_met = gravi_pfits_get_met_wavelength_mean (spectrum_header, met_table);
855
856 /*
857 * Compute the phase of SC and FT. For the SC, we use a guess
858 * of the opd modulation, computed from FT and MET, to unwrap.
859 */
860 cpl_msg_info (cpl_func, "Compute OPD of FT from ellipse");
861
862 cpl_table * ft_table;
863 ft_table = gravi_opds_calibration (spectrumft_table, detectorft_table, NULL);
864 gravi_opds_correct_closures (ft_table, "OPD");
865
866 cpl_msg_info (cpl_func, "Compute OPD of SC from ellipse");
867
868 cpl_table * guess_table;
869 guess_table = gravi_opds_compute_guess (spectrumsc_table, ft_table, vismet_table, dit_sc, lbd_met);
870
871 cpl_table * sc_table;
872 sc_table = gravi_opds_calibration (spectrumsc_table, detectorsc_table, guess_table);
873 gravi_opds_correct_closures (sc_table, "OPD");
874 FREE (cpl_table_delete, guess_table);
875
876 CPLCHECK_MSG ("Cannot calibrate phase");
877
878 /*
879 * Interpolate MET and SC at the sampling frequency of FT
880 * Results are saved in the ft_table table
881 */
882 cpl_msg_info (cpl_func, "Fit MET = a.SC - b.FT + c to get absolute modulation");
883
884 gravi_vis_create_met_ft (ft_table, vismet_table);
885 gravi_vis_create_opdsc_ft (ft_table, sc_table, dit_sc);
886
887 CPLCHECK_MSG ("Cannot resample SC or MET at the FT frequency");
888
889 /*
890 * Compute the scaling coefficients of OPDs by fitting:
891 * OPD_MET_ijt = a.OPD_SC_ijt - b.OPD_FT_ijt + c_ij
892 */
893 cpl_vector * coeff_vector = gravi_opds_fit_opdmet (ft_table, lbd_met);
894
895 CPLCHECK_MSG ("Cannot fit opdmet");
896
897 /* Set the CHI2 of the fit in the QC parameters */
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");
902
903 /* Add the OPD COEFF in QC parameters */
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");
908 }
909
910 /* Set a warning */
911 if (cpl_vector_get (coeff_vector, 2) > 1e-7) {
912 gravi_pfits_add_check (spectrum_header,"Residuals of fit MET=a.SC-b.FT+c are high");
913 }
914
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, "*************************************************");
922 }
923
924 CPLCHECK_MSG ("Cannot set QC parameter");
925
926 /*
927 * Correct opd from overall scaling by fitting MET
928 */
929 double coeff_sc = cpl_vector_get (coeff_vector, GRAVI_SC);
930 cpl_table_multiply_scalar (sc_table, "OPD", coeff_sc);
931
932 double coeff_ft = cpl_vector_get (coeff_vector, GRAVI_FT);
933 cpl_table_multiply_scalar (ft_table, "OPD", coeff_ft);
934
935 FREE (cpl_vector_delete, coeff_vector);
936 CPLCHECK_MSG ("Cannot correct OPDs from scaling coefficients");
937
938 /*
939 * Fill the output
940 */
941// if ((det_type == GRAVI_DET_SC || det_type == GRAVI_DET_ALL))
942// {
943 gravi_data_add_table (spectrum_data, NULL, "OPD_SC", sc_table);
944 CPLCHECK_MSG ("Cannot set OPD_SC table");
945// }
946// else
947// FREE (cpl_table_delete, sc_table);
948//
949//
950// if ((det_type == GRAVI_DET_FT || det_type == GRAVI_DET_ALL))
951// {
952 gravi_data_add_table (spectrum_data, NULL, "OPD_FT", ft_table);
953 CPLCHECK_MSG ("Cannot set OPD_FT table");
954// }
955// else
956// FREE (cpl_table_delete, ft_table);
957
958 gravi_data_add_table (spectrum_data, NULL, GRAVI_OI_VIS_MET_EXT, vismet_table);
959 CPLCHECK_MSG ("Cannot set OI_VIS_MET table");
960
962 return CPL_ERROR_NONE;
963}
964
965/*----------------------------------------------------------------------------*/
989/*----------------------------------------------------------------------------*/
990
991cpl_table * gravi_wave_fibre (cpl_table * spectrum_table,
992 cpl_table * detector_table,
993 cpl_table * opd_table)
994{
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);
999
1000 int nbase = 6;
1001 char name[100];
1002
1003 /* Create the output table */
1004 cpl_table * wave_fibre = cpl_table_new (1);
1005
1006 /* Get the number of wavelength, region, polarisation... */
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);
1010
1011 int npol = (n_region > 24 ? 2 : 1);
1012
1013 /*
1014 * Calibration of each polarization and base
1015 */
1016
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);
1021
1022 /* Get the index of the ABCD. */
1023 int iA = gravi_get_region (detector_table, base, 'A', pol);
1024 int iB = gravi_get_region (detector_table, base, 'B', pol);
1025 int iC = gravi_get_region (detector_table, base, 'C', pol);
1026 int iD = gravi_get_region (detector_table, base, 'D', pol);
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 !!!");
1029 }
1030
1031 /* Sign of this baseline */
1032 double phi_sign = gravi_region_get_base_sign (detector_table, base);
1033
1034 /* Get the OPD into various flarous */
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);
1041 }
1042 CPLCHECK_NUL ("Cannot extract the OPD");
1043
1044 /* Create array to fill the wavelenghts */
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);
1049
1050 /*
1051 * Loop on spectral channel
1052 */
1053 for (cpl_size wave = 0; wave < nwave; wave++) {
1054
1055 cpl_vector * vector_T = NULL, * vector_X, * vector_Y;
1056
1057 /* Define the vector_X = C - A */
1058 vector_X = gravi_table_get_vector (spectrum_table, wave, GRAVI_DATA[iC]);
1059 vector_T = gravi_table_get_vector (spectrum_table, wave, GRAVI_DATA[iA]);
1060 cpl_vector_subtract (vector_X, vector_T);
1061 FREE (cpl_vector_delete, vector_T);
1062
1063 /* Define the vector_Y = D - B */
1064 vector_Y = gravi_table_get_vector (spectrum_table, wave, GRAVI_DATA[iD]);
1065 vector_T = gravi_table_get_vector (spectrum_table, wave, GRAVI_DATA[iB]);
1066 cpl_vector_subtract (vector_Y, vector_T);
1067 FREE (cpl_vector_delete, vector_T);
1068
1069 CPLCHECK_NUL ("Cannot compute vector_X or Y");
1070
1071 /* Compute envelope from OPD for this channel */
1072 cpl_vector * envelope_vector = gravi_compute_envelope (opd_vector, wave, nwave);
1073
1074 /* Compute the phase from the ellipse */
1075 cpl_vector * phase;
1076 phase = gravi_ellipse_phase_create (vector_X, vector_Y,
1077 envelope_vector);
1078 FREE (cpl_vector_delete, vector_X);
1079 FREE (cpl_vector_delete, vector_Y);
1080 FREE (cpl_vector_delete, envelope_vector);
1081
1082 /* If the computation of the ellipse fails, we continue with next wave */
1083 if (phase == NULL) {
1084 cpl_msg_warning (cpl_func, "Cannot compute wave for %lld and base %d", wave, base);
1085 continue;
1086 }
1087
1088 /* Multiply by the sign */
1089 cpl_vector_multiply_scalar (phase, phi_sign);
1090
1091 /* Unwrap phase from the OPD */
1092 double lbd_channel = 1.95e-6 + (2.46e-6 - 1.95e-6) / nwave * wave;
1093 gravi_vector_unwrap_with_guess (phase, opd_vector, CPL_MATH_2PI / lbd_channel);
1094
1095 /* Fit the slope of the phase versus OPD gives the
1096 * wavelength of the spectral element */
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);
1100
1101 /* Compute residuals */
1102 cpl_vector * residuals = cpl_vector_new (nwave);
1103 double rechisq;
1104 cpl_vector_fill_polynomial_fit_residual (residuals, phase, NULL, fit_slope, opd_matrix, &rechisq);
1105
1106 if(PLOT_WAVE_PHASE_VS_OPD && (wave == WAVE_TO_PLOT))
1107 {
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);
1113 CPLCHECK_NUL ("Cannot plot OPD versus phase");
1114 }
1115
1116 CPLCHECK_NUL ("Cannot fit the OPD and phase");
1117
1118 /* Get the slope */
1119 const cpl_size pow_slope = 1;
1120 double slope = cpl_polynomial_get_coeff (fit_slope, &pow_slope);
1121
1122 /* Check slope sign. Should be positive */
1123 if (slope < 0.0 && wave == 0) {
1124 cpl_msg_warning (cpl_func, "Negative wavelength!! "
1125 "Report to DRS team");
1126 }
1127
1128 /* Set the wavelength */
1129 cpl_array_set (wavechi2, wave, sqrt(rechisq));
1130 cpl_array_set (wavelength, wave, CPL_MATH_2PI / fabs(slope));
1131
1132 /* Free memory */
1133 cpl_vector_delete (phase);
1134 cpl_vector_delete (residuals);
1135 cpl_polynomial_delete (fit_slope);
1136
1137 } /* End loop on wave */
1138
1139 cpl_matrix_delete (opd_matrix);
1140 cpl_vector_delete (opd_vector);
1141
1142 /* Set the wavelength array in the output table */
1143 sprintf (name, "BASE_%s_%s", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol,npol));
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);
1148
1149 /* Set the chi2 array in the output table */
1150 sprintf (name, "BASE_%s_%s_CHI2", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol,npol));
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);
1154
1155 } /* End loop on base*/
1156 } /* End loop on polar */
1157
1159 return wave_fibre;
1160}
1161
1162/*----------------------------------------------------------------------------*/
1180/*----------------------------------------------------------------------------*/
1181
1182cpl_table * gravi_wave_fit_2d (cpl_table * wavefibre_table,
1183 cpl_table * detector_table,
1184 gravi_data * wave_param,
1185 cpl_size fullstartx,
1186 int spatial_order,
1187 int spectral_order,
1188 double * rms_residuals)
1189{
1191 cpl_ensure (wavefibre_table, CPL_ERROR_NULL_INPUT, NULL);
1192 cpl_ensure (detector_table, CPL_ERROR_NULL_INPUT, NULL);
1193
1194 int nbase = 6;
1195 char name[100];
1196 *rms_residuals = 0;
1197
1198 /* Get numbers */
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");
1202 CPLCHECK_NUL ("Cannot get data");
1203
1204 /* get the wave param */
1205 cpl_propertylist * wave_param_plist = gravi_data_get_plist(wave_param,GRAVI_PRIMARY_HDR_EXT);
1206
1207 /* Odd index, for SC only */
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);
1212 }
1213
1214 CPLCHECK_NUL ("Cannot buid odd_index");
1215
1216 /* Save the 2D coeficient of each polar */
1217 cpl_polynomial ** coef_poly = cpl_calloc (npol, sizeof (cpl_polynomial*));
1218
1219 /* Loop on polarisation */
1220 for (int pol = 0; pol < npol; pol++) {
1221
1222 /* Prepare the 2D coordinates
1223 * and values to fit */
1224 cpl_vector * coord_X = cpl_vector_new (nbase * nwave);
1225 cpl_vector * coord_Y = cpl_vector_new (nbase * nwave);
1226
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);
1230
1231 /*
1232 * Loop on base and wave to get all
1233 * wavelenght and coordinates
1234 */
1235 for (int base = 0; base < nbase; base ++) {
1236
1237 /* Mean position of this baseline */
1238 int iA = gravi_get_region (detector_table, base, 'A', pol);
1239 int iB = gravi_get_region (detector_table, base, 'B', pol);
1240 int iC = gravi_get_region (detector_table, base, 'C', pol);
1241 int iD = gravi_get_region (detector_table, base, 'D', pol);
1242
1243 /* WAVE_FIBRE data */
1244 sprintf (name, "BASE_%s_%s", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol,npol));
1245 cpl_array * wavelength = cpl_table_get_data_array (wavefibre_table, name)[0];
1246
1247 sprintf (name, "BASE_%s_%s_CHI2", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol,npol));
1248 cpl_array * wavechi2 = cpl_table_get_data_array (wavefibre_table, name)[0];
1249
1250 /* Loop on wave */
1251 for (cpl_size wave = 0; wave < nwave; wave++) {
1252 cpl_size pos = base * nwave + wave;
1253
1254 /* Get the values */
1255 int nv = 0;
1256 double wave_value = cpl_array_get (wavelength, wave, &nv);
1257 double chi2_value = cpl_array_get (wavechi2, wave, &nv);
1258
1259 /* The FT accept all channel for the 2D fit */
1260 cpl_vector_set (all_valid, pos, 1);
1261
1262 if (nwave > 1000) {
1263 /* SC HIGH */
1264 if ((chi2_value > M_PI_4 ||
1265 wave_value < gravi_pfits_get_double_default(wave_param_plist, "ESO OIWAVE HIGH LBD MIN", 2.01e-6) ||
1266 wave_value > gravi_pfits_get_double_default(wave_param_plist, "ESO OIWAVE HIGH LBD MAX", 2.43e-6)))
1267 cpl_vector_set (all_valid, pos, 0);
1268 }
1269 else if (nwave > 100) {
1270 /* SC MEDIUM */
1271 if ((chi2_value > M_PI_4 ||
1272 wave_value < gravi_pfits_get_double_default(wave_param_plist, "ESO OIWAVE MED LBD MIN", 2.01e-6) ||
1273 wave_value > gravi_pfits_get_double_default(wave_param_plist, "ESO OIWAVE MED LBD MAX", 2.43e-6)))
1274 cpl_vector_set (all_valid, pos, 0);
1275 }
1276 else if (nwave > GRAVI_LBD_FTSC) {
1277 /* SC LOW */
1278 if ((chi2_value > M_PI_4 ||
1279 //wave_value < 2.01e-6 ||
1280 wave_value < gravi_pfits_get_double_default(wave_param_plist, "ESO OIWAVE LOW LBD MIN", 1.99e-6) ||
1281 wave_value > gravi_pfits_get_double_default(wave_param_plist, "ESO OIWAVE LOW LBD MAX", 2.5e-6) ||
1282 wave == 0 || wave == nwave-1))
1283 cpl_vector_set (all_valid, pos, 0);
1284 }
1285 else if (nwave == GRAVI_LBD_FTSC) {
1286 /* FT */
1287 if ((chi2_value > M_PI_4 ||
1288 //wave_value < 2.01e-6 ||
1289 wave_value < 1.99e-6 ||
1290 wave_value > 2.5e-6))
1291 cpl_vector_set (all_valid, pos, 0);
1292 }
1293
1294 /* Set the wavelength */
1295 cpl_vector_set (all_wavelength, pos, wave_value);
1296 cpl_vector_set (all_wavechi2, pos, chi2_value);
1297
1298 /* Set the X position as the mean of the 4 regions */
1299 cpl_vector_set (coord_X, pos, (double)(iA + iB + iC + iD) / 4.);
1300
1301 /* Set the Y position. Add a shift of 0.15 pixels
1302 * on odd output for SC detector */
1303 cpl_vector_set (coord_Y, pos, wave + cpl_vector_get (odd_index, wave)*0.15);
1304
1305 } /* End loop on wave */
1306 } /* End loop on base */
1307
1308 CPLCHECK_NUL ("Error get all wavelength");
1309
1310 /*
1311 * Reformat the valid point in vector and matrix
1312 */
1313 cpl_size nvalid = cpl_vector_get_sum (all_valid);
1314
1315 cpl_msg_info (cpl_func, "Remove %lld/%lld bad wave",
1316 nbase*nwave - nvalid, nbase * nwave);
1317
1318 cpl_vector * vector = cpl_vector_new (nvalid);
1319 cpl_matrix * matrix = cpl_matrix_new (2, nvalid);
1320
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));
1326 c++;
1327 }
1328
1329 CPLCHECK_NUL ("Error cropping all wavelength");
1330
1331 /*
1332 * Perform a 2D fit with a polynomial model
1333 * between the position and the wavelength = F(y, x)
1334 */
1335 cpl_size deg2d[2] = {2, 3};
1336 if ( (nwave < 20) && (nwave > 8) ) {deg2d[0] = 2; deg2d[1] = 7;} // FIXME Useless line ?
1337 deg2d[0] = spatial_order;
1338 deg2d[1] = spectral_order;
1339
1340 cpl_msg_info (cpl_func, "Fit a 2d polynomial {%lli..%lli} to the wavelengths map", deg2d[0], deg2d[1]);
1341
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;
1345
1346 CPLCHECK_NUL ("Cannot fit 2D");
1347
1348 /*
1349 * Compute residuals
1350 */
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);
1356 CPLCHECK_NUL ("Cannot compute residuals");
1357
1358
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);
1366
1367 } /* End loop on polarisation */
1368
1369
1370 /*
1371 * Create and fill the interpolated WAVE_DATA table
1372 */
1373
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);
1377
1378 for (cpl_size region = 0 ; region < n_region; region ++) {
1379
1380 int pol = gravi_region_get_pol (detector_table, region);
1381 /* Loop on wave to evaluate the 2D polynome */
1382 for (cpl_size wave = 0; wave < nwave; wave ++) {
1383
1384 /* Evaluate */
1385 cpl_vector_set (pos, 0, region);
1386 cpl_vector_set (pos, 1, wave + cpl_vector_get (odd_index, wave)*0.15);
1387
1388 double result = cpl_polynomial_eval (coef_poly[pol], pos);
1389 cpl_array_set (value, wave, result);
1390
1391 }
1392 /* ensure cresent wavelength */
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);
1398 }
1399
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);
1405 }
1406
1407 /* Add column */
1408 char * data_x = GRAVI_DATA[region];
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);
1411
1412 } /* End loop on regions */
1413
1414 /* Delete allocations */
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);
1419
1421 return wavedata_table;
1422}
1423
1424
1425
1426/*----------------------------------------------------------------------------*/
1427/*
1428 * @brief TBD
1429 */
1430/*----------------------------------------------------------------------------*/
1431
1432cpl_table * gravi_wave_fit_individual (cpl_table * wave_individual_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)
1439{
1440
1442
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);
1449
1450 /* Get the number of wavelength, region, polarisation... */
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);
1455 int nbase = 6;
1456 cpl_size nwave_ref=3000;
1457 if (nwave<10) nwave_ref=600;
1458
1459 CPLCHECK_NUL ("Cannot buid odd_index");
1460
1461 /* Get OPD Table */
1462
1463
1464
1465 /* create temporary variables */
1466
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);
1474
1475 /*initialise arrays */
1476
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)
1479 {
1480 // make matrix todo
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);
1483 }
1484
1485 CPLCHECK_NUL ("Cannot initialize arrays for wavelength fit");
1486
1487
1488 for (cpl_size region = 0 ; region < n_region; region ++) {
1489
1490 int base=gravi_region_get_base (detector_table, region);
1491 char * data_x = GRAVI_DATA[region];
1492
1493 cpl_msg_info_overwritable (cpl_func, "Least square fitting of wavelength for region %s", data_x);
1494 // get data_flux_matrix
1495
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));
1500 }
1501 }
1502
1503
1504 for (cpl_size wave_ref = 0; wave_ref < nwave_ref; wave_ref+=1)
1505 {
1506 // make matrix
1507 double wave_value=cpl_array_get(wave_reference_array,wave_ref,NULL);
1508
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)
1513 {
1514 if (nwave <30)
1515 {
1516 coherence_loss=sin(opd*19500)/(opd*19500);
1517 }
1518 }
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);
1521 }
1522
1523 cpl_matrix * coef_vis = cpl_matrix_solve_normal(vis_to_flux_matrix,data_flux_matrix); // coef_vis is 3xnwave
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); //residuals_fit is nrowxnwave
1527
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);
1535 }
1536
1537 FREE (cpl_matrix_delete, coef_vis);
1538 FREE (cpl_matrix_delete, data_flux_fit);
1539 FREE (cpl_matrix_delete, residuals_fit);
1540
1541 CPLCHECK_NUL ("Cannot do Matrix inversion to calculate optimum wavelength");
1542
1543 }
1544
1545 // get minimum chi2 and amplitude signal for each wave
1546 cpl_size wave_ref=1;
1547 cpl_size discarded=1;
1548 for (cpl_size wave = 0; wave < nwave; wave ++ ) {
1549
1550 cpl_matrix * chi2_extract=cpl_matrix_extract_row(residual_matrix,wave);
1551
1552 cpl_matrix_get_minpos(chi2_extract,&discarded, &wave_ref );
1553
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 ));
1556
1557 cpl_array_set_double (wave_individual_array, wave, wave_value);
1558 cpl_array_set_double (weight_individual_array, wave, weight_value);
1559
1560 FREE (cpl_matrix_delete, chi2_extract);
1561
1562 }
1563
1564 /* Add column */
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);
1571 }
1572
1573 CPLCHECK_NUL ("Cannot get individual wavelength for each pixel");
1574
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);
1582
1583 cpl_msg_info (cpl_func, "Now fitting polynomials on wavelength channels");
1584
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);
1588
1589 // set coordinates
1590 for (cpl_size region = 0 ; region < n_region/ npol; region ++)
1591 {
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);
1598 }
1599
1600 for (int pol = 0; pol < npol; pol++) {
1601
1602 cpl_msg_info (cpl_func, "Looping for polyfit now, with pol: %i",(int) pol);
1603
1604 for (cpl_size wave = 0; wave < nwave; wave ++) {
1605
1606 // get the data for a common wave row
1607 for (cpl_size region = 0 ; region < n_region/ npol; region ++) {
1608
1609 // Get the pointers to the table arrays
1610 char * data_x = GRAVI_DATA[region*npol+pol];
1611
1612
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);
1615
1616
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);
1620
1621 }
1622
1623
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);
1626
1627 // Fit a second order polynomial
1628 cpl_matrix * coeff = cpl_matrix_solve_normal(coef_to_wave2, wavelength2); // 5 x 1
1629 cpl_matrix * wavelength_fitted = cpl_matrix_product_create(coef_to_wave, coeff); // n_region / npol x 1
1630
1631 CPLCHECK_NUL ("Cannot do Matrix inversion to calculate optimum polynomial for wavelength");
1632
1633 // Write result to new wave table
1634 for (cpl_size region = 0 ; region < n_region/ npol; region ++) {
1635
1636 char * data_x = GRAVI_DATA[region*npol+pol];
1637 //cpl_msg_info (cpl_func, "Writing: region %i",region);
1638 cpl_array * wave_array = cpl_table_get_data_array (wave_fitted_table, data_x)[0];
1639
1640 cpl_array_set_double(wave_array,wave,cpl_matrix_get(wavelength_fitted,region,0));
1641
1642 }
1643
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);
1648 }
1649
1650 }
1651 FREE (cpl_matrix_delete ,coef_to_wave);
1652 FREE (cpl_matrix_delete ,coef_to_wave_weight);
1653 FREE (cpl_matrix_delete ,wavelength);
1654
1655 CPLCHECK_NUL ("Cannot fit individual wavelength with 3rd order polynomial");
1656
1657 cpl_msg_info (cpl_func, "Correcting for wavelength error");
1658
1659
1660 /* Correct the computed wavelength from the dispersion */
1661 for (cpl_size region = 0 ; region < n_region; region ++)
1662 {
1663 char * data_x = GRAVI_DATA[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 ++ ) {
1667
1668 double result = cpl_array_get (wavelength_fitted, wave, NULL);
1669 double d_met = (result - LAMBDA_MET) / LAMBDA_MET;
1670 cpl_array_set (wavelength_fitted, wave, result * (n0 + n1*d_met + n2*d_met*d_met));
1671
1672 }
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);
1679 }
1680 }
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);
1687 }
1688 }
1689 }
1690
1691 CPLCHECK_NUL ("Error in correcting for dispersion");
1692
1694 return wave_fitted_table;
1695}
1696
1697
1698
1699
1700/*----------------------------------------------------------------------------*/
1707/*----------------------------------------------------------------------------*/
1708
1709cpl_error_code gravi_wave_correct_dispersion (cpl_table * wave_fibre,
1710 double n0, double n1, double n2)
1711{
1713 cpl_ensure_code (wave_fibre, CPL_ERROR_NULL_INPUT);
1714
1715 int nbase = 6;
1716 char name[100];
1717
1718 /* Get the number of polarisation */
1719 int npol = cpl_table_has_column (wave_fibre, "BASE_12_P") ? 2 : 1;
1720
1721 /* Loop on columns in the table */
1722 for (int pol = 0; pol < npol; pol ++) {
1723 for (int base = 0; base < nbase; base ++) {
1724
1725 /* Get data of this region */
1726 sprintf (name, "BASE_%s_%s", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol, npol));
1727 cpl_array * wavelength = cpl_table_get_data_array (wave_fibre, name)[0];
1728 CPLCHECK_MSG ("Cannot get data");
1729
1730 /* Loop on wave */
1731 cpl_size nwave = cpl_array_get_size (wavelength);
1732 for (cpl_size wave = 0 ; wave < nwave ; wave ++ ) {
1733
1734 /* Correct the computed wavelength from the dispersion */
1735 double result = cpl_array_get (wavelength, wave, NULL);
1736 double d_met = (result - LAMBDA_MET) / LAMBDA_MET;
1737 cpl_array_set (wavelength, wave, result * (n0 + n1*d_met + n2*d_met*d_met));
1738
1739 }
1740 } /* End loop on base */
1741 } /* End loop on pol */
1742
1744 return CPL_ERROR_NONE;
1745}
1746
1747/*----------------------------------------------------------------------------*/
1753/*----------------------------------------------------------------------------*/
1754
1755cpl_error_code gravi_wave_correct_color (gravi_data * vis_data)
1756{
1758 cpl_ensure_code (vis_data, CPL_ERROR_NULL_INPUT);
1759
1760 cpl_propertylist * primary_header = gravi_data_get_header (vis_data);
1761 cpl_propertylist * oiwave_header = NULL;
1762 cpl_table * oiwave_table = NULL;
1763
1764 /* For each type of data SC / FT */
1765 int ntype_data = 2;
1766 for (int type_data = 0; type_data < ntype_data ; type_data ++) {
1767
1768 /* Loop on polarisation */
1769 int npol = gravi_pfits_get_pola_num (primary_header, type_data);
1770 for (int pol = 0 ; pol<npol ; pol++) {
1771 oiwave_table = cpl_table_duplicate ( gravi_data_get_oi_wave ( vis_data, type_data, pol, npol ) );
1772 oiwave_header = cpl_propertylist_duplicate ( gravi_data_get_oi_wave_plist ( vis_data, type_data, pol, npol ) );
1773
1774 /* here you can do what you want on this duplicated oi_wave_table
1775 * to get OI_FLUX table :
1776 * cpl_table * oiflux_table = gravi_data_get_oi_flux(vis_data, type_data, pol, npol)
1777 * */
1778
1779
1780 /* save the new extension with name OI_WAVELENGTH_CORR */
1781 gravi_data_add_table(vis_data, oiwave_header, "OI_WAVELENGTH_CORR", oiwave_table);
1782
1783 CPLCHECK_MSG("Cannot apply color wave correction");
1784 }
1785 /* End loop on polarisation */
1786 }
1787 /* End loop on data_type */
1788
1790 return CPL_ERROR_NONE;
1791}
1792
1793
1794/*----------------------------------------------------------------------------*/
1806/*----------------------------------------------------------------------------*/
1807
1808cpl_imagelist * gravi_wave_test_image (cpl_table * wavedata_table,
1809 cpl_table * wavefibre_table,
1810 cpl_table * profile_table,
1811 cpl_table * detector_table)
1812{
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);
1818
1819 int nv = 0;
1820 char name[100];
1821
1822 cpl_size n_region = cpl_table_get_nrow (detector_table);
1823 int npol = (n_region > 24 ? 2 : 1);
1824
1825 CPLCHECK_NUL ("Cannot get data");
1826
1827 /* Init the output images */
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);
1830
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);
1833
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);
1836
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);
1839
1840 CPLCHECK_NUL ("Cannot prepare output image");
1841
1842 /* Loop on regions */
1843 for (cpl_size region = 0 ; region < n_region; region ++) {
1844
1845 /* Load the profile of this region as an image */
1846 cpl_imagelist * profile_imglist = gravi_imagelist_wrap_column (profile_table, GRAVI_DATA[region]);
1847 cpl_image * profile_image = cpl_imagelist_get (profile_imglist, 0);
1848
1849 CPLCHECK_NUL ("Cannot get data");
1850
1851 /* Sum all profils of all regions */
1852 cpl_image_add (profilesum_image, profile_image);
1853
1854 /*
1855 * Define an image of each region with its WAVE_DATA
1856 */
1857 const cpl_array * wavelength;
1858 wavelength = cpl_table_get_array (wavedata_table, GRAVI_DATA[region], 0);
1859 CPLCHECK_NUL ("Cannot get data");
1860
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));
1866 }
1867 }
1868 CPLCHECK_NUL ("Cannot make image of wave_map");
1869
1870 /*
1871 * Define an image of each region with its WAVE_FIBER
1872 */
1873 int base = gravi_region_get_base (detector_table, region);
1874 int pol = gravi_region_get_pol (detector_table, region);
1875 sprintf (name, "BASE_%s_%s", GRAVI_BASE_NAME[base], GRAVI_POLAR(pol, npol));
1876 wavelength = cpl_table_get_array (wavefibre_table, name, 0);
1877 CPLCHECK_NUL ("Cannot get data");
1878
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));
1884 }
1885 }
1886
1887 gravi_imagelist_unwrap_images (profile_imglist);
1888 CPLCHECK_NUL ("Cannot make image of wave_fibre");
1889 } /* End loop on regions */
1890
1891 /* Create an imagelist */
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);
1896
1898 return testwave_imglist;
1899}
1900
1901/*----------------------------------------------------------------------------*/
1911/*----------------------------------------------------------------------------*/
1912
1913cpl_error_code gravi_wave_qc (gravi_data * wave_map, gravi_data * profile_map)
1914{
1916 cpl_ensure_code (wave_map, CPL_ERROR_NULL_INPUT);
1917
1918 char name[100];
1919
1920 cpl_propertylist * wave_header = gravi_data_get_header (wave_map);
1921
1922 /* Loop on extensions (thus SC/FT) */
1923 for (int type_data = 0; type_data < 2; type_data ++ ) {
1924
1925 /* Get WAVE_FIBRE and IMAGING_DETECTOR tables */
1926 cpl_table * detector_table = gravi_data_get_imaging_detector (wave_map, type_data);
1927 cpl_table * wave_data = gravi_data_get_wave_data (wave_map, type_data);
1928 cpl_size n_region = cpl_table_get_nrow (detector_table);
1929
1930 /* Get the full STARTX */
1931 cpl_propertylist * plist = gravi_data_get_wave_data_plist (wave_map, type_data);
1932 int fullstartx = gravi_pfits_get_fullstartx (plist);
1933
1934 CPLCHECK_MSG ("Cannot get data");
1935
1936 /*
1937 * Compute the min and max wave over all regions
1938 */
1939 double minwave = -1e10;
1940 double maxwave = 1e10;
1941
1942 for (int region = 0; region < n_region; region++) {
1943 const cpl_array * wavelength = cpl_table_get_array (wave_data, GRAVI_DATA[region], 0);
1944
1945 minwave = CPL_MAX (minwave, cpl_array_get_min (wavelength));
1946 maxwave = CPL_MIN (maxwave, cpl_array_get_max (wavelength));
1947 } /* End loop on regions */
1948
1949 cpl_msg_info (cpl_func,"%s = %g [m]", QC_MINWAVE(type_data), minwave);
1950 cpl_msg_info (cpl_func,"%s = %g [m]", QC_MAXWAVE(type_data), maxwave);
1951 cpl_propertylist_update_double (wave_header, QC_MINWAVE(type_data), minwave);
1952 cpl_propertylist_update_double (wave_header, QC_MAXWAVE(type_data), maxwave);
1953 cpl_propertylist_update_double (wave_header, QC_MINWAVE_UM(type_data), minwave * 1e6);
1954 cpl_propertylist_update_double (wave_header, QC_MAXWAVE_UM(type_data), maxwave * 1e6);
1955
1956 CPLCHECK_MSG ("Cannot compute minwave or maxwave");
1957
1958
1959 /*
1960 * Compute the pixel position of QC specified argon wavelength.
1961 * Having two lines allow to check the position and dispersion
1962 */
1963 int qc_reg[3] = {0,23,47};
1964 double qc_wave[2] = {2.099184e-6, 2.313952e-6};
1965
1966 /* Loop on regions */
1967 for (int reg = 0; reg < (n_region > 24 ? 3 : 2); reg++) {
1968 cpl_size region = qc_reg[reg];
1969
1970 const cpl_array * wavelength = cpl_table_get_array (wave_data, GRAVI_DATA[region], 0);
1971 cpl_size nwave = cpl_array_get_size (wavelength);
1972 const double * wave_tab = cpl_array_get_data_double_const (wavelength);
1973 CPLCHECK_MSG ("Cannot get wavelength data");
1974
1975 /* Loop on the two argon lines */
1976 for (int iqc = 0 ; iqc < 2 ; iqc++) {
1977 cpl_size l2 = 0;
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 ++;}
1980
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]);
1983 continue;
1984 }
1985
1986 /* Position on full detector */
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]);
1989
1990 sprintf (name, "ESO QC REFWAVE%i", iqc+1);
1991 cpl_propertylist_update_double (wave_header, name, qc_wave[iqc]);
1992 cpl_propertylist_set_comment (wave_header, name, "[m] value of ref wave");
1993
1994 sprintf (name, "ESO QC REFPOS%i %s%lli", iqc+1, GRAVI_TYPE(type_data),region+1);
1995 cpl_propertylist_update_double (wave_header, name, qc_pos);
1996 cpl_propertylist_set_comment (wave_header, name, "[pix] position of ref wave");
1997
1998 cpl_msg_info (cpl_func, "%s = %f [pix] for %e [m]", name, qc_pos, qc_wave[iqc]);
1999
2000 CPLCHECK_MSG ("Cannot set QC");
2001 } /* End loop on the 2 argon lines */
2002 } /* End loop on 2 or 3 regions */
2003 } /* End loop on SC / FT */
2004
2005
2006 /*
2007 * Create the test image for SC (only used for debug)
2008 */
2009 cpl_table * wavedata_table = gravi_data_get_wave_data (wave_map, GRAVI_SC);
2010 cpl_table * wavefibre_table = gravi_data_get_wave_fibre (wave_map, GRAVI_SC);
2011 cpl_table * profile_table = gravi_data_get_table (profile_map, GRAVI_PROFILE_DATA_EXT);
2012 cpl_table * detector_table = gravi_data_get_imaging_detector (wave_map, GRAVI_SC);
2013
2014 cpl_imagelist * testwave_imglist;
2015 testwave_imglist = gravi_wave_test_image (wavedata_table, wavefibre_table,
2016 profile_table, detector_table);
2017 CPLCHECK_MSG ("Cannot compute TEST_WAVE");
2018
2019 /* Set TEST_WAVE in output */
2020 cpl_propertylist * plist = gravi_data_get_plist (profile_map, GRAVI_PROFILE_DATA_EXT);
2021 gravi_data_add_cube (wave_map, cpl_propertylist_duplicate (plist),
2022 "TEST_WAVE", testwave_imglist);
2023
2024 CPLCHECK_MSG ("Cannot set data");
2025
2027 return CPL_ERROR_NONE;
2028}
2029
2030/*----------------------------------------------------------------------------*/
2058/*----------------------------------------------------------------------------*/
2059
2060cpl_error_code gravi_compute_wave (gravi_data * wave_map,
2061 gravi_data * spectrum_data,
2062 int type_data, const cpl_parameterlist * parlist,
2063 gravi_data * wave_param)
2064{
2066 cpl_ensure_code (wave_map, CPL_ERROR_NULL_INPUT);
2067 cpl_ensure_code (spectrum_data, CPL_ERROR_NULL_INPUT);
2068 cpl_ensure_code (type_data == GRAVI_SC || type_data == GRAVI_FT,
2069 CPL_ERROR_ILLEGAL_INPUT);
2070
2071 /* Get headers */
2072 cpl_propertylist * raw_header = gravi_data_get_header (spectrum_data);
2073 cpl_propertylist * wave_header = gravi_data_get_header (wave_map);
2074
2075 /* Dump full header of wave data */
2076 if (type_data == GRAVI_FT) {
2077 cpl_propertylist_append (wave_header, raw_header);
2078 }
2079
2080
2081 /* Copy IMAGING_DETECTOR in output WAVE */
2082 const char * extname = (type_data == GRAVI_SC) ? GRAVI_IMAGING_DETECTOR_SC_EXT : GRAVI_IMAGING_DETECTOR_FT_EXT;
2083 gravi_data_copy_ext (wave_map, spectrum_data, extname);
2084
2085 /*
2086 * Compute WAVE_FIBRE map
2087 */
2088 cpl_msg_info (cpl_func, "Compute WAVE_FIBRE for %s", GRAVI_TYPE(type_data));
2089
2090 cpl_table * spectrum_table = gravi_data_get_spectrum_data (spectrum_data, type_data);
2091 cpl_table * detector_table = gravi_data_get_imaging_detector (spectrum_data, type_data);
2092 cpl_table * opd_table = gravi_data_get_table (spectrum_data, type_data==GRAVI_SC ? "OPD_SC" : "OPD_FT");
2093
2094 cpl_table * wavefibre_table;
2095 wavefibre_table = gravi_wave_fibre (spectrum_table, detector_table, opd_table);
2096 CPLCHECK_MSG ("Cannot compute the WAVE_FIBRE");
2097
2098 /*
2099 * Correct WAVE_FIBRE from the dispersion model
2100 */
2101 cpl_msg_info (cpl_func, "Correct dispersion in WAVE_FIBRE of %s", GRAVI_TYPE(type_data));
2102
2103 /* Hardcoded values to correct for the fiber dispersion */
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);
2106
2107 cpl_propertylist_update_string (wave_header, "ESO QC WAVE_CORR", "lbd*(N0+N1*(lbd-lbd0)/lbd0+N2*(lbd-lbd0)^2/lbd0^2)");
2108 cpl_propertylist_update_double (wave_header, "ESO QC WAVE_CORR N0", n0);
2109 cpl_propertylist_update_double (wave_header, "ESO QC WAVE_CORR N1", n1);
2110 cpl_propertylist_update_double (wave_header, "ESO QC WAVE_CORR N2", n2);
2111 CPLCHECK_MSG ("Cannot set Keywords");
2112
2113 gravi_wave_correct_dispersion (wavefibre_table, n0, n1, n2);
2114 CPLCHECK_MSG ("Cannot correct dispersion in wave");
2115
2116 /*
2117 * Fit the 2d dispersion WAVE_FIBRE into WAVE_DATA
2118 */
2119 cpl_msg_info (cpl_func,"Fit WAVE_DATA map for %s", GRAVI_TYPE(type_data));
2120
2121 /* Get the fullstartx */
2122 cpl_propertylist * spectrum_plist = gravi_data_get_spectrum_data_plist (spectrum_data, type_data);
2123 int fullstartx = gravi_pfits_get_fullstartx (spectrum_plist);
2124
2125 /* Interpolate table 2D */
2126 cpl_table * wavedata_table;
2127 int spatial_order=2; // default spatial order
2128 int spectral_order=3; // default spectral order
2129 double rms_residuals;
2130 if (type_data == GRAVI_FT && gravi_param_get_bool(parlist, "gravity.calib.force-wave-ft-equal")) {
2131 spatial_order = 0;
2132 cpl_msg_info (cpl_func, "Option force-waveFT-equal applied");
2133 }
2134 // Keep default value
2135 // if (type_data == GRAVI_SC) {
2136 // spectral_order = gravi_param_get_int(parlist, "gravity.calib.wave-spectral-order");
2137 // cpl_msg_info (cpl_func, "Option set_spectral order to %d", spectral_order);
2138 // }
2139
2140 wavedata_table = gravi_wave_fit_2d (wavefibre_table,
2141 detector_table,
2142 wave_param,
2143 fullstartx, spatial_order, spectral_order, &rms_residuals);
2144 cpl_propertylist_update_double (wave_header, QC_RMS_RESIDUALS(type_data), rms_residuals);
2145 cpl_propertylist_update_double (wave_header, QC_RMS_RESIDUALS_UM(type_data), rms_residuals * 1e6);
2146
2147
2148 double rms_fit=cpl_propertylist_get_double (raw_header, QC_PHASECHI2);
2149 cpl_propertylist_update_double (wave_header, QC_CHI2WAVE(type_data),rms_fit*1e9);
2150 cpl_propertylist_set_comment (wave_header, QC_CHI2WAVE(type_data),"[nm]rms a.SC-b.FT+c=MET");
2151
2152 CPLCHECK_MSG ("Cannot fit 2d data");
2153
2154 /*
2155 * New Wavelength interpolation made by sylvestre on January 30 2018
2156 */
2157
2158 if (type_data == GRAVI_SC && !strcmp (gravi_pfits_get_spec_res(raw_header), "LOW"))
2159 {
2160 cpl_msg_info (cpl_func, "Additional Wavelength Fit");
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);
2164
2165 gravi_wave_fit_individual (wave_individual_table,
2166 weight_individual_table,
2167 wave_fitted_table,
2168 opd_table,
2169 spectrum_table,
2170 detector_table,
2171 n0,n1,n2);
2172
2173 cpl_msg_info (cpl_func,"Add tables in wave_map");
2174
2175 gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2176 "WAVE_INDIV_SC", wave_individual_table);
2177 gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2178 "WAVE_WEIGHT_SC", weight_individual_table);
2179 gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2180 "WAVE_FITTED_SC", wavedata_table);
2181
2182 gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2183 GRAVI_WAVE_FIBRE_EXT(type_data), wavefibre_table);
2184 gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2185 GRAVI_WAVE_DATA_EXT(type_data), wave_fitted_table);
2186 } else {
2187 /*
2188 * Add the WAVE_FIBRE and WAVE_DATA table in the wave_map
2189 */
2190 cpl_msg_info (cpl_func,"Add WAVE_FIBRE and WAVE_DATA in wave_map");
2191
2192 gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2193 GRAVI_WAVE_FIBRE_EXT(type_data), wavefibre_table);
2194
2195 gravi_data_add_table (wave_map, cpl_propertylist_duplicate (spectrum_plist),
2196 GRAVI_WAVE_DATA_EXT(type_data), wavedata_table);
2197
2198 }
2199
2200
2201
2202
2203 CPLCHECK_MSG ("Cannot set data");
2204
2206 return CPL_ERROR_NONE;
2207}
2208
2209
#define GRAVI_PRIMARY_HDR_EXT
Definition: gravi-test.h:212
typedefCPL_BEGIN_DECLS struct _gravi_data_ gravi_data
Definition: gravi_data.h:38
#define gravi_data_get_spectrum_data(data, type)
Definition: gravi_data.h:62
#define gravi_data_get_wave_data_plist(data, type)
Definition: gravi_data.h:53
#define gravi_data_get_wave_data(data, type)
Definition: gravi_data.h:52
#define gravi_data_get_spectrum_data_plist(data, type)
Definition: gravi_data.h:63
#define gravi_data_get_header(data)
Definition: gravi_data.h:74
#define gravi_data_get_oi_wave(data, type, pol, npol)
Definition: gravi_data.h:44
#define gravi_data_get_imaging_detector(data, type)
Definition: gravi_data.h:59
#define gravi_data_get_wave_fibre(data, type)
Definition: gravi_data.h:50
#define gravi_data_get_oi_wave_plist(data, type, pol, npol)
Definition: gravi_data.h:68
const cpl_size ntel
#define USE_LINEAR_ELLIPSE
Definition: gravi_ellipse.h:35
gravi_data * wave_data
Definition: gravi_old.c:1994
cpl_propertylist * plist
Definition: gravi_old.c:2000
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"))
cpl_table * wave_data_sc
Definition: gravi_old.c:1997
gravi_vis_create_met_ft(phase_ft, vis_met)
#define GRAVI_OI_VIS_MET_EXT
Definition: gravi_pfits.h:88
#define GRAVI_IMAGING_DETECTOR_FT_EXT
Definition: gravi_pfits.h:82
#define GRAVI_SC
Definition: gravi_pfits.h:165
#define QC_MAXWAVE(type)
Definition: gravi_pfits.h:110
#define QC_MINWAVE_UM(type)
Definition: gravi_pfits.h:111
#define GRAVI_PROFILE_DATA_EXT
Definition: gravi_pfits.h:79
#define GRAVI_WAVE_FIBRE_EXT(type)
Definition: gravi_pfits.h:70
#define QC_MAXWAVE_UM(type)
Definition: gravi_pfits.h:112
#define OPD_COEFF_SIGN(type)
Definition: gravi_pfits.h:115
#define QC_CHI2WAVE(type)
Definition: gravi_pfits.h:108
#define GRAVI_WAVE_DATA_EXT(type)
Definition: gravi_pfits.h:67
#define QC_RMS_RESIDUALS_UM(type)
Definition: gravi_pfits.h:114
#define GRAVI_TYPE(type)
Definition: gravi_pfits.h:167
#define QC_MINWAVE(type)
Definition: gravi_pfits.h:109
#define LAMBDA_MET
Definition: gravi_pfits.h:103
#define GRAVI_IMAGING_DETECTOR_SC_EXT
Definition: gravi_pfits.h:81
#define GRAVI_FT
Definition: gravi_pfits.h:166
#define QC_PHASECHI2
Definition: gravi_pfits.h:107
#define QC_RMS_RESIDUALS(type)
Definition: gravi_pfits.h:113
#define GRAVI_POLAR(pol, npol)
Definition: gravi_utils.h:143
#define gravi_spectrum_get_npol(table)
Definition: gravi_utils.h:92
#define gravi_msg_function_exit(flag)
Definition: gravi_utils.h:85
#define FREE(function, variable)
Definition: gravi_utils.h:69
#define CPLCHECK_NUL(msg)
Definition: gravi_utils.h:48
#define gravi_msg_function_start(flag)
Definition: gravi_utils.h:84
#define CPLCHECK_MSG(msg)
Definition: gravi_utils.h:45
#define GRAVI_LBD_FTSC
Definition: gravi_utils.h:115
#define FREELOOP(function, variable, n)
Definition: gravi_utils.h:72
double gravi_table_get_column_mean(cpl_table *table, const char *name, int base, int nbase)
Definition: gravi_cpl.c:275
cpl_error_code gravi_table_new_column(cpl_table *table, const char *name, const char *unit, cpl_type type)
Definition: gravi_cpl.c:1588
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:2651
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:1689
cpl_error_code gravi_imagelist_unwrap_images(cpl_imagelist *imglist)
Unwrap an imagelist an all its images.
Definition: gravi_cpl.c:1659
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:2167
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:2748
cpl_error_code gravi_table_new_column_array(cpl_table *table, const char *name, const char *unit, cpl_type type, cpl_size size)
Definition: gravi_cpl.c:1610
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:2353
cpl_propertylist * gravi_data_get_plist(gravi_data *self, const char *extname)
Get the propertylist from EXTNAME.
Definition: gravi_data.c:2049
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:2289
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:2096
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:1690
int gravi_param_get_bool(const cpl_parameterlist *parlist, const char *name)
Definition: gravi_dfs.c:1537
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)
Definition: gravi_pfits.c:263
int gravi_pfits_get_fullstartx(const cpl_propertylist *plist)
Definition: gravi_pfits.c:76
cpl_error_code gravi_pfits_add_check(cpl_propertylist *header, const char *msg)
Add a QC.CHECK keyword to the header.
Definition: gravi_pfits.c:1643
const char * gravi_pfits_get_spec_res(const cpl_propertylist *plist)
Definition: gravi_pfits.c:162
double gravi_pfits_get_dit_sc(const cpl_propertylist *plist)
Definition: gravi_pfits.c:664
double gravi_pfits_get_met_wavelength_mean(const cpl_propertylist *plist, cpl_table *met_table)
Definition: gravi_pfits.c:312
double gravi_pfits_get_double_default(const cpl_propertylist *plist, const char *name, double def)
Definition: gravi_pfits.c:1587
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.
Definition: gravi_utils.c:445
int GRAVI_BASE_TEL[6][2]
Definition: gravi_utils.c:56
int gravi_region_get_base(cpl_table *imaging_detector, int region)
Return the base of a region.
Definition: gravi_utils.c:361
char GRAVI_BASE_NAME[6][3]
Definition: gravi_utils.c:57
cpl_size gravi_spectrum_get_nwave(const cpl_table *table)
Definition: gravi_utils.c:1013
char GRAVI_DATA[50][7]
Definition: gravi_utils.c:71
cpl_size gravi_spectrum_get_nregion(const cpl_table *table)
Definition: gravi_utils.c:1018
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
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
int GRAVI_CLO_BASE[4][3]
Definition: gravi_utils.c:68
#define WAVE_TO_PLOT
Definition: gravi_wave.c:79
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_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:458
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:991
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...
Definition: gravi_wave.c:1182
#define PLOT_WAVE_PHASE_VS_OPD
Definition: gravi_wave.c:78
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.
Definition: gravi_wave.c:822
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:353
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.
Definition: gravi_wave.c:2060
cpl_table * gravi_compute_argon_wave(cpl_table *spectrum_table)
Compute a WAVE calibration from the ARGON data (SC only)
Definition: gravi_wave.c:148
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)
Definition: gravi_wave.c:1432
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:1913
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:1709
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:1755
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:1808
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