Go to the documentation of this file.00001
00002
00003
00004
00005 #include "RecoHI/HiEgammaAlgos/interface/HiEgammaSCEnergyCorrectionAlgo.h"
00006 #include "RecoEcal/EgammaCoreTools/interface/SuperClusterShapeAlgo.h"
00007 #include <iostream>
00008 #include <string>
00009 #include <vector>
00010
00011 HiEgammaSCEnergyCorrectionAlgo::HiEgammaSCEnergyCorrectionAlgo(float noise,
00012 reco::CaloCluster::AlgoId theAlgo,
00013 const edm::ParameterSet& pSet,
00014 HiEgammaSCEnergyCorrectionAlgo::VerbosityLevel verbosity
00015 )
00016 {
00017 sigmaElectronicNoise_ = noise;
00018 verbosity_ = verbosity;
00019
00020
00021 p_fBrem_ = pSet.getParameter< std::vector <double> > ("fBremVect");
00022 p_fBremTh_ = pSet.getParameter< std::vector <double> > ("fBremThVect");
00023 p_fEta_ = pSet.getParameter< std::vector <double> > ("fEtaVect");
00024 p_fEtEta_ = pSet.getParameter< std::vector <double> > ("fEtEtaVect");
00025
00026
00027 minR9Barrel_ = pSet.getParameter<double>("minR9Barrel");
00028 minR9Endcap_ = pSet.getParameter<double>("minR9Endcap");
00029
00030
00031 maxR9_ = pSet.getParameter<double>("maxR9");
00032
00033 }
00034
00035
00036
00037 reco::SuperCluster
00038 HiEgammaSCEnergyCorrectionAlgo::applyCorrection(const reco::SuperCluster &cl,
00039 const EcalRecHitCollection &rhc, reco::CaloCluster::AlgoId theAlgo,
00040 const CaloSubdetectorGeometry* geometry,
00041 const CaloTopology *topology, EcalClusterFunctionBaseClass* EnergyCorrection)
00042 {
00043
00044
00045 if (verbosity_ <= pINFO)
00046 {
00047 std::cout << " HiEgammaSCEnergyCorrectionAlgo::applyCorrection" << std::endl;
00048 std::cout << " SC has energy " << cl.energy() << std::endl;
00049 std::cout << " Will correct now.... " << std::endl;
00050 }
00051
00052
00053 reco::CaloClusterPtr seedC = cl.seed();
00054
00055 if (verbosity_ <= pINFO)
00056 {
00057 std::cout << " Seed cluster energy... " << seedC->energy() << std::endl;
00058 }
00059
00060
00061 reco::CaloClusterPtrVector clusters_v;
00062
00063 if (verbosity_ <= pINFO) std::cout << " Constituent cluster energies... ";
00064
00065 for(reco::CaloCluster_iterator cluster = cl.clustersBegin(); cluster != cl.clustersEnd(); cluster ++)
00066 {
00067 clusters_v.push_back(*cluster);
00068 if (verbosity_ <= pINFO) std::cout << (*cluster)->energy() << ", ";
00069 }
00070 if (verbosity_ <= pINFO) std::cout << std::endl;
00071
00072
00073 if (verbosity_ <= pINFO)
00074 {
00075 std::cout << " The seed cluster used algo " << theAlgo;
00076 }
00077
00078
00079
00080 std::vector<std::pair<DetId, float> > const & seedHits = seedC->hitsAndFractions();
00081 EcalSubdetector theBase = EcalSubdetector(seedHits.at(0).first.subdetId());
00082 if (verbosity_ <= pINFO)
00083 {
00084 std::cout << " seed cluster location == " << theBase << std::endl;
00085 }
00086
00087
00088 int nCryGT2Sigma = nCrystalsGT2Sigma(*seedC,rhc);
00089 if (verbosity_ <= pINFO)
00090 {
00091 std::cout << " nCryGT2Sigma " << nCryGT2Sigma << std::endl;
00092 }
00093
00094
00095 float bremsEnergy = cl.energy() - seedC->energy();
00096 if (verbosity_ <= pINFO)
00097 {
00098 std::cout << " bremsEnergy " << bremsEnergy << std::endl;
00099 }
00100
00101
00102
00103 SuperClusterShapeAlgo SCShape(&rhc, geometry);
00104
00105 double phiWidth = 0.;
00106 double etaWidth = 0.;
00107
00108
00109 SCShape.Calculate_Covariances(cl);
00110 phiWidth = SCShape.phiWidth();
00111 etaWidth = SCShape.etaWidth();
00112
00113
00114 float e3x3 = EcalClusterTools::e3x3( *(cl.seed()), &rhc, &(*topology));
00115 float e5x5 = EcalClusterTools::e5x5( *(cl.seed()), &rhc, &(*topology));
00116 float r9 = e3x3 / cl.rawEnergy();
00117
00118
00119
00120
00121 float newEnergy = 0;
00122
00123
00124 if ((r9 < minR9Barrel_&&theBase == EcalBarrel) || (r9 < minR9Endcap_&&theBase == EcalEndcap)) {
00125
00126 newEnergy = (cl.rawEnergy())/ fEta(cl.eta(), theAlgo, theBase) /
00127 fBrem(phiWidth/etaWidth, theAlgo, theBase)
00128 /fEtEta(cl.energy()/cosh(cl.eta()), cl.eta(), theAlgo, theBase);
00129
00130 } else {
00131 if (r9 < maxR9_) {
00132
00133 newEnergy = e5x5 / fEta(cl.eta(), theAlgo, theBase);
00134 } else {
00135
00136 newEnergy = cl.rawEnergy();
00137 }
00138 }
00139
00140
00141 if (verbosity_ <= pINFO)
00142 {
00143 std::cout << " UNCORRECTED SC has energy... " << cl.energy() << std::endl;
00144 std::cout << " CORRECTED SC has energy... " << newEnergy << std::endl;
00145 std::cout << " Size..." <<cl.size() << std::endl;
00146 std::cout << " Seed nCryGT2Sigma Size..." <<nCryGT2Sigma << std::endl;
00147 }
00148
00149 reco::SuperCluster corrCl(newEnergy,
00150 math::XYZPoint(cl.position().X(), cl.position().Y(), cl.position().Z()),
00151 cl.seed(), clusters_v, cl.preshowerEnergy());
00152
00153
00154 corrCl.setFlags(cl.flags());
00155 corrCl.setPhiWidth(phiWidth);
00156 corrCl.setEtaWidth(etaWidth);
00157
00158 return corrCl;
00159 }
00160
00161 float HiEgammaSCEnergyCorrectionAlgo::fEtEta(float et, float eta, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
00162 {
00163 int offset = 0;
00164 float factor;
00165 if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::island)) {
00166 offset = 0;
00167 } else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::island)) {
00168 offset = 7;
00169 }
00170
00171
00172 factor = (p_fEtEta_[0+offset] + p_fEtEta_[1+offset]*sqrt(et));
00173
00174 factor *= (p_fEtEta_[2+offset] + p_fEtEta_[3+offset]*fabs(eta) +
00175 p_fEtEta_[4+offset]*eta*eta + p_fEtEta_[5+offset]*eta*eta*fabs(eta) + + p_fEtEta_[6+offset]*eta*eta*eta*eta);
00176
00177
00178 if (factor< 0.66f ) factor = 0.66f;
00179 if (factor> 1.5f ) factor = 1.5f;
00180
00181 return factor;
00182
00183 }
00184
00185 float HiEgammaSCEnergyCorrectionAlgo::fEta(float eta, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
00186 {
00187 int offset = 0;
00188 float factor;
00189 if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::island)) {
00190 offset = 0;
00191 } else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::island)) {
00192 offset = 3;
00193 }
00194
00195 factor = (p_fEta_[0+offset] + p_fEta_[1+offset]*fabs(eta) + p_fEta_[2+offset]*eta*eta);
00196
00197
00198 if (factor< 0.66f ) factor = 0.66f;
00199 if (factor> 1.5f ) factor = 1.5f;
00200
00201 return factor;
00202 }
00203
00204 float HiEgammaSCEnergyCorrectionAlgo::fBrem(float brem, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
00205 {
00206 int det = 0;
00207 int offset = 0;
00208 float factor;
00209 if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::island)) {
00210 det = 0;
00211 offset = 0;
00212 } else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::island)) {
00213 det = 1;
00214 offset = 6;
00215 }
00216
00217 if (brem < p_fBremTh_[det]) {
00218 factor = p_fBrem_[0+offset] + p_fBrem_[1+offset]*brem + p_fBrem_[2+offset]*brem*brem;
00219 } else {
00220 factor = p_fBrem_[3+offset] + p_fBrem_[4+offset]*brem + p_fBrem_[5+offset]*brem*brem;
00221 };
00222
00223
00224 if (factor< 0.66f ) factor = 0.66f;
00225 if (factor> 1.5f ) factor = 1.5f;
00226
00227 return factor;
00228 }
00229
00230
00231
00232 float HiEgammaSCEnergyCorrectionAlgo::fNCrystals(int nCry, reco::CaloCluster::AlgoId theAlgo, EcalSubdetector theBase) const
00233 {
00234
00235 float x = (float) nCry;
00236 float result =1.f;
00237
00238 if((theBase == EcalBarrel) && (theAlgo == reco::CaloCluster::island)) {
00239 float const p0 = 0.682554f;
00240 float const p1 = 0.0253013f;
00241 float const p2 = -0.0007907f;
00242 float const p3 = 1.166e-5f;
00243 float const p4 = -6.7387e-8f;
00244 if (x < 10.f) x = 10.f;
00245 if (x < 40.f) result = p0 + x*(p1 + x*(p2 + x*(p3 + x*p4))); else result = 1.f;
00246 }
00247
00248 else if((theBase == EcalEndcap) && (theAlgo == reco::CaloCluster::island)) {
00249
00250 float const p0 = 0.712185f;
00251 float const p1 = 0.0273609f;
00252 float const p2 = -0.00103818f;
00253 float const p3 = 2.01828e-05f;
00254 float const p4 = -1.71438e-07f;
00255 if (x < 10.f) x = 10.f;
00256 if (x < 40.f) result = p0 + x*(p1 + x*(p2 + x*(p3 + x*p4))); else result = 1.f;
00257 }
00258
00259 else {
00260 if (verbosity_ <= pINFO)
00261 {
00262 std::cout << "trying to correct unknown cluster!!!" << std::endl;
00263 }
00264 }
00265
00266 if (result > 1.5f) result = 1.5f;
00267 if (result < 0.5f) result = 0.5f;
00268
00269 return result;
00270 }
00271
00272 int HiEgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma(reco::BasicCluster const & seed, EcalRecHitCollection const &rhc) const {
00273
00274
00275
00276 std::vector<std::pair<DetId,float > > const & hits = seed.hitsAndFractions();
00277
00278 if (verbosity_ <= pINFO)
00279 {
00280 std::cout << " HiEgammaSCEnergyCorrectionAlgo::nCrystalsGT2Sigma" << std::endl;
00281 std::cout << " Will calculate number of crystals above 2sigma noise" << std::endl;
00282 std::cout << " sigmaElectronicNoise = " << sigmaElectronicNoise_ << std::endl;
00283 std::cout << " There are " << hits.size() << " recHits" << std::endl;
00284 }
00285
00286 int nCry = 0;
00287 for(std::vector<std::pair<DetId,float > >::const_iterator hit = hits.begin(); hit != hits.end(); ++hit)
00288 {
00289
00290 EcalRecHitCollection::const_iterator aHit = rhc.find((*hit).first);
00291
00292 if(aHit->energy()>2.f*sigmaElectronicNoise_) ++nCry;
00293 }
00294
00295 if (verbosity_ <= pINFO)
00296 {
00297 std::cout << " " << nCry << " of these above 2sigma noise" << std::endl;
00298 }
00299
00300 return nCry;
00301 }
00302