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