00001 #include <iostream>
00002 #include <vector>
00003 #include <memory>
00004
00005
00006 #include "FWCore/Framework/interface/Event.h"
00007 #include "FWCore/Framework/interface/EventSetup.h"
00008 #include "DataFormats/Common/interface/Handle.h"
00009 #include "FWCore/Framework/interface/ESHandle.h"
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/Utilities/interface/Exception.h"
00012
00013 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00014 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00015 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
00016
00017
00018 #include "DataFormats/EgammaTrackReco/interface/TrackCaloClusterAssociation.h"
00019 #include "DataFormats/EgammaCandidates/interface/Photon.h"
00020 #include "DataFormats/EgammaCandidates/interface/PhotonFwd.h"
00021 #include "DataFormats/EgammaCandidates/interface/Conversion.h"
00022 #include "DataFormats/EgammaCandidates/interface/ConversionFwd.h"
00023 #include "DataFormats/TrackReco/interface/Track.h"
00024 #include "DataFormats/TrackReco/interface/TrackExtra.h"
00025 #include "DataFormats/VertexReco/interface/Vertex.h"
00026 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00027
00028 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00029
00030 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00031
00032 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionTrackEcalImpactPoint.h"
00033 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionTrackPairFinder.h"
00034 #include "RecoEgamma/EgammaPhotonAlgos/interface/ConversionVertexFinder.h"
00035 #include "RecoEgamma/EgammaPhotonProducers/interface/ConvertedPhotonProducer.h"
00036
00037 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
00038 #include "RecoTracker/Record/interface/TrackerRecoGeometryRecord.h"
00039
00040 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
00041 #include "TrackingTools/TransientTrack/interface/TrackTransientTrack.h"
00042 #include "TrackingTools/PatternTools/interface/TwoTrackMinimumDistance.h"
00043 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
00044 #include "RecoEgamma/EgammaIsolationAlgos/interface/EgammaTowerIsolation.h"
00045
00046 ConvertedPhotonProducer::ConvertedPhotonProducer(const edm::ParameterSet& config) :
00047
00048 conf_(config),
00049 theTrackPairFinder_(0),
00050 theVertexFinder_(0),
00051 theEcalImpactPositionFinder_(0),
00052 theLikelihoodCalc_(0)
00053
00054
00055 {
00056
00057
00058
00059
00060
00061
00062
00063
00064 bcBarrelCollection_ = conf_.getParameter<edm::InputTag>("bcBarrelCollection");
00065 bcEndcapCollection_ = conf_.getParameter<edm::InputTag>("bcEndcapCollection");
00066
00067 scHybridBarrelProducer_ = conf_.getParameter<edm::InputTag>("scHybridBarrelProducer");
00068 scIslandEndcapProducer_ = conf_.getParameter<edm::InputTag>("scIslandEndcapProducer");
00069
00070 conversionOITrackProducer_ = conf_.getParameter<std::string>("conversionOITrackProducer");
00071 conversionIOTrackProducer_ = conf_.getParameter<std::string>("conversionIOTrackProducer");
00072
00073 outInTrackSCAssociationCollection_ = conf_.getParameter<std::string>("outInTrackSCAssociation");
00074 inOutTrackSCAssociationCollection_ = conf_.getParameter<std::string>("inOutTrackSCAssociation");
00075
00076 algoName_ = conf_.getParameter<std::string>( "AlgorithmName" );
00077
00078 hcalTowers_ = conf_.getParameter<edm::InputTag>("hcalTowers");
00079 hOverEConeSize_ = conf_.getParameter<double>("hOverEConeSize");
00080 maxHOverE_ = conf_.getParameter<double>("maxHOverE");
00081 minSCEt_ = conf_.getParameter<double>("minSCEt");
00082 recoverOneTrackCase_ = conf_.getParameter<bool>( "recoverOneTrackCase" );
00083 dRForConversionRecovery_ = conf_.getParameter<double>("dRForConversionRecovery");
00084 deltaCotCut_ = conf_.getParameter<double>("deltaCotCut");
00085 minApproachDisCut_ = conf_.getParameter<double>("minApproachDisCut");
00086
00087 maxNumOfCandidates_ = conf_.getParameter<int>("maxNumOfCandidates");
00088 risolveAmbiguity_ = conf_.getParameter<bool>("risolveConversionAmbiguity");
00089 likelihoodWeights_= conf_.getParameter<std::string>("MVA_weights_location");
00090
00091
00092
00093 ConvertedPhotonCollection_ = conf_.getParameter<std::string>("convertedPhotonCollection");
00094 CleanedConvertedPhotonCollection_ = conf_.getParameter<std::string>("cleanedConvertedPhotonCollection");
00095
00096
00097
00098 produces< reco::ConversionCollection >(ConvertedPhotonCollection_);
00099 produces< reco::ConversionCollection >(CleanedConvertedPhotonCollection_);
00100
00101
00102 theTrackPairFinder_ = new ConversionTrackPairFinder ();
00103 edm::FileInPath path_mvaWeightFile(likelihoodWeights_.c_str() );
00104 theLikelihoodCalc_ = new ConversionLikelihoodCalculator();
00105 theLikelihoodCalc_->setWeightsFile(path_mvaWeightFile.fullPath().c_str());
00106
00107 theVertexFinder_ = new ConversionVertexFinder ( conf_);
00108
00109
00110
00111 nEvt_=0;
00112
00113
00114 theEcalImpactPositionFinder_ =0;
00115
00116 }
00117
00118 ConvertedPhotonProducer::~ConvertedPhotonProducer() {
00119 delete theTrackPairFinder_;
00120 delete theLikelihoodCalc_;
00121 delete theVertexFinder_;
00122 }
00123
00124
00125
00126 void ConvertedPhotonProducer::beginRun (edm::Run& r, edm::EventSetup const & theEventSetup) {
00127
00128
00129
00130
00131 theEventSetup.get<IdealMagneticFieldRecord>().get(theMF_);
00132
00133
00134 theEcalImpactPositionFinder_ = new ConversionTrackEcalImpactPoint ( &(*theMF_) );
00135
00136
00137 theEventSetup.get<TransientTrackRecord>().get("TransientTrackBuilder",theTransientTrackBuilder_);
00138
00139
00140 }
00141
00142
00143 void ConvertedPhotonProducer::endRun (edm::Run& r, edm::EventSetup const & theEventSetup) {
00144 delete theEcalImpactPositionFinder_;
00145 }
00146
00147
00148 void ConvertedPhotonProducer::endJob () {
00149
00150
00151
00152
00153
00154 }
00155
00156
00157
00158 void ConvertedPhotonProducer::produce(edm::Event& theEvent, const edm::EventSetup& theEventSetup) {
00159
00160 using namespace edm;
00161 nEvt_++;
00162
00163
00164
00165
00166
00167
00168
00169
00170 reco::ConversionCollection outputConvPhotonCollection;
00171 std::auto_ptr<reco::ConversionCollection> outputConvPhotonCollection_p(new reco::ConversionCollection);
00172
00173 reco::ConversionCollection cleanedConversionCollection;
00174 std::auto_ptr<reco::ConversionCollection> cleanedConversionCollection_p(new reco::ConversionCollection);
00175
00176
00177
00178 bool validBarrelSCHandle=true;
00179 edm::Handle<edm::View<reco::CaloCluster> > scBarrelHandle;
00180 theEvent.getByLabel(scHybridBarrelProducer_,scBarrelHandle);
00181 if (!scBarrelHandle.isValid()) {
00182 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<scHybridBarrelProducer_.label();
00183 validBarrelSCHandle=false;
00184 }
00185
00186
00187 bool validEndcapSCHandle=true;
00188 edm::Handle<edm::View<reco::CaloCluster> > scEndcapHandle;
00189 theEvent.getByLabel(scIslandEndcapProducer_,scEndcapHandle);
00190 if (!scEndcapHandle.isValid()) {
00191 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<scIslandEndcapProducer_.label();
00192 validEndcapSCHandle=false;
00193 }
00194
00195
00197 bool validTrackInputs=true;
00198 Handle<reco::TrackCollection> outInTrkHandle;
00199 theEvent.getByLabel(conversionOITrackProducer_, outInTrkHandle);
00200 if (!outInTrkHandle.isValid()) {
00201
00202 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the conversionOITrack " << "\n";
00203 validTrackInputs=false;
00204 }
00205
00206
00207
00209 Handle<reco::TrackCaloClusterPtrAssociation> outInTrkSCAssocHandle;
00210 theEvent.getByLabel( conversionOITrackProducer_, outInTrackSCAssociationCollection_, outInTrkSCAssocHandle);
00211 if (!outInTrkSCAssocHandle.isValid()) {
00212
00213 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product " << outInTrackSCAssociationCollection_.c_str() <<"\n";
00214 validTrackInputs=false;
00215 }
00216
00218 Handle<reco::TrackCollection> inOutTrkHandle;
00219 theEvent.getByLabel(conversionIOTrackProducer_, inOutTrkHandle);
00220 if (!inOutTrkHandle.isValid()) {
00221
00222 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the conversionIOTrack " << "\n";
00223 validTrackInputs=false;
00224 }
00225
00226
00227
00229
00230 Handle<reco::TrackCollection> generalTrkHandle;
00231 if ( recoverOneTrackCase_ ) {
00232 theEvent.getByLabel("generalTracks", generalTrkHandle);
00233 if (!generalTrkHandle.isValid()) {
00234
00235 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the genralTracks " << "\n";
00236 }
00237 }
00238
00240 Handle<reco::TrackCaloClusterPtrAssociation> inOutTrkSCAssocHandle;
00241 theEvent.getByLabel( conversionIOTrackProducer_, inOutTrackSCAssociationCollection_, inOutTrkSCAssocHandle);
00242 if (!inOutTrkSCAssocHandle.isValid()) {
00243
00244 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product " << inOutTrackSCAssociationCollection_.c_str() <<"\n";
00245 validTrackInputs=false;
00246 }
00247
00248
00249
00250
00251
00252 bool validBarrelBCHandle=true;
00253 edm::Handle<edm::View<reco::CaloCluster> > bcBarrelHandle;
00254 theEvent.getByLabel( bcBarrelCollection_, bcBarrelHandle);
00255 if (!bcBarrelHandle.isValid()) {
00256 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<bcBarrelCollection_.label();
00257 validBarrelBCHandle=false;
00258 }
00259
00260
00261
00262 bool validEndcapBCHandle=true;
00263 edm::Handle<edm::View<reco::CaloCluster> > bcEndcapHandle;
00264 theEvent.getByLabel( bcEndcapCollection_, bcEndcapHandle);
00265 if (!bcEndcapHandle.isValid()) {
00266 edm::LogError("ConvertedPhotonProducer") << "Error! Can't get the product "<<bcEndcapCollection_.label();
00267 validEndcapBCHandle=true;
00268 }
00269
00270
00271
00272 Handle<CaloTowerCollection> hcalTowersHandle;
00273 theEvent.getByLabel(hcalTowers_, hcalTowersHandle);
00274
00275
00276 theEventSetup.get<CaloGeometryRecord>().get(theCaloGeom_);
00277
00278
00279 if ( validTrackInputs ) {
00280
00281 std::vector<reco::TransientTrack> t_outInTrk = ( *theTransientTrackBuilder_ ).build(outInTrkHandle );
00282 std::vector<reco::TransientTrack> t_inOutTrk = ( *theTransientTrackBuilder_ ).build(inOutTrkHandle );
00283
00284
00286
00287 std::map<std::vector<reco::TransientTrack>, reco::CaloClusterPtr> allPairs;
00288 allPairs = theTrackPairFinder_->run(t_outInTrk, outInTrkHandle, outInTrkSCAssocHandle, t_inOutTrk, inOutTrkHandle, inOutTrkSCAssocHandle );
00289
00290
00291 buildCollections(theEventSetup, scBarrelHandle, bcBarrelHandle, hcalTowersHandle, generalTrkHandle, allPairs, outputConvPhotonCollection);
00292 buildCollections(theEventSetup, scEndcapHandle, bcEndcapHandle, hcalTowersHandle, generalTrkHandle, allPairs, outputConvPhotonCollection);
00293 }
00294
00295
00296 outputConvPhotonCollection_p->assign(outputConvPhotonCollection.begin(),outputConvPhotonCollection.end());
00297
00298 const edm::OrphanHandle<reco::ConversionCollection> conversionHandle= theEvent.put( outputConvPhotonCollection_p, ConvertedPhotonCollection_);
00299
00300
00301
00302 if ( validBarrelSCHandle) cleanCollections(scBarrelHandle,
00303 conversionHandle,
00304 cleanedConversionCollection);
00305 if ( validEndcapSCHandle) cleanCollections(scEndcapHandle,
00306 conversionHandle,
00307 cleanedConversionCollection);
00308
00309
00310 cleanedConversionCollection_p->assign(cleanedConversionCollection.begin(),cleanedConversionCollection.end());
00311 theEvent.put( cleanedConversionCollection_p, CleanedConvertedPhotonCollection_);
00312
00313
00314 }
00315
00316
00317 void ConvertedPhotonProducer::buildCollections ( edm::EventSetup const & es,
00318 const edm::Handle<edm::View<reco::CaloCluster> > & scHandle,
00319 const edm::Handle<edm::View<reco::CaloCluster> > & bcHandle,
00320 const edm::Handle<CaloTowerCollection> & hcalTowersHandle,
00321 const edm::Handle<reco::TrackCollection> & generalTrkHandle,
00322 std::map<std::vector<reco::TransientTrack>, reco::CaloClusterPtr>& allPairs,
00323 reco::ConversionCollection & outputConvPhotonCollection)
00324
00325 {
00326
00327
00328 reco::Conversion::ConversionAlgorithm algo = reco::Conversion::algoByName(algoName_);
00329
00330 std::vector<reco::TransientTrack> t_generalTrk;
00331 if ( recoverOneTrackCase_ ) t_generalTrk = ( *theTransientTrackBuilder_ ).build(generalTrkHandle );
00332
00333
00334
00335 int myCands=0;
00336 reco::CaloClusterPtrVector scPtrVec;
00337 for (unsigned i = 0; i < scHandle->size(); ++i ) {
00338 reco::CaloClusterPtr aClus= scHandle->ptrAt(i);
00339
00340
00341 if (aClus->energy()/cosh(aClus->eta()) <= minSCEt_) continue;
00342 const reco::CaloCluster* pClus=&(*aClus);
00343 const reco::SuperCluster* sc=dynamic_cast<const reco::SuperCluster*>(pClus);
00344 const CaloTowerCollection* hcalTowersColl = hcalTowersHandle.product();
00345 EgammaTowerIsolation towerIso(hOverEConeSize_,0.,0.,-1,hcalTowersColl) ;
00346 double HoE=towerIso.getTowerESum(sc)/sc->energy();
00347 if (HoE>=maxHOverE_) continue;
00349
00350
00351 std::vector<edm::Ref<reco::TrackCollection> > trackPairRef;
00352 std::vector<math::XYZPoint> trackInnPos;
00353 std::vector<math::XYZVector> trackPin;
00354 std::vector<math::XYZVector> trackPout;
00355 float minAppDist=-99;
00356
00357
00358
00359
00361 const reco::Particle::Point vtx( 0, 0, 0 );
00362 reco::Vertex theConversionVertex;
00363
00364 math::XYZVector direction =aClus->position() - vtx;
00365 math::XYZVector momentum = direction.unit() * aClus->energy();
00366 const reco::Particle::LorentzVector p4(momentum.x(), momentum.y(), momentum.z(), aClus->energy() );
00367
00368 int nFound=0;
00369 if ( allPairs.size() ) {
00370
00371 nFound=0;
00372
00373
00374 for ( std::map<std::vector<reco::TransientTrack>, reco::CaloClusterPtr>::const_iterator iPair= allPairs.begin(); iPair!= allPairs.end(); ++iPair ) {
00375 scPtrVec.clear();
00376
00377 reco::CaloClusterPtr caloPtr=iPair->second;
00378 if ( !( aClus == caloPtr ) ) continue;
00379
00380 scPtrVec.push_back(aClus);
00381 nFound++;
00382
00383 std::vector<math::XYZPoint> trkPositionAtEcal = theEcalImpactPositionFinder_->find( iPair->first, bcHandle );
00384 std::vector<reco::CaloClusterPtr> matchingBC = theEcalImpactPositionFinder_->matchingBC();
00385
00386
00387 minAppDist=-99;
00388 const std::string metname = "ConvertedPhotons|ConvertedPhotonProducer";
00389 if ( (iPair->first).size() > 1 ) {
00390 try{
00391
00392 theVertexFinder_->run(iPair->first, theConversionVertex );
00393
00394
00395 }
00396 catch ( cms::Exception& e ) {
00397
00398 edm::LogWarning(metname) << "cms::Exception caught in ConvertedPhotonProducer::produce\n"
00399 << e.explainSelf();
00400
00401 }
00402
00403
00404
00405
00406
00407
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419
00420
00421
00423 trackPairRef.clear();
00424 trackInnPos.clear();
00425 trackPin.clear();
00426 trackPout.clear();
00427
00428
00429 for ( std::vector<reco::TransientTrack>::const_iterator iTk=(iPair->first).begin(); iTk!= (iPair->first).end(); ++iTk) {
00430
00431
00432 const reco::TrackTransientTrack* ttt = dynamic_cast<const reco::TrackTransientTrack*>(iTk->basicTransientTrack());
00433 reco::TrackRef myTkRef= ttt->persistentTrackRef();
00434
00435
00436 if ( myTkRef->extra().isNonnull() ) {
00437 trackInnPos.push_back( myTkRef->innerPosition());
00438 trackPin.push_back( myTkRef->innerMomentum());
00439 trackPout.push_back( myTkRef->outerMomentum());
00440 }
00441 trackPairRef.push_back(myTkRef);
00442
00443 }
00444
00445
00446
00447
00448
00449 if( theConversionVertex.isValid() ) {
00450
00451
00452
00453 }
00454
00455
00456
00457 minAppDist=calculateMinApproachDistance( trackPairRef[0], trackPairRef[1]);
00458
00459 double like = -999.;
00460 reco::Conversion newCandidate(scPtrVec, trackPairRef, trkPositionAtEcal, theConversionVertex, matchingBC, minAppDist, trackInnPos, trackPin, trackPout, like, algo);
00461
00462 like = theLikelihoodCalc_->calculateLikelihood( newCandidate );
00463
00464 newCandidate.setMVAout(like);
00465 outputConvPhotonCollection.push_back(newCandidate);
00466
00467
00468 myCands++;
00469
00470
00471 } else {
00472
00473
00474
00475
00476
00477 trackPairRef.clear();
00478 trackInnPos.clear();
00479 trackPin.clear();
00480 trackPout.clear();
00481 std::vector<reco::TransientTrack>::const_iterator iTk=(iPair->first).begin();
00482
00483 const reco::TrackTransientTrack* ttt = dynamic_cast<const reco::TrackTransientTrack*>(iTk->basicTransientTrack());
00484 reco::TrackRef myTk= ttt->persistentTrackRef();
00485 if ( myTk->extra().isNonnull() ) {
00486 trackInnPos.push_back( myTk->innerPosition());
00487 trackPin.push_back( myTk->innerMomentum());
00488 trackPout.push_back( myTk->outerMomentum());
00489 }
00490 trackPairRef.push_back(myTk);
00491
00492
00493 if ( recoverOneTrackCase_ ) {
00494 float theta1 = myTk->innerMomentum().Theta();
00495 float dCot=999.;
00496 float dCotTheta=-999.;
00497 reco::TrackRef goodRef;
00498 std::vector<reco::TransientTrack>::const_iterator iGoodGenTran;
00499 for ( std::vector<reco::TransientTrack>::const_iterator iTran= t_generalTrk.begin(); iTran != t_generalTrk.end(); ++iTran) {
00500 const reco::TrackTransientTrack* ttt = dynamic_cast<const reco::TrackTransientTrack*>(iTran->basicTransientTrack());
00501 reco::TrackRef trRef= ttt->persistentTrackRef();
00502 if ( trRef->charge()*myTk->charge() > 0 ) continue;
00503 float dEta = trRef->eta() - myTk->eta();
00504 float dPhi = trRef->phi() - myTk->phi();
00505 if ( sqrt (dEta*dEta + dPhi*dPhi) > dRForConversionRecovery_ ) continue;
00506 float theta2 = trRef->innerMomentum().Theta();
00507 dCotTheta = 1./tan(theta1) - 1./tan(theta2) ;
00508
00509 if ( fabs(dCotTheta) < dCot ) {
00510 dCot = fabs(dCotTheta);
00511 goodRef = trRef;
00512 iGoodGenTran=iTran;
00513 }
00514 }
00515
00516 if ( goodRef.isNonnull() ) {
00517
00518 minAppDist=calculateMinApproachDistance( myTk, goodRef);
00519
00520
00521 if ( fabs(dCotTheta) < deltaCotCut_ && minAppDist > minApproachDisCut_ ) {
00522 trackInnPos.push_back( goodRef->innerPosition());
00523 trackPin.push_back( goodRef->innerMomentum());
00524 trackPout.push_back( goodRef->outerMomentum());
00525 trackPairRef.push_back( goodRef );
00526
00527
00528 std::vector<reco::TransientTrack> mypair;
00529 mypair.push_back(*iTk);
00530 mypair.push_back(*iGoodGenTran);
00531
00532 try{
00533 theVertexFinder_->run(iPair->first, theConversionVertex );
00534
00535 }
00536 catch ( cms::Exception& e ) {
00537
00538 edm::LogWarning(metname) << "cms::Exception caught in ConvertedPhotonProducer::produce\n"
00539 << e.explainSelf();
00540
00541 }
00542 }
00543
00544 }
00545
00546 }
00547 double like = -999.;
00548 reco::Conversion newCandidate(scPtrVec, trackPairRef, trkPositionAtEcal, theConversionVertex, matchingBC, minAppDist, trackInnPos, trackPin, trackPout, like, algo);
00549 like = theLikelihoodCalc_->calculateLikelihood(newCandidate);
00550 newCandidate.setMVAout(like);
00551 outputConvPhotonCollection.push_back(newCandidate);
00552
00553
00554
00555
00556 }
00557
00558
00559
00560
00561 }
00562
00563 }
00564
00565
00566
00567
00568 }
00569
00570
00571
00572
00573
00574 }
00575
00576
00577 void ConvertedPhotonProducer::cleanCollections(const edm::Handle<edm::View<reco::CaloCluster> > & scHandle,
00578 const edm::OrphanHandle<reco::ConversionCollection> & conversionHandle,
00579 reco::ConversionCollection & outputConversionCollection) {
00580
00581
00582 reco::Conversion* newCandidate=0;
00583 for(unsigned int lSC=0; lSC < scHandle->size(); lSC++) {
00584
00585
00586 reco::CaloClusterPtr aClus= scHandle->ptrAt(lSC);
00587
00588
00589 if (aClus->energy()/cosh(aClus->eta()) <= minSCEt_) continue;
00590
00591
00592 if ( conversionHandle.isValid() ) {
00593
00594 if ( risolveAmbiguity_ ) {
00595 std::vector<reco::ConversionRef> bestRef=solveAmbiguity( conversionHandle , aClus);
00596
00597 for ( std::vector<reco::ConversionRef>::iterator iRef=bestRef.begin(); iRef!=bestRef.end(); iRef++ ) {
00598 if ( iRef->isNonnull() ) {
00599 newCandidate= (*iRef)->clone();
00600 outputConversionCollection.push_back(*newCandidate);
00601 delete newCandidate;
00602 }
00603 }
00604
00605 } else {
00606
00607
00608 for( unsigned int icp = 0; icp < conversionHandle->size(); icp++) {
00609 reco::ConversionRef cpRef(reco::ConversionRef(conversionHandle,icp));
00610 if (!( aClus.id() == cpRef->caloCluster()[0].id() && aClus.key() == cpRef->caloCluster()[0].key() )) continue;
00611 if ( !cpRef->isConverted() ) continue;
00612 if ( cpRef->nTracks() <2 ) continue;
00613 newCandidate= (&(*cpRef))->clone();
00614 outputConversionCollection.push_back(*newCandidate);
00615 delete newCandidate;
00616
00617 }
00618
00619 }
00620
00621 }
00622
00623
00624 }
00625 }
00626
00627
00628
00629
00630 std::vector<reco::ConversionRef> ConvertedPhotonProducer::solveAmbiguity(const edm::OrphanHandle<reco::ConversionCollection> & conversionHandle, reco::CaloClusterPtr& scRef) {
00631 std::multimap<double, reco::ConversionRef, std::greater<double> > convMap;
00632
00633 for ( unsigned int icp=0; icp< conversionHandle->size(); icp++) {
00634 reco::ConversionRef cpRef(reco::ConversionRef(conversionHandle,icp));
00635
00636
00637 if (!( scRef.id() == cpRef->caloCluster()[0].id() && scRef.key() == cpRef->caloCluster()[0].key() )) continue;
00638 if ( !cpRef->isConverted() ) continue;
00639 double like = cpRef->MVAout();
00640 if ( cpRef->nTracks() <2 ) continue;
00641
00642 convMap.insert ( std::make_pair(like,cpRef) );
00643
00644 }
00645
00646
00647
00648 std::multimap<double, reco::ConversionRef>::iterator iMap;
00649 std::vector<reco::ConversionRef> bestRefs;
00650 for (iMap=convMap.begin(); iMap!=convMap.end(); iMap++) {
00651
00652 bestRefs.push_back( iMap->second );
00653 if ( int(bestRefs.size()) == maxNumOfCandidates_ ) break;
00654 }
00655
00656
00657 return bestRefs;
00658
00659
00660 }
00661
00662
00663
00664
00665
00666 float ConvertedPhotonProducer::calculateMinApproachDistance ( const reco::TrackRef& track1, const reco::TrackRef& track2) {
00667 float dist=9999.;
00668
00669 double x1, x2, y1, y2;
00670 double xx_1 = track1->innerPosition().x(), yy_1 = track1->innerPosition().y(), zz_1 = track1->innerPosition().z();
00671 double xx_2 = track2->innerPosition().x(), yy_2 = track2->innerPosition().y(), zz_2 = track2->innerPosition().z();
00672 double radius1 = track1->innerMomentum().Rho()/(.3*(theMF_->inTesla(GlobalPoint(xx_1, yy_1, zz_1)).z()))*100;
00673 double radius2 = track2->innerMomentum().Rho()/(.3*(theMF_->inTesla(GlobalPoint(xx_2, yy_2, zz_2)).z()))*100;
00674 getCircleCenter(track1, radius1, x1, y1);
00675 getCircleCenter(track2, radius2, x2, y2);
00676 dist = sqrt((x1-x2)*(x1-x2)+(y1-y2)*(y1-y2)) - radius1 - radius2;
00677
00678 return dist;
00679
00680 }
00681
00682
00683 void ConvertedPhotonProducer::getCircleCenter(const reco::TrackRef& tk, double r, double& x0, double& y0){
00684 double x1, y1, phi;
00685 x1 = tk->innerPosition().x();
00686 y1 = tk->innerPosition().y();
00687 phi = tk->innerMomentum().phi();
00688 const int charge = tk->charge();
00689 x0 = x1 + r*sin(phi)*charge;
00690 y0 = y1 - r*cos(phi)*charge;
00691
00692 }