REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorHitsToSignalProcess.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 
97 
98 #include "TRestDetectorHitsToSignalProcess.h"
99 
100 using namespace std;
101 
103 
108 
122  Initialize();
123 
124  if (LoadConfigFromFile(configFilename) == -1) {
125  LoadDefaultConfig();
126  }
127 
128  if (fReadout == nullptr) {
129  fReadout = new TRestDetectorReadout(configFilename);
130  }
131 }
132 
137  delete fReadout;
138  delete fSignalEvent;
139 }
140 
145  SetName("hitsToSignalProcess-Default");
146  SetTitle("Default config");
147 
148  cout << "Hits to signal metadata not found. Loading default values" << endl;
149 
150  fSampling = 1;
151  fElectricField = 1000;
152  fGasPressure = 10;
153 }
154 
160  SetSectionName(this->ClassName());
161  SetLibraryVersion(LIBRARY_VERSION);
162 
163  fReadout = nullptr;
164  fGas = nullptr;
165 
166  fHitsEvent = nullptr;
167  fSignalEvent = new TRestDetectorSignalEvent();
168 }
169 
175  fGas = GetMetadata<TRestDetectorGas>();
176  if (fGas != nullptr) {
177 #ifndef USE_Garfield
178  RESTError << "A TRestDetectorGas definition was found but REST was not linked to Garfield libraries."
179  << RESTendl;
180  RESTError
181  << "Please, remove the TRestDetectorGas definition, and add gas parameters inside the process "
182  "TRestDetectorHitsToSignalProcess"
183  << RESTendl;
184  if (!fGas->GetError()) {
185  fGas->SetError("REST was not compiled with Garfield.");
186  }
187  if (!this->GetError()) {
188  this->SetError("Attempt to use TRestDetectorGas without Garfield");
189  }
190 #endif
191  if (fGasPressure <= 0) {
192  fGasPressure = fGas->GetPressure();
193  }
194  if (fElectricField <= 0) {
195  fElectricField = fGas->GetElectricField();
196  }
197 
198  fGas->SetPressure(fGasPressure);
199  fGas->SetElectricField(fElectricField);
200 
201  if (fDriftVelocity <= 0) {
202  fDriftVelocity = fGas->GetDriftVelocity();
203  }
204  }
205 
206  if (fDriftVelocity <= 0) {
207  RESTError
208  << "Drift velocity not defined. Please, define it in the process metadata or in TRestDetectorGas"
209  << RESTendl;
210  exit(1);
211  }
212 
213  fReadout = GetMetadata<TRestDetectorReadout>();
214  if (fReadout == nullptr) {
215  if (!this->GetError()) {
216  this->SetError("The readout was not properly initialized.");
217  }
218  }
219 }
220 
225  fHitsEvent = (TRestDetectorHitsEvent*)inputEvent;
226  fSignalEvent->SetEventInfo(fHitsEvent);
227 
228  if (!fReadout) {
229  return nullptr;
230  }
231 
232  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
233  cout << "Number of hits : " << fHitsEvent->GetNumberOfHits() << endl;
234  cout << "--------------------------" << endl;
235  }
236 
237  for (unsigned int hit = 0; hit < fHitsEvent->GetNumberOfHits(); hit++) {
238  Double_t x = fHitsEvent->GetX(hit);
239  Double_t y = fHitsEvent->GetY(hit);
240  Double_t z = fHitsEvent->GetZ(hit);
241  Double_t t = fHitsEvent->GetTime(hit);
242 
243  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Extreme && hit < 20) {
244  cout << "Hit : " << hit << " x : " << x << " y : " << y << " z : " << z << " t : " << t << endl;
245  }
246 
247  for (int p = 0; p < fReadout->GetNumberOfReadoutPlanes(); p++) {
248  const auto [daqId, moduleId, channelId] = fReadout->GetHitsDaqChannelAtReadoutPlane({x, y, z}, p);
249  TRestDetectorReadoutPlane* plane = fReadout->GetReadoutPlane(p);
250 
251  if (daqId >= 0) {
252  auto channel = fReadout->GetReadoutChannelWithDaqID(daqId);
253 
254  const bool isVeto = (channel->GetChannelType() == "veto");
255 
256  Double_t energy = fHitsEvent->GetEnergy(hit);
257  const auto distance = plane->GetDistanceTo({x, y, z});
258  if (distance < 0) {
259  RESTError << "TRestDetectorHitsToSignalProcess: Negative distance to readout plane. "
260  "This should not happen."
261  << RESTendl;
262  exit(1);
263  }
264 
265  auto velocity = fDriftVelocity;
266  if (isVeto) {
267  velocity = REST_Physics::lightSpeed;
268  }
269 
270  Double_t time = t + distance / velocity;
271 
272  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug && hit < 20) {
273  cout << "Module : " << moduleId << " Channel : " << channelId << " daq ID : " << daqId
274  << endl;
275  }
276 
277  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug && hit < 20) {
278  cout << "Energy : " << energy << " time : " << time << endl;
279  }
280  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Extreme && hit < 20) {
281  printf(
282  " TRestDetectorHitsToSignalProcess: x %lf y %lf z %lf energy %lf t %lf "
283  "fDriftVelocity %lf fSampling %lf time %lf\n",
284  x, y, z, energy, t, fDriftVelocity, fSampling, time);
285  }
286 
287  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Extreme) {
288  cout << "Drift velocity : " << fDriftVelocity << " mm/us" << endl;
289  }
290  time = ((Int_t)(time / fSampling)) * fSampling; // now time is in unit "us", but dispersed
291 
292  fSignalEvent->AddChargeToSignal(daqId, time, energy);
293 
294  auto signal = fSignalEvent->GetSignalById(daqId);
295  signal->SetSignalName(channel->GetChannelName());
296  signal->SetSignalType(channel->GetChannelType());
297 
298  } else {
299  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
300  RESTDebug << "TRestDetectorHitsToSignalProcess. Readout channel not find for position ("
301  << x << ", " << y << ", " << z << ")!" << RESTendl;
302  }
303  }
304  }
305  }
306 
307  fSignalEvent->SortSignals();
308 
309  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
310  cout << "TRestDetectorHitsToSignalProcess : Number of signals added : "
311  << fSignalEvent->GetNumberOfSignals() << endl;
312  cout << "TRestDetectorHitsToSignalProcess : Total signals integral : " << fSignalEvent->GetIntegral()
313  << endl;
314  }
315 
316  if (fSignalEvent->GetNumberOfSignals() == 0) {
317  return nullptr;
318  }
319 
320  return fSignalEvent;
321 }
A process to transform a x,y,z coordinate hits into daq identified physical time signals.
void InitProcess() override
Process initialization. This process accesses the information inside TRestGeant4Metadata to identify ...
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
The main processing event function.
void LoadDefaultConfig()
Function to load the default config in absence of RML input.
void Initialize() override
Function to initialize input/output event members and define the section name.
Double_t GetDistanceTo(const TVector3 &pos) const
Returns the perpendicular distance to the readout plane from a given position pos.
A metadata class to generate/store a readout description.
A base class for any REST event.
Definition: TRestEvent.h:38
void SetEventInfo(TRestEvent *eve)
Definition: TRestEvent.cxx:137
@ REST_Debug
+show the defined debug messages
constexpr double lightSpeed
Speed of light in m/s.
Definition: TRestPhysics.h:41