#include <Multi5x5BremRecoveryClusterAlgo.h>
Public Member Functions | |
reco::SuperClusterCollection | makeSuperClusters (reco::CaloClusterPtrVector &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) | |
~Multi5x5BremRecoveryClusterAlgo () | |
Private Member Functions | |
void | makeIslandSuperClusters (reco::CaloClusterPtrVector &clusters_v, double etaRoad, double phiRoad) |
bool | match (reco::CaloClusterPtr seed_p, reco::CaloClusterPtr 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 |
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 |
||
) | [inline] |
Definition at line 29 of file Multi5x5BremRecoveryClusterAlgo.h.
References dynamicPhiRoad_, eb_rdeta_, eb_rdphi_, ec_rdeta_, ec_rdphi_, phiRoadAlgo_, and seedTransverseEnergyThreshold.
{ // e*_rdeta_ and e*_rdphi_ are half the total window // because they correspond to one direction (positive or negative) eb_rdeta_ = eb_sc_road_etasize / 2; eb_rdphi_ = eb_sc_road_phisize / 2; ec_rdeta_ = ec_sc_road_etasize / 2; ec_rdphi_ = ec_sc_road_phisize / 2; seedTransverseEnergyThreshold = theSeedTransverseEnergyThreshold; dynamicPhiRoad_ = dynamicPhiRoad; if (dynamicPhiRoad_) phiRoadAlgo_ = new BremRecoveryPhiRoadAlgo(bremRecoveryPset); }
Multi5x5BremRecoveryClusterAlgo::~Multi5x5BremRecoveryClusterAlgo | ( | ) | [inline] |
Definition at line 52 of file Multi5x5BremRecoveryClusterAlgo.h.
References dynamicPhiRoad_, and phiRoadAlgo_.
{ if (dynamicPhiRoad_) delete phiRoadAlgo_; }
void Multi5x5BremRecoveryClusterAlgo::makeIslandSuperClusters | ( | reco::CaloClusterPtrVector & | clusters_v, |
double | etaRoad, | ||
double | phiRoad | ||
) | [private] |
Definition at line 42 of file Multi5x5BremRecoveryClusterAlgo.cc.
References edm::PtrVector< T >::begin(), edm::PtrVectorBase::clear(), dynamicPhiRoad_, edm::PtrVector< T >::end(), BremRecoveryPhiRoadAlgo::endcapPhiRoad(), relval_parameters_module::energy, spr::find(), LogTrace, match(), phiRoadAlgo_, edm::PtrVector< T >::push_back(), seedTransverseEnergyThreshold, funct::sin(), and superclusters_v.
Referenced by makeSuperClusters().
{ std::vector<DetId> usedSeedDetIds; usedSeedDetIds.clear(); for (reco::CaloCluster_iterator currentSeed = clusters_v.begin(); currentSeed != clusters_v.end(); ++currentSeed) { // check this seed was not already used if (std::find(usedSeedDetIds.begin(), usedSeedDetIds.end(), (*currentSeed)->seed()) != usedSeedDetIds.end()) continue; // Does our highest energy cluster have high enough energy? // changed this to continue from break (to be robust against the order of sorting of the seed clusters) if ((*currentSeed)->energy() * sin((*currentSeed)->position().theta()) < seedTransverseEnergyThreshold) continue; // if yes, make it a seed for a new SuperCluster: double energy = (*currentSeed)->energy(); math::XYZVector position_((*currentSeed)->position().X(), (*currentSeed)->position().Y(), (*currentSeed)->position().Z()); position_ *= energy; usedSeedDetIds.push_back((*currentSeed)->seed()); LogTrace("EcalClusters") << "*****************************"; LogTrace("EcalClusters") << "******NEW SUPERCLUSTER*******"; LogTrace("EcalClusters") << "Seed R = " << (*currentSeed)->position().Rho(); reco::CaloClusterPtrVector constituentClusters; constituentClusters.push_back(*currentSeed); reco::CaloCluster_iterator currentCluster = currentSeed + 1; while (currentCluster != clusters_v.end()) { // if dynamic phi road is enabled then compute the phi road for a cluster // of energy of existing clusters + the candidate cluster. if (dynamicPhiRoad_) phiRoad = phiRoadAlgo_->endcapPhiRoad(energy + (*currentCluster)->energy()); // does the cluster match the phi road for this candidate supercluster if (match(*currentSeed, *currentCluster, etaRoad, phiRoad) && std::find(usedSeedDetIds.begin(), usedSeedDetIds.end(), (*currentCluster)->seed()) == usedSeedDetIds.end()) { // add basic cluster to supercluster constituents constituentClusters.push_back(*currentCluster); energy += (*currentCluster)->energy(); position_ += (*currentCluster)->energy() * math::XYZVector((*currentCluster)->position().X(), (*currentCluster)->position().Y(), (*currentCluster)->position().Z()); // remove cluster from vector of available clusters usedSeedDetIds.push_back((*currentCluster)->seed()); LogTrace("EcalClusters") << "Cluster R = " << (*currentCluster)->position().Rho(); } ++currentCluster; } position_ /= energy; LogTrace("EcalClusters") <<"Final SuperCluster R = " << position_.Rho(); reco::SuperCluster newSuperCluster(energy, math::XYZPoint(position_.X(), position_.Y(), position_.Z()), (*currentSeed), constituentClusters); superclusters_v.push_back(newSuperCluster); LogTrace("EcalClusters") << "created a new supercluster of: "; LogTrace("EcalClusters") << "Energy = " << newSuperCluster.energy(); LogTrace("EcalClusters") << "Position in (R, phi, theta) = (" << newSuperCluster.position().Rho() << ", " << newSuperCluster.position().phi() << ", " << newSuperCluster.position().theta() << ")" ; } clusters_v.clear(); }
reco::SuperClusterCollection Multi5x5BremRecoveryClusterAlgo::makeSuperClusters | ( | reco::CaloClusterPtrVector & | clusters | ) |
Definition at line 5 of file Multi5x5BremRecoveryClusterAlgo.cc.
References edm::PtrVector< T >::begin(), eb_rdeta_, eb_rdphi_, ec_rdeta_, ec_rdphi_, edm::PtrVector< T >::end(), makeIslandSuperClusters(), reco::CaloCluster::multi5x5, edm::PtrVector< T >::push_back(), and superclusters_v.
Referenced by Multi5x5SuperClusterProducer::produceSuperclustersForECALPart().
{ const float etaBorder = 1.479; superclusters_v.clear(); // create vectors of references to clusters of a specific origin... reco::CaloClusterPtrVector islandClustersBarrel_v; reco::CaloClusterPtrVector islandClustersEndCap_v; // ...and populate them: for (reco::CaloCluster_iterator it = clustersCollection.begin(); it != clustersCollection.end(); it++) { reco::CaloClusterPtr cluster_p = *it; if (cluster_p->algo() == reco::CaloCluster::multi5x5) { if (fabs(cluster_p->position().eta()) < etaBorder) { islandClustersBarrel_v.push_back(cluster_p); } else { islandClustersEndCap_v.push_back(cluster_p); } } } // make the superclusters from the Barrel clusters - Island makeIslandSuperClusters(islandClustersBarrel_v, eb_rdeta_, eb_rdphi_); // make the superclusters from the EndCap clusters - Island makeIslandSuperClusters(islandClustersEndCap_v, ec_rdeta_, ec_rdphi_); return superclusters_v; }
bool Multi5x5BremRecoveryClusterAlgo::match | ( | reco::CaloClusterPtr | seed_p, |
reco::CaloClusterPtr | cluster_p, | ||
double | etaRoad, | ||
double | phiRoad | ||
) | [private] |
Definition at line 129 of file Multi5x5BremRecoveryClusterAlgo.cc.
References funct::cos(), and dPhi().
Referenced by makeIslandSuperClusters().
{ math::XYZPoint clusterPosition = cluster_p->position(); math::XYZPoint seedPosition = seed_p->position(); double dPhi = acos(cos(seedPosition.phi() - clusterPosition.phi())); double dEta = fabs(seedPosition.eta() - clusterPosition.eta()); if (dEta > dEtaMax) return false; if (dPhi > dPhiMax) return false; return true; }
bool Multi5x5BremRecoveryClusterAlgo::dynamicPhiRoad_ [private] |
Definition at line 80 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), Multi5x5BremRecoveryClusterAlgo(), and ~Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::eb_rdeta_ [private] |
Definition at line 74 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::eb_rdphi_ [private] |
Definition at line 75 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::ec_rdeta_ [private] |
Definition at line 76 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::ec_rdphi_ [private] |
Definition at line 77 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
Definition at line 81 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), Multi5x5BremRecoveryClusterAlgo(), and ~Multi5x5BremRecoveryClusterAlgo().
double Multi5x5BremRecoveryClusterAlgo::seedTransverseEnergyThreshold [private] |
Definition at line 79 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), and Multi5x5BremRecoveryClusterAlgo().
Definition at line 83 of file Multi5x5BremRecoveryClusterAlgo.h.
Referenced by makeIslandSuperClusters(), and makeSuperClusters().