CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/DQM/DTMonitorModule/src/DTSegmentAnalysisTask.cc

Go to the documentation of this file.
00001 
00002 /*
00003  *  See header file for a description of this class.
00004  *
00005  *  $Date: 2011/11/23 14:23:15 $
00006  *  $Revision: 1.32 $
00007  *  \author G. Cerminara - INFN Torino
00008  *  revised by G. Mila - INFN Torino
00009  */
00010 
00011 #include "DTSegmentAnalysisTask.h"
00012 
00013 // Framework
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/Framework/interface/LuminosityBlock.h"
00017 #include "FWCore/Framework/interface/ESHandle.h"
00018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00019 #include "FWCore/ServiceRegistry/interface/Service.h"
00020 
00021 #include "DQMServices/Core/interface/DQMStore.h"
00022 #include "DQMServices/Core/interface/MonitorElement.h"
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024 
00025 //Geometry
00026 #include "DataFormats/GeometryVector/interface/Pi.h"
00027 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00028 #include "Geometry/DTGeometry/interface/DTTopology.h"
00029 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00030 //RecHit
00031 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00032 
00033 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00034 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00035 
00036 #include "DQM/DTMonitorModule/interface/DTTimeEvolutionHisto.h"
00037 
00038 #include <iterator>
00039 
00040 using namespace edm;
00041 using namespace std;
00042 
00043 DTSegmentAnalysisTask::DTSegmentAnalysisTask(const edm::ParameterSet& pset) : nevents(0) , nEventsInLS(0), hNevtPerLS(0) {
00044 
00045   edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Constructor called!";
00046 
00047   // switch for detailed analysis
00048   detailedAnalysis = pset.getUntrackedParameter<bool>("detailedAnalysis",false);
00049   // the name of the 4D rec hits collection
00050   theRecHits4DLabel = pset.getParameter<string>("recHits4DLabel");
00051   // Get the map of noisy channels
00052   checkNoisyChannels = pset.getUntrackedParameter<bool>("checkNoisyChannels",false);
00053   // # of bins in the time histos
00054   nTimeBins = pset.getUntrackedParameter<int>("nTimeBins",100);
00055   // # of LS per bin in the time histos
00056   nLSTimeBin = pset.getUntrackedParameter<int>("nLSTimeBin",2);
00057   // switch on/off sliding bins in time histos
00058   slideTimeBins = pset.getUntrackedParameter<bool>("slideTimeBins",true);
00059 
00060   // Get the DQM needed services
00061   theDbe = edm::Service<DQMStore>().operator->();
00062   
00063   // top folder for the histograms in DQMStore
00064   topHistoFolder = pset.getUntrackedParameter<string>("topHistoFolder","DT/02-Segments");
00065   // hlt DQM mode
00066   hltDQMMode = pset.getUntrackedParameter<bool>("hltDQMMode",false);
00067 
00068  }
00069 
00070 
00071 DTSegmentAnalysisTask::~DTSegmentAnalysisTask(){
00072     edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Destructor called!";
00073 }
00074 
00075 
00076 void DTSegmentAnalysisTask::beginRun(const Run& run, const edm::EventSetup& context){ 
00077 
00078    // Get the DT Geometry
00079   context.get<MuonGeometryRecord>().get(dtGeom);
00080 
00081   if (!hltDQMMode) {
00082     theDbe->setCurrentFolder("DT/EventInfo/Counters");
00083     nEventMonitor = theDbe->bookFloat("nProcessedEventsSegment");
00084   }
00085     
00086   // loop over all the DT chambers & book the histos
00087   vector<DTChamber*> chambers = dtGeom->chambers();
00088   vector<DTChamber*>::const_iterator ch_it = chambers.begin();
00089   vector<DTChamber*>::const_iterator ch_end = chambers.end();
00090   for (; ch_it != ch_end; ++ch_it) {
00091     bookHistos((*ch_it)->id());
00092   }
00093 
00094   // book sector time-evolution histos
00095   int modeTimeHisto = 0;
00096   if(!slideTimeBins) modeTimeHisto = 1;
00097   for(int wheel = -2; wheel != 3; ++wheel) { // loop over wheels
00098     for(int sector = 1; sector <= 12; ++sector) { // loop over sectors
00099       stringstream wheelstr; wheelstr << wheel; 
00100       stringstream sectorstr; sectorstr << sector;
00101       string sectorHistoName = "NSegmPerEvent_W" + wheelstr.str()
00102         + "_Sec" + sectorstr.str();
00103       string sectorHistoTitle = "# segm. W" + wheelstr.str() + " Sect." + sectorstr.str();
00104 
00105       theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheelstr.str() +
00106                                "/Sector" + sectorstr.str());
00107       histoTimeEvol[wheel][sector] = new DTTimeEvolutionHisto(&(*theDbe),sectorHistoName,sectorHistoTitle,
00108                                                               nTimeBins,nLSTimeBin,slideTimeBins,modeTimeHisto);
00109     }
00110   }
00111 
00112   if(hltDQMMode) theDbe->setCurrentFolder(topHistoFolder);
00113   else theDbe->setCurrentFolder("DT/EventInfo/");
00114 
00115   hNevtPerLS = new DTTimeEvolutionHisto(&(*theDbe),"NevtPerLS","# evt.",nTimeBins,nLSTimeBin,slideTimeBins,2);
00116 
00117 }
00118 
00119 
00120 void DTSegmentAnalysisTask::endJob() {
00121  
00122   edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") <<"[DTSegmentAnalysisTask] endjob called!";
00123 
00124   delete hNevtPerLS;
00125   //theDbe->save("testMonitoring.root");
00126 
00127 }
00128   
00129 
00130 
00131 void DTSegmentAnalysisTask::analyze(const edm::Event& event, const edm::EventSetup& setup) {
00132 
00133   nevents++;
00134   nEventMonitor->Fill(nevents);
00135 
00136   nEventsInLS++;
00137   edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Analyze #Run: " << event.id().run()
00138                                << " #Event: " << event.id().event();
00139   if(!(event.id().event()%1000))
00140     edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Analyze #Run: " << event.id().run()
00141                                  << " #Event: " << event.id().event();
00142   
00143   ESHandle<DTStatusFlag> statusMap;
00144   if(checkNoisyChannels) {
00145     setup.get<DTStatusFlagRcd>().get(statusMap);
00146   } 
00147 
00148 
00149   // -- 4D segment analysis  -----------------------------------------------------
00150   
00151   // Get the 4D segment collection from the event
00152   edm::Handle<DTRecSegment4DCollection> all4DSegments;
00153   event.getByLabel(theRecHits4DLabel, all4DSegments);
00154 
00155   if(!all4DSegments.isValid()) return;
00156   
00157   // Loop over all chambers containing a segment
00158   DTRecSegment4DCollection::id_iterator chamberId;
00159   for (chamberId = all4DSegments->id_begin();
00160        chamberId != all4DSegments->id_end();
00161        ++chamberId){
00162     // Get the range for the corresponding ChamerId
00163     DTRecSegment4DCollection::range  range = all4DSegments->get(*chamberId);
00164 
00165     edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "   Chamber: " << *chamberId << " has " << distance(range.first, range.second)
00166                                  << " 4D segments";
00167 
00168     // Loop over the rechits of this ChamerId
00169     for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
00170          segment4D!=range.second;
00171            ++segment4D){
00172 
00173       //FOR NOISY CHANNELS////////////////////////////////
00174      bool segmNoisy = false;
00175      if(checkNoisyChannels) {
00176        
00177        if((*segment4D).hasPhi()){
00178          const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
00179          vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
00180          map<DTSuperLayerId,vector<DTRecHit1D> > hitsBySLMap; 
00181          for(vector<DTRecHit1D>::const_iterator hit = phiHits.begin();
00182              hit != phiHits.end(); ++hit) {
00183            DTWireId wireId = (*hit).wireId();
00184            
00185            // Check for noisy channels to skip them
00186            bool isNoisy = false;
00187            bool isFEMasked = false;
00188            bool isTDCMasked = false;
00189            bool isTrigMask = false;
00190            bool isDead = false;
00191            bool isNohv = false;
00192            statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
00193            if(isNoisy) {
00194              edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "Wire: " << wireId << " is noisy, skipping!";
00195              segmNoisy = true;
00196            }      
00197          }
00198        }
00199        
00200        if((*segment4D).hasZed()) {
00201          const DTSLRecSegment2D* zSeg = (*segment4D).zSegment();  // zSeg lives in the SL RF
00202          // Check for noisy channels to skip them
00203          vector<DTRecHit1D> zHits = zSeg->specificRecHits();
00204          for(vector<DTRecHit1D>::const_iterator hit = zHits.begin();
00205              hit != zHits.end(); ++hit) {
00206            DTWireId wireId = (*hit).wireId();
00207            bool isNoisy = false;
00208            bool isFEMasked = false;
00209            bool isTDCMasked = false;
00210            bool isTrigMask = false;
00211            bool isDead = false;
00212            bool isNohv = false;
00213            statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
00214            if(isNoisy) {
00215              edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "Wire: " << wireId << " is noisy, skipping!";
00216              segmNoisy = true;
00217            }     
00218          }
00219        } 
00220 
00221      } // end of switch on noisy channels
00222      if (segmNoisy) {
00223        edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")<<"skipping the segment: it contains noisy cells";
00224        continue;
00225      }
00226      //END FOR NOISY CHANNELS////////////////////////////////
00227       
00228      int nHits=0;
00229      if((*segment4D).hasPhi())
00230        nHits = (((*segment4D).phiSegment())->specificRecHits()).size();
00231      if((*segment4D).hasZed()) 
00232        nHits = nHits + ((((*segment4D).zSegment())->specificRecHits()).size());
00233       
00234      fillHistos(*chamberId,
00235                 nHits,
00236                 (*segment4D).chi2()/(*segment4D).degreesOfFreedom());
00237     }
00238   }
00239 
00240   // -----------------------------------------------------------------------------
00241 }
00242   
00243 
00244 // Book a set of histograms for a give chamber
00245 void DTSegmentAnalysisTask::bookHistos(DTChamberId chamberId) {
00246 
00247   edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "   Booking histos for chamber: " << chamberId;
00248 
00249 
00250   // Compose the chamber name
00251   stringstream wheel; wheel << chamberId.wheel();       
00252   stringstream station; station << chamberId.station(); 
00253   stringstream sector; sector << chamberId.sector();
00254   
00255   string chamberHistoName =
00256     "_W" + wheel.str() +
00257     "_St" + station.str() +
00258     "_Sec" + sector.str();
00259   
00260 
00261   for(int wh=-2; wh<=2; wh++){
00262     stringstream wheel; wheel << wh;
00263     theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str());
00264     string histoName =  "numberOfSegments_W" + wheel.str();
00265     summaryHistos[wh] = theDbe->book2D(histoName.c_str(),histoName.c_str(),12,1,13,4,1,5);
00266     summaryHistos[wh]->setAxisTitle("Sector",1);
00267     summaryHistos[wh]->setBinLabel(1,"1",1);
00268     summaryHistos[wh]->setBinLabel(2,"2",1);
00269     summaryHistos[wh]->setBinLabel(3,"3",1);
00270     summaryHistos[wh]->setBinLabel(4,"4",1);
00271     summaryHistos[wh]->setBinLabel(5,"5",1);
00272     summaryHistos[wh]->setBinLabel(6,"6",1);
00273     summaryHistos[wh]->setBinLabel(7,"7",1);
00274     summaryHistos[wh]->setBinLabel(8,"8",1);
00275     summaryHistos[wh]->setBinLabel(9,"9",1);
00276     summaryHistos[wh]->setBinLabel(10,"10",1);
00277     summaryHistos[wh]->setBinLabel(11,"11",1);
00278     summaryHistos[wh]->setBinLabel(12,"12",1);
00279     summaryHistos[wh]->setBinLabel(1,"MB1",2);
00280     summaryHistos[wh]->setBinLabel(2,"MB2",2);
00281     summaryHistos[wh]->setBinLabel(3,"MB3",2);
00282     summaryHistos[wh]->setBinLabel(4,"MB4",2);  
00283   }
00284 
00285 
00286   theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str() +
00287                            "/Sector" + sector.str() +
00288                            "/Station" + station.str());
00289 
00290   // Create the monitor elements
00291   vector<MonitorElement *> histos;
00292   histos.push_back(theDbe->book1D("h4DSegmNHits"+chamberHistoName,
00293                                   "# of hits per segment",
00294                                   16, 0.5, 16.5));
00295   if(detailedAnalysis){
00296     histos.push_back(theDbe->book1D("h4DChi2"+chamberHistoName,
00297                                     "4D Segment reduced Chi2",
00298                                     20, 0, 20));
00299   }
00300   histosPerCh[chamberId] = histos;
00301 }
00302 
00303 
00304 // Fill a set of histograms for a give chamber 
00305 void DTSegmentAnalysisTask::fillHistos(DTChamberId chamberId,
00306                                    int nHits,
00307                                    float chi2) {
00308   int sector = chamberId.sector();
00309   if(chamberId.sector()==13) {
00310     sector = 4;
00311   } else if(chamberId.sector()==14) {
00312      sector = 10;
00313   }
00314 
00315   summaryHistos[chamberId.wheel()]->Fill(sector,chamberId.station());
00316   histoTimeEvol[chamberId.wheel()][sector]->accumulateValueTimeSlot(1);
00317 
00318   vector<MonitorElement *> histos =  histosPerCh[chamberId];                          
00319   histos[0]->Fill(nHits);
00320   if(detailedAnalysis){
00321     histos[1]->Fill(chi2);
00322   }
00323 
00324 }
00325 
00326 
00327 void DTSegmentAnalysisTask::endLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& eSetup) {
00328 
00329   hNevtPerLS->updateTimeSlot(lumiSeg.luminosityBlock(), nEventsInLS);
00330  // book sector time-evolution histos
00331   for(int wheel = -2; wheel != 3; ++wheel) {
00332     for(int sector = 1; sector <= 12; ++sector) {
00333       histoTimeEvol[wheel][sector]->updateTimeSlot(lumiSeg.luminosityBlock(), nEventsInLS);
00334     }
00335   }
00336 }
00337 
00338 
00339 void DTSegmentAnalysisTask::beginLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& eSetup) {
00340   nEventsInLS = 0;
00341 }
00342 
00343