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