CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoParticleFlow/PFProducer/src/PFBlockAlgo.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFProducer/interface/PFBlockAlgo.h"
00002 #include "RecoParticleFlow/PFProducer/interface/Utils.h"
00003 #include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
00004 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
00005 #include "DataFormats/TrackReco/interface/Track.h"
00006 #include "DataFormats/ParticleFlowReco/interface/PFDisplacedVertex.h" // gouzevitch
00007 
00008 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
00009 
00010 #include <stdexcept>
00011 #include "TMath.h"
00012 
00013 using namespace std;
00014 using namespace reco;
00015 
00016 //for debug only 
00017 //#define PFLOW_DEBUG
00018 
00019 const PFBlockAlgo::Mask PFBlockAlgo::dummyMask_;
00020 
00021 PFBlockAlgo::PFBlockAlgo() : 
00022   blocks_( new reco::PFBlockCollection ),
00023   DPtovPtCut_(std::vector<double>(4,static_cast<double>(999.))),
00024   NHitCut_(std::vector<unsigned int>(4,static_cast<unsigned>(0))), 
00025   useIterTracking_(true),
00026   photonSelector_(0),
00027   debug_(false) {}
00028 
00029 
00030 
00031 void PFBlockAlgo::setParameters( std::vector<double>& DPtovPtCut,
00032                                  std::vector<unsigned int>& NHitCut,
00033                                  bool useConvBremPFRecTracks,
00034                                  bool useIterTracking,
00035                                  int nuclearInteractionsPurity,
00036                                  bool useEGPhotons,
00037                                  std::vector<double>& photonSelectionCuts) {
00038   
00039   DPtovPtCut_    = DPtovPtCut;
00040   NHitCut_       = NHitCut;
00041   useIterTracking_ = useIterTracking;
00042   useConvBremPFRecTracks_ = useConvBremPFRecTracks;
00043   nuclearInteractionsPurity_ = nuclearInteractionsPurity;
00044   useEGPhotons_ = useEGPhotons;
00045   // Pt cut; Track iso (constant + slope), Ecal iso (constant + slope), HCAL iso (constant+slope), H/E
00046   if(useEGPhotons_)
00047     photonSelector_ = new PhotonSelectorAlgo(photonSelectionCuts[0],   
00048                                              photonSelectionCuts[1], photonSelectionCuts[2],   
00049                                              photonSelectionCuts[3], photonSelectionCuts[4],    
00050                                              photonSelectionCuts[5], photonSelectionCuts[6],    
00051                                              photonSelectionCuts[7]);
00052 }
00053 
00054 PFBlockAlgo::~PFBlockAlgo() {
00055 
00056 #ifdef PFLOW_DEBUG
00057   if(debug_)
00058     cout<<"~PFBlockAlgo - number of remaining elements: "
00059         <<elements_.size()<<endl;
00060 #endif
00061 
00062   if(photonSelector_) delete photonSelector_;
00063 }
00064 
00065 void 
00066 PFBlockAlgo::findBlocks() {
00067 
00068   //  cout<<"findBlocks : "<<blocks_.get()<<endl;
00069   
00070   // the blocks have not been passed to the event, and need to be cleared
00071   if(blocks_.get() )blocks_->clear();
00072   else 
00073     blocks_.reset( new reco::PFBlockCollection );
00074 
00075   blocks_->reserve(elements_.size());
00076   for(IE ie = elements_.begin(); 
00077       ie != elements_.end();) {
00078     
00079 #ifdef PFLOW_DEBUG
00080     if(debug_) {
00081       cout<<" PFBlockAlgo::findBlocks() ----------------------"<<endl;
00082       cout<<" element "<<**ie<<endl;
00083       cout<<" creating new block"<<endl;
00084     }
00085 #endif
00086     
00087     blocks_->push_back( PFBlock() );
00088     
00089     vector< PFBlockLink > links;
00090 
00091     //    list< IE > used;
00092     ie = associate( elements_.end() , ie, links );
00093 
00094     // build remaining links in current block
00095     packLinks( blocks_->back(), links );
00096   }       
00097 }
00098 
00099 
00100 
00101 
00102 PFBlockAlgo::IE 
00103 PFBlockAlgo::associate( IE last, 
00104                         IE next, 
00105                         vector<PFBlockLink>& links ) {
00106 
00107     
00108 #ifdef PFLOW_DEBUG
00109   if(debug_ ) cout<<"PFBlockAlgo::associate start ----"<<endl;
00110 #endif
00111 
00112   if( last!= elements_.end() ) {
00113     PFBlockLink::Type linktype = PFBlockLink::NONE;
00114     double dist = -1;
00115     PFBlock::LinkTest linktest = PFBlock::LINKTEST_RECHIT;
00116     link( *last, *next, linktype, linktest, dist ); 
00117 
00118    
00119     if(dist<-0.5) {
00120 #ifdef PFLOW_DEBUG
00121       if(debug_ ) cout<<"link failed"<<endl;
00122 #endif
00123       return ++next; // association failed
00124     }
00125     else {
00126       // add next element to the current pflowblock
00127       blocks_->back().addElement( *next );  
00128 
00129       // (*next)->setIndex( blocks_->back()->indexToLastElement() );
00130       
00131       // this is not necessary? 
00132       // next->setPFBlock(this);
00133       
00134       // create a link between next and last
00135       links.push_back( PFBlockLink(linktype, 
00136                                    linktest,
00137                                    dist,
00138                                    (*last)->index(), 
00139                                    (*next)->index() ) );
00140       // not necessary ?
00141       //       next->connect( links_.size()-1 );
00142       //       last->connect( links_.size()-1 );      
00143     }
00144   }
00145   else {
00146     // add next element to this eflowblock
00147 #ifdef PFLOW_DEBUG
00148     if(debug_ ) cout<<"adding to block element "<<(**next)<<endl;
00149 #endif
00150     blocks_->back().addElement( *next );
00151     // (*next)->setIndex( blocks_->back()->indexToLastElement() );   
00152     // next->setPFBlock(this);
00153   }
00154 
00155   // recursive call: associate next and other unused elements 
00156   
00157   //   IE afterNext = next;
00158   //   ++afterNext;
00159   //  cout<<"last "<<**last<<" next "<<**next<<endl;
00160   
00161   for(IE ie = elements_.begin(); 
00162       ie != elements_.end();) {
00163     
00164     if( ie == last || ie == next ) {
00165       ++ie;
00166       continue;
00167     } 
00168     
00169     // *ie already included to a block
00170     if( (*ie)->locked() ) {
00171 #ifdef PFLOW_DEBUG
00172       if(debug_ ) cout<<"element "<<(**ie)<<"already used"<<endl;
00173 #endif
00174       ++ie;
00175       continue;
00176     }    
00177     
00178     
00179 #ifdef PFLOW_DEBUG
00180     if(debug_ ) cout<<"calling associate "<<(**next)<<" & "<<(**ie)<<endl;
00181 #endif
00182     ie = associate(next, ie, links);
00183   }       
00184 
00185 #ifdef PFLOW_DEBUG
00186   if(debug_ ) {
00187     cout<<"**** deleting element "<<endl;
00188     cout<<**next<<endl;
00189   }
00190 #endif
00191   delete *next;
00192 
00193 #ifdef PFLOW_DEBUG
00194   if(debug_ ) {
00195     cout<<"**** removing element "<<endl;
00196   }
00197 #endif
00198 
00199   IE iteratorToNextFreeElement = elements_.erase( next );
00200 
00201 #ifdef PFLOW_DEBUG
00202   if(debug_ ) cout<<"PFBlockAlgo::associate stop ----"<<endl;
00203 #endif
00204 
00205   return iteratorToNextFreeElement;
00206 }
00207 
00208 
00209 
00210 void 
00211 PFBlockAlgo::packLinks( reco::PFBlock& block, 
00212                         const vector<PFBlockLink>& links ) const {
00213   
00214   
00215   const edm::OwnVector< reco::PFBlockElement >& els = block.elements();
00216   
00217   block.bookLinkData();
00218 
00219   //First Loop: update all link data
00220   for( unsigned i1=0; i1<els.size(); i1++ ) {
00221     for( unsigned i2=0; i2<els.size(); i2++ ) {
00222       
00223       // no reflexive link
00224       if( i1==i2 ) continue;
00225       
00226       double dist = -1;
00227       
00228       bool linked = false;
00229       PFBlock::LinkTest linktest 
00230         = PFBlock::LINKTEST_RECHIT; 
00231 
00232       // are these elements already linked ?
00233       // this can be optimized
00234 
00235       for( unsigned il=0; il<links.size(); il++ ) {
00236         if( (links[il].element1() == i1 && 
00237              links[il].element2() == i2) || 
00238             (links[il].element1() == i2 && 
00239              links[il].element2() == i1) ) { // yes
00240           
00241           dist = links[il].dist();
00242           linked = true;
00243 
00244           //modif-beg     
00245           //retrieve type of test used to get distance
00246           linktest = links[il].test();
00247 #ifdef PFLOW_DEBUG
00248           if( debug_ )
00249             cout << "Reading link vector: linktest used=" 
00250                  << linktest 
00251                  << " distance = " << dist 
00252                  << endl; 
00253 #endif
00254           //modif-end
00255           
00256           break;
00257         } 
00258       }
00259       
00260       if(!linked) {
00261         PFBlockLink::Type linktype = PFBlockLink::NONE;
00262         link( & els[i1], & els[i2], linktype, linktest, dist);
00263       }
00264 
00265       //loading link data according to link test used: RECHIT 
00266       //block.setLink( i1, i2, chi2, block.linkData() );
00267 #ifdef PFLOW_DEBUG
00268       if( debug_ )
00269         cout << "Setting link between elements " << i1 << " and " << i2
00270              << " of dist =" << dist << " computed from link test "
00271              << linktest << endl;
00272 #endif
00273       block.setLink( i1, i2, dist, block.linkData(), linktest );
00274     }
00275   }
00276 
00277   // Do not cut the link between the primary track and the clusters. It would be analysed in the PFCandConnector.cc
00278   // checkDisplacedVertexLinks( block );
00279 }
00280 
00281 
00282 
00283 void 
00284 PFBlockAlgo::buildGraph() {
00285   // loop on all blocks and create a big graph
00286 }
00287 
00288 
00289 
00290 void 
00291 PFBlockAlgo::link( const reco::PFBlockElement* el1, 
00292                    const reco::PFBlockElement* el2, 
00293                    PFBlockLink::Type& linktype, 
00294                    reco::PFBlock::LinkTest& linktest,
00295                    double& dist) const {
00296   
00297 
00298 
00299   dist=-1.;
00300   linktest = PFBlock::LINKTEST_RECHIT; //rechit by default 
00301 
00302   PFBlockElement::Type type1 = el1->type();
00303   PFBlockElement::Type type2 = el2->type();
00304 
00305   if( type1==type2 ) {
00306     // cannot link 2 elements of the same type. 
00307     // except if the elements are 2 tracks or 2 ECAL
00308     if( type1!=PFBlockElement::TRACK && type1!=PFBlockElement::GSF &&
00309         type1!=PFBlockElement::ECAL) {
00310       return;
00311     }
00312 
00313     // cannot link two primary tracks  (except if they come from a V0)
00314     if( type1 ==PFBlockElement::TRACK) {
00315       if ( !el1->isLinkedToDisplacedVertex() || !el2->isLinkedToDisplacedVertex()) 
00316       return;
00317     }
00318   }
00319 
00320   linktype = static_cast<PFBlockLink::Type>
00321     ((1<< (type1-1) ) | (1<< (type2-1) ));
00322 
00323   if(debug_ ) std::cout << " PFBlockAlgo links type1 " << type1 << " type2 " << type2 << std::endl;
00324 
00325   PFBlockElement::Type lowType = type1;
00326   PFBlockElement::Type highType = type2;
00327   const PFBlockElement* lowEl = el1;
00328   const PFBlockElement* highEl = el2;
00329   
00330   if(type1>type2) {
00331     lowType = type2;
00332     highType = type1;
00333     lowEl = el2;
00334     highEl = el1;
00335   }
00336   
00337   switch(linktype) {
00338   case PFBlockLink::TRACKandPS1:
00339   case PFBlockLink::TRACKandPS2:
00340     {
00341       //       cout<<"TRACKandPS"<<endl;
00342       PFRecTrackRef trackref = lowEl->trackRefPF();
00343       PFClusterRef  clusterref = highEl->clusterRef();
00344       assert( !trackref.isNull() );
00345       assert( !clusterref.isNull() );
00346       // PJ - 14-May-09 : A link by rechit is needed here !
00347       dist = testTrackAndPS( *trackref, *clusterref );
00348       linktest = PFBlock::LINKTEST_RECHIT;
00349       break;
00350     }
00351     
00352   case PFBlockLink::TRACKandECAL:
00353     {
00354       if(debug_ ) cout<<"TRACKandECAL"<<endl;
00355       PFRecTrackRef trackref = lowEl->trackRefPF();
00356       
00357       if(debug_ ) std::cout << " Track pt " << trackref->trackRef()->pt() << std::endl;
00358       
00359       PFClusterRef  clusterref = highEl->clusterRef();
00360       assert( !trackref.isNull() );
00361       assert( !clusterref.isNull() );
00362       dist = LinkByRecHit::testTrackAndClusterByRecHit( *trackref, *clusterref, false, debug_ );
00363       linktest = PFBlock::LINKTEST_RECHIT;
00364 
00365       if ( debug_ ) { 
00366         if( dist > 0. ) { 
00367           std::cout << " Here a link has been established"
00368                     << " between a track an Ecal with dist  " 
00369                     << dist <<  std::endl;
00370         } else {
00371           std::cout << " No link found " << std::endl;
00372         }
00373       }
00374 
00375       break;
00376     }
00377   case PFBlockLink::TRACKandHCAL:
00378     {
00379       //       cout<<"TRACKandHCAL"<<endl;
00380       PFRecTrackRef trackref = lowEl->trackRefPF();
00381       PFClusterRef  clusterref = highEl->clusterRef();
00382       assert( !trackref.isNull() );
00383       assert( !clusterref.isNull() );
00384       dist = LinkByRecHit::testTrackAndClusterByRecHit( *trackref, *clusterref, false, debug_ );
00385       linktest = PFBlock::LINKTEST_RECHIT;      
00386       break;
00387     }
00388   case PFBlockLink::ECALandHCAL:
00389     {
00390       //       cout<<"ECALandHCAL"<<endl;
00391       PFClusterRef  ecalref = lowEl->clusterRef();
00392       PFClusterRef  hcalref = highEl->clusterRef();
00393       assert( !ecalref.isNull() );
00394       assert( !hcalref.isNull() );
00395       // PJ - 14-May-09 : A link by rechit is needed here !
00396       // dist = testECALAndHCAL( *ecalref, *hcalref );
00397       dist = -1.;
00398       linktest = PFBlock::LINKTEST_RECHIT;
00399       break;
00400     }
00401   case PFBlockLink::PS1andECAL:
00402   case PFBlockLink::PS2andECAL:
00403     {
00404       //       cout<<"PSandECAL"<<endl;
00405       PFClusterRef  psref = lowEl->clusterRef();
00406       PFClusterRef  ecalref = highEl->clusterRef();
00407       assert( !psref.isNull() );
00408       assert( !ecalref.isNull() );
00409       dist = LinkByRecHit::testECALAndPSByRecHit( *ecalref, *psref ,debug_);
00410       linktest = PFBlock::LINKTEST_RECHIT;      
00411       break;
00412     }
00413   case PFBlockLink::PS1andPS2:
00414     {
00415       PFClusterRef  ps1ref = lowEl->clusterRef();
00416       PFClusterRef  ps2ref = highEl->clusterRef();
00417       assert( !ps1ref.isNull() );
00418       assert( !ps2ref.isNull() );
00419       // PJ - 14-May-09 : A link by rechit is needed here !
00420       // dist = testPS1AndPS2( *ps1ref, *ps2ref );
00421       dist = -1.;
00422       linktest = PFBlock::LINKTEST_RECHIT;
00423       break;
00424     }
00425   case PFBlockLink::TRACKandTRACK:
00426     {
00427       if(debug_ ) 
00428         cout<<"TRACKandTRACK"<<endl;
00429       dist = testLinkByVertex(lowEl, highEl);
00430       if(debug_ ) 
00431         std::cout << " PFBlockLink::TRACKandTRACK dist " << dist << std::endl;
00432       linktest = PFBlock::LINKTEST_RECHIT;
00433       break;
00434     }
00435 
00436   case PFBlockLink::ECALandECAL:
00437       {
00438         
00439         PFClusterRef  ecal1ref = lowEl->clusterRef();
00440         PFClusterRef  ecal2ref = highEl->clusterRef();
00441         assert( !ecal1ref.isNull() );
00442         assert( !ecal2ref.isNull() );
00443         if(debug_)
00444           cout << " PFBlockLink::ECALandECAL" << endl;
00445         dist = testLinkBySuperCluster(ecal1ref,ecal2ref);
00446         break;
00447       }
00448 
00449   case PFBlockLink::ECALandGSF:
00450     {
00451       PFClusterRef  clusterref = lowEl->clusterRef();
00452       assert( !clusterref.isNull() );
00453       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00454       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00455       dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, false, debug_ );
00456       linktest = PFBlock::LINKTEST_RECHIT;
00457       
00458       if ( debug_ ) {
00459         if ( dist > 0. ) {
00460           std::cout << " Here a link has been established" 
00461                     << " between a GSF track an Ecal with dist  " 
00462                     << dist <<  std::endl;
00463         } else {
00464           if(debug_ ) std::cout << " No link found " << std::endl;
00465         }
00466       }
00467       break;
00468     }
00469   case PFBlockLink::TRACKandGSF:
00470     {
00471       PFRecTrackRef trackref = lowEl->trackRefPF();
00472       assert( !trackref.isNull() );
00473       const reco::PFBlockElementGsfTrack *  GsfEl =  
00474         dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00475       GsfPFRecTrackRef gsfref = GsfEl->GsftrackRefPF();
00476       reco::TrackRef kftrackref= (*trackref).trackRef();
00477       assert( !gsfref.isNull() );
00478       PFRecTrackRef refkf = (*gsfref).kfPFRecTrackRef();
00479       if(refkf.isNonnull()) {
00480         reco::TrackRef gsftrackref = (*refkf).trackRef();
00481         if (gsftrackref.isNonnull()&&kftrackref.isNonnull()) {
00482           if (kftrackref == gsftrackref) { 
00483             dist = 0.001;
00484           } else { 
00485             dist = -1.;
00486           }
00487         } else { 
00488           dist = -1.;
00489         }
00490       } else {
00491         dist = -1.;
00492       }
00493       
00494       
00495       if(useConvBremPFRecTracks_) {
00496         if(lowEl->isLinkedToDisplacedVertex()){
00497           vector<PFRecTrackRef> pfrectrack_vec = GsfEl->GsftrackRefPF()->convBremPFRecTrackRef();
00498           if(pfrectrack_vec.size() > 0){
00499             for(unsigned int iconv = 0; iconv <  pfrectrack_vec.size(); iconv++) {
00500               if( lowEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)) {
00501                 // use track ref
00502                 if(kftrackref == (*pfrectrack_vec[iconv]).trackRef()) {         
00503                   dist = 0.001;
00504                 }
00505               } 
00506               else{
00507                 // use the track base ref
00508                 reco::TrackBaseRef newTrackBaseRef((*pfrectrack_vec[iconv]).trackRef());
00509                 reco::TrackBaseRef elemTrackBaseRef(kftrackref);              
00510                 if(newTrackBaseRef == elemTrackBaseRef){
00511                   dist = 0.001;
00512                 } 
00513               }
00514             }
00515           }
00516         }
00517       }
00518  
00519 
00520       break;      
00521     }
00522          
00523   case PFBlockLink::GSFandBREM:
00524     {
00525       const reco::PFBlockElementGsfTrack * GsfEl  =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(lowEl);
00526       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00527       GsfPFRecTrackRef gsfref = GsfEl->GsftrackRefPF();
00528       GsfPFRecTrackRef bremref = BremEl->GsftrackRefPF();
00529       assert( !gsfref.isNull() );
00530       assert( !bremref.isNull() );
00531       if (gsfref == bremref)  { 
00532         dist = 0.001;
00533       } else { 
00534         dist = -1.;
00535       }
00536       break;
00537     }
00538   case PFBlockLink::GSFandGSF:
00539     {
00540       const reco::PFBlockElementGsfTrack * lowGsfEl  =  
00541         dynamic_cast<const reco::PFBlockElementGsfTrack*>(lowEl);
00542       const reco::PFBlockElementGsfTrack * highGsfEl  =  
00543         dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00544       
00545       GsfPFRecTrackRef lowgsfref = lowGsfEl->GsftrackRefPF();
00546       GsfPFRecTrackRef highgsfref = highGsfEl->GsftrackRefPF();
00547       assert( !lowgsfref.isNull() );
00548       assert( !highgsfref.isNull() );
00549       
00550       if( (lowGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) == false && 
00551            highGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)) ||
00552           (highGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) == false && 
00553            lowGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV))) {
00554         if(lowgsfref->trackId() == highgsfref->trackId()) {
00555           dist = 0.001;
00556         }
00557         else {
00558           dist = -1.;
00559         }
00560       }
00561       break;
00562     }
00563   case PFBlockLink::ECALandBREM:
00564     {
00565       PFClusterRef  clusterref = lowEl->clusterRef();
00566       assert( !clusterref.isNull() );
00567       const reco::PFBlockElementBrem * BremEl =  
00568         dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00569       const PFRecTrack * myTrack = &(BremEl->trackPF());
00570       /*
00571       double DP = (BremEl->DeltaP())*(-1.);
00572       double SigmaDP = BremEl->SigmaDeltaP();
00573       double SignBremDp = DP/SigmaDP;
00574       */
00575       bool isBrem = true;
00576       dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, isBrem, debug_);
00577       if( debug_ && dist > 0. ) 
00578         std::cout << "ECALandBREM: dist testTrackAndClusterByRecHit " 
00579                   << dist << std::endl;
00580       linktest = PFBlock::LINKTEST_RECHIT;
00581       break;
00582     }
00583   case PFBlockLink::PS1andGSF:
00584   case PFBlockLink::PS2andGSF:
00585     {
00586       PFClusterRef  psref = lowEl->clusterRef();
00587       assert( !psref.isNull() );
00588       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00589       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00590       // PJ - 14-May-09 : A link by rechit is needed here !
00591       dist = testTrackAndPS( *myTrack, *psref );
00592       linktest = PFBlock::LINKTEST_RECHIT;
00593       break;
00594     }
00595   case PFBlockLink::PS1andBREM:
00596   case PFBlockLink::PS2andBREM:
00597     {
00598       PFClusterRef  psref = lowEl->clusterRef();
00599       assert( !psref.isNull() );
00600       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00601       const PFRecTrack * myTrack = &(BremEl->trackPF());
00602       // PJ - 14-May-09 : A link by rechit is needed here !
00603       dist = testTrackAndPS( *myTrack, *psref );
00604       linktest = PFBlock::LINKTEST_RECHIT;
00605       break;
00606     }
00607   case PFBlockLink::HCALandGSF:
00608     {
00609       PFClusterRef  clusterref = lowEl->clusterRef();
00610       assert( !clusterref.isNull() );
00611       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00612       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00613       dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, false, debug_ );
00614       linktest = PFBlock::LINKTEST_RECHIT;
00615       break;
00616     }
00617   case PFBlockLink::HCALandBREM:
00618     {
00619       PFClusterRef  clusterref = lowEl->clusterRef();
00620       assert( !clusterref.isNull() );
00621       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00622       const PFRecTrack * myTrack = &(BremEl->trackPF());
00623       bool isBrem = true;
00624       dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, isBrem, debug_);
00625       break;
00626     }
00627   case PFBlockLink::HFEMandHFHAD:
00628     {
00629       // cout<<"HFEMandHFHAD"<<endl;
00630       PFClusterRef  eref = lowEl->clusterRef();
00631       PFClusterRef  href = highEl->clusterRef();
00632       assert( !eref.isNull() );
00633       assert( !href.isNull() );
00634       dist = LinkByRecHit::testHFEMAndHFHADByRecHit( *eref, *href, debug_ );
00635       linktest = PFBlock::LINKTEST_RECHIT;
00636       break;
00637       
00638     }
00639   case PFBlockLink::SCandECAL:
00640     {
00641       PFClusterRef  clusterref = lowEl->clusterRef();
00642 
00643       assert( !clusterref.isNull() );
00644       
00645       const reco::PFBlockElementSuperCluster * scEl = 
00646         dynamic_cast<const reco::PFBlockElementSuperCluster*>(highEl);
00647       assert (!scEl->superClusterRef().isNull());
00648       dist = testSuperClusterPFCluster(scEl->superClusterRef(),
00649                                        clusterref);
00650       break;
00651     }
00652   default:
00653     dist = -1.;
00654     linktest = PFBlock::LINKTEST_RECHIT;
00655     // cerr<<"link type not implemented yet: 0x"<<hex<<linktype<<dec<<endl;
00656     // assert(0);
00657     return;
00658   }
00659 }
00660 
00661 double
00662 PFBlockAlgo::testTrackAndPS(const PFRecTrack& track, 
00663                             const PFCluster& ps)  const {
00664 
00665 #ifdef PFLOW_DEBUG
00666   //   cout<<"entering testTrackAndPS"<<endl;
00667   // resolution of PS cluster dxdx and dydy from strip pitch and length
00668   double dx=0.;
00669   double dy=0.;
00670   
00671   unsigned layerid =0;
00672   // PS1: vertical strips  PS2: horizontal strips
00673   switch (ps.layer()) {
00674   case PFLayer::PS1:
00675     layerid = reco::PFTrajectoryPoint::PS1;
00676     
00677     // vertical strips in PS1, measure x with pitch precision
00678     dx = resPSpitch_;
00679     dy = resPSlength_; 
00680     break;
00681   case PFLayer::PS2:
00682     layerid = reco::PFTrajectoryPoint::PS2;
00683     // horizontal strips in PS2, measure y with pitch precision
00684     dy = resPSpitch_;
00685     dx = resPSlength_;
00686     break;
00687   default:
00688     break;
00689   }
00690   const reco::PFTrajectoryPoint& atPS
00691     = track.extrapolatedPoint( layerid );  
00692   // did not reach PS, cannot be associated with a cluster.
00693   if( ! atPS.isValid() ) return -1.;   
00694   
00695   double trackx = atPS.position().X();
00696   double tracky = atPS.position().Y();
00697   double trackz = atPS.position().Z(); // MDN jan 09
00698   
00699   // ps position  x, y
00700   double psx = ps.position().X();
00701   double psy = ps.position().Y();
00702   // MDN Jan 09: check that trackz and psz have the same sign
00703   double psz = ps.position().Z();
00704   if( trackz*psz < 0.) return -1.; 
00705   
00706   // double chi2 = (psx-trackx)*(psx-trackx)/(dx*dx + trackresolx*trackresolx)
00707   //  + (psy-tracky)*(psy-tracky)/(dy*dy + trackresoly*trackresoly);
00708 
00709   double dist = std::sqrt( (psx-trackx)*(psx-trackx)
00710                          + (psy-tracky)*(psy-tracky));  
00711   if(debug_) cout<<"testTrackAndPS "<< dist <<" "<<endl;
00712   if(debug_){
00713     cout<<" trackx " << trackx 
00714         <<" tracky " << tracky 
00715         <<" psx "    <<  psx   
00716         <<" psy "    << psy    
00717         << endl;
00718   }
00719 #endif
00720   
00721   // Return -1. as long as no link by rechit is available
00722   return -1.;
00723 }
00724 
00725 double
00726 PFBlockAlgo::testECALAndHCAL(const PFCluster& ecal, 
00727                              const PFCluster& hcal)  const {
00728   
00729   //   cout<<"entering testECALAndHCAL"<<endl;
00730   
00731   /*
00732   double dist = 
00733     computeDist( ecal.positionREP().Eta(),
00734                  ecal.positionREP().Phi(), 
00735                  hcal.positionREP().Eta(), 
00736                  hcal.positionREP().Phi() );
00737   */
00738 
00739 #ifdef PFLOW_DEBUG
00740   if(debug_) cout<<"testECALAndHCAL "<< dist <<" "<<endl;
00741   if(debug_){
00742     cout<<" ecaleta " << ecal.positionREP().Eta()
00743         <<" ecalphi " << ecal.positionREP().Phi()
00744         <<" hcaleta " << hcal.positionREP().Eta()
00745         <<" hcalphi " << hcal.positionREP().Phi()
00746   }
00747 #endif
00748 
00749   // Need to implement a link by RecHit
00750   return -1.;
00751 }
00752 
00753 double
00754 PFBlockAlgo::testLinkBySuperCluster(const PFClusterRef& ecal1, 
00755                                     const PFClusterRef& ecal2)  const {
00756   
00757   //  cout<<"entering testECALAndECAL "<< pfcRefSCMap_.size() << endl;
00758   
00759   double dist = -1;
00760   
00761   // the first one is not in any super cluster
00762   int testindex=pfcSCVec_[ecal1.key()];
00763   if(testindex == -1.) return dist;
00764   //  if(itcheck==pfcRefSCMap_.end()) return dist;
00765   // now retrieve the of PFclusters in this super cluster  
00766 
00767   const std::vector<reco::PFClusterRef> & thePFClusters(scpfcRefs_[testindex]);
00768   
00769   unsigned npf=thePFClusters.size();
00770   for(unsigned i=0;i<npf;++i)
00771     {
00772       if(thePFClusters[i]==ecal2) // yes they are in the same SC 
00773         {
00774           dist=LinkByRecHit::computeDist( ecal1->positionREP().Eta(),
00775                                           ecal1->positionREP().Phi(), 
00776                                           ecal2->positionREP().Eta(), 
00777                                           ecal2->positionREP().Phi() );
00778 //        std::cout << " DETA " << fabs(ecal1->positionREP().Eta()-ecal2->positionREP().Eta()) << std::endl;
00779 //        if(fabs(ecal1->positionREP().Eta()-ecal2->positionREP().Eta())>0.2)
00780 //          {
00781 //            std::cout <<  " Super Cluster " <<  *(superClusters_[testindex]) << std::endl;
00782 //            std::cout <<  " Cluster1 " <<  *ecal1 << std::endl;
00783 //            std::cout <<  " Cluster2 " <<  *ecal2 << std::endl;
00784 //            ClusterClusterMapping::checkOverlap(*ecal1,superClusters_,0.01,true);
00785 //            ClusterClusterMapping::checkOverlap(*ecal2,superClusters_,0.01,true);
00786 //          }
00787           return dist;
00788         }
00789     }
00790   return dist;
00791 }
00792 
00793 
00794 double
00795 PFBlockAlgo::testSuperClusterPFCluster(const SuperClusterRef& ecal1, 
00796                                        const PFClusterRef& ecal2)  const {
00797   
00798   //  cout<<"entering testECALAndECAL "<< pfcRefSCMap_.size() << endl;
00799   
00800   double dist = -1;
00801   
00802   bool overlap=ClusterClusterMapping::overlap(*ecal1,*ecal2);
00803   
00804   if(overlap)   {
00805     dist=LinkByRecHit::computeDist( ecal1->position().eta(),
00806                                     ecal1->position().phi(), 
00807                                     ecal2->positionREP().Eta(), 
00808                                     ecal2->positionREP().Phi() );
00809     return dist;
00810   }
00811   return dist;
00812 }
00813 
00814 
00815 
00816 double
00817 PFBlockAlgo::testPS1AndPS2(const PFCluster& ps1, 
00818                            const PFCluster& ps2)  const {
00819   
00820 #ifdef PFLOW_DEBUG
00821   //   cout<<"entering testPS1AndPS2"<<endl;
00822   
00823   // compute chi2 in y, z using swimming formulae
00824   // y2 = y1 * z2/z1   and x2 = x1 *z2/z1
00825   
00826   // ps position1  x, y, z
00827   double x1 = ps1.position().X();
00828   double y1 = ps1.position().Y();
00829   double z1 = ps1.position().Z();
00830   double x2 = ps2.position().X();
00831   double y2 = ps2.position().Y();
00832   double z2 = ps2.position().Z();
00833   // MDN Bug correction Jan 09: check that z1 and z2 have the same sign!
00834   if (z1*z2<0.) -1.;
00835   // swim to PS2
00836   double scale = z2/z1;
00837   double x1atPS2 = x1*scale;
00838   double y1atPS2 = y1*scale;
00839   // resolution of PS cluster dxdx and dydy from strip pitch and length
00840   // vertical strips in PS1, measure x with pitch precision
00841   double dx1dx1 = resPSpitch_*resPSpitch_*scale*scale;
00842   double dy1dy1 = resPSlength_*resPSlength_*scale*scale;
00843   // horizontal strips in PS2 , measure y with pitch precision
00844   double dy2dy2 = resPSpitch_*resPSpitch_;
00845   double dx2dx2 = resPSlength_*resPSlength_;
00846   
00847   // double chi2 = (x2-x1atPS2)*(x2-x1atPS2)/(dx1dx1 + dx2dx2) 
00848   //  + (y2-y1atPS2)*(y2-y1atPS2)/(dy1dy1 + dy2dy2);
00849   
00850   double dist = std::sqrt( (x2-x1atPS2)*(x2-x1atPS2)
00851                          + (y2-y1atPS2)*(y2-y1atPS2));
00852     
00853   if(debug_) cout<<"testPS1AndPS2 "<<dist<<" "<<endl;
00854   if(debug_){
00855     cout<<" x1atPS2 "<< x1atPS2 << " dx1 "<<resPSpitch_*scale
00856         <<" y1atPS2 "<< y1atPS2 << " dy1 "<<resPSlength_*scale<< endl
00857         <<" x2 " <<x2  << " dx2 "<<resPSlength_
00858         <<" y2 " << y2 << " dy2 "<<resPSpitch_<< endl;
00859   }
00860 #endif
00861 
00862   // Need a link by rechit here
00863   return -1.; 
00864 }
00865 
00866 
00867 
00868 double
00869 PFBlockAlgo::testLinkByVertex( const reco::PFBlockElement* elt1, 
00870                                const reco::PFBlockElement* elt2) const {
00871 
00872   //  cout << "Test link by vertex between" << endl << *elt1 << endl << " and " << endl << *elt2 << endl;
00873 
00874   double result=-1.;
00875 
00876   reco::PFBlockElement::TrackType T_TO_DISP = reco::PFBlockElement::T_TO_DISP;
00877   reco::PFBlockElement::TrackType T_FROM_DISP = reco::PFBlockElement::T_FROM_DISP;
00878   PFDisplacedTrackerVertexRef ni1_TO_DISP = elt1->displacedVertexRef(T_TO_DISP);
00879   PFDisplacedTrackerVertexRef ni2_TO_DISP = elt2->displacedVertexRef(T_TO_DISP);
00880   PFDisplacedTrackerVertexRef ni1_FROM_DISP = elt1->displacedVertexRef(T_FROM_DISP);
00881   PFDisplacedTrackerVertexRef ni2_FROM_DISP = elt2->displacedVertexRef(T_FROM_DISP);
00882   
00883   if( ni1_TO_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
00884     if( ni1_TO_DISP == ni2_FROM_DISP ) { result = 1.0; return result; }
00885 
00886   if( ni1_FROM_DISP.isNonnull() && ni2_TO_DISP.isNonnull())
00887     if( ni1_FROM_DISP == ni2_TO_DISP ) { result = 1.0; return result; }
00888 
00889   if( ni1_FROM_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
00890     if( ni1_FROM_DISP == ni2_FROM_DISP ) { result = 1.0; return result; }
00891     
00892   
00893   if (  elt1->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)  &&
00894              elt2->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)  ) {
00895     
00896     if(debug_ ) std::cout << " testLinkByVertex On Conversions " << std::endl;
00897     
00898     if ( elt1->convRef().isNonnull() && elt2->convRef().isNonnull() ) {
00899       if(debug_ ) std::cout << " PFBlockAlgo.cc testLinkByVertex  Cconversion Refs are non null  " << std::endl;      
00900       if ( elt1->convRef() ==  elt2->convRef() ) {
00901         result=1.0;
00902         if(debug_ ) std::cout << " testLinkByVertex  Cconversion Refs are equal  " << std::endl;      
00903         return result;
00904       }
00905     } 
00906     
00907   }
00908   
00909   if (  elt1->trackType(reco::PFBlockElement::T_FROM_V0)  &&
00910              elt2->trackType(reco::PFBlockElement::T_FROM_V0)  ) {
00911     if(debug_ ) std::cout << " testLinkByVertex On V0 " << std::endl;
00912     if ( elt1->V0Ref().isNonnull() && elt2->V0Ref().isNonnull() ) {
00913       if(debug_ ) std::cout << " PFBlockAlgo.cc testLinkByVertex  V0 Refs are non null  " << std::endl;
00914       if ( elt1->V0Ref() ==  elt2->V0Ref() ) {
00915         result=1.0;
00916         if(debug_ ) std::cout << " testLinkByVertex  V0 Refs are equal  " << std::endl;
00917         return result;
00918       }
00919     }
00920   }
00921 
00922   return result;
00923 }
00924 
00925 
00926 
00927 void 
00928 PFBlockAlgo::checkMaskSize( const reco::PFRecTrackCollection& tracks,
00929                             const reco::GsfPFRecTrackCollection& gsftracks, 
00930                             const reco::PFClusterCollection&  ecals,
00931                             const reco::PFClusterCollection&  hcals,
00932                             const reco::PFClusterCollection&  hfems,
00933                             const reco::PFClusterCollection&  hfhads,
00934                             const reco::PFClusterCollection&  pss, 
00935                             const reco::PhotonCollection&  egphh, 
00936                             const Mask& trackMask, 
00937                             const Mask& gsftrackMask,  
00938                             const Mask& ecalMask,
00939                             const Mask& hcalMask, 
00940                             const Mask& hfemMask,
00941                             const Mask& hfhadMask,                    
00942                             const Mask& psMask,
00943                             const Mask& phMask) const {
00944 
00945   if( !trackMask.empty() && 
00946       trackMask.size() != tracks.size() ) {
00947     string err = "PFBlockAlgo::setInput: ";
00948     err += "The size of the track mask is different ";
00949     err += "from the size of the track vector.";
00950     throw std::length_error( err.c_str() );
00951   }
00952 
00953   if( !gsftrackMask.empty() && 
00954       gsftrackMask.size() != gsftracks.size() ) {
00955     string err = "PFBlockAlgo::setInput: ";
00956     err += "The size of the gsf track mask is different ";
00957     err += "from the size of the gsftrack vector.";
00958     throw std::length_error( err.c_str() );
00959   }
00960 
00961   if( !ecalMask.empty() && 
00962       ecalMask.size() != ecals.size() ) {
00963     string err = "PFBlockAlgo::setInput: ";
00964     err += "The size of the ecal mask is different ";
00965     err += "from the size of the ecal clusters vector.";
00966     throw std::length_error( err.c_str() );
00967   }
00968   
00969   if( !hcalMask.empty() && 
00970       hcalMask.size() != hcals.size() ) {
00971     string err = "PFBlockAlgo::setInput: ";
00972     err += "The size of the hcal mask is different ";
00973     err += "from the size of the hcal clusters vector.";
00974     throw std::length_error( err.c_str() );
00975   }
00976 
00977   if( !hfemMask.empty() && 
00978       hfemMask.size() != hfems.size() ) {
00979     string err = "PFBlockAlgo::setInput: ";
00980     err += "The size of the hfem mask is different ";
00981     err += "from the size of the hfem clusters vector.";
00982     throw std::length_error( err.c_str() );
00983   }
00984 
00985   if( !hfhadMask.empty() && 
00986       hfhadMask.size() != hfhads.size() ) {
00987     string err = "PFBlockAlgo::setInput: ";
00988     err += "The size of the hfhad mask is different ";
00989     err += "from the size of the hfhad clusters vector.";
00990     throw std::length_error( err.c_str() );
00991   }
00992 
00993   if( !psMask.empty() && 
00994       psMask.size() != pss.size() ) {
00995     string err = "PFBlockAlgo::setInput: ";
00996     err += "The size of the ps mask is different ";
00997     err += "from the size of the ps clusters vector.";
00998     throw std::length_error( err.c_str() );
00999   }
01000   
01001     if( !phMask.empty() && 
01002       phMask.size() != egphh.size() ) {
01003     string err = "PFBlockAlgo::setInput: ";
01004     err += "The size of the photon mask is different ";
01005     err += "from the size of the photon vector.";
01006     throw std::length_error( err.c_str() );
01007   }
01008 
01009 }
01010 
01011 
01012 std::ostream& operator<<(std::ostream& out, const PFBlockAlgo& a) {
01013   if(! out) return out;
01014   
01015   out<<"====== Particle Flow Block Algorithm ======= ";
01016   out<<endl;
01017   out<<"number of unassociated elements : "<<a.elements_.size()<<endl;
01018   out<<endl;
01019   
01020   for(PFBlockAlgo::IEC ie = a.elements_.begin(); 
01021       ie != a.elements_.end(); ie++) {
01022     out<<"\t"<<**ie <<endl;
01023   }
01024 
01025   
01026   //   const PFBlockCollection& blocks = a.blocks();
01027 
01028   const std::auto_ptr< reco::PFBlockCollection >& blocks
01029     = a.blocks(); 
01030     
01031   if(!blocks.get() ) {
01032     out<<"blocks already transfered"<<endl;
01033   }
01034   else {
01035     out<<"number of blocks : "<<blocks->size()<<endl;
01036     out<<endl;
01037     
01038     for(PFBlockAlgo::IBC ib=blocks->begin(); 
01039         ib != blocks->end(); ib++) {
01040       out<<(*ib)<<endl;
01041     }
01042   }
01043 
01044   return out;
01045 }
01046 
01047 bool 
01048 PFBlockAlgo::goodPtResolution( const reco::TrackRef& trackref) {
01049 
01050   double P = trackref->p();
01051   double Pt = trackref->pt();
01052   double DPt = trackref->ptError();
01053   unsigned int NHit = trackref->hitPattern().trackerLayersWithMeasurement();
01054   unsigned int NLostHit = trackref->hitPattern().trackerLayersWithoutMeasurement();
01055   unsigned int LostHits = trackref->numberOfLostHits();
01056   double sigmaHad = sqrt(1.20*1.20/P+0.06*0.06) / (1.+LostHits);
01057 
01058   // iteration 1,2,3,4,5 correspond to algo = 1/4,5,6,7,8,9
01059   unsigned int Algo = 0; 
01060   switch (trackref->algo()) {
01061   case TrackBase::ctf:
01062   case TrackBase::iter0:
01063   case TrackBase::iter1:
01064     Algo = 0;
01065     break;
01066   case TrackBase::iter2:
01067     Algo = 1;
01068     break;
01069   case TrackBase::iter3:
01070     Algo = 2;
01071     break;
01072   case TrackBase::iter4:
01073     Algo = 3;
01074     break;
01075   case TrackBase::iter5:
01076     Algo = 4;
01077     break;
01078   default:
01079     Algo = useIterTracking_ ? 5 : 0;
01080     break;
01081   }
01082 
01083   // Protection against 0 momentum tracks
01084   if ( P < 0.05 ) return false;
01085 
01086   // Temporary : Reject all tracking iteration beyond 5th step. 
01087   if ( Algo > 4 ) return false;
01088  
01089   if (debug_) cout << " PFBlockAlgo: PFrecTrack->Track Pt= "
01090                    << Pt << " DPt = " << DPt << endl;
01091   if ( ( DPtovPtCut_[Algo] > 0. && 
01092          DPt/Pt > DPtovPtCut_[Algo]*sigmaHad ) || 
01093        NHit < NHitCut_[Algo] ) { 
01094     // (Algo >= 3 && LostHits != 0) ) {
01095     if (debug_) cout << " PFBlockAlgo: skip badly measured track"
01096                      << ", P = " << P 
01097                      << ", Pt = " << Pt 
01098                      << " DPt = " << DPt 
01099                      << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01100                      << ", Algo = " << Algo
01101                      << endl;
01102     if (debug_) cout << " cut is DPt/Pt < " << DPtovPtCut_[Algo] * sigmaHad << endl;
01103     if (debug_) cout << " cut is NHit >= " << NHitCut_[Algo] << endl;
01104     /*
01105     std::cout << "Track REJECTED : ";
01106     std::cout << ", P = " << P 
01107               << ", Pt = " << Pt 
01108               << " DPt = " << DPt 
01109               << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01110               << ", Algo = " << Algo
01111               << std::endl;
01112     */
01113     return false;
01114   }
01115 
01116   /*
01117   std::cout << "Track Accepted : ";
01118   std::cout << ", P = " << P 
01119        << ", Pt = " << Pt 
01120        << " DPt = " << DPt 
01121        << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01122        << ", Algo = " << Algo
01123        << std::endl;
01124   */
01125   return true;
01126 }
01127 
01128 int
01129 PFBlockAlgo::muAssocToTrack( const reco::TrackRef& trackref,
01130                              const edm::Handle<reco::MuonCollection>& muonh) const {
01131   if(muonh.isValid() ) {
01132     for(unsigned j=0;j<muonh->size(); j++) {
01133       reco::MuonRef muonref( muonh, j );
01134       if (muonref->track().isNonnull()) 
01135         if( muonref->track() == trackref ) return j;
01136     }
01137   }
01138   return -1; // not found
01139 }
01140 
01141 int 
01142 PFBlockAlgo::muAssocToTrack( const reco::TrackRef& trackref,
01143                              const edm::OrphanHandle<reco::MuonCollection>& muonh) const {
01144   if(muonh.isValid() ) {
01145     for(unsigned j=0;j<muonh->size(); j++) {
01146       reco::MuonRef muonref( muonh, j );
01147       if (muonref->track().isNonnull())
01148         if( muonref->track() == trackref ) return j;
01149     }
01150   }
01151   return -1; // not found
01152 }
01153 
01154 
01155 void 
01156 PFBlockAlgo::checkDisplacedVertexLinks( reco::PFBlock& block ) const {
01157   // method which removes link between primary tracks and the clusters
01158   
01159   typedef std::multimap<double, unsigned>::iterator IE;
01160 
01161   const edm::OwnVector< reco::PFBlockElement >& els = block.elements();
01162   // loop on all elements != TRACK
01163   for( unsigned i1=0; i1 != els.size(); ++i1 ) {
01164     if( els[i1].type() == PFBlockElement::TRACK ) continue;
01165     std::multimap<double, unsigned> assocTracks;
01166     // get associated tracks
01167     block.associatedElements( i1,  block.linkData(),
01168                               assocTracks,
01169                               reco::PFBlockElement::TRACK,
01170                               reco::PFBlock::LINKTEST_ALL );
01171     for( IE ie = assocTracks.begin(); ie != assocTracks.end(); ++ie) {
01172       //double   distprim  = ie->first;
01173       unsigned iprim     = ie->second;
01174       // if this track a primary track (T_TO_DISP)
01175       // the new strategy gouzevitch: remove all the links from primary track
01176       if( els[iprim].isPrimary()) {
01177 
01178             block.setLink( i1, iprim, -1, block.linkData(),
01179                            PFBlock::LINKTEST_RECHIT );      
01180       }
01181     } // loop on all associated tracks
01182   } // loop on all elements
01183  
01184 }
01185 
01186 void PFBlockAlgo::fillFromPhoton(const reco::Photon & photon, reco::PFBlockElementSuperCluster * pfbe) {
01187   pfbe->setTrackIso(photon.trkSumPtHollowConeDR04());
01188   pfbe->setEcalIso(photon.ecalRecHitSumEtConeDR04());
01189   pfbe->setHcalIso(photon.hcalTowerSumEtConeDR04());
01190   pfbe->setHoE(photon.hadronicOverEm());
01191   pfbe->setFromPhoton(true);
01192 }