PhysicsTools
PatExamples
bin
JPsiFWLiteAnalyzer.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
13
#include "
DataFormats/FWLite/interface/Handle.h
"
14
#include "
DataFormats/PatCandidates/interface/Muon.h
"
15
#include "
DataFormats/PatCandidates/interface/CompositeCandidate.h
"
16
#include "
FWCore/FWLite/interface/FWLiteEnabler.h
"
17
18
using namespace
std
;
19
20
int
main
(
int
argc
,
char
*
argv
[])
21
{
22
// ----------------------------------------------------------------------
23
// First Part:
24
//
25
// * enable FWLite
26
// * book the histograms of interest
27
// * open the input file
28
// ----------------------------------------------------------------------
29
30
// load framework libraries
31
gSystem->Load(
"libFWCoreFWLite"
);
32
FWLiteEnabler::enable
();
33
34
// open input file (can be located on castor)
35
TFile*
inFile
= TFile::Open(
"file:jpsi.root"
);
36
37
// ----------------------------------------------------------------------
38
// Second Part:
39
//
40
// * loop the events in the input file
41
// * receive the collections of interest via fwlite::Handle
42
// * fill the histograms
43
// * after the loop close the input file
44
// ----------------------------------------------------------------------
45
46
// loop the events
47
unsigned
int
iEvent
=0;
48
fwlite::Event
event
(
inFile
);
49
for
(
event
.toBegin(); !
event
.atEnd(); ++
event
, ++
iEvent
){
50
// break loop after end of file is reached
51
// or after 1000 events have been processed
52
if
(
iEvent
==1000 )
break
;
53
54
// simple event counter
55
if
(
iEvent
>0 &&
iEvent
%1==0){
56
std::cout
<<
" processing event: "
<<
iEvent
<< std::endl;
57
}
58
59
// fwlite::Handle to to jpsi collection
60
fwlite::Handle<std::vector<pat::CompositeCandidate>
> jpsis;
61
jpsis.
getByLabel
(
event
,
"patJPsiCandidates"
);
62
63
// loop jpsi collection and fill histograms
64
for
(
unsigned
i
=0;
i
<jpsis->size(); ++
i
){
65
cout
<<
"jpsi "
<<
i
<<
", mass = "
<< jpsis->at(
i
).mass() <<
", dR = "
<< jpsis->at(
i
).userFloat(
"dR"
) << endl;
66
}
67
}
68
// close input file
69
inFile
->Close();
70
71
// that's it!
72
return
0;
73
}
mps_fire.i
i
Definition:
mps_fire.py:429
std
Definition:
JetResolutionObject.h:76
GCPpyPlots.argv
argv
Definition:
GCPpyPlots.py:3
iEvent
int iEvent
Definition:
GenABIO.cc:224
FWLiteEnabler::enable
static void enable()
enable automatic library loading
Definition:
FWLiteEnabler.cc:46
edmPickEvents.event
event
Definition:
edmPickEvents.py:280
fwlite::Handle
Definition:
Handle.h:39
Handle.h
testHGCalDigi_cfg.inFile
inFile
Definition:
testHGCalDigi_cfg.py:51
FWLiteEnabler.h
CompositeCandidate.h
fwlite::Handle::getByLabel
void getByLabel(const P &iP, const char *iModuleLabel, const char *iProductInstanceLabel=nullptr, const char *iProcessLabel=nullptr)
Definition:
Handle.h:100
dir2webdir.argc
argc
Definition:
dir2webdir.py:39
Muon.h
fwlite::Event
Definition:
Event.h:90
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
event
Definition:
event.py:1
main
int main(int argc, char *argv[])
Definition:
JPsiFWLiteAnalyzer.cc:20
Generated for CMSSW Reference Manual by
1.8.14