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