wavecatcher-analysis
Loading...
Searching...
No Matches
ReadSampic Class Reference

Class containing SAMPIC data loader and Event Builder. More...

#include <ReadSampic.h>

Inheritance diagram for ReadSampic:
Collaboration diagram for ReadSampic:

Classes

struct  HitInfoReduced
 Stores additional information with the waveforms. More...
 

Public Member Functions

void EventBuilder (double, vector< float >={}, vector< int > channels={}, unsigned int=1, bool=false)
 Construct events based on coincident amplitude thresholds within a time window for different channels.
It currently constructs the events from the first waveform above the set amplitude threshold.
ToDo:
.
 
int GetCurrentChannel (int) override
 Get the current channel index for a certain waveform index.
 
int GetCurrentEvent (int) override
 Get the current event index for a certain waveform index.
 
int GetWaveformIndex (int, int) override
 Returns index of a certain event number (if data files are read in parallel threads)
 
TH1FGetwf (int, int, int=1) override
 Helper that returns the waveform histogram for a certain channel number and a certain event number.
 
void PlotWF (int, float=0, float=0)
 Plot any Waveform (can be called before EventBuilder to check data)
 
void ReadFile (string, bool=true, int=0, string="out.root", bool=false, long long=0) override
 Routine to read files created by SAMPIC DAQ.
 
 ReadSampic (int last_bin_file=0, int first_bin_file=0)
 Initializer will call initializer of ReadRun class.
 
bool SkipEvent (int, int=-1) override
 Check if event has skipped flag. Will return true if event flag is false but channel does not exist in data.
 
- Public Member Functions inherited from ReadRun
TH1FBaselineCorrectionResults (int, int, float=-5, float=5, int=200)
 Histograms of the contents of baseline_correction_result.
 
void ChargeCorrelation (float, float, float, float, float, float, int, int, int, bool=false)
 Plot correlation of integrals/amplitudes between two channels.
 
floatChargeList (int, float=20, float=80, float=0, float=300, bool=1)
 Returns array with the individual "charge"/amplitude for all events of one channel.
 
TH1FChargeSpectrum (int, float, float, float=0, float=300, float=-50, float=600, int=750)
 Histogram of the "charge" spectrum for one channel.
 
int CheckBoundsX (int)
 Check if index exists in time of waveforms.
 
void CorrectBaseline (float, float=-999)
 Baseline correction constant window.
 
void CorrectBaseline_function (vector< float > &, float, float, int)
 Helper function called by CorrectBaseline()
 
void CorrectBaselineMin (int=100, double=0.5, int=0, int=0, int=2)
 Baseline correction using minimum sum ( \(\propto\) mean) in range for correction.
 
void CorrectBaselineMin (vector< float >, double, int, int)
 Wrapper for backwards compatibility.
 
void CorrectBaselineMin (vector< float >, double=0, int=2)
 Baseline correction using minimum sum ( \(\propto\) mean) in range for correction.
 
void CorrectBaselineMinSlopeRMS (int=100, bool=false, double=0.5, int=0, int=0, int=2)
 Baseline correction method searching for non-monotonic, rather constant regions.
 
void CorrectBaselineMinSlopeRMS (vector< float >, double, int, int)
 Wrapper for backwards compatibility.
 
void CorrectBaselineMinSlopeRMS (vector< float >, double=0, int=2)
 Baseline correction method searching for non-monotonic, rather constant regions.
 
void FilterAll (double=.3, double=.9, double=.2)
 Filter all waveforms.
 
int GetChannelIndex (int)
 Match channel number (wavecatcher input channel) to channel index.
 
int GetEventIndex (unsigned int)
 Returns index of a certain event number (if data files are read in parallel threads)
 
array< int, 3 > GetIntWindow (const vector< float > &, float, float, float, float, int=0)
 Determine indices for integration window for peaks.
 
array< int, 3 > GetIntWindow (const vector< float > &, int, int, int, int, int=0)
 Determine indices for integration window for peaks with bin indices.
 
array< int, 3 > GetIntWindow (TH1F *, float, float, float, float, int=0)
 Determine indices for integration window for peaks.
 
float GetPeakIntegral (const vector< float > &, float, float, float, float, int=0)
 Calculate the integral around a peak with several options explained in GetIntWindow().
 
