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& ) override;
45  bool withLeptons; // Optionally specify leptons
46  bool withTop; // Optionally include top quarks in the list
47  bool acceptNoDaughters; // Parton with zero daugthers are not considered by default, make it configurable
48  unsigned int skipFirstN; // Default skips first 6 particles, make it configurable
50 };
51 //=========================================================================
52 
54 {
55  produces<reco::GenParticleRefVector>();
56  withLeptons = iConfig.getParameter<bool>("withLeptons");
57  tokenGenParticles_ = consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("src"));
58  if ( iConfig.exists("acceptNoDaughters") ) {
59  acceptNoDaughters = iConfig.getParameter<bool>("acceptNoDaughters");
60  } else {
61  acceptNoDaughters=false;
62  }
63  if ( iConfig.exists("skipFirstN") ) {
64  skipFirstN = iConfig.getParameter<unsigned int>("skipFirstN");
65  } else {
66  skipFirstN=6;
67  }
68  if ( iConfig.exists("withTop") ) {
69  withTop = iConfig.getParameter<bool>("withTop");
70  } else {
71  withTop = false;
72  }
73 }
74 
75 //=========================================================================
76 
78 {
79 }
80 
81 // ------------ method called to produce the data ------------
82 
84 {
85 
86  //edm::Handle <reco::CandidateView> particles;
88  iEvent.getByToken (tokenGenParticles_, particles );
89  edm::LogVerbatim("PartonSelector") << "=== GenParticle size:" << particles->size();
90  int nPart=0;
91 
92  auto_ptr<GenParticleRefVector> thePartons ( new GenParticleRefVector);
93 
94  for (size_t m = 0; m < particles->size(); m++) {
95 
96  // Don't take into account first 6 particles in generator list
97  if (m<skipFirstN) continue;
98 
99  const GenParticle & aParticle = (*particles)[ m ];
100 
101  bool isAParton = false;
102  bool isALepton = false;
103  int flavour = abs(aParticle.pdgId());
104  if(flavour == 1 ||
105  flavour == 2 ||
106  flavour == 3 ||
107  flavour == 4 ||
108  flavour == 5 ||
109  (flavour == 6 && withTop) ||
110  flavour == 21 ) isAParton = true;
111  if(flavour == 11 ||
112  flavour == 12 ||
113  flavour == 13 ||
114  flavour == 14 ||
115  flavour == 15 ||
116  flavour == 16 ) isALepton = true;
117 
118 
119  //Add Partons status 3
120  if( aParticle.status() == 3 && isAParton ) {
121  thePartons->push_back( GenParticleRef( particles, m ) );
122  nPart++;
123  }
124 
125  //Add Partons status 2
126  int nparton_daughters = 0;
127  if( ( aParticle.numberOfDaughters() > 0 || acceptNoDaughters) && isAParton ) {
128 
129  for (unsigned int i=0; i < aParticle.numberOfDaughters(); i++){
130 
131  int daughterFlavour = abs(aParticle.daughter(i)->pdgId());
132  if( (daughterFlavour == 1 || daughterFlavour == 2 || daughterFlavour == 3 ||
133  daughterFlavour == 4 || daughterFlavour == 5 || daughterFlavour == 6 || daughterFlavour == 21)) {
134  nparton_daughters++;
135  }
136 
137  }
138  if(nparton_daughters == 0){
139  nPart++;
140  thePartons->push_back( GenParticleRef( particles, m ) );
141  }
142 
143  }
144 
145  //Add Leptons
146  // Here you have to decide what to do with taus....
147  // Now all leptons, including e and mu from leptonic tau decays, are added
148  if( withLeptons && aParticle.status() == 3 && isALepton ) {
149  thePartons->push_back( GenParticleRef( particles, m ) );
150  nPart++;
151  }
152  }
153 
154  edm::LogVerbatim("PartonSelector") << "=== GenParticle selected:" << nPart;
155  iEvent.put( thePartons );
156 
157 }
158 
160 
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
bool exists(std::string const &parameterName) const
checks if a parameter exists
edm::EDGetTokenT< reco::GenParticleCollection > tokenGenParticles_
virtual int status() const final
status word
int iEvent
Definition: GenABIO.cc:230
unsigned int skipFirstN
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
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) ...
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
PartonSelector(const edm::ParameterSet &)
virtual int pdgId() const =0
PDG identifier.
TString nPart(Int_t part, TString string, TString delimit=";", Bool_t removerest=true)
virtual int pdgId() const final
PDG identifier.
virtual void produce(edm::Event &, const edm::EventSetup &) override
int flavour(const Candidate &part)
Definition: pdgIdUtils.h:31