CMS 3D CMS Logo

ETLUncalibRecHitAlgo.cc
Go to the documentation of this file.
3 
5 
7 public:
11  adcNBits_(conf.getParameter<uint32_t>("adcNbits")),
12  adcSaturation_(conf.getParameter<double>("adcSaturation")),
14  toaLSBToNS_(conf.getParameter<double>("toaLSB_ns")),
15  tofDelay_(conf.getParameter<double>("tofDelay")),
16  timeError_(conf.getParameter<std::string>("timeResolutionInNs")) {}
17 
19  ~ETLUncalibRecHitAlgo() override {}
20 
22  void getEvent(const edm::Event&) final {}
23  void getEventSetup(const edm::EventSetup&) final {}
24 
26  FTLUncalibratedRecHit makeRecHit(const ETLDataFrame& dataFrame) const final;
27 
28 private:
29  const uint32_t adcNBits_;
30  const double adcSaturation_;
31  const double adcLSB_;
32  const double toaLSBToNS_;
33  const double tofDelay_;
35 };
36 
38  constexpr int iSample = 2; //only in-time sample
39  const auto& sample = dataFrame.sample(iSample);
40 
41  const std::array<double, 1> amplitudeV = {{double(sample.data()) * adcLSB_}};
42  // NB: Here amplitudeV is defined as an array in order to be used
43  // below as an input to FormulaEvaluator::evaluate.
44  double time = double(sample.toa()) * toaLSBToNS_ - tofDelay_;
45  unsigned char flag = 0;
46 
47  LogDebug("ETLUncalibRecHit") << "ADC+: set the charge to: " << amplitudeV[0] << ' ' << sample.data() << ' ' << adcLSB_
48  << ' ' << std::endl;
49  LogDebug("ETLUncalibRecHit") << "ADC+: set the time to: " << time << ' ' << sample.toa() << ' ' << toaLSBToNS_ << ' '
50  << std::endl;
51  LogDebug("ETLUncalibRecHit") << "Final uncalibrated amplitude : " << amplitudeV[0] << std::endl;
52 
53  const std::array<double, 1> emptyV = {{0.}};
54  double timeError = timeError_.evaluate(amplitudeV, emptyV);
55 
56  return FTLUncalibratedRecHit(dataFrame.id(),
57  dataFrame.row(),
58  dataFrame.column(),
59  {amplitudeV[0], 0.f},
60  {time, 0.f},
61  timeError,
62  -1.f,
63  -1.f,
64  flag);
65 }
66 
ETLUncalibRecHitAlgo::makeRecHit
FTLUncalibratedRecHit makeRecHit(const ETLDataFrame &dataFrame) const final
make the rec hit
Definition: ETLUncalibRecHitAlgo.cc:37
FTLDataFrameT::id
const D & id() const
det id
Definition: FTLDataFrameT.h:31
MessageLogger.h
ETLUncalibRecHitAlgo::adcLSB_
const double adcLSB_
Definition: ETLUncalibRecHitAlgo.cc:31
simplePhotonAnalyzer_cfi.sample
sample
Definition: simplePhotonAnalyzer_cfi.py:12
ETLUncalibRecHitAlgo::adcSaturation_
const double adcSaturation_
Definition: ETLUncalibRecHitAlgo.cc:30
protons_cff.time
time
Definition: protons_cff.py:39
FTLDataFrameT::column
const int column() const
column
Definition: FTLDataFrameT.h:41
ETLUncalibRecHitAlgo::adcNBits_
const uint32_t adcNBits_
Definition: ETLUncalibRecHitAlgo.cc:29
ETLUncalibRecHitAlgo
Definition: ETLUncalibRecHitAlgo.cc:6
MTDUncalibratedRecHitAlgoBase
Definition: MTDUncalibratedRecHitAlgoBase.h:16
MakerMacros.h
reco::FormulaEvaluator
Definition: FormulaEvaluator.h:67
ETLUncalibRecHitAlgo::~ETLUncalibRecHitAlgo
~ETLUncalibRecHitAlgo() override
Destructor.
Definition: ETLUncalibRecHitAlgo.cc:19
ETLUncalibRecHitAlgo::tofDelay_
const double tofDelay_
Definition: ETLUncalibRecHitAlgo.cc:33
MTDUncalibratedRecHitAlgoBase.h
FTLDataFrameT::sample
const S & sample(int i) const
Definition: FTLDataFrameT.h:57
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
ETLUncalibRecHitAlgo::timeError_
const reco::FormulaEvaluator timeError_
Definition: ETLUncalibRecHitAlgo.cc:34
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:233
edm::ParameterSet
Definition: ParameterSet.h:47
ETLUncalibRecHitAlgo::ETLUncalibRecHitAlgo
ETLUncalibRecHitAlgo(const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
Constructor.
Definition: ETLUncalibRecHitAlgo.cc:9
edmplugin::PluginFactory
Definition: PluginFactory.h:34
FormulaEvaluator.h
edm::EventSetup
Definition: EventSetup.h:58
ETLUncalibRecHitAlgo::toaLSBToNS_
const double toaLSBToNS_
Definition: ETLUncalibRecHitAlgo.cc:32
std
Definition: JetResolutionObject.h:76
FTLDataFrameT
Readout digi for HGC.
Definition: FTLDataFrameT.h:14
reco::FormulaEvaluator::evaluate
double evaluate(V const &iVariables, P const &iParameters) const
Definition: FormulaEvaluator.h:73
ETLUncalibRecHitAlgo::getEvent
void getEvent(const edm::Event &) final
get event and eventsetup information
Definition: ETLUncalibRecHitAlgo.cc:22
FTLDataFrameT::row
const int row() const
row
Definition: FTLDataFrameT.h:36
FTLUncalibratedRecHit
Definition: FTLUncalibratedRecHit.h:7
ETLUncalibRecHitAlgo::getEventSetup
void getEventSetup(const edm::EventSetup &) final
Definition: ETLUncalibRecHitAlgo.cc:23
edm::Event
Definition: Event.h:73
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
RemoveAddSevLevel.flag
flag
Definition: RemoveAddSevLevel.py:116