ERIS Pipeline Reference Manual 1.8.10
eris_ifu_efficiency_response.c
1/* $Id$
2 *
3 * This file is part of the ERIS 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., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
19 */
20
21#ifdef HAVE_CONFIG_H
22#include <config.h>
23#endif
24#include <string.h>
25#include <gsl/gsl_version.h>
26#include <gsl/gsl_multifit.h>
27#include <eris_ifu_efficiency_response.h>
28#include <eris_ifu_star_index.h>
29#include <eris_utils.h>
30#include <eris_ifu_dfs.h>
31#include <eris_ifu_utils.h>
32#include <eris_utils.h>
33#include "eris_pfits.h"
34// TODO: next include is to get a cube saving function. That should be elsewhere
35#include <eris_ifu_resample.h>
36const hdrl_spectrum1D_wave_scale
37global_scale = hdrl_spectrum1D_wave_scale_linear;
38
39
40
41#define KEY_VALUE_HPRO_DID "PRO-1.15"
42#define TELESCOPE_SURFACE 52.8101279
43#define FILE_NAME_SZ 256
44#define KEY_NAME_FILT_NAME "ESO INS3 SPGW ID"
45#define KEY_NAME_FILT_NAME2 "ESO INS FILT1 NAME"
46
47#define RESPONSE_FILENAME "out_response.fits"
48#define DISPERSION_K_DITH 0.000245
49#define DISPERSION_H_DITH 0.000195
50#define DISPERSION_J_DITH 0.000145
51#define DISPERSION_HK_DITH 0.0005
52/*
53 * Band lam_c lam-range Resol Disp
54 * --------------------------------------
55 * J_low 1.25 1.09-1.42 5000 0.00025
56 * H_low 1.66 1.45-1.87 5200 0.000319230769231
57 * K_low 2.21 1.93-2.48 5600 0.000394642857143
58 * ------------------------------------------
59 * J_short 1.18 1.10-1.27 10000 0.000118
60 * J_middle 1.26 1.18-1.35 10000 0.000126
61 * J_long 1.34 1.26-1.43 10000 0.000134
62 * --------------------------------------------
63 * H_short 1.56 1.46-1.67 10000 0.000156
64 * H_middle 1.67 1.56-1.77 10000 0.000167
65 * H_long 1.76 1.66-1.87 10000 0.000176
66 * --------------------------------------------
67 * K_short 2.07 1.93-2.22 11200 0.000207
68 * K_middle 2.20 2.06-2.34 11200 0.000220
69 * K_long 2.33 2.19-2.47 11200 0.000233
70 * -------------------------------------------
71 */
72#define ERIS_GAIN 2.0
73#define SINFONI_GAIN 2.42
74#define DISPERSION_J_low 0.00025
75#define DISPERSION_J_short 0.000118
76#define DISPERSION_J_middle 0.000126
77#define DISPERSION_J_long 0.000134
78
79#define DISPERSION_H_low 0.000319230769231
80#define DISPERSION_H_short 0.000156
81#define DISPERSION_H_middle 0.000167
82#define DISPERSION_H_long 0.000176
83
84#define DISPERSION_K_low 0.000394642857143
85#define DISPERSION_K_short 0.000207
86#define DISPERSION_K_middle 0.000220
87#define DISPERSION_K_long 0.000233
88
89
90typedef struct
91{
92 double *x, *y, *err;
93 int n;
94} eris_poly_data;
95
96static double
97eris_ifu_get_airmass_mean(cpl_propertylist* plist) {
98 double airmass_start = cpl_propertylist_get_double(plist, "ESO TEL AIRM START");
99 double airmass_end = cpl_propertylist_get_double(plist, "ESO TEL AIRM END");
100 double airmass = 0.5 * (airmass_start + airmass_end);
101 return airmass;
102}
103
104static cpl_error_code
105eris_get_wcs_cube(const cpl_propertylist* plist, double *clambda,
106 double *dis, double *cpix, double *cx, double *cy)
107{
108
109 *clambda = eris_pfits_get_crval3(plist);
110 *dis = eris_pfits_get_cd33(plist);
111 *cpix = eris_pfits_get_crpix3(plist);
112 *cx = eris_pfits_get_crpix1(plist);
113 *cy = eris_pfits_get_crpix2(plist);
114
115 eris_check_error_code("eris_get_wcs_cube");
116 return cpl_error_get_code();
117
118}
119/*---------------------------------------------------------------------------*/
126/*---------------------------------------------------------------------------*/
127
128static double
129sinfo_get_dispersion(const char* band)
130{
131 double disp = 0.;
132 if (strcmp(band,"H+K") == 0) {
133 disp = DISPERSION_HK_DITH;
134 } else if (strcmp(band,"K") == 0) {
135 disp = DISPERSION_K_DITH;
136 } else if (strcmp(band,"J") == 0) {
137 disp = DISPERSION_J_DITH;
138 } else if (strcmp(band,"H") == 0) {
139 disp = DISPERSION_H_DITH;
140 }
141 return disp;
142
143}
151static double eris_get_dispersion(const char* band)
152{
153 cpl_ensure(band != NULL, CPL_ERROR_NULL_INPUT, 0);
154 double disp = 0.;
155 if (strcmp(band, "J_low") == 0) {
156 disp = DISPERSION_J_low;
157 } else if (strcmp(band, "J_short") == 0) {
158 disp = DISPERSION_J_short;
159 } else if (strcmp(band, "J_middle") == 0) {
160 disp = DISPERSION_J_middle;
161 } else if (strcmp(band, "J_long") == 0) {
162 disp = DISPERSION_J_long;
163 } else if (strcmp(band, "H_low") == 0) {
164 disp = DISPERSION_J_low;
165 } else if (strcmp(band, "H_short") == 0) {
166 disp = DISPERSION_H_short;
167 } else if (strcmp(band, "H_middle") == 0) {
168 disp = DISPERSION_H_middle;
169 } else if (strcmp(band, "H_long") == 0) {
170 disp = DISPERSION_H_long;
171 } else if (strcmp(band, "K_low") == 0) {
172 disp = DISPERSION_K_low;
173 } else if (strcmp(band, "K_short") == 0) {
174 disp = DISPERSION_K_short;
175 } else if (strcmp(band, "K_middle") == 0) {
176 disp = DISPERSION_K_middle;
177 } else if (strcmp(band, "K_long") == 0) {
178 disp = DISPERSION_K_long;
179 }
180 eris_check_error_code("eris_get_dispersion");
181 return disp;
182
183}
184/*---------------------------------------------------------------------------*/
192/*---------------------------------------------------------------------------*/
193
194static cpl_error_code
195eris_get_band(cpl_frame * ref_frame, char * band)
196{
197 cpl_msg_debug(cpl_func,"ref_frame: %p filename: %s",
198 (const void*)ref_frame, cpl_frame_get_filename(ref_frame));
199
200 cpl_ensure(ref_frame != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
201 cpl_ensure(band != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
202
203 char* ref_file = NULL;
204 cpl_propertylist* plist = NULL;
205
206 ref_file = cpl_strdup(cpl_frame_get_filename(ref_frame)) ;
207
208 if ((cpl_error_code)((plist = cpl_propertylist_load(ref_file, 0)) == NULL)) {
209 cpl_msg_error(cpl_func, "getting header from reference frame %s",
210 ref_file);
211 cpl_propertylist_delete(plist) ;
212 return cpl_error_get_code() ;
213 }
214
215 if (cpl_propertylist_has(plist, KEY_NAME_FILT_NAME)) {
216 strcpy(band, cpl_propertylist_get_string(plist, KEY_NAME_FILT_NAME));
217 /* cpl_msg_info(cpl_func, "%s value is %s", KEY_NAME_FILT_NAME, band); */
218
219 } else {
220 cpl_msg_warning(cpl_func, "keyword %s does not exist",
221 KEY_NAME_FILT_NAME);
222 return cpl_error_get_code() ;
223 }
224
225 cpl_free(ref_file);
226 cpl_propertylist_delete(plist);
227 eris_check_error_code("eris_get_band");
228 return cpl_error_get_code() ;
229}
230
231/*---------------------------------------------------------------------------*/
239/*---------------------------------------------------------------------------*/
240
241static int
242sinfo_get_band(cpl_frame * ref_frame, char * band)
243{
244 cpl_ensure(ref_frame != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
245 cpl_ensure(band != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
246
247 char* ref_file = NULL;
248 cpl_propertylist* plist = NULL;
249
250 ref_file = cpl_strdup(cpl_frame_get_filename(ref_frame)) ;
251 if ((cpl_error_code)((plist = cpl_propertylist_load(ref_file, 0)) == NULL)) {
252 cpl_msg_error(cpl_func, "getting header from reference frame %s",
253 ref_file);
254 cpl_propertylist_delete(plist) ;
255 return cpl_error_get_code();
256 }
257
258 if (cpl_propertylist_has(plist, KEY_NAME_FILT_NAME2)) {
259 strcpy(band, cpl_propertylist_get_string(plist, KEY_NAME_FILT_NAME2));
260 /* cpl_msg_info(cpl_func, "%s value is %s", KEY_NAME_FILT_NAME, band); */
261
262 } else {
263 cpl_msg_warning(cpl_func, "keyword %s does not exist",
264 KEY_NAME_FILT_NAME2);
265 return cpl_error_get_code();
266 }
267
268 cpl_free(ref_file);
269 cpl_propertylist_delete(plist);
270 return cpl_error_get_code();
271}
272cpl_boolean eris_can_extract(cpl_frameset* sof) {
273
274 cpl_boolean has_cube_coadd = CPL_TRUE;
275
276 if (
277 NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_TWK_CUBE_COADD) &&
278 NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_DAR_CUBE_COADD) &&
279 NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_OBJ_CUBE_COADD) &&
280 NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_OBJ_DAR_CUBE_COADD) &&
281 NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_STD_CUBE_COADD) &&
282 NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_STD_FLUX_CUBE_COADD) &&
283 NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_STD_FLUX_CUBE_COADD_NOFLAT) &&
284 NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_PSF_CUBE_COADD)
285 ){
286 cpl_msg_warning(cpl_func,"Provide %s or %s or %s or %s or %s or %s or %s or %s"
287 "to extract..",
288 ERIS_IFU_PRO_JITTER_TWK_CUBE_COADD,
289 ERIS_IFU_PRO_JITTER_DAR_CUBE_COADD,
290 ERIS_IFU_PRO_JITTER_OBJ_CUBE_COADD,
291 ERIS_IFU_PRO_JITTER_OBJ_DAR_CUBE_COADD,
292 ERIS_IFU_PRO_JITTER_STD_CUBE_COADD,
293 ERIS_IFU_PRO_JITTER_STD_FLUX_CUBE_COADD,
294 ERIS_IFU_PRO_JITTER_PSF_CUBE_COADD,
295 ERIS_IFU_PRO_JITTER_STD_FLUX_CUBE_COADD_NOFLAT
296 );
297 has_cube_coadd = CPL_FALSE;
298 }
299 return has_cube_coadd;
300}
301int eris_can_flux_calibrate(cpl_frameset* sof) {
302
303 int has_extcoeff_table=1;
304 int has_response = 1;
305
306
307 if (NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_RESPONSE) ){
308 cpl_msg_warning(cpl_func,"Provide %s to flux calibrate.",
309 ERIS_IFU_PRO_JITTER_RESPONSE);
310 has_response = 0;
311 }
312
313 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_EXTCOEFF_TABLE) ){
314 cpl_msg_warning(cpl_func,"Provide %s to flux calibrate.",
315 ERIS_IFU_CALIB_EXTCOEFF_TABLE);
316 has_extcoeff_table = 0;
317 }
318
319
320 eris_check_error_code("eris_can_flux_calibrate");
321 return has_response *
322 has_extcoeff_table;
323
324}
325
326
327int eris_can_compute_efficiency(cpl_frameset* sof) {
328
329 int has_std_star_spectra=1;
330 int has_extcoeff_table=1;
331 int has_response = 1;
332
333 if (NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_SPECTRUM_NOFLAT)){
334 cpl_msg_warning(cpl_func,"Provide %s to compute efficiency.",
335 ERIS_IFU_PRO_JITTER_SPECTRUM_NOFLAT);
336 has_std_star_spectra=0;
337 }
338
339 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_EFFICIENCY_WINDOWS) ){
340 cpl_msg_warning(cpl_func,"Provide %s to compute efficiency.",
341 ERIS_IFU_CALIB_EFFICIENCY_WINDOWS);
342 has_response = 0;
343 }
344
345 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_EXTCOEFF_TABLE) ){
346 cpl_msg_warning(cpl_func,"Provide %s to compute efficiency.",
347 ERIS_IFU_CALIB_EXTCOEFF_TABLE);
348 has_extcoeff_table = 0;
349 }
350
351
352 eris_check_error_code("eris_can_compute_efficiency");
353 return has_std_star_spectra *
354 has_response *
355 has_extcoeff_table;
356
357}
358
359
360int eris_can_compute_response(cpl_frameset* sof) {
361
362 int has_std_star_spectra=1;
363 int has_flux_std_catalog=1;
364 int has_tel_mod_catalog=1;
365 int has_resp_fit_points_cat=1;
366 int has_extcoeff_table=1;
367 int has_high_abs_regions=1;
368 int has_quality_areas=1;
369 int has_fit_areas=1;
370 int has_resp_windows=1;
371
372 if (NULL == cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_SPECTRUM) ){
373 cpl_msg_warning(cpl_func,"Provide %s to compute response.",
374 ERIS_IFU_PRO_JITTER_SPECTRUM);
375 has_std_star_spectra=0;
376 }
377
378 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_FLUX_STD_CATALOG) ){
379 cpl_msg_warning(cpl_func,"Provide %s to compute response.",
380 ERIS_IFU_CALIB_FLUX_STD_CATALOG);
381 has_flux_std_catalog=0;
382 }
383 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_TELL_MOD_CATALOG) ){
384 cpl_msg_warning(cpl_func,"Provide %s to compute response.",
385 ERIS_IFU_CALIB_TELL_MOD_CATALOG);
386 has_tel_mod_catalog=0;
387 }
388 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_RESP_FIT_POINTS_CATALOG) ){
389 cpl_msg_warning(cpl_func,"Provide %s to compute response.",
390 ERIS_IFU_CALIB_RESP_FIT_POINTS_CATALOG);
391 has_resp_fit_points_cat=0;
392 }
393
394 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_EXTCOEFF_TABLE) ){
395 cpl_msg_warning(cpl_func,"Provide %s to compute response.",
396 ERIS_IFU_CALIB_EXTCOEFF_TABLE);
397 has_extcoeff_table=0;
398 }
399
400 /*
401 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_HIGH_ABS_REGIONS) ){
402 cpl_msg_warning(cpl_func,"Frame %s not provided to compute response.",
403 ERIS_IFU_CALIB_HIGH_ABS_REGIONS);
404 has_high_abs_regions=1;
405 //return 0;
406 }
407 */
408
409 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_QUALITY_AREAS) ){
410 cpl_msg_warning(cpl_func,"Provide %s to compute response.",
411 ERIS_IFU_CALIB_QUALITY_AREAS);
412 has_quality_areas=0;
413 }
414
415 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_FIT_AREAS) ){
416 cpl_msg_warning(cpl_func,"Provide %s to compute response.",
417 ERIS_IFU_CALIB_FIT_AREAS);
418 has_fit_areas=0;
419 }
420
421 if (NULL == cpl_frameset_find(sof, ERIS_IFU_CALIB_RESPONSE_WINDOWS) ){
422 cpl_msg_warning(cpl_func,"Provide %s to compute response.",
423 ERIS_IFU_CALIB_RESPONSE_WINDOWS);
424 has_resp_windows=0;
425 }
426 eris_check_error_code("eris_can_compute_response");
427 return has_std_star_spectra *
428 has_flux_std_catalog *
429 has_tel_mod_catalog *
430 has_resp_fit_points_cat *
431 has_extcoeff_table *
432 has_high_abs_regions *
433 has_quality_areas *
434 has_fit_areas *
435 has_resp_windows;
436
437}
438
439
440
441static cpl_error_code
442eris_eff_qc(cpl_table* tot_eff, cpl_frame* frm_eff_wind_qc, cpl_table** qclog_tbl )
443{
444 const char* ew_name=cpl_frame_get_filename(frm_eff_wind_qc);
445 cpl_table* ew_tab=cpl_table_load(ew_name,1,0);
446 int nrow = cpl_table_get_nrow(ew_tab);
447
448 float* pwmin = cpl_table_get_data_float(ew_tab,"WMIN");
449 float* pwmax = cpl_table_get_data_float(ew_tab,"WMAX");
450 cpl_table* tmp_eff=NULL;
451 cpl_table* sto_eff=NULL;
452 double wmin;
453 double wmax;
454 double eff_avg;
455 double eff_med;
456 double eff_min;
457 double eff_max;
458 double eff_std;
459 char* key_name;
460
461 int nrow_tmp=0;
462 sto_eff = cpl_table_new(0);
463 cpl_table_copy_structure(sto_eff, tot_eff);
464
465 int ncontributes = 0;
466 int nrow_sto = 0;
467
468 for(int i=0;i<nrow;i++) {
469 //cpl_msg_info(cpl_func, "i=%d wmin=%g wmax=%g",i,pwmin[i],pwmax[i]);
470 wmin=pwmin[i];
471 wmax=pwmax[i];
472
473 cpl_table_and_selected_double(tot_eff,"WAVE",CPL_NOT_LESS_THAN,wmin);
474 cpl_table_and_selected_double(tot_eff,"WAVE",CPL_NOT_GREATER_THAN,wmax);
475 tmp_eff=cpl_table_extract_selected(tot_eff);
476 nrow_tmp = cpl_table_get_nrow(tmp_eff);
477 if(nrow_tmp>0) {
478 ncontributes++;
479 cpl_table_insert(sto_eff,tmp_eff, nrow_sto);
480 nrow_sto = cpl_table_get_nrow(sto_eff);
481 eff_avg=cpl_table_get_column_mean(tmp_eff,"EFF");
482 eff_med=cpl_table_get_column_median(tmp_eff,"EFF");
483 eff_min=cpl_table_get_column_min(tmp_eff,"EFF");
484 eff_max=cpl_table_get_column_max(tmp_eff,"EFF");
485 eff_std=cpl_table_get_column_stdev(tmp_eff,"EFF");
486
487
488 key_name = cpl_sprintf("QC EFF WIN%d WLMIN",i);
489 eris_qclog_add_double(*qclog_tbl, key_name, wmin,
490 "[um] Min window wavelength for eff comp");
491 cpl_free(key_name);
492
493 key_name = cpl_sprintf("QC EFF WIN%d WLMAX",i);
494 eris_qclog_add_double(*qclog_tbl, key_name, wmax,
495 "[um] Max window wavelength for eff comp");
496 cpl_free(key_name);
497
498 key_name = cpl_sprintf("QC EFF WIN%d MEAN",i);
499 eris_qclog_add_double(*qclog_tbl, key_name, eff_avg,
500 "Mean efficiency on window");
501 cpl_free(key_name);
502
503 key_name = cpl_sprintf("QC EFF WIN%d MEDIAN",i);
504 eris_qclog_add_double(*qclog_tbl, key_name, eff_med,
505 "Median efficiency on window");
506 cpl_free(key_name);
507
508 key_name = cpl_sprintf("QC EFF WIN%d MIN",i);
509 eris_qclog_add_double(*qclog_tbl, key_name, eff_min,
510 "Min efficiency on window");
511 cpl_free(key_name);
512
513 key_name = cpl_sprintf("QC EFF WIN%d MAX",i);
514 eris_qclog_add_double(*qclog_tbl, key_name, eff_max,
515 "Max efficiency on window");
516 cpl_free(key_name);
517
518 key_name = cpl_sprintf("QC EFF WIN%d RMS",i);
519 eris_qclog_add_double(*qclog_tbl, key_name, eff_std,
520 "RMS efficiency on window");
521 cpl_free(key_name);
522 }
523 cpl_table_select_all(tot_eff);
524 cpl_table_delete(tmp_eff);
525
526 }
527 //cpl_table_save (ew_tab, NULL, NULL, "ew_tbl.fits", CPL_IO_CREATE);
528 wmin=cpl_table_get_column_min(sto_eff,"WAVE");
529 wmax=cpl_table_get_column_max(sto_eff,"WAVE");
530
531 eff_avg=cpl_table_get_column_mean(sto_eff,"EFF");
532 eff_med=cpl_table_get_column_median(sto_eff,"EFF");
533 eff_min=cpl_table_get_column_min(sto_eff,"EFF");
534 eff_max=cpl_table_get_column_max(sto_eff,"EFF");
535 eff_std=cpl_table_get_column_stdev(sto_eff,"EFF");
536
537
538 key_name = cpl_sprintf("QC EFF NWIN");
539 eris_qclog_add_int(*qclog_tbl, key_name, ncontributes,
540 "Number of windows used for eff comp");
541 cpl_free(key_name);
542 key_name = cpl_sprintf("QC EFF WLMIN");
543 eris_qclog_add_double(*qclog_tbl, key_name, wmin,
544 "[um] Min wavelength for eff comp");
545 cpl_free(key_name);
546
547 key_name = cpl_sprintf("QC EFF WLMAX");
548 eris_qclog_add_double(*qclog_tbl, key_name, wmax,
549 "[um] Max wavelength for eff comp");
550 cpl_free(key_name);
551
552 key_name = cpl_sprintf("QC EFF MEAN");
553 eris_qclog_add_double(*qclog_tbl, key_name, eff_avg,
554 "Mean efficiency");
555 cpl_free(key_name);
556
557 key_name = cpl_sprintf("QC EFF MEDIAN");
558 eris_qclog_add_double(*qclog_tbl, key_name, eff_med,
559 "Median efficiency");
560 cpl_free(key_name);
561
562 key_name = cpl_sprintf("QC EFF MIN");
563 eris_qclog_add_double(*qclog_tbl, key_name, eff_min,
564 "Min efficiency");
565 cpl_free(key_name);
566
567 key_name = cpl_sprintf("QC EFF MAX");
568 eris_qclog_add_double(*qclog_tbl, key_name, eff_max,
569 "Max efficiency");
570 cpl_free(key_name);
571
572 key_name = cpl_sprintf("QC EFF RMS");
573 eris_qclog_add_double(*qclog_tbl, key_name, eff_std,
574 "RMS efficiency");
575
576 //cpl_table_save (*qclog_tbl, NULL, NULL, "qclog_tbl.fits", CPL_IO_CREATE);
577
578 cpl_free(key_name);
579 cpl_table_delete(ew_tab);
580 cpl_table_delete(sto_eff);
581 eris_check_error_code("eris_eff_qc");
582 return cpl_error_get_code();
583}
584
585static cpl_table*
586eris_ifu_table_from_sdp_to_normal_format(cpl_frame** frm_sci)
587{
588 cpl_table* sci_tab = NULL;
589 const cpl_array* array_wave = NULL;
590 const char* name_sci = NULL;
591 cpl_propertylist* plist = NULL;
592 name_sci = cpl_frame_get_filename(*frm_sci);
593 plist = cpl_propertylist_load(name_sci,0);
594 sci_tab = cpl_table_load(name_sci, 1, 0);
595 if(NULL != (array_wave = cpl_table_get_array(sci_tab, "WAVE", 0))) {
596 /* convert from SDP table data format to normal table format */
597 const cpl_array* array_flux = NULL;
598 const cpl_array* array_err = NULL;
599 cpl_size nrows = cpl_array_get_size(array_wave);
600
601 cpl_table* spc_tab = cpl_table_new(nrows);
602 cpl_table_new_column(spc_tab, "WAVE", CPL_TYPE_DOUBLE);
603 cpl_table_new_column(spc_tab, "FLUX", CPL_TYPE_DOUBLE);
604 cpl_table_new_column(spc_tab, "ERR", CPL_TYPE_DOUBLE);
605 array_flux = cpl_table_get_array(sci_tab, "TOT_FLUX", 0);
606 array_err = cpl_table_get_array(sci_tab, "ERR", 0);
607
608 cpl_table_copy_data_double(spc_tab, "WAVE",
609 cpl_array_get_data_double_const(array_wave));
610 cpl_table_copy_data_double(spc_tab, "FLUX",
611 cpl_array_get_data_double_const(array_flux));
612 cpl_table_copy_data_double(spc_tab, "ERR",
613 cpl_array_get_data_double_const(array_err));
614
615 cpl_table_save(spc_tab,plist,NULL,"spc_tab.fits", CPL_IO_CREATE);
616 cpl_table_delete(spc_tab);
617 cpl_table_delete(sci_tab);
618
619 sci_tab = cpl_table_load("spc_tab.fits", 1, 0);
620 cpl_frame_set_filename(*frm_sci, "spc_tab.fits");
621
622 }
623 cpl_propertylist_delete(plist);
624 eris_check_error_code("eris_ifu_table_from_sdp_to_normal_format");
625 return sci_tab;
626}
627
628static void
629eris_remove_table_normal_format(void) {
630 char* cmd = cpl_sprintf("rm spc_tab.fits");
631 int status = system(cmd);
632 if(status == -1) {
633 cpl_msg_warning(cpl_func,"call to system failed");
634 }
635 cpl_free(cmd);
636 return;
637}
638/*---------------------------------------------------------------------------*/
652/*---------------------------------------------------------------------------*/
653
654cpl_table*
655eris_efficiency_compute(cpl_frame* frm_sci, cpl_frame* frm_cat,
656 cpl_frame* frm_atmext, cpl_frameset* frameset,
657 const cpl_parameterlist* parlist, const char* pipefile_prefix)
658
659{
660 cpl_ensure(frm_sci != NULL, CPL_ERROR_NULL_INPUT, NULL);
661 cpl_ensure(frm_cat != NULL, CPL_ERROR_NULL_INPUT, NULL);
662 cpl_ensure(frm_atmext != NULL, CPL_ERROR_NULL_INPUT, NULL);
663 cpl_frame* frm_eff_wind_qc = cpl_frameset_find(frameset, ERIS_IFU_CALIB_EFFICIENCY_WINDOWS);
664 cpl_propertylist* plist = NULL;
665
666 cpl_table* tbl_eff = NULL;
667 cpl_table* tbl_ref = NULL;
668 cpl_table* tbl_atmext = NULL;
669 cpl_table* sci_tab = NULL;
670
671
672
673 double dRA = 0;
674 double dDEC = 0;
675 double gain = 0;
676 //double aimprim = 0;
677 double dEpsilon = 0.1;
678 //double um2AA = 1.e4;
679 double um2nm = 1000.;
680 double nm2um = 0.001;
681 double m2tocm2 = 1.e4;
682
683 /* int nrow=0; */
684
685 cpl_boolean is_sinfoni = CPL_FALSE;
686
687 const char* name_sci = NULL;
688 const char* name_atm = NULL;
689
690 name_sci = cpl_frame_get_filename(frm_sci);
691
692 plist = cpl_propertylist_load(name_sci, 0);
693 const char* instrume = cpl_propertylist_get_string(plist, "INSTRUME");
694 //cpl_msg_warning(cpl_func, "fname: %s instrume: %s", name_sci, instrume);
695 if(strcmp(instrume, "SINFONI") == 0) {
696 is_sinfoni = CPL_TRUE;
697 }
698 //sci_tab = cpl_table_load(name_sci, 1, 0);
699
700 dRA = cpl_propertylist_get_double(plist, "RA");
701 dDEC = cpl_propertylist_get_double(plist, "DEC");
702 cpl_frame* frm_sci_dup = cpl_frame_duplicate(frm_sci);
703 sci_tab = eris_ifu_table_from_sdp_to_normal_format(&frm_sci_dup);
704 cpl_frame_delete(frm_sci_dup);
705 double airmass = eris_ifu_get_airmass_mean(plist);
706 if(is_sinfoni) {
707 gain = SINFONI_GAIN;
708 } else {
709 gain = ERIS_GAIN; /* value to be adopted if input spectrum is in ADU */
710 /* gain = 1.0; value to be adopted if input spectrum is in e- (as py script)
711 or input spectrum is in ADU and CUBES are not flux calibrated */
712 //gain = 0.5; /* value to be adopted taking in count that CUBES are not flux
713 //calibrated when making a rectangular spaxel square */
714 }
715 //int biny=1;
716
717 double exptime = 0;
718 if(is_sinfoni) {
719 exptime = cpl_propertylist_get_double(plist, FHDR_S_DIT);
720 } else {
721 exptime = cpl_propertylist_get_double(plist, FHDR_E_DIT);
722 }
723 cpl_msg_debug(cpl_func, "gain=%g airm=%g exptime=%g ra=%g dec=%g dEpsilon=%g",
724 gain, airmass, exptime, dRA, dDEC, dEpsilon);
725
726 //cpl_msg_info(cpl_func, "table sci spectra=%s", name_sci);
727
728 name_atm = cpl_frame_get_filename(frm_atmext);
729 tbl_atmext = cpl_table_load(name_atm, 1, 0);
730 //cpl_msg_info(cpl_func,"frm_atm: %s",cpl_frame_get_filename(frm_atmext));
731 //cpl_msg_info(cpl_func,"frm_cat: %s",cpl_frame_get_filename(frm_cat));
732
733 eris_parse_catalog_std_stars(frm_cat, dRA, dDEC, dEpsilon, &tbl_ref);
734
735 if(tbl_ref == NULL) {
736 cpl_msg_error(cpl_func, "Provide std sar catalog frame");
737 cpl_msg_warning(cpl_func, "STD star not found in catalogue. Skip efficiency computation");
738 cpl_error_set(cpl_func, CPL_ERROR_NONE);
739
740 return NULL;
741 }
742
743 //cpl_table_save(sci_tab, NULL, NULL, "sci.fits", CPL_IO_DEFAULT);
744
745
746 hdrl_spectrum1D * sci_s1d;
747 hdrl_spectrum1D * std_s1d;
748 hdrl_spectrum1D * ext_s1d;
749 hdrl_spectrum1D * eff_s1d;
750 hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
751
752
753 if(is_sinfoni) {
754 sci_s1d = hdrl_spectrum1D_convert_from_table(sci_tab, "counts_bkg",
755 "wavelength", NULL, NULL, scale);
756
757 } else {
758 sci_s1d = hdrl_spectrum1D_convert_from_table(sci_tab, "FLUX",
759 "WAVE", NULL, NULL, scale);
760
761 }
762
763 //hdrl_spectrum1D_save(sci_s1d, "sci_s1d.fits");
764 /* HDRL expects results in nm */
766
767 //hdrl_spectrum1D_save(sci_s1d, "sci_s1d_nm.fits");
768 /* because the reference flux std assumes the flux in units of ergs/cm2/s/A
769 * we need to convert the spectrum to Angstroms and know the size of a bin
770 * in angstroms
771 */
772
773 char band[FILE_NAME_SZ];
774 if(is_sinfoni) {
775 sinfo_get_band(frm_sci, band);
776 } else {
777 eris_get_band(frm_sci, band);
778 }
779
780 cpl_msg_debug(cpl_func, "band: %s", band);
781
782 double conv_factor = 0;
783 double nm2AA = 10.;
784 //double um2AA = um2nm * nm2AA;
785 /* the following is the bin size in um */
786 cpl_size size = cpl_table_get_nrow(sci_tab);
787
788 int status;
789 double cdelt1 =0;
790 if(is_sinfoni) {
791 cdelt1 = cpl_table_get_float(sci_tab, "wavelength", size / 2,
792 &status) - cpl_table_get_float(sci_tab, "wavelength",
793 (size / 2 - 1), &status);
794 } else {
795 cdelt1 = cpl_table_get_double(sci_tab, "WAVE", size / 2,
796 &status) - cpl_table_get_double(sci_tab, "WAVE",
797 (size / 2 - 1), &status);
798 }
799
800 //cpl_msg_info(cpl_func, "cdelt1: %g", cdelt1);
801 double bin_size = 0;
802 if(is_sinfoni) {
803 bin_size = sinfo_get_dispersion(band);
804 } else {
805 bin_size = eris_get_dispersion(band);
806 }
807 //cpl_msg_info(cpl_func, "bin_size: %g", bin_size);
808 bin_size = cdelt1;
809 //cpl_msg_info(cpl_func, "bin_size: %g", bin_size);
810 //bin_size = 0.000394643;
811
812 conv_factor = bin_size * um2nm * nm2AA;
813 hdrl_value operator = {conv_factor, 0.};
814 //cpl_msg_info(cpl_func,"conv factor: %g", conv_factor);
815 hdrl_spectrum1D_div_scalar(sci_s1d, operator);
816 //hdrl_spectrum1D_save(sci_s1d, "sci_s1d_conv_factor.fits");
817 std_s1d = hdrl_spectrum1D_convert_from_table(tbl_ref, "FLUX",
818 "LAMBDA", NULL, NULL, scale);
819 //hdrl_spectrum1D_save(std_s1d, "std_s1d.fits");
820 ext_s1d = hdrl_spectrum1D_convert_from_table(tbl_atmext, "EXTINCTION",
821 "LAMBDA", NULL, NULL, scale);
822 //hdrl_spectrum1D_save(ext_s1d, "ext_s1d.fits");
823
824 hdrl_parameter* eff_pars;
825
826 hdrl_value Ap = {0,0};
827 hdrl_value Am = {airmass,0};
828 hdrl_value G = {gain,0};
829 hdrl_value Tex = {exptime,0};
830
831 /*HDRL expects telescope area in cm2 */
832 hdrl_value Atel = {TELESCOPE_SURFACE * m2tocm2,0};
833 //cpl_msg_info(cpl_func, "Atel: %g", TELESCOPE_SURFACE * m2tocm2);
834 eff_pars = hdrl_efficiency_parameter_create(Ap, Am, G, Tex, Atel);
835
836 eff_s1d = hdrl_efficiency_compute(sci_s1d, std_s1d, ext_s1d, eff_pars);
837 /* HDRL expects results in nm */
839
840 tbl_eff = hdrl_spectrum1D_convert_to_table(eff_s1d, "EFF", "WAVE", NULL,
841 NULL);
842
843 if( tbl_eff != NULL && frm_eff_wind_qc != NULL) {
844 cpl_table* qclog_tbl = eris_qclog_init();
845 eris_eff_qc(tbl_eff, frm_eff_wind_qc, &qclog_tbl);
846 eris_pfits_put_qc(plist, qclog_tbl);
847 cpl_table_delete(qclog_tbl);
848 }
849
850
851 cpl_propertylist_update_string(plist, CPL_DFS_PRO_CATG,
852 ERIS_IFU_PRO_JITTER_EFFICIENCY);
853 cpl_propertylist_erase(plist,"BUNIT");
854
855 char* fname = cpl_sprintf("%s_efficiency.fits", pipefile_prefix);
856 cpl_dfs_save_table(frameset, NULL, parlist, frameset, NULL,
857 tbl_eff, plist, "eris_ifu_stdstar", plist,
858 NULL, PACKAGE "/" PACKAGE_VERSION, fname);
859 cpl_free(fname);
860
861
862 hdrl_spectrum1D_delete(&std_s1d);
863 hdrl_spectrum1D_delete(&sci_s1d);
864 hdrl_spectrum1D_delete(&ext_s1d);
865 hdrl_spectrum1D_delete(&eff_s1d);
866 hdrl_parameter_delete(eff_pars);
867 eris_remove_table_normal_format();
868 cpl_table_delete(sci_tab);
869 cpl_table_delete(tbl_ref);
870 cpl_table_delete(tbl_atmext);
871 cpl_propertylist_delete(plist);
872
873 eris_check_error_code("eris_efficiency_compute");
874 return tbl_eff;
875
876}
877
878
879
880static hdrl_spectrum1D *
881eris_std_cat_frame_to_hdrl_spectrum1D(cpl_frame* frm_std_cat,const double dRA,
882 const double dDEC)
883{
884 cpl_ensure(frm_std_cat != NULL, CPL_ERROR_NULL_INPUT, NULL);
885 hdrl_spectrum1D * std_s1d;
886
887 //const char* name = cpl_frame_get_filename(frm_std_cat);
888 cpl_table* tbl_ref = NULL;
889
890 double dEpsilon=0.1;
891
892 eris_parse_catalog_std_stars(frm_std_cat,dRA,dDEC,dEpsilon,&tbl_ref);
893
894 if(tbl_ref == NULL) {
895 cpl_msg_error(cpl_func, "Provide std sar catalog frame");
896 return NULL;
897 }
898
899 std_s1d=hdrl_spectrum1D_convert_from_table(tbl_ref, "FLUX",
900 "LAMBDA", NULL, NULL, global_scale);
901 cpl_table_delete(tbl_ref);
902 eris_check_error_code("eris_std_cat_frame_to_hdrl_spectrum1D");
903 return std_s1d;
904}
905
906
907
908static hdrl_spectrum1D *
909eris_sci_frame_to_hdrl_spectrum1D(cpl_frame* frm_sci)
910{
911 cpl_ensure(frm_sci != NULL, CPL_ERROR_NULL_INPUT, NULL);
912 hdrl_spectrum1D * sci_s1d;
913 const char* name_sci=NULL;
914 name_sci=cpl_frame_get_filename(frm_sci);
915 cpl_table* sci_tab=NULL;
916 sci_tab=cpl_table_load(name_sci,1,0);
917 cpl_ensure(sci_tab != NULL, CPL_ERROR_NULL_INPUT, NULL);
918
919
920 sci_s1d=hdrl_spectrum1D_convert_from_table(sci_tab, "FLUX",
921 "WAVE", NULL, NULL, global_scale);
922 cpl_table_delete(sci_tab);
923 eris_check_error_code("eris_sci_frame_to_hdrl_spectrum1D");
924 return sci_s1d;
925
926}
927
928static hdrl_spectrum1D *
929eris_atmext_frame_to_hdrl_spectrum1D(cpl_frame* frm_atmext)
930{
931 cpl_ensure(frm_atmext != NULL, CPL_ERROR_NULL_INPUT, NULL);
932 hdrl_spectrum1D * ext_s1d;
933 const char* name_atm = cpl_frame_get_filename(frm_atmext);
934
935 cpl_table* tbl_atmext = cpl_table_load(name_atm,1,0);
936
937
938 ext_s1d=hdrl_spectrum1D_convert_from_table(tbl_atmext, "EXTINCTION",
939 "LAMBDA", NULL, NULL, global_scale);
940 cpl_table_delete(tbl_atmext);
941 eris_check_error_code("eris_atmext_frame_to_hdrl_spectrum1D");
942 return ext_s1d;
943}
944
945static hdrl_spectrum1Dlist *
946eris_get_telluric_models(const cpl_frame * telluric_cat){
947
948 cpl_ensure(telluric_cat != NULL, CPL_ERROR_NULL_INPUT, NULL);
949 const char * cat_name = cpl_frame_get_filename(telluric_cat);
950 const cpl_size next = cpl_frame_get_nextensions(telluric_cat);
951
952 hdrl_spectrum1Dlist * list = hdrl_spectrum1Dlist_new();
953
954 for(cpl_size i = 0; i < next; ++i){
955 cpl_table * tab = cpl_table_load(cat_name, 1 + i, 1);
956
957 hdrl_spectrum1D * s =
959 "flux",
960 "lam",
961 NULL,
962 NULL,
963 global_scale);
964 cpl_table_delete(tab);
965 /*um to nm*/
967 hdrl_spectrum1Dlist_set(list, s, i);
968 }
969 eris_check_error_code("eris_get_telluric_models");
970
971 return list;
972}
973
974static cpl_bivector *
975eris_read_wlen_windows(const cpl_frame * fit_areas_frm){
976
977 cpl_ensure(fit_areas_frm != NULL, CPL_ERROR_NULL_INPUT, NULL);
978 cpl_table * tab = NULL;
979
980 tab = cpl_table_load(cpl_frame_get_filename(fit_areas_frm),1,0);
981
982
983 const cpl_size nrow = cpl_table_get_nrow(tab);
984
985 double * pwmin = cpl_table_unwrap(tab,"LAMBDA_MIN");
986 double * pwmax = cpl_table_unwrap(tab,"LAMBDA_MAX");
987
988 cpl_vector * wmin = cpl_vector_wrap(nrow, pwmin);
989 cpl_vector * wmax = cpl_vector_wrap(nrow, pwmax);
990 cpl_bivector * to_ret = cpl_bivector_wrap_vectors(wmin, wmax);
991 cpl_table_delete(tab);
992
993 eris_check_error_code("eris_read_wlen_windows");
994 return to_ret;
995}
996
997
998static hdrl_parameter*
999eris_response_parameters_calc_par(cpl_propertylist* plist)
1000{
1001 cpl_ensure(plist != NULL, CPL_ERROR_NULL_INPUT, NULL);
1002 hdrl_value Ap = {0, 0};
1003
1004
1005 double airmass = eris_ifu_get_airmass_mean(plist);
1006 hdrl_value Am = {airmass, 0.};
1007
1008
1009
1010 double gain=ERIS_GAIN; // SINFONI had value of SINFONI_GAIN
1011
1012 hdrl_value G= {gain, 0.};
1013 cpl_msg_debug(cpl_func, "resp gain=%g",gain);
1014 double dit = cpl_propertylist_get_double(plist, FHDR_E_DIT);
1015 double exptime = dit;
1016 cpl_msg_debug(cpl_func, "resp exptime=%g",exptime);
1017
1018 hdrl_value Tex= {exptime, 0.};
1019 hdrl_parameter* calc_par = hdrl_response_parameter_create(Ap, Am, G, Tex);
1020 eris_check_error_code("eris_response_parameters_calc_par");
1021 return calc_par;
1022
1023}
1024
1025
1026static cpl_array *
1027eris_read_fit_points(const cpl_frame * frm_fit_points,
1028 const double ra, const double dec, const double ra_dec_tolerance){
1029
1030 cpl_ensure(frm_fit_points != NULL, CPL_ERROR_NULL_INPUT, NULL);
1031 const char * fname = cpl_frame_get_filename(frm_fit_points);
1032
1033 cpl_table * index_table = cpl_table_load(fname, 1, CPL_FALSE);
1034 const cpl_size sz = cpl_table_get_nrow(index_table);
1035
1036 cpl_size selected_ext = -1;
1037
1038 for(cpl_size i = 0; i < sz; ++i){
1039 int rej;
1040 const int ext_id = cpl_table_get_int(index_table, "ext_id", i ,&rej);
1041 const double curr_ra = cpl_table_get(index_table, "ra", i, &rej);
1042 const double curr_dec = cpl_table_get(index_table, "dec", i, &rej);
1043 if ((ext_id > 0) && (fabs(curr_ra - ra) < ra_dec_tolerance) &&
1044 (fabs(curr_dec - dec) < ra_dec_tolerance)){
1045 selected_ext = ext_id;
1046 }
1047 }
1048
1049 cpl_table_delete(index_table);
1050 cpl_ensure(selected_ext >= 0, CPL_ERROR_ILLEGAL_OUTPUT, NULL);
1051
1052
1053 cpl_table * tb = cpl_table_load(fname, selected_ext, CPL_FALSE);
1054
1055 const cpl_size sz_points = cpl_table_get_nrow(tb);
1056 cpl_array * fit_points = cpl_array_new(sz_points, CPL_TYPE_DOUBLE);
1057
1058 for(cpl_size i = 0; i < sz_points; ++i){
1059 int rej = 0;
1060 const double d = cpl_table_get(tb, "LAMBDA", i, &rej);
1061 cpl_array_set(fit_points, i, d);
1062 }
1063
1064 cpl_table_delete(tb);
1065 eris_check_error_code("eris_read_fit_points");
1066 return fit_points;
1067}
1068
1069static inline cpl_size
1070get_before(const double wmin_spectrum, const double wmin_interval, const double step){
1071 const double num = (wmin_spectrum - wmin_interval) / step;
1072 if(num <= 0) return 0;
1073 return (cpl_size) ceil(num);
1074}
1075
1076static inline cpl_size
1077get_after(const double wmax_spectrum, const double wmax_interval, const double step){
1078 const double num = (wmax_interval - wmax_spectrum) / step;
1079 if(num <= 0) return 0;
1080 return (cpl_size) ceil(num);
1081}
1082static inline
1083hdrl_spectrum1D * extend_hdrl_spectrum(const double wmin, const double wmax,
1084 const hdrl_spectrum1D * s){
1085
1086 cpl_ensure(wmin < wmax, CPL_ERROR_ILLEGAL_INPUT, NULL);
1087 cpl_ensure(s != NULL, CPL_ERROR_NULL_INPUT, NULL);
1088
1089 const cpl_array * lambadas_s = hdrl_spectrum1D_get_wavelength(s).wavelength;
1090 const double s_wmin = cpl_array_get_min(lambadas_s);
1091 const double s_wmax = cpl_array_get_max(lambadas_s);
1092
1093 const cpl_size sz_s = cpl_array_get_size(lambadas_s);
1094 const double step = (s_wmax - s_wmin)/sz_s;
1095
1096 cpl_msg_debug(cpl_func, "min=%g max=%g step=%g",s_wmin,s_wmax,step);
1097 cpl_size n_before = get_before(s_wmin, wmin, step);
1098 cpl_size n_after = get_after(s_wmax, wmax, step);
1099
1100 const cpl_size new_size = n_before + n_after + sz_s;
1101 cpl_array * new_lambdas = cpl_array_new(new_size,
1102 cpl_array_get_type(lambadas_s));
1103 hdrl_image * new_flux = hdrl_image_new(new_size, 1);
1104
1105 double lambda = s_wmin;
1106 for(cpl_size i = n_before - 1; i >= 0; i--){
1107 lambda -= step;
1108 cpl_array_set(new_lambdas, i, lambda);
1109 hdrl_image_reject(new_flux, i + 1, 1);
1110 }
1111
1112 lambda = s_wmax;
1113 for(cpl_size i = n_before + sz_s; i < new_size; i++){
1114 lambda += step;
1115 cpl_array_set(new_lambdas, i, lambda);
1116 hdrl_image_reject(new_flux, i + 1, 1);
1117 }
1118
1119 for(cpl_size i = n_before; i < n_before + sz_s; i++){
1120 const cpl_size idx_ori = i - n_before;
1121 int rej = 0;
1122 lambda = hdrl_spectrum1D_get_wavelength_value(s, idx_ori, &rej);
1123 cpl_array_set(new_lambdas, i, lambda);
1124
1125 if(rej) {
1126 hdrl_image_reject(new_flux, i + 1, 1);
1127 continue;
1128 }
1129
1130 hdrl_value val = hdrl_spectrum1D_get_flux_value(s, idx_ori, NULL);
1131 hdrl_image_set_pixel(new_flux, i + 1, 1, val);
1132 }
1133
1134 const hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_get_scale(s);
1135
1136 hdrl_spectrum1D * to_ret =
1138 hdrl_image_get_error(new_flux),
1139 new_lambdas,
1140 scale);
1141
1142
1143 cpl_array_delete(new_lambdas);
1144 hdrl_image_delete(new_flux);
1145
1146 eris_check_error_code("hdrl_response_fit_parameter_create");
1147 return to_ret;
1148}
1149
1150static cpl_error_code
1151eris_resp_qc(cpl_table* resp, cpl_frame* frm_eff_wind_qc, const char* col1,
1152 const char* col2, const double wmin_g,
1153 const double wmax_g, cpl_table** qclog_tbl )
1154{
1155
1156 cpl_ensure(resp != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1157 cpl_ensure(frm_eff_wind_qc != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1158 cpl_ensure(wmin_g < wmax_g, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
1159 cpl_ensure(*qclog_tbl != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1160
1161 const char* ew_name=cpl_frame_get_filename(frm_eff_wind_qc);
1162 cpl_table* ew_tab=cpl_table_load(ew_name,1,0);
1163 int nrow = cpl_table_get_nrow(ew_tab);
1164
1165 float* pwmin = cpl_table_get_data_float(ew_tab,"WMIN");
1166 float* pwmax = cpl_table_get_data_float(ew_tab,"WMAX");
1167 cpl_table* tmp_resp=NULL;
1168 double wmin;
1169 double wmax;
1170 double resp_avg;
1171 double resp_med;
1172 double resp_min;
1173 double resp_max;
1174 double resp_std;
1175 char* key_name;
1176 //double um2nm=1000;
1177 cpl_table* sto_resp = cpl_table_new(0);
1178 cpl_table_copy_structure(sto_resp, resp);
1179 int nrow_tmp=0;
1180 int ncontributes = 0;
1181 int nrow_sto = 0;
1182 for(int i=0;i<nrow;i++) {
1183 cpl_table_select_all(resp);
1184 wmin=pwmin[i];
1185 wmax=pwmax[i];
1186 cpl_table_and_selected_double(resp, col1, CPL_NOT_LESS_THAN, wmin_g);
1187 cpl_table_and_selected_double(resp, col1, CPL_NOT_LESS_THAN, wmin);
1188 cpl_table_and_selected_double(resp, col1, CPL_NOT_GREATER_THAN, wmax);
1189 cpl_table_and_selected_double(resp, col1, CPL_NOT_GREATER_THAN, wmax_g);
1190 tmp_resp=cpl_table_extract_selected(resp);
1191 nrow_tmp = cpl_table_get_nrow(tmp_resp);
1192 if(nrow_tmp>0) {
1193 ncontributes++;
1194 cpl_table_insert(sto_resp,tmp_resp, nrow_sto);
1195 nrow_sto = cpl_table_get_nrow(sto_resp);
1196 resp_avg=cpl_table_get_column_mean(tmp_resp, col2);
1197 resp_med=cpl_table_get_column_median(tmp_resp, col2);
1198 resp_min=cpl_table_get_column_min(tmp_resp, col2);
1199 resp_max=cpl_table_get_column_max(tmp_resp, col2);
1200 resp_std=cpl_table_get_column_stdev(tmp_resp, col2);
1201 /*
1202 cpl_msg_info(cpl_func, "wmin=%g wmax=%g avg=%g med=%g min=%g max=%g std=%g",
1203 wmin, wmax, resp_avg,resp_med,resp_min,resp_max,resp_std);
1204 */
1205 key_name = cpl_sprintf("QC RESP WIN%d WLMIN",i);
1206 eris_qclog_add_double(*qclog_tbl, key_name, wmin,
1207 "[um] Min window wavelength for resp comp");
1208 cpl_free(key_name);
1209
1210 key_name = cpl_sprintf("QC RESP WIN%d WLMAX",i);
1211 eris_qclog_add_double(*qclog_tbl, key_name, wmax,
1212 "[um] Max window wavelength for resp comp");
1213 cpl_free(key_name);
1214
1215
1216 key_name = cpl_sprintf("QC RESP WIN%d MEAN",i);
1217 eris_qclog_add_double(*qclog_tbl, key_name, resp_avg,
1218 "Mean response on window");
1219 cpl_free(key_name);
1220
1221 key_name = cpl_sprintf("QC RESP WIN%d MEDIAN",i);
1222 eris_qclog_add_double(*qclog_tbl, key_name, resp_med,
1223 "Median response on window");
1224 cpl_free(key_name);
1225
1226 key_name = cpl_sprintf("QC RESP WIN%d MIN",i);
1227 eris_qclog_add_double(*qclog_tbl, key_name, resp_min,
1228 "Min response on window");
1229 cpl_free(key_name);
1230
1231 key_name = cpl_sprintf("QC RESP WIN%d MAX",i);
1232 eris_qclog_add_double(*qclog_tbl, key_name, resp_max,
1233 "Max response on window");
1234 cpl_free(key_name);
1235
1236 key_name = cpl_sprintf("QC RESP WIN%d RMS",i);
1237 eris_qclog_add_double(*qclog_tbl, key_name, resp_std,
1238 "RMS response on window");
1239 cpl_free(key_name);
1240 }
1241 cpl_table_select_all(resp);
1242 cpl_table_delete(tmp_resp);
1243
1244 }
1245 cpl_table_select_all(resp);
1246 wmin=cpl_table_get_column_min(sto_resp, col1);
1247 wmax=cpl_table_get_column_max(sto_resp, col1);
1248 cpl_table_and_selected_double(sto_resp, col1, CPL_NOT_LESS_THAN, wmin_g);
1249 cpl_table_and_selected_double(sto_resp, col1, CPL_NOT_GREATER_THAN, wmax_g);
1250 tmp_resp=cpl_table_extract_selected(sto_resp);
1251
1252 resp_avg=cpl_table_get_column_mean(sto_resp, col2);
1253 resp_med=cpl_table_get_column_median(sto_resp, col2);
1254 resp_min=cpl_table_get_column_min(sto_resp, col2);
1255 resp_max=cpl_table_get_column_max(sto_resp, col2);
1256 resp_std=cpl_table_get_column_stdev(sto_resp, col2);
1257 /*
1258 cpl_msg_info(cpl_func,"wmin=%g wmax=%g avg=%g med=%g min=%g max=%g std=%g",
1259 wmin, wmax,resp_avg,resp_med,resp_min,resp_max,resp_std);
1260 */
1261 key_name = cpl_sprintf("QC RESP NWIN");
1262 eris_qclog_add_int(*qclog_tbl, key_name, ncontributes,
1263 "Number of windows used for resp comp");
1264 cpl_free(key_name);
1265 key_name = cpl_sprintf("QC RESP WLMIN");
1266 eris_qclog_add_double(*qclog_tbl, key_name, wmin,
1267 "[um] Min wavelength for resp comp");
1268 cpl_free(key_name);
1269
1270 key_name = cpl_sprintf("QC RESP WLMAX");
1271 eris_qclog_add_double(*qclog_tbl, key_name, wmax,
1272 "[um] Max wavelength for resp comp");
1273 cpl_free(key_name);
1274
1275 key_name = cpl_sprintf("QC RESP MEAN");
1276 eris_qclog_add_double(*qclog_tbl, key_name, resp_avg,
1277 "Mean response");
1278 cpl_free(key_name);
1279
1280 key_name = cpl_sprintf("QC RESP MEDIAN");
1281 eris_qclog_add_double(*qclog_tbl, key_name, resp_med,
1282 "Median response");
1283 cpl_free(key_name);
1284
1285 key_name = cpl_sprintf("QC RESP MIN");
1286 eris_qclog_add_double(*qclog_tbl, key_name, resp_min,
1287 "Min response");
1288 cpl_free(key_name);
1289
1290 key_name = cpl_sprintf("QC RESP MAX");
1291 eris_qclog_add_double(*qclog_tbl, key_name, resp_max,
1292 "Max response");
1293 cpl_free(key_name);
1294
1295 key_name = cpl_sprintf("QC RESP RMS");
1296 eris_qclog_add_double(*qclog_tbl, key_name, resp_std,
1297 "RMS response");
1298 cpl_free(key_name);
1299
1300 cpl_table_delete(tmp_resp);
1301 cpl_table_delete(ew_tab);
1302 cpl_table_delete(sto_resp);
1303 eris_check_error_code("eris_resp_qc");
1304 return cpl_error_get_code();
1305}
1306
1307
1308static void
1309eris_log_pro(char* name_o,
1310 const char* pro_catg,
1311 int frm_type,
1312 cpl_frameset* ref_set,
1313 cpl_frameset** out_set,
1314 cpl_propertylist** plist,
1315 const cpl_parameterlist* parlist,
1316 const char* recid)
1317{
1318 cpl_frame* product_frame = NULL ;
1319 char * pipe_id=NULL;
1320 cpl_errorstate initial_errorstate = cpl_errorstate_get();
1321
1322 pipe_id = cpl_sprintf("%s%s","eris/",PACKAGE_VERSION);
1323 product_frame = cpl_frame_new() ;
1324
1325 cpl_frame_set_filename(product_frame, name_o) ;
1326 cpl_frame_set_tag(product_frame, pro_catg) ;
1327 cpl_frame_set_type(product_frame, frm_type);
1328 cpl_frame_set_group(product_frame, CPL_FRAME_GROUP_PRODUCT);
1329 cpl_frame_set_level(product_frame, CPL_FRAME_LEVEL_FINAL);
1330
1331
1332 if(cpl_dfs_setup_product_header(*plist,product_frame,ref_set,parlist,recid,
1333 pipe_id,KEY_VALUE_HPRO_DID,NULL) != CPL_ERROR_NONE) {
1334 cpl_msg_warning(cpl_func, "Problem in the product DFS-compliance");
1335 cpl_msg_warning(cpl_func, "%s", (char* ) cpl_error_get_message());
1336 cpl_errorstate_dump(initial_errorstate, CPL_FALSE, NULL);
1337 cpl_error_reset();
1338 }
1339
1340 cpl_frameset_insert(*out_set, product_frame);
1341
1342 cpl_free(pipe_id);
1343 eris_check_error_code("eris_log_pro");
1344 return;
1345}
1346
1347
1348
1349static char * eris_new_get_rootname(const char * filename)
1350{
1351 static char path[MAX_STR_SIZE+1];
1352 char * lastdot ;
1353
1354 if (strlen(filename)>MAX_STR_SIZE) return NULL ;
1355 memset(path, 0, MAX_STR_SIZE);
1356 strcpy(path, filename);
1357 lastdot = strrchr(path, '.');
1358 if (lastdot == NULL) return path ;
1359 if ((!strcmp(lastdot, ".fits")) || (!strcmp(lastdot, ".FITS"))
1360 || (!strcmp(lastdot, ".paf")) || (!strcmp(lastdot, ".PAF"))
1361 || (!strcmp(lastdot, ".dat")) || (!strcmp(lastdot, ".DAT"))
1362 || (!strcmp(lastdot, ".fits")) || (!strcmp(lastdot, ".TFITS"))
1363 || (!strcmp(lastdot, ".ascii"))
1364 || (!strcmp(lastdot, ".ASCII")))
1365 {
1366 lastdot[0] = (char)0;
1367 }
1368 return path ;
1369}
1370
1371
1372
1373
1374static void
1375eris_check_name(const char* in, char** ou, /*int type, */char** paf) {
1376
1377 char* name_b;
1378 if (strstr(in, "." ) != NULL ) {
1379 char* tmp = eris_new_get_rootname(in);
1380 name_b= cpl_sprintf("%s",tmp);
1381 } else {
1382 name_b = cpl_sprintf("%s",in);
1383 }
1384 *ou = cpl_sprintf("%s.fits",name_b);
1385 /*
1386 if (type == CPL_FRAME_TYPE_TABLE) {
1387 strcat(*ou,".fits");
1388 } else {
1389 strcat(*ou,".fits");
1390 }
1391 */
1392 *paf = cpl_sprintf("%s.paf",name_b);
1393 cpl_free(name_b);
1394
1395}
1396
1397
1398static int
1399eris_pro_save_tbl(
1400 cpl_table * table,
1401 cpl_frameset * ref,
1402 cpl_frameset * set,
1403 const char * out_file,
1404 const char * pro_catg,
1405 cpl_table * qclog,
1406 const char * recid,
1407 const cpl_parameterlist* parlist)
1408
1409{
1410 char * name_o =NULL;
1411 char * name_p =NULL;
1412 cpl_propertylist * plist=NULL ;
1413 char* ref_file=NULL;
1414 /* Get the reference file */
1415
1416 cpl_frameset_iterator* it = cpl_frameset_iterator_new(ref);
1417 cpl_frame *first_frame = cpl_frameset_iterator_get(it);
1418 ref_file = cpl_strdup(cpl_frame_get_filename(first_frame)) ;
1419
1420 //name_o = cpl_malloc(FILE_NAME_SZ * sizeof(char));
1421 //name_p = cpl_malloc(FILE_NAME_SZ * sizeof(char));
1422 eris_check_name(out_file, &name_o, /*CPL_FRAME_TYPE_TABLE, */&name_p);
1423 cpl_msg_info(cpl_func, "Writing tbl %s pro catg %s" , name_o, pro_catg) ;
1424
1425 /* Get FITS header from reference file */
1426 if ((cpl_error_code)((plist = cpl_propertylist_load(ref_file,0)) == NULL))
1427 {
1428 cpl_msg_error(cpl_func, "getting header from tbl frame %s",ref_file);
1429 cpl_propertylist_delete(plist) ;
1430 cpl_free(ref_file);
1431 cpl_free(name_o);
1432 cpl_free(name_p);
1433 cpl_frameset_iterator_delete(it);
1434 return -1 ;
1435 }
1436
1437 //Remove CHECKSUM & ESO PRO keywords: Why?
1438 //cpl_propertylist_erase_regexp(plist, "CHECKSUM",0);
1439 //cpl_propertylist_erase_regexp(plist, "^ESO PRO .*",0);
1440
1441 /* Add DataFlow keywords and log the saved file in the input frameset */
1442 eris_log_pro(name_o, pro_catg, CPL_FRAME_TYPE_TABLE,
1443 ref,&set,&plist,parlist,recid);
1444
1445 if(qclog != NULL) {
1446 eris_pfits_put_qc(plist, qclog) ;
1447 }
1448 cpl_propertylist_append_string(plist, "PRODCATG", "SCIENCE.SPECTRUM");
1449 /* Save the file */
1450 if (cpl_table_save(table, plist, NULL, name_o, CPL_IO_DEFAULT)
1451 != CPL_ERROR_NONE) {
1452 cpl_msg_error(cpl_func, "Cannot save the product: %s", name_o);
1453 cpl_propertylist_delete(plist) ;
1454 cpl_free(ref_file);
1455 cpl_free(name_o);
1456 cpl_free(name_p);
1457 cpl_frameset_iterator_delete(it);
1458
1459 return -1 ;
1460 }
1461
1462 cpl_propertylist_delete(plist) ;
1463 cpl_msg_indent_less() ;
1464 cpl_free(name_o);
1465 cpl_free(name_p);
1466 cpl_free(ref_file);
1467 cpl_frameset_iterator_delete(it);
1468
1469 eris_check_error_code("eris_pro_save_tbl");
1470 return 0 ;
1471}
1472
1473
1474/*----------------------------------------------------------------------------*/
1486/*----------------------------------------------------------------------------*/
1487
1488static cpl_error_code
1489eris_fit_poly (eris_poly_data *data, const int deg, double *fit_RE,
1490 double *coeffs_RE, double *coeffs_err_RE, double *chisq_RE,
1491 const int print_flag) {
1492
1493 cpl_ensure(data != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1494 cpl_ensure(deg > 1, CPL_ERROR_ILLEGAL_INPUT, CPL_ERROR_ILLEGAL_INPUT);
1495
1496 int i, j, n;
1497 double err, chisq;
1498 gsl_vector *x, *y, *w, *p;
1499 gsl_matrix *X, *covar;
1500 //double *fit;
1501 gsl_multifit_linear_workspace *work;
1502
1503 n = data->n;
1504 x = gsl_vector_alloc (n);
1505 y = gsl_vector_alloc (n);
1506 w = gsl_vector_alloc (n);
1507 p = gsl_vector_alloc (deg);
1508 X = gsl_matrix_alloc (n, deg);
1509 covar = gsl_matrix_alloc (deg, deg);
1510
1511 if (print_flag) {
1512 cpl_msg_info(cpl_func,"n = %d, deg = %d", n, deg);
1513 }
1514
1515 for (i = 0; i < n; i++) {
1516 gsl_vector_set (x, i, data->x[i]);
1517 gsl_vector_set (y, i, data->y[i]);
1518 err = data->err[i];
1519 gsl_vector_set (w, i, 1.0/err/err);
1520
1521 for (j = 0; j < deg; j++)
1522 gsl_matrix_set (X, i, j, gsl_pow_int(gsl_vector_get(x, i), j));
1523 }
1524
1525 if (print_flag) {
1526 cpl_msg_info(cpl_func, "vectors set");
1527 }
1528
1529 work = gsl_multifit_linear_alloc (n, deg);
1530 gsl_multifit_wlinear (X, w, y, p, covar, &chisq, work);
1531 gsl_multifit_linear_free (work);
1532
1533 if (print_flag) {
1534 cpl_msg_info(cpl_func, "Fit done");
1535 }
1536
1537 for (i = 0; i < n; i++) {
1538 fit_RE[i] = 0.;
1539 for (j = 0; j < deg; j++)
1540 fit_RE[i] += gsl_matrix_get (X, i, j) * gsl_vector_get (p, j);
1541 }
1542
1543 chisq = chisq/(n-deg);
1544
1545 for (j = 0; j < deg; j++) {
1546 coeffs_RE[j] = gsl_vector_get (p, j);
1547 coeffs_err_RE[j] = sqrt(gsl_matrix_get (covar, j, j));
1548 }
1549
1550 *chisq_RE = chisq;
1551
1552 gsl_vector_free(x);
1553 gsl_vector_free(y);
1554 gsl_vector_free(w);
1555 gsl_vector_free(p);
1556 gsl_matrix_free(X);
1557 gsl_matrix_free(covar);
1558
1559 eris_check_error_code("eris_fit_poly");
1560 return cpl_error_get_code();
1561}
1562
1563static cpl_table*
1564eris_fit_poly_and_clip_response(const cpl_table* resp_pts,
1565 const cpl_table* resp_smo, const int deg, const double kappa,
1566 const int niter)
1567{
1568
1569 cpl_ensure(resp_pts != NULL, CPL_ERROR_NULL_INPUT, NULL);
1570 cpl_ensure(resp_smo != NULL, CPL_ERROR_NULL_INPUT, NULL);
1571
1572 /* polynomial fit */
1573 eris_poly_data *p_data = (eris_poly_data *)cpl_malloc(sizeof(eris_poly_data));
1574
1575 cpl_size pol_deg = deg;
1576 double chisq = 0;
1577 cpl_table* resp_pts_tmp = cpl_table_duplicate(resp_pts);
1578 double* fit =NULL;
1579 double* coeffs = NULL;
1580 double* coeffs_err = NULL;
1581 double y2, x2;
1582 //double y2err;
1583 cpl_table* xtab = NULL;
1584 for(cpl_size k= 0; k < niter; k++) {
1585 cpl_size nrows = cpl_table_get_nrow(resp_pts_tmp);
1586 cpl_msg_debug(cpl_func,"k: %lld nrows: %lld", k,nrows);
1587 fit = (double *) cpl_calloc (nrows , sizeof(double));
1588 coeffs = (double *) cpl_calloc (pol_deg , sizeof(double));
1589 coeffs_err = (double *) cpl_calloc (pol_deg , sizeof(double));
1590
1591 p_data->x = (double *) cpl_calloc(nrows,sizeof(double));
1592 p_data->y = (double *) cpl_calloc(nrows,sizeof(double));
1593 p_data->err = (double *) cpl_calloc(nrows,sizeof(double));
1594 p_data->n = nrows;
1595
1596 for (cpl_size i = 0; i < nrows; i++) {
1597 p_data->x[i] = cpl_table_get(resp_pts_tmp,"wavelength",i,NULL);
1598 p_data->y[i] = cpl_table_get(resp_pts_tmp,"response_fit_pts",i,NULL);
1599 p_data->err[i] = 1.;
1600 }
1601
1602 eris_fit_poly(p_data, pol_deg, fit, coeffs, coeffs_err, &chisq, 0);
1603
1604 // determine residuals to poly fit
1605 cpl_table_new_column(resp_pts_tmp, "res", CPL_TYPE_DOUBLE);
1606 for (int j = 0; j < nrows; j++) {
1607 x2 = cpl_table_get(resp_pts_tmp, "wavelength", j, NULL);
1608 y2 = 0.;
1609 for (int c = pol_deg-1; c >= 0; c--) {
1610 y2 = y2 * x2 + coeffs[c];
1611 }
1612 cpl_table_set_double(resp_pts_tmp, "res",j, y2-p_data->y[j]);
1613 }
1614 //cpl_table_save (resp_pts_tmp, NULL, NULL, "resp_pts_tmp.fits", CPL_IO_CREATE);
1615
1616 // clip data for which residual is outside 3 sigma from fit
1617 double stddev = cpl_table_get_column_stdev(resp_pts_tmp, "res");
1618 double mean = cpl_table_get_column_mean(resp_pts_tmp, "res");
1619 cpl_table_duplicate_column(resp_pts_tmp,"dist",resp_pts_tmp,"res");
1620 cpl_table_subtract_scalar(resp_pts_tmp,"dist",mean);
1621 cpl_table_abs_column(resp_pts_tmp,"dist");
1622 cpl_msg_debug(cpl_func,"mean: %g, stdev: %g",mean, stddev);
1623 cpl_table_and_selected_double(resp_pts_tmp, "dist", CPL_LESS_THAN,
1624 kappa * stddev);
1625 xtab = cpl_table_extract_selected(resp_pts_tmp);
1626 //cpl_table_save (xtab, NULL, NULL, "xtab.fits", CPL_IO_CREATE);
1627 cpl_table_delete(resp_pts_tmp);
1628 cpl_table_erase_column(xtab,"res");
1629 cpl_table_erase_column(xtab,"dist");
1630 resp_pts_tmp = cpl_table_duplicate(xtab);
1631
1632 cpl_free(coeffs);
1633 cpl_free(coeffs_err);
1634 cpl_free(fit);
1635 cpl_table_delete(xtab);
1636 cpl_free(p_data->x);
1637 cpl_free(p_data->y);
1638 cpl_free(p_data->err);
1639 }
1640
1641 cpl_size next = cpl_table_get_nrow(resp_pts_tmp);
1642
1643 /* now repeat fit after clipping of outliers */
1644
1645 //cpl_free(p_data);
1646 p_data->x = (double *) cpl_calloc(next,sizeof(double));
1647 p_data->y = (double *) cpl_calloc(next,sizeof(double));
1648 p_data->err = (double *) cpl_calloc(next,sizeof(double));
1649 p_data->n = next;
1650
1651 for (cpl_size i = 0; i < next; i++) {
1652 p_data->x[i] = cpl_table_get(resp_pts_tmp,"wavelength",i,NULL);
1653 p_data->y[i] = cpl_table_get(resp_pts_tmp,"response_fit_pts",i,NULL);
1654 p_data->err[i] = 1.;
1655 }
1656 //cpl_free(fit);
1657 //cpl_free(coeffs);
1658 //cpl_free(coeffs_err);
1659
1660 /* repeat fit after kappa-sigma clip of data */
1661 fit = (double *) cpl_calloc (next , sizeof(double));
1662 coeffs = (double *) cpl_calloc (pol_deg , sizeof(double));
1663 coeffs_err = (double *) cpl_calloc (pol_deg , sizeof(double));
1664 eris_fit_poly(p_data, pol_deg, fit, coeffs, coeffs_err, &chisq, 0);
1665
1666 /* now create final result */
1667 cpl_table* resp_poly = cpl_table_duplicate(resp_smo);
1668 cpl_size poly_nval = cpl_table_get_nrow(resp_poly);
1669 cpl_table_name_column(resp_poly, "response_smo", "response_poly");
1670 cpl_table_name_column(resp_poly, "response_smo_error", "response_poly_error");
1671 cpl_table_fill_column_window_double(resp_poly,"response_poly_error",0, poly_nval, 0);
1672 cpl_table_name_column(resp_poly, "response_smo_bpm", "response_poly_bpm");
1673
1674 for (int j = 0; j < poly_nval; j++) {
1675 x2 = cpl_table_get(resp_poly, "wavelength", j, NULL);
1676 //cpl_table_set_double(*s1d_eff_RE, "wavelength", j, x2);
1677 y2 = 0.;
1678 for (int c = pol_deg-1; c >= 0; c--) {
1679 y2 = y2 * x2 + coeffs[c];
1680 }
1681
1682 //for (int c = pol_deg-1; c >= 0; c--) {
1683 // y2err = y2err * x2 + coeffs_err[c];
1684 //}
1685
1686 cpl_table_set_double(resp_poly, "response_poly",j, y2);
1687 //cpl_table_set_double(resp_poly, "response_poly_error",j, y2err);
1688 }
1689
1690 /* free memory */
1691 cpl_free(fit);
1692 cpl_free(coeffs);
1693 cpl_free(coeffs_err);
1694 cpl_table_delete(resp_pts_tmp);
1695 cpl_free(p_data->x);
1696 cpl_free(p_data->y);
1697 cpl_free(p_data->err);
1698 cpl_free(p_data);
1699
1700 eris_check_error_code("eris_fit_poly_and_clip_response");
1701 return resp_poly;
1702}
1703
1704static cpl_error_code
1705eris_get_wave_shift_param(const char* band, double* wmin, double* wmax,
1706 cpl_boolean* enable)
1707{
1708
1709 cpl_msg_info(cpl_func, "observed band: %s", band);
1710 if ((strcmp(band, "J_low") == 0) ||
1711 (strcmp(band, "J_short") == 0) ||
1712 (strcmp(band, "J_middle") == 0) ||
1713 (strcmp(band, "J_long") == 0)) {
1714 //J
1715 *wmin = 7.001; //log scale
1716 *wmax = 7.258;
1717 *enable = CPL_FALSE; // was CPL_TRUE;
1718 } else if ((strcmp(band, "H_low") == 0) ||
1719 (strcmp(band, "H_short") == 0) ||
1720 (strcmp(band, "H_middle") == 0) ||
1721 (strcmp(band, "H_long") == 0)) {
1722 //H
1723 *wmin = 7.268; //log scale
1724 *wmax = 7.531;
1725 *enable = CPL_FALSE;
1726 } else if ((strcmp(band, "K_low") == 0) ||
1727 (strcmp(band, "K_short") == 0) ||
1728 (strcmp(band, "K_middle") == 0) ||
1729 (strcmp(band, "K_long") == 0)) {
1730 //K
1731 *wmin = 7.565; //log scale
1732 *wmax = 7.812;
1733 *enable = CPL_FALSE;
1734 } else if (strcmp(band, "H+K") == 0) {
1735 //HK
1736 *wmin = 7.268; //log scale
1737 *wmax = 7.812;
1738 *enable = CPL_FALSE;
1739 }
1740
1741 eris_check_error_code("eris_get_wave_shift_param");
1742 return cpl_error_get_code();
1743}
1744
1745/*---------------------------------------------------------------------------*/
1760/*---------------------------------------------------------------------------*/
1761
1762
1763cpl_error_code
1764eris_response_compute(const char* plugin_id, const cpl_parameterlist* config,
1765 cpl_frameset* ref_set,cpl_frameset* sof)
1766{
1767
1768 cpl_ensure_code(plugin_id, CPL_ERROR_NULL_INPUT);
1769 cpl_ensure_code(config, CPL_ERROR_NULL_INPUT);
1770 cpl_ensure_code(ref_set, CPL_ERROR_NULL_INPUT);
1771 cpl_ensure_code(sof, CPL_ERROR_NULL_INPUT);
1772
1773 cpl_frame* frm_sci = cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_SPECTRUM);
1774 cpl_frame* frm_std_cat = cpl_frameset_find(sof, ERIS_IFU_CALIB_FLUX_STD_CATALOG);
1775 cpl_frame* frm_atmext = cpl_frameset_find(sof, ERIS_IFU_CALIB_EXTCOEFF_TABLE);
1776 cpl_frame* frm_tell_cat = cpl_frameset_find(sof, ERIS_IFU_CALIB_TELL_MOD_CATALOG);
1777 cpl_frame* frm_resp_fit_points = cpl_frameset_find(sof, ERIS_IFU_CALIB_RESP_FIT_POINTS_CATALOG);
1778 //cpl_frame* frm_high_abs_reg = cpl_frameset_find(sof, ERIS_IFU_CALIB_HIGH_ABS_REGIONS);
1779 cpl_frame* frm_resp_quality_areas = cpl_frameset_find(sof, ERIS_IFU_CALIB_QUALITY_AREAS);
1780 cpl_frame* frm_resp_fit_areas = cpl_frameset_find(sof, ERIS_IFU_CALIB_FIT_AREAS);
1781
1782 cpl_ensure(frm_sci != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1783 cpl_ensure(frm_std_cat != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1784 cpl_ensure(frm_atmext != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1785 cpl_ensure(frm_tell_cat != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1786 cpl_ensure(frm_resp_fit_points != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1787 /* HIGH_ABS_REGIONS is an optional input
1788 cpl_ensure(frm_high_abs_reg != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1789 */
1790 cpl_ensure(frm_resp_quality_areas != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1791 cpl_ensure(frm_resp_fit_areas != NULL, CPL_ERROR_NULL_INPUT, CPL_ERROR_NULL_INPUT);
1792
1793 const char* name_sci = NULL;
1794 cpl_propertylist* plist = NULL;
1795 cpl_table* sci_tab = NULL;
1796 name_sci = cpl_frame_get_filename(frm_sci);
1797 cpl_msg_debug(cpl_func,"fname: %s",name_sci);
1798 plist = cpl_propertylist_load(name_sci,0);
1799
1800 double dRA = cpl_propertylist_get_double(plist, "RA");
1801
1802 double dDEC = cpl_propertylist_get_double(plist, "DEC");
1803 //sci_tab = cpl_table_load(name_sci, 1, 0);
1804
1805
1806
1807 double sci_wmin = 0;
1808 double sci_wmax = 0;
1809 cpl_frame* frm_sci_dup = cpl_frame_duplicate(frm_sci);
1810 //scpl_table_delete(sci_tab);
1811 sci_tab = eris_ifu_table_from_sdp_to_normal_format(&frm_sci_dup);
1812
1813 sci_wmin = cpl_table_get_column_min(sci_tab,"WAVE");
1814 sci_wmax = cpl_table_get_column_max(sci_tab, "WAVE");
1815
1816 double um2nm = 1.e3;
1817
1818 hdrl_spectrum1D * sci_s1d = eris_sci_frame_to_hdrl_spectrum1D(frm_sci_dup);
1819
1820 cpl_frame_delete(frm_sci_dup);
1821 /* HDRL expects results in nm */
1823
1824 /* because the reference flux std assumes the flux in units of ergs/cm2/s/A
1825 * we need to convert the spectrum to Angstroms and know the size of a bin
1826 * in angstroms
1827 */
1828 char band[FILE_NAME_SZ];
1829 eris_get_band(frm_sci,band);
1830
1831 double conv_factor = 0;
1832 double nm2AA = 10.;
1833 /* the following is the bin size in um */
1834 double bin_size=eris_get_dispersion(band);
1835
1836 cpl_size size = cpl_table_get_nrow(sci_tab);
1837 int status;
1838 double cdelt1 = cpl_table_get_double(sci_tab, "WAVE", size / 2, &status) -
1839 cpl_table_get_double(sci_tab, "WAVE", (size / 2 - 1), &status);
1840 cpl_table_delete(sci_tab);
1841 //cpl_msg_info(cpl_func, "cdelt1: %g", cdelt1);
1842 bin_size = cdelt1;
1843 //cpl_msg_info(cpl_func, "bin_size=%g", bin_size);
1844 conv_factor = bin_size * um2nm * nm2AA;
1845
1846 hdrl_spectrum1D_div_scalar(sci_s1d, (hdrl_value){conv_factor,0.});
1847
1848 hdrl_spectrum1D * std_s1d =
1849 eris_std_cat_frame_to_hdrl_spectrum1D(frm_std_cat, dRA, dDEC);
1850
1851 if(std_s1d == NULL) {
1852 cpl_msg_warning(cpl_func,"Ref flux STD star not found in catalogue. Skip response determination");
1853 cpl_error_set(cpl_func,CPL_ERROR_NONE);
1854 return CPL_ERROR_NONE;
1855 }
1856
1857 hdrl_spectrum1D * ext_s1d = eris_atmext_frame_to_hdrl_spectrum1D(frm_atmext);
1858
1859 //hdrl_parameter* eff_pars;
1860 /*HDRL expects telescope area in cm2 */
1861
1862 hdrl_spectrum1Dlist * telluric_models =
1863 eris_get_telluric_models(frm_tell_cat);
1864
1865 cpl_bivector * fit_areas = eris_read_wlen_windows(frm_resp_fit_areas);
1866
1867 cpl_bivector * quality_areas = eris_read_wlen_windows(frm_resp_quality_areas);
1868
1869 /*
1870 cpl_bivector * high_abs_regions = NULL;
1871 if(frm_high_abs_reg!= NULL)
1872 high_abs_regions = eris_read_wlen_windows(frm_high_abs_reg);
1873 */
1874
1875 /* set [parameters */
1876
1877 double wmin = 0;
1878 double wmax = 0;
1879 cpl_msg_info(cpl_func, "observed band: %s", band);
1880 cpl_boolean velocity_enable = CPL_TRUE;
1881 eris_get_wave_shift_param(band, &wmin, &wmax, &velocity_enable);
1882
1883
1884 cpl_msg_debug(cpl_func, "band=%s wmin=%g wmax=%g", band, wmin, wmax);
1885
1886 const hdrl_data_t lmin = (hdrl_data_t)wmin;
1887 const hdrl_data_t lmax = (hdrl_data_t)wmax;
1888
1889
1890 const cpl_boolean xcorr_normalize = (cpl_boolean) CPL_TRUE;
1891 const cpl_boolean shift_in_log_scale = (cpl_boolean) CPL_TRUE;
1892
1893 const hdrl_data_t wguess = 1094.12; // um
1894 const hdrl_data_t range_wmin = 1000;
1895 const hdrl_data_t range_wmax = 1200;
1896 const hdrl_data_t fit_wmin = 1090;
1897 const hdrl_data_t fit_wmax = 1100;
1898 const hdrl_data_t fit_half_win = 1.5;
1899
1900 hdrl_parameter* shift_par = NULL;
1901 if(velocity_enable) {
1903 range_wmin, range_wmax, fit_wmin, fit_wmax, fit_half_win);
1904 }
1905
1906 hdrl_parameter* resp_calc_par = eris_response_parameters_calc_par(plist);
1907
1908 double ra_dec_tolerance = 0.1;
1909 cpl_array * fit_points = eris_read_fit_points(frm_resp_fit_points,
1910 dRA, dDEC, ra_dec_tolerance);
1911
1912 hdrl_parameter * telluric_par = NULL;
1913 if(telluric_models != NULL && quality_areas != NULL && fit_areas != NULL && lmin < lmax) {
1914 const cpl_size xcorr_half_win = 200;
1915 hdrl_data_t xcorr_w_step = 1e-5;
1916 telluric_par = hdrl_response_telluric_evaluation_parameter_create(telluric_models,
1917 xcorr_w_step, xcorr_half_win, xcorr_normalize, shift_in_log_scale, quality_areas,
1918 fit_areas, lmin, lmax);
1919 }
1920
1921 const cpl_size post_smoothing_radius = 11;
1922
1923 //TODO: fill following params
1924 const hdrl_data_t fit_wrange = 4.0;
1925 hdrl_parameter* resp_par =
1926 hdrl_response_fit_parameter_create(post_smoothing_radius, fit_points,
1927 fit_wrange, NULL);
1928
1929 hdrl_spectrum1D * std_s1d_e = extend_hdrl_spectrum(sci_wmin * um2nm,
1930 sci_wmax * um2nm, std_s1d);
1931 hdrl_spectrum1D * ext_s1d_e = extend_hdrl_spectrum(sci_wmin * um2nm,
1932 sci_wmax * um2nm, ext_s1d);
1933
1934 cpl_msg_info(cpl_func,"compute response");
1935 hdrl_response_result * response =
1936 hdrl_response_compute(sci_s1d, std_s1d_e, ext_s1d_e, telluric_par,
1937 shift_par, resp_calc_par, resp_par);
1938 /*
1939 hdrl_spectrum1D_save(sci_s1d, "sci_s1d_resp.fits");
1940 hdrl_spectrum1D_save(std_s1d_e, "std_s1d_e_resp.fits");
1941 hdrl_spectrum1D_save(ext_s1d_e, "ext_s1d_e_resp.fits");
1942 */
1943
1944
1945 hdrl_spectrum1D_delete(&std_s1d_e);
1946 hdrl_spectrum1D_delete(&ext_s1d_e);
1947 /* extends response to same wavelength range as observed sci to flux-cal */
1948
1949
1950 if(cpl_error_get_code() == CPL_ERROR_ILLEGAL_OUTPUT) {
1951 cpl_msg_error(cpl_func,"Problems computing response. Skip this part");
1952 cpl_error_set(cpl_func,CPL_ERROR_NONE);
1953 } else if(response != NULL && cpl_error_get_code() == CPL_ERROR_NONE){
1954 cpl_msg_debug(cpl_func,"size sci = %lld size resp = %lld",
1955 hdrl_spectrum1D_get_size(sci_s1d),
1956 hdrl_spectrum1D_get_size(response->final_response));
1957 const hdrl_spectrum1D * resp_s1d_smo =
1959
1960 const hdrl_spectrum1D * resp_s1d_raw =
1962 //cpl_msg_info(cpl_func,"wmin=%g wmax=%g",sci_wmin,sci_wmax);
1963 //hdrl_spectrum1D_save(resp_s1d_raw, "resp_s1d_raw.fits");
1964 cpl_table* tab=hdrl_spectrum1D_convert_to_table(resp_s1d_smo,
1965 "response", "wavelength", "response_error", "response_bpm");
1966
1967 sci_wmin=cpl_table_get_column_min(tab, "wavelength");
1968 sci_wmax=cpl_table_get_column_max(tab, "wavelength");
1969 //cpl_msg_info(cpl_func, "wmin=%g wmax=%g", sci_wmin, sci_wmax);
1970
1971 cpl_table* resp_smo = NULL;
1972 cpl_table* resp_raw = NULL;
1973 cpl_table* resp_pts = NULL;
1974 resp_smo = hdrl_spectrum1D_convert_to_table(resp_s1d_smo,
1975 "response_smo", "wavelength", "response_smo_error",
1976 "response_smo_bpm");
1977
1978 resp_raw = hdrl_spectrum1D_convert_to_table(resp_s1d_raw,
1979 "response_raw", "wavelength", "response_raw_error",
1980 "response_raw_bpm");
1981
1984 "response_fit_pts", "wavelength", "response_fit_pts_error",
1985 "response_fit_pts_bpm");
1986 cpl_table_delete(tab);
1987
1988
1989 um2nm = 1e3;
1990 cpl_table_divide_scalar(resp_smo, "wavelength", um2nm);
1991 cpl_table_divide_scalar(resp_raw, "wavelength", um2nm);
1992 cpl_table_divide_scalar(resp_pts, "wavelength", um2nm);
1993
1994 cpl_table* qclog_tbl = eris_qclog_init();
1995 cpl_frame* frm_resp_wind_qc = cpl_frameset_find(sof, ERIS_IFU_CALIB_RESPONSE_WINDOWS);
1996
1997 const hdrl_spectrum1D * selected_points = hdrl_response_result_get_selected_response(response);
1998 const cpl_array * waves = hdrl_spectrum1D_get_wavelength(selected_points).wavelength;
1999 const double wmin_g = cpl_array_get_min(waves) / um2nm;
2000 const double wmax_g = cpl_array_get_max(waves) / um2nm;
2001 //cpl_msg_info(cpl_func, "good wavelengths wmin=%g wmax=%g", wmin_g, wmax_g);
2002
2003 /* polynomial fit */
2004 int niter = 3;
2005 int deg = 0;
2006 double kappa = 3.0;
2007 //cpl_parameterlist_dump(config,stdout);
2008 const cpl_parameter* p;
2009 char* param_name;
2010 const char* context = "eris.eris_ifu_stdstar";
2011
2012 eris_print_rec_status(100);
2013 param_name = cpl_sprintf("%s.response-polyfit-deg", context);
2014 p = cpl_parameterlist_find_const(config, param_name);
2015 deg = cpl_parameter_get_int(p);
2016 cpl_free(param_name);
2017 eris_print_rec_status(101);
2018 param_name = cpl_sprintf("%s.response-polyfit-ksigma-kappa", context);
2019 p = cpl_parameterlist_find_const(config, param_name);
2020 kappa = cpl_parameter_get_double(p);
2021 cpl_free(param_name);
2022 eris_print_rec_status(102);
2023 param_name = cpl_sprintf("%s.response-polyfit-ksigma-niter", context);
2024 p = cpl_parameterlist_find_const(config, param_name);
2025 niter = cpl_parameter_get_int(p);
2026 cpl_free(param_name);
2027 eris_print_rec_status(103);
2028
2029 cpl_table* resp_poly = eris_fit_poly_and_clip_response(resp_pts,
2030 resp_smo, deg, kappa, niter);
2031 eris_print_rec_status(104);
2032 if( resp_poly != NULL && frm_resp_wind_qc != NULL) {
2033
2034 eris_resp_qc(resp_poly, frm_resp_wind_qc, "wavelength",
2035 "response_poly", wmin_g, wmax_g, &qclog_tbl);
2036 eris_print_rec_status(105);
2037 }
2038 eris_print_rec_status(106);
2039 char* fname = cpl_sprintf("%s_response.fits",plugin_id);
2040 eris_pro_save_tbl(resp_smo, ref_set, sof, fname,
2041 ERIS_IFU_PRO_JITTER_RESPONSE, qclog_tbl, plugin_id, config);
2042
2043 cpl_table_save(resp_raw, NULL, NULL, fname, CPL_IO_EXTEND);
2044 cpl_table_save(resp_pts, NULL, NULL, fname, CPL_IO_EXTEND);
2045 cpl_table_save(resp_poly, NULL, NULL, fname, CPL_IO_EXTEND);
2046
2047 cpl_free(fname);
2048 cpl_table_delete(resp_poly);
2049 cpl_table_delete(qclog_tbl);
2050 cpl_table_delete(resp_smo);
2051 cpl_table_delete(resp_raw);
2052 cpl_table_delete(resp_pts);
2053
2054 }
2055
2056 cpl_bivector_delete(quality_areas);
2057 cpl_bivector_delete(fit_areas);
2058 //cpl_bivector_delete(high_abs_regions);
2059
2060 hdrl_parameter_delete(telluric_par);
2061 hdrl_parameter_delete(shift_par);
2062 hdrl_parameter_delete(resp_calc_par);
2063 hdrl_parameter_delete(resp_par);
2064
2065 hdrl_spectrum1D_delete(&sci_s1d);
2066 hdrl_spectrum1D_delete(&std_s1d);
2067 hdrl_spectrum1D_delete(&ext_s1d);
2068 cpl_propertylist_delete(plist);
2069 hdrl_spectrum1Dlist_delete(telluric_models);
2070
2071 eris_remove_table_normal_format();
2072 cpl_array_delete(fit_points);
2074 eris_check_error_code("eris_response_compute");
2075 if(cpl_error_get_code() == CPL_ERROR_ILLEGAL_OUTPUT) {
2076 /* occasionally, in long setting, response fail with CPL_ERROR_ILLEGAL_OUTPUT
2077 for the moment we exit.
2078 TODO: we should better handle those residual cases */
2079 cpl_error_set(cpl_func,CPL_ERROR_NONE);
2080 return CPL_ERROR_NONE;
2081 }
2082 return cpl_error_get_code();
2083}
2084
2085
2086//static void hdrl_spectrum1D_save(const hdrl_spectrum1D * s, const char * fname){
2087// if(s == NULL) return;
2088// cpl_table * tb = hdrl_spectrum1D_convert_to_table(s, "FLX", "WLN", "FLX_E",
2089// "FLX_BPM");
2090// cpl_table_save(tb, NULL, NULL, fname, CPL_IO_CREATE);
2091// cpl_table_delete(tb);
2092//}
2093
2094cpl_error_code
2095eris_flux_calibrate_spectra(const char* pipefile_prefix, const char* plugin_id,
2096 const cpl_parameterlist* config, cpl_frameset* ref_set, cpl_frameset* sof){
2097
2098
2099
2100 cpl_ensure_code(pipefile_prefix, CPL_ERROR_NULL_INPUT);
2101 cpl_ensure_code(plugin_id, CPL_ERROR_NULL_INPUT);
2102 cpl_ensure_code(config, CPL_ERROR_NULL_INPUT);
2103 cpl_ensure_code(ref_set, CPL_ERROR_NULL_INPUT);
2104 cpl_ensure_code(sof, CPL_ERROR_NULL_INPUT);
2105
2106 cpl_size ref_set_sz = cpl_frameset_get_size(ref_set);
2107 cpl_size sof_sz = cpl_frameset_get_size(sof);
2108 cpl_ensure_code(ref_set_sz >0, CPL_ERROR_DATA_NOT_FOUND);
2109 cpl_ensure_code(sof_sz >0, CPL_ERROR_DATA_NOT_FOUND);
2110
2111
2112 /* Flux calibrate the spectrum */
2113 const char* fname;
2114 cpl_msg_info(cpl_func,"Flux calibrate the spectrum");
2115 hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
2116 //cpl_frame* response=cpl_frameset_find(sof, PRO_RESPONSE);
2117 cpl_frame* spectra_frm = cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_SPECTRUM);
2118 fname = cpl_frame_get_filename(spectra_frm);
2119
2120
2121 //cpl_table* spectra_tab = cpl_table_load(fname, 1, 0);
2122 cpl_frame* frm_sci_dup = cpl_frame_duplicate(spectra_frm);
2123 cpl_table* spectra_tab = eris_ifu_table_from_sdp_to_normal_format(&frm_sci_dup);
2124 cpl_frame_delete(frm_sci_dup);
2125 char band[FILE_NAME_SZ];
2126 eris_get_band(spectra_frm, band);
2127 double bin_size = eris_get_dispersion(band);
2128 //cpl_msg_info(cpl_func, "bin_size=%g", bin_size);
2129 double um2nm = 1.e3;
2130 double nm2AA = 10;
2131 int status;
2132 double sci_wmin = 0;
2133 double sci_wmax = 0;
2134 double resp_wmin = 0;
2135 double resp_wmax = 0;
2136
2137 cpl_size size = cpl_table_get_nrow(spectra_tab);
2138 cpl_size half_size = size / 2;
2139 double cdelt1 = cpl_table_get_double(spectra_tab, "WAVE", half_size, &status) -
2140 cpl_table_get_double(spectra_tab, "WAVE", (half_size - 1), &status);
2141 //cpl_table_delete(spectra_tab);
2142 bin_size = cdelt1;
2143 double bin_size_scale_factor = bin_size * um2nm * nm2AA;
2144
2145 sci_wmin = cpl_table_get_column_min(spectra_tab,"WAVE");
2146 sci_wmax = cpl_table_get_column_max(spectra_tab,"WAVE");
2147 cpl_msg_debug(cpl_func,"Sci table wmin: %f wmax: %f",sci_wmin, sci_wmax);
2148
2149 cpl_propertylist* plist = cpl_propertylist_load(fname,0);
2150 double dit = cpl_propertylist_get_double(plist, FHDR_E_DIT);
2151
2152 //int ndit = eris_pfits_get_ndit(plist);
2153 double exptime = dit; // for NIR (SINFONI data) we use only DIT
2154 double gain = ERIS_GAIN; // eris_pfits_get_gain(plist);
2155
2156
2157 double airmass = eris_ifu_get_airmass_mean(plist);
2158 cpl_msg_info(cpl_func,"airmass: %g", airmass);
2159 cpl_propertylist_delete(plist);
2160
2161 double factor = 1. / (gain * exptime * bin_size_scale_factor);
2162
2163 /* we load response table. As response is computed in low resolution but
2164 * may be applied also to high resolution bands, we need to trim first the
2165 * response to the same range of the science spectrum, and in case the
2166 * science spectrum is defined over a wider range than the response we also
2167 * trim the observed spectrum so that both the observed spectrum and the
2168 * response have the same wavelength ranges
2169 */
2170 cpl_frame* resp_frm = cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_RESPONSE);
2171 fname = cpl_frame_get_filename(resp_frm);
2172 cpl_size next = cpl_frame_get_nextensions(resp_frm);
2173 cpl_table* resp_tab = NULL;
2174 if(next == 4) {
2175 resp_tab = cpl_table_load(fname, 4, 0);
2176 } else {
2177 resp_tab = cpl_table_load(fname, 1, 0);
2178 }
2179 resp_wmin = cpl_table_get_column_min(resp_tab,"wavelength");
2180 resp_wmax = cpl_table_get_column_max(resp_tab,"wavelength");
2181 cpl_msg_debug(cpl_func,"Response table wmin: %f wmax: %f",resp_wmin, resp_wmax);
2182 //cpl_table_save(resp_tab,NULL,NULL,"resp_tab_before.fits",CPL_IO_DEFAULT);
2183 cpl_table_and_selected_double(resp_tab, "wavelength", CPL_NOT_LESS_THAN, sci_wmin);
2184 cpl_table_and_selected_double(resp_tab, "wavelength", CPL_NOT_GREATER_THAN, sci_wmax);
2185 cpl_table* resp_xtab = cpl_table_extract_selected(resp_tab);
2186 //cpl_table_save(resp_xtab,NULL,NULL,"resp_tab_after.fits",CPL_IO_DEFAULT);
2187
2188 /* trim also observed spectrum to same range as the one of the response */
2189 //cpl_table_save(spectra_tab,NULL,NULL,"sci_tab_before.fits",CPL_IO_DEFAULT);
2190 resp_wmin = cpl_table_get_column_min(resp_xtab,"wavelength");
2191 resp_wmax = cpl_table_get_column_max(resp_xtab,"wavelength");
2192 cpl_msg_info(cpl_func,"Response table after selection wmin: %f wmax: %f",resp_wmin, resp_wmax);
2193 cpl_table_and_selected_double(spectra_tab, "WAVE", CPL_NOT_LESS_THAN, resp_wmin);
2194 cpl_table_and_selected_double(spectra_tab, "WAVE", CPL_NOT_GREATER_THAN, resp_wmax);
2195 cpl_table* spectra_xtab = cpl_table_extract_selected(spectra_tab);
2196 //cpl_table_save(spectra_xtab,NULL,NULL,"sci_tab_after.fits",CPL_IO_DEFAULT);
2197 sci_wmin = cpl_table_get_column_min(spectra_xtab,"WAVE");
2198 sci_wmax = cpl_table_get_column_max(spectra_xtab,"WAVE");
2199 cpl_msg_debug(cpl_func,"Sci table after selection wmin: %f wmax: %f",sci_wmin, sci_wmax);
2200 hdrl_spectrum1D * obs_s1d =
2201 hdrl_spectrum1D_convert_from_table(spectra_xtab, "FLUX",
2202 "WAVE", NULL, NULL, scale);
2203 cpl_table_delete(spectra_xtab);
2204 //hdrl_spectrum1D_save(obs_s1d, "obs_s1d_spct.fits");
2205
2206
2207 //cpl_msg_info(cpl_func,"Response table wmin: %g",cpl_table_get_column_min(resp_tab,"wavelength"));
2208 //cpl_msg_info(cpl_func,"Response table wmax: %g",cpl_table_get_column_max(resp_tab,"wavelength"));
2209 /* As response curve may have a slightly different wavelength span than the observed
2210 spectrum, we extract the relevant common region from the response table */
2211
2212 hdrl_spectrum1D * resp_s1d = NULL;
2213 if(next == 4) {
2214 resp_s1d = hdrl_spectrum1D_convert_from_table(resp_xtab, "response_poly",
2215 "wavelength", "response_poly_error", "response_poly_bpm", scale);
2216 } else {
2217 resp_s1d = hdrl_spectrum1D_convert_from_table(resp_xtab, "response_smo",
2218 "wavelength", "response_smo_error", "response_smo_bpm", scale);
2219 }
2220 //hdrl_spectrum1D_save(resp_s1d, "resp_s1d_spct.fits");
2221 hdrl_spectrum1D * resp_s1d_e = extend_hdrl_spectrum(sci_wmin, sci_wmax, resp_s1d);
2222 hdrl_spectrum1D * obs_s1d_e = extend_hdrl_spectrum(resp_wmin, resp_wmax, obs_s1d);
2223 //hdrl_spectrum1D_save(obs_s1d_e, "obs_s1d_e_spct.fits");
2224
2225 hdrl_parameter * res_par =
2226 hdrl_spectrum1D_resample_interpolate_parameter_create(hdrl_spectrum1D_interp_akima);
2227
2228 hdrl_spectrum1D_wavelength wav = hdrl_spectrum1D_get_wavelength(obs_s1d_e);
2229
2230 hdrl_spectrum1D * resp_s1d_r = hdrl_spectrum1D_resample(resp_s1d_e, &wav, res_par);
2231 //hdrl_spectrum1D_save(resp_s1d_r, "resp_s1d_r_spct.fits");
2232 cpl_msg_debug(cpl_func,"obs_e wmin: %f wmax: %f", cpl_array_get_min(hdrl_spectrum1D_get_wavelength(obs_s1d_e).wavelength), cpl_array_get_max(hdrl_spectrum1D_get_wavelength(obs_s1d_e).wavelength));
2233 cpl_msg_debug(cpl_func,"rsp_r wmin: %f wmax: %f", cpl_array_get_min(hdrl_spectrum1D_get_wavelength(resp_s1d_r).wavelength), cpl_array_get_max(hdrl_spectrum1D_get_wavelength(resp_s1d_r).wavelength));
2234
2235 /* now response and observed spectrum have same range and sampling step
2236 * we can multiply one by the other
2237 */
2238
2239
2240 cpl_frame* atm_ext_frm = cpl_frameset_find(sof, ERIS_IFU_CALIB_EXTCOEFF_TABLE);
2241 fname = cpl_frame_get_filename(atm_ext_frm);
2242 cpl_table* atm_ext_tab = cpl_table_load(fname,1,0);
2243
2244 hdrl_spectrum1D * ext_s1d =
2245 hdrl_spectrum1D_convert_from_table(atm_ext_tab, "EXTINCTION",
2246 "LAMBDA", NULL, NULL, scale);
2247 //hdrl_spectrum1D_save(ext_s1d, "ext_s1d_spct.fits");
2248 hdrl_spectrum1D * ext_s1d_e = extend_hdrl_spectrum(sci_wmin * um2nm, sci_wmax * um2nm, ext_s1d);
2249 //hdrl_spectrum1D_save(ext_s1d_e, "ext_s1d_e.fits");
2250
2251
2252 //hdrl_spectrum1D_save(obs_s1d_e, "obs_s1d_e.fits");
2253 //hdrl_spectrum1D_save(resp_s1d_e, "resp_s1d_e_spct.fits");
2254 hdrl_spectrum1D* flux_s1d = hdrl_spectrum1D_mul_spectrum_create(obs_s1d_e, resp_s1d_r);
2255
2256 //hdrl_spectrum1D_save(flux_s1d, "flux_s1d_0_spct.fits");
2257 hdrl_spectrum1D_mul_scalar(flux_s1d,(hdrl_value){factor,0.});
2258
2259 //hdrl_spectrum1D_save(flux_s1d, "flux_s1d_01_spct.fits");
2260
2261 double* pflux = cpl_image_get_data_double(hdrl_image_get_image(flux_s1d->flux));
2262 int sx = hdrl_image_get_size_x(flux_s1d->flux);
2263
2264
2265 double* pext = cpl_image_get_data_double(hdrl_image_get_image(ext_s1d->flux));
2266 cpl_msg_warning(cpl_func,"sx: %d",sx);
2267 for(int i = 0; i < sx; i++) {
2268 double exp = -0.4 * airmass * pext[i];
2269 pflux[i] *= pow(10., exp);
2270 }
2271 //hdrl_spectrum1D_save(flux_s1d, "flux_s1d_1_spct.fits");
2272
2273 cpl_table* flux_tab = hdrl_spectrum1D_convert_to_table(flux_s1d,
2274 "flux", "wavelength", "flux_error", "flux_bpm");
2275 //cpl_table_save(flux_tab,NULL,NULL,"spectrum_flux_spct.fits",CPL_IO_DEFAULT);
2276
2277 cpl_table* qclog_tbl = eris_qclog_init();
2278 char* file_name = cpl_sprintf("%s_spectrum_fluxcal",pipefile_prefix);
2279 eris_pro_save_tbl(flux_tab, ref_set, sof, file_name,
2280 ERIS_IFU_PRO_JITTER_SPECTRUM_FLUXCAL, qclog_tbl, plugin_id, config);
2281 cpl_free(file_name);
2282 cpl_table_delete(qclog_tbl);
2283 eris_remove_table_normal_format();
2284
2285 hdrl_spectrum1D_delete(&resp_s1d_r);
2286 hdrl_spectrum1D_delete(&obs_s1d_e);
2287 hdrl_spectrum1D_delete(&obs_s1d);
2288 hdrl_spectrum1D_delete(&resp_s1d);
2289 hdrl_spectrum1D_delete(&resp_s1d_e);
2290 hdrl_spectrum1D_delete(&ext_s1d);
2291 hdrl_spectrum1D_delete(&ext_s1d_e);
2292 hdrl_spectrum1D_delete(&flux_s1d);
2293 cpl_table_delete(atm_ext_tab);
2294 cpl_table_delete(flux_tab);
2295 cpl_table_delete(spectra_tab);
2296 cpl_table_delete(resp_tab);
2297 cpl_table_delete(resp_xtab);
2298 hdrl_parameter_delete(res_par);
2299 eris_check_error_code("eris_flux_calibrate_spectra");
2300
2301 return cpl_error_get_code();
2302}
2303
2304//static cpl_error_code
2305//eris_ifu_apply_response_to_cube(hdrl_imagelist* hlist,
2306// hdrl_spectrum1D* resp_s1d_r, hdrl_spectrum1D* resp_s1d_e,
2307// hdrl_spectrum1D * ext_s1d, const double airm, const double factor)
2308//{
2309
2310
2311// //hdrl_spectrum1D_save(ext_s1d, "ext_s1d_cube.fits");
2312// double* presp = cpl_image_get_data_double(hdrl_image_get_image(resp_s1d_r->flux));
2313// double* pext = cpl_image_get_data_double(hdrl_image_get_image(ext_s1d->flux));
2314// cpl_size sx = hdrl_imagelist_get_size_x(hlist);
2315// cpl_size sy = hdrl_imagelist_get_size_y(hlist);
2316// cpl_size sz = hdrl_imagelist_get_size(hlist);
2317
2318// int nx=hdrl_image_get_size_x(resp_s1d_e->flux);
2319// cpl_msg_warning(cpl_func,"nx: %d sz: %lld",nx, sz);
2320// cpl_size kmax = (nx < sz) ? nx : sz;
2321// kmax = sz;
2322// double resp=0;
2323// hdrl_value resp_fct;
2324// hdrl_value factor_value = {factor, 0};
2325
2326// cpl_msg_warning(cpl_func,"kmax: %lld",kmax);
2327// for (cpl_size k = 0; k < kmax; k++) {
2328// double exp = -0.4 * airm * pext[k];
2329
2330// hdrl_image* hima = hdrl_imagelist_get(hlist, k);
2331// hdrl_value atm_ext = {pow(10., exp), 0.};
2332// resp = presp[k];
2333// resp_fct.data = resp;
2334// //if(isfinite(resp) && resp != 0) {
2335// if(isfinite(resp) && resp != 0) {
2336// hdrl_image_mul_scalar(hima, resp_fct);
2337// } else if(resp == 0) {
2338// hdrl_image_mul_scalar(hima, resp_fct);
2339// /* Flag bad pixel in result */
2340// //cpl_msg_warning(cpl_func,"flag bad plane: %lld",k);
2341// for(cpl_size j = 1; j <= sy; j++ ) {
2342// for(cpl_size i = 1; i <= sx; i++ ) {
2343// hdrl_image_reject(hima, i, j);
2344// }
2345// }
2346// } else if(!isfinite(resp)) {
2347// // Flag bad pixel in result
2348// //cpl_msg_warning(cpl_func,"flag bad plane: %lld",k);
2349// for(cpl_size j = 1; j <= sy; j++ ) {
2350// for(cpl_size i = 1; i <= sx; i++ ) {
2351// hdrl_image_reject(hima, i, j);
2352// }
2353// }
2354
2355// }
2356
2357// hdrl_image_mul_scalar(hima, atm_ext);
2358// hdrl_image_mul_scalar(hima, factor_value);
2359// }
2360
2361// eris_check_error_code("eris_ifu_apply_response_to_cube");
2362// return cpl_error_get_code();
2363//}
2364
2365//TODO Here the wmin/wmax should be taken from the cube and not from the extracted spectrum.
2366cpl_error_code
2367eris_flux_calibrate_cube2(const char* procatg, const char* pipefile_prefix,
2368 const char* plugin_id, const cpl_parameterlist* parlist,
2369 cpl_frameset* sof){
2370
2371 cpl_ensure_code(procatg, CPL_ERROR_NULL_INPUT);
2372 cpl_ensure_code(pipefile_prefix, CPL_ERROR_NULL_INPUT);
2373 cpl_ensure_code(plugin_id, CPL_ERROR_NULL_INPUT);
2374
2375 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
2376 cpl_ensure_code(sof, CPL_ERROR_NULL_INPUT);
2377
2378 cpl_size sof_sz = cpl_frameset_get_size(sof);
2379 cpl_ensure_code(sof_sz >0, CPL_ERROR_DATA_NOT_FOUND);
2380
2381
2382
2383
2384 cpl_msg_info(cpl_func,"Flux calibrate the cube %s", procatg);
2385
2386 const char* fname;
2387 hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
2388 /* load input science cube frame */
2389 cpl_frame* coadd_frm = cpl_frameset_find(sof, procatg);
2390 fname = cpl_frame_get_filename(coadd_frm);
2391 cpl_msg_debug(cpl_func, "fname=%s",fname);
2392 cpl_imagelist* cube = cpl_imagelist_load(fname, CPL_TYPE_FLOAT,1);
2393 cpl_imagelist* errs = cpl_imagelist_load(fname, CPL_TYPE_FLOAT,2);
2394 cpl_propertylist* phead = cpl_propertylist_load(fname, 0);
2395 hdrl_imagelist* hlist = hdrl_imagelist_create(cube, errs);
2396 cpl_imagelist_delete(cube);
2397 cpl_imagelist_delete(errs);
2398 double airm = eris_ifu_get_airmass_mean(phead);
2399
2400 /* define bin scale factor */
2401 char band[FILE_NAME_SZ];
2402 eris_get_band(coadd_frm, band);
2403 double bin_size = eris_get_dispersion(band);
2404 cpl_msg_debug(cpl_func, "airm=%g", airm);
2405 double um2nm = 1.e3;
2406 double nm2AA = 10;
2407
2408 double sci_wmin = 0;
2409 double sci_wmax = 0;
2410 double resp_wmin = 0;
2411 double resp_wmax = 0;
2412 cpl_propertylist* data_head = cpl_propertylist_load(fname, 1);
2413 cpl_propertylist* head_wcs = eris_ifu_plist_extract_wcs(data_head);
2414 cpl_propertylist_delete(data_head);
2415 cpl_size size = hdrl_imagelist_get_size(hlist);
2416 double cdelt3 = cpl_propertylist_get_double(head_wcs, "CD3_3");
2417 double crval3 = cpl_propertylist_get_double(head_wcs, "CRVAL3");
2418 bin_size = cdelt3;
2419 double bin_size_scale_factor = bin_size * um2nm * nm2AA;
2420 double dit = cpl_propertylist_get_double(phead, FHDR_E_DIT);
2421 int ndit = cpl_propertylist_get_int(phead, "ESO DET NDIT");
2422 double exptime = dit * ndit;
2423 cpl_msg_debug(cpl_func, "exptime=%g", exptime); //TODO: DIT OR DIT * NDIT??
2424 double gain = ERIS_GAIN; //SINFONI_GAIN; // eris_pfits_get_gain(plist);
2425 cpl_msg_debug(cpl_func, "gain=%g", gain);
2426
2427 //cpl_propertylist_delete(phead);
2428
2429 double factor = 1. / (gain * exptime * bin_size_scale_factor);
2430 sci_wmin = crval3;
2431 sci_wmax = crval3+cdelt3*(size-1);
2432 cpl_msg_warning(cpl_func,"Sci table wmin: %f wmax: %f",sci_wmin, sci_wmax);
2433
2434
2435 /* we load response table. As response is computed in low resolution but
2436 * may be applied also to high resolution bands, we need to trim first the
2437 * response to the same range of the science spectrum, and in case the
2438 * science spectrum is defined over a wider range than the response we also
2439 * trim the observed spectrum so that both the observed spectrum and the
2440 * response have the same wavelength ranges
2441 */
2442 cpl_frame* resp_frm = cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_RESPONSE);
2443 cpl_size next = cpl_frame_get_nextensions(resp_frm);
2444 fname = cpl_frame_get_filename(resp_frm);
2445 cpl_table* resp_tab = NULL;
2446 if(next == 4) {
2447 resp_tab = cpl_table_load(fname, 4, 0);
2448 } else {
2449 resp_tab = cpl_table_load(fname, 1, 0);
2450 }
2451
2452 resp_wmin = cpl_table_get_column_min(resp_tab,"wavelength");
2453 resp_wmax = cpl_table_get_column_max(resp_tab,"wavelength");
2454 cpl_msg_warning(cpl_func,"Response table wmin: %f wmax: %f",resp_wmin, resp_wmax);
2455 //cpl_table_save(resp_tab,NULL,NULL,"resp_tab_before.fits",CPL_IO_DEFAULT);
2456 cpl_table_and_selected_double(resp_tab, "wavelength", CPL_NOT_LESS_THAN, sci_wmin);
2457 cpl_table_and_selected_double(resp_tab, "wavelength", CPL_NOT_GREATER_THAN, sci_wmax);
2458 cpl_table* resp_xtab = cpl_table_extract_selected(resp_tab);
2459 //cpl_table_save(resp_xtab,NULL,NULL,"resp_tab_after.fits",CPL_IO_DEFAULT);
2460
2461 /* trim also observed spectrum to same range as the one of the response */
2462 //cpl_table_save(spectra_tab,NULL,NULL,"sci_tab_before.fits",CPL_IO_DEFAULT);
2463 resp_wmin = cpl_table_get_column_min(resp_xtab,"wavelength");
2464 resp_wmax = cpl_table_get_column_max(resp_xtab,"wavelength");
2465 cpl_msg_info(cpl_func,"Response table after selection wmin: %f wmax: %f",resp_wmin, resp_wmax);
2466
2467 /* FUNCTION TO TRIM CUBE OVER SAME RANGE AS ONE OF RESPONSE */
2468 cpl_size klow_skip = 0;
2469 klow_skip = (resp_wmin > sci_wmin) ? (resp_wmin - sci_wmin) / cdelt3 : 0;
2470 cpl_size kupp_skip = 0;
2471 kupp_skip = (sci_wmax > resp_wmax) ? (sci_wmax - resp_wmax) / cdelt3 : 0;
2472 cpl_msg_info(cpl_func,"klow_skip: %lld kupp_skip: %lld",klow_skip, kupp_skip);
2473 cpl_size z_min = klow_skip;
2474 cpl_size z_max = size - kupp_skip;
2475 cpl_msg_info(cpl_func,"z_min: %lld: z_max: %lld",z_min, z_max);
2476 cpl_imagelist* bpm = cpl_imagelist_new();
2477
2478 for(cpl_size k=0; k< size; k++) {
2479 //cpl_msg_info(cpl_func,"k=%lld",k);
2480 cpl_imagelist_set(bpm, hdrl_image_get_image(hdrl_imagelist_get(hlist,k)), k);
2481 }
2482
2483 //cpl_propertylist* xhead = cpl_propertylist_load(fname, 1);
2484 cpl_propertylist_append_int(head_wcs,"NAXIS3",size);
2485 eris_ifu_cube_trim_nans(z_min, z_max,&hlist, &bpm, head_wcs);
2486 //cpl_imagelist_delete(bpm);
2487 /* updated sci_wmin, sci_wmax */
2488 sci_wmin = crval3;
2489 sci_wmax = crval3+cdelt3*(size-1);
2490 cpl_msg_warning(cpl_func,"Sci table wmin: %f wmax: %f",sci_wmin, sci_wmax);
2491
2492 hdrl_spectrum1D * resp_s1d = NULL;
2493 if(next == 4) {
2494 resp_s1d = hdrl_spectrum1D_convert_from_table(resp_xtab, "response_poly",
2495 "wavelength", "response_poly_error", "response_poly_bpm", scale);
2496 } else {
2497 resp_s1d = hdrl_spectrum1D_convert_from_table(resp_xtab, "response_smo",
2498 "wavelength", "response_smo_error", "response_smo_bpm", scale);
2499 }
2500 cpl_table_delete(resp_xtab);
2501 //hdrl_spectrum1D_save(resp_s1d, "resp_s1d_cube.fits");
2502 hdrl_spectrum1D * resp_s1d_e = extend_hdrl_spectrum(sci_wmin, sci_wmax, resp_s1d);
2503 //hdrl_spectrum1D_save(resp_s1d_e, "resp_s1d_e_cube.fits");
2504
2505 /* creates a wavelength sampling array to re-sample the response properly
2506 * on the same range and with same sampling as the science data cube */
2507 hdrl_parameter * res_par =
2508 hdrl_spectrum1D_resample_interpolate_parameter_create(hdrl_spectrum1D_interp_akima);
2509
2510 cpl_frame* spectra_frm = NULL;
2511
2512 if( NULL == (spectra_frm = cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_SPECTRUM))){
2513 cpl_msg_error(cpl_func,"To flux calibrate a data cube you need to set extract-source=TRUE");
2514 cpl_msg_error(cpl_func,"as the recipe uses the min/max wavelength of the extracted spectrum");
2515 cpl_msg_error(cpl_func,"to determine the range where to flux calibrate the data cube.");
2516 eris_check_error_code("eris_flux_calibrate_cube2");
2517
2518 hdrl_spectrum1D_delete(&resp_s1d);
2519 hdrl_spectrum1D_delete(&resp_s1d_e);
2520 cpl_table_delete(resp_tab);
2521 hdrl_parameter_delete(res_par);
2522 hdrl_imagelist_delete(hlist);
2523 cpl_propertylist_delete(head_wcs);
2524 //cpl_imagelist_delete(bpm);
2525 cpl_propertylist_delete(phead);
2526
2527 return cpl_error_get_code();
2528 }
2529
2530 fname = cpl_frame_get_filename(spectra_frm);
2531
2532 //cpl_table* spectra_tab = cpl_table_load(fname, 1, 0);
2533 cpl_frame* frm_sci_dup = cpl_frame_duplicate(spectra_frm);
2534 cpl_table* spectra_tab = eris_ifu_table_from_sdp_to_normal_format(&frm_sci_dup);
2535 cpl_frame_delete(frm_sci_dup);
2536
2537 cpl_table_and_selected_double(spectra_tab, "WAVE", CPL_NOT_LESS_THAN, sci_wmin);
2538 cpl_table_and_selected_double(spectra_tab, "WAVE", CPL_NOT_GREATER_THAN, sci_wmax);
2539 cpl_table* spectra_xtab = cpl_table_extract_selected(spectra_tab);
2540 cpl_table_delete(spectra_tab);
2541 //cpl_table_save(spectra_xtab,NULL,NULL,"sci_tab_after.fits",CPL_IO_DEFAULT);
2542 sci_wmin = cpl_table_get_column_min(spectra_xtab,"WAVE");
2543 sci_wmax = cpl_table_get_column_max(spectra_xtab,"WAVE");
2544 cpl_msg_debug(cpl_func,"Sci table after selection wmin: %f wmax: %f",sci_wmin, sci_wmax);
2545 hdrl_spectrum1D * obs_s1d = hdrl_spectrum1D_convert_from_table(spectra_xtab, "FLUX",
2546 "WAVE", NULL, NULL, scale);
2547 cpl_table_delete(spectra_xtab);
2548
2549 hdrl_spectrum1D_wavelength wav = hdrl_spectrum1D_get_wavelength(obs_s1d);
2550 hdrl_spectrum1D * resp_s1d_r = hdrl_spectrum1D_resample(resp_s1d_e, &wav, res_par);
2551 //hdrl_spectrum1D_save(resp_s1d_r, "resp_s1d_r_cube.fits");
2552 cpl_msg_debug(cpl_func,"obs_e wmin: %f wmax: %f", sci_wmin, sci_wmax);
2553 //cpl_msg_debug(cpl_func,"rsp_r wmin: %f wmax: %f", cpl_array_get_min(hdrl_spectrum1D_get_wavelength(resp_s1d_r).wavelength), cpl_array_get_max(hdrl_spectrum1D_get_wavelength(resp_s1d_r).wavelength));
2554
2555 /* now response and observed spectrum have same range and sampling step
2556 * we can multiply one by the other
2557 */
2558
2559
2560 cpl_frame* atm_ext_frm = cpl_frameset_find(sof, ERIS_IFU_CALIB_EXTCOEFF_TABLE);
2561 fname = cpl_frame_get_filename(atm_ext_frm);
2562 cpl_table* atm_ext_tab = cpl_table_load(fname,1,0);
2563
2564 hdrl_spectrum1D * ext_s1d =
2565 hdrl_spectrum1D_convert_from_table(atm_ext_tab, "EXTINCTION",
2566 "LAMBDA", NULL, NULL, scale);
2567 //hdrl_spectrum1D_save(ext_s1d, "ext_s1d_cube.fits");
2568 /*
2569 eris_ifu_apply_response_to_cube(hlist, resp_s1d_r, resp_s1d_e, ext_s1d,
2570 airm, factor);
2571 */
2572
2573
2574
2575 double* presp = cpl_image_get_data_double(hdrl_image_get_image(resp_s1d_r->flux));
2576 double* pext = cpl_image_get_data_double(hdrl_image_get_image(ext_s1d->flux));
2577 cpl_size sx = hdrl_imagelist_get_size_x(hlist);
2578 cpl_size sy = hdrl_imagelist_get_size_y(hlist);
2579 cpl_size sz = hdrl_imagelist_get_size(hlist);
2580
2581 int nx=hdrl_image_get_size_x(resp_s1d_e->flux);
2582 cpl_msg_warning(cpl_func,"nx: %d sz: %lld",nx, sz);
2583 cpl_size kmax = (nx < sz) ? nx : sz;
2584 kmax = sz;
2585 double resp=0;
2586 hdrl_value resp_fct;
2587 hdrl_value factor_value = {factor, 0};
2588
2589 cpl_msg_warning(cpl_func,"kmax: %lld",kmax);
2590 for (cpl_size k = 0; k < kmax; k++) {
2591 double exp = -0.4 * airm * pext[k];
2592
2593 hdrl_image* hima = hdrl_imagelist_get(hlist, k);
2594 hdrl_value atm_ext = {pow(10., exp), 0.};
2595 resp = presp[k];
2596 resp_fct.data = resp;
2597 //if(isfinite(resp) && resp != 0) {
2598 if(isfinite(resp) && resp != 0) {
2599 hdrl_image_mul_scalar(hima, resp_fct);
2600 } else if(resp == 0) {
2601 hdrl_image_mul_scalar(hima, resp_fct);
2602 // Flag bad pixel in result
2603 //cpl_msg_warning(cpl_func,"flag bad plane: %lld",k);
2604 for(cpl_size j = 1; j <= sy; j++ ) {
2605 for(cpl_size i = 1; i <= sx; i++ ) {
2606 hdrl_image_reject(hima, i, j);
2607 }
2608 }
2609 } else if(!isfinite(resp)) {
2610 // Flag bad pixel in result
2611 //cpl_msg_warning(cpl_func,"flag bad plane: %lld",k);
2612 for(cpl_size j = 1; j <= sy; j++ ) {
2613 for(cpl_size i = 1; i <= sx; i++ ) {
2614 hdrl_image_reject(hima, i, j);
2615 }
2616 }
2617
2618 }
2619
2620 hdrl_image_mul_scalar(hima, atm_ext);
2621 hdrl_image_mul_scalar(hima, factor_value);
2622 }
2623
2624
2625
2626 /* save results */
2627 char* fcal_cube_ptag = cpl_sprintf("%s_FLUXCAL", procatg);
2628 char* cube_file_name = cpl_sprintf("%s_cube_fluxcal.fits", pipefile_prefix);
2629 //cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
2630 //if(cpl_error_get_code() != CPL_ERROR_NONE) exit(0);
2631 cpl_propertylist_set_string(phead, FHDR_PRO_CATG, fcal_cube_ptag);
2632
2633 hdrl_resample_result* res = cpl_calloc(1, sizeof(hdrl_resample_result));
2634 //res->himlist = hdrl_imagelist_duplicate(hlist);
2635 res->himlist = hlist;
2636 res->header = cpl_propertylist_duplicate(phead);
2637 cpl_propertylist_delete(phead);
2638 cpl_propertylist_append(res->header, head_wcs);
2639 cpl_propertylist_delete(head_wcs);
2640 cpl_propertylist_update_string(res->header, "BUNIT", UNIT_FLUXCAL);
2641
2642
2643 /* collapse data cube */
2644
2645 cpl_propertylist* phead2D = cpl_propertylist_duplicate(res->header);
2646 eris_ifu_plist_erase_wcs3D(phead2D);
2647 char* mean_cube_pctg = cpl_sprintf("%s_%s",fcal_cube_ptag,"MEAN");
2648 const char *filenameSpec ="cube_fluxcal";
2649 char* ima_file_name = cpl_sprintf("%s_%s_mean.fits", pipefile_prefix, filenameSpec);
2650 cpl_propertylist_append_string(phead2D, "PRODCATG", PRODCATG_WHITELIGHT);
2651 cpl_propertylist_update_string(phead2D, "BUNIT", UNIT_FLUXCAL);
2652 hdrl_image* cube_collapsed = NULL;
2653 cpl_image* cube_cmap = NULL;
2654 //TODO: AMO checking why collapsed mean of flux calibrated result is full of NANs.
2655 //hdrl_value fct = {1000000.,0};
2656 //hdrl_imagelist_mul_scalar(res->himlist, fct);
2657 hdrl_imagelist_collapse_mean(res->himlist, &cube_collapsed, &cube_cmap);
2658 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, mean_cube_pctg);
2659 cpl_image* mask_ima = cpl_image_new_from_mask(hdrl_image_get_mask(cube_collapsed));
2660 eris_ifu_save_deq_image(sof, NULL, parlist,sof, NULL,
2661 plugin_id, phead2D, "RADECSYS", ima_file_name,
2662 hdrl_image_get_image(cube_collapsed),
2663 hdrl_image_get_error(cube_collapsed),
2664 rmse, mask_ima, flag16bit, UNIT_FLUXCAL);
2665 cpl_image_delete(mask_ima);
2666 cpl_free(mean_cube_pctg);
2667 cpl_free(ima_file_name);
2668 cpl_propertylist_delete(phead2D);
2669
2670 if(strcmp(plugin_id,"eris_ifu_jitter") == 0) {
2671 char* param_name = cpl_sprintf("%s.crea-phase3", "eris.eris_ifu_jitter");
2672 const cpl_parameter* p = cpl_parameterlist_find_const(parlist, param_name);
2673
2674 cpl_boolean gen_phase3 = cpl_parameter_get_bool(p);
2675 if(gen_phase3) {
2676 cpl_msg_info(cpl_func,"Save flux calibrated COADDED OBJECT cube with phase3 data format");
2677 eris_ifu_resample_save_cube(res, fcal_cube_ptag, plugin_id,
2678 cube_file_name, parlist, sof, CPL_TRUE);
2679 } else {
2680 cpl_msg_info(cpl_func,"Save flux calibrated COADDED OBJECT cube");
2681 eris_ifu_resample_save_cube(res, fcal_cube_ptag, plugin_id,
2682 cube_file_name, parlist, sof, CPL_FALSE);
2683 }
2684 cpl_free(param_name);
2685 } else {
2686 cpl_msg_info(cpl_func,"Save flux calibrated STD STAR cube");
2687 eris_ifu_resample_save_cube(res, fcal_cube_ptag, plugin_id,
2688 cube_file_name, parlist, sof, CPL_FALSE);
2689 }
2690
2691 cpl_free(cube_file_name);
2692
2693 cpl_table* qclog = eris_qclog_init();
2694 hdrl_image_delete(cube_collapsed);
2695 cpl_image_delete(cube_cmap);
2696
2697 cpl_table_delete(qclog);
2698
2699
2700 cpl_free(fcal_cube_ptag);
2701 hdrl_spectrum1D_delete(&resp_s1d);
2702 hdrl_spectrum1D_delete(&ext_s1d);
2703 hdrl_spectrum1D_delete(&resp_s1d_e);
2704 hdrl_spectrum1D_delete(&obs_s1d);
2705 hdrl_spectrum1D_delete(&resp_s1d_r);
2706 cpl_table_delete(atm_ext_tab);
2707 cpl_table_delete(resp_tab);
2708 eris_remove_table_normal_format();
2709 hdrl_parameter_delete(res_par);
2710 hdrl_imagelist_delete(hlist);
2711 //hdrl_resample_result_delete(res);
2712 cpl_free(res);
2713 eris_check_error_code("eris_flux_calibrate_cube2");
2714
2715 return cpl_error_get_code();
2716}
2717
2718//TODO: missing error and bad pixel propagation
2719cpl_error_code
2720eris_flux_calibrate_cube(const char* procatg, const char* pipefile_prefix,
2721 const char* plugin_id, const cpl_parameterlist* parlist,
2722 /*cpl_frameset* ref_set, */cpl_frameset* sof){
2723
2724
2725 cpl_ensure_code(procatg, CPL_ERROR_NULL_INPUT);
2726 cpl_ensure_code(pipefile_prefix, CPL_ERROR_NULL_INPUT);
2727 cpl_ensure_code(plugin_id, CPL_ERROR_NULL_INPUT);
2728
2729 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
2730 cpl_ensure_code(sof, CPL_ERROR_NULL_INPUT);
2731
2732 cpl_size sof_sz = cpl_frameset_get_size(sof);
2733 cpl_ensure_code(sof_sz >0, CPL_ERROR_DATA_NOT_FOUND);
2734
2735 cpl_msg_info(cpl_func,"Flux calibrate the cube %s", procatg);
2736 // Flux calibrate the spectrum
2737 const char* fname;
2738 hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
2739
2740 cpl_frame* coadd_frm = cpl_frameset_find(sof, procatg);
2741 fname = cpl_frame_get_filename(coadd_frm);
2742 cpl_msg_debug(cpl_func, "fname=%s",fname);
2743 cpl_imagelist* cube = cpl_imagelist_load(fname, CPL_TYPE_FLOAT,1);
2744 cpl_imagelist* errs = cpl_imagelist_load(fname, CPL_TYPE_FLOAT,2);
2745 cpl_propertylist* phead = cpl_propertylist_load(fname, 0);
2746 hdrl_imagelist* hlist = hdrl_imagelist_create(cube, errs);
2747 cpl_imagelist_delete(cube);
2748 cpl_imagelist_delete(errs);
2749 double dit = cpl_propertylist_get_double(phead, FHDR_E_DIT);
2750 int ndit = cpl_propertylist_get_int(phead, "ESO DET NDIT");
2751
2752 double centralLambda = 0;
2753 double dis = 0;
2754 double centralpix = 0;
2755 double newcenter_x = 0;
2756 double newcenter_y = 0;
2757
2758 cpl_propertylist* data_head = cpl_propertylist_load(fname, 1);
2759
2760 eris_get_wcs_cube(data_head, &centralLambda, &dis, &centralpix, &newcenter_x,
2761 &newcenter_y);
2762 cpl_propertylist* head_wcs = eris_ifu_plist_extract_wcs(data_head);
2763 cpl_propertylist_delete(data_head);
2764
2765 double exptime = dit * ndit;
2766 cpl_msg_debug(cpl_func, "exptime=%g", exptime); //TODO: DIT OR DIT * NDIT??
2767 double gain = ERIS_GAIN; //SINFONI_GAIN; // eris_pfits_get_gain(plist);
2768 cpl_msg_debug(cpl_func, "gain=%g", gain);
2769
2770 double airm = eris_ifu_get_airmass_mean(phead);
2771 cpl_msg_debug(cpl_func, "airm=%g", airm);
2772 double um2nm = 1000.;
2773 double nm2AA = 10.;
2774 char band[FILE_NAME_SZ];
2775
2776 eris_get_band(coadd_frm,band);
2777
2778 double bin_size = eris_get_dispersion(band);
2779
2780 double bin_size_scale_factor = bin_size * um2nm * nm2AA;
2781 hdrl_value factor = {1. / (gain * exptime * bin_size_scale_factor), 0};
2782 cpl_msg_debug(cpl_func, "factor=%g", factor.data);
2783 cpl_frame* resp_frm = cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_RESPONSE);
2784 fname = cpl_frame_get_filename(resp_frm);
2785 cpl_table* resp_tab = cpl_table_load(fname,1,0);
2786
2787 hdrl_spectrum1D * resp_s1d =
2788 hdrl_spectrum1D_convert_from_table(resp_tab, "response_smo",
2789 "wavelength", "response_smo_error", "response_smo_bpm", scale);
2790 //hdrl_spectrum1D_save(resp_s1d, "resp_s1d_cube.fits");
2791
2792 cpl_frame* atm_ext_frm = cpl_frameset_find(sof, ERIS_IFU_CALIB_EXTCOEFF_TABLE);
2793 fname = cpl_frame_get_filename(atm_ext_frm);
2794 cpl_table* atm_ext_tab = cpl_table_load(fname, 1, 0);
2795
2796 hdrl_spectrum1D * ext_s1d =
2797 hdrl_spectrum1D_convert_from_table(atm_ext_tab, "EXTINCTION",
2798 "LAMBDA", NULL, NULL, scale);
2799 //hdrl_spectrum1D_save(ext_s1d, "ext_s1d_cube.fits");
2800
2801 hdrl_parameter * res_par =
2802 hdrl_spectrum1D_resample_interpolate_parameter_create(hdrl_spectrum1D_interp_akima);
2803
2804 hdrl_spectrum1D_wavelength wav = hdrl_spectrum1D_get_wavelength(resp_s1d);
2805
2806 hdrl_spectrum1D * ext_s1d_e = hdrl_spectrum1D_resample(ext_s1d, &wav, res_par);
2807
2808 hdrl_parameter_delete(res_par);
2809 //BAD
2810
2811 //hdrl_spectrum1D_save(ext_s1d_e, "resp_s1d_e_cube.fits");
2812
2813 int nx=hdrl_image_get_size_x(resp_s1d->flux);
2814 double* pext = cpl_image_get_data_double(hdrl_image_get_image(ext_s1d->flux));
2815
2816 // Flux calibrate the cube
2817 double* presp = cpl_image_get_data_double(hdrl_image_get_image(resp_s1d->flux));
2818
2819
2820
2821 /* TODO:
2822 * here we should also check that hlist and the other arrays are aligned
2823 * in wavelength
2824 */
2825
2826 cpl_size sx = hdrl_imagelist_get_size_x(hlist);
2827 cpl_size sy = hdrl_imagelist_get_size_y(hlist);
2828 cpl_size sz = hdrl_imagelist_get_size(hlist);
2829
2830 cpl_size kmax = (nx < sz) ? nx : sz;
2831 double resp=0;
2832 hdrl_value resp_fct;
2833
2834 for (cpl_size k = 0; k < kmax; k++) {
2835 double exp = -0.4 * airm * pext[k];
2836
2837 hdrl_image* hima = hdrl_imagelist_get(hlist, k);
2838 hdrl_value atm_ext = {pow(10., exp), 0.};
2839 resp = presp[k];
2840 resp_fct.data = resp;
2841 if(isfinite(resp)) {
2842 hdrl_image_mul_scalar(hima, resp_fct);
2843 } else {
2844 /* Flag bad pixel in result */
2845 for(cpl_size j = 1; j <= sy; j++ ) {
2846 for(cpl_size i = 1; i <= sx; i++ ) {
2847 hdrl_image_reject(hima, i, j);
2848 }
2849 }
2850
2851 }
2852
2853 hdrl_image_mul_scalar(hima, atm_ext);
2854 hdrl_image_mul_scalar(hima, factor);
2855 }
2856
2857 char* fcal_cube_ptag = cpl_sprintf("%s_FLUXCAL", procatg);
2858 char* cube_file_name = cpl_sprintf("%s_cube_fluxcal.fits", pipefile_prefix);
2859 cpl_propertylist_set_string(phead, FHDR_PRO_CATG, fcal_cube_ptag);
2860
2861 hdrl_resample_result* res = cpl_calloc(1, sizeof(hdrl_resample_result));
2862 res->himlist = hlist;
2863 res->header = phead;
2864 cpl_propertylist_append(res->header, head_wcs);
2865 cpl_propertylist_delete(head_wcs);
2866 cpl_propertylist_update_string(res->header, "BUNIT", UNIT_FLUXCAL);
2867
2868 /* collapse data cube */
2869 eris_print_rec_status(0);
2870 cpl_propertylist* phead2D = cpl_propertylist_duplicate(res->header);
2871 eris_ifu_plist_erase_wcs3D(phead2D);
2872 char* mean_cube_pctg = cpl_sprintf("%s_%s",fcal_cube_ptag,"MEAN");
2873 const char *filenameSpec ="cube_fluxcal";
2874 char* ima_file_name = cpl_sprintf("%s_%s_mean.fits", pipefile_prefix, filenameSpec);
2875 cpl_propertylist_append_string(phead2D, "PRODCATG", PRODCATG_WHITELIGHT);
2876 cpl_propertylist_update_string(phead2D, "BUNIT", UNIT_FLUXCAL);
2877 hdrl_image* cube_collapsed = NULL;
2878 cpl_image* cube_cmap = NULL;
2879 //TODO: AMO checking why collapsed mean of flux calibrated result is full of NANs.
2880 //hdrl_value fct = {1000000.,0};
2881 //hdrl_imagelist_mul_scalar(res->himlist, fct);
2882 hdrl_imagelist_collapse_mean(res->himlist, &cube_collapsed, &cube_cmap);
2883 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, mean_cube_pctg);
2884 cpl_image* mask_ima = cpl_image_new_from_mask(hdrl_image_get_mask(cube_collapsed));
2885 eris_ifu_save_deq_image(sof, NULL, parlist,sof, NULL,
2886 plugin_id, phead2D, "RADECSYS", ima_file_name,
2887 hdrl_image_get_image(cube_collapsed),
2888 hdrl_image_get_error(cube_collapsed),
2889 rmse, mask_ima, flag16bit, UNIT_FLUXCAL);
2890 cpl_image_delete(mask_ima);
2891 cpl_free(mean_cube_pctg);
2892 cpl_free(ima_file_name);
2893 cpl_propertylist_delete(phead2D);
2894
2895 if(strcmp(plugin_id,"eris_ifu_jitter") == 0) {
2896 char* param_name = cpl_sprintf("%s.crea-phase3", "eris.eris_ifu_jitter");
2897 const cpl_parameter* p = cpl_parameterlist_find_const(parlist, param_name);
2898
2899 cpl_boolean gen_phase3 = cpl_parameter_get_bool(p);
2900 if(gen_phase3) {
2901 cpl_msg_info(cpl_func,"Save flux calibrated COADDED OBJECT cube with phase3 data format");
2902 eris_ifu_resample_save_cube(res, fcal_cube_ptag, plugin_id,
2903 cube_file_name, parlist, sof, CPL_TRUE);
2904 } else {
2905 cpl_msg_info(cpl_func,"Save flux calibrated COADDED OBJECT cube");
2906 eris_ifu_resample_save_cube(res, fcal_cube_ptag, plugin_id,
2907 cube_file_name, parlist, sof, CPL_FALSE);
2908 }
2909 cpl_free(param_name);
2910 } else {
2911 cpl_msg_info(cpl_func,"Save flux calibrated STD STAR cube");
2912 eris_ifu_resample_save_cube(res, fcal_cube_ptag, plugin_id,
2913 cube_file_name, parlist, sof, CPL_FALSE);
2914 }
2915
2916 cpl_free(cube_file_name);
2917
2918 cpl_table* qclog = eris_qclog_init();
2919 hdrl_image_delete(cube_collapsed);
2920 cpl_image_delete(cube_cmap);
2921
2922 cpl_table_delete(qclog);
2923
2924 //hdrl_resample_result_delete(res);
2925 cpl_free(fcal_cube_ptag);
2926 hdrl_spectrum1D_delete(&resp_s1d);
2927 hdrl_spectrum1D_delete(&ext_s1d);
2928 hdrl_spectrum1D_delete(&ext_s1d_e);
2929 cpl_table_delete(atm_ext_tab);
2930 cpl_table_delete(resp_tab);
2931
2932 eris_check_error_code("eris_flux_calibrate_cube");
2933
2934 return cpl_error_get_code();
2935}
cpl_error_code eris_ifu_save_deq_image(cpl_frameset *allframes, cpl_propertylist *header, const cpl_parameterlist *parlist, const cpl_frameset *usedframes, const cpl_frame *inherit, const char *recipe, const cpl_propertylist *applist, const char *remregexp, const char *filename, const cpl_image *image, const cpl_image *error, deqErrorType errorType, const cpl_image *dataQualityMask, deqQualityType dataQualityType, const char *unit)
Write DFS pipeline product with image, error and data qual.
Definition: eris_ifu_dfs.c:300
cpl_error_code eris_ifu_resample_save_cube(hdrl_resample_result *aCube, const char *procatg, const char *recipe, const char *filename, const cpl_parameterlist *parlist, cpl_frameset *frameset, cpl_boolean gen_phase3)
save an cube (imagelist, product of pixel resampling) on a FITS file with proper FITS header
cpl_table * eris_qclog_init(void)
Initialize QC table.
cpl_error_code eris_pfits_put_qc(cpl_propertylist *plist, cpl_table *qclog)
convert table with QC parameter information to a propertylist
cpl_error_code eris_qclog_add_double(cpl_table *table, const char *key_name, const double value, const char *key_help)
add QC double info to table
cpl_error_code eris_qclog_add_int(cpl_table *table, const char *key_name, const int value, const char *key_help)
add QC int info to table
cpl_error_code eris_ifu_cube_trim_nans(const cpl_size zmin, const cpl_size zmax, hdrl_imagelist **iml, cpl_imagelist **bpm, cpl_propertylist *header)
remove planes from input iml and bpm in the range [0, zmin), (zmax,sz]
double eris_pfits_get_crpix2(const cpl_propertylist *plist)
find out the character string associated to the CRPIX2 keyword
Definition: eris_pfits.c:736
double eris_pfits_get_cd33(const cpl_propertylist *plist)
find out the character string associated to the CDELT3 keyword
Definition: eris_pfits.c:686
double eris_pfits_get_crpix1(const cpl_propertylist *plist)
find out the character string associated to the CRPIX1 keyword
Definition: eris_pfits.c:748
double eris_pfits_get_crval3(const cpl_propertylist *plist)
find out the character string associated to the CVRVAL3 keyword
Definition: eris_pfits.c:699
double eris_pfits_get_crpix3(const cpl_propertylist *plist)
find out the character string associated to the CRPIX3 keyword
Definition: eris_pfits.c:723
cpl_error_code eris_check_error_code(const char *func_id)
handle CPL errors
Definition: eris_utils.c:56
hdrl_parameter * hdrl_response_parameter_create(const hdrl_value Ap, const hdrl_value Am, const hdrl_value G, const hdrl_value Tex)
ctor for the hdrl_parameter for response
hdrl_parameter * hdrl_efficiency_parameter_create(const hdrl_value Ap, const hdrl_value Am, const hdrl_value G, const hdrl_value Tex, const hdrl_value Atel)
ctor for the hdrl_parameter for efficiency
hdrl_spectrum1D * hdrl_efficiency_compute(const hdrl_spectrum1D *I_std_arg, const hdrl_spectrum1D *I_std_ref, const hdrl_spectrum1D *E_x, const hdrl_parameter *pars)
efficiency calculation
cpl_error_code hdrl_image_set_pixel(hdrl_image *self, cpl_size xpos, cpl_size ypos, hdrl_value value)
set pixel values of hdrl_image
Definition: hdrl_image.c:594
cpl_error_code hdrl_image_mul_scalar(hdrl_image *self, hdrl_value value)
Elementwise multiplication of an image with a scalar.
cpl_mask * hdrl_image_get_mask(hdrl_image *himg)
get cpl bad pixel mask from image
Definition: hdrl_image.c:157
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
Definition: hdrl_image.c:131
cpl_size hdrl_image_get_size_x(const hdrl_image *self)
return size of X dimension of image
Definition: hdrl_image.c:525
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:105
cpl_error_code hdrl_image_reject(hdrl_image *self, cpl_size xpos, cpl_size ypos)
mark pixel as bad
Definition: hdrl_image.c:427
hdrl_image * hdrl_image_new(cpl_size nx, cpl_size ny)
create new zero filled hdrl image
Definition: hdrl_image.c:311
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
Definition: hdrl_image.c:379
cpl_error_code hdrl_imagelist_collapse_mean(const hdrl_imagelist *himlist, hdrl_image **out, cpl_image **contrib)
Mean collapsing of image list.
cpl_size hdrl_imagelist_get_size_y(const hdrl_imagelist *himlist)
Get number of rows of images in the imagelist.
hdrl_imagelist * hdrl_imagelist_create(cpl_imagelist *imlist, cpl_imagelist *errlist)
Create an hdrl_imagelist out of 2 cpl_imagelist.
void hdrl_imagelist_delete(hdrl_imagelist *himlist)
Free all memory used by a hdrl_imagelist object including the images.
cpl_size hdrl_imagelist_get_size(const hdrl_imagelist *himlist)
Get the number of images in the imagelist.
cpl_size hdrl_imagelist_get_size_x(const hdrl_imagelist *himlist)
Get number of colums of images in the imagelist.
hdrl_image * hdrl_imagelist_get(const hdrl_imagelist *himlist, cpl_size inum)
Get an image from a list of images.
void hdrl_parameter_delete(hdrl_parameter *obj)
shallow delete of a parameter
const hdrl_spectrum1D * hdrl_response_result_get_final_response(const hdrl_response_result *res)
Getter for the final response contained inside the hdrl_response_result.
const hdrl_spectrum1D * hdrl_response_result_get_selected_response(const hdrl_response_result *res)
Getter for the selected response contained inside the hdrl_response_result.
hdrl_parameter * hdrl_response_fit_parameter_create(const cpl_size radius, const cpl_array *fit_points, const hdrl_data_t wrange, const cpl_bivector *high_abs_regions)
ctor for the hdrl_parameter for the final interpolation of the response
void hdrl_response_result_delete(hdrl_response_result *res)
Destructor for hdrl_response_result.
hdrl_response_result * hdrl_response_compute(const hdrl_spectrum1D *obs_s, const hdrl_spectrum1D *ref_s, const hdrl_spectrum1D *E_x, const hdrl_parameter *telluric_par, const hdrl_parameter *velocity_par, const hdrl_parameter *calc_par, const hdrl_parameter *fit_par)
Computation of the response.
hdrl_parameter * hdrl_response_telluric_evaluation_parameter_create(const hdrl_spectrum1Dlist *telluric_models, hdrl_data_t w_step, cpl_size half_win, cpl_boolean normalize, cpl_boolean shift_in_log_scale, const cpl_bivector *quality_areas, const cpl_bivector *fit_areas, hdrl_data_t lmin, hdrl_data_t lmax)
ctor for the hdrl_parameter for the telluric evaluation
const hdrl_spectrum1D * hdrl_response_result_get_raw_response(const hdrl_response_result *res)
Getter for the raw response contained inside the hdrl_response_result.
hdrl_parameter * hdrl_spectrum1D_shift_fit_parameter_create(const hdrl_data_t wguess, const hdrl_data_t range_wmin, const hdrl_data_t range_wmax, const hdrl_data_t fit_wmin, const hdrl_data_t fit_wmax, const hdrl_data_t fit_half_win)
The function create a hdrl_spectrum1D_shift_parameter to be used in hdrl_spectrum1D_compute_shift_fit...
cpl_error_code hdrl_spectrum1D_mul_scalar(hdrl_spectrum1D *self, hdrl_value scalar_operator)
computes the elementwise multiplication of a spectrum by a scalar, the self parameter is modified
cpl_error_code hdrl_spectrum1D_wavelength_mult_scalar_linear(hdrl_spectrum1D *self, hdrl_data_t scale_linear)
computes the elementwise multiplication of the scalar for the wavelength. The scalar is assumed to be...
hdrl_spectrum1D * hdrl_spectrum1D_resample(const hdrl_spectrum1D *self, const hdrl_spectrum1D_wavelength *waves, const hdrl_parameter *par)
resample a hdrl_spectrum1D on the wavelengths contained in waves
hdrl_spectrum1Dlist * hdrl_spectrum1Dlist_new(void)
hdrl_spectrum1Dlist default constructor
hdrl_spectrum1D * hdrl_spectrum1D_mul_spectrum_create(const hdrl_spectrum1D *f1, const hdrl_spectrum1D *f2)
multiply one spectrum by another spectrum
hdrl_parameter * hdrl_spectrum1D_resample_interpolate_parameter_create(const hdrl_spectrum1D_interpolation_method method)
constructor for the hdrl_parameter in the case of interpolation
cpl_size hdrl_spectrum1D_get_size(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter for size
void hdrl_spectrum1D_delete(hdrl_spectrum1D **p_self)
hdrl_spectrum1D destructor
cpl_error_code hdrl_spectrum1Dlist_set(hdrl_spectrum1Dlist *self, hdrl_spectrum1D *s, const cpl_size idx)
hdrl_spectrum1Dlist setter of the i-th element
void hdrl_spectrum1Dlist_delete(hdrl_spectrum1Dlist *l)
hdrl_spectrum1Dlist destructor
hdrl_spectrum1D * hdrl_spectrum1D_convert_from_table(const cpl_table *self, const char *flux_col_name, const char *wavelength_col_name, const char *flux_e_col_name, const char *flux_bpm_col_name, hdrl_spectrum1D_wave_scale scale)
convert a table to a spectrum
hdrl_spectrum1D * hdrl_spectrum1D_create(const cpl_image *arg_flux, const cpl_image *arg_flux_e, const cpl_array *wavelength, hdrl_spectrum1D_wave_scale wave_scale)
hdrl_spectrum1D default constructor
hdrl_spectrum1D_wave_scale hdrl_spectrum1D_get_scale(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter for scale
hdrl_spectrum1D_wavelength hdrl_spectrum1D_get_wavelength(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter for wavelengths
hdrl_data_t hdrl_spectrum1D_get_wavelength_value(const hdrl_spectrum1D *self, int idx, int *rej)
hdrl_spectrum1D getter for a wavelength value
cpl_table * hdrl_spectrum1D_convert_to_table(const hdrl_spectrum1D *self, const char *flux_col_name, const char *wavelength_col_name, const char *flux_e_col_name, const char *flux_bpm_col_name)
converts a spectrum in a table.
cpl_error_code hdrl_spectrum1D_div_scalar(hdrl_spectrum1D *self, hdrl_value scalar_operator)
computes the elementwise division of a spectrum by a scalar, the self parameter is modified
hdrl_value hdrl_spectrum1D_get_flux_value(const hdrl_spectrum1D *self, int idx, int *rej)
hdrl_spectrum1D getter for a flux value