test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PatBasicFWLiteJetUnitTest.cc
Go to the documentation of this file.
1 #include <memory>
2 #include <string>
3 #include <vector>
4 #include <sstream>
5 #include <fstream>
6 #include <iostream>
7 
8 #include <TH1F.h>
9 #include <TROOT.h>
10 #include <TFile.h>
11 #include <TSystem.h>
12 
18 #include "TStopwatch.h"
19 
20 
21 int main(int argc, char* argv[])
22 {
23  // ----------------------------------------------------------------------
24  // First Part:
25  //
26  // * enable FWLite
27  // * book the histograms of interest
28  // * open the input file
29  // ----------------------------------------------------------------------
30 
31  if ( argc < 4 ) return 0;
32 
33  // load framework libraries
34  gSystem->Load( "libFWCoreFWLite" );
36 
37  // book a set of histograms
39  TFileDirectory theDir = fs.mkdir("analyzeBasicPat");
40  TH1F* jetPt_ = theDir.make<TH1F>("jetPt", "pt", 100, 0.,300.);
41  TH1F* jetEta_ = theDir.make<TH1F>("jetEta","eta", 100, -3., 3.);
42  TH1F* jetPhi_ = theDir.make<TH1F>("jetPhi","phi", 100, -5., 5.);
43  TH1F* disc_ = theDir.make<TH1F>("disc", "Discriminant", 100, 0.0, 10.0);
44  TH1F* constituentPt_ = theDir.make<TH1F>("constituentPt", "Constituent pT", 100, 0, 300.0);
45 
46  // open input file (can be located on castor)
47  TFile* inFile = TFile::Open( argv[1] );
48 
49  // ----------------------------------------------------------------------
50  // Second Part:
51  //
52  // * loop the events in the input file
53  // * receive the collections of interest via fwlite::Handle
54  // * fill the histograms
55  // * after the loop close the input file
56  // ----------------------------------------------------------------------
57 
58  // loop the events
59  unsigned int iEvent=0;
60  fwlite::Event ev(inFile);
61  TStopwatch timer;
62  timer.Start();
63 
64  unsigned int nEventsAnalyzed = 0;
65  for(ev.toBegin(); !ev.atEnd(); ++ev, ++iEvent){
66  edm::EventBase const & event = ev;
67 
68  // Handle to the jet collection
70  edm::InputTag jetLabel( argv[3] );
71  event.getByLabel(jetLabel, jets);
72 
73  // loop jet collection and fill histograms
74  for(unsigned i=0; i<jets->size(); ++i){
75  jetPt_ ->Fill( (*jets)[i].pt() );
76  jetEta_->Fill( (*jets)[i].eta() );
77  jetPhi_->Fill( (*jets)[i].phi() );
78  reco::SecondaryVertexTagInfo const * svTagInfos = (*jets)[i].tagInfoSecondaryVertex("secondaryVertex");
79  if ( svTagInfos != 0 ) {
80  if ( svTagInfos->nVertices() > 0 )
81  disc_->Fill( svTagInfos->flightDistance(0).value() );
82  }
83  std::vector<CaloTowerPtr> const & caloConstituents = (*jets)[i].getCaloConstituents();
84  for ( std::vector<CaloTowerPtr>::const_iterator ibegin = caloConstituents.begin(),
85  iend = caloConstituents.end(),
86  iconstituent = ibegin;
87  iconstituent != iend; ++iconstituent ) {
88  constituentPt_->Fill( (*iconstituent)->pt() );
89  }
90  }
91  ++nEventsAnalyzed;
92  }
93  // close input file
94  inFile->Close();
95 
96  timer.Stop();
97 
98  // print some timing statistics
99  Double_t rtime = timer.RealTime();
100  Double_t ctime = timer.CpuTime();
101  printf("Analyzed events: %d \n",nEventsAnalyzed);
102  printf("RealTime=%f seconds, CpuTime=%f seconds\n",rtime,ctime);
103  printf("%4.2f events / RealTime second .\n", (double)nEventsAnalyzed/rtime);
104  printf("%4.2f events / CpuTime second .\n", (double)nEventsAnalyzed/ctime);
105 
106 
107 
108  // ----------------------------------------------------------------------
109  // Third Part:
110  //
111  // * never forget to free the memory of objects you created
112  // ----------------------------------------------------------------------
113 
114  // in this example there is nothing to do
115 
116  // that's it!
117  return 0;
118 }
119 
Measurement1D flightDistance(unsigned int index, bool in2d=false) const
int i
Definition: DBlmapReader.cc:9
bool ev
Event const & toBegin()
Go to the very first Event.
Definition: Event.cc:242
int iEvent
Definition: GenABIO.cc:230
vector< PseudoJet > jets
static void enable()
enable automatic library loading
T * make(const Args &...args) const
make new ROOT object
virtual bool atEnd() const
Definition: Event.cc:285
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
tuple argc
Definition: dir2webdir.py:38
double value() const
Definition: Measurement1D.h:28
Geom::Phi< T > phi() const