Belle II Software  release-05-02-19
SVDClusterCalibrations.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2017 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Giulia Casarosa *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 #include <vxd/dataobjects/VxdID.h>
14 #include <svd/dbobjects/SVDCalibrationsBase.h>
15 #include <svd/dbobjects/SVDCalibrationsScalar.h>
16 #include <svd/dbobjects/SVDClusterCuts.h>
17 #include <svd/dbobjects/SVDHitTimeSelectionFunction.h>
18 #include <framework/database/DBObjPtr.h>
19 #include <framework/logging/Logger.h>
20 #include <string>
21 
22 namespace Belle2 {
32  class SVDClusterCalibrations {
33  public:
34  static std::string name;
35  static std::string time_name;
36  typedef SVDCalibrationsBase< SVDCalibrationsScalar<SVDClusterCuts> >
37  t_payload;
38  typedef SVDCalibrationsBase< SVDCalibrationsScalar<SVDHitTimeSelectionFunction> >
45  {
46  m_aDBObjPtr.addCallback([ this ](const std::string&) -> void {
47  B2INFO("SVDClusterCuts: from now on we are using " <<
48  this->m_aDBObjPtr -> get_uniqueID()); });
49 
50  m_time_aDBObjPtr.addCallback([ this ](const std::string&) -> void {
51  B2INFO("SVDHitTimeSelectionFunction: from now on we are using " <<
52  this->m_time_aDBObjPtr -> get_uniqueID()); });
53  }
54 
68  const Belle2::VxdID& sensorID,
69  const bool& isU, const int& size,
70  const double& raw_error
71  ) const
72  {
73  return m_aDBObjPtr->getReference(sensorID.getLayerNumber(),
74  sensorID.getLadderNumber(),
75  sensorID.getSensorNumber(),
76  m_aDBObjPtr->sideIndex(isU),
77  0).getCorrectedValue(raw_error, size);
78 
79  }
80 
91  inline double getMinSeedSNR(
92  const Belle2::VxdID& sensorID,
93  const bool& isU
94  ) const
95  {
96  return m_aDBObjPtr->getReference(sensorID.getLayerNumber(),
97  sensorID.getLadderNumber(),
98  sensorID.getSensorNumber(),
99  m_aDBObjPtr->sideIndex(isU),
100  0).minSeedSNR;
101 
102  }
103 
114  inline double getMinAdjSNR(
115  const Belle2::VxdID& sensorID,
116  const bool& isU
117  ) const
118  {
119  return m_aDBObjPtr->getReference(sensorID.getLayerNumber(),
120  sensorID.getLadderNumber(),
121  sensorID.getSensorNumber(),
122  m_aDBObjPtr->sideIndex(isU),
123  0 //strip not relevant
124  ).minAdjSNR;
125 
126  }
127 
138  inline double getMinClusterSNR(
139  const Belle2::VxdID& sensorID,
140  const bool& isU
141  ) const
142  {
143  return m_aDBObjPtr->getReference(sensorID.getLayerNumber(),
144  sensorID.getLadderNumber(),
145  sensorID.getSensorNumber(),
146  m_aDBObjPtr->sideIndex(isU),
147  0 // strip not relevant
148  ).minClusterSNR;
149 
150  }
151 
152 
167  inline bool isClusterInTime(
168  const Belle2::VxdID& sensorID,
169  const bool& isU,
170  const double& svdTime, const double& svdTimeError = 0,
171  const double& t0 = 0, const double& t0Error = 0
172  ) const
173  {
174  return m_time_aDBObjPtr->getReference(sensorID.getLayerNumber(),
175  sensorID.getLadderNumber(),
176  sensorID.getSensorNumber(),
177  m_aDBObjPtr->sideIndex(isU),
178  0 //strip not relevant
179  ).isInTime(svdTime, svdTimeError, t0, t0Error);
180 
181  }
182 
193  inline bool areClusterTimesCompatible(
194  const Belle2::VxdID& sensorID,
195  const double& uTime, const double& vTime = 0
196  ) const
197  {
198  return m_time_aDBObjPtr->getReference(sensorID.getLayerNumber(),
199  sensorID.getLadderNumber(),
200  sensorID.getSensorNumber(),
201  m_aDBObjPtr->sideIndex(true), // side not relevant
202  0 // strip not relevant
203  ).areClustersInTime(uTime, vTime);
204  }
205 
206 
218  inline int getTimeSelectionFunction(
219  const Belle2::VxdID& sensorID,
220  const bool& isU
221  ) const
222  {
223  return m_time_aDBObjPtr->getReference(sensorID.getLayerNumber(),
224  sensorID.getLadderNumber(),
225  sensorID.getSensorNumber(),
226  m_aDBObjPtr->sideIndex(isU),
227  0 // strip not relevant
228  ).getFunctionID();
229 
230  }
231 
243  inline float getMinClusterTime(
244  const Belle2::VxdID& sensorID,
245  const bool& isU
246  ) const
247  {
248  return m_time_aDBObjPtr->getReference(sensorID.getLayerNumber(),
249  sensorID.getLadderNumber(),
250  sensorID.getSensorNumber(),
251  m_aDBObjPtr->sideIndex(isU),
252  0 // strip not relevant
253  ).getMinTime();
254 
255  }
256 
257 
259  TString getUniqueID() { return m_aDBObjPtr->get_uniqueID(); }
260 
262  bool isValid() { return m_aDBObjPtr.isValid(); }
263 
264 
265  private:
266 
269  };
271 }
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::SVDClusterCalibrations::time_name
static std::string time_name
name of SVDHitTimeSelectionFunction payload
Definition: SVDClusterCalibrations.h:43
Belle2::SVDClusterCalibrations::areClusterTimesCompatible
bool areClusterTimesCompatible(const Belle2::VxdID &sensorID, const double &uTime, const double &vTime=0) const
Return whether the cluster is estimated to be in time with the event or off-time.
Definition: SVDClusterCalibrations.h:201
Belle2::SVDClusterCalibrations::m_aDBObjPtr
DBObjPtr< t_payload > m_aDBObjPtr
SVDClusterCuts payload.
Definition: SVDClusterCalibrations.h:275
Belle2::SVDClusterCalibrations::getTimeSelectionFunction
int getTimeSelectionFunction(const Belle2::VxdID &sensorID, const bool &isU) const
Return the version of the function used to determine whether the cluster time is acceptable at the SP...
Definition: SVDClusterCalibrations.h:226
Belle2::SVDClusterCalibrations::getMinSeedSNR
double getMinSeedSNR(const Belle2::VxdID &sensorID, const bool &isU) const
Return the minimum SNR for the seed.
Definition: SVDClusterCalibrations.h:99
Belle2::SVDClusterCalibrations::getMinClusterSNR
double getMinClusterSNR(const Belle2::VxdID &sensorID, const bool &isU) const
Return the minimum SNR for the cluster.
Definition: SVDClusterCalibrations.h:146
Belle2::VxdID::getLadderNumber
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:108
Belle2::SVDClusterCalibrations::getUniqueID
TString getUniqueID()
returns the unique ID of the payload
Definition: SVDClusterCalibrations.h:267
Belle2::SVDClusterCalibrations::getMinClusterTime
float getMinClusterTime(const Belle2::VxdID &sensorID, const bool &isU) const
Return the min value of the cluster time to use it for reconstruction.
Definition: SVDClusterCalibrations.h:251
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::SVDClusterCalibrations::t_payload
SVDCalibrationsBase< SVDCalibrationsScalar< SVDClusterCuts > > t_payload
typedef for the of SVDClusterCuts payload of all SVD sensors
Definition: SVDClusterCalibrations.h:45
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::SVDClusterCalibrations::isValid
bool isValid()
returns true if the m_aDBObtPtr is valid in the requested IoV
Definition: SVDClusterCalibrations.h:270
Belle2::SVDClusterCalibrations::SVDClusterCalibrations
SVDClusterCalibrations()
Constructor, no input argument is required.
Definition: SVDClusterCalibrations.h:50
Belle2::SVDClusterCalibrations::getMinAdjSNR
double getMinAdjSNR(const Belle2::VxdID &sensorID, const bool &isU) const
Return the minimum SNR for the adjacent.
Definition: SVDClusterCalibrations.h:122
Belle2::VxdID::getSensorNumber
baseType getSensorNumber() const
Get the sensor id.
Definition: VxdID.h:110
Belle2::VxdID::getLayerNumber
baseType getLayerNumber() const
Get the layer id.
Definition: VxdID.h:106
Belle2::SVDClusterCalibrations::isClusterInTime
bool isClusterInTime(const Belle2::VxdID &sensorID, const bool &isU, const double &svdTime, const double &svdTimeError=0, const double &t0=0, const double &t0Error=0) const
Return whether the cluster is estimated to be in time with the event or off-time.
Definition: SVDClusterCalibrations.h:175
Belle2::SVDClusterCalibrations::m_time_aDBObjPtr
DBObjPtr< t_time_payload > m_time_aDBObjPtr
SVDHitTimeSelectionFunction paylaod.
Definition: SVDClusterCalibrations.h:276
Belle2::SVDClusterCalibrations::t_time_payload
SVDCalibrationsBase< SVDCalibrationsScalar< SVDHitTimeSelectionFunction > > t_time_payload
typedef for the of SVDHitTimeSelectionFunction payload of all SVD sensors
Definition: SVDClusterCalibrations.h:47
Belle2::SVDClusterCalibrations::name
static std::string name
name of SVDClusterCuts payload
Definition: SVDClusterCalibrations.h:42
Belle2::SVDClusterCalibrations::getCorrectedClusterPositionError
double getCorrectedClusterPositionError(const Belle2::VxdID &sensorID, const bool &isU, const int &size, const double &raw_error) const
Return the corrected cluster position error.
Definition: SVDClusterCalibrations.h:75