CMS 3D CMS Logo

ETLUncalibRecHitAlgo.cc
Go to the documentation of this file.
3 
5  public:
8  edm::ConsumesCollector& sumes ) :
10  adcNBits_( conf.getParameter<uint32_t>("adcNbits") ),
11  adcSaturation_( conf.getParameter<double>("adcSaturation") ),
13  toaLSBToNS_( conf.getParameter<double>("toaLSB_ns") ),
14  tofDelay_( conf.getParameter<double>("tofDelay") ),
15  timeError_( conf.getParameter<double>("timeResolutionInNs") )
16  { }
17 
19  ~ETLUncalibRecHitAlgo() override { }
20 
22  void getEvent(const edm::Event&) final {}
23  void getEventSetup(const edm::EventSetup&) final {}
24 
26  FTLUncalibratedRecHit makeRecHit(const ETLDataFrame& dataFrame ) const final;
27 
28  private:
29 
30  const uint32_t adcNBits_;
31  const double adcSaturation_;
32  const double adcLSB_;
33  const double toaLSBToNS_;
34  const double tofDelay_;
35  const double timeError_;
36 
37 };
38 
41  constexpr int iSample=2; //only in-time sample
42  const auto& sample = dataFrame.sample(iSample);
43 
44  double amplitude = double(sample.data()) * adcLSB_;
45  double time = double(sample.toa()) * toaLSBToNS_ - tofDelay_;
46  unsigned char flag = 0;
47 
48  LogDebug("ETLUncalibRecHit") << "ADC+: set the charge to: " << amplitude << ' ' << sample.data()
49  << ' ' << adcLSB_ << ' ' << std::endl;
50  LogDebug("ETLUncalibRecHit") << "ADC+: set the time to: " << time << ' ' << sample.toa()
51  << ' ' << toaLSBToNS_ << ' ' << std::endl;
52  LogDebug("ETLUncalibRecHit") << "Final uncalibrated amplitude : " << amplitude << std::endl;
53 
54  return FTLUncalibratedRecHit( dataFrame.id(), dataFrame.row(), dataFrame.column(),
55  {amplitude, 0.f}, {time, 0.f}, timeError_, flag);
56 }
57 
#define LogDebug(id)
const int row() const
row
Definition: FTLDataFrameT.h:37
FTLUncalibratedRecHit makeRecHit(const ETLDataFrame &dataFrame) const final
make the rec hit
const S & sample(int i) const
Definition: FTLDataFrameT.h:58
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:42
#define DEFINE_EDM_PLUGIN(factory, type, name)
#define constexpr
const D & id() const
det id
Definition: FTLDataFrameT.h:32
void getEventSetup(const edm::EventSetup &) final