CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/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 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00004 
00005 reco::SuperClusterCollection Multi5x5BremRecoveryClusterAlgo::makeSuperClusters(reco::CaloClusterPtrVector & clustersCollection)
00006 {
00007 
00008         const float etaBorder = 1.479;
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::multi5x5) 
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 Multi5x5BremRecoveryClusterAlgo::makeIslandSuperClusters(reco::CaloClusterPtrVector &clusters_v, 
00043                 double etaRoad, double phiRoad)
00044 {
00045 
00046         std::vector<DetId> usedSeedDetIds;
00047         usedSeedDetIds.clear();
00048 
00049         for (reco::CaloCluster_iterator currentSeed = clusters_v.begin(); currentSeed != clusters_v.end(); ++currentSeed)
00050         {
00051 
00052                 // check this seed was not already used
00053                 if (std::find(usedSeedDetIds.begin(), usedSeedDetIds.end(), (*currentSeed)->seed()) 
00054                         != usedSeedDetIds.end()) continue;
00055 
00056                 // Does our highest energy cluster have high enough energy?
00057                 // changed this to continue from break (to be robust against the order of sorting of the seed clusters)
00058                 if ((*currentSeed)->energy() * sin((*currentSeed)->position().theta()) < seedTransverseEnergyThreshold) continue;
00059 
00060                 // if yes, make it a seed for a new SuperCluster:
00061                 double energy = (*currentSeed)->energy();
00062                 math::XYZVector position_((*currentSeed)->position().X(), 
00063                                 (*currentSeed)->position().Y(), 
00064                                 (*currentSeed)->position().Z());
00065                 position_ *= energy;
00066                 usedSeedDetIds.push_back((*currentSeed)->seed());
00067 
00068 LogTrace("EcalClusters") << "*****************************";
00069 LogTrace("EcalClusters") << "******NEW SUPERCLUSTER*******";
00070 LogTrace("EcalClusters") << "Seed R = " << (*currentSeed)->position().Rho();
00071 
00072                 reco::CaloClusterPtrVector constituentClusters;
00073                 constituentClusters.push_back(*currentSeed);
00074                 reco::CaloCluster_iterator currentCluster = currentSeed + 1;
00075 
00076                 while (currentCluster != clusters_v.end())
00077                 {
00078 
00079                         // if dynamic phi road is enabled then compute the phi road for a cluster
00080                         // of energy of existing clusters + the candidate cluster.
00081                         if (dynamicPhiRoad_)
00082                                 phiRoad = phiRoadAlgo_->endcapPhiRoad(energy + (*currentCluster)->energy());
00083 
00084                         // does the cluster match the phi road for this candidate supercluster
00085                         if (match(*currentSeed, *currentCluster, etaRoad, phiRoad) &&
00086                                 std::find(usedSeedDetIds.begin(), usedSeedDetIds.end(), (*currentCluster)->seed()) == usedSeedDetIds.end())
00087                         {
00088 
00089                                 // add basic cluster to supercluster constituents
00090                                 constituentClusters.push_back(*currentCluster);
00091                                 energy   += (*currentCluster)->energy();
00092                                 position_ += (*currentCluster)->energy() * math::XYZVector((*currentCluster)->position().X(),
00093                                                 (*currentCluster)->position().Y(), 
00094                                                 (*currentCluster)->position().Z());
00095 
00096                                 // remove cluster from vector of available clusters
00097                                 usedSeedDetIds.push_back((*currentCluster)->seed());
00098 LogTrace("EcalClusters") << "Cluster R = " << (*currentCluster)->position().Rho();
00099                         }
00100                         ++currentCluster;
00101 
00102                 }
00103 
00104                 position_ /= energy;
00105 
00106 LogTrace("EcalClusters") <<"Final SuperCluster R = " << position_.Rho();
00107 
00108                 reco::SuperCluster newSuperCluster(energy, 
00109                                 math::XYZPoint(position_.X(), position_.Y(), position_.Z()), 
00110                                 (*currentSeed), 
00111                                 constituentClusters);
00112 
00113                 superclusters_v.push_back(newSuperCluster);
00114 LogTrace("EcalClusters") << "created a new supercluster of: ";
00115 LogTrace("EcalClusters") << "Energy = " << newSuperCluster.energy();
00116 LogTrace("EcalClusters") << "Position in (R, phi, theta) = ("
00117                                 << newSuperCluster.position().Rho() << ", "
00118                                 << newSuperCluster.position().phi() << ", "
00119                                 << newSuperCluster.position().theta() << ")" ;
00120 
00121 
00122         }
00123 
00124         clusters_v.clear();
00125 
00126 }
00127 
00128 
00129 bool Multi5x5BremRecoveryClusterAlgo::match(reco::CaloClusterPtr seed_p, 
00130                 reco::CaloClusterPtr cluster_p,
00131                 double dEtaMax, double dPhiMax)
00132 {
00133         math::XYZPoint clusterPosition = cluster_p->position();
00134         math::XYZPoint seedPosition = seed_p->position();
00135 
00136         double dPhi = acos(cos(seedPosition.phi() - clusterPosition.phi()));
00137 
00138         double dEta = fabs(seedPosition.eta() - clusterPosition.eta());
00139 
00140         if (dEta > dEtaMax) return false;
00141         if (dPhi > dPhiMax) return false;
00142 
00143         return true;
00144 }