LDMX Software
trigscint::QIEAnalyzer Class Reference

Public Member Functions

 QIEAnalyzer (const std::string &name, framework::Process &process)
 
void configure (framework::config::Parameters &parameters) override
 Callback for the EventProcessor to configure itself from the given set of parameters.
 
void analyze (const framework::Event &event) override
 Process the event and make histograms or summaries.
 
void onProcessStart () override
 Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.
 
void onProcessEnd () override
 Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.
 
float convertToID (float yVal)
 
- Public Member Functions inherited from framework::Analyzer
 Analyzer (const std::string &name, Process &process)
 Class constructor.
 
virtual void process (Event &event) final
 Processing an event for an Analyzer is calling analyze.
 
virtual void beforeNewRun (ldmx::RunHeader &run_header) final
 Don't allow Analyzers to add parameters to the run header.
 
- Public Member Functions inherited from framework::EventProcessor
 DECLARE_FACTORY (EventProcessor, EventProcessor *, const std::string &, Process &)
 declare that we have a factory for this class
 
 EventProcessor (const std::string &name, Process &process)
 Class constructor.
 
virtual ~EventProcessor ()=default
 Class destructor.
 
virtual void onNewRun (const ldmx::RunHeader &run_header)
 Callback for the EventProcessor to take any necessary action when the run being processed changes.
 
virtual void onFileOpen (EventFile &event_file)
 Callback for the EventProcessor to take any necessary action when a new event input ROOT file is opened.
 
virtual void onFileClose (EventFile &event_file)
 Callback for the EventProcessor to take any necessary action when a event input ROOT file is closed.
 
template<class T >
const T & getCondition (const std::string &condition_name)
 Access a conditions object for the current event.
 
TDirectory * getHistoDirectory ()
 Access/create a directory in the histogram file for this event processor to create histograms and analysis tuples.
 
void setStorageHint (framework::StorageControl::Hint hint)
 Mark the current event as having the given storage control hint from this module_.
 
void setStorageHint (framework::StorageControl::Hint hint, const std::string &purposeString)
 Mark the current event as having the given storage control hint from this module and the given purpose string.
 
int getLogFrequency () const
 Get the current logging frequency from the process.
 
int getRunNumber () const
 Get the run number from the process.
 
std::string getName () const
 Get the processor name.
 
void createHistograms (const std::vector< framework::config::Parameters > &histos)
 Internal function which is used to create histograms passed from the python configuration @parma histos vector of Parameters that configure histograms to create.
 

Private Attributes

std::vector< std::vector< TH1F * > > v_charge_vs_time_
 
std::string input_col_
 
std::string input_pass_name_ {""}
 
std::vector< double > peds_
 
std::vector< double > gain_
 
int start_sample_ {0}
 
int n_ev_ {200}
 
int n_channels_ {16}
 
TH1F * h_out_ [200][16]
 
TH1F * h_pe_ [16]
 
TH2F * h_pe_vs_t_ [16]
 
TH2F * h_ped_subtracted_avg_q_vs_t_ [16]
 
TH2F * h_ped_subtracted_tot_q_vs_ped_ [16]
 
TH2F * h_ped_subtracted_tot_q_vs_n_ [16]
 
TH2F * h_tot_q_vs_ped_ [16]
 
TH2F * h_ped_subtracted_pe_vs_n_ [16]
 
TH2F * h_ped_subtracted_pe_vs_t_ [16]
 
TH2F * h_avg_q_vs_t_ [16]
 
TH2F * h_tdc_fire_chan_vs_event_
 
double y_offset_ {35.}
 
double y_to_id_factor_ {50. / 80.}
 

Additional Inherited Members

- Protected Member Functions inherited from framework::EventProcessor
void abortEvent ()
 Abort the event immediately.
 
- Protected Attributes inherited from framework::EventProcessor
HistogramPool histograms_
 helper object for making and filling histograms
 
NtupleManagerntuple_ {NtupleManager::getInstance()}
 Manager for any ntuples.
 
logging::logger the_log_
 The logger for this EventProcessor.
 

Detailed Description

Definition at line 23 of file QIEAnalyzer.h.

Constructor & Destructor Documentation

◆ QIEAnalyzer()

