Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "FWCore/Framework/interface/stream/EDProducer.h"
0002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0003 #include "FWCore/Framework/interface/Event.h"
0004 #include "FWCore/Framework/interface/EventSetup.h"
0005 
0006 #include "DataFormats/EcalDigi/interface/EcalDigiCollections.h"
0007 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
0008 #include "DataFormats/EcalRecHit/interface/EcalUncalibratedRecHit.h"
0009 
0010 #include "DataFormats/EcalDigi/interface/EcalConstants.h"
0011 
0012 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0013 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0014 
0015 class EcalUncalibRecHitPhase2WeightsProducer : public edm::stream::EDProducer<> {
0016 public:
0017   explicit EcalUncalibRecHitPhase2WeightsProducer(const edm::ParameterSet& ps);
0018   void produce(edm::Event& evt, const edm::EventSetup& es) override;
0019   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0020 
0021 private:
0022   const float tRise_;
0023   const float tFall_;
0024   const std::vector<double> ampWeights_;
0025   const std::vector<double> timeWeights_;
0026 
0027   const edm::EDGetTokenT<EBDigiCollectionPh2> ebDigiCollectionToken_;
0028   const edm::EDPutTokenT<EBUncalibratedRecHitCollection> ebUncalibRecHitCollectionToken_;
0029 };
0030 
0031 void EcalUncalibRecHitPhase2WeightsProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0032   edm::ParameterSetDescription desc;
0033 
0034   desc.add<std::string>("EBhitCollection", "EcalUncalibRecHitsEB");
0035   desc.add<double>("tRise", 0.2);
0036   desc.add<double>("tFall", 2.);
0037   desc.add<std::vector<double>>("ampWeights",
0038                                 {-0.121016,
0039                                  -0.119899,
0040                                  -0.120923,
0041                                  -0.0848959,
0042                                  0.261041,
0043                                  0.509881,
0044                                  0.373591,
0045                                  0.134899,
0046                                  -0.0233605,
0047                                  -0.0913195,
0048                                  -0.112452,
0049                                  -0.118596,
0050                                  -0.121737,
0051                                  -0.121737,
0052                                  -0.121737,
0053                                  -0.121737});
0054   desc.add<std::vector<double>>("timeWeights",
0055                                 {0.429452,
0056                                  0.442762,
0057                                  0.413327,
0058                                  0.858327,
0059                                  4.42324,
0060                                  2.04369,
0061                                  -3.42426,
0062                                  -4.16258,
0063                                  -2.36061,
0064                                  -0.725371,
0065                                  0.0727267,
0066                                  0.326005,
0067                                  0.402035,
0068                                  0.404287,
0069                                  0.434207,
0070                                  0.422775});
0071 
0072   desc.add<edm::InputTag>("BarrelDigis", edm::InputTag("simEcalUnsuppressedDigis", ""));
0073 
0074   descriptions.addWithDefaultLabel(desc);
0075 }
0076 
0077 EcalUncalibRecHitPhase2WeightsProducer::EcalUncalibRecHitPhase2WeightsProducer(const edm::ParameterSet& ps)
0078     : tRise_(ps.getParameter<double>("tRise")),
0079       tFall_(ps.getParameter<double>("tFall")),
0080       ampWeights_(ps.getParameter<std::vector<double>>("ampWeights")),
0081       timeWeights_(ps.getParameter<std::vector<double>>("timeWeights")),
0082       ebDigiCollectionToken_(consumes<EBDigiCollectionPh2>(ps.getParameter<edm::InputTag>("BarrelDigis"))),
0083       ebUncalibRecHitCollectionToken_(
0084           produces<EBUncalibratedRecHitCollection>(ps.getParameter<std::string>("EBhitCollection"))) {}
0085 
0086 void EcalUncalibRecHitPhase2WeightsProducer::produce(edm::Event& evt, const edm::EventSetup& es) {
0087   // retrieve digis
0088   const EBDigiCollectionPh2* pdigis = &evt.get(ebDigiCollectionToken_);
0089 
0090   // prepare output
0091   auto ebUncalibRechits = std::make_unique<EBUncalibratedRecHitCollection>();
0092 
0093   for (auto itdg = pdigis->begin(); itdg != pdigis->end(); ++itdg) {
0094     EBDataFrame digi(*itdg);
0095     EcalDataFrame_Ph2 dataFrame(*itdg);
0096     DetId detId(digi.id());
0097 
0098     bool g1 = false;
0099 
0100     std::vector<float> timetrace;
0101     std::vector<float> adctrace;
0102 
0103     int nSamples = digi.size();
0104 
0105     timetrace.reserve(nSamples);
0106     adctrace.reserve(nSamples);
0107 
0108     float amp = 0;
0109     float t0 = 0;
0110     float gratio;
0111 
0112     for (int sample = 0; sample < nSamples; ++sample) {
0113       EcalLiteDTUSample thisSample = dataFrame[sample];
0114       gratio = ecalPh2::gains[thisSample.gainId()];
0115       adctrace.push_back(thisSample.adc() * gratio);
0116 
0117       amp = amp + adctrace[sample] * ampWeights_[sample];
0118 
0119       t0 = t0 + adctrace[sample] * timeWeights_[sample];
0120 
0121       if (thisSample.gainId() == 1)
0122         g1 = true;
0123 
0124       timetrace.push_back(sample);
0125 
0126     }  // loop on samples
0127 
0128     float amp_e = 1;
0129     float t0_e = 0;
0130 
0131     EcalUncalibratedRecHit rhit(detId, amp, 0., t0, 0., 0);  // rhit(detIt, amp, pedestal, t0, chi2, flags)
0132     rhit.setAmplitudeError(amp_e);
0133     rhit.setJitterError(t0_e);
0134     if (g1)
0135       rhit.setFlagBit(EcalUncalibratedRecHit::kHasSwitchToGain1);
0136 
0137     ebUncalibRechits->push_back(rhit);
0138 
0139   }  // loop on digis
0140 
0141   evt.put(ebUncalibRecHitCollectionToken_, std::move(ebUncalibRechits));
0142 }
0143 
0144 #include "FWCore/Framework/interface/MakerMacros.h"
0145 DEFINE_FWK_MODULE(EcalUncalibRecHitPhase2WeightsProducer);