CR2RE Pipeline Reference Manual 1.6.7
hdrl_spectrum_resample.c
1/*
2 * This file is part of the HDRL
3 * Copyright (C) 2017 European Southern Observatory
4 *
5 * This program is free software; you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation; either version 2 of the License, or
8 * (at your option) any later version.
9 *
10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with this program; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19
20#ifdef HAVE_CONFIG_H
21#include <config.h>
22#endif
23
24/*-----------------------------------------------------------------------------
25 Includes
26 -----------------------------------------------------------------------------*/
27
28#include "hdrl_spectrum_resample.h"
29#include "hdrl_spectrum_defs.h"
30#include "hdrl_parameter.h"
31#include "hdrl_utils.h"
32
33#include <math.h>
34#include <stdlib.h>
35#include <string.h>
36
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>
42
43/*-----------------------------------------------------------------------------
44 Private Functions and Data Structures
45 -----------------------------------------------------------------------------*/
46/*Allocate gsl workspace*/
47static inline
48gsl_bspline_workspace * alloc_workspace(int k, int nCoeff,
49 const double * x_dest, const cpl_size sample_len);
50
51/*Filling gsl structs for fitting*/
52static inline int
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);
57
58/*Private function used for the fit*/
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);
65
66/*Private function used for the windowed fit*/
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);
72
73/* finds the closest element to l in arr*/
74static inline cpl_size
75get_closest_lambda(const double *arr, const cpl_size num_els, const double l);
76
77/*See doc for hdrl_spectrum1D_resample*/
78static inline hdrl_spectrum1D *
79resample_internal(const hdrl_spectrum1D * self,
80 const cpl_array * lambdas_dest, const hdrl_parameter * par);
81
82
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);
87
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);
92
93static inline double
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);
97
98static inline double
99get_stop(const double * v, const cpl_size length, const cpl_size i);
100
101static inline double
102get_start(const double * v, const cpl_size i);
103
104
105/*Parameter used for the interpolation*/
106typedef struct {
107 HDRL_PARAMETER_HEAD;
108 hdrl_spectrum1D_interpolation_method method;
109} hdrl_spectrum1D_resample_interpolate_parameter;
110
111static hdrl_parameter_typeobj
112hdrl_spectrum1D_resample_interpolate_parameter_type = {
113 HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTERPOLATE, /* type */
114 (hdrl_alloc *)&cpl_malloc, /* fp_alloc */
115 (hdrl_free *)&cpl_free, /* fp_free */
116 NULL, /* fp_destroy */
117 sizeof(hdrl_spectrum1D_resample_interpolate_parameter), /* obj_size */
118};
119
120/*Parameter used for the integration*/
121typedef struct {
122 HDRL_PARAMETER_HEAD;
123} hdrl_spectrum1D_resample_integrate_parameter;
124
125static hdrl_parameter_typeobj
126hdrl_spectrum1D_resample_integrate_parameter_type = {
127 HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTEGRATE, /* type */
128 (hdrl_alloc *)&cpl_malloc, /* fp_alloc */
129 (hdrl_free *)&cpl_free, /* fp_free */
130 NULL, /* fp_destroy */
131 sizeof(hdrl_spectrum1D_resample_integrate_parameter), /* obj_size */
132};
133
134
135/*getter for the type of the interpolation*/
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);
142
143 const hdrl_spectrum1D_resample_interpolate_parameter * p =
144 (const hdrl_spectrum1D_resample_interpolate_parameter *)par;
145
146 return p->method;
147}
148
149/*Parameter used for the fit*/
150typedef struct {
151 HDRL_PARAMETER_HEAD;
152 int k;
153 int nCoeff;
154 long window;
155 double factor;
156} hdrl_spectrum1D_resample_fit_parameter;
157
158static hdrl_parameter_typeobj
159hdrl_spectrum1D_resample_fit_parameter_type = {
160 HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT, /* type */
161 (hdrl_alloc *)&cpl_malloc, /* fp_alloc */
162 (hdrl_free *)&cpl_free, /* fp_free */
163 NULL, /* fp_destroy */
164 sizeof(hdrl_spectrum1D_resample_fit_parameter), /* obj_size */
165};
166
167/*getter for the k parameter in fit*/
168int hdrl_spectrum1D_resample_fit_parameter_get_k(const hdrl_parameter * par){
169
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);
174
175 const hdrl_spectrum1D_resample_fit_parameter * p =
176 (const hdrl_spectrum1D_resample_fit_parameter *)par;
177 return p->k;
178}
179
180/*getter for the nCoeff parameter in fit*/
181int hdrl_spectrum1D_resample_fit_parameter_get_nCoeff
182(const hdrl_parameter * par){
183
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);
188
189 const hdrl_spectrum1D_resample_fit_parameter * p =
190 (const hdrl_spectrum1D_resample_fit_parameter *)par;
191 return p->nCoeff;
192}
193
194/*getter for the window parameter in fit*/
195long hdrl_spectrum1D_resample_fit_parameter_get_window
196(const hdrl_parameter * par){
197
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);
202
203 const hdrl_spectrum1D_resample_fit_parameter * p =
204 (const hdrl_spectrum1D_resample_fit_parameter *)par;
205 return p->window;
206}
207
208/*getter for the window parameter in fit*/
209long hdrl_spectrum1D_resample_fit_parameter_get_factor
210(const hdrl_parameter * par){
211
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);
216
217 const hdrl_spectrum1D_resample_fit_parameter * p =
218 (const hdrl_spectrum1D_resample_fit_parameter *)par;
219 return p->factor;
220}
221
226/*-----------------------------------------------------------------------------
227 Functions
228 -----------------------------------------------------------------------------*/
229
230/* ---------------------------------------------------------------------------*/
238/* ---------------------------------------------------------------------------*/
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);
244 p->method = method;
245 return (hdrl_parameter*) p;
246}
247
248/* ---------------------------------------------------------------------------*/
254/* ---------------------------------------------------------------------------*/
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;
260}
261
262
263hdrl_parameter * hdrl_spectrum1D_resample_interpolate_parameter_parse_parlist(
264 const cpl_parameterlist * parlist, const char * prefix){
265
266 cpl_ensure(prefix && parlist, CPL_ERROR_NULL_INPUT, NULL);
267
268 /* Get the Method parameter */
269 char * name = hdrl_join_string(".", 2, prefix, "method");
270
271 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
272 const char * value = cpl_parameter_get_string(par);
273 if (value == NULL) {
274 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
275 "Parameter %s not found", name);
276 cpl_free(name);
277 return NULL;
278 }
279
280 hdrl_spectrum1D_interpolation_method met = hdrl_spectrum1D_interp_linear;
281 /* Switch on the methods */
282 if(!strcmp(value, "LINEAR"))
283 met = hdrl_spectrum1D_interp_linear;
284
285 else if(!strcmp(value, "CSPLINE"))
286 met = hdrl_spectrum1D_interp_cspline;
287
288 else if (!strcmp(value, "AKIMA"))
289 met = hdrl_spectrum1D_interp_akima;
290
291 else{
292 cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
293 "Interpolation method %s not found", value);
294 cpl_free(name);
295 return NULL;
296 }
297
298 cpl_free(name);
300}
301
302cpl_parameterlist *
303hdrl_spectrum1D_resample_interpolate_parameter_create_parlist(const char * base_context,
304 const char * prefix, const char * method_def) {
305
306 cpl_ensure(base_context && prefix, CPL_ERROR_NULL_INPUT, NULL);
307 char * name ;
308 cpl_parameterlist * parlist = cpl_parameterlist_new();
309 cpl_parameter * p ;
310 char * context = hdrl_join_string(".", 2, base_context, prefix);
311
312 /* --prefix.method */
313 name = hdrl_join_string(".", 2, context, "method");
314 p = cpl_parameter_new_enum(name, CPL_TYPE_STRING,
315 "Method used for Spectrum1D interpolation", context,
316 method_def, 3, "LINEAR", "CSPLINE", "AKIMA");
317 cpl_free(name) ;
318 name = hdrl_join_string(".", 2, prefix, "method");
319 cpl_parameter_set_alias(p, CPL_PARAMETER_MODE_CLI, name);
320 cpl_parameter_disable(p, CPL_PARAMETER_MODE_ENV);
321 cpl_free(name) ;
322 cpl_parameterlist_append(parlist, p);
323
324 cpl_free(context);
325 return parlist;
326}
327
328
329/* ---------------------------------------------------------------------------*/
337/* ---------------------------------------------------------------------------*/
339 const int k, const int nCoeff){
340
341 hdrl_spectrum1D_resample_fit_parameter * p
342 = (hdrl_spectrum1D_resample_fit_parameter *)
343 hdrl_parameter_new(&hdrl_spectrum1D_resample_fit_parameter_type);
344 p->k = k;
345 p->nCoeff = nCoeff;
346 p->window = 0;
347 p->factor = 1.0;
348 return (hdrl_parameter*) p;
349}
350
351/* ---------------------------------------------------------------------------*/
365/* ---------------------------------------------------------------------------*/
367 const int k, const int nCoeff, const long window, const double factor){
368
369 cpl_ensure(window > 0, CPL_ERROR_ILLEGAL_INPUT, NULL);
370 cpl_ensure(factor >= 1.0, CPL_ERROR_ILLEGAL_INPUT, NULL);
371
372 hdrl_spectrum1D_resample_fit_parameter * p
373 = (hdrl_spectrum1D_resample_fit_parameter *)
374 hdrl_parameter_new(&hdrl_spectrum1D_resample_fit_parameter_type);
375 p->k = k;
376 p->nCoeff = nCoeff;
377 p->window = window;
378 p->factor = factor;
379 return (hdrl_parameter*) p;
380}
381
382cpl_error_code
383hdrl_resample_parameter_verify(const hdrl_parameter * par){
384
385 cpl_ensure_code(par != NULL, CPL_ERROR_NULL_INPUT);
386
387 hdrl_parameter_enum method = hdrl_parameter_get_parameter_enum(par);
388
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);
393
394 return CPL_ERROR_NONE;
395}
396
397/* ---------------------------------------------------------------------------*/
427/* ---------------------------------------------------------------------------*/
428hdrl_spectrum1D *
429hdrl_spectrum1D_resample(const hdrl_spectrum1D * self,
430 const hdrl_spectrum1D_wavelength* waves,
431 const hdrl_parameter* par){
432
433 cpl_ensure(self != NULL, CPL_ERROR_NULL_INPUT, NULL);
434 cpl_ensure(self-> flux != NULL, CPL_ERROR_NULL_INPUT, NULL);
435
436 cpl_ensure(waves != NULL, CPL_ERROR_NULL_INPUT, NULL);
437 cpl_ensure(waves->wavelength != NULL, CPL_ERROR_NULL_INPUT, NULL);
438
439 cpl_ensure(self->wave_scale == waves->scale, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
440
441 if(hdrl_resample_parameter_verify(par)) return NULL;
442
443 /* if the two spectra are already compatible and we not fitting,
444 * simply return a copy of self */
445 hdrl_spectrum1D_wavelength self_w = hdrl_spectrum1D_get_wavelength(self);
446 if(hdrl_spectrum1D_are_spectra_compatible(&self_w, waves) &&
447 hdrl_parameter_get_parameter_enum(par) != HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_FIT)
448 return hdrl_spectrum1D_duplicate(self);
449
450 return resample_internal(self, waves->wavelength, par);
451}
452
453/* ---------------------------------------------------------------------------*/
482/* ---------------------------------------------------------------------------*/
483hdrl_spectrum1D *
484hdrl_spectrum1D_resample_on_array(const hdrl_spectrum1D * self,
485 const cpl_array * waves,
486 const hdrl_parameter* par){
487
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);
492
493 if(hdrl_resample_parameter_verify(par)) return NULL;
494
495 /* if the two spectra are already compatible and we NOT fitting,
496 * simply return a copy of self */
497 hdrl_spectrum1D_wavelength self_w = hdrl_spectrum1D_get_wavelength(self);
498 if(hdrl_parameter_get_parameter_enum(par) == HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTERPOLATE &&
499 hdrl_spectrum1D_are_wavelengths_compatible(self_w.wavelength, waves))
500 return hdrl_spectrum1D_duplicate(self);
501
502 hdrl_spectrum1D * to_ret = resample_internal(self, waves, par);
503
504 return to_ret;
505}
506
507/*-----------------------------------------------------------------------------
508 Private Functions Implementation
509 -----------------------------------------------------------------------------*/
510
511/*wrapper around the allocation function provided by gsl. It maps our enum in
512 * gls types.*/
513static inline gsl_spline * get_interp_spline
514(const hdrl_spectrum1D_interpolation_method method, const int sample_len){
515
516 switch(method){
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);
523 default:
524 cpl_ensure(0, CPL_ERROR_ILLEGAL_INPUT, NULL);
525 }
526 /*To avoid warning*/
527 return NULL;
528}
529
530/* cleanup for the gsl interpolation structs*/
531static inline void
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);
535 *acc_flux = NULL;
536 }
537
538 if(spline_flux != NULL && *spline_flux != NULL){
539 gsl_spline_free (*spline_flux);
540 *spline_flux = NULL;
541 }
542}
543/*Init of the gls interpolation data structures.
544 * If fails, the memory is cleaned-up */
545static inline cpl_error_code
546init_gsl_interpolate(
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){
550
551 *acc_flux = gsl_interp_accel_alloc ();
552
553 cpl_ensure_code(*acc_flux, CPL_ERROR_UNSPECIFIED);
554
555 *spline_flux = get_interp_spline(method, sample_len);
556
557 if(!*spline_flux){
558 cleanup_gsl_interpolate(acc_flux, spline_flux);
559 cpl_ensure_code(0, CPL_ERROR_UNSPECIFIED);
560 }
561
562 int fail = gsl_spline_init(*spline_flux, x, y, sample_len);
563
564 if(fail){
565 cleanup_gsl_interpolate(acc_flux, spline_flux);
566 cpl_ensure_code(0, CPL_ERROR_UNSPECIFIED);
567 }
568
569 return CPL_ERROR_NONE;
570}
571
572/*Wrapper around gsl_spline_eval. If x is outside the interval provided during
573 * init phase we return NAN. GSL would abort the program.*/
574static inline
575double spline_eval_internal
576(const gsl_spline * spline, double x, gsl_interp_accel * a, cpl_boolean * is_rej){
577
578 *is_rej = x < spline->x[0] || x > spline->x[spline->size - 1];
579
580 if (*is_rej){
581 return NAN;
582 }
583
584 return gsl_spline_eval(spline, x, a);
585}
586
587static inline cpl_size
588get_closest_lambda(const double *arr, const cpl_size num_els, const double l){
589
590 cpl_size best_idx = 0;
591 double smallest_diff = fabs(arr[0] - l);
592
593 for(cpl_size i = 1; i < num_els; i++){
594 const double diff = fabs(arr[i] - l);
595
596 if(diff < smallest_diff){
597 smallest_diff = diff;
598 best_idx = i;
599 }
600
601 /*we assume arr sorted*/
602 if(arr[i] >= l) break;
603 }
604 return best_idx;
605}
606
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){
611
612 cpl_size dest_len = cpl_array_get_size(lambdas_dest);
613
614 gsl_interp_accel * acc_flux = NULL;
615 gsl_spline * spline_flux = NULL;
616
617 cpl_error_code fail =
618 init_gsl_interpolate(x, y, sample_len, method, &acc_flux, &spline_flux);
619
620 cpl_ensure_code(fail == CPL_ERROR_NONE, fail);
621
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);
627 if(is_rejected){
628 cpl_image_reject(dest, i + 1, 1);
629 }
630 else{
631 cpl_image_set(dest, i + 1, 1, val);
632 }
633 }
634 cleanup_gsl_interpolate(&acc_flux, &spline_flux);
635
636 return CPL_ERROR_NONE;
637}
638
639static inline double
640get_min(const double * x, cpl_size length){
641 double best = x[0];
642 for(cpl_size i = 1; i < length; ++i){
643 if(x[i] < best) best = x[i];
644 }
645 return best;
646}
647
648static inline double
649get_max(const double * x, cpl_size length){
650 double best = x[0];
651 for(cpl_size i = 1; i < length; ++i){
652 if(x[i] > best) best = x[i];
653 }
654 return best;
655}
656
657static inline
658gsl_bspline_workspace * alloc_workspace(int k, int nCoeff,
659 const double * x_dest, const cpl_size sample_len){
660
661 const int nKnots = nCoeff + 2 - k;
662
663 gsl_bspline_workspace * ws = gsl_bspline_alloc(k, nKnots);
664
665 double lambda_min = get_min(x_dest, sample_len);
666 double lambda_max = get_max(x_dest, sample_len);
667
668 gsl_bspline_knots_uniform(lambda_min, lambda_max, ws);
669 return ws;
670}
671
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){
677
678 const cpl_size dest_size = cpl_array_get_size(lambdas_dest);
679
680 const cpl_size fit_win = (cpl_size)((double)window * factor);
681 const cpl_size extra_samples_for_fit = (fit_win - window) / 2;
682
683 for(cpl_size dest_start = 0; dest_start < dest_size; dest_start += window){
684
685 const cpl_size dest_stop = CPL_MIN(dest_size - 1, dest_start + window - 1);
686
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);
691
692 /*get corresponding indexes in x_raw (we need to be sure that
693 * min_dest_lambda and max_dest_lambda are inside the interval,
694 * hence the +1 and -1)*/
695 cpl_size raw_start =
696 get_closest_lambda(x_raw, sample_len, min_dest_lambda) - 1;
697 cpl_size raw_end =
698 get_closest_lambda(x_raw, sample_len, max_dest_lambda) + 1;
699
700 /*expand the fit window*/
701 raw_end += extra_samples_for_fit;
702 raw_start -= extra_samples_for_fit;
703
704 raw_end = CPL_MIN(sample_len - 1, raw_end);
705 raw_start = CPL_MAX(0, raw_start);
706
707 const cpl_size win_len = raw_end - raw_start + 1;
708
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);
712
713 if(fail) return fail;
714 }
715
716 return CPL_ERROR_NONE;
717}
718
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){
725
726 cpl_ensure_code(sample_len >= nCoeff, CPL_ERROR_INCOMPATIBLE_INPUT);
727
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);
732
733 int fail = fit_matrixes(x_raw, y_raw, sample_len, nCoeff, ws, B, c, cov);
734
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];
738
739 if(!error)
740 {
741 cpl_size dest_len = cpl_array_get_size(lambdas_dest);
742
743 for (cpl_size i = CPL_MAX(0, lambdas_dest_start);
744 i <= CPL_MIN(dest_len - 1, lambdas_dest_stop); i++)
745 {
746 const double x_dest = cpl_array_get(lambdas_dest, i, NULL);
747
748 /*If outside boundaries reject*/
749 if(x_dest < x_raw_min || x_dest > x_raw_max){
750 cpl_image_reject(flux_dest, i + 1, 1);
751 continue;
752 }
753
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);
757
758 cpl_image_set(flux_dest, i + 1, 1, y_dest);
759 }
760 }
761
762 gsl_matrix_free(cov);
763 gsl_vector_free(B);
764 gsl_vector_free(c);
765 gsl_bspline_free(ws);
766
767 return error;
768}
769
770static inline
771cpl_size count_equals_from_i
772(const double * x, const cpl_size i, const cpl_size sample_len){
773
774 cpl_size idx = i;
775 while( (idx < sample_len - 1) && x[idx] == x[idx + 1])
776 idx++;
777
778 return idx - i + 1;
779}
780
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;
784
785 const double delta = *a_d - *b_d;
786 if(delta > 0.0) return 1;
787 if(delta < 0.0) return -1;
788 return 0;
789}
790
791static inline
792double get_median(double * x, const cpl_size first, const cpl_size n){
793
794 double * v = x + first;
795 qsort(v, n, sizeof(x[0]), compare_d);
796
797 if(n % 2 != 0) return v[n / 2];
798
799 return (v[n / 2] + v[(n - 1)/2]) / 2.0;
800}
801
802cpl_size
803hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median(double * x,
804 double * y1, double * y2, cpl_size sample_len){
805
806 for(cpl_size i = 0; i < sample_len - 1; i++){
807
808 const cpl_size n_equals = count_equals_from_i(x, i, sample_len);
809
810 if(n_equals <= 1) continue;
811
812 y1[i] = get_median(y1, i, n_equals);
813 y2[i] = get_median(y2, i, n_equals);
814
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;
817
818 if(num_bytes > 0){
819 /* using memmove because memory overlaps*/
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);
823 }
824 /* one of the n_equals elements is survived (it is median now), hence the
825 * -1*/
826 sample_len -= n_equals - 1;
827 }
828
829 return sample_len;
830}
831
832
833static inline int
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){
838
839 gsl_matrix* X = gsl_matrix_alloc(sample_len, nCoeff);
840
841 for(cpl_size i = 0; i < sample_len; ++i){
842 int fail = gsl_bspline_eval(x_raw[i], B, ws);
843 if(fail) continue;
844 for(cpl_size j = 0; j < nCoeff; j++){
845 gsl_matrix_set(X, i, j, gsl_vector_get(B, j));
846 }
847 }
848
849 double chisq = 0.0;
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);
853
854 const int error = gsl_multifit_linear(X, &y.vector, c, cov, &chisq, mw);
855
856 gsl_multifit_linear_free(mw);
857 gsl_matrix_free(X);
858
859 return error;
860}
861
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){
867
868 cpl_error_code fail = CPL_ERROR_NONE;
869
870 if(interpolate){
871 hdrl_spectrum1D_interpolation_method int_met =
872 hdrl_spectrum1D_resample_interpolate_parameter_get_method(par);
873
874 fail = fill_cpl_image_with_interpolation(x, y, sample_len, int_met,
875 lambdas_dest, hdrl_image_get_image(flux_dest));
876 } else {
877
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);
882
883 if(window1 == 0){
884 fail = hdrl_spectrum1D_bspline_fit_internal(x, y,
885 sample_len, lambdas_dest, 0,
886 cpl_array_get_size(lambdas_dest) - 1,
887 hdrl_image_get_image(flux_dest), k, nCoeff);
888 } else{
889 fail = hdrl_spectrum1D_fit_windowed_internal(x, y,
890 sample_len, lambdas_dest,
891 hdrl_image_get_image(flux_dest), k,
892 nCoeff, window1, factor);
893 }
894 }
895
896 if(!fail){
897 /* linear interpolate variance */
898 fill_cpl_image_with_interpolation(x, y_var, sample_len,
899 hdrl_spectrum1D_interp_linear,
900 lambdas_dest, hdrl_image_get_error(flux_dest));
901 /* convert variance into error taking the square root*/
902 cpl_image_power(hdrl_image_get_error(flux_dest), 0.5);
903 }
904
905 return fail;
906}
907
908static inline double
909get_start(const double * v, const cpl_size i){
910
911 if(i <= 0){
912 return v[0];
913 }
914
915 const double post = v[i];
916 const double pre = v[i - 1];
917
918 return (post + pre) / 2.0;
919}
920
921static inline double
922get_stop(const double * v, const cpl_size length, const cpl_size i){
923
924 if(i >= length - 1){
925 return v[length - 1];
926 }
927
928 const double post = v[i + 1];
929 const double pre = v[i];
930
931 return (post + pre) / 2.0;
932}
933
934
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){
939 /*edge case, source cover a smaller interval than dest, hence starting and/or
940 * ending bins in destination might be NAN*/
941 const double source_lower_bound = get_start(source, 0);
942 const double source_upper_bound = get_stop(source, source_length, source_length - 1);
943
944 /* destination starts before the start of the first bin*/
945 if(dest_start < source_lower_bound) return CPL_TRUE;
946 /* destination ends after the end of the last bin*/
947 if(dest_stop > source_upper_bound) return CPL_TRUE;
948
949 return CPL_FALSE;
950}
951
952
953static inline double
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){
957
958 double start_source = NAN;
959 double stop_source = NAN;
960 double val = 0.0;
961
962 if(is_destination_outside_source_spectrum(x, sample_len, start_dest, stop_dest))
963 return NAN;
964
965 /*area of destination bin*/
966 const double den = stop_dest - start_dest;
967 *source_idx = CPL_MIN(*source_idx, sample_len - 1);
968
969 do{
970 /*start and stop of source bin*/
971 start_source = get_start(x, *source_idx);
972 stop_source = get_stop(x, sample_len, *source_idx);
973
974 /* There is no intersection, source starts after the end of dest.
975 * We are done with this bin, go back one step so that the next
976 * destination can use the source bin. */
977 if(start_source >= stop_dest) {
978 *source_idx = CPL_MAX(*source_idx - 1, 0);
979 break;
980 }
981
982 /* If stop_source <= start_dest, there is no intersection,
983 * source ends before the start of dest. We can skip this element of source
984 * increase the index so that the next element can be used. We are sure
985 * that there will be an acceptable element because we already checked boundaries
986 * before starting the loop.*/
987 if(stop_source > start_dest){
988 /*Get common area*/
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;
992
993 val += y[*source_idx] * num / den;
994 }
995
996 *source_idx = *source_idx + 1;
997
998 }while(*source_idx < sample_len);
999
1000 return val;
1001}
1002
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){
1007
1008 const cpl_size size_dest = cpl_array_get_size(arg_lambdas_dest);
1009 cpl_bivector * lamb_dest = cpl_bivector_new(size_dest);
1010
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);
1015 }
1016 /*Sort dest lambdas keeping track of the original position*/
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));
1021
1022 for(cpl_size i = 0; i < size_dest; ++i){
1023
1024 const double start_destination = get_start(lambdas_dest_sorted, i);
1025 const double stop_destination = get_stop(lambdas_dest_sorted, size_dest , i);
1026
1027 /*copy index before increase for error computation*/
1028 cpl_size old_idx = source_idx;
1029 const double val = integrate(start_destination, stop_destination,
1030 &source_idx, x, y, sample_len);
1031
1032 const double val_e = sqrt(integrate(start_destination, stop_destination,
1033 &old_idx, x, y_var, sample_len));
1034
1035 const cpl_size dest_idx =
1036 (cpl_size) cpl_vector_get(cpl_bivector_get_y(lamb_dest), i);
1037
1038 if(!isfinite(val) || !isfinite(val_e)) {
1039 hdrl_image_reject(flux_dest, dest_idx + 1, 1);
1040 continue;
1041 }
1042
1043 hdrl_image_set_pixel(flux_dest, dest_idx + 1, 1, (hdrl_value){val, val_e});
1044 }
1045
1046 cpl_bivector_delete(lamb_dest);
1047
1048 return CPL_ERROR_NONE;
1049}
1050
1051static inline hdrl_spectrum1D *
1052resample_internal(const hdrl_spectrum1D * self,
1053 const cpl_array * lambdas_dest, const hdrl_parameter * par){
1054
1055 const cpl_size f_len = hdrl_spectrum1D_get_size(self);
1056 cpl_size sample_len = 0;
1057
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));
1061
1062 const cpl_boolean reject_bad_pix =
1063 hdrl_parameter_get_parameter_enum(par) != HDRL_PARAMETER_SPECTRUM1D_RESAMPLE_INTEGRATE;
1064
1065 for(cpl_size i = 0; i < f_len; i++){
1066
1067 int rej = 0;
1068 const hdrl_value v = hdrl_spectrum1D_get_flux_value(self, i, &rej);
1069 rej = rej || !isfinite((double)v.data) || !isfinite((double)v.error);
1070
1071 if(reject_bad_pix && rej) continue;
1072
1073 y[sample_len] = rej ? NAN : v.data;
1074 /*We interpolate the VARIANCE*/
1075 y_var[sample_len] = rej ? NAN : pow(v.error, 2.0);
1076 x[sample_len] = hdrl_spectrum1D_get_wavelength_value(self, i, NULL);
1077 sample_len++;
1078 }
1079
1080 if(sample_len == 0){
1081 cpl_free(x);
1082 cpl_free(y);
1083 cpl_free(y_var);
1084 cpl_ensure(CPL_FALSE, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
1085 }
1086
1087 /* If wavelengths are not strictly increasing, sort and collapse in case of
1088 * duplicate wavelengths. */
1089 if(!hdrl_is_strictly_monotonic_increasing(x, sample_len)){
1090 /* we need to sort the three arrays so x is strictly increasing*/
1091 hdrl_sort_on_x(x, y, y_var, sample_len, CPL_FALSE);
1092 /* remove duplicated lambdas and substitute flux with median*/
1093 sample_len =
1094 hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1095 (x, y, y_var, sample_len);
1096 }
1097
1098 if(sample_len == 0){
1099 cpl_free(x);
1100 cpl_free(y);
1101 cpl_free(y_var);
1102 cpl_ensure(CPL_FALSE, CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
1103 }
1104
1105 cpl_size dest_len = cpl_array_get_size(lambdas_dest);
1106 hdrl_image * flux_dest = hdrl_image_new(dest_len, 1);
1107
1108 cpl_error_code fail = CPL_ERROR_NONE;
1109
1110 hdrl_parameter_enum method = hdrl_parameter_get_parameter_enum(par);
1111
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);
1117 }
1118 else{
1119 fail = integrate_internal(x, y, y_var, sample_len, lambdas_dest,
1120 flux_dest);
1121 }
1122 cpl_free(x);
1123 cpl_free(y);
1124 cpl_free(y_var);
1125
1126 hdrl_spectrum1D * to_ret = NULL;
1127
1128 if(fail == CPL_ERROR_NONE)
1129 {
1130 const cpl_image * img = hdrl_image_get_image_const(flux_dest);
1131 const cpl_image * img_e = hdrl_image_get_error_const(flux_dest);
1132 to_ret = hdrl_spectrum1D_create(img, img_e, lambdas_dest,
1133 self->wave_scale);
1134 }
1135 hdrl_image_delete(flux_dest);
1136
1137 cpl_ensure(fail == CPL_ERROR_NONE, fail, NULL);
1138
1139 return to_ret;
1140}
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
Definition: hdrl_image.c:594
cpl_image * hdrl_image_get_error(hdrl_image *himg)
get error as cpl image
Definition: hdrl_image.c:131
const cpl_image * hdrl_image_get_error_const(const hdrl_image *himg)
get error as cpl image
Definition: hdrl_image.c:144
cpl_image * hdrl_image_get_image(hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:105
cpl_error_code hdrl_image_reject(hdrl_image *self, cpl_size xpos, cpl_size ypos)
mark pixel as bad
Definition: hdrl_image.c:427
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:118
hdrl_image * hdrl_image_new(cpl_size nx, cpl_size ny)
create new zero filled hdrl image
Definition: hdrl_image.c:311
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
Definition: hdrl_image.c:379
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.
Definition: hdrl_utils.c:767
cpl_boolean hdrl_is_strictly_monotonic_increasing(const double *x, cpl_size l)
returns CPL_TRUE if x is strictly monotonic increasing
Definition: hdrl_utils.c:745
char * hdrl_join_string(const char *sep_, int n,...)
join strings together
Definition: hdrl_utils.c:812