CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoEcal/EgammaClusterAlgos/src/BremRecoveryClusterAlgo.cc

Go to the documentation of this file.
00001 #include "RecoEcal/EgammaClusterAlgos/interface/BremRecoveryClusterAlgo.h"
00002 
00003 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
00004 
00005 reco::SuperClusterCollection BremRecoveryClusterAlgo::makeSuperClusters(reco::CaloClusterPtrVector & clustersCollection)
00006 {
00007   const float etaBorder = 1.479;
00008 
00009   superclusters_v.clear();
00010   
00011   // create vectors of references to clusters of a specific origin...
00012   reco::CaloClusterPtrVector islandClustersBarrel_v;
00013   reco::CaloClusterPtrVector islandClustersEndCap_v;
00014 
00015   // ...and populate them:
00016   for (reco::CaloCluster_iterator it = clustersCollection.begin(); it != clustersCollection.end(); it++)
00017     {
00018       reco::CaloClusterPtr cluster_p = *it;
00019       if (cluster_p->algo() == reco::CaloCluster::island) 
00020       {
00021           if (fabs(cluster_p->position().eta()) < etaBorder)
00022             {
00023               islandClustersBarrel_v.push_back(cluster_p);
00024             }
00025           else
00026             {
00027               islandClustersEndCap_v.push_back(cluster_p);
00028             }
00029       }
00030     }
00031 
00032   // make the superclusters from the Barrel clusters - Island
00033   makeIslandSuperClusters(islandClustersBarrel_v, eb_rdeta_, eb_rdphi_);
00034   // make the superclusters from the EndCap clusters - Island
00035   makeIslandSuperClusters(islandClustersEndCap_v, ec_rdeta_, ec_rdphi_);
00036 
00037   return superclusters_v;
00038 }
00039 
00040 #include "DataFormats/Math/interface/Vector3D.h"
00041 
00042 void BremRecoveryClusterAlgo::makeIslandSuperClusters(reco::CaloClusterPtrVector &clusters_v, 
00043                                                       double etaRoad, double phiRoad)
00044 {
00045 
00046   for ( reco::CaloCluster_iterator currentSeed = clusters_v.begin(); currentSeed != clusters_v.end(); ++currentSeed)
00047     {
00048       // Does our highest energy cluster have high enough energy?
00049       if ((*currentSeed)->energy() * sin((*currentSeed)->position().theta()) < seedTransverseEnergyThreshold) break;
00050 
00051       // if yes, make it a seed for a new SuperCluster:
00052       double energy_ = (*currentSeed)->energy();
00053       math::XYZVector position_((*currentSeed)->position().X(), 
00054                                 (*currentSeed)->position().Y(), 
00055                                 (*currentSeed)->position().Z());
00056       position_ *= energy_;
00057 
00058       if (verbosity <= pINFO)
00059         {
00060           std::cout << "*****************************" << std::endl;
00061           std::cout << "******NEW SUPERCLUSTER*******" << std::endl;
00062           std::cout << "Seed R = " << (*currentSeed)->position().Rho() << std::endl;
00063         }
00064 
00065       // and add the matching clusters:
00066       reco::CaloClusterPtrVector constituentClusters;
00067       constituentClusters.push_back( *currentSeed );
00068       reco::CaloCluster_iterator currentCluster = currentSeed + 1;
00069       while (currentCluster != clusters_v.end())
00070         {
00071           if (match(*currentSeed, *currentCluster, etaRoad, phiRoad))
00072             {
00073               constituentClusters.push_back(*currentCluster);
00074               energy_   += (*currentCluster)->energy();
00075               position_ += (*currentCluster)->energy() * math::XYZVector((*currentCluster)->position().X(), 
00076                                                                          (*currentCluster)->position().Y(), 
00077                                                                          (*currentCluster)->position().Z()); 
00078               if (verbosity <= pINFO) 
00079                 {
00080                   std::cout << "Cluster R = " << (*currentCluster)->position().Rho() << std::endl;
00081                 }
00082 
00083             }
00084           ++currentCluster;
00085         }
00086 
00087       position_ /= energy_;
00088 
00089       if (verbosity <= pINFO)
00090         {
00091           std::cout << "Final SuperCluster R = " << position_.Rho() << std::endl;
00092         }
00093 
00094       reco::SuperCluster newSuperCluster(energy_, 
00095                                          math::XYZPoint(position_.X(), position_.Y(), position_.Z()), 
00096                                          (*currentSeed), 
00097                                          constituentClusters);
00098 
00099       superclusters_v.push_back(newSuperCluster);
00100 
00101       if (verbosity <= pINFO)
00102         {
00103           std::cout << "created a new supercluster of: " << std::endl;
00104           std::cout << "Energy = " << newSuperCluster.energy() << std::endl;
00105           std::cout << "Position in (R, phi, theta) = (" 
00106                     << newSuperCluster.position().Rho() << ", " 
00107                     << newSuperCluster.position().phi() << ", "
00108                     << newSuperCluster.position().theta() << ")" << std::endl;
00109         }
00110     }
00111   clusters_v.clear();
00112 }
00113 
00114 
00115 bool BremRecoveryClusterAlgo::match(reco::CaloClusterPtr seed_p, 
00116                                     reco::CaloClusterPtr cluster_p,
00117                                     double dEtaMax, double dPhiMax)
00118 {
00119   math::XYZPoint clusterPosition = cluster_p->position();
00120   math::XYZPoint seedPosition = seed_p->position();
00121 
00122   double dPhi = acos(cos(seedPosition.phi() - clusterPosition.phi()));
00123  
00124   double dEta = fabs(seedPosition.eta() - clusterPosition.eta());
00125  
00126   if (dEta > dEtaMax) return false;
00127   if (dPhi > dPhiMax) return false;
00128 
00129   return true;
00130 }