27#include "hdrl_spectrumlist.h"
28#include "hdrl_spectrum_resample.h"
29#include "hdrl_imagelist.h"
46cpl_error_code check_getter(
const hdrl_spectrum1Dlist * s,
const cpl_size idx);
49void push_back(hdrl_spectrum1Dlist *self, hdrl_spectrum1D * s);
51static inline hdrl_spectrum1D **
52safe_realloc(hdrl_spectrum1D ** current,
const cpl_size old_capacity,
53 const cpl_size new_capacity);
55static inline hdrl_imagelist *
56create_list(hdrl_spectrum1D ** resampled_spectra,
57 const hdrl_spectrum1Dlist * ori_spectra,
58 const cpl_boolean mark_bpm_in_interpolation);
61delete_spectrum1D_array(hdrl_spectrum1D ** resampled_spectra,
62 const cpl_size num_spectra);
64static inline hdrl_image *
65get_padded_flux(
const hdrl_spectrum1D * resampled_spectrum,
66 const hdrl_spectrum1D * ori_spectrum,
67 const cpl_boolean mark_bpm_in_interpolation);
69static inline cpl_boolean
70are_spectra_valid(
const hdrl_spectrum1Dlist * list);
72static inline cpl_error_code
73get_first_error_code(
const cpl_error_code * cds,
const cpl_size length);
75static inline cpl_boolean
76check_scales_are_same(
const hdrl_spectrum1Dlist * list);
79get_wmin_valid(
const hdrl_spectrum1D * s);
82get_wmax_valid(
const hdrl_spectrum1D * s);
85remove_if_neighbors_are_rejected(hdrl_image * flx,
86 const hdrl_spectrum1D * ori_spectrum,
const cpl_array * wlens);
88static inline hdrl_spectrum1D *
89get_interp_bpm(
const hdrl_spectrum1D * s,
const cpl_array * wlens);
91static inline cpl_boolean
92contains(
const hdrl_spectrum1Dlist * list,
const hdrl_spectrum1D * s);
103 hdrl_spectrum1Dlist * to_ret = cpl_calloc(1,
sizeof(*to_ret));
107 to_ret->spectra = NULL;
124 if(l == NULL)
return NULL;
151 to_ret->spectra = self;
152 to_ret->capacity = to_ret->length = sz;
170 cpl_error_code fail = check_getter(self, idx);
171 cpl_ensure(fail == CPL_ERROR_NONE, fail, NULL);
173 return self->spectra[idx];
188const hdrl_spectrum1D *
191 cpl_error_code fail = check_getter(self, idx);
192 cpl_ensure(fail == CPL_ERROR_NONE, fail, NULL);
194 return self->spectra[idx];
217 hdrl_spectrum1D * s,
const cpl_size idx){
219 cpl_ensure_code(self != NULL, CPL_ERROR_NULL_INPUT);
220 cpl_ensure_code(self->length >= 0, CPL_ERROR_ACCESS_OUT_OF_RANGE);
221 cpl_ensure_code(idx <= self->length, CPL_ERROR_ACCESS_OUT_OF_RANGE);
222 cpl_ensure_code(!contains(self, s), CPL_ERROR_ILLEGAL_INPUT);
224 if(idx == self->length){
226 return CPL_ERROR_NONE;
229 hdrl_spectrum1D * old = self->spectra[idx];
232 self->spectra[idx] = s;
234 return CPL_ERROR_NONE;
256 cpl_error_code fail = check_getter(self, idx);
257 cpl_ensure(fail == CPL_ERROR_NONE, fail, NULL);
259 hdrl_spectrum1D * to_ret = self->spectra[idx];
263 for(cpl_size i = idx; i < sz - 1; ++i){
264 self->spectra[i] = self->spectra[i + 1];
268 const cpl_size new_capacity = self->capacity / 2;
271 self->spectra = safe_realloc(self->spectra, self->capacity, new_capacity);
272 self->capacity = new_capacity;
288 if(l == NULL)
return;
293 cpl_free(l->spectra);
311 cpl_ensure(l != NULL, CPL_ERROR_NULL_INPUT, 0);
351 const hdrl_parameter * stacking_par,
352 const cpl_array * wlengths,
const hdrl_parameter * resample_par,
353 const cpl_boolean mark_bpm_in_interpolation,
354 hdrl_spectrum1D ** result, cpl_image ** contrib,
355 hdrl_imagelist ** resampled_and_aligned_fluxes){
357 cpl_ensure_code(are_spectra_valid(list), CPL_ERROR_NULL_INPUT);
358 cpl_ensure_code(wlengths != NULL, CPL_ERROR_NULL_INPUT);
359 cpl_ensure_code(check_scales_are_same(list), CPL_ERROR_ILLEGAL_INPUT);
360 cpl_ensure_code(result != NULL, CPL_ERROR_NULL_INPUT);
361 cpl_ensure_code(resampled_and_aligned_fluxes != NULL, CPL_ERROR_NULL_INPUT);
366 hdrl_spectrum1D ** resampled_spectra =
367 cpl_calloc(num_spectra,
sizeof(hdrl_spectrum1D*));
369 cpl_ensure_code(num_spectra > 0, CPL_ERROR_ILLEGAL_INPUT);
370 cpl_error_code * errors = cpl_calloc(num_spectra,
sizeof(cpl_error_code));
372 HDRL_OMP(omp parallel
for)
373 for(cpl_size i = 0; i < num_spectra; ++i){
375 resampled_spectra[i] =
377 errors[i] = cpl_error_get_code();
381 cpl_error_code err = get_first_error_code(errors, num_spectra);
385 hdrl_imagelist * stack_list =
386 create_list(resampled_spectra, list, mark_bpm_in_interpolation);
388 hdrl_image * stacked_img = NULL;
391 &stacked_img, contrib);
393 *resampled_and_aligned_fluxes = stack_list;
396 const hdrl_spectrum1D_wave_scale scale =
408 delete_spectrum1D_array(resampled_spectra, num_spectra);
417get_wmin_valid(
const hdrl_spectrum1D * s){
420 double wmin = INFINITY;
422 for(cpl_size i = 0; i < sz; ++i){
426 if(w < wmin) wmin = w;
432get_wmax_valid(
const hdrl_spectrum1D * s){
435 double wmax = -INFINITY;
437 for(cpl_size i = 0; i < sz; ++i){
441 if(w > wmax) wmax = w;
446static inline cpl_image *
447get_img_from_bpm(hdrl_spectrum1D_wavelength * w){
449 return cpl_image_new_from_mask(w->bpm);
451 const cpl_size sz = cpl_array_get_size(w->wavelength);
452 return cpl_image_new(sz, 1, CPL_TYPE_INT);
455static inline hdrl_spectrum1D *
456get_interp_bpm(
const hdrl_spectrum1D * s,
const cpl_array * wlens){
459 cpl_image * flx = get_img_from_bpm(&wlen_ori);
461 wlen_ori.wavelength, wlen_ori.scale);
462 cpl_image_delete(flx);
464 hdrl_parameter * par =
466 hdrl_spectrum1D_interp_linear);
475remove_if_neighbors_are_rejected(hdrl_image * flx,
476 const hdrl_spectrum1D * ori_spectrum,
const cpl_array * wlens){
479 hdrl_spectrum1D * interp_bpm = get_interp_bpm(ori_spectrum, wlens);
483 if(v > HDRL_EPS_DATA){
490static inline hdrl_image *
491get_padded_flux(
const hdrl_spectrum1D * resampled_spectrum,
492 const hdrl_spectrum1D * ori_spectrum,
493 const cpl_boolean mark_bpm_in_interpolation){
495 if(resampled_spectrum == NULL)
return NULL;
497 const double wmin = get_wmin_valid(ori_spectrum);
499 const double wmax = get_wmax_valid(ori_spectrum);
501 if(isinf(wmin) || isinf(wmax))
return NULL;
506 const cpl_array * wlens =
510 const double wlen = cpl_array_get(wlens, i, NULL);
511 if(wlen < wmin || wlen > wmax) {
516 if(mark_bpm_in_interpolation)
517 remove_if_neighbors_are_rejected(flx, ori_spectrum, wlens);
522static inline cpl_boolean
523check_scales_are_same(
const hdrl_spectrum1Dlist * list){
524 if(list == NULL)
return CPL_TRUE;
527 if(sz <= 1)
return CPL_TRUE;
529 const hdrl_spectrum1D_wave_scale scale =
532 for(cpl_size i = 1; i < sz; ++i){
533 const hdrl_spectrum1D_wave_scale scale_this =
536 if(scale_this != scale)
return CPL_FALSE;
544static inline hdrl_imagelist *
545create_list(hdrl_spectrum1D ** resampled_spectra,
546 const hdrl_spectrum1Dlist * ori_spectra,
547 const cpl_boolean mark_bpm_in_interpolation){
550 hdrl_image ** images = cpl_calloc(num_spectra,
sizeof(hdrl_image*));
551 cpl_error_code * errs = cpl_calloc(num_spectra,
sizeof(cpl_error_code));
552 HDRL_OMP(omp parallel
for)
553 for(cpl_size i = 0; i < num_spectra; i++){
554 hdrl_image * img = get_padded_flux(resampled_spectra[i],
556 mark_bpm_in_interpolation);
558 errs[i] = cpl_error_get_code();
561 cpl_error_code fail = get_first_error_code(errs, num_spectra);
564 hdrl_imagelist * list = NULL;
570 for(cpl_size i = 0; i < num_spectra; i++){
571 if(images[i] != NULL)
580cpl_error_code check_getter(
const hdrl_spectrum1Dlist *s,
const cpl_size idx){
582 if(s == NULL)
return CPL_ERROR_NULL_INPUT;
584 if(idx < 0)
return CPL_ERROR_ACCESS_OUT_OF_RANGE;
586 if(idx >= s->length)
return CPL_ERROR_ACCESS_OUT_OF_RANGE;
588 return CPL_ERROR_NONE;
592static inline hdrl_spectrum1D **
593safe_realloc(hdrl_spectrum1D ** current,
const cpl_size old_capacity,
594 const cpl_size new_capacity){
597 if(old_capacity == 0 && new_capacity > 0){
598 return cpl_calloc(new_capacity,
sizeof(hdrl_spectrum1D * ));
602 if(new_capacity == 0){
607 hdrl_spectrum1D ** to_ret =
608 cpl_realloc(current,
sizeof(hdrl_spectrum1D * ) * new_capacity);
610 for(cpl_size i = old_capacity + 1; i < new_capacity; ++i){
619void resize_if_needed(hdrl_spectrum1Dlist *self){
621 if(self->length < self->capacity)
return;
623 const cpl_size new_capacity = self->capacity == 0 ? 1 : 2 * self->capacity;
625 self->spectra = safe_realloc(self->spectra, self->capacity,
627 self->capacity = new_capacity;
632void push_back(hdrl_spectrum1Dlist *self, hdrl_spectrum1D * s){
634 resize_if_needed(self);
635 self->spectra[self->length] = s;
641delete_spectrum1D_array(hdrl_spectrum1D ** resampled_spectra,
642 const cpl_size num_spectra){
644 hdrl_spectrum1Dlist * l =
650static inline cpl_boolean
651are_spectra_valid(
const hdrl_spectrum1Dlist * list){
652 if(list == NULL)
return CPL_FALSE;
655 for(cpl_size i = 0; i < length; ++i){
657 if(s == NULL)
return CPL_FALSE;
662static inline cpl_error_code
663get_first_error_code(
const cpl_error_code * cds,
const cpl_size length){
664 for(cpl_size i = 0; i < length; ++i){
665 if(cds[i])
return cds[i];
667 return CPL_ERROR_NONE;
670static inline cpl_boolean
671contains(
const hdrl_spectrum1Dlist * list,
const hdrl_spectrum1D * s){
673 for(cpl_size i = 0; i < sz; ++i){
675 if(s_in == s)
return CPL_TRUE;
hdrl_image * hdrl_image_duplicate(const hdrl_image *himg)
copy hdrl_image
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl 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
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
cpl_error_code hdrl_imagelist_set(hdrl_imagelist *himlist, hdrl_image *himg, cpl_size pos)
Insert an image into an imagelist.
cpl_error_code hdrl_imagelist_collapse(const hdrl_imagelist *himlist, const hdrl_parameter *param, hdrl_image **out, cpl_image **contrib)
collapsing of image list
hdrl_imagelist * hdrl_imagelist_new(void)
Create an empty imagelist.
void hdrl_parameter_delete(hdrl_parameter *obj)
shallow delete of a parameter
const hdrl_image * hdrl_spectrum1D_get_flux(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter flux
cpl_size hdrl_spectrum1Dlist_get_size(const hdrl_spectrum1Dlist *l)
hdrl_spectrum1Dlist getter for size
hdrl_spectrum1Dlist * hdrl_spectrum1Dlist_new(void)
hdrl_spectrum1Dlist default constructor
hdrl_parameter * hdrl_spectrum1D_resample_interpolate_parameter_create(const hdrl_spectrum1D_interpolation_method method)
constructor for the hdrl_parameter in the case of interpolation
hdrl_spectrum1D * hdrl_spectrum1D_duplicate(const hdrl_spectrum1D *self)
hdrl_spectrum1D copy constructor
cpl_size hdrl_spectrum1D_get_size(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter for size
hdrl_spectrum1D * hdrl_spectrum1Dlist_unset(hdrl_spectrum1Dlist *self, const cpl_size idx)
hdrl_spectrum1Dlist remove of the i-th element
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_spectrum1Dlist_get(hdrl_spectrum1Dlist *self, const cpl_size idx)
hdrl_spectrum1Dlist getter of the i-th element
const hdrl_spectrum1D * hdrl_spectrum1Dlist_get_const(const hdrl_spectrum1Dlist *self, const cpl_size idx)
hdrl_spectrum1Dlist getter of the i-th element
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 * hdrl_spectrum1D_resample_on_array(const hdrl_spectrum1D *self, const cpl_array *waves, const hdrl_parameter *par)
resample a hdrl_spectrum1D on the wavelengths contained in waves
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_spectrum1D * hdrl_spectrum1D_create_error_free(const cpl_image *arg_flux, const cpl_array *wavelength, hdrl_spectrum1D_wave_scale scale)
hdrl_spectrum1D constructor in the case of error-free spectrum (i.e. the error on the flux is zero fo...
hdrl_data_t hdrl_spectrum1D_get_wavelength_value(const hdrl_spectrum1D *self, int idx, int *rej)
hdrl_spectrum1D getter for a wavelength value
hdrl_spectrum1Dlist * hdrl_spectrum1Dlist_duplicate(const hdrl_spectrum1Dlist *l)
hdrl_spectrum1Dlist copy-constructor
hdrl_spectrum1Dlist * hdrl_spectrum1Dlist_wrap(hdrl_spectrum1D **self, const cpl_size sz)
hdrl_spectrum1Dlist wrapper
cpl_error_code hdrl_spectrum1Dlist_collapse(const hdrl_spectrum1Dlist *list, const hdrl_parameter *stacking_par, const cpl_array *wlengths, const hdrl_parameter *resample_par, const cpl_boolean mark_bpm_in_interpolation, hdrl_spectrum1D **result, cpl_image **contrib, hdrl_imagelist **resampled_and_aligned_fluxes)
collapsing a hdrl_spectrum1Dlist. The spectra in list are first resampled on the wavelengths wlengths...
hdrl_value hdrl_spectrum1D_get_flux_value(const hdrl_spectrum1D *self, int idx, int *rej)
hdrl_spectrum1D getter for a flux value