CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Functions
PatCOCExercise.cc File Reference
#include <vector>
#include "TH1.h"
#include "TFile.h"
#include <TROOT.h>
#include <TSystem.h>
#include "DataFormats/Math/interface/deltaR.h"
#include "DataFormats/FWLite/interface/Event.h"
#include "DataFormats/Common/interface/Handle.h"
#include "DataFormats/PatCandidates/interface/Jet.h"
#include "FWCore/ParameterSet/interface/ProcessDesc.h"
#include "FWCore/FWLite/interface/AutoLibraryLoader.h"
#include "PhysicsTools/FWLite/interface/TFileService.h"
#include "FWCore/PythonParameterSet/interface/PythonProcessDesc.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 PatCOCExercise.cc.

References fwlite::Event::atEnd(), edm::PtrVector< T >::begin(), gather_cfg::cout, deltaR(), reco::deltaR(), AutoLibraryLoader::enable(), edm::PtrVector< T >::end(), edm::ParameterSet::getParameter(), iEvent, metsig::jet, fwrapper::jets, TFileDirectory::make(), TFileDirectory::mkdir(), muon::overlap(), analyzePatCOC_cfg::overlaps, PythonProcessDesc::processDesc(), edm::PtrVectorBase::size(), TrackerOfflineValidation_Standalone_cff::TFileService, and fwlite::Event::toBegin().

23  {
24  // ----------------------------------------------------------------------
25  // First Part:
26  //
27  // * enable the AutoLibraryLoader
28  // * book the histograms of interest
29  // * open the input file
30  // ----------------------------------------------------------------------
31 
32  // load framework libraries
33  gSystem->Load( "libFWCoreFWLite" );
35 
36  // only allow one argument for this simple example which should be the
37  // the python cfg file
38  if ( argc < 2 ) {
39  std::cout << "Usage : " << argv[0] << " [parameters.py]" << std::endl;
40  return 0;
41  }
42 
43  // get the python configuration
44  PythonProcessDesc builder(argv[1]);
45  const edm::ParameterSet& fwliteParameters = builder.processDesc()->getProcessPSet()->getParameter<edm::ParameterSet>("FWLiteParams");
46 
47  // now get each parameter
48  std::string input_ ( fwliteParameters.getParameter<std::string >("inputFile" ) );
49  std::string output_ ( fwliteParameters.getParameter<std::string >("outputFile" ) );
50  std::string overlaps_( fwliteParameters.getParameter<std::string >("overlaps") );
51  edm::InputTag jets_ ( fwliteParameters.getParameter<edm::InputTag>("jets" ) );
52 
53  // book a set of histograms
54  fwlite::TFileService fs = fwlite::TFileService(output_.c_str());
55  TFileDirectory theDir = fs.mkdir("analyzePatCOC");
56  TH1F* deltaRElecJet_ = theDir.make<TH1F>("deltaRElecJet" , "#DeltaR (elec, jet)" , 10, 0., 0.5);
57  TH1F* elecOverJet_ = theDir.make<TH1F>("elecOverJet" , "E_{elec}/E_{jet}" , 100, 0., 2.);
58  TH1F* nOverlaps_ = theDir.make<TH1F>("nOverlaps" , "Number of overlaps" , 5, 0., 5.);
59 
60  // open input file (can be located on castor)
61  TFile* inFile = TFile::Open(input_.c_str());
62 
63  // ----------------------------------------------------------------------
64  // Second Part:
65  //
66  // * loop the events in the input file
67  // * receive the collections of interest via fwlite::Handle
68  // * fill the histograms
69  // * after the loop close the input file
70  // ----------------------------------------------------------------------
71 
72  // loop the events
73  unsigned int iEvent=0;
74  fwlite::Event ev(inFile);
75  for(ev.toBegin(); !ev.atEnd(); ++ev, ++iEvent){
76  edm::EventBase const & event = ev;
77 
78  // break loop after end of file is reached
79  // or after 1000 events have been processed
80  if( iEvent==1000 ) break;
81 
82  // simple event counter
83  if(iEvent>0 && iEvent%1==0){
84  std::cout << " processing event: " << iEvent << std::endl;
85  }
86 
87  // handle to jet collection
89  event.getByLabel(jets_, jets);
90 
91  // loop over the jets in the event
92  for( std::vector<pat::Jet>::const_iterator jet = jets->begin(); jet != jets->end(); jet++ ){
93  if(jet->pt()>20 && jet==jets->begin()){
94 
95  if(jet->hasOverlaps(overlaps_)){
96  //get all overlaps
97  const reco::CandidatePtrVector overlaps = jet->overlaps(overlaps_);
98  nOverlaps_->Fill( overlaps.size() );
99  //loop over the overlaps
100  for( reco::CandidatePtrVector::const_iterator overlap = overlaps.begin(); overlap != overlaps.end(); overlap++){
101  float deltaR = reco::deltaR( (*overlap)->eta(), (*overlap)->phi(), jet->eta(), jet->phi() );
102  deltaRElecJet_->Fill( deltaR );
103  elecOverJet_->Fill( (*overlap)->energy()/jet->energy() );
104  }
105  }
106  }
107  }
108  }
109  inFile->Close();
110  return 0;
111 }
T getParameter(std::string const &) const
size_type size() const
Size of the RefVector.
Definition: PtrVectorBase.h:73
double deltaR(double eta1, double phi1, double eta2, double phi2)
Definition: deltaR.h:19
const_iterator begin() const
Definition: PtrVector.h:126
bool overlap(const reco::Muon &muon1, const reco::Muon &muon2, double pullX=1.0, double pullY=1.0, bool checkAdjacentChambers=false)
int iEvent
Definition: GenABIO.cc:243
vector< PseudoJet > jets
const_iterator end() const
Definition: PtrVector.h:131
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
tuple argc
Definition: dir2webdir.py:41
T * make() const
make new ROOT object
static void enable()
enable automatic library loading
tuple cout
Definition: gather_cfg.py:121