CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/RecoEcal/EgammaClusterAlgos/interface/HybridClusterAlgo.h

Go to the documentation of this file.
00001 #ifndef RecoEcal_EgammaClusterAlgos_HybridClusterAlgo_h
00002 #define RecoEcal_EgammaClusterAlgos_HybridClusterAlgo_h
00003 
00004 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00005 #include "DataFormats/EgammaReco/interface/BasicClusterFwd.h"
00006 #include "DataFormats/EgammaReco/interface/SuperCluster.h"
00007 #include "DataFormats/EgammaReco/interface/SuperClusterFwd.h"
00008 #include "RecoCaloTools/Navigation/interface/EcalBarrelNavigator.h"
00009 #include "DataFormats/EcalRecHit/interface/EcalRecHit.h"
00010 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00011 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00012 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00013 #include "Geometry/CaloTopology/interface/EcalBarrelHardcodedTopology.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "RecoEcal/EgammaCoreTools/interface/EcalEtaPhiRegion.h"
00016 #include "RecoEcal/EgammaCoreTools/interface/PositionCalc.h"
00017 #include "RecoEcal/EgammaCoreTools/interface/BremRecoveryPhiRoadAlgo.h"
00018 #include "RecoEcal/EgammaCoreTools/interface/SuperClusterShapeAlgo.h"
00019 
00020 
00021 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00022 
00023 #include <vector>
00024 #include <set>
00025 
00026 class EcalSeverityLevelAlgo;
00027 
00028 //Less than operator for sorting EcalRecHits according to energy.
00029 struct less_mag : public std::binary_function<EcalRecHit, EcalRecHit, bool> {
00030   bool operator()(EcalRecHit x, EcalRecHit y) { return x.energy() > y.energy() ; }
00031 };
00032 
00033 class HybridClusterAlgo
00034 {
00035  private:
00036   //Quick typedef for position calculation.
00037   typedef math::XYZPoint Point;
00038 
00039   //Thresholds for seeds.
00040   double eb_st;
00041 
00042   //Number of steps in phi that the Hybrid algorithm will take
00043   //when clustering.  Remember, uses phiSteps_ in positive direction
00044   //and then phiSteps_ in negative direction.
00045   int phiSteps_;
00046 
00047   // et in 25
00048   double et25(EcalBarrelNavigator &navigator,
00049                 const EcalRecHitCollection *hits,
00050                 const CaloSubdetectorGeometry *geometry);
00051   // ratio Et/e
00052   double e2Et(EcalBarrelNavigator &navigator,
00053                 const EcalRecHitCollection *hits,
00054                 const CaloSubdetectorGeometry *geometry);
00055 
00056 
00057   BremRecoveryPhiRoadAlgo *phiRoadAlgo_;
00058 
00059   //Threshold for basic cluster.
00060   double eThres_;
00061   double eThresA_;
00062   double eThresB_;
00063 
00064   //Threshold for becoming a sub-peak in the supercluster.
00065   double Eseed;
00066 
00067   //Coefficient to increase Eseed as a function of 5x5  
00068   double Xi;
00069 
00070   //Increase Eseed as a function of et_5x5 (othwewise it's e_5x5)  
00071   bool UseEtForXi;
00072 
00073   //Threshold for adding the additional two 'wing' cells to domino. 
00074   double Ewing;
00075 
00076   // do dynamic phi road
00077   bool dynamicPhiRoad_;
00078 
00079   // do dynamic ethres
00080   bool dynamicEThres_;
00081 
00082   //Map of DetId, RecHit relationship.  EcalRecHit knows what DetId it is,
00083   //but DetId doesn't  know what EcalRecHit it is. 
00084   //  std::map<DetId, EcalRecHit>  rechits_m;
00085 
00086   // colection of all rechits
00087   const EcalRecHitCollection *recHits_;
00088 
00089   // topology
00090   EcalBarrelHardcodedTopology *topo_;
00091 
00092   //  SuperClusterShapeAlgo* SCShape_;
00093 
00094   //Set of DetIds that have already been used.
00095   std::set<DetId> useddetids;
00096 
00097   // The vector of seeds:
00098   std::vector<EcalRecHit> seeds;
00099 
00100   //The vector of seed clusters:
00101   std::vector<reco::BasicCluster> seedClus_;
00102 
00103   //Map of basicclusters and what supercluster they will belong to.
00104   std::map<int, std::vector<reco::BasicCluster> > clustered_;
00105 
00106   //algo to calulate position of clusters
00107   PositionCalc posCalculator_;
00108   
00109   // channels not to be used for seeding 
00110   std::vector<int> v_chstatus_; 
00111 
00112   // severity levels to discriminate against
00113   std::vector<int> v_severitylevel_;
00114   float severityRecHitThreshold_;
00115   float severitySpikeThreshold_;
00116 
00117   bool excludeFromCluster_;
00118   std::set<DetId> excludedCrys_;
00119 
00120  public:
00121   
00122   //The default constructor
00123   HybridClusterAlgo(){ }
00124   
00125   //The real constructor
00126   HybridClusterAlgo(double eb_str, 
00127                     int step,
00128                     double ethres,
00129                     double eseed,
00130                     double xi,
00131                     bool useEtForXi,
00132                     double ewing,
00133                     std::vector<int> v_chstatus,
00134                     const PositionCalc& posCalculator,
00135                     bool dynamicEThres = false,
00136                     double eThresA = 0,
00137                     double eThresB = 0.1,
00138                     std::vector<int> severityToExclude=std::vector<int>(999),
00139                     //double severityRecHitThreshold=0.08,
00140                     //int severitySpikeId=1,
00141                     //double severitySpikeThreshold=0,
00142                     bool excludeFromCluster=false
00143                     );
00144 //                    const edm::ParameterSet &bremRecoveryPset,
00145 
00146   // destructor
00147   ~HybridClusterAlgo() 
00148   {
00149      if (dynamicPhiRoad_) delete phiRoadAlgo_;
00150      delete topo_; 
00151     //     delete SCShape_;
00152   } 
00153 
00154   void setDynamicPhiRoad(const edm::ParameterSet &bremRecoveryPset)
00155   {
00156      dynamicPhiRoad_ = true;
00157      phiRoadAlgo_ = new BremRecoveryPhiRoadAlgo(bremRecoveryPset);
00158   }
00159 
00160   //Hand over the map, the geometry, and I'll hand you back clusters.
00161   void makeClusters(const EcalRecHitCollection*,
00162                     const CaloSubdetectorGeometry * geometry,
00163                     reco::BasicClusterCollection &basicClusters,
00164                     const EcalSeverityLevelAlgo * sevLv,
00165                     bool regional = false,
00166                     const std::vector<EcalEtaPhiRegion>& regions = std::vector<EcalEtaPhiRegion>()
00167                     );
00168 
00169   //Make superclusters from the references to the BasicClusters in the event.
00170   reco::SuperClusterCollection makeSuperClusters(const reco::CaloClusterPtrVector&);
00171 
00172   //The routine doing the real work.
00173   void mainSearch(const EcalRecHitCollection* hits, const CaloSubdetectorGeometry * geometry);
00174   
00175   //Make dominos for the hybrid method.
00176   double makeDomino(EcalBarrelNavigator &navigator, std::vector <EcalRecHit> &cells);
00177 
00178 };
00179 
00180 #endif