CMS 3D CMS Logo

ZMassHistogrammer.cc
Go to the documentation of this file.
5 #include "TH1.h"
6 
8 public:
10 private:
11  virtual void analyze(const edm::Event& event, const edm::EventSetup& setup) override;
14  TH1F *h_mZ_, *h_mZMC_;
15 };
16 
24 #include <iostream>
25 
26 using namespace std;
27 using namespace reco;
28 using namespace edm;
29 
31  zToken_(consumes<reco::CandidateView>(pset.getParameter<InputTag>("z"))),
32  genToken_(consumes<reco::CandidateView>(pset.getParameter<InputTag>("gen"))) {
33  cout << ">>> Z Mass constructor" << endl;
35  h_mZ_ = fs->make<TH1F>("ZMass", "Z mass (GeV/c^{2})", 100, 0, 200);
36  h_mZMC_ = fs->make<TH1F>("ZMCMass", "Z MC mass (GeV/c^{2})", 100, 0, 200);
37 }
38 
40  cout << ">>> Z Mass analyze" << endl;
43  event.getByToken(zToken_, z);
44  event.getByToken(genToken_, gen);
45  for(unsigned int i = 0; i < z->size(); ++i) {
46  const Candidate &zCand = (*z)[i];
47  h_mZ_->Fill(zCand.mass());
48  }
49  for(unsigned int i = 0; i < gen->size(); ++i) {
50  const Candidate &genCand = (*gen)[i];
51  if((genCand.pdgId() == 23) && (genCand.status() == 2)) //this is an intermediate Z0
52  cout << ">>> intermediate Z0 found, with " << genCand.numberOfDaughters()
53  << " daughters" << endl;
54  if((genCand.pdgId() == 23)&&(genCand.status() == 3)) { //this is a Z0
55  cout << ">>> Z0 found, with " << genCand.numberOfDaughters()
56  << " daughters" << endl;
57  h_mZMC_->Fill(genCand.mass());
58  if(genCand.numberOfDaughters() == 3) {//Z0 decays in mu+ mu-, the 3rd daughter is the same Z0
59  const Candidate * dauGen0 = genCand.daughter(0);
60  const Candidate * dauGen1 = genCand.daughter(1);
61  const Candidate * dauGen2 = genCand.daughter(2);
62  cout << ">>> daughter MC 0 PDG Id " << dauGen0->pdgId()
63  << ", status " << dauGen0->status()
64  << ", charge " << dauGen0->charge()
65  << endl;
66  cout << ">>> daughter MC 1 PDG Id " << dauGen1->pdgId()
67  << ", status " << dauGen1->status()
68  << ", charge " << dauGen1->charge()
69  << endl;
70  cout << ">>> daughter MC 2 PDG Id " << dauGen2->pdgId()
71  << ", status " << dauGen2->status()
72  << ", charge " << dauGen2->charge() << endl;
73  }
74  }
75  }
76 }
77 
79 
81 
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
size_type size() const
Definition: OwnVector.h:264
virtual void analyze(const edm::Event &event, const edm::EventSetup &setup) override
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
edm::EDGetTokenT< reco::CandidateView > genToken_
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
virtual int status() const =0
status word
virtual int pdgId() const =0
PDG identifier.
edm::EDGetTokenT< reco::CandidateView > zToken_
def gen(fragment, howMuch)
Production test section ####.
ZMassHistogrammer(const edm::ParameterSet &pset)
virtual double mass() const =0
mass
virtual int charge() const =0
electric charge
fixed size matrix
HLT enums.
virtual size_type numberOfDaughters() const =0
number of daughters
Definition: event.py:1