GET library
GETSystemAnalyser Class Reference

#include <GETSystemAnalyser.hh>

Inheritance diagram for GETSystemAnalyser:
GETActarTpcAnalyser

Public Member Functions

Constructors, affectation, destructor
 GETSystemAnalyser (GETSystem *get=NULL)
 
virtual ~GETSystemAnalyser ()
 
General functions
virtual void Reset ()
 
virtual int SetSystem (GETSystem *gptr)
 
u_int GetChannelNumber () const
 
void SetEventsSkip (u_int n)
 
u_int GetEventsSkip ()
 
u_int GetAnalysisEventsCounter ()
 
void SetSelectedData (u_short id)
 
int GetSelectedData () const
 
Main analysis functions, called by GETSystem
virtual void ClearAnalysis ()
 
virtual int AnalyseRunStart ()
 
virtual int AnalyseRunStop ()
 
virtual void ClearEvent ()
 
virtual int AnalyseRawEvent ()
 
virtual int AnalyseCorEvent ()
 
virtual void BadChannelsWarning ()
 
virtual double GetMaximumSignal () const
 
virtual u_int GetMaximumChannel () const
 
Predefined analysis configuration files
virtual int ReadConfigFile (const string &cfg)
 
virtual int ReadConfigFile (FILE *fp)
 
virtual bool ConfigCommand (const GString &code, const GString &args="", FILE *fptr=NULL)
 
Full event data
TH2D * CreateFullEventRawHisto (const string &hname)
 
TH2D * CreateFullEventCorHisto (const string &hname)
 
TH2D * GetFullEventRawHisto () const
 
TH2D * GetFullEventCorHisto () const
 
Channels amplitude analysis
void SetAmplitudeMode (u_short mode, double win=0.)
 
void SetAmplitudeWindow (double win)
 
void SetTimingMode (u_short mode)
 
void SetTimingModeCFD (double frac, double delay)
 
virtual int CalcChannelOutput (u_int ic, double &vmax, double &tmax)
 
virtual int CalcChannelMaxima (u_int ic)
 
void SetMaximaScaling (double fact)
 
void InitAmplitudeData ()
 
void InitMaximaData ()
 
void ResetMaximaData ()
 
void ClearAmplitudeData ()
 
bool IsAmplitudeDataComputed () const
 
bool IsMaximaDataComputed () const
 
virtual void ComputeAmplitudeData ()
 
TH1D * CreateChannelsMaxRawHisto (const string &hname)
 
TH1D * CreateChannelsMaxCorHisto (const string &hname)
 
TH1D * GetChannelsMaxRawHisto () const
 
TH1D * GetChannelsMaxCorHisto () const
 
Analysis of number of data per channel
virtual void InitDataNumberCheck (u_int n0, u_int n1)
 
bool IsDataNumberCheckON () const
 
virtual int ProcessDataNumberCheck ()
 
virtual int GetDataNumberBadCount () const
 
Analysis of data continuity
virtual void InitDataContinuityCheck (double dmax)
 
bool IsDataContinuityCheckON () const
 
virtual bool ContinuityCheck (GETSample &sample)
 
virtual int ProcessDataContinuityCheck ()
 
virtual int GetDataContinuityBadCount () const
 
Noise analysis
virtual void SetNoiseRange (u_int i0=2, u_int i1=509)
 
u_int GetNoiseIndexMin () const
 
u_int GetNoiseIndexMax () const
 
TH1D * GetChannelsNoiseHisto () const
 
bool IsNoiseAnalysisON () const
 
virtual int ProcessNoiseAnalysis ()
 
virtual int FillCumulNoiseHisto (TH2 *hptr)
 
Amplitude resolution analysis
virtual void InitResolutionAnalysis (int n, double amin, double amax)
 
virtual void ClearResolutionAnalysis ()
 
bool IsResolutionAnalysisON () const
 
TH2I * GetResolutionMaxHisto () const
 
TH1D * GetResolutionHisto () const
 
virtual int ProcessResolutionAnalysis ()
 
Baseline fluctuations analysis
virtual void InitBaseFluctuations (u_int i0=2, u_int i1=32, double warn=0.)
 
virtual void ClearBaseFluctuations ()
 
TH1D * GetBaseFluctuationsHisto () const
 
bool IsFluctuationsAnalysisON () const
 
virtual int ProcessFluctuationsAnalysis ()
 
virtual int GetFluctuationsBadCount () const
 
virtual int FillCumulFluctuationsHisto (TH2 *hptr)
 
ROOT related functions
TH1D * CreateChannelsHisto1D (const string &name, const string &ytit="", bool fpn=true)
 
TH2D * CreateChannelsHisto2D (const string &name, int ny, double ymin, double ymax, const string &ytit="", bool fpn=true)
 
TH1I * CreateChannelsHisto1I (const string &name, const string &ytit="", bool fpn=true)
 
TH2I * CreateChannelsHisto2I (const string &name, int ny, double ymin, double ymax, const string &ytit="", bool fpn=true)
 
TH2D * CreateFullEventHisto (const string &name, bool fpn=true)
 
TH2D * CreateFullEventHisto (const string &name, int ny, double ymin, double ymax, bool fpn=true)
 
int FillFullEventHisto (TH2D *hptr, bool reset=true, bool fpn=true, double scale=1.)
 
TH1D * CreateChannelSummaryHisto (const string &name, bool fpn=true)
 
int FillChannelsSummaryHisto (TH1D *hptr, u_short mode=GET::dataMaxSignal, bool fpn=true, bool reset=true)
 
int CumulChannelsSummaryHisto (TH2 *hptr, u_short mode=GET::dataMaxSignal, bool fpn=true)
 
TH1D ** CreateAllChannelsHisto (const string &prefix="GET_", const string &suffix="")
 
int FillAllChannelsHisto (TH1D **htab, bool reset=true)
 
 ClassDef (GETSystemAnalyser, 0)
 

Protected Attributes

double * amax_amplitude
 Array of maximum amplitude values.
 
bool amax_computed
 Whether the maximum amplitude / time is computed.
 
double * amax_time
 Array of maximum amplitude times.
 
