Belle II Software  release-05-02-19
KLMStripEfficiency.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Giacomo De Pietro *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 /* KLM headers. */
14 #include <klm/dataobjects/KLMElementNumbers.h>
15 
16 /* ROOT headers. */
17 #include <TObject.h>
18 
19 /* C++ headers. */
20 #include <cstdint>
21 #include <map>
22 
23 namespace Belle2 {
33  class KLMStripEfficiency: public TObject {
34  public:
35 
40  {
41  }
42 
47  {
48  }
49 
56  void setEfficiency(uint16_t stripId, float efficiency, float efficiencyError = 0.)
57  {
58  m_efficiency.insert(std::pair<uint16_t, float>(stripId, efficiency));
59  m_efficiencyError.insert(std::pair<uint16_t, float>(stripId, efficiencyError));
60  }
61 
72  void setBarrelEfficiency(int section, int sector, int layer, int plane, int strip, float efficiency, float efficiencyError = 0.)
73  {
74  const KLMElementNumbers* elementNumbers = &(KLMElementNumbers::Instance());
75  uint16_t stripId = elementNumbers->channelNumberBKLM(section, sector, layer, plane, strip);
76  setEfficiency(stripId, efficiency, efficiencyError);
77  }
78 
89  void setEndcapEfficiency(int section, int sector, int layer, int plane, int strip, float efficiency, float efficiencyError = 0.)
90  {
91  const KLMElementNumbers* elementNumbers = &(KLMElementNumbers::Instance());
92  uint16_t stripId = elementNumbers->channelNumberEKLM(section, sector, layer, plane, strip);
93  setEfficiency(stripId, efficiency, efficiencyError);
94  }
95 
100  float getEfficiency(uint16_t stripId) const
101  {
102  auto search = m_efficiency.find(stripId);
103  if (search == m_efficiency.end())
104  return std::numeric_limits<float>::quiet_NaN();
105  return search->second;
106  }
107 
116  float getBarrelEfficiency(int section, int sector, int layer, int plane, int strip) const
117  {
118  const KLMElementNumbers* elementNumbers = &(KLMElementNumbers::Instance());
119  uint16_t stripId = elementNumbers->channelNumberBKLM(section, sector, layer, plane, strip);
120  return getEfficiency(stripId);
121  }
122 
131  float getEndcapEfficiency(int section, int sector, int layer, int plane, int strip) const
132  {
133  const KLMElementNumbers* elementNumbers = &(KLMElementNumbers::Instance());
134  uint16_t stripId = elementNumbers->channelNumberEKLM(section, sector, layer, plane, strip);
135  return getEfficiency(stripId);
136  }
137 
142  float getEfficiencyError(uint16_t stripId) const
143  {
144  auto search = m_efficiencyError.find(stripId);
145  if (search == m_efficiencyError.end())
146  return std::numeric_limits<float>::quiet_NaN();
147  return search->second;
148  }
149 
158  float getBarrelEfficiencyError(int section, int sector, int layer, int plane, int strip) const
159  {
160  const KLMElementNumbers* elementNumbers = &(KLMElementNumbers::Instance());
161  uint16_t stripId = elementNumbers->channelNumberBKLM(section, sector, layer, plane, strip);
162  return getEfficiencyError(stripId);
163  }
164 
173  float getEndcapEfficiencyError(int section, int sector, int layer, int plane, int strip) const
174  {
175  const KLMElementNumbers* elementNumbers = &(KLMElementNumbers::Instance());
176  uint16_t stripId = elementNumbers->channelNumberEKLM(section, sector, layer, plane, strip);
177  return getEfficiencyError(stripId);
178  }
179 
180  private:
181 
183  std::map<uint16_t, float> m_efficiency;
184 
186  std::map<uint16_t, float> m_efficiencyError;
187 
190 
191  };
192 
194 }
Belle2::KLMStripEfficiency::getEfficiencyError
float getEfficiencyError(uint16_t stripId) const
Returns error on efficiency of a given KLM strip using directly the stripId.
Definition: KLMStripEfficiency.h:150
Belle2::KLMStripEfficiency::getEfficiency
float getEfficiency(uint16_t stripId) const
Returns efficiency of a given KLM strip using directly the stripId.
Definition: KLMStripEfficiency.h:108
Belle2::KLMStripEfficiency::KLMStripEfficiency
KLMStripEfficiency()
Default constructor.
Definition: KLMStripEfficiency.h:47
Belle2::KLMStripEfficiency::m_efficiency
std::map< uint16_t, float > m_efficiency
KLM strip efficiency.
Definition: KLMStripEfficiency.h:191
Belle2::KLMElementNumbers::Instance
static const KLMElementNumbers & Instance()
Instantiation.
Definition: KLMElementNumbers.cc:31
Belle2::KLMStripEfficiency::getBarrelEfficiencyError
float getBarrelEfficiencyError(int section, int sector, int layer, int plane, int strip) const
Returns error on efficiency of a given BKLM strip using the geometrical infos.
Definition: KLMStripEfficiency.h:166
Belle2::KLMStripEfficiency::ClassDef
ClassDef(KLMStripEfficiency, 1)
Class version.
Belle2::KLMStripEfficiency
DBObject used to store the efficiencies of KLM strips.
Definition: KLMStripEfficiency.h:41
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
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::KLMStripEfficiency::getEndcapEfficiency
float getEndcapEfficiency(int section, int sector, int layer, int plane, int strip) const
Returns efficiency of a given EKLM strip using the geometrical infos.
Definition: KLMStripEfficiency.h:139
Belle2::KLMStripEfficiency::setEfficiency
void setEfficiency(uint16_t stripId, float efficiency, float efficiencyError=0.)
Set efficiency and relative error for a single KLM strip using directly the stripId.
Definition: KLMStripEfficiency.h:64
Belle2::KLMStripEfficiency::getEndcapEfficiencyError
float getEndcapEfficiencyError(int section, int sector, int layer, int plane, int strip) const
Returns error on efficiency of a given EKLM strip using the geometrical infos.
Definition: KLMStripEfficiency.h:181
Belle2::KLMStripEfficiency::getBarrelEfficiency
float getBarrelEfficiency(int section, int sector, int layer, int plane, int strip) const
Returns efficiency of a given BKLM strip using the geometrical infos.
Definition: KLMStripEfficiency.h:124
Belle2::KLMStripEfficiency::~KLMStripEfficiency
~KLMStripEfficiency()
Default destructor.
Definition: KLMStripEfficiency.h:54
Belle2::KLMStripEfficiency::m_efficiencyError
std::map< uint16_t, float > m_efficiencyError
KLM strip efficiency error.
Definition: KLMStripEfficiency.h:194