Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 #include "CommonTools/CandAlgos/interface/ModifyObjectValueBase.h"
0002 #include "DataFormats/Common/interface/ValueMap.h"
0003 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0004 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0005 #include "FWCore/Framework/interface/stream/EDProducer.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011 
0012 class GEDGsfElectronFinalizer : public edm::stream::EDProducer<> {
0013 public:
0014   explicit GEDGsfElectronFinalizer(const edm::ParameterSet&);
0015 
0016   void produce(edm::Event&, const edm::EventSetup&) override;
0017 
0018 private:
0019   const edm::EDGetTokenT<reco::GsfElectronCollection> previousGsfElectrons_;
0020   const edm::EDGetTokenT<reco::PFCandidateCollection> pfCandidates_;
0021   std::vector<edm::EDGetTokenT<edm::ValueMap<float> > > tokenElectronIsoVals_;
0022   std::unique_ptr<ModifyObjectValueBase> gedRegression_;
0023 
0024   const edm::EDPutTokenT<reco::GsfElectronCollection> putToken_;
0025 };
0026 
0027 using edm::InputTag;
0028 using edm::ValueMap;
0029 
0030 using reco::GsfElectronCollection;
0031 
0032 GEDGsfElectronFinalizer::GEDGsfElectronFinalizer(const edm::ParameterSet& cfg)
0033     : previousGsfElectrons_(consumes<GsfElectronCollection>(cfg.getParameter<InputTag>("previousGsfElectronsTag"))),
0034       pfCandidates_(consumes<reco::PFCandidateCollection>(cfg.getParameter<InputTag>("pfCandidatesTag"))),
0035       putToken_{produces<reco::GsfElectronCollection>()} {
0036   edm::ParameterSet pfIsoVals(cfg.getParameter<edm::ParameterSet>("pfIsolationValues"));
0037 
0038   tokenElectronIsoVals_ = {consumes<ValueMap<float> >(pfIsoVals.getParameter<InputTag>("pfSumChargedHadronPt")),
0039                            consumes<ValueMap<float> >(pfIsoVals.getParameter<InputTag>("pfSumPhotonEt")),
0040                            consumes<ValueMap<float> >(pfIsoVals.getParameter<InputTag>("pfSumNeutralHadronEt")),
0041                            consumes<ValueMap<float> >(pfIsoVals.getParameter<InputTag>("pfSumPUPt")),
0042                            consumes<ValueMap<float> >(pfIsoVals.getParameter<InputTag>("pfSumEcalClusterEt")),
0043                            consumes<ValueMap<float> >(pfIsoVals.getParameter<InputTag>("pfSumHcalClusterEt"))};
0044 
0045   if (cfg.existsAs<edm::ParameterSet>("regressionConfig")) {
0046     auto const& iconf = cfg.getParameterSet("regressionConfig");
0047     auto const& mname = iconf.getParameter<std::string>("modifierName");
0048     auto cc = consumesCollector();
0049     gedRegression_ = ModifyObjectValueFactory::get()->create(mname, iconf, cc);
0050   }
0051 }
0052 
0053 void GEDGsfElectronFinalizer::produce(edm::Event& event, const edm::EventSetup& setup) {
0054   // Output collection
0055   reco::GsfElectronCollection outputElectrons;
0056 
0057   if (gedRegression_) {
0058     gedRegression_->setEvent(event);
0059     gedRegression_->setEventContent(setup);
0060   }
0061 
0062   // read input collections
0063   // electrons
0064   auto gedElectronHandle = event.getHandle(previousGsfElectrons_);
0065 
0066   // PFCandidates
0067   auto pfCandidateHandle = event.getHandle(pfCandidates_);
0068   // value maps
0069   std::vector<edm::ValueMap<float> const*> isolationValueMaps(tokenElectronIsoVals_.size());
0070 
0071   for (unsigned i = 0; i < tokenElectronIsoVals_.size(); ++i) {
0072     isolationValueMaps[i] = &event.get(tokenElectronIsoVals_[i]);
0073   }
0074 
0075   // prepare a map of PFCandidates having a valid GsfTrackRef to save time
0076   std::map<reco::GsfTrackRef, const reco::PFCandidate*> gsfPFMap;
0077   for (auto const& pfCand : *pfCandidateHandle) {
0078     // First check that the GsfTrack is non null
0079     if (pfCand.gsfTrackRef().isNonnull()) {
0080       if (abs(pfCand.pdgId()) == 11)  // consider only the electrons
0081         gsfPFMap[pfCand.gsfTrackRef()] = &pfCand;
0082     }
0083   }
0084 
0085   // Now loop on the electrons
0086   unsigned nele = gedElectronHandle->size();
0087   for (unsigned iele = 0; iele < nele; ++iele) {
0088     reco::GsfElectronRef myElectronRef(gedElectronHandle, iele);
0089 
0090     reco::GsfElectron newElectron(*myElectronRef);
0091     reco::GsfElectron::PflowIsolationVariables isoVariables;
0092     isoVariables.sumChargedHadronPt = (*(isolationValueMaps)[0])[myElectronRef];
0093     isoVariables.sumPhotonEt = (*(isolationValueMaps)[1])[myElectronRef];
0094     isoVariables.sumNeutralHadronEt = (*(isolationValueMaps)[2])[myElectronRef];
0095     isoVariables.sumPUPt = (*(isolationValueMaps)[3])[myElectronRef];
0096     isoVariables.sumEcalClusterEt = (*(isolationValueMaps)[4])[myElectronRef];
0097     isoVariables.sumHcalClusterEt = (*(isolationValueMaps)[5])[myElectronRef];
0098 
0099     newElectron.setPfIsolationVariables(isoVariables);
0100 
0101     // now set a status if not already done (in GEDGsfElectronProducer.cc)
0102     //     std::cout << " previous status " << newElectron.mvaOutput().status << std::endl;
0103     if (newElectron.mvaOutput().status <= 0) {
0104       reco::GsfElectron::MvaOutput myMvaOutput(newElectron.mvaOutput());
0105       if (gsfPFMap.find(newElectron.gsfTrack()) != gsfPFMap.end()) {
0106         // it means that there is a PFCandidate with the same GsfTrack
0107         myMvaOutput.status = 3;  //as defined in PFCandidateEGammaExtra.h
0108         //this is currently fully redundant with mvaOutput.stats so candidate for removal
0109         newElectron.setPassPflowPreselection(true);
0110       } else {
0111         myMvaOutput.status = 4;  //
0112         //this is currently fully redundant with mvaOutput.stats so candidate for removal
0113         newElectron.setPassPflowPreselection(false);
0114       }
0115       newElectron.setMvaOutput(myMvaOutput);
0116     }
0117 
0118     if (gedRegression_) {
0119       gedRegression_->modifyObject(newElectron);
0120     }
0121     outputElectrons.push_back(newElectron);
0122   }
0123 
0124   event.emplace(putToken_, std::move(outputElectrons));
0125 }
0126 
0127 #include "FWCore/Framework/interface/MakerMacros.h"
0128 DEFINE_FWK_MODULE(GEDGsfElectronFinalizer);