CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/PhysicsTools/PatExamples/plugins/PatBJetVertexAnalyzer.cc

Go to the documentation of this file.
00001 #include <algorithm>
00002 #include <iostream>
00003 #include <cmath>
00004 #include <vector>
00005 #include <string>
00006 
00007 #include <TH1.h>
00008 #include <TProfile.h>
00009 
00010 #include <Math/VectorUtil.h>
00011 #include <Math/GenVector/PxPyPzE4D.h>
00012 #include <Math/GenVector/PxPyPzM4D.h>
00013 
00014 #include "FWCore/Framework/interface/Frameworkfwd.h"
00015 #include "FWCore/Framework/interface/Event.h"
00016 #include "FWCore/Framework/interface/EventSetup.h"
00017 #include "FWCore/Framework/interface/EDAnalyzer.h"
00018 #include "FWCore/Utilities/interface/InputTag.h"
00019 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00020 #include "FWCore/ServiceRegistry/interface/Service.h"
00021 
00022 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00023 
00024 #include "DataFormats/BTauReco/interface/SecondaryVertexTagInfo.h"
00025 #include "DataFormats/Math/interface/LorentzVector.h"
00026 #include "DataFormats/VertexReco/interface/Vertex.h"
00027 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00028 #include "DataFormats/TrackReco/interface/Track.h"
00029 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00030 #include "DataFormats/PatCandidates/interface/Jet.h"
00031 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
00032 
00033 class PatBJetVertexAnalyzer : public edm::EDAnalyzer  {
00034     public: 
00036         PatBJetVertexAnalyzer(const edm::ParameterSet &params);
00037         ~PatBJetVertexAnalyzer();
00038 
00039         // virtual methods called from base class EDAnalyzer
00040         virtual void beginJob();
00041         virtual void analyze(const edm::Event &event, const edm::EventSetup &es);
00042 
00043     private:
00044         // configuration parameters
00045         edm::InputTag jets_;
00046 
00047         double jetPtCut_;               // minimum (uncorrected) jet energy
00048         double jetEtaCut_;              // maximum |eta| for jet
00049         double maxDeltaR_;              // angle between jet and tracks
00050 
00051         // jet flavour constants
00052 
00053         enum Flavour {
00054                 ALL_JETS = 0,
00055                 UDSG_JETS,
00056                 C_JETS,
00057                 B_JETS,
00058                 NONID_JETS,
00059                 N_JET_TYPES
00060         };
00061 
00062         TH1 * flavours_;
00063 
00064         // one group of plots per jet flavour;
00065         struct Plots {
00066                 TH1 *nVertices;
00067                 TH1 *deltaR, *mass, *dist, *distErr, *distSig;
00068                 TH1 *nTracks, *chi2;
00069         } plots_[N_JET_TYPES];
00070 };
00071 
00072 PatBJetVertexAnalyzer::PatBJetVertexAnalyzer(const edm::ParameterSet &params) :
00073         jets_(params.getParameter<edm::InputTag>("jets")),
00074         jetPtCut_(params.getParameter<double>("jetPtCut")),
00075         jetEtaCut_(params.getParameter<double>("jetEtaCut"))
00076 {
00077 }
00078 
00079 PatBJetVertexAnalyzer::~PatBJetVertexAnalyzer()
00080 {
00081 }
00082 
00083 void PatBJetVertexAnalyzer::beginJob()
00084 {
00085         // retrieve handle to auxiliary service
00086         //  used for storing histograms into ROOT file
00087         edm::Service<TFileService> fs;
00088 
00089         flavours_ = fs->make<TH1F>("flavours", "jet flavours", 5, 0, 5);
00090 
00091         // book histograms for all jet flavours
00092         for(unsigned int i = 0; i < N_JET_TYPES; i++) {
00093                 Plots &plots = plots_[i];
00094                 const char *flavour, *name;
00095 
00096                 switch((Flavour)i) {
00097                     case ALL_JETS:
00098                         flavour = "all jets";
00099                         name = "all";
00100                         break;
00101                     case UDSG_JETS:
00102                         flavour = "light flavour jets";
00103                         name = "udsg";
00104                         break;
00105                     case C_JETS:
00106                         flavour = "charm jets";
00107                         name = "c";
00108                         break;
00109                     case B_JETS:
00110                         flavour = "bottom jets";
00111                         name = "b";
00112                         break;
00113                     default:
00114                         flavour = "unidentified jets";
00115                         name = "ni";
00116                         break;
00117                 }
00118 
00119                 plots.nVertices = fs->make<TH1F>(Form("nVertices_%s", name),
00120                                                  Form("number of secondary vertices in %s", flavour),
00121                                                  5, 0, 5);
00122                 plots.deltaR = fs->make<TH1F>(Form("deltaR_%s", name),
00123                                               Form("\\DeltaR between vertex direction and jet direction in %s", flavour),
00124                                               100, 0, 0.5);
00125                 plots.mass = fs->make<TH1F>(Form("mass_%s", name),
00126                                             Form("vertex mass in %s", flavour),
00127                                             100, 0, 10);
00128                 plots.dist = fs->make<TH1F>(Form("dist_%s", name),
00129                                             Form("Transverse distance between PV and SV in %s", flavour),
00130                                             100, 0, 2);
00131                 plots.distErr = fs->make<TH1F>(Form("distErr_%s", name),
00132                                                Form("Transverse distance error between PV and SV in %s", flavour),
00133                                                100, 0, 0.5);
00134                 plots.distSig = fs->make<TH1F>(Form("distSig_%s", name),
00135                                                Form("Transverse distance significance between PV and SV in %s", flavour),
00136                                                100, 0, 50);
00137                 plots.nTracks = fs->make<TH1F>(Form("nTracks_%s", name),
00138                                                Form("number of tracks at secondary vertex in %s", flavour),
00139                                                20, 0, 20);
00140                 plots.chi2 = fs->make<TH1F>(Form("chi2_%s", name),
00141                                             Form("secondary vertex fit \\chi^{2} in %s", flavour),
00142                                             100, 0, 50);
00143         }
00144 }
00145 
00146 // helper function to sort the tracks by impact parameter significance
00147 
00148 void PatBJetVertexAnalyzer::analyze(const edm::Event &event, const edm::EventSetup &es)
00149 {  
00150         // handle to the jets collection
00151         edm::Handle<pat::JetCollection> jetsHandle;
00152         event.getByLabel(jets_, jetsHandle);
00153 
00154         // now go through all jets
00155         for(pat::JetCollection::const_iterator jet = jetsHandle->begin();
00156             jet != jetsHandle->end(); ++jet) {
00157 
00158                 // only look at jets that pass the pt and eta cut
00159                 if (jet->pt() < jetPtCut_ ||
00160                     std::abs(jet->eta()) > jetEtaCut_)
00161                         continue;
00162 
00163                 Flavour flavour;
00164                 // find out the jet flavour (differs between quark and anti-quark)
00165                 switch(std::abs(jet->partonFlavour())) {
00166                     case 1:
00167                     case 2:
00168                     case 3:
00169                     case 21:
00170                         flavour = UDSG_JETS;
00171                         break;
00172                     case 4:
00173                         flavour = C_JETS;
00174                         break;
00175                     case 5:
00176                         flavour = B_JETS;
00177                         break;
00178                     default:
00179                         flavour = NONID_JETS;
00180                 }
00181 
00182                 // simply count the number of accepted jets
00183                 flavours_->Fill(ALL_JETS);
00184                 flavours_->Fill(flavour);
00185 
00186                 // retrieve the "secondary vertex tag infos"
00187                 // this is the output of the b-tagging reconstruction code
00188                 // and contains secondary vertices in the jets
00189                 const reco::SecondaryVertexTagInfo &svTagInfo =
00190                                         *jet->tagInfoSecondaryVertex();
00191 
00192                 // count the number of secondary vertices
00193                 plots_[ALL_JETS].nVertices->Fill(svTagInfo.nVertices());
00194                 plots_[flavour].nVertices->Fill(svTagInfo.nVertices());
00195 
00196                 // ignore jets without SV from now on
00197                 if (svTagInfo.nVertices() < 1)
00198                         continue;
00199 
00200                 // pick the first secondary vertex (the "best" one)
00201                 const reco::Vertex &sv = svTagInfo.secondaryVertex(0);
00202 
00203                 // and plot number of tracks and chi^2
00204                 plots_[ALL_JETS].nTracks->Fill(sv.tracksSize());
00205                 plots_[flavour].nTracks->Fill(sv.tracksSize());
00206 
00207                 plots_[ALL_JETS].chi2->Fill(sv.chi2());
00208                 plots_[flavour].chi2->Fill(sv.chi2());
00209 
00210                 // the precomputed transverse distance to the primary vertex
00211                 Measurement1D distance = svTagInfo.flightDistance(0, true);
00212 
00213                 plots_[ALL_JETS].dist->Fill(distance.value());
00214                 plots_[flavour].dist->Fill(distance.value());
00215 
00216                 plots_[ALL_JETS].distErr->Fill(distance.error());
00217                 plots_[flavour].distErr->Fill(distance.error());
00218 
00219                 plots_[ALL_JETS].distSig->Fill(distance.significance());
00220                 plots_[flavour].distSig->Fill(distance.significance());
00221 
00222 
00223                 // the precomputed direction with respect to the primary vertex
00224                 GlobalVector dir = svTagInfo.flightDirection(0);
00225 
00226                 // unfortunately CMSSW hsa all kinds of vectors,
00227                 // and sometimes we need to convert them *sigh*
00228                 math::XYZVector dir2(dir.x(), dir.y(), dir.z());
00229 
00230                 // compute a few variables that we are plotting
00231                 double deltaR = ROOT::Math::VectorUtil::DeltaR(
00232                                                 jet->momentum(), dir2);
00233 
00234                 plots_[ALL_JETS].deltaR->Fill(deltaR);
00235                 plots_[flavour].deltaR->Fill(deltaR);
00236 
00237                 // compute the invariant mass from a four-vector sum
00238                 math::XYZTLorentzVector trackFourVectorSum;
00239 
00240                 // loop over all tracks in the vertex
00241                 for(reco::Vertex::trackRef_iterator track = sv.tracks_begin();
00242                     track != sv.tracks_end(); ++track) {
00243                         ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec;
00244                         vec.SetPx((*track)->px());
00245                         vec.SetPy((*track)->py());
00246                         vec.SetPz((*track)->pz());
00247                         vec.SetM(0.13957);      // pion mass
00248                         trackFourVectorSum += vec;
00249                 }
00250 
00251                 // get the invariant mass: sqrt(E² - px² - py² - pz²)
00252                 double vertexMass = trackFourVectorSum.M();
00253 
00254                 plots_[ALL_JETS].mass->Fill(vertexMass);
00255                 plots_[flavour].mass->Fill(vertexMass);
00256         }
00257 }
00258         
00259 #include "FWCore/Framework/interface/MakerMacros.h"
00260 
00261 DEFINE_FWK_MODULE(PatBJetVertexAnalyzer);