REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRawSignalChannelActivityProcess.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
99#include "TRestRawSignalChannelActivityProcess.h"
100
101using namespace std;
102
104
109
114
120 SetSectionName(this->ClassName());
121 SetLibraryVersion(LIBRARY_VERSION);
122
123 fInputEvent = nullptr;
124}
125
136 if (!fReadOnly) {
137 // Create a unique histogram name based on signal type
138 string histogramName = "daqChannelActivityRaw_" + fChannelType;
139
140 // Create the histogram using the unique name
141 fDaqChannelsHisto = new TH1D(histogramName.c_str(), histogramName.c_str(), fDaqChannels,
143 }
144}
145
150 fInputEvent = dynamic_cast<TRestRawSignalEvent*>(inputEvent);
151
152 const auto run = GetRunInfo();
153 if (run != nullptr) {
155 }
156
157 if (fReadoutMetadata == nullptr) {
158 fReadoutMetadata = fInputEvent->GetReadoutMetadata();
159 }
160
161 if (fReadoutMetadata == nullptr && !fChannelType.empty()) {
162 cerr << "TRestRawSignalChannelActivityProcess::ProcessEvent: readout metadata is null, cannot filter "
163 "the process by signal type"
164 << endl;
165 exit(1);
166 }
167
168 for (int s = 0; s < fInputEvent->GetNumberOfSignals(); s++) {
169 const auto signal = fInputEvent->GetSignal(s);
170 if (!fChannelType.empty()) {
171 const auto channelType = fReadoutMetadata->GetTypeForChannelDaqId(signal->GetID());
172 if (fChannelType != channelType) {
173 continue;
174 }
175 }
176 // Adding signal to the channel activity histogram
177 if (!fReadOnly) {
178 Int_t daqChannel = signal->GetID();
179 fDaqChannelsHisto->Fill(daqChannel);
180 }
181 }
182
184 fAnalysisTree->PrintObservables();
185
186 return fInputEvent;
187}
188
195 if (!fReadOnly) {
196 fDaqChannelsHisto->Write();
197 }
198}
TRestAnalysisTree * fAnalysisTree
TRestRun * GetRunInfo() const
Return the pointer of the hosting TRestRun object.
bool fReadOnly
not used, keep for compatibility
A base class for any REST event.
Definition: TRestEvent.h:38
virtual void InitializeReferences(TRestRun *run)
Initialize dynamical references when loading the event from a root file.
Definition: TRestEvent.cxx:175
void SetLibraryVersion(TString version)
Set the library version of this metadata class.
TRestStringOutput::REST_Verbose_Level GetVerboseLevel()
returns the verboselevel in type of REST_Verbose_Level enumerator
void SetSectionName(std::string sName)
set the section name, clear the section content
Int_t fDaqChannels
The number of bins at the daq channels histogram.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
The main processing event function.
Int_t fDaqStartChannel
The first channel at the daq channels histogram.
void InitProcess() override
Process initialization. The ROOT TH1 histograms are created here using the limits defined in the proc...
Int_t fDaqEndChannel
The last channel at the daq channels histogram.
void EndProcess() override
Function to include required actions after all events have been processed. In this process it will ta...
void Initialize() override
Function to initialize input/output event members and define the section name.
TRestRawSignalEvent * fInputEvent
A pointer to the specific TRestRawSignalEvent input.
An event container for time rawdata signals with fixed length.
@ REST_Debug
+show the defined debug messages