CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/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                                              photonSelectionCuts[8],
00053                                              photonSelectionCuts[9],
00054                                              photonSelectionCuts[10]
00055                                              );
00056 
00057 
00058 }
00059 
00060 // Glowinski & Gouzevitch
00061 void PFBlockAlgo::setUseOptimization(bool useKDTreeTrackEcalLinker)
00062 {
00063   useKDTreeTrackEcalLinker_ = useKDTreeTrackEcalLinker;
00064 }
00065 // !Glowinski & Gouzevitch
00066 
00067 
00068 PFBlockAlgo::~PFBlockAlgo() {
00069 
00070 #ifdef PFLOW_DEBUG
00071   if(debug_)
00072     cout<<"~PFBlockAlgo - number of remaining elements: "
00073         <<elements_.size()<<endl;
00074 #endif
00075 
00076   if(photonSelector_) delete photonSelector_;
00077 
00078 }
00079 
00080 void 
00081 PFBlockAlgo::findBlocks() {
00082 
00083   // Glowinski & Gouzevitch
00084   if (useKDTreeTrackEcalLinker_) {
00085     TELinker_.process();
00086     THLinker_.process();
00087     PSELinker_.process();
00088   }
00089   // !Glowinski & Gouzevitch
00090 
00091   // the blocks have not been passed to the event, and need to be cleared
00092   if(blocks_.get() )blocks_->clear();
00093   else 
00094     blocks_.reset( new reco::PFBlockCollection );
00095 
00096   blocks_->reserve(elements_.size());
00097   for(IE ie = elements_.begin(); 
00098       ie != elements_.end();) {
00099     
00100 #ifdef PFLOW_DEBUG
00101     if(debug_) {
00102       cout<<" PFBlockAlgo::findBlocks() ----------------------"<<endl;
00103       cout<<" element "<<**ie<<endl;
00104       cout<<" creating new block"<<endl;
00105     }
00106 #endif
00107 
00108     blocks_->push_back( PFBlock() );
00109     
00110     vector< PFBlockLink > links;
00111 
00112     //    list< IE > used;
00113     ie = associate( elements_.end() , ie, links );
00114 
00115     packLinks( blocks_->back(), links );
00116   }
00117 }
00118 
00119 
00120 
00121 
00122 PFBlockAlgo::IE 
00123 PFBlockAlgo::associate( IE last, 
00124                         IE next, 
00125                         vector<PFBlockLink>& links ) {
00126 
00127     
00128 #ifdef PFLOW_DEBUG
00129   if(debug_ ) cout<<"PFBlockAlgo::associate start ----"<<endl;
00130 #endif
00131 
00132   if( last!= elements_.end() ) {
00133     PFBlockLink::Type linktype = PFBlockLink::NONE;
00134     double dist = -1;
00135     PFBlock::LinkTest linktest = PFBlock::LINKTEST_RECHIT;
00136     link( *last, *next, linktype, linktest, dist ); 
00137 
00138     
00139     if(dist<-0.5) {
00140 #ifdef PFLOW_DEBUG
00141       if(debug_ ) cout<<"link failed"<<endl;
00142 #endif
00143       return ++next; // association failed
00144     }
00145     else {
00146       // add next element to the current pflowblock
00147       blocks_->back().addElement( *next );  
00148 
00149       // (*next)->setIndex( blocks_->back()->indexToLastElement() );
00150       
00151       // this is not necessary? 
00152       // next->setPFBlock(this);
00153       
00154       // create a link between next and last      
00155       links.push_back( PFBlockLink(linktype, 
00156                                    linktest,
00157                                    dist,
00158                                    (*last)->index(), 
00159                                    (*next)->index() ) );
00160       // not necessary ?
00161       //       next->connect( links_.size()-1 );
00162       //       last->connect( links_.size()-1 );      
00163     }
00164   }
00165   else {
00166     // add next element to this eflowblock
00167 #ifdef PFLOW_DEBUG
00168     if(debug_ ) cout<<"adding to block element "<<(**next)<<endl;
00169 #endif
00170     blocks_->back().addElement( *next );
00171     // (*next)->setIndex( blocks_->back()->indexToLastElement() );   
00172     // next->setPFBlock(this);
00173   }
00174 
00175   // recursive call: associate next and other unused elements 
00176   
00177   //   IE afterNext = next;
00178   //   ++afterNext;
00179   //  cout<<"last "<<**last<<" next "<<**next<<endl;
00180   
00181   for(IE ie = elements_.begin(); 
00182       ie != elements_.end();) {
00183     
00184     if( ie == last || ie == next ) {
00185       ++ie;
00186       continue;
00187     } 
00188     bool bTestLink = linkPrefilter(*next, *ie);
00189 
00190     if(!bTestLink) {
00191       ++ie;
00192       continue;
00193     }
00194     
00195     // *ie already included to a block
00196     if( (*ie)->locked() ) {
00197 #ifdef PFLOW_DEBUG
00198       if(debug_ ) cout<<"element "<<(**ie)<<"already used"<<endl;
00199 #endif
00200       ++ie;
00201       continue;
00202     }    
00203     
00204     
00205 #ifdef PFLOW_DEBUG
00206     if(debug_ ) cout<<"calling associate "<<(**next)<<" & "<<(**ie)<<endl;
00207 #endif
00208     ie = associate(next, ie, links);
00209   }       
00210 
00211 #ifdef PFLOW_DEBUG
00212   if(debug_ ) {
00213     cout<<"**** deleting element "<<endl;
00214     cout<<**next<<endl;
00215   }
00216 #endif
00217   delete *next;
00218 
00219 #ifdef PFLOW_DEBUG
00220   if(debug_ ) {
00221     cout<<"**** removing element "<<endl;
00222   }
00223 #endif
00224 
00225   IE iteratorToNextFreeElement = elements_.erase( next );
00226 
00227 #ifdef PFLOW_DEBUG
00228   if(debug_ ) cout<<"PFBlockAlgo::associate stop ----"<<endl;
00229 #endif
00230 
00231   return iteratorToNextFreeElement;
00232 }
00233 
00234 
00235 
00236 void 
00237 PFBlockAlgo::packLinks( reco::PFBlock& block, 
00238                         const vector<PFBlockLink>& links ) const {
00239   
00240   
00241   const edm::OwnVector< reco::PFBlockElement >& els = block.elements();
00242   
00243   block.bookLinkData();
00244 
00245   unsigned elsize = els.size();
00246   unsigned ilStart = 0;
00247   //First Loop: update all link data
00248   for( unsigned i1=0; i1<elsize; ++i1 ) {
00249     for( unsigned i2=0; i2<i1; ++i2 ) {
00250       
00251       // no reflexive link
00252       //if( i1==i2 ) continue;
00253       
00254       double dist = -1;
00255       
00256       bool linked = false;
00257       PFBlock::LinkTest linktest 
00258         = PFBlock::LINKTEST_RECHIT; 
00259 
00260       // are these elements already linked ?
00261       // this can be optimized
00262 
00263       unsigned linksize = links.size();
00264       for( unsigned il = ilStart; il<linksize; ++il ) {
00265         // The following three lines exploits the increasing-element2 ordering of links.
00266         if ( links[il].element2() < i1 ) ilStart = il;
00267         if ( links[il].element2() > i1 ) break;
00268         if( (links[il].element1() == i2 && 
00269              links[il].element2() == i1) ) { // yes
00270           
00271           dist = links[il].dist();
00272           linked = true;
00273 
00274           //modif-beg     
00275           //retrieve type of test used to get distance
00276           linktest = links[il].test();
00277 #ifdef PFLOW_DEBUG
00278           if( debug_ )
00279             cout << "Reading link vector: linktest used=" 
00280                  << linktest 
00281                  << " distance = " << dist 
00282                  << endl; 
00283 #endif
00284           //modif-end
00285           
00286           break;
00287         } 
00288       }
00289       
00290       if(!linked) {
00291         PFBlockLink::Type linktype = PFBlockLink::NONE;
00292         bool bTestLink = linkPrefilter(&els[i1], &els[i2]);
00293         if (bTestLink) link( & els[i1], & els[i2], linktype, linktest, dist);
00294       }
00295 
00296       //loading link data according to link test used: RECHIT 
00297       //block.setLink( i1, i2, chi2, block.linkData() );
00298 #ifdef PFLOW_DEBUG
00299       if( debug_ )
00300         cout << "Setting link between elements " << i1 << " and " << i2
00301              << " of dist =" << dist << " computed from link test "
00302              << linktest << endl;
00303 #endif
00304       block.setLink( i1, i2, dist, block.linkData(), linktest );
00305     }
00306   }
00307 
00308 }
00309 
00310 
00311 
00312 void 
00313 PFBlockAlgo::buildGraph() {
00314   // loop on all blocks and create a big graph
00315 }
00316 
00317 
00318 
00319 void 
00320 PFBlockAlgo::link( const reco::PFBlockElement* el1, 
00321                    const reco::PFBlockElement* el2, 
00322                    PFBlockLink::Type& linktype, 
00323                    reco::PFBlock::LinkTest& linktest,
00324                    double& dist) const {
00325   
00326   // ACHTUNG!!!! If you introduce new links check that they are not desabled in linkPrefilter!!!!
00327 
00328 
00329   dist=-1.;
00330   linktest = PFBlock::LINKTEST_RECHIT; //rechit by default 
00331 
00332   PFBlockElement::Type type1 = el1->type();
00333   PFBlockElement::Type type2 = el2->type();
00334 
00335   linktype = static_cast<PFBlockLink::Type>
00336     ((1<< (type1-1) ) | (1<< (type2-1) ));
00337 
00338   if(debug_ ) std::cout << " PFBlockAlgo links type1 " << type1 << " type2 " << type2 << std::endl;
00339 
00340   const PFBlockElement* lowEl = el1;
00341   const PFBlockElement* highEl = el2;
00342   
00343   if(type1>type2) {
00344     lowEl = el2;
00345     highEl = el1;
00346   }
00347   
00348   switch(linktype) {
00349   // Track and preshower cluster links are not used for now - disable
00350   case PFBlockLink::PS1andECAL:
00351   case PFBlockLink::PS2andECAL:
00352     {
00353       // if(debug_ ) cout<< "PSandECAL" <<endl;
00354       PFClusterRef  psref = lowEl->clusterRef();
00355       PFClusterRef  ecalref = highEl->clusterRef();
00356       assert( !psref.isNull() );
00357       assert( !ecalref.isNull() );
00358 
00359       // Check if the linking has been done using the KDTree algo
00360       // Glowinski & Gouzevitch
00361       if ( useKDTreeTrackEcalLinker_ && lowEl->isMultilinksValide() ) { // KDTree algo
00362         const reco::PFMultilinksType& multilinks = lowEl->getMultilinks();
00363         
00364         double ecalPhi = ecalref->positionREP().Phi();
00365         double ecalEta = ecalref->positionREP().Eta();
00366         
00367         // Check if the link PS/Ecal exist
00368         reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
00369         for (; mlit != multilinks.end(); ++mlit)
00370           if ((mlit->first == ecalPhi) && (mlit->second == ecalEta))
00371             break;
00372         
00373         // If the link exist, we fill dist and linktest. We use old algorithme method.
00374         if (mlit != multilinks.end()){
00375           double xPS = psref->position().X();
00376           double yPS = psref->position().Y();
00377           double xECAL  = ecalref->position().X();
00378           double yECAL  = ecalref->position().Y();
00379 
00380           dist = LinkByRecHit::computeDist(xECAL/1000.,yECAL/1000.,xPS/1000.  ,yPS/1000, false);
00381         }
00382 
00383       } else { //Old algorithm
00384         dist = LinkByRecHit::testECALAndPSByRecHit( *ecalref, *psref ,debug_);
00385       }
00386 
00387       //      linktest = PFBlock::LINKTEST_RECHIT;
00388       
00389       break;
00390     }
00391   case PFBlockLink::TRACKandECAL:
00392     {
00393       if(debug_ ) cout<<"TRACKandECAL"<<endl;
00394 
00395       PFRecTrackRef trackref = lowEl->trackRefPF();
00396       PFClusterRef  clusterref = highEl->clusterRef();
00397       assert( !trackref.isNull() );
00398       assert( !clusterref.isNull() );
00399 
00400       if(debug_ ) std::cout << " Track pt " << trackref->trackRef()->pt() << std::endl;
00401 
00402       // Check if the linking has been done using the KDTree algo
00403       // Glowinski & Gouzevitch
00404       if ( useKDTreeTrackEcalLinker_ && lowEl->isMultilinksValide() ) { //KDTree Algo
00405 
00406         const reco::PFMultilinksType& multilinks = lowEl->getMultilinks();
00407         double ecalphi = clusterref->positionREP().Phi();
00408         double ecaleta = clusterref->positionREP().Eta();
00409 
00410         // Check if the link Track/Ecal exist
00411         reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
00412         for (; mlit != multilinks.end(); ++mlit)
00413           if ((mlit->first == ecalphi) && (mlit->second == ecaleta))
00414             break;
00415 
00416         
00417         // If the link exist, we fill dist and linktest. We use old algorithme method.
00418         if (mlit != multilinks.end()){
00419 
00420 
00421           //Should be something like this :
00422           //      const reco::PFRecTrack& track = *trackref;
00423           //instead of this :
00424           /*
00425           reco::PFRecTrack track (*trackref);
00426           const reco::PFTrajectoryPoint& atECAL_tmp = 
00427             (*trackref).extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax );
00428           if(std::abs(atECAL_tmp.positionREP().Eta())<1E-9 &&
00429              std::abs(atECAL_tmp.positionREP().Phi())<1E-9 &&
00430              atECAL_tmp.positionREP().R()<1E-9) 
00431             track.calculatePositionREP();
00432           */
00433 
00434           const reco::PFTrajectoryPoint& atECAL = 
00435             trackref->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax );
00436           
00437           double tracketa = atECAL.positionREP().Eta();
00438           double trackphi = atECAL.positionREP().Phi();
00439 
00440           dist = LinkByRecHit::computeDist(ecaleta, ecalphi, tracketa, trackphi);
00441         }
00442 
00443       } else {// Old algorithm
00444         if ( trackref->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax ).isValid() )
00445           dist = LinkByRecHit::testTrackAndClusterByRecHit( *trackref, *clusterref, false, debug_ );
00446         else
00447           dist = -1.;
00448       }
00449 
00450       if ( debug_ ) { 
00451         if( dist > 0. ) { 
00452           std::cout << " Here a link has been established"
00453                     << " between a track an Ecal with dist  " 
00454                     << dist <<  std::endl;
00455         } else
00456           std::cout << " No link found " << std::endl;
00457       }
00458       
00459       //     linktest = PFBlock::LINKTEST_RECHIT;
00460       break;
00461     }
00462   case PFBlockLink::TRACKandHCAL:
00463     {
00464       //      if(debug_ ) cout<<"TRACKandHCAL"<<endl;
00465 
00466       PFRecTrackRef trackref = lowEl->trackRefPF();
00467       PFClusterRef  clusterref = highEl->clusterRef();
00468       assert( !trackref.isNull() );
00469       assert( !clusterref.isNull() );
00470 
00471       // Check if the linking has been done using the KDTree algo
00472       // Glowinski & Gouzevitch
00473       if ( useKDTreeTrackEcalLinker_ && highEl->isMultilinksValide() ) { //KDTree Algo
00474         
00475         const reco::PFMultilinksType& multilinks = highEl->getMultilinks();
00476 
00477         /*
00478         reco::PFRecTrack track (*trackref);
00479         const reco::PFTrajectoryPoint& atHCALEntrance_tmp = 
00480           (*trackref).extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance);
00481         if (std::abs(atHCALEntrance_tmp.positionREP().Eta())<1E-9 &&
00482             std::abs(atHCALEntrance_tmp.positionREP().Phi())<1E-9 &&
00483             atHCALEntrance_tmp.positionREP().R()<1E-9)   
00484           track.calculatePositionREP();
00485         */
00486 
00487         const reco::PFTrajectoryPoint& atHCAL = 
00488           trackref->extrapolatedPoint(reco::PFTrajectoryPoint::HCALEntrance);
00489         
00490           
00491         double tracketa = atHCAL.positionREP().Eta();
00492         double trackphi = atHCAL.positionREP().Phi();
00493 
00494         // Check if the link Track/Ecal exist
00495         reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
00496         for (; mlit != multilinks.end(); ++mlit)
00497           if ((mlit->first == trackphi) && (mlit->second == tracketa))
00498             break;
00499 
00500         // If the link exist, we fill dist and linktest. We use old algorithme method.
00501         if (mlit != multilinks.end()){
00502 
00503           const reco::PFTrajectoryPoint& atHCALExit = 
00504             trackref->extrapolatedPoint(reco::PFTrajectoryPoint::HCALExit);
00505           double dHEta = 0.0;
00506           double dHPhi = 0.0;
00507           if (atHCALExit.position().R()>atHCAL.position().R()) {
00508             dHEta = atHCALExit.positionREP().Eta()-atHCAL.positionREP().Eta();
00509             dHPhi = atHCALExit.positionREP().Phi()-atHCAL.positionREP().Phi(); 
00510             if ( dHPhi > M_PI ) dHPhi = dHPhi - 2.*M_PI;
00511             else if ( dHPhi < -M_PI ) dHPhi = dHPhi + 2.*M_PI; 
00512           } else {
00513             std::cout << "Qu'est ce que c'est que ce gag ? " 
00514                       << atHCALExit.position().R() << " is larger than " 
00515                       << atHCAL.position().R() << " !" << std::endl;
00516           }
00517 
00518           tracketa += 0.1 * dHEta;
00519           trackphi += 0.1 * dHPhi;
00520 
00521           double clusterphi = clusterref->positionREP().Phi();
00522           double clustereta = clusterref->positionREP().Eta();
00523           
00524           dist = LinkByRecHit::computeDist(clustereta, clusterphi, tracketa, trackphi);
00525         }       
00526 
00527       } else {// Old algorithm
00528         if ( trackref->extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance ).isValid() )
00529           dist = LinkByRecHit::testTrackAndClusterByRecHit( *trackref, *clusterref, false, debug_ );
00530         else
00531           dist = -1.;
00532       }
00533       
00534       //      linktest = PFBlock::LINKTEST_RECHIT;
00535       break;
00536     }
00537   case PFBlockLink::TRACKandHO:
00538     {
00539       if(debug_ ) cout<<"TRACKandHO"<<endl;
00540 
00541       PFRecTrackRef trackref = lowEl->trackRefPF();
00542       PFClusterRef  clusterref = highEl->clusterRef();
00543       
00544       assert( !trackref.isNull() );
00545       assert( !clusterref.isNull() );
00546 
00547       
00548       // Old algorithm
00549       //      cout<<"TRACKandHO1"<<trackref->pt()<<" "<<trackref->eta()<<" "<<trackref->phi()<<endl;
00550       //Same value is used in PFTrackTransformer::addPoints() for HOLayer, but allow for some rounding precision
00551       if ( lowEl->trackRef()->pt() > 3.00001 && trackref->extrapolatedPoint( reco::PFTrajectoryPoint::HOLayer ).isValid() ) {
00552         //      cout<<"TRACKandHO2"<<endl;
00553         dist = LinkByRecHit::testTrackAndClusterByRecHit( *trackref, *clusterref, false, debug_ );
00554         
00555         //      cout <<"dist TRACKandHO "<<dist<<endl;
00556       } else {
00557         dist = -1.;
00558       }
00559       //      linktest = PFBlock::LINKTEST_RECHIT;
00560       break;
00561     }
00562   case PFBlockLink::ECALandHCAL:
00563     {
00564       //       cout<<"ECALandHCAL"<<endl;
00565       PFClusterRef  ecalref = lowEl->clusterRef();
00566       PFClusterRef  hcalref = highEl->clusterRef();
00567       assert( !ecalref.isNull() );
00568       assert( !hcalref.isNull() );
00569       // PJ - 14-May-09 : A link by rechit is needed here !
00570       dist = testECALAndHCAL( *ecalref, *hcalref );
00571       // dist = -1.;
00572       //     linktest = PFBlock::LINKTEST_RECHIT;
00573       break;
00574     }
00575   case PFBlockLink::HCALandHO:
00576     {
00577       //       cout<<"HCALandH0"<<endl;
00578       PFClusterRef  hcalref = lowEl->clusterRef();
00579       PFClusterRef  horef = highEl->clusterRef();
00580       assert( !hcalref.isNull() );
00581       assert( !horef.isNull() );
00582       // PJ - 14-May-09 : A link by rechit is needed here !
00583       dist = testHCALAndHO( *hcalref, *horef );
00584       // dist = -1.;
00585       //      cout <<"Dist "<<dist<<endl;
00586       //     linktest = PFBlock::LINKTEST_RECHIT;
00587       break;
00588     }
00589   case PFBlockLink::HFEMandHFHAD:
00590     {
00591       // cout<<"HFEMandHFHAD"<<endl;
00592       PFClusterRef  eref = lowEl->clusterRef();
00593       PFClusterRef  href = highEl->clusterRef();
00594       assert( !eref.isNull() );
00595       assert( !href.isNull() );
00596       dist = LinkByRecHit::testHFEMAndHFHADByRecHit( *eref, *href, debug_ );
00597       //    linktest = PFBlock::LINKTEST_RECHIT;
00598       break;      
00599     }
00600 
00601   case PFBlockLink::TRACKandTRACK:
00602     {
00603       if(debug_ ) 
00604         cout<<"TRACKandTRACK"<<endl;
00605       dist = testLinkByVertex(lowEl, highEl);
00606       if(debug_ ) 
00607         std::cout << " PFBlockLink::TRACKandTRACK dist " << dist << std::endl;
00608       //   linktest = PFBlock::LINKTEST_RECHIT;
00609       break;
00610     }
00611 
00612   case PFBlockLink::ECALandECAL:
00613       {
00614         
00615         PFClusterRef  ecal1ref = lowEl->clusterRef();
00616         PFClusterRef  ecal2ref = highEl->clusterRef();
00617         assert( !ecal1ref.isNull() );
00618         assert( !ecal2ref.isNull() );
00619         if(debug_)
00620           cout << " PFBlockLink::ECALandECAL" << endl;
00621         dist = testLinkBySuperCluster(ecal1ref,ecal2ref);
00622         break;
00623       }
00624 
00625   case PFBlockLink::ECALandGSF:
00626     {
00627       PFClusterRef  clusterref = lowEl->clusterRef();
00628       assert( !clusterref.isNull() );
00629       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00630       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00631       if ( myTrack->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax ).isValid() )
00632         dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, false, debug_ );
00633       else
00634         dist = -1.;
00635       //   linktest = PFBlock::LINKTEST_RECHIT;
00636       
00637       if ( debug_ ) {
00638         if ( dist > 0. ) {
00639           std::cout << " Here a link has been established" 
00640                     << " between a GSF track an Ecal with dist  " 
00641                     << dist <<  std::endl;
00642         } else {
00643           if(debug_ ) std::cout << " No link found " << std::endl;
00644         }
00645       }
00646       break;
00647     }
00648   case PFBlockLink::TRACKandGSF:
00649     {
00650       PFRecTrackRef trackref = lowEl->trackRefPF();
00651       assert( !trackref.isNull() );
00652       const reco::PFBlockElementGsfTrack *  GsfEl =  
00653         dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00654       GsfPFRecTrackRef gsfref = GsfEl->GsftrackRefPF();
00655       reco::TrackRef kftrackref= (*trackref).trackRef();
00656       assert( !gsfref.isNull() );
00657       PFRecTrackRef refkf = (*gsfref).kfPFRecTrackRef();
00658       if(refkf.isNonnull()) {
00659         reco::TrackRef gsftrackref = (*refkf).trackRef();
00660         if (gsftrackref.isNonnull()&&kftrackref.isNonnull()) {
00661           if (kftrackref == gsftrackref) { 
00662             dist = 0.001;
00663           } else { 
00664             dist = -1.;
00665           }
00666         } else { 
00667           dist = -1.;
00668         }
00669       } else {
00670         dist = -1.;
00671       }
00672       
00673       
00674       if(useConvBremPFRecTracks_) {
00675         if(lowEl->isLinkedToDisplacedVertex()){
00676           vector<PFRecTrackRef> pfrectrack_vec = GsfEl->GsftrackRefPF()->convBremPFRecTrackRef();
00677           if(pfrectrack_vec.size() > 0){
00678             for(unsigned int iconv = 0; iconv <  pfrectrack_vec.size(); ++iconv) {
00679               if( lowEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)) {
00680                 // use track ref
00681                 if(kftrackref == (*pfrectrack_vec[iconv]).trackRef()) {         
00682                   dist = 0.001;
00683                 }
00684               } 
00685               else{
00686                 // use the track base ref
00687                 reco::TrackBaseRef newTrackBaseRef((*pfrectrack_vec[iconv]).trackRef());
00688                 reco::TrackBaseRef elemTrackBaseRef(kftrackref);              
00689                 if(newTrackBaseRef == elemTrackBaseRef){
00690                   dist = 0.001;
00691                 } 
00692               }
00693             }
00694           }
00695         }
00696       }
00697  
00698 
00699       break;      
00700     }
00701          
00702   case PFBlockLink::GSFandBREM:
00703     {
00704       const reco::PFBlockElementGsfTrack * GsfEl  =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(lowEl);
00705       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00706       GsfPFRecTrackRef gsfref = GsfEl->GsftrackRefPF();
00707       GsfPFRecTrackRef bremref = BremEl->GsftrackRefPF();
00708       assert( !gsfref.isNull() );
00709       assert( !bremref.isNull() );
00710       if (gsfref == bremref)  { 
00711         dist = 0.001;
00712       } else { 
00713         dist = -1.;
00714       }
00715       break;
00716     }
00717   case PFBlockLink::GSFandGSF:
00718     {
00719       const reco::PFBlockElementGsfTrack * lowGsfEl  =  
00720         dynamic_cast<const reco::PFBlockElementGsfTrack*>(lowEl);
00721       const reco::PFBlockElementGsfTrack * highGsfEl  =  
00722         dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00723       
00724       GsfPFRecTrackRef lowgsfref = lowGsfEl->GsftrackRefPF();
00725       GsfPFRecTrackRef highgsfref = highGsfEl->GsftrackRefPF();
00726       assert( !lowgsfref.isNull() );
00727       assert( !highgsfref.isNull() );
00728       
00729       if( (lowGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) == false && 
00730            highGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)) ||
00731           (highGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) == false && 
00732            lowGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV))) {
00733         if(lowgsfref->trackId() == highgsfref->trackId()) {
00734           dist = 0.001;
00735         }
00736         else {
00737           dist = -1.;
00738         }
00739       }
00740       break;
00741     }
00742   case PFBlockLink::ECALandBREM:
00743     {
00744       PFClusterRef  clusterref = lowEl->clusterRef();
00745       assert( !clusterref.isNull() );
00746       const reco::PFBlockElementBrem * BremEl =  
00747         dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00748       const PFRecTrack * myTrack = &(BremEl->trackPF());
00749       /*
00750       double DP = (BremEl->DeltaP())*(-1.);
00751       double SigmaDP = BremEl->SigmaDeltaP();
00752       double SignBremDp = DP/SigmaDP;
00753       */
00754       bool isBrem = true;
00755       if ( myTrack->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax ).isValid() )
00756         dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, isBrem, debug_);
00757       else
00758         dist = -1.;
00759       if( debug_ && dist > 0. ) 
00760         std::cout << "ECALandBREM: dist testTrackAndClusterByRecHit " 
00761                   << dist << std::endl;
00762       //   linktest = PFBlock::LINKTEST_RECHIT;
00763       break;
00764     }
00765   case PFBlockLink::HCALandGSF:
00766     {
00767       PFClusterRef  clusterref = lowEl->clusterRef();
00768       assert( !clusterref.isNull() );
00769       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00770       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00771       if ( myTrack->extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance ).isValid() )
00772         dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, false, debug_ );
00773       else
00774         dist = -1.;
00775 
00776       //    linktest = PFBlock::LINKTEST_RECHIT;
00777       break;
00778     }
00779   case PFBlockLink::HCALandBREM:
00780     {
00781       PFClusterRef  clusterref = lowEl->clusterRef();
00782       assert( !clusterref.isNull() );
00783       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00784       const PFRecTrack * myTrack = &(BremEl->trackPF());
00785       bool isBrem = true;
00786       if ( myTrack->extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance ).isValid() )
00787         dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, isBrem, debug_);
00788       else
00789         dist = -1.;
00790       break;
00791     }
00792   case PFBlockLink::SCandECAL:
00793     {
00794       PFClusterRef  clusterref = lowEl->clusterRef();
00795 
00796       assert( !clusterref.isNull() );
00797       
00798       const reco::PFBlockElementSuperCluster * scEl = 
00799         dynamic_cast<const reco::PFBlockElementSuperCluster*>(highEl);
00800       assert (!scEl->superClusterRef().isNull());
00801       dist = testSuperClusterPFCluster(scEl->superClusterRef(),
00802                                        clusterref);
00803       break;
00804     }
00805     /*
00806   // Links between the two preshower layers are not used for now - disable
00807   case PFBlockLink::PS1andPS2:
00808     {
00809       PFClusterRef  ps1ref = lowEl->clusterRef();
00810       PFClusterRef  ps2ref = highEl->clusterRef();
00811       assert( !ps1ref.isNull() );
00812       assert( !ps2ref.isNull() );
00813       // PJ - 14-May-09 : A link by rechit is needed here !
00814       // dist = testPS1AndPS2( *ps1ref, *ps2ref );
00815       dist = -1.;
00816       linktest = PFBlock::LINKTEST_RECHIT;
00817       break;
00818     }
00819   case PFBlockLink::TRACKandPS1:
00820   case PFBlockLink::TRACKandPS2:
00821     {
00822       //cout<<"TRACKandPS"<<endl;
00823       PFRecTrackRef trackref = lowEl->trackRefPF();
00824       PFClusterRef  clusterref = highEl->clusterRef();
00825       assert( !trackref.isNull() );
00826       assert( !clusterref.isNull() );
00827       // PJ - 14-May-09 : A link by rechit is needed here !
00828       dist = testTrackAndPS( *trackref, *clusterref );
00829       linktest = PFBlock::LINKTEST_RECHIT;
00830       break;
00831     }
00832     // GSF Track/Brem Track and preshower cluster links are not used for now - disable
00833   case PFBlockLink::PS1andGSF:
00834   case PFBlockLink::PS2andGSF:
00835     {
00836       PFClusterRef  psref = lowEl->clusterRef();
00837       assert( !psref.isNull() );
00838       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00839       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00840       // PJ - 14-May-09 : A link by rechit is needed here !
00841       dist = testTrackAndPS( *myTrack, *psref );
00842       linktest = PFBlock::LINKTEST_RECHIT;
00843       break;
00844     }
00845   case PFBlockLink::PS1andBREM:
00846   case PFBlockLink::PS2andBREM:
00847     {
00848       PFClusterRef  psref = lowEl->clusterRef();
00849       assert( !psref.isNull() );
00850       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00851       const PFRecTrack * myTrack = &(BremEl->trackPF());
00852       // PJ - 14-May-09 : A link by rechit is needed here !
00853       dist = testTrackAndPS( *myTrack, *psref );
00854       linktest = PFBlock::LINKTEST_RECHIT;
00855       break;
00856     }
00857     */
00858 
00859   default:
00860     dist = -1.;
00861     //   linktest = PFBlock::LINKTEST_RECHIT;
00862     // cerr<<"link type not implemented yet: 0x"<<hex<<linktype<<dec<<endl;
00863     // assert(0);
00864     return;
00865   }
00866 }
00867 
00868 double
00869 PFBlockAlgo::testTrackAndPS(const PFRecTrack& track, 
00870                             const PFCluster& ps)  const {
00871 
00872 #ifdef PFLOW_DEBUG
00873   //   cout<<"entering testTrackAndPS"<<endl;
00874   // resolution of PS cluster dxdx and dydy from strip pitch and length
00875   double dx=0.;
00876   double dy=0.;
00877   
00878   unsigned layerid =0;
00879   // PS1: vertical strips  PS2: horizontal strips
00880   switch (ps.layer()) {
00881   case PFLayer::PS1:
00882     layerid = reco::PFTrajectoryPoint::PS1;
00883     
00884     // vertical strips in PS1, measure x with pitch precision
00885     dx = resPSpitch_;
00886     dy = resPSlength_; 
00887     break;
00888   case PFLayer::PS2:
00889     layerid = reco::PFTrajectoryPoint::PS2;
00890     // horizontal strips in PS2, measure y with pitch precision
00891     dy = resPSpitch_;
00892     dx = resPSlength_;
00893     break;
00894   default:
00895     break;
00896   }
00897   const reco::PFTrajectoryPoint& atPS
00898     = track.extrapolatedPoint( layerid );  
00899   // did not reach PS, cannot be associated with a cluster.
00900   if( ! atPS.isValid() ) return -1.;   
00901   
00902   double trackx = atPS.position().X();
00903   double tracky = atPS.position().Y();
00904   double trackz = atPS.position().Z(); // MDN jan 09
00905   
00906   // ps position  x, y
00907   double psx = ps.position().X();
00908   double psy = ps.position().Y();
00909   // MDN Jan 09: check that trackz and psz have the same sign
00910   double psz = ps.position().Z();
00911   if( trackz*psz < 0.) return -1.; 
00912   
00913   // double chi2 = (psx-trackx)*(psx-trackx)/(dx*dx + trackresolx*trackresolx)
00914   //  + (psy-tracky)*(psy-tracky)/(dy*dy + trackresoly*trackresoly);
00915 
00916   double dist = std::sqrt( (psx-trackx)*(psx-trackx)
00917                          + (psy-tracky)*(psy-tracky));  
00918   if(debug_) cout<<"testTrackAndPS "<< dist <<" "<<endl;
00919   if(debug_){
00920     cout<<" trackx " << trackx 
00921         <<" tracky " << tracky 
00922         <<" psx "    <<  psx   
00923         <<" psy "    << psy    
00924         << endl;
00925   }
00926 #endif
00927   
00928   // Return -1. as long as no link by rechit is available
00929   return -1.;
00930 }
00931 
00932 double
00933 PFBlockAlgo::testECALAndHCAL(const PFCluster& ecal, 
00934                              const PFCluster& hcal)  const {
00935   
00936   //   cout<<"entering testECALAndHCAL"<<endl;
00937   
00938   double dist = fabs(ecal.positionREP().Eta()) > 2.5 ?
00939     LinkByRecHit::computeDist( ecal.positionREP().Eta(),
00940                                ecal.positionREP().Phi(), 
00941                                hcal.positionREP().Eta(), 
00942                                hcal.positionREP().Phi() )
00943     : 
00944     -1.;
00945 
00946 #ifdef PFLOW_DEBUG
00947   if(debug_) cout<<"testECALAndHCAL "<< dist <<" "<<endl;
00948   if(debug_){
00949     cout<<" ecaleta " << ecal.positionREP().Eta()
00950         <<" ecalphi " << ecal.positionREP().Phi()
00951         <<" hcaleta " << hcal.positionREP().Eta()
00952         <<" hcalphi " << hcal.positionREP().Phi()
00953   }
00954 #endif
00955 
00956   if ( dist < 0.2 ) return dist; 
00957  
00958   // Need to implement a link by RecHit
00959   return -1.;
00960 }
00961 
00962 double
00963 PFBlockAlgo::testHCALAndHO(const PFCluster& hcal, 
00964                              const PFCluster& ho)  const {
00965   
00966   double dist = fabs(hcal.positionREP().Eta()) < 1.5 ?
00967     LinkByRecHit::computeDist( hcal.positionREP().Eta(),
00968                                hcal.positionREP().Phi(), 
00969                                ho.positionREP().Eta(), 
00970                                ho.positionREP().Phi() )
00971     : 
00972     -1.;
00973 
00974 #ifdef PFLOW_DEBUG
00975   if(debug_) cout<<"testHCALAndHO "<< dist <<" "<<endl;
00976   if(debug_){
00977     cout<<" hcaleta " << hcal.positionREP().Eta()
00978         <<" hcalphi " << hcal.positionREP().Phi()
00979         <<" hoeta " << ho.positionREP().Eta()
00980         <<" hophi " << ho.positionREP().Phi()
00981         <<" dist " << dist<<endl;
00982   }
00983 #endif
00984 
00985   if ( dist < 0.20 ) return dist; 
00986  
00987   // Need to implement a link by RecHit
00988   return -1.;
00989 }
00990 
00991 
00992 
00993 double
00994 PFBlockAlgo::testLinkBySuperCluster(const PFClusterRef& ecal1, 
00995                                     const PFClusterRef& ecal2)  const {
00996   
00997   //  cout<<"entering testECALAndECAL "<< pfcRefSCMap_.size() << endl;
00998   
00999   double dist = -1;
01000   
01001   // the first one is not in any super cluster
01002   int testindex=pfcSCVec_[ecal1.key()];
01003   if(testindex == -1.) return dist;
01004   //  if(itcheck==pfcRefSCMap_.end()) return dist;
01005   // now retrieve the of PFclusters in this super cluster  
01006 
01007   const std::vector<reco::PFClusterRef> & thePFClusters(scpfcRefs_[testindex]);
01008   
01009   unsigned npf=thePFClusters.size();
01010   for(unsigned i=0;i<npf;++i)
01011     {
01012       if(thePFClusters[i]==ecal2) // yes they are in the same SC 
01013         {
01014           dist=LinkByRecHit::computeDist( ecal1->positionREP().Eta(),
01015                                           ecal1->positionREP().Phi(), 
01016                                           ecal2->positionREP().Eta(), 
01017                                           ecal2->positionREP().Phi() );
01018 //        std::cout << " DETA " << fabs(ecal1->positionREP().Eta()-ecal2->positionREP().Eta()) << std::endl;
01019 //        if(fabs(ecal1->positionREP().Eta()-ecal2->positionREP().Eta())>0.2)
01020 //          {
01021 //            std::cout <<  " Super Cluster " <<  *(superClusters_[testindex]) << std::endl;
01022 //            std::cout <<  " Cluster1 " <<  *ecal1 << std::endl;
01023 //            std::cout <<  " Cluster2 " <<  *ecal2 << std::endl;
01024 //            ClusterClusterMapping::checkOverlap(*ecal1,superClusters_,0.01,true);
01025 //            ClusterClusterMapping::checkOverlap(*ecal2,superClusters_,0.01,true);
01026 //          }
01027           return dist;
01028         }
01029     }
01030   return dist;
01031 }
01032 
01033 
01034 double
01035 PFBlockAlgo::testSuperClusterPFCluster(const SuperClusterRef& ecal1, 
01036                                        const PFClusterRef& ecal2)  const {
01037   
01038   //  cout<<"entering testECALAndECAL "<< pfcRefSCMap_.size() << endl;
01039   
01040   double dist = -1;
01041   
01042   bool overlap=ClusterClusterMapping::overlap(*ecal1,*ecal2);
01043   
01044   if(overlap)   {
01045     dist=LinkByRecHit::computeDist( ecal1->position().eta(),
01046                                     ecal1->position().phi(), 
01047                                     ecal2->positionREP().Eta(), 
01048                                     ecal2->positionREP().Phi() );
01049     return dist;
01050   }
01051   return dist;
01052 }
01053 
01054 
01055 
01056 double
01057 PFBlockAlgo::testPS1AndPS2(const PFCluster& ps1, 
01058                            const PFCluster& ps2)  const {
01059   
01060 #ifdef PFLOW_DEBUG
01061   //   cout<<"entering testPS1AndPS2"<<endl;
01062   
01063   // compute chi2 in y, z using swimming formulae
01064   // y2 = y1 * z2/z1   and x2 = x1 *z2/z1
01065   
01066   // ps position1  x, y, z
01067   double x1 = ps1.position().X();
01068   double y1 = ps1.position().Y();
01069   double z1 = ps1.position().Z();
01070   double x2 = ps2.position().X();
01071   double y2 = ps2.position().Y();
01072   double z2 = ps2.position().Z();
01073   // MDN Bug correction Jan 09: check that z1 and z2 have the same sign!
01074   if (z1*z2<0.) -1.;
01075   // swim to PS2
01076   double scale = z2/z1;
01077   double x1atPS2 = x1*scale;
01078   double y1atPS2 = y1*scale;
01079   // resolution of PS cluster dxdx and dydy from strip pitch and length
01080   // vertical strips in PS1, measure x with pitch precision
01081   double dx1dx1 = resPSpitch_*resPSpitch_*scale*scale;
01082   double dy1dy1 = resPSlength_*resPSlength_*scale*scale;
01083   // horizontal strips in PS2 , measure y with pitch precision
01084   double dy2dy2 = resPSpitch_*resPSpitch_;
01085   double dx2dx2 = resPSlength_*resPSlength_;
01086   
01087   // double chi2 = (x2-x1atPS2)*(x2-x1atPS2)/(dx1dx1 + dx2dx2) 
01088   //  + (y2-y1atPS2)*(y2-y1atPS2)/(dy1dy1 + dy2dy2);
01089   
01090   double dist = std::sqrt( (x2-x1atPS2)*(x2-x1atPS2)
01091                          + (y2-y1atPS2)*(y2-y1atPS2));
01092     
01093   if(debug_) cout<<"testPS1AndPS2 "<<dist<<" "<<endl;
01094   if(debug_){
01095     cout<<" x1atPS2 "<< x1atPS2 << " dx1 "<<resPSpitch_*scale
01096         <<" y1atPS2 "<< y1atPS2 << " dy1 "<<resPSlength_*scale<< endl
01097         <<" x2 " <<x2  << " dx2 "<<resPSlength_
01098         <<" y2 " << y2 << " dy2 "<<resPSpitch_<< endl;
01099   }
01100 #endif
01101 
01102   // Need a link by rechit here
01103   return -1.; 
01104 }
01105 
01106 
01107 
01108 double
01109 PFBlockAlgo::testLinkByVertex( const reco::PFBlockElement* elt1, 
01110                                const reco::PFBlockElement* elt2) const {
01111 
01112   //  cout << "Test link by vertex between" << endl << *elt1 << endl << " and " << endl << *elt2 << endl;
01113 
01114   double result=-1.;
01115 
01116   reco::PFBlockElement::TrackType T_TO_DISP = reco::PFBlockElement::T_TO_DISP;
01117   reco::PFBlockElement::TrackType T_FROM_DISP = reco::PFBlockElement::T_FROM_DISP;
01118   PFDisplacedTrackerVertexRef ni1_TO_DISP = elt1->displacedVertexRef(T_TO_DISP);
01119   PFDisplacedTrackerVertexRef ni2_TO_DISP = elt2->displacedVertexRef(T_TO_DISP);
01120   PFDisplacedTrackerVertexRef ni1_FROM_DISP = elt1->displacedVertexRef(T_FROM_DISP);
01121   PFDisplacedTrackerVertexRef ni2_FROM_DISP = elt2->displacedVertexRef(T_FROM_DISP);
01122   
01123   if( ni1_TO_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
01124     if( ni1_TO_DISP == ni2_FROM_DISP ) { result = 1.0; return result; }
01125 
01126   if( ni1_FROM_DISP.isNonnull() && ni2_TO_DISP.isNonnull())
01127     if( ni1_FROM_DISP == ni2_TO_DISP ) { result = 1.0; return result; }
01128 
01129   if( ni1_FROM_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
01130     if( ni1_FROM_DISP == ni2_FROM_DISP ) { result = 1.0; return result; }
01131     
01132   
01133   if (  elt1->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)  &&
01134              elt2->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)  ) {
01135     
01136     if(debug_ ) std::cout << " testLinkByVertex On Conversions " << std::endl;
01137     
01138     if ( elt1->convRef().isNonnull() && elt2->convRef().isNonnull() ) {
01139       if(debug_ ) std::cout << " PFBlockAlgo.cc testLinkByVertex  Cconversion Refs are non null  " << std::endl;      
01140       if ( elt1->convRef() ==  elt2->convRef() ) {
01141         result=1.0;
01142         if(debug_ ) std::cout << " testLinkByVertex  Cconversion Refs are equal  " << std::endl;      
01143         return result;
01144       }
01145     } 
01146     
01147   }
01148   
01149   if (  elt1->trackType(reco::PFBlockElement::T_FROM_V0)  &&
01150              elt2->trackType(reco::PFBlockElement::T_FROM_V0)  ) {
01151     if(debug_ ) std::cout << " testLinkByVertex On V0 " << std::endl;
01152     if ( elt1->V0Ref().isNonnull() && elt2->V0Ref().isNonnull() ) {
01153       if(debug_ ) std::cout << " PFBlockAlgo.cc testLinkByVertex  V0 Refs are non null  " << std::endl;
01154       if ( elt1->V0Ref() ==  elt2->V0Ref() ) {
01155         result=1.0;
01156         if(debug_ ) std::cout << " testLinkByVertex  V0 Refs are equal  " << std::endl;
01157         return result;
01158       }
01159     }
01160   }
01161 
01162   return result;
01163 }
01164 
01165 
01166 
01167 void 
01168 PFBlockAlgo::checkMaskSize( const reco::PFRecTrackCollection& tracks,
01169                             const reco::GsfPFRecTrackCollection& gsftracks, 
01170                             const reco::PFClusterCollection&  ecals,
01171                             const reco::PFClusterCollection&  hcals,
01172                             const reco::PFClusterCollection&  hos,
01173                             const reco::PFClusterCollection&  hfems,
01174                             const reco::PFClusterCollection&  hfhads,
01175                             const reco::PFClusterCollection&  pss, 
01176                             const reco::PhotonCollection&  egphh, 
01177                             const Mask& trackMask, 
01178                             const Mask& gsftrackMask,  
01179                             const Mask& ecalMask,
01180                             const Mask& hcalMask,
01181                             const Mask& hoMask, 
01182                             const Mask& hfemMask,
01183                             const Mask& hfhadMask,                    
01184                             const Mask& psMask,
01185                             const Mask& phMask) const {
01186 
01187   if( !trackMask.empty() && 
01188       trackMask.size() != tracks.size() ) {
01189     string err = "PFBlockAlgo::setInput: ";
01190     err += "The size of the track mask is different ";
01191     err += "from the size of the track vector.";
01192     throw std::length_error( err.c_str() );
01193   }
01194 
01195   if( !gsftrackMask.empty() && 
01196       gsftrackMask.size() != gsftracks.size() ) {
01197     string err = "PFBlockAlgo::setInput: ";
01198     err += "The size of the gsf track mask is different ";
01199     err += "from the size of the gsftrack vector.";
01200     throw std::length_error( err.c_str() );
01201   }
01202 
01203   if( !ecalMask.empty() && 
01204       ecalMask.size() != ecals.size() ) {
01205     string err = "PFBlockAlgo::setInput: ";
01206     err += "The size of the ecal mask is different ";
01207     err += "from the size of the ecal clusters vector.";
01208     throw std::length_error( err.c_str() );
01209   }
01210   
01211   if( !hcalMask.empty() && 
01212       hcalMask.size() != hcals.size() ) {
01213     string err = "PFBlockAlgo::setInput: ";
01214     err += "The size of the hcal mask is different ";
01215     err += "from the size of the hcal clusters vector.";
01216     throw std::length_error( err.c_str() );
01217   }
01218 
01219   if( !hoMask.empty() && 
01220       hoMask.size() != hos.size() ) {
01221     string err = "PFBlockAlgo::setInput: ";
01222     err += "The size of the ho mask is different ";
01223     err += "from the size of the ho clusters vector.";
01224     throw std::length_error( err.c_str() );
01225   }
01226 
01227 
01228   if( !hfemMask.empty() && 
01229       hfemMask.size() != hfems.size() ) {
01230     string err = "PFBlockAlgo::setInput: ";
01231     err += "The size of the hfem mask is different ";
01232     err += "from the size of the hfem clusters vector.";
01233     throw std::length_error( err.c_str() );
01234   }
01235 
01236   if( !hfhadMask.empty() && 
01237       hfhadMask.size() != hfhads.size() ) {
01238     string err = "PFBlockAlgo::setInput: ";
01239     err += "The size of the hfhad mask is different ";
01240     err += "from the size of the hfhad clusters vector.";
01241     throw std::length_error( err.c_str() );
01242   }
01243 
01244   if( !psMask.empty() && 
01245       psMask.size() != pss.size() ) {
01246     string err = "PFBlockAlgo::setInput: ";
01247     err += "The size of the ps mask is different ";
01248     err += "from the size of the ps clusters vector.";
01249     throw std::length_error( err.c_str() );
01250   }
01251   
01252     if( !phMask.empty() && 
01253       phMask.size() != egphh.size() ) {
01254     string err = "PFBlockAlgo::setInput: ";
01255     err += "The size of the photon mask is different ";
01256     err += "from the size of the photon vector.";
01257     throw std::length_error( err.c_str() );
01258   }
01259 
01260 }
01261 
01262 
01263 std::ostream& operator<<(std::ostream& out, const PFBlockAlgo& a) {
01264   if(! out) return out;
01265   
01266   out<<"====== Particle Flow Block Algorithm ======= ";
01267   out<<endl;
01268   out<<"number of unassociated elements : "<<a.elements_.size()<<endl;
01269   out<<endl;
01270   
01271   for(PFBlockAlgo::IEC ie = a.elements_.begin(); 
01272       ie != a.elements_.end(); ++ie) {
01273     out<<"\t"<<**ie <<endl;
01274   }
01275 
01276   
01277   //   const PFBlockCollection& blocks = a.blocks();
01278 
01279   const std::auto_ptr< reco::PFBlockCollection >& blocks
01280     = a.blocks(); 
01281     
01282   if(!blocks.get() ) {
01283     out<<"blocks already transfered"<<endl;
01284   }
01285   else {
01286     out<<"number of blocks : "<<blocks->size()<<endl;
01287     out<<endl;
01288     
01289     for(PFBlockAlgo::IBC ib=blocks->begin(); 
01290         ib != blocks->end(); ++ib) {
01291       out<<(*ib)<<endl;
01292     }
01293   }
01294 
01295   return out;
01296 }
01297 
01298 bool 
01299 PFBlockAlgo::goodPtResolution( const reco::TrackRef& trackref) {
01300 
01301   double P = trackref->p();
01302   double Pt = trackref->pt();
01303   double DPt = trackref->ptError();
01304   unsigned int NHit = trackref->hitPattern().trackerLayersWithMeasurement();
01305   unsigned int NLostHit = trackref->hitPattern().trackerLayersWithoutMeasurement();
01306   unsigned int LostHits = trackref->numberOfLostHits();
01307   double sigmaHad = sqrt(1.20*1.20/P+0.06*0.06) / (1.+LostHits);
01308 
01309   // iteration 1,2,3,4,5 correspond to algo = 1/4,5,6,7,8,9
01310   unsigned int Algo = 0; 
01311   switch (trackref->algo()) {
01312   case TrackBase::ctf:
01313   case TrackBase::iter0:
01314   case TrackBase::iter1:
01315   case TrackBase::iter2:
01316     Algo = 0;
01317     break;
01318   case TrackBase::iter3:
01319     Algo = 1;
01320     break;
01321   case TrackBase::iter4:
01322     Algo = 2;
01323     break;
01324   case TrackBase::iter5:
01325     Algo = 3;
01326     break;
01327   case TrackBase::iter6:
01328     Algo = 4;
01329     break;
01330   default:
01331     Algo = useIterTracking_ ? 5 : 0;
01332     break;
01333   }
01334 
01335   // Protection against 0 momentum tracks
01336   if ( P < 0.05 ) return false;
01337 
01338   // Temporary : Reject all tracking iteration beyond 5th step. 
01339   if ( Algo > 4 ) return false;
01340  
01341   if (debug_) cout << " PFBlockAlgo: PFrecTrack->Track Pt= "
01342                    << Pt << " DPt = " << DPt << endl;
01343   if ( ( DPtovPtCut_[Algo] > 0. && 
01344          DPt/Pt > DPtovPtCut_[Algo]*sigmaHad ) || 
01345        NHit < NHitCut_[Algo] ) { 
01346     // (Algo >= 3 && LostHits != 0) ) {
01347     if (debug_) cout << " PFBlockAlgo: skip badly measured track"
01348                      << ", P = " << P 
01349                      << ", Pt = " << Pt 
01350                      << " DPt = " << DPt 
01351                      << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01352                      << ", Algo = " << Algo
01353                      << endl;
01354     if (debug_) cout << " cut is DPt/Pt < " << DPtovPtCut_[Algo] * sigmaHad << endl;
01355     if (debug_) cout << " cut is NHit >= " << NHitCut_[Algo] << endl;
01356     /*
01357     std::cout << "Track REJECTED : ";
01358     std::cout << ", P = " << P 
01359               << ", Pt = " << Pt 
01360               << " DPt = " << DPt 
01361               << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01362               << ", Algo = " << Algo
01363               << std::endl;
01364     */
01365     return false;
01366   }
01367 
01368   /*
01369   std::cout << "Track Accepted : ";
01370   std::cout << ", P = " << P 
01371        << ", Pt = " << Pt 
01372        << " DPt = " << DPt 
01373        << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01374        << ", Algo = " << Algo
01375        << std::endl;
01376   */
01377   return true;
01378 }
01379 
01380 int
01381 PFBlockAlgo::muAssocToTrack( const reco::TrackRef& trackref,
01382                              const edm::Handle<reco::MuonCollection>& muonh) const {
01383   if(muonh.isValid() ) {
01384     for(unsigned j=0;j<muonh->size(); ++j) {
01385       reco::MuonRef muonref( muonh, j );
01386       if (muonref->track().isNonnull()) 
01387         if( muonref->track() == trackref ) return j;
01388     }
01389   }
01390   return -1; // not found
01391 }
01392 
01393 int 
01394 PFBlockAlgo::muAssocToTrack( const reco::TrackRef& trackref,
01395                              const edm::OrphanHandle<reco::MuonCollection>& muonh) const {
01396   if(muonh.isValid() ) {
01397     for(unsigned j=0;j<muonh->size(); ++j) {
01398       reco::MuonRef muonref( muonh, j );
01399       if (muonref->track().isNonnull())
01400         if( muonref->track() == trackref ) return j;
01401     }
01402   }
01403   return -1; // not found
01404 }
01405 
01406 
01407 // This prefilter avoid to call associate when not necessary.
01408 // ACHTUNG!!!! If you introduce new links check that they are not desables here
01409 inline bool
01410 PFBlockAlgo::linkPrefilter(const reco::PFBlockElement* last, const reco::PFBlockElement* next) const {
01411 
01412   PFBlockElement::Type type1 = (last)->type();
01413   PFBlockElement::Type type2 = (next)->type();
01414 
01415   if( type1==type2 ) {
01416     // cannot link 2 elements of the same type. 
01417     // except if the elements are 2 tracks or 2 ECAL
01418     if( type1!=PFBlockElement::TRACK && type1!=PFBlockElement::GSF &&
01419         type1!=PFBlockElement::ECAL) { // && type1!=PFBlockElement::HCAL) {
01420       return false;
01421     }
01422 
01423     if (type1==PFBlockElement::ECAL && bNoSuperclus_) return false;
01424 
01425     // cannot link two primary tracks  (except if they come from a V0)
01426     if( type1 ==PFBlockElement::TRACK) {
01427       if ( !((last)->isLinkedToDisplacedVertex()) || !((next)->isLinkedToDisplacedVertex())) 
01428       return false;
01429     }
01430   }
01431 
01432   if ((type1 == PFBlockElement::PS1 || type1 == PFBlockElement::PS2) && (type2 != PFBlockElement::ECAL)) return false;
01433   if ((type2 == PFBlockElement::PS1 || type2 == PFBlockElement::PS2) && (type1 != PFBlockElement::ECAL)) return false;
01434   if ((type1 == PFBlockElement::HFEM && type2 != PFBlockElement::HFHAD) || (type1 == PFBlockElement::HFHAD && type2 != PFBlockElement::HFEM)) return false;
01435 
01436   if (useKDTreeTrackEcalLinker_){ 
01437   
01438     if ( type1 == PFBlockElement::TRACK && type2 == PFBlockElement::ECAL)
01439       if ( last->isMultilinksValide()  && last->getMultilinks().size()==0 ) return false;
01440     if ( type2 == PFBlockElement::TRACK && type1 == PFBlockElement::ECAL)
01441       if ( next->isMultilinksValide() && next->getMultilinks().size()==0 ) return false;
01442     if ( type1 == PFBlockElement::PS1 || type1 == PFBlockElement::PS2)
01443       if ( last->isMultilinksValide()  && last->getMultilinks().size()==0 ) return false;
01444     if ( type2 == PFBlockElement::PS1 || type2 == PFBlockElement::PS2)
01445       if ( next->isMultilinksValide() && next->getMultilinks().size()==0 ) return false;
01446 
01447   }
01448 
01449   return true;
01450 
01451 }
01452