Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:29

0001 // VI January 2012: needs to be migrated to use cluster directly
0002 
0003 #include "RecoTracker/PixelLowPtUtilities/interface/ClusterShapeTrajectoryFilter.h"
0004 
0005 #include "RecoTracker/PixelLowPtUtilities/interface/ClusterShapeHitFilter.h"
0006 
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "FWCore/Framework/interface/EventSetup.h"
0009 #include "FWCore/Framework/interface/ConsumesCollector.h"
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 
0012 #include "TrackingTools/PatternTools/interface/Trajectory.h"
0013 #include "TrackingTools/PatternTools/interface/TempTrajectory.h"
0014 
0015 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0016 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0017 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2D.h"
0018 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0019 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
0020 #include "DataFormats/SiPixelCluster/interface/SiPixelClusterShapeCache.h"
0021 
0022 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
0023 #include "Geometry/CommonDetUnit/interface/GeomDetType.h"
0024 
0025 #include "MagneticField/Engine/interface/MagneticField.h"
0026 
0027 #include "CondFormats/DataRecord/interface/SiPixelLorentzAngleRcd.h"
0028 #include "CondFormats/DataRecord/interface/SiStripLorentzAngleRcd.h"
0029 
0030 #include "RecoTracker/Record/interface/CkfComponentsRecord.h"
0031 
0032 #include <vector>
0033 
0034 ClusterShapeTrajectoryFilter::ClusterShapeTrajectoryFilter(const edm::ParameterSet& iConfig, edm::ConsumesCollector& iC)
0035     : theCacheToken(iC.consumes<SiPixelClusterShapeCache>(iConfig.getParameter<edm::InputTag>("cacheSrc"))),
0036       theFilterToken(iC.esConsumes(edm::ESInputTag("", "ClusterShapeHitFilter"))),
0037       theFilter(nullptr) {}
0038 
0039 ClusterShapeTrajectoryFilter::~ClusterShapeTrajectoryFilter() {}
0040 
0041 void ClusterShapeTrajectoryFilter::fillPSetDescription(edm::ParameterSetDescription& iDesc) {
0042   iDesc.add<edm::InputTag>("cacheSrc", edm::InputTag("siPixelClusterShapeCache"));
0043 }
0044 
0045 void ClusterShapeTrajectoryFilter::setEvent(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0046   theFilter = &iSetup.getData(theFilterToken);
0047   theCache = &iEvent.get(theCacheToken);
0048 }
0049 
0050 bool ClusterShapeTrajectoryFilter::toBeContinued(Trajectory& trajectory) const {
0051   assert(theCache);
0052   std::vector<TrajectoryMeasurement> tms = trajectory.measurements();
0053 
0054   for (std::vector<TrajectoryMeasurement>::const_iterator tm = tms.begin(); tm != tms.end(); tm++) {
0055     const TrackingRecHit* ttRecHit = &(*((*tm).recHit()));
0056 
0057     if (ttRecHit->isValid()) {
0058       const TrackingRecHit* tRecHit = ttRecHit->hit();
0059 
0060       TrajectoryStateOnSurface ts = (*tm).updatedState();
0061       const GlobalVector gdir = ts.globalDirection();
0062 
0063       if (ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::PixelBarrel ||
0064           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::PixelEndcap ||
0065           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::P1PXB ||
0066           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::P1PXEC ||
0067           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::P2PXB ||
0068           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::P2PXEC) {  // pixel
0069         const SiPixelRecHit* recHit = dynamic_cast<const SiPixelRecHit*>(tRecHit);
0070 
0071         if (recHit != nullptr)
0072           return theFilter->isCompatible(*recHit, gdir, *theCache);
0073       } else if (GeomDetEnumerators::isTrackerStrip(ttRecHit->det()->subDetector())) {  // strip
0074         if (dynamic_cast<const SiStripMatchedRecHit2D*>(tRecHit) != nullptr) {          // glued
0075           const SiStripMatchedRecHit2D* recHit = dynamic_cast<const SiStripMatchedRecHit2D*>(tRecHit);
0076 
0077           if (recHit != nullptr) {
0078             return (theFilter->isCompatible(recHit->monoHit(), gdir) &&
0079                     theFilter->isCompatible(recHit->stereoHit(), gdir));
0080           }
0081         } else {                                                           // single
0082           if (dynamic_cast<const SiStripRecHit2D*>(tRecHit) != nullptr) {  // normal
0083             const SiStripRecHit2D* recHit = dynamic_cast<const SiStripRecHit2D*>(tRecHit);
0084 
0085             if (recHit != nullptr)
0086               return theFilter->isCompatible(*recHit, gdir);
0087           } else {  // projected
0088             const ProjectedSiStripRecHit2D* recHit = dynamic_cast<const ProjectedSiStripRecHit2D*>(tRecHit);
0089 
0090             if (recHit != nullptr)
0091               return theFilter->isCompatible(recHit->originalHit(), gdir);
0092           }
0093         }
0094       }
0095     }
0096   }
0097 
0098   return true;
0099 }
0100 
0101 bool ClusterShapeTrajectoryFilter::toBeContinued(TempTrajectory& trajectory) const {
0102   assert(theCache);
0103   const TempTrajectory::DataContainer& tms = trajectory.measurements();
0104 
0105   for (TempTrajectory::DataContainer::const_iterator tm = tms.rbegin(); tm != tms.rend(); --tm) {
0106     const TrackingRecHit* ttRecHit = &(*((*tm).recHit()));
0107 
0108     if (ttRecHit->isValid()) {
0109       const TrackingRecHit* tRecHit = ttRecHit->hit();
0110 
0111       TrajectoryStateOnSurface ts = (*tm).updatedState();
0112       GlobalVector gdir = ts.globalDirection();
0113 
0114       if (ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::PixelBarrel ||
0115           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::PixelEndcap ||
0116           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::P1PXB ||
0117           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::P1PXEC ||
0118           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::P2PXB ||
0119           ttRecHit->det()->subDetector() == GeomDetEnumerators::SubDetector::P2PXEC) {  // pixel
0120         const SiPixelRecHit* recHit = dynamic_cast<const SiPixelRecHit*>(tRecHit);
0121 
0122         if (recHit != nullptr)
0123           if (!theFilter->isCompatible(*recHit, gdir, *theCache)) {
0124             LogTrace("TrajectFilter") << "  [TrajectFilter] fail pixel";
0125             return false;
0126           }
0127       } else if (GeomDetEnumerators::isTrackerStrip(ttRecHit->det()->subDetector())) {  // strip
0128         if (dynamic_cast<const SiStripMatchedRecHit2D*>(tRecHit) != nullptr) {          // glued
0129           const SiStripMatchedRecHit2D* recHit = dynamic_cast<const SiStripMatchedRecHit2D*>(tRecHit);
0130 
0131           if (recHit != nullptr) {
0132             if (!theFilter->isCompatible(recHit->monoHit(), gdir)) {
0133               LogTrace("TrajectFilter") << "  [TrajectFilter] fail strip matched 1st";
0134               return false;
0135             }
0136 
0137             if (!theFilter->isCompatible(recHit->stereoHit(), gdir)) {
0138               LogTrace("TrajectFilter") << "  [TrajectFilter] fail strip matched 2nd";
0139               return false;
0140             }
0141           }
0142         } else {                                                           // single
0143           if (dynamic_cast<const SiStripRecHit2D*>(tRecHit) != nullptr) {  // normal
0144             const SiStripRecHit2D* recHit = dynamic_cast<const SiStripRecHit2D*>(tRecHit);
0145 
0146             if (recHit != nullptr)
0147               if (!theFilter->isCompatible(*recHit, gdir)) {
0148                 LogTrace("TrajectFilter") << "  [TrajectFilter] fail strip single";
0149                 return false;
0150               }
0151           } else {  // projected
0152             const ProjectedSiStripRecHit2D* recHit = dynamic_cast<const ProjectedSiStripRecHit2D*>(tRecHit);
0153 
0154             if (recHit != nullptr)
0155               if (!theFilter->isCompatible(recHit->originalHit(), gdir)) {
0156                 LogTrace("TrajectFilter") << "  [TrajectFilter] fail strip projected";
0157                 return false;
0158               }
0159           }
0160         }
0161       }
0162     }
0163   }
0164 
0165   return true;
0166 }
0167 
0168 bool ClusterShapeTrajectoryFilter::qualityFilter(const Trajectory& trajectory) const { return true; }
0169 
0170 bool ClusterShapeTrajectoryFilter::qualityFilter(const TempTrajectory& trajectory) const {
0171   TempTrajectory t = trajectory;
0172 
0173   // Check if ok
0174   if (toBeContinued(t))
0175     return true;
0176 
0177   // Should take out last
0178   if (t.measurements().size() <= 3)
0179     return false;
0180 
0181   return true;
0182 }