00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #include <memory>
00013
00014
00015 #include "RecoParticleFlow/PFTracking/interface/PFElecTkProducer.h"
00016 #include "RecoParticleFlow/PFTracking/interface/PFTrackTransformer.h"
00017 #include "RecoParticleFlow/PFTracking/interface/ConvBremPFTrackFinder.h"
00018 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
00019 #include "MagneticField/Engine/interface/MagneticField.h"
00020 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00021 #include "FWCore/Framework/interface/ESHandle.h"
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00024 #include "DataFormats/TrajectorySeed/interface/TrajectorySeed.h"
00025 #include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
00026 #include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
00027 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
00028 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00029 #include "MagneticField/Engine/interface/MagneticField.h"
00030 #include "DataFormats/ParticleFlowReco/interface/PFClusterFwd.h"
00031 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
00032 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
00033 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00034 #include "DataFormats/VertexReco/interface/Vertex.h"
00035 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00036 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedTrackerVertex.h"
00037 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertexFwd.h"
00038 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertex.h"
00039 #include "DataFormats/ParticleFlowReco/interface/PFConversionFwd.h"
00040 #include "DataFormats/ParticleFlowReco/interface/PFConversion.h"
00041 #include "DataFormats/ParticleFlowReco/interface/PFV0Fwd.h"
00042 #include "DataFormats/ParticleFlowReco/interface/PFV0.h"
00043 #include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
00044 #include "RecoParticleFlow/PFClusterTools/interface/ClusterClusterMapping.h"
00045
00046 #include "TMath.h"
00047 using namespace std;
00048 using namespace edm;
00049 using namespace reco;
00050 PFElecTkProducer::PFElecTkProducer(const ParameterSet& iConfig):
00051 conf_(iConfig),
00052 pfTransformer_(0),
00053 convBremFinder_(0)
00054 {
00055 LogInfo("PFElecTkProducer")<<"PFElecTkProducer started";
00056
00057 gsfTrackLabel_ = iConfig.getParameter<InputTag>
00058 ("GsfTrackModuleLabel");
00059
00060 pfTrackLabel_ = iConfig.getParameter<InputTag>
00061 ("PFRecTrackLabel");
00062
00063 primVtxLabel_ = iConfig.getParameter<InputTag>
00064 ("PrimaryVertexLabel");
00065
00066 pfEcalClusters_ = iConfig.getParameter<InputTag>
00067 ("PFEcalClusters");
00068
00069 pfNuclear_ = iConfig.getParameter<InputTag>
00070 ("PFNuclear");
00071
00072 pfConv_ = iConfig.getParameter<InputTag>
00073 ("PFConversions");
00074
00075 pfV0_ = iConfig.getParameter<InputTag>
00076 ("PFV0");
00077
00078 useNuclear_ = iConfig.getParameter<bool>("useNuclear");
00079 useConversions_ = iConfig.getParameter<bool>("useConversions");
00080 useV0_ = iConfig.getParameter<bool>("useV0");
00081 debugGsfCleaning_ = iConfig.getParameter<bool>("debugGsfCleaning");
00082
00083 produces<GsfPFRecTrackCollection>();
00084 produces<GsfPFRecTrackCollection>( "Secondary" ).setBranchAlias( "secondary" );
00085
00086
00087 trajinev_ = iConfig.getParameter<bool>("TrajInEvents");
00088 modemomentum_ = iConfig.getParameter<bool>("ModeMomentum");
00089 applySel_ = iConfig.getParameter<bool>("applyEGSelection");
00090 applyGsfClean_ = iConfig.getParameter<bool>("applyGsfTrackCleaning");
00091 applyAngularGsfClean_ = iConfig.getParameter<bool>("applyAlsoGsfAngularCleaning");
00092 detaCutGsfClean_ = iConfig.getParameter<double>("maxDEtaGsfAngularCleaning");
00093 dphiCutGsfClean_ = iConfig.getParameter<double>("maxDPhiBremTangGsfAngularCleaning");
00094 useFifthStepForTrackDriven_ = iConfig.getParameter<bool>("useFifthStepForTrackerDrivenGsf");
00095 useFifthStepForEcalDriven_ = iConfig.getParameter<bool>("useFifthStepForEcalDrivenGsf");
00096 maxPtConvReco_ = iConfig.getParameter<double>("MaxConvBremRecoPT");
00097 detaGsfSC_ = iConfig.getParameter<double>("MinDEtaGsfSC");
00098 dphiGsfSC_ = iConfig.getParameter<double>("MinDPhiGsfSC");
00099 SCEne_ = iConfig.getParameter<double>("MinSCEnergy");
00100
00101
00102 useConvBremFinder_ = iConfig.getParameter<bool>("useConvBremFinder");
00103 mvaConvBremFinderID_
00104 = iConfig.getParameter<double>("pf_convBremFinderID_mvaCut");
00105
00106 string mvaWeightFileConvBrem
00107 = iConfig.getParameter<string>("pf_convBremFinderID_mvaWeightFile");
00108
00109
00110 if(useConvBremFinder_)
00111 path_mvaWeightFileConvBrem_ = edm::FileInPath ( mvaWeightFileConvBrem.c_str() ).fullPath();
00112
00113 }
00114
00115
00116 PFElecTkProducer::~PFElecTkProducer()
00117 {
00118
00119 delete pfTransformer_;
00120 delete convBremFinder_;
00121 }
00122
00123
00124
00125
00126
00127
00128
00129 void
00130 PFElecTkProducer::produce(Event& iEvent, const EventSetup& iSetup)
00131 {
00132 LogDebug("PFElecTkProducer")<<"START event: "<<iEvent.id().event()
00133 <<" in run "<<iEvent.id().run();
00134
00135
00136 auto_ptr< GsfPFRecTrackCollection >
00137 gsfPFRecTrackCollection(new GsfPFRecTrackCollection);
00138
00139
00140 auto_ptr< GsfPFRecTrackCollection >
00141 gsfPFRecTrackCollectionSecondary(new GsfPFRecTrackCollection);
00142
00143
00144 Handle<GsfTrackCollection> gsftrackscoll;
00145 iEvent.getByLabel(gsfTrackLabel_,gsftrackscoll);
00146
00147
00148 Handle<vector<Trajectory> > TrajectoryCollection;
00149
00150
00151 Handle<PFRecTrackCollection> thePfRecTrackCollection;
00152 iEvent.getByLabel(pfTrackLabel_,thePfRecTrackCollection);
00153 const PFRecTrackCollection& PfRTkColl = *(thePfRecTrackCollection.product());
00154
00155
00156 Handle<PFClusterCollection> theECPfClustCollection;
00157 iEvent.getByLabel(pfEcalClusters_,theECPfClustCollection);
00158 const PFClusterCollection& theEcalClusters = *(theECPfClustCollection.product());
00159
00160
00161 Handle<reco::VertexCollection> thePrimaryVertexColl;
00162 iEvent.getByLabel(primVtxLabel_,thePrimaryVertexColl);
00163
00164
00165
00166
00167 Handle< reco::PFDisplacedTrackerVertexCollection > pfNuclears;
00168 if( useNuclear_ ) {
00169 bool found = iEvent.getByLabel(pfNuclear_, pfNuclears);
00170
00171
00172 if(!found )
00173 LogError("PFElecTkProducer")<<" cannot get PFNuclear : "
00174 << pfNuclear_
00175 << " please set useNuclear=False in RecoParticleFlow/PFTracking/python/pfTrackElec_cfi.py" << endl;
00176 }
00177
00178
00179 Handle< reco::PFConversionCollection > pfConversions;
00180 if( useConversions_ ) {
00181 bool found = iEvent.getByLabel(pfConv_,pfConversions);
00182 if(!found )
00183 LogError("PFElecTkProducer")<<" cannot get PFConversions : "
00184 << pfConv_
00185 << " please set useConversions=False in RecoParticleFlow/PFTracking/python/pfTrackElec_cfi.py" << endl;
00186 }
00187
00188
00189 Handle< reco::PFV0Collection > pfV0;
00190 if( useV0_ ) {
00191 bool found = iEvent.getByLabel(pfV0_, pfV0);
00192
00193 if(!found )
00194 LogError("PFElecTkProducer")<<" cannot get PFV0 : "
00195 << pfV0_
00196 << " please set useV0=False RecoParticleFlow/PFTracking/python/pfTrackElec_cfi.py" << endl;
00197 }
00198
00199
00200
00201 GsfTrackCollection gsftracks = *(gsftrackscoll.product());
00202 vector<Trajectory> tjvec(0);
00203 if (trajinev_){
00204 bool foundTraj = iEvent.getByLabel(gsfTrackLabel_,TrajectoryCollection);
00205 if(!foundTraj)
00206 LogError("PFElecTkProducer")
00207 <<" cannot get Trajectories of : "
00208 << gsfTrackLabel_
00209 << " please set TrajInEvents = False in RecoParticleFlow/PFTracking/python/pfTrackElec_cfi.py" << endl;
00210
00211 tjvec= *(TrajectoryCollection.product());
00212 }
00213
00214
00215 vector<reco::GsfPFRecTrack> selGsfPFRecTracks;
00216 vector<reco::GsfPFRecTrack> primaryGsfPFRecTracks;
00217 std::map<unsigned int, std::vector<reco::GsfPFRecTrack> > GsfPFMap;
00218
00219
00220 for (unsigned int igsf=0; igsf<gsftracks.size();igsf++) {
00221
00222 GsfTrackRef trackRef(gsftrackscoll, igsf);
00223
00224 int kf_ind=FindPfRef(PfRTkColl,gsftracks[igsf],false);
00225
00226 if (kf_ind>=0) {
00227
00228 PFRecTrackRef kf_ref(thePfRecTrackCollection,
00229 kf_ind);
00230
00231
00232 if( useFifthStepForEcalDriven_ == false
00233 || useFifthStepForTrackDriven_ == false) {
00234 bool isFifthStepTrack = isFifthStep(kf_ref);
00235 bool isEcalDriven = true;
00236 bool isTrackerDriven = true;
00237
00238 if (&(*trackRef->seedRef())==0) {
00239 isEcalDriven = false;
00240 isTrackerDriven = false;
00241 }
00242 else {
00243 ElectronSeedRef SeedRef= trackRef->extra()->seedRef().castTo<ElectronSeedRef>();
00244 if(SeedRef->caloCluster().isNull())
00245 isEcalDriven = false;
00246 if(SeedRef->ctfTrack().isNull())
00247 isTrackerDriven = false;
00248 }
00249
00250 if(isFifthStepTrack &&
00251 isEcalDriven &&
00252 isTrackerDriven == false &&
00253 useFifthStepForEcalDriven_ == false) {
00254 continue;
00255 }
00256
00257 if(isFifthStepTrack &&
00258 isTrackerDriven &&
00259 isEcalDriven == false &&
00260 useFifthStepForTrackDriven_ == false) {
00261 continue;
00262 }
00263
00264 if(isFifthStepTrack &&
00265 isTrackerDriven &&
00266 isEcalDriven &&
00267 useFifthStepForTrackDriven_ == false &&
00268 useFifthStepForEcalDriven_ == false) {
00269 continue;
00270 }
00271 }
00272
00273 pftrack_=GsfPFRecTrack( gsftracks[igsf].charge(),
00274 reco::PFRecTrack::GSF,
00275 igsf, trackRef,
00276 kf_ref);
00277 } else {
00278 PFRecTrackRef dummyRef;
00279 pftrack_=GsfPFRecTrack( gsftracks[igsf].charge(),
00280 reco::PFRecTrack::GSF,
00281 igsf, trackRef,
00282 dummyRef);
00283 }
00284
00285
00286 bool validgsfbrem = false;
00287 if(trajinev_) {
00288 validgsfbrem = pfTransformer_->addPointsAndBrems(pftrack_,
00289 gsftracks[igsf],
00290 tjvec[igsf],
00291 modemomentum_);
00292 } else {
00293 validgsfbrem = pfTransformer_->addPointsAndBrems(pftrack_,
00294 gsftracks[igsf],
00295 mtsTransform_);
00296 }
00297
00298 bool passSel = true;
00299 if(applySel_)
00300 passSel = applySelection(gsftracks[igsf]);
00301
00302 if(validgsfbrem && passSel)
00303 selGsfPFRecTracks.push_back(pftrack_);
00304 }
00305
00306
00307 unsigned int count_primary = 0;
00308 if(selGsfPFRecTracks.size() > 0) {
00309 for(unsigned int ipfgsf=0; ipfgsf<selGsfPFRecTracks.size();ipfgsf++) {
00310
00311 vector<unsigned int> secondaries(0);
00312 secondaries.clear();
00313 bool keepGsf = true;
00314
00315 if(applyGsfClean_) {
00316 keepGsf = resolveGsfTracks(selGsfPFRecTracks,ipfgsf,secondaries,theEcalClusters);
00317 }
00318
00319
00320 if(keepGsf == true) {
00321
00322
00323 if(convBremFinder_->foundConvBremPFRecTrack(thePfRecTrackCollection,thePrimaryVertexColl,
00324 pfNuclears,pfConversions,pfV0,
00325 useNuclear_,useConversions_,useV0_,
00326 theEcalClusters,selGsfPFRecTracks[ipfgsf])) {
00327 const vector<PFRecTrackRef>& convBremPFRecTracks(convBremFinder_->getConvBremPFRecTracks());
00328 for(unsigned int ii = 0; ii<convBremPFRecTracks.size(); ii++) {
00329 selGsfPFRecTracks[ipfgsf].addConvBremPFRecTrackRef(convBremPFRecTracks[ii]);
00330 }
00331 }
00332
00333
00334
00335 primaryGsfPFRecTracks.push_back(selGsfPFRecTracks[ipfgsf]);
00336
00337
00338
00339
00340
00341 unsigned int primGsfIndex = selGsfPFRecTracks[ipfgsf].trackId();
00342 vector<reco::GsfPFRecTrack> trueGsfPFRecTracks;
00343 if(secondaries.size() > 0) {
00344
00345 for(unsigned int isecpfgsf=0; isecpfgsf<secondaries.size();isecpfgsf++) {
00346
00347 PFRecTrackRef refsecKF = selGsfPFRecTracks[(secondaries[isecpfgsf])].kfPFRecTrackRef();
00348
00349 unsigned int secGsfIndex = selGsfPFRecTracks[(secondaries[isecpfgsf])].trackId();
00350 GsfTrackRef secGsfRef = selGsfPFRecTracks[(secondaries[isecpfgsf])].gsfTrackRef();
00351
00352 if(refsecKF.isNonnull()) {
00353
00354 secpftrack_= GsfPFRecTrack( gsftracks[secGsfIndex].charge(),
00355 reco::PFRecTrack::GSF,
00356 primGsfIndex, secGsfRef,
00357 refsecKF);
00358 }
00359 else{
00360 PFRecTrackRef dummyRef;
00361
00362 secpftrack_= GsfPFRecTrack( gsftracks[secGsfIndex].charge(),
00363 reco::PFRecTrack::GSF,
00364 primGsfIndex, secGsfRef,
00365 dummyRef);
00366 }
00367
00368 bool validgsfbrem = false;
00369 if(trajinev_) {
00370 validgsfbrem = pfTransformer_->addPointsAndBrems(secpftrack_,
00371 gsftracks[secGsfIndex],
00372 tjvec[secGsfIndex],
00373 modemomentum_);
00374 } else {
00375 validgsfbrem = pfTransformer_->addPointsAndBrems(secpftrack_,
00376 gsftracks[secGsfIndex],
00377 mtsTransform_);
00378 }
00379
00380 if(validgsfbrem) {
00381 gsfPFRecTrackCollectionSecondary->push_back(secpftrack_);
00382 trueGsfPFRecTracks.push_back(secpftrack_);
00383 }
00384 }
00385 }
00386 GsfPFMap.insert(pair<unsigned int,std::vector<reco::GsfPFRecTrack> >(count_primary,trueGsfPFRecTracks));
00387 trueGsfPFRecTracks.clear();
00388 count_primary++;
00389 }
00390 }
00391 }
00392
00393
00394 const edm::OrphanHandle<GsfPFRecTrackCollection> gsfPfRefProd =
00395 iEvent.put(gsfPFRecTrackCollectionSecondary,"Secondary");
00396
00397
00398
00399 createGsfPFRecTrackRef(gsfPfRefProd,primaryGsfPFRecTracks,GsfPFMap);
00400
00401 for(unsigned int iGSF = 0; iGSF<primaryGsfPFRecTracks.size();iGSF++){
00402 gsfPFRecTrackCollection->push_back(primaryGsfPFRecTracks[iGSF]);
00403 }
00404 iEvent.put(gsfPFRecTrackCollection);
00405
00406 selGsfPFRecTracks.clear();
00407 GsfPFMap.clear();
00408 primaryGsfPFRecTracks.clear();
00409 }
00410
00411
00412 void
00413 PFElecTkProducer::createGsfPFRecTrackRef(const edm::OrphanHandle<reco::GsfPFRecTrackCollection>& gsfPfHandle,
00414 std::vector<reco::GsfPFRecTrack>& gsfPFRecTrackPrimary,
00415 const std::map<unsigned int, std::vector<reco::GsfPFRecTrack> >& MapPrimSec) {
00416 unsigned int cgsf=0;
00417 unsigned int csecgsf=0;
00418 for (std::map<unsigned int, std::vector<reco::GsfPFRecTrack> >::const_iterator igsf = MapPrimSec.begin();
00419 igsf != MapPrimSec.end(); igsf++,cgsf++) {
00420 vector<reco::GsfPFRecTrack> SecGsfPF = igsf->second;
00421 for (unsigned int iSecGsf=0; iSecGsf < SecGsfPF.size(); iSecGsf++) {
00422 edm::Ref<reco::GsfPFRecTrackCollection> refgprt(gsfPfHandle,csecgsf);
00423 gsfPFRecTrackPrimary[cgsf].addConvBremGsfPFRecTrackRef(refgprt);
00424 ++csecgsf;
00425 }
00426 }
00427
00428 return;
00429 }
00430
00431 int
00432 PFElecTkProducer::FindPfRef(const reco::PFRecTrackCollection & PfRTkColl,
00433 reco::GsfTrack gsftk,
00434 bool otherColl){
00435
00436
00437 if (&(*gsftk.seedRef())==0) return -1;
00438 ElectronSeedRef ElSeedRef=gsftk.extra()->seedRef().castTo<ElectronSeedRef>();
00439
00440 if (ElSeedRef->ctfTrack().isNull()){
00441 reco::PFRecTrackCollection::const_iterator pft=PfRTkColl.begin();
00442 reco::PFRecTrackCollection::const_iterator pftend=PfRTkColl.end();
00443 unsigned int i_pf=0;
00444 int ibest=-1;
00445 unsigned int ish_max=0;
00446 float dr_min=1000;
00447
00448
00449 for(;pft!=pftend;++pft){
00450 unsigned int ish=0;
00451
00452 float dph= fabs(pft->trackRef()->phi()-gsftk.phi());
00453 if (dph>TMath::Pi()) dph-= TMath::TwoPi();
00454 float det=fabs(pft->trackRef()->eta()-gsftk.eta());
00455 float dr =sqrt(dph*dph+det*det);
00456
00457 trackingRecHit_iterator hhit=
00458 pft->trackRef()->recHitsBegin();
00459 trackingRecHit_iterator hhit_end=
00460 pft->trackRef()->recHitsEnd();
00461
00462
00463
00464 for(;hhit!=hhit_end;++hhit){
00465 if (!(*hhit)->isValid()) continue;
00466 TrajectorySeed::const_iterator hit=
00467 gsftk.seedRef()->recHits().first;
00468 TrajectorySeed::const_iterator hit_end=
00469 gsftk.seedRef()->recHits().second;
00470 for(;hit!=hit_end;++hit){
00471 if (!(hit->isValid())) continue;
00472 if((*hhit)->sharesInput(&*(hit),TrackingRecHit::all)) ish++;
00473
00474
00475 }
00476
00477 }
00478
00479
00480 if ((ish>ish_max)||
00481 ((ish==ish_max)&&(dr<dr_min))){
00482 ish_max=ish;
00483 dr_min=dr;
00484 ibest=i_pf;
00485 }
00486
00487
00488
00489 i_pf++;
00490 }
00491 if (ibest<0) return -1;
00492
00493 if((ish_max==0) || (dr_min>0.05))return -1;
00494 if(otherColl && (ish_max==0)) return -1;
00495 return ibest;
00496 }
00497 else{
00498
00499
00500 reco::PFRecTrackCollection::const_iterator pft=PfRTkColl.begin();
00501 reco::PFRecTrackCollection::const_iterator pftend=PfRTkColl.end();
00502 unsigned int i_pf=0;
00503
00504 for(;pft!=pftend;++pft){
00505
00506 if (pft->trackRef()==ElSeedRef->ctfTrack()){
00507 return i_pf;
00508 }
00509 i_pf++;
00510 }
00511 }
00512 return -1;
00513 }
00514 bool PFElecTkProducer::isFifthStep(reco::PFRecTrackRef pfKfTrack) {
00515
00516 bool isFithStep = false;
00517
00518
00519 TrackRef kfref = pfKfTrack->trackRef();
00520 unsigned int Algo = 0;
00521 switch (kfref->algo()) {
00522 case TrackBase::undefAlgorithm:
00523 case TrackBase::ctf:
00524 case TrackBase::iter0:
00525 case TrackBase::iter1:
00526 Algo = 0;
00527 break;
00528 case TrackBase::iter2:
00529 Algo = 1;
00530 break;
00531 case TrackBase::iter3:
00532 Algo = 2;
00533 break;
00534 case TrackBase::iter4:
00535 Algo = 3;
00536 break;
00537 case TrackBase::iter5:
00538 Algo = 4;
00539 break;
00540 default:
00541 Algo = 5;
00542 break;
00543 }
00544 if ( Algo >= 4 ) {
00545 isFithStep = true;
00546 }
00547
00548 return isFithStep;
00549 }
00550
00551 bool
00552 PFElecTkProducer::applySelection(reco::GsfTrack gsftk) {
00553 if (&(*gsftk.seedRef())==0) return false;
00554 ElectronSeedRef ElSeedRef=gsftk.extra()->seedRef().castTo<ElectronSeedRef>();
00555
00556 bool passCut = false;
00557 if (ElSeedRef->ctfTrack().isNull()){
00558 if(ElSeedRef->caloCluster().isNull()) return passCut;
00559 SuperClusterRef scRef = ElSeedRef->caloCluster().castTo<SuperClusterRef>();
00560
00561 if(scRef.isNonnull()) {
00562 float caloEne = scRef->energy();
00563 float feta = fabs(scRef->eta()-gsftk.etaMode());
00564 float fphi = fabs(scRef->phi()-gsftk.phiMode());
00565 if (fphi>TMath::Pi()) fphi-= TMath::TwoPi();
00566 if(caloEne > SCEne_ && feta < detaGsfSC_ && fabs(fphi) < dphiGsfSC_)
00567 passCut = true;
00568 }
00569 }
00570 else {
00571
00572 passCut = true;
00573 }
00574 return passCut;
00575 }
00576 bool
00577 PFElecTkProducer::resolveGsfTracks(const vector<reco::GsfPFRecTrack> & GsfPFVec,
00578 unsigned int ngsf,
00579 vector<unsigned int> &secondaries,
00580 const reco::PFClusterCollection & theEClus) {
00581 bool debugCleaning = debugGsfCleaning_;
00582 bool n_keepGsf = true;
00583
00584 reco::GsfTrackRef nGsfTrack = GsfPFVec[ngsf].gsfTrackRef();
00585
00586 if (&(*nGsfTrack->seedRef())==0) return false;
00587 ElectronSeedRef nElSeedRef=nGsfTrack->extra()->seedRef().castTo<ElectronSeedRef>();
00588
00589
00590 TrajectoryStateOnSurface inTSOS = mtsTransform_.innerStateOnSurface((*nGsfTrack));
00591 GlobalVector ninnMom;
00592 float nPin = nGsfTrack->pMode();
00593 if(inTSOS.isValid()){
00594 mtsMode_->momentumFromModeCartesian(inTSOS,ninnMom);
00595 nPin = ninnMom.mag();
00596 }
00597
00598 float neta = nGsfTrack->innerMomentum().eta();
00599 float nphi = nGsfTrack->innerMomentum().phi();
00600
00601
00602
00603
00604 if(debugCleaning)
00605 cout << " PFElecTkProducer:: considering track " << nGsfTrack->pt()
00606 << " eta,phi " << nGsfTrack->eta() << ", " << nGsfTrack->phi() << endl;
00607
00608
00609 for (unsigned int igsf=0; igsf< GsfPFVec.size();igsf++) {
00610 if(igsf != ngsf ) {
00611 reco::GsfTrackRef iGsfTrack = GsfPFVec[igsf].gsfTrackRef();
00612
00613 if(debugCleaning)
00614 cout << " PFElecTkProducer:: and comparing with track " << iGsfTrack->pt()
00615 << " eta,phi " << iGsfTrack->eta() << ", " << iGsfTrack->phi() << endl;
00616
00617 float ieta = iGsfTrack->innerMomentum().eta();
00618 float iphi = iGsfTrack->innerMomentum().phi();
00619 float feta = fabs(neta - ieta);
00620 float fphi = fabs(nphi - iphi);
00621 if (fphi>TMath::Pi()) fphi-= TMath::TwoPi();
00622
00623
00624
00625 if(feta < 0.5 && fabs(fphi) < 1.0) {
00626 if(debugCleaning)
00627 cout << " Entering angular superloose preselection " << endl;
00628
00629 TrajectoryStateOnSurface i_inTSOS = mtsTransform_.innerStateOnSurface((*iGsfTrack));
00630 GlobalVector i_innMom;
00631 float iPin = iGsfTrack->pMode();
00632 if(i_inTSOS.isValid()){
00633 mtsMode_->momentumFromModeCartesian(i_inTSOS,i_innMom);
00634 iPin = i_innMom.mag();
00635 }
00636
00637 if (&(*iGsfTrack->seedRef())==0) continue;
00638 ElectronSeedRef iElSeedRef=iGsfTrack->extra()->seedRef().castTo<ElectronSeedRef>();
00639
00640 float SCEnergy = -1.;
00641
00642 bool areBothGsfEcalDriven = false;;
00643 bool isSameSC = isSameEgSC(nElSeedRef,iElSeedRef,areBothGsfEcalDriven,SCEnergy);
00644
00645
00646 if(areBothGsfEcalDriven ) {
00647 if(isSameSC) {
00648 float nEP = SCEnergy/nPin;
00649 float iEP = SCEnergy/iPin;
00650 if(debugCleaning)
00651 cout << " Entering SAME supercluster case "
00652 << " nEP " << nEP
00653 << " iEP " << iEP << endl;
00654
00655
00656
00657
00658
00659
00660
00661 bool isSameLayer = false;
00662 bool iGsfInnermostWithLostHits =
00663 isInnerMostWithLostHits(nGsfTrack,iGsfTrack,isSameLayer);
00664
00665
00666 if(debugCleaning)
00667 cout << " iGsf is InnerMostWithLostHits " << iGsfInnermostWithLostHits
00668 << " isSameLayer " << isSameLayer << endl;
00669
00670 if (iGsfInnermostWithLostHits) {
00671 n_keepGsf = false;
00672 return n_keepGsf;
00673 }
00674 else if(isSameLayer){
00675 if(fabs(iEP-1) < fabs(nEP-1)) {
00676 n_keepGsf = false;
00677 return n_keepGsf;
00678 }
00679 else{
00680 secondaries.push_back(igsf);
00681 }
00682 }
00683 else {
00684
00685 secondaries.push_back(igsf);
00686 }
00687 }
00688 }
00689 else {
00690
00691 float minBremDphi = minTangDist(GsfPFVec[ngsf],GsfPFVec[igsf]);
00692 float nETot = 0.;
00693 float iETot = 0.;
00694 bool isBothGsfTrackerDriven = false;
00695 bool nEcalDriven = false;
00696 bool iEcalDriven = false;
00697 bool isSameScEgPf = isSharingEcalEnergyWithEgSC(GsfPFVec[ngsf],
00698 GsfPFVec[igsf],
00699 nElSeedRef,
00700 iElSeedRef,
00701 theEClus,
00702 isBothGsfTrackerDriven,
00703 nEcalDriven,
00704 iEcalDriven,
00705 nETot,
00706 iETot);
00707
00708
00709 bool isSameLayer = false;
00710 bool iGsfInnermostWithLostHits =
00711 isInnerMostWithLostHits(nGsfTrack,iGsfTrack,isSameLayer);
00712
00713 if(isSameScEgPf) {
00714
00715
00716 if(debugCleaning) {
00717 cout << " Sharing ECAL energy passed "
00718 << " nEtot " << nETot
00719 << " iEtot " << iETot << endl;
00720 if(isBothGsfTrackerDriven)
00721 cout << " Both Track are trackerDriven " << endl;
00722 }
00723
00724
00725 if (iGsfInnermostWithLostHits) {
00726 n_keepGsf = false;
00727 return n_keepGsf;
00728 }
00729 else if(isSameLayer){
00730
00731
00732
00733
00734
00735
00736
00737
00738 if(isBothGsfTrackerDriven == false) {
00739
00740 if(iEcalDriven) {
00741 n_keepGsf = false;
00742 return n_keepGsf;
00743 }
00744 else {
00745 secondaries.push_back(igsf);
00746 }
00747 }
00748 else {
00749
00750
00751
00752 float ETot = -1;
00753 if(nETot != iETot) {
00754 if(nETot > iETot)
00755 ETot = nETot;
00756 else
00757 ETot = iETot;
00758 }
00759 else {
00760 ETot = nETot;
00761 }
00762 float nEP = ETot/nPin;
00763 float iEP = ETot/iPin;
00764
00765
00766 if(debugCleaning)
00767 cout << " nETot " << nETot
00768 << " iETot " << iETot
00769 << " ETot " << ETot << endl
00770 << " nPin " << nPin
00771 << " iPin " << iPin
00772 << " nEP " << nEP
00773 << " iEP " << iEP << endl;
00774
00775
00776 if(fabs(iEP-1) < fabs(nEP-1)) {
00777 n_keepGsf = false;
00778 return n_keepGsf;
00779 }
00780 else{
00781 secondaries.push_back(igsf);
00782 }
00783 }
00784 }
00785 else {
00786 secondaries.push_back(igsf);
00787 }
00788 }
00789 else if(feta < detaCutGsfClean_ && minBremDphi < dphiCutGsfClean_) {
00790
00791 bool secPushedBack = false;
00792 if(nEcalDriven == false && nETot == 0.) {
00793 n_keepGsf = false;
00794 return n_keepGsf;
00795 }
00796 else if(iEcalDriven == false && iETot == 0.) {
00797 secondaries.push_back(igsf);
00798 secPushedBack = true;
00799 }
00800 if(debugCleaning)
00801 cout << " Close Tracks "
00802 << " feta " << feta << " fabs(fphi) " << fabs(fphi)
00803 << " minBremDphi " << minBremDphi
00804 << " nETot " << nETot
00805 << " iETot " << iETot
00806 << " nLostHits " << nGsfTrack->trackerExpectedHitsInner().numberOfLostHits()
00807 << " iLostHits " << iGsfTrack->trackerExpectedHitsInner().numberOfLostHits() << endl;
00808
00809
00810 if(applyAngularGsfClean_) {
00811 if (iGsfInnermostWithLostHits) {
00812 n_keepGsf = false;
00813 return n_keepGsf;
00814 }
00815 else if(isSameLayer == false) {
00816 if(secPushedBack == false)
00817 secondaries.push_back(igsf);
00818 }
00819 }
00820 }
00821 else if(feta < 0.1 && minBremDphi < 0.2){
00822
00823
00824 if(debugCleaning)
00825 cout << " Close Tracks and failed all the conditions "
00826 << " feta " << feta << " fabs(fphi) " << fabs(fphi)
00827 << " minBremDphi " << minBremDphi
00828 << " nETot " << nETot
00829 << " iETot " << iETot
00830 << " nLostHits " << nGsfTrack->trackerExpectedHitsInner().numberOfLostHits()
00831 << " iLostHits " << iGsfTrack->trackerExpectedHitsInner().numberOfLostHits() << endl;
00832
00833 if(nEcalDriven == false && nETot == 0.) {
00834 n_keepGsf = false;
00835 return n_keepGsf;
00836 }
00837
00838 }
00839 }
00840 }
00841 }
00842 }
00843
00844 return n_keepGsf;
00845 }
00846 float
00847 PFElecTkProducer::minTangDist(const reco::GsfPFRecTrack primGsf,
00848 const reco::GsfPFRecTrack secGsf) {
00849
00850 float minDeta = 1000.;
00851 float minDphi = 1000.;
00852
00853
00854 std::vector<reco::PFBrem> primPFBrem = primGsf.PFRecBrem();
00855 std::vector<reco::PFBrem> secPFBrem = secGsf.PFRecBrem();
00856
00857
00858 unsigned int cbrem = 0;
00859 for (unsigned isbrem = 0; isbrem < secPFBrem.size(); isbrem++) {
00860 if(secPFBrem[isbrem].indTrajPoint() == 99) continue;
00861 const reco::PFTrajectoryPoint& atSecECAL
00862 = secPFBrem[isbrem].extrapolatedPoint( reco::PFTrajectoryPoint::ECALEntrance );
00863 if( ! atSecECAL.isValid() ) continue;
00864 float secEta = atSecECAL.positionREP().Eta();
00865 float secPhi = atSecECAL.positionREP().Phi();
00866
00867 unsigned int sbrem = 0;
00868 for(unsigned ipbrem = 0; ipbrem < primPFBrem.size(); ipbrem++) {
00869 if(primPFBrem[ipbrem].indTrajPoint() == 99) continue;
00870 const reco::PFTrajectoryPoint& atPrimECAL
00871 = primPFBrem[ipbrem].extrapolatedPoint( reco::PFTrajectoryPoint::ECALEntrance );
00872 if( ! atPrimECAL.isValid() ) continue;
00873 sbrem++;
00874 if(sbrem <= 3) {
00875 float primEta = atPrimECAL.positionREP().Eta();
00876 float primPhi = atPrimECAL.positionREP().Phi();
00877
00878 float deta = fabs(primEta - secEta);
00879 float dphi = fabs(primPhi - secPhi);
00880 if (dphi>TMath::Pi()) dphi-= TMath::TwoPi();
00881 if(fabs(dphi) < minDphi) {
00882 minDeta = deta;
00883 minDphi = fabs(dphi);
00884 }
00885 }
00886 }
00887
00888
00889 cbrem++;
00890 if(cbrem == 3)
00891 break;
00892 }
00893 return minDphi;
00894 }
00895 bool
00896 PFElecTkProducer::isSameEgSC(const reco::ElectronSeedRef& nSeedRef,
00897 const reco::ElectronSeedRef& iSeedRef,
00898 bool& bothGsfEcalDriven,
00899 float& SCEnergy) {
00900
00901 bool isSameSC = false;
00902
00903 if(nSeedRef->caloCluster().isNonnull() && iSeedRef->caloCluster().isNonnull()) {
00904 SuperClusterRef nscRef = nSeedRef->caloCluster().castTo<SuperClusterRef>();
00905 SuperClusterRef iscRef = iSeedRef->caloCluster().castTo<SuperClusterRef>();
00906
00907 if(nscRef.isNonnull() && iscRef.isNonnull()) {
00908 bothGsfEcalDriven = true;
00909 if(nscRef == iscRef) {
00910 isSameSC = true;
00911
00912 SCEnergy = nscRef->energy();
00913 }
00914 }
00915 }
00916 return isSameSC;
00917 }
00918 bool
00919 PFElecTkProducer::isSharingEcalEnergyWithEgSC(const reco::GsfPFRecTrack& nGsfPFRecTrack,
00920 const reco::GsfPFRecTrack& iGsfPFRecTrack,
00921 const reco::ElectronSeedRef& nSeedRef,
00922 const reco::ElectronSeedRef& iSeedRef,
00923 const reco::PFClusterCollection& theEClus,
00924 bool& bothGsfTrackerDriven,
00925 bool& nEcalDriven,
00926 bool& iEcalDriven,
00927 float& nEnergy,
00928 float& iEnergy) {
00929
00930 bool isSharingEnergy = false;
00931
00932
00933 bool oneEcalDriven = true;
00934 SuperClusterRef scRef;
00935 GsfPFRecTrack gsfPfTrack;
00936
00937 if(nSeedRef->caloCluster().isNonnull()) {
00938 scRef = nSeedRef->caloCluster().castTo<SuperClusterRef>();
00939 nEnergy = scRef->energy();
00940 nEcalDriven = true;
00941 gsfPfTrack = iGsfPFRecTrack;
00942 }
00943 else if(iSeedRef->caloCluster().isNonnull()){
00944 scRef = iSeedRef->caloCluster().castTo<SuperClusterRef>();
00945 iEnergy = scRef->energy();
00946 iEcalDriven = true;
00947 gsfPfTrack = nGsfPFRecTrack;
00948 }
00949 else{
00950 oneEcalDriven = false;
00951 }
00952
00953 if(oneEcalDriven) {
00954
00955
00956 vector<PFCluster> vecPFClusters;
00957 vecPFClusters.clear();
00958
00959 for (PFClusterCollection::const_iterator clus = theEClus.begin();
00960 clus != theEClus.end();
00961 clus++ ) {
00962 PFCluster clust = *clus;
00963 clust.calculatePositionREP();
00964
00965 float deta = fabs(scRef->position().eta() - clust.position().eta());
00966 float dphi = fabs(scRef->position().phi() - clust.position().phi());
00967 if (dphi>TMath::Pi()) dphi-= TMath::TwoPi();
00968
00969
00970
00971 if(deta < 0.5 && fabs(dphi) < 1.0) {
00972 bool foundLink = false;
00973 double distGsf = LinkByRecHit::testTrackAndClusterByRecHit(gsfPfTrack , clust );
00974
00975 if(distGsf > 0.) {
00976 if(nEcalDriven)
00977 iEnergy += clust.energy();
00978 else
00979 nEnergy += clust.energy();
00980 vecPFClusters.push_back(clust);
00981 foundLink = true;
00982 }
00983
00984 if(foundLink == false) {
00985 vector<PFBrem> primPFBrem = gsfPfTrack.PFRecBrem();
00986 for(unsigned ipbrem = 0; ipbrem < primPFBrem.size(); ipbrem++) {
00987 if(primPFBrem[ipbrem].indTrajPoint() == 99) continue;
00988 reco::PFRecTrack pfBremTrack = primPFBrem[ipbrem];
00989 double dist = LinkByRecHit::testTrackAndClusterByRecHit(pfBremTrack , clust, true );
00990 if(dist > 0.) {
00991 if(nEcalDriven)
00992 iEnergy += clust.energy();
00993 else
00994 nEnergy += clust.energy();
00995 vecPFClusters.push_back(clust);
00996 foundLink = true;
00997 }
00998 }
00999 }
01000 }
01001 }
01002 if(vecPFClusters.size() > 0 ) {
01003 for(unsigned int pf = 0; pf < vecPFClusters.size(); pf++) {
01004 bool isCommon = ClusterClusterMapping::overlap(vecPFClusters[pf],*scRef);
01005 if(isCommon) {
01006 isSharingEnergy = true;
01007 }
01008 break;
01009 }
01010 }
01011 }
01012 else {
01013
01014
01015 bothGsfTrackerDriven = true;
01016 vector<PFCluster> nPFCluster;
01017 vector<PFCluster> iPFCluster;
01018
01019 nPFCluster.clear();
01020 iPFCluster.clear();
01021
01022 for (PFClusterCollection::const_iterator clus = theEClus.begin();
01023 clus != theEClus.end();
01024 clus++ ) {
01025 PFCluster clust = *clus;
01026 clust.calculatePositionREP();
01027
01028 float ndeta = fabs(nGsfPFRecTrack.gsfTrackRef()->eta() - clust.position().eta());
01029 float ndphi = fabs(nGsfPFRecTrack.gsfTrackRef()->phi() - clust.position().phi());
01030 if (ndphi>TMath::Pi()) ndphi-= TMath::TwoPi();
01031
01032
01033 if(ndeta < 0.5 && fabs(ndphi) < 1.0) {
01034 bool foundNLink = false;
01035
01036 double distGsf = LinkByRecHit::testTrackAndClusterByRecHit(nGsfPFRecTrack , clust );
01037 if(distGsf > 0.) {
01038 nPFCluster.push_back(clust);
01039 nEnergy += clust.energy();
01040 foundNLink = true;
01041 }
01042 if(foundNLink == false) {
01043 vector<PFBrem> primPFBrem = nGsfPFRecTrack.PFRecBrem();
01044 for(unsigned ipbrem = 0; ipbrem < primPFBrem.size(); ipbrem++) {
01045 if(primPFBrem[ipbrem].indTrajPoint() == 99) continue;
01046 reco::PFRecTrack pfBremTrack = primPFBrem[ipbrem];
01047 if(foundNLink == false) {
01048 double dist = LinkByRecHit::testTrackAndClusterByRecHit(pfBremTrack , clust, true );
01049 if(dist > 0.) {
01050 nPFCluster.push_back(clust);
01051 nEnergy += clust.energy();
01052 foundNLink = true;
01053 }
01054 }
01055 }
01056 }
01057 }
01058
01059 float ideta = fabs(iGsfPFRecTrack.gsfTrackRef()->eta() - clust.position().eta());
01060 float idphi = fabs(iGsfPFRecTrack.gsfTrackRef()->phi() - clust.position().phi());
01061 if (idphi>TMath::Pi()) idphi-= TMath::TwoPi();
01062
01063
01064 if(ideta < 0.5 && fabs(idphi) < 1.0) {
01065 bool foundILink = false;
01066 double dist = LinkByRecHit::testTrackAndClusterByRecHit(iGsfPFRecTrack , clust );
01067 if(dist > 0.) {
01068 iPFCluster.push_back(clust);
01069 iEnergy += clust.energy();
01070 foundILink = true;
01071 }
01072 if(foundILink == false) {
01073 vector<PFBrem> primPFBrem = iGsfPFRecTrack.PFRecBrem();
01074 for(unsigned ipbrem = 0; ipbrem < primPFBrem.size(); ipbrem++) {
01075 if(primPFBrem[ipbrem].indTrajPoint() == 99) continue;
01076 reco::PFRecTrack pfBremTrack = primPFBrem[ipbrem];
01077 if(foundILink == false) {
01078 double dist = LinkByRecHit::testTrackAndClusterByRecHit(pfBremTrack , clust, true );
01079 if(dist > 0.) {
01080 iPFCluster.push_back(clust);
01081 iEnergy += clust.energy();
01082 foundILink = true;
01083 }
01084 }
01085 }
01086 }
01087 }
01088 }
01089
01090
01091 if(nPFCluster.size() > 0 && iPFCluster.size() > 0) {
01092 for(unsigned int npf = 0; npf < nPFCluster.size(); npf++) {
01093 for(unsigned int ipf = 0; ipf < iPFCluster.size(); ipf++) {
01094 bool isCommon = ClusterClusterMapping::overlap(nPFCluster[npf],iPFCluster[ipf]);
01095 if(isCommon) {
01096 isSharingEnergy = true;
01097 break;
01098 }
01099 }
01100 if(isSharingEnergy)
01101 break;
01102 }
01103 }
01104 }
01105
01106 return isSharingEnergy;
01107 }
01108 bool PFElecTkProducer::isInnerMost(const reco::GsfTrackRef& nGsfTrack,
01109 const reco::GsfTrackRef& iGsfTrack,
01110 bool& sameLayer) {
01111
01112
01113
01114
01115 reco::HitPattern gsfHitPattern1 = nGsfTrack->hitPattern();
01116 reco::HitPattern gsfHitPattern2 = iGsfTrack->hitPattern();
01117
01118
01119 int gsfHitCounter1 = 0 ;
01120 trackingRecHit_iterator elHitsIt1 ;
01121 for
01122 ( elHitsIt1 = nGsfTrack->recHitsBegin() ;
01123 elHitsIt1 != nGsfTrack->recHitsEnd() ;
01124 elHitsIt1++, gsfHitCounter1++ )
01125 { if (((**elHitsIt1).isValid())) break ; }
01126
01127 int gsfHitCounter2 = 0 ;
01128 trackingRecHit_iterator elHitsIt2 ;
01129 for
01130 ( elHitsIt2 = iGsfTrack->recHitsBegin() ;
01131 elHitsIt2 != iGsfTrack->recHitsEnd() ;
01132 elHitsIt2++, gsfHitCounter2++ )
01133 { if (((**elHitsIt2).isValid())) break ; }
01134
01135 uint32_t gsfHit1 = gsfHitPattern1.getHitPattern(gsfHitCounter1) ;
01136 uint32_t gsfHit2 = gsfHitPattern2.getHitPattern(gsfHitCounter2) ;
01137
01138
01139 if (gsfHitPattern1.getSubStructure(gsfHit1)!=gsfHitPattern2.getSubStructure(gsfHit2))
01140 {
01141 return (gsfHitPattern2.getSubStructure(gsfHit2)<gsfHitPattern1.getSubStructure(gsfHit1));
01142 }
01143 else if (gsfHitPattern1.getLayer(gsfHit1)!=gsfHitPattern2.getLayer(gsfHit2))
01144 {
01145 return (gsfHitPattern2.getLayer(gsfHit2)<gsfHitPattern1.getLayer(gsfHit1));
01146 }
01147 else
01148 {
01149 sameLayer = true;
01150 return false;
01151 }
01152 }
01153 bool PFElecTkProducer::isInnerMostWithLostHits(const reco::GsfTrackRef& nGsfTrack,
01154 const reco::GsfTrackRef& iGsfTrack,
01155 bool& sameLayer) {
01156
01157
01158 unsigned int nLostHits = nGsfTrack->trackerExpectedHitsInner().numberOfLostHits();
01159 unsigned int iLostHits = iGsfTrack->trackerExpectedHitsInner().numberOfLostHits();
01160
01161 if (nLostHits!=iLostHits) {
01162 return (nLostHits > iLostHits);
01163 }
01164 else {
01165 sameLayer = true;
01166 return false;
01167 }
01168 }
01169
01170
01171
01172
01173 void
01174 PFElecTkProducer::beginRun(edm::Run& run,
01175 const EventSetup& iSetup)
01176 {
01177 ESHandle<MagneticField> magneticField;
01178 iSetup.get<IdealMagneticFieldRecord>().get(magneticField);
01179
01180 ESHandle<TrackerGeometry> tracker;
01181 iSetup.get<TrackerDigiGeometryRecord>().get(tracker);
01182
01183
01184 mtsTransform_ = MultiTrajectoryStateTransform(tracker.product(),magneticField.product());
01185
01186
01187 pfTransformer_= new PFTrackTransformer(math::XYZVector(magneticField->inTesla(GlobalPoint(0,0,0))));
01188
01189
01190 edm::ESHandle<TransientTrackBuilder> builder;
01191 iSetup.get<TransientTrackRecord>().get("TransientTrackBuilder", builder);
01192 TransientTrackBuilder thebuilder = *(builder.product());
01193
01194
01195 if(useConvBremFinder_) {
01196 FILE * fileConvBremID = fopen(path_mvaWeightFileConvBrem_.c_str(), "r");
01197 if (fileConvBremID) {
01198 fclose(fileConvBremID);
01199 }
01200 else {
01201 string err = "PFElecTkProducer: cannot open weight file '";
01202 err += path_mvaWeightFileConvBrem_;
01203 err += "'";
01204 throw invalid_argument( err );
01205 }
01206 }
01207 convBremFinder_ = new ConvBremPFTrackFinder(thebuilder,mvaConvBremFinderID_,path_mvaWeightFileConvBrem_);
01208
01209 }
01210
01211
01212 void
01213 PFElecTkProducer::endRun() {
01214 delete pfTransformer_;
01215 }
01216
01217