REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRawSignalEvent.cxx
1/*************************************************************************
2 * This file is part of the REST software framework. *
3 * *
4 * Copyright (C) 2016 GIFNA/TREX (University of Zaragoza) *
5 * For more information see http://gifna.unizar.es/trex *
6 * *
7 * REST is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * REST is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have a copy of the GNU General Public License along with *
18 * REST in $REST_PATH/LICENSE. *
19 * If not, see http://www.gnu.org/licenses/. *
20 * For the list of contributors see $REST_PATH/CREDITS. *
21 *************************************************************************/
22
63
64#include "TRestRawSignalEvent.h"
65
66#include <TMath.h>
67
68#include "TRestStringHelper.h"
69
70using namespace std;
71
72ClassImp(TRestRawSignalEvent);
73
74TRestRawSignalEvent::TRestRawSignalEvent() {
75 // TRestRawSignalEvent default constructor
76 Initialize();
77}
78
79TRestRawSignalEvent::~TRestRawSignalEvent() {
80 // TRestRawSignalEvent destructor
81}
82
85 fSignal.clear();
86 fPad = nullptr;
87 gr = nullptr;
88
89 fMinValue = numeric_limits<Double_t>::max();
90 fMaxValue = numeric_limits<Double_t>::min();
91 fMinTime = numeric_limits<Double_t>::max();
92 fMaxTime = numeric_limits<Double_t>::min();
93}
94
95void TRestRawSignalEvent::AddSignal(TRestRawSignal& s) {
96 if (signalIDExists(s.GetSignalID())) {
97 cout << "Warning. Signal ID : " << s.GetSignalID()
98 << " already exists. Signal will not be added to signal event" << endl;
99 return;
100 }
101
102 s.CalculateBaseLine(fBaseLineRange.X(), fBaseLineRange.Y());
103 s.SetRange(fRange);
104
105 fSignal.emplace_back(s);
106}
107
108void TRestRawSignalEvent::RemoveSignalWithId(Int_t sId) {
109 Int_t index = GetSignalIndex(sId);
110
111 if (index == -1) {
112 cout << "Warning. Signal ID : " << sId
113 << " does not exist. Signal will not be removed from signal event" << endl;
114 return;
115 }
116
117 fSignal.erase(fSignal.begin() + index);
118}
119
120Int_t TRestRawSignalEvent::GetSignalIndex(Int_t signalID) {
121 for (int i = 0; i < GetNumberOfSignals(); i++)
122 if (fSignal[i].GetSignalID() == signalID) return i;
123 return -1;
124}
125
126Double_t TRestRawSignalEvent::GetIntegral() {
127 Double_t sum = 0;
128
129 for (int i = 0; i < GetNumberOfSignals(); i++) sum += fSignal[i].GetIntegral();
130
131 return sum;
132}
133
137 Double_t sum = 0;
138 for (int i = 0; i < GetNumberOfSignals(); i++) sum += fSignal[i].GetThresholdIntegral();
139 return sum;
140}
141
142TRestRawSignal* TRestRawSignalEvent::GetMaxSignal() {
143 if (GetNumberOfSignals() <= 0) return nullptr;
144
145 Double_t max = fSignal[0].GetIntegral();
146
147 Int_t sId = 0;
148 for (int i = 0; i < GetNumberOfSignals(); i++) {
149 Int_t integ = fSignal[i].GetIntegral();
150 if (max < integ) {
151 max = integ;
152 sId = i;
153 }
154 }
155
156 return &fSignal[sId];
157}
158
159Double_t TRestRawSignalEvent::GetSlopeIntegral() {
160 Double_t sum = 0;
161
162 for (int i = 0; i < GetNumberOfSignals(); i++) sum += fSignal[i].GetSlopeIntegral();
163
164 return sum;
165}
166
167Double_t TRestRawSignalEvent::GetRiseSlope() {
168 Double_t sum = 0;
169
170 Int_t n = 0;
171 for (int i = 0; i < GetNumberOfSignals(); i++) {
172 if (fSignal[i].GetThresholdIntegral() > 0) {
173 sum += fSignal[i].GetSlopeIntegral();
174 n++;
175 }
176 }
177
178 if (n == 0) return 0;
179
180 return sum / n;
181}
182
183Double_t TRestRawSignalEvent::GetRiseTime() {
184 Double_t sum = 0;
185
186 Int_t n = 0;
187 for (int i = 0; i < GetNumberOfSignals(); i++) {
188 if (fSignal[i].GetThresholdIntegral() > 0) {
189 sum += fSignal[i].GetRiseTime();
190 n++;
191 }
192 }
193
194 if (n == 0) return 0;
195
196 return sum / n;
197}
198
199Double_t TRestRawSignalEvent::GetTripleMaxIntegral() {
200 Double_t sum = 0;
201
202 for (int i = 0; i < GetNumberOfSignals(); i++)
203 if (fSignal[i].GetThresholdIntegral() > 0) sum += fSignal[i].GetTripleMaxIntegral();
204
205 return sum;
206}
207
208Double_t TRestRawSignalEvent::GetBaseLineAverage() {
209 Double_t baseLineMean = 0;
210
211 for (int signal = 0; signal < GetNumberOfSignals(); signal++) {
212 Double_t baseline = GetSignal(signal)->GetBaseLine();
213 baseLineMean += baseline;
214 }
215
216 return baseLineMean / GetNumberOfSignals();
217}
218
219Int_t TRestRawSignalEvent::GetLowestWidth(Double_t minPeakAmplitude) {
220 Int_t low = 10000000;
221
222 for (int signal = 0; signal < GetNumberOfSignals(); signal++) {
223 if (GetSignal(signal)->GetMaxPeakValue() > minPeakAmplitude) {
224 Int_t lW = GetSignal(signal)->GetMaxPeakWidth();
225 if (low > lW) low = lW;
226 }
227 }
228
229 return low;
230}
231
232Double_t TRestRawSignalEvent::GetAverageWidth(Double_t minPeakAmplitude) {
233 Double_t avg = 0;
234 Int_t n = 0;
235 for (int signal = 0; signal < GetNumberOfSignals(); signal++) {
236 if (GetSignal(signal)->GetMaxPeakValue() > minPeakAmplitude) {
237 avg += GetSignal(signal)->GetMaxPeakWidth();
238 n++;
239 }
240 }
241
242 if (n == 0)
243 return 0;
244 else
245 return avg / n;
246}
247
248Double_t TRestRawSignalEvent::GetLowAverageWidth(Int_t nSignals, Double_t minPeakAmplitude) {
249 std::vector<Double_t> widths;
250
251 for (int signal = 0; signal < GetNumberOfSignals(); signal++)
252 if (GetSignal(signal)->GetMaxPeakValue() > minPeakAmplitude)
253 widths.push_back(GetSignal(signal)->GetMaxPeakWidth());
254
255 if (widths.size() == 0) return 0;
256
257 std::sort(widths.begin(), widths.end());
258
259 Int_t nMax = nSignals;
260 if (widths.size() < (unsigned int)nSignals) nMax = widths.size();
261
262 Double_t avg = 0;
263 for (int n = 0; n < nMax; n++) avg += widths[n];
264
265 return avg / nSignals;
266}
267
268Double_t TRestRawSignalEvent::GetBaseLineSigmaAverage() {
269 Double_t baseLineSigmaMean = 0;
270
271 for (int signal = 0; signal < GetNumberOfSignals(); signal++) {
272 Double_t baselineSigma = GetSignal(signal)->GetBaseLineSigma();
273 baseLineSigmaMean += baselineSigma;
274 }
275
276 return baseLineSigmaMean / GetNumberOfSignals();
277}
278
284// void TRestRawSignalEvent::subtractBaselines() {
285// for (int signal = 0; signal < GetNumberOfSignals(); signal++)
286// GetSignal(signal)->subtractBaseline();
287//}
288
289void TRestRawSignalEvent::AddChargeToSignal(Int_t signalID, Int_t bin, Short_t value) {
290 Int_t signalIndex = GetSignalIndex(signalID);
291 if (signalIndex == -1) {
292 signalIndex = GetNumberOfSignals();
293
294 TRestRawSignal signal(512); // For the moment we use the default nBins=512
295 signal.SetSignalID(signalID);
296 AddSignal(signal);
297 }
298
299 fSignal[signalIndex].IncreaseBinBy(bin, value);
300}
301
302void TRestRawSignalEvent::PrintEvent() {
304
305 for (int i = 0; i < GetNumberOfSignals(); i++) {
306 cout << "================================================" << endl;
307 cout << "Signal ID : " << fSignal[i].GetSignalID() << endl;
308 cout << "Integral : " << fSignal[i].GetIntegral() << endl;
309 cout << "------------------------------------------------" << endl;
310 fSignal[i].Print();
311 cout << "================================================" << endl;
312 }
313}
314
315// TODO: GetMaxTimeFast, GetMinTimeFast, GetMaxValueFast that return the value
316// of fMinTime, fMaxTime, etc
317void TRestRawSignalEvent::SetMaxAndMin() {
318 fMinValue = numeric_limits<Double_t>::max();
319 fMaxValue = numeric_limits<Double_t>::min();
320 fMinTime = numeric_limits<Double_t>::max();
321 fMaxTime = numeric_limits<Double_t>::min();
322
323 for (int s = 0; s < GetNumberOfSignals(); s++) {
324 if (fMinValue > fSignal[s].GetMinValue()) {
325 fMinValue = fSignal[s].GetMinValue();
326 }
327 if (fMaxValue < fSignal[s].GetMaxValue()) {
328 fMaxValue = fSignal[s].GetMaxValue();
329 }
330 }
331
332 if (GetNumberOfSignals() > 0) {
333 fMaxTime = fSignal[0].GetNumberOfPoints();
334 }
335}
336
337Double_t TRestRawSignalEvent::GetMaxValue() {
338 SetMaxAndMin();
339 return fMaxValue;
340}
341
342Double_t TRestRawSignalEvent::GetMinValue() {
343 SetMaxAndMin();
344 return fMinValue;
345}
346
347Double_t TRestRawSignalEvent::GetMinTime() { return 0; }
348
349Double_t TRestRawSignalEvent::GetMaxTime() {
350 Double_t maxTime = 512;
351
352 if (GetNumberOfSignals() > 0) maxTime = fSignal[0].GetNumberOfPoints();
353
354 return maxTime;
355}
356
357// An auxiliar variable that can be used to store temporal values.
358Double_t fAuxiliar = 0;
359void TRestRawSignalEvent::SetAuxiliar(Double_t aux) { fAuxiliar = aux; }
360auto GetAuxiliar() { return fAuxiliar; }
361
413TPad* TRestRawSignalEvent::DrawEvent(const TString& option) {
414 const int nSignals = GetNumberOfSignals();
415
416 if (fPad != nullptr) {
417 delete fPad;
418 fPad = nullptr;
419 }
420
421 if (nSignals == 0) {
422 cout << "Empty event " << endl;
423 return nullptr;
424 }
425
426 fMinValue = numeric_limits<Double_t>::max();
427 fMaxValue = numeric_limits<Double_t>::min();
428 fMinTime = numeric_limits<Double_t>::max();
429 fMaxTime = numeric_limits<Double_t>::min();
430
431 fPad = new TPad(GetName(), " ", 0, 0, 1, 1);
432 fPad->Draw();
433 fPad->cd();
434 // fPad->DrawFrame(0, GetMinValue() - 1, GetMaxTime() + 1, GetMaxValue() + 1);
435
436 char title[256];
437 vector<TString> optList = Vector_cast<string, TString>(TRestTools::GetOptions((string)option));
438
439 bool thresholdCheck = false;
440 bool baselineCheck = false;
441 bool sRangeID = false;
442 bool printIDs = false;
443
444 double signalTh = 0, pointTh = 0, nOver = 0;
445 int baseLineRangeInit = 0, baseLineRangeEnd = 0;
446 int sRangeInit = 0, sRangeEnd = 0;
447
448 for (const auto& opt : optList) {
449 std::string str = (std::string)opt;
450 // Read threshold option
451 if (str.find("onlyGoodSignals[") != string::npos) {
452 size_t startPos = str.find('[');
453 size_t endPos = str.find(']');
454 TString tmpStr = opt(startPos + 1, endPos - startPos - 1);
455 vector<TString> optList_2 = Vector_cast<string, TString>(Split((string)tmpStr, ","));
456
457 pointTh = StringToDouble((string)optList_2[0]);
458 signalTh = StringToDouble((string)optList_2[1]);
459 nOver = StringToDouble((string)optList_2[2]);
460
461 thresholdCheck = true;
462 }
463
464 // Read baseline option
465 if (str.find("baseLineRange[") != string::npos) {
466 size_t startPos2 = str.find('[');
467 size_t endPos2 = str.find(']');
468 TString tmpStr2 = opt(startPos2 + 1, endPos2 - startPos2 - 1);
469 vector<TString> optList_3 = Vector_cast<string, TString>(Split((string)tmpStr2, ","));
470
471 baseLineRangeInit = StringToInteger((string)optList_3[0]);
472 baseLineRangeEnd = StringToInteger((string)optList_3[1]);
473
474 baselineCheck = true;
475 }
476
477 // Read signal range ID option
478 if (str.find("signalRangeID[") != string::npos || str.find("ids[") != string::npos) {
479 size_t startPos3 = str.find('[');
480 size_t endPos3 = str.find(']');
481 TString tmpStr3 = opt(startPos3 + 1, endPos3 - startPos3 - 1);
482 vector<TString> optList_4;
483 if (str.find(',') != string::npos)
484 optList_4 = Vector_cast<string, TString>(Split((string)tmpStr3, ","));
485 else if (str.find('-') != string::npos)
486 optList_4 = Vector_cast<string, TString>(Split((string)tmpStr3, "-"));
487 else
488 RESTError << "TRestRawSignalEvent::DrawEvent not valid ids format!" << RESTendl;
489
490 sRangeInit = StringToInteger((string)optList_4[0]);
491 sRangeEnd = StringToInteger((string)optList_4[1]);
492
493 sRangeID = true;
494 }
495
496 // Read print ID option
497 if (str.find("printIDs") != string::npos) {
498 printIDs = true;
499 cout << "IDs of printed signals: " << endl;
500 }
501 }
502
503 std::vector<int> signalIDs; // Signal IDs to print
504
506 if ((optList.empty()) || !(isANumber((string)optList[0]))) {
507 // If threshold and baseline options are given
508 if (thresholdCheck && baselineCheck) {
509 RESTDebug << "Draw only good signals with: " << RESTendl;
510 RESTDebug << " Signal threshold: " << signalTh << RESTendl;
511 RESTDebug << " Point threshold: " << pointTh << RESTendl;
512 RESTDebug << " Points over threshold: " << nOver << RESTendl;
513 RESTDebug << " Base line range: (" << baseLineRangeInit << "," << baseLineRangeEnd << ")"
514 << RESTendl;
515
516 for (int n = 0; n < nSignals; n++) {
517 fSignal[n].CalculateBaseLine(baseLineRangeInit, baseLineRangeEnd);
518 fSignal[n].InitializePointsOverThreshold(TVector2(pointTh, signalTh), nOver);
519 if (fSignal[n].GetPointsOverThreshold().size() >= 2) {
520 signalIDs.push_back(fSignal[n].GetID());
521 }
522 }
523 // If no threshold and baseline options are given
524 } else {
525 for (int n = 0; n < nSignals; n++) {
526 signalIDs.push_back(fSignal[n].GetID());
527 }
528 }
529
530 // Remove SIDs which are not in range
531 if (sRangeID) {
532 for (auto it = signalIDs.begin(); it != signalIDs.end();)
533 if (*it >= sRangeInit && *it <= sRangeEnd) {
534 ++it;
535 } else {
536 it = signalIDs.erase(it);
537 }
538 }
539
540 cout << "Number of drawn signals: " << signalIDs.size() << endl;
541
543 } else if (isANumber((string)optList[0])) {
544 string str = (string)optList[0];
545 size_t separation = str.find('-');
546
547 // Signals range //
548 if (separation != string::npos) {
549 TString firstSignal = optList[0](0, separation);
550 TString lastSignal = optList[0](separation + 1, str.size());
551 RESTDebug << "First signal: " << firstSignal << RESTendl;
552 RESTDebug << "Last signal: " << lastSignal << RESTendl;
553
554 if (StringToInteger((string)lastSignal) >= (int)fSignal.size()) {
555 fPad->SetTitle("No Such Signal");
556 cout << "No such last signal" << endl;
557 return fPad;
558 }
559
560 sprintf(title, "Event ID %d", this->GetID());
561
562 if (thresholdCheck && baselineCheck) {
563 RESTDebug << "Draw only good signals with: " << RESTendl;
564 RESTDebug << " Signal threshold: " << signalTh << RESTendl;
565 RESTDebug << " Point threshold: " << pointTh << RESTendl;
566 RESTDebug << " Points over threshold: " << nOver << RESTendl;
567 RESTDebug << " Base line range: (" << baseLineRangeInit << "," << baseLineRangeEnd << ")"
568 << RESTendl;
569
570 for (int n = 0; n < nSignals; n++) {
571 if (n < StringToInteger((string)firstSignal) || n > StringToInteger((string)lastSignal))
572 continue;
573 fSignal[n].CalculateBaseLine(baseLineRangeInit, baseLineRangeEnd);
574 fSignal[n].InitializePointsOverThreshold(TVector2(pointTh, signalTh), nOver);
575 if (fSignal[n].GetPointsOverThreshold().size() >= 2) {
576 signalIDs.push_back(fSignal[n].GetID());
577 }
578 }
579 cout << "Number of good signals in range (" << firstSignal << "," << lastSignal
580 << "): " << signalIDs.size() << endl;
581 // If no threshold and baseline options are given
582 } else {
583 for (int n = StringToInteger((string)firstSignal);
584 n < StringToInteger((string)lastSignal) + 1; n++) {
585 signalIDs.push_back(fSignal[n].GetID());
586 }
587 }
588 // Single signal
589 } else {
590 int signalID = StringToInteger((string)optList[0]);
591 signalIDs.push_back(signalID);
592 }
593 }
594
595 if (signalIDs.empty()) {
596 fPad->SetTitle("No Such Signal");
597 cout << "No signals found" << endl;
598 return fPad;
599 }
600
601 if (printIDs) {
602 cout << "SignalIDs:";
603 auto sortedSignalsIDs = signalIDs;
604 sort(sortedSignalsIDs.begin(), sortedSignalsIDs.end());
605 for (const auto& signalID : sortedSignalsIDs) {
606 cout << " " << signalID;
607 }
608 cout << endl;
609 }
610
611 DrawSignals(fPad, signalIDs);
612
613 return fPad;
614}
615
620void TRestRawSignalEvent::DrawSignals(TPad* pad, const std::vector<Int_t>& signals) {
621 int maxSID = -1;
622 int max = numeric_limits<Short_t>::min();
623 int graphIndex = 1;
624
625 for (const auto& s : signals) {
626 TRestRawSignal* signal = GetSignalById(s);
627 if (!signal) {
628 continue;
629 }
630 TGraph* graph = signal->GetGraph(graphIndex++);
631 const double maxValue = TMath::MaxElement(graph->GetN(), graph->GetY());
632 if (maxValue > max) {
633 max = maxValue;
634 maxSID = s;
635 }
636 }
637
638 RESTDebug << "Max SID " << maxSID << RESTendl;
639
640 if (maxSID == -1) {
641 cout << "No signals ID found" << endl;
642 return;
643 }
644
645 TRestRawSignal* signalMaxID = GetSignalById(maxSID);
646 std::string title = "Event ID " + std::to_string(GetID());
647 if (signals.size() == 1) {
648 title += " Signal ID " + std::to_string(maxSID);
649 }
650
651 signalMaxID->fGraph->SetTitle(title.c_str());
652 signalMaxID->fGraph->GetXaxis()->SetTitle("Time bin");
653 signalMaxID->fGraph->GetYaxis()->SetTitleOffset(1.4);
654 signalMaxID->fGraph->GetYaxis()->SetTitle("Amplitude [a.u.]");
655 pad->Draw();
656 pad->cd();
657 signalMaxID->fGraph->SetLineColor(kBlack);
658 signalMaxID->fGraph->Draw("AL");
659
660 for (const auto& signalID : signals) {
661 if (signalID == maxSID) {
662 continue;
663 }
664 TRestRawSignal* signal = GetSignalById(signalID);
665 pad->cd();
666 signal->fGraph->Draw("L");
667 }
668
669 pad->Update();
670}
671
692TPad* TRestRawSignalEvent::DrawSignal(Int_t signalID, const TString& option) {
693 int nSignals = this->GetNumberOfSignals();
694
695 if (fPad != nullptr) {
696 for (int n = 0; n < nSignals; n++) {
697 delete fSignal[n].fGraph;
698 fSignal[n].fGraph = nullptr;
699 }
700 delete fPad;
701 fPad = nullptr;
702 }
703
704 if (nSignals == 0) {
705 cout << "Empty event " << endl;
706 return nullptr;
707 }
708
709 vector<TString> optList = Vector_cast<string, TString>(TRestTools::GetOptions((string)option));
710
711 double signalTh = 0, pointTh = 0, nOver = 0;
712 int baseLineRangeInit = 0, baseLineRangeEnd = 0;
713
714 for (auto& opt : optList) {
715 string str = (string)opt;
716
717 // Read threshold option
718 size_t goodSigOpt = str.find("goodSignals[");
719
720 if (goodSigOpt != string::npos) {
721 size_t startPos = str.find('[');
722 size_t endPos = str.find(']');
723 TString tmpStr = opt(startPos + 1, endPos - startPos - 1);
724 vector<TString> optList_2 = Vector_cast<string, TString>(Split((string)tmpStr, ","));
725
726 pointTh = StringToDouble((string)optList_2[0]);
727 signalTh = StringToDouble((string)optList_2[1]);
728 nOver = StringToDouble((string)optList_2[2]);
729 }
730
731 // Read baseline option
732 size_t BLOpt = str.find("baseLineRange[");
733
734 if (BLOpt != string::npos) {
735 size_t startPos2 = str.find('[');
736 size_t endPos2 = str.find(']');
737 TString tmpStr2 = opt(startPos2 + 1, endPos2 - startPos2 - 1);
738 vector<TString> optList_3 = Vector_cast<string, TString>(Split((string)tmpStr2, ","));
739
740 baseLineRangeInit = StringToInteger((string)optList_3[0]);
741 baseLineRangeEnd = StringToInteger((string)optList_3[1]);
742 }
743 }
744
745 fPad = new TPad(this->GetName(), " ", 0, 0, 1, 1);
746 fPad->Draw();
747 fPad->cd();
748
749 TGraph* gr = new TGraph();
750
751 TRestRawSignal* signal = this->GetSignalById(signalID);
752 signal->CalculateBaseLine(baseLineRangeInit, baseLineRangeEnd);
753 signal->InitializePointsOverThreshold(TVector2(pointTh, signalTh), nOver);
754
755 RESTInfo << "Drawing signalID. Event ID : " << this->GetID() << " Signal ID : " << signal->GetID()
756 << RESTendl;
757
758 for (int n = 0; n < signal->GetNumberOfPoints(); n++) {
759 gr->SetPoint(n, n, signal->GetData(n));
760 }
761
762 gr->Draw("AC*");
763
764 TGraph* gr2 = new TGraph();
765
766 gr2->SetLineWidth(2);
767 gr2->SetLineColor(2); // Red
768
769 for (int n = baseLineRangeInit; n < baseLineRangeEnd; n++) {
770 gr2->SetPoint(n - baseLineRangeInit, n, signal->GetData(n));
771 }
772
773 gr2->Draw("CP");
774
775 vector<Int_t> pOver = signal->GetPointsOverThreshold();
776
777 TGraph* gr3[5];
778 Int_t nGraphs = 0;
779 gr3[nGraphs] = new TGraph();
780 gr3[nGraphs]->SetLineWidth(2);
781 gr3[nGraphs]->SetLineColor(3);
782 Int_t point = 0;
783 Int_t nPoints = pOver.size();
784 for (int n = 0; n < nPoints; n++) {
785 gr3[nGraphs]->SetPoint(point, pOver[n], signal->GetData(pOver[n]));
786 point++;
787 if (n + 1 < nPoints && pOver[n + 1] - pOver[n] > 1) {
788 gr3[nGraphs]->Draw("CP");
789 nGraphs++;
790 if (nGraphs > 4) cout << "Ngraphs : " << nGraphs << endl;
791 point = 0;
792 gr3[nGraphs] = new TGraph();
793 gr3[nGraphs]->SetLineWidth(2);
794 gr3[nGraphs]->SetLineColor(3); // Green
795 }
796 }
797
798 if (nPoints > 0) {
799 gr3[nGraphs]->Draw("CP");
800 }
801
802 return fPad;
803}
804
805TRestRawReadoutMetadata* TRestRawSignalEvent::GetReadoutMetadata() const {
806 if (fRun == nullptr) {
807 RESTError << "TRestRawSignalEvent::GetReadoutMetadata: fRun is nullptr" << RESTendl;
808 return nullptr;
809 }
810 return dynamic_cast<TRestRawReadoutMetadata*>(fRun->GetMetadataClass("TRestRawReadoutMetadata"));
811}
812
813TRestRawSignalEvent TRestRawSignalEvent::GetSignalEventForType(const string& type) const {
814 return GetSignalEventForTypes({type});
815}
816
817TRestRawSignalEvent TRestRawSignalEvent::GetSignalEventForTypes(
818 const std::set<std::string>& types, const TRestRawReadoutMetadata* readoutMetadata) const {
819 // TODO: verify this works
820 TRestRawSignalEvent signalEvent;
821 signalEvent.SetEventInfo((TRestEvent*)this);
822 auto metadata = readoutMetadata ? readoutMetadata : GetReadoutMetadata();
823 if (metadata == nullptr) {
824 // cerr << "TRestRawSignalEvent::GetSignalEventForTypes: metadata is nullptr" << endl;
825 // exit(1);
826 }
827 for (const auto& signal : fSignal) {
828 if (types.empty() ||
829 types.find(metadata->GetTypeForChannelDaqId(signal.GetSignalID())) != types.end()) {
830 signalEvent.AddSignal(const_cast<TRestRawSignal&>(signal));
831 }
832 }
833 return signalEvent;
834}
A base class for any REST event.
Definition: TRestEvent.h:38
void SetEventInfo(TRestEvent *eve)
Definition: TRestEvent.cxx:137
virtual void PrintEvent() const
Definition: TRestEvent.cxx:187
virtual void Initialize()=0
Definition: TRestEvent.cxx:73
An event container for time rawdata signals with fixed length.
void DrawSignals(TPad *pad, const std::vector< Int_t > &signals)
This method draws selected signal IDs, given by the vector passed as reference.
TPad * DrawEvent(const TString &option="")
This method draws current raw signal event in a TPad.
TPad * DrawSignal(Int_t signalID, const TString &option="")
This method draws selected signalID by ID, with baseline range and points over threshold highlighted.
void AddChargeToSignal(Int_t sgnlID, Int_t bin, Short_t value)
It defines a Short_t array with a physical parameter that evolves in time using a fixed time bin.
Double_t GetBaseLineSigma() const
Int_t GetID() const
Returns the value of signal ID.
TGraph * fGraph
A TGraph pointer used to store the TRestRawSignal drawing.
Double_t GetMaxPeakValue()
It returns the amplitude of the signal maximum, baseline will be corrected if CalculateBaseLine was c...
void CalculateBaseLine(Int_t startBin, Int_t endBin, const std::string &option="")
This method calculates the average and fluctuation of the baseline in the specified range and writes ...
void SetSignalID(Int_t sID)
It sets the id number of the signal.
std::vector< Int_t > GetPointsOverThreshold() const
Returns a std::vector containing the indexes of data points over threshold.
Double_t GetBaseLine() const
void InitializePointsOverThreshold(const TVector2 &thrPar, Int_t nPointsOver, Int_t nPointsFlat=512)
It initializes the fPointsOverThreshold array with the indexes of data points that are found over thr...
Int_t GetMaxPeakWidth()
It returns the temporal width of the peak with maximum amplitude inside the signal.
Double_t GetData(Int_t n) const
It returns the data value of point n including baseline correction.
TGraph * GetGraph(Int_t color=1)
It builds a TGraph object that can be used for drawing.
Int_t GetSignalID() const
Returns the value of signal ID.
void SetRange(const TVector2 &range)
It sets/constrains the range for any calculation.
Int_t GetNumberOfPoints() const
Returns the actual number of points, or size of the signal.
static std::vector< std::string > GetOptions(std::string optionsStr)
Returns all the options in an option string.
Definition: TRestTools.cxx:86
std::vector< std::string > Split(std::string in, std::string separator, bool allowBlankString=false, bool removeWhiteSpaces=false, int startPos=-1)
Split the input string according to the given separator. Returning a vector of fragments.
Double_t StringToDouble(std::string in)
Gets a double from a string.
Int_t StringToInteger(std::string in)
Gets an integer from a string.
Int_t isANumber(std::string in)
Returns 1 only if a valid number is found in the string in. If not it returns 0.