CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TopJetFWLiteAnalyzer.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 #include <cstdlib>
4 #include <sstream>
5 #include <fstream>
6 #include <iostream>
7 
10 
15 
16 
17 int main(int argc, char* argv[])
18 {
19  if( argc<3 ){
20  // -------------------------------------------------
21  std::cerr << "ERROR:: "
22  << "Wrong number of arguments! Please specify:" << std::endl
23  << " * filepath" << std::endl
24  << " * process name" << std::endl;
25  // -------------------------------------------------
26  return -1;
27  }
28 
29  // load framework libraries
30  gSystem->Load( "libFWCoreFWLite" );
32 
33  // set nice style for histograms
34  setNiceStyle();
35 
36  // define some histograms
37  TH1I* noJets = new TH1I("noJets", "N_{Jets}", 10, 0 , 10 );
38  TH1F* ptJets = new TH1F("ptJets", "pt_{Jets}", 100, 0.,300.);
39  TH1F* enJets = new TH1F("enJets", "energy_{Jets}",100, 0.,300.);
40  TH1F* etaJets = new TH1F("etaJets","eta_{Jets}", 100, -3., 3.);
41  TH1F* phiJets = new TH1F("phiJets","phi_{Jets}", 100, -5., 5.);
42 
43  // -------------------------------------------------
44  std::cout << "open file: " << argv[1] << std::endl;
45  // -------------------------------------------------
46  TFile* inFile = TFile::Open(argv[1]);
47  TTree* events_= 0;
48  if( inFile ) inFile->GetObject("Events", events_);
49  if( events_==0 ){
50  // -------------------------------------------------
51  std::cerr << "ERROR:: "
52  << "Unable to retrieve TTree Events!" << std::endl
53  << " Eighter wrong file name or the the tree doesn't exists" << std::endl;
54  // -------------------------------------------------
55  return -1;
56  }
57 
58  // acess branch of elecs
59  char jetsName[50];
60  sprintf(jetsName, "patJets_selectedPatJets__%s.obj", argv[2]);
61  TBranch* jets_ = events_->GetBranch( jetsName ); assert( jets_!=0 );
62 
63  // loop over events and fill histograms
64  std::vector<pat::Jet> jets;
65  int nevt = events_->GetEntries();
66 
67  // -------------------------------------------------
68  std::cout << "start looping " << nevt << " events..." << std::endl;
69  // -------------------------------------------------
70  for(int evt=0; evt<nevt; ++evt){
71  // set branch address
72  jets_->SetAddress( &jets );
73  // get event
74  jets_ ->GetEntry( evt );
75  events_->GetEntry( evt, 0 );
76 
77  // -------------------------------------------------
78  if(evt>0 && !(evt%10)) std::cout << " processing event: " << evt << std::endl;
79  // -------------------------------------------------
80 
81  // fill histograms
82  noJets->Fill(jets.size());
83  for(unsigned idx=0; idx<jets.size(); ++idx){
84  // fill histograms
85  ptJets ->Fill(jets[idx].pt() );
86  enJets ->Fill(jets[idx].energy());
87  etaJets->Fill(jets[idx].eta() );
88  phiJets->Fill(jets[idx].phi() );
89  }
90  }
91  // -------------------------------------------------
92  std::cout << "close file" << std::endl;
93  // -------------------------------------------------
94  inFile->Close();
95 
96  // save histograms to file
97  TFile outFile( "analyzeJets.root", "recreate" );
98  outFile.mkdir("analyzeJet");
99  outFile.cd("analyzeJet");
100  noJets ->Write( );
101  ptJets ->Write( );
102  enJets ->Write( );
103  etaJets->Write( );
104  phiJets->Write( );
105  outFile.Close();
106 
107  // free allocated space
108  delete noJets;
109  delete ptJets;
110  delete enJets;
111  delete etaJets;
112  delete phiJets;
113 
114  return 0;
115 }
void setNiceStyle()
Definition: NiceStyle.cc:3
assert(m_qm.get())
int main(int argc, char **argv)
vector< PseudoJet > jets
static void enable()
enable automatic library loading
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
tuple argc
Definition: dir2webdir.py:38
Geom::Phi< T > phi() const
tuple cout
Definition: gather_cfg.py:121