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