CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/DQM/Physics/src/QcdLowPtDQM.h

Go to the documentation of this file.
00001 // $Id: QcdLowPtDQM.h,v 1.10 2010/03/03 09:32:40 olzem Exp $
00002 
00003 #ifndef QcdLowPtDQM_H
00004 #define QcdLowPtDQM_H
00005 
00006 #include "DataFormats/Common/interface/Handle.h"
00007 #include "DataFormats/GeometryVector/interface/VectorUtil.h"
00008 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00009 #include "FWCore/Framework/interface/Frameworkfwd.h"
00010 #include "FWCore/Framework/interface/EDAnalyzer.h"
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 #include "HLTrigger/HLTcore/interface/HLTConfigProvider.h"
00014 #include <TMath.h>
00015 #include <vector>
00016 
00017 class DQMStore;
00018 class MonitorElement;
00019 class TrackerGeometry;
00020 class TH1F;
00021 class TH2F;
00022 class TH3F;
00023 
00024 class QcdLowPtDQM : public edm::EDAnalyzer 
00025 {
00026   public:
00027     class Pixel {
00028       public:
00029         Pixel(double x=0, double y=0, double z=0, double eta=0, double phi=0, 
00030               double adc=0, double sx=0, double sy=0) : 
00031           x_(x), y_(y), z_(z), rho_(TMath::Sqrt(x_*x_+y_*y_)), 
00032           eta_(eta), phi_(phi), adc_(adc), sizex_(sx), sizey_(sy) {}
00033         Pixel(const GlobalPoint &p, double adc=0, double sx=0, double sy=0) :
00034           x_(p.x()), y_(p.y()), z_(p.z()), rho_(TMath::Sqrt(x_*x_+y_*y_)), 
00035           eta_(p.eta()), phi_(p.phi()), adc_(adc), sizex_(sx), sizey_(sy) {}
00036         double adc()   const { return adc_;   }
00037         double eta()   const { return eta_;   }
00038         double rho()   const { return rho_;   }
00039         double phi()   const { return phi_;   }
00040         double sizex() const { return sizex_; }
00041         double sizey() const { return sizey_; }
00042         double x()     const { return x_;     }
00043         double y()     const { return y_;     }
00044         double z()     const { return z_;     }
00045       protected:    
00046         double x_,y_,z_,rho_,eta_,phi_;
00047         double adc_,sizex_,sizey_;
00048     };
00049     class Tracklet {
00050       public:
00051         Tracklet() : i1_(-1), i2_(-2), deta_(0), dphi_(0) {}
00052         Tracklet(const Pixel &p1, const Pixel &p2) : 
00053           p1_(p1), p2_(p2), i1_(-1), i2_(-1), 
00054           deta_(p1.eta()-p2.eta()), dphi_(Geom::deltaPhi(p1.phi(),p2.phi())) {}
00055         double       deta()  const { return deta_;     }
00056         double       dphi()  const { return dphi_;     }
00057         int          i1()    const { return i1_;       }
00058         int          i2()    const { return i2_;       }
00059         double       eta()   const { return p1_.eta(); }
00060         const Pixel &p1()    const { return p1_;       }
00061         const Pixel &p2()    const { return p2_;       }
00062         void         seti1(int i1) { i1_ = i1;         }      
00063         void         seti2(int i2) { i2_ = i2;         }      
00064       protected:
00065         Pixel  p1_,p2_;
00066         int    i1_, i2_;
00067         double deta_,dphi_;
00068     };
00069     class Vertex {
00070       public:
00071         Vertex(double x=0, double y=0, double z=0, double xs=0, double ys=0, double zs=0, int n=0) :
00072           x_(x), y_(y), z_(z), xs_(xs), ys_(ys), zs_(zs), n_(n) {}
00073         int    n()   const { return n_;  }
00074         double x()   const { return x_;  }
00075         double y()   const { return y_;  }
00076         double z()   const { return z_;  }
00077         double xs()  const { return xs_; }
00078         double ys()  const { return ys_; }
00079         double zs()  const { return zs_; }
00080         void   set(int n, double z, double zs) { n_= n; z_ = z; zs_ = zs; }
00081         void   set(int n, double x,double y,double z, double xs,double ys,double zs) 
00082                  { n_= n; x_ = x; xs_ = xs; y_ = y; ys_ = ys; z_ = z; zs_ = zs; }
00083       protected:    
00084         double x_,y_,z_,xs_,ys_,zs_;
00085         int n_;
00086     };
00087 
00088     QcdLowPtDQM(const edm::ParameterSet &parameters);
00089     virtual ~QcdLowPtDQM();
00090 
00091     void                          analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup);
00092     void                          beginJob(void);
00093     void                          beginLuminosityBlock(const edm::LuminosityBlock &l, 
00094                                                        const edm::EventSetup &iSetup);
00095     void                          beginRun(const edm::Run &r, const edm::EventSetup &iSetup);
00096     void                          endJob(void);
00097     void                          endRun(const edm::Run &r, const edm::EventSetup &iSetup);
00098     void                          endLuminosityBlock(const edm::LuminosityBlock &l, 
00099                                                      const edm::EventSetup &iSetup);
00100 
00101   private:
00102     void                          book1D(std::vector<MonitorElement*> &mes, 
00103                                          const std::string &name, const std::string &title, 
00104                                          int nx, double x1, double x2, bool sumw2=1, bool sbox=1);
00105     void                          book2D(std::vector<MonitorElement*> &mes, 
00106                                          const std::string &name, const std::string &title, 
00107                                          int nx, double x1, double x2, int ny, double y1, double y2,
00108                                          bool sumw2=1, bool sbox=1);
00109     void                          create1D(std::vector<TH1F*> &mes, 
00110                                            const std::string &name, const std::string &title, 
00111                                            int nx, double x1, double x2, bool sumw2=1, bool sbox=1);
00112     void                          create2D(std::vector<TH2F*> &mes, 
00113                                            const std::string &name, const std::string &title, 
00114                                            int nx, double x1, double x2, int ny, double y1, double y2,
00115                                            bool sumw2=1, bool sbox=1);
00116     void                          createHistos();
00117     void                          fill1D(std::vector<TH1F*> &hs, double val, double w=1.);
00118     void                          fill1D(std::vector<MonitorElement*> &mes, double val, double w=1.);
00119     void                          fill2D(std::vector<TH2F*> &hs, 
00120                                          double valx, double valy, double w=1.);
00121     void                          fill2D(std::vector<MonitorElement*> &mes, 
00122                                          double valx, double valy, double w=1.);
00123     void                          fill3D(std::vector<TH3F*> &hs, int gbin, double w=1.);
00124     void                          filldNdeta(const TH3F *AlphaTracklets,
00125                                              const std::vector<TH3F*> &NsigTracklets,
00126                                              const std::vector<TH3F*> &NbkgTracklets,
00127                                              const std::vector<TH1F*> &NEvsPerEta,
00128                                              std::vector<MonitorElement*> &hdNdEtaRawTrkl,
00129                                              std::vector<MonitorElement*> &hdNdEtaSubTrkl,
00130                                              std::vector<MonitorElement*> &hdNdEtaTrklets);
00131     void                          fillHltBits(const edm::Event &iEvent);
00132     void                          fillPixels(const edm::Event &iEvent);
00133     void                          fillPixelClusterInfos(const edm::Event &iEvent, int which=12);
00134     void                          fillPixelClusterInfos(const double vz,
00135                                                         const std::vector<Pixel> &pix, 
00136                                                         std::vector<MonitorElement*> &hClusterYSize,
00137                                                         std::vector<MonitorElement*> &hClusterADC);
00138     void                          fillTracklets(const edm::Event &iEvent, int which=12);
00139     void                          fillTracklets(std::vector<Tracklet> &tracklets, 
00140                                                 const std::vector<Pixel> &pix1, 
00141                                                 const std::vector<Pixel> &pix2,
00142                                                 const Vertex &trackletV);
00143     void                          fillTracklets(const std::vector<Tracklet> &tracklets, 
00144                                                 const std::vector<Pixel> &pixels,
00145                                                 const Vertex &trackletV,
00146                                                 const TH3F *AlphaTracklets,
00147                                                 std::vector<TH3F*> &NsigTracklets,
00148                                                 std::vector<TH3F*> &NbkgTracklets,
00149                                                 std::vector<TH1F*> &eventpereta,
00150                                                 std::vector<MonitorElement*> &detaphi,
00151                                                 std::vector<MonitorElement*> &deta,
00152                                                 std::vector<MonitorElement*> &dphi,
00153                                                 std::vector<MonitorElement*> &etavsvtx);
00154     template <typename TYPE>
00155     void                          getProduct(const std::string name, edm::Handle<TYPE> &prod,
00156                                              const edm::Event &event) const;    
00157     template <typename TYPE>
00158     bool                          getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
00159                                                  const edm::Event &event) const;
00160     void                          print(int level, const char *msg);
00161     void                          print(int level, const std::string &msg)
00162                                     { print(level,msg.c_str()); }
00163     void                          reallyPrint(int level, const char *msg);
00164     void                          trackletVertexUnbinned(const edm::Event &iEvent, int which=12);
00165     void                          trackletVertexUnbinned(std::vector<Pixel> &pix1, 
00166                                                          std::vector<Pixel> &pix2,
00167                                                          Vertex &vtx);
00168     double                        vertexZFromClusters(const std::vector<Pixel> &pix) const;
00169     void                          yieldAlphaHistogram(int which=12);
00170 
00171     std::string                   hltResName_;         //HLT trigger results name
00172     std::vector<std::string>      hltProcNames_;       //HLT process name(s)
00173     std::vector<std::string>      hltTrgNames_;        //HLT trigger name(s)
00174     std::string                   pixelName_;          //pixel reconstructed hits name
00175     std::string                   clusterVtxName_;     //cluster vertex name
00176     double                        ZVCut_;              //Z vertex cut for selected events
00177     double                        ZVEtaRegion_;        //Z vertex eta region
00178     double                        ZVVtxRegion_;        //Z vertex vtx region
00179     double                        dPhiVc_;             //dPhi vertex cut for tracklet based vertex
00180     double                        dZVc_;               //dZ vertex cut for tracklet based vertex
00181     double                        sigEtaCut_;          //signal tracklet eta cut
00182     double                        sigPhiCut_;          //signal tracklet phi cut
00183     double                        bkgEtaCut_;          //bkg tracklet eta cut
00184     double                        bkgPhiCut_;          //bgk tracklet phi cut
00185     int                           verbose_;            //verbosity (0=debug,1=warn,2=error,3=throw)
00186     int                           pixLayers_;          //12 for 12, 13 for 12 and 13, 23 for all 
00187     int                           clusLayers_;         //12 for 12, 13 for 12 and 13, 23 for all 
00188     bool                          useRecHitQ_;         //if true use rec hit quality word
00189     bool                          usePixelQ_;          //if true use pixel hit quality word
00190     std::vector<int>              hltTrgBits_;         //HLT trigger bit(s)
00191     std::vector<bool>             hltTrgDeci_;         //HLT trigger descision(s)
00192     std::vector<std::string>      hltTrgUsedNames_;    //HLT used trigger name(s)
00193     std::string                   hltUsedResName_;     //used HLT trigger results name
00194     std::vector<Pixel>            bpix1_;              //barrel pixels layer 1
00195     std::vector<Pixel>            bpix2_;              //barrel pixels layer 2
00196     std::vector<Pixel>            bpix3_;              //barrel pixels layer 3
00197     std::vector<Tracklet>         btracklets12_;       //barrel tracklets 12
00198     std::vector<Tracklet>         btracklets13_;       //barrel tracklets 13
00199     std::vector<Tracklet>         btracklets23_;       //barrel tracklets 23
00200     Vertex                        trackletV12_;        //reconstructed tracklet vertex 12
00201     Vertex                        trackletV13_;        //reconstructed tracklet vertex 13
00202     Vertex                        trackletV23_;        //reconstructed tracklet vertex 23
00203     std::vector<TH3F*>            NsigTracklets12_;    //number of signal tracklets 12
00204     std::vector<TH3F*>            NbkgTracklets12_;    //number of background tracklets 12
00205     std::vector<TH1F*>            hEvtCountsPerEta12_; //event count per tracklet12 eta-vtx region
00206     TH3F                         *AlphaTracklets12_;   //alpha correction for tracklets 12
00207     std::vector<TH3F*>            NsigTracklets13_;    //number of signal tracklets 13
00208     std::vector<TH3F*>            NbkgTracklets13_;    //number of background tracklets 13
00209     std::vector<TH1F*>            hEvtCountsPerEta13_; //event count per tracklet13 eta-vtx region 
00210     TH3F                         *AlphaTracklets13_;   //alpha correction for tracklets 13
00211     std::vector<TH3F*>            NsigTracklets23_;    //number of signal tracklets 23
00212     std::vector<TH3F*>            NbkgTracklets23_;    //number of background tracklets 23
00213     std::vector<TH1F*>            hEvtCountsPerEta23_; //event count per tracklet23 eta-vtx region
00214     TH3F                         *AlphaTracklets23_;   //alpha correction for tracklets 23
00215     HLTConfigProvider             hltConfig_;
00216     const TrackerGeometry        *tgeo_;               //tracker geometry
00217     DQMStore                     *theDbe_;             //dqm store
00218     MonitorElement               *repSumMap_;          //report summary map
00219     MonitorElement               *repSummary_;         //report summary
00220     MonitorElement               *h2TrigCorr_;         //trigger correlation plot 
00221     std::vector<MonitorElement*>  hNhitsL1_;           //number of hits on layer 1
00222     std::vector<MonitorElement*>  hNhitsL2_;           //number of hits on layer 2
00223     std::vector<MonitorElement*>  hNhitsL3_;           //number of hits on layer 3
00224     std::vector<MonitorElement*>  hNhitsL1z_;          //number of hits on layer 1 (zoomed)
00225     std::vector<MonitorElement*>  hNhitsL2z_;          //number of hits on layer 2 (zoomed)
00226     std::vector<MonitorElement*>  hNhitsL3z_;          //number of hits on layer 3 (zoomed)
00227     std::vector<MonitorElement*>  hdNdEtaHitsL1_;      //dN/dEta of hits on layer 1
00228     std::vector<MonitorElement*>  hdNdEtaHitsL2_;      //dN/dEta of hits on layer 2
00229     std::vector<MonitorElement*>  hdNdEtaHitsL3_;      //dN/dEta of hits on layer 3
00230     std::vector<MonitorElement*>  hdNdPhiHitsL1_;      //dN/dPhi of hits on layer 1
00231     std::vector<MonitorElement*>  hdNdPhiHitsL2_;      //dN/dPhi of hits on layer 2
00232     std::vector<MonitorElement*>  hdNdPhiHitsL3_;      //dN/dPhi of hits on layer 3
00233     std::vector<MonitorElement*>  hTrkVtxZ12_;         //tracklet z vertex 12 histograms
00234     std::vector<MonitorElement*>  hTrkVtxZ13_;         //tracklet z vertex 13 histograms
00235     std::vector<MonitorElement*>  hTrkVtxZ23_;         //tracklet z vertex 23 histograms
00236     std::vector<MonitorElement*>  hRawTrkEtaVtxZ12_;   //tracklet eta vs z vertex 12 histograms
00237     std::vector<MonitorElement*>  hRawTrkEtaVtxZ13_;   //tracklet eta vs z vertex 13 histograms
00238     std::vector<MonitorElement*>  hRawTrkEtaVtxZ23_;   //tracklet eta vs z vertex 23 histograms
00239     std::vector<MonitorElement*>  hTrkRawDetaDphi12_;  //tracklet12 Deta/Dphi distribution
00240     std::vector<MonitorElement*>  hTrkRawDeta12_;      //tracklet12 Deta distribution
00241     std::vector<MonitorElement*>  hTrkRawDphi12_;      //tracklet12 Dphi distribution
00242     std::vector<MonitorElement*>  hTrkRawDetaDphi13_;  //tracklet13 Deta/Dphi distribution
00243     std::vector<MonitorElement*>  hTrkRawDeta13_;      //tracklet13 Deta distribution
00244     std::vector<MonitorElement*>  hTrkRawDphi13_;      //tracklet13 Dphi distribution
00245     std::vector<MonitorElement*>  hTrkRawDetaDphi23_;  //tracklet23 Deta/Dphi distribution
00246     std::vector<MonitorElement*>  hTrkRawDeta23_;      //tracklet23 Deta distribution
00247     std::vector<MonitorElement*>  hTrkRawDphi23_;      //tracklet23 Dphi distribution
00248     std::vector<MonitorElement*>  hdNdEtaRawTrkl12_;   //dN/dEta from raw tracklets 12 
00249     std::vector<MonitorElement*>  hdNdEtaSubTrkl12_;   //dN/dEta from beta tracklets 12 
00250     std::vector<MonitorElement*>  hdNdEtaTrklets12_;   //dN/dEta corrected by alpha 12
00251     std::vector<MonitorElement*>  hdNdEtaRawTrkl13_;   //dN/dEta from raw tracklets 13
00252     std::vector<MonitorElement*>  hdNdEtaSubTrkl13_;   //dN/dEta from beta tracklets 13
00253     std::vector<MonitorElement*>  hdNdEtaTrklets13_;   //dN/dEta corrected by alpha 13
00254     std::vector<MonitorElement*>  hdNdEtaRawTrkl23_;   //dN/dEta from raw tracklets 23
00255     std::vector<MonitorElement*>  hdNdEtaSubTrkl23_;   //dN/dEta from beta tracklets 23
00256     std::vector<MonitorElement*>  hdNdEtaTrklets23_;   //dN/dEta corrected by alpha 23
00257     std::vector<MonitorElement*>  hClusterVertexZ_;    //cluster z vertex histograms
00258     std::vector<MonitorElement*>  hClusterYSize1_;     //cluster y size histograms on layer 1
00259     std::vector<MonitorElement*>  hClusterYSize2_;     //cluster y size histograms on layer 2
00260     std::vector<MonitorElement*>  hClusterYSize3_;     //cluster y size histograms on layer 3
00261     std::vector<MonitorElement*>  hClusterADC1_;       //cluster adc histograms on layer 1
00262     std::vector<MonitorElement*>  hClusterADC2_;       //cluster adc histograms on layer 2
00263     std::vector<MonitorElement*>  hClusterADC3_;       //cluster adc histograms on layer 3
00264 };
00265 
00266 //--------------------------------------------------------------------------------------------------
00267 template <typename TYPE>
00268 inline void QcdLowPtDQM::getProduct(const std::string name, edm::Handle<TYPE> &prod,
00269                                     const edm::Event &event) const
00270 {
00271   // Try to access data collection from EDM file. We check if we really get just one
00272   // product with the given name. If not we throw an exception.
00273 
00274   event.getByLabel(edm::InputTag(name),prod);
00275   if (!prod.isValid()) 
00276     throw edm::Exception(edm::errors::Configuration, "QcdLowPtDQM::GetProduct()\n")
00277       << "Collection with label " << name << " is not valid" <<  std::endl;
00278 }
00279 
00280 //--------------------------------------------------------------------------------------------------
00281 template <typename TYPE>
00282 inline bool QcdLowPtDQM::getProductSafe(const std::string name, edm::Handle<TYPE> &prod,
00283                                         const edm::Event &event) const
00284 {
00285   // Try to safely access data collection from EDM file. We check if we really get just one
00286   // product with the given name. If not, we return false.
00287 
00288   if (name.size()==0)
00289     return false;
00290 
00291   try {
00292     event.getByLabel(edm::InputTag(name),prod);
00293     if (!prod.isValid()) 
00294       return false;
00295   } catch (...) {
00296     return false;
00297   }
00298   return true;
00299 }
00300 
00301 //--------------------------------------------------------------------------------------------------
00302 inline void QcdLowPtDQM::print(int level, const char *msg)
00303 {
00304   // Print out message if it 
00305 
00306   if (level>=verbose_) 
00307     reallyPrint(level,msg);
00308 }
00309 #endif