CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/CalibMuon/DTCalibration/plugins/DTNoiseCalibration.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2011/05/27 12:39:41 $
00005  *  $Revision: 1.15 $
00006  *  \author G. Mila - INFN Torino
00007  */
00008 
00009 #include "DTNoiseCalibration.h"
00010 
00011 // Framework
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/EventSetup.h" 
00014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 
00017 // Geometry
00018 #include "Geometry/DTGeometry/interface/DTLayer.h"
00019 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00020 #include "Geometry/DTGeometry/interface/DTTopology.h"
00021 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00022 
00023 // Digis
00024 #include "DataFormats/DTDigi/interface/DTDigi.h"
00025 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
00026 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00027 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00028 #include "CondFormats/DTObjects/interface/DTReadOutMapping.h"
00029 
00030 // Database
00031 #include "CondFormats/DTObjects/interface/DTTtrig.h"
00032 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
00033 
00034 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
00035 
00036 #include "TH2F.h"
00037 #include "TFile.h"
00038 
00039 using namespace edm;
00040 using namespace std;
00041 
00042 DTNoiseCalibration::DTNoiseCalibration(const edm::ParameterSet& pset):
00043   digiLabel_( pset.getParameter<InputTag>("digiLabel") ),
00044   useTimeWindow_( pset.getParameter<bool>("useTimeWindow") ),
00045   triggerWidth_( pset.getParameter<int>("triggerWidth") ),
00046   timeWindowOffset_( pset.getParameter<int>("timeWindowOffset") ),
00047   maximumNoiseRate_( pset.getParameter<double>("maximumNoiseRate") ),
00048   useAbsoluteRate_( pset.getParameter<bool>("useAbsoluteRate") ), 
00049   readDB_(true), defaultTtrig_(0), 
00050   dbLabel_( pset.getUntrackedParameter<string>("dbLabel", "") ),
00051   //fastAnalysis_( pset.getParameter<bool>("fastAnalysis", true) ),
00052   wireIdWithHisto_( std::vector<DTWireId>() ),
00053   lumiMax_(5000)
00054   {
00055 
00056   // Get the debug parameter for verbose output
00057   //debug = ps.getUntrackedParameter<bool>("debug");
00058   /*// The analysis type
00059   // The wheel & sector interested for the time-dependent analysis
00060   wh = ps.getUntrackedParameter<int>("wheel", 0);
00061   sect = ps.getUntrackedParameter<int>("sector", 6);*/
00062 
00063   if( pset.exists("defaultTtrig") ){
00064      readDB_ = false;
00065      defaultTtrig_ = pset.getParameter<int>("defaultTtrig");
00066   }
00067 
00068   if( pset.exists("cellsWithHisto") ){
00069      vector<string> cellsWithHisto = pset.getParameter<vector<string> >("cellsWithHisto");
00070      for(vector<string>::const_iterator cell = cellsWithHisto.begin(); cell != cellsWithHisto.end(); ++cell){
00071         //FIXME: Use regex to check whether format is right
00072         if( (*cell) != "" && (*cell) != "None"){
00073            stringstream linestr;
00074            int wheel,station,sector,sl,layer,wire;
00075            linestr << (*cell);
00076            linestr >> wheel >> station >> sector >> sl >> layer >> wire;
00077            wireIdWithHisto_.push_back(DTWireId(wheel,station,sector,sl,layer,wire));
00078         }
00079      }
00080   }
00081 
00082   // The root file which will contain the histos
00083   string rootFileName = pset.getUntrackedParameter<string>("rootFileName","noise.root");
00084   rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
00085   rootFile_->cd();
00086 }
00087 
00088 void DTNoiseCalibration::beginJob() {
00089   LogVerbatim("Calibration") << "[DTNoiseCalibration]: Begin job";
00090   
00091   nevents_ = 0;
00092   
00093   TH1::SetDefaultSumw2(true);
00094   int numBin = (triggerWidth_*32/25)/50;
00095   hTDCTriggerWidth_ = new TH1F("TDC_Time_Distribution", "TDC_Time_Distribution", numBin, 0, triggerWidth_*32/25);
00096 }
00097 
00098 void DTNoiseCalibration::beginRun(const edm::Run& run, const edm::EventSetup& setup ) {
00099 
00100   // Get the DT Geometry
00101   setup.get<MuonGeometryRecord>().get(dtGeom_);
00102 
00103   // tTrig 
00104   if( readDB_ ) setup.get<DTTtrigRcd>().get(dbLabel_,tTrigMap_);
00105 
00106   /*
00107   nevents = 0;
00108   counter = 0;
00109 
00110   // TDC time distribution
00111   int numBin = (triggerWidth_*(32/25))/50;
00112   hTDCTriggerWidth = new TH1F("TDC_Time_Distribution", "TDC_Time_Distribution", numBin, 0, triggerWidth_*(32/25));*/
00113 
00114 }
00115 
00116 void DTNoiseCalibration::analyze(const edm::Event& event, const edm::EventSetup& setup){
00117   ++nevents_;
00118   
00119   // Get the digis from the event
00120   Handle<DTDigiCollection> dtdigis;
00121   event.getByLabel(digiLabel_, dtdigis);
00122 
00123   /*TH1F *hOccupancyHisto;
00124   TH2F *hEvtPerWireH;
00125   string Histo2Name;*/
00126 
00127   // Loop over digis
00128   DTDigiCollection::DigiRangeIterator dtLayerId_It;
00129   for (dtLayerId_It=dtdigis->begin(); dtLayerId_It!=dtdigis->end(); ++dtLayerId_It){
00130      // Define time window
00131      float upperLimit = 0.;
00132      if( readDB_ ){
00133         float tTrig,tTrigRMS,kFactor;
00134         DTSuperLayerId slId = ((*dtLayerId_It).first).superlayerId();
00135         int status = tTrigMap_->get( slId, tTrig, tTrigRMS, kFactor, DTTimeUnits::counts );
00136         if(status != 0) throw cms::Exception("DTNoiseCalibration") << "Could not find tTrig entry in DB for" << slId << endl;
00137         upperLimit = tTrig - timeWindowOffset_;
00138      } else {
00139         upperLimit = defaultTtrig_ - timeWindowOffset_;
00140      }
00141 
00142      for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
00143          digiIt!=((*dtLayerId_It).second).second; ++digiIt){
00144 
00145         //Check the TDC trigger width
00146         int tdcTime = (*digiIt).countsTDC();
00147         if( !useTimeWindow_ ){
00148            if( ( ((float)tdcTime*25)/32 ) > triggerWidth_ )
00149               LogError("Calibration") << "Digi has a TDC time (ns) higher than the pre-defined TDC trigger width: " << ((float)tdcTime*25)/32;
00150         }       
00151 
00152         hTDCTriggerWidth_->Fill(tdcTime);
00153 
00154         if( useTimeWindow_ && tdcTime > upperLimit) continue;
00155 
00156         /*LogTrace("Calibration") << "TDC time (ns): " << ((float)tdcTime*25)/32
00157                                 <<" --- trigger width (ns): " << ((float)upperLimit*25)/32;*/
00158 
00159         const DTLayerId dtLId = (*dtLayerId_It).first;
00160         const DTTopology& dtTopo = dtGeom_->layer(dtLId)->specificTopology();
00161         const int firstWire = dtTopo.firstChannel();
00162         const int lastWire = dtTopo.lastChannel();
00163         //const int nWires = dtTopo.channels();
00164         const int nWires = lastWire - firstWire + 1;
00165 
00166         // Book the occupancy histos
00167         if( theHistoOccupancyMap_.find(dtLId) == theHistoOccupancyMap_.end() ){
00168            string histoName = "DigiOccupancy_" + getLayerName(dtLId);
00169            rootFile_->cd();
00170            TH1F* hOccupancyHisto = new TH1F(histoName.c_str(), histoName.c_str(), nWires, firstWire, lastWire+1);
00171            LogTrace("Calibration") << "  Created occupancy Histo: " << hOccupancyHisto->GetName();
00172            theHistoOccupancyMap_[dtLId] = hOccupancyHisto;
00173         }
00174         theHistoOccupancyMap_[dtLId]->Fill((*digiIt).wire());
00175 
00176         const DTWireId wireId(dtLId, (*digiIt).wire());
00177         if( find(wireIdWithHisto_.begin(),wireIdWithHisto_.end(),wireId) != wireIdWithHisto_.end() ){
00178            if( theHistoOccupancyVsLumiMap_.find(wireId) == theHistoOccupancyVsLumiMap_.end() ){
00179               string histoName = "OccupancyVsLumi_" + getChannelName(wireId);
00180               rootFile_->cd();
00181               TH1F* hOccupancyVsLumiHisto = new TH1F(histoName.c_str(), histoName.c_str(), lumiMax_, 0, lumiMax_);
00182               LogTrace("Calibration") << "  Created occupancy histo: " << hOccupancyVsLumiHisto->GetName();
00183               theHistoOccupancyVsLumiMap_[wireId] = hOccupancyVsLumiHisto;
00184            }
00185 
00186            unsigned int lumiSection = event.luminosityBlock();
00187            theHistoOccupancyVsLumiMap_[wireId]->Fill(lumiSection);
00188         }
00189 
00190         /*// Book the digi event plot every 1000 events if the analysis is not "fast" and if is the correct sector
00191         if(!fastAnalysis &&
00192            dtLId.superlayerId().chamberId().wheel()==wh &&
00193            dtLId.superlayerId().chamberId().sector()==sect) {
00194           if(theHistoEvtPerWireMap.find(dtLId) == theHistoEvtPerWireMap.end() ||
00195              (theHistoEvtPerWireMap.find(dtLId) != theHistoEvtPerWireMap.end() &&
00196               skippedPlot[dtLId] != counter)){ 
00197             skippedPlot[dtLId] = counter;
00198             stringstream toAppend; toAppend << counter;
00199             Histo2Name = "DigiPerWirePerEvent_" + getLayerName(dtLId) + "_" + toAppend.str();
00200             theFile->cd();
00201             hEvtPerWireH = new TH2F(Histo2Name.c_str(), Histo2Name.c_str(), 1000,0.5,1000.5,nWires, firstWire, lastWire+1);
00202             if(hEvtPerWireH){
00203               if(debug)
00204                 cout << "  New Histo with the number of digi per evt per wire: " << hEvtPerWireH->GetName() << endl;
00205               theHistoEvtPerWireMap[dtLId]=hEvtPerWireH;
00206             }
00207           }
00208         }*/
00209      }
00210   }
00211     
00212   /*//Fill the plot of the number of digi per event per wire
00213   std::map<int,int > DigiPerWirePerEvent;
00214   // LOOP OVER ALL THE CHAMBERS
00215   vector<DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
00216   vector<DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
00217   for (; ch_it != ch_end; ++ch_it) {
00218     DTChamberId ch = (*ch_it)->id();
00219     vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin(); 
00220     vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
00221     // Loop over the SLs
00222     for(; sl_it != sl_end; ++sl_it) {
00223       DTSuperLayerId sl = (*sl_it)->id();
00224       vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin(); 
00225       vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
00226       // Loop over the Ls
00227       for(; l_it != l_end; ++l_it) {
00228         DTLayerId layerId = (*l_it)->id();
00229         
00230         // Get the number of wires
00231         const DTTopology& dtTopo = dtGeom->layer(layerId)->specificTopology();
00232         const int firstWire = dtTopo.firstChannel();
00233         const int lastWire = dtTopo.lastChannel();
00234           
00235         if (theHistoEvtPerWireMap.find(layerId) != theHistoEvtPerWireMap.end() &&
00236             skippedPlot[layerId] == counter) {
00237           
00238           for (int wire=firstWire; wire<=lastWire; wire++) {
00239             DigiPerWirePerEvent[wire]= 0;
00240           }     
00241           // loop over all the digis of the event
00242           DTDigiCollection::Range layerDigi= dtdigis->get(layerId);
00243           for (DTDigiCollection::const_iterator digi = layerDigi.first;
00244                digi!=layerDigi.second;
00245                ++digi){
00246             if((cosmicRun && (*digi).countsTDC()<upperLimit) || (!cosmicRun))
00247               DigiPerWirePerEvent[(*digi).wire()]+=1;
00248           }
00249           // fill the digi event histo
00250           for (int wire=firstWire; wire<=lastWire; wire++) {
00251             theFile->cd();
00252             int histoEvents = nevents - (counter*1000);
00253             theHistoEvtPerWireMap[layerId]->Fill(histoEvents,wire,DigiPerWirePerEvent[wire]);
00254           }
00255         }
00256       } //Loop Ls
00257     } //Loop SLs
00258   } //Loop chambers
00259   
00260   
00261   if(nevents % 1000 == 0) {
00262     counter++;
00263     // save the digis event plot on file
00264     for(map<DTLayerId,  TH2F* >::const_iterator lHisto = theHistoEvtPerWireMap.begin();
00265         lHisto != theHistoEvtPerWireMap.end();
00266         lHisto++) {
00267       theFile->cd();
00268       if((*lHisto).second)
00269         (*lHisto).second->Write();
00270     }
00271     theHistoEvtPerWireMap.clear();
00272   }*/
00273   
00274 }
00275 
00276 void DTNoiseCalibration::endJob(){
00277 
00278   //LogVerbatim("Calibration") << "[DTNoiseCalibration] endjob called!";
00279   LogVerbatim("Calibration") << "[DTNoiseCalibration] Total number of events analyzed: " << nevents_;
00280 
00281   // Save the TDC digi plot
00282   rootFile_->cd();
00283   hTDCTriggerWidth_->Write();
00284 
00285   for(map<DTWireId, TH1F*>::const_iterator wHisto = theHistoOccupancyVsLumiMap_.begin();
00286       wHisto != theHistoOccupancyVsLumiMap_.end(); ++wHisto) (*wHisto).second->Write();
00287   
00288   // Save on file the occupancy histos and write the list of noisy cells
00289   DTStatusFlag *statusMap = new DTStatusFlag();
00290   for(map<DTLayerId, TH1F*>::const_iterator lHisto = theHistoOccupancyMap_.begin();
00291       lHisto != theHistoOccupancyMap_.end();
00292       ++lHisto){
00293      double triggerWidth_s = 0.;
00294      if( useTimeWindow_ ){
00295         double triggerWidth_ns = 0.;
00296         if( readDB_ ){
00297            float tTrig, tTrigRMS, kFactor;
00298            DTSuperLayerId slId = ((*lHisto).first).superlayerId();
00299            int status = tTrigMap_->get( slId, tTrig, tTrigRMS, kFactor, DTTimeUnits::counts );
00300            if(status != 0) throw cms::Exception("DTNoiseCalibration") << "Could not find tTrig entry in DB for" << slId << endl;
00301            triggerWidth_ns = tTrig - timeWindowOffset_;
00302         } else{
00303            triggerWidth_ns = defaultTtrig_ - timeWindowOffset_;
00304         }
00305         triggerWidth_ns = (triggerWidth_ns*25)/32;
00306         triggerWidth_s = triggerWidth_ns/1e9;
00307      } else{
00308         triggerWidth_s = double(triggerWidth_/1e9);
00309      }
00310      LogTrace("Calibration") << (*lHisto).second->GetName() << " trigger width (s): " << triggerWidth_s;
00311 
00312      double normalization = 1./(nevents_*triggerWidth_s);
00313      if((*lHisto).second){
00314         (*lHisto).second->Scale(normalization);
00315         rootFile_->cd();
00316         (*lHisto).second->Write();
00317         const DTTopology& dtTopo = dtGeom_->layer((*lHisto).first)->specificTopology();
00318         const int firstWire = dtTopo.firstChannel();
00319         const int lastWire = dtTopo.lastChannel();
00320         //const int nWires = dtTopo.channels();
00321         const int nWires = lastWire - firstWire + 1;
00322         // Find average in layer
00323         double averageRate = 0.;  
00324         for(int bin = 1; bin <= (*lHisto).second->GetNbinsX(); ++bin)
00325            averageRate += (*lHisto).second->GetBinContent(bin);
00326 
00327         if(nWires) averageRate /= nWires;  
00328         LogTrace("Calibration") << "  Average rate = " << averageRate;
00329 
00330         for(int i_wire = firstWire; i_wire <= lastWire; ++i_wire){
00331            // From definition of "noisy cell"
00332            int bin = i_wire - firstWire + 1;
00333            double channelRate = (*lHisto).second->GetBinContent(bin);
00334            double rateOffset = (useAbsoluteRate_) ? 0. : averageRate;
00335            if( (channelRate - rateOffset) > maximumNoiseRate_ ){
00336               DTWireId wireID((*lHisto).first, i_wire);
00337               statusMap->setCellNoise(wireID,1);
00338               LogVerbatim("Calibration") << ">>> Channel noisy: " << wireID;
00339            }
00340         }
00341      }
00342   }
00343   LogVerbatim("Calibration") << "Writing noise map object to DB";
00344   string record = "DTStatusFlagRcd";
00345   DTCalibDBUtils::writeToDB<DTStatusFlag>(record, statusMap);
00346 }
00347 
00348 DTNoiseCalibration::~DTNoiseCalibration(){
00349   rootFile_->Close();
00350 }
00351 
00352 string DTNoiseCalibration::getChannelName(const DTWireId& wId) const{
00353   stringstream channelName;
00354   channelName << "Wh" << wId.wheel() << "_St" << wId.station() << "_Sec" << wId.sector()
00355               << "_SL" << wId.superlayer() << "_L" << wId.layer() << "_W"<< wId.wire();
00356 
00357   return channelName.str();
00358 }
00359 
00360 string DTNoiseCalibration::getLayerName(const DTLayerId& lId) const{
00361 
00362   const  DTSuperLayerId dtSLId = lId.superlayerId();
00363   const  DTChamberId dtChId = dtSLId.chamberId(); 
00364   stringstream Layer; Layer << lId.layer();
00365   stringstream superLayer; superLayer << dtSLId.superlayer();
00366   stringstream wheel; wheel << dtChId.wheel();  
00367   stringstream station; station << dtChId.station();    
00368   stringstream sector; sector << dtChId.sector();
00369   
00370   string LayerName = 
00371     "W" + wheel.str()
00372     + "_St" + station.str() 
00373     + "_Sec" + sector.str() 
00374     + "_SL" + superLayer.str()
00375     + "_L" + Layer.str();
00376   
00377   return LayerName;
00378 }
00379 
00380 string DTNoiseCalibration::getSuperLayerName(const DTSuperLayerId& dtSLId) const{
00381 
00382   const  DTChamberId dtChId = dtSLId.chamberId(); 
00383   stringstream superLayer; superLayer << dtSLId.superlayer();
00384   stringstream wheel; wheel << dtChId.wheel();  
00385   stringstream station; station << dtChId.station();    
00386   stringstream sector; sector << dtChId.sector();
00387   
00388   string SuperLayerName = 
00389     "W" + wheel.str()
00390     + "_St" + station.str() 
00391     + "_Sec" + sector.str() 
00392     + "_SL" + superLayer.str();
00393   
00394   return SuperLayerName;
00395 }