float GetPeakIntegral (TH1F *, float, float, float, float, int=0)
 Calculate the integral around a peak with several options explained in GetIntWindow().
 
void GetTimingCFD (float=.3, float=100, float=140, double=0., bool=true, int=2, bool=false, bool=false)
 Determine the timing of the maximum peak with constant fraction discrimination.
 
TH1FGetwf (int)
 Helper that returns the waveform histogram for a certain waveform number number.
 
template<typename T >
Tgetx (double shift=0.)
 Get array of x axis (time of the bin centers) for standard wavecatcher settings.
 
doublegety (int)
 Get array of y values for a certain waveform.
 
doublegety (int, int)
 Get array of y values for a certain waveform.
 
TH1FHis_GetTimingCFD (int, float, float, int=-999)
 Plot results of GetTimingCFD()
 
TH1FHis_GetTimingCFD_diff (vector< int >, vector< int >, float, float, int=-999)
 Plot timing difference between the mean timings of two channel ranges.
 
float IndexToTime (int)
 Convert the bin number of the waveform to the time of the left bin edge
 
void IntegralFilter (vector< float >, vector< bool >, float, float, float=50, float=250, bool=false, bool=false)
 Skip events with threshold on integral.
 
TGraph2DMaxDist (int, float=0, float=300)
 Finds maximum amplitude for a given channel in time window ["from", "to"] and creates 3d map of waveforms ordered by maxima on z axis.
 
virtual int Nevents_good (int=0)
 Number of good events that are not skipped.
 
bool PlotChannel (int)
 Check if a channel index should be plotted according to ReadRun::plot_active_channels.
 
void PlotChannelAverages (bool=false)
 Plot averages only of the good, corrected waveforms for each channel.
 
void PlotChannelSums (bool=false, bool=false, double=0., double=0., int=2)
 Plot sums of all raw waveforms for each channel.
 
void PlotWFHeatmaps (float=-20, float=200, int=880, string="", float=0, EColorPalette=kRainBow)
 Plot stacks of all non-skipped waveforms for all active channels.
 
vector< boolPolarityMap (bool, int)
 Channel map of polarity changes during reading.
For parameters see ReadFile().

 
void Print_GetTimingCFD (float=100, float=140, int=0, int=-999, string="S", bool=true)
 Plot results of GetTimingCFD()
 
void Print_GetTimingCFD_diff (vector< int >, vector< int >, float=100, float=140, int=0, int=-999, float=-999, float=-999, string="RS", bool=true)
 Plot timing difference between the mean timings of two channel ranges.
 
void PrintBaselineCorrectionResults (float=-5, float=5, int=200)
 Print histogram of the baseline correction values for all channels.
 
void PrintChargeSpectrum (float, float, float=0, float=300, float=-50, float=600, int=750, float=0., float=0., int=99, int=0, bool=false)
 Plots the "charge" spectrums of all channels.
 
void PrintChargeSpectrumWF (float, float, float=0, float=300, int=1, float=0., float=0., float=0., float=0.)
 Plot waveforms of all channels for a given event number and add the determined integration windows to the plot.
 
void PrintDCR (float=15, float=85, float=0, float=300, double=3)
 Calculate (SiPM) dark count rate.
 
void PrintMaxDist (float=0, float=300)
 Finds maximum amplitude for a given channel in time window ["from", "to"] and creates 3d map of waveforms ordered by maxima on z axis.
 
void PrintSkippedEvents ()
 Prints a list of all skipped events into the terminal for diagnostics.
 
void PrintTimeDist (float=0, float=300, float=0, float=300, int=100, int=0, float=.3)
 Time distribution of maximum, CFD, or 10% - 90% rise time in a certain time window.
 
void PrintWFProjection (float=0, float=320, float=-50, float=50, int=200)
 Plots waveform projection histograms of all channels.
 
 ReadRun (int last_bin_file=0, int first_bin_file=0)
 Constructor of the class.
 
void SaveChargeLists (float=20, float=80, float=0, float=300, bool=1)
 Saves TGraphs to root file with the individual "charge"/amplitude for all events and all channels.
 
void ShiftAllToAverageCF ()
 This function shifts all waveforms to the average signal starting times for each channel.
 
void SkipEventsPerChannel (vector< float >, float=0, float=0, bool=false)
 Skip events above/below individual thresholds per channel.
 