trigscint::QIEAnalyzer::QIEAnalyzer ( const std::string & name,
framework::Process & process )

Definition at line 12 of file QIEAnalyzer.cxx.

13 : Analyzer(name, process) {}
virtual void process(Event &event) final
Processing an event for an Analyzer is calling analyze.
Analyzer(const std::string &name, Process &process)
Class constructor.

Member Function Documentation

◆ analyze()

void trigscint::QIEAnalyzer::analyze ( const framework::Event & event)
overridevirtual

Process the event and make histograms or summaries.

Parameters
eventThe Event to analyze

Implements framework::Analyzer.

Definition at line 32 of file QIEAnalyzer.cxx.

32 {
33 const auto channels{event.getCollection<trigscint::EventReadout>(
34 input_col_, input_pass_name_)};
35
36 int ev_nb = event.getEventNumber();
37 // while (evNb < 0 ) {
38 // ldmx_log(debug) << "event number = " << evNb << " < 0; incrementing event
39 // number "; evNb++;
40 //}
41 int num_chan = channels.size();
42 ldmx_log(debug) << "in event " << ev_nb << "; num_channels = " << num_chan;
43
44 for (auto chan : channels) {
45 std::vector<float> q = chan.getQ();
46 std::vector<float> q_err = chan.getQError();
47 std::vector<int> tdc = chan.getTDC();
48 // int nTimeSamp = q.size();
49 int bar = chan.getChanID();
50 float q_tot = 0;
51 float q_ped_subtracted_avg = 0;
52 int first_t = start_sample_ - 1;
53 int n_samp_above = 0;
54 int n_samp_above_event_ped = 0;
55 float subtr_pe = 0;
56 float subtr_q = 0;
57 float ped = chan.getPedestal();
58 for (int i_t = 0; i_t < q.size(); i_t++) {
59 ldmx_log(debug) << "in event " << ev_nb << "; channel " << bar
60 << ", got charge[" << i_t << "] = " << q.at(i_t);
61 if (ev_nb < n_ev_ && bar < n_channels_) {
62 // stick within the predefined histogram array
63 // h_out_[evNb][bar]->Fill(iT+start_sample_, q.at(iT));
64 h_out_[ev_nb][bar]->SetBinContent(i_t + start_sample_, q.at(i_t));
65 h_out_[ev_nb][bar]->SetBinError(i_t + start_sample_,
66 fabs(q_err.at(i_t)));
67 if (tdc.at(i_t) < 63) {
68 ldmx_log(info) << "Found fired TDC = " << tdc.at(i_t)
69 << " at time sample " << i_t << " in channel " << bar
70 << " and event " << ev_nb;
71 // for some reason, the style settings are washed out later...
72 h_out_[ev_nb][bar]->SetLineColor(kRed + 1);
73 h_out_[ev_nb][bar]->SetMarkerColor(
74 h_out_[ev_nb][bar]->GetLineColor());
75 h_out_[ev_nb][bar]->SetMarkerSize(0.2);
76
77 if (i_t + start_sample_ > 0)
78 h_tdc_fire_chan_vs_event_->Fill(bar, ev_nb, i_t + start_sample_);
79 else
80 h_tdc_fire_chan_vs_event_->Fill(bar, ev_nb);
81 }
82 } // if within the number of events to plot individually
83 if (q.at(i_t) > 2 * fabs(peds_[bar])) {
84 // integrate all charge well above
85 // ped to convert to a PE count
86 q_tot += q.at(i_t);
87 q_ped_subtracted_avg += q.at(i_t) - chan.getPedestal();
88 // peds_[ bar ];
89 n_samp_above++;
90 ldmx_log(debug) << " above channel overall pedestal: " << q.at(i_t)
91 << " > " << 2 * fabs(peds_[bar]);
92
93 // keep track of first time sample above threshold
94 if (first_t == start_sample_ - 1) first_t = start_sample_ + i_t;
95 } // if above threshold
96 if (q.at(i_t) > ped) {
97 subtr_q += q.at(i_t) - peds_[bar];
98 n_samp_above_event_ped++;
99 ldmx_log(debug) << " above channel event pedestal: " << q.at(i_t)
100 << " > " << ped;
101 } // if above channel event pedestal
102 } // over time samples
103 float pe = q_tot * 6250. / gain_[bar];
104 subtr_pe = subtr_q * 6250. / gain_[bar];
105 h_tot_q_vs_ped_[bar]->Fill(ped, q_tot);
106 h_pe_[bar]->Fill(pe);
107 h_pe_vs_t_[bar]->Fill(first_t, pe);
108 if (n_samp_above > 0) {
109 q_ped_subtracted_avg /= n_samp_above;
110 h_ped_subtracted_avg_q_vs_t_[bar]->Fill(first_t, q_ped_subtracted_avg);
111 h_avg_q_vs_t_[bar]->Fill(first_t, q_tot / n_samp_above);
112 }
113 // if (chan.getPedestal() < 40. ) {
114 // subtrQ = subtrPE/(6250./4.e6); //undo conversion
115 ldmx_log(debug) << "filling qTot histograms";
116 h_ped_subtracted_tot_q_vs_ped_[bar]->Fill(ped, subtr_q);
117 if (ped < 40) // avoid case where we have saturation and a plateau as much
118 // as possible
119 h_ped_subtracted_tot_q_vs_n_[bar]->Fill(n_samp_above_event_ped, subtr_q);
120 h_ped_subtracted_pe_vs_n_[bar]->Fill(n_samp_above_event_ped, subtr_pe);
121 h_ped_subtracted_pe_vs_t_[bar]->Fill(first_t, subtr_pe);
122 ldmx_log(debug) << " done filling qTot histograms";
123 // }
124 } // over channels
125
126 return;
127}
This class represents the linearised QIE output from the trigger scintillator, in charge (fC).

