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 {
25 
26 
27  produces<GenJetCollection>();
28 }
29 
31 
33  const EventSetup& iSetup) const {
34 
36 
37  bool found = iEvent.getByToken( tokenGenParticles_, genParticles);
38 
39  if ( !found ) {
40  std::ostringstream err;
41  err<<" cannot get collection: "
42  <<inputTagGenParticles_<<std::endl;
43  edm::LogError("TauGenJetProducer")<<err.str();
44  throw cms::Exception( "MissingProduct", err.str());
45  }
46 
47  auto pOutVisTaus = std::make_unique<GenJetCollection>();
48 
49  using namespace GenParticlesHelper;
50 
51  GenParticleRefVector allStatus2Taus;
52  findParticles( *genParticles,
53  allStatus2Taus, 15, 2);
54 
55  for ( IGR iTau=allStatus2Taus.begin(); iTau!=allStatus2Taus.end(); ++iTau ) {
56 
57  // look for all status 1 (stable) descendents
58  GenParticleRefVector descendents;
59  findDescendents( *iTau, descendents, 1);
60 
61  // CV: skip status 2 taus that radiate-off a photon
62  // --> have a status 2 tau lepton in the list of descendents
63  GenParticleRefVector status2TauDaughters;
64  findDescendents( *iTau, status2TauDaughters, 2, 15 );
65  if ( status2TauDaughters.size() > 0 ) continue;
66 
67  // loop on descendents, and take all except neutrinos
68  math::XYZTLorentzVector sumVisMom;
69  Particle::Charge charge = 0;
70  Jet::Constituents constituents;
71 
72  if(verbose_)
73  cout<<"tau "<<(*iTau)<<endl;
74 
75  for(IGR igr = descendents.begin();
76  igr!= descendents.end(); ++igr ) {
77 
78  int absPdgId = abs((*igr)->pdgId());
79 
80  // neutrinos
81  if(!includeNeutrinos_ ) {
82  if( absPdgId == 12 ||
83  absPdgId == 14 ||
84  absPdgId == 16 )
85  continue;
86  }
87 
88  if(verbose_)
89  cout<<"\t"<<(*igr)<<endl;
90 
91  charge += (*igr)->charge();
92  sumVisMom += (*igr)->p4();
93 
94  // need to convert the vector of reference to the constituents
95  // to a vector of pointers to build the genjet
96  constituents.push_back( refToPtr( *igr) );
97  }
98 
99  math::XYZPoint vertex;
100  GenJet::Specific specific;
101 
102  GenJet jet( sumVisMom, vertex, specific, constituents);
103 
104  if (charge != (*iTau)->charge() )
105  std::cout<<" charge of Tau: " << (*iTau) << " not equal to charge of sum of charge of all descendents. " << std::endl;
106 
107  jet.setCharge(charge);
108  pOutVisTaus->push_back( jet );
109 
110  }
111  iEvent.put(std::move(pOutVisTaus) );
112 }
113 
114 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
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:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const edm::InputTag inputTagGenParticles_
Input PFCandidates.
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:253
const bool verbose_
verbose ?
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:248
void setCharge(Charge q) final
set electric charge
Definition: LeafCandidate.h:93
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:230
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:24
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
fixed size matrix
HLT enums.
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
TauGenJetProducer(const edm::ParameterSet &)
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
def move(src, dest)
Definition: eostools.py:510
void findDescendents(const reco::GenParticleRef &base, reco::GenParticleRefVector &descendents, int status, int pdgId=0)
find all descendents of a given status and pdgId (recursive)