CMS 3D CMS Logo

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  auto pOutVisTaus = std::make_unique<GenJetCollection>();
33 
34  using namespace GenParticlesHelper;
35 
36  GenParticleRefVector allStatus2Taus;
37  findParticles(*genParticles, allStatus2Taus, 15, 2);
38 
39  for (IGR iTau = allStatus2Taus.begin(); iTau != allStatus2Taus.end(); ++iTau) {
40  // look for all status 1 (stable) descendents
41  GenParticleRefVector descendents;
42  findDescendents(*iTau, descendents, 1);
43 
44  // CV: skip status 2 taus that radiate-off a photon
45  // --> have a status 2 tau lepton in the list of descendents
46  GenParticleRefVector status2TauDaughters;
47  findDescendents(*iTau, status2TauDaughters, 2, 15);
48  if (!status2TauDaughters.empty())
49  continue;
50 
51  // loop on descendents, and take all except neutrinos
52  math::XYZTLorentzVector sumVisMom;
54  Jet::Constituents constituents;
55 
56  if (verbose_)
57  cout << "tau " << (*iTau) << endl;
58 
59  for (IGR igr = descendents.begin(); igr != descendents.end(); ++igr) {
60  int absPdgId = abs((*igr)->pdgId());
61 
62  // neutrinos
63  if (!includeNeutrinos_) {
64  if (absPdgId == 12 || absPdgId == 14 || absPdgId == 16)
65  continue;
66  }
67 
68  if (verbose_)
69  cout << "\t" << (*igr) << endl;
70 
71  charge += (*igr)->charge();
72  sumVisMom += (*igr)->p4();
73 
74  // need to convert the vector of reference to the constituents
75  // to a vector of pointers to build the genjet
76  constituents.push_back(refToPtr(*igr));
77  }
78 
81 
82  GenJet jet(sumVisMom, vertex, specific, constituents);
83 
84  if (charge != (*iTau)->charge())
85  edm::LogError("TauGenJetProducer") << " charge of Tau: " << (*iTau)
86  << " not equal to charge of sum of charge of all descendents.\n"
87  << " Tau's charge: " << (*iTau)->charge() << " sum: " << charge
88  << " # descendents: " << constituents.size() << "\n";
89 
90  jet.setCharge(charge);
91  pOutVisTaus->push_back(jet);
92  }
93  iEvent.put(std::move(pOutVisTaus));
94 }
95 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:99
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_
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
const bool verbose_
verbose ?
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
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
Jets made from MC generator particles.
Definition: GenJet.h:23
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
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
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:228
fixed size matrix
HLT enums.
~TauGenJetProducer() override
TauGenJetProducer(const edm::ParameterSet &)
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:223
def move(src, dest)
Definition: eostools.py:511
void findDescendents(const reco::GenParticleRef &base, reco::GenParticleRefVector &descendents, int status, int pdgId=0)
find all descendents of a given status and pdgId (recursive)