Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-04-02 03:10:56

0001 #include "RecoLocalFastTime/FTLCommonAlgos/interface/MTDUncalibratedRecHitAlgoBase.h"
0002 #include "RecoLocalFastTime/FTLClusterizer/interface/BTLRecHitsErrorEstimatorIM.h"
0003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0004 
0005 #include "CommonTools/Utils/interface/FormulaEvaluator.h"
0006 
0007 class BTLUncalibRecHitAlgo : public BTLUncalibratedRecHitAlgoBase {
0008 public:
0009   /// Constructor
0010   BTLUncalibRecHitAlgo(const edm::ParameterSet& conf, edm::ConsumesCollector& sumes)
0011       : MTDUncalibratedRecHitAlgoBase<BTLDataFrame>(conf, sumes),
0012         adcNBits_(conf.getParameter<uint32_t>("adcNbits")),
0013         adcSaturation_(conf.getParameter<double>("adcSaturation")),
0014         adcLSB_(adcSaturation_ / (1 << adcNBits_)),
0015         toaLSBToNS_(conf.getParameter<double>("toaLSB_ns")),
0016         timeError_(conf.getParameter<std::string>("timeResolutionInNs")),
0017         timeCorr_p0_(conf.getParameter<double>("timeCorr_p0")),
0018         timeCorr_p1_(conf.getParameter<double>("timeCorr_p1")),
0019         timeCorr_p2_(conf.getParameter<double>("timeCorr_p2")),
0020         c_LYSO_(conf.getParameter<double>("c_LYSO")) {}
0021 
0022   /// Destructor
0023   ~BTLUncalibRecHitAlgo() override {}
0024 
0025   /// get event and eventsetup information
0026   void getEvent(const edm::Event&) final {}
0027   void getEventSetup(const edm::EventSetup&) final {}
0028 
0029   /// make the rec hit
0030   FTLUncalibratedRecHit makeRecHit(const BTLDataFrame& dataFrame) const final;
0031 
0032 private:
0033   const uint32_t adcNBits_;
0034   const double adcSaturation_;
0035   const double adcLSB_;
0036   const double toaLSBToNS_;
0037   const reco::FormulaEvaluator timeError_;
0038   const double timeCorr_p0_;
0039   const double timeCorr_p1_;
0040   const double timeCorr_p2_;
0041   const double c_LYSO_;
0042 };
0043 
0044 FTLUncalibratedRecHit BTLUncalibRecHitAlgo::makeRecHit(const BTLDataFrame& dataFrame) const {
0045   // The reconstructed amplitudes and times are saved in a std::pair
0046   //    BTL tile geometry (1 SiPM): only the first value of the amplitude
0047   //                                and time pairs is used.
0048   //    BTL bar geometry (2 SiPMs): both values of the amplitude and
0049   //                                time pairs are filled.
0050 
0051   std::pair<float, float> amplitude(0., 0.);
0052   std::pair<float, float> time(0., 0.);
0053 
0054   unsigned char flag = 0;
0055 
0056   const auto& sampleLeft = dataFrame.sample(0);
0057   const auto& sampleRight = dataFrame.sample(1);
0058 
0059   double nHits = 0.;
0060 
0061   if (sampleLeft.data() > 0) {
0062     amplitude.first = float(sampleLeft.data()) * adcLSB_;
0063     time.first = float(sampleLeft.toa()) * toaLSBToNS_;
0064 
0065     nHits += 1.;
0066 
0067     // Correct the time of the left SiPM for the time-walk
0068     time.first -= timeCorr_p0_ * pow(amplitude.first, timeCorr_p1_) + timeCorr_p2_;
0069     flag |= 0x1;
0070   }
0071 
0072   // --- If available, reconstruct the amplitude and time of the second SiPM
0073   if (sampleRight.data() > 0) {
0074     amplitude.second = sampleRight.data() * adcLSB_;
0075     time.second = sampleRight.toa() * toaLSBToNS_;
0076 
0077     nHits += 1.;
0078 
0079     // Correct the time of the right SiPM for the time-walk
0080     time.second -= timeCorr_p0_ * pow(amplitude.second, timeCorr_p1_) + timeCorr_p2_;
0081     flag |= (0x1 << 1);
0082   }
0083 
0084   // --- Calculate the error on the hit time using the provided parameterization
0085 
0086   const std::array<double, 1> amplitudeV = {{(amplitude.first + amplitude.second) / nHits}};
0087   const std::array<double, 1> emptyV = {{0.}};
0088 
0089   double timeError = (nHits > 0. ? timeError_.evaluate(amplitudeV, emptyV) : -1.);
0090 
0091   // Calculate the position
0092   // Distance from center of bar to hit
0093   float position = 0.5f * (c_LYSO_ * (time.second - time.first));
0094   float positionError = BTLRecHitsErrorEstimatorIM::positionError();
0095 
0096   LogDebug("BTLUncalibRecHit") << "ADC+: set the charge to: (" << amplitude.first << ", " << amplitude.second << ")  ("
0097                                << sampleLeft.data() << ", " << sampleRight.data() << "  " << adcLSB_ << ' '
0098                                << std::endl;
0099   LogDebug("BTLUncalibRecHit") << "TDC+: set the time to: (" << time.first << ", " << time.second << ")  ("
0100                                << sampleLeft.toa() << ", " << sampleRight.toa() << "  " << toaLSBToNS_ << ' '
0101                                << std::endl;
0102 
0103   return FTLUncalibratedRecHit(
0104       dataFrame.id(), dataFrame.row(), dataFrame.column(), amplitude, time, timeError, position, positionError, flag);
0105 }
0106 
0107 #include "FWCore/Framework/interface/MakerMacros.h"
0108 DEFINE_EDM_PLUGIN(BTLUncalibratedRecHitAlgoFactory, BTLUncalibRecHitAlgo, "BTLUncalibRecHitAlgo");