28#include "hdrl_spectrum_resample.h"
29#include "hdrl_spectrum_defs.h"
30#include "hdrl_parameter.h"
31#include "hdrl_utils.h"
37#include <gsl/gsl_blas.h>
38#include <gsl/gsl_errno.h>
39#include <gsl/gsl_spline.h>
40#include <gsl/gsl_bspline.h>
41#include <gsl/gsl_multifit.h>
48gsl_bspline_workspace * alloc_workspace(
int k,
int nCoeff,
49 const double * x_dest,
const cpl_size sample_len);
53fit_matrixes(
const double * x_raw,
const double * y_raw,
54 const cpl_size sample_len,
const int nCoeff,
55 gsl_bspline_workspace * ws, gsl_vector* B,
56 gsl_vector* c, gsl_matrix * cov);
59static inline cpl_error_code
60hdrl_spectrum1D_bspline_fit_internal
61(
const double* x_raw,
const double* y_raw,
const cpl_size sample_len,
62const cpl_array * lambdas_dest,
const cpl_size lambdas_dest_start,
63const cpl_size lambdas_dest_stop, cpl_image * flux_dest,
64const int k,
const int nCoeff);
67static inline cpl_error_code
68hdrl_spectrum1D_fit_windowed_internal
69(
const double* x_raw,
const double* y_raw,
const cpl_size sample_len,
70 const cpl_array * lambdas_dest, cpl_image * flux_dest,
const int k,
71 const int nCoeff,
const long window,
const double factor);
75get_closest_lambda(
const double *arr,
const cpl_size num_els,
const double l);
78static inline hdrl_spectrum1D *
79resample_internal(
const hdrl_spectrum1D * self,
80 const cpl_array * lambdas_dest,
const hdrl_parameter * par);
83static inline cpl_error_code
84integrate_internal(
double * x,
const double * y,
const double * y_var,
85 const cpl_size sample_len,
const cpl_array * arg_lambdas_dest,
86 hdrl_image * flux_dest);
88static inline cpl_boolean
89is_destination_outside_source_spectrum(
const double * source,
90 const cpl_size source_length,
const double dest_start,
91 const double dest_stop);
94integrate(
const double start_dest,
const double stop_dest,
95 cpl_size * source_idx,
const double * x,
const double * y,
96 const cpl_size sample_len);
99get_stop(
const double * v,
const cpl_size length,
const cpl_size i);
102get_start(
const double * v,
const cpl_size i);
108 hdrl_spectrum1D_interpolation_method method;
109} hdrl_spectrum1D_resample_interpolate_parameter;
111static hdrl_parameter_typeobj
112hdrl_spectrum1D_resample_interpolate_parameter_type = {
113 HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTERPOLATE,
114 (hdrl_alloc *)&cpl_malloc,
115 (hdrl_free *)&cpl_free,
117 sizeof(hdrl_spectrum1D_resample_interpolate_parameter),
123} hdrl_spectrum1D_resample_integrate_parameter;
125static hdrl_parameter_typeobj
126hdrl_spectrum1D_resample_integrate_parameter_type = {
127 HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTEGRATE,
128 (hdrl_alloc *)&cpl_malloc,
129 (hdrl_free *)&cpl_free,
131 sizeof(hdrl_spectrum1D_resample_integrate_parameter),
136hdrl_spectrum1D_interpolation_method
137hdrl_spectrum1D_resample_interpolate_parameter_get_method(
const hdrl_parameter * par){
138 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0);
139 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
140 == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTERPOLATE,
141 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
143 const hdrl_spectrum1D_resample_interpolate_parameter * p =
144 (
const hdrl_spectrum1D_resample_interpolate_parameter *)par;
156} hdrl_spectrum1D_resample_fit_parameter;
158static hdrl_parameter_typeobj
159hdrl_spectrum1D_resample_fit_parameter_type = {
160 HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT,
161 (hdrl_alloc *)&cpl_malloc,
162 (hdrl_free *)&cpl_free,
164 sizeof(hdrl_spectrum1D_resample_fit_parameter),
168int hdrl_spectrum1D_resample_fit_parameter_get_k(
const hdrl_parameter * par){
170 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0);
171 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
172 == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT,
173 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
175 const hdrl_spectrum1D_resample_fit_parameter * p =
176 (
const hdrl_spectrum1D_resample_fit_parameter *)par;
181int hdrl_spectrum1D_resample_fit_parameter_get_nCoeff
182(
const hdrl_parameter * par){
184 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0);
185 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
186 == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT,
187 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
189 const hdrl_spectrum1D_resample_fit_parameter * p =
190 (
const hdrl_spectrum1D_resample_fit_parameter *)par;
195long hdrl_spectrum1D_resample_fit_parameter_get_window
196(
const hdrl_parameter * par){
198 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0);
199 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
200 == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT,
201 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
203 const hdrl_spectrum1D_resample_fit_parameter * p =
204 (
const hdrl_spectrum1D_resample_fit_parameter *)par;
209long hdrl_spectrum1D_resample_fit_parameter_get_factor
210(
const hdrl_parameter * par){
212 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0);
213 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
214 == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT,
215 CPL_ERROR_INCOMPATIBLE_INPUT, 0);
217 const hdrl_spectrum1D_resample_fit_parameter * p =
218 (
const hdrl_spectrum1D_resample_fit_parameter *)par;
240(
const hdrl_spectrum1D_interpolation_method method){
241 hdrl_spectrum1D_resample_interpolate_parameter * p
242 = (hdrl_spectrum1D_resample_interpolate_parameter *)
243 hdrl_parameter_new(&hdrl_spectrum1D_resample_interpolate_parameter_type);
245 return (hdrl_parameter*) p;
256 hdrl_spectrum1D_resample_integrate_parameter * p
257 = (hdrl_spectrum1D_resample_integrate_parameter *)
258 hdrl_parameter_new(&hdrl_spectrum1D_resample_integrate_parameter_type);
259 return (hdrl_parameter*) p;
263hdrl_parameter * hdrl_spectrum1D_resample_interpolate_parameter_parse_parlist(
264 const cpl_parameterlist * parlist,
const char * prefix){
266 cpl_ensure(prefix && parlist, CPL_ERROR_NULL_INPUT, NULL);
271 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
272 const char * value = cpl_parameter_get_string(par);
274 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
275 "Parameter %s not found", name);
280 hdrl_spectrum1D_interpolation_method met = hdrl_spectrum1D_interp_linear;
282 if(!strcmp(value,
"LINEAR"))
283 met = hdrl_spectrum1D_interp_linear;
285 else if(!strcmp(value,
"CSPLINE"))
286 met = hdrl_spectrum1D_interp_cspline;
288 else if (!strcmp(value,
"AKIMA"))
289 met = hdrl_spectrum1D_interp_akima;
292 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
293 "Interpolation method %s not found", value);
303hdrl_spectrum1D_resample_interpolate_parameter_create_parlist(
const char * base_context,
304 const char * prefix,
const char * method_def) {
306 cpl_ensure(base_context && prefix, CPL_ERROR_NULL_INPUT, NULL);
308 cpl_parameterlist * parlist = cpl_parameterlist_new();
314 p = cpl_parameter_new_enum(name, CPL_TYPE_STRING,
315 "Method used for Spectrum1D interpolation", context,
316 method_def, 3,
"LINEAR",
"CSPLINE",
"AKIMA");
319 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, name);
320 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
322 cpl_parameterlist_append(parlist, p);
339 const int k,
const int nCoeff){
341 hdrl_spectrum1D_resample_fit_parameter * p
342 = (hdrl_spectrum1D_resample_fit_parameter *)
343 hdrl_parameter_new(&hdrl_spectrum1D_resample_fit_parameter_type);
348 return (hdrl_parameter*) p;
367 const int k,
const int nCoeff,
const long window,
const double factor){
369 cpl_ensure(window > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
370 cpl_ensure(factor >= 1.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
372 hdrl_spectrum1D_resample_fit_parameter * p
373 = (hdrl_spectrum1D_resample_fit_parameter *)
374 hdrl_parameter_new(&hdrl_spectrum1D_resample_fit_parameter_type);
379 return (hdrl_parameter*) p;
383hdrl_resample_parameter_verify(
const hdrl_parameter * par){
385 cpl_ensure_code(par != NULL, CPL_ERROR_NULL_INPUT);
387 hdrl_parameter_enum method = hdrl_parameter_get_parameter_enum(par);
389 cpl_ensure_code(method == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTERPOLATE ||
390 method == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT ||
391 method == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTEGRATE,
392 CPL_ERROR_INCOMPATIBLE_INPUT);
394 return CPL_ERROR_NONE;
430 const hdrl_spectrum1D_wavelength* waves,
431 const hdrl_parameter* par){
433 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
434 cpl_ensure(self-> flux != NULL, CPL_ERROR_NULL_INPUT, NULL);
436 cpl_ensure(waves != NULL, CPL_ERROR_NULL_INPUT, NULL);
437 cpl_ensure(waves->wavelength != NULL, CPL_ERROR_NULL_INPUT, NULL);
439 cpl_ensure(self->wave_scale == waves->scale, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
441 if(hdrl_resample_parameter_verify(par))
return NULL;
447 hdrl_parameter_get_parameter_enum(par) != HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT)
450 return resample_internal(self, waves->wavelength, par);
485 const cpl_array * waves,
486 const hdrl_parameter* par){
488 cpl_ensure(waves != NULL, CPL_ERROR_NULL_INPUT, NULL);
489 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
490 cpl_ensure(self-> flux != NULL, CPL_ERROR_NULL_INPUT, NULL);
491 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, NULL);
493 if(hdrl_resample_parameter_verify(par))
return NULL;
498 if(hdrl_parameter_get_parameter_enum(par) == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTERPOLATE &&
502 hdrl_spectrum1D * to_ret = resample_internal(self, waves, par);
513static inline gsl_spline * get_interp_spline
514(
const hdrl_spectrum1D_interpolation_method method,
const int sample_len){
517 case hdrl_spectrum1D_interp_linear:
518 return gsl_spline_alloc (gsl_interp_linear, sample_len);
519 case hdrl_spectrum1D_interp_cspline:
520 return gsl_spline_alloc (gsl_interp_cspline, sample_len);
521 case hdrl_spectrum1D_interp_akima:
522 return gsl_spline_alloc (gsl_interp_akima, sample_len);
524 cpl_ensure(0, CPL_ERROR_ILLEGAL_INPUT, NULL);
532cleanup_gsl_interpolate(gsl_interp_accel ** acc_flux, gsl_spline ** spline_flux){
533 if(acc_flux != NULL && *acc_flux != NULL){
534 gsl_interp_accel_free (*acc_flux);
538 if(spline_flux != NULL && *spline_flux != NULL){
539 gsl_spline_free (*spline_flux);
545static inline cpl_error_code
547 const double *x,
const double *y,
const cpl_size sample_len,
548 const hdrl_spectrum1D_interpolation_method method,
549 gsl_interp_accel ** acc_flux, gsl_spline ** spline_flux){
551 *acc_flux = gsl_interp_accel_alloc ();
553 cpl_ensure_code(*acc_flux, CPL_ERROR_UNSPECIFIED);
555 *spline_flux = get_interp_spline(method, sample_len);
558 cleanup_gsl_interpolate(acc_flux, spline_flux);
559 cpl_ensure_code(0, CPL_ERROR_UNSPECIFIED);
562 int fail = gsl_spline_init(*spline_flux, x, y, sample_len);
565 cleanup_gsl_interpolate(acc_flux, spline_flux);
566 cpl_ensure_code(0, CPL_ERROR_UNSPECIFIED);
569 return CPL_ERROR_NONE;
575double spline_eval_internal
576(
const gsl_spline * spline,
double x, gsl_interp_accel * a, cpl_boolean * is_rej){
578 *is_rej = x < spline->x[0] || x > spline->x[spline->size - 1];
584 return gsl_spline_eval(spline, x, a);
587static inline cpl_size
588get_closest_lambda(
const double *arr,
const cpl_size num_els,
const double l){
590 cpl_size best_idx = 0;
591 double smallest_diff = fabs(arr[0] - l);
593 for(cpl_size i = 1; i < num_els; i++){
594 const double diff = fabs(arr[i] - l);
596 if(diff < smallest_diff){
597 smallest_diff = diff;
602 if(arr[i] >= l)
break;
607static inline cpl_error_code
608fill_cpl_image_with_interpolation(
const double *x,
const double *y,
609 const cpl_size sample_len,
const hdrl_spectrum1D_interpolation_method method,
610 const cpl_array * lambdas_dest, cpl_image * dest){
612 cpl_size dest_len = cpl_array_get_size(lambdas_dest);
614 gsl_interp_accel * acc_flux = NULL;
615 gsl_spline * spline_flux = NULL;
617 cpl_error_code fail =
618 init_gsl_interpolate(x, y, sample_len, method, &acc_flux, &spline_flux);
620 cpl_ensure_code(fail == CPL_ERROR_NONE, fail);
622 for(cpl_size i = 0; i < dest_len; i++){
623 const double lambda = cpl_array_get(lambdas_dest, i, NULL);
624 cpl_boolean is_rejected = CPL_FALSE;
625 const double val = spline_eval_internal(spline_flux, lambda,
626 acc_flux, &is_rejected);
628 cpl_image_reject(dest, i + 1, 1);
631 cpl_image_set(dest, i + 1, 1, val);
634 cleanup_gsl_interpolate(&acc_flux, &spline_flux);
636 return CPL_ERROR_NONE;
640get_min(
const double * x, cpl_size length){
642 for(cpl_size i = 1; i < length; ++i){
643 if(x[i] < best) best = x[i];
649get_max(
const double * x, cpl_size length){
651 for(cpl_size i = 1; i < length; ++i){
652 if(x[i] > best) best = x[i];
658gsl_bspline_workspace * alloc_workspace(
int k,
int nCoeff,
659 const double * x_dest,
const cpl_size sample_len){
661 const int nKnots = nCoeff + 2 - k;
663 gsl_bspline_workspace * ws = gsl_bspline_alloc(k, nKnots);
665 double lambda_min = get_min(x_dest, sample_len);
666 double lambda_max = get_max(x_dest, sample_len);
668 gsl_bspline_knots_uniform(lambda_min, lambda_max, ws);
672static inline cpl_error_code
673hdrl_spectrum1D_fit_windowed_internal
674(
const double* x_raw,
const double* y_raw,
const cpl_size sample_len,
675 const cpl_array * lambdas_dest, cpl_image * flux_dest,
const int k,
676 const int nCoeff,
const long window,
const double factor){
678 const cpl_size dest_size = cpl_array_get_size(lambdas_dest);
680 const cpl_size fit_win = (cpl_size)((
double)window * factor);
681 const cpl_size extra_samples_for_fit = (fit_win - window) / 2;
683 for(cpl_size dest_start = 0; dest_start < dest_size; dest_start += window){
685 const cpl_size dest_stop = CPL_MIN(dest_size - 1, dest_start + window - 1);
687 const double min_dest_lambda =
688 cpl_array_get(lambdas_dest, dest_start, NULL);
689 const double max_dest_lambda =
690 cpl_array_get(lambdas_dest, dest_stop, NULL);
696 get_closest_lambda(x_raw, sample_len, min_dest_lambda) - 1;
698 get_closest_lambda(x_raw, sample_len, max_dest_lambda) + 1;
701 raw_end += extra_samples_for_fit;
702 raw_start -= extra_samples_for_fit;
704 raw_end = CPL_MIN(sample_len - 1, raw_end);
705 raw_start = CPL_MAX(0, raw_start);
707 const cpl_size win_len = raw_end - raw_start + 1;
709 cpl_error_code fail = hdrl_spectrum1D_bspline_fit_internal(
710 x_raw + raw_start, y_raw + raw_start, win_len,
711 lambdas_dest, dest_start, dest_stop, flux_dest, k, nCoeff);
713 if(fail)
return fail;
716 return CPL_ERROR_NONE;
719static inline cpl_error_code
720hdrl_spectrum1D_bspline_fit_internal
721(
const double* x_raw,
const double* y_raw,
722const cpl_size sample_len,
const cpl_array * lambdas_dest,
723const cpl_size lambdas_dest_start,
const cpl_size lambdas_dest_stop,
724cpl_image * flux_dest,
const int k,
const int nCoeff){
726 cpl_ensure_code(sample_len >= nCoeff, CPL_ERROR_INCOMPATIBLE_INPUT);
728 gsl_vector* B = gsl_vector_alloc(nCoeff);
729 gsl_vector* c = gsl_vector_alloc(nCoeff);
730 gsl_matrix * cov = gsl_matrix_alloc(nCoeff, nCoeff);
731 gsl_bspline_workspace * ws = alloc_workspace(k, nCoeff, x_raw, sample_len);
733 int fail = fit_matrixes(x_raw, y_raw, sample_len, nCoeff, ws, B, c, cov);
735 cpl_error_code error = !fail ? CPL_ERROR_NONE : CPL_ERROR_UNSPECIFIED;
736 const double x_raw_min = x_raw[0];
737 const double x_raw_max = x_raw[sample_len - 1];
741 cpl_size dest_len = cpl_array_get_size(lambdas_dest);
743 for (cpl_size i = CPL_MAX(0, lambdas_dest_start);
744 i <= CPL_MIN(dest_len - 1, lambdas_dest_stop); i++)
746 const double x_dest = cpl_array_get(lambdas_dest, i, NULL);
749 if(x_dest < x_raw_min || x_dest > x_raw_max){
750 cpl_image_reject(flux_dest, i + 1, 1);
754 gsl_bspline_eval(x_dest, B, ws);
755 double y_dest = 0.0, irrelevant = 0.0;
756 gsl_multifit_linear_est(B, c, cov, &y_dest, &irrelevant);
758 cpl_image_set(flux_dest, i + 1, 1, y_dest);
762 gsl_matrix_free(cov);
765 gsl_bspline_free(ws);
771cpl_size count_equals_from_i
772(
const double * x,
const cpl_size i,
const cpl_size sample_len){
775 while( (idx < sample_len - 1) && x[idx] == x[idx + 1])
781static inline int compare_d(
const void * a,
const void * b){
782 const double * a_d = (
const double*)a;
783 const double * b_d = (
const double*)b;
785 const double delta = *a_d - *b_d;
786 if(delta > 0.0)
return 1;
787 if(delta < 0.0)
return -1;
792double get_median(
double * x,
const cpl_size first,
const cpl_size n){
794 double * v = x + first;
795 qsort(v, n,
sizeof(x[0]), compare_d);
797 if(n % 2 != 0)
return v[n / 2];
799 return (v[n / 2] + v[(n - 1)/2]) / 2.0;
803hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median(
double * x,
804 double * y1,
double * y2, cpl_size sample_len){
806 for(cpl_size i = 0; i < sample_len - 1; i++){
808 const cpl_size n_equals = count_equals_from_i(x, i, sample_len);
810 if(n_equals <= 1)
continue;
812 y1[i] = get_median(y1, i, n_equals);
813 y2[i] = get_median(y2, i, n_equals);
815 const cpl_size num_els_to_copy = sample_len - (i + n_equals);
816 const cpl_size num_bytes =
sizeof(double) * num_els_to_copy;
820 memmove( x + i + 1, x + i + n_equals, (
size_t)num_bytes);
821 memmove(y1 + i + 1, y1 + i + n_equals, (
size_t)num_bytes);
822 memmove(y2 + i + 1, y2 + i + n_equals, (
size_t)num_bytes);
826 sample_len -= n_equals - 1;
834fit_matrixes(
const double * x_raw,
const double * y_raw,
835 const cpl_size sample_len,
const int nCoeff,
836 gsl_bspline_workspace * ws, gsl_vector* B,
837 gsl_vector* c, gsl_matrix * cov){
839 gsl_matrix* X = gsl_matrix_alloc(sample_len, nCoeff);
841 for(cpl_size i = 0; i < sample_len; ++i){
842 int fail = gsl_bspline_eval(x_raw[i], B, ws);
844 for(cpl_size j = 0; j < nCoeff; j++){
845 gsl_matrix_set(X, i, j, gsl_vector_get(B, j));
850 gsl_vector_const_view y = gsl_vector_const_view_array(y_raw, sample_len);
851 gsl_multifit_linear_workspace* mw =
852 gsl_multifit_linear_alloc(sample_len, nCoeff);
854 const int error = gsl_multifit_linear(X, &y.vector, c, cov, &chisq, mw);
856 gsl_multifit_linear_free(mw);
862static inline cpl_error_code
863resample_with_interpol_on_variance(
const hdrl_parameter * par,
864 const cpl_boolean interpolate,
const double * x,
const double * y,
865 const double * y_var,
const cpl_size sample_len,
866 const cpl_array * lambdas_dest, hdrl_image * flux_dest){
868 cpl_error_code fail = CPL_ERROR_NONE;
871 hdrl_spectrum1D_interpolation_method int_met =
872 hdrl_spectrum1D_resample_interpolate_parameter_get_method(par);
874 fail = fill_cpl_image_with_interpolation(x, y, sample_len, int_met,
878 int k = hdrl_spectrum1D_resample_fit_parameter_get_k(par);
879 int nCoeff = hdrl_spectrum1D_resample_fit_parameter_get_nCoeff(par);
880 long window1 = hdrl_spectrum1D_resample_fit_parameter_get_window(par);
881 double factor= hdrl_spectrum1D_resample_fit_parameter_get_factor(par);
884 fail = hdrl_spectrum1D_bspline_fit_internal(x, y,
885 sample_len, lambdas_dest, 0,
886 cpl_array_get_size(lambdas_dest) - 1,
889 fail = hdrl_spectrum1D_fit_windowed_internal(x, y,
890 sample_len, lambdas_dest,
892 nCoeff, window1, factor);
898 fill_cpl_image_with_interpolation(x, y_var, sample_len,
899 hdrl_spectrum1D_interp_linear,
909get_start(
const double * v,
const cpl_size i){
915 const double post = v[i];
916 const double pre = v[i - 1];
918 return (post + pre) / 2.0;
922get_stop(
const double * v,
const cpl_size length,
const cpl_size i){
925 return v[length - 1];
928 const double post = v[i + 1];
929 const double pre = v[i];
931 return (post + pre) / 2.0;
935static inline cpl_boolean
936is_destination_outside_source_spectrum(
const double * source,
937 const cpl_size source_length,
const double dest_start,
938 const double dest_stop){
941 const double source_lower_bound = get_start(source, 0);
942 const double source_upper_bound = get_stop(source, source_length, source_length - 1);
945 if(dest_start < source_lower_bound)
return CPL_TRUE;
947 if(dest_stop > source_upper_bound)
return CPL_TRUE;
954integrate(
const double start_dest,
const double stop_dest,
955 cpl_size * source_idx,
const double * x,
const double * y,
956 const cpl_size sample_len){
958 double start_source = NAN;
959 double stop_source = NAN;
962 if(is_destination_outside_source_spectrum(x, sample_len, start_dest, stop_dest))
966 const double den = stop_dest - start_dest;
967 *source_idx = CPL_MIN(*source_idx, sample_len - 1);
971 start_source = get_start(x, *source_idx);
972 stop_source = get_stop(x, sample_len, *source_idx);
977 if(start_source >= stop_dest) {
978 *source_idx = CPL_MAX(*source_idx - 1, 0);
987 if(stop_source > start_dest){
989 const double common_slice_start = CPL_MAX(start_source, start_dest);
990 const double common_slice_stop = CPL_MIN(stop_source, stop_dest);
991 const double num = common_slice_stop - common_slice_start;
993 val += y[*source_idx] * num / den;
996 *source_idx = *source_idx + 1;
998 }
while(*source_idx < sample_len);
1003static inline cpl_error_code
1004integrate_internal(
double * x,
const double * y,
const double * y_var,
1005 const cpl_size sample_len,
const cpl_array * arg_lambdas_dest,
1006 hdrl_image * flux_dest){
1008 const cpl_size size_dest = cpl_array_get_size(arg_lambdas_dest);
1009 cpl_bivector * lamb_dest = cpl_bivector_new(size_dest);
1011 for(cpl_size i = 0; i < size_dest ; ++i){
1012 const double val = cpl_array_get(arg_lambdas_dest, i, NULL);
1013 cpl_vector_set(cpl_bivector_get_x(lamb_dest), i, val);
1014 cpl_vector_set(cpl_bivector_get_y(lamb_dest), i, i);
1017 cpl_bivector_sort(lamb_dest, lamb_dest, CPL_SORT_ASCENDING, CPL_SORT_BY_X);
1018 cpl_size source_idx = 0;
1019 const double * lambdas_dest_sorted =
1020 cpl_vector_get_data_const(cpl_bivector_get_x(lamb_dest));
1022 for(cpl_size i = 0; i < size_dest; ++i){
1024 const double start_destination = get_start(lambdas_dest_sorted, i);
1025 const double stop_destination = get_stop(lambdas_dest_sorted, size_dest , i);
1028 cpl_size old_idx = source_idx;
1029 const double val = integrate(start_destination, stop_destination,
1030 &source_idx, x, y, sample_len);
1032 const double val_e = sqrt(integrate(start_destination, stop_destination,
1033 &old_idx, x, y_var, sample_len));
1035 const cpl_size dest_idx =
1036 (cpl_size) cpl_vector_get(cpl_bivector_get_y(lamb_dest), i);
1038 if(!isfinite(val) || !isfinite(val_e)) {
1046 cpl_bivector_delete(lamb_dest);
1048 return CPL_ERROR_NONE;
1051static inline hdrl_spectrum1D *
1052resample_internal(
const hdrl_spectrum1D * self,
1053 const cpl_array * lambdas_dest,
const hdrl_parameter * par){
1056 cpl_size sample_len = 0;
1058 double * y = cpl_calloc(f_len,
sizeof(
double));
1059 double * y_var = cpl_calloc(f_len,
sizeof(
double));
1060 double * x = cpl_calloc(f_len,
sizeof(
double));
1062 const cpl_boolean reject_bad_pix =
1063 hdrl_parameter_get_parameter_enum(par) != HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTEGRATE;
1065 for(cpl_size i = 0; i < f_len; i++){
1069 rej = rej || !isfinite((
double)v.data) || !isfinite((
double)v.error);
1071 if(reject_bad_pix && rej)
continue;
1073 y[sample_len] = rej ? NAN : v.data;
1075 y_var[sample_len] = rej ? NAN : pow(v.error, 2.0);
1080 if(sample_len == 0){
1084 cpl_ensure(CPL_FALSE, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
1094 hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1095 (x, y, y_var, sample_len);
1098 if(sample_len == 0){
1102 cpl_ensure(CPL_FALSE, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
1105 cpl_size dest_len = cpl_array_get_size(lambdas_dest);
1108 cpl_error_code fail = CPL_ERROR_NONE;
1110 hdrl_parameter_enum method = hdrl_parameter_get_parameter_enum(par);
1112 if(method != HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTEGRATE){
1113 const cpl_boolean is_interpolation =
1114 method == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTERPOLATE;
1115 fail = resample_with_interpol_on_variance(par, is_interpolation, x, y,
1116 y_var, sample_len, lambdas_dest, flux_dest);
1119 fail = integrate_internal(x, y, y_var, sample_len, lambdas_dest,
1126 hdrl_spectrum1D * to_ret = NULL;
1128 if(fail == CPL_ERROR_NONE)
1137 cpl_ensure(fail == CPL_ERROR_NONE, fail, NULL);
cpl_error_code hdrl_image_set_pixel(hdrl_image *self, cpl_size xpos, cpl_size ypos, hdrl_value value)
set pixel values of hdrl_image
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
const cpl_image * hdrl_image_get_error_const(const 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
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
hdrl_image * hdrl_image_new(cpl_size nx, cpl_size ny)
create new zero filled hdrl image
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
hdrl_spectrum1D * hdrl_spectrum1D_resample(const hdrl_spectrum1D *self, const hdrl_spectrum1D_wavelength *waves, const hdrl_parameter *par)
resample a hdrl_spectrum1D on the wavelengths contained in waves
hdrl_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_parameter * hdrl_spectrum1D_resample_fit_parameter_create(const int k, const int nCoeff)
constructor for the hdrl_parameter in the case of interpolation
cpl_boolean hdrl_spectrum1D_are_wavelengths_compatible(const cpl_array *w1, const cpl_array *w2)
checks if two wavelengths array are defined on the same wavelengths.
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_wavelength hdrl_spectrum1D_get_wavelength(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter for wavelengths
hdrl_data_t hdrl_spectrum1D_get_wavelength_value(const hdrl_spectrum1D *self, int idx, int *rej)
hdrl_spectrum1D getter for a wavelength value
hdrl_parameter * hdrl_spectrum1D_resample_fit_windowed_parameter_create(const int k, const int nCoeff, const long window, const double factor)
constructor for the hdrl_parameter in the case of interpolation
hdrl_parameter * hdrl_spectrum1D_resample_integrate_parameter_create(void)
constructor for the hdrl_parameter in the case of integration
cpl_boolean hdrl_spectrum1D_are_spectra_compatible(const hdrl_spectrum1D_wavelength *s1, const hdrl_spectrum1D_wavelength *s2)
checks if two spectrum wavelengths are equal.
hdrl_value hdrl_spectrum1D_get_flux_value(const hdrl_spectrum1D *self, int idx, int *rej)
hdrl_spectrum1D getter for a flux value
void hdrl_sort_on_x(double *x, double *y1, double *y2, const cpl_size sample_len, const cpl_boolean sort_decreasing)
sort in increasing or decreasing order x. Keep aligned with y1 and y2.
cpl_boolean hdrl_is_strictly_monotonic_increasing(const double *x, cpl_size l)
returns CPL_TRUE if x is strictly monotonic increasing
char * hdrl_join_string(const char *sep_, int n,...)
join strings together