CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/SimTracker/TrackHistory/plugins/QualityCutsAnalyzer.cc

Go to the documentation of this file.
00001 #include <algorithm>
00002 #include <cctype>
00003 #include <iomanip>
00004 #include <set>
00005 #include <sstream>
00006 #include <vector>
00007 
00008 #include "TFile.h"
00009 #include "TH1F.h"
00010 
00011 #include "HepPDT/ParticleID.hh"
00012 
00013 // user include files
00014 #include "DataFormats/JetReco/interface/JetTracksAssociation.h"
00015 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00016 #include "DataFormats/VertexReco/interface/Vertex.h"
00017 
00018 #include "FWCore/Framework/interface/Frameworkfwd.h"
00019 #include "FWCore/Framework/interface/EDAnalyzer.h"
00020 #include "FWCore/Framework/interface/ESHandle.h"
00021 #include "FWCore/Framework/interface/MakerMacros.h"
00022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00023 
00024 #include "TrackingTools/IPTools/interface/IPTools.h"
00025 
00026 #include "RecoVertex/PrimaryVertexProducer/interface/PrimaryVertexSorter.h"
00027 
00028 #include "SimTracker/TrackHistory/interface/TrackClassifier.h"
00029 
00030 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00031 
00032 //
00033 // class decleration
00034 //
00035 
00036 class QualityCutsAnalyzer : public edm::EDAnalyzer
00037 {
00038 
00039 public:
00040 
00041     explicit QualityCutsAnalyzer(const edm::ParameterSet&);
00042     ~QualityCutsAnalyzer();
00043 
00044 private:
00045 
00046     virtual void beginJob() ;
00047     virtual void analyze(const edm::Event&, const edm::EventSetup&);
00048     virtual void endJob();
00049 
00050     // Member data
00051 
00052     typedef std::vector<int> vint;
00053     typedef std::vector<std::string> vstring;
00054 
00055     edm::InputTag trackProducer_;
00056     edm::InputTag primaryVertexProducer_;
00057     edm::InputTag jetTracksAssociation_;
00058 
00059     std::string rootFile_;
00060 
00061     int minimumNumberOfHits_, minimumNumberOfPixelHits_;
00062     double minimumTransverseMomentum_, maximumChiSquared_;
00063 
00064     bool useAllQualities_;
00065     reco::TrackBase::TrackQuality trackQuality_;
00066 
00067     void
00068     LoopOverJetTracksAssociation(
00069         const edm::ESHandle<TransientTrackBuilder> &,
00070         const edm::Handle<reco::VertexCollection> &,
00071         const edm::Handle<reco::JetTracksAssociationCollection> &
00072     );
00073 
00074     // Histograms for optimization
00075 
00076     struct histogram_element_t
00077     {
00078         double sdl;  // Signed decay length
00079         double dta;  // Distance to jet axis
00080         double tip;  // Transverse impact parameter
00081         double lip;  // Longitudinal impact parameter
00082         double ips;  // Impact parameter significance.
00083         double pt;   // Transverse momentum
00084         double chi2; // Chi^2
00085         std::size_t hits;      // Number of hits
00086         std::size_t pixelhits; // Number of hits
00087 
00088         histogram_element_t(double d, double a, double t, double l, double i, double p, double c, std::size_t h, std::size_t x)
00089         {
00090             sdl = d;
00091             dta = a;
00092             tip = t;
00093             lip = l;
00094             ips = i;
00095             pt = p;
00096             chi2 = c;
00097             hits = h;
00098             pixelhits = x;
00099         }
00100 
00101         histogram_element_t(const histogram_element_t & orig)
00102         {
00103             sdl = orig.sdl;
00104             dta = orig.dta;
00105             tip = orig.tip;
00106             lip = orig.lip;
00107             ips = orig.ips;
00108             pt = orig.pt;
00109             chi2 = orig.chi2;
00110             hits = orig.hits;
00111             pixelhits = orig.pixelhits;
00112         }
00113     };
00114 
00115     typedef std::vector<std::vector<histogram_element_t> > histogram_data_t;
00116     histogram_data_t histogram_data_;
00117 
00118     class histogram_t
00119     {
00120 
00121         TH1F* sdl;
00122         TH1F* dta;
00123         TH1F* tip;
00124         TH1F* lip;
00125         TH1F* ips;
00126         TH1F* pixelhits;
00127         TH1F* pt_1gev;
00128         TH1F* chi2;
00129         TH1F* hits;
00130 
00131     public:
00132 
00133         histogram_t(const std::string & particleType)
00134         {
00135             std::string name, title;
00136             name = std::string("hits_") + particleType;
00137             title = std::string("Hit distribution for ") + particleType;
00138             hits = new TH1F(name.c_str(), title.c_str(), 19, -0.5, 18.5);
00139 
00140             name = std::string("chi2_") + particleType;
00141             title = std::string("Chi2 distribution for ") + particleType;
00142             chi2 = new TH1F(name.c_str(), title.c_str(), 100, 0., 30.);
00143 
00144             name = std::string("pixelhits_") + particleType;
00145             title = std::string("Pixel hits distribution for ") + particleType;
00146             pixelhits = new TH1F(name.c_str(), title.c_str(), 21, -0.5, 20.5);
00147 
00148             name = std::string("pt_1Gev_") + particleType;
00149             title = std::string("Pt distribution close 1Gev for ") + particleType;
00150             pt_1gev = new TH1F(name.c_str(), title.c_str(), 100, 0., 2.);
00151 
00152             name = std::string("tip_") + particleType;
00153             title = std::string("Transverse impact parameter distribution for ") + particleType;
00154             tip = new TH1F(name.c_str(), title.c_str(), 100, -0.3, 0.3);
00155 
00156             name = std::string("lip_") + particleType;
00157             title = std::string("Longitudinal impact parameter distribution for ") + particleType;
00158             lip = new TH1F(name.c_str(), title.c_str(), 100, -1., 1.);
00159 
00160             name = std::string("ips_") + particleType;
00161             title = std::string("IPS distribution for ") + particleType;
00162             ips = new TH1F(name.c_str(), title.c_str(), 100, -25.0, 25.0);
00163 
00164             name = std::string("sdl_") + particleType;
00165             title = std::string("Decay length distribution for ") + particleType;
00166             sdl = new TH1F(name.c_str(), title.c_str(), 100, -5., 5.);
00167 
00168             name = std::string("dta_") + particleType;
00169             title = std::string("Distance to jet distribution for ") + particleType;
00170             dta = new TH1F(name.c_str(), title.c_str(), 100, 0.0, 0.2);
00171         }
00172 
00173         ~histogram_t()
00174         {
00175             delete hits;
00176             delete chi2;
00177             delete pixelhits;
00178             delete pt_1gev;
00179             delete tip;
00180             delete lip;
00181             delete ips;
00182             delete sdl;
00183             delete dta;
00184         }
00185 
00186         void Fill(const histogram_element_t & data)
00187         {
00188             hits->Fill(data.hits);
00189             chi2->Fill(data.chi2);
00190             pixelhits->Fill(data.pt);
00191             pt_1gev->Fill(data.pt);
00192             ips->Fill(data.ips);
00193             tip->Fill(data.tip);
00194             lip->Fill(data.lip);
00195             sdl->Fill(data.sdl);
00196             dta->Fill(data.dta);
00197         }
00198 
00199         void Write()
00200         {
00201             hits->Write();
00202             chi2->Write();
00203             pixelhits->Write();
00204             pt_1gev->Write();
00205             ips->Write();
00206             tip->Write();
00207             lip->Write();
00208             sdl->Write();
00209             dta->Write();
00210         }
00211     };
00212 
00213     // Track classification.
00214     TrackClassifier classifier_;
00215 
00216 };
00217 
00218 
00219 //
00220 // constructors and destructor
00221 //
00222 QualityCutsAnalyzer::QualityCutsAnalyzer(const edm::ParameterSet& config) : classifier_(config)
00223 {
00224     trackProducer_         = config.getUntrackedParameter<edm::InputTag> ( "trackProducer" );
00225     primaryVertexProducer_ = config.getUntrackedParameter<edm::InputTag> ( "primaryVertexProducer" );
00226     jetTracksAssociation_  = config.getUntrackedParameter<edm::InputTag> ( "jetTracksAssociation" );
00227 
00228     rootFile_ = config.getUntrackedParameter<std::string> ( "rootFile" );
00229 
00230     minimumNumberOfHits_       = config.getUntrackedParameter<int> ( "minimumNumberOfHits" );
00231     minimumNumberOfPixelHits_  = config.getUntrackedParameter<int> ( "minimumNumberOfPixelHits" );
00232     minimumTransverseMomentum_ = config.getUntrackedParameter<double> ( "minimumTransverseMomentum" );
00233     maximumChiSquared_         = config.getUntrackedParameter<double> ( "maximumChiSquared" );
00234 
00235     std::string trackQualityType = config.getUntrackedParameter<std::string>("trackQualityClass"); //used
00236     trackQuality_ =  reco::TrackBase::qualityByName(trackQualityType);
00237     useAllQualities_ = false;
00238 
00239     std::transform(trackQualityType.begin(), trackQualityType.end(), trackQualityType.begin(), (int(*)(int)) std::tolower);
00240     if (trackQualityType == "any")
00241     {
00242         std::cout << "Using any" << std::endl;
00243         useAllQualities_ = true;
00244     }
00245 }
00246 
00247 QualityCutsAnalyzer::~QualityCutsAnalyzer() {}
00248 
00249 //
00250 // member functions
00251 //
00252 
00253 // ------------ method called to for each event  ------------
00254 void
00255 QualityCutsAnalyzer::analyze(const edm::Event& event, const edm::EventSetup& setup)
00256 {
00257     // Track collection
00258     edm::Handle<edm::View<reco::Track> > trackCollection;
00259     event.getByLabel(trackProducer_,trackCollection);
00260     // Primary vertex
00261     edm::Handle<reco::VertexCollection> primaryVertexCollection;
00262     event.getByLabel(primaryVertexProducer_, primaryVertexCollection);
00263     // Jet to tracks associator
00264     edm::Handle<reco::JetTracksAssociationCollection> jetTracks;
00265     event.getByLabel(jetTracksAssociation_, jetTracks);
00266     // Trasient track builder
00267     edm::ESHandle<TransientTrackBuilder> TTbuilder;
00268     setup.get<TransientTrackRecord>().get("TransientTrackBuilder", TTbuilder);
00269 
00270     // Setting up event information for the track categories.
00271     classifier_.newEvent(event, setup);
00272 
00273     LoopOverJetTracksAssociation(
00274         TTbuilder,
00275         primaryVertexCollection,
00276         jetTracks
00277     );
00278 }
00279 
00280 
00281 // ------------ method called once each job just before starting event loop  ------------
00282 void
00283 QualityCutsAnalyzer::beginJob()
00284 {
00285     histogram_data_.resize(6);
00286 }
00287 
00288 
00289 // ------------ method called once each job just after ending the event loop  ------------
00290 void
00291 QualityCutsAnalyzer::endJob()
00292 {
00293     TFile file(rootFile_.c_str(), "RECREATE");
00294     file.cd();
00295 
00296     // saving the histograms
00297     for (std::size_t i=0; i<6; i++)
00298     {
00299         std::string particle;
00300         if (i == 0)
00301             particle = std::string("B_tracks");
00302         else if (i == 1)
00303             particle = std::string("C_tracks");
00304         else if (i == 2)
00305             particle = std::string("nonB_tracks");
00306         else if (i == 3)
00307             particle = std::string("displaced_tracks");
00308         else if (i == 4)
00309             particle = std::string("bad_tracks");
00310         else
00311             particle = std::string("fake_tracks");
00312 
00313         histogram_t histogram(particle);
00314 
00315         for (std::size_t j=0; j<histogram_data_[i].size(); j++)
00316             histogram.Fill(histogram_data_[i][j]);
00317 
00318         histogram.Write();
00319     }
00320 
00321     file.Flush();
00322 }
00323 
00324 
00325 void
00326 QualityCutsAnalyzer::LoopOverJetTracksAssociation(
00327     const edm::ESHandle<TransientTrackBuilder> & TTbuilder,
00328     const edm::Handle<reco::VertexCollection> & primaryVertexProducer_,
00329     const edm::Handle<reco::JetTracksAssociationCollection> & jetTracksAssociation
00330 )
00331 {
00332     const TransientTrackBuilder * bproduct = TTbuilder.product();
00333 
00334     // getting the primary vertex
00335     // use first pv of the collection
00336     reco::Vertex pv;
00337 
00338     if (primaryVertexProducer_->size() != 0)
00339     {
00340         PrimaryVertexSorter pvs;
00341         std::vector<reco::Vertex> sortedList = pvs.sortedList(*(primaryVertexProducer_.product()));
00342         pv = (sortedList.front());
00343     }
00344     else
00345     { // create a dummy PV
00346         // cout << "NO PV FOUND" << endl;
00347         reco::Vertex::Error e;
00348         e(0,0)=0.0015*0.0015;
00349         e(1,1)=0.0015*0.0015;
00350         e(2,2)=15.*15.;
00351         reco::Vertex::Point p(0,0,0);
00352         pv = reco::Vertex(p,e,1,1,1);
00353     }
00354 
00355     reco::JetTracksAssociationCollection::const_iterator it = jetTracksAssociation->begin();
00356 
00357     int i=0;
00358 
00359     for (; it != jetTracksAssociation->end(); it++, i++)
00360     {
00361         // get jetTracks object
00362         reco::JetTracksAssociationRef jetTracks(jetTracksAssociation, i);
00363 
00364         double pvZ = pv.z();
00365         GlobalVector direction(jetTracks->first->px(), jetTracks->first->py(), jetTracks->first->pz());
00366 
00367         // get the tracks associated to the jet
00368         reco::TrackRefVector tracks = jetTracks->second;
00369         for (std::size_t index = 0; index < tracks.size(); index++)
00370         {
00371             edm::RefToBase<reco::Track> track(tracks[index]);
00372 
00373             double pt = tracks[index]->pt();
00374             double chi2 = tracks[index]->normalizedChi2();
00375             int hits = tracks[index]->hitPattern().numberOfValidHits();
00376             int pixelHits = tracks[index]->hitPattern().numberOfValidPixelHits();
00377 
00378             if (
00379                 hits < minimumNumberOfHits_ ||
00380                 pixelHits < minimumNumberOfPixelHits_ ||
00381                 pt < minimumTransverseMomentum_ ||
00382                 chi2 >  maximumChiSquared_ ||
00383                 (!useAllQualities_ && !tracks[index]->quality(trackQuality_))
00384             ) continue;
00385 
00386             const reco::TransientTrack transientTrack = bproduct->build(&(*tracks[index]));
00387             double dta = - IPTools::jetTrackDistance(transientTrack, direction, pv).second.value();
00388             double sdl = IPTools::signedDecayLength3D(transientTrack, direction, pv).second.value();
00389             double ips = IPTools::signedImpactParameter3D(transientTrack, direction, pv).second.value();
00390             double d0 = IPTools::signedTransverseImpactParameter(transientTrack, direction, pv).second.value();
00391             double dz = tracks[index]->dz() - pvZ;
00392 
00393             // Classify the reco track;
00394             classifier_.evaluate( edm::RefToBase<reco::Track>(tracks[index]) );
00395 
00396             // Check for the different categories
00397             if ( classifier_.is(TrackClassifier::Fake) )
00398                 histogram_data_[5].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
00399             else if ( classifier_.is(TrackClassifier::BWeakDecay) )
00400                 histogram_data_[0].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
00401             else if ( classifier_.is(TrackClassifier::Bad) )
00402                 histogram_data_[4].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
00403             else if (
00404                 !classifier_.is(TrackClassifier::CWeakDecay) &&
00405                 !classifier_.is(TrackClassifier::PrimaryVertex)
00406             )
00407                 histogram_data_[3].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
00408             else if ( classifier_.is(TrackClassifier::CWeakDecay) )
00409                 histogram_data_[1].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
00410             else
00411                 histogram_data_[2].push_back(histogram_element_t(sdl, dta, d0, dz, ips, pt, chi2, hits, pixelHits));
00412 
00413         }
00414     }
00415 }
00416 
00417 DEFINE_FWK_MODULE(QualityCutsAnalyzer);