CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
BTLUncalibRecHitAlgo Class Reference
Inheritance diagram for BTLUncalibRecHitAlgo:
MTDUncalibratedRecHitAlgoBase< DataFrame >

Public Member Functions

 BTLUncalibRecHitAlgo (const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
 Constructor. More...
 
void getEvent (const edm::Event &) final
 get event and eventsetup information More...
 
void getEventSetup (const edm::EventSetup &) final
 
FTLUncalibratedRecHit makeRecHit (const BTLDataFrame &dataFrame) const final
 make the rec hit More...
 
 ~BTLUncalibRecHitAlgo () override
 Destructor. More...
 
- Public Member Functions inherited from MTDUncalibratedRecHitAlgoBase< DataFrame >
virtual FTLUncalibratedRecHit makeRecHit (const DataFrame &dataFrame) const =0
 make the rec hit More...
 
 MTDUncalibratedRecHitAlgoBase (const edm::ParameterSet &conf, edm::ConsumesCollector &sumes)
 Constructor. More...
 
const std::string & name () const
 
virtual ~MTDUncalibratedRecHitAlgoBase ()
 Destructor. More...
 

Private Attributes

const double adcLSB_
 
const uint32_t adcNBits_
 
const double adcSaturation_
 
const double c_LYSO_
 
const double timeCorr_p0_
 
const double timeCorr_p1_
 
const double timeCorr_p2_
 
const reco::FormulaEvaluator timeError_
 
const double toaLSBToNS_
 

Detailed Description

Definition at line 8 of file BTLUncalibRecHitAlgo.cc.

Constructor & Destructor Documentation

◆ BTLUncalibRecHitAlgo()

BTLUncalibRecHitAlgo::BTLUncalibRecHitAlgo ( const edm::ParameterSet conf,
edm::ConsumesCollector sumes 
)
inline

Constructor.

Definition at line 11 of file BTLUncalibRecHitAlgo.cc.

13  adcNBits_(conf.getParameter<uint32_t>("adcNbits")),
14  adcSaturation_(conf.getParameter<double>("adcSaturation")),
16  toaLSBToNS_(conf.getParameter<double>("toaLSB_ns")),
17  timeError_(conf.getParameter<std::string>("timeResolutionInNs")),
18  timeCorr_p0_(conf.getParameter<double>("timeCorr_p0")),
19  timeCorr_p1_(conf.getParameter<double>("timeCorr_p1")),
20  timeCorr_p2_(conf.getParameter<double>("timeCorr_p2")),
21  c_LYSO_(conf.getParameter<double>("c_LYSO")) {}
const reco::FormulaEvaluator timeError_

◆ ~BTLUncalibRecHitAlgo()

BTLUncalibRecHitAlgo::~BTLUncalibRecHitAlgo ( )
inlineoverride

Destructor.

Definition at line 24 of file BTLUncalibRecHitAlgo.cc.

24 {}

Member Function Documentation

◆ getEvent()

void BTLUncalibRecHitAlgo::getEvent ( const edm::Event )
inlinefinalvirtual

get event and eventsetup information

Implements MTDUncalibratedRecHitAlgoBase< DataFrame >.

Definition at line 27 of file BTLUncalibRecHitAlgo.cc.

27 {}

◆ getEventSetup()

void BTLUncalibRecHitAlgo::getEventSetup ( const edm::EventSetup )
inlinefinalvirtual

Implements MTDUncalibratedRecHitAlgoBase< DataFrame >.

Definition at line 28 of file BTLUncalibRecHitAlgo.cc.

28 {}

◆ makeRecHit()

FTLUncalibratedRecHit BTLUncalibRecHitAlgo::makeRecHit ( const BTLDataFrame dataFrame) const
final

make the rec hit

Definition at line 45 of file BTLUncalibRecHitAlgo.cc.

References adcLSB_, CustomPhysics_cfi::amplitude, c_LYSO_, FTLDataFrameT< D, S, DECODE >::column(), reco::FormulaEvaluator::evaluate(), RemoveAddSevLevel::flag, nano_mu_digi_cff::float, FTLDataFrameT< D, S, DECODE >::id(), LogDebug, nHits, position, BTLRecHitsErrorEstimatorIM::positionError(), funct::pow(), FTLDataFrameT< D, S, DECODE >::row(), FTLDataFrameT< D, S, DECODE >::sample(), hcalRecHitTable_cff::time, timeCorr_p0_, timeCorr_p1_, timeCorr_p2_, timeError_, and toaLSBToNS_.

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& sampleRight = dataFrame.sample(0);
58  const auto& sampleLeft = dataFrame.sample(1);
59 
60  double nHits = 0.;
61 
62  LogDebug("BTLUncalibRecHit") << "Original input time t1,t2 " << float(sampleRight.toa()) * toaLSBToNS_ << ", "
63  << float(sampleLeft.toa()) * toaLSBToNS_ << std::endl;
64 
65  if (sampleRight.data() > 0) {
66  amplitude.first = float(sampleRight.data()) * adcLSB_;
67  time.first = float(sampleRight.toa()) * toaLSBToNS_;
68 
69  nHits += 1.;
70 
71  // Correct the time of the left SiPM for the time-walk
73  flag |= 0x1;
74  }
75 
76  // --- If available, reconstruct the amplitude and time of the second SiPM
77  if (sampleLeft.data() > 0) {
78  amplitude.second = float(sampleLeft.data()) * adcLSB_;
79  time.second = float(sampleLeft.toa()) * toaLSBToNS_;
80 
81  nHits += 1.;
82 
83  // Correct the time of the right SiPM for the time-walk
84  time.second -= timeCorr_p0_ * pow(amplitude.second, timeCorr_p1_) + timeCorr_p2_;
85  flag |= (0x1 << 1);
86  }
87 
88  // --- Calculate the error on the hit time using the provided parameterization
89 
90  const std::array<double, 1> amplitudeV = {{(amplitude.first + amplitude.second) / nHits}};
91  const std::array<double, 1> emptyV = {{0.}};
92 
93  double timeError = (nHits > 0. ? timeError_.evaluate(amplitudeV, emptyV) : -1.);
94 
95  // Calculate the position
96  // Distance from center of bar to hit
97 
98  float position = 0.5f * (c_LYSO_ * (time.second - time.first));
99  float positionError = BTLRecHitsErrorEstimatorIM::positionError();
100 
101  LogDebug("BTLUncalibRecHit") << "DetId: " << dataFrame.id().rawId() << " x position = " << position << " +/- "
102  << positionError;
103  LogDebug("BTLUncalibRecHit") << "ADC+: set the charge to: (" << amplitude.first << ", " << amplitude.second << ") ("
104  << sampleRight.data() << ", " << sampleLeft.data() << ") " << adcLSB_ << ' '
105  << std::endl;
106  LogDebug("BTLUncalibRecHit") << "TDC+: set the time to: (" << time.first << ", " << time.second << ") ("
107  << sampleRight.toa() << ", " << sampleLeft.toa() << ") " << toaLSBToNS_ << ' '
108  << std::endl;
109 
110  return FTLUncalibratedRecHit(
111  dataFrame.id(), dataFrame.row(), dataFrame.column(), amplitude, time, timeError, position, positionError, flag);
112 }
double evaluate(V const &iVariables, P const &iParameters) const
const D & id() const
det id
Definition: FTLDataFrameT.h:30
const S & sample(int i) const
Definition: FTLDataFrameT.h:56
const reco::FormulaEvaluator timeError_
const int column() const
column
Definition: FTLDataFrameT.h:40
static int position[264][3]
Definition: ReadPGInfo.cc:289
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:35
#define LogDebug(id)

