CMS 3D CMS Logo

OptTOA.cc

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

Generated on Tue Jun 9 17:47:59 2009 for CMSSW by  doxygen 1.5.4