CR2RE Pipeline Reference Manual 1.6.2
hdrl_collapse-test.c
1/*
2 * This file is part of the HDRL
3 * Copyright (C) 2013,2014 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_combine.h"
29#include "hdrl_types.h"
30
31#include <cpl.h>
32
33#include <math.h>
34#include <stdlib.h>
35
36#ifndef ARRAY_LEN
37#define ARRAY_LEN(a) sizeof((a))/sizeof((a)[0])
38#endif
39
40
41/*----------------------------------------------------------------------------*/
45/*----------------------------------------------------------------------------*/
46typedef struct {
47 double v1;
48 double v2;
49} pair;
50
51/* cpl_test_abs that handles nan, rej is the rej variable of cpl_image_get to
52 * test if NAN is marked as a bad pixel or -1 if not applicable */
53#define hdrl_test_abs(a, b, tol, rej) \
54 if (isnan(a)) { \
55 cpl_test_eq(isnan(a), isnan(b)); \
56 cpl_test(rej == -1 || rej != 0); \
57 } \
58 else { \
59 cpl_test_abs(a, b, tol); \
60 }
61
62
63/* cpl_image_get that really returns the value regardless of bad or not */
64static double
65_hdrl_image_get(cpl_image * img, cpl_size x, cpl_size y, int * rej)
66{
67 double v = cpl_image_get(img, x, y, rej);
68 if (*rej) {
69 size_t nx = cpl_image_get_size_x(img);
70 if (cpl_image_get_type(img) == CPL_TYPE_DOUBLE) {
71 v = cpl_image_get_data_double(img)[(y - 1) * nx + (x - 1)];
72 }
73 else {
74 v = cpl_image_get_data_float(img)[(y - 1) * nx + (x - 1)];
75 }
76 }
77 return v;
78}
79
80/* ---------------------------------------------------------------------------*/
91/* ---------------------------------------------------------------------------*/
92static void
93prep_l2v_input(cpl_imagelist * data,
94 cpl_imagelist * errs,
95 cpl_size x,
96 cpl_size y,
97 cpl_imagelist * vl,
98 cpl_imagelist * el)
99{
100 int d;
101 size_t n = cpl_imagelist_get_size(data);
102 cpl_image * vimg = cpl_image_new(1, n, HDRL_TYPE_DATA);
103 cpl_image * verr = cpl_image_new(1, n, HDRL_TYPE_ERROR);
104 for (size_t i = 0; i < n; i++) {
105 cpl_image * img = cpl_imagelist_get(data, i);
106 cpl_image_set(vimg, 1, i + 1, cpl_image_get(img, x, y, &d));
107 if (d) {
108 cpl_image_reject(vimg, 1, i + 1);
109 }
110 img = cpl_imagelist_get(errs, i);
111 cpl_image_set(verr, 1, i + 1, cpl_image_get(img, x, y, &d));
112 if (d) {
113 cpl_image_reject(verr, 1, i + 1);
114 }
115 }
116 cpl_imagelist_set(vl, vimg, 0);
117 cpl_imagelist_set(el, verr, 0);
118}
119
120
121/* ---------------------------------------------------------------------------*/
135/* ---------------------------------------------------------------------------*/
136void test_l2i_and_l2v(cpl_imagelist * data,
137 cpl_imagelist * errs,
138 hdrl_collapse_imagelist_to_image_t * l2im,
139 hdrl_collapse_imagelist_to_vector_t * l2vm,
140 cpl_size x,
141 cpl_size y,
142 pair v,
143 pair e,
144 cpl_image * expcontrib)
145{
146 int d;
147 /* test list to image */
148 cpl_image * outimg, *outerr, *contrib;
149 hdrl_imagelist_combine(data, errs, l2im, &outimg, &outerr, &contrib);
150
151 hdrl_test_abs(_hdrl_image_get(outimg, x, y, &d), v.v1, v.v2, d);
152 hdrl_test_abs(_hdrl_image_get(outerr, x, y, &d), e.v1, e.v2, d);
153 cpl_test_image_abs(contrib, expcontrib, 0);
154 cpl_image_delete(outimg);
155 cpl_image_delete(outerr);
156 cpl_image_delete(contrib);
157
158 /* map the list to image input into a list to vector input that should give
159 * the same result */
160 cpl_imagelist * vl = cpl_imagelist_new();
161 cpl_imagelist * el = cpl_imagelist_new();
162 prep_l2v_input(data, errs, x, y, vl, el);
163
164 /* test list to vector */
165 cpl_vector * voutimg, *vouterr;
166 cpl_array * acontrib;
167 hdrl_collapse_imagelist_to_vector_call(l2vm, vl, el,
168 &voutimg, &vouterr, &acontrib,
169 NULL);
170 cpl_imagelist_delete(vl);
171 cpl_imagelist_delete(el);
172
173 hdrl_test_abs(cpl_vector_get(voutimg, 0), v.v1, v.v2, -1);
174 hdrl_test_abs(cpl_vector_get(vouterr, 0), e.v1, e.v2, -1);
175
176 cpl_test_abs(cpl_array_get_int(acontrib, 0, NULL),
177 cpl_image_get(expcontrib, x, y, &d), 0);
178
179 cpl_vector_delete(voutimg);
180 cpl_vector_delete(vouterr);
181 cpl_array_delete(acontrib);
182 hdrl_collapse_imagelist_to_image_delete(l2im);
183 hdrl_collapse_imagelist_to_vector_delete(l2vm);
184}
185
186void test_parameters(void)
187{
188 hdrl_parameter * hpar = hdrl_collapse_mean_parameter_create();
189 cpl_test(hdrl_collapse_parameter_is_mean(hpar));
192 cpl_test(hdrl_collapse_parameter_is_median(hpar));
201 cpl_test(hdrl_collapse_parameter_is_minmax(hpar));
203 hpar = hdrl_collapse_mode_parameter_create(100., 200., 1., HDRL_MODE_FIT, 100);
204 cpl_test(hdrl_collapse_parameter_is_mode(hpar));
206
207 cpl_test(hdrl_collapse_parameter_is_mean(HDRL_COLLAPSE_MEAN));
208 cpl_test(!hdrl_collapse_parameter_is_mean(HDRL_COLLAPSE_MEDIAN));
209 hdrl_parameter_delete((hdrl_parameter*)HDRL_COLLAPSE_MEAN);
210
211 cpl_test(hdrl_collapse_parameter_is_median(HDRL_COLLAPSE_MEDIAN));
212 hdrl_parameter_delete((hdrl_parameter*)HDRL_COLLAPSE_MEDIAN);
213
215 HDRL_COLLAPSE_WEIGHTED_MEAN));
216 hdrl_parameter_delete((hdrl_parameter*)HDRL_COLLAPSE_WEIGHTED_MEAN);
217}
218
219
220void test_parlist(void)
221{
222 /* parameter parsing smoketest */
223 hdrl_parameter * hpar;
224 hdrl_parameter * sigclip_def =
226 hdrl_parameter * minmax_def =
228 hdrl_parameter * mode_def =
229 hdrl_collapse_mode_parameter_create(10., 1., 0., HDRL_MODE_MEDIAN, 0);
230
231 cpl_test_error(CPL_ERROR_NONE);
232 cpl_parameterlist *parlist = hdrl_collapse_parameter_create_parlist(
233 "RECIPE", "collapse", "UNKNOWN", sigclip_def, minmax_def, mode_def) ;
234 cpl_test_error(CPL_ERROR_NONE);
235
236 cpl_parameterlist *parlist1 = hdrl_collapse_parameter_create_parlist(
237 "RECIPE", "collapse", "MEAN", sigclip_def, minmax_def, mode_def) ;
238 cpl_test_error(CPL_ERROR_NONE);
239
240 cpl_parameterlist *parlist2 = hdrl_collapse_parameter_create_parlist(
241 "RECIPE", "collapse", "WEIGHTED_MEAN", sigclip_def, minmax_def, mode_def) ;
242 cpl_test_error(CPL_ERROR_NONE);
243
244 cpl_parameterlist *parlist3 = hdrl_collapse_parameter_create_parlist(
245 "RECIPE", "collapse", "MEDIAN", sigclip_def, minmax_def, mode_def) ;
246 cpl_test_error(CPL_ERROR_NONE);
247
248 cpl_parameterlist *parlist4 = hdrl_collapse_parameter_create_parlist(
249 "RECIPE", "collapse", "SIGCLIP", sigclip_def, minmax_def, mode_def) ;
250 cpl_test_error(CPL_ERROR_NONE);
251
252 cpl_parameterlist *parlist5 = hdrl_collapse_parameter_create_parlist(
253 "RECIPE", "collapse", "MINMAX", sigclip_def, minmax_def, mode_def) ;
254 cpl_test_error(CPL_ERROR_NONE);
255
256 cpl_parameterlist *parlist6 = hdrl_collapse_parameter_create_parlist(
257 "RECIPE", "collapse", "MODE", sigclip_def, minmax_def, mode_def) ;
258 cpl_test_error(CPL_ERROR_NONE);
259
260
261/*
262 Befor adding the mode method:
263 cpl_test_eq(cpl_parameterlist_get_size(parlist4), 6);
264*/
265 cpl_test_eq(cpl_parameterlist_get_size(parlist5), 6 + 5);
266
267 hdrl_parameter_delete(sigclip_def);
268 hdrl_parameter_delete(minmax_def);
269 hdrl_parameter_delete(mode_def);
270
271
272 /* UNKNOWN */
273/* hpar = hdrl_collapse_parameter_parse_parlist(parlist, "RECIPE.invalid");
274 cpl_test_null(hpar);
275 cpl_test_error(CPL_ERROR_DATA_NOT_FOUND);
276 hpar = hdrl_collapse_parameter_parse_parlist(parlist, "RECIPE.collapse");
277 cpl_test_error(CPL_ERROR_NONE);
278 hdrl_parameter_delete(hpar);
279*/
280 /* MEAN */
281 hpar = hdrl_collapse_parameter_parse_parlist(parlist1, "RECIPE.invalid");
282 cpl_test_null(hpar);
283 cpl_test_error(CPL_ERROR_DATA_NOT_FOUND);
284 hpar = hdrl_collapse_parameter_parse_parlist(parlist1, "RECIPE.collapse");
285 cpl_test_error(CPL_ERROR_NONE);
287
288
289 /* WEIGHTED_MEAN */
290 hpar = hdrl_collapse_parameter_parse_parlist(parlist2, "RECIPE.invalid");
291 cpl_test_null(hpar);
292 cpl_test_error(CPL_ERROR_DATA_NOT_FOUND);
293 hpar = hdrl_collapse_parameter_parse_parlist(parlist2, "RECIPE.collapse");
294 cpl_test_error(CPL_ERROR_NONE);
296
297
298 /* MEDIAN */
299 hpar = hdrl_collapse_parameter_parse_parlist(parlist3, "RECIPE.invalid");
300 cpl_test_null(hpar);
301 cpl_test_error(CPL_ERROR_DATA_NOT_FOUND);
302 hpar = hdrl_collapse_parameter_parse_parlist(parlist3, "RECIPE.collapse");
303 cpl_test_error(CPL_ERROR_NONE);
305
306
307 /* SIGCLIP */
308 hpar = hdrl_collapse_parameter_parse_parlist(parlist4, "RECIPE.invalid");
309 cpl_test_null(hpar);
310 cpl_test_error(CPL_ERROR_DATA_NOT_FOUND);
311 hpar = hdrl_collapse_parameter_parse_parlist(parlist4, "RECIPE.collapse");
312 cpl_test_error(CPL_ERROR_NONE);
313
315 cpl_test(!hdrl_collapse_parameter_is_median(hpar));
318 cpl_test_eq(hdrl_collapse_sigclip_parameter_get_niter(hpar), 5);
320
321
322 /* MINMAX */
323 hpar = hdrl_collapse_parameter_parse_parlist(parlist5, "RECIPE.invalid");
324 cpl_test_null(hpar);
325 cpl_test_error(CPL_ERROR_DATA_NOT_FOUND);
326 hpar = hdrl_collapse_parameter_parse_parlist(parlist5, "RECIPE.collapse");
327 cpl_test_error(CPL_ERROR_NONE);
328
329 cpl_test(hdrl_collapse_parameter_is_minmax(hpar));
330 cpl_test(!hdrl_collapse_parameter_is_median(hpar));
331 cpl_test_eq(hdrl_collapse_minmax_parameter_get_nlow(hpar), 1.);
332 cpl_test_eq(hdrl_collapse_minmax_parameter_get_nhigh(hpar), 2.);
334
335 /* MODE */
336 hpar = hdrl_collapse_parameter_parse_parlist(parlist6, "RECIPE.invalid");
337 cpl_test_null(hpar);
338 cpl_test_error(CPL_ERROR_DATA_NOT_FOUND);
339 hpar = hdrl_collapse_parameter_parse_parlist(parlist6, "RECIPE.collapse");
340 cpl_test_error(CPL_ERROR_NONE);
341
342 cpl_test(hdrl_collapse_parameter_is_mode(hpar));
343 cpl_test(!hdrl_collapse_parameter_is_median(hpar));
344 cpl_test_eq(hdrl_collapse_mode_parameter_get_histo_min(hpar), 10.);
345 cpl_test_eq(hdrl_collapse_mode_parameter_get_histo_max(hpar), 1.);
346 cpl_test_eq(hdrl_collapse_mode_parameter_get_bin_size(hpar), 0);
347 cpl_test_eq(hdrl_collapse_mode_parameter_get_method(hpar), HDRL_MODE_MEDIAN);
350
351
352 /* Clean up */
353 cpl_parameterlist_delete(parlist );
354 cpl_parameterlist_delete(parlist1);
355 cpl_parameterlist_delete(parlist2);
356 cpl_parameterlist_delete(parlist3);
357 cpl_parameterlist_delete(parlist4);
358 cpl_parameterlist_delete(parlist5);
359 cpl_parameterlist_delete(parlist6);
360
361
363 cpl_test(hdrl_collapse_parameter_is_mean(hpar));
366 cpl_test(hdrl_collapse_parameter_is_median(hpar));
371}
372
373
374void test_eout(void)
375{
376 cpl_size n = 40;
377
378
379 /* Method with vector mean */
380 hdrl_collapse_imagelist_to_vector_t * vMethod1;
381 vMethod1 = hdrl_collapse_imagelist_to_vector_mean();
382
383 void *eout1 = hdrl_collapse_imagelist_to_vector_create_eout(vMethod1, n);
384
385 void *eout2 = NULL;
386 hdrl_collapse_imagelist_to_vector_move_eout(vMethod1, eout1, eout2, n);
387
388 hdrl_collapse_imagelist_to_vector_delete_eout(vMethod1, eout2);
389 hdrl_collapse_imagelist_to_vector_delete_eout(vMethod1, eout1);
390 hdrl_collapse_imagelist_to_vector_delete(vMethod1);
391
392
393 /* Method with vector minmax */
394 hdrl_collapse_imagelist_to_vector_t *vMethod2;
395 vMethod2 = hdrl_collapse_imagelist_to_vector_minmax(3., 3.);
396
397 void *eout3 = hdrl_collapse_imagelist_to_vector_create_eout(vMethod2, n);
398
399 hdrl_collapse_imagelist_to_vector_delete_eout(vMethod2, eout3);
400 hdrl_collapse_imagelist_to_vector_delete(vMethod2);
401
402
403 /* Method with vector sigclip */
404 hdrl_collapse_imagelist_to_vector_t *vMethod3;
405 vMethod3 = hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3);
406
407 void *eout4 = hdrl_collapse_imagelist_to_vector_create_eout(vMethod3, n);
408 void *eout5 = hdrl_collapse_imagelist_to_vector_create_eout(vMethod3, n);
409
410 hdrl_collapse_imagelist_to_vector_move_eout(vMethod3, NULL, eout5, n);
411 cpl_test_error(CPL_ERROR_NULL_INPUT);
412 hdrl_collapse_imagelist_to_vector_move_eout(vMethod3, eout4, NULL, n);
413 cpl_test_error(CPL_ERROR_NULL_INPUT);
414 hdrl_collapse_imagelist_to_vector_move_eout(vMethod3, eout4, eout4, n);
415 cpl_test_error(CPL_ERROR_ACCESS_OUT_OF_RANGE);
416 hdrl_collapse_imagelist_to_vector_move_eout(vMethod3, eout4, eout5, n);
417 cpl_test_error(CPL_ERROR_ACCESS_OUT_OF_RANGE);
418 hdrl_collapse_imagelist_to_vector_move_eout(vMethod3, eout5, eout4, n);
419 cpl_test_error(CPL_ERROR_ACCESS_OUT_OF_RANGE);
420 hdrl_collapse_imagelist_to_vector_move_eout(vMethod3, eout5, eout5, n);
421 cpl_test_error(CPL_ERROR_ACCESS_OUT_OF_RANGE);
422
423 hdrl_collapse_imagelist_to_vector_delete_eout(vMethod3, eout4);
424 hdrl_collapse_imagelist_to_vector_delete_eout(vMethod3, eout5);
425 hdrl_collapse_imagelist_to_vector_delete(vMethod3);
426
427
428 /* Method with image */
429 hdrl_collapse_imagelist_to_image_t *vMethod4;
430 vMethod4 = hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3);
431
432 cpl_image *img = cpl_image_new(10, 10, CPL_TYPE_DOUBLE);
433 void *eout6 = hdrl_collapse_imagelist_to_image_create_eout(vMethod4, img);
434 cpl_image_delete(img);
435
436 hdrl_collapse_imagelist_to_image_delete_eout(vMethod4, eout6);
437 hdrl_collapse_imagelist_to_image_delete(vMethod4);
438}
439
440
441void test_results(void)
442{
443
444#define TST_FREE \
445 cpl_image_delete(outimg); outimg = NULL; \
446 cpl_image_delete(outerr); outerr = NULL; \
447 cpl_image_delete(contrib); contrib = NULL; \
448 cpl_vector_delete(voutimg); voutimg = NULL; \
449 cpl_vector_delete(vouterr); vouterr = NULL; \
450 cpl_array_delete(acontrib); acontrib = NULL; \
451 hdrl_collapse_imagelist_to_image_delete(method);
452
453 /* create data, value 5., error +-1 */
454 cpl_imagelist * data = cpl_imagelist_new();
455 cpl_imagelist * errs = cpl_imagelist_new();
456 cpl_size nz = 5;
457 cpl_size nx = 40;
458 cpl_size ny = 37;
459 cpl_image * img = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
460 cpl_image * err = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
461 cpl_image_add_scalar(img, 5.);
462 cpl_image_add_scalar(err, 2.);
463
464 /* create expected results (err / sqrt(nz) for mean) */
465 cpl_image * expect_err = cpl_image_duplicate(err);
466 cpl_image_divide_scalar(expect_err, sqrt(nz));
467 cpl_image * expect_contrib = cpl_image_new(nx, ny, CPL_TYPE_INT);
468 cpl_image_add_scalar(expect_contrib, nz);
469 for (cpl_size i = 0; i < nz; i++) {
470 cpl_imagelist_set(data, cpl_image_duplicate(img),
471 cpl_imagelist_get_size(data));
472 cpl_imagelist_set(errs, cpl_image_duplicate(err),
473 cpl_imagelist_get_size(errs));
474 }
475 cpl_image * outimg, * outerr, * contrib;
476
477 cpl_vector * expect_vimg = cpl_vector_new(nz);
478 cpl_vector * expect_verr = cpl_vector_new(nz);
479 cpl_array * expect_acontrib = cpl_array_new(nz, CPL_TYPE_INT);
480 cpl_vector * voutimg, * vouterr;
481 cpl_array * acontrib;
482 cpl_vector_fill(expect_vimg, 5.);
483 cpl_vector_fill(expect_verr, 2. / sqrt(nx * ny));
484 cpl_array_fill_window_int(expect_acontrib, 0, nz, nx * ny);
485
486 /* test reductions on uniform error cases */
487 {
488 /* mean */
489 hdrl_collapse_imagelist_to_image_t * method =
490 hdrl_collapse_imagelist_to_image_mean();
491 hdrl_collapse_imagelist_to_image_call(method, data, errs,
492 &outimg, &outerr, &contrib,
493 NULL);
494
495 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA);
496 cpl_test_image_abs(outerr, expect_err, HDRL_EPS_ERROR);
497 cpl_test_image_abs(contrib, expect_contrib, 0);
498
499 hdrl_collapse_imagelist_to_vector_t * vmethod =
500 hdrl_collapse_imagelist_to_vector_mean();
501 hdrl_collapse_imagelist_to_vector_call(vmethod, data, errs,
502 &voutimg, &vouterr, &acontrib,
503 NULL);
504 cpl_test_vector_abs(voutimg, expect_vimg, HDRL_EPS_DATA);
505 cpl_test_vector_abs(vouterr, expect_verr, HDRL_EPS_ERROR);
506 cpl_test_array_abs(acontrib, expect_acontrib, 0);
507 hdrl_collapse_imagelist_to_vector_delete(vmethod);
508 TST_FREE;
509
510 /* sigclip */
511 method = hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3);
512 hdrl_collapse_imagelist_to_image_call(method, data, errs,
513 &outimg, &outerr, &contrib,
514 NULL);
515
516 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA);
517 cpl_test_image_abs(outerr, expect_err, HDRL_EPS_ERROR);
518 cpl_test_image_abs(contrib, expect_contrib, 0);
519
520
521 vmethod = hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3);
522 hdrl_collapse_imagelist_to_vector_call(vmethod, data, errs,
523 &voutimg, &vouterr, &acontrib,
524 NULL);
525 cpl_test_vector_abs(voutimg, expect_vimg, HDRL_EPS_DATA);
526 cpl_test_vector_abs(vouterr, expect_verr, HDRL_EPS_ERROR);
527 cpl_test_array_abs(acontrib, expect_acontrib, 0);
528 hdrl_collapse_imagelist_to_vector_delete(vmethod);
529 TST_FREE;
530
531 /* minmax */
532 method = hdrl_collapse_imagelist_to_image_minmax(0., 0.);
533 hdrl_collapse_imagelist_to_image_call(method, data, errs,
534 &outimg, &outerr, &contrib,
535 NULL);
536
537 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA);
538 cpl_test_image_abs(outerr, expect_err, HDRL_EPS_ERROR);
539 cpl_test_image_abs(contrib, expect_contrib, 0);
540
541
542 vmethod = hdrl_collapse_imagelist_to_vector_minmax(0., 0.);
543 hdrl_collapse_imagelist_to_vector_call(vmethod, data, errs,
544 &voutimg, &vouterr, &acontrib,
545 NULL);
546 cpl_test_vector_abs(voutimg, expect_vimg, HDRL_EPS_DATA);
547 cpl_test_vector_abs(vouterr, expect_verr, HDRL_EPS_ERROR);
548 cpl_test_array_abs(acontrib, expect_acontrib, 0);
549 hdrl_collapse_imagelist_to_vector_delete(vmethod);
550 TST_FREE;
551
552 /* weighted mean */
553 method = hdrl_collapse_imagelist_to_image_weighted_mean();
554 hdrl_collapse_imagelist_to_image_call(method, data, errs,
555 &outimg, &outerr, &contrib,
556 NULL);
557
558 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA);
559 cpl_test_image_abs(outerr, expect_err, HDRL_EPS_ERROR);
560 cpl_test_image_abs(contrib, expect_contrib, 0);
561 TST_FREE;
562
563 /* median */
564 cpl_image_multiply_scalar(expect_err, sqrt(CPL_MATH_PI_2));
565 method = hdrl_collapse_imagelist_to_image_median();
566 hdrl_collapse_imagelist_to_image_call(method, data, errs,
567 &outimg, &outerr, &contrib,
568 NULL);
569
570 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA);
571 cpl_test_image_abs(outerr, expect_err, HDRL_EPS_ERROR);
572 cpl_test_image_abs(contrib, expect_contrib, 0);
573
574 cpl_vector_multiply_scalar(expect_verr, sqrt(CPL_MATH_PI_2));
575 vmethod = hdrl_collapse_imagelist_to_vector_median();
576 hdrl_collapse_imagelist_to_vector_call(vmethod, data, errs,
577 &voutimg, &vouterr, &acontrib,
578 NULL);
579 cpl_test_vector_abs(voutimg, expect_vimg, HDRL_EPS_DATA);
580 cpl_test_vector_abs(vouterr, expect_verr, HDRL_EPS_ERROR);
581 cpl_test_array_abs(acontrib, expect_acontrib, 0);
582 hdrl_collapse_imagelist_to_vector_delete(vmethod);
583 TST_FREE;
584 }
585 /* test non uniform error cases */
586 {
587 double v[] = {1, 2, 1, 3, 2};
588 double e[] = {0.5, 0.7, 0.1, 1.0, 0.01};
589 for (cpl_size i = 0; i < nz; i++) {
590 cpl_image * tmp = cpl_imagelist_get(data, i);
591 cpl_image_set(tmp, 1, 1, v[i]);
592 tmp = cpl_imagelist_get(errs, i);
593 cpl_image_set(tmp, 1, 1, e[i]);
594 }
595 hdrl_collapse_imagelist_to_image_t * method =
596 hdrl_collapse_imagelist_to_image_mean();
597 hdrl_collapse_imagelist_to_vector_t * vmethod =
598 hdrl_collapse_imagelist_to_vector_mean();
599 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
600 (pair){1.8, HDRL_EPS_DATA},
601 (pair){0.26458269028793246, HDRL_EPS_ERROR},
602 expect_contrib);
603
604 method = hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3);
605 vmethod = hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3);
606 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
607 (pair){1.8, HDRL_EPS_DATA},
608 (pair){0.26458269028793246, HDRL_EPS_ERROR},
609 expect_contrib);
610
611 method = hdrl_collapse_imagelist_to_image_minmax(1, 1);
612 vmethod = hdrl_collapse_imagelist_to_vector_minmax(1, 1);
613 cpl_image * expect_contrib_minmax =
614 cpl_image_subtract_scalar_create(expect_contrib, 2);
615 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
616 (pair){5. / 3., HDRL_EPS_DATA},
617 (pair){sqrt(0.1 * 0.1 + 0.7 * 0.7 + 0.01 * 0.01) / 3.,
618 HDRL_EPS_ERROR},
619 expect_contrib_minmax);
620 cpl_image_delete(expect_contrib_minmax);
621
622 method = hdrl_collapse_imagelist_to_image_weighted_mean();
623 vmethod = hdrl_collapse_imagelist_to_vector_weighted_mean();
624 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
625 (pair){1.9898090843925733, HDRL_EPS_ERROR},
626 (pair){0.0099469054598625289, HDRL_EPS_ERROR},
627 expect_contrib);
628 }
629 /* test non uniform error cases with rejected values */
630 {
631 double v[] = {1, 2, 1, 3, 2};
632 double e[] = {0.5, 0.7, 0.1, 1.0, 0.01};
633 for (cpl_size i = 0; i < nz; i++) {
634 cpl_image * tmp = cpl_imagelist_get(data, i);
635 cpl_image_set(tmp, 1, 1, v[i]);
636 if (i == 3) {
637 cpl_image_reject(tmp, 1, 1);
638 }
639 tmp = cpl_imagelist_get(errs, i);
640 cpl_image_set(tmp, 1, 1, e[i]);
641 if (i == 3) {
642 cpl_image_reject(tmp, 1, 1);
643 }
644 }
645 cpl_image_delete(expect_contrib);
646 expect_contrib = cpl_image_new_from_accepted(data);
647
648 hdrl_collapse_imagelist_to_image_t * method =
649 hdrl_collapse_imagelist_to_image_mean();
650 hdrl_collapse_imagelist_to_vector_t * vmethod =
651 hdrl_collapse_imagelist_to_vector_mean();
652 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
653 (pair){1.5, HDRL_EPS_DATA},
654 (pair){0.21652078422174625, HDRL_EPS_ERROR},
655 expect_contrib);
656
657 method = hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3);
658 vmethod = hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3);
659 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
660 (pair){1.5, HDRL_EPS_DATA},
661 (pair){0.21652078422174625, HDRL_EPS_ERROR},
662 expect_contrib);
663
664 method = hdrl_collapse_imagelist_to_image_minmax(1, 1);
665 vmethod = hdrl_collapse_imagelist_to_vector_minmax(1, 1);
666 cpl_image * expect_contrib_minmax =
667 cpl_image_subtract_scalar_create(expect_contrib, 2);
668 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
669 (pair){3. / 2., HDRL_EPS_DATA},
670 (pair){sqrt(0.1 * 0.1 + 0.01 * 0.01) / 2.,
671 HDRL_EPS_ERROR},
672 expect_contrib_minmax);
673 cpl_image_delete(expect_contrib_minmax);
674
675 method = hdrl_collapse_imagelist_to_image_weighted_mean();
676 vmethod = hdrl_collapse_imagelist_to_vector_weighted_mean();
677 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
678 (pair){1.9897091252756485, HDRL_EPS_ERROR},
679 (pair){0.0099473975744101273, HDRL_EPS_ERROR},
680 expect_contrib);
681 }
682 /* test all values rejected in line, list to image */
683 {
684 double v[] = {1, 2, 1, 3, 2};
685 double e[] = {0.5, 0.7, 0.1, 1.0, 0.01};
686 for (cpl_size i = 0; i < nz; i++) {
687 cpl_image * tmp = cpl_imagelist_get(data, i);
688 cpl_image_set(tmp, 1, 1, v[i]);
689 cpl_image_reject(tmp, 1, 1);
690 tmp = cpl_imagelist_get(errs, i);
691 cpl_image_set(tmp, 1, 1, e[i]);
692 cpl_image_reject(tmp, 1, 1);
693 }
694 cpl_image_delete(expect_contrib);
695 expect_contrib = cpl_image_new_from_accepted(data);
696
697 hdrl_collapse_imagelist_to_image_t * img_meth[] = {
698 hdrl_collapse_imagelist_to_image_mean(),
699 hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3),
700 hdrl_collapse_imagelist_to_image_weighted_mean(),
701 hdrl_collapse_imagelist_to_image_median(),
702 };
703 hdrl_collapse_imagelist_to_vector_t * vec_meth[] = {
704 hdrl_collapse_imagelist_to_vector_mean(),
705 hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3),
706 hdrl_collapse_imagelist_to_vector_weighted_mean(),
707 hdrl_collapse_imagelist_to_vector_median(),
708 };
709
710 for (size_t i = 0; i < ARRAY_LEN(img_meth); i++) {
711 test_l2i_and_l2v(data, errs, img_meth[i], vec_meth[i], 1, 1,
712 (pair){NAN, HDRL_EPS_DATA},
713 (pair){NAN, HDRL_EPS_ERROR},
714 expect_contrib);
715 }
716
717 hdrl_collapse_imagelist_to_image_t * method =
718 hdrl_collapse_imagelist_to_image_minmax(1., 1.);
719 hdrl_collapse_imagelist_to_vector_t * vmethod =
720 hdrl_collapse_imagelist_to_vector_minmax(1., 1.);
721 cpl_image * expect_contrib_minmax =
722 cpl_image_subtract_scalar_create(expect_contrib, 2);
723 cpl_image_set(expect_contrib_minmax, 1, 1, 0);
724 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
725 (pair){NAN, HDRL_EPS_DATA},
726 (pair){NAN, HDRL_EPS_ERROR},
727 expect_contrib_minmax);
728 cpl_image_delete(expect_contrib_minmax);
729 }
730 /* test all values rejected, list to image */
731 {
732 for (cpl_size i = 0; i < nz; i++) {
733 cpl_image * tmp = cpl_imagelist_get(data, i);
734 cpl_image_accept_all(tmp);
735 cpl_mask_not(cpl_image_get_bpm(tmp));
736 }
737 cpl_image_delete(expect_contrib);
738 expect_contrib = cpl_image_new_from_accepted(data);
739
740 hdrl_collapse_imagelist_to_image_t * img_meth[] = {
741 hdrl_collapse_imagelist_to_image_mean(),
742 hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3),
743 hdrl_collapse_imagelist_to_image_weighted_mean(),
744 hdrl_collapse_imagelist_to_image_median(),
745 hdrl_collapse_imagelist_to_image_minmax(1., 1.),
746 };
747 hdrl_collapse_imagelist_to_vector_t * vec_meth[] = {
748 hdrl_collapse_imagelist_to_vector_mean(),
749 hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3),
750 hdrl_collapse_imagelist_to_vector_weighted_mean(),
751 hdrl_collapse_imagelist_to_vector_median(),
752 hdrl_collapse_imagelist_to_vector_minmax(1., 1.),
753 };
754
755 for (size_t i = 0; i < ARRAY_LEN(img_meth); i++) {
756 cpl_image * outimg_loc, *outerr_loc, *contrib_loc;
757 hdrl_imagelist_combine(data, errs, img_meth[i], &outimg_loc, &outerr_loc, &contrib_loc);
758 /* unlike cpl, don't emit error */
759 cpl_test_error(CPL_ERROR_NONE);
760
761 cpl_test_image_abs(contrib_loc, expect_contrib, 0);
762 cpl_test_eq(cpl_image_count_rejected(outimg_loc), nx *ny);
763 cpl_test_eq(cpl_image_count_rejected(outerr_loc), nx *ny);
764 cpl_image_delete(outimg_loc);
765 cpl_image_delete(outerr_loc);
766 cpl_image_delete(contrib_loc);
767
768 /* also check vector variant */
769 test_l2i_and_l2v(data, errs, img_meth[i], vec_meth[i], 1, 1,
770 (pair){NAN, HDRL_EPS_DATA},
771 (pair){NAN, HDRL_EPS_ERROR},
772 expect_contrib);
773 cpl_test_error(CPL_ERROR_NONE);
774 }
775 }
776 /* test median error propagation with rejects, only makes sense on uniform
777 * errors as the scaling relies on gaussian errors */
778 {
779 double v[] = {1, 2, 1, 3, 2};
780 double e[] = {1., 1., 1., 1., 1.};
781 int d;
782 for (cpl_size i = 0; i < nz; i++) {
783 cpl_image * tmp = cpl_imagelist_get(data, i);
784 cpl_image_set(tmp, 1, 1, v[i]);
785 cpl_image_set(tmp, 2, 2, v[i]);
786 if (i > 1) {
787 cpl_image_reject(tmp, 1, 1);
788 }
789 tmp = cpl_imagelist_get(errs, i);
790 cpl_image_set(tmp, 1, 1, e[i]);
791 cpl_image_set(tmp, 2, 2, e[i]);
792 if (i > 1) {
793 cpl_image_reject(tmp, 1, 1);
794 }
795 }
796 cpl_image_delete(expect_contrib);
797 expect_contrib = cpl_image_new_from_accepted(data);
798
799 hdrl_collapse_imagelist_to_image_t * method =
800 hdrl_collapse_imagelist_to_image_median();
801 hdrl_imagelist_combine(data, errs, method, &outimg, &outerr, &contrib);
802
803 /* contrib > 2 -> sqrt(nz * pi / 2) error scaling */
804 cpl_test_abs(cpl_image_get(outimg, 2, 2, &d), 2., HDRL_EPS_DATA);
805 cpl_test_abs(cpl_image_get(outerr, 2, 2, &d),
806 1. / sqrt(nz) * sqrt(CPL_MATH_PI_2), HDRL_EPS_ERROR);
807 /* contrib <= 2 -> median is a mean, no scaling */
808 cpl_test_abs(cpl_image_get(outimg, 1, 1, &d), 1.5, HDRL_EPS_DATA);
809 cpl_test_abs(cpl_image_get(outerr, 1, 1, &d), 1. / sqrt(2.), HDRL_EPS_ERROR);
810 cpl_test_image_abs(contrib, expect_contrib, 0);
811
812 cpl_imagelist * vl = cpl_imagelist_new();
813 cpl_imagelist * el = cpl_imagelist_new();
814 prep_l2v_input(data, errs, 1, 1, vl, el);
815 hdrl_collapse_imagelist_to_vector_t * vmethod =
816 hdrl_collapse_imagelist_to_vector_median();
817 hdrl_collapse_imagelist_to_vector_call(vmethod, vl, el,
818 &voutimg, &vouterr, &acontrib,
819 NULL);
820
821 cpl_test_abs(cpl_vector_get(voutimg, 0), 1.5, HDRL_EPS_DATA);
822 cpl_test_abs(cpl_vector_get(vouterr, 0), 1. / sqrt(2.), HDRL_EPS_ERROR);
823
824 cpl_test_abs(cpl_array_get_int(acontrib, 0, NULL), 2, 0);
825 cpl_vector_delete(voutimg);
826 cpl_vector_delete(vouterr);
827 cpl_array_delete(acontrib);
828 cpl_imagelist_empty(vl);
829 cpl_imagelist_empty(el);
830
831 prep_l2v_input(data, errs, 2, 2, vl, el);
832 hdrl_collapse_imagelist_to_vector_call(vmethod, vl, el,
833 &voutimg, &vouterr, &acontrib,
834 NULL);
835
836 cpl_test_abs(cpl_vector_get(voutimg, 0), 2., HDRL_EPS_DATA);
837 cpl_test_abs(cpl_vector_get(vouterr, 0),
838 1. / sqrt(nz) * sqrt(CPL_MATH_PI_2), HDRL_EPS_ERROR);
839
840 cpl_test_abs(cpl_array_get_int(acontrib, 0, NULL), 5, 0);
841 cpl_imagelist_delete(vl);
842 cpl_imagelist_delete(el);
843 hdrl_collapse_imagelist_to_vector_delete(vmethod);
844
845 TST_FREE;
846 }
847
848 cpl_imagelist_delete(data);
849 cpl_imagelist_delete(errs);
850 cpl_image_delete(expect_err);
851 cpl_image_delete(expect_contrib);
852 cpl_image_delete(img);
853 cpl_image_delete(err);
854 cpl_vector_delete(expect_vimg);
855 cpl_vector_delete(expect_verr);
856 cpl_array_delete(expect_acontrib);
857
858#undef TST_FREE
859}
860
861void test_results_mode(void)
862{
863#define TST_FREE_MODE \
864 cpl_image_delete(outimg); outimg = NULL; \
865 cpl_image_delete(outerr); outerr = NULL; \
866 cpl_image_delete(contrib); contrib = NULL; \
867 cpl_vector_delete(voutimg); voutimg = NULL; \
868 cpl_vector_delete(vouterr); vouterr = NULL; \
869 cpl_array_delete(acontrib); acontrib = NULL; \
870 hdrl_collapse_imagelist_to_image_delete(method);
871
872 /* create data, value 5., error +-1 */
873 cpl_imagelist * data = cpl_imagelist_new();
874 cpl_imagelist * errs = cpl_imagelist_new();
875 cpl_size nz = 10;
876 cpl_size nx = 2;
877 cpl_size ny = 2;
878 cpl_image * img = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
879 cpl_image * err = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
880 cpl_image_add_scalar(img, 3.5);
881 cpl_image_add_scalar(err, 1.5);
882
883 /* create expected results (err / sqrt(nz) for mean) */
884 cpl_image * expect_err = cpl_image_duplicate(err);
885 cpl_image_divide_scalar(expect_err, sqrt(nz));
886 cpl_image * expect_contrib = cpl_image_new(nx, ny, CPL_TYPE_INT);
887 cpl_image_add_scalar(expect_contrib, nz);
888 for (cpl_size i = 0; i < nz; i++) {
889 cpl_imagelist_set(data, cpl_image_duplicate(img),
890 cpl_imagelist_get_size(data));
891 cpl_imagelist_set(errs, cpl_image_duplicate(err),
892 cpl_imagelist_get_size(errs));
893 }
894 cpl_image * outimg, * outerr, * contrib;
895
896 cpl_vector * expect_vimg = cpl_vector_new(nz);
897 cpl_vector * expect_verr = cpl_vector_new(nz);
898 cpl_array * expect_acontrib = cpl_array_new(nz, CPL_TYPE_INT);
899 cpl_vector * voutimg = NULL, * vouterr = NULL;
900 cpl_array * acontrib = NULL;
901 cpl_vector_fill(expect_vimg, 3.5);
902 cpl_vector_fill(expect_verr, 1.5 / sqrt(nx * ny));
903 cpl_array_fill_window_int(expect_acontrib, 0, nz, nx * ny);
904
905 /* test reductions on uniform error cases */
906 {
907 /* mode MEDIAN */
908
909 cpl_image * mode_expect_err = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
910 cpl_vector * mode_expect_verr = cpl_vector_new(nz);
911 cpl_vector_fill(mode_expect_verr, 0.);
912
913 hdrl_collapse_imagelist_to_image_t * method =
914 hdrl_collapse_imagelist_to_image_mode(5, 5, 0, HDRL_MODE_MEDIAN, 20);
915 hdrl_collapse_imagelist_to_image_call(method, data, errs,
916 &outimg, &outerr, &contrib,
917 NULL);
918
919 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA);
920 cpl_test_image_abs(outerr, mode_expect_err, HDRL_EPS_ERROR);
921 cpl_test_image_abs(contrib, expect_contrib, 0);
922 TST_FREE_MODE;
923
924 method = hdrl_collapse_imagelist_to_image_mode(3, 4, 1, HDRL_MODE_MEDIAN, 20);
925 hdrl_collapse_imagelist_to_image_call(method, data, errs,
926 &outimg, &outerr, &contrib,
927 NULL);
928
929 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA);
930 cpl_test_image_abs(outerr, mode_expect_err, HDRL_EPS_ERROR);
931 cpl_test_image_abs(contrib, expect_contrib, 0);
932 cpl_image_delete(mode_expect_err);
933
934
935 hdrl_collapse_imagelist_to_vector_t * vmethod =
936 hdrl_collapse_imagelist_to_vector_mode(3., 4. ,1., HDRL_MODE_MEDIAN, 20);
937 hdrl_collapse_imagelist_to_vector_call(vmethod, data, errs,
938 &voutimg, &vouterr, &acontrib,
939 NULL);
940 cpl_test_vector_abs(voutimg, expect_vimg, HDRL_EPS_DATA);
941 cpl_test_vector_abs(vouterr, mode_expect_verr, HDRL_EPS_ERROR);
942 cpl_test_array_abs(acontrib, expect_acontrib, 0);
943 hdrl_collapse_imagelist_to_vector_delete(vmethod);
944 cpl_vector_delete(mode_expect_verr);
945 TST_FREE_MODE;
946
947 /* mode WEIGHTED */
948
949 mode_expect_err = cpl_image_new(nx, ny, CPL_TYPE_DOUBLE);
950 mode_expect_verr = cpl_vector_new(nz);
951 cpl_vector_fill(mode_expect_verr, 0.);
952 method = hdrl_collapse_imagelist_to_image_mode(5, 5, 0, HDRL_MODE_WEIGHTED, 20);
953 hdrl_collapse_imagelist_to_image_call(method, data, errs,
954 &outimg, &outerr, &contrib,
955 NULL);
956
957 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA * 2);
958 cpl_test_image_abs(outerr, mode_expect_err, HDRL_EPS_ERROR);
959 cpl_test_image_abs(contrib, expect_contrib, 0);
960 TST_FREE_MODE;
961
962 method = hdrl_collapse_imagelist_to_image_mode(3., 4., 1, HDRL_MODE_WEIGHTED, 20);
963 hdrl_collapse_imagelist_to_image_call(method, data, errs,
964 &outimg, &outerr, &contrib,
965 NULL);
966
967 //cpl_image_add_scalar(img, 0.5);
968 cpl_test_image_abs(outimg, img, HDRL_EPS_DATA);
969 //cpl_image_subtract_scalar(img, 0.5);
970
971 cpl_test_image_abs(outerr, mode_expect_err, HDRL_EPS_ERROR);
972 cpl_test_image_abs(contrib, expect_contrib, 0);
973 cpl_image_delete(mode_expect_err);
974
975
976 vmethod = hdrl_collapse_imagelist_to_vector_mode(3., 4. ,1., HDRL_MODE_WEIGHTED, 20);
977 hdrl_collapse_imagelist_to_vector_call(vmethod, data, errs,
978 &voutimg, &vouterr, &acontrib,
979 NULL);
980 //cpl_vector_add_scalar(expect_vimg, 0.5);
981 cpl_test_vector_abs(voutimg, expect_vimg, HDRL_EPS_DATA);
982 //cpl_vector_subtract_scalar(expect_vimg, 0.5);
983 cpl_test_vector_abs(vouterr, mode_expect_verr, HDRL_EPS_ERROR);
984 cpl_test_array_abs(acontrib, expect_acontrib, 0);
985 hdrl_collapse_imagelist_to_vector_delete(vmethod);
986 cpl_vector_delete(mode_expect_verr);
987 TST_FREE_MODE;
988
989 /* mode FIT */
990
991 method = hdrl_collapse_imagelist_to_image_mode(5, 5, 0, HDRL_MODE_FIT, 20);
992 hdrl_collapse_imagelist_to_image_call(method, data, errs,
993 &outimg, &outerr, &contrib,
994 NULL);
995
996 /* You can not make a fit with one bin only */
997 cpl_test_eq(cpl_image_count_rejected(outimg), nx * ny); /* All bad */
998 cpl_test_eq(cpl_image_count_rejected(outerr), nx * ny); /* All bad */
999 cpl_test_eq(cpl_image_get_sqflux(contrib), 0); /* All Zero */
1000 TST_FREE_MODE;
1001
1002 method = hdrl_collapse_imagelist_to_image_mode(2, 5, 1, HDRL_MODE_FIT, 20);
1003 hdrl_collapse_imagelist_to_image_call(method, data, errs,
1004 &outimg, &outerr, &contrib,
1005 NULL);
1006
1007 cpl_test_eq(cpl_image_count_rejected(outimg), 0); /* All good */
1008 cpl_test_eq(cpl_image_count_rejected(outerr), 0); /* All good */
1009 cpl_test_eq(cpl_image_get_sqflux(contrib), 400.); /* All 4 */
1010
1011 vmethod = hdrl_collapse_imagelist_to_vector_mode(2., 5. ,1., HDRL_MODE_FIT, 20);
1012 hdrl_collapse_imagelist_to_vector_call(vmethod, data, errs,
1013 &voutimg, &vouterr, &acontrib,
1014 NULL);
1015
1016 /* You can not make a fit with one bin only */
1017 cpl_test_eq(cpl_image_count_rejected(outimg), 0); /* All good */
1018 cpl_test_eq(cpl_image_count_rejected(outerr), 0); /* All good */
1019 cpl_test_eq(cpl_image_get_sqflux(contrib), 400.); /* 4 * 10*10 */
1020 hdrl_collapse_imagelist_to_vector_delete(vmethod);
1021 TST_FREE_MODE;
1022
1023 }
1024 /* test non uniform error cases */
1025 {
1026 double v[] = {1.5, 2.5, 2.5, 3.5, 3.5, 3.5, 3.5, 4.5, 4.5, 5.5};
1027 double e[] = {0, 0, 0, 0, 0, 0, 0, 0, 0, 0}; /* Not used */
1028 for (cpl_size i = 0; i < nz; i++) {
1029 cpl_image * tmp = cpl_imagelist_get(data, i);
1030 cpl_image_set(tmp, 1, 1, v[i]);
1031 tmp = cpl_imagelist_get(errs, i);
1032 cpl_image_set(tmp, 1, 1, e[i]);
1033 }
1034
1035 /* MODE MEDIAN */
1036
1037 /* histo_max, histo_min and binsize is determined internally by the
1038 * algorithm */
1039 hdrl_collapse_imagelist_to_image_t * method =
1040 hdrl_collapse_imagelist_to_image_mode(0, 0, 0, HDRL_MODE_MEDIAN, 0);
1041 hdrl_collapse_imagelist_to_vector_t * vmethod =
1042 hdrl_collapse_imagelist_to_vector_mode(0, 0, 0, HDRL_MODE_MEDIAN, 0);
1043
1044 /* The code that deives automatically the bines gives:
1045 Histogram bin size: 4.80337 min: -0.901685 max: 8.70506 number of bins: 2
1046 Therefore the final vector in the innermost bin is:
1047 [1.5, 2.5, 2.5, 3.5, 3.5, 3.5, 3.5,]
1048 The error is the standard deviation of the final vector (using cpl) */
1049
1050
1051 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1052 (pair){3.5, HDRL_EPS_DATA},
1053 (pair){0.7867957924694431, HDRL_EPS_ERROR},
1054 expect_contrib);
1055
1056 /* histo_max, histo_min and binsize is given by the user */
1057 method = hdrl_collapse_imagelist_to_image_mode(3, 4, 1, HDRL_MODE_MEDIAN, 0);
1058 vmethod = hdrl_collapse_imagelist_to_vector_mode(3, 4, 1, HDRL_MODE_MEDIAN, 0);
1059
1060 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1061 (pair){3.5, HDRL_EPS_DATA},
1062 (pair){0., HDRL_EPS_ERROR},
1063 expect_contrib);
1064
1065 /* MODE WEIGHTED */
1066
1067 /* histo_max, and histo_min are determined internally by the
1068 * algorithm */
1069 method =
1070 hdrl_collapse_imagelist_to_image_mode(0, 0, 1, HDRL_MODE_WEIGHTED, 0);
1071 vmethod =
1072 hdrl_collapse_imagelist_to_vector_mode(0, 0, 1, HDRL_MODE_WEIGHTED, 0);
1073
1074 /* The code that deives automatically the bines gives:
1075 Histogram bin size: 4.80337 min: -0.901685 max: 8.70506 number of bins: 2
1076 Therefore the final vector in the innermost bin is:
1077 [1.5, 2.5, 2.5, 3.5, 3.5, 3.5, 3.5,]
1078 The error is the standard deviation of the final vector (using cpl) */
1079
1080 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1081 (pair){3.5, HDRL_EPS_DATA},
1082 (pair){0.4330127018922193, HDRL_EPS_ERROR},
1083 expect_contrib);
1084
1085
1086
1087 /* histo_max, histo_min and binsize is given by the user */
1088 method = hdrl_collapse_imagelist_to_image_mode(3, 4, 1, HDRL_MODE_WEIGHTED, 0);
1089 vmethod = hdrl_collapse_imagelist_to_vector_mode(3, 4, 1, HDRL_MODE_WEIGHTED, 0);
1090
1091 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1092 (pair){3.5, HDRL_EPS_DATA},
1093 (pair){0.1767766952966369, HDRL_EPS_ERROR},
1094 expect_contrib);
1095
1096
1097 /* MODE FIT */
1098
1099 /* histo_max, histo_min and binsize is given by the user */
1100 method = hdrl_collapse_imagelist_to_image_mode(2, 5, 1, HDRL_MODE_FIT, 0);
1101 vmethod = hdrl_collapse_imagelist_to_vector_mode(2, 5, 1, HDRL_MODE_FIT, 0);
1102
1103 test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1104 (pair){3.6, HDRL_EPS_DATA * 10},
1105 (pair){0.2683281572999749, HDRL_EPS_ERROR * 10.},
1106 expect_contrib);
1107
1108 }
1109
1110
1111// /* test non uniform error cases with rejected values */
1112// if (0){
1113// double v[] = {1, 2, 1, 3, 2};
1114// double e[] = {0.5, 0.7, 0.1, 1.0, 0.01};
1115// for (cpl_size i = 0; i < nz; i++) {
1116// cpl_image * tmp = cpl_imagelist_get(data, i);
1117// cpl_image_set(tmp, 1, 1, v[i]);
1118// if (i == 3) {
1119// cpl_image_reject(tmp, 1, 1);
1120// }
1121// tmp = cpl_imagelist_get(errs, i);
1122// cpl_image_set(tmp, 1, 1, e[i]);
1123// if (i == 3) {
1124// cpl_image_reject(tmp, 1, 1);
1125// }
1126// }
1127// cpl_image_delete(expect_contrib);
1128// expect_contrib = cpl_image_new_from_accepted(data);
1129//
1130// hdrl_collapse_imagelist_to_image_t * method =
1131// hdrl_collapse_imagelist_to_image_mean();
1132// hdrl_collapse_imagelist_to_vector_t * vmethod =
1133// hdrl_collapse_imagelist_to_vector_mean();
1134// test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1135// (pair){1.5, HDRL_EPS_DATA},
1136// (pair){0.21652078422174625, HDRL_EPS_ERROR},
1137// expect_contrib);
1138//
1139// method = hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3);
1140// vmethod = hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3);
1141// test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1142// (pair){1.5, HDRL_EPS_DATA},
1143// (pair){0.21652078422174625, HDRL_EPS_ERROR},
1144// expect_contrib);
1145//
1146// method = hdrl_collapse_imagelist_to_image_minmax(1, 1);
1147// vmethod = hdrl_collapse_imagelist_to_vector_minmax(1, 1);
1148// cpl_image * expect_contrib_minmax =
1149// cpl_image_subtract_scalar_create(expect_contrib, 2);
1150// test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1151// (pair){3. / 2., HDRL_EPS_DATA},
1152// (pair){sqrt(0.1 * 0.1 + 0.01 * 0.01) / 2.,
1153// HDRL_EPS_ERROR},
1154// expect_contrib_minmax);
1155// cpl_image_delete(expect_contrib_minmax);
1156//
1157// method = hdrl_collapse_imagelist_to_image_weighted_mean();
1158// vmethod = hdrl_collapse_imagelist_to_vector_weighted_mean();
1159// test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1160// (pair){1.9897091252756485, HDRL_EPS_ERROR},
1161// (pair){0.0099473975744101273, HDRL_EPS_ERROR},
1162// expect_contrib);
1163// }
1164// /* test all values rejected in line, list to image */
1165// if (0){
1166// double v[] = {1, 2, 1, 3, 2};
1167// double e[] = {0.5, 0.7, 0.1, 1.0, 0.01};
1168// for (cpl_size i = 0; i < nz; i++) {
1169// cpl_image * tmp = cpl_imagelist_get(data, i);
1170// cpl_image_set(tmp, 1, 1, v[i]);
1171// cpl_image_reject(tmp, 1, 1);
1172// tmp = cpl_imagelist_get(errs, i);
1173// cpl_image_set(tmp, 1, 1, e[i]);
1174// cpl_image_reject(tmp, 1, 1);
1175// }
1176// cpl_image_delete(expect_contrib);
1177// expect_contrib = cpl_image_new_from_accepted(data);
1178//
1179// hdrl_collapse_imagelist_to_image_t * img_meth[] = {
1180// hdrl_collapse_imagelist_to_image_mean(),
1181// hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3),
1182// hdrl_collapse_imagelist_to_image_weighted_mean(),
1183// hdrl_collapse_imagelist_to_image_median(),
1184// };
1185// hdrl_collapse_imagelist_to_vector_t * vec_meth[] = {
1186// hdrl_collapse_imagelist_to_vector_mean(),
1187// hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3),
1188// hdrl_collapse_imagelist_to_vector_weighted_mean(),
1189// hdrl_collapse_imagelist_to_vector_median(),
1190// };
1191//
1192// for (size_t i = 0; i < ARRAY_LEN(img_meth); i++) {
1193// test_l2i_and_l2v(data, errs, img_meth[i], vec_meth[i], 1, 1,
1194// (pair){NAN, HDRL_EPS_DATA},
1195// (pair){NAN, HDRL_EPS_ERROR},
1196// expect_contrib);
1197// }
1198//
1199// hdrl_collapse_imagelist_to_image_t * method =
1200// hdrl_collapse_imagelist_to_image_minmax(1., 1.);
1201// hdrl_collapse_imagelist_to_vector_t * vmethod =
1202// hdrl_collapse_imagelist_to_vector_minmax(1., 1.);
1203// cpl_image * expect_contrib_minmax =
1204// cpl_image_subtract_scalar_create(expect_contrib, 2);
1205// cpl_image_set(expect_contrib_minmax, 1, 1, 0);
1206// test_l2i_and_l2v(data, errs, method, vmethod, 1, 1,
1207// (pair){NAN, HDRL_EPS_DATA},
1208// (pair){NAN, HDRL_EPS_ERROR},
1209// expect_contrib_minmax);
1210// cpl_image_delete(expect_contrib_minmax);
1211// }
1212// /* test all values rejected, list to image */
1213// if (0){
1214// for (cpl_size i = 0; i < nz; i++) {
1215// cpl_image * tmp = cpl_imagelist_get(data, i);
1216// cpl_image_accept_all(tmp);
1217// cpl_mask_not(cpl_image_get_bpm(tmp));
1218// }
1219// cpl_image_delete(expect_contrib);
1220// expect_contrib = cpl_image_new_from_accepted(data);
1221//
1222// hdrl_collapse_imagelist_to_image_t * img_meth[] = {
1223// hdrl_collapse_imagelist_to_image_mean(),
1224// hdrl_collapse_imagelist_to_image_sigclip(3., 3., 3),
1225// hdrl_collapse_imagelist_to_image_weighted_mean(),
1226// hdrl_collapse_imagelist_to_image_median(),
1227// hdrl_collapse_imagelist_to_image_minmax(1., 1.),
1228// };
1229// hdrl_collapse_imagelist_to_vector_t * vec_meth[] = {
1230// hdrl_collapse_imagelist_to_vector_mean(),
1231// hdrl_collapse_imagelist_to_vector_sigclip(3., 3., 3),
1232// hdrl_collapse_imagelist_to_vector_weighted_mean(),
1233// hdrl_collapse_imagelist_to_vector_median(),
1234// hdrl_collapse_imagelist_to_vector_minmax(1., 1.),
1235// };
1236//
1237// for (size_t i = 0; i < ARRAY_LEN(img_meth); i++) {
1238// cpl_image * outimg_loc, *outerr_loc, *contrib_loc;
1239// hdrl_imagelist_combine(data, errs, img_meth[i], &outimg_loc, &outerr_loc, &contrib_loc);
1240// /* unlike cpl, don't emit error */
1241// cpl_test_error(CPL_ERROR_NONE);
1242//
1243// cpl_test_image_abs(contrib_loc, expect_contrib, 0);
1244// cpl_test_eq(cpl_image_count_rejected(outimg_loc), nx *ny);
1245// cpl_test_eq(cpl_image_count_rejected(outerr_loc), nx *ny);
1246// cpl_image_delete(outimg_loc);
1247// cpl_image_delete(outerr_loc);
1248// cpl_image_delete(contrib_loc);
1249//
1250// /* also check vector variant */
1251// test_l2i_and_l2v(data, errs, img_meth[i], vec_meth[i], 1, 1,
1252// (pair){NAN, HDRL_EPS_DATA},
1253// (pair){NAN, HDRL_EPS_ERROR},
1254// expect_contrib);
1255// cpl_test_error(CPL_ERROR_NONE);
1256// }
1257// }
1258// /* test median error propagation with rejects, only makes sense on uniform
1259// * errors as the scaling relies on gaussian errors */
1260// if (0){
1261// double v[] = {1, 2, 1, 3, 2};
1262// double e[] = {1., 1., 1., 1., 1.};
1263// int d;
1264// for (cpl_size i = 0; i < nz; i++) {
1265// cpl_image * tmp = cpl_imagelist_get(data, i);
1266// cpl_image_set(tmp, 1, 1, v[i]);
1267// cpl_image_set(tmp, 2, 2, v[i]);
1268// if (i > 1) {
1269// cpl_image_reject(tmp, 1, 1);
1270// }
1271// tmp = cpl_imagelist_get(errs, i);
1272// cpl_image_set(tmp, 1, 1, e[i]);
1273// cpl_image_set(tmp, 2, 2, e[i]);
1274// if (i > 1) {
1275// cpl_image_reject(tmp, 1, 1);
1276// }
1277// }
1278// cpl_image_delete(expect_contrib);
1279// expect_contrib = cpl_image_new_from_accepted(data);
1280//
1281// hdrl_collapse_imagelist_to_image_t * method =
1282// hdrl_collapse_imagelist_to_image_median();
1283// hdrl_imagelist_combine(data, errs, method, &outimg, &outerr, &contrib);
1284//
1285// /* contrib > 2 -> sqrt(nz * pi / 2) error scaling */
1286// cpl_test_abs(cpl_image_get(outimg, 2, 2, &d), 2., HDRL_EPS_DATA);
1287// cpl_test_abs(cpl_image_get(outerr, 2, 2, &d),
1288// 1. / sqrt(nz) * sqrt(CPL_MATH_PI_2), HDRL_EPS_ERROR);
1289// /* contrib <= 2 -> median is a mean, no scaling */
1290// cpl_test_abs(cpl_image_get(outimg, 1, 1, &d), 1.5, HDRL_EPS_DATA);
1291// cpl_test_abs(cpl_image_get(outerr, 1, 1, &d), 1. / sqrt(2.), HDRL_EPS_ERROR);
1292// cpl_test_image_abs(contrib, expect_contrib, 0);
1293//
1294// cpl_imagelist * vl = cpl_imagelist_new();
1295// cpl_imagelist * el = cpl_imagelist_new();
1296// prep_l2v_input(data, errs, 1, 1, vl, el);
1297// hdrl_collapse_imagelist_to_vector_t * vmethod =
1298// hdrl_collapse_imagelist_to_vector_median();
1299// hdrl_collapse_imagelist_to_vector_call(vmethod, vl, el,
1300// &voutimg, &vouterr, &acontrib,
1301// NULL);
1302//
1303// cpl_test_abs(cpl_vector_get(voutimg, 0), 1.5, HDRL_EPS_DATA);
1304// cpl_test_abs(cpl_vector_get(vouterr, 0), 1. / sqrt(2.), HDRL_EPS_ERROR);
1305//
1306// cpl_test_abs(cpl_array_get_int(acontrib, 0, NULL), 2, 0);
1307// cpl_vector_delete(voutimg);
1308// cpl_vector_delete(vouterr);
1309// cpl_array_delete(acontrib);
1310// cpl_imagelist_empty(vl);
1311// cpl_imagelist_empty(el);
1312//
1313// prep_l2v_input(data, errs, 2, 2, vl, el);
1314// hdrl_collapse_imagelist_to_vector_call(vmethod, vl, el,
1315// &voutimg, &vouterr, &acontrib,
1316// NULL);
1317//
1318// cpl_test_abs(cpl_vector_get(voutimg, 0), 2., HDRL_EPS_DATA);
1319// cpl_test_abs(cpl_vector_get(vouterr, 0),
1320// 1. / sqrt(nz) * sqrt(CPL_MATH_PI_2), HDRL_EPS_ERROR);
1321//
1322// cpl_test_abs(cpl_array_get_int(acontrib, 0, NULL), 5, 0);
1323// cpl_imagelist_delete(vl);
1324// cpl_imagelist_delete(el);
1325// hdrl_collapse_imagelist_to_vector_delete(vmethod);
1326//
1327// TST_FREE_MODE_MODE;
1328// }
1329
1330 cpl_imagelist_delete(data);
1331 cpl_imagelist_delete(errs);
1332 cpl_image_delete(expect_err);
1333 cpl_image_delete(expect_contrib);
1334 cpl_image_delete(img);
1335 cpl_image_delete(err);
1336 cpl_vector_delete(expect_vimg);
1337 cpl_vector_delete(expect_verr);
1338 cpl_array_delete(expect_acontrib);
1339
1340#undef TST_FREE_MODE
1341}
1342
1343
1344/*----------------------------------------------------------------------------*/
1348/*----------------------------------------------------------------------------*/
1349int main(void)
1350{
1351 cpl_test_init(PACKAGE_BUGREPORT, CPL_MSG_WARNING);
1352
1353 test_parameters();
1354 test_parlist();
1355 test_eout();
1356 test_results();
1357 test_results_mode();
1358 return cpl_test_end(0);
1359}
1360
1361
1362
double hdrl_collapse_mode_parameter_get_bin_size(const hdrl_parameter *p)
get size of the histogram bins
double hdrl_collapse_mode_parameter_get_histo_min(const hdrl_parameter *p)
get min value
cpl_boolean hdrl_collapse_parameter_is_weighted_mean(const hdrl_parameter *self)
check if parameter is a weighted mean parameter
double hdrl_collapse_sigclip_parameter_get_kappa_low(const hdrl_parameter *p)
get low kappa
hdrl_parameter * hdrl_collapse_mean_parameter_create(void)
create a parameter object for mean
hdrl_parameter * hdrl_collapse_sigclip_parameter_create(double kappa_low, double kappa_high, int niter)
create a parameter object for sigclipped mean
hdrl_parameter * hdrl_collapse_weighted_mean_parameter_create(void)
create a parameter object for weighted mean
cpl_boolean hdrl_collapse_parameter_is_mean(const hdrl_parameter *self)
check if parameter is a mean parameter
cpl_boolean hdrl_collapse_parameter_is_median(const hdrl_parameter *self)
check if parameter is a median parameter
double hdrl_collapse_mode_parameter_get_histo_max(const hdrl_parameter *p)
get high value
int hdrl_collapse_sigclip_parameter_get_niter(const hdrl_parameter *p)
get maximum number of clipping iterations
cpl_size hdrl_collapse_mode_parameter_get_error_niter(const hdrl_parameter *p)
get the error type of the mode
cpl_boolean hdrl_collapse_parameter_is_minmax(const hdrl_parameter *self)
check if parameter is a minmax mean parameter
cpl_boolean hdrl_collapse_parameter_is_mode(const hdrl_parameter *self)
check if parameter is a mode parameter
cpl_boolean hdrl_collapse_parameter_is_sigclip(const hdrl_parameter *self)
check if parameter is a sigclip mean parameter
double hdrl_collapse_sigclip_parameter_get_kappa_high(const hdrl_parameter *p)
get high kappa
double hdrl_collapse_minmax_parameter_get_nlow(const hdrl_parameter *p)
get low value
hdrl_parameter * hdrl_collapse_median_parameter_create(void)
create a parameter object for median
hdrl_parameter * hdrl_collapse_mode_parameter_create(double histo_min, double histo_max, double bin_size, hdrl_mode_type mode_method, cpl_size error_niter)
create a parameter object for the mode
cpl_parameterlist * hdrl_collapse_parameter_create_parlist(const char *base_context, const char *prefix, const char *method_def, hdrl_parameter *sigclip_def, hdrl_parameter *minmax_def, hdrl_parameter *mode_def)
Create parameters for the collapse.
hdrl_mode_type hdrl_collapse_mode_parameter_get_method(const hdrl_parameter *p)
get the mode determination method
hdrl_parameter * hdrl_collapse_minmax_parameter_create(double nlow, double nhigh)
create a parameter object for min-max rejected mean
hdrl_parameter * hdrl_collapse_parameter_parse_parlist(const cpl_parameterlist *parlist, const char *prefix)
parse parameterlist for imagelist reduction method
double hdrl_collapse_minmax_parameter_get_nhigh(const hdrl_parameter *p)
get high value
void hdrl_parameter_delete(hdrl_parameter *obj)
shallow delete of a parameter