CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/DQM/L1TMonitor/src/L1TDTTPG.cc

Go to the documentation of this file.
00001 /*
00002  * \file L1TDTTPG.cc
00003  *
00004  * $Date: 2009/11/19 14:38:34 $
00005  * $Revision: 1.21 $
00006  * \author J. Berryhill
00007  *
00008  * $Log: L1TDTTPG.cc,v $
00009  * Revision 1.21  2009/11/19 14:38:34  puigh
00010  * modify beginJob
00011  *
00012  * Revision 1.20  2008/03/20 19:38:25  berryhil
00013  *
00014  *
00015  * organized message logger
00016  *
00017  * Revision 1.19  2008/03/14 20:35:46  berryhil
00018  *
00019  *
00020  * stripped out obsolete parameter settings
00021  *
00022  * rpc tpg restored with correct dn access and dbe handling
00023  *
00024  * Revision 1.18  2008/03/12 17:24:24  berryhil
00025  *
00026  *
00027  * eliminated log files, truncated HCALTPGXana histo output
00028  *
00029  * Revision 1.17  2008/03/10 09:29:52  lorenzo
00030  * added MEs
00031  *
00032  * Revision 1.16  2008/03/01 00:40:00  lat
00033  * DQM core migration.
00034  *
00035  * $Log: L1TDTTPG.cc,v $
00036  * Revision 1.21  2009/11/19 14:38:34  puigh
00037  * modify beginJob
00038  *
00039  * Revision 1.20  2008/03/20 19:38:25  berryhil
00040  *
00041  *
00042  * organized message logger
00043  *
00044  * Revision 1.19  2008/03/14 20:35:46  berryhil
00045  *
00046  *
00047  * stripped out obsolete parameter settings
00048  *
00049  * rpc tpg restored with correct dn access and dbe handling
00050  *
00051  * Revision 1.18  2008/03/12 17:24:24  berryhil
00052  *
00053  *
00054  * eliminated log files, truncated HCALTPGXana histo output
00055  *
00056  * Revision 1.17  2008/03/10 09:29:52  lorenzo
00057  * added MEs
00058  *
00059  * Revision 1.15  2008/01/22 18:56:01  muzaffar
00060  * include cleanup. Only for cc/cpp files
00061  *
00062  * Revision 1.14  2007/12/21 17:41:20  berryhil
00063  *
00064  *
00065  * try/catch removal
00066  *
00067  * Revision 1.13  2007/11/19 15:08:22  lorenzo
00068  * changed top folder name
00069  *
00070  * Revision 1.12  2007/08/15 18:56:25  berryhil
00071  *
00072  *
00073  * split histograms by bx; add Maiken's bx classifier plots
00074  *
00075  * Revision 1.11  2007/07/26 09:37:09  berryhil
00076  *
00077  *
00078  * set verbose false for all modules
00079  * set verbose fix for DTTPG tracks
00080  *
00081  * Revision 1.10  2007/07/25 09:03:58  berryhil
00082  *
00083  *
00084  * conform to DTTFFEDReader input tag.... for now
00085  *
00086  * Revision 1.9  2007/07/12 16:06:18  wittich
00087  * add simple phi output track histograms.
00088  * note that the label of this class is different than others
00089  * from the DTFFReader creates.
00090  *
00091  */
00092 #include <iostream>
00093 #include <sstream>
00094 #include <fstream>
00095 #include <vector>
00096 
00097 #include "FWCore/ServiceRegistry/interface/Service.h"
00098 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00099 
00100 #include "DQM/L1TMonitor/interface/L1TDTTPG.h"
00101 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhContainer.h"
00102 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambPhDigi.h"
00103 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThContainer.h"
00104 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTChambThDigi.h"
00105 #include "DataFormats/L1DTTrackFinder/interface/L1MuDTTrackContainer.h"
00106 #include "DQMServices/Core/interface/DQMStore.h"
00107 
00108 using namespace std;
00109 using namespace edm;
00110 
00111 L1TDTTPG::L1TDTTPG(const ParameterSet& ps)
00112   : dttpgSource_( ps.getParameter< InputTag >("dttpgSource") )
00113 {
00114 
00115   // verbosity switch
00116   verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
00117 
00118   if(verbose_) cout << "L1TDTTPG: constructor...." << endl;
00119 
00120 
00121   dbe = NULL;
00122   if ( ps.getUntrackedParameter<bool>("DQMStore", false) ) 
00123     {
00124       dbe = Service<DQMStore>().operator->();
00125       dbe->setVerbose(0);
00126     }
00127   
00128   outputFile_ = ps.getUntrackedParameter<string>("outputFile", "");
00129   if ( outputFile_.size() != 0 ) {
00130     cout << "L1T Monitoring histograms will be saved to " << outputFile_.c_str() << endl;
00131   }
00132 
00133   bool disable = ps.getUntrackedParameter<bool>("disableROOToutput", false);
00134   if(disable){
00135     outputFile_="";
00136   }
00137 
00138 
00139   if ( dbe !=NULL ) {
00140     dbe->setCurrentFolder("L1T/L1TDTTPG");
00141   }
00142 
00143 
00144 }
00145 
00146 L1TDTTPG::~L1TDTTPG()
00147 {
00148 }
00149 
00150 void L1TDTTPG::beginJob(void)
00151 {
00152 
00153   nev_ = 0;
00154 
00155   // get hold of back-end interface
00156   DQMStore* dbe = 0;
00157   dbe = Service<DQMStore>().operator->();
00158   if ( dbe ) {
00159     dbe->setCurrentFolder("L1T/L1TDTTPG");
00160     dbe->rmdir("L1T/L1TDTTPG");
00161   }
00162 
00163 
00164   if ( dbe ) 
00165     {
00166       dbe->setCurrentFolder("L1T/L1TDTTPG");
00167 
00168 
00169       //hist1[0]
00170       dttpgphbx[0] = dbe->book1D("BxEncoding_PHI",
00171                                  "Bunch encoding DTTF Phi",11,0,11);
00172       //hist1[1]
00173       dttpgphbx[1] = dbe->book1D("BxEncoding_OUT",
00174                                  "Bunch encoding DTTF Output",11,0,11);
00175 
00176       //hist1[2]
00177       dttpgphbx[2] = dbe->book1D("NumberOfSegmentsPHI_BunchNeg1",
00178                                  "Number of segments for bunch -1 Dttf Phi",
00179                                  20,0,20);
00180       //hist1[3]
00181       dttpgphbx[3] = dbe->book1D("NumberOfSegmentsPHI_Bunch0",
00182                                  "Number of segments for bunch 0 Dttf Phi",
00183                                  20,0,20);
00184       //hist1[4]
00185       dttpgphbx[4] = dbe->book1D("NumberOfSegmentsPHI_Bunch1",
00186                                  "Number of segments for bunch 1 Dttf Phi",
00187                                  20,0,20);
00188 
00189       //hist1[5]
00190       dttpgphbx[5] = dbe->book1D("NumberOfSegmentsOUT_BunchNeg1",
00191                                  "Number of segments for bunch -1 Dttf Output",
00192                                  20,0,20);
00193       //hist1[6] 
00194       dttpgphbx[6] = dbe->book1D("NumberOfSegmentsOUT_Bunch0",
00195                                  "Number of segments for bunch 0 Dttf Output",
00196                                  20,0,20);
00197       //hist1[7]
00198       dttpgphbx[7] = dbe->book1D("NumberOfSegmentsOUT_Bunch1",
00199                                  "Number of segments for bunch 1 Dttf Output",
00200                                  20,0,20);
00201 
00202       for(int i=0;i<2;i++){
00203         dttpgphbx[i]->setBinLabel(1,"None");
00204         dttpgphbx[i]->setBinLabel(3,"Only bx=-1");
00205         dttpgphbx[i]->setBinLabel(4,"Only bx= 0");
00206         dttpgphbx[i]->setBinLabel(5,"Only bx=+1");
00207         dttpgphbx[i]->setBinLabel(7,"Bx=-1,0");
00208         dttpgphbx[i]->setBinLabel(8,"Bx=-1,1");
00209         dttpgphbx[i]->setBinLabel(9,"Bx= 0,1");
00210         dttpgphbx[i]->setBinLabel(11,"All bx");
00211       }
00212    
00213       dttpgphbxcomp = dbe->book2D("BxEncoding_PHI_OUT",
00214                                   "Bunch encoding: DTTF Phi vs. Output",
00215                                   11,0,11,11,0,11);
00216       dttpgphbxcomp->setAxisTitle("DTTF (output)",1);
00217       dttpgphbxcomp->setAxisTitle("PHI-TF",2);
00218       for(int i=1;i<=2;i++){
00219         dttpgphbxcomp->setBinLabel(1,"None",i);
00220         dttpgphbxcomp->setBinLabel(3,"Only bx=-1",i);
00221         dttpgphbxcomp->setBinLabel(4,"Only bx= 0",i);
00222         dttpgphbxcomp->setBinLabel(5,"Only bx=+1",i);
00223         dttpgphbxcomp->setBinLabel(7,"Bx=-1,0",i);
00224         dttpgphbxcomp->setBinLabel(8,"Bx=-1,1",i);
00225         dttpgphbxcomp->setBinLabel(9,"Bx= 0,1",i);
00226         dttpgphbxcomp->setBinLabel(11,"All bx",i);
00227       }
00228 
00229       dttpgphntrack = dbe->book1D("DT_TPG_phi_ntrack", 
00230                                   "DT TPG phi ntrack", 20, -0.5, 19.5 ) ;  
00231       dttpgthntrack = dbe->book1D("DT_TPG_theta_ntrack", 
00232                                   "DT TPG theta ntrack", 20, -0.5, 19.5 ) ;  
00233 
00234       for (int ibx=0 ; ibx<=2; ibx++) {
00235         
00236         ostringstream bxnum;
00237         bxnum << ibx-1;
00238         string bxn;
00239         if (ibx<2)
00240           bxn = bxnum.str();
00241         else
00242           bxn = "+" + bxnum.str();
00243         
00244         // Phi
00245         dttpgphwheel[ibx] = dbe->book1D("DT_TPG_phi_wheel_number_"+bxn, 
00246                                             "DT TPG phi wheel number "+bxn, 5, -2.5, 2.5 ) ;  
00247         dttpgphsector[ibx] = dbe->book1D("DT_TPG_phi_sector_number_"+bxn, 
00248                                          "DT TPG phi sector number "+bxn, 12, -0.5, 11.5 );  
00249         dttpgphstation[ibx] = dbe->book1D("DT_TPG_phi_station_number_"+bxn, 
00250                                           "DT TPG phi station number "+bxn, 5, 0.5, 4.5 ) ;
00251 //      dttpgphphi[ibx] = dbe->book1D("DT_TPG_phi_"+bxn, 
00252 //                                    "DT TPG phi "+bxn, 100, -2100., 2100. ) ;  
00253 //      dttpgphphiB[ibx] = dbe->book1D("DT_TPG_phiB_"+bxn, 
00254 //                                     "DT TPG phiB "+bxn, 100, -550., 550. ) ;  
00255         dttpgphquality[ibx] = dbe->book1D("DT_TPG_phi_quality_"+bxn, 
00256                                           "DT TPG phi quality "+bxn, 8, -0.5, 7.5 ) ;  
00257         dttpgphts2tag[ibx] = dbe->book1D("DT_TPG_phi_Ts2Tag_"+bxn, 
00258                                          "DT TPG phi Ts2Tag "+bxn, 2, -0.5, 1.5 ) ;  
00259 //      dttpgphbxcnt[ibx] = dbe->book1D("DT_TPG_phi_BxCnt_"+bxn, 
00260 //                                      "DT TPG phi BxCnt "+bxn, 10, -0.5, 9.5 ) ;  
00261         dttpgphmapbx[ibx] = dbe->book2D("DT_TPG_phi_map_bx"+bxn,
00262                                       "Map of triggers per station (BX="+bxn+")",20,1,21,12,0,12);
00263         setMapPhLabel(dttpgphmapbx[ibx]);
00264 
00265         //Theta
00266         dttpgthbx[ibx] = dbe->book1D("DT_TPG_theta_bx_"+bxn, 
00267                                      "DT TPG theta bx "+bxn, 50, -24.5, 24.5 ) ;  
00268         dttpgthwheel[ibx] = dbe->book1D("DT_TPG_theta_wheel_number_"+bxn, 
00269                                         "DT TPG theta wheel number "+bxn, 5, -2.5, 2.5 ) ;  
00270         dttpgthsector[ibx] = dbe->book1D("DT_TPG_theta_sector_number_"+bxn, 
00271                                          "DT TPG theta sector number "+bxn, 12, -0.5, 11.5 ) ;  
00272         dttpgthstation[ibx] = dbe->book1D("DT_TPG_theta_station_number_"+bxn, 
00273                                           "DT TPG theta station number "+bxn, 5, -0.5, 4.5 ) ;  
00274         dttpgththeta[ibx] = dbe->book1D("DT_TPG_theta_"+bxn, 
00275                                         "DT TPG theta "+bxn, 20, -0.5, 19.5 ) ;  
00276         dttpgthquality[ibx] = dbe->book1D("DT_TPG_theta_quality_"+bxn, 
00277                                           "DT TPG theta quality "+bxn, 8, -0.5, 7.5 ) ;  
00278         dttpgthmapbx[ibx] = dbe->book2D("DT_TPG_theta_map_bx_"+bxn,
00279                                         "Map of triggers per station (BX="+bxn+")",15,1,16,12,0,12);
00280         setMapThLabel(dttpgthmapbx[ibx]);
00281 
00282         // Phi output
00283         dttf_p_phi[ibx] = dbe->book1D("dttf_p_phi_"+bxn, "dttf phi output #phi "+bxn, 256, 
00284                                       -0.5, 255.5);
00285         dttf_p_qual[ibx] = dbe->book1D("dttf_p_qual_"+bxn, "dttf phi output qual "+bxn, 8, -0.5, 7.5);
00286         dttf_p_q[ibx] = dbe->book1D("dttf_p_q_"+bxn, "dttf phi output q "+bxn, 2, -0.5, 1.5);
00287         dttf_p_pt[ibx] = dbe->book1D("dttf_p_pt_"+bxn, "dttf phi output p_{t} "+bxn, 32, -0.5, 31.5);
00288       
00289       }
00290 
00291       dttpgphmap = dbe->book2D("DT_TPG_phi_map",
00292                                "Map of triggers per station",20,1,21,12,0,12);
00293       dttpgphmapcorr = dbe->book2D("DT_TPG_phi_map_corr",
00294                                    "Map of correlated triggers per station",20,1,21,12,0,12);
00295       dttpgphmap2nd = dbe->book2D("DT_TPG_phi_map_2nd",
00296                                   "Map of second tracks per station",20,1,21,12,0,12);
00297       dttpgphbestmap = dbe->book2D("DT_TPG_phi_best_map",
00298                                    "Map of best triggers per station",20,1,21,12,0,12);
00299       dttpgphbestmapcorr = dbe->book2D("DT_TPG_phi_best_map_corr",
00300                                        "Map of correlated best triggers per station",20,1,21,12,0,12);
00301       setMapPhLabel(dttpgphmap);
00302       setMapPhLabel(dttpgphmapcorr);
00303       setMapPhLabel(dttpgphmap2nd);
00304       setMapPhLabel(dttpgphbestmap);
00305       setMapPhLabel(dttpgphbestmapcorr);
00306       
00307 
00308 
00309       dttpgthmap = dbe->book2D("DT_TPG_theta_map",
00310                                "Map of triggers per station",15,1,16,12,0,12);
00311       dttpgthmaph = dbe->book2D("DT_TPG_theta_map_h",
00312                                 "Map of H quality triggers per station",15,1,16,12,0,12);
00313       dttpgthbestmap = dbe->book2D("DT_TPG_theta_best_map",
00314                                    "Map of besttriggers per station",15,1,16,12,0,12);
00315       dttpgthbestmaph = dbe->book2D("DT_TPG_theta_best_map_h",
00316                                     "Map of H quality best triggers per station",15,1,16,12,0,12);
00317       setMapThLabel(dttpgthmap);
00318       setMapThLabel(dttpgthmaph);
00319       setMapThLabel(dttpgthbestmap);
00320       setMapThLabel(dttpgthbestmaph);
00321 
00322 
00323     }  
00324 }
00325 
00326 
00327 void L1TDTTPG::endJob(void)
00328 {
00329   if(verbose_) cout << "L1TDTTPG: end job...." << endl;
00330   LogInfo("EndJob") << "analyzed " << nev_ << " events"; 
00331 
00332   if ( outputFile_.size() != 0  && dbe ) dbe->save(outputFile_);
00333 
00334   return;
00335 }
00336 
00337 void L1TDTTPG::analyze(const Event& e, const EventSetup& c)
00338 {
00339 
00340   nev_++; 
00341   if(verbose_) cout << "L1TDTTPG: analyze...." << endl;
00342 
00343   edm::Handle<L1MuDTChambPhContainer > myL1MuDTChambPhContainer;  
00344   e.getByLabel(dttpgSource_,myL1MuDTChambPhContainer);
00345   
00346   if (!myL1MuDTChambPhContainer.isValid()) {
00347     edm::LogInfo("DataNotFound") << "can't find L1MuDTChambPhContainer with label "
00348                              << dttpgSource_.label() ;
00349     return;
00350   }
00351   L1MuDTChambPhContainer::Phi_Container *myPhContainer =  
00352     myL1MuDTChambPhContainer->getContainer();
00353 
00354   edm::Handle<L1MuDTChambThContainer > myL1MuDTChambThContainer;  
00355   e.getByLabel(dttpgSource_,myL1MuDTChambThContainer);
00356   
00357   if (!myL1MuDTChambThContainer.isValid()) {
00358     edm::LogInfo("DataNotFound") << "can't find L1MuDTChambThContainer with label "
00359                              << dttpgSource_.label() ;
00360     edm::LogInfo("DataNotFound") << "if this fails try to add DATA to the process name." ;
00361 
00362     return;
00363   }
00364   L1MuDTChambThContainer::The_Container* myThContainer =  
00365     myL1MuDTChambThContainer->getContainer();
00366 
00367   int ndttpgphtrack = 0;
00368   int ndttpgthtrack = 0; 
00369   int NumberOfSegmentsPhi[3]={0,0,0};
00370   
00371   for( L1MuDTChambPhContainer::Phi_Container::const_iterator 
00372          DTPhDigiItr =  myPhContainer->begin() ;
00373        DTPhDigiItr != myPhContainer->end() ;
00374        ++DTPhDigiItr ) 
00375     {           
00376       int bx = DTPhDigiItr->bxNum() - DTPhDigiItr->Ts2Tag();
00377       if(bx == -1)
00378        NumberOfSegmentsPhi[0]++;
00379       if(bx == 0)
00380        NumberOfSegmentsPhi[1]++;
00381       if(bx == 1)
00382        NumberOfSegmentsPhi[2]++;   
00383     }
00384    /*Fill Histos for Segment counter for each bx separately */
00385 
00386   for(int k=0;k<3;k++){
00387      dttpgphbx[k+2]->Fill(NumberOfSegmentsPhi[k]);
00388    }
00389    int bxCounterDttfPhi=0; // = no. of BX's with non-zero data
00390    for (int k=0;k<3;k++){
00391      if (NumberOfSegmentsPhi[k]>0)
00392        bxCounterDttfPhi++;
00393    }
00394 
00395    /* the BX "code" */
00396 
00397    int bxCodePhi=0;
00398    if(bxCounterDttfPhi==0){
00399      bxCodePhi=0;
00400    }else if(bxCounterDttfPhi==1){
00401      for(int k=0;k<3;k++){
00402        if(NumberOfSegmentsPhi[k]>0)
00403         bxCodePhi=k+2;
00404      }
00405    }else if(bxCounterDttfPhi==2){
00406      for(int k=0;k<3;k++){
00407        if(NumberOfSegmentsPhi[k]==0)
00408         bxCodePhi=8-k;
00409      }
00410    }else if(bxCounterDttfPhi==3){
00411      bxCodePhi=10;
00412    }
00413 
00414    //The bx analyzer histo
00415    dttpgphbx[0]->Fill(bxCodePhi);
00416 
00417 
00418    const L1MuDTChambPhDigi* bestPhQualMap[5][12][4];
00419    memset(bestPhQualMap,0,240*sizeof(L1MuDTChambPhDigi*));
00420 
00421    for( L1MuDTChambPhContainer::Phi_Container::const_iterator 
00422          DTPhDigiItr =  myPhContainer->begin() ;
00423        DTPhDigiItr != myPhContainer->end() ;
00424        ++DTPhDigiItr ) 
00425     {           
00426 
00427       ndttpgphtrack++;
00428 
00429       int bxindex = DTPhDigiItr->bxNum() - DTPhDigiItr->Ts2Tag() + 1;
00430 
00431       dttpgphwheel[bxindex]->Fill(DTPhDigiItr->whNum());
00432       if (verbose_)
00433         {
00434           cout << "DTTPG phi wheel number " << DTPhDigiItr->whNum() << endl;
00435         }
00436       dttpgphstation[bxindex]->Fill(DTPhDigiItr->stNum());
00437       if (verbose_)
00438         {   
00439           cout << "DTTPG phi station number " << DTPhDigiItr->stNum() << endl;
00440         }
00441       dttpgphsector[bxindex]->Fill(DTPhDigiItr->scNum());
00442       if (verbose_)
00443         {
00444           cout << "DTTPG phi sector number " << DTPhDigiItr->scNum() << endl;
00445         }
00446 //       dttpgphphi[bxindex]->Fill(DTPhDigiItr->phi());
00447 //       if (verbose_)
00448 //      {
00449 //        cout << "DTTPG phi phi " << DTPhDigiItr->phi() << endl;
00450 //      }
00451 //       dttpgphphiB[bxindex]->Fill(DTPhDigiItr->phiB());
00452 //       if (verbose_)
00453 //      {
00454 //        cout << "DTTPG phi phiB " << DTPhDigiItr->phiB() << endl;
00455 //      }
00456       dttpgphquality[bxindex]->Fill(DTPhDigiItr->code());
00457       if (verbose_)
00458         {
00459           cout << "DTTPG phi quality " << DTPhDigiItr->code() << endl;
00460         }
00461       dttpgphts2tag[bxindex]->Fill(DTPhDigiItr->Ts2Tag());
00462       if (verbose_)
00463         {
00464           cout << "DTTPG phi ts2tag " << DTPhDigiItr->Ts2Tag() << endl;
00465         }
00466 //       dttpgphbxcnt[bxindex]->Fill(DTPhDigiItr->BxCnt());
00467 //       if (verbose_)
00468 //      {
00469 //        cout << "DTTPG phi bxcnt " << DTPhDigiItr->BxCnt() << endl;
00470 //      }
00471     
00472       int ypos = DTPhDigiItr->scNum();
00473       int xpos = DTPhDigiItr->stNum()+4*(DTPhDigiItr->whNum()+2);
00474       dttpgphmap->Fill(xpos,ypos);
00475       if (DTPhDigiItr->Ts2Tag())
00476         dttpgphmap2nd->Fill(xpos,ypos);
00477       dttpgphmapbx[bxindex]->Fill(xpos,ypos);
00478       if (DTPhDigiItr->code()>3)
00479         dttpgphmapcorr->Fill(xpos,ypos);
00480 
00481       if (bestPhQualMap[DTPhDigiItr->whNum()+2][ DTPhDigiItr->scNum()][DTPhDigiItr->stNum()-1]==0 ||
00482           bestPhQualMap[DTPhDigiItr->whNum()+2][ DTPhDigiItr->scNum()][DTPhDigiItr->stNum()-1]->code()<DTPhDigiItr->code())
00483         {
00484           bestPhQualMap[DTPhDigiItr->whNum()+2][ DTPhDigiItr->scNum()][DTPhDigiItr->stNum()-1]=&(*DTPhDigiItr);
00485         }
00486 
00487     }
00488 
00489    for (int iwh=0; iwh<5; iwh++){
00490      for (int isec=0; isec<12; isec++){
00491        for (int ist=0; ist<4; ist++){
00492          if (bestPhQualMap[iwh][isec][ist]){
00493            int xpos = iwh*4+ist+1;
00494            dttpgphbestmap->Fill(xpos,isec);
00495            if(bestPhQualMap[iwh][isec][ist]->code()>3)
00496              dttpgphbestmapcorr->Fill(xpos,isec);
00497          }
00498        }
00499      }
00500    }
00501 
00502 
00503    int bestThQualMap[5][12][3];
00504    memset(bestThQualMap,0,180*sizeof(int));
00505    //for( vector<L1MuDTChambThDigi>::const_iterator 
00506    for( L1MuDTChambThContainer::The_Container::const_iterator 
00507          DTThDigiItr =  myThContainer->begin() ;
00508        DTThDigiItr != myThContainer->end() ;
00509        ++DTThDigiItr ) 
00510     {                           
00511       ndttpgthtrack++;
00512 
00513       int bxindex = DTThDigiItr->bxNum() + 1;
00514 
00515       dttpgthwheel[bxindex]->Fill(DTThDigiItr->whNum());
00516       if (verbose_)
00517         {
00518           cout << "DTTPG theta wheel number " << DTThDigiItr->whNum() << endl;
00519         }
00520       dttpgthstation[bxindex]->Fill(DTThDigiItr->stNum());
00521       if (verbose_)
00522         {   
00523           cout << "DTTPG theta station number " << DTThDigiItr->stNum() << endl;
00524         }
00525       dttpgthsector[bxindex]->Fill(DTThDigiItr->scNum());
00526       if (verbose_)
00527         {
00528           cout << "DTTPG theta sector number " << DTThDigiItr->scNum() << endl;
00529         }
00530       dttpgthbx[bxindex]->Fill(DTThDigiItr->bxNum());
00531       if (verbose_)
00532         {
00533           cout << "DTTPG theta bx number " << DTThDigiItr->bxNum() << endl;
00534         }
00535       int thcode[7]= {0,0,0,0,0,0,0};
00536       for (int j = 0; j < 7; j++)
00537         {
00538           dttpgththeta[bxindex]->Fill(DTThDigiItr->position(j));
00539           if (verbose_)
00540             {
00541               cout << "DTTPG theta position " << DTThDigiItr->position(j) << endl;
00542             }
00543           thcode[j]=DTThDigiItr->code(j);
00544           dttpgthquality[bxindex]->Fill(thcode[j]);
00545           if (verbose_)
00546             {
00547               cout << "DTTPG theta quality " << DTThDigiItr->code(j) << endl;
00548             }
00549         }
00550       
00551       int ypos = DTThDigiItr->scNum();
00552       int xpos = DTThDigiItr->stNum()+4*(DTThDigiItr->whNum()+2);
00553       int bestqual=0;
00554       dttpgthmap->Fill(xpos,ypos);
00555       dttpgthmapbx[bxindex]->Fill(xpos,ypos);
00556       for (int pos = 0; pos < 7; pos++){
00557         if (thcode[pos]>bestqual)
00558           bestqual=thcode[pos];
00559         if(thcode[pos]==2)
00560           dttpgthmaph->Fill(xpos,ypos);
00561       }
00562 
00563       if (bestThQualMap[DTThDigiItr->whNum()+2][ DTThDigiItr->scNum()][DTThDigiItr->stNum()-1] < bestqual)
00564         {
00565           bestThQualMap[DTThDigiItr->whNum()+2][ DTThDigiItr->scNum()][DTThDigiItr->stNum()-1]=bestqual;
00566         }
00567     }
00568 
00569    for (int iwh=0; iwh<5; iwh++){
00570      for (int isec=0; isec<12; isec++){
00571        for (int ist=0; ist<3; ist++){
00572          if (bestThQualMap[iwh][isec][ist]){
00573            int xpos = iwh*4+ist+1;
00574            dttpgthbestmap->Fill(xpos,isec);
00575            if(bestThQualMap[iwh][isec][ist]==2)
00576              dttpgthbestmaph->Fill(xpos,isec);
00577          }
00578        }
00579      }
00580    }
00581 
00582 
00583   dttpgphntrack->Fill(ndttpgphtrack);
00584   if (verbose_)
00585     {
00586       cout << "DTTPG phi ntrack " << ndttpgphtrack << endl;
00587     }
00588   dttpgthntrack->Fill(ndttpgthtrack);
00589   if (verbose_) {
00590     cout << "DTTPG theta ntrack " << ndttpgthtrack << endl;
00591   }
00592 
00593   edm::Handle<L1MuDTTrackContainer > myL1MuDTTrackContainer;
00594 
00595   
00596     std::string trstring;
00597     trstring = dttpgSource_.label()+":"+"DATA"+":"+dttpgSource_.process();
00598     edm::InputTag trInputTag(trstring);
00599     e.getByLabel(trInputTag,myL1MuDTTrackContainer);
00600 
00601   if (!myL1MuDTTrackContainer.isValid()) {
00602     edm::LogInfo("DataNotFound") << "can't find L1MuDTTrackContainer with label "
00603                                << dttpgSource_.label() ;
00604     return;
00605   }
00606 
00607   L1MuDTTrackContainer::TrackContainer *t =  myL1MuDTTrackContainer->getContainer();
00608 
00609 
00610 
00611   int NumberOfSegmentsOut[3]={0,0,0};
00612   for ( L1MuDTTrackContainer::TrackContainer::const_iterator i 
00613           = t->begin(); i != t->end(); ++i ) {
00614     if(i->bx() ==-1)
00615        NumberOfSegmentsOut[0]++;
00616     if(i->bx() ==0)
00617        NumberOfSegmentsOut[1]++;
00618     if(i->bx() ==1)
00619        NumberOfSegmentsOut[2]++;
00620   }
00621 
00622 
00623    /*Fill Histos for Segment counter*/
00624    for(int k=0;k<3;k++){
00625      dttpgphbx[k+5]->Fill(NumberOfSegmentsOut[k]);
00626    }
00627 
00628    /*Bunch assigments*/
00629 
00630    int bxCounterDttfOut=0;
00631    for (int k=0;k<3;k++){
00632    if (NumberOfSegmentsOut[k]>0)
00633      bxCounterDttfOut++;
00634    }
00635 
00636    int bxCodeOut=0;
00637    if(bxCounterDttfOut==0){
00638      bxCodeOut=0;
00639    }else if(bxCounterDttfOut==1){
00640      for(int k=0;k<3;k++){
00641        if(NumberOfSegmentsOut[k]>0)
00642          bxCodeOut=k+2;
00643      }
00644    }else if(bxCounterDttfOut==2){
00645      for(int k=0;k<3;k++){
00646        if(NumberOfSegmentsOut[k]==0)
00647          bxCodeOut=8-k;
00648      }
00649    }else if(bxCounterDttfOut==3){
00650      bxCodeOut=10;
00651    }
00652 
00653    //The bx analyzer histo
00654    dttpgphbx[1]->Fill(bxCodeOut);
00655 
00656    /*End Dttf Output Bunch analysis*/
00657 
00658    // the 2-DIM histo with phi.input vs. output
00659    dttpgphbxcomp->Fill(bxCodePhi,bxCodeOut);
00660 
00661 
00662   for ( L1MuDTTrackContainer::TrackContainer::const_iterator i 
00663           = t->begin(); i != t->end(); ++i ) {
00664     if ( verbose_ ) {
00665       std::cout << "bx = " << i->bx() 
00666                 << std::endl;
00667       std::cout << "quality (packed) = " << i->quality_packed() 
00668                 << std::endl;
00669       std::cout << "pt      (packed) = " << i->pt_packed() 
00670                 << std::endl;
00671       std::cout << "phi     (packed) = " << i->phi_packed() 
00672                 << std::endl;
00673       std::cout << "charge  (packed) = " << i->charge_packed() 
00674                 << std::endl;
00675     }
00676 
00677 
00678     int bxindex = i->bx() + 1;
00679     dttf_p_phi[bxindex]->Fill(i->phi_packed());
00680     dttf_p_qual[bxindex]->Fill(i->quality_packed());
00681     dttf_p_pt[bxindex]->Fill(i->pt_packed());
00682     dttf_p_q[bxindex]->Fill(i->charge_packed());
00683   }
00684     
00685 }
00686 
00687 void L1TDTTPG::setMapPhLabel(MonitorElement *me)
00688 {
00689 
00690   me->setAxisTitle("DTTF Sector",2);
00691       for(int i=0;i<5;i++){
00692         ostringstream wheel;
00693         wheel << i-2;
00694         me->setBinLabel(1+i*4,"Wheel "+ wheel.str(),1);
00695       }
00696   
00697 }
00698 
00699 void L1TDTTPG::setMapThLabel(MonitorElement *me)
00700 {
00701 
00702   me->setAxisTitle("DTTF Sector",2);
00703       for(int i=0;i<5;i++){
00704         ostringstream wheel;
00705         wheel << i-2;
00706         me->setBinLabel(1+i*3,"Wheel "+ wheel.str(),1);
00707       }
00708   
00709 }