CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/DQM/DTMonitorModule/src/DTLocalTriggerTask.cc

Go to the documentation of this file.
00001 /*
00002  * \file DTLocalTriggerTask.cc
00003  * 
00004  * $Date: 2012/09/24 16:08:07 $
00005  * $Revision: 1.40 $
00006  * \author M. Zanetti - INFN Padova
00007  *
00008 */
00009 
00010 #include "DQM/DTMonitorModule/src/DTLocalTriggerTask.h"
00011 
00012 // Framework
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 
00015 // DT trigger
00016 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
00017 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
00018 #include "DataFormats/LTCDigi/interface/LTCDigi.h"
00019 #include "DQM/DTMonitorModule/interface/DTTrigGeomUtils.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/DTTopology.h"
00027 
00028 //Root
00029 #include"TH1.h"
00030 #include"TAxis.h"
00031 
00032 #include <sstream>
00033 #include <iostream>
00034 #include <fstream>
00035 
00036 
00037 using namespace edm;
00038 using namespace std;
00039 
00040 DTLocalTriggerTask::DTLocalTriggerTask(const edm::ParameterSet& ps) :
00041   trigGeomUtils(0),
00042   isLocalRun(ps.getUntrackedParameter<bool>("localrun", true))
00043  {
00044   if (!isLocalRun) {
00045     ltcDigiCollectionTag = ps.getParameter<edm::InputTag>("ltcDigiCollectionTag");
00046   }
00047   
00048   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerTask") << "[DTLocalTriggerTask]: Constructor"<<endl;
00049 
00050   tpMode           = ps.getUntrackedParameter<bool>("testPulseMode", false);
00051   detailedAnalysis = ps.getUntrackedParameter<bool>("detailedAnalysis", false);
00052   doDCCTheta       = ps.getUntrackedParameter<bool>("enableDCCTheta", false);
00053   
00054   if (tpMode) {
00055     baseFolderDCC = "DT/11-LocalTriggerTP-DCC/";
00056     baseFolderDDU = "DT/12-LocalTriggerTP-DDU/";
00057   }
00058   else {
00059     baseFolderDCC = "DT/03-LocalTrigger-DCC/";
00060     baseFolderDDU = "DT/04-LocalTrigger-DDU/";
00061   }
00062 
00063   parameters = ps;
00064   
00065   dbe = edm::Service<DQMStore>().operator->();
00066 
00067 }
00068 
00069 
00070 DTLocalTriggerTask::~DTLocalTriggerTask() {
00071 
00072   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerTask") << "[DTLocalTriggerTask]: analyzed " << nevents << " events" << endl;
00073   if (trigGeomUtils) { delete trigGeomUtils; }
00074 
00075 }
00076 
00077 
00078 void DTLocalTriggerTask::beginJob(){
00079  
00080   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerTask") << "[DTLocalTriggerTask]: BeginJob" << endl;
00081 
00082   nevents = 0;
00083 
00084 }
00085 
00086 void DTLocalTriggerTask::beginRun(const edm::Run& run, const edm::EventSetup& context) {
00087 
00088   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerTask") << "[DTLocalTriggerTask]: BeginRun" << endl;   
00089 
00090   context.get<MuonGeometryRecord>().get(muonGeom);
00091   trigGeomUtils = new DTTrigGeomUtils(muonGeom);
00092 
00093   if(parameters.getUntrackedParameter<bool>("staticBooking", true)) {  // Static histo booking
00094    
00095     vector<string> trigSources;
00096     if(parameters.getUntrackedParameter<bool>("localrun", true)) {
00097       trigSources.push_back("");
00098     }
00099     else {
00100       trigSources.push_back("_DTonly");
00101       trigSources.push_back("_NoDT");
00102       trigSources.push_back("_DTalso");
00103     }
00104     vector<string>::const_iterator trigSrcIt  = trigSources.begin();
00105     vector<string>::const_iterator trigSrcEnd = trigSources.end();
00106    
00107     if(parameters.getUntrackedParameter<bool>("process_dcc", true)) {
00108       bookBarrelHistos("DCC_ErrorsChamberID");
00109     }
00110 
00111     if (tpMode) {
00112       for (int stat=1;stat<5;++stat){
00113         for (int wh=-2;wh<3;++wh){
00114           for (int sect=1;sect<13;++sect){
00115             DTChamberId dtChId(wh,stat,sect);
00116             
00117             if (parameters.getUntrackedParameter<bool>("process_dcc", true)){ // DCC data
00118               bookHistos(dtChId,"LocalTriggerPhi","DCC_BXvsQual"+(*trigSrcIt));
00119               bookHistos(dtChId,"LocalTriggerPhi","DCC_QualvsPhirad"+(*trigSrcIt));
00120             }
00121             
00122             if (parameters.getUntrackedParameter<bool>("process_ros", true)){ // DDU data             
00123               bookHistos(dtChId,"LocalTriggerPhi","DDU_BXvsQual"+(*trigSrcIt));
00124             }
00125             
00126           }
00127         }
00128       } // end of loop
00129     }
00130     else {
00131       for (;trigSrcIt!=trigSrcEnd;++trigSrcIt){
00132         for (int wh=-2;wh<3;++wh){
00133           if (parameters.getUntrackedParameter<bool>("process_dcc", true) &&
00134               parameters.getUntrackedParameter<bool>("process_ros", true)){ // DCC+DDU data
00135             bookWheelHistos(wh,"COM_BXDiff"+(*trigSrcIt));
00136           }
00137           for (int sect=1;sect<13;++sect){
00138             for (int stat=1;stat<5;++stat){
00139               DTChamberId dtChId(wh,stat,sect);
00140               if (parameters.getUntrackedParameter<bool>("process_dcc", true)){ // DCC data
00141               
00142                 bookHistos(dtChId,"LocalTriggerPhi","DCC_BXvsQual"+(*trigSrcIt));
00143                 if (detailedAnalysis) {
00144                   bookHistos(dtChId,"LocalTriggerPhi","DCC_QualvsPhirad"+(*trigSrcIt));
00145                   bookHistos(dtChId,"LocalTriggerPhi","DCC_QualvsPhibend"+(*trigSrcIt));
00146                 }
00147                 bookHistos(dtChId,"LocalTriggerPhi","DCC_Flag1stvsQual"+(*trigSrcIt));
00148                 bookHistos(dtChId,"LocalTriggerPhi","DCC_BestQual"+(*trigSrcIt));
00149                 if (stat!=4 && doDCCTheta){                                                  
00150                   bookHistos(dtChId,"LocalTriggerTheta","DCC_PositionvsBX"+(*trigSrcIt));
00151 //                bookHistos(dtChId,"LocalTriggerTheta","DCC_PositionvsQual"+(*trigSrcIt));  // DCC theta quality not available!        
00152 //                bookHistos(dtChId,"LocalTriggerTheta","DCC_ThetaBXvsQual"+(*trigSrcIt));
00153 //                bookHistos(dtChId,"LocalTriggerTheta","DCC_ThetaBestQual"+(*trigSrcIt));
00154                 }
00155               
00156                 if (parameters.getUntrackedParameter<bool>("process_seg", true)){ // DCC + Segemnt
00157                   bookHistos(dtChId,"Segment","DCC_PhitkvsPhitrig"+(*trigSrcIt));
00158                   bookHistos(dtChId,"Segment","DCC_PhibtkvsPhibtrig"+(*trigSrcIt));
00159                   bookHistos(dtChId,"Segment","DCC_PhiResidual"+(*trigSrcIt));
00160                   bookHistos(dtChId,"Segment","DCC_PhiResidualvsLUTPhi"+(*trigSrcIt));
00161                   bookHistos(dtChId,"Segment","DCC_PhibResidual"+(*trigSrcIt));
00162                   bookHistos(dtChId,"Segment","DCC_HitstkvsQualtrig"+(*trigSrcIt));
00163                   bookHistos(dtChId,"Segment","DCC_TrackPosvsAngle"+(*trigSrcIt));
00164                   bookHistos(dtChId,"Segment","DCC_TrackPosvsAngleandTrig"+(*trigSrcIt));
00165                   bookHistos(dtChId,"Segment","DCC_TrackPosvsAngleandTrigHHHL"+(*trigSrcIt));
00166                   if(stat!=4){
00167                     bookHistos(dtChId,"Segment","DCC_TrackThetaPosvsAngle"+(*trigSrcIt)); // theta view
00168                     bookHistos(dtChId,"Segment","DCC_TrackThetaPosvsAngleandTrig"+(*trigSrcIt));
00169 //                  bookHistos(dtChId,"Segment","DCC_TrackThetaPosvsAngleandTrigH"+(*trigSrcIt));     // DCC theta quality not available!
00170                   }
00171                 }
00172               
00173               }
00174             
00175               if (parameters.getUntrackedParameter<bool>("process_ros", true)){ // DDU data
00176               
00177                 bookHistos(dtChId,"LocalTriggerPhi","DDU_BXvsQual"+(*trigSrcIt));
00178                 bookHistos(dtChId,"LocalTriggerPhi","DDU_Flag1stvsQual"+(*trigSrcIt));
00179                 bookHistos(dtChId,"LocalTriggerPhi","DDU_BestQual"+(*trigSrcIt));
00180                 if(stat!=4){                                                    // theta view
00181                   bookHistos(dtChId,"LocalTriggerTheta","DDU_ThetaBXvsQual"+(*trigSrcIt));
00182                   bookHistos(dtChId,"LocalTriggerTheta","DDU_ThetaBestQual"+(*trigSrcIt));
00183                 }
00184               
00185                 if (parameters.getUntrackedParameter<bool>("process_seg", true)){ // DDU + Segment
00186                   bookHistos(dtChId,"Segment","DDU_HitstkvsQualtrig"+(*trigSrcIt));
00187                   bookHistos(dtChId,"Segment","DDU_TrackPosvsAngle"+(*trigSrcIt));
00188                   bookHistos(dtChId,"Segment","DDU_TrackPosvsAngleandTrig"+(*trigSrcIt));
00189                   bookHistos(dtChId,"Segment","DDU_TrackPosvsAngleandTrigHHHL"+(*trigSrcIt));
00190                   if(stat!=4){
00191                     bookHistos(dtChId,"Segment","DDU_TrackThetaPosvsAngle"+(*trigSrcIt)); // theta view
00192                     bookHistos(dtChId,"Segment","DDU_TrackThetaPosvsAngleandTrig"+(*trigSrcIt));
00193                     bookHistos(dtChId,"Segment","DDU_TrackThetaPosvsAngleandTrigH"+(*trigSrcIt));
00194                   }
00195                 }
00196               
00197               }
00198             
00199               if (parameters.getUntrackedParameter<bool>("process_dcc", true) &&
00200                   parameters.getUntrackedParameter<bool>("process_ros", true)){ // DCC+DDU data
00201                 bookHistos(dtChId,"LocalTriggerPhi","COM_QualDDUvsQualDCC"+(*trigSrcIt));
00202               }
00203             
00204             }
00205           }
00206           for (int sect=13;sect<15;++sect){
00207             DTChamberId dtChId(wh,4,sect);
00208             if (parameters.getUntrackedParameter<bool>("process_dcc", true) &&
00209                 parameters.getUntrackedParameter<bool>("process_seg", true)){ // DCC+SEG LUTs data
00210               bookHistos(dtChId,"Segment","DCC_PhitkvsPhitrig"+(*trigSrcIt));
00211               bookHistos(dtChId,"Segment","DCC_PhibtkvsPhibtrig"+(*trigSrcIt));
00212               bookHistos(dtChId,"Segment","DCC_PhiResidual"+(*trigSrcIt));
00213               bookHistos(dtChId,"Segment","DCC_PhiResidualvsLUTPhi"+(*trigSrcIt));
00214               bookHistos(dtChId,"Segment","DCC_PhibResidual"+(*trigSrcIt));
00215             }
00216           }     
00217         }                 
00218       }// end of loop
00219     }
00220     
00221   }
00222 
00223 }
00224 
00225 
00226 void DTLocalTriggerTask::beginLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
00227   
00228   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerTask") << "[DTLocalTriggerTask]: Begin of LS transition" << endl;
00229   
00230   if(lumiSeg.id().luminosityBlock()%parameters.getUntrackedParameter<int>("ResetCycle", 3) == 0) {
00231     for(map<uint32_t, map<string, MonitorElement*> > ::const_iterator histo = digiHistos.begin();
00232         histo != digiHistos.end();
00233         histo++) {
00234       for(map<string, MonitorElement*> ::const_iterator ht = (*histo).second.begin();
00235           ht != (*histo).second.end();
00236           ht++) {
00237         (*ht).second->Reset();
00238       }
00239     }
00240   }
00241   
00242 }
00243 
00244 
00245 void DTLocalTriggerTask::endJob(){
00246 
00247   LogVerbatim("DTDQM|DTMonitorModule|DTLocalTriggerTask") << "[DTLocalTriggerTask]: analyzed " << nevents << " events" << endl;
00248   if (parameters.getUntrackedParameter<bool>("process_dcc", true)) dbe->rmdir(topFolder(0)); //DDU folder
00249   if (parameters.getUntrackedParameter<bool>("process_ros", true)) dbe->rmdir(topFolder(1)); //DCC folder
00250 
00251 }
00252 
00253 
00254 void DTLocalTriggerTask::analyze(const edm::Event& e, const edm::EventSetup& c){
00255   
00256   string dcc_label  = parameters.getUntrackedParameter<string>("dcc_label", "dttpgprod");
00257   string ros_label  = parameters.getUntrackedParameter<string>("ros_label", "dtunpacker");
00258   string seg_label  = parameters.getUntrackedParameter<string>("seg_label", "dt4DSegments");
00259 
00260   if (!nevents){
00261 
00262     edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh;
00263     e.getByLabel(dcc_label, l1DTTPGPh);
00264     edm::Handle<L1MuDTChambThContainer> l1DTTPGTh;
00265     e.getByLabel(dcc_label, l1DTTPGTh);
00266     useDCC = (l1DTTPGPh.isValid() || l1DTTPGTh.isValid()) && parameters.getUntrackedParameter<bool>("process_dcc", true) ;
00267 
00268     Handle<DTLocalTriggerCollection> l1DDUTrigs;
00269     e.getByLabel(ros_label,l1DDUTrigs);
00270     useDDU = l1DDUTrigs.isValid() && parameters.getUntrackedParameter<bool>("process_ros", true) ;
00271 
00272     Handle<DTRecSegment4DCollection> all4DSegments;
00273     e.getByLabel(seg_label, all4DSegments);  
00274     useSEG = all4DSegments.isValid() && parameters.getUntrackedParameter<bool>("process_seg", true) ;
00275     
00276   }
00277 
00278   nevents++;
00279     
00280   triggerSource(e);  
00281 
00282   if ( useDCC ) {   
00283     edm::Handle<L1MuDTChambPhContainer> l1DTTPGPh;
00284     e.getByLabel(dcc_label,l1DTTPGPh);
00285     vector<L1MuDTChambPhDigi>*  l1PhTrig = l1DTTPGPh->getContainer();
00286 
00287     edm::Handle<L1MuDTChambThContainer> l1DTTPGTh;
00288     e.getByLabel(dcc_label,l1DTTPGTh);
00289     vector<L1MuDTChambThDigi>*  l1ThTrig = l1DTTPGTh->getContainer();
00290 
00291     runDCCAnalysis(l1PhTrig,l1ThTrig);
00292   }  
00293   if ( useDDU ) {
00294     Handle<DTLocalTriggerCollection> l1DDUTrigs;
00295     e.getByLabel(ros_label,l1DDUTrigs);
00296 
00297     runDDUAnalysis(l1DDUTrigs);
00298   }
00299   if ( !tpMode && useSEG ) {
00300     Handle<DTRecSegment4DCollection> segments4D;
00301     e.getByLabel(seg_label, segments4D);  
00302     
00303     runSegmentAnalysis(segments4D);
00304   } 
00305   if ( !tpMode && useDCC && useDDU ) {
00306     runDDUvsDCCAnalysis(trigsrc);
00307   }
00308 
00309 }
00310 
00311 
00312 void DTLocalTriggerTask::bookBarrelHistos(string histoTag) {
00313 
00314   bool isDCC = histoTag.substr(0,3) == "DCC";
00315   dbe->setCurrentFolder(topFolder(isDCC));
00316   if (histoTag == "DCC_ErrorsChamberID") {
00317     dcc_IDDataErrorPlot = dbe->book1D(histoTag.c_str(),"DCC Data ID Error",5,-2,3);
00318     dcc_IDDataErrorPlot->setAxisTitle("wheel",1);
00319   }
00320   
00321   return;
00322 
00323 }
00324 
00325 void DTLocalTriggerTask::bookHistos(const DTChamberId& dtCh, string folder, string histoTag) {
00326   
00327   int wh=dtCh.wheel();          
00328   int sc=dtCh.sector(); 
00329   stringstream wheel; wheel << wh;      
00330   stringstream station; station << dtCh.station();      
00331   stringstream sector; sector << sc;    
00332 
00333   double minBX=0;
00334   double maxBX=0;
00335   int  rangeBX=0;
00336 
00337   string histoType = histoTag.substr(4,histoTag.find("_",4)-4);
00338   bool isDCC = histoTag.substr(0,3) == "DCC"; 
00339 
00340   dbe->setCurrentFolder(topFolder(isDCC) + "Wheel" + wheel.str() +
00341                         "/Sector" + sector.str() +
00342                         "/Station" + station.str() + "/" + folder);
00343 
00344   string histoName = histoTag + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
00345     
00346   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerTask") << "[DTLocalTriggerTask]: booking " << topFolder(isDCC) << "Wheel" << wheel.str()
00347                                                        << "/Sector" << sector.str()
00348                                                        << "/Station"<< station.str() << "/" << folder << "/" << histoName << endl;
00349     
00350   if (histoType.find("BX") != string::npos){
00351     if (histoTag.substr(0,3) == "DCC"){
00352       minBX= parameters.getUntrackedParameter<int>("minBXDCC",-2) - 0.5;
00353       maxBX= parameters.getUntrackedParameter<int>("maxBXDCC",2) + 0.5;
00354     }
00355     else {
00356       minBX= parameters.getUntrackedParameter<int>("minBXDDU",0) - 0.5;
00357       maxBX= parameters.getUntrackedParameter<int>("maxBXDDU",20) + 0.5;
00358     }
00359     rangeBX = (int)(maxBX-minBX);
00360   }
00361     
00362   if ( folder == "LocalTriggerPhi") {
00363       
00364     if( histoType == "BXvsQual" ){
00365       (digiHistos[dtCh.rawId()])[histoTag] = 
00366         dbe->book2D(histoName,"BX vs trigger quality",7,-0.5,6.5,rangeBX,minBX,maxBX);
00367       setQLabels((digiHistos[dtCh.rawId()])[histoTag],1);
00368       return ;
00369     }
00370     if( histoType == "BestQual" ){
00371       (digiHistos[dtCh.rawId()])[histoTag] = 
00372         dbe->book1D(histoName,"Trigger quality of best primitives",7,-0.5,6.5);
00373       setQLabels((digiHistos[dtCh.rawId()])[histoTag],1);
00374       return ;
00375     }
00376     if( histoType == "QualvsPhirad" ){
00377       (digiHistos[dtCh.rawId()])[histoTag] = 
00378         dbe->book2D(histoName,"Trigger quality vs local position",100,-500.,500.,7,-0.5,6.5);
00379       setQLabels((digiHistos[dtCh.rawId()])[histoTag],2);
00380       return ;
00381     }
00382     if( histoType == "QualvsPhibend" ) { 
00383       (digiHistos[dtCh.rawId()])[histoTag] = 
00384         dbe->book2D(histoName,"Trigger quality vs local direction",200,-40.,40.,7,-0.5,6.5);
00385       setQLabels((digiHistos[dtCh.rawId()])[histoTag],2);      
00386       return ;
00387     }
00388     if( histoType == "Flag1stvsQual" ) {
00389       (digiHistos[dtCh.rawId()])[histoTag] = 
00390         dbe->book2D(histoName,"1st/2nd trig flag vs quality",7,-0.5,6.5,2,-0.5,1.5);
00391       setQLabels((digiHistos[dtCh.rawId()])[histoTag],1);
00392       return ;
00393     }
00394     if( histoType == "QualDDUvsQualDCC" ){
00395       (digiHistos[dtCh.rawId()])[histoTag] = 
00396         dbe->book2D(histoName,"DDU quality vs DCC quality",8,-1.5,6.5,8,-1.5,6.5);
00397       setQLabels((digiHistos[dtCh.rawId()])[histoTag],1);
00398       setQLabels((digiHistos[dtCh.rawId()])[histoTag],2);      
00399       return ;
00400     }    
00401       
00402   }
00403   else if ( folder == "LocalTriggerTheta")   {
00404       
00405     if( histoType == "PositionvsBX" ) {
00406       (digiHistos[dtCh.rawId()])[histoTag] = 
00407         dbe->book2D(histoName,"Theta trigger position vs BX",rangeBX,minBX,maxBX,7,-0.5,6.5);
00408       return ;
00409     }
00410     if( histoType == "PositionvsQual" ) {
00411       (digiHistos[dtCh.rawId()])[histoTag] = 
00412         dbe->book2D(histoName,"Theta trigger position vs quality",7,-0.5,6.5,7,-0.5,6.5);
00413       setQLabels((digiHistos[dtCh.rawId()])[histoTag],1);
00414       return ;
00415     }  
00416     if( histoType == "ThetaBXvsQual" ) {
00417       (digiHistos[dtCh.rawId()])[histoTag] = 
00418         dbe->book2D(histoName,"BX vs trigger quality",7,-0.5,6.5,rangeBX,minBX,maxBX);
00419       setQLabels((digiHistos[dtCh.rawId()])[histoTag],1);
00420     }
00421     if( histoType == "ThetaBestQual" ){
00422       (digiHistos[dtCh.rawId()])[histoTag] = 
00423         dbe->book1D(histoName,"Trigger quality of best primitives (theta)",7,-0.5,6.5);
00424       setQLabels((digiHistos[dtCh.rawId()])[histoTag],1);
00425       return ;
00426     }
00427 
00428   }
00429   else if ( folder == "Segment")   {
00430       
00431     if( histoType.find("TrackThetaPosvsAngle" ) == 0 ) {
00432 
00433       string histoLabel = "Position vs Angle (theta)";
00434       if (histoType.find("andTrigH") != string::npos) histoLabel += " for H triggers";
00435       else if (histoType.find("andTrig") != string::npos) histoLabel += " for triggers";
00436 
00437       float min,max;
00438       int nbins;
00439       trigGeomUtils->thetaRange(dtCh,min,max,nbins);
00440       (digiHistos[dtCh.rawId()])[histoTag] = 
00441         dbe->book2D(histoName,histoLabel,16,-40.,40.,nbins,min,max);
00442       return ;
00443     }
00444     if( histoType.find("TrackPosvsAngle") == 0 ){
00445 
00446       float min,max;
00447       int nbins;
00448       trigGeomUtils->phiRange(dtCh,min,max,nbins);
00449 
00450       string histoLabel = "Position vs Angle (phi)";
00451       if (histoType.find("andTrigHHHL")  != string::npos) histoLabel += " for HH/HL triggers";
00452       else if (histoType.find("andTrig") != string::npos) histoLabel += " for triggers";
00453 
00454       (digiHistos[dtCh.rawId()])[histoTag] = 
00455         dbe->book2D(histoName,histoLabel,16,-40.,40.,nbins,min,max);
00456       return ;
00457     }
00458     if( histoType == "PhitkvsPhitrig" ){ 
00459       (digiHistos[dtCh.rawId()])[histoTag] = 
00460         dbe->book2D(histoName,"Local position: segment vs trigger",100,-500.,500.,100,-500.,500.);
00461       return ;
00462     }
00463     if( histoType == "PhibtkvsPhibtrig" ){ 
00464       (digiHistos[dtCh.rawId()])[histoTag] = 
00465         dbe->book2D(histoName,"Local direction : segment vs trigger",200,-40.,40.,200,-40.,40.);
00466       return ;
00467     }
00468     if( histoType == "PhiResidual" ){ 
00469       (digiHistos[dtCh.rawId()])[histoTag] = 
00470         dbe->book1D(histoName,"Trigger local position - Segment local position (correlated triggers)",400,-10.,10.);
00471       return ;
00472     }
00473     if( histoType == "PhibResidual" ){ 
00474       (digiHistos[dtCh.rawId()])[histoTag] = 
00475         dbe->book1D(histoName,"Trigger local direction - Segment local direction (correlated triggers)",500,-10.,10.);
00476       return ;
00477     }
00478     if( histoType == "HitstkvsQualtrig" ){ 
00479       (digiHistos[dtCh.rawId()])[histoTag] = 
00480         dbe->book2D(histoName,"Segment hits (phi) vs trigger quality",7,-0.5,6.5,10,0.5,10.5);
00481       setQLabels((digiHistos[dtCh.rawId()])[histoTag],1);
00482       return ;
00483     }
00484 
00485   }
00486 
00487 }
00488 
00489 void DTLocalTriggerTask::bookWheelHistos(int wh, string histoTag) {
00490   
00491   stringstream wheel; wheel << wh;      
00492 
00493   string histoType = histoTag.substr(4,histoTag.find("_",4)-4);
00494   bool isDCC = histoTag.substr(0,3) == "DCC"; 
00495 
00496   dbe->setCurrentFolder(topFolder(isDCC) + "Wheel" + wheel.str() + "/");
00497 
00498   string histoName = histoTag + "_W" + wheel.str();
00499     
00500   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerTask") << "[DTLocalTriggerTask]: booking " << topFolder(isDCC) 
00501                                                        << "Wheel" << wheel.str() << "/" << histoName << endl;
00502     
00503   if( histoType.find("BXDiff") != string::npos ){
00504     MonitorElement *me = dbe->bookProfile2D(histoName,"DDU-DCC BX Difference",12,1,13,4,1,5,0.,20.);
00505     me->setAxisTitle("Sector",1);
00506     me->setAxisTitle("station",2);
00507     (wheelHistos[wh])[histoTag] = me;
00508     return ;
00509   }
00510 
00511 }
00512 
00513 void DTLocalTriggerTask::runDCCAnalysis( std::vector<L1MuDTChambPhDigi>* phTrigs, 
00514                                          std::vector<L1MuDTChambThDigi>* thTrigs ){
00515 
00516   string histoType ;
00517   string histoTag ;
00518 
00519   // define best quality trigger segment (phi and theta)  
00520   // in any station start from 1 and zero is kept empty
00521   for (int i=0;i<5;++i){
00522     for (int j=0;j<6;++j){
00523       for (int k=0;k<13;++k){
00524         phcode_best[j][i][k] = -1;
00525         thcode_best[j][i][k] = -1;
00526       }
00527     }
00528   }
00529 
00530   vector<L1MuDTChambPhDigi>::const_iterator iph  = phTrigs->begin();
00531   vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end();
00532   for(; iph !=iphe ; ++iph) {
00533 
00534     int phwheel = iph->whNum();
00535     int phsec   = iph->scNum() + 1; // SM The track finder goes from 0 to 11. I need them from 1 to 12 !!!!!
00536     int phst    = iph->stNum();
00537     int phbx    = iph->bxNum();
00538     int phcode  = iph->code();
00539     int phi1st  = iph->Ts2Tag();
00540 
00541     // FIXME: workaround for DCC data with station ID 
00542     if(phst == 0) {
00543       dcc_IDDataErrorPlot->Fill(phwheel);
00544       continue;
00545     }
00546 
00547     if(phcode>phcode_best[phwheel+3][phst][phsec] && phcode<7) {
00548       phcode_best[phwheel+3][phst][phsec]=phcode; 
00549       iphbest[phwheel+3][phst][phsec] = &(*iph);
00550     }
00551       
00552     DTChamberId dtChId(phwheel,phst,phsec);
00553       
00554     float x     = trigGeomUtils->trigPos(&(*iph));
00555     float angle = trigGeomUtils->trigDir(&(*iph));
00556     uint32_t indexCh = dtChId.rawId();
00557     
00558     map<string, MonitorElement*> &innerME = digiHistos[indexCh];
00559     if (innerME.find("DCC_BXvsQual"+trigsrc) == innerME.end()){
00560       if (tpMode) {
00561         bookHistos(dtChId,"LocalTriggerPhi","DCC_BXvsQual"+trigsrc);
00562         bookHistos(dtChId,"LocalTriggerPhi","DCC_QualvsPhirad"+trigsrc);
00563       }
00564       else {
00565         bookHistos(dtChId,"LocalTriggerPhi","DCC_BXvsQual"+trigsrc);
00566         bookHistos(dtChId,"LocalTriggerPhi","DCC_Flag1stvsQual"+trigsrc);
00567         if (detailedAnalysis) {
00568           bookHistos(dtChId,"LocalTriggerPhi","DCC_QualvsPhirad"+trigsrc);
00569           bookHistos(dtChId,"LocalTriggerPhi","DCC_QualvsPhibend"+trigsrc);
00570         }
00571       }
00572     }
00573 
00574     if (tpMode) {
00575       innerME.find("DCC_BXvsQual"+trigsrc)->second->Fill(phcode,phbx-phi1st);    // SM BX vs Qual Phi view (1st tracks) 
00576       innerME.find("DCC_QualvsPhirad"+trigsrc)->second->Fill(x,phcode);          // SM Qual vs radial angle Phi view
00577     }
00578     else {
00579       innerME.find("DCC_BXvsQual"+trigsrc)->second->Fill(phcode,phbx-phi1st);    // SM BX vs Qual Phi view (1st tracks) 
00580       innerME.find("DCC_Flag1stvsQual"+trigsrc)->second->Fill(phcode,phi1st);    // SM Qual 1st/2nd track flag Phi view
00581       if (detailedAnalysis) {
00582         innerME.find("DCC_QualvsPhirad"+trigsrc)->second->Fill(x,phcode);          // SM Qual vs radial angle Phi view
00583         innerME.find("DCC_QualvsPhibend"+trigsrc)->second->Fill(angle,phcode);     // SM Qual vs bending Phi view
00584       }
00585     }
00586       
00587   } 
00588 
00589   if (doDCCTheta) {
00590     int thcode[7];
00591     vector<L1MuDTChambThDigi>::const_iterator ith  = thTrigs->begin();
00592     vector<L1MuDTChambThDigi>::const_iterator ithe = thTrigs->end();
00593     for(; ith != ithe; ++ith) {
00594       int thwheel = ith->whNum();
00595       int thsec   = ith->scNum() + 1; // SM The track finder goes from 0 to 11. I need them from 1 to 12 !!!!!
00596       int thst    = ith->stNum();
00597       int thbx    = ith->bxNum();
00598       
00599       for (int pos=0; pos<7; pos++) {
00600         thcode[pos] = ith->code(pos);
00601 
00602         if(thcode[pos]>thcode_best[thwheel+3][thst][thsec] ) {
00603           thcode_best[thwheel+3][thst][thsec]=thcode[pos]; 
00604           ithbest[thwheel+3][thst][thsec] = &(*ith);
00605         }
00606       } 
00607       
00608       DTChamberId dtChId(thwheel,thst,thsec);
00609       uint32_t indexCh = dtChId.rawId();   
00610 
00611       map<string, MonitorElement*> &innerME = digiHistos[indexCh];
00612       if (innerME.find("DCC_PositionvsBX"+trigsrc) == innerME.end()){
00613         bookHistos(dtChId,"LocalTriggerTheta","DCC_PositionvsBX"+trigsrc);
00614 //      if (detailedAnalysis) {
00615 //        bookHistos(dtChId,"LocalTriggerTheta","DCC_PositionvsBX"+trigsrc);
00616 //        bookHistos(dtChId,"LocalTriggerTheta","DCC_PositionvsQual"+trigsrc);
00617 //      }
00618       }
00619 
00620       for (int pos=0; pos<7; pos++) { //SM fill position for non zero position bit in theta view
00621         if(thcode[pos]>0){
00622           innerME.find("DCC_PositionvsBX"+trigsrc)->second->Fill(thbx,pos);          // SM BX vs Position Theta view
00623 //        if (detailedAnalysis) {
00624 //          int thqual = (thcode[pos]/2)*2+1;
00625 //          innerME.find("DCC_ThetaBXvsQual"+trigsrc)->second->Fill(thqual,thbx);      // SM BX vs Code Theta view
00626 //          innerME.find("DCC_PositionvsQual"+trigsrc)->second->Fill(thqual,pos);      // SM Code vs Position Theta view
00627 //        }
00628         }
00629       }
00630     }
00631   }
00632 
00633   
00634   // Fill Quality plots with best DCC triggers in phi & theta
00635   if (!tpMode) {
00636     for (int st=1;st<5;++st){
00637       for (int wh=-2;wh<3;++wh){
00638         for (int sc=1;sc<13;++sc){
00639           if (phcode_best[wh+3][st][sc]>-1 && phcode_best[wh+3][st][sc]<7){
00640             DTChamberId id(wh,st,sc);
00641             uint32_t indexCh = id.rawId();
00642             map<string, MonitorElement*> &innerME = digiHistos[indexCh];
00643             if (innerME.find("DCC_BestQual"+trigsrc) == innerME.end())
00644               bookHistos(id,"LocalTriggerPhi","DCC_BestQual"+trigsrc);
00645             innerME.find("DCC_BestQual"+trigsrc)->second->Fill(phcode_best[wh+3][st][sc]);  // Best Qual Trigger Phi view
00646           }
00647 //        if (thcode_best[wh+3][st][sc]>0){
00648 //          DTChamberId id(wh,st,sc);
00649 //          uint32_t indexCh = id.rawId();
00650 //          map<string, MonitorElement*> &innerME = digiHistos[indexCh];
00651 //          if (innerME.find("DCC_ThetaBestQual"+trigsrc) == innerME.end())
00652 //            bookHistos(id,"LocalTriggerTheta","DCC_ThetaBestQual"+trigsrc);
00653 //          innerME.find("DCC_ThetaBestQual"+trigsrc)->second->Fill(thcode_best[wh+3][st][sc]); // Best Qual Trigger Theta view
00654 //        }
00655         }
00656       }
00657     }
00658   }
00659 
00660 }
00661 
00662 void DTLocalTriggerTask::runDDUAnalysis(Handle<DTLocalTriggerCollection>& trigsDDU){
00663     
00664   DTLocalTriggerCollection::DigiRangeIterator detUnitIt;
00665     
00666   for (int i=0;i<5;++i){
00667     for (int j=0;j<6;++j){
00668       for (int k=0;k<13;++k){
00669         dduphcode_best[j][i][k] = -1;
00670         dduthcode_best[j][i][k] = -1;
00671       }
00672     }
00673   }
00674 
00675   for (detUnitIt=trigsDDU->begin();
00676        detUnitIt!=trigsDDU->end();
00677        ++detUnitIt){
00678       
00679     const DTChamberId& id = (*detUnitIt).first;
00680     const DTLocalTriggerCollection::Range& range = (*detUnitIt).second;
00681     uint32_t indexCh = id.rawId();
00682     map<string, MonitorElement*> &innerME = digiHistos[indexCh];
00683       
00684     int wh = id.wheel();
00685     int sec = id.sector();
00686     int st = id.station();
00687 
00688     for (DTLocalTriggerCollection::const_iterator trigIt = range.first;
00689          trigIt!=range.second;
00690          ++trigIt){
00691         
00692       int bx = trigIt->bx();
00693       int quality = trigIt->quality();
00694       int thqual = trigIt->trTheta();
00695       int flag1st = trigIt->secondTrack() ? 1 : 0;
00696 
00697       // check if SC data exist: fill for any trigger
00698       if( quality>-1 && quality<7 ) {     // it is a phi trigger
00699         
00700         if(quality>dduphcode_best[wh+3][st][sec]) {
00701           dduphcode_best[wh+3][st][sec]=quality;
00702           iphbestddu[wh+3][st][sec] = &(*trigIt);
00703         }
00704         
00705         if (innerME.find("DDU_BXvsQual"+trigsrc) == innerME.end()){
00706           bookHistos(id,"LocalTriggerPhi","DDU_BXvsQual"+trigsrc);
00707           bookHistos(id,"LocalTriggerPhi","DDU_Flag1stvsQual"+trigsrc);
00708         }
00709 
00710         if(tpMode) {      
00711           innerME.find("DDU_BXvsQual"+trigsrc)->second->Fill(quality,bx-flag1st);     // SM BX vs Qual Phi view 
00712         }
00713         else {
00714           innerME.find("DDU_BXvsQual"+trigsrc)->second->Fill(quality,bx-flag1st);     // SM BX vs Qual Phi view 
00715           innerME.find("DDU_Flag1stvsQual"+trigsrc)->second->Fill(quality,flag1st); // SM Quality vs 1st/2nd track flag Phi view
00716         }
00717       }
00718       if( thqual>0 && !tpMode ) {  // it is a theta trigger
00719         
00720         if(thqual>dduthcode_best[wh+3][st][sec] ) {
00721           dduthcode_best[wh+3][st][sec]=thqual; 
00722         }
00723         if (innerME.find("DDU_ThetaBXvsQual"+trigsrc) == innerME.end())
00724           bookHistos(id,"LocalTriggerTheta","DDU_ThetaBXvsQual"+trigsrc);
00725         innerME.find("DDU_ThetaBXvsQual"+trigsrc)->second->Fill(thqual,bx);     // SM BX vs Qual Theta view
00726         
00727       }
00728     }
00729     
00730     // Fill Quality plots with best ddu triggers in phi & theta
00731     if (!tpMode) {
00732       if (dduphcode_best[wh+3][st][sec]>-1 &&
00733           dduphcode_best[wh+3][st][sec]<7){
00734         if (innerME.find("DDU_BestQual"+trigsrc) == innerME.end())
00735           bookHistos(id,"LocalTriggerPhi","DDU_BestQual"+trigsrc);
00736         innerME.find("DDU_BestQual"+trigsrc)->second->Fill(dduphcode_best[wh+3][st][sec]);  // Best Qual Trigger Phi view
00737       }
00738       if (dduthcode_best[wh+3][st][sec]>0){
00739         if (innerME.find("DDU_ThetaBestQual"+trigsrc) == innerME.end())
00740           bookHistos(id,"LocalTriggerTheta","DDU_ThetaBestQual"+trigsrc);
00741         innerME.find("DDU_ThetaBestQual"+trigsrc)->second->Fill(dduthcode_best[wh+3][st][sec]); // Best Qual Trigger Theta view
00742       }  
00743     }
00744   }
00745   
00746 }
00747 
00748 
00749 void DTLocalTriggerTask::runSegmentAnalysis(Handle<DTRecSegment4DCollection>& segments4D){    
00750 
00751   DTRecSegment4DCollection::const_iterator track;
00752 
00753   // Find best tracks & good tracks
00754   memset(track_ok,false,450*sizeof(bool));
00755 
00756   DTRecSegment4DCollection::id_iterator chamberId;
00757   vector<const DTRecSegment4D*> best4DSegments;
00758 
00759   // Preliminary loop finds best 4D Segment and high quality ones
00760   for (chamberId = segments4D->id_begin(); chamberId != segments4D->id_end(); ++chamberId){
00761 
00762     DTRecSegment4DCollection::range  range = segments4D->get(*chamberId);
00763     const DTRecSegment4D* tmpBest=0;
00764     int tmpdof = 0;
00765     int dof = 0;
00766 
00767     for ( track = range.first; track != range.second; ++track){
00768 
00769       if( (*track).hasPhi() ) {
00770         
00771         dof = (*track).phiSegment()->degreesOfFreedom();
00772         if ( dof>tmpdof ){
00773           tmpBest = &(*track);
00774           tmpdof = dof;
00775         
00776           int wheel = (*track).chamberId().wheel();
00777           int sector = (*track).chamberId().sector();
00778           int station = (*track).chamberId().station();
00779           if (sector==13){
00780             sector=4;
00781           }
00782           else if (sector==14){
00783             sector=10;
00784           }
00785           track_ok[wheel+3][station][sector] = (!track_ok[wheel+3][station][sector] && dof>=2);
00786         }
00787 
00788       }
00789     }
00790     if (tmpBest) best4DSegments.push_back(tmpBest);
00791   }
00792     
00793   vector<const DTRecSegment4D*>::const_iterator btrack;
00794 
00795   for ( btrack = best4DSegments.begin(); btrack != best4DSegments.end(); ++btrack ){
00796 
00797     if( (*btrack)->hasPhi() ) { // Phi component
00798         
00799       int wheel    = (*btrack)->chamberId().wheel();
00800       int station  = (*btrack)->chamberId().station();
00801       int sector   = (*btrack)->chamberId().sector();
00802       int scsector = 0;
00803       float x_track, y_track, x_angle, y_angle;
00804       trigGeomUtils->computeSCCoordinates((*btrack),scsector,x_track,x_angle,y_track,y_angle);
00805       int nHitsPhi = (*btrack)->phiSegment()->degreesOfFreedom()+2;
00806 
00807       DTChamberId dtChId(wheel,station,sector);      // get chamber for LUTs histograms (Sectors 1 to 14)
00808       uint32_t indexCh = dtChId.rawId(); 
00809       map<string, MonitorElement*> &innerMECh = digiHistos[indexCh];
00810       
00811       DTChamberId dtChIdSC = DTChamberId(wheel,station,scsector);  // get chamber for histograms SC granularity (sectors 1 to 12)
00812       indexCh = dtChIdSC.rawId(); 
00813       map<string, MonitorElement*> &innerME = digiHistos[indexCh];
00814 
00815       if (useDDU && 
00816           dduphcode_best[wheel+3][station][scsector] > -1 && 
00817           dduphcode_best[wheel+3][station][scsector] < 7 ) {
00818 
00819         // SM hits of the track vs quality of the trigger
00820         if (innerME.find("DDU_HitstkvsQualtrig"+trigsrc) == innerME.end())
00821           bookHistos(dtChIdSC,"Segment","DDU_HitstkvsQualtrig"+trigsrc);
00822         innerME.find("DDU_HitstkvsQualtrig"+trigsrc)->second->Fill(dduphcode_best[wheel+3][station][scsector],nHitsPhi);
00823 
00824       }
00825         
00826       if (useDCC &&
00827           phcode_best[wheel+3][station][scsector] > -1 && 
00828           phcode_best[wheel+3][station][scsector] < 7 ) {
00829         
00830         if (innerME.find("DCC_HitstkvsQualtrig"+trigsrc) == innerME.end()){
00831           bookHistos(dtChIdSC,"Segment","DCC_HitstkvsQualtrig"+trigsrc);
00832         }
00833         innerME.find("DCC_HitstkvsQualtrig"+trigsrc)->second->Fill(phcode_best[wheel+3][station][scsector],nHitsPhi);
00834      
00835         if (phcode_best[wheel+3][station][scsector]>3 && nHitsPhi>=7){
00836           
00837           float x_trigger     = trigGeomUtils->trigPos(iphbest[wheel+3][station][scsector]);
00838           float angle_trigger = trigGeomUtils->trigDir(iphbest[wheel+3][station][scsector]);
00839           trigGeomUtils->trigToSeg(station,x_trigger,x_angle);
00840             
00841           if (innerMECh.find("DCC_PhibResidual"+trigsrc) == innerMECh.end()){
00842             bookHistos(dtChId,"Segment","DCC_PhiResidual"+trigsrc);
00843             bookHistos(dtChId,"Segment","DCC_PhibResidual"+trigsrc);
00844             bookHistos(dtChId,"Segment","DCC_PhitkvsPhitrig"+trigsrc);
00845             bookHistos(dtChId,"Segment","DCC_PhibtkvsPhibtrig"+trigsrc);
00846             bookHistos(dtChId,"Segment","DCC_HitstkvsQualtrig"+trigsrc);
00847           }
00848 
00849           innerMECh.find("DCC_PhitkvsPhitrig"+trigsrc)->second->Fill(x_trigger,x_track);
00850           innerMECh.find("DCC_PhibtkvsPhibtrig"+trigsrc)->second->Fill(angle_trigger,x_angle);
00851           innerMECh.find("DCC_PhiResidual"+trigsrc)->second->Fill(x_trigger-x_track);
00852           innerMECh.find("DCC_PhibResidual"+trigsrc)->second->Fill(angle_trigger-x_angle);
00853         }
00854 
00855 
00856 
00857       }
00858 
00859       
00860       if (useDCC) {
00861           
00862         // check for triggers elsewhere in the sector
00863         bool trigFlagDCC =false;
00864         for (int ist=1; ist<5; ist++){
00865           if (ist!=station &&
00866               phcode_best[wheel+3][ist][scsector]>=2 && 
00867               phcode_best[wheel+3][ist][scsector]<7 &&
00868               track_ok[wheel+3][ist][scsector]==true){
00869             trigFlagDCC = true;
00870             break;
00871           }
00872         }
00873           
00874         if (trigFlagDCC && fabs(x_angle)<40. && nHitsPhi>=7){
00875 
00876           if (innerME.find("DCC_TrackPosvsAngle"+trigsrc) == innerME.end()){
00877             bookHistos(dtChIdSC,"Segment","DCC_TrackPosvsAngle"+trigsrc);
00878             bookHistos(dtChIdSC,"Segment","DCC_TrackPosvsAngleandTrig"+trigsrc);
00879             bookHistos(dtChIdSC,"Segment","DCC_TrackPosvsAngleandTrigHHHL"+trigsrc);
00880           }
00881             
00882           // position vs angle of track for reconstruced tracks (denom. for trigger efficiency)
00883           innerME.find("DCC_TrackPosvsAngle"+trigsrc)->second->Fill(x_angle,x_track);     
00884           if (phcode_best[wheel+3][station][scsector] >= 2 && phcode_best[wheel+3][station][scsector] < 7) {
00885             innerME.find("DCC_TrackPosvsAngleandTrig"+trigsrc)->second->Fill(x_angle,x_track);    
00886             if (phcode_best[wheel+3][station][scsector] > 4){  //HH & HL Triggers
00887               innerME.find("DCC_TrackPosvsAngleandTrigHHHL"+trigsrc)->second->Fill(x_angle,x_track);      
00888             }         
00889           }
00890 
00891         }
00892           
00893         if ((*btrack)->hasZed() && trigFlagDCC && fabs(y_angle)<40. && (*btrack)->zSegment()->degreesOfFreedom()>=1){
00894 
00895           if (innerME.find("DCC_TrackThetaPosvsAngle"+trigsrc) == innerME.end()){
00896             bookHistos(dtChIdSC,"Segment","DCC_TrackThetaPosvsAngle"+trigsrc);
00897             bookHistos(dtChIdSC,"Segment","DCC_TrackThetaPosvsAngleandTrig"+trigsrc);
00898 //          bookHistos(dtChIdSC,"Segment","DCC_TrackThetaPosvsAngleandTrigH"+trigsrc);          no DCC theta qual info
00899           }
00900 
00901           // position va angle of track for reconstruced tracks (denom. for trigger efficiency) along theta direction
00902           innerME.find("DCC_TrackThetaPosvsAngle"+trigsrc)->second->Fill(y_angle,y_track);
00903           if (thcode_best[wheel+3][station][scsector] > 0) {            
00904             innerME.find("DCC_TrackThetaPosvsAngleandTrig"+trigsrc)->second->Fill(y_angle,y_track);
00905 //          if (thcode_best[wheel+3][station][scsector] == 3) {
00906 //            innerME.find("DCC_TrackThetaPosvsAngleH"+trigsrc)->second->Fill(y_angle,y_track);
00907 //          }            
00908           }
00909 
00910         }  
00911       }
00912 
00913       if (useDDU) {
00914           
00915         // check for triggers elsewhere in the sector
00916         bool trigFlagDDU =false;
00917         for (int ist=1; ist<5; ist++){
00918           if (ist!=station &&
00919               dduphcode_best[wheel+3][ist][scsector]>=2 && 
00920               dduphcode_best[wheel+3][ist][scsector]<7 &&
00921               track_ok[wheel+3][ist][scsector]==true){
00922             trigFlagDDU = true;
00923             break;
00924           }
00925         }
00926 
00927         if (trigFlagDDU && fabs(x_angle)<40. && nHitsPhi>=7){
00928 
00929           if (innerME.find("DDU_TrackPosvsAngle"+trigsrc) == innerME.end()){
00930             bookHistos(dtChIdSC,"Segment","DDU_TrackPosvsAngle"+trigsrc);
00931             bookHistos(dtChIdSC,"Segment","DDU_TrackPosvsAngleandTrig"+trigsrc);
00932             bookHistos(dtChIdSC,"Segment","DDU_TrackPosvsAngleandTrigHHHL"+trigsrc);
00933           }
00934 
00935           // position vs angle of track for reconstruced tracks (denom. for trigger efficiency)
00936           innerME.find("DDU_TrackPosvsAngle"+trigsrc)->second->Fill(x_angle,x_track);
00937           if (dduphcode_best[wheel+3][station][scsector] >= 2 && dduphcode_best[wheel+3][station][scsector] < 7) {
00938           innerME.find("DDU_TrackPosvsAngleandTrig"+trigsrc)->second->Fill(x_angle,x_track);
00939             if (dduphcode_best[wheel+3][station][scsector] > 4){  //HH & HL Triggers
00940               innerME.find("DDU_TrackPosvsAngleandTrigHHHL"+trigsrc)->second->Fill(x_angle,x_track);      
00941             }         
00942           }
00943 
00944         }
00945           
00946         if ((*btrack)->hasZed() && trigFlagDDU && fabs(y_angle)<40. && (*btrack)->zSegment()->degreesOfFreedom()>=1){
00947 
00948           if (innerME.find("DDU_TrackThetaPosvsAngle"+trigsrc) == innerME.end()){
00949             bookHistos(dtChIdSC,"Segment","DDU_TrackThetaPosvsAngle"+trigsrc);
00950             bookHistos(dtChIdSC,"Segment","DDU_TrackThetaPosvsAngleandTrig"+trigsrc);
00951             bookHistos(dtChIdSC,"Segment","DDU_TrackThetaPosvsAngleandTrigH"+trigsrc);
00952           }
00953 
00954           // position va angle of track for reconstruced tracks (denom. for trigger efficiency) along theta direction
00955           innerME.find("DDU_TrackThetaPosvsAngle"+trigsrc)->second->Fill(y_angle,y_track);
00956           if (dduthcode_best[wheel+3][station][scsector] > 0) {         
00957             innerME.find("DDU_TrackThetaPosvsAngleandTrig"+trigsrc)->second->Fill(y_angle,y_track);
00958             if (dduthcode_best[wheel+3][station][scsector] == 3) {
00959               innerME.find("DDU_TrackThetaPosvsAngleandTrigH"+trigsrc)->second->Fill(y_angle,y_track);
00960             }            
00961           }
00962 
00963         }  
00964       }
00965     }
00966   } 
00967 
00968 }
00969 
00970 
00971 void DTLocalTriggerTask::runDDUvsDCCAnalysis(string& trigsrc){
00972 
00973   string histoType ;
00974   string histoTag ;
00975 
00976   for (int st=1;st<5;++st){
00977     for (int wh=-2;wh<3;++wh){
00978       for (int sc=1;sc<13;++sc){
00979         if ( (phcode_best[wh+3][st][sc]>-1 && phcode_best[wh+3][st][sc]<7) ||
00980              (dduphcode_best[wh+3][st][sc]>-1 && dduphcode_best[wh+3][st][sc]<7) ){
00981           DTChamberId id(wh,st,sc);
00982           uint32_t indexCh = id.rawId();
00983           map<string, MonitorElement*> &innerME = digiHistos[indexCh];
00984           if (innerME.find("COM_QualDDUvsQualDCC"+trigsrc) == innerME.end())
00985             bookHistos(id,"LocalTriggerPhi","COM_QualDDUvsQualDCC"+trigsrc);
00986           innerME.find("COM_QualDDUvsQualDCC"+trigsrc)->second->Fill(phcode_best[wh+3][st][sc],dduphcode_best[wh+3][st][sc]);
00987           if ( (phcode_best[wh+3][st][sc]>-1 && phcode_best[wh+3][st][sc]<7) &&
00988                (dduphcode_best[wh+3][st][sc]>-1 && dduphcode_best[wh+3][st][sc]<7) ){
00989             int bxDDU = iphbestddu[wh+3][st][sc]->bx() - iphbestddu[wh+3][st][sc]->secondTrack();
00990             int bxDCC = iphbest[wh+3][st][sc]->bxNum() - iphbest[wh+3][st][sc]->Ts2Tag();
00991             (wheelHistos[wh]).find("COM_BXDiff"+trigsrc)->second->Fill(sc,st,bxDDU-bxDCC);
00992           }
00993         }
00994       }
00995     }
00996   }
00997 
00998 }
00999 
01000 void DTLocalTriggerTask::setQLabels(MonitorElement* me, short int iaxis){
01001 
01002   TH1* histo = me->getTH1();
01003   if (!histo) return;
01004   
01005   TAxis* axis=0;
01006   if (iaxis==1) {
01007     axis=histo->GetXaxis();
01008   }
01009   else if(iaxis==2) {
01010     axis=histo->GetYaxis();
01011   }
01012   if (!axis) return;
01013 
01014   string labels[7] = {"LI","LO","HI","HO","LL","HL","HH"};
01015   int istart = axis->GetXmin()<-1 ? 2 : 1;
01016   for (int i=0;i<7;i++) {
01017     axis->SetBinLabel(i+istart,labels[i].c_str());
01018   }
01019 
01020 }
01021 
01022 void DTLocalTriggerTask::triggerSource(const edm::Event& e) {
01023   
01024   
01025   if (!isLocalRun){
01026     
01027     Handle<LTCDigiCollection> ltcdigis;
01028     e.getByLabel(ltcDigiCollectionTag, ltcdigis);
01029     
01030     for (std::vector<LTCDigi>::const_iterator ltc_it = ltcdigis->begin(); ltc_it != ltcdigis->end(); ltc_it++){
01031       
01032       size_t otherTriggerSum=0;
01033       for (size_t i = 1; i < 6; i++) {
01034         otherTriggerSum += size_t((*ltc_it).HasTriggered(i));
01035       }
01036       if ((*ltc_it).HasTriggered(0) && otherTriggerSum == 0) 
01037         trigsrc = "_DTonly";
01038       else if (!(*ltc_it).HasTriggered(0))
01039         trigsrc = "_NoDT";
01040       else if ((*ltc_it).HasTriggered(0) && otherTriggerSum > 0)
01041         trigsrc = "_DTalso";
01042       
01043     }
01044     return;
01045   }
01046 
01047   trigsrc = "";
01048   return;
01049 
01050 }