CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_2_9_HLT1_bphpatch4/src/RecoEcal/EgammaClusterAlgos/src/Multi5x5BremRecoveryClusterAlgo.cc

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