Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:25:11

0001 // Author: Felice Pantaleo,Marco Rovere - felice.pantaleo@cern.ch, marco.rovere@cern.ch
0002 // Date: 11/2018
0003 
0004 #ifndef __RecoHGCal_TICL_HGCDoublet_H__
0005 #define __RecoHGCal_TICL_HGCDoublet_H__
0006 
0007 #include <cmath>
0008 #include <vector>
0009 
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
0012 #include "DataFormats/HGCalReco/interface/TICLSeedingRegion.h"
0013 
0014 class HGCDoublet {
0015 public:
0016   using HGCntuplet = std::vector<unsigned int>;
0017 
0018   HGCDoublet(const int innerClusterId,
0019              const int outerClusterId,
0020              const int doubletId,
0021              const std::vector<reco::CaloCluster> *layerClusters,
0022              const int seedIndex,
0023              bool areSiblingClusters = false)
0024       : layerClusters_(layerClusters),
0025         theDoubletId_(doubletId),
0026         innerClusterId_(innerClusterId),
0027         outerClusterId_(outerClusterId),
0028         innerR_((*layerClusters)[innerClusterId].position().r()),
0029         outerR_((*layerClusters)[outerClusterId].position().r()),
0030         innerX_((*layerClusters)[innerClusterId].x()),
0031         outerX_((*layerClusters)[outerClusterId].x()),
0032         innerY_((*layerClusters)[innerClusterId].y()),
0033         outerY_((*layerClusters)[outerClusterId].y()),
0034         innerZ_((*layerClusters)[innerClusterId].z()),
0035         outerZ_((*layerClusters)[outerClusterId].z()),
0036         seedIndex_(seedIndex),
0037         alreadyVisited_(false),
0038         areSiblingClusters_(areSiblingClusters) {}
0039 
0040   double innerX() const { return innerX_; }
0041 
0042   double outerX() const { return outerX_; }
0043 
0044   double innerY() const { return innerY_; }
0045 
0046   double outerY() const { return outerY_; }
0047 
0048   double innerZ() const { return innerZ_; }
0049 
0050   double outerZ() const { return outerZ_; }
0051 
0052   double innerR() const { return innerR_; }
0053 
0054   double outerR() const { return outerR_; }
0055 
0056   int seedIndex() const { return seedIndex_; }
0057 
0058   int innerClusterId() const { return innerClusterId_; }
0059 
0060   int outerClusterId() const { return outerClusterId_; }
0061 
0062   bool areSiblingClusters() const { return areSiblingClusters_; }
0063 
0064   void tagAsOuterNeighbor(unsigned int otherDoublet) { outerNeighbors_.push_back(otherDoublet); }
0065 
0066   void tagAsInnerNeighbor(unsigned int otherDoublet) { innerNeighbors_.push_back(otherDoublet); }
0067 
0068   bool checkCompatibilityAndTag(std::vector<HGCDoublet> &allDoublets,
0069                                 const std::vector<int> &innerDoublets,
0070                                 const GlobalVector &refDir,
0071                                 float minCosTheta,
0072                                 float minCosPointing = 1.,
0073                                 bool debug = false);
0074 
0075   int areAligned(double xi,
0076                  double yi,
0077                  double zi,
0078                  double xo,
0079                  double yo,
0080                  double zo,
0081                  float minCosTheta,
0082                  float minCosPointing,
0083                  const GlobalVector &refDir,
0084                  bool debug = false) const;
0085 
0086   void findNtuplets(std::vector<HGCDoublet> &allDoublets,
0087                     HGCntuplet &tmpNtuplet,
0088                     int seedIndex,
0089                     const bool outInDFS,
0090                     const unsigned int outInHops,
0091                     const unsigned int maxOutInHops,
0092                     std::vector<std::pair<unsigned int, unsigned int> > &outInToVisit);
0093 
0094   void setVisited(bool visited) { alreadyVisited_ = visited; }
0095 
0096 private:
0097   const std::vector<reco::CaloCluster> *layerClusters_;
0098   std::vector<int> outerNeighbors_;
0099   std::vector<int> innerNeighbors_;
0100 
0101   const int theDoubletId_;
0102   const int innerClusterId_;
0103   const int outerClusterId_;
0104 
0105   const double innerR_;
0106   const double outerR_;
0107   const double innerX_;
0108   const double outerX_;
0109   const double innerY_;
0110   const double outerY_;
0111   const double innerZ_;
0112   const double outerZ_;
0113   int seedIndex_;
0114   bool alreadyVisited_;
0115   bool areSiblingClusters_;
0116 };
0117 
0118 #endif /*HGCDoublet_H_ */