Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:23:23

0001 /// This is only to make sure that our FWLite tools also compile with gcc
0002 /// that usually spots errors in a much more readable way
0003 
0004 #include "FWCore/FWLite/interface/FWLiteEnabler.h"
0005 #include "PhysicsTools/FWLite/interface/Scanner.h"
0006 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0007 #include "DataFormats/TrackReco/interface/Track.h"
0008 
0009 #include <TROOT.h>
0010 #include <TFile.h>
0011 #include <TCanvas.h>
0012 #include <TStyle.h>
0013 #include <iostream>
0014 
0015 int main(int argc, char *argv[]) {
0016   if (argc != 2) {
0017     std::cerr << "usage: " << argv[0] << " cmssw_reco_file.root" << std::endl;
0018     return 2;
0019   }
0020 
0021   FWLiteEnabler::enable();
0022   gROOT->SetStyle("Plain");
0023   gStyle->SetPalette(1);
0024   gStyle->SetHistMinimumZero(1);
0025 
0026   TFile *f = TFile::Open(argv[1]);
0027 
0028   fwlite::Event ev(f);
0029   fwlite::Scanner<std::vector<reco::Track> > sc(&ev, "generalTracks");
0030 
0031   TCanvas *c1 = new TCanvas("c1", "c1");
0032   sc.scan("pt:eta", "quality('highPurity')", 2);
0033 
0034   sc.setMaxEvents(200);
0035   double c = sc.count("quality('highPurity')");
0036   std::cout << "Found " << c << " highPurity tracks." << std::endl;
0037   double ce = sc.countEvents();
0038   std::cout << "Found " << ce << " events." << std::endl;
0039 
0040   TH1 *heta = sc.draw("eta");
0041   heta->Sumw2();
0042   heta->Scale(1.0 / sc.countEvents());
0043   c1->Print("eta.compiled.png");
0044   heta->SetLineColor(kBlue);
0045   TH1 *hp = sc.draw("eta", "quality('highPurity')", "SAME");
0046   hp->Sumw2();
0047   hp->Scale(1.0 / sc.countEvents());
0048   c1->Print("eta_and_hp.compiled.png");
0049 
0050   TH1 *hb = sc.draw("pt", "abs(eta)<1 && pt < 5", "", "hbarrel");
0051   TH1 *he = sc.draw("pt", "abs(eta)>1 && pt < 5", "", "hendcaps");
0052   hb->Draw();
0053   c1->Print("barrel_1.compiled.png");
0054   gROOT->FindObject("hbarrel")->Draw();
0055   c1->Print("barrel_2.compiled.png");
0056   he->Draw();
0057   c1->Print("endcaps_1.compiled.png");
0058   gROOT->FindObject("hendcaps")->Draw();
0059   c1->Print("endcaps_2.compiled.png");
0060 
0061   hb = sc.draw("pt", "abs(eta)<1 && pt < 5", "NORM");
0062   hb->SetLineColor(2);
0063   c1->Print("normalized_pt.compiled.png");
0064   he = sc.draw("pt", "abs(eta)>1 && pt < 5", "NORM SAME");
0065   he->SetLineColor(4);
0066   c1->Print("normalized_pts.compiled.png");
0067 
0068   sc.draw("eta", 5, -2.5, 2.5);
0069   c1->Print("eta_binned.compiled.png");
0070 
0071   double etabins[4] = {-2.5, -1, 1, 2.5};
0072   sc.draw("eta", 3, etabins, "quality('loose')");
0073   c1->Print("eta_specialbins.compiled.png");
0074 
0075   TProfile *p = sc.drawProf("eta", 5, -2.5, 2.5, "pt", "quality('loose')");
0076   p->SetLineColor(kRed);
0077   c1->Print("eta_prof_pt.compiled.png");
0078 
0079   sc.drawProf("eta", "pt", "quality('highPurity')", "SAME");
0080   c1->Print("eta_prof_pt_same.compiled.png");
0081 
0082   sc.draw2D("pt", 5, 0, 5, "eta", 4, -2.5, 2.5, "quality('highPurity')", "COLZ");
0083   c1->Print("pt_eta_2d_manual.compiled.png");
0084   sc.draw2D("pt", "eta", "quality('highPurity') && pt <= 5", "COLZ");
0085   c1->Print("pt_eta_2d_auto.compiled.png");
0086 
0087   TGraph *g = sc.drawGraph("pt", "eta", "quality('highPurity') && pt <= 5");
0088   g->SetMarkerStyle(8);
0089   c1->Print("pt_eta_2d_graph.compiled.png");
0090 
0091   // Now make dN/deta only for events that have two tracks with |eta|<1, pt > 500 MeV
0092   fwlite::ObjectCountSelector<std::vector<reco::Track> > *ntracks;
0093   ntracks =
0094       new fwlite::ObjectCountSelector<std::vector<reco::Track> >("generalTracks", "", "", "pt > 0.5 && abs(eta)<1", 2);
0095   sc.addEventSelector(ntracks);
0096   TH1 *heta2 = sc.draw("eta", 5, -2.5, 2.5);
0097   heta2->Sumw2();
0098   heta2->Scale(1.0 / sc.countEvents());
0099   heta2->SetLineColor(4);
0100   ntracks->setMin(0);
0101   heta = sc.draw("eta", "", "SAME");
0102   heta->Sumw2();
0103   heta->Scale(1.0 / sc.countEvents());
0104   c1->Print("eta_twotracks.compiled.png");
0105 
0106   sc.setMaxEvents(20);
0107   RooDataSet *ds = sc.fillDataSet(
0108       "pt:eta:@hits=hitPattern.numberOfValidHits", "@highPurity=quality('highPurity'):@highPt=pt>2", "pt > 0.5");
0109   ds->Print("v");
0110   delete ds;
0111 
0112   return 0;
0113 }