CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/DQM/DTMonitorModule/src/DTLocalTriggerSynchTask.cc

Go to the documentation of this file.
00001 
00002 /*
00003  * \file DTLocalTriggerSynchTask.cc
00004  * 
00005  * $Date: 2010/04/19 13:59:00 $
00006  * $Revision: 1.7 $
00007  * \author C. Battilana - CIEMAT
00008  *
00009 */
00010 
00011 #include "DQM/DTMonitorModule/src/DTLocalTriggerSynchTask.h"
00012 
00013 // Framework
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Utilities/interface/InputTag.h"
00016 
00017 // DT trigger
00018 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
00019 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
00020 
00021 // Geometry
00022 #include "DataFormats/GeometryVector/interface/Pi.h"
00023 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00024 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00025 #include "Geometry/DTGeometry/interface/DTLayer.h"
00026 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
00027 #include "Geometry/DTGeometry/interface/DTTopology.h"
00028 
00029 // tTrigs
00030 #include "CalibMuon/DTDigiSync/interface/DTTTrigSyncFactory.h"
00031 #include "CalibMuon/DTDigiSync/interface/DTTTrigBaseSync.h"
00032 
00033 // DT Digi
00034 #include <DataFormats/DTDigi/interface/DTDigi.h>
00035 #include <DataFormats/DTDigi/interface/DTDigiCollection.h>
00036 
00037 
00038 //Root
00039 #include"TH1.h"
00040 #include"TAxis.h"
00041 
00042 #include <sstream>
00043 #include <iostream>
00044 #include <fstream>
00045 
00046 
00047 using namespace edm;
00048 using namespace std;
00049 
00050 DTLocalTriggerSynchTask::DTLocalTriggerSynchTask(const edm::ParameterSet& ps) : nevents(0) {
00051   
00052   edm::LogVerbatim ("DTLocalTriggerSynchTask")  << "[DTLocalTriggerSynchTask]: Constructor" << endl;
00053   parameters = ps;
00054 
00055 }
00056 
00057 
00058 DTLocalTriggerSynchTask::~DTLocalTriggerSynchTask() {
00059 
00060   edm::LogVerbatim ("DTLocalTriggerSynchTask")  << "[DTLocalTriggerSynchTask]: analyzed " << nevents << " events" << endl;
00061 
00062 }
00063 
00064 
00065 void DTLocalTriggerSynchTask::beginJob(){
00066 
00067   edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: BeginJob" << endl;
00068 
00069   bxTime        = parameters.getParameter<double>("bxTimeInterval");   // CB move this to static const or DB
00070   rangeInBX     = parameters.getParameter<bool>("rangeWithinBX");
00071   nBXLow        = parameters.getParameter<int>("nBXLow");
00072   nBXHigh       = parameters.getParameter<int>("nBXHigh");
00073   angleRange    = parameters.getParameter<double>("angleRange");
00074   minHitsPhi    = parameters.getParameter<int>("minHitsPhi");
00075   baseDirectory = parameters.getParameter<string>("baseDir");
00076 
00077   dbe = edm::Service<DQMStore>().operator->();
00078   dbe->setCurrentFolder(baseDir());
00079   dbe->bookFloat("BXTimeSpacing")->Fill(bxTime);
00080   
00081 }
00082 
00083 void DTLocalTriggerSynchTask::beginRun(const Run& run, const EventSetup& context) {
00084 
00085   edm::LogVerbatim ("DTLocalTriggerSynchTask") <<"[DTLocalTriggerSynchTask]: Begin Run"<<endl;
00086 
00087   context.get<MuonGeometryRecord>().get(muonGeom);
00088   tTrigSync = DTTTrigSyncFactory::get()->create(parameters.getParameter<std::string>("tTrigMode"),
00089                                                 parameters.getParameter<edm::ParameterSet>("tTrigModeConfig"));
00090   tTrigSync->setES(context);
00091 
00092 
00093   std::vector<DTChamber*>::const_iterator chambIt  = muonGeom->chambers().begin();
00094   std::vector<DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end();
00095 
00096   for (; chambIt!=chambEnd; ++chambIt) { 
00097     bookHistos((*chambIt)->id());
00098     triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL1"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(),1,1,2)));
00099     triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL3"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(),3,1,2)));
00100   } 
00101   
00102 }
00103 
00104 
00105 void DTLocalTriggerSynchTask::endJob() {
00106 
00107   edm::LogVerbatim ("DTLocalTriggerSynchTask")  << "[DTLocalTriggerSynchTask]: analyzed " << nevents << " events" << endl;
00108   dbe->rmdir(baseDir());
00109 
00110 }
00111 
00112 
00113 void DTLocalTriggerSynchTask::analyze(const edm::Event& event, const edm::EventSetup& context){
00114 
00115   nevents++;
00116 
00117   InputTag inputTagDCC  = parameters.getParameter<edm::InputTag>("DCCInputTag");
00118   InputTag inputTagDDU  = parameters.getParameter<edm::InputTag>("DDUInputTag");
00119   InputTag inputTagSEG  = parameters.getParameter<edm::InputTag>("SEGInputTag");
00120 
00121   for (int i=0;i<5;++i){
00122     for (int j=0;j<6;++j){
00123       for (int k=0;k<13;++k){
00124         phCodeBestDCC[j][i][k] = -1;
00125         phCodeBestDDU[j][i][k] = -1;
00126       }
00127     }
00128   }
00129 
00130   // Get best DCC triggers
00131   edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh;
00132   event.getByLabel(inputTagDCC,l1DTTPGPh);
00133   vector<L1MuDTChambPhDigi>*  phTrigs = l1DTTPGPh->getContainer();
00134   
00135   vector<L1MuDTChambPhDigi>::const_iterator iph  = phTrigs->begin();
00136   vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end();
00137   for(; iph !=iphe ; ++iph) {
00138     
00139     int phwheel = iph->whNum();
00140     int phsec   = iph->scNum() + 1; // DTTF[0-11] -> DT[1-12] Sector Numbering
00141     int phst    = iph->stNum();
00142     int phcode  = iph->code();
00143 
00144     if(phcode>phCodeBestDCC[phwheel+3][phst][phsec] && phcode<7) {
00145       phCodeBestDCC[phwheel+3][phst][phsec]=phcode; 
00146     }
00147     
00148   }
00149 
00150   // Get best DDU triggers
00151   Handle<DTLocalTriggerCollection> trigsDDU;
00152   event.getByLabel(inputTagDDU,trigsDDU);
00153   DTLocalTriggerCollection::DigiRangeIterator detUnitIt;
00154 
00155   for (detUnitIt=trigsDDU->begin();detUnitIt!=trigsDDU->end();++detUnitIt){
00156       
00157     const DTChamberId& id = (*detUnitIt).first;
00158     const DTLocalTriggerCollection::Range& range = (*detUnitIt).second;
00159 
00160     int wh = id.wheel();
00161     int sec = id.sector();
00162     int st = id.station();
00163 
00164     for (DTLocalTriggerCollection::const_iterator trigIt = range.first; trigIt!=range.second;++trigIt){
00165         
00166       int quality = trigIt->quality();
00167       
00168       if(quality>-1 && quality<7 &&
00169          quality>phCodeBestDDU[wh+3][st][sec]) {
00170         phCodeBestDDU[wh+3][st][sec]=quality;
00171       }
00172     }
00173   }
00174   
00175   //Get best segments (highest number of phi hits)
00176   vector<const DTRecSegment4D*> bestSegments4D;
00177   Handle<DTRecSegment4DCollection> segments4D;
00178   event.getByLabel(inputTagSEG, segments4D);
00179   DTRecSegment4DCollection::const_iterator track;
00180   DTRecSegment4DCollection::id_iterator chambIdIt;
00181 
00182   for (chambIdIt = segments4D->id_begin(); chambIdIt != segments4D->id_end(); ++chambIdIt){
00183     
00184     DTRecSegment4DCollection::range  range = segments4D->get(*chambIdIt);
00185     const DTRecSegment4D* best=0;
00186     int hitsBest = 0;
00187     int hits = 0;
00188 
00189     for ( track = range.first; track != range.second; ++track){
00190       if( (*track).hasPhi() ) {
00191         hits = (*track).phiSegment()->degreesOfFreedom()+2;
00192         if ( hits>hitsBest ){
00193           best = &(*track);
00194           hitsBest = hits;
00195         }
00196       }
00197     }
00198     if (best) {
00199       bestSegments4D.push_back(best);
00200     }
00201   }
00202     
00203   
00204   // Filling histos
00205   vector<const DTRecSegment4D*>::const_iterator bestSegIt  = bestSegments4D.begin();
00206   vector<const DTRecSegment4D*>::const_iterator bestSegEnd = bestSegments4D.end();
00207   for (; bestSegIt!=bestSegEnd; ++bestSegIt ){
00208 
00209     float dir = atan((*bestSegIt)->localDirection().x()/ (*bestSegIt)->localDirection().z())*180/Geom::pi(); // CB cerca un modo migliore x farlo
00210     const DTRecSegment2D* seg2D = (*bestSegIt)->phiSegment();
00211     int nHitsPhi = seg2D->degreesOfFreedom()+2; 
00212     DTChamberId chambId = (*bestSegIt)->chamberId();
00213     map<string, MonitorElement*> &innerME = triggerHistos[chambId.rawId()];
00214     
00215     if (fabs(dir)<angleRange && 
00216         nHitsPhi>=minHitsPhi && 
00217         seg2D->ist0Valid()){
00218       
00219       float t0seg = (*bestSegIt)->phiSegment()->t0();
00220       float tTrig = (tTrigSync->offset(DTWireId(chambId,1,1,2)) + tTrigSync->offset(DTWireId(chambId,3,1,2)) )/2;
00221       float time = tTrig+t0seg;
00222       float htime = rangeInBX ? time-int(time/bxTime)*bxTime : time-int(tTrig/bxTime)*bxTime;
00223       
00224       int wheel   = chambId.wheel();
00225       int sector  = chambId.sector();
00226       int station = chambId.station();
00227       int scsector = sector>12 ? sector==13 ? 4 : 10 : sector;
00228 
00229       int qualDCC = phCodeBestDCC[wheel+3][station][scsector];
00230       int qualDDU = phCodeBestDDU[wheel+3][station][scsector];
00231 
00232       if (fabs(t0seg)>0.01) {
00233         innerME.find("SEG_TrackCrossingTime")->second->Fill(htime);
00234         if ( qualDCC>=0 ) innerME.find("DCC_TrackCrossingTimeAll")->second->Fill(htime);          
00235         if ( qualDCC==6 ) innerME.find("DCC_TrackCrossingTimeHH")->second->Fill(htime);
00236         if ( qualDDU>=0 ) innerME.find("DDU_TrackCrossingTimeAll")->second->Fill(htime);          
00237         if ( qualDDU==6 ) innerME.find("DDU_TrackCrossingTimeHH")->second->Fill(htime);
00238       }
00239 
00240     }
00241   }
00242 
00243 }
00244 
00245 void DTLocalTriggerSynchTask::bookHistos(const DTChamberId& dtChId) {
00246   
00247   stringstream wheel; wheel << dtChId.wheel();  
00248   stringstream station; station << dtChId.station();    
00249   stringstream sector; sector << dtChId.sector();
00250   uint32_t chRawId = dtChId.rawId();
00251   
00252   dbe->setCurrentFolder(baseDir() + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() );
00253   
00254 
00255   string histoTag[5] = { "SEG_TrackCrossingTime", "DCC_TrackCrossingTimeAll", "DCC_TrackCrossingTimeHH", "DDU_TrackCrossingTimeAll", "DDU_TrackCrossingTimeHH" };
00256  
00257   float min = rangeInBX ?      0 : nBXLow*bxTime;
00258   float max = rangeInBX ? bxTime : nBXHigh*bxTime;
00259   int nbins = static_cast<int>(ceil( rangeInBX ? bxTime : (nBXHigh-nBXLow)*bxTime));
00260 
00261   for (int iHisto=0;iHisto<5;++iHisto) {
00262     string histoName = histoTag[iHisto] + (rangeInBX ? "InBX" : "") + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
00263     edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: booking " 
00264                                                  << baseDir() + "/Wheel" << wheel.str()
00265                                             << "/Sector" << sector.str()
00266                                             << "/Station"<< station.str() 
00267                                             << "/" << histoName << endl;
00268 
00269     triggerHistos[chRawId][histoTag[iHisto]] = dbe->book1D(histoName.c_str(),"Track time distribution",nbins,min,max);
00270   }
00271 
00272   string floatTag[2] = { "tTrig_SL1", "tTrig_SL3" };
00273 
00274   for (int iFloat=0;iFloat<2;++iFloat) { 
00275     string floatName = floatTag[iFloat] + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
00276     triggerHistos[chRawId][floatTag[iFloat]] = dbe->bookFloat(floatName);
00277   }
00278 
00279 }