Member Data Documentation

◆ adcLSB_

const double BTLUncalibRecHitAlgo::adcLSB_
private

Definition at line 36 of file BTLUncalibRecHitAlgo.cc.

Referenced by makeRecHit().

◆ adcNBits_

const uint32_t BTLUncalibRecHitAlgo::adcNBits_
private

Definition at line 34 of file BTLUncalibRecHitAlgo.cc.

◆ adcSaturation_

const double BTLUncalibRecHitAlgo::adcSaturation_
private

Definition at line 35 of file BTLUncalibRecHitAlgo.cc.

◆ c_LYSO_

const double BTLUncalibRecHitAlgo::c_LYSO_
private

Definition at line 42 of file BTLUncalibRecHitAlgo.cc.

Referenced by makeRecHit().

◆ timeCorr_p0_

const double BTLUncalibRecHitAlgo::timeCorr_p0_
private

Definition at line 39 of file BTLUncalibRecHitAlgo.cc.

Referenced by makeRecHit().

◆ timeCorr_p1_

const double BTLUncalibRecHitAlgo::timeCorr_p1_
private

Definition at line 40 of file BTLUncalibRecHitAlgo.cc.

Referenced by makeRecHit().

◆ timeCorr_p2_

const double BTLUncalibRecHitAlgo::timeCorr_p2_
private

Definition at line 41 of file BTLUncalibRecHitAlgo.cc.

Referenced by makeRecHit().

◆ timeError_

const reco::FormulaEvaluator BTLUncalibRecHitAlgo::timeError_
private

Definition at line 38 of file BTLUncalibRecHitAlgo.cc.

Referenced by makeRecHit().

◆ toaLSBToNS_

const double BTLUncalibRecHitAlgo::toaLSBToNS_
private

Definition at line 37 of file BTLUncalibRecHitAlgo.cc.

Referenced by makeRecHit().