CMS 3D CMS Logo

RecoTauCommonUtilities.cc
Go to the documentation of this file.
2 
8 
10 
11 #include <algorithm>
12 
13 typedef std::vector<reco::CandidatePtr> CandPtrs;
14 typedef CandPtrs::iterator CandIter;
15 
16 namespace reco {
17  namespace tau {
18 
19  namespace {
20  // Re-implemented from PFCandidate.cc
21  int translateTypeToAbsPdgId(int type) {
22  switch (type) {
24  return 211; // pi+
26  return 11;
28  return 13;
30  return 22;
32  return 130; // K_L0
34  return 1; // dummy pdg code
36  return 2; // dummy pdg code
38  default:
39  return 0;
40  }
41  }
42  } // namespace
43  std::vector<CandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>::const_iterator& piZerosBegin,
44  const std::vector<RecoTauPiZero>::const_iterator& piZerosEnd) {
45  std::vector<CandidatePtr> output;
46 
47  for (std::vector<RecoTauPiZero>::const_iterator piZero = piZerosBegin; piZero != piZerosEnd; ++piZero) {
48  for (size_t iDaughter = 0; iDaughter < piZero->numberOfDaughters(); ++iDaughter) {
49  output.push_back(CandidatePtr(piZero->daughterPtr(iDaughter)));
50  }
51  }
52  return output;
53  }
54 
55  std::vector<CandidatePtr> flattenPiZeros(const std::vector<RecoTauPiZero>& piZeros) {
56  return flattenPiZeros(piZeros.begin(), piZeros.end());
57  }
58 
59  std::vector<reco::CandidatePtr> pfCandidates(const reco::Jet& jet, int particleId, bool sort) {
60  return pfCandidatesByPdgId(jet, translateTypeToAbsPdgId(particleId), sort);
61  }
62 
63  std::vector<CandidatePtr> pfCandidates(const Jet& jet, const std::vector<int>& particleIds, bool sort) {
64  std::vector<int> pdgIds;
65  pdgIds.reserve(particleIds.size());
66  for (auto particleId : particleIds)
67  pdgIds.push_back(translateTypeToAbsPdgId(particleId));
69  }
70 
71  std::vector<reco::CandidatePtr> pfCandidatesByPdgId(const reco::Jet& jet, int pdgId, bool sort) {
72  CandPtrs pfCands = jet.daughterPtrVector();
73  CandPtrs selectedPFCands = filterPFCandidates(pfCands.begin(), pfCands.end(), pdgId, sort);
74  return selectedPFCands;
75  }
76 
77  std::vector<reco::CandidatePtr> pfCandidatesByPdgId(const reco::Jet& jet,
78  const std::vector<int>& pdgIds,
79  bool sort) {
80  const CandPtrs& pfCands = jet.daughterPtrVector();
82  // Get each desired candidate type, unsorted for now
83  for (std::vector<int>::const_iterator pdgId = pdgIds.begin(); pdgId != pdgIds.end(); ++pdgId) {
84  CandPtrs&& selectedPFCands = filterPFCandidates(pfCands.begin(), pfCands.end(), *pdgId, false);
85  output.insert(output.end(), selectedPFCands.begin(), selectedPFCands.end());
86  }
87  if (sort)
89  return output;
90  }
91 
92  std::vector<reco::CandidatePtr> pfGammas(const reco::Jet& jet, bool sort) { return pfCandidates(jet, 22, sort); }
93 
94  std::vector<reco::CandidatePtr> pfChargedCands(const reco::Jet& jet, bool sort) {
95  const CandPtrs& pfCands = jet.daughterPtrVector();
97  CandPtrs&& mus = filterPFCandidates(pfCands.begin(), pfCands.end(), 13, false);
98  CandPtrs&& es = filterPFCandidates(pfCands.begin(), pfCands.end(), 11, false);
99  CandPtrs&& chs = filterPFCandidates(pfCands.begin(), pfCands.end(), 211, false);
100  output.reserve(mus.size() + es.size() + chs.size());
101  output.insert(output.end(), mus.begin(), mus.end());
102  output.insert(output.end(), es.begin(), es.end());
103  output.insert(output.end(), chs.begin(), chs.end());
104  if (sort)
105  std::sort(output.begin(), output.end(), SortPFCandsDescendingPt());
106  return output;
107  }
108 
110  const reco::PFCandidate* pfCand = dynamic_cast<const reco::PFCandidate*>(part);
111  if (pfCand)
112  return pfCand->positionAtECALEntrance();
113 
116  RawParticle(math::XYZTLorentzVector(part->px(), part->py(), part->pz(), part->energy()),
117  math::XYZTLorentzVector(part->vertex().x(), part->vertex().y(), part->vertex().z(), 0.)),
118  part->charge(),
119  0.,
120  0.,
121  bField);
122  theParticle.propagateToEcalEntrance(false);
123  if (theParticle.getSuccess() != 0) {
124  pos = math::XYZPointF(theParticle.particle().vertex());
125  }
126  return pos;
127  }
128 
129  } // namespace tau
130 } // namespace reco
reco::tau::pfChargedCands
std::vector< CandidatePtr > pfChargedCands(const Jet &jet, bool sort=true)
Extract all non-neutral candidates from a PFJet.
Definition: RecoTauCommonUtilities.cc:94
RawParticle
Definition: RawParticle.h:37
reco::Jet
Base class for all types of Jets.
Definition: Jet.h:20
BaseParticlePropagator::propagateToEcalEntrance
bool propagateToEcalEntrance(bool first=true)
Definition: BaseParticlePropagator.cc:450
BaseParticlePropagator::particle
RawParticle const & particle() const
The particle being propagated.
Definition: BaseParticlePropagator.h:164
reco::PFCandidate::e
Definition: PFCandidate.h:47
metsig::tau
Definition: SignAlgoResolutions.h:49
PFCandidate.h
convertSQLitetoXML_cfg.output
output
Definition: convertSQLitetoXML_cfg.py:72
RawParticle::vertex
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:320
pos
Definition: PixelAliasList.h:18
reco::PFCandidate::positionAtECALEntrance
const math::XYZPointF & positionAtECALEntrance() const
Definition: PFCandidate.h:363
reco::tau::atECALEntrance
math::XYZPointF atECALEntrance(const reco::Candidate *part, double bField)
Definition: RecoTauCommonUtilities.cc:109
reco::PFCandidate::h
Definition: PFCandidate.h:46
BaseParticlePropagator.h
reco::PFCandidate::h_HF
Definition: PFCandidate.h:51
Jet.h
BaseParticlePropagator::getSuccess
int getSuccess() const
Has propagation been performed and was barrel or endcap reached ?
Definition: BaseParticlePropagator.h:296
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
RawParticle::charge
double charge() const
get the MEASURED charge
Definition: RawParticle.h:294
reco::tau::flattenPiZeros
std::vector< CandidatePtr > flattenPiZeros(const std::vector< RecoTauPiZero >::const_iterator &, const std::vector< RecoTauPiZero >::const_iterator &)
Flatten a list of pi zeros into a list of there constituent PFCandidates.
Definition: RecoTauCommonUtilities.cc:43
BaseParticlePropagator
Definition: BaseParticlePropagator.h:82
reco::PFCandidate::X
Definition: PFCandidate.h:45
CandidateFwd.h
part
part
Definition: HCALResponse.h:20
reco::PFCandidate::mu
Definition: PFCandidate.h:48
Jet
Definition: Jet.py:1
RecoTauCommonUtilities.h
reco::tau::SortPFCandsDescendingPt
Definition: RecoTauCommonUtilities.h:33
Vertex.h
reco::tau::pfGammas
std::vector< CandidatePtr > pfGammas(const Jet &jet, bool sort=true)
Extract all pfGammas from a PFJet.
Definition: RecoTauCommonUtilities.cc:92
type
type
Definition: SiPixelVCal_PayloadInspector.cc:37
jetUpdater_cfi.sort
sort
Definition: jetUpdater_cfi.py:29
reco::PFCandidate::egamma_HF
Definition: PFCandidate.h:52
CosmicGenFilterHelix_cfi.pdgIds
pdgIds
Definition: CosmicGenFilterHelix_cfi.py:5
CandPtrs
std::vector< reco::CandidatePtr > CandPtrs
Definition: RecoTauCommonUtilities.cc:13
reco::tau::filterPFCandidates
std::vector< CandidatePtr > filterPFCandidates(const Iterator &begin, const Iterator &end, int pdgId, bool sort=true)
Definition: RecoTauCommonUtilities.h:41
reco::PFCandidate::gamma
Definition: PFCandidate.h:49
EgammaValidation_cff.pdgId
pdgId
Definition: EgammaValidation_cff.py:118
reco::Candidate
Definition: Candidate.h:27
reco::tau::pfCandidatesByPdgId
std::vector< CandidatePtr > pfCandidatesByPdgId(const Jet &jet, int pdgId, bool sort=true)
Definition: RecoTauCommonUtilities.cc:71
Calorimetry_cff.bField
bField
Definition: Calorimetry_cff.py:284
math::XYZTLorentzVector
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
metsig::jet
Definition: SignAlgoResolutions.h:47
reco::PFCandidate::h0
Definition: PFCandidate.h:50
reco::PFCandidate
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
Candidate.h
reco::CandidatePtr
edm::Ptr< Candidate > CandidatePtr
persistent reference to an object in a collection of Candidate objects
Definition: CandidateFwd.h:25
reco::tau::pfCandidates
std::vector< CandidatePtr > pfCandidates(const Jet &jet, int particleId, bool sort=true)
Definition: RecoTauCommonUtilities.cc:59
CandIter
CandPtrs::iterator CandIter
Definition: RecoTauCommonUtilities.cc:14
math::XYZPointF
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10