CMS 3D CMS Logo

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 <TMath.h>
38 
39 #include <string>
40 #include <iostream>
41 #include <iomanip>
42 
44 {
45  public:
48  void produce(edm::Event&, const edm::EventSetup&) override;
49 
50  private:
51  typedef edm::AssociationMap<edm::OneToMany<std::vector<reco::PFJet>, std::vector<reco::PFCandidate>, unsigned int> > JetToPFCandidateAssociation;
52 
54 
58 
60 };
61 
63  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
64 {
65  srcSubjets_ = consumes<JetView>(cfg.getParameter<edm::InputTag>("subjetSrc"));
66  srcPFCandidates_ = consumes<reco::PFCandidateCollection>(cfg.getParameter<edm::InputTag>("pfCandidateSrc"));
67 
68  verbosity_ = ( cfg.exists("verbosity") ) ?
69  cfg.getParameter<int>("verbosity") : 0;
70 
71  produces<reco::PFJetCollection>();
72  produces<JetToPFCandidateAssociation>("pfCandAssocMapForIsolation");
73  //produces<JetToPFCandidateAssociation>("pfCandAssocMapForIsoDepositVetos");
74 }
75 
76 namespace
77 {
78  typedef std::vector<std::unordered_set<uint32_t> > JetToConstitMap;
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 JetToConstitMap::value_type& constitmap, const reco::Jet::Constituents& jetConstituents, double /*dRmatch*/, bool invert)
163  {
164  //const double dRmatch2 = dRmatch*dRmatch; // comment out for now in case someone needs a dR-based search again
165  auto const & collection_cand = (*pfCandidates);
166  std::vector<reco::PFCandidateRef> pfCandidates_exclJetConstituents;
167  size_t numPFCandidates = pfCandidates->size();
168  for ( size_t pfCandidateIdx = 0; pfCandidateIdx < numPFCandidates; ++pfCandidateIdx ) {
169  if(!(deltaR2(collection_cand[pfCandidateIdx], jet)<1.0)) continue;
170  bool isJetConstituent = constitmap.count(pfCandidateIdx);
171  if ( !(isJetConstituent^invert) ) {
172  reco::PFCandidateRef pfCandidate(pfCandidates, pfCandidateIdx);
173  pfCandidates_exclJetConstituents.push_back(pfCandidate);
174  }
175  }
176  return pfCandidates_exclJetConstituents;
177  }
178 
179  void printJetConstituents(const std::string& label, const reco::Jet::Constituents& jetConstituents)
180  {
181  std::cout << "#" << label << " = " << jetConstituents.size() << ":" << std::endl;
182  int idx = 0;
183  for ( reco::Jet::Constituents::const_iterator jetConstituent = jetConstituents.begin();
184  jetConstituent != jetConstituents.end(); ++jetConstituent ) {
185  std::cout << " jetConstituent #" << idx << ": Pt = " << (*jetConstituent)->pt() << ", eta = " << (*jetConstituent)->eta() << ", phi = " << (*jetConstituent)->phi() << std::endl;
186  ++idx;
187  }
188  }
189 }
190 
192 {
193  if ( verbosity_ >= 1 ) {
194  std::cout << "<BoostedTauSeedsProducer::produce (moduleLabel = " << moduleLabel_ << ")>:" << std::endl;
195  }
196 
197  edm::Handle<JetView> subjets;
198  evt.getByToken(srcSubjets_, subjets);
199  if ( verbosity_ >= 1 ) {
200  std::cout << "#subjets = " << subjets->size() << std::endl;
201  }
202  assert((subjets->size() % 2) == 0); // CV: ensure that subjets come in pairs
203 
205  evt.getByToken(srcPFCandidates_, pfCandidates);
206  if ( verbosity_ >= 1 ) {
207  std::cout << "#pfCandidates = " << pfCandidates->size() << std::endl;
208  }
209 
210  auto selectedSubjets = std::make_unique<reco::PFJetCollection>();
212 
213  auto selectedSubjetPFCandidateAssociationForIsolation = std::make_unique<JetToPFCandidateAssociation>(&evt.productGetter());
214  //auto selectedSubjetPFCandidateAssociationForIsoDepositVetos = std::make_unique<JetToPFCandidateAssociation>(&evt.productGetter());
215 
216  // cache for jet->pfcandidate
217  JetToConstitMap constitmap(subjets->size());
218 
219  // fill constituents map
220  const auto& thesubjets = *subjets;
221  for( unsigned i = 0; i < thesubjets.size(); ++i ) {
222  for ( unsigned j = 0; j < thesubjets[i].numberOfDaughters(); ++j ) {
223  constitmap[i].emplace(thesubjets[i].daughterPtr(j).key());
224  }
225  }
226 
227  for ( size_t idx = 0; idx < (subjets->size() / 2); ++idx ) {
228  const reco::Jet* subjet1 = &subjets->at(2*idx);
229  const reco::Jet* subjet2 = &subjets->at(2*idx + 1);
230  assert(subjet1 && subjet2);
231  if ( verbosity_ >= 1 ) {
232  std::cout << "processing jet #" << idx << ":" << std::endl;
233  std::cout << " subjet1: Pt = " << subjet1->pt() << ", eta = " << subjet1->eta() << ", phi = " << subjet1->phi() << ", mass = " << subjet1->mass()
234  << " (#constituents = " << subjet1->nConstituents() << ", area = " << subjet1->jetArea() << ")" << std::endl;
235  std::cout << " subjet2: Pt = " << subjet2->pt() << ", eta = " << subjet2->eta() << ", phi = " << subjet2->phi() << ", mass = " << subjet2->mass()
236  << " (#constituents = " << subjet2->nConstituents() << ", area = " << subjet2->jetArea() << ")" << std::endl;
237  }
238 
239  if ( !(subjet1->nConstituents() >= 1 && subjet1->pt() > 1. &&
240  subjet2->nConstituents() >= 1 && subjet2->pt() > 1.) ) continue; // CV: skip pathological cases
241 
242  // find PFCandidate constituents of each subjet
243  reco::Jet::Constituents subjetConstituents1;
244  getJetConstituents(*subjet1, subjetConstituents1);
245  reco::Jet::Constituents subjetConstituents2;
246  getJetConstituents(*subjet2, subjetConstituents2);
247  if ( verbosity_ >= 1 ) {
248  printJetConstituents("subjetConstituents1", subjetConstituents1);
249  printJetConstituents("subjetConstituents2", subjetConstituents2);
250  }
251 
252  selectedSubjets->push_back(convertToPFJet(*subjet1, subjetConstituents1));
253  edm::Ref<reco::PFJetCollection> subjetRef1(selectedSubjetRefProd, selectedSubjets->size() - 1);
254  selectedSubjets->push_back(convertToPFJet(*subjet2, subjetConstituents2));
255  edm::Ref<reco::PFJetCollection> subjetRef2(selectedSubjetRefProd, selectedSubjets->size() - 1);
256 
257  // find all PFCandidates that are not constituents of the **other** subjet
258  std::vector<reco::PFCandidateRef> pfCandidatesNotInSubjet1 = getPFCandidates_exclJetConstituents(*subjet1, pfCandidates, constitmap[2*idx], subjetConstituents2, 1.e-4, false);
259  std::vector<reco::PFCandidateRef> pfCandidatesNotInSubjet2 = getPFCandidates_exclJetConstituents(*subjet2, pfCandidates, constitmap[2*idx+1], subjetConstituents1, 1.e-4, false);
260  if ( verbosity_ >= 1 ) {
261  std::cout << "#pfCandidatesNotInSubjet1 = " << pfCandidatesNotInSubjet1.size() << std::endl;
262  std::cout << "#pfCandidatesNotInSubjet2 = " << pfCandidatesNotInSubjet2.size() << std::endl;
263  }
264 
265  // build JetToPFCandidateAssociation
266  // (key = subjet, value = collection of PFCandidates that are not constituents of subjet)
267  for(auto const& pfCandidate : pfCandidatesNotInSubjet1 ) {
268  selectedSubjetPFCandidateAssociationForIsolation->insert(subjetRef1, pfCandidate);
269  }
270  for(auto const& pfCandidate : pfCandidatesNotInSubjet2 ) {
271  selectedSubjetPFCandidateAssociationForIsolation->insert(subjetRef2, pfCandidate);
272  }
273  }
274 
275  evt.put(std::move(selectedSubjets));
276  evt.put(std::move(selectedSubjetPFCandidateAssociationForIsolation), "pfCandAssocMapForIsolation");
277 }
278 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
double eta() const final
momentum pseudorapidity
int mChargedMultiplicity
Definition: PFJet.h:75
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Base class for all types of Jets.
Definition: Jet.h:20
EDProductGetter const & productGetter() const
Definition: Event.cc:93
bool exists(std::string const &parameterName) const
checks if a parameter exists
std::vector< Constituent > Constituents
Definition: Jet.h:23
size_type size() const
double pt() const final
transverse momentum
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 &)
virtual void setJetArea(float fArea)
set jet area
Definition: Jet.h:103
char const * label
float mChargedEmEnergy
Definition: PFJet.h:72
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::View< reco::Jet > JetView
float mNeutralEmEnergy
Definition: PFJet.h:74
edm::EDGetTokenT< reco::PFCandidateCollection > srcPFCandidates_
const Point & vertex() const override
vertex position (overwritten by PF...)
double energy() const final
energy
const LorentzVector & p4() const final
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
RefProd< PROD > getRefBeforePut()
Definition: Event.h:150
edm::AssociationMap< edm::OneToMany< std::vector< reco::PFJet >, std::vector< reco::PFCandidate >, unsigned int > > JetToPFCandidateAssociation
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
int mNeutralMultiplicity
Definition: PFJet.h:76
std::vector< PFJet > PFJetCollection
collection of PFJet objects
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:40
virtual float jetArea() const
get jet area
Definition: Jet.h:105
virtual int nConstituents() const
of constituents
Definition: Jet.h:65
const_reference at(size_type pos) const
edm::EDGetTokenT< JetView > srcSubjets_
void produce(edm::Event &, const edm::EventSetup &) override
virtual ParticleType particleId() const
Definition: PFCandidate.h:374
double phi() const final
momentum azimuthal angle
def move(src, dest)
Definition: eostools.py:511
double mass() const final
mass