Belle II Software  release-05-02-19
KLMStripEfficiencyAlgorithm.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Vitaliy Popov, Dmytro Minchenko *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 /* Own header. */
12 #include <klm/calibration/KLMStripEfficiencyAlgorithm.h>
13 
14 /* KLM headers. */
15 #include <klm/dataobjects/KLMChannelIndex.h>
16 
17 /* ROOT headers. */
18 #include <TFile.h>
19 #include <TH1F.h>
20 
21 using namespace Belle2;
22 
24 {
25 }
26 
28 {
29  m_AchievedPrecision = results.m_AchievedPrecision;
30  m_MatchedDigits = results.m_MatchedDigits;
31  m_ExtHits = results.m_ExtHits;
33  m_Efficiency = new float[nPlanes];
34  m_ExtHitsPlane = new int[nPlanes];
35  std::memcpy(m_Efficiency, results.m_Efficiency, nPlanes * sizeof(float));
36  std::memcpy(m_ExtHitsPlane, results.m_ExtHitsPlane, nPlanes * sizeof(int));
37 }
38 
40 {
41  if (m_Efficiency != nullptr)
42  delete m_Efficiency;
43  if (m_ExtHitsPlane != nullptr)
44  delete m_ExtHitsPlane;
45 }
46 
48  CalibrationAlgorithm("KLMStripEfficiencyCollector"),
49  m_ElementNumbers(&(KLMElementNumbers::Instance())),
50  m_PlaneArrayIndex(&(KLMPlaneArrayIndex::Instance())),
51  m_StripEfficiency(nullptr)
52 {
53  int nPlanes = m_PlaneArrayIndex->getNElements();
54  m_Results.m_Efficiency = new float[nPlanes];
55  m_Results.m_ExtHitsPlane = new int[nPlanes];
56 }
57 
59 {
60  if (m_StripEfficiency != nullptr)
61  delete m_StripEfficiency;
62 }
63 
65 {
66  /* Get collected results. */
67  int nPlanes = m_PlaneArrayIndex->getNElements();
68  TH1F* efficiencyHistogram =
69  new TH1F("plane_efficiency", "KLM plane efficiency",
70  nPlanes, -0.5, double(nPlanes) - 0.5);
71  std::shared_ptr<TH1F> matchedDigitsInPlane;
72  matchedDigitsInPlane = getObjectPtr<TH1F>("matchedDigitsInPlane");
73  m_Results.m_MatchedDigits = matchedDigitsInPlane->Integral();
74  std::shared_ptr<TH1F> allExtHitsInPlane;
75  allExtHitsInPlane = getObjectPtr<TH1F>("allExtHitsInPlane");
76  m_Results.m_ExtHits = allExtHitsInPlane->Integral();
77  matchedDigitsInPlane.get()->Sumw2();
78  allExtHitsInPlane.get()->Sumw2();
79  efficiencyHistogram->Divide(matchedDigitsInPlane.get(),
80  allExtHitsInPlane.get(), 1, 1, "B");
81  for (int i = 0; i < nPlanes; ++i) {
82  m_Results.m_Efficiency[i] = efficiencyHistogram->GetBinContent(i + 1);
83  m_Results.m_ExtHitsPlane[i] = allExtHitsInPlane->GetBinContent(i + 1);
84  }
85  /* Check whether the amount of data is sufficient. */
86  bool notEnoughData = false;
89  for (KLMChannelIndex& klmPlane : klmPlanes) {
90  uint16_t plane = klmPlane.getKLMPlaneNumber();
91  uint16_t planeIndex = m_PlaneArrayIndex->getIndex(plane);
92  int extHits = allExtHitsInPlane->GetBinContent(planeIndex + 1);
93  float efficiencyError = efficiencyHistogram->GetBinError(planeIndex + 1);
94  if (efficiencyError > m_Results.m_AchievedPrecision)
95  m_Results.m_AchievedPrecision = efficiencyError;
96  /*
97  * No hits is not considered as "not enough data", because this can
98  * happen in case KLM is excluded.
99  */
100  switch (m_CalibrationStage) {
102  if (extHits != 0 && extHits < m_MinimalExtHits)
103  notEnoughData = true;
104  break;
106  if (efficiencyError > m_RequestedPrecision)
107  notEnoughData = true;
108  break;
109  }
110  }
111  /*
112  * Fill the payload. A new object is created, because saveCalibration()
113  * stores a pointer to KLMStripEfficiency, and it is necessary to save
114  * the payloads to commit them at the end of calibration.
115  */
117  (!notEnoughData || m_ForcedCalibration)) {
119  KLMChannelIndex klmChannels;
120  for (KLMChannelIndex& klmChannel : klmChannels) {
121  int subdetector = klmChannel.getSubdetector();
122  int section = klmChannel.getSection();
123  int sector = klmChannel.getSector();
124  int layer = klmChannel.getLayer();
125  int plane = klmChannel.getPlane();
126  int strip = klmChannel.getStrip();
127  uint16_t planeKLM = 0;
128  if (subdetector == KLMElementNumbers::c_BKLM) {
129  planeKLM = m_ElementNumbers->planeNumberBKLM(
130  section, sector, layer, plane);
131  } else {
132  planeKLM = m_ElementNumbers->planeNumberEKLM(
133  section, sector, layer, plane);
134  }
135  uint16_t planeIndex = m_PlaneArrayIndex->getIndex(planeKLM);
136  float efficiency = efficiencyHistogram->GetBinContent(planeIndex + 1);
137  float efficiencyError = efficiencyHistogram->GetBinError(planeIndex + 1);
138  /* Fill the efficiency for this strip. */
139  if (subdetector == KLMElementNumbers::c_BKLM) {
141  section, sector, layer, plane, strip, efficiency, efficiencyError);
142  } else {
144  section, sector, layer, plane, strip, efficiency, efficiencyError);
145  }
146  }
147  saveCalibration(m_StripEfficiency, "KLMStripEfficiency");
148  }
149  /* Write histograms to output file. */
150  TFile* outputFile = new TFile(m_OutputFileName.c_str(), "recreate");
151  outputFile->cd();
152  matchedDigitsInPlane.get()->Write();
153  allExtHitsInPlane.get()->Write();
154  efficiencyHistogram->Write();
155  delete efficiencyHistogram;
156  delete outputFile;
157  /* Set output status. */
158  if (notEnoughData && !m_ForcedCalibration)
161 }
162 
164  float* efficiency) const
165 {
166  const int nPlanes = KLMPlaneArrayIndex::Instance().getNElements();
167  int newPlanes = 0;
168  for (int i = 0; i < nPlanes; ++i) {
169  if (m_Efficiency[i] > 0 && efficiency[i] == 0)
170  newPlanes++;
171  }
172  return newPlanes;
173 }
174 
176  int* extHitsPlane) const
177 {
178  const int nPlanes = KLMPlaneArrayIndex::Instance().getNElements();
179  int newPlanes = 0;
180  for (int i = 0; i < nPlanes; ++i) {
181  if (m_ExtHitsPlane[i] > 0 && extHitsPlane[i] == 0)
182  newPlanes++;
183  }
184  return newPlanes;
185 }
Belle2::KLMStripEfficiencyAlgorithm::m_PlaneArrayIndex
const KLMPlaneArrayIndex * m_PlaneArrayIndex
Plane array index.
Definition: KLMStripEfficiencyAlgorithm.h:253
Belle2::KLMStripEfficiencyAlgorithm::Results::m_ExtHits
int m_ExtHits
Number of ExtHits (overall).
Definition: KLMStripEfficiencyAlgorithm.h:144
Belle2::KLMStripEfficiencyAlgorithm::m_OutputFileName
std::string m_OutputFileName
Output root file.
Definition: KLMStripEfficiencyAlgorithm.h:235
Belle2::CalibrationAlgorithm::saveCalibration
void saveCalibration(TClonesArray *data, const std::string &name)
Store DBArray payload with given name with default IOV.
Definition: CalibrationAlgorithm.cc:290
Belle2::KLMElementNumbers::planeNumberBKLM
uint16_t planeNumberBKLM(int section, int sector, int layer, int plane) const
Get plane number for BKLM.
Definition: KLMElementNumbers.cc:128
Belle2::KLMStripEfficiencyAlgorithm::Results::~Results
~Results()
Destructor,.
Definition: KLMStripEfficiencyAlgorithm.cc:39
Belle2::KLMStripEfficiencyAlgorithm::m_MinimalExtHits
int m_MinimalExtHits
Minimal number of ExtHits per plane.
Definition: KLMStripEfficiencyAlgorithm.h:244
Belle2::KLMPlaneArrayIndex
KLM plane array index.
Definition: KLMPlaneArrayIndex.h:33
Belle2::CalibrationAlgorithm::c_OK
@ c_OK
Finished successfuly =0 in Python.
Definition: CalibrationAlgorithm.h:51
Belle2::KLMStripEfficiencyAlgorithm::Results::m_MatchedDigits
int m_MatchedDigits
Number of matched digits.
Definition: KLMStripEfficiencyAlgorithm.h:138
Belle2::KLMChannelIndex::getSubdetector
int getSubdetector() const
Get subdetector.
Definition: KLMChannelIndex.h:116
Belle2::KLMStripEfficiencyAlgorithm::Results::m_ExtHitsPlane
int * m_ExtHitsPlane
Number of ExtHits per plane.
Definition: KLMStripEfficiencyAlgorithm.h:147
Belle2::KLMStripEfficiencyAlgorithm::m_RequestedPrecision
float m_RequestedPrecision
Requested precision of efficiency measurement.
Definition: KLMStripEfficiencyAlgorithm.h:247
Belle2::KLMStripEfficiency
DBObject used to store the efficiencies of KLM strips.
Definition: KLMStripEfficiency.h:41
Belle2::KLMChannelIndex::c_IndexLevelPlane
@ c_IndexLevelPlane
Plane.
Definition: KLMChannelIndex.h:55
Belle2::KLMPlaneArrayIndex::Instance
static const KLMPlaneArrayIndex & Instance()
Instantiation.
Definition: KLMPlaneArrayIndex.cc:28
Belle2::KLMStripEfficiency::setEndcapEfficiency
void setEndcapEfficiency(int section, int sector, int layer, int plane, int strip, float efficiency, float efficiencyError=0.)
Set efficiency and relative error for a single EKLM strip using the geometrical infos.
Definition: KLMStripEfficiency.h:97
Belle2::KLMStripEfficiencyAlgorithm::Results::Results
Results()
Constructor.
Definition: KLMStripEfficiencyAlgorithm.cc:23
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::KLMStripEfficiency::setBarrelEfficiency
void setBarrelEfficiency(int section, int sector, int layer, int plane, int strip, float efficiency, float efficiencyError=0.)
Set efficiency and relative error for a single BKLM strip using the geometrical infos.
Definition: KLMStripEfficiency.h:80
Belle2::KLMStripEfficiencyAlgorithm::m_ForcedCalibration
bool m_ForcedCalibration
Whether the calibration is forced.
Definition: KLMStripEfficiencyAlgorithm.h:238
Belle2::KLMElementNumbers::c_BKLM
@ c_BKLM
BKLM.
Definition: KLMElementNumbers.h:47
Belle2::KLMStripEfficiencyAlgorithm::m_ElementNumbers
const KLMElementNumbers * m_ElementNumbers
Element numbers.
Definition: KLMStripEfficiencyAlgorithm.h:250
Belle2::KLMStripEfficiencyAlgorithm::Results::newMeasuredPlanes
int newMeasuredPlanes(float *efficiency) const
Get number of new measured planes.
Definition: KLMStripEfficiencyAlgorithm.cc:163
Belle2::KLMStripEfficiencyAlgorithm::Results::m_Efficiency
float * m_Efficiency
Efficiency.
Definition: KLMStripEfficiencyAlgorithm.h:141
Belle2::KLMStripEfficiencyAlgorithm::m_CalibrationStage
enum CalibrationStage m_CalibrationStage
Calibration stage.
Definition: KLMStripEfficiencyAlgorithm.h:241
Belle2::KLMElementArrayIndex::getNElements
uint16_t getNElements() const
Get number of elements.
Definition: KLMElementArrayIndex.h:65
Belle2::CalibrationAlgorithm::EResult
EResult
The result of calibration.
Definition: CalibrationAlgorithm.h:50
Belle2::KLMStripEfficiencyAlgorithm::c_EfficiencyMeasurement
@ c_EfficiencyMeasurement
Measurement.
Definition: KLMStripEfficiencyAlgorithm.h:51
Belle2::KLMElementNumbers::planeNumberEKLM
uint16_t planeNumberEKLM(int section, int sector, int layer, int plane) const
Get channel number for EKLM.
Definition: KLMElementNumbers.cc:136
Belle2::KLMStripEfficiencyAlgorithm::~KLMStripEfficiencyAlgorithm
~KLMStripEfficiencyAlgorithm()
Destructor.
Definition: KLMStripEfficiencyAlgorithm.cc:58
Belle2::CalibrationAlgorithm::c_NotEnoughData
@ c_NotEnoughData
Needs more data =2 in Python.
Definition: CalibrationAlgorithm.h:53
Belle2::KLMChannelIndex
KLM channel index.
Definition: KLMChannelIndex.h:33
Belle2::CalibrationAlgorithm
Base class for calibration algorithms.
Definition: CalibrationAlgorithm.h:47
Belle2::KLMStripEfficiencyAlgorithm::calibrate
CalibrationAlgorithm::EResult calibrate() override
Calibration.
Definition: KLMStripEfficiencyAlgorithm.cc:64
Belle2::KLMStripEfficiencyAlgorithm::Results
Calibration results or supplementary results calculated from the input data.
Definition: KLMStripEfficiencyAlgorithm.h:59
Belle2::KLMStripEfficiencyAlgorithm::Results::newExtHitsPlanes
int newExtHitsPlanes(int *extHitsPlane) const
Get number of new measured planes with ExtHits.
Definition: KLMStripEfficiencyAlgorithm.cc:175
Belle2::KLMElementNumbers
KLM element numbers.
Definition: KLMElementNumbers.h:37
Belle2::KLMStripEfficiencyAlgorithm::m_StripEfficiency
KLMStripEfficiency * m_StripEfficiency
Efficiency data object.
Definition: KLMStripEfficiencyAlgorithm.h:259
Belle2::KLMStripEfficiencyAlgorithm::Results::m_AchievedPrecision
float m_AchievedPrecision
Achieved precision of efficiency measurement.
Definition: KLMStripEfficiencyAlgorithm.h:135
Belle2::KLMStripEfficiencyAlgorithm::KLMStripEfficiencyAlgorithm
KLMStripEfficiencyAlgorithm()
Constructor.
Definition: KLMStripEfficiencyAlgorithm.cc:47
Belle2::KLMStripEfficiencyAlgorithm::c_MeasurablePlaneCheck
@ c_MeasurablePlaneCheck
Check of set of planes with determined efficiency.
Definition: KLMStripEfficiencyAlgorithm.h:48
Belle2::KLMElementArrayIndex::getIndex
uint16_t getIndex(uint16_t number) const
Get element index.
Definition: KLMElementArrayIndex.cc:54
Belle2::KLMStripEfficiencyAlgorithm::m_Results
Results m_Results
Calibration results.
Definition: KLMStripEfficiencyAlgorithm.h:256