CMS 3D CMS Logo

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 public:
39  ~PartonSelector() override;
40 
41 private:
42  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
43  bool withLeptons; // Optionally specify leptons
44  bool withTop; // Optionally include top quarks in the list
45  bool acceptNoDaughters; // Parton with zero daugthers are not considered by default, make it configurable
46  unsigned int skipFirstN; // Default skips first 6 particles, make it configurable
48 };
49 //=========================================================================
50 
52  produces<reco::GenParticleRefVector>();
53  withLeptons = iConfig.getParameter<bool>("withLeptons");
54  tokenGenParticles_ = consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("src"));
55  if (iConfig.exists("acceptNoDaughters")) {
56  acceptNoDaughters = iConfig.getParameter<bool>("acceptNoDaughters");
57  } else {
58  acceptNoDaughters = false;
59  }
60  if (iConfig.exists("skipFirstN")) {
61  skipFirstN = iConfig.getParameter<unsigned int>("skipFirstN");
62  } else {
63  skipFirstN = 6;
64  }
65  if (iConfig.exists("withTop")) {
66  withTop = iConfig.getParameter<bool>("withTop");
67  } else {
68  withTop = false;
69  }
70 }
71 
72 //=========================================================================
73 
75 
76 // ------------ method called to produce the data ------------
77 
79  //edm::Handle <reco::CandidateView> particles;
81  iEvent.getByToken(tokenGenParticles_, particles);
82  edm::LogVerbatim("PartonSelector") << "=== GenParticle size:" << particles->size();
83  int nPart = 0;
84 
85  auto thePartons = std::make_unique<GenParticleRefVector>();
86 
87  for (size_t m = 0; m < particles->size(); m++) {
88  // Don't take into account first 6 particles in generator list
89  if (m < skipFirstN)
90  continue;
91 
92  const GenParticle& aParticle = (*particles)[m];
93 
94  bool isAParton = false;
95  bool isALepton = false;
96  int flavour = abs(aParticle.pdgId());
97  if (flavour == 1 || flavour == 2 || flavour == 3 || flavour == 4 || flavour == 5 || (flavour == 6 && withTop) ||
98  flavour == 21)
99  isAParton = true;
100  if (flavour == 11 || flavour == 12 || flavour == 13 || flavour == 14 || flavour == 15 || flavour == 16)
101  isALepton = true;
102 
103  //Add Partons status 3
104  if (aParticle.status() == 3 && isAParton) {
105  thePartons->push_back(GenParticleRef(particles, m));
106  nPart++;
107  }
108 
109  //Add Partons status 2
110  int nparton_daughters = 0;
111  if ((aParticle.numberOfDaughters() > 0 || acceptNoDaughters) && isAParton) {
112  for (unsigned int i = 0; i < aParticle.numberOfDaughters(); i++) {
113  int daughterFlavour = abs(aParticle.daughter(i)->pdgId());
114  if ((daughterFlavour == 1 || daughterFlavour == 2 || daughterFlavour == 3 || daughterFlavour == 4 ||
115  daughterFlavour == 5 || daughterFlavour == 6 || daughterFlavour == 21)) {
116  nparton_daughters++;
117  }
118  }
119  if (nparton_daughters == 0) {
120  nPart++;
121  thePartons->push_back(GenParticleRef(particles, m));
122  }
123  }
124 
125  //Add Leptons
126  // Here you have to decide what to do with taus....
127  // Now all leptons, including e and mu from leptonic tau decays, are added
128  if (withLeptons && aParticle.status() == 3 && isALepton) {
129  thePartons->push_back(GenParticleRef(particles, m));
130  nPart++;
131  }
132  }
133 
134  edm::LogVerbatim("PartonSelector") << "=== GenParticle selected:" << nPart;
135  iEvent.put(std::move(thePartons));
136 }
137 
PartonSelector::~PartonSelector
~PartonSelector() override
Definition: PartonSelector.cc:74
edm::StreamID
Definition: StreamID.h:30
PartonSelector::withLeptons
bool withLeptons
Definition: PartonSelector.cc:43
GenHFHadronMatcher_cff.flavour
flavour
Definition: GenHFHadronMatcher_cff.py:8
Handle.h
mps_fire.i
i
Definition: mps_fire.py:428
reco::CompositeRefCandidateT::numberOfDaughters
size_t numberOfDaughters() const override
number of daughters
MessageLogger.h
ESHandle.h
reco::LeafCandidate::status
int status() const final
status word
Definition: LeafCandidate.h:180
edm::EDGetTokenT< reco::GenParticleCollection >
edm
HLT enums.
Definition: AlignableModifier.h:19
GenParticleRef
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
Definition: MultiTrackValidator.cc:43
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::Handle< reco::GenParticleCollection >
GenParticle
Definition: GenParticle.py:1
ecalTrigSettings_cff.particles
particles
Definition: ecalTrigSettings_cff.py:11
PartonSelector
Definition: PartonSelector.cc:36
GenParticle.h
PartonSelector::PartonSelector
PartonSelector(const edm::ParameterSet &)
Definition: PartonSelector.cc:51
MakerMacros.h
DEFINE_FWK_MODULE
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
visualization-live-secondInstance_cfg.m
m
Definition: visualization-live-secondInstance_cfg.py:72
GenParticleFwd.h
PartonSelector::skipFirstN
unsigned int skipFirstN
Definition: PartonSelector.cc:46
edm::global::EDProducer
Definition: EDProducer.h:32
nPart
TString nPart(Int_t part, TString string, TString delimit=";", Bool_t removerest=true)
edm::ParameterSet::exists
bool exists(std::string const &parameterName) const
checks if a parameter exists
Definition: ParameterSet.cc:681
edm::ParameterSet
Definition: ParameterSet.h:47
hltBtagJetMCTools_cff.withLeptons
withLeptons
Definition: hltBtagJetMCTools_cff.py:7
Event.h
reco::LeafCandidate::pdgId
int pdgId() const final
PDG identifier.
Definition: LeafCandidate.h:176
PartonSelector::tokenGenParticles_
edm::EDGetTokenT< reco::GenParticleCollection > tokenGenParticles_
Definition: PartonSelector.cc:47
iEvent
int iEvent
Definition: GenABIO.cc:224
edm::EventSetup
Definition: EventSetup.h:57
reco::Candidate::pdgId
virtual int pdgId() const =0
PDG identifier.
InputTag.h
eostools.move
def move(src, dest)
Definition: eostools.py:511
std
Definition: JetResolutionObject.h:76
Ref.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
EventSetup.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
reco::CompositeRefCandidateT::daughter
const Candidate * daughter(size_type) const override
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode)
PartonSelector::acceptNoDaughters
bool acceptNoDaughters
Definition: PartonSelector.cc:45
PartonSelector::withTop
bool withTop
Definition: PartonSelector.cc:44
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ParameterSet.h
PartonSelector::produce
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Definition: PartonSelector.cc:78
EDProducer.h
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15