CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/PhysicsTools/PatExamples/plugins/PatBJetTrackAnalyzer.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 
00012 #include "FWCore/Framework/interface/Frameworkfwd.h"
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Framework/interface/EDAnalyzer.h"
00016 #include "FWCore/Utilities/interface/InputTag.h"
00017 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00018 #include "FWCore/ServiceRegistry/interface/Service.h"
00019 
00020 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00021 
00022 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00023 #include "DataFormats/VertexReco/interface/Vertex.h"
00024 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00025 #include "DataFormats/TrackReco/interface/Track.h"
00026 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00027 #include "DataFormats/PatCandidates/interface/Jet.h"
00028 #include "DataFormats/GeometryCommonDetAlgo/interface/Measurement1D.h"
00029 
00030 class PatBJetTrackAnalyzer : public edm::EDAnalyzer  {
00031     public: 
00033         PatBJetTrackAnalyzer(const edm::ParameterSet &params);
00034         ~PatBJetTrackAnalyzer();
00035 
00036         // virtual methods called from base class EDAnalyzer
00037         virtual void beginJob();
00038         virtual void analyze(const edm::Event &event, const edm::EventSetup &es);
00039 
00040     private:
00041         // configuration parameters
00042         edm::InputTag jets_;
00043         edm::InputTag tracks_;
00044         edm::InputTag beamSpot_;
00045         edm::InputTag primaryVertices_;
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         double minPt_;                  // track pt quality cut
00052         unsigned int minPixelHits_;     // minimum number of pixel hits
00053         unsigned int minTotalHits_;     // minimum number of total hits
00054 
00055         unsigned int nThTrack_;         // n-th hightest track to choose
00056 
00057         // jet flavour constants
00058 
00059         enum Flavour {
00060                 ALL_JETS = 0,
00061                 UDSG_JETS,
00062                 C_JETS,
00063                 B_JETS,
00064                 NONID_JETS,
00065                 N_JET_TYPES
00066         };
00067 
00068         TH1 *flavours_;
00069 
00070         // one group of plots per jet flavour;
00071         struct Plots {
00072                 TH1 *allIP, *allIPErr, *allIPSig;
00073                 TH1 *trackIP, *trackIPErr, *trackIPSig;
00074                 TH1 *negativeIP, *negativeIPErr, *negativeIPSig;
00075                 TH1 *nTracks, *allDeltaR;
00076         } plots_[N_JET_TYPES];
00077 };
00078 
00079 PatBJetTrackAnalyzer::PatBJetTrackAnalyzer(const edm::ParameterSet &params) :
00080         jets_(params.getParameter<edm::InputTag>("jets")),
00081         tracks_(params.getParameter<edm::InputTag>("tracks")),
00082         beamSpot_(params.getParameter<edm::InputTag>("beamSpot")),
00083         primaryVertices_(params.getParameter<edm::InputTag>("primaryVertices")),
00084         jetPtCut_(params.getParameter<double>("jetPtCut")),
00085         jetEtaCut_(params.getParameter<double>("jetEtaCut")),
00086         maxDeltaR_(params.getParameter<double>("maxDeltaR")),
00087         minPt_(params.getParameter<double>("minPt")),
00088         minPixelHits_(params.getParameter<unsigned int>("minPixelHits")),
00089         minTotalHits_(params.getParameter<unsigned int>("minTotalHits")),
00090         nThTrack_(params.getParameter<unsigned int>("nThTrack"))
00091 {
00092 }
00093 
00094 PatBJetTrackAnalyzer::~PatBJetTrackAnalyzer()
00095 {
00096 }
00097 
00098 void PatBJetTrackAnalyzer::beginJob()
00099 {
00100         // retrieve handle to auxiliary service
00101         //  used for storing histograms into ROOT file
00102         edm::Service<TFileService> fs;
00103 
00104         flavours_ = fs->make<TH1F>("flavours", "jet flavours", 5, 0, 5);
00105 
00106         // book histograms for all jet flavours
00107         for(unsigned int i = 0; i < N_JET_TYPES; i++) {
00108                 Plots &plots = plots_[i];
00109                 const char *flavour, *name;
00110 
00111                 switch((Flavour)i) {
00112                     case ALL_JETS:
00113                         flavour = "all jets";
00114                         name = "all";
00115                         break;
00116                     case UDSG_JETS:
00117                         flavour = "light flavour jets";
00118                         name = "udsg";
00119                         break;
00120                     case C_JETS:
00121                         flavour = "charm jets";
00122                         name = "c";
00123                         break;
00124                     case B_JETS:
00125                         flavour = "bottom jets";
00126                         name = "b";
00127                         break;
00128                     default:
00129                         flavour = "unidentified jets";
00130                         name = "ni";
00131                         break;
00132                 }
00133 
00134                 plots.allIP = fs->make<TH1F>(Form("allIP_%s", name),
00135                                              Form("signed IP for all tracks in %s", flavour),
00136                                              100, -0.1, 0.2);
00137                 plots.allIPErr = fs->make<TH1F>(Form("allIPErr_%s", name),
00138                                                 Form("error of signed IP for all tracks in %s", flavour),
00139                                                 100, 0, 0.05);
00140                 plots.allIPSig = fs->make<TH1F>(Form("allIPSig_%s", name),
00141                                                 Form("signed IP significance for all tracks in %s", flavour),
00142                                                 100, -10, 20);
00143 
00144                 plots.trackIP = fs->make<TH1F>(Form("trackIP_%s", name),
00145                                                Form("signed IP for selected positive track in %s", flavour),
00146                                                100, -0.1, 0.2);
00147                 plots.trackIPErr = fs->make<TH1F>(Form("trackIPErr_%s", name),
00148                                                   Form("error of signed IP for selected positive track in %s", flavour),
00149                                                   100, 0, 0.05);
00150                 plots.trackIPSig = fs->make<TH1F>(Form("trackIPSig_%s", name),
00151                                                   Form("signed IP significance for selected positive track in %s", flavour),
00152                                                   100, -10, 20);
00153 
00154                 plots.negativeIP = fs->make<TH1F>(Form("negativeIP_%s", name),
00155                                                   Form("signed IP for selected negative track in %s", flavour),
00156                                                   100, -0.2, 0.1);
00157                 plots.negativeIPErr = fs->make<TH1F>(Form("negativeIPErr_%s", name),
00158                                                      Form("error of signed IP for selected negative track in %s", flavour),
00159                                                      100, 0, 0.05);
00160                 plots.negativeIPSig = fs->make<TH1F>(Form("negativeIPSig_%s", name),
00161                                                      Form("signed IP significance for selected negative track in %s", flavour),
00162                                                      100, -20, 10);
00163 
00164                 plots.nTracks = fs->make<TH1F>(Form("nTracks_%s", name),
00165                                                Form("number of usable tracks in %s", flavour),
00166                                                30, 0, 30);
00167                 plots.allDeltaR = fs->make<TH1F>(Form("allDeltaR_%s", name),
00168                                                  Form("\\DeltaR between track and %s", flavour),
00169                                                  100, 0, 1);
00170         }
00171 }
00172 
00173 // helper function to sort the tracks by impact parameter significance
00174 
00175 static bool significanceHigher(const Measurement1D &meas1,
00176                                const Measurement1D &meas2)
00177 { return meas1.significance() > meas2.significance(); }
00178 
00179 void PatBJetTrackAnalyzer::analyze(const edm::Event &event, const edm::EventSetup &es)
00180 {  
00181         // handle to the primary vertex collection
00182         edm::Handle<reco::VertexCollection> pvHandle;
00183         event.getByLabel(primaryVertices_, pvHandle);
00184 
00185         // handle to the tracks collection
00186         edm::Handle<reco::TrackCollection> tracksHandle;
00187         event.getByLabel(tracks_, tracksHandle);
00188 
00189         // handle to the jets collection
00190         edm::Handle<pat::JetCollection> jetsHandle;
00191         event.getByLabel(jets_, jetsHandle);
00192 
00193         // handle to the beam spot
00194         edm::Handle<reco::BeamSpot> beamSpot;
00195         event.getByLabel(beamSpot_, beamSpot);
00196 
00197         // rare case of no reconstructed primary vertex
00198         if (pvHandle->empty())
00199                 return;
00200 
00201         // extract the position of the (most probable) reconstructed vertex
00202         math::XYZPoint pv = (*pvHandle)[0].position();
00203 
00204         // now go through all jets
00205         for(pat::JetCollection::const_iterator jet = jetsHandle->begin();
00206             jet != jetsHandle->end(); ++jet) {
00207 
00208                 // only look at jets that pass the pt and eta cut
00209                 if (jet->pt() < jetPtCut_ ||
00210                     std::abs(jet->eta()) > jetEtaCut_)
00211                         continue;
00212 
00213                 Flavour flavour;
00214                 // find out the jet flavour (differs between quark and anti-quark)
00215                 switch(std::abs(jet->partonFlavour())) {
00216                     case 1:
00217                     case 2:
00218                     case 3:
00219                     case 21:
00220                         flavour = UDSG_JETS;
00221                         break;
00222                     case 4:
00223                         flavour = C_JETS;
00224                         break;
00225                     case 5:
00226                         flavour = B_JETS;
00227                         break;
00228                     default:
00229                         flavour = NONID_JETS;
00230                 }
00231         
00232                 // simply count the number of accepted jets
00233                 flavours_->Fill(ALL_JETS);
00234                 flavours_->Fill(flavour);
00235 
00236                 // this vector will contain IP value / error pairs
00237                 std::vector<Measurement1D> ipValErr;
00238 
00239                 // Note: PAT is also able to store associated tracks
00240                 //       within the jet object, so we don't have to do the
00241                 //       matching ourselves
00242                 // (see ->associatedTracks() method)
00243                 // However, using this we can't play with the DeltaR cone
00244                 // withour rerunning the PAT producer
00245 
00246                 // now loop through all tracks
00247                 for(reco::TrackCollection::const_iterator track = tracksHandle->begin();
00248                     track != tracksHandle->end(); ++track) {
00249 
00250                         // check the quality criteria
00251                         if (track->pt() < minPt_ ||
00252                             track->hitPattern().numberOfValidHits() < (int)minTotalHits_ ||
00253                             track->hitPattern().numberOfValidPixelHits() < (int)minPixelHits_)
00254                                 continue;
00255 
00256                         // check the Delta R between jet axis and track
00257                         // (Delta_R^2 = Delta_Eta^2 + Delta_Phi^2)
00258                         double deltaR = ROOT::Math::VectorUtil::DeltaR(
00259                                         jet->momentum(), track->momentum());
00260 
00261                         plots_[ALL_JETS].allDeltaR->Fill(deltaR);
00262                         plots_[flavour].allDeltaR->Fill(deltaR);
00263 
00264                         // only look at tracks in jet cone
00265                         if (deltaR > maxDeltaR_)
00266                                 continue;
00267 
00268                         // What follows here is an approximation!
00269                         //
00270                         // The dxy() method of the tracks does a linear
00271                         // extrapolation from the track reference position
00272                         // given as the closest point to the beam spot
00273                         // with respect to the given vertex.
00274                         // Since we are using primary vertices, this
00275                         // approximation works well
00276                         //
00277                         // In order to get better results, the
00278                         // "TransientTrack" and specialised methods have
00279                         // to be used.
00280                         // Or look at the "impactParameterTagInfos",
00281                         // which contains the precomputed information
00282                         // from the official b-tagging algorithms
00283                         //
00284                         // see ->tagInfoTrackIP() method
00285 
00286                         double ipError = track->dxyError();
00287                         double ipValue = std::abs(track->dxy(pv));
00288 
00289                         // in order to compute the sign, we check if
00290                         // the point of closest approach to the vertex
00291                         // is in front or behind the vertex.
00292                         // Again, we a linear approximation
00293                         // 
00294                         // dot product between reference point and jet axis
00295 
00296                         math::XYZVector closestPoint = track->referencePoint() - beamSpot->position();
00297                         // only interested in transverse component, z -> 0
00298                         closestPoint.SetZ(0.);
00299                         double sign = closestPoint.Dot(jet->momentum());
00300 
00301                         if (sign < 0)
00302                                 ipValue = -ipValue;
00303 
00304                         ipValErr.push_back(Measurement1D(ipValue, ipError));
00305                 }
00306 
00307                 // now order all tracks by significance (highest first)
00308                 std::sort(ipValErr.begin(), ipValErr.end(), significanceHigher);
00309 
00310                 plots_[ALL_JETS].nTracks->Fill(ipValErr.size());
00311                 plots_[flavour].nTracks->Fill(ipValErr.size());
00312 
00313                 // plot all tracks
00314 
00315                 for(std::vector<Measurement1D>::const_iterator iter = ipValErr.begin();
00316                     iter != ipValErr.end(); ++iter) {
00317                         plots_[ALL_JETS].allIP->Fill(iter->value());
00318                         plots_[flavour].allIP->Fill(iter->value());
00319 
00320                         plots_[ALL_JETS].allIPErr->Fill(iter->error());
00321                         plots_[flavour].allIPErr->Fill(iter->error());
00322 
00323                         // significance (is really just value / error)
00324                         plots_[ALL_JETS].allIPSig->Fill(iter->significance());
00325                         plots_[flavour].allIPSig->Fill(iter->significance());
00326                 }
00327 
00328                 // check if we have enough tracks to fulfill the
00329                 // n-th track requirement
00330                 if (ipValErr.size() < nThTrack_)
00331                         continue;
00332 
00333                 // pick the n-th highest track
00334                 const Measurement1D *trk = &ipValErr[nThTrack_ - 1];
00335 
00336                 plots_[ALL_JETS].trackIP->Fill(trk->value());
00337                 plots_[flavour].trackIP->Fill(trk->value());
00338 
00339                 plots_[ALL_JETS].trackIPErr->Fill(trk->error());
00340                 plots_[flavour].trackIPErr->Fill(trk->error());
00341 
00342                 plots_[ALL_JETS].trackIPSig->Fill(trk->significance());
00343                 plots_[flavour].trackIPSig->Fill(trk->significance());
00344 
00345                 // here we define a "negative tagger", i.e. we take
00346                 // the track with the n-lowest signed IP
00347                 // (i.e. preferrably select tracks that appear to become
00348                 //  from "behind" the primary vertex, supposedly mismeasured
00349                 //  tracks really coming from the primary vertex, and
00350                 //  the contamination of displaced tracks should be small)
00351                 trk = &ipValErr[ipValErr.size() - nThTrack_];
00352 
00353                 plots_[ALL_JETS].negativeIP->Fill(trk->value());
00354                 plots_[flavour].negativeIP->Fill(trk->value());
00355 
00356                 plots_[ALL_JETS].negativeIPErr->Fill(trk->error());
00357                 plots_[flavour].negativeIPErr->Fill(trk->error());
00358 
00359                 plots_[ALL_JETS].negativeIPSig->Fill(trk->significance());
00360                 plots_[flavour].negativeIPSig->Fill(trk->significance());
00361         }
00362 }
00363         
00364 #include "FWCore/Framework/interface/MakerMacros.h"
00365 
00366 DEFINE_FWK_MODULE(PatBJetTrackAnalyzer);