CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
BoostedTauSeedsProducer.cc
Go to the documentation of this file.
1 
2 /* =============================================================================
3  * Filename: BoostedTauSeedsProducer.cc
4  *
5  * Description: Take the two subjets found by CMSBoostedTauSeedingAlgorithm
6  * and add the data-formats for
7  * o seeding the reconstruction of hadronic taus (PFJets, collection of PFCandidates not within subjet)
8  * o computation of electron and muon isolation Pt-sums, excluding the particles in the other subjet (collection of PFCandidates within subjet, used to define "Vetos")
9  *
10  * Created: 10/22/2013 16:05:00
11  *
12  * Authors: Christian Veelken (LLR)
13  *
14  * =============================================================================
15  */
16 
24 
36 
37 #include <boost/foreach.hpp>
38 
39 #include <TMath.h>
40 
41 #include <string>
42 #include <iostream>
43 #include <iomanip>
44 
46 {
47  public:
50  virtual void produce(edm::Event&, const edm::EventSetup&);
51 
52  private:
53  typedef edm::AssociationMap<edm::OneToMany<std::vector<reco::PFJet>, std::vector<reco::PFCandidate>, unsigned int> > JetToPFCandidateAssociation;
54 
56 
60 
62 };
63 
65  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
66 {
67  srcSubjets_ = consumes<JetView>(cfg.getParameter<edm::InputTag>("subjetSrc"));
68  srcPFCandidates_ = consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("pfCandidateSrc"));
69 
70  verbosity_ = ( cfg.exists("verbosity") ) ?
71  cfg.getParameter<int>("verbosity") : 0;
72 
73  produces<reco::PFJetCollection>();
74  produces<JetToPFCandidateAssociation>("pfCandAssocMapForIsolation");
75  //produces<JetToPFCandidateAssociation>("pfCandAssocMapForIsoDepositVetos");
76 }
77 
78 namespace
79 {
80  reco::PFJet convertToPFJet(const reco::Jet& jet, const reco::Jet::Constituents& jetConstituents)
81  {
82  // CV: code for filling pfJetSpecific objects taken from
83  // RecoParticleFlow/PFRootEvent/src/JetMaker.cc
84  double chargedHadronEnergy = 0.;
85  double neutralHadronEnergy = 0.;
86  double chargedEmEnergy = 0.;
87  double neutralEmEnergy = 0.;
88  double chargedMuEnergy = 0.;
89  int chargedMultiplicity = 0;
90  int neutralMultiplicity = 0;
91  int muonMultiplicity = 0;
92  for ( reco::Jet::Constituents::const_iterator jetConstituent = jetConstituents.begin();
93  jetConstituent != jetConstituents.end(); ++jetConstituent ) {
94  const reco::PFCandidate* pfCandidate = dynamic_cast<const reco::PFCandidate*>(jetConstituent->get());
95  if ( pfCandidate ) {
96  switch ( pfCandidate->particleId() ) {
97  case reco::PFCandidate::h : // charged hadron
98  chargedHadronEnergy += pfCandidate->energy();
99  ++chargedMultiplicity;
100  break;
101  case reco::PFCandidate::e : // electron
102  chargedEmEnergy += pfCandidate->energy();
103  ++chargedMultiplicity;
104  break;
105  case reco::PFCandidate::mu : // muon
106  chargedMuEnergy += pfCandidate->energy();
107  ++chargedMultiplicity;
109  break;
110  case reco::PFCandidate::gamma : // photon
111  case reco::PFCandidate::egamma_HF : // electromagnetic in HF
112  neutralEmEnergy += pfCandidate->energy();
113  ++neutralMultiplicity;
114  break;
115  case reco::PFCandidate::h0 : // neutral hadron
116  case reco::PFCandidate::h_HF : // hadron in HF
117  neutralHadronEnergy += pfCandidate->energy();
118  ++neutralMultiplicity;
119  break;
120  default:
121  edm::LogWarning("convertToPFJet")
122  << "PFCandidate: Pt = " << pfCandidate->pt() << ", eta = " << pfCandidate->eta() << ", phi = " << pfCandidate->phi()
123  << " has invalid particleID = " << pfCandidate->particleId() << " !!" << std::endl;
124  break;
125  }
126  } else {
127  edm::LogWarning("convertToPFJet")
128  << "Jet constituent: Pt = " << pfCandidate->pt() << ", eta = " << pfCandidate->eta() << ", phi = " << pfCandidate->phi()
129  << " is not of type PFCandidate !!" << std::endl;
130  }
131  }
132  reco::PFJet::Specific pfJetSpecific;
133  pfJetSpecific.mChargedHadronEnergy = chargedHadronEnergy;
134  pfJetSpecific.mNeutralHadronEnergy = neutralHadronEnergy;
135  pfJetSpecific.mChargedEmEnergy = chargedEmEnergy;
136  pfJetSpecific.mChargedMuEnergy = chargedMuEnergy;
137  pfJetSpecific.mNeutralEmEnergy = neutralEmEnergy;
138  pfJetSpecific.mChargedMultiplicity = chargedMultiplicity;
139  pfJetSpecific.mNeutralMultiplicity = neutralMultiplicity;
140  pfJetSpecific.mMuonMultiplicity = muonMultiplicity;
141 
142  reco::PFJet pfJet(jet.p4(), jet.vertex(), pfJetSpecific, jetConstituents);
143  pfJet.setJetArea(jet.jetArea());
144 
145  return pfJet;
146  }
147 
148  void getJetConstituents(const reco::Jet& jet, reco::Jet::Constituents& jet_and_subjetConstituents)
149  {
150  reco::Jet::Constituents jetConstituents = jet.getJetConstituents();
151  for ( reco::Jet::Constituents::const_iterator jetConstituent = jetConstituents.begin();
152  jetConstituent != jetConstituents.end(); ++jetConstituent ) {
153  const reco::Jet* subjet = dynamic_cast<const reco::Jet*>(jetConstituent->get());
154  if ( subjet ) {
155  getJetConstituents(*subjet, jet_and_subjetConstituents);
156  } else {
157  jet_and_subjetConstituents.push_back(*jetConstituent);
158  }
159  }
160  }
161 
162  std::vector<reco::PFCandidateRef> getPFCandidates_exclJetConstituents(const reco::Jet& jet, const edm::Handle<reco::PFCandidateCollection>& pfCandidates, const reco::Jet::Constituents& jetConstituents, double dRmatch, bool invert)
163  {
164  auto const & collection_cand = (*pfCandidates);
165  std::vector<reco::PFCandidateRef> pfCandidates_exclJetConstituents;
166  size_t numPFCandidates = pfCandidates->size();
167  for ( size_t pfCandidateIdx = 0; pfCandidateIdx < numPFCandidates; ++pfCandidateIdx ) {
168  if(!(deltaR(collection_cand[pfCandidateIdx].p4(), jet.p4())<1.0)) continue;
169  bool isJetConstituent = false;
170  for ( reco::Jet::Constituents::const_iterator jetConstituent = jetConstituents.begin();
171  jetConstituent != jetConstituents.end(); ++jetConstituent ) {
172  double dR = deltaR(collection_cand[pfCandidateIdx].p4(), (*jetConstituent)->p4());
173  if ( dR < dRmatch ) {
174  isJetConstituent = true;
175  break;
176  }
177  }
178  if ( !(isJetConstituent^invert) ) {
179  reco::PFCandidateRef pfCandidate(pfCandidates, pfCandidateIdx);
180  pfCandidates_exclJetConstituents.push_back(pfCandidate);
181  }
182  }
183  return pfCandidates_exclJetConstituents;
184  }
185 
186  void printJetConstituents(const std::string& label, const reco::Jet::Constituents& jetConstituents)
187  {
188  std::cout << "#" << label << " = " << jetConstituents.size() << ":" << std::endl;
189  int idx = 0;
190  for ( reco::Jet::Constituents::const_iterator jetConstituent = jetConstituents.begin();
191  jetConstituent != jetConstituents.end(); ++jetConstituent ) {
192  std::cout << " jetConstituent #" << idx << ": Pt = " << (*jetConstituent)->pt() << ", eta = " << (*jetConstituent)->eta() << ", phi = " << (*jetConstituent)->phi() << std::endl;
193  ++idx;
194  }
195  }
196 }
197 
199 {
200  if ( verbosity_ >= 1 ) {
201  std::cout << "<BoostedTauSeedsProducer::produce (moduleLabel = " << moduleLabel_ << ")>:" << std::endl;
202  }
203 
204  edm::Handle<JetView> subjets;
205  evt.getByToken(srcSubjets_, subjets);
206  if ( verbosity_ >= 1 ) {
207  std::cout << "#subjets = " << subjets->size() << std::endl;
208  }
209  assert((subjets->size() % 2) == 0); // CV: ensure that subjets come in pairs
210 
212  evt.getByToken(srcPFCandidates_, pfCandidates);
213  if ( verbosity_ >= 1 ) {
214  std::cout << "#pfCandidates = " << pfCandidates->size() << std::endl;
215  }
216 
217  std::auto_ptr<reco::PFJetCollection> selectedSubjets(new reco::PFJetCollection());
219 
220  std::auto_ptr<JetToPFCandidateAssociation> selectedSubjetPFCandidateAssociationForIsolation(new JetToPFCandidateAssociation(&evt.productGetter()));
221  //std::auto_ptr<JetToPFCandidateAssociation> selectedSubjetPFCandidateAssociationForIsoDepositVetos(new JetToPFCandidateAssociation(&evt.productGetter()));
222 
223  for ( size_t idx = 0; idx < (subjets->size() / 2); ++idx ) {
224  const reco::Jet* subjet1 = &subjets->at(2*idx);
225  const reco::Jet* subjet2 = &subjets->at(2*idx + 1);
226  assert(subjet1 && subjet2);
227  if ( verbosity_ >= 1 ) {
228  std::cout << "processing jet #" << idx << ":" << std::endl;
229  std::cout << " subjet1: Pt = " << subjet1->pt() << ", eta = " << subjet1->eta() << ", phi = " << subjet1->phi() << ", mass = " << subjet1->mass()
230  << " (#constituents = " << subjet1->nConstituents() << ", area = " << subjet1->jetArea() << ")" << std::endl;
231  std::cout << " subjet2: Pt = " << subjet2->pt() << ", eta = " << subjet2->eta() << ", phi = " << subjet2->phi() << ", mass = " << subjet2->mass()
232  << " (#constituents = " << subjet2->nConstituents() << ", area = " << subjet2->jetArea() << ")" << std::endl;
233  }
234 
235  if ( !(subjet1->nConstituents() >= 1 && subjet1->pt() > 1. &&
236  subjet2->nConstituents() >= 1 && subjet2->pt() > 1.) ) continue; // CV: skip pathological cases
237 
238  // find PFCandidate constituents of each subjet
239  reco::Jet::Constituents subjetConstituents1;
240  getJetConstituents(*subjet1, subjetConstituents1);
241  reco::Jet::Constituents subjetConstituents2;
242  getJetConstituents(*subjet2, subjetConstituents2);
243  if ( verbosity_ >= 1 ) {
244  printJetConstituents("subjetConstituents1", subjetConstituents1);
245  printJetConstituents("subjetConstituents2", subjetConstituents2);
246  }
247 
248  selectedSubjets->push_back(convertToPFJet(*subjet1, subjetConstituents1));
249  edm::Ref<reco::PFJetCollection> subjetRef1(selectedSubjetRefProd, selectedSubjets->size() - 1);
250  selectedSubjets->push_back(convertToPFJet(*subjet2, subjetConstituents2));
251  edm::Ref<reco::PFJetCollection> subjetRef2(selectedSubjetRefProd, selectedSubjets->size() - 1);
252 
253  // find all PFCandidates that are not constituents of the **other** subjet
254  std::vector<reco::PFCandidateRef> pfCandidatesNotInSubjet1 = getPFCandidates_exclJetConstituents(*subjet1, pfCandidates, subjetConstituents2, 1.e-4, false);
255  std::vector<reco::PFCandidateRef> pfCandidatesNotInSubjet2 = getPFCandidates_exclJetConstituents(*subjet2, pfCandidates, subjetConstituents1, 1.e-4, false);
256  if ( verbosity_ >= 1 ) {
257  std::cout << "#pfCandidatesNotInSubjet1 = " << pfCandidatesNotInSubjet1.size() << std::endl;
258  std::cout << "#pfCandidatesNotInSubjet2 = " << pfCandidatesNotInSubjet2.size() << std::endl;
259  }
260 
261  // build JetToPFCandidateAssociation
262  // (key = subjet, value = collection of PFCandidates that are not constituents of subjet)
263  BOOST_FOREACH( const reco::PFCandidateRef& pfCandidate, pfCandidatesNotInSubjet1 ) {
264  selectedSubjetPFCandidateAssociationForIsolation->insert(subjetRef1, pfCandidate);
265  }
266  BOOST_FOREACH( const reco::PFCandidateRef& pfCandidate, pfCandidatesNotInSubjet2 ) {
267  selectedSubjetPFCandidateAssociationForIsolation->insert(subjetRef2, pfCandidate);
268  }
269  }
270 
271  evt.put(selectedSubjets);
272  evt.put(selectedSubjetPFCandidateAssociationForIsolation, "pfCandAssocMapForIsolation");
273 }
274 
T getParameter(std::string const &) const
tuple cfg
Definition: looper.py:293
virtual double energy() const final
energy
int mChargedMultiplicity
Definition: PFJet.h:75
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
virtual const Point & vertex() const
vertex position (overwritten by PF...)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Base class for all types of Jets.
Definition: Jet.h:20
EDProductGetter const & productGetter() const
Definition: Event.cc:55
assert(m_qm.get())
bool exists(std::string const &parameterName) const
checks if a parameter exists
virtual double phi() const final
momentum azimuthal angle
std::vector< Constituent > Constituents
Definition: Jet.h:23
float mNeutralHadronEnergy
Definition: PFJet.h:55
float mChargedMuEnergy
Definition: PFJet.h:73
virtual Constituents getJetConstituents() const
list of constituents
Jets made from PFObjects.
Definition: PFJet.h:21
float mChargedHadronEnergy
Definition: PFJet.h:54
BoostedTauSeedsProducer(const edm::ParameterSet &)
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:103
float mChargedEmEnergy
Definition: PFJet.h:72
virtual void produce(edm::Event &, const edm::EventSetup &)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
double p4[4]
Definition: TauolaWrapper.h:92
edm::View< reco::Jet > JetView
float mNeutralEmEnergy
Definition: PFJet.h:74
edm::EDGetTokenT< reco::PFCandidateCollection > srcPFCandidates_
virtual double mass() const final
mass
RefProd< PROD > getRefBeforePut()
Definition: Event.h:141
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
edm::AssociationMap< edm::OneToMany< std::vector< reco::PFJet >, std::vector< reco::PFCandidate >, unsigned int > > JetToPFCandidateAssociation
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
int mNeutralMultiplicity
Definition: PFJet.h:76
std::vector< PFJet > PFJetCollection
collection of PFJet objects
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:39
virtual float jetArea() const
get jet area
Definition: Jet.h:105
virtual int nConstituents() const
of constituents
Definition: Jet.h:65
edm::EDGetTokenT< JetView > srcSubjets_
tuple cout
Definition: gather_cfg.py:145
virtual double eta() const final
momentum pseudorapidity
virtual ParticleType particleId() const
Definition: PFCandidate.h:373
moduleLabel_(iConfig.getParameter< string >("@module_label"))
virtual const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
virtual double pt() const final
transverse momentum