Belle II Software  release-08-01-10
PXDDataMCGainCalibrationAlgorithm.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 #include <pxd/calibration/PXDDataMCGainCalibrationAlgorithm.h>
10 #include <pxd/dbobjects/PXDClusterChargeMapPar.h>
11 #include <pxd/dbobjects/PXDGainMapPar.h>
12 
13 #include <string>
14 #include <algorithm>
15 #include <map>
16 #include <set>
17 
18 #include <sstream>
19 #include <iostream>
20 
21 #include <boost/format.hpp>
22 #include <cmath>
23 
24 
25 //ROOT
26 #include <TRandom.h>
27 #include <TF1.h>
28 #include <TH2I.h>
29 #include <TCanvas.h>
30 
31 using namespace std;
32 using boost::format;
33 using namespace Belle2;
34 
35 
36 // Anonymous namespace for data objects used by PXDDataMCGainCalibrationAlgorithm class
37 namespace {
38 
40  int m_signal;
42  int m_run;
44  int m_exp;
45 
47  void getNumberOfBins(const std::shared_ptr<TH1I>& histo_ptr, unsigned short& nBinsU, unsigned short& nBinsV)
48  {
49  set<unsigned short> uBinSet;
50  set<unsigned short> vBinSet;
51 
52  // Loop over all bins of input histo
53  for (auto histoBin = 1; histoBin <= histo_ptr->GetNbinsX(); histoBin++) {
54  // The bin label contains the vxdid, uBin and vBin
55  string label = histo_ptr->GetXaxis()->GetBinLabel(histoBin);
56 
57  // Parse label string format to read sensorID, uBin and vBin
58  istringstream stream(label);
59  string token;
60  getline(stream, token, '_');
61  getline(stream, token, '_');
62  unsigned short uBin = std::stoi(token);
63 
64  getline(stream, token, '_');
65  unsigned short vBin = std::stoi(token);
66 
67  uBinSet.insert(uBin);
68  vBinSet.insert(vBin);
69  }
70 
71  if (uBinSet.empty() || vBinSet.empty()) {
72  B2FATAL("Not able to determine the grid size. Something is wrong with collected data.");
73  } else {
74  nBinsU = *uBinSet.rbegin() + 1;
75  nBinsV = *vBinSet.rbegin() + 1;
76  }
77  }
78 
80  unsigned short getNumberOfSensors(const std::shared_ptr<TH1I>& histo_ptr)
81  {
82  set<unsigned short> sensorSet;
83 
84  // Loop over all bins of input histo
85  for (auto histoBin = 1; histoBin <= histo_ptr->GetNbinsX(); histoBin++) {
86  // The bin label contains the vxdid, uBin and vBin
87  string label = histo_ptr->GetXaxis()->GetBinLabel(histoBin);
88 
89  // Parse label string format to read sensorID, uBin and vBin
90  istringstream stream(label);
91  string token;
92  getline(stream, token, '_');
93  VxdID sensorID(token);
94  sensorSet.insert(sensorID.getID());
95  }
96 
97  return sensorSet.size();
98  }
99 
100 }
101 
102 
103 PXDDataMCGainCalibrationAlgorithm::PXDDataMCGainCalibrationAlgorithm():
104  CalibrationAlgorithm("PXDClusterChargeCollector"),
105  minClusters(5000), noiseSigma(0.6), safetyFactor(2.0), forceContinue(false), strategy(0),
106  doCalibration(false), useChargeHistogram(false), chargePayloadName("PXDClusterChargeMapPar")
107 {
109  " ------------------------------ PXDDataMCGainCalibrationAlgorithm -----------------------------------\n"
110  " \n"
111  " Algorithm for estimating pxd median/MPV cluster charges for different position on sensor in ADU, and optionally calculate gain \n"
112  " ------------------------------------------------------------------------------------------------\n"
113  );
114 
115 }
116 
117 
119 {
120 
121  // Get counter histogram
122  auto cluster_counter = getObjectPtr<TH1I>("PXDClusterCounter");
123  if (!cluster_counter) {
124  B2INFO("Not enough Data: cluster counter does not exist ");
125  return c_NotEnoughData;
126  }
127 
128  // Extract number of sensors from counter histograms
129  auto nSensors = getNumberOfSensors(cluster_counter);
130 
131  // Extract the number of grid bins from counter histograms
132  unsigned short nBinsU = 0;
133  unsigned short nBinsV = 0;
134  getNumberOfBins(cluster_counter, nBinsU, nBinsV);
135 
136  // Check that we have collected enough Data
137  if (cluster_counter->GetEntries() < int(safetyFactor * minClusters * nSensors * nBinsU * nBinsV)) {
138  if (not forceContinue) {
139  B2INFO("Not enough Data: Only " << cluster_counter->GetEntries() << " hits were collected but " << int(safetyFactor * minClusters*
140  nSensors * nBinsU * nBinsV) << " needed!");
141  return c_NotEnoughData;
142  } else {
143  B2INFO("Continue despite low statistics: Only " << cluster_counter->GetEntries() << " hits were collected but" << int(
145  nSensors * nBinsU * nBinsV) << " would be desirable!");
146  }
147  }
148 
149  B2INFO("Start calibration using a " << nBinsU << "x" << nBinsV << " grid per sensor.");
150 
151  // This is the PXD charge calibration payload for conditions DB
152  PXDClusterChargeMapPar* chargeMapPar = new PXDClusterChargeMapPar(nBinsU, nBinsV);
153  // This is the PXD gain correction payload for conditions DB
154  PXDGainMapPar* gainMapPar = new PXDGainMapPar(nBinsU, nBinsV);
155  set<VxdID> pxdSensors;
156 
157  // Read back existing DB payloads
158  PXDClusterChargeMapPar* chargeMapPtr = nullptr;
159  PXDGainMapPar* gainMapPtr = nullptr;
160  auto dbtree = getObjectPtr<TTree>("dbtree");
161  dbtree->SetBranchAddress("run", &m_run);
162  dbtree->SetBranchAddress("exp", &m_exp);
163  dbtree->SetBranchAddress("chargeMap", &chargeMapPtr);
164  dbtree->SetBranchAddress("gainMap", &gainMapPtr);
165 
166  // Loop over all bins of input histo
167  for (auto histoBin = 1; histoBin <= cluster_counter->GetNbinsX(); histoBin++) {
168 
169  // The bin label contains the vxdid, uBin and vBin
170  string label = cluster_counter->GetXaxis()->GetBinLabel(histoBin);
171 
172  // Parse label string format to read sensorID, uBin and vBin
173  istringstream stream(label);
174  string token;
175  getline(stream, token, '_');
176  VxdID sensorID(token);
177 
178  getline(stream, token, '_');
179  unsigned short uBin = std::stoi(token);
180 
181  getline(stream, token, '_');
182  unsigned short vBin = std::stoi(token);
183 
184  // Read back the counters for number of collected clusters
185  int numberOfDataHits = cluster_counter->GetBinContent(histoBin);
186 
187  // Only perform fitting, when enough data is available
188  if (numberOfDataHits >= minClusters) {
189 
190  B2INFO("start EstimateCharge");
191  // Compute the cluster charge or gain for the part of PXD
192  auto Charge = EstimateCharge(sensorID, uBin, vBin, histoBin);
193 
194  // Store the charge or gain
195  if (!doCalibration) {
196 
197  B2INFO("EstimateCharge: sensor " << sensorID.getID() << " U " << uBin << " V " << vBin << " Charge " << Charge);
198  chargeMapPar->setContent(sensorID.getID(), uBin, vBin, Charge);
199 
200  } else {
201 
202  auto Gain = 0.0;
203 
204  if (Charge <= 0.0) {
205  B2WARNING("Retrieved negative charge for data for sensor=" << sensorID << " uBin=" << uBin << " vBin=" << vBin <<
206  ". Set gain to default value (=1.0).");
207  Gain = 1.0;
208  } else {
209  dbtree->GetEntry(0);
210  double mcCharge = chargeMapPtr->getContent(sensorID.getID(), uBin, vBin);
211  //GetChargeFromDB(sensorID, uBin, vBin, dbtree);
212  if (mcCharge <= 0.0) {
213  B2WARNING("Retrieved negative charge for MC from DB for sensor=" << sensorID << " uBin=" << uBin << " vBin=" << vBin <<
214  ". Set gain to default value (=1.0).");
215  Gain = 1.0;
216  } else {
217  Gain = Charge / mcCharge;
218  B2INFO("Estimated Gain: sensor " << sensorID.getID() << " U " << uBin << " V " << vBin << " Gain " << Gain << " = " << Charge <<
219  " / " << mcCharge);
220  }
221 
222  }
223 
224  gainMapPar->setContent(sensorID.getID(), uBin, vBin, Gain);
225 
226  }
227 
228  } else {
229 
230  B2WARNING(label << ": Number of data hits too small for fitting (" << numberOfDataHits << " < " << minClusters <<
231  "). Use default value of 0 for charge, 1 for gain.");
232  if (!doCalibration) chargeMapPar->setContent(sensorID.getID(), uBin, vBin, 0.0);
233  else gainMapPar->setContent(sensorID.getID(), uBin, vBin, 1.0);
234 
235  }
236 
237  pxdSensors.insert(sensorID);
238 
239  }
240 
241  chargeMapPtr = nullptr;
242  gainMapPtr = nullptr;
243 
244  // Save the charge map to database. Note that this will set the database object name to the same as the collector but you
245  // are free to change it.
246  if (!doCalibration) {
247  saveCalibration(chargeMapPar, chargePayloadName);
248 
249  B2INFO("PXD Cluster Charge Calibration Successful");
250  return c_OK;
251 
252  } else {
253 
254  // Post processing of gain map. It is possible that the gain
255  // computation failed on some parts. Here, we replace default
256  // values (1.0) by local averages of neighboring sensor parts.
257 
258  for (const auto& sensorID : pxdSensors) {
259  for (unsigned short vBin = 0; vBin < nBinsV; ++vBin) {
260  float meanGain = 0;
261  unsigned short nGood = 0;
262  unsigned short nBad = 0;
263  for (unsigned short uBin = 0; uBin < nBinsU; ++uBin) {
264  auto gain = gainMapPar->getContent(sensorID.getID(), uBin, vBin);
265  // Filter default gains
266  if (gain != 1.0) {
267  nGood += 1;
268  meanGain += gain;
269  } else {
270  nBad += 1;
271  }
272  }
273  B2RESULT("Gain calibration on sensor=" << sensorID << " and vBin=" << vBin << " was successful on " << nGood << "/" << nBinsU <<
274  " uBins.");
275 
276  // Check if we can repair bad calibrations with a local avarage
277  if (nGood > 0 && nBad > 0) {
278  meanGain /= nGood;
279  for (unsigned short uBin = 0; uBin < nBinsU; ++uBin) {
280  auto gain = gainMapPar->getContent(sensorID.getID(), uBin, vBin);
281  if (gain == 1.0) {
282  gainMapPar->setContent(sensorID.getID(), uBin, vBin, meanGain);
283  B2RESULT("Gain calibration on sensor=" << sensorID << ", vBin=" << vBin << " uBin " << uBin << ": Replace default gain wih average "
284  << meanGain);
285  }
286  }
287  }
288  }
289  }
290 
291  // Save the gain map to database. Note that this will set the database object name to the same as the collector but you
292  // are free to change it.
293  saveCalibration(gainMapPar, "PXDGainMapPar");
294 
295  B2INFO("PXD Gain Calibration Successful");
296  return c_OK;
297 
298  }
299 
300 }
301 
302 
303 double PXDDataMCGainCalibrationAlgorithm::EstimateCharge(VxdID sensorID, unsigned short uBin, unsigned short vBin,
304  unsigned short histoBin)
305 {
306 
307  double charge = -1.;
308 
309  if (strategy < 0 || strategy > 2) {
310  B2FATAL("strategy unavailable, use 0 for medians, 1 for landau fit and 2 for mean!");
311  }
312 
313  // Construct a tree name for requested part of PXD
314  auto layerNumber = sensorID.getLayerNumber();
315  auto ladderNumber = sensorID.getLadderNumber();
316  auto sensorNumber = sensorID.getSensorNumber();
317 
318  if (!useChargeHistogram) {
319  const string treename = str(format("tree_%1%_%2%_%3%_%4%_%5%") % layerNumber % ladderNumber % sensorNumber % uBin % vBin);
320  // Vector with cluster signals from collected data
321  vector<double> signals;
322  // Fill data_signal vector from input data
323  auto tree = getObjectPtr<TTree>(treename);
324  tree->SetBranchAddress("signal", &m_signal);
325 
326  // Loop over tree
327  const auto nEntries = tree->GetEntries();
328  int incr(1);
329  if (int(nEntries / minClusters) > 2) incr = int(2 * nEntries / minClusters);
330  for (int i = 0; i < nEntries; i += incr) {
331  tree->GetEntry(i);
332  double noise = gRandom->Gaus(0.0, noiseSigma);
333  signals.push_back(m_signal + noise); //qyliu: why we introduce noise simulation here?
334  }
335 
336  if (strategy == 0) {
337  double median = CalculateMedian(signals);
338  B2INFO("EstimateCharge: sensor " << sensorID.getID() << "(" << layerNumber << "," << ladderNumber << "," << sensorNumber
339  << ") U " << uBin << " V " << vBin
340  << " Charge " << median);
341  charge = median;
342  } else if (strategy == 1) {
343  double median = CalculateMedian(signals);
344  double landaumpv = FitLandau(signals);
345  double diff = (landaumpv - median);
346  double difff = 0.;
347  if (landaumpv > 0.) difff = (landaumpv - median) / landaumpv;
348  B2INFO("EstimateCharge: sensor " << sensorID.getID() << "(" << layerNumber << "," << ladderNumber << "," << sensorNumber
349  << ") U " << uBin << " V " << vBin
350  << " Charge " << landaumpv << " Median " << median << " diff = " << diff << "/" << difff);
351  charge = landaumpv; //FitLandau(signals);
352  } else if (strategy == 2) {
353  double mean = 0;
354  for (auto& each : signals)
355  mean += each;
356  if (signals.size() > 0) mean = mean / signals.size();
357  charge = mean;
358  }
359 
360  } else {
361 
362  auto cluster_counter = getObjectPtr<TH2I>("PXDClusterCharge");
363  TH1D* hist_signals = cluster_counter->ProjectionY("proj", histoBin, histoBin);
364 
365  if (strategy == 0) {
366  double median = CalculateMedian(hist_signals);
367  B2INFO("EstimateCharge: sensor " << sensorID.getID() << "(" << layerNumber << "," << ladderNumber << "," << sensorNumber
368  << ") U " << uBin << " V " << vBin
369  << " Charge " << median);
370  charge = median;
371  } else if (strategy == 1) {
372  double median = CalculateMedian(hist_signals);
373  double landaumpv = FitLandau(hist_signals);
374  double diff = (landaumpv - median);
375  double difff = 0.;
376  if (landaumpv > 0.) difff = (landaumpv - median) / landaumpv;
377  B2INFO("EstimateCharge: sensor " << sensorID.getID() << "(" << layerNumber << "," << ladderNumber << "," << sensorNumber
378  << ") U " << uBin << " V " << vBin
379  << " Charge " << landaumpv << " Median " << median << " diff = " << diff << "/" << difff);
380 
381  //return landaumpv; //FitLandau(signals);
382  charge = landaumpv;
383  } else if (strategy == 2) {
384  charge = hist_signals->GetMean();
385  }
386 
387  delete hist_signals;
388  }
389 
390  return charge;
391 
392 }
393 
395 {
396  auto size = signals.size();
397 
398  if (size == 0) {
399  return 0.0; // Undefined, really.
400  } else {
401  sort(signals.begin(), signals.end());
402  if (size % 2 == 0) {
403  return (signals[size / 2 - 1] + signals[size / 2]) / 2;
404  } else {
405  return signals[size / 2];
406  }
407  }
408 }
409 
411 {
412  auto size = hist_signals->GetEntries();
413 
414  if (size == 0) return 0.0; // Undefined.
415 
416  int sum = 0;
417  for (int ibin = 0; ibin < hist_signals->GetNbinsX(); ++ibin) {
418  sum += hist_signals->GetBinContent(ibin + 1);
419  if (sum > size / 2) {
420  return hist_signals->GetBinLowEdge(ibin + 1);
421  }
422  }
423 
424  B2WARNING("Could not find median! using default value 0.0!");
425  return 0.0;
426 }
427 
428 double PXDDataMCGainCalibrationAlgorithm::FitLandau(vector<double>& signals)
429 {
430  auto size = signals.size();
431  if (size == 0) return 0.0; // Undefined, really.
432 
433  // get max and min values of vector
434  int max = *max_element(signals.begin(), signals.end());
435  int min = *min_element(signals.begin(), signals.end());
436  // make even bins
437  if ((max - min) % 2) max--;
438  int nbin = max - min;
439 
440  // create histogram to hold signals and fill it
441  TH1D* hist_signals = new TH1D("", "", nbin, min, max);
442  for (auto it = signals.begin(); it != signals.end(); ++it) {
443  hist_signals->Fill(*it);
444  }
445 
446  // create fit function
447  TF1* landau = new TF1("landau", "TMath::Landau(x,[0],[1])*[2]", min, max);
448  landau->SetParNames("MPV", "sigma", "scale");
449  landau->SetParameters(35, 8, 1000);
450  landau->SetParLimits(0, 0., 80.);
451 
452  //do fit and get results, fit range restricted to exclude low charge peak
453  float fitmin(min);
454  float fitmax(350);
455  Int_t status = hist_signals->Fit("landau", "Lq", "", fitmin, fitmax);
456  double MPV = landau->GetParameter("MPV");
457 
458  B2INFO("Fit result: " << status << " MPV " << MPV << " sigma " << landau->GetParameter("sigma")
459  << " scale " << landau->GetParameter("scale") << " chi2 " << landau->GetChisquare());
460 
461  // clean up
462  delete hist_signals;
463  delete landau;
464 
465  // check fit status
466  if (status == 0) return MPV;
467  else {
468  B2WARNING("Fit failed! using default value 0.0!");
469  return 0.0;
470  }
471 }
472 
474 {
475  auto size = hist_signals->GetEntries();
476  if (size == 0) return 0.0; // Undefined.
477 
478  int max = hist_signals->GetBinLowEdge(hist_signals->GetNbinsX() + 1);
479  int min = hist_signals->GetBinLowEdge(1);
480 
481  // create fit function
482  TF1* landau = new TF1("landau", "TMath::Landau(x,[0],[1])*[2]", min, max);
483  landau->SetParNames("MPV", "sigma", "scale");
484  landau->SetParameters(35, 8, 1000);
485  landau->SetParLimits(0, 0., 80.);
486 
487  // do fit and get results, fit range restricted to exclude low charge peak
488  float fitmin(min);
489  float fitmax(250);
490  Int_t status = hist_signals->Fit("landau", "Lq", "", fitmin, fitmax);
491  double MPV = landau->GetParameter("MPV");
492 
493  B2INFO("Fit result: " << status << " MPV " << MPV << " sigma " << landau->GetParameter("sigma")
494  << " scale " << landau->GetParameter("scale") << " chi2 " << landau->GetChisquare());
495 
496  // clean up
497  delete landau;
498 
499  // check fit status
500  if (status == 0) return MPV;
501  else {
502  B2WARNING("Fit failed! using default value 0.0!");
503  return 0.0;
504  }
505 }
506 
Base class for calibration algorithms.
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
void setDescription(const std::string &description)
Set algorithm description (in constructor)
EResult
The result of calibration.
@ c_OK
Finished successfuly =0 in Python.
@ c_NotEnoughData
Needs more data =2 in Python.
The payload class for PXD cluster charge calibrations.
float getContent(unsigned short sensorID, unsigned short globalID) const
Get content.
void setContent(unsigned short sensorID, unsigned short globalID, float value)
Set map content.
bool doCalibration
flag to perform full calibration or only esitmate charge: False: only estimate charge,...
double FitLandau(std::vector< double > &signals)
calculate MPV of unsorted signal vector using a Landau fit
double EstimateCharge(VxdID sensorID, unsigned short uBin, unsigned short vBin, unsigned short histoBin)
Estimate median charge form collected clusters on part of PXD.
int minClusters
Minimum number of collected clusters for estimating median charge.
std::string chargePayloadName
Payload name for Cluster Charge.
float noiseSigma
Artificial noise sigma for smearing cluster charge.
bool useChargeHistogram
Flag to use histogram as charge input.
virtual EResult calibrate() override
Run algo on data.
int strategy
strategy to used for gain calibration, 0 for medians, 1 for landau fit
float safetyFactor
Safety factor for determining whether the collected number of clusters is enough.
double CalculateMedian(std::vector< double > &signals)
Calculate a median from unsorted signal vector. The input vector gets sorted.
bool forceContinue
Force continue in low statistics runs instead of returning c_NotEnoughData.
The payload class for PXD gain corrections.
Definition: PXDGainMapPar.h:43
float getContent(unsigned short sensorID, unsigned short globalID) const
Get content.
void setContent(unsigned short sensorID, unsigned short globalID, float value)
Set map content.
Definition: PXDGainMapPar.h:68
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:33
baseType getID() const
Get the unique id.
Definition: VxdID.h:94
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:100
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:98
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:96
unsigned short getNumberOfSensors(const std::shared_ptr< TH1I > &histo_ptr)
Helper function to extract number of sensors from counter histogram labels.
void getNumberOfBins(const std::shared_ptr< TH1I > &histo_ptr, unsigned short &nBinsU, unsigned short &nBinsV)
Helper function to extract number of bins along u side and v side from counter histogram labels.
Abstract base class for different kinds of events.