REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
TRestDetectorSignalToHitsProcess.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 
126 #include "TRestDetectorSignalToHitsProcess.h"
127 
128 #include "TRestDetectorSetup.h"
129 
130 using namespace std;
131 
133 
138 
152  Initialize();
153 
154  if (LoadConfigFromFile(configFilename) == -1) LoadDefaultConfig();
155 }
156 
160 void TRestDetectorSignalToHitsProcess::LoadDefaultConfig() { SetTitle("Default config"); }
161 
166 
172  SetSectionName(this->ClassName());
173  SetLibraryVersion(LIBRARY_VERSION);
174 
175  fHitsEvent = new TRestDetectorHitsEvent();
176  fSignalEvent = nullptr;
177 
178  fGas = nullptr;
179  fReadout = nullptr;
180 }
181 
188  fGas = GetMetadata<TRestDetectorGas>();
189  if (fGas != nullptr) {
190 #ifndef USE_Garfield
191  RESTError << "A TRestDetectorGas definition was found but REST was not linked to Garfield libraries."
192  << RESTendl;
193  RESTError
194  << "Please, remove the TRestDetectorGas definition, and add gas parameters inside the process "
195  "TRestDetectorSignalToHitsProcess"
196  << RESTendl;
197  if (!fGas->GetError()) fGas->SetError("REST was not compiled with Garfield.");
198  if (!this->GetError()) this->SetError("Attempt to use TRestDetectorGas without Garfield");
199 #endif
200 
201  if (fGasPressure <= 0) fGasPressure = fGas->GetPressure();
202  if (fElectricField <= 0) fElectricField = fGas->GetElectricField();
203 
204  fGas->SetPressure(fGasPressure);
205  fGas->SetElectricField(fElectricField);
206 
207  if (fDriftVelocity <= 0) fDriftVelocity = fGas->GetDriftVelocity();
208  } else {
209  if (fDriftVelocity < 0) {
210  if (!this->GetError()) this->SetError("Drift velocity is negative.");
211  }
212  }
213 
214  fReadout = GetMetadata<TRestDetectorReadout>();
215 
216  if (fReadout == nullptr) {
217  if (!this->GetError()) this->SetError("The readout was not properly initialized.");
218  }
219 }
220 
225  fSignalEvent = (TRestDetectorSignalEvent*)inputEvent;
226 
227  if (!fReadout) return nullptr;
228 
229  fHitsEvent->SetID(fSignalEvent->GetID());
230  fHitsEvent->SetSubID(fSignalEvent->GetSubID());
231  fHitsEvent->SetTimeStamp(fSignalEvent->GetTimeStamp());
232  fHitsEvent->SetSubEventTag(fSignalEvent->GetSubEventTag());
233 
234  RESTDebug << "TRestDetectorSignalToHitsProcess. Event id : " << fHitsEvent->GetID() << RESTendl;
235  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Extreme) fSignalEvent->PrintEvent();
236 
237  Int_t numberOfSignals = fSignalEvent->GetNumberOfSignals();
238 
239  if (numberOfSignals == 0) return nullptr;
240 
241  Int_t planeID, readoutChannel = -1, readoutModule;
242  for (int i = 0; i < numberOfSignals; i++) {
243  TRestDetectorSignal* sgnl = fSignalEvent->GetSignal(i);
244  Int_t signalID = sgnl->GetSignalID();
245 
246  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug)
247  cout << "Searching readout coordinates for signal ID : " << signalID << endl;
248 
249  fReadout->GetPlaneModuleChannel(signalID, planeID, readoutModule, readoutChannel);
250 
251  if (readoutChannel == -1) {
252  // cout << "REST Warning : Readout channel not found for daq ID : " << signalID << endl;
253  continue;
254  }
256 
257  TRestDetectorReadoutPlane* plane = fReadout->GetReadoutPlaneWithID(planeID);
258 
259  // For the moment this will only be valid for a TPC with its axis (field
260  // direction) being in z
261  Double_t fieldZDirection = plane->GetNormal().Z();
262  Double_t zPosition = plane->GetPosition().Z();
263 
264  Double_t x = plane->GetX(readoutModule, readoutChannel);
265  Double_t y = plane->GetY(readoutModule, readoutChannel);
266 
267  REST_HitType type = XYZ;
268  TRestDetectorReadoutModule* mod = plane->GetModuleByID(readoutModule);
269  if (TMath::IsNaN(x)) {
270  x = mod->GetPlaneCoordinates(TVector2(mod->GetSize().X() / 2, mod->GetSize().Y() / 2)).X();
271  type = YZ;
272  } else if (TMath::IsNaN(y)) {
273  y = mod->GetPlaneCoordinates(TVector2(mod->GetSize().X() / 2, mod->GetSize().Y() / 2)).Y();
274  type = XZ;
275  }
276 
277  if (fMethod == "onlyMax") {
278  Double_t hitTime = sgnl->GetMaxPeakTime();
279  Double_t distanceToPlane = hitTime * fDriftVelocity;
280 
281  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug)
282  cout << "Distance to plane : " << distanceToPlane << endl;
283 
284  Double_t z = zPosition + fieldZDirection * distanceToPlane;
285 
286  Double_t energy = sgnl->GetMaxPeakValue();
287 
288  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug)
289  cout << "Adding hit. Time : " << hitTime << " x : " << x << " y : " << y << " z : " << z
290  << " Energy : " << energy << endl;
291 
292  fHitsEvent->AddHit(x, y, z, energy, 0, type);
293  } else if (fMethod == "tripleMax") {
294  Int_t bin = sgnl->GetMaxIndex();
295  int binprev = (bin - 1) < 0 ? bin : bin - 1;
296  int binnext = (bin + 1) > sgnl->GetNumberOfPoints() - 1 ? bin : bin + 1;
297 
298  Double_t hitTime = sgnl->GetTime(bin);
299  Double_t energy = sgnl->GetData(bin);
300 
301  Double_t distanceToPlane = hitTime * fDriftVelocity;
302  Double_t z = zPosition + fieldZDirection * distanceToPlane;
303 
304  fHitsEvent->AddHit(x, y, z, energy, 0, type);
305 
306  hitTime = sgnl->GetTime(binprev);
307  energy = sgnl->GetData(binprev);
308 
309  distanceToPlane = hitTime * fDriftVelocity;
310  z = zPosition + fieldZDirection * distanceToPlane;
311 
312  fHitsEvent->AddHit(x, y, z, energy, 0, type);
313 
314  hitTime = sgnl->GetTime(binnext);
315  energy = sgnl->GetData(binnext);
316 
317  distanceToPlane = hitTime * fDriftVelocity;
318  z = zPosition + fieldZDirection * distanceToPlane;
319 
320  fHitsEvent->AddHit(x, y, z, energy, 0, type);
321 
322  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
323  cout << "Distance to plane : " << distanceToPlane << endl;
324  cout << "Adding hit. Time : " << hitTime << " x : " << x << " y : " << y << " z : " << z
325  << " Energy : " << energy << endl;
326  }
327  } else if (fMethod == "tripleMaxAverage") {
328  Int_t bin = sgnl->GetMaxIndex();
329  int binprev = (bin - 1) < 0 ? bin : bin - 1;
330  int binnext = (bin + 1) > sgnl->GetNumberOfPoints() - 1 ? bin : bin + 1;
331 
332  Double_t hitTime = sgnl->GetTime(bin);
333  Double_t energy1 = sgnl->GetData(bin);
334 
335  Double_t distanceToPlane = hitTime * fDriftVelocity;
336  Double_t z1 = zPosition + fieldZDirection * distanceToPlane;
337 
338  hitTime = sgnl->GetTime(binprev);
339  Double_t energy2 = sgnl->GetData(binprev);
340 
341  distanceToPlane = hitTime * fDriftVelocity;
342  Double_t z2 = zPosition + fieldZDirection * distanceToPlane;
343 
344  hitTime = sgnl->GetTime(binnext);
345  Double_t energy3 = sgnl->GetData(binnext);
346 
347  distanceToPlane = hitTime * fDriftVelocity;
348  Double_t z3 = zPosition + fieldZDirection * distanceToPlane;
349 
350  Double_t eTot = energy1 + energy2 + energy3;
351 
352  Double_t zAvg = ((z1 * energy1) + (z2 * energy2) + (z3 * energy3)) / eTot;
353  // Double_t zAvg = (z1 + z2 + z3) / 3.0;
354  Double_t eAvg = eTot / 3.0;
355 
356  fHitsEvent->AddHit(x, y, zAvg, eAvg, 0, type);
357 
358  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
359  cout << "Distance to plane: " << distanceToPlane << endl;
360  cout << "Adding hit. Time : " << hitTime << " x : " << x << " y : " << y << " z : " << zAvg
361  << " Energy : " << eAvg << endl;
362  cout << "z1, z2, z3 = " << z1 << ", " << z2 << ", " << z3 << endl;
363  cout << "E1, E2, E3 = " << energy1 << ", " << energy2 << ", " << energy3 << endl;
364  }
365 
366  } else if (fMethod == "gaussFit") {
367  TVector2 gaussFit = sgnl->GetMaxGauss();
368 
369  Double_t hitTime = 0;
370  Double_t z = -1.0;
371  if (gaussFit.X() >= 0.0) {
372  hitTime = gaussFit.X();
373  Double_t distanceToPlane = hitTime * fDriftVelocity;
374  z = zPosition + fieldZDirection * distanceToPlane;
375  }
376  Double_t energy = -1.0;
377  if (gaussFit.Y() >= 0.0) {
378  energy = gaussFit.Y();
379  }
380 
381  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
382  cout << "Signal event : " << sgnl->GetSignalID()
383  << "--------------------------------------------------------" << endl;
384  cout << "GausFit : time bin " << gaussFit.X() << " and energy : " << gaussFit.Y() << endl;
385  cout << "Signal to hit info : zPosition : " << zPosition
386  << "; fieldZDirection : " << fieldZDirection << " and driftV : " << fDriftVelocity
387  << endl;
388  cout << "Adding hit. Time : " << hitTime << " x : " << x << " y : " << y << " z : " << z
389  << " Energy : " << energy << endl;
390  }
391 
392  fHitsEvent->AddHit(x, y, z, energy, 0, type);
393 
394  } else if (fMethod == "landauFit") {
395  TVector2 landauFit = sgnl->GetMaxLandau();
396 
397  Double_t hitTime = 0;
398  Double_t z = -1.0;
399  if (landauFit.X() >= 0.0) {
400  hitTime = landauFit.X();
401  Double_t distanceToPlane = hitTime * fDriftVelocity;
402  z = zPosition + fieldZDirection * distanceToPlane;
403  }
404  Double_t energy = -1.0;
405  if (landauFit.Y() >= 0.0) {
406  energy = landauFit.Y();
407  }
408 
409  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
410  cout << "Signal event : " << sgnl->GetSignalID()
411  << "--------------------------------------------------------" << endl;
412  cout << "landauFit : time bin " << landauFit.X() << " and energy : " << landauFit.Y() << endl;
413  cout << "Signal to hit info : zPosition : " << zPosition
414  << "; fieldZDirection : " << fieldZDirection << " and driftV : " << fDriftVelocity
415  << endl;
416  cout << "Adding hit. Time : " << hitTime << " x : " << x << " y : " << y << " z : " << z
417  << " Energy : " << energy << endl;
418  }
419 
420  fHitsEvent->AddHit(x, y, z, energy, 0, type);
421 
422  } else if (fMethod == "agetFit") {
423  TVector2 agetFit = sgnl->GetMaxAget();
424 
425  Double_t hitTime = 0;
426  Double_t z = -1.0;
427  if (agetFit.X() >= 0.0) {
428  hitTime = agetFit.X();
429  Double_t distanceToPlane = hitTime * fDriftVelocity;
430  z = zPosition + fieldZDirection * distanceToPlane;
431  }
432  Double_t energy = -1.0;
433  if (agetFit.Y() >= 0.0) {
434  energy = agetFit.Y();
435  }
436 
437  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
438  cout << "Signal event : " << sgnl->GetSignalID()
439  << "--------------------------------------------------------" << endl;
440  cout << "agetFit : time bin " << agetFit.X() << " and energy : " << agetFit.Y() << endl;
441  cout << "Signal to hit info : zPosition : " << zPosition
442  << "; fieldZDirection : " << fieldZDirection << " and driftV : " << fDriftVelocity
443  << endl;
444  cout << "Adding hit. Time : " << hitTime << " x : " << x << " y : " << y << " z : " << z
445  << " Energy : " << energy << endl;
446  }
447 
448  fHitsEvent->AddHit(x, y, z, energy, 0, type);
449 
450  } else if (fMethod == "qCenter") {
451  Double_t energy_signal = 0;
452  Double_t distanceToPlane = 0;
453 
454  for (int j = 0; j < sgnl->GetNumberOfPoints(); j++) {
455  Double_t energy_point = sgnl->GetData(j);
456  energy_signal += energy_point;
457  distanceToPlane += sgnl->GetTime(j) * fDriftVelocity * energy_point;
458  }
459  Double_t energy = energy_signal / sgnl->GetNumberOfPoints();
460 
461  Double_t z = zPosition + fieldZDirection * (distanceToPlane / energy_signal);
462  fHitsEvent->AddHit(x, y, z, energy, 0, type);
463  } else if (fMethod == "all") {
464  for (int j = 0; j < sgnl->GetNumberOfPoints(); j++) {
465  Double_t energy = sgnl->GetData(j);
466 
467  Double_t distanceToPlane = sgnl->GetTime(j) * fDriftVelocity;
468 
469  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug) {
470  cout << "Time : " << sgnl->GetTime(j) << " Drift velocity : " << fDriftVelocity << endl;
471  cout << "Distance to plane : " << distanceToPlane << endl;
472  }
473 
474  Double_t z = zPosition + fieldZDirection * distanceToPlane;
475 
476  if (GetVerboseLevel() >= TRestStringOutput::REST_Verbose_Level::REST_Debug)
477  cout << "Adding hit. Time : " << sgnl->GetTime(j) << " x : " << x << " y : " << y
478  << " z : " << z << endl;
479 
480  fHitsEvent->AddHit(x, y, z, energy, 0, type);
481  }
482  } else if (fMethod == "intwindow") {
483  Int_t nPoints = sgnl->GetNumberOfPoints();
484  std::map<int, std::pair<int, double> > windowMap;
485  for (int j = 0; j < nPoints; j++) {
486  int index = sgnl->GetTime(j) / fIntWindow;
487  auto it = windowMap.find(index);
488  if (it != windowMap.end()) {
489  it->second.first++;
490  it->second.second += sgnl->GetData(j);
491  } else {
492  windowMap[index] = std::make_pair(1, sgnl->GetData(j));
493  }
494  }
495 
496  for (const auto& [index, pair] : windowMap) {
497  Double_t hitTime = index * fIntWindow + fIntWindow / 2.;
498  Double_t energy = pair.second / pair.first;
499  if (energy < fThreshold) continue;
500  RESTDebug << "TimeBin " << index << " Time " << hitTime << " Charge: " << energy
501  << " Thr: " << (fThreshold) << RESTendl;
502  Double_t distanceToPlane = hitTime * fDriftVelocity;
503  Double_t z = zPosition + fieldZDirection * distanceToPlane;
504 
505  RESTDebug << "Time : " << hitTime << " Drift velocity : " << fDriftVelocity
506  << "\nDistance to plane : " << distanceToPlane << RESTendl;
507  RESTDebug << "Adding hit. Time : " << hitTime << " x : " << x << " y : " << y << " z : " << z
508  << " type " << type << RESTendl;
509 
510  fHitsEvent->AddHit(x, y, z, energy, 0, type);
511  }
512  } else {
513  string errMsg = "The method " + (string)fMethod + " is not implemented!";
514  SetError(errMsg);
515  }
516  }
517 
518  RESTDebug << "TRestDetectorSignalToHitsProcess. Hits added : " << fHitsEvent->GetNumberOfHits()
519  << RESTendl;
520  RESTDebug << "TRestDetectorSignalToHitsProcess. Hits total energy : " << fHitsEvent->GetTotalEnergy()
521  << RESTendl;
522 
523  if (this->GetVerboseLevel() == TRestStringOutput::REST_Verbose_Level::REST_Debug) {
524  fHitsEvent->PrintEvent(30);
525  } else if (this->GetVerboseLevel() == TRestStringOutput::REST_Verbose_Level::REST_Extreme) {
526  fHitsEvent->PrintEvent(-1);
527  }
528 
529  if (fHitsEvent->GetNumberOfHits() <= 0) {
530  string errMsg = "Last event id: " + IntegerToString(fHitsEvent->GetID()) +
531  ". Failed to find readout positions in channel to hit conversion.";
532  SetWarning(errMsg);
533  return nullptr;
534  }
535 
536  return fHitsEvent;
537 }
TVector2 GetSize() const
Returns the module size (x, y) in mm.
TVector2 GetPlaneCoordinates(const TVector2 &p)
TRestDetectorReadoutModule * GetModuleByID(Int_t modID)
Returns a pointer to a module using its internal module id.
TVector3 GetNormal() const
Returns a TVector3 with a vector normal to the readout plane.
TVector3 GetPosition() const
Returns a TVector3 with the readout plane position.
Double_t GetX(Int_t modID, Int_t chID)
Returns the X coordinate of a given channel in a given module using their internal module and channel...
Double_t GetY(Int_t modID, Int_t chID)
Returns the Y coordinate of a given channel in a given module using their internal module and channel...
A process to transform a daq channel and physical time to spatial coordinates.
void InitProcess() override
Process initialization. Observable names are interpreted and auxiliar observable members,...
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
The main processing event function.
void Initialize() override
Function to initialize input/output event members and define the section name.
void LoadDefaultConfig()
Function to load the default config in absence of RML input.
A base class for any REST event.
Definition: TRestEvent.h:38
@ REST_Debug
+show the defined debug messages
std::string IntegerToString(Int_t n, std::string format="%d")
Gets a string from an integer.