CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Multi5x5BremRecoveryClusterAlgo.h
Go to the documentation of this file.
1 #ifndef RecoEcal_EgammaClusterAlgos_Multi5x5BremRecoveryClusterAlgo_h_
2 #define RecoEcal_EgammaClusterAlgos_Multi5x5BremRecoveryClusterAlgo_h_
3 
10 
12 
13 #include <vector>
14 
15 
16 /*
17  The Multi5x5BremRecoveryClusterAlgo class encapsulates the functionality needed
18  to perform the SuperClustering.
19 
20  WARNING: This code assumes that the BasicClusters
21  from the event are sorted by energy
22 */
23 
25 {
26  public:
27 
28  enum VerbosityLevel { pDEBUG = 0, pWARNING = 1, pINFO = 2, pERROR = 3 };
29 
31  double eb_sc_road_etasize = 0.06, // Search window in eta - Barrel
32  double eb_sc_road_phisize = 0.80, // Search window in phi - Barrel
33  double ec_sc_road_etasize = 0.14, // Search window in eta - Endcap
34  double ec_sc_road_phisize = 0.40, // Search window in eta - Endcap
35  bool dynamicPhiRoad = true,
36  double theSeedTransverseEnergyThreshold = 0.40,
37  VerbosityLevel the_verbosity = pERROR
38  )
39  {
40  // e*_rdeta_ and e*_rdphi_ are half the total window
41  // because they correspond to one direction (positive or negative)
42  eb_rdeta_ = eb_sc_road_etasize / 2;
43  eb_rdphi_ = eb_sc_road_phisize / 2;
44  ec_rdeta_ = ec_sc_road_etasize / 2;
45  ec_rdphi_ = ec_sc_road_phisize / 2;
46 
47  seedTransverseEnergyThreshold = theSeedTransverseEnergyThreshold;
48  dynamicPhiRoad_ = dynamicPhiRoad;
49  if (dynamicPhiRoad_) phiRoadAlgo_ = new BremRecoveryPhiRoadAlgo(bremRecoveryPset);
50 
51  verbosity = the_verbosity;
52  }
53 
54  // destructor
56  {
57  if (dynamicPhiRoad_) delete phiRoadAlgo_;
58  }
59 
60  void setVerbosity(VerbosityLevel the_verbosity)
61  {
62  verbosity = the_verbosity;
63  }
64 
65  // the method called from outside to do the SuperClustering - returns a vector of SCs:
67 
68  private:
69 
70  // make superclusters out of clusters produced by the Island algorithm:
72  double etaRoad, double phiRoad);
73 
74  // return true if the cluster is within the search phi-eta window of the seed
75  bool match(reco::CaloClusterPtr seed_p,
76  reco::CaloClusterPtr cluster_p,
77  double etaRoad, double phiRoad);
78 
79  //
80 
82 
83  double eb_rdeta_;
84  double eb_rdphi_;
85  double ec_rdeta_;
86  double ec_rdphi_;
87 
91 
93 
94 };
95 
96 #endif
void makeIslandSuperClusters(reco::CaloClusterPtrVector &clusters_v, double etaRoad, double phiRoad)
reco::SuperClusterCollection makeSuperClusters(reco::CaloClusterPtrVector &clusters)
std::vector< SuperCluster > SuperClusterCollection
collection of SuperCluser objectr
void setVerbosity(VerbosityLevel the_verbosity)
reco::SuperClusterCollection superclusters_v
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)
bool match(reco::CaloClusterPtr seed_p, reco::CaloClusterPtr cluster_p, double etaRoad, double phiRoad)