Belle II Software  release-08-01-10
PXDClusterPositionCalibrationAlgorithm.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/PXDClusterPositionCalibrationAlgorithm.h>
10 
11 #include <string>
12 #include <tuple>
13 #include "TH2F.h"
14 #include "TMatrixDSym.h"
15 #include "TVectorD.h"
16 #include "TMatrixDSymEigen.h"
17 
18 #include <boost/format.hpp>
19 #include <cmath>
20 
21 using namespace std;
22 using boost::format;
23 using namespace Belle2;
24 
25 PXDClusterPositionCalibrationAlgorithm::PXDClusterPositionCalibrationAlgorithm():
26  CalibrationAlgorithm("PXDClusterPositionCollector"),
27  minClusterForShapeLikelyhood(500), minClusterForPositionOffset(2000), maxEtaBins(10),
28  // Branches from TTree
29  m_clusterEta(0), m_positionOffsetU(0), m_positionOffsetV(0), m_sizeV(0),
30  m_pitchV(0), m_clusterKind(0)
31 {
33  " -------------------------- PXDClusterPositionCalibrationAlgorithm ----------------------\n"
34  " \n"
35  " Algorithm for estimating cluster position offsets and shape likelyhoods. \n"
36  " ----------------------------------------------------------------------------------------\n"
37  );
38 }
39 
41 {
42 
43  // Read back the V pitch of all cluster kinds in source data
44  // This avoids relying on VXD::GeoCache.
45  auto pitchtree = getObjectPtr<TTree>("pitchtree");
46  pitchtree->SetBranchAddress("PitchV", &m_pitchV);
47  pitchtree->SetBranchAddress("ClusterKind", &m_clusterKind);
48 
49  for (int i = 0; i < pitchtree->GetEntries(); ++i) {
50  pitchtree->GetEntry(i);
52  }
53 
54  // Buffer temporary payloads for shape calibration for all
55  // clusterkinds and angle bins
56  typedef tuple<int, int, int> bufferkey_t;
57  typedef pair<PXDClusterShapeIndexPar, PXDClusterShapeClassifierPar> buffervalue_t;
58  map<bufferkey_t, buffervalue_t> localCalibrationMap;
59 
60  for (auto clusterKind : clusterKinds) {
61 
62  B2INFO("Start calibration of clusterkind=" << clusterKind << " ...");
63 
64  string gridname = str(format("GridKind_%1%") % clusterKind);
65  auto grid = getObjectPtr<TH2F>(gridname);
66 
67  for (auto uBin = 1; uBin <= grid->GetXaxis()->GetNbins(); uBin++) {
68  for (auto vBin = 1; vBin <= grid->GetYaxis()->GetNbins(); vBin++) {
69 
70  // Bin is centered around angles
71  auto thetaU = grid->GetXaxis()->GetBinCenter(uBin);
72  auto thetaV = grid->GetYaxis()->GetBinCenter(vBin);
73 
74  if (thetaV < 0) {
75  B2INFO("Skip training estimator on thetaU=" << thetaU << ", thetaV=" << thetaV);
76  continue;
77  } else {
78  B2INFO("Start training estimator on thetaU=" << thetaU << ", thetaV=" << thetaV);
79  }
80 
81  string treename = str(format("tree_%1%_%2%_%3%") % clusterKind % uBin % vBin);
82 
83  auto localShapeIndexer = PXDClusterShapeIndexPar();
84  auto localShapeClassifier = PXDClusterShapeClassifierPar();
85  createShapeClassifier(treename, &localShapeClassifier, &localShapeIndexer);
86 
87  bufferkey_t key = std::make_tuple(clusterKind, uBin, vBin);
88  localCalibrationMap[key] = buffervalue_t(localShapeIndexer, localShapeClassifier);
89  }
90  }
91  }
92 
93  // Create a ShapeIndexer payload
95 
96  int globalShapeIndex = 0;
97  for (auto it = m_shapeSet.begin(); it != m_shapeSet.end(); ++it) {
98  shapeIndexer->addShape(*it, globalShapeIndex);
99  globalShapeIndex++;
100  }
101 
102  B2INFO("Number of cluster shapes is " << globalShapeIndex);
103 
104  // Save the cluster shape index table to database.
105  saveCalibration(shapeIndexer, "PXDClusterShapeIndexPar");
106 
107  // Create position estimator
109 
110  for (auto clusterKind : clusterKinds) {
111 
112  string gridname = str(format("GridKind_%1%") % clusterKind);
113  auto grid = getObjectPtr<TH2F>(gridname);
114 
115  positionEstimator->addGrid(clusterKind, *grid);
116 
117  for (auto uBin = 1; uBin <= grid->GetXaxis()->GetNbins(); uBin++) {
118  for (auto vBin = 1; vBin <= grid->GetYaxis()->GetNbins(); vBin++) {
119 
120  // Bin is centered around angles
121  auto thetaU = grid->GetXaxis()->GetBinCenter(uBin);
122  auto thetaV = grid->GetYaxis()->GetBinCenter(vBin);
123 
124  if (thetaV < 0) {
125  // We skipped this part in training before
126  continue;
127  }
128 
129  // Find the local calibration results
130  auto iter = localCalibrationMap.find(std::make_tuple(clusterKind, uBin, vBin));
131  auto localShapeIndexer = iter->second.first;
132  auto localShapeClassifer = iter->second.second;
133 
134  // Require that all shape classifiers use a common shape indexer
135  auto shapeClassifier = localToGlobal(&localShapeClassifer, &localShapeIndexer, shapeIndexer);
136 
137  // Mirror the shape classifier along v
138  auto mirror_vBin = grid->GetYaxis()->FindBin(-thetaV);
139  auto mirroredClassifier = mirrorShapeClassifier(&shapeClassifier, shapeIndexer, clusterKind);
140 
141  // and fill into position estimator payload
142  B2INFO("Add shape classifier for angles thetaU=" << thetaU << ", thetaV=" << thetaV << ", clusterkind=" << clusterKind);
143  positionEstimator->setShapeClassifier(shapeClassifier, uBin, vBin, clusterKind);
144  B2INFO("Add mirrored shape classifier for angles thetaU=" << thetaU << ", thetaV=" << -thetaV << ", clusterkind=" << clusterKind);
145  positionEstimator->setShapeClassifier(mirroredClassifier, uBin, mirror_vBin, clusterKind);
146  }
147  }
148  }
149 
150  // Save the cluster positions to database.
151  saveCalibration(positionEstimator, "PXDClusterPositionEstimatorPar");
152 
153  B2INFO("PXDClusterPosition Calibration Successful");
154  return c_OK;
155 }
156 
157 
159  shapeClassifier, PXDClusterShapeIndexPar* shapeIndexer, int clusterKind)
160 {
161  // Create a mirrored shape classifier
162  auto mirroredShapeClassifier = PXDClusterShapeClassifierPar();
163 
164  // Mirror the shape likelyhood map
165  auto shapeLikelyhoodMap = shapeClassifier->getShapeLikelyhoodMap();
166  for (auto indexAndValue : shapeLikelyhoodMap) {
167  // Compute the mirrored shape index
168  auto shapeIndex = indexAndValue.first;
169  auto shapeName = shapeIndexer->getShapeName(shapeIndex);
170  auto mirroredName = m_mirrorMap[shapeName];
171  auto mirroredIndex = shapeIndexer->getShapeIndex(mirroredName);
172  // Store the result
173  mirroredShapeClassifier.addShapeLikelyhood(mirroredIndex, indexAndValue.second);
174  }
175 
176  // Mirror the offset related maps
177  auto offsetMap = shapeClassifier->getOffsetMap();
178  auto percentileMap = shapeClassifier->getPercentilesMap();
179  auto likelyhoodMap = shapeClassifier->getLikelyhoodMap();
180  for (auto indexAndValue : offsetMap) {
181  // Compute the mirrored shape index
182  auto shapeIndex = indexAndValue.first;
183  auto shapeName = shapeIndexer->getShapeName(shapeIndex);
184  auto mirroredName = m_mirrorMap[shapeName];
185  auto mirroredIndex = shapeIndexer->getShapeIndex(mirroredName);
186 
187  mirroredShapeClassifier.addShape(mirroredIndex);
188 
189  int etaBin = 0;
190  for (auto offset : indexAndValue.second) {
191  // Copy over percentile
192  auto percentile = percentileMap[shapeIndex][etaBin];
193  mirroredShapeClassifier.addEtaPercentile(mirroredIndex, percentile);
194  // Copy over likelyhood
195  auto likelyhood = likelyhoodMap[shapeIndex][etaBin];
196  mirroredShapeClassifier.addEtaLikelyhood(mirroredIndex, likelyhood);
197  // Mirror the offset: v offset shifts and covariance swaps sign
198  double shift = (m_sizeMap[shapeName] - 1) * m_pitchMap[clusterKind];
199  auto mirroredOffset = PXDClusterOffsetPar(offset.getU(), shift - offset.getV(), offset.getUSigma2(), offset.getVSigma2(),
200  -offset.getUVCovariance());
201  mirroredShapeClassifier.addEtaOffset(mirroredIndex, mirroredOffset);
202  etaBin++;
203  }
204  }
205 
206  return mirroredShapeClassifier;
207 }
208 
210  PXDClusterShapeIndexPar* shapeIndexer, PXDClusterShapeIndexPar* globalShapeIndexer)
211 {
212  // Create a shape classifier using global shape indices
213  auto globalShapeClassifier = PXDClusterShapeClassifierPar();
214 
215  // Re-index the the shape likelyhood map
216  auto shapeLikelyhoodMap = shapeClassifier->getShapeLikelyhoodMap();
217  for (auto indexAndValue : shapeLikelyhoodMap) {
218  // Compute the global shape index
219  auto shapeIndex = indexAndValue.first;
220  auto shapeName = shapeIndexer->getShapeName(shapeIndex);
221  auto globalIndex = globalShapeIndexer->getShapeIndex(shapeName);
222  // Store the result
223  globalShapeClassifier.addShapeLikelyhood(globalIndex, indexAndValue.second);
224  }
225 
226  // Re-index the offset related maps
227  auto offsetMap = shapeClassifier->getOffsetMap();
228  auto percentileMap = shapeClassifier->getPercentilesMap();
229  auto likelyhoodMap = shapeClassifier->getLikelyhoodMap();
230  for (auto indexAndValue : offsetMap) {
231  // Compute the global shape index
232  auto shapeIndex = indexAndValue.first;
233  auto shapeName = shapeIndexer->getShapeName(shapeIndex);
234  auto globalIndex = globalShapeIndexer->getShapeIndex(shapeName);
235 
236  globalShapeClassifier.addShape(globalIndex);
237 
238  int etaBin = 0;
239  for (auto offset : indexAndValue.second) {
240  // Copy over percentile
241  auto percentile = percentileMap[shapeIndex][etaBin];
242  globalShapeClassifier.addEtaPercentile(globalIndex, percentile);
243  // Copy over likelyhood
244  auto likelyhood = likelyhoodMap[shapeIndex][etaBin];
245  globalShapeClassifier.addEtaLikelyhood(globalIndex, likelyhood);
246  // Copy over offset
247  globalShapeClassifier.addEtaOffset(globalIndex, offset);
248  etaBin++;
249  }
250  }
251  return globalShapeClassifier;
252 }
253 
255  PXDClusterShapeClassifierPar* shapeClassifier, PXDClusterShapeIndexPar* shapeIndexer)
256 {
257 
258  auto tree = getObjectPtr<TTree>(treename);
259 
260  const auto nEntries = tree->GetEntries();
261  B2INFO("Number of clusters is " << nEntries);
262 
263  string* shapeNamePtr = &m_shapeName;
264  string* mirroredShapeNamePtr = &m_mirroredShapeName;
265 
266  tree->SetBranchAddress("ShapeName", &shapeNamePtr);
267  tree->SetBranchAddress("MirroredShapeName", &mirroredShapeNamePtr);
268  tree->SetBranchAddress("ClusterEta", &m_clusterEta);
269  tree->SetBranchAddress("OffsetU", &m_positionOffsetU);
270  tree->SetBranchAddress("OffsetV", &m_positionOffsetV);
271  tree->SetBranchAddress("SizeV", &m_sizeV);
272 
273  // Vector to enumerate all shapes by unique name and count their
274  // occurence in training data.
275  vector< pair<string, float> > shapeList;
276 
277  for (int i = 0; i < nEntries; ++i) {
278  tree->GetEntry(i);
279 
280  string shapeName = m_shapeName;
281 
282  auto it = std::find_if(shapeList.begin(), shapeList.end(),
283  [&](const pair<string, float>& element) { return element.first == shapeName;});
284 
285  //Shape name exists in vector
286  if (it != shapeList.end()) {
287  //increment key in map
288  it->second++;
289  }
290  //Shape name does not exist
291  else {
292  //Not found, insert in vector
293  shapeList.push_back(pair<string, int>(shapeName, 1));
294  // Remember the relation between name of a shape and name of mirrored shape
295  m_mirrorMap[shapeName] = m_mirroredShapeName;
296  // Remember the relation between name of a shape and its size
297  m_sizeMap[shapeName] = m_sizeV;
298  }
299  }
300 
301  // Loop over shapeList to select shapes for
302  // next calibration step
303 
304  // Vector with eta histograms for selected shapes
305  vector< pair<string, TH1D> > etaHistos;
306 
307  // Index for enumerating selected shapes
308  int tmpIndex = 0;
309 
310  // Coverage of position offsets on training data
311  double coverage = 0.0;
312 
313  for (auto iter : shapeList) {
314  auto name = iter.first;
315  auto counter = iter.second;
316 
317  double likelyhood = counter / nEntries;
318 
319  if (counter >= minClusterForShapeLikelyhood) {
320  //B2INFO("Adding shape " << name << " with index " << tmpIndex << " and shape likelyhood " << 100*likelyhood << "% and count " << counter);
321  shapeIndexer->addShape(name, tmpIndex);
322  shapeClassifier->addShapeLikelyhood(tmpIndex, likelyhood);
323 
324  // Add name of shape to global (all clusterkinds + all angle bins) shape set
325  m_shapeSet.insert(name);
326  // Add name of mirrored shape as well
327  m_shapeSet.insert(m_mirrorMap[name]);
328  // Increment the index
329  tmpIndex++;
330  }
331 
332  if (counter >= minClusterForPositionOffset) {
333  coverage += likelyhood;
334  string etaname = str(format("eta_%1%") % name);
335 
336  // Single pixel case: Eta value is cluster charge
337  if (name == "SD0.0") {
338  TH1D etaHisto(etaname.c_str(), etaname.c_str(), 255, 0, 255);
339  etaHisto.SetDirectory(0);
340  etaHistos.push_back(pair<string, TH1D>(name, etaHisto));
341  }
342  // Multipixel case: Eta value is ratio head/(tail+head) of charges (to be less gain sensitive)
343  else {
344  TH1D etaHisto(etaname.c_str(), etaname.c_str(), 301, 0, 1);
345  etaHisto.SetDirectory(0);
346  etaHistos.push_back(pair<string, TH1D>(name, etaHisto));
347  }
348  }
349  }
350 
351  B2INFO("Offset coverage is " << 100 * coverage << "%");
352 
353  // Loop over the tree is to fill the eta histograms for
354  // selected shapes.
355 
356  for (int i = 0; i < nEntries; ++i) {
357  tree->GetEntry(i);
358 
359  string shapeName = m_shapeName;
360  auto it = std::find_if(etaHistos.begin(), etaHistos.end(),
361  [&](const pair<string, TH1D>& element) { return element.first == shapeName;});
362  //Item exists in map
363  if (it != etaHistos.end()) {
364  // increment key in map
365  it->second.Fill(m_clusterEta);
366  }
367  }
368 
369  // Vector for offset histograms stored by offset shape name and eta bin
370  vector< pair< string, vector<TH2D> > > offsetHistosVec;
371 
372  for (auto iter : etaHistos) {
373  auto name = iter.first;
374  auto& histo = iter.second;
375  int nClusters = histo.GetEntries();
376 
377  // Add shape for offset correction
378  int shapeIndex = shapeIndexer->getShapeIndex(name);
379  shapeClassifier->addShape(shapeIndex);
380 
381  // Try to split clusters into n bins with minClusterForPositionOffset clusters
382  int nEtaBins = std::max(int(nClusters / minClusterForPositionOffset), 1);
383  nEtaBins = std::min(nEtaBins, maxEtaBins);
384 
385  //B2INFO("SHAPE NAME:" << name << " WITH BINS " << nEtaBins);
386 
387  vector< TH2D > offsetHistos;
388 
389  for (int i = 0; i < nEtaBins; i++) {
390  // Position where to compute the quantiles in [0,1]
391  double xq = double(i) / nEtaBins;
392  // Double to contain the quantile
393  double yq = 0;
394  histo.GetQuantiles(1, &yq, &xq);
395  //B2INFO(" Quantile at xq =" << xq << " is yq=" << yq);
396  shapeClassifier->addEtaPercentile(shapeIndex, yq);
397 
398  string offsetname = str(format("offset_%1%_%2%") % name % i);
399  TH2D offsetHisto(offsetname.c_str(), offsetname.c_str(), 1, 0, 1, 1, 0, 1);
400  offsetHisto.SetDirectory(0);
401  offsetHisto.StatOverflows();
402  offsetHistos.push_back(offsetHisto);
403 
404  }
405  offsetHistosVec.push_back(pair< string, vector<TH2D> >(name, offsetHistos));
406  }
407 
408  // Loop over the tree is to fill offset histograms
409 
410  for (int i = 0; i < nEntries; ++i) {
411  tree->GetEntry(i);
412 
413  string shapeName = m_shapeName;
414  auto it = std::find_if(offsetHistosVec.begin(), offsetHistosVec.end(),
415  [&](const pair<string, vector<TH2D>>& element) { return element.first == shapeName;});
416  //Item exists in map
417  if (it != offsetHistosVec.end()) {
418  int shapeIndex = shapeIndexer->getShapeIndex(shapeName);
419  int etaBin = shapeClassifier->getEtaIndex(shapeIndex, m_clusterEta);
420  it->second.at(etaBin).Fill(m_positionOffsetU, m_positionOffsetV);
421  }
422  }
423 
424  // Compute the moments of the offset histograms and finalize the shape classifier object
425 
426  // Loop over shape names
427  for (auto iter : offsetHistosVec) {
428  auto name = iter.first;
429  auto& histovec = iter.second;
430 
431  int shapeIndex = shapeIndexer->getShapeIndex(name);
432 
433  // Loop over eta bins
434  for (auto& histo : histovec) {
435  // Compute offset moments
436  double etaLikelyhood = double(histo.GetEntries()) / nEntries;
437  double offsetU = histo.GetMean(1);
438  double offsetV = histo.GetMean(2);
439  double covUV = histo.GetCovariance();
440  double covU = pow(histo.GetRMS(1), 2);
441  double covV = pow(histo.GetRMS(2), 2);
442 
443  B2INFO("Name " << name << ", posU=" << offsetU << ", posV=" << offsetV << ", covU=" << covU << ", covV=" << covV << ", covUV=" <<
444  covUV);
445 
446  TMatrixDSym HitCov(2);
447  HitCov(0, 0) = covU;
448  HitCov(1, 0) = covUV;
449  HitCov(0, 1) = covUV;
450  HitCov(1, 1) = covV;
451 
452  TMatrixDSymEigen HitCovE(HitCov);
453  TVectorD eigenval = HitCovE.GetEigenValues();
454  if (eigenval(0) <= 0 || eigenval(1) <= 0) {
455  B2ERROR("Estimated covariance matrix not positive definite.");
456  }
457 
458  auto offset = PXDClusterOffsetPar(offsetU, offsetV, covU, covV, covUV);
459  shapeClassifier->addEtaLikelyhood(shapeIndex, etaLikelyhood);
460  shapeClassifier->addEtaOffset(shapeIndex, offset);
461  }
462  }
463 
464  B2INFO("Added shape classifier with coverage " << 100 * coverage << "% on training data sample.");
465 
466  return;
467 }
468 
469 
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.
The class for PXD cluster position offset payload.
void createShapeClassifier(std::string treename, PXDClusterShapeClassifierPar *shapeClassifier, PXDClusterShapeIndexPar *shapeIndexer)
Returns a new classifier and index trained on cluster tree.
std::string m_mirroredShapeName
Name of mirrored cluster shape.
std::vector< int > clusterKinds
Vector of clusterkinds to calibrate.
PXDClusterShapeClassifierPar mirrorShapeClassifier(PXDClusterShapeClassifierPar *shapeClassifier, PXDClusterShapeIndexPar *shapeIndexer, int clusterKind)
Returns a mirrored version of shape classifier.
PXDClusterShapeClassifierPar localToGlobal(PXDClusterShapeClassifierPar *localShapeClassifier, PXDClusterShapeIndexPar *localShapeIndexer, PXDClusterShapeIndexPar *globalShapeIndexer)
Returns a shape classifier using global shape indices instead of local ones.
std::map< std::string, int > m_sizeMap
Helper needed to map the name of a shape to the V size of the cluster.
int maxEtaBins
Maximum number of eta bins for estimating cluster position offsets.
int minClusterForPositionOffset
Minimum number of collected clusters for estimating cluster position offsets.
int minClusterForShapeLikelyhood
Minimum number of collected clusters for estimating shape likelyhood.
std::map< int, float > m_pitchMap
Helper needed to map the clusterkind to the V pitch of the sensor.
std::set< std::string > m_shapeSet
Set of unique shape names.
std::map< std::string, std::string > m_mirrorMap
Helper needed to map the name of a shape to the name of the mirrored shape.
The class for PXD cluster position lookup table payload.
void setShapeClassifier(const PXDClusterShapeClassifierPar &classifier, int uBin, int vBin, int clusterkind)
Set shape classifier.
void addGrid(int clusterkind, const TH2F &grid)
Add grid for clusterkind.
The class for PXD cluster shape classifier payload.
const std::map< int, std::vector< float > > & getPercentilesMap() const
Return percentiles map for position correction.
void addShapeLikelyhood(int shape_index, float likelyhood)
Add shape likelyhood.
const std::map< int, float > & getShapeLikelyhoodMap() const
Return shape likelyhood map
void addEtaLikelyhood(int shape_index, float likelyhood)
Add eta likelyhood to shape for position correction.
void addEtaPercentile(int shape_index, float percentile)
Add eta percentile to shape for position correction.
unsigned int getEtaIndex(int shape_index, float eta) const
Get eta index for position correction.
const std::map< int, std::vector< float > > & getLikelyhoodMap() const
Return likelyhood map for position correction.
void addShape(int shape_index)
Add shape for position correction.
const std::map< int, std::vector< PXDClusterOffsetPar > > & getOffsetMap() const
Return offset map for position correction.
void addEtaOffset(int shape_index, PXDClusterOffsetPar &offset)
Add offset to shape for position correction.
The class for PXD cluster shape index payload.
int getShapeIndex(const std::string &name) const
Returns shape index from name.
void addShape(const std::string &name, int index)
Add shape with name and index
const std::string & getShapeName(int index) const
Returns shape name from index.
Abstract base class for different kinds of events.