File indexing completed on 2024-09-07 04:37:18
0001
0002 #include "FWCore/Framework/interface/Frameworkfwd.h"
0003 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0004 #include "FWCore/Framework/interface/ESHandle.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/Framework/interface/EventSetup.h"
0007 #include "FWCore/Framework/interface/MakerMacros.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0009
0010 #include "DataFormats/Candidate/interface/Candidate.h"
0011 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0012 #include "DataFormats/JetReco/interface/Jet.h"
0013 #include "DataFormats/JetReco/interface/JetCollection.h"
0014 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0015 #include "DataFormats/JetReco/interface/JetFloatAssociation.h"
0016
0017 #include "DataFormats/Common/interface/View.h"
0018 #include "DataFormats/Common/interface/Ref.h"
0019 #include "DataFormats/Common/interface/getRef.h"
0020
0021 #include "DataFormats/JetMatching/interface/JetFlavour.h"
0022 #include "DataFormats/JetMatching/interface/JetFlavourMatching.h"
0023 #include "DataFormats/JetMatching/interface/MatchedPartons.h"
0024 #include "DataFormats/JetMatching/interface/JetMatchedPartons.h"
0025
0026
0027 #include <memory>
0028 #include <string>
0029 #include <iostream>
0030 #include <vector>
0031 #include <Math/VectorUtil.h>
0032 #include <TMath.h>
0033
0034 class printJetFlavour : public edm::one::EDAnalyzer<> {
0035 public:
0036 explicit printJetFlavour(const edm::ParameterSet&);
0037 ~printJetFlavour() {}
0038 void analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup);
0039
0040 private:
0041 edm::InputTag sourcePartons_;
0042 edm::EDGetTokenT<reco::JetMatchedPartonsCollection> sourceByReferToken_;
0043 edm::EDGetTokenT<reco::JetFlavourMatchingCollection> sourceByValueToken_;
0044 edm::Handle<reco::JetMatchedPartonsCollection> theTagByRef;
0045 edm::Handle<reco::JetFlavourMatchingCollection> theTagByValue;
0046 };
0047
0048 using namespace std;
0049 using namespace reco;
0050 using namespace edm;
0051 using namespace ROOT::Math::VectorUtil;
0052
0053 printJetFlavour::printJetFlavour(const edm::ParameterSet& iConfig) {
0054 sourceByReferToken_ = consumes<reco::JetMatchedPartonsCollection>(iConfig.getParameter<InputTag>("srcByReference"));
0055 sourceByValueToken_ = consumes<reco::JetFlavourMatchingCollection>(iConfig.getParameter<InputTag>("srcByValue"));
0056 }
0057
0058 void printJetFlavour::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0059 cout << "[printJetFlavour] analysing event " << iEvent.id() << endl;
0060 try {
0061 iEvent.getByToken(sourceByReferToken_, theTagByRef);
0062 iEvent.getByToken(sourceByValueToken_, theTagByValue);
0063 } catch (std::exception& ce) {
0064 cerr << "[printJetFlavour] caught std::exception " << ce.what() << endl;
0065 return;
0066 }
0067
0068 cout << "-------------------- Jet Flavour by Ref From Partons--------------" << endl;
0069 for (JetMatchedPartonsCollection::const_iterator j = theTagByRef->begin(); j != theTagByRef->end(); j++) {
0070 const Jet* aJet = (*j).first.get();
0071 const MatchedPartons aMatch = (*j).second;
0072 printf("[printJetFlavour] (pt,eta,phi) jet = %7.2f %6.3f %6.3f \n", aJet->et(), aJet->eta(), aJet->phi());
0073 const GenParticleRef theHeaviest = aMatch.heaviest();
0074 if (theHeaviest.isNonnull()) {
0075 float dist = DeltaR(aJet->p4(), theHeaviest.get()->p4());
0076 cout << setprecision(2) << setw(6) << fixed
0077 << " theHeaviest flav (pt,eta,phi)=" << theHeaviest.get()->pdgId() << " ("
0078 << theHeaviest.get()->et() << "," << theHeaviest.get()->eta() << "," << theHeaviest.get()->phi()
0079 << ") Dr=" << dist << endl;
0080 }
0081 const GenParticleRef theNearest2 = aMatch.nearest_status2();
0082 if (theNearest2.isNonnull()) {
0083 float dist = DeltaR(aJet->p4(), theNearest2.get()->p4());
0084 cout << " theNearest Stat2 flav (pt,eta,phi)=" << theNearest2.get()->pdgId() << " ("
0085 << theNearest2.get()->et() << "," << theNearest2.get()->eta() << "," << theNearest2.get()->phi()
0086 << ") Dr=" << dist << endl;
0087 }
0088 const GenParticleRef theNearest3 = aMatch.nearest_status3();
0089 if (theNearest3.isNonnull()) {
0090 float dist = DeltaR(aJet->p4(), theNearest3.get()->p4());
0091 cout << " theNearest Stat3 flav (pt,eta,phi)=" << theNearest3.get()->pdgId() << " ("
0092 << theNearest3.get()->et() << "," << theNearest3.get()->eta() << "," << theNearest3.get()->phi()
0093 << ") Dr=" << dist << endl;
0094 }
0095 const GenParticleRef thePhyDef = aMatch.physicsDefinitionParton();
0096 if (thePhyDef.isNonnull()) {
0097 float dist = DeltaR(aJet->p4(), thePhyDef.get()->p4());
0098 cout << " thePhysDefinition flav (pt,eta,phi)=" << thePhyDef.get()->pdgId() << " ("
0099 << thePhyDef.get()->et() << "," << thePhyDef.get()->eta() << "," << thePhyDef.get()->phi() << ") Dr=" << dist
0100 << endl;
0101 }
0102 const GenParticleRef theAlgDef = aMatch.algoDefinitionParton();
0103 if (theAlgDef.isNonnull()) {
0104 float dist = DeltaR(aJet->p4(), theAlgDef.get()->p4());
0105 cout << " theAlgoDefinition flav (pt,eta,phi)=" << theAlgDef.get()->pdgId() << " ("
0106 << theAlgDef.get()->et() << "," << theAlgDef.get()->eta() << "," << theAlgDef.get()->phi() << ") Dr=" << dist
0107 << endl;
0108 }
0109 }
0110
0111 cout << "-------------------- Jet Flavour by Value ------------------------" << endl;
0112 for (JetFlavourMatchingCollection::const_iterator j = theTagByValue->begin(); j != theTagByValue->end(); j++) {
0113 RefToBase<Jet> aJet = (*j).first;
0114 const JetFlavour aFlav = (*j).second;
0115
0116 printf("[printJetFlavour] (pt,eta,phi) jet = %7.2f %6.3f %6.3f | parton = %7.2f %6.3f %6.3f | %4d\n",
0117 aJet.get()->et(),
0118 aJet.get()->eta(),
0119 aJet.get()->phi(),
0120 aFlav.getLorentzVector().pt(),
0121 aFlav.getLorentzVector().eta(),
0122 aFlav.getLorentzVector().phi(),
0123 aFlav.getFlavour());
0124 }
0125 }
0126
0127 DEFINE_FWK_MODULE(printJetFlavour);