CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauJetRegionProducer.cc
Go to the documentation of this file.
1 /*
2  * RecoTauJetRegionProducer
3  *
4  * Given a set of PFJets, make new jets with the same p4 but collect all the
5  * PFCandidates from a cone of a given size into the constituents.
6  *
7  * Author: Evan K. Friis, UC Davis
8  *
9  */
10 
11 #include <boost/bind.hpp>
12 
18 
21 
28 
29 #include <string>
30 #include <iostream>
31 
33 {
34  public:
36  typedef edm::AssociationMap<edm::OneToMany<std::vector<reco::PFJet>, std::vector<reco::PFCandidate>, unsigned int> > JetToPFCandidateAssociation;
37  explicit RecoTauJetRegionProducer(const edm::ParameterSet& pset);
39 
40  void produce(edm::Event& evt, const edm::EventSetup& es) override;
41 
42  private:
44 
48 
52 
53  double deltaR2_;
54 
56 };
57 
59  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
60 {
62  pfCandSrc_ = cfg.getParameter<edm::InputTag>("pfCandSrc");
63  pfCandAssocMapSrc_ = cfg.getParameter<edm::InputTag>("pfCandAssocMapSrc");
64 
65  pf_token = consumes<reco::PFCandidateCollection>(pfCandSrc_);
66  Jets_token = consumes<reco::CandidateView>(inputJets_);
67  pfCandAssocMap_token = consumes<JetToPFCandidateAssociation>(pfCandAssocMapSrc_);
68 
69  double deltaR = cfg.getParameter<double>("deltaR");
70  deltaR2_ = deltaR*deltaR;
71 
72  verbosity_ = ( cfg.exists("verbosity") ) ?
73  cfg.getParameter<int>("verbosity") : 0;
74 
75  produces<reco::PFJetCollection>("jets");
76  produces<PFJetMatchMap>();
77 }
78 
80 {
81  if ( verbosity_ ) {
82  std::cout << "<RecoTauJetRegionProducer::produce (moduleLabel = " << moduleLabel_ << ")>:" << std::endl;
83  std::cout << " inputJets = " << inputJets_ << std::endl;
84  std::cout << " pfCandSrc = " << pfCandSrc_ << std::endl;
85  std::cout << " pfCandAssocMapSrc_ = " << pfCandAssocMapSrc_ << std::endl;
86  }
87 
89  evt.getByToken(pf_token, pfCandsHandle);
90 
91  // Build Ptrs for all the PFCandidates
92  typedef edm::Ptr<reco::PFCandidate> PFCandPtr;
93  std::vector<PFCandPtr> pfCands;
94  pfCands.reserve(pfCandsHandle->size());
95  for ( size_t icand = 0; icand < pfCandsHandle->size(); ++icand ) {
96  pfCands.push_back(PFCandPtr(pfCandsHandle, icand));
97  }
98 
99  // Get the jets
101  evt.getByToken(Jets_token, jetView);
102  // Convert to a vector of PFJetRefs
103  reco::PFJetRefVector jets = reco::tau::castView<reco::PFJetRefVector>(jetView);
104  size_t nJets = jets.size();
105 
106  // Get the association map matching jets to PFCandidates
107  // (needed for recinstruction of boosted taus)
109  if ( pfCandAssocMapSrc_.label() != "" ) {
110  evt.getByToken(pfCandAssocMap_token, jetToPFCandMap);
111  }
112 
113  // Get the original product, so we can match against it - otherwise the
114  // indices don't match up.
115  edm::ProductID originalId = jets.id();
117  size_t nOriginalJets = 0;
118  // We have to make sure that we have some selected jets, otherwise we don't
119  // actually have a valid product ID to the original jets.
120  if ( nJets ) {
121  try {
122  evt.get(originalId, originalJets);
123  } catch(const cms::Exception &e) {
124  edm::LogError("MissingOriginalCollection")
125  << "Can't get the original jets that made: " << inputJets_
126  << " that have product ID: " << originalId
127  << " from the event!!";
128  throw e;
129  }
130  nOriginalJets = originalJets->size();
131  }
132 
133  std::auto_ptr<reco::PFJetCollection> newJets(new reco::PFJetCollection);
134 
135  // Keep track of the indices of the current jet and the old (original) jet
136  // -1 indicates no match.
137  std::vector<int> matchInfo(nOriginalJets, -1);
138  newJets->reserve(nJets);
139  for ( size_t ijet = 0; ijet < nJets; ++ijet ) {
140  // Get a ref to jet
141  reco::PFJetRef jetRef = jets[ijet];
142  // Make an initial copy.
143  newJets->emplace_back(*jetRef);
144  reco::PFJet& newJet = newJets->back();
145  // Clear out all the constituents
146  newJet.clearDaughters();
147  // Loop over all the PFCands
148  for ( std::vector<PFCandPtr>::const_iterator pfCand = pfCands.begin();
149  pfCand != pfCands.end(); ++pfCand ) {
150  bool isMappedToJet = false;
151  if ( pfCandAssocMapSrc_.label() != "" ) {
152  edm::RefVector<reco::PFCandidateCollection> pfCandsMappedToJet = (*jetToPFCandMap)[jetRef];
153  for ( edm::RefVector<reco::PFCandidateCollection>::const_iterator pfCandMappedToJet = pfCandsMappedToJet.begin();
154  pfCandMappedToJet != pfCandsMappedToJet.end(); ++pfCandMappedToJet ) {
155  if ( reco::deltaR2(**pfCandMappedToJet, **pfCand) < 1.e-8 ) {
156  isMappedToJet = true;
157  break;
158  }
159  }
160  } else {
161  isMappedToJet = true;
162  }
163  if ( reco::deltaR2(*jetRef, **pfCand) < deltaR2_ && isMappedToJet ) newJet.addDaughter(*pfCand);
164  }
165  if ( verbosity_ ) {
166  std::cout << "jet #" << ijet << ": Pt = " << jetRef->pt() << ", eta = " << jetRef->eta() << ", phi = " << jetRef->eta() << ","
167  << " mass = " << jetRef->mass() << ", area = " << jetRef->jetArea() << std::endl;
168  std::vector<reco::PFCandidatePtr> jetConstituents = newJet.getPFConstituents();
169  int idx = 0;
170  for ( std::vector<reco::PFCandidatePtr>::const_iterator jetConstituent = jetConstituents.begin();
171  jetConstituent != jetConstituents.end(); ++jetConstituent ) {
172  std::cout << " constituent #" << idx << ": Pt = " << (*jetConstituent)->pt() << ", eta = " << (*jetConstituent)->eta() << ", phi = " << (*jetConstituent)->phi() << std::endl;
173  }
174  }
175  // Match the index of the jet we just made to the index into the original
176  // collection.
177  matchInfo[jetRef.key()] = ijet;
178  }
179 
180  // Put our new jets into the event
181  edm::OrphanHandle<reco::PFJetCollection> newJetsInEvent = evt.put(newJets, "jets");
182 
183  // Create a matching between original jets -> extra collection
184  std::auto_ptr<PFJetMatchMap> matching(new PFJetMatchMap(newJetsInEvent));
185  if ( nJets ) {
186  PFJetMatchMap::Filler filler(*matching);
187  filler.insert(originalJets, matchInfo.begin(), matchInfo.end());
188  filler.fill();
189  }
190  evt.put(matching);
191 }
192 
T getParameter(std::string const &) const
tuple cfg
Definition: looper.py:259
RecoTauJetRegionProducer(const edm::ParameterSet &pset)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
edm::Association< reco::PFJetCollection > PFJetMatchMap
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< reco::PFCandidateCollection > pf_token
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
bool exists(std::string const &parameterName) const
checks if a parameter exists
key_type key() const
Accessor for product key.
Definition: Ref.h:266
const_iterator end() const
Termination of iteration.
Definition: RefVector.h:249
void produce(edm::Event &evt, const edm::EventSetup &es) override
Jets made from PFObjects.
Definition: PFJet.h:21
const_iterator begin() const
Initialize an iterator over the RefVector.
Definition: RefVector.h:244
ProductID id() const
Accessor for product ID.
Definition: RefVector.h:104
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
edm::EDGetTokenT< reco::CandidateView > Jets_token
vector< PseudoJet > jets
edm::EDGetTokenT< JetToPFCandidateAssociation > pfCandAssocMap_token
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:318
double deltaR2(const T1 &t1, const T2 &t2)
Definition: deltaR.h:36
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
void clearDaughters()
clear daughter references
std::string const & label() const
Definition: InputTag.h:42
std::vector< PFJet > PFJetCollection
collection of PFJet objects
edm::AssociationMap< edm::OneToMany< std::vector< reco::PFJet >, std::vector< reco::PFCandidate >, unsigned int > > JetToPFCandidateAssociation
void addDaughter(const CandidatePtr &)
add a daughter via a reference
virtual std::vector< reco::PFCandidatePtr > getPFConstituents() const
get all constituents
Definition: PFJet.cc:52
size_type size() const
Size of the RefVector.
Definition: RefVector.h:89
tuple cout
Definition: gather_cfg.py:121
moduleLabel_(iConfig.getParameter< string >("@module_label"))