Belle II Software  release-05-02-19
SVDCoGTimeCalibrations.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/SVDCoGCalibrationFunction.h>
17 #include <framework/database/DBObjPtr.h>
18 #include <framework/logging/Logger.h>
19 #include <string>
20 
21 namespace Belle2 {
32  class SVDCoGTimeCalibrations {
33  public:
34  static std::string name;
35  typedef SVDCalibrationsBase< SVDCalibrationsScalar< SVDCoGCalibrationFunction > >
36  t_payload;
40  {
41  m_aDBObjPtr.addCallback([ this ](const std::string&) -> void {
42  B2INFO("SVDCoGTimeCalibrations: from now on we are using " <<
43  this->m_aDBObjPtr -> get_uniqueID()); });
44  }
45 
58  inline double getCorrectedTime(
59  const Belle2::VxdID& sensorID,
60  const bool& isU, const unsigned short& strip,
61  const double& raw_time,
62  const int& bin
63  ) const
64  {
65  return m_aDBObjPtr->getReference(sensorID.getLayerNumber(),
66  sensorID.getLadderNumber(),
67  sensorID.getSensorNumber(),
68  m_aDBObjPtr->sideIndex(isU),
69  strip).calibratedValue(raw_time, bin);
70 
71  }
72 
87  inline double getCorrectedTimeError(
88  const Belle2::VxdID& sensorID,
89  const bool& isU, const unsigned short& strip,
90  const double& raw_time,
91  const double& raw_timeErr,
92  const int& bin
93  ) const
94  {
95  return m_aDBObjPtr->get(sensorID.getLayerNumber(),
96  sensorID.getLadderNumber(),
97  sensorID.getSensorNumber(),
98  m_aDBObjPtr->sideIndex(isU),
99  strip).calibratedValueError(raw_time, raw_timeErr, bin);
100 
101  }
102 
104  TString getUniqueID() { return m_aDBObjPtr->get_uniqueID(); }
105 
107  bool isValid() { return m_aDBObjPtr.isValid(); }
108 
109 
110  private:
111 
113  };
115 }
Belle2::SVDCoGTimeCalibrations::SVDCoGTimeCalibrations
SVDCoGTimeCalibrations()
Constructor, no input argument is required.
Definition: SVDCoGTimeCalibrations.h:47
Belle2::SVDCoGTimeCalibrations::name
static std::string name
name of the SVDCoGCalibrationFunction payload
Definition: SVDCoGTimeCalibrations.h:42
Belle2::VxdID
Class to uniquely identify a any structure of the PXD and SVD.
Definition: VxdID.h:43
Belle2::SVDCoGTimeCalibrations::m_aDBObjPtr
DBObjPtr< t_payload > m_aDBObjPtr
SVDCoGCalibrationFunction payload.
Definition: SVDCoGTimeCalibrations.h:120
Belle2::SVDCoGTimeCalibrations::getUniqueID
TString getUniqueID()
returns the unique ID of the payload
Definition: SVDCoGTimeCalibrations.h:112
Belle2::VxdID::getLadderNumber
baseType getLadderNumber() const
Get the ladder id.
Definition: VxdID.h:108
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::SVDCoGTimeCalibrations::getCorrectedTimeError
double getCorrectedTimeError(const Belle2::VxdID &sensorID, const bool &isU, const unsigned short &strip, const double &raw_time, const double &raw_timeErr, const int &bin) const
Return the strip time error, given the raw strip time, and tje raw time error.
Definition: SVDCoGTimeCalibrations.h:95
Belle2::SVDCoGTimeCalibrations::t_payload
SVDCalibrationsBase< SVDCalibrationsScalar< SVDCoGCalibrationFunction > > t_payload
typedef for the SVDCoGCalibrationFunction payload of all SVD sensors
Definition: SVDCoGTimeCalibrations.h:44
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::SVDCoGTimeCalibrations::getCorrectedTime
double getCorrectedTime(const Belle2::VxdID &sensorID, const bool &isU, const unsigned short &strip, const double &raw_time, const int &bin) const
Return the strip time, given the raw strip time.
Definition: SVDCoGTimeCalibrations.h:66
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::SVDCoGTimeCalibrations::isValid
bool isValid()
returns true if the m_aDBObtPtr is valid in the requested IoV
Definition: SVDCoGTimeCalibrations.h:115