Belle II Software  release-05-01-25
TOPGainEfficiencyCalculatorModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2017 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Maeda Yosuke, Okuto Rikuya *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 // own include
12 #include <top/modules/TOPGainEfficiencyMonitor/TOPGainEfficiencyCalculatorModule.h>
13 
14 // standard libraries
15 #include <iostream>
16 #include <sstream>
17 #include <iomanip>
18 
19 // ROOT
20 #include <TROOT.h>
21 #include <TObject.h>
22 #include <TFile.h>
23 #include <TF1.h>
24 #include <TCanvas.h>
25 #include <TStyle.h>
26 #include <TAxis.h>
27 #include <TLine.h>
28 #include <TArrow.h>
29 #include <TLatex.h>
30 #include <TMath.h>
31 
32 namespace Belle2 {
38  //-----------------------------------------------------------------
39  // Register the Module
40  //-----------------------------------------------------------------
41 
42  REG_MODULE(TOPGainEfficiencyCalculator)
43 
44 
45  //-----------------------------------------------------------------
46  // Implementation
47  //-----------------------------------------------------------------
48 
50  {
51  // Set description()
52  setDescription("Calculate pixel gain and efficiency for a given PMT from 2D histogram of hit timing and pulse charge, "
53  "created by TOPLaserHitSelectorModule.");
54 
55  // Add parameters
56  addParam("inputFile", m_inputFile, "input file name containing 2D histogram", std::string(""));
57  addParam("outputPDFFile", m_outputPDFFile, "output PDF file to store plots", std::string(""));
58  addParam("targetSlotId", m_targetSlotId, "TOP module ID in slot number (1-16)", (short)0);
59  addParam("targetPmtId", m_targetPmtId, "PMT number (1-32)", (short)0);
60  addParam("targetPmtChId", m_targetPmtChId, "PMT channel number (1-16)", (short) - 1);
61  addParam("hvDiff", m_hvDiff, "HV difference from nominal HV value", (short)0);
62  addParam("fitHalfWidth", m_fitHalfWidth, "half fit width for direct laser hit peak in [ns] unit", (float)1.0);
63  addParam("threshold", m_threshold,
64  "pulse height (or integrated charge) threshold in fitting its distribution and calculating efficiency", (float)100.);
65  addParam("p0HeightIntegral", m_p0HeightIntegral,
66  "Parameter from p0 + x*p1 function fitting height-integral distribtion.", (float) - 50.0);
67  addParam("p1HeightIntegral", m_p1HeightIntegral,
68  "Parameter from p0 + x*p1 function fitting height-integral distribtion.", (float)6.0);
69  addParam("fracFit", m_fracFit, "fraction of events to be used in fitting. "
70  "An upper limit of a fit range is given to cover this fraction of events. "
71  "Set negative value to calculate the fraction to exclude only 10 events in tail.", (float)(-1)); //,(float)0.99);
72  addParam("initialP0", m_initialP0, "initial value of the fit parameter p0 divided by histogram entries."
73  "Set negative value to calculate from histogram inforamtion automatically.", (float)(0.0001)); //,(float)1e-6);
74  addParam("initialP1", m_initialP1, "initial value of the fit parameter p1."
75  "Set negative value to calculate from histogram inforamtion automatically.", (float)(1.0)); //,(float)1.0);
76  addParam("initialP2", m_initialP2, "initial value of the fit parameter p2."
77  "Set negative value to calculate from histogram inforamtion automatically.", (float)(1.0)); //,(float)1.0);
78  addParam("initialX0", m_initialX0, "initial value of the fit parameter x0 divided by histogram bin width."
79  "Set negative value to calculate from histogram inforamtion automatically.", (float)(100)); //, (float)100.);
80  addParam("pedestalSigma", m_pedestalSigma, "sigma of pedestal width", (float)10.);
81  addParam("fitoption", m_fitoption, "fit option likelihood: default chisquare: R", std::string("L"));
82 
83  }
84 
86 
88  {
89  REG_HISTOGRAM;
90  }
91 
93  {
94  m_tree = new TTree("tree", "TTree for gain/efficiency monitor summary");
95  m_branch[0].push_back(m_tree->Branch("slotId", &m_targetSlotId, "slotId/S"));
96  m_branch[0].push_back(m_tree->Branch("pmtId", &m_targetPmtId, "pmtId/S"));
97  m_branch[0].push_back(m_tree->Branch("pixelId", &m_pixelId, "pixelId/S"));
98  m_branch[0].push_back(m_tree->Branch("pmtChId", &m_pmtChId, "pmtChId/S"));
99  m_branch[0].push_back(m_tree->Branch("hvDiff", &m_hvDiff, "hvDiff/S"));
100  m_branch[0].push_back(m_tree->Branch("threshold", &m_threshold, "threshold/F"));
101  m_branch[0].push_back(m_tree->Branch("thresholdForIntegral", &m_thresholdForIntegral, "thresholdForIntegral/F"));
102 
103  m_branch[1].push_back(m_tree->Branch("nCalPulse", &m_nCalPulse, "nCalPulse/I"));
104  m_branch[1].push_back(m_tree->Branch("hitTimingForGain", &m_hitTiming, "hitTimingForGain/F"));
105  m_branch[1].push_back(m_tree->Branch("hitTimingSigmaForGain", &m_hitTimingSigma, "hitTimingSigmaForGain/F"));
106  m_branch[1].push_back(m_tree->Branch("nEntriesForGain", &m_nEntries, "nEntriesForGain/I"));
107  m_branch[1].push_back(m_tree->Branch("nOverflowEventsForGain", &m_nOverflowEvents, "nOverflowEventsForGain/I"));
108  m_branch[1].push_back(m_tree->Branch("meanPulseHeightForGain", &m_meanPulseHeight, "meanPulseHeightForGain/F"));
109  m_branch[1].push_back(m_tree->Branch("meanPulseHeightErrorForGain", &m_meanPulseHeightError, "meanPulseHeightErrorForGain/F"));
110  m_branch[1].push_back(m_tree->Branch("fitMaxForGain", &m_fitMax, "fitMaxForGain/F"));
111  m_branch[1].push_back(m_tree->Branch("gain", &m_gain, "gain/F"));
112  m_branch[1].push_back(m_tree->Branch("efficiency", &m_efficiency, "efficiency/F"));
113  m_branch[1].push_back(m_tree->Branch("p0ForGain", &m_p0, "p0ForGain/F"));
114  m_branch[1].push_back(m_tree->Branch("p1ForGain", &m_p1, "p1ForGain/F"));
115  m_branch[1].push_back(m_tree->Branch("p2ForGain", &m_p2, "p2ForGain/F"));
116  m_branch[1].push_back(m_tree->Branch("x0ForGain", &m_x0, "x0ForGain/F"));
117  m_branch[1].push_back(m_tree->Branch("p0ErrorForGain", &m_p0Error, "p0ErrorForGain/F"));
118  m_branch[1].push_back(m_tree->Branch("p1ErrorForGain", &m_p1Error, "p1ErrorForGain/F"));
119  m_branch[1].push_back(m_tree->Branch("p2ErrorForGain", &m_p2Error, "p2ErrorForGain/F"));
120  m_branch[1].push_back(m_tree->Branch("x0ErrorForGain", &m_x0Error, "x0ErrorForGain/F"));
121  m_branch[1].push_back(m_tree->Branch("chisquareForGain", &m_chisquare, "chisquareForGain/F"));
122  m_branch[1].push_back(m_tree->Branch("ndfForGain", &m_ndf, "ndfForGain/I"));
123  m_branch[1].push_back(m_tree->Branch("funcFullRangeIntegralForGain", &m_funcFullRangeIntegral, "funcFullRangeIntegralForGain/F"));
124  m_branch[1].push_back(m_tree->Branch("funcFitRangeIntegralForGain", &m_funcFitRangeIntegral, "funcFitRangeIntegralForGain/F"));
125  m_branch[1].push_back(m_tree->Branch("histoFitRangeIntegralForGain", &m_histoFitRangeIntegral, "histoFitRangeIntegralForGain/F"));
126  m_branch[1].push_back(m_tree->Branch("histoMeanAboveThreForGain", &m_histoMeanAboveThre, "histoMeanAboveThreForGain/F"));
127 
128  m_branch[2].push_back(m_tree->Branch("hitTimingForEff", &m_hitTiming, "hitTimingForEff/F"));
129  m_branch[2].push_back(m_tree->Branch("hitTimingSigmaForEff", &m_hitTimingSigma, "hitTimingSigmaForEff/F"));
130  m_branch[2].push_back(m_tree->Branch("nEntriesForEff", &m_nEntries, "nEntriesForEff/I"));
131  m_branch[2].push_back(m_tree->Branch("nOverflowEventsForEff", &m_nOverflowEvents, "nOverflowEventsForEff/I"));
132  m_branch[2].push_back(m_tree->Branch("meanPulseHeightForEff", &m_meanPulseHeight, "meanPulseHeightForEff/F"));
133  m_branch[2].push_back(m_tree->Branch("meanPulseHeightErrorForEff", &m_meanPulseHeightError, "meanPulseHeightErrorForEff/F"));
134  m_branch[2].push_back(m_tree->Branch("histoMeanAboveThreForEff", &m_histoMeanAboveThre, "histoMeanAboveThreForEff/F"));
135 
136  m_branch[3].push_back(m_tree->Branch("meanIntegral", &m_meanPulseHeight, "meanIntegral/F"));
137  m_branch[3].push_back(m_tree->Branch("meanIntegralError", &m_meanPulseHeightError, "meanIntegralError/F"));
138  m_branch[3].push_back(m_tree->Branch("nOverflowEventsUseIntegral", &m_nOverflowEvents, "nOverflowEventsUseIntegral/I"));
139  m_branch[3].push_back(m_tree->Branch("fitMaxUseIntegral", &m_fitMax, "fitMaxUseIntegral/F"));
140  m_branch[3].push_back(m_tree->Branch("gainUseIntegral", &m_gain, "gainUseIntegral/F"));
141  m_branch[3].push_back(m_tree->Branch("efficiencyUseIntegral", &m_efficiency, "efficiencyUseIntegral/F"));
142  m_branch[3].push_back(m_tree->Branch("p0UseIntegral", &m_p0, "p0UseIntegral/F"));
143  m_branch[3].push_back(m_tree->Branch("p1UseIntegral", &m_p1, "p1UseIntegral/F"));
144  m_branch[3].push_back(m_tree->Branch("p2UseIntegral", &m_p2, "p2UseIntegral/F"));
145  m_branch[3].push_back(m_tree->Branch("x0UseIntegral", &m_x0, "x0UseIntegral/F"));
146  m_branch[3].push_back(m_tree->Branch("p0ErrorUseIntegral", &m_p0Error, "p0ErrorUseIntegral/F"));
147  m_branch[3].push_back(m_tree->Branch("p1ErrorUseIntegral", &m_p1Error, "p1ErrorUseIntegral/F"));
148  m_branch[3].push_back(m_tree->Branch("p2ErrorUseIntegral", &m_p2Error, "p2ErrorUseIntegral/F"));
149  m_branch[3].push_back(m_tree->Branch("x0ErrorUseIntegral", &m_x0Error, "x0ErrorUseIntegral/F"));
150  m_branch[3].push_back(m_tree->Branch("chisquareUseIntegral", &m_chisquare, "chisquareUseIntegral/F"));
151  m_branch[3].push_back(m_tree->Branch("ndfUseIntegral", &m_ndf, "ndfUseIntegral/I"));
152  m_branch[3].push_back(m_tree->Branch("funcFullRangeIntegralUseIntegral", &m_funcFullRangeIntegral,
153  "funcFullRangeIntegralUseIntegral/F"));
154  m_branch[3].push_back(m_tree->Branch("funcFitRangeIntegralUseIntegral", &m_funcFitRangeIntegral,
155  "funcFitRangeIntegralUseIntegral/F"));
156  m_branch[3].push_back(m_tree->Branch("histoFitRangeIntegralUseIntegral", &m_histoFitRangeIntegral,
157  "histoFitRangeIntegralUseIntegral/F"));
158  m_branch[3].push_back(m_tree->Branch("histoMeanAboveThreUseIntegral", &m_histoMeanAboveThre, "histoMeanAboveThreUseIntegral/F"));
159 
160  }
161 
163  {
164  }
165 
167  {
168  }
169 
171  {
172  }
173 
174 
176  {
177  //first, check validity of input parameters
178  bool areGoodParameters = true;
179  if (m_targetSlotId < 1 || m_targetSlotId > c_NModule) {
180  B2ERROR("TOPGainEfficiencyCalculator : invalid slotID : " << m_targetSlotId);
181  areGoodParameters = false;
182  }
183  if (m_targetPmtId < 1 || m_targetPmtId > c_NPMTPerModule) {
184  B2ERROR("TOPGainEfficiencyCalculator : invalid PMT ID : " << m_targetPmtId);
185  areGoodParameters = false;
186  }
187  if (m_pedestalSigma < 1e-6) {
188  B2ERROR("TOPGainEfficiencyCalculator : pedestal sigma must be non-zero positive value");
189  areGoodParameters = false;
190  }
191 
192  //do not proceed to the main process unless all the input parameters are not valid
193  if (areGoodParameters) {
194  if (m_outputPDFFile.size() == 0) {
195  if (m_inputFile.rfind(".") == std::string::npos)
197  else
198  m_outputPDFFile = m_inputFile.substr(0, m_inputFile.rfind("."));
199  }
200 
201  //gain run using height distribution
202  LoadHistograms("Height_gain");
203  FitHistograms(c_LoadForFitHeight);
204  DrawResult("Height_gain", c_LoadForFitHeight);
205 
206  //efficiency run
207  LoadHistograms("Height_efficiency");
208  FitHistograms(c_LoadHitRateHeight);
209  DrawResult("Height_efficiency", c_LoadHitRateHeight);
210 
211  //gain run using integral distribution
212  LoadHistograms("Integral_gain");
213  FitHistograms(c_LoadForFitIntegral);
214  DrawResult("Integral_gain", c_LoadForFitIntegral);
215 
216  }
217  }
218 
219 
220  void TOPGainEfficiencyCalculatorModule::LoadHistograms(const std::string& histotype)
221  {
222 
223  TFile* f = new TFile(m_inputFile.c_str());
224  if (!f->IsOpen()) {
225  B2ERROR("TOPGainEfficiencyCalculator : fail to open input file \"" << m_inputFile << "\"");
226  return;
227  }
228 
229  for (int iHisto = 0 ; iHisto < c_NChannelPerPMT ; iHisto++) {
230  if (m_targetPmtChId != -1 && iHisto + 1 != m_targetPmtChId) continue;
231  std::ostringstream pixelstr;
232  pixelstr << histotype << "_"
233  << "s" << std::setw(2) << std::setfill('0') << m_targetSlotId << "_PMT"
234  << std::setw(2) << std::setfill('0') << m_targetPmtId
235  << "_" << std::setw(2) << std::setfill('0') << (iHisto + 1);
236  std::ostringstream hname;
237  hname << "hTime" << pixelstr.str();
238 
239  //first get 2D histogram from a given input (=an output file of TOPLaserHitSelector)
240  m_timeChargeHistogram[iHisto] = (TH2F*)f->Get(hname.str().c_str());
241  TH2F* h2D = m_timeChargeHistogram[iHisto];
242  if (!h2D) continue;
243 
244  //create a projection histogram along the x-axis and fit the distribution (hit timing) to get direct laser hit timing
245  std::ostringstream hnameProj[2];
246  hnameProj[0] << "hTime_" << pixelstr.str();
247  hnameProj[1] << "hCharge_" << pixelstr.str();
248  TH1D* hTime = (TH1D*)h2D->ProjectionX(hnameProj[0].str().c_str());
249  m_timeHistogram[iHisto] = hTime;
250  double peakTime = FindPeakForSmallerXThan(hTime, 0);
251  //double peakTime = hTime->GetXaxis()->GetBinCenter(hTime->GetMaximumBin());
252  double fitMin = peakTime - m_fitHalfWidth;
253  double fitMax = peakTime + m_fitHalfWidth;
254  TF1* funcLaser = new TF1(std::string(std::string("func_") + hnameProj[1].str()).c_str(),
255  "gaus(0)", fitMin, fitMax);
256  funcLaser->SetParameters(hTime->GetBinContent(hTime->GetXaxis()->FindBin(peakTime)), peakTime, m_fitHalfWidth);
257  funcLaser->SetParLimits(1, fitMin, fitMax);
258  hTime->Fit(funcLaser, "Q", "", fitMin, fitMax);
259  //if (funcLaser->GetNDF() < 1) continue;
260  m_funcForLaser[iHisto] = funcLaser;
261 
262  //if the fitting is successful, create y-projection histogram with timing cut
263  m_hitTiming = funcLaser->GetParameter(1);
264  int binNumMin = hTime->GetXaxis()->FindBin(m_hitTiming - 2 * m_fitHalfWidth);
265  int binNumMax = hTime->GetXaxis()->FindBin(m_hitTiming + 2 * m_fitHalfWidth);
266  TH1D* hCharge = (TH1D*)h2D->ProjectionY(hnameProj[1].str().c_str(),
267  binNumMin, binNumMax);
268  m_chargeHistogram[iHisto] = hCharge;
269  }
270 
271  m_nCalPulseHistogram = (TH1F*)f->Get("hNCalPulse");
273  B2WARNING("TOPGainEfficiencyCalculator : no histogram for the number of events with calibration pulses identified in the given input file");
275  return;
276  }
277 
278  void TOPGainEfficiencyCalculatorModule::FitHistograms(EHistogramType LoadHisto)
279  {
280  float threshold = m_threshold;
281  int globalAsicId = 0;
282  if (LoadHisto == c_LoadForFitIntegral || LoadHisto == c_LoadHitRateIntegral) threshold = m_thresholdForIntegral;
283 
284  for (int iHisto = 0 ; iHisto < c_NChannelPerPMT ; iHisto++) {
285  if (m_targetPmtChId != -1 && iHisto + 1 != m_targetPmtChId) continue;
286 
287  m_pixelId = ((m_targetPmtId - 1) % c_NPMTPerRow) * c_NChannelPerPMTRow
288  + ((m_targetPmtId - 1) / c_NPMTPerRow) * c_NPixelPerModule / 2
289  + (iHisto / c_NChannelPerPMTRow) * c_NPixelPerRow + (iHisto % c_NChannelPerPMTRow) + 1;
290  m_pmtChId = (iHisto + 1);
291  globalAsicId = ((m_targetSlotId - 1) * c_NPixelPerModule + (m_pixelId - 1)) / c_NChannelPerAsic;
292  if (LoadHisto == c_LoadForFitHeight) {
293  for (auto itr = m_branch[0].begin(); itr != m_branch[0].end(); ++itr) {
294  (*itr)->Fill();
295  }
296  }
297 
298  TH1D* hCharge = m_chargeHistogram[iHisto];
299  if (!hCharge) { DummyFillBranch(LoadHisto); continue;}
300 
301  std::cout << " ***** fitting charge distribution for " << hCharge->GetName() << " *****" << std::endl;
302  int nBins = hCharge->GetXaxis()->GetNbins();
303  double binWidth = hCharge->GetXaxis()->GetBinUpEdge(1) - hCharge->GetXaxis()->GetBinLowEdge(1);
304  double histoMax = hCharge->GetXaxis()->GetBinUpEdge(nBins);
305  m_fitMax = threshold;
306  double wholeIntegral = hCharge->Integral(0, hCharge->GetXaxis()->GetNbins() + 1);
307  double fitRangeFraction = (m_fracFit > 0 ? m_fracFit : 1. - 10. / wholeIntegral);
308  while (hCharge->Integral(0, hCharge->GetXaxis()->FindBin(m_fitMax - 0.01 * binWidth)) / wholeIntegral < fitRangeFraction)
309  m_fitMax += binWidth;
310  if (m_fitMax < threshold + c_NParameterGainFit * binWidth) {
311  B2WARNING("TOPGainEfficiencyCalculator : no enough entries for fitting at slot"
312  << std::setw(2) << std::setfill('0') << m_targetSlotId << ", PMT"
313  << std::setw(2) << std::setfill('0') << m_targetPmtId << ", Ch"
314  << std::setw(2) << std::setfill('0') << m_pmtChId);
315  DummyFillBranch(LoadHisto); continue;
316  }
317 
318  std::ostringstream fname;
319  fname << "func_" << (iHisto + 1);
320  TObject* object = gROOT->FindObject(fname.str().c_str());
321  if (object) delete object;
322  TF1* func = new TF1(fname.str().c_str(), TOPGainFunc, threshold, m_fitMax, c_NParameterGainFit);
323  double initGain = TMath::Max(hCharge->GetMean(), 26.1) - 25;
324  double initP1 = TMath::Min(4.0, TMath::Max(10000.*TMath::Power(initGain - 25, -2), 0.01));
325  double initP2 = TMath::Min(0.8 + 0.005 * TMath::Power(initP1, -3), 4.);
326  double initX0 = TMath::Max(initGain * 2 - 150, 10.);
327  //if (initP1 > initP2) initX0 = initX0 / 10.
328  double initP1overP2 = initP1 / initP2;
329  double initP0 = hCharge->GetBinContent(hCharge->GetMaximumBin())
330  / (TMath::Power(initP1overP2, initP1overP2) * TMath::Exp(-initP1overP2)) / binWidth;
331  if (m_initialX0 < 0)func->SetParameter(3, initX0);
332  else if (LoadHisto == c_LoadForFitHeight)
333  func->SetParameter(3, 150 + 0.5 * hCharge->GetMean());
334  else if (LoadHisto == c_LoadForFitIntegral)
335  func->SetParameter(3, 1000 + 0.5 * hCharge->GetMean());
336  if (m_initialP2 < 0)func->SetParameter(2, initP2);
337  else func->SetParameter(2, m_initialP2);
338  if (m_initialP1 < 0)func->SetParameter(1, initP1);
339  else func->SetParameter(1, m_initialP1);
340  if (m_initialP0 < 0)func->SetParameter(0, initP0);
341  else func->SetParameter(0, m_initialP0 * hCharge->GetEntries()*binWidth);
342 
343  func->FixParameter(4, 0);
344  func->FixParameter(5, m_pedestalSigma);
345  func->SetParName(0, "#it{p}_{0}");
346  func->SetParName(1, "#it{p}_{1}");
347  func->SetParName(2, "#it{p}_{2}");
348  func->SetParName(3, "#it{x}_{0}");
349  func->SetParName(4, "pedestal");
350  func->SetParName(5, "pedestal #sigma");
351  func->SetParLimits(0, 1e-8, 1e8);
352  func->SetParLimits(1, 1e-8, 10);
353  func->SetParLimits(2, 1e-8, 10);
354  func->SetParLimits(3, 1e-8, 1e8);
355  func->SetLineColor(2);
356  func->SetLineWidth(1);
357  TF1* funcFull = NULL;
358  if (LoadHisto == c_LoadForFitHeight or LoadHisto == c_LoadForFitIntegral) {
359  hCharge->Fit(func, m_fitoption.c_str(), "", threshold , m_fitMax);
360 
361  if (func->GetNDF() < 2) { DummyFillBranch(LoadHisto); continue;}
362 
363  double funcFullMax = histoMax * 2;
364  funcFull = new TF1((fname.str() + "_full").c_str(), TOPGainFunc, (-1)*func->GetParameter(5), funcFullMax, c_NParameterGainFit);
365  for (int iPar = 0 ; iPar < c_NParameterGainFit ; iPar++)
366  funcFull->SetParameter(iPar, func->GetParameter(iPar));
367  funcFull->SetLineColor(3);
368  funcFull->SetLineWidth(2);
369 
370  double totalWeight = 0;
371  double weightedIntegral = 0;
372  double x = (-1) * func->GetParameter(5);
373  while (x < funcFullMax) {
374  double funcVal = funcFull->Eval(x);
375  totalWeight += funcVal;
376  weightedIntegral += funcVal * x;
377  x += binWidth / 5.;
378  }
379 
380  //fill results to the output TTree
381  m_gain = weightedIntegral / totalWeight;
382  m_efficiency = funcFull->Integral(threshold, funcFullMax) / funcFull->Integral((-1) * func->GetParameter(5), funcFullMax);
383  m_p0 = func->GetParameter(0);
384  m_p1 = func->GetParameter(1);
385  m_p2 = func->GetParameter(2);
386  m_x0 = func->GetParameter(3);
387  m_p0Error = func->GetParError(0);
388  m_p1Error = func->GetParError(1);
389  m_p2Error = func->GetParError(2);
390  m_x0Error = func->GetParError(3);
391  m_chisquare = func->GetChisquare();
392  m_ndf = func->GetNDF();
393  m_funcFullRangeIntegral = funcFull->Integral((-1) * func->GetParameter(5), funcFullMax) / binWidth;
394  m_funcFitRangeIntegral = funcFull->Integral(threshold, m_fitMax) / binWidth;
395  } else std::cout << "*****fitting is skipped***** " << std::endl;
396  int threBin = hCharge->GetXaxis()->FindBin(threshold + 0.01 * binWidth);
397  int fitMaxBin = hCharge->GetXaxis()->FindBin(m_fitMax - 0.01 * binWidth);
398  m_histoFitRangeIntegral = hCharge->Integral(threBin, fitMaxBin);
399 
401  for (int iBin = threBin ; iBin < nBins + 1 ; iBin++) {
402  m_histoMeanAboveThre += (hCharge->GetBinContent(iBin) * hCharge->GetXaxis()->GetBinCenter(iBin));
403  }
404  m_histoMeanAboveThre /= hCharge->Integral(threBin, nBins);
405  m_nEntries = hCharge->GetEntries();
406  m_nCalPulse = (m_nCalPulseHistogram ? m_nCalPulseHistogram->GetBinContent(globalAsicId + 1) : -1);
407  m_nOverflowEvents = TMath::FloorNint(hCharge->GetBinContent(nBins + 1));
408  m_meanPulseHeight = hCharge->GetMean();
409  m_meanPulseHeightError = hCharge->GetMeanError();
410  m_hitTiming = 0;
411  m_hitTimingSigma = -1;
412 
413  TF1* funcLaser = m_funcForLaser[iHisto];
414  if (m_timeHistogram[iHisto] && funcLaser) {
415  m_hitTiming = funcLaser->GetParameter(1);
416  m_hitTimingSigma = funcLaser->GetParameter(2);
417  }
418 
419  m_funcForFitRange[iHisto] = func;
420  m_funcForFullRange[iHisto] = funcFull;
421 
422  for (auto itr = m_branch[LoadHisto].begin(); itr != m_branch[LoadHisto].end(); ++itr) {
423  (*itr)->Fill();
424  }
425 
426  std::cout << std::endl;
427  }
428  if (m_targetPmtChId == -1) m_tree->SetEntries(c_NChannelPerPMT);
429  else m_tree->SetEntries(1);
430 
431  return;
432  }
433 
434  void TOPGainEfficiencyCalculatorModule::DummyFillBranch(EHistogramType LoadHisto)
435  {
436  m_fitMax = -1;
437  m_hitTiming = -1;
438  m_hitTimingSigma = -1;
439  m_nEntries = -1;
440  m_nCalPulse = -1;
441  m_nOverflowEvents = -1;
444  m_gain = -1;
445  m_efficiency = -1;
446  m_p0 = -1;
447  m_p1 = -1;
448  m_p2 = -1;
449  m_x0 = -1;
450  m_p0Error = -1;
451  m_p1Error = -1;
452  m_p2Error = -1;
453  m_x0Error = -1;
454  m_chisquare = -1;
455  m_ndf = -1;
460 
461  for (auto itr = m_branch[LoadHisto].begin(); itr != m_branch[LoadHisto].end(); ++itr) {
462  (*itr)->Fill();
463  }
464  }
465 
466 
467  void TOPGainEfficiencyCalculatorModule::DrawResult(const std::string& histotype, EHistogramType LoadHisto)
468  {
469  std::ostringstream pdfFilename;
470  pdfFilename << m_outputPDFFile << "_" << histotype;
471  if (m_targetPmtChId != -1)
472  pdfFilename << "_" << "ch" << std::setw(2) << std::setfill('0') << m_targetPmtChId;
473  pdfFilename << ".pdf";
474 
475  gStyle->SetFrameFillStyle(0);
476  gStyle->SetFillStyle(0);
477  gStyle->SetStatStyle(0);
478  gStyle->SetOptStat(112210);
479  gStyle->SetOptFit(1110);
480  TCanvas* canvas = new TCanvas();
481  canvas->SetFillStyle(0);
482  canvas->Print((pdfFilename.str() + "[").c_str());
483 
484  TLine* line = new TLine();
485  line->SetLineWidth(1);
486  line->SetLineStyle(1);
487  line->SetLineColor(4);
488  TArrow* arrow = new TArrow();
489  arrow->SetLineWidth(1);
490  arrow->SetLineStyle(1);
491  arrow->SetLineColor(3);
492  TLatex* latex = new TLatex();
493  latex->SetNDC();
494  latex->SetTextFont(22);
495  latex->SetTextSize(0.05);
496  latex->SetTextAlign(32);
497  TObject* object;
498 
499  float threshold;
500  if (LoadHisto == c_LoadForFitIntegral) threshold = m_thresholdForIntegral;
501  else threshold = m_threshold;
502 
503  for (int iHisto = 0 ; iHisto < c_NChannelPerPMT ; iHisto++) {
504 
505  if ((iHisto % c_NChannelPerPage) == 0) {
506  canvas->Clear();
507  canvas->Divide(c_NPlotsPerChannel, c_NChannelPerPage);
508  }
509 
510  //2D (time vs pulse charge) histogram
511  canvas->cd(3 * (iHisto % c_NChannelPerPage) + 1);
512  gPad->SetFrameFillStyle(0);
513  gPad->SetFillStyle(0);
514  TH2F* h2D = m_timeChargeHistogram[iHisto];
515  if (h2D) {
516  h2D->Draw("colz");
517  h2D->GetXaxis()->SetTitle("hit timing [ns]");
518  h2D->GetYaxis()->SetTitle("pulse charge [ADC count]");
519  }
520 
521  //timing histogram
522  canvas->cd(c_NPlotsPerChannel * (iHisto % c_NChannelPerPage) + 2);
523  gPad->SetFrameFillStyle(0);
524  gPad->SetFillStyle(0);
525  TH1D* hTime = m_timeHistogram[iHisto];
526  if (hTime) {
527  gPad->SetLogy();
528  hTime->Draw();
529  hTime->SetLineColor(1);
530  hTime->GetXaxis()->SetTitle("hit timing [ns]");
531  float binWidth = hTime->GetXaxis()->GetBinUpEdge(1) - hTime->GetXaxis()->GetBinLowEdge(1);
532  std::ostringstream ytitle;
533  ytitle << "Entries [/(" << binWidth << " ns)]";
534  hTime->GetYaxis()->SetTitle(ytitle.str().c_str());
535 
536  TF1* funcLaser = m_funcForLaser[iHisto];
537  if (funcLaser) {
538  double charge = funcLaser->GetParameter(0);
539  double peakTime = funcLaser->GetParameter(1);
540  float xMin = hTime->GetXaxis()->GetBinLowEdge(hTime->GetXaxis()->FindBin(peakTime - 2 * m_fitHalfWidth));
541  float xMax = hTime->GetXaxis()->GetBinUpEdge(hTime->GetXaxis()->FindBin(peakTime + 2 * m_fitHalfWidth));
542  line->DrawLine(xMin, 0.5, xMin, charge * 2.);
543  line->DrawLine(xMax, 0.5, xMax, charge * 2.);
544  arrow->DrawArrow(xMin, charge * 1.5, xMax, charge * 1.5, 0.01, "<>");
545  }
546  }
547 
548  //charge histogram with fit result (after timing cut)
549  canvas->cd(c_NPlotsPerChannel * (iHisto % c_NChannelPerPage) + 3);
550  gPad->SetFrameFillStyle(0);
551  gPad->SetFillStyle(0);
552  TH1D* hCharge = m_chargeHistogram[iHisto];
553  if (hCharge) {
554  gPad->SetLogy();
555  hCharge->Draw();
556  hCharge->SetLineColor(1);
557  hCharge->GetXaxis()->SetTitle("charge [ADC counts]");
558  float binWidth = hCharge->GetXaxis()->GetBinUpEdge(1) - hCharge->GetXaxis()->GetBinLowEdge(1);
559  std::ostringstream ytitle;
560  ytitle << "Entries [/(" << binWidth << " ADC counts)]";
561  hCharge->GetYaxis()->SetTitle(ytitle.str().c_str());
562 
563  if (m_funcForFitRange[iHisto] && m_funcForFullRange[iHisto]) {
564  m_funcForFullRange[iHisto]->Draw("same");
565  m_funcForFitRange[iHisto]->Draw("same");
566  double charge = hCharge->GetBinContent(hCharge->GetMaximumBin());
567  line->DrawLine(threshold, 0.5, threshold, charge * 2.);
568 
569  if ((object = gROOT->FindObject("dummy"))) delete object;
570  std::ostringstream cut;
571  cut << "pmtChId==" << (iHisto + 1);
572  long nEntries = 0;
573  std::ostringstream summarystr[2];
574  if (LoadHisto == c_LoadForFitHeight) {
575  nEntries = m_tree->Project("dummy", "gain:efficiency", cut.str().c_str());
576  if (nEntries == 1) {
577  summarystr[0] << "gain = " << std::setiosflags(std::ios::fixed) << std::setprecision(1)
578  << m_tree->GetV1()[0];
579  latex->DrawLatex(0.875, 0.34, summarystr[0].str().c_str());
580  summarystr[1] << "efficiency = " << std::setiosflags(std::ios::fixed) << std::setprecision(1)
581  << (m_tree->GetV2()[0] * 100) << " %";
582  latex->DrawLatex(0.875, 0.29, summarystr[1].str().c_str());
583  }
584  } else if (LoadHisto == c_LoadForFitIntegral) {
585  nEntries = m_tree->Project("dummy", "gainUseIntegral:efficiencyUseIntegral", cut.str().c_str());
586  if (nEntries == 1) {
587  summarystr[0] << "gain = " << std::setiosflags(std::ios::fixed) << std::setprecision(1)
588  << m_tree->GetV1()[0];
589  latex->DrawLatex(0.875, 0.34, summarystr[0].str().c_str());
590  summarystr[1] << "efficiency = " << std::setiosflags(std::ios::fixed) << std::setprecision(1)
591  << (m_tree->GetV2()[0] * 100) << " %";
592  latex->DrawLatex(0.875, 0.29, summarystr[1].str().c_str());
593  }
594  }
595 
596  if (nEntries > 1) {
597  B2WARNING("TOPGainEfficiencyCalculator : mutliple entries with the same channel ID ("
598  << m_pmtChId << ") in the output TTree");
599  }
600  }
601  }
602 
603  if (((iHisto + 1) % c_NChannelPerPage) == 0)
604  canvas->Print(pdfFilename.str().c_str());
605  }
606  for (int iHisto = (c_NChannelPerPMT - 1) % c_NChannelPerPage + 1 ; iHisto < c_NChannelPerPage ; iHisto++) {
607  for (int iPad = 0 ; iPad < c_NPlotsPerChannel ; iPad++) {
608  canvas->cd(c_NPlotsPerChannel * (iHisto % c_NChannelPerPage) + iPad + 1);
609  gPad->SetFrameFillStyle(0);
610  gPad->SetFillStyle(0);
611  }
612  if (((iHisto + 1) % c_NChannelPerPage) == 0)
613  canvas->Print(pdfFilename.str().c_str());
614  }
615 
616  canvas->Print((pdfFilename.str() + "]").c_str());
617 
618  delete latex;
619  delete arrow;
620  delete line;
621  delete canvas;
622  if ((object = gROOT->FindObject("dummy"))) delete object;
623 
624  return;
625  }
626 
627  double TOPGainEfficiencyCalculatorModule::TOPGainFunc(double* var, double* par)
628  {
629 
630  double pedestal = par[4];
631  double pedestalWidth = par[5];
632  double x = (var[0] - pedestal);
633 
634  double output = 0;
635  double step = pedestalWidth / 100.;
636  double t = TMath::Max(step, x - pedestalWidth * 10);
637  //double t = TMath::Max( g_xStep, x ); //for test
638  while (t < x + pedestalWidth * 10) {
639  output += TMath::Gaus(x, t, pedestalWidth) * TMath::Power(t / par[3], par[1])
640  * TMath::Exp((-1) * TMath::Power(t / par[3], par[2]));
641  t += step;
642  }
643  // TF1* func = new TF1( "func", "[0]*pow(x-[4],[1])*exp(-pow(x-[4],[2])/[3])",
644 
645  return par[0] * output * step;
646  }
647 
648  double TOPGainEfficiencyCalculatorModule::FindPeakForSmallerXThan(TH1* histo, double xmax)
649  {
650  double histo_xmax = histo->GetXaxis()->GetBinUpEdge(histo->GetXaxis()->GetNbins() - 1);
651  if (xmax > histo_xmax) xmax = histo_xmax;
652 
653  int iBin = 1;
654  double peakPos = histo->GetXaxis()->GetBinCenter(iBin);
655  double peakCharge = histo->GetBinContent(iBin);
656  while (true) {
657  iBin++;
658  double x = histo->GetXaxis()->GetBinCenter(iBin);
659  if (x > xmax) break;
660 
661  double binEntry = histo->GetBinContent(iBin);
662  if (binEntry > peakCharge) {
663  peakPos = x;
664  peakCharge = binEntry;
665  }
666  }
667 
668  return peakPos;
669  }
670 
671 
672 
674 } // end Belle2 namespace
Belle2::EvtPDLUtil::charge
double charge(int pdgCode)
Returns electric charge of a particle with given pdg code.
Definition: EvtPDLUtil.cc:46
Belle2::TOPGainEfficiencyCalculatorModule::m_p2
float m_p2
fit result of p2
Definition: TOPGainEfficiencyCalculatorModule.h:198
Belle2::TOPGainEfficiencyCalculatorModule::m_chargeHistogram
TH1D * m_chargeHistogram[c_NChannelPerPMT]
pulse charge distribution, extracted from m_timeChargeHistogram as a projection along its y-axis with...
Definition: TOPGainEfficiencyCalculatorModule.h:152
Belle2::TOPGainEfficiencyCalculatorModule::m_funcForLaser
TF1 * m_funcForLaser[c_NChannelPerPMT]
array of TF1 pointer to store fit function for hit timing distribution
Definition: TOPGainEfficiencyCalculatorModule.h:155
Belle2::TOPGainEfficiencyCalculatorModule::m_threshold
float m_threshold
pulse charge threshold, which defines lower limit of fit region and efficiency calculation
Definition: TOPGainEfficiencyCalculatorModule.h:172
Belle2::TOPGainEfficiencyCalculatorModule::m_funcFullRangeIntegral
float m_funcFullRangeIntegral
integral of fit function for its full range
Definition: TOPGainEfficiencyCalculatorModule.h:206
Belle2::TOPGainEfficiencyCalculatorModule::DrawResult
void DrawResult(const std::string &histotype, EHistogramType LoadHisto)
Draw results of gain/efficiency calculation for each channel to a given output file.
Definition: TOPGainEfficiencyCalculatorModule.cc:475
Belle2::TOPGainEfficiencyCalculatorModule::m_targetPmtChId
short m_targetPmtChId
PMT channel ID.
Definition: TOPGainEfficiencyCalculatorModule.h:168
Belle2::TOPGainEfficiencyCalculatorModule::m_hitTimingSigma
float m_hitTimingSigma
Gaussian fit sigma for a peak of laser direct photons in hit timing distribution.
Definition: TOPGainEfficiencyCalculatorModule.h:188
Belle2::TOPGainEfficiencyCalculatorModule::m_thresholdForIntegral
float m_thresholdForIntegral
pulse integral threshold, which defines lower limit of fit region and efficiency calculation
Definition: TOPGainEfficiencyCalculatorModule.h:173
Belle2::TOPGainEfficiencyCalculatorModule::m_hitTiming
float m_hitTiming
timing of laser direct photon hits, given by Gaussian fit mean
Definition: TOPGainEfficiencyCalculatorModule.h:187
Belle2::TOPGainEfficiencyCalculatorModule::m_targetPmtId
short m_targetPmtId
PMT ID.
Definition: TOPGainEfficiencyCalculatorModule.h:167
Belle2::TOPGainEfficiencyCalculatorModule::m_nCalPulseHistogram
TH1F * m_nCalPulseHistogram
histogram to store the number of events with calibration pulse(s) identified for each asic (1,...
Definition: TOPGainEfficiencyCalculatorModule.h:153
Belle2::Module::setDescription
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:216
Belle2::TOPGainEfficiencyCalculatorModule::m_initialP1
float m_initialP1
initial value of the fit parameter p1
Definition: TOPGainEfficiencyCalculatorModule.h:180
Belle2::TOPGainEfficiencyCalculatorModule::m_meanPulseHeight
float m_meanPulseHeight
histogram mean of pulse height distribution
Definition: TOPGainEfficiencyCalculatorModule.h:192
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::TOPGainEfficiencyCalculatorModule::m_outputPDFFile
std::string m_outputPDFFile
output PDF file to store plots of 2D histogram, timing, and charge distribution for each channel
Definition: TOPGainEfficiencyCalculatorModule.h:160
Belle2::TOPGainEfficiencyCalculatorModule::m_fitHalfWidth
float m_fitHalfWidth
half fit width for direct laser hit peak in [ns] unit
Definition: TOPGainEfficiencyCalculatorModule.h:171
Belle2::TOPGainEfficiencyCalculatorModule::m_fracFit
float m_fracFit
fraction of events which are covered by an area [0,m_fitMax]
Definition: TOPGainEfficiencyCalculatorModule.h:178
Belle2::TOPGainEfficiencyCalculatorModule::event
virtual void event() override
This will be empty as the all the processes are done in beginRun() function thus input file can be a ...
Definition: TOPGainEfficiencyCalculatorModule.cc:174
Belle2::TOPGainEfficiencyCalculatorModule::terminate
virtual void terminate() override
Termination action.
Definition: TOPGainEfficiencyCalculatorModule.cc:183
Belle2::TOPGainEfficiencyCalculatorModule::m_p1Error
float m_p1Error
fit error of p1
Definition: TOPGainEfficiencyCalculatorModule.h:201
Belle2::TOPGainEfficiencyCalculatorModule::m_nEntries
int m_nEntries
entries of pulse charge distribution
Definition: TOPGainEfficiencyCalculatorModule.h:189
Belle2::TOPGainEfficiencyCalculatorModule::m_p1
float m_p1
fit result of p1
Definition: TOPGainEfficiencyCalculatorModule.h:197
Belle2::TOPGainEfficiencyCalculatorModule::m_histoMeanAboveThre
float m_histoMeanAboveThre
mean of histogram above threshold, ignore overflow bin
Definition: TOPGainEfficiencyCalculatorModule.h:209
Belle2::TOPGainEfficiencyCalculatorModule::m_p1HeightIntegral
float m_p1HeightIntegral
Parameter from p0 + x*p1 function that fits height-integral distribution.
Definition: TOPGainEfficiencyCalculatorModule.h:176
Belle2::TOPGainEfficiencyCalculatorModule::beginRun
virtual void beginRun() override
The main processes, fitting charge distribution and calculating gain/efficiency, are done in this fun...
Definition: TOPGainEfficiencyCalculatorModule.cc:170
Belle2::TOPGainEfficiencyCalculatorModule::m_tree
TTree * m_tree
ntuple to store summary
Definition: TOPGainEfficiencyCalculatorModule.h:147
Belle2::TOPGainEfficiencyCalculatorModule::m_nCalPulse
int m_nCalPulse
the number of events with calibration pulse(s) identified
Definition: TOPGainEfficiencyCalculatorModule.h:190
Belle2::TOPGainEfficiencyCalculatorModule::m_efficiency
float m_efficiency
calculated efficiency from fitting of pulse charge distribution
Definition: TOPGainEfficiencyCalculatorModule.h:195
Belle2::TOPGainEfficiencyCalculatorModule::m_branch
std::vector< TBranch * > m_branch[4]
ntuple to store summary of gain using height distribution.
Definition: TOPGainEfficiencyCalculatorModule.h:148
Belle2::TOPGainEfficiencyCalculatorModule::m_inputFile
std::string m_inputFile
input file containing timing vs charge 2D histograms (output of TOPLaserHitSelector)
Definition: TOPGainEfficiencyCalculatorModule.h:159
Belle2::TOPGainEfficiencyCalculatorModule::m_nOverflowEvents
int m_nOverflowEvents
the number of events outside histogram range
Definition: TOPGainEfficiencyCalculatorModule.h:191
Belle2::TOPGainEfficiencyCalculatorModule::m_initialX0
float m_initialX0
initial value of the fit parameter x0
Definition: TOPGainEfficiencyCalculatorModule.h:182
Belle2::TOPGainEfficiencyCalculatorModule::FitHistograms
void FitHistograms(EHistogramType LoadHisto)
Fit charge (or integrated charged) distribution to calculate gain and efficiency for each channel.
Definition: TOPGainEfficiencyCalculatorModule.cc:286
Belle2::TOPGainEfficiencyCalculatorModule::m_p0HeightIntegral
float m_p0HeightIntegral
Parameter from p0 + x*p1 function that fits height-integral distribution.
Definition: TOPGainEfficiencyCalculatorModule.h:175
Belle2::TOPGainEfficiencyCalculatorModule::endRun
virtual void endRun() override
Draw plots to show fitting results for each channel and save them into a given PDF file (outputPDFFil...
Definition: TOPGainEfficiencyCalculatorModule.cc:178
Belle2::TOPGainEfficiencyCalculatorModule::TOPGainEfficiencyCalculatorModule
TOPGainEfficiencyCalculatorModule()
Constructor.
Definition: TOPGainEfficiencyCalculatorModule.cc:57
Belle2::TOPGainEfficiencyCalculatorModule::m_histoFitRangeIntegral
float m_histoFitRangeIntegral
integral of histogram for a range [threshold, fitMax]
Definition: TOPGainEfficiencyCalculatorModule.h:208
Belle2::TOPGainEfficiencyCalculatorModule::m_x0
float m_x0
fit result of x0
Definition: TOPGainEfficiencyCalculatorModule.h:199
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::TOPGainEfficiencyCalculatorModule::m_pmtChId
short m_pmtChId
PMT channel ID.
Definition: TOPGainEfficiencyCalculatorModule.h:186
Belle2::TOPGainEfficiencyCalculatorModule::LoadHistograms
void LoadHistograms(const std::string &histotype)
Load 2D histograms from a given input file (output of TOPLaserHitSelector) and create timing and char...
Definition: TOPGainEfficiencyCalculatorModule.cc:228
Belle2::TOPGainEfficiencyCalculatorModule::m_p2Error
float m_p2Error
fit error of p2
Definition: TOPGainEfficiencyCalculatorModule.h:202
Belle2::TOPGainEfficiencyCalculatorModule::m_pedestalSigma
float m_pedestalSigma
sigma of pedestal
Definition: TOPGainEfficiencyCalculatorModule.h:183
Belle2::TOPGainEfficiencyCalculatorModule::m_meanPulseHeightError
float m_meanPulseHeightError
histogram mean error of pulse height distribution
Definition: TOPGainEfficiencyCalculatorModule.h:193
Belle2::TOPGainEfficiencyCalculatorModule::TOPGainFunc
static double TOPGainFunc(double *var, double *par)
Fit function of pulse charge (or charnge) distribution for channel(pixel)-by-channel gain extraction,...
Definition: TOPGainEfficiencyCalculatorModule.cc:635
Belle2::TOPGainEfficiencyCalculatorModule::m_initialP0
float m_initialP0
initial value of the fit parameter p0
Definition: TOPGainEfficiencyCalculatorModule.h:179
Belle2::TOPGainEfficiencyCalculatorModule::FindPeakForSmallerXThan
static double FindPeakForSmallerXThan(TH1 *histo, double xmax=0)
Find peak and return its position for a limited range of x (x smaller than the given value (xmax))
Definition: TOPGainEfficiencyCalculatorModule.cc:656
Belle2::TOPGainEfficiencyCalculatorModule::m_initialP2
float m_initialP2
initial value of the fit parameter p2
Definition: TOPGainEfficiencyCalculatorModule.h:181
Belle2::TOPGainEfficiencyCalculatorModule::m_x0Error
float m_x0Error
fit error of x0
Definition: TOPGainEfficiencyCalculatorModule.h:203
Belle2::TOPGainEfficiencyCalculatorModule::m_hvDiff
short m_hvDiff
HV difference from nominal HV value.
Definition: TOPGainEfficiencyCalculatorModule.h:169
Belle2::TOPGainEfficiencyCalculatorModule::m_timeHistogram
TH1D * m_timeHistogram[c_NChannelPerPMT]
hit timing distribution, extracted from m_timeChargeHistogram as a projection along its x-axis.
Definition: TOPGainEfficiencyCalculatorModule.h:151
Belle2::TOPGainEfficiencyCalculatorModule::DummyFillBranch
void DummyFillBranch(EHistogramType LoadHisto)
Fill Dummy for Branch.
Definition: TOPGainEfficiencyCalculatorModule.cc:442
Belle2::TOPGainEfficiencyCalculatorModule::m_p0Error
float m_p0Error
fit error of p0
Definition: TOPGainEfficiencyCalculatorModule.h:200
Belle2::TOPGainEfficiencyCalculatorModule::initialize
virtual void initialize() override
Load time vs charge 2D histogram from a given input file (paramter "inputFile") and prepare hit timin...
Definition: TOPGainEfficiencyCalculatorModule.cc:95
Belle2::TOPGainEfficiencyCalculatorModule::m_pixelId
short m_pixelId
pixel ID, calculated from PMT ID and PMT channel ID
Definition: TOPGainEfficiencyCalculatorModule.h:185
Belle2::Module::addParam
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:562
Belle2::TOPGainEfficiencyCalculatorModule::m_funcFitRangeIntegral
float m_funcFitRangeIntegral
integral of fit function for a range [threshold, fitMax]
Definition: TOPGainEfficiencyCalculatorModule.h:207
Belle2::TOPGainEfficiencyCalculatorModule::defineHisto
virtual void defineHisto() override
Define TTree branches to store fit results for each channel This TTree is saved in an output file giv...
Definition: TOPGainEfficiencyCalculatorModule.cc:100
Belle2::TOPGainEfficiencyCalculatorModule::m_gain
float m_gain
calculated gain from fitting of pulse charge distribution
Definition: TOPGainEfficiencyCalculatorModule.h:194
Belle2::TOPGainEfficiencyCalculatorModule::m_p0
float m_p0
fit result of p0
Definition: TOPGainEfficiencyCalculatorModule.h:196
Belle2::TOPGainEfficiencyCalculatorModule::m_fitMax
float m_fitMax
upper limit of fit region for pulse charge distribution, determined based on m_fracFit value
Definition: TOPGainEfficiencyCalculatorModule.h:177
Belle2::TOPGainEfficiencyCalculatorModule::m_chisquare
float m_chisquare
chi2 of fitting
Definition: TOPGainEfficiencyCalculatorModule.h:204
Belle2::TOPGainEfficiencyCalculatorModule::m_fitoption
std::string m_fitoption
charge histograms fitting option.
Definition: TOPGainEfficiencyCalculatorModule.h:163
Belle2::TOPGainEfficiencyCalculatorModule::m_funcForFitRange
TF1 * m_funcForFitRange[c_NChannelPerPMT]
array of TF1 pointer to store fit function for pulse charge distribution, defined only for fit region
Definition: TOPGainEfficiencyCalculatorModule.h:156
Belle2::TOPGainEfficiencyCalculatorModule::m_targetSlotId
short m_targetSlotId
slot ID
Definition: TOPGainEfficiencyCalculatorModule.h:166
Belle2::TOPGainEfficiencyCalculatorModule::~TOPGainEfficiencyCalculatorModule
virtual ~TOPGainEfficiencyCalculatorModule()
Destructor.
Definition: TOPGainEfficiencyCalculatorModule.cc:93
Belle2::TOPGainEfficiencyCalculatorModule::m_funcForFullRange
TF1 * m_funcForFullRange[c_NChannelPerPMT]
array of TF1 pointer to store fit function for pulse charge distribution, defined only for full range...
Definition: TOPGainEfficiencyCalculatorModule.h:157
Belle2::TOPGainEfficiencyCalculatorModule::m_ndf
int m_ndf
NDF of fitting.
Definition: TOPGainEfficiencyCalculatorModule.h:205
Belle2::TOPGainEfficiencyCalculatorModule::m_timeChargeHistogram
TH2F * m_timeChargeHistogram[c_NChannelPerPMT]
2D histogram of hit timing and pulse charge (or charge), taken from an output file of TOPLaserHitSele...
Definition: TOPGainEfficiencyCalculatorModule.h:150