CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10/src/DQM/Physics/src/QcdUeDQM.h

Go to the documentation of this file.
00001 
00002 /*
00003     This is the DQM code for UE physics plots
00004     11/12/2009 Sunil Bansal
00005 */
00006 
00007 #ifndef QcdUeDQM_H
00008 #define QcdUeDQM_H
00009 
00010 #include "DataFormats/Common/interface/Handle.h"
00011 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
00012 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00013 #include "FWCore/Framework/interface/Frameworkfwd.h"
00014 #include "FWCore/Framework/interface/EDAnalyzer.h"
00015 #include "FWCore/Framework/interface/Event.h"
00016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00017 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfo.h"
00018 #include "AnalysisDataFormats/TrackInfo/interface/TrackInfoTrackAssociation.h"
00019 #include "DataFormats/TrackReco/interface/Track.h"
00020 #include "DataFormats/Candidate/interface/CandidateFwd.h"
00021 #include "DataFormats/Candidate/interface/CandMatchMap.h"
00022 #include "DataFormats/JetReco/interface/CaloJetCollection.h"
00023 #include "DataFormats/JetReco/interface/TrackJetCollection.h"
00024 #include "DataFormats/VertexReco/interface/Vertex.h"
00025 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00026 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h" 
00027 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00028 #include "DataFormats/Common/interface/TriggerResults.h"
00029 #include "FWCore/Common/interface/TriggerNames.h"
00030 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
00031 #include "DataFormats/HLTReco/interface/TriggerEvent.h"
00032 #include <TMath.h>
00033 #include <vector>
00034 #define PI 3.141592654
00035 class DQMStore;
00036 class MonitorElement;
00037 class TrackerGeometry;
00038 class TH1F;
00039 class TH2F;
00040 class TH3F;
00041 class TProfile;
00042 
00043 class PtSorter {
00044 public:
00045   template <class T> bool operator() ( const T& a, const T& b ) {
00046     return ( a.pt() > b.pt() );
00047   }
00048 };
00049 
00050 
00051 class QcdUeDQM : public edm::EDAnalyzer 
00052 {
00053   public:
00054 
00055     QcdUeDQM(const edm::ParameterSet &parameters);
00056     virtual ~QcdUeDQM();
00057 
00058     void                          analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup);
00059     void                          beginJob(void);
00060     void                          beginLuminosityBlock(const edm::LuminosityBlock &l, 
00061                                                        const edm::EventSetup &iSetup);
00062     void                          beginRun(const edm::Run &r, const edm::EventSetup &iSetup);
00063     void                          endJob(void);
00064     void                          endRun(const edm::Run &r, const edm::EventSetup &iSetup);
00065     void                          endLuminosityBlock(const edm::LuminosityBlock &l, 
00066                                                      const edm::EventSetup &iSetup);
00067 
00068   private:
00069     bool isHltConfigSuccessful_; // to prevent processing in case of problems
00070 
00071     void                          book1D(std::vector<MonitorElement*> &mes, 
00072                                          const std::string &name, const std::string &title, 
00073                                          int nx, double x1, double x2, bool sumw2=1, bool sbox=1);
00074     void                          book2D(std::vector<MonitorElement*> &mes, 
00075                                          const std::string &name, const std::string &title, 
00076                                          int nx, double x1, double x2, int ny, double y1, double y2,
00077                                          bool sumw2=1, bool sbox=1);
00078     void                          bookProfile(std::vector<MonitorElement*> &mes, 
00079                                          const std::string &name, const std::string &title, 
00080                                          int nx, double x1, double x2,  double y1, double y2,
00081                                          bool sumw2=1, bool sbox=1); 
00082     void                          create1D(std::vector<TH1F*> &mes, 
00083                                            const std::string &name, const std::string &title, 
00084                                            int nx, double x1, double x2, bool sumw2=1, bool sbox=1);
00085     void                          create2D(std::vector<TH2F*> &mes, 
00086                                            const std::string &name, const std::string &title, 
00087                                            int nx, double x1, double x2, int ny, double y1, double y2,
00088                                            bool sumw2=1, bool sbox=1);
00089     void                          createProfile(std::vector<TProfile*> &mes, 
00090                                            const std::string &name, const std::string &title, 
00091                                            int nx, double x1, double x2,  double y1, double y2,
00092                                            bool sumw2=1, bool sbox=1);  
00093     void                          createHistos();
00094     void                          fill1D(std::vector<TH1F*> &hs, double val, double w=1.);
00095     void                          fill1D(std::vector<MonitorElement*> &mes, double val, double w=1.);
00096     void                          fill2D(std::vector<TH2F*> &hs, 
00097                                          double valx, double valy, double w=1.);
00098     void                          fill2D(std::vector<MonitorElement*> &mes, 
00099                                          double valx, double valy, double w=1.);
00100     void                          fillProfile(std::vector<TProfile*> &hs, 
00101                                          double valx, double valy, double w=1.);
00102     void                          fillProfile(std::vector<MonitorElement*> &mes, 
00103                                         double valx, double valy, double w=1.);
00104     void                          fill3D(std::vector<TH3F*> &hs, int gbin, double w=1.);
00105     void                          setLabel1D(std::vector<MonitorElement*> &mes);
00106    
00107 
00108     bool                          trackSelection(const reco::Track &trk, const reco::BeamSpot* bs, const reco::Vertex& vtx, int sizevtx);
00109     void                          fillHltBits(const edm::Event &iEvent,const edm::EventSetup &iSetup);
00110     
00111     bool                          fillVtxPlots(const reco::BeamSpot* bs, const edm::Handle< reco::VertexCollection > vtxColl);
00112     void                          fillpTMaxRelated(const std::vector<const reco::Track *>  &track);
00113 
00114     void                          fillChargedJetSpectra( const edm::Handle<reco::TrackJetCollection> trackJets);
00115   
00116     //void                          fillCaloJetSpectra(const edm::Handle<reco::CaloJetCollection> caloJets);
00117 
00118     void                          fillUE_with_ChargedJets(const std::vector<const reco::Track *>  &track,  const edm::Handle<reco::TrackJetCollection> &trackJets);
00119     //void                          fillUE_with_CaloJets(const std::vector<const reco::Track *>  &track, const edm::Handle<reco::CaloJetCollection> &caloJets);
00120     void                          fillUE_with_MaxpTtrack(const std::vector<const reco::Track *>  &track ); 
00121 
00122    
00123     template <typename TYPE>
00124     void                          getProduct(const std::string name, edm::Handle<TYPE> &prod,
00125                                              const edm::Event &event) const;    
00126     template <typename TYPE>
00127     bool                          getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
00128                                                  const edm::Event &event) const;
00129 
00130     HLTConfigProvider hltConfig;   
00131 
00132     std::string                   hltResName_;         //HLT trigger results name
00133     std::vector<std::string>      hltProcNames_;       //HLT process name(s)
00134     std::vector<std::string>      hltTrgNames_;        //HLT trigger name(s)
00135     
00136    
00137     std::vector<int>              hltTrgBits_;         //HLT trigger bit(s)
00138     std::vector<bool>             hltTrgDeci_;         //HLT trigger descision(s)
00139     std::vector<std::string>      hltTrgUsedNames_;    //HLT used trigger name(s)
00140     std::string                   hltUsedResName_;     //used HLT trigger results name 
00141     int                           verbose_;            //verbosity (0=debug,1=warn,2=error,3=throw)
00142     const TrackerGeometry        *tgeo_;               //tracker geometry
00143     DQMStore                     *theDbe_;             //dqm store
00144     MonitorElement               *repSumMap_;          //report summary map
00145     MonitorElement               *repSummary_;         //report summary
00146     MonitorElement               *h2TrigCorr_;         //trigger correlation plot
00147 
00148     std::vector<MonitorElement*> hNevts_;
00149     std::vector<MonitorElement*> hNtrackerLayer_;
00150     std::vector<MonitorElement*> hNtrackerPixelLayer_;
00151     std::vector<MonitorElement*> hNtrackerStripPixelLayer_;
00152     std::vector<MonitorElement*> hRatioPtErrorPt_;
00153     std::vector<MonitorElement*> hTrkPt_;
00154     std::vector<MonitorElement*> hTrkEta_;  
00155     std::vector<MonitorElement*> hTrkPhi_;
00156     std::vector<MonitorElement*> hNgoodTrk_;
00157     std::vector<MonitorElement*> hGoodTrkPt500_;
00158     std::vector<MonitorElement*> hGoodTrkEta500_;
00159     std::vector<MonitorElement*> hGoodTrkPhi500_; 
00160     std::vector<MonitorElement*> hGoodTrkPt900_;
00161     std::vector<MonitorElement*> hGoodTrkEta900_;
00162     std::vector<MonitorElement*> hGoodTrkPhi900_;
00163     std::vector<MonitorElement*> hRatioDxySigmaDxyBS_;
00164     std::vector<MonitorElement*> hRatioDxySigmaDxyPV_;
00165     std::vector<MonitorElement*> hRatioDzSigmaDzBS_;
00166     std::vector<MonitorElement*> hRatioDzSigmaDzPV_;
00167     std::vector<MonitorElement*> hTrkChi2_;
00168     std::vector<MonitorElement*> hTrkNdof_; 
00169 
00170 
00171 
00172     std::vector<MonitorElement*>  hNvertices_;           // Number of vertices
00173     std::vector<MonitorElement*>  hVertex_z_;            // z-position of vertex
00174     std::vector<MonitorElement*>  hVertex_y_;
00175     std::vector<MonitorElement*>  hVertex_x_;
00176     std::vector<MonitorElement*>  hVertex_ndof_;
00177     std::vector<MonitorElement*>  hVertex_rho_;
00178     std::vector<MonitorElement*>  hVertex_z_bs_; 
00179 
00180 
00181     std::vector<MonitorElement*>  hBeamSpot_z_;            // z-position of vertex
00182     std::vector<MonitorElement*>  hBeamSpot_y_;
00183     std::vector<MonitorElement*>  hBeamSpot_x_;
00184 
00185 
00186     std::vector<MonitorElement*>  hLeadingTrack_pTSpectrum_;      //pt spectrum of leading track
00187     std::vector<MonitorElement*>  hLeadingTrack_etaSpectrum_;      //eta spectrum of leading track
00188     std::vector<MonitorElement*>  hLeadingTrack_phiSpectrum_;      //phi spectrum of leading track
00189 
00190     
00191    
00192    
00193 
00194     std::vector<MonitorElement*>   hChargedJetMulti_;         // Number of charged jets
00195     std::vector<MonitorElement*>   hChargedJetConstituent_;         // Number of constituent of charged jets 
00196     std::vector<MonitorElement*>   hLeadingChargedJet_pTSpectrum_;         // pT spectrum of charged jets 
00197 
00198     /*std::vector<MonitorElement*>   hCaloJetMulti_;         // Number of calo jets
00199     std::vector<MonitorElement*>   hCaloJetConstituent_;         // Number of constituent of calo jets
00200     std::vector<MonitorElement*>   hLeadingCaloJet_pTSpectrum_;         // pT spectrum of calo jets
00201     std::vector<MonitorElement*>  hLeadingCaloJet_etaSpectrum_;      //eta spectrum of leading calo jet
00202     std::vector<MonitorElement*>  hLeadingCaloJet_phiSpectrum_;      //phi spectrum of leading calo jet
00203    */
00204     std::vector<MonitorElement*>  hdPhi_maxpTTrack_tracks_;  // delta phi between leading track and tracks
00205     //std::vector<MonitorElement*>  hdPhi_caloJet_tracks_;  // delta phi between leading calo jet and tracks  
00206     std::vector<MonitorElement*>  hdPhi_chargedJet_tracks_;  // delta phi between leading charged jet and tracks
00207    
00208 
00209     std::vector<MonitorElement*>  hLeadingChargedJet_etaSpectrum_;      //eta spectrum of leading charged jet
00210     std::vector<MonitorElement*>  hLeadingChargedJet_phiSpectrum_;      //phi spectrum of leading charged jet
00211   
00212     std::vector<MonitorElement*>  hdNdEtadPhi_pTMax_Toward500_;    // number of tracks in toward region of leadin track (pT > 500 MeV)
00213     std::vector<MonitorElement*>  hdNdEtadPhi_pTMax_Transverse500_;    // number of tracks in transverse region of leadin track (pT > 500 MeV)  
00214     std::vector<MonitorElement*>  hdNdEtadPhi_pTMax_Away500_;    // number of tracks in away region of leadin track (pT > 500 MeV)
00215    /* std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Toward500_;    // number of tracks in toward region of leadin calo Jet (pT > 500 MeV)
00216     std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Transverse500_;    // number of tracks in transverse region of leadin calo Jet (pT > 500 MeV)  
00217     std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Away500_;    // number of tracks in away region of leadin calo Jet (pT > 500 MeV)
00218 */
00219     std::vector<MonitorElement*>  hdNdEtadPhi_trackJet_Toward500_;    // number of tracks in toward region of leadin calo Jet (pT > 500 MeV)
00220     std::vector<MonitorElement*>  hdNdEtadPhi_trackJet_Transverse500_;    // number of tracks in transverse region of leadin calo Jet (pT > 500 MeV)  
00221     std::vector<MonitorElement*>  hdNdEtadPhi_trackJet_Away500_;    // number of tracks in away region of leadin calo Jet (pT > 500 MeV) 
00222 
00223 
00224     
00225     std::vector<MonitorElement*>  hpTSumdEtadPhi_pTMax_Toward500_;    // pT sum of tracks in toward region of leadin track (pT > 500 MeV)
00226     std::vector<MonitorElement*>  hpTSumdEtadPhi_pTMax_Transverse500_;    // pT sum of tracks in transverse region of leadin track (pT > 500 MeV)  
00227     std::vector<MonitorElement*>  hpTSumdEtadPhi_pTMax_Away500_;    // pT sum of tracks in away region of leadin track (pT > 500 MeV)
00228   /*  std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Toward500_;    // pT sum of tracks in toward region of leadin calo Jet (pT > 500 MeV)
00229     std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Transverse500_;    // pT sum of tracks in transverse region of leadin calo Jet (pT > 500 MeV)  
00230     std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Away500_;    // pT sum of tracks in away region of leadin calo Jet (pT > 500 MeV)
00231 */
00232     std::vector<MonitorElement*>  hpTSumdEtadPhi_trackJet_Toward500_;    // pT sum of tracks in toward region of leadin calo Jet (pT > 500 MeV)
00233     std::vector<MonitorElement*>  hpTSumdEtadPhi_trackJet_Transverse500_;    // pT sum of tracks in transverse region of leadin calo Jet (pT > 500 MeV)  
00234     std::vector<MonitorElement*>  hpTSumdEtadPhi_trackJet_Away500_;    // pT sum of tracks in away region of leadin calo Jet (pT > 500 MeV)
00235    
00236   
00237 
00238     std::vector<MonitorElement*>  hdNdEtadPhi_pTMax_Toward900_;    // number of tracks in toward region of leadin track (pT > 900 MeV)
00239     std::vector<MonitorElement*>  hdNdEtadPhi_pTMax_Transverse900_;    // number of tracks in transverse region of leadin track (pT > 900 MeV)  
00240     std::vector<MonitorElement*>  hdNdEtadPhi_pTMax_Away900_;    // number of tracks in away region of leadin track (pT > 900 MeV)
00241   /*  std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Toward900_;    // number of tracks in toward region of leadin calo Jet (pT > 900 MeV)
00242     std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Transverse900_;    // number of tracks in transverse region of leadin calo Jet (pT > 900 MeV)  
00243     std::vector<MonitorElement*>  hdNdEtadPhi_caloJet_Away900_;    // number of tracks in away region of leadin calo Jet (pT > 900 MeV)
00244 */
00245     std::vector<MonitorElement*>  hdNdEtadPhi_trackJet_Toward900_;    // number of tracks in toward region of leadin calo Jet (pT > 900 MeV)
00246     std::vector<MonitorElement*>  hdNdEtadPhi_trackJet_Transverse900_;    // number of tracks in transverse region of leadin calo Jet (pT > 900 MeV)  
00247     std::vector<MonitorElement*>  hdNdEtadPhi_trackJet_Away900_;    // number of tracks in away region of leadin calo Jet (pT > 900 MeV) 
00248 
00249 
00250 
00251     std::vector<MonitorElement*>  hpTSumdEtadPhi_pTMax_Toward900_;    // pT sum of tracks in toward region of leadin track (pT > 900 MeV)
00252     std::vector<MonitorElement*>  hpTSumdEtadPhi_pTMax_Transverse900_;    // pT sum of tracks in transverse region of leadin track (pT > 900 MeV)  
00253     std::vector<MonitorElement*>  hpTSumdEtadPhi_pTMax_Away900_;    // pT sum of tracks in away region of leadin track (pT > 900 MeV)
00254   /*  std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Toward900_;    // pT sum of tracks in toward region of leadin calo Jet (pT > 900 MeV)
00255     std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Transverse900_;    // pT sum of tracks in transverse region of leadin calo Jet (pT > 900 MeV)  
00256     std::vector<MonitorElement*>  hpTSumdEtadPhi_caloJet_Away900_;    // pT sum of tracks in away region of leadin calo Jet (pT > 900 MeV)
00257 */
00258     std::vector<MonitorElement*>  hpTSumdEtadPhi_trackJet_Toward900_;    // pT sum of tracks in toward region of leadin calo Jet (pT > 900 MeV)
00259     std::vector<MonitorElement*>  hpTSumdEtadPhi_trackJet_Transverse900_;    // pT sum of tracks in transverse region of leadin calo Jet (pT > 900 MeV)  
00260     std::vector<MonitorElement*>  hpTSumdEtadPhi_trackJet_Away900_;    // pT sum of tracks in away region of leadin calo Jet (pT > 900 MeV)  
00261 
00262  
00263    
00264     double ptMin_;
00265     double minRapidity_;
00266     double maxRapidity_;
00267     double tip_;
00268     double lip_;
00269     double diffvtxbs_;
00270     double ptErr_pt_;
00271     double vtxntk_;
00272     int minHit_;
00273     double pxlLayerMinCut_;
00274     bool requirePIX1_;
00275     int min3DHit_;
00276     double maxChi2_;
00277     bool bsuse_;
00278     bool allowTriplets_;
00279     double bsPos_;
00280     edm::InputTag caloJetLabel_;
00281     edm::InputTag chargedJetLabel_;
00282     edm::InputTag trackLabel_;
00283     edm::InputTag vtxLabel_;
00284     edm::InputTag bsLabel_;
00285     std::vector<reco::TrackBase::TrackQuality> quality_; 
00286     std::vector<reco::TrackBase::TrackAlgorithm> algorithm_;
00287     typedef std::vector<const reco::Track *> container;
00288     container selected_;
00289     reco::Vertex vtx1;
00290     
00291 };
00292 
00293 //--------------------------------------------------------------------------------------------------
00294 template <typename TYPE>
00295 inline void QcdUeDQM::getProduct(const std::string name, edm::Handle<TYPE> &prod,
00296                                     const edm::Event &event) const
00297 {
00298   // Try to access data collection from EDM file. We check if we really get just one
00299   // product with the given name. If not we throw an exception.
00300 
00301   event.getByLabel(edm::InputTag(name),prod);
00302   if (!prod.isValid()) 
00303     throw edm::Exception(edm::errors::Configuration, "QcdUeDQM::GetProduct()\n")
00304       << "Collection with label " << name << " is not valid" <<  std::endl;
00305 }
00306 
00307 //--------------------------------------------------------------------------------------------------
00308 template <typename TYPE>
00309 inline bool QcdUeDQM::getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
00310                                         const edm::Event &event) const
00311 {
00312   // Try to safely access data collection from EDM file. We check if we really get just one
00313   // product with the given name. If not, we return false.
00314 
00315   if (name.size()==0)
00316     return false;
00317 
00318   try {
00319     event.getByLabel(edm::InputTag(name),prod);
00320     if (!prod.isValid()) 
00321       return false;
00322   } catch (...) {
00323     return false;
00324   }
00325   return true;
00326 }
00327 
00328 //--------------------------------------------------------------------------------------------------
00329 #endif