CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ZMCHistogrammer.cc
Go to the documentation of this file.
3 #include "TH1.h"
4 
6 public:
8 private:
9  virtual void analyze(const edm::Event& event, const edm::EventSetup& setup) override;
15  TH1F *h_invmMuMu_;
18  //TH1F *h_mZ2vs3MC_, *h_ptZ2vs3MC_, *h_phiZ2vs3MC_, *h_thetaZ2vs3MC_, *h_etaZ2vs3MC_, *h_rapidityZ2vs3MC_;
22  bool isMCatNLO_;
23 };
24 
31 
32 #include "HepMC/WeightContainer.h"
33 #include "HepMC/GenEvent.h"
39 #include <cmath>
40 #include <iostream>
41 
42 using namespace std;
43 using namespace reco;
44 using namespace edm;
45 
47  z_(pset.getParameter<InputTag>("z")),
48  gen_(pset.getParameter<InputTag>("gen")),
49  match_(pset.getParameter<InputTag>("match")),
50  hepMCProductTag_(pset.getParameter<InputTag>("hepMCProductTag")),
51  nbinsMass_(pset.getUntrackedParameter<unsigned int>("nbinsMass")),
52  nbinsPt_(pset.getUntrackedParameter<unsigned int>("nbinsPt")),
53  nbinsAng_(pset.getUntrackedParameter<unsigned int>("nbinsAng")),
54  nbinsMassRes_(pset.getUntrackedParameter<unsigned int>("nbinsMassRes")),
55  massMax_(pset.getUntrackedParameter<double>("massMax")),
56  ptMax_(pset.getUntrackedParameter<double>("ptMax")),
57  angMax_(pset.getUntrackedParameter<double>("angMax")),
58  massResMax_(pset.getUntrackedParameter<double>("massResMax")),
59  isMCatNLO_(pset.getUntrackedParameter<bool>("isMCatNLO")) {
60  cout << ">>> Z Histogrammer constructor" << endl;
62  TFileDirectory ZHisto = fs->mkdir( "ZRecoHisto" );
63  TFileDirectory ZMCHisto = fs->mkdir( "ZMCHisto" );
64  TFileDirectory ZResHisto = fs->mkdir( "ZResolutionHisto" );
65  //TFileDirectory Z2vs3MCHisto = fs->mkdir( "Z2vs3MCHisto" );
66  h_nZ_ = ZHisto.make<TH1F>("ZNumber", "number of Z particles", 11, -0.5, 10.5);
67  h_mZ_ = ZHisto.make<TH1F>("ZMass", "Z mass (GeV/c^{2})", nbinsMass_, 0, massMax_);
68  h_ptZ_ = ZHisto.make<TH1F>("ZPt", "Z p_{t} (GeV/c)", nbinsPt_, 0, ptMax_);
69  h_phiZ_ = ZHisto.make<TH1F>("ZPhi", "Z #phi", nbinsAng_, -angMax_, angMax_);
70  h_thetaZ_ = ZHisto.make<TH1F>("Ztheta", "Z #theta", nbinsAng_, 0, angMax_);
71  h_etaZ_ = ZHisto.make<TH1F>("ZEta", "Z #eta", nbinsAng_, -angMax_, angMax_);
72  h_rapidityZ_ = ZHisto.make<TH1F>("ZRapidity", "Z rapidity", nbinsAng_, -angMax_, angMax_);
73  h_invmMuMu_ = ZHisto.make<TH1F>("MuMuMass", "#mu #mu invariant mass",nbinsMass_, 0, massMax_);
74  h_weight_histo = ZHisto.make<TH1F>("weight_histo","weight_histo",20,-10,10);
75  h_nZMC_ = ZMCHisto.make<TH1F>("ZMCNumber", "number of Z MC particles", 11, -0.5, 10.5);
76  h_mZMC_ = ZMCHisto.make<TH1F>("ZMCMass", "Z MC mass (GeV/c^{2})", nbinsMass_, 0, massMax_);
77  h_ptZMC_ = ZMCHisto.make<TH1F>("ZMCPt", "Z MC p_{t} (GeV/c)", nbinsPt_, 0, ptMax_);
78  h_phiZMC_ = ZMCHisto.make<TH1F>("ZMCPhi", "Z MC #phi", nbinsAng_, -angMax_, angMax_);
79  h_thetaZMC_ = ZMCHisto.make<TH1F>("ZMCTheta", "Z MC #theta", nbinsAng_, 0, angMax_);
80  h_etaZMC_ = ZMCHisto.make<TH1F>("ZMCEta", "Z MC #eta", nbinsAng_, -angMax_, angMax_);
81  h_rapidityZMC_ = ZMCHisto.make<TH1F>("ZMCRapidity", "Z MC rapidity", nbinsAng_, -angMax_, angMax_);
82  h_invmMuMuMC_ = ZMCHisto.make<TH1F>("MuMuMCMass", "#mu #mu MC invariant mass",
83  nbinsMass_, 0, massMax_);
84  /*
85  h_mZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCMass", "Z MC st 2 vs st 3 mass (GeV/c^{2})",
86  nbinsMassRes_, -massResMax_, massResMax_);
87  h_ptZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCPt", "Z MC st 2 vs st 3 p_{t} (GeV/c)",
88  nbinsPt_, -ptMax_, ptMax_);
89  h_phiZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCPhi", "Z MC st 2 vs st 3 #phi",
90  nbinsAng_, -angMax_, angMax_);
91  h_thetaZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCTheta", "Z MC st 2 vs st 3 #theta",
92  nbinsAng_, -angMax_, angMax_);
93  h_etaZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCEta", "Z MC st 2 vs st 3 #eta",
94  nbinsAng_, -angMax_, angMax_);
95  h_rapidityZ2vs3MC_ = Z2vs3MCHisto.make<TH1F>("Z2vs3MCRapidity", "Z MC st 2 vs st 3 rapidity",
96  nbinsAng_, -angMax_, angMax_);
97  */
98  h_mResZ_ = ZResHisto.make<TH1F>("ZMassResolution", "Z mass Resolution (GeV/c^{2})",
100  h_ptResZ_ = ZResHisto.make<TH1F>("ZPtResolution", "Z p_{t} Resolution (GeV/c)",
101  nbinsPt_, -ptMax_, ptMax_);
102  h_phiResZ_ = ZResHisto.make<TH1F>("ZPhiResolution", "Z #phi Resolution",
104  h_thetaResZ_ = ZResHisto.make<TH1F>("ZThetaResolution", "Z #theta Resolution",
106  h_etaResZ_ = ZResHisto.make<TH1F>("ZEtaResolution", "Z #eta Resolution",
108  h_rapidityResZ_ = ZResHisto.make<TH1F>("ZRapidityResolution", "Z rapidity Resolution",
110  h_mResZMuMu_ = ZResHisto.make<TH1F>("ZToMuMuRecoMassResolution",
111  "Z Reco vs matched final state #mu #mu mass Difference (GeV/c^{2})",
113  h_mRatioZMuMu_ = ZResHisto.make<TH1F>("ZToMuMuRecoMassRatio",
114  "Z Reco vs matched final state #mu #mu mass Ratio",
115  4000, 0, 2);
116  h_mResZMuMuMC_ = ZResHisto.make<TH1F>("ZToMuMuMCMassResolution",
117  "Z vs final state #mu #mu MC mass Difference (GeV/c^{2})",
118  nbinsMassRes_/2 + 1, -2*massResMax_/nbinsMassRes_, massResMax_);
119  h_mRatioZMuMuMC_ = ZResHisto.make<TH1F>("ZToMuMuMCMassRatio",
120  "Z vs final state #mu #mu MC mass Ratio",
121  2002, 0.999, 2);
122 }
123 
125  cout << ">>> Z Histogrammer analyze" << endl;
129  event.getByLabel(z_, z);
130  event.getByLabel(gen_, gen);
131  event.getByLabel(match_, match);
132  h_nZ_->Fill(z->size());
133 
134 // get HepMC::GenEvent ...
135  Handle<HepMCProduct> evt_h;
136  event.getByLabel(hepMCProductTag_, evt_h);
137  HepMC::GenEvent* evt = new HepMC::GenEvent(*(evt_h->GetEvent()));
138 
139 
140  // get weight and fill it to histogram
141  HepMC::WeightContainer weights = evt->weights();
142  double weight = weights.front();
143  if(!weight) weight=1.;
144  h_weight_histo->Fill(weight);
145 
146  if(isMCatNLO_) {
147  weight > 0 ? weight=1. : weight=-1.;
148  }
149 
150  for(unsigned int i = 0; i < z->size(); ++i) {
151  const Candidate &zCand = (*z)[i];
152  h_mZ_->Fill(zCand.mass(),weight );
153  h_ptZ_->Fill(zCand.pt(),weight);
154  h_phiZ_->Fill(zCand.phi(),weight);
155  h_thetaZ_->Fill(zCand.theta(),weight);
156  h_etaZ_->Fill(zCand.eta(),weight);
157  h_rapidityZ_->Fill(zCand.rapidity(),weight);
158  CandidateBaseRef zCandRef = z->refAt(i);
159 
160  GenParticleRef zMCMatch = (*match)[i];
161  if(zMCMatch.isNonnull() && zMCMatch->pdgId()==23) {
162  h_mResZ_->Fill(zCandRef->mass() - zMCMatch->mass());
163  h_ptResZ_->Fill(zCandRef->pt() - zMCMatch->pt());
164  h_phiResZ_->Fill(zCandRef->phi() - zMCMatch->phi());
165  h_thetaResZ_->Fill(zCandRef->theta() - zMCMatch->theta());
166  h_etaResZ_->Fill(zCandRef->eta() - zMCMatch->eta());
167  h_rapidityResZ_->Fill(zCandRef->rapidity() - zMCMatch->rapidity());
168  const Candidate * dau0 = zMCMatch->daughter(0);
169  const Candidate * dau1 = zMCMatch->daughter(1);
170  for(unsigned int i0 = 0; i0 < dau0->numberOfDaughters(); ++i0) {
171  const Candidate * ddau0 = dau0->daughter(i0);
172  if(abs(ddau0->pdgId())==13 && ddau0->status()==1) {
173  dau0 = ddau0; break;
174  }
175  }
176  for(unsigned int i1 = 0; i1 < dau1->numberOfDaughters(); ++i1) {
177  const Candidate * ddau1 = dau1->daughter(i1);
178  if(abs(ddau1->pdgId())==13 && ddau1->status()==1) {
179  dau1 = ddau1; break;
180  }
181  }
182  assert(abs(dau0->pdgId())==13 && dau0->status()==1);
183  assert(abs(dau1->pdgId())==13 && dau1->status()==1);
184  double invMass = (dau0->p4()+dau1->p4()).mass();
185  h_invmMuMu_->Fill(invMass,weight);
186  h_mResZMuMu_->Fill(zCand.mass() - invMass);
187  h_mRatioZMuMu_->Fill(zCand.mass()/invMass);
188  }
189  }
190  h_nZMC_->Fill(gen->size());
191  for(unsigned int i = 0; i < gen->size(); ++i) {
192  const Candidate &genCand = (*gen)[i];
193  if((genCand.pdgId() == 23) && (genCand.status() == 2)) //this is an intermediate Z0
194  cout << ">>> intermediate Z0 found, with " << genCand.numberOfDaughters()
195  << " daughters" << endl;
196  if((genCand.pdgId() == 23)&&(genCand.status() == 3)) { //this is a Z0
197  cout << ">>> Z0 found, with " << genCand.numberOfDaughters()
198  << " daughters" << endl;
199  h_mZMC_->Fill(genCand.mass(),weight);
200  h_ptZMC_->Fill(genCand.pt(),weight);
201  h_phiZMC_->Fill(genCand.phi(),weight);
202  h_thetaZMC_->Fill(genCand.theta(),weight);
203  h_etaZMC_->Fill(genCand.eta(),weight);
204  h_rapidityZMC_->Fill(genCand.rapidity(),weight);
205  Particle::LorentzVector pZ(0, 0, 0, 0);
206  int nMu = 0;
207  for(unsigned int j = 0; j < genCand.numberOfDaughters(); ++j) {
208  const Candidate *dauGen = genCand.daughter(j);
209  /*
210  if((dauGen->pdgId() == 23) && (dauGen->status() == 2)) {
211  h_mZ2vs3MC_->Fill(genCand.mass() - dauGen->mass());
212  h_ptZ2vs3MC_->Fill(genCand.pt() - dauGen->pt());
213  h_phiZ2vs3MC_->Fill(genCand.phi() - dauGen->phi());
214  h_thetaZ2vs3MC_->Fill(genCand.theta() - dauGen->theta());
215  h_etaZ2vs3MC_->Fill(genCand.eta() - dauGen->eta());
216  h_rapidityZ2vs3MC_->Fill(genCand.rapidity() - dauGen->rapidity());
217  }
218  */
219  if((abs(dauGen->pdgId()) == 13) && (dauGen->numberOfDaughters() != 0)) {
220  //we are looking for photons of final state radiation
221  cout << ">>> The muon " << j
222  << " has " << dauGen->numberOfDaughters() << " daughters" <<endl;
223  for(unsigned int k = 0; k < dauGen->numberOfDaughters(); ++k) {
224  const Candidate * dauMuGen = dauGen->daughter(k);
225  cout << ">>> Mu " << j
226  << " daughter MC " << k
227  << " PDG Id " << dauMuGen->pdgId()
228  << ", status " << dauMuGen->status()
229  << ", charge " << dauMuGen->charge()
230  << endl;
231  if(abs(dauMuGen->pdgId()) == 13 && dauMuGen->status() ==1) {
232  pZ += dauMuGen->p4();
233  nMu ++;
234  }
235  }
236  }
237  }
238  assert(nMu == 2);
239  double mZ = pZ.mass();
240  h_invmMuMuMC_->Fill(mZ,weight);
241  h_mResZMuMuMC_->Fill(genCand.mass() - mZ);
242  h_mRatioZMuMuMC_->Fill(genCand.mass()/mZ);
243  }
244  }
245 }
246 
248 
250 
int i
Definition: DBlmapReader.cc:9
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
virtual void analyze(const edm::Event &event, const edm::EventSetup &setup) override
virtual float mass() const =0
mass
edm::InputTag z_
virtual float eta() const =0
momentum pseudorapidity
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
ZMCHistogrammer(const edm::ParameterSet &pset)
virtual int status() const =0
status word
virtual float phi() const =0
momentum azimuthal angle
virtual double rapidity() const =0
rapidity
float float float z
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:250
virtual size_type numberOfDaughters() const =0
number of daughters
unsigned int nbinsAng_
virtual float pt() const =0
transverse momentum
edm::InputTag hepMCProductTag_
virtual double theta() const =0
momentum polar angle
edm::InputTag gen_
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
edm::InputTag match_
virtual int charge() const =0
electric charge
T * make(const Args &...args) const
make new ROOT object
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
int k[5][pyjets_maxn]
virtual int pdgId() const =0
PDG identifier.
unsigned int nbinsMass_
TFileDirectory mkdir(const std::string &dir, const std::string &descr="")
create a new subdirectory
Definition: TFileService.h:69
tuple cout
Definition: gather_cfg.py:121
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
int weight
Definition: histoStyle.py:50
unsigned int nbinsPt_
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Particle.h:27
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
unsigned int nbinsMassRes_