GRAVI Pipeline Reference Manual 1.9.3
Loading...
Searching...
No Matches
gravi_calib.c
Go to the documentation of this file.
1/* $Id: gravi_calib.c,v 1.10 2012/03/23 15:10:40 nazouaoui Exp $
2 *
3 * This file is part of the GRAVI Pipeline
4 * Copyright (C) 2002,2003 European Southern Observatory
5 *
6 * This program is free software; you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation; either version 2 of the License, or
9 * (at your option) any later version.
10 *
11 * This program is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program; if not, write to the Free Software
18 * Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
19 */
20
37/*
38 * History :
39 * 11/01/2019 : Move global parameter from gravi_calib.h
40 * initialize sigma with 0
41 * add brackets to 'for' statements
42 */
43/*----------------------------------------------------------------------------
44 DEBUG
45 -----------------------------------------------------------------------------*/
46
47#define INFO_DEBUG 0
48/*-----------------------------------------------------------------------------
49 Includes
50 -----------------------------------------------------------------------------*/
51
52#ifdef HAVE_CONFIG_H
53#include <config.h>
54#endif
55
56#include <cpl.h>
57#include <string.h>
58#include <stdio.h>
59#include <time.h>
60#include <math.h>
61#include <complex.h>
62
63#include "gravi_data.h"
64#include "gravi_dfs.h"
65#include "gravi_pfits.h"
66#include "gravi_cpl.h"
67
68#include "gravi_utils.h"
69
70#include "gravi_preproc.h"
71#include "gravi_calib.h"
72
73#include "gravi_pca.h"
74
75/*-----------------------------------------------------------------------------
76 ekw 11/01/2019 : Global parameter from gravi_calib.h
77 -----------------------------------------------------------------------------*/
78
79static double met_Sep_2016[64] = {-0.000338233 ,-0.000282107 ,0.000345582 ,0.000381979 ,0.000326764 ,-0.000225664 ,5.35555e-05 ,0.00060938 ,-0.000205619 ,-0.000265577 ,0.00017356 ,0.000265687 ,0.000295987 ,3.05255e-05 ,-0.000498407 ,0.000232075 ,0.00502058 ,0.002045925 ,-0.000110657 ,-0.000403592 ,-0.000065043 ,-0.000433645 ,-0.000626545 ,3.29765e-05 ,-0.001694525 ,-0.0011684 ,6.2795e-06 ,-0.001584115 ,-0.000735553 ,-0.000868538 ,-0.000985087 ,-0.0008204 ,0.000768577 ,0.000848342 ,-0.000134943 ,-0.000385157 ,-8.68822e-05 ,-0.000757366 ,0.000446051 ,-0.000231723 ,0.000790425 ,-0.000638897 ,0.000503496 ,-7.78205e-05 ,0.000287366 ,0.000243789 ,0.000083288 ,-0.000125138 ,-0.000147337 ,4.14224e-05 ,0.000123082 ,-0.00117179 ,4.54785e-05 ,-0.000186707 ,0.000682836 ,0.00090649 ,0.000357256 ,-0.002133845 ,-0.00151895 ,-0.00150048 ,-0.00266423 ,-0.0030716 ,0.000599228 ,-0.001078583};
80static double met_Mar_2017[64] = {0.000481474 ,0.000678173 ,0.00022232 ,9.04735e-05 ,0.001030283 ,0.000400478 ,0.000199579 ,0.000610902 ,0.000870061 ,0.001148845 ,0.000706636 ,0.000491999 ,0.00093158 ,0.001224225 ,0.000652115 ,0.001117025 ,0.00481452 ,0.002359475 ,0.000641491 ,0.000126103 ,0.000158341 ,0.001455785 ,0.000227113 ,0.000366087 ,-0.000372672 ,-0.000814455 ,0.000403834 ,-0.00072791 ,-0.000422227 ,0.000413887 ,-0.000024651 ,0.000106683 ,0.00092596 ,0.000327427 ,0.000775269 ,0.000906505 ,0.000108337 ,-0.000214467 ,0.001249965 ,0.000694693 ,0.000718101 ,0.00083926 ,0.00138818 ,0.00131215 ,0.001113065 ,0.00134113 ,0.000972572 ,0.00073247 ,-3.35943e-05 ,0.000312545 ,0.000365923 ,0.000510906 ,0.000955084 ,0.000808904 ,0.000403238 ,-0.000150186 ,0.000200673 ,-0.00078668 ,-0.00078133 ,-0.00039955 ,-0.00226333 ,-0.00177287 ,0.000435995 ,-0.000397403};
81static double met_Jun_2017[64] = {0.00029923 ,5.42185e-05 ,-0.000023823 ,0.000187444 ,0.000650675 ,0.000135629 ,-0.000047364 ,0.000479489 ,0.000463565 ,0.000536066 ,0.000528847 ,0.000874895 ,0.000702853 ,0.000438988 ,0.000183642 ,0.000536044 ,0.004829985 ,0.00215068 ,5.94134e-05 ,-0.000164468 ,-4.89301e-05 ,0.000400028 ,-0.000277333 ,-0.000169422 ,-0.00136967 ,-0.000998661 ,0.000244959 ,-0.00123726 ,-0.000404182 ,-0.000511517 ,-0.000398515 ,-0.000314167 ,0.000546602 ,0.000547136 ,0.000388259 ,0.000108443 ,-0.000266431 ,-0.000734324 ,0.000946866 ,0.0001005 ,0.000739434 ,0.000530113 ,0.000821013 ,0.000748435 ,0.00094666 ,0.00102861 ,-0.00013052 ,0.000576223 ,-0.000657692 ,-0.000293389 ,0.000146337 ,-0.000338176 ,0.000241386 ,0.00005168 ,0.000495891 ,-1.76273e-05 ,8.84135e-05 ,-0.00153857 ,-0.000850848 ,-0.000986058 ,-0.002700175 ,-0.002620785 ,0.000440797 ,-0.00063843};
82static double met_Jul_2017[64] = {8.19405e-05 ,-0.000061716 ,0.000067827 ,0.000234624 ,0.000521161 ,0.000031451 ,0.000211139 ,0.000578157 ,0.000316435 ,0.000431433 ,0.000431228 ,0.000761259 ,0.000576742 ,0.000299923 ,0.00034315 ,0.000561323 ,0.004834095 ,0.00208259 ,-0.000155192 ,-0.000120373 ,-0.000328756 ,0.000265791 ,-0.000249398 ,-0.000310416 ,-0.001047503 ,-0.000884526 ,0.000431645 ,-0.000890023 ,-0.000330709 ,-0.000540866 ,-0.0002067 ,-0.000334498 ,0.00035945 ,0.000683828 ,0.000464533 ,0.000166613 ,-0.000276387 ,-0.00088093 ,0.001283446 ,0.00016521 ,0.000578248 ,0.000334565 ,0.000762716 ,0.000521964 ,0.000757276 ,0.000959739 ,0.000035001 ,0.000442429 ,-0.0006928 ,-0.000105863 ,-0.00006137 ,-0.000454825 ,-0.000157343 ,-5.68635e-05 ,0.000296553 ,-3.80711e-05 ,0.000384616 ,-0.00149983 ,-0.000650183 ,-0.00076473 ,-0.002494945 ,-0.00246317 ,0.000794264 ,-0.000337915};
83static double met_Aug_2017[64] = {0.000225917 ,0.00033629 ,1.70956e-05 ,0.000178987 ,0.00103744 ,0.000275415 ,0.000107582 ,0.000678724 ,0.000402821 ,0.000713535 ,0.000518528 ,0.000607008 ,0.000506748 ,0.000640514 ,0.000300349 ,0.000522693 ,0.005585285 ,0.00305222 ,0.000771553 ,0.001006485 ,0.000618719 ,0.00130769 ,0.000419028 ,0.000491215 ,-0.000815277 ,-0.000646205 ,0.000640702 ,-0.000881775 ,-0.000200994 ,-0.000258582 ,-0.000465663 ,-0.00004838 ,0.000591627 ,0.000615608 ,0.000197913 ,0.000352664 ,-0.000113565 ,-0.000637229 ,0.001241975 ,-6.659e-06 ,0.000649339 ,0.000641595 ,0.000673658 ,0.000657426 ,0.001263897 ,0.001258565 ,0.000050901 ,0.000318844 ,-9.23866e-05 ,0.00065074 ,0.000793493 ,0.000600095 ,0.001049984 ,0.000858337 ,0.000786737 ,0.000737051 ,0.000576774 ,-0.001384295 ,-0.000642989 ,-0.0008569 ,-0.00208747 ,-0.00247888 ,0.00108834 ,-0.000190064};
84
85/*-----------------------------------------------------------------------------
86 Private prototypes
87 -----------------------------------------------------------------------------*/
88
89cpl_error_code gravi_fit_profile (cpl_vector * values_x0,
90 cpl_vector * values_y0,
91 cpl_vector * values_sigma,
92 cpl_image * mean_img,
93 int ref_x, int ref_y,
94 int size_profile,
95 const char * resolution);
96
97cpl_image * gravi_create_profile_image (cpl_image * mean_img,
98 cpl_vector * values_x0,
99 cpl_vector * values_y0,
100 cpl_vector * values_sigma,
101 cpl_size iy_min, cpl_size iy_max,
102 const char * resolution);
103
104/*-----------------------------------------------------------------------------
105 Functions code
106 -----------------------------------------------------------------------------*/
107
108/*----------------------------------------------------------------------------*/
123/*----------------------------------------------------------------------------*/
124
126{
128 cpl_ensure (raw_data, CPL_ERROR_NULL_INPUT, NULL);
129
130 /* Create the output DARK or SKY map */
131 gravi_data * dark_map = gravi_data_new (0);
132
133 /* Dump full header of RAW data */
134 cpl_propertylist * dark_header = gravi_data_get_header (dark_map);
135 cpl_propertylist * raw_header = gravi_data_get_header (raw_data);
136 cpl_propertylist_append (dark_header, raw_header);
137
138 /* Check if this is a SKY or a DARK */
139 const char * dpr_type = gravi_pfits_get_dpr_type (raw_header);
140 int isSky = strstr (dpr_type, "SKY")?1:0;
141
142 /* The dark file must contains all the shutter close */
143 if ( isSky==0 && !gravi_data_check_shutter_closed (raw_data) ) {
144 gravi_pfits_add_check (dark_header, "DARK has some shutter OPEN !!");
145 }
146
147 /* The sky file must contains all the shutter open */
148 if ( isSky==1 && !gravi_data_check_shutter_open (raw_data) ) {
149 gravi_pfits_add_check (dark_header, "SKY has some shutter CLOSED !!");
150 }
151
152 /*
153 * Compute the SC DARK
154 */
156 cpl_msg_warning (cpl_func,"The DARK data has no IMAGING_DATA_SC");
157 }
158 else
159 {
160 cpl_msg_info (cpl_func, "Computing the %s of SC",isSky?"SKY":"DARK");
161
162 /* Copy IMAGING_DETECTOR into product */
164
165 /* Load the IMAGING_DATA table or image list */
166 cpl_imagelist * imglist = gravi_data_get_cube (raw_data, GRAVI_IMAGING_DATA_SC_EXT);
167
168 /* Compute the median image of the imagelist */
169 cpl_image * median_img = cpl_imagelist_collapse_create (imglist);
170 CPLCHECK_NUL ("Cannot compute the median dark");
171
172 /* Compute STD. We should see if we use sigma-clipping
173 * or not for the collapse, it may change the bad pixel detection */
174 cpl_msg_info (cpl_func,"Compute std with imglist");
175 cpl_imagelist * temp_imglist = cpl_imagelist_duplicate (imglist);
176 cpl_imagelist_subtract_image (temp_imglist, median_img);
177 cpl_imagelist_power (temp_imglist, 2.0);
178 cpl_image * stdev_img = cpl_imagelist_collapse_create (temp_imglist);
179 FREE (cpl_imagelist_delete, temp_imglist);
180 cpl_image_power (stdev_img, 0.5);
181 CPLCHECK_NUL ("Cannot compute the STD of the DARK");
182
183 /* Compute the QC parameters RMS, MEDIAN, ZERO.NB */
184 cpl_msg_info (cpl_func, "Compute QC parameters");
185 double mean_qc = cpl_image_get_median (median_img);
186 double darkrms = cpl_image_get_median (stdev_img);
187 cpl_propertylist_update_double (dark_header, isSky?QC_MEANSKY_SC:QC_MEANDARK_SC, mean_qc);
188 cpl_propertylist_update_double (dark_header, isSky?QC_SKYRMS_SC:QC_DARKRMS_SC, darkrms);
190 {
191 cpl_imagelist * acq_imglist = gravi_data_get_cube (raw_data, GRAVI_IMAGING_DATA_ACQ_EXT);
192
193 size_t acq_dark_zero_count = 0;
194 for(int i=0; i < cpl_imagelist_get_size(acq_imglist); i++)
195 {
196 cpl_mask * acq_zero_mask = cpl_mask_threshold_image_create(cpl_imagelist_get(acq_imglist, i), -FLT_MIN, FLT_MIN);
197 acq_dark_zero_count+= cpl_mask_count(acq_zero_mask);
198 cpl_mask_delete(acq_zero_mask);
199 }
200 cpl_propertylist_update_double (dark_header, QC_ACQ_ZERO_NB, acq_dark_zero_count);
201 }
202
203 /* Verbose */
204 cpl_msg_info (cpl_func, "QC_MEDIAN%s_SC = %e",isSky?"SKY":"DARK", mean_qc);
205 cpl_msg_info (cpl_func, "QC_%sRMS_SC = %e",isSky?"SKY":"DARK", darkrms);
206
207 /* Put the data in the output table : dark_map */
208 cpl_propertylist * img_plist = gravi_data_get_plist (raw_data, GRAVI_IMAGING_DATA_SC_EXT);
209 img_plist = cpl_propertylist_duplicate (img_plist);
210 gravi_data_add_img (dark_map, img_plist, GRAVI_IMAGING_DATA_SC_EXT, median_img);
211
212 img_plist = cpl_propertylist_duplicate (img_plist);
213 gravi_data_add_img (dark_map, img_plist, GRAVI_IMAGING_ERR_SC_EXT, stdev_img);
214 CPLCHECK_NUL ("Cannot set the SC data");
215
216 } /* End SC case */
217
218 /*
219 * Compute the FT DARK
220 */
222 cpl_msg_warning (cpl_func,"The DARK data has no IMAGING_DATA_FT");
223 }
224 else
225 {
226 cpl_msg_info(cpl_func, "Computing the %s of FT",isSky?"SKY":"DARK");
227
228 /* Copy IMAGING_DETECTOR into product */
230
231 /* Load the IMAGING_DATA table as DOUBLE */
232 cpl_msg_info (cpl_func,"Load data");
233 cpl_table * table_ft = gravi_data_get_table (raw_data, GRAVI_IMAGING_DATA_FT_EXT);
234 cpl_table_cast_column (table_ft, "PIX", "PIX", CPL_TYPE_DOUBLE);
235 cpl_imagelist * imglist = gravi_imagelist_wrap_column (table_ft, "PIX");
236 CPLCHECK_NUL ("Cannot load the FT data");
237
238 /* Compute the median image of the imagelist */
239 cpl_msg_info (cpl_func,"Compute mean and median");
240 cpl_image * median_img = cpl_imagelist_collapse_median_create (imglist);
241 cpl_image * mean_img = cpl_imagelist_collapse_create (imglist);
242 CPLCHECK_NUL ("Cannot compute the MEAN dark");
243
244 /* Compute the std of each pixels */
245 cpl_msg_info (cpl_func,"Compute std");
246 cpl_imagelist_subtract_image (imglist, mean_img);
247 cpl_imagelist_power (imglist, 2.0);
248 cpl_image * stdev_img = cpl_imagelist_collapse_create (imglist);
249 cpl_image_power (stdev_img, 0.5);
250
251 /* We power back the data, because we are working in-place */
252 cpl_imagelist_power (imglist, 0.5);
253 cpl_imagelist_add_image (imglist, mean_img);
254
255 /* Delete data */
256 cpl_msg_info (cpl_func,"Delete data");
258 FREE (cpl_image_delete, mean_img);
259 CPLCHECK_NUL ("Cannot compute the STD of the DARK");
260
261 /* Compute the QC parameters RMS and MEDIAN */
262 cpl_msg_info (cpl_func, "Compute QC parameters");
263 double mean_qc = cpl_image_get_mean (median_img);
264 double darkrms = cpl_image_get_mean (stdev_img);
265 cpl_propertylist_update_double (dark_header, isSky?QC_MEANSKY_FT:QC_MEANDARK_FT, mean_qc);
266 cpl_propertylist_update_double (dark_header, isSky?QC_SKYRMS_FT:QC_DARKRMS_FT, darkrms);
267 CPLCHECK_NUL ("Cannot compute the QC");
268
269 /* Verbose */
270 cpl_msg_info (cpl_func, "QC_MEDIAN%s_FT = %e",isSky?"SKY":"DARK", mean_qc);
271 cpl_msg_info (cpl_func, "QC_%sRMS_FT = %e",isSky?"SKY":"DARK", darkrms);
272
273 /* Create the output DARK table, with a single row */
274 cpl_table * median_table = cpl_table_extract (table_ft, 0, 1);
275 cpl_array * median_array = gravi_array_wrap_image (median_img);
276 cpl_table_set_array (median_table, "PIX", 0, median_array);
277 FREE (cpl_array_unwrap, median_array);
278 FREE (cpl_image_delete, median_img);
279 CPLCHECK_NUL("Cannot set median in table");
280
281 /* Put median dark in the output gravi_data */
282 gravi_data_add_table (dark_map, NULL, GRAVI_IMAGING_DATA_FT_EXT, median_table);
283 CPLCHECK_NUL("Cannot save median in gravi_data");
284
285 /* Create the output DARK RMS table, with a single row */
286 cpl_table * stdev_table = cpl_table_extract (table_ft, 0, 1);
287 cpl_array * stdev_array = gravi_array_wrap_image (stdev_img);
288 cpl_table_set_array (stdev_table, "PIX", 0, stdev_array);
289 FREE (cpl_array_unwrap, stdev_array);
290 FREE (cpl_image_delete, stdev_img);
291 CPLCHECK_NUL("Cannot set rms in table");
292
293 /* Put median dark in the output gravi_data */
294 gravi_data_add_table (dark_map, NULL, GRAVI_IMAGING_ERR_FT_EXT, stdev_table);
295 CPLCHECK_NUL("Cannot save median in gravi_data");
296
297 } /* End case FT */
298
299 /*
300 * Compute the METROLOGY DARK
301 */
302 if (( isSky==1 )||(!gravi_data_has_extension (raw_data, GRAVI_METROLOGY_EXT))) {
303
304 if ( isSky==0 )
305 cpl_msg_warning (cpl_func,"The DARK data has no METROLOGY");
306
307 }
308 else
309 {
310
311 cpl_msg_info(cpl_func, "Computing the %s of METROLOGY",isSky?"SKY":"DARK");
312
313 /* Load the IMAGING_DATA table as DOUBLE */
314 cpl_msg_info (cpl_func,"Load data");
315 cpl_table * table_met = gravi_data_get_table (raw_data, GRAVI_METROLOGY_EXT);
316 cpl_table_cast_column (table_met, "VOLT", "VOLT", CPL_TYPE_DOUBLE);
317 cpl_imagelist * imglist = gravi_imagelist_wrap_column (table_met, "VOLT");
318 CPLCHECK_NUL ("Cannot load the VOLT data");
319
320 /* Compute the median image of the imagelist */
321 cpl_msg_info (cpl_func,"Compute mean and median");
322 cpl_image * median_img = cpl_imagelist_collapse_median_create (imglist);
323 cpl_image * mean_img = cpl_imagelist_collapse_create (imglist);
324 CPLCHECK_NUL ("Cannot compute the MEAN dark");
325
326 /* Compute the std of each pixels */
327 cpl_msg_info (cpl_func,"Compute std");
328 cpl_imagelist_subtract_image (imglist, mean_img);
329 cpl_imagelist_power (imglist, 2.0);
330 cpl_image * stdev_img = cpl_imagelist_collapse_create (imglist);
331 cpl_image_power (stdev_img, 0.5);
332
333 /* We power back the data, because we are working in-place */
334 cpl_imagelist_power (imglist, 0.5);
335 cpl_imagelist_add_image (imglist, mean_img);
336
337 /* Delete data */
338 cpl_msg_info (cpl_func,"Delete data");
340 FREE (cpl_image_delete, mean_img);
341 CPLCHECK_NUL ("Cannot compute the STD of the DARK");
342
343 /* Compute the QC parameters RMS and MEDIAN */
344 cpl_msg_info (cpl_func, "Compute QC parameters");
345 double mean_qc = cpl_image_get_mean_window (median_img,1,1,64,1);
346 double darkrms = cpl_image_get_mean_window (stdev_img,1,1,64,1);
347 double darkmin = cpl_image_get_min_window (median_img,1,1,64,1);
348 double darkmax = cpl_image_get_max_window (median_img,1,1,64,1);
349 if (darkmin>0) darkmin=0;
350 if (darkmax<0) darkmax=0;
351 cpl_propertylist_update_double (dark_header, QC_MEANDARK_MET, mean_qc);
352 cpl_propertylist_update_double (dark_header, QC_DARKRMS_MET, darkrms);
353 cpl_propertylist_update_double (dark_header, QC_DARKRANGE_MET, darkmax-darkmin);
354 CPLCHECK_NUL ("Cannot compute the QC");
355
356 /* Verbose */
357 cpl_msg_info (cpl_func, "QC_MEDIAN%s_MET = %e","DARK", mean_qc);
358 cpl_msg_info (cpl_func, "QC_%sRMS_MET = %e","DARK", darkrms);
359 cpl_msg_info (cpl_func, "QC_%sRANGE_MET = %e","DARK", darkmax-darkmin);
360
361 /* Create the output DARK table, with a single row */
362 cpl_table * median_table = cpl_table_extract (table_met, 0, 1);
363 cpl_array * median_array = gravi_array_wrap_image (median_img);
364
365 /* Put to zero the dark of the fiber coupler unit */
366 for (cpl_size diode = 64; diode < 80; diode++){
367 cpl_array_set(median_array, diode, 0);
368 }
369
370 /* Check on time */
371 double time_mjd_obs= cpl_propertylist_get_double (raw_header, "MJD-OBS");
372
373 /* If the data is too old, we are using an old dataset */
374 for (cpl_size diode = 0; diode < 64; diode++)
375 {
376 if (time_mjd_obs < 57747) cpl_array_set(median_array, diode, met_Sep_2016[diode]); /* 25 decembre 2016 */
377 else if (time_mjd_obs < 57851) cpl_array_set(median_array, diode, met_Mar_2017[diode]); /* 8 April */
378 else if (time_mjd_obs < 57924) cpl_array_set(median_array, diode, met_Jun_2017[diode]); /* 20 Juin */
379 else if (time_mjd_obs < 57990) cpl_array_set(median_array, diode, met_Jul_2017[diode]); /* 25 Juillet */
380 else if (time_mjd_obs < 58208.01) cpl_array_set(median_array, diode, met_Aug_2017[diode]);
381 }
382
383 cpl_table_set_array (median_table, "VOLT", 0, median_array);
384 FREE (cpl_array_unwrap, median_array);
385 FREE (cpl_image_delete, median_img);
386 CPLCHECK_NUL("Cannot set median in table");
387
388 /* Put median dark in the output gravi_data */
389 gravi_data_add_table (dark_map, NULL, GRAVI_METROLOGY_EXT, median_table);
390 CPLCHECK_NUL("Cannot save median in gravi_data");
391
392 /* Create the output DARK RMS table, with a single row */
393 cpl_table * stdev_table = cpl_table_extract (table_met, 0, 1);
394 cpl_array * stdev_array = gravi_array_wrap_image (stdev_img);
395 cpl_table_set_array (stdev_table, "VOLT", 0, stdev_array);
396 FREE (cpl_array_unwrap, stdev_array);
397 FREE (cpl_image_delete, stdev_img);
398 CPLCHECK_NUL("Cannot set rms in table");
399
400 /* Put median dark in the output gravi_data */
401 gravi_data_add_table (dark_map, NULL, GRAVI_METROLOGY_ERR_EXT, stdev_table);
402 CPLCHECK_NUL("Cannot save median in gravi_data");
403
404 } /* End case METROLOGY */
405
406 /*
407 * Compute the ACQ DARK
408 */
410 cpl_msg_warning (cpl_func,"The DARK data has no IMAGING_DATA_ACQ");
411 }
412 else
413 {
414 cpl_msg_info (cpl_func, "Computing the %s of ACQ",isSky?"SKY":"DARK");
415
416 /* Load the IMAGING_DATA table or image list */
417 cpl_imagelist * imglist = gravi_data_get_cube (raw_data, GRAVI_IMAGING_DATA_ACQ_EXT);
418
419 /* Compute the median image of the imagelist */
420 cpl_image * median_img = cpl_imagelist_collapse_median_create (imglist);
421 CPLCHECK_NUL ("Cannot compute the median dark");
422
423 /* Compute the QC parameters RMS and MEDIAN */
424 cpl_msg_info (cpl_func, "Compute QC parameters");
425 double mean_qc = cpl_image_get_median (median_img);
426
427 /* Verbose */
428 cpl_msg_info (cpl_func, "QC_MEDIAN%s_ACQ = %e",isSky?"SKY":"DARK", mean_qc);
429 cpl_propertylist_update_double (dark_header, isSky?"ESO QC MEDIANSKY ACQ":"ESO QC MEDIANDARK ACQ",
430 mean_qc);
431
432 /* Put the data in the output table : dark_map */
433 cpl_propertylist * img_plist = gravi_data_get_plist (raw_data, GRAVI_IMAGING_DATA_ACQ_EXT);
434 img_plist = cpl_propertylist_duplicate (img_plist);
435 gravi_data_add_img (dark_map, img_plist, GRAVI_IMAGING_DATA_ACQ_EXT, median_img);
436 CPLCHECK_NUL("Cannot save median in gravi_data");
437 }
438
440 return dark_map;
441}
442
443/*---------------------------------------------------------------------------*/
455/*---------------------------------------------------------------------------*/
456
457gravi_data * gravi_average_dark (gravi_data ** data, cpl_size ndata)
458{
460 cpl_ensure (data, CPL_ERROR_NULL_INPUT, NULL);
461
462 gravi_msg_warning ("FIXME","Averaging DARK or SKY is experimental");
463
464 gravi_data * output_data = gravi_data_duplicate (data[0]);
465
466 /* Average the IMAGING_DATA of SC */
467 cpl_msg_info (cpl_func,"Average IMAGING_DATA of SC");
468
469 cpl_image * darksc_image = gravi_data_get_img (output_data, GRAVI_IMAGING_DATA_SC_EXT);
470 for (int file = 1; file < ndata; file++) {
471 cpl_image_add (darksc_image, gravi_data_get_img (data[file], GRAVI_IMAGING_DATA_SC_EXT));
472 }
473 cpl_image_divide_scalar (darksc_image, ndata);
474 CPLCHECK_NUL ("Cannot average DARK of SC");
475
476 /* Average the IMAGING_ERR of SC */
477 cpl_msg_info (cpl_func,"Average IMAGING_ERR of SC");
478
479 cpl_image * stdevsc_image = gravi_data_get_img (output_data, GRAVI_IMAGING_ERR_SC_EXT);
480 for (int file = 1; file < ndata; file++) {
481 cpl_image_add (stdevsc_image, gravi_data_get_img (data[file], GRAVI_IMAGING_ERR_SC_EXT));
482 }
483 cpl_image_divide_scalar (stdevsc_image, ndata);
484 CPLCHECK_NUL ("Cannot average DARKERR of SC");
485
486 /* Average the IMAGING_DATA of FT */
487 cpl_msg_info (cpl_func,"Average IMAGING_DATA of FT");
488
489 cpl_array * darkft_array;
490 darkft_array = cpl_table_get_data_array (gravi_data_get_table (output_data, GRAVI_IMAGING_DATA_FT_EXT), "PIX")[0];
491 for (int file = 1; file < ndata; file++) {
492 cpl_array_add (darkft_array, cpl_table_get_data_array (gravi_data_get_table (data[file], GRAVI_IMAGING_DATA_FT_EXT), "PIX")[0]);
493 }
494 cpl_array_divide_scalar (darkft_array, ndata);
495 CPLCHECK_NUL ("Cannot average DARK of FT");
496
497 /* Average the IMAGING_ERR of FT */
498 cpl_msg_info (cpl_func,"Average IMAGING_ERR of FT");
499
500 cpl_array * errft_array;
501 errft_array = cpl_table_get_data_array (gravi_data_get_table (output_data, GRAVI_IMAGING_ERR_FT_EXT), "PIX")[0];
502 for (int file = 1; file < ndata; file++) {
503 cpl_array_add (errft_array, cpl_table_get_data_array (gravi_data_get_table (data[file], GRAVI_IMAGING_ERR_FT_EXT), "PIX")[0]);
504 }
505 cpl_array_divide_scalar (errft_array, ndata);
506 CPLCHECK_NUL ("Cannot average DARK of FT");
507
509 return output_data;
510}
511
512
513/*----------------------------------------------------------------------------*/
534/*----------------------------------------------------------------------------*/
535
536cpl_error_code gravi_fit_profile (cpl_vector * values_x0,
537 cpl_vector * values_y0,
538 cpl_vector * values_sigma,
539 cpl_image * mean_img,
540 int ref_x, int ref_y,
541 int size_profile,
542 const char * resolution)
543{
545 cpl_ensure_code (values_x0, CPL_ERROR_NULL_INPUT);
546 cpl_ensure_code (values_y0, CPL_ERROR_NULL_INPUT);
547 cpl_ensure_code (values_sigma, CPL_ERROR_NULL_INPUT);
548 cpl_ensure_code (mean_img, CPL_ERROR_NULL_INPUT);
549 cpl_ensure_code (resolution, CPL_ERROR_NULL_INPUT);
550
551 /* Compute middle of size_profile */
552 int middle = floor (size_profile / 2);
553
554 /* Get the image dimensions */
555 int nx = cpl_image_get_size_x (mean_img);
556 int ny = cpl_image_get_size_y (mean_img);
557
558 cpl_ensure_code (ref_x > 0 && ref_x < nx, CPL_ERROR_ILLEGAL_INPUT);
559 cpl_ensure_code (ref_y > 0 && ref_y < ny, CPL_ERROR_ILLEGAL_INPUT);
560
561 /* y0: best fit profile position (y0_ is of previous channel)
562 * sigma: best fit profile width (sigma_ is of previous channel) */
563 double y0, y0_ = ref_y;
564 double sigma, sigma_ = 0;
565
566 /* Loop on spectral channels. Here fit all the
567 * points under the reference point */
568 int coord_x = ref_x;
569 while (coord_x <= nx) {
570
571 /* Extract the corresponding column in the 2D image */
572 cpl_vector * line_data;
573 line_data = cpl_vector_new_from_image_column (mean_img,
574 coord_x);
575
576 /* Construction of the vector of references points y_data and
577 * its values z_data for the gaussian fit */
578 cpl_vector * y_data, * z_data;
579 y_data = cpl_vector_new (size_profile);
580
581 if (! (strcmp(resolution, "LOW") && strcmp(resolution, "MED")) ) {
582 /* Case LOW and MED */
583 for (cpl_size i = 0; i < size_profile; i++){
584 cpl_vector_set (y_data, i, ref_y - size_profile/2 + i);
585 }
586 z_data = cpl_vector_extract (line_data,
587 ref_y - size_profile/2,
588 ref_y - size_profile/2+size_profile-1, 1);
589 }
590 else {
591 /* Case HIGH */
592 if ((ref_y - middle) <= 0 ) {
593 for (cpl_size i = 0; i < size_profile; i++){
594 cpl_vector_set (y_data, i, i);
595 }
596 z_data = cpl_vector_extract (line_data, 0,
597 size_profile - 1, 1);
598 }
599 else if ((ref_y - middle + size_profile) >= ny ){
600 for (cpl_size i = 0; i < size_profile; i++){
601 cpl_vector_set (y_data, i, ny - size_profile + i);
602 }
603 z_data = cpl_vector_extract (line_data,
604 ny - size_profile,
605 ny - 1, 1);
606 }
607 else {
608 for (cpl_size i = 0; i < size_profile; i++){
609 cpl_vector_set(y_data, i, ref_y + (i - middle));
610 }
611 z_data = cpl_vector_extract (line_data, ref_y - middle,
612 ref_y - middle + size_profile - 1, 1);
613 }
614 }
615
616 /* Fit a Gaussian to the extracted data */
617 cpl_errorstate prestate = cpl_errorstate_get();
618 double area, mse, offset = 0;
619 cpl_vector_fit_gaussian (y_data, NULL, z_data, NULL,
620 CPL_FIT_ALL, &y0, &sigma, &area,
621 &offset, &mse, NULL, NULL);
622
623 /* If the fit fail, we keep the value of the
624 * previous spectral channel */
625 if (cpl_error_get_code() == CPL_ERROR_CONTINUE){
626 cpl_errorstate_set (prestate);
627 if (coord_x == ref_x) {
628 y0 = ref_y;
629 sigma = 1;
630 }else{
631 y0 = y0_;
632 sigma = sigma_;
633 }
634 cpl_msg_warning (cpl_func, "Cannot fit profile of channel %d, new y0=%e", coord_x, y0);
635 }
636
637 CPLCHECK_MSG ("Error during the gaussian fit");
638
639 /* If more than 1 pixel shift, we keep the value of the
640 * previous spectral channel */
641 if ((fabs(y0 - y0_) >= 1) && (coord_x != ref_x)) {
642 cpl_msg_warning (cpl_func, "Too much difference of channel %d with previous (%e pixels)", coord_x, y0 - y0_);
643 y0 = y0_;
644 sigma = sigma_;
645 }
646
647 /* In HIGH Resolution, we update ref_y in order to
648 * follow the curvature when extracting the data */
649 if (! (strcmp(resolution, "HIGH"))) {
650 ref_y = floor(y0);
651 }
652
653 /* Record x0 and sigma to compare to next channel */
654 y0_ = y0;
655 sigma_ = sigma;
656
657 /* Save best fit y0 and sigma into output vectors */
658 cpl_vector_set (values_x0, (coord_x - 1), coord_x - 1);
659 cpl_vector_set (values_y0, (coord_x - 1), y0);
660 cpl_vector_set (values_sigma, (coord_x - 1), sigma);
661
662 /* Increment spectral channel number */
663 coord_x ++;
664
665 FREE (cpl_vector_delete, line_data);
666 FREE (cpl_vector_delete, z_data);
667 FREE (cpl_vector_delete, y_data);
668
669 } /* End fit the column ref_x -> nx */
670
671
672 /* Reset coord_x to the middle of the spectra */
673 coord_x = ref_x - 1;
674
675 /* In HIGH, we use the best-fit position for this channel */
676 if (! (strcmp(resolution, "HIGH"))) {
677 ref_y = cpl_vector_get (values_y0, (ref_x - 1));
678 }
679
680
681 /* Loop on spectral channels. Here fit all the
682 * points under the reference point */
683 while (coord_x > 0) {
684
685 /* Extract the corresponding column in the 2D image */
686 cpl_vector * line_data;
687 line_data = cpl_vector_new_from_image_column (mean_img,
688 coord_x);
689
690 /* Construction of the vector of references points y_data and
691 * its values z_data for the gaussian fit */
692 cpl_vector * y_data, * z_data;
693 y_data = cpl_vector_new(size_profile);
694
695 if (! (strcmp(resolution, "LOW") && strcmp(resolution, "MED")) ) {
696 /* Case LOW and MED */
697 for (cpl_size i = 0; i < size_profile; i++){
698 cpl_vector_set(y_data, i, ref_y - size_profile/2 + i);
699 }
700 z_data = cpl_vector_extract (line_data,
701 ref_y - size_profile/2,
702 ref_y - size_profile/2 +size_profile-1, 1);
703 }
704 else {
705 /* Case HIGH */
706 if ((ref_y - middle) <= 0 ){
707 for (cpl_size i = 0; i < size_profile; i++){
708 cpl_vector_set(y_data, i, i);
709 }
710 z_data = cpl_vector_extract (line_data, 0,
711 size_profile - 1, 1);
712 }
713 else if ((ref_y - middle + size_profile) >= ny ){
714 for (cpl_size i = 0; i < size_profile; i++){
715 cpl_vector_set(y_data, i, ny - size_profile + i);
716 }
717 z_data = cpl_vector_extract (line_data,
718 ny - size_profile,
719 ny - 1, 1);
720 }
721 else{
722 for (cpl_size i = 0; i < size_profile; i++){
723 cpl_vector_set(y_data, i, ref_y + (i - middle));
724 }
725 z_data = cpl_vector_extract (line_data, ref_y - middle,
726 ref_y - middle + size_profile - 1, 1);
727 }
728 }
729
730 /* Fit a Gaussian to the extracted data */
731 cpl_errorstate prestate = cpl_errorstate_get();
732 double area, mse, offset = 0;
733 cpl_vector_fit_gaussian (y_data, NULL, z_data, NULL,
734 CPL_FIT_ALL, &y0, &sigma, &area,
735 &offset, &mse, NULL, NULL);
736
737 /* If the fit fail, we keep the value of the
738 * previous spectral channel */
739 if (cpl_error_get_code() == CPL_ERROR_CONTINUE){
740 cpl_errorstate_set (prestate);
741 if (coord_x == ref_x-1){
742 y0 = ref_y;
743 sigma = 1;
744 } else {
745 y0 = y0_;
746 sigma = sigma_;
747 }
748 cpl_msg_warning (cpl_func, "Cannot fit profile of channel %d, new y0=%e", coord_x, y0);
749 }
750
751 CPLCHECK_MSG ("Error during the gaussian fit");
752
753 /* If more than 1 pixel shift, we keep the value of the
754 * previous spectral channel */
755 if ((fabs(y0 - y0_) >= 1) && (coord_x != ref_x-1)) {
756 cpl_msg_warning (cpl_func, "Too much difference of channel %d with previous (%e pixels)", coord_x, y0 - y0_);
757 y0 = y0_;
758 sigma = sigma_;
759 }
760
761 /* In HIGH Resolution, we update ref_y in order to
762 * follow the curvature when extracting the data */
763 if (! (strcmp(resolution, "HIGH"))) {
764 ref_y = floor(y0);
765 }
766
767 /* Record x0 and sigma to compare to next channel */
768 y0_ = y0;
769 sigma_ = sigma;
770
771 /* Save best fit y0 and sigma into output vectors */
772 cpl_vector_set (values_x0, (coord_x - 1), coord_x);
773 cpl_vector_set (values_y0, (coord_x - 1), y0 );
774 cpl_vector_set (values_sigma, (coord_x - 1), sigma);
775
776 /* Increment spectral channel number */
777 coord_x --;
778
779 FREE (cpl_vector_delete, line_data);
780 FREE (cpl_vector_delete, z_data);
781 FREE (cpl_vector_delete, y_data);
782
783 } /* End fit the column ref_x -> 0 */
784
785
787 return CPL_ERROR_NONE;
788}
789
790/*----------------------------------------------------------------------------*/
812/*----------------------------------------------------------------------------*/
813
814cpl_image * gravi_create_profile_image (cpl_image * mean_img,
815 cpl_vector * values_x0,
816 cpl_vector * values_y0,
817 cpl_vector * values_sigma,
818 cpl_size iy_min,
819 cpl_size iy_max,
820 const char * mode)
821{
822 int nv = 0;
824 cpl_ensure (values_x0, CPL_ERROR_NULL_INPUT, NULL);
825 cpl_ensure (values_y0, CPL_ERROR_NULL_INPUT, NULL);
826 cpl_ensure (values_sigma, CPL_ERROR_NULL_INPUT, NULL);
827 cpl_ensure (mean_img, CPL_ERROR_NULL_INPUT, NULL);
828 cpl_ensure (mode, CPL_ERROR_NULL_INPUT, NULL);
829
830 /* Get the image dimensions */
831 cpl_size nx = cpl_image_get_size_x (mean_img);
832 cpl_size ny = cpl_image_get_size_y (mean_img);
833
834 cpl_msg_info (cpl_func, "iy_min = %lld, iy_max = %lld, ny = %lld",
835 iy_min, iy_max, ny);
836
837 cpl_ensure (iy_min >= 0 && iy_min < ny, CPL_ERROR_ILLEGAL_INPUT, NULL);
838 cpl_ensure (iy_max >= 0 && iy_max < ny, CPL_ERROR_ILLEGAL_INPUT, NULL);
839
840 /* Filter the profile params with a polynomial fit in spectral direction
841 * FIXME: could be a running median instead of fit, surely more stable.
842 * However the values contains some 'accidents' or modulations that
843 * may not be real and are efficiently removed by polynomial fit. */
844 cpl_matrix * matrix = cpl_matrix_wrap (1, nx, cpl_vector_get_data(values_x0));
845 cpl_size power;
846
847 /* Fit the y0 -- apply a median filter first */
848 cpl_polynomial * y0_poly = cpl_polynomial_new(1);
849 power = 4;
850 cpl_vector * temp_median = cpl_vector_filter_median_create(values_y0,3);
851 cpl_polynomial_fit (y0_poly, matrix, NULL, temp_median, NULL,
852 CPL_FALSE, NULL, &power);
853 cpl_vector_delete (temp_median);
854 CPLCHECK_NUL ("Cannot fit the y0");
855
856 /* Fit the sigma */
857 cpl_polynomial * sigma_poly = cpl_polynomial_new(1);
858 power = 5;
859 cpl_polynomial_fit (sigma_poly, matrix, NULL, values_sigma,
860 NULL, CPL_FALSE, NULL, &power);
861 CPLCHECK_NUL ("Cannot fit the sigma");
862
863 /* Compute the new y0 and sigma from these fits */
864 cpl_vector * valuesfit_y0 = cpl_vector_new(nx);
865 cpl_vector * valuesfit_sig = cpl_vector_new(nx);
866
867 for (cpl_size ix = 0; ix < nx; ix++){
868 double result;
869 result = cpl_polynomial_eval_1d (y0_poly, cpl_vector_get(values_x0, ix), NULL);
870 cpl_vector_set(valuesfit_y0, ix, result);
871 result = cpl_polynomial_eval_1d (sigma_poly, cpl_vector_get(values_x0, ix), NULL);
872 cpl_vector_set(valuesfit_sig, ix, result);
873 }
874 FREE (cpl_polynomial_delete, y0_poly);
875 FREE (cpl_polynomial_delete, sigma_poly);
876 cpl_matrix_unwrap (matrix);
877
878 /*
879 * Allocate image profile
880 */
881
882 cpl_image * region_img;
883 region_img = cpl_image_new (nx, ny, CPL_TYPE_DOUBLE);
884 cpl_image_fill_window (region_img, 1,1,nx,ny, 0.0);
885
886 /* Loop on spectral direction */
887 for (cpl_size ix = 0; ix < nx; ix++){
888
889 double sum_flux = 0, sum_flux2 = 0;
890
891 /* Loop on spatial direction */
892 for (cpl_size iy = iy_min; iy <= iy_max; iy++ ){
893
894 double result;
895
896 /* We use the measured profile */
897 if (!strcmp (mode, "PROFILE")) {
898 result = cpl_image_get (mean_img, ix+1, iy+1, &nv);
899 }
900 /* We use a fixed 5 pixel boxcar profile */
901 else if (!strcmp (mode, "BOX")) {
902 result = ( (fabs(iy - cpl_vector_get (valuesfit_y0, ix)) < 3 ) ? 1.0 : 0.0);
903 }
904 /* We use a Gaussian profile */
905 else if (!strcmp (mode, "GAUSS")) {
906 result = (iy - cpl_vector_get (valuesfit_y0, ix)) /
907 cpl_vector_get (valuesfit_sig, ix);
908 result = exp( - pow (result, 2) / 2);
909 } else {
910 cpl_msg_error (cpl_func, "BUG, report to DRS team");
911 return NULL;
912 }
913
914 /* Fill the profile image of this region */
915 cpl_image_set (region_img, ix + 1, iy + 1, result);
916
917 /* Compute normalization coefficients */
918 sum_flux += result;
919 sum_flux2 += result * result;
920 } /* End loop on spatial direction */
921
922 /* Keep only effective part of the profile
923 * Force pixel <1e-7 to zero */
924 for (cpl_size iy = iy_min; iy <= iy_max; iy++ ) {
925 double img_j = cpl_image_get (region_img, ix+1, iy+1, &nv);
926 if (img_j / sum_flux < 1e-7)
927 cpl_image_set (region_img, ix+1, iy+1, 0.0);
928 }
929
930 /* When we use a complex profile, we have to normalize
931 * it to make it flux conservative at extraction */
932 if (!strcmp (mode, "PROFILE") || !strcmp (mode, "GAUSS") ) {
933
934 double sum_flux = 0.0;
935 double sum_flux2 = 0.0;
936 for (cpl_size iy = iy_min; iy <= iy_max; iy++ ) {
937 double img_j = cpl_image_get (region_img, ix+1, iy+1, &nv);
938 sum_flux += img_j;
939 sum_flux2 += img_j * img_j;
940 }
941
942 for (cpl_size iy = iy_min; iy <= iy_max; iy++ ) {
943 double img_j = cpl_image_get (region_img, ix+1, iy+1, &nv);
944 cpl_image_set (region_img, ix+1, iy+1, img_j * sum_flux / sum_flux2);
945 }
946 }
947
948 } /* End loop on columns = spectral direction */
949
950 /* Deallocated of the fitted variables */
951 FREE (cpl_vector_delete, valuesfit_y0);
952 FREE (cpl_vector_delete, valuesfit_sig);
953
955 return region_img;
956}
957
958/*----------------------------------------------------------------------------*/
992/*----------------------------------------------------------------------------*/
993
995 gravi_data * dark_map, gravi_data * bad_map,
996 int nflat, const cpl_parameterlist * params)
997{
998 /* Verbose */
1000 cpl_ensure (flats_data, CPL_ERROR_NULL_INPUT, NULL);
1001 cpl_ensure (dark_map, CPL_ERROR_NULL_INPUT, NULL);
1002
1003 /* Ensure all beam open once */
1004 cpl_ensure (gravi_data_check_shutter_beam (flats_data, nflat),
1005 CPL_ERROR_ILLEGAL_INPUT, NULL);
1006
1007 /* Construct the product */
1008 gravi_data * out_data = gravi_data_new (0);
1009 cpl_propertylist * out_header = gravi_data_get_header (out_data);
1010
1011 /* Dump full header of first RAW data */
1012 gravi_data * raw0 = flats_data[0];
1013 cpl_propertylist_append (out_header, gravi_data_get_header (flats_data[0]));
1014
1015 /* Copy IMAGING_DETECTOR into product */
1018
1019 /* Create the FLUX array*/
1020 cpl_array * array_flux_FT = cpl_array_new(nflat, CPL_TYPE_DOUBLE);
1021 cpl_array * array_flux_SC = cpl_array_new(nflat, CPL_TYPE_DOUBLE);
1022
1023 /*
1024 * (1) Compute the FLAT of the FT
1025 */
1027 cpl_msg_warning (cpl_func,"The FLAT data has no IMAGING_DATA_FT");
1028 }
1029 else
1030 {
1031 cpl_msg_info (cpl_func, "Computing the FLAT of FT");
1032
1033 /* Get the dark of FT */
1034 cpl_table * darkft_table;
1035 darkft_table = gravi_data_get_table (dark_map, GRAVI_IMAGING_DATA_FT_EXT);
1036 CPLCHECK_NUL ("Cannot get DARK data of FT");
1037
1038 /* Get the DARK as an image */
1039 cpl_image * darkft_img;
1040 darkft_img = gravi_image_from_column (darkft_table, "PIX", 0);
1041
1042 /* Collapse each FLAT file */
1043 cpl_imagelist *imglist_ft = cpl_imagelist_new ();
1044 for (int file = 0; file < nflat; file++) {
1045
1046 cpl_table * dataft_table;
1047 dataft_table = gravi_data_get_table (flats_data[file], GRAVI_IMAGING_DATA_FT_EXT);
1048 CPLCHECK_NUL ("Cannot get data");
1049
1050 cpl_imagelist * imglistft_tmp;
1051 imglistft_tmp = gravi_imagelist_wrap_column (dataft_table, "PIX");
1052
1053 /* Remove DARK */
1054 cpl_imagelist_subtract_image (imglistft_tmp, darkft_img);
1055
1056 /* Collapse the DITs of this FLAT */
1057 cpl_imagelist_set (imglist_ft, cpl_imagelist_collapse_create (imglistft_tmp), file);
1058 gravi_imagelist_unwrap_images (imglistft_tmp);
1059
1060 /* compute the flux of this flat*/
1061 cpl_array_set(array_flux_FT, file, cpl_image_get_flux(cpl_imagelist_get(imglist_ft, file)));
1062 }
1063
1064 /* Collapse the FLAT files for all telescopes together */
1065 cpl_image * flatft_img = cpl_imagelist_collapse_create (imglist_ft);
1066 cpl_imagelist_delete (imglist_ft);
1067
1068 /* Create the flat_table */
1069 cpl_table * flatft_table = cpl_table_extract (gravi_data_get_table (flats_data[nflat-1],
1071
1072 cpl_array * flatft_array = gravi_array_wrap_image (flatft_img);
1073 cpl_table_set_array (flatft_table, "PIX", 0, flatft_array);
1074
1075 /* Remove median image and array */
1076 FREE (cpl_array_unwrap, flatft_array);
1077 FREE (cpl_image_delete, flatft_img);
1078 FREE (cpl_image_delete, darkft_img);
1079
1080 /* Set the FLAT as IMAGING_DATA into product */
1081 gravi_data_add_table (out_data, NULL, GRAVI_IMAGING_DATA_FT_EXT, flatft_table);
1082
1083 } /* End FLAT of FT */
1084
1085 /*
1086 * (2) General data used by FLAT and PROFILE for SC
1087 */
1088
1089 /* Get that the detector_table */
1090 cpl_table * detector_table;
1091 detector_table = gravi_data_get_table (flats_data[0], GRAVI_IMAGING_DETECTOR_SC_EXT);
1092 cpl_ensure (detector_table, CPL_ERROR_ILLEGAL_INPUT, NULL);
1093
1094 /* Get that the header of first file and of dark */
1095 cpl_propertylist * dark_header = gravi_data_get_header (dark_map);
1096 cpl_propertylist * flat0_header = gravi_data_get_header (flats_data[0]);
1097
1098 /* Get necessary information */
1099 int det_startx = gravi_pfits_get_window_start (flat0_header);
1100 const char * resolution = gravi_pfits_get_resolution (flat0_header);
1101 const char * pola_mode = gravi_pfits_get_pola_mode(flat0_header, GRAVI_SC);
1102 int nb_region = cpl_table_get_nrow (detector_table);
1103
1104 /* Get the DARK and BAD data */
1105 cpl_image * dark_img, * bad_img;
1106 dark_img = gravi_data_get_img (dark_map, GRAVI_IMAGING_DATA_SC_EXT);
1107 bad_img = gravi_data_get_img (bad_map, GRAVI_IMAGING_DATA_SC_EXT);
1108 cpl_ensure (dark_img && bad_img, CPL_ERROR_ILLEGAL_INPUT, NULL);
1109
1110 /*
1111 * (3) Compute the FLAT for the SC
1112 */
1113 cpl_msg_info (cpl_func, "Computing the FLAT of SC");
1114
1115 /* Imagelist to store the collapsed FLATs with all telescopes */
1116 cpl_imagelist * temp_imglist = cpl_imagelist_new ();
1117
1118 for (int file = 0; file < nflat; file++){
1119
1120 /* Extract necessary parameters and construct the output table */
1121 cpl_imagelist * data_imglist;
1122 data_imglist = gravi_data_get_cube (flats_data[file], GRAVI_IMAGING_DATA_SC_EXT);
1123
1124 /* Extract data with DARK, BADPIX in [ADU] */
1125 data_imglist = cpl_imagelist_duplicate (data_imglist);
1126 cpl_imagelist_subtract_image (data_imglist, dark_img);
1127 gravi_remove_badpixel_sc (data_imglist, bad_img);
1128
1129 /* Collapse the DITs of this FLAT */
1130 cpl_image * collapsed_img = cpl_imagelist_collapse_create (data_imglist);
1131 FREE (cpl_imagelist_delete, data_imglist);
1132
1133 /* Save this FLAT in the imagelist to collapse them */
1134 cpl_imagelist_set (temp_imglist, collapsed_img,
1135 cpl_imagelist_get_size (temp_imglist));
1136
1137 /* compute the flux of this flat*/
1138 cpl_array_set(array_flux_SC, file, cpl_image_get_flux(collapsed_img));
1139
1140 CPLCHECK_NUL ("Error");
1141 } /* End loop on FLATs*/
1142
1143 /* Collapse the FLATs for all telescopes together */
1144 cpl_image * allflat_img;
1145 allflat_img = cpl_imagelist_collapse_create (temp_imglist);
1146 FREE (cpl_imagelist_delete, temp_imglist);
1147
1148 CPLCHECK_NUL ("Cannot collapse FLATs");
1149
1150 /* Get the extension (illumated part) of the combined FLAT */
1151 int * ext_dim = gravi_image_extract_dimension (allflat_img);
1152 int fullstartx = ext_dim[0] + det_startx - 2;
1153
1154 /* Crop the FLAT into these new dimensions
1155 * FIXME: would be better to no crop FLAT */
1156 cpl_image * flatsc_img;
1157 flatsc_img = cpl_image_extract (allflat_img, ext_dim[0], 1,
1158 ext_dim[0] + ext_dim[1] - 1,
1159 cpl_image_get_size_y (allflat_img));
1160 FREE (cpl_image_delete, allflat_img);
1161
1162 /* Set the cropped FLAT of SC in output data
1163 * STARTX is in FITS convention (start 1) while
1164 * FULLSTARTX seems in C convention (start 0) */
1165 cpl_propertylist * flat_plist = cpl_propertylist_new ();
1166 cpl_propertylist_update_int (flat_plist, PROFILE_FULLSTARTX, fullstartx);
1167 cpl_propertylist_update_int (flat_plist, PROFILE_STARTX, ext_dim[0]);
1168 cpl_propertylist_update_int (flat_plist, PROFILE_NX, ext_dim[1]);
1169
1170 gravi_data_add_img (out_data, flat_plist, GRAVI_IMAGING_DATA_SC_EXT, flatsc_img);
1171 CPLCHECK_NUL ("Cannot set FLAT");
1172
1173
1174 /*
1175 * (4) Compute the PROFILE for the SC
1176 */
1177 cpl_msg_info (cpl_func, "Computing the PROFILE of SC");
1178
1179
1180 /* List of image to store each FLAT */
1181 cpl_image ** median_img = cpl_calloc (4, sizeof(cpl_image *));
1182
1183 /* kernel for median filter. Use more
1184 * pixels (~15) in HIGH, to skip the cluster */
1185 cpl_size size_kernel = !strcmp (resolution, "HIGH") ? 15 : 5;
1186 cpl_msg_info (cpl_func, "Median filtering over %lld spectral pixels", size_kernel);
1187 cpl_mask * kernel = cpl_mask_new (size_kernel, 1);
1188 cpl_mask_not (kernel);
1189
1190 for (int file = 0; file < nflat; file++){
1191 /* Get the primary header of each file */
1192 cpl_propertylist * flat_header;
1193 flat_header = gravi_data_get_header (flats_data[file]);
1194
1195 /* Extract necessary parameters and construct the output table */
1196 cpl_imagelist * data_imglist;
1197 data_imglist = gravi_data_get_cube (flats_data[file], GRAVI_IMAGING_DATA_SC_EXT);
1198 cpl_ensure (flat_header && data_imglist, CPL_ERROR_NULL_INPUT, NULL);
1199
1200 /* Extract data with DARK, BADPIX in [ADU] */
1201 data_imglist = cpl_imagelist_duplicate (data_imglist);
1202 cpl_imagelist_subtract_image (data_imglist, dark_img);
1203 gravi_remove_badpixel_sc (data_imglist, bad_img);
1204
1205 /* Collapse the DITs of this FLAT */
1206 cpl_image * collapsed_img = cpl_imagelist_collapse_create (data_imglist);
1207
1208 FREE (cpl_imagelist_delete, data_imglist);
1209
1210 /* Create a filtered version of this FLAT */
1211 cpl_image * filtered_img = cpl_image_duplicate (collapsed_img);
1212 cpl_image_filter_mask (filtered_img, collapsed_img, kernel,
1213 CPL_FILTER_MEDIAN, CPL_BORDER_FILTER);
1214 FREE (cpl_image_delete, collapsed_img);
1215
1216 /* Crop it */
1217 cpl_image * crop_img;
1218 crop_img = cpl_image_extract (filtered_img, ext_dim[0], 1,
1219 ext_dim[0] + ext_dim[1] - 1,
1220 cpl_image_get_size_y (filtered_img));
1221 FREE (cpl_image_delete, filtered_img);
1222
1223 /* Save this filtered version in the median_img[]
1224 * According to the beam open */
1225 int id = gravi_get_shutter_id (flat_header);
1226 median_img[id] = crop_img;
1227
1228 CPLCHECK_NUL ("Error");
1229 } /* End loop on FLATs*/
1230
1231 FREE (cpl_mask_delete, kernel);
1232
1233 /* Get the image dimensions after crop */
1234 int nx = cpl_image_get_size_x (median_img[0]);
1235 int ny = cpl_image_get_size_y (median_img[0]);
1236
1237 /* Get the profile width parameter */
1238 int profile_width = gravi_param_get_int (params, "gravity.calib.profile-width");
1239 cpl_ensure (profile_width > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1240
1241 /* Compute the size_profile. For LOW and MEDIUM it is automatic.
1242 * For HIGH it comes from the option */
1243 int n_darkline;
1244 cpl_propertylist * flat_header;
1245 int window_mode;
1246 flat_header = gravi_data_get_header (flats_data[0]);
1247 int size_profile;
1248 /* if the windowing mode if after the 05/12/2016 */
1249 if ( cpl_propertylist_get_float(flat_header, "MJD-OBS") > 57728 ){
1250 n_darkline= 0;
1251 window_mode = 1;
1252 cpl_msg_info(cpl_func, "Windowing after MJD = 57728");
1253 if (! (strcmp(resolution, "LOW") && strcmp(resolution, "MED")) && !strcmp(pola_mode, "COMB")) {
1254 size_profile = ny/nb_region;
1255 cpl_msg_info(cpl_func, "Use a computed size_profile of %d", size_profile);
1256 }
1257 else{
1258 size_profile = (fmod(profile_width, 2) == 0) ? profile_width + 1 : profile_width;
1259 cpl_msg_info(cpl_func, "Use a given size profile of %d", size_profile);
1260 }
1261 }
1262 /* else keep old windowing for backward compatibility */
1263 else {
1264 n_darkline= 1;
1265 window_mode = 0;
1266 if (! (strcmp(resolution, "LOW") && strcmp(resolution, "MED"))) {
1267 size_profile = (ny-(nb_region+1)*n_darkline)/nb_region;
1268 cpl_msg_info(cpl_func, "Use a computed size_profile of %d", size_profile);
1269 }
1270 else{
1271 size_profile = (fmod(profile_width, 2) == 0) ? profile_width + 1 : profile_width;
1272 cpl_msg_info(cpl_func, "Use a given size profile of %d", size_profile);
1273 }
1274 }
1275
1276 /* Construction of the PROFILE_DATA table */
1277 cpl_table * profile_table = cpl_table_new (1);
1278
1279 /* Create the DATA# columns with their dimension */
1280 cpl_array * dimension = cpl_array_new (2, CPL_TYPE_INT);
1281 cpl_array_set (dimension, 0, nx);
1282 cpl_array_set (dimension, 1, ny);
1283 for (int region = 0; region < nb_region; region ++){
1284 const char * data = GRAVI_DATA[region];
1285 cpl_table_new_column_array (profile_table, data,
1286 CPL_TYPE_DOUBLE, nx * ny);
1287 cpl_table_set_column_dimensions (profile_table, data, dimension);
1288 }
1289 FREE (cpl_array_delete, dimension);
1290
1291 /* Construction of the PROFILE_PARAMS table */
1292 cpl_table * params_table = cpl_table_new (nb_region);
1293 cpl_table_new_column_array (params_table, "CENTERY", CPL_TYPE_DOUBLE, nx);
1294 cpl_table_new_column_array (params_table, "WIDTH", CPL_TYPE_DOUBLE, nx);
1295
1296
1297 /* Construction of a map used to zero the badpixels
1298 * Mask is: good pixels = 1 ; bad pixels = 0 */
1299 cpl_image * mask_img = NULL;
1300
1301 if ( !gravi_param_get_bool (params,
1302 "gravity.calib.force-badpix-to-zero") )
1303 {
1304 cpl_msg_info (cpl_func, "Bad pixels are *not*"
1305 "forced to zero in profiles");
1306 cpl_image_fill_window (mask_img, 1, 1, nx, ny, 1.0);
1307 }
1308 else {
1309 cpl_msg_info (cpl_func, "Bad pixels are "
1310 "forced to zero in profiles");
1311 mask_img = cpl_image_extract (bad_img, ext_dim[0], 1,
1312 ext_dim[0] + ext_dim[1] - 1, ny);
1313 cpl_image_threshold (mask_img, 0.5, 0.5, 1.0, 0.0);
1314 }
1315
1316 /* Update mask from FLAT itself in LOW,
1317 * FIXME: to decide what we zero exactly. */
1318 if ( !strcmp (resolution, "LOW") )
1319 {
1320 cpl_msg_info (cpl_func, "Pixels with low FLAT values"
1321 " are forced to zero in profiles.");
1322
1323 double threshold = cpl_propertylist_get_double (dark_header, QC_DARKRMS_SC);
1324 cpl_image * mask2_img = cpl_image_duplicate (flatsc_img);
1325 cpl_image_threshold (mask2_img, threshold, threshold, 0.0, 1.0);
1326 cpl_image_multiply (mask_img, mask2_img);
1327 FREE (cpl_image_delete, mask2_img);
1328 }
1329
1330 /* Define the mode */
1331 const char * mode = gravi_param_get_string_default (params,
1332 "gravity.calib.profile-mode", "AUTO");
1333
1334 if (!strcmp(mode, "AUTO")) {
1335 if (!strcmp(resolution, "LOW")) mode = "PROFILE";
1336 if (!strcmp(resolution, "MED")) mode = "PROFILE";
1337 if (!strcmp(resolution, "HIGH")) mode = "BOX";
1338 }
1339
1340 cpl_msg_info (cpl_func, "Profile computed with mode: %s (%s)", mode, resolution);
1341
1342 /* Loop on regions */
1343 for (int region = 0; region < nb_region ; region++) {
1344
1345 int tel_1 = gravi_region_get_tel (detector_table, region, 0);
1346 int tel_2 = gravi_region_get_tel (detector_table, region, 1);
1347 CPLCHECK_NUL ("Cannot get the telescope from region");
1348
1349 /* Compute the mean image between the first
1350 * and second telescope */
1351 cpl_image * mean_img;
1352 mean_img = cpl_image_add_create (median_img[tel_1],
1353 median_img[tel_2]);
1354 cpl_image_divide_scalar (mean_img, 2);
1355
1356 /* Get the reference coordinates of this region. That is
1357 * a point where the spectra is supposed to be for sure
1358 * In unit of the full detector window */
1359 int ref0_x = gravi_table_get_value (detector_table, "CENTER", region, 0);
1360 int ref0_y = gravi_table_get_value (detector_table, "CENTER", region, 1);
1361
1362 /* Convert the reference coordinates in unit of the cropped data
1363 * Actually we force ref_x and ref_y in LOW and MED. */
1364 int ref_x, ref_y;
1365
1366 /* if the windowing mode if after the 05/12/2016 */
1367 if ( window_mode == 1 ){
1368 if (! (strcmp(resolution, "LOW") && strcmp(resolution, "MED")) && !strcmp(pola_mode, "COMB") ) {
1369 ref_x = nx / 2;
1370 ref_y = (region+1)*(size_profile)-size_profile/2;
1371 } else {
1372 //ref_x = ref0_x - det_startx - (1 + ext_dim[0]);
1373 ref_x = ref0_x - (1 + ext_dim[0]);
1374 ref_y = ref0_y;
1375 }
1376 }
1377 else {
1378 if (! (strcmp(resolution, "LOW") && strcmp(resolution, "MED")) ) {
1379 ref_x = nx / 2;
1380 ref_y = (region+1)*(size_profile+n_darkline)-size_profile/2;
1381 } else {
1382 //ref_x = ref0_x - det_startx - (1 + ext_dim[0]);
1383 ref_x = ref0_x - (1 + ext_dim[0]);
1384 ref_y = ref0_y;
1385 }
1386 }
1387
1388 /*
1389 * Fit Gaussian to profile of each spectral channel
1390 */
1391 cpl_vector * values_x0 = cpl_vector_new (nx);
1392 cpl_vector * values_y0 = cpl_vector_new (nx);
1393 cpl_vector * values_sigma = cpl_vector_new (nx);
1394
1395 gravi_fit_profile (values_x0, values_y0, values_sigma,
1396 mean_img, ref_x, ref_y, size_profile,
1397 resolution);
1398 CPLCHECK_NUL ("Cannot fit data into profile params");
1399
1400 /* Fill the PROFILE_PARAM table for this region */
1401 cpl_array * values_arr;
1402 values_arr = cpl_array_wrap_double (cpl_vector_get_data (values_y0), nx);
1403 cpl_table_set_array (params_table, "CENTERY", region, values_arr);
1404 cpl_array_unwrap (values_arr);
1405
1406 values_arr = cpl_array_wrap_double (cpl_vector_get_data(values_sigma), nx);
1407 cpl_table_set_array (params_table, "WIDTH", region, values_arr);
1408 cpl_array_unwrap (values_arr);
1409
1410
1411 /*
1412 * Create and compute the profile image of this region
1413 */
1414
1415 /* Define the spatial pixels over which we fill the profile */
1416 cpl_size iy_min = 0, iy_max = ny-1;
1417 if (! (strcmp(resolution, "LOW") && strcmp(resolution, "MED")) ){
1418 iy_min = ref_y - size_profile/2;
1419 iy_max = ref_y + size_profile/2-1;
1420 }
1421
1422 /* Fill the profile image */
1423 cpl_image * region_img;
1424 region_img = gravi_create_profile_image (mean_img, values_x0,
1425 values_y0, values_sigma,
1426 iy_min, iy_max, mode);
1427 FREE (cpl_vector_delete, values_x0);
1428 FREE (cpl_vector_delete, values_y0);
1429 FREE (cpl_vector_delete, values_sigma);
1430 CPLCHECK_NUL ("Cannot build profile image");
1431
1432 /* Force some pixel to zero with computed mask */
1433 cpl_image_multiply (region_img, mask_img);
1434
1435 /* Fill the profile_table for this region */
1436 cpl_array * array = gravi_array_wrap_image (region_img);
1437 cpl_table_set_array (profile_table, GRAVI_DATA[region], 0, array);
1438 FREE (cpl_array_unwrap, array);
1439 FREE (cpl_image_delete, region_img);
1440
1441 /* Free the mean image of the 2 flats of this region */
1442 FREE (cpl_image_delete, mean_img);
1443
1444 } /* End loop on regions */
1445
1446 /* Save the PROFILE_DATA table */
1447 cpl_propertylist * profile_plist = cpl_propertylist_new ();
1448 cpl_propertylist_copy_property (profile_plist, flat_plist, PROFILE_FULLSTARTX);
1449 cpl_propertylist_copy_property (profile_plist, flat_plist, PROFILE_STARTX);
1450 cpl_propertylist_copy_property (profile_plist, flat_plist, PROFILE_NX);
1451 gravi_data_add_table (out_data, profile_plist,
1452 GRAVI_PROFILE_DATA_EXT, profile_table);
1453
1454 /* Save the PROFILE_PARAMS table */
1455 cpl_propertylist * params_plist = cpl_propertylist_duplicate (profile_plist);
1456 gravi_data_add_table (out_data, params_plist,
1457 GRAVI_PROFILE_PARAMS_EXT, params_table);
1458
1459
1460 /*
1461 * (5) Add the QC parameter of the lateral positioning of the first region
1462 */
1463
1464 cpl_propertylist * main_header = gravi_data_get_header (out_data);
1465 double qc_value;
1466 char qc_name[100];
1467 cpl_size idx;
1468
1469 for (int file = 0; file < nflat; file++){
1470 sprintf (qc_name, "ESO QC FLATFLUX FT%i", file+1);
1471 qc_value = cpl_array_get(array_flux_FT, file, NULL);
1472 cpl_propertylist_append_double (main_header, qc_name, qc_value);
1473 cpl_propertylist_set_comment (main_header, qc_name, "[ADU] Total flux per DIT");
1474 cpl_msg_info (cpl_func,"%s = %f [ADU]", qc_name, qc_value);
1475 }
1476 for (int file = 0; file < nflat; file++){
1477 sprintf (qc_name, "ESO QC FLATFLUX SC%i", file+1);
1478 qc_value = cpl_array_get(array_flux_SC, file, NULL);
1479 cpl_propertylist_append_double (main_header, qc_name, qc_value);
1480 cpl_propertylist_set_comment (main_header, qc_name, "[ADU] Total flux per DIT");
1481 cpl_msg_info (cpl_func,"%s = %f [ADU]", qc_name, qc_value);
1482 }
1483
1484 for (int reg = 0; reg < nb_region; reg += 12) {
1485
1486 /* Median lateral position */
1487 qc_value = cpl_array_get_median (cpl_table_get_array (params_table, "CENTERY", reg));
1488 qc_value = floor (qc_value * 1e6) * 1e-6;
1489 sprintf (qc_name, "ESO QC PROFILE_CENTER SC%i MED", reg+1);
1490 cpl_propertylist_append_double (main_header, qc_name, qc_value);
1491 cpl_propertylist_set_comment (main_header, qc_name, "[pixel] position of region");
1492 cpl_msg_info (cpl_func,"%s = %f [pixel]", qc_name, qc_value);
1493
1494 /* Median width */
1495 qc_value = cpl_array_get_median (cpl_table_get_array (params_table, "WIDTH", reg));
1496 qc_value = floor (qc_value * 1e6) * 1e-6;
1497 sprintf (qc_name, "ESO QC PROFILE_WIDTH SC%i MED", reg+1);
1498 cpl_propertylist_append_double (main_header, qc_name, qc_value);
1499 cpl_propertylist_set_comment (main_header, qc_name, "[pixel] width of region");
1500 cpl_msg_info (cpl_func,"%s = %f [pixel]", qc_name, qc_value);
1501
1502 /* Lateral position in the left */
1503 idx = 1 * ext_dim[1] / 6;
1504 qc_value = cpl_array_get (cpl_table_get_array (params_table, "CENTERY", reg), idx, NULL);
1505 qc_value = floor (qc_value * 1e6) * 1e-6;
1506 sprintf (qc_name, "ESO QC PROFILE_CENTER SC%i LEFT", reg+1);
1507 cpl_propertylist_append_double (main_header, qc_name, qc_value);
1508 cpl_propertylist_set_comment (main_header, qc_name, "[pixel] at STARTX+NX/6");
1509 cpl_msg_info (cpl_func,"%s = %f [pixel] for x=%lld", qc_name, qc_value, idx);
1510
1511 /* Lateral position in the right */
1512 idx = 5 * ext_dim[1] / 6;
1513 qc_value = cpl_array_get (cpl_table_get_array (params_table, "CENTERY", reg), idx, NULL);
1514 qc_value = floor (qc_value * 1e6) * 1e-6;
1515 sprintf (qc_name, "ESO QC PROFILE_CENTER SC%i RIGHT", reg+1);
1516 cpl_propertylist_append_double (main_header, qc_name, qc_value);
1517 cpl_propertylist_set_comment (main_header, qc_name, "[pixel] at STARTX+5*NX/6");
1518 cpl_msg_info (cpl_func,"%s = %f [pixel] for x=%lld", qc_name, qc_value, idx);
1519
1520 CPLCHECK_NUL ("Cannot compute QC parameters");
1521 }
1522
1523 /* Deallocation of the variables */
1524 FREE (cpl_image_delete, mask_img);
1525 FREELOOP (cpl_image_delete, median_img, 4);
1526 FREE (cpl_free, ext_dim);
1527 cpl_array_delete(array_flux_FT);
1528 cpl_array_delete(array_flux_SC);
1529
1530 /* Verbose */
1532 return out_data;
1533}
1534
1535/*----------------------------------------------------------------------------*/
1549/*----------------------------------------------------------------------------*/
1550
1551cpl_propertylist * gravi_compute_gain (gravi_data ** flats_data,
1552 int nrawgain,
1553 gravi_data * dark_map)
1554{
1555 int nv;
1556 const cpl_size maxdeg = 1, mindeg = 0;
1557 const cpl_size slope_deg = 1;
1558
1559 /* Verbose */
1561 cpl_ensure (flats_data, CPL_ERROR_NULL_INPUT, NULL);
1562 cpl_ensure (dark_map, CPL_ERROR_NULL_INPUT, NULL);
1563
1564 /* Ensure all beam open once */
1565 cpl_ensure (gravi_data_check_shutter_beam (flats_data, nrawgain),
1566 CPL_ERROR_ILLEGAL_INPUT, NULL);
1567
1568 /* Get a sequence of 4 single-shutter open files
1569 * FIXME: remove this part isn't it ?? */
1570 cpl_msg_info (cpl_func, "Search for 4-shutter sequence");
1571
1572 gravi_data ** gain_file = cpl_calloc(4,sizeof(gravi_data*));
1573
1574 for (int i = 0; i < nrawgain; i++){
1575 cpl_propertylist * flat_header = gravi_data_get_header (flats_data[i]);
1576
1577 if (gravi_check_shutter (flat_header, 1,0,0,0)) {
1578 gain_file[0] = (flats_data[i]);
1579 }
1580 if (gravi_check_shutter (flat_header, 0,1,0,0)) {
1581 gain_file[1] = (flats_data[i]);
1582 }
1583 if (gravi_check_shutter (flat_header, 0,0,1,0)) {
1584 gain_file[2] = (flats_data[i]);
1585 }
1586 if (gravi_check_shutter (flat_header, 0,0,0,1)) {
1587 gain_file[3] = (flats_data[i]);
1588 }
1589 }
1590
1591 /*
1592 * Create the output. It is a plist only since we
1593 * don't create a gravi_data
1594 */
1595 cpl_propertylist * output_plist = cpl_propertylist_new();
1596
1597
1598 /*
1599 *
1600 * Compute the gain for SC
1601 *
1602 */
1604 cpl_msg_warning (cpl_func,"The FLAT data has no IMAGING_DATA_SC");
1605 }
1606 else
1607 {
1608 cpl_msg_info (cpl_func, "Computing the gain of SC");
1609
1610 /* Get the size of the image */
1611 cpl_image * image = gravi_data_get_img (gain_file[0], GRAVI_IMAGING_DATA_SC_EXT);
1612 cpl_size size = cpl_image_get_size_x (image) * cpl_image_get_size_y (image);
1613
1614 /* Build matrix and vector to fit var = f(mean) */
1615 cpl_vector * vector_var = cpl_vector_new (4*size);
1616 cpl_matrix * matrix_mean = cpl_matrix_new (1, 4*size);
1617
1618 /* Get a pointer to the dark image */
1619 cpl_msg_info (cpl_func,"DARK and BADPIX of SC are not used");
1620
1621 /* Loop on files to fill the variance and mean vectors */
1622 for (int file = 0; file < 4; file++) {
1623
1624 /* Get the imagelist and property lists */
1625 cpl_imagelist * data_imglist;
1626 data_imglist = gravi_data_get_cube (gain_file[file], GRAVI_IMAGING_DATA_SC_EXT);
1627
1628 /* Compute <image> and <image^2> */
1629 cpl_imagelist * imglist1 = cpl_imagelist_new ();
1630 cpl_imagelist * imglist2 = cpl_imagelist_new ();
1631 for (cpl_size i = 0; i < cpl_imagelist_get_size (data_imglist); i ++){
1632 image = cpl_image_cast (cpl_imagelist_get (data_imglist, i), CPL_TYPE_DOUBLE);
1633 cpl_imagelist_set (imglist1, image, i);
1634 cpl_imagelist_set (imglist2, cpl_image_power_create (image, 2), i);
1635 }
1636 cpl_image * image1_mean = cpl_imagelist_collapse_create (imglist1);
1637 cpl_image * image2_mean = cpl_imagelist_collapse_create (imglist2);
1638 FREE (cpl_imagelist_delete, imglist1);
1639 FREE (cpl_imagelist_delete, imglist2);
1640
1641 /* Loop on pixels to fill the vector and matrix structures
1642 * used latter in the PTC fit */
1643 cpl_size nx = cpl_image_get_size_x (image1_mean);
1644 cpl_size ny = cpl_image_get_size_y (image1_mean);
1645 for (cpl_size i = 0; i < nx; i ++) {
1646 for (cpl_size j = 0; j < ny; j ++) {
1647 cpl_vector_set (vector_var, file * size + (j + ny*i),
1648 cpl_image_get (image2_mean, i+1, j+1, &nv) -
1649 pow (cpl_image_get (image1_mean, i+1, j+1, &nv), 2));
1650 cpl_matrix_set (matrix_mean, 0, file * size + (j + ny*i),
1651 cpl_image_get (image1_mean, i+1, j+1, &nv));
1652 }
1653 }
1654 /* End loop on pixels */
1655
1656 /* Deallocation of variables */
1657 FREE (cpl_image_delete, image1_mean);
1658 FREE (cpl_image_delete, image2_mean);
1659 }
1660 /* End loop on files */
1661
1662 /* Fit the curve variance versus the median,
1663 thus the output gains are in ADU/e */
1664 cpl_polynomial * fit_slope = cpl_polynomial_new (1);
1665 cpl_polynomial_fit (fit_slope, matrix_mean, NULL, vector_var, NULL,
1666 CPL_FALSE, &mindeg, &maxdeg);
1667
1668 /* Get the slope */
1669 double slope = cpl_polynomial_get_coeff (fit_slope, &slope_deg);
1670 cpl_msg_info (cpl_func, "mean gain SC = %.4f [adu/e-] Mean gain of detector", slope);
1671 cpl_propertylist_append_double (output_plist, QC_MEANGAIN_SC, slope);
1672 cpl_propertylist_set_comment (output_plist, QC_MEANGAIN_SC, "[adu/e-] Mean gain of SC detector" );
1673
1674 /* Delete */
1675 cpl_vector_delete (vector_var);
1676 cpl_matrix_delete (matrix_mean);
1677 cpl_polynomial_delete (fit_slope);
1678 }
1679 /* End SC case */
1680
1681 /*
1682 *
1683 * Compute the gain for FT
1684 *
1685 */
1687 cpl_msg_warning (cpl_func,"The FLAT data has no IMAGING_DATA_FT");
1688 }
1689 else
1690 {
1691 cpl_msg_info (cpl_func, "Computing the gain of FT");
1692
1693 /* Get DARK image */
1694 cpl_table * dark_table;
1695 dark_table = gravi_data_get_table (dark_map, GRAVI_IMAGING_DATA_FT_EXT);
1696 cpl_image * dark_img = gravi_image_from_column (dark_table, "PIX", 0);
1697
1698 /* Get the size of the image */
1699 cpl_size size = cpl_table_get_column_depth (dark_table, "PIX");
1700
1701 /* Build matrix and vector to fit var = f(mean) */
1702 cpl_vector * vector_var = cpl_vector_new (4 * size);
1703 cpl_matrix * matrix_mean = cpl_matrix_new (1, 4 * size);
1704
1705 /* Compute the gain for each file */
1706 for (int file = 0; file < 4; file++) {
1707
1708 /* Get the tables */
1709 cpl_table * data_table;
1710 data_table = gravi_data_get_table (gain_file[file], GRAVI_IMAGING_DATA_FT_EXT);
1711 cpl_table_cast_column (data_table, "PIX", "PIX", CPL_TYPE_DOUBLE);
1712 cpl_imagelist * data_imglist = gravi_imagelist_from_column (data_table,"PIX");
1713 CPLCHECK_NUL ("Cannot get data");
1714
1715 /* Remove DARK */
1716 cpl_imagelist_subtract_image (data_imglist, dark_img);
1717
1718 /* Compute MEAN */
1719 cpl_image * mean_img = cpl_imagelist_collapse_create (data_imglist);
1720
1721 /* Compute VARIANCE */
1722 cpl_imagelist_subtract_image (data_imglist, mean_img);
1723 cpl_imagelist_power (data_imglist, 2.0);
1724 cpl_image * var_img = cpl_imagelist_collapse_create (data_imglist);
1725
1726 FREE (cpl_imagelist_delete, data_imglist);
1727
1728 /* Loop on pixels to fill the vector and matrix structures
1729 * used latter in the PTC fit */
1730 cpl_size nx = cpl_image_get_size_x (mean_img);
1731 cpl_size ny = cpl_image_get_size_y (mean_img);
1732 for (cpl_size i = 0; i < nx; i ++) {
1733 for (cpl_size j = 0; j < ny; j ++) {
1734 cpl_vector_set (vector_var, file * size + (j + ny*i),
1735 cpl_image_get (var_img, i+1, j+1, &nv));
1736 cpl_matrix_set (matrix_mean, 0, file * size + (j + ny*i),
1737 cpl_image_get (mean_img, i+1, j+1, &nv));
1738 }
1739 }
1740 /* End loop on pixels */
1741
1742 /* Deallocation of variables */
1743 FREE (cpl_image_delete, var_img);
1744 FREE (cpl_image_delete, mean_img);
1745 } /* End loop on files */
1746
1747 /* Fit the curve variance versus the median,
1748 thus the output gains are in ADU/e */
1749 cpl_polynomial * fit_slope = cpl_polynomial_new (1);
1750 cpl_polynomial_fit (fit_slope, matrix_mean, NULL, vector_var, NULL,
1751 CPL_FALSE, &mindeg, &maxdeg);
1752
1753 /* Get the slope */
1754 double slope = cpl_polynomial_get_coeff (fit_slope, &slope_deg);
1755 cpl_msg_info (cpl_func, "mean gain FT = %.4f [adu/e-] Mean gain of detector", slope);
1756 cpl_propertylist_append_double (output_plist, QC_MEANGAIN_FT, slope);
1757 cpl_propertylist_set_comment (output_plist, QC_MEANGAIN_FT, "[adu/e-] Mean gain of FT detector");
1758
1759 /* Delete */
1760 FREE (cpl_vector_delete, vector_var);
1761 FREE (cpl_matrix_delete, matrix_mean);
1762 FREE (cpl_polynomial_delete, fit_slope);
1763 FREE (cpl_image_delete, dark_img);
1764
1765 } /* End FT case*/
1766
1767 cpl_free (gain_file);
1768
1769 /* Verbose */
1771 return output_plist;
1772}
1773
1774/*----------------------------------------------------------------------------*/
1804/*----------------------------------------------------------------------------*/
1805
1807 gravi_data ** flats_data,
1808 int nflat,
1809 const cpl_parameterlist * params)
1810{
1811 /* Verbose */
1813 cpl_ensure (dark_map, CPL_ERROR_NULL_INPUT, NULL);
1814 cpl_ensure (params, CPL_ERROR_NULL_INPUT, NULL);
1815 cpl_ensure (nflat==4 || nflat==0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1816
1817 /* Construction of the bad pixels map */
1818 gravi_data * bad_map = gravi_data_new(0);
1819 cpl_propertylist * bad_header = gravi_data_get_header (bad_map);
1820
1821 /* Dump full header into product */
1822 cpl_propertylist * dark_header = gravi_data_get_header (dark_map);
1823 cpl_propertylist_append (bad_header, dark_header);
1824
1825 /*
1826 * Compute bad pixels of FT
1827 */
1828
1830 cpl_msg_warning (cpl_func,"The DARK map has no IMAGING_DATA_FT");
1831 }
1832 else
1833 {
1834 cpl_msg_info (cpl_func,"Compute BADPIXEL of FT");
1835
1836 int bad_pix_number_A = gravi_param_get_int (params, "gravity.calib.bad-pixel-A-ft")-1;
1837 int bad_pix_number_B = gravi_param_get_int (params, "gravity.calib.bad-pixel-B-ft")-1;
1838
1839 /* Copy necessary tables */
1841
1842 /* This is the FT */
1843 cpl_table * dark_table = gravi_data_get_table (dark_map, GRAVI_IMAGING_DATA_FT_EXT);
1844 cpl_table * std_table = gravi_data_get_table (dark_map, GRAVI_IMAGING_ERR_FT_EXT);
1845
1846 /* Verify it has one single row */
1847 cpl_size n_row = cpl_table_get_nrow (dark_table);
1848 cpl_ensure (n_row == 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1849 cpl_ensure (n_row == cpl_table_get_nrow (std_table), CPL_ERROR_ILLEGAL_INPUT, NULL);
1850
1851 /* Get the rms factor for dark bad pixel threshold */
1852 int bad_dark_factor = gravi_param_get_int (params, "gravity.calib.bad-dark-threshold");
1853
1854 CPLCHECK_NUL ("Cannot get data");
1855
1856 /* Get the dark rms and the mean dark */
1857 double dark_rms = cpl_propertylist_get_double (dark_header, QC_DARKRMS_FT);
1858 double dark_mean = cpl_propertylist_get_double (dark_header, QC_MEANDARK_FT);
1859 CPLCHECK_NUL ("Cannot get QC");
1860
1861
1862 /* Compute the specified range are declared as bad pixels */
1863 double range_max = dark_mean + bad_dark_factor * dark_rms;
1864 cpl_msg_info (cpl_func,"FT threshold on mean value = %f",range_max);
1865
1866 /* Get the min std to declare pixel as bad */
1867 double std_min = 1e-3;
1868
1869 /* Get the max std to declare pixel as bad */
1870 cpl_size npix = cpl_table_get_column_depth (dark_table, "PIX");
1871 cpl_array * std_array = cpl_table_get_data_array (std_table, "PIX")[0];
1872 cpl_vector * std_vector = cpl_vector_new (npix);
1873 for (cpl_size pix = 0; pix < npix; pix++) {
1874 cpl_vector_set (std_vector, pix, cpl_array_get (std_array, pix, NULL));
1875 }
1876 cpl_vector_sort (std_vector, CPL_SORT_ASCENDING);
1877 /* use 21/24 percentil to include the rms of the metrology*/
1878 cpl_size percentil = (int) npix*0.875;
1879 double std_max = cpl_vector_get (std_vector, percentil)*(1+bad_dark_factor/8.75);
1880 cpl_msg_info (cpl_func,"FT threshold on std value = %f",std_max);
1881
1882 /* Compute the number of bad pixel */
1883 int count_bp_dark = 0;
1884
1885 /* Create output table */
1886 cpl_table * bad_table = cpl_table_extract (dark_table, 0, 1);
1887
1888 /* Loop on pixels */
1889 cpl_array * dark_array = cpl_table_get_data_array (dark_table, "PIX")[0];
1890 cpl_array * bad_array = cpl_table_get_data_array (bad_table, "PIX")[0];
1891
1892 for (cpl_size pix = 0; pix < npix; pix++) {
1893 cpl_array_set (bad_array, pix, 0);
1894 int pixel_removed=0;
1895 /* flag designated pixels */
1896 if ((pix == bad_pix_number_A)||(pix == bad_pix_number_B))
1897 {
1898 cpl_array_set (bad_array, pix, BADPIX_DARK);
1899 pixel_removed=1;
1900 cpl_msg_info(cpl_func,"Detected a bad FT pixel at position %lli based on reduction option", pix+1);
1901 }
1902 /* flag on the mean value of the dark */
1903 if (cpl_array_get (dark_array, pix, NULL) > range_max) {
1904 cpl_array_set (bad_array, pix, BADPIX_DARK);
1905 pixel_removed=1;
1906 cpl_msg_info(cpl_func,"Detected a bad FT pixel at position %lli based on its mean flux: %f ADU", pix+1, cpl_array_get (dark_array, pix, NULL));
1907 }
1908 /* flag on the std value (min and max) of the dark */
1909 if (cpl_array_get (std_array, pix, NULL) < std_min) {
1910 cpl_array_set (bad_array, pix, BADPIX_DARK);
1911 pixel_removed=1;
1912 cpl_msg_warning(cpl_func,"Detected a bad FT pixel at position %lli based on its low rms: %f ADU", pix+1, cpl_array_get (std_array, pix, NULL));
1913 } else if (cpl_array_get (std_array, pix, NULL) > std_max) {
1914 cpl_array_set (bad_array, pix, BADPIX_DARK);
1915 pixel_removed=1;
1916 cpl_msg_warning(cpl_func,"Detected a bad FT pixel at position %lli based on its high rms: %f ADU", pix+1, cpl_array_get (std_array, pix, NULL));
1917 }
1918 if (pixel_removed == 1) count_bp_dark ++;
1919 }
1920
1921 /* check if low flux flag is on. If yes, set to remove low flux pixels */
1922 if (!nflat && gravi_param_get_bool(params, "gravity.calib.lowflux-pixels-ft")) {
1923 cpl_msg_warning (cpl_func, "Option to remove low flux pixels applied, but no FLATs were provided");
1924 }
1925 else if (gravi_param_get_bool(params, "gravity.calib.lowflux-pixels-ft")) {
1926 cpl_ensure (flats_data, CPL_ERROR_NULL_INPUT, NULL);
1927 cpl_msg_info (cpl_func,"FLATs used to remove low flux values on FT");
1928
1929 cpl_array * flat_array_sum = cpl_array_new (npix, CPL_TYPE_DOUBLE);
1930 cpl_array_fill_window (flat_array_sum, 0, npix, 0);
1931
1932 int count_bp_lowflux = 0;
1933
1934 for (int f = 0; f<nflat; f++) {
1935 cpl_table * flat = gravi_data_get_table (flats_data[f], GRAVI_IMAGING_DATA_FT_EXT);
1936 cpl_size n_row = cpl_table_get_nrow (flat);
1937 for (int row = 0; row<n_row; row++) {
1938 cpl_array * flat_array = cpl_table_get_data_array (flat, "PIX")[row];
1939 cpl_array_add (flat_array_sum, flat_array);
1940 cpl_array_subtract (flat_array_sum, dark_array);
1941 }
1942 }
1943 for (cpl_size pix = 0; pix < npix-1; pix+=2) {
1944 /* remove pixel if the flat flux is lower than 1/3 of the flux on the adjacent pixel*/
1945 double pix1=cpl_array_get (flat_array_sum, pix, NULL);
1946 double pix2=cpl_array_get (flat_array_sum, pix+1, NULL);
1947
1948 /* make sure the pixel flux are positive, otherwise, scale then up */
1949 if (pix1<0) {
1950 pix2-=pix1;
1951 pix1-=pix1;
1952 }
1953 if (pix2<0) {
1954 pix2-=pix2;
1955 pix1-=pix2;
1956 }
1957
1958 if ((cpl_array_get (bad_array, pix, NULL)==0)&&(cpl_array_get (bad_array, pix+1, NULL)==0))
1959 {
1960 if (pix1 < pix2 / 3)
1961 {
1962 cpl_array_set (bad_array, pix, BADPIX_DARK);
1963 count_bp_lowflux ++;
1964 }
1965 else if (pix2 < pix1 / 3)
1966 {
1967 cpl_array_set (bad_array, pix+1, BADPIX_DARK);
1968 count_bp_lowflux ++;
1969 }
1970 }
1971 }
1972 cpl_msg_info (cpl_func,"removed %.3f percents of FT pixels because of low flux",(100.0 * count_bp_lowflux)/ npix);
1973
1974 cpl_array_delete(flat_array_sum);
1975
1976 CPLCHECK_NUL ("Cannot use FT flats");
1977 } /* End case FLAT provided */
1978
1979 /* Set QC parameter */
1980 cpl_propertylist_append_int (bad_header, QC_BADPIX_FT, count_bp_dark);
1981 cpl_msg_info (cpl_func, "QC_BADPIX_FT = %d", count_bp_dark);
1982
1983 /* Set the badpixel map of FT in output data */
1984 gravi_data_add_table (bad_map, NULL, GRAVI_IMAGING_DATA_FT_EXT, bad_table);
1985
1986 CPLCHECK_NUL ("Cannot get the FT bad pixel mask");
1987 FREE (cpl_vector_delete, std_vector);
1988 } /* End FT case */
1989
1990 /*
1991 * Compute bad pixels of SC
1992 */
1993
1995 cpl_msg_warning (cpl_func,"The DARK map has no IMAGING_DATA_SC");
1996 }
1997 else
1998 {
1999 cpl_msg_info (cpl_func,"Compute BADPIXEL of SC");
2000
2001 /* Copy necessary tables */
2003
2004 /* This is the SC */
2005 cpl_image * dark_img = gravi_data_get_img (dark_map, GRAVI_IMAGING_DATA_SC_EXT);
2006 cpl_image * std_img = gravi_data_get_img (dark_map, GRAVI_IMAGING_ERR_SC_EXT);
2007 cpl_size nx = cpl_image_get_size_x (dark_img);
2008 cpl_size ny = cpl_image_get_size_y (dark_img);
2009 CPLCHECK_NUL ("Cannot get the SC data");
2010
2011 /* Compute an image with only the high-frequency of dark */
2012 cpl_image * darkhf_img = cpl_image_cast (dark_img, CPL_TYPE_DOUBLE);
2013 cpl_mask * kernel = cpl_mask_new (9, 9);
2014 cpl_mask_not (kernel);
2015 cpl_image_filter_mask (darkhf_img, dark_img, kernel, CPL_FILTER_MEDIAN, CPL_BORDER_FILTER);
2016 FREE (cpl_mask_delete, kernel);
2017
2018 cpl_image_subtract (darkhf_img, dark_img);
2019 cpl_image_multiply_scalar (darkhf_img, -1.0);
2020 CPLCHECK_NUL ("Cannot create darkhf");
2021
2022 /* Get the rms factor for dark bad pixel threshold */
2023 int bad_dark_factor = gravi_param_get_int (params, "gravity.calib.bad-dark-threshold");
2024
2025 /* Accepted range for DARK mean */
2026 double dark_rms = cpl_propertylist_get_double (dark_header, QC_DARKRMS_SC);
2027 double dark_mean = cpl_image_get_mean (darkhf_img);
2028 double dark_max = dark_mean + 2*bad_dark_factor * dark_rms;
2029 double dark_min = dark_mean - 2*bad_dark_factor * dark_rms;
2030
2031 /* Accepted range for dark STD */
2032 double dark_rms_std = cpl_image_get_stdev (std_img);
2033 double std_max = bad_dark_factor * dark_rms_std;
2034 double std_min = 0.05 * dark_rms_std;
2035
2036 /* Cannot detect bad pix on LOW mode */
2037 if (nx < 60) {
2038 cpl_msg_warning (cpl_func, "Don't detect SC bapixels in mode LOW");
2039 dark_max = 50000; //7; can't find a good param for slpit and combined
2040 dark_min = -5000;
2041 std_max *= 5;
2042 std_min *= 5; // FIXME: weird ??
2043 }
2044
2045 /* If we provide FLATs, we also use them. flat_img
2046 * will be filled with a filtered version of the FLAT */
2047 cpl_image * flat_img = NULL, * flatmed_img = NULL, * flatmask_img = NULL;
2048
2049 if (!nflat) {
2050 cpl_msg_info (cpl_func,"No FLATs provided, detect only DARK");
2051 }
2052 else {
2053 cpl_ensure (flats_data, CPL_ERROR_NULL_INPUT, NULL);
2054 cpl_msg_info (cpl_func,"FLATs used to detect badpix on SC");
2055
2056 /* Init co-add FLAT image */
2057 flat_img = cpl_image_new (nx,ny,CPL_TYPE_DOUBLE);
2058 cpl_image_fill_window (flat_img, 1, 1, nx, ny, 0.0);
2059
2060 /* Co-add FLATs, remove a self-estimate bias value
2061 * (the normal bias-pixels maybe wrong because of defocus),
2062 * assuming only half of the detector is illuminated */
2063 for (int f = 0; f<nflat; f++) {
2064 cpl_image * img = cpl_imagelist_collapse_median_create (
2065 gravi_data_get_cube (flats_data[f],
2067 cpl_image_subtract (img, dark_img);
2068
2069 double bias = gravi_image_get_quantile (img, 0.25);
2070 cpl_image_subtract_scalar (img, bias);
2071 cpl_image_add (flat_img, img);
2072 cpl_image_delete (img);
2073 CPLCHECK_NUL ("Cannot add flats");
2074 } /* End co-add flats */
2075
2076 /* Create kernel for horizontal median filtering
2077 * Note that the large cluster is ~11 pixel wide */
2078 cpl_size kernel_x = 5;
2079 if (nx > 100) kernel_x = 11;
2080 if (nx > 1000) kernel_x = 31;
2081 cpl_msg_info (cpl_func,"Kernel of (%lld,%i) pixels for median filtering", kernel_x, 1);
2082 cpl_mask * kernel = cpl_mask_new (kernel_x, 1);
2083 cpl_mask_not (kernel);
2084
2085 /* Run the median filter */
2086 flatmed_img = cpl_image_duplicate (flat_img);
2087 cpl_image_filter_mask (flatmed_img, flat_img, kernel, CPL_FILTER_MEDIAN, CPL_BORDER_FILTER);
2088 FREE (cpl_mask_delete, kernel);
2089
2090 /* Compute a map of illuminated pixels at 100 [adu] */
2091 double threshold = 100.0;
2092 flatmask_img = cpl_image_duplicate (flatmed_img);
2093 cpl_image_threshold (flatmask_img, threshold, threshold, 0, 1.0);
2094
2095 double fraction = cpl_image_get_flux (flatmask_img) / (nx*ny) * 100.0;
2096 cpl_msg_info (cpl_func, "Fraction of detector illuminated = %.2f%% (>%.1f adu)",
2097 fraction, threshold);
2098
2099 /* Add this image of illuminated part in BAD map */
2101 flatmask_img);
2102
2103 } /* End case FLAT provided */
2104
2105
2106 /* Compute the number of bad pixel */
2107 int count_bp_dark = 0;
2108 int count_bp_rms = 0;
2109 int count_bp_flat = 0;
2110 int count_bp = 0;
2111
2112 /* Create the badpixel image */
2113 cpl_image * bad_img = cpl_image_new (nx, ny, CPL_TYPE_INT);
2114
2115 /* Loop on pixels */
2116 for (cpl_size i = 0; i < nx; i++){
2117 for (cpl_size j = 0; j < ny; j++){
2118 int nv = 0, is_bad = 0, flag = 0;
2119
2120 /* Tag the bad pixel on its mean value ... */
2121 double dij = cpl_image_get (darkhf_img, i + 1, j + 1, &nv);
2122 if ( (dij > dark_max)||
2123 (dij < dark_min)) {
2124 flag += BADPIX_DARK;
2125 count_bp_dark ++;
2126 is_bad = 1;
2127 }
2128 /* ... and its variance */
2129 double stdij = cpl_image_get (std_img, i + 1, j + 1, &nv);
2130 if ( (stdij > std_max) ||
2131 (stdij < std_min) ) {
2132 flag += BADPIX_RMS;
2133 count_bp_rms ++;
2134 is_bad = 1;
2135 }
2136 /* ... and the flat */
2137 if (flat_img) {
2138 double flat = cpl_image_get (flat_img, i + 1, j + 1, &nv);
2139 double med = cpl_image_get (flatmed_img, i + 1, j + 1, &nv);
2140 double mask = cpl_image_get (flatmask_img, i + 1, j + 1, &nv);
2141 if ( flat < 0.5 * med && mask && i>4 && i<nx-4 ) {
2142 flag += BADPIX_FLAT;
2143 count_bp_flat ++;
2144 is_bad = 1;
2145 } }
2146
2147 /* Set this tag */
2148 count_bp += is_bad;
2149 cpl_image_set (bad_img, i + 1, j + 1, flag);
2150
2151 CPLCHECK_NUL ("Cannot compute bad pixel");
2152 }
2153 } /* End loop on pixels */
2154
2155 /* Set the badpixel map of SC in IMAGING_DATA_SC */
2156 gravi_data_add_img (bad_map, NULL, GRAVI_IMAGING_DATA_SC_EXT, bad_img);
2157
2158 /* Update QC parameter of main header */
2159 cpl_propertylist_append_int (bad_header, QC_BADPIX_SC, count_bp);
2160 cpl_msg_info (cpl_func, "QC_BADPIX_SC (total) = %d (%.2f%%)",
2161 count_bp, (100.0 * count_bp) / (nx*ny));
2162 cpl_propertylist_append_int (bad_header, QC_BADPIX_DARK_SC, count_bp_dark);
2163 cpl_msg_info (cpl_func, "QC_BADPIX_DARK_SC = %d", count_bp_dark);
2164 cpl_propertylist_append_int (bad_header, QC_BADPIX_RMS_SC, count_bp_rms);
2165 cpl_msg_info (cpl_func, "QC_BADPIX_RMS_SC = %d", count_bp_rms);
2166 cpl_propertylist_append_int (bad_header, QC_BADPIX_FLAT_SC, count_bp_flat);
2167 cpl_msg_info (cpl_func, "QC_BADPIX_FLAT_SC = %d", count_bp_flat);
2168
2169 FREE (cpl_image_delete, darkhf_img);
2170 FREE (cpl_image_delete, flat_img);
2171 FREE (cpl_image_delete, flatmed_img);
2172 } /* End SC case*/
2173
2174 /*
2175 * Compute bad pixels of ACQ
2176 */
2177
2179 cpl_msg_warning (cpl_func,"The DARK map has no IMAGING_DATA_ACQ");
2180 }
2181 else
2182 {
2183 cpl_msg_info (cpl_func,"Compute BADPIXEL of ACQ");
2184
2185 /* This is the SC */
2186 cpl_image * dark_img = gravi_data_get_img (dark_map, GRAVI_IMAGING_DATA_ACQ_EXT);
2187 CPLCHECK_NUL ("Cannot get dark");
2188
2189 /* Compute an image with only the high-frequency of dark */
2190 cpl_image * darkhf_img = cpl_image_cast (dark_img, CPL_TYPE_DOUBLE);
2191 cpl_mask * kernel = cpl_mask_new (21, 21);
2192 cpl_mask_not (kernel);
2193 cpl_image_filter_mask (darkhf_img, dark_img, kernel, CPL_FILTER_MEDIAN, CPL_BORDER_FILTER);
2194 FREE (cpl_mask_delete, kernel);
2195
2196 cpl_image_subtract (darkhf_img, dark_img);
2197 cpl_image_multiply_scalar (darkhf_img, -1.0);
2198 cpl_image_abs (darkhf_img);
2199 CPLCHECK_NUL ("Cannot create darkhf");
2200
2201 /* Get the STD of the dark_hf */
2202 double dark_std = cpl_image_get_median (darkhf_img);
2203 cpl_msg_info (cpl_func, "DARK_ACQ_STD = %e", dark_std);
2204
2205 /* Detect bad pixel as >5 STD */
2206 double threshold = 20 * dark_std + 1e-10;
2207 cpl_image_threshold (darkhf_img, threshold, threshold, 0, 1);
2208 cpl_image * bad_img = cpl_image_cast (darkhf_img, CPL_TYPE_INT);
2209
2210 /* Remove pixel in the middle of a square pattern */
2211 for (cpl_size x=2; x<cpl_image_get_size_x (bad_img);x++) {
2212 for (cpl_size y=2; y<cpl_image_get_size_y (bad_img);y++) {
2213 int nv = 0;
2214 if (cpl_image_get (bad_img,x,y-1,&nv) &&
2215 cpl_image_get (bad_img,x,y+1,&nv) &&
2216 cpl_image_get (bad_img,x-1,y,&nv) &&
2217 cpl_image_get (bad_img,x+1,y,&nv)) {
2218 cpl_image_set (bad_img,x,y, 1);
2219 }
2220 }
2221 }
2222
2223 /* Set QC parameter */
2224 cpl_size count_bp = 0;
2225 int nv = 0;
2226
2227 /* EKW 10/01/2019 Attention, I added brackets to avoid . [-Wmisleading-indentat */
2228 for (cpl_size x=0; x<cpl_image_get_size_x (bad_img);x++) {
2229 for (cpl_size y=0; y<cpl_image_get_size_y (bad_img);y++) {
2230 if (cpl_image_get (bad_img,x+1,y+1,&nv)) count_bp++;
2231 }
2232 }
2233
2234 cpl_propertylist_append_int (bad_header, "ESO QC BADPIX ACQ", count_bp);
2235
2236 /* Set the badpixel map of ACQ in IMAGING_DATA_ACQ */
2237 gravi_data_add_img (bad_map, NULL, GRAVI_IMAGING_DATA_ACQ_EXT, bad_img);
2238
2239 FREE (cpl_image_delete, darkhf_img);
2240 }
2241
2242 /* Verbose */
2244 return bad_map;
2245}
2246
2247
2248/*----------------------------------------------------------------------------*/
2269/*----------------------------------------------------------------------------*/
2270
2272 gravi_data ** flats_data,
2273 int nflat,
2274 const cpl_parameterlist * params)
2275{
2276 /* Verbose */
2278 cpl_ensure (dark_map, CPL_ERROR_NULL_INPUT, NULL);
2279 cpl_ensure (params, CPL_ERROR_NULL_INPUT, NULL);
2280 cpl_ensure (nflat==4 || nflat==0, CPL_ERROR_ILLEGAL_INPUT, NULL);
2281
2282 /* Construction of the bad pixels map */
2283 gravi_data * biasmask_map = gravi_data_new(0);
2284 cpl_propertylist * biasmask_header = gravi_data_get_header (biasmask_map);
2285
2286 /* Dump full header into product */
2287 cpl_propertylist * dark_header = gravi_data_get_header (dark_map);
2288 cpl_propertylist_append (biasmask_header, dark_header);
2289
2290 /* Copy necessary tables */
2291 gravi_data_copy_ext (biasmask_map, dark_map, GRAVI_IMAGING_DETECTOR_SC_EXT);
2292
2293 /* This is the SC */
2294 cpl_image * dark_img = gravi_data_get_img (dark_map, GRAVI_IMAGING_DATA_SC_EXT);
2295 cpl_size nx = cpl_image_get_size_x (dark_img);
2296 cpl_size ny = cpl_image_get_size_y (dark_img);
2297 CPLCHECK_NUL ("Cannot get the SC data");
2298
2299 /* Init co-add FLAT image */
2300 cpl_image * flat_img = cpl_image_new (nx,ny,CPL_TYPE_DOUBLE);
2301 cpl_image_fill_window (flat_img, 1, 1, nx, ny, 0.0);
2302
2303 /* Co-add FLATs, remove a self-estimate bias value
2304 * (the normal bias-pixels maybe wrong because of defocus),
2305 * assuming only half of the detector is illuminated */
2306 for (int f = 0; f<nflat; f++) {
2307 cpl_image * img = cpl_imagelist_collapse_median_create (
2308 gravi_data_get_cube (flats_data[f],
2310 cpl_image_subtract (img, dark_img);
2311
2312 double bias = gravi_image_get_quantile (img, 0.25);
2313 cpl_image_subtract_scalar (img, bias);
2314 cpl_image_add (flat_img, img);
2315 cpl_image_delete (img);
2316 CPLCHECK_NUL ("Cannot add flats");
2317 } /* End co-add flats */
2318
2319 /* Create kernel for horizontal median filtering
2320 * Note that the large cluster is ~11 pixel wide */
2321 cpl_size kernel_x = 5;
2322 if (nx > 100) kernel_x = 11;
2323 if (nx > 1000) kernel_x = 31;
2324 cpl_msg_info (cpl_func,"Kernel of (%lld,%i) pixels for median filtering", kernel_x, 1);
2325 cpl_mask * kernel = cpl_mask_new (kernel_x, 1);
2326 cpl_mask_not (kernel);
2327
2328 /* Run the median filter */
2329 cpl_image * flatmed_img = cpl_image_duplicate (flat_img);
2330 cpl_image_filter_mask (flatmed_img, flat_img, kernel, CPL_FILTER_MEDIAN, CPL_BORDER_FILTER);
2331 FREE (cpl_mask_delete, kernel);
2332
2333 /* Compute a map of illuminated pixels at 100 [adu] */
2334 double threshold = 100.0;
2335 cpl_image * biasmask_img = cpl_image_duplicate (flatmed_img);
2336 cpl_image_threshold (biasmask_img, threshold, threshold, 1.0, 0.0);
2337
2338 double fraction = cpl_image_get_flux (biasmask_img) / (nx*ny) * 100.0;
2339 cpl_msg_info (cpl_func, "Fraction of detector in bias mask = %.2f%% (<%.1f adu)",
2340 fraction, threshold);
2341
2342 /* Add this image of illuminated part in BAD map */
2343 gravi_data_add_img (biasmask_map, NULL, GRAVI_BIAS_MASK_SC_EXT,
2344 biasmask_img);
2345
2346 FREE (cpl_image_delete, flat_img);
2347 FREE (cpl_image_delete, flatmed_img);
2348
2349 /* Verbose */
2351 return biasmask_map;
2352}
2353
2354/*----------------------------------------------------------------------------*/
2374/*----------------------------------------------------------------------------*/
2375
2377 const cpl_parameterlist * params)
2378{
2379 /* Verbose */
2381 cpl_ensure (data, CPL_ERROR_NULL_INPUT, NULL);
2382 cpl_ensure (params, CPL_ERROR_NULL_INPUT, NULL);
2383
2384 /* Construction of the piezo tf product */
2385 gravi_data * piezo_tf = gravi_data_new(0);
2386 cpl_propertylist * piezotf_header = gravi_data_get_header (piezo_tf);
2387
2388 /* Dump full header into product */
2389 cpl_propertylist * data_header = gravi_data_get_header (data);
2390 cpl_propertylist_append (piezotf_header, data_header);
2391
2392 /* Copy necessary tables */
2393 gravi_data_copy_ext (piezo_tf, data, GRAVI_OPDC_EXT);
2394
2395 /* Get the OPDC table */
2396 cpl_table * opdc = gravi_data_get_table (piezo_tf, GRAVI_OPDC_EXT);
2397
2398 /* Get the necessary variables */
2399 int ntel = 4;
2400 int nresp = 5;
2401 int nsmooth = 201;
2402 char qc_name[100];
2403 char name[100];
2404 double phase;
2405 int base1,base2,base3,sign1,sign2,sign3;
2406 int ndit_small, dit_matrix;
2407 double twoPi = 2.0 * CPL_MATH_PI;
2408
2409
2410 /* DO THE COMPUTATION*/
2411
2412 // Determine the actuator to use, and associated gain, based on ESO.DPR.TYPE
2413
2414 char * actuator;
2415 const char * dpr_type = cpl_propertylist_get_string(data_header, "ESO DPR TYPE");
2416 if (strcmp(dpr_type, "PIEZOTF") == 0)
2417 actuator = cpl_strdup("PIEZO_DL_OFFSET");
2418 else if (strcmp(dpr_type, "VLTITF") == 0)
2419 actuator = cpl_strdup("VLTI_DL_OFFSET");
2420 else
2421 {
2422 cpl_msg_error(cpl_func, "ESO DPR TYPE = %s is not supported!!!", dpr_type);
2423 return NULL;
2424 }
2425 cpl_msg_info(cpl_func, "Using %s actuator (DPR.TYPE=%s)", actuator, dpr_type);
2426
2427 // Read data array from OPDC table
2428
2429 cpl_size ndit = cpl_table_get_nrow (opdc);
2430 ndit_small=ndit-nsmooth-(nresp+1);
2431 cpl_msg_info (cpl_func, "Preparing matrix inversion with NDIT = %lld",ndit);
2432 cpl_array** opd = cpl_table_get_data_array (opdc, "OPD");
2433 cpl_array** piezo = cpl_table_get_data_array (opdc, actuator);
2434 CPLCHECK_NUL ("Cannot read the OPDC data");
2435
2436 // create temporary arrays
2437
2438 cpl_array * phase_array = cpl_array_new(GRAVI_NBASE, CPL_TYPE_DOUBLE);
2439 cpl_array * piezo_array = cpl_array_new(ntel, CPL_TYPE_DOUBLE);
2440 cpl_matrix * phase_matrix = cpl_matrix_new (ndit_small*GRAVI_NBASE,1);
2441 cpl_matrix * piezo_matrix = cpl_matrix_new (ndit_small*GRAVI_NBASE, nresp*ntel);
2442 cpl_matrix * piezo_header_resp = cpl_matrix_new (nresp*ntel,1);
2443
2444 // Unwrap phase of OPD
2445 for (cpl_size dit = 0 ; dit < ndit-1 ; dit ++)
2446 for (cpl_size base = 0 ; base < GRAVI_NBASE; base ++)
2447 {
2448 phase=cpl_array_get(opd[dit+1],base, NULL)-cpl_array_get(opd[dit],base, NULL);
2449 phase-= twoPi * floor( phase / twoPi );
2450 if (phase > CPL_MATH_PI) phase -= twoPi;
2451 cpl_array_set(opd[dit+1],base, phase+cpl_array_get(opd[dit],base, NULL));
2452 }
2453
2454 // filter OPDs, Commands and create Matrixes
2455
2456 for (cpl_size dit = 0 ; dit < ndit-nsmooth; dit ++)
2457 {
2458 // Filter OPDs
2459 cpl_array_fill_window (phase_array, 0, GRAVI_NBASE, 0.0);
2460 for (cpl_size smooth = 0 ; smooth < nsmooth; smooth ++)
2461 cpl_array_add( phase_array, opd[smooth+dit] );
2462 cpl_array_multiply_scalar (phase_array, -1.0 / nsmooth);
2463 cpl_array_add( phase_array, opd[ (int) (nsmooth/2+dit) ] );
2464
2465 // Filter PIEZO commands
2466 cpl_array_fill_window_double (piezo_array, 0, ntel, 0.0);
2467 for (cpl_size smooth = 0 ; smooth < nsmooth; smooth ++)
2468 cpl_array_add( piezo_array, piezo[smooth+dit] );
2469 cpl_array_multiply_scalar (piezo_array, -1.0 / nsmooth);
2470 cpl_array_add( piezo_array, piezo[ (int) (nsmooth/2+dit) ] );
2471
2472 // Store values into phase matrix for SVD inversion
2473 if (( dit-(nresp+1) >= 0 )&( dit-(nresp+1) < ndit_small))
2474 for (cpl_size base = 0 ; base < GRAVI_NBASE; base ++)
2475 {
2476 int dit_matrix =dit-(nresp+1)+base*ndit_small;
2477 cpl_matrix_set(phase_matrix,dit_matrix,0,cpl_array_get(phase_array, base, NULL));
2478 }
2479
2480 // Store values into piezo matrix for SVD inversion
2481 for (cpl_size tel = 0 ; tel < ntel; tel ++)
2482 {
2483 switch (tel)
2484 {
2485 default:
2486 case 0:
2487 base1=0;
2488 sign1=-1;
2489 base2=1;
2490 sign2=-1;
2491 base3=2;
2492 sign3=-1;
2493 break;
2494 case 1:
2495 base1=0;
2496 sign1=1;
2497 base2=3;
2498 sign2=-1;
2499 base3=4;
2500 sign3=-1;
2501 break;
2502 case 2:
2503 base1=1;
2504 sign1=1;
2505 base2=3;
2506 sign2=1;
2507 base3=5;
2508 sign3=-1;
2509 break;
2510 case 3:
2511 base1=2;
2512 sign1=1;
2513 base2=4;
2514 sign2=1;
2515 base3=5;
2516 sign3=1;
2517 break;
2518 }
2519
2520 for (cpl_size resp = 0 ; resp < nresp; resp ++)
2521 if (( dit-(nresp+1)+(1+resp) >= 0 )&( dit-(nresp+1)+(1+resp) < ndit_small))
2522 {
2523 dit_matrix =dit-(nresp+1)+(1+resp)+base1*ndit_small;
2524 cpl_matrix_set(piezo_matrix,dit_matrix,resp*ntel+tel,sign1*cpl_array_get(piezo_array, tel, NULL));
2525 dit_matrix =dit-(nresp+1)+(1+resp)+base2*ndit_small;
2526 cpl_matrix_set(piezo_matrix,dit_matrix,resp*ntel+tel,sign2*cpl_array_get(piezo_array, tel, NULL));
2527 dit_matrix =dit-(nresp+1)+(1+resp)+base3*ndit_small;
2528 cpl_matrix_set(piezo_matrix,dit_matrix,resp*ntel+tel,sign3*cpl_array_get(piezo_array, tel, NULL));
2529 }
2530 }
2531 }
2532 CPLCHECK_NUL ("Cannot create matrix for SVD inversion");
2533
2534 // resolve SVD
2535 cpl_msg_info (cpl_func, "Doing SVD inversion" );
2536 cpl_matrix * piezo_resp = cpl_matrix_solve_normal(piezo_matrix,phase_matrix); // coef_vis is 20x1
2537 cpl_matrix * residuals_fit = cpl_matrix_product_create(piezo_matrix,piezo_resp);
2538 cpl_matrix_subtract (residuals_fit,phase_matrix); //residuals_fit is (ndit-5)*GRAVI_NBASE
2539 CPLCHECK_NUL ("Failed to do SVD inversion");
2540
2541 // Scale DL response to arbitrary PIEZO response of 17.4 rad/Volt
2542 if (strcmp(dpr_type, "VLTITF") == 0)
2543 {
2544 for (cpl_size tel = 0 ; tel < ntel; tel ++)
2545 {
2546 double gain = 0.0;
2547 for (cpl_size resp = 0 ; resp < nresp; resp ++)
2548 gain += cpl_matrix_get (piezo_resp, resp*ntel+tel, 0);
2549 gain /= 17.4;
2550 for (cpl_size resp = 0 ; resp < nresp; resp ++)
2551 cpl_matrix_set (piezo_resp, resp*ntel+tel, 0, cpl_matrix_get (piezo_resp, resp*ntel+tel, 0) / gain );
2552 }
2553 }
2554
2555 // Get piezo response from header
2556
2557 for (cpl_size tel = 0 ; tel < ntel; tel ++)
2558 for (cpl_size resp = 0 ; resp < nresp; resp ++)
2559 {
2560 sprintf (name, "ESO FT KAL P%lld_RESP%lld", tel+1, resp+1);
2561 if (cpl_propertylist_has (piezotf_header, name))
2562 cpl_matrix_set( piezo_header_resp, resp*ntel+ tel, 0, cpl_propertylist_get_double (piezotf_header, name));
2563 else
2564 cpl_matrix_set( piezo_header_resp, resp*ntel+ tel, 0, 0.0);
2565 }
2566
2567 // output QC parameters
2568
2569 sprintf (qc_name, "ESO QC FT KAL P_FIT");
2570 cpl_propertylist_update_double (piezotf_header, qc_name, cpl_matrix_get_stdev( residuals_fit ) );
2571 cpl_propertylist_set_comment (piezotf_header, qc_name, "Fitting standard deviation [rad]");
2572 cpl_msg_info (cpl_func, "Fit standard deviation = %e [rad]", cpl_matrix_get_stdev( residuals_fit ) );
2573
2574
2575 sprintf (name, "ESO FT RATE");
2576 double sampling = cpl_propertylist_get_double (piezotf_header, name);
2577 char ft_name[100];
2578
2579 for (cpl_size tel = 0 ; tel < ntel; tel ++)
2580 {
2581 for (cpl_size resp = 0 ; resp < nresp; resp ++)
2582 {
2583 sprintf (qc_name, "ESO QC FT KAL P%lld_RESP%lld", tel+1, resp+1);
2584 cpl_propertylist_update_double (piezotf_header, qc_name, cpl_matrix_get( piezo_resp, resp*ntel+ tel,0 ) );
2585 cpl_propertylist_set_comment (piezotf_header, qc_name, "Computed Kalman response");
2586
2587 cpl_msg_info (cpl_func, "QC FT KAL P%lld_RESP%lld = %5.5g [rad/Volts]", tel+1, resp+1, cpl_matrix_get( piezo_resp, resp*ntel+ tel,0 ));
2588
2589 sprintf (qc_name, "ESO QC FT KAL P%lld_RESP%lld DIFF", tel+1, resp+1);
2590 sprintf (ft_name, "ESO FT KAL P%lld_RESP%lld", tel+1, resp+1);
2591 double ft_kal_diff = cpl_matrix_get( piezo_resp, resp*ntel+ tel,0 ) -
2592 cpl_propertylist_get_double (data_header, ft_name);
2593 cpl_propertylist_update_double (piezotf_header, qc_name, ft_kal_diff);
2594 cpl_propertylist_set_comment (piezotf_header, qc_name, "Comp-Nom Kalman resp");
2595 }
2596
2597 // Get gain in rad/Volts
2598
2599 double QC_gain=0.0;
2600 for (cpl_size resp = 0 ; resp < nresp; resp ++)
2601 QC_gain+=cpl_matrix_get( piezo_resp, resp*ntel + tel,0 );
2602
2603 sprintf (qc_name, "ESO QC FT KAL P%lld_GAIN", tel+1);
2604 cpl_propertylist_update_double (piezotf_header, qc_name, QC_gain );
2605 cpl_propertylist_set_comment (piezotf_header, qc_name, "Open loop gain [rad/Volts]");
2606
2607 // Get pur delay in milliseconds
2608
2609 double QC_delay=0.0;
2610 for (cpl_size resp = 0 ; resp < nresp; resp ++)
2611 QC_delay+=(resp+1)*cpl_matrix_get( piezo_resp, resp*ntel + tel,0 )*sampling;
2612 QC_delay/=QC_gain;
2613
2614 sprintf (qc_name, "ESO QC FT KAL P%lld_DELAY", tel+1);
2615 cpl_propertylist_update_double (piezotf_header, qc_name, QC_delay );
2616 cpl_propertylist_set_comment (piezotf_header, qc_name, "Open loop latency [ms]");
2617
2618 // Get standard deviation
2619 double QC_std=0.0;
2620 for (cpl_size resp = 0 ; resp < nresp; resp ++)
2621 QC_std+=(cpl_matrix_get( piezo_resp, resp*ntel + tel,0 ) - cpl_matrix_get( piezo_header_resp, resp*ntel+ tel, 0))
2622 * (cpl_matrix_get( piezo_resp, resp*ntel + tel,0 ) - cpl_matrix_get( piezo_header_resp, resp*ntel+ tel, 0));
2623
2624 sprintf (qc_name, "ESO QC FT KAL P%lld_STDEV", tel+1);
2625 cpl_propertylist_update_double (piezotf_header, qc_name, sqrt(QC_std)/sqrt(nresp) );
2626 cpl_propertylist_set_comment (piezotf_header, qc_name, "Stdev of RTC [radians]");
2627
2628 }
2629 CPLCHECK_NUL ("Failed to generate and store QC parameters");
2630
2631 // delete disposable arrays
2632 cpl_array_delete(phase_array);
2633 cpl_array_delete(piezo_array);
2634 cpl_matrix_delete(phase_matrix);
2635 cpl_matrix_delete(piezo_matrix);
2636 cpl_matrix_delete(piezo_header_resp);
2637 cpl_matrix_delete(piezo_resp);
2638 cpl_matrix_delete(residuals_fit);
2639 cpl_free(actuator);
2640
2641 /* Verbose */
2643 return piezo_tf;
2644}
2645
2646
2647/*---------------------------------------------------------------------------*/
2659/*---------------------------------------------------------------------------*/
2660
2661cpl_error_code gravi_remove_cosmicrays_sc (cpl_imagelist * imglist_sc)
2662{
2664 cpl_ensure_code (imglist_sc, CPL_ERROR_NULL_INPUT);
2665
2666 cpl_image * img;
2667 double * img_ptr;
2668 cpl_binary * bpm_ptr;
2669
2670 const cpl_size nrow = cpl_imagelist_get_size (imglist_sc);
2671 img = cpl_imagelist_get (imglist_sc, 0);
2672 const cpl_size nx = cpl_image_get_size_x (img);
2673 const cpl_size ny = cpl_image_get_size_y (img);
2674
2675 /* default clip is 10 sigma, but to be increased in case of small number of images */
2676 double clip_thresh = 5.;
2677
2678 if (nrow <= 32) clip_thresh = 10;
2679 if (nrow <= 16) clip_thresh = 15;
2680 if (nrow <= 12) clip_thresh = 20;
2681 if (nrow <= 8) clip_thresh = 30;
2682 if (nrow <= 4) clip_thresh = 40;
2683 if (nrow < 4) clip_thresh = 50;
2684
2685 cpl_msg_info (cpl_func, "Number of images is %lld => cosmic ray detection threshold set to %f sigma", nrow, clip_thresh);
2686
2687 /* Declare arrays and get pointers for speed */
2688
2689 /* Median and SD across images */
2690 cpl_array * med_val = cpl_array_new (nx, CPL_TYPE_DOUBLE);
2691 double * med_val_ptr = cpl_array_get_data_double (med_val);
2692 cpl_array * std_val = cpl_array_new (nx, CPL_TYPE_DOUBLE);
2693 double * std_val_ptr = cpl_array_get_data_double (std_val);
2694
2695 /* Pixel values at fixed (x,y) for each frame */
2696 cpl_vector * row_val = cpl_vector_new (nrow);
2697 double * row_ptr = cpl_vector_get_data (row_val);
2698
2699 /* Indices of valid and CR pixels */
2700 cpl_array * good_x = cpl_array_new (nx, CPL_TYPE_DOUBLE);
2701 double * good_x_ptr = cpl_array_get_data_double (good_x);
2702 cpl_array * good_y = cpl_array_new (nx, CPL_TYPE_DOUBLE);
2703 double * good_y_ptr = cpl_array_get_data_double (good_y);
2704 cpl_array * CR_x = cpl_array_new (nx, CPL_TYPE_DOUBLE);
2705 double * CR_x_ptr = cpl_array_get_data_double (CR_x);
2706
2707 /* Work arrays for interpolation of CR pixels */
2708 /* Allocated to maximum possible length, actual number of valid values will vary */
2709 cpl_vector * xref = cpl_vector_new(nx + 2);
2710 double * xref_ptr = cpl_vector_get_data(xref);
2711 cpl_vector * yref = cpl_vector_new(nx + 2);
2712 double * yref_ptr = cpl_vector_get_data(yref);
2713 cpl_vector * xout = cpl_vector_new(nx);
2714 double * xout_ptr = cpl_vector_get_data(xout);
2715 cpl_vector * yout = cpl_vector_new(nx);
2716 double * yout_ptr = cpl_vector_get_data(yout);
2717
2718 /* Merged work arrays for bivector interpolation */
2719 cpl_bivector * fref = cpl_bivector_wrap_vectors (xref, yref);
2720 cpl_bivector * fout = cpl_bivector_wrap_vectors (xout, yout);
2721
2722 /* count number of CR pixels for log message */
2723 cpl_vector * cCR_vector = cpl_vector_new (nrow);
2724 cpl_vector_fill (cCR_vector, 0.0);
2725 double *cCR_vector_ptr = cpl_vector_get_data(cCR_vector);
2726
2727 /* macro for indexing 2d array */
2728 #define IMAGE_IDX(_nx, _x, _y) (_x + _nx * _y)
2729
2730 /* loop through all image rows of the image */
2731 for (cpl_size k = 0; k < ny; k++) {
2732 /* loop through all pixels in the image row */
2733 for (cpl_size i = 0; i < nx; i++) {
2734 /* find pixel value across all images */
2735 for (cpl_size row = 0; row < nrow; row++) {
2736 img_ptr = cpl_image_get_data (cpl_imagelist_get (imglist_sc, row));
2737 cpl_size idx = IMAGE_IDX (nx, i, k);
2738 row_ptr[row] = img_ptr[idx];
2739 } /* End image loop */
2740
2741 /* calculate pixel mean and standard deviation */
2742 double median = cpl_vector_get_median (row_val);
2743 med_val_ptr[i] = median;
2744 cpl_vector_subtract_scalar (row_val, median);
2745 cpl_vector_multiply (row_val, row_val);
2746 cpl_vector_sqrt (row_val); /* = abs (val) */
2747 std_val_ptr[i] = cpl_vector_get_median (row_val) * CPL_MATH_STD_MAD;
2748 } /* End pixel loop */
2749
2750 /* loop through the values for each image and identify the outliers */
2751 /* this is done per image per row */
2752
2753 for (cpl_size row = 0; row < nrow; row++) {
2754 cpl_size nGood = 0, nCR = 0;
2755
2756 /* load image */
2757 img = cpl_imagelist_get (imglist_sc, row);
2758 img_ptr = cpl_image_get_data_double (img);
2759 bpm_ptr = cpl_mask_get_data (cpl_image_get_bpm (img));
2760
2761 /* Separate good pixels from those with CRs */
2762 for (cpl_size i = 0; i < nx; i++) {
2763 cpl_size idx = IMAGE_IDX (nx, i, k);
2764 cpl_binary nv = bpm_ptr[idx];
2765 double val = img_ptr[idx];
2766 if ((nv == 1) || (val > med_val_ptr[i] + clip_thresh * std_val_ptr[i])) {
2767 CR_x_ptr[nCR++] = i;
2768 } else {
2769 good_x_ptr[nGood] = i;
2770 good_y_ptr[nGood++] = val;
2771 }
2772 } /* End column loop */
2773
2774 /* add counter of CR to vector */
2775 cCR_vector_ptr[row] += nCR;
2776 //cpl_msg_warning (cpl_func,"TT Cosmic rays detected: %lli", nCR);
2777
2778 /* interpolate CR affected pixels */
2779 if (nCR > 0) {
2780 /* clear vectors */
2781 cpl_vector_fill(xref, 0.0);
2782 cpl_vector_fill(yref, 0.0);
2783 cpl_vector_fill(xout, 0.0);
2784 cpl_vector_fill(yout, 0.0);
2785
2786 /* fill vectors with good and CR pixels*/
2787 for (cpl_size i = 0; i < nGood; i++) {
2788 xref_ptr[i+1] = good_x_ptr[i];
2789 yref_ptr[i+1] = good_y_ptr[i];
2790 } /* End nGood loop */
2791
2792 /* Fix the non-extrapolation inability of cpl_bivector_interpolate_linear */
2793 xref_ptr[0] = 0;
2794 yref_ptr[0] = good_y_ptr[0];
2795 for (int m = nGood+1; m < nx + 2; m++) {
2796 xref_ptr[m] = nx;
2797 yref_ptr[m] = good_y_ptr[nGood-1];
2798 }
2799
2800 for (cpl_size i = 0; i < nCR; i++) {
2801 xout_ptr[i] = CR_x_ptr[i];
2802 } /* End nCR loop */
2803
2804 /* interpolate CR positions */
2805 cpl_bivector_interpolate_linear (fout, fref);
2806 CPLCHECK_MSG ("Cannot interpolate CR pixels!");
2807
2808 /* replace CR pixels in image with interpolated values */
2809 for (cpl_size i = 0; i < nCR; i++) {
2810 cpl_size idx = IMAGE_IDX (nx, xout_ptr[i], k);
2811 img_ptr[idx] = yout_ptr[i];
2812 bpm_ptr[idx] = CPL_BINARY_1;
2813 }
2814 } /* End IF nCR*/
2815 } /* End loop through all images */
2816
2817 } /* End loop through all image rows */
2818
2819 /* print the number of flagged pixels to the log */
2820 double percentage_CR_perpixel = cpl_vector_get_mean (cCR_vector) * 100 / (nx * ny);
2821 if (percentage_CR_perpixel > 1)
2822 cpl_msg_warning (cpl_func,"Cosmic rays flagged on %g percents of the pixels", percentage_CR_perpixel);
2823 else
2824 cpl_msg_info (cpl_func,"Cosmic rays flagged on %g percent of the pixels", percentage_CR_perpixel);
2825 for (cpl_size row = 0; row < nrow; row++)
2826 if (cpl_vector_get (cCR_vector, row) > (nx * ny) / 1000)
2827 cpl_msg_warning (cpl_func,"Cosmic rays detected on image %lli: %g", row+1, cpl_vector_get (cCR_vector, row));
2828
2829 /* Delete temporary arrays */
2830 FREE (cpl_array_delete, med_val);
2831 FREE (cpl_array_delete, std_val);
2832 FREE (cpl_vector_delete, row_val);
2833
2834 FREE (cpl_array_delete, good_x);
2835 FREE (cpl_array_delete, good_y);
2836 FREE (cpl_array_delete, CR_x);
2837 FREE (cpl_vector_delete, cCR_vector);
2838
2839 FREE(cpl_vector_delete, xref);
2840 FREE(cpl_vector_delete, yref);
2841 FREE(cpl_vector_delete, xout);
2842 FREE(cpl_vector_delete, yout);
2843
2844 FREE (cpl_bivector_unwrap_vectors, fref);
2845 FREE (cpl_bivector_unwrap_vectors, fout);
2846
2848 return CPL_ERROR_NONE;
2849#undef IMAGE_IDX
2850}
2851
2852/*----------------------------------------------------------------------------*/
2853
2862cpl_array *gravi_filter_array_median(const cpl_array *arr, int size)
2863{
2864 cpl_size n = cpl_array_get_size(arr);
2865 cpl_array *result_arr = cpl_array_new(n, CPL_TYPE_DOUBLE);
2866 cpl_array *window_arr = cpl_array_new(size, CPL_TYPE_DOUBLE);
2867 const double *source_data = cpl_array_get_data_double_const(arr);
2868
2869 for (int i = 0; i < n; i++) { // over wl
2870 // Use reflecting boundary conditions
2871 // d c b a | a b c d | d c b a
2872 for (int w = 0; w < size; w++) { // over wl, with bcs
2873 int v = i + w - size / 2;
2874 v = (v < 0) ? -(v + 1) : ((v > n - 1) ? 2 * n - 1 - v : v);
2875 cpl_array_set_double(window_arr, w, source_data[v]);
2876 }
2877 cpl_array_set_double(result_arr, i, cpl_array_get_median(window_arr));
2878 }
2879 cpl_array_delete(window_arr);
2880 return result_arr;
2881}
2882
2891double gravi_calc_sigmaclipped_stddev(const cpl_array *arr, double nstd)
2892{
2893 cpl_ensure(arr, CPL_ERROR_NULL_INPUT, 0.0);
2894 cpl_ensure(nstd > 0.0, CPL_ERROR_ILLEGAL_INPUT, 0.0);
2895
2896 int n = cpl_array_get_size(arr);
2897 cpl_array *work_arr = cpl_array_duplicate(arr);
2898
2899 int n_changed = 0, iterations = 0;
2900 const int maxiters = 5;
2901 do {
2902 n_changed = 0;
2903 double centre = cpl_array_get_median(work_arr);
2904 double std = cpl_array_get_stdev(work_arr);
2905
2906 double lower_bound = centre - nstd * std;
2907 double upper_bound = centre + nstd * std;
2908
2909 for (int i = 0; i < n; i++) {
2910 int invalid;
2911 double v = cpl_array_get_double(work_arr, i, &invalid);
2912 if (!invalid && (v < lower_bound || v > upper_bound)) {
2913 n_changed++;
2914 cpl_array_set_invalid(work_arr, i);
2915 }
2916 }
2917 if (++iterations > maxiters) break;
2918 } while(n_changed > 0);
2919 double std = cpl_array_get_stdev(work_arr);
2920 cpl_array_delete(work_arr);
2921 return std;
2922}
2923
2944/*----------------------------------------------------------------------------*/
2946 int naccept,
2947 const cpl_parameterlist * params)
2948{
2949 /* Verbose */
2951 cpl_ensure(data, CPL_ERROR_NULL_INPUT, NULL);
2952 cpl_ensure(naccept, CPL_ERROR_DATA_NOT_FOUND, NULL);
2953 cpl_ensure(params, CPL_ERROR_NULL_INPUT, NULL);
2954
2955 int npol, nwave = 0;
2956
2957 int num_components = cpl_parameter_get_int(
2958 cpl_parameterlist_find_const(params, "gravity.calib.pca-components"));
2959
2960 int median_filter_size = cpl_parameter_get_int(
2961 cpl_parameterlist_find_const(params, "gravity.calib.pca-clean-size"));
2962
2963 double median_filter_nstd = cpl_parameter_get_double(
2964 cpl_parameterlist_find_const(params, "gravity.calib.pca-clean-nstd"));
2965
2966 const char *fit_type = cpl_parameter_get_string(
2967 cpl_parameterlist_find_const(params, "gravity.calib.pca-fit-type"));
2968
2969 int fit_degree = cpl_parameter_get_int(
2970 cpl_parameterlist_find_const(params, "gravity.calib.pca-fit-degree"));
2971
2972 cpl_boolean save_residuals = cpl_parameter_get_bool(
2973 cpl_parameterlist_find_const(params, "gravity.calib.pca-save-residuals"));
2974
2975 cpl_propertylist *hdr = NULL, *wave_plist = NULL;
2976 cpl_table *wave_table = NULL, *vis_table = NULL, *pca_table = NULL;
2977
2978 const char *telescope = NULL, *pola_mode = NULL, *spec_res = NULL;
2979 char colname_vis[100], colname_pca[100], colname_pca_fit[100], colname_mean[100], colname_residual[100];
2980
2981 /* Get header data */
2982 hdr = gravi_data_get_header(data[0]);
2983 telescope = cpl_propertylist_get_string(hdr, "TELESCOP");
2984 pola_mode = gravi_pfits_get_pola_mode(hdr, GRAVI_SC);
2986 spec_res = gravi_pfits_get_spec_res(hdr);
2987
2988 /* Get size of wavelength axis */
2989 wave_plist = gravi_data_get_oi_wave_plist(data[0], GRAVI_SC, 0, npol);
2990 nwave = cpl_propertylist_get_int(wave_plist, "NWAVE");
2991
2992 /* Get wavelength and cast to double */
2993 wave_table = gravi_data_get_oi_wave(data[0], GRAVI_SC, 0, npol);
2994 cpl_table_cast_column(wave_table, "EFF_WAVE", "EFF_WAVE", CPL_TYPE_DOUBLE);
2995 cpl_table_multiply_scalar(wave_table, "EFF_WAVE", 1.0e6);
2996 cpl_array *wave_arr = cpl_array_wrap_double(cpl_table_get_data_double(wave_table, "EFF_WAVE"), nwave);
2997
2998 /* Prepare table for collated visphi */
2999 vis_table = cpl_table_new(naccept); // One row per result - contains data which is different for every input frame
3000 for(int i = 0; i < GRAVI_NBASE; i++) {
3001 for(int j = 0; j < npol; j++) {
3002 sprintf(colname_vis, "PCA_VISPHI_BASE%d_POL%d", i, j);
3003 cpl_table_new_column_array(vis_table, colname_vis, CPL_TYPE_DOUBLE, nwave);
3004 }
3005 }
3006
3007 /* Collate visphi and apply median smoothing */
3008 for (int n = 0; n < naccept; n++) {
3009 hdr = gravi_data_get_header(data[n]);
3010
3011 for (int j = 0; j < npol; j++) {
3012 cpl_table *vis_tmp = gravi_data_get_oi_vis(data[n], GRAVI_SC, j, npol);
3013 for (int i = 0; i < GRAVI_NBASE; i++) {
3014 /* Copy the visphi data */
3015 cpl_array *vis_arr = cpl_array_duplicate(cpl_table_get_array(vis_tmp, "VISPHI", i));
3016
3017 /* Calculate median-smoothed visphi */
3018 cpl_array *vis_filt = gravi_filter_array_median(vis_arr, median_filter_size);
3019
3020 /* Calculate standard deviation using sigma-clipping */
3021 cpl_array *vis_centred = cpl_array_duplicate(vis_arr);
3022 cpl_array_subtract(vis_centred, vis_filt); // vis_centred = arr - med_filt
3023 double std = gravi_calc_sigmaclipped_stddev(vis_centred, 3.0);
3024
3025 /* Replace outlier data with the filtered median */
3026 for (int k = 0; k < nwave; k++) {
3027 if (fabs(cpl_array_get(vis_centred, k, NULL)) > median_filter_nstd * std) {
3028 cpl_array_set(vis_arr, k, cpl_array_get(vis_filt, k, NULL));
3029 }
3030 }
3031
3032 /* Store data with outliers removed */
3033 sprintf(colname_vis, "PCA_VISPHI_BASE%d_POL%d", i, j);
3034 cpl_table_set_array(vis_table, colname_vis, n, vis_arr);
3035
3036 cpl_array_delete(vis_filt);
3037 cpl_array_delete(vis_centred);
3038 cpl_array_delete(vis_arr);
3039 }
3040 }
3041 }
3042
3043 /* Calculate PCA decomposition */
3044 pca_table = cpl_table_new(nwave); // One row per wavelength - contains aggregate data
3045 gravi_pca_result **pca_decomps = cpl_malloc(GRAVI_NBASE * npol * sizeof(gravi_pca_result *));
3046
3047 for (int i = 0; i < GRAVI_NBASE; i++) {
3048 for (int j = 0; j < npol; j++) {
3049 sprintf(colname_vis, "PCA_VISPHI_BASE%d_POL%d", i, j);
3050
3051 /* Copy data into matrix of size (naccept, nwave) */
3052 cpl_matrix *vis_tmp = cpl_matrix_new(naccept, nwave);
3053 for (int n = 0; n < naccept; n++) {
3054 const cpl_array *vis_arr_tmp = cpl_table_get_array(vis_table, colname_vis, n);
3055 for (int w = 0; w < nwave; w++)
3056 cpl_matrix_set(vis_tmp, n, w, cpl_array_get(vis_arr_tmp, w, NULL));
3057 }
3058
3059 /* Perform decomposition */
3060 gravi_pca_result *decomp = gravi_pca_create_result(vis_tmp, /* mask */ NULL);
3062 gravi_pca_set_component_signs(decomp, wave_arr, num_components);
3063 pca_decomps[i + GRAVI_NBASE * j] = decomp;
3064
3065 /* Extract mean values */
3066 cpl_vector *mean = gravi_pca_get_component(pca_decomps[i + GRAVI_NBASE * j], 0);
3067 sprintf(colname_mean, "PCA_MEAN_BASE%d_POL%d", i, j);
3068 cpl_table_new_column(pca_table, colname_mean, CPL_TYPE_DOUBLE);
3069 cpl_table_copy_data_double(pca_table, colname_mean, cpl_vector_get_data_const(mean));
3070
3071 /* Extract components */
3072 for (int c = 0; c < num_components; c++) {
3073 sprintf(colname_pca, "PCA_C%d_BASE%d_POL%d", c+1, i, j);
3074 cpl_vector *cv = gravi_pca_get_component(decomp, c+1);
3075 cpl_table_new_column(pca_table, colname_pca, CPL_TYPE_DOUBLE);
3076 cpl_table_copy_data_double(pca_table, colname_pca, cpl_vector_get_data_const(cv));
3077 cpl_vector_delete(cv);
3078 }
3079
3080 cpl_vector_delete(mean);
3081 cpl_matrix_delete(vis_tmp);
3082 }
3083 }
3084
3085 /* Calculate median-averaged PCA components */
3087 (const gravi_pca_result**)pca_decomps, GRAVI_NBASE * npol, num_components);
3088 CPLCHECK_NUL("Failed to compute median averaged PCA components");
3089
3090 /* Store median component values */
3091 for (int c = 0; c < num_components; c++) {
3092 sprintf(colname_pca, "PCA_C%d_MEDIAN", c+1);
3093 cpl_vector *cmed = gravi_pca_get_component_median(pca_model, c+1);
3094 cpl_table_new_column(pca_table, colname_pca, CPL_TYPE_DOUBLE);
3095 cpl_table_copy_data_double(pca_table, colname_pca, cpl_vector_get_data(cmed));
3096 cpl_vector_delete(cmed);
3097 }
3098
3099 /* Compute fit to median-averaged components (either polynomial or B-spline) */
3100 cpl_error_code err;
3101 if (strcmp(fit_type, "POLYNOMIAL") == 0) {
3102 err = gravi_pca_fit_components_polynomial(pca_model, wave_arr, fit_degree, num_components);
3103 if (err) cpl_error_set(cpl_func, err);
3104 } else {
3105 err = gravi_pca_fit_components_bspline(pca_model, wave_arr, fit_degree, num_components);
3106 if (err) cpl_error_set(cpl_func, err);
3107 }
3108 CPLCHECK_NUL("Failed to fit median PCA components");
3109
3110 /* Store fit values */
3111 for (int c = 0; c < num_components; c++) {
3112 cpl_vector *fit_vals = gravi_pca_get_component_fit(pca_model, c+1);
3113 sprintf(colname_pca_fit, "PCA_C%d_FIT", c+1);
3114 cpl_table_new_column(pca_table, colname_pca_fit, CPL_TYPE_DOUBLE);
3115 cpl_table_copy_data_double(pca_table, colname_pca_fit, cpl_vector_get_data_const(fit_vals));
3116 cpl_vector_delete(fit_vals);
3117 }
3118
3119 /* Construct noise-free model from linear combination of components and fit to data */
3120 for (int i = 0; i < GRAVI_NBASE; i++) {
3121 for (int j = 0; j < npol; j++) {
3122 err = gravi_pca_fit_model(pca_decomps[i + GRAVI_NBASE * j], pca_model,
3123 /* fit_mean_subtracted */ CPL_FALSE, /* verbose */ CPL_FALSE);
3124 if (err) cpl_error_set(cpl_func, err);
3125 CPLCHECK_NUL("Failed to fit noise-free model to data");
3126
3127 /* Store residual */
3128 sprintf(colname_residual, "PCA_RESID_BASE%d_POL%d", i, j);
3129 cpl_table_new_column_array(vis_table, colname_residual, CPL_TYPE_DOUBLE, nwave);
3130 cpl_matrix *resid = gravi_pca_get_data_residual(pca_decomps[i + GRAVI_NBASE * j]);
3131 for (int n = 0; n < naccept; n++) {
3132 cpl_matrix *resid_row = cpl_matrix_extract_row(resid, n);
3133 cpl_array *resid_row_arr = cpl_array_wrap_double(cpl_matrix_get_data(resid_row), nwave);
3134 cpl_table_set_array(vis_table, colname_residual, n, resid_row_arr);
3135 cpl_array_unwrap(resid_row_arr);
3136 cpl_matrix_delete(resid_row);
3137 }
3138 cpl_matrix_delete(resid);
3139 }
3140 }
3141
3142 /* "Refine" the mean i.e. adjust it to be the residuals */
3143 for (int i = 0; i < GRAVI_NBASE; i++) {
3144 for (int j = 0; j < npol; j++) {
3145 /* Copy residual into matrix of size (naccept, nwave) */
3146 sprintf(colname_residual, "PCA_RESID_BASE%d_POL%d", i, j);
3147 cpl_matrix *mat_tmp = cpl_matrix_new(naccept, nwave);
3148 for (int n = 0; n < naccept; n++) {
3149 const cpl_array *residual_arr = cpl_table_get_array(vis_table, colname_residual, n);
3150 for (int w = 0; w < nwave; w++) {
3151 cpl_matrix_set(mat_tmp, n, w, cpl_array_get(residual_arr, w, NULL));
3152 }
3153 }
3154
3155 /* Refine mean and re-centre data accordingly */
3156 gravi_pca_refine_mean(pca_decomps[i + GRAVI_NBASE * j], mat_tmp);
3157
3158 /* Store refined mean values */
3159 cpl_vector *mean = gravi_pca_get_component(pca_decomps[i + GRAVI_NBASE * j], 0);
3160 sprintf(colname_mean, "PCA_FLAT_MEAN_BASE%d_POL%d", i, j);
3161 cpl_table_new_column(pca_table, colname_mean, CPL_TYPE_DOUBLE);
3162 cpl_table_copy_data_double(pca_table, colname_mean, cpl_vector_get_data_const(mean));
3163
3164 cpl_matrix_delete(mat_tmp);
3165 }
3166 }
3167
3168 /* Now fit again using the mean-subtracted data */
3169 for (int i = 0; i < GRAVI_NBASE; i++) {
3170 for (int j = 0; j < npol; j++) {
3171 err = gravi_pca_fit_model(pca_decomps[i + GRAVI_NBASE * j], pca_model,
3172 /* fit_mean_subtracted */ CPL_TRUE, /* verbose */ CPL_FALSE);
3173 if (err) cpl_error_set(cpl_func, err);
3174 CPLCHECK_NUL("Failed to fit noise-free model to data");
3175
3176 /* Store residual */
3177 sprintf(colname_residual, "PCA_FLAT_RESID_BASE%d_POL%d", i, j);
3178 cpl_table_new_column_array(vis_table, colname_residual, CPL_TYPE_DOUBLE, nwave);
3179 cpl_matrix *resid = gravi_pca_get_data_residual(pca_decomps[i + GRAVI_NBASE * j]);
3180 for (int n = 0; n < naccept; n++) {
3181 cpl_matrix *resid_row = cpl_matrix_extract_row(resid, n);
3182 cpl_array *resid_row_arr = cpl_array_wrap_double(cpl_matrix_get_data(resid_row), nwave);
3183 cpl_table_set_array(vis_table, colname_residual, n, resid_row_arr);
3184 cpl_array_unwrap(resid_row_arr);
3185 cpl_matrix_delete(resid_row);
3186 }
3187 cpl_matrix_delete(resid);
3188 }
3189 }
3190
3191 /* Freeing the initial decompositions */
3192 for (int i = 0; i < GRAVI_NBASE * npol; i++)
3193 FREE(gravi_pca_result_delete, pca_decomps[i]);
3194 FREE(gravi_pca_model_delete, pca_model);
3195
3196 /* Create the output table */
3197 cpl_propertylist *header = cpl_propertylist_new();
3198 cpl_propertylist_append_string(header, "TELESCOP", telescope);
3199 cpl_propertylist_append_string(header, "ESO INS POLA MODE", pola_mode);
3200 cpl_propertylist_append_string(header, "ESO INS SPEC RES", spec_res);
3201 cpl_propertylist_append_int(header, "NWAVE", nwave);
3202 cpl_propertylist_append_int(header, "PCA NCOMP", num_components);
3203
3204 cpl_table *pca_result = cpl_table_new(nwave);
3205 cpl_table_duplicate_column(pca_result, "EFF_WAVE", wave_table, "EFF_WAVE");
3206 cpl_table_multiply_scalar(pca_result, "EFF_WAVE", 1.0e-6);
3207 cpl_table_set_column_unit(pca_result, "EFF_WAVE", "m");
3208
3209 for (int i = 0; i < GRAVI_NBASE; i++) {
3210 for (int j = 0; j < npol; j++) {
3211 sprintf(colname_mean, "PCA_MEAN_BASE%d_POL%d", i, j);
3212 cpl_table_duplicate_column(pca_result, colname_mean, pca_table, colname_mean);
3213
3214 // sprintf(colname_mean, "PCA_FLAT_MEAN_BASE%d_POL%d", i, j);
3215 // cpl_table_duplicate_column(pca_result, colname_mean, pca_table, colname_mean);
3216
3217 for (int c = 0; c < num_components; c++) {
3218 sprintf(colname_pca, "PCA_C%d_BASE%d_POL%d", c+1, i, j);
3219 cpl_table_duplicate_column(pca_result, colname_pca, pca_table, colname_pca);
3220 }
3221 }
3222 }
3223
3224 for (int c = 0; c < num_components; c++) {
3225 // sprintf(colname_pca, "PCA_C%d_MEDIAN", c+1);
3226 // cpl_table_duplicate_column(pca_result, colname_pca, pca_table, colname_pca);
3227
3228 sprintf(colname_pca, "PCA_C%d_FIT", c+1);
3229 cpl_table_duplicate_column(pca_result, colname_pca, pca_table, colname_pca);
3230 }
3231
3232 gravi_data *pca = gravi_data_new(0);
3233 gravi_data_add_table(pca, header, GRAVI_PCA_EXT, pca_result);
3234
3235 if (save_residuals) {
3236 /* Create optional extra table with filtered visphi and residuals */
3237 cpl_table *pca_resid = cpl_table_new(naccept);
3238 cpl_propertylist *header_resid = cpl_propertylist_duplicate(header);
3239
3240 for (int i = 0; i < GRAVI_NBASE; i++) {
3241 for (int j = 0; j < npol; j++) {
3242 sprintf(colname_vis, "PCA_VISPHI_BASE%d_POL%d", i, j);
3243 cpl_table_duplicate_column(pca_resid, colname_vis, vis_table, colname_vis);
3244
3245 sprintf(colname_residual, "PCA_RESID_BASE%d_POL%d", i, j);
3246 cpl_table_duplicate_column(pca_resid, colname_residual, vis_table, colname_residual);
3247
3248 sprintf(colname_residual, "PCA_FLAT_RESID_BASE%d_POL%d", i, j);
3249 cpl_table_duplicate_column(pca_resid, colname_residual, vis_table, colname_residual);
3250 }
3251 }
3252 gravi_data_add_table(pca, header_resid, GRAVI_PCA_RESID_EXT, pca_resid);
3253 }
3254
3255 cpl_array_unwrap(wave_arr);
3256 FREE(cpl_table_delete, pca_table);
3257 FREE(cpl_table_delete, vis_table);
3258
3259 /* Verbose */
3261 return pca;
3262}
3263
3264/*----------------------------------------------------------------------------*/
3265
3278/*----------------------------------------------------------------------------*/
3279cpl_error_code gravi_flatten_vis (gravi_data * vis_data,
3280 gravi_data * calib_data)
3281{
3282 /* Verbose */
3284 cpl_ensure_code(vis_data, CPL_ERROR_NULL_INPUT);
3285 cpl_ensure_code(calib_data, CPL_ERROR_NULL_INPUT);
3286
3287 int ncomp, nwave, nwave_obs, npol;
3288
3289 cpl_propertylist *calib_header = NULL, *vis_header = NULL, *obs_wave_plist = NULL;
3290 cpl_table *calib_table = NULL;
3291 const char *telescope = NULL, *pola_mode = NULL, *spec_res = NULL;
3292 char colname_calib[100], colname_vis[100], colname_mask[100], colname_vis_fit[100], colname_vis_resid[100];
3293 double vis_mjd;
3294
3295 vis_header = gravi_data_get_header(vis_data);
3296 telescope = cpl_propertylist_get_string(vis_header, "TELESCOP");
3297 pola_mode = gravi_pfits_get_pola_mode(vis_header, GRAVI_SC);
3298 npol = gravi_pfits_get_pola_num(vis_header, GRAVI_SC);
3299 spec_res = gravi_pfits_get_spec_res(vis_header);
3300 vis_mjd = cpl_propertylist_get_double(vis_header, "MJD-OBS");
3301
3302 /* Check calibrator for mode and epoch compatibility */
3303 calib_table = gravi_data_get_table(calib_data, GRAVI_PCA_EXT);
3304 if (!calib_table)
3305 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "No VISPHI flattening calibration");
3306
3307 calib_header = gravi_data_get_plist(calib_data, GRAVI_PCA_EXT);
3308 if (vis_mjd < cpl_propertylist_get_double(calib_header, "PCA EPOCH BEGIN"))
3309 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "VISPHI flattening calibration is too new for data");
3310 else if (vis_mjd > cpl_propertylist_get_double(calib_header, "PCA EPOCH END"))
3311 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "VISPHI flattening calibration is too old for data");
3312
3313 if (strcmp(gravi_pfits_get_spec_res(calib_header), spec_res))
3314 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "Calibration resolution does not match data");
3315 if (strcmp(gravi_pfits_get_pola_mode(calib_header, GRAVI_SC), pola_mode))
3316 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "Calibration polarisation mode does not match data");
3317 if (strcmp(cpl_propertylist_get_string(calib_header, "TELESCOP"), telescope))
3318 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND, "Calibration telescope does not match data");
3319
3320 ncomp = cpl_propertylist_get_int(calib_header, "PCA NCOMP");
3321 nwave = cpl_propertylist_get_int(calib_header, "NWAVE");
3322
3323 /* Copy calibration data into model */
3324 cpl_matrix *tmp_matrix = cpl_matrix_new(ncomp, nwave);
3325 for (int c = 0; c < ncomp; c++) {
3326 sprintf(colname_calib, "PCA_C%d_FIT", c + 1);
3327 for (int w = 0; w < nwave; w++) {
3328 double val = cpl_table_get(calib_table, colname_calib, w, NULL);
3329 cpl_matrix_set(tmp_matrix, c, w, val);
3330 }
3331 }
3332 gravi_pca_model *calib_model = gravi_pca_load_model(tmp_matrix);
3333 cpl_matrix_delete(tmp_matrix);
3334
3335 /* Check the observations for compatibility */
3336 obs_wave_plist = gravi_data_get_oi_wave_plist(vis_data, GRAVI_SC, 0, npol);
3337 nwave_obs = cpl_propertylist_get_int(obs_wave_plist, "NWAVE");
3338 if (nwave != nwave_obs) {
3339 return cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
3340 "Input file wavelength axis does not match the calibrator");
3341 }
3342
3343 /* Prepare table for collated visphi */
3344 cpl_table *vis_table = cpl_table_new(nwave);
3345 for(int i = 0; i < GRAVI_NBASE; i++) {
3346 for(int j = 0; j < npol; j++) {
3347 sprintf(colname_vis, "PCA_VISPHI_BASE%d_POL%d", i, j);
3348 cpl_table_new_column(vis_table, colname_vis, CPL_TYPE_DOUBLE);
3349 sprintf(colname_mask, "PCA_MASK_BASE%d_POL%d", i, j);
3350 cpl_table_new_column(vis_table, colname_mask, CPL_TYPE_INT);
3351 }
3352 }
3353
3354 /* Collate visphi */
3355 for (int j = 0; j < npol; j++) {
3356 cpl_table *vis_tmp = gravi_data_get_oi_vis(vis_data, GRAVI_SC, j, npol);
3357 for (int i = 0; i < GRAVI_NBASE; i++) {
3358 /* Copy the visphi data */
3359 const cpl_array *arr_tmp = cpl_table_get_array(vis_tmp, "VISPHI", i);
3360 sprintf(colname_vis, "PCA_VISPHI_BASE%d_POL%d", i, j);
3361 cpl_table_copy_data_double(vis_table, colname_vis, cpl_array_get_data_double_const(arr_tmp));
3362
3363 /* Copy the flag data */
3364 arr_tmp = cpl_table_get_array(vis_tmp, "FLAG", i);
3365 sprintf(colname_mask, "PCA_MASK_BASE%d_POL%d", i, j);
3366 cpl_table_copy_data_int(vis_table, colname_mask, cpl_array_get_data_int_const(arr_tmp));
3367 cpl_table_cast_column(vis_table, colname_mask, colname_mask, CPL_TYPE_DOUBLE);
3368 }
3369 }
3370
3371 /* Construct the PCA object */
3375 gravi_pca_result **results = cpl_malloc(GRAVI_NBASE * npol * sizeof(gravi_pca_result *));
3376 for(int i = 0; i < GRAVI_NBASE; i++) {
3377 for(int j = 0; j < npol; j++) {
3378 sprintf(colname_vis, "PCA_VISPHI_BASE%d_POL%d", i, j);
3379 sprintf(colname_mask, "PCA_MASK_BASE%d_POL%d", i, j);
3380 sprintf(colname_calib, "PCA_MEAN_BASE%d_POL%d", i, j);
3381
3382 /* Copy data into matrix of size (n, nwave) */
3383 cpl_matrix *vis_tmp = cpl_matrix_wrap(1, nwave, cpl_table_get_data_double(vis_table, colname_vis));
3384 cpl_matrix *vis_mean = cpl_matrix_wrap(1, nwave, cpl_table_get_data_double(calib_table, colname_calib));
3385 cpl_matrix_subtract(vis_tmp, vis_mean);
3386 cpl_matrix *mask_tmp = cpl_matrix_wrap(1, nwave, cpl_table_get_data_double(vis_table, colname_mask));
3387 results[i + GRAVI_NBASE * j] = gravi_pca_create_result(vis_tmp, mask_tmp);
3388 cpl_matrix_unwrap(vis_tmp);
3389 cpl_matrix_unwrap(mask_tmp);
3390 cpl_matrix_unwrap(vis_mean);
3391 }
3392 }
3393
3394 /* Fit the calibrated PCA components to the data */
3395 for(int i = 0; i < GRAVI_NBASE; i++) {
3396 for(int j = 0; j < npol; j++) {
3397 cpl_error_code err = gravi_pca_fit_model(results[i + GRAVI_NBASE * j], calib_model,
3398 /* fit_mean_subtracted */ CPL_FALSE, /* verbose */ CPL_FALSE);
3399 if (err) cpl_error_set(cpl_func, err);
3400 CPLCHECK_MSG("Failed to fit PCA model to data");
3401
3402 /* Store model */
3403 sprintf(colname_vis_fit, "PCA_VISPHI_MODEL_BASE%d_POL%d", i, j);
3404 cpl_table_new_column(vis_table, colname_vis_fit, CPL_TYPE_DOUBLE);
3405 cpl_matrix *model_eval = gravi_pca_get_data_fit(results[i + GRAVI_NBASE * j]);
3406
3407 sprintf(colname_calib, "PCA_MEAN_BASE%d_POL%d", i, j);
3408 cpl_matrix *vis_mean = cpl_matrix_wrap(1, nwave, cpl_table_get_data_double(calib_table, colname_calib));
3409 cpl_matrix_add(model_eval, vis_mean);
3410 cpl_table_copy_data_double(vis_table, colname_vis_fit, cpl_matrix_get_data_const(model_eval));
3411
3412 cpl_matrix_delete(model_eval);
3413 cpl_matrix_unwrap(vis_mean);
3414
3415 /* Store residual */
3416 sprintf(colname_vis_resid, "PCA_VISPHI_RESID_BASE%d_POL%d", i, j);
3417 cpl_table_new_column(vis_table, colname_vis_resid, CPL_TYPE_DOUBLE);
3418 cpl_matrix *resid = gravi_pca_get_data_residual(results[i + GRAVI_NBASE * j]);
3419 cpl_table_copy_data_double(vis_table, colname_vis_resid, cpl_matrix_get_data_const(resid));
3420 cpl_matrix_delete(resid);
3421 }
3422 }
3423
3424 /* Generate the output table */
3425 for (int j = 0; j < npol; j++) {
3426 cpl_table *vis_tmp = gravi_data_get_oi_vis(vis_data, GRAVI_SC, j, npol);
3427 cpl_table_new_column_array(vis_tmp, "VISPHI_MODEL", CPL_TYPE_DOUBLE, nwave);
3428 cpl_table_new_column_array(vis_tmp, "VISPHI_FLAT", CPL_TYPE_DOUBLE, nwave);
3429 for (int i = 0; i < GRAVI_NBASE; i++) {
3430 /* Copy the PCA model */
3431 sprintf(colname_vis_fit, "PCA_VISPHI_MODEL_BASE%d_POL%d", i, j);
3432 cpl_array *arr_tmp = cpl_array_wrap_double(
3433 cpl_table_get_data_double(vis_table, colname_vis_fit), nwave);
3434 cpl_table_set_array(vis_tmp, "VISPHI_MODEL", i, arr_tmp);
3435 cpl_array_unwrap(arr_tmp);
3436
3437 /* Copy the flattened visphi data */
3438 sprintf(colname_vis_resid, "PCA_VISPHI_RESID_BASE%d_POL%d", i, j);
3439 arr_tmp = cpl_array_wrap_double(
3440 cpl_table_get_data_double(vis_table, colname_vis_resid), nwave);
3441 cpl_table_set_array(vis_tmp, "VISPHI_FLAT", i, arr_tmp);
3442 cpl_array_unwrap(arr_tmp);
3443 }
3444 }
3445 CPLCHECK_MSG("Failed to add calibrated VISPHI to table");
3446
3447 FREE(gravi_pca_model_delete, calib_model);
3448 FREE(cpl_table_delete, vis_table);
3450
3451 /* Verbose */
3453 return CPL_ERROR_NONE;
3454}
3455
#define IMAGE_IDX(_nx, _x, _y)
#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_header(data)
Definition: gravi_data.h:75
#define gravi_data_get_oi_wave(data, type, pol, npol)
Definition: gravi_data.h:45
#define gravi_data_get_oi_vis(data, type, pol, npol)
Definition: gravi_data.h:46
#define gravi_data_get_oi_wave_plist(data, type, pol, npol)
Definition: gravi_data.h:69
#define gravi_data_get_img(data, ext)
Definition: gravi_data.h:79
const cpl_size ntel
cpl_image_multiply(flat_profiled, profile_crop)
cpl_propertylist * header
Definition: gravi_old.c:2004
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)
cpl_error_code gravi_pca_fit_components_polynomial(gravi_pca_model *self, const cpl_array *wave, int degree, int ncomp)
Fit polynomial model to each of a set of median-averaged PCA components.
Definition: gravi_pca.c:587
cpl_error_code gravi_pca_refine_mean(gravi_pca_result *self, const cpl_matrix *residual)
Override the mean component calculated from the PCA decomposition. In conjunction with the fit_mean_s...
Definition: gravi_pca.c:652
cpl_matrix * gravi_pca_get_data_residual(const gravi_pca_result *self)
Get residual (data - model).
Definition: gravi_pca.c:956
cpl_error_code gravi_pca_set_component_signs(gravi_pca_result *self, const cpl_array *wave, const int ncomp)
Impose sign convention on PCA components.
Definition: gravi_pca.c:393
gravi_pca_model * gravi_pca_create_model(const gravi_pca_result **results, int nres, int ncomp)
Compute median values of PCA components over a set of decomposition results.
Definition: gravi_pca.c:423
void gravi_pca_result_delete(gravi_pca_result *self)
Deallocate a gravi_pca_result object.
Definition: gravi_pca.c:99
gravi_pca_result * gravi_pca_create_result(const cpl_matrix *data, const cpl_matrix *mask)
Construct a new gravi_pca_result object from a matrix of visphi data.
Definition: gravi_pca.c:219
cpl_vector * gravi_pca_get_component(const gravi_pca_result *self, int component)
Get components from PCA decomposition.
Definition: gravi_pca.c:907
cpl_matrix * gravi_pca_get_data_fit(const gravi_pca_result *self)
Get noise-free model.
Definition: gravi_pca.c:935
cpl_vector * gravi_pca_get_component_median(const gravi_pca_model *self, int component)
Get median-averaged component from a set of PCA decompositions.
Definition: gravi_pca.c:694
cpl_error_code gravi_pca_decomp_matrix_svd(gravi_pca_result *self)
Perform PCA decomposition by calculating singular value decomposition.
Definition: gravi_pca.c:283
void gravi_pca_model_delete(gravi_pca_model *self)
Deallocate a gravi_pca_model object.
Definition: gravi_pca.c:118
cpl_error_code gravi_pca_fit_components_bspline(gravi_pca_model *self, const cpl_array *wave, int degree, int ncomp)
Fit B-spline model to each of a set of median-averaged PCA components.
Definition: gravi_pca.c:511
cpl_vector * gravi_pca_get_component_fit(const gravi_pca_model *self, int component)
Get fit to components from PCA decomposition.
Definition: gravi_pca.c:716
gravi_pca_model * gravi_pca_load_model(const cpl_matrix *components)
Create PCA model from existing set of components.
Definition: gravi_pca.c:483
cpl_error_code gravi_pca_fit_model(gravi_pca_result *self, const gravi_pca_model *model, cpl_boolean fit_mean_subtracted, cpl_boolean verbose)
Fit model formed from linear combination of PCA components to data.
Definition: gravi_pca.c:792
#define QC_MEANDARK_SC
Definition: gravi_pfits.h:118
#define QC_MEANDARK_MET
Definition: gravi_pfits.h:128
#define QC_BADPIX_FT
Definition: gravi_pfits.h:144
#define GRAVI_PROFILE_PARAMS_EXT
Definition: gravi_pfits.h:80
#define GRAVI_IMAGING_DETECTOR_FT_EXT
Definition: gravi_pfits.h:82
#define QC_ACQ_ZERO_NB
Definition: gravi_pfits.h:135
#define GRAVI_OPDC_EXT
Definition: gravi_pfits.h:62
#define GRAVI_SC
Definition: gravi_pfits.h:165
#define QC_DARKRMS_MET
Definition: gravi_pfits.h:129
#define QC_DARKRMS_SC
Definition: gravi_pfits.h:119
#define QC_SKYRMS_SC
Definition: gravi_pfits.h:121
#define GRAVI_PROFILE_DATA_EXT
Definition: gravi_pfits.h:79
#define QC_MEANGAIN_FT
Definition: gravi_pfits.h:143
#define PROFILE_FULLSTARTX
Definition: gravi_pfits.h:146
#define GRAVI_IMAGING_DATA_FT_EXT
Definition: gravi_pfits.h:43
#define GRAVI_BIAS_MASK_SC_EXT
Definition: gravi_pfits.h:58
#define GRAVI_IMAGING_ERR_FT_EXT
Definition: gravi_pfits.h:46
#define PROFILE_NX
Definition: gravi_pfits.h:147
#define QC_BADPIX_RMS_SC
Definition: gravi_pfits.h:141
#define QC_MEANSKY_SC
Definition: gravi_pfits.h:120
#define QC_DARKRANGE_MET
Definition: gravi_pfits.h:130
#define GRAVI_PCA_EXT
Definition: gravi_pfits.h:76
#define GRAVI_IMAGING_DATA_SC_EXT
Definition: gravi_pfits.h:44
#define QC_BADPIX_SC
Definition: gravi_pfits.h:139
#define QC_MEANDARK_FT
Definition: gravi_pfits.h:123
#define GRAVI_IMAGING_DATA_ACQ_EXT
Definition: gravi_pfits.h:41
#define PROFILE_STARTX
Definition: gravi_pfits.h:145
#define QC_MEANSKY_FT
Definition: gravi_pfits.h:125
#define QC_DARKRMS_FT
Definition: gravi_pfits.h:124
#define GRAVI_IMAGING_ERR_SC_EXT
Definition: gravi_pfits.h:45
#define QC_BADPIX_FLAT_SC
Definition: gravi_pfits.h:142
#define GRAVI_PCA_RESID_EXT
Definition: gravi_pfits.h:77
#define QC_SKYRMS_FT
Definition: gravi_pfits.h:126
#define GRAVI_IMAGING_DETECTOR_SC_EXT
Definition: gravi_pfits.h:81
#define QC_BADPIX_DARK_SC
Definition: gravi_pfits.h:140
#define GRAVI_METROLOGY_EXT
Definition: gravi_pfits.h:60
#define QC_MEANGAIN_SC
Definition: gravi_pfits.h:138
#define GRAVI_IMAGING_MASK_SC_EXT
Definition: gravi_pfits.h:47
#define GRAVI_METROLOGY_ERR_EXT
Definition: gravi_pfits.h:61
#define gravi_msg_function_exit(flag)
Definition: gravi_utils.h:85
#define gravi_data_check_shutter_closed(data)
Definition: gravi_utils.h:90
#define FREE(function, variable)
Definition: gravi_utils.h:69
#define CPLCHECK_NUL(msg)
Definition: gravi_utils.h:48
#define gravi_msg_function_start(flag)
Definition: gravi_utils.h:84
#define CPLCHECK_MSG(msg)
Definition: gravi_utils.h:45
#define gravi_data_check_shutter_open(data)
Definition: gravi_utils.h:91
#define FREELOOP(function, variable, n)
Definition: gravi_utils.h:72
#define GRAVI_NBASE
Definition: gravi_utils.h:105
static double met_Jun_2017[64]
Definition: gravi_calib.c:81
gravi_data * gravi_compute_piezotf(gravi_data *data, const cpl_parameterlist *params)
Create piezo transfer function for Kalman Calibration & monitoring.
Definition: gravi_calib.c:2376
static double met_Mar_2017[64]
Definition: gravi_calib.c:80
cpl_error_code gravi_flatten_vis(gravi_data *vis_data, gravi_data *calib_data)
Use PCA model to flatten observed visphi. The flattened data are added to the existing VIS table.
Definition: gravi_calib.c:3279
cpl_propertylist * gravi_compute_gain(gravi_data **flats_data, int nrawgain, gravi_data *dark_map)
Compute mean detector gain.
Definition: gravi_calib.c:1551
static double met_Aug_2017[64]
Definition: gravi_calib.c:83
gravi_data * gravi_compute_profile(gravi_data **flats_data, gravi_data *dark_map, gravi_data *bad_map, int nflat, const cpl_parameterlist *params)
Computes the spatial profile of each spectrum for optimal extraction purpose.
Definition: gravi_calib.c:994
cpl_image * gravi_create_profile_image(cpl_image *mean_img, cpl_vector *values_x0, cpl_vector *values_y0, cpl_vector *values_sigma, cpl_size iy_min, cpl_size iy_max, const char *resolution)
Create the profile image from image and/or fitted params.
Definition: gravi_calib.c:814
double gravi_calc_sigmaclipped_stddev(const cpl_array *arr, double nstd)
Perform sigma-clipping on input array and return std of clipped array.
Definition: gravi_calib.c:2891
gravi_data * gravi_compute_pca(gravi_data **data, int naccept, const cpl_parameterlist *params)
Fit model for visphi flattening using PCA.
Definition: gravi_calib.c:2945
gravi_data * gravi_average_dark(gravi_data **data, cpl_size ndata)
Average several DARK calibration map.
Definition: gravi_calib.c:457
gravi_data * gravi_compute_dark(gravi_data *raw_data)
Compute the DARK calibration map.
Definition: gravi_calib.c:125
static double met_Sep_2016[64]
Definition: gravi_calib.c:79
static double met_Jul_2017[64]
Definition: gravi_calib.c:82
cpl_array * gravi_filter_array_median(const cpl_array *arr, int size)
Apply median filter on 1-D array.
Definition: gravi_calib.c:2862
cpl_error_code gravi_remove_cosmicrays_sc(cpl_imagelist *imglist_sc)
Remove cosmic rays via filtering through images.
Definition: gravi_calib.c:2661
gravi_data * gravi_compute_badpix(gravi_data *dark_map, gravi_data **flats_data, int nflat, const cpl_parameterlist *params)
Identify the bad pixels in the DARK map and create the BAD map.
Definition: gravi_calib.c:1806
gravi_data * gravi_compute_biasmask(gravi_data *dark_map, gravi_data **flats_data, int nflat, const cpl_parameterlist *params)
Create BIASMASK for SC from raw FLATs and raw DARK.
Definition: gravi_calib.c:2271
cpl_error_code gravi_fit_profile(cpl_vector *values_x0, cpl_vector *values_y0, cpl_vector *values_sigma, cpl_image *mean_img, int ref_x, int ref_y, int size_profile, const char *resolution)
Fit the profile parameters in an image.
Definition: gravi_calib.c:536
cpl_image * gravi_image_from_column(cpl_table *table_data, const char *data_x, cpl_size row)
Create an image from a column array in table.
Definition: gravi_cpl.c:2052
cpl_imagelist * gravi_imagelist_from_column(cpl_table *table_data, const char *data_x)
Create an imagelist from a column array in table.
Definition: gravi_cpl.c:2084
double gravi_image_get_quantile(const cpl_image *img, double thr)
Compute the quantile of an image.
Definition: gravi_cpl.c:2477
cpl_array * gravi_array_wrap_image(cpl_image *img)
Wrap the data of na image into an array.
Definition: gravi_cpl.c:2141
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
gravi_data * gravi_data_duplicate(const gravi_data *self)
Create a copy of the gravi data.
Definition: gravi_data.c:250
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_error_code gravi_data_add_img(gravi_data *self, cpl_propertylist *plist, const char *extname, cpl_image *image)
Add an IMAGE (single image) extension in gravi_data.
Definition: gravi_data.c:2394
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(const cpl_parameterlist *parlist, const char *name)
Definition: gravi_dfs.c:1524
const char * gravi_param_get_string_default(const cpl_parameterlist *parlist, const char *name, const char *def)
Definition: gravi_dfs.c:1499
int gravi_pfits_get_pola_num(const cpl_propertylist *plist, int type_data)
Definition: gravi_pfits.c:263
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
const char * gravi_pfits_get_resolution(const cpl_propertylist *plist)
Definition: gravi_pfits.c:155
const char * gravi_pfits_get_dpr_type(const cpl_propertylist *plist)
Definition: gravi_pfits.c:148
int gravi_pfits_get_window_start(const cpl_propertylist *plist)
Definition: gravi_pfits.c:92
cpl_error_code gravi_remove_badpixel_sc(cpl_imagelist *imglist_sc, cpl_image *bad_img)
Remove the badpixel of the SC.
cpl_error_code gravi_msg_warning(const char *component, const char *msg)
Definition: gravi_utils.c:127
int gravi_check_shutter(cpl_propertylist *header, int t0, int t1, int t2, int t3)
Definition: gravi_utils.c:943
char GRAVI_DATA[50][7]
Definition: gravi_utils.c:71
int * gravi_image_extract_dimension(cpl_image *img_profile)
Compute startx and nx of the illuminated part of the image.
Definition: gravi_utils.c:727
int gravi_get_shutter_id(cpl_propertylist *header)
Definition: gravi_utils.c:956
int gravi_data_check_shutter_beam(gravi_data **datas, int nb_datas)
Definition: gravi_utils.c:986
int gravi_region_get_tel(cpl_table *imaging_detector, int region, int beam)
Return the telescope id (0,1,2,3) in a beam of a region.
Definition: gravi_utils.c:512
Type to hold average (median) components obtained from a set of PCA decompositions and/or best-fit mo...
Definition: gravi_pca.c:87
Type to hold results of a PCA decomposition.
Definition: gravi_pca.c:64