Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 /*
0003     This is the DQM code for UE physics plots
0004     11/12/2009 Sunil Bansal
0005 */
0006 
0007 #ifndef QcdUeDQM_H
0008 #define QcdUeDQM_H
0009 
0010 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0011 #include "DataFormats/Common/interface/Handle.h"
0012 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
0013 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
0014 #include "FWCore/Framework/interface/Frameworkfwd.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfo.h"
0018 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfoTrackAssociation.h"
0019 #include "DataFormats/TrackReco/interface/Track.h"
0020 #include "DataFormats/Candidate/interface/CandidateFwd.h"
0021 #include "DataFormats/Candidate/interface/CandMatchMap.h"
0022 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
0023 #include "DataFormats/JetReco/interface/TrackJetCollection.h"
0024 #include "DataFormats/VertexReco/interface/Vertex.h"
0025 #include "DataFormats/VertexReco/interface/VertexFwd.h"
0026 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
0027 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0028 #include "DataFormats/Common/interface/TriggerResults.h"
0029 #include "FWCore/Common/interface/TriggerNames.h"
0030 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0031 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
0032 #include "DQMServices/Core/interface/DQMStore.h"
0033 
0034 #include <TMath.h>
0035 #include <vector>
0036 
0037 #define PI 3.141592654
0038 
0039 class TrackerGeometry;
0040 class TH1F;
0041 class TH2F;
0042 class TH3F;
0043 class TProfile;
0044 
0045 class PtSorter {
0046 public:
0047   template <class T>
0048   bool operator()(const T &a, const T &b) {
0049     return (a.pt() > b.pt());
0050   }
0051 };
0052 
0053 class QcdUeDQM : public DQMEDAnalyzer {
0054 public:
0055   QcdUeDQM(const edm::ParameterSet &parameters);
0056   ~QcdUeDQM() override;
0057   void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override;
0058   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0059   void analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) override;
0060 
0061 private:
0062   bool isHltConfigSuccessful_;  // to prevent processing in case of problems
0063 
0064   void book1D(DQMStore::IBooker &,
0065               std::vector<MonitorElement *> &mes,
0066               const std::string &name,
0067               const std::string &title,
0068               int nx,
0069               double x1,
0070               double x2,
0071               bool sumw2 = true,
0072               bool sbox = true);
0073   void bookProfile(DQMStore::IBooker &,
0074                    std::vector<MonitorElement *> &mes,
0075                    const std::string &name,
0076                    const std::string &title,
0077                    int nx,
0078                    double x1,
0079                    double x2,
0080                    double y1,
0081                    double y2,
0082                    bool sumw2 = true,
0083                    bool sbox = true);
0084   void fill1D(std::vector<TH1F *> &hs, double val, double w = 1.);
0085   void fill1D(std::vector<MonitorElement *> &mes, double val, double w = 1.);
0086   void fill2D(std::vector<TH2F *> &hs, double valx, double valy, double w = 1.);
0087   void fill2D(std::vector<MonitorElement *> &mes, double valx, double valy, double w = 1.);
0088   void fillProfile(std::vector<TProfile *> &hs, double valx, double valy, double w = 1.);
0089   void fillProfile(std::vector<MonitorElement *> &mes, double valx, double valy, double w = 1.);
0090   void fill3D(std::vector<TH3F *> &hs, int gbin, double w = 1.);
0091   void setLabel1D(std::vector<MonitorElement *> &mes);
0092 
0093   bool trackSelection(const reco::Track &trk, const reco::BeamSpot *bs, const reco::Vertex &vtx, int sizevtx);
0094   void fillHltBits(const edm::Event &iEvent, const edm::EventSetup &iSetup);
0095   bool fillVtxPlots(const reco::BeamSpot *bs, const edm::Handle<reco::VertexCollection> vtxColl);
0096   void fillpTMaxRelated(const std::vector<const reco::Track *> &track);
0097   void fillChargedJetSpectra(const edm::Handle<reco::TrackJetCollection> trackJets);
0098   void fillUE_with_ChargedJets(const std::vector<const reco::Track *> &track,
0099                                const edm::Handle<reco::TrackJetCollection> &trackJets);
0100   void fillUE_with_MaxpTtrack(const std::vector<const reco::Track *> &track);
0101 
0102   template <typename TYPE>
0103   void getProduct(const std::string name, edm::Handle<TYPE> &prod, const edm::Event &event) const;
0104   template <typename TYPE>
0105   bool getProductSafe(const std::string name, edm::Handle<TYPE> &prod, const edm::Event &event) const;
0106 
0107   HLTConfigProvider hltConfig;
0108 
0109   std::string hltResName_;                 // HLT trigger results name
0110   std::vector<std::string> hltProcNames_;  // HLT process name(s)
0111   std::vector<std::string> hltTrgNames_;   // HLT trigger name(s)
0112 
0113   std::vector<int> hltTrgBits_;               // HLT trigger bit(s)
0114   std::vector<bool> hltTrgDeci_;              // HLT trigger descision(s)
0115   std::vector<std::string> hltTrgUsedNames_;  // HLT used trigger name(s)
0116   std::string hltUsedResName_;                // used HLT trigger results name
0117   int verbose_;                               // verbosity (0=debug,1=warn,2=error,3=throw)
0118   const TrackerGeometry *tgeo_;               // tracker geometry
0119   MonitorElement *repSumMap_;                 // report summary map
0120   MonitorElement *repSummary_;                // report summary
0121   MonitorElement *h2TrigCorr_;                // trigger correlation plot
0122 
0123   std::vector<MonitorElement *> hNevts_;
0124   std::vector<MonitorElement *> hNtrackerLayer_;
0125   std::vector<MonitorElement *> hNtrackerPixelLayer_;
0126   std::vector<MonitorElement *> hNtrackerStripPixelLayer_;
0127   std::vector<MonitorElement *> hRatioPtErrorPt_;
0128   std::vector<MonitorElement *> hTrkPt_;
0129   std::vector<MonitorElement *> hTrkEta_;
0130   std::vector<MonitorElement *> hTrkPhi_;
0131   std::vector<MonitorElement *> hNgoodTrk_;
0132   std::vector<MonitorElement *> hGoodTrkPt500_;
0133   std::vector<MonitorElement *> hGoodTrkEta500_;
0134   std::vector<MonitorElement *> hGoodTrkPhi500_;
0135   std::vector<MonitorElement *> hGoodTrkPt900_;
0136   std::vector<MonitorElement *> hGoodTrkEta900_;
0137   std::vector<MonitorElement *> hGoodTrkPhi900_;
0138   std::vector<MonitorElement *> hRatioDxySigmaDxyBS_;
0139   std::vector<MonitorElement *> hRatioDxySigmaDxyPV_;
0140   std::vector<MonitorElement *> hRatioDzSigmaDzBS_;
0141   std::vector<MonitorElement *> hRatioDzSigmaDzPV_;
0142   std::vector<MonitorElement *> hTrkChi2_;
0143   std::vector<MonitorElement *> hTrkNdof_;
0144 
0145   std::vector<MonitorElement *> hNvertices_;  // Number of vertices
0146   std::vector<MonitorElement *> hVertex_z_;   // z-position of vertex
0147   std::vector<MonitorElement *> hVertex_y_;
0148   std::vector<MonitorElement *> hVertex_x_;
0149   std::vector<MonitorElement *> hVertex_ndof_;
0150   std::vector<MonitorElement *> hVertex_rho_;
0151   std::vector<MonitorElement *> hVertex_z_bs_;
0152 
0153   std::vector<MonitorElement *> hBeamSpot_z_;  // z-position of vertex
0154   std::vector<MonitorElement *> hBeamSpot_y_;
0155   std::vector<MonitorElement *> hBeamSpot_x_;
0156 
0157   std::vector<MonitorElement *> hLeadingTrack_pTSpectrum_;   // pt spectrum of
0158                                                              // leading track
0159   std::vector<MonitorElement *> hLeadingTrack_etaSpectrum_;  // eta spectrum of
0160                                                              // leading track
0161   std::vector<MonitorElement *> hLeadingTrack_phiSpectrum_;  // phi spectrum of
0162                                                              // leading track
0163 
0164   std::vector<MonitorElement *> hChargedJetMulti_;               // Number of charged jets
0165   std::vector<MonitorElement *> hChargedJetConstituent_;         // Number of
0166                                                                  // constituent of
0167                                                                  // charged jets
0168   std::vector<MonitorElement *> hLeadingChargedJet_pTSpectrum_;  // pT spectrum
0169                                                                  // of charged
0170                                                                  // jets
0171 
0172   /*std::vector<MonitorElement*>   hCaloJetMulti_;         // Number of calo
0173   jets
0174   std::vector<MonitorElement*>   hCaloJetConstituent_;         // Number of
0175   constituent of calo jets
0176   std::vector<MonitorElement*>   hLeadingCaloJet_pTSpectrum_;         // pT
0177   spectrum of calo jets
0178   std::vector<MonitorElement*>  hLeadingCaloJet_etaSpectrum_;      //eta
0179   spectrum of leading calo jet
0180   std::vector<MonitorElement*>  hLeadingCaloJet_phiSpectrum_;      //phi
0181   spectrum of leading calo jet
0182  */
0183   std::vector<MonitorElement *> hdPhi_maxpTTrack_tracks_;  // delta phi between
0184                                                            // leading track and
0185                                                            // tracks
0186   // std::vector<MonitorElement*>  hdPhi_caloJet_tracks_;  // delta phi between
0187   // leading calo jet and tracks
0188   std::vector<MonitorElement *> hdPhi_chargedJet_tracks_;  // delta phi between
0189                                                            // leading charged
0190                                                            // jet and tracks
0191 
0192   std::vector<MonitorElement *> hLeadingChargedJet_etaSpectrum_;  // eta
0193                                                                   // spectrum of
0194                                                                   // leading
0195                                                                   // charged jet
0196   std::vector<MonitorElement *> hLeadingChargedJet_phiSpectrum_;  // phi
0197                                                                   // spectrum of
0198                                                                   // leading
0199                                                                   // charged jet
0200 
0201   std::vector<MonitorElement *> hdNdEtadPhi_pTMax_Toward500_;      // number of
0202                                                                    // tracks in
0203                                                                    // toward region
0204                                                                    // of leadin
0205                                                                    // track (pT >
0206                                                                    // 500 MeV)
0207   std::vector<MonitorElement *> hdNdEtadPhi_pTMax_Transverse500_;  // number of
0208                                                                    // tracks in
0209                                                                    // transverse
0210                                                                    // region of
0211                                                                    // leadin
0212                                                                    // track (pT
0213                                                                    // > 500 MeV)
0214   std::vector<MonitorElement *> hdNdEtadPhi_pTMax_Away500_;        // number of tracks
0215                                                                    // in away region
0216                                                                    // of leadin track
0217                                                                    // (pT > 500 MeV)
0218   /* std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Toward500_;    // number
0219    of tracks in toward region of leadin calo Jet (pT > 500 MeV)
0220    std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Transverse500_;    //
0221    number of tracks in transverse region of leadin calo Jet (pT > 500 MeV)
0222    std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Away500_;    // number of
0223    tracks in away region of leadin calo Jet (pT > 500 MeV)
0224 */
0225   std::vector<MonitorElement *> hdNdEtadPhi_trackJet_Toward500_;      // number of
0226                                                                       // tracks in
0227                                                                       // toward
0228                                                                       // region of
0229                                                                       // leadin calo
0230                                                                       // Jet (pT >
0231                                                                       // 500 MeV)
0232   std::vector<MonitorElement *> hdNdEtadPhi_trackJet_Transverse500_;  // number of tracks in transverse
0233                                                                       // region of leadin calo Jet (pT >
0234                                                                       // 500 MeV)
0235   std::vector<MonitorElement *> hdNdEtadPhi_trackJet_Away500_;        // number of
0236                                                                       // tracks in
0237                                                                       // away region
0238                                                                       // of leadin
0239                                                                       // calo Jet (pT
0240                                                                       // > 500 MeV)
0241 
0242   std::vector<MonitorElement *> hpTSumdEtadPhi_pTMax_Toward500_;      // pT sum of
0243                                                                       // tracks in
0244                                                                       // toward
0245                                                                       // region of
0246                                                                       // leadin
0247                                                                       // track (pT >
0248                                                                       // 500 MeV)
0249   std::vector<MonitorElement *> hpTSumdEtadPhi_pTMax_Transverse500_;  // pT sum of tracks in transverse
0250                                                                       // region of leadin track (pT > 500
0251                                                                       // MeV)
0252   std::vector<MonitorElement *> hpTSumdEtadPhi_pTMax_Away500_;        // pT sum of
0253                                                                       // tracks in
0254                                                                       // away region
0255                                                                       // of leadin
0256                                                                       // track (pT >
0257                                                                       // 500 MeV)
0258   /*  std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Toward500_;    // pT
0259     sum of tracks in toward region of leadin calo Jet (pT > 500 MeV)
0260     std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Transverse500_;    //
0261     pT sum of tracks in transverse region of leadin calo Jet (pT > 500 MeV)
0262     std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Away500_;    // pT sum
0263     of tracks in away region of leadin calo Jet (pT > 500 MeV)
0264 */
0265   std::vector<MonitorElement *> hpTSumdEtadPhi_trackJet_Toward500_;      // pT sum
0266                                                                          // of
0267                                                                          // tracks
0268                                                                          // in
0269                                                                          // toward
0270                                                                          // region
0271                                                                          // of
0272                                                                          // leadin
0273                                                                          // calo Jet
0274                                                                          // (pT >
0275                                                                          // 500 MeV)
0276   std::vector<MonitorElement *> hpTSumdEtadPhi_trackJet_Transverse500_;  // pT sum of tracks in transverse
0277                                                                          // region of leadin calo Jet (pT
0278                                                                          // > 500 MeV)
0279   std::vector<MonitorElement *> hpTSumdEtadPhi_trackJet_Away500_;        // pT sum of
0280                                                                          // tracks in
0281                                                                          // away
0282                                                                          // region of
0283                                                                          // leadin
0284                                                                          // calo Jet
0285                                                                          // (pT > 500
0286                                                                          // MeV)
0287 
0288   std::vector<MonitorElement *> hdNdEtadPhi_pTMax_Toward900_;      // number of
0289                                                                    // tracks in
0290                                                                    // toward region
0291                                                                    // of leadin
0292                                                                    // track (pT >
0293                                                                    // 900 MeV)
0294   std::vector<MonitorElement *> hdNdEtadPhi_pTMax_Transverse900_;  // number of
0295                                                                    // tracks in
0296                                                                    // transverse
0297                                                                    // region of
0298                                                                    // leadin
0299                                                                    // track (pT
0300                                                                    // > 900 MeV)
0301   std::vector<MonitorElement *> hdNdEtadPhi_pTMax_Away900_;        // number of tracks
0302                                                                    // in away region
0303                                                                    // of leadin track
0304                                                                    // (pT > 900 MeV)
0305   /*  std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Toward900_;    // number
0306     of tracks in toward region of leadin calo Jet (pT > 900 MeV)
0307     std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Transverse900_;    //
0308     number of tracks in transverse region of leadin calo Jet (pT > 900 MeV)
0309     std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Away900_;    // number of
0310     tracks in away region of leadin calo Jet (pT > 900 MeV)
0311 */
0312   std::vector<MonitorElement *> hdNdEtadPhi_trackJet_Toward900_;      // number of
0313                                                                       // tracks in
0314                                                                       // toward
0315                                                                       // region of
0316                                                                       // leadin calo
0317                                                                       // Jet (pT >
0318                                                                       // 900 MeV)
0319   std::vector<MonitorElement *> hdNdEtadPhi_trackJet_Transverse900_;  // number of tracks in transverse
0320                                                                       // region of leadin calo Jet (pT >
0321                                                                       // 900 MeV)
0322   std::vector<MonitorElement *> hdNdEtadPhi_trackJet_Away900_;        // number of
0323                                                                       // tracks in
0324                                                                       // away region
0325                                                                       // of leadin
0326                                                                       // calo Jet (pT
0327                                                                       // > 900 MeV)
0328 
0329   std::vector<MonitorElement *> hpTSumdEtadPhi_pTMax_Toward900_;      // pT sum of
0330                                                                       // tracks in
0331                                                                       // toward
0332                                                                       // region of
0333                                                                       // leadin
0334                                                                       // track (pT >
0335                                                                       // 900 MeV)
0336   std::vector<MonitorElement *> hpTSumdEtadPhi_pTMax_Transverse900_;  // pT sum of tracks in transverse
0337                                                                       // region of leadin track (pT > 900
0338                                                                       // MeV)
0339   std::vector<MonitorElement *> hpTSumdEtadPhi_pTMax_Away900_;        // pT sum of
0340                                                                       // tracks in
0341                                                                       // away region
0342                                                                       // of leadin
0343                                                                       // track (pT >
0344                                                                       // 900 MeV)
0345   /*  std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Toward900_;    // pT
0346     sum of tracks in toward region of leadin calo Jet (pT > 900 MeV)
0347     std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Transverse900_;    //
0348     pT sum of tracks in transverse region of leadin calo Jet (pT > 900 MeV)
0349     std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Away900_;    // pT sum
0350     of tracks in away region of leadin calo Jet (pT > 900 MeV)
0351 */
0352   std::vector<MonitorElement *> hpTSumdEtadPhi_trackJet_Toward900_;      // pT sum
0353                                                                          // of
0354                                                                          // tracks
0355                                                                          // in
0356                                                                          // toward
0357                                                                          // region
0358                                                                          // of
0359                                                                          // leadin
0360                                                                          // calo Jet
0361                                                                          // (pT >
0362                                                                          // 900 MeV)
0363   std::vector<MonitorElement *> hpTSumdEtadPhi_trackJet_Transverse900_;  // pT sum of tracks in transverse
0364                                                                          // region of leadin calo Jet (pT
0365                                                                          // > 900 MeV)
0366   std::vector<MonitorElement *> hpTSumdEtadPhi_trackJet_Away900_;        // pT sum of
0367                                                                          // tracks in
0368                                                                          // away
0369                                                                          // region of
0370                                                                          // leadin
0371                                                                          // calo Jet
0372                                                                          // (pT > 900
0373                                                                          // MeV)
0374 
0375   double ptMin_;
0376   double minRapidity_;
0377   double maxRapidity_;
0378   double tip_;
0379   double lip_;
0380   double diffvtxbs_;
0381   double ptErr_pt_;
0382   double vtxntk_;
0383   int minHit_;
0384   double pxlLayerMinCut_;
0385   bool requirePIX1_;
0386   int min3DHit_;
0387   double maxChi2_;
0388   bool bsuse_;
0389   bool allowTriplets_;
0390   double bsPos_;
0391   edm::EDGetTokenT<reco::CaloJetCollection> caloJetLabel_;
0392   edm::EDGetTokenT<reco::TrackJetCollection> chargedJetLabel_;
0393   edm::EDGetTokenT<reco::TrackCollection> trackLabel_;
0394   edm::EDGetTokenT<reco::VertexCollection> vtxLabel_;
0395   edm::EDGetTokenT<reco::BeamSpot> bsLabel_;
0396   std::vector<reco::TrackBase::TrackQuality> quality_;
0397   std::vector<reco::TrackBase::TrackAlgorithm> algorithm_;
0398   typedef std::vector<const reco::Track *> container;
0399   container selected_;
0400   reco::Vertex vtx1;
0401 };
0402 
0403 //--------------------------------------------------------------------------------------------------
0404 template <typename TYPE>
0405 inline void QcdUeDQM::getProduct(const std::string name, edm::Handle<TYPE> &prod, const edm::Event &event) const {
0406   // Try to access data collection from EDM file. We check if we really get just
0407   // one
0408   // product with the given name. If not we throw an exception.
0409 
0410   event.getByLabel(edm::InputTag(name), prod);
0411   if (!prod.isValid())
0412     throw edm::Exception(edm::errors::Configuration, "QcdUeDQM::GetProduct()\n")
0413         << "Collection with label " << name << " is not valid" << std::endl;
0414 }
0415 
0416 //--------------------------------------------------------------------------------------------------
0417 template <typename TYPE>
0418 inline bool QcdUeDQM::getProductSafe(const std::string name, edm::Handle<TYPE> &prod, const edm::Event &event) const {
0419   // Try to safely access data collection from EDM file. We check if we really
0420   // get just one
0421   // product with the given name. If not, we return false.
0422 
0423   if (name.empty())
0424     return false;
0425 
0426   try {
0427     event.getByLabel(edm::InputTag(name), prod);
0428     if (!prod.isValid())
0429       return false;
0430   } catch (...) {
0431     return false;
0432   }
0433   return true;
0434 }
0435 
0436 //--------------------------------------------------------------------------------------------------
0437 #endif
0438 
0439 // Local Variables:
0440 // show-trailing-whitespace: t
0441 // truncate-lines: t
0442 // End: