CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/DQM/TrackingMonitor/src/VertexMonitor.cc

Go to the documentation of this file.
00001 
00002 /*
00003  *  See header file for a description of this class.
00004  *
00005  *  $Date: 2012/03/29 17:21:28 $
00006  *  $Revision: 1.2 $
00007  *  \author:  Mia Tosi,40 3-B32,+41227671609 
00008  */
00009 
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 #include "FWCore/Utilities/interface/InputTag.h"
00012 #include "DQMServices/Core/interface/DQMStore.h"
00013 
00014 #include "DQM/TrackingMonitor/interface/VertexMonitor.h"
00015 
00016 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00017 
00018 #include "DQM/TrackingMonitor/interface/GetLumi.h"
00019 
00020 #include "TMath.h"
00021 
00022 
00023 VertexMonitor::VertexMonitor(const edm::ParameterSet& iConfig, edm::InputTag primaryVertexInputTag, edm::InputTag selectedPrimaryVertexInputTag, std::string pvLabel)
00024     : conf_( iConfig )
00025     , primaryVertexInputTag_         ( primaryVertexInputTag )
00026     , selectedPrimaryVertexInputTag_ ( selectedPrimaryVertexInputTag )
00027     , label_                         ( pvLabel )
00028     , NumberOfPVtx(NULL)
00029     , NumberOfPVtxVsBXlumi(NULL)
00030     , NumberOfPVtxVsGoodPVtx(NULL)
00031     , NumberOfGoodPVtx(NULL)
00032     , NumberOfGoodPVtxVsBXlumi(NULL)
00033     , FractionOfGoodPVtx(NULL)
00034     , FractionOfGoodPVtxVsBXlumi(NULL)
00035     , FractionOfGoodPVtxVsGoodPVtx(NULL)
00036     , FractionOfGoodPVtxVsPVtx(NULL)
00037     , NumberOfBADndofPVtx(NULL)
00038     , NumberOfBADndofPVtxVsBXlumi(NULL)
00039     , NumberOfBADndofPVtxVsGoodPVtx(NULL)
00040     , GoodPVtxSumPt(NULL)
00041     , GoodPVtxSumPtVsBXlumi(NULL)
00042     , GoodPVtxSumPtVsGoodPVtx(NULL)
00043     , GoodPVtxNumberOfTracks(NULL)
00044     , GoodPVtxNumberOfTracksVsBXlumi(NULL)
00045     , GoodPVtxNumberOfTracksVsGoodPVtx(NULL)
00046     , GoodPVtxNumberOfTracksVsGoodPVtxNdof(NULL)
00047     , GoodPVtxChi2oNDFVsGoodPVtx(NULL)
00048     , GoodPVtxChi2oNDFVsBXlumi(NULL)
00049     , GoodPVtxChi2ProbVsGoodPVtx(NULL)
00050     , GoodPVtxChi2ProbVsBXlumi(NULL)
00051     , doAllPlots_       ( conf_.getParameter<bool>("doAllPlots")        )
00052     , doPlotsVsBXlumi_  ( conf_.getParameter<bool>("doPlotsVsBXlumi")   )
00053     , doPlotsVsGoodPVtx_( conf_.getParameter<bool>("doPlotsVsGoodPVtx") )
00054 
00055 {
00056    //now do what ever initialization is needed
00057   if ( doPlotsVsBXlumi_ )
00058     lumiDetails_ = new GetLumi( iConfig.getParameter<edm::ParameterSet>("BXlumiSetup") );
00059 
00060 }
00061 
00062 
00063 VertexMonitor::~VertexMonitor()
00064 {
00065  
00066    // do anything here that needs to be done at desctruction time
00067    // (e.g. close files, deallocate resources etc.)
00068   //  if (lumiDetails_) delete lumiDetails_;
00069 }
00070 
00071 
00072 //
00073 // member functions
00074 //
00075 
00076 // -- Analyse
00077 // ------------ method called for each event  ------------
00078 // ------------------------------------------------------- //
00079 void
00080 VertexMonitor::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00081 {
00082 
00083   double bxlumi = 0.;
00084   if ( doPlotsVsBXlumi_ )
00085     bxlumi = lumiDetails_->getValue(iEvent);
00086   std::cout << "bxlumi : " << bxlumi << std::endl;
00087 
00088   size_t totalNumPV = 0;
00089   size_t totalNumBADndofPV = 0;
00090   edm::Handle< reco::VertexCollection > pvHandle;
00091   iEvent.getByLabel(primaryVertexInputTag_, pvHandle );
00092   if ( pvHandle.isValid() )
00093     {
00094       totalNumPV = pvHandle->size();
00095       std::cout << "totalNumPV : " << totalNumPV << std::endl;
00096       for (reco::VertexCollection::const_iterator pv = pvHandle->begin();
00097            pv != pvHandle->end(); ++pv) {
00098         //--- count pv w/ ndof < 4 
00099         if (pv->ndof() <  4.) totalNumBADndofPV++;
00100       }
00101     } else return;
00102   NumberOfPVtx        -> Fill( totalNumPV        );
00103   NumberOfBADndofPVtx -> Fill( totalNumBADndofPV );
00104   if ( doPlotsVsBXlumi_ ) {
00105     NumberOfPVtxVsBXlumi        -> Fill( bxlumi, totalNumPV        );
00106     NumberOfBADndofPVtxVsBXlumi -> Fill( bxlumi, totalNumBADndofPV );
00107   }
00108   
00109   size_t totalNumGoodPV = 0;
00110   edm::Handle< reco::VertexCollection > selpvHandle;
00111   iEvent.getByLabel(selectedPrimaryVertexInputTag_, selpvHandle );
00112   if ( selpvHandle.isValid() )
00113     totalNumGoodPV = selpvHandle->size();
00114   else return;
00115   std::cout << "totalNumGoodPV: " << totalNumGoodPV << std::endl;
00116   if ( doPlotsVsGoodPVtx_ ) {
00117     NumberOfPVtxVsGoodPVtx        -> Fill( totalNumGoodPV, totalNumPV        );
00118     NumberOfBADndofPVtxVsGoodPVtx -> Fill( totalNumGoodPV, totalNumBADndofPV );
00119   }
00120 
00121   double fracGoodPV = double(totalNumGoodPV)/double(totalNumPV);
00122   std::cout << "fracGoodPV: " << fracGoodPV << std::endl;
00123 
00124   NumberOfGoodPVtx    -> Fill( totalNumGoodPV    );
00125   FractionOfGoodPVtx  -> Fill( fracGoodPV        );
00126   if ( doPlotsVsBXlumi_ ) {
00127     NumberOfGoodPVtxVsBXlumi    -> Fill( bxlumi, totalNumGoodPV    );
00128     FractionOfGoodPVtxVsBXlumi  -> Fill( bxlumi, fracGoodPV        );
00129   }
00130   if ( doPlotsVsGoodPVtx_ ) {
00131     FractionOfGoodPVtxVsGoodPVtx  -> Fill( totalNumGoodPV, fracGoodPV        );
00132     FractionOfGoodPVtxVsPVtx      -> Fill( totalNumPV,     fracGoodPV        );
00133   }
00134 
00135   if ( selpvHandle->size() ) {
00136     double sumpt    = 0;
00137     size_t ntracks  = 0;
00138     double chi2ndf  = 0.; 
00139     double chi2prob = 0.;
00140 
00141     if (!selpvHandle->at(0).isFake()) {
00142       
00143       reco::Vertex pv = selpvHandle->at(0);
00144       
00145       ntracks  = pv.tracksSize();
00146       chi2ndf  = pv.normalizedChi2();
00147       chi2prob = TMath::Prob(pv.chi2(),(int)pv.ndof());
00148       
00149       for (reco::Vertex::trackRef_iterator itrk = pv.tracks_begin();
00150            itrk != pv.tracks_end(); ++itrk) {
00151         double pt = (**itrk).pt();
00152         sumpt += pt*pt;
00153       }
00154       GoodPVtxSumPt           -> Fill( sumpt   );
00155       GoodPVtxNumberOfTracks  -> Fill( ntracks );
00156 
00157       if ( doPlotsVsBXlumi_ ) {
00158         GoodPVtxSumPtVsBXlumi           -> Fill( bxlumi, sumpt    );
00159         GoodPVtxNumberOfTracksVsBXlumi  -> Fill( bxlumi, ntracks  );
00160         GoodPVtxChi2oNDFVsBXlumi        -> Fill( bxlumi, chi2ndf  );
00161         GoodPVtxChi2ProbVsBXlumi        -> Fill( bxlumi, chi2prob );
00162       }
00163       if ( doPlotsVsGoodPVtx_ ) {
00164         GoodPVtxSumPtVsGoodPVtx          -> Fill( totalNumGoodPV, sumpt    );
00165         GoodPVtxNumberOfTracksVsGoodPVtx -> Fill( totalNumGoodPV, ntracks  );
00166         GoodPVtxChi2oNDFVsGoodPVtx       -> Fill( totalNumGoodPV, chi2ndf  );
00167         GoodPVtxChi2ProbVsGoodPVtx       -> Fill( totalNumGoodPV, chi2prob );
00168       }
00169     }
00170   }
00171 }
00172 
00173 
00174 // ------------ method called once each job just before starting event loop  ------------
00175 void 
00176 VertexMonitor::beginJob(DQMStore * dqmStore_)
00177 {
00178     // parameters from the configuration
00179     std::string MEFolderName   = conf_.getParameter<std::string>("PVFolderName"); 
00180 
00181     // get binning from the configuration
00182     int    GoodPVtxBin   = conf_.getParameter<int>("GoodPVtxBin");
00183     double GoodPVtxMin   = conf_.getParameter<double>("GoodPVtxMin");
00184     double GoodPVtxMax   = conf_.getParameter<double>("GoodPVtxMax");
00185 
00186     // book histo
00187     // ----------------------//
00188     dqmStore_->setCurrentFolder(MEFolderName+"/"+label_);
00189 
00190     histname = "NumberOfPVtx_" + label_;
00191     NumberOfPVtx = dqmStore_->book1D(histname,histname, GoodPVtxBin,GoodPVtxMin,GoodPVtxMax);
00192     NumberOfPVtx->setAxisTitle("Number of PV",1);
00193     NumberOfPVtx->setAxisTitle("Number of Events",2);
00194     
00195     histname = "NumberOfGoodPVtx_" + label_;
00196     NumberOfGoodPVtx = dqmStore_->book1D(histname,histname, GoodPVtxBin,GoodPVtxMin,GoodPVtxMax);
00197     NumberOfGoodPVtx->setAxisTitle("Number of Good PV",1);
00198     NumberOfGoodPVtx->setAxisTitle("Number of Events",2);
00199     
00200     histname = "FractionOfGoodPVtx_" + label_;
00201     FractionOfGoodPVtx = dqmStore_->book1D(histname,histname, 100,0.,1.);
00202     FractionOfGoodPVtx->setAxisTitle("fraction of Good PV",1);
00203     FractionOfGoodPVtx->setAxisTitle("Number of Events",2);
00204     
00205     histname = "NumberOfBADndofPVtx_" + label_;
00206     NumberOfBADndofPVtx = dqmStore_->book1D(histname,histname,GoodPVtxBin,GoodPVtxMin,GoodPVtxMax);
00207     NumberOfBADndofPVtx->setAxisTitle("Number of BADndof #PV",1);
00208     NumberOfBADndofPVtx->setAxisTitle("Number of Events",2);
00209     
00210     histname = "GoodPVtxSumPt_" + label_;
00211     GoodPVtxSumPt = dqmStore_->book1D(histname,histname,100,0.,500.);
00212     GoodPVtxSumPt->setAxisTitle("primary vertex #Sum p_{T}^{2} [GeV^{2}/c^{2}]",1);
00213     GoodPVtxSumPt->setAxisTitle("Number of events",2);
00214 
00215     histname = "GoodPVtxNumberOfTracks_" + label_;
00216     GoodPVtxNumberOfTracks = dqmStore_->book1D(histname,histname,100,0.,100.);
00217     GoodPVtxNumberOfTracks->setAxisTitle("primary vertex number of tracks",1);
00218     GoodPVtxNumberOfTracks->setAxisTitle("Number of events",2);
00219 
00220     if ( doPlotsVsBXlumi_ ) {
00221       // get binning from the configuration
00222       edm::ParameterSet BXlumiParameters = conf_.getParameter<edm::ParameterSet>("BXlumiSetup");
00223       int    BXlumiBin   = BXlumiParameters.getParameter<int>("BXlumiBin");
00224       double BXlumiMin   = BXlumiParameters.getParameter<double>("BXlumiMin");
00225       double BXlumiMax   = BXlumiParameters.getParameter<double>("BXlumiMax");
00226 
00227       dqmStore_->setCurrentFolder(MEFolderName+"/"+label_+"/PUmonitoring/");
00228 
00229       histname = "NumberOfPVtxVsBXlumi_" + label_;
00230       NumberOfPVtxVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax,0.,60.,"");
00231       NumberOfPVtxVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00232       NumberOfPVtxVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00233       NumberOfPVtxVsBXlumi->setAxisTitle("Mean number of PV",2);
00234 
00235       histname = "NumberOfGoodPVtxVsBXlumi_" + label_;
00236       NumberOfGoodPVtxVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax,0.,60.,"");
00237       NumberOfGoodPVtxVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00238       NumberOfGoodPVtxVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00239       NumberOfGoodPVtxVsBXlumi->setAxisTitle("Mean number of PV",2);
00240 
00241       histname = "FractionOfGoodPVtxVsBXlumi_" + label_;
00242       FractionOfGoodPVtxVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax,0.,1.5,"");
00243       FractionOfGoodPVtxVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00244       FractionOfGoodPVtxVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00245       FractionOfGoodPVtxVsBXlumi->setAxisTitle("Mean number of PV",2);
00246 
00247       histname = "NumberOfBADndofPVtxVsBXlumi_" + label_;
00248       NumberOfBADndofPVtxVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax,0.,60.,"");
00249       NumberOfBADndofPVtxVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00250       NumberOfBADndofPVtxVsBXlumi->setAxisTitle("BADndof #PV",1);
00251       NumberOfBADndofPVtxVsBXlumi->setAxisTitle("Number of Events",2);
00252     
00253       histname = "GoodPVtxSumPtVsBXlumi_" + label_;
00254       GoodPVtxSumPtVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax,0.,500.,"");
00255       GoodPVtxSumPtVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00256       GoodPVtxSumPtVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00257       GoodPVtxSumPtVsBXlumi->setAxisTitle("Mean pv #Sum p_{T}^{2} [GeV^{2}/c]^{2}",2);
00258       
00259       histname = "GoodPVtxNumberOfTracksVsBXlumi_" + label_;
00260       GoodPVtxNumberOfTracksVsBXlumi = dqmStore_->bookProfile(histname,histname, BXlumiBin,BXlumiMin,BXlumiMax,0.,100.,"");
00261       GoodPVtxNumberOfTracksVsBXlumi->getTH1()->SetBit(TH1::kCanRebin);
00262       GoodPVtxNumberOfTracksVsBXlumi->setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00263       GoodPVtxNumberOfTracksVsBXlumi->setAxisTitle("Mean pv number of tracks",2);
00264 
00265       // get binning from the configuration
00266       double Chi2NDFMin = conf_.getParameter<double>("Chi2NDFMin");
00267       double Chi2NDFMax = conf_.getParameter<double>("Chi2NDFMax");
00268 
00269       double Chi2ProbMin  = conf_.getParameter<double>("Chi2ProbMin");
00270       double Chi2ProbMax  = conf_.getParameter<double>("Chi2ProbMax");
00271 
00272       histname = "Chi2oNDFVsBXlumi_" + label_;
00273       Chi2oNDFVsBXlumi  = dqmStore_->bookProfile(histname,histname,BXlumiBin, BXlumiMin,BXlumiMax,Chi2NDFMin,Chi2NDFMax,"");
00274       Chi2oNDFVsBXlumi -> getTH1()->SetBit(TH1::kCanRebin);
00275       Chi2oNDFVsBXlumi -> setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00276       Chi2oNDFVsBXlumi -> setAxisTitle("Mean #chi^{2}/ndof",2);
00277 
00278       histname = "Chi2ProbVsBXlumi_" + label_;
00279       Chi2ProbVsBXlumi  = dqmStore_->bookProfile(histname,histname,BXlumiBin, BXlumiMin,BXlumiMax,Chi2ProbMin,Chi2ProbMax,"");
00280       Chi2ProbVsBXlumi -> getTH1()->SetBit(TH1::kCanRebin);
00281       Chi2ProbVsBXlumi -> setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00282       Chi2ProbVsBXlumi -> setAxisTitle("Mean #chi^{2}/prob",2);
00283 
00284       histname = "GoodPVtxChi2oNDFVsBXlumi_" + label_;
00285       GoodPVtxChi2oNDFVsBXlumi  = dqmStore_->bookProfile(histname,histname,BXlumiBin, BXlumiMin,BXlumiMax,Chi2NDFMin,Chi2NDFMax,"");
00286       GoodPVtxChi2oNDFVsBXlumi -> getTH1()->SetBit(TH1::kCanRebin);
00287       GoodPVtxChi2oNDFVsBXlumi -> setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00288       GoodPVtxChi2oNDFVsBXlumi -> setAxisTitle("Mean PV #chi^{2}/ndof",2);
00289 
00290       histname = "GoodPVtxChi2ProbVsBXlumi_" + label_;
00291       GoodPVtxChi2ProbVsBXlumi  = dqmStore_->bookProfile(histname,histname,BXlumiBin, BXlumiMin,BXlumiMax,Chi2ProbMin,Chi2ProbMax,"");
00292       GoodPVtxChi2ProbVsBXlumi -> getTH1()->SetBit(TH1::kCanRebin);
00293       GoodPVtxChi2ProbVsBXlumi -> setAxisTitle("lumi BX [10^{30}Hzcm^{-2}]",1);
00294       GoodPVtxChi2ProbVsBXlumi -> setAxisTitle("Mean PV #chi^{2}/prob",2);
00295     }
00296 
00297     if ( doPlotsVsGoodPVtx_ ) {
00298 
00299       dqmStore_->setCurrentFolder(MEFolderName+"/"+label_+"/PUmonitoring/VsGoodPVtx");
00300 
00301       histname = "NumberOfPVtxVsGoodPVtx_" + label_;
00302       NumberOfPVtxVsGoodPVtx = dqmStore_->bookProfile(histname,histname, GoodPVtxBin,GoodPVtxMin,GoodPVtxMax,0.,60.,"");
00303       NumberOfPVtxVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00304       NumberOfPVtxVsGoodPVtx->setAxisTitle("Number of Good PV",1);
00305       NumberOfPVtxVsGoodPVtx->setAxisTitle("Mean number of PV",2);
00306 
00307       histname = "FractionOfGoodPVtxVsGoodPVtx_" + label_;
00308       FractionOfGoodPVtxVsGoodPVtx = dqmStore_->bookProfile(histname,histname, GoodPVtxBin,GoodPVtxMin,GoodPVtxMax,0.,60.,"");
00309       FractionOfGoodPVtxVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00310       FractionOfGoodPVtxVsGoodPVtx->setAxisTitle("Number of Good PV",1);
00311       FractionOfGoodPVtxVsGoodPVtx->setAxisTitle("Mean fraction of Good PV",2);
00312 
00313       histname = "FractionOfGoodPVtxVsPVtx_" + label_;
00314       FractionOfGoodPVtxVsPVtx = dqmStore_->bookProfile(histname,histname, GoodPVtxBin,GoodPVtxMin,GoodPVtxMax,0.,60.,"");
00315       FractionOfGoodPVtxVsPVtx->getTH1()->SetBit(TH1::kCanRebin);
00316       FractionOfGoodPVtxVsPVtx->setAxisTitle("Number of Good PV",1);
00317       FractionOfGoodPVtxVsPVtx->setAxisTitle("Mean number of Good PV",2);
00318 
00319       histname = "NumberOfBADndofPVtxVsGoodPVtx_" + label_;
00320       NumberOfBADndofPVtxVsGoodPVtx = dqmStore_->bookProfile(histname,histname, GoodPVtxBin,GoodPVtxMin,GoodPVtxMax,0.,60.,"");
00321       NumberOfBADndofPVtxVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00322       NumberOfBADndofPVtxVsGoodPVtx->setAxisTitle("Number of Good PV",1);
00323       NumberOfBADndofPVtxVsGoodPVtx->setAxisTitle("Mean Number of BAD PV",2);
00324     
00325       histname = "GoodPVtxSumPtVsGoodPVtx_" + label_;
00326       GoodPVtxSumPtVsGoodPVtx = dqmStore_->bookProfile(histname,histname, GoodPVtxBin,GoodPVtxMin,GoodPVtxMax,0.,500.,"");
00327       GoodPVtxSumPtVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00328       GoodPVtxSumPtVsGoodPVtx->setAxisTitle("Number of Good PV",1);
00329       GoodPVtxSumPtVsGoodPVtx->setAxisTitle("Mean pv #Sum p_{T}^{2} [GeV^{2}/c]^{2}",2);
00330       
00331       histname = "GoodPVtxNumberOfTracksVsGoodPVtx_" + label_;
00332       GoodPVtxNumberOfTracksVsGoodPVtx = dqmStore_->bookProfile(histname,histname, GoodPVtxBin,GoodPVtxMin,GoodPVtxMax,0.,100.,"");
00333       GoodPVtxNumberOfTracksVsGoodPVtx->getTH1()->SetBit(TH1::kCanRebin);
00334       GoodPVtxNumberOfTracksVsGoodPVtx->setAxisTitle("Number of Good PV",1);
00335       GoodPVtxNumberOfTracksVsGoodPVtx->setAxisTitle("Mean pv number of tracks",2);
00336       
00337       // get binning from the configuration
00338       double Chi2NDFMin = conf_.getParameter<double>("Chi2NDFMin");
00339       double Chi2NDFMax = conf_.getParameter<double>("Chi2NDFMax");
00340 
00341       double Chi2ProbMin  = conf_.getParameter<double>("Chi2ProbMin");
00342       double Chi2ProbMax  = conf_.getParameter<double>("Chi2ProbMax");
00343 
00344       histname = "GoodPVtxChi2oNDFVsGoodPVtx_" + label_;
00345       GoodPVtxChi2oNDFVsGoodPVtx  = dqmStore_->bookProfile(histname,histname,GoodPVtxBin, GoodPVtxMin,GoodPVtxMax,Chi2NDFMin,Chi2NDFMax,"");
00346       GoodPVtxChi2oNDFVsGoodPVtx -> getTH1()->SetBit(TH1::kCanRebin);
00347       GoodPVtxChi2oNDFVsGoodPVtx -> setAxisTitle("Number of Good PV",1);
00348       GoodPVtxChi2oNDFVsGoodPVtx -> setAxisTitle("Mean PV #chi^{2}/ndof",2);
00349 
00350       histname = "GoodPVtxChi2ProbVsGoodPVtx_" + label_;
00351       GoodPVtxChi2ProbVsGoodPVtx  = dqmStore_->bookProfile(histname,histname,GoodPVtxBin, GoodPVtxMin,GoodPVtxMax,Chi2ProbMin,Chi2ProbMax,"");
00352       GoodPVtxChi2ProbVsGoodPVtx -> getTH1()->SetBit(TH1::kCanRebin);
00353       GoodPVtxChi2ProbVsGoodPVtx -> setAxisTitle("Number of Good PV",1);
00354       GoodPVtxChi2ProbVsGoodPVtx -> setAxisTitle("Mean PV #chi^{2}/prob",2);
00355 
00356     }
00357 }
00358  
00359 
00360 void 
00361 VertexMonitor::beginLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
00362 {
00363 }
00364 
00365 // ------------ method called when ending the processing of a luminosity block  ------------
00366 void 
00367 VertexMonitor::endLuminosityBlock(edm::LuminosityBlock const&, edm::EventSetup const&)
00368 {
00369 }
00370 
00371 // ------------ method fills 'descriptions' with the allowed parameters for the module  ------------
00372 void
00373 VertexMonitor::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
00374   //The following says we do not know what parameters are allowed so do no validation
00375   // Please change this to state exactly what you do use, even if it is no parameters
00376   edm::ParameterSetDescription desc;
00377   desc.setUnknown();
00378   descriptions.addDefault(desc);
00379 }