CMS 3D CMS Logo

ETLUncalibRecHitAlgo.cc
Go to the documentation of this file.
3 
5 
7 public:
11  adcNBits_(conf.getParameter<uint32_t>("adcNbits")),
12  adcSaturation_(conf.getParameter<double>("adcSaturation")),
14  toaLSBToNS_(conf.getParameter<double>("toaLSB_ns")),
15  tofDelay_(conf.getParameter<double>("tofDelay")),
16  timeError_(conf.getParameter<std::string>("timeResolutionInNs")),
17  timeCorr_p0_(conf.getParameter<double>("timeCorr_p0")),
18  timeCorr_p1_(conf.getParameter<double>("timeCorr_p1")),
19  timeCorr_p2_(conf.getParameter<double>("timeCorr_p2")),
20  timeCorr_p3_(conf.getParameter<double>("timeCorr_p3")) {}
22  ~ETLUncalibRecHitAlgo() override {}
23 
25  void getEvent(const edm::Event&) final {}
26  void getEventSetup(const edm::EventSetup&) final {}
27 
29  FTLUncalibratedRecHit makeRecHit(const ETLDataFrame& dataFrame) const final;
30 
31 private:
32  const uint32_t adcNBits_;
33  const double adcSaturation_;
34  const double adcLSB_;
35  const double toaLSBToNS_;
36  const double tofDelay_;
38  const double timeCorr_p0_;
39  const double timeCorr_p1_;
40  const double timeCorr_p2_;
41  const double timeCorr_p3_;
42 };
43 
45  constexpr int iSample = 2; //only in-time sample
46  const auto& sample = dataFrame.sample(iSample);
47 
48  double time = double(sample.toa()) * toaLSBToNS_ - tofDelay_;
49  double time_over_threshold = double(sample.tot()) * toaLSBToNS_;
50  const std::array<double, 1> time_over_threshold_V = {{time_over_threshold}};
51 
52  unsigned char flag = 0;
53 
54  LogDebug("ETLUncalibRecHit") << "ADC+: set the charge to: " << time_over_threshold << ' ' << sample.tot() << ' '
55  << toaLSBToNS_;
56 
57  if (time_over_threshold == 0) {
58  LogDebug("ETLUncalibRecHit") << "ADC+: set the time to: " << time << ' ' << sample.toa() << ' ' << toaLSBToNS_;
59 
60  } else {
61  // Time-walk correction for toa
62  double timeWalkCorr = timeCorr_p0_ + timeCorr_p1_ * time_over_threshold +
63  timeCorr_p2_ * time_over_threshold * time_over_threshold +
64  timeCorr_p3_ * time_over_threshold * time_over_threshold * time_over_threshold;
65 
66  time -= timeWalkCorr;
67 
68  LogDebug("ETLUncalibRecHit") << "ADC+: set the time to: " << time << ' ' << sample.toa() << ' ' << toaLSBToNS_
69  << " .Timewalk correction: " << timeWalkCorr;
70  }
71 
72  LogDebug("ETLUncalibRecHit") << "Final uncalibrated time_over_threshold: " << time_over_threshold;
73 
74  const std::array<double, 1> emptyV = {{0.}};
75 
76  double timeError = timeError_.evaluate(time_over_threshold_V, emptyV);
77 
78  return FTLUncalibratedRecHit(dataFrame.id(),
79  dataFrame.row(),
80  dataFrame.column(),
81  {time_over_threshold_V[0], 0.f},
82  {time, 0.f},
83  timeError,
84  -1.f,
85  -1.f,
86  flag);
87 }
double evaluate(V const &iVariables, P const &iParameters) const
FTLUncalibratedRecHit makeRecHit(const ETLDataFrame &dataFrame) const final
make the rec hit
const reco::FormulaEvaluator timeError_
ETLUncalibRecHitAlgo(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Constructor.
const D & id() const
det id
Definition: FTLDataFrameT.h:30
~ETLUncalibRecHitAlgo() override
Destructor.
const S & sample(int i) const
Definition: FTLDataFrameT.h:56
const int column() const
column
Definition: FTLDataFrameT.h:40
void getEvent(const edm::Event &) final
get event and eventsetup information
Readout digi for HGC.
Definition: FTLDataFrameT.h:14
#define DEFINE_EDM_PLUGIN(factory, type, name)
const int row() const
row
Definition: FTLDataFrameT.h:35
#define LogDebug(id)
void getEventSetup(const edm::EventSetup &) final