CMS 3D CMS Logo

ETLUncalibRecHitAlgo.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  tofDelay_(conf.getParameter<double>("tofDelay")),
14  timeError_(conf.getParameter<double>("timeResolutionInNs")) {}
15 
17  ~ETLUncalibRecHitAlgo() override {}
18 
20  void getEvent(const edm::Event&) final {}
21  void getEventSetup(const edm::EventSetup&) final {}
22 
24  FTLUncalibratedRecHit makeRecHit(const ETLDataFrame& dataFrame) const final;
25 
26 private:
27  const uint32_t adcNBits_;
28  const double adcSaturation_;
29  const double adcLSB_;
30  const double toaLSBToNS_;
31  const double tofDelay_;
32  const double timeError_;
33 };
34 
36  constexpr int iSample = 2; //only in-time sample
37  const auto& sample = dataFrame.sample(iSample);
38 
39  double amplitude = double(sample.data()) * adcLSB_;
40  double time = double(sample.toa()) * toaLSBToNS_ - tofDelay_;
41  unsigned char flag = 0;
42 
43  LogDebug("ETLUncalibRecHit") << "ADC+: set the charge to: " << amplitude << ' ' << sample.data() << ' ' << adcLSB_
44  << ' ' << std::endl;
45  LogDebug("ETLUncalibRecHit") << "ADC+: set the time to: " << time << ' ' << sample.toa() << ' ' << toaLSBToNS_ << ' '
46  << std::endl;
47  LogDebug("ETLUncalibRecHit") << "Final uncalibrated amplitude : " << amplitude << std::endl;
48 
49  return FTLUncalibratedRecHit(
50  dataFrame.id(), dataFrame.row(), dataFrame.column(), {amplitude, 0.f}, {time, 0.f}, timeError_, flag);
51 }
52 
#define LogDebug(id)
const int row() const
row
Definition: FTLDataFrameT.h:36
FTLUncalibratedRecHit makeRecHit(const ETLDataFrame &dataFrame) const final
make the rec hit
const S & sample(int i) const
Definition: FTLDataFrameT.h:57
ETLUncalibRecHitAlgo(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Constructor.
~ETLUncalibRecHitAlgo() override
Destructor.
void getEvent(const edm::Event &) final
get event and eventsetup information
Readout digi for HGC.
Definition: FTLDataFrameT.h:14
const int column() const
column
Definition: FTLDataFrameT.h:41
#define DEFINE_EDM_PLUGIN(factory, type, name)
#define constexpr
const D & id() const
det id
Definition: FTLDataFrameT.h:31
void getEventSetup(const edm::EventSetup &) final