test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Types | Private Attributes
reco::tau::RecoTauBuilderConePlugin Class Reference
Inheritance diagram for reco::tau::RecoTauBuilderConePlugin:
reco::tau::RecoTauBuilderPlugin reco::tau::RecoTauEventHolderPlugin reco::tau::RecoTauNamedPlugin

Public Member Functions

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
 
 RecoTauBuilderConePlugin (const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
 
 ~RecoTauBuilderConePlugin ()
 
- Public Member Functions inherited from reco::tau::RecoTauBuilderPlugin
virtual void beginEvent ()
 
const edm::Handle
< PFCandidateCollection > & 
getPFCands () const
 Hack to be able to convert Ptrs to Refs. More...
 
reco::VertexRef primaryVertex (const reco::PFJetRef &jet) const
 Get primary vertex associated to this jet. More...
 
 RecoTauBuilderPlugin (const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
 
virtual ~RecoTauBuilderPlugin ()
 
- Public Member Functions inherited from reco::tau::RecoTauEventHolderPlugin
const edm::Eventevt () const
 
edm::Eventevt ()
 
const edm::EventSetupevtSetup () const
 
 RecoTauEventHolderPlugin (const edm::ParameterSet &pset)
 
void setup (edm::Event &, const edm::EventSetup &)
 
virtual ~RecoTauEventHolderPlugin ()
 
- Public Member Functions inherited from reco::tau::RecoTauNamedPlugin
const std::string & name () const
 
 RecoTauNamedPlugin (const edm::ParameterSet &pset)
 
virtual ~RecoTauNamedPlugin ()
 

Private Types

typedef StringObjectFunction
< reco::PFJet
JetFunc
 

Private Attributes

JetFunc isoConeChargedHadrons_
 
JetFunc isoConeNeutralHadrons_
 
JetFunc isoConePiZeros_
 
double leadObjecPtThreshold_
 
JetFunc matchingCone_
 
int maxSignalConeChargedHadrons_
 
RecoTauQualityCuts qcuts_
 
JetFunc signalConeChargedHadrons_
 
JetFunc signalConeNeutralHadrons_
 
JetFunc signalConePiZeros_
 
bool usePFLeptonsAsChargedHadrons_
 

Additional Inherited Members

- Public Types inherited from reco::tau::RecoTauBuilderPlugin
typedef boost::ptr_vector
< reco::PFTau
output_type
 
typedef std::auto_ptr
< output_type
return_type
 

Detailed Description

Definition at line 34 of file RecoTauBuilderConePlugin.cc.

Member Typedef Documentation

Definition at line 51 of file RecoTauBuilderConePlugin.cc.

Constructor & Destructor Documentation

reco::tau::RecoTauBuilderConePlugin::RecoTauBuilderConePlugin ( const edm::ParameterSet pset,
edm::ConsumesCollector &&  iC 
)
explicit

Definition at line 66 of file RecoTauBuilderConePlugin.cc.

69  "qualityCuts").getParameterSet("signalQualityCuts")),
70  usePFLeptonsAsChargedHadrons_(pset.getParameter<bool>("usePFLeptons")),
71  leadObjecPtThreshold_(pset.getParameter<double>("leadObjectPt")),
72  matchingCone_(pset.getParameter<std::string>("matchingCone")),
74  "signalConeChargedHadrons")),
76  pset.getParameter<std::string>("isoConeChargedHadrons")),
78  pset.getParameter<std::string>("signalConePiZeros")),
80  pset.getParameter<std::string>("isoConePiZeros")),
82  pset.getParameter<std::string>("signalConeNeutralHadrons")),
84  pset.getParameter<std::string>("isoConeNeutralHadrons")),
86  pset.getParameter<int>("maxSignalConeChargedHadrons"))
87 {}
T getParameter(std::string const &) const
def move
Definition: eostools.py:510
ParameterSet const & getParameterSet(std::string const &) const
RecoTauBuilderPlugin(const edm::ParameterSet &pset, edm::ConsumesCollector &&iC)
reco::tau::RecoTauBuilderConePlugin::~RecoTauBuilderConePlugin ( )
inline

Definition at line 37 of file RecoTauBuilderConePlugin.cc.

37 {}

Member Function Documentation

