Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:31:00

0001 #ifndef SimTracker_TrackAssociation_trackHitsToClusterRefs_h
0002 #define SimTracker_TrackAssociation_trackHitsToClusterRefs_h
0003 
0004 #include "DataFormats/SiPixelDetId/interface/PixelSubdetector.h"
0005 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
0006 #include "DataFormats/TrackerRecHit2D/interface/Phase2TrackerRecHit1D.h"
0007 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0008 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0009 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
0010 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0011 #include "DataFormats/TrackerRecHit2D/interface/VectorHit.h"
0012 #include "DataFormats/TrackingRecHit/interface/TrackingRecHitFwd.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 
0015 namespace track_associator {
0016   inline const TrackingRecHit *getHitFromIter(trackingRecHit_iterator iter) { return &(**iter); }
0017 
0018   inline const TrackingRecHit *getHitFromIter(TrackingRecHitCollection::const_iterator iter) { return &(*iter); }
0019 
0020   template <typename iter>
0021   std::vector<OmniClusterRef> hitsToClusterRefs(iter begin, iter end) {
0022     std::vector<OmniClusterRef> returnValue;
0023     for (iter iRecHit = begin; iRecHit != end; ++iRecHit) {
0024       const TrackingRecHit *rhit = getHitFromIter(iRecHit);
0025       if (trackerHitRTTI::isFromDet(*rhit)) {
0026         int subdetid = rhit->geographicalId().subdetId();
0027         if (subdetid == PixelSubdetector::PixelBarrel || subdetid == PixelSubdetector::PixelEndcap) {
0028           const SiPixelRecHit *pRHit = dynamic_cast<const SiPixelRecHit *>(rhit);
0029           if (pRHit && !pRHit->cluster().isNonnull())
0030             edm::LogError("TrackAssociator") << ">>> RecHit does not have an associated cluster!"
0031                                              << " file: " << __FILE__ << " line: " << __LINE__;
0032           returnValue.push_back(pRHit->omniClusterRef());
0033         } else if (subdetid == SiStripDetId::TIB || subdetid == SiStripDetId::TOB || subdetid == SiStripDetId::TID ||
0034                    subdetid == SiStripDetId::TEC) {
0035           const std::type_info &tid = typeid(*rhit);
0036           if (tid == typeid(SiStripMatchedRecHit2D)) {
0037             const SiStripMatchedRecHit2D *sMatchedRHit = dynamic_cast<const SiStripMatchedRecHit2D *>(rhit);
0038             if (!sMatchedRHit->monoHit().cluster().isNonnull() || !sMatchedRHit->stereoHit().cluster().isNonnull())
0039               edm::LogError("TrackAssociator") << ">>> RecHit does not have an associated cluster!"
0040                                                << " file: " << __FILE__ << " line: " << __LINE__;
0041             returnValue.push_back(sMatchedRHit->monoClusterRef());
0042             returnValue.push_back(sMatchedRHit->stereoClusterRef());
0043           } else if (tid == typeid(SiStripRecHit2D)) {
0044             const SiStripRecHit2D *sRHit = dynamic_cast<const SiStripRecHit2D *>(rhit);
0045             if (!sRHit->cluster().isNonnull())
0046               edm::LogError("TrackAssociator") << ">>> RecHit does not have an associated cluster!"
0047                                                << " file: " << __FILE__ << " line: " << __LINE__;
0048             returnValue.push_back(sRHit->omniClusterRef());
0049           } else if (tid == typeid(SiStripRecHit1D)) {
0050             const SiStripRecHit1D *sRHit = dynamic_cast<const SiStripRecHit1D *>(rhit);
0051             if (!sRHit->cluster().isNonnull())
0052               edm::LogError("TrackAssociator") << ">>> RecHit does not have an associated cluster!"
0053                                                << " file: " << __FILE__ << " line: " << __LINE__;
0054             returnValue.push_back(sRHit->omniClusterRef());
0055           } else if (tid == typeid(Phase2TrackerRecHit1D)) {
0056             const Phase2TrackerRecHit1D *ph2Hit = dynamic_cast<const Phase2TrackerRecHit1D *>(rhit);
0057             if (!ph2Hit->cluster().isNonnull())
0058               edm::LogError("TrackAssociator") << ">>> RecHit does not have an associated cluster!"
0059                                                << " file: " << __FILE__ << " line: " << __LINE__;
0060             returnValue.push_back(ph2Hit->omniClusterRef());
0061           } else if (tid == typeid(VectorHit)) {
0062             const VectorHit *vectorHit = dynamic_cast<const VectorHit *>(rhit);
0063             if (!vectorHit->cluster().isNonnull())
0064               edm::LogError("TrackAssociator") << ">>> RecHit does not have an associated cluster!"
0065                                                << " file: " << __FILE__ << " line: " << __LINE__;
0066             returnValue.push_back(vectorHit->firstClusterRef());
0067 
0068           } else {
0069             auto const &thit = static_cast<BaseTrackerRecHit const &>(*rhit);
0070             if (thit.isProjected()) {
0071             } else {
0072               edm::LogError("TrackAssociator") << ">>> getMatchedClusters: TrackingRecHit not associated to "
0073                                                   "any SiStripCluster! subdetid = "
0074                                                << subdetid;
0075             }
0076           }
0077         } else {
0078           edm::LogError("TrackAssociator") << ">>> getMatchedClusters: TrackingRecHit not associated to any "
0079                                               "cluster! subdetid = "
0080                                            << subdetid;
0081         }
0082       }
0083     }
0084     return returnValue;
0085   }
0086 }  // namespace track_associator
0087 
0088 #endif