CMS 3D CMS Logo

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