u_short amplitude_mode
 Signal amplitude analysis mode.
 
double amplitude_window
 Signal amplitude analysis window size (for fit or average)
 
u_int ch_max_index
 Channel with maximum amplitude.
 
double ch_max_signal
 Maximum channel amplitude value.
 
bool * channel_bad
 Whether the channel contains bad data.
 
bool channel_bad_flag
 Whether there are bad channels.
 
int data_cont_bad
 Number of bad channels in analysis of the data continuity.
 
double data_cont_max
 Maximum variation for continuity.
 
bool data_cont_on
 Whether the analysis of data discontinuity is ON.
 
int data_num_bad
 Number of bad channels in analysis of the data per channel.
 
int data_num_max
 Upper value for bad range of data / channel.
 
int data_num_min
 Lower value for bad range of data / channel.
 
bool data_num_on
 Whether the analysis of number of data per channel is ON.
 
u_int events_count
 Number of events since analysis cleard.
 
u_int events_skip
 Number of first events to skip for analysis.
 
int fluct_bad
 Number of bad channels in fluctuations analysis.
 
int * fluct_cnt
 Number of events for baseline check.
 
u_int fluct_imax
 Upper limit for baseline fluctuations analysis.
 
u_int fluct_imin
 Lower limit for baseline fluctuations analysis.
 
bool fluct_on
 Whether baseline fluctuations analysis is ON.
 
double * fluct_sum
 Cumulated values for baseline check.
 
double fluct_warn
 Baseline fluctuation value to issue a warning.
 
GETSystemget_ptr
 Pointer to the GET system object.
 
TH1I * h_channels_data_num
 Pointer to the channels summary histogram with number of raw data per channel.
 
TH1D * h_channels_fluct
 Pointer to the channels baseline fluctuations analysis summary average histogram.
 
TH1D * h_channels_max_cor
 Pointer to the channels summary histogram with maximum corrected output values.
 
TH1D * h_channels_max_raw
 Pointer to the channels summary histogram with maximum raw output values.
 
TH1D * h_channels_rms
 Pointer to the channels RMS analysis summary histogram.
 
TH2D * h_full_event_cor
 Pointer to the histogram with full data after corrections (if any)
 
TH2D * h_full_event_raw
 Pointer to the histogram with full raw data (before processing)
 
TH2I * h_resol_ampl
 Pointer to the channels versus amplitude maximum histogram.
 
TH1D * h_resol_rms
 Pointer to the channels amplitude variations histogram.
 
GETHitMaskhit_mask
 Pointer to the hit mask (same as base class hit mask, but known as GETHitMaskXY)
 
double ** mmax_amplitude
 Maxima values arrays for each channel.
 
bool mmax_analysis
 Whether multiple maxima analysis is requested.
 
bool mmax_computed
 Whether multiple maxima analysis is performed.
 
RExtremaFinder mmax_finder
 Multiple maxima analyser.
 
u_int mmax_max
 Maximum number of maxima per channel.
 
u_int * mmax_num
 Number of maxima for each channel.
 
double mmax_scale
 Scaling factor for amplitude.
 
double ** mmax_time
 Maxima times arrays for each channel.
 
int resol_dim
 Dimension of amplitude axis.
 
double resol_max
 Upper limit of amplitude axis.
 
double resol_min
 Lower limit of amplitude axis.
 
bool resol_on
 Whether amplitude resolution analysis is ON.
 
double resol_step
 Lower limit of amplitude axis.
 
u_int rms_imax
 Upper limit for RMS analysis.
 
u_int rms_imin
 Lower limit for RMS analysis.
 
bool rms_on
 Whether noise analysis is ON.
 
u_short selected_data
 Information on the channel data to analyse: GET::signalOut (default), GET::signalTst or GET::signalRec.
 
double timing_cfd_delay
 CFD offset for timing analysis.
 
double timing_cfd_fraction
 CFD fraction for timing analysis.
 
u_short timing_mode
 Signal timing analysis mode.
 

Private Member Functions

 GObject (GETSystemAnalyser)
 

Friends

class GETSystem
 

Detailed Description

This is a base class for various analysis of experimental data from the GETSystem objects. When a GETSystemAnalyser object is associated to a GETSystem object, the analysis functions from this class are called by the GETSystem specific functions:

In addition, the function GETSystemAnalyser::ClearAnalysis may be called when all analysis data should be resetted. Note that this function is never called by the GETSystem class.

The GETSystemAnalyser::ClearEvent function is automatically called by the GETSystem object when its clears its own event data.

This class may be derived in order to perform specific analysis, eventually from a derived class of GETSystem: see for example the GETActarDemAnalyser for the GETActarDem system class.

Note
When a GETSystemAnalyser (or derived) class object is used to perform events analysis while the events are not read using the GETSystem (or derived) object event processing functions, then the analyser functions (GETSystemAnalyser::AnalyseRunStart, GETSystemAnalyser::AnalyseRunStop, GETSystemAnalyser::AnalyseRawEvent and GETSystemAnalyser::AnalyseCorEvent) are not called automatically.
This for example the case when events are read from simulation and imported in a GETSystem object. In this case, the functions must be explicitely called for the analysis to be performed and the histograms to be filled properly.

The GETSystemAnalyser class can perform various types of analysis, depending on what has been defined in terms of analysis parameters and histograms.

The analysis available for the base class are:

In addition, it is possible to read the analysis parameters from a text file: see functions GETSystemAnalyser::ReadConfigFile and GETSystemAnalyser::ConfigCommand.

Constructor & Destructor Documentation

GETSystemAnalyser::~GETSystemAnalyser ( )
virtual

Destructor. The histograms associated to the analysis are deleted.

References channel_bad, and Reset().

Member Function Documentation

int GETSystemAnalyser::AnalyseCorEvent ( )
virtual

Function called when the GETSystem has read an event, after applying corrections to the output samples. The function does the following actions:

  • update the hit mask;
  • fill automatic full corrected event (all channels) histogram if defined;
  • compute channels maximum amplitude and time for corrected data;
  • fill automatic corrected data channels maximum amplitude histogram if defined;
  • perfom some event checks if defined:
    • resolution analysis (for pulser runs)
    • noise analysis

