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  * $Id $
10  */
11 
12 #include <vector>
13 #include <algorithm>
14 
19 
22 
27 
29 
30 namespace reco { namespace tau {
31 
33  public:
34  explicit RecoTauBuilderConePlugin(const edm::ParameterSet& pset);
36  // Build a tau from a jet
38  const std::vector<RecoTauPiZero>& piZeros,
39  const std::vector<PFCandidatePtr>& regionalExtras) const;
40  private:
44  /* String function to extract values from PFJets */
46  // Cone defintions
54 };
55 
56 // ctor - initialize all of our variables
58  const edm::ParameterSet& pset):RecoTauBuilderPlugin(pset),
59  qcuts_(pset.getParameterSet(
60  "qualityCuts").getParameterSet("signalQualityCuts")),
61  usePFLeptonsAsChargedHadrons_(pset.getParameter<bool>("usePFLeptons")),
62  leadObjecPtThreshold_(pset.getParameter<double>("leadObjectPt")),
63  matchingCone_(pset.getParameter<std::string>("matchingCone")),
64  signalConeChargedHadrons_(pset.getParameter<std::string>(
65  "signalConeChargedHadrons")),
66  isoConeChargedHadrons_(
67  pset.getParameter<std::string>("isoConeChargedHadrons")),
68  signalConePiZeros_(
69  pset.getParameter<std::string>("signalConePiZeros")),
70  isoConePiZeros_(
71  pset.getParameter<std::string>("isoConePiZeros")),
72  signalConeNeutralHadrons_(
73  pset.getParameter<std::string>("signalConeNeutralHadrons")),
74  isoConeNeutralHadrons_(
75  pset.getParameter<std::string>("isoConeNeutralHadrons")) {}
76 
78  const reco::PFJetRef& jet,
79  const std::vector<RecoTauPiZero>& piZeros,
80  const std::vector<PFCandidatePtr>& regionalExtras) const {
81  // Get access to our cone tools
82  using namespace cone;
83  // Define output. We only produce one tau per jet for the cone algo.
85 
86  // Our tau builder - the true indicates to automatically copy gamma candidates
87  // from the pizeros.
88  RecoTauConstructor tau(jet, getPFCands(), true);
89  // Setup our quality cuts to use the current vertex (supplied by base class)
91 
92  typedef std::vector<PFCandidatePtr> PFCandPtrs;
93 
94  // Get the PF Charged hadrons + quality cuts
95  PFCandPtrs pfchs;
98  } else {
99  // Check if we want to include electrons in muons in "charged hadron"
100  // collection. This is the preferred behavior, as the PF lepton selections
101  // are very loose.
102  pfchs = qcuts_.filterRefs(pfChargedCands(*jet));
103  }
104 
105  // Get the PF gammas
108  // Neutral hadrons
109  PFCandPtrs pfnhs = qcuts_.filterRefs(
111 
112  // All the extra junk
113  PFCandPtrs regionalJunk = qcuts_.filterRefs(regionalExtras);
114 
115  /***********************************************
116  ****** Lead Candidate Finding **********
117  ***********************************************/
118 
119  // Define our matching cone and filters
120  double matchingCone = matchingCone_(*jet);
121  PFCandPtrDRFilter matchingConeFilter(jet->p4(), 0, matchingCone);
122 
123  // Find the maximum PFCharged hadron in the matching cone. The call to
124  // PFCandidates always a sorted list, so we can just take the first if it
125  // if it exists.
126  PFCandidatePtr leadPFCH;
127  PFCandPtrs::iterator leadPFCH_iter =
128  std::find_if(pfchs.begin(), pfchs.end(), matchingConeFilter);
129 
130  if (leadPFCH_iter != pfchs.end()) {
131  leadPFCH = *leadPFCH_iter;
132  // Set leading candidate
133  tau.setleadPFChargedHadrCand(leadPFCH);
134  } else {
135  // If there is no leading charged candidate at all, return nothing - the
136  // producer class that owns the plugin will build a null tau if desired.
137  return output.release();
138  }
139 
140  // Find the leading neutral candidate
141  PFCandidatePtr leadPFGamma;
142  PFCandPtrs::iterator leadPFGamma_iter =
143  std::find_if(pfGammas.begin(), pfGammas.end(), matchingConeFilter);
144 
145  if (leadPFGamma_iter != pfGammas.end()) {
146  leadPFGamma = *leadPFGamma_iter;
147  // Set leading neutral candidate
148  tau.setleadPFNeutralCand(leadPFGamma);
149  }
150 
152  // Always use the leadPFCH if it is above our threshold
153  if (leadPFCH.isNonnull() && leadPFCH->pt() > leadObjecPtThreshold_) {
154  leadPFCand = leadPFCH;
155  } else if (leadPFGamma.isNonnull() &&
156  leadPFGamma->pt() > leadObjecPtThreshold_) {
157  // Otherwise use the lead Gamma if it is above threshold
158  leadPFCand = leadPFGamma;
159  } else {
160  // If both are too low PT, just take the charged one
161  leadPFCand = leadPFCH;
162  }
163 
165 
166  // Our cone axis is defined about the lead charged hadron
167  reco::Candidate::LorentzVector coneAxis = leadPFCH->p4();
168 
169  /***********************************************
170  ****** Cone Construction **********
171  ***********************************************/
172 
173  // Define the signal and isolation cone sizes for this jet and build filters
174  // to select elements in the given DeltaR regions
175 
176  PFCandPtrDRFilter signalConePFCHFilter(
177  coneAxis, -0.1, signalConeChargedHadrons_(*jet));
178  PFCandPtrDRFilter signalConePFNHFilter(
179  coneAxis, -0.1, signalConeNeutralHadrons_(*jet));
180  PiZeroDRFilter signalConePiZeroFilter(
181  coneAxis, -0.1, signalConePiZeros_(*jet));
182 
183  PFCandPtrDRFilter isoConePFCHFilter(
184  coneAxis, signalConeChargedHadrons_(*jet), isoConeChargedHadrons_(*jet));
185  PFCandPtrDRFilter isoConePFGammaFilter(
186  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
187  PFCandPtrDRFilter isoConePFNHFilter(
188  coneAxis, signalConeNeutralHadrons_(*jet), isoConeNeutralHadrons_(*jet));
189  PiZeroDRFilter isoConePiZeroFilter(
190  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
191 
192  // Additionally make predicates to select the different PF object types
193  // of the regional junk objects to add to the iso cone.
195  PFCandPtrDRFilter> RegionalJunkConeAndIdFilter;
196 
198  pfchCandSelector(reco::PFCandidate::h);
200  pfgammaCandSelector(reco::PFCandidate::gamma);
202  pfnhCandSelector(reco::PFCandidate::h0);
203 
204  // Predicate to select the regional junk in the iso cone by PF id
205  RegionalJunkConeAndIdFilter pfChargedJunk(
206  pfchCandSelector, // select charged stuff from junk
207  isoConePFCHFilter // only take those in iso cone
208  );
209 
210  RegionalJunkConeAndIdFilter pfGammaJunk(
211  pfgammaCandSelector, // select gammas from junk
212  isoConePFGammaFilter // only take those in iso cone
213  );
214 
215  RegionalJunkConeAndIdFilter pfNeutralJunk(
216  pfnhCandSelector, // select neutral stuff from junk
217  isoConePFNHFilter // select stuff in iso cone
218  );
219 
220  // Build filter iterators select the signal charged stuff.
221  PFCandPtrDRFilterIter signalPFCHs_begin(
222  signalConePFCHFilter, pfchs.begin(), pfchs.end());
223  PFCandPtrDRFilterIter signalPFCHs_end(
224  signalConePFCHFilter, pfchs.end(), pfchs.end());
225 
226  // Cross clean pi zero content using signal cone charged hadron constituents.
228  signalPFCHs_begin, signalPFCHs_end);
229  std::vector<reco::RecoTauPiZero> cleanPiZeros = xCleaner(piZeros);
230 
231  // For the rest of the constituents, we need to filter any constituents that
232  // are already contained in the pizeros (i.e. electrons)
233  xclean::CrossCleanPtrs pfCandXCleaner(cleanPiZeros);
234 
235  // Build signal charged hadrons
238  signalPFCHs_begin, signalPFCHs_end);
239 
242  boost::make_filter_iterator(
243  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner),
244  pfnhs.begin(), pfnhs.end()),
245  boost::make_filter_iterator(
246  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner),
247  pfnhs.end(), pfnhs.end()));
248 
249  // Build signal PiZeros
251  PiZeroDRFilterIter(signalConePiZeroFilter,
252  cleanPiZeros.begin(), cleanPiZeros.end()),
253  PiZeroDRFilterIter(signalConePiZeroFilter,
254  cleanPiZeros.end(), cleanPiZeros.end()));
255 
256  // Build isolation charged hadrons
259  boost::make_filter_iterator(
260  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner),
261  pfchs.begin(), pfchs.end()),
262  boost::make_filter_iterator(
263  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner),
264  pfchs.end(), pfchs.end()));
265 
266  // Add all the stuff in the isolation cone that wasn't in the jet constituents
269  boost::make_filter_iterator(
270  pfChargedJunk, regionalJunk.begin(), regionalJunk.end()),
271  boost::make_filter_iterator(
272  pfChargedJunk, regionalJunk.end(), regionalJunk.end())
273  );
274 
275  // Build isolation neutral hadrons
278  boost::make_filter_iterator(
279  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner),
280  pfnhs.begin(), pfnhs.end()),
281  boost::make_filter_iterator(
282  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner),
283  pfnhs.end(), pfnhs.end()));
284 
285  // Add regional stuff not in jet
288  boost::make_filter_iterator(
289  pfNeutralJunk, regionalJunk.begin(), regionalJunk.end()),
290  boost::make_filter_iterator(
291  pfNeutralJunk, regionalJunk.end(), regionalJunk.end())
292  );
293 
294  // Build isolation PiZeros
296  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.begin(),
297  cleanPiZeros.end()),
298  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.end(),
299  cleanPiZeros.end()));
300 
301  // Add regional stuff not in jet
304  boost::make_filter_iterator(
305  pfGammaJunk, regionalJunk.begin(), regionalJunk.end()),
306  boost::make_filter_iterator(
307  pfGammaJunk, regionalJunk.end(), regionalJunk.end())
308  );
309 
310  // Put our built tau in the output - 'false' indicates don't build the
311  // leading candidates, we already did that explicitly above.
312 
313  std::auto_ptr<reco::PFTau> tauPtr = tau.get(false);
314 
315  // Set event vertex position for tau
316  reco::VertexRef primaryVertexRef = primaryVertex(jet);
317  if ( primaryVertexRef.isNonnull() )
318  tauPtr->setVertex(primaryVertexRef->position());
319 
320  output.push_back(tauPtr);
321  return output.release();
322 }
323 }} // end namespace reco::tauk
324 
328  "RecoTauBuilderConePlugin");
const edm::Handle< PFCandidateCollection > & getPFCands() const
Hack to be able to convert Ptrs to Refs.
return_type operator()(const reco::PFJetRef &jet, const std::vector< RecoTauPiZero > &piZeros, const std::vector< PFCandidatePtr > &regionalExtras) const
std::auto_ptr< reco::PFTau > get(bool setupLeadingCandidates=true)
Coll filterRefs(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:48
boost::ptr_vector< reco::PFTau > output_type
std::auto_ptr< output_type > return_type
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 pi zeros 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)
StringObjectFunction< reco::PFJet > JetFunc
RecoTauBuilderConePlugin(const edm::ParameterSet &pset)
DeltaRFilter< RecoTauPiZero > PiZeroDRFilter
Definition: ConeTools.h:46
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:38
void setleadPFChargedHadrCand(const T &cand)
Set leading PFChargedHadron candidate.
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:44
#define DEFINE_EDM_PLUGIN(factory, type, name)
DeltaRPtrFilter< PFCandidatePtr > PFCandPtrDRFilter
Definition: ConeTools.h:42
Transform a pizero to remove given candidates.
PredicateAND< P1, P2 > makePredicateAND(const P1 &p1, const P2 &p2)