CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TauGenJetProducer.cc
Go to the documentation of this file.
2 
4 
8 
12 
14 
15 using namespace std;
16 using namespace edm;
17 using namespace reco;
18 
20  : inputTagGenParticles_(iConfig.getParameter<InputTag>("GenParticles")),
21  tokenGenParticles_(consumes<GenParticleCollection>(inputTagGenParticles_)),
22  includeNeutrinos_(iConfig.getParameter<bool>("includeNeutrinos")),
23  verbose_(iConfig.getUntrackedParameter<bool>("verbose", false)) {
24  produces<GenJetCollection>();
25 }
26 
28 
31 
32  bool found = iEvent.getByToken(tokenGenParticles_, genParticles);
33 
34  if (!found) {
35  std::ostringstream err;
36  err << " cannot get collection: " << inputTagGenParticles_ << std::endl;
37  edm::LogError("TauGenJetProducer") << err.str();
38  throw cms::Exception("MissingProduct", err.str());
39  }
40 
41  auto pOutVisTaus = std::make_unique<GenJetCollection>();
42 
43  using namespace GenParticlesHelper;
44 
45  GenParticleRefVector allStatus2Taus;
46  findParticles(*genParticles, allStatus2Taus, 15, 2);
47 
48  for (IGR iTau = allStatus2Taus.begin(); iTau != allStatus2Taus.end(); ++iTau) {
49  // look for all status 1 (stable) descendents
50  GenParticleRefVector descendents;
51  findDescendents(*iTau, descendents, 1);
52 
53  // CV: skip status 2 taus that radiate-off a photon
54  // --> have a status 2 tau lepton in the list of descendents
55  GenParticleRefVector status2TauDaughters;
56  findDescendents(*iTau, status2TauDaughters, 2, 15);
57  if (!status2TauDaughters.empty())
58  continue;
59 
60  // loop on descendents, and take all except neutrinos
61  math::XYZTLorentzVector sumVisMom;
63  Jet::Constituents constituents;
64 
65  if (verbose_)
66  cout << "tau " << (*iTau) << endl;
67 
68  for (IGR igr = descendents.begin(); igr != descendents.end(); ++igr) {
69  int absPdgId = abs((*igr)->pdgId());
70 
71  // neutrinos
72  if (!includeNeutrinos_) {
73  if (absPdgId == 12 || absPdgId == 14 || absPdgId == 16)
74  continue;
75  }
76 
77  if (verbose_)
78  cout << "\t" << (*igr) << endl;
79 
80  charge += (*igr)->charge();
81  sumVisMom += (*igr)->p4();
82 
83  // need to convert the vector of reference to the constituents
84  // to a vector of pointers to build the genjet
85  constituents.push_back(refToPtr(*igr));
86  }
87 
88  math::XYZPoint vertex;
90 
91  GenJet jet(sumVisMom, vertex, specific, constituents);
92 
93  if (charge != (*iTau)->charge())
94  std::cout << " charge of Tau: " << (*iTau) << " not equal to charge of sum of charge of all descendents. "
95  << std::endl;
96 
97  jet.setCharge(charge);
98  pOutVisTaus->push_back(jet);
99  }
100  iEvent.put(std::move(pOutVisTaus));
101 }
102 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
dictionary specific
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
Ptr< typename C::value_type > refToPtr(Ref< C, typename C::value_type, refhelper::FindUsingAdvance< C, typename C::value_type > > const &ref)
Definition: RefToPtr.h:18
void findParticles(const reco::GenParticleCollection &sourceParticles, reco::GenParticleRefVector &particleRefs, int pdgId, int status)
find all particles of a given pdgId and status
const edm::EDGetTokenT< reco::GenParticleCollection > tokenGenParticles_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::InputTag inputTagGenParticles_
Input PFCandidates.
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
Log< level::Error, false > LogError
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
const bool verbose_
verbose ?
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
void setCharge(Charge q) final
set electric charge
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:224
const bool includeNeutrinos_
if yes, neutrinos will be included, for debug purposes
def move
Definition: eostools.py:511
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:23
builds a GenJet from the visible daughters of each status 2 tau in the event.
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
~TauGenJetProducer() override
TauGenJetProducer(const edm::ParameterSet &)
tuple cout
Definition: gather_cfg.py:144
void findDescendents(const reco::GenParticleRef &base, reco::GenParticleRefVector &descendents, int status, int pdgId=0)
find all descendents of a given status and pdgId (recursive)