Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:08

0001 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
0002 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0003 #include "TrackingTools/Records/interface/TransientRecHitRecord.h"
0004 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0005 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0006 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
0007 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
0008 #include "RecoTracker/TransientTrackingRecHit/interface/TkTransientTrackingRecHitBuilder.h"
0009 #include "RecoTracker/TkTrackingRegions/interface/RectangularEtaPhiTrackingRegion.h"
0010 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0011 #include "MagneticField/Engine/interface/MagneticField.h"
0012 
0013 #include "DataFormats/Candidate/interface/Candidate.h"
0014 #include "DataFormats/Math/interface/deltaPhi.h"
0015 
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "FWCore/Utilities/interface/InputTag.h"
0018 #include "DQMServices/Core/interface/DQMStore.h"
0019 #include "DQM/TrackingMonitor/interface/TrackBuildingAnalyzer.h"
0020 #include <string>
0021 #include "TMath.h"
0022 
0023 TrackBuildingAnalyzer::TrackBuildingAnalyzer(const edm::ParameterSet& iConfig)
0024     : doAllPlots(iConfig.getParameter<bool>("doAllPlots")),
0025       doAllSeedPlots(iConfig.getParameter<bool>("doSeedParameterHistos")),
0026       doTCPlots(iConfig.getParameter<bool>("doTrackCandHistos")),
0027       doAllTCPlots(iConfig.getParameter<bool>("doAllTrackCandHistos")),
0028       doPT(iConfig.getParameter<bool>("doSeedPTHisto")),
0029       doETA(iConfig.getParameter<bool>("doSeedETAHisto")),
0030       doPHI(iConfig.getParameter<bool>("doSeedPHIHisto")),
0031       doPHIVsETA(iConfig.getParameter<bool>("doSeedPHIVsETAHisto")),
0032       doTheta(iConfig.getParameter<bool>("doSeedThetaHisto")),
0033       doQ(iConfig.getParameter<bool>("doSeedQHisto")),
0034       doDxy(iConfig.getParameter<bool>("doSeedDxyHisto")),
0035       doDz(iConfig.getParameter<bool>("doSeedDzHisto")),
0036       doNRecHits(iConfig.getParameter<bool>("doSeedNRecHitsHisto")),
0037       doProfPHI(iConfig.getParameter<bool>("doSeedNVsPhiProf")),
0038       doProfETA(iConfig.getParameter<bool>("doSeedNVsEtaProf")),
0039       doStopSource(iConfig.getParameter<bool>("doStopSource")),
0040       doMVAPlots(iConfig.getParameter<bool>("doMVAPlots")),
0041       doRegionPlots(iConfig.getParameter<bool>("doRegionPlots")),
0042       doRegionCandidatePlots(iConfig.getParameter<bool>("doRegionCandidatePlots")) {}
0043 
0044 void TrackBuildingAnalyzer::initHisto(DQMStore::IBooker& ibooker, const edm::ParameterSet& iConfig) {
0045   // parameters from the configuration
0046   std::string AlgoName = iConfig.getParameter<std::string>("AlgoName");
0047   std::string MEFolderName = iConfig.getParameter<std::string>("FolderName");
0048 
0049   //  std::cout << "[TrackBuildingAnalyzer::beginRun] AlgoName: " << AlgoName << std::endl;
0050 
0051   // use the AlgoName and Quality Name
0052   const std::string& CatagoryName = AlgoName;
0053 
0054   // get binning from the configuration
0055   int TrackPtBin = iConfig.getParameter<int>("TrackPtBin");
0056   double TrackPtMin = iConfig.getParameter<double>("TrackPtMin");
0057   double TrackPtMax = iConfig.getParameter<double>("TrackPtMax");
0058 
0059   int PhiBin = iConfig.getParameter<int>("PhiBin");
0060   double PhiMin = iConfig.getParameter<double>("PhiMin");
0061   double PhiMax = iConfig.getParameter<double>("PhiMax");
0062   phiBinWidth = PhiBin > 0 ? (PhiMax - PhiMin) / PhiBin : 0.;
0063 
0064   int EtaBin = iConfig.getParameter<int>("EtaBin");
0065   double EtaMin = iConfig.getParameter<double>("EtaMin");
0066   double EtaMax = iConfig.getParameter<double>("EtaMax");
0067   etaBinWidth = EtaBin > 0 ? (EtaMax - EtaMin) / EtaBin : 0.;
0068 
0069   int ThetaBin = iConfig.getParameter<int>("ThetaBin");
0070   double ThetaMin = iConfig.getParameter<double>("ThetaMin");
0071   double ThetaMax = iConfig.getParameter<double>("ThetaMax");
0072 
0073   int TrackQBin = iConfig.getParameter<int>("TrackQBin");
0074   double TrackQMin = iConfig.getParameter<double>("TrackQMin");
0075   double TrackQMax = iConfig.getParameter<double>("TrackQMax");
0076 
0077   int SeedDxyBin = iConfig.getParameter<int>("SeedDxyBin");
0078   double SeedDxyMin = iConfig.getParameter<double>("SeedDxyMin");
0079   double SeedDxyMax = iConfig.getParameter<double>("SeedDxyMax");
0080 
0081   int SeedDzBin = iConfig.getParameter<int>("SeedDzBin");
0082   double SeedDzMin = iConfig.getParameter<double>("SeedDzMin");
0083   double SeedDzMax = iConfig.getParameter<double>("SeedDzMax");
0084 
0085   int SeedHitBin = iConfig.getParameter<int>("SeedHitBin");
0086   double SeedHitMin = iConfig.getParameter<double>("SeedHitMin");
0087   double SeedHitMax = iConfig.getParameter<double>("SeedHitMax");
0088 
0089   int TCDxyBin = iConfig.getParameter<int>("TCDxyBin");
0090   double TCDxyMin = iConfig.getParameter<double>("TCDxyMin");
0091   double TCDxyMax = iConfig.getParameter<double>("TCDxyMax");
0092 
0093   int TCDzBin = iConfig.getParameter<int>("TCDzBin");
0094   double TCDzMin = iConfig.getParameter<double>("TCDzMin");
0095   double TCDzMax = iConfig.getParameter<double>("TCDzMax");
0096 
0097   int TCHitBin = iConfig.getParameter<int>("TCHitBin");
0098   double TCHitMin = iConfig.getParameter<double>("TCHitMin");
0099   double TCHitMax = iConfig.getParameter<double>("TCHitMax");
0100 
0101   int MVABin = iConfig.getParameter<int>("MVABin");
0102   double MVAMin = iConfig.getParameter<double>("MVAMin");
0103   double MVAMax = iConfig.getParameter<double>("MVAMax");
0104 
0105   edm::InputTag seedProducer = iConfig.getParameter<edm::InputTag>("SeedProducer");
0106   edm::InputTag tcProducer = iConfig.getParameter<edm::InputTag>("TCProducer");
0107   std::vector<std::string> mvaProducers = iConfig.getParameter<std::vector<std::string>>("MVAProducers");
0108   edm::InputTag regionProducer = iConfig.getParameter<edm::InputTag>("RegionProducer");
0109 
0110   //    if (doAllPlots){doAllSeedPlots=true; doTCPlots=true;}
0111 
0112   ibooker.setCurrentFolder(MEFolderName);
0113 
0114   // book the Seed histograms
0115   // ---------------------------------------------------------------------------------//
0116   //  std::cout << "[TrackBuildingAnalyzer::beginRun] MEFolderName: " << MEFolderName << std::endl;
0117   ibooker.setCurrentFolder(MEFolderName + "/TrackBuilding");
0118 
0119   if (doAllSeedPlots || doPT) {
0120     histname = "SeedPt_" + seedProducer.label() + "_";
0121     SeedPt = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, TrackPtBin, TrackPtMin, TrackPtMax);
0122     SeedPt->setAxisTitle("Seed p_{T} (GeV/c)", 1);
0123     SeedPt->setAxisTitle("Number of Seeds", 2);
0124   }
0125 
0126   if (doAllSeedPlots || doETA) {
0127     histname = "SeedEta_" + seedProducer.label() + "_";
0128     SeedEta = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax);
0129     SeedEta->setAxisTitle("Seed #eta", 1);
0130     SeedEta->setAxisTitle("Number of Seeds", 2);
0131   }
0132 
0133   if (doAllSeedPlots || doPHI) {
0134     histname = "SeedPhi_" + seedProducer.label() + "_";
0135     SeedPhi = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, PhiBin, PhiMin, PhiMax);
0136     SeedPhi->setAxisTitle("Seed #phi", 1);
0137     SeedPhi->setAxisTitle("Number of Seed", 2);
0138   }
0139 
0140   if (doAllSeedPlots || doPHIVsETA) {
0141     histname = "SeedPhiVsEta_" + seedProducer.label() + "_";
0142     SeedPhiVsEta = ibooker.book2D(
0143         histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax);
0144     SeedPhiVsEta->setAxisTitle("Seed #eta", 1);
0145     SeedPhiVsEta->setAxisTitle("Seed #phi", 2);
0146   }
0147 
0148   if (doAllSeedPlots || doTheta) {
0149     histname = "SeedTheta_" + seedProducer.label() + "_";
0150     SeedTheta = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, ThetaBin, ThetaMin, ThetaMax);
0151     SeedTheta->setAxisTitle("Seed #theta", 1);
0152     SeedTheta->setAxisTitle("Number of Seeds", 2);
0153   }
0154 
0155   if (doAllSeedPlots || doQ) {
0156     histname = "SeedQ_" + seedProducer.label() + "_";
0157     SeedQ = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, TrackQBin, TrackQMin, TrackQMax);
0158     SeedQ->setAxisTitle("Seed Charge", 1);
0159     SeedQ->setAxisTitle("Number of Seeds", 2);
0160   }
0161 
0162   if (doAllSeedPlots || doDxy) {
0163     histname = "SeedDxy_" + seedProducer.label() + "_";
0164     SeedDxy = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, SeedDxyBin, SeedDxyMin, SeedDxyMax);
0165     SeedDxy->setAxisTitle("Seed d_{xy} (cm)", 1);
0166     SeedDxy->setAxisTitle("Number of Seeds", 2);
0167   }
0168 
0169   if (doAllSeedPlots || doDz) {
0170     histname = "SeedDz_" + seedProducer.label() + "_";
0171     SeedDz = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, SeedDzBin, SeedDzMin, SeedDzMax);
0172     SeedDz->setAxisTitle("Seed d_{z} (cm)", 1);
0173     SeedDz->setAxisTitle("Number of Seeds", 2);
0174   }
0175 
0176   if (doAllSeedPlots || doNRecHits) {
0177     histname = "NumberOfRecHitsPerSeed_" + seedProducer.label() + "_";
0178     NumberOfRecHitsPerSeed =
0179         ibooker.book1D(histname + CatagoryName, histname + CatagoryName, SeedHitBin, SeedHitMin, SeedHitMax);
0180     NumberOfRecHitsPerSeed->setAxisTitle("Number of RecHits per Seed", 1);
0181     NumberOfRecHitsPerSeed->setAxisTitle("Number of Seeds", 2);
0182   }
0183 
0184   if (doAllSeedPlots || doProfPHI) {
0185     histname = "NumberOfRecHitsPerSeedVsPhiProfile_" + seedProducer.label() + "_";
0186     NumberOfRecHitsPerSeedVsPhiProfile = ibooker.bookProfile(histname + CatagoryName,
0187                                                              histname + CatagoryName,
0188                                                              PhiBin,
0189                                                              PhiMin,
0190                                                              PhiMax,
0191                                                              SeedHitBin,
0192                                                              SeedHitMin,
0193                                                              SeedHitMax,
0194                                                              "s");
0195     NumberOfRecHitsPerSeedVsPhiProfile->setAxisTitle("Seed #phi", 1);
0196     NumberOfRecHitsPerSeedVsPhiProfile->setAxisTitle("Number of RecHits of each Seed", 2);
0197   }
0198 
0199   if (doAllSeedPlots || doProfETA) {
0200     histname = "NumberOfRecHitsPerSeedVsEtaProfile_" + seedProducer.label() + "_";
0201     NumberOfRecHitsPerSeedVsEtaProfile = ibooker.bookProfile(histname + CatagoryName,
0202                                                              histname + CatagoryName,
0203                                                              EtaBin,
0204                                                              EtaMin,
0205                                                              EtaMax,
0206                                                              SeedHitBin,
0207                                                              SeedHitMin,
0208                                                              SeedHitMax,
0209                                                              "s");
0210     NumberOfRecHitsPerSeedVsEtaProfile->setAxisTitle("Seed #eta", 1);
0211     NumberOfRecHitsPerSeedVsEtaProfile->setAxisTitle("Number of RecHits of each Seed", 2);
0212   }
0213 
0214   if (doRegionPlots) {
0215     if (doAllSeedPlots || doETA) {
0216       histname = "TrackingRegionEta_" + seedProducer.label() + "_";
0217       TrackingRegionEta = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax);
0218       TrackingRegionEta->setAxisTitle("TrackingRegion-covered #eta", 1);
0219       TrackingRegionEta->setAxisTitle("Number of TrackingRegions", 2);
0220     }
0221 
0222     if (doAllSeedPlots || doPHI) {
0223       histname = "TrackingRegionPhi_" + seedProducer.label() + "_";
0224       TrackingRegionPhi = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, PhiBin, PhiMin, PhiMax);
0225       TrackingRegionPhi->setAxisTitle("TrackingRegion-covered #phi", 1);
0226       TrackingRegionPhi->setAxisTitle("Number of TrackingRegions", 2);
0227     }
0228 
0229     if (doAllSeedPlots || doPHIVsETA) {
0230       histname = "TrackingRegionPhiVsEta_" + seedProducer.label() + "_";
0231       TrackingRegionPhiVsEta = ibooker.book2D(
0232           histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax);
0233       TrackingRegionPhiVsEta->setAxisTitle("TrackingRegion-covered #eta", 1);
0234       TrackingRegionPhiVsEta->setAxisTitle("TrackingRegion-covered #phi", 2);
0235     }
0236 
0237     if (doRegionCandidatePlots) {
0238       if (doAllSeedPlots || doPT) {
0239         auto ptBin = iConfig.getParameter<int>("RegionCandidatePtBin");
0240         auto ptMin = iConfig.getParameter<double>("RegionCandidatePtMin");
0241         auto ptMax = iConfig.getParameter<double>("RegionCandidatePtMax");
0242 
0243         histname = "TrackingRegionCandidatePt_" + seedProducer.label() + "_";
0244         TrackingRegionCandidatePt =
0245             ibooker.book1D(histname + CatagoryName, histname + CatagoryName, ptBin, ptMin, ptMax);
0246         TrackingRegionCandidatePt->setAxisTitle("TrackingRegion Candidate p_{T} (GeV/c)", 1);
0247         TrackingRegionCandidatePt->setAxisTitle("Number of TrackingRegion Candidates", 2);
0248       }
0249 
0250       if (doAllSeedPlots || doETA) {
0251         histname = "TrackingRegionCandidateEta_" + seedProducer.label() + "_";
0252         TrackingRegionCandidateEta =
0253             ibooker.book1D(histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax);
0254         TrackingRegionCandidateEta->setAxisTitle("TrackingRegion Candidate #eta", 1);
0255         TrackingRegionCandidateEta->setAxisTitle("Number of TrackingRegion Candidates", 2);
0256       }
0257 
0258       if (doAllSeedPlots || doPHI) {
0259         histname = "TrackingRegionCandidatePhi_" + seedProducer.label() + "_";
0260         TrackingRegionCandidatePhi =
0261             ibooker.book1D(histname + CatagoryName, histname + CatagoryName, PhiBin, PhiMin, PhiMax);
0262         TrackingRegionCandidatePhi->setAxisTitle("TrackingRegion Candidate #phi", 1);
0263         TrackingRegionCandidatePhi->setAxisTitle("Number of TrackingRegion Candidates", 2);
0264       }
0265 
0266       if (doAllSeedPlots || doPHIVsETA) {
0267         histname = "TrackingRegionCandidatePhiVsEta_" + seedProducer.label() + "_";
0268         TrackingRegionCandidatePhiVsEta = ibooker.book2D(
0269             histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax);
0270         TrackingRegionCandidatePhiVsEta->setAxisTitle("TrackingRegion Candidate #eta", 1);
0271         TrackingRegionCandidatePhiVsEta->setAxisTitle("TrackingRegion Candidate #phi", 2);
0272       }
0273     }
0274   }
0275 
0276   if (doAllSeedPlots || doStopSource) {
0277     const auto stopReasonSize = static_cast<unsigned int>(SeedStopReason::SIZE);
0278 
0279     const auto candsBin = iConfig.getParameter<int>("SeedCandBin");
0280     const auto candsMin = iConfig.getParameter<double>("SeedCandMin");
0281     const auto candsMax = iConfig.getParameter<double>("SeedCandMax");
0282 
0283     histname = "SeedStoppingSource_" + seedProducer.label() + "_";
0284     seedStoppingSource =
0285         ibooker.book1D(histname + CatagoryName, histname + CatagoryName, stopReasonSize, 0., stopReasonSize);
0286     seedStoppingSource->setAxisTitle("Stopping reason", 1);
0287     seedStoppingSource->setAxisTitle("Number of seeds", 2);
0288 
0289     histname = "SeedStoppingSourceVsPhi_" + seedProducer.label() + "_";
0290     seedStoppingSourceVsPhi =
0291         ibooker.bookProfile(histname + CatagoryName, histname + CatagoryName, PhiBin, PhiMin, PhiMax, 2, 0., 2.);
0292     seedStoppingSourceVsPhi->setAxisTitle("seed #phi", 1);
0293     seedStoppingSourceVsPhi->setAxisTitle("fraction stopped", 2);
0294 
0295     histname = "SeedStoppingSourceVsEta_" + seedProducer.label() + "_";
0296     seedStoppingSourceVsEta =
0297         ibooker.bookProfile(histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax, 2, 0., 2.);
0298     seedStoppingSourceVsEta->setAxisTitle("seed #eta", 1);
0299     seedStoppingSourceVsEta->setAxisTitle("fraction stopped", 2);
0300 
0301     histname = "NumberOfTrajCandsPerSeed_" + seedProducer.label() + "_";
0302     numberOfTrajCandsPerSeed =
0303         ibooker.book1D(histname + CatagoryName, histname + CatagoryName, candsBin, candsMin, candsMax);
0304     numberOfTrajCandsPerSeed->setAxisTitle("Number of Trajectory Candidate for each Seed", 1);
0305     numberOfTrajCandsPerSeed->setAxisTitle("Number of Seeds", 2);
0306 
0307     histname = "NumberOfTrajCandsPerSeedVsPhi_" + seedProducer.label() + "_";
0308     numberOfTrajCandsPerSeedVsPhi = ibooker.bookProfile(
0309         histname + CatagoryName, histname + CatagoryName, PhiBin, PhiMin, PhiMax, candsBin, candsMin, candsMax);
0310     numberOfTrajCandsPerSeedVsPhi->setAxisTitle("Seed #phi", 1);
0311     numberOfTrajCandsPerSeedVsPhi->setAxisTitle("Number of Trajectory Candidates for each Seed", 2);
0312 
0313     histname = "NumberOfTrajCandsPerSeedVsEta_" + seedProducer.label() + "_";
0314     numberOfTrajCandsPerSeedVsEta = ibooker.bookProfile(
0315         histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax, candsBin, candsMin, candsMax);
0316     numberOfTrajCandsPerSeedVsEta->setAxisTitle("Seed #eta", 1);
0317     numberOfTrajCandsPerSeedVsEta->setAxisTitle("Number of Trajectory Candidates for each Seed", 2);
0318 
0319     histname = "SeedStoppingSourceVsNumberOfTrajCandsPerSeed_" + seedProducer.label() + "_";
0320     seedStoppingSourceVsNumberOfTrajCandsPerSeed =
0321         ibooker.bookProfile(histname + CatagoryName, histname + CatagoryName, candsBin, candsMin, candsMax, 2, 0., 2.);
0322     seedStoppingSourceVsNumberOfTrajCandsPerSeed->setAxisTitle("Number of Trajectory Candidates for each Seed", 1);
0323     seedStoppingSourceVsNumberOfTrajCandsPerSeed->setAxisTitle("fraction stopped", 2);
0324 
0325     for (unsigned int i = 0; i < stopReasonSize; ++i) {
0326       seedStoppingSource->setBinLabel(i + 1, SeedStopReasonName::SeedStopReasonName[i], 1);
0327     }
0328   }
0329 
0330   if (doAllTCPlots || doStopSource) {
0331     // DataFormats/TrackReco/interface/TrajectoryStopReasons.h
0332     size_t StopReasonNameSize = static_cast<size_t>(StopReason::SIZE);
0333 
0334     histname = "StoppingSource_" + seedProducer.label() + "_";
0335     stoppingSource = ibooker.book1D(
0336         histname + CatagoryName, histname + CatagoryName, StopReasonNameSize, 0., double(StopReasonNameSize));
0337     stoppingSource->setAxisTitle("stopping reason", 1);
0338     stoppingSource->setAxisTitle("Number of Tracks", 2);
0339 
0340     histname = "StoppingSourceVSeta_" + seedProducer.label() + "_";
0341     stoppingSourceVSeta =
0342         ibooker.bookProfile(histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax, 2, 0., 2.);
0343     stoppingSourceVSeta->setAxisTitle("track #eta", 1);
0344     stoppingSourceVSeta->setAxisTitle("fraction stopped", 2);
0345 
0346     histname = "StoppingSourceVSphi_" + seedProducer.label() + "_";
0347     stoppingSourceVSphi =
0348         ibooker.bookProfile(histname + CatagoryName, histname + CatagoryName, PhiBin, PhiMin, PhiMax, 2, 0., 2.);
0349     stoppingSourceVSphi->setAxisTitle("track #phi", 1);
0350     stoppingSourceVSphi->setAxisTitle("fraction stopped", 2);
0351 
0352     for (size_t ibin = 0; ibin < StopReasonNameSize; ibin++) {
0353       stoppingSource->setBinLabel(ibin + 1, StopReasonName::StopReasonName[ibin], 1);
0354     }
0355   }
0356 
0357   // book the TrackCandidate histograms
0358   // ---------------------------------------------------------------------------------//
0359 
0360   if (doTCPlots) {
0361     ibooker.setCurrentFolder(MEFolderName + "/TrackBuilding");
0362 
0363     histname = "TrackCandPt_" + tcProducer.label() + "_";
0364     TrackCandPt = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, TrackPtBin, TrackPtMin, TrackPtMax);
0365     TrackCandPt->setAxisTitle("Track Candidate p_{T} (GeV/c)", 1);
0366     TrackCandPt->setAxisTitle("Number of Track Candidates", 2);
0367 
0368     histname = "TrackCandEta_" + tcProducer.label() + "_";
0369     TrackCandEta = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax);
0370     TrackCandEta->setAxisTitle("Track Candidate #eta", 1);
0371     TrackCandEta->setAxisTitle("Number of Track Candidates", 2);
0372 
0373     histname = "TrackCandPhi_" + tcProducer.label() + "_";
0374     TrackCandPhi = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, PhiBin, PhiMin, PhiMax);
0375     TrackCandPhi->setAxisTitle("Track Candidate #phi", 1);
0376     TrackCandPhi->setAxisTitle("Number of Track Candidates", 2);
0377 
0378     if (doTheta) {
0379       histname = "TrackCandTheta_" + tcProducer.label() + "_";
0380       TrackCandTheta = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, ThetaBin, ThetaMin, ThetaMax);
0381       TrackCandTheta->setAxisTitle("Track Candidate #theta", 1);
0382       TrackCandTheta->setAxisTitle("Number of Track Candidates", 2);
0383     }
0384 
0385     if (doAllTCPlots) {
0386       histname = "TrackCandQ_" + tcProducer.label() + "_";
0387       TrackCandQ = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, TrackQBin, TrackQMin, TrackQMax);
0388       TrackCandQ->setAxisTitle("Track Candidate Charge", 1);
0389       TrackCandQ->setAxisTitle("Number of Track Candidates", 2);
0390 
0391       histname = "TrackCandDxy_" + tcProducer.label() + "_";
0392       TrackCandDxy = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, TCDxyBin, TCDxyMin, TCDxyMax);
0393       TrackCandDxy->setAxisTitle("Track Candidate d_{xy} (cm)", 1);
0394       TrackCandDxy->setAxisTitle("Number of Track Candidates", 2);
0395 
0396       histname = "TrackCandDz_" + tcProducer.label() + "_";
0397       TrackCandDz = ibooker.book1D(histname + CatagoryName, histname + CatagoryName, TCDzBin, TCDzMin, TCDzMax);
0398       TrackCandDz->setAxisTitle("Track Candidate d_{z} (cm)", 1);
0399       TrackCandDz->setAxisTitle("Number of Track Candidates", 2);
0400 
0401       histname = "NumberOfRecHitsPerTrackCand_" + tcProducer.label() + "_";
0402       NumberOfRecHitsPerTrackCand =
0403           ibooker.book1D(histname + CatagoryName, histname + CatagoryName, TCHitBin, TCHitMin, TCHitMax);
0404       NumberOfRecHitsPerTrackCand->setAxisTitle("Number of RecHits per Track Candidate", 1);
0405       NumberOfRecHitsPerTrackCand->setAxisTitle("Number of Track Candidates", 2);
0406 
0407       histname = "NumberOfRecHitsPerTrackCandVsPhiProfile_" + tcProducer.label() + "_";
0408       NumberOfRecHitsPerTrackCandVsPhiProfile = ibooker.bookProfile(
0409           histname + CatagoryName, histname + CatagoryName, PhiBin, PhiMin, PhiMax, TCHitBin, TCHitMin, TCHitMax, "s");
0410       NumberOfRecHitsPerTrackCandVsPhiProfile->setAxisTitle("Track Candidate #phi", 1);
0411       NumberOfRecHitsPerTrackCandVsPhiProfile->setAxisTitle("Number of RecHits of each Track Candidate", 2);
0412 
0413       histname = "NumberOfRecHitsPerTrackCandVsEtaProfile_" + tcProducer.label() + "_";
0414       NumberOfRecHitsPerTrackCandVsEtaProfile = ibooker.bookProfile(
0415           histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax, TCHitBin, TCHitMin, TCHitMax, "s");
0416       NumberOfRecHitsPerTrackCandVsEtaProfile->setAxisTitle("Track Candidate #eta", 1);
0417       NumberOfRecHitsPerTrackCandVsEtaProfile->setAxisTitle("Number of RecHits of each Track Candidate", 2);
0418     }
0419 
0420     histname = "TrackCandPhiVsEta_" + tcProducer.label() + "_";
0421     TrackCandPhiVsEta = ibooker.book2D(
0422         histname + CatagoryName, histname + CatagoryName, EtaBin, EtaMin, EtaMax, PhiBin, PhiMin, PhiMax);
0423     TrackCandPhiVsEta->setAxisTitle("Track Candidate #eta", 1);
0424     TrackCandPhiVsEta->setAxisTitle("Track Candidate #phi", 2);
0425 
0426     if (doAllTCPlots || doMVAPlots) {
0427       for (size_t i = 1, end = mvaProducers.size(); i <= end; ++i) {
0428         auto num = std::to_string(i);
0429         std::string pfix;
0430 
0431         if (i == 1) {
0432           trackMVAsHP.push_back(nullptr);
0433           trackMVAsHPVsPtProfile.push_back(nullptr);
0434           trackMVAsHPVsEtaProfile.push_back(nullptr);
0435         } else {
0436           pfix = " (not loose-selected)";
0437           std::string pfix2 = " (not HP-selected)";
0438           histname = "TrackMVA" + num + "HP_" + tcProducer.label() + "_";
0439           trackMVAsHP.push_back(
0440               ibooker.book1D(histname + CatagoryName, histname + CatagoryName + pfix2, MVABin, MVAMin, MVAMax));
0441           trackMVAsHP.back()->setAxisTitle("Track selection MVA" + num, 1);
0442           trackMVAsHP.back()->setAxisTitle("Number of tracks", 2);
0443 
0444           histname = "TrackMVA" + num + "HPVsPtProfile_" + tcProducer.label() + "_";
0445           trackMVAsHPVsPtProfile.push_back(ibooker.bookProfile(histname + CatagoryName,
0446                                                                histname + CatagoryName + pfix2,
0447                                                                TrackPtBin,
0448                                                                TrackPtMin,
0449                                                                TrackPtMax,
0450                                                                MVABin,
0451                                                                MVAMin,
0452                                                                MVAMax));
0453           trackMVAsHPVsPtProfile.back()->setAxisTitle("Track p_{T} (GeV/c)", 1);
0454           trackMVAsHPVsPtProfile.back()->setAxisTitle("Track selection MVA" + num, 2);
0455 
0456           histname = "TrackMVA" + num + "HPVsEtaProfile_" + tcProducer.label() + "_";
0457           trackMVAsHPVsEtaProfile.push_back(ibooker.bookProfile(
0458               histname + CatagoryName, histname + CatagoryName + pfix2, EtaBin, EtaMin, EtaMax, MVABin, MVAMin, MVAMax));
0459           trackMVAsHPVsEtaProfile.back()->setAxisTitle("Track #eta", 1);
0460           trackMVAsHPVsEtaProfile.back()->setAxisTitle("Track selection MVA" + num, 2);
0461         }
0462 
0463         histname = "TrackMVA" + num + "_" + tcProducer.label() + "_";
0464         trackMVAs.push_back(
0465             ibooker.book1D(histname + CatagoryName, histname + CatagoryName + pfix, MVABin, MVAMin, MVAMax));
0466         trackMVAs.back()->setAxisTitle("Track selection MVA" + num, 1);
0467         trackMVAs.back()->setAxisTitle("Number of tracks", 2);
0468 
0469         histname = "TrackMVA" + num + "VsPtProfile_" + tcProducer.label() + "_";
0470         trackMVAsVsPtProfile.push_back(ibooker.bookProfile(histname + CatagoryName,
0471                                                            histname + CatagoryName + pfix,
0472                                                            TrackPtBin,
0473                                                            TrackPtMin,
0474                                                            TrackPtMax,
0475                                                            MVABin,
0476                                                            MVAMin,
0477                                                            MVAMax));
0478         trackMVAsVsPtProfile.back()->setAxisTitle("Track p_{T} (GeV/c)", 1);
0479         trackMVAsVsPtProfile.back()->setAxisTitle("Track selection MVA" + num, 2);
0480 
0481         histname = "TrackMVA" + num + "VsEtaProfile_" + tcProducer.label() + "_";
0482         trackMVAsVsEtaProfile.push_back(ibooker.bookProfile(
0483             histname + CatagoryName, histname + CatagoryName + pfix, EtaBin, EtaMin, EtaMax, MVABin, MVAMin, MVAMax));
0484         trackMVAsVsEtaProfile.back()->setAxisTitle("Track #eta", 1);
0485         trackMVAsVsEtaProfile.back()->setAxisTitle("Track selection MVA" + num, 2);
0486       }
0487     }
0488   }
0489 }
0490 
0491 // -- Analyse
0492 // ---------------------------------------------------------------------------------//
0493 void TrackBuildingAnalyzer::analyze(const edm::Event& iEvent,
0494                                     const edm::EventSetup& iSetup,
0495                                     const TrajectorySeed& candidate,
0496                                     const SeedStopInfo& stopInfo,
0497                                     const reco::BeamSpot& bs,
0498                                     const MagneticField& theMF,
0499                                     const TransientTrackingRecHitBuilder& theTTRHBuilder) {
0500   TSCBLBuilderNoMaterial tscblBuilder;
0501 
0502   //get parameters and errors from the candidate state
0503   auto const& theG = ((TkTransientTrackingRecHitBuilder const*)(&theTTRHBuilder))->geometry();
0504   auto const& candSS = candidate.startingState();
0505   TrajectoryStateOnSurface state =
0506       trajectoryStateTransform::transientState(candSS, &(theG->idToDet(candSS.detId())->surface()), &theMF);
0507   TrajectoryStateClosestToBeamLine tsAtClosestApproachSeed =
0508       tscblBuilder(*state.freeState(), bs);  //as in TrackProducerAlgorithm
0509   if (!(tsAtClosestApproachSeed.isValid())) {
0510     edm::LogVerbatim("TrackBuilding") << "TrajectoryStateClosestToBeamLine not valid";
0511     return;
0512   }
0513   GlobalPoint v0 = tsAtClosestApproachSeed.trackStateAtPCA().position();
0514   GlobalVector p = tsAtClosestApproachSeed.trackStateAtPCA().momentum();
0515   GlobalPoint v(v0.x() - bs.x0(), v0.y() - bs.y0(), v0.z() - bs.z0());
0516 
0517   double pt = sqrt(state.globalMomentum().perp2());
0518   double eta = state.globalPosition().eta();
0519   double phi = state.globalPosition().phi();
0520   double theta = state.globalPosition().theta();
0521   //double pm           = sqrt(state.globalMomentum().mag2());
0522   //double pz           = state.globalMomentum().z();
0523   //double qoverp       = tsAtClosestApproachSeed.trackStateAtPCA().charge()/p.mag();
0524   //double theta        = p.theta();
0525   //double lambda       = M_PI/2-p.theta();
0526   double numberOfHits = candidate.nHits();
0527   double dxy = (-v.x() * sin(p.phi()) + v.y() * cos(p.phi()));
0528   double dz = v.z() - (v.x() * p.x() + v.y() * p.y()) / p.perp() * p.z() / p.perp();
0529 
0530   // fill the ME's
0531   if (doAllSeedPlots || doQ)
0532     SeedQ->Fill(state.charge());
0533   if (doAllSeedPlots || doPT)
0534     SeedPt->Fill(pt);
0535   if (doAllSeedPlots || doETA)
0536     SeedEta->Fill(eta);
0537   if (doAllSeedPlots || doPHI)
0538     SeedPhi->Fill(phi);
0539   if (doAllSeedPlots || doPHIVsETA)
0540     SeedPhiVsEta->Fill(eta, phi);
0541   if (doAllSeedPlots || doTheta)
0542     SeedTheta->Fill(theta);
0543   if (doAllSeedPlots || doDxy)
0544     SeedDxy->Fill(dxy);
0545   if (doAllSeedPlots || doDz)
0546     SeedDz->Fill(dz);
0547   if (doAllSeedPlots || doNRecHits)
0548     NumberOfRecHitsPerSeed->Fill(numberOfHits);
0549   if (doAllSeedPlots || doProfETA)
0550     NumberOfRecHitsPerSeedVsEtaProfile->Fill(eta, numberOfHits);
0551   if (doAllSeedPlots || doProfPHI)
0552     NumberOfRecHitsPerSeedVsPhiProfile->Fill(phi, numberOfHits);
0553   if (doAllSeedPlots || doStopSource) {
0554     const double stopped = stopInfo.stopReason() == SeedStopReason::NOT_STOPPED ? 0. : 1.;
0555     seedStoppingSource->Fill(stopInfo.stopReasonUC());
0556     seedStoppingSourceVsPhi->Fill(phi, stopped);
0557     seedStoppingSourceVsEta->Fill(eta, stopped);
0558 
0559     const auto ncands = stopInfo.candidatesPerSeed();
0560     numberOfTrajCandsPerSeed->Fill(ncands);
0561     numberOfTrajCandsPerSeedVsPhi->Fill(phi, ncands);
0562     numberOfTrajCandsPerSeedVsEta->Fill(eta, ncands);
0563 
0564     seedStoppingSourceVsNumberOfTrajCandsPerSeed->Fill(ncands, stopped);
0565   }
0566 }
0567 
0568 // -- Analyse
0569 // ---------------------------------------------------------------------------------//
0570 void TrackBuildingAnalyzer::analyze(const edm::Event& iEvent,
0571                                     const edm::EventSetup& iSetup,
0572                                     const TrackCandidate& candidate,
0573                                     const reco::BeamSpot& bs,
0574                                     const MagneticField& theMF,
0575                                     const TransientTrackingRecHitBuilder& theTTRHBuilder) {
0576   TSCBLBuilderNoMaterial tscblBuilder;
0577 
0578   //get parameters and errors from the candidate state
0579   auto const& theG = ((TkTransientTrackingRecHitBuilder const*)(&theTTRHBuilder))->geometry();
0580   auto const& candSS = candidate.trajectoryStateOnDet();
0581   TrajectoryStateOnSurface state =
0582       trajectoryStateTransform::transientState(candSS, &(theG->idToDet(candSS.detId())->surface()), &theMF);
0583   TrajectoryStateClosestToBeamLine tsAtClosestApproachTrackCand =
0584       tscblBuilder(*state.freeState(), bs);  //as in TrackProducerAlgorithm
0585   if (!(tsAtClosestApproachTrackCand.isValid())) {
0586     edm::LogVerbatim("TrackBuilding") << "TrajectoryStateClosestToBeamLine not valid";
0587     return;
0588   }
0589   GlobalPoint v0 = tsAtClosestApproachTrackCand.trackStateAtPCA().position();
0590   GlobalVector p = tsAtClosestApproachTrackCand.trackStateAtPCA().momentum();
0591   GlobalPoint v(v0.x() - bs.x0(), v0.y() - bs.y0(), v0.z() - bs.z0());
0592 
0593   double pt = sqrt(state.globalMomentum().perp2());
0594   double eta = state.globalPosition().eta();
0595   double phi = state.globalPosition().phi();
0596   double theta = state.globalPosition().theta();
0597   //double pm           = sqrt(state.globalMomentum().mag2());
0598   //double pz           = state.globalMomentum().z();
0599   //double qoverp       = tsAtClosestApproachTrackCand.trackStateAtPCA().charge()/p.mag();
0600   //double theta        = p.theta();
0601   //double lambda       = M_PI/2-p.theta();
0602   double numberOfHits = candidate.nRecHits();
0603   double dxy = (-v.x() * sin(p.phi()) + v.y() * cos(p.phi()));
0604 
0605   double dz = v.z() - (v.x() * p.x() + v.y() * p.y()) / p.perp() * p.z() / p.perp();
0606 
0607   if (doAllTCPlots || doStopSource) {
0608     // stopping source
0609     int max = stoppingSource->getNbinsX();
0610     double stop = candidate.stopReason() > max ? double(max - 1) : static_cast<double>(candidate.stopReason());
0611     double stopped = int(StopReason::NOT_STOPPED) == candidate.stopReason() ? 0. : 1.;
0612     stoppingSource->Fill(stop);
0613     stoppingSourceVSeta->Fill(eta, stopped);
0614     stoppingSourceVSphi->Fill(phi, stopped);
0615   }
0616 
0617   if (doTCPlots) {
0618     // fill the ME's
0619     if (doAllTCPlots)
0620       TrackCandQ->Fill(state.charge());
0621     TrackCandPt->Fill(pt);
0622     TrackCandEta->Fill(eta);
0623     TrackCandPhi->Fill(phi);
0624     TrackCandPhiVsEta->Fill(eta, phi);
0625     if (doTheta)
0626       TrackCandTheta->Fill(theta);
0627     if (doAllTCPlots)
0628       TrackCandDxy->Fill(dxy);
0629     if (doAllTCPlots)
0630       TrackCandDz->Fill(dz);
0631     if (doAllTCPlots)
0632       NumberOfRecHitsPerTrackCand->Fill(numberOfHits);
0633     if (doAllTCPlots)
0634       NumberOfRecHitsPerTrackCandVsEtaProfile->Fill(eta, numberOfHits);
0635     if (doAllTCPlots)
0636       NumberOfRecHitsPerTrackCandVsPhiProfile->Fill(phi, numberOfHits);
0637   }
0638 }
0639 
0640 namespace {
0641   bool trackSelected(unsigned char mask, unsigned char qual) { return mask & 1 << qual; }
0642 }  // namespace
0643 void TrackBuildingAnalyzer::analyze(const edm::View<reco::Track>& trackCollection,
0644                                     const std::vector<const MVACollection*>& mvaCollections,
0645                                     const std::vector<const QualityMaskCollection*>& qualityMaskCollections) {
0646   if (!(doAllTCPlots || doMVAPlots))
0647     return;
0648   if (trackCollection.empty())
0649     return;
0650 
0651   const auto ntracks = trackCollection.size();
0652   const auto nmva = mvaCollections.size();
0653   for (const auto mva : mvaCollections) {
0654     if (mva->size() != ntracks) {
0655       edm::LogError("LogicError") << "TrackBuildingAnalyzer: Incompatible size of MVACollection, " << mva->size()
0656                                   << " differs from the size of the track collection " << ntracks;
0657       return;
0658     }
0659   }
0660   for (const auto qual : qualityMaskCollections) {
0661     if (qual->size() != ntracks) {
0662       edm::LogError("LogicError") << "TrackBuildingAnalyzer: Incompatible size of QualityMaskCollection, "
0663                                   << qual->size() << " differs from the size of the track collection " << ntracks;
0664       return;
0665     }
0666   }
0667 
0668   for (size_t iTrack = 0; iTrack < ntracks; ++iTrack) {
0669     // Fill MVA1 histos with all tracks, MVA2 histos only with tracks
0670     // not selected by MVA1 etc
0671     bool selectedLoose = false;
0672     bool selectedHP = false;
0673 
0674     const auto pt = trackCollection[iTrack].pt();
0675     const auto eta = trackCollection[iTrack].eta();
0676 
0677     for (size_t iMVA = 0; iMVA < nmva; ++iMVA) {
0678       const auto mva = (*(mvaCollections[iMVA]))[iTrack];
0679       if (!selectedLoose) {
0680         trackMVAs[iMVA]->Fill(mva);
0681         trackMVAsVsPtProfile[iMVA]->Fill(pt, mva);
0682         trackMVAsVsEtaProfile[iMVA]->Fill(eta, mva);
0683       }
0684       if (iMVA >= 1 && !selectedHP) {
0685         trackMVAsHP[iMVA]->Fill(mva);
0686         trackMVAsHPVsPtProfile[iMVA]->Fill(pt, mva);
0687         trackMVAsHPVsEtaProfile[iMVA]->Fill(eta, mva);
0688       }
0689 
0690       const auto qual = (*(qualityMaskCollections)[iMVA])[iTrack];
0691       selectedLoose |= trackSelected(qual, reco::TrackBase::loose);
0692       selectedHP |= trackSelected(qual, reco::TrackBase::highPurity);
0693 
0694       if (selectedLoose && selectedHP)
0695         break;
0696     }
0697   }
0698 }
0699 
0700 void TrackBuildingAnalyzer::analyze(const reco::CandidateView& regionCandidates) {
0701   if (!doRegionPlots || !doRegionCandidatePlots)
0702     return;
0703 
0704   for (const auto& candidate : regionCandidates) {
0705     const auto eta = candidate.eta();
0706     const auto phi = candidate.phi();
0707     if (doAllSeedPlots || doPT)
0708       TrackingRegionCandidatePt->Fill(candidate.pt());
0709     if (doAllSeedPlots || doETA)
0710       TrackingRegionCandidateEta->Fill(eta);
0711     if (doAllSeedPlots || doPHI)
0712       TrackingRegionCandidatePhi->Fill(phi);
0713     if (doAllSeedPlots || doPHIVsETA)
0714       TrackingRegionCandidatePhiVsEta->Fill(eta, phi);
0715   }
0716 }
0717 
0718 void TrackBuildingAnalyzer::analyze(const std::vector<std::unique_ptr<TrackingRegion>>& regions) {
0719   analyzeRegions(regions);
0720 }
0721 void TrackBuildingAnalyzer::analyze(const TrackingRegionsSeedingLayerSets& regions) { analyzeRegions(regions); }
0722 
0723 namespace {
0724   const TrackingRegion* regionPtr(const std::unique_ptr<TrackingRegion>& region) { return region.get(); }
0725   const TrackingRegion* regionPtr(const TrackingRegionsSeedingLayerSets::RegionLayers& regionLayers) {
0726     return &(regionLayers.region());
0727   }
0728 }  // namespace
0729 
0730 template <typename T>
0731 void TrackBuildingAnalyzer::analyzeRegions(const T& regions) {
0732   if (!doRegionPlots && etaBinWidth <= 0. && phiBinWidth <= 0.)
0733     return;
0734 
0735   for (const auto& tmp : regions) {
0736     if (const auto* etaPhiRegion = dynamic_cast<const RectangularEtaPhiTrackingRegion*>(regionPtr(tmp))) {
0737       const auto& etaRange = etaPhiRegion->etaRange();
0738       const auto& phiMargin = etaPhiRegion->phiMargin();
0739 
0740       const auto etaMin = etaRange.min();
0741       const auto etaMax = etaRange.max();
0742 
0743       const auto phiMin = etaPhiRegion->phiDirection() - phiMargin.left();
0744       const auto phiMax = etaPhiRegion->phiDirection() + phiMargin.right();
0745 
0746       if (doAllSeedPlots || doETA) {
0747         for (auto eta = etaMin; eta < etaMax; eta += etaBinWidth) {
0748           TrackingRegionEta->Fill(eta);
0749         }
0750       }
0751       if (doAllSeedPlots || doPHI) {
0752         for (auto phi = phiMin; phi < phiMax; phi += phiBinWidth) {
0753           TrackingRegionPhi->Fill(reco::reduceRange(phi));
0754         }
0755       }
0756       if (doAllSeedPlots || doPHIVsETA) {
0757         for (auto phi = phiMin; phi < phiMax; phi += phiBinWidth) {
0758           const auto reducedPhi = reco::reduceRange(phi);
0759           for (auto eta = etaMin; eta < etaMax; eta += etaBinWidth) {
0760             TrackingRegionPhiVsEta->Fill(eta, reducedPhi);
0761           }
0762         }
0763       }
0764     }
0765   }
0766 }