Reimplemented in GETActarTpcAnalyser.

References amax_computed, ComputeAmplitudeData(), GET::dataMaxSignal, events_count, events_skip, FillChannelsSummaryHisto(), FillFullEventHisto(), h_channels_max_cor, h_full_event_cor, hit_mask, GETHitMask::IsMaskSet(), ProcessNoiseAnalysis(), ProcessResolutionAnalysis(), resol_on, rms_on, GET::signalOut, and GETHitMask::UpdateMask().

Referenced by GETSystem::AnalyseCorEvent().

int GETSystemAnalyser::AnalyseRawEvent ( )
virtual

Function called when the GETSystem has read an event, before applying corrections to the output samples. The function does the following actions:

  • fill automatic full raw event (all channels) histogram if defined;
  • compute channels maximum amplitude and time for raw data;
  • fill automatic raw data channels maximum amplitude histogram if defined;
  • perfom some event checks if defined:
    • number of data per channel,
    • signal (dis)continuity check,
    • fluctuations analysis

Reimplemented in GETActarTpcAnalyser.

References amax_computed, ComputeAmplitudeData(), data_cont_on, data_num_on, GET::dataMaxSignal, events_count, events_skip, FillChannelsSummaryHisto(), FillFullEventHisto(), fluct_on, h_channels_max_raw, h_full_event_raw, mmax_analysis, ProcessDataContinuityCheck(), ProcessDataNumberCheck(), ProcessFluctuationsAnalysis(), and GET::signalOut.

Referenced by GETActarTpcAnalyser::AnalyseRawEvent(), and GETSystem::AnalyseRawEvent().

int GETSystemAnalyser::AnalyseRunStart ( )
virtual

Function called when the GETSystem starts processing a new file. Note that the ClearAnalysis function is not automatically called in this function.

Referenced by GETSystem::OpenRunFile(), and GETSystem::OpenRunSerie().

int GETSystemAnalyser::AnalyseRunStop ( )
virtual

Function called when the GETSystem stops processing a run file.

Referenced by GETSystem::CloseRunFile().

void GETSystemAnalyser::BadChannelsWarning ( )
virtual

Function sending a detailed warning message indicating where bad data have been detected in event analysis.

References channel_bad, channel_bad_flag, get_ptr, GETCoBo::GetChildrenNumber(), GETAsAd::GetChildrenNumber(), GETAGet::GetChildrenNumber(), and GETSystem::GetChildrenNumber().

Referenced by GETActarEventReader::AnalysisEvent().

int GETSystemAnalyser::CalcChannelMaxima ( u_int  ic)
virtual

Function performing the maxima search for a channel. Only the maxima above the channel threshold are considered. The function returns 0 if no error occured.

The amplitude of maxima corresponds to the integral amplitude of the signal around the maxima. The scaling factor GETSystemAnalyser::mmax_scale can be used to set a scaling factor (for example to convert from integral to signal maximum value, if the response function is known).

Parameters
icindex of channel to analyse

References get_ptr, GETSystem::GetChannel(), GETHitMask::GetThreshold(), hit_mask, GETHitMask::IsValidHit(), mmax_amplitude, mmax_finder, mmax_max, mmax_num, mmax_scale, mmax_time, GETChannel::OutSample(), GETChannel::RecSample(), selected_data, GET::signalRec, GET::signalTst, and GETChannel::TstSample().

Referenced by ComputeAmplitudeData().

int GETSystemAnalyser::CalcChannelOutput ( u_int  ic,
double &  vmax,
double &  tmax 
)
virtual

Function computing the amplitude and the timing for a channel. They are estimated according to the calculation modes defined for the analyser (see SetAmplitudeMode function).

Note that the calculations that requires to analyse data are used only:

The function returns 0 if no error occured.

Parameters
icindex of channel to analyse
vmaxreturned value for the amplitude
tmaxreturned value for the timing

References amplitude_mode, amplitude_window, GET::amplitudeLocalAvg, GET::amplitudeLocalFitP2, GET::amplitudeMaxSignal, get_ptr, GETSystem::GetChannel(), GETSystem::GetHitMask(), GETObject::GetSampleDim(), GETObject::GetTimeStep(), GETTimingCFD(), hit_mask, GETChannel::IsFPNChannel(), GETHitMask::IsValidHit(), GETChannel::OutSample(), GETChannel::RecSample(), selected_data, GET::signalRec, GET::signalTst, timing_cfd_delay, timing_cfd_fraction, timing_mode, GET::timingCFD, GET::timingMaxSignal, and GETChannel::TstSample().

Referenced by ComputeAmplitudeData().

GETSystemAnalyser::ClassDef ( GETSystemAnalyser  ,
 
)

for use within ROOT.

void GETSystemAnalyser::ClearAmplitudeData ( )

Function that sets the maximum amplitudes and times to 0. For multiple maxima analysis, the number of maxima is set to 0 for each channel.

References amax_amplitude, amax_computed, amax_time, ch_max_index, ch_max_signal, GetChannelNumber(), mmax_computed, and mmax_num.

Referenced by ClearEvent(), ComputeAmplitudeData(), InitAmplitudeData(), and SetSelectedData().

void GETSystemAnalyser::ClearAnalysis ( )
virtual

Function that clears all current analysis results: histograms, baseline fluctuations,...

Reimplemented in GETActarTpcAnalyser.

References ClearBaseFluctuations(), ClearEvent(), ClearResolutionAnalysis(), events_count, fluct_on, and resol_on.

Referenced by GETActarTpcAnalyser::ClearAnalysis().

void GETSystemAnalyser::ClearBaseFluctuations ( )
virtual

Function that clears the baseline fluctuations analysis. The sum and event counts per channel are set to 0.

References fluct_bad, fluct_cnt, fluct_sum, get_ptr, GETSystem::GetChannelCount(), and h_channels_fluct.

Referenced by ClearAnalysis(), and InitBaseFluctuations().

