10 #include <top/modules/TOPGainEfficiencyMonitor/TOPLaserHitSelectorModule.h>
13 #include <framework/datastore/StoreArray.h>
16 #include <top/dataobjects/TOPDigit.h>
47 setDescription(
"TOP pixel-by-pixel gain analysis - first step : create hit timing-pulse charge histogram");
50 m_timeHistogramBinning.push_back(140);
51 m_timeHistogramBinning.push_back(-55);
52 m_timeHistogramBinning.push_back(15);
53 m_chargeHistogramBinning.push_back(125);
54 m_chargeHistogramBinning.push_back(-50);
55 m_chargeHistogramBinning.push_back(2450);
56 m_chargeHistogramBinning.push_back(150);
57 m_chargeHistogramBinning.push_back(-500);
58 m_chargeHistogramBinning.push_back(29500);
60 addParam(
"timeHistogramBinning", m_timeHistogramBinning,
61 "histogram binning (the number of bins, lower limit, upper limit) for hit timing distribution (should be integer value)",
62 m_timeHistogramBinning);
63 addParam(
"chargeHistogramBinning", m_chargeHistogramBinning,
64 "histogram binning (the number of bins, lower limit, upper limit) for pulse charge distribution (should be integer value)",
65 m_chargeHistogramBinning);
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.",
72 addParam(
"minHeightFirstCalPulse", m_calibrationPulseThreshold1,
"pulse height threshold for the first cal. pulse",
74 addParam(
"minHeightSecondCalPulse", m_calibrationPulseThreshold2,
"pulse height threshold for the second cal. pulse",
76 addParam(
"nominalDeltaT", m_calibrationPulseInterval,
"nominal DeltaT (time interval of the double calibration pulses) [ns]",
78 addParam(
"nominalDeltaTRange", m_calibrationPulseIntervalRange,
79 "acceptable DeltaT range from the nominal value before calibration [ns]",
81 addParam(
"windowSelect", m_windowSelect,
82 "select window number (All=0, Odd=2, Even=1)",
84 addParam(
"includePrimaryChargeShare", m_includePrimaryChargeShare,
85 "set ture when you require without primary chargeshare cut for making 2D histogram",
87 addParam(
"includeAllChargeShare", m_includeAllChargeShare,
88 "set ture when you require without all chargeshare cut for making 2D histogram",
92 TOPLaserHitSelectorModule::~TOPLaserHitSelectorModule() {}
94 void TOPLaserHitSelectorModule::initialize()
102 void TOPLaserHitSelectorModule::defineHisto()
104 for (
int iPixel = 0 ; iPixel < c_NPixelPerModule * c_NModule ; iPixel++) {
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;
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;
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(),
122 m_timeHistogramBinning[0], m_timeHistogramBinning[1], m_timeHistogramBinning[2],
123 m_chargeHistogramBinning[0], m_chargeHistogramBinning[1],
124 m_chargeHistogramBinning[2]);
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(),
131 m_timeHistogramBinning[0], m_timeHistogramBinning[1], m_timeHistogramBinning[2],
132 m_chargeHistogramBinning[3], m_chargeHistogramBinning[4],
133 m_chargeHistogramBinning[5]);
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(),
140 m_timeHistogramBinning[0], m_timeHistogramBinning[1], m_timeHistogramBinning[2],
141 m_chargeHistogramBinning[0], m_chargeHistogramBinning[1],
142 m_chargeHistogramBinning[2]);
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);
150 void TOPLaserHitSelectorModule::beginRun()
154 void TOPLaserHitSelectorModule::event()
158 std::map<short, float> refTimingMap;
159 std::map<short, std::vector<hitInfo_t> >
163 for (
const auto& digit : digits) {
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() });
175 for (
const auto& calPulse : calPulsesMap) {
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) {
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;
191 if (maxPulseHeight > m_calibrationPulseThreshold1)
192 calPulseTiming = maxHeightTiming;
194 for (
unsigned iVec = 0 ; iVec < nCalPulseCandidates ; iVec++) {
195 for (
unsigned jVec = 0 ; jVec < nCalPulseCandidates ; jVec++) {
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) {
204 if (refTimingMap.count(globalAsicId) == 0 || refTimingMap[globalAsicId] > vec[iVec].m_time)
205 calPulseTiming = vec[iVec].m_time;
211 if (calPulseTiming < 9999998) {
212 refTimingMap[globalAsicId] = calPulseTiming;
213 m_nCalPulseHistogram->Fill(globalAsicId);
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;
224 if (digit.getHitQuality() == TOPDigit::c_Junk
225 || digit.getHitQuality() == TOPDigit::c_CrossTalk)
continue;
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;
233 if (m_includeAllChargeShare) {
234 m_TimeHeightHistogramForHitRate[globalPixelId]->Fill(hitTime, pulseHeight);
235 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
236 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
240 if (!digit.isSecondaryChargeShare()) {
241 m_TimeHeightHistogramForHitRate[globalPixelId]->Fill(hitTime, pulseHeight);
242 if (m_includePrimaryChargeShare) {
243 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
244 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
249 if (digit.isSecondaryChargeShare() || digit.isPrimaryChargeShare())
continue;
250 m_TimeIntegralHistogramForFit[globalPixelId]->Fill(hitTime, Integral);
251 m_TimeHeightHistogramForFit[globalPixelId]->Fill(hitTime, pulseHeight);
256 void TOPLaserHitSelectorModule::endRun()
260 void TOPLaserHitSelectorModule::terminate()
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Accessor to arrays stored in the data store.
Module for pixel-by-pixel gain/efficiency analysis.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.
structure to hold hit information, used in double cal.