GRAVI Pipeline Reference Manual  1.2.3
gravi_preproc.c
1 /* $Id: gravi_preproc.c,v 1.12 2011/04/31 06: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 
33 /*
34  * History
35  * 11/01/2019 Fix Warnings , unused parameter : profile_header
36  * brackets of if statements
37  */
38 
39 /*-----------------------------------------------------------------------------
40  Includes
41  -----------------------------------------------------------------------------*/
42 
43 #ifdef HAVE_CONFIG_H
44 #include <config.h>
45 #endif
46 
47 #include <cpl.h>
48 #include <stdio.h>
49 #include <string.h>
50 #include <math.h>
51 #include <time.h>
52 #include <complex.h>
53 
54 #include "gravi_data.h"
55 #include "gravi_calib.h"
56 #include "gravi_dfs.h"
57 #include "gravi_pfits.h"
58 #include "gravi_cpl.h"
59 
60 #include "gravi_utils.h"
61 #include "gravi_preproc.h"
62 
63 /*-----------------------------------------------------------------------------
64  Private prototypes
65  -----------------------------------------------------------------------------*/
66 cpl_table * gravi_table_ft_format (cpl_table * table_ft, cpl_table * sky_table_std,
67  cpl_table * sky_table_avg, cpl_table * badft_table,
68  int n_region, double gain);
69 
70 cpl_table * gravi_imglist_sc_collapse (cpl_table * profile_table,
71  cpl_imagelist * raw_imglist,
72  cpl_imagelist * rawVar_imglist,
73  cpl_size startx);
74 
75 cpl_error_code gravi_interpolate_spectrum_table (cpl_table * spectrum_table,
76  cpl_table * wave_table,
77  cpl_table ** oiwave_tables,
78  cpl_table * detector_table,
79  cpl_table * specflat_table,
80  const cpl_parameterlist * parlist,
81  int type_data);
82 
83 /*-----------------------------------------------------------------------------
84  Function code
85  -----------------------------------------------------------------------------*/
86 
87 /*---------------------------------------------------------------------------*/
97 /*---------------------------------------------------------------------------*/
98 
99 int gravi_pixel_is_good (cpl_image * bad_img, int x, int y)
100 {
101  int bad;
102  int nv;
103 
104  bad = cpl_image_get (bad_img, x, y, &nv);
105  if ((bad & BADPIX_DARK) || (bad & BADPIX_RMS) || (bad & BADPIX_FLAT))
106  return 0;
107  else
108  return 1;
109 }
110 
111 
112 /*---------------------------------------------------------------------------*/
123 /*---------------------------------------------------------------------------*/
124 
125 cpl_error_code gravi_remove_badpixel_sc (cpl_imagelist * imglist_sc, cpl_image * bad_img)
126 {
127  gravi_msg_function_start(0);
128  cpl_ensure_code (imglist_sc, CPL_ERROR_NULL_INPUT);
129  cpl_ensure_code (bad_img, CPL_ERROR_NULL_INPUT);
130 
131  cpl_image * img;
132  int middel = 3, size = 7, nv, comp, badpix_comp;
133  cpl_vector * x = NULL;
134 
135  cpl_size nrow = cpl_imagelist_get_size (imglist_sc);
136  cpl_size nx = cpl_image_get_size_x (bad_img);
137  cpl_size ny = cpl_image_get_size_y (bad_img);
138 
139  /* Check consistency of BADPIX and DATA */
140  if ((nx != cpl_image_get_size_x (cpl_imagelist_get (imglist_sc, 0))) ||
141  (ny != cpl_image_get_size_y (cpl_imagelist_get (imglist_sc, 0))) ){
142  return cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
143  "The image lists have not the same size");
144  }
145 
146  /* Compute the amount of badpixels */
147  cpl_size nbad = 0;
148  for (cpl_size k = 0; k < ny; k++)
149  for (cpl_size i = 0; i < nx; i++)
150  if (cpl_image_get (bad_img, i+1, k+1, &nv) != 0) nbad++;
151 
152  /* Check the fraction of badpixels */
153  if (nbad > 0.25 * nx*ny) {
154  return cpl_error_set_message (cpl_func, CPL_ERROR_ILLEGAL_INPUT,
155  "Too many bad pixels (more than 25 percent)");
156  }
157 
158  /* Loop on the pixels of the image */
159  for (cpl_size k = 0; k < ny; k++){
160  for (cpl_size i = 0; i < nx; i++){
161 
162  /* This is a bad pixel */
163  if (cpl_image_get (bad_img, i+1, k+1, &nv) != 0) {
164 
165  /* Loop on the images of the RAW imagelist */
166  for (cpl_size row = 0; row < nrow; row ++){
167 
168  img = cpl_imagelist_get (imglist_sc, row);
169  x = cpl_vector_new (size-1);
170 
171  if ((i - middel) < 0 ){
172  comp = 0;
173  badpix_comp = 1;
174  for (cpl_size j = 0; j < size; j++){
175  if (j != i) {
176  while (cpl_image_get (bad_img, i + badpix_comp, k + 1, &nv) != 0){
177  badpix_comp ++;
178  }
179  cpl_vector_set(x, comp, cpl_image_get (img, i + badpix_comp, k + 1, &nv));
180  comp ++ ;
181  }
182  }
183  CPLCHECK_MSG("Fail 1");
184  }
185  else if ((i + middel) >= nx){
186  comp = 0;
187  badpix_comp = 1;
188  for (cpl_size j = 0; (j < size); j++){
189  if (j != (nx - i)){
190  while (cpl_image_get (bad_img, nx - (j + 1) + badpix_comp, k + 1, &nv) != 0){
191 
192  badpix_comp --;
193  }
194  cpl_vector_set(x, comp, cpl_image_get (img, nx - (j + 1) + badpix_comp, k + 1 , &nv));
195  comp ++;
196  }
197  }
198  CPLCHECK_MSG("Fail 2");
199  }
200  else{
201  comp = 0;
202  badpix_comp = 1;
203  int test1 = 0, test2 = 0;
204  for (cpl_size j = 0; j < size; j++){
205  if (j != middel){
206  if (i + j - middel + badpix_comp <= 1){
207  while (cpl_image_get (bad_img, i + j - middel + badpix_comp, k+1, &nv) != 0){
208  badpix_comp ++;
209  test1 = 1;
210  }
211  }
212  else if (i + j - middel + badpix_comp >= nx){
213  while (cpl_image_get (bad_img, i + j - middel + badpix_comp, k+1, &nv) != 0){
214  badpix_comp --;
215  test2 = 1;
216  }
217  }
218  else {
219  if ((test1 == 0) && (test2 == 0)){
220  if (j > middel)
221  badpix_comp ++;
222  else
223  badpix_comp --;
224  }
225  else if (test1)
226  badpix_comp ++;
227  else
228  badpix_comp --;
229  }
230  cpl_vector_set(x, comp, cpl_image_get (img, i + j - middel + badpix_comp, k+1, &nv));
231  comp ++;
232  }
233  }
234  CPLCHECK_MSG("Fail 3");
235  }
236 
237  double mean = cpl_vector_get_median (x);
238  cpl_vector_delete (x);
239  cpl_image_set (img, i + 1, k + 1, mean);
240 
241  CPLCHECK_MSG("Fail 4");
242 
243  } /* End loop on row */
244  } /* End case this is a bad pixel */
245  } /* End loop k */
246  } /* End loop k */
247 
248  gravi_msg_function_exit(0);
249  return CPL_ERROR_NONE;
250 }
251 
252 /*---------------------------------------------------------------------------*/
269 /*---------------------------------------------------------------------------*/
270 
271 cpl_table * gravi_table_ft_format (cpl_table * pix_table,
272  cpl_table * skystd_table,
273  cpl_table * skyavg_table,
274  cpl_table * badft_table,
275  int n_region, double gain)
276 {
277  /* Verbose */
278  gravi_msg_function_start(0);
279  cpl_ensure (pix_table, CPL_ERROR_NULL_INPUT, NULL);
280 
281  /* Get the number of frames */
282  cpl_size nrow = cpl_table_get_nrow (pix_table);
283  cpl_table * output_table = cpl_table_new (nrow);
284 
285  /* Get PIX of data */
286  cpl_array ** arr_data = cpl_table_get_data_array (pix_table, "PIX");
287  CPLCHECK_NUL ("Cannot get data");
288 
289  /* Get PIX of badft */
290  cpl_array ** arr_badft = cpl_table_get_data_array (badft_table, "PIX");
291  CPLCHECK_NUL ("Cannot get data");
292 
293  /* Get the dimensions for the descramling */
294  cpl_size npix = cpl_table_get_column_depth (pix_table, "PIX");
295  cpl_size sizex = cpl_table_get_column_dimension(pix_table, "PIX", 0);
296  cpl_size sizey = cpl_table_get_column_dimension(pix_table, "PIX", 1);
297 
298  int npol = n_region > 24 ? 2 : 1;
299  cpl_size ny = sizey / npol;
300 
301  /* Keep the ny_out to 5 for now to keep the same size.
302  * could be ny if the pipeline accept 6 pixels spectra */
303  cpl_size ny_out = 5;
304 
305  /* Number of RAW pixels */
306  cpl_size n_output = 24;
307  cpl_size nx = sizex / n_output;
308  cpl_msg_info (cpl_func, "Descramble %lld x %lld as %lld outputs x %lld pixels x %lld channels x %i pol",
309  sizex, sizey, n_output, nx, ny, npol);
310 
311  /* Copy the column TIME */
312  if (cpl_table_has_column (pix_table, "TIME")) {
313  cpl_table_duplicate_column (output_table, "TIME", pix_table, "TIME");
314  CPLCHECK_NUL ("Cannot get TIME data");
315  }
316 
317  /* Get pointer to the sky mean [e] out of the loop */
318  double * pSky = cpl_calloc (sizex * sizey, sizeof(double));
319  if (skyavg_table) {
320  for (cpl_size pix = 0; pix < npix; pix++) {
321  pSky[pix] = gravi_table_get_value (skyavg_table, "PIX", 0, pix) / gain;
322  }
323  }
324  CPLCHECK_NUL ("Cannot get the sky data");
325 
326  /* Get pointer to the sky variance [e^2] out of the loop */
327  double * pSkyVar = cpl_calloc (sizex * sizey, sizeof(double));
328  if (skystd_table) {
329  for (cpl_size pix = 0; pix < npix; pix++) {
330  pSkyVar[pix] = gravi_table_get_value (skystd_table, "PIX", 0, pix);
331  pSkyVar[pix] = pow (pSkyVar[pix] / gain, 2.0);
332  }
333  }
334  CPLCHECK_NUL ("Cannot get sky variance");
335 
336 
337  /* Loop on regions */
338  for (cpl_size region = 0; region < n_output; region++) {
339 
340  /* Loop on polarisation */
341  for (int pol = 0; pol < npol; pol ++) {
342 
343  /* Verbose every 6 regions */
344  if ( !(region+pol) || !((region*npol+pol+1)%6) )
345  cpl_msg_info_overwritable (cpl_func,
346  "Extract region of FT %lld over %lld (fast-no-cpl)",
347  (region*npol+pol+1), n_output*npol);
348 
349  /* Create DATA column */
350  const char * data = GRAVI_DATA[region*npol + pol];
351  cpl_table_new_column_array (output_table, data, CPL_TYPE_DOUBLE, ny_out);
352  cpl_array ** tData = cpl_table_get_data_array (output_table, data);
353  CPLCHECK_NUL ("Cannot create DATA column");
354 
355  /* Create DATAERR column */
356  const char * dataerr = GRAVI_DATAERR[region*npol + pol];
357  cpl_table_new_column_array (output_table, dataerr, CPL_TYPE_DOUBLE, ny_out);
358  cpl_array ** tDataErr = cpl_table_get_data_array (output_table, dataerr);
359  CPLCHECK_NUL ("Cannot create DATAERR column");
360 
361  /* Compute the flux in in [e] by loop on row and spectral direction
362  * - Desinterlace the FT frames
363  * - data = data - mean_sky
364  * - var_data = data - mean_sky + var_sky
365  * If 2 pixels per element
366  * - data = data1 - mean_sky1 + data2 - mean_sky2
367  * - var_data = data + var_sky1 + var_sky2
368  * Only considere the pixel that are not flag as BAD
369  * Allocating the arrays wrap data takes most of
370  * the time of the extract (80%) */
371 
372  for (cpl_size row = 0; row < nrow; row ++) {
373  double *pData = cpl_malloc (ny_out * sizeof(double));
374  double *pDataErr = cpl_malloc (ny_out * sizeof(double));
375  for (cpl_size j = 0; j < ny_out; j++) {
376  long idx = sizex * (j + ny*pol) + region*nx;
377  double value = 0;
378  double var = 0;
379  if (cpl_array_get (arr_badft[0], idx, NULL) == 0) { // if not bad
380  value += cpl_array_get (arr_data[row], idx, NULL) / gain - pSky[idx];
381  var += pSkyVar[idx];
382  }
383  if (nx>1 && cpl_array_get (arr_badft[0], idx+1, NULL) == 0) { // if not bad and if there is 2 pixels
384  value += cpl_array_get (arr_data[row], idx+1, NULL) / gain - pSky[idx+1];
385  var += pSkyVar[idx+1];
386  }
387  pData[j] = value;
388  pDataErr[j] = sqrt (CPL_MAX (value,0.0) + var);
389  }
390  tData[row] = cpl_array_wrap_double (pData, ny_out);
391  tDataErr[row] = cpl_array_wrap_double (pDataErr, ny_out);
392 
393  CPLCHECK_NUL ("Cannot extract and wrap data");
394  } /* End loop on rows */
395 
396  CPLCHECK_NUL ("Cannot compute region");
397  }
398  /* End loop on polarisation */
399  }
400  /* End loop on regions */
401 
402  FREE (cpl_free, pSky);
403  FREE (cpl_free, pSkyVar);
404 
405  CPLCHECK_NUL ("Cannot format ft");
406 
407  gravi_msg_function_exit(0);
408  return output_table;
409 }
410 
411 /*----------------------------------------------------------------------------*/
434 /*----------------------------------------------------------------------------*/
435 
436 cpl_table * gravi_imglist_sc_collapse (cpl_table * profile_table,
437  cpl_imagelist * raw_imglist,
438  cpl_imagelist * rawVar_imglist,
439  cpl_size startx)
440 {
441  int nv;
442  gravi_msg_function_start(1);
443 
444  cpl_ensure (profile_table, CPL_ERROR_NULL_INPUT, NULL);
445  cpl_ensure (raw_imglist, CPL_ERROR_NULL_INPUT, NULL);
446  cpl_ensure (rawVar_imglist, CPL_ERROR_NULL_INPUT, NULL);
447 
448  /* Get data */
449  cpl_size n_region = gravi_spectrum_get_nregion (profile_table);
450  cpl_size n_row = cpl_imagelist_get_size (raw_imglist);
451  cpl_ensure (n_row > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
452  cpl_ensure (n_region==24 || n_region==48, CPL_ERROR_ILLEGAL_INPUT, NULL);
453 
454  /* Ensure the input images are of type DOUBLE */
455  cpl_type type0 = cpl_image_get_type (cpl_imagelist_get (raw_imglist, 0));
456  cpl_type type1 = cpl_image_get_type (cpl_imagelist_get (rawVar_imglist, 0));
457  cpl_ensure (type0 == CPL_TYPE_DOUBLE, CPL_ERROR_ILLEGAL_INPUT, NULL);
458  cpl_ensure (type1 == CPL_TYPE_DOUBLE, CPL_ERROR_ILLEGAL_INPUT, NULL);
459 
460  /* Get output tables */
461  cpl_table * spectrum_table = cpl_table_new (n_row);
462 
463  /* Get the profile dimension */
464  cpl_size nx = cpl_table_get_column_dimension (profile_table,"DATA1",0);
465  cpl_size ny = cpl_table_get_column_dimension (profile_table,"DATA1",1);
466  CPLCHECK_NUL ("Cannot get profile dimension");
467 
468  /* Loop on regions (output of beam combiner) */
469  for (cpl_size region = 0; region < n_region; region++){
470 
471  /* Verbose every 6 regions */
472  if ( !region || !((region+1)%6) )
473  cpl_msg_info_overwritable(cpl_func, "Extract region of SC %lld over %lld",region+1,n_region);
474 
475  /* Construction of the DATA column */
476  const char * regionNameData = GRAVI_DATA[region];
477  cpl_table_new_column_array (spectrum_table, regionNameData, CPL_TYPE_DOUBLE, nx);
478 
479  /* Construction of the DATAERR column */
480  const char * regionNameErr = GRAVI_DATAERR[region];
481  cpl_table_new_column_array (spectrum_table, regionNameErr, CPL_TYPE_DOUBLE, nx);
482 
483  /* Get the profile of this region as a 2D image */
484  cpl_imagelist * profile_imglist = gravi_imagelist_wrap_column (profile_table, GRAVI_DATA[region]);
485  cpl_image * profile_img = cpl_imagelist_get (profile_imglist, 0);
486  CPLCHECK_NUL ("Cannot get data");
487 
488  /* Get the bound of the profile of this region.
489  * x = spectral, y = spatial */
490  int xmin = 1, xmax = nx;
491  int ymin = ny, ymax = 0;
492  for (cpl_size jy = 1 ; jy <= ny ; jy++ ) {
493  for (cpl_size ix = 1 ; ix <= nx ; ix++ ) {
494  if ( !cpl_image_get (profile_img, ix, jy, &nv) ) continue;
495  if ( jy<ymin ) ymin=jy;
496  if ( jy>ymax ) ymax=jy;
497  }
498  }
499  CPLCHECK_NUL ("Cannot get profile limits");
500 
501  /* Dump the found limits */
502  cpl_msg_debug (cpl_func, "Found limits: x=[%4d, %4d] y=[%4d, %4d] (nx=%lld,ny=%lld, FITS convention)",
503  xmin,xmax,ymin,ymax,nx,ny);
504 
505  /* Extract a cropped version of the profile */
506  cpl_image * profile_crop = cpl_image_extract (profile_img, xmin, ymin, xmax, ymax);
507  CPLCHECK_NUL ("Cannot extract profile");
508 
509  /* Get the pointer to rows, to avoid calling this into the loop */
510  cpl_array ** tData = cpl_table_get_data_array (spectrum_table, regionNameData);
511  cpl_array ** tDataErr = cpl_table_get_data_array (spectrum_table, regionNameErr);
512 
513  /* Loop on frame */
514  for (cpl_size row = 0; row < n_row; row ++) {
515 
516  /* Extracted flux in [e] for this frame
517  * rawFlux = < image * profile > */
518  cpl_image * rawFlux_profiled = cpl_image_extract (cpl_imagelist_get (raw_imglist, row),
519  xmin+startx-1, ymin, xmax+startx-1, ymax);
520  cpl_image_multiply (rawFlux_profiled, profile_crop);
521 
522  cpl_image *rawFlux = cpl_image_collapse_create (rawFlux_profiled,0);
523  cpl_image_delete (rawFlux_profiled);
524  CPLCHECK_NUL ("Cannot collapse flux");
525 
526  /* Extracted variance for this frame in [e]
527  * rawVar = < variance * profile^2 > */
528  cpl_image * rawVar_profiled = cpl_image_extract (cpl_imagelist_get (rawVar_imglist, row),
529  xmin+startx-1, ymin, xmax+startx-1, ymax);
530  cpl_image_multiply (rawVar_profiled, profile_crop);
531  cpl_image_multiply (rawVar_profiled, profile_crop);
532 
533  cpl_image *rawErr = cpl_image_collapse_create (rawVar_profiled,0);
534  cpl_image_delete (rawVar_profiled);
535  cpl_image_threshold (rawErr, 0.0, DBL_MAX, 0.0, DBL_MAX);
536  cpl_image_power (rawErr, 0.5);
537  CPLCHECK_NUL ("Cannot collapse variance");
538 
539  /* Fill the output table for the given region and frame : flux in [e] */
540  tData[row] = cpl_array_wrap_double (cpl_image_get_data_double (rawFlux), nx);
541 
542  /* Fill the output table for the given region and frame : error in [e] */
543  tDataErr[row] = cpl_array_wrap_double (cpl_image_get_data_double (rawErr), nx);
544 
545  /* Delete tmp images */
546  FREE (cpl_image_unwrap, rawFlux);
547  FREE (cpl_image_unwrap, rawErr);
548  }
549 
550  /* Delete all arrays of the loop */
551  FREE (gravi_imagelist_unwrap_images, profile_imglist);
552  FREE (cpl_image_delete, profile_crop);
553  }
554  /* End loop on region */
555 
556  gravi_msg_function_exit(1);
557  return spectrum_table;
558 }
559 
560 /*----------------------------------------------------------------------------*/
588 /*----------------------------------------------------------------------------*/
589 
590 gravi_data * gravi_extract_spectrum (gravi_data * raw_data,
591  gravi_data * profile_map,
592  gravi_data * dark_map,
593  gravi_data * bad_map,
594  gravi_data * sky_map,
595  const cpl_parameterlist * parlist,
596  enum gravi_detector_type det_type)
597 {
598  gravi_msg_function_start(1);
599  cpl_ensure (raw_data, CPL_ERROR_NULL_INPUT, NULL);
600  cpl_ensure (profile_map, CPL_ERROR_NULL_INPUT, NULL);
601  cpl_ensure (bad_map, CPL_ERROR_NULL_INPUT, NULL);
602  cpl_ensure (dark_map || sky_map, CPL_ERROR_NULL_INPUT, NULL);
603  cpl_ensure (parlist, CPL_ERROR_NULL_INPUT, NULL);
604 
605  /* Get header of input gravi_data */
606  cpl_propertylist * raw_header = gravi_data_get_header (raw_data);
607  /* cpl_propertylist * profile_header = gravi_data_get_header (profile_map); */
608 
609  /* Create output gravi_data */
610  gravi_data * spectrum_data = gravi_data_new(0);
611  cpl_propertylist * spectrum_header = gravi_data_get_header (spectrum_data);
612 
613  /* Dump all header from RAW product */
614  cpl_propertylist_append (spectrum_header, raw_header);
615 
616  if ((det_type == GRAVI_DET_FT || det_type == GRAVI_DET_ALL) &&
617  gravi_data_has_detector (raw_data, GRAVI_FT)) {
618 
619  /*
620  * Compute the flux and variance of each spectral element of FT
621  * The spectrum are already extracted from spatial direction.
622  */
623 
624  /* Dupplicate necessary extension in the output spectrum_data */
625  gravi_data_copy_ext (spectrum_data, raw_data, GRAVI_IMAGING_DETECTOR_FT_EXT);
626  CPLCHECK_NUL ("Cannot duplicate necessary extensions");
627 
628  /* Get the data */
629  cpl_table * imaging_data_ft = gravi_data_get_table (raw_data, GRAVI_IMAGING_DATA_FT_EXT);
630  cpl_size n_region = cpl_table_get_nrow (gravi_data_get_table (raw_data,GRAVI_IMAGING_DETECTOR_FT_EXT));
631 
632  /* Get the background mean and std. In case we have both SKY and DARK we still use
633  * only the SKY because it is surely the closest in time (taken at night) and with the same setup
634  * FIXME: could be interesting to check the setup and time distance... */
635  cpl_table * skyavg_table, * skystd_table;
636  if (sky_map != NULL) {
637  cpl_msg_info (cpl_func, "Extract FT spectra with SKY as background and variance");
638  skyavg_table = gravi_data_get_table (sky_map, GRAVI_IMAGING_DATA_FT_EXT);
639  skystd_table = gravi_data_get_table (sky_map, GRAVI_IMAGING_ERR_FT_EXT);
640  }
641  else if (dark_map != NULL) {
642  if ( !gravi_data_is_internal(raw_data) )
643  gravi_pfits_add_check (spectrum_header, "Extract FT spectra with DARK as background and variance");
644  else
645  cpl_msg_info (cpl_func, "Extract FT spectra with DARK as background and variance");
646  skyavg_table = gravi_data_get_table (dark_map, GRAVI_IMAGING_DATA_FT_EXT);
647  skystd_table = gravi_data_get_table (dark_map, GRAVI_IMAGING_ERR_FT_EXT);
648  }
649 
650  /* Get the FT gain in [ADU/e] */
651  double gain_ft = gravi_pfits_get_ft_gain (raw_header);
652  cpl_propertylist_update_double (spectrum_header, "ESO QC USEDGAIN FT", gain_ft);
653  cpl_propertylist_set_comment (spectrum_header, "ESO QC USEDGAIN FT", "[ADU/e-] value used for reduction");
654 
655  CPLCHECK_NUL("Problem while getting the tables");
656 
657  /* Get the badpixel image */
658  cpl_table * badft_table = gravi_data_get_table(bad_map, GRAVI_IMAGING_DATA_FT_EXT);
659 
660  /* Convert PIX column to DATA# and DATAERR# */
661  cpl_table * spectrum_ft;
662  spectrum_ft = gravi_table_ft_format (imaging_data_ft, skystd_table, skyavg_table, badft_table, n_region, gain_ft);
663  CPLCHECK_NUL ("Cannot format FT data");
664 
665  /* Set units */
666  for (cpl_size reg=0; reg<n_region; reg++) {
667  cpl_table_set_column_unit (spectrum_ft, GRAVI_DATA[reg], "e");
668  cpl_table_set_column_unit (spectrum_ft, GRAVI_DATAERR[reg], "e");
669  }
670 
671  /* Set the data in output array */
672  cpl_propertylist * spectrum_plist = cpl_propertylist_new ();
673  int nx = gravi_spectrum_get_nwave (spectrum_ft);
674  cpl_propertylist_update_int (spectrum_plist, PROFILE_NX, nx);
675  cpl_propertylist_update_int (spectrum_plist, PROFILE_STARTX, 1);
676  cpl_propertylist_update_int (spectrum_plist, PROFILE_FULLSTARTX, 0);
677 
678  gravi_data_add_table (spectrum_data, spectrum_plist, GRAVI_SPECTRUM_DATA_FT_EXT, spectrum_ft);
679  }
680 
681 
682  if ((det_type == GRAVI_DET_SC || det_type == GRAVI_DET_ALL) &&
683  gravi_data_has_detector (raw_data, GRAVI_SC)) {
684 
685  /*
686  * Compute the flux and variance of
687  * each spectral element of SC
688  */
689 
690  /* Dupplicate necessary extension in the output spectrum_data */
691  gravi_data_copy_ext (spectrum_data, raw_data, GRAVI_IMAGING_DETECTOR_SC_EXT);
692  CPLCHECK_NUL ("Cannot duplicate necessary extensions");
693 
694  /* Get the profile data */
695  cpl_table * profile_table = gravi_data_get_table (profile_map, GRAVI_PROFILE_DATA_EXT);
696  cpl_ensure (profile_table, CPL_ERROR_ILLEGAL_INPUT, NULL);
697 
698  /* Get the data */
699  cpl_imagelist * imaging_data = gravi_data_get_cube (raw_data, GRAVI_IMAGING_DATA_SC_EXT);
700  cpl_propertylist * profile_plist = gravi_data_get_plist (profile_map, GRAVI_PROFILE_DATA_EXT);
701  CPLCHECK_NUL ("Cannot get data");
702 
703  /* Get the SC gain in [ADU/e] */
704  double gain_sc = gravi_pfits_get_sc_gain (raw_header);
705  cpl_propertylist_update_double (spectrum_header, "ESO QC USEDGAIN SC", gain_sc);
706  cpl_propertylist_set_comment (spectrum_header, "ESO QC USEDGAIN SC", "[ADU/e-] value used for reduction");
707 
708  /* Get the badpixel image */
709  cpl_image * badpix_img = gravi_data_get_img (bad_map, GRAVI_IMAGING_DATA_SC_EXT);
710 
711  /* Estimate the stellar flux in [e].
712  * Give priority to SKY for estimating the background */
713  cpl_image * skyavg_img;
714  if (sky_map != NULL){
715  cpl_msg_info (cpl_func,"Extract SC spectra with SKY as background");
716  skyavg_img = gravi_data_get_img (sky_map, GRAVI_IMAGING_DATA_SC_EXT);
717  }
718  else if (dark_map != NULL) {
719  if ( !gravi_data_is_internal(raw_data) )
720  gravi_pfits_add_check (spectrum_header, "Extract SC spectra with DARK as background");
721  else
722  cpl_msg_info (cpl_func,"Extract SC spectra with DARK as background");
723  skyavg_img = gravi_data_get_img (dark_map, GRAVI_IMAGING_DATA_SC_EXT);
724  }
725 
726  /* Estimate the total variance in [e^2] composed of background variance
727  * and additional photonic variance. Give priority to DARK for estimating
728  * the background variance */
729  cpl_image * darkavg_img;
730  cpl_image * darkstd_img;
731  if (dark_map != NULL) {
732  cpl_msg_info (cpl_func,"Extract SC photonic variance with DARK as background and variance");
733  darkavg_img = gravi_data_get_img (dark_map, GRAVI_IMAGING_DATA_SC_EXT);
734  darkstd_img = gravi_data_get_img (dark_map, GRAVI_IMAGING_ERR_SC_EXT);
735  }
736  else if (sky_map != NULL) {
737  gravi_pfits_add_check (spectrum_header, "Extract SC photonic variance with SKY as background and variance");
738  darkavg_img = gravi_data_get_img (sky_map, GRAVI_IMAGING_DATA_SC_EXT);
739  darkstd_img = gravi_data_get_img (sky_map, GRAVI_IMAGING_ERR_SC_EXT);
740  }
741 
742  /* Build the dimension for the raw_imglist, that shall match the profile dimesion */
743  cpl_size startx = gravi_pfits_get_startx (profile_plist);
744  CPLCHECK_NUL ("The coordinate dimensions of the new window is missing");
745 
746  /* raw = (DATA - SKY) / gain [e] */
747  cpl_msg_info (cpl_func,"Compute flux image");
748  cpl_imagelist * raw_imglist;
749  raw_imglist = cpl_imagelist_duplicate (imaging_data);
750  cpl_imagelist_subtract_image (raw_imglist, skyavg_img);
751  cpl_imagelist_divide_scalar (raw_imglist, gain_sc);
752  gravi_remove_badpixel_sc (raw_imglist, badpix_img);
753  CPLCHECK_NUL ("Cannot extract the data");
754 
755  /* rawVar = (DATA - DARK) / gain [e^2] */
756  cpl_msg_info (cpl_func,"Compute variance image");
757  cpl_imagelist * rawVar_imglist;
758  rawVar_imglist = cpl_imagelist_duplicate (imaging_data);
759  cpl_imagelist_subtract_image (rawVar_imglist, darkavg_img);
760  cpl_imagelist_divide_scalar (rawVar_imglist, gain_sc);
761  CPLCHECK_NUL ("Cannot extract the data");
762 
763  /* rawVar += (std(DARK) / gain) ** 2 [e^2] */
764  cpl_image * darkvar_img = cpl_image_power_create (darkstd_img, 2.0);
765  cpl_image_divide_scalar (darkvar_img, gain_sc * gain_sc);
766  cpl_imagelist_add_image (rawVar_imglist, darkvar_img);
767  FREE (cpl_image_delete, darkvar_img);
768  CPLCHECK_NUL ("Cannot extract the data");
769 
770  /* Collapse images with spectrum and create
771  * spectrum_table for SC */
772  cpl_table * spectrum_table;
773  spectrum_table = gravi_imglist_sc_collapse (profile_table, raw_imglist,
774  rawVar_imglist, startx);
775  CPLCHECK_NUL ("Cannot collapse the spectrum");
776 
777  cpl_size n_row = cpl_table_get_nrow (spectrum_table);
778  cpl_size n_region = cpl_table_get_ncol (profile_table);
779 
780  /* Set units */
781  for (cpl_size reg=0; reg<n_region; reg++) {
782  cpl_table_set_column_unit (spectrum_table, GRAVI_DATA[reg], "e");
783  cpl_table_set_column_unit (spectrum_table, GRAVI_DATAERR[reg], "e");
784  }
785 
786  /* Get the user-specified DIT shift if any */
787  int ditshift = gravi_param_get_int_default (parlist, "gravity.preproc.ditshift-sc", 0);
788  cpl_msg_info (cpl_func, "ditshift = %i", ditshift);
789 
790  if (ditshift !=0 )
791  gravi_msg_warning ("CRITICAL","DITSHIFT is not set to zero... are you sure??");
792 
793  /* Compute the time of the middle of the SC DIT, in [us]
794  * with respect to the RMN start PCR.ACQ.START */
795  cpl_table_new_column (spectrum_table, "TIME", CPL_TYPE_INT);
796  for (cpl_size row = 0; row < n_row; row ++) {
797  cpl_table_set_int (spectrum_table, "TIME", row,
798  gravi_pfits_get_time_sc (raw_header, row + ditshift));
799  }
800 
801  /* Check if the SC profil extraction is flux conservative */
802  double full_flux_img = gravi_imagelist_get_flux (raw_imglist);
803  double full_flux_reg = gravi_spectrum_get_flux (spectrum_table);
804 
805  cpl_msg_info (cpl_func, "Total flux in REGIONs: %.2e [e], in IMGs:%.2e [e] (ratio=%.5f)",
806  full_flux_reg,full_flux_img,full_flux_reg/full_flux_img);
807 
808  cpl_propertylist_update_double (spectrum_header, "ESO QC TRANS PROFILE SC", full_flux_reg/full_flux_img);
809  cpl_propertylist_set_comment (spectrum_header, "ESO QC TRANS PROFILE SC", "[e/e] at profile extraction");
810 
811  CPLCHECK_NUL ("Cannot verify if extraction was flux-conservative");
812 
813  /* Set the SPECTRUM_DATA_SC tables to the gravi_data */
814  cpl_propertylist * spectrum_plist = cpl_propertylist_new ();
815  cpl_propertylist_copy_property (spectrum_plist, profile_plist, PROFILE_FULLSTARTX);
816  cpl_propertylist_copy_property (spectrum_plist, profile_plist, PROFILE_STARTX);
817  cpl_propertylist_copy_property (spectrum_plist, profile_plist, PROFILE_NX);
818 
819  gravi_data_add_table (spectrum_data, spectrum_plist,
820  GRAVI_SPECTRUM_DATA_SC_EXT, spectrum_table);
821 
822 
823  /* Apply the same extraction to the flat
824  * FIXME: no units since it was not in e- */
825  cpl_msg_info (cpl_func, "Extract the FLAT");
826 
827  cpl_table * spectrumflat_table;
828  cpl_imagelist * flat_imglist = gravi_data_get_cube (profile_map, GRAVI_IMAGING_DATA_SC_EXT);
829  spectrumflat_table = gravi_imglist_sc_collapse (profile_table, flat_imglist,
830  flat_imglist, 1);
831  CPLCHECK_NUL ("Cannot collapse the flat");
832 
833  /* Set the SPECTRUMFLAT_DATA_SC tables to the gravi_data */
834  cpl_propertylist * spectrumflat_plist = cpl_propertylist_new ();
835  gravi_data_add_table (spectrum_data, spectrumflat_plist,
836  GRAVI_SPECTRUMFLAT_DATA_SC_EXT, spectrumflat_table);
837 
838  /* Delete variable */
839  FREE (cpl_imagelist_delete, raw_imglist);
840  FREE (cpl_imagelist_delete, rawVar_imglist);
841  } /* End SC case */
842 
843  gravi_msg_function_exit(1);
844  return spectrum_data;
845 }
846 
847 
848 
849 /*----------------------------------------------------------------------------*/
858 /*----------------------------------------------------------------------------*/
859 
860 cpl_error_code gravi_subtract_met_dark (gravi_data * preproc_data,
861  gravi_data * dark_map)
862 {
863 
864  gravi_msg_function_start(1);
865  cpl_ensure_code (dark_map , CPL_ERROR_NULL_INPUT);
866 
867  if (!gravi_data_has_extension (dark_map, GRAVI_METROLOGY_EXT)) {
868  cpl_msg_warning (cpl_func,"The DARK data has no METROLOGY extension");
869  }
870  else
871  {
872  cpl_table * met_data_table = gravi_data_get_table (preproc_data, GRAVI_METROLOGY_EXT);
873  cpl_table * met_dark_table = gravi_data_get_table (dark_map, GRAVI_METROLOGY_EXT);
874  CPLCHECK_MSG("Problem while getting the tables");
875 
876  /* Get pointer to the data */
877  const char * data_x="VOLT";
878  cpl_size nrow = cpl_table_get_nrow (met_data_table);
879  cpl_array ** met_data_array = cpl_table_get_data_array (met_data_table, data_x);
880  cpl_ensure_code (met_data_array, CPL_ERROR_ILLEGAL_INPUT);
881  cpl_array ** met_dark_array = cpl_table_get_data_array (met_dark_table, data_x);
882  cpl_ensure_code (met_dark_array, CPL_ERROR_ILLEGAL_INPUT);
883 
884  /* Subtract data */
885 
886  for (cpl_size j = 0; j < nrow ; j++)
887  {
888  cpl_array_subtract(met_data_array[j],met_dark_array[0]);
889  }
890 
891  cpl_msg_info (cpl_func,"METROLOGY Volts have been subtracted from dark");
892  }
893 
894  gravi_msg_function_exit(1);
895  return CPL_ERROR_NONE;
896 
897 }
898 
899 /*----------------------------------------------------------------------------*/
918 /*----------------------------------------------------------------------------*/
919 
920 cpl_error_code gravi_interpolate_spectrum_table (cpl_table * spectrum_table,
921  cpl_table * wave_table,
922  cpl_table ** oiwave_tables,
923  cpl_table * detector_table,
924  cpl_table * specflat_table,
925  const cpl_parameterlist * parlist,
926  int type_data)
927 {
928  gravi_msg_function_start(1);
929 
930  /* Check needed data */
931  cpl_ensure_code (spectrum_table, CPL_ERROR_NULL_INPUT);
932  cpl_ensure_code (wave_table, CPL_ERROR_NULL_INPUT);
933  cpl_ensure_code (oiwave_tables, CPL_ERROR_NULL_INPUT);
934  cpl_ensure_code (oiwave_tables[0], CPL_ERROR_NULL_INPUT);
935 
936  /* Sizes */
937  cpl_size nb_region = gravi_spectrum_get_nregion (spectrum_table);
938  cpl_size nb_pol = gravi_spectrum_get_npol (spectrum_table);
939  cpl_size nb_wave = gravi_spectrum_get_nwave (spectrum_table);
940  cpl_size nb_row = cpl_table_get_nrow (spectrum_table);
941  cpl_size nb_oiwave = cpl_table_get_nrow (oiwave_tables[0]);
942 
943  /* Ensure the two OI_WAVELENGTH have same size */
944  if (nb_pol > 1) {
945  cpl_ensure_code (oiwave_tables[1], CPL_ERROR_NULL_INPUT);
946  cpl_ensure_code (cpl_table_get_nrow (oiwave_tables[0]) ==
947  cpl_table_get_nrow (oiwave_tables[1]),
948  CPL_ERROR_ILLEGAL_INPUT);
949  }
950 
951  /* Only 24 and 48 regions are accepted */
952  cpl_ensure_code (nb_region == 24 || nb_region == 48,
953  CPL_ERROR_ILLEGAL_INPUT);
954 
955  /* So far oversampling is not supported */
956  cpl_msg_info (cpl_func,"nb_oiwave=%lld, nb_wave=%lld", nb_oiwave, nb_wave);
957  cpl_ensure_code (nb_oiwave <= nb_wave, CPL_ERROR_ILLEGAL_INPUT);
958 
959  /* Verbose about the target effective map */
960  if (specflat_table)
961  cpl_msg_info (cpl_func, "Change the target wavelength because of FLAT");
962  else
963  cpl_msg_info (cpl_func, "Don't change the target wavelenght because of FLAT");
964 
965 
966  /* Allocate memory for the temporary target wavelength (xout)
967  * Allocate memory for the temporary output (yout) */
968  double * xout = cpl_malloc (nb_oiwave * sizeof(double));
969  double * yout = cpl_malloc (nb_oiwave * sizeof(double));
970  cpl_size * id = cpl_malloc (nb_oiwave * sizeof(cpl_size));
971 
972  double * outputData = cpl_malloc (nb_oiwave * sizeof(double));
973  double * outputErr = cpl_malloc (nb_oiwave * sizeof(double));
974 
975  /* Loop on polarisations */
976  for (int pol = 0; pol < nb_pol; pol++) {
977 
978  /* Get max and min of the new wavelength table */
979  double oiwave_min = cpl_table_get_column_min (oiwave_tables[pol], "EFF_WAVE");
980  double oiwave_max = cpl_table_get_column_max (oiwave_tables[pol], "EFF_WAVE");
981 
982  cpl_msg_info(cpl_func,"Reinterpolate pol %i: oiwave_min=%g, oiwave_max=%g [um]",
983  pol, oiwave_min*1e6, oiwave_max*1e6);
984 
985  /* Loop on regions -- skip regions not in requested polarisation */
986  for (cpl_size reg = 0; reg < nb_region; reg++){
987  if (gravi_region_get_pol (detector_table, reg) != pol) continue;
988 
989  /* Verbose every regions */
990 
991  /* Name and dimension of this region */
992  const char * data_x = GRAVI_DATA[reg];
993  const char * data_errx = GRAVI_DATAERR[reg];
994 
995  /* Check size of SPECTRUM of this region */
996  int nbd = cpl_table_get_column_depth (spectrum_table, data_x);
997  cpl_ensure_code (nbd == nb_wave, CPL_ERROR_ILLEGAL_INPUT);
998 
999  /* Check size of WAVE of this region */
1000  int nbw = cpl_table_get_column_depth (wave_table, data_x);
1001  cpl_ensure_code (nbw == nb_wave, CPL_ERROR_ILLEGAL_INPUT);
1002 
1003  /* Ensure we have no extrapolation and have increasing wavelengths */
1004  double wave_first = gravi_table_get_value (wave_table, data_x, 0, 0);
1005  double wave_last = gravi_table_get_value (wave_table, data_x, 0, nb_wave-1);
1006  cpl_ensure_code (wave_first < wave_last, CPL_ERROR_ILLEGAL_INPUT);
1007  //cpl_ensure_code (oiwave_min > wave_first, CPL_ERROR_ILLEGAL_INPUT);
1008  //cpl_ensure_code (oiwave_max < wave_last, CPL_ERROR_ILLEGAL_INPUT);
1009 
1010 // /* Get directly the pointer to arrays to speed-up */
1011 // cpl_array ** inputData = cpl_table_get_data_array (spectrum_table, data_x);
1012 // cpl_array ** inputErr = cpl_table_get_data_array (spectrum_table, data_errx);
1013 // CPLCHECK_MSG("Cannot get data");
1014 
1015 
1016  /* Get directly the pointer to arrays to speed-up */
1017  cpl_array ** inputData = cpl_table_get_data_array (spectrum_table, data_x);
1018  cpl_array ** inputErr = cpl_table_get_data_array (spectrum_table, data_errx);
1019  const cpl_array * flat = NULL;
1020  if (specflat_table){
1021  flat = cpl_table_get_array (specflat_table, data_x, 0);
1022  }
1023  cpl_array * wave = cpl_table_get_data_array(wave_table, data_x)[0];
1024  CPLCHECK_MSG("Cannot get data");
1025 
1026  /* Verify the data are double since we will get pointer to them */
1027  cpl_ensure_code ((cpl_array_get_type (inputData[0]) == CPL_TYPE_DOUBLE) &&
1028  (cpl_array_get_type (inputErr[0]) == CPL_TYPE_DOUBLE) &&
1029  (cpl_array_get_type (wave) == CPL_TYPE_DOUBLE),
1030  CPL_ERROR_ILLEGAL_INPUT);
1031 
1032 
1033  /* Get the current wavelength as pointer */
1034  double * wave_ref = cpl_array_get_data_double (wave);
1035  CPLCHECK_MSG("Cannot get data");
1036 
1037 
1038  /* Check the Option
1039  * Go to 3 pixels interp if :
1040  * 1- if LOW ie. nwave < 50 (old: Option set)
1041  * 2- Flat is given in input
1042  * 3- Not FT data
1043  */
1044  //if ( (gravi_param_get_bool (parlist,"gravity.preproc.interp-3pixels") == TRUE)
1045  // && (flat != NULL) && ( nb_oiwave !=5 ))
1046  if ( ( nb_oiwave <= 50 ) && ( type_data == GRAVI_SC) && (flat != NULL) ) {
1047  cpl_msg_info_overwritable (cpl_func, "Reinterpolate (3 pixels) region %lld "
1048  "over %lld (%lld->%lld channels)",
1049  reg+1,nb_region,nb_wave,nb_oiwave);
1050  /* Init the index and weights for the interpolation such that
1051  * yout[iw] = inputData[id[iw]-1] * weight_m[iw]
1052  * + inputData[id[iw]] * weight_m[iw]
1053  * + inputData[id[iw]+1] * weight_p[iw]
1054  * yerr[iw] = sqrt( (inputErr[id[iw]-1] * weight_m[iw])^2
1055  * + (inputErr[id[iw]] * weight_m[iw])^2
1056  * + (inputErr[id[iw]+1] * weight_p[iw])^2 )
1057  * We here assume monotonicity and no extrapolation */
1058  double * weight_m = cpl_malloc (nb_oiwave * sizeof(double));
1059  double * weight_i = cpl_malloc (nb_oiwave * sizeof(double));
1060  double * weight_p = cpl_malloc (nb_oiwave * sizeof(double));
1061 
1062  for (cpl_size iw = 0 ; iw < nb_oiwave ; iw ++)
1063  {
1064  // xo is the wavenumber of the targeted wavelength (of indice iw)
1065  double l0=cpl_table_get (oiwave_tables[pol], "EFF_WAVE", iw, NULL);
1066 
1067  // id[iw] is the index of the closest wavenumber to xo
1068  // id[iw] must be above 1, and below nb_wave-1
1069  // wave is assumed to be in increasing order
1070 
1071  if (iw == 0) id[iw] = 1;
1072  else id[iw] = id[iw-1];
1073  while ( (1./l0 < 0.5/wave_ref[id[iw]]+0.5/wave_ref[id[iw]+1]) && ( id[iw] < nb_wave-2 ) ){
1074  id[iw]++;
1075  }
1076 
1077  // compute the wavenumbers (normalized to x0)
1078  // FE: declaration should be further up
1079  double xm=l0/wave_ref[id[iw]-1]-1.;
1080  double xi=l0/wave_ref[id[iw]]-1.;
1081  double xp=l0/wave_ref[id[iw]+1]-1.;
1082  double xT=(xp-xm)/2;
1083  double Fm = 1.;
1084  double Fi = 1.;
1085  double Fp = 1.;
1086  double norm= 1.;
1087 
1088  if (specflat_table)
1089  {
1090  Fm=cpl_array_get (flat, id[iw]-1, NULL);
1091  Fi=cpl_array_get (flat, id[iw], NULL);
1092  Fp=cpl_array_get (flat, id[iw]+1, NULL);
1093  }
1094 
1095  // Compute the weigths (initialized at zero)
1096  weight_m[iw]=0;
1097  weight_i[iw]=0;
1098  weight_p[iw]=0;
1099 
1100  // only compute the weights if the three flats are not null (bad pixel)
1101  // and if xo is below (xm+xi)/2
1102  // and if xo is above (xi+xp)/2
1103  if ( (Fi>1e-4) && (Fp>1e-4) && (Fm>1e-4) && (xm-xi>1e-4) && (xi-xp>1e-4) )
1104  {
1105 
1106  /* The three equations
1107  * 1-> fi = 1 => conservation
1108  * 2-> fm*xm + fi*xi + fp*xp = 0 => First moment
1109  * 3-> fm*xm^3 + fi*xi^3 + fp*xp^3 = 0 => Third moment
1110  are now approximated by:
1111  * 1-> fi = 1 => conservation
1112  * 2-> fm*xm + fi*xi + fp*xp = 0 => First moment
1113  * 3-> fm*xm*(-xT-xi)^2 + fi*xi^3 + fp*xp*(xT-xi)^2 = 0 => Third moment
1114  * so the solution is simpler, and do not diverge: */
1115 
1116  weight_p[iw]= (xT-2*xi)/(4*xp) ;
1117  weight_i[iw]= 1. ;
1118  weight_m[iw]= -(xT+2*xi)/(4*xm);
1119 
1120  /* Acounting for flat
1121  * W/=Flat
1122  */
1123  weight_p[iw]/=Fp;
1124  weight_i[iw]/=Fi;
1125  weight_m[iw]/=Fm;
1126 
1127  // normalize the weight to conserve the variance
1128  //norm = sqrt( pow(weight_m[iw],2) + pow(weight_i[iw],2) + pow(weight_p[iw],2));
1129  /* FE: the original 2 pixel interpolation applies normalization to preserve photon noise) */
1130  norm = weight_m[iw] + weight_i[iw] + weight_p[iw];
1131  /* SL: changed to simple sum(weights) =1 */
1132  weight_p[iw]/=norm;
1133  weight_i[iw]/=norm;
1134  weight_m[iw]/=norm;
1135 
1136  }
1137 
1138  /* FE: for debugging puposes*/
1139  // cpl_msg_info(cpl_func," iw: %i l0: %g id[iw]: %i wm: %g wi: %g wp: %g Fm: %g Fi: %g Fp: %g xm: %g xi: %g xp: %g norm: %g\n",
1140  // iw, l0, id[iw], weight_m[iw], weight_i[iw], weight_p[iw], Fm, Fi, Fp, xm, xi, xp,norm);
1141 
1142 
1143  }// End loop on nb_oiwave
1144 
1145  /* Loop on frames */
1146  for (cpl_size j = 0; j < nb_row; j++){
1147  double * yref;
1148 
1149  /*
1150  * Compute the fluxes
1151  */
1152  yref = cpl_array_get_data_double (inputData[j]);
1153  /* loop on wavelength */
1154  for (cpl_size iw=0 ; iw < nb_oiwave ; iw ++) {
1155  outputData[iw]=yref[id[iw]-1]*weight_m[iw] + yref[id[iw]]*weight_i[iw] + yref[id[iw]+1]*weight_p[iw];
1156  }
1157  /* Copy output array into input array */
1158  for (cpl_size iw = 0 ; iw < nb_oiwave ; iw ++){
1159  yref[iw]=outputData[iw];
1160  }
1161 
1162  /*
1163  * Compute the noise (sqrt of variance)
1164  */
1165  yref = cpl_array_get_data_double (inputErr[j]);
1166  /* loop on wavelength */
1167  for (cpl_size iw=0 ; iw < nb_oiwave ; iw ++) {
1168  outputErr[iw]=sqrt(pow(yref[id[iw]-1]*weight_m[iw],2) + pow(yref[id[iw]]*weight_i[iw],2) + pow(yref[id[iw]+1]*weight_p[iw],2));
1169  }
1170  /* Copy output array into input array */
1171  for (cpl_size iw = 0 ; iw < nb_oiwave ; iw ++){
1172  yref[iw]=outputErr[iw];
1173  }
1174  } /* end loop on frames */
1175 
1176  /* Free the weight for interpolation */
1177  FREE (cpl_free, weight_m);
1178  FREE (cpl_free, weight_i);
1179  FREE (cpl_free, weight_p);
1180 
1181 
1182  } /* End of interpolation on 3 pixels */
1183  else
1184  {
1185  /* starting interpolation 2 pixels */
1186  cpl_msg_info_overwritable (cpl_func, "Reinterpolate region %lld "
1187  "over %lld (%lld->%lld channels)",
1188  reg+1,nb_region,nb_wave,nb_oiwave);
1189  /* starting interpolation 2 pixels */
1190  /* Init the index and weights for the interpolation such that
1191  * yout[iw] = inputData[id[iw]-1] * weight_m[iw]
1192  * + inputData[id[iw]] * weight_p[iw]
1193  * yerr[iw] = sqrt( (inputErr[id[iw]-1] * weight_m[iw])^2
1194  * + (inputErr[id[iw]] * weight_p[iw])^2 )
1195  * We here assume monotonicity and no extrapolation */
1196  double * weight_m = cpl_malloc (nb_oiwave * sizeof(double));
1197  double * weight_p = cpl_malloc (nb_oiwave * sizeof(double));
1198 
1199  for (cpl_size iw = 0 ; iw < nb_oiwave ; iw ++)
1200  {
1201  // xo is the wavenumber of the targeted wavelength (of indice iw)
1202  double l0=cpl_table_get (oiwave_tables[pol], "EFF_WAVE", iw, NULL);
1203 
1204  // id[iw] is the index of the closest wavenumber to xo
1205  // id[iw] must be above 1, and below nb_wave-1
1206  // wave is assumed to be in increasing order
1207 
1208  if (iw == 0) id[iw] = 1;
1209  else id[iw] = id[iw-1];
1210  while ( (1./l0 < 1./wave_ref[id[iw]]) && ( id[iw] < nb_wave-1 ) ){
1211  id[iw]++;
1212  }
1213 
1214  // compute the wavenumbers (normalized to x0)
1215  // FE: declaration should be further up
1216  double xm=l0/wave_ref[id[iw]-1]-1.;
1217  double xp=l0/wave_ref[id[iw]]-1.;
1218  double Fm = 1.;
1219  double Fp = 1.;
1220  double norm = 1.;
1221 
1222 
1223  if (specflat_table)
1224  {
1225  Fm=cpl_array_get (flat, id[iw]-1, NULL);
1226  Fp=cpl_array_get (flat, id[iw], NULL);
1227  }
1228 
1229  // Compute the weigths (initialized at zero)
1230  weight_m[iw]=0;
1231  weight_p[iw]=0;
1232 
1233  // only compute the weights if the three flats are not null (bad pixel)
1234  // and if xo is below (xm+xi)/2
1235  // and if xo is above (xi+xp)/2
1236  if ( (Fp>1e-4) && (Fm>1e-4) && (xm>=0.0) && (xp<=0.0) )
1237  {
1238 
1239  /* The three equations
1240  * 1-> fm + fp = xm - xp => conservation
1241  * 2-> fm*xm + fp*xp = 0 => First moment */
1242 
1243  weight_p[iw]= xm ;
1244  weight_m[iw]= -xp ;
1245 
1246  /* Acounting for flat
1247  * W/=Flat
1248  */
1249  weight_p[iw]/=Fp;
1250  weight_m[iw]/=Fm;
1251 
1252  // normalize the weight to conserve the variance
1253  //norm = sqrt( pow(weight_m[iw],2) + pow(weight_p[iw],2));
1254  norm = weight_m[iw]+ weight_p[iw];
1255  /* FE: the original 2 pixel interpolation applies normalization to preserve photon noise) */
1256  /* SL: changed to simple sum(weights) =1 */
1257  weight_p[iw]/=norm;
1258  weight_m[iw]/=norm;
1259 
1260  }
1261 
1262  /* FE: for debugging puposes
1263  cpl_msg_info(cpl_func," reg: %i iw: %i l0: %g id[iw]: %i wm: %g wp: %g Fm: %g Fp: %g xm: %g xp: %g \n",
1264  iw, iw, l0, id[iw], weight_m[iw], weight_p[iw], Fm, Fp, xm, xp);*/
1265 
1266 
1267 
1268  }// End loop on nb_oiwave
1269 
1270  /* Loop on frames */
1271  for (cpl_size j = 0; j < nb_row; j++){
1272  double * yref;
1273 
1274  /*
1275  * Compute the fluxes
1276  */
1277  yref = cpl_array_get_data_double (inputData[j]);
1278  /* loop on wavelength */
1279  for (cpl_size iw=0 ; iw < nb_oiwave ; iw ++) {
1280  outputData[iw]=yref[id[iw]-1]*weight_m[iw] + yref[id[iw]]*weight_p[iw];
1281  }
1282  /* Copy output array into input array */
1283  for (cpl_size iw = 0 ; iw < nb_oiwave ; iw ++){
1284  yref[iw]=outputData[iw];
1285  }
1286 
1287  /*
1288  * Compute the noise (sqrt of variance)
1289  */
1290  yref = cpl_array_get_data_double (inputErr[j]);
1291  /* loop on wavelength */
1292  for (cpl_size iw=0 ; iw < nb_oiwave ; iw ++) {
1293  outputErr[iw]=sqrt(pow(yref[id[iw]-1]*weight_m[iw],2) + pow(yref[id[iw]]*weight_p[iw],2));
1294  }
1295  /* Copy output array into input array */
1296  for (cpl_size iw = 0 ; iw < nb_oiwave ; iw ++){
1297  yref[iw]=outputErr[iw];
1298  }
1299  } /* end loop on frames */
1300 
1301  /* Free the weight for interpolation */
1302  FREE (cpl_free, weight_m);
1303  FREE (cpl_free, weight_p);
1304 
1305  /* end of interpolation on 2 pixels */
1306  }
1307 
1308  /* Modify the depth of the region (remove useless pixels) */
1309  if (nb_wave > nb_oiwave) {
1310  cpl_msg_debug (cpl_func,"Modify depth of column %s (%lld->%lld)", data_x, nb_wave, nb_oiwave);
1311  cpl_table_set_column_depth (spectrum_table, data_x, nb_oiwave);
1312  cpl_table_set_column_depth (spectrum_table, data_errx, nb_oiwave);
1313  CPLCHECK_MSG ("Cannot change column depth");
1314  }
1315 
1316  } /* End loop on regions i */
1317 
1318  } /* End loop on polarisation */
1319 
1320  /* Desallocate the modified wavelength array */
1321  FREE (cpl_free, yout);
1322  FREE (cpl_free, xout);
1323 
1324  /* Free the weight for interpolation */
1325  FREE (cpl_free, id);
1326  FREE (cpl_free, outputData);
1327  FREE (cpl_free, outputErr);
1328 
1329 
1330  gravi_msg_function_exit(1);
1331  return CPL_ERROR_NONE;
1332 }
1333 
1334 
1335 /*----------------------------------------------------------------------------*/
1356 /*----------------------------------------------------------------------------*/
1357 
1358 cpl_error_code gravi_align_spectrum (gravi_data * spectrum_data,
1359  gravi_data * wave_map,
1360  gravi_data * p2vm_map,
1361  enum gravi_detector_type det_type,
1362  const cpl_parameterlist * parlist)
1363 {
1364  gravi_msg_function_start(1);
1365 
1366 
1367  /* Check needed data */
1368  cpl_ensure_code (spectrum_data, CPL_ERROR_NULL_INPUT);
1369  cpl_ensure_code (wave_map, CPL_ERROR_NULL_INPUT);
1370  cpl_ensure_code (p2vm_map, CPL_ERROR_NULL_INPUT);
1371 
1372  /* Duplicate extensions in the spectrum data from p2vm_map
1373  * FIXME: shall copy only required from the type FT/SC */
1374  gravi_data_copy_ext (spectrum_data, p2vm_map, GRAVI_OI_WAVELENGTH_EXT);
1375  CPLCHECK_MSG ("Cannot copy OI_WAVELENGTH extention(s)");
1376 
1377  /* Loop on FT/SC */
1378  int init_type_data = 1;
1379  int end_type_data = 1;
1380  if(det_type == GRAVI_DET_SC || det_type == GRAVI_DET_ALL) {
1381  init_type_data = 0;
1382  }
1383  if(det_type == GRAVI_DET_FT || det_type == GRAVI_DET_ALL) {
1384  end_type_data = 1;
1385  }
1386  for (int type_data = init_type_data; type_data <= end_type_data; type_data++ ) {
1387 
1388  /* Check if SPECTRUM data exists */
1389  if (!gravi_data_has_spectrum (spectrum_data, type_data)) {
1390  cpl_msg_info (cpl_func,"No data for %s, skip it", GRAVI_TYPE(type_data));
1391  continue;
1392  }
1393 
1394  /* Get a timer */
1395  cpl_msg_info (cpl_func, "Re-interpolate the %s",GRAVI_TYPE(type_data));
1396 
1397  /* Get all the input tables */
1398  cpl_table * spectrum_table, * wave_table, * detector_table;
1399  spectrum_table = gravi_data_get_spectrum_data (spectrum_data, type_data);
1400  wave_table = gravi_data_get_wave_data (wave_map, type_data);
1401  detector_table = gravi_data_get_imaging_detector (spectrum_data, type_data);
1402 
1403  cpl_table * specflat_table;
1404  if (type_data == GRAVI_SC)
1405  specflat_table = gravi_data_get_table (spectrum_data, GRAVI_SPECTRUMFLAT_DATA_SC_EXT);
1406  else
1407  specflat_table = NULL;
1408 
1409  CPLCHECK_MSG ("Cannot load the extention data");
1410 
1411  /* Measure total flux before interpolation */
1412  double full_flux_reg0 = gravi_spectrum_get_flux (spectrum_table);
1413 
1414  /* Get the OI_WAVE table of each polarisation */
1415  int npol = gravi_spectrum_get_npol (wave_table);
1416  cpl_table ** oiwave_tables = gravi_data_get_oiwave_tables (p2vm_map, type_data, npol);
1417 
1418  /* Re-interpolate the spectrum into OI_WAVELENGTH
1419  * Run interpolation -- in-place to save enormous time */
1420  gravi_interpolate_spectrum_table (spectrum_table,
1421  wave_table,
1422  oiwave_tables,
1423  detector_table,
1424  specflat_table,
1425  parlist, type_data);
1426 
1427  FREE (cpl_free, oiwave_tables);
1428  CPLCHECK_MSG ("Cannot interpolate");
1429 
1430  /* Measure total flux after interpolation */
1431  double full_flux_reg1 = gravi_spectrum_get_flux (spectrum_table);
1432 
1433  /* Replace (FULLSTARTX,STARTX,NX) by NWAVE in the table plist */
1434  cpl_propertylist * spectrum_plist = gravi_data_get_spectrum_data_plist (spectrum_data, type_data);
1435  int nwave = gravi_spectrum_get_nwave (spectrum_table);
1436  cpl_propertylist_erase (spectrum_plist, PROFILE_FULLSTARTX);
1437  cpl_propertylist_erase (spectrum_plist, PROFILE_STARTX);
1438  cpl_propertylist_erase (spectrum_plist, PROFILE_NX);
1439  cpl_propertylist_update_int (spectrum_plist, "NWAVE", nwave);
1440  CPLCHECK_MSG ("Cannot replace STARTX by NWAVE");
1441 
1442  /* Check if the interpolation is flux conservative */
1443  cpl_msg_info (cpl_func, "Total flux in REGION1s: %.2e [e], in REGION0s:%.2e [e] (ratio=%.5f)",
1444  full_flux_reg1,full_flux_reg0,full_flux_reg1/full_flux_reg0);
1445 
1446  char qc_name[80];
1447  cpl_propertylist * spectrum_header = gravi_data_get_header (spectrum_data);
1448  sprintf (qc_name, "ESO QC TRANS INTERP %s",GRAVI_TYPE(type_data));
1449  cpl_propertylist_update_double (spectrum_header, qc_name, full_flux_reg1/full_flux_reg0);
1450  cpl_propertylist_set_comment (spectrum_header, qc_name, "[e/e] at interpolation");
1451  CPLCHECK_MSG ("Cannot set QC parameters");
1452 
1453  } /* End loop on SC / FT */
1454 
1455  gravi_msg_function_exit(1);
1456  return CPL_ERROR_NONE;
1457 }
1458 
gravi_data * gravi_data_new(int nb_ext)
Create an empty gravi_data.
Definition: gravi_data.c:102
double gravi_pfits_get_time_sc(const cpl_propertylist *header, cpl_size row)
Time of the middle of the SC exposure row in [us], counted from PRC.ACQ.START.
Definition: gravi_pfits.c:590
int gravi_region_get_pol(cpl_table *imaging_detector, int region)
Return the polarisation id of a region.
Definition: gravi_utils.c:445
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_interpolate_spectrum_table(cpl_table *spectrum_table, cpl_table *wave_table, cpl_table **oiwave_tables, cpl_table *detector_table, cpl_table *specflat_table, const cpl_parameterlist *parlist, int type_data)
Re-interpolate in-place a spectrum table.
int gravi_data_has_extension(gravi_data *raw_calib, const char *ext_name)
Check if data has extension with given EXTNAME.
Definition: gravi_data.c:1443
int gravi_pfits_get_startx(const cpl_propertylist *plist)
find out the name of the propertylist
Definition: gravi_pfits.c:66
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
double gravi_pfits_get_sc_gain(const cpl_propertylist *plist)
SC system gain in [ADU/e].
Definition: gravi_pfits.c:766
double gravi_imagelist_get_flux(const cpl_imagelist *imglist)
Return the total flux in imagelist.
Definition: gravi_utils.c:1117
cpl_table ** gravi_data_get_oiwave_tables(gravi_data *data, int type_data, int npol)
Get pointer to the OI_WAVELENGTH tables of both polarisations.
Definition: gravi_data.c:2253
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
double gravi_pfits_get_ft_gain(const cpl_propertylist *plist)
FT system gain in [ADU/e].
Definition: gravi_pfits.c:727
cpl_error_code gravi_remove_badpixel_sc(cpl_imagelist *imglist_sc, cpl_image *bad_img)
Remove the badpixel of the SC.
int gravi_pixel_is_good(cpl_image *bad_img, int x, int y)
Check if the pixel in the BADPIX map is a good pixel.
Definition: gravi_preproc.c:99
cpl_table * gravi_table_ft_format(cpl_table *table_ft, cpl_table *sky_table_std, cpl_table *sky_table_avg, cpl_table *badft_table, int n_region, double gain)
Extract FT spectrum from PIX column.
double gravi_spectrum_get_flux(const cpl_table *table)
Return the total flux in DATA# regions.
Definition: gravi_utils.c:1042
cpl_error_code gravi_imagelist_unwrap_images(cpl_imagelist *imglist)
Unwrap an imagelist an all its images.
Definition: gravi_cpl.c:1496
cpl_propertylist * gravi_data_get_plist(gravi_data *self, const char *extname)
Get the propertylist from EXTNAME.
Definition: gravi_data.c:1684
cpl_imagelist * gravi_data_get_cube(gravi_data *self, const char *extname)
Return a pointer on an IMAGE extension by its EXTNAME.
Definition: gravi_data.c:1751
cpl_error_code gravi_align_spectrum(gravi_data *spectrum_data, gravi_data *wave_map, gravi_data *p2vm_map, enum gravi_detector_type det_type, const cpl_parameterlist *parlist)
Regrid the regions into a common wavelength (in-place)
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_table * gravi_imglist_sc_collapse(cpl_table *profile_table, cpl_imagelist *raw_imglist, cpl_imagelist *rawVar_imglist, cpl_size startx)
Extract the SC spectrum with profile.
cpl_error_code gravi_subtract_met_dark(gravi_data *preproc_data, gravi_data *dark_map)
Substract metrology dark.
gravi_data * gravi_extract_spectrum(gravi_data *raw_data, gravi_data *profile_map, gravi_data *dark_map, gravi_data *bad_map, gravi_data *sky_map, const cpl_parameterlist *parlist, enum gravi_detector_type det_type)
Create the SPECTRUM gravi_data with extracted spectrum per region.
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