126 #include "TRestDetectorSignalToHitsProcess.h"
128 #include "TRestDetectorSetup.h"
154 if (LoadConfigFromFile(configFilename) == -1) LoadDefaultConfig();
172 SetSectionName(this->ClassName());
173 SetLibraryVersion(LIBRARY_VERSION);
176 fSignalEvent =
nullptr;
188 fGas = GetMetadata<TRestDetectorGas>();
189 if (fGas !=
nullptr) {
191 RESTError <<
"A TRestDetectorGas definition was found but REST was not linked to Garfield libraries."
194 <<
"Please, remove the TRestDetectorGas definition, and add gas parameters inside the process "
195 "TRestDetectorSignalToHitsProcess"
197 if (!fGas->GetError()) fGas->SetError(
"REST was not compiled with Garfield.");
198 if (!this->GetError()) this->SetError(
"Attempt to use TRestDetectorGas without Garfield");
201 if (fGasPressure <= 0) fGasPressure = fGas->GetPressure();
202 if (fElectricField <= 0) fElectricField = fGas->GetElectricField();
204 fGas->SetPressure(fGasPressure);
205 fGas->SetElectricField(fElectricField);
207 if (fDriftVelocity <= 0) fDriftVelocity = fGas->GetDriftVelocity();
209 if (fDriftVelocity < 0) {
210 if (!this->GetError()) this->SetError(
"Drift velocity is negative.");
214 fReadout = GetMetadata<TRestDetectorReadout>();
216 if (fReadout ==
nullptr) {
217 if (!this->GetError()) this->SetError(
"The readout was not properly initialized.");
227 if (!fReadout)
return nullptr;
229 fHitsEvent->SetID(fSignalEvent->GetID());
230 fHitsEvent->SetSubID(fSignalEvent->GetSubID());
231 fHitsEvent->SetTimeStamp(fSignalEvent->GetTimeStamp());
232 fHitsEvent->SetSubEventTag(fSignalEvent->GetSubEventTag());
234 RESTDebug <<
"TRestDetectorSignalToHitsProcess. Event id : " << fHitsEvent->GetID() << RESTendl;
237 Int_t numberOfSignals = fSignalEvent->GetNumberOfSignals();
239 if (numberOfSignals == 0)
return nullptr;
241 Int_t planeID, readoutChannel = -1, readoutModule;
242 for (
int i = 0; i < numberOfSignals; i++) {
244 Int_t signalID = sgnl->GetSignalID();
247 cout <<
"Searching readout coordinates for signal ID : " << signalID << endl;
249 fReadout->GetPlaneModuleChannel(signalID, planeID, readoutModule, readoutChannel);
251 if (readoutChannel == -1) {
261 Double_t fieldZDirection = plane->
GetNormal().Z();
264 Double_t x = plane->
GetX(readoutModule, readoutChannel);
265 Double_t y = plane->
GetY(readoutModule, readoutChannel);
267 REST_HitType type = XYZ;
269 if (TMath::IsNaN(x)) {
272 }
else if (TMath::IsNaN(y)) {
277 if (fMethod ==
"onlyMax") {
278 Double_t hitTime = sgnl->GetMaxPeakTime();
279 Double_t distanceToPlane = hitTime * fDriftVelocity;
282 cout <<
"Distance to plane : " << distanceToPlane << endl;
284 Double_t z = zPosition + fieldZDirection * distanceToPlane;
286 Double_t energy = sgnl->GetMaxPeakValue();
289 cout <<
"Adding hit. Time : " << hitTime <<
" x : " << x <<
" y : " << y <<
" z : " << z
290 <<
" Energy : " << energy << endl;
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;
298 Double_t hitTime = sgnl->GetTime(bin);
299 Double_t energy = sgnl->GetData(bin);
301 Double_t distanceToPlane = hitTime * fDriftVelocity;
302 Double_t z = zPosition + fieldZDirection * distanceToPlane;
304 fHitsEvent->AddHit(x, y, z, energy, 0, type);
306 hitTime = sgnl->GetTime(binprev);
307 energy = sgnl->GetData(binprev);
309 distanceToPlane = hitTime * fDriftVelocity;
310 z = zPosition + fieldZDirection * distanceToPlane;
312 fHitsEvent->AddHit(x, y, z, energy, 0, type);
314 hitTime = sgnl->GetTime(binnext);
315 energy = sgnl->GetData(binnext);
317 distanceToPlane = hitTime * fDriftVelocity;
318 z = zPosition + fieldZDirection * distanceToPlane;
320 fHitsEvent->AddHit(x, y, z, energy, 0, type);
323 cout <<
"Distance to plane : " << distanceToPlane << endl;
324 cout <<
"Adding hit. Time : " << hitTime <<
" x : " << x <<
" y : " << y <<
" z : " << z
325 <<
" Energy : " << energy << endl;
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;
332 Double_t hitTime = sgnl->GetTime(bin);
333 Double_t energy1 = sgnl->GetData(bin);
335 Double_t distanceToPlane = hitTime * fDriftVelocity;
336 Double_t z1 = zPosition + fieldZDirection * distanceToPlane;
338 hitTime = sgnl->GetTime(binprev);
339 Double_t energy2 = sgnl->GetData(binprev);
341 distanceToPlane = hitTime * fDriftVelocity;
342 Double_t z2 = zPosition + fieldZDirection * distanceToPlane;
344 hitTime = sgnl->GetTime(binnext);
345 Double_t energy3 = sgnl->GetData(binnext);
347 distanceToPlane = hitTime * fDriftVelocity;
348 Double_t z3 = zPosition + fieldZDirection * distanceToPlane;
350 Double_t eTot = energy1 + energy2 + energy3;
352 Double_t zAvg = ((z1 * energy1) + (z2 * energy2) + (z3 * energy3)) / eTot;
354 Double_t eAvg = eTot / 3.0;
356 fHitsEvent->AddHit(x, y, zAvg, eAvg, 0, type);
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;
366 }
else if (fMethod ==
"gaussFit") {
367 TVector2 gaussFit = sgnl->GetMaxGauss();
369 Double_t hitTime = 0;
371 if (gaussFit.X() >= 0.0) {
372 hitTime = gaussFit.X();
373 Double_t distanceToPlane = hitTime * fDriftVelocity;
374 z = zPosition + fieldZDirection * distanceToPlane;
376 Double_t energy = -1.0;
377 if (gaussFit.Y() >= 0.0) {
378 energy = gaussFit.Y();
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
388 cout <<
"Adding hit. Time : " << hitTime <<
" x : " << x <<
" y : " << y <<
" z : " << z
389 <<
" Energy : " << energy << endl;
392 fHitsEvent->AddHit(x, y, z, energy, 0, type);
394 }
else if (fMethod ==
"landauFit") {
395 TVector2 landauFit = sgnl->GetMaxLandau();
397 Double_t hitTime = 0;
399 if (landauFit.X() >= 0.0) {
400 hitTime = landauFit.X();
401 Double_t distanceToPlane = hitTime * fDriftVelocity;
402 z = zPosition + fieldZDirection * distanceToPlane;
404 Double_t energy = -1.0;
405 if (landauFit.Y() >= 0.0) {
406 energy = landauFit.Y();
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
416 cout <<
"Adding hit. Time : " << hitTime <<
" x : " << x <<
" y : " << y <<
" z : " << z
417 <<
" Energy : " << energy << endl;
420 fHitsEvent->AddHit(x, y, z, energy, 0, type);
422 }
else if (fMethod ==
"agetFit") {
423 TVector2 agetFit = sgnl->GetMaxAget();
425 Double_t hitTime = 0;
427 if (agetFit.X() >= 0.0) {
428 hitTime = agetFit.X();
429 Double_t distanceToPlane = hitTime * fDriftVelocity;
430 z = zPosition + fieldZDirection * distanceToPlane;
432 Double_t energy = -1.0;
433 if (agetFit.Y() >= 0.0) {
434 energy = agetFit.Y();
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
444 cout <<
"Adding hit. Time : " << hitTime <<
" x : " << x <<
" y : " << y <<
" z : " << z
445 <<
" Energy : " << energy << endl;
448 fHitsEvent->AddHit(x, y, z, energy, 0, type);
450 }
else if (fMethod ==
"qCenter") {
451 Double_t energy_signal = 0;
452 Double_t distanceToPlane = 0;
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;
459 Double_t energy = energy_signal / sgnl->GetNumberOfPoints();
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);
467 Double_t distanceToPlane = sgnl->GetTime(j) * fDriftVelocity;
470 cout <<
"Time : " << sgnl->GetTime(j) <<
" Drift velocity : " << fDriftVelocity << endl;
471 cout <<
"Distance to plane : " << distanceToPlane << endl;
474 Double_t z = zPosition + fieldZDirection * distanceToPlane;
477 cout <<
"Adding hit. Time : " << sgnl->GetTime(j) <<
" x : " << x <<
" y : " << y
478 <<
" z : " << z << endl;
480 fHitsEvent->AddHit(x, y, z, energy, 0, type);
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()) {
490 it->second.second += sgnl->GetData(j);
492 windowMap[index] = std::make_pair(1, sgnl->GetData(j));
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;
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;
510 fHitsEvent->AddHit(x, y, z, energy, 0, type);
513 string errMsg =
"The method " + (string)fMethod +
" is not implemented!";
518 RESTDebug <<
"TRestDetectorSignalToHitsProcess. Hits added : " << fHitsEvent->GetNumberOfHits()
520 RESTDebug <<
"TRestDetectorSignalToHitsProcess. Hits total energy : " << fHitsEvent->GetTotalEnergy()
524 fHitsEvent->PrintEvent(30);
526 fHitsEvent->PrintEvent(-1);
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.";
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.
~TRestDetectorSignalToHitsProcess()
Default destructor.
TRestDetectorSignalToHitsProcess()
Default constructor.
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.
@ REST_Extreme
show everything
@ REST_Debug
+show the defined debug messages
std::string IntegerToString(Int_t n, std::string format="%d")
Gets a string from an integer.