REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRawToDetectorSignalProcess.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 
121 
122 #include "TRestRawToDetectorSignalProcess.h"
123 
124 using namespace std;
125 
127 
132 
137 
143  SetSectionName(this->ClassName());
144  SetLibraryVersion(LIBRARY_VERSION);
145 
146  fInputSignalEvent = nullptr;
147  fOutputSignalEvent = new TRestDetectorSignalEvent();
148 }
149 
154  fInputSignalEvent = (TRestRawSignalEvent*)inputEvent;
155 
156  Int_t rejectedSignal = 0;
157 
158  if (fZeroSuppression) {
159  fInputSignalEvent->SetBaseLineRange(fBaseLineRange);
160  fInputSignalEvent->SetRange(fIntegralRange);
161  }
162 
163  for (int n = 0; n < fInputSignalEvent->GetNumberOfSignals(); n++) {
164  TRestDetectorSignal signal;
165  TRestRawSignal* rawSignal = fInputSignalEvent->GetSignal(n);
166  signal.SetID(rawSignal->GetID());
167 
168  if (fZeroSuppression) {
169  ZeroSuppresion(rawSignal, signal);
170  } else {
171  for (int p = 0; p < int(rawSignal->GetNumberOfPoints()); p++) {
172  if (rawSignal->GetData(p) > fThreshold) {
173  signal.NewPoint(fTriggerStarts + fSampling * p, fGain * rawSignal->GetData(p));
174  }
175  }
176  }
177 
178  if (signal.GetNumberOfPoints() > 0) {
179  fOutputSignalEvent->AddSignal(signal);
180  } else {
181  rejectedSignal++;
182  }
183  }
184 
185  SetObservableValue("NSignalsRejected", rejectedSignal);
186 
187  if (fOutputSignalEvent->GetNumberOfSignals() <= 0) {
188  return nullptr;
189  }
190 
191  return fOutputSignalEvent;
192 }
193 
194 void TRestRawToDetectorSignalProcess::ZeroSuppresion(TRestRawSignal* rawSignal, TRestDetectorSignal& signal) {
195  rawSignal->InitializePointsOverThreshold(TVector2(fPointThreshold, fSignalThreshold),
196  fNPointsOverThreshold, 512);
197 
198  std::vector<Int_t> pOver = rawSignal->GetPointsOverThreshold();
199  for (unsigned int n = 0; n < pOver.size(); n++) {
200  int j = pOver[n];
201  signal.NewPoint(fTriggerStarts + fSampling * j, fGain * rawSignal->GetData(j));
202  }
203 }
A base class for any REST event.
Definition: TRestEvent.h:38
An event container for time rawdata signals with fixed length.
It defines a Short_t array with a physical parameter that evolves in time using a fixed time bin.
Int_t GetID() const
Returns the value of signal ID.
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...
std::vector< Int_t > GetPointsOverThreshold() const
Returns a std::vector containing the indexes of data points over threshold.
Double_t GetData(Int_t n) const
It returns the data value of point n including baseline correction.
Int_t GetNumberOfPoints() const
Returns the actual number of points, or size of the signal.
A process to convert a TRestRawSignalEvent into a TRestDetectorSignalEvent.
void Initialize() override
Function to initialize input/output event members and define the section name.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
The main processing event function.