Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:42

0001 #include "HeterogeneousCore/SonicTriton/interface/TritonEDProducer.h"
0002 #include "HeterogeneousCore/SonicTriton/interface/TritonData.h"
0003 
0004 #include "FWCore/Framework/interface/stream/EDProducer.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/MakerMacros.h"
0007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0008 
0009 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0010 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
0011 #include "DataFormats/Common/interface/ValueMap.h"
0012 
0013 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0014 
0015 #include "RecoEcal/EgammaClusterAlgos/interface/SCEnergyCorrectorDRN.h"
0016 
0017 #include <sstream>
0018 #include <string>
0019 #include <vector>
0020 #include <random>
0021 
0022 /*
0023  * SCEnergyCorrectorDRNProducer
0024  *
0025  * Simple producer to generate a set of corrected superclusters with the DRN regression
0026  * Based on RecoEcal/EgammaClusterProducers/SCEnergyCorrectorProducer by S. Harper (RAL/CERN)
0027  *
0028  * Author: Simon Rothman (UMN, MIT)
0029  *
0030  */
0031 
0032 namespace {
0033   float sigmoid(float x) { return 1.0f / (1.0f + exp(-x)); }
0034 
0035   float logcorrection(float x) {
0036     static float ln2 = log(2);
0037     return ln2 * 2 * (sigmoid(x) - 0.5);
0038   }
0039 
0040   float correction(float x) { return exp(-logcorrection(x)); }
0041 }  // namespace
0042 
0043 class SCEnergyCorrectorDRNProducer : public TritonEDProducer<> {
0044 public:
0045   explicit SCEnergyCorrectorDRNProducer(const edm::ParameterSet& iConfig);
0046 
0047   void beginLuminosityBlock(const edm::LuminosityBlock& iLumi, const edm::EventSetup& iSetup) override;
0048 
0049   void acquire(edm::Event const& iEvent, edm::EventSetup const& iSetup, Input& input) override;
0050   void produce(edm::Event& iEvent, const edm::EventSetup& iSetup, Output const& iOutput) override;
0051 
0052   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0053 
0054 private:
0055   SCEnergyCorrectorDRN energyCorrector_;
0056   edm::EDGetTokenT<reco::SuperClusterCollection> inputSCToken_;
0057 };
0058 
0059 SCEnergyCorrectorDRNProducer::SCEnergyCorrectorDRNProducer(const edm::ParameterSet& iConfig)
0060     : TritonEDProducer<>(iConfig),
0061       energyCorrector_(iConfig.getParameterSet("correctorCfg"), consumesCollector()),
0062       inputSCToken_(consumes<reco::SuperClusterCollection>(iConfig.getParameter<edm::InputTag>("inputSCs"))) {
0063   produces<reco::SuperClusterCollection>();
0064 }
0065 
0066 void SCEnergyCorrectorDRNProducer::beginLuminosityBlock(const edm::LuminosityBlock& iLumi,
0067                                                         const edm::EventSetup& iSetup) {
0068   energyCorrector_.setEventSetup(iSetup);
0069 }
0070 
0071 void SCEnergyCorrectorDRNProducer::acquire(edm::Event const& iEvent, edm::EventSetup const& iSetup, Input& iInput) {
0072   const auto& inputSCs = iEvent.get(inputSCToken_);
0073 
0074   if (inputSCs.empty()) {
0075     client_->setBatchSize(0);
0076     return;
0077   } else {
0078     client_->setBatchSize(1);
0079   }
0080 
0081   energyCorrector_.setEvent(iEvent);
0082   energyCorrector_.makeInput(iEvent, iInput, inputSCs);
0083 }
0084 
0085 void SCEnergyCorrectorDRNProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup, Output const& iOutput) {
0086   const auto& inputSCs = iEvent.get(inputSCToken_);
0087   if (inputSCs.empty())
0088     return;
0089 
0090   const auto& serverout = energyCorrector_.getOutput(iOutput);
0091 
0092   auto corrSCs = std::make_unique<reco::SuperClusterCollection>();
0093   unsigned i = 0;
0094   for (const auto& inputSC : inputSCs) {
0095     float corrEn = correction(serverout[0][0 + 6 * i]) * inputSC.rawEnergy();
0096     corrSCs->push_back(inputSC);
0097     corrSCs->back().setEnergy(corrEn);
0098     corrSCs->back().setCorrectedEnergy(corrEn);
0099     ++i;
0100   }
0101 
0102   auto scHandle = iEvent.put(std::move(corrSCs));
0103 }
0104 
0105 void SCEnergyCorrectorDRNProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0106   edm::ParameterSetDescription desc;
0107   desc.add<edm::ParameterSetDescription>("correctorCfg", SCEnergyCorrectorDRN::makePSetDescription());
0108   TritonClient::fillPSetDescription(desc);
0109   desc.add<edm::InputTag>("inputSCs", edm::InputTag("particleFlowSuperClusterECAL"));
0110   descriptions.add("scEnergyCorrectorDRNProducer", desc);
0111 }
0112 
0113 DEFINE_FWK_MODULE(SCEnergyCorrectorDRNProducer);