CMS 3D CMS Logo

Functions
FWLiteHistograms.cc File Reference
#include <memory>
#include <string>
#include <vector>
#include <sstream>
#include <fstream>
#include <iostream>
#include <TH1F.h>
#include <TROOT.h>
#include <TFile.h>
#include <TSystem.h>
#include "DataFormats/FWLite/interface/Event.h"
#include "DataFormats/Common/interface/Handle.h"
#include "FWCore/FWLite/interface/FWLiteEnabler.h"
#include "DataFormats/MuonReco/interface/Muon.h"
#include "DataFormats/PatCandidates/interface/Muon.h"
#include "PhysicsTools/FWLite/interface/TFileService.h"
#include "PhysicsTools/FWLite/interface/CommandLineParser.h"

Go to the source code of this file.

Functions

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

Function Documentation

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

Definition at line 23 of file FWLiteHistograms.cc.

References fwlite::Event::atEnd(), gather_cfg::cout, dir, FWLiteEnabler::enable(), ev, optutl::VariableMapCont::integerValue(), TFileDirectory::make(), ResonanceBuilder::mass, TFileDirectory::mkdir(), nanoDQM_cfi::Muon, extraflags_cff::muons, optutl::CommandLineParser::parseArguments(), writedatasetfile::parser, AlCaHLTBitMon_QueryRunRegistry::string, optutl::VariableMapCont::stringValue(), optutl::VariableMapCont::stringVector(), align_cfg::TFileService, and fwlite::Event::toBegin().

24 {
25  // define what muon you are using; this is necessary as FWLite is not
26  // capable of reading edm::Views
27  using reco::Muon;
28 
29  // ----------------------------------------------------------------------
30  // First Part:
31  //
32  // * enable FWLite
33  // * book the histograms of interest
34  // * open the input file
35  // ----------------------------------------------------------------------
36 
37  // load framework libraries
38  gSystem->Load( "libFWCoreFWLite" );
40 
41  // initialize command line parser
42  optutl::CommandLineParser parser ("Analyze FWLite Histograms");
43 
44  // set defaults
45  parser.integerValue ("maxEvents" ) = 1000;
46  parser.integerValue ("outputEvery") = 10;
47  parser.stringValue ("outputFile" ) = "analyzeFWLiteHistograms.root";
48 
49  // parse arguments
50  parser.parseArguments (argc, argv);
51  int maxEvents_ = parser.integerValue("maxEvents");
52  unsigned int outputEvery_ = parser.integerValue("outputEvery");
53  std::string outputFile_ = parser.stringValue("outputFile");
54  std::vector<std::string> inputFiles_ = parser.stringVector("inputFiles");
55 
56  // book a set of histograms
58  TFileDirectory dir = fs.mkdir("analyzeBasicPat");
59  TH1F* muonPt_ = dir.make<TH1F>("muonPt" , "pt" , 100, 0., 300.);
60  TH1F* muonEta_ = dir.make<TH1F>("muonEta" , "eta" , 100, -3., 3.);
61  TH1F* muonPhi_ = dir.make<TH1F>("muonPhi" , "phi" , 100, -5., 5.);
62  TH1F* mumuMass_= dir.make<TH1F>("mumuMass", "mass", 90, 30., 120.);
63 
64  // loop the events
65  int ievt=0;
66  for(unsigned int iFile=0; iFile<inputFiles_.size(); ++iFile){
67  // open input file (can be located on castor)
68  TFile* inFile = TFile::Open(inputFiles_[iFile].c_str());
69  if( inFile ){
70  // ----------------------------------------------------------------------
71  // Second Part:
72  //
73  // * loop the events in the input file
74  // * receive the collections of interest via fwlite::Handle
75  // * fill the histograms
76  // * after the loop close the input file
77  // ----------------------------------------------------------------------
78  fwlite::Event ev(inFile);
79  for(ev.toBegin(); !ev.atEnd(); ++ev, ++ievt){
80  edm::EventBase const & event = ev;
81  // break loop if maximal number of events is reached
82  if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
83  // simple event counter
84  if(outputEvery_!=0 ? (ievt>0 && ievt%outputEvery_==0) : false)
85  std::cout << " processing event: " << ievt << std::endl;
86 
87  // Handle to the muon collection
89  event.getByLabel(std::string("muons"), muons);
90 
91  // loop muon collection and fill histograms
92  for(std::vector<Muon>::const_iterator mu1=muons->begin(); mu1!=muons->end(); ++mu1){
93  muonPt_ ->Fill( mu1->pt () );
94  muonEta_->Fill( mu1->eta() );
95  muonPhi_->Fill( mu1->phi() );
96  if( mu1->pt()>20 && fabs(mu1->eta())<2.1 ){
97  for(std::vector<Muon>::const_iterator mu2=muons->begin(); mu2!=muons->end(); ++mu2){
98  if(mu2>mu1){ // prevent double conting
99  if( mu1->charge()*mu2->charge()<0 ){ // check only muon pairs of unequal charge
100  if( mu2->pt()>20 && fabs(mu2->eta())<2.1 ){
101  mumuMass_->Fill( (mu1->p4()+mu2->p4()).mass() );
102  }
103  }
104  }
105  }
106  }
107  }
108  }
109  // close input file
110  inFile->Close();
111  }
112  // break loop if maximal number of events is reached:
113  // this has to be done twice to stop the file loop as well
114  if(maxEvents_>0 ? ievt+1>maxEvents_ : false) break;
115  }
116  return 0;
117 }
bool ev
static void enable()
enable automatic library loading
T * make(const Args &...args) const
make new ROOT object
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
dbl *** dir
Definition: mlp_gen.cc:35