CMS 3D CMS Logo

BTLUncalibRecHitAlgo.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  timeError_( conf.getParameter<double>("timeResolutionInNs") ),
15  timeCorr_p0_( conf.getParameter<double>("timeCorr_p0") ),
16  timeCorr_p1_( conf.getParameter<double>("timeCorr_p1") ),
17  timeCorr_p2_( conf.getParameter<double>("timeCorr_p2") )
18  { }
19 
21  ~BTLUncalibRecHitAlgo() override { }
22 
24  void getEvent(const edm::Event&) final {}
25  void getEventSetup(const edm::EventSetup&) final {}
26 
28  FTLUncalibratedRecHit makeRecHit(const BTLDataFrame& dataFrame ) const final;
29 
30  private:
31 
32  const uint32_t adcNBits_;
33  const double adcSaturation_;
34  const double adcLSB_;
35  const double toaLSBToNS_;
36  const double timeError_;
37  const double timeCorr_p0_;
38  const double timeCorr_p1_;
39  const double timeCorr_p2_;
40 
41 };
42 
45 
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& sampleLeft = dataFrame.sample(0);
58  const auto& sampleRight = dataFrame.sample(1);
59 
60  if ( sampleLeft.data() > 0 ) {
61 
62  amplitude.first = float(sampleLeft.data()) * adcLSB_;
63  time.first = float(sampleLeft.toa()) * toaLSBToNS_;
64 
65  // Correct the time of the left SiPM for the time-walk
66  time.first -= timeCorr_p0_*pow(amplitude.first,timeCorr_p1_) + timeCorr_p2_;
67  flag |= 0x1;
68 
69  }
70 
71  // --- If available, reconstruct the amplitude and time of the second SiPM
72  if ( sampleRight.data() > 0 ) {
73 
74  amplitude.second = sampleRight.data() * adcLSB_;
75  time.second = sampleRight.toa() * toaLSBToNS_;
76 
77  // Correct the time of the right SiPM for the time-walk
78  time.second -= timeCorr_p0_*pow(amplitude.second,timeCorr_p1_) + timeCorr_p2_;
79  flag |= (0x1 << 1);
80 
81  }
82 
83  LogDebug("BTLUncalibRecHit") << "ADC+: set the charge to: (" << amplitude.first << ", "
84  << amplitude.second << ") ("
85  << sampleLeft.data() << ", " << sampleRight.data()
86  << " " << adcLSB_ << ' ' << std::endl;
87  LogDebug("BTLUncalibRecHit") << "TDC+: set the time to: (" << time.first << ", "
88  << time.second << ") ("
89  << sampleLeft.toa() << ", " << sampleRight.toa()
90  << " " << toaLSBToNS_ << ' ' << std::endl;
91 
92  return FTLUncalibratedRecHit( dataFrame.id(), dataFrame.row(), dataFrame.column(),
94 }
95 
#define LogDebug(id)
const int row() const
row
Definition: FTLDataFrameT.h:37
~BTLUncalibRecHitAlgo() override
Destructor.
BTLUncalibRecHitAlgo(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Constructor.
const S & sample(int i) const
Definition: FTLDataFrameT.h:58
FTLUncalibratedRecHit makeRecHit(const BTLDataFrame &dataFrame) const final
make the rec hit
void getEvent(const edm::Event &) final
get event and eventsetup information
Readout digi for HGC.
Definition: FTLDataFrameT.h:14
void getEventSetup(const edm::EventSetup &) final
const int column() const
column
Definition: FTLDataFrameT.h:42
#define DEFINE_EDM_PLUGIN(factory, type, name)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
const D & id() const
det id
Definition: FTLDataFrameT.h:32