void SkipEventsTimeDiffCut (int, int, double, double, bool=false)
 Skip events where the time difference between two channels is outside of specified range.
 
void SmoothAll (double, int)
 Smoothing all waveforms which are not skipped (for testing, careful when using for analysis!)
 
void SmoothAll (double, string="Gaus")
 Smoothing all waveforms which are not skipped (for testing, careful when using for analysis!)
 
TH1FTimeDist (int, float=0, float=300, float=0, float=300, int=100, int=0, float=.3)
 Time distribution of maximum, CFD, or 10% - 90% rise time in a certain time window.
 
int TimeToIndex (float)
 Convert time to the bin number of the waveform.
 
void UnskipAll ()
 Sets skip_event flag to false for all events, removing any previous cuts.
 
TH2FWFHeatmapChannel (int, float=-20, float=200, int=880)
 2D histogram of all non-skipped waveforms for one channel
 
TH1FWFProjectionChannel (int, int=0, int=1024, float=-50, float=50, int=200)
 Waveform projections for one channel.
 
virtual ~ReadRun ()
 Destructor.
 

Public Attributes

vector< vector< int > > ch_nr_event_storage
 Stores all channel numbers assigned to events.
 
vector< doubleevent_time_stamps
 Stores all event times.
 
vector< HitInfoReducedhitInfo
 Stores additional info about wf.
 
vector< vector< int > > wf_nr_event_storage
 Stores all waveform numbers assigned to events.
 
- Public Attributes inherited from ReadRun
vector< intactive_channels
 Stores the numbers of the active channels.
 
vector< vector< float > > amplValuessum
 Collects sums of all waveforms for each channel.
 
vector< vector< float > > baseline_correction_result
 Stores baseline correction results for CorrectBaseline() and related functions.
 
int binNumber
 Number of bins (usually 1024 samples per waveform).
 
float DAQ_factor = 250. / 4096.
 DAQ conversion factor for wavecatcher output to mV.
 
string data_path
 Path to data.
 
bool discard_original_eventnr = false
 Can be used to discard the original event numbering of the data.
 
int end_read_at_channel = -1
 See ReadRun::start_read_at_channel.
 
vector< unsigned inteventnr_storage
 Events will be stored here in the order they have been read.
 
int FirstBinFileToRead
 Number first of .bin file to be read in.
 
vector< TFitResultPtrfit_results
 Stores the fit results of PrintChargeSpectrum() for all channels and all function calls in ascending order.
 
int LastBinFileToRead
 Number of last .bin file to be read in.
 
int maxNWF = 1e7
 Maximum possible number of waveforms in data: number of active channels x number of events.
 
vector< intmaxSumBin
 Stores bin numbers where the sum of waveforms have their maximum.
 
vector< floatmean_integral
 Stores the mean integral/lightyield from PrintChargeSpectrum() for all channels.
 
int nchannels
 Number of active channels in data.
 
int nChannelsWC = 64
 Wavecatcher hardware max. number of channels (reduce/increase if not using the 64 channel crate)
 
int nevents = 0
 Number of triggered events in data.
 
int nwf
 Total number of waveforms read from data: number of active channels x number of events.
 
vector< intplot_active_channels
 Stores the numbers of the active channels which should be plotted.
 
vector< vector< float > > PrintChargeSpectrum_cal = vector(nChannelsWC, vector<float>(2, 1))
 Calibration values to normalize charge spectrum to number of photoelectrons Chennels must be ordered as in plot_active_channels.
The first entry must be the gain and the second the position of the pedestal.
 
vector< floatPrintChargeSpectrum_pars
 Starting values of the fit parameters for PrintChargeSpectrum()
 
TFileroot_out
 Stores results of analysis.
 
vector< vector< float > > rundata
 Stores waveforms.
 
vector< boolskip_event
 Stores the event numbers which should be skipped in the analysis.
 
float SP = .3125
 Sampling: ns per bin of data, sampling rate 3.2 GS/s -> 0.3125 ns.
 
float SP_inv = 1 / .3125
 1/SP
 
int start_read_at_channel = -1
 Do analysis only for limited range of channels to reduce memory usage.
 
vector< intswitch_polarity_for_channels
 Stores the channel number where the polarity should be inverted. Example use to switch polarity for channel 0 and channel 14: mymeas.switch_polarity_for_channels = {0, 14};. Needs to be called before ReadFile().
 
