REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRawSignal.h
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
23#ifndef RestCore_TRestRawSignal
24#define RestCore_TRestRawSignal
25
26#include <TGraph.h>
27#include <TRandom.h>
28#include <TString.h>
29#include <TVector2.h>
30
31#include <iostream>
32#include <string>
33#include <tuple>
34#include <vector>
35
38 private:
40
41 void CalculateBaseLineSigmaSD(Int_t startBin, Int_t endBin);
42
43 void CalculateBaseLineSigmaIQR(Int_t startBin, Int_t endBin);
44
45 void CalculateBaseLineSigmaExcludeOutliers(Int_t startBin, Int_t endBin);
46
47 std::vector<Float_t> GetSignalSmoothed_ExcludeOutliers(Int_t averagingPoints);
48
49 protected:
51 Int_t fSignalID;
52
54 std::vector<Short_t> fSignalData;
55
56 Bool_t fShowWarnings = true;
57
59 UInt_t fSeed = gRandom->GetSeed();
60
61 public:
63 TGraph* fGraph;
64
66 std::vector<Int_t> fPointsOverThreshold;
67
69 Double_t fThresholdIntegral = -1;
70
73
76
78 Double_t fBaseLine = 0;
79
81 Double_t fBaseLineSigma = 0;
82
84 TVector2 fRange = TVector2(0, 0);
85
87 inline Int_t GetSignalID() const { return fSignalID; }
88
90 inline Int_t GetID() const { return fSignalID; }
91
93 inline Int_t GetNumberOfPoints() const { return fSignalData.size(); }
94
96 inline std::vector<Int_t> GetPointsOverThreshold() const { return fPointsOverThreshold; }
97
99 inline Double_t GetMaxValue() { return GetMaxPeakValue(); }
100
102 inline Double_t GetMinValue() { return GetMinPeakValue(); }
103
105 inline Int_t GetHeadPoints() const { return fHeadPoints; }
106
108 inline Int_t GetTailPoints() const { return fTailPoints; }
109
112 inline Double_t GetBaseLine() const { return fBaseLine; }
113
116 inline Double_t GetBaseLineSigma() const { return fBaseLineSigma; }
117
119 inline TVector2 GetRange() const { return fRange; }
120
122 inline Bool_t isBaseLineInitialized() const { return !(fBaseLineSigma == 0 && fBaseLine == 0); }
123
124 Double_t GetData(Int_t n) const;
125
126 Double_t GetRawData(Int_t n) const;
127
128 Short_t operator[](Int_t n);
129
131 inline void SetSignalID(Int_t sID) { fSignalID = sID; }
132
134 inline void SetID(Int_t sID) { fSignalID = sID; }
135
137 inline void SetHeadPoints(Int_t p) { fHeadPoints = p; }
138
140 inline void SetTailPoints(Int_t p) { fTailPoints = p; }
141
143 inline void SetRange(const TVector2& range) { fRange = range; }
144
145 inline void SetRangeToMax() { fRange = TVector2(0, GetNumberOfPoints()); }
146
148 inline void SetRange(Int_t from, Int_t to) { fRange = TVector2(from, to); }
149
150 void Reset();
151
152 void Initialize();
153
154 void AddPoint(Short_t);
155
156 void AddPoint(Double_t);
157
158 void IncreaseBinBy(Int_t bin, Double_t data);
159
160 void InitializePointsOverThreshold(const TVector2& thrPar, Int_t nPointsOver, Int_t nPointsFlat = 512);
161
162 UInt_t GetSeed() const { return fSeed; }
163
164 Double_t GetIntegral();
165
166 Double_t GetIntegralInRange(Int_t startBin, Int_t endBin);
167
168 Double_t GetThresholdIntegral();
169
170 Double_t GetTripleMaxIntegral();
171
172 Double_t GetSlopeIntegral();
173
174 Double_t GetRiseSlope();
175
176 Int_t GetRiseTime();
177
178 Double_t GetAverageInRange(Int_t startBin, Int_t endBin);
179
180 Int_t GetMaxPeakWidth();
181
182 Double_t GetMaxPeakValue();
183
184 Int_t GetMaxPeakBin();
185
186 Double_t GetMinPeakValue();
187
188 Int_t GetMinPeakBin();
189
190 Bool_t IsADCSaturation(int Nflat = 3);
191
192 void GetDifferentialSignal(TRestRawSignal* diffSignal, Int_t smearPoints);
193
194 void GetSignalSmoothed(TRestRawSignal* smoothedSignal, Int_t averagingPoints);
195
196 std::vector<Float_t> GetSignalSmoothed(Int_t averagingPoints, std::string option = "");
197
198 void GetWhiteNoiseSignal(TRestRawSignal* noiseSignal, Double_t noiseLevel = 1.);
199
200 void CalculateBaseLineMean(Int_t startBin, Int_t endBin);
201
202 void CalculateBaseLineMedian(Int_t startBin, Int_t endBin);
203
204 void CalculateBaseLineMedianExcludeOutliers(Int_t startBin, Int_t endBin);
205
206 void CalculateBaseLine(Int_t startBin, Int_t endBin, const std::string& option = "");
207
208 void GetBaseLineCorrected(TRestRawSignal* smoothedSignal, Int_t averagingPoints);
209
210 void AddOffset(Short_t offset);
211
212 void SignalAddition(const TRestRawSignal& signal);
213
214 void Scale(Double_t value);
215
216 void WriteSignalToTextFile(const TString& filename);
217
218 void Print() const;
219
220 void SetSeed(UInt_t seed) { fSeed = seed; }
221
222 TGraph* GetGraph(Int_t color = 1);
223
227 std::vector<std::tuple<double, UShort_t, double>> GetPeaks(double threshold, UShort_t distance = 5,
228 double signalBaseLine = 0.0) const;
229 std::vector<std::tuple<double, UShort_t, double>> GetPeaksVeto(double threshold, UShort_t distance = 5,
230 double signalBaseLine = 0.0) const;
231
233 TRestRawSignal(Int_t nBins);
235
236 ClassDef(TRestRawSignal, 3);
237};
238#endif
It defines a Short_t array with a physical parameter that evolves in time using a fixed time bin.
Int_t GetMinPeakBin()
It returns the bin at which the minimum peak amplitude happens.
Double_t GetBaseLineSigma() const
void WriteSignalToTextFile(const TString &filename)
This method dumps to a text file the data inside fSignalData.
Double_t GetThresholdIntegral()
It returns the integral of points identified over threshold. InitializePointsOverThreshold should hav...
Int_t GetID() const
Returns the value of signal ID.
void SetID(Int_t sID)
It sets the id number of the signal.
~TRestRawSignal()
Default destructor.
TGraph * fGraph
A TGraph pointer used to store the TRestRawSignal drawing.
Double_t fBaseLine
This baseline value will be subtracted from GetData for any raw signal observable calculation.
void GetSignalSmoothed(TRestRawSignal *smoothedSignal, Int_t averagingPoints)
It smoothes the existing signal and places it at the signal pointer given by argument.
Double_t GetRiseSlope()
It returns the rate of change or slope from the points that have been identified over threshlold on t...
std::vector< Short_t > fSignalData
Vector with the data of the signal.
Double_t GetMaxValue()
Returns the maximum value found in the data points. It includes baseline correction.
void Print() const
It prints the signal data on screen.
Double_t fThresholdIntegral
It stores the integral value obtained from the points identified over threshold.
Int_t fSignalID
An integer value used to attribute a unique identification number to the signal.
Double_t GetMaxPeakValue()
It returns the amplitude of the signal maximum, baseline will be corrected if CalculateBaseLine was c...
void CalculateBaseLineSigmaExcludeOutliers(Int_t startBin, Int_t endBin)
This method is called by CalculateBaseLine with the "OUTLIERS"-option to determine the value of the b...
std::vector< Float_t > GetSignalSmoothed_ExcludeOutliers(Int_t averagingPoints)
It smooths the existing signal and returns it in a vector of Float_t values. This method excludes poi...
Double_t GetMinPeakValue()
It returns the amplitude of the signal minimum, 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 ...
Int_t GetTailPoints() const
Returns the number of tail points used on points over threshold definition.
void CalculateThresholdIntegral()
This method will be called each time InitializePointsOverThreshold is called to re-define the value o...
void Initialize()
Initialization of TRestRawSignal members.
void CalculateBaseLineMean(Int_t startBin, Int_t endBin)
This method is called by CalculateBaseLine and is used to determine the value of the baseline as aver...
TVector2 GetRange() const
Returns the range defined by user.
void CalculateBaseLineSigmaSD(Int_t startBin, Int_t endBin)
This method is called by CalculateBaseLine to determine the value of the baseline fluctuation as its ...
void SetSignalID(Int_t sID)
It sets the id number of the signal.
TRestRawSignal()
Default constructor.
std::vector< Int_t > GetPointsOverThreshold() const
Returns a std::vector containing the indexes of data points over threshold.
Double_t GetBaseLine() const
Short_t operator[](Int_t n)
It overloads the operator [] so that we can retrieve a particular point n in the form rawSignal[n].
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...
Double_t GetIntegralInRange(Int_t startBin, Int_t endBin)
It returns the integral of points found in the specific range given by (startBin,endBin).
void Reset()
Initializes the existing signal data and sets it to zero while keeping the array size.
void SignalAddition(const TRestRawSignal &signal)
This method adds the signal provided by argument to the existing signal.
Int_t fTailPoints
It defines the number of points to include after point over threshold definition. NOT implemented.
void AddOffset(Short_t offset)
This method adds an offset to the signal data.
Int_t GetMaxPeakWidth()
It returns the temporal width of the peak with maximum amplitude inside the signal.
Int_t GetRiseTime()
It returns the time required from the signal to reach the maximum. InitializePointsOverThreshold shou...
void CalculateBaseLineMedianExcludeOutliers(Int_t startBin, Int_t endBin)
This method is called by CalculateBaseLine with the "OUTLIERS"-option and is used to determine the va...
Double_t GetAverageInRange(Int_t startBin, Int_t endBin)
It returns the average of the points found in the range (startBin, endBin)
void SetRange(Int_t from, Int_t to)
It sets/constrains the range for any calculation.
void IncreaseBinBy(Int_t bin, Double_t data)
It adds the content of data to fSignalData[bin].
std::vector< std::tuple< double, UShort_t, double > > GetPeaks(double threshold, UShort_t distance=5, double signalBaseLine=0.0) const
void GetDifferentialSignal(TRestRawSignal *diffSignal, Int_t smearPoints)
It calculates the differential signal of the existing signal and it will place at the signal pointer ...
void CalculateBaseLineMedian(Int_t startBin, Int_t endBin)
This method is called by CalculateBaseLine with the "ROBUST"-option and is used to determine the valu...
std::vector< Int_t > fPointsOverThreshold
A std::vector containing the index of points that are identified over threshold.
Double_t GetRawData(Int_t n) const
It returns the original data value of point n without baseline correction.
void SetTailPoints(Int_t p)
It sets the number of tail points.
Double_t GetData(Int_t n) const
It returns the data value of point n including baseline correction.
Double_t GetTripleMaxIntegral()
It returns the integral calculated using the maximum signal amplitude and its neightbour points.
void Scale(Double_t value)
This method scales the signal by a given value.
Double_t GetIntegral()
It returns the integral of points found in the region defined by fRange. If fRange was not defined th...
void AddPoint(Short_t)
Adds a new point to the end of the signal data array.
TGraph * GetGraph(Int_t color=1)
It builds a TGraph object that can be used for drawing.
Double_t GetMinValue()
Returns the lowest value found in the data points. It includes baseline correction.
Int_t GetHeadPoints() const
Returns the number of head points used on points over threshold definition.
UInt_t fSeed
Seed used for random number generation.
Double_t GetSlopeIntegral()
It returns the integral of points identified over threshold found in the first positive rise of the s...
Double_t fBaseLineSigma
The baseline fluctuation calculated as the standard deviation of the baseline.
void GetWhiteNoiseSignal(TRestRawSignal *noiseSignal, Double_t noiseLevel=1.)
It calculates an arbitrary Gaussian noise placing it at the signal pointer given by argument....
Bool_t IsADCSaturation(int Nflat=3)
It returns whether the signal has ADC saturation.
Int_t fHeadPoints
It defines the number of points to include before point over threshold definition....
Int_t GetMaxPeakBin()
It returns the bin at which the maximum peak amplitude happens.
Int_t GetSignalID() const
Returns the value of signal ID.
void GetBaseLineCorrected(TRestRawSignal *smoothedSignal, Int_t averagingPoints)
It applies the moving average filter (GetSignalSmoothed) to the signal, which is then subtracted from...
TVector2 fRange
Any signal calculation will be restricted to the following range definition.
void SetHeadPoints(Int_t p)
It sets the number of head points.
void SetRange(const TVector2 &range)
It sets/constrains the range for any calculation.
void CalculateBaseLineSigmaIQR(Int_t startBin, Int_t endBin)
This method is called by CalculateBaseLine with the "ROBUST"-option to determine the value of the bas...
Int_t GetNumberOfPoints() const
Returns the actual number of points, or size of the signal.
Bool_t isBaseLineInitialized() const
Returns false if the baseline and its baseline fluctuation was not initialized.