ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/HbbAnalysis/src/JetFlavour.cc
Revision: 1.1
Committed: Thu May 28 15:12:45 2009 UTC (15 years, 11 months ago) by amagnan
Content type: text/plain
Branch: MAIN
Log Message:
add JetFlavour class + selectors

File Contents

# User Rev Content
1 amagnan 1.1 #include <iostream>
2     #include <fstream>
3    
4     #include "Math/VectorUtil.h"
5    
6     #include "UserCode/HbbAnalysis/interface/JetFlavour.hh"
7    
8     namespace HbbAnalysis {
9    
10     void JetFlavour::Initialise(TFileDirectory & aDir, unsigned short aDebug, unsigned int aFlavour){
11    
12     debug_ = aDebug;
13     flavour_ = aFlavour;
14    
15     partons_.clear();
16     flav_.clear();
17     eid_.clear();
18     muid_.clear();
19    
20     p_nQuarks = aDir.make<TH1F>("p_nQuarks",";n_{q};N_{entries}",10,0,10);
21     p_nAntiQuarks = aDir.make<TH1F>("p_nAntiQuarks",";n_{#bar{q}};N_{entries}",10,0,10);
22     p_partonFlavour = aDir.make<TH1F>("p_partonFlavour",";pdgId of initial parton;N_{entries}",28,-6,22);
23     p_parentOfStableFlavour = aDir.make<TH1F>("p_parentOfStableFlavour",";pdgId of parent of stable particle;N_{entries}",12000,-6000,6000);
24     p_stableFlavour = aDir.make<TH1F>("p_stableFlavour",";pdgId of stable particle;N_{entries}",2000,-1000,1000);
25     p_numberOfParentsOfStable = aDir.make<TH1F>("p_numberOfParentsOfStable",";number of parents of stable particle;N_{entries}",50,0,50);
26    
27    
28     nTot_ = 0;
29     nLight_ = 0;
30     nGluon_ = 0;
31     nCharm_ = 0;
32     nBottom_ = 0;
33     nTotal_ = 0;
34     nHeavy_ = 0;
35     nSemiHeavy_ = 0;
36    
37     nZEvts_ = 0;
38     nWEvts_ = 0;
39    
40     if (debug_) std::cout << " --- JetFlavour::Initialise success ! " << std::endl;
41    
42     }
43    
44     unsigned int JetFlavour::fillPartons(const edm::Handle<reco::GenParticleCollection> & aGenParticles) {//genMatching
45    
46     nTot_++;
47    
48     partons_.clear();
49     flav_.clear();
50     eid_.clear();
51     muid_.clear();
52     partons_.reserve(10);
53     eid_.reserve(10);
54     muid_.reserve(10);
55     flav_.reserve(10);
56    
57     unsigned int n_quarks = 0;
58     unsigned int n_antiquarks = 0;
59    
60     bool isZ = false;
61     bool isW = false;
62    
63     unsigned int lNMC = aGenParticles->size();
64    
65     if (debug_) std::cout << "* nMC = " << lNMC << std::endl;
66    
67     //loop on MC part to find b and bbar initial particles.
68     for (unsigned int imc(0); imc < lNMC; imc++){//loop on MC particles
69     const reco::GenParticle & lPart = (*aGenParticles).at(imc);
70     //break loop if W events
71     //for QCD, it should never fullfill these criteria anyway....
72     if ( lPart.status() == 3 && lPart.numberOfMothers() > 1 &&
73     ( abs(lPart.pdgId()) == 24 ||
74     abs(lPart.pdgId()) == 12 ||
75     abs(lPart.pdgId()) == 14 ||
76     abs(lPart.pdgId()) == 16 ) ) {
77     isW = true;
78     nWEvts_++;
79     if (!isZ) break;
80     }
81     //record when Z event
82     if ( lPart.status() == 3 && lPart.numberOfMothers() > 1 &&
83     ( abs(lPart.pdgId()) == 23 ||
84     lPart.pdgId() == 11 ||
85     lPart.pdgId() == 13 ||
86     lPart.pdgId() == 15 ) ) {
87     isZ = true;
88     nZEvts_++;
89     }
90     //looking for the final-states partons, which will have 2 parents being the 2 initial states partons.
91     if (lPart.status() == 3 && lPart.numberOfMothers() > 1 &&
92     ( ( flavour_ == 45 && (abs(lPart.pdgId()) == 4 || abs(lPart.pdgId()) == 5) ) ||
93     ( flavour_ > 3 && flavour_ < 22 && abs(lPart.pdgId()) == flavour_) ||
94     ( flavour_ == 3 && ( (abs(lPart.pdgId()) <= flavour_ ||
95     lPart.pdgId() == 21)
96     )
97     )
98     )
99     ) {
100     if (debug_ > 0) {
101     std::cout << "** initial parton #" << imc ;
102     printParticle(lPart);
103     }
104     partons_.push_back(imc);
105     if (lPart.pdgId() > 0) n_quarks++;
106     else n_antiquarks++;
107     p_partonFlavour->Fill(lPart.pdgId());
108     }
109     else if ((debug_ > 0 && imc < 15) || debug_ > 1) {
110     std::cout << "** not selected: #" << imc;
111     printParticle(lPart);
112     }
113    
114     }//loop on MC particles
115    
116     //if (isW == isZ) std::cout << "*** WARNING *** isW = " << isW << ", isZ = " << isZ << std::endl;
117     if ((flavour_ > 3 && flavour_ != 21) && !isZ) return 0;
118    
119     if (debug_) {
120     std::cout << "nQuarks = " << n_quarks << ", n_antiQuarks = " << n_antiquarks << std::endl;
121     std::cout << "nPartons = " << partons_.size() << " with indices " << std::endl;
122     for (unsigned int ip(0); ip<partons_.size();ip++){
123     std::cout << "#" << partons_.at(ip) << ": " ;
124     printParticle(parton(aGenParticles,ip));
125     }
126     }
127    
128     const unsigned int nPartons = partons_.size();
129    
130     if (nPartons == 0) {
131     nHeavy_++;
132     return 0;
133     }
134     if (nPartons == 1) nSemiHeavy_++;
135    
136     nTotal_ += nPartons;
137     assert (nPartons == n_quarks+n_antiquarks);
138    
139     if (flavour_ > 3) {
140     assert(n_quarks >= 1 && n_antiquarks >= 1);
141     assert (nPartons == 2);
142     }
143    
144     p_nQuarks->Fill(n_quarks);
145     p_nAntiQuarks->Fill(n_antiquarks);
146    
147    
148     //loop on MC to find partons daughters
149     unsigned int *nDaugh = new unsigned int[nPartons];
150     //vectors of indices of parents/daughters
151     //indices with respect to initial genParticle collection (always read in the same order)
152     std::vector<unsigned int> *parents = new std::vector<unsigned int>[nPartons];
153     std::vector<unsigned int> *daughters = new std::vector<unsigned int>[nPartons];
154    
155     unsigned int *flav = new unsigned int[nPartons];
156    
157     //indices of stable electron/muon, if any
158     unsigned int *eid = new unsigned int[nPartons];
159     unsigned int *muid = new unsigned int[nPartons];
160    
161    
162     for (unsigned int ip(0); ip<nPartons; ip++){//loop on partons
163     nDaugh[ip] = 0;
164     flav[ip] = 1;
165     if (flavour_ == 3 || flavour_ == 21){
166     if (parton(aGenParticles,ip).pdgId() == 21) flav[ip] = 3;
167     else flav[ip] = 2;
168     }
169     eid[ip] = 0;
170     muid[ip] = 0;
171     daughters[ip].clear();
172     daughters[ip].reserve(250);
173     parents[ip].clear();
174     parents[ip].reserve(250);
175     parents[ip].push_back(partons_.at(ip));
176    
177     for (unsigned int pmc = 0; pmc < parents[ip].size(); pmc++){//loop "interactively" on parents
178     unsigned int lParIndex = parents[ip].at(pmc);
179     if (debug_ > 2) std::cout << "** Looking for parent : #" << lParIndex << " (id " << (*aGenParticles).at(lParIndex).pdgId() << "), vector size = " << parents[ip].size() << std::endl;
180    
181     //counter to check a match is found
182     unsigned int lNNotFound = 0;
183     unsigned int lNTot = 0;
184    
185     for (unsigned int imc(0); imc < lNMC; imc++){//loop on parts
186     const reco::GenParticle & pgen = (*aGenParticles).at(imc);
187    
188     for (unsigned int im(0); im<pgen.numberOfMothers(); im++){//loop on parents
189     //recursive loop to find all parents
190     if ( compareParticles(pgen.mother(im),(*aGenParticles).at(lParIndex)) ){//if par=parton
191     lNTot++;
192     if (pgen.status() != 1) {
193     nDaugh[ip]++;
194     if (debug_ > 2) {
195     std::cout << "*** b daughter #" << nDaugh[ip] << " : " << imc ;
196     printParticle(pgen);
197     }
198     //to recursively account for every possible parent until the last stable particles...
199     if (std::find(parents[ip].begin(),parents[ip].end(),imc)==parents[ip].end()) parents[ip].push_back(imc);
200     }
201     else {
202     nDaugh[ip]++;
203     if (debug_ > 2) {
204     std::cout << "*** b stable daughter #" << nDaugh[ip] << " : " << imc;
205     printParticle(pgen);
206     }
207     if (std::find(daughters[ip].begin(),daughters[ip].end(),imc)==daughters[ip].end()) {
208     daughters[ip].push_back(imc);
209     p_stableFlavour->Fill(pgen.pdgId());
210     p_numberOfParentsOfStable->Fill(pgen.numberOfMothers());
211     p_parentOfStableFlavour->Fill(pgen.mother(0)->pdgId());
212     if (fabs(pgen.pdgId()) == 11) {
213     assert(pgen.numberOfMothers() == 1);
214     }
215     else if (fabs(pgen.pdgId()) == 13) {
216     assert(pgen.numberOfMothers() == 1);
217     }
218     }
219     }
220     }//if par=parton
221     else {
222     lNNotFound++;
223     lNTot++;
224     }
225     }//loop on parents
226    
227     }//loop on parts
228     if (lNNotFound == lNTot ) {
229     std::cout << " --- WARNING: no particle found for parton #" << lParIndex ;
230     printParticle((*aGenParticles).at(lParIndex));
231     }
232    
233     }//loop on parents
234    
235     //loop on daughters to find e/mu coming from B and D mesons only
236     //and associate flav=2 if b,c decay into e, 3 if mu
237    
238     if (flavour_ > 3 && flavour_ != 21){//for b and c jets
239     if (abs(parton(aGenParticles,ip).pdgId()) == 4) nCharm_++;
240     else nBottom_++;
241     for (unsigned int dmc = 0; dmc < daughters[ip].size(); dmc++){//loop on daughters
242     const reco::GenParticle & pgen = (*aGenParticles).at(daughters[ip].at(dmc));
243     if (fabs(pgen.pdgId()) == 11 && pgen.numberOfMothers() > 1) std::cout << "***WARNING*** b Electron has more than one parent !! : " << pgen.numberOfMothers() << std::endl;
244     if (fabs(pgen.pdgId()) == 13 && pgen.numberOfMothers() > 1) std::cout << "***WARNING*** b Muon has more than one parent !! : " << pgen.numberOfMothers() << std::endl;
245     if (fabs(pgen.pdgId()) == 11 && pgen.numberOfMothers() == 1 && ((fabs(pgen.mother(0)->pdgId()) > 400 && fabs(pgen.mother(0)->pdgId()) < 600) || fabs(pgen.mother(0)->pdgId()) > 4000))
246     {
247     flav[ip] = 2;
248     eid[ip] = daughters[ip].at(dmc);
249    
250     break;
251     }
252     if (fabs(pgen.pdgId()) == 13 && pgen.numberOfMothers() == 1 && ((fabs(pgen.mother(0)->pdgId()) > 400 && fabs(pgen.mother(0)->pdgId()) < 600) || fabs(pgen.mother(0)->pdgId()) > 4000))
253     {
254     flav[ip] = 3;
255     muid[ip] = daughters[ip].at(dmc);
256     break;
257     }
258     }//loop on daughters
259    
260     flav_.push_back(flav[ip]);
261     eid_.push_back(eid[ip]);
262     muid_.push_back(muid[ip]);
263    
264     }//for b and c jets
265    
266    
267     if (parents[ip].size() > 250) std::cout << "WARNING ! Initial default size of parents not sufficient = " << parents[ip].size() << std::endl;
268     if (daughters[ip].size() > 250) std::cout << "WARNING ! Initial default size of daughters not sufficient = " << daughters[ip].size() << std::endl;
269    
270    
271     }//loop on partons
272    
273     //reject QCD if B or C in the decay chain
274     if (flavour_ == 3 || flavour_ == 21){//for light jets
275     int pp = 0;
276     for (std::vector<unsigned int>::iterator lPart = partons_.begin(); lPart!= partons_.end(); lPart++){//loop on partons
277     bool isHeavy = false;
278     for (unsigned int dmc = 0; dmc < daughters[pp].size(); dmc++){//loop on daughters
279     const reco::GenParticle & pgen = (*aGenParticles).at(daughters[pp].at(dmc));
280    
281     if (abs(pgen.pdgId()) == 4 || abs(pgen.pdgId()) == 5) {
282     //if (debug_ > 2)
283     std::cout << "*** Found daughter with pdgId = " << pgen.pdgId() << " in decay chain of parton #" << pp << ", deleting from partons vector and going to next parton." << std::endl;
284     if (abs(pgen.pdgId()) == 4) nCharm_++;
285     else nBottom_++;
286     isHeavy = true;
287     partons_.erase(lPart);
288     break;
289     }
290     }//loop on daughters
291     if (!isHeavy) {
292     flav_.push_back(flav[pp]);
293     //to ensure all vectors have the same size...
294     eid_.push_back(0);
295     muid_.push_back(0);
296    
297     if (flav[pp] == 3) nGluon_++;
298     else nLight_++;
299     }
300     pp++;
301    
302     }//loop on partons
303    
304     }//for light jets
305    
306     delete[] nDaugh;
307     delete[] parents;
308     delete[] daughters;
309     delete[] flav;
310     delete[] eid;
311     delete[] muid;
312    
313     return nPartons;
314    
315    
316     }//fillPartons
317    
318    
319     std::pair<int,unsigned int> JetFlavour::partonMatchingGenJet(const pat::Jet & aJet, const edm::Handle<reco::GenParticleCollection> & aGenParticles, double aDR){//parton match with gen jet
320    
321     int aIndex = -1;
322     double deltaRMin = 1000;
323    
324     if (!aJet.genJet()) return std::pair<int,unsigned int>(aIndex,0);
325    
326     for (unsigned int lPart(0); lPart<nPartons(); lPart++){
327     double lDR = ROOT::Math::VectorUtil::DeltaR(parton(aGenParticles,lPart).p4(),aJet.genJet()->p4());
328     if (lDR < deltaRMin) {
329     deltaRMin = lDR;
330     aIndex = lPart;
331     }
332     }
333    
334     if (deltaRMin < aDR) return std::pair<int,unsigned int>(aIndex,partonFlavour(aIndex));
335     else return std::pair<int,unsigned int>(aIndex,0);
336    
337     }//parton match with gen jet
338    
339     unsigned int JetFlavour::leptonMatchingGenJet(const pat::Jet & aJet, const edm::Handle<reco::GenParticleCollection> & aGenParticles, double aDR){//lepton match with gen jet
340    
341     int aIndex = -1;
342     double deltaRMin = 1000;
343    
344     if (!aJet.genJet()) return 0;
345    
346     if (flavour_ > 3 && flavour_ != 21) {
347     for (unsigned int lPart(0); lPart<nPartons(); lPart++){
348     if (eid_.at(lPart)){
349     double lDR = ROOT::Math::VectorUtil::DeltaR(stableElectron(aGenParticles,lPart).p4(),aJet.genJet()->p4());
350     if (lDR < deltaRMin) {
351     deltaRMin = lDR;
352     aIndex = lPart;
353     }
354     }
355     else if (muid_.at(lPart)){
356     double lDR = ROOT::Math::VectorUtil::DeltaR(stableMuon(aGenParticles,lPart).p4(),aJet.genJet()->p4());
357     if (lDR < deltaRMin) {
358     deltaRMin = lDR;
359     aIndex = lPart;
360     }
361     }
362     if (eid_.at(lPart) && muid_.at(lPart)) {
363     std::cerr << " ---- Problem !! found electron and muon for parton ";
364     printParticle(parton(aGenParticles,lPart));
365     }
366     }
367    
368     if (deltaRMin < aDR) return partonFlavour(aIndex);
369     else return 0;
370     }
371     else return 0;
372    
373     }//lepton match with gen jet
374    
375     unsigned int JetFlavour::partonMatchingRecoJet(const pat::Jet & aJet, const edm::Handle<reco::GenParticleCollection> & aGenParticles, double aDR){//parton match with reco jet
376    
377     int aIndex = -1;
378     double deltaRMin = 1000;
379    
380     for (unsigned int lPart(0); lPart<nPartons(); lPart++){
381     double lDR = ROOT::Math::VectorUtil::DeltaR(parton(aGenParticles,lPart).p4(),aJet.p4());
382     if (lDR < deltaRMin) {
383     deltaRMin = lDR;
384     aIndex = lPart;
385     }
386     }
387    
388     if (deltaRMin < aDR) return partonFlavour(aIndex);
389     else return 0;
390    
391     }//parton match with reco jet
392    
393     unsigned int JetFlavour::leptonMatchingRecoJet(const pat::Jet & aJet, const edm::Handle<reco::GenParticleCollection> & aGenParticles, double aDR){//lepton match with reco jet
394    
395     int aIndex = -1;
396     double deltaRMin = 1000;
397    
398     if (flavour_ > 3 && flavour_ != 21) {
399     for (unsigned int lPart(0); lPart<nPartons(); lPart++){
400     if (eid_.at(lPart)){
401     double lDR = ROOT::Math::VectorUtil::DeltaR(stableElectron(aGenParticles,lPart).p4(),aJet.p4());
402     if (lDR < deltaRMin) {
403     deltaRMin = lDR;
404     aIndex = lPart;
405     }
406     }
407     else if (muid_.at(lPart)){
408     double lDR = ROOT::Math::VectorUtil::DeltaR(stableMuon(aGenParticles,lPart).p4(),aJet.p4());
409     if (lDR < deltaRMin) {
410     deltaRMin = lDR;
411     aIndex = lPart;
412     }
413     }
414     if (eid_.at(lPart) && muid_.at(lPart)) {
415     std::cerr << " ---- Problem !! found electron and muon for parton ";
416     printParticle(parton(aGenParticles,lPart));
417     }
418     }
419    
420     if (deltaRMin < aDR) return partonFlavour(aIndex);
421     else return 0;
422     }
423     else return 0;
424    
425     }//lepton match with reco jet
426    
427     bool JetFlavour::compareParticles(const reco::Candidate * aMum, const reco::GenParticle & aP)
428     {
429    
430     if (!aMum) return false;
431     if ( aMum->pdgId() == aP.pdgId() &&
432     aMum->status() == aP.status() &&
433     aMum->pt() >= aP.pt()-0.001 && aMum->pt() <= aP.pt()+0.001 &&
434     ROOT::Math::VectorUtil::DeltaR(aMum->p4(),aP.p4()) < 0.001 )
435     return true;
436    
437     return false;
438     }
439    
440     void JetFlavour::printParticle(const reco::GenParticle & aP) const{
441     std::cout << " HistosJets::printParticle: " << std::endl <<
442     " ---- pid = " << aP.pdgId() << std::endl <<
443     " ---- status = " << aP.status() << std::endl <<
444     " ---- charge = " << aP.charge() << std::endl <<
445     " ---- pt = " << aP.pt() << std::endl <<
446     //" ---- eta = " << eta_ << std::endl <<
447     //" ---- phi = " << phi_ << std::endl <<
448     //" ---- vertex = " << vx_ << " " << vy_ << " " << vz_ << std::endl <<
449     " ---- parents : " << aP.numberOfMothers() << ". IDs,status = " ;
450     for (unsigned int ip(0); ip<aP.numberOfMothers(); ip++){
451     std::cout << aP.mother(ip)->pdgId() << "," << aP.mother(ip)->status() << " ";
452     }
453     std::cout << std::endl;
454    
455     }
456    
457     void JetFlavour::printSummary() {
458    
459     std::cout <<
460     "=====================================================================" << std::endl <<
461     "======== Summary: " << std::endl <<
462     "======== Looking for partons with flavour = " << flavour_ << std::endl <<
463     "=====================================================================" << std::endl <<
464     "**** total = " << nTot_ << " events," << std::endl <<
465     "**** Z events = " << nZEvts_ << ", W events = " << nWEvts_ << std::endl <<
466     "**** discarded " << nHeavy_ << " Z events with no selected partons," << std::endl <<
467     "**** selected " << nSemiHeavy_ << " Z events with only 1 selected parton." << std::endl <<
468     "=====================================================================" << std::endl <<
469     "**** Selected (Z) events composed of : " << nTotal_ << " quarks, " << std::endl <<
470     "**** --------------------------- " << nLight_ << " u,d,s " << std::endl <<
471     "**** --------------------------- " << nGluon_ << " gluon " << std::endl <<
472     "**** --------------------------- " << nCharm_ << " c " << std::endl <<
473     "**** --------------------------- " << nBottom_ << " b " << std::endl <<
474     "=====================================================================" << std::endl;
475    
476     }
477    
478    
479    
480     }//namespace
481    
482