◆ configure()

void trigscint::QIEAnalyzer::configure ( framework::config::Parameters & parameters)
overridevirtual

Callback for the EventProcessor to configure itself from the given set of parameters.

The parameters a processor has access to are the member variables of the python class in the sequence that has className equal to the EventProcessor class name.

For an example, look at MyProcessor.

Parameters
parametersParameters for configuration.

Reimplemented from framework::EventProcessor.

Definition at line 15 of file QIEAnalyzer.cxx.

15 {
16 input_col_ = parameters.get<std::string>("inputCollection");
17 input_pass_name_ = parameters.get<std::string>("inputPassName");
18 peds_ = parameters.get<std::vector<double> >("pedestals");
19 gain_ = parameters.get<std::vector<double> >("gain");
20 start_sample_ = parameters.get<int>("startSample");
21
22 ldmx_log(trace) << "In configure(), got parameters "
23 << "\n\t inputCollection = " << input_col_
24 << "\n\t inputPassName = " << input_pass_name_
25 << "\n\t startSample = " << start_sample_
26 << "\n\t pedestals[0] = " << peds_[0]
27 << "\n\t gain[0] = " << gain_[0] << "\t.";
28
29 return;
30}
const T & get(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:78

References framework::config::Parameters::get().

◆ convertToID()

float trigscint::QIEAnalyzer::convertToID ( float yVal)
inline

Definition at line 37 of file QIEAnalyzer.h.

37{ return (yVal + y_offset_) * y_to_id_factor_; }

◆ onProcessEnd()

void trigscint::QIEAnalyzer::onProcessEnd ( )
overridevirtual

Callback for the EventProcessor to take any necessary action when the processing of events finishes, such as calculating job-summary quantities.

Reimplemented from framework::EventProcessor.

Definition at line 213 of file QIEAnalyzer.cxx.

213{ return; }

◆ onProcessStart()

void trigscint::QIEAnalyzer::onProcessStart ( )
overridevirtual

Callback for the EventProcessor to take any necessary action when the processing of events starts, such as creating histograms.

Reimplemented from framework::EventProcessor.

Definition at line 129 of file QIEAnalyzer.cxx.

129 {
130 ldmx_log(trace)
131 << "\n\n Process starts! My analyzer should do something -- like "
132 "print this \n\n";
134
135 int n_time_samp = 68; // 40
136 int p_emax = 100;
137 int n_p_ebins = 5 * p_emax;
138 float qmax = p_emax / (6250. / 4.e6);
139 float qmin = -10;
140 int n_qbins = (qmax - qmin) / 4;
141
142 ldmx_log(debug) << "Setting up histograms... ";
143
144 for (int i_b = 0; i_b < n_channels_; i_b++) {
145 h_pe_[i_b] = new TH1F(Form("h_pe_chan%i", i_b), Form(";PE, chan%i", i_b),
146 n_p_ebins, 0, p_emax);
147 h_pe_vs_t_[i_b] = new TH2F(
148 Form("h_pe_vs_t__chan%i", i_b),
149 Form(";First time sample above summing threshold;PE, chan%i", i_b),
150 n_time_samp + 1, -1.5, n_time_samp - 0.5, n_p_ebins, 0, p_emax);
151 h_ped_subtracted_avg_q_vs_t_[i_b] = new TH2F(
152 Form("hPedSubtrAvgQvsT_chan%i", i_b),
153 Form(";First time sample above threshold;Pedestal subtracted average "
154 "Q, chan%i [fC]",
155 i_b),
156 n_time_samp + 1, -1.5, n_time_samp - 0.5, n_qbins / 10, qmin,
157 qmax / 10.);
158 h_ped_subtracted_tot_q_vs_ped_[i_b] =
159 new TH2F(Form("hPedSubtrTotQvsPed_chan%i", i_b),
160 Form(";Channel event pedestal [fC];Event pedestal subtracted "
161 "total Q, chan%i [fC]",
162 i_b),
163 1010, qmin, 1000, 10010, -10,
164 10000); // nQbins/2,Qmin,Qmax/5., nQbins,Qmin,2*Qmax);
165 h_ped_subtracted_tot_q_vs_n_[i_b] =
166 new TH2F(Form("hPedSubtrTotQvsN_chan%i", i_b),
167 Form(";Number of time samples added; Event pedestal "
168 "subtracted total Q, chan%i [fC]",
169 i_b),
170 n_time_samp + 1, -1.5, n_time_samp - 0.5, 10010, -10, 10000);
171 h_tot_q_vs_ped_[i_b] = new TH2F(
172 Form("h_tot_q_vs_ped__chan%i", i_b),
173 Form(";Channel event pedestal [fC];Event total Q, chan%i [fC]", i_b),
174 1010, qmin, 1000, 10010, -10,
175 10000); // nQbins/2,Qmin,Qmax/5., nQbins,Qmin,2*Qmax);
176 h_ped_subtracted_pe_vs_n_[i_b] = new TH2F(
177 Form("hPedSubtrPEvsN_chan%i", i_b),
178 Form(";Number of time samples above threshold;Pedestal "
179 "subtracted PE, chan%i [fC]",
180 i_b),
181 n_time_samp + 1, -1.5, n_time_samp - 0.5, n_p_ebins, 0, p_emax);
182 h_ped_subtracted_pe_vs_t_[i_b] = new TH2F(
183 Form("hPedSubtrPEvsT_chan%i", i_b),
184 Form(";First time sample above threshold;Pedestal subtracted "
185 "PE, chan%i [fC]",
186 i_b),
187 n_time_samp + 1, -1.5, n_time_samp - 0.5, n_p_ebins, 0, p_emax);
188 h_avg_q_vs_t_[i_b] = new TH2F(
189 Form("h_avg_q_vs_t_chan%i", i_b),
190 Form(";First time sample above threshold;Average Q, chan%i [fC]", i_b),
191 n_time_samp + 1, -1.5, n_time_samp - 0.5, n_qbins / 10, qmin,
192 qmax / 10);
193 }
194
195 for (int i_e = 0; i_e < n_ev_; i_e++) {
196 for (int i_b = 0; i_b < n_channels_; i_b++) {
197 h_out_[i_e][i_b] =
198 new TH1F(Form("hCharge_chan%i_ev%i", i_b, i_e),
199 Form(";time sample; Q, channel %i, event %i [fC]", i_b, i_e),
200 n_time_samp, -0.5, n_time_samp - 0.5);
201 }
202 }
203
204 h_tdc_fire_chan_vs_event_ = new TH2F(
205 "h_tdc_fire_chan_vs_event", ";channel with TDC < 63;event number",
206 n_channels_, -0.5, n_channels_ - 0.5, n_ev_, 0, n_ev_);
207
208 ldmx_log(debug) << "done setting up histograms";
209
210 return;
211}
TDirectory * getHistoDirectory()
Access/create a directory in the histogram file for this event processor to create histograms and ana...

Member Data Documentation

◆ gain_

std::vector<double> trigscint::QIEAnalyzer::gain_
private

Definition at line 46 of file QIEAnalyzer.h.

◆ h_avg_q_vs_t_

TH2F* trigscint::QIEAnalyzer::h_avg_q_vs_t_[16]
private

Definition at line 64 of file QIEAnalyzer.h.

◆ h_out_

TH1F* trigscint::QIEAnalyzer::h_out_[200][16]
private

Definition at line 55 of file QIEAnalyzer.h.

◆ h_pe_

TH1F* trigscint::QIEAnalyzer::h_pe_[16]
private

Definition at line 56 of file QIEAnalyzer.h.

◆ h_pe_vs_t_

TH2F* trigscint::QIEAnalyzer::h_pe_vs_t_[16]
private

Definition at line 57 of file QIEAnalyzer.h.

◆ h_ped_subtracted_avg_q_vs_t_

TH2F* trigscint::QIEAnalyzer::h_ped_subtracted_avg_q_vs_t_[16]
private

Definition at line 58 of file QIEAnalyzer.h.

◆ h_ped_subtracted_pe_vs_n_

TH2F* trigscint::QIEAnalyzer::h_ped_subtracted_pe_vs_n_[16]
private

Definition at line 62 of file QIEAnalyzer.h.

◆ h_ped_subtracted_pe_vs_t_

TH2F* trigscint::QIEAnalyzer::h_ped_subtracted_pe_vs_t_[16]
private

Definition at line 63 of file QIEAnalyzer.h.

◆ h_ped_subtracted_tot_q_vs_n_

TH2F* trigscint::QIEAnalyzer::h_ped_subtracted_tot_q_vs_n_[16]
private

Definition at line 60 of file QIEAnalyzer.h.

◆ h_ped_subtracted_tot_q_vs_ped_

TH2F* trigscint::QIEAnalyzer::h_ped_subtracted_tot_q_vs_ped_[16]
private

Definition at line 59 of file QIEAnalyzer.h.

◆ h_tdc_fire_chan_vs_event_

TH2F* trigscint::QIEAnalyzer::h_tdc_fire_chan_vs_event_
private

Definition at line 66 of file QIEAnalyzer.h.

◆ h_tot_q_vs_ped_

TH2F* trigscint::QIEAnalyzer::h_tot_q_vs_ped_[16]
private

Definition at line 61 of file QIEAnalyzer.h.

◆ input_col_

std::string trigscint::QIEAnalyzer::input_col_
private

Definition at line 43 of file QIEAnalyzer.h.

◆ input_pass_name_

std::string trigscint::QIEAnalyzer::input_pass_name_ {""}
private

Definition at line 44 of file QIEAnalyzer.h.

44{""};

◆ n_channels_

int trigscint::QIEAnalyzer::n_channels_ {16}
private

Definition at line 51 of file QIEAnalyzer.h.

51{16};

◆ n_ev_

int trigscint::QIEAnalyzer::n_ev_ {200}
private

Definition at line 50 of file QIEAnalyzer.h.

50{200};

◆ peds_

std::vector<double> trigscint::QIEAnalyzer::peds_
private

Definition at line 45 of file QIEAnalyzer.h.

◆ start_sample_

int trigscint::QIEAnalyzer::start_sample_ {0}
private

Definition at line 47 of file QIEAnalyzer.h.

47{0};

◆ v_charge_vs_time_

std::vector<std::vector<TH1F*> > trigscint::QIEAnalyzer::v_charge_vs_time_
private

Definition at line 40 of file QIEAnalyzer.h.

◆ y_offset_

double trigscint::QIEAnalyzer::y_offset_ {35.}
private

Definition at line 67 of file QIEAnalyzer.h.

67{35.};

◆ y_to_id_factor_

double trigscint::QIEAnalyzer::y_to_id_factor_ {50. / 80.}
private

Definition at line 68 of file QIEAnalyzer.h.

68{50. / 80.};

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