Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 // File: TestAssociator.cc
0002 // Author:  P. Azzi
0003 // Creation Date:  PA May 2006 Initial version.
0004 //                 Pixel RecHits added by V.Chiochia - 18/5/06
0005 // 25/9/17 (W.T.Ford) Add Phase 2 Outer Tracker, common template function
0006 //
0007 //--------------------------------------------
0008 #include <memory>
0009 #include <string>
0010 #include <iostream>
0011 
0012 #include "SimTracker/TrackerHitAssociation/test/TestAssociator.h"
0013 
0014 //--- for Geometry:
0015 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0016 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0017 #include "Geometry/TrackerNumberingBuilder/interface/GeometricDet.h"
0018 #include "Geometry/CommonDetUnit/interface/PixelGeomDetType.h"
0019 #include "Geometry/TrackerGeometryBuilder/interface/StripGeomDetType.h"
0020 #include "DataFormats/DetId/interface/DetId.h"
0021 #include "DataFormats/GeometryVector/interface/LocalPoint.h"
0022 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0024 
0025 //--- for RecHits
0026 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHitCollection.h"
0027 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
0028 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
0029 #include "DataFormats/TrackerRecHit2D/interface/Phase2TrackerRecHit1D.h"
0030 
0031 //--- for SimHits
0032 #include "SimDataFormats/TrackingHit/interface/PSimHit.h"
0033 #include "SimDataFormats/TrackingHit/interface/PSimHitContainer.h"
0034 
0035 using namespace std;
0036 using namespace edm;
0037 
0038 void TestAssociator::analyze(const edm::Event& e, const edm::EventSetup& es) {
0039   using namespace edm;
0040   int pixelcounter = 0;
0041   int stripcounter = 0;
0042 
0043   edm::LogVerbatim("TrackAssociator") << " === TestAssociator ";
0044 
0045   // Get inputs
0046   edm::Handle<SiStripMatchedRecHit2DCollection> rechitsmatched;
0047   edm::Handle<SiStripRecHit2DCollection> rechitsrphi;
0048   edm::Handle<SiStripRecHit2DCollection> rechitsstereo;
0049   edm::Handle<SiPixelRecHitCollection> pixelrechits;
0050   edm::Handle<Phase2TrackerRecHit1DCollectionNew> phase2rechits;
0051 
0052   // Construct the associator object
0053   TrackerHitAssociator associate(e, trackerHitAssociatorConfig_);
0054 
0055   // Process each RecHit collection in turn
0056   if (doPixel_) {
0057     e.getByToken(siPixelRecHitsToken, pixelrechits);
0058     if (pixelrechits.isValid())
0059       printRechitSimhit(pixelrechits, "Pixel        ", pixelcounter, associate);
0060   }
0061   if (doStrip_) {
0062     if (useOTph2_) {
0063       e.getByToken(siPhase2RecHitsToken, phase2rechits);
0064       if (phase2rechits.isValid())
0065         printRechitSimhit(phase2rechits, "Phase 2 OT   ", stripcounter, associate);
0066     } else {
0067       e.getByToken(rphiRecHitToken, rechitsrphi);
0068       if (rechitsrphi.isValid())
0069         printRechitSimhit(rechitsrphi, "Strip rphi   ", stripcounter, associate);
0070       e.getByToken(stereoRecHitToken, rechitsstereo);
0071       if (rechitsstereo.isValid())
0072         printRechitSimhit(rechitsstereo, "Strip stereo ", stripcounter, associate);
0073       e.getByToken(matchedRecHitToken, rechitsmatched);
0074       if (rechitsmatched.isValid())
0075         printRechitSimhit(rechitsmatched, "Strip matched", stripcounter, associate);
0076     }
0077   }
0078   if (!doPixel_ && !doStrip_)
0079     throw edm::Exception(errors::Configuration, "Strip and pixel association disabled");
0080 
0081   edm::LogVerbatim("TrackAssociator") << " === TestAssociator end\n ";
0082 }
0083 
0084 template <typename rechitType>
0085 void TestAssociator::printRechitSimhit(const edm::Handle<edmNew::DetSetVector<rechitType>> rechitCollection,
0086                                        const char* rechitName,
0087                                        int hitCounter,
0088                                        TrackerHitAssociator& associate) const {
0089   std::vector<PSimHit> matched;
0090   // Loop over sensors with detected rechits of type rechitType
0091   for (auto const& theDetSet : *rechitCollection) {
0092     DetId detid = theDetSet.detId();
0093     uint32_t myid = detid.rawId();
0094     // Loop over the RecHits in this sensor
0095     for (auto const& rechit : theDetSet) {
0096       hitCounter++;
0097       edm::LogVerbatim("TrackAssociator") << hitCounter << ") " << rechitName << " RecHit subDet, DetId "
0098                                           << detid.subdetId() << ", " << myid << " Pos = " << rechit.localPosition();
0099       bool isPixel = false;
0100       float mindist = 999999;
0101       float dist, distx, disty;
0102       PSimHit closest;
0103       // Find the vector of SimHits matching this RecHit
0104       matched.clear();
0105       matched = associate.associateHit(rechit);
0106       if (!matched.empty()) {
0107         // Print out the SimHit positions and residuals
0108         for (auto const& m : matched) {
0109           edm::LogVerbatim("TrackAssociator")
0110               << " simtrack ID = " << m.trackId() << "                            Simhit Pos = " << m.localPosition();
0111           // Seek the smallest residual
0112           if (const SiPixelRecHit* dummy = dynamic_cast<const SiPixelRecHit*>(&rechit)) {
0113             isPixel = true;
0114             dist = (rechit.localPosition() - m.localPosition()).mag();  // pixels measure 2 dimensions
0115           } else {
0116             isPixel = false;
0117             dist = fabs(rechit.localPosition().x() - m.localPosition().x());
0118           }
0119           if (dist < mindist) {
0120             mindist = dist;
0121             closest = m;
0122           }
0123         }
0124         std::ostringstream st1;
0125         st1 << " Closest Simhit = " << closest.localPosition();
0126         if (isPixel) {
0127           distx = fabs(rechit.localPosition().x() - closest.localPosition().x());
0128           disty = fabs(rechit.localPosition().y() - closest.localPosition().y());
0129           st1 << ", diff(x,y) = (" << distx << ", " << disty << ")";
0130         }
0131         edm::LogVerbatim("TrackAssociator") << st1.str() << ", |diff| = " << mindist;
0132       }
0133     }
0134   }  // end loop on detSets
0135 }
0136 
0137 //---------------
0138 // Constructor --
0139 //---------------
0140 
0141 TestAssociator::TestAssociator(edm::ParameterSet const& conf)
0142     : trackerHitAssociatorConfig_(conf, consumesCollector()),
0143       doPixel_(conf.getParameter<bool>("associatePixel")),
0144       doStrip_(conf.getParameter<bool>("associateStrip")),
0145       useOTph2_(conf.getParameter<bool>("usePhase2Tracker")) {
0146   matchedRecHitToken =
0147       consumes<edmNew::DetSetVector<SiStripMatchedRecHit2D>>(conf.getParameter<edm::InputTag>("matchedRecHit"));
0148   rphiRecHitToken = consumes<edmNew::DetSetVector<SiStripRecHit2D>>(conf.getParameter<edm::InputTag>("rphiRecHit"));
0149   stereoRecHitToken = consumes<edmNew::DetSetVector<SiStripRecHit2D>>(conf.getParameter<edm::InputTag>("stereoRecHit"));
0150   siPixelRecHitsToken =
0151       consumes<edmNew::DetSetVector<SiPixelRecHit>>(conf.getParameter<edm::InputTag>("siPixelRecHits"));
0152   siPhase2RecHitsToken =
0153       consumes<edmNew::DetSetVector<Phase2TrackerRecHit1D>>(conf.getParameter<edm::InputTag>("siPhase2RecHits"));
0154 }
0155 
0156 //---------------
0157 // Destructor --
0158 //---------------
0159 
0160 TestAssociator::~TestAssociator() {}