CMS 3D CMS Logo

DTNoiseCalibration.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2008/08/19 13:33:00 $
00005  *  $Revision: 1.7 $
00006  *  \author G. Mila - INFN Torino
00007  */
00008 
00009 
00010 #include "CalibMuon/DTCalibration/plugins/DTNoiseCalibration.h"
00011 
00012 
00013 // Framework
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include <FWCore/Framework/interface/EventSetup.h>
00017 
00018 // Geometry
00019 #include "Geometry/DTGeometry/interface/DTLayer.h"
00020 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00021 #include "Geometry/DTGeometry/interface/DTTopology.h"
00022 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00023 
00024 // Digis
00025 #include <DataFormats/DTDigi/interface/DTDigi.h>
00026 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
00027 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00028 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00029 #include "CondFormats/DTObjects/interface/DTReadOutMapping.h"
00030 
00031 // Database
00032 #include <CondFormats/DTObjects/interface/DTTtrig.h>
00033 #include <CondFormats/DataRecord/interface/DTTtrigRcd.h>
00034 
00035 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
00036 
00037 #include "TH2F.h"
00038 #include "TFile.h"
00039 
00040 using namespace edm;
00041 using namespace std;
00042 
00043 
00044 DTNoiseCalibration::DTNoiseCalibration(const edm::ParameterSet& ps){
00045 
00046   cout << "[DTNoiseCalibration]: Constructor" <<endl;
00047 
00048   // Get the debug parameter for verbose output
00049   debug = ps.getUntrackedParameter<bool>("debug");
00050 
00051   // The analysis type
00052   fastAnalysis = ps.getUntrackedParameter<bool>("fastAnalysis", true);
00053   // The wheel & sector interested for the time-dependent analysis
00054   wh = ps.getUntrackedParameter<int>("wheel", 0);
00055   sect = ps.getUntrackedParameter<int>("sector", 6);
00056 
00057   // The trigger mode
00058   cosmicRun = ps.getUntrackedParameter<bool>("cosmicRun", false);
00059 
00060   // The trigger width (if noise run)
00061   TriggerWidth = ps.getUntrackedParameter<int>("TriggerWidth");
00062 
00063   // The root file which will contain the histos
00064   string rootFileName = ps.getUntrackedParameter<string>("rootFileName");
00065   theFile = new TFile(rootFileName.c_str(), "RECREATE");
00066   theFile->cd();
00067 
00068   parameters=ps;
00069 
00070 }
00071 
00072 
00073 void DTNoiseCalibration::beginJob(const edm::EventSetup& context){
00074 
00075   cout <<"[DTNoiseCalibration]: BeginJob"<<endl; 
00076   nevents = 0;
00077   counter = 0;
00078 
00079   // Get the DT Geometry
00080   context.get<MuonGeometryRecord>().get(dtGeom);
00081 
00082   // tTrig 
00083   if (parameters.getUntrackedParameter<bool>("readDB", true)) 
00084     context.get<DTTtrigRcd>().get(tTrigMap);
00085 
00086   // TDC time distribution
00087   int numBin = (TriggerWidth*(32/25))/50;
00088   hTDCTriggerWidth = new TH1F("TDC_Time_Distribution", "TDC_Time_Distribution", numBin, 0, TriggerWidth*(32/25));
00089 
00090 }
00091 
00092 
00093 void DTNoiseCalibration::analyze(const edm::Event& e, const edm::EventSetup& context){
00094 
00095   nevents++;
00096   if(debug)
00097     cout<<"nevents: "<<nevents<<endl;
00098   
00099   // Get the digis from the event
00100   edm::Handle<DTDigiCollection> dtdigis;
00101   e.getByLabel("dtunpacker", dtdigis);
00102 
00103   TH1F *hOccupancyHisto;
00104   TH2F *hEvtPerWireH;
00105   string Histo2Name;
00106 
00107   // LOOP OVER ALL THE DIGIS OF THE EVENT
00108   DTDigiCollection::DigiRangeIterator dtLayerId_It;
00109   for (dtLayerId_It=dtdigis->begin(); dtLayerId_It!=dtdigis->end(); ++dtLayerId_It){
00110     for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
00111          digiIt!=((*dtLayerId_It).second).second; ++digiIt){
00112 
00113       //Check the TDC trigger width
00114       int tdcTime = (*digiIt).countsTDC();
00115       if(!cosmicRun){   
00116         if(debug)
00117           cout<<"tdcTime (ns): "<<(tdcTime*25)/32<<endl;
00118         if(((tdcTime*25)/32)>TriggerWidth){
00119           cout<<"***Error*** : your digi has a tdcTime (ns) higher than the TDC trigger width :"<<(tdcTime*25)/32<<endl;
00120           abort();
00121         }
00122       } 
00123 
00124       if((!fastAnalysis &&
00125          (*dtLayerId_It).first.superlayerId().chamberId().wheel()==wh &&
00126          (*dtLayerId_It).first.superlayerId().chamberId().sector()==sect) ||
00127          fastAnalysis)
00128         hTDCTriggerWidth->Fill(tdcTime);
00129 
00130       // Set the window of interest if the run is triggered by cosmics
00131       if ( parameters.getUntrackedParameter<bool>("readDB", true) ) {
00132         tTrigMap->slTtrig( ((*dtLayerId_It).first).superlayerId(), tTrig, tTrigRMS);
00133         upperLimit = tTrig-500;
00134       }
00135       else { 
00136         tTrig = parameters.getUntrackedParameter<int>("defaultTtrig", 4000);
00137         upperLimit = tTrig-500;
00138       }
00139         
00140       if((cosmicRun && (*digiIt).countsTDC()<upperLimit) || (!cosmicRun) ){
00141 
00142         if(debug && cosmicRun)
00143           cout<<"tdcTime (ns): "<<((*digiIt).countsTDC()*25)/32<<" --- TriggerWidth (ns): "<<(upperLimit*25)/32<<endl;
00144 
00145         // Get the number of wires
00146         const  DTLayerId dtLId = (*dtLayerId_It).first;
00147         const DTTopology& dtTopo = dtGeom->layer(dtLId)->specificTopology();
00148         const int nWires = dtTopo.channels();
00149         const int firstWire = dtTopo.firstChannel();
00150         const int lastWire = dtTopo.lastChannel();
00151         
00152         // book the occupancy histos
00153         theFile->cd();
00154         if((!fastAnalysis &&
00155            dtLId.superlayerId().chamberId().wheel()==wh &&
00156            dtLId.superlayerId().chamberId().sector()==sect) ||
00157            fastAnalysis){
00158           hOccupancyHisto = theHistoOccupancyMap[dtLId];
00159           if(hOccupancyHisto == 0) {
00160             string HistoName = "DigiOccupancy_" + getLayerName(dtLId);
00161             theFile->cd();
00162             hOccupancyHisto = new TH1F(HistoName.c_str(), HistoName.c_str(), nWires, firstWire, lastWire+1);
00163             if(debug)
00164               cout << "  New Occupancy Histo: " << hOccupancyHisto->GetName() << endl;
00165             theHistoOccupancyMap[dtLId] = hOccupancyHisto;
00166           }
00167           hOccupancyHisto->Fill((*digiIt).wire());
00168         }
00169 
00170         // book the digi event plot every 1000 events if the analysis is not "fast" and if is the correct sector
00171         if(!fastAnalysis &&
00172            dtLId.superlayerId().chamberId().wheel()==wh &&
00173            dtLId.superlayerId().chamberId().sector()==sect) {
00174           if(theHistoEvtPerWireMap.find(dtLId) == theHistoEvtPerWireMap.end() ||
00175              (theHistoEvtPerWireMap.find(dtLId) != theHistoEvtPerWireMap.end() &&
00176               skippedPlot[dtLId] != counter)){ 
00177             skippedPlot[dtLId] = counter;
00178             stringstream toAppend; toAppend << counter;
00179             Histo2Name = "DigiPerWirePerEvent_" + getLayerName(dtLId) + "_" + toAppend.str();
00180             theFile->cd();
00181             hEvtPerWireH = new TH2F(Histo2Name.c_str(), Histo2Name.c_str(), 1000,0.5,1000.5,nWires, firstWire, lastWire+1);
00182             if(hEvtPerWireH){
00183               if(debug)
00184                 cout << "  New Histo with the number of digi per evt per wire: " << hEvtPerWireH->GetName() << endl;
00185               theHistoEvtPerWireMap[dtLId]=hEvtPerWireH;
00186             }
00187           }
00188         }
00189       }
00190     }
00191   }
00192     
00193   //Fill the plot of the number of digi per event per wire
00194   std::map<int,int > DigiPerWirePerEvent;
00195   // LOOP OVER ALL THE CHAMBERS
00196   vector<DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
00197   vector<DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
00198   for (; ch_it != ch_end; ++ch_it) {
00199     DTChamberId ch = (*ch_it)->id();
00200     vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin(); 
00201     vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
00202     // Loop over the SLs
00203     for(; sl_it != sl_end; ++sl_it) {
00204       DTSuperLayerId sl = (*sl_it)->id();
00205       vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin(); 
00206       vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
00207       // Loop over the Ls
00208       for(; l_it != l_end; ++l_it) {
00209         DTLayerId layerId = (*l_it)->id();
00210         
00211         // Get the number of wires
00212         const DTTopology& dtTopo = dtGeom->layer(layerId)->specificTopology();
00213         const int firstWire = dtTopo.firstChannel();
00214         const int lastWire = dtTopo.lastChannel();
00215           
00216         if (theHistoEvtPerWireMap.find(layerId) != theHistoEvtPerWireMap.end() &&
00217             skippedPlot[layerId] == counter) {
00218           
00219           for (int wire=firstWire; wire<=lastWire; wire++) {
00220             DigiPerWirePerEvent[wire]= 0;
00221           }     
00222           // loop over all the digis of the event
00223           DTDigiCollection::Range layerDigi= dtdigis->get(layerId);
00224           for (DTDigiCollection::const_iterator digi = layerDigi.first;
00225                digi!=layerDigi.second;
00226                ++digi){
00227             if((cosmicRun && (*digi).countsTDC()<upperLimit) || (!cosmicRun))
00228               DigiPerWirePerEvent[(*digi).wire()]+=1;
00229           }
00230           // fill the digi event histo
00231           for (int wire=firstWire; wire<=lastWire; wire++) {
00232             theFile->cd();
00233             int histoEvents = nevents - (counter*1000);
00234             theHistoEvtPerWireMap[layerId]->Fill(histoEvents,wire,DigiPerWirePerEvent[wire]);
00235           }
00236         }
00237       } //Loop Ls
00238     } //Loop SLs
00239   } //Loop chambers
00240   
00241   
00242   if(nevents % 1000 == 0) {
00243     counter++;
00244     // save the digis event plot on file
00245     for(map<DTLayerId,  TH2F* >::const_iterator lHisto = theHistoEvtPerWireMap.begin();
00246         lHisto != theHistoEvtPerWireMap.end();
00247         lHisto++) {
00248       theFile->cd();
00249       if((*lHisto).second)
00250         (*lHisto).second->Write();
00251     }
00252     theHistoEvtPerWireMap.clear();
00253   }
00254   
00255 }
00256 
00257 
00258 void DTNoiseCalibration::endJob(){
00259 
00260   cout << "[DTNoiseCalibration] endjob called!" <<endl;
00261 
00262   // save the TDC digi plot
00263   theFile->cd();
00264   hTDCTriggerWidth->Write();
00265 
00266   // save on file the occupancy histo and write the list of noisy cells
00267   double TriggerWidth_s=0;
00268   DTStatusFlag *statusMap = new DTStatusFlag();
00269   for(map<DTLayerId, TH1F*>::const_iterator lHisto = theHistoOccupancyMap.begin();
00270       lHisto != theHistoOccupancyMap.end();
00271       lHisto++) {
00272     if(cosmicRun){
00273       if ( parameters.getUntrackedParameter<bool>("readDB", true) ) 
00274         tTrigMap->slTtrig( ((*lHisto).first).superlayerId(), tTrig, tTrigRMS); 
00275       else tTrig = parameters.getUntrackedParameter<int>("defaultTtrig", 4000);
00276       double TriggerWidth_ns = ((tTrig-500)*25)/32;
00277       TriggerWidth_s = TriggerWidth_ns/1e9;
00278     }
00279     if(!cosmicRun)
00280       TriggerWidth_s = double(TriggerWidth/1e9);
00281     if(debug)
00282       cout<<"TriggerWidth (s): "<<TriggerWidth_s<<"  TotEvents: "<<nevents<<endl;
00283     double normalization = 1/double(nevents*TriggerWidth_s);
00284     if((*lHisto).second){
00285       (*lHisto).second->Scale(normalization);
00286       theFile->cd();
00287       (*lHisto).second->Write();
00288       const DTTopology& dtTopo = dtGeom->layer((*lHisto).first)->specificTopology();
00289       const int firstWire = dtTopo.firstChannel();
00290       const int lastWire = dtTopo.lastChannel();
00291       for(int bin=firstWire; bin<=lastWire; bin++){
00292         //from definition of "noisy cell"
00293         if((*lHisto).second->GetBinContent(bin)>500){
00294           DTWireId wireID((*lHisto).first, bin);
00295           statusMap->setCellNoise(wireID,1);
00296         }
00297       }
00298     }
00299   }
00300   cout << "Writing Noise Map object to DB!" << endl;
00301   string record = "DTStatusFlagRcd";
00302   DTCalibDBUtils::writeToDB<DTStatusFlag>(record, statusMap);
00303 
00304  
00305   /*
00306   //save the digi event plot per SuperLayer 
00307   bool histo=false;
00308   map<DTSuperLayerId, vector<int> > maxPerSuperLayer;
00309   int numPlot = (nevents/1000);
00310   int num=0;
00311   // loop over the numPlot 
00312   for(int i=0; i<numPlot; i++){  
00313     vector<DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
00314     vector<DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
00315     // Loop over the chambers
00316     for (; ch_it != ch_end; ++ch_it) {
00317       DTChamberId ch = (*ch_it)->id();
00318       vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin(); 
00319       vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
00320       // Loop over the SLs
00321       for(; sl_it != sl_end; ++sl_it) {
00322         DTSuperLayerId sl = (*sl_it)->id();
00323         vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin(); 
00324         vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
00325         double dummy = pow(10.,10.);
00326         maxPerSuperLayer[sl].push_back(0);
00327         // Loop over the Ls
00328         for(; l_it != l_end; ++l_it) {
00329           DTLayerId layerId = (*l_it)->id();
00330 
00331           if (theHistoEvtPerWireMap.find(layerId) != theHistoEvtPerWireMap.end() &&
00332               theHistoEvtPerWireMap[layerId].size() > i){
00333             if (theHistoEvtPerWireMap[layerId][i]->GetMaximum(dummy)>maxPerSuperLayer[sl][i])
00334               maxPerSuperLayer[sl][i] = theHistoEvtPerWireMap[layerId][i]->GetMaximum(dummy);
00335           }
00336         }
00337       } // loop over SLs
00338     } // loop over chambers
00339   } // loop over numPlot
00340 
00341   // loop over the numPlot 
00342   for(int i=0; i<numPlot; i++){
00343     vector<DTChamber*>::const_iterator chamber_it = dtGeom->chambers().begin();
00344     vector<DTChamber*>::const_iterator chamber_end = dtGeom->chambers().end();
00345     // Loop over the chambers
00346     for (; chamber_it != chamber_end; ++chamber_it) {
00347       DTChamberId ch = (*chamber_it)->id();
00348       vector<const DTSuperLayer*>::const_iterator sl_it = (*chamber_it)->superLayers().begin(); 
00349       vector<const DTSuperLayer*>::const_iterator sl_end = (*chamber_it)->superLayers().end();
00350       // Loop over the SLs
00351       for(; sl_it != sl_end; ++sl_it) {
00352         DTSuperLayerId sl = (*sl_it)->id();
00353         vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin(); 
00354         vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
00355 
00356         stringstream num; num << i;
00357         string canvasName = "c" + getSuperLayerName(sl) + "_" + num.str();
00358         TCanvas c1(canvasName.c_str(),canvasName.c_str(),600,780);
00359         TLegend *leg=new TLegend(0.5,0.6,0.7,0.8);
00360         for(; l_it != l_end; ++l_it) {
00361           DTLayerId layerId = (*l_it)->id();
00362 
00363           if (theHistoEvtPerWireMap.find(layerId) != theHistoEvtPerWireMap.end() &&
00364             theHistoEvtPerWireMap[layerId].size() > i){
00365 
00366             string TitleName = "DigiPerWirePerEvent_" + getSuperLayerName(sl) + "_" + num.str();
00367             theHistoEvtPerWireMap[layerId][i]->SetTitle(TitleName.c_str());
00368             stringstream layer; layer << layerId.layer();       
00369             string legendHisto = "layer " + layer.str();
00370             leg->AddEntry(theHistoEvtPerWireMap[layerId][i],legendHisto.c_str(),"L");
00371             theHistoEvtPerWireMap[layerId][i]->SetMaximum(maxPerSuperLayer[sl][i]);
00372             if(histo==false)
00373               theHistoEvtPerWireMap[layerId][i]->Draw("lego");
00374             else
00375               theHistoEvtPerWireMap[layerId][i]->Draw("same , lego");
00376             theHistoEvtPerWireMap[layerId][i]->SetLineColor(layerId.layer());
00377             histo=true;
00378           }
00379         } // loop over Ls
00380         if(histo){
00381           leg->Draw("same");
00382           theFile->cd();
00383           c1.Write();
00384         }
00385         histo=false;
00386       } // loop over SLs
00387     } // loop over chambers 
00388   } // loop over numPlot
00389   */
00390 
00391 }
00392 
00393 
00394 
00395 DTNoiseCalibration::~DTNoiseCalibration(){
00396 
00397   cout << "DTNoiseCalibration: analyzed " << nevents << " events" <<endl;
00398   theFile->Close();
00399 
00400 }
00401 
00402 
00403 
00404 string DTNoiseCalibration::getLayerName(const DTLayerId& lId) const {
00405 
00406   const  DTSuperLayerId dtSLId = lId.superlayerId();
00407   const  DTChamberId dtChId = dtSLId.chamberId(); 
00408   stringstream Layer; Layer << lId.layer();
00409   stringstream superLayer; superLayer << dtSLId.superlayer();
00410   stringstream wheel; wheel << dtChId.wheel();  
00411   stringstream station; station << dtChId.station();    
00412   stringstream sector; sector << dtChId.sector();
00413   
00414   string LayerName = 
00415     "W" + wheel.str()
00416     + "_St" + station.str() 
00417     + "_Sec" + sector.str() 
00418     + "_SL" + superLayer.str()
00419     + "_L" + Layer.str();
00420   
00421   return LayerName;
00422 
00423 }
00424 
00425 
00426 string DTNoiseCalibration::getSuperLayerName(const DTSuperLayerId& dtSLId) const {
00427 
00428   const  DTChamberId dtChId = dtSLId.chamberId(); 
00429   stringstream superLayer; superLayer << dtSLId.superlayer();
00430   stringstream wheel; wheel << dtChId.wheel();  
00431   stringstream station; station << dtChId.station();    
00432   stringstream sector; sector << dtChId.sector();
00433   
00434   string SuperLayerName = 
00435     "W" + wheel.str()
00436     + "_St" + station.str() 
00437     + "_Sec" + sector.str() 
00438     + "_SL" + superLayer.str();
00439   
00440   return SuperLayerName;
00441 
00442 }
00443 
00444 
00445 
00446   

Generated on Tue Jun 9 17:25:28 2009 for CMSSW by  doxygen 1.5.4