CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/RecoEcal/EgammaClusterAlgos/src/HybridClusterAlgo.cc

Go to the documentation of this file.
00001 #include "RecoEcal/EgammaClusterAlgos/interface/HybridClusterAlgo.h"
00002 #include "RecoCaloTools/Navigation/interface/EcalBarrelNavigator.h"
00003 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00004 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00005 #include <iostream>
00006 #include <map>
00007 #include <vector>
00008 #include <set>
00009 #include "RecoEcal/EgammaCoreTools/interface/ClusterEtLess.h"
00010 
00011 //The real constructor
00012 HybridClusterAlgo::HybridClusterAlgo(double eb_str, 
00013                                      int step, 
00014                                      double ethres,
00015                                      double eseed,
00016                                      double ewing,
00017                                      std::vector<int> v_chstatus,
00018                                      const PositionCalc& posCalculator,
00019                                      DebugLevel debugLevel,
00020                                      bool dynamicEThres,
00021                                      double eThresA,
00022                                      double eThresB,
00023                                      std::vector<int> severityToExclude,
00024                                      double severityRecHitThreshold,
00025                                      int severitySpikeId,
00026                                      double severitySpikeThreshold,
00027                                      bool excludeFromCluster
00028                                      ) :
00029   
00030   eb_st(eb_str), phiSteps_(step), 
00031   eThres_(ethres), eThresA_(eThresA), eThresB_(eThresB),
00032   Eseed(eseed),  Ewing(ewing), 
00033   dynamicEThres_(dynamicEThres), debugLevel_(debugLevel),
00034   v_chstatus_(v_chstatus), v_severitylevel_(severityToExclude),severityRecHitThreshold_(severityRecHitThreshold), severitySpikeThreshold_(severitySpikeThreshold), excludeFromCluster_(excludeFromCluster)
00035   
00036   
00037 {
00038   spId_ = EcalSeverityLevelAlgo::SpikeId(severitySpikeId);
00039   dynamicPhiRoad_ = false;
00040   if ( debugLevel_ == pDEBUG ) {
00041     //std::cout << "dynamicEThres: " << dynamicEThres_ 
00042     //          << " : A,B " << eThresA_ << ", " << eThresB_ << std::endl;
00043   }
00044   
00045   //if (dynamicPhiRoad_) phiRoadAlgo_ = new BremRecoveryPhiRoadAlgo(bremRecoveryPset);
00046   posCalculator_ = posCalculator;
00047   topo_ = new EcalBarrelHardcodedTopology();
00048   
00049   std::sort( v_chstatus_.begin(), v_chstatus_.end() );
00050   
00051   
00052 }
00053 
00054 // Return a vector of clusters from a collection of EcalRecHits:
00055 void HybridClusterAlgo::makeClusters(const EcalRecHitCollection*recColl, 
00056                                      const CaloSubdetectorGeometry*geometry,
00057                                      reco::BasicClusterCollection &basicClusters,
00058                                      bool regional,
00059                                      const std::vector<EcalEtaPhiRegion>& regions,
00060                                      const EcalChannelStatus*chStatus
00061                                      )
00062 {
00063   //clear vector of seeds
00064   seeds.clear();
00065   //clear flagged crystals
00066   excludedCrys_.clear();
00067   //clear map of supercluster/basiccluster association
00068   clustered_.clear();
00069   //clear set of used detids
00070   useddetids.clear();
00071   //clear vector of seed clusters
00072   seedClus_.clear();
00073   //Pass in a pointer to the collection.
00074   recHits_ = recColl;
00075   
00076   //
00077   //  SCShape_ = new SuperClusterShapeAlgo(recHits_, geometry);
00078   
00079   if ( debugLevel_ == pDEBUG ) {
00080     std::cout << "Cleared vectors, starting clusterization..." << std::endl;
00081     std::cout << "Purple monkey aardvark." << std::endl;
00082   }
00083   
00084   int nregions=0;
00085   if(regional) nregions=regions.size();
00086   
00087   if(!regional || nregions) {
00088     
00089     EcalRecHitCollection::const_iterator it;
00090     
00091     for (it = recHits_->begin(); it != recHits_->end(); it++){
00092       
00093       //Make the vector of seeds that we're going to use.
00094       //One of the few places position is used, needed for ET calculation.    
00095       const CaloCellGeometry *this_cell = (*geometry).getGeometry(it->id());
00096       GlobalPoint position = this_cell->getPosition();
00097       
00098       
00099       // Require that RecHit is within clustering region in case
00100       // of regional reconstruction
00101       bool withinRegion = false;
00102       if (regional) {
00103         std::vector<EcalEtaPhiRegion>::const_iterator region;
00104         for (region=regions.begin(); region!=regions.end(); region++) {
00105           if (region->inRegion(position)) {
00106             withinRegion =  true;
00107             break;
00108           }
00109         }
00110       }
00111       
00112       if (!regional || withinRegion) {
00113 
00114         //Must pass seed threshold
00115         float ET = it->energy() * sin(position.theta());
00116 
00117         if (ET > eb_st) {
00118           if (debugLevel_==pDEBUG) std::cout << "Seed crystal: " << std::endl;
00119           // avoid seeding for anomalous channels (recoFlag based)
00120           uint32_t rhFlag = (*it).recoFlag();
00121           if (debugLevel_==pDEBUG) std::cout << "rhFlag: " << rhFlag << std::endl;
00122           std::vector<int>::const_iterator vit = std::find( v_chstatus_.begin(), v_chstatus_.end(), rhFlag );
00123           if ( vit != v_chstatus_.end() ){
00124             if (excludeFromCluster_)
00125               excludedCrys_.insert(it->id());
00126             continue; // the recHit has to be excluded from seeding
00127           }
00128           int severityFlag =  EcalSeverityLevelAlgo::severityLevel( it->id(), 
00129                                                                   (*recHits_), 
00130                                                                   (*chStatus),
00131                                                                   severityRecHitThreshold_,
00132                                                                   spId_,
00133                                                                   severitySpikeThreshold_);
00134           std::vector<int>::const_iterator sit = std::find( v_severitylevel_.begin(), v_severitylevel_.end(), severityFlag);
00135           if (debugLevel_ == pDEBUG){
00136             std::cout << "found flag: " << severityFlag << std::endl;
00137           }
00138           
00139           
00140           if (sit!=v_severitylevel_.end()){ 
00141             if (excludeFromCluster_)
00142               excludedCrys_.insert(it->id());
00143             continue;
00144           }
00145           seeds.push_back(*it);
00146           if ( debugLevel_ == pDEBUG ){
00147             std::cout << "Seed ET: " << ET << std::endl;
00148             std::cout << "Seed E: " << it->energy() << std::endl;
00149           }
00150         }
00151       }
00152     }
00153     
00154   }
00155   
00156   
00157   //Yay sorting.
00158   if ( debugLevel_ == pDEBUG )
00159     std::cout << "Built vector of seeds, about to sort them...";
00160   
00161   //Needs three argument sort with seed comparison operator
00162   sort(seeds.begin(), seeds.end(), less_mag());
00163   
00164   if ( debugLevel_ == pDEBUG )
00165     std::cout << "done" << std::endl;
00166   
00167   //Now to do the work.
00168   if ( debugLevel_ ==pDEBUG ) 
00169     std::cout << "About to call mainSearch...";
00170   mainSearch(recColl,geometry);
00171   if ( debugLevel_ == pDEBUG ) 
00172     std::cout << "done" << std::endl;
00173   
00174   //Hand the basicclusters back to the producer.  It has to 
00175   //put them in the event.  Then we can make superclusters.
00176   std::map<int, reco::BasicClusterCollection>::iterator bic; 
00177   for (bic= clustered_.begin();bic!=clustered_.end();bic++){
00178     reco::BasicClusterCollection bl = bic->second;
00179     for (int j=0;j<int(bl.size());++j){
00180       basicClusters.push_back(bl[j]);
00181     }
00182   }
00183   
00184   //Yay more sorting.
00185   sort(basicClusters.rbegin(), basicClusters.rend(), ClusterEtLess() );
00186   //Done!
00187   if ( debugLevel_ == pDEBUG )
00188     std::cout << "returning to producer. " << std::endl;
00189 }
00190 
00191 
00192 void HybridClusterAlgo::mainSearch(const EcalRecHitCollection* hits, const CaloSubdetectorGeometry*geometry)
00193 {
00194 
00195         if ( debugLevel_ ==pDEBUG ) {
00196                 std::cout << "HybridClusterAlgo Algorithm - looking for clusters" << std::endl;
00197                 std::cout << "Found the following clusters:" << std::endl;
00198         }
00199 
00200         // Loop over seeds:
00201         std::vector<EcalRecHit>::iterator it;
00202         int clustercounter=0;
00203 
00204         for (it = seeds.begin(); it != seeds.end(); it++){
00205                 std::vector <reco::BasicCluster> thisseedClusters;
00206                 DetId itID = it->id();
00207 
00208                 // make sure the current seed has not been used/will not be used in the future:
00209                 std::set<DetId>::iterator seed_in_rechits_it = useddetids.find(itID);
00210 
00211                 if (seed_in_rechits_it != useddetids.end()) continue;
00212                 //If this seed is already used, then don't use it again.
00213 
00214                 // output some info on the hit:
00215                 if ( debugLevel_ == pDEBUG ){
00216                         std::cout << "*****************************************************" << std::endl;
00217                         std::cout << "Seed of energy E = " << it->energy() << " @ " << EBDetId(itID) 
00218                                 << std::endl;
00219                         std::cout << "*****************************************************" << std::endl;
00220                 }
00221 
00222                 //Make a navigator, and set it to the seed cell.
00223                 EcalBarrelNavigator navigator(itID, topo_);
00224 
00225                 //Now use the navigator to start building dominoes.
00226 
00227                 //Walking positive in phi:
00228                 std::vector <double> dominoEnergyPhiPlus;  //Here I will store the results of the domino sums
00229                 std::vector <std::vector <EcalRecHit> > dominoCellsPhiPlus; //These are the actual EcalRecHit for dominos.
00230 
00231                 //Walking negative in phi
00232                 std::vector <double> dominoEnergyPhiMinus;  //Here I will store the results of the domino sums
00233                 std::vector <std::vector <EcalRecHit> > dominoCellsPhiMinus; //These are the actual EcalRecHit for dominos.
00234 
00235                 //The two sets together.
00236                 std::vector <double> dominoEnergy;  //Here I will store the results of the domino sums
00237                 std::vector <std::vector <EcalRecHit> > dominoCells; //These are the actual EcalRecHit for dominos.
00238 
00239                 //First, the domino about the seed:
00240                 std::vector <EcalRecHit> initialdomino;
00241                 double e_init = makeDomino(navigator, initialdomino);
00242         if (e_init < Eseed) continue;
00243 
00244                 if ( debugLevel_ == pDEBUG )
00245                 {
00246                         std::cout << "Made initial domino" << std::endl;
00247                 }
00248 
00249                 //
00250                 // compute the phi road length
00251                 double phiSteps;
00252                 if (dynamicPhiRoad_ && e_init > 0)
00253                 {
00254                         double et5x5 = et25(navigator, hits, geometry);
00255                         phiSteps = phiRoadAlgo_->barrelPhiRoad(et5x5);
00256                         navigator.home();
00257                 } else phiSteps = phiSteps_;
00258 
00259                 //Positive phi steps.
00260                 for (int i=0;i<phiSteps;++i){
00261                         //remember, this always increments the current position of the navigator.
00262                         DetId centerD = navigator.north();
00263                         if (centerD.null())
00264                                 continue;
00265                         if ( debugLevel_ == pDEBUG )
00266                         {
00267                                 std::cout << "Step ++" << i << " @ " << EBDetId(centerD) << std::endl;
00268                         }
00269                         EcalBarrelNavigator dominoNav(centerD, topo_);
00270 
00271                         //Go get the new domino.
00272                         std::vector <EcalRecHit> dcells;
00273                         double etemp = makeDomino(dominoNav, dcells);
00274 
00275                         //save this information
00276                         dominoEnergyPhiPlus.push_back(etemp);
00277                         dominoCellsPhiPlus.push_back(dcells);
00278                 }
00279 
00280                 if ( debugLevel_ == pDEBUG )
00281                         std::cout << "Got positive dominos" << std::endl;
00282                 //return to initial position
00283                 navigator.home();
00284 
00285                 //Negative phi steps.
00286                 for (int i=0;i<phiSteps;++i){
00287                         //remember, this always decrements the current position of the navigator.
00288                         DetId centerD = navigator.south();
00289                         if (centerD.null())
00290                                 continue;
00291 
00292                         if ( debugLevel_ == pDEBUG )
00293                         {
00294                                 std::cout << "Step --" << i << " @ " << EBDetId(centerD) << std::endl;
00295                         }
00296                         EcalBarrelNavigator dominoNav(centerD, topo_);
00297 
00298                         //Go get the new domino.
00299                         std::vector <EcalRecHit> dcells;
00300                         double etemp = makeDomino(dominoNav, dcells);
00301 
00302                         //save this information
00303                         dominoEnergyPhiMinus.push_back(etemp);
00304                         dominoCellsPhiMinus.push_back(dcells);
00305                 }
00306 
00307                 if ( debugLevel_ == pDEBUG )
00308                         std::cout << "Got negative dominos: " << std::endl;
00309 
00310                 //Assemble this information:
00311                 for (int i=int(dominoEnergyPhiMinus.size())-1;i >= 0;--i){
00312                         dominoEnergy.push_back(dominoEnergyPhiMinus[i]);
00313                         dominoCells.push_back(dominoCellsPhiMinus[i]);
00314                 }
00315                 dominoEnergy.push_back(e_init);
00316                 dominoCells.push_back(initialdomino);
00317                 for (int i=0;i<int(dominoEnergyPhiPlus.size());++i){
00318                         dominoEnergy.push_back(dominoEnergyPhiPlus[i]);
00319                         dominoCells.push_back(dominoCellsPhiPlus[i]);
00320                 }
00321 
00322                 //Ok, now I have all I need in order to go ahead and make clusters.
00323                 if ( debugLevel_ == pDEBUG ){
00324                         std::cout << "Dumping domino energies: " << std::endl;
00325                         for (int i=0;i<int(dominoEnergy.size());++i){
00326                                 std::cout << "Domino: " << i << " E: " << dominoEnergy[i] << std::endl;
00327                         }
00328                 }
00329 
00330 
00331                 //Identify the peaks in this set of dominos:
00332                 //Peak==a domino whose energy is greater than the two adjacent dominos.
00333                 //thus a peak in the local sense.
00334                 std::vector <int> PeakIndex;
00335                 for (int i=1;i<int(dominoEnergy.size())-1;++i){
00336                         if (dominoEnergy[i] > dominoEnergy[i-1]
00337                                         && dominoEnergy[i] >= dominoEnergy[i+1]
00338                                         && dominoEnergy[i] > Eseed){
00339                                 PeakIndex.push_back(i);
00340                         }
00341                 }
00342 
00343                 if ( debugLevel_ == pDEBUG )
00344                         std::cout << "Found: " << PeakIndex.size() << " peaks." << std::endl;
00345 
00346                 //Order these peaks by energy:
00347                 for (int i=0;i<int(PeakIndex.size());++i){
00348                         for (int j=0;j<int(PeakIndex.size())-1;++j){
00349                                 if (dominoEnergy[PeakIndex[j]] < dominoEnergy[PeakIndex[j+1]]){
00350                                         int ihold = PeakIndex[j+1];
00351                                         PeakIndex[j+1] = PeakIndex[j];
00352                                         PeakIndex[j] = ihold;
00353                                 }
00354                         }
00355                 }
00356 
00357                 std::vector<int> OwnerShip;
00358                 std::vector<double> LumpEnergy;
00359                 for (int i=0;i<int(dominoEnergy.size());++i) OwnerShip.push_back(-1);
00360 
00361                 //Loop over peaks.  
00362                 double eThres = eThres_;
00363                 double e5x5 = 0.0;
00364                 for (int i = 0; i < int(PeakIndex.size()); ++i)
00365                 {
00366 
00367                         int idxPeak = PeakIndex[i];
00368                         OwnerShip[idxPeak] = i;
00369                         double lump = dominoEnergy[idxPeak];
00370 
00371                         // compute eThres for this peak
00372                         // if set to dynamic (otherwise uncanged from
00373                         // fixed setting
00374                         if (dynamicEThres_) {
00375 
00376                                 //std::cout << "i : " << i << " idxPeak " << idxPeak << std::endl;
00377                                 //std::cout << "    the dominoEnergy.size() = " << dominoEnergy.size() << std::endl;
00378                                 // compute e5x5 for this seed crystal
00379                                 //std::cout << "idxPeak, phiSteps " << idxPeak << ", " << phiSteps << std::endl;
00380                                 e5x5 = lump;
00381                                 //std::cout << "lump " << e5x5 << std::endl;
00382                                 if ((idxPeak + 1) < (int)dominoEnergy.size()) e5x5 += dominoEnergy[idxPeak + 1];
00383                                 //std::cout << "+1 " << e5x5 << std::endl;
00384                                 if ((idxPeak + 2) < (int)dominoEnergy.size()) e5x5 += dominoEnergy[idxPeak + 2];
00385                                 //std::cout << "+2 " << e5x5 << std::endl;
00386                                 if ((idxPeak - 1) > 0) e5x5 += dominoEnergy[idxPeak - 1];
00387                                 //std::cout << "-1 " << e5x5 << std::endl;
00388                                 if ((idxPeak - 2) > 0) e5x5 += dominoEnergy[idxPeak - 2];
00389                                 //std::cout << "-2 " << e5x5 << std::endl;
00390                                 // compute eThres
00391                                 eThres = (eThresA_ * e5x5) + eThresB_;   
00392                                 //std::cout << eThres << std::endl;
00393                                 //std::cout << std::endl;
00394                         }
00395 
00396                         //Loop over adjacent dominos at higher phi
00397                         for (int j=idxPeak+1;j<int(dominoEnergy.size());++j){
00398                                 if (OwnerShip[j]==-1 && 
00399                                                 dominoEnergy[j] > eThres
00400                                                 && dominoEnergy[j] <= dominoEnergy[j-1]){
00401                                         OwnerShip[j]= i;
00402                                         lump+=dominoEnergy[j];
00403                                 }
00404                                 else{
00405                                         break;
00406                                 }
00407                         }
00408                         //loop over adjacent dominos at lower phi.  Sum up energy of lumps.
00409                         for (int j=idxPeak-1;j>=0;--j){
00410                                 if (OwnerShip[j]==-1 && 
00411                                                 dominoEnergy[j] > eThres
00412                                                 && dominoEnergy[j] <= dominoEnergy[j+1]){
00413                                         OwnerShip[j]= i;
00414                                         lump+=dominoEnergy[j];
00415                                 }
00416                                 else{
00417                                         break;
00418                                 }
00419                         }
00420                         LumpEnergy.push_back(lump);
00421                 }
00422 
00423                 //Make the basic clusters:
00424                 for (int i=0;i<int(PeakIndex.size());++i){
00425                         bool HasSeedCrystal = false;
00426                         //One cluster for each peak.
00427                         std::vector<EcalRecHit> recHits;
00428                         std::vector< std::pair<DetId, float> > dets;
00429                         int nhits=0;
00430                         for (int j=0;j<int(dominoEnergy.size());++j){   
00431                                 if (OwnerShip[j] == i){
00432                                         std::vector <EcalRecHit> temp = dominoCells[j];
00433                                         for (int k=0;k<int(temp.size());++k){
00434                                                 dets.push_back( std::pair<DetId, float>(temp[k].id(), 1.) ); // by default energy fractions are 1
00435                                                 if (temp[k].id()==itID)
00436                                                         HasSeedCrystal = true;
00437                                                 recHits.push_back(temp[k]);
00438                                                 nhits++;
00439                                         }
00440                                 }  
00441                         }
00442                         if ( debugLevel_ == pDEBUG ){
00443                                 std::cout << "Adding a cluster with: " << nhits << std::endl;
00444                                 std::cout << "total E: " << LumpEnergy[i] << std::endl;
00445                                 std::cout << "total dets: " << dets.size() << std::endl;
00446                         }
00447 
00448                         //Get Calorimeter position
00449                         Point pos = posCalculator_.Calculate_Location(dets,hits,geometry);
00450 
00451                         //double totChi2=0;
00452                         //double totE=0;
00453                         std::vector<std::pair<DetId, float> > usedHits;
00454                         for (int blarg=0;blarg<int(recHits.size());++blarg){
00455                                 //totChi2 +=0;
00456                                 //totE+=recHits[blarg].energy();
00457                                 usedHits.push_back(std::make_pair<DetId, float>(recHits[blarg].id(), 1.0));
00458                                 useddetids.insert(recHits[blarg].id());
00459                         }
00460 
00461                         //if (totE>0)
00462                         //totChi2/=totE;
00463 
00464                         if (HasSeedCrystal) {
00465                                 // note that this "basiccluster" has the seed crystal of the hyrbid, so record it
00466                                 seedClus_.push_back(reco::BasicCluster(LumpEnergy[i], pos, 
00467                                         reco::CaloID(reco::CaloID::DET_ECAL_BARREL), usedHits, 
00468                                         reco::CaloCluster::hybrid, itID));
00469                                 // and also add to the vector of clusters that will be used in constructing
00470                                 // the supercluster
00471                                 thisseedClusters.push_back(reco::BasicCluster(LumpEnergy[i], pos,                                                               reco::CaloID(reco::CaloID::DET_ECAL_BARREL), 
00472                                      usedHits, reco::CaloCluster::hybrid, itID));
00473                         }
00474                         else {
00475                                 // note that if this "basiccluster" is not the one that seeded the hybrid, 
00476                                 // the seed crystal is unset in this entry in the vector of clusters that will
00477                                 // be used in constructing the super cluster
00478                                 thisseedClusters.push_back(reco::BasicCluster(LumpEnergy[i], pos,                                                               reco::CaloID(reco::CaloID::DET_ECAL_BARREL), 
00479                                      usedHits, reco::CaloCluster::hybrid));
00480                         }
00481                 }
00482 
00483 
00484                 // Make association so that superclusters can be made later.
00485                 // but only if some BasicClusters have been found...
00486                 if (thisseedClusters.size() > 0) 
00487                 {
00488                         clustered_.insert(std::make_pair(clustercounter, thisseedClusters));
00489                         clustercounter++;
00490                 }
00491         }//Seed loop
00492 
00493 }
00494 
00495 reco::SuperClusterCollection HybridClusterAlgo::makeSuperClusters(const reco::CaloClusterPtrVector& clustersCollection)
00496 {
00497         //Here's what we'll return.
00498         reco::SuperClusterCollection SCcoll;
00499 
00500         //Here's our map iterator that gives us the appropriate association.
00501         std::map<int, reco::BasicClusterCollection>::iterator mapit;
00502         for (mapit = clustered_.begin();mapit!=clustered_.end();mapit++){
00503 
00504                 reco::CaloClusterPtrVector thissc;
00505                 reco::CaloClusterPtr seed;//This is not really a seed, but I need to tell SuperCluster something.
00506                 //So I choose the highest energy basiccluster in the SuperCluster.
00507 
00508                 std::vector <reco::BasicCluster> thiscoll = mapit->second; //This is the set of BasicClusters in this
00509                 //SuperCluster
00510 
00511                 double ClusterE = 0; //Sum of cluster energies for supercluster.
00512                 //Holders for position of this supercluster.
00513                 double posX=0;
00514                 double posY=0;
00515                 double posZ=0;
00516 
00517                 //Loop over this set of basic clusters, find their references, and add them to the
00518                 //supercluster.  This could be somehow more efficient.
00519 
00520                 for (int i=0;i<int(thiscoll.size());++i){
00521                         reco::BasicCluster thisclus = thiscoll[i]; //The Cluster in question.
00522                         for (int j=0;j<int(clustersCollection.size());++j){
00523                                 //Find the appropriate cluster from the list of references
00524                                 reco::BasicCluster cluster_p = *clustersCollection[j];
00525                                 if (thisclus== cluster_p){ //Comparison based on energy right now.
00526                                         thissc.push_back(clustersCollection[j]);
00527                                         bool isSeed = false;
00528                                         for (int qu=0;qu<int(seedClus_.size());++qu){
00529                                                 if (cluster_p == seedClus_[qu])
00530                                                         isSeed = true;
00531                                         }
00532                                         if (isSeed) seed = clustersCollection[j];
00533 
00534                                         ClusterE += cluster_p.energy();
00535                                         posX += cluster_p.energy() * cluster_p.position().X();
00536                                         posY += cluster_p.energy() * cluster_p.position().Y();
00537                                         posZ += cluster_p.energy() * cluster_p.position().Z();
00538 
00539                                 }
00540                         }//End loop over finding references.
00541                 }//End loop over clusters.
00542 
00543                 posX /= ClusterE;
00544                 posY /= ClusterE;
00545                 posZ /= ClusterE;
00546 
00547                 /* //This part is moved to EgammaSCEnergyCorrectionAlgo
00548                    double preshowerE = 0.;
00549                    double phiWidth = 0.;
00550                    double etaWidth = 0.;
00551                 //Calculate phiWidth & etaWidth for SuperClusters
00552                 reco::SuperCluster suCl(ClusterE, math::XYZPoint(posX, posY, posZ), seed, thissc, preshowerE, phiWidth, etaWidth);
00553                 SCShape_->Calculate_Covariances(suCl);
00554                 phiWidth = SCShape_->phiWidth();
00555                 etaWidth = SCShape_->etaWidth();
00556                 //Assign phiWidth & etaWidth to SuperCluster as data members
00557                 suCl.setPhiWidth(phiWidth);
00558                 suCl.setEtaWidth(etaWidth);
00559                  */
00560 
00561                 reco::SuperCluster suCl(ClusterE, math::XYZPoint(posX, posY, posZ), seed, thissc);
00562                 SCcoll.push_back(suCl);
00563 
00564                 if ( debugLevel_ == pDEBUG ){
00565                         std::cout << "Super cluster sum: " << ClusterE << std::endl;
00566                         std::cout << "Made supercluster with energy E: " << suCl.energy() << std::endl;
00567                 }
00568         }//end loop over map
00569         sort(SCcoll.rbegin(), SCcoll.rend(), ClusterEtLess());
00570         return SCcoll;
00571 }
00572 
00573 double HybridClusterAlgo::makeDomino(EcalBarrelNavigator &navigator, std::vector <EcalRecHit> &cells)
00574 {
00575         //At the beginning of this function, the navigator starts at the middle of the domino,
00576         //and that's where EcalBarrelNavigator::home() should send it.
00577         //Walk one crystal in eta to either side of the initial point.  Sum the three cell energies.
00578         //If the resultant energy is larger than Ewing, then go an additional cell to either side.
00579         //Returns:  Total domino energy.  Also, stores the cells used to create domino in the vector.
00580         cells.clear();
00581         double Etot = 0;
00582 
00583         //Ready?  Get the starting cell.
00584         DetId center = navigator.pos();
00585         EcalRecHitCollection::const_iterator center_it = recHits_->find(center);
00586 
00587         if (center_it!=recHits_->end()){
00588                 EcalRecHit SeedHit = *center_it;
00589                 if (useddetids.find(center) == useddetids.end() && excludedCrys_.find(center)==excludedCrys_.end()){ 
00590                         Etot += SeedHit.energy();
00591                         cells.push_back(SeedHit);
00592                 }
00593         }
00594         //One step upwards in Ieta:
00595         DetId ieta1 = navigator.west();
00596         EcalRecHitCollection::const_iterator eta1_it = recHits_->find(ieta1);
00597         if (eta1_it !=recHits_->end()){
00598                 EcalRecHit UpEta = *eta1_it;
00599                 if (useddetids.find(ieta1) == useddetids.end() && excludedCrys_.find(ieta1)==excludedCrys_.end()){
00600                         Etot+=UpEta.energy();
00601                         cells.push_back(UpEta);
00602                 }
00603         }
00604 
00605         //Go back to the middle.
00606         navigator.home();
00607 
00608         //One step downwards in Ieta:
00609         DetId ieta2 = navigator.east();
00610         EcalRecHitCollection::const_iterator eta2_it = recHits_->find(ieta2);
00611         if (eta2_it !=recHits_->end()){
00612                 EcalRecHit DownEta = *eta2_it;
00613                 if (useddetids.find(ieta2)==useddetids.end() && excludedCrys_.find(ieta2)==excludedCrys_.end()){
00614                         Etot+=DownEta.energy();
00615                         cells.push_back(DownEta);
00616                 }
00617         }
00618 
00619         //Now check the energy.  If smaller than Ewing, then we're done.  If greater than Ewing, we have to
00620         //add two additional cells, the 'wings'
00621         if (Etot < Ewing) {
00622                 navigator.home(); //Needed even here!!
00623                 return Etot;  //Done!  Not adding 'wings'.
00624         }
00625 
00626         //Add the extra 'wing' cells.  Remember, we haven't sent the navigator home,
00627         //we're still on the DownEta cell.
00628         if (ieta2 != DetId(0)){
00629                 DetId ieta3 = navigator.east(); //Take another step downward.
00630                 EcalRecHitCollection::const_iterator eta3_it = recHits_->find(ieta3);
00631                 if (eta3_it != recHits_->end()){
00632                         EcalRecHit DownEta2 = *eta3_it;
00633                         if (useddetids.find(ieta3)==useddetids.end() && excludedCrys_.find(ieta3)==excludedCrys_.end()){
00634                                 Etot+=DownEta2.energy();
00635                                 cells.push_back(DownEta2);
00636                         }
00637                 }
00638         }
00639         //Now send the navigator home.
00640         navigator.home();
00641         navigator.west(); //Now you're on eta1_it
00642         if (ieta1 != DetId(0)){
00643                 DetId ieta4 = navigator.west(); //Take another step upward.
00644                 EcalRecHitCollection::const_iterator eta4_it = recHits_->find(ieta4);
00645                 if (eta4_it != recHits_->end()){
00646                         EcalRecHit UpEta2 = *eta4_it;
00647                         if (useddetids.find(ieta4) == useddetids.end() && excludedCrys_.find(ieta4)==excludedCrys_.end()){
00648                                 Etot+=UpEta2.energy();
00649                                 cells.push_back(UpEta2);
00650                         }
00651                 }
00652         }
00653         navigator.home();
00654         return Etot;
00655 }
00656 
00657 double HybridClusterAlgo::et25(EcalBarrelNavigator &navigator, 
00658                 const EcalRecHitCollection *hits, 
00659                 const CaloSubdetectorGeometry *geometry)
00660 {
00661 
00662         DetId thisDet;
00663         std::vector< std::pair<DetId, float> > dets;
00664         dets.clear();
00665         EcalRecHitCollection::const_iterator hit;
00666         double energySum = 0.0;
00667 
00668         for (int dx = -2; dx < 3; ++dx)
00669         {
00670                 for (int dy = -2; dy < 3; ++ dy)
00671                 {
00672                         //std::cout << "dx, dy " << dx << ", " << dy << std::endl;
00673                         thisDet = navigator.offsetBy(dx, dy);
00674                         navigator.home();
00675 
00676                         if (thisDet != DetId(0))
00677                         {
00678                                 hit = recHits_->find(thisDet);
00679                                 if (hit != recHits_->end()) 
00680                                 {
00681                                         dets.push_back( std::pair<DetId, float>(thisDet, 1.) ); // by default hit energy fraction is set to 1
00682                                         energySum += hit->energy();
00683                                 }
00684                         }
00685                 }
00686         }
00687 
00688         // convert it to ET
00689         //std::cout << "dets.size(), energySum: " << dets.size() << ", " << energySum << std::endl;
00690         Point pos = posCalculator_.Calculate_Location(dets, hits, geometry);
00691         double et = energySum/cosh(pos.eta());
00692         return et;
00693 
00694 }
00695