float tCutEndg
 End of time window for baseline correction when using ReadRun::Using_BaselineCorrection_in_file_loop.
 
float tCutg
 Start of time window for baseline correction when using ReadRun::Using_BaselineCorrection_in_file_loop.
 
vector< TFitResultPtrtiming_fit_results
 Stores the fit results of Print_GetTimingCFD() for all channels.
 
vector< vector< float > > timing_results
 Matrix to store timing of peaks from GetTimingCFD()
 
bool Using_BaselineCorrection_in_file_loop = false
 Set true for baseline correction during data reading. Needs to be called before ReadFile().
 

Protected Member Functions

void checkData (bool isBaselineCorrection=false) const override
 Primitive check to see if data has been loaded.
 

Protected Attributes

bool eventsBuilt = false
 Check if the data is grouped in events (SAMPIC needs manual event building)
 
- Protected Attributes inherited from ReadRun
int PlotChannelAverages_cnt
 Index for multiple executions of the same plotting function.
 
int PlotWFHeatmaps_cnt
 Index for multiple executions of the same plotting function.
 
int PrintChargeSpectrum_cnt
 Index for multiple executions of the same plotting function.
 
int PrintWFProjection_cnt
 Index for multiple executions of the same plotting function.
 

Additional Inherited Members

- Static Public Member Functions inherited from ReadRun
static pair< float, boolLinearInterpolation (float, float, float, float, float, bool=false)
 Simple linear interpolation for x.
 

Detailed Description

Class containing SAMPIC data loader and Event Builder.

Definition at line 8 of file ReadSampic.h.

Constructor & Destructor Documentation

◆ ReadSampic()

ReadSampic::ReadSampic ( int  last_bin_file = 0,
int  first_bin_file = 0 
)
inline

Initializer will call initializer of ReadRun class.

Parameters
last_bin_fileNumber of last .bin file to be read in.
Set it to =>1 in order to constrain the number of .bin files to be read from the target folder.
Intended for quick tests on a fraction of the full dataset or for batch reading if combined with min_no_of_bin_files_to_read.
first_bin_fileNumber first of .bin file to be read in.
Can be used to batch read large datasets in chunks of files.

Definition at line 65 of file ReadSampic.h.

Member Function Documentation

◆ checkData()

void ReadSampic::checkData ( bool  isBaselineCorrection = false) const
inlineoverrideprotectedvirtual

Primitive check to see if data has been loaded.

Reimplemented from ReadRun.

Definition at line 43 of file ReadSampic.h.

◆ EventBuilder()

void ReadSampic::EventBuilder ( double  coincidence_time_window,
vector< float thresholds = {},
vector< int channels = {},
unsigned int  min_coincidences = 1,
bool  shift_relative = false 
)

Construct events based on coincident amplitude thresholds within a time window for different channels.
It currently constructs the events from the first waveform above the set amplitude threshold.
ToDo:
.

  • Add option to create event from most coincidences, not from first
  • Add option for integral threshold? (Would require parameters for integration window -> separate function?)
    Parameters
    coincidence_time_windowEvent time window in ns. All hits within time window will be added to event.
    Starts at the beginning of the first channel above threshold and will include all channels which start recording before the end of the coincidence time window. This means the peak might be after the end of the coincidence time window. If you want to include only waveforms with the signal inside the time window you would need to subtract the length of a waveform (depends on SAMPIC settings) from the time window.
    Must be shorter than the dead time of SAMPIC (<1000 ns).
    thresholdsVector of thresholds for all channels in ascending order (lowest channel number in data is the first entry).
    If you specify only a single value {val}, it will be used for all channels. The amplitude must be larger than the specified threshold.
    channelsVector of channels to be used. If left empty {} all channels will be used. Provide the actual channel numbers from SAMPIC (e. g. {10, 28}).
    min_coincidencesMinimum number of channels above threshold. Events will be constructed if at least min_coincidences channels are above their threshold.
    shift_relativeShift waveforms relative to first timestamp in event
    NOT WORKING YET

Definition at line 228 of file ReadSampic.cc.

◆ GetCurrentChannel()

int ReadSampic::GetCurrentChannel ( int  waveform_index)
overridevirtual

