Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*
0002  *  See header file for a description of this class.
0003  **
0004  *  \author Suchandra Dutta , Giorgia Mila
0005  */
0006 
0007 #include "DataFormats/TrackReco/interface/Track.h"
0008 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0009 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0010 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0011 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0012 
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/Utilities/interface/InputTag.h"
0015 #include "FWCore/Utilities/interface/transform.h"
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 #include "DQM/TrackingMonitor/interface/TrackBuildingAnalyzer.h"
0018 #include "DQM/TrackingMonitor/interface/TrackAnalyzer.h"
0019 #include "DQM/TrackingMonitor/interface/VertexMonitor.h"
0020 #include "DQM/TrackingMonitor/interface/TrackingMonitor.h"
0021 
0022 #include "FWCore/ParameterSet/interface/Registry.h"
0023 
0024 #include "TrackingTools/TrajectoryState/interface/TrajectoryStateTransform.h"
0025 #include "TrackingTools/PatternTools/interface/TSCBLBuilderNoMaterial.h"
0026 #include "TrackingTools/PatternTools/interface/TSCPBuilderNoMaterial.h"
0027 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0028 #include "CommonTools/TriggerUtils/interface/GenericTriggerEventFlag.h"
0029 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0030 
0031 #include "DataFormats/Candidate/interface/Candidate.h"
0032 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0033 
0034 #include "DQM/TrackingMonitor/interface/GetLumi.h"
0035 
0036 #include <string>
0037 
0038 #ifdef VI_DEBUG
0039 #define COUT(x) std::cout << x << ' '
0040 #else
0041 #define COUT(x) LogDebug(x)
0042 #endif
0043 
0044 // TrackingMonitor
0045 // ----------------------------------------------------------------------------------//
0046 
0047 TrackingMonitor::TrackingMonitor(const edm::ParameterSet& iConfig)
0048     : confID_(iConfig.id()),
0049       theTrackBuildingAnalyzer(new TrackBuildingAnalyzer(iConfig)),
0050       NumberOfTracks(nullptr),
0051       NumberOfMeanRecHitsPerTrack(nullptr),
0052       NumberOfMeanLayersPerTrack(nullptr)
0053       //    , NumberOfGoodTracks(NULL)
0054       ,
0055       FractionOfGoodTracks(nullptr),
0056       NumberOfTrackingRegions(nullptr),
0057       NumberOfSeeds(nullptr),
0058       NumberOfSeeds_lumiFlag(nullptr),
0059       NumberOfTrackCandidates(nullptr),
0060       FractionCandidatesOverSeeds(nullptr)
0061       //    , NumberOfGoodTrkVsClus(NULL)
0062       ,
0063       NumberEventsOfVsLS(nullptr),
0064       NumberOfTracksVsLS(nullptr)
0065       //    , NumberOfGoodTracksVsLS(NULL)
0066       ,
0067       GoodTracksFractionVsLS(nullptr)
0068       //    , GoodTracksNumberOfRecHitsPerTrackVsLS(NULL)
0069       // ADD by Mia for PU monitoring
0070       // vertex plots to be moved in ad hoc class
0071       ,
0072       NumberOfGoodPVtxVsLS(nullptr),
0073       NumberOfGoodPVtxWO0VsLS(nullptr),
0074       NumberEventsOfVsBX(nullptr),
0075       NumberOfTracksVsBX(nullptr),
0076       GoodTracksFractionVsBX(nullptr),
0077       NumberOfRecHitsPerTrackVsBX(nullptr),
0078       NumberOfGoodPVtxVsBX(nullptr),
0079       NumberOfGoodPVtxWO0VsBX(nullptr),
0080       NumberOfTracksVsBXlumi(nullptr),
0081       NumberOfTracksVsGoodPVtx(nullptr),
0082       NumberOfTracksVsPUPVtx(nullptr),
0083       NumberEventsOfVsGoodPVtx(nullptr),
0084       GoodTracksFractionVsGoodPVtx(nullptr),
0085       NumberOfRecHitsPerTrackVsGoodPVtx(nullptr),
0086       NumberOfPVtxVsGoodPVtx(nullptr),
0087       NumberOfPixelClustersVsGoodPVtx(nullptr),
0088       NumberOfStripClustersVsGoodPVtx(nullptr),
0089       NumberEventsOfVsLUMI(nullptr),
0090       NumberOfTracksVsLUMI(nullptr),
0091       GoodTracksFractionVsLUMI(nullptr),
0092       NumberOfRecHitsPerTrackVsLUMI(nullptr),
0093       NumberOfGoodPVtxVsLUMI(nullptr),
0094       NumberOfGoodPVtxWO0VsLUMI(nullptr),
0095       NumberOfPixelClustersVsLUMI(nullptr),
0096       NumberOfStripClustersVsLUMI(nullptr),
0097       NumberOfTracks_lumiFlag(nullptr)
0098       //    , NumberOfGoodTracks_lumiFlag(NULL)
0099 
0100       ,
0101       builderName(iConfig.getParameter<std::string>("TTRHBuilder")),
0102       doTrackerSpecific_(iConfig.getParameter<bool>("doTrackerSpecific")),
0103       doLumiAnalysis(iConfig.getParameter<bool>("doLumiAnalysis")),
0104       doProfilesVsLS_(iConfig.getParameter<bool>("doProfilesVsLS")),
0105       doAllSeedPlots(iConfig.getParameter<bool>("doSeedParameterHistos")),
0106       doAllPlots(iConfig.getParameter<bool>("doAllPlots")),
0107       doGeneralPropertiesPlots_(iConfig.getParameter<bool>("doGeneralPropertiesPlots")),
0108       doHitPropertiesPlots_(iConfig.getParameter<bool>("doHitPropertiesPlots")),
0109       doTkCandPlots(iConfig.getParameter<bool>("doTrackCandHistos")),
0110       doPUmonitoring_(iConfig.getParameter<bool>("doPUmonitoring")),
0111       genTriggerEventFlag_(new GenericTriggerEventFlag(
0112           iConfig.getParameter<edm::ParameterSet>("genericTriggerEventPSet"), consumesCollector(), *this)),
0113       numSelection_(iConfig.getParameter<std::string>("numCut")),
0114       denSelection_(iConfig.getParameter<std::string>("denCut")),
0115       pvNDOF_(iConfig.getParameter<int>("pvNDOF")),
0116       forceSCAL_(iConfig.getParameter<bool>("forceSCAL")) {
0117   edm::ConsumesCollector c{consumesCollector()};
0118   theTrackAnalyzer = new tadqm::TrackAnalyzer(iConfig, c);
0119 
0120   // input tags for collections from the configuration
0121   bsSrc_ = iConfig.getParameter<edm::InputTag>("beamSpot");
0122   pvSrc_ = iConfig.getParameter<edm::InputTag>("primaryVertex");
0123   bsSrcToken_ = consumes<reco::BeamSpot>(bsSrc_);
0124   pvSrcToken_ = mayConsume<reco::VertexCollection>(pvSrc_);
0125 
0126   lumiscalersToken_ = consumes<LumiScalersCollection>(iConfig.getParameter<edm::InputTag>("scal"));
0127   metaDataToken_ = consumes<OnlineLuminosityRecord>(iConfig.getParameter<edm::InputTag>("metadata"));
0128 
0129   edm::InputTag alltrackProducer = iConfig.getParameter<edm::InputTag>("allTrackProducer");
0130   edm::InputTag trackProducer = iConfig.getParameter<edm::InputTag>("TrackProducer");
0131   edm::InputTag tcProducer = iConfig.getParameter<edm::InputTag>("TCProducer");
0132   edm::InputTag seedProducer = iConfig.getParameter<edm::InputTag>("SeedProducer");
0133   allTrackToken_ = consumes<edm::View<reco::Track> >(alltrackProducer);
0134   trackToken_ = consumes<edm::View<reco::Track> >(trackProducer);
0135   trackCandidateToken_ = consumes<TrackCandidateCollection>(tcProducer);
0136   seedToken_ = consumes<edm::View<TrajectorySeed> >(seedProducer);
0137   seedStopInfoToken_ = consumes<std::vector<SeedStopInfo> >(tcProducer);
0138 
0139   doMVAPlots = iConfig.getParameter<bool>("doMVAPlots");
0140   if (doMVAPlots) {
0141     mvaQualityTokens_ = edm::vector_transform(
0142         iConfig.getParameter<std::vector<std::string> >("MVAProducers"), [&](const std::string& tag) {
0143           return std::make_tuple(consumes<MVACollection>(edm::InputTag(tag, "MVAValues")),
0144                                  consumes<QualityMaskCollection>(edm::InputTag(tag, "QualityMasks")));
0145         });
0146     mvaTrackToken_ = consumes<edm::View<reco::Track> >(iConfig.getParameter<edm::InputTag>("TrackProducerForMVA"));
0147   }
0148 
0149   doRegionPlots = iConfig.getParameter<bool>("doRegionPlots");
0150   doRegionCandidatePlots = iConfig.getParameter<bool>("doRegionCandidatePlots");
0151   if (doRegionPlots) {
0152     const auto& regionTag = iConfig.getParameter<edm::InputTag>("RegionProducer");
0153     if (!regionTag.label().empty()) {
0154       regionToken_ = consumes(regionTag);
0155     }
0156     const auto& regionLayersTag = iConfig.getParameter<edm::InputTag>("RegionSeedingLayersProducer");
0157     if (!regionLayersTag.label().empty()) {
0158       if (!regionToken_.isUninitialized()) {
0159         throw cms::Exception("Configuration")
0160             << "Only one of 'RegionProducer' and 'RegionSeedingLayersProducer' can be non-empty, now both are.";
0161       }
0162       regionLayerSetsToken_ = consumes<TrackingRegionsSeedingLayerSets>(regionLayersTag);
0163     } else if (regionToken_.isUninitialized()) {
0164       throw cms::Exception("Configuration") << "With doRegionPlots=True either 'RegionProducer' or "
0165                                                "'RegionSeedingLayersProducer' must be non-empty, now both are empty.";
0166     }
0167 
0168     if (doRegionCandidatePlots) {
0169       regionCandidateToken_ = consumes<reco::CandidateView>(iConfig.getParameter<edm::InputTag>("RegionCandidates"));
0170     }
0171   }
0172 
0173   edm::InputTag stripClusterInputTag_ = iConfig.getParameter<edm::InputTag>("stripCluster");
0174   edm::InputTag pixelClusterInputTag_ = iConfig.getParameter<edm::InputTag>("pixelCluster");
0175   stripClustersToken_ = mayConsume<edmNew::DetSetVector<SiStripCluster> >(stripClusterInputTag_);
0176   pixelClustersToken_ = mayConsume<edmNew::DetSetVector<SiPixelCluster> >(pixelClusterInputTag_);
0177 
0178   runTrackBuildingAnalyzerForSeed =
0179       (doAllSeedPlots || iConfig.getParameter<bool>("doSeedPTHisto") || iConfig.getParameter<bool>("doSeedETAHisto") ||
0180        iConfig.getParameter<bool>("doSeedPHIHisto") || iConfig.getParameter<bool>("doSeedPHIVsETAHisto") ||
0181        iConfig.getParameter<bool>("doSeedThetaHisto") || iConfig.getParameter<bool>("doSeedQHisto") ||
0182        iConfig.getParameter<bool>("doSeedDxyHisto") || iConfig.getParameter<bool>("doSeedDzHisto") ||
0183        iConfig.getParameter<bool>("doSeedNRecHitsHisto") || iConfig.getParameter<bool>("doSeedNVsPhiProf") ||
0184        iConfig.getParameter<bool>("doSeedNVsEtaProf"));
0185 
0186   if (doTkCandPlots || doAllSeedPlots || runTrackBuildingAnalyzerForSeed) {
0187     magneticFieldToken_ = esConsumes<MagneticField, IdealMagneticFieldRecord>();
0188     transientTrackingRecHitBuilderToken_ =
0189         esConsumes<TransientTrackingRecHitBuilder, TransientRecHitRecord>(edm::ESInputTag("", builderName));
0190   }
0191 
0192   doFractionPlot_ = true;
0193   if (alltrackProducer.label() == trackProducer.label())
0194     doFractionPlot_ = false;
0195 
0196   Quality_ = iConfig.getParameter<std::string>("Quality");
0197   AlgoName_ = iConfig.getParameter<std::string>("AlgoName");
0198 
0199   // get flag from the configuration
0200   doPlotsVsBXlumi_ = iConfig.getParameter<bool>("doPlotsVsBXlumi");
0201   if (doPlotsVsBXlumi_)
0202     theLumiDetails_ = new GetLumi(iConfig.getParameter<edm::ParameterSet>("BXlumiSetup"), c);
0203   doPlotsVsGoodPVtx_ = iConfig.getParameter<bool>("doPlotsVsGoodPVtx");
0204   doPlotsVsLUMI_ = iConfig.getParameter<bool>("doPlotsVsLUMI");
0205   doPlotsVsBX_ = iConfig.getParameter<bool>("doPlotsVsBX");
0206 
0207   if (doPUmonitoring_) {
0208     std::vector<edm::InputTag> primaryVertexInputTags =
0209         iConfig.getParameter<std::vector<edm::InputTag> >("primaryVertexInputTags");
0210     std::vector<edm::InputTag> selPrimaryVertexInputTags =
0211         iConfig.getParameter<std::vector<edm::InputTag> >("selPrimaryVertexInputTags");
0212     std::vector<std::string> pvLabels = iConfig.getParameter<std::vector<std::string> >("pvLabels");
0213 
0214     if (primaryVertexInputTags.size() == pvLabels.size() and
0215         primaryVertexInputTags.size() == selPrimaryVertexInputTags.size()) {
0216       for (size_t i = 0; i < primaryVertexInputTags.size(); i++) {
0217         edm::InputTag iPVinputTag = primaryVertexInputTags[i];
0218         edm::InputTag iSelPVinputTag = selPrimaryVertexInputTags[i];
0219         std::string iPVlabel = pvLabels[i];
0220 
0221         theVertexMonitor.push_back(new VertexMonitor(iConfig, iPVinputTag, iSelPVinputTag, iPVlabel, c));
0222       }
0223     }
0224   }
0225 }
0226 
0227 TrackingMonitor::~TrackingMonitor() {
0228   if (theTrackAnalyzer)
0229     delete theTrackAnalyzer;
0230   if (theTrackBuildingAnalyzer)
0231     delete theTrackBuildingAnalyzer;
0232   if (doPUmonitoring_)
0233     for (size_t i = 0; i < theVertexMonitor.size(); i++)
0234       if (theVertexMonitor[i])
0235         delete theVertexMonitor[i];
0236   if (genTriggerEventFlag_)
0237     delete genTriggerEventFlag_;
0238 }
0239 
0240 void TrackingMonitor::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& iRun, edm::EventSetup const& iSetup) {
0241   // parameters from the configuration
0242   auto const* conf = edm::pset::Registry::instance()->getMapped(confID_);
0243   assert(conf != nullptr);
0244   std::string Quality = conf->getParameter<std::string>("Quality");
0245   std::string AlgoName = conf->getParameter<std::string>("AlgoName");
0246   MEFolderName = conf->getParameter<std::string>("FolderName");
0247   std::string Folder = MEFolderName.substr(0, 2);
0248 
0249   // test for the Quality veriable validity
0250   if (!Quality_.empty()) {
0251     if (Quality_ != "highPurity" && Quality_ != "tight" && Quality_ != "loose") {
0252       edm::LogWarning("TrackingMonitor") << "Qualty Name is invalid, using no quality criterea by default";
0253       Quality_ = "";
0254     }
0255   }
0256 
0257   // use the AlgoName and Quality Name
0258   std::string CategoryName = !Quality_.empty() ? AlgoName_ + "_" + Quality_ : AlgoName_;
0259 
0260   // get binning from the configuration
0261   int TKNoBin = conf->getParameter<int>("TkSizeBin");
0262   double TKNoMin = conf->getParameter<double>("TkSizeMin");
0263   double TKNoMax = conf->getParameter<double>("TkSizeMax");
0264 
0265   int TCNoBin = conf->getParameter<int>("TCSizeBin");
0266   double TCNoMin = conf->getParameter<double>("TCSizeMin");
0267   double TCNoMax = conf->getParameter<double>("TCSizeMax");
0268 
0269   int TKNoSeedBin = conf->getParameter<int>("TkSeedSizeBin");
0270   double TKNoSeedMin = conf->getParameter<double>("TkSeedSizeMin");
0271   double TKNoSeedMax = conf->getParameter<double>("TkSeedSizeMax");
0272 
0273   //int RecHitBin = conf->getParameter<int>("RecHitBin");
0274   double RecHitMin = conf->getParameter<double>("RecHitMin");
0275   double RecHitMax = conf->getParameter<double>("RecHitMax");
0276 
0277   int MeanHitBin = conf->getParameter<int>("MeanHitBin");
0278   double MeanHitMin = conf->getParameter<double>("MeanHitMin");
0279   double MeanHitMax = conf->getParameter<double>("MeanHitMax");
0280 
0281   int MeanLayBin = conf->getParameter<int>("MeanLayBin");
0282   double MeanLayMin = conf->getParameter<double>("MeanLayMin");
0283   double MeanLayMax = conf->getParameter<double>("MeanLayMax");
0284 
0285   int LSBin = conf->getParameter<int>("LSBin");
0286   int LSMin = conf->getParameter<double>("LSMin");
0287   int LSMax = conf->getParameter<double>("LSMax");
0288 
0289   std::string StateName = conf->getParameter<std::string>("MeasurementState");
0290   if (StateName != "OuterSurface" && StateName != "InnerSurface" && StateName != "ImpactPoint" &&
0291       StateName != "default" && StateName != "All") {
0292     // print warning
0293     edm::LogWarning("TrackingMonitor") << "State Name is invalid, using 'ImpactPoint' by default";
0294   }
0295 
0296   ibooker.setCurrentFolder(MEFolderName);
0297 
0298   // book the General Property histograms
0299   // ---------------------------------------------------------------------------------//
0300 
0301   if (doGeneralPropertiesPlots_ || doAllPlots) {
0302     ibooker.setCurrentFolder(MEFolderName + "/GeneralProperties");
0303 
0304     histname = "NumberOfTracks_" + CategoryName;
0305     // MODIFY by Mia in order to cope w/ high multiplicity
0306     NumberOfTracks = ibooker.book1D(histname, histname, TKNoBin, TKNoMin, TKNoMax);
0307     NumberOfTracks->setAxisTitle("Number of Tracks per Event", 1);
0308     NumberOfTracks->setAxisTitle("Number of Events", 2);
0309 
0310     if (Folder == "Tr") {
0311       histname = "NumberOfTracks_PUvtx_" + CategoryName;
0312       NumberOfTracks_PUvtx = ibooker.book1D(histname, histname, TKNoBin, TKNoMin, TKNoMax);
0313       NumberOfTracks_PUvtx->setAxisTitle("Number of Tracks per Event (matched a PU vertex)", 1);
0314       NumberOfTracks_PUvtx->setAxisTitle("Number of Events", 2);
0315 
0316       histname = "NumberofTracks_Hardvtx_" + CategoryName;
0317       NumberofTracks_Hardvtx =
0318           ibooker.book1D(histname, histname, TKNoBin / 10, TKNoMin, ((TKNoMax - TKNoMin) / 10) - 0.5);
0319       NumberofTracks_Hardvtx->setAxisTitle("Number of Tracks per Event (matched main vertex)", 1);
0320       NumberofTracks_Hardvtx->setAxisTitle("Number of Events", 2);
0321 
0322       histname = "NumberofTracks_Hardvtx_PUvtx_" + CategoryName;
0323       NumberofTracks_Hardvtx_PUvtx = ibooker.book1D(histname, histname, 2, 0., 2.);
0324       NumberofTracks_Hardvtx_PUvtx->setAxisTitle("Number of Tracks per PU/Hard vertex", 1);
0325       NumberofTracks_Hardvtx_PUvtx->setAxisTitle("Number of Tracks", 2);
0326       NumberofTracks_Hardvtx_PUvtx->setBinLabel(1, "PU_Vertex");
0327       NumberofTracks_Hardvtx_PUvtx->setBinLabel(2, "Hard_Vertex");
0328     }
0329 
0330     histname = "NumberOfMeanRecHitsPerTrack_" + CategoryName;
0331     NumberOfMeanRecHitsPerTrack = ibooker.book1D(histname, histname, MeanHitBin, MeanHitMin, MeanHitMax);
0332     NumberOfMeanRecHitsPerTrack->setAxisTitle("Mean number of valid RecHits per Track", 1);
0333     NumberOfMeanRecHitsPerTrack->setAxisTitle("Entries", 2);
0334 
0335     histname = "NumberOfMeanLayersPerTrack_" + CategoryName;
0336     NumberOfMeanLayersPerTrack = ibooker.book1D(histname, histname, MeanLayBin, MeanLayMin, MeanLayMax);
0337     NumberOfMeanLayersPerTrack->setAxisTitle("Mean number of Layers per Track", 1);
0338     NumberOfMeanLayersPerTrack->setAxisTitle("Entries", 2);
0339 
0340     if (doFractionPlot_) {
0341       histname = "FractionOfGoodTracks_" + CategoryName;
0342       FractionOfGoodTracks = ibooker.book1D(histname, histname, 101, -0.005, 1.005);
0343       FractionOfGoodTracks->setAxisTitle("Fraction of Tracks (w.r.t. generalTracks)", 1);
0344       FractionOfGoodTracks->setAxisTitle("Entries", 2);
0345     }
0346   }
0347 
0348   if (doLumiAnalysis) {
0349     // add by Mia in order to deal with LS transitions
0350     ibooker.setCurrentFolder(MEFolderName + "/LSanalysis");
0351     auto scope = DQMStore::IBooker::UseLumiScope(ibooker);
0352 
0353     histname = "NumberOfTracks_lumiFlag_" + CategoryName;
0354     NumberOfTracks_lumiFlag = ibooker.book1D(histname, histname, TKNoBin, TKNoMin, TKNoMax);
0355     NumberOfTracks_lumiFlag->setAxisTitle("Number of Tracks per Event", 1);
0356     NumberOfTracks_lumiFlag->setAxisTitle("Number of Events", 2);
0357   }
0358 
0359   // book profile plots vs LS :
0360   //---------------------------
0361 
0362   if (doProfilesVsLS_ || doAllPlots) {
0363     ibooker.setCurrentFolder(MEFolderName + "/GeneralProperties");
0364 
0365     histname = "NumberOfTracksVsLS_" + CategoryName;
0366     NumberOfTracksVsLS = ibooker.bookProfile(histname, histname, LSBin, LSMin, LSMax, TKNoMin, TKNoMax, "");
0367     NumberOfTracksVsLS->getTH1()->SetCanExtend(TH1::kAllAxes);
0368     NumberOfTracksVsLS->setAxisTitle("#Lumi section", 1);
0369     NumberOfTracksVsLS->setAxisTitle("Number of  Tracks", 2);
0370 
0371     histname = "NumberOfRecHitsPerTrackVsLS_" + CategoryName;
0372     NumberOfRecHitsPerTrackVsLS =
0373         ibooker.bookProfile(histname, histname, LSBin, LSMin, LSMax, RecHitMin, RecHitMax * 5, "");
0374     NumberOfRecHitsPerTrackVsLS->getTH1()->SetCanExtend(TH1::kAllAxes);
0375     NumberOfRecHitsPerTrackVsLS->setAxisTitle("#Lumi section", 1);
0376     NumberOfRecHitsPerTrackVsLS->setAxisTitle("Mean number of Valid RecHits per track", 2);
0377 
0378     histname = "NumberEventsVsLS_" + CategoryName;
0379     NumberEventsOfVsLS = ibooker.book1D(histname, histname, LSBin, LSMin, LSMax);
0380     NumberEventsOfVsLS->getTH1()->SetCanExtend(TH1::kAllAxes);
0381     NumberEventsOfVsLS->setAxisTitle("#Lumi section", 1);
0382     NumberEventsOfVsLS->setAxisTitle("Number of events", 2);
0383 
0384     edm::ParameterSet ParametersGoodPVtx = conf->getParameter<edm::ParameterSet>("GoodPVtx");
0385 
0386     double GoodPVtxMin = ParametersGoodPVtx.getParameter<double>("GoodPVtxMin");
0387     double GoodPVtxMax = ParametersGoodPVtx.getParameter<double>("GoodPVtxMax");
0388 
0389     histname = "NumberOfGoodPVtxVsLS_" + CategoryName;
0390     NumberOfGoodPVtxVsLS =
0391         ibooker.bookProfile(histname, histname, LSBin, LSMin, LSMax, GoodPVtxMin, 3. * GoodPVtxMax, "");
0392     NumberOfGoodPVtxVsLS->getTH1()->SetCanExtend(TH1::kAllAxes);
0393     NumberOfGoodPVtxVsLS->setAxisTitle("#Lumi section", 1);
0394     NumberOfGoodPVtxVsLS->setAxisTitle("Mean number of good PV", 2);
0395 
0396     histname = "NumberOfGoodPVtxWO0VsLS_" + CategoryName;
0397     NumberOfGoodPVtxWO0VsLS =
0398         ibooker.bookProfile(histname, histname, LSBin, LSMin, LSMax, GoodPVtxMin, 3. * GoodPVtxMax, "");
0399     NumberOfGoodPVtxWO0VsLS->setAxisTitle("#Lumi section", 1);
0400     NumberOfGoodPVtxWO0VsLS->setAxisTitle("Mean number of good PV", 2);
0401 
0402     if (doFractionPlot_) {
0403       histname = "GoodTracksFractionVsLS_" + CategoryName;
0404       GoodTracksFractionVsLS = ibooker.bookProfile(histname, histname, LSBin, LSMin, LSMax, 0, 1.1, "");
0405       GoodTracksFractionVsLS->getTH1()->SetCanExtend(TH1::kAllAxes);
0406       GoodTracksFractionVsLS->setAxisTitle("#Lumi section", 1);
0407       GoodTracksFractionVsLS->setAxisTitle("Fraction of Good Tracks", 2);
0408     }
0409 
0410     if (doPlotsVsBX_ || doAllPlots) {
0411       ibooker.setCurrentFolder(MEFolderName + "/BXanalysis");
0412       int BXBin = 3564;
0413       double BXMin = 0.5;
0414       double BXMax = 3564.5;
0415 
0416       histname = "NumberEventsVsBX_" + CategoryName;
0417       NumberEventsOfVsBX = ibooker.book1D(histname, histname, BXBin, BXMin, BXMax);
0418       NumberEventsOfVsBX->setAxisTitle("BX", 1);
0419       NumberEventsOfVsBX->setAxisTitle("Number of events", 2);
0420 
0421       histname = "NumberOfTracksVsBX_" + CategoryName;
0422       NumberOfTracksVsBX = ibooker.bookProfile(histname, histname, BXBin, BXMin, BXMax, TKNoMin, TKNoMax * 3., "");
0423       NumberOfTracksVsBX->setAxisTitle("BX", 1);
0424       NumberOfTracksVsBX->setAxisTitle("Number of  Tracks", 2);
0425 
0426       histname = "NumberOfRecHitsPerTrackVsBX_" + CategoryName;
0427       NumberOfRecHitsPerTrackVsBX =
0428           ibooker.bookProfile(histname, histname, BXBin, BXMin, BXMax, RecHitMin, RecHitMax * 5, "");
0429       NumberOfRecHitsPerTrackVsBX->setAxisTitle("BX", 1);
0430       NumberOfRecHitsPerTrackVsBX->setAxisTitle("Mean number of Valid RecHits per track", 2);
0431 
0432       histname = "NumberOfGoodPVtxVsBX_" + CategoryName;
0433       NumberOfGoodPVtxVsBX =
0434           ibooker.bookProfile(histname, histname, BXBin, BXMin, BXMax, GoodPVtxMin, 3. * GoodPVtxMax, "");
0435       NumberOfGoodPVtxVsBX->setAxisTitle("BX", 1);
0436       NumberOfGoodPVtxVsBX->setAxisTitle("Mean number of good PV", 2);
0437 
0438       histname = "NumberOfGoodPVtxWO0VsBX_" + CategoryName;
0439       NumberOfGoodPVtxWO0VsBX =
0440           ibooker.bookProfile(histname, histname, BXBin, BXMin, BXMax, GoodPVtxMin, 3. * GoodPVtxMax, "");
0441       NumberOfGoodPVtxWO0VsBX->setAxisTitle("BX", 1);
0442       NumberOfGoodPVtxWO0VsBX->setAxisTitle("Mean number of good PV", 2);
0443 
0444       if (doFractionPlot_) {
0445         histname = "GoodTracksFractionVsBX_" + CategoryName;
0446         GoodTracksFractionVsBX = ibooker.bookProfile(histname, histname, BXBin, BXMin, BXMax, 0, 1.1, "");
0447         GoodTracksFractionVsBX->setAxisTitle("BX", 1);
0448         GoodTracksFractionVsBX->setAxisTitle("Fraction of Good Tracks", 2);
0449       }
0450     }
0451   }
0452 
0453   // book PU monitoring plots :
0454   //---------------------------
0455 
0456   if (doPUmonitoring_) {
0457     for (size_t i = 0; i < theVertexMonitor.size(); i++)
0458       theVertexMonitor[i]->initHisto(ibooker);
0459   }
0460 
0461   if (doPlotsVsGoodPVtx_) {
0462     ibooker.setCurrentFolder(MEFolderName + "/PUmonitoring");
0463     // get binning from the configuration
0464     int PVBin = conf->getParameter<int>("PVBin");
0465     float PVMin = conf->getParameter<double>("PVMin");
0466     float PVMax = conf->getParameter<double>("PVMax");
0467 
0468     histname = "NumberOfTracksVsGoodPVtx";
0469     NumberOfTracksVsGoodPVtx = ibooker.bookProfile(histname, histname, PVBin, PVMin, PVMax, TKNoMin, TKNoMax * 5., "");
0470     NumberOfTracksVsGoodPVtx->setAxisTitle("Number of PV", 1);
0471     NumberOfTracksVsGoodPVtx->setAxisTitle("Mean number of Tracks per Event", 2);
0472 
0473     histname = "NumberOfTracksVsPUPVtx";
0474     NumberOfTracksVsPUPVtx = ibooker.bookProfile(histname, histname, PVBin, PVMin, PVMax, 0., TKNoMax * 5., "");
0475     NumberOfTracksVsPUPVtx->setAxisTitle("Number of PU", 1);
0476     NumberOfTracksVsPUPVtx->setAxisTitle("Mean number of Tracks per PUvtx", 2);
0477 
0478     histname = "NumberEventsVsGoodPVtx";
0479     NumberEventsOfVsGoodPVtx = ibooker.book1D(histname, histname, PVBin, PVMin, PVMax);
0480     NumberEventsOfVsGoodPVtx->setAxisTitle("Number of good PV (PU)", 1);
0481     NumberEventsOfVsGoodPVtx->setAxisTitle("Number of events", 2);
0482 
0483     if (doFractionPlot_) {
0484       histname = "GoodTracksFractionVsGoodPVtx";
0485       GoodTracksFractionVsGoodPVtx = ibooker.bookProfile(histname, histname, PVBin, PVMin, PVMax, 0., 1.1, "");
0486       GoodTracksFractionVsGoodPVtx->setAxisTitle("Number of good PV (PU)", 1);
0487       GoodTracksFractionVsGoodPVtx->setAxisTitle("Mean fraction of good tracks", 2);
0488     }
0489 
0490     histname = "NumberOfRecHitsPerTrackVsGoodPVtx";
0491     NumberOfRecHitsPerTrackVsGoodPVtx = ibooker.bookProfile(histname, histname, PVBin, PVMin, PVMax, 0., 200., "");
0492     NumberOfRecHitsPerTrackVsGoodPVtx->setAxisTitle("Number of good PV (PU)", 1);
0493     NumberOfRecHitsPerTrackVsGoodPVtx->setAxisTitle("Mean number of valid rechits per Tracks", 2);
0494 
0495     histname = "NumberOfPVtxVsGoodPVtx";
0496     NumberOfPVtxVsGoodPVtx = ibooker.bookProfile(histname, histname, PVBin, PVMin, PVMax, 0., 3. * PVMax, "");
0497     NumberOfPVtxVsGoodPVtx->setAxisTitle("Number of good PV (PU)", 1);
0498     NumberOfPVtxVsGoodPVtx->setAxisTitle("Mean number of vertices", 2);
0499 
0500     double NClusPxMin = conf->getParameter<double>("NClusPxMin");
0501     double NClusPxMax = conf->getParameter<double>("NClusPxMax");
0502     histname = "NumberOfPixelClustersVsGoodPVtx";
0503     NumberOfPixelClustersVsGoodPVtx =
0504         ibooker.bookProfile(histname, histname, PVBin, PVMin, PVMax, NClusPxMin, 3. * NClusPxMax, "");
0505     NumberOfPixelClustersVsGoodPVtx->setAxisTitle("Number of good PV (PU)", 1);
0506     NumberOfPixelClustersVsGoodPVtx->setAxisTitle("Mean number of pixel clusters", 2);
0507 
0508     double NClusStrMin = conf->getParameter<double>("NClusStrMin");
0509     double NClusStrMax = conf->getParameter<double>("NClusStrMax");
0510     histname = "NumberOfStripClustersVsGoodPVtx";
0511     NumberOfStripClustersVsGoodPVtx =
0512         ibooker.bookProfile(histname, histname, PVBin, PVMin, PVMax, NClusStrMin, 3. * NClusStrMax, "");
0513     NumberOfStripClustersVsGoodPVtx->setAxisTitle("Number of good PV (PU)", 1);
0514     NumberOfStripClustersVsGoodPVtx->setAxisTitle("Mean number of strip clusters", 2);
0515   }
0516 
0517   bool isMC = ((iRun.runAuxiliary().run() == 1) && (iRun.runAuxiliary().beginTime().value() == 1));
0518 
0519   if ((doPlotsVsLUMI_ || doAllPlots) && !isMC) {
0520     ibooker.setCurrentFolder(MEFolderName + "/LUMIanalysis");
0521     int LUMIBin = conf->getParameter<int>("LUMIBin");
0522     float LUMIMin = conf->getParameter<double>("LUMIMin");
0523     float LUMIMax = conf->getParameter<double>("LUMIMax");
0524 
0525     histname = "NumberEventsVsLUMI";
0526     NumberEventsOfVsLUMI = ibooker.book1D(histname, histname, LUMIBin, LUMIMin, LUMIMax);
0527     NumberEventsOfVsLUMI->setAxisTitle("online lumi [1e30 Hz cm^{-2}]", 1);
0528     NumberEventsOfVsLUMI->setAxisTitle("Number of events", 2);
0529 
0530     histname = "NumberOfTracksVsLUMI";
0531     NumberOfTracksVsLUMI = ibooker.bookProfile(histname, histname, LUMIBin, LUMIMin, LUMIMax, TKNoMin, TKNoMax * 3, "");
0532     NumberOfTracksVsLUMI->setAxisTitle("online lumi [1e30 Hz cm^{-2}]", 1);
0533     NumberOfTracksVsLUMI->setAxisTitle("Mean number of vertices", 2);
0534 
0535     if (doFractionPlot_) {
0536       histname = "GoodTracksFractionVsLUMI";
0537       GoodTracksFractionVsLUMI = ibooker.bookProfile(histname, histname, LUMIBin, LUMIMin, LUMIMax, 0., 1.1, "");
0538       GoodTracksFractionVsLUMI->setAxisTitle("online lumi [1e30 Hz cm^{-2}]", 1);
0539       GoodTracksFractionVsLUMI->setAxisTitle("Mean number of vertices", 2);
0540     }
0541 
0542     histname = "NumberOfRecHitsPerTrackVsLUMI";
0543     NumberOfRecHitsPerTrackVsLUMI =
0544         ibooker.bookProfile(histname, histname, LUMIBin, LUMIMin, LUMIMax, RecHitMin, RecHitMax * 5, "");
0545     NumberOfRecHitsPerTrackVsLUMI->setAxisTitle("online lumi [1e30 Hz cm^{-2}]", 1);
0546     NumberOfRecHitsPerTrackVsLUMI->setAxisTitle("Mean number of vertices", 2);
0547 
0548     double PVMin = conf->getParameter<double>("PVMin");
0549     double PVMax = conf->getParameter<double>("PVMax");
0550 
0551     histname = "NumberOfGoodPVtxVsLUMI";
0552     NumberOfGoodPVtxVsLUMI = ibooker.bookProfile(histname, histname, LUMIBin, LUMIMin, LUMIMax, PVMin, 3. * PVMax, "");
0553     NumberOfGoodPVtxVsLUMI->setAxisTitle("online lumi [1e30 Hz cm^{-2}]", 1);
0554     NumberOfGoodPVtxVsLUMI->setAxisTitle("Mean number of vertices", 2);
0555 
0556     histname = "NumberOfGoodPVtxWO0VsLUMI";
0557     NumberOfGoodPVtxWO0VsLUMI =
0558         ibooker.bookProfile(histname, histname, LUMIBin, LUMIMin, LUMIMax, PVMin, 3. * PVMax, "");
0559     NumberOfGoodPVtxWO0VsLUMI->setAxisTitle("online lumi [1e30 Hz cm^{-2}]", 1);
0560     NumberOfGoodPVtxWO0VsLUMI->setAxisTitle("Mean number of vertices", 2);
0561 
0562     double NClusPxMin = conf->getParameter<double>("NClusPxMin");
0563     double NClusPxMax = conf->getParameter<double>("NClusPxMax");
0564     histname = "NumberOfPixelClustersVsGoodPVtx";
0565     NumberOfPixelClustersVsLUMI =
0566         ibooker.bookProfile(histname, histname, LUMIBin, LUMIMin, LUMIMax, NClusPxMin, 3. * NClusPxMax, "");
0567     NumberOfPixelClustersVsLUMI->setAxisTitle("online lumi [1e30 Hz cm^{-2}]", 1);
0568     NumberOfPixelClustersVsLUMI->setAxisTitle("Mean number of pixel clusters", 2);
0569 
0570     double NClusStrMin = conf->getParameter<double>("NClusStrMin");
0571     double NClusStrMax = conf->getParameter<double>("NClusStrMax");
0572     histname = "NumberOfStripClustersVsLUMI";
0573     NumberOfStripClustersVsLUMI =
0574         ibooker.bookProfile(histname, histname, LUMIBin, LUMIMin, LUMIMax, NClusStrMin, 3. * NClusStrMax, "");
0575     NumberOfStripClustersVsLUMI->setAxisTitle("online lumi [1e30 Hz cm^{-2}]", 1);
0576     NumberOfStripClustersVsLUMI->setAxisTitle("Mean number of strip clusters", 2);
0577   }
0578 
0579   if (doPlotsVsBXlumi_) {
0580     ibooker.setCurrentFolder(MEFolderName + "/PUmonitoring");
0581     // get binning from the configuration
0582     edm::ParameterSet BXlumiParameters = conf->getParameter<edm::ParameterSet>("BXlumiSetup");
0583     int BXlumiBin = BXlumiParameters.getParameter<int>("BXlumiBin");
0584     double BXlumiMin = BXlumiParameters.getParameter<double>("BXlumiMin");
0585     double BXlumiMax = BXlumiParameters.getParameter<double>("BXlumiMax");
0586 
0587     histname = "NumberOfTracksVsBXlumi_" + CategoryName;
0588     NumberOfTracksVsBXlumi =
0589         ibooker.bookProfile(histname, histname, BXlumiBin, BXlumiMin, BXlumiMax, TKNoMin, 3. * TKNoMax, "");
0590     NumberOfTracksVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]", 1);
0591     NumberOfTracksVsBXlumi->setAxisTitle("Mean number of Tracks", 2);
0592   }
0593 
0594   if (doLumiAnalysis) {
0595     auto scope = DQMStore::IBooker::UseLumiScope(ibooker);
0596     theTrackAnalyzer->initHisto(ibooker, iSetup, *conf);
0597   } else {
0598     theTrackAnalyzer->initHisto(ibooker, iSetup, *conf);
0599   }
0600 
0601   // book the Seed Property histograms
0602   // ---------------------------------------------------------------------------------//
0603 
0604   ibooker.setCurrentFolder(MEFolderName + "/TrackBuilding");
0605 
0606   doSeedNumberPlot = conf->getParameter<bool>("doSeedNumberHisto");
0607   doSeedLumiAnalysis_ = conf->getParameter<bool>("doSeedLumiAnalysis");
0608   doSeedVsClusterPlot = conf->getParameter<bool>("doSeedVsClusterHisto");
0609 
0610   edm::InputTag seedProducer = conf->getParameter<edm::InputTag>("SeedProducer");
0611 
0612   if (doAllSeedPlots || doSeedNumberPlot) {
0613     ibooker.setCurrentFolder(MEFolderName + "/TrackBuilding");
0614     histname = "NumberOfSeeds_" + seedProducer.label() + "_" + CategoryName;
0615     NumberOfSeeds = ibooker.book1D(histname, histname, TKNoSeedBin, TKNoSeedMin, TKNoSeedMax);
0616     NumberOfSeeds->setAxisTitle("Number of Seeds per Event", 1);
0617     NumberOfSeeds->setAxisTitle("Number of Events", 2);
0618 
0619     if (doSeedLumiAnalysis_) {
0620       ibooker.setCurrentFolder(MEFolderName + "/LSanalysis");
0621       auto scope = DQMStore::IBooker::UseLumiScope(ibooker);
0622       histname = "NumberOfSeeds_lumiFlag_" + seedProducer.label() + "_" + CategoryName;
0623       NumberOfSeeds_lumiFlag = ibooker.book1D(histname, histname, TKNoSeedBin, TKNoSeedMin, TKNoSeedMax);
0624       NumberOfSeeds_lumiFlag->setAxisTitle("Number of Seeds per Event", 1);
0625       NumberOfSeeds_lumiFlag->setAxisTitle("Number of Events", 2);
0626     }
0627   }
0628 
0629   if (doAllSeedPlots || doSeedVsClusterPlot) {
0630     ibooker.setCurrentFolder(MEFolderName + "/TrackBuilding");
0631 
0632     ClusterLabels = conf->getParameter<std::vector<std::string> >("ClusterLabels");
0633 
0634     std::vector<double> histoMin, histoMax;
0635     std::vector<int> histoBin;  //these vectors are for max min and nbins in histograms
0636 
0637     int NClusPxBin = conf->getParameter<int>("NClusPxBin");
0638     double NClusPxMin = conf->getParameter<double>("NClusPxMin");
0639     double NClusPxMax = conf->getParameter<double>("NClusPxMax");
0640 
0641     int NClusStrBin = conf->getParameter<int>("NClusStrBin");
0642     double NClusStrMin = conf->getParameter<double>("NClusStrMin");
0643     double NClusStrMax = conf->getParameter<double>("NClusStrMax");
0644 
0645     setMaxMinBin(
0646         histoMin, histoMax, histoBin, NClusStrMin, NClusStrMax, NClusStrBin, NClusPxMin, NClusPxMax, NClusPxBin);
0647 
0648     for (uint i = 0; i < ClusterLabels.size(); i++) {
0649       histname = "SeedsVsClusters_" + seedProducer.label() + "_Vs_" + ClusterLabels[i] + "_" + CategoryName;
0650       SeedsVsClusters.push_back(dynamic_cast<MonitorElement*>(ibooker.book2D(
0651           histname, histname, histoBin[i], histoMin[i], histoMax[i], TKNoSeedBin, TKNoSeedMin, TKNoSeedMax)));
0652       SeedsVsClusters[i]->setAxisTitle("Number of Clusters", 1);
0653       SeedsVsClusters[i]->setAxisTitle("Number of Seeds", 2);
0654       SeedsVsClusters[i]->getTH2F()->SetCanExtend(TH1::kAllAxes);
0655     }
0656   }
0657 
0658   if (doRegionPlots) {
0659     ibooker.setCurrentFolder(MEFolderName + "/TrackBuilding");
0660 
0661     int regionBin = conf->getParameter<int>("RegionSizeBin");
0662     double regionMin = conf->getParameter<double>("RegionSizeMin");
0663     double regionMax = conf->getParameter<double>("RegionSizeMax");
0664 
0665     histname = "TrackingRegionsNumberOf_" + seedProducer.label() + "_" + CategoryName;
0666     NumberOfTrackingRegions = ibooker.book1D(histname, histname, regionBin, regionMin, regionMax);
0667     NumberOfTrackingRegions->setAxisTitle("Number of TrackingRegions per Event", 1);
0668     NumberOfTrackingRegions->setAxisTitle("Number of Events", 2);
0669   }
0670 
0671   if (doTkCandPlots) {
0672     ibooker.setCurrentFolder(MEFolderName + "/TrackBuilding");
0673 
0674     edm::InputTag tcProducer = conf->getParameter<edm::InputTag>("TCProducer");
0675 
0676     histname = "NumberOfTrackCandidates_" + tcProducer.label() + "_" + CategoryName;
0677     NumberOfTrackCandidates = ibooker.book1D(histname, histname, TCNoBin, TCNoMin, TCNoMax);
0678     NumberOfTrackCandidates->setAxisTitle("Number of Track Candidates per Event", 1);
0679     NumberOfTrackCandidates->setAxisTitle("Number of Event", 2);
0680 
0681     histname = "FractionOfCandOverSeeds_" + tcProducer.label() + "_" + CategoryName;
0682     FractionCandidatesOverSeeds = ibooker.book1D(histname, histname, 101, 0., 1.01);
0683     FractionCandidatesOverSeeds->setAxisTitle("Number of Track Candidates / Number of Seeds per Event", 1);
0684     FractionCandidatesOverSeeds->setAxisTitle("Number of Event", 2);
0685   }
0686 
0687   theTrackBuildingAnalyzer->initHisto(ibooker, *conf);
0688 
0689   if (doTrackerSpecific_ || doAllPlots) {
0690     ClusterLabels = conf->getParameter<std::vector<std::string> >("ClusterLabels");
0691 
0692     std::vector<double> histoMin, histoMax;
0693     std::vector<int> histoBin;  //these vectors are for max min and nbins in histograms
0694 
0695     int NClusStrBin = conf->getParameter<int>("NClusStrBin");
0696     double NClusStrMin = conf->getParameter<double>("NClusStrMin");
0697     double NClusStrMax = conf->getParameter<double>("NClusStrMax");
0698 
0699     int NClusPxBin = conf->getParameter<int>("NClusPxBin");
0700     double NClusPxMin = conf->getParameter<double>("NClusPxMin");
0701     double NClusPxMax = conf->getParameter<double>("NClusPxMax");
0702 
0703     edm::ParameterSet ParametersNTrk2D = conf->getParameter<edm::ParameterSet>("NTrk2D");
0704     int NTrk2DBin = ParametersNTrk2D.getParameter<int>("NTrk2DBin");
0705     double NTrk2DMin = ParametersNTrk2D.getParameter<double>("NTrk2DMin");
0706     double NTrk2DMax = ParametersNTrk2D.getParameter<double>("NTrk2DMax");
0707 
0708     setMaxMinBin(
0709         histoMin, histoMax, histoBin, NClusStrMin, NClusStrMax, NClusStrBin, NClusPxMin, NClusPxMax, NClusPxBin);
0710 
0711     ibooker.setCurrentFolder(MEFolderName + "/HitProperties");
0712 
0713     for (uint i = 0; i < ClusterLabels.size(); i++) {
0714       ibooker.setCurrentFolder(MEFolderName + "/HitProperties");
0715       histname = "TracksVs" + ClusterLabels[i] + "Cluster_" + CategoryName;
0716       NumberOfTrkVsClusters.push_back(dynamic_cast<MonitorElement*>(
0717           ibooker.book2D(histname, histname, histoBin[i], histoMin[i], histoMax[i], NTrk2DBin, NTrk2DMin, NTrk2DMax)));
0718       std::string title = "Number of " + ClusterLabels[i] + " Clusters";
0719       if (ClusterLabels[i] == "Tot")
0720         title = "# of Clusters in (Pixel+Strip) Detectors";
0721       NumberOfTrkVsClusters[i]->setAxisTitle(title, 1);
0722       NumberOfTrkVsClusters[i]->setAxisTitle("Number of Tracks", 2);
0723       NumberOfTrkVsClusters[i]->getTH1()->SetCanExtend(TH1::kXaxis);
0724     }
0725   }
0726 
0727   // Initialize the GenericTriggerEventFlag
0728   if (genTriggerEventFlag_->on())
0729     genTriggerEventFlag_->initRun(iRun, iSetup);
0730 }
0731 
0732 /*
0733 // -- BeginRun
0734 //---------------------------------------------------------------------------------//
0735 void TrackingMonitor::beginRun(const edm::Run& iRun, const edm::EventSetup& iSetup)
0736 {
0737   
0738 }
0739 */
0740 
0741 // -- Analyse
0742 // ---------------------------------------------------------------------------------//
0743 void TrackingMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0744   // Filter out events if Trigger Filtering is requested
0745   if (genTriggerEventFlag_->on() && !genTriggerEventFlag_->accept(iEvent, iSetup))
0746     return;
0747   auto const* conf = edm::pset::Registry::instance()->getMapped(confID_);
0748   MEFolderName = conf->getParameter<std::string>("FolderName");
0749   std::string Folder = MEFolderName.substr(0, 2);
0750   float lumi = -1.;
0751   if (forceSCAL_) {
0752     edm::Handle<LumiScalersCollection> lumiScalers = iEvent.getHandle(lumiscalersToken_);
0753     if (lumiScalers.isValid() && !lumiScalers->empty()) {
0754       LumiScalersCollection::const_iterator scalit = lumiScalers->begin();
0755       lumi = scalit->instantLumi();
0756     }
0757   } else {
0758     edm::Handle<OnlineLuminosityRecord> metaData = iEvent.getHandle(metaDataToken_);
0759     if (metaData.isValid())
0760       lumi = metaData->instLumi();
0761   }
0762 
0763   if ((doPlotsVsLUMI_ || doAllPlots) && iEvent.isRealData())
0764     NumberEventsOfVsLUMI->Fill(lumi);
0765 
0766   //  Analyse the tracks
0767   //  if the collection is empty, do not fill anything
0768   // ---------------------------------------------------------------------------------//
0769 
0770   size_t bx = iEvent.bunchCrossing();
0771   if (doPlotsVsBX_ || doAllPlots)
0772     NumberEventsOfVsBX->Fill(bx);
0773 
0774   // get the track collection
0775   edm::Handle<edm::View<reco::Track> > trackHandle = iEvent.getHandle(trackToken_);
0776 
0777   int numberOfTracks_den = 0;
0778   edm::Handle<edm::View<reco::Track> > allTrackHandle = iEvent.getHandle(allTrackToken_);
0779   if (allTrackHandle.isValid()) {
0780     for (edm::View<reco::Track>::const_iterator track = allTrackHandle->begin(); track != allTrackHandle->end();
0781          ++track) {
0782       if (denSelection_(*track))
0783         numberOfTracks_den++;
0784     }
0785   }
0786 
0787   edm::Handle<reco::VertexCollection> pvHandle = iEvent.getHandle(pvSrcToken_);
0788   reco::Vertex const* pv0 = nullptr;
0789   if (pvHandle.isValid()) {
0790     pv0 = &pvHandle->front();
0791     //--- pv fake (the pv collection should have size==1 and the pv==beam spot)
0792     if (pv0->isFake() ||
0793         pv0->tracksSize() == 0
0794         // definition of goodOfflinePrimaryVertex
0795         || pv0->ndof() < pvNDOF_ || pv0->z() > 24.)
0796       pv0 = nullptr;
0797   }
0798 
0799   if (trackHandle.isValid()) {
0800     int numberOfTracks = trackHandle->size();
0801     int numberOfTracks_num = 0;
0802     int numberOfTracks_pv0 = 0;
0803 
0804     const edm::View<reco::Track>& trackCollection = *trackHandle;
0805     // calculate the mean # rechits and layers
0806     int totalRecHits = 0, totalLayers = 0;
0807 
0808     theTrackAnalyzer->setNumberOfGoodVertices(iEvent);
0809     theTrackAnalyzer->setBX(iEvent);
0810     theTrackAnalyzer->setLumi(iEvent, iSetup);
0811     for (edm::View<reco::Track>::const_iterator track = trackCollection.begin(); track != trackCollection.end();
0812          ++track) {
0813       if (doPlotsVsBX_ || doAllPlots)
0814         NumberOfRecHitsPerTrackVsBX->Fill(bx, track->numberOfValidHits());
0815       if (numSelection_(*track)) {
0816         numberOfTracks_num++;
0817         if (pv0 && std::abs(track->dz(pv0->position())) < 0.15)
0818           ++numberOfTracks_pv0;
0819       }
0820 
0821       if (doProfilesVsLS_ || doAllPlots)
0822         NumberOfRecHitsPerTrackVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()),
0823                                           track->numberOfValidHits());
0824 
0825       if ((doPlotsVsLUMI_ || doAllPlots) && iEvent.isRealData())
0826         NumberOfRecHitsPerTrackVsLUMI->Fill(lumi, track->numberOfValidHits());
0827 
0828       totalRecHits += track->numberOfValidHits();
0829       totalLayers += track->hitPattern().trackerLayersWithMeasurement();
0830 
0831       // do analysis per track
0832       theTrackAnalyzer->analyze(iEvent, iSetup, *track);
0833     }
0834 
0835     double frac = -1.;
0836     //      if (numberOfAllTracks > 0) frac = static_cast<double>(numberOfTracks)/static_cast<double>(numberOfAllTracks);
0837     if (numberOfTracks_den > 0)
0838       frac = static_cast<double>(numberOfTracks_num) / static_cast<double>(numberOfTracks_den);
0839 
0840     if (doGeneralPropertiesPlots_ || doAllPlots) {
0841       NumberOfTracks->Fill(double(numberOfTracks));
0842 
0843       if (Folder == "Tr") {
0844         NumberofTracks_Hardvtx->Fill(double(numberOfTracks_pv0));
0845         NumberOfTracks_PUvtx->Fill(double(numberOfTracks - numberOfTracks_pv0));
0846         NumberofTracks_Hardvtx_PUvtx->Fill(0.5, double(numberOfTracks_pv0));
0847         NumberofTracks_Hardvtx_PUvtx->Fill(1.5, double(numberOfTracks - numberOfTracks_pv0));
0848       }
0849 
0850       if (doPlotsVsBX_ || doAllPlots)
0851         NumberOfTracksVsBX->Fill(bx, numberOfTracks);
0852       if ((doPlotsVsLUMI_ || doAllPlots) && iEvent.isRealData())
0853         NumberOfTracksVsLUMI->Fill(lumi, numberOfTracks);
0854       if (doFractionPlot_) {
0855         FractionOfGoodTracks->Fill(frac);
0856 
0857         if (doFractionPlot_) {
0858           if (doPlotsVsBX_ || doAllPlots)
0859             GoodTracksFractionVsBX->Fill(bx, frac);
0860           if ((doPlotsVsLUMI_ || doAllPlots) && iEvent.isRealData())
0861             GoodTracksFractionVsLUMI->Fill(lumi, frac);
0862         }
0863       }
0864       if (numberOfTracks > 0) {
0865         double meanRecHits = static_cast<double>(totalRecHits) / static_cast<double>(numberOfTracks);
0866         double meanLayers = static_cast<double>(totalLayers) / static_cast<double>(numberOfTracks);
0867         NumberOfMeanRecHitsPerTrack->Fill(meanRecHits);
0868         NumberOfMeanLayersPerTrack->Fill(meanLayers);
0869       }
0870     }
0871 
0872     if (doProfilesVsLS_ || doAllPlots) {
0873       float nLS = static_cast<double>(iEvent.id().luminosityBlock());
0874       NumberEventsOfVsLS->Fill(nLS);
0875       NumberOfTracksVsLS->Fill(nLS, numberOfTracks);
0876       if (doFractionPlot_)
0877         GoodTracksFractionVsLS->Fill(nLS, frac);
0878     }
0879 
0880     if (doLumiAnalysis) {
0881       NumberOfTracks_lumiFlag->Fill(numberOfTracks);
0882     }
0883 
0884     //  Analyse the Track Building variables
0885     //  if the collection is empty, do not fill anything
0886     // ---------------------------------------------------------------------------------//
0887 
0888     // fill the TrackCandidate info
0889     if (doTkCandPlots) {
0890       // magnetic field
0891       MagneticField const& theMF = iSetup.getData(magneticFieldToken_);
0892 
0893       // get the candidate collection
0894       edm::Handle<TrackCandidateCollection> theTCHandle = iEvent.getHandle(trackCandidateToken_);
0895       const TrackCandidateCollection& theTCCollection = *theTCHandle;
0896 
0897       if (theTCHandle.isValid()) {
0898         // get the beam spot
0899         edm::Handle<reco::BeamSpot> recoBeamSpotHandle = iEvent.getHandle(bsSrcToken_);
0900         const reco::BeamSpot& bs = *recoBeamSpotHandle;
0901 
0902         NumberOfTrackCandidates->Fill(theTCCollection.size());
0903 
0904         // get the seed collection
0905         edm::Handle<edm::View<TrajectorySeed> > seedHandle = iEvent.getHandle(seedToken_);
0906         const edm::View<TrajectorySeed>& seedCollection = *seedHandle;
0907         if (seedHandle.isValid() && !seedCollection.empty())
0908           FractionCandidatesOverSeeds->Fill(double(theTCCollection.size()) / double(seedCollection.size()));
0909 
0910         TransientTrackingRecHitBuilder const& theTTRHBuilder = iSetup.getData(transientTrackingRecHitBuilderToken_);
0911         for (TrackCandidateCollection::const_iterator cand = theTCCollection.begin(); cand != theTCCollection.end();
0912              ++cand) {
0913           theTrackBuildingAnalyzer->analyze(iEvent, iSetup, *cand, bs, theMF, theTTRHBuilder);
0914         }
0915       } else {
0916         edm::LogWarning("TrackingMonitor") << "No Track Candidates in the event.  Not filling associated histograms";
0917       }
0918 
0919       if (doMVAPlots) {
0920         // Get MVA and quality mask collections
0921         std::vector<const MVACollection*> mvaCollections;
0922         std::vector<const QualityMaskCollection*> qualityMaskCollections;
0923 
0924         edm::Handle<edm::View<reco::Track> > htracks = iEvent.getHandle(mvaTrackToken_);
0925 
0926         edm::Handle<MVACollection> hmva;
0927         edm::Handle<QualityMaskCollection> hqual;
0928         for (const auto& tokenTpl : mvaQualityTokens_) {
0929           hmva = iEvent.getHandle(std::get<0>(tokenTpl));
0930           hqual = iEvent.getHandle(std::get<1>(tokenTpl));
0931           mvaCollections.push_back(hmva.product());
0932           qualityMaskCollections.push_back(hqual.product());
0933         }
0934         theTrackBuildingAnalyzer->analyze(*htracks, mvaCollections, qualityMaskCollections);
0935       }
0936     }
0937 
0938     //plots for trajectory seeds
0939 
0940     if (doAllSeedPlots || doSeedNumberPlot || doSeedVsClusterPlot || runTrackBuildingAnalyzerForSeed) {
0941       // get the seed collection
0942       edm::Handle<edm::View<TrajectorySeed> > seedHandle = iEvent.getHandle(seedToken_);
0943 
0944       // fill the seed info
0945       if (seedHandle.isValid()) {
0946         const auto& seedCollection = *seedHandle;
0947 
0948         if (doAllSeedPlots || doSeedNumberPlot) {
0949           NumberOfSeeds->Fill(seedCollection.size());
0950           if (doSeedLumiAnalysis_)
0951             NumberOfSeeds_lumiFlag->Fill(seedCollection.size());
0952         }
0953 
0954         if (doAllSeedPlots || doSeedVsClusterPlot) {
0955           std::vector<int> NClus;
0956           setNclus(iEvent, NClus);
0957           for (uint i = 0; i < ClusterLabels.size(); i++) {
0958             SeedsVsClusters[i]->Fill(NClus[i], seedCollection.size());
0959           }
0960         }
0961 
0962         if (doAllSeedPlots || runTrackBuildingAnalyzerForSeed) {
0963           edm::Handle<std::vector<SeedStopInfo> > stopHandle = iEvent.getHandle(seedStopInfoToken_);
0964           const auto& seedStopInfo = *stopHandle;
0965 
0966           if (seedStopInfo.size() == seedCollection.size()) {
0967             //here duplication of mag field and be informations is needed to allow seed and track cand histos to be independent
0968             // magnetic field
0969             MagneticField const& theMF = iSetup.getData(magneticFieldToken_);
0970 
0971             // get the beam spot
0972             edm::Handle<reco::BeamSpot> recoBeamSpotHandle = iEvent.getHandle(bsSrcToken_);
0973             const reco::BeamSpot& bs = *recoBeamSpotHandle;
0974 
0975             TransientTrackingRecHitBuilder const& theTTRHBuilder = iSetup.getData(transientTrackingRecHitBuilderToken_);
0976             for (size_t i = 0; i < seedCollection.size(); ++i) {
0977               theTrackBuildingAnalyzer->analyze(
0978                   iEvent, iSetup, seedCollection[i], seedStopInfo[i], bs, theMF, theTTRHBuilder);
0979             }
0980           } else {
0981             edm::LogWarning("TrackingMonitor")
0982                 << "Seed collection size (" << seedCollection.size()
0983                 << ") differs from seed stop info collection size (" << seedStopInfo.size()
0984                 << "). This is a sign of inconsistency in the configuration. Not filling associated histograms.";
0985           }
0986         }
0987 
0988       } else {
0989         edm::LogWarning("TrackingMonitor") << "No Trajectory seeds in the event.  Not filling associated histograms";
0990       }
0991     }
0992 
0993     // plots for tracking regions
0994     if (doRegionPlots) {
0995       if (!regionToken_.isUninitialized()) {
0996         auto const& regions = iEvent.get(regionToken_);
0997         NumberOfTrackingRegions->Fill(regions.size());
0998 
0999         theTrackBuildingAnalyzer->analyze(regions);
1000       } else if (!regionLayerSetsToken_.isUninitialized()) {
1001         edm::Handle<TrackingRegionsSeedingLayerSets> hregions = iEvent.getHandle(regionLayerSetsToken_);
1002         const auto& regions = *hregions;
1003         NumberOfTrackingRegions->Fill(regions.regionsSize());
1004 
1005         theTrackBuildingAnalyzer->analyze(regions);
1006       }
1007 
1008       if (doRegionCandidatePlots) {
1009         edm::Handle<reco::CandidateView> hcandidates = iEvent.getHandle(regionCandidateToken_);
1010         theTrackBuildingAnalyzer->analyze(*hcandidates);
1011       }
1012     }
1013 
1014     if (doTrackerSpecific_ || doAllPlots) {
1015       std::vector<int> NClus;
1016       setNclus(iEvent, NClus);
1017       for (uint i = 0; i < ClusterLabels.size(); i++) {
1018         NumberOfTrkVsClusters[i]->Fill(NClus[i], numberOfTracks);
1019       }
1020     }
1021 
1022     if (doPUmonitoring_) {
1023       // do vertex monitoring
1024       for (size_t i = 0; i < theVertexMonitor.size(); i++)
1025         theVertexMonitor[i]->analyze(iEvent, iSetup);
1026     }
1027     if (doPlotsVsGoodPVtx_) {
1028       size_t totalNumGoodPV = 0;
1029       if (pvHandle.isValid()) {
1030         for (reco::VertexCollection::const_iterator pv = pvHandle->begin(); pv != pvHandle->end(); ++pv) {
1031           //--- pv fake (the pv collection should have size==1 and the pv==beam spot)
1032           if (pv->isFake() || pv->tracksSize() == 0)
1033             continue;
1034 
1035           // definition of goodOfflinePrimaryVertex
1036           if (pv->ndof() < pvNDOF_ || pv->z() > 24.)
1037             continue;
1038           totalNumGoodPV++;
1039         }
1040 
1041         NumberEventsOfVsGoodPVtx->Fill(float(totalNumGoodPV));
1042         NumberOfTracksVsGoodPVtx->Fill(float(totalNumGoodPV), numberOfTracks);
1043         if (totalNumGoodPV > 1)
1044           NumberOfTracksVsPUPVtx->Fill(totalNumGoodPV - 1,
1045                                        double(numberOfTracks - numberOfTracks_pv0) / double(totalNumGoodPV - 1));
1046         NumberOfPVtxVsGoodPVtx->Fill(float(totalNumGoodPV), pvHandle->size());
1047 
1048         for (edm::View<reco::Track>::const_iterator track = trackCollection.begin(); track != trackCollection.end();
1049              ++track) {
1050           NumberOfRecHitsPerTrackVsGoodPVtx->Fill(float(totalNumGoodPV), track->numberOfValidHits());
1051         }
1052 
1053         if (doProfilesVsLS_ || doAllPlots)
1054           NumberOfGoodPVtxVsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()), totalNumGoodPV);
1055         if (doPlotsVsBX_ || doAllPlots)
1056           NumberOfGoodPVtxVsBX->Fill(bx, float(totalNumGoodPV));
1057 
1058         if (doFractionPlot_)
1059           GoodTracksFractionVsGoodPVtx->Fill(float(totalNumGoodPV), frac);
1060 
1061         if ((doPlotsVsLUMI_ || doAllPlots) && iEvent.isRealData())
1062           NumberOfGoodPVtxVsLUMI->Fill(lumi, float(totalNumGoodPV));
1063       }
1064 
1065       std::vector<int> NClus;
1066       setNclus(iEvent, NClus);
1067       std::ostringstream ss;
1068       ss << "VI stat " << totalNumGoodPV << ' ' << numberOfTracks;
1069       for (uint i = 0; i < ClusterLabels.size(); i++) {
1070         ss << ' ' << NClus[i];
1071         if ((doPlotsVsLUMI_ || doAllPlots) && iEvent.isRealData()) {
1072           if (ClusterLabels[i] == "Pix")
1073             NumberOfPixelClustersVsLUMI->Fill(lumi, NClus[i]);
1074           if (ClusterLabels[i] == "Strip")
1075             NumberOfStripClustersVsLUMI->Fill(lumi, NClus[i]);
1076         }
1077         if (ClusterLabels[i] == "Pix")
1078           NumberOfPixelClustersVsGoodPVtx->Fill(float(totalNumGoodPV), NClus[i]);
1079         if (ClusterLabels[i] == "Strip")
1080           NumberOfStripClustersVsGoodPVtx->Fill(float(totalNumGoodPV), NClus[i]);
1081       }
1082       COUT(MEFolderName) << ss.str() << std::endl;
1083       if (doPlotsVsBXlumi_) {
1084         double bxlumi = theLumiDetails_->getValue(iEvent);
1085         NumberOfTracksVsBXlumi->Fill(bxlumi, numberOfTracks);
1086       }
1087 
1088       if (doProfilesVsLS_ || doAllPlots)
1089         if (totalNumGoodPV != 0)
1090           NumberOfGoodPVtxWO0VsLS->Fill(static_cast<double>(iEvent.id().luminosityBlock()), float(totalNumGoodPV));
1091       if (doPlotsVsBX_ || doAllPlots)
1092         if (totalNumGoodPV != 0)
1093           NumberOfGoodPVtxWO0VsBX->Fill(bx, float(totalNumGoodPV));
1094       if ((doPlotsVsLUMI_ || doAllPlots) && iEvent.isRealData())
1095         if (totalNumGoodPV != 0)
1096           NumberOfGoodPVtxWO0VsLUMI->Fill(lumi, float(totalNumGoodPV));
1097 
1098     }  // PU monitoring
1099 
1100   }  // trackHandle is valid
1101 }
1102 
1103 void TrackingMonitor::setMaxMinBin(std::vector<double>& arrayMin,
1104                                    std::vector<double>& arrayMax,
1105                                    std::vector<int>& arrayBin,
1106                                    double smin,
1107                                    double smax,
1108                                    int sbin,
1109                                    double pmin,
1110                                    double pmax,
1111                                    int pbin) {
1112   arrayMin.resize(ClusterLabels.size());
1113   arrayMax.resize(ClusterLabels.size());
1114   arrayBin.resize(ClusterLabels.size());
1115 
1116   for (uint i = 0; i < ClusterLabels.size(); ++i) {
1117     if (ClusterLabels[i] == "Pix") {
1118       arrayMin[i] = pmin;
1119       arrayMax[i] = pmax;
1120       arrayBin[i] = pbin;
1121     } else if (ClusterLabels[i] == "Strip") {
1122       arrayMin[i] = smin;
1123       arrayMax[i] = smax;
1124       arrayBin[i] = sbin;
1125     } else if (ClusterLabels[i] == "Tot") {
1126       arrayMin[i] = smin;
1127       arrayMax[i] = smax + pmax;
1128       arrayBin[i] = sbin;
1129     } else {
1130       edm::LogWarning("TrackingMonitor") << "Cluster Label " << ClusterLabels[i]
1131                                          << " not defined, using strip parameters ";
1132       arrayMin[i] = smin;
1133       arrayMax[i] = smax;
1134       arrayBin[i] = sbin;
1135     }
1136   }
1137 }
1138 
1139 void TrackingMonitor::setNclus(const edm::Event& iEvent, std::vector<int>& arrayNclus) {
1140   int ncluster_pix = -1;
1141   int ncluster_strip = -1;
1142 
1143   edm::Handle<edmNew::DetSetVector<SiStripCluster> > strip_clusters = iEvent.getHandle(stripClustersToken_);
1144   edm::Handle<edmNew::DetSetVector<SiPixelCluster> > pixel_clusters = iEvent.getHandle(pixelClustersToken_);
1145 
1146   if (strip_clusters.isValid() && pixel_clusters.isValid()) {
1147     ncluster_pix = (*pixel_clusters).dataSize();
1148     ncluster_strip = (*strip_clusters).dataSize();
1149   }
1150 
1151   arrayNclus.resize(ClusterLabels.size());
1152   for (uint i = 0; i < ClusterLabels.size(); ++i) {
1153     if (ClusterLabels[i] == "Pix")
1154       arrayNclus[i] = ncluster_pix;
1155     else if (ClusterLabels[i] == "Strip")
1156       arrayNclus[i] = ncluster_strip;
1157     else if (ClusterLabels[i] == "Tot")
1158       arrayNclus[i] = ncluster_pix + ncluster_strip;
1159     else {
1160       edm::LogWarning("TrackingMonitor") << "Cluster Label " << ClusterLabels[i]
1161                                          << " not defined using stri parametrs ";
1162       arrayNclus[i] = ncluster_strip;
1163     }
1164   }
1165 }
1166 DEFINE_FWK_MODULE(TrackingMonitor);