RecoTauBuilderConePlugin::return_type reco::tau::RecoTauBuilderConePlugin::operator() ( const reco::PFJetRef ,
const std::vector< reco::PFRecoTauChargedHadron > &  ,
const std::vector< RecoTauPiZero > &  ,
const std::vector< PFCandidatePtr > &   
) const
overridevirtual

Construct one or more PFTaus from the a PFJet and its asscociated reconstructed PiZeros and regional extras i.e. objects in a 0.8 cone about the jet

Implements reco::tau::RecoTauBuilderPlugin.

Definition at line 110 of file RecoTauBuilderConePlugin.cc.

References reco::tau::RecoTauConstructor::addPFCands(), reco::tau::RecoTauConstructor::addPiZeros(), reco::tau::RecoTauQualityCuts::filterCandRefs(), reco::PFCandidate::gamma, reco::tau::RecoTauConstructor::get(), reco::tau::RecoTauBuilderPlugin::getPFCands(), reco::PFCandidate::h, reco::PFCandidate::h0, isoConeChargedHadrons_, isoConeNeutralHadrons_, isoConePiZeros_, reco::tau::RecoTauConstructor::kChargedHadron, reco::tau::RecoTauConstructor::kGamma, reco::tau::RecoTauConstructor::kIsolation, reco::tau::RecoTauConstructor::kNeutralHadron, reco::tau::RecoTauConstructor::kSignal, leadObjecPtThreshold_, reco::tau::leadPFCand(), reco::tau::xclean::makePredicateAND(), matchingCone_, maxSignalConeChargedHadrons_, convertSQLitetoXML_cfg::output, reco::tau::pfCandidates(), reco::tau::pfChargedCands(), reco::tau::pfGammas(), reco::tau::RecoTauBuilderPlugin::primaryVertex(), qcuts_, reco::tau::RecoTauConstructor::setleadPFCand(), reco::tau::RecoTauConstructor::setleadPFChargedHadrCand(), reco::tau::RecoTauConstructor::setleadPFNeutralCand(), reco::tau::RecoTauQualityCuts::setPV(), signalConeChargedHadrons_, signalConeNeutralHadrons_, signalConePiZeros_, metsig::tau, and usePFLeptonsAsChargedHadrons_.

