CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/RecoLocalCalo/HcalRecAlgos/src/HcalHF_PETalgorithm.cc

Go to the documentation of this file.
00001 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalHF_PETalgorithm.h"
00002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00003 #include "Geometry/HcalTowerAlgo/src/HcalHardcodeGeometryData.h" // for eta bounds
00004 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalCaloFlagLabels.h"
00005 #include "CondFormats/HcalObjects/interface/HcalChannelQuality.h"
00006 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputer.h"
00007 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalSeverityLevelComputerRcd.h"
00008 
00009 #include <algorithm> // for "max"
00010 #include <cmath>
00011 #include <iostream>
00012 
00013 HcalHF_PETalgorithm::HcalHF_PETalgorithm()
00014 {
00015   // no params given; revert to old 'algo 1' fixed settings
00016   short_R.clear();
00017   short_R.push_back(0.995);
00018   short_R_29.clear();
00019   short_R_29.push_back(0.995);
00020   short_Energy_Thresh.clear();
00021   short_Energy_Thresh.push_back(50);
00022   short_ET_Thresh.clear();
00023   short_ET_Thresh.push_back(0);
00024   
00025   long_R.clear();
00026   long_R.push_back(0.995);
00027   long_R_29.clear();
00028   long_R_29.push_back(0.995);
00029   long_Energy_Thresh.clear();
00030   long_Energy_Thresh.push_back(50);
00031   long_ET_Thresh.clear();
00032   long_ET_Thresh.push_back(0);
00033   flagsToSkip_=0;
00034 }
00035 
00036 HcalHF_PETalgorithm::HcalHF_PETalgorithm(std::vector<double> shortR, 
00037                                          std::vector<double> shortEnergyParams, 
00038                                          std::vector<double> shortETParams, 
00039                                          std::vector<double> longR, 
00040                                          std::vector<double> longEnergyParams, 
00041                                          std::vector<double> longETParams,
00042                                          int flagsToSkip,
00043                                          std::vector<double> shortR29,
00044                                          std::vector<double> longR29)
00045 {
00046   // R is parameterized depending on the energy of the cells, so just store the parameters here
00047   short_R    = shortR;
00048   long_R     = longR;
00049   short_R_29 = shortR29;
00050   long_R_29  = longR29;
00051 
00052   // Energy and ET cuts are ieta-dependent, and only need to be calculated once!
00053   short_Energy_Thresh.clear();
00054   short_ET_Thresh.clear();
00055   long_Energy_Thresh.clear();
00056   long_ET_Thresh.clear();
00057 
00058   //separate short, long cuts provided for each ieta 
00059   short_Energy_Thresh=shortEnergyParams;
00060   short_ET_Thresh=shortETParams;
00061   long_Energy_Thresh=longEnergyParams;
00062   long_ET_Thresh=longETParams;
00063 
00064   flagsToSkip_=flagsToSkip;
00065 }
00066 
00067 HcalHF_PETalgorithm::~HcalHF_PETalgorithm(){}
00068 
00069 
00070 
00071 void HcalHF_PETalgorithm::HFSetFlagFromPET(HFRecHit& hf,
00072                                            HFRecHitCollection& rec,
00073                                            HcalChannelQuality* myqual,
00074                                            const HcalSeverityLevelComputer* mySeverity)
00075 {
00076   /*  Set the HFLongShort flag by comparing the ratio |L-S|/|L+S|.  Channels must first pass energy and ET cuts, and channels whose partners are known to be dead are skipped, since those channels can never satisfy the ratio cut.  */
00077 
00078   int ieta=hf.id().ieta(); // get coordinates of rechit being checked
00079   int depth=hf.id().depth();
00080   int iphi=hf.id().iphi();
00081   double fEta = 0.5*(theHFEtaBounds[abs(ieta)-29] + theHFEtaBounds[abs(ieta)-28]); // calculate eta as average of eta values at ieta boundaries
00082   double energy=hf.energy();
00083   double ET = energy/fabs(cosh(fEta));
00084 
00085   // Step 1:  Check energy and ET cuts
00086   double ETthresh=0, Energythresh=0; // set ET, energy thresholds
00087 
00088   if (depth==1)  // set thresholds for long fibers
00089     {
00090       Energythresh = long_Energy_Thresh[abs(ieta)-29];
00091       ETthresh     = long_ET_Thresh[abs(ieta)-29];
00092     }
00093   else if (depth==2) // short fibers
00094     {
00095       Energythresh = short_Energy_Thresh[abs(ieta)-29];
00096       ETthresh     = short_ET_Thresh[abs(ieta)-29];
00097     }
00098   
00099   if (energy<Energythresh || ET < ETthresh)
00100     {
00101       hf.setFlagField(0, HcalCaloFlagLabels::HFLongShort); // shouldn't be necessary, but set bit to 0 just to be sure
00102       hf.setFlagField(0, HcalCaloFlagLabels::HFPET);
00103       return;
00104     }
00105 
00106   // Step 2:  Get partner info, check if partner is excluded from rechits already
00107   HcalDetId partner(HcalForward, ieta, iphi, 3-depth); //  if depth=1, 3-depth=2, and vice versa
00108   DetId detpartner=DetId(partner);
00109   const HcalChannelStatus* partnerstatus=myqual->getValues(detpartner.rawId());
00110 
00111   // Don't set the bit if the partner has been dropped from the rechit collection ('dead' or 'off')
00112   if (mySeverity->dropChannel(partnerstatus->getValue() ) )
00113     {
00114       hf.setFlagField(0, HcalCaloFlagLabels::HFLongShort); // shouldn't be necessary, but set bit to 0 just to be sure
00115       hf.setFlagField(0, HcalCaloFlagLabels::HFPET);
00116       return;
00117     }
00118 
00119   // Step 3:  Compute ratio
00120   double Ratio=1;
00121   HFRecHitCollection::const_iterator part=rec.find(partner);
00122   if (part!=rec.end())
00123     {
00124       if ((flagsToSkip_&part->flags())==0)
00125         Ratio=(energy-part->energy())/(energy+part->energy());
00126     }
00127 
00128   double RatioThresh=0;
00129   // Allow for the ratio cut to be parameterized in terms of energy
00130   if (abs(ieta)==29)
00131     {
00132       if (depth==1) RatioThresh=CalcThreshold(energy,long_R_29);
00133       else if (depth==2) RatioThresh=CalcThreshold(energy,short_R_29);
00134     }
00135   else
00136     {
00137       if (depth==1) RatioThresh=CalcThreshold(energy,long_R);
00138       else if (depth==2) RatioThresh=CalcThreshold(energy,short_R);
00139     }
00140   if (Ratio<=RatioThresh)
00141     {
00142       hf.setFlagField(0, HcalCaloFlagLabels::HFLongShort); // shouldn't be necessary, but set bit to 0 just to be sure
00143       hf.setFlagField(0, HcalCaloFlagLabels::HFPET);
00144       return;
00145     }
00146   // Made it this far -- ratio is > threshold, and cell should be flagged!
00147   // 'HFLongShort' is overall topological flag, and 'HFPET' is the flag for this
00148   // specific test
00149   hf.setFlagField(1, HcalCaloFlagLabels::HFLongShort);
00150   hf.setFlagField(1, HcalCaloFlagLabels::HFPET);
00151 }//void HcalHF_PETalgorithm::HFSetFlagFromPET
00152 
00153 
00154 
00155 double HcalHF_PETalgorithm::CalcThreshold(double abs_x,std::vector<double> params)
00156 {
00157   /* CalcEnergyThreshold calculates the polynomial [0]+[1]*x + [2]*x^2 + ....,
00158      where x is an integer provided by the first argument (int double abs_x),
00159      and [0],[1],[2] is a vector of doubles provided by the second (std::vector<double> params).
00160      The output of the polynomial calculation (threshold) is returned by the function.
00161   */
00162   double threshold=0;
00163   for (std::vector<double>::size_type i=0;i<params.size();++i)
00164     {
00165       threshold+=params[i]*pow(abs_x, (int)i);
00166     }
00167   return threshold;
00168 } //double HcalHF_PETalgorithm::CalcThreshold(double abs_x,std::vector<double> params)
00169