CMS 3D CMS Logo

DTChamberEfficiencyTest.cc

Go to the documentation of this file.
00001 
00002 
00003 /*
00004  *  See header file for a description of this class.
00005  *
00006  *  $Date: 2008/10/22 09:40:34 $
00007  *  $Revision: 1.14 $
00008  *  \author G. Mila - INFN Torino
00009  */
00010 
00011 
00012 #include <DQM/DTMonitorClient/src/DTChamberEfficiencyTest.h>
00013 #include "DQMServices/Core/interface/MonitorElement.h"
00014 #include "DQMServices/Core/interface/DQMStore.h"
00015 
00016 // Framework
00017 #include <FWCore/Framework/interface/EventSetup.h>
00018 
00019 
00020 // Geometry
00021 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00022 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00023 
00024 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00025 
00026 #include <stdio.h>
00027 #include <sstream>
00028 #include <math.h>
00029 
00030 
00031 using namespace edm;
00032 using namespace std;
00033 
00034 
00035 
00036 DTChamberEfficiencyTest::DTChamberEfficiencyTest(const edm::ParameterSet& ps){
00037 
00038   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "[DTChamberEfficiencyTest]: Constructor";
00039 
00040   parameters = ps;
00041 
00042   dbe = edm::Service<DQMStore>().operator->();
00043 
00044   prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1);
00045 
00046 }
00047 
00048 
00049 
00050 DTChamberEfficiencyTest::~DTChamberEfficiencyTest(){
00051 
00052   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "DTChamberEfficiencyTest: analyzed " << nevents << " events";
00053 
00054 }
00055 
00056 
00057 void DTChamberEfficiencyTest::beginJob(const edm::EventSetup& context){
00058 
00059   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") <<"[DTChamberEfficiencyTest]: BeginJob"; 
00060 
00061   nevents = 0;
00062 
00063   // Get the geometry
00064   context.get<MuonGeometryRecord>().get(muonGeom);
00065 
00066 }
00067 
00068 
00069 void DTChamberEfficiencyTest::beginLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& context) {
00070 
00071   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") <<"[DTChamberEfficiencyTest]: Begin of LS transition";
00072 
00073   // Get the run number
00074   run = lumiSeg.run();
00075 
00076 }
00077 
00078 
00079 void DTChamberEfficiencyTest::beginRun(const edm::Run& run, const edm::EventSetup& setup){
00080   
00081   // Get the DT Geometry
00082   setup.get<MuonGeometryRecord>().get(muonGeom);
00083 
00084   // Loop over all the chambers
00085   vector<DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
00086   vector<DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
00087   for (; ch_it != ch_end; ++ch_it) {
00088     // histo booking
00089     bookHistos((*ch_it)->id());
00090   }
00091 
00092   //summary histo booking
00093   bookHistos();
00094 }
00095 
00096 void DTChamberEfficiencyTest::analyze(const edm::Event& e, const edm::EventSetup& context){
00097 
00098   nevents++;
00099   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "[DTChamberEfficiencyTest]: "<<nevents<<" events";
00100 
00101 }
00102 
00103 
00104 void DTChamberEfficiencyTest::endLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& context) {
00105   
00106   // counts number of updats (online mode) or number of events (standalone mode)
00107   //nevents++;
00108   // if running in standalone perform diagnostic only after a reasonalbe amount of events
00109   //if ( parameters.getUntrackedParameter<bool>("runningStandalone", false) && 
00110   //     nevents%parameters.getUntrackedParameter<int>("diagnosticPrescale", 1000) != 0 ) return;  
00111   //edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "[DTChamberEfficiencyTest]: "<<nevents<<" updates";
00112   
00113   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") <<"[DTChamberEfficiencyTest]: End of LS transition, performing the DQM client operation";
00114 
00115   // counts number of lumiSegs 
00116   nLumiSegs = lumiSeg.id().luminosityBlock();
00117 
00118   // prescale factor
00119   if ( nLumiSegs%prescaleFactor != 0 ) return;
00120 
00121   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") <<"[DTChamberEfficiencyTest]: "<<nLumiSegs<<" updates";
00122   
00123 
00124   vector<DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
00125   vector<DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
00126 
00127   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "[DTChamberEfficiencyTest]: ChamberEfficiency tests results"; 
00128   
00129   // Loop over the chambers
00130   for (; ch_it != ch_end; ++ch_it) {
00131     DTChamberId chID = (*ch_it)->id();
00132     
00133     stringstream wheel; wheel << chID.wheel();
00134     stringstream station; station << chID.station();
00135     stringstream sector; sector << chID.sector();
00136     
00137     string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
00138     
00139     // Get the ME produced by EfficiencyTask Source
00140     MonitorElement * GoodSegDen_histo = dbe->get(getMEName("hEffGoodSegVsPosDen", chID));       
00141     MonitorElement * GoodCloseSegNum_histo = dbe->get(getMEName("hEffGoodCloseSegVsPosNum", chID));
00142     
00143     // ME -> TH1F
00144     if(GoodSegDen_histo && GoodCloseSegNum_histo) {       
00145       TH2F * GoodSegDen_histo_root = GoodSegDen_histo->getTH2F();
00146       TH2F * GoodCloseSegNum_histo_root = GoodCloseSegNum_histo->getTH2F();
00147         
00148       int lastBinX=(*GoodSegDen_histo_root).GetNbinsX();
00149       TH1D* proxN=GoodCloseSegNum_histo_root->ProjectionX();
00150       TH1D* proxD=GoodSegDen_histo_root->ProjectionX();
00151 
00152       int lastBinY=(*GoodSegDen_histo_root).GetNbinsY();
00153       TH1D* proyN=GoodCloseSegNum_histo_root->ProjectionY();
00154       TH1D* proyD=GoodSegDen_histo_root->ProjectionY();
00155           
00156       for(int xBin=1; xBin<=lastBinX; xBin++) {
00157         if(proxD->GetBinContent(xBin)!=0){
00158           float Xefficiency = proxN->GetBinContent(xBin) / proxD->GetBinContent(xBin);
00159           xEfficiencyHistos.find(HistoName)->second->setBinContent(xBin, Xefficiency);
00160         }
00161 
00162         for(int yBin=1; yBin<=lastBinY; yBin++) {
00163           if(GoodSegDen_histo_root->GetBinContent(xBin, yBin)!=0){
00164             float XvsYefficiency = GoodCloseSegNum_histo_root->GetBinContent(xBin, yBin) / GoodSegDen_histo_root->GetBinContent(xBin, yBin);
00165             xVSyEffHistos.find(HistoName)->second->setBinContent(xBin, yBin, XvsYefficiency);
00166           }
00167         }
00168             
00169       }
00170           
00171       for(int yBin=1; yBin<=lastBinY; yBin++) {
00172         if(proyD->GetBinContent(yBin)!=0){
00173           float Yefficiency = proyN->GetBinContent(yBin) / proyD->GetBinContent(yBin);
00174           yEfficiencyHistos.find(HistoName)->second->setBinContent(yBin, Yefficiency);
00175         }
00176       }
00177     }
00178   } // loop on chambers
00179   
00180   
00181   // ChamberEfficiency test on X axis
00182   string XEfficiencyCriterionName = parameters.getUntrackedParameter<string>("XEfficiencyTestName","ChEfficiencyInRangeX"); 
00183   for(map<string, MonitorElement*>::const_iterator hXEff = xEfficiencyHistos.begin();
00184       hXEff != xEfficiencyHistos.end();
00185       hXEff++) {
00186     const QReport * theXEfficiencyQReport = (*hXEff).second->getQReport(XEfficiencyCriterionName);
00187     if(theXEfficiencyQReport) {
00188       vector<dqm::me_util::Channel> badChannels = theXEfficiencyQReport->getBadChannels();
00189       for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); 
00190            channel != badChannels.end(); channel++) {
00191         edm::LogError ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "Chamber : " << (*hXEff).first << " Bad XChamberEfficiency channels: "<<(*channel).getBin()<<"  Contents : "<<(*channel).getContents();
00192       }
00193       // FIXME: getMessage() sometimes returns and invalid string (null pointer inside QReport data member)
00194       // edm::LogWarning ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "-------- Chamber : "<<(*hXEff).first<<"  "<<theXEfficiencyQReport->getMessage()<<" ------- "<<theXEfficiencyQReport->getStatus();
00195     }
00196   }
00197 
00198 
00199   // ChamberEfficiency test on Y axis
00200   string YEfficiencyCriterionName = parameters.getUntrackedParameter<string>("YEfficiencyTestName","ChEfficiencyInRangeY"); 
00201   for(map<string, MonitorElement*>::const_iterator hYEff = yEfficiencyHistos.begin();
00202       hYEff != yEfficiencyHistos.end();
00203       hYEff++) {
00204     const QReport * theYEfficiencyQReport = (*hYEff).second->getQReport(YEfficiencyCriterionName);
00205     if(theYEfficiencyQReport) {
00206       vector<dqm::me_util::Channel> badChannels = theYEfficiencyQReport->getBadChannels();
00207       for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); 
00208            channel != badChannels.end(); channel++) {
00209         edm::LogError ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "Chamber : " << (*hYEff).first <<" Bad YChamberEfficiency channels: "<<(*channel).getBin()<<"  Contents : "<<(*channel).getContents();
00210       }
00211       // FIXME: getMessage() sometimes returns and invalid string (null pointer inside QReport data member)
00212       // edm::LogWarning ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "-------- Chamber : "<<(*hYEff).first<<"  "<<theYEfficiencyQReport->getMessage()<<" ------- "<<theYEfficiencyQReport->getStatus();
00213     }
00214   }
00215 
00216   //Fill the report summary histos
00217   for(int wh=-2; wh<=2; wh++){
00218     for(int sec=1; sec<=12; sec++){
00219       for(int st=1; st<=4; st++){
00220 
00221         summaryHistos[wh]->Fill(sec,st,1);
00222       
00223       }
00224     }
00225   }
00226 
00227 }
00228 
00229 
00230 
00231 void DTChamberEfficiencyTest::endJob(){
00232 
00233   edm::LogVerbatim ("DTDQM|DTMonitorClient|DTChamberEfficiencyTest") << "[DTChamberEfficiencyTest] endjob called!";
00234 
00235 }
00236 
00237 
00238 
00239 
00240 string DTChamberEfficiencyTest::getMEName(string histoTag, const DTChamberId & chID) {
00241 
00242   stringstream wheel; wheel << chID.wheel();
00243   stringstream station; station << chID.station();
00244   stringstream sector; sector << chID.sector();
00245  
00246   string folderRoot = parameters.getUntrackedParameter<string>("folderRoot", "Collector/FU0/");
00247   string folderName = 
00248     folderRoot + "DT/01-DTChamberEfficiency/Task/Wheel" +  wheel.str() +
00249     "/Sector" + sector.str() +
00250     "/Station" + station.str() + "/";
00251 
00252   string histoname = folderName + histoTag  
00253     + "_W" + wheel.str() 
00254     + "_St" + station.str() 
00255     + "_Sec" + sector.str();
00256   
00257   return histoname;
00258   
00259 }
00260 
00261 
00262 void DTChamberEfficiencyTest::bookHistos(const DTChamberId & chId) {
00263 
00264   stringstream wheel; wheel << chId.wheel();
00265   stringstream station; station << chId.station();      
00266   stringstream sector; sector << chId.sector();
00267 
00268   string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
00269   string xEfficiencyHistoName =  "xEfficiency_" + HistoName; 
00270   string yEfficiencyHistoName =  "yEfficiency_" + HistoName; 
00271   string xVSyEffHistoName =  "xVSyEff_" + HistoName; 
00272 
00273   dbe->setCurrentFolder("DT/01-DTChamberEfficiency/Wheel" + wheel.str() +
00274                         "/Sector" + sector.str() +
00275                         "/Station" + station.str());
00276 
00277   xEfficiencyHistos[HistoName] = dbe->book1D(xEfficiencyHistoName.c_str(),xEfficiencyHistoName.c_str(),25,-250.,250.);
00278   yEfficiencyHistos[HistoName] = dbe->book1D(yEfficiencyHistoName.c_str(),yEfficiencyHistoName.c_str(),25,-250.,250.);
00279   xVSyEffHistos[HistoName] = dbe->book2D(xVSyEffHistoName.c_str(),xVSyEffHistoName.c_str(),25,-250.,250., 25,-250.,250.);
00280 
00281 }
00282 
00283 
00284 void DTChamberEfficiencyTest::bookHistos() {
00285 
00286   for(int wh=-2; wh<=2; wh++){
00287     stringstream wheel; wheel << wh;
00288     string histoName =  "chEfficiencySummary_W" + wheel.str();
00289     dbe->setCurrentFolder("DT/01-DTChamberEfficiency");
00290     summaryHistos[wh] = dbe->book2D(histoName.c_str(),histoName.c_str(),12,1,13,4,1,5);
00291     summaryHistos[wh]->setAxisTitle("Sector",1);
00292     summaryHistos[wh]->setBinLabel(1,"MB1",2);
00293     summaryHistos[wh]->setBinLabel(2,"MB2",2);
00294     summaryHistos[wh]->setBinLabel(3,"MB3",2);
00295     summaryHistos[wh]->setBinLabel(4,"MB4",2);
00296   }
00297 
00298 }

Generated on Tue Jun 9 17:32:34 2009 for CMSSW by  doxygen 1.5.4