36#include "hdrl_spectrum_resample.h"
37#include "irplib_wavecal.h"
38#include "irplib_framelist.h"
40#define SBRM_UNDERSCORE_MACRO
43#include "visir_spectro.h"
44#include "visir_utils.h"
45#include "visir_pfits.h"
46#include "visir_inputs.h"
47#include "visir_parameter.h"
48#include "visir_spc_distortion.h"
58#include <gsl/gsl_fit.h>
68#define RESET SBRM_RESET
70#define ABORT SBRM_ABORT
71#define CLEANUP SBRM_CLEANUP
77#define HOLE(type) &(type){0}
79#define skip_if_error_present() skip_if(0)
81#define MSG_WARN(...) cpl_msg_warning(cpl_func, __VA_ARGS__)
82#define MSG_INFO(...) cpl_msg_info(cpl_func, __VA_ARGS__)
83#define MSG_ERR(...) cpl_msg_error(cpl_func, __VA_ARGS__)
84#define MSG_DBG(...) cpl_msg_debug(cpl_func, __VA_ARGS__)
104 const cpl_vector * vsymm;
108 const cpl_bivector * lines;
110 const cpl_bivector * tqeff;
112} visir_spectrum_model;
114typedef cpl_bivector * (extract_func)(cpl_image *,
int,
int, cpl_propertylist *,
115 cpl_image **,
const visir_spc_config *,
116 const visir_apdefs *,
const bool,
124visir_polynomial_shift_1d_from_correlation(cpl_polynomial *,
126 irplib_base_spectrum_model *,
129 const cpl_polynomial *,
130 irplib_base_spectrum_model *),
131 int,
int, cpl_boolean,
134static cpl_error_code visir_spectro_refine(cpl_polynomial *,
136 visir_spectrum_model *,
137 const cpl_polynomial *,
138 int, cpl_boolean, visir_spc_resol,
139 double *, cpl_boolean *,
double *);
141static cpl_error_code visir_spectro_fill(cpl_vector *,
const cpl_polynomial *,
142 irplib_base_spectrum_model *);
144static extract_func visir_spc_oldex;
145static extract_func visir_spc_newex;
146static extract_func visir_spc_extract;
148static cpl_error_code visir_spc_emission(cpl_bivector *,
const cpl_vector *,
149 const cpl_bivector *,
150 const cpl_bivector *,
151 const cpl_vector *,
double);
153static cpl_polynomial * visir_spc_phys_disp(
int,
double, visir_spc_resol,
int,
155static cpl_polynomial * visir_spc_phys_lrp(
void);
156static double visir_spc_get_dispersion(
const cpl_polynomial *,
double);
157static cpl_error_code visir_vector_convolve_symm(cpl_vector *,
159static cpl_vector * cpl_spc_convolve_init(
int,
double,
double,
int);
161static cpl_error_code visir_spectro_qclist_wcal(cpl_propertylist *,
164 const cpl_polynomial *,
165 const cpl_polynomial *);
167static void * visir_spectro_qclist_obs(cpl_propertylist *,
double,
double);
169static const double N_upper = 13.4e-6;
170static const double whechelle = 35.8/2;
173#define VISIR_XC_LEN 50
175#ifndef VISIR_XC_SUBSEARCH
176#define VISIR_XC_SUBSEARCH 100
179#ifndef VISIR_SPECTRO_SIGMA
180#define VISIR_SPECTRO_SIGMA 3.0
197static const char * pn(
const int oo)
200 const char * sign = oo ? (oo > 0 ?
"+" :
"-") :
"";
201 snprintf(buf,
sizeof(buf),
"%s%d", sign, abs(oo));
205cpl_error_code visir_spc_extract_order(cpl_image ** order,
206 cpl_image ** comorder,
207 int * lcol,
int * rcol,
208 const cpl_image * combined,
209 const cpl_image * imhcycle,
211 const visir_spc_config * cfg,
212 const cpl_boolean do_ech,
219 VISIR_PARAM_REJLEFT);
221 VISIR_PARAM_REJRIGHT);
224 cpl_msg_debug(cpl_func,
"extracting order, wlen=%f, do_ech=%d, jcol1=%d, "
225 "jcol2=%d", wlen, do_ech, jcol1, jcol2);
228 skip_if (visir_spc_echelle_limit(&icol1, &icol2, wlen, cfg, 1,
229 cpl_image_get_size_y(combined),
233 icol2 = cpl_image_get_size_x(imhcycle);
238 cpl_msg_info(cpl_func,
"Ignoring %d leftmost columns from %d to %d",
239 jcol1, icol1, icol1 + jcol1);
243 cpl_msg_info(cpl_func,
"Ignoring %d rightmost columns from %d to %d",
244 jcol2, icol2 - jcol2, icol2);
249 cpl_msg_info(cpl_func,
"Ignoring %d leftmost columns", jcol1);
253 cpl_msg_info(cpl_func,
"Ignoring %d rightmost columns", jcol2);
258 if (icol1 != 1 || icol2 != cpl_image_get_size_x(imhcycle)) {
259 *order = visir_spc_column_extract(imhcycle, icol1, icol2, cfg->plot);
260 skip_if_error_present();
262 *comorder = visir_spc_column_extract(combined, icol1, icol2, cfg->plot);
263 skip_if_error_present();
266 *order = cpl_image_duplicate(imhcycle);
267 *comorder = cpl_image_duplicate(combined);
275 return cpl_error_get_code();
295visir_spc_resol visir_spc_get_res_wl(
const irplib_framelist * rawframes,
296 double * pwlen,
double * pslitw,
297 double * ptemp,
double * pfwhm,
300 cpl_errorstate cleanstate = cpl_errorstate_get();
302 visir_spc_resol resol = VISIR_SPC_R_ERR;
303 char ptmp[IRPLIB_FITS_STRLEN+1];
304 double wl, spx, pfov = 0.127;
306 cpl_boolean need_temp = ptemp != NULL;
310 cpl_ensure(rawframes != NULL, CPL_ERROR_NULL_INPUT, VISIR_SPC_R_ERR);
311 cpl_ensure(pwlen != NULL, CPL_ERROR_NULL_INPUT, VISIR_SPC_R_ERR);
312 cpl_ensure(pslitw != NULL, CPL_ERROR_NULL_INPUT, VISIR_SPC_R_ERR);
313 cpl_ensure(pfwhm != NULL, CPL_ERROR_NULL_INPUT, VISIR_SPC_R_ERR);
315 n = irplib_framelist_get_size(rawframes);
317 cpl_ensure(n > 0, CPL_ERROR_DATA_NOT_FOUND, VISIR_SPC_R_ERR);
320 skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_PIXSPACE,
321 CPL_TYPE_DOUBLE, CPL_TRUE, 1e-6));
324 skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_DOUBLE_SLITWIDTH,
325 CPL_TYPE_DOUBLE, CPL_FALSE, 0.0));
327 skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_STRING_RESOL,
328 CPL_TYPE_STRING, CPL_TRUE, 0.0));
330 skip_if(irplib_framelist_contains(rawframes, VISIR_PFITS_STRING_SLITNAME,
331 CPL_TYPE_STRING, CPL_TRUE, 0.0));
333 for (
int i = 0; i < n; i++) {
334 const cpl_propertylist * plist;
335 const char * filename =
336 cpl_frame_get_filename(irplib_framelist_get_const(rawframes, i));
338 double wl_tmp, sl_tmp, spx_tmp, pfov_tmp;
341 cpl_ensure(!cpl_error_get_code(), CPL_ERROR_DATA_NOT_FOUND,
344 cpl_ensure(filename != NULL, CPL_ERROR_DATA_NOT_FOUND,
347 plist = irplib_framelist_get_propertylist_const(rawframes, i);
349 cpl_ensure(plist != NULL, CPL_ERROR_DATA_NOT_FOUND, VISIR_SPC_R_ERR);
352 if (wl_tmp <= 0.0 || !cpl_errorstate_is_equal(cleanstate)) {
353 irplib_error_recover(cleanstate,
"Missing or invalid FITS card");
354 wl_tmp = VISIR_SPC_LRP_CWLEN;
357 if (pfits == NULL || !cpl_errorstate_is_equal(cleanstate)) {
358 irplib_error_recover(cleanstate,
"Missing or invalid FITS card");
359 pfits = VISIR_SPC_LRP_NAME;
367 if (pfov_tmp <= 0.) {
368 cpl_errorstate_set(cleanstate);
369 cpl_msg_warning(cpl_func, VISIR_PFITS_STRING_PIXSCALE
370 " not set, falling back to 0.127");
375 cpl_ensure(!cpl_error_get_code(), CPL_ERROR_DATA_NOT_FOUND,
380 visir_optmod ins_settings;
393 strncpy(ptmp, pfits, IRPLIB_FITS_STRLEN);
394 ptmp[IRPLIB_FITS_STRLEN] =
'\0';
396 cpl_msg_info(cpl_func,
"RESOL [" VISIR_SPC_LRP_NAME
"|LR|MR|HRS|HRG]"
397 " and WLEN [m] (%d frames): %s %g", n, ptmp, *pwlen);
400 cpl_msg_error(cpl_func,
"Pixel Spacing (%g) in %s is non-"
401 "positive", spx, filename);
402 cpl_ensure(0, CPL_ERROR_ILLEGAL_INPUT, VISIR_SPC_R_ERR);
406 cpl_msg_error(cpl_func,
"Slit Width (%g) in %s is non-positive",
408 cpl_ensure(0, CPL_ERROR_ILLEGAL_INPUT, VISIR_SPC_R_ERR);
411 cpl_msg_info(cpl_func,
"Slit Width [pixel] and Pixel Spacing [m]: "
412 "%g %g", *pslitw, spx);
414 if (!strcmp(VISIR_SPC_LRP_NAME, ptmp)) {
415 resol = VISIR_SPC_R_LRP;
416 }
else if (!strcmp(
"LR", ptmp)) {
417 resol = VISIR_SPC_R_LR;
418 }
else if (!strcmp(
"MR", ptmp)) {
419 resol = VISIR_SPC_R_MR;
420 }
else if (!strcmp(
"HRS", ptmp)) {
421 resol = VISIR_SPC_R_HR;
422 }
else if (!strcmp(
"HRG", ptmp)) {
423 resol = VISIR_SPC_R_GHR;
425 cpl_msg_error(cpl_func,
"Unsupported resolution (%s) in %s",
427 cpl_ensure(0, CPL_ERROR_UNSUPPORTED_MODE, VISIR_SPC_R_ERR);
430 if (resol != VISIR_SPC_R_LRP) {
432 skip_if(irplib_framelist_contains(rawframes,
433 VISIR_PFITS_DOUBLE_WLEN,
434 CPL_TYPE_DOUBLE, CPL_TRUE,
438 if (visir_spc_optmod_init(resol, *pwlen, &ins_settings, is_aqu)) {
439 cpl_msg_error(cpl_func,
"Resolution %s does not support "
440 "Central Wavelength [m]: %g", ptmp, *pwlen);
441 cpl_ensure(0, CPL_ERROR_INCOMPATIBLE_INPUT, VISIR_SPC_R_ERR);
444 cpl_msg_info(cpl_func,
"The %s-Spectral Resolution at %gm: %g",
446 visir_spc_optmod_resolution(&ins_settings));
447 cpl_msg_info(cpl_func,
"The %s-Linear Dispersion at %gm [pixel/m]: "
449 visir_spc_optmod_dispersion(&ins_settings));
451 *pfwhm = *pwlen * visir_spc_optmod_dispersion(&ins_settings)
452 / visir_spc_optmod_resolution(&ins_settings);
454 cpl_msg_info(cpl_func,
"The %s-FWHM at %gm [pixel]: %g",
455 ptmp, *pwlen, *pfwhm);
457 if (fabs(sl-sl_tmp) > 1e-3) {
458 cpl_msg_error(cpl_func,
"Inconsistent slit width (%g <=>"
459 " %g) in %s (%d of %d)",
460 sl, sl_tmp, filename, i+1, n);
461 cpl_ensure(0, CPL_ERROR_INCOMPATIBLE_INPUT, VISIR_SPC_R_ERR);
463 if (fabs(pfov-pfov_tmp) > 1e-4) {
464 cpl_msg_error(cpl_func,
"Inconsistent pfov (%g <=>"
465 " %g) in %s (%d of %d)",
466 pfov, pfov_tmp, filename, i+1, n);
467 cpl_ensure(0, CPL_ERROR_INCOMPATIBLE_INPUT, VISIR_SPC_R_ERR);
473 if (cpl_error_get_code()) {
474 visir_error_reset(
"Could not get FITS key");
475 }
else if ((-20 < temp) && (temp < 60)) {
477 need_temp = CPL_FALSE;
485 cpl_msg_warning(cpl_func,
"No FITS-files specify the M1 temperature, "
493 cpl_msg_info(cpl_func,
"The M1 temperature [Kelvin]: %g", *ptemp);
524cpl_error_code visir_vector_resample(cpl_vector * self,
525 const cpl_vector * xbounds,
526 const cpl_bivector * source)
529 const cpl_vector * xsource = cpl_bivector_get_x_const(source);
530 const cpl_vector * ysource = cpl_bivector_get_y_const(source);
532 const double * pxsource = cpl_vector_get_data_const(xsource);
533 const double * pysource = cpl_vector_get_data_const(ysource);
534 const double * pxbounds = cpl_vector_get_data_const(xbounds);
537 cpl_vector * ybounds = cpl_vector_new(cpl_vector_get_size(xbounds));
538 IRPLIB_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual)
539 cpl_bivector * boundary = cpl_bivector_wrap_vectors((cpl_vector*)xbounds,
541 IRPLIB_DIAG_PRAGMA_POP
542 double * pybounds = cpl_vector_get_data(ybounds);
544 double * pself = cpl_vector_get_data(self);
545 const int npix = cpl_vector_get_size(self);
550 cpl_ensure_code(cpl_bivector_get_size(boundary) == npix + 1,
551 CPL_ERROR_ILLEGAL_INPUT);
553 skip_if_error_present();
555 itt = cpl_vector_find(xsource, pxbounds[0]);
557 skip_if_error_present();
559 skip_if (cpl_bivector_interpolate_linear(boundary, source));
563 while (pxsource[itt] < pxbounds[0]) itt++;
565 for (i=0; i < npix; i++) {
570 double xlow = pxbounds[i];
571 double x = pxsource[itt];
573 if (x > pxbounds[i+1]) x = pxbounds[i+1];
575 pself[i] = pybounds[i] * (x - xlow);
578 while (pxsource[itt] < pxbounds[i+1]) {
579 const double xprev = x;
581 if (x > pxbounds[i+1]) x = pxbounds[i+1];
582 pself[i] += pysource[itt] * (x - xlow);
588 pself[i] += pybounds[i+1] * (pxbounds[i+1] - xlow);
592 pself[i] /= 2 * (pxbounds[i+1] - pxbounds[i]);
599 cpl_vector_delete(ybounds);
600 cpl_bivector_unwrap_vectors(boundary);
602 return cpl_error_get_code();
634void * visir_spc_extract_wcal(
const cpl_image * combined,
635 const cpl_image * hcycle,
636 const int lcol,
const int rcol,
637 const double wlen,
const double slitw,
638 const double temp,
const double fwhm,
639 const visir_spc_resol resol,
640 const visir_spc_config * cfg,
641 const char * spc_cal_lines,
642 const char * spc_cal_qeff,
644 const visir_apdefs * aps,
645 const cpl_size ncomb,
const bool rev,
646 cpl_table ** pspc_table,
647 cpl_image ** pweight2d,
648 cpl_propertylist * qclist)
651 const int npix = cpl_image_get_size_y(combined);
653 if (!pspc_table)
return ABORT(CPL_ERROR_NULL_INPUT);
654 if (!pweight2d)
return ABORT(CPL_ERROR_NULL_INPUT);
657 if (npix < 1)
return ABORT(CPL_ERROR_ILLEGAL_INPUT);
658 if (npix != cpl_image_get_size_y(hcycle))
return ABORT(
659 CPL_ERROR_ILLEGAL_INPUT,
"Sky frame does not have same size as the "
660 "object frame. %d vs %d pixels", (
int)cpl_image_get_size_y(hcycle),
665 visir_spc_wavecal _(hcycle, qclist, wlen, slitw, temp, fwhm, resol,
666 cfg, spc_cal_lines, spc_cal_qeff, pspc_table,
670 SET(flipped, cpl_image) = cpl_image_cast _(combined, CPL_TYPE_DOUBLE);
674 SET(spc_n_err, cpl_bivector) = visir_spc_extract _(
675 flipped->o, lcol, rcol, qclist, pweight2d, cfg, aps, rev, ncomb);
679 if (visir_str_par_is_empty(cfg->respcal)) {
680 cpl_table_new_column _(*pspc_table,
"SPC_EXTRACTED", CPL_TYPE_DOUBLE);
681 cpl_table_new_column _(*pspc_table,
"SPC_ERROR", CPL_TYPE_DOUBLE);
683 cpl_table_copy_data_double _(*pspc_table,
"SPC_EXTRACTED",
684 cpl_bivector_get_x_data(spc_n_err->o));
685 cpl_table_copy_data_double _(*pspc_table,
"SPC_ERROR",
686 cpl_bivector_get_y_data(spc_n_err->o));
689 MSG_INFO(
"Applying response calibration...");
692 SET(calib, cpl_table) = cpl_table_load _(cfg->respcal, 1, 0);
693 cpl_size nrows = cpl_table_get_nrow _(calib->o);
694 double * pwlen = cpl_table_get_data_double _(calib->o,
"WLEN");
695 double * pflux = cpl_table_get_data_double _(calib->o,
"FLUX");
696 double * perr = cpl_table_get_data_double(calib->o,
"ERR");
697 if (cpl_error_get_code()) WARN(
"ERR column missing from %s, continuing "
698 "with zero error", cfg->respcal);
699 SET(wlen, cpl_array, p) = cpl_array_wrap_double _(pwlen, nrows);
700 cpl_array_multiply_scalar _(wlen->o, 1e-6);
701 SET(flux, cpl_image, p) = cpl_image_wrap_double _(nrows, 1, pflux);
702 SET(err, cpl_image, p) = !perr ? NULL : cpl_image_wrap_double _(
704 SET(response, hdrl_spectrum1D, w) = hdrl_spectrum1D_create _(
705 flux->o, err->o, wlen->o, hdrl_spectrum1D_wave_scale_linear);
708 nrows = cpl_table_get_nrow _(*pspc_table);
709 pwlen = cpl_table_get_data_double _(*pspc_table,
"WLEN");
710 pflux = cpl_bivector_get_x_data _(spc_n_err->o);
711 perr = cpl_bivector_get_y_data _(spc_n_err->o);
712 RESET(wlen) = cpl_array_wrap_double _(pwlen, nrows);
713 RESET(flux) = cpl_image_wrap_double _(nrows, 1, pflux);
714 RESET(err) = cpl_image_wrap_double _(nrows, 1, perr);
715 SET(spectrum, hdrl_spectrum1D, w) = hdrl_spectrum1D_create _(
716 flux->o, err->o, wlen->o, hdrl_spectrum1D_wave_scale_linear);
719 const hdrl_spectrum1D_wavelength spec_wav =
720 hdrl_spectrum1D_get_wavelength _(spectrum->o);
721 SET(params, hdrl_parameter) =
722 hdrl_spectrum1D_resample_interpolate_parameter_create _(
723 hdrl_spectrum1D_interp_linear);
724 SET(result, hdrl_spectrum1D, w) = hdrl_spectrum1D_resample _(
725 response->o, &spec_wav, params->o);
728 hdrl_spectrum1D_div_spectrum _(spectrum->o, result->o);
729 hdrl_spectrum1D_append_to_table _(
730 spectrum->o, *pspc_table,
"SPC_EXTRACTED", NULL,
"SPC_ERROR", NULL);
733 cpl_table_set_column_unit _(*pspc_table,
"SPC_EXTRACTED",
"ADU/s");
734 cpl_table_set_column_unit _(*pspc_table,
"SPC_ERROR",
"ADU/s");
737 visir_table_plot(
"set grid;set xlabel 'Wavelength [m]';",
738 "t 'Extracted Spectrum' w linespoints",
739 "", *pspc_table,
"WLEN",
"SPC_EXTRACTED");
740 visir_table_plot(
"set grid;set xlabel 'Wavelength [m]';",
741 "t 'Error on Extracted Spectrum' w linespoints",
742 "", *pspc_table,
"WLEN",
"SPC_ERROR");
775cpl_error_code visir_spc_wavecal(
const cpl_image * hcycle,
776 cpl_propertylist * qclist,
777 double wlen,
double slitw,
778 double temp,
double fwhm,
779 visir_spc_resol resol,
780 const visir_spc_config * cfg,
781 const char * linefile,
782 const char * qefffile,
783 cpl_table ** pspc_table,
788 cpl_polynomial * phdisp = NULL;
790 cpl_polynomial * xcdisp = NULL;
792 visir_spectrum_model mymodel;
793 cpl_vector * wlvals = NULL;
794 cpl_vector * spmodel = NULL;
796 cpl_bivector * emission = NULL;
797 cpl_vector * boundary = NULL;
799 cpl_bivector * temiss = NULL;
800 cpl_bivector * tqeff = NULL;
802 cpl_image * corrected = NULL;
804 cpl_image * xc_image = NULL;
805 cpl_vector * xc_vector = NULL;
807 cpl_vector * vsymm = NULL;
808 cpl_vector * vxc = NULL;
810 const int npix = cpl_image_get_size_y(hcycle);
814 double qcxc = -1.0, qcsubdelta = 0.;
816 const cpl_size i0 = 0;
817 const cpl_size i1 = 1;
818 cpl_boolean didshift = CPL_FALSE;
821 cpl_ensure_code(!cpl_error_get_code(), cpl_error_get_code());
822 cpl_ensure_code(pspc_table, CPL_ERROR_NULL_INPUT);
823 cpl_ensure_code(npix > 0, CPL_ERROR_ILLEGAL_INPUT);
827 corrected = cpl_image_cast(hcycle, CPL_TYPE_DOUBLE);
828 skip_if_error_present();
830 hc_min = cpl_image_get_min(corrected);
831 skip_if_error_present();
832 cpl_msg_info(cpl_func,
"Half-cycle image [%d X %d] has minimum intensity: %g",
833 (
int)cpl_image_get_size_x(hcycle), npix, hc_min);
835 cpl_msg_warning(cpl_func,
"Thresholding negative intensities in half-"
836 "cycle image: %g", hc_min);
837 skip_if (cpl_image_threshold(corrected, 0.0, DBL_MAX, 0.0, DBL_MAX));
838 }
else if (hc_min > 0) {
839 skip_if (cpl_image_subtract_scalar(corrected, hc_min));
842 xc_image = cpl_image_duplicate(corrected);
845 cpl_image_delete(corrected);
846 corrected = cpl_image_collapse_create(xc_image, 1);
847 cpl_image_delete(xc_image);
848 xc_image = corrected;
851 skip_if(cpl_image_divide_scalar(xc_image, npix));
853 xc_vector = cpl_vector_wrap(npix, cpl_image_get_data(xc_image));
855 skip_if_error_present();
858 phdisp = visir_spc_phys_lrp();
859 cpl_msg_info(cpl_func,
"Central Dispersion (physical model) [pixel/m]: %g",
860 1.0/visir_spc_get_dispersion(phdisp, npix/2.0 + 0.5));
861 cpl_msg_info(cpl_func,
"Central Wavelength (physical model) [m]: %g",
862 cpl_polynomial_eval_1d(phdisp, npix/2.0 + 0.5, NULL));
863 cpl_msg_info(cpl_func,
"First Wavelength (physical model) [m]: %g",
864 cpl_polynomial_eval_1d(phdisp, 1.0, NULL));
865 cpl_msg_info(cpl_func,
"Last Wavelength (physical model) [m]: %g",
866 cpl_polynomial_eval_1d(phdisp, 1024, NULL));
867 cpl_polynomial_dump(phdisp, stdout);
868 cpl_polynomial_delete(phdisp);
871 phdisp = visir_spc_phys_disp(npix, wlen, resol, cfg->orderoffset, is_aqu);
872 skip_if_error_present();
874 if (cpl_polynomial_get_degree(phdisp) == 2) {
875 const cpl_size i2 = 2;
876 cpl_msg_info(cpl_func,
"Dispersion polynomial of physical model:"
877 " %gmum + ipix * %gmum/pixel + ipix^2 * (%g)mum/pixel^2 "
878 "[ipix = 1, 2, ..., %d]",
879 cpl_polynomial_get_coeff(phdisp, &i0) * 1e6,
880 cpl_polynomial_get_coeff(phdisp, &i1) * 1e6,
881 cpl_polynomial_get_coeff(phdisp, &i2) * 1e6,
885 cpl_msg_info(cpl_func,
"Dispersion polynomial of physical model:"
886 " %gmum + ipix * %gmum/pixel [ipix = 1, 2, ..., %d]",
887 cpl_polynomial_get_coeff(phdisp, &i0) * 1e6,
888 cpl_polynomial_get_coeff(phdisp, &i1) * 1e6, npix);
891 temiss = visir_bivector_load_fits(linefile,
"Wavelength",
"Emission", 1);
892 any_if (
"Could not load file with Emission Lines");
894 tqeff = visir_bivector_load_fits(qefffile,
"Wavelength",
"Efficiency",
896 any_if(
"Could not load file with Quantum-Efficiencies");
899 visir_bivector_plot(
"set grid;set xlabel 'Wavelength [m]';",
"t '"
900 "Quantum Efficiency' w linespoints",
"", tqeff);
903 vsymm = cpl_spc_convolve_init(npix, slitw, fwhm, cfg->plot);
905 skip_if (vsymm == NULL);
907 vxc = cpl_vector_new(1);
908 xcdisp = cpl_polynomial_new(1);
910 mymodel.lines = temiss;
911 mymodel.tqeff = tqeff;
912 mymodel.vsymm = vsymm;
918 skip_if(visir_spectro_refine(xcdisp, xc_vector, &mymodel, phdisp,
919 VISIR_XC_LEN, cfg->plot, resol,
920 &qcxc, &didshift, &qcsubdelta));
923 if (fabs(qcsubdelta) >= VISIR_XC_LEN) {
924 cpl_msg_warning(cpl_func,
"Cross-correlation (%g pixel shift): %g",
927 cpl_msg_info(cpl_func,
"Cross-correlation (%g pixel shift): %g",
932 cpl_msg_info(cpl_func,
"Dispersion polynomial from cross-correlation: "
933 "%gm + ipix * %gm/pixel [ipix = 1, 2, ..., %d]",
934 cpl_polynomial_get_coeff(xcdisp, &i0),
935 cpl_polynomial_get_coeff(xcdisp, &i1), npix);
937 cpl_msg_info(cpl_func,
"New Central Wavelength [m]: %g",
938 cpl_polynomial_eval_1d(xcdisp, 0.5*npix+0.5, NULL));
940 *pspc_table = cpl_table_new(npix);
941 skip_if_error_present();
944 wlvals = cpl_vector_new(npix);
945 bug_if (cpl_table_wrap_double(*pspc_table, cpl_vector_get_data(wlvals),
948 skip_if (cpl_vector_fill_polynomial(wlvals, xcdisp, 1.0, 1.0));
951 spmodel = cpl_vector_new(npix);
952 bug_if (cpl_table_wrap_double(*pspc_table, cpl_vector_get_data(spmodel),
954 skip_if (visir_spectro_fill(spmodel, phdisp,
955 (irplib_base_spectrum_model*)&mymodel));
958 (void)cpl_vector_unwrap(spmodel);
959 spmodel = cpl_vector_new(npix);
960 bug_if (cpl_table_wrap_double(*pspc_table, cpl_vector_get_data(spmodel),
963 skip_if (visir_spectro_fill(spmodel, xcdisp,
964 (irplib_base_spectrum_model*)&mymodel));
966 bug_if (cpl_table_wrap_double(*pspc_table,
967 cpl_image_get_data_double(xc_image),
969 (void)cpl_image_unwrap(xc_image);
973 (void)cpl_vector_unwrap(spmodel);
974 spmodel = cpl_vector_new(npix);
975 bug_if (cpl_table_wrap_double(*pspc_table, cpl_vector_get_data(spmodel),
978 boundary = cpl_vector_new(npix + 1);
979 skip_if (cpl_vector_fill_polynomial(boundary, xcdisp, 0.5, 1.0));
980 skip_if (visir_vector_resample(spmodel, boundary, temiss));
982 bug_if (cpl_table_set_column_unit(*pspc_table,
"WLEN",
"m"));
983 bug_if (cpl_table_set_column_unit(*pspc_table,
"SPC_MODEL_PH",
985 bug_if (cpl_table_set_column_unit(*pspc_table,
"SPC_MODEL_XC",
987 bug_if (cpl_table_set_column_unit(*pspc_table,
"SPC_SKY",
"ADU/s"));
993 if (resol != VISIR_SPC_R_LRP && cpl_vector_get(wlvals, 0) < N_upper &&
994 N_upper < cpl_vector_get(wlvals, npix-1))
995 cpl_msg_warning(cpl_func,
"Spectrum goes above N-band (%gm). Wavelength"
996 " Calibration may be entirely inaccurate", N_upper);
998 bug_if(visir_spectro_qclist_wcal(qclist, npix, qcxc, didshift, qcsubdelta,
1002 cpl_bivector * plot = cpl_bivector_wrap_vectors(wlvals, xc_vector);
1004 visir_bivector_plot(
"set grid;set xlabel 'Wavelength [m]';",
"t 'Spec"
1005 "trum from Half-cycle' w linespoints",
"", plot);
1006 cpl_bivector_unwrap_vectors(plot);
1008 visir_table_plot(
"set grid;set xlabel 'Wavelength [m]';",
1009 "t 'Calibrated Model Spectrum' w linespoints",
1010 "", *pspc_table,
"WLEN",
"SPC_MODEL_XC");
1013 visir_table_plot(
"set grid;set xlabel 'Wavelength [m]';",
1014 "t 'Physical Model Spectrum' w linespoints",
1015 "", *pspc_table,
"WLEN",
"SPC_MODEL_PH");
1017 if (resol != VISIR_SPC_R_LRP) {
1020 emission = cpl_bivector_new(2 * npix);
1022 cpl_vector_delete(boundary);
1023 boundary = cpl_vector_new(2 * npix + 1);
1025 cpl_vector_fill_polynomial(cpl_bivector_get_x(emission),
1026 phdisp, -0.5*npix, 1);
1027 cpl_vector_fill_polynomial(boundary, phdisp, -0.5*(npix+1), 1);
1030 visir_spc_emission(emission, boundary, temiss, tqeff, vsymm, temp);
1031 cpl_vector_delete(boundary);
1034 visir_bivector_plot(
"set grid;set xlabel 'Wavelength [m]';",
1035 "t 'Extended Model Spectrum' w linespoints",
1042 (void)cpl_vector_unwrap(wlvals);
1043 (void)cpl_vector_unwrap(spmodel);
1044 cpl_polynomial_delete(phdisp);
1045 cpl_polynomial_delete(xcdisp);
1046 cpl_image_delete(xc_image);
1047 cpl_vector_delete(vsymm);
1048 cpl_image_delete(corrected);
1049 cpl_bivector_delete(temiss);
1050 cpl_bivector_delete(tqeff);
1051 cpl_vector_delete(boundary);
1052 cpl_bivector_delete(emission);
1053 (void)cpl_vector_unwrap(xc_vector);
1054 cpl_vector_delete(vxc);
1056 return cpl_error_get_code();
1078cpl_error_code visir_spc_echelle_limit(
int * pcol1,
int * pcol2,
double wlen,
1079 const visir_spc_config * cfg,
1080 int icolmin,
int icolmax,
1084 visir_optmod ins_settings;
1091 cpl_ensure_code(wlen > 0, CPL_ERROR_ILLEGAL_INPUT);
1092 cpl_ensure_code(pcol1, CPL_ERROR_NULL_INPUT);
1093 cpl_ensure_code(pcol2, CPL_ERROR_NULL_INPUT);
1094 cpl_ensure_code(icolmin > 0, CPL_ERROR_ILLEGAL_INPUT);
1095 cpl_ensure_code(icolmax >= icolmin, CPL_ERROR_ILLEGAL_INPUT);
1097 cpl_ensure_code(cfg->orderoffset >= -4, CPL_ERROR_ILLEGAL_INPUT);
1098 cpl_ensure_code(cfg->orderoffset <= 4, CPL_ERROR_ILLEGAL_INPUT);
1100 error = visir_spc_optmod_init(VISIR_SPC_R_GHR, wlen, &ins_settings, is_aqu);
1102 MSG_ERR(
"HRG Optical model initialization (%p) failed: %d (%g)",
1103 (
void*)&ins_settings, error, wlen);
1104 cpl_ensure_code(0, CPL_ERROR_ILLEGAL_INPUT);
1106 order = cfg->orderoffset + visir_spc_optmod_get_echelle_order(&ins_settings);
1109 cpl_ensure_code(order > 0, CPL_ERROR_ILLEGAL_INPUT);
1110 cpl_ensure_code(order <= 18, CPL_ERROR_ILLEGAL_INPUT);
1112 wleni = visir_spc_optmod_echelle(&ins_settings, wlen, cfg->orderoffset );
1114 echpos = visir_spc_optmod_cross_dispersion(&ins_settings, wleni);
1115 if (echpos <= 0 || echpos >= icolmax) {
1116 MSG_ERR(
"Echelle order %2d: offset %s: location out of range [%d;%d]: "
1117 "%g", order, pn(cfg->orderoffset), icolmin, icolmax, echpos);
1118 cpl_ensure_code(0, CPL_ERROR_DATA_NOT_FOUND);
1121 *pcol1 = ceil(echpos - whechelle);
1122 *pcol2 = echpos + whechelle;
1124 if (*pcol1 < icolmin) *pcol1 = icolmin;
1125 if (*pcol2 > icolmax) *pcol2 = icolmax;
1127 MSG_INFO(
"Echelle order %2d: offset %s: at col %g [%d; %d]", order,
1128 pn(cfg->orderoffset), echpos, *pcol1, *pcol2);
1131 char * label = cpl_sprintf(
"ESO DRS APGUI OFFS%d", order);
1132 cpl_propertylist_update_int(cfg->phu, label, cfg->orderoffset);
1134 label = cpl_sprintf(
"ESO DRS APGUI WLEN%d", order);
1135 cpl_propertylist_update_double(cfg->phu, label, wleni);
1137 label = cpl_sprintf(
"ESO DRS APGUI CPIX%d", order);
1138 cpl_propertylist_update_double(cfg->phu, label, echpos);
1140 label = cpl_sprintf(
"ESO DRS APGUI LPIX%d", order);
1141 cpl_propertylist_update_int(cfg->phu, label, *pcol1);
1143 label = cpl_sprintf(
"ESO DRS APGUI RPIX%d", order);
1144 cpl_propertylist_update_int(cfg->phu, label, *pcol2);
1148 return cpl_error_get_code();
1166cpl_image * visir_spc_column_extract(
const cpl_image * self,
int icol1,
1167 int icol2,
int doplot)
1170 cpl_image * band = NULL;
1171 cpl_image * spatial = NULL;
1172 const int nrow = cpl_image_get_size_y(self);
1173 const int ncol = cpl_image_get_size_x(self);
1175 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
1176 cpl_ensure(icol1 > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1177 cpl_ensure(icol2 >= icol1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1179 cpl_ensure(ncol >= icol2, CPL_ERROR_ILLEGAL_INPUT, NULL);
1181 band = cpl_image_extract(self, icol1, 1, icol2, nrow);
1182 skip_if_error_present();
1185 visir_image_plot(
"",
"t 'The full-width image'",
"", self);
1189 spatial = cpl_image_collapse_create(self, 0);
1190 skip_if_error_present();
1191 skip_if (cpl_image_divide_scalar(spatial, nrow));
1193 visir_image_row_plot(
"set grid;",
"t 'Spectral direction "
1194 "collapsed' w linespoints",
"",
1201 cpl_image_delete(spatial);
1202 if (cpl_error_get_code() && band != NULL) {
1203 cpl_image_delete(band);
1226cpl_error_code visir_spectro_qc(cpl_propertylist * qclist,
1227 cpl_propertylist * paflist,
1228 cpl_boolean drop_wcs,
1229 const irplib_framelist * rawframes,
1230 const char * regcopy,
1231 const char * regcopypaf)
1234 const cpl_propertylist * reflist
1235 = irplib_framelist_get_propertylist_const(rawframes, 0);
1239 bug_if (visir_qc_append_capa(qclist, rawframes));
1241 if (regcopy != NULL)
1242 bug_if (cpl_propertylist_copy_property_regexp(qclist, reflist,
1245 if (regcopypaf != NULL)
1246 bug_if (cpl_propertylist_copy_property_regexp(paflist, reflist,
1249 bug_if (cpl_propertylist_append(paflist, qclist));
1252 cpl_propertylist * pcopy = cpl_propertylist_new();
1253 const cpl_error_code error
1254 = cpl_propertylist_copy_property_regexp(pcopy, reflist,
"^("
1255 IRPLIB_PFITS_WCS_REGEXP
1257 if (!error && cpl_propertylist_get_size(pcopy) > 0) {
1258 cpl_msg_warning(cpl_func,
"Combined image will have no WCS "
1261 cpl_propertylist_delete(pcopy);
1264 bug_if(cpl_propertylist_copy_property_regexp(qclist, reflist,
"^("
1265 IRPLIB_PFITS_WCS_REGEXP
1271 return cpl_error_get_code();
1293static cpl_error_code visir_spectro_qclist_wcal(cpl_propertylist * self,
1294 int npix,
double xc,
1295 cpl_boolean didshift,
1297 const cpl_polynomial * phdisp,
1298 const cpl_polynomial * xcdisp)
1301 const cpl_size phdegree = cpl_polynomial_get_degree(phdisp);
1302 const cpl_size xcdegree = cpl_polynomial_get_degree(xcdisp);
1304 const double phdisp0 = cpl_polynomial_eval_1d(phdisp, 1.0, NULL);
1305 const double xcdisp0 = cpl_polynomial_eval_1d(xcdisp, 1.0, NULL);
1307 const double xcwlen = cpl_polynomial_eval_1d(xcdisp, 0.5*(
double)npix+0.5,
1309 const double phcdisp = visir_spc_get_dispersion(phdisp, npix/2.0 + 0.5);
1310 const double xccdisp = visir_spc_get_dispersion(xcdisp, npix/2.0 + 0.5);
1315 skip_if (phdegree < 1);
1316 skip_if (xcdegree < 1);
1318 cpl_msg_info(cpl_func,
"Central Dispersion (physical model) [m/pixel]: %g",
1320 cpl_msg_info(cpl_func,
"Central Dispersion (calibrated) [m/pixel]: %g",
1323 bug_if (cpl_propertylist_append_double(self,
"ESO QC XC", xc));
1326 bug_if (cpl_propertylist_append_double(self,
"ESO QC XCSHIFT",
1329 bug_if (cpl_propertylist_append_int(self,
"ESO QC PHDEGREE", phdegree));
1330 bug_if (cpl_propertylist_append_double(self,
"ESO QC PHDISPX0", phdisp0));
1331 for (i = 1; i <= phdegree; i++) {
1332 const double coeff = cpl_polynomial_get_coeff(phdisp, &i);
1333 char * label = cpl_sprintf(
"ESO QC PHDISPX%d", (
int)i);
1335 bug_if (cpl_propertylist_append_double(self, label, coeff));
1339 bug_if (cpl_propertylist_append_double(self,
"ESO QC XCWLEN", xcwlen));
1341 bug_if (cpl_propertylist_append_int(self,
"ESO QC XCDEGREE", xcdegree));
1342 bug_if (cpl_propertylist_append_double(self,
"ESO QC XCDISPX0", xcdisp0));
1344 for (i = 1; i <= xcdegree; i++) {
1345 const double coeff = cpl_polynomial_get_coeff(xcdisp, &i);
1346 char * label = cpl_sprintf(
"ESO QC XCDISPX%d", (
int)i);
1348 bug_if (cpl_propertylist_append_double(self, label, coeff));
1354 return cpl_error_get_code();
1373static void * visir_spectro_qclist_obs(cpl_propertylist * self,
double xfwhm,
1378 cpl_propertylist_append_double _(self,
"ESO QC XFWHM", xfwhm);
1379 cpl_propertylist_append_double _(self,
"ESO QC XCENTROI", xcentro);
1399static cpl_error_code visir_vector_convolve_symm(cpl_vector * self,
1400 const cpl_vector * vsymm)
1403 const int npix = cpl_vector_get_size(self);
1404 const int ihwidth = cpl_vector_get_size(vsymm) - 1;
1405 cpl_vector * raw = cpl_vector_duplicate(self);
1406 double * pself= cpl_vector_get_data(self);
1407 double * praw = cpl_vector_get_data(raw);
1408 const double * psymm = cpl_vector_get_data_const(vsymm);
1413 skip_if_error_present();
1416 skip_if (ihwidth >= npix);
1419 for (i = 0; i < ihwidth; i++) {
1420 pself[i] = praw[i] * psymm[0];
1421 for (j = 1; j <= ihwidth; j++) {
1422 const int k = i-j < 0 ? 0 : i-j;
1423 pself[i] += (praw[k]+praw[i+j]) * psymm[j];
1428 for (i = ihwidth; i < npix-ihwidth; i++) {
1429 pself[i] = praw[i] * psymm[0];
1430 for (j = 1; j <= ihwidth; j++)
1431 pself[i] += (praw[i-j]+praw[i+j]) * psymm[j];
1434 for (i = npix-ihwidth; i < npix; i++) {
1435 pself[i] = praw[i] * psymm[0];
1436 for (j = 1; j <= ihwidth; j++) {
1437 const int k = i+j > npix-1 ? npix - 1 : i+j;
1438 pself[i] += (praw[k]+praw[i-j]) * psymm[j];
1445 cpl_vector_delete(raw);
1447 return cpl_error_get_code();
1472cpl_image * visir_spc_flip(
const cpl_image * image,
double wlen,
1473 visir_spc_resol resol, visir_data_type dtype,
1476 cpl_image * flipped = cpl_image_cast(image, CPL_TYPE_DOUBLE);
1477 visir_optmod ins_settings;
1478 if (is_flipped) *is_flipped =
false;
1480 skip_if_error_present();
1482 if ((resol == VISIR_SPC_R_HR || resol == VISIR_SPC_R_GHR) &&
1483 visir_spc_optmod_init(resol, wlen, &ins_settings,
1484 visir_data_is_aqu(dtype))) {
1485 visir_error_set(CPL_ERROR_ILLEGAL_INPUT);
1492 if (visir_data_is_aqu(dtype)) {
1493 skip_if (cpl_image_turn(flipped, 1));
1494 if ((resol == VISIR_SPC_R_HR || resol == VISIR_SPC_R_GHR) &&
1495 visir_spc_optmod_side_is_A(&ins_settings) == 0) {
1496 cpl_msg_info(cpl_func,
"Flipping image");
1497 skip_if (cpl_image_flip(flipped, 0));
1498 if (is_flipped) *is_flipped =
true;
1502 else if ((resol != VISIR_SPC_R_HR && resol != VISIR_SPC_R_GHR) ||
1503 visir_spc_optmod_side_is_A(&ins_settings) > 0) {
1505 cpl_msg_info(cpl_func,
"Flipping image");
1507 skip_if (cpl_image_flip(flipped, 0));
1508 if (is_flipped) *is_flipped =
true;
1513 if (cpl_error_get_code() && flipped) {
1514 cpl_image_delete(flipped);
1539static cpl_polynomial * visir_spc_phys_disp(
int npix,
double wlen,
1540 visir_spc_resol resol,
int ioffset,
1544 cpl_polynomial * phdisp = NULL;
1545 visir_optmod ins_settings;
1551 const cpl_size i1 = 1;
1552 const cpl_size i0 = 0;
1555 cpl_ensure(resol, CPL_ERROR_ILLEGAL_INPUT, NULL);
1556 cpl_ensure(wlen > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1557 cpl_ensure(npix > 1, CPL_ERROR_ILLEGAL_INPUT, NULL);
1564 cpl_ensure(!visir_spc_optmod_init(resol, wlen, &ins_settings, is_aqu),
1565 CPL_ERROR_ILLEGAL_INPUT, NULL);
1570 dwl = visir_spc_optmod_wlen(&ins_settings, &wlen0, &wlen1);
1572 cpl_ensure(dwl >= 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1577 if (fabs(dwl) > 2*wlen*DBL_EPSILON) cpl_msg_warning(cpl_func,
"Too large res"
1578 "idual in Scan-Angle determination [meps]: %g", dwl/DBL_EPSILON/wlen);
1580 if ((resol == VISIR_SPC_R_HR || resol == VISIR_SPC_R_GHR) &&
1581 !visir_spc_optmod_side_is_A(&ins_settings)) {
1582 const double swap = wlen1;
1586 cpl_ensure(wlen1 > wlen0, CPL_ERROR_ILLEGAL_INPUT, NULL);
1588 if (resol == VISIR_SPC_R_LRP) {
1589 phdisp = visir_spc_phys_lrp();
1593 phdisp = cpl_polynomial_new(1);
1596 disp = (wlen1-wlen0)/(npix-1);
1598 skip_if_error_present();
1600 skip_if (cpl_polynomial_set_coeff(phdisp, &i1, disp));
1602 skip_if (cpl_polynomial_set_coeff(phdisp, &i0, wlen0-disp));
1605 if ((resol == VISIR_SPC_R_HR || resol == VISIR_SPC_R_GHR) &&
1606 !visir_spc_optmod_side_is_A(&ins_settings)) {
1607 cpl_msg_info(cpl_func,
"HR B-side WLMin, WLMax, Disp: %g %g %g", wlen0,
1608 wlen1, cpl_polynomial_get_coeff(phdisp, &i1));
1610 cpl_msg_info(cpl_func,
"WLMin, WLMax, Disp: %g %g %g", wlen0, wlen1,
1611 cpl_polynomial_get_coeff(phdisp, &i1));
1614 if (resol == VISIR_SPC_R_GHR && ioffset != 0) {
1617 const double dispi = visir_spc_optmod_echelle(&ins_settings,
1618 cpl_polynomial_get_coeff(phdisp, &i1), ioffset);
1619 const double wlen0i= visir_spc_optmod_echelle(&ins_settings,
1620 cpl_polynomial_get_coeff(phdisp, &i0), ioffset);
1622 skip_if (cpl_polynomial_set_coeff(phdisp, &i1, dispi));
1624 skip_if (cpl_polynomial_set_coeff(phdisp, &i0, wlen0i));
1626 cpl_msg_info(cpl_func,
"WLc relative error(%d): %g", ioffset,
1627 (wlen0i - cpl_polynomial_eval_1d(phdisp, 1, NULL))/wlen0i);
1633 if (cpl_error_get_code() && phdisp != NULL) {
1634 cpl_polynomial_delete(phdisp);
1659cpl_bivector * visir_bivector_load_fits(
const char * file,
1660 const char * labelx,
1661 const char * labely,
1665 cpl_bivector * result = NULL;
1666 cpl_table * table = NULL;
1667 cpl_propertylist * extlist = NULL;
1668 cpl_vector * xwrapper;
1669 cpl_vector * ywrapper;
1677 bug_if (extnum < 1);
1679 next = cpl_fits_count_extensions(file);
1680 any_if(
"Could not load FITS table from (extension %d in) file: %s",
1681 extnum, file ? file :
"<NULL>");
1683 skip_if_lt(next, extnum,
"extensions in file: %s", file);
1685 table = cpl_table_load(file, extnum, 0);
1686 any_if (
"Could not load FITS table from extension %d of %d in file: %s",
1687 extnum, next, file ? file :
"<NULL>");
1689 extlist = cpl_propertylist_load_regexp(file, extnum,
"EXTNAME", 0);
1690 if (cpl_propertylist_has(extlist,
"EXTNAME")) {
1691 const char * extname = cpl_propertylist_get_string(extlist,
"EXTNAME");
1692 sext = cpl_sprintf(
" (EXTNAME=%s)", extname);
1695 nlines = cpl_table_get_nrow(table);
1696 skip_if_lt(nlines, 2,
"rows in table from extension %d%s of %d "
1697 "in %s", extnum, sext, next, file);
1699 prowx = cpl_table_get_data_double(table, labelx);
1700 any_if(
"Table from extension %d%s of %d in %s has no column %s",
1701 extnum, sext, next, file, labelx);
1703 prowy = cpl_table_get_data_double(table, labely);
1704 any_if(
"Table from extension %d%s of %d in %s has no column %s",
1705 extnum, sext, next, file, labely);
1707 xwrapper = cpl_vector_wrap(nlines, prowx);
1708 ywrapper = cpl_vector_wrap(nlines, prowy);
1710 result = cpl_bivector_wrap_vectors(xwrapper, ywrapper);
1711 cpl_table_unwrap(table, labelx);
1712 cpl_table_unwrap(table, labely);
1714 cpl_msg_info(cpl_func,
"Read %d rows from extension %d%s of %d "
1715 "in %s [%g;%g]", nlines, extnum, sext, next, file,
1716 cpl_vector_get(xwrapper, 0),
1717 cpl_vector_get(ywrapper, nlines-1));
1722 cpl_table_delete(table);
1723 cpl_propertylist_delete(extlist);
1725 if (result && cpl_error_get_code()) {
1726 cpl_bivector_delete(result);
1763static cpl_error_code visir_spc_emission(cpl_bivector * emission,
1764 const cpl_vector * boundary,
1765 const cpl_bivector * temiss,
1766 const cpl_bivector * tqeff,
1767 const cpl_vector * vsymm,
1770 cpl_bivector * tqeffi = NULL;
1771 cpl_vector * planck = NULL;
1772 const int npix = cpl_bivector_get_size(emission);
1775 bug_if(emission == NULL);
1776 bug_if(boundary == NULL);
1777 bug_if(temiss == NULL);
1778 bug_if(tqeff == NULL);
1783 skip_if(cpl_vector_get_size(boundary) != npix + 1);
1785 planck = cpl_vector_new(npix);
1786 skip_if_error_present();
1790 cpl_photom_fill_blackbody(planck, CPL_UNIT_ENERGYRADIANCE,
1791 cpl_bivector_get_x(emission),
1792 CPL_UNIT_LENGTH, 253);
1794 skip_if (visir_vector_resample(cpl_bivector_get_y(emission),
1798 skip_if (visir_vector_convolve_symm(cpl_bivector_get_y(emission),
1801 skip_if (cpl_vector_multiply(cpl_bivector_get_y(emission), planck));
1805 cpl_photom_fill_blackbody(planck, CPL_UNIT_ENERGYRADIANCE,
1806 cpl_bivector_get_x(emission),
1807 CPL_UNIT_LENGTH, temp);
1810 skip_if (cpl_vector_multiply_scalar(planck, 0.12));
1813 skip_if (cpl_vector_add(cpl_bivector_get_y(emission), planck));
1816 tqeffi = cpl_bivector_duplicate(emission);
1817 skip_if (cpl_bivector_interpolate_linear(tqeffi, tqeff));
1819 skip_if (cpl_vector_multiply(cpl_bivector_get_y(emission),
1820 cpl_bivector_get_y(tqeffi)));
1824 cpl_bivector_delete(tqeffi);
1825 cpl_vector_delete(planck);
1827 return cpl_error_get_code();
1855static cpl_vector * cpl_spc_convolve_init(
int maxlen,
double slitw,
1856 double fwhm,
int doplot)
1859 const double sigma = fwhm * CPL_MATH_SIG_FWHM;
1860 const int ihtophat = (int)slitw/2;
1861 const int gausshlen = 1 + 5 * sigma + ihtophat < maxlen/2
1862 ? 1 + 5 * sigma + ihtophat : maxlen/2 - 1;
1864 const int convolen = 1 + 10 * sigma + 8*ihtophat;
1865 cpl_vector * self = cpl_vector_new(gausshlen);
1866 cpl_vector * tophat = cpl_vector_new(convolen);
1870 cpl_image * iself = cpl_image_wrap_double(gausshlen, 1,
1871 cpl_vector_get_data(self));
1874 skip_if_error_present();
1876 skip_if( slitw <= 0.0);
1877 skip_if( fwhm <= 0.0);
1878 skip_if( convolen < 2 * gausshlen);
1881 skip_if (cpl_image_fill_gaussian(iself, 1.0, 1.0, CPL_MATH_SQRT2PI,
1884 if (doplot > 2) visir_vector_plot(
"set grid;",
"t 'Right Half of Gaussian' "
1885 "w linespoints",
"", self);
1888 skip_if( cpl_vector_fill(tophat, 0.0));
1890 for (i = convolen/2-ihtophat; i < 1+convolen/2+ihtophat; i++)
1891 skip_if (cpl_vector_set(tophat, i, 1.0/(1.0+2.0*ihtophat)));
1894 skip_if (visir_vector_convolve_symm(tophat, self));
1896 if (doplot > 2) visir_vector_plot(
"set grid;",
"t 'Full Width Convolution' "
1897 "w linespoints",
"", tophat);
1902 memcpy(cpl_vector_get_data(self),
1903 cpl_vector_get_data(tophat) + convolen/2,
1904 sizeof(
double)*gausshlen);
1907 for (i = 0 ; i < gausshlen; i++)
1908 skip_if (cpl_vector_set(self, i, cpl_vector_get(tophat,
1912 skip_if_error_present();
1914 cpl_msg_info(cpl_func,
"Convolving Model Spectrum, Gauss-sigma=%g, "
1915 "Tophat-width=%d, Truncation-Error=%g with width=%d", sigma,
1917 cpl_vector_get(self,gausshlen-1)/cpl_vector_get(self,0),
1920 if (doplot > 1) visir_vector_plot(
"set grid;",
"t 'Right Half of Convolution"
1921 "' w linespoints",
"", self);
1925 cpl_vector_delete(tophat);
1926 cpl_image_unwrap(iself);
1928 if (cpl_error_get_code()) {
1929 cpl_vector_delete(self);
1938static cpl_error_code
1939fit_gaussians(
const cpl_image * flipped,
const cpl_vector * error,
1940 cpl_size icollo, cpl_size icolhi,
1941 cpl_propertylist * qclist)
1943 cpl_size nrow = cpl_image_get_size_y(flipped);
1944 cpl_size ncol = cpl_image_get_size_x(flipped);
1945 icollo = CX_MAX(1, icollo);
1946 icolhi = CX_MIN(ncol, icolhi);
1947 cpl_errorstate cleanstate = cpl_errorstate_get();
1949 double sigs_err = 0.;
1951 double peaks_err = 0.;
1953 for (cpl_size row = 0; row < nrow; row++) {
1954 const cpl_binary * dmask = cpl_image_get_bpm_const(flipped) ?
1955 cpl_mask_get_data_const(cpl_image_get_bpm_const(flipped)) : NULL;
1956 const double *dflipped = cpl_image_get_data_double_const(flipped);
1957 double * dx = cpl_malloc(ncol *
sizeof(*dx));
1958 double * dy = cpl_malloc(ncol *
sizeof(*dy));
1959 double * dye = cpl_malloc(ncol *
sizeof(*dye));
1964 for (cpl_size i = icollo; i <= icolhi; i++) {
1965 if (dmask == NULL || !dmask[row * ncol + i]) {
1967 dy[n] = dflipped[row * ncol + (i - 1)];
1968 dye[n] = cpl_vector_get(error, (i - 1));
1973 x = cpl_vector_wrap(n, dx);
1974 y = cpl_vector_wrap(n, dy);
1975 ye = cpl_vector_wrap(n, dye);
1976 double x0, sigma, sigma_err, peak, peak_err;
1977 fit_1d_gauss(x, y, ye, &x0, NULL, &peak, &peak_err, &sigma, &sigma_err);
1978 if (cpl_error_get_code() != CPL_ERROR_NONE) {
1979 cpl_msg_debug(cpl_func,
"FIT row %lld failed", row);
1980 cpl_errorstate_set(cleanstate);
1983 sigs[nmeas] = sigma;
1984 peaks[nmeas] = peak;
1985 sigs_err += sigma * sigma;
1986 peaks_err += peak * peak;
1988 cpl_msg_debug(cpl_func,
"FIT row %lld x %g sig %g +- %g "
1990 row, x0, sigma, sigma_err, peak, peak_err);
1992 cpl_vector_delete(x);
1993 cpl_vector_delete(y);
1994 cpl_vector_delete(ye);
2002 cpl_vector * sigv = cpl_vector_wrap(nmeas, sigs);
2003 cpl_vector * peakv = cpl_vector_wrap(nmeas, peaks);
2004 double medsigma = cpl_vector_get_median(sigv);
2005 double medsigma_err = sqrt(sigs_err) * sqrt(CPL_MATH_PI_2) / nmeas;
2006 double medpeak = cpl_vector_get_median(peakv);
2007 double medpeak_err = sqrt(peaks_err) * sqrt(CPL_MATH_PI_2) / nmeas;
2008 cpl_msg_info(cpl_func,
"Median FWHM of spectrum: %g +- %g, Peak %g +- %g",
2009 medsigma, medsigma_err, medpeak, medpeak_err);
2010 cpl_propertylist_append_double(qclist,
"ESO QC GAUSSFIT FWHM",
2012 cpl_propertylist_set_comment(qclist,
"ESO QC GAUSSFIT FWHM",
"[pix]");
2013 cpl_propertylist_append_double(qclist,
"ESO QC GAUSSFIT FWHM_ERR",
2014 medsigma_err * 2.355);
2015 cpl_propertylist_append_double(qclist,
"ESO QC GAUSSFIT PEAK", medpeak);
2016 cpl_propertylist_set_comment(qclist,
"ESO QC GAUSSFIT PEAK",
"[adu/s]");
2017 cpl_propertylist_append_double(qclist,
"ESO QC GAUSSFIT PEAK_ERR",
2019 cpl_vector_unwrap(sigv);
2020 cpl_vector_unwrap(peakv);
2022 return cpl_error_get_code();
2033static cpl_error_code
2034add_qc_background_sigma(
const cpl_image * flipped, cpl_propertylist * qclist)
2038 if (cpl_image_get_size_y(flipped) > VISIR_AQU_APPROX_WLEN13) {
2039 lly = VISIR_AQU_APPROX_WLEN8;
2040 ury = VISIR_AQU_APPROX_WLEN13;
2044 ury = cpl_image_get_size_y(flipped);
2047 cpl_image * cutimg =
2048 cpl_image_extract(flipped, 1, lly, cpl_image_get_size_x(flipped), ury);
2051 double bkgmad, bkgmed;
2052 bkgmed = cpl_image_get_mad(cutimg, &bkgmad);
2053 for (
size_t i = 0; i < 3; i++) {
2055 cpl_mask_threshold_image_create(cutimg,
2056 bkgmed - bkgmad * CPL_MATH_STD_MAD * 3,
2057 bkgmed + bkgmad * CPL_MATH_STD_MAD * 3);
2059 cpl_image_reject_from_mask(cutimg, rej);
2060 cpl_mask_delete(rej);
2061 bkgmed = cpl_image_get_mad(cutimg, &bkgmad);
2064 cpl_propertylist_append_double(qclist,
"ESO QC BACKGD SIGMA",
2065 bkgmad * CPL_MATH_STD_MAD);
2066 cpl_propertylist_set_comment(qclist,
"ESO QC BACKGD SIGMA",
2067 "[adu/s] background corrected");
2068 cpl_image_delete(cutimg);
2070 return cpl_error_get_code();
2080static double * visir_bkg_linfit(
const cpl_image * row) {
2081 INIT(3, rv,
double, v, cpl_free) = NULL;
2083 const int n = cpl_image_get_size_x _(row);
2084 const int ngood = n - cpl_image_count_rejected _(row);
2087 SET(x,
double, v, cpl_free) = cpl_malloc(ngood *
sizeof(
double));
2088 SET(y,
double, v, cpl_free) = cpl_malloc(ngood *
sizeof(
double));
2089 for (
int bad, i = 0, j = 0; i < n; ++i) {
2090 const double cand = cpl_image_get _(row, i+1, 1, &bad);
2097 double c0, c1, cov00, cov01, cov11;
2098 gsl_fit_linear(x->o, 1, y->o, 1, ngood, &c0, &c1, &cov00, &cov01, &cov11,
2102 RESET(rv) = cpl_malloc(n *
sizeof(
double));
2103 for (
int i = 0; i < n; ++i) gsl_fit_linear_est(
2104 i, c0, c1, cov00, cov01, cov11, rv->o + i, HOLE(
double));
2129int visir_norm_coord(
const bool rev,
const float coord,
const int lcol,
2130 const int rcol,
const visir_apdefs * aps)
2132 const int x = coord < 0 ? -coord : coord;
2133 if (copysign(1.0, coord) > 0.0)
return
2134 rev ? rcol - aps->limits[x].l + 1 : aps->limits[x].r - lcol + 1;
2136 rev ? rcol - aps->limits[x].r + 1 : aps->limits[x].l - lcol + 1;
2148static cpl_image * visir_image_filter_median(
const cpl_image * image,
2149 const int nx,
const int ny)
2151 INIT(1, filtered, cpl_image) = cpl_image_duplicate _(image);
2152 if (nx < 2)
return CLEANUP();
2155 return ABORT(CPL_ERROR_ILLEGAL_INPUT,
"kernel size must be odd");
2157 const cpl_size xsz = cpl_image_get_size_x _(image);
2158 const cpl_size ysz = cpl_image_get_size_y _(image);
2159 const cpl_type type = cpl_image_get_type _(image);
2160 RESET(filtered) = cpl_image_new _(xsz, ysz, type);
2162 SET(kernel, cpl_mask) = cpl_mask_new _(nx, ny);
2163 cpl_mask_not _(kernel->o);
2164 cpl_image_filter_mask _(filtered->o, image, kernel->o, CPL_FILTER_MEDIAN,
2186static void * visir_extraction(
const cpl_image * insci,
const cpl_image * invar,
2187 const cpl_image * insky, cpl_vector * outext,
2188 cpl_vector * outsky, cpl_vector * outerr,
2189 cpl_image * outwgt,
const visir_spc_config * cfg,
2190 const int method,
const int ncomb,
const int beg,
2195 if (!insci || !invar || !insky || !cfg ||
2196 !outext || !outsky || !outerr || !outwgt)
2197 return ABORT(CPL_ERROR_NULL_INPUT);
2199 const int specLen = cpl_image_get_size_x _(insci);
2200 const int numRows = cpl_image_get_size_y _(insci);
2201 if (beg < 0 || beg >= numRows || end < 0 || end >= numRows || beg > end)
2202 return ABORT(CPL_ERROR_INCOMPATIBLE_INPUT);
2204 double * oext_d = cpl_vector_get_data _(outext);
2205 double * osky_d = cpl_vector_get_data _(outsky);
2206 double * oerr_d = cpl_vector_get_data _(outerr);
2207 double * owgt_d = cpl_image_get_data _(outwgt);
2209 const double *
const isci_d = cpl_image_get_data_const _(insci);
2210 const double *
const isky_d = cpl_image_get_data_const _(insky);
2211 const double *
const ivar_d = invar ? cpl_image_get_data_const(invar) : NULL;
2216 if (method && abs(end-beg) > cfg->ox_kernel + 2 ) {
2217 SET(smoothed, cpl_image) = visir_image_filter_median _(
2218 insci, cfg->ox_kernel, cfg->ox_kernel);
2219 double * smoo_d = cpl_image_get_data _(smoothed->o);
2220 for (
int i = 0; i < specLen; i++) {
2222 for (
int j = beg; j <= end; j++)
2223 oext_d[i] += smoo_d[i + j * specLen];
2227 for (
int i = 0; i < specLen; i++) {
2229 for (
int j = beg; j <= end; j++)
2230 oext_d[i] += isci_d[i + j * specLen];
2236 SET(buf,
double, v, cpl_free) = cpl_calloc(specLen,
sizeof(
double));
2237 for (
int iter = 0; iter < cfg->ox_niters; iter++) {
2240 for (
int i = 0; i < specLen; i++) {
2241 for (
int j = beg; j <= end; j++) {
2242 const int index = i + j * specLen;
2244 owgt_d[index] = fabs(oext_d[i]) > 0.00001
2245 ? isci_d[index] / oext_d[i] : 0.0;
2249 SET(wrap, cpl_vector, p) = NULL;
2250 for (
int j = beg; j <= end; j++) {
2254 for (
int i = 0; i < specLen - cfg->ox_smooth; i++) {
2255 RESET(wrap) = cpl_vector_wrap _(cfg->ox_smooth, owgt_d +
2257 double value = cpl_vector_get_median_const _(wrap->o);
2258 if (value < 0) value = 0.0;
2259 buf->o[i + cfg->ox_smooth / 2] = value;
2263 RESET(wrap) = cpl_vector_wrap _(cfg->ox_smooth / 2, owgt_d +
2265 double value = cpl_vector_get_mean _(wrap->o);
2266 if (value < 0) value = 0.0;
2267 for (
int i = 0; i < cfg->ox_smooth / 2; i++) buf->o[i] = value;
2270 RESET(wrap) = cpl_vector_wrap _(cfg->ox_smooth / 2, owgt_d +
2271 specLen - cfg->ox_smooth / 2 + j * specLen);
2272 value = cpl_vector_get_mean _(wrap->o);
2273 if (value < 0) value = 0.0;
2274 for (
int i = 0; i < cfg->ox_smooth / 2; i++)
2275 buf->o[i + specLen - cfg->ox_smooth / 2] = value;
2276 for (
int i = 0; i < specLen; i++)
2277 owgt_d[i + j * specLen] = buf->o[i];
2281 for (
int i = 0; i < specLen; i++) {
2283 for (
int j = beg; j <= end; j++)
2284 value += owgt_d[i + j * specLen];
2285 if (value > 0.00001)
2286 for (
int j = beg; j <= end; j++)
2287 owgt_d[i + j * specLen] /= value;
2289 for (
int j = beg; j <= end; j++)
2290 owgt_d[i + j * specLen] = 0.0;
2295 for (
int i = 0; i < specLen; i++) {
2296 double sumSci, sumSky, sumWgt, sumProf, sumVar;
2297 sumSci = sumSky = sumWgt = sumProf = sumVar = 0.0;
2298 for (
int j = beg; j <= end; j++) {
2299 const int index = i + j * specLen;
2306 double var = cfg->ron * cfg->ron + fabs(oext_d[i] *
2307 owgt_d[index] + isky_d[index]) / cfg->gain;
2312 double value = isci_d[index] - oext_d[i] * owgt_d[index];
2313 if (fabs(value) / sqrt(var) < cfg->ox_sigma) {
2314 const double weight = 1000000 * owgt_d[index] / var;
2315 sumSci += weight * isci_d[index];
2316 sumSky += weight * isky_d[index];
2317 sumWgt += weight * owgt_d[index];
2318 sumProf += owgt_d[index];
2323 if (ivar_d) sumVar += weight * weight * ivar_d[index];
2327 if (sumWgt > 0.00001) {
2328 oext_d[i] = sumSci / sumWgt;
2329 osky_d[i] = sumSky / sumWgt;
2332 oerr_d[i] = sqrt(sumVar / sumWgt / sumWgt);
2336 oerr_d[i] = 1000 * sqrt(sumProf / sumWgt);
2351 for (
int i = 0; i < specLen; i++) {
2353 for (
int j = beg; j <= end; j++)
2354 osky_d[i] += isky_d[i + j * specLen];
2358 for (
int i = 0; i < specLen; i++) {
2362 for (
int j = beg; j <= end; j++)
2363 oerr_d[i] += ivar_d[i + j * specLen];
2364 oerr_d[i] = sqrt(oerr_d[i]);
2367 oerr_d[i] = sqrt(cfg->ron * cfg->ron + fabs(oext_d[i] +
2368 osky_d[i]) / cfg->gain);
2375static cpl_bivector * visir_spc_extract(cpl_image * flipped,
2377 cpl_propertylist * qclist,
2378 cpl_image ** pweight2d,
2379 const visir_spc_config * cfg,
2380 const visir_apdefs * aps,
2381 const bool rev,
const cpl_size ncomb)
2383 extract_func * meth = aps->ident < 0 ? visir_spc_oldex : visir_spc_newex;
2384 return meth(flipped, lcol, rcol, qclist, pweight2d, cfg, aps, rev, ncomb);
2400static cpl_bivector * visir_spc_newex(cpl_image * flipped,
2402 cpl_propertylist * qclist,
2403 cpl_image ** pweight2d,
2404 const visir_spc_config * cfg,
2405 const visir_apdefs * aps,
2406 const bool rev,
const cpl_size ncomb)
2408 INIT(20, rv, cpl_bivector) = NULL;
2410 if (!flipped || !qclist || !cfg || !aps)
return ABORT(CPL_ERROR_NULL_INPUT);
2411 const int ncol = cpl_image_get_size_x _(flipped);
2412 const int specLen = cpl_image_get_size_y _(flipped);
2413 const int oo = cfg->orderoffset;
2414 const size_t xn = cfg->extract;
2416 if (aps->nlimits < 1)
return ABORT(CPL_ERROR_INCOMPATIBLE_INPUT);
2417 if (!pweight2d || *pweight2d)
return ABORT(CPL_ERROR_INCOMPATIBLE_INPUT);
2418 if (ncol != rcol-lcol+1)
return ABORT(CPL_ERROR_INCOMPATIBLE_INPUT);
2419 MSG_DBG(
":%s:%ld: > [%d;%d] ([%d;%d]) <", pn(oo), xn, 1, ncol, rcol, lcol);
2423 SET(key,
char, v, cpl_free) = cpl_sprintf(
"ESO DRS APDEF%d", aps->ident);
2424 cpl_propertylist_append_string _(cfg->phu, key->o, line->o);
2427 add_qc_background_sigma _(flipped, qclist);
2429 const bool apex = aps->extract_method == VISIR_EXTRACT_METHOD_APERTURE;
2431 SET(bkg, cpl_image) = NULL;
2432 SET(diff, cpl_image) = cpl_image_duplicate _(flipped);
2433 if (!cfg->bkgcorrect) {
2434 const cpl_type type = cpl_image_get_type _(flipped);
2435 RESET(bkg) = cpl_image_new _(ncol, specLen, type);
2436 MSG_WARN(
"Sky subtraction is not enabled: extraction results may be "
2441 RESET(bkg) = cpl_image_duplicate _(flipped);
2445 int beg = rev ? aps->nlimits - 1 : 1;
2446 int end = rev ? 0 : aps->nlimits;
2447 const int inc = rev ? -1 : 1;
2450 for (
int a = beg; a != end; a += inc) {
2451 int l = visir_norm_coord(rev, -a, lcol, rcol, aps);
2452 int r = visir_norm_coord(rev, +a, lcol, rcol, aps);
2453 MSG_DBG(
":%s:%ld: [%d;%d] ([%d;%d])", pn(oo), xn, l, r,
2454 aps->limits[a].r, aps->limits[a].l);
2457 if ((1 <= lp && lp <= ncol) || (1 <= l && l <= ncol)) {
2458 const int trunc_lp = lp < 1 ? 1 : lp;
2459 const int trunc_l = l > ncol ? ncol : l;
2460 if (trunc_lp <= trunc_l) {
2461 cpl_image_fill_window _(bkg->o, trunc_lp, 1, trunc_l,
2462 specLen, -INFINITY);
2463 MSG_DBG(
":%s:%ld: [%d;%d] rejected", pn(oo), xn, trunc_lp,
2469 if (1 <= lp && lp <= ncol) {
2470 cpl_image_fill_window _(bkg->o, lp, 1, ncol, specLen, -INFINITY);
2471 MSG_DBG(
":%s:%ld: [%d;%d] rejected", pn(oo), xn, lp, ncol);
2473 cpl_image_reject_value _(bkg->o, CPL_VALUE_MINUSINF);
2476 double (*method)(
const cpl_image *) = cpl_image_get_median;
2477 if (apex && aps->sky_method == VISIR_SKY_METHOD_AVERAGE)
2478 method = cpl_image_get_mean;
2479 const bool linear = apex && aps->sky_method == VISIR_SKY_METHOD_LINFIT;
2482 for (
int r = 0; r < specLen; ++r) {
2483 SET(row, cpl_image) = cpl_image_extract _(bkg->o, 1, r+1, ncol, r+1);
2485 SET(levels,
double, v, cpl_free) = visir_bkg_linfit _(row->o);
2486 for (cpl_size c = 0; c < ncol; ++c) {
2487 cpl_image_set _(bkg->o, c+1, r+1, levels->o[c]); }
2490 double level = method _(row->o);
2491 cpl_image_fill_window _(bkg->o, 1, r+1, ncol, r+1, level);
2500 cpl_image_subtract _(diff->o, bkg->o);
2504 SET(var, cpl_image) = cpl_image_duplicate _(flipped);
2505 cpl_image_abs _(var->o);
2506 cpl_image_divide_scalar _(var->o, cfg->gain);
2507 cpl_image_add_scalar _(var->o, cfg->ron * cfg->ron);
2510 cpl_image_turn _(diff->o, 1);
2511 cpl_image_turn _(var->o, 1);
2512 cpl_image_turn _(bkg->o, 1);
2515 SET(spc, cpl_vector) = cpl_vector_new _(specLen);
2516 SET(sky, cpl_vector) = cpl_vector_new _(specLen);
2517 SET(err, cpl_vector) = cpl_vector_new _(specLen);
2518 SET(wgt, cpl_image) = cpl_image_new _(specLen, ncol, CPL_TYPE_DOUBLE);
2519 const int beg = ncol - visir_norm_coord(rev, +0.0, lcol, rcol, aps);
2520 const int end = ncol - visir_norm_coord(rev, -0.0, lcol, rcol, aps);
2521 visir_extraction _(diff->o, var->o, bkg->o, spc->o, sky->o, err->o, wgt->o,
2522 cfg, apex ? 0 : 1, ncomb, beg, end);
2524 cpl_image_turn _(wgt->o, -1);
2525 RESET(rv) = cpl_bivector_wrap_vectors _(spc->o, err->o);
2526 *pweight2d = YANK(wgt, cpl_image);
2527 YANK(spc); YANK(err);
2553static cpl_bivector * visir_spc_oldex(cpl_image * flipped,
2555 cpl_propertylist * qclist,
2556 cpl_image ** pweight2d,
2557 const visir_spc_config * cfg,
2558 const visir_apdefs * aps,
2559 const bool rev,
const cpl_size ncomb)
2561 INIT(20, rv, cpl_bivector) = NULL;
2563 if (!flipped || !qclist || !cfg || !aps)
return ABORT(CPL_ERROR_NULL_INPUT);
2564 const int ncol = cpl_image_get_size_x _(flipped);
2565 const int nrow = cpl_image_get_size_y _(flipped);
2566 const int oo = cfg->orderoffset;
2567 const size_t xn = cfg->extract;
2571 const double sigma = VISIR_SPECTRO_SIGMA;
2573 if (sigma <= 0.0)
return ABORT(CPL_ERROR_UNSUPPORTED_MODE);
2574 if (aps->nlimits < 1)
return ABORT(CPL_ERROR_INCOMPATIBLE_INPUT);
2575 if (!pweight2d || *pweight2d)
return ABORT(CPL_ERROR_INCOMPATIBLE_INPUT);
2576 if (ncol != rcol-lcol+1)
return ABORT(CPL_ERROR_INCOMPATIBLE_INPUT);
2577 MSG_DBG(
":%s:%ld: > [%d;%d] ([%d;%d]) <", pn(oo), xn, 1, ncol, rcol, lcol);
2581 SET(key,
char, v, cpl_free) = cpl_sprintf(
"ESO DRS APDEF%d", aps->ident);
2582 cpl_propertylist_append_string _(cfg->phu, key->o, line->o);
2585 add_qc_background_sigma _(flipped, qclist);
2588 SET(orig, cpl_image) = cpl_image_duplicate _(flipped);
2594 if (cfg->bkgcorrect) {
2597 SET(work, cpl_image) = cpl_image_duplicate _(flipped);
2599 for (
int r = 0; r < nrow; ++r) {
2600 SET(row, cpl_image) = cpl_image_extract _(work->o, 1, r+1, ncol, r+1);
2601 double level = cpl_image_get_median _(row->o);
2602 cpl_image_fill_window _(work->o, 1, r+1, ncol, r+1, level);
2606 cpl_image_subtract _(flipped, work->o);
2610 const int is_echelle = ncol <= 2 * (whechelle + 1);
2615 double mean = cpl_image_get_mean _(flipped);
2616 MSG_INFO(
"Combined image has mean: %g", mean);
2618 SET(col, cpl_vector) = cpl_vector_new _(nrow);
2621 double * pweight = cpl_image_get_data_double _(flipped);
2622 for (
int r=0; r < nrow; r++, pweight += ncol) {
2625 SET(imrow, cpl_image, p) = cpl_image_wrap_double _(ncol, 1, pweight);
2628 mean = cpl_image_get_mean _(imrow->o);
2629 cpl_vector_set _(col->o, r, mean);
2632 cpl_image_subtract_scalar _(imrow->o, mean);
2637 if (cfg->plot > 1) visir_vector_plot(
2638 "set grid;",
"t 'Estimated Background' w linespoints",
"", col->o);
2644 SET(spatial, cpl_image) = cpl_image_collapse_create _(flipped, 0);
2645 cpl_image_divide_scalar _(spatial->o, nrow);
2649 SET(iweight, cpl_image) = cpl_image_duplicate _(spatial->o);
2650 cpl_image_normalise _(iweight->o, CPL_NORM_ABSFLUX);
2651 const double sqflux = cpl_image_get_sqflux _(iweight->o);
2652 const double weight_2norm = sqrt(sqflux);
2653 MSG_INFO(
"2-norm of weights: %g", weight_2norm);
2655 if (cfg->plot > 1) visir_image_row_plot(
2656 "set grid;",
"t 'Cleaned, normalized combined image with spectral "
2657 "direction averaged' w linespoints",
"", iweight->o, 1, 1, 1);
2661 const double sp_median = cpl_image_get_median _(spatial->o);
2662 double stdev2d = visir_img_phot_sigma_clip _(flipped);
2663 stdev2d /= sqrt(nrow);
2664 MSG_INFO(
"spatial median %g and stdev %g", sp_median, stdev2d);
2668 SET(binary, cpl_mask) = cpl_mask_threshold_image_create _(spatial->o,
2669 sp_median - sigma * stdev2d, sp_median + sigma * stdev2d);
2670 int mspix = cpl_mask_count _(binary->o);
2671 if (mspix == ncol)
return ABORT(CPL_ERROR_DATA_NOT_FOUND,
"%d spatial "
2672 "weights too noisy. sigma=%g. stdev2d=%g. Spatial median=%g", ncol,
2673 sigma, stdev2d, sp_median);
2674 MSG_INFO(
"Pixels of noise (%g +/- %g*%g): %d", sp_median, stdev2d, sigma,
2676 cpl_image_reject_from_mask _(spatial->o, binary->o);
2682 cpl_image_get_maxpos _(spatial->o, &ifwhm, HOLE(cpl_size));
2683 const double max = cpl_image_get _(spatial->o, ifwhm, 1, &rejected);
2684 if (rejected)
return ABORT(CPL_ERROR_ILLEGAL_OUTPUT);
2685 if (max <= 0.0)
return ABORT(CPL_ERROR_DATA_NOT_FOUND,
"Cannot compute "
2686 "FWHM on a collapsed spectrum with a non-positive maximum: %g (at "
2687 "i=%lld)", max, ifwhm);
2689 if (cfg->plot > 1) {
2690 visir_image_col_plot(
"",
"t 'Most intense column' w linespoints",
2691 "", flipped, ifwhm, ifwhm, 1);
2692 visir_image_row_plot(
"set grid;",
"t 'Combined image with spectral "
2693 "direction collapsed' w linespoints",
2694 "", spatial->o, 1, 1, 1);
2699 int ilnoise, ihnoise;
2701 for (ilnoise = ifwhm; ilnoise > 0 &&
2702 !cpl_image_is_rejected(spatial->o, ilnoise, 1); ilnoise--);
2704 for (ihnoise = ifwhm; ihnoise <= ncol &&
2705 !cpl_image_is_rejected(spatial->o, ihnoise, 1); ihnoise++);
2707 if (!ilnoise) ilnoise = 1;
2708 if (ihnoise > ncol) ihnoise = ncol;
2711 const double xcentro = cpl_image_get_centroid_x_window _(spatial->o,
2712 ilnoise, 1, ihnoise, 1);
2714 cpl_image_get_fwhm _(spatial->o, ifwhm, 1, &xfwhm, HOLE(
double));
2715 visir_spectro_qclist_obs _(qclist, xfwhm, xcentro);
2716 MSG_INFO(
"Spatial FWHM(%d:%lld:%d:%g): %g", ilnoise, ifwhm, ihnoise,
2724 int iright = ncol - 5;
2726 if (ileft > xcentro - xfwhm * 2) ileft = xcentro - xfwhm * 2;
2727 if (iright < xcentro + xfwhm * 2) iright = xcentro + xfwhm * 2;
2728 MSG_INFO(
"HRG pixels of noise: [1 %d] [%d %d]", ileft, iright, ncol);
2732 cpl_mask_xor _(binary->o, binary->o);
2736 cpl_binary * pbin = cpl_mask_get_data _(binary->o);
2737 for (
int i = 0; i < ncol; i++) pbin[i] = CPL_BINARY_0;
2738 for (
int i = 0; i < ileft; i++) pbin[i] = CPL_BINARY_1;
2739 for (
int i = iright; i < ncol; i++) pbin[i] = CPL_BINARY_1;
2741 mspix = cpl_mask_count _(binary->o);
2742 MSG_INFO(
"Pixels of noise (post-echelle refinement): %d", mspix);
2744 if (mspix < 2)
return ABORT(CPL_ERROR_DATA_NOT_FOUND,
"Cannot estimate "
2745 "spectrum noise with just %d pixels of noise", mspix);
2749 SET(locnoise, cpl_image) = cpl_image_new_from_mask _(binary->o);
2753 SET(error, cpl_vector) = cpl_vector_new _(nrow);
2754 for (
int r = 0; r < nrow; r++) {
2758 SET(imrow, cpl_image) = cpl_image_extract _(flipped, 1, r+1, ncol, r+1);
2762 SET(objs, cpl_apertures) = cpl_apertures_new_from_image _(imrow->o,
2767 double stdev1d = cpl_apertures_get_stdev _(objs->o, 1);
2772 double npp = weight_2norm * stdev1d;
2776 cpl_vector_set _(error->o, r, npp);
2780 fit_gaussians _(flipped, error->o, ifwhm - 20, ifwhm + 20, qclist);
2788 cpl_vector * spectrum = NULL;
2789 cpl_vector * divisor = NULL;
2790 for (
int c = 1; c <= ncol; c++) {
2792 SET(ntor, cpl_vector) = cpl_vector_new_from_image_column _(flipped, c);
2793 SET(dtor, cpl_vector) = NULL;
2796 const double weight = cpl_image_get _(iweight->o, c, 1, &rejected);
2797 if (rejected)
return ABORT(CPL_ERROR_DATA_NOT_FOUND);
2799 if (weight == 0)
continue;
2800 cpl_vector_multiply_scalar _(ntor->o, weight);
2804 cpl_vector_add _(spectrum, ntor->o);
2806 spectrum = YANK(ntor, cpl_vector);
2809 cpl_vector_add _(divisor, dtor->o);
2811 divisor = YANK(dtor, cpl_vector);
2813 if (!spectrum)
return ABORT(CPL_ERROR_ILLEGAL_OUTPUT);
2815 if (divisor) { cpl_vector_divide _(spectrum, divisor); }
2816 double min = cpl_vector_get_min _(spectrum);
2817 if (min < 0) MSG_WARN(
"Extracted spectrum has negative intensity: %g", min);
2820 *pweight2d = cpl_image_new _(ncol, nrow, CPL_TYPE_DOUBLE);
2821 for (
int r=1; r <= nrow; r++) {
2822 cpl_image_copy _(*pweight2d, iweight->o, 1, r); }
2824 visir_image_plot(
"",
"t 'The weight map'",
"", *pweight2d);
2826 RESET(rv) = cpl_bivector_wrap_vectors _(spectrum, error->o);
2829 visir_bivector_plot(
"",
"t 'error versus spectrum'",
"", rv->o);
2857static cpl_error_code visir_spectro_fill(cpl_vector * self,
2858 const cpl_polynomial * disp,
2859 irplib_base_spectrum_model * model)
2862 visir_spectrum_model * mymodel = (visir_spectrum_model*)model;
2863 const cpl_size npix = cpl_vector_get_size(self);
2865 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
2866 cpl_ensure_code(model, CPL_ERROR_NULL_INPUT);
2867 cpl_ensure_code(disp, CPL_ERROR_NULL_INPUT);
2869 cpl_vector * wavelength = cpl_vector_new(npix);
2870 cpl_bivector * emission = cpl_bivector_wrap_vectors(wavelength, self);
2871 cpl_vector * boundary = cpl_vector_new(npix + 1);
2875 skip_if (cpl_vector_fill_polynomial(cpl_bivector_get_x(emission),
2877 skip_if (cpl_vector_fill_polynomial(boundary, disp, 0.5, 1));
2880 skip_if (visir_spc_emission(emission, boundary, mymodel->lines,
2881 mymodel->tqeff, mymodel->vsymm,
2885 cpl_bivector_unwrap_vectors(emission);
2886 cpl_vector_delete(wavelength);
2887 cpl_vector_delete(boundary);
2889 return cpl_error_get_code();
2910static cpl_error_code visir_spectro_refine(cpl_polynomial * self,
2911 const cpl_vector * xc_vector,
2912 visir_spectrum_model * pmymodel,
2913 const cpl_polynomial * phdisp,
2914 int hsize, cpl_boolean doplot,
2915 visir_spc_resol resol,
2917 cpl_boolean * pdidshift,
2920 const int subres = VISIR_XC_SUBSEARCH;
2921 cpl_polynomial * shifted = NULL;
2922#ifdef VISIR_SPC_CAL_HIGH
2923 const int fitdeg = 2;
2924 double pixstep = 0.5;
2925 double pixtol = 1e-5;
2926 const int maxite = fitdeg * 200;
2929 const int clines = (int)(cpl_bivector_get_size(pmymodel->lines) *
2930 cpl_vector_get_size(xc_vector));
2931 cpl_errorstate prestate = cpl_errorstate_get();
2934 cpl_ensure_code(self, CPL_ERROR_NULL_INPUT);
2935 cpl_ensure_code(xc_vector, CPL_ERROR_NULL_INPUT);
2936 cpl_ensure_code(pmymodel, CPL_ERROR_NULL_INPUT);
2937 cpl_ensure_code(phdisp, CPL_ERROR_NULL_INPUT);
2938 cpl_ensure_code(pxc, CPL_ERROR_NULL_INPUT);
2939 cpl_ensure_code(pdidshift, CPL_ERROR_NULL_INPUT);
2940 cpl_ensure_code(pdelta, CPL_ERROR_NULL_INPUT);
2942 skip_if(cpl_polynomial_copy(self, phdisp));
2944#ifdef VISIR_SPC_CAL_HIGH
2945 if (irplib_polynomial_find_1d_from_correlation_all
2946 (self, fitdeg, xc_vector, 1, clines,
2947 (irplib_base_spectrum_model*)pmymodel,
2948 visir_spectro_fill, pixtol, pixstep,
2949 hsize, maxite, maxfail, maxcont, doplot, pxc) || *pxc <= 0.0) {
2951 irplib_error_recover(prestate,
"Could not optimize %d "
2952 "coefficients, trying shifting", fitdeg);
2953 skip_if(cpl_polynomial_copy(self, phdisp));
2955 skip_if(visir_polynomial_shift_1d_from_correlation
2956 (self, xc_vector, (irplib_base_spectrum_model*) pmymodel,
2957 visir_spectro_fill, hsize, subres, doplot, pxc, pdelta));
2958 *pdidshift = CPL_TRUE;
2961 shifted = cpl_polynomial_duplicate(self);
2963 if (irplib_polynomial_find_1d_from_correlation_all
2964 (self, fitdeg, xc_vector, 1, clines,
2965 (irplib_base_spectrum_model*)pmymodel,
2966 visir_spectro_fill, pixtol, pixstep,
2967 hsize, maxite, maxfail, maxcont, doplot, pxc) || *pxc <= 0.0) {
2969 irplib_error_recover(prestate,
"Could not re-optimize %d "
2970 "coefficients, keeping shifted", fitdeg);
2971 skip_if(cpl_polynomial_copy(self, shifted));
2976 cpl_size clow = 0, chigh = 0;
2978 if (resol == VISIR_SPC_R_LRP) {
2982 cpl_vector * xc_vector_cut = cpl_vector_extract(xc_vector, clow,
2983 cpl_vector_get_size(xc_vector)
2985 cpl_polynomial_shift_1d(self, 0, clow);
2986 skip_if(visir_polynomial_shift_1d_from_correlation
2987 (self, xc_vector_cut, (irplib_base_spectrum_model*) pmymodel,
2988 visir_spectro_fill, hsize, subres, doplot, pxc, pdelta));
2989 cpl_polynomial_shift_1d(self, 0, -clow);
2990 cpl_vector_delete(xc_vector_cut);
2991 *pdidshift = CPL_TRUE;
2994 error_if (*pxc <= 0.0, CPL_ERROR_DATA_NOT_FOUND,
"Atmospheric and Model "
2995 "Spectra have non-positive cross-correlation (%g pixel shift): "
2996 "%g", *pdelta, *pxc);
3000 cpl_polynomial_delete(shifted);
3002 return cpl_error_get_code();
3030static cpl_error_code
3031visir_polynomial_shift_1d_from_correlation(cpl_polynomial * self,
3032 const cpl_vector * obs,
3033 irplib_base_spectrum_model * model,
3034 cpl_error_code (*filler)
3036 const cpl_polynomial *,
3037 irplib_base_spectrum_model *),
3041 double * pxc,
double *pshift)
3043 const int nobs = cpl_vector_get_size(obs);
3044 cpl_polynomial * cand = NULL;
3045 cpl_bivector * xcplot = NULL;
3046 double * xcplotx = NULL;
3047 double * xcploty = NULL;
3048 cpl_vector * mspec1d = NULL;
3050 double bestxc = -1.0;
3051 double bestdelta = -1.0;
3055 cpl_ensure_code(self != NULL, CPL_ERROR_NULL_INPUT);
3056 cpl_ensure_code(obs != NULL, CPL_ERROR_NULL_INPUT);
3057 cpl_ensure_code(model != NULL, CPL_ERROR_NULL_INPUT);
3058 cpl_ensure_code(filler != NULL, CPL_ERROR_NULL_INPUT);
3059 cpl_ensure_code(subres > 0, CPL_ERROR_ILLEGAL_INPUT);
3060 cpl_ensure_code(hsize > 0, CPL_ERROR_ILLEGAL_INPUT);
3062 cand = cpl_polynomial_new(1);
3063 mspec1d = cpl_vector_new(2 * hsize + nobs);
3064 vxc = cpl_vector_new(2 * hsize + 1);
3066 xcplot = cpl_bivector_new(subres * (2 * hsize + 1));
3067 xcplotx = cpl_bivector_get_x_data(xcplot);
3068 xcploty = cpl_bivector_get_y_data(xcplot);
3072 for (i = 0; i < subres; i++) {
3073 const double delta = i / (double)subres;
3077 bug_if (cpl_polynomial_copy(cand, self));
3078 bug_if (cpl_polynomial_shift_1d(cand, 0, delta - hsize));
3080 skip_if (filler(mspec1d, cand, model));
3082 ixc = cpl_vector_correlate(vxc, mspec1d, obs);
3083 xc = cpl_vector_get(vxc, ixc);
3087 bestxxc = ixc - hsize;
3088 bestdelta = delta + bestxxc;
3089 cpl_msg_debug(cpl_func,
"Shifting %g = %d + %g pixels (XC=%g)",
3090 bestdelta, bestxxc, delta, bestxc);
3094 for (j = 0; j <= 2 * hsize; j++) {
3095 const double xcj = cpl_vector_get(vxc, j);
3096 xcplotx[i + j * subres] = (double)(j - hsize) + delta;
3097 xcploty[i + j * subres] = xcj;
3102#ifdef IRPLIB_SPC_DUMP
3104 irplib_polynomial_dump_corr_step(self, vxc,
"Shift");
3107 skip_if(cpl_polynomial_shift_1d(self, 0, bestdelta));
3110 cpl_vector_set_size(vxc, 1);
3111 cpl_vector_set_size(mspec1d, nobs);
3112 skip_if (filler(mspec1d, self, model));
3113 bug_if(cpl_vector_correlate(vxc, mspec1d, obs));
3116 char * title = cpl_sprintf(
"t 'Cross-correlation of %d-pixel spectrum "
3117 "(max=%.4g at %g pixel)' w points", nobs,
3118 cpl_vector_get(vxc, 0), bestdelta);
3120 cpl_plot_bivector(
"set grid;set xlabel 'Offset [pixel]';set ylabel "
3121 "'Cross-correlation';", title,
"", xcplot);
3124 irplib_plot_spectrum_and_model(obs, self, model, filler);
3127 cpl_msg_info(cpl_func,
"Shifting %g = %d + %g pixels (XC: %g <=> %g)",
3128 bestdelta, bestxxc, bestdelta - (
double)bestxxc,
3129 cpl_vector_get(vxc, 0), bestxc);
3131 if (pxc != NULL) *pxc = cpl_vector_get(vxc, 0);
3132 if (pshift != NULL) *pshift = bestdelta;
3136 cpl_vector_delete(mspec1d);
3137 cpl_polynomial_delete(cand);
3138 cpl_vector_delete(vxc);
3139 cpl_bivector_delete(xcplot);
3141 return cpl_error_get_code();
3151static cpl_polynomial * visir_spc_phys_lrp(
void)
3153 const double xval[] = {161, 307, 336, 449, 491, 518, 623, 760, 795, 839};
3154 const double yval[] = {8.22e-6, 9.50e-06, 9.660e-06, 10.5e-06, 10.82e-6,
3155 11.e-06, 11.7e-06, 12.54e-06, 12.76e-06,
3158 const cpl_size maxdeg1d = 2;
3160 cpl_polynomial * self = cpl_polynomial_new(1);
3161 const cpl_boolean sampsym = CPL_FALSE;
3162 const size_t nvals =
sizeof(xval)/
sizeof(*xval);
3164 IRPLIB_DIAG_PRAGMA_PUSH_IGN(-Wcast-qual)
3165 cpl_matrix * xmatrix = cpl_matrix_wrap(1, nvals, (
double*)xval);
3166 cpl_vector * yvector = cpl_vector_wrap(nvals, (
double*)yval);
3167 IRPLIB_DIAG_PRAGMA_POP;
3168 cpl_vector * fitres = cpl_vector_new(nvals);
3170 const cpl_error_code error = cpl_polynomial_fit(self, xmatrix, &sampsym,
3172 CPL_FALSE, NULL, &maxdeg1d)
3173 || cpl_vector_fill_polynomial_fit_residual(fitres, yvector, NULL, self,
3176 const double mse = cpl_vector_product(fitres, fitres) / (double)nvals;
3178 (void)cpl_matrix_unwrap(xmatrix);
3179 (void)cpl_vector_unwrap(yvector);
3180 cpl_vector_delete(fitres);
3183 cpl_error_set_where(cpl_func);
3184 cpl_polynomial_delete(self);
3187 cpl_msg_info(cpl_func,
"Fitted %d degree 1D-polynomial to %u "
3188 "wavelengths with a root mean square error [m]: %g",
3189 (
int)maxdeg1d, (
unsigned)nvals, sqrt(mse));
3203static double visir_spc_get_dispersion(
const cpl_polynomial * self,
double xval)
3206 cpl_errorstate prestate = cpl_errorstate_get();
3209 (void)cpl_polynomial_eval_1d(self, xval, &disp);
3211 if (!cpl_errorstate_is_equal(prestate)) {
3212 (void)cpl_error_set_where(cpl_func);
int visir_parameterlist_get_int(const cpl_parameterlist *self, const char *recipe, visir_parameter bitmask)
Retrieve the value of a VISIR integer parameter.
double visir_pfits_get_temp(const cpl_propertylist *self)
The telescope (M1) temperature [Celcius].
double visir_pfits_get_slitwidth(const cpl_propertylist *self)
The slit width in Arcseconds.
double visir_pfits_get_pixscale(const cpl_propertylist *self)
The pixel scale.
double visir_pfits_get_pixspace(const cpl_propertylist *self)
The pixel spacing.
double visir_pfits_get_wlen(const cpl_propertylist *self)
The central wavelength.
const char * visir_pfits_get_resol(const cpl_propertylist *self)
The spectral resolution.