LDMX Software
TestBeamHitAnalyzer.cxx
Go to the documentation of this file.
1
9
10namespace trigscint {
11
12TestBeamHitAnalyzer::TestBeamHitAnalyzer(const std::string& name,
13 framework::Process& process)
14 : Analyzer(name, process) {}
15
16void TestBeamHitAnalyzer::configure(framework::config::Parameters& parameters) {
17 inputCol_ = parameters.getParameter<std::string>("inputCollection");
18 inputPassName_ = parameters.getParameter<std::string>("inputPassName");
19 peds_ = parameters.getParameter<std::vector<double> >("pedestals");
20 startSample_ = parameters.getParameter<int>("startSample");
21
22 std::cout << " [ TestBeamHitAnalyzer ] In configure(), got parameters "
23 << "\n\t inputCollection = " << inputCol_
24 << "\n\t inputPassName = " << inputPassName_
25 << "\n\t startSample = " << startSample_
26 << "\n\t pedestals[0] = " << peds_[0] << "\t." << std::endl;
27
28 return;
29}
30
31void TestBeamHitAnalyzer::analyze(const framework::Event& event) {
32 const auto channels{
33 event.getCollection<trigscint::TestBeamHit>(inputCol_, inputPassName_)};
34
35 int evNb = event.getEventNumber();
36 // int nChan = channels.size();
37 int leadBar = -1;
38 int subleadBar = -1;
39 float peLead = -1;
40 float peSublead = -1;
41 bool existsIntermediatePE = false;
42
43 for (auto chan : channels) {
44 int bar = chan.getBarID();
45 float PE = chan.getPE();
46 if (evNb < nEv &&
47 bar < nChannels) { // stick within the predefined histogram array
48 hEvDisp->Fill(evNb, bar, PE);
49 } // if within event display range
50 if (PE > peLead) {
51 peSublead = peLead;
52 subleadBar = leadBar;
53 peLead = PE;
54 leadBar = bar;
55 } else if (PE > peSublead) { // need a specific check, bars not sorted in
56 // PE so leadPE might be found before or after
57 peSublead = PE;
58 subleadBar = bar;
59 }
60 hPE[bar]->Fill(PE);
61 if (chan.getQualityFlag() == 0 && bar < 12 && 15 < PE && PE < 40)
62 existsIntermediatePE = true;
63 } // over channels
64
65 if (existsIntermediatePE && fillNb < nEv) {
66 for (auto chan : channels) {
67 int bar = chan.getBarID();
68 if (bar < 12) { // stick within the predefined histogram array
69 float PE = chan.getPE();
70 hEvDispPE->Fill(fillNb, bar, PE);
71 }
72 } // if within event display range
73 fillNb++;
74 } // if PE range triggers writing these event displays
75
76 if (subleadBar == -1) {
77 subleadBar = leadBar;
78 peSublead = peLead;
79 }
80 hPEVsDelta[leadBar]->Fill(leadBar - subleadBar, peLead);
81 hDeltaPEVsDelta[leadBar]->Fill(leadBar - subleadBar, peLead - peSublead);
82
83 // if ( (subleadBar%2) == (leadBar%2) ) // in same layer (or even, hit).
84 // skip if we're not seeing a max across layers
85 // return; -- on the other hand this is evident from the plot: delta is even
86
87 hPEmaxVsDelta->Fill(leadBar - subleadBar, peLead);
88
89 return;
90}
91
92void TestBeamHitAnalyzer::onProcessStart() {
93 std::cout << "\n\n Process starts! My analyzer should do something -- like "
94 "print this \n\n"
95 << std::endl;
96
97 getHistoDirectory();
98
99 int nTimeSamp = 40;
100 int PEmax = 400;
101 int nPEbins = 2 * PEmax;
102 // float Qmax = PEmax / (6250. / 4.e6);
103 // float Qmin = -10;
104 // int nQbins = (Qmax - Qmin) / 4;
105
106 for (int iB = 0; iB < nChannels; iB++) {
107 hPE[iB] = new TH1F(Form("hPE_chan%i", iB), Form(";PE, chan%i", iB), nPEbins,
108 0, PEmax);
109 hPEinClusters[iB] = new TH1F(Form("hPEinClusters_chan%i", iB),
110 Form(";PE, chan%i", iB), nPEbins, 0, PEmax);
111 hPEVsDelta[iB] = new TH2F(Form("hPEVsDelta_chan%i", iB),
112 Form(";#Delta_{barID};PE, chan%i has max PE", iB),
113 nChannels + 1, -nChannels / 2 - 0.5,
114 nChannels / 2 + 0.5, nPEbins, 0, PEmax);
115 hDeltaPEVsDelta[iB] = new TH2F(
116 Form("hDeltaPEVsDelta_chan%i", iB),
117 Form(";#Delta_{barID};#Delta_PE, chan%i has max PE", iB), nChannels + 1,
118 -nChannels / 2 - 0.5, nChannels / 2 + 0.5, nPEbins, 0, PEmax);
119 }
120
121 // make event displays for events where there are channels with weird event
122 // hit PE counts
123 for (int iE = 0; iE < nEv; iE++) {
124 for (int iB = 0; iB < nChannels; iB++) {
125 hOut[iE][iB] =
126 new TH1F(Form("hCharge_chan%i_nv%i", iB, iE),
127 Form(";time sample; Q, chan %i, ev %i [fC]", iB, iE),
128 nTimeSamp, -0.5, nTimeSamp - 0.5);
129 }
130 }
131
132 hPEmaxVsDelta =
133 new TH2F("hPEmaxVsDelta", ";#Delta_{barID};PE, max hit", nChannels,
134 -nChannels / 2, nChannels / 2, nPEbins, 0, PEmax);
135 hEvDisp = new TH2F(Form("hEvDisp_ev%i", nEv), ";Event number; Bar ID; PE",
136 nEv, 0.5, nEv + 0.5, nChannels, -0.5, nChannels - 0.5);
137 hEvDispPE = new TH2F("hEvDispPEcut", ";Event number; Bar ID; PE", nEv, 0.5,
138 nEv + 0.5, nChannels, -0.5, nChannels - 0.5);
139
140 fillNb = 0;
141
142 return;
143}
144
145void TestBeamHitAnalyzer::onProcessEnd() { return; }
146
147} // namespace trigscint
148
149DECLARE_ANALYZER_NS(trigscint, TestBeamHitAnalyzer)
#define DECLARE_ANALYZER_NS(NS, CLASS)
Macro which allows the framework to construct an analyzer given its name during configuration.
Implements an event buffer system for storing event data.
Definition Event.h:41
Class which represents the process under execution.
Definition Process.h:36
Class encapsulating parameters for configuring a processor.
Definition Parameters.h:27
T getParameter(const std::string &name) const
Retrieve the parameter of the given name.
Definition Parameters.h:89
int getBarID() const
Get the bar ID.
This class represents the linearised QIE output from the trigger scintillator, in charge (fC).
Definition TestBeamHit.h:24