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"
36#include <eris_ifu_resample.h>
37const hdrl_spectrum1D_wave_scale
38global_scale = hdrl_spectrum1D_wave_scale_linear;
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"
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
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
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
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
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);
106eris_get_wcs_cube(
const cpl_propertylist* plist,
double *clambda,
107 double *dis,
double *cpix,
double *cx,
double *cy)
117 return cpl_error_get_code();
130sinfo_get_dispersion(
const char* band)
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;
152static double eris_get_dispersion(
const char* band)
154 cpl_ensure(band != NULL, CPL_ERROR_NULL_INPUT, 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;
196eris_get_band(cpl_frame * ref_frame,
char * band)
198 cpl_msg_debug(cpl_func,
"ref_frame: %p filename: %s",
199 (
const void*)ref_frame, cpl_frame_get_filename(ref_frame));
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);
204 char* ref_file = NULL;
205 cpl_propertylist* plist = NULL;
207 ref_file = cpl_strdup(cpl_frame_get_filename(ref_frame)) ;
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",
212 cpl_propertylist_delete(plist) ;
213 return cpl_error_get_code() ;
216 if (cpl_propertylist_has(plist, KEY_NAME_FILT_NAME)) {
217 strcpy(band, cpl_propertylist_get_string(plist, KEY_NAME_FILT_NAME));
221 cpl_msg_warning(cpl_func,
"keyword %s does not exist",
223 return cpl_error_get_code() ;
227 cpl_propertylist_delete(plist);
229 return cpl_error_get_code() ;
243sinfo_get_band(cpl_frame * ref_frame,
char * band)
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);
248 char* ref_file = NULL;
249 cpl_propertylist* plist = NULL;
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",
255 cpl_propertylist_delete(plist) ;
256 return cpl_error_get_code();
259 if (cpl_propertylist_has(plist, KEY_NAME_FILT_NAME2)) {
260 strcpy(band, cpl_propertylist_get_string(plist, KEY_NAME_FILT_NAME2));
264 cpl_msg_warning(cpl_func,
"keyword %s does not exist",
265 KEY_NAME_FILT_NAME2);
266 return cpl_error_get_code();
270 cpl_propertylist_delete(plist);
271 return cpl_error_get_code();
273cpl_boolean eris_can_extract(cpl_frameset* sof) {
275 cpl_boolean has_cube_coadd = CPL_TRUE;
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)
287 cpl_msg_warning(cpl_func,
"Provide %s or %s or %s or %s or %s or %s or %s or %s"
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
298 has_cube_coadd = CPL_FALSE;
300 return has_cube_coadd;
302int eris_can_flux_calibrate(cpl_frameset* sof) {
304 int has_extcoeff_table=1;
305 int has_response = 1;
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);
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;
322 return has_response *
328int eris_can_compute_efficiency(cpl_frameset* sof) {
330 int has_std_star_spectra=1;
331 int has_extcoeff_table=1;
332 int has_response = 1;
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;
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);
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;
354 return has_std_star_spectra *
361int eris_can_compute_response(cpl_frameset* sof) {
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;
371 int has_resp_windows=1;
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;
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;
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;
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;
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;
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);
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);
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);
428 return has_std_star_spectra *
429 has_flux_std_catalog *
430 has_tel_mod_catalog *
431 has_resp_fit_points_cat *
433 has_high_abs_regions *
443eris_eff_qc(cpl_table* tot_eff, cpl_frame* frm_eff_wind_qc, cpl_table** qclog_tbl )
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);
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;
463 sto_eff = cpl_table_new(0);
464 cpl_table_copy_structure(sto_eff, tot_eff);
466 int ncontributes = 0;
469 for(
int i=0;i<nrow;i++) {
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);
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");
489 key_name = cpl_sprintf(
"QC EFF WIN%d WLMIN",i);
491 "[um] Min window wavelength for eff comp");
494 key_name = cpl_sprintf(
"QC EFF WIN%d WLMAX",i);
496 "[um] Max window wavelength for eff comp");
499 key_name = cpl_sprintf(
"QC EFF WIN%d MEAN",i);
501 "Mean efficiency on window");
504 key_name = cpl_sprintf(
"QC EFF WIN%d MEDIAN",i);
506 "Median efficiency on window");
509 key_name = cpl_sprintf(
"QC EFF WIN%d MIN",i);
511 "Min efficiency on window");
514 key_name = cpl_sprintf(
"QC EFF WIN%d MAX",i);
516 "Max efficiency on window");
519 key_name = cpl_sprintf(
"QC EFF WIN%d RMS",i);
521 "RMS efficiency on window");
524 cpl_table_select_all(tot_eff);
525 cpl_table_delete(tmp_eff);
529 wmin=cpl_table_get_column_min(sto_eff,
"WAVE");
530 wmax=cpl_table_get_column_max(sto_eff,
"WAVE");
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");
539 key_name = cpl_sprintf(
"QC EFF NWIN");
541 "Number of windows used for eff comp");
543 key_name = cpl_sprintf(
"QC EFF WLMIN");
545 "[um] Min wavelength for eff comp");
548 key_name = cpl_sprintf(
"QC EFF WLMAX");
550 "[um] Max wavelength for eff comp");
553 key_name = cpl_sprintf(
"QC EFF MEAN");
558 key_name = cpl_sprintf(
"QC EFF MEDIAN");
560 "Median efficiency");
563 key_name = cpl_sprintf(
"QC EFF MIN");
568 key_name = cpl_sprintf(
"QC EFF MAX");
573 key_name = cpl_sprintf(
"QC EFF RMS");
580 cpl_table_delete(ew_tab);
581 cpl_table_delete(sto_eff);
583 return cpl_error_get_code();
587eris_ifu_table_from_sdp_to_normal_format(cpl_frame** frm_sci)
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))) {
598 const cpl_array* array_flux = NULL;
599 const cpl_array* array_err = NULL;
600 cpl_size nrows = cpl_array_get_size(array_wave);
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);
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));
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);
620 sci_tab = cpl_table_load(
"spc_tab.fits", 1, 0);
621 cpl_frame_set_filename(*frm_sci,
"spc_tab.fits");
624 cpl_propertylist_delete(plist);
630eris_remove_table_normal_format(
void) {
631 char* cmd = cpl_sprintf(
"rm spc_tab.fits");
632 int status = system(cmd);
634 cpl_msg_warning(cpl_func,
"call to system failed");
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)
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;
667 cpl_table* tbl_eff = NULL;
668 cpl_table* tbl_ref = NULL;
669 cpl_table* tbl_atmext = NULL;
670 cpl_table* sci_tab = NULL;
678 double dEpsilon = 0.1;
680 double um2nm = 1000.;
681 double nm2um = 0.001;
682 double m2tocm2 = 1.e4;
686 cpl_boolean is_sinfoni = CPL_FALSE;
688 const char* name_sci = NULL;
689 const char* name_atm = NULL;
691 name_sci = cpl_frame_get_filename(frm_sci);
693 plist = cpl_propertylist_load(name_sci, 0);
694 const char* instrume = cpl_propertylist_get_string(plist,
"INSTRUME");
696 if(strcmp(instrume,
"SINFONI") == 0) {
697 is_sinfoni = CPL_TRUE;
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);
720 exptime = cpl_propertylist_get_double(plist, FHDR_S_DIT);
722 exptime = cpl_propertylist_get_double(plist, FHDR_E_DIT);
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);
729 name_atm = cpl_frame_get_filename(frm_atmext);
730 tbl_atmext = cpl_table_load(name_atm, 1, 0);
734 eris_parse_catalog_std_stars(frm_cat, dRA, dDEC, dEpsilon, &tbl_ref);
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);
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;
756 "wavelength", NULL, NULL, scale);
760 "WAVE", NULL, NULL, scale);
774 char band[FILE_NAME_SZ];
776 sinfo_get_band(frm_sci, band);
778 eris_get_band(frm_sci, band);
781 cpl_msg_debug(cpl_func,
"band: %s", band);
783 double conv_factor = 0;
787 cpl_size size = cpl_table_get_nrow(sci_tab);
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);
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);
804 bin_size = sinfo_get_dispersion(band);
806 bin_size = eris_get_dispersion(band);
813 conv_factor = bin_size * um2nm * nm2AA;
814 hdrl_value
operator = {conv_factor, 0.};
819 "LAMBDA", NULL, NULL, scale);
822 "LAMBDA", NULL, NULL, scale);
825 hdrl_parameter* eff_pars;
827 hdrl_value Ap = {0,0};
828 hdrl_value Am = {airmass,0};
829 hdrl_value G = {gain,0};
830 hdrl_value Tex = {exptime,0};
833 hdrl_value Atel = {TELESCOPE_SURFACE * m2tocm2,0};
844 if( tbl_eff != NULL && frm_eff_wind_qc != NULL) {
846 eris_eff_qc(tbl_eff, frm_eff_wind_qc, &qclog_tbl);
848 cpl_table_delete(qclog_tbl);
852 cpl_propertylist_update_string(plist, CPL_DFS_PRO_CATG,
853 ERIS_IFU_PRO_JITTER_EFFICIENCY);
854 cpl_propertylist_erase(plist,
"BUNIT");
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);
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);
881static hdrl_spectrum1D *
882eris_std_cat_frame_to_hdrl_spectrum1D(cpl_frame* frm_std_cat,
const double dRA,
885 cpl_ensure(frm_std_cat != NULL, CPL_ERROR_NULL_INPUT, NULL);
886 hdrl_spectrum1D * std_s1d;
889 cpl_table* tbl_ref = NULL;
893 eris_parse_catalog_std_stars(frm_std_cat,dRA,dDEC,dEpsilon,&tbl_ref);
895 if(tbl_ref == NULL) {
896 cpl_msg_error(cpl_func,
"Provide std sar catalog frame");
901 "LAMBDA", NULL, NULL, global_scale);
902 cpl_table_delete(tbl_ref);
909static hdrl_spectrum1D *
910eris_sci_frame_to_hdrl_spectrum1D(cpl_frame* frm_sci)
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);
922 "WAVE", NULL, NULL, global_scale);
923 cpl_table_delete(sci_tab);
929static hdrl_spectrum1D *
930eris_atmext_frame_to_hdrl_spectrum1D(cpl_frame* frm_atmext)
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);
936 cpl_table* tbl_atmext = cpl_table_load(name_atm,1,0);
940 "LAMBDA", NULL, NULL, global_scale);
941 cpl_table_delete(tbl_atmext);
946static hdrl_spectrum1Dlist *
947eris_get_telluric_models(
const cpl_frame * telluric_cat){
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);
955 for(cpl_size i = 0; i < next; ++i){
956 cpl_table * tab = cpl_table_load(cat_name, 1 + i, 1);
958 hdrl_spectrum1D * s =
965 cpl_table_delete(tab);
976eris_read_wlen_windows(
const cpl_frame * fit_areas_frm){
978 cpl_ensure(fit_areas_frm != NULL, CPL_ERROR_NULL_INPUT, NULL);
979 cpl_table * tab = NULL;
981 tab = cpl_table_load(cpl_frame_get_filename(fit_areas_frm),1,0);
984 const cpl_size nrow = cpl_table_get_nrow(tab);
986 double * pwmin = cpl_table_unwrap(tab,
"LAMBDA_MIN");
987 double * pwmax = cpl_table_unwrap(tab,
"LAMBDA_MAX");
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);
999static hdrl_parameter*
1000eris_response_parameters_calc_par(cpl_propertylist* plist)
1002 cpl_ensure(plist != NULL, CPL_ERROR_NULL_INPUT, NULL);
1003 hdrl_value Ap = {0, 0};
1006 double airmass = eris_ifu_get_airmass_mean(plist);
1007 hdrl_value Am = {airmass, 0.};
1011 double gain=ERIS_GAIN;
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);
1019 hdrl_value Tex= {exptime, 0.};
1028eris_read_fit_points(
const cpl_frame * frm_fit_points,
1029 const double ra,
const double dec,
const double ra_dec_tolerance){
1031 cpl_ensure(frm_fit_points != NULL, CPL_ERROR_NULL_INPUT, NULL);
1032 const char * fname = cpl_frame_get_filename(frm_fit_points);
1034 cpl_table * index_table = cpl_table_load(fname, 1, CPL_FALSE);
1035 const cpl_size sz = cpl_table_get_nrow(index_table);
1037 cpl_size selected_ext = -1;
1039 for(cpl_size i = 0; i < sz; ++i){
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;
1050 cpl_table_delete(index_table);
1051 cpl_ensure(selected_ext >= 0, CPL_ERROR_ILLEGAL_OUTPUT, NULL);
1054 cpl_table * tb = cpl_table_load(fname, selected_ext, CPL_FALSE);
1056 const cpl_size sz_points = cpl_table_get_nrow(tb);
1057 cpl_array * fit_points = cpl_array_new(sz_points, CPL_TYPE_DOUBLE);
1059 for(cpl_size i = 0; i < sz_points; ++i){
1061 const double d = cpl_table_get(tb,
"LAMBDA", i, &rej);
1062 cpl_array_set(fit_points, i, d);
1065 cpl_table_delete(tb);
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);
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);
1084hdrl_spectrum1D * extend_hdrl_spectrum(
const double wmin,
const double wmax,
1085 const hdrl_spectrum1D * s){
1087 cpl_ensure(wmin < wmax, CPL_ERROR_ILLEGAL_INPUT, NULL);
1088 cpl_ensure(s != NULL, CPL_ERROR_NULL_INPUT, NULL);
1091 const double s_wmin = cpl_array_get_min(lambadas_s);
1092 const double s_wmax = cpl_array_get_max(lambadas_s);
1094 const cpl_size sz_s = cpl_array_get_size(lambadas_s);
1095 const double step = (s_wmax - s_wmin)/sz_s;
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);
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));
1106 double lambda = s_wmin;
1107 for(cpl_size i = n_before - 1; i >= 0; i--){
1109 cpl_array_set(new_lambdas, i, lambda);
1114 for(cpl_size i = n_before + sz_s; i < new_size; i++){
1116 cpl_array_set(new_lambdas, i, lambda);
1120 for(cpl_size i = n_before; i < n_before + sz_s; i++){
1121 const cpl_size idx_ori = i - n_before;
1124 cpl_array_set(new_lambdas, i, lambda);
1137 hdrl_spectrum1D * to_ret =
1144 cpl_array_delete(new_lambdas);
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 )
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);
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);
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;
1178 cpl_table* sto_resp = cpl_table_new(0);
1179 cpl_table_copy_structure(sto_resp, resp);
1181 int ncontributes = 0;
1183 for(
int i=0;i<nrow;i++) {
1184 cpl_table_select_all(resp);
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);
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);
1206 key_name = cpl_sprintf(
"QC RESP WIN%d WLMIN",i);
1208 "[um] Min window wavelength for resp comp");
1211 key_name = cpl_sprintf(
"QC RESP WIN%d WLMAX",i);
1213 "[um] Max window wavelength for resp comp");
1217 key_name = cpl_sprintf(
"QC RESP WIN%d MEAN",i);
1219 "Mean response on window");
1222 key_name = cpl_sprintf(
"QC RESP WIN%d MEDIAN",i);
1224 "Median response on window");
1227 key_name = cpl_sprintf(
"QC RESP WIN%d MIN",i);
1229 "Min response on window");
1232 key_name = cpl_sprintf(
"QC RESP WIN%d MAX",i);
1234 "Max response on window");
1237 key_name = cpl_sprintf(
"QC RESP WIN%d RMS",i);
1239 "RMS response on window");
1242 cpl_table_select_all(resp);
1243 cpl_table_delete(tmp_resp);
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);
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);
1262 key_name = cpl_sprintf(
"QC RESP NWIN");
1264 "Number of windows used for resp comp");
1266 key_name = cpl_sprintf(
"QC RESP WLMIN");
1268 "[um] Min wavelength for resp comp");
1271 key_name = cpl_sprintf(
"QC RESP WLMAX");
1273 "[um] Max wavelength for resp comp");
1276 key_name = cpl_sprintf(
"QC RESP MEAN");
1281 key_name = cpl_sprintf(
"QC RESP MEDIAN");
1286 key_name = cpl_sprintf(
"QC RESP MIN");
1291 key_name = cpl_sprintf(
"QC RESP MAX");
1296 key_name = cpl_sprintf(
"QC RESP RMS");
1301 cpl_table_delete(tmp_resp);
1302 cpl_table_delete(ew_tab);
1303 cpl_table_delete(sto_resp);
1305 return cpl_error_get_code();
1310eris_log_pro(
char* name_o,
1311 const char* pro_catg,
1313 cpl_frameset* ref_set,
1314 cpl_frameset** out_set,
1315 cpl_propertylist** plist,
1316 const cpl_parameterlist* parlist,
1319 cpl_frame* product_frame = NULL ;
1320 char * pipe_id=NULL;
1321 cpl_errorstate initial_errorstate = cpl_errorstate_get();
1323 pipe_id = cpl_sprintf(
"%s%s",
"eris/",PACKAGE_VERSION);
1324 product_frame = cpl_frame_new() ;
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);
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);
1341 cpl_frameset_insert(*out_set, product_frame);
1350static char * eris_new_get_rootname(
const char * filename)
1352 static char path[MAX_STR_SIZE+1];
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")))
1367 lastdot[0] = (char)0;
1376eris_check_name(
const char* in,
char** ou,
char** paf) {
1379 if (strstr(in,
"." ) != NULL ) {
1380 char* tmp = eris_new_get_rootname(in);
1381 name_b= cpl_sprintf(
"%s",tmp);
1383 name_b = cpl_sprintf(
"%s",in);
1385 *ou = cpl_sprintf(
"%s.fits",name_b);
1393 *paf = cpl_sprintf(
"%s.paf",name_b);
1404 const char * out_file,
1405 const char * pro_catg,
1408 const cpl_parameterlist* parlist)
1411 char * name_o =NULL;
1412 char * name_p =NULL;
1413 cpl_propertylist * plist=NULL ;
1414 char* ref_file=NULL;
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)) ;
1423 eris_check_name(out_file, &name_o, &name_p);
1424 cpl_msg_info(cpl_func,
"Writing tbl %s pro catg %s" , name_o, pro_catg) ;
1427 if ((cpl_error_code)((plist = cpl_propertylist_load(ref_file,0)) == NULL))
1429 cpl_msg_error(cpl_func,
"getting header from tbl frame %s",ref_file);
1430 cpl_propertylist_delete(plist) ;
1434 cpl_frameset_iterator_delete(it);
1443 eris_log_pro(name_o, pro_catg, CPL_FRAME_TYPE_TABLE,
1444 ref,&set,&plist,parlist,recid);
1449 cpl_propertylist_append_string(plist,
"PRODCATG",
"SCIENCE.SPECTRUM");
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) ;
1458 cpl_frameset_iterator_delete(it);
1463 cpl_propertylist_delete(plist) ;
1464 cpl_msg_indent_less() ;
1468 cpl_frameset_iterator_delete(it);
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) {
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);
1499 gsl_vector *x, *y, *w, *p;
1500 gsl_matrix *X, *covar;
1502 gsl_multifit_linear_workspace *work;
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);
1513 cpl_msg_info(cpl_func,
"n = %d, deg = %d", n, deg);
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]);
1520 gsl_vector_set (w, i, 1.0/err/err);
1522 for (j = 0; j < deg; j++)
1523 gsl_matrix_set (X, i, j, gsl_pow_int(gsl_vector_get(x, i), j));
1527 cpl_msg_info(cpl_func,
"vectors set");
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);
1535 cpl_msg_info(cpl_func,
"Fit done");
1538 for (i = 0; i < n; i++) {
1540 for (j = 0; j < deg; j++)
1541 fit_RE[i] += gsl_matrix_get (X, i, j) * gsl_vector_get (p, j);
1544 chisq = chisq/(n-deg);
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));
1558 gsl_matrix_free(covar);
1561 return cpl_error_get_code();
1565eris_fit_poly_and_clip_response(
const cpl_table* resp_pts,
1566 const cpl_table* resp_smo,
const int deg,
const double kappa,
1570 cpl_ensure(resp_pts != NULL, CPL_ERROR_NULL_INPUT, NULL);
1571 cpl_ensure(resp_smo != NULL, CPL_ERROR_NULL_INPUT, NULL);
1574 eris_poly_data *p_data = (eris_poly_data *)cpl_malloc(
sizeof(eris_poly_data));
1576 cpl_size pol_deg = deg;
1578 cpl_table* resp_pts_tmp = cpl_table_duplicate(resp_pts);
1580 double* coeffs = NULL;
1581 double* coeffs_err = NULL;
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));
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));
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.;
1603 eris_fit_poly(p_data, pol_deg, fit, coeffs, coeffs_err, &chisq, 0);
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);
1610 for (
int c = pol_deg-1; c >= 0; c--) {
1611 y2 = y2 * x2 + coeffs[c];
1613 cpl_table_set_double(resp_pts_tmp,
"res",j, y2-p_data->y[j]);
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,
1626 xtab = cpl_table_extract_selected(resp_pts_tmp);
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);
1634 cpl_free(coeffs_err);
1636 cpl_table_delete(xtab);
1637 cpl_free(p_data->x);
1638 cpl_free(p_data->y);
1639 cpl_free(p_data->err);
1642 cpl_size next = cpl_table_get_nrow(resp_pts_tmp);
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));
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.;
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);
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");
1676 for (
int j = 0; j < poly_nval; j++) {
1677 x2 = cpl_table_get(resp_poly,
"wavelength", j, NULL);
1680 for (
int c = pol_deg-1; c >= 0; c--) {
1681 y2 = y2 * x2 + coeffs[c];
1690 cpl_table_set_double(resp_poly,
"response_poly",j, y2);
1692 cpl_table_set_double(resp_poly,
"response_poly_error",j, 0.1*y2 );
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);
1709static cpl_error_code
1710eris_get_wave_shift_param(
const char* band,
double* wmin,
double* wmax,
1711 cpl_boolean* enable)
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)) {
1722 *enable = CPL_FALSE;
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)) {
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)) {
1738 *enable = CPL_FALSE;
1739 }
else if (strcmp(band,
"H+K") == 0) {
1743 *enable = CPL_FALSE;
1747 return cpl_error_get_code();
1769eris_response_compute(
const char* plugin_id,
const cpl_parameterlist* config,
1770 cpl_frameset* ref_set,cpl_frameset* sof)
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);
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);
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);
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);
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);
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);
1805 double dRA = cpl_propertylist_get_double(plist,
"RA");
1807 double dDEC = cpl_propertylist_get_double(plist,
"DEC");
1812 double sci_wmin = 0;
1813 double sci_wmax = 0;
1814 cpl_frame* frm_sci_dup = cpl_frame_duplicate(frm_sci);
1816 sci_tab = eris_ifu_table_from_sdp_to_normal_format(&frm_sci_dup);
1818 sci_wmin = cpl_table_get_column_min(sci_tab,
"WAVE");
1819 sci_wmax = cpl_table_get_column_max(sci_tab,
"WAVE");
1821 double um2nm = 1.e3;
1823 hdrl_spectrum1D * sci_s1d = eris_sci_frame_to_hdrl_spectrum1D(frm_sci_dup);
1825 cpl_frame_delete(frm_sci_dup);
1833 char band[FILE_NAME_SZ];
1834 eris_get_band(frm_sci,band);
1836 double conv_factor = 0;
1839 double bin_size=eris_get_dispersion(band);
1841 cpl_size size = cpl_table_get_nrow(sci_tab);
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);
1849 conv_factor = bin_size * um2nm * nm2AA;
1853 hdrl_spectrum1D * std_s1d =
1854 eris_std_cat_frame_to_hdrl_spectrum1D(frm_std_cat, dRA, dDEC);
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;
1862 hdrl_spectrum1D * ext_s1d = eris_atmext_frame_to_hdrl_spectrum1D(frm_atmext);
1867 hdrl_spectrum1Dlist * telluric_models =
1868 eris_get_telluric_models(frm_tell_cat);
1870 cpl_bivector * fit_areas = eris_read_wlen_windows(frm_resp_fit_areas);
1872 cpl_bivector * quality_areas = eris_read_wlen_windows(frm_resp_quality_areas);
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);
1889 cpl_msg_debug(cpl_func,
"band=%s wmin=%g wmax=%g", band, wmin, wmax);
1891 const hdrl_data_t lmin = (hdrl_data_t)wmin;
1892 const hdrl_data_t lmax = (hdrl_data_t)wmax;
1895 const cpl_boolean xcorr_normalize = (cpl_boolean) CPL_TRUE;
1896 const cpl_boolean shift_in_log_scale = (cpl_boolean) CPL_TRUE;
1898 const hdrl_data_t wguess = 1094.12;
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;
1905 hdrl_parameter* shift_par = NULL;
1906 if(velocity_enable) {
1908 range_wmin, range_wmax, fit_wmin, fit_wmax, fit_half_win);
1911 hdrl_parameter* resp_calc_par = eris_response_parameters_calc_par(plist);
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);
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;
1922 xcorr_w_step, xcorr_half_win, xcorr_normalize, shift_in_log_scale, quality_areas,
1923 fit_areas, lmin, lmax);
1926 const cpl_size post_smoothing_radius = 11;
1929 const hdrl_data_t fit_wrange = 4.0;
1930 hdrl_parameter* resp_par =
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);
1939 cpl_msg_info(cpl_func,
"compute response");
1940 hdrl_response_result * response =
1942 shift_par, resp_calc_par, resp_par);
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",
1962 const hdrl_spectrum1D * resp_s1d_smo =
1965 const hdrl_spectrum1D * resp_s1d_raw =
1970 "response",
"wavelength",
"response_error",
"response_bpm");
1972 sci_wmin=cpl_table_get_column_min(tab,
"wavelength");
1973 sci_wmax=cpl_table_get_column_max(tab,
"wavelength");
1976 cpl_table* resp_smo = NULL;
1977 cpl_table* resp_raw = NULL;
1978 cpl_table* resp_pts = NULL;
1980 "response_smo",
"wavelength",
"response_smo_error",
1981 "response_smo_bpm");
1984 "response_raw",
"wavelength",
"response_raw_error",
1985 "response_raw_bpm");
1989 "response_fit_pts",
"wavelength",
"response_fit_pts_error",
1990 "response_fit_pts_bpm");
1991 cpl_table_delete(tab);
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);
2000 cpl_frame* frm_resp_wind_qc = cpl_frameset_find(sof, ERIS_IFU_CALIB_RESPONSE_WINDOWS);
2004 const double wmin_g = cpl_array_get_min(waves) / um2nm;
2005 const double wmax_g = cpl_array_get_max(waves) / um2nm;
2013 const cpl_parameter* p;
2015 const char* context =
"eris.eris_ifu_stdstar";
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);
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) {
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);
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);
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);
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);
2061 cpl_bivector_delete(quality_areas);
2062 cpl_bivector_delete(fit_areas);
2073 cpl_propertylist_delete(plist);
2076 eris_remove_table_normal_format();
2077 cpl_array_delete(fit_points);
2080 if(cpl_error_get_code() == CPL_ERROR_ILLEGAL_OUTPUT) {
2084 cpl_error_set(cpl_func,CPL_ERROR_NONE);
2085 return CPL_ERROR_NONE;
2087 return cpl_error_get_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){
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);
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);
2119 cpl_msg_info(cpl_func,
"Flux calibrate the spectrum");
2120 hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
2122 cpl_frame* spectra_frm = cpl_frameset_find(sof, ERIS_IFU_PRO_JITTER_SPECTRUM);
2123 fname = cpl_frame_get_filename(spectra_frm);
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);
2134 double um2nm = 1.e3;
2137 double sci_wmin = 0;
2138 double sci_wmax = 0;
2139 double resp_wmin = 0;
2140 double resp_wmax = 0;
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);
2148 double bin_size_scale_factor = bin_size * um2nm * nm2AA;
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);
2154 cpl_propertylist* plist = cpl_propertylist_load(fname,0);
2155 double dit = cpl_propertylist_get_double(plist, FHDR_E_DIT);
2158 double exptime = dit;
2159 double gain = ERIS_GAIN;
2162 double airmass = eris_ifu_get_airmass_mean(plist);
2163 cpl_msg_info(cpl_func,
"airmass: %g", airmass);
2164 cpl_propertylist_delete(plist);
2166 double factor = 1. / (gain * exptime * bin_size_scale_factor);
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;
2180 resp_tab = cpl_table_load(fname, 4, 0);
2182 resp_tab = cpl_table_load(fname, 1, 0);
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);
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);
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);
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 =
2207 "WAVE",
"ERR", NULL, scale);
2208 cpl_table_delete(spectra_xtab);
2217 hdrl_spectrum1D * resp_s1d = NULL;
2220 "wavelength",
"response_poly_error",
"response_poly_bpm", scale);
2223 "wavelength",
"response_smo_error",
"response_smo_bpm", scale);
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);
2230 hdrl_parameter * res_par =
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);
2249 hdrl_spectrum1D * ext_s1d =
2251 "LAMBDA", NULL, NULL, scale);
2253 hdrl_spectrum1D * ext_s1d_e = extend_hdrl_spectrum(sci_wmin * um2nm, sci_wmax * um2nm, ext_s1d);
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);
2277 perr[i] *= pow(10., exp);
2282 "flux",
"wavelength",
"flux_error", NULL);
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();
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);
2310 return cpl_error_get_code();
2376eris_flux_calibrate_cube2(
const char* procatg,
const char* pipefile_prefix,
2377 const char* plugin_id,
const cpl_parameterlist* parlist,
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);
2384 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
2385 cpl_ensure_code(sof, CPL_ERROR_NULL_INPUT);
2387 cpl_size sof_sz = cpl_frameset_get_size(sof);
2388 cpl_ensure_code(sof_sz >0, CPL_ERROR_DATA_NOT_FOUND);
2393 cpl_msg_info(cpl_func,
"Flux calibrate the cube %s", procatg);
2396 hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
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);
2405 cpl_imagelist_delete(cube);
2406 cpl_imagelist_delete(errs);
2407 double airm = eris_ifu_get_airmass_mean(phead);
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;
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);
2425 double cdelt3 = cpl_propertylist_get_double(head_wcs,
"CD3_3");
2426 double crval3 = cpl_propertylist_get_double(head_wcs,
"CRVAL3");
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);
2433 double gain = ERIS_GAIN;
2434 cpl_msg_debug(cpl_func,
"gain=%g", gain);
2438 double factor = 1. / (gain * exptime * bin_size_scale_factor);
2440 sci_wmax = crval3+cdelt3*(size-1);
2441 cpl_msg_warning(cpl_func,
"Sci table wmin: %f wmax: %f",sci_wmin, sci_wmax);
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;
2456 resp_tab = cpl_table_load(fname, 4, 0);
2458 resp_tab = cpl_table_load(fname, 1, 0);
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);
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);
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);
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();
2487 for(cpl_size k=0; k< size; k++) {
2493 cpl_propertylist_append_int(head_wcs,
"NAXIS3",size);
2498 sci_wmax = crval3+cdelt3*(size-1);
2499 cpl_msg_warning(cpl_func,
"Sci table wmin: %f wmax: %f",sci_wmin, sci_wmax);
2501 hdrl_spectrum1D * resp_s1d = NULL;
2504 "wavelength",
"response_poly_error",
"response_poly_bpm", scale);
2507 "wavelength",
"response_smo_error",
"response_smo_bpm", scale);
2509 cpl_table_delete(resp_xtab);
2511 hdrl_spectrum1D * resp_s1d_e = extend_hdrl_spectrum(sci_wmin, sci_wmax, resp_s1d);
2516 hdrl_parameter * res_par =
2519 cpl_frame* spectra_frm = NULL;
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.");
2529 cpl_table_delete(resp_tab);
2532 cpl_propertylist_delete(head_wcs);
2534 cpl_propertylist_delete(phead);
2536 return cpl_error_get_code();
2539 fname = cpl_frame_get_filename(spectra_frm);
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);
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);
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);
2555 "WAVE", NULL, NULL, scale);
2556 cpl_table_delete(spectra_xtab);
2561 cpl_msg_debug(cpl_func,
"obs_e wmin: %f wmax: %f", sci_wmin, sci_wmax);
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);
2573 hdrl_spectrum1D * ext_s1d =
2575 "LAMBDA", NULL, NULL, scale);
2591 cpl_msg_warning(cpl_func,
"nx: %d sz: %lld",nx, sz);
2592 cpl_size kmax = (nx < sz) ? nx : sz;
2595 hdrl_value resp_fct;
2596 hdrl_value factor_value = {factor, 0};
2598 cpl_msg_warning(cpl_func,
"kmax: %lld",kmax);
2600 for (cpl_size k = 0; k < kmax; k++) {
2601 double exp = -0.4 * airm * pext[k];
2604 hdrl_value atm_ext = {pow(10., exp), 0.};
2608 resp_fct.data = resp;
2610 if(isfinite(resp) && resp != 0) {
2612 }
else if(resp == 0) {
2616 for(cpl_size j = 1; j <= sy; j++ ) {
2617 for(cpl_size i = 1; i <= sx; i++ ) {
2621 }
else if(!isfinite(resp)) {
2624 for(cpl_size j = 1; j <= sy; j++ ) {
2625 for(cpl_size i = 1; i <= sx; i++ ) {
2639 char* fcal_cube_ptag = cpl_sprintf(
"%s_FLUXCAL", procatg);
2640 char* cube_file_name = cpl_sprintf(
"%s_cube_fluxcal.fits", pipefile_prefix);
2643 cpl_propertylist_set_string(phead, FHDR_PRO_CATG, fcal_cube_ptag);
2645 hdrl_resample_result* res = cpl_calloc(1,
sizeof(hdrl_resample_result));
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);
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;
2670 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, mean_cube_pctg);
2673 plugin_id, phead2D,
"RADECSYS", ima_file_name,
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);
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);
2686 cpl_boolean gen_phase3 = cpl_parameter_get_bool(p);
2688 cpl_msg_info(cpl_func,
"Save flux calibrated COADDED OBJECT cube with phase3 data format");
2690 cube_file_name, parlist, sof, CPL_TRUE);
2692 cpl_msg_info(cpl_func,
"Save flux calibrated COADDED OBJECT cube");
2694 cube_file_name, parlist, sof, CPL_FALSE);
2696 cpl_free(param_name);
2698 cpl_msg_info(cpl_func,
"Save flux calibrated STD STAR cube");
2700 cube_file_name, parlist, sof, CPL_FALSE);
2703 cpl_free(cube_file_name);
2707 cpl_image_delete(cube_cmap);
2709 cpl_table_delete(qclog);
2712 cpl_free(fcal_cube_ptag);
2718 cpl_table_delete(atm_ext_tab);
2719 cpl_table_delete(resp_tab);
2720 eris_remove_table_normal_format();
2727 return cpl_error_get_code();
2732eris_flux_calibrate_cube(
const char* procatg,
const char* pipefile_prefix,
2733 const char* plugin_id,
const cpl_parameterlist* parlist,
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);
2741 cpl_ensure_code(parlist, CPL_ERROR_NULL_INPUT);
2742 cpl_ensure_code(sof, CPL_ERROR_NULL_INPUT);
2744 cpl_size sof_sz = cpl_frameset_get_size(sof);
2745 cpl_ensure_code(sof_sz >0, CPL_ERROR_DATA_NOT_FOUND);
2747 cpl_msg_info(cpl_func,
"Flux calibrate the cube %s", procatg);
2750 hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
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);
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");
2764 double centralLambda = 0;
2766 double centralpix = 0;
2767 double newcenter_x = 0;
2768 double newcenter_y = 0;
2770 cpl_propertylist* data_head = cpl_propertylist_load(fname, 1);
2772 eris_get_wcs_cube(data_head, ¢ralLambda, &dis, ¢ralpix, &newcenter_x,
2774 cpl_propertylist* head_wcs = eris_ifu_plist_extract_wcs(data_head);
2775 cpl_propertylist_delete(data_head);
2777 double exptime = dit * ndit;
2778 cpl_msg_debug(cpl_func,
"exptime=%g", exptime);
2779 double gain = ERIS_GAIN;
2780 cpl_msg_debug(cpl_func,
"gain=%g", gain);
2782 double airm = eris_ifu_get_airmass_mean(phead);
2783 cpl_msg_debug(cpl_func,
"airm=%g", airm);
2784 double um2nm = 1000.;
2786 char band[FILE_NAME_SZ];
2788 eris_get_band(coadd_frm,band);
2790 double bin_size = eris_get_dispersion(band);
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);
2799 hdrl_spectrum1D * resp_s1d =
2801 "wavelength",
"response_smo_error",
"response_smo_bpm", scale);
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);
2808 hdrl_spectrum1D * ext_s1d =
2810 "LAMBDA", NULL, NULL, scale);
2813 hdrl_parameter * res_par =
2842 cpl_size kmax = (nx < sz) ? nx : sz;
2844 hdrl_value resp_fct;
2846 for (cpl_size k = 0; k < kmax; k++) {
2847 double exp = -0.4 * airm * pext[k];
2850 hdrl_value atm_ext = {pow(10., exp), 0.};
2852 resp_fct.data = resp;
2853 if(isfinite(resp)) {
2857 for(cpl_size j = 1; j <= sy; j++ ) {
2858 for(cpl_size i = 1; i <= sx; i++ ) {
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);
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);
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;
2895 cpl_propertylist_update_string(phead2D, CPL_DFS_PRO_CATG, mean_cube_pctg);
2898 plugin_id, phead2D,
"RADECSYS", ima_file_name,
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);
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);
2911 cpl_boolean gen_phase3 = cpl_parameter_get_bool(p);
2913 cpl_msg_info(cpl_func,
"Save flux calibrated COADDED OBJECT cube with phase3 data format");
2915 cube_file_name, parlist, sof, CPL_TRUE);
2917 cpl_msg_info(cpl_func,
"Save flux calibrated COADDED OBJECT cube");
2919 cube_file_name, parlist, sof, CPL_FALSE);
2921 cpl_free(param_name);
2923 cpl_msg_info(cpl_func,
"Save flux calibrated STD STAR cube");
2925 cube_file_name, parlist, sof, CPL_FALSE);
2928 cpl_free(cube_file_name);
2932 cpl_image_delete(cube_cmap);
2934 cpl_table_delete(qclog);
2937 cpl_free(fcal_cube_ptag);
2941 cpl_table_delete(atm_ext_tab);
2942 cpl_table_delete(resp_tab);
2946 return cpl_error_get_code();
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.
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
double eris_pfits_get_cd33(const cpl_propertylist *plist)
find out the character string associated to the CDELT3 keyword
double eris_pfits_get_crpix1(const cpl_propertylist *plist)
find out the character string associated to the CRPIX1 keyword
double eris_pfits_get_crval3(const cpl_propertylist *plist)
find out the character string associated to the CVRVAL3 keyword
double eris_pfits_get_crpix3(const cpl_propertylist *plist)
find out the character string associated to the CRPIX3 keyword
cpl_error_code eris_check_error_code(const char *func_id)
handle CPL errors
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
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
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
cpl_size hdrl_image_get_size_x(const hdrl_image *self)
return size of X dimension of image
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
cpl_error_code hdrl_image_reject(hdrl_image *self, cpl_size xpos, cpl_size ypos)
mark pixel as bad
hdrl_image * hdrl_image_new(cpl_size nx, cpl_size ny)
create new zero filled hdrl image
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
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