114  {
115  //std::cout << "<RecoTauBuilderConePlugin::operator()>:" << std::endl;
116  //std::cout << " jet: Pt = " << jet->pt() << ", eta = " << jet->eta() << ", phi = " << jet->phi() << std::endl;
117 
118  // Get access to our cone tools
119  using namespace cone;
120  // Define output. We only produce one tau per jet for the cone algo.
122 
123  // Our tau builder - the true indicates to automatically copy gamma candidates
124  // from the pizeros.
125  RecoTauConstructor tau(jet, getPFCands(), true);
126  // Setup our quality cuts to use the current vertex (supplied by base class)
128 
129  typedef std::vector<PFCandidatePtr> PFCandPtrs;
130 
131  // Get the PF Charged hadrons + quality cuts
132  PFCandPtrs pfchs;
135  } else {
136  // Check if we want to include electrons in muons in "charged hadron"
137  // collection. This is the preferred behavior, as the PF lepton selections
138  // are very loose.
140  }
141 
142  // CV: sort collection of PF Charged hadrons by descending Pt
143  std::sort(pfchs.begin(), pfchs.end(), SortPFCandsDescendingPt());
144 
145  // Get the PF gammas
148  // Neutral hadrons
151 
152  // All the extra junk
153  PFCandPtrs regionalJunk = qcuts_.filterCandRefs(regionalExtras);
154 
155  /***********************************************
156  ****** Lead Candidate Finding **********
157  ***********************************************/
158 
159  // Define our matching cone and filters
160  double matchingCone = matchingCone_(*jet);
161  PFCandPtrDRFilter matchingConeFilter(jet->p4(), 0, matchingCone);
162 
163  // Find the maximum PFCharged hadron in the matching cone. The call to
164  // PFCandidates always a sorted list, so we can just take the first if it
165  // if it exists.
166  PFCandidatePtr leadPFCH;
167  PFCandPtrs::iterator leadPFCH_iter =
168  std::find_if(pfchs.begin(), pfchs.end(), matchingConeFilter);
169 
170  if (leadPFCH_iter != pfchs.end()) {
171  leadPFCH = *leadPFCH_iter;
172  // Set leading candidate
173  tau.setleadPFChargedHadrCand(leadPFCH);
174  } else {
175  // If there is no leading charged candidate at all, return nothing - the
176  // producer class that owns the plugin will build a null tau if desired.
177  return output.release();
178  }
179 
180  // Find the leading neutral candidate
181  PFCandidatePtr leadPFGamma;
182  PFCandPtrs::iterator leadPFGamma_iter =
183  std::find_if(pfGammas.begin(), pfGammas.end(), matchingConeFilter);
184 
185  if (leadPFGamma_iter != pfGammas.end()) {
186  leadPFGamma = *leadPFGamma_iter;
187  // Set leading neutral candidate
188  tau.setleadPFNeutralCand(leadPFGamma);
189  }
190 
192  // Always use the leadPFCH if it is above our threshold
193  if (leadPFCH.isNonnull() && leadPFCH->pt() > leadObjecPtThreshold_) {
194  leadPFCand = leadPFCH;
195  } else if (leadPFGamma.isNonnull() &&
196  leadPFGamma->pt() > leadObjecPtThreshold_) {
197  // Otherwise use the lead Gamma if it is above threshold
198  leadPFCand = leadPFGamma;
199  } else {
200  // If both are too low PT, just take the charged one
201  leadPFCand = leadPFCH;
202  }
203 
204  tau.setleadPFCand(leadPFCand);
205 
206  // Our cone axis is defined about the lead charged hadron
207  reco::Candidate::LorentzVector coneAxis = leadPFCH->p4();
208 
209  /***********************************************
210  ****** Cone Construction **********
211  ***********************************************/
212 
213  // Define the signal and isolation cone sizes for this jet and build filters
214  // to select elements in the given DeltaR regions
215 
216  PFCandPtrDRFilter signalConePFCHFilter(
217  coneAxis, -0.1, signalConeChargedHadrons_(*jet));
218  PFCandPtrDRFilter signalConePFNHFilter(
219  coneAxis, -0.1, signalConeNeutralHadrons_(*jet));
220  PiZeroDRFilter signalConePiZeroFilter(
221  coneAxis, -0.1, signalConePiZeros_(*jet));
222 
223  PFCandPtrDRFilter isoConePFCHFilter(
225  PFCandPtrDRFilter isoConePFGammaFilter(
226  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
227  PFCandPtrDRFilter isoConePFNHFilter(
229  PiZeroDRFilter isoConePiZeroFilter(
230  coneAxis, signalConePiZeros_(*jet), isoConePiZeros_(*jet));
231 
232  // Additionally make predicates to select the different PF object types
233  // of the regional junk objects to add to the iso cone.
234  typedef xclean::PredicateAND<xclean::FilterPFCandByParticleId, PFCandPtrDRFilter> RegionalJunkConeAndIdFilter;
235 
236  xclean::FilterPFCandByParticleId pfchCandSelector(reco::PFCandidate::h);
237  xclean::FilterPFCandByParticleId pfgammaCandSelector(reco::PFCandidate::gamma);
238  xclean::FilterPFCandByParticleId pfnhCandSelector(reco::PFCandidate::h0);
239 
240  // Predicate to select the regional junk in the iso cone by PF id
241  RegionalJunkConeAndIdFilter pfChargedJunk(
242  pfchCandSelector, // select charged stuff from junk
243  isoConePFCHFilter // only take those in iso cone
244  );
245 
246  RegionalJunkConeAndIdFilter pfGammaJunk(
247  pfgammaCandSelector, // select gammas from junk
248  isoConePFGammaFilter // only take those in iso cone
249  );
250 
251  RegionalJunkConeAndIdFilter pfNeutralJunk(
252  pfnhCandSelector, // select neutral stuff from junk
253  isoConePFNHFilter // select stuff in iso cone
254  );
255 
256  // Build filter iterators select the signal charged stuff.
257  PFCandPtrDRFilterIter signalPFCHCands_begin(
258  signalConePFCHFilter, pfchs.begin(), pfchs.end());
259  PFCandPtrDRFilterIter signalPFCHCands_end(
260  signalConePFCHFilter, pfchs.end(), pfchs.end());
261  PFCandPtrs signalPFCHs;
262  int numSignalPFCHs = 0;
263  PFCandPtrs isolationPFCHs;
264  int numIsolationPFCHs = 0;
265  for ( PFCandPtrDRFilterIter iter = signalPFCHCands_begin; iter != signalPFCHCands_end; ++iter ) {
266  if ( numSignalPFCHs < maxSignalConeChargedHadrons_ || maxSignalConeChargedHadrons_ == -1 ) {
267  //std::cout << "adding signalPFCH #" << numSignalPFCHs << ": Pt = " << (*iter)->pt() << ", eta = " << (*iter)->eta() << ", phi = " << (*iter)->phi() << std::endl;
268  signalPFCHs.push_back(*iter);
269  ++numSignalPFCHs;
270  } else {
271  //std::cout << "maxSignalConeChargedHadrons reached"
272  // << " --> adding isolationPFCH #" << numIsolationPFCHs << ": Pt = " << (*iter)->pt() << ", eta = " << (*iter)->eta() << ", phi = " << (*iter)->phi() << std::endl;
273  isolationPFCHs.push_back(*iter);
274  ++numIsolationPFCHs;
275  }
276  }
277  PFCandPtrs::const_iterator signalPFCHs_begin = signalPFCHs.begin();
278  PFCandPtrs::const_iterator signalPFCHs_end = signalPFCHs.end();
279 
280  // Cross clean pi zero content using signal cone charged hadron constituents.
281  xclean::CrossCleanPiZeros<PFCandPtrDRFilterIter> piZeroXCleaner(
282  signalPFCHCands_begin, signalPFCHCands_end);
283  std::vector<reco::RecoTauPiZero> cleanPiZeros = piZeroXCleaner(piZeros);
284 
285  // For the rest of the constituents, we need to filter any constituents that
286  // are already contained in the pizeros (i.e. electrons)
287  xclean::CrossCleanPtrs<PiZeroList::const_iterator> pfCandXCleaner(cleanPiZeros.begin(), cleanPiZeros.end());
288 
289  auto isolationPFCHCands_begin(
290  boost::make_filter_iterator(
291  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner),
292  pfchs.begin(), pfchs.end()));
293  auto isolationPFCHCands_end(
294  boost::make_filter_iterator(
295  xclean::makePredicateAND(isoConePFCHFilter, pfCandXCleaner),
296  pfchs.end(), pfchs.end()));
297  for ( auto iter = isolationPFCHCands_begin; iter != isolationPFCHCands_end; ++iter ) {
298  //std::cout << "adding isolationPFCH #" << numIsolationPFCHs << ": Pt = " << (*iter)->pt() << ", eta = " << (*iter)->eta() << ", phi = " << (*iter)->phi() << std::endl;
299  isolationPFCHs.push_back(*iter);
300  ++numIsolationPFCHs;
301  }
302  PFCandPtrs::const_iterator isolationPFCHs_begin = isolationPFCHs.begin();
303  PFCandPtrs::const_iterator isolationPFCHs_end = isolationPFCHs.end();
304 
305  // Build signal charged hadrons
306  tau.addPFCands(RecoTauConstructor::kSignal,
308  signalPFCHs_begin, signalPFCHs_end);
309 
310  tau.addPFCands(RecoTauConstructor::kSignal,
312  boost::make_filter_iterator(
313  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner),
314  pfnhs.begin(), pfnhs.end()),
315  boost::make_filter_iterator(
316  xclean::makePredicateAND(signalConePFNHFilter, pfCandXCleaner),
317  pfnhs.end(), pfnhs.end()));
318 
319  // Build signal PiZeros
320  tau.addPiZeros(RecoTauConstructor::kSignal,
321  PiZeroDRFilterIter(signalConePiZeroFilter,
322  cleanPiZeros.begin(), cleanPiZeros.end()),
323  PiZeroDRFilterIter(signalConePiZeroFilter,
324  cleanPiZeros.end(), cleanPiZeros.end()));
325 
326  // Build isolation charged hadrons
329  isolationPFCHs_begin, isolationPFCHs_end);
330 
331  // Add all the stuff in the isolation cone that wasn't in the jet constituents
334  boost::make_filter_iterator(
335  pfChargedJunk, regionalJunk.begin(), regionalJunk.end()),
336  boost::make_filter_iterator(
337  pfChargedJunk, regionalJunk.end(), regionalJunk.end())
338  );
339 
340  // Build isolation neutral hadrons
343  boost::make_filter_iterator(
344  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner),
345  pfnhs.begin(), pfnhs.end()),
346  boost::make_filter_iterator(
347  xclean::makePredicateAND(isoConePFNHFilter, pfCandXCleaner),
348  pfnhs.end(), pfnhs.end()));
349 
350  // Add regional stuff not in jet
353  boost::make_filter_iterator(
354  pfNeutralJunk, regionalJunk.begin(), regionalJunk.end()),
355  boost::make_filter_iterator(
356  pfNeutralJunk, regionalJunk.end(), regionalJunk.end())
357  );
358 
359  // Build isolation PiZeros
361  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.begin(),
362  cleanPiZeros.end()),
363  PiZeroDRFilterIter(isoConePiZeroFilter, cleanPiZeros.end(),
364  cleanPiZeros.end()));
365 
366  // Add regional stuff not in jet
369  boost::make_filter_iterator(
370  pfGammaJunk, regionalJunk.begin(), regionalJunk.end()),
371  boost::make_filter_iterator(
372  pfGammaJunk, regionalJunk.end(), regionalJunk.end())
373  );
374 
375  // Put our built tau in the output - 'false' indicates don't build the
376  // leading candidates, we already did that explicitly above.
377 
378  std::auto_ptr<reco::PFTau> tauPtr = tau.get(false);
379 
380  // Set event vertex position for tau
381  reco::VertexRef primaryVertexRef = primaryVertex(jet);
382  if ( primaryVertexRef.isNonnull() )
383  tauPtr->setVertex(primaryVertexRef->position());
384 
385  output.push_back(tauPtr);
386  return output.release();
387 }
const edm::Handle< PFCandidateCollection > & getPFCands() const
Hack to be able to convert Ptrs to Refs.
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.
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::vector< PFCandidatePtr > pfCandidates(const PFJet &jet, int particleId, bool sort=true)
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.
InputIterator leadPFCand(InputIterator begin, InputIterator end)
DeltaRFilter< RecoTauPiZero > PiZeroDRFilter
Definition: ConeTools.h:56
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
std::vector< PFCandidatePtr > pfChargedCands(const PFJet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
boost::filter_iterator< PFCandPtrDRFilter, std::vector< PFCandidatePtr >::const_iterator > PFCandPtrDRFilterIter
Definition: ConeTools.h:50
edm::Ptr< PFCandidate > PFCandidatePtr
persistent Ptr to a PFCandidate
DeltaRPtrFilter< PFCandidatePtr > PFCandPtrDRFilter
Definition: ConeTools.h:48
PredicateAND< P1, P2 > makePredicateAND(const P1 &p1, const P2 &p2)

Member Data Documentation

JetFunc reco::tau::RecoTauBuilderConePlugin::isoConeChargedHadrons_
private

Definition at line 56 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

JetFunc reco::tau::RecoTauBuilderConePlugin::isoConeNeutralHadrons_
private

Definition at line 60 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

JetFunc reco::tau::RecoTauBuilderConePlugin::isoConePiZeros_
private

Definition at line 58 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

double reco::tau::RecoTauBuilderConePlugin::leadObjecPtThreshold_
private

Definition at line 48 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

JetFunc reco::tau::RecoTauBuilderConePlugin::matchingCone_
private

Definition at line 54 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

int reco::tau::RecoTauBuilderConePlugin::maxSignalConeChargedHadrons_
private

Definition at line 62 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

RecoTauQualityCuts reco::tau::RecoTauBuilderConePlugin::qcuts_
private

Definition at line 44 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

JetFunc reco::tau::RecoTauBuilderConePlugin::signalConeChargedHadrons_
private

Definition at line 55 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

JetFunc reco::tau::RecoTauBuilderConePlugin::signalConeNeutralHadrons_
private

Definition at line 59 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

JetFunc reco::tau::RecoTauBuilderConePlugin::signalConePiZeros_
private

Definition at line 57 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().

bool reco::tau::RecoTauBuilderConePlugin::usePFLeptonsAsChargedHadrons_
private

Definition at line 46 of file RecoTauBuilderConePlugin.cc.

Referenced by operator()().