Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-06-22 02:23:34

0001 // system includes
0002 #include <map>
0003 #include <set>
0004 #include <string>
0005 #include <vector>
0006 
0007 // user includes
0008 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0009 #include "DQMServices/Core/interface/DQMStore.h"
0010 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0011 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0012 #include "DataFormats/Common/interface/Handle.h"
0013 #include "DataFormats/TrackCandidate/interface/TrajectoryStopReasons.h"
0014 #include "DataFormats/DetId/interface/DetId.h"
0015 #include "DataFormats/JetReco/interface/PFJet.h"
0016 #include "DataFormats/TrackReco/interface/Track.h"
0017 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0018 #include "DataFormats/TrackerRecHit2D/interface/BaseTrackerRecHit.h"
0019 #include "DataFormats/TrackerRecHit2D/interface/OmniClusterRef.h"
0020 #include "DataFormats/TrackerRecHit2D/interface/SiStripMatchedRecHit2D.h"
0021 #include "DataFormats/VertexReco/interface/Vertex.h"
0022 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0023 #include "DataFormats/Math/interface/deltaR.h"
0024 #include "FWCore/Framework/interface/ESHandle.h"
0025 #include "FWCore/Framework/interface/Event.h"
0026 #include "FWCore/Framework/interface/Frameworkfwd.h"
0027 #include "FWCore/Framework/interface/MakerMacros.h"
0028 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0029 #include "FWCore/Utilities/interface/InputTag.h"
0030 #include "Geometry/CommonDetUnit/interface/GluedGeomDet.h"
0031 #include "Geometry/Records/interface/TrackerDigiGeometryRecord.h"
0032 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
0033 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0034 #include "RecoLocalTracker/SiStripClusterizer/interface/SiStripClusterInfo.h"
0035 #include "SimDataFormats/PileupSummaryInfo/interface/PileupSummaryInfo.h"
0036 #include "TrackingTools/IPTools/interface/IPTools.h"
0037 #include "TrackingTools/Records/interface/TransientTrackRecord.h"
0038 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0039 #include "TrackingTools/TransientTrack/interface/TransientTrackBuilder.h"
0040 #include "RecoVertex/VertexPrimitives/interface/TransientVertex.h"
0041 #include "RecoVertex/KalmanVertexFit/interface/KalmanVertexFitter.h"
0042 #include "RecoVertex/VertexTools/interface/VertexDistance3D.h"
0043 
0044 // ROOT includes
0045 #include "TFile.h"
0046 #include "TH1.h"
0047 #include "TLorentzVector.h"
0048 #include "TMath.h"
0049 #include "TPRegexp.h"
0050 
0051 using MVACollection = std::vector<float>;
0052 using QualityMaskCollection = std::vector<unsigned char>;
0053 
0054 class StandaloneTrackMonitor : public DQMEDAnalyzer {
0055 public:
0056   StandaloneTrackMonitor(const edm::ParameterSet&);
0057   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0058 
0059 protected:
0060   void analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) override;
0061   void addClusterToMap(uint32_t detid, const SiStripCluster* cluster);
0062   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0063   void dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) override;
0064   void processClusters(edm::Event const& iEvent,
0065                        edm::EventSetup const& iSetup,
0066                        const TrackerGeometry& tkGeom,
0067                        double wfac = 1);
0068   void processHit(const TrackingRecHit& recHit,
0069                   edm::EventSetup const& iSetup,
0070                   const TrackerGeometry& tkGeom,
0071                   double wfac = 1);
0072 
0073 private:
0074   const reco::Vertex* findClosestVertex(const TransientVertex aTransVtx, const reco::VertexCollection* vertices)
0075       const;  //same as in /DQMOffline/Alignment/src/DiMuonVertexSelector.cc
0076   const std::string moduleName_;
0077   const std::string folderName_;
0078 
0079   const bool isRECO_;
0080 
0081   SiStripClusterInfo siStripClusterInfo_;
0082 
0083   const edm::InputTag trackTag_;
0084   const edm::InputTag bsTag_;
0085   const edm::InputTag vertexTag_;
0086   const edm::InputTag puSummaryTag_;
0087   const edm::InputTag clusterTag_;
0088   const edm::InputTag jetsTag_;
0089   const edm::EDGetTokenT<reco::TrackCollection> trackToken_;
0090   const edm::EDGetTokenT<reco::BeamSpot> bsToken_;
0091   const edm::EDGetTokenT<reco::VertexCollection> vertexToken_;
0092   const edm::EDGetTokenT<std::vector<PileupSummaryInfo> > puSummaryToken_;
0093   const edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster> > clusterToken_;
0094   const edm::EDGetTokenT<std::vector<reco::PFJet> > jetsToken_;
0095   // track MVA
0096   const std::string trackQuality_;
0097   const bool doPUCorrection_;
0098   const bool doTrackCorrection_;
0099   const bool isMC_;
0100   const bool haveAllHistograms_;
0101   const std::string puScaleFactorFile_;
0102   const std::string trackScaleFactorFile_;
0103   const std::vector<std::string> mvaProducers_;
0104   const edm::InputTag mvaTrackTag_;
0105   const edm::EDGetTokenT<edm::View<reco::Track> > mvaTrackToken_;
0106   const edm::InputTag tcProducer_;
0107   const std::string algoName_;
0108 
0109   int nevt = 0;
0110   int chi2it = 0, chi2itGt = 0, chi2itLt = 0;
0111   const bool verbose_;
0112   std::vector<std::tuple<edm::EDGetTokenT<MVACollection>, edm::EDGetTokenT<QualityMaskCollection> > > mvaQualityTokens_;
0113   const edm::ESGetToken<TrackerGeometry, TrackerDigiGeometryRecord> geomToken_;
0114   const edm::ESGetToken<TransientTrackBuilder, TransientTrackRecord> transTrackToken_;
0115   const edm::ParameterSet TrackEtaHistoPar_;
0116   const edm::ParameterSet TrackPhiHistoPar_;
0117   const edm::ParameterSet TrackPtHistoPar_;
0118   const TrackerGeometry* tkGeom_ = nullptr;
0119 
0120   MonitorElement* trackEtaH_;
0121   MonitorElement* trackEtaerrH_;
0122   MonitorElement* trackCosThetaH_;
0123   MonitorElement* trackThetaerrH_;
0124   MonitorElement* trackPhiH_;
0125   MonitorElement* trackPhierrH_;
0126   MonitorElement* trackPH_;
0127   MonitorElement* trackPtH_;
0128   MonitorElement* trackPt_ZoomH_;
0129   MonitorElement* trackPtUpto2GeVH_;
0130   MonitorElement* trackPtOver10GeVH_;
0131   MonitorElement* trackPterrH_;
0132   MonitorElement* trackqOverpH_;
0133   MonitorElement* trackqOverperrH_;
0134   MonitorElement* trackChargeH_;
0135   MonitorElement* trackChi2H_;
0136   MonitorElement* tracknDOFH_;
0137   MonitorElement* trackChi2ProbH_;
0138   MonitorElement* trackChi2oNDFH_;
0139   MonitorElement* trackd0H_;
0140   MonitorElement* trackChi2bynDOFH_;
0141   MonitorElement* trackalgoH_;
0142   MonitorElement* trackorigalgoH_;
0143   MonitorElement* trackStoppingSourceH_;
0144 
0145   MonitorElement* DistanceOfClosestApproachToPVH_;
0146   MonitorElement* DistanceOfClosestApproachToPVZoomedH_;
0147   MonitorElement* DistanceOfClosestApproachToPVVsPhiH_;
0148   MonitorElement* xPointOfClosestApproachVsZ0wrtPVH_;
0149   MonitorElement* yPointOfClosestApproachVsZ0wrtPVH_;
0150   MonitorElement* trackDeltaRwrtClosestTrack_;
0151 
0152   MonitorElement* ip2dToPVH_;
0153   MonitorElement* iperr2dToPVH_;
0154 
0155   MonitorElement* ip3dToBSH_;
0156   MonitorElement* iperr3dToBSH_;
0157   MonitorElement* iperr3dToBSWtH_;
0158   MonitorElement* iperr2dToBSH_;
0159   MonitorElement* iperr2dToPVWtH_;
0160 
0161   MonitorElement* ip2dToBSH_;
0162   MonitorElement* sip2dToBSH_;
0163 
0164   MonitorElement* ip3dToPVH_;
0165   MonitorElement* iperr3dToPVH_;
0166   MonitorElement* iperr3dToPVWtH_;
0167   MonitorElement* sip3dToPVH_;
0168   MonitorElement* sip3dToBSH_;
0169   MonitorElement* sip2dToPVH_;
0170   MonitorElement* sip2dToPVWtH_;
0171   MonitorElement* sipDxyToPVH_;
0172   MonitorElement* sipDzToPVH_;
0173 
0174   MonitorElement* ip3dToPV2validpixelhitsH_;
0175   MonitorElement* ip3dToBS2validpixelhitsH_;
0176   MonitorElement* iperr3dToPV2validpixelhitsH_;
0177   MonitorElement* iperr3dToBS2validpixelhitsH_;
0178   MonitorElement* sip3dToPV2validpixelhitsH_;
0179   MonitorElement* sip3dToBS2validpixelhitsH_;
0180 
0181   MonitorElement* nallHitsH_;
0182   MonitorElement* ntrackerHitsH_;
0183 
0184   MonitorElement* nvalidTrackerHitsH_;
0185   MonitorElement* nvalidPixelHitsH_;
0186   MonitorElement* nvalidPixelBHitsH_;
0187   MonitorElement* nvalidPixelEHitsH_;
0188   MonitorElement* nvalidStripHitsH_;
0189   MonitorElement* nvalidTIBHitsH_;
0190   MonitorElement* nvalidTOBHitsH_;
0191   MonitorElement* nvalidTIDHitsH_;
0192   MonitorElement* nvalidTECHitsH_;
0193 
0194   MonitorElement* nlostTrackerHitsH_;
0195   MonitorElement* nlostPixelHitsH_;
0196   MonitorElement* nlostPixelBHitsH_;
0197   MonitorElement* nlostPixelEHitsH_;
0198   MonitorElement* nlostStripHitsH_;
0199   MonitorElement* nlostTIBHitsH_;
0200   MonitorElement* nlostTOBHitsH_;
0201   MonitorElement* nlostTIDHitsH_;
0202   MonitorElement* nlostTECHitsH_;
0203 
0204   MonitorElement* nMissingInnerHitBH_;
0205   MonitorElement* nMissingInnerHitEH_;
0206   MonitorElement* nMissingOuterHitBH_;
0207   MonitorElement* nMissingOuterHitEH_;
0208 
0209   MonitorElement* trkLayerwithMeasurementH_;
0210   MonitorElement* pixelLayerwithMeasurementH_;
0211   MonitorElement* pixelBLayerwithMeasurementH_;
0212   MonitorElement* pixelELayerwithMeasurementH_;
0213   MonitorElement* stripLayerwithMeasurementH_;
0214   MonitorElement* stripTIBLayerwithMeasurementH_;
0215   MonitorElement* stripTOBLayerwithMeasurementH_;
0216   MonitorElement* stripTIDLayerwithMeasurementH_;
0217   MonitorElement* stripTECLayerwithMeasurementH_;
0218 
0219   MonitorElement* nlostHitsH_;
0220   MonitorElement* nMissingExpectedInnerHitsH_;
0221   MonitorElement* nMissingExpectedOuterHitsH_;
0222 
0223   MonitorElement* beamSpotXYposH_;
0224   MonitorElement* beamSpotXYposerrH_;
0225   MonitorElement* beamSpotZposH_;
0226   MonitorElement* beamSpotZposerrH_;
0227 
0228   MonitorElement* vertexXposH_;
0229   MonitorElement* vertexYposH_;
0230   MonitorElement* vertexZposH_;
0231   MonitorElement* nVertexH_;
0232   MonitorElement* nVtxH_;
0233 
0234   MonitorElement* nPixBarrelH_;
0235   MonitorElement* nPixEndcapH_;
0236   MonitorElement* nStripTIBH_;
0237   MonitorElement* nStripTOBH_;
0238   MonitorElement* nStripTECH_;
0239   MonitorElement* nStripTIDH_;
0240   MonitorElement* nTracksH_;
0241 
0242   // MC only
0243   MonitorElement* bunchCrossingH_;
0244   MonitorElement* nPUH_;
0245   MonitorElement* trueNIntH_;
0246 
0247   // Exclusive Quantities
0248   MonitorElement* nLostHitByLayerH_;
0249   MonitorElement* nLostHitByLayerPixH_;
0250   MonitorElement* nLostHitByLayerStripH_;
0251   MonitorElement* nLostHitsVspTH_;
0252   MonitorElement* nLostHitsVsEtaH_;
0253   MonitorElement* nLostHitsVsCosThetaH_;
0254   MonitorElement* nLostHitsVsPhiH_;
0255   MonitorElement* nLostHitsVsIterationH_;
0256 
0257   MonitorElement* nHitsTIBSVsEtaH_;
0258   MonitorElement* nHitsTOBSVsEtaH_;
0259   MonitorElement* nHitsTECSVsEtaH_;
0260   MonitorElement* nHitsTIDSVsEtaH_;
0261   MonitorElement* nHitsStripSVsEtaH_;
0262 
0263   MonitorElement* nHitsTIBDVsEtaH_;
0264   MonitorElement* nHitsTOBDVsEtaH_;
0265   MonitorElement* nHitsTECDVsEtaH_;
0266   MonitorElement* nHitsTIDDVsEtaH_;
0267   MonitorElement* nHitsStripDVsEtaH_;
0268 
0269   MonitorElement* nValidHitsVspTH_;
0270   MonitorElement* nValidHitsVsnVtxH_;
0271   MonitorElement* nValidHitsVsEtaH_;
0272   MonitorElement* nValidHitsVsCosThetaH_;
0273   MonitorElement* nValidHitsVsPhiH_;
0274 
0275   MonitorElement* nValidHitsPixVsEtaH_;
0276   MonitorElement* nValidHitsPixBVsEtaH_;
0277   MonitorElement* nValidHitsPixEVsEtaH_;
0278   MonitorElement* nValidHitsStripVsEtaH_;
0279   MonitorElement* nValidHitsTIBVsEtaH_;
0280   MonitorElement* nValidHitsTOBVsEtaH_;
0281   MonitorElement* nValidHitsTECVsEtaH_;
0282   MonitorElement* nValidHitsTIDVsEtaH_;
0283 
0284   MonitorElement* nValidHitsPixVsPhiH_;
0285   MonitorElement* nValidHitsPixBVsPhiH_;
0286   MonitorElement* nValidHitsPixEVsPhiH_;
0287   MonitorElement* nValidHitsStripVsPhiH_;
0288   MonitorElement* nValidHitsTIBVsPhiH_;
0289   MonitorElement* nValidHitsTOBVsPhiH_;
0290   MonitorElement* nValidHitsTECVsPhiH_;
0291   MonitorElement* nValidHitsTIDVsPhiH_;
0292 
0293   MonitorElement* nLostHitsPixVsEtaH_;
0294   MonitorElement* nLostHitsPixBVsEtaH_;
0295   MonitorElement* nLostHitsPixEVsEtaH_;
0296   MonitorElement* nLostHitsStripVsEtaH_;
0297   MonitorElement* nLostHitsTIBVsEtaH_;
0298   MonitorElement* nLostHitsTOBVsEtaH_;
0299   MonitorElement* nLostHitsTECVsEtaH_;
0300   MonitorElement* nLostHitsTIDVsEtaH_;
0301 
0302   MonitorElement* nLostHitsPixVsIterationH_;
0303   MonitorElement* nLostHitsPixBVsIterationH_;
0304   MonitorElement* nLostHitsPixEVsIterationH_;
0305   MonitorElement* nLostHitsStripVsIterationH_;
0306   MonitorElement* nLostHitsTIBVsIterationH_;
0307   MonitorElement* nLostHitsTOBVsIterationH_;
0308   MonitorElement* nLostHitsTECVsIterationH_;
0309   MonitorElement* nLostHitsTIDVsIterationH_;
0310 
0311   MonitorElement* nLostHitsPixVsPhiH_;
0312   MonitorElement* nLostHitsPixBVsPhiH_;
0313   MonitorElement* nLostHitsPixEVsPhiH_;
0314   MonitorElement* nLostHitsStripVsPhiH_;
0315   MonitorElement* nLostHitsTIBVsPhiH_;
0316   MonitorElement* nLostHitsTOBVsPhiH_;
0317   MonitorElement* nLostHitsTECVsPhiH_;
0318   MonitorElement* nLostHitsTIDVsPhiH_;
0319 
0320   MonitorElement* trackChi2oNDFVsEtaH_;
0321   MonitorElement* trackChi2oNDFVsPhiH_;
0322   MonitorElement* trackChi2probVsEtaH_;
0323   MonitorElement* trackChi2probVsPhiH_;
0324 
0325   MonitorElement* trackIperr3dVsEtaH_;
0326   MonitorElement* trackIperr3dVsChi2probH_;
0327 
0328   MonitorElement* trackSip2dVsEtaH_;
0329 
0330   MonitorElement* trackIperr3dVsEta2DH_;
0331   MonitorElement* trackIperr3dVsChi2prob2DH_;
0332 
0333   MonitorElement* trackSip2dVsEta2DH_;
0334   MonitorElement* trackSip2dVsChi2prob2DH_;
0335 
0336   MonitorElement* hOnTrkClusChargeThinH_;
0337   MonitorElement* hOnTrkClusWidthThinH_;
0338   MonitorElement* hOnTrkClusChargeThickH_;
0339   MonitorElement* hOnTrkClusWidthThickH_;
0340 
0341   MonitorElement* hOffTrkClusChargeThinH_;
0342   MonitorElement* hOffTrkClusWidthThinH_;
0343   MonitorElement* hOffTrkClusChargeThickH_;
0344   MonitorElement* hOffTrkClusWidthThickH_;
0345 
0346   std::vector<MonitorElement*> trackMVAs;
0347   std::vector<MonitorElement*> trackMVAsHP;
0348   std::vector<MonitorElement*> trackMVAsVsPtProfile;
0349   std::vector<MonitorElement*> trackMVAsHPVsPtProfile;
0350   std::vector<MonitorElement*> trackMVAsVsEtaProfile;
0351   std::vector<MonitorElement*> trackMVAsHPVsEtaProfile;
0352 
0353   MonitorElement* nJet_;
0354   MonitorElement* Jet_pt_;
0355   MonitorElement* Jet_eta_;
0356   MonitorElement* Jet_energy_;
0357   MonitorElement* Jet_chargedMultiplicity_;
0358 
0359   MonitorElement* Zpt_;
0360   MonitorElement* ZInvMass_;
0361 
0362   MonitorElement* cosPhi3DdileptonH_;
0363 
0364   std::vector<int> lumivec1;
0365   std::vector<int> lumivec2;
0366   std::vector<float> vpu_;
0367   std::vector<float> vtrack_;
0368   std::map<uint32_t, std::set<const SiStripCluster*> > clusterMap_;
0369 };
0370 
0371 void StandaloneTrackMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0372   edm::ParameterSetDescription desc;
0373   desc.addUntracked<std::string>("moduleName", "StandaloneTrackMonitor");
0374   desc.addUntracked<std::string>("folderName", "highPurityTracks");
0375   desc.addUntracked<bool>("isRECO", false);
0376   desc.addUntracked<edm::InputTag>("trackInputTag", edm::InputTag("generalTracks"));
0377   desc.addUntracked<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"));
0378   desc.addUntracked<edm::InputTag>("vertexTag", edm::InputTag("offlinePrimaryVertices"));
0379   desc.addUntracked<edm::InputTag>("puTag", edm::InputTag("addPileupInfo"));
0380   desc.addUntracked<edm::InputTag>("clusterTag", edm::InputTag("siStripClusters"));
0381   desc.addUntracked<edm::InputTag>("PFJetsCollection", edm::InputTag("ak4PFJetsCHS"));
0382 
0383   desc.addUntracked<std::string>("trackQuality", "highPurity");
0384   desc.addUntracked<bool>("doPUCorrection", false);
0385   desc.addUntracked<bool>("doTrackCorrection", false);
0386   desc.addUntracked<bool>("isMC", false);
0387   desc.addUntracked<bool>("haveAllHistograms", false);
0388   desc.addUntracked<std::string>("puScaleFactorFile", "PileupScaleFactor.root");
0389   desc.addUntracked<std::string>("trackScaleFactorFile", "PileupScaleFactor.root");
0390   desc.addUntracked<std::vector<std::string> >("MVAProducers");
0391   desc.addUntracked<edm::InputTag>("TrackProducerForMVA");
0392 
0393   desc.addUntracked<edm::InputTag>("TCProducer");
0394   desc.addUntracked<std::string>("AlgoName");
0395   desc.addUntracked<bool>("verbose", false);
0396 
0397   {
0398     edm::ParameterSetDescription TrackEtaHistoPar;
0399     TrackEtaHistoPar.add<int>("Xbins", 60);
0400     TrackEtaHistoPar.add<double>("Xmin", -3.0);
0401     TrackEtaHistoPar.add<double>("Xmax", 3.0);
0402     desc.add<edm::ParameterSetDescription>("trackEtaH", TrackEtaHistoPar);
0403   }
0404 
0405   {
0406     edm::ParameterSetDescription TrackPtHistoPar;
0407     TrackPtHistoPar.add<int>("Xbins", 60);
0408     TrackPtHistoPar.add<double>("Xmin", 0.0);
0409     TrackPtHistoPar.add<double>("Xmax", 100.0);
0410     desc.add<edm::ParameterSetDescription>("trackPtH", TrackPtHistoPar);
0411   }
0412 
0413   {
0414     edm::ParameterSetDescription TrackPhiHistoPar;
0415     TrackPhiHistoPar.add<int>("Xbins", 100);
0416     TrackPhiHistoPar.add<double>("Xmin", -M_PI);
0417     TrackPhiHistoPar.add<double>("Xmax", M_PI);
0418     desc.add<edm::ParameterSetDescription>("trackPhiH", TrackPhiHistoPar);
0419   }
0420 
0421   descriptions.add("standaloneTrackMonitorDefault", desc);
0422 }
0423 
0424 // -----------------------------
0425 // constructors and destructor
0426 // -----------------------------
0427 StandaloneTrackMonitor::StandaloneTrackMonitor(const edm::ParameterSet& ps)
0428     : moduleName_(ps.getUntrackedParameter<std::string>("moduleName", "StandaloneTrackMonitor")),
0429       folderName_(ps.getUntrackedParameter<std::string>("folderName", "highPurityTracks")),
0430       isRECO_(ps.getUntrackedParameter<bool>("isRECO", false)),
0431       siStripClusterInfo_(consumesCollector()),
0432       trackTag_(ps.getUntrackedParameter<edm::InputTag>("trackInputTag", edm::InputTag("generalTracks"))),
0433       bsTag_(ps.getUntrackedParameter<edm::InputTag>("offlineBeamSpot", edm::InputTag("offlineBeamSpot"))),
0434       vertexTag_(ps.getUntrackedParameter<edm::InputTag>("vertexTag", edm::InputTag("offlinePrimaryVertices"))),
0435       puSummaryTag_(ps.getUntrackedParameter<edm::InputTag>("puTag", edm::InputTag("addPileupInfo"))),
0436       clusterTag_(ps.getUntrackedParameter<edm::InputTag>("clusterTag", edm::InputTag("siStripClusters"))),
0437       jetsTag_(ps.getUntrackedParameter<edm::InputTag>("PFJetsCollection", edm::InputTag("ak4PFJetsCHS"))),
0438       trackToken_(consumes<reco::TrackCollection>(trackTag_)),
0439       bsToken_(consumes<reco::BeamSpot>(bsTag_)),
0440       vertexToken_(consumes<reco::VertexCollection>(vertexTag_)),
0441       puSummaryToken_(consumes<std::vector<PileupSummaryInfo> >(puSummaryTag_)),
0442       clusterToken_(consumes<edmNew::DetSetVector<SiStripCluster> >(clusterTag_)),
0443       jetsToken_(consumes<std::vector<reco::PFJet> >(jetsTag_)),
0444       trackQuality_(ps.getUntrackedParameter<std::string>("trackQuality", "highPurity")),
0445       doPUCorrection_(ps.getUntrackedParameter<bool>("doPUCorrection", false)),
0446       doTrackCorrection_(ps.getUntrackedParameter<bool>("doTrackCorrection", false)),
0447       isMC_(ps.getUntrackedParameter<bool>("isMC", false)),
0448       haveAllHistograms_(ps.getUntrackedParameter<bool>("haveAllHistograms", false)),
0449       puScaleFactorFile_(ps.getUntrackedParameter<std::string>("puScaleFactorFile", "PileupScaleFactor.root")),
0450       trackScaleFactorFile_(ps.getUntrackedParameter<std::string>("trackScaleFactorFile", "PileupScaleFactor.root")),
0451       mvaProducers_(ps.getUntrackedParameter<std::vector<std::string> >("MVAProducers")),
0452       mvaTrackTag_(ps.getUntrackedParameter<edm::InputTag>("TrackProducerForMVA")),
0453       mvaTrackToken_(consumes<edm::View<reco::Track> >(mvaTrackTag_)),
0454       tcProducer_(ps.getUntrackedParameter<edm::InputTag>("TCProducer")),
0455       algoName_(ps.getUntrackedParameter<std::string>("AlgoName")),
0456       verbose_(ps.getUntrackedParameter<bool>("verbose", false)),
0457       geomToken_(esConsumes<TrackerGeometry, TrackerDigiGeometryRecord, edm::Transition::BeginRun>()),
0458       transTrackToken_(esConsumes<TransientTrackBuilder, TransientTrackRecord, edm::Transition::Event>(
0459           edm::ESInputTag{"", "TransientTrackBuilder"})),
0460       TrackEtaHistoPar_(ps.getParameter<edm::ParameterSet>("trackEtaH")),
0461       TrackPhiHistoPar_(ps.getParameter<edm::ParameterSet>("trackPhiH")),
0462       TrackPtHistoPar_(ps.getParameter<edm::ParameterSet>("trackPtH")) {
0463   for (const auto& v : mvaProducers_) {
0464     mvaQualityTokens_.push_back(std::make_tuple(consumes<MVACollection>(edm::InputTag(v, "MVAValues")),
0465                                                 consumes<QualityMaskCollection>(edm::InputTag(v, "QualityMasks"))));
0466   }
0467 
0468   trackEtaH_ = nullptr;
0469   trackEtaerrH_ = nullptr;
0470   trackCosThetaH_ = nullptr;
0471   trackThetaerrH_ = nullptr;
0472   trackPhiH_ = nullptr;
0473   trackPhierrH_ = nullptr;
0474   trackPH_ = nullptr;
0475   trackPtH_ = nullptr;
0476   trackPtUpto2GeVH_ = nullptr;
0477   trackPtOver10GeVH_ = nullptr;
0478   trackPterrH_ = nullptr;
0479   trackqOverpH_ = nullptr;
0480   trackqOverperrH_ = nullptr;
0481   trackChargeH_ = nullptr;
0482   nlostHitsH_ = nullptr;
0483   nvalidTrackerHitsH_ = nullptr;
0484   nvalidPixelHitsH_ = nullptr;
0485   nvalidStripHitsH_ = nullptr;
0486   trkLayerwithMeasurementH_ = nullptr;
0487   pixelLayerwithMeasurementH_ = nullptr;
0488   stripLayerwithMeasurementH_ = nullptr;
0489   beamSpotXYposH_ = nullptr;
0490   beamSpotXYposerrH_ = nullptr;
0491   beamSpotZposH_ = nullptr;
0492   beamSpotZposerrH_ = nullptr;
0493   trackChi2H_ = nullptr;
0494   tracknDOFH_ = nullptr;
0495   trackd0H_ = nullptr;
0496   trackChi2bynDOFH_ = nullptr;
0497   vertexXposH_ = nullptr;
0498   vertexYposH_ = nullptr;
0499   vertexZposH_ = nullptr;
0500 
0501   nPixBarrelH_ = nullptr;
0502   nPixEndcapH_ = nullptr;
0503   nStripTIBH_ = nullptr;
0504   nStripTOBH_ = nullptr;
0505   nStripTECH_ = nullptr;
0506   nStripTIDH_ = nullptr;
0507 
0508   // for MC only
0509   nVtxH_ = nullptr;
0510   nVertexH_ = nullptr;
0511   bunchCrossingH_ = nullptr;
0512   nPUH_ = nullptr;
0513   trueNIntH_ = nullptr;
0514 
0515   nLostHitsVspTH_ = nullptr;
0516   nLostHitsVsEtaH_ = nullptr;
0517   nLostHitsVsCosThetaH_ = nullptr;
0518   nLostHitsVsPhiH_ = nullptr;
0519   nLostHitsVsIterationH_ = nullptr;
0520 
0521   nHitsTIBSVsEtaH_ = nullptr;
0522   nHitsTOBSVsEtaH_ = nullptr;
0523   nHitsTECSVsEtaH_ = nullptr;
0524   nHitsTIDSVsEtaH_ = nullptr;
0525   nHitsStripSVsEtaH_ = nullptr;
0526 
0527   nHitsTIBDVsEtaH_ = nullptr;
0528   nHitsTOBDVsEtaH_ = nullptr;
0529   nHitsTECDVsEtaH_ = nullptr;
0530   nHitsTIDDVsEtaH_ = nullptr;
0531   nHitsStripDVsEtaH_ = nullptr;
0532 
0533   nValidHitsVspTH_ = nullptr;
0534   nValidHitsVsEtaH_ = nullptr;
0535   nValidHitsVsCosThetaH_ = nullptr;
0536   nValidHitsVsPhiH_ = nullptr;
0537   nValidHitsVsnVtxH_ = nullptr;
0538 
0539   nValidHitsPixVsEtaH_ = nullptr;
0540   nValidHitsPixBVsEtaH_ = nullptr;
0541   nValidHitsPixEVsEtaH_ = nullptr;
0542   nValidHitsStripVsEtaH_ = nullptr;
0543   nValidHitsTIBVsEtaH_ = nullptr;
0544   nValidHitsTOBVsEtaH_ = nullptr;
0545   nValidHitsTECVsEtaH_ = nullptr;
0546   nValidHitsTIDVsEtaH_ = nullptr;
0547 
0548   nValidHitsPixVsPhiH_ = nullptr;
0549   nValidHitsPixBVsPhiH_ = nullptr;
0550   nValidHitsPixEVsPhiH_ = nullptr;
0551   nValidHitsStripVsPhiH_ = nullptr;
0552   nValidHitsTIBVsPhiH_ = nullptr;
0553   nValidHitsTOBVsPhiH_ = nullptr;
0554   nValidHitsTECVsPhiH_ = nullptr;
0555   nValidHitsTIDVsPhiH_ = nullptr;
0556 
0557   nLostHitsPixVsEtaH_ = nullptr;
0558   nLostHitsPixBVsEtaH_ = nullptr;
0559   nLostHitsPixEVsEtaH_ = nullptr;
0560   nLostHitsStripVsEtaH_ = nullptr;
0561   nLostHitsTIBVsEtaH_ = nullptr;
0562   nLostHitsTOBVsEtaH_ = nullptr;
0563   nLostHitsTECVsEtaH_ = nullptr;
0564   nLostHitsTIDVsEtaH_ = nullptr;
0565 
0566   nLostHitsPixVsPhiH_ = nullptr;
0567   nLostHitsPixBVsPhiH_ = nullptr;
0568   nLostHitsPixEVsPhiH_ = nullptr;
0569   nLostHitsStripVsPhiH_ = nullptr;
0570   nLostHitsTIBVsPhiH_ = nullptr;
0571   nLostHitsTOBVsPhiH_ = nullptr;
0572   nLostHitsTECVsPhiH_ = nullptr;
0573   nLostHitsTIDVsPhiH_ = nullptr;
0574 
0575   nLostHitsPixVsIterationH_ = nullptr;
0576   nLostHitsPixBVsIterationH_ = nullptr;
0577   nLostHitsPixEVsIterationH_ = nullptr;
0578   nLostHitsStripVsIterationH_ = nullptr;
0579   nLostHitsTIBVsIterationH_ = nullptr;
0580   nLostHitsTOBVsIterationH_ = nullptr;
0581   nLostHitsTECVsIterationH_ = nullptr;
0582   nLostHitsTIDVsIterationH_ = nullptr;
0583 
0584   hOnTrkClusChargeThinH_ = nullptr;
0585   hOnTrkClusWidthThinH_ = nullptr;
0586   hOnTrkClusChargeThickH_ = nullptr;
0587   hOnTrkClusWidthThickH_ = nullptr;
0588 
0589   hOffTrkClusChargeThinH_ = nullptr;
0590   hOffTrkClusWidthThinH_ = nullptr;
0591   hOffTrkClusChargeThickH_ = nullptr;
0592   hOffTrkClusWidthThickH_ = nullptr;
0593 
0594   cosPhi3DdileptonH_ = nullptr;
0595   ip3dToPV2validpixelhitsH_ = nullptr;
0596   ip3dToBS2validpixelhitsH_ = nullptr;
0597   iperr3dToPV2validpixelhitsH_ = nullptr;
0598   iperr3dToBS2validpixelhitsH_ = nullptr;
0599   sip3dToPV2validpixelhitsH_ = nullptr;
0600   sip3dToBS2validpixelhitsH_ = nullptr;
0601 
0602   // Read pileup weight factors
0603 
0604   if (isMC_ && doPUCorrection_ && doTrackCorrection_) {
0605     throw std::runtime_error("if isMC is true, only one of doPUCorrection and doTrackCorrection can be true");
0606   }
0607 
0608   if (isMC_ && doPUCorrection_) {
0609     vpu_.clear();
0610     TFile* f1 = TFile::Open(puScaleFactorFile_.c_str());
0611     TH1F* h1 = dynamic_cast<TH1F*>(f1->Get("pileupweight"));
0612     for (int i = 1; i <= h1->GetNbinsX(); ++i)
0613       vpu_.push_back(h1->GetBinContent(i));
0614     f1->Close();
0615   }
0616 
0617   if (isMC_ && doTrackCorrection_) {
0618     vtrack_.clear();
0619     TFile* f1 = TFile::Open(trackScaleFactorFile_.c_str());
0620     TH1F* h1 = dynamic_cast<TH1F*>(f1->Get("trackweight"));
0621     for (int i = 1; i <= h1->GetNbinsX(); ++i)
0622       vtrack_.push_back(h1->GetBinContent(i));
0623     f1->Close();
0624   }
0625 }
0626 
0627 void StandaloneTrackMonitor::dqmBeginRun(const edm::Run& iRun, const edm::EventSetup& iSetup) {
0628   tkGeom_ = &(iSetup.getData(geomToken_));
0629 }
0630 
0631 void StandaloneTrackMonitor::bookHistograms(DQMStore::IBooker& ibook,
0632                                             edm::Run const& iRun,
0633                                             edm::EventSetup const& iSetup) {
0634   std::string currentFolder = moduleName_ + "/" + folderName_;
0635   ibook.setCurrentFolder(currentFolder);
0636 
0637   // The following are common with the official tool
0638   if (haveAllHistograms_) {
0639     trackEtaH_ = ibook.book1DD("trackEta",
0640                                "Track Eta",
0641                                TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0642                                TrackEtaHistoPar_.getParameter<double>("Xmin"),
0643                                TrackEtaHistoPar_.getParameter<double>("Xmax"));
0644 
0645     trackEtaerrH_ = ibook.book1DD("trackEtaerr", "Track Eta Error", 50, 0.0, 1.0);
0646     trackCosThetaH_ = ibook.book1DD("trackCosTheta", "Track Cos(Theta)", 50, -1.0, 1.0);
0647     trackThetaerrH_ = ibook.book1DD("trackThetaerr", "Track Theta Error", 50, 0.0, 1.0);
0648     trackPhiH_ = ibook.book1DD("trackPhi", "Track Phi", 70, -3.5, 3.5);
0649     trackPhierrH_ = ibook.book1DD("trackPhierr", "Track Phi Error", 50, 0.0, 1.0);
0650 
0651     trackPH_ = ibook.book1DD("trackP", "Track 4-momentum", 50, 0.0, 10.0);
0652     trackPtH_ = ibook.book1DD("trackPt",
0653                               "Track Pt",
0654                               TrackPtHistoPar_.getParameter<int32_t>("Xbins"),
0655                               TrackPtHistoPar_.getParameter<double>("Xmin"),
0656                               TrackPtHistoPar_.getParameter<double>("Xmax"));
0657     trackPt_ZoomH_ = ibook.book1DD("trackPt_Zoom", "Track Pt", 100, 60, 70);
0658     trackPterrH_ = ibook.book1DD("trackPterr", "Track Pt Error", 100, 0.0, 100.0);
0659     trackqOverpH_ = ibook.book1DD("trackqOverp", "q Over p", 100, -0.5, 0.5);
0660     trackqOverperrH_ = ibook.book1DD("trackqOverperr", "q Over p Error", 50, 0.0, 0.5);
0661     trackChargeH_ = ibook.book1DD("trackCharge", "Track Charge", 3, -1.5, 1.5);
0662     trackChi2H_ = ibook.book1DD("trackChi2", "Chi2", 100, 0.0, 100.0);
0663     tracknDOFH_ = ibook.book1DD("tracknDOF", "nDOF", 100, 0.0, 100.0);
0664     trackChi2ProbH_ = ibook.book1DD("trackChi2Prob", "Chi2prob", 50, 0.0, 1.0);
0665     trackChi2oNDFH_ = ibook.book1DD("trackChi2oNDF", "Chi2oNDF", 100, 0.0, 100.0);
0666     trackd0H_ = ibook.book1DD("trackd0", "Track d0", 100, -1, 1);
0667     trackChi2bynDOFH_ = ibook.book1DD("trackChi2bynDOF", "Chi2 Over nDOF", 100, 0.0, 10.0);
0668     trackalgoH_ = ibook.book1DD("trackalgo", "Track Algo", 46, 0.0, 46.0);
0669     trackorigalgoH_ = ibook.book1DD("trackorigalgo", "Track Original Algo", 46, 0.0, 46.0);
0670 
0671     for (size_t ibin = 0; ibin < reco::TrackBase::algoSize - 1; ibin++) {
0672       trackalgoH_->setBinLabel(ibin + 1, reco::TrackBase::algoNames[ibin], 1);
0673       trackorigalgoH_->setBinLabel(ibin + 1, reco::TrackBase::algoNames[ibin], 1);
0674     }
0675 
0676     trackStoppingSourceH_ = ibook.book1DD("trackstoppingsource", "Track Stopping Source", 12, 0.0, 12.0);
0677 
0678     // DataFormats/TrackCandidate/interface/TrajectoryStopReasons.h
0679     size_t StopReasonNameSize = static_cast<size_t>(StopReason::SIZE);
0680     for (unsigned int i = 0; i < StopReasonNameSize; ++i) {
0681       trackStoppingSourceH_->setBinLabel(i + 1, StopReasonName::StopReasonName[i], 1);
0682     }
0683 
0684     DistanceOfClosestApproachToPVH_ =
0685         ibook.book1DD("DistanceOfClosestApproachToPV", "DistanceOfClosestApproachToPV", 1000, -1.0, 1.0);
0686     DistanceOfClosestApproachToPVZoomedH_ =
0687         ibook.book1DD("DistanceOfClosestApproachToPVZoomed", "DistanceOfClosestApproachToPV", 1000, -0.1, 0.1);
0688     DistanceOfClosestApproachToPVVsPhiH_ = ibook.bookProfile(
0689         "DistanceOfClosestApproachToPVVsPhi", "DistanceOfClosestApproachToPVVsPhi", 100, -3.5, 3.5, 0.0, 0.0, "g");
0690     xPointOfClosestApproachVsZ0wrtPVH_ = ibook.bookProfile(
0691         "xPointOfClosestApproachVsZ0wrtPV", "xPointOfClosestApproachVsZ0wrtPV", 120, -60, 60, 0.0, 0.0, "g");
0692     yPointOfClosestApproachVsZ0wrtPVH_ = ibook.bookProfile(
0693         "yPointOfClosestApproachVsZ0wrtPV", "yPointOfClosestApproachVsZ0wrtPV", 120, -60, 60, 0.0, 0.0, "g");
0694     trackDeltaRwrtClosestTrack_ =
0695         ibook.book1DD("trackDeltaRwrtClosestTrack", "min#DeltaR(considered track,other tracks)", 500, 0, 10);
0696 
0697     ip2dToPVH_ = ibook.book1DD("ip2dToPV", "IP in 2d To PV", 1000, -1.0, 1.0);
0698     unsigned int niperrbins = 100;
0699     float iperrbinning[niperrbins + 1];
0700     float miniperr = 0.0001, maxiperr = 5;
0701     iperrbinning[0] = miniperr;
0702     for (unsigned int i = 1; i != niperrbins + 1; i++) {
0703       iperrbinning[i] = iperrbinning[i - 1] * pow(maxiperr / miniperr, 1. / niperrbins);
0704     }
0705     iperr2dToPVH_ = ibook.book1DD("iperr2dToPV", "IP error in 2d To PV", niperrbins, iperrbinning);
0706     iperr2dToPVWtH_ = ibook.book1DD("iperr2dToPVWt", "IP error in 2d To PV", niperrbins, iperrbinning);
0707 
0708     ip3dToPVH_ = ibook.book1DD("ip3dToPV", "IP in 3d To PV", 200, -20, 20);
0709     ip3dToBSH_ = ibook.book1DD("ip3dToBS", "IP in 3d To BS", 200, -20, 20);
0710     iperr3dToPVH_ = ibook.book1DD("iperr3dToPV", "IP error in 3d To PV", niperrbins, iperrbinning);
0711     iperr3dToBSH_ = ibook.book1DD("iperr3dToBS", "IP error in 3d To BS", niperrbins, iperrbinning);
0712     sip3dToPVH_ = ibook.book1DD("sip3dToPV", "IP significance in 3d To PV", 200, -10, 10);
0713     sip3dToBSH_ = ibook.book1DD("sip3dToBS", "IP significance in 3d To BS", 200, -10, 10);
0714 
0715     ip3dToPV2validpixelhitsH_ =
0716         ibook.book1DD("ip3dToPV2validpixelhits", "IP in 3d To PV (nValidPixelHits>2)", 200, -0.20, 0.20);
0717     ip3dToBS2validpixelhitsH_ =
0718         ibook.book1DD("ip3dToBS2validpixelhits", "IP in 3d To BS (nValidPixelHits>2)", 200, -0.20, 0.20);
0719     iperr3dToPV2validpixelhitsH_ = ibook.book1DD(
0720         "iperr3dToPV2validpixelhits", "IP error in 3d To PV (nValidPixelHits>2)", niperrbins, iperrbinning);
0721     iperr3dToBS2validpixelhitsH_ = ibook.book1DD(
0722         "iperr3dToBS2validpixelhits", "IP error in 3d To BS (nValidPixelHits>2)", niperrbins, iperrbinning);
0723     sip3dToPV2validpixelhitsH_ =
0724         ibook.book1DD("sip3dToPV2validpixelhits", "IP significance in 3d To PV (nValidPixelHits>2)", 200, -10, 10);
0725     sip3dToBS2validpixelhitsH_ =
0726         ibook.book1DD("sip3dToBS2validpixelhits", "IP significance in 3d To BS (nValidPixelHits>2)", 200, -10, 10);
0727 
0728     ip2dToBSH_ = ibook.book1DD("ip2dToBS", "IP in 2d To BS", 1000, -1., 1.);  //Beamspot
0729     iperr2dToBSH_ = ibook.book1DD("iperr2dToBS", "IP error in 2d To BS", niperrbins, iperrbinning);
0730     sip2dToBSH_ = ibook.book1DD("sip2dToBS", "IP significance in 2d To BS", 200, -10, 10);
0731 
0732     iperr3dToPVWtH_ = ibook.book1DD("iperr3dToPVWt", "IP error in 3d To PV", niperrbins, iperrbinning);
0733     sip2dToPVH_ = ibook.book1DD("sip2dToPV", "IP significance in 2d To PV", 200, -10, 10);
0734 
0735     sip2dToPVWtH_ = ibook.book1DD("sip2dToPVWt", "IP significance in 2d To PV", 200, -10, 10);
0736     sipDxyToPVH_ = ibook.book1DD("sipDxyToPV", "IP significance in dxy To PV", 100, -10, 10);
0737     sipDzToPVH_ = ibook.book1DD("sipDzToPV", "IP significance in dz To PV", 100, -10, 10);
0738 
0739     nallHitsH_ = ibook.book1DD("nallHits", "No. of All Hits", 60, -0.5, 59.5);
0740     ntrackerHitsH_ = ibook.book1DD("ntrackerHits", "No. of Tracker Hits", 60, -0.5, 59.5);
0741 
0742     nvalidTrackerHitsH_ = ibook.book1DD("nvalidTrackerhits", "No. of Valid Tracker Hits", 47, -0.5, 46.5);
0743     nvalidPixelHitsH_ = ibook.book1DD("nvalidPixelHits", "No. of Valid Hits in Pixel", 8, -0.5, 7.5);
0744     nvalidPixelBHitsH_ = ibook.book1DD("nvalidPixelBarrelHits", "No. of Valid Hits in Pixel Barrel", 6, -0.5, 5.5);
0745     nvalidPixelEHitsH_ = ibook.book1DD("nvalidPixelEndcapHits", "No. of Valid Hits in Pixel Endcap", 6, -0.5, 6.5);
0746     nvalidStripHitsH_ = ibook.book1DD("nvalidStripHits", "No. of Valid Hits in Strip", 36, -0.5, 35.5);
0747     nvalidTIBHitsH_ = ibook.book1DD("nvalidTIBHits", "No. of Valid Hits in Strip TIB", 6, -0.5, 5.5);
0748     nvalidTOBHitsH_ = ibook.book1DD("nvalidTOBHits", "No. of Valid Hits in Strip TOB", 11, -0.5, 10.5);
0749     nvalidTIDHitsH_ = ibook.book1DD("nvalidTIDHits", "No. of Valid Hits in Strip TID", 6, -0.5, 5.5);
0750     nvalidTECHitsH_ = ibook.book1DD("nvalidTECHits", "No. of Valid Hits in Strip TEC", 11, -0.5, 10.5);
0751 
0752     nlostTrackerHitsH_ = ibook.book1DD("nlostTrackerhits", "No. of Lost Tracker Hits", 15, -0.5, 14.5);
0753     nlostPixelHitsH_ = ibook.book1DD("nlostPixelHits", "No. of Lost Hits in Pixel", 8, -0.5, 7.5);
0754     nlostPixelBHitsH_ = ibook.book1DD("nlostPixelBarrelHits", "No. of Lost Hits in Pixel Barrel", 5, -0.5, 4.5);
0755     nlostPixelEHitsH_ = ibook.book1DD("nlostPixelEndcapHits", "No. of Lost Hits in Pixel Endcap", 4, -0.5, 3.5);
0756     nlostStripHitsH_ = ibook.book1DD("nlostStripHits", "No. of Lost Hits in Strip", 10, -0.5, 9.5);
0757     nlostTIBHitsH_ = ibook.book1DD("nlostTIBHits", "No. of Lost Hits in Strip TIB", 5, -0.5, 4.5);
0758     nlostTOBHitsH_ = ibook.book1DD("nlostTOBHits", "No. of Lost Hits in Strip TOB", 10, -0.5, 9.5);
0759     nlostTIDHitsH_ = ibook.book1DD("nlostTIDHits", "No. of Lost Hits in Strip TID", 5, -0.5, 4.5);
0760     nlostTECHitsH_ = ibook.book1DD("nlostTECHits", "No. of Lost Hits in Strip TEC", 10, -0.5, 9.5);
0761 
0762     trkLayerwithMeasurementH_ = ibook.book1DD("trkLayerwithMeasurement", "No. of Layers per Track", 20, 0.0, 20.0);
0763     pixelLayerwithMeasurementH_ =
0764         ibook.book1DD("pixelLayerwithMeasurement", "No. of Pixel Layers per Track", 10, 0.0, 10.0);
0765     pixelBLayerwithMeasurementH_ =
0766         ibook.book1DD("pixelBLayerwithMeasurement", "No. of Pixel Barrel Layers per Track", 5, 0.0, 5.0);
0767     pixelELayerwithMeasurementH_ =
0768         ibook.book1DD("pixelELayerwithMeasurement", "No. of Pixel Endcap Layers per Track", 5, 0.0, 5.0);
0769     stripLayerwithMeasurementH_ =
0770         ibook.book1DD("stripLayerwithMeasurement", "No. of Strip Layers per Track", 20, 0.0, 20.0);
0771     stripTIBLayerwithMeasurementH_ =
0772         ibook.book1DD("stripTIBLayerwithMeasurement", "No. of Strip TIB Layers per Track", 10, 0.0, 10.0);
0773     stripTOBLayerwithMeasurementH_ =
0774         ibook.book1DD("stripTOBLayerwithMeasurement", "No. of Strip TOB Layers per Track", 10, 0.0, 10.0);
0775     stripTIDLayerwithMeasurementH_ =
0776         ibook.book1DD("stripTIDLayerwithMeasurement", "No. of Strip TID Layers per Track", 5, 0.0, 5.0);
0777     stripTECLayerwithMeasurementH_ =
0778         ibook.book1DD("stripTECLayerwithMeasurement", "No. of Strip TEC Layers per Track", 15, 0.0, 15.0);
0779 
0780     nlostHitsH_ = ibook.book1DD("nlostHits", "No. of Lost Hits", 10, -0.5, 9.5);
0781     nMissingExpectedInnerHitsH_ =
0782         ibook.book1DD("nMissingExpectedInnerHits", "No. of Missing Expected Inner Hits", 10, -0.5, 9.5);
0783     nMissingExpectedOuterHitsH_ =
0784         ibook.book1DD("nMissingExpectedOuterHits", "No. of Missing Expected Outer Hits", 10, -0.5, 9.5);
0785 
0786     beamSpotXYposH_ = ibook.book1DD("beamSpotXYpos", "d_{xy} w.r.t beam spot", 100, -1.0, 1.0);
0787     beamSpotXYposerrH_ = ibook.book1DD("beamSpotXYposerr", "error on d_{xy} w.r.t beam spot", 50, 0.0, 0.25);
0788     beamSpotZposH_ = ibook.book1DD("beamSpotZpos", "d_{z} w.r.t beam spot", 100, -20.0, 20.0);
0789     beamSpotZposerrH_ = ibook.book1DD("beamSpotZposerr", "error on d_{z} w.r.t. beam spot", 50, 0.0, 0.25);
0790 
0791     vertexXposH_ = ibook.book1DD("vertexXpos", "Vertex X position", 100, -0.1, 0.1);
0792     vertexYposH_ = ibook.book1DD("vertexYpos", "Vertex Y position", 200, -0.1, 0.1);
0793     vertexZposH_ = ibook.book1DD("vertexZpos", "Vertex Z position", 100, -20.0, 20.0);
0794     nVertexH_ = ibook.book1DD("nVertex", "# of vertices", 120, -0.5, 119.5);
0795     nVtxH_ = ibook.book1DD("nVtx", "# of vtxs", 120, -0.5, 119.5);
0796 
0797     nMissingInnerHitBH_ = ibook.book1DD("nMissingInnerHitB", "No. missing inner hit per Track in Barrel", 6, -0.5, 5.5);
0798     nMissingInnerHitEH_ = ibook.book1DD("nMissingInnerHitE", "No. missing inner hit per Track in Endcap", 6, -0.5, 5.5);
0799     nMissingOuterHitBH_ =
0800         ibook.book1DD("nMissingOuterHitB", "No. missing outer hit per Track in Barrel", 11, -0.5, 10.5);
0801     nMissingOuterHitEH_ =
0802         ibook.book1DD("nMissingOuterHitE", "No. missing outer hit per Track in Endcap", 11, -0.5, 10.5);
0803     nPixBarrelH_ = ibook.book1DD("nHitPixelBarrel", "No. of hits in Pixel Barrel per Track", 20, 0, 20.0);
0804     nPixEndcapH_ = ibook.book1DD("nHitPixelEndcap", "No. of hits in Pixel Endcap per Track", 20, 0, 20.0);
0805     nStripTIBH_ = ibook.book1DD("nHitStripTIB", "No. of hits in Strip TIB per Track", 30, 0, 30.0);
0806     nStripTOBH_ = ibook.book1DD("nHitStripTOB", "No. of hits in Strip TOB per Track", 30, 0, 30.0);
0807     nStripTECH_ = ibook.book1DD("nHitStripTEC", "No. of hits in Strip TEC per Track", 30, 0, 30.0);
0808     nStripTIDH_ = ibook.book1DD("nHitStripTID", "No. of hits in Strip TID per Tracks", 30, 0, 30.0);
0809     nTracksH_ = ibook.book1DD("nTracks", "No. of Tracks", 1200, -0.5, 1199.5);
0810     nJet_ = ibook.book1DD("nJet", "Number of Jets", 101, -0.5, 100.5);
0811     Jet_pt_ = ibook.book1DD("Jet_pt", "Jet p_{T}", 200, 0., 200.);
0812     Jet_eta_ = ibook.book1DD("Jet_eta", "Jet #eta", 100, -5.2, 5.2);
0813     Jet_energy_ = ibook.book1DD("Jet_energy", "Jet Energy", 200, 0., 200.);
0814     Jet_chargedMultiplicity_ =
0815         ibook.book1DD("Jet_chargedMultiplicity", "Jet charged Hadron Multiplicity", 201, -0.5, 200.5);
0816     Zpt_ = ibook.book1DD("Zpt", "Z-boson transverse momentum", 100, 0, 100);
0817     ZInvMass_ = ibook.book1DD("ZInvMass", "m_{ll}", 120, 75, 105);
0818     cosPhi3DdileptonH_ = ibook.book1DD("cosPhi3Ddilepton", "cos#Phi_{3D,ll}", 202, -1.01, 1.01);
0819   }
0820   if (isMC_) {
0821     bunchCrossingH_ = ibook.book1DD("bunchCrossing", "Bunch Crossing", 60, 0, 60.0);
0822     nPUH_ = ibook.book1DD("nPU", "No of Pileup", 100, 0, 100.0);
0823     trueNIntH_ = ibook.book1DD("trueNInt", "True no of Interactions", 100, 0, 100.0);
0824   }
0825   // Exclusive histograms
0826 
0827   nLostHitByLayerH_ = ibook.book1DD("nLostHitByLayer", "No. of Lost Hit per Layer", 29, 0.5, 29.5);
0828 
0829   nLostHitByLayerPixH_ =
0830       ibook.book1DD("nLostHitByLayerPix", "No. of Lost Hit per Layer for Pixel detector", 7, 0.5, 7.5);
0831 
0832   nLostHitByLayerStripH_ =
0833       ibook.book1DD("nLostHitByLayerStrip", "No. of Lost Hit per Layer for SiStrip detector", 22, 0.5, 22.5);
0834 
0835   nLostHitsVspTH_ = ibook.bookProfile("nLostHitsVspT",
0836                                       "Number of Lost Hits Vs pT",
0837                                       TrackPtHistoPar_.getParameter<int32_t>("Xbins"),
0838                                       TrackPtHistoPar_.getParameter<double>("Xmin"),
0839                                       TrackPtHistoPar_.getParameter<double>("Xmax"),
0840                                       0.0,
0841                                       0.0,
0842                                       "g");
0843   nLostHitsVsEtaH_ = ibook.bookProfile("nLostHitsVsEta",
0844                                        "Number of Lost Hits Vs Eta",
0845                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0846                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0847                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0848                                        0.0,
0849                                        0.0,
0850                                        "g");
0851   nLostHitsVsCosThetaH_ =
0852       ibook.bookProfile("nLostHitsVsCosTheta", "Number of Lost Hits Vs Cos(Theta)", 50, -1.0, 1.0, 0.0, 0.0, "g");
0853   nLostHitsVsPhiH_ = ibook.bookProfile("nLostHitsVsPhi", "Number of Lost Hits Vs Phi", 100, -3.5, 3.5, 0.0, 0.0, "g");
0854   nLostHitsVsIterationH_ =
0855       ibook.bookProfile("nLostHitsVsIteration", "Number of Lost Hits Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
0856 
0857   nHitsTIBSVsEtaH_ = ibook.bookProfile("nHitsTIBSVsEta",
0858                                        "Number of Hits in TIB Vs Eta (Single-sided)",
0859                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0860                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0861                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0862                                        0.0,
0863                                        0.0,
0864                                        "g");
0865   nHitsTOBSVsEtaH_ = ibook.bookProfile("nHitsTOBSVsEta",
0866                                        "Number of Hits in TOB Vs Eta (Single-sided)",
0867                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0868                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0869                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0870                                        0.0,
0871                                        0.0,
0872                                        "g");
0873   nHitsTECSVsEtaH_ = ibook.bookProfile("nHitsTECSVsEta",
0874                                        "Number of Hits in TEC Vs Eta (Single-sided)",
0875                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0876                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0877                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0878                                        0.0,
0879                                        0.0,
0880                                        "g");
0881   nHitsTIDSVsEtaH_ = ibook.bookProfile("nHitsTIDSVsEta",
0882                                        "Number of Hits in TID Vs Eta (Single-sided)",
0883                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0884                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0885                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0886                                        0.0,
0887                                        0.0,
0888                                        "g");
0889 
0890   nHitsStripSVsEtaH_ = ibook.bookProfile("nHitsStripSVsEta",
0891                                          "Number of Strip Hits Vs Eta (Single-sided)",
0892                                          TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0893                                          TrackEtaHistoPar_.getParameter<double>("Xmin"),
0894                                          TrackEtaHistoPar_.getParameter<double>("Xmax"),
0895                                          0.0,
0896                                          0.0,
0897                                          "g");
0898 
0899   nHitsTIBDVsEtaH_ = ibook.bookProfile("nHitsTIBDVsEta",
0900                                        "Number of Hits in TIB Vs Eta (Double-sided)",
0901                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0902                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0903                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0904                                        0.0,
0905                                        0.0,
0906                                        "g");
0907   nHitsTOBDVsEtaH_ = ibook.bookProfile("nHitsTOBDVsEta",
0908                                        "Number of Hits in TOB Vs Eta (Double-sided)",
0909                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0910                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0911                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0912                                        0.0,
0913                                        0.0,
0914                                        "g");
0915   nHitsTECDVsEtaH_ = ibook.bookProfile("nHitsTECDVsEta",
0916                                        "Number of Hits in TEC Vs Eta (Double-sided)",
0917                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0918                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0919                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0920                                        0.0,
0921                                        0.0,
0922                                        "g");
0923   nHitsTIDDVsEtaH_ = ibook.bookProfile("nHitsTIDDVsEta",
0924                                        "Number of Hits in TID Vs Eta (Double-sided)",
0925                                        TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0926                                        TrackEtaHistoPar_.getParameter<double>("Xmin"),
0927                                        TrackEtaHistoPar_.getParameter<double>("Xmax"),
0928                                        0.0,
0929                                        0.0,
0930                                        "g");
0931   nHitsStripDVsEtaH_ = ibook.bookProfile("nHitsStripDVsEta",
0932                                          "Number of Strip Hits Vs Eta (Double-sided)",
0933                                          TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0934                                          TrackEtaHistoPar_.getParameter<double>("Xmin"),
0935                                          TrackEtaHistoPar_.getParameter<double>("Xmax"),
0936                                          0.0,
0937                                          0.0,
0938                                          "g");
0939 
0940   nValidHitsVspTH_ = ibook.bookProfile("nValidHitsVspT",
0941                                        "Number of Valid Hits Vs pT",
0942                                        TrackPtHistoPar_.getParameter<int32_t>("Xbins"),
0943                                        TrackPtHistoPar_.getParameter<double>("Xmin"),
0944                                        TrackPtHistoPar_.getParameter<double>("Xmax"),
0945                                        0.0,
0946                                        0.0,
0947                                        "g");
0948   nValidHitsVsnVtxH_ =
0949       ibook.bookProfile("nValidHitsVsnVtx", "Number of Valid Hits Vs Number of Vertex", 100, 0., 100., 0.0, 0.0, "g");
0950   nValidHitsVsEtaH_ = ibook.bookProfile("nValidHitsVsEta",
0951                                         "Number of Hits Vs Eta",
0952                                         TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0953                                         TrackEtaHistoPar_.getParameter<double>("Xmin"),
0954                                         TrackEtaHistoPar_.getParameter<double>("Xmax"),
0955                                         0.0,
0956                                         0.0,
0957                                         "g");
0958 
0959   nValidHitsVsCosThetaH_ =
0960       ibook.bookProfile("nValidHitsVsCosTheta", "Number of Valid Hits Vs Cos(Theta)", 50, -1.0, 1.0, 0.0, 0.0, "g");
0961   nValidHitsVsPhiH_ =
0962       ibook.bookProfile("nValidHitsVsPhi", "Number of Valid Hits Vs Phi", 100, -3.5, 3.5, 0.0, 0.0, "g");
0963 
0964   nValidHitsPixVsEtaH_ = ibook.bookProfile("nValidHitsPixVsEta",
0965                                            "Number of Valid Hits in Pixel Vs Eta",
0966                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0967                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
0968                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
0969                                            0.0,
0970                                            0.0,
0971                                            "g");
0972   nValidHitsPixBVsEtaH_ = ibook.bookProfile("nValidHitsPixBVsEta",
0973                                             "Number of Valid Hits in Pixel Barrel Vs Eta",
0974                                             TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0975                                             TrackEtaHistoPar_.getParameter<double>("Xmin"),
0976                                             TrackEtaHistoPar_.getParameter<double>("Xmax"),
0977                                             0.0,
0978                                             0.0,
0979                                             "g");
0980   nValidHitsPixEVsEtaH_ = ibook.bookProfile("nValidHitsPixEVsEta",
0981                                             "Number of Valid Hits in Pixel Endcap Vs Eta",
0982                                             TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0983                                             TrackEtaHistoPar_.getParameter<double>("Xmin"),
0984                                             TrackEtaHistoPar_.getParameter<double>("Xmax"),
0985                                             0.0,
0986                                             0.0,
0987                                             "g");
0988   nValidHitsStripVsEtaH_ = ibook.bookProfile("nValidHitsStripVsEta",
0989                                              "Number of Valid Hits in SiStrip Vs Eta",
0990                                              TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0991                                              TrackEtaHistoPar_.getParameter<double>("Xmin"),
0992                                              TrackEtaHistoPar_.getParameter<double>("Xmax"),
0993                                              0.0,
0994                                              0.0,
0995                                              "g");
0996   nValidHitsTIBVsEtaH_ = ibook.bookProfile("nValidHitsTIBVsEta",
0997                                            "Number of Valid Hits in TIB Vs Eta",
0998                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
0999                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1000                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1001                                            0.0,
1002                                            0.0,
1003                                            "g");
1004   nValidHitsTOBVsEtaH_ = ibook.bookProfile("nValidHitsTOBVsEta",
1005                                            "Number of Valid Hits in TOB Vs Eta",
1006                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1007                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1008                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1009                                            0.0,
1010                                            0.0,
1011                                            "g");
1012   nValidHitsTECVsEtaH_ = ibook.bookProfile("nValidHitsTECVsEta",
1013                                            "Number of Valid Hits in TEC Vs Eta",
1014                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1015                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1016                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1017                                            0.0,
1018                                            0.0,
1019                                            "g");
1020   nValidHitsTIDVsEtaH_ = ibook.bookProfile("nValidHitsTIDVsEta",
1021                                            "Number of Valid Hits in TID Vs Eta",
1022                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1023                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1024                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1025                                            0.0,
1026                                            0.0,
1027                                            "g");
1028 
1029   nValidHitsPixVsPhiH_ = ibook.bookProfile("nValidHitsPixVsPhi",
1030                                            "Number of Valid Hits in Pixel Vs Phi",
1031                                            TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1032                                            TrackPhiHistoPar_.getParameter<double>("Xmin"),
1033                                            TrackPhiHistoPar_.getParameter<double>("Xmax"),
1034                                            0.0,
1035                                            0.0,
1036                                            "g");
1037   nValidHitsPixBVsPhiH_ = ibook.bookProfile("nValidHitsPixBVsPhi",
1038                                             "Number of Valid Hits in Pixel Barrel Vs Phi",
1039                                             TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1040                                             TrackPhiHistoPar_.getParameter<double>("Xmin"),
1041                                             TrackPhiHistoPar_.getParameter<double>("Xmax"),
1042                                             0.0,
1043                                             0.0,
1044                                             "g");
1045   nValidHitsPixEVsPhiH_ = ibook.bookProfile("nValidHitsPixEVsPhi",
1046                                             "Number of Valid Hits in Pixel Endcap Vs Phi",
1047                                             TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1048                                             TrackEtaHistoPar_.getParameter<double>("Xmin"),
1049                                             TrackEtaHistoPar_.getParameter<double>("Xmax"),
1050                                             0.0,
1051                                             0.0,
1052                                             "g");
1053   nValidHitsStripVsPhiH_ = ibook.bookProfile("nValidHitsStripVsPhi",
1054                                              "Number of Valid Hits in SiStrip Vs Phi",
1055                                              TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1056                                              TrackPhiHistoPar_.getParameter<double>("Xmin"),
1057                                              TrackPhiHistoPar_.getParameter<double>("Xmax"),
1058                                              0.0,
1059                                              0.0,
1060                                              "g");
1061   nValidHitsTIBVsPhiH_ = ibook.bookProfile("nValidHitsTIBVsPhi",
1062                                            "Number of Valid Hits in TIB Vs Phi",
1063                                            TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1064                                            TrackPhiHistoPar_.getParameter<double>("Xmin"),
1065                                            TrackPhiHistoPar_.getParameter<double>("Xmax"),
1066                                            0.0,
1067                                            0.0,
1068                                            "g");
1069   nValidHitsTOBVsPhiH_ = ibook.bookProfile("nValidHitsTOBVsPhi",
1070                                            "Number of Valid Hits in TOB Vs Phi",
1071                                            TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1072                                            TrackPhiHistoPar_.getParameter<double>("Xmin"),
1073                                            TrackPhiHistoPar_.getParameter<double>("Xmax"),
1074                                            0.0,
1075                                            0.0,
1076                                            "g");
1077   nValidHitsTECVsPhiH_ = ibook.bookProfile("nValidHitsTECVsPhi",
1078                                            "Number of Valid Hits in TEC Vs Phi",
1079                                            TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1080                                            TrackPhiHistoPar_.getParameter<double>("Xmin"),
1081                                            TrackPhiHistoPar_.getParameter<double>("Xmax"),
1082                                            0.0,
1083                                            0.0,
1084                                            "g");
1085   nValidHitsTIDVsPhiH_ = ibook.bookProfile("nValidHitsTIDVsPhi",
1086                                            "Number of Valid Hits in TID Vs Phi",
1087                                            TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1088                                            TrackPhiHistoPar_.getParameter<double>("Xmin"),
1089                                            TrackPhiHistoPar_.getParameter<double>("Xmax"),
1090                                            0.0,
1091                                            0.0,
1092                                            "g");
1093 
1094   nLostHitsPixVsEtaH_ = ibook.bookProfile("nLostHitsPixVsEta",
1095                                           "Number of Lost Hits in Pixel Vs Eta",
1096                                           TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1097                                           TrackEtaHistoPar_.getParameter<double>("Xmin"),
1098                                           TrackEtaHistoPar_.getParameter<double>("Xmax"),
1099                                           0.0,
1100                                           0.0,
1101                                           "g");
1102   nLostHitsPixBVsEtaH_ = ibook.bookProfile("nLostHitsPixBVsEta",
1103                                            "Number of Lost Hits in Pixel Barrel Vs Eta",
1104                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1105                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1106                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1107                                            0.0,
1108                                            0.0,
1109                                            "g");
1110   nLostHitsPixEVsEtaH_ = ibook.bookProfile("nLostHitsPixEVsEta",
1111                                            "Number of Lost Hits in Pixel Endcap Vs Eta",
1112                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1113                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1114                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1115                                            0.0,
1116                                            0.0,
1117                                            "g");
1118   nLostHitsStripVsEtaH_ = ibook.bookProfile("nLostHitsStripVsEta",
1119                                             "Number of Lost Hits in SiStrip Vs Eta",
1120                                             TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1121                                             TrackEtaHistoPar_.getParameter<double>("Xmin"),
1122                                             TrackEtaHistoPar_.getParameter<double>("Xmax"),
1123                                             0.0,
1124                                             0.0,
1125                                             "g");
1126   nLostHitsTIBVsEtaH_ = ibook.bookProfile("nLostHitsTIBVsEta",
1127                                           "Number of Lost Hits in TIB Vs Eta",
1128                                           TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1129                                           TrackEtaHistoPar_.getParameter<double>("Xmin"),
1130                                           TrackEtaHistoPar_.getParameter<double>("Xmax"),
1131                                           0.0,
1132                                           0.0,
1133                                           "g");
1134   nLostHitsTOBVsEtaH_ = ibook.bookProfile("nLostHitsTOBVsEta",
1135                                           "Number of Lost Hits in TOB Vs Eta",
1136                                           TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1137                                           TrackEtaHistoPar_.getParameter<double>("Xmin"),
1138                                           TrackEtaHistoPar_.getParameter<double>("Xmax"),
1139                                           0.0,
1140                                           0.0,
1141                                           "g");
1142   nLostHitsTECVsEtaH_ = ibook.bookProfile("nLostHitsTECVsEta",
1143                                           "Number of Lost Hits in TEC Vs Eta",
1144                                           TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1145                                           TrackEtaHistoPar_.getParameter<double>("Xmin"),
1146                                           TrackEtaHistoPar_.getParameter<double>("Xmax"),
1147                                           0.0,
1148                                           0.0,
1149                                           "g");
1150   nLostHitsTIDVsEtaH_ = ibook.bookProfile("nLostHitsTIDVsEta",
1151                                           "Number of Lost Hits in TID Vs Eta",
1152                                           TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1153                                           TrackEtaHistoPar_.getParameter<double>("Xmin"),
1154                                           TrackEtaHistoPar_.getParameter<double>("Xmax"),
1155                                           0.0,
1156                                           0.0,
1157                                           "g");
1158 
1159   nLostHitsPixVsPhiH_ = ibook.bookProfile("nLostHitsPixVsPhi",
1160                                           "Number of Lost Hits in Pixel Vs Phi",
1161                                           TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1162                                           TrackPhiHistoPar_.getParameter<double>("Xmin"),
1163                                           TrackPhiHistoPar_.getParameter<double>("Xmax"),
1164                                           0.0,
1165                                           0.0,
1166                                           "g");
1167   nLostHitsPixBVsPhiH_ = ibook.bookProfile("nLostHitsPixBVsPhi",
1168                                            "Number of Lost Hits in Pixel Barrel Vs Phi",
1169                                            TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1170                                            TrackPhiHistoPar_.getParameter<double>("Xmin"),
1171                                            TrackPhiHistoPar_.getParameter<double>("Xmax"),
1172                                            0.0,
1173                                            0.0,
1174                                            "g");
1175   nLostHitsPixEVsPhiH_ = ibook.bookProfile("nLostHitsPixEVsPhi",
1176                                            "Number of Lost Hits in Pixel Endcap Vs Phi",
1177                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1178                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1179                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1180                                            0.0,
1181                                            0.0,
1182                                            "g");
1183   nLostHitsStripVsPhiH_ = ibook.bookProfile("nLostHitsStripVsPhi",
1184                                             "Number of Lost Hits in SiStrip Vs Phi",
1185                                             TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1186                                             TrackPhiHistoPar_.getParameter<double>("Xmin"),
1187                                             TrackPhiHistoPar_.getParameter<double>("Xmax"),
1188                                             0.0,
1189                                             0.0,
1190                                             "g");
1191   nLostHitsTIBVsPhiH_ = ibook.bookProfile("nLostHitsTIBVsPhi",
1192                                           "Number of Lost Hits in TIB Vs Phi",
1193                                           TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1194                                           TrackPhiHistoPar_.getParameter<double>("Xmin"),
1195                                           TrackPhiHistoPar_.getParameter<double>("Xmax"),
1196                                           0.0,
1197                                           0.0,
1198                                           "g");
1199   nLostHitsTOBVsPhiH_ = ibook.bookProfile("nLostHitsTOBVsPhi",
1200                                           "Number of Lost Hits in TOB Vs Phi",
1201                                           TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1202                                           TrackPhiHistoPar_.getParameter<double>("Xmin"),
1203                                           TrackPhiHistoPar_.getParameter<double>("Xmax"),
1204                                           0.0,
1205                                           0.0,
1206                                           "g");
1207   nLostHitsTECVsPhiH_ = ibook.bookProfile("nLostHitsTECVsPhi",
1208                                           "Number of Lost Hits in TEC Vs Phi",
1209                                           TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1210                                           TrackPhiHistoPar_.getParameter<double>("Xmin"),
1211                                           TrackPhiHistoPar_.getParameter<double>("Xmax"),
1212                                           0.0,
1213                                           0.0,
1214                                           "g");
1215   nLostHitsTIDVsPhiH_ = ibook.bookProfile("nLostHitsTIDVsPhi",
1216                                           "Number of Lost Hits in TID Vs Phi",
1217                                           TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1218                                           TrackPhiHistoPar_.getParameter<double>("Xmin"),
1219                                           TrackPhiHistoPar_.getParameter<double>("Xmax"),
1220                                           0.0,
1221                                           0.0,
1222                                           "g");
1223 
1224   nLostHitsPixVsIterationH_ = ibook.bookProfile(
1225       "nLostHitsPixVsIteration", "Number of Lost Hits in Pixel Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
1226   nLostHitsPixBVsIterationH_ = ibook.bookProfile(
1227       "nLostHitsPixBVsIteration", "Number of Lost Hits in Pixel Barrel Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
1228   nLostHitsPixEVsIterationH_ = ibook.bookProfile(
1229       "nLostHitsPixEVsIteration", "Number of Lost Hits in Pixel Endcap Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
1230   nLostHitsStripVsIterationH_ = ibook.bookProfile(
1231       "nLostHitsStripVsIteration", "Number of Lost Hits in SiStrip Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
1232   nLostHitsTIBVsIterationH_ = ibook.bookProfile(
1233       "nLostHitsTIBVsIteration", "Number of Lost Hits in TIB Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
1234   nLostHitsTOBVsIterationH_ = ibook.bookProfile(
1235       "nLostHitsTOBVsIteration", "Number of Lost Hits in TOB Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
1236   nLostHitsTECVsIterationH_ = ibook.bookProfile(
1237       "nLostHitsTECVsIteration", "Number of Lost Hits in TEC Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
1238   nLostHitsTIDVsIterationH_ = ibook.bookProfile(
1239       "nLostHitsTIDVsIteration", "Number of Lost Hits in TID Vs Iteration", 47, -0.5, 46.5, 0.0, 0.0, "g");
1240 
1241   trackChi2oNDFVsEtaH_ = ibook.bookProfile("trackChi2oNDFVsEta",
1242                                            "chi2/ndof of Tracks Vs Eta",
1243                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1244                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1245                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1246                                            0.0,
1247                                            10.0,
1248                                            "g");
1249   trackChi2oNDFVsPhiH_ = ibook.bookProfile("trackChi2oNDFVsPhi",
1250                                            "chi2/ndof of Tracks Vs Phi",
1251                                            TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1252                                            TrackPhiHistoPar_.getParameter<double>("Xmin"),
1253                                            TrackPhiHistoPar_.getParameter<double>("Xmax"),
1254                                            0.0,
1255                                            10.0,
1256                                            "g");
1257 
1258   trackChi2probVsEtaH_ = ibook.bookProfile("trackChi2probVsEta",
1259                                            "chi2 probability of Tracks Vs Eta",
1260                                            TrackEtaHistoPar_.getParameter<int32_t>("Xbins"),
1261                                            TrackEtaHistoPar_.getParameter<double>("Xmin"),
1262                                            TrackEtaHistoPar_.getParameter<double>("Xmax"),
1263                                            0.0,
1264                                            1.0,
1265                                            "g");
1266   trackChi2probVsPhiH_ = ibook.bookProfile("trackChi2probVsPhi",
1267                                            "chi2 probability of Tracks Vs Phi",
1268                                            TrackPhiHistoPar_.getParameter<int32_t>("Xbins"),
1269                                            TrackPhiHistoPar_.getParameter<double>("Xmin"),
1270                                            TrackPhiHistoPar_.getParameter<double>("Xmax"),
1271                                            0.0,
1272                                            1.0,
1273                                            "g");
1274 
1275   trackIperr3dVsEtaH_ =
1276       ibook.bookProfile("trackIperr3dVsEta", "ip3d error of Tracks Vs Eta", 80, -4., 4., 0.0, 0.0, "g");
1277 
1278   trackSip2dVsEtaH_ = ibook.bookProfile("trackSip2dVsEta", "sip2d of Tracks Vs Eta", 80, -4., 4., 0.0, 0.0, "g");
1279 
1280   trackIperr3dVsEta2DH_ =
1281       ibook.book2D("trackIperr3dVsEta2D", "ip3d error of Tracks Vs Eta 2d", 80, -4., 4., 100, 0., 5.);
1282   trackIperr3dVsChi2prob2DH_ =
1283       ibook.book2D("trackIperr3dVsChi2prob2D", "ip3d error of Tracks Vs chi2prob 2d", 50, 0., 1., 100, 0., 5.);
1284   trackSip2dVsEta2DH_ = ibook.book2D("trackSip2dVsEta2D", "sip2d of Tracks Vs Eta 2d", 80, -4., 4., 200, -10., 10.);
1285   trackSip2dVsChi2prob2DH_ =
1286       ibook.book2D("trackSip2dVsChi2prob2D", "sip2d of Tracks Vs chi2prob 2d", 50, 0., 1., 200, -10., 10.);
1287 
1288   // On and off-track cluster properties
1289   hOnTrkClusChargeThinH_ = ibook.book1DD("hOnTrkClusChargeThin", "On-track Cluster Charge (Thin Sensor)", 100, 0, 1000);
1290   hOnTrkClusWidthThinH_ = ibook.book1DD("hOnTrkClusWidthThin", "On-track Cluster Width (Thin Sensor)", 20, -0.5, 19.5);
1291   hOnTrkClusChargeThickH_ =
1292       ibook.book1DD("hOnTrkClusChargeThick", "On-track Cluster Charge (Thick Sensor)", 100, 0, 1000);
1293   hOnTrkClusWidthThickH_ =
1294       ibook.book1DD("hOnTrkClusWidthThick", "On-track Cluster Width (Thick Sensor)", 20, -0.5, 19.5);
1295 
1296   hOffTrkClusChargeThinH_ =
1297       ibook.book1DD("hOffTrkClusChargeThin", "Off-track Cluster Charge (Thin Sensor)", 100, 0, 1000);
1298   hOffTrkClusWidthThinH_ =
1299       ibook.book1DD("hOffTrkClusWidthThin", "Off-track Cluster Width (Thin Sensor)", 20, -0.5, 19.5);
1300   hOffTrkClusChargeThickH_ =
1301       ibook.book1DD("hOffTrkClusChargeThick", "Off-track Cluster Charge (Thick Sensor)", 100, 0, 1000);
1302   hOffTrkClusWidthThickH_ =
1303       ibook.book1DD("hOffTrkClusWidthThick", "Off-track Cluster Width (Thick Sensor)", 20, -0.5, 19.5);
1304 }
1305 void StandaloneTrackMonitor::analyze(edm::Event const& iEvent, edm::EventSetup const& iSetup) {
1306   if (verbose_)
1307     edm::LogInfo("StandaloneTrackMonitor") << "Begin StandaloneTrackMonitor" << std::endl;
1308 
1309   nevt++;
1310 
1311   // Get event setup (to get global transformation)
1312   const TrackerGeometry& tkGeom = (*tkGeom_);
1313 
1314   // Primary vertex collection
1315   edm::Handle<reco::VertexCollection> vertexColl;
1316   iEvent.getByToken(vertexToken_, vertexColl);
1317   if (!vertexColl.isValid()) {
1318     edm::LogError("DqmTrackStudy") << "Error! Failed to get reco::Vertex Collection, " << vertexTag_;
1319   }
1320   if (vertexColl->empty()) {
1321     edm::LogError("StandalonaTrackMonitor") << "No good vertex in the event!!" << std::endl;
1322     return;
1323   }
1324   const reco::Vertex& pv = (*vertexColl)[0];
1325 
1326   // Beam spot
1327   edm::Handle<reco::BeamSpot> beamSpot;
1328   iEvent.getByToken(bsToken_, beamSpot);
1329   if (!beamSpot.isValid())
1330     edm::LogError("StandalonaTrackMonitor") << "Beamspot for input tag: " << bsTag_ << " not found!!";
1331 
1332   // Track collection
1333   edm::Handle<reco::TrackCollection> tracks;
1334   iEvent.getByToken(trackToken_, tracks);
1335   if (!tracks.isValid())
1336     edm::LogError("StandalonaTrackMonitor") << "TrackCollection for input tag: " << trackTag_ << " not found!!";
1337 
1338   // Access PU information
1339   double wfac = 1.0;  // for data
1340   if (!iEvent.isRealData()) {
1341     edm::Handle<std::vector<PileupSummaryInfo> > PupInfo;
1342     iEvent.getByToken(puSummaryToken_, PupInfo);
1343 
1344     if (verbose_)
1345       edm::LogInfo("StandaloneTrackMonitor") << "nPUColl = " << PupInfo->size();
1346     if (PupInfo.isValid()) {
1347       for (auto const& v : *PupInfo) {
1348         int bx = v.getBunchCrossing();
1349         if (bunchCrossingH_)
1350           bunchCrossingH_->Fill(bx);
1351         if (bx == 0) {
1352           if (nPUH_)
1353             nPUH_->Fill(v.getPU_NumInteractions());
1354           int ntrueInt = v.getTrueNumInteractions();
1355           int nVertex = (vertexColl.isValid() ? vertexColl->size() : 0);
1356           if (trueNIntH_)
1357             trueNIntH_->Fill(ntrueInt);
1358           if (doPUCorrection_) {
1359             if (nVertex > -1 && nVertex < int(vpu_.size()))
1360               wfac = vpu_.at(nVertex);
1361             else
1362               wfac = 0.0;
1363           }
1364         }
1365       }
1366     } else
1367       edm::LogError("StandalonaTrackMonitor") << "PUSummary for input tag: " << puSummaryTag_ << " not found!!";
1368     if (doTrackCorrection_) {
1369       int ntrack = 0;
1370       for (auto const& track : *tracks) {
1371         if (!track.quality(reco::Track::qualityByName(trackQuality_)))
1372           continue;
1373         ++ntrack;
1374       }
1375       if (ntrack > -1 && ntrack < int(vtrack_.size()))
1376         wfac = vtrack_.at(ntrack);
1377       else
1378         wfac = 0.0;
1379     }
1380   }
1381   if (verbose_)
1382     edm::LogInfo("StandaloneTrackMonitor") << "PU reweight factor = " << wfac;
1383 
1384   if (haveAllHistograms_) {
1385     int nvtx = (vertexColl.isValid() ? vertexColl->size() : 0);
1386     nVertexH_->Fill(nvtx, wfac);
1387     nVtxH_->Fill(nvtx);
1388   }
1389 
1390   // Get MVA and quality mask collections
1391   int ntracks = 0;
1392   std::vector<TLorentzVector> list;
1393 
1394   const TransientTrackBuilder* theB = &iSetup.getData(transTrackToken_);
1395   TransientVertex mumuTransientVtx;
1396   std::vector<reco::TransientTrack> tks;
1397 
1398   if (tracks.isValid()) {
1399     if (verbose_)
1400       edm::LogInfo("StandaloneTrackMonitor") << "Total # of Tracks: " << tracks->size();
1401     reco::Track::TrackQuality quality = reco::Track::qualityByName(trackQuality_);
1402 
1403     for (auto const& track : *tracks) {
1404       if (!track.quality(quality))
1405         continue;
1406       ++ntracks;
1407 
1408       reco::TransientTrack trajectory = theB->build(track);
1409       tks.push_back(trajectory);
1410 
1411       double eta = track.eta();
1412       double theta = track.theta();
1413       double phi = track.phi();
1414       double pt = track.pt();
1415 
1416       const reco::HitPattern& hitp = track.hitPattern();
1417       double nAllHits = hitp.numberOfAllHits(reco::HitPattern::TRACK_HITS);
1418       double nAllTrackerHits = hitp.numberOfAllTrackerHits(reco::HitPattern::TRACK_HITS);
1419 
1420       double trackdeltaR = std::numeric_limits<double>::max();
1421       double muonmass = 0.1056583745;
1422       TLorentzVector trackpmu;
1423       trackpmu.SetPtEtaPhiM(track.pt(), track.eta(), track.phi(), muonmass);
1424       for (auto const& TRACK : *tracks) {
1425         if (&track == &TRACK)
1426           continue;
1427         double tmpdr = reco::deltaR(track.eta(), track.phi(), TRACK.eta(), TRACK.phi());
1428         if (tmpdr < trackdeltaR)
1429           trackdeltaR = tmpdr;
1430       }
1431 
1432       list.push_back(trackpmu);
1433 
1434       double nValidTrackerHits = hitp.numberOfValidTrackerHits();
1435       double nValidPixelHits = hitp.numberOfValidPixelHits();
1436       double nValidPixelBHits = hitp.numberOfValidPixelBarrelHits();
1437       double nValidPixelEHits = hitp.numberOfValidPixelEndcapHits();
1438       double nValidStripHits = hitp.numberOfValidStripHits();
1439       double nValidTIBHits = hitp.numberOfValidStripTIBHits();
1440       double nValidTOBHits = hitp.numberOfValidStripTOBHits();
1441       double nValidTIDHits = hitp.numberOfValidStripTIDHits();
1442       double nValidTECHits = hitp.numberOfValidStripTECHits();
1443 
1444       int missingInnerHit = hitp.numberOfAllHits(reco::HitPattern::MISSING_INNER_HITS);
1445       int missingOuterHit = hitp.numberOfAllHits(reco::HitPattern::MISSING_OUTER_HITS);
1446 
1447       nValidHitsVspTH_->Fill(pt, nValidTrackerHits);
1448       nValidHitsVsEtaH_->Fill(eta, nValidTrackerHits);
1449       nValidHitsVsCosThetaH_->Fill(std::cos(theta), nValidTrackerHits);
1450       nValidHitsVsPhiH_->Fill(phi, nValidTrackerHits);
1451       nValidHitsVsnVtxH_->Fill(vertexColl->size(), nValidTrackerHits);
1452 
1453       nValidHitsPixVsEtaH_->Fill(eta, nValidPixelHits);
1454       nValidHitsPixBVsEtaH_->Fill(eta, nValidPixelBHits);
1455       nValidHitsPixEVsEtaH_->Fill(eta, nValidPixelEHits);
1456       nValidHitsStripVsEtaH_->Fill(eta, nValidStripHits);
1457       nValidHitsTIBVsEtaH_->Fill(eta, nValidTIBHits);
1458       nValidHitsTOBVsEtaH_->Fill(eta, nValidTOBHits);
1459       nValidHitsTECVsEtaH_->Fill(eta, nValidTECHits);
1460       nValidHitsTIDVsEtaH_->Fill(eta, nValidTIDHits);
1461 
1462       nValidHitsPixVsPhiH_->Fill(phi, nValidPixelHits);
1463       nValidHitsPixBVsPhiH_->Fill(phi, nValidPixelBHits);
1464       nValidHitsPixEVsPhiH_->Fill(phi, nValidPixelEHits);
1465       nValidHitsStripVsPhiH_->Fill(phi, nValidStripHits);
1466       nValidHitsTIBVsPhiH_->Fill(phi, nValidTIBHits);
1467       nValidHitsTOBVsPhiH_->Fill(phi, nValidTOBHits);
1468       nValidHitsTECVsPhiH_->Fill(phi, nValidTECHits);
1469       nValidHitsTIDVsPhiH_->Fill(phi, nValidTIDHits);
1470 
1471       int nLostHits = track.numberOfLostHits();
1472       int nMissingExpectedInnerHits = hitp.numberOfLostHits(reco::HitPattern::MISSING_INNER_HITS);
1473       int nMissingExpectedOuterHits = hitp.numberOfLostHits(reco::HitPattern::MISSING_OUTER_HITS);
1474       int nLostTrackerHits = hitp.numberOfLostTrackerHits(reco::HitPattern::TRACK_HITS);
1475       int nLostPixHits = hitp.numberOfLostPixelHits(reco::HitPattern::TRACK_HITS);
1476       int nLostPixBHits = hitp.numberOfLostPixelBarrelHits(reco::HitPattern::TRACK_HITS);
1477       int nLostPixEHits = hitp.numberOfLostPixelEndcapHits(reco::HitPattern::TRACK_HITS);
1478       int nLostStripHits = hitp.numberOfLostStripHits(reco::HitPattern::TRACK_HITS);
1479       int nLostStripTIBHits = hitp.numberOfLostStripTIBHits(reco::HitPattern::TRACK_HITS);
1480       int nLostStripTIDHits = hitp.numberOfLostStripTIDHits(reco::HitPattern::TRACK_HITS);
1481       int nLostStripTOBHits = hitp.numberOfLostStripTOBHits(reco::HitPattern::TRACK_HITS);
1482       int nLostStripTECHits = hitp.numberOfLostStripTECHits(reco::HitPattern::TRACK_HITS);
1483       int nIteration = track.originalAlgo();
1484 
1485       nLostHitsVspTH_->Fill(pt, nLostTrackerHits);
1486       nLostHitsVsEtaH_->Fill(eta, nLostTrackerHits);
1487       nLostHitsVsCosThetaH_->Fill(std::cos(theta), nLostTrackerHits);
1488       nLostHitsVsPhiH_->Fill(phi, nLostTrackerHits);
1489       nLostHitsVsIterationH_->Fill(nIteration, nLostTrackerHits);
1490 
1491       nLostHitsPixVsEtaH_->Fill(eta, nLostPixHits);
1492       nLostHitsPixBVsEtaH_->Fill(eta, nLostPixBHits);
1493       nLostHitsPixEVsEtaH_->Fill(eta, nLostPixEHits);
1494       nLostHitsStripVsEtaH_->Fill(eta, nLostStripHits);
1495       nLostHitsTIBVsEtaH_->Fill(eta, nLostStripTIBHits);
1496       nLostHitsTOBVsEtaH_->Fill(eta, nLostStripTOBHits);
1497       nLostHitsTECVsEtaH_->Fill(eta, nLostStripTECHits);
1498       nLostHitsTIDVsEtaH_->Fill(eta, nLostStripTIDHits);
1499 
1500       nLostHitsPixVsPhiH_->Fill(phi, nLostPixHits);
1501       nLostHitsPixBVsPhiH_->Fill(phi, nLostPixBHits);
1502       nLostHitsPixEVsPhiH_->Fill(phi, nLostPixEHits);
1503       nLostHitsStripVsPhiH_->Fill(phi, nLostStripHits);
1504       nLostHitsTIBVsPhiH_->Fill(phi, nLostStripTIBHits);
1505       nLostHitsTOBVsPhiH_->Fill(phi, nLostStripTOBHits);
1506       nLostHitsTECVsPhiH_->Fill(phi, nLostStripTECHits);
1507       nLostHitsTIDVsPhiH_->Fill(phi, nLostStripTIDHits);
1508 
1509       nLostHitsPixVsIterationH_->Fill(nIteration, nLostPixHits);
1510       nLostHitsPixBVsIterationH_->Fill(nIteration, nLostPixBHits);
1511       nLostHitsPixEVsIterationH_->Fill(nIteration, nLostPixEHits);
1512       nLostHitsStripVsIterationH_->Fill(nIteration, nLostStripHits);
1513       nLostHitsTIBVsIterationH_->Fill(nIteration, nLostStripTIBHits);
1514       nLostHitsTOBVsIterationH_->Fill(nIteration, nLostStripTOBHits);
1515       nLostHitsTECVsIterationH_->Fill(nIteration, nLostStripTECHits);
1516       nLostHitsTIDVsIterationH_->Fill(nIteration, nLostStripTIDHits);
1517 
1518       if (abs(eta) <= 1.4) {
1519         nMissingInnerHitBH_->Fill(missingInnerHit, wfac);
1520         nMissingOuterHitBH_->Fill(missingOuterHit, wfac);
1521       } else {
1522         nMissingInnerHitEH_->Fill(missingInnerHit, wfac);
1523         nMissingOuterHitEH_->Fill(missingOuterHit, wfac);
1524       }
1525 
1526       for (int i = 0; i < hitp.numberOfAllHits(reco::HitPattern::TRACK_HITS); i++) {
1527         uint32_t hit = hitp.getHitPattern(reco::HitPattern::TRACK_HITS, i);
1528         if (hitp.missingHitFilter(hit)) {
1529           double losthitBylayer = -1.0;
1530           double losthitBylayerPix = -1.0;
1531           double losthitBylayerStrip = -1.0;
1532           int layer = hitp.getLayer(hit);
1533           if (hitp.pixelBarrelHitFilter(hit)) {
1534             losthitBylayer = layer;
1535             losthitBylayerPix = layer;
1536           } else if (hitp.pixelEndcapHitFilter(hit)) {
1537             losthitBylayer = layer + 4;
1538             losthitBylayerPix = layer + 4;
1539           } else if (hitp.stripTIBHitFilter(hit)) {
1540             losthitBylayer = layer + 7;
1541             losthitBylayerStrip = layer;
1542           } else if (hitp.stripTIDHitFilter(hit)) {
1543             losthitBylayer = layer + 11;
1544             losthitBylayerStrip = layer + 4;
1545           } else if (hitp.stripTOBHitFilter(hit)) {
1546             losthitBylayer = layer + 14;
1547             losthitBylayerStrip = layer + 7;
1548           } else if (hitp.stripTECHitFilter(hit)) {
1549             losthitBylayer = layer + 20;
1550             losthitBylayerStrip = layer + 13;
1551           }
1552           if (losthitBylayer > -1)
1553             nLostHitByLayerH_->Fill(losthitBylayer, wfac);
1554           if (losthitBylayerPix > -1)
1555             nLostHitByLayerPixH_->Fill(losthitBylayerPix, wfac);
1556           if (losthitBylayerStrip > -1)
1557             nLostHitByLayerStripH_->Fill(losthitBylayerStrip, wfac);
1558         }
1559       }
1560 
1561       if (haveAllHistograms_) {
1562         double etaError = track.etaError();
1563         double thetaError = track.thetaError();
1564         double phiError = track.phiError();
1565         double p = track.p();
1566         double ptError = track.ptError();
1567         double qoverp = track.qoverp();
1568         double qoverpError = track.qoverpError();
1569         double charge = track.charge();
1570 
1571         double dxy = track.dxy(beamSpot->position());
1572         double dxyError = track.dxyError();
1573         double dz = track.dz(beamSpot->position());
1574         double dzError = track.dzError();
1575 
1576         double trkd0 = track.d0();
1577         double chi2 = track.chi2();
1578         double ndof = track.ndof();
1579         double chi2prob = TMath::Prob(track.chi2(), (int)track.ndof());
1580         double chi2oNDF = track.normalizedChi2();
1581         double vx = track.vx();
1582         double vy = track.vy();
1583         double vz = track.vz();
1584         //// algorithm
1585         unsigned int track_algo = track.algo();
1586         unsigned int track_origalgo = track.originalAlgo();
1587         // stopping source
1588         int ssmax = trackStoppingSourceH_->getNbinsX();
1589         double stop = track.stopReason() > ssmax ? double(ssmax - 1) : static_cast<double>(track.stopReason());
1590         double distanceOfClosestApproachToPV = track.dxy(pv.position());
1591         double xPointOfClosestApproachwrtPV = track.vx() - pv.position().x();
1592         double yPointOfClosestApproachwrtPV = track.vy() - pv.position().y();
1593         double positionZ0 = track.dz(pv.position());
1594 
1595         reco::TransientTrack transTrack = iSetup.get<TransientTrackRecord>().get(transTrackToken_).build(track);
1596 
1597         double ip3dToPV = 0, iperr3dToPV = 0, sip3dToPV = 0, sip2dToPV = 0;
1598         GlobalVector dir(track.px(), track.py(), track.pz());
1599         std::pair<bool, Measurement1D> ip3d = IPTools::signedImpactParameter3D(transTrack, dir, pv);
1600         std::pair<bool, Measurement1D> ip2d = IPTools::signedTransverseImpactParameter(transTrack, dir, pv);
1601         if (ip3d.first) {
1602           sip3dToPV = ip3d.second.value() / ip3d.second.error();
1603           ip3dToPV = ip3d.second.value();
1604           iperr3dToPV = ip3d.second.error();
1605         }
1606 
1607         double ip3dToBS = 0, iperr3dToBS = 0, sip3dToBS = 0, sip2dToBS = 0;
1608         reco::Vertex beamspotvertex((*beamSpot).position(), (*beamSpot).covariance3D());
1609         std::pair<bool, Measurement1D> ip3dbs = IPTools::signedImpactParameter3D(transTrack, dir, beamspotvertex);
1610         std::pair<bool, Measurement1D> ip2dbs =
1611             IPTools::signedTransverseImpactParameter(transTrack, dir, beamspotvertex);
1612         if (ip3dbs.first) {
1613           sip3dToBS = ip3dbs.second.value() / ip3dbs.second.error();
1614           ip3dToBS = ip3dbs.second.value();
1615           iperr3dToBS = ip3dbs.second.error();
1616         }
1617 
1618         double ip2dToPV = 0, iperr2dToPV = 0;
1619         if (ip2d.first) {
1620           ip2dToPV = ip2d.second.value();
1621           iperr2dToPV = ip2d.second.error();
1622           sip2dToPV = ip2d.second.value() / ip2d.second.error();
1623         }
1624 
1625         double ip2dToBS = 0, iperr2dToBS = 0;
1626         if (ip2dbs.first) {
1627           ip2dToBS = ip2dbs.second.value();
1628           iperr2dToBS = ip2dbs.second.error();
1629           sip2dToBS = ip2d.second.value() / ip2d.second.error();
1630         }
1631 
1632         if (ip2d.first)
1633           sip2dToPV = ip2d.second.value() / ip2d.second.error();
1634         double sipDxyToPV = track.dxy(pv.position()) / track.dxyError();
1635         double sipDzToPV = track.dz(pv.position()) / track.dzError();
1636 
1637         // Fill the histograms
1638         trackDeltaRwrtClosestTrack_->Fill(trackdeltaR, wfac);
1639         trackEtaH_->Fill(eta, wfac);
1640         trackEtaerrH_->Fill(etaError, wfac);
1641         trackCosThetaH_->Fill(std::cos(theta), wfac);
1642         trackThetaerrH_->Fill(thetaError, wfac);
1643         trackPhiH_->Fill(phi, wfac);
1644         trackPhierrH_->Fill(phiError, wfac);
1645         trackPH_->Fill(p, wfac);
1646         trackPtH_->Fill(pt, wfac);
1647         trackPt_ZoomH_->Fill(pt, wfac);
1648         trackPterrH_->Fill(ptError, wfac);
1649         trackqOverpH_->Fill(qoverp, wfac);
1650         trackqOverperrH_->Fill(qoverpError, wfac);
1651         trackChargeH_->Fill(charge, wfac);
1652         trackChi2H_->Fill(chi2, wfac);
1653         trackChi2ProbH_->Fill(chi2prob, wfac);
1654         trackChi2oNDFH_->Fill(chi2oNDF, wfac);
1655         trackd0H_->Fill(trkd0, wfac);
1656         tracknDOFH_->Fill(ndof, wfac);
1657         trackChi2bynDOFH_->Fill(chi2 / ndof, wfac);
1658         trackalgoH_->Fill(track_algo, wfac);
1659         trackorigalgoH_->Fill(track_origalgo, wfac);
1660         trackStoppingSourceH_->Fill(stop, wfac);
1661         trackChi2oNDFVsEtaH_->Fill(eta, chi2oNDF);
1662         trackChi2oNDFVsPhiH_->Fill(phi, chi2oNDF);
1663         trackChi2probVsEtaH_->Fill(eta, chi2prob);
1664         trackChi2probVsPhiH_->Fill(phi, chi2prob);
1665 
1666         nlostHitsH_->Fill(nLostHits, wfac);
1667         nMissingExpectedInnerHitsH_->Fill(nMissingExpectedInnerHits, wfac);
1668         nMissingExpectedOuterHitsH_->Fill(nMissingExpectedOuterHits, wfac);
1669         nlostTrackerHitsH_->Fill(nLostTrackerHits, wfac);
1670 
1671         beamSpotXYposH_->Fill(dxy, wfac);
1672         beamSpotXYposerrH_->Fill(dxyError, wfac);
1673         beamSpotZposH_->Fill(dz, wfac);
1674         beamSpotZposerrH_->Fill(dzError, wfac);
1675 
1676         vertexXposH_->Fill(vx, wfac);
1677         vertexYposH_->Fill(vy, wfac);
1678         vertexZposH_->Fill(vz, wfac);
1679 
1680         int nPixBarrel = 0, nPixEndcap = 0, nStripTIB = 0, nStripTOB = 0, nStripTEC = 0, nStripTID = 0;
1681         if (isRECO_) {
1682           for (auto it = track.recHitsBegin(); it != track.recHitsEnd(); ++it) {
1683             const TrackingRecHit& hit = (**it);
1684             if (hit.isValid()) {
1685               if (hit.geographicalId().det() == DetId::Tracker) {
1686                 int subdetId = hit.geographicalId().subdetId();
1687                 if (subdetId == PixelSubdetector::PixelBarrel)
1688                   ++nPixBarrel;
1689                 else if (subdetId == PixelSubdetector::PixelEndcap)
1690                   ++nPixEndcap;
1691                 else if (subdetId == StripSubdetector::TIB)
1692                   ++nStripTIB;
1693                 else if (subdetId == StripSubdetector::TOB)
1694                   ++nStripTOB;
1695                 else if (subdetId == StripSubdetector::TEC)
1696                   ++nStripTEC;
1697                 else if (subdetId == StripSubdetector::TID)
1698                   ++nStripTID;
1699 
1700                 // Find on-track clusters
1701                 processHit(hit, iSetup, tkGeom, wfac);
1702               }
1703             }
1704           }
1705         }
1706         if (verbose_)
1707           edm::LogInfo("StandaloneTrackMonitor")
1708               << " >>> HITs: nPixBarrel: " << nPixBarrel << " nPixEndcap: " << nPixEndcap << " nStripTIB: " << nStripTIB
1709               << " nStripTOB: " << nStripTOB << " nStripTEC: " << nStripTEC << " nStripTID: " << nStripTID;
1710         if (haveAllHistograms_) {
1711           nPixBarrelH_->Fill(nPixBarrel, wfac);
1712           nPixEndcapH_->Fill(nPixEndcap, wfac);
1713           nStripTIBH_->Fill(nStripTIB, wfac);
1714           nStripTOBH_->Fill(nStripTOB, wfac);
1715           nStripTECH_->Fill(nStripTEC, wfac);
1716           nStripTIDH_->Fill(nStripTID, wfac);
1717         }
1718 
1719         DistanceOfClosestApproachToPVH_->Fill(distanceOfClosestApproachToPV, wfac);
1720         DistanceOfClosestApproachToPVZoomedH_->Fill(distanceOfClosestApproachToPV, wfac);
1721         DistanceOfClosestApproachToPVVsPhiH_->Fill(phi, distanceOfClosestApproachToPV);
1722         xPointOfClosestApproachVsZ0wrtPVH_->Fill(positionZ0, xPointOfClosestApproachwrtPV);
1723         yPointOfClosestApproachVsZ0wrtPVH_->Fill(positionZ0, yPointOfClosestApproachwrtPV);
1724 
1725         ip3dToPVH_->Fill(ip3dToPV, wfac);
1726         iperr3dToPVH_->Fill(iperr3dToPV, wfac);
1727         ip3dToBSH_->Fill(ip3dToBS, wfac);
1728         iperr3dToBSH_->Fill(iperr3dToBS, wfac);
1729         ip2dToPVH_->Fill(ip2dToPV, wfac);
1730         iperr2dToPVH_->Fill(iperr2dToPV, wfac);
1731         iperr2dToPVWtH_->Fill(iperr2dToPV, wfac);
1732         ip2dToBSH_->Fill(ip2dToBS, wfac);
1733         iperr2dToBSH_->Fill(iperr2dToBS, wfac);
1734 
1735         iperr3dToPVWtH_->Fill(iperr3dToPV, wfac);
1736         sip3dToPVH_->Fill(sip3dToPV, wfac);
1737         sip2dToPVH_->Fill(sip2dToPV, wfac);
1738         sip3dToBSH_->Fill(sip3dToBS, wfac);
1739         sip2dToBSH_->Fill(sip2dToBS, wfac);
1740         sip2dToPVWtH_->Fill(sip2dToPV, wfac);
1741         sipDxyToPVH_->Fill(sipDxyToPV, wfac);
1742         sipDzToPVH_->Fill(sipDzToPV, wfac);
1743 
1744         if (nValidPixelHits >= 2) {
1745           ip3dToPV2validpixelhitsH_->Fill(ip3dToPV, wfac);
1746           iperr3dToPV2validpixelhitsH_->Fill(iperr3dToPV, wfac);
1747           ip3dToBS2validpixelhitsH_->Fill(ip3dToBS, wfac);
1748           iperr3dToBS2validpixelhitsH_->Fill(iperr3dToBS, wfac);
1749           sip3dToPV2validpixelhitsH_->Fill(sip3dToPV, wfac);
1750           sip3dToBS2validpixelhitsH_->Fill(sip3dToBS, wfac);
1751         }
1752 
1753         trackIperr3dVsEta2DH_->Fill(eta, iperr3dToPV, wfac);
1754         trackSip2dVsEta2DH_->Fill(eta, sip2dToPV, wfac);
1755         trackIperr3dVsChi2prob2DH_->Fill(chi2prob, iperr3dToPV, wfac);
1756         trackSip2dVsChi2prob2DH_->Fill(chi2prob, sip2dToPV, wfac);
1757 
1758         trackIperr3dVsEtaH_->Fill(eta, iperr3dToPV);
1759         trackSip2dVsEtaH_->Fill(eta, sip2dToPV);
1760 
1761         double trackerLayersWithMeasurement = hitp.trackerLayersWithMeasurement();
1762         double pixelLayersWithMeasurement = hitp.pixelLayersWithMeasurement();
1763         double pixelBLayersWithMeasurement = hitp.pixelBarrelLayersWithMeasurement();
1764         double pixelELayersWithMeasurement = hitp.pixelEndcapLayersWithMeasurement();
1765         double stripLayersWithMeasurement = hitp.stripLayersWithMeasurement();
1766         double stripTIBLayersWithMeasurement = hitp.stripTIBLayersWithMeasurement();
1767         double stripTOBLayersWithMeasurement = hitp.stripTOBLayersWithMeasurement();
1768         double stripTIDLayersWithMeasurement = hitp.stripTIDLayersWithMeasurement();
1769         double stripTECLayersWithMeasurement = hitp.stripTECLayersWithMeasurement();
1770 
1771         trkLayerwithMeasurementH_->Fill(trackerLayersWithMeasurement, wfac);
1772         pixelLayerwithMeasurementH_->Fill(pixelLayersWithMeasurement, wfac);
1773         pixelBLayerwithMeasurementH_->Fill(pixelBLayersWithMeasurement, wfac);
1774         pixelELayerwithMeasurementH_->Fill(pixelELayersWithMeasurement, wfac);
1775         stripLayerwithMeasurementH_->Fill(stripLayersWithMeasurement, wfac);
1776         stripTIBLayerwithMeasurementH_->Fill(stripTIBLayersWithMeasurement, wfac);
1777         stripTOBLayerwithMeasurementH_->Fill(stripTOBLayersWithMeasurement, wfac);
1778         stripTIDLayerwithMeasurementH_->Fill(stripTIDLayersWithMeasurement, wfac);
1779         stripTECLayerwithMeasurementH_->Fill(stripTECLayersWithMeasurement, wfac);
1780 
1781         nallHitsH_->Fill(nAllHits, wfac);
1782         ntrackerHitsH_->Fill(nAllTrackerHits, wfac);
1783         nvalidTrackerHitsH_->Fill(nValidTrackerHits, wfac);
1784         nvalidPixelHitsH_->Fill(nValidPixelHits, wfac);
1785         nvalidPixelBHitsH_->Fill(nValidPixelBHits, wfac);
1786         nvalidPixelEHitsH_->Fill(nValidPixelEHits, wfac);
1787         nvalidStripHitsH_->Fill(nValidStripHits, wfac);
1788         nvalidTIBHitsH_->Fill(nValidTIBHits, wfac);
1789         nvalidTOBHitsH_->Fill(nValidTOBHits, wfac);
1790         nvalidTIDHitsH_->Fill(nValidTIDHits, wfac);
1791         nvalidTECHitsH_->Fill(nValidTECHits, wfac);
1792 
1793         nlostTrackerHitsH_->Fill(nLostTrackerHits, wfac);
1794         nlostPixelHitsH_->Fill(nLostPixHits, wfac);
1795         nlostPixelBHitsH_->Fill(nLostPixBHits, wfac);
1796         nlostPixelEHitsH_->Fill(nLostPixEHits, wfac);
1797         nlostStripHitsH_->Fill(nLostStripHits, wfac);
1798         nlostTIBHitsH_->Fill(nLostStripTIBHits, wfac);
1799         nlostTOBHitsH_->Fill(nLostStripTOBHits, wfac);
1800         nlostTIDHitsH_->Fill(nLostStripTIDHits, wfac);
1801         nlostTECHitsH_->Fill(nLostStripTECHits, wfac);
1802       }
1803     }
1804   } else {
1805     edm::LogError("StandaloneTrackMonitor") << "Error! Failed to get reco::Track collection, " << trackTag_;
1806   }
1807 
1808   if (haveAllHistograms_) {
1809     nTracksH_->Fill(ntracks, wfac);
1810     edm::Handle<std::vector<reco::PFJet> > jetsColl;
1811     iEvent.getByToken(jetsToken_, jetsColl);
1812     nJet_->Fill(jetsColl->size());
1813 
1814     for (auto const& jet : *jetsColl) {
1815       Jet_pt_->Fill(jet.pt(), wfac);
1816       Jet_eta_->Fill(jet.eta(), wfac);
1817       Jet_energy_->Fill(jet.energy(), wfac);
1818       Jet_chargedMultiplicity_->Fill(jet.chargedHadronMultiplicity(), wfac);
1819     }
1820     if (list.size() >= 2) {
1821       TLorentzVector Zp = list[0] + list[1];
1822       Zpt_->Fill(Zp.Pt(), wfac);
1823       ZInvMass_->Fill(Zp.Mag(), wfac);
1824       KalmanVertexFitter kalman(true);
1825       mumuTransientVtx = kalman.vertex(tks);
1826       if (mumuTransientVtx.isValid()) {
1827         const reco::Vertex* theClosestVertex;
1828         edm::Handle<reco::VertexCollection> vertexHandle = iEvent.getHandle(vertexToken_);
1829         if (vertexHandle.isValid()) {
1830           const reco::VertexCollection* vertices = vertexHandle.product();
1831           theClosestVertex = this->findClosestVertex(mumuTransientVtx, vertices);
1832           reco::Vertex theMainVtx;
1833           if (theClosestVertex == nullptr) {
1834             theMainVtx = vertexHandle.product()->front();
1835           } else {
1836             theMainVtx = *theClosestVertex;
1837           }
1838           if (theMainVtx.isValid()) {
1839             const auto& mainVertexPos = theMainVtx.position();
1840             const auto& myVertexPos = mumuTransientVtx.position();
1841 
1842             const double deltaX = myVertexPos.x() - mainVertexPos.x();
1843             const double deltaY = myVertexPos.y() - mainVertexPos.y();
1844             const double deltaZ = myVertexPos.z() - mainVertexPos.z();
1845 
1846             const double deltaNorm = sqrt(deltaX * deltaX + deltaY * deltaY + deltaZ * deltaZ);
1847             const double zpNorm = sqrt(Zp.Px() * Zp.Px() + Zp.Py() * Zp.Py() + Zp.Pz() * Zp.Pz());
1848 
1849             const double cosPhi3D = (Zp.Px() * deltaX + Zp.Py() * deltaY + Zp.Pz() * deltaZ) / (zpNorm * deltaNorm);
1850             cosPhi3DdileptonH_->Fill(cosPhi3D, wfac);
1851           }
1852         }
1853       }
1854     }
1855   }
1856 
1857   // off track cluster properties (only on RECO data-tier)
1858   if (isRECO_) {
1859     processClusters(iEvent, iSetup, tkGeom, wfac);
1860   }
1861 
1862   if (verbose_)
1863     edm::LogInfo("StandaloneTrackMonitor") << "Ends StandaloneTrackMonitor successfully";
1864 }
1865 void StandaloneTrackMonitor::processClusters(edm::Event const& iEvent,
1866                                              edm::EventSetup const& iSetup,
1867                                              const TrackerGeometry& tkGeom,
1868                                              double wfac) {
1869   // SiStripClusters
1870   edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusterHandle;
1871   iEvent.getByToken(clusterToken_, clusterHandle);
1872 
1873   if (clusterHandle.isValid()) {
1874     // Loop on Dets
1875     for (edmNew::DetSetVector<SiStripCluster>::const_iterator dsvit = clusterHandle->begin();
1876          dsvit != clusterHandle->end();
1877          ++dsvit) {
1878       uint32_t detId = dsvit->id();
1879       std::map<uint32_t, std::set<const SiStripCluster*> >::iterator jt = clusterMap_.find(detId);
1880       bool detid_found = (jt != clusterMap_.end()) ? true : false;
1881 
1882       // Loop on Clusters
1883       for (edmNew::DetSet<SiStripCluster>::const_iterator clusit = dsvit->begin(); clusit != dsvit->end(); ++clusit) {
1884         if (detid_found) {
1885           std::set<const SiStripCluster*>& s = jt->second;
1886           if (s.find(&*clusit) != s.end())
1887             continue;
1888         }
1889 
1890         siStripClusterInfo_.setCluster(*clusit, detId);
1891         float charge = siStripClusterInfo_.charge();
1892         float width = siStripClusterInfo_.width();
1893 
1894         const GeomDetUnit* detUnit = tkGeom.idToDetUnit(detId);
1895         float thickness = detUnit->surface().bounds().thickness();  // unit cm
1896         if (thickness > 0.035) {
1897           hOffTrkClusChargeThickH_->Fill(charge, wfac);
1898           hOffTrkClusWidthThickH_->Fill(width, wfac);
1899         } else {
1900           hOffTrkClusChargeThinH_->Fill(charge, wfac);
1901           hOffTrkClusWidthThinH_->Fill(width, wfac);
1902         }
1903       }
1904     }
1905   } else {
1906     edm::LogError("StandaloneTrackMonitor") << "ClusterCollection " << clusterTag_ << " not valid!!" << std::endl;
1907   }
1908 }
1909 void StandaloneTrackMonitor::processHit(const TrackingRecHit& recHit,
1910                                         edm::EventSetup const& iSetup,
1911                                         const TrackerGeometry& tkGeom,
1912                                         double wfac) {
1913   uint32_t detid = recHit.geographicalId();
1914   const GeomDetUnit* detUnit = tkGeom.idToDetUnit(detid);
1915   float thickness = detUnit->surface().bounds().thickness();  // unit cm
1916 
1917   auto const& thit = static_cast<BaseTrackerRecHit const&>(recHit);
1918   if (!thit.isValid())
1919     return;
1920 
1921   auto const& clus = thit.firstClusterRef();
1922   if (!clus.isValid())
1923     return;
1924   if (!clus.isStrip())
1925     return;
1926 
1927   if (thit.isMatched()) {
1928     const SiStripMatchedRecHit2D& matchedHit = dynamic_cast<const SiStripMatchedRecHit2D&>(recHit);
1929 
1930     auto& clusterM = matchedHit.monoCluster();
1931     siStripClusterInfo_.setCluster(clusterM, detid);
1932 
1933     if (thickness > 0.035) {
1934       hOnTrkClusChargeThickH_->Fill(siStripClusterInfo_.charge(), wfac);
1935       hOnTrkClusWidthThickH_->Fill(siStripClusterInfo_.width(), wfac);
1936     } else {
1937       hOnTrkClusChargeThinH_->Fill(siStripClusterInfo_.charge(), wfac);
1938       hOnTrkClusWidthThinH_->Fill(siStripClusterInfo_.width(), wfac);
1939     }
1940     addClusterToMap(detid, &clusterM);
1941 
1942     auto& clusterS = matchedHit.stereoCluster();
1943 
1944     siStripClusterInfo_.setCluster(clusterS, detid);
1945     if (thickness > 0.035) {
1946       hOnTrkClusChargeThickH_->Fill(siStripClusterInfo_.charge(), wfac);
1947       hOnTrkClusWidthThickH_->Fill(siStripClusterInfo_.width(), wfac);
1948     } else {
1949       hOnTrkClusChargeThinH_->Fill(siStripClusterInfo_.charge(), wfac);
1950       hOnTrkClusWidthThinH_->Fill(siStripClusterInfo_.width(), wfac);
1951     }
1952     addClusterToMap(detid, &clusterS);
1953   } else {
1954     auto& cluster = clus.stripCluster();
1955     siStripClusterInfo_.setCluster(cluster, detid);
1956     if (thickness > 0.035) {
1957       hOnTrkClusChargeThickH_->Fill(siStripClusterInfo_.charge(), wfac);
1958       hOnTrkClusWidthThickH_->Fill(siStripClusterInfo_.width(), wfac);
1959     } else {
1960       hOnTrkClusChargeThinH_->Fill(siStripClusterInfo_.charge(), wfac);
1961       hOnTrkClusWidthThinH_->Fill(siStripClusterInfo_.width(), wfac);
1962     }
1963 
1964     addClusterToMap(detid, &cluster);
1965   }
1966 }
1967 void StandaloneTrackMonitor::addClusterToMap(uint32_t detid, const SiStripCluster* cluster) {
1968   std::map<uint32_t, std::set<const SiStripCluster*> >::iterator it = clusterMap_.find(detid);
1969   if (it == clusterMap_.end()) {
1970     std::set<const SiStripCluster*> s;
1971     s.insert(cluster);
1972     clusterMap_.insert(std::pair<uint32_t, std::set<const SiStripCluster*> >(detid, s));
1973   } else {
1974     std::set<const SiStripCluster*>& s = it->second;
1975     s.insert(cluster);
1976   }
1977 }
1978 
1979 const reco::Vertex* StandaloneTrackMonitor::findClosestVertex(
1980     const TransientVertex aTransVtx,
1981     const reco::VertexCollection* vertices) const {  //same as in /DQMOffline/Alignment/src/DiMuonVertexSelector.cc
1982   reco::Vertex* defaultVtx = nullptr;
1983 
1984   if (!aTransVtx.isValid())
1985     return defaultVtx;
1986 
1987   // find the closest vertex to the secondary vertex in 3D
1988   VertexDistance3D vertTool3D;
1989   float minD = 9999.;
1990   int closestVtxIndex = 0;
1991   int counter = 0;
1992   for (const auto& vtx : *vertices) {
1993     double dist3D = vertTool3D.distance(aTransVtx, vtx).value();
1994     if (dist3D < minD) {
1995       minD = dist3D;
1996       closestVtxIndex = counter;
1997     }
1998     counter++;
1999   }
2000 
2001   if ((*vertices).at(closestVtxIndex).isValid()) {
2002     return &(vertices->at(closestVtxIndex));
2003   } else {
2004     return defaultVtx;
2005   }
2006 }
2007 
2008 // Define this as a plug-in
2009 #include "FWCore/Framework/interface/MakerMacros.h"
2010 DEFINE_FWK_MODULE(StandaloneTrackMonitor);