40#include "irplib_framelist.h"
42#include "visir_spc_photom.h"
44#include "visir_parameter.h"
45#include "visir_utils.h"
46#include "visir_pfits.h"
47#include "visir_spc_distortion.h"
53static cpl_bivector * visir_spc_phot_model_from_cat(
const char *,
double,
55static int visir_spc_phot_model_rebin(
const cpl_bivector *, cpl_table *);
57static cpl_error_code visir_spc_sens_stat(
const visir_spc_config *,
58 const visir_spc_resol,
59 double *,
double *,
double *,
62static cpl_error_code visir_spc_phot_qc(cpl_propertylist *, cpl_propertylist *,
63 const char *,
double,
double,
double,
77#ifndef VISIR_WAVELEN_QC_MIN
78#define VISIR_WAVELEN_QC_MIN 8e-6
81#ifndef VISIR_WAVELEN_QC_MAX
82#define VISIR_WAVELEN_QC_MAX 13e-6
86#ifndef VISIR_WAVELEN_QC_IGNORE
87#define VISIR_WAVELEN_QC_IGNORE 3
111 const visir_spc_config * cfg,
112 const cpl_propertylist * plist,
113 const char * star_cat,
114 cpl_image ** pweight2d,
115 cpl_propertylist * qclist,
116 cpl_table * spc_table,
117 const visir_spc_resol resol,
118 const char * dit_key)
121 cpl_errorstate cleanstate = cpl_errorstate_get();
122 cpl_bivector * flux_model = NULL;
126 double sens_median = -1.0;
127 double sens_mean = -1.0;
128 double sens_stdev = -1.0;
129 double exptime, factor, exectime = -1.;
133 const char * star_name;
138 skip_if(irplib_framelist_contains(rawframes, dit_key,
139 CPL_TYPE_DOUBLE, CPL_FALSE, 0.0));
140 if(irplib_framelist_contains(rawframes, dit_key,
141 CPL_TYPE_DOUBLE, CPL_TRUE, 1e-5)) {
144 visir_error_reset(
"DIT differs by more than %g", 1e-5);
148 skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_RA,
149 CPL_TYPE_DOUBLE, CPL_FALSE, 0.0));
152 skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_DEC,
153 CPL_TYPE_DOUBLE, CPL_TRUE, 1.0));
155 skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_INT_CHOP_NCYCLES,
156 CPL_TYPE_INT, CPL_TRUE, 0.0));
158 skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_INT_NDIT,
159 CPL_TYPE_INT, CPL_TRUE, 0.0));
161 if (irplib_framelist_contains(rawframes, VISIR_PFITS_STRING_STARNAME,
162 CPL_TYPE_STRING, CPL_TRUE, 0.0)) {
163 visir_error_reset(
"Rawframe(s) missing standard star name");
166 if (cpl_propertylist_has(plist,
"ESO QC EXPTIME"))
167 exptime = cpl_propertylist_get_double(plist,
"ESO QC EXPTIME");
169 const int nnod = irplib_framelist_get_size(rawframes);
170 exptime = visir_utils_get_exptime(nnod, plist);
172 cpl_msg_info(cpl_func,
"Exposure time: %.3f", exptime);
173 if (cpl_propertylist_has(plist,
"ESO QC EXECTIME"))
174 exectime = cpl_propertylist_get_double(plist,
"ESO QC EXECTIME");
176 skip_if (exptime <= 0.0);
178 factor = sqrt(exptime)/6.0;
182 if (star_name == NULL) visir_error_reset(
"Could not get standard star name");
185 if ((flux_model=visir_spc_phot_model_from_cat(star_cat, ra, dec)) == NULL) {
186 cpl_msg_error(cpl_func,
"Cannot retrieve the flux model from the cat.");
191 if (visir_spc_phot_model_rebin(flux_model, spc_table) == -1) {
192 cpl_msg_error(cpl_func,
"Cannot rebin the flux model");
195 cpl_bivector_delete(flux_model);
201 cpl_table_duplicate_column(spc_table,
"SENSITIVITY", spc_table,
203 cpl_table_set_column_unit(spc_table,
"SENSITIVITY",
"mJy");
205 ext_array = cpl_table_get_data_double(spc_table,
"SENSITIVITY");
206 spc_array = cpl_table_get_data_double(spc_table,
"SPC_EXTRACTED");
207 err_array = cpl_table_get_data_double(spc_table,
"SPC_ERROR");
209 for (cpl_size i = 0; i < cpl_table_get_nrow(spc_table); i++) {
210 ext_array[i] *= factor * err_array[i];
211 if (ext_array[i] < 0 || spc_array[i] <= 0) {
212 cpl_msg_warning(cpl_func,
"Setting non-physical sensitivity in row %d "
213 "to 0: %g/%g", (
int)i+1, ext_array[i], spc_array[i]);
216 ext_array[i] /= spc_array[i];
220 skip_if(visir_spc_sens_stat(cfg, resol, &sens_median, &sens_mean,
221 &sens_stdev, spc_table));
224 skip_if (cpl_table_erase_column(spc_table,
"SPC_EMISSIVITY"));
226 skip_if (visir_spc_phot_qc(qclist, cfg->phu, star_name, exptime, exectime,
227 sens_median, sens_mean, sens_stdev));
231 visir_table_plot(
"",
"t 'Extracted spectrum' w lines",
"", spc_table,
232 "WLEN",
"SPC_EXTRACTED");
233 visir_table_plot(
"",
"t 'Extracted spectrum error' w lines",
"",
234 spc_table,
"WLEN",
"SPC_ERROR");
235 visir_table_plot(
"",
"t 'Standard star model' w lines",
"", spc_table,
236 "WLEN",
"STD_STAR_MODEL");
237 visir_table_plot(
"set grid;",
"t 'Sensitivity (mJy)' w lines",
"",
238 spc_table,
"WLEN",
"SENSITIVITY");
239 visir_image_plot(
"",
"t 'The weight map'",
"", *pweight2d);
243 cpl_bivector_delete(flux_model);
245 return cpl_error_get_code();
272 const irplib_framelist * rawframes,
273 const visir_spc_config * pconfig,
274 const char * star_cat,
275 const char * spc_cal_lines,
276 const char * spc_cal_qeff,
277 cpl_image ** pweight2d,
278 cpl_propertylist * qclist,
284 const visir_spc_resol resol,
285 const char * dit_key)
287 const cpl_propertylist * plist;
288 const cpl_frame * frm = NULL;
289 cpl_imagelist * hcycle = NULL;
290 cpl_image * imhcycle = NULL;
291 cpl_image * flipped = NULL;
292 cpl_image * comorder = NULL;
293 cpl_image * order = NULL;
294 cpl_table * spc_table = NULL;
295 visir_data_type dtype;
299 skip_if (rawframes == NULL);
302 plist = irplib_framelist_get_propertylist_const(rawframes, 0);
303 frm = irplib_framelist_get_const(rawframes, 0);
305 skip_if (visir_get_data_type(frm, plist, &dtype, NULL));
311 skip_if (pconfig->do_fixcombi &&
312 visir_spc_det_fix(combined, 1, CPL_TRUE,
321 imhcycle = cpl_imagelist_unset(hcycle, 0);
324 skip_if (visir_spc_det_fix(&imhcycle, 1, CPL_FALSE,
334 flipped = visir_spc_flip(*combined, wlen, resol, dtype, &rev);
336 cpl_image_delete(*combined);
341 flipped = visir_spc_flip(imhcycle, wlen, resol, dtype, NULL);
343 cpl_image_delete(imhcycle);
347 const int ncomb = visir_get_ncombine(rawframes);
350 int lcol = -1, rcol = -1;
351 skip_if(visir_spc_extract_order(&order, &comorder, &lcol, &rcol, *combined,
352 imhcycle, wlen, pconfig, do_ech,
353 visir_data_is_aqu(dtype)));
354 skip_if(visir_qc_append_background(qclist, rawframes,
358 int const ident = -1;
359 visir_apdefs *
const aps = visir_apdefs_new(1, ident,
'O', 0);
360 aps->limits[0] = (visir_aplimits){ lcol, rcol };
363 skip_if (visir_spc_extract_wcal(comorder, order, lcol, rcol, wlen, slitw,
364 temp, fwhm, resol, pconfig, spc_cal_lines,
365 spc_cal_qeff, visir_data_is_aqu(dtype), aps,
366 ncomb, rev, &spc_table, pweight2d, qclist));
370 qclist, spc_table, resol, dit_key);
372 cpl_imagelist_delete(hcycle);
378 cpl_imagelist_delete(hcycle);
380 cpl_image_delete(flipped);
381 cpl_image_delete(imhcycle);
382 cpl_image_delete(order);
383 cpl_image_delete(comorder);
385 if (cpl_error_get_code()) {
386 cpl_table_delete(spc_table);
388 cpl_image_delete(*pweight2d);
410static cpl_bivector * visir_spc_phot_model_from_cat(
const char * star_cat,
414 const double max_radius = VISIR_STAR_MAX_RADIUS;
415 const cpl_table * catalog = NULL;
418 const cpl_array * da;
419 cpl_vector * v_ra = NULL;
420 cpl_vector * v_dec = NULL;
422 cpl_bivector * model = NULL;
423 cpl_vector * model_x = NULL;
424 cpl_vector * model_y = NULL;
425 const double conv_mJy = 3.33E8;
433 catalog = cpl_table_load(star_cat, 1, 0);
434 if (catalog == NULL) {
435 cpl_msg_error(cpl_func,
"Could not open the star catalog: %s",
436 star_cat ? star_cat :
"<NULL>");
440 nb_stars = cpl_table_get_nrow(catalog);
442 skip_if(nb_stars < 1);
445 dd = cpl_table_get_data_double_const(catalog,
"RA");
446 skip_if (dd == NULL);
447 IRPLIB_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual);
448 v_ra = cpl_vector_wrap(nb_stars, (
double*)dd);
449 IRPLIB_DIAG_PRAGMA_POP;
450 bug_if( v_ra == NULL);
452 dd = cpl_table_get_data_double_const(catalog,
"DEC");
453 skip_if (dd == NULL);
454 IRPLIB_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual);
455 v_dec = cpl_vector_wrap(nb_stars, (
double*)dd);
456 IRPLIB_DIAG_PRAGMA_POP;
457 bug_if( v_dec == NULL);
460 min_dist_ind = visir_star_find(v_ra, v_dec, ra, dec, max_radius, &star_dist);
462 skip_if (min_dist_ind < 0);
464 cpl_msg_info(cpl_func,
"The standard star closest to (RA,DEC)=(%g,%g) is "
465 "no. %d, '%s' at (RA,DEC)=(%g,%g) with the distance [degree]: "
466 "%g", ra, dec, 1+min_dist_ind,
467 cpl_table_get_string(catalog,
"STARS", min_dist_ind),
468 cpl_table_get_double(catalog,
"RA", min_dist_ind, NULL),
469 cpl_table_get_double(catalog,
"DEC", min_dist_ind, NULL),
473 da = cpl_table_get_array(catalog,
"WAVELENGTHS", min_dist_ind);
474 skip_if (da == NULL);
475 dd = cpl_array_get_data_double_const(da);
476 skip_if (dd == NULL);
478 nb_vals = cpl_array_get_size(da);
480 model_x = cpl_vector_new(nb_vals);
481 memcpy(cpl_vector_get_data(model_x), dd, nb_vals *
sizeof(
double));
484 da = cpl_table_get_array(catalog,
"MODEL_FLUX", min_dist_ind);
485 skip_if (da == NULL);
486 dd = cpl_array_get_data_double_const(da);
487 skip_if (dd == NULL);
489 skip_if (nb_vals != cpl_array_get_size(da));
491 model_y = cpl_vector_new(nb_vals);
492 memcpy(cpl_vector_get_data(model_y), dd, nb_vals *
sizeof(
double));
495 bug_if (cpl_vector_multiply_scalar(model_y, conv_mJy));
496 bug_if (cpl_vector_multiply(model_y, model_x));
497 bug_if (cpl_vector_multiply(model_y, model_x));
500 bug_if (cpl_vector_multiply_scalar(model_x, 1e-6));
502 model = cpl_bivector_wrap_vectors(model_x, model_y);
504 bug_if (model == NULL);
508 (void)cpl_vector_unwrap(v_ra);
509 (void)cpl_vector_unwrap(v_dec);
510 IRPLIB_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual);
511 cpl_table_delete((cpl_table*)catalog);
512 IRPLIB_DIAG_PRAGMA_POP;
514 if (cpl_error_get_code()) {
516 cpl_vector_delete(model_x);
517 cpl_vector_delete(model_y);
530static int visir_spc_phot_model_rebin(
531 const cpl_bivector * model,
532 cpl_table * spec_tab)
534 cpl_vector * bounds = NULL;
535 cpl_vector * spec = NULL;
537 const int nrow = cpl_table_get_nrow(spec_tab);
544 bounds = cpl_vector_new(nrow + 1);
545 for (i=1 ; i<cpl_vector_get_size(bounds) - 1 ; i++) {
546 bin_pos = (cpl_table_get(spec_tab,
"WLEN", i-1, NULL) +
547 cpl_table_get(spec_tab,
"WLEN", i, NULL)) / 2.0;
548 cpl_vector_set(bounds, i, bin_pos);
550 bin_pos = cpl_table_get(spec_tab,
"WLEN", 0, NULL) -
551 ((cpl_table_get(spec_tab,
"WLEN", 1, NULL) -
552 cpl_table_get(spec_tab,
"WLEN", 0, NULL)) / 2.0);
553 cpl_vector_set(bounds, 0, bin_pos);
555 cpl_table_get(spec_tab,
"WLEN", nrow-1, NULL) +
556 ((cpl_table_get(spec_tab,
"WLEN", 1, NULL) -
557 cpl_table_get(spec_tab,
"WLEN", 0, NULL)) / 2.0);
558 cpl_vector_set(bounds, cpl_vector_get_size(bounds)-1, bin_pos);
561 spec = cpl_vector_new(nrow);
564 if (visir_vector_resample(spec, bounds, model)) {
565 cpl_msg_error(cpl_func,
"Cannot rebin the spectrum");
570 cpl_table_new_column(spec_tab,
"STD_STAR_MODEL", CPL_TYPE_DOUBLE);
571 cpl_table_set_column_unit(spec_tab,
"STD_STAR_MODEL",
"mJy");
572 for (i=0 ; i<nrow ; i++) {
573 cpl_table_set_double(spec_tab,
"STD_STAR_MODEL", i,
574 cpl_vector_get(spec, i));
579 cpl_vector_delete(spec);
580 cpl_vector_delete(bounds);
582 return cpl_error_get_code();
596static cpl_error_code visir_spc_wl_start_end(
const cpl_table * self,
597 cpl_size * pstart, cpl_size * pend)
599 const cpl_size npix = cpl_table_get_nrow(self);
600 cpl_vector * vwave = NULL;
604 IRPLIB_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual);
605 vwave = cpl_vector_wrap(npix, (
double*)
606 cpl_table_get_data_double_const(self,
"WLEN"));
607 IRPLIB_DIAG_PRAGMA_POP;
609 assert(VISIR_WAVELEN_QC_MIN < VISIR_WAVELEN_QC_MAX);
611 cpl_size l = visir_lower_bound(vwave, VISIR_WAVELEN_QC_MIN);
612 cpl_size u = visir_upper_bound(vwave, VISIR_WAVELEN_QC_MAX);
614 l = CPL_MAX(l, VISIR_WAVELEN_QC_IGNORE);
615 u = CPL_MIN(u, npix - VISIR_WAVELEN_QC_IGNORE);
618 (void)cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
619 "For QC the calibrated wavelength range "
620 "(%g, %g) does not intersect with the "
621 "required (%g, %g) [m]",
622 cpl_vector_get(vwave, l),
623 cpl_vector_get(vwave, u),
624 VISIR_WAVELEN_QC_MIN, VISIR_WAVELEN_QC_MAX);
630 cpl_vector_unwrap(vwave);
632 return cpl_error_get_code();
648static cpl_error_code visir_spc_sens_stat(
const visir_spc_config * pconfig,
649 const visir_spc_resol resol,
650 double * psens_median,
652 double * psens_stdev,
653 const cpl_table * spc_table)
656 cpl_vector * sens_ignore = NULL;
657 cpl_vector * sens_mini = NULL;
658 double emis_min, emis_max;
659 cpl_size npix_start = 0;
660 cpl_size npix_end = 0;
666 skip_if(pconfig == NULL);
668 if (resol == VISIR_SPC_R_LRP || resol == VISIR_SPC_R_LR) {
669 skip_if(visir_spc_wl_start_end(spc_table, &npix_start, &npix_end));
671 npix_start = VISIR_WAVELEN_QC_IGNORE;
672 npix_end = cpl_table_get_nrow(spc_table) - VISIR_WAVELEN_QC_IGNORE;
675 nuse = npix_end - npix_start;
677 cpl_msg_info(cpl_func,
"QC computed from %lld spectral values in the range "
678 "(%g, %g) [m]", (
long long)nuse, VISIR_WAVELEN_QC_MIN,
679 VISIR_WAVELEN_QC_MAX);
685 IRPLIB_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual);
686 sens_ignore = cpl_vector_wrap( nuse, (
double*)
687 cpl_table_get_data_double_const(spc_table,
688 "SENSITIVITY") + npix_start);
689 IRPLIB_DIAG_PRAGMA_POP;
693 visir_vector_plot(
"set grid;",
"t 'Truncated Sensitivity (mJy)' w lines",
696 sens_mini = cpl_vector_duplicate(sens_ignore);
699 emis_min = cpl_table_get_column_min(spc_table,
"SPC_EMISSIVITY");
700 emis_max = cpl_table_get_column_max(spc_table,
"SPC_EMISSIVITY");
703 emis_max = emis_min + pconfig->phot_emis_tol * (emis_max - emis_min);
706 for (cpl_size i=0; i < nuse; i++) {
707 const double emis = cpl_table_get(spc_table,
"SPC_EMISSIVITY",
708 i + npix_start, NULL);
712 if (emis > emis_max)
continue;
715 skip_if(cpl_vector_set(sens_mini, iok, cpl_vector_get(sens_mini, i)));
723 skip_if(cpl_vector_set_size(sens_mini, iok));
725 *psens_mean = cpl_vector_get_mean(sens_mini);
728 cpl_msg_warning(cpl_func,
"Sensitivity computed on only 1 wavelength "
729 "with emissivity %f", emis_max);
732 cpl_msg_info(cpl_func,
"Sensitivity computed on %d wavelengths with "
733 "emissivity at most %f", iok, emis_max);
735 *psens_stdev = cpl_vector_get_stdev(sens_mini);
739 *psens_median = cpl_vector_get_median(sens_mini);
743 cpl_vector_unwrap(sens_ignore);
744 cpl_vector_delete(sens_mini);
746 return cpl_error_get_code();
766static cpl_error_code visir_spc_phot_qc(cpl_propertylist * qclist,
767 cpl_propertylist * phu_qclist,
768 const char * star_name,
776 bug_if (cpl_propertylist_append_double(phu_qclist,
"ESO QC EXPTIME", exptime));
778 bug_if (cpl_propertylist_append_double(phu_qclist,
"ESO QC EXECTIME",
781 bug_if (cpl_propertylist_append_double(qclist,
"ESO QC SENS MEDIAN",
783 bug_if (cpl_propertylist_append_double(qclist,
"ESO QC SENS MEAN", sens_mean));
784 bug_if (cpl_propertylist_append_double(qclist,
"ESO QC SENS STDEV",
787 bug_if (cpl_propertylist_append_string(phu_qclist,
"ESO QC STARNAME",
788 star_name ? star_name :
""));
792 return cpl_error_get_code();
double visir_pfits_get_ra(const cpl_propertylist *self)
The RA.
double visir_pfits_get_dec(const cpl_propertylist *self)
The DEC.
const char * visir_pfits_get_starname(const cpl_propertylist *self)
The std star name.
cpl_error_code visir_spc_phot_sensit(const irplib_framelist *rawframes, const visir_spc_config *cfg, const cpl_propertylist *plist, const char *star_cat, cpl_image **pweight2d, cpl_propertylist *qclist, cpl_table *spc_table, const visir_spc_resol resol, const char *dit_key)
Compute the sensitivity from an extracted spectrum.
cpl_table * visir_spc_phot_sensit_from_image(cpl_image **combined, const irplib_framelist *rawframes, const visir_spc_config *pconfig, const char *star_cat, const char *spc_cal_lines, const char *spc_cal_qeff, cpl_image **pweight2d, cpl_propertylist *qclist, cpl_boolean do_ech, double wlen, double slitw, double temp, double fwhm, const visir_spc_resol resol, const char *dit_key)
Extract spectrum from image and compute sensitivity.