wavecatcher-analysis
Loading...
Searching...
No Matches
ReadSampic.h
Go to the documentation of this file.
1
2#ifndef _ReadSampic
3#define _ReadSampic
4
5#include "ReadRun.h"
6// #include "utils/Helpers.h"
7
8class ReadSampic : public virtual ReadRun
9{
10private:
12 double coincidence_time;
13
14 #pragma pack(1)
15 typedef struct {
16 int HitNumber;
17 unsigned char Channel;
18 double FirstSampleTimeStamp;
19 unsigned short RawToTValue;
20 float TOTValue;
21 float Time;
22 float Baseline;
23 float Amplitude;
24 unsigned char FirstCellIndex;
25 unsigned char WaveformSize;
26 } HitStructInfoForWaveformAndMeasurements_t;
27
28 typedef struct {
29 int HitNumber;
30 unsigned char Channel;
31 double FirstSampleTimeStamp;
32 short RawToTValue;
33 float TOTValue;
34 unsigned char FirstCellIndex;
35 unsigned char WaveformSize;
36 } HitStructInfoForWaveformOnly_t;
37 #pragma pack()
38
39protected:
41 bool eventsBuilt = false;
42
43 void checkData(bool isBaselineCorrection = false) const override {
44 if (hitInfo.empty()) {
45 throw runtime_error(
46 "Error: No data has been loaded yet!\n \
47 Please call ReadFile() before calling any functions which manipulate data.\n \
48 Aborting execution."
49 );
50 }
51 else if (!eventsBuilt && !isBaselineCorrection) {
52 throw runtime_error(
53 "Error: Sampic::EventBuilder() must be called before any analysis can be performed."
54 );
55 }
56 }
57
58public:
66
67 void ReadFile(string, bool = true, int = 0, string = "out.root", bool = false, long long = 0) override;
68
69 void PlotWF(int, float = 0, float = 0);
70
71 void EventBuilder(double, vector<float> = {}, vector<int> channels = {}, unsigned int = 1, bool = false);
72
73 TH1F* Getwf(int, int, int = 1) override; // channel, eventnr, color
74 int GetWaveformIndex(int, int) override;
75 int GetCurrentChannel(int) override;
76 int GetCurrentEvent(int) override;
77
78 bool SkipEvent(int, int = -1) override;
79
80 // ---------------- variables ----------------
81 #pragma pack(1)
87 int EventNumber = -1;
88 float Max = 0;
89 float Min = 0;
90 bool IsEventCandidate = false;
91 bool IsEvent = false;
92 };
93 #pragma pack()
94
103};
104#endif
int start_read_at_channel
Do analysis only for limited range of channels to reduce memory usage.
Definition ReadRun.h:277
Class containing SAMPIC data loader and Event Builder.
Definition ReadSampic.h:9
int GetCurrentEvent(int) override
Get the current event index for a certain waveform index.
vector< HitInfoReduced > hitInfo
Stores additional info about wf.
Definition ReadSampic.h:96
void PlotWF(int, float=0, float=0)
Plot any Waveform (can be called before EventBuilder to check data)
ReadSampic(int last_bin_file=0, int first_bin_file=0)
Initializer will call initializer of ReadRun class.
Definition ReadSampic.h:65
int GetWaveformIndex(int, int) override
Returns index of a certain event number (if data files are read in parallel threads)
void checkData(bool isBaselineCorrection=false) const override
Primitive check to see if data has been loaded.
Definition ReadSampic.h:43
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...
void ReadFile(string, bool=true, int=0, string="out.root", bool=false, long long=0) override
Routine to read files created by SAMPIC DAQ.
Definition ReadSampic.cc:23
TH1F * Getwf(int, int, int=1) override
Helper that returns the waveform histogram for a certain channel number and a certain event number.
vector< vector< int > > wf_nr_event_storage
Stores all waveform numbers assigned to events.
Definition ReadSampic.h:98
bool eventsBuilt
Check if the data is grouped in events (SAMPIC needs manual event building)
Definition ReadSampic.h:41
int GetCurrentChannel(int) override
Get the current channel index for a certain waveform index.
vector< double > event_time_stamps
Stores all event times.
Definition ReadSampic.h:102
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 i...
vector< vector< int > > ch_nr_event_storage
Stores all channel numbers assigned to events.
Definition ReadSampic.h:100
Stores additional information with the waveforms.
Definition ReadSampic.h:83