00001 #include "RecoEgamma/EgammaHFProducers/interface/HFClusterAlgo.h"
00002 #include <sstream>
00003 #include <iostream>
00004 #include <algorithm>
00005 #include <list>
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 #include "RecoLocalCalo/HcalRecAlgos/interface/HcalCaloFlagLabels.h"
00008
00009 using namespace std;
00010 using namespace reco;
00019 HFClusterAlgo::HFClusterAlgo()
00020 {
00021 m_isMC=true;
00022
00023 }
00024
00025 class CompareHFCompleteHitET {
00026 public:
00027 bool operator()(const HFClusterAlgo::HFCompleteHit& h1,const HFClusterAlgo::HFCompleteHit& h2) const {
00028 return h1.et>h2.et;
00029 }
00030 };
00031
00032 class CompareHFCore {
00033 public:
00034 bool operator()(const double c1,const double c2) const {
00035 return c1>c2;
00036 }
00037 };
00038
00039 static int indexByEta(HcalDetId id) {
00040 return (id.zside()>0)?(id.ietaAbs()-29+13):(41-id.ietaAbs());
00041 }
00042
00043 static const double MCMaterialCorrections[] = { 1.000,1.000,1.105,0.970,0.965,0.975,0.956,0.958,0.981,1.005,0.986,1.086,1.000,
00044 1.000,1.086,0.986,1.005,0.981,0.958,0.956,0.975,0.965,0.970,1.105,1.000,1.000 };
00045
00046 static const double ZplusMC2010Corrections[] = { 1.0, 1.0, 1.075, 0.902, 0.983, 0.992, 0.948, 0.941, 0.975, 0.990, 0.988, 1.086, 1.0,
00047 1.0, 1.126, 0.959, 0.979, 0.989, 0.959, 0.901, 0.937, 0.962, 0.967, 1.006, 1.0, 1.0};
00048
00049
00050
00051 void HFClusterAlgo::setup(double minTowerEnergy, double seedThreshold,double maximumSL,double maximumRenergy,
00052 bool usePMTflag,bool usePulseflag, bool forcePulseFlagMC, int correctionSet){
00053 m_seedThreshold=seedThreshold;
00054 m_minTowerEnergy=minTowerEnergy;
00055 m_maximumSL=maximumSL;
00056 m_usePMTFlag=usePMTflag;
00057 m_usePulseFlag=usePulseflag;
00058 m_forcePulseFlagMC=forcePulseFlagMC;
00059 m_maximumRenergy=maximumRenergy;
00060 m_correctionSet = correctionSet;
00061
00062 for(int ii=0;ii<13;ii++){
00063 m_cutByEta.push_back(-1);
00064 }
00065
00066
00067 for (int ii=0; ii<13*2; ii++)
00068 m_correctionByEta.push_back(1.0);
00069
00070 if (m_correctionSet==1) {
00071 for (int ii=0; ii<13*2; ii++)
00072 m_correctionByEta[ii]=MCMaterialCorrections[ii];
00073 }
00074 if (m_correctionSet==2) {
00075 for (int ii=0; ii<13*2; ii++)
00076 m_correctionByEta[ii]=ZplusMC2010Corrections[ii];
00077 }
00078
00079 }
00080
00082 void HFClusterAlgo::clusterize(const HFRecHitCollection& hf,
00083 const CaloGeometry& geom,
00084 HFEMClusterShapeCollection& clusterShapes,
00085 SuperClusterCollection& SuperClusters) {
00086
00087 std::vector<HFCompleteHit> protoseeds, seeds;
00088 HFRecHitCollection::const_iterator j,j2;
00089 std::vector<HFCompleteHit>::iterator i;
00090 std::vector<HFCompleteHit>::iterator k;
00091 int dP, dE, PWrap;
00092 bool isok=true;
00093 HFEMClusterShape clusShp;
00094
00095 SuperCluster Sclus;
00096 bool doCluster=false;
00097
00098 for (j=hf.begin(); j!= hf.end(); j++) {
00099 const int aieta=j->id().ietaAbs();
00100 int iz=(aieta-29);
00101
00102 if (j->id().depth()!=1) continue;
00103 if (aieta==40 || aieta==41 || aieta==29) continue;
00104
00105
00106 if (iz<0 || iz>12) {
00107 edm::LogWarning("HFClusterAlgo") << "Strange invalid HF hit: " << j->id();
00108 continue;
00109 }
00110
00111 if (m_cutByEta[iz]<0) {
00112 double eta=geom.getPosition(j->id()).eta();
00113 m_cutByEta[iz]=m_seedThreshold*cosh(eta);
00114 }
00115 double elong=j->energy()*m_correctionByEta[indexByEta(j->id())];
00116 if (elong>m_cutByEta[iz]) {
00117 j2=hf.find(HcalDetId(HcalForward,j->id().ieta(),j->id().iphi(),2));
00118 double eshort=(j2==hf.end())?(0):(j2->energy());
00119 if (j2!=hf.end())
00120 eshort*=m_correctionByEta[indexByEta(j2->id())];
00121 if (((elong-eshort)/(elong+eshort))>m_maximumSL) continue;
00122
00123
00124 if(isPMTHit(*j)) continue;
00125
00126 HFCompleteHit ahit;
00127 double eta=geom.getPosition(j->id()).eta();
00128 ahit.id=j->id();
00129 ahit.energy=elong;
00130 ahit.et=ahit.energy/cosh(eta);
00131 protoseeds.push_back(ahit);
00132 }
00133 }
00134
00135 if(!protoseeds.empty()){
00136 std::sort(protoseeds.begin(), protoseeds.end(), CompareHFCompleteHitET());
00137 for (i=protoseeds.begin(); i!= protoseeds.end(); i++) {
00138 isok=true;
00139 doCluster=false;
00140
00141 if ( (i==protoseeds.begin()) && (isok) ) {
00142 doCluster=true;
00143 }else {
00144
00145 for (k=seeds.begin(); isok && k!=seeds.end(); k++) {
00146
00147 for (dE=-2; dE<=2; dE++)
00148 for (dP=-4;dP<=4; dP+=2) {
00149 PWrap=k->id.iphi()+dP;
00150 if (PWrap<0)
00151 PWrap+=72;
00152 if (PWrap>72)
00153 PWrap-=72;
00154
00155 if ( (i->id.iphi()==PWrap) && (i->id.ieta()==k->id.ieta()+dE))
00156 isok = false;
00157 }
00158 }
00159 if (isok) {
00160 doCluster=true;
00161 }
00162 }
00163 if (doCluster) {
00164 seeds.push_back(*i);
00165
00166 bool clusterOk=makeCluster( i->id(),hf, geom,clusShp,Sclus);
00167 if (clusterOk) {
00168 clusterShapes.push_back(clusShp);
00169 SuperClusters.push_back(Sclus);
00170 }
00171
00172 }
00173 }
00174 }
00175 }
00176
00177
00178 bool HFClusterAlgo::makeCluster(const HcalDetId& seedid,
00179 const HFRecHitCollection& hf,
00180 const CaloGeometry& geom,
00181 HFEMClusterShape& clusShp,
00182 SuperCluster& Sclus) {
00183
00184
00185 double w=0;
00186 double wgt=0;
00187 double w_e=0;
00188 double w_x=0;
00189 double w_y=0;
00190 double w_z=0;
00191 double wp_e=0;
00192 double e_e=0;
00193 double e_ep=0;
00194
00195 double l_3=0;
00196 double s_3=0;
00197 double l_5=0;
00198 double s_5=0;
00199 double l_1=0;
00200 double s_1=0;
00201 int de, dp, phiWrap;
00202 double l_1e=0;
00203 GlobalPoint sp=geom.getPosition(seedid);
00204 std::vector<double> coreCanid;
00205 std::vector<double>::const_iterator ci;
00206 HFRecHitCollection::const_iterator i,is,il;
00207 std::vector<DetId> usedHits;
00208
00209 HFRecHitCollection::const_iterator si;
00210 HcalDetId sid(HcalForward,seedid.ieta(),seedid.iphi(),1);
00211 si=hf.find(sid);
00212
00213 bool clusterOk=true;
00214
00215
00216
00217 bool edge_type1=seedid.ietaAbs()==39 && (seedid.iphi()%4)==3;
00218
00219 double e_seed=si->energy()*m_correctionByEta[indexByEta(si->id())];
00220
00221 for (de=-2; de<=2; de++)
00222 for (dp=-4;dp<=4; dp+=2) {
00223 phiWrap=seedid.iphi()+dp;
00224 if (phiWrap<0)
00225 phiWrap+=72;
00226 if (phiWrap>72)
00227 phiWrap-=72;
00228
00229
00230
00231 if (edge_type1 && de==seedid.zside()) {
00232 if (dp==-2) {
00233 phiWrap-=2;
00234 if (phiWrap<0)
00235 phiWrap+=72;
00236 }
00237 else if (dp==-4) {
00238 continue;
00239 }
00240 }
00241
00242 HcalDetId idl(HcalForward,seedid.ieta()+de,phiWrap,1);
00243 HcalDetId ids(HcalForward,seedid.ieta()+de,phiWrap,2);
00244
00245
00246 il=hf.find(idl);
00247 is=hf.find(ids);
00248
00249
00250
00251
00252 double e_long=1.0;
00253 double e_short=0.0;
00254 if (il!=hf.end()) e_long=il->energy()*m_correctionByEta[indexByEta(il->id())];
00255 if (e_long <= m_minTowerEnergy) e_long=0.0;
00256 if (is!=hf.end()) e_short=is->energy()*m_correctionByEta[indexByEta(is->id())];
00257 if (e_short <= m_minTowerEnergy) e_short=0.0;
00258 double eRatio=(e_long-e_short)/std::max(1.0,(e_long+e_short));
00259
00260
00261 if ((abs(eRatio) > m_maximumSL)&&(std::max(e_long,e_short) > m_maximumRenergy)) {
00262 if (dp==0 && de==0) clusterOk=false;
00263 continue;
00264 }
00265
00266 if((il!=hf.end())&&(isPMTHit(*il))){
00267 if (dp==0 && de==0) clusterOk=false;
00268 continue;
00269 }
00270
00271
00272
00273
00274
00275
00276
00277
00278
00279
00280
00281
00282
00283
00284
00285
00286
00287
00288
00289
00290
00291 if (e_long > m_minTowerEnergy && il!=hf.end()) {
00292
00293
00294 usedHits.push_back(idl.rawId());
00295
00296 l_5+=e_long;
00297
00298 if ((de>-2)&&(de<2)&&(dp>-4)&&(dp<4)) {
00299 l_3+=e_long;
00300 }
00301
00302 if ((dp==0)&&(de==0)) {
00303 l_1=e_long;
00304 }
00305
00306
00307 if ((de>-2)&&(de<2)&&(dp>-4)&&(dp<4)&&(e_long>(.5*e_seed))) {
00308 coreCanid.push_back(e_long);
00309 }
00310
00311
00312 GlobalPoint p=geom.getPosition(idl);
00313
00314 double d_p = p.phi()-sp.phi();
00315 while (d_p < -M_PI)
00316 d_p+=2*M_PI;
00317 while (d_p > M_PI)
00318 d_p-=2*M_PI;
00319 double d_e = p.eta()-sp.eta();
00320
00321 wgt=log((e_long));
00322 if (wgt>0){
00323 w+=wgt;
00324 w_e+=(d_e)*wgt;
00325 wp_e+=(d_p)*wgt;
00326 e_e+=d_e;
00327 e_ep+=d_p;
00328
00329 w_x+=(p.x())*wgt;
00330 w_y+=(p.y())*wgt;
00331 w_z+=(p.z())*wgt;
00332 }
00333 } else {
00334 if (dp==0 && de==0) clusterOk=false;
00335 }
00336
00337 if (e_short > m_minTowerEnergy && is!=hf.end()) {
00338
00339 usedHits.push_back(ids.rawId());
00340
00341 s_5+=e_short;
00342
00343 if ((de>-2)&&(de<2)&&(dp>-4)&&(dp<4)) {
00344 s_3+=e_short;
00345 }
00346
00347 if ((dp==0)&&(de==0)) {
00348 s_1=e_short;
00349 }
00350 }
00351 }
00352
00353
00354 if (!clusterOk) return false;
00355
00356
00357 std::sort(coreCanid.begin(), coreCanid.end(), CompareHFCore());
00358 for (ci=coreCanid.begin();ci!=coreCanid.end();ci++){
00359 if(ci==coreCanid.begin()){
00360 l_1e=*ci;
00361 }else if (*ci>.5*l_1e){
00362 l_1e+=*ci;
00363 }
00364 }
00365
00366 double z_=w_z/w;
00367 double x_=w_x/w;
00368 double y_=w_y/w;
00369
00370
00371 double eta=w_e/w+sp.eta();
00372
00373 double phi=(wp_e/w)+sp.phi();
00374
00375 while (phi < -M_PI)
00376 phi+=2*M_PI;
00377 while (phi > M_PI)
00378 phi-=2*M_PI;
00379
00380
00381 static const double HFEtaBounds[14] = {2.853, 2.964, 3.139, 3.314, 3.489, 3.664, 3.839, 4.013, 4.191, 4.363, 4.538, 4.716, 4.889, 5.191};
00382 double RcellEta = fabs(eta);
00383 double Cphi = (phi>0.)?(fmod((phi),0.087*2)/(0.087*2)):((fmod((phi),0.087*2)/(0.087*2))+1.0);
00384 double Rbin = -1.0;
00385 for (int icell = 0; icell < 12; icell++ ){
00386 if ( (RcellEta>HFEtaBounds[icell]) && (RcellEta<HFEtaBounds[icell+1]) )
00387 Rbin = (RcellEta - HFEtaBounds[icell])/(HFEtaBounds[icell+1] - HFEtaBounds[icell]);
00388 }
00389 double Ceta=Rbin;
00390
00391 while (phi< -M_PI)
00392 phi+=2*M_PI;
00393 while (phi > M_PI)
00394 phi-=2*M_PI;
00395
00396
00397 math::XYZPoint xyzclus(x_,y_,z_);
00398
00399
00400 HFEMClusterShape myClusShp(l_1, s_1, l_3, s_3, l_5,s_5, l_1e,Ceta, Cphi,seedid);
00401 clusShp = myClusShp;
00402
00403 SuperCluster MySclus(l_3,xyzclus);
00404 Sclus=MySclus;
00405
00406 return clusterOk;
00407
00408 }
00409 bool HFClusterAlgo::isPMTHit(const HFRecHit& hfr){
00410
00411 bool pmthit=false;
00412
00413 if((hfr.flagField(HcalCaloFlagLabels::HFLongShort))&&(m_usePMTFlag)) pmthit=true;
00414 if (!(m_isMC && !m_forcePulseFlagMC))
00415 if((hfr.flagField(HcalCaloFlagLabels::HFDigiTime))&&(m_usePulseFlag)) pmthit=true;
00416
00417 return pmthit;
00418
00419
00420 }