CMS 3D CMS Logo

MTDUncalibRecHitAlgo.cc
Go to the documentation of this file.
3 
5  public:
8  edm::ConsumesCollector& sumes ) :
9  MTDUncalibratedRecHitAlgoBase( conf, sumes ) {
10  uint32_t BTL_nBits = conf.getParameter<uint32_t>("BTLadcNbits");
11  double BTL_saturation = conf.getParameter<double>("BTLadcSaturation");
12  BTL_adcLSB_ = BTL_saturation/(1<<BTL_nBits);
13  BTL_toaLSBToNS_ = conf.getParameter<double>("BTLtoaLSB_ns");
14  BTL_timeError_ = conf.getParameter<double>("BTLtimeResolutionInNs");
15 
16  uint32_t ETL_nBits = conf.getParameter<uint32_t>("ETLadcNbits");
17  double ETL_saturation = conf.getParameter<double>("ETLadcSaturation");
18  ETL_adcLSB_ = ETL_saturation/(1<<ETL_nBits);
19  ETL_toaLSBToNS_ = conf.getParameter<double>("ETLtoaLSB_ns");
20  ETL_timeError_ = conf.getParameter<double>("ETLtimeResolutionInNs");
21  }
22 
24  ~MTDUncalibRecHitAlgo() 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  FTLUncalibratedRecHit makeRecHit(const ETLDataFrame& dataFrame ) const final;
33 
34  private:
37 };
38 
41  constexpr int iSample=2; //only in-time sample
42  const auto& sample = dataFrame.sample(iSample);
43 
44  double amplitude = double(sample.data()) * BTL_adcLSB_;
45  double time = double(sample.toa()) * BTL_toaLSBToNS_;
46 
47  unsigned char flag = 0;
48 
49  LogDebug("MTDUncalibRecHit") << "ADC+: set the charge to: " << amplitude << ' ' << sample.data()
50  << ' ' << BTL_adcLSB_ << ' ' << std::endl;
51  LogDebug("MTDUncalibRecHit") << "ADC+: set the time to: " << time << ' ' << sample.toa()
52  << ' ' << BTL_toaLSBToNS_ << ' ' << std::endl;
53  LogDebug("MTDUncalibRecHit") << "Final uncalibrated amplitude : " << amplitude << std::endl;
54 
55  return FTLUncalibratedRecHit( dataFrame.id(), amplitude, time, BTL_timeError_, flag);
56 }
59  constexpr int iSample=2; //only in-time sample
60  const auto& sample = dataFrame.sample(iSample);
61 
62  double amplitude = double(sample.data()) * ETL_adcLSB_;
63  double time = double(sample.toa()) * ETL_toaLSBToNS_;
64  unsigned char flag = 0;
65 
66  LogDebug("MTDUncalibRecHit") << "ADC+: set the charge to: " << amplitude << ' ' << sample.data()
67  << ' ' << ETL_adcLSB_ << ' ' << std::endl;
68  LogDebug("MTDUncalibRecHit") << "ADC+: set the time to: " << time << ' ' << sample.toa()
69  << ' ' << ETL_toaLSBToNS_ << ' ' << std::endl;
70  LogDebug("MTDUncalibRecHit") << "Final uncalibrated amplitude : " << amplitude << std::endl;
71 
72  return FTLUncalibratedRecHit( dataFrame.id(), amplitude, time, ETL_timeError_, flag);
73 }
74 
#define LogDebug(id)
T getParameter(std::string const &) const
MTDUncalibRecHitAlgo(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Constructor.
const D & id() const
det id
Definition: FTLDataFrameT.h:32
#define constexpr
void getEvent(const edm::Event &) final
get event and eventsetup information
const S & sample(int i) const
Definition: FTLDataFrameT.h:48
~MTDUncalibRecHitAlgo() override
Destructor.
void getEventSetup(const edm::EventSetup &) final
FTLUncalibratedRecHit makeRecHit(const BTLDataFrame &dataFrame) const final
make the rec hit
#define DEFINE_EDM_PLUGIN(factory, type, name)