REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Loading...
Searching...
No Matches
TRestDetectorSignal.cxx
1
21// dec 2015:
22//
23// Javier Galan
25
26#include "TRestDetectorSignal.h"
27
28#include <TCanvas.h>
29#include <TF1.h>
30#include <TFitResult.h>
31#include <TH1.h>
32#include <TMath.h>
33#include <TRandom3.h>
34
35#include <limits>
36
37using namespace std;
38
39ClassImp(TRestDetectorSignal);
40
41TRestDetectorSignal::TRestDetectorSignal() {
42 // TRestDetectorSignal default constructor
43 fGraph = nullptr;
44 fSignalID = -1;
45 fSignalTime.clear();
46 fSignalCharge.clear();
47
48 fPointsOverThreshold.clear();
49}
50
51TRestDetectorSignal::~TRestDetectorSignal() {
52 // TRestDetectorSignal destructor
53}
54
55void TRestDetectorSignal::NewPoint(Float_t time, Float_t data) {
56 fSignalTime.push_back(time);
57 fSignalCharge.push_back(data);
58}
59
64void TRestDetectorSignal::IncreaseAmplitude(Double_t t, Double_t d) {
65 TVector2 p(t, d);
67}
68
76 Int_t index = GetTimeIndex(p.X());
77 Float_t x = p.X();
78 Float_t y = p.Y();
79
80 if (index >= 0) {
81 fSignalTime[index] = x;
82 fSignalCharge[index] += y;
83 } else {
84 fSignalTime.push_back(x);
85 fSignalCharge.push_back(y);
86 }
87}
88
97 Int_t index = GetTimeIndex(p.X());
98 Float_t x = p.X();
99 Float_t y = p.Y();
100
101 if (index >= 0) {
102 fSignalTime[index] = x;
103 fSignalCharge[index] = y;
104 } else {
105 fSignalTime.push_back(x);
106 fSignalCharge.push_back(y);
107 }
108}
109
117void TRestDetectorSignal::SetPoint(Double_t t, Double_t d) {
118 TVector2 p(t, d);
119 SetPoint(p);
120}
121
126void TRestDetectorSignal::SetPoint(Int_t index, Double_t t, Double_t d) {
127 fSignalTime[index] = t;
128 fSignalCharge[index] = d;
129}
130
131Double_t TRestDetectorSignal::GetIntegral(Int_t startBin, Int_t endBin) const {
132 if (startBin < 0) {
133 startBin = 0;
134 }
135 if (endBin <= 0 || endBin > GetNumberOfPoints()) {
136 endBin = GetNumberOfPoints();
137 }
138
139 Double_t sum = 0;
140 for (int i = startBin; i < endBin; i++) {
141 sum += GetData(i);
142 }
143
144 return sum;
145}
146
147void TRestDetectorSignal::Normalize(Double_t scale) {
148 Double_t sum = GetIntegral();
149
150 for (int i = 0; i < GetNumberOfPoints(); i++) fSignalCharge[i] = scale * GetData(i) / sum;
151}
152
153Double_t TRestDetectorSignal::GetIntegralWithTime(Double_t startTime, Double_t endTime) const {
154 Double_t sum = 0;
155 for (int i = 0; i < GetNumberOfPoints(); i++) {
156 if (GetTime(i) >= startTime && GetTime(i) < endTime) {
157 sum += GetData(i);
158 }
159 }
160
161 return sum;
162}
163
164Double_t TRestDetectorSignal::GetMaxPeakWithTime(Double_t startTime, Double_t endTime) {
165 Double_t max = std::numeric_limits<double>::min();
166
167 for (int i = 0; i < GetNumberOfPoints(); i++)
168 if (GetTime(i) >= startTime && GetTime(i) < endTime) {
169 if (this->GetData(i) > max) max = GetData(i);
170 }
171
172 return max;
173}
174
175/* {{{
176Double_t TRestDetectorSignal::GetIntegralWithThreshold(Int_t from, Int_t to, Int_t startBaseline, Int_t
177endBaseline,
178 Double_t nSigmas, Int_t nPointsOverThreshold,
179 Double_t nMinSigmas) {
180 if (startBaseline < 0) startBaseline = 0;
181 if (endBaseline <= 0 || endBaseline > GetNumberOfPoints()) endBaseline = GetNumberOfPoints();
182
183 Double_t baseLine = GetBaseLine(startBaseline, endBaseline);
184
185 Double_t pointThreshold = nSigmas * GetBaseLineSigma(startBaseline, endBaseline);
186 Double_t signalThreshold = nMinSigmas * GetBaseLineSigma(startBaseline, endBaseline);
187
188 return GetIntegralWithThreshold(from, to, baseLine, pointThreshold, nPointsOverThreshold,
189 signalThreshold);
190}
191
192Double_t TRestDetectorSignal::GetIntegralWithThreshold(Int_t from, Int_t to, Double_t baseline,
193 Double_t pointThreshold, Int_t nPointsOverThreshold,
194 Double_t signalThreshold) {
195 Double_t sum = 0;
196 Int_t nPoints = 0;
197 fPointsOverThreshold.clear();
198
199 if (to > GetNumberOfPoints()) to = GetNumberOfPoints();
200
201 Float_t maxValue = 0;
202 for (int i = from; i < to; i++) {
203 if (GetData(i) > baseline + pointThreshold) {
204 if (GetData(i) > maxValue) maxValue = GetData(i);
205 nPoints++;
206 } else {
207 if (nPoints >= nPointsOverThreshold) {
208 Double_t sig = GetStandardDeviation(i - nPoints, i);
209 if (sig > signalThreshold) {
210 for (int j = i - nPoints; j < i; j++) {
211 sum += this->GetData(j);
212 fPointsOverThreshold.push_back(j);
213 }
214 }
215 }
216 nPoints = 0;
217 maxValue = 0;
218 }
219 }
220
221 if (nPoints >= nPointsOverThreshold) {
222 Double_t sig = GetStandardDeviation(to - nPoints, to);
223 if (sig > signalThreshold) {
224 for (int j = to - nPoints; j < to; j++) {
225 sum += this->GetData(j);
226 fPointsOverThreshold.push_back(j);
227 }
228 }
229 }
230
231 return sum;
232}
233}}} */
234
235Double_t TRestDetectorSignal::GetAverage(Int_t start, Int_t end) {
236 this->Sort();
237
238 if (end == 0) end = this->GetNumberOfPoints();
239
240 Double_t sum = 0;
241 for (int i = start; i <= end; i++) {
242 sum += this->GetData(i);
243 }
244 return sum / (end - start + 1);
245}
246
247Int_t TRestDetectorSignal::GetMaxPeakWidth() {
248 this->Sort();
249
250 Int_t mIndex = this->GetMaxIndex();
251 Double_t maxValue = this->GetData(mIndex);
252
253 Double_t value = maxValue;
254 Int_t rightIndex = mIndex;
255 while (value > maxValue / 2) {
256 value = this->GetData(rightIndex);
257 rightIndex++;
258 }
259 Int_t leftIndex = mIndex;
260 value = maxValue;
261 while (value > maxValue / 2) {
262 value = this->GetData(leftIndex);
263 leftIndex--;
264 }
265
266 return rightIndex - leftIndex;
267}
268
269Double_t TRestDetectorSignal::GetMaxPeakValue() { return GetData(GetMaxIndex()); }
270
271Int_t TRestDetectorSignal::GetMaxIndex(Int_t from, Int_t to) {
272 Double_t max = std::numeric_limits<double>::min();
273 Int_t index = 0;
274
275 if (from < 0) from = 0;
276 if (to > GetNumberOfPoints()) to = GetNumberOfPoints();
277
278 if (to == 0) to = GetNumberOfPoints();
279
280 for (int i = from; i < to; i++) {
281 if (this->GetData(i) > max) {
282 max = GetData(i);
283 index = i;
284 }
285 }
286
287 return index;
288}
289
290// z position by gaussian fit
291
292TVector2
293TRestDetectorSignal::GetMaxGauss() // returns a 2vector with the time of the peak time in us and the energy
294{
295 Int_t maxRaw = GetMaxIndex(); // The bin where the maximum of the raw signal is found
296 Double_t maxRawTime =
297 GetTime(maxRaw); // The time of the bin where the maximum of the raw signal is found
298 Double_t energy = 0, time = 0;
299 Double_t lowerLimit = maxRawTime - 0.2; // us
300 Double_t upperLimit = maxRawTime + 0.4; // us
301
302 TF1* gaus = new TF1("gaus", "gaus", lowerLimit, upperLimit);
303 TH1F* h1 = new TH1F("h1", "h1", 1000, 0,
304 10); // Histogram to store the signal. For now the number of bins is fixed.
305
306 // copying the signal peak to a histogram
307 for (int i = 0; i < GetNumberOfPoints(); i++) {
308 h1->Fill(GetTime(i), GetData(i));
309 }
310 /*
311 TCanvas* c = new TCanvas("c", "Signal fit", 200, 10, 1280, 720);
312 h1->GetXaxis()->SetTitle("Time (us)");
313 h1->GetYaxis()->SetTitle("Amplitude");
314 h1->Draw();
315 */
316
317 TFitResultPtr fitResult =
318 h1->Fit(gaus, "QNRS"); // Q = quiet, no info in screen; N = no plot; R = fit in the function range; S
319 // = save and return the fit result
320
321 if (fitResult->IsValid()) {
322 energy = gaus->GetParameter(0);
323 time = gaus->GetParameter(1);
324 } else {
325 // the fit failed, return -1 to indicate failure
326 energy = -1;
327 time = -1;
328 cout << endl
329 << "WARNING: bad fit to signal with ID " << GetID() << " with maximum at time = " << maxRawTime
330 << " ns "
331 << "\n"
332 << "Failed fit parameters = " << gaus->GetParameter(0) << " || " << gaus->GetParameter(1)
333 << " || " << gaus->GetParameter(2) << "\n"
334 << "Assigned fit parameters : energy = " << energy << ", time = " << time << endl;
335 /*
336 TCanvas* c2 = new TCanvas("c2", "Signal fit", 200, 10, 1280, 720);
337 h1->Draw();
338 c2->Update();
339 getchar();
340 delete c2;
341 */
342 }
343
344 TVector2 fitParam(time, energy);
345
346 delete h1;
347 delete gaus;
348
349 return fitParam;
350}
351
352// z position by landau fit
353
354TVector2
355TRestDetectorSignal::GetMaxLandau() // returns a 2vector with the time of the peak time in us and the energy
356{
357 Int_t maxRaw = GetMaxIndex(); // The bin where the maximum of the raw signal is found
358 Double_t maxRawTime =
359 GetTime(maxRaw); // The time of the bin where the maximum of the raw signal is found
360 Double_t energy = 0, time = 0;
361 Double_t lowerLimit = maxRawTime - 0.2; // us
362 Double_t upperLimit = maxRawTime + 0.4; // us
363
364 TF1* landau = new TF1("landau", "landau", lowerLimit, upperLimit);
365 TH1F* h1 = new TH1F("h1", "h1", 1000, 0,
366 10); // Histogram to store the signal. For now the number of bins is fixed.
367
368 // copying the signal peak to a histogram
369 for (int i = 0; i < GetNumberOfPoints(); i++) {
370 h1->Fill(GetTime(i), GetData(i));
371 }
372
373 TFitResultPtr fitResult =
374 h1->Fit(landau, "QNRS"); // Q = quiet, no info in screen; N = no plot; R = fit in the function range;
375 // S = save and return the fit result
376 if (fitResult->IsValid()) {
377 energy = landau->GetParameter(0);
378 time = landau->GetParameter(1);
379 } else {
380 // the fit failed, return -1 to indicate failure
381 energy = -1;
382 time = -1;
383 cout << endl
384 << "WARNING: bad fit to signal with ID " << GetID() << " with maximum at time = " << maxRawTime
385 << " ns "
386 << "\n"
387 << "Failed fit parameters = " << landau->GetParameter(0) << " || " << landau->GetParameter(1)
388 << " || " << landau->GetParameter(2) << "\n"
389 << "Assigned fit parameters : energy = " << energy << ", time = " << time << endl;
390 /*
391 TCanvas* c2 = new TCanvas("c2", "Signal fit", 200, 10, 1280, 720);
392 h1->Draw();
393 c2->Update();
394 getchar();
395 delete c2;
396 */
397 }
398
399 TVector2 fitParam(time, energy);
400
401 delete h1;
402 delete landau;
403
404 return fitParam;
405}
406
407// z position by aget fit
408
409Double_t agetResponseFunction(Double_t* x, Double_t* par) { // x contains as many elements as the number of
410 // dimensions (in this case 1, i.e. x[0]), and
411 // par contains as many elements as the number of free parameters in my function.
412
413 Double_t arg =
414 (x[0] - par[1] + 1.1664) /
415 par[2]; // 1.1664 is the x value where the maximum of the base function (i.e. without parameters) is.
416 Double_t f = par[0] / 0.0440895 * exp(-3 * (arg)) * (arg) * (arg) *
417 (arg)*sin(arg); // to rescale the Y axis and get amplitude.
418 return f;
419}
420
421TVector2
422TRestDetectorSignal::GetMaxAget() // returns a 2vector with the time of the peak time in us and the energy
423{
424 Int_t maxRaw = GetMaxIndex(); // The bin where the maximum of the raw signal is found
425 Double_t maxRawTime =
426 GetTime(maxRaw); // The time of the bin where the maximum of the raw signal is found
427 Double_t energy = 0, time = 0;
428 // The intervals below are small because otherwise the function doesn't fit anymore.
429 Double_t lowerLimit = maxRawTime - 0.2; // us
430 Double_t upperLimit = maxRawTime + 0.7; // us
431
432 TF1* aget = new TF1("aget", agetResponseFunction, lowerLimit, upperLimit, 3); //
433 TH1F* h1 = new TH1F("h1", "h1", 1000, 0,
434 10); // Histogram to store the signal. For now the number of bins is fixed.
435 aget->SetParameters(500, maxRawTime, 1.2);
436
437 // copying the signal peak to a histogram
438 for (int i = 0; i < GetNumberOfPoints(); i++) {
439 h1->Fill(GetTime(i), GetData(i));
440 }
441
442 TFitResultPtr fitResult =
443 h1->Fit(aget, "QNRS"); // Q = quiet, no info in screen; N = no plot; R = fit in
444 // the function range; S = save and return the fit result
445
446 if (fitResult->IsValid()) {
447 energy = aget->GetParameter(0);
448 time = aget->GetParameter(1);
449 } else {
450 // the fit failed, return -1 to indicate failure
451 energy = -1;
452 time = -1;
453 cout << endl
454 << "WARNING: bad fit to signal with ID " << GetID() << " with maximum at time = " << maxRawTime
455 << " ns "
456 << "\n"
457 << "Failed fit parameters = " << aget->GetParameter(0) << " || " << aget->GetParameter(1)
458 << " || " << aget->GetParameter(2) << "\n"
459 << "Assigned fit parameters : energy = " << energy << ", time = " << time << endl;
460 }
461
462 TVector2 fitParam(time, energy);
463
464 delete h1;
465 delete aget;
466
467 return fitParam;
468}
469Double_t TRestDetectorSignal::GetMaxPeakTime(Int_t from, Int_t to) { return GetTime(GetMaxIndex(from, to)); }
470
471Double_t TRestDetectorSignal::GetMinPeakValue() { return GetData(GetMinIndex()); }
472
473Int_t TRestDetectorSignal::GetMinIndex() {
474 Double_t min = numeric_limits<double>::max();
475 Int_t index = 0;
476
477 for (int i = 0; i < GetNumberOfPoints(); i++) {
478 if (this->GetData(i) < min) {
479 min = GetData(i);
480 index = i;
481 }
482 }
483
484 return index;
485}
486
487Double_t TRestDetectorSignal::GetMinTime() const {
488 if (GetNumberOfPoints() == 0) {
489 return 0;
490 }
491 Double_t minTime = numeric_limits<float>::max();
492 for (int n = 0; n < GetNumberOfPoints(); n++) {
493 if (minTime > fSignalTime[n]) {
494 minTime = fSignalTime[n];
495 }
496 }
497
498 return minTime;
499}
500
501Double_t TRestDetectorSignal::GetMaxTime() const {
502 if (GetNumberOfPoints() == 0) {
503 return 0;
504 }
505 Double_t maxTime = numeric_limits<float>::min();
506 for (int n = 0; n < GetNumberOfPoints(); n++) {
507 if (maxTime < fSignalTime[n]) {
508 maxTime = fSignalTime[n];
509 }
510 }
511 return maxTime;
512}
513
514Int_t TRestDetectorSignal::GetTimeIndex(Double_t t) {
515 Float_t time = t;
516
517 for (int n = 0; n < GetNumberOfPoints(); n++)
518 if (time == fSignalTime[n]) return n;
519 return -1;
520}
521
522Bool_t TRestDetectorSignal::isSorted() {
523 for (int i = 0; i < GetNumberOfPoints() - 1; i++) {
524 if (GetTime(i + 1) < GetTime(i)) return false;
525 }
526 return true;
527}
528
529void TRestDetectorSignal::Sort() {
530 while (!isSorted()) {
531 for (int i = 0; i < GetNumberOfPoints(); i++) {
532 for (int j = i; j < GetNumberOfPoints(); j++) {
533 if (GetTime(i) > GetTime(j)) {
534 iter_swap(fSignalTime.begin() + i, fSignalTime.begin() + j);
535 iter_swap(fSignalCharge.begin() + i, fSignalCharge.begin() + j);
536 }
537 }
538 }
539 }
540}
541
542void TRestDetectorSignal::GetDifferentialSignal(TRestDetectorSignal* diffSignal, Int_t smearPoints) {
543 this->Sort();
544
545 for (int i = 0; i < smearPoints; i++) diffSignal->IncreaseAmplitude(GetTime(i), 0);
546
547 for (int i = smearPoints; i < this->GetNumberOfPoints() - smearPoints; i++) {
548 Double_t value = (this->GetData(i + smearPoints) - GetData(i - smearPoints)) /
549 (GetTime(i + smearPoints) - GetTime(i - smearPoints));
550 Double_t time = (GetTime(i + smearPoints) + GetTime(i - smearPoints)) / 2.;
551
552 diffSignal->IncreaseAmplitude(time, value);
553 }
554
555 for (int i = GetNumberOfPoints() - smearPoints; i < GetNumberOfPoints(); i++) {
556 diffSignal->IncreaseAmplitude(GetTime(i), 0);
557 }
558}
559
560void TRestDetectorSignal::GetSignalDelayed(TRestDetectorSignal* delayedSignal, Int_t delay) {
561 this->Sort();
562
563 for (int i = 0; i < delay; i++) delayedSignal->IncreaseAmplitude(GetTime(i), GetData(i));
564
565 for (int i = delay; i < GetNumberOfPoints(); i++) {
566 delayedSignal->IncreaseAmplitude(GetTime(i), GetData(i - delay));
567 }
568}
569
570void TRestDetectorSignal::GetSignalSmoothed(TRestDetectorSignal* smthSignal, Int_t averagingPoints) {
571 this->Sort();
572
573 averagingPoints = (averagingPoints / 2) * 2 + 1; // make it odd >= averagingPoints
574
575 Double_t sum = GetIntegral(0, averagingPoints);
576 for (int i = 0; i <= averagingPoints / 2; i++)
577 smthSignal->IncreaseAmplitude(GetTime(i), sum / averagingPoints);
578
579 for (int i = averagingPoints / 2 + 1; i < GetNumberOfPoints() - averagingPoints / 2; i++) {
580 sum -= this->GetData(i - (averagingPoints / 2 + 1));
581 sum += this->GetData(i + averagingPoints / 2);
582 smthSignal->IncreaseAmplitude(this->GetTime(i), sum / averagingPoints);
583 }
584
585 for (int i = GetNumberOfPoints() - averagingPoints / 2; i < GetNumberOfPoints(); i++)
586 smthSignal->IncreaseAmplitude(GetTime(i), sum / averagingPoints);
587}
588
589Double_t TRestDetectorSignal::GetBaseLine(Int_t startBin, Int_t endBin) {
590 if (endBin - startBin <= 0) return 0.;
591
592 Double_t baseLine = 0;
593 for (int i = startBin; i < endBin; i++) baseLine += fSignalCharge[i];
594
595 return baseLine / (endBin - startBin);
596}
597
598Double_t TRestDetectorSignal::GetStandardDeviation(Int_t startBin, Int_t endBin) {
599 Double_t bL = GetBaseLine(startBin, endBin);
600 return GetBaseLineSigma(startBin, endBin, bL);
601}
602
603Double_t TRestDetectorSignal::GetBaseLineSigma(Int_t startBin, Int_t endBin, Double_t baseline) {
604 Double_t bL = baseline;
605 if (bL == 0) bL = GetBaseLine(startBin, endBin);
606
607 Double_t baseLineSigma = 0;
608 for (int i = startBin; i < endBin; i++)
609 baseLineSigma += (bL - fSignalCharge[i]) * (bL - fSignalCharge[i]);
610
611 return TMath::Sqrt(baseLineSigma / (endBin - startBin));
612}
613
614Double_t TRestDetectorSignal::SubstractBaseline(Int_t startBin, Int_t endBin) {
615 Double_t bL = GetBaseLine(startBin, endBin);
616
617 AddOffset(-bL);
618
619 return bL;
620}
621
622void TRestDetectorSignal::AddOffset(Double_t offset) {
623 for (int i = 0; i < GetNumberOfPoints(); i++) fSignalCharge[i] = fSignalCharge[i] + offset;
624}
625
626void TRestDetectorSignal::MultiplySignalBy(Double_t factor) {
627 for (int i = 0; i < GetNumberOfPoints(); i++) fSignalCharge[i] = factor * fSignalCharge[i];
628}
629
630void TRestDetectorSignal::ExponentialConvolution(Double_t fromTime, Double_t decayTime, Double_t offset) {
631 for (int i = 0; i < GetNumberOfPoints(); i++) {
632 if (fSignalTime[i] > fromTime)
633 fSignalCharge[i] =
634 (fSignalCharge[i] - offset) * exp(-(fSignalTime[i] - fromTime) / decayTime) + offset;
635 }
636}
637
638void TRestDetectorSignal::SignalAddition(TRestDetectorSignal* inSignal) {
639 if (this->GetNumberOfPoints() != inSignal->GetNumberOfPoints()) {
640 cout << "ERROR : I cannot add two signals with different number of points" << endl;
641 return;
642 }
643
644 Int_t badSignalTimes = 0;
645
646 for (int i = 0; i < GetNumberOfPoints(); i++)
647 if (GetTime(i) != inSignal->GetTime(i)) {
648 cout << "Time : " << GetTime(i) << " != " << inSignal->GetTime(i) << endl;
649 badSignalTimes++;
650 }
651
652 if (badSignalTimes) {
653 cout << "ERROR : The times of signal addition must be the same" << endl;
654 return;
655 }
656
657 for (int i = 0; i < GetNumberOfPoints(); i++) fSignalCharge[i] += inSignal->GetData(i);
658}
659
660void TRestDetectorSignal::AddGaussianSignal(Double_t amp, Double_t sigma, Double_t time, Int_t N,
661 Double_t fromTime, Double_t toTime) {
662 for (int i = 0; i < N; i++) {
663 Double_t tme = fromTime + (double)i / (N - 1) * (toTime - fromTime);
664
665 Double_t dta = 300 + amp * TMath::Exp(-0.5 * (tme - time) * (tme - time) / sigma / sigma);
666
667 cout << "T : " << tme << " D : " << dta << endl;
668 IncreaseAmplitude(tme, dta);
669 }
670}
671
672void TRestDetectorSignal::GetWhiteNoiseSignal(TRestDetectorSignal* noiseSignal, Double_t noiseLevel) {
673 this->Sort();
674
675 for (int i = 0; i < GetNumberOfPoints(); i++) {
676 TRandom3* fRandom = new TRandom3(0);
677
678 noiseSignal->IncreaseAmplitude(GetTime(i), GetData(i) + fRandom->Gaus(0, noiseLevel));
679
680 delete fRandom;
681 }
682}
683
684void TRestDetectorSignal::GetSignalGaussianConvolution(TRestDetectorSignal* convSignal, Double_t sigma,
685 Int_t nSigmas) {
686 this->Sort();
687
688 Int_t nPoints = GetMaxTime() - GetMinTime();
689 TF1* fGaus = new TF1("fGauss", "exp(-0.5*((x-[1])/[2])**2)", -nPoints, nPoints);
690 sigma = sigma * 1000.; // conversion to nanoseconds
691 fGaus->SetParameter(2, sigma); // the width of the gaussian is set
692
693 Double_t totChargeInitial = 0.;
694 Double_t totChargeFinal = 0.;
695
696 Double_t sum;
697
698 // We calculate the charge of the event before convolution
699 for (int i = 0; i < GetNumberOfPoints(); i++) totChargeInitial += fSignalCharge[i];
700
701 // The gaussian convolution of the initial signal is performed
702 for (int i = GetMinTime() - nSigmas * sigma; i < GetMaxTime() + nSigmas * sigma; i++) {
703 for (int j = 0; j < GetNumberOfPoints(); j++) {
704 if (TMath::Abs(i - GetTime(j)) > nSigmas * sigma) continue;
705 if (TMath::Abs(i - GetTime(j)) > nSigmas * sigma && i < GetTime(j)) break;
706
707 fGaus->SetParameter(1, GetTime(j));
708 sum = fSignalCharge[j] / TMath::Sqrt(2. * TMath::Pi()) / sigma * fGaus->Integral(i, i + 1);
709
710 convSignal->IncreaseAmplitude(i, sum);
711 totChargeFinal += sum;
712 }
713 }
714
715 cout << "Initial charge of the pulse " << totChargeInitial << endl;
716 cout << "Final charge of the pulse " << totChargeFinal << endl;
717}
718
719void TRestDetectorSignal::WriteSignalToTextFile(const TString& filename) {
720 FILE* fff = fopen(filename.Data(), "w");
721 for (int i = 0; i < GetNumberOfPoints(); i++) {
722 fprintf(fff, "%e\t%e\n", GetTime(i), GetData(i));
723 }
724 fclose(fff);
725}
726
727void TRestDetectorSignal::Print() const {
728 cout << "Signal ID : " << GetSignalID() << endl;
729 cout << "Integral : " << GetIntegral() << endl;
730 if (!GetSignalName().empty()) {
731 cout << "Name: " << GetSignalName() << endl;
732 }
733 if (!GetSignalType().empty()) {
734 cout << "Type: " << GetSignalType() << endl;
735 }
736
737 cout << "------------------------------------------------" << endl;
738 for (int i = 0; i < GetNumberOfPoints(); i++) {
739 cout << "Time : " << GetTime(i) << " Charge : " << GetData(i) << endl;
740 }
741
742 cout << "================================================" << endl;
743}
744
745TGraph* TRestDetectorSignal::GetGraph(Int_t color) {
746 if (fGraph != nullptr) {
747 delete fGraph;
748 fGraph = nullptr;
749 }
750
751 fGraph = new TGraph();
752
753 // cout << "Signal ID " << this->GetSignalID( ) << " points " <<
754 // this->GetNumberOfPoints() << endl;
755
756 fGraph->SetLineWidth(2);
757 fGraph->SetLineColor(color);
758 fGraph->SetMarkerStyle(7);
759
760 int points = 0;
761 for (int n = 0; n < GetNumberOfPoints(); n++) {
762 fGraph->SetPoint(points, GetTime(n), GetData(n));
763 points++;
764 }
765
766 return fGraph;
767}
void IncreaseAmplitude(TVector2 p)
If the point already exists inside the detector signal event, the amplitude value will be added to th...
void SetPoint(TVector2 p)
If the point already exists inside the detector signal event, it will be overwritten....