00001 #include <TFile.h>
00002 #include "EgammaAnalysis/ElectronTools/interface/PFIsolationEstimator.h"
00003 #include <cmath>
00004 #include "DataFormats/Math/interface/deltaR.h"
00005 using namespace std;
00006
00007 #ifndef STANDALONE
00008 #include "DataFormats/TrackReco/interface/Track.h"
00009 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
00010 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00011 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00012 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00013 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00014 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
00015 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
00016 #include "DataFormats/Common/interface/RefToPtr.h"
00017 #include "DataFormats/VertexReco/interface/Vertex.h"
00018 #include "RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h"
00019 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00020 #include "TrackingTools/IPTools/interface/IPTools.h"
00021
00022 #endif
00023
00024 using namespace reco;
00025
00026
00027
00028
00029 PFIsolationEstimator::PFIsolationEstimator() :
00030 fisInitialized(kFALSE)
00031 {
00032
00033 }
00034
00035
00036
00037
00038 PFIsolationEstimator::~PFIsolationEstimator()
00039 {
00040
00041 }
00042
00043
00044 void PFIsolationEstimator::initialize( Bool_t bApplyVeto, int iParticleType ) {
00045
00046 setParticleType(iParticleType);
00047
00048
00049 checkClosestZVertex = kTRUE;
00050
00051
00052 setApplyVeto(bApplyVeto);
00053
00054 setDeltaRVetoBarrelPhotons();
00055 setDeltaRVetoBarrelNeutrals();
00056 setDeltaRVetoBarrelCharged();
00057 setDeltaRVetoEndcapPhotons();
00058 setDeltaRVetoEndcapNeutrals();
00059 setDeltaRVetoEndcapCharged();
00060
00061
00062 setRectangleDeltaPhiVetoBarrelPhotons();
00063 setRectangleDeltaPhiVetoBarrelNeutrals();
00064 setRectangleDeltaPhiVetoBarrelCharged();
00065 setRectangleDeltaPhiVetoEndcapPhotons();
00066 setRectangleDeltaPhiVetoEndcapNeutrals();
00067 setRectangleDeltaPhiVetoEndcapCharged();
00068
00069
00070 setRectangleDeltaEtaVetoBarrelPhotons();
00071 setRectangleDeltaEtaVetoBarrelNeutrals();
00072 setRectangleDeltaEtaVetoBarrelCharged();
00073 setRectangleDeltaEtaVetoEndcapPhotons();
00074 setRectangleDeltaEtaVetoEndcapNeutrals();
00075 setRectangleDeltaEtaVetoEndcapCharged();
00076
00077
00078 if(bApplyVeto && iParticleType==kElectron){
00079
00080
00081 setDeltaRVetoBarrel(kTRUE);
00082 setDeltaRVetoEndcap(kTRUE);
00083 setRectangleVetoBarrel(kFALSE);
00084 setRectangleVetoEndcap(kFALSE);
00085 setApplyDzDxyVeto(kFALSE);
00086 setApplyPFPUVeto(kTRUE);
00087 setApplyMissHitPhVeto(kTRUE);
00088
00089 setUseCrystalSize(kFALSE);
00090
00091
00092
00093
00094 setDeltaRVetoEndcapPhotons(0.08);
00095 setDeltaRVetoEndcapCharged(0.015);
00096
00097
00098 setConeSize(0.4);
00099
00100
00101 }else{
00102
00103 setApplyDzDxyVeto(kTRUE);
00104 setApplyPFPUVeto(kTRUE);
00105 setApplyMissHitPhVeto(kFALSE);
00106 setDeltaRVetoBarrel(kTRUE);
00107 setDeltaRVetoEndcap(kTRUE);
00108 setRectangleVetoBarrel(kTRUE);
00109 setRectangleVetoEndcap(kFALSE);
00110 setUseCrystalSize(kTRUE);
00111 setConeSize(0.3);
00112
00113 setDeltaRVetoBarrelPhotons(-1);
00114 setDeltaRVetoBarrelNeutrals(-1);
00115 setDeltaRVetoBarrelCharged(0.02);
00116 setRectangleDeltaPhiVetoBarrelPhotons(1.);
00117 setRectangleDeltaPhiVetoBarrelNeutrals(-1);
00118 setRectangleDeltaPhiVetoBarrelCharged(-1);
00119 setRectangleDeltaEtaVetoBarrelPhotons(0.015);
00120 setRectangleDeltaEtaVetoBarrelNeutrals(-1);
00121 setRectangleDeltaEtaVetoBarrelCharged(-1);
00122
00123 setDeltaRVetoEndcapPhotons(0.07);
00124 setDeltaRVetoEndcapNeutrals(-1);
00125 setDeltaRVetoEndcapCharged(0.02);
00126 setRectangleDeltaPhiVetoEndcapPhotons(-1);
00127 setRectangleDeltaPhiVetoEndcapNeutrals(-1);
00128 setRectangleDeltaPhiVetoEndcapCharged(-1);
00129 setRectangleDeltaEtaVetoEndcapPhotons(-1);
00130 setRectangleDeltaEtaVetoEndcapNeutrals(-1);
00131 setRectangleDeltaEtaVetoEndcapCharged(-1);
00132 setNumberOfCrystalEndcapPhotons(4.);
00133 }
00134
00135
00136 }
00137
00138
00139 void PFIsolationEstimator::initializeElectronIsolation( Bool_t bApplyVeto){
00140 initialize(bApplyVeto,kElectron);
00141 initializeRings(1, fConeSize);
00142
00143
00144
00145
00146
00147
00148
00149
00150
00151
00152 }
00153
00154
00155 void PFIsolationEstimator::initializePhotonIsolation( Bool_t bApplyVeto){
00156 initialize(bApplyVeto,kPhoton);
00157 initializeRings(1, fConeSize);
00158 }
00159
00160
00161
00162 void PFIsolationEstimator::initializeElectronIsolationInRings( Bool_t bApplyVeto, int iNumberOfRings, float fRingSize ){
00163 initialize(bApplyVeto,kElectron);
00164 initializeRings(iNumberOfRings, fRingSize);
00165 }
00166
00167
00168 void PFIsolationEstimator::initializePhotonIsolationInRings( Bool_t bApplyVeto, int iNumberOfRings, float fRingSize ){
00169 initialize(bApplyVeto,kPhoton);
00170 initializeRings(iNumberOfRings, fRingSize);
00171 }
00172
00173
00174
00175 void PFIsolationEstimator::initializeRings(int iNumberOfRings, float fRingSize){
00176
00177 setRingSize(fRingSize);
00178 setNumbersOfRings(iNumberOfRings);
00179
00180 fIsolationInRings.clear();
00181 for(int isoBin =0;isoBin<iNumberOfRings;isoBin++){
00182 float fTemp = 0.0;
00183 fIsolationInRings.push_back(fTemp);
00184
00185 float fTempPhoton = 0.0;
00186 fIsolationInRingsPhoton.push_back(fTempPhoton);
00187
00188 float fTempNeutral = 0.0;
00189 fIsolationInRingsNeutral.push_back(fTempNeutral);
00190
00191 float fTempCharged = 0.0;
00192 fIsolationInRingsCharged.push_back(fTempCharged);
00193
00194 float fTempChargedAll = 0.0;
00195 fIsolationInRingsChargedAll.push_back(fTempChargedAll);
00196
00197 }
00198
00199 fConeSize = fRingSize * (float)iNumberOfRings;
00200
00201 }
00202
00203
00204
00205 float PFIsolationEstimator::fGetIsolation(const reco::PFCandidate * pfCandidate, const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices) {
00206
00207 fGetIsolationInRings( pfCandidate, pfParticlesColl, vtx, vertices);
00208 refSC = SuperClusterRef();
00209 fIsolation = fIsolationInRings[0];
00210
00211 return fIsolation;
00212 }
00213
00214
00215
00216 vector<float > PFIsolationEstimator::fGetIsolationInRings(const reco::PFCandidate * pfCandidate, const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices) {
00217
00218 int isoBin;
00219
00220
00221 for(isoBin =0;isoBin<iNumberOfRings;isoBin++){
00222 fIsolationInRings[isoBin]=0.;
00223 fIsolationInRingsPhoton[isoBin]=0.;
00224 fIsolationInRingsNeutral[isoBin]=0.;
00225 fIsolationInRingsCharged[isoBin]=0.;
00226 fIsolationInRingsChargedAll[isoBin]=0.;
00227 }
00228
00229
00230
00231 fEta = pfCandidate->eta();
00232 fPhi = pfCandidate->phi();
00233 fPt = pfCandidate->pt();
00234 fVx = pfCandidate->vx();
00235 fVy = pfCandidate->vy();
00236 fVz = pfCandidate->vz();
00237
00238 pivotInBarrel = fabs(pfCandidate->positionAtECALEntrance().eta())<1.479;
00239
00240 for(unsigned iPF=0; iPF<pfParticlesColl->size(); iPF++) {
00241
00242 const reco::PFCandidate& pfParticle= (*pfParticlesColl)[iPF];
00243
00244 if(&pfParticle==(pfCandidate))
00245 continue;
00246
00247 if(pfParticle.pdgId()==22){
00248
00249 if(isPhotonParticleVetoed( &pfParticle)>=0.){
00250 isoBin = (int)(fDeltaR/fRingSize);
00251 fIsolationInRingsPhoton[isoBin] = fIsolationInRingsPhoton[isoBin] + pfParticle.pt();
00252 }
00253
00254 }else if(abs(pfParticle.pdgId())==130){
00255
00256 if(isNeutralParticleVetoed( &pfParticle)>=0.){
00257 isoBin = (int)(fDeltaR/fRingSize);
00258 fIsolationInRingsNeutral[isoBin] = fIsolationInRingsNeutral[isoBin] + pfParticle.pt();
00259 }
00260
00261
00262
00263 }else if(abs(pfParticle.pdgId()) == 211){
00264 if(isChargedParticleVetoed( &pfParticle, vtx, vertices)>=0.){
00265 isoBin = (int)(fDeltaR/fRingSize);
00266 fIsolationInRingsCharged[isoBin] = fIsolationInRingsCharged[isoBin] + pfParticle.pt();
00267 }
00268
00269 }
00270 }
00271
00272
00273 for(int isoBin =0;isoBin<iNumberOfRings;isoBin++){
00274 fIsolationInRings[isoBin]= fIsolationInRingsPhoton[isoBin]+ fIsolationInRingsNeutral[isoBin] + fIsolationInRingsCharged[isoBin];
00275 }
00276
00277 return fIsolationInRings;
00278 }
00279
00280
00281
00282 float PFIsolationEstimator::fGetIsolation(const reco::Photon * photon, const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices) {
00283
00284 fGetIsolationInRings( photon, pfParticlesColl, vtx, vertices);
00285 fIsolation = fIsolationInRings[0];
00286
00287 return fIsolation;
00288 }
00289
00290
00291
00292 vector<float > PFIsolationEstimator::fGetIsolationInRings(const reco::Photon * photon, const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices) {
00293
00294 int isoBin;
00295
00296 for(isoBin =0;isoBin<iNumberOfRings;isoBin++){
00297 fIsolationInRings[isoBin]=0.;
00298 fIsolationInRingsPhoton[isoBin]=0.;
00299 fIsolationInRingsNeutral[isoBin]=0.;
00300 fIsolationInRingsCharged[isoBin]=0.;
00301 fIsolationInRingsChargedAll[isoBin]=0.;
00302 }
00303
00304 iMissHits = 0;
00305
00306 refSC = photon->superCluster();
00307 pivotInBarrel = fabs((refSC->position().eta()))<1.479;
00308
00309 for(unsigned iPF=0; iPF<pfParticlesColl->size(); iPF++) {
00310
00311 const reco::PFCandidate& pfParticle= (*pfParticlesColl)[iPF];
00312
00313 if (pfParticle.superClusterRef().isNonnull() &&
00314 photon->superCluster().isNonnull() &&
00315 pfParticle.superClusterRef() == photon->superCluster())
00316 continue;
00317
00318 if(pfParticle.pdgId()==22){
00319
00320
00321 math::XYZVector direction = math::XYZVector(photon->superCluster()->x() - pfParticle.vx(),
00322 photon->superCluster()->y() - pfParticle.vy(),
00323 photon->superCluster()->z() - pfParticle.vz());
00324
00325 fEta = direction.Eta();
00326 fPhi = direction.Phi();
00327 fVx = pfParticle.vx();
00328 fVy = pfParticle.vy();
00329 fVz = pfParticle.vz();
00330
00331 if(isPhotonParticleVetoed(&pfParticle)>=0.){
00332 isoBin = (int)(fDeltaR/fRingSize);
00333 fIsolationInRingsPhoton[isoBin] = fIsolationInRingsPhoton[isoBin] + pfParticle.pt();
00334 }
00335
00336 }else if(abs(pfParticle.pdgId())==130){
00337
00338
00339 math::XYZVector direction = math::XYZVector(photon->superCluster()->x() - pfParticle.vx(),
00340 photon->superCluster()->y() - pfParticle.vy(),
00341 photon->superCluster()->z() - pfParticle.vz());
00342
00343 fEta = direction.Eta();
00344 fPhi = direction.Phi();
00345 fVx = pfParticle.vx();
00346 fVy = pfParticle.vy();
00347 fVz = pfParticle.vz();
00348
00349 if(isNeutralParticleVetoed( &pfParticle)>=0.){
00350 isoBin = (int)(fDeltaR/fRingSize);
00351 fIsolationInRingsNeutral[isoBin] = fIsolationInRingsNeutral[isoBin] + pfParticle.pt();
00352 }
00353
00354
00355 }else if(abs(pfParticle.pdgId()) == 211){
00356
00357
00358 math::XYZVector direction = math::XYZVector(photon->superCluster()->x() - (*vtx).x(),
00359 photon->superCluster()->y() - (*vtx).y(),
00360 photon->superCluster()->z() - (*vtx).z());
00361
00362 fEta = direction.Eta();
00363 fPhi = direction.Phi();
00364 fVx = (*vtx).x();
00365 fVy = (*vtx).y();
00366 fVz = (*vtx).z();
00367
00368 if(isChargedParticleVetoed( &pfParticle, vtx, vertices)>=0.){
00369 isoBin = (int)(fDeltaR/fRingSize);
00370 fIsolationInRingsCharged[isoBin] = fIsolationInRingsCharged[isoBin] + pfParticle.pt();
00371 }
00372
00373 }
00374 }
00375
00376
00377 for(int isoBin =0;isoBin<iNumberOfRings;isoBin++){
00378 fIsolationInRings[isoBin]= fIsolationInRingsPhoton[isoBin]+ fIsolationInRingsNeutral[isoBin] + fIsolationInRingsCharged[isoBin];
00379 }
00380
00381 return fIsolationInRings;
00382 }
00383
00384
00385
00386
00387 float PFIsolationEstimator::fGetIsolation(const reco::GsfElectron * electron, const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices) {
00388
00389 fGetIsolationInRings( electron, pfParticlesColl, vtx, vertices);
00390 fIsolation = fIsolationInRings[0];
00391
00392 return fIsolation;
00393 }
00394
00395
00396
00397 vector<float > PFIsolationEstimator::fGetIsolationInRings(const reco::GsfElectron * electron, const reco::PFCandidateCollection* pfParticlesColl,reco::VertexRef vtx, edm::Handle< reco::VertexCollection > vertices) {
00398
00399 int isoBin;
00400
00401 for(isoBin =0;isoBin<iNumberOfRings;isoBin++){
00402 fIsolationInRings[isoBin]=0.;
00403 fIsolationInRingsPhoton[isoBin]=0.;
00404 fIsolationInRingsNeutral[isoBin]=0.;
00405 fIsolationInRingsCharged[isoBin]=0.;
00406 fIsolationInRingsChargedAll[isoBin]=0.;
00407 }
00408
00409
00410
00411
00412 fEta = electron->eta();
00413 fPhi = electron->phi();
00414 fPt = electron->pt();
00415 fVx = electron->vx();
00416 fVy = electron->vy();
00417 fVz = electron->vz();
00418 iMissHits = electron->gsfTrack()->trackerExpectedHitsInner().numberOfHits();
00419
00420
00421 refSC = electron->superCluster();
00422 pivotInBarrel = fabs((refSC->position().eta()))<1.479;
00423
00424 for(unsigned iPF=0; iPF<pfParticlesColl->size(); iPF++) {
00425
00426 const reco::PFCandidate& pfParticle= (*pfParticlesColl)[iPF];
00427
00428
00429 if(pfParticle.pdgId()==22){
00430
00431 if(isPhotonParticleVetoed(&pfParticle)>=0.){
00432 isoBin = (int)(fDeltaR/fRingSize);
00433 fIsolationInRingsPhoton[isoBin] = fIsolationInRingsPhoton[isoBin] + pfParticle.pt();
00434
00435 }
00436
00437 }else if(abs(pfParticle.pdgId())==130){
00438
00439 if(isNeutralParticleVetoed( &pfParticle)>=0.){
00440 isoBin = (int)(fDeltaR/fRingSize);
00441 fIsolationInRingsNeutral[isoBin] = fIsolationInRingsNeutral[isoBin] + pfParticle.pt();
00442 }
00443
00444
00445 }else if(abs(pfParticle.pdgId()) == 211){
00446 if(isChargedParticleVetoed( &pfParticle, vtx, vertices)>=0.){
00447 isoBin = (int)(fDeltaR/fRingSize);
00448
00449 fIsolationInRingsCharged[isoBin] = fIsolationInRingsCharged[isoBin] + pfParticle.pt();
00450 }
00451
00452 }
00453 }
00454
00455
00456 for(int isoBin =0;isoBin<iNumberOfRings;isoBin++){
00457 fIsolationInRings[isoBin]= fIsolationInRingsPhoton[isoBin]+ fIsolationInRingsNeutral[isoBin] + fIsolationInRingsCharged[isoBin];
00458 }
00459
00460 return fIsolationInRings;
00461 }
00462
00463
00464
00465 float PFIsolationEstimator::isPhotonParticleVetoed( const reco::PFCandidate* pfIsoCand ){
00466
00467
00468 fDeltaR = deltaR(fEta,fPhi,pfIsoCand->eta(),pfIsoCand->phi());
00469
00470 if(fDeltaR > fConeSize)
00471 return -999.;
00472
00473 fDeltaPhi = deltaPhi(fPhi,pfIsoCand->phi());
00474 fDeltaEta = fEta-pfIsoCand->eta();
00475
00476 if(!bApplyVeto)
00477 return fDeltaR;
00478
00479
00480
00481
00482 if(bApplyMissHitPhVeto) {
00483 if(iMissHits > 0)
00484 if(pfIsoCand->mva_nothing_gamma() > 0.99) {
00485 if(pfIsoCand->superClusterRef().isNonnull() && refSC.isNonnull()) {
00486 if(pfIsoCand->superClusterRef() == refSC)
00487 return -999.;
00488 }
00489 }
00490 }
00491
00492 if(pivotInBarrel){
00493 if(bDeltaRVetoBarrel){
00494 if(fDeltaR < fDeltaRVetoBarrelPhotons)
00495 return -999.;
00496 }
00497
00498 if(bRectangleVetoBarrel){
00499 if(abs(fDeltaEta) < fRectangleDeltaEtaVetoBarrelPhotons && abs(fDeltaPhi) < fRectangleDeltaPhiVetoBarrelPhotons){
00500 return -999.;
00501 }
00502 }
00503 }else{
00504 if (bUseCrystalSize == true) {
00505 fDeltaRVetoEndcapPhotons = 0.00864*fabs(sinh(refSC->position().eta()))*fNumberOfCrystalEndcapPhotons;
00506 }
00507
00508 if(bDeltaRVetoEndcap){
00509 if(fDeltaR < fDeltaRVetoEndcapPhotons)
00510 return -999.;
00511 }
00512 if(bRectangleVetoEndcap){
00513 if(abs(fDeltaEta) < fRectangleDeltaEtaVetoEndcapPhotons && abs(fDeltaPhi) < fRectangleDeltaPhiVetoEndcapPhotons){
00514 return -999.;
00515 }
00516 }
00517 }
00518
00519 return fDeltaR;
00520 }
00521
00522
00523 float PFIsolationEstimator::isNeutralParticleVetoed( const reco::PFCandidate* pfIsoCand ){
00524
00525 fDeltaR = deltaR(fEta,fPhi,pfIsoCand->eta(),pfIsoCand->phi());
00526
00527 if(fDeltaR > fConeSize)
00528 return -999;
00529
00530 fDeltaPhi = deltaPhi(fPhi,pfIsoCand->phi());
00531 fDeltaEta = fEta-pfIsoCand->eta();
00532
00533 if(!bApplyVeto)
00534 return fDeltaR;
00535
00536
00537
00538 if(pivotInBarrel){
00539 if(!bDeltaRVetoBarrel&&!bRectangleVetoBarrel){
00540 return fDeltaR;
00541 }
00542
00543 if(bDeltaRVetoBarrel){
00544 if(fDeltaR < fDeltaRVetoBarrelNeutrals)
00545 return -999.;
00546 }
00547 if(bRectangleVetoBarrel){
00548 if(abs(fDeltaEta) < fRectangleDeltaEtaVetoBarrelNeutrals && abs(fDeltaPhi) < fRectangleDeltaPhiVetoBarrelNeutrals){
00549 return -999.;
00550 }
00551 }
00552
00553 }else{
00554 if(!bDeltaRVetoEndcap&&!bRectangleVetoEndcap){
00555 return fDeltaR;
00556 }
00557 if(bDeltaRVetoEndcap){
00558 if(fDeltaR < fDeltaRVetoEndcapNeutrals)
00559 return -999.;
00560 }
00561 if(bRectangleVetoEndcap){
00562 if(abs(fDeltaEta) < fRectangleDeltaEtaVetoEndcapNeutrals && abs(fDeltaPhi) < fRectangleDeltaPhiVetoEndcapNeutrals){
00563 return -999.;
00564 }
00565 }
00566 }
00567
00568 return fDeltaR;
00569 }
00570
00571
00572
00573 float PFIsolationEstimator::isChargedParticleVetoed(const reco::PFCandidate* pfIsoCand, edm::Handle< reco::VertexCollection > vertices ){
00574
00575
00576 return -999;
00577 }
00578
00579
00580 float PFIsolationEstimator::isChargedParticleVetoed(const reco::PFCandidate* pfIsoCand,reco::VertexRef vtxMain, edm::Handle< reco::VertexCollection > vertices ){
00581
00582 VertexRef vtx = chargedHadronVertex(vertices, *pfIsoCand );
00583 if(vtx.isNull())
00584 return -999.;
00585
00586
00587
00588 float fVtxMainZ = (*vtxMain).z();
00589
00590 if(bApplyPFPUVeto) {
00591 if(vtx != vtxMain)
00592 return -999.;
00593 }
00594
00595
00596 if(bApplyDzDxyVeto) {
00597 if(iParticleType==kPhoton){
00598
00599 float dz = fabs( pfIsoCand->trackRef()->dz( (*vtxMain).position() ) );
00600 if (dz > 0.2)
00601 return -999.;
00602
00603 double dxy = pfIsoCand->trackRef()->dxy( (*vtxMain).position() );
00604 if (fabs(dxy) > 0.1)
00605 return -999.;
00606
00607
00608
00609
00610
00611
00612
00613
00614
00615
00616
00617
00618 }else{
00619
00620
00621 float dz = fabs(vtx->z() - fVtxMainZ);
00622 if (dz > 1.)
00623 return -999.;
00624
00625 double dxy = ( -(vtx->x() - fVx)*pfIsoCand->py() + (vtx->y() - fVy)*pfIsoCand->px()) / pfIsoCand->pt();
00626 if(fabs(dxy) > 0.1)
00627 return -999.;
00628 }
00629 }
00630
00631 fDeltaR = deltaR(pfIsoCand->eta(),pfIsoCand->phi(),fEta,fPhi);
00632
00633 if(fDeltaR > fConeSize)
00634 return -999.;
00635
00636 fDeltaPhi = deltaPhi(fPhi,pfIsoCand->phi());
00637 fDeltaEta = fEta-pfIsoCand->eta();
00638
00639
00640
00641
00642
00643
00644
00645 if(!bApplyVeto)
00646 return fDeltaR;
00647
00648
00649
00650 if(pivotInBarrel){
00651 if(!bDeltaRVetoBarrel&&!bRectangleVetoBarrel){
00652 return fDeltaR;
00653 }
00654
00655 if(bDeltaRVetoBarrel){
00656 if(fDeltaR < fDeltaRVetoBarrelCharged)
00657 return -999.;
00658 }
00659 if(bRectangleVetoBarrel){
00660 if(abs(fDeltaEta) < fRectangleDeltaEtaVetoBarrelCharged && abs(fDeltaPhi) < fRectangleDeltaPhiVetoBarrelCharged){
00661 return -999.;
00662 }
00663 }
00664
00665 }else{
00666 if(!bDeltaRVetoEndcap&&!bRectangleVetoEndcap){
00667 return fDeltaR;
00668 }
00669 if(bDeltaRVetoEndcap){
00670 if(fDeltaR < fDeltaRVetoEndcapCharged)
00671 return -999.;
00672 }
00673 if(bRectangleVetoEndcap){
00674 if(abs(fDeltaEta) < fRectangleDeltaEtaVetoEndcapCharged && abs(fDeltaPhi) < fRectangleDeltaPhiVetoEndcapCharged){
00675 return -999.;
00676 }
00677 }
00678 }
00679
00680
00681
00682
00683 return fDeltaR;
00684 }
00685
00686
00687
00688 VertexRef PFIsolationEstimator::chargedHadronVertex( edm::Handle< reco::VertexCollection > verticesColl, const reco::PFCandidate& pfcand ){
00689
00690
00691
00692 reco::TrackBaseRef trackBaseRef( pfcand.trackRef() );
00693
00694 size_t iVertex = 0;
00695 unsigned index=0;
00696 unsigned nFoundVertex = 0;
00697
00698 float bestweight=0;
00699
00700 const reco::VertexCollection& vertices = *(verticesColl.product());
00701
00702 for( reco::VertexCollection::const_iterator iv=vertices.begin(); iv!=vertices.end(); ++iv, ++index) {
00703
00704 const reco::Vertex& vtx = *iv;
00705
00706
00707 for(reco::Vertex::trackRef_iterator iTrack=vtx.tracks_begin();iTrack!=vtx.tracks_end(); ++iTrack) {
00708
00709 const reco::TrackBaseRef& baseRef = *iTrack;
00710
00711
00712
00713 if(baseRef == trackBaseRef ) {
00714 float w = vtx.trackWeight(baseRef);
00715
00716 if (w > bestweight){
00717 bestweight=w;
00718 iVertex=index;
00719 nFoundVertex++;
00720 }
00721 }
00722 }
00723
00724 }
00725
00726
00727
00728 if (nFoundVertex>0){
00729 if (nFoundVertex!=1)
00730 edm::LogWarning("TrackOnTwoVertex")<<"a track is shared by at least two verteces. Used to be an assert";
00731 return VertexRef( verticesColl, iVertex);
00732 }
00733
00734
00735
00736 if ( checkClosestZVertex ) {
00737
00738 double dzmin = 10000.;
00739 double ztrack = pfcand.vertex().z();
00740 bool foundVertex = false;
00741 index = 0;
00742 for( reco::VertexCollection::const_iterator iv=vertices.begin(); iv!=vertices.end(); ++iv, ++index) {
00743
00744 double dz = fabs(ztrack - iv->z());
00745 if(dz<dzmin) {
00746 dzmin = dz;
00747 iVertex = index;
00748 foundVertex = true;
00749 }
00750 }
00751
00752 if( foundVertex )
00753 return VertexRef( verticesColl, iVertex);
00754
00755 }
00756
00757 return VertexRef( );
00758 }
00759
00760
00761
00762 int PFIsolationEstimator::matchPFObject(const reco::Photon* photon, const reco::PFCandidateCollection* Candidates ){
00763
00764 Int_t iMatch = -1;
00765
00766 int i=0;
00767 for(reco::PFCandidateCollection::const_iterator iPF=Candidates->begin();iPF !=Candidates->end();iPF++){
00768 const reco::PFCandidate& pfParticle = (*iPF);
00769
00770 if((((pfParticle.pdgId()==22 ) || TMath::Abs(pfParticle.pdgId())==11) )){
00771
00772 if(pfParticle.superClusterRef()==photon->superCluster())
00773 iMatch= i;
00774
00775 }
00776
00777 i++;
00778 }
00779
00780
00781
00782
00783
00784
00785
00786
00787
00788
00789
00790
00791
00792
00793
00794
00795
00796
00797
00798
00799
00800 return iMatch;
00801
00802 }
00803
00804
00805
00806
00807
00808 int PFIsolationEstimator::matchPFObject(const reco::GsfElectron* electron, const reco::PFCandidateCollection* Candidates ){
00809
00810 Int_t iMatch = -1;
00811
00812 int i=0;
00813 for(reco::PFCandidateCollection::const_iterator iPF=Candidates->begin();iPF !=Candidates->end();iPF++){
00814 const reco::PFCandidate& pfParticle = (*iPF);
00815
00816 if((((pfParticle.pdgId()==22 ) || TMath::Abs(pfParticle.pdgId())==11) )){
00817
00818 if(pfParticle.superClusterRef()==electron->superCluster())
00819 iMatch= i;
00820
00821 }
00822
00823 i++;
00824 }
00825
00826 if(iMatch == -1){
00827 i=0;
00828 float fPt = -1;
00829 for(reco::PFCandidateCollection::const_iterator iPF=Candidates->begin();iPF !=Candidates->end();iPF++){
00830 const reco::PFCandidate& pfParticle = (*iPF);
00831 if((((pfParticle.pdgId()==22 ) || TMath::Abs(pfParticle.pdgId())==11) )){
00832 if(pfParticle.pt()>fPt){
00833 fDeltaR = deltaR(pfParticle.eta(),pfParticle.phi(),electron->eta(),electron->phi());
00834 if(fDeltaR<0.1){
00835 iMatch = i;
00836 fPt = pfParticle.pt();
00837 }
00838 }
00839 }
00840 i++;
00841 }
00842 }
00843
00844 return iMatch;
00845
00846 }
00847