void GETSystemAnalyser::ClearEvent ( )
virtual
void GETSystemAnalyser::ClearResolutionAnalysis ( )
virtual

Function that clears the resolution analysis. The corresponding histograms are resetted.

References h_resol_ampl, and h_resol_rms.

Referenced by ClearAnalysis().

bool GETSystemAnalyser::ConfigCommand ( const GString &  code,
const GString &  args = "",
FILE *  fptr = NULL 
)
virtual

Function that sets a configuration parameter from a command line.

The valid command line in the class are:

  • analysis configuration:
    • INCL <cfg_file>: include another configuration file
    • THR <val>: set global channels (software) threshold
    • MMAX <scale>: perform multiple maxima analysis
    • AMPL <mode> [info]: set channel amplitude calculation mode:
      • mode = MAX: maximum signal amplitude
      • mode = AVG: average signal amplitude, the additional info is the time window for averaging
      • mode = P2: quadratic fit around maximum, the additional info is the time window for fitting
    • TIME <mode> [info]: set channel timing calculation mode
      • mode = MAX: maximum signal time
      • mode = CFD: use of a CFD algorithm, the additional information is the fraction and the delay for the CFD (positive values)
    • SIGNAL <info>: select the analysed signal: [ TST | OUT | REC ] (default is OUT)
  • output signals correction processes:
    • FPN <val>: set FPN corrections (1, 2, 4 or 0 for none)
    • SMOOTH <width>: add a smoothing (Gaussian) filter
    • CAL <file>: add calibration from file
    • BL_CORR <baseline_file>: set baseline correction from file
    • AUTO_BL <nb_lo> <off_lo> <nb_hi> <off_hi>: set linear (low & high time buckets) automatic baseline correction
    • AUTO_BL_LOW <nb_lo> <off_lo> : set constant (low time buckets) automatic baseline correction
    • AUTO_BL_HIGH <nb_hi> <off_hi> : set constant (high time buckets) automatic baseline correction
  • predefined analysis:
    • NOISE <i0> <i1>: set the noise (RMS) analysis between time buckets i0 and i1 (included)
    • BASE_CHK <i0> <i1> <limit>: set baseline fluctuation analysis between time buckets i0 and i1 (included), that should stay below the limit
    • READ_CHK <i0> <i1>: check the number of data per channel, that should be outside the range [i0,i1]
    • DELTA_CHK <diff>: check the continuity between time buckets, by setting a maximum difference between contiguous buckets
Parameters
codeparameter(s) identifier
argsadditional arguments
fptrfile pointer (in case more information should be read)

Reimplemented in GETActarTpcAnalyser.

References GET::amplitudeLocalAvg, GET::amplitudeLocalFitP2, GET::amplitudeMaxSignal, GET::CoBoCorrectFPN, GET::CoBoCorrectFPN1, GET::CoBoCorrectFPN2, GET::CoBoCorrectFPN4, get_ptr, GETSystem::GetChannelCount(), GETSystem::GetHitMask(), GETObject::GetSampleDim(), GetSelectedData(), GETObject::GetTimeStep(), InitAmplitudeData(), InitBaseFluctuations(), InitDataContinuityCheck(), InitDataNumberCheck(), InitResolutionAnalysis(), mmax_analysis, ReadConfigFile(), SetAmplitudeMode(), GETSystem::SetAutoBaseline(), GETSystem::SetAutoBaselineHigh(), GETSystem::SetAutoBaselineLow(), GETSystem::SetBaselineCorrection(), GETSystem::SetCalibProcess(), SetEventsSkip(), SetNoiseRange(), GETSystem::SetOption(), GETSystem::SetOutputSmoothing(), SetSelectedData(), GETHitMask::SetThreshold(), SetTimingMode(), SetTimingModeCFD(), GET::signalOut, GET::signalRec, GET::signalTst, GET::timingCFD, and GET::timingMaxSignal.

Referenced by GETActarTpcAnalyser::ConfigCommand(), and ReadConfigFile().

bool GETSystemAnalyser::ContinuityCheck ( GETSample sample)
virtual

Function returning true if no discontinuity is found in the data. A discontinuity is identified as data variation higher than the maximaum accepted value compared to both the next and the previous data.

Parameters
sampledata sample

References data_cont_max, and data_cont_on.

Referenced by ProcessDataContinuityCheck().

TH1D ** GETSystemAnalyser::CreateAllChannelsHisto ( const string &  prefix = "GET_",
const string &  suffix = "" 
)

Function allocating an array of histogram pointers and the corresponding histograms for each channel.

Each histogram in the array and the array itself must be deleted by the calling program after use.

Parameters
prefixprefix for channels histograms name
suffixsuffix for channels histograms name

References get_ptr, GETSystem::GetChannel(), GETSystem::GetChannelCount(), GETObject::GetFullId(), and GETChannel::OutSample().

TH1D * GETSystemAnalyser::CreateChannelsHisto1D ( const string &  name,
const string &  ytit = "",
bool  fpn = true 
)

Create a 1D histogram with channels indexes as X-axis.

Parameters
namehistogram name
ytittitle for the Y-axis
fpnwhether the FPN channels are included

References get_ptr, GETSystem::GetChannelCount(), and GETObject::GetSignalChannelCount().

Referenced by CreateChannelSummaryHisto(), and InitResolutionAnalysis().

TH1I * GETSystemAnalyser::CreateChannelsHisto1I ( const string &  name,
const string &  ytit = "",
bool  fpn = true 
)

Create a 1D histogram with channels indexes as X-axis.

Parameters
namehistogram name
ytittitle for the Y-axis
fpnwhether the FPN channels are included

References get_ptr, GETSystem::GetChannelCount(), and GETObject::GetSignalChannelCount().

TH2D * GETSystemAnalyser::CreateChannelsHisto2D ( const string &  name,
int  ny,
double  ymin,
double  ymax,
const string &  ytit = "",
bool  fpn = true 
)

Create a 2D histogram with channels indexes as X-axis.

Parameters
namename of the histogram
nynumber of bins for Y dimension
yminlower limit Y axis
ymaxupper limit Y axis
ytittitle for the Y-axis
fpnwhether the FPN channels are included

