Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-04-22 22:55:13

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