CR2RE Pipeline Reference Manual 1.6.10
hdrl_spectrum1d-test.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.h"
29#include "hdrl_spectrumlist.h"
30#include "hdrl_spectrum_resample.h"
31#include "hdrl_DER_SNR.h"
32
33#include <cpl.h>
34#include <math.h>
35#include <time.h>
36#include <stdlib.h>
37
38/*-----------------------------------------------------------------------------
39 Define
40 -----------------------------------------------------------------------------*/
41
42#define HDRL_DELTA_COMPARE_VALUE CPL_MAX(HDRL_EPS_DATA, HDRL_EPS_ERROR) * 1.
43#define HDRL_DELTA_COMPARE_VALUE_ABS CPL_MAX(HDRL_EPS_DATA, HDRL_EPS_ERROR) * 4.
44
45#define LENGTH_ARRAY(a) ((cpl_size)(sizeof((a))/sizeof((a[0]))))
46
47/*-----------------------------------------------------------------------------
48 Function prototipes
49 -----------------------------------------------------------------------------*/
50
51static inline cpl_image * get_random_1d_img(cpl_size length,
52 double min, double max, cpl_type type);
53
54static inline void set_1d_bpm(cpl_image * img);
55
56static inline cpl_array * get_wavelength(cpl_size length, cpl_type type);
57
58static inline cpl_boolean are_hdrl_eq(const hdrl_image* flux_compound,
59 const cpl_image * flux, const cpl_image * flux_e);
60
61static inline cpl_error_code get_error_code_and_reset();
62
63static inline double rand_0_to_1();
64
65static inline hdrl_spectrum1D *
66get_random_spectrum(int length, hdrl_spectrum1D_wave_scale scale);
67
68typedef hdrl_spectrum1D * (* operate_spectra_create)
69 (const hdrl_spectrum1D * self,
70 const hdrl_spectrum1D * other);
71
72typedef cpl_error_code (* operate_spectra)(hdrl_spectrum1D * self,
73 const hdrl_spectrum1D * other);
74
75static inline void test_error_create_func(const hdrl_spectrum1D * s1,
76 const hdrl_spectrum1D * s2, operate_spectra_create f);
77
78static inline void
79test_calc_creat_error(operate_spectra_create f);
80
81static inline void test_calc_error(operate_spectra f);
82
83static inline hdrl_spectrum1D * get_spectrum1D_sin_shuffled
84(cpl_size sz, int start, cpl_boolean add_peak, cpl_array ** unshuffled_lambda);
85
86typedef struct{
87 double lambda_min;
88 double lambda_max;
89 cpl_boolean is_internal;
90}sel_window;
91
92hdrl_spectrum1D * select_window(hdrl_spectrum1D * s, const sel_window w){
93
94 cpl_bivector * vec = cpl_bivector_new(1);
95
96 cpl_vector_set(cpl_bivector_get_x(vec), 0, w.lambda_min);
97 cpl_vector_set(cpl_bivector_get_y(vec), 0, w.lambda_max);
98
99 hdrl_spectrum1D * to_ret = hdrl_spectrum1D_select_wavelengths(s, vec, w.is_internal);
100 cpl_bivector_delete(vec);
101
102 return to_ret;
103}
104
105/*----------------------------------------------------------------------------*/
109/*----------------------------------------------------------------------------*/
110
111void test_spectrum1D_constructor(cpl_type type){
112
113 const cpl_size sz = 40;
114 cpl_image * spectrum1d = get_random_1d_img(sz, 0.0f, 128.0f, type);
115 cpl_image * spectrum1d_error = get_random_1d_img(sz, 0.0f, 1.0f,type);
116 set_1d_bpm(spectrum1d);
117 cpl_array * wavelengths = get_wavelength(sz,type);
118
119 /* test when passing error array */
120 hdrl_spectrum1D * spec;
121
123 spectrum1d, spectrum1d_error,
124 wavelengths, hdrl_spectrum1D_wave_scale_linear);
125
126 cpl_test_eq(hdrl_spectrum1D_get_size(spec), 40);
127 cpl_boolean are_equal = are_hdrl_eq(
128 hdrl_spectrum1D_get_flux(spec), spectrum1d, spectrum1d_error
129 );
130 cpl_test_eq(are_equal, CPL_TRUE);
131
133 cpl_test_null(spec);
134
135 const cpl_size wn = 10;
136
137 /* test when using DER SNR error*/
139 spectrum1d, wn,
140 wavelengths, hdrl_spectrum1D_wave_scale_linear);
141 cpl_test_nonnull(spec);
142
143 /*some casting to make sure that we are always using hdrl_data_t types*/
144 cpl_image * flux_casted = cpl_image_cast(spectrum1d, HDRL_TYPE_DATA);
145 const hdrl_data_t * flux =
146 (const hdrl_data_t *) cpl_image_get_data_const(flux_casted);
147 const cpl_mask *msk_in = cpl_image_get_bpm_const(flux_casted);
148 const cpl_binary *msk_bn = cpl_mask_get_data_const(msk_in);
149
150 cpl_image * noise = estimate_noise_DER_SNR(flux, msk_bn, wavelengths,
151 sz, wn);
152
153 cpl_image_delete(flux_casted);
154
155 cpl_mask * msk = cpl_image_unset_bpm(noise);
156 cpl_mask_delete(cpl_image_set_bpm(spectrum1d, msk));
157
158 {
159 hdrl_spectrum1D * spec2 = hdrl_spectrum1D_create(
160 spectrum1d, noise,
161 wavelengths, hdrl_spectrum1D_wave_scale_linear);
162
163 const cpl_image * flux2 =
165
166 const cpl_image * flux2_e =
168
169 cpl_boolean are_equal2 =
170 are_hdrl_eq(hdrl_spectrum1D_get_flux(spec), flux2, flux2_e);
171 cpl_test_eq(are_equal2, CPL_TRUE);
172
174 cpl_test_null(spec2);
175 }
176
178 cpl_test_null(spec);
179
181 spectrum1d, wavelengths, hdrl_spectrum1D_wave_scale_linear);
182 cpl_test_nonnull(spec);
183
184 cpl_size not_rej = 0;
185 for(cpl_size i = 0; i < sz; ++i){
186 int rej;
187 const hdrl_value v = hdrl_spectrum1D_get_flux_value(spec, i, &rej);
188 if(rej) continue;
189 cpl_test_eq(v.error, 0.0);
190 not_rej ++;
191 }
192 cpl_test(not_rej > 0);
193
194 cpl_array_delete(wavelengths);
195 cpl_image_delete(spectrum1d);
196 cpl_image_delete(noise);
197 cpl_image_delete(spectrum1d_error);
198
200 cpl_test_null(spec);
201
202 /* deleting a nullptr should be NOOP */
204}
205
206/*----------------------------------------------------------------------------*/
211/*----------------------------------------------------------------------------*/
212
213hdrl_value test_val(hdrl_data_t lambda){
214 return (hdrl_value){lambda * 2.0, lambda *3.0};
215}
216
217void test_spectrum1D_constructor_analytical(void){
218
219 cpl_size sz = 10;
220 cpl_array * wav = cpl_array_new(sz, CPL_TYPE_DOUBLE);
221
222 for(cpl_size i = 0; i < sz; i++){
223 cpl_array_set(wav, i, (1 + i) * 10.0);
224 }
225
226 hdrl_spectrum1D * spec = hdrl_spectrum1D_create_analytic
227 (test_val, wav, hdrl_spectrum1D_wave_scale_linear);
228
229 for(cpl_size i = 0; i < sz; i++){
230
231 hdrl_value v = hdrl_spectrum1D_get_flux_value(spec, i, NULL);
232
233 cpl_test_abs(v.data, (1 + i) * 20.0, 1e-3);
234 cpl_test_abs(v.error, (1 + i) * 30.0, 1e-3);
235
236 hdrl_data_t ws = hdrl_spectrum1D_get_wavelength_value(spec, i, NULL);
237 hdrl_data_t ww = cpl_array_get(wav, i, NULL);
238 cpl_test_abs(ws, ww, 1e-3);
239 }
240
241 cpl_test_eq(hdrl_spectrum1D_wave_scale_linear,
243
244 cpl_array_delete(wav);
246}
247
248/*----------------------------------------------------------------------------*/
253/*----------------------------------------------------------------------------*/
254
255void test_spectrum1D_constructor_error(void){
256
257 cpl_image * spectrum1d_40 =
258 get_random_1d_img(40, 0.0f, 128.0f, CPL_TYPE_DOUBLE);
259 cpl_image * spectrum1d_error_40 =
260 get_random_1d_img(40, 0.0f, 1.0f,CPL_TYPE_DOUBLE);
261 set_1d_bpm(spectrum1d_40);
262 cpl_array * wavelengths_40 = get_wavelength(40,CPL_TYPE_DOUBLE);
263
264 cpl_image * spectrum1d_42 =
265 get_random_1d_img(42, 0.0f, 128.0f, CPL_TYPE_DOUBLE);
266 cpl_image * spectrum1d_error_42 =
267 get_random_1d_img(42, 0.0f, 1.0f,CPL_TYPE_DOUBLE);
268 set_1d_bpm(spectrum1d_42);
269 cpl_array * wavelengths_42 = get_wavelength(42,CPL_TYPE_DOUBLE);
270
271 hdrl_spectrum1D * spec = hdrl_spectrum1D_create(NULL, NULL,
272 wavelengths_40,
273 hdrl_spectrum1D_wave_scale_linear);
274
275 cpl_test_null(spec);
276 cpl_test_eq(get_error_code_and_reset(), CPL_ERROR_NULL_INPUT);
277
279
280 spec = hdrl_spectrum1D_create(spectrum1d_40, NULL,
281 NULL, hdrl_spectrum1D_wave_scale_linear);
282
283 cpl_test_null(spec);
284 cpl_test_eq(get_error_code_and_reset(), CPL_ERROR_NULL_INPUT);
285
287
288 spec = hdrl_spectrum1D_create(spectrum1d_40, spectrum1d_error_42,
289 wavelengths_40, hdrl_spectrum1D_wave_scale_linear);
290
291 cpl_test_null(spec);
292 cpl_test_eq(get_error_code_and_reset(), CPL_ERROR_INCOMPATIBLE_INPUT);
293
294 spec = hdrl_spectrum1D_create(spectrum1d_42, spectrum1d_error_40,
295 wavelengths_40, hdrl_spectrum1D_wave_scale_linear);
296
297 cpl_test_null(spec);
298 cpl_test_eq(get_error_code_and_reset(), CPL_ERROR_INCOMPATIBLE_INPUT);
299
300 spec = hdrl_spectrum1D_create(spectrum1d_40, spectrum1d_error_40,
301 wavelengths_42, hdrl_spectrum1D_wave_scale_linear);
302
303 cpl_test_null(spec);
304 cpl_test_eq(get_error_code_and_reset(), CPL_ERROR_INCOMPATIBLE_INPUT);
305
306 double el0 = cpl_array_get(wavelengths_40, 0, NULL);
307 cpl_array_set(wavelengths_40, 1, el0);
308
309 spec = hdrl_spectrum1D_create(spectrum1d_40, spectrum1d_error_40,
310 wavelengths_40, hdrl_spectrum1D_wave_scale_linear);
311 cpl_test_nonnull(spec);
312 cpl_test_eq(get_error_code_and_reset(), CPL_ERROR_NONE);
313
315
316 cpl_array_delete(wavelengths_40);
317 cpl_image_delete(spectrum1d_40);
318 cpl_image_delete(spectrum1d_error_40);
319
320 cpl_array_delete(wavelengths_42);
321 cpl_image_delete(spectrum1d_42);
322 cpl_image_delete(spectrum1d_error_42);
323}
324
325
326/*----------------------------------------------------------------------------*/
330/*----------------------------------------------------------------------------*/
331
332void test_spectrum1D_duplication(void){
333
334 cpl_image * spectrum1d =
335 get_random_1d_img(140, 0.0f, 128.0f, CPL_TYPE_DOUBLE);
336 cpl_image * spectrum1d_error =
337 get_random_1d_img(140, 0.0f, 1.0f,CPL_TYPE_DOUBLE);
338 set_1d_bpm(spectrum1d);
339 cpl_array * wavelengths = get_wavelength(140, CPL_TYPE_DOUBLE);
340
341 hdrl_spectrum1D * spec = hdrl_spectrum1D_create(
342 spectrum1d, NULL,
343 wavelengths, hdrl_spectrum1D_wave_scale_linear);
344
345 hdrl_spectrum1D * spec_copy = hdrl_spectrum1D_duplicate(spec);
346
348 hdrl_spectrum1D_delete(&spec_copy);
349
350 /*Issues with NULL?*/
351 hdrl_spectrum1D * should_be_null = hdrl_spectrum1D_duplicate(NULL);
352 cpl_test_null(should_be_null);
353
354 cpl_array_delete(wavelengths);
355 cpl_image_delete(spectrum1d);
356 cpl_image_delete(spectrum1d_error);
357}
358
359/*----------------------------------------------------------------------------*/
363/*----------------------------------------------------------------------------*/
364
365void test_spectrum1D_calculation_scalar(void){
366
367 cpl_image * spectrum1d =
368 get_random_1d_img(40, 0.0f, 128.0f, CPL_TYPE_DOUBLE);
369 cpl_image * spectrum1d_error =
370 get_random_1d_img(40, 0.0f, 1.0f, CPL_TYPE_DOUBLE);
371 set_1d_bpm(spectrum1d);
372 cpl_array * wavelengths = get_wavelength(40, CPL_TYPE_DOUBLE);
373
374 cpl_image_set(spectrum1d, 3, 1, 5.0);
375 cpl_image_set(spectrum1d_error, 3, 1, 2.1);
376
377 hdrl_spectrum1D * spec = hdrl_spectrum1D_create(
378 spectrum1d, spectrum1d_error,
379 wavelengths, hdrl_spectrum1D_wave_scale_linear);
380
381 hdrl_value vs = {1.5, 0.3};
382
383 hdrl_spectrum1D * spec_copy = hdrl_spectrum1D_mul_scalar_create(spec, vs);
385
386 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).data,
387 7.5, 1e-3);
388 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).error,
389 3.4889, 1e-3);
390 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).data,
391 7.5, 1e-3);
392 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).error,
393 3.4889, 1e-3);
394 hdrl_spectrum1D_delete(&spec_copy);
395
396 spec_copy = hdrl_spectrum1D_div_scalar_create(spec, vs);
398
399 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).data,
400 5.0, 1e-3);
401 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).error,
402 2.53179, 1e-3);
403 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy , 2, NULL).data,
404 5.0, 1e-3);
405 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy , 2, NULL).error,
406 2.53179, 1e-3);
407 hdrl_spectrum1D_delete(&spec_copy);
408
409 spec_copy = hdrl_spectrum1D_add_scalar_create(spec, vs);
411
412 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).data,
413 6.5, 1e-3);
414 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).error,
415 2.54951, 1e-3);
416 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).data,
417 6.5, 1e-3);
418 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).error,
419 2.54951, 1e-3);
420 hdrl_spectrum1D_delete(&spec_copy);
421
422 spec_copy = hdrl_spectrum1D_sub_scalar_create(spec, vs);
424
425 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).data,
426 5.0, 1e-3);
427 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).error,
428 2.5671, 1e-3);
429 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).data,
430 5.0, 1e-3);
431 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).error,
432 2.5671, 1e-3);
433 hdrl_spectrum1D_delete(&spec_copy);
434
435 vs.data = 2.0;
436 spec_copy = hdrl_spectrum1D_pow_scalar_create(spec, vs);
438
439 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).data,
440 pow(5.0, 2.0), 1e-3);
441 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).error,
442 28.3673, 1e-3);
443 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).data,
444 pow(5.0, 2.0), 1e-3);
445 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).error,
446 28.3673, 1e-3);
447 hdrl_spectrum1D_delete(&spec_copy);
448
450 spec = hdrl_spectrum1D_create(spectrum1d, spectrum1d_error,
451 wavelengths, hdrl_spectrum1D_wave_scale_linear);
452
453 vs.data = 2.0;
454 vs.error = 0.2;
455 spec_copy = hdrl_spectrum1D_exp_scalar_create(spec, vs);
457 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).data,
458 pow(2.0, 5.0), 1e-3);
459 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec, 2, NULL).error,
460 49.25087754, 1e-3);
461 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).data,
462 pow(2.0, 5.0), 1e-3);
463 cpl_test_abs(hdrl_spectrum1D_get_flux_value(spec_copy, 2, NULL).error,
464 49.25087754, 1e-3);
465 hdrl_spectrum1D_delete(&spec_copy);
466
467 cpl_array_delete(wavelengths);
468 cpl_image_delete(spectrum1d);
469 cpl_image_delete(spectrum1d_error);
471
472
473 /* Issues with NULL */
474 hdrl_value v = {100.0, 5.0};
475 cpl_test_null(hdrl_spectrum1D_sub_scalar_create(NULL, v));
476 cpl_test_null(hdrl_spectrum1D_add_scalar_create(NULL, v));
477 cpl_test_null(hdrl_spectrum1D_div_scalar_create(NULL, v));
478 cpl_test_null(hdrl_spectrum1D_mul_scalar_create(NULL, v));
479 cpl_test_null(hdrl_spectrum1D_pow_scalar_create(NULL, v));
480 cpl_test_null(hdrl_spectrum1D_exp_scalar_create(NULL, v));
481
482 cpl_test_eq(hdrl_spectrum1D_sub_scalar(NULL, v), CPL_ERROR_NONE);
483 cpl_test_eq(hdrl_spectrum1D_add_scalar(NULL, v), CPL_ERROR_NONE);
484 cpl_test_eq(hdrl_spectrum1D_div_scalar(NULL, v), CPL_ERROR_NONE);
485 cpl_test_eq(hdrl_spectrum1D_mul_scalar(NULL, v), CPL_ERROR_NONE);
486 cpl_test_eq(hdrl_spectrum1D_pow_scalar(NULL, v), CPL_ERROR_NONE);
487 cpl_test_eq(hdrl_spectrum1D_exp_scalar(NULL, v), CPL_ERROR_NONE);
488}
489
490/*----------------------------------------------------------------------------*/
494/*----------------------------------------------------------------------------*/
495
496void test_spectrum1D_calculation(void){
497
498 cpl_image * spectrum1d1 =
499 get_random_1d_img(40, 1.0f, 128.0f, CPL_TYPE_DOUBLE);
500 cpl_image * spectrum1d2 =
501 get_random_1d_img(40, 1.0f, 128.0f, CPL_TYPE_DOUBLE);
502
503 cpl_image * spectrum1d_error1 =
504 get_random_1d_img(40, 0.5f, 2.0f, CPL_TYPE_DOUBLE);
505 cpl_image * spectrum1d_error2 =
506 get_random_1d_img(40, 0.5f, 2.0f, CPL_TYPE_DOUBLE);
507 set_1d_bpm(spectrum1d1);
508 cpl_array * wavelengths = get_wavelength(40, CPL_TYPE_DOUBLE);
509
510 cpl_image_set(spectrum1d1, 15, 1, 8.0);
511 cpl_image_set(spectrum1d2, 15, 1, 4.0);
512 cpl_image_set(spectrum1d_error1, 15, 1, 2.0);
513 cpl_image_set(spectrum1d_error2, 15, 1, 1.0);
514
515 hdrl_spectrum1D * s1 = hdrl_spectrum1D_create(
516 spectrum1d1, spectrum1d_error1,
517 wavelengths, hdrl_spectrum1D_wave_scale_linear);
518
519 hdrl_spectrum1D * s2 = hdrl_spectrum1D_create(
520 spectrum1d2, spectrum1d_error2,
521 wavelengths, hdrl_spectrum1D_wave_scale_linear);
522
523 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s1, 14, NULL).data, 8.0, 1e-3);
524 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s2, 14, NULL).data, 4.0, 1e-3);
525 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s1, 14, NULL).error, 2.0, 1e-3);
526 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s2, 14, NULL).error, 1.0, 1e-3);
527
528 hdrl_spectrum1D * s3 = hdrl_spectrum1D_div_spectrum_create(s1, s2);
529 hdrl_spectrum1D * s4 = hdrl_spectrum1D_duplicate(s1);
531 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s3, 14, NULL).data, 2.0, 1e-3);
532 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s3, 14, NULL).error,
533 0.707107, 1e-3);
534 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s4, 14, NULL).data, 2.0,1e-3);
535 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s4, 14, NULL).error,
536 0.707107, 1e-3);
539
543 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s3, 14, NULL).data,
544 32.0, 1e-3);
545 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s3, 14, NULL).error,
546 11.3137, 1e-3);
547 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s4, 14, NULL).data,
548 32.0, 1e-3);
549 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s4, 14, NULL).error,
550 11.3137, 1e-3);
553
554
558 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s3, 14, NULL).data, 4.0, 1e-3);
559 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s3, 14, NULL).error, 2.23607,
560 1e-3);
561 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s4, 14, NULL).data, 4.0, 1e-3);
562 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s4, 14, NULL).error, 2.23607,
563 1e-3);
566
570 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s3, 14, NULL).data, 12.0, 1e-3);
571 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s3, 14, NULL).error, 2.23607,
572 1e-3);
573 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s4, 14, NULL).data, 12.0, 1e-3);
574 cpl_test_abs(hdrl_spectrum1D_get_flux_value(s4, 14, NULL).error, 2.23607,
575 1e-3);
578
581
582 cpl_array_delete(wavelengths);
583 cpl_image_delete(spectrum1d1);
584 cpl_image_delete(spectrum1d2);
585 cpl_image_delete(spectrum1d_error1);
586 cpl_image_delete(spectrum1d_error2);
587}
588
589/*----------------------------------------------------------------------------*/
594/*----------------------------------------------------------------------------*/
595
596void test_spectrum1D_calculation_error(void){
597
598 test_calc_creat_error(hdrl_spectrum1D_div_spectrum_create);
599 test_calc_creat_error(hdrl_spectrum1D_add_spectrum_create);
600 test_calc_creat_error(hdrl_spectrum1D_sub_spectrum_create);
601 test_calc_creat_error(hdrl_spectrum1D_mul_spectrum_create);
602
603 test_calc_error(hdrl_spectrum1D_div_spectrum);
604 test_calc_error(hdrl_spectrum1D_add_spectrum);
605 test_calc_error(hdrl_spectrum1D_sub_spectrum);
606 test_calc_error(hdrl_spectrum1D_mul_spectrum);
607}
608
609/*----------------------------------------------------------------------------*/
614/*----------------------------------------------------------------------------*/
615
616void test_spectrum1D_conversion_wavelength_scale(void){
617
618 cpl_image * spectrum1d =
619 get_random_1d_img(40, 1.0f, 128.0f, CPL_TYPE_DOUBLE);
620
621 cpl_image * spectrum1d_error1 =
622 get_random_1d_img(40, 0.5f, 2.0f, CPL_TYPE_DOUBLE);
623
624 set_1d_bpm(spectrum1d);
625 cpl_array * wavelengths = get_wavelength(40, CPL_TYPE_DOUBLE);
626
627 hdrl_spectrum1D * sp = hdrl_spectrum1D_create(
628 spectrum1d, spectrum1d_error1,
629 wavelengths, hdrl_spectrum1D_wave_scale_linear);
630
631 hdrl_data_t w1 = cpl_array_get(wavelengths, 4, NULL);
632 hdrl_data_t w2 = cpl_array_get(wavelengths, 6, NULL);
633
634 hdrl_spectrum1D * s_lg =
636 hdrl_spectrum1D_wavelength w_lg =
638
639 cpl_test_eq(w_lg.scale, hdrl_spectrum1D_wave_scale_log);
640
641 cpl_test_abs(log(w1), cpl_array_get(w_lg.wavelength, 4, NULL), 1e-3);
642 cpl_test_abs(log(w2), cpl_array_get(w_lg.wavelength, 6, NULL), 1e-3);
643
644 hdrl_spectrum1D * s_lin =
646 hdrl_spectrum1D_wavelength w_lin = hdrl_spectrum1D_get_wavelength(s_lin);
647 cpl_test_eq(w_lin.scale, hdrl_spectrum1D_wave_scale_linear);
648
649 cpl_test_abs(w1, cpl_array_get(w_lin.wavelength, 4, NULL), 1e-3);
650 cpl_test_abs(w2, cpl_array_get(w_lin.wavelength, 6, NULL), 1e-3);
651
652 /*conversions on functions already in the right scale should behave like a
653 * duplicate*/
654
655 hdrl_spectrum1D * s_lg2 =
657 hdrl_spectrum1D_wavelength w_lg2 = hdrl_spectrum1D_get_wavelength(s_lg);
658 cpl_test_eq(w_lg2.scale, hdrl_spectrum1D_wave_scale_log);
659 cpl_test_noneq_ptr(s_lg2, s_lg);
660
661 hdrl_spectrum1D * s_lin2 =
663 hdrl_spectrum1D_wavelength w_lin2 = hdrl_spectrum1D_get_wavelength(s_lin);
664 cpl_test_eq(w_lin2.scale, hdrl_spectrum1D_wave_scale_linear);
665 cpl_test_noneq_ptr(s_lin2, s_lin);
666
670 hdrl_spectrum1D_delete(&s_lin2);
671
672 /* Mutator functions */
673
674 /* linear -> log */
675 cpl_error_code e = hdrl_spectrum1D_wavelength_convert_to_log(sp);
676 cpl_test_eq(e, CPL_ERROR_NONE);
677
678 hdrl_spectrum1D_wavelength sp_w_lg = hdrl_spectrum1D_get_wavelength(sp);
679
680 cpl_test_eq(sp_w_lg.scale, hdrl_spectrum1D_wave_scale_log);
681 cpl_test_abs(log(w1), cpl_array_get(sp_w_lg.wavelength, 4, NULL), 1e-3);
682 cpl_test_abs(log(w2), cpl_array_get(sp_w_lg.wavelength, 6, NULL), 1e-3);
683
684 /* log -> log NOOP*/
686 cpl_test_eq(e, CPL_ERROR_NONE);
687
688 sp_w_lg = hdrl_spectrum1D_get_wavelength(sp);
689
690 cpl_test_eq(sp_w_lg.scale, hdrl_spectrum1D_wave_scale_log);
691 cpl_test_abs(log(w1), cpl_array_get(sp_w_lg.wavelength, 4, NULL), 1e-3);
692 cpl_test_abs(log(w2), cpl_array_get(sp_w_lg.wavelength, 6, NULL), 1e-3);
693
694 /* log -> linear */
696 cpl_test_eq(e, CPL_ERROR_NONE);
697
698 hdrl_spectrum1D_wavelength sp_w_lin = hdrl_spectrum1D_get_wavelength(sp);
699 cpl_test_eq(sp_w_lin.scale, hdrl_spectrum1D_wave_scale_linear);
700 cpl_test_abs(w1, cpl_array_get(sp_w_lin.wavelength, 4, NULL), 1e-3);
701 cpl_test_abs(w2, cpl_array_get(sp_w_lin.wavelength, 6, NULL), 1e-3);
702
703 /* linear -> linear NOOP*/
705 cpl_test_eq(e, CPL_ERROR_NONE);
706
707 sp_w_lin = hdrl_spectrum1D_get_wavelength(sp);
708 cpl_test_eq(sp_w_lin.scale, hdrl_spectrum1D_wave_scale_linear);
709 cpl_test_abs(w1, cpl_array_get(sp_w_lin.wavelength, 4, NULL), 1e-3);
710 cpl_test_abs(w2, cpl_array_get(sp_w_lin.wavelength, 6, NULL), 1e-3);
711
713 cpl_array_delete(wavelengths);
714 cpl_image_delete(spectrum1d_error1);
715 cpl_image_delete(spectrum1d);
716
717 /* Test for null */
718
720 CPL_ERROR_NONE);
722 CPL_ERROR_NONE);
723
726}
727
728/*----------------------------------------------------------------------------*/
733/*----------------------------------------------------------------------------*/
734
735void test_spectrum1D_mul_wavelength(void){
736
737 cpl_image * spectrum1d =
738 get_random_1d_img(40, 1.0f, 128.0f, CPL_TYPE_DOUBLE);
739
740 cpl_image * spectrum1d_error1 =
741 get_random_1d_img(40, 0.5f, 2.0f, CPL_TYPE_DOUBLE);
742
743 set_1d_bpm(spectrum1d);
744 cpl_array * wavelengths = get_wavelength(40, CPL_TYPE_DOUBLE);
745
746 hdrl_spectrum1D * sp = hdrl_spectrum1D_create(
747 spectrum1d, spectrum1d_error1,
748 wavelengths, hdrl_spectrum1D_wave_scale_linear);
749
750 hdrl_data_t w1 = cpl_array_get(wavelengths, 4, NULL);
751 hdrl_data_t w2 = cpl_array_get(wavelengths, 6, NULL);
752
753 /* mutator functions */
754
755 cpl_error_code e = hdrl_spectrum1D_wavelength_mult_scalar_linear(sp, 1e3);
756 cpl_test_eq(e, CPL_ERROR_NONE);
757
758 hdrl_spectrum1D_wavelength sp_mul_1e3 = hdrl_spectrum1D_get_wavelength(sp);
759
760 cpl_test_abs(w1 * 1e3, cpl_array_get(sp_mul_1e3.wavelength, 4, NULL), 1e-3);
761 cpl_test_abs(w2 * 1e3, cpl_array_get(sp_mul_1e3.wavelength, 6, NULL), 1e-3);
762
764 cpl_test_eq(e, CPL_ERROR_NONE);
765
767
769 cpl_test_eq(e, CPL_ERROR_NONE);
770
771 hdrl_spectrum1D_wavelength sp_mul_log = hdrl_spectrum1D_get_wavelength(sp);
772
773 cpl_test_abs(log(w1 * 1e3), cpl_array_get(sp_mul_log.wavelength, 4, NULL),
774 1e-3);
775 cpl_test_abs(log(w2 * 1e3), cpl_array_get(sp_mul_log.wavelength, 6, NULL),
776 1e-3);
777
778 /* non mutator functions */
780
782 spectrum1d, spectrum1d_error1,
783 wavelengths, hdrl_spectrum1D_wave_scale_linear);
784
785 hdrl_spectrum1D * sp_lin =
788 hdrl_spectrum1D * sp_log =
791
792 sp_mul_1e3 = hdrl_spectrum1D_get_wavelength(sp_lin);
793 sp_mul_log = hdrl_spectrum1D_get_wavelength(sp_log);
794
795 cpl_test_abs(w1 * 1e-4, cpl_array_get(sp_mul_1e3.wavelength, 4, NULL), 1e-6);
796 cpl_test_abs(w2 * 1e-4, cpl_array_get(sp_mul_1e3.wavelength, 6, NULL), 1e-6);
797
798 cpl_test_abs(
799 log(w1 * 1e-4), cpl_array_get(sp_mul_log.wavelength, 4, NULL), 1e-6);
800 cpl_test_abs(
801 log(w2 * 1e-4), cpl_array_get(sp_mul_log.wavelength, 6, NULL), 1e-6);
802
803 cpl_test_noneq(hdrl_spectrum1D_wavelength_mult_scalar_linear(sp_log, -2.0),
804 CPL_ERROR_NONE);
805 cpl_test_noneq(get_error_code_and_reset(), CPL_ERROR_NONE);
806
807 cpl_test_null(
809 cpl_test_noneq(get_error_code_and_reset(), CPL_ERROR_NONE);
810
811 hdrl_spectrum1D_delete(&sp_log);
812 hdrl_spectrum1D_delete(&sp_lin);
813 cpl_array_delete(wavelengths);
814 cpl_image_delete(spectrum1d_error1);
815 cpl_image_delete(spectrum1d);
816
818 CPL_ERROR_NONE);
820 3.0));
821}
822
823/*----------------------------------------------------------------------------*/
827/*----------------------------------------------------------------------------*/
828
829void test_spectrum1D_shift_wavelength(void){
830
831 cpl_image * spectrum1d =
832 get_random_1d_img(40, 1.0f, 128.0f, CPL_TYPE_DOUBLE);
833
834 cpl_image * spectrum1d_error1 =
835 get_random_1d_img(40, 0.5f, 2.0f, CPL_TYPE_DOUBLE);
836
837 set_1d_bpm(spectrum1d);
838 cpl_array * wavelengths = get_wavelength(40, CPL_TYPE_DOUBLE);
839
840 hdrl_data_t w1 = cpl_array_get(wavelengths, 4, NULL);
841 hdrl_data_t w2 = cpl_array_get(wavelengths, 6, NULL);
842
843 hdrl_spectrum1D * sp1 = hdrl_spectrum1D_create(
844 spectrum1d, spectrum1d_error1,
845 wavelengths, hdrl_spectrum1D_wave_scale_linear);
846
847 hdrl_spectrum1D * sp2 = hdrl_spectrum1D_create(
848 spectrum1d, spectrum1d_error1,
849 wavelengths, hdrl_spectrum1D_wave_scale_log);
850
851 /* shift works in the same way both for linear and log scale */
852 cpl_error_code err = hdrl_spectrum1D_wavelength_shift(sp1, 3.0);
853 cpl_test_eq(err, CPL_ERROR_NONE);
854
855 err = hdrl_spectrum1D_wavelength_shift(sp2, -3.0);
856 cpl_test_eq(err, CPL_ERROR_NONE);
857
858 hdrl_spectrum1D_wavelength sp_p_3 = hdrl_spectrum1D_get_wavelength(sp1);
859 hdrl_spectrum1D_wavelength sp_m_3 = hdrl_spectrum1D_get_wavelength(sp2);
860
861 cpl_test_abs(w1 + 3.0, cpl_array_get(sp_p_3.wavelength, 4, NULL), 1e-3);
862 cpl_test_abs(w2 + 3.0, cpl_array_get(sp_p_3.wavelength, 6, NULL), 1e-3);
863
864 cpl_test_abs(w1 - 3.0, cpl_array_get(sp_m_3.wavelength, 4, NULL), 1e-3);
865 cpl_test_abs(w2 - 3.0, cpl_array_get(sp_m_3.wavelength, 6, NULL), 1e-3);
866
867 hdrl_spectrum1D* sp_new_1 =
869 hdrl_spectrum1D* sp_new_2 =
871
874
875 hdrl_spectrum1D_wavelength sp_n_1_w =
877 hdrl_spectrum1D_wavelength sp_n_2_w =
879
880 cpl_test_abs(w1, cpl_array_get(sp_n_1_w.wavelength, 4, NULL), 1e-3);
881 cpl_test_abs(w2, cpl_array_get(sp_n_1_w.wavelength, 6, NULL), 1e-3);
882 cpl_test_abs(w1, cpl_array_get(sp_n_2_w.wavelength, 4, NULL), 1e-3);
883 cpl_test_abs(w2, cpl_array_get(sp_n_2_w.wavelength, 6, NULL), 1e-3);
884
885 hdrl_spectrum1D_delete(&sp_new_1);
886 hdrl_spectrum1D_delete(&sp_new_2);
887
888 /* Issues with NULL? */
889 err = hdrl_spectrum1D_wavelength_shift(NULL, 3.0);
890 cpl_test_eq(err, CPL_ERROR_NONE);
891 hdrl_spectrum1D * should_be_null =
893 cpl_test_null(should_be_null);
894
895 cpl_array_delete(wavelengths);
896 cpl_image_delete(spectrum1d_error1);
897 cpl_image_delete(spectrum1d);
898}
899
900/*----------------------------------------------------------------------------*/
904/*----------------------------------------------------------------------------*/
905void test_spectrum1D_table_conversion(void){
906
907 const cpl_size sz_ori = 17;
908 cpl_array * unshuffled_lambda = NULL;
909 hdrl_spectrum1D * sp1 = get_spectrum1D_sin_shuffled(sz_ori, 2, CPL_TRUE,
910 &unshuffled_lambda);
911 cpl_array_delete(unshuffled_lambda);
912
913 cpl_table * tab = hdrl_spectrum1D_convert_to_table(sp1, "flux",
914 "lambdas", "flux_e", "flux_bpm");
915
916 /*ad rejected pixels*/
917 cpl_table_set_int(tab, "flux_bpm", 0, 1);
918 cpl_table_set_int(tab, "flux_bpm", sz_ori - 1, 1);
919
920 hdrl_spectrum1D * sp2 = hdrl_spectrum1D_convert_from_table
921 (tab, "flux", "lambdas", "flux_e", "flux_bpm",
922 hdrl_spectrum1D_wave_scale_linear);
923
924 const cpl_image * flux1 =
926 const cpl_image * flux1_e =
928
929 const cpl_array * lambdas1 = hdrl_spectrum1D_get_wavelength(sp1).wavelength;
930
931 const cpl_image * flux2 =
933 const cpl_image * flux2_e =
935 const cpl_array * lambdas2 = hdrl_spectrum1D_get_wavelength(sp2).wavelength;
936
937 const cpl_size sz_x = cpl_image_get_size_x(flux1);
938 const cpl_size sz_y = cpl_image_get_size_y(flux1);
939 const cpl_size sz = cpl_array_get_size(lambdas1);
940
941 cpl_test_eq(sz_x, cpl_image_get_size_x(flux2));
942 cpl_test_eq(sz_y, cpl_image_get_size_y(flux2));
943 cpl_test_eq(sz, cpl_array_get_size(lambdas2));
944
945 for(cpl_size i = 1; i < sz - 1; i++){
946
947 int rej1 = 0;
948 int rej2 = 0;
949 const double flx1 = cpl_image_get(flux1, i + 1, 1, &rej1);
950 const double flx2 = cpl_image_get(flux2, i + 1, 1, &rej2);
951
952
953 cpl_test_abs(flx1, flx2, 1e-3);
954 cpl_test_eq(rej1, rej2);
955 rej1 = rej2 = 0;
956
957 const double flx1_e = cpl_image_get(flux1_e, i + 1, 1, &rej1);
958 const double flx2_e = cpl_image_get(flux2_e, i + 1, 1, &rej2);
959
960 cpl_test_abs(flx1_e, flx2_e, 1e-3);
961 cpl_test_eq(rej1, rej2);
962 rej1 = rej2 = 0;
963
964
965 const double l1 = cpl_array_get(lambdas1, i, &rej1);
966 const double l2 = cpl_array_get(lambdas2, i, &rej2);
967
968 cpl_test_abs(l1, l2, 1e-3);
969 cpl_test_eq(rej1, rej2);
970 }
971
972 cpl_test(cpl_image_is_rejected(flux2, 1, 1));
973 cpl_test(cpl_image_is_rejected(flux2, sz_ori, 1));
974
975
977 "flux2", "lambdas2", NULL, "flux2_bpm");
978
979 for(cpl_size i = 0; i < sz; i++){
980
981 int rej = 0;
982
983 const int bpm1 = cpl_table_get_int(tab, "flux_bpm", i, &rej);
984 const int bpm2 = cpl_table_get_int(tab, "flux2_bpm", i, &rej);
985 cpl_test_eq(bpm1, bpm2);
986
987 if(bpm1) continue;
988
989 const double flx1 = cpl_table_get(tab, "flux", i, &rej);
990 const double flx2 = cpl_table_get(tab, "flux2", i, &rej);
991 cpl_test_abs(flx1, flx2, 1e-3);
992
993 const double l1 = cpl_table_get(tab, "lambdas", i, &rej);
994 const double l2 = cpl_table_get(tab, "lambdas2", i, &rej);
995 cpl_test_abs(l1, l2, 1e-3);
996 }
997
1000 (tab, "flux", "lambdas", NULL, NULL,
1001 hdrl_spectrum1D_wave_scale_linear);
1002
1003 /* Check save */
1004 const char* filename = "check_spectrum1D.fits";
1005 hdrl_spectrum1D_save(sp2, filename);
1006 remove(filename);
1007
1008 for(cpl_size i = 0; i < hdrl_spectrum1D_get_size(sp2); ++i){
1009 int rej;
1010
1011 const double fx1 = hdrl_spectrum1D_get_flux_value(sp1, i, &rej).data;
1012
1013 const double fx2 = hdrl_spectrum1D_get_flux_value(sp2, i, &rej).data;
1014 const double fx2_e = hdrl_spectrum1D_get_flux_value(sp2, i, &rej).error;
1015
1016 cpl_test_eq(fx2_e, 0.0);
1017 cpl_test_abs(fx1, fx2, 1e-3);
1018 }
1019
1021 cpl_table_delete(tab);
1022
1023 cpl_table * f_only =
1024 hdrl_spectrum1D_convert_to_table(sp1, "flux", NULL, NULL, NULL);
1025 cpl_table * l_only =
1026 hdrl_spectrum1D_convert_to_table(sp1, NULL, "wav", NULL, NULL);
1027 cpl_table * f_and_e =
1028 hdrl_spectrum1D_convert_to_table(sp1, "flux", NULL, "error", NULL);
1029
1030 for(cpl_size i = 0; i < hdrl_spectrum1D_get_size(sp1); ++i){
1031 const double f = hdrl_spectrum1D_get_flux_value(sp1, i, NULL).data;
1032 const double e = hdrl_spectrum1D_get_flux_value(sp1, i, NULL).error;
1033 const double w = hdrl_spectrum1D_get_wavelength_value(sp1, i, NULL);
1034
1035 const double f1 = cpl_table_get(f_only, "flux", i, NULL);
1036 const double f2 = cpl_table_get(f_and_e, "flux", i, NULL);
1037 cpl_test_eq(f, f1);
1038 cpl_test_eq(f, f2);
1039
1040 const double w1 = cpl_table_get(l_only, "wav", i, NULL);
1041 cpl_test_eq(w, w1);
1042
1043 const double e1 = cpl_table_get(f_and_e, "error", i, NULL);
1044 cpl_test_eq(e, e1);
1045 }
1046
1047
1049 cpl_table_delete(f_only);
1050 cpl_table_delete(l_only);
1051 cpl_table_delete(f_and_e);
1052}
1053/*----------------------------------------------------------------------------*/
1057/*----------------------------------------------------------------------------*/
1058
1059void test_spectrum1D_resample_spectrum(cpl_boolean add_peak){
1060
1061 /* Initialize variables */
1062 cpl_array * unshuffled_lambda = NULL;
1063 hdrl_spectrum1D * sp = get_spectrum1D_sin_shuffled(17, 2, add_peak, &unshuffled_lambda);
1064
1065 hdrl_spectrum1D_wavelength wl = hdrl_spectrum1D_get_wavelength(sp);
1066 cpl_size sz = hdrl_spectrum1D_get_size(sp);
1067 cpl_array * new_lambda = cpl_array_new(sz, HDRL_TYPE_DATA);
1068 for(cpl_size i = 0; i < sz; i++){
1069 double d = cpl_array_get(unshuffled_lambda, i, NULL)
1070 + cpl_array_get(unshuffled_lambda, CPL_MIN(sz - 1, i + 1), NULL);
1071 cpl_array_set(new_lambda, i, d/2);
1072 }
1073
1074 wl.wavelength = new_lambda;
1075
1076
1077 /* Tests variables */
1078 int rej;
1079 hdrl_parameter *pars;
1080 hdrl_spectrum1D *resampled_lambda;
1081
1082
1083 /* Test 1 */
1085 resampled_lambda = hdrl_spectrum1D_resample(sp, &wl, pars);
1086
1087 const double data2_fit = add_peak ? 116.368 : 209.577;
1088 const double data3_fit = add_peak ? 303.376 : 199.524;
1089
1090 cpl_test_abs(hdrl_spectrum1D_get_flux_value(resampled_lambda, 2, &rej).data, data2_fit, 1e-3);
1091 cpl_test_abs(hdrl_spectrum1D_get_flux_value(resampled_lambda, 3, &rej).data, data3_fit, 1e-3);
1092
1094 hdrl_spectrum1D_delete(&resampled_lambda);
1095
1096
1097 /* Test 2 */
1098 pars = hdrl_spectrum1D_resample_interpolate_parameter_create(hdrl_spectrum1D_interp_akima);
1099 resampled_lambda = hdrl_spectrum1D_resample(sp, &wl, pars);
1100
1101 const double data2_interp = add_peak ? 208.699 : 209.65;
1102 const double data3_interp = add_peak ? 247.949 : 199.585;
1103
1104 cpl_test_abs(hdrl_spectrum1D_get_flux_value(resampled_lambda, 2, &rej).data, data2_interp, 1e-3);
1105 cpl_test_abs(hdrl_spectrum1D_get_flux_value(resampled_lambda, 3, &rej).data, data3_interp, 1e-3);
1106
1108 hdrl_spectrum1D_delete(&resampled_lambda);
1109
1110
1111 /* Test 3 */
1113 resampled_lambda = hdrl_spectrum1D_resample(sp, &wl, pars);
1114
1115 const double data2_integrate = 207.878;
1116 const double data3_integrate = add_peak ? 245.443 : 197.992;
1117
1118 cpl_test_abs(hdrl_spectrum1D_get_flux_value(resampled_lambda, 2, &rej).data, data2_integrate, 1e-3);
1119 cpl_test_abs(hdrl_spectrum1D_get_flux_value(resampled_lambda, 3, &rej).data, data3_integrate, 1e-3);
1120
1122 hdrl_spectrum1D_delete(&resampled_lambda);
1123
1124
1125 /* Clean up */
1127 cpl_array_delete(new_lambda);
1128 cpl_array_delete(unshuffled_lambda);
1129}
1130
1131/*----------------------------------------------------------------------------*/
1136/*----------------------------------------------------------------------------*/
1137
1138void test_spectrum1D_resample_spectrum_private_funcs(void){
1139
1140 double testValue1 = 2.1;
1141 double testValue2 = 3.5;
1142 double testValue3 = 5.5;
1143
1144 {
1145 double x[] = {3, 2.1, 5.5, 8.7, 3.3, 5.6, 2.1};
1146 double y1[] = {11, 88, -22, 56, 4, 22, 23};
1147 double y2[] = { 2, 55, 2, 27, 23, 1, 5};
1148
1149 double x_sorted[] = {2.1, 2.1, 3, 3.3, 5.5, 5.6, 8.7};
1150 double y1_sorted[] = {88, 23, 11, 4, -22, 22, 56};
1151 double y2_sorted[] = {55, 5, 2, 23, 2, 1, 27};
1152
1153 const cpl_size l = LENGTH_ARRAY(x);
1154
1155 hdrl_sort_on_x(x, y1, y2, l, CPL_FALSE);
1156
1157 cpl_test_eq(x[0], testValue1);
1158 cpl_test_eq(x[1], testValue1);
1159
1160 /*test the case were x is duplicated.*/
1161 cpl_test(y1[1] != y1[0]);
1162 cpl_test(y2[1] != y2[0]);
1163
1164 for(cpl_size i = 2; i < l; ++i){
1165 cpl_test_eq(x[i], x_sorted[i]);
1166 cpl_test_eq(y1[i], y1_sorted[i]);
1167 cpl_test_eq(y2[i], y2_sorted[i]);
1168 }
1169 }
1170
1171 {
1172 double x[] = {3, 2.1, 5.5, 8.7, 3.3, 5.6, 2.1};
1173 double y2[] = { 2, 55, 2, 27, 23, 1, 5};
1174
1175 double x_sorted[] = {2.1, 2.1, 3, 3.3, 5.5, 5.6, 8.7};
1176 double y2_sorted[] = {55, 5, 2, 23, 2, 1, 27};
1177
1178 const cpl_size l = LENGTH_ARRAY(x);
1179
1180 hdrl_sort_on_x(x, NULL, y2, l, CPL_FALSE);
1181
1182 cpl_test_eq(x[0], testValue1);
1183 cpl_test_eq(x[1], testValue1);
1184
1185 cpl_test(y2[1] != y2[0]);
1186
1187 for(cpl_size i = 2; i < l; ++i){
1188 cpl_test_eq(x[i], x_sorted[i]);
1189 cpl_test_eq(y2[i], y2_sorted[i]);
1190 }
1191 }
1192
1193 {
1194 double x[] = {3, 2.1, 5.5, 8.7, 3.3, 5.6, 2.1};
1195 double y1[] = {11, 88, -22, 56, 4, 22, 23};
1196
1197 double x_sorted[] = {2.1, 2.1, 3, 3.3, 5.5, 5.6, 8.7};
1198 double y1_sorted[] = {88, 23, 11, 4, -22, 22, 56};
1199
1200 const cpl_size l = LENGTH_ARRAY(x);
1201
1202 hdrl_sort_on_x(x, y1, NULL, l, CPL_FALSE);
1203
1204 cpl_test_eq(x[0], testValue1);
1205 cpl_test_eq(x[1], testValue1);
1206
1207 /*test the case were x is duplicated.*/
1208 cpl_test(y1[1] != y1[0]);
1209
1210 for(cpl_size i = 2; i < l; ++i){
1211 cpl_test_eq(x[i], x_sorted[i]);
1212 cpl_test_eq(y1[i], y1_sorted[i]);
1213 }
1214 }
1215
1216 /* test duplicates search and median */
1217 {
1218 /* edge case, all x are equal, even number of samples*/
1219 double x[] = {1,1,1,1,1};
1220 double y1[] = {5,4,3,2,5};
1221 double y2[] = {8,7,5,2,6};
1222
1223 cpl_size l = LENGTH_ARRAY(x);
1224
1225 l = hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1226 (x, y1, y2, l);
1227
1228 cpl_test_eq(l, 1);
1229 cpl_test_eq(x[0], 1);
1230 cpl_test_eq(y1[0], 4);
1231 cpl_test_eq(y2[0], 6);
1232 }
1233
1234 {
1235 /* edge case, all x are equal, odd number of samples*/
1236 double x[] = {1,1,1,1,1,1};
1237 double y1[] = {5,4,3,2,5,2.5};
1238 double y2[] = {8,7,5,2,6,4.6};
1239
1240 cpl_size l = LENGTH_ARRAY(x);
1241
1242 l = hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1243 (x, y1, y2, l);
1244
1245 cpl_test_eq(l, 1);
1246 cpl_test_eq(x[0], 1);
1247 cpl_test_eq(y1[0], testValue2);
1248 cpl_test_eq(y2[0], testValue3);
1249 }
1250
1251 /* test duplicates search and median */
1252 {
1253 /* edge case, all x are equal, except the first, even number of samples*/
1254 double x[] = {1,2,2,2,2,2};
1255 double y1[] = {55,5,4,3,2,5};
1256 double y2[] = {88,8,7,5,2,6};
1257
1258 cpl_size l = LENGTH_ARRAY(x);
1259
1260 l = hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1261 (x, y1, y2, l);
1262
1263 cpl_test_eq(l, 2);
1264 cpl_test_eq(x[1], 2);
1265 cpl_test_eq(y1[1], 4);
1266 cpl_test_eq(y2[1], 6);
1267 cpl_test_eq(x[0], 1);
1268 cpl_test_eq(y1[0], 55);
1269 cpl_test_eq(y2[0], 88);
1270 }
1271
1272 {
1273 /* edge case, all x are equal, except the last, odd number of samples*/
1274 double x[] = {1, 1, 1, 1, 1, 1, 8};
1275 double y1[] = {5, 4, 3, 2, 5, 2.5, 77};
1276 double y2[] = {8, 7, 5, 2, 6, 4.6, 96};
1277
1278 cpl_size l = LENGTH_ARRAY(x);
1279
1280 l = hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1281 (x, y1, y2, l);
1282
1283 cpl_test_eq(l, 2);
1284 cpl_test_eq(x[0] , 1);
1285 cpl_test_eq(y1[0], testValue2);
1286 cpl_test_eq(y2[0], testValue3);
1287
1288 cpl_test_eq(x[1] , 8);
1289 cpl_test_eq(y1[1], 77);
1290 cpl_test_eq(y2[1], 96);
1291 }
1292
1293 /*more realist case with hunks of equal elements*/
1294 {
1295 {
1296 double x[] = {1, 2, 2, 3, 3, 3, 5, 6, 7, 7, 8, 9, 10, 10, 10, 11};
1297 double y1[] = {4, 3, 7, 8, 9, 4, 3, 7, 2, 4, 5, 2, 8, 7, 1, 12};
1298 double y2[] = {3, 6, 7, 8, 4, 5, 8, 3, 5, 1, 3, 8, 44, 33, 55, 45};
1299
1300 double x_f[] = {1, 2, 3, 5, 6, 7, 8, 9, 10, 11};
1301 double y1_f[] = {4, 5, 8, 3, 7, 3, 5, 2, 7, 12};
1302 double y2_f[] = {3, 6.5, 5, 8, 3, 3, 3, 8, 44, 45};
1303
1304 cpl_size l = LENGTH_ARRAY(x);
1305
1306 l = hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1307 (x, y1, y2, l);
1308
1309 cpl_test_eq(l, 10);
1310
1311 for(cpl_size i = 0; i < l; ++i){
1312 cpl_test_eq(x[i], x_f[i]);
1313 cpl_test_eq(y1[i], y1_f[i]);
1314 cpl_test_eq(y2[i], y2_f[i]);
1315 }
1316 }
1317 }
1318
1319 /*more realist case with hunks of equal elements, one hunk at the end*/
1320 {
1321 {
1322 double x[] = {1, 2, 2, 3, 3, 3, 5, 6, 7, 7, 8, 9, 10, 10, 10, 11, 11};
1323 double y1[] = {4, 3, 7, 8, 9, 4, 3, 7, 2, 4, 5, 2, 8, 7, 1, 12, 2};
1324 double y2[] = {3, 6, 7, 8, 4, 5, 8, 3, 5, 1, 3, 8, 44, 33, 55, 45, 5};
1325
1326 double x_f[] = {1, 2, 3, 5, 6, 7, 8, 9, 10, 11};
1327 double y1_f[] = {4, 5, 8, 3, 7, 3, 5, 2, 7, 7};
1328 double y2_f[] = {3, 6.5, 5, 8, 3, 3, 3, 8, 44, 25};
1329
1330 cpl_size l = LENGTH_ARRAY(x);
1331
1332 l = hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1333 (x, y1, y2, l);
1334
1335 cpl_test_eq(l, 10);
1336
1337 for(cpl_size i = 0; i < l; ++i){
1338 cpl_test_eq(x[i], x_f[i]);
1339 cpl_test_eq(y1[i], y1_f[i]);
1340 cpl_test_eq(y2[i], y2_f[i]);
1341 }
1342 }
1343 }
1344
1345 /*more realist case with hunks of equal elements, one hunk at the beginning*/
1346 {
1347 {
1348 double x[] = {1, 1, 2, 2, 3, 3, 3, 5, 6, 7, 7, 8, 9, 10, 10, 10, 11};
1349 double y1[] = {5, 7, 3, 7, 8, 9, 4, 3, 7, 2, 4, 5, 2, 8, 7, 1, 12};
1350 double y2[] = {1, 3, 6, 7, 8, 4, 5, 8, 3, 5, 1, 3, 8, 44, 33, 55, 45};
1351
1352 double x_f[] = {1, 2, 3, 5, 6, 7, 8, 9, 10, 11};
1353 double y1_f[] = {6, 5, 8, 3, 7, 3, 5, 2, 7, 12};
1354 double y2_f[] = {2, 6.5, 5, 8, 3, 3, 3, 8, 44, 45};
1355
1356 cpl_size l = LENGTH_ARRAY(x);
1357
1358 l = hdrl_spectrum1D_resample_filter_dups_and_substitute_with_median
1359 (x, y1, y2, l);
1360
1361 cpl_test_eq(l, 10);
1362
1363 for(cpl_size i = 0; i < l; ++i){
1364 cpl_test_eq(x[i], x_f[i]);
1365 cpl_test_eq(y1[i], y1_f[i]);
1366 cpl_test_eq(y2[i], y2_f[i]);
1367 }
1368 }
1369 }
1370}
1371
1372/*----------------------------------------------------------------------------*/
1377/*----------------------------------------------------------------------------*/
1378void test_spectrum1D_resample_spectrum_interpolation_error_test(void){
1379
1380 double y[] = { 0, 1, 2, 1, 0, -1, -2, -1, 0, 1, 2, 1, 0, -1};
1381 double y_e[] = {.1, .2, .3, .4, .5, .6, .7, .8, .9, 1, 1.1, 1.2, 1.3, 1.4};
1382 double x[] = { 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14};
1383
1384 double x_r[] =
1385 {1.1, 2.2, 3.3, 4.4, 5.5, 6.6, 7.7, 8.8, 9.9, 10.9, 11.1, 12.2, 13.9};
1386
1387 cpl_size closer_idx[] =
1388 {0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12};
1389
1390 const cpl_size l = LENGTH_ARRAY(y);
1391 const cpl_size l2 = LENGTH_ARRAY(x_r);
1392 cpl_test_eq(LENGTH_ARRAY(y), LENGTH_ARRAY(y_e));
1393 cpl_test_eq(LENGTH_ARRAY(y), LENGTH_ARRAY(x));
1394
1395 cpl_image * flux = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1396 cpl_image * flux_e = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1397 cpl_array * wavelengths = cpl_array_new(l, CPL_TYPE_DOUBLE);
1398 cpl_array * wavelengths_resamp = cpl_array_new(l2, CPL_TYPE_DOUBLE);
1399
1400 for(cpl_size i = 0; i < l; ++i){
1401 cpl_image_set(flux, i + 1, 1, y[i]);
1402 cpl_image_set(flux_e, i + 1, 1, y_e[i]);
1403 cpl_array_set(wavelengths, i, x[i]);
1404 }
1405
1406 for(cpl_size i = 0; i < l2; ++i){
1407 cpl_array_set(wavelengths_resamp, i, x_r[i]);
1408 }
1409
1410 hdrl_spectrum1D * sp1 = hdrl_spectrum1D_create(flux, flux_e, wavelengths,
1411 hdrl_spectrum1D_wave_scale_linear);
1412
1413 cpl_array_delete(wavelengths);
1414 cpl_image_delete(flux);
1415 cpl_image_delete(flux_e);
1416
1417 hdrl_spectrum1D_wavelength wl = hdrl_spectrum1D_get_wavelength(sp1);
1418
1419 wl.wavelength = wavelengths_resamp;
1420 wl.bpm = NULL;
1421
1423 (hdrl_spectrum1D_interp_akima);
1424
1425 hdrl_spectrum1D * sp2 = hdrl_spectrum1D_resample(sp1, &wl, pars);
1426
1427 hdrl_parameter_delete(pars); pars = NULL;
1428
1429 /*check that for each position in x_r[i] we linearly interpolate the noise*/
1430 for(cpl_size i = 0; i < l2; ++i){
1431
1432 hdrl_error_t err = hdrl_spectrum1D_get_flux_value(sp2, i, NULL).error;
1433 const cpl_size idx = closer_idx[i];
1434 hdrl_error_t err_ori = pow(y_e[idx], 2.0) * fabs(x[idx + 1] - x_r[i]) +
1435 pow(y_e[idx + 1], 2.0) * fabs(x[idx] - x_r[i]);
1436
1437 err_ori = sqrt(err_ori);
1438
1439 cpl_test_abs(err, err_ori, HDRL_DELTA_COMPARE_VALUE_ABS);
1440 hdrl_data_t w = hdrl_spectrum1D_get_wavelength_value(sp2, i, NULL);
1441 cpl_test_abs(w, x_r[i], HDRL_DELTA_COMPARE_VALUE_ABS);
1442 }
1443
1446 cpl_array_delete(wavelengths_resamp);
1447}
1448
1449/*----------------------------------------------------------------------------*/
1454/*----------------------------------------------------------------------------*/
1455void test_spectrum1D_resample_spectrum_fit_error_test_error_interpol(void){
1456
1457 double y[] = { 0, 1, 2, 1, 0, -1, -2, -1, 0, 1, 2, 1, 0, -1};
1458 double y_e[] = {.1, .2, .3, .4, .5, .6, .7, .8, .9, 1, 1.1, 1.2, 1.3, 1.4};
1459 double x[] = { 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14};
1460
1461 double x_r[] =
1462 {1.1, 2.2, 3.3, 4.4, 5.5, 6.6, 7.7, 8.8, 9.9, 10.9, 11.1, 12.2, 13.9};
1463
1464 cpl_size closer_idx[] =
1465 {0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12};
1466
1467 const cpl_size l = LENGTH_ARRAY(y);
1468 const cpl_size l2 = LENGTH_ARRAY(x_r);
1469 cpl_test_eq(LENGTH_ARRAY(y), LENGTH_ARRAY(y_e));
1470 cpl_test_eq(LENGTH_ARRAY(y), LENGTH_ARRAY(x));
1471
1472 cpl_image * flux = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1473 cpl_image * flux_e = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1474 cpl_array * wavelengths = cpl_array_new(l, CPL_TYPE_DOUBLE);
1475 cpl_array * wavelengths_resamp = cpl_array_new(l2, CPL_TYPE_DOUBLE);
1476
1477 for(cpl_size i = 0; i < l; ++i){
1478 cpl_image_set(flux, i + 1, 1, y[i]);
1479 cpl_image_set(flux_e, i + 1, 1, y_e[i]);
1480 cpl_array_set(wavelengths, i, x[i]);
1481 }
1482
1483 for(cpl_size i = 0; i < l2; ++i){
1484 cpl_array_set(wavelengths_resamp, i, x_r[i]);
1485 }
1486
1487 hdrl_spectrum1D * sp1 = hdrl_spectrum1D_create(flux, flux_e, wavelengths,
1488 hdrl_spectrum1D_wave_scale_linear);
1489
1490 cpl_array_delete(wavelengths);
1491 cpl_image_delete(flux);
1492 cpl_image_delete(flux_e);
1493
1494 hdrl_spectrum1D_wavelength wl = hdrl_spectrum1D_get_wavelength(sp1);
1495
1496 wl.wavelength = wavelengths_resamp;
1497 wl.bpm = NULL;
1498
1499 hdrl_parameter * pars = hdrl_spectrum1D_resample_fit_parameter_create(2, 5);
1500
1501
1502 hdrl_spectrum1D * sp2 = hdrl_spectrum1D_resample(sp1, &wl, pars);
1503
1504 hdrl_parameter_delete(pars); pars = NULL;
1505
1506 /*check that for each position in x_r[i] we linearly interpolate the noise*/
1507 for(cpl_size i = 0; i < l2; ++i){
1508
1509 hdrl_error_t err = hdrl_spectrum1D_get_flux_value(sp2, i, NULL).error;
1510 const cpl_size idx = closer_idx[i];
1511 hdrl_error_t err_ori = pow(y_e[idx], 2.0) * fabs(x[idx + 1] - x_r[i]) +
1512 pow(y_e[idx + 1], 2.0) * fabs(x[idx] - x_r[i]);
1513
1514 err_ori = sqrt(err_ori);
1515
1516 cpl_test_abs(err, err_ori, HDRL_DELTA_COMPARE_VALUE_ABS);
1517 hdrl_data_t w = hdrl_spectrum1D_get_wavelength_value(sp2, i, NULL);
1518 cpl_test_abs(w, x_r[i], HDRL_DELTA_COMPARE_VALUE_ABS);
1519 }
1520
1523 cpl_array_delete(wavelengths_resamp);
1524}
1525
1526
1527/*----------------------------------------------------------------------------*/
1532/*----------------------------------------------------------------------------*/
1533void test_spectrum1D_resample_spectrum_bpm(cpl_boolean interpolate){
1534
1535 double x[] = { 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14,
1536 15};
1537
1538 double x_r[] =
1539 {1.1, 2.2, 3.3, 4.4, 5.5, 6.6, 7.7, 8.8, 9.9, 10.9, 11.1, 12.2, 13.9, 14.1};
1540
1541 const cpl_size l = LENGTH_ARRAY(x);
1542 const cpl_size l2 = LENGTH_ARRAY(x_r);
1543
1544 cpl_image * flux = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1545 cpl_array * wavelengths = cpl_array_new(l, CPL_TYPE_DOUBLE);
1546 cpl_array * wavelengths_resamp = cpl_array_new(l2, CPL_TYPE_DOUBLE);
1547
1548 /* set all the bad pixels to 10 and all the others to 0. */
1549 for(cpl_size i = 0; i < l; ++i){
1550 if(i % 2 == 0)
1551 cpl_image_set(flux, i + 1, 1, 10);
1552 cpl_array_set(wavelengths, i, x[i]);
1553 }
1554
1555 for(cpl_size i = 0; i < l; i+=2)
1556 cpl_image_reject(flux, i + 1, 1);
1557
1558 for(cpl_size i = 0; i < l2; ++i){
1559 cpl_array_set(wavelengths_resamp, i, x_r[i]);
1560 }
1561
1562
1563 hdrl_spectrum1D * sp1 = hdrl_spectrum1D_create_error_free(flux, wavelengths,
1564 hdrl_spectrum1D_wave_scale_linear);
1565
1566 cpl_array_delete(wavelengths);
1567 cpl_image_delete(flux);
1568
1569 hdrl_spectrum1D_wavelength wl = hdrl_spectrum1D_get_wavelength(sp1);
1570
1571 wl.wavelength = wavelengths_resamp;
1572 wl.bpm = NULL;
1573
1574 hdrl_parameter * pars = NULL;
1575 if(interpolate)
1577 (hdrl_spectrum1D_interp_akima);
1578 else
1580
1581 hdrl_spectrum1D * sp2 = hdrl_spectrum1D_resample(sp1, &wl, pars);
1582
1583 hdrl_parameter_delete(pars); pars = NULL;
1584
1585 /* the first element in the interpolated must be rejected since the element
1586 * a wlength 1 is bad and we mark as bad every pixel outside the interval of the
1587 * source spectra. I.e. we truncate*/
1588 int rej = 0;
1589 hdrl_spectrum1D_get_flux_value(sp2, 0, &rej);
1590 cpl_test_eq(rej, 1);
1591
1592 rej = 0;
1593 hdrl_spectrum1D_get_flux_value(sp2, l2 - 1, &rej);
1594 cpl_test_eq(rej, 1);
1595
1596 /*check that all the pixels are zero, the non zero pixels of sp1 should have
1597 * not contributed to the interpolation, because there were all bad*/
1598 for(cpl_size i = 1; i < l2 - 1; ++i){
1599
1600 hdrl_data_t data = hdrl_spectrum1D_get_flux_value(sp2, i, NULL).data;
1601
1602 cpl_test_rel(data, 0.0, 1e-6);
1603 }
1604
1607 cpl_array_delete(wavelengths_resamp);
1608}
1609
1610/*----------------------------------------------------------------------------*/
1618/*----------------------------------------------------------------------------*/
1619void test_spectrum1D_resample_spectrum_fit_error_test_shift(cpl_boolean is_error_free){
1620
1621 const hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
1622 double y[] = { 0, 1, 2, 1, 0, -1, -2, -1, 0, 1, 2, 1, 0, -1};
1623 double y_e[] = {.1, .2, .3, .2, .1, .2, .3, .2, .1, .2, .3, .2, .1, .2};
1624 double x[] = { 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14};
1625
1626 const cpl_size l = LENGTH_ARRAY(x);
1627
1628 cpl_image * flux = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1629 cpl_image * flux_e = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1630 cpl_array * wavelengths = cpl_array_new(l, CPL_TYPE_DOUBLE);
1631
1632
1633 for(cpl_size i = 0; i < l; ++i){
1634 cpl_image_set(flux, i + 1, 1, y[i]);
1635 if(!is_error_free)
1636 cpl_image_set(flux_e, i + 1, 1, y_e[i]);
1637 cpl_array_set(wavelengths, i, x[i]);
1638 }
1639
1640 cpl_array * wavelengths_resampled1 = cpl_array_new(l - 1, CPL_TYPE_DOUBLE);
1641 cpl_array * wavelengths_resampled2 = cpl_array_new(l - 2, CPL_TYPE_DOUBLE);
1642
1643 for(cpl_size i = 0; i < l - 1; ++i){
1644 const double d = cpl_array_get(wavelengths, i, NULL) + .2;
1645 cpl_array_set(wavelengths_resampled1, i, d);
1646 }
1647
1648 for(cpl_size i = 1; i < l - 1; ++i){
1649 const double d = cpl_array_get(wavelengths, i, NULL);
1650 cpl_array_set(wavelengths_resampled2, i - 1, d);
1651 }
1652
1653 hdrl_spectrum1D * sp1 = hdrl_spectrum1D_create(flux, flux_e,
1654 wavelengths, scale);
1655 cpl_test_nonnull(sp1);
1656
1657 const hdrl_spectrum1D_wavelength wl1 = {wavelengths_resampled1, NULL, scale};
1658 const hdrl_spectrum1D_wavelength wl2 = {wavelengths_resampled2, NULL, scale};
1659
1660 hdrl_parameter * pars =
1662
1663 hdrl_spectrum1D * sp2 = hdrl_spectrum1D_resample(sp1, &wl1, pars);
1664 cpl_test_nonnull(sp2);
1665
1666 hdrl_spectrum1D * sp3 = hdrl_spectrum1D_resample(sp2, &wl2, pars);
1667 cpl_test_nonnull(sp3);
1668
1669 for(cpl_size i = 1; i < l - 2; ++i){
1670 const hdrl_value v1 =
1671 hdrl_spectrum1D_get_flux_value(sp1, i, NULL);
1672 const hdrl_value v3 =
1673 hdrl_spectrum1D_get_flux_value(sp3, i - 1, NULL);
1674
1675 cpl_test_eq(v1.error, v3.error);
1676 cpl_test_abs(v1.data, v3.data, .5);
1677 }
1678
1680 cpl_image_delete(flux);
1681 cpl_image_delete(flux_e);
1682 cpl_array_delete(wavelengths);
1683 cpl_array_delete(wavelengths_resampled1);
1684 cpl_array_delete(wavelengths_resampled2);
1685
1689}
1690
1691static inline double
1692func(const double t)
1693{
1694 double x = sin(10.0 * t);
1695 return exp(x*x*x);
1696}
1697
1698/*----------------------------------------------------------------------------*/
1702/*----------------------------------------------------------------------------*/
1703void test_spectrum1D_resample_spectrum_fit_windowed(void){
1704
1705 const cpl_size nblocks = 5e2;
1706 const cpl_size length = 100e3;
1707 const cpl_size window = length / nblocks;
1708
1709 const double dt = 1.0 / (length - 1);
1710
1711 const hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
1712
1713 hdrl_parameter * pars_no_win =
1715
1716 hdrl_parameter * pars_win =
1718 (4, 6, window, 1.2);
1719
1720 cpl_array * lambdas = cpl_array_new(length, CPL_TYPE_DOUBLE);
1721 cpl_array * lambdas_resampled = cpl_array_new(length - 2, CPL_TYPE_DOUBLE);
1722
1723 cpl_image * flux = cpl_image_new(length, 1, CPL_TYPE_DOUBLE);
1724 cpl_image * flux_real = cpl_image_new(length - 2, 1, CPL_TYPE_DOUBLE);
1725
1726 for(cpl_size i = 0; i < length; ++i){
1727 const double d = i * dt;
1728 cpl_array_set(lambdas, i, d);
1729 cpl_image_set(flux, i + 1, 1, func(d));
1730 }
1731
1732 for(cpl_size i = 0; i < length - 2; ++i){
1733 const double d = (i + .5) * dt;
1734 cpl_array_set(lambdas_resampled, i, d);
1735 cpl_image_set(flux_real, i + 1, 1, func(d));
1736 }
1737
1738 hdrl_spectrum1D_wavelength wav = {lambdas_resampled, NULL, scale};
1739
1740 hdrl_spectrum1D * sp_source =
1741 hdrl_spectrum1D_create_error_free(flux, lambdas, scale);
1742
1743 hdrl_spectrum1D * sp_ideal_resampled =
1745 (flux_real, lambdas_resampled, scale);
1746
1747 hdrl_spectrum1D * sp_win_resampled = hdrl_spectrum1D_resample
1748 (sp_source, &wav, pars_win);
1749
1750 hdrl_spectrum1D * sp_resampled = hdrl_spectrum1D_resample
1751 (sp_source, &wav, pars_no_win);
1752
1753 for(cpl_size i = 0; i < length - 2; ++i){
1754 hdrl_data_t v_ideal = hdrl_spectrum1D_get_flux_value
1755 (sp_ideal_resampled, i, NULL).data;
1756 hdrl_data_t v_win = hdrl_spectrum1D_get_flux_value
1757 (sp_win_resampled, i, NULL).data;
1758 hdrl_data_t v_no_win = hdrl_spectrum1D_get_flux_value
1759 (sp_resampled, i, NULL).data;
1760
1761 cpl_test_rel(v_ideal, v_win, 1.5e-1);
1762 cpl_test_rel(v_ideal, v_no_win, 1.5e-1);
1763 }
1764
1765 hdrl_parameter_delete(pars_no_win);
1766 hdrl_parameter_delete(pars_win);
1767
1768 cpl_array_delete(lambdas);
1769 cpl_array_delete(lambdas_resampled);
1770
1771 cpl_image_delete(flux);
1772 cpl_image_delete(flux_real);
1773 hdrl_spectrum1D_delete(&sp_ideal_resampled);
1774 hdrl_spectrum1D_delete(&sp_win_resampled);
1775 hdrl_spectrum1D_delete(&sp_resampled);
1776 hdrl_spectrum1D_delete(&sp_source);
1777}
1778
1779/*----------------------------------------------------------------------------*/
1783/*----------------------------------------------------------------------------*/
1784void test_spectrum1D_wavelength_select(void){
1785
1786 const hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
1787 double y[] = { 0, 1, 2, 1, 0, -1, -2, -1, 0, 1, 2, 1, 0, -1};
1788 double y_e[] = {.1, .2, .3, .2, .1, .2, .3, .2, .1, .2, .3, .2, .1, .2};
1789 double x[] = { 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14};
1790
1791 const cpl_size l = LENGTH_ARRAY(x);
1792
1793 cpl_image * flux = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1794 cpl_image * flux_e = cpl_image_new(l, 1, CPL_TYPE_DOUBLE);
1795 cpl_array * wavelengths = cpl_array_new(l, CPL_TYPE_DOUBLE);
1796
1797
1798 for(cpl_size i = 0; i < l; ++i){
1799 cpl_image_set(flux, i + 1, 1, y[i]);
1800 cpl_image_set(flux_e, i + 1, 1, y_e[i]);
1801 cpl_array_set(wavelengths, i, x[i]);
1802 }
1803
1804 cpl_image_reject(flux, 5, 1);
1805
1806 hdrl_spectrum1D * sp1 =
1807 hdrl_spectrum1D_create(flux, flux_e, wavelengths, scale);
1808
1809 hdrl_spectrum1D * sp2 = select_window(sp1, (sel_window){3, 10, CPL_TRUE});
1810
1811 cpl_test_eq(8, hdrl_spectrum1D_get_size(sp2));
1812 {
1813 int rej = 0;
1814 hdrl_spectrum1D_get_flux_value(sp2, 2, &rej);
1815 cpl_test_eq(rej, 1);
1816 }
1817 const cpl_size l2 = hdrl_spectrum1D_get_size(sp2);
1818 for(cpl_size i = 0; i < l2; ++i){
1819
1820 if(i == 2) continue;
1821
1822 int rej1_1 = 0;
1823 int rej2_1 = 0;
1824
1825 const hdrl_data_t w1 =
1826 hdrl_spectrum1D_get_wavelength_value(sp1, i + 2, &rej1_1);
1827 const hdrl_data_t w2 =
1828 hdrl_spectrum1D_get_wavelength_value(sp2, i, &rej2_1);
1829
1830 cpl_test_eq(rej1_1, 0);
1831 cpl_test_eq(rej2_1, 0);
1832
1833 cpl_test_eq(w1, w2);
1834
1835 int rej1_2 = 0;
1836 int rej2_2 = 0;
1837 const hdrl_value s1 = hdrl_spectrum1D_get_flux_value(sp1, i + 2, &rej1_2);
1838 const hdrl_value s2 = hdrl_spectrum1D_get_flux_value(sp2, i, &rej2_2);
1839
1840 cpl_test_eq(rej1_2, 0);
1841 cpl_test_eq(rej2_2, 0);
1842
1843 cpl_test_rel(s1.data, s2.data, HDRL_DELTA_COMPARE_VALUE);
1844 cpl_test_rel(s1.error, s2.error, HDRL_DELTA_COMPARE_VALUE);
1845 }
1846
1848
1849 hdrl_spectrum1D * sp3 = select_window(sp1, (sel_window){3, 10, CPL_FALSE});
1850
1851 cpl_test_eq(6, hdrl_spectrum1D_get_size(sp3));
1852
1853 cpl_size idxes[] = {0, 1, 10, 11, 12, 13};
1854
1855 for(cpl_size i = 0; i < 6; ++i){
1856
1857 int rej1_1 = 0;
1858 int rej2_1 = 0;
1859
1860 const hdrl_data_t w1 =
1861 hdrl_spectrum1D_get_wavelength_value(sp1, idxes[i], &rej1_1);
1862 const hdrl_data_t w2 =
1863 hdrl_spectrum1D_get_wavelength_value(sp3, i, &rej2_1);
1864
1865 cpl_test_eq(rej1_1, 0);
1866 cpl_test_eq(rej2_1, 0);
1867
1868 cpl_test_eq(w1, w2);
1869
1870 int rej1_2 = 0;
1871 int rej2_2 = 0;
1872 const hdrl_value s1 =
1873 hdrl_spectrum1D_get_flux_value(sp1, idxes[i], &rej1_2);
1874 const hdrl_value s2 = hdrl_spectrum1D_get_flux_value(sp3, i, &rej2_2);
1875
1876 cpl_test_eq(rej1_2, 0);
1877 cpl_test_eq(rej2_2, 0);
1878
1879 cpl_test_rel(s1.data, s2.data, HDRL_DELTA_COMPARE_VALUE);
1880 cpl_test_rel(s1.error, s2.error, HDRL_DELTA_COMPARE_VALUE);
1881 }
1882
1885 cpl_image_delete(flux);
1886 cpl_image_delete(flux_e);
1887 cpl_array_delete(wavelengths);
1888
1889 /* Check for the function that rejects some elements inside the spectrum */
1890 const cpl_size sz = 10;
1891
1892 hdrl_spectrum1D * sp =
1893 get_spectrum1D_sin_shuffled(sz, 3, CPL_FALSE, NULL);
1894
1895 cpl_array * arr = cpl_array_new(sz, CPL_TYPE_INT);
1896
1897 for(cpl_size i = 0; i < sz; i++){
1898 cpl_array_set(arr, i, i % 2);
1899 }
1900
1901 hdrl_spectrum1D * sp_r1 =
1903
1904 for(cpl_size i = 0; i < sz; i++){
1905 int rej = 0;
1906 hdrl_value d = hdrl_spectrum1D_get_flux_value(sp_r1, i, &rej);
1907
1908 if(i % 2 == 1){
1909 cpl_test(rej);
1910 }
1911 else{
1912 hdrl_value d1 = hdrl_spectrum1D_get_flux_value(sp, i, NULL);
1913 cpl_test_rel(d.data, d1.data, HDRL_EPS_DATA);
1914 cpl_test_rel(d.error, d1.error, HDRL_EPS_DATA);
1915 }
1916 }
1917
1918 cpl_array_delete(arr);
1920 hdrl_spectrum1D_delete(&sp_r1);
1921}
1922
1923/*----------------------------------------------------------------------------*/
1927/*----------------------------------------------------------------------------*/
1928void test_spectrum1D_test_uniformly_sampled(void){
1929
1930 const hdrl_spectrum1D_wave_scale scale = hdrl_spectrum1D_wave_scale_linear;
1931
1932 const cpl_size sz = 100;
1933 hdrl_spectrum1D * sp =
1934 get_spectrum1D_sin_shuffled(sz, 3, CPL_FALSE, NULL);
1935
1936 double delta = 0.0;
1937 cpl_boolean is_OK = hdrl_spectrum1D_is_uniformly_sampled(sp, &delta);
1938
1939 cpl_test(is_OK);
1940 cpl_test_abs(delta, 2 * CPL_MATH_PI / sz, HDRL_EPS_DATA);
1941
1943
1944 cpl_array * arr = cpl_array_new(sz, CPL_TYPE_DOUBLE);
1945 cpl_image * flx = cpl_image_new(sz, 1, CPL_TYPE_DOUBLE);
1946 for(cpl_size i = 0; i < sz; ++i){
1947 cpl_array_set(arr, i, i + 1);
1948 cpl_image_set(flx, i + 1, 1, .1);
1949 }
1950
1951 cpl_array_set(arr, 4, 5.1);
1952 sp = hdrl_spectrum1D_create_error_free(flx, arr, scale);
1953 is_OK = hdrl_spectrum1D_is_uniformly_sampled(sp, &delta);
1954 cpl_test(!is_OK);
1956
1957 cpl_array_set(arr, 4, 5);
1958 sp = hdrl_spectrum1D_create_error_free(flx, arr, scale);
1959 is_OK = hdrl_spectrum1D_is_uniformly_sampled(sp, &delta);
1960 cpl_test(is_OK);
1961 cpl_test_rel(delta, 1.0, HDRL_EPS_DATA);
1963
1964 cpl_image_delete(flx);
1965 cpl_array_delete(arr);
1966}
1967
1968/*----------------------------------------------------------------------------*/
1972/*----------------------------------------------------------------------------*/
1973void test_spectrum1Dlist(void){
1974
1975 {
1976 hdrl_spectrum1Dlist * list = hdrl_spectrum1Dlist_new();
1978 }
1979
1980 const cpl_size sz = 40;
1981 {
1982 hdrl_spectrum1Dlist * list = hdrl_spectrum1Dlist_new();
1983
1984 hdrl_spectrum1D * s1 =
1985 get_spectrum1D_sin_shuffled(sz, 4, CPL_TRUE, NULL);
1986 hdrl_spectrum1D * s2 =
1987 get_spectrum1D_sin_shuffled(sz, 4, CPL_TRUE, NULL);
1988 hdrl_spectrum1D * s3 =
1989 get_spectrum1D_sin_shuffled(sz, 4, CPL_TRUE, NULL);
1990
1991 hdrl_spectrum1Dlist_set(list, s1, 0);
1992 hdrl_spectrum1Dlist_set(list, s2, 1);
1993 hdrl_spectrum1Dlist_set(list, s3, 2);
1994
1995 cpl_test_eq(hdrl_spectrum1Dlist_get_size(list), 3);
1996 cpl_test_eq(list->capacity, 4);
1997
1998 hdrl_spectrum1D * s22 = hdrl_spectrum1Dlist_get(list, 1);
1999 const hdrl_spectrum1D * s22_c = hdrl_spectrum1Dlist_get_const(list, 1);
2000
2001 cpl_test_eq_ptr(s2, s22);
2002 cpl_test_eq_ptr(s2, s22_c);
2003
2004 hdrl_spectrum1D * s22_u = hdrl_spectrum1Dlist_unset(list, 1);
2005 cpl_test_eq_ptr(s22_u, s22_c);
2006
2007 cpl_test_eq(hdrl_spectrum1Dlist_get_size(list), 2);
2008
2009 hdrl_spectrum1D_delete(&s22_u);
2010
2011 cpl_size i = 0;
2012 while(hdrl_spectrum1Dlist_get_size(list) > 0) {
2013 hdrl_spectrum1D * s = hdrl_spectrum1Dlist_unset(list, 0);
2014
2015 if(i == 0)
2016 cpl_test_eq_ptr(s, s1);
2017 else
2018 cpl_test_eq_ptr(s, s3);
2019
2021 i++;
2022 }
2023
2024 cpl_test_eq(hdrl_spectrum1Dlist_get_size(list), 0);
2025 cpl_test_eq(list->capacity, 0);
2026 cpl_test_eq_ptr(list->spectra, NULL);
2027
2029 }
2030}
2031
2032static inline cpl_array *
2033get_waves(double start_wave, cpl_size sz, double step){
2034
2035 cpl_array * to_ret = cpl_array_new(sz, HDRL_TYPE_DATA);
2036
2037 for(cpl_size i = 0; i < sz; ++i){
2038 cpl_array_set(to_ret, i, start_wave + step * i);
2039 }
2040
2041 return to_ret;
2042}
2043
2044static inline hdrl_spectrum1D *
2045generate_stair_spectrum(int start, int stop, double start_wave, double step_wave){
2046
2047 const cpl_size sz = stop - start + 1;
2048 cpl_array * wave = get_waves(start_wave, sz, step_wave);
2049 cpl_image * flx = cpl_image_new(sz, 1, CPL_TYPE_DOUBLE);
2050 cpl_image * flx_e = cpl_image_new(sz, 1, CPL_TYPE_DOUBLE);
2051
2052
2053 for(cpl_size i = 0; i < sz; ++i){
2054 const double f = start + i;
2055 cpl_image_set(flx, i + 1, 1, f);
2056 cpl_image_set(flx_e, i + 1, 1, f / 10.0);
2057 }
2058
2059 hdrl_spectrum1D * to_ret = hdrl_spectrum1D_create(flx, flx_e, wave,
2060 hdrl_spectrum1D_wave_scale_linear);
2061
2062 cpl_image_delete(flx);
2063 cpl_image_delete(flx_e);
2064 cpl_array_delete(wave);
2065
2066
2067 return to_ret;
2068}
2069
2070cpl_boolean contains(const cpl_array * arr, cpl_size idx){
2071 for(cpl_size i = 0; i < cpl_array_get_size(arr); ++i){
2072 const cpl_size this_idx = (cpl_size)cpl_array_get(arr, i, NULL);
2073 if(this_idx == idx) return CPL_TRUE;
2074 }
2075 return CPL_FALSE;
2076}
2077
2078static inline hdrl_spectrum1D *
2079generate_bad_stair_spectrum(int start, int stop, double start_wave,
2080 double step_wave, cpl_array * bad_idxes){
2081
2082 const cpl_size sz = stop - start + 1;
2083 cpl_array * wave = get_waves(start_wave, sz, step_wave);
2084 cpl_image * flx = cpl_image_new(sz, 1, CPL_TYPE_DOUBLE);
2085 cpl_image * flx_e = cpl_image_new(sz, 1, CPL_TYPE_DOUBLE);
2086
2087
2088 for(cpl_size i = 0; i < sz; ++i){
2089
2090 if(contains(bad_idxes, i))
2091 {
2092 cpl_image_reject(flx, i + 1, 1);
2093 cpl_image_reject(flx_e, i + 1, 1);
2094 continue;
2095 }
2096
2097 const double f = start + i;
2098 cpl_image_set(flx, i + 1, 1, f);
2099 cpl_image_set(flx_e, i + 1, 1, f / 10.0);
2100 }
2101
2102 hdrl_spectrum1D * to_ret = hdrl_spectrum1D_create(flx, flx_e, wave,
2103 hdrl_spectrum1D_wave_scale_linear);
2104
2105 cpl_image_delete(flx);
2106 cpl_image_delete(flx_e);
2107 cpl_array_delete(wave);
2108
2109
2110 return to_ret;
2111}
2112
2113
2114static inline hdrl_spectrum1D *
2115generate_stair_spectrum_shuffled(int start, int stop, double start_wave,
2116 double step_wave, cpl_size * idxes){
2117
2118 const cpl_size sz = stop - start + 1;
2119 cpl_array * wave = get_waves(start_wave, sz, step_wave);
2120 cpl_array * wave_s = cpl_array_duplicate(wave);
2121
2122 cpl_image * flx = cpl_image_new(sz, 1, CPL_TYPE_DOUBLE);
2123 cpl_image * flx_e = cpl_image_new(sz, 1, CPL_TYPE_DOUBLE);
2124
2125 for(cpl_size i = 0; i < sz; ++i){
2126 const double f = start + i;
2127 const cpl_size dest_idx = idxes[i];
2128 cpl_image_set(flx, dest_idx + 1, 1, f);
2129 cpl_image_set(flx_e, dest_idx + 1, 1, f / 10.0);
2130
2131 const double w = cpl_array_get(wave, i, NULL);
2132 cpl_array_set(wave_s, dest_idx, w);
2133 }
2134
2135 hdrl_spectrum1D * to_ret = hdrl_spectrum1D_create(flx, flx_e, wave_s,
2136 hdrl_spectrum1D_wave_scale_linear);
2137
2138 cpl_image_delete(flx);
2139 cpl_image_delete(flx_e);
2140 cpl_array_delete(wave);
2141 cpl_array_delete(wave_s);
2142
2143 return to_ret;
2144}
2145
2146void test1(void){
2148
2149 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 8, 20.0, 2.0);
2150
2151 cpl_array * wavs_integrate = get_waves(21, 9, 1);
2152 hdrl_spectrum1D * integrated_s =
2153 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2154
2155 /* first element: the initial bin is asimmetric, tested as
2156 * special case*/
2157 {
2158 int rej = 0;
2159 hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, 0, &rej);
2160 cpl_test_eq(rej, 0);
2161 cpl_test_rel(f.data, 2.0, HDRL_DELTA_COMPARE_VALUE);
2162 cpl_test_rel(f.error, f.data / 10.0, HDRL_DELTA_COMPARE_VALUE);
2163 }
2164
2165 const cpl_size int_size = cpl_array_get_size(wavs_integrate);
2166
2167 /* last element: the final bin is asimmetric, tested as
2168 * special case*/
2169 {
2170 int rej = 0;
2171 hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, int_size - 1, &rej);
2172 cpl_test_eq(rej, 0);
2173 cpl_test_rel(f.data, 5.0, HDRL_DELTA_COMPARE_VALUE);
2174 cpl_test_rel(f.error, f.data / 10.0, HDRL_DELTA_COMPARE_VALUE);
2175 }
2176
2177 double src_flx = 2.0;
2178 for(cpl_size i = 1; i < int_size - 1; ++i){
2179
2180 int rej = 0;
2181 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2182
2183 cpl_test_eq(rej, 0);
2184
2185 if((i % 2) == 1){
2186 /*odd elements: they are in the center of the source bin*/
2187 cpl_test_rel(f.data, src_flx, HDRL_DELTA_COMPARE_VALUE);
2188 cpl_test_rel(f.error, src_flx / 10.0, HDRL_DELTA_COMPARE_VALUE);
2189
2190 }
2191 else{
2192 /*even elements: they are between bins*/
2193 const double el = src_flx + .5;
2194 cpl_test_rel(f.data, el, HDRL_DELTA_COMPARE_VALUE);
2195
2196 const double el_e = pow(src_flx, 2.0) + pow((src_flx + 1.0), 2.0);
2197 cpl_test_rel(f.error, sqrt(el_e) / (10.0 * sqrt(2.0)), HDRL_DELTA_COMPARE_VALUE);
2198 src_flx++;
2199 }
2200 }
2201
2203 hdrl_spectrum1D_delete(&ori_s);
2204 hdrl_spectrum1D_delete(&integrated_s);
2205 cpl_array_delete(wavs_integrate);
2206}
2207
2208void test2(void){
2209
2211
2212 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 8, 20.0, 2.0);
2213
2214 cpl_array * wavs_integrate = get_waves(20, 15, 1);
2215 hdrl_spectrum1D * integrated_s =
2216 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2217
2218 const cpl_size int_size = cpl_array_get_size(wavs_integrate);
2219
2220 double src_flx = 1.0;
2221 for(cpl_size i = 0; i < int_size; ++i){
2222
2223 int rej = 0;
2224 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2225
2226 cpl_test_eq(rej, 0);
2227
2228 if((i % 2) == 0){
2229 /*odd elements: they are in the center of the source bin*/
2230 cpl_test_rel(f.data, src_flx, HDRL_DELTA_COMPARE_VALUE);
2231 cpl_test_rel(f.error, src_flx / 10.0, HDRL_DELTA_COMPARE_VALUE);
2232
2233 }
2234 else{
2235 /*even elements: they are between bins*/
2236 const double el = src_flx + .5;
2237 cpl_test_rel(f.data, el, HDRL_DELTA_COMPARE_VALUE);
2238
2239 const double el_e = pow(src_flx, 2.0) + pow((src_flx + 1.0), 2.0);
2240 cpl_test_rel(f.error, sqrt(el_e) / (10.0 * sqrt(2.0)), HDRL_DELTA_COMPARE_VALUE);
2241 src_flx++;
2242 }
2243 }
2244
2246 hdrl_spectrum1D_delete(&ori_s);
2247 hdrl_spectrum1D_delete(&integrated_s);
2248 cpl_array_delete(wavs_integrate);
2249}
2250
2251void test3(void){
2252
2254
2255 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 8, 20.0, 2.0);
2256
2257 cpl_array * wavs_integrate = get_waves(19, 17, 1);
2258 hdrl_spectrum1D * integrated_s =
2259 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2260
2261 const cpl_size int_size = cpl_array_get_size(wavs_integrate);
2262
2263 /*2 bins are rejected, 1 because is completely outside, the other half outside*/
2264 for(cpl_size i = 0; i < 2; ++i){
2265 int rej = 0;
2266 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2267 cpl_test_eq(rej, 1);
2268 cpl_test(isnan(f.data));
2269 cpl_test(isnan(f.error));
2270 }
2271
2272 /*2 bins are rejected, 1 because is completely outside, the other half outside*/
2273 for(cpl_size i = int_size - 2; i < int_size; ++i){
2274 int rej = 0;
2275 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2276 cpl_test_eq(rej, 1);
2277 cpl_test(isnan(f.data));
2278 cpl_test(isnan(f.error));
2279 }
2280
2281 double src_flx = 1.0;
2282 for(cpl_size i = 2; i < int_size - 2; ++i){
2283
2284 int rej = 0;
2285 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2286
2287 cpl_test_eq(rej, 0);
2288
2289 if((i % 2) == 1){
2290 /*odd elements: they are in the center of the source bin*/
2291 cpl_test_rel(f.data, src_flx, HDRL_DELTA_COMPARE_VALUE);
2292 cpl_test_rel(f.error, src_flx / 10.0, HDRL_DELTA_COMPARE_VALUE);
2293
2294 }
2295 else{
2296 /*even elements: they are between bins*/
2297 const double el = src_flx + .5;
2298 cpl_test_rel(f.data, el, HDRL_DELTA_COMPARE_VALUE);
2299
2300 const double el_e = pow(src_flx, 2.0) + pow((src_flx + 1.0), 2.0);
2301 cpl_test_rel(f.error, sqrt(el_e) / (10.0 * sqrt(2.0)), HDRL_DELTA_COMPARE_VALUE);
2302 src_flx++;
2303 }
2304 }
2305
2307 hdrl_spectrum1D_delete(&ori_s);
2308 hdrl_spectrum1D_delete(&integrated_s);
2309 cpl_array_delete(wavs_integrate);
2310}
2311
2312static inline cpl_array *
2313shuffle(const cpl_array * ori, const cpl_size * idxes){
2314 cpl_array * to_ret = cpl_array_duplicate(ori);
2315 cpl_size sz = cpl_array_get_size(ori);
2316 for(cpl_size i = 0; i < sz; ++i){
2317 const double w = cpl_array_get(ori, i, NULL);
2318 cpl_array_set(to_ret, idxes[i], w);
2319 }
2320 return to_ret;
2321}
2322
2323void test4(void){
2325
2326 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 8, 20.0, 2.0);
2327
2328 hdrl_spectrum1D * ori_s_shuffled = generate_stair_spectrum_shuffled(1, 8, 20.0, 2.0,
2329 (cpl_size[]){3, 2, 1, 4, 7, 6, 0, 5});
2330
2331 cpl_array * wavs_integrate = get_waves(21, 9, 1);
2332
2333 cpl_size shuffles[] = {1, 3, 5, 0, 8, 7, 6, 2, 4};
2334 cpl_array * wavs_integrate_shuffle = shuffle(wavs_integrate,
2335 shuffles);
2336
2337 hdrl_spectrum1D * integrated_s_shuffled =
2338 hdrl_spectrum1D_resample_on_array(ori_s_shuffled, wavs_integrate_shuffle, par);
2339
2340 hdrl_spectrum1D * integrated_s=
2341 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2342
2343 const cpl_size sz = hdrl_spectrum1D_get_size(integrated_s_shuffled);
2344 cpl_test_eq(sz, cpl_array_get_size(wavs_integrate));
2345 cpl_test_eq(sz, hdrl_spectrum1D_get_size(integrated_s));
2346
2347 for(cpl_size i = 0; i < sz; ++i){
2348 const double wav_s = hdrl_spectrum1D_get_wavelength_value(integrated_s,
2349 i, NULL);
2350 const double wav = cpl_array_get(wavs_integrate, i, NULL);
2351 cpl_test_rel(wav_s, wav, HDRL_DELTA_COMPARE_VALUE);
2352 }
2353
2354 for(cpl_size i = 0; i < sz; ++i){
2355 const double wav_s = hdrl_spectrum1D_get_wavelength_value(integrated_s_shuffled,
2356 i, NULL);
2357 const double wav = cpl_array_get(wavs_integrate_shuffle, i, NULL);
2358 cpl_test_rel(wav_s, wav, HDRL_DELTA_COMPARE_VALUE);
2359 }
2360
2361 for(cpl_size i = 0; i < sz; ++i){
2362 int rej = 0, rej_shuffled = 0;
2363
2364 const hdrl_value flx = hdrl_spectrum1D_get_flux_value(integrated_s,
2365 i, &rej);
2366 const hdrl_value flx_shuffled = hdrl_spectrum1D_get_flux_value(integrated_s_shuffled,
2367 shuffles[i], &rej_shuffled);
2368
2369 const hdrl_data_t wav = hdrl_spectrum1D_get_wavelength_value(integrated_s,
2370 i, NULL);
2371
2372 const hdrl_data_t wav_shuffled = hdrl_spectrum1D_get_wavelength_value(integrated_s_shuffled,
2373 shuffles[i], NULL);
2374
2375 cpl_test_rel(wav, wav_shuffled, HDRL_DELTA_COMPARE_VALUE);
2376 cpl_test_eq(rej, rej_shuffled);
2377 cpl_test_rel(flx.data, flx_shuffled.data, HDRL_DELTA_COMPARE_VALUE);
2378 cpl_test_rel(flx.error, flx_shuffled.error, HDRL_DELTA_COMPARE_VALUE);
2379 }
2380
2382 hdrl_spectrum1D_delete(&ori_s);
2383 hdrl_spectrum1D_delete(&ori_s_shuffled);
2384 hdrl_spectrum1D_delete(&integrated_s);
2385 hdrl_spectrum1D_delete(&integrated_s_shuffled);
2386 cpl_array_delete(wavs_integrate);
2387 cpl_array_delete(wavs_integrate_shuffle);
2388}
2389
2390void test5(void){
2391
2393
2394 cpl_array * bad_idxes = cpl_array_new(4, CPL_TYPE_INT);
2395
2396 cpl_array_set(bad_idxes, 0, 0);
2397 cpl_array_set(bad_idxes, 1, 7);
2398 cpl_array_set(bad_idxes, 2, 2);
2399 cpl_array_set(bad_idxes, 3, 5);
2400
2401 hdrl_spectrum1D * ori_s = generate_bad_stair_spectrum(1, 8, 20.0, 2.0 , bad_idxes);
2402
2403 cpl_array * wavs_integrate = get_waves(19, 17, 1);
2404 hdrl_spectrum1D * integrated_s =
2405 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2406
2407 const cpl_size int_size = cpl_array_get_size(wavs_integrate);
2408
2409 /*3 bins are rejected, 1 because is completely outside, the other half outside,
2410 * the third in a bad pixel*/
2411 for(cpl_size i = 0; i < 3; ++i){
2412 int rej = 0;
2413 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2414 cpl_test_eq(rej, 1);
2415 cpl_test(isnan(f.data));
2416 cpl_test(isnan(f.error));
2417 }
2418
2419 /*3 bins are rejected, 1 because is completely outside, the other half outside,
2420 * the third in a bad pixel*/
2421 for(cpl_size i = int_size - 3; i < int_size; ++i){
2422 int rej = 0;
2423 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2424 cpl_test_eq(rej, 1);
2425 cpl_test(isnan(f.data));
2426 cpl_test(isnan(f.error));
2427 }
2428
2429 double src_flx = 2.0;
2430 for(cpl_size i = 3; i < int_size - 3; ++i){
2431
2432 /*rejected because inside a bad pixel*/
2433 if((i >=4 && i <= 6) || (i >= 10 && i <= 12)){
2434 int rej = 0;
2435 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2436 cpl_test_eq(rej, 1);
2437 cpl_test(isnan(f.data));
2438 cpl_test(isnan(f.error));
2439
2440 if((i % 2) == 0) src_flx++;
2441
2442 continue;
2443 }
2444
2445 int rej = 0;
2446 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2447
2448 cpl_test_eq(rej, 0);
2449
2450 if((i % 2) == 1){
2451 /*odd elements: they are in the center of the source bin*/
2452 cpl_test_rel(f.data, src_flx, HDRL_DELTA_COMPARE_VALUE);
2453 cpl_test_rel(f.error, src_flx / 10.0, HDRL_DELTA_COMPARE_VALUE);
2454
2455 }
2456 else{
2457 /*even elements: they are between bins*/
2458 const double el = src_flx + .5;
2459 cpl_test_rel(f.data, el, HDRL_DELTA_COMPARE_VALUE);
2460
2461 const double el_e = pow(src_flx, 2.0) + pow((src_flx + 1.0), 2.0);
2462 cpl_test_rel(f.error, sqrt(el_e) / (10.0 * sqrt(2.0)), HDRL_DELTA_COMPARE_VALUE);
2463 src_flx++;
2464 }
2465 }
2466
2468 hdrl_spectrum1D_delete(&ori_s);
2469 hdrl_spectrum1D_delete(&integrated_s);
2470 cpl_array_delete(wavs_integrate);
2471 cpl_array_delete(bad_idxes);
2472}
2473
2474void test6(void){
2475
2477
2478 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 8, 20.0, 2.0);
2479
2480 cpl_array * wavs_integrate = get_waves(20.5, 15, 1);
2481 hdrl_spectrum1D * integrated_s =
2482 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2483
2484 double src_flx = 1.0;
2485
2486 for(cpl_size i = 0; i < cpl_array_get_size(wavs_integrate) - 1; ++i){
2487 int rej = 0;
2488 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2489 cpl_test_eq(rej, 0);
2490 cpl_test_rel(f.data, src_flx, HDRL_DELTA_COMPARE_VALUE);
2491 cpl_test_rel(f.error, src_flx / 10.0, HDRL_DELTA_COMPARE_VALUE);
2492 src_flx = src_flx + ((i+1) % 2);
2493 }
2494
2495 int rej = 0;
2496 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s,
2497 cpl_array_get_size(wavs_integrate) - 1, &rej);
2498 cpl_test_eq(rej, 1);
2499 cpl_test(isnan(f.data));
2500 cpl_test(isnan(f.error));
2501
2503 hdrl_spectrum1D_delete(&ori_s);
2504 hdrl_spectrum1D_delete(&integrated_s);
2505 cpl_array_delete(wavs_integrate);
2506}
2507
2508void test7(void){
2509
2511
2512 cpl_array * bads = cpl_array_new(3, CPL_TYPE_INT);
2513
2514 cpl_array_set(bads, 0, 0);
2515 cpl_array_set(bads, 1, 7);
2516 cpl_array_set(bads, 2, 4);
2517
2518 hdrl_spectrum1D * ori_s = generate_bad_stair_spectrum(1, 8, 20.0, 2.0, bads);
2519
2520 cpl_array * wavs_integrate = get_waves(20.5, 15, 1);
2521 hdrl_spectrum1D * integrated_s =
2522 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2523
2524 double src_flx = 1.0;
2525
2526 for(cpl_size i = 0; i < cpl_array_get_size(wavs_integrate); ++i){
2527 int rej = 0;
2528 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2529
2530 if(i == 0 || i == 7 || i == 8 || i >= 13){
2531 cpl_test_eq(rej, 1);
2532 cpl_test(isnan(f.data));
2533 cpl_test(isnan(f.error));
2534 }
2535 else{
2536 cpl_test_eq(rej, 0);
2537 cpl_test_rel(f.data, src_flx, HDRL_DELTA_COMPARE_VALUE);
2538 cpl_test_rel(f.error, src_flx / 10.0, HDRL_DELTA_COMPARE_VALUE);
2539 }
2540 src_flx = src_flx + ((i+1) % 2);
2541 }
2542
2544 hdrl_spectrum1D_delete(&ori_s);
2545 hdrl_spectrum1D_delete(&integrated_s);
2546 cpl_array_delete(wavs_integrate);
2547 cpl_array_delete(bads);
2548}
2549
2550
2551void test8(void){
2553
2554 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 8, 20.0, 1.0);
2555
2556 cpl_array * wavs_integrate = get_waves(21, 3, 2);
2557 hdrl_spectrum1D * integrated_s =
2558 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2559
2560 double res_flx[] = {2.5, 4, 5.5};
2561 double res_flx_e[] = {sqrt(6.5), sqrt(16.5), sqrt(30.5)};
2562
2563 for(cpl_size i = 0; i < 3; ++i){
2564
2565 int rej = 0;
2566 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2567
2568 cpl_test_eq(rej, 0);
2569 cpl_test_rel(f.data, res_flx[i], HDRL_DELTA_COMPARE_VALUE);
2570 cpl_test_rel(f.error, res_flx_e[i] / 10.0, HDRL_DELTA_COMPARE_VALUE);
2571 }
2572
2574 hdrl_spectrum1D_delete(&ori_s);
2575 hdrl_spectrum1D_delete(&integrated_s);
2576 cpl_array_delete(wavs_integrate);
2577}
2578
2579void test9(void){
2581
2582 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 8, 20.0, 1.0);
2583
2584 cpl_array * wavs_integrate = get_waves(19, 5, 2);
2585 hdrl_spectrum1D * integrated_s =
2586 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2587
2588 {
2589 int rej = 0;
2590 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, 0, &rej);
2591 cpl_test(rej);
2592 cpl_test(isnan(f.data));
2593 cpl_test(isnan(f.error));
2594 }
2595
2596 double res_flx[] = {2, 4, 6, 7.5};
2597 for(cpl_size i = 1; i <= 4; ++i){
2598
2599 int rej = 0;
2600 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2601
2602 cpl_test_eq(rej, 0);
2603 cpl_test_rel(f.data, res_flx[i - 1], HDRL_DELTA_COMPARE_VALUE);
2604 }
2605
2607 hdrl_spectrum1D_delete(&ori_s);
2608 hdrl_spectrum1D_delete(&integrated_s);
2609 cpl_array_delete(wavs_integrate);
2610}
2611
2612
2613void test10(void){
2615
2616 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 8, 20.0, 1.0);
2617
2618 cpl_array * wavs_integrate = get_waves(19, 6, 2);
2619 hdrl_spectrum1D * integrated_s =
2620 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2621 {
2622 int rej = 0;
2623 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, 0, &rej);
2624 cpl_test(rej);
2625 cpl_test(isnan(f.data));
2626 cpl_test(isnan(f.error));
2627 }
2628
2629 double res_flx[] = {2, 4, 6};
2630 for(cpl_size i = 1; i <= 3; ++i){
2631
2632 int rej = 0;
2633 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2634
2635 cpl_test_eq(rej, 0);
2636 cpl_test_rel(f.data, res_flx[i - 1], HDRL_DELTA_COMPARE_VALUE);
2637 }
2638
2639 for(cpl_size i = 4; i < 6; ++i){
2640 int rej = 0;
2641 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2642 cpl_test(rej);
2643 cpl_test(isnan(f.data));
2644 cpl_test(isnan(f.error));
2645 }
2646
2648 hdrl_spectrum1D_delete(&ori_s);
2649 hdrl_spectrum1D_delete(&integrated_s);
2650 cpl_array_delete(wavs_integrate);
2651}
2652
2653
2654void test11(void){
2656
2657 cpl_array * bads = cpl_array_new(1, CPL_TYPE_INT);
2658
2659 cpl_array_set(bads, 0, 3);
2660
2661 hdrl_spectrum1D * ori_s = generate_bad_stair_spectrum(1, 8, 20.0, 1.0, bads);
2662
2663 cpl_array * wavs_integrate = get_waves(19, 6, 2);
2664 hdrl_spectrum1D * integrated_s =
2665 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2666 {
2667 int rej = 0;
2668 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, 0, &rej);
2669 cpl_test(rej);
2670 cpl_test(isnan(f.data));
2671 cpl_test(isnan(f.error));
2672 }
2673
2674 double res_flx[] = {2, 4, 6};
2675 for(cpl_size i = 1; i <= 3; ++i){
2676
2677 int rej = 0;
2678 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2679
2680 if(i == 2){
2681 cpl_test(rej);
2682 cpl_test(isnan(f.data));
2683 cpl_test(isnan(f.error));
2684 continue;
2685 }
2686
2687 cpl_test_eq(rej, 0);
2688 cpl_test_rel(f.data, res_flx[i - 1], HDRL_DELTA_COMPARE_VALUE);
2689 }
2690
2691 for(cpl_size i = 4; i < 6; ++i){
2692 int rej = 0;
2693 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2694 cpl_test(rej);
2695 cpl_test(isnan(f.data));
2696 cpl_test(isnan(f.error));
2697 }
2698
2700 hdrl_spectrum1D_delete(&ori_s);
2701 hdrl_spectrum1D_delete(&integrated_s);
2702 cpl_array_delete(wavs_integrate);
2703 cpl_array_delete(bads);
2704}
2705
2706void test12(void){
2708
2709 cpl_array * bads = cpl_array_new(1, CPL_TYPE_INT);
2710
2711 cpl_array_set(bads, 0, 4);
2712
2713 hdrl_spectrum1D * ori_s = generate_bad_stair_spectrum(1, 8, 20.0, 1.0, bads);
2714
2715 cpl_array * wavs_integrate = get_waves(19, 6, 2);
2716 hdrl_spectrum1D * integrated_s =
2717 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2718
2719 {
2720 int rej = 0;
2721 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, 0, &rej);
2722 cpl_test(rej);
2723 cpl_test(isnan(f.data));
2724 cpl_test(isnan(f.error));
2725 }
2726
2727 double res_flx[] = {2, 4, 6};
2728 for(cpl_size i = 1; i <= 3; ++i){
2729
2730 int rej = 0;
2731 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2732
2733 if(i == 2 || i == 3){
2734 cpl_test(rej);
2735 cpl_test(isnan(f.data));
2736 cpl_test(isnan(f.error));
2737 continue;
2738 }
2739
2740 cpl_test_eq(rej, 0);
2741 cpl_test_rel(f.data, res_flx[i - 1], HDRL_DELTA_COMPARE_VALUE);
2742 }
2743
2744 for(cpl_size i = 4; i < 6; ++i){
2745 int rej = 0;
2746 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2747 cpl_test(rej);
2748 cpl_test(isnan(f.data));
2749 cpl_test(isnan(f.error));
2750 }
2751
2753 hdrl_spectrum1D_delete(&ori_s);
2754 hdrl_spectrum1D_delete(&integrated_s);
2755 cpl_array_delete(wavs_integrate);
2756 cpl_array_delete(bads);
2757}
2758
2759
2760void test13(void){
2762
2763 cpl_array * bads = cpl_array_new(2, CPL_TYPE_INT);
2764
2765 cpl_array_set(bads, 0, 0);
2766 cpl_array_set(bads, 1, 7);
2767
2768 hdrl_spectrum1D * ori_s = generate_bad_stair_spectrum(1, 8, 20.0, 1.0, bads);
2769
2770 cpl_array * wavs_integrate = get_waves(19, 6, 2);
2771 hdrl_spectrum1D * integrated_s =
2772 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2773 {
2774 int rej = 0;
2775 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, 0, &rej);
2776 cpl_test(rej);
2777 cpl_test(isnan(f.data));
2778 cpl_test(isnan(f.error));
2779 }
2780
2781 double res_flx[] = {2, 4, 6};
2782 for(cpl_size i = 0; i <= 3; ++i){
2783
2784 int rej = 0;
2785 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2786
2787 if(i == 0 || i == 1){
2788 cpl_test(rej);
2789 cpl_test(isnan(f.data));
2790 cpl_test(isnan(f.error));
2791 continue;
2792 }
2793
2794 cpl_test_eq(rej, 0);
2795 cpl_test_rel(f.data, res_flx[i - 1], HDRL_DELTA_COMPARE_VALUE);
2796 }
2797
2798 for(cpl_size i = 4; i < 6; ++i){
2799 int rej = 0;
2800 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2801 cpl_test(rej);
2802 cpl_test(isnan(f.data));
2803 cpl_test(isnan(f.error));
2804 }
2805
2807 hdrl_spectrum1D_delete(&ori_s);
2808 hdrl_spectrum1D_delete(&integrated_s);
2809 cpl_array_delete(wavs_integrate);
2810 cpl_array_delete(bads);
2811}
2812
2813void test14(void){
2814
2816
2817 cpl_array * bads = cpl_array_new(2, CPL_TYPE_INT);
2818
2819 cpl_array_set(bads, 0, 0);
2820 cpl_array_set(bads, 1, 4);
2821
2822 hdrl_spectrum1D * ori_s = generate_bad_stair_spectrum(1, 9, 20.5, 1.0, bads);
2823
2824 cpl_array * wavs_integrate = get_waves(21, 5, 2);
2825 hdrl_spectrum1D * integrated_s =
2826 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2827
2828
2829 double res_flx[] = {2, 3.5, 4.5, 7.5};
2830
2831 for(cpl_size i = 0; i < cpl_array_get_size(wavs_integrate) - 1; ++i){
2832 int rej = 0;
2833 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2834
2835 if(i == 2){
2836 cpl_test(rej);
2837 cpl_test(isnan(f.data));
2838 cpl_test(isnan(f.error));
2839 continue;
2840 }
2841
2842 cpl_test_eq(rej, 0);
2843 cpl_test_rel(f.data, res_flx[i], HDRL_DELTA_COMPARE_VALUE);
2844 }
2845
2846 int rej = 0;
2847 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s,
2848 cpl_array_get_size(wavs_integrate) - 1, &rej);
2849 cpl_test(rej);
2850 cpl_test(isnan(f.data));
2851 cpl_test(isnan(f.error));
2852
2854 hdrl_spectrum1D_delete(&ori_s);
2855 hdrl_spectrum1D_delete(&integrated_s);
2856 cpl_array_delete(wavs_integrate);
2857 cpl_array_delete(bads);
2858}
2859
2860void test15(void){
2861
2863
2864 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 3, 20.5, 1.0);
2865
2866 double wa_mx = hdrl_spectrum1D_get_wavelength_value(ori_s, 2, NULL);
2867 cpl_test_rel(wa_mx, 22.5, HDRL_DELTA_COMPARE_VALUE);
2868
2869 cpl_array * wavs_integrate = get_waves(22.5, 5, 2);
2870 hdrl_spectrum1D * integrated_s =
2871 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2872
2873 for(cpl_size i = 0; i < cpl_array_get_size(wavs_integrate); ++i){
2874 int rej = 0;
2875 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2876 cpl_test(rej);
2877 cpl_test(isnan(f.data));
2878 cpl_test(isnan(f.error));
2879 }
2880
2882 hdrl_spectrum1D_delete(&ori_s);
2883 hdrl_spectrum1D_delete(&integrated_s);
2884 cpl_array_delete(wavs_integrate);
2885}
2886
2887void test16(void){
2888
2890
2891 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 3, 20.5, 1.0);
2892
2893 cpl_array * wavs_integrate = get_waves(14.5, 4, 2);
2894
2895 double wa_mx = cpl_array_get(wavs_integrate, 3, NULL);
2896 cpl_test_rel(wa_mx, 20.5, HDRL_DELTA_COMPARE_VALUE);
2897
2898 hdrl_spectrum1D * integrated_s =
2899 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2900
2901 for(cpl_size i = 0; i < cpl_array_get_size(wavs_integrate); ++i){
2902 int rej = 0;
2903 const hdrl_value f = hdrl_spectrum1D_get_flux_value(integrated_s, i, &rej);
2904 cpl_test(rej);
2905 cpl_test(isnan(f.data));
2906 cpl_test(isnan(f.error));
2907 }
2908
2910 hdrl_spectrum1D_delete(&ori_s);
2911 hdrl_spectrum1D_delete(&integrated_s);
2912 cpl_array_delete(wavs_integrate);
2913}
2914
2915void test17(void){
2916
2918 cpl_array * bads = cpl_array_new(3, CPL_TYPE_INT);
2919
2920 cpl_array_set(bads, 0, 0);
2921 cpl_array_set(bads, 1, 4);
2922 cpl_array_set(bads, 2, 7);
2923
2924 hdrl_spectrum1D * ori_s = generate_bad_stair_spectrum(1, 8, 20.5, 1.0, bads);
2925
2926 const cpl_array * wavs_integrate =
2927 hdrl_spectrum1D_get_wavelength(ori_s).wavelength;
2928
2929 hdrl_spectrum1D * integrated_s =
2930 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2931
2932 for(cpl_size i = 0; i < cpl_array_get_size(wavs_integrate); ++i){
2933 int ori_rej = 0;
2934 const hdrl_value ori_flx = hdrl_spectrum1D_get_flux_value(ori_s, i, &ori_rej);
2935
2936 int int_rej = 0;
2937 const hdrl_value int_flx = hdrl_spectrum1D_get_flux_value(integrated_s, i, &int_rej);
2938
2939 cpl_test_eq(ori_rej, int_rej);
2940 if(ori_rej) continue;
2941
2942 cpl_test_rel(ori_flx.data, int_flx.data, HDRL_DELTA_COMPARE_VALUE);
2943 cpl_test_rel(ori_flx.error, int_flx.error, HDRL_DELTA_COMPARE_VALUE);
2944 }
2945
2947 hdrl_spectrum1D_delete(&ori_s);
2948 hdrl_spectrum1D_delete(&integrated_s);
2949 cpl_array_delete(bads);
2950}
2951
2952double get_bin_size(const hdrl_spectrum1D * s, const cpl_size i){
2953
2954 if(i == 0){
2955 const double w0 = hdrl_spectrum1D_get_wavelength_value(s, i, NULL);
2956 const double w1 = hdrl_spectrum1D_get_wavelength_value(s, i + 1, NULL);
2957 return (w0 + w1) / 2 - w0;
2958 }
2959
2960 const cpl_size sz = hdrl_spectrum1D_get_size(s);
2961 if(i == sz - 1){
2962 const double w0 = hdrl_spectrum1D_get_wavelength_value(s, i - 1, NULL);
2963 const double w1 = hdrl_spectrum1D_get_wavelength_value(s, i, NULL);
2964 return (w0 + w1) / 2 - w0;
2965 }
2966
2967 const double w0 = hdrl_spectrum1D_get_wavelength_value(s, i - 1, NULL);
2968 const double w1 = hdrl_spectrum1D_get_wavelength_value(s, i + 1, NULL);
2969 return (w0 + w1) / 2 - w0;
2970}
2971
2972double calc_total_flux(const hdrl_spectrum1D * s){
2973 double flux = 0;
2974 const cpl_size sz = hdrl_spectrum1D_get_size(s);
2975 for(cpl_size i = 0; i < sz; ++i)
2976 flux += hdrl_spectrum1D_get_flux_value(s, i, NULL).data * get_bin_size(s, i);
2977
2978 return flux;
2979}
2980
2981void test18(void){
2982
2983
2985
2986 hdrl_spectrum1D * ori_s = generate_stair_spectrum(1, 4, 20.0, 2.0);
2987
2988 cpl_array * wavs_integrate = get_waves(20, 3, 3);
2989
2990 hdrl_spectrum1D * integrated_s =
2991 hdrl_spectrum1D_resample_on_array(ori_s, wavs_integrate, par);
2992
2993 cpl_test_rel(calc_total_flux(ori_s), 15.0, HDRL_DELTA_COMPARE_VALUE);
2994
2995 cpl_test_rel(calc_total_flux(ori_s), calc_total_flux(integrated_s), HDRL_DELTA_COMPARE_VALUE);
2996
2998 hdrl_spectrum1D_delete(&ori_s);
2999 hdrl_spectrum1D_delete(&integrated_s);
3000 cpl_array_delete(wavs_integrate);
3001}
3002
3003/*----------------------------------------------------------------------------*/
3007/*----------------------------------------------------------------------------*/
3008void test_spectrum1D_resample_spectrum_integrate(void){
3009
3010 /*tests 1 to 6 test the case of upsampling doing integration*/
3011 /*test when the destination spectrum starts after the source and ends before
3012 * the source end.*/
3013 test1();
3014
3015 /*test when source and destination points start and stop at the same bin*/
3016 test2();
3017
3018 /*test when destination points cover a wider range of the source*/
3019 test3();
3020
3021 /*like test1 but the output and input values are shuffled*/
3022 test4();
3023
3024 /*like test3 but we have bad pixels inside (and the edges) of the spectrum*/
3025 test5();
3026
3027 /* Test for rebinning, the dest bins split exactly source bin in half,
3028 * except the last one which is outside*/
3029 test6();
3030
3031 /*Test like test6 but with bad pixels*/
3032 test7();
3033
3034 /*The following test test the case of dest bins bigger than source*/
3035
3036 /*Test when the destination spectrum starts after the source and ends before
3037 * the source end.*/
3038 test8();
3039
3040 /*test when source and destination points stop at the same bin*/
3041 test9();
3042
3043 /*test when destination points cover a wider range of the source*/
3044 test10();
3045
3046 /*like test 10 but we have a bad pixel*/
3047 test11();
3048 test12();
3049 test13();
3050
3051 /* Test for rebinning, the dest bins split exactly source bin in half,
3052 * except the last one which is outside, with bad pixels*/
3053 test14();
3054
3055 /*Source and destination have empty intersections*/
3056 test15();
3057 test16();
3058
3059 /*Source and destination with same bins*/
3060 test17();
3061
3062 /*Flux conservation*/
3063 test18();
3064
3065 cpl_test_eq(cpl_error_get_code(), CPL_ERROR_NONE);
3066}
3067
3068void test_parlist(void)
3069{
3070 /* parameter parsing smoketest */
3071 const char *base_context = "RECIPE";
3072 const char *prefix = "test";
3073 const char *method_def = "LINEAR";
3074
3075 cpl_parameterlist *pos;
3076 pos = hdrl_spectrum1D_resample_interpolate_parameter_create_parlist(base_context, prefix, method_def);
3077 cpl_test_error(CPL_ERROR_NONE);
3078
3079 hdrl_parameter *hpar;
3080 hpar = hdrl_spectrum1D_resample_interpolate_parameter_parse_parlist(
3081 pos, "RECIPE.test");
3082 cpl_test_error(CPL_ERROR_NONE);
3083
3085
3086 cpl_parameterlist_delete(pos);
3087}
3088
3089/*----------------------------------------------------------------------------*/
3093/*----------------------------------------------------------------------------*/
3094int main(void)
3095{
3096 cpl_test_init(PACKAGE_BUGREPORT, CPL_MSG_WARNING);
3097
3098 srand(500);
3099
3100 test_parlist();
3101
3102 test_spectrum1D_constructor(CPL_TYPE_DOUBLE);
3103 test_spectrum1D_constructor(CPL_TYPE_FLOAT);
3104 test_spectrum1D_constructor_error();
3105 test_spectrum1D_constructor_analytical();
3106 test_spectrum1D_duplication();
3107
3108
3109 test_spectrum1D_calculation_scalar();
3110
3111 test_spectrum1D_calculation();
3112 test_spectrum1D_calculation_error();
3113
3114 test_spectrum1D_conversion_wavelength_scale();
3115 test_spectrum1D_mul_wavelength();
3116 test_spectrum1D_shift_wavelength();
3117
3118 test_spectrum1D_wavelength_select();
3119
3120 test_spectrum1D_resample_spectrum(CPL_TRUE);
3121 test_spectrum1D_resample_spectrum(CPL_FALSE);
3122 test_spectrum1D_resample_spectrum_private_funcs();
3123 test_spectrum1D_resample_spectrum_bpm(CPL_TRUE);
3124 test_spectrum1D_resample_spectrum_bpm(CPL_FALSE);
3125 test_spectrum1D_resample_spectrum_interpolation_error_test();
3126 test_spectrum1D_resample_spectrum_fit_error_test_shift(CPL_TRUE);
3127 test_spectrum1D_resample_spectrum_fit_error_test_shift(CPL_FALSE);
3128 test_spectrum1D_resample_spectrum_fit_error_test_error_interpol();
3129 test_spectrum1D_resample_spectrum_fit_windowed();
3130
3131 test_spectrum1D_resample_spectrum_integrate();
3132
3133 test_spectrum1D_table_conversion();
3134
3135 test_spectrum1D_test_uniformly_sampled();
3136
3137 test_spectrum1Dlist();
3138
3139 return cpl_test_end(0);
3140}
3141
3142/*----------------------------------------------------------------------------*/
3146/*----------------------------------------------------------------------------*/
3147
3148static inline cpl_image * get_random_1d_img(cpl_size length,
3149 double min, double max, cpl_type type){
3150 cpl_image * to_ret = cpl_image_new(length, 1, type);
3151
3152 for(cpl_size i = 0; i < length; i++){
3153 double d = rand_0_to_1() * (max - min) + min;
3154 cpl_image_set(to_ret, i + 1, 1, d);
3155 }
3156 return to_ret;
3157}
3158
3159static inline void set_1d_bpm(cpl_image * img){
3160
3161 int sz_x = cpl_image_get_size_x(img);
3162 int sz_y = cpl_image_get_size_y(img);
3163
3164 cpl_mask * msk = cpl_image_get_bpm(img);
3165
3166 for(cpl_size x = 0; x < sz_x; x++){
3167 for(cpl_size y = 0; y < sz_y; y++){
3168 cpl_mask_set(msk, x + 1, y + 1, rand_0_to_1() > .5);
3169 }
3170 }
3171
3172}
3173
3174static inline cpl_array * get_wavelength(cpl_size length, cpl_type type){
3175
3176 cpl_array * to_ret = cpl_array_new(length, type);
3177 double d = rand_0_to_1();
3178
3179 for(cpl_size i = 0; i < length; i++){
3180 cpl_array_set(to_ret, i, d);
3181 d += 1.0 + rand_0_to_1();
3182 }
3183 return to_ret;
3184}
3185
3186static inline cpl_boolean are_cpl_img_eq(
3187 const cpl_image *im1, const cpl_image *im2){
3188
3189 if (cpl_image_get_size_x(im1) != cpl_image_get_size_x(im2) ||
3190 cpl_image_get_size_y(im1) != cpl_image_get_size_y(im2) ){
3191 return CPL_FALSE;
3192 }
3193
3194 cpl_size sz_x = cpl_image_get_size_x(im1);
3195 cpl_size sz_y = cpl_image_get_size_y(im1);
3196
3197 int rej1 = 0;
3198 int rej2 = 0;
3199
3200 for(cpl_size x = 0; x < sz_x; x++){
3201 for(cpl_size y = 0; y < sz_y; y++){
3202 double px1 = cpl_image_get(im1, x + 1, y + 1, &rej1);
3203 double px2 = cpl_image_get(im2, x + 1, y + 1, &rej2);
3204
3205 if(px1 != px2 || rej1 != rej2) {
3206 return CPL_FALSE;
3207 }
3208 }
3209 }
3210
3211 return CPL_TRUE;
3212}
3213
3214static inline cpl_boolean are_hdrl_eq(const hdrl_image* flux_compound,
3215 const cpl_image * flux, const cpl_image * flux_e){
3216
3217 if (cpl_image_get_size_x(flux) != cpl_image_get_size_x(flux_e) ||
3218 cpl_image_get_size_y(flux) != cpl_image_get_size_y(flux_e) ||
3219 hdrl_image_get_size_x(flux_compound) != cpl_image_get_size_x(flux) ||
3220 hdrl_image_get_size_y(flux_compound) != cpl_image_get_size_y(flux) ){
3221 return CPL_FALSE;
3222 }
3223
3224 const cpl_image * flux_hdrl = hdrl_image_get_image_const(flux_compound);
3225 const cpl_image * flux_e_hdrl = hdrl_image_get_error_const(flux_compound);
3226
3227 hdrl_image * hdrl_img = hdrl_image_create(flux, flux_e);
3228
3229 cpl_boolean is_success =
3230 are_cpl_img_eq(hdrl_image_get_image_const(hdrl_img), flux_hdrl);
3231
3232 is_success &= are_cpl_img_eq(hdrl_image_get_error_const(hdrl_img),
3233 flux_e_hdrl);
3234
3235 hdrl_image_delete(hdrl_img); hdrl_img = NULL;
3236
3237 return is_success;
3238}
3239
3240static inline cpl_error_code get_error_code_and_reset(void){
3241
3242 cpl_error_code err = cpl_error_get_code();
3243 cpl_error_reset();
3244 return err;
3245}
3246
3247static inline double rand_0_to_1(void){
3248 double r = rand();
3249 r /= RAND_MAX;
3250 return r;
3251}
3252
3253static inline hdrl_spectrum1D *
3254get_random_spectrum(int length, hdrl_spectrum1D_wave_scale scale){
3255
3256 cpl_image * spectrum1d =
3257 get_random_1d_img(length, 1.0f, 128.0f, CPL_TYPE_DOUBLE);
3258 cpl_image * spectrum1d_error =
3259 get_random_1d_img(length, 0.5f, 2.0f, CPL_TYPE_DOUBLE);
3260
3261 set_1d_bpm(spectrum1d);
3262 cpl_array * wavelengths = get_wavelength(length, CPL_TYPE_DOUBLE);
3263
3264 hdrl_spectrum1D * s1 = hdrl_spectrum1D_create(
3265 spectrum1d, spectrum1d_error,
3266 wavelengths, scale);
3267
3268 cpl_array_delete(wavelengths);
3269 cpl_image_delete(spectrum1d);
3270 cpl_image_delete(spectrum1d_error);
3271
3272 return s1;
3273}
3274
3275
3276static inline void
3277test_error_create_func(const hdrl_spectrum1D * s1, const hdrl_spectrum1D * s2,
3278 operate_spectra_create f){
3279
3280 hdrl_spectrum1D * res = f(s1, s2);
3281 cpl_test_null(res);
3282 const cpl_error_code cd = get_error_code_and_reset();
3283 cpl_test_noneq(cd, CPL_ERROR_NONE);
3284}
3285
3286static inline void
3287test_error_mutate_func(hdrl_spectrum1D * s1, const hdrl_spectrum1D * s2,
3288 operate_spectra f){
3289
3290 const cpl_error_code res = f(s1, s2);
3291 cpl_test_noneq(res, CPL_ERROR_NONE);
3292
3293 const cpl_error_code cd = get_error_code_and_reset();
3294 cpl_test_noneq(cd, CPL_ERROR_NONE);
3295}
3296
3297static inline void
3298test_calc_creat_error(operate_spectra_create f){
3299
3300 hdrl_spectrum1D * spec_l40_linear =
3301 get_random_spectrum(40, hdrl_spectrum1D_wave_scale_linear);
3302 hdrl_spectrum1D * spec_l40_log =
3303 get_random_spectrum(40, hdrl_spectrum1D_wave_scale_log);
3304
3305 hdrl_spectrum1D * spec_l41_linear =
3306 get_random_spectrum(41, hdrl_spectrum1D_wave_scale_linear);
3307 hdrl_spectrum1D * spec_l41_log =
3308 get_random_spectrum(41, hdrl_spectrum1D_wave_scale_log);
3309
3310 test_error_create_func
3311 (spec_l40_linear, spec_l40_log, f);
3312 test_error_create_func
3313 (spec_l41_linear, spec_l40_linear, f);
3314 test_error_create_func
3315 (spec_l40_log, spec_l40_linear, f);
3316 test_error_create_func
3317 (spec_l40_linear, spec_l41_linear, f);
3318
3319 test_error_create_func
3320 (NULL, spec_l40_log, f);
3321 test_error_create_func
3322 (spec_l41_linear, NULL, f);
3323 test_error_create_func(NULL, NULL, f);
3324
3325
3326 hdrl_spectrum1D_delete(&spec_l40_linear);
3327 hdrl_spectrum1D_delete(&spec_l41_linear);
3328
3329 hdrl_spectrum1D_delete(&spec_l40_log);
3330 hdrl_spectrum1D_delete(&spec_l41_log);
3331}
3332
3333
3334static inline void test_calc_error(operate_spectra f){
3335
3336 hdrl_spectrum1D * spec_l40_linear =
3337 get_random_spectrum(40, hdrl_spectrum1D_wave_scale_linear);
3338 hdrl_spectrum1D * spec_l40_log =
3339 get_random_spectrum(40, hdrl_spectrum1D_wave_scale_log);
3340
3341 hdrl_spectrum1D * spec_l41_linear =
3342 get_random_spectrum(41, hdrl_spectrum1D_wave_scale_linear);
3343 hdrl_spectrum1D * spec_l41_log =
3344 get_random_spectrum(41, hdrl_spectrum1D_wave_scale_log);
3345
3346 test_error_mutate_func
3347 (spec_l40_linear, spec_l40_log, f);
3348 test_error_mutate_func
3349 (spec_l41_linear, spec_l40_linear, f);
3350 test_error_mutate_func
3351 (spec_l40_log, spec_l40_linear, f);
3352 test_error_mutate_func
3353 (spec_l40_linear, spec_l41_linear, f);
3354
3355 test_error_mutate_func
3356 (NULL, spec_l40_log, f);
3357 test_error_mutate_func
3358 (spec_l41_linear, NULL, f);
3359 test_error_mutate_func(NULL, NULL, f);
3360
3361 hdrl_spectrum1D_delete(&spec_l40_linear);
3362 hdrl_spectrum1D_delete(&spec_l41_linear);
3363
3364 hdrl_spectrum1D_delete(&spec_l40_log);
3365 hdrl_spectrum1D_delete(&spec_l41_log);
3366}
3367
3368static inline hdrl_spectrum1D * get_spectrum1D_sin_shuffled(cpl_size sz, int start,
3369 cpl_boolean add_peak, cpl_array ** unshuffled_lambda){
3370
3371 static const double peak = 100.0;
3372 const double delta = 2 * CPL_MATH_PI / sz;
3373
3374 cpl_array * lambda = cpl_array_new(sz, HDRL_TYPE_DATA);
3375 cpl_image * flux = cpl_image_new(sz, 1, HDRL_TYPE_DATA);
3376
3377 for(cpl_size i = 0; i < sz; i++){
3378 double l = delta * (i + start);
3379 double f = fabs(peak * (sin(l) + 1.1));
3380
3381 if(i == 4 && add_peak) f *= 1.5;
3382
3383 cpl_array_set(lambda, i, l);
3384 cpl_image_set(flux, i + 1, 1, f);
3385 }
3386
3387 if(unshuffled_lambda)
3388 *unshuffled_lambda = cpl_array_duplicate(lambda);
3389
3390 /* scramble array */
3391 for(cpl_size i1 = 0; i1 < sz; i1++){
3392 int rej;
3393 double l1 = cpl_array_get(lambda, i1, &rej);
3394 double f1 = cpl_image_get(flux, i1 + 1, 1, &rej);
3395
3396 cpl_size i2 = (cpl_size)(rand_0_to_1() * (sz - 1));
3397
3398 double l2 = cpl_array_get(lambda, i2, &rej);
3399 double f2 = cpl_image_get(flux, i2 + 1, 1, &rej);
3400
3401 cpl_array_set(lambda, i1, l2);
3402 cpl_image_set(flux, i1 + 1, 1, f2);
3403
3404 cpl_array_set(lambda, i2, l1);
3405 cpl_image_set(flux, i2 + 1, 1, f1);
3406 }
3407
3408
3409 hdrl_spectrum1D * sp =
3411 (flux, 10, lambda, hdrl_spectrum1D_wave_scale_linear);
3412
3413 cpl_test_nonnull(sp);
3414 cpl_test_eq(get_error_code_and_reset(), CPL_ERROR_NONE);
3415
3416
3417 cpl_array_delete(lambda);
3418 cpl_image_delete(flux);
3419
3420 return sp;
3421}
cpl_size hdrl_image_get_size_y(const hdrl_image *self)
return size of Y dimension of image
Definition: hdrl_image.c:540
cpl_size hdrl_image_get_size_x(const hdrl_image *self)
return size of X dimension of image
Definition: hdrl_image.c:525
const cpl_image * hdrl_image_get_error_const(const hdrl_image *himg)
get error as cpl image
Definition: hdrl_image.c:144
hdrl_image * hdrl_image_create(const cpl_image *image, const cpl_image *error)
create a new hdrl_image from to existing images by copying them
Definition: hdrl_image.c:295
const cpl_image * hdrl_image_get_image_const(const hdrl_image *himg)
get data as cpl image
Definition: hdrl_image.c:118
void hdrl_image_delete(hdrl_image *himg)
delete hdrl_image
Definition: hdrl_image.c:379
void hdrl_parameter_delete(hdrl_parameter *obj)
shallow delete of a parameter
hdrl_spectrum1D * hdrl_spectrum1D_pow_scalar_create(const hdrl_spectrum1D *self, hdrl_value scalar_operator)
subtract a scalar from a spectrum
const hdrl_image * hdrl_spectrum1D_get_flux(const hdrl_spectrum1D *self)
hdrl_spectrum1D getter flux
cpl_error_code hdrl_spectrum1D_mul_scalar(hdrl_spectrum1D *self, hdrl_value scalar_operator)
computes the elementwise multiplication of a spectrum by a scalar, the self parameter is modified
cpl_error_code hdrl_spectrum1D_wavelength_mult_scalar_linear(hdrl_spectrum1D *self, hdrl_data_t scale_linear)
computes the elementwise multiplication of the scalar for the wavelength. The scalar is assumed to be...
cpl_error_code hdrl_spectrum1D_wavelength_convert_to_linear(hdrl_spectrum1D *self)
converts the wavelength scale to linear.
hdrl_spectrum1D * hdrl_spectrum1D_sub_spectrum_create(const hdrl_spectrum1D *f1, const hdrl_spectrum1D *f2)
subtract two spectra
cpl_error_code hdrl_spectrum1D_wavelength_shift(hdrl_spectrum1D *self, hdrl_data_t shift)
computes the elementwise shift of the wavelength by the shift parameter. The self parameter is modifi...
cpl_size hdrl_spectrum1Dlist_get_size(const hdrl_spectrum1Dlist *l)
hdrl_spectrum1Dlist getter for size
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
cpl_error_code hdrl_spectrum1D_div_spectrum(hdrl_spectrum1D *self, const hdrl_spectrum1D *other)
divide one spectrum by another spectrum
hdrl_spectrum1Dlist * hdrl_spectrum1Dlist_new(void)
hdrl_spectrum1Dlist default constructor
cpl_error_code hdrl_spectrum1D_pow_scalar(hdrl_spectrum1D *self, hdrl_value scalar_operator)
computes the elementwise power of the flux to the scalar, the self parameter is modified
hdrl_spectrum1D * hdrl_spectrum1D_mul_spectrum_create(const hdrl_spectrum1D *f1, const hdrl_spectrum1D *f2)
multiply one spectrum by another spectrum
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_spectrum1D_mul_spectrum(hdrl_spectrum1D *self, const hdrl_spectrum1D *other)
multiply one spectrum by another spectrum
cpl_error_code hdrl_spectrum1D_wavelength_convert_to_log(hdrl_spectrum1D *self)
converts the wavelength scale to log. If the spectrum is already in log scale nothing is done.
hdrl_spectrum1D * hdrl_spectrum1D_div_spectrum_create(const hdrl_spectrum1D *num, const hdrl_spectrum1D *den)
divide one spectrum by another spectrum
cpl_error_code hdrl_spectrum1D_exp_scalar(hdrl_spectrum1D *self, hdrl_value scalar_operator)
computes the elementwise power of the scalar to the flux, the self parameter is modified
hdrl_spectrum1D * hdrl_spectrum1D_mul_scalar_create(const hdrl_spectrum1D *self, hdrl_value scalar_operator)
multiply a spectrum by a scalar
cpl_error_code hdrl_spectrum1D_sub_spectrum(hdrl_spectrum1D *self, const hdrl_spectrum1D *other)
subtract two spectra
cpl_error_code hdrl_spectrum1Dlist_set(hdrl_spectrum1Dlist *self, hdrl_spectrum1D *s, const cpl_size idx)
hdrl_spectrum1Dlist setter of the i-th element
hdrl_parameter * hdrl_spectrum1D_resample_fit_parameter_create(const int k, const int nCoeff)
constructor for the hdrl_parameter in the case of interpolation
void hdrl_spectrum1Dlist_delete(hdrl_spectrum1Dlist *l)
hdrl_spectrum1Dlist destructor
hdrl_spectrum1D * hdrl_spectrum1D_wavelength_shift_create(const hdrl_spectrum1D *self, hdrl_data_t shift)
computes the elementwise shift of the wavelength by the shift parameter.
hdrl_spectrum1D * hdrl_spectrum1D_div_scalar_create(const hdrl_spectrum1D *self, hdrl_value scalar_operator)
divide a spectrum by a scalar
hdrl_spectrum1D * hdrl_spectrum1D_create_error_DER_SNR(const cpl_image *arg_flux, cpl_size half_window, const cpl_array *wavelength, hdrl_spectrum1D_wave_scale scale)
hdrl_spectrum1D constructor when no error information is available, in this case we use DER_SNR to es...
hdrl_spectrum1D * hdrl_spectrum1Dlist_get(hdrl_spectrum1Dlist *self, const cpl_size idx)
hdrl_spectrum1Dlist getter of the i-th element
hdrl_spectrum1D * hdrl_spectrum1D_convert_from_table(const cpl_table *self, const char *flux_col_name, const char *wavelength_col_name, const char *flux_e_col_name, const char *flux_bpm_col_name, hdrl_spectrum1D_wave_scale scale)
convert a table to a spectrum
hdrl_spectrum1D * hdrl_spectrum1D_wavelength_convert_to_log_create(const hdrl_spectrum1D *self)
converts the wavelength scale to log. It returns a modified version of self. self is not modified....
const hdrl_spectrum1D * hdrl_spectrum1Dlist_get_const(const hdrl_spectrum1Dlist *self, const cpl_size idx)
hdrl_spectrum1Dlist getter of the i-th element
cpl_error_code hdrl_spectrum1D_sub_scalar(hdrl_spectrum1D *self, hdrl_value scalar_operator)
computes the elementwise subtraction of a spectrum by a scalar, the self parameter is modified
cpl_boolean hdrl_spectrum1D_is_uniformly_sampled(const hdrl_spectrum1D *self, double *bin)
checks if the spectrum is defined on uniformly sampled 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
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 * 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
cpl_error_code hdrl_spectrum1D_add_spectrum(hdrl_spectrum1D *self, const hdrl_spectrum1D *other)
sum two spectra
hdrl_spectrum1D * hdrl_spectrum1D_exp_scalar_create(const hdrl_spectrum1D *self, hdrl_value scalar_operator)
subtract a scalar from a spectrum
hdrl_spectrum1D * hdrl_spectrum1D_sub_scalar_create(const hdrl_spectrum1D *self, hdrl_value scalar_operator)
subtract a scalar from a spectrum
hdrl_spectrum1D * hdrl_spectrum1D_wavelength_mult_scalar_linear_create(const hdrl_spectrum1D *self, hdrl_data_t scale_linear)
computes the elementwise multiplication of the scalar for the wavelength. The scalar is assumed to be...
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_spectrum1D * hdrl_spectrum1D_add_scalar_create(const hdrl_spectrum1D *self, hdrl_value scalar_operator)
add a scalar to a spectrum
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_spectrum1D * hdrl_spectrum1D_create_analytic(calculate_analytic_spectrum_point func, const cpl_array *wavelength, hdrl_spectrum1D_wave_scale scale)
hdrl_spectrum1D constructor in the case of a spectrum defined by an analytical function
cpl_error_code hdrl_spectrum1D_append_to_table(const hdrl_spectrum1D *self, cpl_table *dest, const char *flux_col_name, const char *wavelength_col_name, const char *flux_e_col_name, const char *flux_bpm_col_name)
append a spectrum to a table.
hdrl_spectrum1D * hdrl_spectrum1D_wavelength_convert_to_linear_create(const hdrl_spectrum1D *self)
converts the wavelength scale to linear.
cpl_table * hdrl_spectrum1D_convert_to_table(const hdrl_spectrum1D *self, const char *flux_col_name, const char *wavelength_col_name, const char *flux_e_col_name, const char *flux_bpm_col_name)
converts a spectrum in a table.
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...
hdrl_spectrum1D * hdrl_spectrum1D_add_spectrum_create(const hdrl_spectrum1D *f1, const hdrl_spectrum1D *f2)
sum one spectrum to another spectrum
cpl_error_code hdrl_spectrum1D_div_scalar(hdrl_spectrum1D *self, hdrl_value scalar_operator)
computes the elementwise division of a spectrum by a scalar, the self parameter is modified
cpl_image * estimate_noise_DER_SNR(const hdrl_data_t *flux_in, const cpl_binary *msk_in, const cpl_array *wavelengths, const cpl_size length, const cpl_size half_window)
For every pixel in position i in img_arg, the function estimates the noise using the pixels in the wi...
Definition: hdrl_DER_SNR.c:160
hdrl_parameter * hdrl_spectrum1D_resample_integrate_parameter_create(void)
constructor for the hdrl_parameter in the case of integration
hdrl_spectrum1D * hdrl_spectrum1D_reject_pixels(const hdrl_spectrum1D *self, const cpl_array *bad_samples)
For every i-th element in bad_samples having value CPL_TRUE, the i-th pixel in the 1D spectrum is mar...
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