CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoHI/HiEgammaAlgos/src/TxCalculator.cc

Go to the documentation of this file.
00001 // ROOT includes
00002 #include <Math/VectorUtil.h>
00003 
00004 #include "RecoHI/HiEgammaAlgos/interface/TxCalculator.h"
00005 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00006 
00007 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00008 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00009 
00010 #include "DataFormats/Common/interface/Handle.h"
00011 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00012 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00013 #include "DataFormats/Math/interface/Vector3D.h"
00014 
00015 #include "FWCore/ServiceRegistry/interface/Service.h"
00016 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00017 
00018 using namespace edm;
00019 using namespace reco;
00020 using namespace std;
00021 using namespace ROOT::Math::VectorUtil; 
00022 
00023 
00024 TxCalculator::TxCalculator (const edm::Event &iEvent, const edm::EventSetup &iSetup, edm::InputTag trackLabel)
00025 {
00026    iEvent.getByLabel(trackLabel, recCollection); 
00027    edm::Service<edm::RandomNumberGenerator> rng;
00028    if ( ! rng.isAvailable()) {
00029       throw cms::Exception("Configuration")
00030          << "XXXXXXX requires the RandomNumberGeneratorService\n"
00031          "which is not present in the configuration file.  You must add the service\n"
00032          "in the configuration file or remove the modules that require it.";
00033    }
00034    CLHEP::HepRandomEngine& engine = rng->getEngine();
00035    theDice = new CLHEP::RandFlat(engine, 0, 1);
00036    
00037 } 
00038 
00039 
00040 double TxCalculator::getJurassicArea( double r1, double r2, double width) {
00041    
00042    float theta1 = asin( width / r1);
00043    float theta2 = asin( width / r2);
00044    float theA   = sqrt ( r1*r1 + r2*r2 - 2 * r1 * r2 * cos ( theta1 - theta2) );
00045    float area1 =  0.5 * r1*r1 * ( 3.141592 - 2 * theta1 )   ;
00046    float area2 =  0.5 * r2*r2 * ( 3.141592 - 2 * theta2 )   ;
00047    float area3 =  width * theA;
00048    float finalArea = 2 * ( area1 - area2 - area3);
00049    return finalArea;
00050 }
00051 
00052 
00053 double TxCalculator::getMPT( double ptCut     ,   double etaCut  )
00054 {
00055    using namespace edm;
00056    using namespace reco;
00057    
00058    double sumpx(0), sumpy(0);
00059    
00060    for(reco::TrackCollection::const_iterator
00061           recTrack = recCollection->begin(); recTrack!= recCollection->end(); recTrack++)
00062       {
00063          double pt = recTrack->pt();
00064          double eta = recTrack->eta();
00065          
00066          if(pt < ptCut ) 
00067             continue;
00068          if ( fabs( eta) > etaCut ) 
00069             continue;
00070          
00071          double pxTemp = recTrack->px();
00072          double pyTemp = recTrack->py();
00073          
00074          sumpx = sumpx + pxTemp;
00075          sumpy = sumpy + pyTemp;
00076          //      cout << " pt  = " << recTrack->pt() <<  "    and px = " << pxTemp << " and  py = " << pyTemp << endl;
00077       }
00078    //   cout << " square = " << sumpx*sumpx + sumpy*sumpy << endl;
00079    double theMPT = sqrt(sumpx*sumpx + sumpy*sumpy) ;
00080    //  cout << " mpt    = "<< theMPT << endl;
00081    
00082    return theMPT;
00083 }
00084 
00085 
00086 double TxCalculator::getTx(const reco::Photon cluster, double x, double threshold, double innerDR, double effRatio)
00087 {
00088 
00089    using namespace edm;
00090    using namespace reco;
00091 
00092    
00093    
00094    double SClusterEta = cluster.eta();
00095    double SClusterPhi = cluster.phi();
00096    double TotalPt = 0;
00097    
00098    for(reco::TrackCollection::const_iterator
00099           recTrack = recCollection->begin(); recTrack!= recCollection->end(); recTrack++)
00100       {
00101          double diceNum = theDice->fire();
00102          if ( (effRatio < 1 ) &&  ( diceNum > effRatio))
00103             continue;
00104          
00105          double pt = recTrack->pt();
00106          double eta2 = recTrack->eta();
00107          double phi2 = recTrack->phi();
00108          
00109          if(dRDistance(SClusterEta,SClusterPhi,eta2,phi2) >= 0.1 * x)
00110             continue;
00111          if(dRDistance(SClusterEta,SClusterPhi,eta2,phi2) < innerDR)
00112             continue;
00113          if(pt > threshold)
00114             TotalPt = TotalPt + pt;
00115       }
00116    
00117    return TotalPt;
00118 }
00119 
00120 
00121 
00122 
00123 
00124 double TxCalculator::getCTx(const reco::Photon cluster, double x, double threshold, double innerDR,double effRatio)
00125 {
00126    using namespace edm;
00127    using namespace reco;
00128 
00129    double SClusterEta = cluster.eta();
00130    double SClusterPhi = cluster.phi();
00131    double TotalPt = 0;
00132 
00133    TotalPt = 0;
00134 
00135    for(reco::TrackCollection::const_iterator
00136           recTrack = recCollection->begin(); recTrack!= recCollection->end(); recTrack++)
00137       {
00138          double diceNum = theDice->fire();
00139          if ( (effRatio < 1 ) &&  ( diceNum > effRatio))
00140             continue;
00141          
00142          
00143          double pt = recTrack->pt();
00144          double eta2 = recTrack->eta();
00145          double phi2 = recTrack->phi();
00146          double dEta = fabs(eta2-SClusterEta);
00147          
00148       if(dEta >= 0.1 * x)
00149          continue;
00150       if(dRDistance(SClusterEta,SClusterPhi,eta2,phi2) < innerDR)
00151          continue;
00152       
00153       if(pt > threshold)
00154          TotalPt = TotalPt + pt;
00155    }
00156    
00157    double Tx = getTx(cluster,x,threshold,innerDR,effRatio);
00158    double CTx = Tx - TotalPt / 40.0 * x;
00159 
00160    return CTx;
00161 }
00162 
00163 
00164 
00165 double TxCalculator::getJt(const reco::Photon cluster, double r1, double r2, double jWidth, double threshold)
00166 {
00167 
00168    using namespace edm;
00169    using namespace reco;
00170 
00171 
00172    double SClusterEta = cluster.eta();
00173    double SClusterPhi = cluster.phi();
00174    double TotalPt = 0;
00175 
00176    for(reco::TrackCollection::const_iterator
00177           recTrack = recCollection->begin(); recTrack!= recCollection->end(); recTrack++)
00178       {
00179          double pt = recTrack->pt();
00180          double eta = recTrack->eta();
00181          double phi = recTrack->phi();
00182          double dEta = fabs(eta-SClusterEta);
00183          double dPhi = phi-SClusterPhi;
00184          if ( dPhi < -PI )    dPhi = dPhi + 2*PI ;
00185          if ( dPhi >  PI )    dPhi = dPhi - 2*PI ;
00186          if ( fabs(dPhi) >PI )   cout << " error!!! dphi > 2pi   : " << dPhi << endl;
00187          double dR = sqrt(dEta*dEta+dPhi*dPhi);
00188          // Jurassic Cone /////                                                                                                     
00189          if ( dR > r1 ) continue;
00190          if ( dR < r2 ) continue;
00191          if ( fabs(dEta) <  jWidth)  continue;
00192          // stupid bug if ( fabs(dPhi) >  jWidth)  continue;
00194          
00195          if(pt > threshold)
00196             TotalPt = TotalPt + pt;
00197       }
00198 
00199    return TotalPt;
00200 }
00201 
00202 
00203 double TxCalculator::getJct(const reco::Photon cluster, double r1, double r2, double jWidth, double threshold)
00204 {
00205 
00206    using namespace edm;
00207    using namespace reco;
00208 
00209 
00210    double SClusterEta = cluster.eta();
00211    double SClusterPhi = cluster.phi();
00212    double TotalPt = 0;
00213 
00214    for(reco::TrackCollection::const_iterator
00215           recTrack = recCollection->begin(); recTrack!= recCollection->end(); recTrack++)
00216       {
00217          double pt = recTrack->pt();
00218          double eta = recTrack->eta();
00219          double phi = recTrack->phi();
00220          double dEta = fabs(eta-SClusterEta);
00221          double dPhi = phi-SClusterPhi;
00222          if ( dPhi < -PI )    dPhi = dPhi + 2*PI ;
00223          if ( dPhi >  PI )    dPhi = dPhi - 2*PI ;
00224          if ( fabs(dPhi) >PI )   cout << " error!!! dphi > 2pi   : " << dPhi << endl;
00225          //         double dR = sqrt(dEta*dEta+dPhi*dPhi);
00226          
00227 
00229          if ( fabs(dEta) > r1 ) continue;
00230          if ( fabs(dPhi) <r1 ) continue;
00232          
00233          if(pt > threshold)
00234             TotalPt = TotalPt + pt;
00235       }
00236    
00237    double areaStrip = 4*PI*r1 -  4*r1*r1;
00238    double areaJura  = getJurassicArea(r1,r2, jWidth) ;
00239    double theCJ     = getJt(cluster,r1,  r2, jWidth, threshold);
00240 
00241    double theCCJ   = theCJ - TotalPt * areaJura / areaStrip ;
00242    return theCCJ;
00243 
00244 }
00245