References get_ptr, GETSystem::GetChannelCount(), and GETObject::GetSignalChannelCount().

Referenced by CreateFullEventHisto().

TH2I * GETSystemAnalyser::CreateChannelsHisto2I ( const string &  name,
int  ny,
double  ymin,
double  ymax,
const string &  ytit = "",
bool  fpn = true 
)

Create a 2D histogram with channels indexes as X-axis.

Parameters
namename of the histogram
nynumber of bins for Y dimension
yminlower limit Y axis
ymaxupper limit Y axis
ytittitle for the Y-axis
fpnwhether the FPN channels are included

References get_ptr, GETSystem::GetChannelCount(), and GETObject::GetSignalChannelCount().

Referenced by InitResolutionAnalysis().

TH1D * GETSystemAnalyser::CreateChannelsMaxCorHisto ( const string &  hname)

Function that defines the histogram for channel maximum amplitude of corrected output data. This histogram is filled automatically by the AnalyseRawEvent function.

Parameters
hnamename for the histogram

References CreateChannelSummaryHisto(), and h_channels_max_cor.

TH1D * GETSystemAnalyser::CreateChannelsMaxRawHisto ( const string &  hname)

Function that defines the histogram for channel maximum amplitude of uncorrected output data. This histogram is filled automatically by the AnalyseRawEvent function.

Parameters
hnamename for the histogram

References CreateChannelSummaryHisto(), and h_channels_max_raw.

TH1D * GETSystemAnalyser::CreateChannelSummaryHisto ( const string &  name,
bool  fpn = true 
)
inline

Create an histogram to display a 1D channel summary, where the X dimension is the channel number.

Parameters
namename of the histogram
fpnwhether the FPN channels are included

References CreateChannelsHisto1D().

Referenced by CreateChannelsMaxCorHisto(), and CreateChannelsMaxRawHisto().

TH2D * GETSystemAnalyser::CreateFullEventCorHisto ( const string &  hname)

Function that defines the histogram for full event corrected output data. This histogram is filled automatically by the AnalyseRawEvent function.

Parameters
hnamename for the histogram

References CreateFullEventHisto(), and h_full_event_cor.

Referenced by GETActarEventReader::PresetAnalysis().

TH2D * GETSystemAnalyser::CreateFullEventHisto ( const string &  name,
bool  fpn = true 
)

Create an histogram to display a full event. The X dimension is the channel number and the Y dimension is the channel content (time axis of channels).

The histogram must be deleted by calling program after use.

The Y dimension is set to the ADC range.

Parameters
namename of the histogram
fpnwhether the FPN channels are included

References get_ptr, GETObject::GetTotalChannel(), and GETChannel::OutSample().

Referenced by CreateFullEventCorHisto(), and CreateFullEventRawHisto().

TH2D * GETSystemAnalyser::CreateFullEventHisto ( const string &  name,
int  ny,
double  ymin,
double  ymax,
bool  fpn = true 
)
inline

Create an histogram to display a full event. The X dimension is the channel number and the Y dimension is the channel content. The Y range is set argument values.

The histogram must be deleted by calling program after use.

Parameters
namename of the histogram
nynumber of bins for Y dimension
yminlower limit Y axis
ymaxupper limit Y axis
fpnwhether the FPN channels are included

References CreateChannelsHisto2D().

TH2D * GETSystemAnalyser::CreateFullEventRawHisto ( const string &  hname)

Function that defines the histogram for full event uncorrected output data. This histogram is filled automatically by the AnalyseRawEvent function.

Parameters
hnamename for the histogram

References CreateFullEventHisto(), and h_full_event_raw.

Referenced by GETActarEventReader::PresetAnalysis().

int GETSystemAnalyser::CumulChannelsSummaryHisto ( TH2 *  hptr,
u_short  mode = GET::dataMaxSignal,
bool  fpn = true 
)

Function filling a 2D histogram (i.e. from CreateFullEventHisto function) with data from current event channels (output, reconstructed or test sample, according to selected data). One value per channel is incremented depending on mode argument. The histogram is not resetted before filling.

Parameters
hptrhistogram pointer
modedata selected for filling: GET::dataMaxSignal (default), GET::dataIntegral or GET::dataTime
fpnwhether the FPN channels are included

References amax_amplitude, amax_time, ComputeAmplitudeData(), GET::dataIntegral, GET::dataTime, get_ptr, GETSystem::GetChannel(), GETSystem::GetChannelCount(), GETChannel::GetReadoutFlag(), GETObject::GetSignalChannelCount(), IsAmplitudeDataComputed(), GETChannel::IsSignalChannel(), GETChannel::OutSample(), GETChannel::RecSample(), selected_data, GET::signalRec, GET::signalTst, and GETChannel::TstSample().

Referenced by ProcessResolutionAnalysis().

int GETSystemAnalyser::FillAllChannelsHisto ( TH1D **  htab,
bool  reset = true 
)

Function filling a set of single histogram with channels signal (output, reconstructed or test sample, according to selected data). The histogram pointers must have been stored in the htab array (i.e. from CreateAllChannelsHisto function). The histograms are resetted before filling. The function returns 0 if no error occured, a positive value in case of error, and a negative value if some histograms were not filled.

Parameters
htabarray of histogram pointers
resetwhether the histograms are resetted before filling

References get_ptr, GETSystem::GetChannel(), GETSystem::GetChannelCount(), GETObject::GetSampleDim(), GETChannel::OutSample(), GETChannel::RecSample(), selected_data, GET::signalRec, GET::signalTst, and GETChannel::TstSample().

int GETSystemAnalyser::FillChannelsSummaryHisto ( TH1D *  hptr,
u_short  mode = GET::dataMaxSignal,
bool  reset = true,
bool  fpn = true 
)

Function filling a 2D histogram (i.e. from CreateFullEventHisto function) with data from current event channels (output, reconstructed or test sample, according to selected data). One value per channel is incremented depending on mode argument. The histogram is not resetted before filling.

Parameters
hptrhistogram pointer
modedata selected for filling: GET::dataMaxSignal (default), GET::dataIntegral or GET::dataTime
fpnwhether the FPN channels are included
resetwhether the histogram is resetted before filling

