Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-11-17 23:29:39

0001 #ifndef DQMOffline_RecoB_PrimaryVertexMonitor_H
0002 #define DQMOffline_RecoB_PrimaryVertexMonitor_H
0003 
0004 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0005 #include "DQMServices/Core/interface/DQMStore.h"
0006 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0007 #include "DataFormats/Common/interface/Association.h"
0008 #include "DataFormats/Common/interface/ValueMap.h"
0009 #include "DataFormats/VertexReco/interface/Vertex.h"
0010 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/Frameworkfwd.h"
0013 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/Utilities/interface/EDGetToken.h"
0016 #include "FWCore/Utilities/interface/InputTag.h"
0017 
0018 /** \class PrimaryVertexMonitor
0019  *
0020  *
0021  */
0022 
0023 namespace pvMonitor {
0024   // same logic used for the MTV:
0025   // cf https://github.com/cms-sw/cmssw/blob/master/Validation/RecoTrack/src/MTVHistoProducerAlgoForTracker.cc
0026   typedef dqm::reco::DQMStore DQMStore;
0027 
0028   void setBinLog(TAxis *axis) {
0029     int bins = axis->GetNbins();
0030     float from = axis->GetXmin();
0031     float to = axis->GetXmax();
0032     float width = (to - from) / bins;
0033     std::vector<float> new_bins(bins + 1, 0);
0034     for (int i = 0; i <= bins; i++) {
0035       new_bins[i] = TMath::Power(10, from + i * width);
0036     }
0037     axis->Set(bins, new_bins.data());
0038   }
0039 
0040   void setBinLogX(TH1 *h) {
0041     TAxis *axis = h->GetXaxis();
0042     setBinLog(axis);
0043   }
0044   void setBinLogY(TH1 *h) {
0045     TAxis *axis = h->GetYaxis();
0046     setBinLog(axis);
0047   }
0048 
0049   template <typename... Args>
0050   dqm::reco::MonitorElement *makeProfileIfLog(DQMStore::IBooker &ibook, bool logx, bool logy, Args &&...args) {
0051     auto prof = std::make_unique<TProfile>(std::forward<Args>(args)...);
0052     if (logx)
0053       setBinLogX(prof.get());
0054     if (logy)
0055       setBinLogY(prof.get());
0056     const auto &name = prof->GetName();
0057     return ibook.bookProfile(name, prof.release());
0058   }
0059 
0060   template <typename... Args>
0061   dqm::reco::MonitorElement *makeTH1IfLog(DQMStore::IBooker &ibook, bool logx, bool logy, Args &&...args) {
0062     auto h1 = std::make_unique<TH1F>(std::forward<Args>(args)...);
0063     if (logx)
0064       setBinLogX(h1.get());
0065     if (logy)
0066       setBinLogY(h1.get());
0067     const auto &name = h1->GetName();
0068     return ibook.book1D(name, h1.release());
0069   }
0070 
0071 }  // namespace pvMonitor
0072 
0073 class PrimaryVertexMonitor : public DQMEDAnalyzer {
0074 public:
0075   explicit PrimaryVertexMonitor(const edm::ParameterSet &pSet);
0076   ~PrimaryVertexMonitor() override = default;
0077 
0078   static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0079 
0080   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0081   void analyze(const edm::Event &, const edm::EventSetup &) override;
0082 
0083   struct IPMonitoring {
0084     std::string varname_;
0085     float pTcut_;
0086     dqm::reco::MonitorElement *IP_, *IPErr_, *IPPull_;
0087     dqm::reco::MonitorElement *IPVsPhi_, *IPVsEta_, *IPVsPt_;
0088     dqm::reco::MonitorElement *IPErrVsPhi_, *IPErrVsEta_, *IPErrVsPt_;
0089     dqm::reco::MonitorElement *IPVsEtaVsPhi_, *IPErrVsEtaVsPhi_;
0090 
0091     void bookIPMonitor(DQMStore::IBooker &, const edm::ParameterSet &);
0092 
0093   private:
0094     int PhiBin_, EtaBin_, PtBin_;
0095     double PhiMin_, PhiMax_, EtaMin_, EtaMax_, PtMin_, PtMax_;
0096   };
0097 
0098 private:
0099   void pvTracksPlots(const reco::Vertex &v);
0100   void vertexPlots(const reco::Vertex &v, const reco::BeamSpot &beamSpot, int i);
0101 
0102   // event data
0103 
0104   const edm::InputTag vertexInputTag_;
0105   const edm::InputTag beamSpotInputTag_;
0106   const edm::EDGetTokenT<reco::VertexCollection> vertexToken_;
0107   using VertexScore = edm::ValueMap<float>;
0108   const edm::EDGetTokenT<VertexScore> scoreToken_;
0109   const edm::EDGetTokenT<reco::BeamSpot> beamspotToken_;
0110 
0111   // configuration
0112 
0113   const edm::ParameterSet conf_;
0114   const std::string dqmLabel;
0115   const std::string TopFolderName_;
0116   const std::string AlignmentLabel_;
0117   const int ndof_;
0118   const bool useHPfoAlignmentPlots_;
0119   bool errorPrinted_;
0120 
0121   static constexpr int cmToUm = 10000;
0122 
0123   // the histos
0124   MonitorElement *nbvtx, *nbgvtx, *nbtksinvtx[2], *trksWeight[2], *score[2];
0125   MonitorElement *tt[2];
0126   MonitorElement *xrec[2], *yrec[2], *zrec[2], *xDiff[2], *yDiff[2], *xerr[2], *yerr[2], *zerr[2];
0127   MonitorElement *xerrVsTrks[2], *yerrVsTrks[2], *zerrVsTrks[2];
0128   MonitorElement *ntracksVsZ[2];
0129   MonitorElement *vtxchi2[2], *vtxndf[2], *vtxprob[2], *nans[2];
0130   MonitorElement *type[2];
0131   MonitorElement *bsX, *bsY, *bsZ, *bsSigmaZ, *bsDxdz, *bsDydz, *bsBeamWidthX, *bsBeamWidthY, *bsType;
0132 
0133   MonitorElement *sumpt, *ntracks, *weight, *chi2ndf, *chi2prob;
0134   MonitorElement *trackpt;
0135   MonitorElement *phi_pt1, *eta_pt1;
0136   MonitorElement *phi_pt10, *eta_pt10;
0137 
0138   MonitorElement *dxy2;
0139 
0140   // IP monitoring structs
0141   IPMonitoring dxy_pt1;
0142   IPMonitoring dxy_pt10;
0143 
0144   IPMonitoring dz_pt1;
0145   IPMonitoring dz_pt10;
0146 };
0147 
0148 #endif