CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch13/src/RecoParticleFlow/PFClusterTools/src/LinkByRecHit.cc

Go to the documentation of this file.
00001 #include "RecoParticleFlow/PFClusterTools/interface/LinkByRecHit.h"
00002 #include "DataFormats/Math/interface/normalizedPhi.h"
00003 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
00004 #include "TMath.h"
00005 
00006 // to enable debugs
00007 //#define PFLOW_DEBUG
00008 
00009 double 
00010 LinkByRecHit::testTrackAndClusterByRecHit ( const reco::PFRecTrack& track, 
00011                                             const reco::PFCluster&  cluster,
00012                                             bool isBrem,
00013                                             bool debug)  {
00014   
00015 #ifdef PFLOW_DEBUG
00016   if( debug ) 
00017     std::cout<<"entering test link by rechit function"<<std::endl;
00018 #endif
00019 
00020   //cluster position
00021   double clustereta  = cluster.positionREP().Eta();
00022   double clusterphi  = cluster.positionREP().Phi();
00023   double clusterX    = cluster.position().X();
00024   double clusterY    = cluster.position().Y();
00025   double clusterZ    = cluster.position().Z();
00026 
00027   bool barrel = false;
00028   bool hcal = false;
00029   double distance = 999999.9;
00030 
00031   //track extrapolation
00032   const reco::PFTrajectoryPoint& atVertex 
00033     = track.extrapolatedPoint( reco::PFTrajectoryPoint::ClosestApproach );
00034   const reco::PFTrajectoryPoint& atECAL 
00035     = track.extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax );
00036 
00037 
00038   //track at calo's
00039   double tracketa = 999999.9;
00040   double trackphi = 999999.9;
00041   double track_X  = 999999.9;
00042   double track_Y  = 999999.9;
00043   double track_Z  = 999999.9;
00044   double dHEta = 0.;
00045   double dHPhi = 0.;
00046 
00047   // Quantities at vertex
00048   double trackPt = isBrem ? 999. : sqrt(atVertex.momentum().Vect().Perp2());
00049   // double trackEta = isBrem ? 999. : atVertex.momentum().Vect().Eta();
00050 
00051 
00052   switch (cluster.layer()) {
00053   case PFLayer::ECAL_BARREL: barrel = true;
00054   case PFLayer::ECAL_ENDCAP:
00055 #ifdef PFLOW_DEBUG
00056     if( debug )
00057       std::cout << "Fetching Ecal Resolution Maps"
00058            << std::endl;
00059 #endif
00060     // did not reach ecal, cannot be associated with a cluster.
00061     if( ! atECAL.isValid() ) return -1.;   
00062     
00063     tracketa = atECAL.positionREP().Eta();
00064     trackphi = atECAL.positionREP().Phi();
00065     track_X  = atECAL.position().X();
00066     track_Y  = atECAL.position().Y();
00067     track_Z  = atECAL.position().Z();
00068 
00069     distance 
00070       = std::sqrt( (track_X-clusterX)*(track_X-clusterX)
00071                   +(track_Y-clusterY)*(track_Y-clusterY)
00072                   +(track_Z-clusterZ)*(track_Z-clusterZ)
00073                    );
00074                                    
00075     break;
00076    
00077   case PFLayer::HCAL_BARREL1: barrel = true; 
00078   case PFLayer::HCAL_ENDCAP:  
00079 #ifdef PFLOW_DEBUG
00080     if( debug )
00081       std::cout << "Fetching Hcal Resolution Maps"
00082            << std::endl;
00083 #endif
00084     if( isBrem ) {  
00085       return  -1.;
00086     } else { 
00087       hcal=true;
00088       const reco::PFTrajectoryPoint& atHCAL 
00089         = track.extrapolatedPoint( reco::PFTrajectoryPoint::HCALEntrance );
00090       const reco::PFTrajectoryPoint& atHCALExit 
00091         = track.extrapolatedPoint( reco::PFTrajectoryPoint::HCALExit );
00092       // did not reach hcal, cannot be associated with a cluster.
00093       if( ! atHCAL.isValid() ) return -1.;   
00094       
00095       // The link is computed between 0 and ~1 interaction length in HCAL
00096       dHEta = atHCALExit.positionREP().Eta()-atHCAL.positionREP().Eta();
00097       dHPhi = atHCALExit.positionREP().Phi()-atHCAL.positionREP().Phi(); 
00098       if ( dHPhi > M_PI ) dHPhi = dHPhi - 2.*M_PI;
00099       else if ( dHPhi < -M_PI ) dHPhi = dHPhi + 2.*M_PI; 
00100       tracketa = atHCAL.positionREP().Eta() + 0.1*dHEta;
00101       trackphi = atHCAL.positionREP().Phi() + 0.1*dHPhi;
00102       track_X  = atHCAL.position().X();
00103       track_Y  = atHCAL.position().Y();
00104       track_Z  = atHCAL.position().Z();
00105       distance 
00106         = -std::sqrt( (track_X-clusterX)*(track_X-clusterX)
00107                      +(track_Y-clusterY)*(track_Y-clusterY)
00108                      +(track_Z-clusterZ)*(track_Z-clusterZ)
00109                      );
00110                                    
00111     }
00112     break;
00113   case PFLayer::PS1:
00114   case PFLayer::PS2:
00115     //Note Alex: Nothing implemented for the
00116     //PreShower (No resolution maps yet)
00117 #ifdef PFLOW_DEBUG
00118     if( debug )
00119       std::cout << "No link by rechit possible for pre-shower yet!"
00120            << std::endl;
00121 #endif
00122     return -1.;
00123   default:
00124     return -1.;
00125   }
00126 
00127 
00128   // Check that, if the cluster is in the endcap, 
00129   // 0) the track indeed points to the endcap at vertex (DISABLED)
00130   // 1) the track extrapolation is in the endcap too !
00131   // 2) the track is in the same end-cap !
00132   // PJ - 10-May-09
00133   if ( !barrel ) { 
00134     // if ( fabs(trackEta) < 1.0 ) return -1; 
00135     if ( !hcal && fabs(track_Z) < 300. ) return -1.;
00136     if ( track_Z * clusterZ < 0. ) return -1.;
00137   }
00138   // Check that, if the cluster is in the barrel, 
00139   // 1) the track is in the barrel too !
00140   if ( barrel ) 
00141     if ( !hcal && fabs(track_Z) > 300. ) return -1.;
00142 
00143   // Finally check that, if the track points to the central barrel (|eta| < 1), 
00144   // it cannot be linked to a cluster in Endcaps (avoid low pt loopers)
00145 
00146 
00147   double dist = LinkByRecHit::computeDist( clustereta, clusterphi, 
00148                                          tracketa, trackphi);
00149   
00150 #ifdef PFLOW_DEBUG
00151   if(debug) std::cout<<"test link by rechit "<< dist <<" "<<std::endl;
00152   if(debug){
00153     std::cout<<" clustereta "  << clustereta 
00154         <<" clusterphi "  << clusterphi 
00155         <<" tracketa " << tracketa
00156         <<" trackphi " << trackphi << std::endl;
00157   }
00158 #endif
00159   
00160   //Testing if Track can be linked by rechit to a cluster.
00161   //A cluster can be linked to a track if the extrapolated position 
00162   //of the track to the ECAL ShowerMax/HCAL entrance falls within 
00163   //the boundaries of any cell that belongs to this cluster.
00164 
00165   const std::vector< reco::PFRecHitFraction >& 
00166     fracs = cluster.recHitFractions();
00167   
00168   bool linkedbyrechit = false;
00169   //loop rechits
00170   for(unsigned int rhit = 0; rhit < fracs.size(); ++rhit){
00171 
00172     const reco::PFRecHitRef& rh = fracs[rhit].recHitRef();
00173     double fraction = fracs[rhit].fraction();
00174     if(fraction < 1E-4) continue;
00175     if(rh.isNull()) continue;
00176     
00177     //getting rechit center position
00178     const reco::PFRecHit& rechit_cluster = *rh;
00179     const math::XYZPoint& posxyz 
00180       = rechit_cluster.position();
00181     const reco::PFRecHit::REPPoint& posrep 
00182       = rechit_cluster.positionREP();
00183     
00184     //getting rechit corners
00185     const std::vector< math::XYZPoint >& 
00186       cornersxyz = rechit_cluster.getCornersXYZ();
00187     const std::vector<reco::PFRecHit::REPPoint>& corners = 
00188       rechit_cluster.getCornersREP();
00189     assert(corners.size() == 4);
00190     
00191     if( barrel || hcal ){ // barrel case matching in eta/phi 
00192                           // (and HCAL endcap too!)
00193       
00194       //rechit size determination 
00195       // blown up by 50% (HCAL) to 100% (ECAL) to include cracks & gaps
00196       // also blown up to account for multiple scattering at low pt.
00197       double rhsizeEta 
00198         = fabs(corners[0].Eta() - corners[2].Eta());
00199       double rhsizePhi 
00200         = fabs(corners[0].Phi() - corners[2].Phi());
00201       if ( rhsizePhi > M_PI ) rhsizePhi = 2.*M_PI - rhsizePhi;
00202       if ( hcal ) { 
00203         rhsizeEta = rhsizeEta * (1.50 + 0.5/fracs.size()) + 0.2*fabs(dHEta);
00204         rhsizePhi = rhsizePhi * (1.50 + 0.5/fracs.size()) + 0.2*fabs(dHPhi); 
00205         
00206       } else { 
00207         rhsizeEta *= 2.00 + 1.0/fracs.size()/std::min(1.,trackPt/2.);
00208         rhsizePhi *= 2.00 + 1.0/fracs.size()/std::min(1.,trackPt/2.); 
00209       }
00210       
00211 #ifdef PFLOW_DEBUG
00212       if( debug ) {
00213         std::cout << rhit         << " Hcal RecHit=" 
00214              << posrep.Eta() << " " 
00215              << posrep.Phi() << " "
00216              << rechit_cluster.energy() 
00217              << std::endl; 
00218         for ( unsigned jc=0; jc<4; ++jc ) 
00219           std::cout<<"corners "<<jc<<" "<<corners[jc].Eta()
00220               <<" "<<corners[jc].Phi()<<std::endl;
00221         
00222         std::cout << "RecHit SizeEta=" << rhsizeEta
00223              << " SizePhi=" << rhsizePhi << std::endl;
00224       }
00225 #endif
00226       
00227       //distance track-rechit center
00228       // const math::XYZPoint& posxyz 
00229       // = rechit_cluster.position();
00230       double deta = fabs(posrep.Eta() - tracketa);
00231       double dphi = fabs(posrep.Phi() - trackphi);
00232       if ( dphi > M_PI ) dphi = 2.*M_PI - dphi;
00233       
00234 #ifdef PFLOW_DEBUG
00235       if( debug ){
00236         std::cout << "distance=" 
00237              << deta << " " 
00238              << dphi << " ";
00239         if(deta < (rhsizeEta/2.) && dphi < (rhsizePhi/2.))
00240           std::cout << " link here !" << std::endl;
00241         else std::cout << std::endl;
00242       }
00243 #endif
00244       
00245       if(deta < (rhsizeEta/2.) && dphi < (rhsizePhi/2.)){ 
00246         linkedbyrechit = true;
00247         break;
00248       }
00249     }
00250     else { //ECAL & PS endcap case, matching in X,Y
00251       
00252 #ifdef PFLOW_DEBUG
00253       if( debug ){
00254         const math::XYZPoint& posxyz 
00255           = rechit_cluster.position();
00256         
00257         std::cout << "RH " << posxyz.X()
00258              << " "   << posxyz.Y()
00259              << std::endl;
00260         
00261         std::cout << "TRACK " << track_X
00262              << " "      << track_Y
00263              << std::endl;
00264       }
00265 #endif
00266       
00267       double x[5];
00268       double y[5];
00269       
00270       for ( unsigned jc=0; jc<4; ++jc ) {
00271         math::XYZPoint cornerposxyz = cornersxyz[jc];
00272         x[jc] = cornerposxyz.X() + (cornerposxyz.X()-posxyz.X())
00273           * (1.00+0.50/fracs.size()/std::min(1.,trackPt/2.));
00274         y[jc] = cornerposxyz.Y() + (cornerposxyz.Y()-posxyz.Y())
00275           * (1.00+0.50/fracs.size()/std::min(1.,trackPt/2.));
00276         
00277 #ifdef PFLOW_DEBUG
00278         if( debug ){
00279           std::cout<<"corners "<<jc
00280               << " " << cornerposxyz.X()
00281               << " " << cornerposxyz.Y()
00282               << std::endl;
00283         }
00284 #endif
00285       }//loop corners
00286       
00287       //need to close the polygon in order to
00288       //use the TMath::IsInside fonction from root lib
00289       x[4] = x[0];
00290       y[4] = y[0];
00291       
00292       //Check if the extrapolation point of the track falls 
00293       //within the rechit boundaries
00294       bool isinside = TMath::IsInside(track_X,
00295                                       track_Y,
00296                                       5,x,y);
00297       
00298       if( isinside ){
00299         linkedbyrechit = true;
00300         break;
00301       }
00302     }//
00303     
00304   }//loop rechits
00305   
00306   if( linkedbyrechit ) {
00307 #ifdef PFLOW_DEBUG
00308     if( debug ) 
00309       std::cout << "Track and Cluster LINKED BY RECHIT" << std::endl;
00310 #endif
00311     /*    
00312     //if ( distance > 40. || distance < -100. ) 
00313     double clusterr = std::sqrt(clusterX*clusterX+clusterY*clusterY);
00314     double trackr = std::sqrt(track_X*track_X+track_Y*track_Y);
00315     if ( distance > 40. ) 
00316     std::cout << "Distance = " << distance 
00317     << ", Barrel/Hcal/Brem ? " << barrel << " " << hcal << " " << isBrem << std::endl
00318     << " Cluster " << clusterr << " " << clusterZ << " " << clusterphi << " " << clustereta << std::endl
00319     << " Track   " << trackr << " " << track_Z << " " << trackphi << " " << tracketa << std::endl;
00320     if ( !barrel && fabs(trackEta) < 1.0 ) { 
00321       double clusterr = std::sqrt(clusterX*clusterX+clusterY*clusterY);
00322       double trackr = std::sqrt(track_X*track_X+track_Y*track_Y);
00323       std::cout << "TrackEta/Pt = " << trackEta << " " << trackPt << ", distance = " << distance << std::endl 
00324                 << ", Barrel/Hcal/Brem ? " << barrel << " " << hcal << " " << isBrem << std::endl
00325                 << " Cluster " << clusterr << " " << clusterZ << " " << clusterphi << " " << clustereta << std::endl
00326                 << " Track   " << trackr << " " << track_Z << " " << trackphi << " " << tracketa << " " << trackEta << " " << trackPt << std::endl;
00327     } 
00328     */
00329     return dist;
00330   } else {
00331     return -1.;
00332   }
00333 
00334 }
00335 
00336 
00337 
00338 double
00339 LinkByRecHit::testECALAndPSByRecHit( const reco::PFCluster& clusterECAL, 
00340                                      const reco::PFCluster& clusterPS,
00341                                      bool debug)  {
00342 
00343 // 0.19 <-> strip_pitch
00344 // 6.1  <-> strip_length
00345   static double resPSpitch = 0.19/sqrt(12.);
00346   static double resPSlength = 6.1/sqrt(12.);
00347 
00348   // Check that clusterECAL is in ECAL endcap and that clusterPS is a preshower cluster
00349   if ( clusterECAL.layer() != PFLayer::ECAL_ENDCAP ||
00350        ( clusterPS.layer() != PFLayer::PS1 && 
00351          clusterPS.layer() != PFLayer::PS2 ) ) return -1.;
00352 
00353 #ifdef PFLOW_DEBUG
00354   if( debug ) 
00355     std::cout<<"entering test link by rechit function for ECAL and PS"<<std::endl;
00356 #endif
00357 
00358   //ECAL cluster position
00359   double zECAL  = clusterECAL.position().Z();
00360   double xECAL  = clusterECAL.position().X();
00361   double yECAL  = clusterECAL.position().Y();
00362 
00363   // PS cluster position, extrapolated to ECAL
00364   double zPS = clusterPS.position().Z();
00365   double xPS = clusterPS.position().X(); //* zECAL/zPS;
00366   double yPS = clusterPS.position().Y(); //* zECAL/zPS;
00367 // MDN jan09 : check that zEcal and zPs have the same sign
00368         if (zECAL*zPS <0.) return -1.;
00369   double deltaX = 0.;
00370   double deltaY = 0.;
00371   double sqr12 = std::sqrt(12.);
00372   switch (clusterPS.layer()) {
00373   case PFLayer::PS1:
00374     // vertical strips, measure x with pitch precision
00375     deltaX = resPSpitch * sqr12;
00376     deltaY = resPSlength * sqr12;
00377     break;
00378   case PFLayer::PS2:
00379     // horizontal strips, measure y with pitch precision
00380     deltaY = resPSpitch * sqr12;
00381     deltaX = resPSlength * sqr12;
00382     break;
00383   default:
00384     break;
00385   }
00386 
00387   // Get the rechits
00388   const std::vector< reco::PFRecHitFraction >&  fracs = clusterECAL.recHitFractions();
00389   bool linkedbyrechit = false;
00390   //loop rechits
00391   for(unsigned int rhit = 0; rhit < fracs.size(); ++rhit){
00392 
00393     const reco::PFRecHitRef& rh = fracs[rhit].recHitRef();
00394     double fraction = fracs[rhit].fraction();
00395     if(fraction < 1E-4) continue;
00396     if(rh.isNull()) continue;
00397 
00398     //getting rechit center position
00399     const reco::PFRecHit& rechit_cluster = *rh;
00400     
00401     //getting rechit corners
00402     const std::vector< math::XYZPoint >&  corners = rechit_cluster.getCornersXYZ();
00403     assert(corners.size() == 4);
00404     
00405     const math::XYZPoint& posxyz = rechit_cluster.position() * zPS/zECAL;
00406 #ifdef PFLOW_DEBUG
00407     if( debug ){
00408       std::cout << "Ecal rechit " << posxyz.X() << " "   << posxyz.Y() << std::endl;
00409       std::cout << "PS cluster  " << xPS << " " << yPS << std::endl;
00410     }
00411 #endif
00412     
00413     double x[5];
00414     double y[5];
00415     for ( unsigned jc=0; jc<4; ++jc ) {
00416       // corner position projected onto the preshower
00417       math::XYZPoint cornerpos = corners[jc] * zPS/zECAL;
00418       // Inflate the size by the size of the PS strips, and by 5% to include ECAL cracks.
00419       x[jc] = cornerpos.X() + (cornerpos.X()-posxyz.X()) * (0.05 +1.0/fabs((cornerpos.X()-posxyz.X()))*deltaX/2.);
00420       y[jc] = cornerpos.Y() + (cornerpos.Y()-posxyz.Y()) * (0.05 +1.0/fabs((cornerpos.Y()-posxyz.Y()))*deltaY/2.);
00421       
00422 #ifdef PFLOW_DEBUG
00423       if( debug ){
00424         std::cout<<"corners "<<jc
00425             << " " << cornerpos.X() << " " << x[jc] 
00426             << " " << cornerpos.Y() << " " << y[jc]
00427             << std::endl;
00428       }
00429 #endif
00430     }//loop corners
00431     
00432     //need to close the polygon in order to
00433     //use the TMath::IsInside fonction from root lib
00434     x[4] = x[0];
00435     y[4] = y[0];
00436     
00437     //Check if the extrapolation point of the track falls 
00438     //within the rechit boundaries
00439     bool isinside = TMath::IsInside(xPS,yPS,5,x,y);
00440       
00441     if( isinside ){
00442       linkedbyrechit = true;
00443       break;
00444     }
00445 
00446   }//loop rechits
00447   
00448   if( linkedbyrechit ) {
00449     if( debug ) std::cout << "Cluster PS and Cluster ECAL LINKED BY RECHIT" << std::endl;
00450     double dist = computeDist( xECAL/1000.,yECAL/1000.,
00451                                xPS/1000.  ,yPS/1000);    
00452     return dist;
00453   } else { 
00454     return -1.;
00455   }
00456 
00457 }
00458 
00459 
00460 
00461 double 
00462 LinkByRecHit::testHFEMAndHFHADByRecHit(const reco::PFCluster& clusterHFEM, 
00463                                       const reco::PFCluster& clusterHFHAD,
00464                                       bool debug) {
00465   
00466   math::XYZPoint posxyzEM = clusterHFEM.position();
00467   math::XYZPoint posxyzHAD = clusterHFHAD.position();
00468 
00469   double dX = posxyzEM.X()-posxyzHAD.X();
00470   double dY = posxyzEM.Y()-posxyzHAD.Y();
00471   double sameZ = posxyzEM.Z()*posxyzHAD.Z();
00472 
00473   if(sameZ<0) return -1.;
00474 
00475   double dist2 = dX*dX + dY*dY; 
00476 
00477   if( dist2 < 0.1 ) {
00478     // less than one mm
00479     double dist = sqrt( dist2 );
00480     return dist;;
00481   }
00482   else 
00483     return -1.;
00484 
00485 }
00486 
00487 double
00488 LinkByRecHit::computeDist( double eta1, double phi1, 
00489                           double eta2, double phi2 )  {
00490   
00491   double phicor = normalizedPhi(phi1 - phi2);
00492   
00493   // double chi2 =  
00494   //  (eta1 - eta2)*(eta1 - eta2) / ( reta1*reta1+ reta2*reta2 ) +
00495   //  phicor*phicor / ( rphi1*rphi1+ rphi2*rphi2 );
00496 
00497   double dist = std::sqrt( (eta1 - eta2)*(eta1 - eta2) 
00498                           + phicor*phicor);
00499 
00500   return dist;
00501 
00502 }