Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-09-20 02:50:17

0001 // system includes
0002 #include <memory>
0003 
0004 // user includes
0005 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0006 #include "DataFormats/Common/interface/Handle.h"
0007 #include "DataFormats/TrackReco/interface/HitPattern.h"
0008 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0009 #include "FWCore/Framework/interface/Event.h"
0010 #include "FWCore/Framework/interface/Frameworkfwd.h"
0011 #include "FWCore/Framework/interface/global/EDProducer.h"
0012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0013 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 #include "FWCore/Utilities/interface/StreamID.h"
0018 
0019 // Electron selector
0020 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0021 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
0022 #include "DataFormats/EgammaCandidates/interface/Electron.h"
0023 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0024 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0025 
0026 // Muon Selector
0027 #include "DataFormats/MuonReco/interface/Muon.h"
0028 #include "DataFormats/MuonReco/interface/MuonFwd.h"
0029 #include "DataFormats/MuonReco/interface/MuonSelectors.h"
0030 
0031 // Jet Selector
0032 #include "DataFormats/JetReco/interface/PFJet.h"
0033 #include "DataFormats/JetReco/interface/PFJetCollection.h"
0034 #include "DataFormats/JetReco/interface/PileupJetIdentifier.h"
0035 #include "DataFormats/PatCandidates/interface/Jet.h"
0036 
0037 // b-jet Selector
0038 #include "DataFormats/BTauReco/interface/JetTag.h"
0039 
0040 // Met Selector
0041 #include "DataFormats/METReco/interface/PFMET.h"
0042 #include "DataFormats/METReco/interface/PFMETCollection.h"
0043 
0044 // ROOT includes
0045 #include "TLorentzVector.h"
0046 
0047 class TtbarTrackProducer : public edm::global::EDProducer<> {
0048 public:
0049   explicit TtbarTrackProducer(const edm::ParameterSet&);
0050   ~TtbarTrackProducer() override = default;
0051   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0052 
0053   void produce(edm::StreamID streamID, edm::Event& iEvent, edm::EventSetup const& iSetup) const override;
0054 
0055 private:
0056   // ----------member data ---------------------------
0057 
0058   const edm::InputTag electronTag_;
0059   const edm::InputTag jetsTag_;
0060   const edm::InputTag bjetsTag_;
0061   const edm::InputTag pfmetTag_;
0062   const edm::InputTag muonTag_;
0063   const edm::InputTag bsTag_;
0064   const edm::EDGetTokenT<reco::GsfElectronCollection> electronToken_;
0065   const edm::EDGetTokenT<reco::PFJetCollection> jetsToken_;
0066   const edm::EDGetTokenT<reco::JetTagCollection> bjetsToken_;
0067   const edm::EDGetTokenT<reco::PFMETCollection> pfmetToken_;
0068   const edm::EDGetTokenT<reco::MuonCollection> muonToken_;
0069   const edm::EDGetTokenT<reco::BeamSpot> bsToken_;
0070 
0071   const double maxEtaEle_;
0072   const double maxEtaMu_;
0073   const double minPt_;
0074   const double maxDeltaPhiInEB_;
0075   const double maxDeltaEtaInEB_;
0076   const double maxHOEEB_;
0077   const double maxSigmaiEiEEB_;
0078   const double maxDeltaPhiInEE_;
0079   const double maxDeltaEtaInEE_;
0080   const double maxHOEEE_;
0081   const double maxSigmaiEiEEE_;
0082 
0083   const double minChambers_;
0084 
0085   const double minEta_Jets_;
0086 
0087   const double btagFactor_;
0088 
0089   const double maxNormChi2_;
0090   const double maxD0_;
0091   const double maxDz_;
0092   const int minPixelHits_;
0093   const int minStripHits_;
0094   const double maxIsoEle_;
0095   const double maxIsoMu_;
0096   const double minPtHighestMu_;
0097   const double minPtHighestEle_;
0098   const double minPtHighest_Jets_;
0099   const double minPt_Jets_;
0100   const double minInvMass_;
0101   const double maxInvMass_;
0102   const double minMet_;
0103   const double maxMet_;
0104   const double minWmass_;
0105   const double maxWmass_;
0106 };
0107 
0108 void TtbarTrackProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0109   edm::ParameterSetDescription desc;
0110   desc.addUntracked<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"));
0111   desc.addUntracked<edm::InputTag>("jetsInputTag", edm::InputTag("ak4PFJetsCHS"));
0112   desc.addUntracked<edm::InputTag>("bjetsInputTag", edm::InputTag("pfDeepCSVJetTags", "probb"));
0113   desc.addUntracked<edm::InputTag>("pfmetTag", edm::InputTag("pfMet"));
0114   desc.addUntracked<edm::InputTag>("muonInputTag", edm::InputTag("muons"));
0115   desc.addUntracked<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"));
0116   desc.addUntracked<double>("maxEtaEle", 2.4);
0117   desc.addUntracked<double>("maxEtaMu", 2.1);
0118   desc.addUntracked<double>("minPt", 5);
0119 
0120   // for Electron only
0121   desc.addUntracked<double>("maxDeltaPhiInEB", .15);
0122   desc.addUntracked<double>("maxDeltaEtaInEB", .007);
0123   desc.addUntracked<double>("maxHOEEB", .12);
0124   desc.addUntracked<double>("maxSigmaiEiEEB", .01);
0125   desc.addUntracked<double>("maxDeltaPhiInEE", .1);
0126   desc.addUntracked<double>("maxDeltaEtaInEE", .009);
0127   desc.addUntracked<double>("maxHOEEB_", .10);
0128   desc.addUntracked<double>("maxSigmaiEiEEE", .03);
0129 
0130   // for Muon only
0131   desc.addUntracked<uint32_t>("minChambers", 2);
0132 
0133   // for Jets only
0134   desc.addUntracked<double>("minEta_Jets", 3.0);
0135 
0136   // for b-tag only
0137   desc.addUntracked<double>("btagFactor", 0.6);
0138 
0139   desc.addUntracked<double>("maxNormChi2", 10);
0140   desc.addUntracked<double>("maxD0", 0.02);
0141   desc.addUntracked<double>("maxDz", 20.);
0142   desc.addUntracked<uint32_t>("minPixelHits", 1);
0143   desc.addUntracked<uint32_t>("minStripHits", 8);
0144   desc.addUntracked<double>("maxIsoEle", 0.5);
0145   desc.addUntracked<double>("maxIsoMu", 0.3);
0146   desc.addUntracked<double>("minPtHighestMu", 24);
0147   desc.addUntracked<double>("minPtHighestEle", 32);
0148   desc.addUntracked<double>("minPtHighest_Jets", 30);
0149   desc.addUntracked<double>("minPt_Jets", 20);
0150   desc.addUntracked<double>("minInvMass", 140);
0151   desc.addUntracked<double>("maxInvMass", 200);
0152   desc.addUntracked<double>("minMet", 50);
0153   desc.addUntracked<double>("maxMet", 80);
0154   desc.addUntracked<double>("minWmass", 50);
0155   desc.addUntracked<double>("maxWmass", 130);
0156   descriptions.addWithDefaultLabel(desc);
0157 }
0158 
0159 using namespace std;
0160 using namespace edm;
0161 
0162 TtbarTrackProducer::TtbarTrackProducer(const edm::ParameterSet& ps)
0163     : electronTag_(ps.getUntrackedParameter<edm::InputTag>("electronInputTag", edm::InputTag("gedGsfElectrons"))),
0164       jetsTag_(ps.getUntrackedParameter<edm::InputTag>("jetsInputTag", edm::InputTag("ak4PFJetsCHS"))),
0165 
0166       bjetsTag_(ps.getUntrackedParameter<edm::InputTag>("bjetsInputTag", edm::InputTag("pfDeepCSVJetTags", "probb"))),
0167       pfmetTag_(ps.getUntrackedParameter<edm::InputTag>("pfmetTag", edm::InputTag("pfMet"))),  //("pfMetT1T2Txy"))),
0168       muonTag_(ps.getUntrackedParameter<edm::InputTag>("muonInputTag", edm::InputTag("muons"))),
0169       bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
0170       electronToken_(consumes<reco::GsfElectronCollection>(electronTag_)),
0171       jetsToken_(consumes<reco::PFJetCollection>(jetsTag_)),
0172       bjetsToken_(consumes<reco::JetTagCollection>(bjetsTag_)),
0173       pfmetToken_(consumes<reco::PFMETCollection>(pfmetTag_)),
0174       muonToken_(consumes<reco::MuonCollection>(muonTag_)),
0175       bsToken_(consumes<reco::BeamSpot>(bsTag_)),
0176       maxEtaEle_(ps.getUntrackedParameter<double>("maxEtaEle", 2.4)),
0177       maxEtaMu_(ps.getUntrackedParameter<double>("maxEtaMu", 2.1)),
0178       minPt_(ps.getUntrackedParameter<double>("minPt", 5)),
0179 
0180       // for Electron only
0181       maxDeltaPhiInEB_(ps.getUntrackedParameter<double>("maxDeltaPhiInEB", .15)),
0182       maxDeltaEtaInEB_(ps.getUntrackedParameter<double>("maxDeltaEtaInEB", .007)),
0183       maxHOEEB_(ps.getUntrackedParameter<double>("maxHOEEB", .12)),
0184       maxSigmaiEiEEB_(ps.getUntrackedParameter<double>("maxSigmaiEiEEB", .01)),
0185       maxDeltaPhiInEE_(ps.getUntrackedParameter<double>("maxDeltaPhiInEE", .1)),
0186       maxDeltaEtaInEE_(ps.getUntrackedParameter<double>("maxDeltaEtaInEE", .009)),
0187       maxHOEEE_(ps.getUntrackedParameter<double>("maxHOEEB_", .10)),
0188       maxSigmaiEiEEE_(ps.getUntrackedParameter<double>("maxSigmaiEiEEE", .03)),
0189 
0190       // for Muon only
0191       minChambers_(ps.getUntrackedParameter<uint32_t>("minChambers", 2)),
0192 
0193       // for Jets only
0194       minEta_Jets_(ps.getUntrackedParameter<double>("minEta_Jets", 3.0)),
0195 
0196       // for b-tag only
0197       btagFactor_(ps.getUntrackedParameter<double>("btagFactor", 0.6)),
0198 
0199       maxNormChi2_(ps.getUntrackedParameter<double>("maxNormChi2", 10)),
0200       maxD0_(ps.getUntrackedParameter<double>("maxD0", 0.02)),
0201       maxDz_(ps.getUntrackedParameter<double>("maxDz", 20.)),
0202       minPixelHits_(ps.getUntrackedParameter<uint32_t>("minPixelHits", 1)),
0203       minStripHits_(ps.getUntrackedParameter<uint32_t>("minStripHits", 8)),
0204       maxIsoEle_(ps.getUntrackedParameter<double>("maxIsoEle", 0.5)),
0205       maxIsoMu_(ps.getUntrackedParameter<double>("maxIsoMu", 0.3)),
0206       minPtHighestMu_(ps.getUntrackedParameter<double>("minPtHighestMu", 24)),
0207       minPtHighestEle_(ps.getUntrackedParameter<double>("minPtHighestEle", 32)),
0208       minPtHighest_Jets_(ps.getUntrackedParameter<double>("minPtHighest_Jets", 30)),
0209       minPt_Jets_(ps.getUntrackedParameter<double>("minPt_Jets", 20)),
0210       minInvMass_(ps.getUntrackedParameter<double>("minInvMass", 140)),
0211       maxInvMass_(ps.getUntrackedParameter<double>("maxInvMass", 200)),
0212       minMet_(ps.getUntrackedParameter<double>("minMet", 50)),
0213       maxMet_(ps.getUntrackedParameter<double>("maxMet", 80)),
0214       minWmass_(ps.getUntrackedParameter<double>("minWmass", 50)),
0215       maxWmass_(ps.getUntrackedParameter<double>("maxWmass", 130)) {
0216   produces<reco::TrackCollection>("");
0217 }
0218 
0219 void TtbarTrackProducer::produce(edm::StreamID streamID, edm::Event& iEvent, edm::EventSetup const& iSetup) const {
0220   std::unique_ptr<reco::TrackCollection> outputTColl(new reco::TrackCollection());
0221 
0222   // beamspot
0223   edm::Handle<reco::BeamSpot> beamSpot;
0224   iEvent.getByToken(bsToken_, beamSpot);
0225 
0226   // Read Electron Collection
0227   edm::Handle<reco::GsfElectronCollection> electronColl;
0228   iEvent.getByToken(electronToken_, electronColl);
0229   std::vector<TLorentzVector> list_ele;
0230   std::vector<int> chrgeList_ele;
0231 
0232   // for Track jet collections
0233   edm::Handle<reco::PFJetCollection> jetColl;
0234   iEvent.getByToken(jetsToken_, jetColl);
0235   reco::PFJetCollection Jets;
0236   std::vector<TLorentzVector> list_jets;
0237 
0238   if (jetColl.isValid()) {
0239     for (const auto& jets : *jetColl) {
0240       if (jets.pt() < minPt_Jets_)
0241         continue;
0242       if (std::fabs(jets.eta()) < minEta_Jets_)
0243         continue;
0244       TLorentzVector lv_jets;  // lv_bJets;
0245       lv_jets.SetPtEtaPhiE(jets.pt(), jets.eta(), jets.phi(), jets.energy());
0246       list_jets.push_back(lv_jets);
0247       Jets.push_back(jets);
0248     }
0249   }
0250 
0251   edm::Handle<reco::JetTagCollection> bTagHandle;
0252   iEvent.getByToken(bjetsToken_, bTagHandle);
0253   const reco::JetTagCollection& bTags = *(bTagHandle.product());
0254   std::vector<TLorentzVector> list_bjets;
0255 
0256   if (!bTags.empty()) {
0257     for (unsigned bj = 0; bj != bTags.size(); ++bj) {
0258       TLorentzVector lv_bjets;
0259       lv_bjets.SetPtEtaPhiE(
0260           bTags[bj].first->pt(), bTags[bj].first->eta(), bTags[bj].first->phi(), bTags[bj].first->energy());
0261       if (bTags[bj].second > btagFactor_)
0262         list_bjets.push_back(lv_bjets);
0263     }
0264   }
0265 
0266   for (unsigned int i = 0; i != Jets.size(); i++) {
0267     reco::TrackRefVector vector = Jets[i].getTrackRefs();
0268     for (unsigned int j = 0; j != vector.size(); j++) {
0269       outputTColl->push_back(*vector[j]);
0270     }
0271   }
0272 
0273   iEvent.put(std::move(outputTColl));
0274 }
0275 
0276 // Define this as a plug-in
0277 #include "FWCore/Framework/interface/MakerMacros.h"
0278 DEFINE_FWK_MODULE(TtbarTrackProducer);