#include <RecoEcal/EgammaClusterAlgos/interface/Multi5x5BremRecoveryClusterAlgo.h>
Public Types | |
enum | VerbosityLevel { pDEBUG = 0, pWARNING = 1, pINFO = 2, pERROR = 3 } |
Public Member Functions | |
reco::SuperClusterCollection | makeSuperClusters (reco::BasicClusterRefVector &clusters) |
Multi5x5BremRecoveryClusterAlgo (const edm::ParameterSet &bremRecoveryPset, double eb_sc_road_etasize=0.06, double eb_sc_road_phisize=0.80, double ec_sc_road_etasize=0.14, double ec_sc_road_phisize=0.40, bool dynamicPhiRoad=true, double theSeedTransverseEnergyThreshold=0.40, VerbosityLevel the_verbosity=pERROR) | |
void | setVerbosity (VerbosityLevel the_verbosity) |
~Multi5x5BremRecoveryClusterAlgo () | |
Private Member Functions | |
void | makeIslandSuperClusters (reco::BasicClusterRefVector &clusters_v, double etaRoad, double phiRoad) |
bool | match (reco::BasicClusterRef seed_p, reco::BasicClusterRef cluster_p, double etaRoad, double phiRoad) |
Private Attributes | |
bool | dynamicPhiRoad_ |
double | eb_rdeta_ |
double | eb_rdphi_ |
double | ec_rdeta_ |
double | ec_rdphi_ |
BremRecoveryPhiRoadAlgo * | phiRoadAlgo_ |
double | seedTransverseEnergyThreshold |
reco::SuperClusterCollection | superclusters_v |
VerbosityLevel | verbosity |
Definition at line 24 of file Multi5x5BremRecoveryClusterAlgo.h.
Multi5x5BremRecoveryClusterAlgo::Multi5x5BremRecoveryClusterAlgo | ( | const edm::ParameterSet & | bremRecoveryPset, | |
double | eb_sc_road_etasize = 0.06 , |
|||
double | eb_sc_road_phisize = 0.80 , |
|||
double | ec_sc_road_etasize = 0.14 , |
|||
double | ec_sc_road_phisize = 0.40 , |
|||
bool | dynamicPhiRoad = true , |
|||
double | theSeedTransverseEnergyThreshold = 0.40 , |
|||
VerbosityLevel | the_verbosity = pERROR | |||
) | [inline] |
Definition at line 30 of file Multi5x5BremRecoveryClusterAlgo.h.
References dynamicPhiRoad_, eb_rdeta_, eb_rdphi_, ec_rdeta_, ec_rdphi_, phiRoadAlgo_, seedTransverseEnergyThreshold, and verbosity.
00039 { 00040 // e*_rdeta_ and e*_rdphi_ are half the total window 00041 // because they correspond to one direction (positive or negative) 00042 eb_rdeta_ = eb_sc_road_etasize / 2; 00043 eb_rdphi_ = eb_sc_road_phisize / 2; 00044 ec_rdeta_ = ec_sc_road_etasize / 2; 00045 ec_rdphi_ = ec_sc_road_phisize / 2; 00046 00047 seedTransverseEnergyThreshold = theSeedTransverseEnergyThreshold; 00048 dynamicPhiRoad_ = dynamicPhiRoad; 00049 if (dynamicPhiRoad_) phiRoadAlgo_ = new BremRecoveryPhiRoadAlgo(bremRecoveryPset); 00050 00051 verbosity = the_verbosity; 00052 }
Multi5x5BremRecoveryClusterAlgo::~Multi5x5BremRecoveryClusterAlgo | ( | ) | [inline] |
Definition at line 55 of file Multi5x5BremRecoveryClusterAlgo.h.
References dynamicPhiRoad_, and phiRoadAlgo_.
00056 { 00057 if (dynamicPhiRoad_) delete phiRoadAlgo_; 00058 }
void Multi5x5BremRecoveryClusterAlgo::makeIslandSuperClusters | ( | reco::BasicClusterRefVector & | clusters_v, | |
double | etaRoad, | |||
double | phiRoad | |||
) | [private] |
Definition at line 42 of file Multi5x5BremRecoveryClusterAlgo.cc.
References edm::RefVector< C, T, F >::begin(), GenMuonPlsPt100GeV_cfg::cout, dynamicPhiRoad_, edm::RefVector< C, T, F >::empty(), edm::RefVector< C, T, F >::end(), BremRecoveryPhiRoadAlgo::endcapPhiRoad(), lat::endl(), relval_parameters_module::energy, edm::RefVector< C, T, F >::erase(), match(), phiRoadAlgo_, pINFO, edm::RefVector< C, T, F >::push_back(), seedTransverseEnergyThreshold, funct::sin(), superclusters_v, and verbosity.
Referenced by makeSuperClusters().
00044 { 00045 00046 reco::basicCluster_iterator currentSeed; 00047 for (currentSeed = clusters_v.begin(); !clusters_v.empty(); clusters_v.erase(currentSeed)) 00048 { 00049 // Does our highest energy cluster have high enough energy? 00050 // changed this to continue from break (to be robust against the order of sorting of the seed clusters) 00051 if ((*currentSeed)->energy() * sin((*currentSeed)->position().theta()) < seedTransverseEnergyThreshold) continue; 00052 00053 // if yes, make it a seed for a new SuperCluster: 00054 double energy = (*currentSeed)->energy(); 00055 math::XYZVector position_((*currentSeed)->position().X(), 00056 (*currentSeed)->position().Y(), 00057 (*currentSeed)->position().Z()); 00058 position_ *= energy; 00059 00060 00061 if (verbosity <= pINFO) 00062 { 00063 std::cout << "*****************************" << std::endl; 00064 std::cout << "******NEW SUPERCLUSTER*******" << std::endl; 00065 std::cout << "Seed R = " << (*currentSeed)->position().Rho() << std::endl; 00066 } 00067 00068 // and add the matching clusters: 00069 reco::BasicClusterRefVector constituentClusters; 00070 constituentClusters.push_back(*currentSeed); 00071 reco::basicCluster_iterator currentCluster = currentSeed + 1; 00072 while (currentCluster != clusters_v.end()) 00073 { 00074 00075 // if dynamic phi road is enabled then compute the phi road for a cluster 00076 // of energy of existing clusters + the candidate cluster. 00077 if (dynamicPhiRoad_) { 00078 phiRoad = phiRoadAlgo_->endcapPhiRoad(energy + (*currentCluster)->energy()); 00079 } 00080 00081 if (match(*currentSeed, *currentCluster, etaRoad, phiRoad)) 00082 { 00083 constituentClusters.push_back(*currentCluster); 00084 energy += (*currentCluster)->energy(); 00085 position_ += (*currentCluster)->energy() * math::XYZVector((*currentCluster)->position().X(), 00086 (*currentCluster)->position().Y(), 00087 (*currentCluster)->position().Z()); 00088 if (verbosity <= pINFO) 00089 { 00090 std::cout << "Cluster R = " << (*currentCluster)->position().Rho() << std::endl; 00091 } 00092 00093 clusters_v.erase(currentCluster); 00094 } 00095 else currentCluster++; 00096 } 00097 00098 position_ /= energy; 00099 00100 if (verbosity <= pINFO) 00101 { 00102 std::cout << "Final SuperCluster R = " << position_.Rho() << std::endl; 00103 } 00104 00105 reco::SuperCluster newSuperCluster(energy, 00106 math::XYZPoint(position_.X(), position_.Y(), position_.Z()), 00107 (*currentSeed), 00108 constituentClusters); 00109 00110 superclusters_v.push_back(newSuperCluster); 00111 00112 if (verbosity <= pINFO) 00113 { 00114 std::cout << "created a new supercluster of: " << std::endl; 00115 std::cout << "Energy = " << newSuperCluster.energy() << std::endl; 00116 std::cout << "Position in (R, phi, theta) = (" 00117 << newSuperCluster.position().Rho() << ", " 00118 << newSuperCluster.position().phi() << ", " 00119 << newSuperCluster.position().theta() << ")" << std::endl; 00120 } 00121 } 00122 currentSeed = clusters_v.end(); 00123 }
reco::SuperClusterCollection Multi5x5BremRecoveryClusterAlgo::makeSuperClusters | ( | reco::BasicClusterRefVector & | clusters | ) |
Definition at line 4 of file Multi5x5BremRecoveryClusterAlgo.cc.
References edm::RefVector< C, T, F >::begin(), eb_rdeta_, eb_rdphi_, ec_rdeta_, ec_rdphi_, edm::RefVector< C, T, F >::end(), reco::island, it, makeIslandSuperClusters(), edm::RefVector< C, T, F >::push_back(), and superclusters_v.
Referenced by Multi5x5SuperClusterProducer::produceSuperclustersForECALPart().
00005 { 00006 const float etaBorder = 1.479; 00007 00008 superclusters_v.clear(); 00009 00010 // create vectors of references to clusters of a specific origin... 00011 reco::BasicClusterRefVector islandClustersBarrel_v; 00012 reco::BasicClusterRefVector islandClustersEndCap_v; 00013 00014 // ...and populate them: 00015 reco::basicCluster_iterator it; 00016 for (it = clustersCollection.begin(); it != clustersCollection.end(); it++) 00017 { 00018 reco::BasicClusterRef cluster_p = *it; 00019 if (cluster_p->algo() == reco::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 }
bool Multi5x5BremRecoveryClusterAlgo::match | ( | reco::BasicClusterRef | seed_p, | |
reco::BasicClusterRef | cluster_p, | |||
double | etaRoad, | |||
double | phiRoad | |||
) | [private] |
Definition at line 126 of file Multi5x5BremRecoveryClusterAlgo.cc.
References funct::cos(), and dPhi().
Referenced by makeIslandSuperClusters().
00129 { 00130 math::XYZPoint clusterPosition = cluster_p->position(); 00131 math::XYZPoint seedPosition = seed_p->position(); 00132 00133 double dPhi = acos(cos(seedPosition.phi() - clusterPosition.phi())); 00134 00135 double dEta = fabs(seedPosition.eta() - clusterPosition.eta()); 00136 00137 if (dEta > dEtaMax) return false; 00138 if (dPhi > dPhiMax) return false; 00139 00140 return true; 00141 }
void Multi5x5BremRecoveryClusterAlgo::setVerbosity | ( | VerbosityLevel | the_verbosity | ) | [inline] |
Definition at line 60 of file Multi5x5BremRecoveryClusterAlgo.h.
References verbosity.
00061 { 00062 verbosity = the_verbosity; 00063 }
Definition at line 89 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), Multi5x5BremRecoveryClusterAlgo(), and ~Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::eb_rdeta_ [private] |
Definition at line 83 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::eb_rdphi_ [private] |
Definition at line 84 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::ec_rdeta_ [private] |
Definition at line 85 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::ec_rdphi_ [private] |
Definition at line 86 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
Definition at line 90 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), Multi5x5BremRecoveryClusterAlgo(), and ~Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::seedTransverseEnergyThreshold [private] |
Definition at line 88 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
Definition at line 92 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), and makeSuperClusters().
Definition at line 81 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), Multi5x5BremRecoveryClusterAlgo(), and setVerbosity().