CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
RecoTauBuilderConePlugin.cc
Go to the documentation of this file.
1 /*
2  * RecoTauBuilderConePlugin
3  *
4  * Build a PFTau using cones defined in DeltaR.
5  *
6  * Original Authors: Ludovic Houchu, Simone Gennai
7  * Modifications: Evan K. Friis
8  *
9  */
10 
11 #include <vector>
12 #include <algorithm>
13 
18 
21 
27 
29 
30 namespace reco { namespace tau {
31 
32 typedef std::vector<RecoTauPiZero> PiZeroList;
33 
35  public:
38  // Build a tau from a jet
40  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons,
41  const std::vector<RecoTauPiZero>& piZeros,
42  const std::vector<PFCandidatePtr>& regionalExtras) const override;
43  private:
45 
47 
49 
50  /* String function to extract values from PFJets */
52 
53  // Cone defintions
61 };
62 
63 // ctor - initialize all of our variables
65  const edm::ParameterSet& pset, edm::ConsumesCollector &&iC):RecoTauBuilderPlugin(pset,std::move(iC)),
66  qcuts_(pset.getParameterSet(
67  "qualityCuts").getParameterSet("signalQualityCuts")),
68  usePFLeptonsAsChargedHadrons_(pset.getParameter<bool>("usePFLeptons")),
69  leadObjecPtThreshold_(pset.getParameter<double>("leadObjectPt")),
70  matchingCone_(pset.getParameter<std::string>("matchingCone")),
71  signalConeChargedHadrons_(pset.getParameter<std::string>(
72  "signalConeChargedHadrons")),
73  isoConeChargedHadrons_(
74  pset.getParameter<std::string>("isoConeChargedHadrons")),
75  signalConePiZeros_(
76  pset.getParameter<std::string>("signalConePiZeros")),
77  isoConePiZeros_(
78  pset.getParameter<std::string>("isoConePiZeros")),
79  signalConeNeutralHadrons_(
80  pset.getParameter<std::string>("signalConeNeutralHadrons")),
81  isoConeNeutralHadrons_(
82  pset.getParameter<std::string>("isoConeNeutralHadrons"))
83 {}
84 
85 namespace xclean
86 {
87  // define template specialization for cross-cleaning
88  template<>
90  {
91  // Get the list of objects we need to clean
92  for ( cone::PFCandPtrDRFilterIter signalTrack = signalTracksBegin; signalTrack != signalTracksEnd; ++signalTrack ) {
93  toRemove_.insert(reco::CandidatePtr(*signalTrack));
94  }
95  }
96 
97  template<>
98  inline void CrossCleanPtrs<PiZeroList::const_iterator>::initialize(const PiZeroList::const_iterator& piZerosBegin, const PiZeroList::const_iterator& piZerosEnd)
99  {
100  BOOST_FOREACH( const PFCandidatePtr &ptr, flattenPiZeros(piZerosBegin, piZerosEnd) ) {
101  toRemove_.insert(CandidatePtr(ptr));
102  }
103  }
104 }
105 
107  const reco::PFJetRef& jet,
108  const std::vector<reco::PFRecoTauChargedHadron>& chargedHadrons,
109  const std::vector<RecoTauPiZero>& piZeros,
110  const std::vector<PFCandidatePtr>& regionalExtras) const {
111  // Get access to our cone tools
112  using namespace cone;
113  // Define output. We only produce one tau per jet for the cone algo.
115 
116  // Our tau builder - the true indicates to automatically copy gamma candidates
117  // from the pizeros.
118  RecoTauConstructor tau(jet, getPFCands(), true);
119  // Setup our quality cuts to use the current vertex (supplied by base class)
120  qcuts_.setPV(primaryVertex(jet));
121 
122  typedef std::vector<PFCandidatePtr> PFCandPtrs;
123 
124  // Get the PF Charged hadrons + quality cuts
125  PFCandPtrs pfchs;
128  } else {
129  // Check if we want to include electrons in muons in "charged hadron"
130  // collection. This is the preferred behavior, as the PF lepton selections
131  // are very loose.
132  pfchs = qcuts_.filterCandRefs(pfChargedCands(*jet));
133  }
134 
135  // Get the PF gammas
138  // Neutral hadrons
141 
142  // All the extra junk
143  PFCandPtrs regionalJunk = qcuts_.filterCandRefs(regionalExtras);
144 
145  /***********************************************
146  ****** Lead Candidate Finding **********
147  ***********************************************/
148 
149  // Define our matching cone and filters
150  double matchingCone = matchingCone_(*jet);
151  PFCandPtrDRFilter matchingConeFilter(jet->p4(), 0, matchingCone);
152 
153  // Find the maximum PFCharged hadron in the matching cone. The call to
154  // PFCandidates always a sorted list, so we can just take the first if it
155  // if it exists.
156  PFCandidatePtr leadPFCH;
157  PFCandPtrs::iterator leadPFCH_iter =
158  std::find_if(pfchs.begin(), pfchs.end(), matchingConeFilter);
159 
160  if (leadPFCH_iter != pfchs.end()) {
161  leadPFCH = *leadPFCH_iter;
162  // Set leading candidate
163  tau.setleadPFChargedHadrCand(leadPFCH);
164  } else {
165  // If there is no leading charged candidate at all, return nothing - the
166  // producer class that owns the plugin will build a null tau if desired.
167  return output.release();
168  }
169 
170  // Find the leading neutral candidate
171  PFCandidatePtr leadPFGamma;
172  PFCandPtrs::iterator leadPFGamma_iter =
173  std::find_if(pfGammas.begin(), pfGammas.end(), matchingConeFilter);
174 
175  if (leadPFGamma_iter != pfGammas.end()) {
176  leadPFGamma = *leadPFGamma_iter;
177  // Set leading neutral candidate
178  tau.setleadPFNeutralCand(leadPFGamma);
179  }
180 
182  // Always use the leadPFCH if it is above our threshold
183  if (leadPFCH.isNonnull() && leadPFCH->pt() > leadObjecPtThreshold_) {
184  leadPFCand = leadPFCH;
185  } else if (leadPFGamma.isNonnull() &&
186  leadPFGamma->pt() > leadObjecPtThreshold_) {
187  // Otherwise use the lead Gamma if it is above threshold
188  leadPFCand = leadPFGamma;
189  } else {
190  // If both are too low PT, just take the charged one
191  leadPFCand = leadPFCH;
192  }
193 
195 
196  // Our cone axis is defined about the lead charged hadron
197  reco::Candidate::LorentzVector coneAxis = leadPFCH->p4();
198 
199  /***********************************************
200  ****** Cone Construction **********
201  ***********************************************/
202 
203  // Define the signal and isolation cone sizes for this jet and build filters
204  // to select elements in the given DeltaR regions
205 
206  PFCandPtrDRFilter signalConePFCHFilter(
207  coneAxis, -0.1, signalConeChargedHadrons_(*jet));
208  PFCandPtrDRFilter signalConePFNHFilter(
209  coneAxis, -0.1, signalConeNeutralHadrons_(*jet));
210  PiZeroDRFilter signalConePiZeroFilter(
211  coneAxis, -0.1, signalConePiZeros_(*jet));
212 
213  PFCandPtrDRFilter isoConePFCHFilter(
214  coneAxis, signalConeChargedHadrons_(*jet), isoConeChargedHadrons_(*jet));
215  PFCandPtrDRFilter isoConePFGammaFilter(
216  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
217  PFCandPtrDRFilter isoConePFNHFilter(
218  coneAxis, signalConeNeutralHadrons_(*jet), isoConeNeutralHadrons_(*jet));
219  PiZeroDRFilter isoConePiZeroFilter(
220  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
221 
222  // Additionally make predicates to select the different PF object types
223  // of the regional junk objects to add to the iso cone.
225 
229 
230  // Predicate to select the regional junk in the iso cone by PF id
231  RegionalJunkConeAndIdFilter pfChargedJunk(
232  pfchCandSelector, // select charged stuff from junk
233  isoConePFCHFilter // only take those in iso cone
234  );
235 
236  RegionalJunkConeAndIdFilter pfGammaJunk(
237  pfgammaCandSelector, // select gammas from junk
238  isoConePFGammaFilter // only take those in iso cone
239  );
240 
241  RegionalJunkConeAndIdFilter pfNeutralJunk(
242  pfnhCandSelector, // select neutral stuff from junk
243  isoConePFNHFilter // select stuff in iso cone
244  );
245 
246  // Build filter iterators select the signal charged stuff.
247  PFCandPtrDRFilterIter signalPFCHs_begin(
248  signalConePFCHFilter, pfchs.begin(), pfchs.end());
249  PFCandPtrDRFilterIter signalPFCHs_end(
250  signalConePFCHFilter, pfchs.end(), pfchs.end());
251 
252  // Cross clean pi zero content using signal cone charged hadron constituents.
254  signalPFCHs_begin, signalPFCHs_end);
255  std::vector<reco::RecoTauPiZero> cleanPiZeros = piZeroXCleaner(piZeros);
256 
257  // For the rest of the constituents, we need to filter any constituents that
258  // are already contained in the pizeros (i.e. electrons)
259  xclean::CrossCleanPtrs<PiZeroList::const_iterator> pfCandXCleaner(cleanPiZeros.begin(), cleanPiZeros.end());
260 
261  // Build signal charged hadrons
264  signalPFCHs_begin, signalPFCHs_end);
265 
268  boost::make_filter_iterator(
269  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner),
270  pfnhs.begin(), pfnhs.end()),
271  boost::make_filter_iterator(
272  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner),
273  pfnhs.end(), pfnhs.end()));
274 
275  // Build signal PiZeros
277  PiZeroDRFilterIter(signalConePiZeroFilter,
278  cleanPiZeros.begin(), cleanPiZeros.end()),
279  PiZeroDRFilterIter(signalConePiZeroFilter,
280  cleanPiZeros.end(), cleanPiZeros.end()));
281 
282  // Build isolation charged hadrons
285  boost::make_filter_iterator(
286  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner),
287  pfchs.begin(), pfchs.end()),
288  boost::make_filter_iterator(
289  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner),
290  pfchs.end(), pfchs.end()));
291 
292  // Add all the stuff in the isolation cone that wasn't in the jet constituents
295  boost::make_filter_iterator(
296  pfChargedJunk, regionalJunk.begin(), regionalJunk.end()),
297  boost::make_filter_iterator(
298  pfChargedJunk, regionalJunk.end(), regionalJunk.end())
299  );
300 
301  // Build isolation neutral hadrons
304  boost::make_filter_iterator(
305  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner),
306  pfnhs.begin(), pfnhs.end()),
307  boost::make_filter_iterator(
308  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner),
309  pfnhs.end(), pfnhs.end()));
310 
311  // Add regional stuff not in jet
314  boost::make_filter_iterator(
315  pfNeutralJunk, regionalJunk.begin(), regionalJunk.end()),
316  boost::make_filter_iterator(
317  pfNeutralJunk, regionalJunk.end(), regionalJunk.end())
318  );
319 
320  // Build isolation PiZeros
322  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.begin(),
323  cleanPiZeros.end()),
324  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.end(),
325  cleanPiZeros.end()));
326 
327  // Add regional stuff not in jet
330  boost::make_filter_iterator(
331  pfGammaJunk, regionalJunk.begin(), regionalJunk.end()),
332  boost::make_filter_iterator(
333  pfGammaJunk, regionalJunk.end(), regionalJunk.end())
334  );
335 
336  // Put our built tau in the output - 'false' indicates don't build the
337  // leading candidates, we already did that explicitly above.
338 
339  std::auto_ptr<reco::PFTau> tauPtr = tau.get(false);
340 
341  // Set event vertex position for tau
342  reco::VertexRef primaryVertexRef = primaryVertex(jet);
343  if ( primaryVertexRef.isNonnull() )
344  tauPtr->setVertex(primaryVertexRef->position());
345 
346  output.push_back(tauPtr);
347  return output.release();
348 }
349 }} // end namespace reco::tauk
350 
354  "RecoTauBuilderConePlugin");
const edm::Handle< PFCandidateCollection > & getPFCands() const
Hack to be able to convert Ptrs to Refs.
std::auto_ptr< reco::PFTau > get(bool setupLeadingCandidates=true)
Coll filterCandRefs(const Coll &refcoll, bool invert=false) const
Filter a ref vector of PFCandidates.
std::vector< PFCandidatePtr > pfGammas(const PFJet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
ParameterSet const & getParameterSet(ParameterSetID const &id)
std::vector< reco::PFCandidatePtr > PFCandPtrs
boost::filter_iterator< PiZeroDRFilter, std::vector< RecoTauPiZero >::const_iterator > PiZeroDRFilterIter
Definition: ConeTools.h:58
boost::ptr_vector< reco::PFTau > output_type
std::auto_ptr< output_type > return_type
RecoTauBuilderConePlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
std::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
void addPiZeros(Region region, const InputIterator &begin, const InputIterator &end)
Add a list of pizeros to the input collection.
reco::VertexRef primaryVertex(const reco::PFJetRef &jet) const
Get primary vertex associated to this jet.
void setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
void setleadPFCand(const T &cand)
Set leading PF candidate.
void addPFCands(Region region, ParticleType type, const InputIterator &begin, const InputIterator &end)
InputIterator leadPFCand(InputIterator begin, InputIterator end)
void initialize(const PtrIter &chargedHadronsBegin, const PtrIter &chargedHadronsEnd)
StringObjectFunction< reco::PFJet > JetFunc
std::vector< PFCandidatePtr > flattenPiZeros(const std::vector< RecoTauPiZero >::const_iterator &, const std::vector< RecoTauPiZero >::const_iterator &)
Flatten a list of pi zeros into a list of there constituent PFCandidates.
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
DeltaRFilter< RecoTauPiZero > PiZeroDRFilter
Definition: ConeTools.h:56
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:41
void setleadPFChargedHadrCand(const T &cand)
Set leading PFChargedHadron candidate.
return_type operator()(const reco::PFJetRef &jet, const std::vector< reco::PFRecoTauChargedHadron > &chargedHadrons, const std::vector< RecoTauPiZero > &piZeros, const std::vector< PFCandidatePtr > &regionalExtras) const override
std::vector< PFCandidatePtr > pfChargedCands(const PFJet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
void setleadPFNeutralCand(const T &cand)
Set leading PFGamma candidate.
boost::filter_iterator< PFCandPtrDRFilter, std::vector< PFCandidatePtr >::const_iterator > PFCandPtrDRFilterIter
Definition: ConeTools.h:50
#define DEFINE_EDM_PLUGIN(factory, type, name)
void initialize(const PtrIter &particlesBegin, const PtrIter &particlesEnd)
DeltaRPtrFilter< PFCandidatePtr > PFCandPtrDRFilter
Definition: ConeTools.h:48
std::vector< RecoTauPiZero > PiZeroList
Transform a pizero to remove given candidates.
PredicateAND< P1, P2 > makePredicateAND(const P1 &p1, const P2 &p2)