CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/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  {
00040 
00041   LogVerbatim(metname_) << "[DTtTrigDBValidation] Constructor called!";
00042 
00043   // Get the DQM needed services
00044   dbe_ = edm::Service<DQMStore>().operator->();
00045   dbe_->setCurrentFolder("DT/DtCalib/TTrigDBValidation");
00046 
00047   outputMEsInRootFile_ = false;
00048   if( pset.exists("OutputFileName") ){
00049      outputMEsInRootFile_ = true;
00050      outputFileName_ = pset.getParameter<std::string>("OutputFileName");
00051   }
00052 }
00053 
00054 DTtTrigDBValidation::~DTtTrigDBValidation(){}
00055 
00056 void DTtTrigDBValidation::beginRun(const edm::Run& run, const EventSetup& setup) {
00057 
00058   LogVerbatim(metname_) << "[DTtTrigDBValidation] Parameters initialization";
00059  
00060   ESHandle<DTTtrig> tTrig_Ref;
00061   setup.get<DTTtrigRcd>().get(labelDBRef_, tTrig_Ref);
00062   const DTTtrig* DTTtrigRefMap = &*tTrig_Ref;
00063   LogVerbatim(metname_) << "[DTtTrigDBValidation] reference Ttrig version: " << tTrig_Ref->version();
00064 
00065   ESHandle<DTTtrig> tTrig;
00066   setup.get<DTTtrigRcd>().get(labelDB_, tTrig);
00067   const DTTtrig* DTTtrigMap = &*tTrig;
00068   LogVerbatim(metname_) << "[DTtTrigDBValidation] Ttrig to validate version: " << tTrig->version();
00069 
00070   //book&reset the summary histos
00071   for(int wheel=-2; wheel<=2; wheel++){
00072     bookHistos(wheel);
00073     tTrigDiffWheel_[wheel]->Reset();
00074   }
00075 
00076   // Get the geometry
00077   setup.get<MuonGeometryRecord>().get(dtGeom_);
00078 
00079   // Loop over Ref DB entries
00080   for(DTTtrig::const_iterator it = DTTtrigRefMap->begin();
00081       it != DTTtrigRefMap->end(); ++it) {
00082     DTSuperLayerId slId((*it).first.wheelId,
00083                         (*it).first.stationId,
00084                         (*it).first.sectorId,
00085                         (*it).first.slId);
00086     float tTrigMean;
00087     float tTrigRms;
00088     float kFactor;
00089     DTTtrigRefMap->get(slId, tTrigMean, tTrigRms, kFactor, DTTimeUnits::ns);
00090     float tTrigCorr = tTrigMean + kFactor*tTrigRms;
00091     LogTrace(metname_)<< "Ref Superlayer: " <<  slId << "\n"
00092                      << " Ttrig mean (ns): " << tTrigMean
00093                      << " Ttrig rms (ns): " << tTrigRms
00094                      << " Ttrig k-Factor: " << kFactor
00095                      << " Ttrig value (ns): " << tTrigCorr;
00096                      
00097     //tTrigRefMap[slId] = std::pair<float,float>(tTrigmean,tTrigrms);
00098     tTrigRefMap_[slId] = pair<float,float>(tTrigCorr,tTrigRms);
00099   }
00100 
00101   // Loop over Ref DB entries
00102   for(DTTtrig::const_iterator it = DTTtrigMap->begin();
00103       it != DTTtrigMap->end(); ++it) {
00104     DTSuperLayerId slId((*it).first.wheelId,
00105                         (*it).first.stationId,
00106                         (*it).first.sectorId,
00107                         (*it).first.slId);
00108     float tTrigMean;
00109     float tTrigRms;
00110     float kFactor;
00111     DTTtrigMap->get(slId, tTrigMean, tTrigRms, kFactor, DTTimeUnits::ns);
00112     float tTrigCorr = tTrigMean + kFactor*tTrigRms;
00113     LogTrace(metname_)<< "Superlayer: " <<  slId << "\n"
00114                      << " Ttrig mean (ns): " << tTrigMean
00115                      << " Ttrig rms (ns): " << tTrigRms
00116                      << " Ttrig k-Factor: " << kFactor
00117                      << " Ttrig value (ns): " << tTrigCorr;
00118 
00119     //tTrigMap[slId] = std::pair<float,float>(tTrigmean,tTrigrms);
00120     tTrigMap_[slId] = pair<float,float>(tTrigCorr,tTrigRms);
00121   }
00122 
00123   for(map<DTSuperLayerId, pair<float,float> >::const_iterator it = tTrigRefMap_.begin();
00124       it != tTrigRefMap_.end(); ++it) {  
00125       if(tTrigMap_.find((*it).first) == tTrigMap_.end()) continue;
00126 
00127       // compute the difference
00128       float difference = tTrigMap_[(*it).first].first - (*it).second.first;
00129 
00130       //book histo
00131       int wheel = (*it).first.chamberId().wheel();
00132       int sector = (*it).first.chamberId().sector();    
00133       if(tTrigDiffHistos_.find(make_pair(wheel,sector)) == tTrigDiffHistos_.end()) bookHistos(wheel,sector);
00134                         
00135       LogTrace(metname_) << "Filling histos for super-layer: " << (*it).first << "  difference: " << difference;
00136  
00137       // Fill the test histos
00138       int entry = -1;
00139       int station = (*it).first.chamberId().station();  
00140       if(station == 1) entry=0;
00141       if(station == 2) entry=3;
00142       if(station == 3) entry=6;
00143       if(station == 4) entry=9;
00144 
00145       int slBin = entry + (*it).first.superLayer();
00146       if(slBin == 12) slBin=11; 
00147 
00148       tTrigDiffHistos_[make_pair(wheel,sector)]->setBinContent(slBin, difference);
00149       if(abs(difference) < lowerLimit_){
00150         tTrigDiffWheel_[wheel]->setBinContent(slBin,sector,1);
00151       }else if(abs(difference) < higherLimit_){
00152         tTrigDiffWheel_[wheel]->setBinContent(slBin,sector,2);
00153       }else{
00154         tTrigDiffWheel_[wheel]->setBinContent(slBin,sector,3);
00155       }
00156         
00157 
00158   } // Loop over the tTrig map reference
00159   
00160 }
00161 
00162 void DTtTrigDBValidation::endRun(edm::Run const& run, edm::EventSetup const& setup) {}
00163 
00164 void DTtTrigDBValidation::endJob(){
00165 
00166   if(outputMEsInRootFile_){
00167      // write the histos on a file
00168      dbe_->save(outputFileName_);
00169   }
00170 }
00171 
00172 void DTtTrigDBValidation::bookHistos(int wheel, int sector) {
00173 
00174   LogTrace(metname_) << "   Booking histos for Wheel, Sector: " << wheel << ", " << sector;
00175 
00176   // Compose the chamber name
00177   stringstream str_wheel; str_wheel << wheel;
00178   stringstream str_sector; str_sector << sector;
00179 
00180   string lHistoName = "_W" + str_wheel.str() + "_Sec" + str_sector.str();
00181 
00182   dbe_->setCurrentFolder("DT/DtCalib/TTrigDBValidation/Wheel" + str_wheel.str());
00183 
00184   // Create the monitor elements
00185   MonitorElement * hDifference;
00186   hDifference = dbe_->book1D("TTrigDifference"+lHistoName, "difference between the two tTrig values",11,0,11);
00187 
00188   pair<int,int> mypair(wheel,sector);
00189   tTrigDiffHistos_[mypair] = hDifference;
00190 
00191   (tTrigDiffHistos_[mypair])->setBinLabel(1,"MB1_SL1",1);
00192   (tTrigDiffHistos_[mypair])->setBinLabel(2,"MB1_SL2",1);
00193   (tTrigDiffHistos_[mypair])->setBinLabel(3,"MB1_SL3",1);
00194   (tTrigDiffHistos_[mypair])->setBinLabel(4,"MB2_SL1",1);
00195   (tTrigDiffHistos_[mypair])->setBinLabel(5,"MB2_SL2",1);
00196   (tTrigDiffHistos_[mypair])->setBinLabel(6,"MB2_SL3",1);
00197   (tTrigDiffHistos_[mypair])->setBinLabel(7,"MB3_SL1",1);
00198   (tTrigDiffHistos_[mypair])->setBinLabel(8,"MB3_SL2",1);
00199   (tTrigDiffHistos_[mypair])->setBinLabel(9,"MB3_SL3",1);
00200   (tTrigDiffHistos_[mypair])->setBinLabel(10,"MB4_SL1",1);
00201   (tTrigDiffHistos_[mypair])->setBinLabel(11,"MB4_SL3",1);
00202 
00203 }
00204 
00205 // Book the summary histos
00206 void DTtTrigDBValidation::bookHistos(int wheel) {
00207 
00208   stringstream wh; wh << wheel;
00209 
00210   dbe_->setCurrentFolder("DT/DtCalib/TTrigDBValidation");
00211   tTrigDiffWheel_[wheel] = dbe_->book2D("TTrigDifference_W"+wh.str(), "W"+wh.str()+": summary of tTrig differences",11,1,12,14,1,15);
00212   tTrigDiffWheel_[wheel]->setBinLabel(1,"MB1_SL1",1);
00213   tTrigDiffWheel_[wheel]->setBinLabel(2,"MB1_SL2",1);
00214   tTrigDiffWheel_[wheel]->setBinLabel(3,"MB1_SL3",1);
00215   tTrigDiffWheel_[wheel]->setBinLabel(4,"MB2_SL1",1);
00216   tTrigDiffWheel_[wheel]->setBinLabel(5,"MB2_SL2",1);
00217   tTrigDiffWheel_[wheel]->setBinLabel(6,"MB2_SL3",1);
00218   tTrigDiffWheel_[wheel]->setBinLabel(7,"MB3_SL1",1);
00219   tTrigDiffWheel_[wheel]->setBinLabel(8,"MB3_SL2",1);
00220   tTrigDiffWheel_[wheel]->setBinLabel(9,"MB3_SL3",1);
00221   tTrigDiffWheel_[wheel]->setBinLabel(10,"MB4_SL1",1);
00222   tTrigDiffWheel_[wheel]->setBinLabel(11,"MB4_SL3",1);
00223 
00224 }
00225 
00226 
00227 int DTtTrigDBValidation::stationFromBin(int bin) const {
00228   return (int) (bin /3.1)+1;
00229 }
00230  
00231 int DTtTrigDBValidation::slFromBin(int bin) const {
00232   int ret = bin%3;
00233   if(ret == 0 || bin == 11) ret = 3;
00234   
00235   return ret;
00236 }