00001 #include <iostream>
00002
00003 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
00004 #include "RecoEcal/EgammaCoreTools/interface/ClusterShapeAlgo.h"
00005 #include "RecoCaloTools/Navigation/interface/CaloNavigator.h"
00006 #include "Geometry/EcalAlgo/interface/EcalBarrelGeometry.h"
00007 #include "Geometry/CaloTopology/interface/EcalBarrelTopology.h"
00008 #include "Geometry/CaloTopology/interface/EcalEndcapTopology.h"
00009 #include "Geometry/CaloTopology/interface/EcalPreshowerTopology.h"
00010
00011 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00012 #include "DataFormats/Math/interface/Point3D.h"
00013 #include "DataFormats/Math/interface/Vector3D.h"
00014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00015
00016 #include "CLHEP/Geometry/Transform3D.h"
00017 #include "CLHEP/Units/SystemOfUnits.h"
00018
00019 ClusterShapeAlgo::ClusterShapeAlgo(const std::map<std::string,double> & passedParameterMap) : parameterMap_(passedParameterMap) {}
00020
00021 reco::ClusterShape ClusterShapeAlgo::Calculate(const reco::BasicCluster &passedCluster,
00022 const EcalRecHitCollection *hits,
00023 const CaloSubdetectorGeometry * geometry,
00024 const CaloSubdetectorTopology* topology)
00025 {
00026 Calculate_TopEnergy(passedCluster,hits);
00027 Calculate_2ndEnergy(passedCluster,hits);
00028 Create_Map(hits,topology);
00029 Calculate_e2x2();
00030 Calculate_e3x2();
00031 Calculate_e3x3();
00032 Calculate_e4x4();
00033 Calculate_e5x5();
00034 Calculate_e2x5Right();
00035 Calculate_e2x5Left();
00036 Calculate_e2x5Top();
00037 Calculate_e2x5Bottom();
00038 Calculate_Covariances(passedCluster,hits,geometry);
00039 Calculate_BarrelBasketEnergyFraction(passedCluster,hits, Eta, geometry);
00040 Calculate_BarrelBasketEnergyFraction(passedCluster,hits, Phi, geometry);
00041 Calculate_EnergyDepTopology (passedCluster,hits,geometry,true) ;
00042 Calculate_lat(passedCluster);
00043 Calculate_ComplexZernikeMoments(passedCluster);
00044
00045 return reco::ClusterShape(covEtaEta_, covEtaPhi_, covPhiPhi_, eMax_, eMaxId_,
00046 e2nd_, e2ndId_, e2x2_, e3x2_, e3x3_,e4x4_, e5x5_,
00047 e2x5Right_, e2x5Left_, e2x5Top_, e2x5Bottom_,
00048 e3x2Ratio_, lat_, etaLat_, phiLat_, A20_, A42_,
00049 energyBasketFractionEta_, energyBasketFractionPhi_);
00050 }
00051
00052 void ClusterShapeAlgo::Calculate_TopEnergy(const reco::BasicCluster &passedCluster,const EcalRecHitCollection *hits)
00053 {
00054 double eMax=0;
00055 DetId eMaxId(0);
00056
00057 std::vector<DetId> clusterDetIds = passedCluster.getHitsByDetId();
00058 std::vector<DetId>::iterator posCurrent;
00059
00060 EcalRecHit testEcalRecHit;
00061
00062 for(posCurrent = clusterDetIds.begin(); posCurrent != clusterDetIds.end(); posCurrent++)
00063 {
00064 if ((*posCurrent != DetId(0)) && (hits->find(*posCurrent) != hits->end()))
00065 {
00066 EcalRecHitCollection::const_iterator itt = hits->find(*posCurrent);
00067 testEcalRecHit = *itt;
00068
00069 if(testEcalRecHit.energy() > eMax)
00070 {
00071 eMax = testEcalRecHit.energy();
00072 eMaxId = testEcalRecHit.id();
00073 }
00074 }
00075 }
00076
00077 eMax_ = eMax;
00078 eMaxId_ = eMaxId;
00079 }
00080
00081 void ClusterShapeAlgo::Calculate_2ndEnergy(const reco::BasicCluster &passedCluster,const EcalRecHitCollection *hits)
00082 {
00083 double e2nd=0;
00084 DetId e2ndId(0);
00085
00086 std::vector<DetId> clusterDetIds = passedCluster.getHitsByDetId();
00087 std::vector<DetId>::iterator posCurrent;
00088
00089 EcalRecHit testEcalRecHit;
00090
00091 for(posCurrent = clusterDetIds.begin(); posCurrent != clusterDetIds.end(); posCurrent++)
00092 {
00093 if ((*posCurrent != DetId(0)) && (hits->find(*posCurrent) != hits->end()))
00094 {
00095 EcalRecHitCollection::const_iterator itt = hits->find(*posCurrent);
00096 testEcalRecHit = *itt;
00097
00098 if(testEcalRecHit.energy() > e2nd && testEcalRecHit.id() != eMaxId_)
00099 {
00100 e2nd = testEcalRecHit.energy();
00101 e2ndId = testEcalRecHit.id();
00102 }
00103 }
00104 }
00105
00106 e2nd_ = e2nd;
00107 e2ndId_ = e2ndId;
00108 }
00109
00110 void ClusterShapeAlgo::Create_Map(const EcalRecHitCollection *hits,const CaloSubdetectorTopology* topology)
00111 {
00112 EcalRecHit tempEcalRecHit;
00113 CaloNavigator<DetId> posCurrent = CaloNavigator<DetId>(eMaxId_,topology );
00114
00115 for(int x = 0; x < 5; x++)
00116 for(int y = 0; y < 5; y++)
00117 {
00118 posCurrent.home();
00119 posCurrent.offsetBy(-2+x,-2+y);
00120
00121 if((*posCurrent != DetId(0)) && (hits->find(*posCurrent) != hits->end()))
00122 {
00123 EcalRecHitCollection::const_iterator itt = hits->find(*posCurrent);
00124 tempEcalRecHit = *itt;
00125 energyMap_[y][x] = std::make_pair(tempEcalRecHit.id(),tempEcalRecHit.energy());
00126 }
00127 else
00128 energyMap_[y][x] = std::make_pair(DetId(0), 0);
00129 }
00130 }
00131
00132 void ClusterShapeAlgo::Calculate_e2x2()
00133 {
00134 double e2x2Max = 0;
00135 double e2x2Test = 0;
00136
00137 int deltaX=0, deltaY=0;
00138
00139 for(int corner = 0; corner < 4; corner++)
00140 {
00141 switch(corner)
00142 {
00143 case 0: deltaX = -1; deltaY = -1; break;
00144 case 1: deltaX = -1; deltaY = 1; break;
00145 case 2: deltaX = 1; deltaY = -1; break;
00146 case 3: deltaX = 1; deltaY = 1; break;
00147 }
00148
00149 e2x2Test = energyMap_[2][2].second;
00150 e2x2Test += energyMap_[2+deltaY][2].second;
00151 e2x2Test += energyMap_[2][2+deltaX].second;
00152 e2x2Test += energyMap_[2+deltaY][2+deltaX].second;
00153
00154 if(e2x2Test > e2x2Max)
00155 {
00156 e2x2Max = e2x2Test;
00157 e2x2_Diagonal_X_ = 2+deltaX;
00158 e2x2_Diagonal_Y_ = 2+deltaY;
00159 }
00160 }
00161
00162 e2x2_ = e2x2Max;
00163
00164 }
00165
00166 void ClusterShapeAlgo::Calculate_e3x2()
00167 {
00168 double e3x2 = 0.0;
00169 double e3x2Ratio = 0.0, e3x2RatioNumerator = 0.0, e3x2RatioDenominator = 0.0;
00170
00171 int e2ndX = 2, e2ndY = 2;
00172 int deltaY = 0, deltaX = 0;
00173
00174 double nextEnergy = -999;
00175 int nextEneryDirection = -1;
00176
00177 for(int cardinalDirection = 0; cardinalDirection < 4; cardinalDirection++)
00178 {
00179 switch(cardinalDirection)
00180 {
00181 case 0: deltaX = -1; deltaY = 0; break;
00182 case 1: deltaX = 1; deltaY = 0; break;
00183 case 2: deltaX = 0; deltaY = -1; break;
00184 case 3: deltaX = 0; deltaY = 1; break;
00185 }
00186
00187 if(energyMap_[2+deltaY][2+deltaX].second >= nextEnergy)
00188 {
00189 nextEnergy = energyMap_[2+deltaY][2+deltaX].second;
00190 nextEneryDirection = cardinalDirection;
00191
00192 e2ndX = 2+deltaX;
00193 e2ndY = 2+deltaY;
00194 }
00195 }
00196
00197 switch(nextEneryDirection)
00198 {
00199 case 0: ;
00200 case 1: deltaX = 0; deltaY = 1; break;
00201 case 2: ;
00202 case 3: deltaX = 1; deltaY = 0; break;
00203 }
00204
00205 for(int sign = -1; sign <= 1; sign++)
00206 e3x2 += (energyMap_[2+deltaY*sign][2+deltaX*sign].second + energyMap_[e2ndY+deltaY*sign][e2ndX+deltaX*sign].second);
00207
00208 e3x2RatioNumerator = energyMap_[e2ndY+deltaY][e2ndX+deltaX].second + energyMap_[e2ndY-deltaY][e2ndX-deltaX].second;
00209 e3x2RatioDenominator = 0.5 + energyMap_[2+deltaY][2+deltaX].second + energyMap_[2-deltaY][2-deltaX].second;
00210 e3x2Ratio = e3x2RatioNumerator / e3x2RatioDenominator;
00211
00212 e3x2_ = e3x2;
00213 e3x2Ratio_ = e3x2Ratio;
00214 }
00215
00216 void ClusterShapeAlgo::Calculate_e3x3()
00217 {
00218 double e3x3=0;
00219
00220 for(int i = 1; i <= 3; i++)
00221 for(int j = 1; j <= 3; j++)
00222 e3x3 += energyMap_[j][i].second;
00223
00224 e3x3_ = e3x3;
00225
00226 }
00227
00228 void ClusterShapeAlgo::Calculate_e4x4()
00229 {
00230 double e4x4=0;
00231
00232 int startX=-1, startY=-1;
00233
00234 switch(e2x2_Diagonal_X_)
00235 {
00236 case 1: startX = 0; break;
00237 case 3: startX = 1; break;
00238 }
00239
00240 switch(e2x2_Diagonal_Y_)
00241 {
00242 case 1: startY = 0; break;
00243 case 3: startY = 1; break;
00244 }
00245
00246 for(int i = startX; i <= startX+3; i++)
00247 for(int j = startY; j <= startY+3; j++)
00248 e4x4 += energyMap_[j][i].second;
00249
00250 e4x4_ = e4x4;
00251 }
00252
00253 void ClusterShapeAlgo::Calculate_e5x5()
00254 {
00255 double e5x5=0;
00256
00257 for(int i = 0; i <= 4; i++)
00258 for(int j = 0; j <= 4; j++)
00259 e5x5 += energyMap_[j][i].second;
00260
00261 e5x5_ = e5x5;
00262
00263 }
00264
00265 void ClusterShapeAlgo::Calculate_e2x5Right()
00266 {
00267 double e2x5R=0.0;
00268 for(int i = 0; i <= 4; i++){
00269 for(int j = 0; j <= 4; j++){
00270 if(j>2){e2x5R +=energyMap_[i][j].second;}
00271 }
00272 }
00273 e2x5Right_=e2x5R;
00274 }
00275
00276 void ClusterShapeAlgo::Calculate_e2x5Left()
00277 {
00278 double e2x5L=0.0;
00279 for(int i = 0; i <= 4; i++){
00280 for(int j = 0; j <= 4; j++){
00281 if(j<2){e2x5L +=energyMap_[i][j].second;}
00282 }
00283 }
00284 e2x5Left_=e2x5L;
00285 }
00286
00287 void ClusterShapeAlgo::Calculate_e2x5Bottom()
00288 {
00289 double e2x5B=0.0;
00290 for(int i = 0; i <= 4; i++){
00291 for(int j = 0; j <= 4; j++){
00292
00293 if(i>2){e2x5B +=energyMap_[i][j].second;}
00294 }
00295 }
00296 e2x5Bottom_=e2x5B;
00297 }
00298
00299 void ClusterShapeAlgo::Calculate_e2x5Top()
00300 {
00301 double e2x5T=0.0;
00302 for(int i = 0; i <= 4; i++){
00303 for(int j = 0; j <= 4; j++){
00304
00305 if(i<2){e2x5T +=energyMap_[i][j].second;}
00306 }
00307 }
00308 e2x5Top_=e2x5T;
00309 }
00310
00311 void ClusterShapeAlgo::Calculate_Covariances(const reco::BasicCluster &passedCluster, const EcalRecHitCollection* hits,
00312 const CaloSubdetectorGeometry* geometry)
00313 {
00314 if (e5x5_ > 0.)
00315 {
00316 double w0_ = parameterMap_.find("W0")->second;
00317
00318
00319
00320 math::XYZVector meanPosition(0.0, 0.0, 0.0);
00321 for (int i = 0; i < 5; ++i)
00322 {
00323 for (int j = 0; j < 5; ++j)
00324 {
00325 DetId id = energyMap_[i][j].first;
00326 if (id != DetId(0))
00327 {
00328 GlobalPoint positionGP = geometry->getGeometry(id)->getPosition();
00329 math::XYZVector position(positionGP.x(),positionGP.y(),positionGP.z());
00330 meanPosition = meanPosition + energyMap_[i][j].second * position;
00331 }
00332 }
00333 }
00334
00335 meanPosition /= e5x5_;
00336
00337
00338 double numeratorEtaEta = 0;
00339 double numeratorEtaPhi = 0;
00340 double numeratorPhiPhi = 0;
00341 double denominator = 0;
00342
00343 for (int i = 0; i < 5; ++i)
00344 {
00345 for (int j = 0; j < 5; ++j)
00346 {
00347 DetId id = energyMap_[i][j].first;
00348 if (id != DetId(0))
00349 {
00350 GlobalPoint position = geometry->getGeometry(id)->getPosition();
00351
00352 double dPhi = position.phi() - meanPosition.phi();
00353 if (dPhi > + Geom::pi()) { dPhi = Geom::twoPi() - dPhi; }
00354 if (dPhi < - Geom::pi()) { dPhi = Geom::twoPi() + dPhi; }
00355
00356 double dEta = position.eta() - meanPosition.eta();
00357 double w = 0.;
00358 if ( energyMap_[i][j].second > 0.)
00359 w = std::max(0.0, w0_ + log( energyMap_[i][j].second / e5x5_));
00360
00361 denominator += w;
00362 numeratorEtaEta += w * dEta * dEta;
00363 numeratorEtaPhi += w * dEta * dPhi;
00364 numeratorPhiPhi += w * dPhi * dPhi;
00365 }
00366 }
00367 }
00368
00369 covEtaEta_ = numeratorEtaEta / denominator;
00370 covEtaPhi_ = numeratorEtaPhi / denominator;
00371 covPhiPhi_ = numeratorPhiPhi / denominator;
00372 }
00373 else
00374 {
00375
00376
00377 covEtaEta_ = 0;
00378 covEtaPhi_ = 0;
00379 covPhiPhi_ = 0;
00380 }
00381 }
00382
00383 void ClusterShapeAlgo::Calculate_BarrelBasketEnergyFraction(const reco::BasicCluster &passedCluster,
00384 const EcalRecHitCollection *hits,
00385 const int EtaPhi,
00386 const CaloSubdetectorGeometry* geometry)
00387 {
00388 if( (hits!=0) && ( ((*hits)[0]).id().subdetId() != EcalBarrel ) ) {
00389
00390 return;
00391 }
00392
00393 std::map<int,double> indexedBasketEnergy;
00394 std::vector<DetId> clusterDetIds = passedCluster.getHitsByDetId();
00395 const EcalBarrelGeometry* subDetGeometry = (const EcalBarrelGeometry*) geometry;
00396
00397 for(std::vector<DetId>::iterator posCurrent = clusterDetIds.begin(); posCurrent != clusterDetIds.end(); posCurrent++)
00398 {
00399 int basketIndex = 999;
00400
00401 if(EtaPhi == Eta) {
00402 int unsignedIEta = abs(EBDetId(*posCurrent).ieta());
00403 std::vector<int> etaBasketSize = subDetGeometry->getEtaBaskets();
00404
00405 for(unsigned int i = 0; i < etaBasketSize.size(); i++) {
00406 unsignedIEta -= etaBasketSize[i];
00407 if(unsignedIEta - 1 < 0)
00408 {
00409 basketIndex = i;
00410 break;
00411 }
00412 }
00413 basketIndex = (basketIndex+1)*(EBDetId(*posCurrent).ieta() > 0 ? 1 : -1);
00414
00415 } else if(EtaPhi == Phi) {
00416 int halfNumBasketInPhi = (EBDetId::MAX_IPHI - EBDetId::MIN_IPHI + 1)/subDetGeometry->getBasketSizeInPhi()/2;
00417
00418 basketIndex = (EBDetId(*posCurrent).iphi() - 1)/subDetGeometry->getBasketSizeInPhi()
00419 - (EBDetId(clusterDetIds[0]).iphi() - 1)/subDetGeometry->getBasketSizeInPhi();
00420
00421 if(basketIndex >= halfNumBasketInPhi) basketIndex -= 2*halfNumBasketInPhi;
00422 else if(basketIndex < -1*halfNumBasketInPhi) basketIndex += 2*halfNumBasketInPhi;
00423
00424 } else throw(std::runtime_error("\n\nOh No! Calculate_BarrelBasketEnergyFraction called on invalid index.\n\n"));
00425
00426 indexedBasketEnergy[basketIndex] += (hits->find(*posCurrent))->energy();
00427 }
00428
00429 std::vector<double> energyFraction;
00430 for(std::map<int,double>::iterator posCurrent = indexedBasketEnergy.begin(); posCurrent != indexedBasketEnergy.end(); posCurrent++)
00431 {
00432 energyFraction.push_back(indexedBasketEnergy[posCurrent->first]/passedCluster.energy());
00433 }
00434
00435 switch(EtaPhi)
00436 {
00437 case Eta: energyBasketFractionEta_ = energyFraction; break;
00438 case Phi: energyBasketFractionPhi_ = energyFraction; break;
00439 }
00440
00441 }
00442
00443 void ClusterShapeAlgo::Calculate_lat(const reco::BasicCluster &passedCluster) {
00444
00445 double r,redmoment=0;
00446 double phiRedmoment = 0 ;
00447 double etaRedmoment = 0 ;
00448 int n,n1,n2,tmp;
00449 int clusterSize=energyDistribution_.size();
00450 if (clusterSize<3) {
00451 etaLat_ = 0.0 ;
00452 lat_ = 0.0;
00453 return;
00454 }
00455
00456 n1=0; n2=1;
00457 if (energyDistribution_[1].deposited_energy >
00458 energyDistribution_[0].deposited_energy)
00459 {
00460 tmp=n2; n2=n1; n1=tmp;
00461 }
00462 for (int i=2; i<clusterSize; i++) {
00463 n=i;
00464 if (energyDistribution_[i].deposited_energy >
00465 energyDistribution_[n1].deposited_energy)
00466 {
00467 tmp = n2;
00468 n2 = n1; n1 = i; n=tmp;
00469 } else {
00470 if (energyDistribution_[i].deposited_energy >
00471 energyDistribution_[n2].deposited_energy)
00472 {
00473 tmp=n2; n2=i; n=tmp;
00474 }
00475 }
00476
00477 r = energyDistribution_[n].r;
00478 redmoment += r*r* energyDistribution_[n].deposited_energy;
00479 double rphi = r * cos (energyDistribution_[n].phi) ;
00480 phiRedmoment += rphi * rphi * energyDistribution_[n].deposited_energy;
00481 double reta = r * sin (energyDistribution_[n].phi) ;
00482 etaRedmoment += reta * reta * energyDistribution_[n].deposited_energy;
00483 }
00484 double e1 = energyDistribution_[n1].deposited_energy;
00485 double e2 = energyDistribution_[n2].deposited_energy;
00486
00487 lat_ = redmoment/(redmoment+2.19*2.19*(e1+e2));
00488 phiLat_ = phiRedmoment/(phiRedmoment+2.19*2.19*(e1+e2));
00489 etaLat_ = etaRedmoment/(etaRedmoment+2.19*2.19*(e1+e2));
00490 }
00491
00492 void ClusterShapeAlgo::Calculate_ComplexZernikeMoments(const reco::BasicCluster &passedCluster) {
00493
00494
00495 A20_ = absZernikeMoment(passedCluster,2,0);
00496 A42_ = absZernikeMoment(passedCluster,4,2);
00497 }
00498
00499 double ClusterShapeAlgo::absZernikeMoment(const reco::BasicCluster &passedCluster,
00500 int n, int m, double R0) {
00501
00502 if ((m>n) || ((n-m)%2 != 0) || (n<0) || (m<0)) return -1;
00503
00504
00505
00506 if ((n>20) || (R0<=2.19)) return -1;
00507 if (n<=5) return fast_AbsZernikeMoment(passedCluster,n,m,R0);
00508 else return calc_AbsZernikeMoment(passedCluster,n,m,R0);
00509 }
00510
00511 double ClusterShapeAlgo::f00(double r) { return 1; }
00512
00513 double ClusterShapeAlgo::f11(double r) { return r; }
00514
00515 double ClusterShapeAlgo::f20(double r) { return 2.0*r*r-1.0; }
00516
00517 double ClusterShapeAlgo::f22(double r) { return r*r; }
00518
00519 double ClusterShapeAlgo::f31(double r) { return 3.0*r*r*r - 2.0*r; }
00520
00521 double ClusterShapeAlgo::f33(double r) { return r*r*r; }
00522
00523 double ClusterShapeAlgo::f40(double r) { return 6.0*r*r*r*r-6.0*r*r+1.0; }
00524
00525 double ClusterShapeAlgo::f42(double r) { return 4.0*r*r*r*r-3.0*r*r; }
00526
00527 double ClusterShapeAlgo::f44(double r) { return r*r*r*r; }
00528
00529 double ClusterShapeAlgo::f51(double r) { return 10.0*pow(r,5)-12.0*pow(r,3)+3.0*r; }
00530
00531 double ClusterShapeAlgo::f53(double r) { return 5.0*pow(r,5) - 4.0*pow(r,3); }
00532
00533 double ClusterShapeAlgo::f55(double r) { return pow(r,5); }
00534
00535 double ClusterShapeAlgo::fast_AbsZernikeMoment(const reco::BasicCluster &passedCluster,
00536 int n, int m, double R0) {
00537 double r,ph,e,Re=0,Im=0,result;
00538 double TotalEnergy = passedCluster.energy();
00539 int index = (n/2)*(n/2)+(n/2)+m;
00540 int clusterSize=energyDistribution_.size();
00541 if(clusterSize<3) return 0.0;
00542
00543 for (int i=0; i<clusterSize; i++)
00544 {
00545 r = energyDistribution_[i].r / R0;
00546 if (r<1) {
00547 fcn_.clear();
00548 Calculate_Polynomials(r);
00549 ph = (energyDistribution_[i]).phi;
00550 e = energyDistribution_[i].deposited_energy;
00551 Re = Re + e/TotalEnergy * fcn_[index] * cos( (double) m * ph);
00552 Im = Im - e/TotalEnergy * fcn_[index] * sin( (double) m * ph);
00553 }
00554 }
00555 result = sqrt(Re*Re+Im*Im);
00556
00557 return result;
00558 }
00559
00560 double ClusterShapeAlgo::calc_AbsZernikeMoment(const reco::BasicCluster &passedCluster,
00561 int n, int m, double R0) {
00562 double r,ph,e,Re=0,Im=0,f_nm,result;
00563 double TotalEnergy = passedCluster.energy();
00564 std::vector<DetId> clusterDetIds = passedCluster.getHitsByDetId();
00565 int clusterSize=energyDistribution_.size();
00566 if(clusterSize<3) return 0.0;
00567
00568 for (int i=0; i<clusterSize; i++)
00569 {
00570 r = energyDistribution_[i].r / R0;
00571 if (r<1) {
00572 ph = (energyDistribution_[i]).phi;
00573 e = energyDistribution_[i].deposited_energy;
00574 f_nm=0;
00575 for (int s=0; s<=(n-m)/2; s++) {
00576 if (s%2==0)
00577 {
00578 f_nm = f_nm + factorial(n-s)*pow(r,(double) (n-2*s))/(factorial(s)*factorial((n+m)/2-s)*factorial((n-m)/2-s));
00579 }else {
00580 f_nm = f_nm - factorial(n-s)*pow(r,(double) (n-2*s))/(factorial(s)*factorial((n+m)/2-s)*factorial((n-m)/2-s));
00581 }
00582 }
00583 Re = Re + e/TotalEnergy * f_nm * cos( (double) m*ph);
00584 Im = Im - e/TotalEnergy * f_nm * sin( (double) m*ph);
00585 }
00586 }
00587 result = sqrt(Re*Re+Im*Im);
00588
00589 return result;
00590 }
00591
00592 void ClusterShapeAlgo::Calculate_EnergyDepTopology (const reco::BasicCluster &passedCluster,
00593 const EcalRecHitCollection *hits,
00594 const CaloSubdetectorGeometry* geometry,
00595 bool logW) {
00596
00597 energyDistribution_.clear();
00598
00599
00600
00601 CLHEP::Hep3Vector clVect(passedCluster.position().x(),
00602 passedCluster.position().y(),
00603 passedCluster.position().z());
00604 CLHEP::Hep3Vector clDir(clVect);
00605 clDir*=1.0/clDir.mag();
00606
00607 CLHEP::Hep3Vector theta_axis(clDir.y(),-clDir.x(),0.0);
00608 theta_axis *= 1.0/theta_axis.mag();
00609 Hep3Vector phi_axis = theta_axis.cross(clDir);
00610
00611 std::vector<DetId> clusterDetIds = passedCluster.getHitsByDetId();
00612
00613 EcalClusterEnergyDeposition clEdep;
00614 EcalRecHit testEcalRecHit;
00615 std::vector<DetId>::iterator posCurrent;
00616
00617 for(posCurrent=clusterDetIds.begin(); posCurrent!=clusterDetIds.end(); ++posCurrent) {
00618 EcalRecHitCollection::const_iterator itt = hits->find(*posCurrent);
00619 testEcalRecHit=*itt;
00620
00621 if((*posCurrent != DetId(0)) && (hits->find(*posCurrent) != hits->end())) {
00622 clEdep.deposited_energy = testEcalRecHit.energy();
00623
00624
00625 if(logW) {
00626 double w0_ = parameterMap_.find("W0")->second;
00627
00628 double weight = std::max(0.0, w0_ + log(fabs(clEdep.deposited_energy)/passedCluster.energy()) );
00629 if(weight==0) {
00630 LogDebug("ClusterShapeAlgo") << "Crystal has insufficient energy: E = "
00631 << clEdep.deposited_energy << " GeV; skipping... ";
00632 continue;
00633 }
00634 else LogDebug("ClusterShapeAlgo") << "===> got crystal. Energy = " << clEdep.deposited_energy << " GeV. ";
00635 }
00636 DetId id_ = *posCurrent;
00637 const CaloCellGeometry *this_cell = geometry->getGeometry(id_);
00638 GlobalPoint cellPos = this_cell->getPosition();
00639 CLHEP::Hep3Vector gblPos (cellPos.x(),cellPos.y(),cellPos.z());
00640
00641 CLHEP::Hep3Vector diff = gblPos - clVect;
00642
00643
00644
00645 CLHEP::Hep3Vector DigiVect = diff - diff.dot(clDir)*clDir;
00646 clEdep.r = DigiVect.mag();
00647 LogDebug("ClusterShapeAlgo") << "E = " << clEdep.deposited_energy
00648 << "\tdiff = " << diff.mag()
00649 << "\tr = " << clEdep.r;
00650 clEdep.phi = DigiVect.angle(theta_axis);
00651 if(DigiVect.dot(phi_axis)<0) clEdep.phi = 2*M_PI - clEdep.phi;
00652 energyDistribution_.push_back(clEdep);
00653 }
00654 }
00655 }
00656
00657 void ClusterShapeAlgo::Calculate_Polynomials(double rho) {
00658 fcn_.push_back(f00(rho));
00659 fcn_.push_back(f11(rho));
00660 fcn_.push_back(f20(rho));
00661 fcn_.push_back(f31(rho));
00662 fcn_.push_back(f22(rho));
00663 fcn_.push_back(f33(rho));
00664 fcn_.push_back(f40(rho));
00665 fcn_.push_back(f51(rho));
00666 fcn_.push_back(f42(rho));
00667 fcn_.push_back(f53(rho));
00668 fcn_.push_back(f44(rho));
00669 fcn_.push_back(f55(rho));
00670 }
00671
00672 double ClusterShapeAlgo::factorial(int n) const {
00673 double res=1.0;
00674 for(int i=2; i<=n; i++) res*=(double) i;
00675 return res;
00676 }