CMS 3D CMS Logo

ZMassHistogrammer.cc
Go to the documentation of this file.
5 #include "TH1.h"
6 
8 public:
10 
11 private:
12  void analyze(const edm::Event& event, const edm::EventSetup& setup) override;
15  TH1F *h_mZ_, *h_mZMC_;
16 };
17 
25 #include <iostream>
26 
27 using namespace std;
28 using namespace reco;
29 using namespace edm;
30 
32  : zToken_(consumes<reco::CandidateView>(pset.getParameter<InputTag>("z"))),
33  genToken_(consumes<reco::CandidateView>(pset.getParameter<InputTag>("gen"))) {
34  cout << ">>> Z Mass constructor" << endl;
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 
41  cout << ">>> Z Mass analyze" << endl;
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 
73 
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:355
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:29
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:36
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:500
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