CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoEcal/EgammaCoreTools/src/ClusterShapeAlgo.cc

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