CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/DQM/DTMonitorClient/src/DTDeadChannelTest.cc

Go to the documentation of this file.
00001 
00002 /*
00003  *  See header file for a description of this class.
00004  *
00005  *  $Date: 2010/01/05 10:15:46 $
00006  *  $Revision: 1.15 $
00007  *  \author G. Mila - INFN Torino
00008  */
00009 
00010 
00011 #include <DQM/DTMonitorClient/src/DTDeadChannelTest.h>
00012 
00013 // Framework
00014 #include <FWCore/Framework/interface/EventSetup.h>
00015 
00016 
00017 // Geometry
00018 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00019 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00020 #include "Geometry/DTGeometry/interface/DTLayer.h"
00021 #include "Geometry/DTGeometry/interface/DTTopology.h"
00022 
00023 #include <CondFormats/DTObjects/interface/DTTtrig.h>
00024 #include <CondFormats/DataRecord/interface/DTTtrigRcd.h>
00025 
00026 #include "DQMServices/Core/interface/DQMStore.h"
00027 #include "DQMServices/Core/interface/MonitorElement.h"
00028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00029 
00030 #include <stdio.h>
00031 #include <sstream>
00032 #include <math.h>
00033 
00034 
00035 using namespace edm;
00036 using namespace std;
00037 
00038 
00039 DTDeadChannelTest::DTDeadChannelTest(const edm::ParameterSet& ps){
00040  
00041   edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: Constructor";
00042 
00043   parameters = ps;
00044 
00045   dbe = edm::Service<DQMStore>().operator->();
00046 
00047   prescaleFactor = parameters.getUntrackedParameter<int>("diagnosticPrescale", 1);
00048 
00049 }
00050 
00051 DTDeadChannelTest::~DTDeadChannelTest(){
00052 
00053   edm::LogVerbatim ("deadChannel") << "DTDeadChannelTest: analyzed " << nevents << " events";
00054 
00055 }
00056 
00057 
00058 void DTDeadChannelTest::beginJob(){
00059 
00060   edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: BeginJob";
00061 
00062   nevents = 0;
00063 
00064 }
00065 
00066 void DTDeadChannelTest::beginRun(Run const& run, EventSetup const& context) {
00067 
00068   edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: BeginRun";
00069 
00070   // Get the geometry
00071   context.get<MuonGeometryRecord>().get(muonGeom);
00072 
00073 }
00074 
00075 
00076 void DTDeadChannelTest::beginLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& context) {
00077 
00078   edm::LogVerbatim ("deadChannel") <<"[DTDeadChannelTest]: Begin of LS transition";
00079 
00080   // Get the run number
00081   run = lumiSeg.run();
00082 
00083 }
00084 
00085 
00086 
00087 void DTDeadChannelTest::analyze(const edm::Event& e, const edm::EventSetup& context){
00088 
00089   nevents++;
00090   edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: "<<nevents<<" events";
00091 
00092 }
00093 
00094 
00095 
00096 void DTDeadChannelTest::endLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& context) {
00097   
00098   // counts number of updats (online mode) or number of events (standalone mode)
00099   //nevents++;
00100   // if running in standalone perform diagnostic only after a reasonalbe amount of events
00101   //if ( parameters.getUntrackedParameter<bool>("runningStandalone", false) && 
00102   //     nevents%parameters.getUntrackedParameter<int>("diagnosticPrescale", 1000) != 0 ) return;
00103   //edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: "<<nevents<<" updates";
00104 
00105 
00106   edm::LogVerbatim ("deadChannel") <<"[DTDeadChannelTest]: End of LS transition, performing the DQM client operation";
00107 
00108   // counts number of lumiSegs 
00109   nLumiSegs = lumiSeg.id().luminosityBlock();
00110 
00111   // prescale factor
00112   if ( nLumiSegs%prescaleFactor != 0 ) return;
00113 
00114   edm::LogVerbatim ("deadChannel") <<"[DTDeadChannelTest]: "<<nLumiSegs<<" updates";
00115 
00116 
00117   vector<DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
00118   vector<DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
00119 
00120   edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest]: Occupancy tests results";
00121 
00122   // Loop over the chambers
00123   for (; ch_it != ch_end; ++ch_it) {
00124     DTChamberId chID = (*ch_it)->id();
00125     vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin(); 
00126     vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
00127 
00128     stringstream wheel; wheel << chID.wheel();
00129     stringstream station; station << chID.station();
00130     stringstream sector; sector << chID.sector();
00131     
00132     context.get<DTTtrigRcd>().get(tTrigMap);
00133 
00134     string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str(); 
00135 
00136     // Get the ME produced by DigiTask Source
00137     MonitorElement * noise_histo = dbe->get(getMEName("OccupancyNoise_perCh", chID));   
00138     MonitorElement * hitInTime_histo = dbe->get(getMEName("OccupancyInTimeHits_perCh", chID));
00139 
00140     // ME -> TH2F
00141     if(noise_histo && hitInTime_histo) {          
00142       TH2F * noise_histo_root = noise_histo->getTH2F();
00143       TH2F * hitInTime_histo_root = hitInTime_histo->getTH2F();
00144 
00145       // Loop over the SuperLayers
00146       for(; sl_it != sl_end; ++sl_it) {
00147         DTSuperLayerId slID = (*sl_it)->id();
00148         vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin();
00149         vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
00150             
00151         // ttrig and rms are counts
00152         float tTrig, tTrigRMS, kFactor;
00153         tTrigMap->get(slID, tTrig, tTrigRMS, kFactor, DTTimeUnits::counts);
00154       
00155         // Loop over the layers
00156         for(; l_it != l_end; ++l_it) {
00157           DTLayerId lID = (*l_it)->id();
00158 
00159           //Parameters to fill histos
00160           stringstream superLayer; superLayer << slID.superlayer();
00161           stringstream layer; layer << lID.layer();
00162           string HistoNameTest = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() +  "_SL" + superLayer.str() +  "_L" + layer.str();
00163 
00164           const int firstWire = muonGeom->layer(lID)->specificTopology().firstChannel();
00165           const int lastWire = muonGeom->layer(lID)->specificTopology().lastChannel();
00166 
00167           int entry=-1;
00168           if(slID.superlayer() == 1) entry=0;
00169           if(slID.superlayer() == 2) entry=4;
00170           if(slID.superlayer() == 3) entry=8;
00171           int YBinNumber = entry+lID.layer();
00172               
00173 
00174           // Loop over the TH2F bin and fill the ME to be used for the Quality Test
00175           for(int bin=firstWire; bin <= lastWire; bin++) {
00176             if (OccupancyDiffHistos.find(HistoNameTest) == OccupancyDiffHistos.end()) bookHistos(lID, firstWire, lastWire);
00177             // tMax default value
00178             float tMax = 450.0;
00179 
00180             float difference = (hitInTime_histo_root->GetBinContent(bin, YBinNumber) / tMax) 
00181                                - (noise_histo_root->GetBinContent(bin, YBinNumber) / tTrig);
00182             OccupancyDiffHistos.find(HistoNameTest)->second->setBinContent(bin, difference);
00183           }
00184         } // loop on layers
00185       } // loop on superlayers
00186     }
00187   } // loop on chambers
00188 
00189   // Occupancy Difference test 
00190   string OccupancyDiffCriterionName = parameters.getUntrackedParameter<string>("OccupancyDiffTestName","OccupancyDiffInRange"); 
00191   for(map<string, MonitorElement*>::const_iterator hOccDiff = OccupancyDiffHistos.begin();
00192       hOccDiff != OccupancyDiffHistos.end();
00193       hOccDiff++) {
00194     const QReport * theOccupancyDiffQReport = (*hOccDiff).second->getQReport(OccupancyDiffCriterionName);
00195     if(theOccupancyDiffQReport) {
00196       vector<dqm::me_util::Channel> badChannels = theOccupancyDiffQReport->getBadChannels();
00197       for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); 
00198            channel != badChannels.end(); channel++) {
00199         edm::LogError ("deadChannel") << "Layer : "<<(*hOccDiff).first<<" Bad occupancy difference channels: "<<(*channel).getBin()<<" Contents : "<<(*channel).getContents();
00200       }
00201       // FIXME: getMessage() sometimes returns and invalid string (null pointer inside QReport data member)
00202       // edm::LogWarning("deadChannel")<< "-------- Layer : "<<(*hOccDiff).first<<"  "<<theOccupancyDiffQReport->getMessage()<<" ------- "<<theOccupancyDiffQReport->getStatus(); 
00203     }
00204   }
00205 
00206 }
00207 
00208 
00209 void DTDeadChannelTest::endJob(){
00210 
00211   edm::LogVerbatim ("deadChannel") << "[DTDeadChannelTest] endjob called!";
00212 
00213   dbe->rmdir("DT/Tests/DTDeadChannel");
00214 
00215 }
00216 
00217 string DTDeadChannelTest::getMEName(string histoTag, const DTChamberId & chId) {
00218 
00219   stringstream wheel; wheel << chId.wheel();
00220   stringstream station; station << chId.station();
00221   stringstream sector; sector << chId.sector();
00222 
00223   string folderRoot = parameters.getUntrackedParameter<string>("folderRoot", "Collector/FU0/");
00224   string folderName = 
00225     folderRoot + "DT/DTDigiTask/Wheel" +  wheel.str() +
00226     "/Station" + station.str() +
00227     "/Sector" + sector.str() + 
00228     "/Occupancies" + "/";
00229   
00230   string histoname = folderName + histoTag  
00231     + "_W" + wheel.str() 
00232     + "_St" + station.str() 
00233     + "_Sec" + sector.str();
00234   
00235   return histoname;
00236   
00237 }
00238 
00239 
00240 void DTDeadChannelTest::bookHistos(const DTLayerId & lId, int firstWire, int lastWire) {
00241 
00242   stringstream wheel; wheel << lId.superlayerId().wheel();
00243   stringstream station; station << lId.superlayerId().station();        
00244   stringstream sector; sector << lId.superlayerId().sector();
00245   stringstream superLayer; superLayer << lId.superlayerId().superlayer();
00246   stringstream layer; layer << lId.layer();
00247 
00248   string HistoName = "W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() +  "_SL" + superLayer.str() +  "_L" + layer.str();
00249   string OccupancyDiffHistoName =  "OccupancyDiff_" + HistoName; 
00250 
00251   dbe->setCurrentFolder("DT/Tests/DTDeadChannel/Wheel" + wheel.str() +
00252                            "/Station" + station.str() +
00253                            "/Sector" + sector.str());
00254 
00255   OccupancyDiffHistos[HistoName] = dbe->book1D(OccupancyDiffHistoName.c_str(),OccupancyDiffHistoName.c_str(),lastWire-firstWire+1, firstWire-0.5, lastWire+0.5);
00256 
00257 }