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 {
32
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 true when you require without primary chargeshare cut for making 2D histogram",
86 (bool)false);
87 addParam("includeAllChargeShare", m_includeAllChargeShare,
88 "set true when you require without all chargeshare cut for making 2D histogram",
89 (bool)false);
90 }
91
92
94 {
95 REG_HISTOGRAM;
96
98 digits.isRequired();
99 }
100
102 {
103 for (int iPixel = 0 ; iPixel < c_NPixelPerModule * c_NModule ; iPixel++) {
104
105 short slotId = (iPixel / c_NPixelPerModule) + 1;
106 short pixelId = (iPixel % c_NPixelPerModule) + 1;
107 short pmtId = ((pixelId - 1) % c_NPixelPerRow) / c_NChannelPerPMTRow + c_NPMTPerRow * ((pixelId - 1) / (c_NPixelPerModule / 2)) + 1;
108 short pmtChId = (pixelId - 1) % c_NChannelPerPMTRow + c_NChannelPerPMTRow * ((pixelId - 1) %
109 (c_NPixelPerModule / 2) / c_NPixelPerRow) + 1;
110
111 std::ostringstream pixelstr;
112 pixelstr << "s" << std::setw(2) << std::setfill('0') << slotId << "_PMT"
113 << std::setw(2) << std::setfill('0') << pmtId
114 << "_" << std::setw(2) << std::setfill('0') << pmtChId;
115
116 std::ostringstream hnameForgain;
117 hnameForgain << "hTimeHeight_gain_" << pixelstr.str();
118 std::ostringstream htitleForgain;
119 htitleForgain << "2D distribution of hit timing and pulse height for Gain" << pixelstr.str();
120 m_TimeHeightHistogramForFit[iPixel] = new TH2F(hnameForgain.str().c_str(), htitleForgain.str().c_str(),
124
125 std::ostringstream hnameForIntegral;
126 hnameForIntegral << "hTimeIntegral_gain_" << pixelstr.str();
127 std::ostringstream htitleForIntegral;
128 htitleForIntegral << "2D distribution of hit timing and integral for Gain" << pixelstr.str();
129 m_TimeIntegralHistogramForFit[iPixel] = new TH2F(hnameForIntegral.str().c_str(), htitleForIntegral.str().c_str(),
133
134 std::ostringstream hnameForeff;
135 hnameForeff << "hTimeHeight_efficiency_" << pixelstr.str();
136 std::ostringstream htitleForeff;
137 htitleForeff << "2D distribution of hit timing and pulse height for efficiency" << pixelstr.str();
138 m_TimeHeightHistogramForHitRate[iPixel] = new TH2F(hnameForeff.str().c_str(), htitleForeff.str().c_str(),
142 }
143
144 const short nAsic = c_NPixelPerModule / c_NChannelPerAsic * c_NChannelPerPMT;
145 m_nCalPulseHistogram = new TH1F("hNCalPulse", "number of calibration pulses identificed for each asic",
146 nAsic, -0.5, nAsic - 0.5);
147 }
148
150 {
152
153 std::map<short, float> refTimingMap;//map of pixel id and time
154 std::map<short, std::vector<hitInfo_t> >
155 calPulsesMap;//map of pixel id and hit information (timing and pulse charge) for cal. pulse candidetes
156
157 //first, identify cal. pulse
158 for (const auto& digit : digits) {
159
160 short slotId = digit.getModuleID();
161 short pixelId = digit.getPixelID();
162 short globalPixelId = (slotId - 1) * c_NPixelPerModule + pixelId - 1;
163 if (digit.getHitQuality() == TOPDigit::c_CalPulse) {
164 calPulsesMap[globalPixelId].push_back((hitInfo_t) { (float)digit.getTime(), (float)digit.getPulseHeight() });
165 }
166 }// for digit pair
167
168
169 //cal. pulse timing calculation
170 for (const auto& calPulse : calPulsesMap) {
171
172 short globalPixelId = calPulse.first;
173 short globalAsicId = globalPixelId / c_NChannelPerAsic;
174 std::vector<hitInfo_t> vec = calPulse.second;
175 double calPulseTiming = 9999999;
176 unsigned nCalPulseCandidates = vec.size();
177 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)
178 double maxPulseHeight = -1;
179 double maxHeightTiming = 9999999;
180 for (unsigned iVec = 0 ; iVec < nCalPulseCandidates ; iVec++) {
181 if (maxPulseHeight < vec[iVec].m_height) {
182 maxPulseHeight = vec[iVec].m_height;
183 maxHeightTiming = vec[iVec].m_time;
184 }
185 }
186 if (maxPulseHeight > m_calibrationPulseThreshold1)
187 calPulseTiming = maxHeightTiming;
188 } else { //try to identify both of double cal. pulses when specified by option (default)
189 for (unsigned iVec = 0 ; iVec < nCalPulseCandidates ; iVec++) {
190 for (unsigned jVec = 0 ; jVec < nCalPulseCandidates ; jVec++) {
191
192 if (iVec == jVec || vec[iVec].m_time > vec[jVec].m_time) continue;
193 if (vec[iVec].m_height > m_calibrationPulseThreshold1
194 && vec[jVec].m_height > m_calibrationPulseThreshold2
195 && TMath::Abs(vec[jVec].m_time - vec[iVec].m_time - m_calibrationPulseInterval) < m_calibrationPulseIntervalRange) {
196
197 //in case multiple candidates of double cal. pulses are found,
198 //choose a pair with the earliest 1st cal. pulse timing
199 if (refTimingMap.count(globalAsicId) == 0 || refTimingMap[globalAsicId] > vec[iVec].m_time)
200 calPulseTiming = vec[iVec].m_time;
201 }
202 }//for(jVec)
203 }//for(iVec)
204 }//if(m_useDoublePulse)
205
206 if (calPulseTiming < 9999998) { //when cal. pulse(s) are properly identified
207 refTimingMap[globalAsicId] = calPulseTiming;
208 m_nCalPulseHistogram->Fill(globalAsicId);
209 }
210 }//for(pair)
211
212 //calculate hit timing with respect to cal. pulse and fill hit info. histogram
213 for (const auto& digit : digits) {
214 short slotId = digit.getModuleID();
215 short pixelId = digit.getPixelID();
216 short globalPixelId = (slotId - 1) * c_NPixelPerModule + pixelId - 1;
217 short globalAsicId = globalPixelId / c_NChannelPerAsic;
218
219 if (digit.getHitQuality() == TOPDigit::c_Junk
220 || digit.getHitQuality() == TOPDigit::c_CrossTalk) continue;
221
222 float hitTime = digit.getTime() - refTimingMap[globalAsicId];
223 float pulseHeight = digit.getPulseHeight();
224 float Integral = digit.getIntegral();
225 short windowNumberOfHit = (short)(digit.getFirstWindow()) + (short)(digit.getRawTime() / 64);
226 if (m_windowSelect && windowNumberOfHit % 2 != (m_windowSelect - 1)) continue;
227
229 m_TimeHeightHistogramForHitRate[globalPixelId]->Fill(hitTime, pulseHeight);
230 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
231 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
232 continue;
233 }
234
235 if (!digit.isSecondaryChargeShare()) {
236 m_TimeHeightHistogramForHitRate[globalPixelId]->Fill(hitTime, pulseHeight);
238 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
239 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
240 continue;
241 }
242 }
243
244 if (digit.isSecondaryChargeShare() || digit.isPrimaryChargeShare()) continue;
245 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
246 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
247 }
248
249 }
250
252} // end Belle2 namespace
HistoModule()
Constructor.
Definition HistoModule.h:32
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 second 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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
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.