CMS 3D CMS Logo

PFBlockBasedIsolation.cc
Go to the documentation of this file.
7 
8 #include <cmath>
9 
10 //--------------------------------------------------------------------------------------------------
11 
13  // Default Constructor.
14 }
15 
16 //--------------------------------------------------------------------------------------------------
18 
19 void PFBlockBasedIsolation::setup(const edm::ParameterSet& conf) { coneSize_ = conf.getParameter<double>("coneSize"); }
20 
21 std::vector<reco::PFCandidateRef> PFBlockBasedIsolation::calculate(
23  const reco::PFCandidateRef pfEGCand,
24  const edm::Handle<reco::PFCandidateCollection> pfCandidateHandle) {
25  std::vector<reco::PFCandidateRef> myVec;
26 
27  math::XYZVector candidateMomentum(p4.px(), p4.py(), p4.pz());
28  math::XYZVector candidateDirection = candidateMomentum.Unit();
29 
30  const reco::PFCandidate::ElementsInBlocks& theElementsInpfEGcand = (*pfEGCand).elementsInBlocks();
31  reco::PFCandidate::ElementsInBlocks::const_iterator ieg = theElementsInpfEGcand.begin();
32  const reco::PFBlockRef egblock = ieg->first;
33 
34  unsigned nObj = pfCandidateHandle->size();
35  for (unsigned int lCand = 0; lCand < nObj; lCand++) {
36  reco::PFCandidateRef pfCandRef(reco::PFCandidateRef(pfCandidateHandle, lCand));
37 
38  float dR = 0.0;
39  if (coneSize_ < 10.0) {
40  dR = deltaR(candidateDirection.Eta(), candidateDirection.Phi(), pfCandRef->eta(), pfCandRef->phi());
41  if (dR > coneSize_)
42  continue;
43  }
44 
45  const reco::PFCandidate::ElementsInBlocks& theElementsInPFcand = pfCandRef->elementsInBlocks();
46 
47  bool elementFound = false;
48  for (reco::PFCandidate::ElementsInBlocks::const_iterator ipf = theElementsInPFcand.begin();
49  ipf < theElementsInPFcand.end();
50  ++ipf) {
51  if (ipf->first == egblock && !elementFound) {
52  for (ieg = theElementsInpfEGcand.begin(); ieg < theElementsInpfEGcand.end(); ++ieg) {
53  if (ipf->second == ieg->second && !elementFound) {
54  if (elementPassesCleaning(pfCandRef, pfEGCand)) {
55  myVec.push_back(pfCandRef);
56  elementFound = true;
57  }
58  }
59  }
60  }
61  }
62  }
63 
64  return myVec;
65 }
66 
68  const reco::PFCandidateRef& pfEGCand) {
69  if (pfCand->particleId() == reco::PFCandidate::h)
70  return passesCleaningChargedHadron(pfCand, pfEGCand);
71  else if (pfCand->particleId() == reco::PFCandidate::h0)
72  return passesCleaningNeutralHadron(pfCand, pfEGCand);
73  else if (pfCand->particleId() == reco::PFCandidate::gamma)
74  return passesCleaningPhoton(pfCand, pfEGCand);
75  else
76  return true; //doesnt really matter here as if its not a photon,neutral or charged it wont be included in isolation
77 }
78 
79 //currently the record of which candidates came from the charged hadron is acceptable, no further cleaning is needed
81  const reco::PFCandidateRef& pfEGCand) {
82  return true;
83 }
84 
85 //neutral hadrons are not part of the PF E/gamma reco, therefore they cant currently come from an electron/photon and so should be rejected
86 //but we still think there may be some useful info here and given we can easily
87 //fix this at AOD level, we will auto accept them for now and clean later
89  const reco::PFCandidateRef& pfEGCand) {
90  return true;
91 }
92 
93 //the highest et ECAL element of the photon must match to the electron superclusters or one of its sub clusters
95  const reco::PFCandidateRef& pfEGCand) {
96  bool passesCleaning = false;
97  const reco::PFBlockElementCluster* ecalClusWithMaxEt = getHighestEtECALCluster(*pfCand);
98  if (ecalClusWithMaxEt) {
99  if (ecalClusWithMaxEt->superClusterRef().isNonnull() &&
100  ecalClusWithMaxEt->superClusterRef()->seed()->seed() ==
101  pfEGCand->superClusterRef()
102  ->seed()
103  ->seed()) { //being sure to match, some concerned about different collections, shouldnt be but to be safe
104  passesCleaning = true;
105  } else {
106  for (auto cluster : pfEGCand->superClusterRef()->clusters()) {
107  //the PF clusters there are in two different collections so cant reference match
108  //but we can match on the seed id, no clusters can share a seed so if the seeds are
109  //equal, it must be the same cluster
110  if (ecalClusWithMaxEt->clusterRef()->seed() == cluster->seed()) {
111  passesCleaning = true;
112  }
113  } //end of loop over clusters
114  }
115  } //end of null check for highest ecal cluster
116  return passesCleaning;
117 }
118 
120  float maxECALEt = -1;
121  const reco::PFBlockElement* maxEtECALCluster = nullptr;
122  const reco::PFCandidate::ElementsInBlocks& elementsInPFCand = pfCand.elementsInBlocks();
123  for (auto& elemIndx : elementsInPFCand) {
124  const reco::PFBlockElement* elem =
125  elemIndx.second < elemIndx.first->elements().size() ? &elemIndx.first->elements()[elemIndx.second] : nullptr;
126  if (elem && elem->type() == reco::PFBlockElement::ECAL && elem->clusterRef()->pt() > maxECALEt) {
127  maxECALEt = elem->clusterRef()->pt();
128  maxEtECALCluster = elem;
129  }
130  }
131  return dynamic_cast<const reco::PFBlockElementCluster*>(maxEtECALCluster);
132 }
Abstract base class for a PFBlock element (track, cluster...)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
ROOT::Math::LorentzVector< ROOT::Math::PxPyPzE4D< double > > XYZTLorentzVectorD
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:14
bool elementPassesCleaning(const reco::PFCandidateRef &pfCand, const reco::PFCandidateRef &pfEGCand)
std::vector< reco::PFCandidateRef > calculate(math::XYZTLorentzVectorD p4, const reco::PFCandidateRef pfEGCand, const edm::Handle< reco::PFCandidateCollection > pfCandidateHandle)
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:232
const SuperClusterRef & superClusterRef() const
std::vector< ElementInBlock > ElementsInBlocks
Definition: PFCandidate.h:404
const ElementsInBlocks & elementsInBlocks() const
Definition: PFCandidate.cc:665
bool passesCleaningNeutralHadron(const reco::PFCandidateRef &pfCand, const reco::PFCandidateRef &pfEGCand)
void setup(const edm::ParameterSet &conf)
virtual const PFClusterRef & clusterRef() const
bool passesCleaningPhoton(const reco::PFCandidateRef &pfCand, const reco::PFCandidateRef &pfEGCand)
const PFClusterRef & clusterRef() const override
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
Particle reconstructed by the particle flow algorithm.
Definition: PFCandidate.h:41
const reco::PFBlockElementCluster * getHighestEtECALCluster(const reco::PFCandidate &pfCand)
bool passesCleaningChargedHadron(const reco::PFCandidateRef &pfCand, const reco::PFCandidateRef &pfEGCand)