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