CMS 3D CMS Logo

PFRecoTauChargedHadronProducer.cc
Go to the documentation of this file.
1 /*
2  * PFRecoTauChargedHadronProducer
3  *
4  * Author: Christian Veelken, LLR
5  *
6  * Associates reconstructed ChargedHadrons to PFJets. The ChargedHadrons are built using one
7  * or more RecoTauBuilder plugins. Any overlaps (ChargedHadrons sharing tracks)
8  * are removed, with the best ChargedHadron candidates taken. The 'best' are defined
9  * via the input list of PFRecoTauChargedHadronQualityPlugins, which form a
10  * lexicograpical ranking.
11  *
12  */
13 
18 
21 
26 
37 
39 
40 #include <boost/ptr_container/ptr_vector.hpp>
41 #include <boost/ptr_container/ptr_list.hpp>
42 
43 #include <string>
44 #include <vector>
45 #include <list>
46 #include <set>
47 #include <algorithm>
48 #include <functional>
49 #include <cmath>
50 
52 {
53 public:
56 
59  void produce(edm::Event& evt, const edm::EventSetup& es) override;
60  template <typename T>
61  void print(const T& chargedHadrons);
62 
63  private:
64  typedef boost::ptr_vector<Builder> builderList;
65  typedef boost::ptr_vector<Ranker> rankerList;
66  typedef boost::ptr_vector<reco::PFRecoTauChargedHadron> ChargedHadronVector;
67  typedef boost::ptr_list<reco::PFRecoTauChargedHadron> ChargedHadronList;
68 
70 
72 
73  // input jet collection
76  double minJetPt_;
77  double maxJetAbsEta_;
78 
79  // plugins for building and ranking ChargedHadron candidates
80  builderList builders_;
81  rankerList rankers_;
82 
83  std::auto_ptr<ChargedHadronPredicate> predicate_;
84 
85  // output selector
86  std::auto_ptr<StringCutObjectSelector<reco::PFRecoTauChargedHadron> > outputSelector_;
87 
88  // flag to enable/disable debug print-out
90 };
91 
93  : moduleLabel_(cfg.getParameter<std::string>("@module_label"))
94 {
95  srcJets_ = cfg.getParameter<edm::InputTag>("jetSrc");
96  Jets_token = consumes<reco::CandidateView>(srcJets_);
97  minJetPt_ = ( cfg.exists("minJetPt") ) ? cfg.getParameter<double>("minJetPt") : -1.0;
98  maxJetAbsEta_ = ( cfg.exists("maxJetAbsEta") ) ? cfg.getParameter<double>("maxJetAbsEta") : 99.0;
99  verbosity_ = ( cfg.exists("verbosity") ) ?
100  cfg.getParameter<int>("verbosity") : 0;
101 
102  // get set of ChargedHadron builder plugins
103  edm::VParameterSet psets_builders = cfg.getParameter<edm::VParameterSet>("builders");
104  for ( edm::VParameterSet::const_iterator pset = psets_builders.begin();
105  pset != psets_builders.end(); ++pset ) {
106  std::string pluginType = pset->getParameter<std::string>("plugin");
107  edm::ParameterSet pset_modified = (*pset);
108  pset_modified.addParameter<int>("verbosity", verbosity_);
109  builders_.push_back(PFRecoTauChargedHadronBuilderPluginFactory::get()->create(pluginType, pset_modified, consumesCollector()));
110  }
111 
112  // get set of plugins for ranking ChargedHadrons in quality
113  edm::VParameterSet psets_rankers = cfg.getParameter<edm::VParameterSet>("ranking");
114  for ( edm::VParameterSet::const_iterator pset = psets_rankers.begin();
115  pset != psets_rankers.end(); ++pset ) {
116  std::string pluginType = pset->getParameter<std::string>("plugin");
117  edm::ParameterSet pset_modified = (*pset);
118  pset_modified.addParameter<int>("verbosity", verbosity_);
119  rankers_.push_back(PFRecoTauChargedHadronQualityPluginFactory::get()->create(pluginType, pset_modified));
120  }
121 
122  // build the sorting predicate
123  predicate_ = std::auto_ptr<ChargedHadronPredicate>(new ChargedHadronPredicate(rankers_));
124 
125  // check if we want to apply a final output selection
126  if ( cfg.exists("outputSelection") ) {
127  std::string selection = cfg.getParameter<std::string>("outputSelection");
128  if ( !selection.empty() ) {
130  }
131  }
132 
133  produces<reco::PFJetChargedHadronAssociation>();
134 }
135 
137 {
138  if ( verbosity_ ) {
139  edm::LogPrint("PFRecoTauChHProducer")<< "<PFRecoTauChargedHadronProducer::produce>:" ;
140  edm::LogPrint("PFRecoTauChHProducer")<< " moduleLabel = " << moduleLabel_ ;
141  }
142 
143  // give each of our plugins a chance at doing something with the edm::Event
144  for(auto& builder : builders_ ) {
145  builder.setup(evt, es);
146  }
147 
148  // get a view of our jets via the base candidates
150  evt.getByToken(Jets_token, jets);
151 
152  // convert the view to a RefVector of actual PFJets
153  reco::PFJetRefVector pfJets = reco::tau::castView<reco::PFJetRefVector>(jets);
154 
155  // make our association
156  std::unique_ptr<reco::PFJetChargedHadronAssociation> pfJetChargedHadronAssociations;
157 
158 
159  if ( !pfJets.empty() ) {
160  edm::Handle<reco::PFJetCollection> pfJetCollectionHandle;
161  evt.get(pfJets.id(), pfJetCollectionHandle);
162  pfJetChargedHadronAssociations = std::make_unique<reco::PFJetChargedHadronAssociation>(reco::PFJetRefProd(pfJetCollectionHandle));
163  } else {
164  pfJetChargedHadronAssociations = std::make_unique<reco::PFJetChargedHadronAssociation>();
165  }
166 
167  // loop over our jets
168  for(auto const& pfJet : pfJets ) {
169 
170  if(pfJet->pt() - minJetPt_ < 1e-5) continue;
171  if(std::abs(pfJet->eta()) - maxJetAbsEta_ > -1e-5) continue;
172 
173  // build global list of ChargedHadron candidates for each jet
174  ChargedHadronList uncleanedChargedHadrons;
175 
176  // merge candidates reconstructed by all desired algorithm plugins
177  for(auto const& builder : builders_ ) {
178  try {
179  ChargedHadronVector result(builder(*pfJet));
180  if ( verbosity_ ) {
181  edm::LogPrint("PFRecoTauChHProducer")<< "result of builder = " << builder.name() << ":" ;
182  print(result);
183  }
184  uncleanedChargedHadrons.transfer(uncleanedChargedHadrons.end(), result);
185  } catch ( cms::Exception& exception ) {
186  edm::LogError("BuilderPluginException")
187  << "Exception caught in builder plugin " << builder.name()
188  << ", rethrowing" << std::endl;
189  throw exception;
190  }
191  }
192 
193  // rank the candidates according to our quality plugins
194  uncleanedChargedHadrons.sort(*predicate_);
195 
196  // define collection of cleaned ChargedHadrons;
197  std::vector<reco::PFRecoTauChargedHadron> cleanedChargedHadrons;
198 
199  // keep track of neutral PFCandidates, charged PFCandidates and tracks "used" by ChargedHadron candidates in the clean collection
200  typedef std::pair<double, double> etaPhiPair;
201  std::list<etaPhiPair> tracksInCleanCollection;
202  std::set<reco::PFCandidatePtr> neutralPFCandsInCleanCollection;
203 
204  while ( !uncleanedChargedHadrons.empty() ) {
205 
206  // get next best ChargedHadron candidate
207  std::auto_ptr<reco::PFRecoTauChargedHadron> nextChargedHadron(uncleanedChargedHadrons.pop_front().release());
208  if ( verbosity_ ) {
209  edm::LogPrint("PFRecoTauChHProducer")<< "processing nextChargedHadron:" ;
210  edm::LogPrint("PFRecoTauChHProducer")<< (*nextChargedHadron);
211  }
212 
213  // discard candidates which fail final output selection
214  if ( !(*outputSelector_)(*nextChargedHadron) ) continue;
215 
216  const reco::Track* track = nullptr;
217  if ( nextChargedHadron->getChargedPFCandidate().isNonnull() ) {
218  const reco::PFCandidatePtr& chargedPFCand = nextChargedHadron->getChargedPFCandidate();
219  if ( chargedPFCand->trackRef().isNonnull() ) track = chargedPFCand->trackRef().get();
220  else if ( chargedPFCand->muonRef().isNonnull() && chargedPFCand->muonRef()->innerTrack().isNonnull() ) track = chargedPFCand->muonRef()->innerTrack().get();
221  else if ( chargedPFCand->muonRef().isNonnull() && chargedPFCand->muonRef()->globalTrack().isNonnull() ) track = chargedPFCand->muonRef()->globalTrack().get();
222  else if ( chargedPFCand->muonRef().isNonnull() && chargedPFCand->muonRef()->outerTrack().isNonnull() ) track = chargedPFCand->muonRef()->outerTrack().get();
223  else if ( chargedPFCand->gsfTrackRef().isNonnull() ) track = chargedPFCand->gsfTrackRef().get();
224  }
225  if ( nextChargedHadron->getTrack().isNonnull() && !track ) {
226  track = nextChargedHadron->getTrack().get();
227  }
228 
229  // discard candidate in case its track is "used" by any ChargedHadron in the clean collection
230  bool isTrack_overlap = false;
231  if ( track ) {
232  double track_eta = track->eta();
233  double track_phi = track->phi();
234  for ( std::list<etaPhiPair>::const_iterator trackInCleanCollection = tracksInCleanCollection.begin();
235  trackInCleanCollection != tracksInCleanCollection.end(); ++trackInCleanCollection ) {
236  double dR = deltaR(track_eta, track_phi, trackInCleanCollection->first, trackInCleanCollection->second);
237  if ( dR < 1.e-4 ) isTrack_overlap = true;
238  }
239  }
240  if ( verbosity_ ) {
241  edm::LogPrint("PFRecoTauChHProducer")<< "isTrack_overlap = " << isTrack_overlap ;
242  }
243  if ( isTrack_overlap ) continue;
244 
245  // discard ChargedHadron candidates without track in case they are close to neutral PFCandidates "used" by ChargedHadron candidates in the clean collection
246  bool isNeutralPFCand_overlap = false;
247  if ( nextChargedHadron->algoIs(reco::PFRecoTauChargedHadron::kPFNeutralHadron) ) {
248  for ( std::set<reco::PFCandidatePtr>::const_iterator neutralPFCandInCleanCollection = neutralPFCandsInCleanCollection.begin();
249  neutralPFCandInCleanCollection != neutralPFCandsInCleanCollection.end(); ++neutralPFCandInCleanCollection ) {
250  if ( (*neutralPFCandInCleanCollection) == nextChargedHadron->getChargedPFCandidate() ) isNeutralPFCand_overlap = true;
251  }
252  }
253  if ( verbosity_ ) {
254  edm::LogPrint("PFRecoTauChHProducer")<< "isNeutralPFCand_overlap = " << isNeutralPFCand_overlap ;
255  }
256  if ( isNeutralPFCand_overlap ) continue;
257 
258  // find neutral PFCandidates that are not "used" by any ChargedHadron in the clean collection
259  std::vector<reco::PFCandidatePtr> uniqueNeutralPFCands;
260  std::set_difference(nextChargedHadron->getNeutralPFCandidates().begin(),
261  nextChargedHadron->getNeutralPFCandidates().end(),
262  neutralPFCandsInCleanCollection.begin(),
263  neutralPFCandsInCleanCollection.end(),
264  std::back_inserter(uniqueNeutralPFCands));
265 
266  if ( uniqueNeutralPFCands.size() == nextChargedHadron->getNeutralPFCandidates().size() ) { // all neutral PFCandidates are unique, add ChargedHadron candidate to clean collection
267  if ( track ) tracksInCleanCollection.push_back(std::make_pair(track->eta(), track->phi()));
268  neutralPFCandsInCleanCollection.insert(nextChargedHadron->getNeutralPFCandidates().begin(), nextChargedHadron->getNeutralPFCandidates().end());
269  if ( verbosity_ ) {
270  edm::LogPrint("PFRecoTauChHProducer")<< "--> adding nextChargedHadron to output collection." ;
271  }
272  cleanedChargedHadrons.push_back(*nextChargedHadron);
273  } else { // remove overlapping neutral PFCandidates, reevaluate ranking criterion and process ChargedHadron candidate again
274  nextChargedHadron->neutralPFCandidates_.clear();
275  for(auto const& neutralPFCand : uniqueNeutralPFCands ) {
276  nextChargedHadron->neutralPFCandidates_.push_back(neutralPFCand);
277  }
278  // update ChargedHadron four-momentum
279  reco::tau::setChargedHadronP4(*nextChargedHadron);
280  // reinsert ChargedHadron candidate into list of uncleaned candidates,
281  // at position according to new rank
282  ChargedHadronList::iterator insertionPoint = std::lower_bound(uncleanedChargedHadrons.begin(), uncleanedChargedHadrons.end(), *nextChargedHadron, *predicate_);
283  if ( verbosity_ ) {
284  edm::LogPrint("PFRecoTauChHProducer")<< "--> removing non-unique neutral PFCandidates and reinserting nextChargedHadron in uncleaned collection." ;
285  }
286  uncleanedChargedHadrons.insert(insertionPoint, nextChargedHadron);
287  }
288  }
289 
290  if ( verbosity_ ) {
291  print(cleanedChargedHadrons);
292  }
293 
294  // add ChargedHadron-to-jet association
295  pfJetChargedHadronAssociations->setValue(pfJet.key(), cleanedChargedHadrons);
296  }
297 
298  evt.put(std::move(pfJetChargedHadronAssociations));
299 }
300 
301 template <typename T>
302 void PFRecoTauChargedHadronProducer::print(const T& chargedHadrons)
303 {
304  for ( typename T::const_iterator chargedHadron = chargedHadrons.begin();
305  chargedHadron != chargedHadrons.end(); ++chargedHadron ) {
306  edm::LogPrint("PFRecoTauChHProducer") << (*chargedHadron);
307  edm::LogPrint("PFRecoTauChHProducer") << "Rankers:" ;
308  for ( rankerList::const_iterator ranker = rankers_.begin();
309  ranker != rankers_.end(); ++ranker) {
310  const unsigned width = 25;
311  edm::LogPrint("PFRecoTauChHProducer") << " " << std::setiosflags(std::ios::left) << std::setw(width) << ranker->name()
312  << " " << std::resetiosflags(std::ios::left) << std::setprecision(3) << (*ranker)(*chargedHadron) << std::endl;
313  }
314  }
315 }
316 
318 
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
def create(alignables, pedeDump, additionalData, outputFile, config)
reco::tau::RecoTauLexicographicalRanking< rankerList, reco::PFRecoTauChargedHadron > ChargedHadronPredicate
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
PFRecoTauChargedHadronProducer(const edm::ParameterSet &cfg)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< ParameterSet > VParameterSet
Definition: ParameterSet.h:33
T const * get() const
Returns C++ pointer to the item.
Definition: Ptr.h:159
selection
main part
Definition: corrVsCorr.py:99
bool exists(std::string const &parameterName) const
checks if a parameter exists
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:678
void produce(edm::Event &evt, const edm::EventSetup &es) override
bool empty() const
Is the RefVector empty.
Definition: RefVector.h:104
edm::EDGetTokenT< reco::CandidateView > Jets_token
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:684
reco::tau::PFRecoTauChargedHadronBuilderPlugin Builder
ProductID id() const
Accessor for product ID.
Definition: RefVector.h:122
boost::ptr_vector< reco::PFRecoTauChargedHadron > ChargedHadronVector
std::auto_ptr< ChargedHadronPredicate > predicate_
vector< PseudoJet > jets
void addParameter(std::string const &name, T const &value)
Definition: ParameterSet.h:143
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool get(ProductID const &oid, Handle< PROD > &result) const
Definition: Event.h:356
edm::RefProd< PFJetCollection > PFJetRefProd
bool isNonnull() const
Checks for non-null.
Definition: Ptr.h:168
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
reco::tau::PFRecoTauChargedHadronQualityPlugin Ranker
std::auto_ptr< StringCutObjectSelector< reco::PFRecoTauChargedHadron > > outputSelector_
boost::ptr_list< reco::PFRecoTauChargedHadron > ChargedHadronList
void setChargedHadronP4(reco::PFRecoTauChargedHadron &chargedHadron, double scaleFactor_neutralPFCands=1.0)
long double T
def move(src, dest)
Definition: eostools.py:511
T get(const Candidate &c)
Definition: component.h:55