CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/RecoLocalCalo/HcalRecAlgos/src/HcalHF_S9S1algorithm.cc

Go to the documentation of this file.
00001 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalHF_S9S1algorithm.h"
00002 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalCaloFlagLabels.h"
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 #include "Geometry/HcalTowerAlgo/src/HcalHardcodeGeometryData.h" // for eta bounds
00005 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
00006 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputerRcd.h"
00007 #include "CondFormats/HcalObjects/interface/HcalChannelQuality.h"
00008 
00009 #include <algorithm> // for "max"
00010 #include <cmath>
00011 #include <iostream>
00012 
00013 HcalHF_S9S1algorithm::HcalHF_S9S1algorithm()
00014 { 
00015   // Default settings:  Energy > 50 GeV, slope = 0, ET = 0
00016   std::vector<double> blank;
00017   blank.clear();
00018   blank.push_back(0);
00019   std::vector<double> EnergyDefault;
00020   EnergyDefault.clear();
00021   EnergyDefault.push_back(50);
00022 
00023   // Thresholds only need to be computed once, not every event!
00024   LongSlopes.clear();
00025   ShortSlopes.clear();
00026   for (int i=29;i<=41;++i)
00027     {
00028       LongSlopes.push_back(0);
00029       ShortSlopes.push_back(0);
00030     }
00031   LongEnergyThreshold.clear();
00032   LongETThreshold.clear();
00033   ShortEnergyThreshold.clear();
00034   ShortETThreshold.clear();
00035   for (int i=29;i<=41;++i)
00036     {
00037       LongEnergyThreshold.push_back(EnergyDefault[0]);
00038       LongETThreshold.push_back(blank[0]);
00039       ShortEnergyThreshold.push_back(EnergyDefault[0]);
00040       ShortETThreshold.push_back(blank[0]);
00041     }
00042   HcalAcceptSeverityLevel_=0;
00043   isS8S1_=false; // S8S1 is almost the same as S9S1
00044 }
00045 
00046 
00047 HcalHF_S9S1algorithm::HcalHF_S9S1algorithm(const std::vector<double>& short_optimumSlope, 
00048                                            const std::vector<double>& short_Energy, 
00049                                            const std::vector<double>& short_ET, 
00050                                            const std::vector<double>& long_optimumSlope, 
00051                                            const std::vector<double>& long_Energy, 
00052                                            const std::vector<double>& long_ET,
00053                                            int HcalAcceptSeverityLevel,
00054                                            bool isS8S1)
00055 
00056 {
00057   // Constructor in the case where all parameters are provided by the user
00058 
00059   // Thresholds only need to be computed once, not every event!
00060 
00061   LongSlopes=long_optimumSlope;
00062   ShortSlopes=short_optimumSlope;
00063   
00064   while (LongSlopes.size()<13)
00065     LongSlopes.push_back(0); // should be unnecessary, but include this protection to avoid crashes
00066   while (ShortSlopes.size()<13)
00067     ShortSlopes.push_back(0);
00068 
00069   // Get long, short energy thresholds (different threshold for each |ieta|)
00070   LongEnergyThreshold.clear();
00071   LongETThreshold.clear();
00072   ShortEnergyThreshold.clear();
00073   ShortETThreshold.clear();
00074   LongEnergyThreshold=long_Energy;
00075   LongETThreshold=long_ET;
00076   ShortEnergyThreshold=short_Energy;
00077   ShortETThreshold=short_ET;
00078 
00079   HcalAcceptSeverityLevel_=HcalAcceptSeverityLevel;
00080   isS8S1_=isS8S1;
00081 } // HcalHF_S9S1algorithm constructor with parameters
00082 
00083 HcalHF_S9S1algorithm::~HcalHF_S9S1algorithm(){}
00084 
00085 
00086 void HcalHF_S9S1algorithm::HFSetFlagFromS9S1(HFRecHit& hf,
00087                                              HFRecHitCollection& rec,
00088                                              const HcalChannelQuality* myqual,
00089                                              const HcalSeverityLevelComputer* mySeverity)
00090 
00091 {
00092   int ieta=hf.id().ieta(); // get coordinates of rechit being checked
00093   int depth=hf.id().depth();
00094   int iphi=hf.id().iphi();
00095   double fEta = 0.5*(theHFEtaBounds[abs(ieta)-29] + theHFEtaBounds[abs(ieta)-28]); // calculate eta as average of eta values at ieta boundaries
00096   double energy=hf.energy();
00097   double ET = energy/fabs(cosh(fEta));
00098 
00099   // Step 1:  Check eta-dependent energy and ET thresholds -- same as PET algorithm
00100   double ETthresh=0, Energythresh=0; // set ET, energy thresholds
00101   if (depth==1)  // set thresholds for long fibers
00102     {
00103       Energythresh = LongEnergyThreshold[abs(ieta)-29];
00104       ETthresh     = LongETThreshold[abs(ieta)-29];
00105     }
00106   else if (depth==2) // short fibers
00107     {
00108       Energythresh = ShortEnergyThreshold[abs(ieta)-29];
00109       ETthresh     = ShortETThreshold[abs(ieta)-29];
00110     }
00111   if (energy<Energythresh || ET < ETthresh)
00112     return;
00113 
00114   // Step 1A:
00115   // Check that EL<ES when evaluating short fibers  (S8S1 check only)
00116   if (depth==2 && abs(ieta)>29 && isS8S1_)
00117     {
00118       double EL=0;
00119       // look for long partner
00120       HcalDetId neighbor(HcalForward, ieta,iphi,1);
00121       HFRecHitCollection::const_iterator neigh=rec.find(neighbor);
00122       if (neigh!=rec.end())
00123         EL=neigh->energy();
00124       
00125       if (EL>=energy)
00126         return;
00127     }
00128 
00129   // Step 2:  Find all neighbors, and calculate S9/S1
00130   double S9S1=0;
00131   int testphi=-99;
00132 
00133   // Part A:  Check fixed iphi, and vary ieta
00134   for (int d=1;d<=2;++d) // depth loop
00135     {
00136       for (int i=ieta-1;i<=ieta+1;++i) // ieta loop
00137         {
00138           testphi=iphi;
00139           // Special case when ieta=39, since ieta=40 only has phi values at 3,7,11,...
00140           // phi=3 covers 3,4,5,6
00141           if (abs(ieta)==39 && abs(i)>39 && testphi%4==1)
00142             testphi-=2;
00143           while (testphi<0) testphi+=72;
00144           if (i==ieta)
00145             if (d==depth || isS8S1_==true) continue;  // don't add the cell itself; don't count neighbor in same ieta-phi if S8S1 test enabled
00146 
00147           // Look to see if neighbor is in rechit collection
00148           HcalDetId neighbor(HcalForward, i,testphi,d);
00149           HFRecHitCollection::const_iterator neigh=rec.find(neighbor);
00150           // require that neighbor exists, and that it doesn't have a prior flag already set
00151           if (neigh!=rec.end())
00152             {
00153               const uint32_t chanstat = myqual->getValues(neighbor)->getValue();
00154               int SeverityLevel=mySeverity->getSeverityLevel(neighbor, neigh->flags(),
00155                                                              chanstat);
00156               if (SeverityLevel<=HcalAcceptSeverityLevel_)
00157                 S9S1+=neigh->energy();
00158             }
00159         }
00160     }
00161 
00162   // Part B: Fix ieta, and loop over iphi.  A bit more tricky, because of iphi wraparound and different segmentation at 40, 41
00163   
00164   int phiseg=2; // 10 degree segmentation for most of HF (1 iphi unit = 5 degrees)
00165   if (abs(ieta)>39) phiseg=4; // 20 degree segmentation for |ieta|>39
00166   for (int d=1;d<=2;++d)
00167     {
00168       for (int i=iphi-phiseg;i<=iphi+phiseg;i+=phiseg)
00169         {
00170           if (i==iphi) continue;  // don't add the cell itself, or its depthwise partner (which is already counted above)
00171           testphi=i;
00172           // Our own modular function, since default produces results -1%72 = -1
00173           while (testphi<0) testphi+=72;
00174           while (testphi>72) testphi-=72;
00175           // Look to see if neighbor is in rechit collection
00176           HcalDetId neighbor(HcalForward, ieta,testphi,d);
00177           HFRecHitCollection::const_iterator neigh=rec.find(neighbor);
00178           if (neigh!=rec.end())
00179             {
00180               const uint32_t chanstat = myqual->getValues(neighbor)->getValue();
00181               int SeverityLevel=mySeverity->getSeverityLevel(neighbor, neigh->flags(),
00182                                                              chanstat);
00183               if (SeverityLevel<=HcalAcceptSeverityLevel_)
00184                 S9S1+=neigh->energy();
00185             }
00186         }
00187     }
00188   
00189   if (abs(ieta)==40) // add extra cells for 39/40 boundary due to increased phi size at ieta=40.
00190     {
00191       for (int d=1;d<=2;++d) // add cells from both depths!
00192         {
00193           HcalDetId neighbor(HcalForward, 39*abs(ieta)/ieta,(iphi+2)%72,d);  
00194           HFRecHitCollection::const_iterator neigh=rec.find(neighbor);
00195           if (neigh!=rec.end())
00196             {
00197               const uint32_t chanstat = myqual->getValues(neighbor)->getValue();
00198               int SeverityLevel=mySeverity->getSeverityLevel(neighbor, neigh->flags(),
00199                                                              chanstat);
00200               if (SeverityLevel<=HcalAcceptSeverityLevel_)
00201                 S9S1+=neigh->energy();
00202             }
00203 
00204         }
00205     }
00206     
00207   // So far, S9S1 is the sum of the neighbors; divide to form ratio
00208   S9S1/=energy;
00209 
00210   // Now compare to threshold
00211   double slope=0;
00212   if (depth==1) slope = LongSlopes[abs(ieta)-29];
00213   else if (depth==2) slope=ShortSlopes[abs(ieta)-29];
00214   double intercept = 0;
00215   if (depth==1) intercept = LongEnergyThreshold[abs(ieta)-29];
00216   else if (depth==2)  intercept = ShortEnergyThreshold[abs(ieta)-29];
00217 
00218   // S9S1 cut has the form [0] + [1]*log[E];  S9S1 value should be above this line
00219   double S9S1cut = 0;
00220   // Protection in case intercept or energy are ever less than 0.  Do we have some other default value of S9S1cut we'd like touse in this case?
00221   if (intercept>0 && energy>0)  
00222     S9S1cut=-1.*slope*log(intercept) + slope*log(energy);
00223   if (S9S1<S9S1cut)
00224     {
00225       // Only set HFS8S1Ratio if S8/S1 ratio test fails
00226       if (isS8S1_==true)
00227         hf.setFlagField(1,HcalCaloFlagLabels::HFS8S1Ratio);
00228       // *Always* set the HFLongShort bit if either S8S1 or S9S1 fail
00229       hf.setFlagField(1,HcalCaloFlagLabels::HFLongShort);
00230     }
00231   return;
00232 } // void HcalHF_S9S1algorithm::HFSetFlagFromS9S1
00233 
00234 
00235 
00236 double HcalHF_S9S1algorithm::CalcSlope(int abs_ieta, const std::vector<double>& params)
00237 {
00238   /* CalcSlope calculates the polynomial [0]+[1]*x + [2]*x^2 + ....,
00239      where x is an integer provided by the first argument (int abs_ieta),
00240      and [0],[1],[2] is a vector of doubles provided by the second (std::vector<double> params).
00241      The output of the polynomial calculation (threshold) is returned by the function.
00242      This function should no longer be needed, since we pass slopes for all ietas into the function via the parameter set.
00243   */
00244   double threshold=0;
00245   for (std::vector<double>::size_type i=0;i<params.size();++i)
00246     {
00247       threshold+=params[i]*pow(static_cast<double>(abs_ieta), (int)i);
00248     }
00249   return threshold;
00250 } // HcalHF_S9S1algorithm::CalcRThreshold(int abs_ieta, std::vector<double> params)
00251   
00252 
00253 
00254 
00255 double HcalHF_S9S1algorithm::CalcEnergyThreshold(double abs_energy,const std::vector<double>& params)
00256 {
00257   /* CalcEnergyThreshold calculates the polynomial [0]+[1]*x + [2]*x^2 + ....,
00258      where x is an integer provided by the first argument (int abs_ieta),
00259      and [0],[1],[2] is a vector of doubles provided by the second (std::vector<double> params).
00260      The output of the polynomial calculation (threshold) is returned by the function.
00261   */
00262   double threshold=0;
00263   for (std::vector<double>::size_type i=0;i<params.size();++i)
00264     {
00265       threshold+=params[i]*pow(abs_energy, (int)i);
00266     }
00267   return threshold;
00268 } //double HcalHF_S9S1algorithm::CalcEnergyThreshold(double abs_energy,std::vector<double> params)
00269