CMS 3D CMS Logo

BTLUncalibRecHitAlgo.cc
Go to the documentation of this file.
4 
6 
8 public:
12  adcNBits_(conf.getParameter<uint32_t>("adcNbits")),
13  adcSaturation_(conf.getParameter<double>("adcSaturation")),
15  toaLSBToNS_(conf.getParameter<double>("toaLSB_ns")),
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  c_LYSO_(conf.getParameter<double>("c_LYSO")) {}
21 
23  ~BTLUncalibRecHitAlgo() override {}
24 
26  void getEvent(const edm::Event&) final {}
27  void getEventSetup(const edm::EventSetup&) final {}
28 
30  FTLUncalibratedRecHit makeRecHit(const BTLDataFrame& dataFrame) const final;
31 
32 private:
33  const uint32_t adcNBits_;
34  const double adcSaturation_;
35  const double adcLSB_;
36  const double toaLSBToNS_;
38  const double timeCorr_p0_;
39  const double timeCorr_p1_;
40  const double timeCorr_p2_;
41  const double c_LYSO_;
42 };
43 
45  // The reconstructed amplitudes and times are saved in a std::pair
46  // BTL tile geometry (1 SiPM): only the first value of the amplitude
47  // and time pairs is used.
48  // BTL bar geometry (2 SiPMs): both values of the amplitude and
49  // time pairs are filled.
50 
51  std::pair<float, float> amplitude(0., 0.);
52  std::pair<float, float> time(0., 0.);
53 
54  unsigned char flag = 0;
55 
56  const auto& sampleRight = dataFrame.sample(0);
57  const auto& sampleLeft = dataFrame.sample(1);
58 
59  double nHits = 0.;
60 
61  if (sampleRight.data() > 0) {
62  amplitude.first = float(sampleRight.data()) * adcLSB_;
63  time.first = float(sampleRight.toa()) * toaLSBToNS_;
64 
65  nHits += 1.;
66 
67  // Correct the time of the left SiPM for the time-walk
69  flag |= 0x1;
70  }
71 
72  // --- If available, reconstruct the amplitude and time of the second SiPM
73  if (sampleLeft.data() > 0) {
74  amplitude.second = float(sampleLeft.data()) * adcLSB_;
75  time.second = float(sampleLeft.toa()) * toaLSBToNS_;
76 
77  nHits += 1.;
78 
79  // Correct the time of the right SiPM for the time-walk
80  time.second -= timeCorr_p0_ * pow(amplitude.second, timeCorr_p1_) + timeCorr_p2_;
81  flag |= (0x1 << 1);
82  }
83 
84  // --- Calculate the error on the hit time using the provided parameterization
85 
86  const std::array<double, 1> amplitudeV = {{(amplitude.first + amplitude.second) / nHits}};
87  const std::array<double, 1> emptyV = {{0.}};
88 
89  double timeError = (nHits > 0. ? timeError_.evaluate(amplitudeV, emptyV) : -1.);
90 
91  // Calculate the position
92  // Distance from center of bar to hit
93  float position = 0.5f * (c_LYSO_ * (time.second - time.first));
94  float positionError = BTLRecHitsErrorEstimatorIM::positionError();
95 
96  LogDebug("BTLUncalibRecHit") << "ADC+: set the charge to: (" << amplitude.first << ", " << amplitude.second << ") ("
97  << sampleRight.data() << ", " << sampleLeft.data() << " " << adcLSB_ << ' '
98  << std::endl;
99  LogDebug("BTLUncalibRecHit") << "TDC+: set the time to: (" << time.first << ", " << time.second << ") ("
100  << sampleRight.toa() << ", " << sampleLeft.toa() << " " << toaLSBToNS_ << ' '
101  << std::endl;
102 
103  return FTLUncalibratedRecHit(
104  dataFrame.id(), dataFrame.row(), dataFrame.column(), amplitude, time, timeError, position, positionError, flag);
105 }
106 
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)