References amax_amplitude, amax_time, ComputeAmplitudeData(), GET::dataIntegral, GET::dataTime, get_ptr, GETSystem::GetChannel(), GETSystem::GetChannelCount(), GETObject::GetSignalChannelCount(), IsAmplitudeDataComputed(), GETChannel::IsSignalChannel(), GETChannel::OutSample(), GETChannel::RecSample(), selected_data, GET::signalRec, GET::signalTst, and GETChannel::TstSample().

Referenced by AnalyseCorEvent(), and AnalyseRawEvent().

int GETSystemAnalyser::FillCumulFluctuationsHisto ( TH2 *  hptr)
virtual

Fills the argument 2D histogram with fluctuations analysis result.

Parameters
hptrhistogram pointer (2D)

References fluct_on, GetChannelNumber(), and h_channels_fluct.

int GETSystemAnalyser::FillCumulNoiseHisto ( TH2 *  hptr)
virtual

Fills the argument 2D histogram with noise analysis result.

Parameters
hptrhistogram pointer (2D)

References GetChannelNumber(), h_channels_rms, and rms_on.

int GETSystemAnalyser::FillFullEventHisto ( TH2D *  hptr,
bool  reset = true,
bool  fpn = true,
double  scale = 1. 
)

Function filling a 2D histogram (i.e. from CreateFullEventHisto function) with current event (output, reconstructed or test sample, according to selected data). The histogram is resetted before filling.

Parameters
hptrhistogram pointer
resetwhether the histograms are resetted before filling
fpnwhether the FPN channels are included
scalescaling of signal before filling (this may be used if the histogram binning is different than the data samples)

References get_ptr, GETSystem::GetChannel(), GETSystem::GetChannelCount(), GETObject::GetSignalChannelCount(), GETObject::GetTimeStep(), GETChannel::IsSignalChannel(), GETChannel::OutSample(), GETChannel::RecSample(), selected_data, GET::signalRec, GET::signalTst, and GETChannel::TstSample().

Referenced by AnalyseCorEvent(), and AnalyseRawEvent().

u_int GETSystemAnalyser::GetAnalysisEventsCounter ( )
inline

Return the number of events since the last analysis clear.

References events_count.

TH1D * GETSystemAnalyser::GetBaseFluctuationsHisto ( ) const
inline

Return the pointer to the corresponding histogram.

References h_channels_fluct.

Referenced by GETActarEventReader::AnalysisEvent().

TH1D * GETSystemAnalyser::GetChannelsMaxCorHisto ( ) const
inline

Return the pointer to the corresponding histogram (corrected output).

References h_channels_max_cor.

TH1D * GETSystemAnalyser::GetChannelsMaxRawHisto ( ) const
inline

Return the pointer to the corresponding histogram (raw output).

References h_channels_max_raw.

TH1D * GETSystemAnalyser::GetChannelsNoiseHisto ( ) const
inline

Return the pointer to the channels noise histogram of current event.

References h_channels_rms.

Referenced by GETActarEventReader::AnalysisEvent().

int GETSystemAnalyser::GetDataContinuityBadCount ( ) const
inlinevirtual

Return the number of channels for which data continuity is bad.

References data_cont_bad.

Referenced by GETActarEventReader::AnalysisEvent().

int GETSystemAnalyser::GetDataNumberBadCount ( ) const
inlinevirtual

Return the number of channels for which number of data is in the rejection range.

References data_num_bad.

Referenced by GETActarEventReader::AnalysisEvent().

u_int GETSystemAnalyser::GetEventsSkip ( )
inline

Return the number of events to skip at the beginning of analysis.

References events_skip.

int GETSystemAnalyser::GetFluctuationsBadCount ( ) const
inlinevirtual

Return the number of channels for which the baseline difference with the average is over the warning threshold.

References fluct_bad.

Referenced by GETActarEventReader::AnalysisEvent().

TH2D * GETSystemAnalyser::GetFullEventCorHisto ( ) const
inline

Return the pointer to the corresponding histogram (corrected output).

References h_full_event_raw.

TH2D * GETSystemAnalyser::GetFullEventRawHisto ( ) const
inline

Return the pointer to the corresponding histogram (raw output).

References h_full_event_raw.

u_int GETSystemAnalyser::GetMaximumChannel ( ) const
inlinevirtual

Return the index of the signal channel with maximum amplitude value for the selected data (input/test or reconstructed). This value is stored when the event is processed by the analyser, which make this function faster than the GETSystem::GetMaximumChannel function.

References ch_max_index.

Referenced by GETActarEventReader::AnalysisEvent().

double GETSystemAnalyser::GetMaximumSignal ( ) const
inlinevirtual

Return the maximum channel amplitude value for the selected data (input/test or reconstructed). This value is stored when the event is processed by the analyser, which make this function faster than the GETSystem::GetMaximumChannel function.

References ch_max_signal.

Referenced by GETActarEventReader::AnalysisEvent().

u_int GETSystemAnalyser::GetNoiseIndexMax ( ) const
inline

Return the upper time bucket index of noise analysis. (the index is included in the analysis).

References rms_imax.

u_int GETSystemAnalyser::GetNoiseIndexMin ( ) const
inline

Return the lower time bucket index of noise analysis (the index is included in the analysis).

References rms_imin.

TH1D * GETSystemAnalyser::GetResolutionHisto ( ) const
inline

Return the pointer to the channels resolution histogram.

References h_resol_rms.

Referenced by GETActarEventReader::AnalysisEvent().

TH2I * GETSystemAnalyser::GetResolutionMaxHisto ( ) const
inline

Return the pointer to the channel versus maximum amplitude histogram.

References h_resol_ampl.

int GETSystemAnalyser::GetSelectedData ( ) const
inline

Return the analysed data: output (default), input/test or reconstructed signals.

References selected_data.

Referenced by GETActarEventReader::AnalysisEvent(), ConfigCommand(), and GETActarEventReader::DefinePadSpectrum().

GETSystemAnalyser::GObject ( GETSystemAnalyser  )
private

