CMS 3D CMS Logo

BTLUncalibRecHitAlgo.cc
Go to the documentation of this file.
2 
5 
7 
9 public:
13  adcNBits_(conf.getParameter<uint32_t>("adcNbits")),
14  adcSaturation_(conf.getParameter<double>("adcSaturation")),
16  toaLSBToNS_(conf.getParameter<double>("toaLSB_ns")),
17  timeError_(conf.getParameter<std::string>("timeResolutionInNs")),
18  timeCorr_p0_(conf.getParameter<double>("timeCorr_p0")),
19  timeCorr_p1_(conf.getParameter<double>("timeCorr_p1")),
20  timeCorr_p2_(conf.getParameter<double>("timeCorr_p2")),
21  c_LYSO_(conf.getParameter<double>("c_LYSO")) {}
22 
24  ~BTLUncalibRecHitAlgo() override {}
25 
27  void getEvent(const edm::Event&) final {}
28  void getEventSetup(const edm::EventSetup&) final {}
29 
31  FTLUncalibratedRecHit makeRecHit(const BTLDataFrame& dataFrame) const final;
32 
33 private:
34  const uint32_t adcNBits_;
35  const double adcSaturation_;
36  const double adcLSB_;
37  const double toaLSBToNS_;
39  const double timeCorr_p0_;
40  const double timeCorr_p1_;
41  const double timeCorr_p2_;
42  const double c_LYSO_;
43 };
44 
46  // The reconstructed amplitudes and times are saved in a std::pair
47  // BTL tile geometry (1 SiPM): only the first value of the amplitude
48  // and time pairs is used.
49  // BTL bar geometry (2 SiPMs): both values of the amplitude and
50  // time pairs are filled.
51 
52  std::pair<float, float> amplitude(0., 0.);
53  std::pair<float, float> time(0., 0.);
54 
55  unsigned char flag = 0;
56 
57  const auto& sampleRight = dataFrame.sample(0);
58  const auto& sampleLeft = dataFrame.sample(1);
59 
60  double nHits = 0.;
61 
62  LogDebug("BTLUncalibRecHit") << "Original input time t1,t2 " << float(sampleRight.toa()) * toaLSBToNS_ << ", "
63  << float(sampleLeft.toa()) * toaLSBToNS_ << std::endl;
64 
65  if (sampleRight.data() > 0) {
66  amplitude.first = float(sampleRight.data()) * adcLSB_;
67  time.first = float(sampleRight.toa()) * toaLSBToNS_;
68 
69  nHits += 1.;
70 
71  // Correct the time of the left SiPM for the time-walk
73  flag |= 0x1;
74  }
75 
76  // --- If available, reconstruct the amplitude and time of the second SiPM
77  if (sampleLeft.data() > 0) {
78  amplitude.second = float(sampleLeft.data()) * adcLSB_;
79  time.second = float(sampleLeft.toa()) * toaLSBToNS_;
80 
81  nHits += 1.;
82 
83  // Correct the time of the right SiPM for the time-walk
84  time.second -= timeCorr_p0_ * pow(amplitude.second, timeCorr_p1_) + timeCorr_p2_;
85  flag |= (0x1 << 1);
86  }
87 
88  // --- Calculate the error on the hit time using the provided parameterization
89 
90  const std::array<double, 1> amplitudeV = {{(amplitude.first + amplitude.second) / nHits}};
91  const std::array<double, 1> emptyV = {{0.}};
92 
93  double timeError = (nHits > 0. ? timeError_.evaluate(amplitudeV, emptyV) : -1.);
94 
95  // Calculate the position
96  // Distance from center of bar to hit
97 
98  float position = 0.5f * (c_LYSO_ * (time.second - time.first));
99  float positionError = BTLRecHitsErrorEstimatorIM::positionError();
100 
101  LogDebug("BTLUncalibRecHit") << "DetId: " << dataFrame.id().rawId() << " x position = " << position << " +/- "
102  << positionError;
103  LogDebug("BTLUncalibRecHit") << "ADC+: set the charge to: (" << amplitude.first << ", " << amplitude.second << ") ("
104  << sampleRight.data() << ", " << sampleLeft.data() << ") " << adcLSB_ << ' '
105  << std::endl;
106  LogDebug("BTLUncalibRecHit") << "TDC+: set the time to: (" << time.first << ", " << time.second << ") ("
107  << sampleRight.toa() << ", " << sampleLeft.toa() << ") " << toaLSBToNS_ << ' '
108  << std::endl;
109 
110  return FTLUncalibratedRecHit(
111  dataFrame.id(), dataFrame.row(), dataFrame.column(), amplitude, time, timeError, position, positionError, flag);
112 }
113 
double evaluate(V const &iVariables, P const &iParameters) const
~BTLUncalibRecHitAlgo() override
Destructor.
BTLUncalibRecHitAlgo(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Constructor.
const D & id() const
det id
Definition: FTLDataFrameT.h:31
void getEvent(const edm::Event &) final
get event and eventsetup information
const S & sample(int i) const
Definition: FTLDataFrameT.h:57
const reco::FormulaEvaluator timeError_
const int column() const
column
Definition: FTLDataFrameT.h:41
FTLUncalibratedRecHit makeRecHit(const BTLDataFrame &dataFrame) const final
make the rec hit
static int position[264][3]
Definition: ReadPGInfo.cc:289
Readout digi for HGC.
Definition: FTLDataFrameT.h:14
void getEventSetup(const edm::EventSetup &) final
#define DEFINE_EDM_PLUGIN(factory, type, name)
TupleMultiplicity< TrackerTraits > const *__restrict__ uint32_t nHits
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
const int row() const
row
Definition: FTLDataFrameT.h:36
#define LogDebug(id)