CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_1/src/Validation/RPCRecHits/src/RPCRecHitValidClient.cc

Go to the documentation of this file.
00001 #include "Validation/RPCRecHits/interface/RPCRecHitValidClient.h"
00002 
00003 #include "FWCore/Framework/interface/MakerMacros.h"
00004 
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00007 
00008 using namespace std;
00009 
00010 typedef MonitorElement* MEP;
00011 
00012 RPCRecHitValidClient::RPCRecHitValidClient(const edm::ParameterSet& pset)
00013 {
00014   subDir_ = pset.getParameter<std::string>("subDir");
00015 }
00016 
00017 void RPCRecHitValidClient::endRun(const edm::Run& run, const edm::EventSetup& eventSetup)
00018 {
00019   DQMStore* dbe = edm::Service<DQMStore>().operator->();
00020   if ( !dbe ) return;
00021 
00022   dbe->setCurrentFolder(subDir_);
00023   MEP me_rollEfficiencyBarrel_eff = dbe->book1D("RollEfficiencyBarrel_eff", "Roll efficiency in Barrel;Efficiency [%]", 50+2, -2, 100+2);
00024   MEP me_rollEfficiencyEndcap_eff = dbe->book1D("RollEfficiencyEndcap_eff", "Roll efficiency in Endcap;Efficiency [%]", 50+2, -2, 100+2);
00025   MEP me_rollEfficiencyStatCutOffBarrel_eff = dbe->book1D("RollEfficiencyCutOffBarrel_eff", "Roll efficiency in Barrel without low stat chamber;Efficiency [%]", 50+2, -2, 100+2);
00026   MEP me_rollEfficiencyStatCutOffEndcap_eff = dbe->book1D("RollEfficiencyCutOffEndcap_eff", "Roll efficiency in Endcap without low stat chamber;Efficiency [%]", 50+2, -2, 100+2);
00027 
00028   const double maxNoise = 1e-7;
00029   MEP me_rollNoiseBarrel_noise = dbe->book1D("RollNoiseBarrel_noise", "Roll noise in Barrel;Noise level [Event^{-1}cm^{-2}]", 25+2, -maxNoise/25, maxNoise+maxNoise/25);
00030   MEP me_rollNoiseEndcap_noise = dbe->book1D("RollNoiseEndcap_noise", "Roll noise in Endcap;Noise level [Event^{-1}cm^{-2}]", 25+2, -maxNoise/25, maxNoise+maxNoise/25);
00031 
00032   MEP me_matchOccupancyBarrel_detId = dbe->get(subDir_+"/Occupancy/MatchOccupancyBarrel_detId");
00033   MEP me_matchOccupancyEndcap_detId = dbe->get(subDir_+"/Occupancy/MatchOccupancyEndcap_detId");
00034   MEP me_refOccupancyBarrel_detId = dbe->get(subDir_+"/Occupancy/RefOccupancyBarrel_detId");
00035   MEP me_refOccupancyEndcap_detId = dbe->get(subDir_+"/Occupancy/RefOccupancyEndcap_detId");
00036 
00037   if ( me_matchOccupancyBarrel_detId and me_refOccupancyBarrel_detId )
00038   {
00039     TH1* h_matchOccupancyBarrel_detId = me_matchOccupancyBarrel_detId->getTH1();
00040     TH1* h_refOccupancyBarrel_detId = me_refOccupancyBarrel_detId->getTH1();
00041 
00042     for ( int bin = 1, nBin = h_matchOccupancyBarrel_detId->GetNbinsX(); bin <= nBin; ++bin )
00043     {
00044       const double nRec = h_matchOccupancyBarrel_detId->GetBinContent(bin);
00045       const double nRef = h_refOccupancyBarrel_detId->GetBinContent(bin);
00046 
00047       const double eff = nRef ? nRec/nRef*100 : -1;
00048 
00049       me_rollEfficiencyBarrel_eff->Fill(eff);
00050       if ( nRef >= 20 ) me_rollEfficiencyStatCutOffBarrel_eff->Fill(eff);
00051     }
00052   }
00053 
00054   if ( me_matchOccupancyEndcap_detId and me_refOccupancyEndcap_detId )
00055   {
00056     TH1* h_matchOccupancyEndcap_detId = me_matchOccupancyEndcap_detId->getTH1();
00057     TH1* h_refOccupancyEndcap_detId = me_refOccupancyEndcap_detId->getTH1();
00058 
00059     for ( int bin = 1, nBin = h_matchOccupancyEndcap_detId->GetNbinsX(); bin <= nBin; ++bin )
00060     {
00061       const double nRec = h_matchOccupancyEndcap_detId->GetBinContent(bin);
00062       const double nRef = h_refOccupancyEndcap_detId->GetBinContent(bin);
00063 
00064       const double eff = nRef ? nRec/nRef*100 : -1;
00065 
00066       me_rollEfficiencyEndcap_eff->Fill(eff);
00067       if ( nRef >= 20 ) me_rollEfficiencyStatCutOffEndcap_eff->Fill(eff);
00068     }
00069   }
00070 
00071   MEP me_eventCount = dbe->get(subDir_+"/Occupancy/EventCount");
00072   const double nEvent = me_eventCount ? me_eventCount->getTH1()->GetBinContent(1) : 1;
00073   MEP me_noiseOccupancyBarrel_detId = dbe->get(subDir_+"/Occupancy/NoiseOccupancyBarrel_detId");
00074   MEP me_rollAreaBarrel_detId = dbe->get(subDir_+"/Occupancy/RollAreaBarrel_detId");
00075   if ( me_noiseOccupancyBarrel_detId and me_rollAreaBarrel_detId )
00076   {
00077     TH1* h_noiseOccupancyBarrel_detId = me_noiseOccupancyBarrel_detId->getTH1();
00078     TH1* h_rollAreaBarrel_detId = me_rollAreaBarrel_detId->getTH1();
00079 
00080     for ( int bin = 1, nBin = h_noiseOccupancyBarrel_detId->GetNbinsX(); bin <= nBin; ++bin )
00081     {
00082       const double noiseCount = h_noiseOccupancyBarrel_detId->GetBinContent(bin);
00083       const double area = h_rollAreaBarrel_detId->GetBinContent(bin);
00084       const double noiseLevel = area > 0 ? noiseCount/area/nEvent : 0;
00085       if ( noiseLevel == 0. ) me_rollNoiseBarrel_noise->Fill(-maxNoise/50); // Fill underflow bin if noise is exactly zero
00086       else me_rollNoiseBarrel_noise->Fill(std::min(noiseLevel, maxNoise));
00087     }
00088   }
00089 
00090   MEP me_noiseOccupancyEndcap_detId = dbe->get(subDir_+"/Occupancy/NoiseOccupancyEndcap_detId");
00091   MEP me_rollAreaEndcap_detId = dbe->get(subDir_+"/Occupancy/RollAreaEndcap_detId");
00092   if ( me_noiseOccupancyEndcap_detId and me_rollAreaEndcap_detId )
00093   {
00094     TH1* h_noiseOccupancyEndcap_detId = me_noiseOccupancyEndcap_detId->getTH1();
00095     TH1* h_rollAreaEndcap_detId = me_rollAreaEndcap_detId->getTH1();
00096 
00097     for ( int bin = 1, nBin = h_noiseOccupancyEndcap_detId->GetNbinsX(); bin <= nBin; ++bin )
00098     {
00099       const double noiseCount = h_noiseOccupancyEndcap_detId->GetBinContent(bin);
00100       const double area = h_rollAreaEndcap_detId->GetBinContent(bin);
00101       const double noiseLevel = area > 0 ? noiseCount/area/nEvent : 0;
00102       if ( noiseLevel == 0 ) me_rollNoiseEndcap_noise->Fill(-maxNoise/50); // Fill underflow bin if noise if exactly zero
00103       else me_rollNoiseEndcap_noise->Fill(std::min(noiseLevel, maxNoise));
00104     }
00105   }
00106 
00107 }
00108 
00109 DEFINE_FWK_MODULE(RPCRecHitValidClient);
00110