Main Page
Namespaces
Namespace List
Namespace Members
All
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Functions
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Variables
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Typedefs
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Enumerations
a
b
c
d
e
f
g
h
i
j
k
l
m
o
p
q
r
s
t
u
v
w
z
Enumerator
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Classes
Class List
Class Index
Class Hierarchy
Class Members
All
:
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Functions
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
~
Variables
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Typedefs
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Enumerations
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
Enumerator
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
y
z
Properties
_
a
d
e
f
l
m
o
p
s
t
u
v
Related Functions
:
_
a
b
c
d
e
f
g
h
i
j
k
l
m
n
o
p
q
r
s
t
u
v
w
x
Package Documentation
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
ElectroWeakAnalysis
ZMuMu
plugins
ZMassHistogrammer.cc
Go to the documentation of this file.
1
#include "
FWCore/Framework/interface/EDAnalyzer.h
"
2
#include "
FWCore/Utilities/interface/InputTag.h
"
3
#include "
DataFormats/Candidate/interface/Candidate.h
"
4
#include "
DataFormats/Candidate/interface/CandidateFwd.h
"
5
#include "TH1.h"
6
7
class
ZMassHistogrammer
:
public
edm::EDAnalyzer
{
8
public
:
9
ZMassHistogrammer
(
const
edm::ParameterSet
&
pset
);
10
11
private
:
12
void
analyze
(
const
edm::Event
&
event
,
const
edm::EventSetup
&
setup
)
override
;
13
edm::EDGetTokenT<reco::CandidateView>
zToken_
;
14
edm::EDGetTokenT<reco::CandidateView>
genToken_
;
15
TH1F *
h_mZ_
, *
h_mZMC_
;
16
};
17
18
#include "
DataFormats/HepMCCandidate/interface/GenParticle.h
"
19
#include "
DataFormats/HepMCCandidate/interface/GenParticleFwd.h
"
20
#include "
FWCore/Framework/interface/Event.h
"
21
#include "
DataFormats/Common/interface/Handle.h
"
22
#include "
FWCore/ParameterSet/interface/ParameterSet.h
"
23
#include "
FWCore/ServiceRegistry/interface/Service.h
"
24
#include "
CommonTools/UtilAlgos/interface/TFileService.h
"
25
#include <iostream>
26
27
using namespace
std
;
28
using namespace
reco
;
29
using namespace
edm
;
30
31
ZMassHistogrammer::ZMassHistogrammer
(
const
ParameterSet
&
pset
)
32
: zToken_(consumes<
reco
::
CandidateView
>(
pset
.getParameter<
InputTag
>(
"z"
))),
33
genToken_(consumes<
reco
::
CandidateView
>(
pset
.getParameter<
InputTag
>(
"gen"
))) {
34
cout
<<
">>> Z Mass constructor"
<< endl;
35
Service<TFileService>
fs;
36
h_mZ_
= fs->
make
<TH1F>(
"ZMass"
,
"Z mass (GeV/c^{2})"
, 100, 0, 200);
37
h_mZMC_
= fs->
make
<TH1F>(
"ZMCMass"
,
"Z MC mass (GeV/c^{2})"
, 100, 0, 200);
38
}
39
40
void
ZMassHistogrammer::analyze
(
const
edm::Event
&
event
,
const
edm::EventSetup
&
setup
) {
41
cout
<<
">>> Z Mass analyze"
<< endl;
42
Handle<CandidateCollection>
z
;
43
Handle<CandidateCollection>
gen
;
44
event
.getByToken(
zToken_
,
z
);
45
event
.getByToken(
genToken_
,
gen
);
46
for
(
unsigned
int
i
= 0;
i
<
z
->size(); ++
i
) {
47
const
Candidate
& zCand = (*z)[
i
];
48
h_mZ_
->Fill(zCand.
mass
());
49
}
50
for
(
unsigned
int
i
= 0;
i
<
gen
->size(); ++
i
) {
51
const
Candidate
& genCand = (*gen)[
i
];
52
if
((genCand.
pdgId
() == 23) && (genCand.
status
() == 2))
//this is an intermediate Z0
53
cout
<<
">>> intermediate Z0 found, with "
<< genCand.
numberOfDaughters
() <<
" daughters"
<< endl;
54
if
((genCand.
pdgId
() == 23) && (genCand.
status
() == 3)) {
//this is a Z0
55
cout
<<
">>> Z0 found, with "
<< genCand.
numberOfDaughters
() <<
" daughters"
<< endl;
56
h_mZMC_
->Fill(genCand.
mass
());
57
if
(genCand.
numberOfDaughters
() == 3) {
//Z0 decays in mu+ mu-, the 3rd daughter is the same Z0
58
const
Candidate
* dauGen0 = genCand.
daughter
(0);
59
const
Candidate
* dauGen1 = genCand.
daughter
(1);
60
const
Candidate
* dauGen2 = genCand.
daughter
(2);
61
cout
<<
">>> daughter MC 0 PDG Id "
<< dauGen0->
pdgId
() <<
", status "
<< dauGen0->
status
() <<
", charge "
62
<< dauGen0->
charge
() << endl;
63
cout
<<
">>> daughter MC 1 PDG Id "
<< dauGen1->
pdgId
() <<
", status "
<< dauGen1->
status
() <<
", charge "
64
<< dauGen1->
charge
() << endl;
65
cout
<<
">>> daughter MC 2 PDG Id "
<< dauGen2->
pdgId
() <<
", status "
<< dauGen2->
status
() <<
", charge "
66
<< dauGen2->
charge
() << endl;
67
}
68
}
69
}
70
}
71
72
#include "
FWCore/Framework/interface/MakerMacros.h
"
73
74
DEFINE_FWK_MODULE
(
ZMassHistogrammer
);
reco::Candidate::daughter
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode)
Handle.h
mps_fire.i
i
Definition:
mps_fire.py:428
ZMassHistogrammer::analyze
void analyze(const edm::Event &event, const edm::EventSetup &setup) override
Definition:
ZMassHistogrammer.cc:40
reco::Candidate::mass
virtual double mass() const =0
mass
edm::EDGetTokenT
Definition:
EDGetToken.h:33
edm
HLT enums.
Definition:
AlignableModifier.h:19
gather_cfg.cout
cout
Definition:
gather_cfg.py:144
ZMassHistogrammer::genToken_
edm::EDGetTokenT< reco::CandidateView > genToken_
Definition:
ZMassHistogrammer.cc:14
reco::Candidate::status
virtual int status() const =0
status word
EDAnalyzer.h
reco
fixed size matrix
Definition:
AlignmentAlgorithmBase.h:45
ZMassHistogrammer
Definition:
ZMassHistogrammer.cc:7
edm::Handle
Definition:
AssociativeIterator.h:50
singleTopDQM_cfi.setup
setup
Definition:
singleTopDQM_cfi.py:37
edm::EDAnalyzer
Definition:
EDAnalyzer.h:28
GenParticle.h
CandidateFwd.h
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition:
MakerMacros.h:16
reco::Candidate::charge
virtual int charge() const =0
electric charge
Service.h
GenParticleFwd.h
DDAxes::z
reco::Candidate::numberOfDaughters
virtual size_type numberOfDaughters() const =0
number of daughters
gen
Definition:
PythiaDecays.h:13
TFileService.h
edm::View
Definition:
CaloClusterFwd.h:14
ZMassHistogrammer::zToken_
edm::EDGetTokenT< reco::CandidateView > zToken_
Definition:
ZMassHistogrammer.cc:13
ZMassHistogrammer::h_mZMC_
TH1F * h_mZMC_
Definition:
ZMassHistogrammer.cc:15
edm::ParameterSet
Definition:
ParameterSet.h:47
Event.h
edm::Service< TFileService >
edm::EventSetup
Definition:
EventSetup.h:57
reco::Candidate::pdgId
virtual int pdgId() const =0
PDG identifier.
InputTag.h
reco::Candidate
Definition:
Candidate.h:27
std
Definition:
JetResolutionObject.h:76
ZMassHistogrammer::ZMassHistogrammer
ZMassHistogrammer(const edm::ParameterSet &pset)
Definition:
ZMassHistogrammer.cc:31
relval_steps.gen
def gen(fragment, howMuch)
Production test section ####.
Definition:
relval_steps.py:509
Candidate.h
ParameterSet.h
event
Definition:
event.py:1
edm::Event
Definition:
Event.h:73
edm::InputTag
Definition:
InputTag.h:15
TFileService::make
T * make(const Args &... args) const
make new ROOT object
Definition:
TFileService.h:64
muonDTDigis_cfi.pset
pset
Definition:
muonDTDigis_cfi.py:27
ZMassHistogrammer::h_mZ_
TH1F * h_mZ_
Definition:
ZMassHistogrammer.cc:15
Generated for CMSSW Reference Manual by
1.8.16