-
Notifications
You must be signed in to change notification settings - Fork 1
/
fit.c
1347 lines (1261 loc) · 55.4 KB
/
fit.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
Jaakko's Backscattering Simulator (JaBS)
Copyright (C) 2021 - 2023 Jaakko Julin
This program is free software; you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation; either version 2 of the License, or
(at your option) any later version.
See LICENSE.txt for the full license.
Some parts of this source file under different license, see below!
*/
/*
* This file is based partially based around the GSL examples for nonlinear least-squares fitting. See e.g.
* https://www.gnu.org/software/gsl/doc/html/nls.html . Please note that GSL (https://www.gnu.org/software/gsl/) is
* is distributed under the terms of the GNU General Public License (GPL).
*/
#ifndef _GNU_SOURCE
#define _GNU_SOURCE
#endif
#include <stdio.h>
#include "win_compat.h"
#include <string.h>
#include <assert.h>
#include <gsl/gsl_matrix.h>
#include <gsl/gsl_blas.h>
#include "jabs_debug.h"
#include "defaults.h"
#include "generic.h"
#include "jabs.h"
#include "spectrum.h"
#include "message.h"
#include "gsl_inline.h"
#include "histogram.h"
#include "fit.h"
#ifdef _OPENMP
#include <omp.h>
#endif
int fit_detector(const jibal *jibal, const fit_data_det *fdd, const simulation *sim, result_spectra *spectra, gsl_vector *f) {
detector *det = sim_det(sim, fdd->i_det);
if(!det) {
jabs_message(MSG_ERROR, "No detector set.\n");
return EXIT_FAILURE;
}
detector_update(det); /* non-const pointer inside const... not great */
sim_workspace *ws = sim_workspace_init(jibal, sim, det);
if(!ws) {
jabs_message(MSG_ERROR, "Could not initialize workspace of detector %s.\n", det->name);
return EXIT_FAILURE;
}
if(simulate_with_ds(ws)) {
jabs_message(MSG_ERROR, "Simulation of detector %s spectrum failed!\n", det->name);
return EXIT_FAILURE;
}
fit_data_det_residual_vector_set(fdd, ws->histo_sum, f);
#if 0
if(iter_call == 1) { /* TODO: give spectra (and f_iter) as pointers, set these if pointer given */
gsl_vector_memcpy(fdd->f_iter, f); /* Copy residual vector, this can be used on later calls if fdd is not active */
}
#endif
if(spectra) {
result_spectra_free(spectra);
fit_data_spectra_copy_to_spectra_from_ws(spectra, det, fdd->exp, ws);
}
sim_workspace_free(ws);
return EXIT_SUCCESS;
}
void fit_deriv_prepare_jspaces(const gsl_vector *x, fit_data *fit) {
jacobian_space *space = fit->jspace;
for(size_t j = 0; j < fit->fit_params->n_active; j++) { /* Make copies of sim (shallow) with deep copies of detector and sample. This way we can parallelize. */
struct jacobian_space *spc = &space[j];
spc->n_spectra_calculated = 0;
fit_variable *var = spc->var;
double xj = gsl_vector_get(x, j);
double delta = fit->h_df * fabs(xj);
if(delta == 0.0) {
delta = fit->h_df; /* TODO: this is what GSL does, but it doesn't always work */
}
spc->delta_inv = 1.0/delta;
*(var->value) = (xj + delta) * var->value_orig; /* Perturb */
spc->sim = *fit->sim; /* Shallow copy! */
switch(var->type) { /* Deepen the copy based on the variable */
case FIT_VARIABLE_SAMPLE:
spc->sim.sample = sample_from_sample_model(fit->sm);
break;
case FIT_VARIABLE_DETECTOR:
spc->sim.det = spc->det;
spc->sim.det[var->i_det] = detector_clone(fit->sim->det[var->i_det]);
break;
default:
break;
}
*(var->value) = (xj) * var->value_orig; /* Unperturb */
}
}
void fit_deriv_cleanup_jspaces(fit_data *fit) {
jacobian_space *space = fit->jspace;
for(size_t j = 0; j < fit->fit_params->n_active; j++) {
struct jacobian_space *spc = &space[j];
fit_variable *var = spc->var;
fit->stats.n_spectra_iter += spc->n_spectra_calculated;
switch(var->type) {
case FIT_VARIABLE_SAMPLE:
sample_free(spc->sim.sample);
break;
case FIT_VARIABLE_DETECTOR:
detector_free(spc->sim.det[var->i_det]);
break;
default:
break;
}
}
}
int fit_deriv_function(const gsl_vector *x, void *params, gsl_matrix *J) {
struct fit_data *fit = (struct fit_data *) params;
assert(fit->stats.iter_call >= 1); /* This function relies on data that was computed when iter_call == 1 */
assert(fit->fdf->p == fit->fit_params->n_active);
gsl_matrix_set_zero(J); /* We might only set parts of the Jacobian, the rest of the elements of the matrix should be zero. */
fit_deriv_prepare_jspaces(x, fit);
double start = jabs_clock();
volatile int error = FALSE;
const int n = (int) fit->fit_params->n_active;
int j;
#pragma omp parallel default(none) shared(fit, J, error, n)
#pragma omp for schedule(dynamic)
for(j = 0; j < n; j++) {
//fprintf(stderr, "Thread id %i got %zu.\n", omp_get_thread_num(), j);
struct jacobian_space *spc = &fit->jspace[j];
fit_variable *var = spc->var;
size_t i_start, i_stop;
if(var->type == FIT_VARIABLE_DETECTOR) {
i_start = var->i_det;
i_stop = var->i_det;
} else { /* All detectors */
i_start = 0;
i_stop = fit->sim->n_det - 1;
}
DEBUGMSG("Variable %s (i_v = %zu, j = %i) Jacobian. i_det = [%zu, %zu]", var->name, var->i_v, j, i_start, i_stop);
for(size_t i_det = i_start; i_det <= i_stop; i_det++) {
fit_data_det *fdd = &fit->fdd[i_det];
gsl_vector_view v = gsl_vector_subvector(spc->f_param, fdd->f_offset, fdd->n_ch); /* The ROIs of this detector are a subvector of f_param (all channels in fit) */
DEBUGMSG("Detector %zu (FDD %p), %zu ranges, offset: %zu, len: %zu", i_det, (void *) fdd, fdd->n_ranges, fdd->f_offset, fdd->n_ch);
if(fdd->n_ranges == 0) {
continue;
}
if(fit_detector(fit->jibal, fdd, &spc->sim, NULL, &v.vector)) {
error = TRUE;
break;
}
spc->n_spectra_calculated++;
for(size_t i = 0; i < fdd->n_ch; i++) {
double fnext = jabs_gsl_vector_get(&v.vector, i);
double fi = jabs_gsl_vector_get(fit->f_iter, fdd->f_offset + i);
jabs_gsl_matrix_set(J, i + fdd->f_offset, j, (fnext - fi) * spc->delta_inv);
}
}
if(error) {
continue;
}
}
fit_deriv_cleanup_jspaces(fit);
double end = jabs_clock();
fit->stats.cputime_iter += (end - start);
return error ? GSL_FAILURE : GSL_SUCCESS;
}
int fit_function(const gsl_vector *x, void *params, gsl_vector *f) {
struct fit_data *fit = (struct fit_data *) params;
fit->stats.iter_call++;
if(fit->stats.iter_call == 1) { /* Clear stored spectra from previous iteration */
for(size_t i = 0; i < fit->n_det_spectra; i++) {
result_spectra_free(&fit->spectra[i]);
}
}
if(fit_parameters_set_from_vector(fit, x)) { /* This also sets some helper numbers and arrays */
return GSL_FAILURE;
}
DEBUGMSG("Fit iteration %zu call %zu. Size of vector x: %zu, f: %zu. %zu active fit parameters.",
fit->stats.iter, fit->stats.iter_call, x->size, f->size, fit->fit_params->n_active);
if(sim_sanity_check(fit->sim) == GSL_FAILURE) {
fit->stats.error = FIT_ERROR_SANITY;
return GSL_FAILURE;
}
sample_free(fit->sim->sample); /* TODO: this is only necessary on first call of iter and when sample specific parameters are active */
fit->sim->sample = sample_from_sample_model(fit->sm);
if(!fit->sim->sample) {
fit->stats.error = FIT_ERROR_SANITY;
return GSL_FAILURE;
}
double start = jabs_clock();
volatile int error = FALSE;
if(fit->sim->n_det == 1) { /* Skip OpenMP if only one workspace */
fit_data_det *fdd = &fit->fdd[0];
result_spectra *spectra = (fit->stats.iter_call == 1 ? &fit->spectra[fdd->i_det] : NULL);
error = fit_detector(fit->jibal, fdd, fit->sim, spectra, f);
} else {
int i;
const int n = (int) fit->sim->n_det;
#pragma omp parallel default(none) shared(fit, error, f, n)
#pragma omp for
for(i = 0; i < n; i++) {
#ifdef _OPENMP
DEBUGMSG("Thread id %i got %i.", omp_get_thread_num(), i);
#endif
fit_data_det *fdd = &fit->fdd[i];
if(fdd->n_ranges == 0) { /* This will NOT simulate those detectors that don't participate in fit! */
continue;
}
result_spectra *spectra = (fit->stats.iter_call == 1 ? &fit->spectra[fdd->i_det] : NULL); /* Pass spectra pointer on first call */
gsl_vector_view f_det = gsl_vector_subvector(f, fdd->f_offset, fdd->n_ch);
if(fit_detector(fit->jibal, fdd, fit->sim, spectra, &f_det.vector)) {
error = TRUE;
}
}
}
double end = jabs_clock();
DEBUGMSG("Fit iteration %zu call %zu simulation done.", fit->stats.iter, fit->stats.iter_call);
if(error) {
jabs_message(MSG_ERROR, "Error in simulating (during fit).");
return EXIT_FAILURE;
}
fit->stats.cputime_iter += (end - start);
fit->stats.n_evals_iter++;
fit->stats.n_spectra_iter += fit->sim->n_det;
if(fit->stats.iter_call == 1) {
gsl_vector_memcpy(fit->f_iter, f); /* Store residual vector on first call, this acts as baseline for everything we do later */
}
//gsl_vector_memcpy(f, fit->f);
DEBUGSTR("Successful fit function call.");
return GSL_SUCCESS;
}
int fit_sanity_check(const fit_data *fit) {
if(sample_model_sanity_check(fit->sm)) {
return GSL_FAILURE;
}
if(sim_sanity_check(fit->sim)) {
return GSL_FAILURE;
}
return GSL_SUCCESS;
}
int fit_parameters_set_from_vector(struct fit_data *fit, const gsl_vector *x) {
DEBUGMSG("Fit iteration has %zu active parameters from a total of %zu possible.", fit->fit_params->n_active, fit->fit_params->n)
for(size_t i = 0; i < fit->fit_params->n; i++) {
fit_variable *var = &fit->fit_params->vars[i];
if(!var->active) {
continue;
}
if(!isfinite(*(var->value))) {
DEBUGMSG("Fit iteration %zu, call %zu, fit variable %zu (%s) is not finite.", fit->stats.iter, fit->stats.iter_call, var->i_v, var->name);
fit->stats.error = FIT_ERROR_SANITY;
return GSL_FAILURE;
}
*(var->value) = gsl_vector_get(x, var->i_v) * var->value_orig;
}
return GSL_SUCCESS;
}
void fit_iter_stats_update(struct fit_data *fit_data, const gsl_multifit_nlinear_workspace *w) {
gsl_vector *f = gsl_multifit_nlinear_residual(w);
/* compute reciprocal condition number of J(x) */
gsl_multifit_nlinear_rcond(&fit_data->stats.rcond, w);
gsl_blas_ddot(f, f, &fit_data->stats.chisq);
fit_data->stats.norm = gsl_blas_dnrm2(f);
fit_data->stats.chisq_dof = fit_data->stats.chisq / fit_data->dof;
fit_data->stats.n_evals += fit_data->stats.n_evals_iter;
fit_data->stats.n_spectra += fit_data->stats.n_spectra_iter;
fit_data->stats.cputime_cumul += fit_data->stats.cputime_iter;
}
void fit_iter_stats_print(const struct fit_stats *stats) {
jabs_message(MSG_INFO, "%4zu | %10.2e | %14.7e | %12.7lf | %7zu | %7zu | %10.3lf | %13.1lf |\n",
stats->iter, 1.0 / stats->rcond, stats->norm,
stats->chisq_dof, stats->n_evals, stats->n_spectra,
stats->cputime_cumul, 1000.0 * stats->cputime_iter / stats->n_spectra_iter);
}
void fit_stats_print(const struct fit_stats *stats, jabs_msg_level msg_level) {
if(stats->chisq_dof > 0.0) {
jabs_message(msg_level, "Final chisq/dof = %.7lf\n", stats->chisq_dof);
}
}
int fit_data_fit_range_add(struct fit_data *fit_data, const struct roi *range) { /* Makes a deep copy */
if(range->low == 0 && range->high == 0) {
return EXIT_FAILURE;
}
if(range->high < range->low) {
return EXIT_FAILURE;
}
fit_data->n_fit_ranges++;
fit_data->fit_ranges = realloc(fit_data->fit_ranges, fit_data->n_fit_ranges * sizeof(roi));
if(!fit_data->fit_ranges) {
fit_data->n_fit_ranges = 0;
return EXIT_FAILURE;
}
fit_data->fit_ranges[fit_data->n_fit_ranges - 1] = *range;
qsort(fit_data->fit_ranges, fit_data->n_fit_ranges, sizeof(roi), fit_range_compare); /* Sorting by detector is important for fit residual vector to work */
return EXIT_SUCCESS;
}
void fit_data_fit_ranges_free(struct fit_data *fit_data) {
if(!fit_data)
return;
free(fit_data->fit_ranges);
fit_data->fit_ranges = NULL;
fit_data->n_fit_ranges = 0;
}
void fit_data_det_residual_vector_set(const fit_data_det *fdd, const jabs_histogram *histo_sum, gsl_vector *f) { /* Sets residual vector (f) based on histo_sum */
size_t i_vec = 0;
if(fdd->n_ranges == 0) {
fprintf(stderr, "No ranges!\n");
}
for(size_t i_range = 0; i_range < fdd->n_ranges; i_range++) {
const roi *range = &fdd->ranges[i_range];
for(size_t i = range->low; i <= range->high; i++) {
if(i >= histo_sum->n) { /* Outside range of simulated spectrum (simulated is zero) */
jabs_gsl_vector_set(f, i_vec, fdd->exp->bin[i]);
} else {
jabs_gsl_vector_set(f, i_vec, fdd->exp->bin[i] - histo_sum->bin[i]);
}
i_vec++;
}
}
DEBUGMSG("Set %zu elements of vector from %zu ranges.", i_vec, fdd->n_ranges);
}
fit_data *fit_data_new(const jibal *jibal, simulation *sim) {
struct fit_data *fit = calloc(1, sizeof(struct fit_data));
fit->jibal = jibal;
fit->sim = sim;
fit_data_exp_alloc(fit);
fit_data_defaults(fit);
return fit;
}
void fit_data_defaults(fit_data *f) {
f->n_iters_max = FIT_ITERS_MAX;
f->xtol = FIT_XTOL;
f->chisq_tol = FIT_CHISQ_TOL;
f->chisq_fast_tol = FIT_FAST_CHISQ_TOL;
f->phase_start = FIT_PHASE_FAST;
f->phase_stop = FIT_PHASE_SLOW;
}
int fit_data_jspace_init(fit_data *fit, size_t n_channels_in_fit) {
if(!fit || !fit->fit_params || !fit->fit_params->n_active) {
return EXIT_FAILURE;
}
fit->jspace = calloc(fit->fit_params->n_active, sizeof(jacobian_space)); /* Space for each parameter for Jacobian calculation */
for(size_t j = 0; j < fit->fit_params->n_active; j++) {
jacobian_space *spc = &fit->jspace[j];
spc->var = fit_params_find_active(fit->fit_params, j);
assert(spc->var);
if(spc->var->type == FIT_VARIABLE_DETECTOR) {
spc->det = calloc(fit->sim->n_det, sizeof(detector *)); /* Array for detector pointers. This overrides array in sim, so it needs to be an array, although we only use one element of it. */
} else {
spc->det = NULL;
}
spc->f_param = gsl_vector_alloc(n_channels_in_fit); /* Residuals for this parameter */
}
return EXIT_SUCCESS;
}
void fit_data_jspace_free(fit_data *fit) {
if(!fit) {
return;
}
for(size_t j = 0; j < fit->fit_params->n_active; j++) {
jacobian_space *spc = &fit->jspace[j];
free(spc->det);
gsl_vector_free(spc->f_param);
}
free(fit->jspace);
fit->jspace = NULL;
}
void fit_data_free(fit_data *fit) {
if(!fit)
return;
fit_data_reset(fit);
fit_data_exp_free(fit);
fit_data_fdd_free(fit);
for(size_t i = 0; i < fit->n_det_spectra; i++) {
result_spectra_free(&fit->spectra[i]);
}
free(fit->spectra);
free(fit);
}
void fit_data_reset(fit_data *fit) {
if(!fit) {
return;
}
fit_data_fit_ranges_free(fit);
fit_params_free(fit->fit_params);
fit->fit_params = NULL;
if(fit->covar) {
gsl_matrix_free(fit->covar);
fit->covar = NULL;
}
fit_data_spectra_free(fit);
}
void fit_data_exp_reset(fit_data *fit) {
for(size_t i = 0; i < fit->n_exp; i++) {
jabs_histogram_free(fit->exp[i]);
fit->exp[i] = NULL;
}
}
void fit_data_roi_print(const struct fit_data *fit_data, const struct roi *roi) {
if(!fit_data) {
return;
}
jabs_histogram *histo_sim = fit_data_histo_sum(fit_data, roi->i_det);
jabs_histogram *histo_exp = fit_data_exp(fit_data, roi->i_det);
jabs_histogram *histo_ref = fit_data->ref;
size_t n_sim = jabs_histogram_channels_in_range(histo_sim, roi->low, roi->high);
size_t n_exp = jabs_histogram_channels_in_range(histo_exp, roi->low, roi->high);
size_t n_ref = jabs_histogram_channels_in_range(histo_ref, roi->low, roi->high);
double sim_cts = jabs_histogram_roi(histo_sim, roi->low, roi->high);
double exp_cts = jabs_histogram_roi(histo_exp, roi->low, roi->high);
double ref_cts = jabs_histogram_roi(histo_ref, roi->low, roi->high);
const detector *det = sim_det(fit_data->sim, roi->i_det);
jabs_message(MSG_INFO, " low = %12zu\n", roi->low);
jabs_message(MSG_INFO, " high = %12zu\n", roi->high);
if(det) {
jabs_message(MSG_INFO, " E_low = %12.3lf keV (low energy edge of bin)\n", detector_calibrated(det, JIBAL_ANY_Z, roi->low) / C_KEV);
jabs_message(MSG_INFO, " E_high = %12.3lf keV (high energy edge of bin)\n", detector_calibrated(det, JIBAL_ANY_Z, roi->high + 1) / C_KEV);
}
if(histo_sim) {
jabs_message(MSG_INFO, " n_sim = %12zu (number of channels)\n", n_sim);
jabs_message(MSG_INFO, " sim = %12.8g (number of counts)\n", sim_cts);
}
if(histo_exp) {
jabs_message(MSG_INFO, " n_exp = %12zu (number of channels)\n", n_exp);
jabs_message(MSG_INFO, " exp = %12.8g (number of counts)\n", exp_cts);
jabs_message(MSG_INFO, " sqrt(exp) = %12.5lf\n", sqrt(exp_cts));
}
if(histo_ref) {
jabs_message(MSG_INFO, " n_ref = %12zu (number of channels)\n", n_ref);
jabs_message(MSG_INFO, " ref = %12.8g (number of counts)\n", ref_cts);
jabs_message(MSG_INFO, " sqrt(ref) = %12.5lf\n", sqrt(ref_cts));
}
if(histo_sim && histo_exp && sim_cts > 0 && exp_cts > 0) {
jabs_message(MSG_INFO, " exp-sim = %12.8g\n", exp_cts - sim_cts);
jabs_message(MSG_INFO, " sim/exp = %12.5lf\n", sim_cts / exp_cts);
jabs_message(MSG_INFO, " exp/sim = %12.5lf\n", exp_cts / sim_cts);
jabs_message(MSG_INFO, " 1/sqrt(exp) = %12.5lf%%\n", 100.0 / sqrt(exp_cts));
jabs_message(MSG_INFO, "(exp-sim)/exp = %12.5lf%%\n", 100.0 * (exp_cts - sim_cts) / exp_cts);
}
if(histo_sim && histo_ref && sim_cts > 0 && ref_cts > 0) {
jabs_message(MSG_INFO, " ref-sim = %12g\n", ref_cts - sim_cts);
jabs_message(MSG_INFO, " sim/ref = %12.5lf\n", sim_cts / ref_cts);
jabs_message(MSG_INFO, " ref/sim = %12.5lf\n", ref_cts / sim_cts);
jabs_message(MSG_INFO, " 1/sqrt(ref) = %12.5lf%%\n", 100.0 / sqrt(ref_cts));
jabs_message(MSG_INFO, "(ref-sim)/ref = %12.5lf%%\n", 100.0 * (ref_cts - sim_cts) / ref_cts);
}
}
jabs_histogram *fit_data_exp(const fit_data *fit, size_t i_det) {
if(!fit || !fit->exp)
return NULL;
if(i_det >= fit->sim->n_det)
return NULL;
assert(fit->n_exp == fit->sim->n_det);
return fit->exp[i_det];
}
fit_params *fit_params_all(fit_data *fit) {
simulation *sim = fit->sim;
sample_model *sm = fit->sm;
if(!sim)
return NULL;
size_t param_name_max_len = 256; /* Laziness. We use a fixed size temporary string. snprintf is used, so no overflows should occur, but very long names may be truncated. */
char *param_name = malloc(sizeof(char) * param_name_max_len);
fit_params *params = fit_params_new();
fit_params_add_parameter(params, FIT_VARIABLE_BEAM, &sim->fluence, "fluence", "", 1.0, 0); /* This must be the first parameter always, as there is a speedup in the fit routine */
fit_params_add_parameter(params, FIT_VARIABLE_GEOMETRY, &sim->sample_theta, "alpha", "deg", C_DEG, 0);
fit_params_add_parameter(params, FIT_VARIABLE_BEAM, &sim->beam_E, "energy", "keV", C_KEV, 0);
if(sm) {
for(size_t i_range = 0; i_range < sm->n_ranges; i_range++) {
sample_range *r = &(sm->ranges[i_range]);
size_t range_index = i_range + 1; /* Human readable indexing */
if(r->x < 0.0) { /* Don't fit anything if thickness is zero (negative shouldn't be possible) */
continue;
}
if(r->rough.model != ROUGHNESS_FILE) { /* Layer thickness is not a parameter with arbitrary roughness from a file */
snprintf(param_name, param_name_max_len, "thick%zu", range_index);
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, &(r->x), param_name, "tfu", C_TFU, 0);
}
snprintf(param_name, param_name_max_len, "yield%zu", range_index);
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, &(r->yield), param_name, "", 1.0, 0);
snprintf(param_name, param_name_max_len, "yield_slope%zu", range_index);
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, &(r->yield_slope), param_name, "", 1.0, 0);
if(i_range == sm->n_ranges - 1) { /* Last range, add channeling "aliases" (=yield corrections) */
snprintf(param_name, param_name_max_len, "channeling");
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, &(r->yield), param_name, "", 1.0, 0);
snprintf(param_name, param_name_max_len, "channeling_slope");
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, &(r->yield_slope), param_name, "", 1.0, 0);
}
snprintf(param_name, param_name_max_len, "bragg%zu", range_index);
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, &(r->bragg), param_name, "", 1.0, 0);
snprintf(param_name, param_name_max_len, "stragg%zu", range_index);
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, &(r->stragg), param_name, "", 1.0, 0);
if(r->rough.model == ROUGHNESS_GAMMA && r->rough.x > 0.0) {
snprintf(param_name, param_name_max_len, "rough%zu", range_index);
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, &(r->rough.x), param_name, "tfu", C_TFU, 0);
}
for(size_t i_mat = 0; i_mat < sm->n_materials; i_mat++) {
if(*sample_model_conc_bin(sm, i_range, i_mat) < CONC_TOLERANCE) /* Don't add fit variables for negative, zero or very low concentrations. */
continue;
snprintf(param_name, param_name_max_len, "conc%zu_%s", range_index, sm->materials[i_mat]->name);
fit_params_add_parameter(params, FIT_VARIABLE_SAMPLE, sample_model_conc_bin(sm, i_range, i_mat), param_name, "%", C_PERCENT, 0);
}
}
}
for(size_t i_det = 0; i_det < sim->n_det; i_det++) {
detector *det = sim_det(sim, i_det);
char *det_name = NULL;
if(asprintf(&det_name, "det%zu_", i_det + 1) < 0) {
return NULL;
}
snprintf(param_name, param_name_max_len, "%ssolid", det_name);
fit_params_add_parameter(params, FIT_VARIABLE_DETECTOR, &det->solid, param_name, "msr", C_MSR, i_det);
snprintf(param_name, param_name_max_len, "%stheta", det_name);
fit_params_add_parameter(params, FIT_VARIABLE_DETECTOR, &det->theta, param_name, "deg", C_DEG, i_det);
snprintf(param_name, param_name_max_len, "%sphi", det_name);
fit_params_add_parameter(params, FIT_VARIABLE_DETECTOR, &det->phi, param_name, "deg", C_DEG, i_det);
for(int Z = JIBAL_ANY_Z; Z <= det->cal_Z_max; Z++) {
calibration *c = detector_get_calibration(det, Z);
if(Z != JIBAL_ANY_Z && c == det->calibration) /* No Z-specific calibration */
continue;
assert(c);
size_t n = calibration_get_number_of_params(c);
for(int i = CALIBRATION_PARAM_RESOLUTION_SLOPE; i < (int) n; i++) {
if(i < 0 && calibration_get_param(c, i) == 0.0) {
DEBUGMSG("Fit variable %i omitted, since it is zero.\n", i);
continue;
}
char *calib_param_name = calibration_param_name(c->type, i);
snprintf(param_name, param_name_max_len, "%scalib%s%s_%s",
det_name,
(Z == JIBAL_ANY_Z) ? "" : "_",
(Z == JIBAL_ANY_Z) ? "" : jibal_element_name(fit->jibal->elements, Z),
calib_param_name);
free(calib_param_name);
fit_params_add_parameter(params, FIT_VARIABLE_DETECTOR, calibration_get_param_ref(c, i), param_name, "keV", C_KEV, i_det);
}
}
free(det_name);
}
free(param_name);
return params;
}
int fit_data_fdd_init(fit_data *fit) {
if(!fit) {
return EXIT_FAILURE;
}
if(fit->fdd) {
free(fit->fdd);
}
fit->fdd = calloc(fit->sim->n_det, sizeof(fit_data_det));
for(size_t i_roi = 0; i_roi < fit->n_fit_ranges; i_roi++) { /* Loop over rois, fdds, count rois and n of channels */
const roi *roi = &fit->fit_ranges[i_roi];
for(size_t i_det = 0; i_det < fit->sim->n_det; i_det++) {
fit_data_det *fdd = &fit->fdd[i_det];
if(i_det != roi->i_det) { /* This ROI is not for this fdd */
continue;
}
fdd->n_ranges++;
fdd->n_ch += (roi->high - roi->low + 1);
}
}
for(size_t i_det = 0; i_det < fit->sim->n_det; i_det++) {
fit_data_det *fdd = &fit->fdd[i_det];
fdd->exp = fit_data_exp(fit, i_det);
fdd->i_det = i_det;
fdd->det = sim_det(fit->sim, i_det);
fdd->f_iter = gsl_vector_alloc(fdd->n_ch);
fdd->ranges = calloc(fdd->n_ranges, sizeof(roi));
if(fdd->n_ranges == 0) {
DEBUGMSG("Detector %zu no ranges, length: %zu", i_det + 1, fdd->n_ch);
}
size_t i = 0; /* Index of fdd roi */
size_t i_vec = 0; /* fit residual vector index at beginning of roi */
for(size_t i_roi = 0; i_roi < fit->n_fit_ranges; i_roi++) {
const roi *roi = &fit->fit_ranges[i_roi];
size_t l = (roi->high - roi->low + 1);
if(i_det != roi->i_det) {
i_vec += l;
continue;
}
if(i == 0) {
fdd->f_offset = i_vec;
DEBUGMSG("FDD %zu, %zu ROIs, First ROI is %zu/%zu of all, subvector offset %zu, length %zu)", i_det, fdd->n_ranges, i_roi, fit->n_fit_ranges, fdd->f_offset, fdd->n_ch);
}
fdd->ranges[i] = *roi;
i++;
i_vec += l;
}
}
return EXIT_SUCCESS;
}
void fit_data_fdd_free(fit_data *fit) {
if(!fit || !fit->fdd) {
return;
}
for(size_t i_fdd = 0; i_fdd < fit->sim->n_det; i_fdd++) {
fit_data_det *fdd = &fit->fdd[i_fdd];
if(!fdd) {
continue;
}
gsl_vector_free(fdd->f_iter);
free(fdd->ranges);
}
free(fit->fdd);
fit->fdd = NULL;
}
void fit_data_exp_alloc(fit_data *fit) {
if(!fit) {
return;
}
size_t n_alloc = fit->sim->n_det;
if(fit->n_exp == 0) { /* This could be handled by realloc too, but this is an easy way to get null pointers as contents. */
fit->exp = calloc(n_alloc, sizeof(jabs_histogram *));
} else {
fit->exp = realloc(fit->exp, sizeof(jabs_histogram *) * n_alloc);
for(size_t i = fit->n_exp; i < n_alloc; i++) { /* Reset newly allocated space */
fit->exp[i] = NULL;
}
}
fit->n_exp = n_alloc;
}
void fit_data_exp_free(fit_data *fit) {
if(!fit->exp)
return;
fit_data_exp_reset(fit);
free(fit->exp);
fit->exp = NULL;
fit->n_exp = 0;
}
int fit_data_load_exp(struct fit_data *fit, size_t i_det, const char *filename) {
jabs_histogram *h = spectrum_read_detector(filename, sim_det(fit->sim, i_det));
if(!h) {
jabs_message(MSG_ERROR, "Reading spectrum from file \"%s\" was not successful.\n", filename);
return EXIT_FAILURE;
}
if(fit->exp[i_det]) {
jabs_histogram_free(fit->exp[i_det]);
}
fit->exp[i_det] = h;
return EXIT_SUCCESS;
}
int fit_data_set_sample_model(fit_data *fit, sample_model *sm_new) {
int sanity = sample_model_sanity_check(sm_new);
if(sanity) {
return sanity;
}
sample_model_free(fit->sm);
fit->sm = sm_new;
if(fit->sim->n_reactions > 0) {
jabs_message(MSG_WARNING, "Reactions were reset automatically, since the sample was changed.\n");
sim_reactions_free(fit->sim);
}
return EXIT_SUCCESS;
}
jabs_histogram *fit_data_histo_sum(const fit_data *fit, size_t i_det) {
if(!fit)
return NULL;
if(i_det >= fit->n_det_spectra)
return NULL;
if(RESULT_SPECTRA_SIMULATED > fit->spectra[i_det].n_spectra)
return NULL;
return result_spectra_simulated_histo(&fit->spectra[i_det]);
}
void fit_data_spectra_copy_to_spectra_from_ws(result_spectra *spectra, const detector *det, const jabs_histogram *exp, const sim_workspace *ws) {
spectra->n_spectra = ws->n_reactions + RESULT_SPECTRA_N_FIXED; /* Simulated, experimental, confidence limits (2) + reaction spectra */
spectra->s = calloc(spectra->n_spectra, sizeof(result_spectrum));
result_spectrum_set(&spectra->s[RESULT_SPECTRA_SIMULATED], ws->histo_sum, "Simulated", NULL, REACTION_NONE);
result_spectrum_set(&spectra->s[RESULT_SPECTRA_EXPERIMENTAL], exp, "Experimental", NULL, REACTION_NONE);
calibration_apply_to_histogram(det->calibration, result_spectra_experimental_histo(spectra));
DEBUGMSG("Copying %zu spectra (%zu reactions) from ws to spectra structure.", spectra->n_spectra, ws->n_reactions);
for(size_t i = 0; i < ws->n_reactions; i++) {
const reaction *r = ws->reactions[i]->r;
result_spectrum *s = &spectra->s[RESULT_SPECTRA_REACTION_SPECTRUM(i)];
s->histo = jabs_histogram_clone(ws->reactions[i]->histo);
if(asprintf(&s->name, "%s (%s)", r->target->name, reaction_type_to_string(r->type)) < 0) {
s->name = NULL;
}
s->target_isotope = r->target;
s->type = r->type;
}
}
int fit_data_spectra_alloc(fit_data *fit) {
fit_data_spectra_free(fit);
fit->spectra = calloc(fit->sim->n_det, sizeof(result_spectra));
if(!fit->spectra) {
fit->n_det_spectra = 0;
return EXIT_FAILURE;
}
fit->n_det_spectra = fit->sim->n_det;
return EXIT_SUCCESS;
}
void fit_data_spectra_free(fit_data *fit) {
if(!fit) {
return;
}
for(size_t i = 0; i < fit->n_det_spectra; i++) {
result_spectra_free(&fit->spectra[i]);
}
free(fit->spectra);
fit->spectra = NULL;
fit->n_det_spectra = 0;
}
int fit_data_add_det(struct fit_data *fit, detector *det) {
if(!fit || !det)
return EXIT_FAILURE;
if(sim_det_add(fit->sim, det)) {
return EXIT_FAILURE;
}
fit_data_exp_alloc(fit); /* Number of detectors in sim changed, let these guys know it too */
return EXIT_SUCCESS;
}
fit_data_det *fit_data_fdd(const fit_data *fit, size_t i_det) {
if(!fit|| !fit->fdd)
return NULL;
if(i_det >= fit->sim->n_det)
return NULL;
return &fit->fdd[i_det];
}
size_t fit_data_ranges_calculate_number_of_channels(const struct fit_data *fit_data) {
size_t sum = 0;
for(size_t i = 0; i < fit_data->n_fit_ranges; i++) {
roi *r = &fit_data->fit_ranges[i];
detector *det = sim_det(fit_data->sim, r->i_det);
if(!det) {
continue;
}
if(r->high >= det->channels) { /* Limited by detector. TODO: ideally ROIs should not be outside detector channels.. */
sum += (det->channels - r->low);
} else {
sum += (r->high - r->low) + 1;
}
}
return sum;
}
struct fit_stats fit_stats_init(void) {
struct fit_stats s;
s.n_evals = 0;
s.n_evals_iter = 0;
s.n_spectra = 0;
s.n_spectra_iter = 0;
s.cputime_cumul = 0.0;
s.cputime_iter = 0.0;
s.chisq0 = 0.0;
s.chisq = 0.0;
s.chisq_dof = 0.0;
s.rcond = 0.0;
s.iter = 0;
s.iter_call = 0;
s.error = FIT_ERROR_NONE;
return s;
}
void fit_data_print(const fit_data *fit, jabs_msg_level msg_level) {
if(!fit) {
return;
}
if(fit->n_fit_ranges == 0) {
jabs_message(MSG_ERROR, "No fit ranges.\n");
return;
}
jabs_message(msg_level, "Fit ROI # | detector name | low ch | high ch | exp counts | sim counts\n");
for(size_t i = 0; i < fit->n_fit_ranges; i++) {
roi *range = &fit->fit_ranges[i];
double exp = jabs_histogram_roi(fit_data_exp(fit, range->i_det), range->low, range->high);
double sim = jabs_histogram_roi(fit_data_histo_sum(fit, range->i_det), range->low, range->high);
jabs_message(msg_level, " %8zu | %13s | %7zu | %7zu | %10g | %10g\n",
i + 1, detector_name(sim_det(fit->sim, range->i_det)), range->low, range->high,
exp, sim);
}
jabs_message(msg_level, "\nFit %zu ranges with %zu channels total.\n", fit->n_fit_ranges, fit_data_ranges_calculate_number_of_channels(fit));
}
int jabs_test_delta(const gsl_vector *dx, const gsl_vector *x, double epsabs, double epsrel) { /* test_delta() copied from GSL convergence.c and modified */
int ok = TRUE;
DEBUGVERBOSEMSG("Test deltas to x->size=%zu\n", x->size);
for(size_t i = 0; i < x->size; i++) {
double xi = gsl_vector_get(x, i);
double dxi = gsl_vector_get(dx, i);
double tolerance = epsabs + epsrel * fabs(xi);
double rel = fabs(dxi) / tolerance; /* "How many times over the acceptable tolerance are we */
DEBUGVERBOSEMSG("Test delta: i %zu, xi %g, dxi %g, tolerance %g, rel %g\n", i, xi, dxi, tolerance, rel);
if(rel >= 1.0) {
DEBUGVERBOSEMSG("Fails because %g > 1.0.\n", rel);
ok = FALSE;
break;
}
}
if(ok)
return GSL_SUCCESS;
return GSL_CONTINUE;
}
int multifit_nlinear_print_jacobian(const gsl_multifit_nlinear_workspace *w, const char *filename) {
gsl_matrix *J = gsl_multifit_nlinear_jac(w);
FILE *f = fopen_file_or_stream(filename, "w");
if(!f) {
return EXIT_FAILURE;
}
for(size_t i = 0; i < J->size1; i++) {
for(size_t j = 0; j < J->size2; j++) {
double val = gsl_matrix_get(J, i, j);
fprintf(f, " %12e", val);
}
fprintf(f, "\n");
}
fclose_file_or_stream(f);
return EXIT_SUCCESS;
}
int jabs_gsl_multifit_nlinear_driver(const size_t maxiter, const double xtol, const double chisq_tol, struct fit_data *fit_data, gsl_multifit_nlinear_workspace *w) {
int status = 0;
size_t iter;
double chisq_dof_old;
jabs_message(MSG_INFO, "iter | cond(J) | |f(x)| | chisq/dof | evals | spectra | time cumul | time/spectrum |\n");
jabs_message(MSG_INFO, " | | | | cumul. | cumul. | s | ms |\n");
for(iter = 0; iter <= maxiter; iter++) {
fit_data->stats.iter_call = 0;
fit_data->stats.iter = iter;
if(iter) {
chisq_dof_old = fit_data->stats.chisq_dof;
fit_data->stats.cputime_iter = 0.0;
fit_data->stats.n_evals_iter = 0;
fit_data->stats.n_spectra_iter = 0;
status = gsl_multifit_nlinear_iterate(w);
DEBUGMSG("Iteration status %i (%s)", status, gsl_strerror(status));
}
if(fit_data->stats.error) {
return fit_data->stats.error;
}
if(status == GSL_ENOPROG && iter == 1) {
return FIT_ERROR_NO_PROGRESS;
}
fit_iter_stats_update(fit_data, w);
fit_iter_stats_print(&fit_data->stats);
#ifdef DEBUG
if(fit_data->stats.phase > FIT_PHASE_FAST) {
char *jacobian_filename;
asprintf(&jacobian_filename, "jacobian_iter%zu.dat", fit_data->stats.iter);
if(jacobian_filename) {
multifit_nlinear_print_jacobian(w, jacobian_filename);
free(jacobian_filename);
}
}
#endif
if(fit_data->fit_iter_callback) {
if(fit_data->fit_iter_callback(fit_data->stats)) {
return FIT_ERROR_ABORTED;
}
}
if(iter == 0)
continue;
/* test for convergence */
status = jabs_test_delta(w->dx, w->x, xtol * xtol, xtol);
if(status == GSL_SUCCESS) {
return FIT_SUCCESS_DELTA;
}
double chisq_change = 1.0 - fit_data->stats.chisq_dof / chisq_dof_old;
if(fit_data->stats.chisq_dof > chisq_dof_old) {
jabs_message(MSG_WARNING, "Chisq increased, this probably shouldn't happen.\n");
}
if(chisq_change < chisq_tol) {
return FIT_SUCCESS_CHISQ;
}
}
return FIT_ERROR_MAXITER;
}
void fit_report_results(const fit_data *fit, const gsl_multifit_nlinear_workspace *w, const gsl_multifit_nlinear_fdf *fdf) {
jabs_message(MSG_INFO, "summary from method '%s/%s'\n", gsl_multifit_nlinear_name(w), gsl_multifit_nlinear_trs_name(w));
jabs_message(MSG_INFO, "number of iterations: %zu\n", gsl_multifit_nlinear_niter(w));
jabs_message(MSG_INFO, "function evaluations: %zu\n", fit->stats.n_evals);
#ifdef DEBUG
jabs_message(MSG_INFO, "function evaluations (GSL): %zu\n", fdf->nevalf);
#endif
jabs_message(MSG_INFO, "Jacobian evaluations: %zu\n", fdf->nevaldf);
jabs_message(MSG_INFO, "number of spectra simulated: %zu\n", fit->stats.n_spectra);
jabs_message(MSG_INFO, "reason for stopping: %s\n", fit_error_str(fit->stats.error));
jabs_message(MSG_INFO, "initial |f(x)| = %f\n", sqrt(fit->stats.chisq0));
jabs_message(MSG_INFO, "final |f(x)| = %f\n", sqrt(fit->stats.chisq));
}
void fit_correlation_print(const gsl_matrix *covar, jabs_msg_level msg_level) {
jabs_message(msg_level, "\nCorrelation coefficients (sigma_ij/(sigma_i*sigma_j)) matrix:\n | ");
for(size_t i = 0; i < covar->size1; i++) {
jabs_message(msg_level, " %4zu ", i + 1);
}
jabs_message(msg_level, "\n");
for(size_t i = 0; i < covar->size1; i++) {
jabs_message(msg_level, "%6zu | ", i + 1);
for(size_t j = 0; j <= i && j < covar->size2; j++) {
jabs_message(msg_level, " %6.3f", gsl_matrix_get(covar, i, j) / sqrt(gsl_matrix_get(covar, i, i) * gsl_matrix_get(covar, j, j)));
}
jabs_message(msg_level, "\n");
}
}
int fit_uncertainty_spectra(const fit_data *fit, const gsl_matrix *J, const gsl_matrix *covar, const gsl_vector *f, const gsl_vector *w, const char *filename) {
FILE *f_errvec = NULL;
if(filename) {
f_errvec = fopen(filename, "w");
if(!f_errvec) {
return EXIT_FAILURE;
}
}
gsl_matrix *JC = gsl_matrix_alloc(J->size1, covar->size2); /* Product of J*C */
gsl_matrix *JCJT = gsl_matrix_alloc(J->size1, J->size1); /* Product J*C*JT */
gsl_blas_dgemm(CblasNoTrans, CblasNoTrans, 1.0, J, covar, 0.0, JC); /* TODO: we can probably optimize these matrix multiplications, since eventually we are only interested about the diagonal */
gsl_blas_dgemm(CblasNoTrans, CblasTrans, 1.0, JC, J, 0.0, JCJT);
gsl_vector_view err_vec = gsl_matrix_diagonal(JCJT);
size_t i_vec = 0;
jabs_histogram **h_minus = calloc(fit->n_det_spectra, sizeof(jabs_histogram *));
jabs_histogram **h_plus = calloc(fit->n_det_spectra, sizeof(jabs_histogram *));
for(size_t i_det = 0; i_det < fit->n_det_spectra; i_det++) {
for(size_t i_roi = 0; i_roi < fit->n_fit_ranges; i_roi++) {
const roi *roi = &fit->fit_ranges[i_roi];
const result_spectra *s = &(fit->spectra[i_det]);
const jabs_histogram *sim = result_spectra_simulated_histo(s);
if(!sim) {
continue;
}
if(roi->i_det == i_det) {
h_minus[i_det] = jabs_histogram_alloc(sim->n);
h_plus[i_det] = jabs_histogram_alloc(sim->n);
DEBUGMSG("Allocated uncertainty histograms for detector %zu/%zu, they have %zu channels. Pointers %p and %p.", i_det, fit->n_det_spectra, sim->n, (void *)h_minus[i_det], (void *)h_plus[i_det]);
break; /* Only once per detector */
}
}
}
for(size_t i_roi = 0; i_roi < fit->n_fit_ranges; i_roi++) {
const roi *roi = &fit->fit_ranges[i_roi];
const result_spectra *s = &(fit->spectra[roi->i_det]);
const jabs_histogram *sim = result_spectra_simulated_histo(s);
const jabs_histogram *exp = result_spectra_experimental_histo(s);
const jabs_histogram *h_uncertainty_low = h_minus[roi->i_det];
const jabs_histogram *h_uncertainty_high = h_plus[roi->i_det];
for(size_t ch = roi->low; ch <= roi->high; ch++) {
double sim_counts = sim ? jabs_histogram_get(sim, ch) : 0.0;
double exp_counts = exp ? jabs_histogram_get(exp, ch) : 0.0;