CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_9_patch3/src/DQM/DTMonitorModule/src/DTTriggerEfficiencyTask.cc

Go to the documentation of this file.
00001 /*
00002  * \file DTTriggerEfficiencyTask.cc
00003  * 
00004  * $Date: 2011/10/21 18:10:23 $
00005  * $Revision: 1.9 $
00006  * \author C.Battilana - CIEMAT
00007  *
00008  */
00009 
00010 #include "DQM/DTMonitorModule/src/DTTriggerEfficiencyTask.h"
00011 
00012 // Framework
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 
00015 // DT trigger
00016 #include "DQM/DTMonitorModule/interface/DTTrigGeomUtils.h"
00017 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
00018 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
00019 
00020 // Geometry
00021 #include "DataFormats/GeometryVector/interface/Pi.h"
00022 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00023 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00024 #include "Geometry/DTGeometry/interface/DTLayer.h"
00025 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
00026 #include "Geometry/DTGeometry/interface/DTTopology.h"
00027 
00028 // DT Digi
00029 #include <DataFormats/DTDigi/interface/DTDigi.h>
00030 #include <DataFormats/DTDigi/interface/DTDigiCollection.h>
00031 
00032 // Muon tracks
00033 #include <DataFormats/MuonReco/interface/Muon.h>
00034 #include <DataFormats/MuonReco/interface/MuonFwd.h>
00035 
00036 //Root
00037 #include"TH1.h"
00038 #include"TAxis.h"
00039 
00040 #include <sstream>
00041 #include <iostream>
00042 #include <fstream>
00043 
00044 
00045 using namespace edm;
00046 using namespace std;
00047 
00048 DTTriggerEfficiencyTask::DTTriggerEfficiencyTask(const edm::ParameterSet& ps) : trigGeomUtils(0) {
00049 
00050   LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask")  << "[DTTriggerEfficiencyTask]: Constructor" << endl;
00051 
00052   parameters = ps;
00053   dbe = edm::Service<DQMStore>().operator->();
00054 
00055 }
00056 
00057 
00058 DTTriggerEfficiencyTask::~DTTriggerEfficiencyTask() {
00059 
00060   LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask")  << "[DTTriggerEfficiencyTask]: analyzed " << nevents << " events" << endl;
00061 
00062 }
00063 
00064 
00065 void DTTriggerEfficiencyTask::beginJob(){
00066 
00067   LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") << "[DTTriggerEfficiencyTask]: BeginJob" << endl;
00068 
00069   inputTagMuons = parameters.getUntrackedParameter<edm::InputTag>("inputTagMuons");
00070 
00071   inputTagDCC = parameters.getUntrackedParameter<edm::InputTag>("inputTagDCC");
00072   inputTagDDU = parameters.getUntrackedParameter<edm::InputTag>("inputTagDDU");
00073   inputTagSEG = parameters.getUntrackedParameter<edm::InputTag>("inputTagSEG");
00074   inputTagGMT = parameters.getUntrackedParameter<edm::InputTag>("inputTagGMT");
00075 
00076   SegmArbitration = parameters.getUntrackedParameter<std::string>("SegmArbitration");
00077 
00078   detailedPlots = parameters.getUntrackedParameter<bool>("detailedAnalysis");
00079   processDCC = parameters.getUntrackedParameter<bool>("processDCC");
00080   processDDU = parameters.getUntrackedParameter<bool>("processDDU");
00081   minBXDDU = parameters.getUntrackedParameter<int>("minBXDDU");
00082   maxBXDDU = parameters.getUntrackedParameter<int>("maxBXDDU");
00083 
00084   nMinHitsPhi = parameters.getUntrackedParameter<int>("nMinHitsPhi");
00085   phiAccRange = parameters.getUntrackedParameter<double>("phiAccRange");
00086 
00087   if (processDCC) processTags.push_back("DCC");
00088   if (processDDU) processTags.push_back("DDU");
00089 
00090 }
00091 
00092 void DTTriggerEfficiencyTask::beginRun(const edm::Run& run, const edm::EventSetup& context){
00093 
00094   LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") << "[DTTriggerEfficiencyTask]: BeginRun" << endl;
00095 
00096   nevents = 0;
00097 
00098   for (int wh=-2;wh<=2;++wh){
00099     vector<string>::const_iterator tagIt  = processTags.begin();
00100     vector<string>::const_iterator tagEnd = processTags.end();
00101     for (; tagIt!=tagEnd; ++tagIt) {
00102       bookWheelHistos(wh,(*tagIt),"Task");
00103       if (detailedPlots) {
00104         for (int stat=1;stat<=4;++stat){
00105           for (int sect=1;sect<=12;++sect){
00106             bookChamberHistos(DTChamberId(wh,stat,sect),(*tagIt),"Segment");
00107           }
00108         }
00109       }
00110     }
00111   }
00112 
00113   context.get<MuonGeometryRecord>().get(muonGeom);
00114   trigGeomUtils = new DTTrigGeomUtils(muonGeom);
00115 
00116 }
00117 
00118 void DTTriggerEfficiencyTask::beginLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
00119 
00120   LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") <<"[DTTriggerEfficiencyTask]: Begin of LS transition"<<endl;
00121 
00122 }
00123 
00124 
00125 void DTTriggerEfficiencyTask::endJob(){
00126 
00127   LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask")  << "[DTTriggerEfficiencyTask]: analyzed " << nevents << " events" << endl;
00128 
00129 }
00130 
00131 
00132 void DTTriggerEfficiencyTask::analyze(const edm::Event& e, const edm::EventSetup& c){
00133 
00134   nevents++;
00135 
00136   if (!hasRPCTriggers(e)) { return; }
00137 
00138   map<DTChamberId,const L1MuDTChambPhDigi*> phBestDCC;
00139   map<DTChamberId,const DTLocalTrigger*>    phBestDDU;
00140 
00141   // Getting best DCC Stuff
00142   edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh;
00143   e.getByLabel(inputTagDCC,l1DTTPGPh);
00144   vector<L1MuDTChambPhDigi>*  phTrigs = l1DTTPGPh->getContainer();
00145 
00146   vector<L1MuDTChambPhDigi>::const_iterator iph  = phTrigs->begin();
00147   vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end();
00148   for(; iph !=iphe ; ++iph) {
00149 
00150     int phwheel = iph->whNum();
00151     int phsec   = iph->scNum() + 1; // DTTF numbering [0:11] -> DT numbering [1:12]
00152     int phst    = iph->stNum();
00153     int phcode  = iph->code();
00154 
00155     DTChamberId chId(phwheel,phst,phsec);
00156 
00157     if( phcode < 7 && (phBestDCC.find(chId) == phBestDCC.end() || 
00158           phcode>phBestDCC[chId]->code()) ) phBestDCC[chId] = &(*iph);
00159   }
00160 
00161   //Getting Best DDU Stuff
00162   Handle<DTLocalTriggerCollection> trigsDDU;
00163   e.getByLabel(inputTagDDU,trigsDDU);
00164   DTLocalTriggerCollection::DigiRangeIterator detUnitIt;
00165 
00166   for (detUnitIt=trigsDDU->begin();detUnitIt!=trigsDDU->end();++detUnitIt){
00167 
00168     const DTChamberId& id = (*detUnitIt).first;
00169     const DTLocalTriggerCollection::Range& range = (*detUnitIt).second;
00170 
00171     DTLocalTriggerCollection::const_iterator trigIt  = range.first;
00172     DTLocalTriggerCollection::const_iterator trigEnd = range.second;
00173     for (; trigIt!= trigEnd;++trigIt){
00174       int quality = trigIt->quality();
00175       if(quality>-1 && quality<7 &&
00176           (phBestDDU.find(id) == phBestDDU.end() || 
00177            quality>phBestDDU[id]->quality()) ) phBestDDU[id] = &(*trigIt);
00178     }
00179 
00180   }
00181 
00182   //Getting Best Segments
00183   vector<const DTRecSegment4D*> best4DSegments;
00184 
00185   Handle<reco::MuonCollection> muons;
00186   e.getByLabel(inputTagMuons, muons);
00187   reco::MuonCollection::const_iterator mu;
00188 
00189   for( mu = muons->begin(); mu != muons->end(); ++mu ) {
00190 
00191     // Make sure that is standalone muon
00192     if( !((*mu).isStandAloneMuon()) ) {continue;}
00193 
00194     // Get the chambers compatible with the muon
00195     const vector<reco::MuonChamberMatch> matchedChambers = (*mu).matches(); 
00196     vector<reco::MuonChamberMatch>::const_iterator chamber;
00197 
00198     for( chamber = matchedChambers.begin(); chamber != matchedChambers.end(); ++chamber ) {
00199 
00200       // look only in DTs
00201       if( chamber->detector() != MuonSubdetId::DT ) {continue;}
00202 
00203       // Get the matched segments in the chamber
00204       const vector<reco::MuonSegmentMatch> matchedSegments = (*chamber).segmentMatches; 
00205       vector<reco::MuonSegmentMatch>::const_iterator segment;
00206 
00207       for( segment = matchedSegments.begin(); segment != matchedSegments.end(); ++segment ) {
00208 
00209         edm::Ref<DTRecSegment4DCollection> dtSegment = segment->dtSegmentRef;
00210 
00211         // Segment Arbitration
00212         if( SegmArbitration == "SegmentArbitration" 
00213             && !((*segment).isMask(reco::MuonSegmentMatch::BestInChamberByDR)) ) {continue;}
00214 
00215         if( SegmArbitration == "SegmentAndTrackArbitration" 
00216             && (!((*segment).isMask(reco::MuonSegmentMatch::BestInChamberByDR)) ||
00217               !((*segment).isMask(reco::MuonSegmentMatch::BelongsToTrackByDR))) ) {continue;} 
00218 
00219         if( SegmArbitration == "SegmentAndTrackArbitrationCleaned" 
00220             && (!((*segment).isMask(reco::MuonSegmentMatch::BestInChamberByDR))  ||
00221               !((*segment).isMask(reco::MuonSegmentMatch::BelongsToTrackByDR)) ||
00222               !((*segment).isMask(reco::MuonSegmentMatch::BelongsToTrackByCleaning))) ) {continue;} 
00223 
00224 
00225         if( (*dtSegment).hasPhi() ) {
00226           best4DSegments.push_back(&(*dtSegment));
00227         }
00228 
00229       }// end loop on matched segments
00230     }// end loop on compatible chambers
00231   }// end loop on muons
00232 
00233   // Plot filling
00234   vector<const DTRecSegment4D*>::const_iterator btrack;
00235   for ( btrack = best4DSegments.begin(); btrack != best4DSegments.end(); ++btrack ){
00236 
00237     int wheel    = (*btrack)->chamberId().wheel();
00238     int station  = (*btrack)->chamberId().station();
00239     int scsector = 0;
00240     float x, xdir, y, ydir;
00241     trigGeomUtils->computeSCCoordinates((*btrack),scsector,x,xdir,y,ydir);
00242     int nHitsPhi = (*btrack)->phiSegment()->degreesOfFreedom()+2;       
00243     DTChamberId dtChId(wheel,station,scsector);
00244     uint32_t indexCh = dtChId.rawId(); 
00245     map<string, MonitorElement*> &innerChME = chamberHistos[indexCh];
00246     map<string, MonitorElement*> &innerWhME = wheelHistos[wheel];
00247 
00248     if (fabs(xdir)<phiAccRange && nHitsPhi>=nMinHitsPhi){
00249 
00250       vector<string>::const_iterator tagIt  = processTags.begin();
00251       vector<string>::const_iterator tagEnd = processTags.end();
00252 
00253       for (; tagIt!=tagEnd; ++tagIt) { 
00254 
00255         int qual   = (*tagIt) == "DCC" ?
00256           phBestDCC.find(dtChId) != phBestDCC.end() ? phBestDCC[dtChId]->code() : -1 :
00257           phBestDDU.find(dtChId) != phBestDDU.end() ? phBestDDU[dtChId]->quality() : -1;
00258 
00259         innerWhME.find((*tagIt) + "_TrigEffDenum")->second->Fill(scsector,station);
00260         if ( qual>=0 && qual<7 ) {
00261           innerWhME.find((*tagIt) + "_TrigEffNum")->second->Fill(scsector,station);
00262           if ( qual>=4 ) {
00263             innerWhME.find((*tagIt) + "_TrigEffCorrNum")->second->Fill(scsector,station);
00264           }
00265         }
00266         if (detailedPlots) {
00267           innerChME.find((*tagIt) + "_TrackPosvsAngle")->second->Fill(xdir,x);
00268           if ( qual>=0 && qual<7 ) {
00269             innerChME.find((*tagIt) + "_TrackPosvsAngleAnyQual")->second->Fill(xdir,x);
00270             if ( qual>=4 ) {
00271               innerChME.find((*tagIt) + "_TrackPosvsAngleCorr")->second->Fill(xdir,x);  
00272             }
00273           }
00274         }
00275       }
00276     }
00277   }
00278 
00279 }
00280 
00281 bool DTTriggerEfficiencyTask::hasRPCTriggers(const edm::Event& e) {
00282 
00283   edm::Handle<L1MuGMTReadoutCollection> gmtrc; 
00284   e.getByLabel(inputTagGMT,gmtrc);
00285 
00286   std::vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords();
00287   std::vector<L1MuGMTReadoutRecord>::const_iterator igmtrr = gmt_records.begin();
00288   std::vector<L1MuGMTReadoutRecord>::const_iterator egmtrr = gmt_records.end();
00289   for(; igmtrr!=egmtrr; igmtrr++) {
00290 
00291     std::vector<L1MuGMTExtendedCand> candsGMT = igmtrr->getGMTCands();
00292     std::vector<L1MuGMTExtendedCand>::const_iterator candGMTIt   = candsGMT.begin();
00293     std::vector<L1MuGMTExtendedCand>::const_iterator candGMTEnd  = candsGMT.end();
00294 
00295     for(; candGMTIt!=candGMTEnd; ++candGMTIt){
00296       if(!candGMTIt->empty()) {
00297         int quality = candGMTIt->quality();
00298         if(candGMTIt->bx()==0 && 
00299             (quality == 5 || quality == 7)){
00300           return true;
00301         }
00302       }
00303     }
00304   }
00305 
00306   return false;
00307 
00308 }
00309 
00310 void DTTriggerEfficiencyTask::bookChamberHistos(const DTChamberId& dtCh, string histoType, string folder) {
00311 
00312   int wh = dtCh.wheel();                
00313   int sc = dtCh.sector();       
00314   int st = dtCh.station();
00315   stringstream wheel; wheel << wh;      
00316   stringstream station; station << st;  
00317   stringstream sector; sector << sc;    
00318 
00319   string hwFolder      = topFolder(histoType); 
00320   string bookingFolder = hwFolder + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() + "/" + folder;
00321   string histoTag      = "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
00322 
00323   dbe->setCurrentFolder(bookingFolder);
00324 
00325   LogTrace ("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") 
00326     << "[DTTriggerEfficiencyTask]: booking histos in " << bookingFolder << endl;
00327 
00328   float min, max;
00329   int nbins;
00330   trigGeomUtils->phiRange(dtCh,min,max,nbins,20);
00331 
00332   string histoName = histoType + "_TrackPosvsAngle" +  histoTag;
00333   string histoLabel = "Position vs Angle (phi)";
00334   (chamberHistos[dtCh.rawId()])[histoType + "_TrackPosvsAngle"] = 
00335     dbe->book2D(histoName,histoLabel,12,-30.,30.,nbins,min,max);
00336 
00337   histoName = histoType + "_TrackPosvsAngleAnyQual" +  histoTag;
00338   histoLabel = "Position vs Angle (phi) for any qual triggers";
00339   (chamberHistos[dtCh.rawId()])[histoType + "_TrackPosvsAngleAnyQual"] = 
00340     dbe->book2D(histoName,histoLabel,12,-30.,30.,nbins,min,max);
00341 
00342   histoName = histoType + "_TrackPosvsAngleCorr" +  histoTag;
00343   histoLabel = "Position vs Angle (phi) for correlated triggers";
00344   (chamberHistos[dtCh.rawId()])[histoType + "_TrackPosvsAngleCorr"] = 
00345     dbe->book2D(histoName,histoLabel,12,-30.,30.,nbins,min,max);
00346 
00347 }
00348 
00349 void DTTriggerEfficiencyTask::bookWheelHistos(int wheel,string hTag,string folder) {  
00350 
00351   stringstream wh; wh << wheel;
00352   string basedir;  
00353   if (hTag.find("Summary") != string::npos ) {
00354     basedir = topFolder(hTag);   //Book summary histo outside folder directory
00355   } else {
00356     basedir = topFolder(hTag) + folder + "/" ;
00357   }
00358 
00359   dbe->setCurrentFolder(basedir);
00360 
00361   string hTagName = "_W" + wh.str();
00362 
00363   LogTrace("DTDQM|DTMonitorModule|DTTriggerEfficiencyTask") 
00364     << "[DTTriggerEfficiencyTask]: booking histos in "<< basedir << endl;  
00365 
00366   string hName = hTag + "_TrigEffDenum" + hTagName;
00367   MonitorElement* me = dbe->book2D(hName.c_str(),hName.c_str(),12,1,13,4,1,5);
00368 
00369   me->setBinLabel(1,"MB1",2);
00370   me->setBinLabel(2,"MB2",2);
00371   me->setBinLabel(3,"MB3",2);
00372   me->setBinLabel(4,"MB4",2);
00373   me->setAxisTitle("Sector",1);
00374 
00375   wheelHistos[wheel][hTag + "_TrigEffDenum"] = me;
00376 
00377   hName = hTag + "_TrigEffNum" + hTagName;
00378   me = dbe->book2D(hName.c_str(),hName.c_str(),12,1,13,4,1,5);
00379 
00380   me->setBinLabel(1,"MB1",2);
00381   me->setBinLabel(2,"MB2",2);
00382   me->setBinLabel(3,"MB3",2);
00383   me->setBinLabel(4,"MB4",2);
00384   me->setAxisTitle("Sector",1);
00385 
00386   wheelHistos[wheel][hTag + "_TrigEffNum"] = me;
00387 
00388   hName = hTag + "_TrigEffCorrNum" + hTagName;
00389   me = dbe->book2D(hName.c_str(),hName.c_str(),12,1,13,4,1,5);
00390 
00391   me->setBinLabel(1,"MB1",2);
00392   me->setBinLabel(2,"MB2",2);
00393   me->setBinLabel(3,"MB3",2);
00394   me->setBinLabel(4,"MB4",2);
00395   me->setAxisTitle("Sector",1);
00396 
00397   wheelHistos[wheel][hTag + "_TrigEffCorrNum"] = me;
00398 
00399   return;
00400 
00401 }
00402