Macro from GCpp library that defines the following functions:

  • ClassName(): return the real class name of the object
  • StaticClassName(): return the used class name of the object (that may be a base class of the real object).
  • IsInstanceOf<T>(): return true if the current object is an instance of the template class name argument
  • Clone(): return an allocated copy of the object.
void GETSystemAnalyser::InitAmplitudeData ( )

Function allocating the arrays for maximum amplitude and time estimates for all channels.

References amax_amplitude, amax_time, ClearAmplitudeData(), GetChannelNumber(), InitMaximaData(), mmax_analysis, mmax_finder, and ResetMaximaData().

Referenced by ConfigCommand(), GETActarTpcAnalyser::SetSystem(), and SetSystem().

void GETSystemAnalyser::InitBaseFluctuations ( u_int  i0 = 2,
u_int  i1 = 32,
double  warn = 0. 
)
virtual

Function initializing the baseline fluctuation analysis. If the histograms are not defined, they are allocated, else only the analysis condition are resetted.

Parameters
i0lower time bucket index of baseline analysis range
i1upper time bucket index of baseline analysis range
warnfluctuation level threshold to issue a warning

References ClearBaseFluctuations(), fluct_bad, fluct_cnt, fluct_imax, fluct_imin, fluct_on, fluct_sum, fluct_warn, get_ptr, GETSystem::GetChannelCount(), GETObject::GetSampleDim(), and h_channels_fluct.

Referenced by ConfigCommand().

void GETSystemAnalyser::InitDataContinuityCheck ( double  dmax = 0.)
virtual

Set the analysis of the data continuity. Set the argument to 0 to switch the analysis OFF.

Parameters
dmaxMaximum in a channel sample variation for continuity

References data_cont_bad, data_cont_max, data_cont_on, and get_ptr.

Referenced by ConfigCommand().

void GETSystemAnalyser::InitDataNumberCheck ( u_int  n0,
u_int  n1 
)
virtual

Set the analysis of the number of data read per channel. If the number of data read is between n1 and n2 (included), a warning is issued.

For example, in full readout mode, without zero suppression, the limits should be 0 and 511 since for each channel, 512 data may be read.

In partial readout mode, either 0 or 512 data shoulde be read for each channel, and the limits 1 and 511 may be defined

Set both limits to 0 to turn off the analysis

Parameters
n0lower limit of bad number of read data per channel
n1upper limit of bad number of read data per channel

References data_num_bad, data_num_max, data_num_min, data_num_on, and get_ptr.

Referenced by ConfigCommand().

void GETSystemAnalyser::InitMaximaData ( )

Function allocating the arrays for multiple maxima estimates for all channels.

References get_ptr, GetChannelNumber(), GETObject::GetSampleDim(), mmax_amplitude, mmax_max, mmax_num, mmax_time, and ResetMaximaData().

Referenced by InitAmplitudeData().

void GETSystemAnalyser::InitResolutionAnalysis ( int  n,
double  amin,
double  amax 
)
virtual

Function that defines the maximum amplitude resolution analysis.

A 2D histogram is created to store maximum amplitudes versus the channel number. If this histogram existed before, it is deleted and rect=reated in order to adapt the new dimensions

A 1D histogram is created (if not already existing) to get the resolution for each channel. This histogram is not recreated, since the channel number may not change.

Parameters
namplitude axis dimension
aminlower limit of amplitude axis
amaxupper limit of amplitude axis

References CreateChannelsHisto1D(), CreateChannelsHisto2I(), get_ptr, h_resol_ampl, h_resol_rms, resol_dim, resol_max, resol_min, resol_on, and resol_step.

Referenced by ConfigCommand().

bool GETSystemAnalyser::IsDataContinuityCheckON ( ) const
inline

Return true if the check of data continuity is set.

References data_cont_on.

bool GETSystemAnalyser::IsDataNumberCheckON ( ) const
inline

Return true if the check of data per channel is set.

References data_num_on.

bool GETSystemAnalyser::IsFluctuationsAnalysisON ( ) const
inline

Return true if the baseline fluctuations analysis is defined.

References fluct_on.

bool GETSystemAnalyser::IsMaximaDataComputed ( ) const
inline

Function that returns true is the channels maxima have been computed.

References mmax_computed.

bool GETSystemAnalyser::IsNoiseAnalysisON ( ) const
inline

Return true if the noise analysis is defined.

References rms_on.

Referenced by GETActarEventReader::PresetAnalysis().

bool GETSystemAnalyser::IsResolutionAnalysisON ( ) const
inline

Return true if the maximum amplitude resolution analysis is defined.

References resol_on.

int GETSystemAnalyser::ProcessDataContinuityCheck ( )
virtual

Perform the analysis of data continuity check. The function returns the number of bad channels, for which the number some variations passes the threshold.

References channel_bad, channel_bad_flag, ContinuityCheck(), data_cont_bad, data_cont_on, data_num_bad, get_ptr, GETSystem::GetChannel(), GetChannelNumber(), GETSystem::GetEventCount(), and GETChannel::OutSample().

Referenced by AnalyseRawEvent().

int GETSystemAnalyser::ProcessDataNumberCheck ( )
virtual

Perform the analysis of number of data per channel. The function returns the number of bad channels, for which the number of data is between the preset limits (included).

References channel_bad, channel_bad_flag, data_num_bad, data_num_max, data_num_min, data_num_on, get_ptr, GETSystem::GetChannel(), GetChannelNumber(), GETSystem::GetEventCount(), and GETChannel::GetReadoutCount().

Referenced by AnalyseRawEvent().

int GETSystemAnalyser::ProcessFluctuationsAnalysis ( )
virtual

Perform the baseline fluctuations analysis: the mean sample value in the analysis time range is compared to the average of previous events. The function returns the number of bad channels.

References channel_bad, channel_bad_flag, GETSystem::ChannelHasData(), fluct_bad, fluct_cnt, fluct_imax, fluct_imin, fluct_on, fluct_sum, fluct_warn, get_ptr, GETSystem::GetChannel(), GetChannelNumber(), GETSystem::GetEventCount(), h_channels_fluct, and GETChannel::OutSample().

