00001 #include "RecoLocalCalo/EcalRecAlgos/interface/EcalSeverityLevelAlgo.h"
00002
00003 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004
00005 int EcalSeverityLevelAlgo::severityLevel( const DetId id,
00006 const EcalRecHitCollection & recHits,
00007 const EcalChannelStatus & chStatus,
00008 float recHitEtThreshold,
00009 SpikeId spId,
00010 float spIdThreshold,
00011 float recHitEnergyThresholdForTiming,
00012 float recHitEnergyThresholdForEE,
00013 float spIdThresholdIEta85
00014 )
00015 {
00016
00017
00018 uint16_t dbStatus = retrieveDBStatus( id, chStatus );
00019
00020 EcalRecHitCollection::const_iterator it = recHits.find( id );
00021 if ( it == recHits.end() ) {
00022
00023
00024 if ( dbStatus >= 10 ) {
00025 return kBad;
00026 } else if ( dbStatus > 0 && dbStatus < 10 ) {
00027
00028 return kProblematic;
00029 } else if ( dbStatus == 0 ) {
00030
00031 return kGood;
00032 }
00033 } else {
00034
00035
00036
00037 if ( id.subdetId() == EcalBarrel ) {
00038 if ( abs(((EBDetId)id).ieta()) == 85 && spId == kSwissCrossBordersIncluded && spikeFromNeighbours(id, recHits, recHitEtThreshold, spId) > spIdThresholdIEta85 ) return kWeird;
00039 if ( spikeFromNeighbours(id, recHits, recHitEtThreshold, spId) > spIdThreshold ) return kWeird;
00040 }
00041
00042 if ( id.subdetId() == EcalBarrel && spikeFromTiming(*it, recHitEnergyThresholdForTiming) ) return kTime;
00043
00044
00045 float re = 0;
00046 if ( id.subdetId() == EcalEndcap && spId >= kSwissCross && (re = recHitE(id, recHits)) > 0 && ( 1-swissCross(id, recHits, recHitEnergyThresholdForEE, spId) < 0.02*log(re/4.) ) ) return kWeird;
00047
00048
00049 return severityLevel( *it, chStatus );
00050 }
00051 return kGood;
00052 }
00053
00054 int EcalSeverityLevelAlgo::severityLevel( const EcalRecHit &recHit,
00055 const EcalChannelStatus &chStatus )
00056 {
00057
00058
00059 uint32_t rhFlag = recHit.recoFlag();
00060 uint16_t dbStatus = retrieveDBStatus( recHit.id(), chStatus );
00061 return severityLevel( rhFlag, dbStatus );
00062 }
00063
00064 int EcalSeverityLevelAlgo::severityLevel( uint32_t rhFlag, uint16_t chStatus )
00065 {
00066
00067 if ( rhFlag == EcalRecHit::kPoorReco
00068 || rhFlag == EcalRecHit::kOutOfTime
00069 || rhFlag == EcalRecHit::kNoisy
00070 || rhFlag == EcalRecHit::kPoorCalib
00071 || rhFlag == EcalRecHit::kFaultyHardware
00072 ) {
00073
00074 return kProblematic;
00075 } else if ( rhFlag == EcalRecHit::kLeadingEdgeRecovered
00076 || rhFlag == EcalRecHit::kNeighboursRecovered
00077 || rhFlag == EcalRecHit::kTowerRecovered
00078 ) {
00079
00080 return kRecovered;
00081 } else if ( rhFlag == EcalRecHit::kDead
00082 || rhFlag == EcalRecHit::kSaturated
00083
00084 || rhFlag == EcalRecHit::kFakeNeighbours
00085 || rhFlag == EcalRecHit::kKilled ) {
00086
00087
00088 return kBad;
00089 }
00090
00091 return kGood;
00092 }
00093
00094 uint16_t EcalSeverityLevelAlgo::retrieveDBStatus( const DetId id, const EcalChannelStatus &chStatus )
00095 {
00096 EcalChannelStatus::const_iterator chIt = chStatus.find( id );
00097 uint16_t dbStatus = 0;
00098 if ( chIt != chStatus.end() ) {
00099 dbStatus = chIt->getStatusCode();
00100 } else {
00101 edm::LogError("EcalSeverityLevelError") << "No channel status found for xtal "
00102 << id.rawId()
00103 << "! something wrong with EcalChannelStatus in your DB? ";
00104 }
00105 return dbStatus;
00106 }
00107
00108 float EcalSeverityLevelAlgo::spikeFromNeighbours( const DetId id,
00109 const EcalRecHitCollection & recHits,
00110 float recHitThreshold,
00111 SpikeId spId
00112 )
00113 {
00114 switch( spId ) {
00115 case kE1OverE9:
00116 return E1OverE9( id, recHits, recHitThreshold );
00117 break;
00118 case kSwissCross:
00119 return swissCross( id, recHits, recHitThreshold , true);
00120 break;
00121 case kSwissCrossBordersIncluded:
00122 return swissCross( id, recHits, recHitThreshold , false);
00123 break;
00124 default:
00125 edm::LogInfo("EcalSeverityLevelAlgo") << "Algorithm number " << spId
00126 << " not known. Please check the enum in EcalSeverityLevelAlgo.h";
00127 break;
00128
00129 }
00130 return 0;
00131 }
00132
00133 float EcalSeverityLevelAlgo::E1OverE9( const DetId id, const EcalRecHitCollection & recHits, float recHitEtThreshold )
00134 {
00135
00136 if ( id.subdetId() == EcalBarrel ) {
00137
00138 if ( recHitApproxEt( id, recHits ) < recHitEtThreshold ) return 0;
00139 EBDetId ebId( id );
00140 float s9 = 0;
00141 for ( int deta = -1; deta <= +1; ++deta ) {
00142 for ( int dphi = -1; dphi <= +1; ++dphi ) {
00143 s9 += recHitE( id, recHits, deta, dphi );
00144 }
00145 }
00146 return recHitE(id, recHits) / s9;
00147 } else if( id.subdetId() == EcalEndcap ) {
00148
00149 if ( recHitApproxEt( id, recHits ) < recHitEtThreshold ) return 0;
00150 EEDetId eeId( id );
00151 float s9 = 0;
00152 for ( int dx = -1; dx <= +1; ++dx ) {
00153 for ( int dy = -1; dy <= +1; ++dy ) {
00154 s9 += recHitE( id, recHits, dx, dy );
00155 }
00156 }
00157 return recHitE(id, recHits) / s9;
00158
00159 }
00160 return 0;
00161 }
00162
00163 float EcalSeverityLevelAlgo::swissCross( const DetId id, const EcalRecHitCollection & recHits, float recHitThreshold , bool avoidIeta85)
00164 {
00165
00166 if ( id.subdetId() == EcalBarrel ) {
00167 EBDetId ebId( id );
00168
00169
00170
00171
00172 if ( abs(ebId.ieta())==85 && avoidIeta85) return 0;
00173
00174 if ( recHitApproxEt( id, recHits ) < recHitThreshold ) return 0;
00175 float s4 = 0;
00176 float e1 = recHitE( id, recHits );
00177
00178 if ( e1 == 0 ) return 0;
00179 s4 += recHitE( id, recHits, 1, 0 );
00180 s4 += recHitE( id, recHits, -1, 0 );
00181 s4 += recHitE( id, recHits, 0, 1 );
00182 s4 += recHitE( id, recHits, 0, -1 );
00183 return 1 - s4 / e1;
00184 } else if ( id.subdetId() == EcalEndcap ) {
00185 EEDetId eeId( id );
00186
00187 float e1 = recHitE( id, recHits );
00188 if ( e1 < recHitThreshold ) return 0;
00189 float s4 = 0;
00190
00191 if ( e1 == 0 ) return 0;
00192 s4 += recHitE( id, recHits, 1, 0 );
00193 s4 += recHitE( id, recHits, -1, 0 );
00194 s4 += recHitE( id, recHits, 0, 1 );
00195 s4 += recHitE( id, recHits, 0, -1 );
00196 return 1 - s4 / e1;
00197 }
00198 return 0;
00199 }
00200
00201 float EcalSeverityLevelAlgo::recHitE( const DetId id, const EcalRecHitCollection & recHits,
00202 int di, int dj )
00203 {
00204
00205
00206
00207 DetId nid;
00208 if( id.subdetId() == EcalBarrel) nid = EBDetId::offsetBy( id, di, dj );
00209 else if( id.subdetId() == EcalEndcap) nid = EEDetId::offsetBy( id, di, dj );
00210
00211 return ( nid == DetId(0) ? 0 : recHitE( nid, recHits ) );
00212 }
00213
00214 float EcalSeverityLevelAlgo::recHitE( const DetId id, const EcalRecHitCollection &recHits )
00215 {
00216 if ( id == DetId(0) ) {
00217 return 0;
00218 } else {
00219 EcalRecHitCollection::const_iterator it = recHits.find( id );
00220 if ( it != recHits.end() ) return (*it).energy();
00221 }
00222 return 0;
00223 }
00224
00225
00226 float EcalSeverityLevelAlgo::recHitApproxEt( const DetId id, const EcalRecHitCollection &recHits )
00227 {
00228
00229 if ( id.subdetId() == EcalBarrel ) {
00230 return recHitE( id, recHits ) / cosh( EBDetId::approxEta( id ) );
00231 }
00232 return 0;
00233 }
00234
00235
00236 bool EcalSeverityLevelAlgo::spikeFromTiming( const EcalRecHit &recHit, float recHitEnergyThreshold)
00237 {
00238 if ( recHit.energy() < recHitEnergyThreshold ) return false;
00239 if ( recHit.recoFlag() == EcalRecHit::kOutOfTime ) return true;
00240 return false;
00241 }
00242
00243
00244
00245 float EcalSeverityLevelAlgo::E2overE9( const DetId id, const EcalRecHitCollection & recHits,
00246 float recHitEtThreshold, float recHitEtThreshold2 ,
00247 bool avoidIeta85, bool KillSecondHit)
00248 {
00249
00250
00251
00252
00253
00254
00255
00256
00257
00258
00259
00260
00261
00262
00263
00264
00265
00266
00267
00268
00269 if ( id.subdetId() == EcalBarrel ) {
00270
00271 EBDetId ebId( id );
00272
00273
00274 if ( abs(ebId.ieta())==85 && avoidIeta85) return 0;
00275
00276
00277
00278
00279 float e1 = recHitE( id, recHits );
00280 float ete1=recHitApproxEt( id, recHits );
00281
00282
00283
00284
00285 if (ete1 < std::min(recHitEtThreshold,recHitEtThreshold2) ) return 0;
00286
00287 if (ete1 < recHitEtThreshold && !KillSecondHit ) return 0;
00288
00289
00290 float e2=-1;
00291 float ete2=0;
00292 float s9 = 0;
00293
00294
00295 int e2eta=0;
00296 int e2phi=0;
00297
00298
00299
00300 for ( int deta = -1; deta <= +1; ++deta ) {
00301 for ( int dphi = -1; dphi <= +1; ++dphi ) {
00302
00303
00304
00305 float etmp=recHitE( id, recHits, deta, dphi );
00306 s9 += etmp;
00307
00308 EBDetId idtmp=EBDetId::offsetBy(id,deta,dphi);
00309 float eapproxet=recHitApproxEt( idtmp, recHits );
00310
00311
00312 if (etmp>e2 && eapproxet>recHitEtThreshold2 && !(deta==0 && dphi==0)) {
00313
00314 e2=etmp;
00315 ete2=eapproxet;
00316 e2eta=deta;
00317 e2phi=dphi;
00318
00319 }
00320
00321 }
00322 }
00323
00324 if ( e1 == 0 ) return 0;
00325
00326
00327 if ( e2 == -1 ) return 0;
00328
00329
00330
00331 float e2nd=e1+e2;
00332 float e2e9=0;
00333
00334 if (s9!=0) e2e9=e2nd/s9;
00335
00336
00337
00338
00339 if (e1 > e2 && ete1>recHitEtThreshold) return e2e9;
00340
00341
00342
00343 if ( e2 > e1 ) {
00344
00345
00346
00347
00348
00349
00350 if (!KillSecondHit || ete2<recHitEtThreshold || ete1<recHitEtThreshold2) {
00351
00352 return 0;
00353
00354 }
00355
00356
00357 else {
00358
00359
00360
00361 float s92nd=0;
00362
00363 float e2nd_prime=0;
00364 int e2prime_eta=0;
00365 int e2prime_phi=0;
00366
00367 EBDetId secondid=EBDetId::offsetBy(id,e2eta,e2phi);
00368
00369
00370 for ( int deta = -1; deta <= +1; ++deta ) {
00371 for ( int dphi = -1; dphi <= +1; ++dphi ) {
00372
00373
00374
00375 float etmp=recHitE( secondid, recHits, deta, dphi );
00376 s92nd += etmp;
00377
00378 if (etmp>e2nd_prime && !(deta==0 && dphi==0)) {
00379 e2nd_prime=etmp;
00380 e2prime_eta=deta;
00381 e2prime_phi=dphi;
00382 }
00383
00384 }
00385 }
00386
00387
00388
00389 if (!(e2prime_eta==-e2eta && e2prime_phi==-e2phi))
00390 {
00391 return 0;
00392 }
00393
00394
00395
00396 float e2e9_2=0;
00397 if (s92nd!=0) e2e9_2=e2nd/s92nd;
00398
00399
00400
00401 return e2e9_2;
00402
00403
00404 }
00405
00406 }
00407
00408
00409 } else if ( id.subdetId() == EcalEndcap ) {
00410
00411 return 0;
00412 }
00413 return 0;
00414 }
00415