CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 {
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  std::auto_ptr<GenJetCollection>
48  pOutVisTaus(new GenJetCollection());
49 
50  using namespace GenParticlesHelper;
51 
52  GenParticleRefVector allStatus2Taus;
53  findParticles( *genParticles,
54  allStatus2Taus, 15, 2);
55 
56  for ( IGR iTau=allStatus2Taus.begin(); iTau!=allStatus2Taus.end(); ++iTau ) {
57 
58  // look for all status 1 (stable) descendents
59  GenParticleRefVector descendents;
60  findDescendents( *iTau, descendents, 1);
61 
62  // CV: skip status 2 taus that radiate-off a photon
63  // --> have a status 2 tau lepton in the list of descendents
64  GenParticleRefVector status2TauDaughters;
65  findDescendents( *iTau, status2TauDaughters, 2, 15 );
66  if ( status2TauDaughters.size() > 0 ) continue;
67 
68  // loop on descendents, and take all except neutrinos
69  math::XYZTLorentzVector sumVisMom;
70  Particle::Charge charge = 0;
71  Jet::Constituents constituents;
72 
73  if(verbose_)
74  cout<<"tau "<<(*iTau)<<endl;
75 
76  for(IGR igr = descendents.begin();
77  igr!= descendents.end(); ++igr ) {
78 
79  int absPdgId = abs((*igr)->pdgId());
80 
81  // neutrinos
82  if(!includeNeutrinos_ ) {
83  if( absPdgId == 12 ||
84  absPdgId == 14 ||
85  absPdgId == 16 )
86  continue;
87  }
88 
89  if(verbose_)
90  cout<<"\t"<<(*igr)<<endl;
91 
92  charge += (*igr)->charge();
93  sumVisMom += (*igr)->p4();
94 
95  // need to convert the vector of reference to the constituents
96  // to a vector of pointers to build the genjet
97  constituents.push_back( refToPtr( *igr) );
98  }
99 
100  math::XYZPoint vertex;
102 
103  GenJet jet( sumVisMom, vertex, specific, constituents);
104 
105  if (charge != (*iTau)->charge() )
106  std::cout<<" charge of Tau: " << (*iTau) << " not equal to charge of sum of charge of all descendents. " << std::endl;
107 
108  jet.setCharge(charge);
109  pOutVisTaus->push_back( jet );
110 
111  }
112  iEvent.put( pOutVisTaus );
113 }
114 
115 
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
dictionary specific
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:462
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
const edm::InputTag inputTagGenParticles_
Input PFCandidates.
virtual void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
std::vector< GenJet > GenJetCollection
collection of GenJet objects
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
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
const bool includeNeutrinos_
if yes, neutrinos will be included, for debug purposes
virtual void setCharge(Charge q) final
set electric charge
Definition: LeafCandidate.h:93
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
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
TauGenJetProducer(const edm::ParameterSet &)
tuple cout
Definition: gather_cfg.py:145
volatile std::atomic< bool > shutdown_flag false
void findDescendents(const reco::GenParticleRef &base, reco::GenParticleRefVector &descendents, int status, int pdgId=0)
find all descendents of a given status and pdgId (recursive)