CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch9/src/DQMOffline/CalibMuon/src/DTtTrigDBValidation.cc

Go to the documentation of this file.
00001 #include "DQMOffline/CalibMuon/interface/DTtTrigDBValidation.h"
00002 
00003 // Framework
00004 #include "FWCore/Framework/interface/Event.h"
00005 #include "FWCore/Framework/interface/EventSetup.h"
00006 #include "FWCore/Framework/interface/ESHandle.h"
00007 #include "FWCore/ServiceRegistry/interface/Service.h"
00008 
00009 #include "DQMServices/Core/interface/DQMStore.h"
00010 #include "DQMServices/Core/interface/MonitorElement.h"
00011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00012 
00013 // Geometry
00014 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00015 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00016 #include "Geometry/DTGeometry/interface/DTTopology.h"
00017 
00018 // DataFormats
00019 #include "DataFormats/MuonDetId/interface/DTSuperLayerId.h"
00020 
00021 // tTrig
00022 #include "CondFormats/DTObjects/interface/DTTtrig.h"
00023 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
00024 
00025 #include <stdio.h>
00026 #include <sstream>
00027 #include <math.h>
00028 #include "TFile.h"
00029 
00030 using namespace edm;
00031 using namespace std;
00032 
00033 DTtTrigDBValidation::DTtTrigDBValidation(const ParameterSet& pset):
00034    metname("tTrigdbValidation"),
00035    labelDBRef(pset.getParameter<string>("labelDBRef")),
00036    labelDB(pset.getParameter<string>("labelDB")),
00037    lowerLimit(pset.getUntrackedParameter<int>("lowerLimit",1)),
00038    higherLimit(pset.getUntrackedParameter<int>("higherLimit",3)),
00039    outputMEsInRootFile(pset.getUntrackedParameter<bool>("OutputMEsInRootFile",false)),
00040    outputFileName(pset.getUntrackedParameter<string>("OutputFileName","tTrigDBMonitoring.root"))
00041  {
00042 
00043   LogVerbatim(metname) << "[DTtTrigDBValidation] Constructor called!";
00044 
00045   // Get the DQM needed services
00046   dbe = edm::Service<DQMStore>().operator->();
00047   dbe->setCurrentFolder("DT/DTDBValidation");
00048 }
00049 
00050 DTtTrigDBValidation::~DTtTrigDBValidation(){}
00051 
00052 void DTtTrigDBValidation::beginRun(const edm::Run& run, const EventSetup& setup) {
00053 
00054   LogVerbatim(metname) << "[DTtTrigDBValidation] Parameters initialization";
00055  
00056   ESHandle<DTTtrig> tTrig_Ref;
00057   setup.get<DTTtrigRcd>().get(labelDBRef, tTrig_Ref);
00058   const DTTtrig* DTTtrigRefMap = &*tTrig_Ref;
00059   LogVerbatim(metname) << "[DTtTrigDBValidation] reference Ttrig version: " << tTrig_Ref->version();
00060 
00061   ESHandle<DTTtrig> tTrig;
00062   setup.get<DTTtrigRcd>().get(labelDB, tTrig);
00063   const DTTtrig* DTTtrigMap = &*tTrig;
00064   LogVerbatim(metname) << "[DTtTrigDBValidation] Ttrig to validate version: " << tTrig->version();
00065 
00066   //book&reset the summary histos
00067   for(int wheel=-2; wheel<=2; wheel++){
00068     bookHistos(wheel);
00069     tTrigDiffWheel[wheel]->Reset();
00070   }
00071 
00072   // Get the geometry
00073   setup.get<MuonGeometryRecord>().get(dtGeom);
00074 
00075   // Loop over Ref DB entries
00076   for(DTTtrig::const_iterator it = DTTtrigRefMap->begin();
00077       it != DTTtrigRefMap->end(); ++it) {
00078     DTSuperLayerId slId((*it).first.wheelId,
00079                         (*it).first.stationId,
00080                         (*it).first.sectorId,
00081                         (*it).first.slId);
00082     float tTrigMean;
00083     float tTrigRms;
00084     float kFactor;
00085     DTTtrigRefMap->get(slId, tTrigMean, tTrigRms, kFactor, DTTimeUnits::ns);
00086     float tTrigCorr = tTrigMean + kFactor*tTrigRms;
00087     LogTrace(metname)<< "Ref Superlayer: " <<  slId << "\n"
00088                      << " Ttrig mean (ns): " << tTrigMean
00089                      << " Ttrig rms (ns): " << tTrigRms
00090                      << " Ttrig k-Factor: " << kFactor
00091                      << " Ttrig value (ns): " << tTrigCorr;
00092                      
00093     //tTrigRefMap[slId] = std::pair<float,float>(tTrigmean,tTrigrms);
00094     tTrigRefMap[slId] = pair<float,float>(tTrigCorr,tTrigRms);
00095   }
00096 
00097   // Loop over Ref DB entries
00098   for(DTTtrig::const_iterator it = DTTtrigMap->begin();
00099       it != DTTtrigMap->end(); ++it) {
00100     DTSuperLayerId slId((*it).first.wheelId,
00101                         (*it).first.stationId,
00102                         (*it).first.sectorId,
00103                         (*it).first.slId);
00104     float tTrigMean;
00105     float tTrigRms;
00106     float kFactor;
00107     DTTtrigMap->get(slId, tTrigMean, tTrigRms, kFactor, DTTimeUnits::ns);
00108     float tTrigCorr = tTrigMean + kFactor*tTrigRms;
00109     LogTrace(metname)<< "Superlayer: " <<  slId << "\n"
00110                      << " Ttrig mean (ns): " << tTrigMean
00111                      << " Ttrig rms (ns): " << tTrigRms
00112                      << " Ttrig k-Factor: " << kFactor
00113                      << " Ttrig value (ns): " << tTrigCorr;
00114 
00115     //tTrigMap[slId] = std::pair<float,float>(tTrigmean,tTrigrms);
00116     tTrigMap[slId] = pair<float,float>(tTrigCorr,tTrigRms);
00117   }
00118 
00119   for(map<DTSuperLayerId, pair<float,float> >::const_iterator it = tTrigRefMap.begin();
00120       it != tTrigRefMap.end(); ++it) {  
00121       if(tTrigMap.find((*it).first) == tTrigMap.end()) continue;
00122 
00123       // compute the difference
00124       float difference = tTrigMap[(*it).first].first - (*it).second.first;
00125 
00126       //book histo
00127       int wheel = (*it).first.chamberId().wheel();
00128       int sector = (*it).first.chamberId().sector();    
00129       if(tTrigDiffHistos.find(make_pair(wheel,sector)) == tTrigDiffHistos.end()) bookHistos(wheel,sector);
00130                         
00131       LogTrace(metname) << "Filling histos for super-layer: " << (*it).first << "  difference: " << difference;
00132  
00133       // Fill the test histos
00134       int entry = -1;
00135       int station = (*it).first.chamberId().station();  
00136       if(station == 1) entry=0;
00137       if(station == 2) entry=3;
00138       if(station == 3) entry=6;
00139       if(station == 4) entry=9;
00140 
00141       int slBin = entry + (*it).first.superLayer();
00142       if(slBin == 12) slBin=11; 
00143 
00144       tTrigDiffHistos[make_pair(wheel,sector)]->setBinContent(slBin, difference);
00145       if(abs(difference)<lowerLimit){
00146         tTrigDiffWheel[wheel]->setBinContent(slBin,sector,1);
00147       }else if(abs(difference)<higherLimit){
00148         tTrigDiffWheel[wheel]->setBinContent(slBin,sector,2);
00149       }else{
00150         tTrigDiffWheel[wheel]->setBinContent(slBin,sector,3);
00151       }
00152         
00153 
00154   } // Loop over the tTrig map reference
00155   
00156 }
00157 
00158 void DTtTrigDBValidation::endRun(edm::Run const& run, edm::EventSetup const& setup) {}
00159 
00160 void DTtTrigDBValidation::endJob(){
00161 
00162   if(outputMEsInRootFile){
00163      // write the histos on a file
00164      dbe->save(outputFileName);
00165   }
00166 }
00167 
00168 void DTtTrigDBValidation::bookHistos(int wheel, int sector) {
00169 
00170   LogTrace(metname) << "   Booking histos for Wheel, Sector: " << wheel << ", " << sector;
00171 
00172   // Compose the chamber name
00173   stringstream str_wheel; str_wheel << wheel;
00174   stringstream str_sector; str_sector << sector;
00175 
00176   string lHistoName = "_W" + str_wheel.str() + "_Sec" + str_sector.str();
00177 
00178   dbe->setCurrentFolder("DT/tTrigValidation/Wheel" + str_wheel.str());
00179 
00180   // Create the monitor elements
00181   MonitorElement * hDifference;
00182   hDifference = dbe->book1D("htTrigDifference"+lHistoName, "difference between the two tTrig values",11,0,11);
00183 
00184   pair<int,int> mypair(wheel,sector);
00185   tTrigDiffHistos[mypair] = hDifference;
00186 
00187   (tTrigDiffHistos[mypair])->setBinLabel(1,"MB1_SL1",1);
00188   (tTrigDiffHistos[mypair])->setBinLabel(2,"MB1_SL2",1);
00189   (tTrigDiffHistos[mypair])->setBinLabel(3,"MB1_SL3",1);
00190   (tTrigDiffHistos[mypair])->setBinLabel(4,"MB2_SL1",1);
00191   (tTrigDiffHistos[mypair])->setBinLabel(5,"MB2_SL2",1);
00192   (tTrigDiffHistos[mypair])->setBinLabel(6,"MB2_SL3",1);
00193   (tTrigDiffHistos[mypair])->setBinLabel(7,"MB3_SL1",1);
00194   (tTrigDiffHistos[mypair])->setBinLabel(8,"MB3_SL2",1);
00195   (tTrigDiffHistos[mypair])->setBinLabel(9,"MB3_SL3",1);
00196   (tTrigDiffHistos[mypair])->setBinLabel(10,"MB4_SL1",1);
00197   (tTrigDiffHistos[mypair])->setBinLabel(11,"MB4_SL3",1);
00198 }
00199 
00200 // Book the summary histos
00201 void DTtTrigDBValidation::bookHistos(int wheel) {
00202 
00203   stringstream wh; wh << wheel;
00204 
00205   dbe->setCurrentFolder("DT/tTrigValidation/Wheel" + wh.str());
00206   tTrigDiffWheel[wheel] = dbe->book2D("htTrigDifference_W"+wh.str(), "W"+wh.str()+": summary of tTrig differences",11,1,12,14,1,15);
00207   tTrigDiffWheel[wheel]->setBinLabel(1,"MB1_SL1",1);
00208   tTrigDiffWheel[wheel]->setBinLabel(2,"MB1_SL2",1);
00209   tTrigDiffWheel[wheel]->setBinLabel(3,"MB1_SL3",1);
00210   tTrigDiffWheel[wheel]->setBinLabel(4,"MB2_SL1",1);
00211   tTrigDiffWheel[wheel]->setBinLabel(5,"MB2_SL2",1);
00212   tTrigDiffWheel[wheel]->setBinLabel(6,"MB2_SL3",1);
00213   tTrigDiffWheel[wheel]->setBinLabel(7,"MB3_SL1",1);
00214   tTrigDiffWheel[wheel]->setBinLabel(8,"MB3_SL2",1);
00215   tTrigDiffWheel[wheel]->setBinLabel(9,"MB3_SL3",1);
00216   tTrigDiffWheel[wheel]->setBinLabel(10,"MB4_SL1",1);
00217   tTrigDiffWheel[wheel]->setBinLabel(11,"MB4_SL3",1);
00218 
00219 
00220 }
00221 
00222 
00223 int DTtTrigDBValidation::stationFromBin(int bin) const {
00224   return (int) (bin /3.1)+1;
00225 }
00226  
00227 int DTtTrigDBValidation::slFromBin(int bin) const {
00228   int ret = bin%3;
00229   if(ret == 0 || bin == 11) ret = 3;
00230   
00231   return ret;
00232 }