CMS 3D CMS Logo

BTLUncalibRecHitAlgo.cc
Go to the documentation of this file.
3 
5 public:
9  adcNBits_(conf.getParameter<uint32_t>("adcNbits")),
10  adcSaturation_(conf.getParameter<double>("adcSaturation")),
12  toaLSBToNS_(conf.getParameter<double>("toaLSB_ns")),
13  timeError_(conf.getParameter<double>("timeResolutionInNs")),
14  timeCorr_p0_(conf.getParameter<double>("timeCorr_p0")),
15  timeCorr_p1_(conf.getParameter<double>("timeCorr_p1")),
16  timeCorr_p2_(conf.getParameter<double>("timeCorr_p2")) {}
17 
19  ~BTLUncalibRecHitAlgo() override {}
20 
22  void getEvent(const edm::Event&) final {}
23  void getEventSetup(const edm::EventSetup&) final {}
24 
26  FTLUncalibratedRecHit makeRecHit(const BTLDataFrame& dataFrame) const final;
27 
28 private:
29  const uint32_t adcNBits_;
30  const double adcSaturation_;
31  const double adcLSB_;
32  const double toaLSBToNS_;
33  const double timeError_;
34  const double timeCorr_p0_;
35  const double timeCorr_p1_;
36  const double timeCorr_p2_;
37 };
38 
40  // The reconstructed amplitudes and times are saved in a std::pair
41  // BTL tile geometry (1 SiPM): only the first value of the amplitude
42  // and time pairs is used.
43  // BTL bar geometry (2 SiPMs): both values of the amplitude and
44  // time pairs are filled.
45 
46  std::pair<float, float> amplitude(0., 0.);
47  std::pair<float, float> time(0., 0.);
48 
49  unsigned char flag = 0;
50 
51  const auto& sampleLeft = dataFrame.sample(0);
52  const auto& sampleRight = dataFrame.sample(1);
53 
54  if (sampleLeft.data() > 0) {
55  amplitude.first = float(sampleLeft.data()) * adcLSB_;
56  time.first = float(sampleLeft.toa()) * toaLSBToNS_;
57 
58  // Correct the time of the left SiPM for the time-walk
60  flag |= 0x1;
61  }
62 
63  // --- If available, reconstruct the amplitude and time of the second SiPM
64  if (sampleRight.data() > 0) {
65  amplitude.second = sampleRight.data() * adcLSB_;
66  time.second = sampleRight.toa() * toaLSBToNS_;
67 
68  // Correct the time of the right SiPM for the time-walk
69  time.second -= timeCorr_p0_ * pow(amplitude.second, timeCorr_p1_) + timeCorr_p2_;
70  flag |= (0x1 << 1);
71  }
72 
73  LogDebug("BTLUncalibRecHit") << "ADC+: set the charge to: (" << amplitude.first << ", " << amplitude.second << ") ("
74  << sampleLeft.data() << ", " << sampleRight.data() << " " << adcLSB_ << ' '
75  << std::endl;
76  LogDebug("BTLUncalibRecHit") << "TDC+: set the time to: (" << time.first << ", " << time.second << ") ("
77  << sampleLeft.toa() << ", " << sampleRight.toa() << " " << toaLSBToNS_ << ' '
78  << std::endl;
79 
80  return FTLUncalibratedRecHit(dataFrame.id(), dataFrame.row(), dataFrame.column(), amplitude, time, timeError_, flag);
81 }
82 
FTLDataFrameT::id
const D & id() const
det id
Definition: FTLDataFrameT.h:31
BTLUncalibRecHitAlgo::adcNBits_
const uint32_t adcNBits_
Definition: BTLUncalibRecHitAlgo.cc:29
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:34
BTLUncalibRecHitAlgo::getEvent
void getEvent(const edm::Event &) final
get event and eventsetup information
Definition: BTLUncalibRecHitAlgo.cc:22
BTLUncalibRecHitAlgo::timeCorr_p2_
const double timeCorr_p2_
Definition: BTLUncalibRecHitAlgo.cc:36
BTLUncalibRecHitAlgo::timeCorr_p1_
const double timeCorr_p1_
Definition: BTLUncalibRecHitAlgo.cc:35
FTLDataFrameT::column
const int column() const
column
Definition: FTLDataFrameT.h:41
BTLUncalibRecHitAlgo::getEventSetup
void getEventSetup(const edm::EventSetup &) final
Definition: BTLUncalibRecHitAlgo.cc:23
BTLUncalibRecHitAlgo::adcLSB_
const double adcLSB_
Definition: BTLUncalibRecHitAlgo.cc:31
MTDUncalibratedRecHitAlgoBase
Definition: MTDUncalibratedRecHitAlgoBase.h:16
MakerMacros.h
BTLUncalibRecHitAlgo
Definition: BTLUncalibRecHitAlgo.cc:4
BTLUncalibRecHitAlgo::makeRecHit
FTLUncalibratedRecHit makeRecHit(const BTLDataFrame &dataFrame) const final
make the rec hit
Definition: BTLUncalibRecHitAlgo.cc:39
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:7
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
edm::ParameterSet
Definition: ParameterSet.h:36
edmplugin::PluginFactory
Definition: PluginFactory.h:34
BTLUncalibRecHitAlgo::adcSaturation_
const double adcSaturation_
Definition: BTLUncalibRecHitAlgo.cc:30
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:30
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:19
edm::ConsumesCollector
Definition: ConsumesCollector.h:39
RemoveAddSevLevel.flag
flag
Definition: RemoveAddSevLevel.py:116
BTLUncalibRecHitAlgo::timeError_
const double timeError_
Definition: BTLUncalibRecHitAlgo.cc:33
BTLUncalibRecHitAlgo::toaLSBToNS_
const double toaLSBToNS_
Definition: BTLUncalibRecHitAlgo.cc:32