Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:56

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