Belle II Software development
TOPLaserHitSelectorModule.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9// own include
10#include <top/modules/TOPGainEfficiencyMonitor/TOPLaserHitSelectorModule.h>
11
12// framework - DataStore
13#include <framework/datastore/StoreArray.h>
14
15// Hit classes
16#include <top/dataobjects/TOPDigit.h>
17
18// standard libraries
19#include <iostream>
20#include <sstream>
21#include <iomanip>
22
23#include <map>
24
25using namespace std;
26
27namespace Belle2 {
33 //-----------------------------------------------------------------
35 //-----------------------------------------------------------------
36
37 REG_MODULE(TOPLaserHitSelector);
38
39
40 //-----------------------------------------------------------------
41 // Implementation
42 //-----------------------------------------------------------------
43
45 {
46 // Set description()
47 setDescription("TOP pixel-by-pixel gain analysis - first step : create hit timing-pulse charge histogram");
48
49 // Add parameters
50 m_timeHistogramBinning.push_back(140);//140 for time axis
51 m_timeHistogramBinning.push_back(-55);//-25
52 m_timeHistogramBinning.push_back(15);//45
53 m_chargeHistogramBinning.push_back(125);//for height
54 m_chargeHistogramBinning.push_back(-50);//
55 m_chargeHistogramBinning.push_back(2450);//
56 m_chargeHistogramBinning.push_back(150);//for integral
57 m_chargeHistogramBinning.push_back(-500);//
58 m_chargeHistogramBinning.push_back(29500);//
59
60 addParam("timeHistogramBinning", m_timeHistogramBinning,
61 "histogram binning (the number of bins, lower limit, upper limit) for hit timing distribution (should be integer value)",
63 addParam("chargeHistogramBinning", m_chargeHistogramBinning,
64 "histogram binning (the number of bins, lower limit, upper limit) for pulse charge distribution (should be integer value)",
66
67 // addParam("calibrationChannel", m_calibrationChannel, "asic channel # where calibration pulses routed",
68 // (unsigned)0);
69 addParam("useDoublePulse", m_useDoublePulse,
70 "set true when you require both of double calibration pulses for reference timing. You need to enable offline feature extraction.",
71 (bool)true);
72 addParam("minHeightFirstCalPulse", m_calibrationPulseThreshold1, "pulse height threshold for the first cal. pulse",
73 (float)300.);
74 addParam("minHeightSecondCalPulse", m_calibrationPulseThreshold2, "pulse height threshold for the second cal. pulse",
75 (float)100.);
76 addParam("nominalDeltaT", m_calibrationPulseInterval, "nominal DeltaT (time interval of the double calibration pulses) [ns]",
77 (float)25.5);
78 addParam("nominalDeltaTRange", m_calibrationPulseIntervalRange,
79 "acceptable DeltaT range from the nominal value before calibration [ns]",
80 (float)2.);
81 addParam("windowSelect", m_windowSelect,
82 "select window number (All=0, Odd=2, Even=1)",
83 0);
84 addParam("includePrimaryChargeShare", m_includePrimaryChargeShare,
85 "set ture when you require without primary chargeshare cut for making 2D histogram",
86 (bool)false);
87 addParam("includeAllChargeShare", m_includeAllChargeShare,
88 "set ture when you require without all chargeshare cut for making 2D histogram",
89 (bool)false);
90 }
91
93
95 {
96 REG_HISTOGRAM;
97
99 digits.isRequired();
100 }
101
103 {
104 for (int iPixel = 0 ; iPixel < c_NPixelPerModule * c_NModule ; iPixel++) {
105
106 short slotId = (iPixel / c_NPixelPerModule) + 1;
107 short pixelId = (iPixel % c_NPixelPerModule) + 1;
108 short pmtId = ((pixelId - 1) % c_NPixelPerRow) / c_NChannelPerPMTRow + c_NPMTPerRow * ((pixelId - 1) / (c_NPixelPerModule / 2)) + 1;
109 short pmtChId = (pixelId - 1) % c_NChannelPerPMTRow + c_NChannelPerPMTRow * ((pixelId - 1) %
110 (c_NPixelPerModule / 2) / c_NPixelPerRow) + 1;
111
112 std::ostringstream pixelstr;
113 pixelstr << "s" << std::setw(2) << std::setfill('0') << slotId << "_PMT"
114 << std::setw(2) << std::setfill('0') << pmtId
115 << "_" << std::setw(2) << std::setfill('0') << pmtChId;
116
117 std::ostringstream hnameForgain;
118 hnameForgain << "hTimeHeight_gain_" << pixelstr.str();
119 std::ostringstream htitleForgain;
120 htitleForgain << "2D distribution of hit timing and pulse height for Gain" << pixelstr.str();
121 m_TimeHeightHistogramForFit[iPixel] = new TH2F(hnameForgain.str().c_str(), htitleForgain.str().c_str(),
125
126 std::ostringstream hnameForIntegral;
127 hnameForIntegral << "hTimeIntegral_gain_" << pixelstr.str();
128 std::ostringstream htitleForIntegral;
129 htitleForIntegral << "2D distribution of hit timing and integral for Gain" << pixelstr.str();
130 m_TimeIntegralHistogramForFit[iPixel] = new TH2F(hnameForIntegral.str().c_str(), htitleForIntegral.str().c_str(),
134
135 std::ostringstream hnameForeff;
136 hnameForeff << "hTimeHeight_efficiency_" << pixelstr.str();
137 std::ostringstream htitleForeff;
138 htitleForeff << "2D distribution of hit timing and pulse height for efficiency" << pixelstr.str();
139 m_TimeHeightHistogramForHitRate[iPixel] = new TH2F(hnameForeff.str().c_str(), htitleForeff.str().c_str(),
143 }
144
145 const short nAsic = c_NPixelPerModule / c_NChannelPerAsic * c_NChannelPerPMT;
146 m_nCalPulseHistogram = new TH1F("hNCalPulse", "number of calibration pulses identificed for each asic",
147 nAsic, -0.5, nAsic - 0.5);
148 }
149
151 {
152 }
153
155 {
157
158 std::map<short, float> refTimingMap;//map of pixel id and time
159 std::map<short, std::vector<hitInfo_t> >
160 calPulsesMap;//map of pixel id and hit information (timing and pulse charge) for cal. pulse candidetes
161
162 //first, identify cal. pulse
163 for (const auto& digit : digits) {
164
165 short slotId = digit.getModuleID();
166 short pixelId = digit.getPixelID();
167 short globalPixelId = (slotId - 1) * c_NPixelPerModule + pixelId - 1;
168 if (digit.getHitQuality() == TOPDigit::c_CalPulse) {
169 calPulsesMap[globalPixelId].push_back((hitInfo_t) { (float)digit.getTime(), (float)digit.getPulseHeight() });
170 }
171 }// for digit pair
172
173
174 //cal. pulse timing calculation
175 for (const auto& calPulse : calPulsesMap) {
176
177 short globalPixelId = calPulse.first;
178 short globalAsicId = globalPixelId / c_NChannelPerAsic;
179 std::vector<hitInfo_t> vec = calPulse.second;
180 double calPulseTiming = 9999999;
181 unsigned nCalPulseCandidates = vec.size();
182 if (!m_useDoublePulse) { //try to find the first cal. pulse in case that both of double cal. pulses are not required (choose the largest hit)
183 double maxPulseHeight = -1;
184 double maxHeightTiming = 9999999;
185 for (unsigned iVec = 0 ; iVec < nCalPulseCandidates ; iVec++) {
186 if (maxPulseHeight < vec[iVec].m_height) {
187 maxPulseHeight = vec[iVec].m_height;
188 maxHeightTiming = vec[iVec].m_time;
189 }
190 }
191 if (maxPulseHeight > m_calibrationPulseThreshold1)
192 calPulseTiming = maxHeightTiming;
193 } else { //try to identify both of double cal. pulses when specified by option (default)
194 for (unsigned iVec = 0 ; iVec < nCalPulseCandidates ; iVec++) {
195 for (unsigned jVec = 0 ; jVec < nCalPulseCandidates ; jVec++) {
196
197 if (iVec == jVec || vec[iVec].m_time > vec[jVec].m_time) continue;
198 if (vec[iVec].m_height > m_calibrationPulseThreshold1
199 && vec[jVec].m_height > m_calibrationPulseThreshold2
200 && TMath::Abs(vec[jVec].m_time - vec[iVec].m_time - m_calibrationPulseInterval) < m_calibrationPulseIntervalRange) {
201
202 //in case multiple candidates of double cal. pulses are found,
203 //choose a pair with the earliest 1st cal. pulse timing
204 if (refTimingMap.count(globalAsicId) == 0 || refTimingMap[globalAsicId] > vec[iVec].m_time)
205 calPulseTiming = vec[iVec].m_time;
206 }
207 }//for(jVec)
208 }//for(iVec)
209 }//if(m_useDoublePulse)
210
211 if (calPulseTiming < 9999998) { //when cal. pulse(s) are properly identified
212 refTimingMap[globalAsicId] = calPulseTiming;
213 m_nCalPulseHistogram->Fill(globalAsicId);
214 }
215 }//for(pair)
216
217 //calculate hit timing with respect to cal. pulse and fill hit info. histogram
218 for (const auto& digit : digits) {
219 short slotId = digit.getModuleID();
220 short pixelId = digit.getPixelID();
221 short globalPixelId = (slotId - 1) * c_NPixelPerModule + pixelId - 1;
222 short globalAsicId = globalPixelId / c_NChannelPerAsic;
223
224 if (digit.getHitQuality() == TOPDigit::c_Junk
225 || digit.getHitQuality() == TOPDigit::c_CrossTalk) continue;
226
227 float hitTime = digit.getTime() - refTimingMap[globalAsicId];
228 float pulseHeight = digit.getPulseHeight();
229 float Integral = digit.getIntegral();
230 short windowNumberOfHit = (short)(digit.getFirstWindow()) + (short)(digit.getRawTime() / 64);
231 if (m_windowSelect && windowNumberOfHit % 2 != (m_windowSelect - 1)) continue;
232
234 m_TimeHeightHistogramForHitRate[globalPixelId]->Fill(hitTime, pulseHeight);
235 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
236 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
237 continue;
238 }
239
240 if (!digit.isSecondaryChargeShare()) {
241 m_TimeHeightHistogramForHitRate[globalPixelId]->Fill(hitTime, pulseHeight);
243 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
244 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
245 continue;
246 }
247 }
248
249 if (digit.isSecondaryChargeShare() || digit.isPrimaryChargeShare()) continue;
250 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
251 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
252 }
253
254 }
255
257 {
258 }
259
261 {
262 }
263
265} // end Belle2 namespace
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
bool m_useDoublePulse
set true when you require both of double calibration pulses for reference timing
float m_calibrationPulseInterval
nominal DeltaT value (time interval of two calibration signals) in a unit of ns
float m_calibrationPulseIntervalRange
tolerable shift of DeltaT from its nominal before calibration in a unit of ns
float m_calibrationPulseThreshold2
minimum pulse height for the secon calibration pulse to be qualified as calibration signals
std::vector< int > m_timeHistogramBinning
histogram binning of hit timing distribution, in the order of number of bins, lower limit,...
TH2F * m_TimeHeightHistogramForFit[c_NPixelPerModule *c_NModule]
array of histogram pointer to 2D histogram of hit timing vs pulse height distribution for each pixel ...
std::vector< int > m_chargeHistogramBinning
histogram binning of pulse height distribution, in the order of number of bins, lower limit,...
bool m_includePrimaryChargeShare
set true when you require without chargeshare cut for making 2D histogram
TH2F * m_TimeHeightHistogramForHitRate[c_NPixelPerModule *c_NModule]
array of histogram pointer to 2D histogram of hit timing vs pulse height distribution for each pixel ...
TH2F * m_TimeIntegralHistogramForFit[c_NPixelPerModule *c_NModule]
array of histogram pointer to 2D histogram of hit timing vs integral distribution for each pixel (all...
float m_calibrationPulseThreshold1
minimum pulse height for the first calibration pulse to be qualified as calibration signals
bool m_includeAllChargeShare
set true when you require without chargeshare cut for making 2D histogram
int m_windowSelect
select window number is [All=0, Odd=2, Even=1]
TH1F * m_nCalPulseHistogram
histogram to store the number of events with calibration pulse(s) identified for each asic (1,...
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void endRun() override
End-of-run action.
virtual void terminate() override
Termination action.
virtual void beginRun() override
Called when entering a new run.
virtual void defineHisto() override
create timing-height 2D histograms for all 8192 pixels
Abstract base class for different kinds of events.
STL namespace.
structure to hold hit information, used in double cal.