Referenced by AnalyseRawEvent().

int GETSystemAnalyser::ProcessNoiseAnalysis ( )
virtual

Compute the noise (RMS) for all channels int the time range for noise analysis. The results are stored in the corresponding histogram.

Reimplemented in GETActarTpcAnalyser.

References GETSystem::ChannelHasData(), get_ptr, GETSystem::GetChannel(), GetChannelNumber(), h_channels_rms, GETChannel::OutSample(), rms_imax, rms_imin, and rms_on.

Referenced by AnalyseCorEvent(), and GETActarTpcAnalyser::ProcessNoiseAnalysis().

int GETSystemAnalyser::ProcessResolutionAnalysis ( )
virtual

Function processing current event in the maximum amplitude resolution analysis. The maximum amplitude must have been previously computed (for example this is done by GETSystemAnalyser::ComputeAmplitudeData, that is called in the GETSystemAnalyser::AnalyseCorEvent function).

The function first add the current event to the 2D amplitude summary, then it computes the RMS for all channels, from the 2D histogram analysis, and stores the result in the resolution histogram.

References CumulChannelsSummaryHisto(), GET::dataMaxSignal, get_ptr, GETSystem::GetChannelCount(), h_resol_ampl, h_resol_rms, resol_dim, resol_on, and resol_step.

Referenced by GETActarTpcAnalyser::AnalyseCorEvent(), and AnalyseCorEvent().

int GETSystemAnalyser::ReadConfigFile ( const string &  cfg)
virtual

Function that reads setup parameters for GET system.

The function returns 0 if no error occured, a positive value in case of file error and a negative value if some invalid commands are found.

Parameters
cfgconfiguration file name

Referenced by ConfigCommand().

int GETSystemAnalyser::ReadConfigFile ( FILE *  fp)
virtual

Function that reads setup parameters for GET system.

The function returns 0 if no error occured, a positive value in case of file error and a negative value if some invalid commands are found.

Parameters
fpsetup file pointer

References ConfigCommand().

void GETSystemAnalyser::ResetMaximaData ( )

Function deleting the arrays for multiple maxima estimates for all channels.

References mmax_amplitude, mmax_computed, mmax_num, and mmax_time.

Referenced by InitAmplitudeData(), and InitMaximaData().

void GETSystemAnalyser::SetAmplitudeMode ( u_short  mode,
double  win = 0. 
)

This function defines the procedure for the calculation of the channels amplitude analysis. The valid modes are:

References amplitude_mode, amplitude_window, GET::amplitudeLocalFitP2, GET::amplitudeMaxSignal, and SetAmplitudeWindow().

Referenced by ConfigCommand().

void GETSystemAnalyser::SetAmplitudeWindow ( double  win)

This function amplitude analysis window around the maximum for local maximum analysis modes.

Parameters
winassociated analysis time window for local analysis

References amplitude_window.

Referenced by SetAmplitudeMode().

void GETSystemAnalyser::SetEventsSkip ( u_int  n)
inline

Set the number of events to skip at the beginning of analysis: they are not analysed at all. An call the the GETSystemAnalyser::ClearAnalysis function resets the analysis event counter, and next n events will be skipped again.

Parameters
nnumber of events to skip

References events_skip.

Referenced by ConfigCommand().

void GETSystemAnalyser::SetMaximaScaling ( double  fact)
inline

Function that sets the scaling factor for maxima analysis. The maxima are basically the integral of signal around the maxima. For example, for Gauss peaks, a factor 1/(sqrt(2*pi)*sigma) would then give amplitudes instead of integral signal for multiple maxima analysis.

Parameters
factscaling factor

References mmax_scale.

void GETSystemAnalyser::SetNoiseRange ( u_int  i0 = 2,
u_int  i1 = 509 
)
virtual

Defines or modify the limits (time bucket indexes) for the noise analysis. If not yet defined, the corresponding histogram is created and the noise analysis is turned ON.

If both limits are set to 0, the analysis is turned OFF (but the histogram is not deleted).

Parameters
i0lower time bucket for the analysis range
i1upper time bucket for the analysis range

Reimplemented in GETActarTpcAnalyser.

References get_ptr, GETSystem::GetChannelCount(), GETObject::GetSampleDim(), h_channels_rms, rms_imax, rms_imin, and rms_on.

Referenced by ConfigCommand(), and GETActarTpcAnalyser::SetNoiseRange().

void GETSystemAnalyser::SetSelectedData ( u_short  id)
inline

Define the data to be analysed: output (default), input/test or reconstructed signals.

Parameters
idselected data (GET::signalOut, GET::signalTst or GET::signalRec)

References ClearAmplitudeData(), hit_mask, GETObject::id, selected_data, and GETHitMask::SetDataInfo().

Referenced by ConfigCommand().

int GETSystemAnalyser::SetSystem ( GETSystem gptr)
virtual

Define the GET system from which the data will be analysed. Note that all previously defined information is resetted. The GET system must be initialized (number of channels defined).

Parameters
gptrPointer to the GET system

Reimplemented in GETActarTpcAnalyser.

References channel_bad, channel_bad_flag, get_ptr, GETSystem::GetHitMask(), GETObject::GetTotalChannelCount(), hit_mask, InitAmplitudeData(), and Reset().

Referenced by GETSystem::SetAnalyser(), and GETActarTpcAnalyser::SetSystem().

void GETSystemAnalyser::SetTimingMode ( u_short  mode)

This function defines the procedure for the calculation of the channels timing analysis. The valid modes are:

References timing_mode, GET::timingCFD, and GET::timingMaxSignal.

Referenced by ConfigCommand().

void GETSystemAnalyser::SetTimingModeCFD ( double  frac,
double  delay 
)

This function selects the CFD procedure for the calculation of the channels timing analysis and sets the associated parameters.

Parameters
fracconstant fraction value
delayfraction signal delay

References timing_cfd_delay, timing_cfd_fraction, timing_mode, and GET::timingCFD.

Referenced by ConfigCommand().


The documentation for this class was generated from the following files: