CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Tau.h
Go to the documentation of this file.
1 //
2 //
3 
4 #ifndef DataFormats_PatCandidates_Tau_h
5 #define DataFormats_PatCandidates_Tau_h
6 
29 
34 
36 // Define typedefs for convenience
37 namespace pat {
38  class Tau;
39  typedef std::vector<Tau> TauCollection;
42 }
43 
44 namespace reco {
46  std::ostream& operator<<(std::ostream& out, const pat::Tau& obj);
47  // bool sortByPt(const CandidatePtrVector::const_iterator &lhs, const CandidatePtrVector::const_iterator &rhs) { return (*lhs)->pt() < (*rhs)->pt(); }
48 }
49 
50 // Class definition
51 namespace pat {
52 
53 
54  class PATTauSlimmer;
55 
56  class Tau : public Lepton<reco::BaseTau> {
60  friend class PATTauProducer;
61 
62  public:
63 
64  typedef std::pair<std::string, float> IdPair;
65 
67  Tau();
69  Tau(const reco::BaseTau & aTau);
71  Tau(const edm::RefToBase<reco::BaseTau> & aTauRef);
73  Tau(const edm::Ptr<reco::BaseTau> & aTauRef);
75  virtual ~Tau();
76 
78  virtual Tau * clone() const { return new Tau(*this); }
79 
80  // ---- methods for content embedding ----
82  const reco::TrackRefVector & isolationTracks() const;
84  reco::TrackRef leadTrack() const;
86  const reco::TrackRefVector & signalTracks() const;
88  void embedIsolationTracks();
90  void embedLeadTrack();
92  void embedSignalTracks();
95  void embedLeadPFCand();
101  void embedSignalPFCands();
109  void embedIsolationPFCands();
116 
117  // ---- matched GenJet methods ----
119  const reco::GenJet * genJet() const;
121  void setGenJet(const reco::GenJetRef & ref);
122 
123  // ---- CaloTau accessors (getters only) ----
125  bool isCaloTau() const { return !caloSpecific_.empty(); }
127  const pat::tau::TauCaloSpecific & caloSpecific() const ;
154  float maximumHCALhitEt() const { return caloSpecific().maximumHCALhitEt_; }
155 
156  // ---- PFTau accessors (getters only) ----
158  bool isPFTau() const { return !pfSpecific_.empty(); }
160  const pat::tau::TauPFSpecific & pfSpecific() const ;
161  const pat::tau::TauPFEssential & pfEssential() const;
164  const reco::PFJetRef & pfJetRef() const { return pfSpecific().pfJetRef_; }
179  const reco::PFCandidatePtr leadPFCand() const;
182  const std::vector<reco::PFCandidatePtr>& signalPFCands() const;
185  const std::vector<reco::PFCandidatePtr>& signalPFChargedHadrCands() const;
188  const std::vector<reco::PFCandidatePtr>& signalPFNeutrHadrCands() const;
191  const std::vector<reco::PFCandidatePtr>& signalPFGammaCands() const;
194  const std::vector<reco::PFRecoTauChargedHadron> & signalTauChargedHadronCandidates() const;
197  const std::vector<reco::RecoTauPiZero> & signalPiZeroCandidates() const;
200  const std::vector<reco::PFCandidatePtr>& isolationPFCands() const;
203  const std::vector<reco::PFCandidatePtr>& isolationPFChargedHadrCands() const;
206  const std::vector<reco::PFCandidatePtr>& isolationPFNeutrHadrCands() const;
209  const std::vector<reco::PFCandidatePtr>& isolationPFGammaCands() const;
212  const std::vector<reco::PFRecoTauChargedHadron> & isolationTauChargedHadronCandidates() const;
215  const std::vector<reco::RecoTauPiZero> & isolationPiZeroCandidates() const;
227  float emFraction() const { return pfSpecific().emFraction_; }
230  float hcalTotOverPLead() const { return pfSpecific().hcalTotOverPLead_; }
233  float hcalMaxOverPLead() const { return pfSpecific().hcalMaxOverPLead_; }
236  float hcal3x3OverPLead() const { return pfSpecific().hcal3x3OverPLead_; }
254  float caloComp() const { return pfSpecific().caloComp_; }
257  float segComp() const { return pfSpecific().segComp_; }
260  bool muonDecision() const { return pfSpecific().muonDecision_; }
261 
266  const reco::CandidatePtr leadNeutralCand() const;
268  const reco::CandidatePtr leadCand() const;
271  bool ExistSignalCands() const;
272  bool ExistIsolationCands() const;
295 
303 
306  size_t numberOfSourceCandidatePtrs() const ;
309 
310 
315  float dxy() const { return pfEssential().dxy_; }
316  float dxy_error() const { return pfEssential().dxy_error_; }
317  float dxy_Sig() const;
318  const reco::VertexRef& primaryVertex() const { return pfEssential().pv_; }
321  bool hasSecondaryVertex() const { return pfEssential().hasSV_; }
323  float flightLengthSig() const { return pfEssential().flightLengthSig_; }
325  const reco::VertexRef& secondaryVertex() const { return pfEssential().sv_; }
328 
332  float etaetaMoment() const;
333  float phiphiMoment() const;
334  float etaphiMoment() const;
335 
337  int decayMode() const { return pfEssential().decayMode_; }
339  void setDecayMode(int);
340 
341  // ---- methods for tau ID ----
347  float tauID(const std::string & name) const;
348  float tauID(const char* name ) const {return tauID( std::string(name) );}
350  bool isTauIDAvailable(const std::string & name) const;
353  const std::vector<IdPair> & tauIDs() const { return tauIDs_; }
356  void setTauIDs(const std::vector<IdPair> & ids) { tauIDs_ = ids; }
357 
359  friend std::ostream& reco::operator<<(std::ostream& out, const Tau& obj);
360 
363  const std::vector<std::string> availableJECSets() const;
364  // returns the available JEC Levels for a given jecSet
365  const std::vector<std::string> availableJECLevels(const int& set = 0) const;
366  // returns the available JEC Levels for a given jecSet
367  const std::vector<std::string> availableJECLevels(const std::string& set) const { return availableJECLevels(jecSet(set)); };
370  bool jecSetsAvailable() const { return !jec_.empty(); }
373  bool jecSetAvailable(const std::string& set) const {return (jecSet(set) >= 0); };
376  bool jecSetAvailable(const unsigned int& set) const {return (set < jec_.size()); };
379  return currentJECSet_<jec_.size() ? jec_.at(currentJECSet_).jecSet() : std::string("ERROR");
380  }
383  return currentJECSet_<jec_.size() ? jec_.at(currentJECSet_).jecLevel(currentJECLevel_) : std::string("ERROR");
384  }
387  float jecFactor(const std::string& level, const std::string& set = "") const;
390  float jecFactor(const unsigned int& level, const unsigned int& set = 0) const;
393  Tau correctedTauJet(const std::string& level, const std::string& set = "") const;
396  Tau correctedTauJet(const unsigned int& level, const unsigned int& set = 0) const;
399  const LorentzVector& correctedP4(const std::string& level, const std::string& set = "") const {
400  return correctedTauJet(level, set).p4();
401  }
404  const LorentzVector& correctedP4(const unsigned int& level, const unsigned int& set = 0) const {
405  return correctedTauJet(level, set).p4();
406  }
407 
408 
409  friend class PATTauSlimmer;
410 
411  protected:
412 
415  int jecSet(const std::string& label) const;
417  void currentJECSet(const unsigned int& set) { currentJECSet_=set; };
419  void currentJECLevel(const unsigned int& level) { currentJECLevel_=level; };
421  void addJECFactors(const TauJetCorrFactors& jec) {jec_.push_back(jec); };
423  void initializeJEC(unsigned int level, const unsigned int set = 0);
424 
425 
426  private:
427  // ---- for content embedding ----
429  std::vector<reco::Track> isolationTracks_;
432  std::vector<reco::Track> leadTrack_;
434  std::vector<reco::Track> signalTracks_;
436  // specific for PFTau
437  std::vector<reco::PFCandidate> leadPFCand_;
439  std::vector<reco::PFCandidate> leadPFChargedHadrCand_;
441  std::vector<reco::PFCandidate> leadPFNeutralCand_;
443 
444  std::vector<reco::PFCandidate> signalPFCands_;
447  std::vector<reco::PFCandidate> signalPFChargedHadrCands_;
450  std::vector<reco::PFCandidate> signalPFNeutralHadrCands_;
453  std::vector<reco::PFCandidate> signalPFGammaCands_;
456  std::vector<reco::PFCandidate> isolationPFCands_;
459  std::vector<reco::PFCandidate> isolationPFChargedHadrCands_;
462  std::vector<reco::PFCandidate> isolationPFNeutralHadrCands_;
465  std::vector<reco::PFCandidate> isolationPFGammaCands_;
468 
469  // ---- matched GenJet holder ----
470  std::vector<reco::GenJet> genJet_;
471 
472  // ---- tau ID's holder ----
473  std::vector<IdPair> tauIDs_;
474 
475  // ---- PFTau specific variables ----
477  std::vector<pat::tau::TauPFSpecific> pfSpecific_;
478 
479  // ---- CaloTau specific variables ----
481  std::vector<pat::tau::TauCaloSpecific> caloSpecific_;
482 
483  // ---- energy scale correction factors ----
484  // energy scale correction factors; the string carries a potential label if
485  // more then one set of correction factors is embedded. The label corresponds
486  // to the label of the jetCorrFactors module that has been embedded.
487  std::vector<pat::TauJetCorrFactors> jec_;
488  // currently applied set of jet energy correction factors (i.e. the index in
489  // jetEnergyCorrections_)
490  unsigned int currentJECSet_;
491  // currently applied jet energy correction level
492  unsigned int currentJECLevel_;
493 
494  // ---- references to packed pf candidates -----
498 
502 
503  // -- essential info to keep
504 
505  std::vector<pat::tau::TauPFEssential> pfEssential_;
506 
507 
508  };
509 }
510 
511 #endif
bool embeddedLeadPFNeutralCand_
Definition: Tau.h:442
virtual Tau * clone() const
required reimplementation of the Candidate&#39;s clone method
Definition: Tau.h:78
int i
Definition: DBlmapReader.cc:9
const std::vector< IdPair > & tauIDs() const
Definition: Tau.h:353
const LorentzVector & correctedP4(const std::string &level, const std::string &set="") const
Definition: Tau.h:399
edm::AtomicPtrCache< reco::TrackRefVector > signalTracksTransientRefVector_
Definition: Tau.h:435
const std::vector< reco::PFCandidatePtr > & isolationPFCands() const
const std::vector< std::string > availableJECSets() const
float hcalTotOverPLead() const
Definition: Tau.h:230
void currentJECLevel(const unsigned int &level)
update the current JEC level; used by correctedJet
Definition: Tau.h:419
const LorentzVector & correctedP4(const unsigned int &level, const unsigned int &set=0) const
Definition: Tau.h:404
const pat::tau::TauCaloSpecific & caloSpecific() const
return CaloTau info or throw exception &#39;not CaloTau&#39;
reco::CandidatePtrVector isolationGammaCandPtrs_
Definition: Tau.h:501
reco::CaloTauTagInfoRef CaloTauTagInfoRef_
pat::tau::TauPFEssential::CovMatrix flightLengthCov() const
const reco::CandidatePtr leadCand() const
return the PFCandidate if available (reference or embedded), or the PackedPFCandidate on miniAOD ...
const reco::TrackRefVector & isolationTracks() const
override the reco::BaseTau::isolationTracks method, to access the internal storage of the isolation t...
reco::CandidatePtrVector isolationChargedHadrCands() const
void embedIsolationPFCands()
method to store the isolation candidates internally
float caloComp() const
Definition: Tau.h:254
size_t size_type
Definition: Candidate.h:30
bool embeddedSignalPFChargedHadrCands_
Definition: Tau.h:448
float isolationTracksPtSum() const
Definition: Tau.h:148
std::vector< reco::PFCandidate > leadPFNeutralCand_
Definition: Tau.h:441
void embedIsolationPFGammaCands()
method to store the isolation gamma candidates internally
reco::TrackRef electronPreIDTrack_
Definition: TauPFSpecific.h:54
edm::AtomicPtrCache< std::vector< reco::PFCandidatePtr > > isolationPFChargedHadrCandsTransientPtrs_
Definition: Tau.h:461
float dxy() const
Definition: Tau.h:315
std::vector< reco::PFCandidate > signalPFNeutralHadrCands_
Definition: Tau.h:450
std::vector< reco::PFCandidate > signalPFCands_
Definition: Tau.h:444
size_t numberOfSourceCandidatePtrs() const
edm::AtomicPtrCache< std::vector< reco::PFCandidatePtr > > signalPFNeutralHadrCandsTransientPtrs_
Definition: Tau.h:452
bool embeddedIsolationPFChargedHadrCands_
Definition: Tau.h:460
bool jecSetAvailable(const std::string &set) const
Definition: Tau.h:373
bool embeddedIsolationPFCands_
Definition: Tau.h:457
const std::vector< reco::PFCandidatePtr > & isolationPFChargedHadrCands() const
float hcal3x3OverPLead() const
Definition: Tau.h:236
bool jecSetAvailable(const unsigned int &set) const
Definition: Tau.h:376
Slimmer of PAT Taus.
const pat::tau::TauPFEssential & pfEssential() const
const reco::PFCandidatePtr leadPFChargedHadrCand() const
reco::CaloTauTagInfoRef caloTauTagInfoRef() const
Definition: Tau.h:130
void embedLeadTrack()
method to store the leading track internally
const pat::tau::TauPFEssential::Point & secondaryVertexPos() const
Definition: Tau.h:326
void embedSignalPFChargedHadrCands()
method to store the signal charged hadrons candidates internally
std::vector< reco::PFCandidate > isolationPFGammaCands_
Definition: Tau.h:465
std::vector< reco::PFCandidate > isolationPFNeutralHadrCands_
Definition: Tau.h:462
edm::AtomicPtrCache< std::vector< reco::PFCandidatePtr > > isolationPFGammaCandsTransientPtrs_
Definition: Tau.h:467
float TracksInvariantMass() const
Definition: Tau.h:145
std::vector< Tau > TauCollection
Definition: Tau.h:38
const reco::CandidatePtr leadNeutralCand() const
return the PFCandidate if available (reference or embedded), or the PackedPFCandidate on miniAOD ...
const std::vector< reco::PFRecoTauChargedHadron > & signalTauChargedHadronCandidates() const
std::vector< reco::GenJet > genJet_
Definition: Tau.h:470
void setDecayMode(int)
set decay mode
reco::CandidatePtr sourceCandidatePtr(size_type i) const
get the source candidate pointer with index i
const reco::TrackRef & electronPreIDTrack() const
Definition: Tau.h:245
float leadTrackHCAL3x3hottesthitDEta() const
Definition: Tau.h:139
edm::Ref< TauCollection > TauRef
Definition: Tau.h:40
const std::vector< std::string > availableJECLevels(const std::string &set) const
Definition: Tau.h:367
void setTauIDs(const std::vector< IdPair > &ids)
Definition: Tau.h:356
edm::AtomicPtrCache< std::vector< reco::PFCandidatePtr > > signalPFCandsTransientPtrs_
Definition: Tau.h:446
float tauID(const std::string &name) const
bool embeddedIsolationPFNeutralHadrCands_
Definition: Tau.h:463
reco::CandidatePtrVector signalCands() const
float leadPFChargedHadrCandsignedSipt_
Definition: TauPFSpecific.h:29
bool embeddedLeadTrack_
Definition: Tau.h:431
float isolationPFGammaCandsEtSum() const
Definition: Tau.h:221
edm::AtomicPtrCache< std::vector< reco::PFCandidatePtr > > signalPFGammaCandsTransientPtrs_
Definition: Tau.h:455
int jecSet(const std::string &label) const
float phiphiMoment() const
edm::RefVector< TauCollection > TauRefVector
Definition: Tau.h:41
float maximumHCALPFClusterEt() const
Definition: Tau.h:224
std::pair< std::string, float > IdPair
Definition: Tau.h:64
std::vector< reco::PFCandidate > leadPFCand_
Definition: Tau.h:437
reco::CandidatePtrVector signalGammaCands() const
const pat::tau::TauPFEssential::CovMatrix & primaryVertexCov() const
Definition: Tau.h:320
const std::vector< reco::PFRecoTauChargedHadron > & isolationTauChargedHadronCandidates() const
bool embeddedSignalPFGammaCands_
Definition: Tau.h:454
reco::CandidatePtrVector isolationNeutrHadrCands() const
float maximumHCALhitEt() const
Definition: Tau.h:154
unsigned int currentJECLevel_
Definition: Tau.h:492
edm::AtomicPtrCache< reco::TrackRefVector > isolationTracksTransientRefVector_
Definition: Tau.h:430
const std::vector< reco::PFCandidatePtr > & signalPFCands() const
reco::PFRecoTauChargedHadronRef leadTauChargedHadronCandidate() const
float hcalMaxOverPLead() const
Definition: Tau.h:233
void setIsolationGammaCands(const reco::CandidatePtrVector &ptrs)
Definition: Tau.h:302
void embedSignalPFCands()
method to store the signal candidates internally
void currentJECSet(const unsigned int &set)
update the current JEC set; used by correctedJet
Definition: Tau.h:417
float isolationPFChargedHadrCandsPtSum() const
Definition: Tau.h:218
reco::CandidatePtrVector isolationChargedHadrCandPtrs_
Definition: Tau.h:499
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:71
std::vector< reco::Track > signalTracks_
Definition: Tau.h:434
std::vector< IdPair > tauIDs_
Definition: Tau.h:473
float signalTracksInvariantMass() const
Definition: Tau.h:142
edm::AtomicPtrCache< std::vector< reco::PFCandidatePtr > > signalPFChargedHadrCandsTransientPtrs_
Definition: Tau.h:449
bool embeddedIsolationPFGammaCands_
Definition: Tau.h:466
math::XYZVectorF Vector
const reco::VertexRef & primaryVertex() const
Definition: Tau.h:318
float etaetaMoment() const
std::vector< reco::PFCandidate > signalPFChargedHadrCands_
Definition: Tau.h:447
Tau()
default constructor
void setIsolationNeutralHadrCands(const reco::CandidatePtrVector &ptrs)
Definition: Tau.h:301
const pat::tau::TauPFEssential::Vector & flightLength() const
Definition: Tau.h:322
Analysis-level lepton class.
Definition: Lepton.h:30
void setIsolationChargedHadrCands(const reco::CandidatePtrVector &ptrs)
Definition: Tau.h:300
float isolationECALhitsEtSum() const
Definition: Tau.h:151
Tau correctedTauJet(const std::string &level, const std::string &set="") const
float jecFactor(const std::string &level, const std::string &set="") const
reco::CandidatePtrVector isolationNeutralHadrCandPtrs_
Definition: Tau.h:500
Jets made from MC generator particles.
Definition: GenJet.h:24
bool embeddedSignalTracks_
Definition: Tau.h:433
const reco::VertexRef & secondaryVertex() const
Definition: Tau.h:325
void setSignalChargedHadrCands(const reco::CandidatePtrVector &ptrs)
setters for the PtrVectors (for miniAOD)
Definition: Tau.h:297
bool embeddedLeadPFCand_
Definition: Tau.h:438
bool hasSecondaryVertex() const
Definition: Tau.h:321
reco::CandidatePtrVector signalNeutrHadrCands() const
const reco::CandidatePtr leadChargedHadrCand() const
std::vector< pat::tau::TauPFEssential > pfEssential_
Definition: Tau.h:505
void embedIsolationPFNeutralHadrCands()
method to store the isolation neutral hadrons candidates internally
const std::vector< reco::RecoTauPiZero > & signalPiZeroCandidates() const
void embedLeadPFChargedHadrCand()
method to store the leading charged hadron candidate internally
std::string currentJECLevel() const
return the name of the current step of jet energy corrections
Definition: Tau.h:382
const reco::PFCandidatePtr leadPFCand() const
Analysis-level tau class.
Definition: Tau.h:56
const pat::tau::TauPFEssential::CovMatrix & secondaryVertexCov() const
Definition: Tau.h:327
bool embeddedLeadPFChargedHadrCand_
Definition: Tau.h:440
void embedIsolationPFChargedHadrCands()
method to store the isolation charged hadrons candidates internally
bool embeddedSignalPFNeutralHadrCands_
Definition: Tau.h:451
std::vector< reco::PFCandidate > isolationPFCands_
Definition: Tau.h:456
std::vector< reco::Track > isolationTracks_
Definition: Tau.h:429
tuple out
Definition: dbtoconf.py:99
void setGenJet(const reco::GenJetRef &ref)
set the matched GenJet
bool isTauIDAvailable(const std::string &name) const
Returns true if a specific ID is available in this pat::Tau.
const std::vector< reco::PFCandidatePtr > & signalPFChargedHadrCands() const
float isolationPFChargedHadrCandsPtSum_
Definition: TauPFSpecific.h:44
void embedLeadPFNeutralCand()
method to store the leading neutral candidate internally
const pat::tau::TauPFEssential::Point & dxy_PCA() const
Definition: Tau.h:314
float emFraction() const
Definition: Tau.h:227
std::vector< pat::tau::TauCaloSpecific > caloSpecific_
holder for CaloTau info, or empty vector if PFTau
Definition: Tau.h:481
void setSignalGammaCands(const reco::CandidatePtrVector &ptrs)
Definition: Tau.h:299
float ecalStripSumEOverPLead() const
Definition: Tau.h:239
bool isCaloTau() const
Returns true if this pat::Tau was made from a reco::CaloTau.
Definition: Tau.h:125
const std::vector< reco::PFCandidatePtr > & signalPFNeutrHadrCands() const
float leadTrackHCAL3x3hitsEtSum() const
Definition: Tau.h:136
const std::vector< reco::PFCandidatePtr > & isolationPFGammaCands() const
reco::CandidatePtrVector signalGammaCandPtrs_
Definition: Tau.h:497
int decayMode() const
reconstructed tau decay mode (specific to PFTau)
Definition: Tau.h:337
const reco::PFJetRef & pfJetRef() const
Definition: Tau.h:164
float bremsRecoveryEOverPLead() const
Definition: Tau.h:242
std::vector< pat::TauJetCorrFactors > jec_
Definition: Tau.h:487
float etaphiMoment() const
edm::AtomicPtrCache< std::vector< reco::PFCandidatePtr > > isolationPFNeutralHadrCandsTransientPtrs_
Definition: Tau.h:464
reco::PFJetRef pfJetRef_
Definition: TauPFSpecific.h:27
std::vector< reco::PFCandidate > signalPFGammaCands_
Definition: Tau.h:453
float leadPFChargedHadrCandsignedSipt() const
Definition: Tau.h:173
std::vector< reco::Track > leadTrack_
Definition: Tau.h:432
std::vector< reco::PFCandidate > isolationPFChargedHadrCands_
Definition: Tau.h:459
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
bool jecSetsAvailable() const
Definition: Tau.h:370
void initializeJEC(unsigned int level, const unsigned int set=0)
initialize the jet to a given JEC level during creation starting from Uncorrected ...
const pat::tau::TauPFSpecific & pfSpecific() const
return PFTau info or throw exception &#39;not PFTau&#39;
void embedSignalPFNeutralHadrCands()
method to store the signal neutral hadrons candidates internally
reco::CandidatePtrVector signalNeutralHadrCandPtrs_
Definition: Tau.h:496
const reco::GenJet * genJet() const
return matched GenJet, built from the visible particles of a generated tau
float flightLengthSig() const
Definition: Tau.h:323
void addJECFactors(const TauJetCorrFactors &jec)
add more sets of energy correction factors
Definition: Tau.h:421
bool embeddedSignalPFCands_
Definition: Tau.h:445
float dxy_Sig() const
float tauID(const char *name) const
Definition: Tau.h:348
const reco::TrackRefVector & signalTracks() const
override the reco::BaseTau::signalTracks method, to access the internal storage of the signal tracks ...
reco::TrackRef leadTrack() const
override the reco::BaseTau::leadTrack method, to access the internal storage of the leading track ...
float segComp() const
Definition: Tau.h:257
void embedSignalTracks()
method to store the signal tracks internally
reco::CandidatePtrVector isolationGammaCands() const
bool isPFTau() const
Returns true if this pat::Tau was made from a reco::PFTau.
Definition: Tau.h:158
const std::vector< std::string > availableJECLevels(const int &set=0) const
bool embeddedIsolationTracks_
Definition: Tau.h:428
float leadTracksignedSipt() const
Definition: Tau.h:133
reco::CandidatePtrVector signalChargedHadrCands() const
tuple level
Definition: testEve_cfg.py:34
void embedSignalPFGammaCands()
method to store the signal gamma candidates internally
const std::vector< reco::RecoTauPiZero > & isolationPiZeroCandidates() const
edm::AtomicPtrCache< std::vector< reco::PFCandidatePtr > > isolationPFCandsTransientPtrs_
Definition: Tau.h:458
bool ExistIsolationCands() const
std::vector< reco::PFCandidate > leadPFChargedHadrCand_
Definition: Tau.h:439
std::string currentJECSet() const
returns the label of the current set of jet energy corrections
Definition: Tau.h:378
bool ExistSignalCands() const
const std::vector< reco::PFCandidatePtr > & signalPFGammaCands() const
virtual const LorentzVector & p4() const
four-momentum Lorentz vector
Definition: LeafCandidate.h:99
reco::CandidatePtrVector isolationCands() const
void embedIsolationTracks()
method to store the isolation tracks internally
virtual ~Tau()
destructor
reco::Candidate::LorentzVector p4Jet() const
std::vector< pat::tau::TauPFSpecific > pfSpecific_
holder for PFTau info, or empty vector if CaloTau
Definition: Tau.h:477
const std::vector< reco::PFCandidatePtr > & isolationPFNeutrHadrCands() const
float electronPreIDOutput() const
Definition: Tau.h:248
reco::CandidatePtrVector signalChargedHadrCandPtrs_
Definition: Tau.h:495
unsigned int currentJECSet_
Definition: Tau.h:490
bool electronPreIDDecision() const
Definition: Tau.h:251
void setSignalNeutralHadrCands(const reco::CandidatePtrVector &ptrs)
Definition: Tau.h:298
bool muonDecision() const
Definition: Tau.h:260
Produces pat::Tau&#39;s.
const pat::tau::TauPFEssential::Point & primaryVertexPos() const
Definition: Tau.h:319
void embedLeadPFCand()
float dxy_error() const
Definition: Tau.h:316
math::ErrorF< 3 >::type CovMatrix
const reco::PFCandidatePtr leadPFNeutralCand() const