27#include "hdrl_spectrum_shift.h"
28#include "hdrl_parameter.h"
29#include "hdrl_spectrum_resample.h"
45static inline hdrl_data_t
46hdrl_spectrum1D_compute_min_fit(
const hdrl_spectrum1D * s,
47 const hdrl_parameter * shift_par);
49static inline cpl_bivector * get_win(hdrl_data_t wmin, hdrl_data_t wmax);
52static inline cpl_array *
53convert_to_sorted_array(
const hdrl_spectrum1D * s1);
55static inline hdrl_data_t
56hdrl_shift_fit_parameter_get_fit_half_win(
const hdrl_parameter * par);
57static inline hdrl_data_t
58hdrl_shift_fit_parameter_get_fit_wmax(
const hdrl_parameter * par);
59static inline hdrl_data_t
60hdrl_shift_fit_parameter_get_fit_wmin(
const hdrl_parameter * par);
61static inline hdrl_data_t
62hdrl_shift_fit_parameter_get_range_wmin(
const hdrl_parameter * par);
63static inline hdrl_data_t
64hdrl_shift_fit_parameter_get_range_wmax(
const hdrl_parameter * par);
65static inline hdrl_data_t
66hdrl_shift_fit_parameter_get_wguess(
const hdrl_parameter * par);
68static inline hdrl_spectrum1D *
69get_polyfit_for_slope(
const int degree,
const hdrl_spectrum1D * s,
70 const cpl_array * wlengths);
72static inline hdrl_spectrum1D *
73hdrl_spectrum1D_fit(
const hdrl_spectrum1D * obs,
const int degree,
74 const hdrl_data_t wmin,
const hdrl_data_t wmax);
93hdrl_xcorrelation_result *
95 const hdrl_spectrum1D * s2,
97 const cpl_boolean normalize){
100 cpl_ensure(s1 != NULL, CPL_ERROR_NULL_INPUT, NULL);
101 cpl_ensure(s2 != NULL, CPL_ERROR_NULL_INPUT, NULL);
108 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
113 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
115 cpl_array * f1 = convert_to_sorted_array(s1);
116 cpl_array * f2 = convert_to_sorted_array(s2);
119 half_win, normalize, bin, 0.0005);
121 cpl_array_delete(f1);
122 cpl_array_delete(f2);
131 hdrl_data_t range_wmin;
132 hdrl_data_t range_wmax;
133 hdrl_data_t fit_wmin;
134 hdrl_data_t fit_wmax;
135 hdrl_data_t fit_half_win;
136} hdrl_spectrum1D_shift_parameter;
138static hdrl_parameter_typeobj
139hdrl_shift_fit_parameters_type = {
140 HDRL_PARAMETER_SPECTRUM1D_SHIFT,
141 (hdrl_alloc *)&cpl_malloc,
142 (hdrl_free *)&cpl_free,
144 sizeof(hdrl_spectrum1D_shift_parameter),
163 const hdrl_data_t range_wmin,
const hdrl_data_t range_wmax,
164 const hdrl_data_t fit_wmin,
const hdrl_data_t fit_wmax,
165 const hdrl_data_t fit_half_win){
167 hdrl_spectrum1D_shift_parameter * p
168 = (hdrl_spectrum1D_shift_parameter *)
169 hdrl_parameter_new(&hdrl_shift_fit_parameters_type);
171 p->fit_half_win = fit_half_win;
172 p->fit_wmax = fit_wmax;
173 p->fit_wmin = fit_wmin;
174 p->range_wmax = range_wmax;
175 p->range_wmin = range_wmin;
178 return (hdrl_parameter*) p;
199 const hdrl_parameter * par){
201 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0.0);
203 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
204 == HDRL_PARAMETER_SPECTRUM1D_SHIFT, CPL_ERROR_ILLEGAL_INPUT, 0.0);
206 cpl_ensure(hdrl_shift_fit_parameter_get_range_wmin(par)
207 < hdrl_shift_fit_parameter_get_range_wmax(par), CPL_ERROR_ILLEGAL_INPUT, 0.0);
209 cpl_ensure(hdrl_shift_fit_parameter_get_fit_wmin(par)
210 < hdrl_shift_fit_parameter_get_fit_wmax(par), CPL_ERROR_ILLEGAL_INPUT, 0.0);
212 cpl_ensure(hdrl_shift_fit_parameter_get_range_wmin(par)
213 < hdrl_shift_fit_parameter_get_fit_wmin(par), CPL_ERROR_ILLEGAL_INPUT, 0.0);
215 cpl_ensure(hdrl_shift_fit_parameter_get_range_wmax(par)
216 > hdrl_shift_fit_parameter_get_fit_wmax(par), CPL_ERROR_ILLEGAL_INPUT, 0.0);
218 cpl_bivector * win = get_win(hdrl_shift_fit_parameter_get_range_wmin(par),
219 hdrl_shift_fit_parameter_get_range_wmax(par));
224 const int degree = 4;
226 hdrl_spectrum1D * obs_fitted = hdrl_spectrum1D_fit(obs_sel, degree,
227 hdrl_shift_fit_parameter_get_fit_wmin(par),
228 hdrl_shift_fit_parameter_get_fit_wmax(par));
230 cpl_ensure(obs_fitted != NULL, CPL_ERROR_ILLEGAL_OUTPUT, 0.0);
239 hdrl_data_t min_wlen = hdrl_spectrum1D_compute_min_fit(obs_sel, par);
242 cpl_bivector_delete(win);
243 const hdrl_data_t wguess = hdrl_shift_fit_parameter_get_wguess(par);
244 const hdrl_data_t offset = (min_wlen - wguess)/wguess;
253cpl_size convert_to_matrix_and_vector(cpl_matrix ** x, cpl_vector ** values,
254 const hdrl_spectrum1D * s){
260 double * x_vals = cpl_calloc(sz,
sizeof(
double));
261 double * p_values = cpl_calloc(sz,
sizeof(
double));
263 cpl_size real_sz = 0;
264 for(cpl_size i = 0; i < sz; ++i){
270 p_values[real_sz] = f;
282 *values = cpl_vector_wrap(real_sz, p_values);
283 *x = cpl_matrix_wrap(1, real_sz, x_vals);
287static inline cpl_polynomial *
288polynomial_fit_1d_create(
289 const hdrl_spectrum1D * s,
292 cpl_polynomial * fit1d = cpl_polynomial_new(1);
295 cpl_size loc_deg=(cpl_size)degree;
296 cpl_matrix * samppos = NULL;
297 cpl_vector * values = NULL;
299 const cpl_size x_size = convert_to_matrix_and_vector(&samppos, &values, s);
301 cpl_ensure(x_size > 0, CPL_ERROR_ILLEGAL_OUTPUT, NULL);
303 cpl_vector * fitresidual = cpl_vector_new(x_size);
305 cpl_polynomial_fit(fit1d, samppos, NULL, values, NULL,
306 CPL_FALSE, NULL, &loc_deg);
307 cpl_ensure(!cpl_error_get_code(), cpl_error_get_code(), NULL);
309 if ( x_size > (degree + 1) ) {
310 cpl_vector_fill_polynomial_fit_residual(fitresidual, values, NULL, fit1d,
312 cpl_ensure(!cpl_error_get_code(), cpl_error_get_code(), NULL);
315 cpl_matrix_delete(samppos);
316 cpl_vector_delete(fitresidual);
317 cpl_vector_delete(values);
322static inline hdrl_spectrum1D *
323get_polyfit_for_slope(
const int degree,
const hdrl_spectrum1D * s,
324 const cpl_array * wlengths){
328 cpl_polynomial* pfit = polynomial_fit_1d_create(s, degree);
329 cpl_ensure(pfit != NULL, CPL_ERROR_ILLEGAL_OUTPUT, NULL);
331 const cpl_size sz = cpl_array_get_size(wlengths);
332 cpl_image * new_flux = cpl_image_new(sz, 1, HDRL_TYPE_DATA);
333 for(cpl_size i = 0; i < sz; ++i){
334 double v = cpl_polynomial_eval_1d( pfit, cpl_array_get(wlengths, i, NULL), NULL);
335 cpl_image_set(new_flux, i + 1, 1, v);
340 cpl_polynomial_delete(pfit);
341 cpl_image_delete(new_flux);
347static inline hdrl_spectrum1D *
349 const hdrl_spectrum1D * obs,
const int degree,
const hdrl_data_t wmin,
350 const hdrl_data_t wmax){
352 cpl_bivector * win = get_win(wmin, wmax);
356 cpl_ensure(obs_sel != NULL, CPL_ERROR_ILLEGAL_OUTPUT, NULL);
358 hdrl_spectrum1D * obs_fitted = get_polyfit_for_slope(degree, obs_sel,
361 cpl_ensure(obs_fitted != NULL, CPL_ERROR_ILLEGAL_OUTPUT, NULL);
363 cpl_bivector_delete(win);
368static inline hdrl_data_t
369hdrl_spectrum1D_compute_min_fit(
const hdrl_spectrum1D * s,
370 const hdrl_parameter * shift_par){
372 const hdrl_data_t wguess = hdrl_shift_fit_parameter_get_wguess(shift_par);
373 const hdrl_data_t fit_half_win =
374 hdrl_shift_fit_parameter_get_fit_half_win(shift_par);
375 cpl_bivector * win = get_win(wguess - fit_half_win, wguess + fit_half_win);
377 hdrl_spectrum1D * s_core =
382 hdrl_spectrum1D * s_resampled =
383 get_polyfit_for_slope(4, s_core, waves.wavelength);
385 cpl_bivector_delete(win);
389 cpl_size x = 0, y = 0;
392 const hdrl_data_t real_min =
399static inline cpl_bivector * get_win(hdrl_data_t wmin, hdrl_data_t wmax){
400 cpl_bivector * v = cpl_bivector_new(1);
402 cpl_vector_set(cpl_bivector_get_x(v), 0, wmin);
403 cpl_vector_set(cpl_bivector_get_y(v), 0, wmax);
410static inline cpl_array *
411convert_to_sorted_array(
const hdrl_spectrum1D * s1){
414 double * flx = cpl_calloc(sz,
sizeof(
double));
415 double * wav = cpl_calloc(sz,
sizeof(
double));
416 double * is_rej = cpl_calloc(sz,
sizeof(
double));
417 for(cpl_size i = 0; i < sz; ++i){
428 cpl_array * to_ret = cpl_array_wrap_double(flx, sz);
430 for(cpl_size i = 0; i < sz; ++i){
431 if(fabs(is_rej[i]) < 1e-4)
continue;
433 cpl_array_set_invalid(to_ret, i);
441static inline hdrl_data_t
442hdrl_shift_fit_parameter_get_fit_half_win(
const hdrl_parameter * par){
443 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0.0);
444 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
445 == HDRL_PARAMETER_SPECTRUM1D_SHIFT, CPL_ERROR_ILLEGAL_INPUT, 0.0);
447 const hdrl_spectrum1D_shift_parameter * p =
448 (
const hdrl_spectrum1D_shift_parameter *)par;
449 return p->fit_half_win;
452static inline hdrl_data_t
453hdrl_shift_fit_parameter_get_fit_wmax(
const hdrl_parameter * par){
455 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0.0);
456 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
457 == HDRL_PARAMETER_SPECTRUM1D_SHIFT, CPL_ERROR_ILLEGAL_INPUT, 0.0);
459 const hdrl_spectrum1D_shift_parameter * p =
460 (
const hdrl_spectrum1D_shift_parameter *)par;
464static inline hdrl_data_t
465hdrl_shift_fit_parameter_get_fit_wmin(
const hdrl_parameter * par){
467 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0.0);
468 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
469 == HDRL_PARAMETER_SPECTRUM1D_SHIFT, CPL_ERROR_ILLEGAL_INPUT, 0.0);
471 const hdrl_spectrum1D_shift_parameter * p =
472 (
const hdrl_spectrum1D_shift_parameter *)par;
476static inline hdrl_data_t
477hdrl_shift_fit_parameter_get_range_wmin(
const hdrl_parameter * par){
479 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0.0);
480 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
481 == HDRL_PARAMETER_SPECTRUM1D_SHIFT, CPL_ERROR_ILLEGAL_INPUT, 0.0);
483 const hdrl_spectrum1D_shift_parameter * p =
484 (
const hdrl_spectrum1D_shift_parameter *)par;
485 return p->range_wmin;
488static inline hdrl_data_t
489hdrl_shift_fit_parameter_get_range_wmax(
const hdrl_parameter * par){
491 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0.0);
492 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
493 == HDRL_PARAMETER_SPECTRUM1D_SHIFT, CPL_ERROR_ILLEGAL_INPUT, 0.0);
495 const hdrl_spectrum1D_shift_parameter * p =
496 (
const hdrl_spectrum1D_shift_parameter *)par;
497 return p->range_wmax;
500static inline hdrl_data_t
501hdrl_shift_fit_parameter_get_wguess(
const hdrl_parameter * par){
503 cpl_ensure(par != NULL, CPL_ERROR_NULL_INPUT, 0.0);
504 cpl_ensure(hdrl_parameter_get_parameter_enum(par)
505 == HDRL_PARAMETER_SPECTRUM1D_SHIFT, CPL_ERROR_ILLEGAL_INPUT, 0.0);
507 const hdrl_spectrum1D_shift_parameter * p =
508 (
const hdrl_spectrum1D_shift_parameter *)par;
hdrl_xcorrelation_result * hdrl_compute_offset_gaussian(const cpl_array *arr1, const cpl_array *arr2, const cpl_size half_win, const cpl_boolean normalize, const double bin, const double wrange)
Calculate gaussian fit on cross-correlation, does a second fitting for refinement.
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
hdrl_parameter * hdrl_spectrum1D_shift_fit_parameter_create(const hdrl_data_t wguess, const hdrl_data_t range_wmin, const hdrl_data_t range_wmax, const hdrl_data_t fit_wmin, const hdrl_data_t fit_wmax, const hdrl_data_t fit_half_win)
The function create a hdrl_spectrum1D_shift_parameter to be used in hdrl_spectrum1D_compute_shift_fit...
const hdrl_image * hdrl_spectrum1D_get_flux(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter flux
hdrl_data_t hdrl_spectrum1D_compute_shift_fit(const hdrl_spectrum1D *obs, const hdrl_parameter *par)
The function compute the shift due to radial velocity. If wguess is the reference line and wfound is ...
cpl_error_code hdrl_spectrum1D_div_spectrum(hdrl_spectrum1D *self, const hdrl_spectrum1D *other)
divide one spectrum by another spectrum
cpl_size hdrl_spectrum1D_get_size(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter for size
void hdrl_spectrum1D_delete(hdrl_spectrum1D **p_self)
hdrl_spectrum1D destructor
hdrl_xcorrelation_result * hdrl_spectrum1D_compute_shift_xcorrelation(const hdrl_spectrum1D *s1, const hdrl_spectrum1D *s2, cpl_size half_win, const cpl_boolean normalize)
The function computes the shift between the two spectra.
cpl_boolean hdrl_spectrum1D_is_uniformly_sampled(const hdrl_spectrum1D *self, double *bin)
checks if the spectrum is defined on uniformly sampled wavelengths.
cpl_error_code hdrl_spectrum1D_add_scalar(hdrl_spectrum1D *self, hdrl_value scalar_operator)
computes the elementwise addition of a spectrum by a scalar, the self parameter is modified
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_spectrum1D * hdrl_spectrum1D_select_wavelengths(const hdrl_spectrum1D *self, const cpl_bivector *windows, const cpl_boolean is_internal)
the function selects or discards flux values according to whether the value of the corresponding wave...
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.