GRSISort "v4.0.0.5"
An extension of the ROOT analysis Framework
Loading...
Searching...
No Matches
TPulseAnalyzer.cxx
Go to the documentation of this file.
1#include <utility>
2
3#include "TPulseAnalyzer.h"
4
9
10TPulseAnalyzer::TPulseAnalyzer(const TFragment& fragment, double noise_fac)
11{
12 Clear();
13 SetData(fragment, noise_fac);
14}
15
16TPulseAnalyzer::TPulseAnalyzer(const std::vector<Short_t>& wave, double noise_fac, std::string name)
17 : fName(std::move(name))
18{
19 Clear();
20 SetData(wave, noise_fac);
21}
22
24{
25 delete cWpar;
26 delete spar;
27 delete shpar;
28}
29
30void TPulseAnalyzer::Clear(Option_t*)
31{
32 SetCsI(false);
33 set = false;
34 cN = 0;
35 FILTER = 8;
36 T0RANGE = 8;
37 LARGECHISQ = 1E111;
38 EPS = 0.001;
39
40 lineq_dim = 0;
41 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
42 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
43 memset(lineq_solution, 0, sizeof(lineq_solution)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
44 memset(copy_matrix, 0, sizeof(copy_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
45}
46
47void TPulseAnalyzer::SetData(const TFragment& fragment, double noise_fac)
48{
49 if(fragment.HasWave()) {
50 SetData(*fragment.GetWaveform(), noise_fac);
51 }
52}
53
54// NOLINTBEGIN(cppcoreguidelines-narrowing-conversions)
55void TPulseAnalyzer::SetData(const std::vector<Short_t>& wave, double noise_fac)
56{
57 SetCsI(false);
58 cN = wave.size();
59 if(cN > 0) {
60 cWavebuffer = wave;
61 set = true;
62 if(noise_fac > 0) {
63 FILTER = 8 * noise_fac;
64 T0RANGE = 8 * noise_fac;
65 }
66 }
67}
68
69////////////////////////////////////////
70// Linear Equation Solver
71////////////////////////////////////////
72
73// "Very efficient" apparently, written by Kris S,
74// Solve the currently stored n dimentional linear eqaution
76{
77 memcpy(copy_matrix, lineq_matrix, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
78 long double w = determinant(lineq_dim);
79 if(w == 0.) {
80 return 0;
81 }
82 for(int i = 0; i < lineq_dim; i++) {
83 memcpy(copy_matrix, lineq_matrix, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
84 memcpy(copy_matrix[i], lineq_vector, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
86 }
87 return 1;
88}
89
90// solve the determinant of the currently stored copy_matrix for dimentions m
92{
93 if(m == 1) {
94 return copy_matrix[0][0];
95 }
96 int sign = 1;
97 if(copy_matrix[m - 1][m - 1] == 0.) {
98 int j = m - 1;
99 while(copy_matrix[m - 1][j] == 0 && j >= 0) {
100 j--;
101 }
102 if(j < 0) {
103 return 0.;
104 }
105 for(int i = 0; i < m; i++) {
106 long double s = copy_matrix[i][m - 1];
107 copy_matrix[i][m - 1] = copy_matrix[i][j];
108 copy_matrix[i][j] = s;
109 }
110 sign *= -1;
111 }
112 for(int j = m - 2; j >= 0; j--) {
113 for(int i = 0; i < m; i++) {
114 copy_matrix[i][j] -= copy_matrix[i][m - 1] / copy_matrix[m - 1][m - 1] * copy_matrix[m - 1][j];
115 }
116 }
117 return copy_matrix[m - 1][m - 1] * sign * determinant(m - 1);
118}
119
120////////////////////////////////////////
121// Waveform Fits Functions
122////////////////////////////////////////
123
124int TPulseAnalyzer::fit_smooth_parabola(int low, int high, double x0, ParPar* pp)
125{
126 memset(pp, 0, sizeof(ParPar));
127 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
128 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
129 lineq_dim = 2;
130 double chisq = 0.;
131 int ndf = 0;
132 int k = static_cast<int>(rint(x0));
133
134 for(int i = low; i < k; i++) {
135 lineq_matrix[0][0] += 1;
136 lineq_vector[0] += cWavebuffer[i];
137 ndf++;
138 chisq += cWavebuffer[i] * cWavebuffer[i];
139 }
140
141 for(int i = k; i < high; i++) {
142 double x = (i - x0) * (i - x0);
143 lineq_matrix[0][0] += 1;
144 lineq_matrix[0][1] += x;
145 lineq_matrix[1][1] += x * x;
146 lineq_vector[0] += cWavebuffer[i];
147 lineq_vector[1] += cWavebuffer[i] * x;
148 ndf++;
149 chisq += cWavebuffer[i] * cWavebuffer[i];
150 }
151 lineq_matrix[1][0] = lineq_matrix[0][1];
152
153 if(solve_lin_eq() == 0) {
154 pp->chisq = BADCHISQ_MAT;
155 return -1;
156 }
157 chisq -= lineq_vector[0] * lineq_solution[0];
158 chisq -= lineq_vector[1] * lineq_solution[1];
159
160 pp->constant = lineq_solution[0];
161 pp->linear = 0.;
162 pp->quadratic = lineq_solution[1];
163 pp->chisq = chisq;
164 pp->ndf = ndf;
165
166 return 1;
167}
168
169////////////////////////////////////////
170// RF Fit Run Functions
171////////////////////////////////////////
172
174{
175 if(!set || cN < 10) {
176 return -1;
177 }
178 delete spar;
179 spar = new SinPar;
180
181 spar->t0 = -1;
182
183 return 5 * get_sin_par(T);
184}
185
186////////////////////////////////////////
187// Waveform Time Fit Run Functions
188////////////////////////////////////////
189
190// Overall function which determins limits and fits the 3 trial functions
192{
193 if(!set || cN < 10) {
194 return -1;
195 }
196
197 delete cWpar;
198 cWpar = new WaveFormPar;
199 cWpar->t0 = -1;
200
201 std::array<double, 3> chisq;
202 std::array<WaveFormPar, 3> w;
203
204 cWpar->baseline_range = T0RANGE; // default only 8 samples!
205 get_baseline();
206 get_tmax();
207
208 if(good_baseline() == 0) {
209 return -1;
210 }
211
212 get_t30();
213 get_t50();
214 cWpar->thigh = cWpar->t50;
215
216 for(double& i : chisq) {
217 i = LARGECHISQ;
218 }
219
220 size_t swp = sizeof(WaveFormPar);
221 chisq[0] = get_smooth_T0();
222 memcpy(w.data(), cWpar, swp);
223 chisq[1] = get_parabolic_T0();
224 memcpy(&w[1], cWpar, swp);
225 chisq[2] = get_linear_T0();
226 memcpy(&w[2], cWpar, swp);
227
228 double chimin = LARGECHISQ;
229 int imin = 0;
230
231 for(int i = 0; i < 3; i++) {
232 if(chisq[i] < chimin && chisq[i] > 0) {
233 chimin = chisq[i];
234 imin = i;
235 }
236 }
237
238 if(imin < 2) {
239 memcpy(cWpar, &w[imin], swp);
240 }
241
243 return cWpar->t0;
244}
245/*================================================================*/
246
247/*================================================================*/
248int TPulseAnalyzer::fit_parabola(int low, int high, ParPar* pp)
249{
250 memset(pp, 0, sizeof(ParPar));
251 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
252 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
253 lineq_dim = 3;
254 double chisq = 0.;
255 int ndf = 0;
256 for(int i = low; i < high; i++) {
257 lineq_matrix[0][0] += 1;
258 lineq_matrix[0][1] += i;
259 lineq_matrix[0][2] += i * i;
260 lineq_matrix[1][2] += i * i * i;
261 lineq_matrix[2][2] += i * i * i * i;
262 lineq_vector[0] += cWavebuffer[i];
263 lineq_vector[1] += cWavebuffer[i] * i;
264 lineq_vector[2] += cWavebuffer[i] * i * i;
265 ndf++;
266 chisq += cWavebuffer[i] * cWavebuffer[i];
267 }
268 lineq_matrix[1][0] = lineq_matrix[0][1];
269 lineq_matrix[1][1] = lineq_matrix[0][2];
270 lineq_matrix[2][0] = lineq_matrix[0][2];
271 lineq_matrix[2][1] = lineq_matrix[1][2];
272
273 if(solve_lin_eq() == 0) {
274 pp->chisq = BADCHISQ_MAT;
275 return -1;
276 }
277 chisq -= lineq_vector[0] * lineq_solution[0];
278 chisq -= lineq_vector[1] * lineq_solution[1];
279 chisq -= lineq_vector[2] * lineq_solution[2];
280 pp->constant = lineq_solution[0];
281 pp->linear = lineq_solution[1];
282 pp->quadratic = lineq_solution[2];
283 pp->chisq = chisq;
284 pp->ndf = ndf;
285 return 1;
286}
287
288/*================================================================*/
289int TPulseAnalyzer::fit_line(int low, int high, LinePar* lp)
290{
291 memset(lp, 0, sizeof(LinePar));
292 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
293 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
294 lineq_dim = 2;
295 double chisq = 0.;
296 int ndf = 0;
297 for(int i = low; i < high; i++) {
298 lineq_matrix[0][0] += 1;
299 lineq_matrix[0][1] += i;
300 lineq_matrix[1][1] += i * i;
301 lineq_vector[0] += cWavebuffer[i];
302 lineq_vector[1] += cWavebuffer[i] * i;
303 ndf++;
304 chisq += cWavebuffer[i] * cWavebuffer[i];
305 }
306 lineq_matrix[1][0] = lineq_matrix[0][1];
307
308 if(solve_lin_eq() == 0) {
309 lp->chisq = BADCHISQ_MAT;
310 return -1;
311 }
312 chisq -= lineq_vector[0] * lineq_solution[0];
313 chisq -= lineq_vector[1] * lineq_solution[1];
314 lp->slope = lineq_solution[1];
315 lp->intercept = lineq_solution[0];
316 lp->chisq = chisq;
317 lp->ndf = ndf;
318 return 1;
319}
320/*================================================================*/
321
322/*======================================================*/
323
325{
326 LinePar lp = {};
327 LinePar lpl = {};
328 int k = 0;
329 double chit = 0.;
330 double chitmin = 0.;
331 double b = 0.;
332 double c = 0.;
333 double t = 0.;
334
335 chitmin = LARGECHISQ;
336
337 for(k = T0RANGE / 2; k < cWpar->thigh - T0RANGE / 8; k++) {
338 // fit line to the baseline
339 fit_line(0, k, &lp);
340
341 // fit line to the risetime
342 fit_line(k, cWpar->thigh, &lpl);
343
344 chit = lp.chisq + lpl.chisq;
345
346 if(chit < chitmin) {
347 chitmin = chit;
348 cWpar->b0 = lp.intercept;
349 cWpar->b1 = lp.slope;
350 cWpar->s0 = lpl.intercept;
351 cWpar->s1 = lpl.slope;
352 cWpar->s2 = 0.;
353 }
354 } // end of the loop over k
355 b = cWpar->s1 - cWpar->b1;
356 c = cWpar->s0 - cWpar->b0;
357 t = -c / b;
358
359 cWpar->t0 = -1;
360 cWpar->temin = 0;
362 if(t < cN && t > 0) {
363 cWpar->t0 = t;
364 cWpar->temin = static_cast<int>(rint(cWpar->t0)) - 2;
365 cWpar->temax = static_cast<int>(rint(cWpar->t0)) + 2;
366 return (chitmin / (cWpar->thigh - 5));
367 }
368
369 return BADCHISQ_LIN_T0;
370}
371/*================================================================*/
373{
374 ParPar pp;
375 ParPar ppmin;
376 int k = 0;
377 int kmin = 0;
378 double chit = 0.;
379 double chitmin = 0.;
380 double c = 0.;
381 double t = 0.;
382
383 memset(&ppmin, 0, sizeof(ParPar));
384
385 chitmin = LARGECHISQ;
386 kmin = 0;
387 // corse search first
388 for(k = T0RANGE / 2; k < cWpar->thigh - T0RANGE / 2; k++) {
389 fit_smooth_parabola(0, cWpar->thigh, static_cast<double>(k), &pp);
390
391 chit = pp.chisq;
392 if(chit < chitmin) {
393 chitmin = chit;
394 kmin = k;
395 }
396 } // end of the corse search loop over k
397 c = kmin;
398
399 chitmin = LARGECHISQ;
400 // fine search next
401 for(t = kmin - 1; t < kmin + 1; t += 0.1) {
402 fit_smooth_parabola(0, cWpar->thigh, t, &pp);
403 chit = pp.chisq;
404 if(chit < chitmin) {
405 memcpy(&ppmin, &pp, sizeof(ParPar));
406 chitmin = chit;
407 c = t;
408 }
409 } // end of the fine search loop over k
410
411 memcpy(&pp, &ppmin, sizeof(ParPar));
412 t = c;
413 cWpar->s0 = pp.constant + pp.quadratic * t * t;
414 cWpar->s1 = -2. * pp.quadratic * t;
415 cWpar->s2 = pp.quadratic;
416 cWpar->b0 = pp.constant;
417 cWpar->b1 = 0.;
418
419 cWpar->t0 = -1;
420 cWpar->temin = 0;
422 if(t < cN && t > 0) {
423 cWpar->t0 = t;
424 cWpar->temin = static_cast<int>(rint(cWpar->t0)) - 2;
425 cWpar->temax = static_cast<int>(rint(cWpar->t0)) + 2;
426 return (chitmin / (cWpar->thigh - 2));
427 }
428 return BADCHISQ_SMOOTH_T0;
429}
430/*================================================================*/
432{
433 LinePar lp;
434 ParPar pp;
435 double chit = 0.;
436 double chitmin = 0.;
437 double a = 0.;
438 double b = 0.;
439 double c = 0.;
440 double d = 0.;
441 double t = 0.;
442
443 chitmin = LARGECHISQ;
444 for(int k = T0RANGE / 2; k < cWpar->thigh - T0RANGE / 2; k++) {
445 // fit line to the baseline
446 fit_line(0, k, &lp);
447
448 // fit parabola to the risetime
449 fit_parabola(k, cWpar->thigh, &pp);
450
451 chit = lp.chisq + pp.chisq;
452
453 if(chit < chitmin) {
454 chitmin = chit;
455 cWpar->b0 = lp.intercept;
456 cWpar->b1 = lp.slope;
457 cWpar->s0 = pp.constant;
458 cWpar->s1 = pp.linear;
459 cWpar->s2 = pp.quadratic;
460 }
461 } // end loop through k
462
463 a = cWpar->s2;
464 b = cWpar->s1 - cWpar->b1;
465 c = cWpar->s0 - cWpar->b0;
466 d = b * b - 4 * a * c;
467
468 if(a == 0.) {
469 t = -c / b;
470 } else {
471 if(d >= 0) {
472 if(d == 0.) {
473 t = -0.5 * b / a;
474 } else {
475 d = sqrt(d);
476 t = 0.5 * (-b + d) / a;
477 }
478 } else {
479 return BADCHISQ_PAR_T0;
480 }
481 }
482
483 cWpar->t0 = -1;
484 cWpar->temin = 0;
486 if(t < cN && t > 0) {
487 cWpar->t0 = t;
488 cWpar->temin = static_cast<int>(rint(cWpar->t0)) - 2;
489 cWpar->temax = static_cast<int>(rint(cWpar->t0)) + 2;
490 return (chitmin / (cWpar->thigh - 5));
491 }
492 return BADCHISQ_PAR_T0;
493}
494
495// Measure the baseline and standard deviation of the waveform, over the tick range specified by
496// cWpar->baseline_range
498{
499 cWpar->baseline = 0.;
500 cWpar->baselineStDev = 0.;
501
502 // error if waveform length cN is shorter than baseline range
503 if(cN < cWpar->baseline_range) {
504 std::cout << "Baseline range (" << cWpar->baseline_range << ") larger than waveform length!" << std::endl;
505 std::cout << "Terminating program" << std::endl;
506 exit(0);
507 }
508
509 for(int i = 0; i < cWpar->baseline_range; i++) {
512 }
513
518 cWpar->bflag = true; // flag after establishing baseline
519}
520
521/*======================================================*/
522
523// Measure the baseline and standard deviation up to tick cWpar->t0 after a fit
525{
526 cWpar->baselinefin = 0.;
528 int tb = cWpar->t0; // t0 non integer, result always too small before.
529 if(tb > T0RANGE + 10) {
530 tb -= 10;
531 }
532
533 // error if waveform length cN is shorter than baseline range
534 if(cN > tb && tb > 0) {
535 for(int i = 0; i < tb; i++) {
538 }
539
540 cWpar->baselineStDevfin /= tb;
541 cWpar->baselinefin /= tb;
543 cWpar->baselineStDevfin = sqrt(std::abs(cWpar->baselineStDevfin));
544 }
545}
546
547/*======================================================*/
548
549// Find the maximum of the wavefunction, smoothed with a moving average filter
551{
552 int i = 0;
553 int j = 0;
554 int sum = 0;
555 int D = FILTER / 2;
556
557 cWpar->max = cWavebuffer[0];
558 cWpar->tmax = 0;
559
560 // applies the filter to the waveform
561 // cout<<" "<<cWpar->tmax<<" "<< cWpar->max<<flush;
562 for(i = D; i < cN - D; i++) {
563 sum = 0;
564 for(j = i - D; j < i + D; j++) {
565 sum += cWavebuffer[j];
566 }
567 sum /= FILTER; // the value of the filtered waveform at i
568 if(sum > cWpar->max) {
569 // if the value of the filtered waveform at i is larger than the current maximum, max=value and tmax = i
570 cWpar->max = sum;
571 cWpar->tmax = i;
572 }
573 }
574 cWpar->mflag = 1; // flag after finding tmax
575}
576
577/*===========================================================*/
578double TPulseAnalyzer::get_tfrac(double frac, double fraclow, double frachigh)
579{
580 if(!cWpar->bflag) {
581 std::cout << "Baseline not deterimned for the tfraction" << std::endl;
582 exit(1);
583 }
584
585 if(cWpar->mflag != 1) {
586 std::cout << "Maximum not deterimned for the tfraction" << std::endl;
587 exit(1);
588 }
589
590 int t = cWpar->tmax;
591
592 double f = cWpar->baseline + frac * (cWpar->max - cWpar->baseline);
593 double flow = cWpar->baseline + fraclow * (cWpar->max - cWpar->baseline);
594 double fhigh = cWpar->baseline + frachigh * (cWpar->max - cWpar->baseline);
595
596 while(cWavebuffer[t] > f) {
597 t--;
598 if(t <= 4) {
599 break;
600 }
601 }
602 int imin = t;
603 while(cWavebuffer[imin] > flow) {
604 imin--;
605 if(imin <= 1) {
606 break;
607 }
608 }
609
610 int imax = t;
611
612 while(cWavebuffer[imax] < fhigh) {
613 imax++;
614 if(imax >= cN - 1) {
615 break;
616 }
617 }
618
619 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
620 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
621 lineq_dim = 3;
622
623 int i = imax - imin;
624 int64_t a = i;
625 lineq_matrix[0][0] = a + 1;
626 lineq_matrix[0][1] = 0.5 * a;
627 lineq_matrix[2][0] = a / 6.;
628 lineq_matrix[2][2] = -a / 30.;
629 a *= i;
630 lineq_matrix[0][1] += 0.5 * a;
631 lineq_matrix[2][0] += 0.5 * a;
632 lineq_matrix[2][1] = 0.25 * a;
633 a *= i;
634 lineq_matrix[2][0] += a / 3.;
635 lineq_matrix[2][1] += 0.5 * a;
636 lineq_matrix[2][2] += a / 3.;
637 a *= i;
638 lineq_matrix[2][1] += 0.25 * a;
639 lineq_matrix[2][2] += 0.5 * a;
640 a *= i;
641 lineq_matrix[2][2] += 0.2 * a;
642
643 lineq_matrix[1][0] = lineq_matrix[0][1];
644 lineq_matrix[1][1] = lineq_matrix[2][0];
645 lineq_matrix[0][2] = lineq_matrix[2][0];
646 lineq_matrix[1][2] = lineq_matrix[2][1];
647
648 for(i = 0; i < lineq_dim; i++) {
649 lineq_vector[i] = 0;
650 }
651
652 for(i = imin; i < imax + 1; i++) {
653 a = i - imin;
654 lineq_vector[0] += cWavebuffer[i];
655 lineq_vector[1] += cWavebuffer[i] * a;
656 lineq_vector[2] += cWavebuffer[i] * a * a;
657 }
658
659 if(solve_lin_eq() == 0) {
660 return -4;
661 }
662 double p = lineq_solution[0] - f;
663 double q = lineq_solution[1];
664 double r = lineq_solution[2];
665
666 if(r != 0) {
667 double d = q * q - 4 * r * p;
668 if(d < 0) {
669 return -5;
670 }
671 f = -q + sqrt(d);
672 f *= 0.5;
673 f /= r;
674 f += imin;
675 return f;
676 }
677 if(q != 0) {
678 f = -p / q;
679 return f;
680 }
681 return -6;
682
683 return -7;
684}
685
686/* ==================================================== */
688{
689 int t = get_tfrac(0.5, 0.3, 0.8);
690 if((t > 0) && (t < MAX_SAMPLES)) {
691 cWpar->t50_flag = 1;
692 cWpar->t50 = t;
693 } else {
694 cWpar->t50_flag = -1;
695 cWpar->t50 = -1;
696 }
697}
698/* ==================================================== */
700{
701 int t = get_tfrac(0.9, 0.8, 0.98);
702
703 if((t > 0) && (t < MAX_SAMPLES)) {
704 cWpar->t90_flag = 1;
705 cWpar->t90 = t;
706 } else {
707 cWpar->t90_flag = -1;
708 cWpar->t90 = -1;
709 }
710}
711/*===========================================================*/
713{
714 int t = get_tfrac(0.1, 0.05, 0.2);
715
716 if((t > 0) && (t < MAX_SAMPLES)) {
717 cWpar->t10_flag = 1;
718 cWpar->t10 = t;
719 } else {
720 cWpar->t10_flag = -1;
721 cWpar->t10 = -1;
722 }
723}
724/*===========================================================*/
726{
727 int t = get_tfrac(0.3, 0.15, 0.45);
728 if((t > 0) && (t < MAX_SAMPLES)) {
729 cWpar->t30_flag = 1;
730 cWpar->t30 = t;
731 } else {
732 cWpar->t30_flag = -1;
733 cWpar->t30 = -1;
734 }
735}
736
738{
739 double s = 0.;
740 double sn = 0.;
741 double snm = 0.;
742 double s2 = 0.;
743 double s2n = 0.;
744 double s2nm = 0.;
745 double c = 0.;
746 double cn = 0.;
747 double cnm = 0.;
748 double c2 = 0.;
749 double c2n = 0.;
750 double c2nm = 0.;
751 double w = 0.;
752 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
753 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
754 lineq_dim = 3;
755
756 w = 2 * TMath::Pi() / T;
757
758 s = sin(w);
759 sn = sin(cN * w);
760 snm = sin((cN - 1) * w);
761 s2 = sin(2 * w);
762 s2n = sin(2 * cN * w);
763 s2nm = sin(2 * (cN - 1) * w);
764
765 c = cos(w);
766 cn = cos(cN * w);
767 cnm = cos((cN - 1) * w);
768 c2 = cos(2 * w);
769 c2n = cos(2 * cN * w);
770 c2nm = cos(2 * (cN - 1) * w);
771
772 lineq_matrix[0][0] = 0.5 * cN - 0.25 * (1 - c2 - c2n + c2nm) / (1 - c2);
773 lineq_matrix[0][1] = 0.25 * (s2 + s2nm - s2n) / (1 - c2);
774 lineq_matrix[1][0] = lineq_matrix[0][1];
775 lineq_matrix[0][2] = 0.5 * (s + snm - sn) / (1 - c);
776 lineq_matrix[2][0] = lineq_matrix[0][2];
777 lineq_matrix[1][1] = 0.5 * cN + 0.25 * (1 - c2 - c2n + c2nm) / (1 - c2);
778 lineq_matrix[1][2] = 0.5 * (1 - c - cn + cnm) / (1 - c);
779 lineq_matrix[2][1] = lineq_matrix[1][2];
780 lineq_matrix[2][2] = cN;
781
782 for(int i = 0; i < lineq_dim; i++) {
783 lineq_vector[i] = 0;
784 }
785
786 for(int i = 0; i < cN; i++) {
787 lineq_vector[0] += cWavebuffer[i] * sin(w * i);
788 lineq_vector[1] += cWavebuffer[i] * cos(w * i);
789 lineq_vector[2] += cWavebuffer[i];
790 }
791 if(solve_lin_eq() == 0) {
792 return 0;
793 }
795 spar->C = lineq_solution[2];
796
797 s = -lineq_solution[1] / spar->A;
798 c = lineq_solution[0] / spar->A;
799
800 if(s >= 0) {
801 spar->t0 = acos(c) * T / (2 * TMath::Pi());
802 } else {
803 spar->t0 = (1 - acos(c) / (2 * TMath::Pi())) * T;
804 }
805
806 return spar->t0;
807}
808
809/*======================================================*/
810// void TPulseAnalyzer::get_sig2noise(){ if(cN==0)return;
811// if(set){
812// get_baseline();
813// get_tmax();
814// cWpar->sig2noise=(cWpar->max-cWpar->baseline)/cWpar->baselineStDev;
815// }
816// return;
817// }
818/*======================================================*/
820{
821 if(set && (cWpar != nullptr)) {
822 if(cWpar->t0 > 0) {
824 }
825 }
826 return -1;
827}
828
830{
831 if(set && (cWpar != nullptr)) {
832 if(cWpar->tmax < T0RANGE) {
833 return 0;
834 }
835 if((cWpar->max - cWpar->baseline) < (cWpar->baseline - cWavebuffer[0]) * 10) {
836 return 0;
837 }
839 return 0;
840 }
841 return 1;
842 }
843 return 0;
844}
845
846//=====================================================//
847// CsI functions:
848//=====================================================//
849
851{
852 if(CsIIsSet()) {
853 return shpar->chisq;
854 }
855 return -1;
856}
857
859{
860 if(CsIIsSet()) {
861 return shpar->type;
862 }
863 return -1;
864}
865
867{
868
869 if(CsIIsSet()) {
870 return shpar->t[0];
871 }
872 if(!set || cN < 10) {
873 return -1.;
874 }
875
876 delete cWpar;
877 cWpar = new WaveFormPar;
878 delete shpar;
879 shpar = new ShapePar;
880
882 if(cWpar->teflag == 1) {
883 // good exclusion zone
884 int tmpchisq = GetCsIShape();
885 if(tmpchisq >= 0) {
886 SetCsI();
887 return shpar->t[0];
888 }
889 }
890
891 return -1.0;
892}
893
895{
896 // std::cout<<"Fitting CsI PID"<<std::endl;
897 if(CsIIsSet()) {
898 double f = shpar->am[2];
899 double s = shpar->am[3];
900 double r = s / f * 100;
901
902 return r;
903 }
904 if(!set || cN < 10) {
905 return -1.;
906 }
907
908 delete cWpar;
909 cWpar = new WaveFormPar;
910 delete shpar;
911 shpar = new ShapePar;
912
913 shpar->t[1] = 4510;
914 shpar->t[2] = 64.3;
915 shpar->t[3] = 380.0;
916
918
919 if(!cWpar->bflag) {
920 shpar->type = -2; // type for failed exclusion zone determination
921 return BAD_BASELINE_RANGE;
922 }
923 if(cWpar->mflag == 0) {
924 shpar->type = -2; // type for failed exclusion zone determination
925 return BAD_MAX;
926 }
927 if(cWpar->teflag == 0) {
928 shpar->type = -2; // type for failed exclusion zone determination
929 return BAD_EXCLUSION_ZONE;
930 }
931
932 int tmpchisq = GetCsIShape();
933 if(tmpchisq >= 0) {
934 double f = shpar->am[2];
935 double s = shpar->am[3];
936 double r = s / f * 100;
937
938 SetCsI();
939
940 return r;
941 }
942
943 return BAD_PID;
944}
945
947{
948
949 std::array<double, 4> chisq; // chisq array for fit types: 0=2 comp, 1=fast, 2=slow, 3=gamma on PIN
950
951 int imin = 0;
952
953 // calculate ndf assuming two component (4 parameter) fit
954 int ndf = -4;
955 for(int i = 0; i < cWpar->temin; i++) {
956 ndf++;
957 }
958 for(int i = cWpar->temax; i < cN; i++) {
959 ndf++;
960 }
961
962 // std::cout<<"ndf 4 parameters: "<<ndf<<std::endl;
963
964 // initialize chisq to large value and set up waveform parameters
965 for(int i = 0; i < 4; i++) {
966 chisq[i] = LARGECHISQ;
967 csiTestShpar[i] = new ShapePar;
968 csiTestWpar[i] = new WaveFormPar;
969 memcpy(csiTestShpar[i], shpar, sizeof(ShapePar));
970 memcpy(csiTestWpar[i], cWpar, sizeof(WaveFormPar));
971 }
972
973 // two component fit
975 chisq[0] = csiTestShpar[0]->chisq / ndf;
976
977 // for 3 parameter fits, ndf is one larger
978 ndf++;
979
980 // std::cout<<"ndf 3 parameters: "<<ndf<<std::endl;
981
982 // fast only for high Z recoils
984 chisq[1] = csiTestShpar[1]->chisq / ndf;
985
986 // slow only for gamma on CsI
987 csiTestShpar[2]->t[2] = shpar->t[3];
988 csiTestShpar[2]->t[3] = shpar->t[2];
990 chisq[2] = csiTestShpar[2]->chisq / ndf;
991
992 // gamma on PIN
993 csiTestShpar[3]->t[2] = shpar->t[4];
994 csiTestShpar[3]->t[4] = shpar->t[2];
996 chisq[3] = csiTestShpar[3]->chisq / ndf;
997
998 /*std::cout<<"0=two comp 1=fast 2=slow 3=gamma on PIN"<<std::endl;
999 for(i = 0; i < 4; i++)
1000 std::cout<<"chisq "<<csiTestShpar[i]->chisq<<", chisq["<<i<<"]/ndf "<<chisq[i]<<std::endl;*/
1001
1002 // find minimum chisq
1003 imin = -1;
1004 double chimin = LARGECHISQ;
1005 for(int i = 0; i < 4; i++) {
1006 if((chisq[i] < chimin) && (chisq[i] > 0)) {
1007 chimin = chisq[i];
1008 imin = i;
1009 }
1010 }
1011
1012 // std::cout<<"minimum chisq["<<imin<<"]/ndf "<<chisq[imin]<<", t0min "<<csiTestShpar[imin]->t[0]<<std::endl;
1013
1014 switch(imin) {
1015 case 3:
1016 // gamma on PIN fit type
1017 memcpy(shpar, csiTestShpar[imin], sizeof(ShapePar));
1018 memcpy(cWpar, csiTestWpar[imin], sizeof(WaveFormPar));
1019 shpar->t[2] = csiTestShpar[imin]->t[4];
1020 shpar->am[2] = csiTestShpar[imin]->am[4];
1021 shpar->t[4] = csiTestShpar[imin]->t[2];
1022 shpar->am[4] = csiTestShpar[imin]->am[2];
1023 shpar->type = imin + 1; // gamma on PIN type
1024 break;
1025 case 2:
1026 // slow component only fit type
1027 memcpy(shpar, csiTestShpar[imin], sizeof(ShapePar));
1028 memcpy(cWpar, csiTestWpar[imin], sizeof(WaveFormPar));
1029 shpar->t[2] = csiTestShpar[imin]->t[3];
1030 shpar->am[2] = csiTestShpar[imin]->am[3];
1031 shpar->t[3] = csiTestShpar[imin]->t[2];
1032 shpar->am[3] = csiTestShpar[imin]->am[2];
1033 shpar->type = imin + 1; // slow only type
1034 break;
1035 case 1:
1036 // fast component only fit type
1037 memcpy(shpar, csiTestShpar[imin], sizeof(ShapePar));
1038 memcpy(cWpar, csiTestWpar[imin], sizeof(WaveFormPar));
1039 shpar->type = imin + 1; // fast only type
1040 break;
1041 case 0:
1042 // two component fit type
1043 memcpy(shpar, csiTestShpar[imin], sizeof(ShapePar));
1044 memcpy(cWpar, csiTestWpar[imin], sizeof(WaveFormPar));
1045 shpar->type = imin + 1; // two component type
1046 break;
1047 default:
1048 // failed fit
1049 shpar->type = -1; // fit failure type
1050 shpar->chisq = BADCHISQ_FAIL_DIRECT; // set it here so it still frees memory on bad fits
1051 break;
1052 }
1053
1054 // free memory allocated for the fit
1055 for(int i = 0; i < 4; i++) {
1056 delete csiTestShpar[i];
1057 delete csiTestWpar[i];
1058 }
1059
1060 return shpar->chisq; // generic chisq return statement for all types
1061}
1062
1064{
1065 long double sum = 0.;
1066 long double tau = 0.;
1067 long double tau_i = 0.;
1068 long double tau_j = 0.;
1069 int p = 0;
1070 int q = 0;
1071 int d = 0;
1072
1073 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1074 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1075 memset(lineq_solution, 0, sizeof(lineq_solution)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1076
1077 /* q is the low limit of the signal section */
1078 q = wpar->temax;
1079 if(q >= cN || q <= 0) {
1080 par->chisq = -1;
1081 return -1.;
1082 }
1083
1084 /* p is the high limit of the baseline section */
1085 p = wpar->temin;
1086 if(p >= cN || p <= 0) {
1087 par->chisq = -1;
1088 return -1.;
1089 }
1090 lineq_dim = dim;
1091
1092 // initialize amplitudes to 0
1093 for(long double& i : par->am) { // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1094 i = 0.;
1095 }
1096
1097 d = cN;
1098
1099 // initialize chi square 0 and ndf = n-k to -k where k=dim
1100 par->chisq = 0.;
1101 par->ndf = -lineq_dim;
1102
1103 /**************************************************************************
1104 linearized chi square fit is Mu = v where M is a data matrix
1105 u, v are vectors; u is the parameter vector (solution)
1106 note that in this formulation, chisq_min = y_i^2-sum(u_iv_i)
1107 **************************************************************************/
1108
1109 // create matrix for linearized fit
1110 for(int i = 1; i < lineq_dim; i++) {
1111 tau = GetCsITau(i, par);
1112 tau_i = tau;
1113 sum = -(static_cast<double>(q)) / tau + log(1. - exp(-(static_cast<double>(d - q)) / tau));
1114 sum -= log(1. - exp(-1. / tau));
1115 lineq_matrix[i][0] = exp(sum);
1116 lineq_matrix[0][i] = exp(sum);
1117
1118 tau /= 2.;
1119 sum = -(static_cast<double>(q)) / tau + log(1. - exp(-(static_cast<double>(d - q)) / tau));
1120 sum -= log(1. - exp(-1. / tau));
1121 lineq_matrix[i][i] = exp(sum);
1122
1123 for(int j = i + 1; j < lineq_dim; j++) {
1124 tau_j = GetCsITau(j, par);
1125 tau = (tau_i * tau_j) / (tau_i + tau_j);
1126 sum = -(static_cast<double>(q)) / tau + log(1. - exp(-(static_cast<double>(d - q)) / tau));
1127 sum -= log(1. - exp(-1. / tau));
1128 lineq_matrix[i][j] = exp(sum);
1129 lineq_matrix[j][i] = exp(sum);
1130 }
1131 }
1132
1133 lineq_vector[0] = 0;
1134 lineq_matrix[0][0] = 0;
1135
1136 for(int j = q; j < cN; j++) {
1137 lineq_vector[0] += cWavebuffer[j];
1138 lineq_matrix[0][0] += 1;
1139 par->chisq += cWavebuffer[j] * cWavebuffer[j];
1140 par->ndf += 1;
1141 }
1142
1143 if(lineq_dim >= cN) {
1144 par->chisq = -1;
1145 return -1.;
1146 }
1147
1148 for(int i = 1; i < lineq_dim; i++) {
1149 tau = GetCsITau(i, par);
1150 lineq_vector[i] = 0;
1151 for(int j = q; j < cN; j++) {
1152 lineq_vector[i] += cWavebuffer[j] * exp(-(static_cast<double>(j)) / tau);
1153 }
1154 }
1155
1156 for(int j = 0; j < p; j++) {
1157 lineq_vector[0] += cWavebuffer[j];
1158 lineq_matrix[0][0] += 1;
1159 par->chisq += cWavebuffer[j] * cWavebuffer[j];
1160 par->ndf += 1;
1161 }
1162
1163 // solve the matrix equation Mu = v -> u = M^(-1)v where M^(-1) is the inverse
1164 // of M. note this has no solution if M is not invertable!
1165
1166 // error if the matrix cannot be inverted
1167 if(solve_lin_eq() == 0) {
1168 par->chisq = BADCHISQ_MAT;
1169 par->ndf = 1;
1170
1171 return BADCHISQ_MAT;
1172 } // else try and find t0 and calculate amplitudes
1173
1174 // see the function comments for find_t0 for details
1175
1176 par->t[0] = GetCsIt0(par, wpar);
1177
1178 // if t0 is less than 0, return a T0FAIL
1179 if(par->t[0] <= 0) {
1180 par->chisq = BADCHISQ_T0;
1181 par->ndf = 1;
1182 return BADCHISQ_T0;
1183 }
1184
1185 // calculate amplitudes
1186 par->am[0] = lineq_solution[0];
1187
1188 for(int i = 1; i < lineq_dim; i++) {
1189 tau = GetCsITau(i, par);
1190 par->am[i] = lineq_solution[i] * exp(-par->t[0] / tau);
1191 }
1192 // done calculating amplitudes
1193
1194 for(int i = 0; i < lineq_dim; i++) {
1195 par->chisq -= lineq_solution[i] * lineq_vector[i];
1196 }
1197
1198 if(par->chisq < 0) {
1199 par->chisq = BADCHISQ_NEG;
1200 par->ndf = 1;
1201 return BADCHISQ_NEG;
1202 }
1203
1204 for(int i = 2; i < lineq_dim; i++) {
1205 par->am[i] *= -1;
1206 }
1207
1208 par->type = dim - 2;
1209
1210 // return BADCHISQ_AMPL if a component amplitude is less than 0
1211 //(apart from the baseline which can be negative)
1212 for(int i = 1; i < lineq_dim; i++) {
1213 if(par->am[i] < 0) {
1214 par->chisq = BADCHISQ_AMPL;
1215 par->ndf = 1;
1216 return BADCHISQ_AMPL;
1217 }
1218 }
1219
1220 // std::cout<<"chisq from FitCsIShape: "<<par->chisq<<std::endl;
1221
1222 return par->chisq;
1223}
1224
1226{
1227 int i = 0;
1228 int j = 0;
1229 int D = FILTER / 2; // filter half width
1230
1231 // initilize the fit parameters for the risetime to 0 for safety
1232 cWpar->afit = 0.;
1233 cWpar->bfit = 0.;
1234
1235 // make sure the baseline is established prior to finding the exclusion zone
1237 get_baseline();
1238
1239 // Here we determine the x position temax of the upper limit for the exclusion zone.
1240 // find tmax and define baselineMax
1241 get_tmax();
1242
1243 // If tmax is established, continue.
1244 if(cWpar->mflag == 1) {
1246
1247 // Starting at tmax and working backwards along the waveform get the value of the filtered waveform at i and
1248 // when
1249 // the value of the filtered waveform goes below baselineMax, set the upper limit of the exclusion zone temax =
1250 // i.
1251 // The exclusion zone cannot be defined in the area of the waveform used to calculate the baseline.
1252 for(i = cWpar->tmax; i > cWpar->baseline_range; i--) {
1253 double sum = 0.;
1254 for(j = i - D; j < i + D; j++) {
1255 sum += cWavebuffer[j];
1256 }
1257 sum /= FILTER;
1258 if(sum < cWpar->baselineMax) {
1259 cWpar->temax = i;
1260 cWpar->teflag = 1;
1261 break;
1262 }
1263 }
1264
1265 if(cWpar->temax > cWpar->tmax || cWpar->temax < cWpar->baseline_range || cWpar->temax < 0) {
1266 cWpar->teflag = 0;
1267 }
1268 // End of the determination of the upper limit of the exclusion zone.
1269
1270 // Here we determine the x position of the lower limit for the exclusion zone
1271 /***** Fitting the risetime *****/
1272 if(cWpar->teflag == 1) {
1273 // Set baselineMin
1275
1276 // Here we fit a line y=ax+b from temax to temax + 3*FILTER and find the intersection with baselineMin. The
1277 // x
1278 // coordinate of this intersection becomes temin.
1279 // Matrix for the fit
1280 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1281 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1282 memset(lineq_solution, 0, sizeof(lineq_solution)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1283
1284 lineq_dim = 2;
1285 for(i = cWpar->temax; i <= cWpar->temax + 3 * FILTER; i++) {
1286 lineq_matrix[0][0] += 1;
1287 lineq_matrix[0][1] += i;
1288 lineq_matrix[1][1] += i * i;
1289
1290 lineq_vector[0] += cWavebuffer[i];
1291
1292 lineq_vector[1] += cWavebuffer[i] * i;
1293 }
1294 lineq_matrix[1][0] = lineq_matrix[0][1];
1295
1296 // solve_lin_eq returns 0 if the determinant of the matrix is 0 the system is unsolvable. If there is no
1297 // solution, set temin to the upper limit of the baseline range.
1298 if(solve_lin_eq() == 0) {
1300 } else {
1303 // solve baselineMin = afit*x + bfit for x to find the crossing point. If the crossing point is outside
1304 // the
1305 // acceptable range, set temin to the upper limit of the baseline range.
1306 cWpar->temin = static_cast<int>(floor((cWpar->baselineMin - cWpar->bfit) / cWpar->afit));
1309 }
1310 if(cWpar->temin > cWpar->temax) {
1312 }
1313 }
1314 }
1315 // End of the determination of the lower limit of the exclusion zone.
1316 }
1317}
1318
1320{
1321
1322 if(i == 1) {
1323 return par->t[1];
1324 }
1325
1326 if(i >= 2) {
1327 if(i < NSHAPE) {
1328 return par->t[i] * par->t[1] / (par->t[i] + par->t[1]);
1329 }
1330 }
1331
1332 return -1.;
1333}
1334
1336{
1337 /*************************************************************************
1338 This function calculates t0 given the matrix solution from FitCsIShape.
1339 In this case, the fit function is written as follows:
1340
1341 f(t) = C + (Af+As)*exp(t0/tRC)*exp(-t/tRC) - Af*exp(t0/tF')*exp(-t/tF')
1342 - As*exp(t0/tS')*exp(-t/tS')
1343
1344 This can be re-written as:
1345
1346 f(t) = C' + alpha*exp(-t/tRC) + beta*exp(-t/tF') + gamma*exp(-t/tS')
1347
1348Where:
1349C = C'
1350alpha = (Af+As)*exp(t0/tRC)
1351beta = -Af*exp(t0/tF')
1352gamma = -As*exp(t0/tS')
1353
1354Ignoring the constant, we have:
1355
1356f'(t0) = alpha*exp(-t0/tRC) + beta*exp(-t0/tF') + gamma*exp(-t0/tS') = 0
1357
1358For t<t0, f'(t)< 0, and for t>t0, f'(t)>0. This function finds the
1359intersection of f'(t) and 0 by linear interpolation from these endpoints.
1360 *************************************************************************/
1361
1362 double ta = wpar->baseline_range;
1363 // ta=wpar->temin;
1364 double fa = 0.;
1365
1366 // t0 must be between the baseline and the max
1367 // calculates fit value (no constant) at the end of the baseline range
1368 // this is the t<t0 point
1369 for(int i = 1; i < lineq_dim; i++) {
1370 double tau = GetCsITau(i, par);
1371 // getc(stdin);
1372 fa += lineq_solution[i] * exp(-ta / tau);
1373 }
1374
1375 double tb = wpar->tmax;
1376 // tb=wpar->temax;
1377 double fb = 0.;
1378
1379 // calculates fit value (no constant) at tmax
1380 // this is the t>t0 point
1381 for(int i = 1; i < lineq_dim; i++) {
1382 double tau = GetCsITau(i, par);
1383 fb += lineq_solution[i] * exp(-tb / tau);
1384 }
1385
1386 double delta = 1.;
1387
1388 double tc = 0.; // need this later to set wpar->t0
1389
1390 if((fa < 0) && (fb > 0)) {
1391 // keep the interpolation going until you get below epsilon
1392 /* |f(t0) - 0| = |f(t0)|< epsilon */
1393 while(delta > EPS) {
1394 double slope = -fa / (fb - fa); // interpolation slope for dependent variable t
1395
1396 //"reasonable" interpolation slopes
1397 if(slope > 0.99) {
1398 slope = 0.99;
1399 }
1400
1401 if(slope < 0.01) {
1402 slope = 0.01;
1403 }
1404 // its pretty harmless computationally
1405
1406 // tc is the estimate for t0
1407 tc = ta + slope * (tb - ta);
1408 double fc = 0.;
1409 for(int i = 1; i < lineq_dim; i++) {
1410 double tau = GetCsITau(i, par);
1411 fc += lineq_solution[i] * exp(-tc / tau);
1412 }
1413
1414 // really should have this, just to be safe
1415 if(fc == 0) {
1416 return tc;
1417 }
1418
1419 if(fc > 0) {
1420 tb = tc;
1421 fb = fc;
1422 } else {
1423 ta = tc;
1424 fa = fc;
1425 }
1426 delta = std::abs(fc);
1427 }
1428 } else {
1429 return -1;
1430 }
1431
1432 // set wpar->t0 here
1433 wpar->t0 = tc;
1434
1435 return tc;
1436}
1437
1439{
1440 if(!IsSet()) { return; }
1441 delete cWpar;
1442 cWpar = new WaveFormPar;
1443 cWpar->baseline_range = T0RANGE; // default only 8 samples! but can be increased with a multiplier in TPulseAnalyzer constructor
1444 cWpar->t90_flag = 0;
1445 cWpar->t50_flag = 0;
1446 cWpar->t10_flag = 0;
1447 cWpar->mflag = 0;
1448 cWpar->bflag = false;
1449 cWpar->t0 = 0;
1450 cWpar->baselinefin = 0;
1451
1452 get_baseline(); // Takes a small sample baseline
1453 get_tmax(); // Does a filtered max search
1454
1455 if(!(cWpar->mflag == 1 && cWpar->bflag)) { return; }
1456
1457 if(cWpar->tmax > cN) { cWpar->tmax = cN - 1; }
1458
1459 double amp = cWpar->max - cWpar->baseline;
1460 double y9 = cWpar->baseline + amp * .9;
1461 double y5 = cWpar->baseline + amp * .5;
1462 double y1 = cWpar->baseline + amp * .1;
1463
1464 for(int t = cWpar->tmax; t > 0; t--) {
1465 if(cWavebuffer[t] < y5) {
1466 cWpar->t50_flag = 1;
1467 cWpar->t50 = t + 0.5;
1468 break;
1469 }
1470 }
1471
1472 if(cWpar->t50_flag != 1) { return; }
1473
1474 for(int t = cWpar->t50; t < cWpar->tmax; t++) {
1475 if(cWavebuffer[t] > y9) {
1476 cWpar->t90_flag = 1;
1477 cWpar->t90 = t - 0.5;
1478 break;
1479 }
1480 }
1481
1482 for(int t = cWpar->t50; t > 0; t--) {
1483 if(cWavebuffer[t] < y1) {
1484 cWpar->t10_flag = 1;
1485 cWpar->t10 = t + 0.5;
1486 break;
1487 }
1488 }
1489
1490 if(cWpar->t10_flag != 1) { return; }
1491
1492 double t0 = cWpar->t50 - ((cWpar->t50 - cWpar->t10) * 1.2);
1493 if(cWpar->t90_flag == 1) {
1494 t0 += cWpar->t90 - ((cWpar->t90 - cWpar->t10) * 1.125);
1495 t0 *= 0.5;
1496 }
1497 if(t0 < 0.) { t0 = 0.; }
1498
1499 // std::cout<<std::endl<<t0<<std::flush;
1500 cWpar->t0 = t0;
1501 get_baseline_fin(); // baseline all the way up to t0
1502}
1503
1504bool TPulseAnalyzer::SiliShapePrepare(double tauDecay, double tauRise)
1505{
1506 if(IsSet()) {
1507 // double t0=fit_newT0();//fits the T0 in the SFU way with my added bit at the end for a nice baseline calc
1508 // double t0=cWpar->t0;
1509 // int exclusion=t0+3;
1510 // if(t0<1){//if the fit_newT0() failed
1511 // exclusion=10;
1512 // if(abs(baseline-cWavebuffer[0])>100)baseline=cWavebuffer[0];
1513 // }
1514
1515 // New simplified guesses because fit_newT0 was taking 1000% longer
1516 GetQuickPara();
1517 cWpar->amplitude = 0;
1518 cWpar->tauDecay = tauDecay;
1519 cWpar->tauRise = tauRise;
1520 cWpar->bflag = false; // baseline
1521 cWpar->baseamp = 0;
1522 cWpar->basefreq = 0;
1523 cWpar->basephase = 0;
1524 cWpar->osciflag = 0;
1525
1526 if(cWpar->t10_flag == 0) { return false; }
1527
1528 // GetQuickPara() Returns values that are spurious if the baseline is missing or <<T0RANGE
1529 if(cWpar->t0 < cWpar->baseline_range) { // Is there is no clear baseline baseline
1530 if(!(cWpar->baselineStDevfin / (cWpar->max - cWpar->baselinefin) < 0.035)) { // Strict (previously 0.05) limit determined from data
1531 return false;
1532 }
1533 }
1534 cWpar->bflag = true;
1535 return true;
1536 }
1537 return false;
1538}
1539
1540bool TPulseAnalyzer::GetSiliShape(double tauDecay, double tauRise)
1541{
1542 if(IsSet()) {
1543
1544 if(!SiliShapePrepare(tauDecay, tauRise)) { return false; }
1545
1546 int exclusion = cWpar->t10;
1547 double baseline = cWpar->baselinefin;
1548
1549 /**************************************************************************
1550 // Parametes for this function
1551 //fShpar->t0 t0 (time where signal starts, calculated)
1552 //fShpar->tau[0]) decay (provided)
1553 //fShpar->tau[1]) rise (provided)
1554 //fShpar->tau[2]) slow (not used)
1555 //fShpar->tau[3]) diode (not used)
1556 //fShpar->am[0]) baseline (provided)
1557 //fShpar->am[1]) fast (amplitude, calculated )
1558 //fShpar->am[2]) slow (not used)
1559 //fShpar->am[3]) diode (not used)
1560
1561 linearized chi square fit is Mu = v where M is a data matrix
1562 u, v are vectors; u is the parameter vector (solution)
1563 note that in this formulation, chisq_min = y_i^2-sum(u_iv_i)
1564 **************************************************************************/
1565
1566 // cout<<baseline<<" "<<exclusion<<endl ;// ResetShapeAmplitudes();/ fShpar->am[0] = baseline ;
1567
1568 lineq_dim = 2;
1569 memset(lineq_matrix, 0, sizeof(lineq_matrix)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1570 memset(lineq_vector, 0, sizeof(lineq_vector)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1571 memset(lineq_solution, 0, sizeof(lineq_solution)); // NOLINT(cppcoreguidelines-pro-bounds-array-to-pointer-decay)
1572
1573 if(exclusion >= cN) { return false; }
1574 if(lineq_dim >= cN) { return false; }
1575
1576 // setting M[0,0] V[0] V[1]
1577 for(int j = exclusion; j < cN; j++) {
1578 // vector
1579 if((cWavebuffer[j] - baseline) < 0) {
1580 exclusion++;
1581 continue;
1582 } // this is crucial for oscillations
1583 // if (j%10==0) cout<<j<<" "<< cWavebuffer[j]<<" "<<cWavebuffer[j]-baseline<<endl ;
1584 double signal = log(cWavebuffer[j] - baseline) + j / tauDecay; // sum of Y_i where Y_i = (y_i - baseline)*exp(-t_i/tauDecay)
1585 lineq_vector[0] += exp(signal);
1586 lineq_vector[1] -= exp(signal - j / tauRise); // sum of Y_i*X_i
1587 // Matrix
1588 lineq_matrix[0][0] += 1;
1589 }
1590
1591 // create matrix for linearized fit
1592 // setting elements M[0,1] M[1,0] M[1,1]
1593 long double sum = -static_cast<double>(exclusion) / tauRise + log(1. - exp(-(static_cast<double>(cN - exclusion)) / tauRise));
1594 sum -= log(1. - exp(-1. / tauRise)); // finishing the geometric sequence sum
1595 lineq_matrix[1][0] = -exp(sum);
1596 lineq_matrix[0][1] = -exp(sum);
1597
1598 double tauRise_2 = tauRise / 2.;
1599 sum = -static_cast<double>(exclusion) / tauRise_2 + log(1. - exp(-(static_cast<double>(cN - exclusion)) / tauRise_2));
1600 sum -= log(1. - exp(-1. / tauRise_2));
1601 lineq_matrix[1][1] = exp(sum);
1602
1603 // cout<<lineq_matrix[0][0]<<" "<<lineq_matrix[0][1]<<" ---------- "<<lineq_vector[0] <<endl ;
1604 // cout<<lineq_matrix[1][0]<<" "<<lineq_matrix[1][1]<<" ---------- "<<lineq_vector[1] <<endl ;
1605
1606 // solve the matrix equation Mu = v -> u = M^(-1)v where M^(-1) is the inverse
1607 // of M. note this has no solution if M is not invertable!
1608
1609 // error if the matrix cannot be inverted
1610 if(solve_lin_eq() == 0) {
1611 return false;
1612 }
1613 // calculate amplitudes
1614 double beta = lineq_solution[0];
1615 double alpha = lineq_solution[1];
1616
1617 double dom = exp(((log(alpha) - log(beta)) * tauRise) / tauDecay);
1618 if(dom > 0 || dom < 0) { cWpar->amplitude = beta / dom; }
1619
1620 double tt = (log(alpha) - log(beta)) * tauRise;
1621 if(tt > 0) { cWpar->t0 = tt; }
1622
1623 return true;
1624 }
1625 return false;
1626}
1627
1628// Significantly slower and should only be used in non-sorting analysis of poor waveform
1629// Needs initial estimates even if fitting those parameters
1630// Setting basefreq>0 opens a very experimental/bad mode
1631bool TPulseAnalyzer::GetSiliShapeTF1(double tauDecay, double tauRise, double baseline, double basefreq)
1632{
1633 TGraph* h = GetWaveGraph();
1634 if(h != nullptr) { // Graph better than hist for stats and simplicity
1635 SiliShapePrepare(tauDecay, tauRise);
1636 TF1 g = Getsilifit();
1637
1638 // g.SetParameter(0,cWpar->t0);
1639 // g.SetParameter(1,cWpar->tauDecay);
1640 // g.SetParameter(2,cWpar->tauRise);
1641 // g.SetParameter(3,cWpar->baselinefin);
1642 // g.SetParameter(4,cWpar->amplitude);
1643
1644 // Currently constrained for positive waveforms
1645 double r = cWpar->max - cWavebuffer[0];
1646 g.SetParameter(4, r * 1.05);
1647
1648 if(cWpar->bflag) { // Have reasonable T0 & base, fit the shape
1649 g.SetParLimits(0, cWpar->t0 * 0.5, cWpar->t0 * 1.5);
1650 g.SetParLimits(1, tauDecay * 0.3, tauDecay * 3.0);
1651 g.SetParLimits(2, tauRise * 0.3, tauRise * 1.5);
1652 g.FixParameter(3, cWpar->baselinefin);
1653 g.SetParLimits(4, r * 0.5, r * 2.0);
1654 } else { // Have no T0 or base, FIX the shape
1655 g.SetParLimits(0, -cN, cWpar->baseline_range);
1656 g.FixParameter(1, tauDecay);
1657 g.FixParameter(2, tauRise);
1658 g.FixParameter(3, baseline);
1659 // g.SetParameter(3,baseline);
1660 // g.SetParLimits(3,baseline-300,baseline+300);
1661 g.SetParLimits(4, r * 0.1, r * 10.0);
1662 }
1663
1664 if(basefreq > 0) {
1665 cWpar->osciflag = 1;
1666 g.ReleaseParameter(5);
1667 g.ReleaseParameter(6);
1668 g.ReleaseParameter(7);
1669 g.SetParameter(5, basefreq);
1670 g.SetParameter(6, 0.5);
1671 g.SetParameter(7, r * 0.1);
1672 g.SetParLimits(5, basefreq * 0.5, basefreq * 2.0);
1673 g.SetParLimits(6, 0, 1);
1674 g.SetParLimits(7, 0, r * 10.0);
1675 }
1676
1677 int res = h->Fit(&g, "QN");
1678 delete h;
1679
1680 if(res == 0) {
1681 cWpar->t0 = g.GetParameter(0);
1682 cWpar->tauDecay = g.GetParameter(1);
1683 cWpar->tauRise = g.GetParameter(2);
1684 cWpar->baselinefin = g.GetParameter(3);
1685 cWpar->amplitude = g.GetParameter(4);
1686 cWpar->basefreq = g.GetParameter(5);
1687 cWpar->basephase = g.GetParameter(6);
1688 cWpar->baseamp = g.GetParameter(7);
1689
1690 return true;
1691 }
1692 }
1693 return false;
1694}
1695
1696double TPulseAnalyzer::SiLiFitFunction(double* i, double* p) // NOLINT(readability-non-const-parameter)
1697{
1698 // p[0]-p[2] are t0, RC, Tau
1699 // p[3]-p[4] are baseline, A0
1700 // p[5]-p[7] are osci freq,phase,amp
1701
1702 double x = i[0] - p[0];
1703
1704 double s = p[3];
1705 if(x > 0) { s += p[4] * (1 - exp(-x / p[2])) * exp(-x / p[1]); }
1706 if(p[7] > 0) { s += p[7] * sin((p[6] + i[0] / p[5]) * 2 * TMath::Pi()); }
1707
1708 return s;
1709}
1710
1712{
1713 if(set && cWpar != nullptr) {
1714 std::ostringstream name;
1715 name << "Fit" << fNameIter;
1716 ++fNameIter;
1717 TF1 g(name.str().c_str(), SiLiFitFunction, 0, cN, 8);
1718
1719 g.SetParameter(0, cWpar->t0);
1720 g.SetParameter(1, cWpar->tauDecay);
1721 g.SetParameter(2, cWpar->tauRise);
1722 g.SetParameter(3, cWpar->baselinefin);
1723 g.SetParameter(4, cWpar->amplitude);
1724 g.FixParameter(5, cWpar->basefreq);
1725 g.FixParameter(6, cWpar->basephase);
1726 g.FixParameter(7, cWpar->baseamp);
1727
1728 g.SetLineColor(kRed);
1729
1730 return g;
1731 }
1732
1733 return {};
1734}
1735
1737{
1738 double Smirnov = 0;
1739 if(set && cWpar != nullptr) {
1740 double gsum = 0.;
1741 double wsum = 0.;
1742 TF1 g = Getsilifit();
1743 for(Int_t i = 0; i < cN; i++) {
1744 wsum += abs(cWavebuffer[i]);
1745 gsum += abs(g.Eval(i + 0.5));
1746 Smirnov += abs(wsum - gsum);
1747 }
1748 }
1749 return Smirnov;
1750}
1751
1753{
1754 if(!set) { return; }
1755 DrawWave();
1756 if(cWpar != nullptr) {
1757 Getsilifit().DrawCopy("same");
1758 std::cout << "t0:\t" << cWpar->t0 << ", A:\t" << cWpar->amplitude << std::endl;
1759 }
1760}
1761
1763{
1764 TH1I* h = GetWaveHist();
1765 if(h != nullptr) {
1766 h->DrawCopy();
1767 delete h;
1768 }
1769}
1770
1773{
1774 if(cN == 0 || !set) {
1775 return nullptr;
1776 }
1777 std::ostringstream name;
1778 name << "WaveformHist" << fNameIter;
1779 ++fNameIter; // Avoid naming conflicts with TNamed
1780 TH1I* h = new TH1I(name.str().c_str(), name.str().c_str(), cN, -0.5, cN - 0.5); // midpoint should be the value, else time is off
1781 for(Int_t i = 0; i < cN; i++) {
1782 h->SetBinContent(i + 1, cWavebuffer[i]);
1783 }
1784 return h;
1785}
1786
1788{
1789 if(cN == 0 || !set) {
1790 return nullptr;
1791 }
1792 auto* g = new TGraph();
1793 for(Int_t i = 0; i < cN; i++) {
1794 g->SetPoint(i, i, cWavebuffer[i]);
1795 }
1796 return g;
1797}
1798
1800{
1801 if(cN == 0 || !set) {
1802 return;
1803 }
1804
1805 DrawWave();
1806
1807 if(spar != nullptr) {
1808 TF1 f("fit", "[2] + [0]*sin(6.283185307 * (x - [1])/(2*8.48409))", 0, cN);
1809
1810 f.SetParameter(0, spar->A);
1811 f.SetParameter(1, spar->t0);
1812 f.SetParameter(2, spar->C);
1813
1814 f.DrawCopy("same");
1815
1816 std::cout << "t0:\t" << spar->t0 << ", A:\t" << spar->A << ", O:\t" << spar->C << std::endl;
1817 }
1818}
1819
1821{
1822
1823 if(cN == 0 || !set) {
1824 return;
1825 }
1826
1827 DrawWave();
1828
1829 if(cWpar != nullptr) {
1830 TF1 g("fit", "[0]+[1]*x", 0, cWpar->temax);
1831 TF1 f("fit", "[0]+[1]*x+[2]*x*x", cWpar->temin, cWpar->thigh);
1832
1833 g.SetParameter(0, cWpar->b0);
1834 g.SetParameter(1, cWpar->b1);
1835 g.SetLineColor(kRed);
1836
1837 f.SetParameter(0, cWpar->s0);
1838 f.SetParameter(1, cWpar->s1);
1839 f.SetParameter(2, cWpar->s2);
1840 f.SetLineColor(8);
1841
1842 f.DrawCopy("same");
1843 g.DrawCopy("same");
1844
1845 std::cout << "t0:\t" << cWpar->t0 << std::endl;
1846 }
1847}
1848
1850{
1851
1852 if(cN == 0 || !set) {
1853 return;
1854 }
1855
1856 DrawWave();
1857 if(cWpar != nullptr) {
1858 TF1 f("base", "[0]", 0, cWpar->baseline_range);
1859 TF1 g("basemin", "[0]", cWpar->temin, cWpar->temax);
1860 TF1 h("basemax", "[0]", cWpar->temin, cWpar->temax);
1861 TF1 r("risetime", "[0]*x+[1]", cWpar->temin, cWpar->temax + 3 * FILTER);
1862
1863 std::cout << "Baseline:\t" << cWpar->baseline << std::endl;
1864 std::cout << "Zero crossing:\t" << cWpar->t0 << std::endl;
1865
1866 f.SetParameter(0, cWpar->baseline);
1867 f.SetLineColor(kGreen);
1868
1869 g.SetParameter(0, cWpar->baselineMin);
1870 g.SetLineColor(kBlue);
1871
1872 h.SetParameter(0, cWpar->baselineMax);
1873 h.SetLineColor(kBlack);
1874
1875 r.SetParameter(0, cWpar->afit);
1876 r.SetParameter(1, cWpar->bfit);
1877 r.SetLineColor(kRed);
1878
1879 f.DrawCopy("same");
1880 g.DrawCopy("same");
1881 h.DrawCopy("same");
1882 r.DrawCopy("same");
1883 }
1884}
1885
1887{
1888
1889 if(cN == 0 || !set || !CsIIsSet()) {
1890 return;
1891 }
1892
1893 DrawWave();
1894 if(cWpar != nullptr) {
1895 TF1 shape("shape", TGRSIFunctions::CsIFitFunction, 0, cN, 9);
1896
1897 shape.SetParameter(0, shpar->t[0]);
1898 shape.SetParameter(1, shpar->t[1]);
1899 shape.SetParameter(2, shpar->t[2]);
1900 shape.SetParameter(3, shpar->t[3]);
1901 shape.SetParameter(4, shpar->t[4]);
1902 shape.SetParameter(5, shpar->am[0]);
1903 shape.SetParameter(6, shpar->am[2]);
1904 shape.SetParameter(7, shpar->am[3]);
1905 shape.SetParameter(8, shpar->am[4]);
1906 shape.SetLineColor(kRed);
1907
1908 std::cout << "t0:\t" << shpar->t[0] << ",\ttRC:\t" << shpar->t[1] << ",\ttF:\t" << shpar->t[2] << ",\ttS:\t" << shpar->t[3] << ",\tTGamma:\t" << shpar->t[4] << std::endl;
1909 std::cout << "Baseline:\t" << shpar->am[0] << ",\tFast:\t" << shpar->am[2] << ",\tSlow:\t" << shpar->am[3] << ",\tGamma:\t" << shpar->am[4] << std::endl;
1910
1911 shape.DrawCopy("same");
1912 }
1913}
1914// NOLINTEND(cppcoreguidelines-narrowing-conversions)
1915
1916/*======================================================*/
1918{
1919 std::cout << "== Currently established waveform parameters ============" << std::endl;
1920 std::cout << "baseline : " << std::setw(10) << cWpar->baseline << std::endl;
1921 std::cout << "baseline st. dev.: " << std::setw(10) << cWpar->baselineStDev << std::endl;
1922 std::cout << "max : " << std::setw(10) << cWpar->max << std::endl;
1923 std::cout << "tmax : " << std::setw(10) << cWpar->tmax << std::endl;
1924 std::cout << "temin : " << std::setw(10) << static_cast<double>(cWpar->temin) << std::endl;
1925 std::cout << "temax : " << std::setw(10) << static_cast<double>(cWpar->temax) << std::endl;
1926 std::cout << "t0 : " << std::setw(10) << cWpar->t0 << std::endl;
1927}
const Double_t s
const std::vector< Short_t > * GetWaveform() const
!
virtual bool HasWave() const
!
long double lineq_solution[20]
std::vector< Short_t > cWavebuffer
double fit_rf(double=2 *8.48409)
long double determinant(int)
bool CsIIsSet() const
static const int BADCHISQ_PAR_T0
static const int BADCHISQ_LIN_T0
bool GetSiliShapeTF1(double tauDecay, double tauRise, double baseline, double basefreq=0)
static const int CSI_BASELINE_RANGE
static const int BAD_EXCLUSION_ZONE
double GetCsIt0(ShapePar *, WaveFormPar *)
double get_tfrac(double, double, double)
static const int NSHAPE
void SetCsI(bool option=true)
static const int BADCHISQ_FAIL_DIRECT
void Clear(Option_t *opt="")
bool IsSet() const
double GetCsITau(int, ShapePar *)
static const int BAD_BASELINE_RANGE
static double SiLiFitFunction(double *i, double *p)
long double copy_matrix[20][20]
ShapePar * csiTestShpar[4]
static const int BADCHISQ_SMOOTH_T0
long double lineq_vector[20]
bool SiliShapePrepare(double tauDecay, double tauRise)
int fit_parabola(int, int, ParPar *)
int fit_line(int, int, LinePar *)
WaveFormPar * cWpar
static int fNameIter
static const int BADCHISQ_AMPL
int fit_smooth_parabola(int, int, double, ParPar *)
static const int MAX_SAMPLES
static const int BADCHISQ_MAT
static const int BAD_MAX
WaveFormPar * csiTestWpar[4]
int16_t good_baseline()
long double lineq_matrix[20][20]
virtual ~TPulseAnalyzer()
int FitCsIShape(int, ShapePar *, WaveFormPar *)
double get_parabolic_T0()
static const int NOISE_LEVEL_CSI
static const int BADCHISQ_NEG
static const int BADCHISQ_T0
void SetData(const TFragment &fragment, double=0)
TGraph * GetWaveGraph()
bool GetSiliShape(double tauDecay, double tauRise)
double get_sin_par(double)
static const int BAD_PID
Double_t CsIFitFunction(Double_t *time, Double_t *par)