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  constexpr int iSample=2; //only in-time sample
46  const auto& sample = dataFrame.sample(iSample);
47 
48  double amplitude = double(sample.data()) * adcLSB_;
49  double time = double(sample.toa()) * toaLSBToNS_;
50 
51  // --- Correct the time for the time-walk and the constant delays
52  if ( amplitude > 0. )
53  time -= timeCorr_p0_*pow(amplitude,timeCorr_p1_) + timeCorr_p2_;
54 
55  unsigned char flag = 0;
56 
57  LogDebug("BTLUncalibRecHit") << "ADC+: set the charge to: " << amplitude << ' ' << sample.data()
58  << ' ' << adcLSB_ << ' ' << std::endl;
59  LogDebug("BTLUncalibRecHit") << "ADC+: set the time to: " << time << ' ' << sample.toa()
60  << ' ' << toaLSBToNS_ << ' ' << std::endl;
61  LogDebug("BTLUncalibRecHit") << "Final uncalibrated amplitude : " << amplitude << std::endl;
62 
63  return FTLUncalibratedRecHit( dataFrame.id(), amplitude, time, timeError_, flag);
64 }
65 
#define LogDebug(id)
const D & id() const
det id
Definition: FTLDataFrameT.h:32
~BTLUncalibRecHitAlgo() override
Destructor.
BTLUncalibRecHitAlgo(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Constructor.
#define constexpr
FTLUncalibratedRecHit makeRecHit(const BTLDataFrame &dataFrame) const final
make the rec hit
void getEvent(const edm::Event &) final
get event and eventsetup information
const S & sample(int i) const
Definition: FTLDataFrameT.h:48
void getEventSetup(const edm::EventSetup &) final
#define DEFINE_EDM_PLUGIN(factory, type, name)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40