CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoMET/METAlgorithms/src/GlobalHaloAlgo.cc

Go to the documentation of this file.
00001 #include "RecoMET/METAlgorithms/interface/GlobalHaloAlgo.h"
00002 /*
00003   [class]:  GlobalHaloAlgo
00004   [authors]: R. Remington, The University of Florida
00005   [description]: See GlobalHaloAlgo.h
00006   [date]: October 15, 2009
00007 */
00008 using namespace std;
00009 using namespace edm;
00010 using namespace reco;
00011 
00012 int Phi_To_HcaliPhi(float phi) 
00013 {
00014   phi = phi < 0 ? phi + 2.*TMath::Pi() : phi ;
00015   float phi_degrees = phi * (360.) / ( 2. * TMath::Pi() ) ;
00016   int iPhi = (int) ( ( phi_degrees/5. ) + 1.);
00017    
00018   return iPhi < 73 ? iPhi : 73 ;
00019 }
00020 
00021 int Phi_To_EcaliPhi(float phi) 
00022 {
00023   phi = phi < 0 ? phi + 2.*TMath::Pi() : phi ;
00024   float phi_degrees = phi * (360.) / ( 2. * TMath::Pi() ) ;
00025   int iPhi = (int) (  phi_degrees  + 1.);
00026    
00027   return iPhi < 361 ? iPhi : 360 ;
00028 }
00029 
00030 GlobalHaloAlgo::GlobalHaloAlgo()
00031 {
00032   // Defaults are "loose"
00033   Ecal_R_Min = 110.;   // Tight: 200.
00034   Ecal_R_Max = 330.;   // Tight: 250. 
00035   Hcal_R_Min = 110.;   // Tight: 220.
00036   Hcal_R_Max = 490.;   // Tight: 350.
00037   
00038 }
00039 
00040 reco::GlobalHaloData GlobalHaloAlgo::Calculate(const CaloGeometry& TheCaloGeometry, const CSCGeometry& TheCSCGeometry, const reco::CaloMET& TheCaloMET, edm::Handle< edm::View<Candidate> >& TheCaloTowers, edm::Handle<CSCSegmentCollection>& TheCSCSegments, edm::Handle<CSCRecHit2DCollection>& TheCSCRecHits,  const CSCHaloData& TheCSCHaloData, const EcalHaloData& TheEcalHaloData, const HcalHaloData& TheHcalHaloData)
00041 {
00042   
00043   GlobalHaloData TheGlobalHaloData;
00044   float METOverSumEt = TheCaloMET.sumEt() ? TheCaloMET.pt() / TheCaloMET.sumEt() : 0 ;
00045   TheGlobalHaloData.SetMETOverSumEt(METOverSumEt);
00046 
00047   //int EcalOverlapping_CSCRecHits[73];
00048   //int EcalOverlapping_CSCSegments[73];
00049 
00050   int EcalOverlapping_CSCRecHits[361];
00051   int EcalOverlapping_CSCSegments[361];
00052   int HcalOverlapping_CSCRecHits[73];
00053   int HcalOverlapping_CSCSegments[73];
00054 
00055   if( TheCSCSegments.isValid() )
00056     {
00057       for(CSCSegmentCollection::const_iterator iSegment = TheCSCSegments->begin(); iSegment != TheCSCSegments->end(); iSegment++) 
00058         {
00059           bool EcalOverlap[361];
00060           bool HcalOverlap[73];
00061           for( int i = 0 ; i < 361 ; i++ ) 
00062             {
00063               EcalOverlap[i] = false;
00064               if( i < 73 ) HcalOverlap[i] = false;
00065             }
00066           
00067           std::vector<CSCRecHit2D> Hits = iSegment->specificRecHits() ;
00068           for(std::vector<CSCRecHit2D>::iterator iHit = Hits.begin() ; iHit != Hits.end(); iHit++ )
00069             {
00070               DetId TheDetUnitId(iHit->geographicalId());
00071               if( TheDetUnitId.det() != DetId::Muon ) continue;
00072               if( TheDetUnitId.subdetId() != MuonSubdetId::CSC ) continue;
00073 
00074               const GeomDetUnit *TheUnit = TheCSCGeometry.idToDetUnit(TheDetUnitId);
00075               LocalPoint TheLocalPosition = iHit->localPosition();  
00076               const BoundPlane& TheSurface = TheUnit->surface();
00077               const GlobalPoint TheGlobalPosition = TheSurface.toGlobal(TheLocalPosition);
00078 
00079               int Hcal_iphi = Phi_To_HcaliPhi( TheGlobalPosition.phi() ) ;
00080               int Ecal_iphi = Phi_To_EcaliPhi( TheGlobalPosition.phi() ) ;
00081               float x = TheGlobalPosition.x(); 
00082               float y = TheGlobalPosition.y();
00083               
00084               float r = TMath::Sqrt( x*x + y*y);
00085               
00086               if( r < Ecal_R_Max && r > Ecal_R_Min )
00087                 EcalOverlap[Ecal_iphi] = true;
00088               if( r < Hcal_R_Max && r > Hcal_R_Max ) 
00089                 HcalOverlap[Hcal_iphi] = true;
00090             }
00091           for( int i = 0 ; i < 361 ; i++ ) 
00092             {
00093               if( EcalOverlap[i] )  EcalOverlapping_CSCSegments[i]++;
00094               if( i < 73 && HcalOverlap[i] )
00095                 HcalOverlapping_CSCSegments[i]++;
00096             }
00097         } 
00098     }
00099   if( TheCSCRecHits.isValid() )
00100     {
00101       for(CSCRecHit2DCollection::const_iterator iCSCRecHit = TheCSCRecHits->begin();   iCSCRecHit != TheCSCRecHits->end(); iCSCRecHit++ )
00102         {
00103           
00104           DetId TheDetUnitId(iCSCRecHit->geographicalId());
00105           if( TheDetUnitId.det() != DetId::Muon ) continue;
00106           if( TheDetUnitId.subdetId() != MuonSubdetId::CSC ) continue;
00107           
00108           const GeomDetUnit *TheUnit = TheCSCGeometry.idToDetUnit(TheDetUnitId);
00109           LocalPoint TheLocalPosition = iCSCRecHit->localPosition();  
00110           const BoundPlane& TheSurface = TheUnit->surface();
00111           const GlobalPoint TheGlobalPosition = TheSurface.toGlobal(TheLocalPosition);
00112           
00113           int Hcaliphi = Phi_To_HcaliPhi( TheGlobalPosition.phi() ) ;
00114           int Ecaliphi = Phi_To_EcaliPhi( TheGlobalPosition.phi() ) ;
00115           float x = TheGlobalPosition.x(); 
00116           float y = TheGlobalPosition.y();
00117           
00118           float r = TMath::Sqrt(x*x + y*y);
00119           
00120           if( r < Ecal_R_Max && r > Ecal_R_Min )
00121             EcalOverlapping_CSCRecHits[Ecaliphi] ++;
00122           if( r < Hcal_R_Max && r > Hcal_R_Max ) 
00123             HcalOverlapping_CSCRecHits[Hcaliphi] ++ ;
00124         }
00125     }  
00126 
00127   // In development....
00128   // Get Ecal Wedges
00129   std::vector<PhiWedge> EcalWedges = TheEcalHaloData.GetPhiWedges();
00130   
00131   // Get Hcal Wedges
00132   std::vector<PhiWedge> HcalWedges = TheHcalHaloData.GetPhiWedges();
00133   
00134   //Get Ref to CSC Tracks
00135   //edm::RefVector<reco::TrackCollection> TheCSCTracks = TheCSCHaloData.GetTracks();
00136   //for(unsigned int i = 0 ; i < TheCSCTracks.size() ; i++ )
00137   //edm::Ref<reco::TrackCollection> iTrack( TheCSCTracks, i );
00138   
00139   // Get global positions of central most rechit of CSC Halo tracks
00140   std::vector<GlobalPoint> TheGlobalPositions = TheCSCHaloData.GetCSCTrackImpactPositions();
00141 
00142   // Container to store Ecal/Hcal iPhi values matched to impact point of CSC tracks
00143   std::vector<int> vEcaliPhi, vHcaliPhi;
00144 
00145   // Keep track of number of calo pointing CSC halo tracks that do not match to Phi wedges
00146   int N_Unmatched_Tracks = 0;  
00147   
00148   for( std::vector<GlobalPoint>::iterator Pos = TheGlobalPositions.begin() ; Pos != TheGlobalPositions.end() ; Pos ++ ) 
00149     {
00150       // Calculate global phi coordinate for central most rechit in the track
00151       float global_phi = Pos->phi();
00152       float global_r = TMath::Sqrt(Pos->x()*Pos->x() + Pos->y()*Pos->y());
00153       
00154       // Convert global phi to iPhi
00155       int global_EcaliPhi = Phi_To_EcaliPhi( global_phi );
00156       int global_HcaliPhi = Phi_To_HcaliPhi( global_phi );
00157       bool MATCHED = false;
00158       
00159       //Loop over Ecal Phi Wedges 
00160       for( std::vector<PhiWedge>::iterator iWedge = EcalWedges.begin() ; iWedge != EcalWedges.end() ; iWedge++ )
00161         {
00162           if( (TMath::Abs( global_EcaliPhi - iWedge->iPhi() ) <= 5 ) && (global_r >  Ecal_R_Min && global_r < Ecal_R_Max ) )
00163             {
00164               bool StoreWedge = true;
00165               for( unsigned int i = 0 ; i< vEcaliPhi.size() ; i++ ) if ( vEcaliPhi[i] == iWedge->iPhi() ) StoreWedge = false;
00166               
00167               if( StoreWedge ) 
00168                 {
00169                   PhiWedge NewWedge(*iWedge);
00170                   NewWedge.SetOverlappingCSCSegments( EcalOverlapping_CSCSegments[iWedge->iPhi()] );
00171                   NewWedge.SetOverlappingCSCRecHits( EcalOverlapping_CSCRecHits[iWedge->iPhi()] );
00172                   vEcaliPhi.push_back( iWedge->iPhi() );
00173                   TheGlobalHaloData.GetMatchedEcalPhiWedges().push_back( NewWedge );
00174                 }
00175               MATCHED = true;
00176             }
00177         }
00178       //Loop over Hcal Phi Wedges 
00179       for( std::vector<PhiWedge>::iterator iWedge = HcalWedges.begin() ; iWedge != HcalWedges.end() ; iWedge++ )
00180         {
00181           if(  (TMath::Abs( global_HcaliPhi - iWedge->iPhi() ) <=  2 ) && (global_r > Hcal_R_Min && global_r < Hcal_R_Max ) )
00182             {
00183               bool StoreWedge  = true;
00184               for( unsigned int i = 0 ; i < vHcaliPhi.size() ; i++ ) if(  vHcaliPhi[i] == iWedge->iPhi() ) StoreWedge = false;
00185               
00186               if( StoreWedge ) 
00187                 {
00188                   vHcaliPhi.push_back( iWedge->iPhi() ) ;
00189                   PhiWedge NewWedge(*iWedge);
00190                   NewWedge.SetOverlappingCSCSegments( HcalOverlapping_CSCSegments[iWedge->iPhi()] );
00191                   NewWedge.SetOverlappingCSCRecHits(  HcalOverlapping_CSCRecHits[iWedge->iPhi()] );               
00192                   PhiWedge wedge(*iWedge);
00193                   TheGlobalHaloData.GetMatchedHcalPhiWedges().push_back( NewWedge ) ; 
00194                 }
00195               MATCHED = true;
00196             }
00197         }
00198       if( !MATCHED ) N_Unmatched_Tracks ++;
00199     }
00200   
00201   // Corrections to MEx, MEy
00202   float dMEx = 0.; 
00203   float dMEy = 0.;
00204   // Loop over calotowers and correct the MET for the towers that lie in the trajectory of the CSC Halo Tracks
00205   for( edm::View<Candidate>::const_iterator iCandidate = TheCaloTowers->begin() ; iCandidate != TheCaloTowers->end() ; iCandidate++ )
00206     {
00207       const Candidate* c = &(*iCandidate);
00208       if ( c )
00209         {
00210           const CaloTower* iTower = dynamic_cast<const CaloTower*> (c);
00211           if( iTower->et() < TowerEtThreshold ) continue;
00212           if( abs(iTower->ieta()) > 24 )  continue;   // not in barrel/endcap
00213           int iphi = iTower->iphi();
00214           for( unsigned int x = 0 ; x < vEcaliPhi.size() ; x++ )
00215             {
00216               if( iphi == vEcaliPhi[x] ) 
00217                 {
00218                   dMEx += ( TMath::Cos(iTower->phi())*iTower->emEt() );
00219                   dMEy += ( TMath::Sin(iTower->phi())*iTower->emEt() );
00220                 }
00221             }
00222           for( unsigned int x = 0 ; x < vHcaliPhi.size() ; x++ )
00223             {
00224               if( iphi == vHcaliPhi[x] ) 
00225                 {
00226                   dMEx += ( TMath::Cos(iTower->phi() )*iTower->hadEt() ) ;
00227                   dMEy += ( TMath::Sin(iTower->phi() )*iTower->hadEt() ) ;
00228                 }
00229             }
00230         }
00231     }
00232   
00233   TheGlobalHaloData.SetMETCorrections(dMEx, dMEy);
00234   return TheGlobalHaloData;
00235 }