CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/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   PFBlockElement::Type lowType = type1;
00341   PFBlockElement::Type highType = type2;
00342   const PFBlockElement* lowEl = el1;
00343   const PFBlockElement* highEl = el2;
00344   
00345   if(type1>type2) {
00346     lowType = type2;
00347     highType = type1;
00348     lowEl = el2;
00349     highEl = el1;
00350   }
00351   
00352   switch(linktype) {
00353   // Track and preshower cluster links are not used for now - disable
00354   case PFBlockLink::PS1andECAL:
00355   case PFBlockLink::PS2andECAL:
00356     {
00357       // if(debug_ ) cout<< "PSandECAL" <<endl;
00358       PFClusterRef  psref = lowEl->clusterRef();
00359       PFClusterRef  ecalref = highEl->clusterRef();
00360       assert( !psref.isNull() );
00361       assert( !ecalref.isNull() );
00362 
00363       // Check if the linking has been done using the KDTree algo
00364       // Glowinski & Gouzevitch
00365       if ( useKDTreeTrackEcalLinker_ && lowEl->isMultilinksValide() ) { // KDTree algo
00366         const reco::PFMultilinksType& multilinks = lowEl->getMultilinks();
00367         
00368         double ecalPhi = ecalref->positionREP().Phi();
00369         double ecalEta = ecalref->positionREP().Eta();
00370         
00371         // Check if the link PS/Ecal exist
00372         reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
00373         for (; mlit != multilinks.end(); ++mlit)
00374           if ((mlit->first == ecalPhi) && (mlit->second == ecalEta))
00375             break;
00376         
00377         // If the link exist, we fill dist and linktest. We use old algorithme method.
00378         if (mlit != multilinks.end()){
00379           double xPS = psref->position().X();
00380           double yPS = psref->position().Y();
00381           double xECAL  = ecalref->position().X();
00382           double yECAL  = ecalref->position().Y();
00383 
00384           dist = LinkByRecHit::computeDist(xECAL/1000.,yECAL/1000.,xPS/1000.  ,yPS/1000, false);
00385         }
00386 
00387       } else { //Old algorithm
00388         dist = LinkByRecHit::testECALAndPSByRecHit( *ecalref, *psref ,debug_);
00389       }
00390 
00391       //      linktest = PFBlock::LINKTEST_RECHIT;
00392       
00393       break;
00394     }
00395   case PFBlockLink::TRACKandECAL:
00396     {
00397       if(debug_ ) cout<<"TRACKandECAL"<<endl;
00398 
00399       PFRecTrackRef trackref = lowEl->trackRefPF();
00400       PFClusterRef  clusterref = highEl->clusterRef();
00401       assert( !trackref.isNull() );
00402       assert( !clusterref.isNull() );
00403 
00404       if(debug_ ) std::cout << " Track pt " << trackref->trackRef()->pt() << std::endl;
00405 
00406       // Check if the linking has been done using the KDTree algo
00407       // Glowinski & Gouzevitch
00408       if ( useKDTreeTrackEcalLinker_ && lowEl->isMultilinksValide() ) { //KDTree Algo
00409 
00410         const reco::PFMultilinksType& multilinks = lowEl->getMultilinks();
00411         double ecalphi = clusterref->positionREP().Phi();
00412         double ecaleta = clusterref->positionREP().Eta();
00413 
00414         // Check if the link Track/Ecal exist
00415         reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
00416         for (; mlit != multilinks.end(); ++mlit)
00417           if ((mlit->first == ecalphi) && (mlit->second == ecaleta))
00418             break;
00419 
00420         
00421         // If the link exist, we fill dist and linktest. We use old algorithme method.
00422         if (mlit != multilinks.end()){
00423 
00424 
00425           //Should be something like this :
00426           //      const reco::PFRecTrack& track = *trackref;
00427           //instead of this :
00428           /*
00429           reco::PFRecTrack track (*trackref);
00430           const reco::PFTrajectoryPoint& atECAL_tmp = 
00431             (*trackref).extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax );
00432           if(std::abs(atECAL_tmp.positionREP().Eta())<1E-9 &&
00433              std::abs(atECAL_tmp.positionREP().Phi())<1E-9 &&
00434              atECAL_tmp.positionREP().R()<1E-9) 
00435             track.calculatePositionREP();
00436           */
00437 
00438           const reco::PFTrajectoryPoint& atECAL = 
00439             trackref->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax );
00440           
00441           double tracketa = atECAL.positionREP().Eta();
00442           double trackphi = atECAL.positionREP().Phi();
00443 
00444           dist = LinkByRecHit::computeDist(ecaleta, ecalphi, tracketa, trackphi);
00445         }
00446 
00447       } else {// Old algorithm
00448         if ( trackref->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax ).isValid() )
00449           dist = LinkByRecHit::testTrackAndClusterByRecHit( *trackref, *clusterref, false, debug_ );
00450         else
00451           dist = -1.;
00452       }
00453 
00454       if ( debug_ ) { 
00455         if( dist > 0. ) { 
00456           std::cout << " Here a link has been established"
00457                     << " between a track an Ecal with dist  " 
00458                     << dist <<  std::endl;
00459         } else
00460           std::cout << " No link found " << std::endl;
00461       }
00462       
00463       //     linktest = PFBlock::LINKTEST_RECHIT;
00464       break;
00465     }
00466   case PFBlockLink::TRACKandHCAL:
00467     {
00468       //      if(debug_ ) cout<<"TRACKandHCAL"<<endl;
00469 
00470       PFRecTrackRef trackref = lowEl->trackRefPF();
00471       PFClusterRef  clusterref = highEl->clusterRef();
00472       assert( !trackref.isNull() );
00473       assert( !clusterref.isNull() );
00474 
00475       // Check if the linking has been done using the KDTree algo
00476       // Glowinski & Gouzevitch
00477       if ( useKDTreeTrackEcalLinker_ && highEl->isMultilinksValide() ) { //KDTree Algo
00478         
00479         const reco::PFMultilinksType& multilinks = highEl->getMultilinks();
00480 
00481         /*
00482         reco::PFRecTrack track (*trackref);
00483         const reco::PFTrajectoryPoint& atHCALEntrance_tmp = 
00484           (*trackref).extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance);
00485         if (std::abs(atHCALEntrance_tmp.positionREP().Eta())<1E-9 &&
00486             std::abs(atHCALEntrance_tmp.positionREP().Phi())<1E-9 &&
00487             atHCALEntrance_tmp.positionREP().R()<1E-9)   
00488           track.calculatePositionREP();
00489         */
00490 
00491         const reco::PFTrajectoryPoint& atHCAL = 
00492           trackref->extrapolatedPoint(reco::PFTrajectoryPoint::HCALEntrance);
00493         
00494           
00495         double tracketa = atHCAL.positionREP().Eta();
00496         double trackphi = atHCAL.positionREP().Phi();
00497 
00498         // Check if the link Track/Ecal exist
00499         reco::PFMultilinksType::const_iterator mlit = multilinks.begin();
00500         for (; mlit != multilinks.end(); ++mlit)
00501           if ((mlit->first == trackphi) && (mlit->second == tracketa))
00502             break;
00503 
00504         // If the link exist, we fill dist and linktest. We use old algorithme method.
00505         if (mlit != multilinks.end()){
00506 
00507           const reco::PFTrajectoryPoint& atHCALExit = 
00508             trackref->extrapolatedPoint(reco::PFTrajectoryPoint::HCALExit);
00509           double dHEta = atHCALExit.positionREP().Eta()-atHCAL.positionREP().Eta();
00510           double dHPhi = atHCALExit.positionREP().Phi()-atHCAL.positionREP().Phi(); 
00511           if ( dHPhi > M_PI ) dHPhi = dHPhi - 2.*M_PI;
00512           else if ( dHPhi < -M_PI ) dHPhi = dHPhi + 2.*M_PI; 
00513           tracketa += 0.1 * dHEta;
00514           trackphi += 0.1 * dHPhi;
00515 
00516           double clusterphi = clusterref->positionREP().Phi();
00517           double clustereta = clusterref->positionREP().Eta();
00518           
00519           dist = LinkByRecHit::computeDist(clustereta, clusterphi, tracketa, trackphi);
00520         }       
00521 
00522       } else {// Old algorithm
00523         if ( trackref->extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance ).isValid() )
00524           dist = LinkByRecHit::testTrackAndClusterByRecHit( *trackref, *clusterref, false, debug_ );
00525         else
00526           dist = -1.;
00527       }
00528 
00529       //      linktest = PFBlock::LINKTEST_RECHIT;
00530       break;
00531     }
00532   case PFBlockLink::ECALandHCAL:
00533     {
00534       //       cout<<"ECALandHCAL"<<endl;
00535       PFClusterRef  ecalref = lowEl->clusterRef();
00536       PFClusterRef  hcalref = highEl->clusterRef();
00537       assert( !ecalref.isNull() );
00538       assert( !hcalref.isNull() );
00539       // PJ - 14-May-09 : A link by rechit is needed here !
00540       dist = testECALAndHCAL( *ecalref, *hcalref );
00541       // dist = -1.;
00542       //     linktest = PFBlock::LINKTEST_RECHIT;
00543       break;
00544     }
00545   case PFBlockLink::HFEMandHFHAD:
00546     {
00547       // cout<<"HFEMandHFHAD"<<endl;
00548       PFClusterRef  eref = lowEl->clusterRef();
00549       PFClusterRef  href = highEl->clusterRef();
00550       assert( !eref.isNull() );
00551       assert( !href.isNull() );
00552       dist = LinkByRecHit::testHFEMAndHFHADByRecHit( *eref, *href, debug_ );
00553       //    linktest = PFBlock::LINKTEST_RECHIT;
00554       break;      
00555     }
00556 
00557   case PFBlockLink::TRACKandTRACK:
00558     {
00559       if(debug_ ) 
00560         cout<<"TRACKandTRACK"<<endl;
00561       dist = testLinkByVertex(lowEl, highEl);
00562       if(debug_ ) 
00563         std::cout << " PFBlockLink::TRACKandTRACK dist " << dist << std::endl;
00564       //   linktest = PFBlock::LINKTEST_RECHIT;
00565       break;
00566     }
00567 
00568   case PFBlockLink::ECALandECAL:
00569       {
00570         
00571         PFClusterRef  ecal1ref = lowEl->clusterRef();
00572         PFClusterRef  ecal2ref = highEl->clusterRef();
00573         assert( !ecal1ref.isNull() );
00574         assert( !ecal2ref.isNull() );
00575         if(debug_)
00576           cout << " PFBlockLink::ECALandECAL" << endl;
00577         dist = testLinkBySuperCluster(ecal1ref,ecal2ref);
00578         break;
00579       }
00580 
00581   case PFBlockLink::ECALandGSF:
00582     {
00583       PFClusterRef  clusterref = lowEl->clusterRef();
00584       assert( !clusterref.isNull() );
00585       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00586       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00587       if ( myTrack->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax ).isValid() )
00588         dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, false, debug_ );
00589       else
00590         dist = -1.;
00591       //   linktest = PFBlock::LINKTEST_RECHIT;
00592       
00593       if ( debug_ ) {
00594         if ( dist > 0. ) {
00595           std::cout << " Here a link has been established" 
00596                     << " between a GSF track an Ecal with dist  " 
00597                     << dist <<  std::endl;
00598         } else {
00599           if(debug_ ) std::cout << " No link found " << std::endl;
00600         }
00601       }
00602       break;
00603     }
00604   case PFBlockLink::TRACKandGSF:
00605     {
00606       PFRecTrackRef trackref = lowEl->trackRefPF();
00607       assert( !trackref.isNull() );
00608       const reco::PFBlockElementGsfTrack *  GsfEl =  
00609         dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00610       GsfPFRecTrackRef gsfref = GsfEl->GsftrackRefPF();
00611       reco::TrackRef kftrackref= (*trackref).trackRef();
00612       assert( !gsfref.isNull() );
00613       PFRecTrackRef refkf = (*gsfref).kfPFRecTrackRef();
00614       if(refkf.isNonnull()) {
00615         reco::TrackRef gsftrackref = (*refkf).trackRef();
00616         if (gsftrackref.isNonnull()&&kftrackref.isNonnull()) {
00617           if (kftrackref == gsftrackref) { 
00618             dist = 0.001;
00619           } else { 
00620             dist = -1.;
00621           }
00622         } else { 
00623           dist = -1.;
00624         }
00625       } else {
00626         dist = -1.;
00627       }
00628       
00629       
00630       if(useConvBremPFRecTracks_) {
00631         if(lowEl->isLinkedToDisplacedVertex()){
00632           vector<PFRecTrackRef> pfrectrack_vec = GsfEl->GsftrackRefPF()->convBremPFRecTrackRef();
00633           if(pfrectrack_vec.size() > 0){
00634             for(unsigned int iconv = 0; iconv <  pfrectrack_vec.size(); ++iconv) {
00635               if( lowEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)) {
00636                 // use track ref
00637                 if(kftrackref == (*pfrectrack_vec[iconv]).trackRef()) {         
00638                   dist = 0.001;
00639                 }
00640               } 
00641               else{
00642                 // use the track base ref
00643                 reco::TrackBaseRef newTrackBaseRef((*pfrectrack_vec[iconv]).trackRef());
00644                 reco::TrackBaseRef elemTrackBaseRef(kftrackref);              
00645                 if(newTrackBaseRef == elemTrackBaseRef){
00646                   dist = 0.001;
00647                 } 
00648               }
00649             }
00650           }
00651         }
00652       }
00653  
00654 
00655       break;      
00656     }
00657          
00658   case PFBlockLink::GSFandBREM:
00659     {
00660       const reco::PFBlockElementGsfTrack * GsfEl  =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(lowEl);
00661       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00662       GsfPFRecTrackRef gsfref = GsfEl->GsftrackRefPF();
00663       GsfPFRecTrackRef bremref = BremEl->GsftrackRefPF();
00664       assert( !gsfref.isNull() );
00665       assert( !bremref.isNull() );
00666       if (gsfref == bremref)  { 
00667         dist = 0.001;
00668       } else { 
00669         dist = -1.;
00670       }
00671       break;
00672     }
00673   case PFBlockLink::GSFandGSF:
00674     {
00675       const reco::PFBlockElementGsfTrack * lowGsfEl  =  
00676         dynamic_cast<const reco::PFBlockElementGsfTrack*>(lowEl);
00677       const reco::PFBlockElementGsfTrack * highGsfEl  =  
00678         dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00679       
00680       GsfPFRecTrackRef lowgsfref = lowGsfEl->GsftrackRefPF();
00681       GsfPFRecTrackRef highgsfref = highGsfEl->GsftrackRefPF();
00682       assert( !lowgsfref.isNull() );
00683       assert( !highgsfref.isNull() );
00684       
00685       if( (lowGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) == false && 
00686            highGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)) ||
00687           (highGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV) == false && 
00688            lowGsfEl->trackType(reco::PFBlockElement::T_FROM_GAMMACONV))) {
00689         if(lowgsfref->trackId() == highgsfref->trackId()) {
00690           dist = 0.001;
00691         }
00692         else {
00693           dist = -1.;
00694         }
00695       }
00696       break;
00697     }
00698   case PFBlockLink::ECALandBREM:
00699     {
00700       PFClusterRef  clusterref = lowEl->clusterRef();
00701       assert( !clusterref.isNull() );
00702       const reco::PFBlockElementBrem * BremEl =  
00703         dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00704       const PFRecTrack * myTrack = &(BremEl->trackPF());
00705       /*
00706       double DP = (BremEl->DeltaP())*(-1.);
00707       double SigmaDP = BremEl->SigmaDeltaP();
00708       double SignBremDp = DP/SigmaDP;
00709       */
00710       bool isBrem = true;
00711       if ( myTrack->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax ).isValid() )
00712         dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, isBrem, debug_);
00713       else
00714         dist = -1.;
00715       if( debug_ && dist > 0. ) 
00716         std::cout << "ECALandBREM: dist testTrackAndClusterByRecHit " 
00717                   << dist << std::endl;
00718       //   linktest = PFBlock::LINKTEST_RECHIT;
00719       break;
00720     }
00721   case PFBlockLink::HCALandGSF:
00722     {
00723       PFClusterRef  clusterref = lowEl->clusterRef();
00724       assert( !clusterref.isNull() );
00725       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00726       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00727       if ( myTrack->extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance ).isValid() )
00728         dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, false, debug_ );
00729       else
00730         dist = -1.;
00731 
00732       //    linktest = PFBlock::LINKTEST_RECHIT;
00733       break;
00734     }
00735   case PFBlockLink::HCALandBREM:
00736     {
00737       PFClusterRef  clusterref = lowEl->clusterRef();
00738       assert( !clusterref.isNull() );
00739       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00740       const PFRecTrack * myTrack = &(BremEl->trackPF());
00741       bool isBrem = true;
00742       if ( myTrack->extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance ).isValid() )
00743         dist = LinkByRecHit::testTrackAndClusterByRecHit( *myTrack, *clusterref, isBrem, debug_);
00744       else
00745         dist = -1.;
00746       break;
00747     }
00748   case PFBlockLink::SCandECAL:
00749     {
00750       PFClusterRef  clusterref = lowEl->clusterRef();
00751 
00752       assert( !clusterref.isNull() );
00753       
00754       const reco::PFBlockElementSuperCluster * scEl = 
00755         dynamic_cast<const reco::PFBlockElementSuperCluster*>(highEl);
00756       assert (!scEl->superClusterRef().isNull());
00757       dist = testSuperClusterPFCluster(scEl->superClusterRef(),
00758                                        clusterref);
00759       break;
00760     }
00761     /*
00762   // Links between the two preshower layers are not used for now - disable
00763   case PFBlockLink::PS1andPS2:
00764     {
00765       PFClusterRef  ps1ref = lowEl->clusterRef();
00766       PFClusterRef  ps2ref = highEl->clusterRef();
00767       assert( !ps1ref.isNull() );
00768       assert( !ps2ref.isNull() );
00769       // PJ - 14-May-09 : A link by rechit is needed here !
00770       // dist = testPS1AndPS2( *ps1ref, *ps2ref );
00771       dist = -1.;
00772       linktest = PFBlock::LINKTEST_RECHIT;
00773       break;
00774     }
00775   case PFBlockLink::TRACKandPS1:
00776   case PFBlockLink::TRACKandPS2:
00777     {
00778       //cout<<"TRACKandPS"<<endl;
00779       PFRecTrackRef trackref = lowEl->trackRefPF();
00780       PFClusterRef  clusterref = highEl->clusterRef();
00781       assert( !trackref.isNull() );
00782       assert( !clusterref.isNull() );
00783       // PJ - 14-May-09 : A link by rechit is needed here !
00784       dist = testTrackAndPS( *trackref, *clusterref );
00785       linktest = PFBlock::LINKTEST_RECHIT;
00786       break;
00787     }
00788     // GSF Track/Brem Track and preshower cluster links are not used for now - disable
00789   case PFBlockLink::PS1andGSF:
00790   case PFBlockLink::PS2andGSF:
00791     {
00792       PFClusterRef  psref = lowEl->clusterRef();
00793       assert( !psref.isNull() );
00794       const reco::PFBlockElementGsfTrack *  GsfEl =  dynamic_cast<const reco::PFBlockElementGsfTrack*>(highEl);
00795       const PFRecTrack * myTrack =  &(GsfEl->GsftrackPF());
00796       // PJ - 14-May-09 : A link by rechit is needed here !
00797       dist = testTrackAndPS( *myTrack, *psref );
00798       linktest = PFBlock::LINKTEST_RECHIT;
00799       break;
00800     }
00801   case PFBlockLink::PS1andBREM:
00802   case PFBlockLink::PS2andBREM:
00803     {
00804       PFClusterRef  psref = lowEl->clusterRef();
00805       assert( !psref.isNull() );
00806       const reco::PFBlockElementBrem * BremEl =  dynamic_cast<const reco::PFBlockElementBrem*>(highEl);
00807       const PFRecTrack * myTrack = &(BremEl->trackPF());
00808       // PJ - 14-May-09 : A link by rechit is needed here !
00809       dist = testTrackAndPS( *myTrack, *psref );
00810       linktest = PFBlock::LINKTEST_RECHIT;
00811       break;
00812     }
00813     */
00814 
00815   default:
00816     dist = -1.;
00817     //   linktest = PFBlock::LINKTEST_RECHIT;
00818     // cerr<<"link type not implemented yet: 0x"<<hex<<linktype<<dec<<endl;
00819     // assert(0);
00820     return;
00821   }
00822 }
00823 
00824 double
00825 PFBlockAlgo::testTrackAndPS(const PFRecTrack& track, 
00826                             const PFCluster& ps)  const {
00827 
00828 #ifdef PFLOW_DEBUG
00829   //   cout<<"entering testTrackAndPS"<<endl;
00830   // resolution of PS cluster dxdx and dydy from strip pitch and length
00831   double dx=0.;
00832   double dy=0.;
00833   
00834   unsigned layerid =0;
00835   // PS1: vertical strips  PS2: horizontal strips
00836   switch (ps.layer()) {
00837   case PFLayer::PS1:
00838     layerid = reco::PFTrajectoryPoint::PS1;
00839     
00840     // vertical strips in PS1, measure x with pitch precision
00841     dx = resPSpitch_;
00842     dy = resPSlength_; 
00843     break;
00844   case PFLayer::PS2:
00845     layerid = reco::PFTrajectoryPoint::PS2;
00846     // horizontal strips in PS2, measure y with pitch precision
00847     dy = resPSpitch_;
00848     dx = resPSlength_;
00849     break;
00850   default:
00851     break;
00852   }
00853   const reco::PFTrajectoryPoint& atPS
00854     = track.extrapolatedPoint( layerid );  
00855   // did not reach PS, cannot be associated with a cluster.
00856   if( ! atPS.isValid() ) return -1.;   
00857   
00858   double trackx = atPS.position().X();
00859   double tracky = atPS.position().Y();
00860   double trackz = atPS.position().Z(); // MDN jan 09
00861   
00862   // ps position  x, y
00863   double psx = ps.position().X();
00864   double psy = ps.position().Y();
00865   // MDN Jan 09: check that trackz and psz have the same sign
00866   double psz = ps.position().Z();
00867   if( trackz*psz < 0.) return -1.; 
00868   
00869   // double chi2 = (psx-trackx)*(psx-trackx)/(dx*dx + trackresolx*trackresolx)
00870   //  + (psy-tracky)*(psy-tracky)/(dy*dy + trackresoly*trackresoly);
00871 
00872   double dist = std::sqrt( (psx-trackx)*(psx-trackx)
00873                          + (psy-tracky)*(psy-tracky));  
00874   if(debug_) cout<<"testTrackAndPS "<< dist <<" "<<endl;
00875   if(debug_){
00876     cout<<" trackx " << trackx 
00877         <<" tracky " << tracky 
00878         <<" psx "    <<  psx   
00879         <<" psy "    << psy    
00880         << endl;
00881   }
00882 #endif
00883   
00884   // Return -1. as long as no link by rechit is available
00885   return -1.;
00886 }
00887 
00888 double
00889 PFBlockAlgo::testECALAndHCAL(const PFCluster& ecal, 
00890                              const PFCluster& hcal)  const {
00891   
00892   //   cout<<"entering testECALAndHCAL"<<endl;
00893   
00894   double dist = fabs(ecal.positionREP().Eta()) > 2.5 ?
00895     LinkByRecHit::computeDist( ecal.positionREP().Eta(),
00896                                ecal.positionREP().Phi(), 
00897                                hcal.positionREP().Eta(), 
00898                                hcal.positionREP().Phi() )
00899     : 
00900     -1.;
00901 
00902 #ifdef PFLOW_DEBUG
00903   if(debug_) cout<<"testECALAndHCAL "<< dist <<" "<<endl;
00904   if(debug_){
00905     cout<<" ecaleta " << ecal.positionREP().Eta()
00906         <<" ecalphi " << ecal.positionREP().Phi()
00907         <<" hcaleta " << hcal.positionREP().Eta()
00908         <<" hcalphi " << hcal.positionREP().Phi()
00909   }
00910 #endif
00911 
00912   if ( dist < 0.2 ) return dist; 
00913  
00914   // Need to implement a link by RecHit
00915   return -1.;
00916 }
00917 
00918 double
00919 PFBlockAlgo::testLinkBySuperCluster(const PFClusterRef& ecal1, 
00920                                     const PFClusterRef& ecal2)  const {
00921   
00922   //  cout<<"entering testECALAndECAL "<< pfcRefSCMap_.size() << endl;
00923   
00924   double dist = -1;
00925   
00926   // the first one is not in any super cluster
00927   int testindex=pfcSCVec_[ecal1.key()];
00928   if(testindex == -1.) return dist;
00929   //  if(itcheck==pfcRefSCMap_.end()) return dist;
00930   // now retrieve the of PFclusters in this super cluster  
00931 
00932   const std::vector<reco::PFClusterRef> & thePFClusters(scpfcRefs_[testindex]);
00933   
00934   unsigned npf=thePFClusters.size();
00935   for(unsigned i=0;i<npf;++i)
00936     {
00937       if(thePFClusters[i]==ecal2) // yes they are in the same SC 
00938         {
00939           dist=LinkByRecHit::computeDist( ecal1->positionREP().Eta(),
00940                                           ecal1->positionREP().Phi(), 
00941                                           ecal2->positionREP().Eta(), 
00942                                           ecal2->positionREP().Phi() );
00943 //        std::cout << " DETA " << fabs(ecal1->positionREP().Eta()-ecal2->positionREP().Eta()) << std::endl;
00944 //        if(fabs(ecal1->positionREP().Eta()-ecal2->positionREP().Eta())>0.2)
00945 //          {
00946 //            std::cout <<  " Super Cluster " <<  *(superClusters_[testindex]) << std::endl;
00947 //            std::cout <<  " Cluster1 " <<  *ecal1 << std::endl;
00948 //            std::cout <<  " Cluster2 " <<  *ecal2 << std::endl;
00949 //            ClusterClusterMapping::checkOverlap(*ecal1,superClusters_,0.01,true);
00950 //            ClusterClusterMapping::checkOverlap(*ecal2,superClusters_,0.01,true);
00951 //          }
00952           return dist;
00953         }
00954     }
00955   return dist;
00956 }
00957 
00958 
00959 double
00960 PFBlockAlgo::testSuperClusterPFCluster(const SuperClusterRef& ecal1, 
00961                                        const PFClusterRef& ecal2)  const {
00962   
00963   //  cout<<"entering testECALAndECAL "<< pfcRefSCMap_.size() << endl;
00964   
00965   double dist = -1;
00966   
00967   bool overlap=ClusterClusterMapping::overlap(*ecal1,*ecal2);
00968   
00969   if(overlap)   {
00970     dist=LinkByRecHit::computeDist( ecal1->position().eta(),
00971                                     ecal1->position().phi(), 
00972                                     ecal2->positionREP().Eta(), 
00973                                     ecal2->positionREP().Phi() );
00974     return dist;
00975   }
00976   return dist;
00977 }
00978 
00979 
00980 
00981 double
00982 PFBlockAlgo::testPS1AndPS2(const PFCluster& ps1, 
00983                            const PFCluster& ps2)  const {
00984   
00985 #ifdef PFLOW_DEBUG
00986   //   cout<<"entering testPS1AndPS2"<<endl;
00987   
00988   // compute chi2 in y, z using swimming formulae
00989   // y2 = y1 * z2/z1   and x2 = x1 *z2/z1
00990   
00991   // ps position1  x, y, z
00992   double x1 = ps1.position().X();
00993   double y1 = ps1.position().Y();
00994   double z1 = ps1.position().Z();
00995   double x2 = ps2.position().X();
00996   double y2 = ps2.position().Y();
00997   double z2 = ps2.position().Z();
00998   // MDN Bug correction Jan 09: check that z1 and z2 have the same sign!
00999   if (z1*z2<0.) -1.;
01000   // swim to PS2
01001   double scale = z2/z1;
01002   double x1atPS2 = x1*scale;
01003   double y1atPS2 = y1*scale;
01004   // resolution of PS cluster dxdx and dydy from strip pitch and length
01005   // vertical strips in PS1, measure x with pitch precision
01006   double dx1dx1 = resPSpitch_*resPSpitch_*scale*scale;
01007   double dy1dy1 = resPSlength_*resPSlength_*scale*scale;
01008   // horizontal strips in PS2 , measure y with pitch precision
01009   double dy2dy2 = resPSpitch_*resPSpitch_;
01010   double dx2dx2 = resPSlength_*resPSlength_;
01011   
01012   // double chi2 = (x2-x1atPS2)*(x2-x1atPS2)/(dx1dx1 + dx2dx2) 
01013   //  + (y2-y1atPS2)*(y2-y1atPS2)/(dy1dy1 + dy2dy2);
01014   
01015   double dist = std::sqrt( (x2-x1atPS2)*(x2-x1atPS2)
01016                          + (y2-y1atPS2)*(y2-y1atPS2));
01017     
01018   if(debug_) cout<<"testPS1AndPS2 "<<dist<<" "<<endl;
01019   if(debug_){
01020     cout<<" x1atPS2 "<< x1atPS2 << " dx1 "<<resPSpitch_*scale
01021         <<" y1atPS2 "<< y1atPS2 << " dy1 "<<resPSlength_*scale<< endl
01022         <<" x2 " <<x2  << " dx2 "<<resPSlength_
01023         <<" y2 " << y2 << " dy2 "<<resPSpitch_<< endl;
01024   }
01025 #endif
01026 
01027   // Need a link by rechit here
01028   return -1.; 
01029 }
01030 
01031 
01032 
01033 double
01034 PFBlockAlgo::testLinkByVertex( const reco::PFBlockElement* elt1, 
01035                                const reco::PFBlockElement* elt2) const {
01036 
01037   //  cout << "Test link by vertex between" << endl << *elt1 << endl << " and " << endl << *elt2 << endl;
01038 
01039   double result=-1.;
01040 
01041   reco::PFBlockElement::TrackType T_TO_DISP = reco::PFBlockElement::T_TO_DISP;
01042   reco::PFBlockElement::TrackType T_FROM_DISP = reco::PFBlockElement::T_FROM_DISP;
01043   PFDisplacedTrackerVertexRef ni1_TO_DISP = elt1->displacedVertexRef(T_TO_DISP);
01044   PFDisplacedTrackerVertexRef ni2_TO_DISP = elt2->displacedVertexRef(T_TO_DISP);
01045   PFDisplacedTrackerVertexRef ni1_FROM_DISP = elt1->displacedVertexRef(T_FROM_DISP);
01046   PFDisplacedTrackerVertexRef ni2_FROM_DISP = elt2->displacedVertexRef(T_FROM_DISP);
01047   
01048   if( ni1_TO_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
01049     if( ni1_TO_DISP == ni2_FROM_DISP ) { result = 1.0; return result; }
01050 
01051   if( ni1_FROM_DISP.isNonnull() && ni2_TO_DISP.isNonnull())
01052     if( ni1_FROM_DISP == ni2_TO_DISP ) { result = 1.0; return result; }
01053 
01054   if( ni1_FROM_DISP.isNonnull() && ni2_FROM_DISP.isNonnull())
01055     if( ni1_FROM_DISP == ni2_FROM_DISP ) { result = 1.0; return result; }
01056     
01057   
01058   if (  elt1->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)  &&
01059              elt2->trackType(reco::PFBlockElement::T_FROM_GAMMACONV)  ) {
01060     
01061     if(debug_ ) std::cout << " testLinkByVertex On Conversions " << std::endl;
01062     
01063     if ( elt1->convRef().isNonnull() && elt2->convRef().isNonnull() ) {
01064       if(debug_ ) std::cout << " PFBlockAlgo.cc testLinkByVertex  Cconversion Refs are non null  " << std::endl;      
01065       if ( elt1->convRef() ==  elt2->convRef() ) {
01066         result=1.0;
01067         if(debug_ ) std::cout << " testLinkByVertex  Cconversion Refs are equal  " << std::endl;      
01068         return result;
01069       }
01070     } 
01071     
01072   }
01073   
01074   if (  elt1->trackType(reco::PFBlockElement::T_FROM_V0)  &&
01075              elt2->trackType(reco::PFBlockElement::T_FROM_V0)  ) {
01076     if(debug_ ) std::cout << " testLinkByVertex On V0 " << std::endl;
01077     if ( elt1->V0Ref().isNonnull() && elt2->V0Ref().isNonnull() ) {
01078       if(debug_ ) std::cout << " PFBlockAlgo.cc testLinkByVertex  V0 Refs are non null  " << std::endl;
01079       if ( elt1->V0Ref() ==  elt2->V0Ref() ) {
01080         result=1.0;
01081         if(debug_ ) std::cout << " testLinkByVertex  V0 Refs are equal  " << std::endl;
01082         return result;
01083       }
01084     }
01085   }
01086 
01087   return result;
01088 }
01089 
01090 
01091 
01092 void 
01093 PFBlockAlgo::checkMaskSize( const reco::PFRecTrackCollection& tracks,
01094                             const reco::GsfPFRecTrackCollection& gsftracks, 
01095                             const reco::PFClusterCollection&  ecals,
01096                             const reco::PFClusterCollection&  hcals,
01097                             const reco::PFClusterCollection&  hfems,
01098                             const reco::PFClusterCollection&  hfhads,
01099                             const reco::PFClusterCollection&  pss, 
01100                             const reco::PhotonCollection&  egphh, 
01101                             const Mask& trackMask, 
01102                             const Mask& gsftrackMask,  
01103                             const Mask& ecalMask,
01104                             const Mask& hcalMask, 
01105                             const Mask& hfemMask,
01106                             const Mask& hfhadMask,                    
01107                             const Mask& psMask,
01108                             const Mask& phMask) const {
01109 
01110   if( !trackMask.empty() && 
01111       trackMask.size() != tracks.size() ) {
01112     string err = "PFBlockAlgo::setInput: ";
01113     err += "The size of the track mask is different ";
01114     err += "from the size of the track vector.";
01115     throw std::length_error( err.c_str() );
01116   }
01117 
01118   if( !gsftrackMask.empty() && 
01119       gsftrackMask.size() != gsftracks.size() ) {
01120     string err = "PFBlockAlgo::setInput: ";
01121     err += "The size of the gsf track mask is different ";
01122     err += "from the size of the gsftrack vector.";
01123     throw std::length_error( err.c_str() );
01124   }
01125 
01126   if( !ecalMask.empty() && 
01127       ecalMask.size() != ecals.size() ) {
01128     string err = "PFBlockAlgo::setInput: ";
01129     err += "The size of the ecal mask is different ";
01130     err += "from the size of the ecal clusters vector.";
01131     throw std::length_error( err.c_str() );
01132   }
01133   
01134   if( !hcalMask.empty() && 
01135       hcalMask.size() != hcals.size() ) {
01136     string err = "PFBlockAlgo::setInput: ";
01137     err += "The size of the hcal mask is different ";
01138     err += "from the size of the hcal clusters vector.";
01139     throw std::length_error( err.c_str() );
01140   }
01141 
01142   if( !hfemMask.empty() && 
01143       hfemMask.size() != hfems.size() ) {
01144     string err = "PFBlockAlgo::setInput: ";
01145     err += "The size of the hfem mask is different ";
01146     err += "from the size of the hfem clusters vector.";
01147     throw std::length_error( err.c_str() );
01148   }
01149 
01150   if( !hfhadMask.empty() && 
01151       hfhadMask.size() != hfhads.size() ) {
01152     string err = "PFBlockAlgo::setInput: ";
01153     err += "The size of the hfhad mask is different ";
01154     err += "from the size of the hfhad clusters vector.";
01155     throw std::length_error( err.c_str() );
01156   }
01157 
01158   if( !psMask.empty() && 
01159       psMask.size() != pss.size() ) {
01160     string err = "PFBlockAlgo::setInput: ";
01161     err += "The size of the ps mask is different ";
01162     err += "from the size of the ps clusters vector.";
01163     throw std::length_error( err.c_str() );
01164   }
01165   
01166     if( !phMask.empty() && 
01167       phMask.size() != egphh.size() ) {
01168     string err = "PFBlockAlgo::setInput: ";
01169     err += "The size of the photon mask is different ";
01170     err += "from the size of the photon vector.";
01171     throw std::length_error( err.c_str() );
01172   }
01173 
01174 }
01175 
01176 
01177 std::ostream& operator<<(std::ostream& out, const PFBlockAlgo& a) {
01178   if(! out) return out;
01179   
01180   out<<"====== Particle Flow Block Algorithm ======= ";
01181   out<<endl;
01182   out<<"number of unassociated elements : "<<a.elements_.size()<<endl;
01183   out<<endl;
01184   
01185   for(PFBlockAlgo::IEC ie = a.elements_.begin(); 
01186       ie != a.elements_.end(); ++ie) {
01187     out<<"\t"<<**ie <<endl;
01188   }
01189 
01190   
01191   //   const PFBlockCollection& blocks = a.blocks();
01192 
01193   const std::auto_ptr< reco::PFBlockCollection >& blocks
01194     = a.blocks(); 
01195     
01196   if(!blocks.get() ) {
01197     out<<"blocks already transfered"<<endl;
01198   }
01199   else {
01200     out<<"number of blocks : "<<blocks->size()<<endl;
01201     out<<endl;
01202     
01203     for(PFBlockAlgo::IBC ib=blocks->begin(); 
01204         ib != blocks->end(); ++ib) {
01205       out<<(*ib)<<endl;
01206     }
01207   }
01208 
01209   return out;
01210 }
01211 
01212 bool 
01213 PFBlockAlgo::goodPtResolution( const reco::TrackRef& trackref) {
01214 
01215   double P = trackref->p();
01216   double Pt = trackref->pt();
01217   double DPt = trackref->ptError();
01218   unsigned int NHit = trackref->hitPattern().trackerLayersWithMeasurement();
01219   unsigned int NLostHit = trackref->hitPattern().trackerLayersWithoutMeasurement();
01220   unsigned int LostHits = trackref->numberOfLostHits();
01221   double sigmaHad = sqrt(1.20*1.20/P+0.06*0.06) / (1.+LostHits);
01222 
01223   // iteration 1,2,3,4,5 correspond to algo = 1/4,5,6,7,8,9
01224   unsigned int Algo = 0; 
01225   switch (trackref->algo()) {
01226   case TrackBase::ctf:
01227   case TrackBase::iter0:
01228   case TrackBase::iter1:
01229   case TrackBase::iter2:
01230     Algo = 0;
01231     break;
01232   case TrackBase::iter3:
01233     Algo = 1;
01234     break;
01235   case TrackBase::iter4:
01236     Algo = 2;
01237     break;
01238   case TrackBase::iter5:
01239     Algo = 3;
01240     break;
01241   case TrackBase::iter6:
01242     Algo = 4;
01243     break;
01244   default:
01245     Algo = useIterTracking_ ? 5 : 0;
01246     break;
01247   }
01248 
01249   // Protection against 0 momentum tracks
01250   if ( P < 0.05 ) return false;
01251 
01252   // Temporary : Reject all tracking iteration beyond 5th step. 
01253   if ( Algo > 4 ) return false;
01254  
01255   if (debug_) cout << " PFBlockAlgo: PFrecTrack->Track Pt= "
01256                    << Pt << " DPt = " << DPt << endl;
01257   if ( ( DPtovPtCut_[Algo] > 0. && 
01258          DPt/Pt > DPtovPtCut_[Algo]*sigmaHad ) || 
01259        NHit < NHitCut_[Algo] ) { 
01260     // (Algo >= 3 && LostHits != 0) ) {
01261     if (debug_) cout << " PFBlockAlgo: skip badly measured track"
01262                      << ", P = " << P 
01263                      << ", Pt = " << Pt 
01264                      << " DPt = " << DPt 
01265                      << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01266                      << ", Algo = " << Algo
01267                      << endl;
01268     if (debug_) cout << " cut is DPt/Pt < " << DPtovPtCut_[Algo] * sigmaHad << endl;
01269     if (debug_) cout << " cut is NHit >= " << NHitCut_[Algo] << endl;
01270     /*
01271     std::cout << "Track REJECTED : ";
01272     std::cout << ", P = " << P 
01273               << ", Pt = " << Pt 
01274               << " DPt = " << DPt 
01275               << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01276               << ", Algo = " << Algo
01277               << std::endl;
01278     */
01279     return false;
01280   }
01281 
01282   /*
01283   std::cout << "Track Accepted : ";
01284   std::cout << ", P = " << P 
01285        << ", Pt = " << Pt 
01286        << " DPt = " << DPt 
01287        << ", N(hits) = " << NHit << " (Lost : " << LostHits << "/" << NLostHit << ")"
01288        << ", Algo = " << Algo
01289        << std::endl;
01290   */
01291   return true;
01292 }
01293 
01294 int
01295 PFBlockAlgo::muAssocToTrack( const reco::TrackRef& trackref,
01296                              const edm::Handle<reco::MuonCollection>& muonh) const {
01297   if(muonh.isValid() ) {
01298     for(unsigned j=0;j<muonh->size(); ++j) {
01299       reco::MuonRef muonref( muonh, j );
01300       if (muonref->track().isNonnull()) 
01301         if( muonref->track() == trackref ) return j;
01302     }
01303   }
01304   return -1; // not found
01305 }
01306 
01307 int 
01308 PFBlockAlgo::muAssocToTrack( const reco::TrackRef& trackref,
01309                              const edm::OrphanHandle<reco::MuonCollection>& muonh) const {
01310   if(muonh.isValid() ) {
01311     for(unsigned j=0;j<muonh->size(); ++j) {
01312       reco::MuonRef muonref( muonh, j );
01313       if (muonref->track().isNonnull())
01314         if( muonref->track() == trackref ) return j;
01315     }
01316   }
01317   return -1; // not found
01318 }
01319 
01320 
01321 // This prefilter avoid to call associate when not necessary.
01322 // ACHTUNG!!!! If you introduce new links check that they are not desables here
01323 inline bool
01324 PFBlockAlgo::linkPrefilter(const reco::PFBlockElement* last, const reco::PFBlockElement* next) const {
01325 
01326   PFBlockElement::Type type1 = (last)->type();
01327   PFBlockElement::Type type2 = (next)->type();
01328 
01329   if( type1==type2 ) {
01330     // cannot link 2 elements of the same type. 
01331     // except if the elements are 2 tracks or 2 ECAL
01332     if( type1!=PFBlockElement::TRACK && type1!=PFBlockElement::GSF &&
01333         type1!=PFBlockElement::ECAL) {
01334       return false;
01335     }
01336 
01337     if (type1==PFBlockElement::ECAL && bNoSuperclus_) return false;
01338 
01339     // cannot link two primary tracks  (except if they come from a V0)
01340     if( type1 ==PFBlockElement::TRACK) {
01341       if ( !((last)->isLinkedToDisplacedVertex()) || !((next)->isLinkedToDisplacedVertex())) 
01342       return false;
01343     }
01344   }
01345 
01346   if ((type1 == PFBlockElement::PS1 || type1 == PFBlockElement::PS2) && (type2 != PFBlockElement::ECAL)) return false;
01347   if ((type2 == PFBlockElement::PS1 || type2 == PFBlockElement::PS2) && (type1 != PFBlockElement::ECAL)) return false;
01348   if ((type1 == PFBlockElement::HFEM && type2 != PFBlockElement::HFHAD) || (type1 == PFBlockElement::HFHAD && type2 != PFBlockElement::HFEM)) return false;
01349 
01350   if (useKDTreeTrackEcalLinker_){ 
01351   
01352     if ( type1 == PFBlockElement::TRACK && type2 == PFBlockElement::ECAL)
01353       if ( last->isMultilinksValide()  && last->getMultilinks().size()==0 ) return false;
01354     if ( type2 == PFBlockElement::TRACK && type1 == PFBlockElement::ECAL)
01355       if ( next->isMultilinksValide() && next->getMultilinks().size()==0 ) return false;
01356     if ( type1 == PFBlockElement::PS1 || type1 == PFBlockElement::PS2)
01357       if ( last->isMultilinksValide()  && last->getMultilinks().size()==0 ) return false;
01358     if ( type2 == PFBlockElement::PS1 || type2 == PFBlockElement::PS2)
01359       if ( next->isMultilinksValide() && next->getMultilinks().size()==0 ) return false;
01360 
01361   }
01362 
01363   return true;
01364 
01365 }
01366