Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:57:03

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