Get the current channel index for a certain waveform index.

Parameters
waveform_index
Returns
Current channel index

Reimplemented from ReadRun.

Definition at line 399 of file ReadSampic.cc.

◆ GetCurrentEvent()

int ReadSampic::GetCurrentEvent ( int  waveform_index)
overridevirtual

Get the current event index for a certain waveform index.

Parameters
waveform_index
Returns
Current event index. If the waveform is not assigned to an event it will return -1.

Reimplemented from ReadRun.

Definition at line 406 of file ReadSampic.cc.

◆ GetWaveformIndex()

int ReadSampic::GetWaveformIndex ( int  eventnr,
int  channel_index 
)
overridevirtual

Returns index of a certain event number (if data files are read in parallel threads)

Parameters
eventnrEvent number as stored in the data.
channel_indexChannel index as stored in the data.
Returns
Corresponding waveform number in the internal data structure.

Reimplemented from ReadRun.

Definition at line 384 of file ReadSampic.cc.

◆ Getwf()

TH1F * ReadSampic::Getwf ( int  channelnr,
int  eventnr,
int  color = 1 
)
overridevirtual

Helper that returns the waveform histogram for a certain channel number and a certain event number.

Parameters
channelnrChannel number index (not the actual channel number)
eventnrEvent number
colorChoose color of histogram
Returns
Waveform histogram

Reimplemented from ReadRun.

Definition at line 348 of file ReadSampic.cc.

◆ PlotWF()

void ReadSampic::PlotWF ( int  wfNumber,
float  ymin = 0,
float  ymax = 0 
)

Plot any Waveform (can be called before EventBuilder to check data)

Parameters
wfNumberNumber of the waveform (hit number)
yminscale of y axis
ymaxscale of y axis

Definition at line 194 of file ReadSampic.cc.

◆ ReadFile()

void ReadSampic::ReadFile ( string  path,
bool  change_polarity = true,
int  change_sign_from_to_ch_num = 0,
string  out_file_name = "out.root",
bool  debug = false,
long long  max_nfs_to_read = 0 
)
overridevirtual

Routine to read files created by SAMPIC DAQ.

Reads the data and can already do simple data manipulation during reading:

Parameters
pathPath to the data. All files in this folder containing .bin in the file name will be read in.
change_polaritySet true to change polarity (sign) of certain channels (see change_sign_from_to_ch_num below). You can also define a list of channels where the polarity should be switched with switch_polarity_for_channels.
change_sign_from_to_ch_numAll channels \( \geq \) change_sign_from_to_ch_num will be inverted if change_polarity is true.
If negative number all channels \( \leq \) abs(change_sign_from_to_ch_num) will be inverted if change_polarity is true.
out_file_nameName of the .root file which stores the results, e. g. results.root.
debugSet true to increase the verbosity.
max_nfs_to_readMaximum number of events to be read per file. For quick testing of analysis on subset of the data.

Reimplemented from ReadRun.

Definition at line 23 of file ReadSampic.cc.

◆ SkipEvent()

bool ReadSampic::SkipEvent ( int  event_index,
int  channel_index = -1 
)
overridevirtual

Check if event has skipped flag. Will return true if event flag is false but channel does not exist in data.

Parameters
event_indexIndex of the event
channel_indexIndex of the channel

Reimplemented from ReadRun.

Definition at line 413 of file ReadSampic.cc.

Member Data Documentation

◆ ch_nr_event_storage

vector<vector<int> > ReadSampic::ch_nr_event_storage

Stores all channel numbers assigned to events.

Definition at line 100 of file ReadSampic.h.

◆ event_time_stamps

vector<double> ReadSampic::event_time_stamps

Stores all event times.

Definition at line 102 of file ReadSampic.h.

◆ eventsBuilt

bool ReadSampic::eventsBuilt = false
protected

Check if the data is grouped in events (SAMPIC needs manual event building)

Definition at line 41 of file ReadSampic.h.

◆ hitInfo

vector<HitInfoReduced> ReadSampic::hitInfo

Stores additional info about wf.

Definition at line 96 of file ReadSampic.h.

◆ wf_nr_event_storage

vector<vector<int> > ReadSampic::wf_nr_event_storage

Stores all waveform numbers assigned to events.

Definition at line 98 of file ReadSampic.h.


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