CMS 3D CMS Logo

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