CMS 3D CMS Logo

BTLUncalibRecHitAlgo.cc
Go to the documentation of this file.
4 
6 public:
10  adcNBits_(conf.getParameter<uint32_t>("adcNbits")),
11  adcSaturation_(conf.getParameter<double>("adcSaturation")),
13  toaLSBToNS_(conf.getParameter<double>("toaLSB_ns")),
14  timeError_(conf.getParameter<double>("timeResolutionInNs")),
15  timeCorr_p0_(conf.getParameter<double>("timeCorr_p0")),
16  timeCorr_p1_(conf.getParameter<double>("timeCorr_p1")),
17  timeCorr_p2_(conf.getParameter<double>("timeCorr_p2")),
18  c_LYSO_(conf.getParameter<double>("c_LYSO")) {}
19 
21  ~BTLUncalibRecHitAlgo() override {}
22 
24  void getEvent(const edm::Event&) final {}
25  void getEventSetup(const edm::EventSetup&) final {}
26 
28  FTLUncalibratedRecHit makeRecHit(const BTLDataFrame& dataFrame) const final;
29 
30 private:
31  const uint32_t adcNBits_;
32  const double adcSaturation_;
33  const double adcLSB_;
34  const double toaLSBToNS_;
35  const double timeError_;
36  const double timeCorr_p0_;
37  const double timeCorr_p1_;
38  const double timeCorr_p2_;
39  const double c_LYSO_;
40 };
41 
43  // The reconstructed amplitudes and times are saved in a std::pair
44  // BTL tile geometry (1 SiPM): only the first value of the amplitude
45  // and time pairs is used.
46  // BTL bar geometry (2 SiPMs): both values of the amplitude and
47  // time pairs are filled.
48 
49  std::pair<float, float> amplitude(0., 0.);
50  std::pair<float, float> time(0., 0.);
51 
52  unsigned char flag = 0;
53 
54  const auto& sampleLeft = dataFrame.sample(0);
55  const auto& sampleRight = dataFrame.sample(1);
56 
57  if (sampleLeft.data() > 0) {
58  amplitude.first = float(sampleLeft.data()) * adcLSB_;
59  time.first = float(sampleLeft.toa()) * toaLSBToNS_;
60 
61  // Correct the time of the left SiPM for the time-walk
63  flag |= 0x1;
64  }
65 
66  // --- If available, reconstruct the amplitude and time of the second SiPM
67  if (sampleRight.data() > 0) {
68  amplitude.second = sampleRight.data() * adcLSB_;
69  time.second = sampleRight.toa() * toaLSBToNS_;
70 
71  // Correct the time of the right SiPM for the time-walk
72  time.second -= timeCorr_p0_ * pow(amplitude.second, timeCorr_p1_) + timeCorr_p2_;
73  flag |= (0x1 << 1);
74  }
75 
76  // Calculate the position
77  // Distance from center of bar to hit
78  float position = 0.5f * (c_LYSO_ * (time.second - time.first));
79  float positionError = BTLRecHitsErrorEstimatorIM::positionError();
80 
81  LogDebug("BTLUncalibRecHit") << "ADC+: set the charge to: (" << amplitude.first << ", " << amplitude.second << ") ("
82  << sampleLeft.data() << ", " << sampleRight.data() << " " << adcLSB_ << ' '
83  << std::endl;
84  LogDebug("BTLUncalibRecHit") << "TDC+: set the time to: (" << time.first << ", " << time.second << ") ("
85  << sampleLeft.toa() << ", " << sampleRight.toa() << " " << toaLSBToNS_ << ' '
86  << std::endl;
87 
88  return FTLUncalibratedRecHit(
89  dataFrame.id(), dataFrame.row(), dataFrame.column(), amplitude, time, timeError_, position, positionError, flag);
90 }
91 
FTLDataFrameT::id
const D & id() const
det id
Definition: FTLDataFrameT.h:31
BTLUncalibRecHitAlgo::adcNBits_
const uint32_t adcNBits_
Definition: BTLUncalibRecHitAlgo.cc:31
CustomPhysics_cfi.amplitude
amplitude
Definition: CustomPhysics_cfi.py:12
MessageLogger.h
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
BTLUncalibRecHitAlgo::timeCorr_p0_
const double timeCorr_p0_
Definition: BTLUncalibRecHitAlgo.cc:36
BTLUncalibRecHitAlgo::getEvent
void getEvent(const edm::Event &) final
get event and eventsetup information
Definition: BTLUncalibRecHitAlgo.cc:24
BTLUncalibRecHitAlgo::timeCorr_p2_
const double timeCorr_p2_
Definition: BTLUncalibRecHitAlgo.cc:38
BTLUncalibRecHitAlgo::timeCorr_p1_
const double timeCorr_p1_
Definition: BTLUncalibRecHitAlgo.cc:37
FTLDataFrameT::column
const int column() const
column
Definition: FTLDataFrameT.h:41
BTLUncalibRecHitAlgo::getEventSetup
void getEventSetup(const edm::EventSetup &) final
Definition: BTLUncalibRecHitAlgo.cc:25
BTLUncalibRecHitAlgo::adcLSB_
const double adcLSB_
Definition: BTLUncalibRecHitAlgo.cc:33
MTDUncalibratedRecHitAlgoBase
Definition: MTDUncalibratedRecHitAlgoBase.h:16
MakerMacros.h
BTLUncalibRecHitAlgo
Definition: BTLUncalibRecHitAlgo.cc:5
BTLUncalibRecHitAlgo::makeRecHit
FTLUncalibratedRecHit makeRecHit(const BTLDataFrame &dataFrame) const final
make the rec hit
Definition: BTLUncalibRecHitAlgo.cc:42
MTDUncalibratedRecHitAlgoBase.h
FTLDataFrameT::sample
const S & sample(int i) const
Definition: FTLDataFrameT.h:57
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
BTLUncalibRecHitAlgo::BTLUncalibRecHitAlgo
BTLUncalibRecHitAlgo(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Constructor.
Definition: BTLUncalibRecHitAlgo.cc:8
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:223
edm::ParameterSet
Definition: ParameterSet.h:47
edmplugin::PluginFactory
Definition: PluginFactory.h:34
position
static int position[264][3]
Definition: ReadPGInfo.cc:289
BTLUncalibRecHitAlgo::adcSaturation_
const double adcSaturation_
Definition: BTLUncalibRecHitAlgo.cc:32
edm::EventSetup
Definition: EventSetup.h:57
FTLDataFrameT
Readout digi for HGC.
Definition: FTLDataFrameT.h:14
FTLDataFrameT::row
const int row() const
row
Definition: FTLDataFrameT.h:36
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
FTLUncalibratedRecHit
Definition: FTLUncalibratedRecHit.h:7
ntuplemaker.time
time
Definition: ntuplemaker.py:310
edm::Event
Definition: Event.h:73
BTLUncalibRecHitAlgo::~BTLUncalibRecHitAlgo
~BTLUncalibRecHitAlgo() override
Destructor.
Definition: BTLUncalibRecHitAlgo.cc:21
BTLRecHitsErrorEstimatorIM::positionError
static float positionError()
Definition: BTLRecHitsErrorEstimatorIM.h:42
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
BTLRecHitsErrorEstimatorIM.h
BTLUncalibRecHitAlgo::c_LYSO_
const double c_LYSO_
Definition: BTLUncalibRecHitAlgo.cc:39
RemoveAddSevLevel.flag
flag
Definition: RemoveAddSevLevel.py:116
BTLUncalibRecHitAlgo::timeError_
const double timeError_
Definition: BTLUncalibRecHitAlgo.cc:35
BTLUncalibRecHitAlgo::toaLSBToNS_
const double toaLSBToNS_
Definition: BTLUncalibRecHitAlgo.cc:34