CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PartonSelector.cc
Go to the documentation of this file.
1 //
2 // Select the partons status 2 and 3 for MC Jet Flavour
3 // Author: Attilio
4 // Date: 10.10.2007
5 //
6 
7 //=======================================================================
8 
9 // user include files
13 
19 
22 //#include "DataFormats/Candidate/interface/CandidateFwd.h"
23 //#include "DataFormats/HepMCCandidate/interface/GenParticleCandidate.h"
26 
27 #include <memory>
28 #include <string>
29 #include <iostream>
30 #include <vector>
31 
32 using namespace std;
33 using namespace reco;
34 using namespace edm;
35 
37 {
38  public:
40  ~PartonSelector();
41 
42  private:
43 
44  virtual void produce(edm::Event&, const edm::EventSetup& );
45  bool withLeptons; // Optionally specify leptons
46  bool withTop; // Optionally include top quarks in the list
47  edm::InputTag inputTagGenParticles_; // input collection
48 };
49 //=========================================================================
50 
52 {
53  produces<reco::GenParticleRefVector>();
54  withLeptons = iConfig.getParameter<bool>("withLeptons");
55  inputTagGenParticles_ = iConfig.getParameter<edm::InputTag>("src");
56  if ( iConfig.exists("withTop") ) {
57  withTop = iConfig.getParameter<bool>("withTop");
58  } else {
59  withTop = false;
60  }
61 }
62 
63 //=========================================================================
64 
66 {
67 }
68 
69 // ------------ method called to produce the data ------------
70 
72 {
73 
74  //edm::Handle <reco::CandidateView> particles;
76  iEvent.getByLabel (inputTagGenParticles_, particles );
77  edm::LogVerbatim("PartonSelector") << "=== GenParticle size:" << particles->size();
78  int nPart=0;
79 
80  auto_ptr<GenParticleRefVector> thePartons ( new GenParticleRefVector);
81 
82  for (size_t m = 0; m < particles->size(); m++) {
83 
84  // Don't take into account first 6 particles in generator list
85  if (m<6) continue;
86 
87  const GenParticle & aParticle = (*particles)[ m ];
88 
89  bool isAParton = false;
90  bool isALepton = false;
91  int flavour = abs(aParticle.pdgId());
92  if(flavour == 1 ||
93  flavour == 2 ||
94  flavour == 3 ||
95  flavour == 4 ||
96  flavour == 5 ||
97  (flavour == 6 && withTop) ||
98  flavour == 21 ) isAParton = true;
99  if(flavour == 11 ||
100  flavour == 12 ||
101  flavour == 13 ||
102  flavour == 14 ||
103  flavour == 15 ||
104  flavour == 16 ) isALepton = true;
105 
106 
107  //Add Partons status 3
108  if( aParticle.status() == 3 && isAParton ) {
109  thePartons->push_back( GenParticleRef( particles, m ) );
110  nPart++;
111  }
112 
113  //Add Partons status 2
114  int nparton_daughters = 0;
115  if( aParticle.numberOfDaughters() > 0 && isAParton ) {
116 
117  for (unsigned int i=0; i < aParticle.numberOfDaughters(); i++){
118 
119  int daughterFlavour = abs(aParticle.daughter(i)->pdgId());
120  if( (daughterFlavour == 1 || daughterFlavour == 2 || daughterFlavour == 3 ||
121  daughterFlavour == 4 || daughterFlavour == 5 || daughterFlavour == 6 || daughterFlavour == 21)) {
122  nparton_daughters++;
123  }
124 
125  }
126  if(nparton_daughters == 0){
127  nPart++;
128  thePartons->push_back( GenParticleRef( particles, m ) );
129  }
130 
131  }
132 
133  //Add Leptons
134  // Here you have to decide what to do with taus....
135  // Now all leptons, including e and mu from leptonic tau decays, are added
136  if( withLeptons && aParticle.status() == 3 && isALepton ) {
137  thePartons->push_back( GenParticleRef( particles, m ) );
138  nPart++;
139  }
140  }
141 
142  edm::LogVerbatim("PartonSelector") << "=== GenParticle selected:" << nPart;
143  iEvent.put( thePartons );
144 
145 }
146 
148 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual int pdgId() const
PDG identifier.
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual int status() const
status word
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
bool exists(std::string const &parameterName) const
checks if a parameter exists
#define abs(x)
Definition: mlp_lapack.h:159
edm::InputTag inputTagGenParticles_
int iEvent
Definition: GenABIO.cc:243
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
virtual size_t numberOfDaughters() const
number of daughters
virtual const Candidate * daughter(size_type) const
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
PartonSelector(const edm::ParameterSet &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
virtual int pdgId() const =0
PDG identifier.
virtual void produce(edm::Event &, const edm::EventSetup &)
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31