CMS 3D CMS Logo

HPSPFTau.cc
Go to the documentation of this file.
3 
4 // default constructor
6  : tauType_(kUndefined),
7  sumChargedIso_(0.),
8  sumNeutralIso_(0.),
9  sumCombinedIso_(0.),
10  sumChargedIsoPileup_(0.),
11  passTightIso_(false),
12  passMediumIso_(false),
13  passLooseIso_(false),
14  passVLooseIso_(false),
15  passTightRelIso_(false),
16  passMediumRelIso_(false),
17  passLooseRelIso_(false),
18  passVLooseRelIso_(false) {}
19 
20 // destructor
22 
23 // print to stream
24 ostream& operator<<(ostream& os, const l1t::HPSPFTau& l1PFTau) {
25  os << "pT = " << l1PFTau.pt() << ", eta = " << l1PFTau.eta() << ", phi = " << l1PFTau.phi()
26  << " (type = " << l1PFTau.tauType() << ")" << std::endl;
27  os << "lead. ChargedPFCand:" << std::endl;
28  if (l1PFTau.leadChargedPFCand().isNonnull()) {
29  printPFCand(os, *l1PFTau.leadChargedPFCand(), l1PFTau.primaryVertex());
30  } else {
31  os << " N/A" << std::endl;
32  }
33  os << "seed:";
34  if (l1PFTau.isChargedPFCandSeeded()) {
35  os << " chargedPFCand";
36  } else if (l1PFTau.isJetSeeded()) {
37  os << " CaloJet";
38  } else {
39  cms::Exception ex("InconsistentTau");
40  ex.addContext("Calling HPSPFTau::operator <<");
41  ex.addAdditionalInfo("This tau is not seed by either a chargedPFCand or a PFJet!");
42  throw ex;
43  }
44  os << std::endl;
45  os << "signalPFCands:" << std::endl;
46  for (const auto& l1PFCand : l1PFTau.signalAllL1PFCandidates()) {
47  printPFCand(os, *l1PFCand, l1PFTau.primaryVertex());
48  }
49  os << "stripPFCands:" << std::endl;
50  for (const auto& l1PFCand : l1PFTau.stripAllL1PFCandidates()) {
51  printPFCand(os, *l1PFCand, l1PFTau.primaryVertex());
52  }
53  os << "strip pT = " << l1PFTau.stripP4().pt() << std::endl;
54  os << "isolationPFCands:" << std::endl;
55  for (const auto& l1PFCand : l1PFTau.isoAllL1PFCandidates()) {
56  printPFCand(os, *l1PFCand, l1PFTau.primaryVertex());
57  }
58  os << "isolation pT-sum: charged = " << l1PFTau.sumChargedIso() << ", neutral = " << l1PFTau.sumNeutralIso()
59  << " (charged from pileup = " << l1PFTau.sumChargedIsoPileup() << ")" << std::endl;
60  return os;
61 }
62 
63 void printPFCand(ostream& os, const l1t::PFCandidate& l1PFCand, const l1t::TkPrimaryVertexRef& primaryVertex) {
64  float primaryVertexZ = (primaryVertex.isNonnull()) ? primaryVertex->zvertex() : 0.;
65  printPFCand(os, l1PFCand, primaryVertexZ);
66 }
67 
68 void printPFCand(ostream& os, const l1t::PFCandidate& l1PFCand, float primaryVertexZ) {
69  std::string typeString;
70  if (l1PFCand.id() == l1t::PFCandidate::ChargedHadron)
71  typeString = "PFChargedHadron";
72  else if (l1PFCand.id() == l1t::PFCandidate::Electron)
73  typeString = "PFElectron";
74  else if (l1PFCand.id() == l1t::PFCandidate::NeutralHadron)
75  typeString = "PFNeutralHadron";
76  else if (l1PFCand.id() == l1t::PFCandidate::Photon)
77  typeString = "PFPhoton";
78  else if (l1PFCand.id() == l1t::PFCandidate::Muon)
79  typeString = "PFMuon";
80  else
81  typeString = "N/A";
82  os << " " << typeString << " with pT = " << l1PFCand.pt() << ", eta = " << l1PFCand.eta()
83  << ", phi = " << l1PFCand.phi() << ","
84  << " mass = " << l1PFCand.mass() << ", charge = " << l1PFCand.charge();
85  if (l1PFCand.charge() != 0 && primaryVertexZ != 0.) {
86  os << " (dz = " << std::fabs(l1PFCand.pfTrack()->vertex().z() - primaryVertexZ) << ")";
87  }
88  os << std::endl;
89 }
const l1t::PFCandidateRefVector & isoAllL1PFCandidates() const
Definition: HPSPFTau.h:41
const l1t::PFCandidateRefVector & signalAllL1PFCandidates() const
Definition: HPSPFTau.h:30
double pt() const final
transverse momentum
float sumChargedIsoPileup() const
Definition: HPSPFTau.h:69
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
void printPFCand(ostream &os, const l1t::PFCandidate &l1PFCand, const l1t::TkPrimaryVertexRef &primaryVertex)
Definition: HPSPFTau.cc:63
float sumChargedIso() const
Definition: HPSPFTau.h:66
~HPSPFTau() override
destructor
Definition: HPSPFTau.cc:21
void addAdditionalInfo(std::string const &info)
Definition: Exception.cc:169
const l1t::TkPrimaryVertexRef & primaryVertex() const
Definition: HPSPFTau.h:55
std::ostream & operator<<(std::ostream &os, const l1t::CaloParamsHelper &p)
HPSPFTau()
default constructor
Definition: HPSPFTau.cc:5
void addContext(std::string const &context)
Definition: Exception.cc:165
float sumNeutralIso() const
Definition: HPSPFTau.h:67
bool isJetSeeded() const
Definition: HPSPFTau.h:24
double mass() const final
mass
const l1t::PFCandidateRefVector & stripAllL1PFCandidates() const
Definition: HPSPFTau.h:37
const reco::Particle::LorentzVector & stripP4() const
Definition: HPSPFTau.h:60
const PFTrackRef & pfTrack() const
Definition: PFCandidate.h:36
bool isChargedPFCandSeeded() const
accessor functions for reco level quantities
Definition: HPSPFTau.h:23
const l1t::PFCandidateRef & leadChargedPFCand() const
Definition: HPSPFTau.h:28
Kind tauType() const
Definition: HPSPFTau.h:58
double phi() const final
momentum azimuthal angle
primaryVertex
hltOfflineBeamSpot for HLTMON
ParticleType id() const
Definition: PFCandidate.h:34
int charge() const final
electric charge
double eta() const final
momentum pseudorapidity