CMS 3D CMS Logo

RecoTauGenericJetRegionProducer.cc
Go to the documentation of this file.
1 /*
2  * RecoTauJetRegionProducer
3  *
4  * Given a set of Jets, make new jets with the same p4 but collect all the
5  * particle-flow candidates from a cone of a given size into the constituents.
6  *
7  * Author: Evan K. Friis, UC Davis
8  *
9  */
10 
18 
21 
28 
31 
32 #include <string>
33 #include <iostream>
34 
35 template <class JetType, class CandType>
37 public:
39  typedef edm::AssociationMap<edm::OneToMany<std::vector<JetType>, std::vector<CandType>, unsigned int> >
43 
44  void produce(edm::Event& evt, const edm::EventSetup& es) override;
45 
46  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
47  static void fillDescriptionsBase(edm::ConfigurationDescriptions& descriptions, const std::string& name);
48 
49 private:
51 
55 
59 
60  double minJetPt_;
61  double maxJetAbsEta_;
62  double deltaR2_;
63 
65 };
66 
67 template <class JetType, class CandType>
69  : moduleLabel_(cfg.getParameter<std::string>("@module_label")) {
70  inputJets_ = cfg.getParameter<edm::InputTag>("src");
71  pfCandSrc_ = cfg.getParameter<edm::InputTag>("pfCandSrc");
72  pfCandAssocMapSrc_ = cfg.getParameter<edm::InputTag>("pfCandAssocMapSrc");
73 
74  pf_token = consumes<std::vector<CandType> >(pfCandSrc_);
75  Jets_token = consumes<reco::CandidateView>(inputJets_);
76  pfCandAssocMap_token = consumes<JetToCandidateAssociation>(pfCandAssocMapSrc_);
77 
78  double deltaR = cfg.getParameter<double>("deltaR");
80  minJetPt_ = cfg.getParameter<double>("minJetPt");
81  maxJetAbsEta_ = cfg.getParameter<double>("maxJetAbsEta");
82 
83  verbosity_ = cfg.getParameter<int>("verbosity");
84 
85  produces<std::vector<JetType> >("jets");
86  produces<JetMatchMap>();
87 }
88 
89 template <class JetType, class CandType>
91  if (verbosity_) {
92  std::cout << "<RecoTauJetRegionProducer::produce (moduleLabel = " << moduleLabel_ << ")>:" << std::endl;
93  std::cout << " inputJets = " << inputJets_ << std::endl;
94  std::cout << " pfCandSrc = " << pfCandSrc_ << std::endl;
95  std::cout << " pfCandAssocMapSrc_ = " << pfCandAssocMapSrc_ << std::endl;
96  }
97 
98  edm::Handle<std::vector<CandType> > pfCandsHandle;
99  evt.getByToken(pf_token, pfCandsHandle);
100 
101  // Build Ptrs for all the PFCandidates
102  typedef edm::Ptr<CandType> CandPtr;
103  std::vector<CandPtr> pfCands;
104  pfCands.reserve(pfCandsHandle->size());
105  for (size_t icand = 0; icand < pfCandsHandle->size(); ++icand) {
106  pfCands.push_back(CandPtr(pfCandsHandle, icand));
107  }
108 
109  // Get the jets
111  evt.getByToken(Jets_token, jetView);
112  // Convert to a vector of JetRefs
113  edm::RefVector<std::vector<JetType> > jets = reco::tau::castView<edm::RefVector<std::vector<JetType> > >(jetView);
114  size_t nJets = jets.size();
115 
116  // Get the association map matching jets to Candidates
117  // (needed for reconstruction of boosted taus)
119  std::vector<std::unordered_set<unsigned> > fastJetToPFCandMap;
120  if (!pfCandAssocMapSrc_.label().empty()) {
121  evt.getByToken(pfCandAssocMap_token, jetToPFCandMap);
122  fastJetToPFCandMap.resize(nJets);
123  for (size_t ijet = 0; ijet < nJets; ++ijet) {
124  // Get a ref to jet
125  const edm::Ref<std::vector<JetType> >& jetRef = jets[ijet];
126  const auto& pfCandsMappedToJet = (*jetToPFCandMap)[jetRef];
127  for (const auto& pfCandMappedToJet : pfCandsMappedToJet) {
128  fastJetToPFCandMap[ijet].emplace(pfCandMappedToJet.key());
129  }
130  }
131  }
132 
133  // Get the original product, so we can match against it - otherwise the
134  // indices don't match up.
135  edm::ProductID originalId = jets.id();
136  edm::Handle<std::vector<JetType> > originalJets;
137  size_t nOriginalJets = 0;
138  // We have to make sure that we have some selected jets, otherwise we don't
139  // actually have a valid product ID to the original jets.
140  if (nJets) {
141  try {
142  evt.get(originalId, originalJets);
143  } catch (const cms::Exception& e) {
144  edm::LogError("MissingOriginalCollection") << "Can't get the original jets that made: " << inputJets_
145  << " that have product ID: " << originalId << " from the event!!";
146  throw e;
147  }
148  nOriginalJets = originalJets->size();
149  }
150 
151  auto newJets = std::make_unique<std::vector<JetType> >();
152 
153  // Keep track of the indices of the current jet and the old (original) jet
154  // -1 indicates no match.
155  std::vector<int> matchInfo(nOriginalJets, -1);
156  newJets->reserve(nJets);
157  size_t nNewJets = 0;
158  for (size_t ijet = 0; ijet < nJets; ++ijet) {
159  // Get a ref to jet
160  const edm::Ref<std::vector<JetType> >& jetRef = jets[ijet];
161  if (jetRef->pt() - minJetPt_ < 1e-5)
162  continue;
163  if (std::abs(jetRef->eta()) - maxJetAbsEta_ > -1e-5)
164  continue;
165  // Make an initial copy.
166  newJets->emplace_back(*jetRef);
167  JetType& newJet = newJets->back();
168  // Clear out all the constituents
169  newJet.clearDaughters();
170  // Loop over all the PFCands
171  for (const auto& pfCand : pfCands) {
172  bool isMappedToJet = false;
173  if (jetToPFCandMap.isValid()) {
174  auto temp = jetToPFCandMap->find(jetRef);
175  if (temp == jetToPFCandMap->end()) {
176  edm::LogWarning("WeirdCandidateMap") << "Candidate map for jet " << jetRef.key() << " is empty!";
177  continue;
178  }
179  isMappedToJet = fastJetToPFCandMap[ijet].count(pfCand.key());
180  } else {
181  isMappedToJet = true;
182  }
183  if (reco::deltaR2(*jetRef, *pfCand) < deltaR2_ && isMappedToJet)
184  newJet.addDaughter(pfCand);
185  }
186  if (verbosity_) {
187  std::cout << "jet #" << ijet << ": Pt = " << jetRef->pt() << ", eta = " << jetRef->eta()
188  << ", phi = " << jetRef->eta() << ","
189  << " mass = " << jetRef->mass() << ", area = " << jetRef->jetArea() << std::endl;
190  auto jetConstituents = newJet.daughterPtrVector();
191  int idx = 0;
192  for (const auto& jetConstituent : jetConstituents) {
193  std::cout << " constituent #" << idx << ": Pt = " << jetConstituent->pt() << ", eta = " << jetConstituent->eta()
194  << ", phi = " << jetConstituent->phi() << std::endl;
195  ++idx;
196  }
197  }
198  // Match the index of the jet we just made to the index into the original
199  // collection.
200  //matchInfo[jetRef.key()] = ijet;
201  matchInfo[jetRef.key()] = nNewJets;
202  nNewJets++;
203  }
204 
205  // Put our new jets into the event
206  edm::OrphanHandle<std::vector<JetType> > newJetsInEvent = evt.put(std::move(newJets), "jets");
207 
208  // Create a matching between original jets -> extra collection
209  auto matching = (nJets != 0)
210  ? std::make_unique<JetMatchMap>(
212  : std::make_unique<JetMatchMap>();
213  for (size_t ijet = 0; ijet < nJets; ++ijet) {
214  matching->insert(edm::RefToBase<reco::Jet>(jets[ijet]),
215  edm::RefToBase<reco::Jet>(edm::Ref<std::vector<JetType> >(newJetsInEvent, matchInfo[ijet])));
216  }
217  evt.put(std::move(matching));
218 }
219 
220 template <class JetType, class CandType>
222  edm::ConfigurationDescriptions& descriptions, const std::string& name) {
223  // RecoTauGenericJetRegionProducer
225  desc.add<edm::InputTag>("src", edm::InputTag("ak4PFJets"));
226  desc.add<double>("deltaR", 0.8);
227  desc.add<edm::InputTag>("pfCandAssocMapSrc", edm::InputTag(""));
228  desc.add<int>("verbosity", 0);
229  desc.add<double>("maxJetAbsEta", 2.5);
230  desc.add<double>("minJetPt", 14.0);
231  desc.add<edm::InputTag>("pfCandSrc", edm::InputTag("particleFlow"));
232  descriptions.add(name, desc);
233 }
234 
235 template <>
237  edm::ConfigurationDescriptions& descriptions) {
238  // RecoTauGenericJetRegionProducer
239  RecoTauGenericJetRegionProducer::fillDescriptionsBase(descriptions, "RecoTauJetRegionProducer");
240 }
241 
242 template <>
244  edm::ConfigurationDescriptions& descriptions) {
245  // RecoTauGenericJetRegionProducer
246  RecoTauGenericJetRegionProducer::fillDescriptionsBase(descriptions, "RecoTauPatJetRegionProducer");
247 }
248 
251 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:344
static void fillDescriptionsBase(edm::ConfigurationDescriptions &descriptions, const std::string &name)
edm::EDGetTokenT< JetToCandidateAssociation > pfCandAssocMap_token
RecoTauGenericJetRegionProducer< reco::PFJet, reco::PFCandidate > RecoTauJetRegionProducer
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
RecoTauGenericJetRegionProducer< pat::Jet, pat::PackedCandidate > RecoTauPatJetRegionProducer
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
void produce(edm::Event &evt, const edm::EventSetup &es) override
Log< level::Error, false > LogError
static constexpr int nJets
key_type key() const
Accessor for product key.
Definition: Ref.h:250
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
edm::EDGetTokenT< reco::CandidateView > Jets_token
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::AssociationMap< edm::OneToMany< std::vector< JetType >, std::vector< CandType >, unsigned int > > JetToCandidateAssociation
constexpr auto deltaR2(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:16
edm::AssociationMap< edm::OneToOne< reco::JetView, reco::JetView > > JetMatchMap
void add(std::string const &label, ParameterSetDescription const &psetDescription)
bool isValid() const
Definition: HandleBase.h:70
RecoTauGenericJetRegionProducer(const edm::ParameterSet &pset)
edm::EDGetTokenT< std::vector< CandType > > pf_token
RefToBaseProd< T > makeRefToBaseProdFrom(RefToBase< T > const &iRef, Event const &iEvent)
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511