CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions
TopElecFWLiteAnalyzer.cc File Reference
#include <memory>
#include <string>
#include <cstdlib>
#include <sstream>
#include <fstream>
#include <iostream>
#include "FWCore/FWLite/interface/FWLiteEnabler.h"
#include "DataFormats/PatCandidates/interface/Electron.h"
#include "TopQuarkAnalysis/Examples/bin/NiceStyle.cc"
#include "TopQuarkAnalysis/Examples/interface/RootSystem.h"
#include "TopQuarkAnalysis/Examples/interface/RootHistograms.h"
#include "TopQuarkAnalysis/Examples/interface/RootPostScript.h"

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

int main ( int  argc,
char *  argv[] 
)

CALO JETS


PF JETS

Definition at line 17 of file TopElecFWLiteAnalyzer.cc.

References assert(), ecal_dqm_sourceclient-live_cfg::cerr, gather_cfg::cout, FWLiteEnabler::enable(), relval_parameters_module::energy, eta, customizeTrackingMonitorSeedNumber::idx, nevt, GetRecoTauVFromDQM_MC_cff::outFile, phi(), EnergyCorrector::pt, and setNiceStyle().

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* noElecs = new TH1I("noElecs", "N_{Elecs}", 10, 0 , 10 );
38  TH1F* ptElecs = new TH1F("ptElecs", "pt_{Elecs}", 100, 0.,300.);
39  TH1F* enElecs = new TH1F("enElecs", "energy_{Elecs}",100, 0.,300.);
40  TH1F* etaElecs = new TH1F("etaElecs","eta_{Elecs}", 100, -3., 3.);
41  TH1F* phiElecs = new TH1F("phiElecs","phi_{Elecs}", 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 elecName[50];
60  sprintf(elecName, "patElectrons_selectedPatElectrons__%s.obj", argv[2]);
61  TBranch* elecs_ = events_->GetBranch( elecName ); assert( elecs_!=0 );
62 
63  // loop over events and fill histograms
64  std::vector<pat::Electron> elecs;
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  elecs_->SetAddress( &elecs );
73  // get event
74  elecs_ ->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  noElecs->Fill(elecs.size());
83  for(unsigned idx=0; idx<elecs.size(); ++idx){
84  // fill histograms
85  ptElecs ->Fill(elecs[idx].pt() );
86  enElecs ->Fill(elecs[idx].energy());
87  etaElecs->Fill(elecs[idx].eta() );
88  phiElecs->Fill(elecs[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( "analyzeElecs.root", "recreate" );
98  outFile.mkdir("analyzeElec");
99  outFile.cd("analyzeElec");
100  noElecs ->Write( );
101  ptElecs ->Write( );
102  enElecs ->Write( );
103  etaElecs->Write( );
104  phiElecs->Write( );
105  outFile.Close();
106 
107  // free allocated space
108  delete noElecs;
109  delete ptElecs;
110  delete enElecs;
111  delete etaElecs;
112  delete phiElecs;
113 
114  return 0;
115 }
void setNiceStyle()
Definition: NiceStyle.cc:3
assert(m_qm.get())
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