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 
20 
24 
26 
27 namespace reco { namespace tau {
28 
30  public:
33  // Build a tau from a jet
35  const std::vector<RecoTauPiZero>& piZeros) const;
36  private:
40  /* String function to extract values from PFJets */
42  // Cone defintions
50 };
51 
52 // ctor - initialize all of our variables
55  qcuts_(pset.getParameter<edm::ParameterSet>("qualityCuts")),
56  usePFLeptonsAsChargedHadrons_(pset.getParameter<bool>("usePFLeptons")),
57  leadObjecPtThreshold_(pset.getParameter<double>("leadObjectPt")),
58  matchingCone_(pset.getParameter<std::string>("matchingCone")),
59  signalConeChargedHadrons_(pset.getParameter<std::string>(
60  "signalConeChargedHadrons")),
61  isoConeChargedHadrons_(
62  pset.getParameter<std::string>("isoConeChargedHadrons")),
63  signalConePiZeros_(
64  pset.getParameter<std::string>("signalConePiZeros")),
65  isoConePiZeros_(
66  pset.getParameter<std::string>("isoConePiZeros")),
67  signalConeNeutralHadrons_(
68  pset.getParameter<std::string>("signalConeNeutralHadrons")),
69  isoConeNeutralHadrons_(
70  pset.getParameter<std::string>("isoConeNeutralHadrons"))
71 {}
72 
74  const reco::PFJetRef& jet,
75  const std::vector<RecoTauPiZero>& piZeros) const {
76  // Get access to our cone tools
77  using namespace cone;
78  // Define output. We only produce one tau per jet for the cone algo.
80 
81  // Our tau builder - the true indicates to automatically copy gamma candidates
82  // from the pizeros.
83  RecoTauConstructor tau(jet, getPFCands(), true);
84  // Setup our quality cuts to use the current vertex (supplied by base class)
86 
87  typedef std::vector<PFCandidatePtr> PFCandPtrs;
88 
89  // Get the PF Charged hadrons + quality cuts
90  PFCandPtrs pfchs;
93  } else {
94  // Check if we want to include electrons in muons in "charged hadron"
95  // collection. This is the preferred behavior, as the PF lepton selections
96  // are very loose.
97  pfchs = qcuts_.filterRefs(pfChargedCands(*jet));
98  }
99 
100  // Get the PF gammas
103  // Neutral hadrons
104  PFCandPtrs pfnhs = qcuts_.filterRefs(
106 
107  /***********************************************
108  ****** Lead Candidate Finding **********
109  ***********************************************/
110 
111  // Define our matching cone and filters
112  double matchingCone = matchingCone_(*jet);
113  PFCandPtrDRFilter matchingConeFilter(jet->p4(), 0, matchingCone);
114 
115  // Find the maximum PFCharged hadron in the matching cone. The call to
116  // PFCandidates always a sorted list, so we can just take the first if it
117  // if it exists.
118  PFCandidatePtr leadPFCH;
119  PFCandPtrs::iterator leadPFCH_iter =
120  std::find_if(pfchs.begin(), pfchs.end(), matchingConeFilter);
121 
122  if (leadPFCH_iter != pfchs.end()) {
123  leadPFCH = *leadPFCH_iter;
124  // Set leading candidate
125  tau.setleadPFChargedHadrCand(leadPFCH);
126  } else {
127  // If there is no leading charged candidate at all, return nothing - the
128  // producer class that owns the plugin will build a null tau if desired.
129  return output.release();
130  }
131 
132  // Find the leading neutral candidate
133  PFCandidatePtr leadPFGamma;
134  PFCandPtrs::iterator leadPFGamma_iter =
135  std::find_if(pfGammas.begin(), pfGammas.end(), matchingConeFilter);
136 
137  if (leadPFGamma_iter != pfGammas.end()) {
138  leadPFGamma = *leadPFGamma_iter;
139  // Set leading neutral candidate
140  tau.setleadPFNeutralCand(leadPFGamma);
141  }
142 
144  // Always use the leadPFCH if it is above our threshold
145  if (leadPFCH.isNonnull() && leadPFCH->pt() > leadObjecPtThreshold_) {
146  leadPFCand = leadPFCH;
147  } else if (leadPFGamma.isNonnull() &&
148  leadPFGamma->pt() > leadObjecPtThreshold_) {
149  // Otherwise use the lead Gamma if it is above threshold
150  leadPFCand = leadPFGamma;
151  } else {
152  // If both are too low PT, just take the charged one
153  leadPFCand = leadPFCH;
154  }
155 
157 
158  // Our cone axis is defined about the lead charged hadron
159  reco::Candidate::LorentzVector coneAxis = leadPFCH->p4();
160 
161  /***********************************************
162  ****** Cone Construction **********
163  ***********************************************/
164 
165  // Define the signal and isolation cone sizes for this jet and build filters
166  // to select elements in the given DeltaR regions
167 
168  PFCandPtrDRFilter signalConePFCHFilter(
169  coneAxis, -0.1, signalConeChargedHadrons_(*jet));
170  PFCandPtrDRFilter signalConePFNHFilter(
171  coneAxis, -0.1, signalConeNeutralHadrons_(*jet));
172  PiZeroDRFilter signalConePiZeroFilter(
173  coneAxis, -0.1, signalConePiZeros_(*jet));
174 
175  PFCandPtrDRFilter isoConePFCHFilter(
176  coneAxis, signalConeChargedHadrons_(*jet), isoConeChargedHadrons_(*jet));
177  PFCandPtrDRFilter isoConePFNHFilter(
178  coneAxis, signalConeNeutralHadrons_(*jet), isoConeNeutralHadrons_(*jet));
179  PiZeroDRFilter isoConePiZeroFilter(
180  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
181 
182  // Build signal charged hadrons
185  PFCandPtrDRFilterIter(signalConePFCHFilter, pfchs.begin(),
186  pfchs.end()),
187  PFCandPtrDRFilterIter(signalConePFCHFilter, pfchs.end(),
188  pfchs.end()));
189 
190  // Build signal neutral hadrons
193  PFCandPtrDRFilterIter(signalConePFNHFilter, pfnhs.begin(),
194  pfnhs.end()),
195  PFCandPtrDRFilterIter(signalConePFNHFilter, pfnhs.end(),
196  pfnhs.end()));
197 
198  // Build signal PiZeros
200  PiZeroDRFilterIter(signalConePiZeroFilter,
201  piZeros.begin(), piZeros.end()),
202  PiZeroDRFilterIter(signalConePiZeroFilter,
203  piZeros.end(), piZeros.end()));
204 
205  // Build isolation charged hadrons
208  PFCandPtrDRFilterIter(isoConePFCHFilter, pfchs.begin(),
209  pfchs.end()),
210  PFCandPtrDRFilterIter(isoConePFCHFilter, pfchs.end(),
211  pfchs.end()));
212 
213  // Build isolation neutral hadrons
216  PFCandPtrDRFilterIter(isoConePFNHFilter, pfnhs.begin(),
217  pfnhs.end()),
218  PFCandPtrDRFilterIter(isoConePFNHFilter, pfnhs.end(),
219  pfnhs.end()));
220 
221  // Build isolation PiZeros
223  PiZeroDRFilterIter(isoConePiZeroFilter, piZeros.begin(),
224  piZeros.end()),
225  PiZeroDRFilterIter(isoConePiZeroFilter, piZeros.end(),
226  piZeros.end()));
227 
228  // Put our built tau in the output - 'false' indicates don't build the
229  // leading candidtes, we already did that explicitly above.
230 
231  output.push_back(tau.get(false));
232  return output.release();
233 }
234 }} // end namespace reco::tauk
235 
239  "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)
void addPiZeros(Region region, InputIterator begin, InputIterator end)
Add a list of pi zeros to the input collection.
std::vector< PFCandidatePtr > pfGammas(const PFJet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
return_type operator()(const reco::PFJetRef &jet, const std::vector< RecoTauPiZero > &piZeros) const
std::vector< reco::PFCandidatePtr > PFCandPtrs
boost::filter_iterator< PiZeroDRFilter, std::vector< RecoTauPiZero >::const_iterator > PiZeroDRFilterIter
Definition: ConeTools.h:47
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 setPV(const reco::VertexRef &vtx) const
Update the primary vertex.
Coll filterRefs(const Coll &refcoll) const
Filter a ref vector of PFCandidates.
tuple pset
Definition: CrabTask.py:85
void setleadPFCand(const T &cand)
Set leading PF candidate.
InputIterator leadPFCand(InputIterator begin, InputIterator end)
StringObjectFunction< reco::PFJet > JetFunc
RecoTauBuilderConePlugin(const edm::ParameterSet &pset)
DeltaRFilter< RecoTauPiZero > PiZeroDRFilter
Definition: ConeTools.h:45
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:39
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.
void addPFCands(Region region, ParticleType type, InputIterator begin, InputIterator end)
boost::filter_iterator< PFCandPtrDRFilter, std::vector< PFCandidatePtr >::const_iterator > PFCandPtrDRFilterIter
Definition: ConeTools.h:43
#define DEFINE_EDM_PLUGIN(factory, type, name)
const reco::VertexRef & primaryVertex() const
Get primary vertex associated to this event.
DeltaRPtrFilter< PFCandidatePtr > PFCandPtrDRFilter
Definition: ConeTools.h:41