CR2RE Pipeline Reference Manual 1.6.10
hdrl_mode.c
1/*
2 * hdrl_mode.c
3 *
4 * Created on: Mar 1, 2021
5 * Author: amodigli, agabasch
6 */
7
8/*
9 * This file is part of the HDRL
10 * Copyright (C) 2021 European Southern Observatory
11 *
12 * This program is free software; you can redistribute it and/or modify
13 * it under the terms of the GNU General Public License as published by
14 * the Free Software Foundation; either version 2 of the License, or
15 * (at your option) any later version.
16 *
17 * This program is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
20 * GNU General Public License for more details.
21 *
22 * You should have received a copy of the GNU General Public License
23 * along with this program; if not, write to the Free Software
24 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 */
26
27#ifdef HAVE_CONFIG_H
28#include <config.h>
29#endif
30
31/*-----------------------------------------------------------------------------
32 Includes
33-----------------------------------------------------------------------------*/
34#include "hdrl_mode.h"
35#include "hdrl_sigclip.h"
36#include "hdrl_utils.h"
37#include "hdrl_collapse.h"
38#include "hdrl_random.h"
39
40#ifndef _OPENMP
41#define omp_get_max_threads() 1
42#define omp_get_thread_num() 0
43#else
44#include <omp.h>
45#endif
46
47#include <cpl.h>
48//#include <cpl_vector.h>
49#include <string.h>
50#include <math.h>
51#include <time.h>
52#include <gsl/gsl_multifit.h>
53#include <gsl/gsl_vector.h>
54#include <gsl/gsl_histogram.h>
55#include <gsl/gsl_poly.h>
56
57static const char * method_to_string(const hdrl_mode_type method)
58{
59 switch (method) {
60 case HDRL_MODE_MEDIAN:
61 return "MEDIAN";
62 break;
63 case HDRL_MODE_WEIGHTED:
64 return "WEIGHTED";
65 break;
66 case HDRL_MODE_FIT:
67 return "FIT";
68 break;
69 default :
70 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
71 "mode method unknown");
72 return "";
73 break;
74 }
75}
76
77/* ---------------------------------------------------------------------------*/
95/* ---------------------------------------------------------------------------*/
96cpl_parameterlist * hdrl_mode_parameter_create_parlist(
97 const char *base_context,
98 const char *prefix,
99 const hdrl_parameter *defaults)
100{
101 cpl_ensure(base_context && prefix && defaults,
102 CPL_ERROR_NULL_INPUT, NULL);
103
104 cpl_ensure(hdrl_collapse_parameter_is_mode(defaults),
105 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
106
107 cpl_parameterlist * parlist = cpl_parameterlist_new();
108
109 /* --prefix.histo-min */
110 hdrl_setup_vparameter(parlist, prefix, ".", "",
111 "histo-min", base_context,
112 "Minimum pixel value to accept for mode computation",
113 CPL_TYPE_DOUBLE,
115
116 /* --prefix.histo-max */
117 hdrl_setup_vparameter(parlist, prefix, ".", "",
118 "histo-max", base_context,
119 "Maximum pixel value to accept for mode computation",
120 CPL_TYPE_DOUBLE,
122
123 /* --prefix.bin-size */
124 hdrl_setup_vparameter(parlist, prefix, ".", "",
125 "bin-size", base_context,
126 "Binsize of the histogram",
127 CPL_TYPE_DOUBLE,
129
130 /* --prefix.method */
131 char * name ;
132 char * context = hdrl_join_string(".", 2, base_context, prefix);
133
134 hdrl_mode_type method = hdrl_collapse_mode_parameter_get_method(defaults);
135 const char * method_def = method_to_string(method);
136 name = hdrl_join_string(".", 2, context, "method");
137 cpl_free(context);
138 cpl_parameter * par = cpl_parameter_new_enum(name, CPL_TYPE_STRING,
139 "Mode method (algorithm) to use",
140 base_context, method_def, 3, "MEDIAN", "WEIGHTED", "FIT");
141 cpl_free(name);
142 name = hdrl_join_string(".", 2, prefix, "method");
143 cpl_parameter_set_alias(par, CPL_PARAMETER_MODE_CLI, name);
144 cpl_parameter_disable(par, CPL_PARAMETER_MODE_ENV);
145 cpl_free(name);
146 cpl_parameterlist_append(parlist, par);
147
148/*
149 --prefix.method
150 hdrl_setup_vparameter(parlist, prefix, ".", "",
151 "method", base_context,
152 "Mode method",
153 CPL_TYPE_INT,
154 hdrl_collapse_mode_parameter_get_method(defaults));
155*/
156
157 /* --prefix.error-niter */
158 hdrl_setup_vparameter(parlist, prefix, ".", "",
159 "error-niter", base_context,
160 "Iterations to compute the mode error",
161 CPL_TYPE_INT,
163
164 if (cpl_error_get_code()) {
165 cpl_parameterlist_delete(parlist);
166 return NULL;
167 }
168
169 return parlist;
170}
171/* ---------------------------------------------------------------------------*/
192/* ---------------------------------------------------------------------------*/
193cpl_error_code hdrl_mode_parameter_parse_parlist(
194 const cpl_parameterlist * parlist,
195 const char * prefix,
196 double * histo_min,
197 double * histo_max,
198 double * bin_size,
199 hdrl_mode_type * method,
200 cpl_size * error_niter)
201{
202 cpl_ensure_code(prefix && parlist, CPL_ERROR_NULL_INPUT);
203 char * name;
204
205 if (histo_min) {
206 name = hdrl_join_string(".", 2, prefix, "mode.histo-min");
207 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
208 *histo_min = cpl_parameter_get_double(par);
209 cpl_free(name);
210 }
211 if (histo_max) {
212 name = hdrl_join_string(".", 2, prefix, "mode.histo-max");
213 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
214 *histo_max = cpl_parameter_get_double(par);
215 cpl_free(name);
216 }
217
218 if (bin_size) {
219 name = hdrl_join_string(".", 2, prefix, "mode.bin-size");
220 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
221 *bin_size = cpl_parameter_get_double(par);
222 cpl_free(name);
223 }
224
225
226 /* --method */
227 if (method) {
228 name = hdrl_join_string(".", 2, prefix, "mode.method");
229 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name) ;
230 const char * tmp_str = cpl_parameter_get_string(par);
231 if (tmp_str == NULL) {
232 cpl_free(name);
233 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
234 "Parameter mode.method not found");
235 }
236 if(!strcmp(tmp_str, "MEDIAN")) {
237 *method = HDRL_MODE_MEDIAN;
238 } else if(!strcmp(tmp_str, "WEIGHTED")) {
239 *method = HDRL_MODE_WEIGHTED;
240 } else if(!strcmp(tmp_str, "FIT")) {
241 *method = HDRL_MODE_FIT;
242 }
243 cpl_free(name) ;
244
245 }
246
247/*
248 if (method) {
249 name = hdrl_join_string(".", 2, prefix, "mode.method");
250 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
251 *method = cpl_parameter_get_int(par);
252 cpl_free(name);
253 }
254*/
255
256 if (error_niter) {
257 name = hdrl_join_string(".", 2, prefix, "mode.error-niter");
258 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
259 *error_niter = cpl_parameter_get_int(par);
260 cpl_free(name);
261 }
262
263 if (cpl_error_get_code()) {
264 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
265 "Error while parsing parameterlist "
266 "with prefix %s", prefix);
267 }
268
269 return CPL_ERROR_NONE;
270}
271
284static double
285hdrl_mode_compute_binsize(cpl_vector* vec)
286{
287
288 cpl_size size = cpl_vector_get_size(vec);
289 double mad = 0;
290 hcpl_vector_get_mad_window(vec, 1, size, &mad);
291 double std_deviation_from_mad = mad * CPL_MATH_STD_MAD;
292 /*
293 cpl_msg_debug(cpl_func, "mad: %g standard deviation from mad: %16.10f",
294 mad, std_deviation_from_mad);
295 */
296
297 double binsize = 2. * 3.49 * std_deviation_from_mad / pow(size, 1. / 3. );
298
299 if(binsize <= 0){
300 binsize = nextafter(0,1.0);
301 }
302
303 return binsize ;
304
305}
313static cpl_size
314hdrl_mode_get_nbin(const double min, const double max, const double size)
315{
316
317 return (cpl_size)floor( (max - min) / size ) + 1;
318
319}
320
333static gsl_histogram *
334hdrl_mode_histogram(const cpl_vector* vec, const double histogram_min,
335 const double histogram_max, const cpl_size nbin)
336{
337
338 cpl_error_ensure(nbin > 0, CPL_ERROR_ILLEGAL_INPUT,
339 return NULL, "Number of bins must be > 0");
340
341 cpl_error_ensure(histogram_max > histogram_min, CPL_ERROR_ILLEGAL_INPUT,
342 return NULL, "histo_max must be larger than histo_min");
343
344 gsl_histogram * histogram = gsl_histogram_alloc(nbin);
345
346 gsl_histogram_set_ranges_uniform (histogram, histogram_min, histogram_max);
347
348 const cpl_size size = cpl_vector_get_size(vec);
349 const double* data = cpl_vector_get_data_const(vec);
350 for(cpl_size i = 0; i < size; i++) {
351
352 gsl_histogram_increment(histogram, data[i]);
353
354 }
355
356 return histogram;
357}
358/*-------------------------------------------------------------------------*/
379static cpl_table *
380hdrl_mode_histogram_to_table(gsl_histogram* histogram,
381 const double histogram_min,
382 const double histogram_step, const cpl_size nbin){
383
384 cpl_table* table = cpl_table_new(nbin);
385
386 cpl_table_new_column(table, "BIN", CPL_TYPE_DOUBLE);
387 cpl_table_new_column(table, "INTERVAL_LOWER", CPL_TYPE_DOUBLE);
388 cpl_table_new_column(table, "INTERVAL_UPPER", CPL_TYPE_DOUBLE);
389 cpl_table_new_column(table, "COUNTS", CPL_TYPE_DOUBLE);
390
391 cpl_table_fill_column_window(table, "BIN", 0, nbin, 0);
392 cpl_table_fill_column_window(table, "INTERVAL_LOWER", 0, nbin, 0);
393 cpl_table_fill_column_window(table, "INTERVAL_UPPER", 0, nbin, 0);
394 cpl_table_fill_column_window(table, "COUNTS", 0, nbin, 0);
395
396 double* phb = cpl_table_get_data_double(table, "BIN");
397 double* phl = cpl_table_get_data_double(table, "INTERVAL_LOWER");
398 double* phu = cpl_table_get_data_double(table, "INTERVAL_UPPER");
399 double* phc = cpl_table_get_data_double(table, "COUNTS");
400
401 //cpl_msg_debug(cpl_func,"histogram bin step: %16.10f", histo_step);
402 for(cpl_size i = 0; i < nbin; i++) {
403 phb[i] = (double)i;
404 phl[i] = histogram_min + phb[i] * histogram_step;
405 phu[i] = phl[i] + histogram_step;
406 phc[i] = histogram->bin[i];
407 }
408 return table;
409}
410/*----------------------------------------------------------------------------*/
429/*----------------------------------------------------------------------------*/
430static gsl_matrix *
431hdrl_gsl_fit_poly (const double* data_x,
432 const double* data_y,
433 const double* errs_y,
434 const size_t n_sampling_points,
435 const cpl_size poly_degree,
436 double *fit,
437 double *coeffs_val,
438 double *coeffs_err,
439 double *chisq) {
440
441 cpl_size i, j;
442 double err, chisq_loc;
443 gsl_vector *x = NULL, *y = NULL, *w = NULL, *p = NULL;
444 gsl_matrix *X = NULL, *covar = NULL;
445
446 //double *fit;
447 gsl_multifit_linear_workspace *work = NULL;
448
449 x = gsl_vector_alloc (n_sampling_points);
450 y = gsl_vector_alloc (n_sampling_points);
451 w = gsl_vector_alloc (n_sampling_points);
452 p = gsl_vector_alloc (poly_degree);
453 X = gsl_matrix_alloc (n_sampling_points, poly_degree);
454 covar = gsl_matrix_alloc (poly_degree, poly_degree);
455
456 for (i = 0; i < (cpl_size) n_sampling_points; i++) {
457 gsl_vector_set (x, i, data_x[i]);
458 gsl_vector_set (y, i, data_y[i]);
459 err = errs_y[i];
460
461 gsl_vector_set (w, i, 1.0 / err / err);
462
463 for (j = 0; j < poly_degree; j++)
464 gsl_matrix_set (X, i, j, gsl_pow_int(gsl_vector_get(x, i), j));
465 }
466
467 work = gsl_multifit_linear_alloc (n_sampling_points, poly_degree);
468 gsl_multifit_wlinear (X, w, y, p, covar, &chisq_loc, work);
469 gsl_multifit_linear_free (work);
470
471 for (i = 0; i < (cpl_size) n_sampling_points; i++) {
472
473 fit[i] = 0.;
474 for (j = 0; j < poly_degree; j++)
475 fit[i] += gsl_matrix_get (X, i, j) * gsl_vector_get (p, j);
476 }
477
478 *chisq = chisq_loc/(n_sampling_points - poly_degree);
479 /*
480 * NOTE we scale covar arbitrarily by chisq_loc to get results as
481 * python code by Lodo
482 */
483 for (j = 0; j < poly_degree; j++) {
484 gsl_matrix_set (covar, j, j, gsl_matrix_get (covar, j, j)*chisq_loc);
485 coeffs_val[j] = gsl_vector_get (p, j);
486 coeffs_err[j] = sqrt(gsl_matrix_get (covar, j, j));
487 }
488
489 /*
490 cpl_msg_debug(cpl_func,"GSL poly fit c0: %16.8f c1: %16.8f, c2: %16.8f",
491 coeffs_val[0], coeffs_val[1], coeffs_val[2]);
492 cpl_msg_debug(cpl_func,"GSL poly fit e0: %16.8f e1: %16.8f, e2: %16.8f",
493 coeffs_err[0], coeffs_err[1], coeffs_err[2]);
494 cpl_msg_debug(cpl_func,"GSL poly fit chi2: %16.8f", *chisq);
495 */
496
497 gsl_vector_free(x);
498 gsl_vector_free(y);
499 gsl_vector_free(w);
500 gsl_vector_free(p);
501 gsl_matrix_free(X);
502
503
504 return covar;
505}
506
507/*----------------------------------------------------------------------------*/
520static cpl_vector*
521hdrl_mode_vector_trim(const cpl_vector* vec, const double min,
522 const double max) {
523 cpl_size size = cpl_vector_get_size(vec);
524 cpl_error_ensure(size > 0, CPL_ERROR_ILLEGAL_INPUT,
525 return NULL, "vector size must be > 0");
526 cpl_vector * vec_trim = cpl_vector_new(size);
527 const double * pvec = cpl_vector_get_data_const(vec);
528 double * pvec_trim = cpl_vector_get_data(vec_trim);
529
530 cpl_size index = 0;
531
532 for (cpl_size i = 0; i < size; i++) {
533 if (pvec[i] >= min && pvec[i] <= max) {
534 pvec_trim[index] = pvec[i];
535 index++;
536 } else {
537 continue;
538 }
539 }
540
541 if (index > 0) {
542 /* return resized vector */
543 cpl_vector_set_size(vec_trim, index);
544 return vec_trim;
545 } else {
546 /* return NULL */
547 cpl_vector_delete(vec_trim);
548 return NULL;
549 }
550}
551
552/*----------------------------------------------------------------------------*/
570/*----------------------------------------------------------------------------*/
571static cpl_error_code hdrl_mode_median(
572 const cpl_vector * vec,
573 const double histo_min,
574 const double histo_max,
575 const cpl_size nbin,
576 const cpl_size error_niter,
577 double * mode_median,
578 double * mode_median_err)
579{
580
581
582 cpl_error_ensure(vec != NULL, CPL_ERROR_NULL_INPUT,
583 return CPL_ERROR_NULL_INPUT, "Null input vector data");
584
585 gsl_histogram * histogram = hdrl_mode_histogram(vec, histo_min, histo_max,
586 nbin);
587
588 cpl_error_ensure(histogram != NULL, CPL_ERROR_NULL_INPUT,
589 return CPL_ERROR_NULL_INPUT, "Histogram can not be created");
590
591 /* for debug purposes
592 cpl_table* t = hdrl_mode_histogram_to_table(h, histo_min, nbin);
593 cpl_table_save(t,NULL,NULL,"debug.fits",CPL_IO_CREATE);
594 */
595
596 cpl_size histogram_bin_max = gsl_histogram_max_bin(histogram);
597 //cpl_msg_debug(cpl_func, "histogram bin max: %lld", histogram_bin_max);
598
599 /* get median and error associated to the values in the histogram max bin */
600
601 /* to compute the mode we need to determine the median of the values in the
602 * histogram bin max. The error is the standard deviation of these values.
603 * */
604 double lower = 0.;
605 double upper = 0.;
606 gsl_histogram_get_range(histogram, histogram_bin_max, &lower, &upper);
607
608 cpl_vector * vec_final = hdrl_mode_vector_trim(vec, lower, upper);
609 *mode_median = cpl_vector_get_median(vec_final);
610
611 if(error_niter == 0) {
612 *mode_median_err = cpl_vector_get_stdev(vec_final);
613 cpl_msg_debug(cpl_func, "(method median) computed mode: %g, associated error: %g",
614 *mode_median, *mode_median_err);
615 } else {
616 /* Derive only the analytical error, if requested */
617 *mode_median_err = 0;
618 }
619
620 /* free memory */
621 gsl_histogram_free(histogram);
622 cpl_vector_delete(vec_final);
623
624 return cpl_error_get_code();
625}
626/*----------------------------------------------------------------------------*/
646/*----------------------------------------------------------------------------*/
647static cpl_error_code hdrl_mode_weight(
648 const cpl_vector * vec,
649 const double histo_min,
650 const double histo_max,
651 const double bin_size,
652 const cpl_size nbin,
653 const cpl_size error_niter,
654 double * mode_weight,
655 double * mode_weight_err)
656{
657
658 cpl_error_ensure(vec != NULL, CPL_ERROR_NULL_INPUT,
659 return CPL_ERROR_NULL_INPUT, "Null input vector data");
660
661 gsl_histogram * histogram = hdrl_mode_histogram(vec, histo_min, histo_max,
662 nbin);
663
664 cpl_error_ensure(histogram != NULL, CPL_ERROR_NULL_INPUT,
665 return CPL_ERROR_NULL_INPUT, "Histogram can not be created");
666
667 cpl_table* table = hdrl_mode_histogram_to_table(histogram, histo_min,
668 bin_size, nbin);
669 //cpl_table_save(table, NULL, NULL, "histo_tab.fits", CPL_IO_DEFAULT);
670
671 double histogram_max = gsl_histogram_max_val(histogram);
672 //cpl_msg_debug(cpl_func, "histogram value at max: %g", histogram_max);
673 cpl_size histogram_bin_max = gsl_histogram_max_bin(histogram);
674 //cpl_msg_debug(cpl_func, "histogram bin max: %lld", histogram_bin_max);
675 //cpl_msg_debug(cpl_func, "histogram bin size: %ld", gsl_histogram_bins(histogram));
676 if(histogram_bin_max > 0 &&
677 (histogram_bin_max < ((cpl_size)gsl_histogram_bins(histogram) - 1))) {
678 cpl_msg_debug(cpl_func, "histogram (bin_max-1) value: %16.8g",
679 gsl_histogram_get(histogram,histogram_bin_max-1));
680 cpl_msg_debug(cpl_func, "histogram (bin_max+1) value: %16.8g",
681 gsl_histogram_get(histogram,histogram_bin_max+1));
682 }
683
684 /* determine the mean of the values that lead to the histogram bin max */
685 double lower = 0;
686 double upper = 0;
687 gsl_histogram_get_range(histogram, histogram_bin_max, &lower, &upper);
688 //cpl_msg_debug(cpl_func, "histogram range lower: %16.8g upper: %16.8g", lower, upper);
689 cpl_table_and_selected_double(table,"COUNTS",CPL_EQUAL_TO, histogram_max);
690 cpl_table* extract = cpl_table_extract_selected(table);
691
692 //cpl_table_save(extract, NULL, NULL, "histo_tab_extr1.fits", CPL_IO_DEFAULT);
693
694 double mean = cpl_table_get_column_mean(extract, "INTERVAL_LOWER");
695
696 cpl_table_delete(extract);
697
698 cpl_size max_pos = 0;
699 cpl_table_get_column_maxpos(table, "INTERVAL_LOWER", &max_pos);
700
701 cpl_table_delete(table);
702
703 double freq1 = histogram_max;
704 double level1 = mean;
705
706 double freq0, freq2;
707 //cpl_msg_debug(cpl_func, "nbin: %lld histogram_bin_max: %lld", nbin,histogram_bin_max);
708 if(histogram_bin_max + 1 <= ((cpl_size)nbin - 1)) {
709 freq2 = gsl_histogram_get(histogram, histogram_bin_max + 1);
710 } else {
711 freq2 = 0;
712 }
713
714
715 if(histogram_bin_max > 0) {
716 freq0 = gsl_histogram_get(histogram, histogram_bin_max - 1);
717 } else {
718 freq0 = 0;
719 }
720
721 /*
722 cpl_msg_debug(cpl_func, "Frequencies: f0: %16.8g f1: %16.10g f2: %16.8g l1: %16.10g",
723 freq0, freq1, freq2, level1);
724 */
725
726 double diff1 = freq1 - freq0;
727 double diff2 = freq1 - freq2;
728 double factor = diff1 / (diff1 + diff2);
729 if( factor == 0 || isnan(factor) ) {
730 factor = 0.5;
731 }
732 /* compute mode */
733 *mode_weight = level1 + bin_size * factor;
734
735 /* compute error associated to mode */
736 if (error_niter == 0) {
737 double dd1 = sqrt(freq0 + freq1);
738 double dd2 = sqrt(freq1 + freq2);
739 double dw2=1;
740
741 double term1 = diff2 * dd1 / pow((diff1 + diff2), 2);
742 double term2 = diff1 * dd2 / pow((diff1 + diff2), 2);
743 dw2 = pow(term1, 2.) + pow(term2, 2.);
744 *mode_weight_err = bin_size * sqrt(dw2);
745 /*
746 cpl_msg_debug(cpl_func,"fct: %16.8g m: %16.8g, dd1: %16.8g dd2: %16.8g",
747 factor,*mode_weight,dd1,dd2);
748 cpl_msg_debug(cpl_func,"dw2: %16.8g e: %16.8g", dw2, *mode_weight_err);
749 */
750
751 } else {
752 /* Derive only the analytical error, if requested */
753 *mode_weight_err = 0;
754 }
755
756 cpl_msg_debug(cpl_func, "(method weight) computed mode: %16.10g error: %16.10g",
757 *mode_weight, *mode_weight_err);
758
759 /* free memory */
760 gsl_histogram_free(histogram);
761
762 return cpl_error_get_code();
763}
772static double
773hdrl_mode_fit_analytical_error(const double* coeffs_val,
774 const double* coeffs_err, gsl_matrix* covar,
775 const double scale_factor){
776
777
778 double a2 = coeffs_val[2];
779 double a1 = coeffs_val[1];
780 //double a0 = coeffs_val[0];
781
782 /*NOTE we scale arbirtarily by scale_factor = chi2 / dof, where dof is the
783 * degree of freedom to get results as python code by Lodo
784 */
785 double cvar_matrix = gsl_matrix_get (covar, 2, 1) * scale_factor;
786
787 double da2 = coeffs_err[2];
788 double da1 = coeffs_err[1];
789 //double da0 = coeffs_err[0];
790
791 double cvr_term = 2 * (-1. / (2 * a2)) * (a1 / (2 * a2 * a2)) * cvar_matrix;
792
793 double add1 = (da1 / (2 * a2));
794 add1 *= add1;
795 double a_square = a2 * a2;
796 double add2 = (a1 * da2 / (2 * a_square ));
797 add2 *= add2;
798 double err2 = add1 + add2;
799
800 err2 += cvr_term;
801
802 return sqrt(err2);;
803}
804
805/*----------------------------------------------------------------------------*/
824/*----------------------------------------------------------------------------*/
825static cpl_error_code hdrl_mode_fit(
826 const cpl_vector * vec,
827 const double histo_min,
828 const double histo_max,
829 const double bin_size,
830 const cpl_size nbin,
831 const cpl_size error_niter,
832 double * mode_fit,
833 double * mode_fit_err)
834{
835
836 cpl_error_ensure(vec != NULL, CPL_ERROR_NULL_INPUT,
837 return CPL_ERROR_NULL_INPUT, "Null input vector data");
838
839 /* We use semi_fit_window = 2 as the fit polynomial degree used is 2
840 * We cannot use a larger value to prevent asymmetries (for asymmetric
841 * distributions) that cannot be fit well by a parabola */
842 cpl_size semi_fit_window = 2;
843
844 double hmin = histo_min;
845 double hmax = histo_max;
846
847 gsl_histogram * h = hdrl_mode_histogram(vec, hmin, hmax, nbin);
848
849 cpl_error_ensure(h != NULL, CPL_ERROR_NULL_INPUT,
850 return CPL_ERROR_NULL_INPUT, "Histogram can not be created");
851
852 cpl_size histogram_bin_max = (cpl_size) gsl_histogram_max_bin(h);
853 cpl_size histo_nbins = (cpl_size) gsl_histogram_bins(h);
854
855 if(histogram_bin_max > 0) {
856 cpl_msg_debug(cpl_func, "histogram (bin_max-1) value: %16.8g",
857 gsl_histogram_get(h, histogram_bin_max - 1));
858 }
859 if(histogram_bin_max < (histo_nbins-1)) {
860 cpl_msg_debug(cpl_func, "histogram (bin_max+1) value: %16.8g",
861 gsl_histogram_get(h, histogram_bin_max + 1));
862 }
863
864 double lower = 0.;
865 double upper = 0.;
866 gsl_histogram_get_range(h, gsl_histogram_max_bin(h), &lower, &upper);
867 double value_at_max = lower;
868
869 /* check if one has enough points (at least 3) to be able to do a polynomial
870 * fit. If too less points, switch to weight mode
871 */
872 if( histo_nbins < 3 ) {
873 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
874 cpl_msg_error(cpl_func,"Cannot do polynomial fit with less than 3 points.");
875 gsl_histogram_free(h);
876 return cpl_error_get_code();
877 } else {
878 cpl_size mn = histogram_bin_max - semi_fit_window;
879 cpl_size mx = histogram_bin_max + semi_fit_window;
880 mn = histogram_bin_max - semi_fit_window;
881 mn = (mn > 0) ? mn : 0;
882 mx = histogram_bin_max + semi_fit_window;
883 mx = (mx < histo_nbins) ? mx : histo_nbins-1; // Cannot go till histo_nbins else after for loop seg fault
884
885 /* fit a parabola */
886 const cpl_size good_points = (mx - mn + 1);
887 const cpl_size win_size = 2 * semi_fit_window + 1;
888 const cpl_size sz = (good_points < win_size) ? good_points : win_size;
889
890 double* x_vals = cpl_calloc(sz, sizeof(double));
891 double* y_vals = cpl_calloc(sz, sizeof(double));
892 double* y_errs = cpl_calloc(sz, sizeof(double));
893 cpl_size j = 0;
894 /* we assume uniform error as we do not know how to compute the error
895 * associated to each bin size
896 */
897 for(cpl_size i = mn; i <= mx; i++) {
898 double lower_local = 0.;
899 double upper_local = 0.;
900 gsl_histogram_get_range(h, i, &lower_local, &upper_local);
901 x_vals[j] = lower_local;
902 y_vals[j] = gsl_histogram_get(h, i);
903 y_errs[j] = 1;
904 j++;
905 }
906
907 /* do the polynomial fit, using GSL */
908 cpl_size deg = 2;
909 double chi2;
910 double* coeffs_val = (double *) cpl_calloc (sz, sizeof(double));
911 double* coeffs_err = (double *) cpl_calloc (sz, sizeof(double));
912 double* fit = (double *) cpl_calloc (sz, sizeof(double));
913 gsl_matrix* covar =
914 hdrl_gsl_fit_poly (x_vals, y_vals, y_errs, sz, deg+1, fit, coeffs_val,
915 coeffs_err, &chi2);
916
917 double m = -coeffs_val[1] / 2. / coeffs_val[2];
918 double value_at_mode = gsl_poly_eval(coeffs_val, sz, m);
919
920 *mode_fit = m + bin_size / 2.;
921
922 /* check that we are not at an edge and that we found a true maximum
923 * (that means the point distribution had a proper shape)
924 */
925 double y_min_val = gsl_poly_eval(coeffs_val, sz,x_vals[0]);
926 double y_max_val = gsl_poly_eval(coeffs_val, sz,x_vals[sz - 1]);
927 double max_parab = (y_max_val > y_min_val ) ? y_max_val : y_min_val;
928 cpl_boolean is_not_supported_fit = CPL_FALSE;
929 if( (fabs(value_at_max - m) > bin_size / 2.) ||
930 (value_at_mode < max_parab ) ) {
931 /* if we are near an edge point or the point distribution is not
932 * proper it is safer to use the weight method to get the mode
933 */
934 if( fabs(value_at_max - m) > bin_size / 2. ) {
935 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
936 cpl_msg_error(cpl_func, "Max too close to point distribution edge: abs(value_at_max+bin_size/2-m) > bin_size");
937 is_not_supported_fit = CPL_TRUE;
938 }
939
940 if(value_at_mode < max_parab) {
941 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
942 cpl_msg_error(cpl_func, "Value at mode is not a valid maximum: value_at_mode < np.max(parab)");
943 is_not_supported_fit = CPL_TRUE;
944 }
945 if(is_not_supported_fit == CPL_TRUE) {
946 gsl_matrix_free(covar);
947 gsl_histogram_free(h);
948 cpl_free(fit);
949 cpl_free(coeffs_val);
950 cpl_free(coeffs_err);
951 cpl_free(y_errs);
952 cpl_free(x_vals);
953 cpl_free(y_vals);
954 return cpl_error_get_code();
955 }
956
957 } else {
958 if (error_niter == 0) {
959 /* As requested, compute the analytical error, */
960 double dof = sz - (deg +1);
961 double scale_factor = (chi2 / dof);
962
963 *mode_fit_err =
964 hdrl_mode_fit_analytical_error(coeffs_val, coeffs_err, covar,
965 scale_factor);
966 } else {
967 /* As requested, do not compute the analytical error */
968 *mode_fit_err = 0;
969 }
970
971 }
972 /* check if problem occurred during poly fit or in error computation*/
973 if(!isfinite(*mode_fit_err) || !isfinite(*mode_fit)) {
974 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT);
975 *mode_fit_err = NAN;
976 *mode_fit = NAN;
977 }
978 cpl_msg_debug(cpl_func, "(method fit) computed mode: %16.10g err: %16.10g ",
979 *mode_fit, *mode_fit_err);
980
981 gsl_matrix_free(covar);
982 cpl_free(fit);
983 cpl_free(coeffs_val);
984 cpl_free(coeffs_err);
985 cpl_free(x_vals);
986 cpl_free(y_vals);
987 cpl_free(y_errs);
988
989 }
990
991 /* free memory */
992 gsl_histogram_free(h);
993
994 return cpl_error_get_code();
995}
996
997/*----------------------------------------------------------------------------*/
1072/*----------------------------------------------------------------------------*/
1073cpl_error_code hdrl_mode_clip_image(
1074 const cpl_image * source,
1075 const double histo_min,
1076 const double histo_max,
1077 const double bin_size,
1078 const hdrl_mode_type method,
1079 const cpl_size error_niter,
1080 double * mode,
1081 double * mode_error,
1082 cpl_size * naccepted)
1083{
1084 cpl_vector * vec_source = NULL;
1085
1086 /* Check Entries */
1087 cpl_error_ensure(source != NULL, CPL_ERROR_NULL_INPUT,
1088 return CPL_ERROR_NULL_INPUT, "Null input source image!");
1089
1090 cpl_image* data =(cpl_image*) source;
1091
1092 /* compress images to vectors excluding the bad pixels */
1093 vec_source = hdrl_image_to_vector(data, cpl_image_get_bpm_const(source));
1094
1095 if (vec_source != NULL) {
1096 /* compute mode */
1097 hdrl_mode_clip(vec_source, histo_min, histo_max, bin_size, method,
1098 error_niter, mode, mode_error, naccepted);
1099
1100 if(error_niter > 0) {
1101 /* Calculate the error using the bootstrap technique */
1102 hdrl_mode_bootstrap(vec_source, histo_min, histo_max, bin_size,
1103 method, error_niter, mode_error);
1104 }
1105
1106 if (CPL_ERROR_NONE != cpl_error_get_code()) {
1107 /* if error occurred exit after cleaning memory */
1108 cpl_vector_delete(vec_source);
1109 return cpl_error_get_code();
1110 }
1111
1112 } else {
1113 /* no good pixels */
1114 *mode = NAN;
1115 *mode_error = NAN;
1116 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
1117 }
1118
1119 cpl_vector_delete(vec_source);
1120 return cpl_error_get_code();
1121}
1122
1123
1124/*----------------------------------------------------------------------------*/
1143/*----------------------------------------------------------------------------*/
1144
1145cpl_error_code hdrl_mode_clip(
1146 cpl_vector * vec,
1147 const double histo_min,
1148 const double histo_max,
1149 const double bin_size,
1150 const hdrl_mode_type method,
1151 const cpl_size error_niter,
1152 double * mode,
1153 double * mode_error,
1154 cpl_size * naccepted)
1155{
1156 cpl_error_ensure(vec != NULL, CPL_ERROR_NULL_INPUT,
1157 return CPL_ERROR_NULL_INPUT, "Null input source image!");
1158
1159 *naccepted = 0; /* In case an error happens - if not, re-set it later */
1160
1161 cpl_vector * loc_vec = NULL;
1162 double bsize = bin_size;
1163
1164 /* if binsize <=0 : the routine derives the value from a formula */
1165 if (bsize <= DBL_EPSILON) {
1166 bsize = hdrl_mode_compute_binsize(vec);
1167 }
1168
1169 double hmin = histo_min;
1170 double hmax = histo_max;
1171 cpl_size nbin = 0;
1172 cpl_vector* trim_vec = NULL;
1173 /* if histo_min >= histo_max: the routine derives the two values from the
1174 * data (and also the number of histogram bins, nbin)
1175 */
1176 if ( histo_min >= histo_max) {
1177
1178 trim_vec = cpl_vector_duplicate(vec);
1179 hmin = cpl_vector_get_min(vec) - bsize / 2;
1180 hmax = cpl_vector_get_max(vec) + bsize / 2;
1181 /* the following two are required to have same results as Lodo's code
1182 * as implemented in the unit tests
1183 */
1184 nbin = hdrl_mode_get_nbin(hmin, hmax, bsize);
1185 hmax = hmin + (nbin * bsize);
1186 /* Fallback solution for the case there is only one single value */
1187 if(hmin == hmax) {
1188 hmin = nextafter(hmin, hmin - FLT_EPSILON);
1189 hmax = nextafter(hmax, hmax + FLT_EPSILON);
1190 bsize = nextafter(0.0, 1.0);
1191 nbin = 1 ;
1192 }
1193
1194
1195 } else {
1196 nbin = hdrl_mode_get_nbin(hmin, hmax, bsize);
1197 trim_vec = hdrl_mode_vector_trim(vec, hmin, hmax);
1198 if( hmin + (nbin * bsize) >= hmax) {
1199 /* adjust hmax to use gsl_histogram without changing bsize */
1200 //cpl_msg_warning(cpl_func,"Strange");
1201 hmax = hmin + (nbin * bsize);
1202 }
1203 }
1204 cpl_msg_debug(cpl_func,"Histogram bin size: %g min: %g max: %g number of bins: %lld",
1205 bsize, hmin, hmax, nbin);
1206
1207 cpl_error_ensure(trim_vec != NULL, CPL_ERROR_NULL_INPUT,
1208 return CPL_ERROR_NULL_INPUT, "No data for mode "
1209 "computation. Try to change mode parameters ... ");
1210
1211 loc_vec = trim_vec;
1212
1213
1214 /* method == enum with mode_fit, mode_weight, mode_median */
1215 switch(method) {
1216
1217 case HDRL_MODE_FIT:
1218
1219 if (CPL_ERROR_NONE != hdrl_mode_fit(loc_vec, hmin, hmax, bsize, nbin, error_niter, mode,
1220 mode_error) ) {
1221 //cpl_error_set(cpl_func,CPL_ERROR_ILLEGAL_INPUT);
1222 cpl_msg_error(cpl_func,"Mode computation failed using method fit. Try method weight or median.");
1223 }
1224 break;
1225
1226 case HDRL_MODE_WEIGHTED:
1227
1228 if (CPL_ERROR_NONE != hdrl_mode_weight(loc_vec, hmin, hmax, bsize, nbin,
1229 error_niter, mode, mode_error) ) {
1230 //cpl_error_set(cpl_func,CPL_ERROR_UNSUPPORTED_MODE);
1231 cpl_msg_error(cpl_func,"Mode computation failed using method weight. Try method fit or median.");
1232
1233 }
1234 break;
1235
1236 case HDRL_MODE_MEDIAN:
1237
1238 if (CPL_ERROR_NONE != hdrl_mode_median(loc_vec, hmin, hmax, nbin,
1239 error_niter, mode, mode_error) ) {
1240 //cpl_error_set(cpl_func,CPL_ERROR_UNSUPPORTED_MODE);
1241 cpl_msg_error(cpl_func,"Mode computation failed using method median. Try method fit or weight");
1242
1243 }
1244
1245 break;
1246
1247 default:
1248 cpl_msg_error(cpl_func,"Unsupported mode method. Supported methods are: fit, weight, median");
1249 return CPL_ERROR_UNSUPPORTED_MODE;
1250
1251 }
1252
1253 *naccepted = cpl_vector_get_size(vec);
1254 cpl_vector_delete(trim_vec);
1255
1256 return cpl_error_get_code();
1257}
1258
1259/* ---------------------------------------------------------------------------*/
1278/* ---------------------------------------------------------------------------*/
1279
1280cpl_error_code
1281hdrl_mode_bootstrap(
1282 const cpl_vector * vec,
1283 const double histo_min,
1284 const double histo_max,
1285 const double bin_size,
1286 const hdrl_mode_type method,
1287 const cpl_size error_niter,
1288 double * mode_error)
1289{
1290
1291
1292 /* To ensure new random numbers each time the compiled program is called,
1293 * call srand(), e.g. srand(time(NULL)) */
1294
1295 /* Predefine the state for the threads as they get manipulated in the thread*/
1296 hdrl_random_state ** pstate = cpl_calloc(omp_get_max_threads(),
1297 sizeof(hdrl_random_state *));
1298
1299 for(cpl_size i = 0 ; i < (cpl_size)omp_get_max_threads(); i++){
1300 uint64_t seed[2] = {(uint64_t)rand(), (uint64_t)rand()};
1301 pstate[i] = hdrl_random_state_new(1, seed);
1302 }
1303
1304 cpl_size vec_size = cpl_vector_get_size(vec);
1305 const double * pvec = cpl_vector_get_data_const(vec);
1306 cpl_image * ima_mode = cpl_image_new(1, error_niter, CPL_TYPE_DOUBLE);
1307 double * pima_mode = cpl_image_get_data_double(ima_mode);
1308
1309 /* In case the loop gets parallelized use cpl_image_get_bpm outside the loop
1310 * once to get a single bad pixel mask for the image */
1311 cpl_binary * pima_mode_bpm = cpl_mask_get_data(cpl_image_get_bpm(ima_mode));
1312
1313 HDRL_OMP(omp parallel for)
1314 for (cpl_size i = 0; i < error_niter; i++) {
1315 cpl_error_code err = CPL_ERROR_NONE;
1316 cpl_vector * vec_simul = cpl_vector_new(vec_size);
1317 double mode = 0.;
1318 double mode_err = 0.;
1319 cpl_size naccepted = 0;
1320 double * pvec_simul = cpl_vector_get_data(vec_simul);
1321
1322 /* simulate the new vector by bootstraping */
1323 for (cpl_size j = 0; j < vec_size; j++) {
1324 pvec_simul[j] = pvec[(cpl_size)hdrl_random_uniform_int64(pstate[omp_get_thread_num()], 0, vec_size - 1)];
1325 }
1326
1327 err = hdrl_mode_clip(vec_simul, histo_min, histo_max, bin_size, method,
1328 -1, &mode, &mode_err, &naccepted);
1329 cpl_vector_delete(vec_simul);
1330
1331 if (err == CPL_ERROR_NONE) {
1332 pima_mode[i] = mode;
1333 pima_mode_bpm[i] = CPL_BINARY_0; // may be overdoing ....
1334 } else {
1335 pima_mode[i] = NAN;
1336 pima_mode_bpm[i] = CPL_BINARY_1;
1337 cpl_error_reset();
1338 }
1339 }
1340
1341 /* Return the mad scaled standard deviation */
1342 //double mad = 0;
1343 //cpl_image_get_mad(ima_mode, &mad);
1344 //*mode_error = mad * CPL_MATH_STD_MAD;
1345
1346 /* Return the pure standard deviation as the mode on the various simulations
1347 * is not very smoothly distributed - the MAD based mode depends thus more on
1348 * the number of iterations */
1349 *mode_error = cpl_image_get_stdev(ima_mode);
1350
1351 /* clean memory */
1352 cpl_image_delete(ima_mode);
1353
1354 for(cpl_size i = 0 ; i < (cpl_size)omp_get_max_threads(); i++){
1355 hdrl_random_state_delete(pstate[i]);
1356 }
1357 cpl_free(pstate);
1358
1359 return cpl_error_get_code();
1360}
1361
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
double hdrl_collapse_mode_parameter_get_histo_max(const hdrl_parameter *p)
get high value
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_mode(const hdrl_parameter *self)
check if parameter is a mode parameter
hdrl_mode_type hdrl_collapse_mode_parameter_get_method(const hdrl_parameter *p)
get the mode determination method
char * hdrl_join_string(const char *sep_, int n,...)
join strings together
Definition: hdrl_utils.c:812