REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestRawFeminosRootToSignalProcess.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
64
65#include "TRestRawFeminosRootToSignalProcess.h"
66
67using namespace std;
68
70
71TRestRawFeminosRootToSignalProcess::TRestRawFeminosRootToSignalProcess() { Initialize(); }
72
73TRestRawFeminosRootToSignalProcess::TRestRawFeminosRootToSignalProcess(const char* configFilename) {
74 Initialize();
75}
76
77TRestRawFeminosRootToSignalProcess::~TRestRawFeminosRootToSignalProcess() {}
78
80 delete fSignalEvent;
81 fSignalEvent = new TRestRawSignalEvent();
82
83 fIsExternal = true; // We need this in order to prevent error since we are not reading a rest root file
84 fSingleThreadOnly = true;
85}
86
88 // print input aqs file
89 const auto inputFilename = fRunInfo->GetInputFileName(0);
90
91 // verify it is a root file reading the last 5 characters
92 if (inputFilename.substr(inputFilename.size() - 5) != ".root") {
93 cerr << "TRestRawFeminosRootToSignalProcess::InitProcess: Input file is not a root file" << endl;
94 exit(1);
95 }
96
97 fInputFile = TFile::Open(inputFilename.c_str(), "READ");
98 if (!fInputFile) {
99 cerr << "TRestRawFeminosRootToSignalProcess::InitProcess: Error opening input file" << endl;
100 exit(1);
101 }
102 fInputFile->cd();
103
104 fInputRunTree = fInputFile->Get<TTree>("run");
105 fInputEventTree = fInputFile->Get<TTree>("events");
106
107 if (!fInputRunTree || !fInputEventTree) {
108 cerr << "TRestRawFeminosRootToSignalProcess::InitProcess: Error opening input trees" << endl;
109 exit(1);
110 }
111
112 fRunInfo->SetFeminosDaqTotalEvents(fInputEventTree->GetEntries());
113
114 fInputEventTree->SetBranchAddress("timestamp", &fInputEventTreeTimestamp);
115 fInputEventTree->SetBranchAddress("signal_ids", &fInputEventTreeSignalIds);
116 fInputEventTree->SetBranchAddress("signal_values", &fInputEventTreeSignalValues);
117}
118
120 fSignalEvent->Initialize();
121
122 fInputEventTree->GetEntry(fInputTreeEntry);
123
124 fSignalEvent->SetID(fInputTreeEntry);
125
126 // TODO: double check this is correct
127 fSignalEvent->SetTime(fInputEventTreeTimestamp / 1000, fInputEventTreeTimestamp % 1000);
128
129 for (size_t i = 0; i < fInputEventTreeSignalIds->size(); i++) {
130 auto signal = TRestRawSignal();
131 signal.Initialize();
132
133 const auto id = fInputEventTreeSignalIds->at(i);
134 signal.SetID(id);
135
136 for (int j = 0; j < 512; j++) {
137 signal.AddPoint(short(fInputEventTreeSignalValues->at(i * 512 + j)));
138 }
139
140 fSignalEvent->AddSignal(signal);
141 }
142
143 fInputTreeEntry += 1;
144
145 if (fInputTreeEntry >= fInputEventTree->GetEntries()) {
146 return nullptr;
147 }
148
149 return fSignalEvent;
150}
bool fIsExternal
It defines if the process reads event data from an external source.
TRestRun * fRunInfo
< Pointer to TRestRun object where to find metadata.
A base class for any REST event.
Definition: TRestEvent.h:38
void SetTime(Double_t time)
Definition: TRestEvent.cxx:85
void InitProcess() override
To be executed at the beginning of the run (outside event loop)
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
Process one event.
void Initialize() override
Making default settings.
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.