REST-for-Physics  v2.3
Rare Event Searches ToolKit for Physics
Loading...
Searching...
No Matches
TRestRawSignalAnalysisProcess.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
249
250#include "TRestRawSignalAnalysisProcess.h"
251
252using namespace std;
253
255
260
265
271 SetSectionName(this->ClassName());
272 SetLibraryVersion(LIBRARY_VERSION);
273 fSignalEvent = nullptr;
274}
275
280 if (fSignalsRange.X() != -1 && fSignalsRange.Y() != -1) {
281 fRangeEnabled = true;
282 }
283}
284
287 const auto filterType = GetParameter("channelType", "");
288 if (!filterType.empty()) {
289 fChannelTypes.insert(filterType);
290 }
291}
292
297 fSignalEvent = (TRestRawSignalEvent*)inputEvent;
299 auto event = fSignalEvent->GetSignalEventForTypes(fChannelTypes, fReadoutMetadata);
300
301 // we save some complex typed analysis result
302 map<int, Double_t> baseline;
303 map<int, Double_t> baselinesigma;
304 map<int, Double_t> ampsgn_maxmethod;
305 map<int, Double_t> ampsgn_intmethod;
306 map<int, int> risetime;
307 map<int, int> peak_time;
308 map<int, int> npointsot;
309 vector<int> saturatedchnId;
310
311 baseline.clear();
312 baselinesigma.clear();
313 ampsgn_maxmethod.clear();
314 ampsgn_intmethod.clear();
315 risetime.clear();
316 npointsot.clear();
317
318 Int_t nGoodSignals = 0;
319
322 // This will affect the calculation of observables, but not the stored
323 // TRestRawSignal data.
325 event.SetRange(fIntegralRange);
326
327 for (int s = 0; s < event.GetNumberOfSignals(); s++) {
328 TRestRawSignal* sgnl = event.GetSignal(s);
329
333
334 if (fRangeEnabled && (sgnl->GetID() < fSignalsRange.X() || sgnl->GetID() > fSignalsRange.Y())) {
335 continue;
336 }
337
338 // We do not want that signals that are not identified as such contribute to
339 // define our observables
340 // nkx: we still need to store all the signals in baseline/rise time maps in
341 // case for noise analysis
342 // if (sgnl->GetPointsOverThreshold().size() < 2) continue;
343 if (sgnl->GetPointsOverThreshold().size() >= 2) nGoodSignals++;
344
345 // Now TRestRawSignal returns directly baseline subtracted values
346 baseline[sgnl->GetID()] = sgnl->GetBaseLine();
347 baselinesigma[sgnl->GetID()] = sgnl->GetBaseLineSigma();
348 ampsgn_intmethod[sgnl->GetID()] = sgnl->GetThresholdIntegral();
349 ampsgn_maxmethod[sgnl->GetID()] = sgnl->GetMaxPeakValue();
350 risetime[sgnl->GetID()] = sgnl->GetRiseTime();
351 peak_time[sgnl->GetID()] = sgnl->GetMaxPeakBin();
352 npointsot[sgnl->GetID()] = sgnl->GetPointsOverThreshold().size();
353 if (sgnl->IsADCSaturation()) saturatedchnId.push_back(sgnl->GetID());
354 }
355
356 SetObservableValue("pointsoverthres_map", npointsot);
357 SetObservableValue("risetime_map", risetime);
358 SetObservableValue("peak_time_map", peak_time);
359 SetObservableValue("baseline_map", baseline);
360 SetObservableValue("baselinesigma_map", baselinesigma);
361 SetObservableValue("max_amplitude_map", ampsgn_maxmethod);
362 SetObservableValue("thr_integral_map", ampsgn_intmethod);
363 SetObservableValue("SaturatedChannelID", saturatedchnId);
364
365 Double_t baseLineMean = event.GetBaseLineAverage();
366 SetObservableValue("BaseLineMean", baseLineMean);
367
368 Double_t baseLineSigma = event.GetBaseLineSigmaAverage();
369 SetObservableValue("BaseLineSigmaMean", baseLineSigma);
370
371 Double_t timeDelay = event.GetMaxTime() - event.GetMinTime();
372 SetObservableValue("TimeBinsLength", timeDelay);
373
374 Int_t nSignals = event.GetNumberOfSignals();
375 SetObservableValue("NumberOfSignals", nSignals);
376 SetObservableValue("NumberOfGoodSignals", nGoodSignals);
377
378 // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
379 // SubstractBaselines
380 // After this the signal gets zero-ed, for the following analysis
381 // Keep in mind, to add raw signal analysis, we must write code at before
382 // This is where most of the problems occur
383 // !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
384 // Javier: I believe we should not substract baseline in the analysis process
385 // then ...
386 // ... of course we need to consider baseline substraction for each
387 // observable. TRestRawSignal methods
388 // should do that internally. I have updated that to be like that, but we need
389 // to be with open eyes for
390 // some period.
391 // Baseline substraction will always happen when we transfer a TRestRawSignal
392 // to TRestDetectorSignal
393 //
394 // We do not substract baselines then now, as it was done before
395 //
396 // fSignalEvent->SubstractBaselines(fBaseLineRange.X(), fBaseLineRange.Y());
397 //
398 // Methods in TRestRawSignal have been updated to consider baseline.
399 // TRestRawSignal now implements that internally. We need to define the
400 // baseline range, and the range
401 // where calculations take place. All we need is to call at some point to the
402 // following methods.
403 //
404 // TRestRawSignalEvent::SetBaseLineRange and TRestRawSignalEvent::SetRange.
405 //
406 // Then, if any method accepts a different range it will be given in the
407 // method name,
408 // for example: GetIntegralInRange( Int_t startBin, Int_t endBin );
409 //
410
411 Double_t fullIntegral = event.GetIntegral();
412 SetObservableValue("FullIntegral", fullIntegral);
413
414 Double_t thrIntegral = event.GetThresholdIntegral();
415 SetObservableValue("ThresholdIntegral", thrIntegral);
416
417 Double_t riseSlope = event.GetRiseSlope();
418 SetObservableValue("RiseSlopeAvg", riseSlope);
419
420 Double_t slopeIntegral = event.GetSlopeIntegral();
421 SetObservableValue("SlopeIntegral", slopeIntegral);
422
423 Double_t rateOfChange = riseSlope / slopeIntegral;
424 if (slopeIntegral == 0) rateOfChange = 0;
425 SetObservableValue("RateOfChangeAvg", rateOfChange);
426
427 Double_t riseTime = event.GetRiseTime();
428 SetObservableValue("RiseTimeAvg", riseTime);
429
430 Double_t tripleMaxIntegral = event.GetTripleMaxIntegral();
431 SetObservableValue("TripleMaxIntegral", tripleMaxIntegral);
432
433 Double_t integralRatio = (fullIntegral - thrIntegral) / (fullIntegral + thrIntegral);
434 SetObservableValue("IntegralBalance", integralRatio);
435
436 Double_t maxValue = 0;
437 Double_t minValue = 1.e6;
438 Double_t maxValueIntegral = 0;
439 Double_t minDownValue = 1.e6;
440
441 Double_t minPeakTime = 1000; // TODO substitute this for something better
442 Double_t maxPeakTime = 0;
443 Double_t peakTimeAverage = 0;
444
445 for (int s = 0; s < event.GetNumberOfSignals(); s++) {
446 TRestRawSignal* sgnl = event.GetSignal(s);
447
448 if (fRangeEnabled && (sgnl->GetID() < fSignalsRange.X() || sgnl->GetID() > fSignalsRange.Y()))
449 continue;
450
451 if (sgnl->GetPointsOverThreshold().size() > 1) {
452 Double_t value = event.GetSignal(s)->GetMaxValue();
453 maxValueIntegral += value;
454
455 if (value > maxValue) maxValue = value;
456 if (value < minValue) minValue = value;
457
458 Double_t peakBin = sgnl->GetMaxPeakBin();
459 peakTimeAverage += peakBin;
460
461 if (minPeakTime > peakBin) minPeakTime = peakBin;
462 if (maxPeakTime < peakBin) maxPeakTime = peakBin;
463 }
464 Double_t mindownvalue = event.GetSignal(s)->GetMinValue();
465 if (mindownvalue < minDownValue) {
466 minDownValue = mindownvalue;
467 }
468 }
469
470 if (nGoodSignals > 0) peakTimeAverage /= nGoodSignals;
471
472 Double_t ampIntRatio = thrIntegral / maxValueIntegral;
473 if (maxValueIntegral == 0) ampIntRatio = 0;
474
475 SetObservableValue("AmplitudeIntegralRatio", ampIntRatio);
476 SetObservableValue("MinPeakAmplitude", minValue);
477 SetObservableValue("MaxPeakAmplitude", maxValue);
478 SetObservableValue("PeakAmplitudeIntegral", maxValueIntegral);
479
480 SetObservableValue("MinEventValue", minDownValue);
481
482 Double_t amplitudeRatio = maxValueIntegral / maxValue;
483 if (maxValue == 0) amplitudeRatio = 0;
484
485 SetObservableValue("AmplitudeRatio", amplitudeRatio);
486 SetObservableValue("MaxPeakTime", maxPeakTime);
487 SetObservableValue("MinPeakTime", minPeakTime);
488
489 Double_t peakTimeDelay = maxPeakTime - minPeakTime;
490
491 SetObservableValue("MaxPeakTimeDelay", peakTimeDelay);
492 SetObservableValue("AveragePeakTime", peakTimeAverage);
493
495 for (const auto& i : fObservablesDefined) {
496 fAnalysisTree->PrintObservable(i.second);
497 }
498 }
499
500 // If cut condition matches the event will be not registered.
501 if (ApplyCut()) {
502 return nullptr;
503 }
504
505 return fSignalEvent;
506}
virtual void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
TRestAnalysisTree * fAnalysisTree
TRestRun * GetRunInfo() const
Return the pointer of the hosting TRestRun object.
std::map< std::string, int > fObservablesDefined
Stores the list of all the appeared process observables in the code.
void SetObservableValue(const std::string &name, const T &value)
Set observable value for AnalysisTree.
A base class for any REST event.
Definition TRestEvent.h:38
virtual void InitializeReferences(TRestRun *run)
Initialize dynamical references when loading the event from a root file.
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 SetSectionName(std::string sName)
set the section name, clear the section content
std::string GetParameter(std::string parName, TiXmlElement *e, TString defaultValue=PARAMETER_NOT_FOUND_STR)
Returns the value for the parameter named parName in the given section.
An analysis process to extract valuable information from a TRestRawSignalEvent.
void InitProcess() override
Process initialization.
TVector2 fSignalsRange
It defines the signals id range where analysis is applied.
TRestEvent * ProcessEvent(TRestEvent *inputEvent) override
The main processing event function.
Double_t fPointThreshold
The number of sigmas over baseline fluctuations to identify a point over threshold.
std::string fBaseLineOption
Option for calculation of baseline parameters, can be set to "ROBUST".
void Initialize() override
Function to initialize input/output event members and define the section name.
Double_t fSignalThreshold
A parameter to define a minimum signal fluctuation. Measured in sigmas.
Bool_t fRangeEnabled
Just a flag to quickly determine if we have to apply the range filter.
Int_t fPointsOverThreshold
The minimum number of points over threshold to identify a signal as such.
TRestRawSignalEvent * fSignalEvent
A pointer to the specific TRestRawSignalEvent input.
TVector2 fIntegralRange
The range where the observables will be calculated.
TVector2 fBaseLineRange
The range where the baseline range will be calculated.
void InitFromConfigFile() override
To make settings from rml file. This method must be implemented in the derived class.
An event container for time rawdata signals with fixed length.
void SetBaseLineRange(TVector2 blRange, std::string option="")
It defines a Short_t array with a physical parameter that evolves in time using a fixed time bin.
Double_t GetBaseLineSigma() const
Double_t GetThresholdIntegral()
It returns the integral of points identified over threshold. InitializePointsOverThreshold should hav...
Int_t GetID() const
Returns the value of signal ID.
Double_t GetMaxValue()
Returns the maximum value found in the data points. It includes baseline correction.
Double_t GetMaxPeakValue()
It returns the amplitude of the signal maximum, baseline will be corrected if CalculateBaseLine was c...
std::vector< Int_t > GetPointsOverThreshold() const
Returns a std::vector containing the indexes of data points over threshold.
Double_t GetBaseLine() const
void InitializePointsOverThreshold(const TVector2 &thrPar, Int_t nPointsOver, Int_t nPointsFlat=512)
It initializes the fPointsOverThreshold array with the indexes of data points that are found over thr...
Int_t GetRiseTime()
It returns the time required from the signal to reach the maximum. InitializePointsOverThreshold shou...
Bool_t IsADCSaturation(int Nflat=3)
It returns whether the signal has ADC saturation.
Int_t GetMaxPeakBin()
It returns the bin at which the maximum peak amplitude happens.
@ REST_Debug
+show the defined debug messages