Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-10-25 09:58:57

0001 #ifndef RecoEcal_EgammaClusterAlgos_Multi5x5BremRecoveryClusterAlgo_h_
0002 #define RecoEcal_EgammaClusterAlgos_Multi5x5BremRecoveryClusterAlgo_h_
0003 
0004 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0005 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
0006 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
0007 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
0008 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
0009 #include "RecoEcal/EgammaCoreTools/interface/BremRecoveryPhiRoadAlgo.h"
0010 
0011 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0012 
0013 #include <vector>
0014 
0015 /*
0016   The Multi5x5BremRecoveryClusterAlgo class encapsulates the functionality needed
0017   to perform the SuperClustering.
0018   
0019   WARNING: This code assumes that the BasicClusters 
0020   from the event are sorted by energy
0021 */
0022 
0023 class Multi5x5BremRecoveryClusterAlgo {
0024 public:
0025   Multi5x5BremRecoveryClusterAlgo(const edm::ParameterSet &bremRecoveryPset,
0026                                   double eb_sc_road_etasize = 0.06,  // Search window in eta - Barrel
0027                                   double eb_sc_road_phisize = 0.80,  // Search window in phi - Barrel
0028                                   double ec_sc_road_etasize = 0.14,  // Search window in eta - Endcap
0029                                   double ec_sc_road_phisize = 0.40,  // Search window in eta - Endcap
0030                                   bool dynamicPhiRoad = true,
0031                                   double theSeedTransverseEnergyThreshold = 0.40) {
0032     // e*_rdeta_ and e*_rdphi_ are half the total window
0033     // because they correspond to one direction (positive or negative)
0034     eb_rdeta_ = eb_sc_road_etasize / 2;
0035     eb_rdphi_ = eb_sc_road_phisize / 2;
0036     ec_rdeta_ = ec_sc_road_etasize / 2;
0037     ec_rdphi_ = ec_sc_road_phisize / 2;
0038 
0039     seedTransverseEnergyThreshold = theSeedTransverseEnergyThreshold;
0040     dynamicPhiRoad_ = dynamicPhiRoad;
0041     if (dynamicPhiRoad_)
0042       phiRoadAlgo_ = new BremRecoveryPhiRoadAlgo(bremRecoveryPset);
0043   }
0044 
0045   // destructor
0046   ~Multi5x5BremRecoveryClusterAlgo() {
0047     if (dynamicPhiRoad_)
0048       delete phiRoadAlgo_;
0049   }
0050 
0051   // the method called from outside to do the SuperClustering - returns a vector of SCs:
0052   reco::SuperClusterCollection makeSuperClusters(const reco::CaloClusterPtrVector &clusters);
0053 
0054 private:
0055   // make superclusters out of clusters produced by the Island algorithm:
0056   void makeIslandSuperClusters(reco::CaloClusterPtrVector &clusters_v, double etaRoad, double phiRoad);
0057 
0058   double eb_rdeta_;
0059   double eb_rdphi_;
0060   double ec_rdeta_;
0061   double ec_rdphi_;
0062 
0063   double seedTransverseEnergyThreshold;
0064   bool dynamicPhiRoad_;
0065   BremRecoveryPhiRoadAlgo *phiRoadAlgo_;
0066 
0067   reco::SuperClusterCollection superclusters_v;
0068 };
0069 
0070 #endif