Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-01-17 23:57:20

0001 // -*- C++ -*-
0002 //
0003 // Package:    Alignment/OfflineValidation
0004 // Class:      GeneralPurposeTrackAnalyzer
0005 //
0006 /**\class GeneralPurposeTrackAnalyzer GeneralPurposeTrackAnalyzer.cc Alignment/OfflineValidation/plugins/GeneralPurposeTrackAnalyzer.cc
0007 
0008 */
0009 //
0010 // Original Author:  Marco Musich
0011 //         Created:  Mon, 13 Jun 2016 15:07:11 GMT
0012 //
0013 //
0014 
0015 // ROOT includes files
0016 
0017 #include "TMath.h"
0018 #include "TFile.h"
0019 #include "TH1D.h"
0020 #include "TH1I.h"
0021 #include "TH2D.h"
0022 #include "TProfile.h"
0023 
0024 // system includes files
0025 
0026 #include <cmath>
0027 #include <fstream>
0028 #include <iostream>
0029 #include <map>
0030 #include <string>
0031 #include <vector>
0032 #include <boost/range/adaptor/indexed.hpp>
0033 
0034 // user include files
0035 #include "CommonTools/TrackerMap/interface/TrackerMap.h"
0036 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0037 #include "CondCore/SiPixelPlugins/interface/SiPixelPayloadInspectorHelper.h"
0038 #include "CondFormats/AlignmentRecord/interface/GlobalPositionRcd.h"
0039 #include "CondFormats/DataRecord/interface/SiPixelFedCablingMapRcd.h"
0040 #include "CondFormats/DataRecord/interface/SiStripCondDataRecords.h"
0041 #include "CondFormats/SiPixelObjects/interface/SiPixelFedCablingMap.h"
0042 #include "CondFormats/SiStripObjects/interface/SiStripLatency.h"
0043 #include "DQM/SiPixelPhase1Common/interface/SiPixelCoordinates.h"
0044 #include "DQM/TrackerRemapper/interface/Phase1PixelMaps.h"
0045 #include "DQM/TrackerRemapper/interface/Phase1PixelROCMaps.h"
0046 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0047 #include "DataFormats/Common/interface/TriggerResults.h"
0048 #include "DataFormats/DetId/interface/DetId.h"
0049 #include "DataFormats/GeometrySurface/interface/LocalError.h"
0050 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0051 #include "DataFormats/SiStripDetId/interface/SiStripDetId.h"
0052 #include "DataFormats/TrackReco/interface/HitPattern.h"
0053 #include "DataFormats/TrackReco/interface/Track.h"
0054 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0055 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
0056 #include "DataFormats/TrackerRecHit2D/interface/ProjectedSiStripRecHit2D.h"
0057 #include "DataFormats/TrackerRecHit2D/interface/SiPixelRecHit.h"
0058 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0059 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2DCollection.h"
0060 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit1D.h"
0061 #include "DataFormats/TrackerRecHit2D/interface/SiStripRecHit2DCollection.h"
0062 #include "DataFormats/TrackerRecHit2D/interface/SiTrackerMultiRecHit.h"
0063 #include "DataFormats/TrackingRecHit/interface/TrackingRecHit.h"
0064 #include "DataFormats/VertexReco/interface/Vertex.h"
0065 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0066 #include "FWCore/Common/interface/TriggerNames.h"
0067 #include "FWCore/Framework/interface/Event.h"
0068 #include "FWCore/Framework/interface/EventSetup.h"
0069 #include "FWCore/Framework/interface/Frameworkfwd.h"
0070 #include "FWCore/Framework/interface/MakerMacros.h"
0071 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0072 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0073 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0074 #include "FWCore/ServiceRegistry/interface/Service.h"
0075 #include "Geometry/Records/interface/IdealGeometryRecord.h"
0076 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0077 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0078 #include "MagneticField/Engine/interface/MagneticField.h"
0079 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
0080 #include "TrackingTools/PatternTools/interface/Trajectory.h"
0081 #include "TrackingTools/TrackFitters/interface/TrajectoryStateCombiner.h"
0082 
0083 // toggle to enable debugging
0084 #define DEBUG 0
0085 
0086 const int kBPIX = PixelSubdetector::PixelBarrel;
0087 const int kFPIX = PixelSubdetector::PixelEndcap;
0088 
0089 class GeneralPurposeTrackAnalyzer : public edm::one::EDAnalyzer<edm::one::WatchRuns, edm::one::SharedResources> {
0090 public:
0091   GeneralPurposeTrackAnalyzer(const edm::ParameterSet &pset)
0092       : geomToken_(esConsumes()),
0093         magFieldToken_(esConsumes<edm::Transition::BeginRun>()),
0094         geomTokenBR_(esConsumes<edm::Transition::BeginRun>()),
0095         trackerTopologyTokenBR_(esConsumes<edm::Transition::BeginRun>()),
0096         siPixelFedCablingMapTokenBR_(esConsumes<edm::Transition::BeginRun>()) {
0097     doLatencyAnalysis_ = pset.getParameter<bool>("doLatencyAnalysis");
0098     if (doLatencyAnalysis_) {
0099       latencyToken_ = esConsumes<edm::Transition::BeginRun>();
0100     }
0101 
0102     usesResource(TFileService::kSharedResource);
0103 
0104     TkTag_ = pset.getParameter<edm::InputTag>("TkTag");
0105     theTrackCollectionToken_ = consumes<reco::TrackCollection>(TkTag_);
0106 
0107     TriggerResultsTag_ = pset.getParameter<edm::InputTag>("TriggerResultsTag");
0108     hltresultsToken_ = consumes<edm::TriggerResults>(TriggerResultsTag_);
0109 
0110     BeamSpotTag_ = pset.getParameter<edm::InputTag>("BeamSpotTag");
0111     beamspotToken_ = consumes<reco::BeamSpot>(BeamSpotTag_);
0112 
0113     VerticesTag_ = pset.getParameter<edm::InputTag>("VerticesTag");
0114     vertexToken_ = consumes<reco::VertexCollection>(VerticesTag_);
0115 
0116     isCosmics_ = pset.getParameter<bool>("isCosmics");
0117 
0118     pmap = std::make_unique<TrackerMap>("Pixel");
0119     pmap->onlyPixel(true);
0120     pmap->setTitle("Pixel Hit entries");
0121     pmap->setPalette(1);
0122 
0123     tmap = std::make_unique<TrackerMap>("Strip");
0124     tmap->setTitle("Strip Hit entries");
0125     tmap->setPalette(1);
0126 
0127     pixelmap = std::make_unique<Phase1PixelMaps>("COLZ0 L");
0128     pixelmap->bookBarrelHistograms("entriesBarrel", "# hits", "# pixel hits");
0129     pixelmap->bookForwardHistograms("entriesForward", "# hits", "# pixel hits");
0130 
0131     pixelrocsmap_ = std::make_unique<Phase1PixelROCMaps>("");
0132   }
0133 
0134   ~GeneralPurposeTrackAnalyzer() override = default;
0135 
0136   static void fillDescriptions(edm::ConfigurationDescriptions &);
0137 
0138   template <class OBJECT_TYPE>
0139   int index(const std::vector<OBJECT_TYPE *> &vec, const TString &name) {
0140     for (const auto &iter : vec | boost::adaptors::indexed(0)) {
0141       if (iter.value() && iter.value()->GetName() == name) {
0142         return iter.index();
0143       }
0144     }
0145     edm::LogError("GeneralPurposeTrackAnalyzer") << "@SUB=GeneralPurposeTrackAnalyzer::index"
0146                                                  << " could not find " << name;
0147     return -1;
0148   }
0149 
0150   template <typename T, typename... Args>
0151   T *book(const Args &...args) const {
0152     T *t = fs->make<T>(args...);
0153     return t;
0154   }
0155 
0156 private:
0157   // tokens for the event setup
0158   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0159   const edm::ESGetToken<MagneticField, IdealMagneticFieldRecord> magFieldToken_;
0160   edm::ESGetToken<SiStripLatency, SiStripLatencyRcd> latencyToken_;
0161 
0162   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomTokenBR_;
0163   const edm::ESGetToken<TrackerTopology, TrackerTopologyRcd> trackerTopologyTokenBR_;
0164   const edm::ESGetToken<SiPixelFedCablingMap, SiPixelFedCablingMapRcd> siPixelFedCablingMapTokenBR_;
0165 
0166   edm::ESHandle<MagneticField> magneticField_;
0167 
0168   SiPixelCoordinates coord_;
0169 
0170   edm::Service<TFileService> fs;
0171 
0172   std::unique_ptr<TrackerMap> tmap;
0173   std::unique_ptr<TrackerMap> pmap;
0174 
0175   std::unique_ptr<Phase1PixelMaps> pixelmap;
0176   std::unique_ptr<Phase1PixelROCMaps> pixelrocsmap_;
0177 
0178   TH1D *hchi2ndof;
0179   TH1D *hNtrk;
0180   TH1D *hNtrkZoom;
0181   TH1I *htrkQuality;
0182   TH1I *htrkAlgo;
0183   TH1I *htrkOriAlgo;
0184   TH1D *hNhighPurity;
0185   TH1D *hP;
0186   TH1D *hPPlus;
0187   TH1D *hPMinus;
0188   TH1D *hPt;
0189   TH1D *hMinPt;
0190   TH1D *hPtPlus;
0191   TH1D *hPtMinus;
0192   TH1D *hHit;
0193   TH1D *hHit2D;
0194 
0195   TH1D *hHitCountVsXBPix;
0196   TH1D *hHitCountVsXFPix;
0197   TH1D *hHitCountVsYBPix;
0198   TH1D *hHitCountVsYFPix;
0199   TH1D *hHitCountVsZBPix;
0200   TH1D *hHitCountVsZFPix;
0201 
0202   TH1D *hHitCountVsThetaBPix;
0203   TH1D *hHitCountVsPhiBPix;
0204 
0205   TH1D *hHitCountVsThetaFPix;
0206   TH1D *hHitCountVsPhiFPix;
0207 
0208   TH1D *hHitPlus;
0209   TH1D *hHitMinus;
0210 
0211   TH1D *hPhp;
0212   TH1D *hPthp;
0213   TH1D *hHithp;
0214   TH1D *hEtahp;
0215   TH1D *hPhihp;
0216   TH1D *hchi2ndofhp;
0217   TH1D *hchi2Probhp;
0218 
0219   TH1D *hCharge;
0220   TH1D *hQoverP;
0221   TH1D *hQoverPZoom;
0222   TH1D *hEta;
0223   TH1D *hEtaPlus;
0224   TH1D *hEtaMinus;
0225   TH1D *hPhi;
0226   TH1D *hPhiBarrel;
0227   TH1D *hPhiOverlapPlus;
0228   TH1D *hPhiOverlapMinus;
0229   TH1D *hPhiEndcapPlus;
0230   TH1D *hPhiEndcapMinus;
0231   TH1D *hPhiPlus;
0232   TH1D *hPhiMinus;
0233 
0234   TH1D *hDeltaPhi;
0235   TH1D *hDeltaEta;
0236   TH1D *hDeltaR;
0237 
0238   TH1D *hvx;
0239   TH1D *hvy;
0240   TH1D *hvz;
0241   TH1D *hd0;
0242   TH1D *hdz;
0243   TH1D *hdxy;
0244 
0245   TH2D *hd0PVvsphi;
0246   TH2D *hd0PVvseta;
0247   TH2D *hd0PVvspt;
0248 
0249   TH2D *hd0vsphi;
0250   TH2D *hd0vseta;
0251   TH2D *hd0vspt;
0252 
0253   TH1D *hnhpxb;
0254   TH1D *hnhpxe;
0255   TH1D *hnhTIB;
0256   TH1D *hnhTID;
0257   TH1D *hnhTOB;
0258   TH1D *hnhTEC;
0259 
0260   TH1D *hHitComposition;
0261 
0262   TProfile *pNBpixHitsVsVx;
0263   TProfile *pNBpixHitsVsVy;
0264   TProfile *pNBpixHitsVsVz;
0265 
0266   TH1D *hMultCand;
0267 
0268   TH1D *hdxyBS;
0269   TH1D *hd0BS;
0270   TH1D *hdzBS;
0271   TH1D *hdxyPV;
0272   TH1D *hd0PV;
0273   TH1D *hdzPV;
0274   TH1D *hrun;
0275   TH1D *hlumi;
0276 
0277   std::vector<TH1 *> vTrackHistos_;
0278   std::vector<TH1 *> vTrackProfiles_;
0279   std::vector<TH1 *> vTrack2DHistos_;
0280 
0281   TH1D *tksByTrigger_;
0282   TH1D *evtsByTrigger_;
0283 
0284   TH1D *modeByRun_;
0285   TH1D *fieldByRun_;
0286 
0287   int ievt;
0288   int itrks;
0289   int mode;
0290   float etaMax_;
0291   SiPixelPI::phase phase_;
0292 
0293   const TrackerGeometry *trackerGeometry_;
0294 
0295   edm::InputTag TkTag_;
0296   edm::InputTag TriggerResultsTag_;
0297   edm::InputTag BeamSpotTag_;
0298   edm::InputTag VerticesTag_;
0299 
0300   bool isCosmics_;
0301   bool doLatencyAnalysis_;
0302 
0303   edm::EDGetTokenT<reco::TrackCollection> theTrackCollectionToken_;
0304   edm::EDGetTokenT<edm::TriggerResults> hltresultsToken_;
0305   edm::EDGetTokenT<reco::BeamSpot> beamspotToken_;
0306   edm::EDGetTokenT<reco::VertexCollection> vertexToken_;
0307 
0308   std::map<std::string, std::pair<int, int> > triggerMap_;
0309   std::map<int, std::pair<int, float> > conditionsMap_;
0310   std::map<int, std::pair<int, int> > runInfoMap_;
0311 
0312   //*************************************************************
0313   void analyze(const edm::Event &event, const edm::EventSetup &setup) override
0314   //*************************************************************
0315   {
0316     ievt++;
0317 
0318     // geometry setup
0319     const TrackerGeometry *theGeometry = &setup.getData(geomToken_);
0320 
0321     edm::Handle<reco::TrackCollection> trackCollection = event.getHandle(theTrackCollectionToken_);
0322     const reco::TrackCollection tC = *(trackCollection.product());
0323     itrks += tC.size();
0324 
0325     runInfoMap_[event.run()].first += 1;
0326     runInfoMap_[event.run()].second += tC.size();
0327 
0328     if (DEBUG) {
0329       edm::LogInfo("GeneralPurposeTrackAnalyzer") << "Reconstructed " << tC.size() << " tracks" << std::endl;
0330     }
0331     //int iCounter=0;
0332 
0333     edm::Handle<edm::TriggerResults> hltresults = event.getHandle(hltresultsToken_);
0334     if (hltresults.isValid()) {
0335       const edm::TriggerNames &triggerNames_ = event.triggerNames(*hltresults);
0336       int ntrigs = hltresults->size();
0337       //const vector<std::string> &triggernames = triggerNames_.triggerNames();
0338 
0339       for (int itrig = 0; itrig != ntrigs; ++itrig) {
0340         const std::string &trigName = triggerNames_.triggerName(itrig);
0341         bool accept = hltresults->accept(itrig);
0342         if (accept == 1) {
0343           if (DEBUG) {
0344             edm::LogInfo("GeneralPurposeTrackAnalyzer")
0345                 << trigName << " " << accept << " ,track size: " << tC.size() << std::endl;
0346           }
0347           triggerMap_[trigName].first += 1;
0348           triggerMap_[trigName].second += tC.size();
0349           // triggerInfo.push_back(pair <std::string, int> (trigName, accept));
0350         }
0351       }
0352     }
0353 
0354     hrun->Fill(event.run());
0355     hlumi->Fill(event.luminosityBlock());
0356 
0357     int nHighPurityTracks = 0;
0358 
0359     for (auto track = tC.cbegin(); track != tC.cend(); track++) {
0360       unsigned int nHit2D = 0;
0361       std::bitset<16> rocsToMask;
0362       for (auto iHit = track->recHitsBegin(); iHit != track->recHitsEnd(); ++iHit) {
0363         if (this->isHit2D(**iHit)) {
0364           ++nHit2D;
0365         }
0366         // rest the ROCs for the map
0367         rocsToMask.reset();
0368         const DetId &detId = (*iHit)->geographicalId();
0369         const GeomDet *geomDet(theGeometry->idToDet(detId));
0370 
0371         const SiPixelRecHit *pixhit = dynamic_cast<const SiPixelRecHit *>(*iHit);
0372 
0373         if (pixhit) {
0374           if (pixhit->isValid()) {
0375             unsigned int subid = detId.subdetId();
0376             int detid_db = detId.rawId();
0377 
0378             // get the cluster
0379             auto clustp = pixhit->cluster();
0380 
0381             if (clustp.isNull())
0382               continue;
0383             auto const &cluster = *clustp;
0384             int row = cluster.x() - 0.5, col = cluster.y() - 0.5;
0385             int rocId = coord_.roc(detId, std::make_pair(row, col));
0386 
0387             if (phase_ == SiPixelPI::phase::zero) {
0388               pmap->fill(detid_db, 1);
0389             } else if (phase_ == SiPixelPI::phase::one) {
0390               rocsToMask.set(rocId);
0391               pixelrocsmap_->fillSelectedRocs(detid_db, rocsToMask, 1);
0392 
0393               if (subid == PixelSubdetector::PixelBarrel) {
0394                 pixelmap->fillBarrelBin("entriesBarrel", detid_db, 1);
0395               } else {
0396                 pixelmap->fillForwardBin("entriesForward", detid_db, 1);
0397               }
0398             }
0399 
0400             LocalPoint lp = (*iHit)->localPosition();
0401             //LocalError le = (*iHit)->localPositionError();
0402 
0403             GlobalPoint GP = geomDet->surface().toGlobal(lp);
0404 
0405             if ((subid == PixelSubdetector::PixelBarrel) || (subid == PixelSubdetector::PixelEndcap)) {
0406               // 1 = PXB, 2 = PXF
0407               if (subid == PixelSubdetector::PixelBarrel) {
0408                 hHitCountVsThetaBPix->Fill(GP.theta());
0409                 hHitCountVsPhiBPix->Fill(GP.phi());
0410 
0411                 hHitCountVsZBPix->Fill(GP.z());
0412                 hHitCountVsXBPix->Fill(GP.x());
0413                 hHitCountVsYBPix->Fill(GP.y());
0414 
0415               } else if (subid == PixelSubdetector::PixelEndcap) {
0416                 hHitCountVsThetaFPix->Fill(GP.theta());
0417                 hHitCountVsPhiFPix->Fill(GP.phi());
0418 
0419                 hHitCountVsZFPix->Fill(GP.z());
0420                 hHitCountVsXFPix->Fill(GP.x());
0421                 hHitCountVsYFPix->Fill(GP.y());
0422               }
0423             }
0424           }
0425         } else {
0426           if ((*iHit)->isValid() && phase_ != SiPixelPI::phase::two) {
0427             tmap->fill(detId.rawId(), 1);
0428           }
0429         }
0430       }
0431 
0432       hHit2D->Fill(nHit2D);
0433       hHit->Fill(track->numberOfValidHits());
0434       hnhpxb->Fill(track->hitPattern().numberOfValidPixelBarrelHits());
0435       hnhpxe->Fill(track->hitPattern().numberOfValidPixelEndcapHits());
0436       hnhTIB->Fill(track->hitPattern().numberOfValidStripTIBHits());
0437       hnhTID->Fill(track->hitPattern().numberOfValidStripTIDHits());
0438       hnhTOB->Fill(track->hitPattern().numberOfValidStripTOBHits());
0439       hnhTEC->Fill(track->hitPattern().numberOfValidStripTECHits());
0440 
0441       // fill hit composition histogram
0442       if (track->hitPattern().numberOfValidPixelBarrelHits() != 0) {
0443         hHitComposition->Fill(0., track->hitPattern().numberOfValidPixelBarrelHits());
0444 
0445         pNBpixHitsVsVx->Fill(track->vx(), track->hitPattern().numberOfValidPixelBarrelHits());
0446         pNBpixHitsVsVy->Fill(track->vy(), track->hitPattern().numberOfValidPixelBarrelHits());
0447         pNBpixHitsVsVz->Fill(track->vz(), track->hitPattern().numberOfValidPixelBarrelHits());
0448       }
0449       if (track->hitPattern().numberOfValidPixelEndcapHits() != 0) {
0450         hHitComposition->Fill(1., track->hitPattern().numberOfValidPixelEndcapHits());
0451       }
0452       if (track->hitPattern().numberOfValidStripTIBHits() != 0) {
0453         hHitComposition->Fill(2., track->hitPattern().numberOfValidStripTIBHits());
0454       }
0455       if (track->hitPattern().numberOfValidStripTIDHits() != 0) {
0456         hHitComposition->Fill(3., track->hitPattern().numberOfValidStripTIDHits());
0457       }
0458       if (track->hitPattern().numberOfValidStripTOBHits() != 0) {
0459         hHitComposition->Fill(4., track->hitPattern().numberOfValidStripTOBHits());
0460       }
0461       if (track->hitPattern().numberOfValidStripTECHits() != 0) {
0462         hHitComposition->Fill(5., track->hitPattern().numberOfValidStripTECHits());
0463       }
0464 
0465       hCharge->Fill(track->charge());
0466       hQoverP->Fill(track->qoverp());
0467       hQoverPZoom->Fill(track->qoverp());
0468       hPt->Fill(track->pt());
0469       hP->Fill(track->p());
0470       hchi2ndof->Fill(track->normalizedChi2());
0471       hEta->Fill(track->eta());
0472       hPhi->Fill(track->phi());
0473 
0474       if (fabs(track->eta()) < 0.8) {
0475         hPhiBarrel->Fill(track->phi());
0476       }
0477       if (track->eta() > 0.8 && track->eta() < 1.4) {
0478         hPhiOverlapPlus->Fill(track->phi());
0479       }
0480       if (track->eta() < -0.8 && track->eta() > -1.4) {
0481         hPhiOverlapMinus->Fill(track->phi());
0482       }
0483       if (track->eta() > 1.4) {
0484         hPhiEndcapPlus->Fill(track->phi());
0485       }
0486       if (track->eta() < -1.4) {
0487         hPhiEndcapMinus->Fill(track->phi());
0488       }
0489 
0490       hd0->Fill(track->d0());
0491       hdz->Fill(track->dz());
0492       hdxy->Fill(track->dxy());
0493       hvx->Fill(track->vx());
0494       hvy->Fill(track->vy());
0495       hvz->Fill(track->vz());
0496 
0497       htrkAlgo->Fill(static_cast<int>(track->algo()));
0498       htrkOriAlgo->Fill(static_cast<int>(track->originalAlgo()));
0499 
0500       int myquality = -99;
0501       if (track->quality(reco::TrackBase::undefQuality)) {
0502         myquality = -1;
0503         htrkQuality->Fill(myquality);
0504       }
0505       if (track->quality(reco::TrackBase::loose)) {
0506         myquality = 0;
0507         htrkQuality->Fill(myquality);
0508       }
0509       if (track->quality(reco::TrackBase::tight)) {
0510         myquality = 1;
0511         htrkQuality->Fill(myquality);
0512       }
0513       if (track->quality(reco::TrackBase::highPurity) && (!isCosmics_)) {
0514         myquality = 2;
0515         htrkQuality->Fill(myquality);
0516         hPhp->Fill(track->p());
0517         hPthp->Fill(track->pt());
0518         hHithp->Fill(track->numberOfValidHits());
0519         hEtahp->Fill(track->eta());
0520         hPhihp->Fill(track->phi());
0521         hchi2ndofhp->Fill(track->normalizedChi2());
0522         hchi2Probhp->Fill(TMath::Prob(track->chi2(), track->ndof()));
0523         nHighPurityTracks++;
0524       }
0525       if (track->quality(reco::TrackBase::confirmed)) {
0526         myquality = 3;
0527         htrkQuality->Fill(myquality);
0528       }
0529       if (track->quality(reco::TrackBase::goodIterative)) {
0530         myquality = 4;
0531         htrkQuality->Fill(myquality);
0532       }
0533 
0534       // Fill 1D track histos
0535       static const int etaindex = this->index(vTrackHistos_, "h_tracketa");
0536       vTrackHistos_[etaindex]->Fill(track->eta());
0537       static const int phiindex = this->index(vTrackHistos_, "h_trackphi");
0538       vTrackHistos_[phiindex]->Fill(track->phi());
0539       static const int numOfValidHitsindex = this->index(vTrackHistos_, "h_trackNumberOfValidHits");
0540       vTrackHistos_[numOfValidHitsindex]->Fill(track->numberOfValidHits());
0541       static const int numOfLostHitsindex = this->index(vTrackHistos_, "h_trackNumberOfLostHits");
0542       vTrackHistos_[numOfLostHitsindex]->Fill(track->numberOfLostHits());
0543 
0544       GlobalPoint gPoint(track->vx(), track->vy(), track->vz());
0545       double theLocalMagFieldInInverseGeV = magneticField_->inInverseGeV(gPoint).z();
0546       double kappa = -track->charge() * theLocalMagFieldInInverseGeV / track->pt();
0547 
0548       static const int kappaindex = this->index(vTrackHistos_, "h_curvature");
0549       vTrackHistos_[kappaindex]->Fill(kappa);
0550       static const int kappaposindex = this->index(vTrackHistos_, "h_curvature_pos");
0551       if (track->charge() > 0)
0552         vTrackHistos_[kappaposindex]->Fill(fabs(kappa));
0553       static const int kappanegindex = this->index(vTrackHistos_, "h_curvature_neg");
0554       if (track->charge() < 0)
0555         vTrackHistos_[kappanegindex]->Fill(fabs(kappa));
0556 
0557       double chi2Prob = TMath::Prob(track->chi2(), track->ndof());
0558       double normchi2 = track->normalizedChi2();
0559 
0560       static const int normchi2index = this->index(vTrackHistos_, "h_normchi2");
0561       vTrackHistos_[normchi2index]->Fill(normchi2);
0562       static const int chi2index = this->index(vTrackHistos_, "h_chi2");
0563       vTrackHistos_[chi2index]->Fill(track->chi2());
0564       static const int chi2Probindex = this->index(vTrackHistos_, "h_chi2Prob");
0565       vTrackHistos_[chi2Probindex]->Fill(chi2Prob);
0566       static const int ptindex = this->index(vTrackHistos_, "h_pt");
0567       static const int pt2index = this->index(vTrackHistos_, "h_ptrebin");
0568       vTrackHistos_[ptindex]->Fill(track->pt());
0569       vTrackHistos_[pt2index]->Fill(track->pt());
0570       if (track->ptError() != 0.) {
0571         static const int ptResolutionindex = this->index(vTrackHistos_, "h_ptResolution");
0572         vTrackHistos_[ptResolutionindex]->Fill(track->ptError() / track->pt());
0573       }
0574       // Fill track profiles
0575       static const int d0phiindex = this->index(vTrackProfiles_, "p_d0_vs_phi");
0576       vTrackProfiles_[d0phiindex]->Fill(track->phi(), track->d0());
0577       static const int dzphiindex = this->index(vTrackProfiles_, "p_dz_vs_phi");
0578       vTrackProfiles_[dzphiindex]->Fill(track->phi(), track->dz());
0579       static const int d0etaindex = this->index(vTrackProfiles_, "p_d0_vs_eta");
0580       vTrackProfiles_[d0etaindex]->Fill(track->eta(), track->d0());
0581       static const int dzetaindex = this->index(vTrackProfiles_, "p_dz_vs_eta");
0582       vTrackProfiles_[dzetaindex]->Fill(track->eta(), track->dz());
0583       static const int chiProbphiindex = this->index(vTrackProfiles_, "p_chi2Prob_vs_phi");
0584       vTrackProfiles_[chiProbphiindex]->Fill(track->phi(), chi2Prob);
0585       static const int chiProbabsd0index = this->index(vTrackProfiles_, "p_chi2Prob_vs_d0");
0586       vTrackProfiles_[chiProbabsd0index]->Fill(fabs(track->d0()), chi2Prob);
0587       static const int chiProbabsdzindex = this->index(vTrackProfiles_, "p_chi2Prob_vs_dz");
0588       vTrackProfiles_[chiProbabsdzindex]->Fill(track->dz(), chi2Prob);
0589       static const int chiphiindex = this->index(vTrackProfiles_, "p_chi2_vs_phi");
0590       vTrackProfiles_[chiphiindex]->Fill(track->phi(), track->chi2());
0591       static const int normchiphiindex = this->index(vTrackProfiles_, "p_normchi2_vs_phi");
0592       vTrackProfiles_[normchiphiindex]->Fill(track->phi(), normchi2);
0593       static const int chietaindex = this->index(vTrackProfiles_, "p_chi2_vs_eta");
0594       vTrackProfiles_[chietaindex]->Fill(track->eta(), track->chi2());
0595       static const int normchiptindex = this->index(vTrackProfiles_, "p_normchi2_vs_pt");
0596       vTrackProfiles_[normchiptindex]->Fill(track->pt(), normchi2);
0597       static const int normchipindex = this->index(vTrackProfiles_, "p_normchi2_vs_p");
0598       vTrackProfiles_[normchipindex]->Fill(track->p(), normchi2);
0599       static const int chiProbetaindex = this->index(vTrackProfiles_, "p_chi2Prob_vs_eta");
0600       vTrackProfiles_[chiProbetaindex]->Fill(track->eta(), chi2Prob);
0601       static const int normchietaindex = this->index(vTrackProfiles_, "p_normchi2_vs_eta");
0602       vTrackProfiles_[normchietaindex]->Fill(track->eta(), normchi2);
0603       static const int kappaphiindex = this->index(vTrackProfiles_, "p_kappa_vs_phi");
0604       vTrackProfiles_[kappaphiindex]->Fill(track->phi(), kappa);
0605       static const int kappaetaindex = this->index(vTrackProfiles_, "p_kappa_vs_eta");
0606       vTrackProfiles_[kappaetaindex]->Fill(track->eta(), kappa);
0607       static const int ptResphiindex = this->index(vTrackProfiles_, "p_ptResolution_vs_phi");
0608       vTrackProfiles_[ptResphiindex]->Fill(track->phi(), track->ptError() / track->pt());
0609       static const int ptResetaindex = this->index(vTrackProfiles_, "p_ptResolution_vs_eta");
0610       vTrackProfiles_[ptResetaindex]->Fill(track->eta(), track->ptError() / track->pt());
0611 
0612       // Fill 2D track histos
0613       static const int etaphiindex_2d = this->index(vTrack2DHistos_, "h2_phi_vs_eta");
0614       vTrack2DHistos_[etaphiindex_2d]->Fill(track->eta(), track->phi());
0615       static const int d0phiindex_2d = this->index(vTrack2DHistos_, "h2_d0_vs_phi");
0616       vTrack2DHistos_[d0phiindex_2d]->Fill(track->phi(), track->d0());
0617       static const int dzphiindex_2d = this->index(vTrack2DHistos_, "h2_dz_vs_phi");
0618       vTrack2DHistos_[dzphiindex_2d]->Fill(track->phi(), track->dz());
0619       static const int d0etaindex_2d = this->index(vTrack2DHistos_, "h2_d0_vs_eta");
0620       vTrack2DHistos_[d0etaindex_2d]->Fill(track->eta(), track->d0());
0621       static const int dzetaindex_2d = this->index(vTrack2DHistos_, "h2_dz_vs_eta");
0622       vTrack2DHistos_[dzetaindex_2d]->Fill(track->eta(), track->dz());
0623       static const int chiphiindex_2d = this->index(vTrack2DHistos_, "h2_chi2_vs_phi");
0624       vTrack2DHistos_[chiphiindex_2d]->Fill(track->phi(), track->chi2());
0625       static const int chiProbphiindex_2d = this->index(vTrack2DHistos_, "h2_chi2Prob_vs_phi");
0626       vTrack2DHistos_[chiProbphiindex_2d]->Fill(track->phi(), chi2Prob);
0627       static const int chiProbabsd0index_2d = this->index(vTrack2DHistos_, "h2_chi2Prob_vs_d0");
0628       vTrack2DHistos_[chiProbabsd0index_2d]->Fill(fabs(track->d0()), chi2Prob);
0629       static const int normchiphiindex_2d = this->index(vTrack2DHistos_, "h2_normchi2_vs_phi");
0630       vTrack2DHistos_[normchiphiindex_2d]->Fill(track->phi(), normchi2);
0631       static const int chietaindex_2d = this->index(vTrack2DHistos_, "h2_chi2_vs_eta");
0632       vTrack2DHistos_[chietaindex_2d]->Fill(track->eta(), track->chi2());
0633       static const int chiProbetaindex_2d = this->index(vTrack2DHistos_, "h2_chi2Prob_vs_eta");
0634       vTrack2DHistos_[chiProbetaindex_2d]->Fill(track->eta(), chi2Prob);
0635       static const int normchietaindex_2d = this->index(vTrack2DHistos_, "h2_normchi2_vs_eta");
0636       vTrack2DHistos_[normchietaindex_2d]->Fill(track->eta(), normchi2);
0637       static const int kappaphiindex_2d = this->index(vTrack2DHistos_, "h2_kappa_vs_phi");
0638       vTrack2DHistos_[kappaphiindex_2d]->Fill(track->phi(), kappa);
0639       static const int kappaetaindex_2d = this->index(vTrack2DHistos_, "h2_kappa_vs_eta");
0640       vTrack2DHistos_[kappaetaindex_2d]->Fill(track->eta(), kappa);
0641       static const int normchi2kappa_2d = this->index(vTrack2DHistos_, "h2_normchi2_vs_kappa");
0642       vTrack2DHistos_[normchi2kappa_2d]->Fill(normchi2, kappa);
0643 
0644       //dxy with respect to the beamspot
0645       reco::BeamSpot beamSpot;
0646       edm::Handle<reco::BeamSpot> beamSpotHandle = event.getHandle(beamspotToken_);
0647       if (beamSpotHandle.isValid()) {
0648         beamSpot = *beamSpotHandle;
0649         math::XYZPoint point(beamSpot.x0(), beamSpot.y0(), beamSpot.z0());
0650         double dxy = track->dxy(point);
0651         double dz = track->dz(point);
0652         hdxyBS->Fill(dxy);
0653         hd0BS->Fill(-dxy);
0654         hdzBS->Fill(dz);
0655       }
0656 
0657       //dxy with respect to the primary vertex
0658       reco::Vertex pvtx;
0659       edm::Handle<reco::VertexCollection> vertexHandle = event.getHandle(vertexToken_);
0660       double mindxy = 100.;
0661       double dz = 100;
0662       if (vertexHandle.isValid() && !isCosmics_) {
0663         for (auto pvtx = vertexHandle->cbegin(); pvtx != vertexHandle->cend(); ++pvtx) {
0664           math::XYZPoint mypoint(pvtx->x(), pvtx->y(), pvtx->z());
0665           if (abs(mindxy) > abs(track->dxy(mypoint))) {
0666             mindxy = track->dxy(mypoint);
0667             dz = track->dz(mypoint);
0668           }
0669         }
0670 
0671         hdxyPV->Fill(mindxy);
0672         hd0PV->Fill(-mindxy);
0673         hdzPV->Fill(dz);
0674 
0675         hd0PVvsphi->Fill(track->phi(), -mindxy);
0676         hd0PVvseta->Fill(track->eta(), -mindxy);
0677         hd0PVvspt->Fill(track->pt(), -mindxy);
0678 
0679       } else {
0680         hdxyPV->Fill(100);
0681         hd0PV->Fill(100);
0682         hdzPV->Fill(100);
0683 
0684         hd0vsphi->Fill(track->phi(), -track->dxy());
0685         hd0vseta->Fill(track->eta(), -track->dxy());
0686         hd0vspt->Fill(track->pt(), -track->dxy());
0687       }
0688 
0689       if (DEBUG) {
0690         edm::LogInfo("GeneralPurposeTrackAnalyzer") << "end of track loop" << std::endl;
0691       }
0692     }
0693 
0694     hNtrk->Fill(tC.size());
0695     hNtrkZoom->Fill(tC.size());
0696     hNhighPurity->Fill(nHighPurityTracks);
0697 
0698     if (DEBUG) {
0699       edm::LogInfo("GeneralPurposeTrackAnalyzer") << "end of analysis" << std::endl;
0700     }
0701   }
0702 
0703   //*************************************************************
0704   void beginRun(edm::Run const &run, edm::EventSetup const &setup) override
0705   //*************************************************************
0706   {
0707     // initialize runInfoMap_
0708     runInfoMap_[run.run()].first = 0;
0709     runInfoMap_[run.run()].second = 0;
0710 
0711     // Magnetic Field setup
0712     magneticField_ = setup.getHandle(magFieldToken_);
0713     float B_ = magneticField_.product()->inTesla(GlobalPoint(0, 0, 0)).mag();
0714 
0715     if (DEBUG) {
0716       edm::LogInfo("GeneralPurposeTrackAnalyzer")
0717           << "run number:" << run.run() << " magnetic field: " << B_ << " [T]" << std::endl;
0718     }
0719 
0720     const TrackerGeometry *trackerGeometry = &setup.getData(geomTokenBR_);
0721     if (trackerGeometry->isThere(GeomDetEnumerators::P2PXB) || trackerGeometry->isThere(GeomDetEnumerators::P2PXEC)) {
0722       phase_ = SiPixelPI::phase::two;
0723     } else if (trackerGeometry->isThere(GeomDetEnumerators::P1PXB) ||
0724                trackerGeometry->isThere(GeomDetEnumerators::P1PXEC)) {
0725       phase_ = SiPixelPI::phase::one;
0726     } else {
0727       phase_ = SiPixelPI::phase::zero;
0728     }
0729 
0730     // if it's a phase-2 geometry there are no phase-1 conditions
0731     if (phase_ == SiPixelPI::phase::two) {
0732       mode = 0;
0733     } else {
0734       if (doLatencyAnalysis_) {
0735         //SiStrip Latency
0736         const SiStripLatency *apvlat = &setup.getData(latencyToken_);
0737         if (apvlat->singleReadOutMode() == 1) {
0738           mode = 1;  // peak mode
0739         } else if (apvlat->singleReadOutMode() == 0) {
0740           mode = -1;  // deco mode
0741         }
0742       } else {
0743         mode = 0;
0744       }
0745     }
0746 
0747     conditionsMap_[run.run()].first = mode;
0748     conditionsMap_[run.run()].second = B_;
0749 
0750     // init the sipixel coordinates
0751     const TrackerTopology *trackerTopology = &setup.getData(trackerTopologyTokenBR_);
0752     const SiPixelFedCablingMap *siPixelFedCablingMap = &setup.getData(siPixelFedCablingMapTokenBR_);
0753 
0754     // Pixel Phase-1 helper class
0755     coord_.init(trackerTopology, trackerGeometry, siPixelFedCablingMap);
0756   }
0757 
0758   //*************************************************************
0759   void endRun(edm::Run const &, edm::EventSetup const &) override {}
0760   //*************************************************************
0761 
0762   //*************************************************************
0763   void beginJob() override
0764   //*************************************************************
0765   {
0766     if (DEBUG) {
0767       edm::LogInfo("GeneralPurposeTrackAnalyzer") << __LINE__ << std::endl;
0768     }
0769 
0770     TH1D::SetDefaultSumw2(kTRUE);
0771     etaMax_ = 3.0;
0772 
0773     ievt = 0;
0774     itrks = 0;
0775 
0776     hrun = book<TH1D>("h_run", "run", 100000, 230000, 240000);
0777     hlumi = book<TH1D>("h_lumi", "lumi", 1000, 0, 1000);
0778 
0779     hchi2ndof = book<TH1D>("h_chi2ndof", "chi2/ndf;#chi^{2}/ndf;tracks", 100, 0, 5.);
0780     hCharge = book<TH1D>("h_charge", "charge;Charge of the track;tracks", 5, -2.5, 2.5);
0781     hNtrk = book<TH1D>("h_Ntrk", "ntracks;Number of Tracks;events", 200, 0., 200.);
0782     hNtrkZoom = book<TH1D>("h_NtrkZoom", "Number of tracks; number of tracks;events", 10, 0., 10.);
0783     hNhighPurity =
0784         book<TH1D>("h_NhighPurity", "n. high purity tracks;Number of high purity tracks;events", 200, 0., 200.);
0785 
0786     htrkAlgo = book<TH1I>("h_trkAlgo",
0787                           "tracking step;iterative tracking step;tracks",
0788                           reco::TrackBase::algoSize,
0789                           0.,
0790                           double(reco::TrackBase::algoSize));
0791 
0792     htrkOriAlgo = book<TH1I>("h_trkOriAlgo",
0793                              "original tracking step;original iterative tracking step;tracks",
0794                              reco::TrackBase::algoSize,
0795                              0.,
0796                              double(reco::TrackBase::algoSize));
0797 
0798     for (size_t ibin = 0; ibin < reco::TrackBase::algoSize - 1; ibin++) {
0799       htrkAlgo->GetXaxis()->SetBinLabel(ibin + 1, (reco::TrackBase::algoNames[ibin]).c_str());
0800       htrkOriAlgo->GetXaxis()->SetBinLabel(ibin + 1, (reco::TrackBase::algoNames[ibin]).c_str());
0801     }
0802 
0803     htrkQuality = book<TH1I>("h_trkQuality", "track quality;track quality;tracks", 6, -1, 5);
0804     std::string qualities[7] = {"undef", "loose", "tight", "highPurity", "confirmed", "goodIterative"};
0805     for (int nbin = 1; nbin <= htrkQuality->GetNbinsX(); nbin++) {
0806       htrkQuality->GetXaxis()->SetBinLabel(nbin, (qualities[nbin - 1]).c_str());
0807     }
0808 
0809     hP = book<TH1D>("h_P", "Momentum;track momentum [GeV];tracks", 100, 0., 100.);
0810     hQoverP = book<TH1D>("h_qoverp", "Track q/p; track q/p [GeV^{-1}];tracks", 100, -1., 1.);
0811     hQoverPZoom = book<TH1D>("h_qoverpZoom", "Track q/p; track q/p [GeV^{-1}];tracks", 100, -0.1, 0.1);
0812     hPt = book<TH1D>("h_Pt", "Transverse Momentum;track p_{T} [GeV];tracks", 100, 0., 100.);
0813     hHit = book<TH1D>("h_nHits", "Number of hits;track n. hits;tracks", 50, -0.5, 49.5);
0814     hHit2D = book<TH1D>("h_nHit2D", "Number of 2D hits; number of 2D hits;tracks", 20, 0, 20);
0815 
0816     hHitCountVsZBPix = book<TH1D>("h_HitCountVsZBpix", "Number of BPix hits vs z;hit global z;hits", 60, -30, 30);
0817     hHitCountVsZFPix = book<TH1D>("h_HitCountVsZFpix", "Number of FPix hits vs z;hit global z;hits", 100, -100, 100);
0818 
0819     hHitCountVsXBPix = book<TH1D>("h_HitCountVsXBpix", "Number of BPix hits vs x;hit global x;hits", 20, -20, 20);
0820     hHitCountVsXFPix = book<TH1D>("h_HitCountVsXFpix", "Number of FPix hits vs x;hit global x;hits", 20, -20, 20);
0821 
0822     hHitCountVsYBPix = book<TH1D>("h_HitCountVsYBpix", "Number of BPix hits vs y;hit global y;hits", 20, -20, 20);
0823     hHitCountVsYFPix = book<TH1D>("h_HitCountVsYFpix", "Number of FPix hits vs y;hit global y;hits", 20, -20, 20);
0824 
0825     hHitCountVsThetaBPix =
0826         book<TH1D>("h_HitCountVsThetaBpix", "Number of BPix hits vs #theta;hit global #theta;hits", 20, 0., M_PI);
0827     hHitCountVsPhiBPix =
0828         book<TH1D>("h_HitCountVsPhiBpix", "Number of BPix hits vs #phi;hit global #phi;hits", 20, -M_PI, M_PI);
0829 
0830     hHitCountVsThetaFPix =
0831         book<TH1D>("h_HitCountVsThetaFpix", "Number of FPix hits vs #theta;hit global #theta;hits", 40, 0., M_PI);
0832     hHitCountVsPhiFPix =
0833         book<TH1D>("h_HitCountVsPhiFpix", "Number of FPix hits vs #phi;hit global #phi;hits", 20, -M_PI, M_PI);
0834 
0835     hEta = book<TH1D>("h_Eta", "Track pseudorapidity; track #eta;tracks", 100, -etaMax_, etaMax_);
0836     hPhi = book<TH1D>("h_Phi", "Track azimuth; track #phi;tracks", 100, -M_PI, M_PI);
0837 
0838     hPhiBarrel = book<TH1D>("h_PhiBarrel", "hPhiBarrel (0<|#eta|<0.8);track #Phi;tracks", 100, -M_PI, M_PI);
0839     hPhiOverlapPlus =
0840         book<TH1D>("h_PhiOverlapPlus", "hPhiOverlapPlus (0.8<#eta<1.4);track #phi;tracks", 100, -M_PI, M_PI);
0841     hPhiOverlapMinus =
0842         book<TH1D>("h_PhiOverlapMinus", "hPhiOverlapMinus (-1.4<#eta<-0.8);track #phi;tracks", 100, -M_PI, M_PI);
0843     hPhiEndcapPlus = book<TH1D>("h_PhiEndcapPlus", "hPhiEndcapPlus (#eta>1.4);track #phi;track", 100, -M_PI, M_PI);
0844     hPhiEndcapMinus = book<TH1D>("h_PhiEndcapMinus", "hPhiEndcapMinus (#eta<1.4);track #phi;tracks", 100, -M_PI, M_PI);
0845 
0846     if (!isCosmics_) {
0847       hPhp = book<TH1D>("h_P_hp", "Momentum (high purity);track momentum [GeV];tracks", 100, 0., 100.);
0848       hPthp = book<TH1D>("h_Pt_hp", "Transverse Momentum (high purity);track p_{T} [GeV];tracks", 100, 0., 100.);
0849       hHithp = book<TH1D>("h_nHit_hp", "Number of hits (high purity);track n. hits;tracks", 30, 0, 30);
0850       hEtahp = book<TH1D>("h_Eta_hp", "Track pseudorapidity (high purity); track #eta;tracks", 100, -etaMax_, etaMax_);
0851       hPhihp = book<TH1D>("h_Phi_hp", "Track azimuth (high purity); track #phi;tracks", 100, -M_PI, M_PI);
0852       hchi2ndofhp = book<TH1D>("h_chi2ndof_hp", "chi2/ndf (high purity);#chi^{2}/ndf;tracks", 100, 0, 5.);
0853       hchi2Probhp = book<TH1D>(
0854           "h_chi2_Prob_hp", "#chi^{2} probability (high purity);#chi^{2}prob_{Track};Number of Tracks", 100, 0.0, 1.);
0855 
0856       hvx = book<TH1D>("h_vx", "Track v_{x} ; track v_{x} [cm];tracks", 100, -1.5, 1.5);
0857       hvy = book<TH1D>("h_vy", "Track v_{y} ; track v_{y} [cm];tracks", 100, -1.5, 1.5);
0858       hvz = book<TH1D>("h_vz", "Track v_{z} ; track v_{z} [cm];tracks", 100, -20., 20.);
0859       hd0 = book<TH1D>("h_d0", "Track d_{0} ; track d_{0} [cm];tracks", 100, -1., 1.);
0860       hdxy = book<TH1D>("h_dxy", "Track d_{xy}; track d_{xy} [cm]; tracks", 100, -0.5, 0.5);
0861       hdz = book<TH1D>("h_dz", "Track d_{z} ; track d_{z} [cm]; tracks", 100, -20, 20);
0862 
0863       hd0PVvsphi =
0864           book<TH2D>("h2_d0PVvsphi", "hd0PVvsphi;track #phi;track d_{0}(PV) [cm]", 160, -M_PI, M_PI, 100, -1., 1.);
0865       hd0PVvseta =
0866           book<TH2D>("h2_d0PVvseta", "hdPV0vseta;track #eta;track d_{0}(PV) [cm]", 160, -2.5, 2.5, 100, -1., 1.);
0867       hd0PVvspt = book<TH2D>("h2_d0PVvspt", "hdPV0vspt;track p_{T};d_{0}(PV) [cm]", 50, 0., 100., 100, -1, 1.);
0868 
0869       hdxyBS = book<TH1D>("h_dxyBS", "hdxyBS; track d_{xy}(BS) [cm];tracks", 100, -0.1, 0.1);
0870       hd0BS = book<TH1D>("h_d0BS", "hd0BS ; track d_{0}(BS) [cm];tracks", 100, -0.1, 0.1);
0871       hdzBS = book<TH1D>("h_dzBS", "hdzBS ; track d_{z}(BS) [cm];tracks", 100, -12, 12);
0872       hdxyPV = book<TH1D>("h_dxyPV", "hdxyPV; track d_{xy}(PV) [cm];tracks", 100, -0.1, 0.1);
0873       hd0PV = book<TH1D>("h_d0PV", "hd0PV ; track d_{0}(PV) [cm];tracks", 100, -0.15, 0.15);
0874       hdzPV = book<TH1D>("h_dzPV", "hdzPV ; track d_{z}(PV) [cm];tracks", 100, -0.1, 0.1);
0875 
0876       hnhTIB = book<TH1D>("h_nHitTIB", "nhTIB;# hits in TIB; tracks", 20, 0., 20.);
0877       hnhTID = book<TH1D>("h_nHitTID", "nhTID;# hits in TID; tracks", 20, 0., 20.);
0878       hnhTOB = book<TH1D>("h_nHitTOB", "nhTOB;# hits in TOB; tracks", 20, 0., 20.);
0879       hnhTEC = book<TH1D>("h_nHitTEC", "nhTEC;# hits in TEC; tracks", 20, 0., 20.);
0880 
0881     } else {
0882       hvx = book<TH1D>("h_vx", "Track v_{x};track v_{x} [cm];tracks", 100, -100., 100.);
0883       hvy = book<TH1D>("h_vy", "Track v_{y};track v_{y} [cm];tracks", 100, -100., 100.);
0884       hvz = book<TH1D>("h_vz", "Track v_{z};track v_{z} [cm];track", 100, -100., 100.);
0885       hd0 = book<TH1D>("h_d0", "Track d_{0};track d_{0} [cm];track", 100, -100., 100.);
0886       hdxy = book<TH1D>("h_dxy", "Track d_{xy};track d_{xy} [cm];tracks", 100, -100, 100);
0887       hdz = book<TH1D>("h_dz", "Track d_{z};track d_{z} [cm];tracks", 100, -200, 200);
0888 
0889       hd0vsphi = book<TH2D>(
0890           "h2_d0vsphi", "Track d_{0} vs #phi; track #phi;track d_{0} [cm]", 160, -3.20, 3.20, 100, -100., 100.);
0891       hd0vseta = book<TH2D>(
0892           "h2_d0vseta", "Track d_{0} vs #eta; track #eta;track d_{0} [cm]", 160, -3.20, 3.20, 100, -100., 100.);
0893       hd0vspt =
0894           book<TH2D>("h2_d0vspt", "Track d_{0} vs p_{T};track p_{T};track d_{0} [cm]", 50, 0., 100., 100, -100, 100);
0895 
0896       hdxyBS = book<TH1D>("h_dxyBS", "Track d_{xy}(BS);d_{xy}(BS) [cm];tracks", 100, -100., 100.);
0897       hd0BS = book<TH1D>("h_d0BS", "Track d_{0}(BS);d_{0}(BS) [cm];tracks", 100, -100., 100.);
0898       hdzBS = book<TH1D>("h_dzBS", "Track d_{z}(BS);d_{z}(BS) [cm];tracks", 100, -100., 100.);
0899       hdxyPV = book<TH1D>("h_dxyPV", "Track d_{xy}(PV); d_{xy}(PV) [cm];tracks", 100, -100., 100.);
0900       hd0PV = book<TH1D>("h_d0PV", "Track d_{0}(PV); d_{0}(PV) [cm];tracks", 100, -100., 100.);
0901       hdzPV = book<TH1D>("h_dzPV", "Track d_{z}(PV); d_{z}(PV) [cm];tracks", 100, -100., 100.);
0902 
0903       hnhTIB = book<TH1D>("h_nHitTIB", "nhTIB;# hits in TIB; tracks", 30, 0., 30.);
0904       hnhTID = book<TH1D>("h_nHitTID", "nhTID;# hits in TID; tracks", 30, 0., 30.);
0905       hnhTOB = book<TH1D>("h_nHitTOB", "nhTOB;# hits in TOB; tracks", 30, 0., 30.);
0906       hnhTEC = book<TH1D>("h_nHitTEC", "nhTEC;# hits in TEC; tracks", 30, 0., 30.);
0907     }
0908 
0909     hnhpxb = book<TH1D>("h_nHitPXB", "nhpxb;# hits in Pixel Barrel; tracks", 10, 0., 10.);
0910     hnhpxe = book<TH1D>("h_nHitPXE", "nhpxe;# hits in Pixel Endcap; tracks", 10, 0., 10.);
0911 
0912     hHitComposition = book<TH1D>("h_hitcomposition", "track hit composition;;# hits", 6, -0.5, 5.5);
0913 
0914     pNBpixHitsVsVx =
0915         book<TProfile>("p_NpixHits_vs_Vx", "n. Barrel Pixel hits vs. v_{x};v_{x} (cm);n. BPix hits", 20, -20, 20);
0916 
0917     pNBpixHitsVsVy =
0918         book<TProfile>("p_NpixHits_vs_Vy", "n. Barrel Pixel hits vs. v_{y};v_{y} (cm);n. BPix hits", 20, -20, 20);
0919 
0920     pNBpixHitsVsVz =
0921         book<TProfile>("p_NpixHits_vs_Vz", "n. Barrel Pixel hits vs. v_{z};v_{z} (cm);n. BPix hits", 20, -100, 100);
0922 
0923     std::string dets[6] = {"PXB", "PXF", "TIB", "TID", "TOB", "TEC"};
0924 
0925     for (int i = 1; i <= hHitComposition->GetNbinsX(); i++) {
0926       hHitComposition->GetXaxis()->SetBinLabel(i, (dets[i - 1]).c_str());
0927     }
0928 
0929     // vector of track histograms
0930 
0931     // clang-format off
0932 
0933     vTrackHistos_.push_back(book<TH1F>("h_tracketa", "Track #eta;#eta_{Track};Number of Tracks", 90, -etaMax_, etaMax_));
0934     vTrackHistos_.push_back(book<TH1F>("h_trackphi", "Track #phi;#phi_{Track};Number of Tracks", 90, -M_PI, M_PI));
0935     vTrackHistos_.push_back(book<TH1F>("h_trackNumberOfValidHits", "Track # of valid hits;# of valid hits _{Track};Number of Tracks", 40, 0., 40.));
0936     vTrackHistos_.push_back(book<TH1F>("h_trackNumberOfLostHits", "Track # of lost hits;# of lost hits _{Track};Number of Tracks", 10, 0., 10.));
0937     vTrackHistos_.push_back(book<TH1F>("h_curvature", "Curvature #kappa;#kappa_{Track};Number of Tracks", 100, -.05, .05));
0938     vTrackHistos_.push_back(book<TH1F>("h_curvature_pos", "Curvature |#kappa| Positive Tracks;|#kappa_{pos Track}|;Number of Tracks", 100, .0, .05));
0939     vTrackHistos_.push_back(book<TH1F>("h_curvature_neg", "Curvature |#kappa| Negative Tracks;|#kappa_{neg Track}|;Number of Tracks", 100, .0, .05));
0940     vTrackHistos_.push_back(book<TH1F>("h_diff_curvature","Curvature |#kappa| Tracks Difference;|#kappa_{Track}|;# Pos Tracks - # Neg Tracks",100, .0, .05));
0941     vTrackHistos_.push_back(book<TH1F>("h_chi2", "Track #chi^{2};#chi^{2}_{Track};Number of Tracks", 500, -0.01, 500.));
0942     vTrackHistos_.push_back(book<TH1F>("h_chi2Prob", "#chi^{2} probability;Track Prob(#chi^{2},ndof);Number of Tracks", 100, 0.0, 1.));
0943     vTrackHistos_.push_back(book<TH1F>("h_normchi2", "#chi^{2}/ndof;#chi^{2}/ndof;Number of Tracks", 100, -0.01, 10.));
0944     //variable binning for chi2/ndof vs. pT
0945     double xBins[19] = {0., 0.15, 0.5, 1., 1.5, 2., 2.5, 3., 3.5, 4., 4.5, 5., 7., 10., 15., 25., 40., 100., 200.};
0946     vTrackHistos_.push_back(book<TH1F>("h_pt", "Track p_{T};p_{T}^{track} [GeV];Number of Tracks", 250, 0., 250));
0947     vTrackHistos_.push_back(book<TH1F>("h_ptrebin", "Track p_{T};p_{T}^{track} [GeV];Number of Tracks", 18, xBins));
0948     vTrackHistos_.push_back(book<TH1F>("h_ptResolution", "#delta_{p_{T}}/p_{T}^{track};#delta_{p_{T}}/p_{T}^{track};Number of Tracks", 100, 0., 0.5));
0949 
0950     vTrackProfiles_.push_back(book<TProfile>("p_d0_vs_phi", "Transverse Impact Parameter vs. #phi;#phi_{Track};#LT d_{0} #GT [cm]", 100, -M_PI, M_PI));
0951     vTrackProfiles_.push_back(book<TProfile>("p_dz_vs_phi", "Longitudinal Impact Parameter vs. #phi;#phi_{Track};#LT d_{z} #GT [cm]", 100, -M_PI, M_PI));
0952     vTrackProfiles_.push_back(book<TProfile>("p_d0_vs_eta", "Transverse Impact Parameter vs. #eta;#eta_{Track};#LT d_{0} #GT [cm]", 100, -etaMax_, etaMax_));
0953     vTrackProfiles_.push_back(book<TProfile>("p_dz_vs_eta", "Longitudinal Impact Parameter vs. #eta;#eta_{Track};#LT d_{z} #GT [cm]", 100, -etaMax_, etaMax_));
0954     vTrackProfiles_.push_back(book<TProfile>("p_chi2_vs_phi", "#chi^{2} vs. #phi;#phi_{Track};#LT #chi^{2} #GT", 100, -M_PI, M_PI));
0955     vTrackProfiles_.push_back(book<TProfile>("p_chi2Prob_vs_phi", "#chi^{2} probablility vs. #phi;#phi_{Track};#LT #chi^{2} probability#GT", 100, -M_PI, M_PI));
0956     vTrackProfiles_.push_back(book<TProfile>("p_chi2Prob_vs_d0", "#chi^{2} probablility vs. |d_{0}|;|d_{0}|[cm];#LT #chi^{2} probability#GT", 100, 0, 80));
0957     vTrackProfiles_.push_back(book<TProfile>("p_chi2Prob_vs_dz", "#chi^{2} probablility vs. dz;d_{z} [cm];#LT #chi^{2} probability#GT", 100, -30, 30));
0958     vTrackProfiles_.push_back(book<TProfile>("p_normchi2_vs_phi", "#chi^{2}/ndof vs. #phi;#phi_{Track};#LT #chi^{2}/ndof #GT", 100, -M_PI, M_PI));
0959     vTrackProfiles_.push_back(book<TProfile>("p_chi2_vs_eta", "#chi^{2} vs. #eta;#eta_{Track};#LT #chi^{2} #GT", 100, -etaMax_, etaMax_));
0960     vTrackProfiles_.push_back(book<TProfile>("p_normchi2_vs_pt", "norm #chi^{2} vs. p_{T}_{Track}; p_{T}_{Track};#LT #chi^{2}/ndof #GT", 18, xBins));
0961     vTrackProfiles_.push_back(book<TProfile>("p_normchi2_vs_p", "#chi^{2}/ndof vs. p_{Track};p_{Track};#LT #chi^{2}/ndof #GT", 18, xBins));
0962     vTrackProfiles_.push_back(book<TProfile>("p_chi2Prob_vs_eta","#chi^{2} probability vs. #eta;#eta_{Track};#LT #chi^{2} probability #GT",100,-etaMax_,etaMax_));
0963     vTrackProfiles_.push_back(book<TProfile>("p_normchi2_vs_eta", "#chi^{2}/ndof vs. #eta;#eta_{Track};#LT #chi^{2}/ndof #GT", 100, -etaMax_, etaMax_));
0964     vTrackProfiles_.push_back(book<TProfile>("p_kappa_vs_phi", "#kappa vs. #phi;#phi_{Track};#kappa", 100, -M_PI, M_PI));
0965     vTrackProfiles_.push_back(book<TProfile>("p_kappa_vs_eta", "#kappa vs. #eta;#eta_{Track};#kappa", 100, -etaMax_, etaMax_));
0966     vTrackProfiles_.push_back(book<TProfile>("p_ptResolution_vs_phi","#delta_{p_{T}}/p_{T}^{track};#phi^{track};#delta_{p_{T}}/p_{T}^{track}",100,-M_PI,M_PI));
0967     vTrackProfiles_.push_back(book<TProfile>("p_ptResolution_vs_eta","#delta_{p_{T}}/p_{T}^{track};#eta^{track};#delta_{p_{T}}/p_{T}^{track}",100,-etaMax_,etaMax_));
0968 
0969     vTrack2DHistos_.push_back(book<TH2F>("h2_phi_vs_eta", "Track #phi vs. #eta;#eta_{Track};#phi_{Track}",50, -etaMax_, etaMax_, 50, -M_PI, M_PI));
0970     vTrack2DHistos_.push_back(book<TH2F>("h2_d0_vs_phi", "Transverse Impact Parameter vs. #phi;#phi_{Track};d_{0} [cm]", 100, -M_PI, M_PI, 100, -1., 1.));
0971     vTrack2DHistos_.push_back(book<TH2F>("h2_dz_vs_phi", "Longitudinal Impact Parameter vs. #phi;#phi_{Track};d_{z} [cm]",100, -M_PI, M_PI, 100, -100., 100.));
0972     vTrack2DHistos_.push_back(book<TH2F>("h2_d0_vs_eta", "Transverse Impact Parameter vs. #eta;#eta_{Track};d_{0} [cm]", 100, -etaMax_, etaMax_, 100, -1., 1.));
0973     vTrack2DHistos_.push_back(book<TH2F>("h2_dz_vs_eta", "Longitudinal Impact Parameter vs. #eta;#eta_{Track};d_{z} [cm]", 100, -etaMax_, etaMax_, 100, -100., 100.));
0974     vTrack2DHistos_.push_back(book<TH2F>("h2_chi2_vs_phi", "#chi^{2} vs. #phi;#phi_{Track};#chi^{2}", 100, -M_PI, M_PI, 500, 0., 500.));
0975     vTrack2DHistos_.push_back(book<TH2F>("h2_chi2Prob_vs_phi", "#chi^{2} probability vs. #phi;#phi_{Track};#chi^{2} probability", 100, -M_PI, M_PI, 100, 0., 1.));
0976     vTrack2DHistos_.push_back(book<TH2F>("h2_chi2Prob_vs_d0", "#chi^{2} probability vs. |d_{0}|;|d_{0}| [cm];#chi^{2} probability", 100, 0, 80, 100, 0., 1.));
0977     vTrack2DHistos_.push_back(book<TH2F>("h2_normchi2_vs_phi", "#chi^{2}/ndof vs. #phi;#phi_{Track};#chi^{2}/ndof", 100, -M_PI, M_PI, 100, 0., 10.));
0978     vTrack2DHistos_.push_back(book<TH2F>("h2_chi2_vs_eta", "#chi^{2} vs. #eta;#eta_{Track};#chi^{2}", 100, -etaMax_, etaMax_, 500, 0., 500.));
0979     vTrack2DHistos_.push_back(book<TH2F>("h2_chi2Prob_vs_eta", "#chi^{2} probaility vs. #eta;#eta_{Track};#chi^{2} probability", 100, -etaMax_, etaMax_, 100, 0., 1.));
0980     vTrack2DHistos_.push_back(book<TH2F>("h2_normchi2_vs_eta", "#chi^{2}/ndof vs. #eta;#eta_{Track};#chi^{2}/ndof", 100, -etaMax_, etaMax_, 100, 0., 10.));
0981     vTrack2DHistos_.push_back(book<TH2F>("h2_kappa_vs_phi", "#kappa vs. #phi;#phi_{Track};#kappa", 100, -M_PI, M_PI, 100, .0, .05));
0982     vTrack2DHistos_.push_back(book<TH2F>("h2_kappa_vs_eta", "#kappa vs. #eta;#eta_{Track};#kappa", 100, -etaMax_, etaMax_, 100, .0, .05));
0983     vTrack2DHistos_.push_back(book<TH2F>("h2_normchi2_vs_kappa", "#kappa vs. #chi^{2}/ndof;#chi^{2}/ndof;#kappa", 100, 0., 10, 100, -.03, .03));
0984 
0985     // clang-format on
0986 
0987   }  //beginJob
0988 
0989   //*************************************************************
0990   void endJob() override
0991   //*************************************************************
0992   {
0993     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
0994     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "Events run in total: " << ievt << std::endl;
0995     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "n. tracks: " << itrks << std::endl;
0996     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
0997 
0998     int nFiringTriggers = triggerMap_.size();
0999     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "firing triggers: " << nFiringTriggers << std::endl;
1000     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
1001 
1002     tksByTrigger_ =
1003         book<TH1D>("tksByTrigger", "tracks by HLT path;;% of # traks", nFiringTriggers, -0.5, nFiringTriggers - 0.5);
1004     evtsByTrigger_ =
1005         book<TH1D>("evtsByTrigger", "events by HLT path;;% of # events", nFiringTriggers, -0.5, nFiringTriggers - 0.5);
1006 
1007     int i = 0;
1008     for (const auto &it : triggerMap_) {
1009       i++;
1010 
1011       double trkpercent = ((it.second).second) * 100. / double(itrks);
1012       double evtpercent = ((it.second).first) * 100. / double(ievt);
1013 
1014       std::cout.precision(4);
1015 
1016       edm::LogPrint("GeneralPurposeTrackAnalyzer")
1017           << "HLT path: " << std::setw(60) << std::left << it.first << " | events firing: " << std::right
1018           << std::setw(8) << (it.second).first << " (" << std::setw(8) << std::fixed << std::setprecision(4)
1019           << evtpercent << "%)"
1020           << " | tracks collected: " << std::setw(10) << (it.second).second << " (" << std::setw(8) << std::fixed
1021           << std::setprecision(4) << trkpercent << "%)";
1022 
1023       tksByTrigger_->SetBinContent(i, trkpercent);
1024       tksByTrigger_->GetXaxis()->SetBinLabel(i, (it.first).c_str());
1025 
1026       evtsByTrigger_->SetBinContent(i, evtpercent);
1027       evtsByTrigger_->GetXaxis()->SetBinLabel(i, (it.first).c_str());
1028     }
1029 
1030     int nRuns = conditionsMap_.size();
1031 
1032     std::vector<int> theRuns_;
1033     for (const auto &it : conditionsMap_) {
1034       theRuns_.push_back(it.first);
1035     }
1036 
1037     sort(theRuns_.begin(), theRuns_.end());
1038     int runRange = theRuns_.back() - theRuns_.front() + 1;
1039 
1040     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
1041     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "first run: " << theRuns_.front() << std::endl;
1042     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "last run:  " << theRuns_.back() << std::endl;
1043     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "considered runs: " << nRuns << std::endl;
1044     edm::LogPrint("GeneralPurposeTrackAnalyzer") << "*******************************" << std::endl;
1045 
1046     modeByRun_ = book<TH1D>("modeByRun",
1047                             "Strip APV mode by run number;;APV mode (-1=deco,+1=peak)",
1048                             runRange,
1049                             theRuns_.front() - 0.5,
1050                             theRuns_.back() + 0.5);
1051 
1052     fieldByRun_ = book<TH1D>("fieldByRun",
1053                              "CMS B-field intensity by run number;;B-field intensity [T]",
1054                              runRange,
1055                              theRuns_.front() - 0.5,
1056                              theRuns_.back() + 0.5);
1057 
1058     for (const auto &the_r : theRuns_) {
1059       if (conditionsMap_.find(the_r)->second.first != 0) {
1060         edm::LogPrint("GeneralPurposeTrackAnalyzer")
1061             << "run:" << the_r << " | isPeak: " << std::setw(4) << conditionsMap_.find(the_r)->second.first
1062             << "| B-field: " << conditionsMap_.find(the_r)->second.second << " [T]"
1063             << "| events: " << std::setw(10) << runInfoMap_.find(the_r)->second.first << ", tracks " << std::setw(10)
1064             << runInfoMap_.find(the_r)->second.second << std::endl;
1065       }
1066 
1067       modeByRun_->SetBinContent((the_r - theRuns_.front()) + 1, conditionsMap_.find(the_r)->second.first);
1068       fieldByRun_->SetBinContent((the_r - theRuns_.front()) + 1, conditionsMap_.find(the_r)->second.second);
1069       modeByRun_->GetXaxis()->SetBinLabel((the_r - theRuns_.front()) + 1, std::to_string(the_r).c_str());
1070       fieldByRun_->GetXaxis()->SetBinLabel((the_r - theRuns_.front()) + 1, std::to_string(the_r).c_str());
1071     }
1072 
1073     if (phase_ < SiPixelPI::phase::two) {
1074       if (phase_ == SiPixelPI::phase::zero) {
1075         pmap->save(true, 0, 0, "PixelHitMap.pdf", 600, 800);
1076         pmap->save(true, 0, 0, "PixelHitMap.png", 500, 750);
1077       }
1078 
1079       tmap->save(true, 0, 0, "StripHitMap.pdf");
1080       tmap->save(true, 0, 0, "StripHitMap.png");
1081 
1082       gStyle->SetPalette(kRainBow);
1083       pixelmap->beautifyAllHistograms();
1084 
1085       TCanvas cB("CanvBarrel", "CanvBarrel", 1200, 1000);
1086       pixelmap->drawBarrelMaps("entriesBarrel", cB);
1087       cB.SaveAs("pixelBarrelEntries.png");
1088 
1089       TCanvas cF("CanvForward", "CanvForward", 1600, 1000);
1090       pixelmap->drawForwardMaps("entriesForward", cF);
1091       cF.SaveAs("pixelForwardEntries.png");
1092 
1093       TCanvas cRocs = TCanvas("cRocs", "cRocs", 1200, 1600);
1094       pixelrocsmap_->drawMaps(cRocs, "Pixel on-track clusters occupancy");
1095       cRocs.SaveAs("Phase1PixelROCMaps_fullROCs.png");
1096     }
1097   }
1098 
1099   //*************************************************************
1100   bool isHit2D(const TrackingRecHit &hit)
1101   //*************************************************************
1102   {
1103     bool countStereoHitAs2D_ = true;
1104     // we count SiStrip stereo modules as 2D if selected via countStereoHitAs2D_
1105     // (since they provide theta information)
1106     if (!hit.isValid() ||
1107         (hit.dimension() < 2 && !countStereoHitAs2D_ && !dynamic_cast<const SiStripRecHit1D *>(&hit))) {
1108       return false;  // real RecHit1D - but SiStripRecHit1D depends on countStereoHitAs2D_
1109     } else {
1110       const DetId detId(hit.geographicalId());
1111       if (detId.det() == DetId::Tracker) {
1112         if (detId.subdetId() == kBPIX || detId.subdetId() == kFPIX) {
1113           return true;  // pixel is always 2D
1114         } else {        // should be SiStrip now
1115           const SiStripDetId stripId(detId);
1116           if (stripId.stereo())
1117             return countStereoHitAs2D_;  // stereo modules
1118           else if (dynamic_cast<const SiStripRecHit1D *>(&hit) || dynamic_cast<const SiStripRecHit2D *>(&hit))
1119             return false;  // rphi modules hit
1120           //the following two are not used any more since ages...
1121           else if (dynamic_cast<const SiStripMatchedRecHit2D *>(&hit))
1122             return true;  // matched is 2D
1123           else if (dynamic_cast<const ProjectedSiStripRecHit2D *>(&hit)) {
1124             const ProjectedSiStripRecHit2D *pH = static_cast<const ProjectedSiStripRecHit2D *>(&hit);
1125             return (countStereoHitAs2D_ && this->isHit2D(pH->originalHit()));  // depends on original...
1126           } else {
1127             edm::LogError("UnknownType") << "@SUB=GeneralPurposeTrackAnalyzer::isHit2D"
1128                                          << "Tracker hit not in pixel, neither SiStripRecHit[12]D nor "
1129                                          << "SiStripMatchedRecHit2D nor ProjectedSiStripRecHit2D.";
1130             return false;
1131           }
1132         }
1133       } else {  // not tracker??
1134         edm::LogWarning("DetectorMismatch") << "@SUB=GeneralPurposeTrackAnalyzer::isHit2D"
1135                                             << "Hit not in tracker with 'official' dimension >=2.";
1136         return true;  // dimension() >= 2 so accept that...
1137       }
1138     }
1139     // never reached...
1140   }
1141 };
1142 
1143 //*************************************************************
1144 void GeneralPurposeTrackAnalyzer::fillDescriptions(edm::ConfigurationDescriptions &descriptions)
1145 //*************************************************************
1146 {
1147   edm::ParameterSetDescription desc;
1148   desc.setComment("Generic track analyzer to check ALCARECO sample quantities");
1149   desc.add<edm::InputTag>("TkTag", edm::InputTag("generalTracks"));
1150   desc.add<edm::InputTag>("TriggerResultsTag", edm::InputTag("TriggerResults", "", "HLT"));
1151   desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("offlineBeamSpot"));
1152   desc.add<edm::InputTag>("VerticesTag", edm::InputTag("offlinePrimaryVertices"));
1153   desc.add<bool>("isCosmics", false);
1154   desc.add<bool>("doLatencyAnalysis", true);
1155   descriptions.addWithDefaultLabel(desc);
1156 }
1157 
1158 DEFINE_FWK_MODULE(GeneralPurposeTrackAnalyzer);