35#include "hdrl_sigclip.h"
36#include "hdrl_utils.h"
37#include "hdrl_collapse.h"
38#include "hdrl_random.h"
41#define omp_get_max_threads() 1
42#define omp_get_thread_num() 0
52#include <gsl/gsl_multifit.h>
53#include <gsl/gsl_vector.h>
54#include <gsl/gsl_histogram.h>
55#include <gsl/gsl_poly.h>
57static const char * method_to_string(
const hdrl_mode_type method)
60 case HDRL_MODE_MEDIAN:
63 case HDRL_MODE_WEIGHTED:
70 cpl_error_set_message(cpl_func, CPL_ERROR_ILLEGAL_INPUT,
71 "mode method unknown");
96cpl_parameterlist * hdrl_mode_parameter_create_parlist(
97 const char *base_context,
99 const hdrl_parameter *defaults)
101 cpl_ensure(base_context && prefix && defaults,
102 CPL_ERROR_NULL_INPUT, NULL);
105 CPL_ERROR_INCOMPATIBLE_INPUT, NULL);
107 cpl_parameterlist * parlist = cpl_parameterlist_new();
110 hdrl_setup_vparameter(parlist, prefix,
".",
"",
111 "histo-min", base_context,
112 "Minimum pixel value to accept for mode computation",
117 hdrl_setup_vparameter(parlist, prefix,
".",
"",
118 "histo-max", base_context,
119 "Maximum pixel value to accept for mode computation",
124 hdrl_setup_vparameter(parlist, prefix,
".",
"",
125 "bin-size", base_context,
126 "Binsize of the histogram",
135 const char * method_def = method_to_string(method);
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");
143 cpl_parameter_set_alias(par, CPL_PARAMETER_MODE_CLI, name);
144 cpl_parameter_disable(par, CPL_PARAMETER_MODE_ENV);
146 cpl_parameterlist_append(parlist, par);
158 hdrl_setup_vparameter(parlist, prefix,
".",
"",
159 "error-niter", base_context,
160 "Iterations to compute the mode error",
164 if (cpl_error_get_code()) {
165 cpl_parameterlist_delete(parlist);
193cpl_error_code hdrl_mode_parameter_parse_parlist(
194 const cpl_parameterlist * parlist,
199 hdrl_mode_type * method,
200 cpl_size * error_niter)
202 cpl_ensure_code(prefix && parlist, CPL_ERROR_NULL_INPUT);
207 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
208 *histo_min = cpl_parameter_get_double(par);
213 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
214 *histo_max = cpl_parameter_get_double(par);
220 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
221 *bin_size = cpl_parameter_get_double(par);
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) {
233 return cpl_error_set_message(cpl_func, CPL_ERROR_DATA_NOT_FOUND,
234 "Parameter mode.method not found");
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;
258 const cpl_parameter * par = cpl_parameterlist_find_const(parlist, name);
259 *error_niter = cpl_parameter_get_int(par);
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);
269 return CPL_ERROR_NONE;
285hdrl_mode_compute_binsize(cpl_vector* vec)
288 cpl_size size = cpl_vector_get_size(vec);
290 hcpl_vector_get_mad_window(vec, 1, size, &mad);
291 double std_deviation_from_mad = mad * CPL_MATH_STD_MAD;
297 double binsize = 2. * 3.49 * std_deviation_from_mad / pow(size, 1. / 3. );
300 binsize = nextafter(0,1.0);
314hdrl_mode_get_nbin(
const double min,
const double max,
const double size)
317 return (cpl_size)floor( (max - min) / size ) + 1;
333static gsl_histogram *
334hdrl_mode_histogram(
const cpl_vector* vec,
const double histogram_min,
335 const double histogram_max,
const cpl_size nbin)
338 cpl_error_ensure(nbin > 0, CPL_ERROR_ILLEGAL_INPUT,
339 return NULL,
"Number of bins must be > 0");
341 cpl_error_ensure(histogram_max > histogram_min, CPL_ERROR_ILLEGAL_INPUT,
342 return NULL,
"histo_max must be larger than histo_min");
344 gsl_histogram * histogram = gsl_histogram_alloc(nbin);
346 gsl_histogram_set_ranges_uniform (histogram, histogram_min, histogram_max);
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++) {
352 gsl_histogram_increment(histogram, data[i]);
380hdrl_mode_histogram_to_table(gsl_histogram* histogram,
381 const double histogram_min,
382 const double histogram_step,
const cpl_size nbin){
384 cpl_table* table = cpl_table_new(nbin);
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);
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);
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");
402 for(cpl_size i = 0; i < nbin; i++) {
404 phl[i] = histogram_min + phb[i] * histogram_step;
405 phu[i] = phl[i] + histogram_step;
406 phc[i] = histogram->bin[i];
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,
442 double err, chisq_loc;
443 gsl_vector *x = NULL, *y = NULL, *w = NULL, *p = NULL;
444 gsl_matrix *X = NULL, *covar = NULL;
447 gsl_multifit_linear_workspace *work = NULL;
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);
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]);
461 gsl_vector_set (w, i, 1.0 / err / err);
463 for (j = 0; j < poly_degree; j++)
464 gsl_matrix_set (X, i, j, gsl_pow_int(gsl_vector_get(x, i), j));
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);
471 for (i = 0; i < (cpl_size) n_sampling_points; i++) {
474 for (j = 0; j < poly_degree; j++)
475 fit[i] += gsl_matrix_get (X, i, j) * gsl_vector_get (p, j);
478 *chisq = chisq_loc/(n_sampling_points - poly_degree);
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));
521hdrl_mode_vector_trim(
const cpl_vector* vec,
const double min,
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);
532 for (cpl_size i = 0; i < size; i++) {
533 if (pvec[i] >= min && pvec[i] <= max) {
534 pvec_trim[index] = pvec[i];
543 cpl_vector_set_size(vec_trim, index);
547 cpl_vector_delete(vec_trim);
571static cpl_error_code hdrl_mode_median(
572 const cpl_vector * vec,
573 const double histo_min,
574 const double histo_max,
576 const cpl_size error_niter,
577 double * mode_median,
578 double * mode_median_err)
582 cpl_error_ensure(vec != NULL, CPL_ERROR_NULL_INPUT,
583 return CPL_ERROR_NULL_INPUT,
"Null input vector data");
585 gsl_histogram * histogram = hdrl_mode_histogram(vec, histo_min, histo_max,
588 cpl_error_ensure(histogram != NULL, CPL_ERROR_NULL_INPUT,
589 return CPL_ERROR_NULL_INPUT,
"Histogram can not be created");
596 cpl_size histogram_bin_max = gsl_histogram_max_bin(histogram);
606 gsl_histogram_get_range(histogram, histogram_bin_max, &lower, &upper);
608 cpl_vector * vec_final = hdrl_mode_vector_trim(vec, lower, upper);
609 *mode_median = cpl_vector_get_median(vec_final);
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);
617 *mode_median_err = 0;
621 gsl_histogram_free(histogram);
622 cpl_vector_delete(vec_final);
624 return cpl_error_get_code();
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,
653 const cpl_size error_niter,
654 double * mode_weight,
655 double * mode_weight_err)
658 cpl_error_ensure(vec != NULL, CPL_ERROR_NULL_INPUT,
659 return CPL_ERROR_NULL_INPUT,
"Null input vector data");
661 gsl_histogram * histogram = hdrl_mode_histogram(vec, histo_min, histo_max,
664 cpl_error_ensure(histogram != NULL, CPL_ERROR_NULL_INPUT,
665 return CPL_ERROR_NULL_INPUT,
"Histogram can not be created");
667 cpl_table* table = hdrl_mode_histogram_to_table(histogram, histo_min,
671 double histogram_max = gsl_histogram_max_val(histogram);
673 cpl_size histogram_bin_max = gsl_histogram_max_bin(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));
687 gsl_histogram_get_range(histogram, histogram_bin_max, &lower, &upper);
689 cpl_table_and_selected_double(table,
"COUNTS",CPL_EQUAL_TO, histogram_max);
690 cpl_table* extract = cpl_table_extract_selected(table);
694 double mean = cpl_table_get_column_mean(extract,
"INTERVAL_LOWER");
696 cpl_table_delete(extract);
698 cpl_size max_pos = 0;
699 cpl_table_get_column_maxpos(table,
"INTERVAL_LOWER", &max_pos);
701 cpl_table_delete(table);
703 double freq1 = histogram_max;
704 double level1 = mean;
708 if(histogram_bin_max + 1 <= ((cpl_size)nbin - 1)) {
709 freq2 = gsl_histogram_get(histogram, histogram_bin_max + 1);
715 if(histogram_bin_max > 0) {
716 freq0 = gsl_histogram_get(histogram, histogram_bin_max - 1);
726 double diff1 = freq1 - freq0;
727 double diff2 = freq1 - freq2;
728 double factor = diff1 / (diff1 + diff2);
729 if( factor == 0 || isnan(factor) ) {
733 *mode_weight = level1 + bin_size * factor;
736 if (error_niter == 0) {
737 double dd1 = sqrt(freq0 + freq1);
738 double dd2 = sqrt(freq1 + freq2);
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);
753 *mode_weight_err = 0;
756 cpl_msg_debug(cpl_func,
"(method weight) computed mode: %16.10g error: %16.10g",
757 *mode_weight, *mode_weight_err);
760 gsl_histogram_free(histogram);
762 return cpl_error_get_code();
773hdrl_mode_fit_analytical_error(
const double* coeffs_val,
774 const double* coeffs_err, gsl_matrix* covar,
775 const double scale_factor){
778 double a2 = coeffs_val[2];
779 double a1 = coeffs_val[1];
785 double cvar_matrix = gsl_matrix_get (covar, 2, 1) * scale_factor;
787 double da2 = coeffs_err[2];
788 double da1 = coeffs_err[1];
791 double cvr_term = 2 * (-1. / (2 * a2)) * (a1 / (2 * a2 * a2)) * cvar_matrix;
793 double add1 = (da1 / (2 * a2));
795 double a_square = a2 * a2;
796 double add2 = (a1 * da2 / (2 * a_square ));
798 double err2 = add1 + add2;
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,
831 const cpl_size error_niter,
833 double * mode_fit_err)
836 cpl_error_ensure(vec != NULL, CPL_ERROR_NULL_INPUT,
837 return CPL_ERROR_NULL_INPUT,
"Null input vector data");
842 cpl_size semi_fit_window = 2;
844 double hmin = histo_min;
845 double hmax = histo_max;
847 gsl_histogram * h = hdrl_mode_histogram(vec, hmin, hmax, nbin);
849 cpl_error_ensure(h != NULL, CPL_ERROR_NULL_INPUT,
850 return CPL_ERROR_NULL_INPUT,
"Histogram can not be created");
852 cpl_size histogram_bin_max = (cpl_size) gsl_histogram_max_bin(h);
853 cpl_size histo_nbins = (cpl_size) gsl_histogram_bins(h);
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));
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));
866 gsl_histogram_get_range(h, gsl_histogram_max_bin(h), &lower, &upper);
867 double value_at_max = lower;
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();
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;
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;
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));
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);
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));
914 hdrl_gsl_fit_poly (x_vals, y_vals, y_errs, sz, deg+1, fit, coeffs_val,
917 double m = -coeffs_val[1] / 2. / coeffs_val[2];
918 double value_at_mode = gsl_poly_eval(coeffs_val, sz, m);
920 *mode_fit = m + bin_size / 2.;
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 ) ) {
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;
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;
945 if(is_not_supported_fit == CPL_TRUE) {
946 gsl_matrix_free(covar);
947 gsl_histogram_free(h);
949 cpl_free(coeffs_val);
950 cpl_free(coeffs_err);
954 return cpl_error_get_code();
958 if (error_niter == 0) {
960 double dof = sz - (deg +1);
961 double scale_factor = (chi2 / dof);
964 hdrl_mode_fit_analytical_error(coeffs_val, coeffs_err, covar,
973 if(!isfinite(*mode_fit_err) || !isfinite(*mode_fit)) {
974 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_OUTPUT);
978 cpl_msg_debug(cpl_func,
"(method fit) computed mode: %16.10g err: %16.10g ",
979 *mode_fit, *mode_fit_err);
981 gsl_matrix_free(covar);
983 cpl_free(coeffs_val);
984 cpl_free(coeffs_err);
992 gsl_histogram_free(h);
994 return cpl_error_get_code();
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,
1081 double * mode_error,
1082 cpl_size * naccepted)
1084 cpl_vector * vec_source = NULL;
1087 cpl_error_ensure(source != NULL, CPL_ERROR_NULL_INPUT,
1088 return CPL_ERROR_NULL_INPUT,
"Null input source image!");
1090 cpl_image* data =(cpl_image*) source;
1093 vec_source = hdrl_image_to_vector(data, cpl_image_get_bpm_const(source));
1095 if (vec_source != NULL) {
1097 hdrl_mode_clip(vec_source, histo_min, histo_max, bin_size, method,
1098 error_niter, mode, mode_error, naccepted);
1100 if(error_niter > 0) {
1102 hdrl_mode_bootstrap(vec_source, histo_min, histo_max, bin_size,
1103 method, error_niter, mode_error);
1106 if (CPL_ERROR_NONE != cpl_error_get_code()) {
1108 cpl_vector_delete(vec_source);
1109 return cpl_error_get_code();
1116 cpl_error_set(cpl_func, CPL_ERROR_ILLEGAL_INPUT);
1119 cpl_vector_delete(vec_source);
1120 return cpl_error_get_code();
1145cpl_error_code hdrl_mode_clip(
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,
1153 double * mode_error,
1154 cpl_size * naccepted)
1156 cpl_error_ensure(vec != NULL, CPL_ERROR_NULL_INPUT,
1157 return CPL_ERROR_NULL_INPUT,
"Null input source image!");
1161 cpl_vector * loc_vec = NULL;
1162 double bsize = bin_size;
1165 if (bsize <= DBL_EPSILON) {
1166 bsize = hdrl_mode_compute_binsize(vec);
1169 double hmin = histo_min;
1170 double hmax = histo_max;
1172 cpl_vector* trim_vec = NULL;
1176 if ( histo_min >= histo_max) {
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;
1184 nbin = hdrl_mode_get_nbin(hmin, hmax, bsize);
1185 hmax = hmin + (nbin * bsize);
1188 hmin = nextafter(hmin, hmin - FLT_EPSILON);
1189 hmax = nextafter(hmax, hmax + FLT_EPSILON);
1190 bsize = nextafter(0.0, 1.0);
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) {
1201 hmax = hmin + (nbin * bsize);
1204 cpl_msg_debug(cpl_func,
"Histogram bin size: %g min: %g max: %g number of bins: %lld",
1205 bsize, hmin, hmax, nbin);
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 ... ");
1219 if (CPL_ERROR_NONE != hdrl_mode_fit(loc_vec, hmin, hmax, bsize, nbin, error_niter, mode,
1222 cpl_msg_error(cpl_func,
"Mode computation failed using method fit. Try method weight or median.");
1226 case HDRL_MODE_WEIGHTED:
1228 if (CPL_ERROR_NONE != hdrl_mode_weight(loc_vec, hmin, hmax, bsize, nbin,
1229 error_niter, mode, mode_error) ) {
1231 cpl_msg_error(cpl_func,
"Mode computation failed using method weight. Try method fit or median.");
1236 case HDRL_MODE_MEDIAN:
1238 if (CPL_ERROR_NONE != hdrl_mode_median(loc_vec, hmin, hmax, nbin,
1239 error_niter, mode, mode_error) ) {
1241 cpl_msg_error(cpl_func,
"Mode computation failed using method median. Try method fit or weight");
1248 cpl_msg_error(cpl_func,
"Unsupported mode method. Supported methods are: fit, weight, median");
1249 return CPL_ERROR_UNSUPPORTED_MODE;
1253 *naccepted = cpl_vector_get_size(vec);
1254 cpl_vector_delete(trim_vec);
1256 return cpl_error_get_code();
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)
1296 hdrl_random_state ** pstate = cpl_calloc(omp_get_max_threads(),
1297 sizeof(hdrl_random_state *));
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);
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);
1311 cpl_binary * pima_mode_bpm = cpl_mask_get_data(cpl_image_get_bpm(ima_mode));
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);
1318 double mode_err = 0.;
1319 cpl_size naccepted = 0;
1320 double * pvec_simul = cpl_vector_get_data(vec_simul);
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)];
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);
1331 if (err == CPL_ERROR_NONE) {
1332 pima_mode[i] = mode;
1333 pima_mode_bpm[i] = CPL_BINARY_0;
1336 pima_mode_bpm[i] = CPL_BINARY_1;
1349 *mode_error = cpl_image_get_stdev(ima_mode);
1352 cpl_image_delete(ima_mode);
1354 for(cpl_size i = 0 ; i < (cpl_size)omp_get_max_threads(); i++){
1355 hdrl_random_state_delete(pstate[i]);
1359 return cpl_error_get_code();
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