CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_9/src/DQM/TrigXMonitor/src/L1Scalers.cc

Go to the documentation of this file.
00001 // $Id: L1Scalers.cc,v 1.27 2011/12/31 12:26:03 eulisse Exp $
00002 #include <iostream>
00003 
00004 
00005 // FW
00006 #include "FWCore/Framework/interface/Event.h"
00007 
00008 #include "FWCore/ServiceRegistry/interface/Service.h"
00009 
00010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00011 
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 
00014 // L1
00015 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00016 
00017 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuRegionalCand.h"
00018 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTCand.h"
00019 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTExtendedCand.h"
00020 #include "DataFormats/L1GlobalMuonTrigger/interface/L1MuGMTReadoutCollection.h"
00021 
00022 
00023 #include "FWCore/Framework/interface/LuminosityBlock.h"
00024 
00025 
00026 #include "DataFormats/Scalers/interface/L1TriggerScalers.h"
00027 #include "DataFormats/Scalers/interface/L1TriggerRates.h"
00028 #include "DataFormats/Scalers/interface/LumiScalers.h"
00029 #include "DQM/TrigXMonitor/interface/L1Scalers.h"
00030 #include "DataFormats/Common/interface/Handle.h"
00031 #include "DQMServices/Core/interface/DQMStore.h"
00032 
00033 // HACK START
00034 #include "DataFormats/FEDRawData/interface/FEDRawDataCollection.h"
00035 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00036 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00037 // HACK END
00038 
00039 using namespace edm;
00040 
00041 
00042 
00043 L1Scalers::L1Scalers(const edm::ParameterSet &ps):
00044   dbe_(0), nev_(0),
00045   verbose_(ps.getUntrackedParameter < bool > ("verbose", false)),
00046   l1GtDataSource_( ps.getParameter< edm::InputTag >("l1GtData")),
00047   denomIsTech_(ps.getUntrackedParameter <bool> ("denomIsTech", true)),
00048   denomBit_(ps.getUntrackedParameter <unsigned int> ("denomBit", 40)),
00049   tfIsTech_(ps.getUntrackedParameter <bool> ("tfIsTech", true)),
00050   tfBit_(ps.getUntrackedParameter <unsigned int> ("tfBit", 41)),
00051   algoSelected_(ps.getUntrackedParameter<std::vector<unsigned int> >("algoMonitorBits", std::vector<unsigned int>())),
00052   techSelected_(ps.getUntrackedParameter<std::vector<unsigned int> >("techMonitorBits", std::vector<unsigned int>())),
00053   folderName_( ps.getUntrackedParameter< std::string>("dqmFolder", 
00054                                           std::string("L1T/L1Scalers_EvF"))),
00055   l1scalers_(0),
00056   l1techScalers_(0),
00057   l1Correlations_(0),
00058   bxNum_(0),
00059   l1scalersBx_(0),
00060   l1techScalersBx_(0),
00061 //   pixFedSizeBx_(0),
00062 //   hfEnergyMaxTowerBx_(0),
00063   nLumiBlock_(0),
00064   l1AlgoCounter_(0),
00065   l1TtCounter_(0),
00066 //   pixFedSize_(0),
00067 //   hfEnergy_(0),
00068   fedStart_(ps.getUntrackedParameter<unsigned int>("firstFED", 0)),
00069   fedStop_(ps.getUntrackedParameter<unsigned int>("lastFED", 931)), 
00070   rateAlgoCounter_(0),                                    
00071   rateTtCounter_(0),                                      
00072   fedRawCollection_(ps.getParameter<edm::InputTag>("fedRawData")),
00073   maskedList_(ps.getUntrackedParameter<std::vector<int> >("maskedChannels", std::vector<int>())), //this is using the ashed index
00074   HcalRecHitCollection_(ps.getParameter<edm::InputTag>("HFRecHitCollection"))
00075 {
00076   LogDebug("Status") << "constructor" ;
00077 } 
00078 
00079 
00080 
00081 
00082 void L1Scalers::beginJob(void)
00083 {
00084   LogDebug("Status") << "L1Scalers::beginJob()...";
00085 
00086   dbe_ = Service<DQMStore>().operator->();
00087   if (dbe_ ) {
00088     dbe_->setVerbose(0);
00089     dbe_->setCurrentFolder(folderName_);
00090 
00091 
00092     l1scalers_ = dbe_->book1D("l1AlgoBits", "L1 Algorithm Bits",
00093                               128, -0.5, 127.5);
00094     l1scalersBx_ = dbe_->book2D("l1AlgoBits_Vs_Bx", "L1 Algorithm Bits vs "
00095                                 "Bunch Number",
00096                                 3600, -0.5, 3599.5,
00097                                 128, -0.5, 127.5);
00098     l1Correlations_ = dbe_->book2D("l1Correlations", "L1 Algorithm Bits " 
00099                                     "Correlations",
00100                                    128, -0.5, 127.5,
00101                                    128, -0.5, 127.5);
00102     l1techScalers_ = dbe_->book1D("l1TechBits", "L1 Tech. Trigger Bits",
00103                                   64, -0.5, 63.5);
00104     l1techScalersBx_ = dbe_->book2D("l1TechBits_Vs_Bx", "L1 Technical "
00105                                     "Trigger "
00106                                     "Bits vs Bunch Number",
00107                                     3600, -0.5, 3599.5, 64, -0.5, 63.5);
00108 //     pixFedSizeBx_ = dbe_->book2D("pixFedSize_Vs_Bx", "Size of Pixel FED data vs "
00109 //                              "Bunch Number",
00110 //                              3600, -0.5, 3599.5,
00111 //                              200, 0., 20000.);
00112 //     hfEnergyMaxTowerBx_ = dbe_->book2D("hfEnergyMaxTower_Vs_Bx", "HF Energy Max Tower vs "
00113 //                              "Bunch Number",
00114 //                              3600, -0.5, 3599.5,
00115 //                              100, 0., 500.);
00116     bxNum_ = dbe_->book1D("bxNum", "Bunch number from GTFE",
00117                           3600, -0.5, 3599.5);
00118 
00119     nLumiBlock_ = dbe_->bookInt("nLumiBlock");
00120 
00121 //  l1 total rate
00122     
00123     l1AlgoCounter_ = dbe_->bookInt("l1AlgoCounter");
00124     l1TtCounter_ = dbe_->bookInt("l1TtCounter");
00125 
00126 
00127     //int maxNbins = 200;
00128     //int maxLumi = 2000;
00129 
00130     //timing plots
00131     std::stringstream sdenom;
00132     if(denomIsTech_) sdenom << "tech" ;
00133     else sdenom << "algo" ;
00134 
00135     dbe_->setCurrentFolder(folderName_ + "/Synch");
00136     algoBxDiff_.clear();
00137     algoBxDiff_.clear();
00138     algoBxDiffLumi_.clear();
00139     techBxDiffLumi_.clear();
00140     for(uint ibit = 0; ibit < algoSelected_.size(); ibit++){
00141       std::stringstream ss;
00142       ss << algoSelected_[ibit] << "_" << sdenom.str() << denomBit_;
00143       algoBxDiff_.push_back(dbe_->book1D("BX_diff_algo"+ ss.str(),"BX_diff_algo"+ ss.str(),9,-4,5));
00144       algoBxDiffLumi_.push_back(dbe_->book2D("BX_diffvslumi_algo"+ ss.str(),"BX_diff_algo"+ss.str(),MAX_LUMI_BIN,-0.5,double(MAX_LUMI_SEG)-0.5,9,-4,5));
00145       //algoBxDiffLumi_[ibit]->setAxisTitle("Lumi Section", 1);
00146     }
00147     for(uint ibit = 0; ibit < techSelected_.size(); ibit++){
00148       std::stringstream ss;
00149       ss << techSelected_[ibit] << "_" << sdenom.str() << denomBit_;
00150       techBxDiff_.push_back(dbe_->book1D("BX_diff_tech"+ ss.str(),"BX_diff_tech"+ ss.str(),9,-4,5));
00151       techBxDiffLumi_.push_back(dbe_->book2D("BX_diffvslumi_tech"+ ss.str(),"BX_diff_tech"+ss.str(),MAX_LUMI_BIN,-0.5,double(MAX_LUMI_SEG)-0.5,9,-4,5));
00152       //techBxDiffLumi_[ibit]->setAxisTitle("Lumi Section", 1);
00153     }
00154 
00155     //GMT timing plots
00156     std::stringstream ss1;
00157     ss1 << "_" << sdenom.str() << denomBit_;
00158     dtBxDiff_ = dbe_->book1D("BX_diff_DT" + ss1.str(),"BX_diff_DT" + ss1.str(),9,-4,5);
00159     dtBxDiffLumi_ = dbe_->book2D("BX_diffvslumi_DT" + ss1.str(),"BX_diffvslumi_DT" + ss1.str(),MAX_LUMI_BIN,-0.5,double(MAX_LUMI_SEG)-0.5,9,-4,5);
00160     cscBxDiff_ = dbe_->book1D("BX_diff_CSC" + ss1.str(),"BX_diff_CSC" + ss1.str(),9,-4,5);
00161     cscBxDiffLumi_ = dbe_->book2D("BX_diffvslumi_CSC" + ss1.str(),"BX_diffvslumi_CSC" + ss1.str(),MAX_LUMI_BIN,-0.5,double(MAX_LUMI_SEG)-0.5,9,-4,5);
00162     rpcbBxDiff_ = dbe_->book1D("BX_diff_RPCb" + ss1.str(),"BX_diff_RPCb" + ss1.str(),9,-4,5);
00163     rpcbBxDiffLumi_ = dbe_->book2D("BX_diffvslumi_RPCb" + ss1.str(),"BX_diffvslumi_RPCb" + ss1.str(),MAX_LUMI_BIN,-0.5,double(MAX_LUMI_SEG)-0.5,9,-4,5);
00164     rpcfBxDiff_ = dbe_->book1D("BX_diff_RPCf" + ss1.str(),"BX_diff_RPCf" + ss1.str(),9,-4,5);
00165     rpcfBxDiffLumi_ = dbe_->book2D("BX_diffvslumi_RPCf" + ss1.str(),"BX_diffvslumi_RPCf" + ss1.str(),MAX_LUMI_BIN,-0.5,double(MAX_LUMI_SEG)-0.5,9,-4,5);
00166 
00167 
00168     // early triggers
00169 //     pixFedSize_ = dbe_->book1D("pixFedSize", "Size of Pixel FED data",
00170 //                             200, 0., 20000.);
00171 //     hfEnergy_   = dbe_->book1D("hfEnergy", "HF energy",
00172 //                             100, 0., 500.);
00173 
00174   }
00175   
00176   
00177   return;
00178 }
00179 
00180 void L1Scalers::endJob(void)
00181 {
00182 }
00183 
00184 void L1Scalers::analyze(const edm::Event &e, const edm::EventSetup &iSetup)
00185 {
00186   nev_++;
00187 
00188   LogDebug("Status") << "L1Scalers::analyze  event " << nev_ ;
00189 
00190   // int myGTFEbx = -1;
00191   // get Global Trigger decision and the decision word
00192   // these are locally derived
00193   edm::Handle<L1GlobalTriggerReadoutRecord> gtRecord;
00194   bool t = e.getByLabel(l1GtDataSource_,gtRecord);
00195 
00196   if ( ! t ) {
00197     LogDebug("Product") << "can't find L1GlobalTriggerReadoutRecord "
00198                         << "with label " << l1GtDataSource_.label() ;
00199   }
00200   else {
00201 
00202     // DEBUG
00203     //gtRecord->print(std::cout);
00204     // DEBUG
00205     
00206     L1GtfeWord gtfeWord = gtRecord->gtfeWord();
00207     int gtfeBx = gtfeWord.bxNr();
00208     bxNum_->Fill(gtfeBx);
00209     // myGTFEbx = gtfeBx;
00210     
00211     bool tfBitGood = false;
00212 
00213     // First, the default
00214     // vector of bool
00215     for(int iebx=0; iebx<=4; iebx++) {
00216 
00217       //Algorithm Bits
00218       DecisionWord gtDecisionWord = gtRecord->decisionWord(iebx-2);
00219       //    DecisionWord gtDecisionWord = gtRecord->decisionWord();
00220       if ( ! gtDecisionWord.empty() ) { // if board not there this is zero
00221         // loop over dec. bit to get total rate (no overlap)
00222         for ( uint i = 0; i < gtDecisionWord.size(); ++i ) {
00223           if ( gtDecisionWord[i] ) {
00224             rateAlgoCounter_++;
00225             l1AlgoCounter_->Fill(rateAlgoCounter_);
00226             break;
00227           }
00228         }
00229         // loop over decision bits
00230         for ( uint i = 0; i < gtDecisionWord.size(); ++i ) {
00231           if ( gtDecisionWord[i] ) {
00232             l1scalers_->Fill(i);
00233             l1scalersBx_->Fill(gtfeBx-2+iebx,i);
00234             for ( uint j = i + 1; j < gtDecisionWord.size(); ++j ) {
00235               if ( gtDecisionWord[j] ) {
00236                 l1Correlations_->Fill(i,j);
00237                 l1Correlations_->Fill(j,i);
00238               }
00239             }
00240           }
00241         }
00242       }
00243  
00244 
00245       // loop over technical triggers
00246       // vector of bool again. 
00247       TechnicalTriggerWord tw = gtRecord->technicalTriggerWord(iebx-2);
00248       //    TechnicalTriggerWord tw = gtRecord->technicalTriggerWord();
00249       if ( ! tw.empty() ) {
00250         // loop over dec. bit to get total rate (no overlap)
00251         for ( uint i = 0; i < tw.size(); ++i ) {
00252           if ( tw[i] ) {
00253             rateTtCounter_++;
00254             l1TtCounter_->Fill(rateTtCounter_);
00255             break;
00256           }
00257         }        
00258         for ( uint i = 0; i < tw.size(); ++i ) {
00259           if ( tw[i] ) {
00260             l1techScalers_->Fill(i);
00261             l1techScalersBx_->Fill(gtfeBx-2+iebx, i);
00262           }
00263         } 
00264 
00265         // check if bit used to filter timing plots fired in this event 
00266         // (anywhere in the bx window)
00267         if(tfIsTech_){
00268           if(tfBit_ < tw.size()){
00269             if( tw[tfBit_] ) tfBitGood = true;
00270           }
00271         }
00272       } // ! tw.empty
00273 
00274     }//bx
00275 
00276 
00277     //timing plots
00278     earliestDenom_ = 9;
00279     earliestAlgo_.clear();
00280     earliestTech_.clear();
00281     for(uint i=0; i < techSelected_.size(); i++) earliestTech_.push_back(9);
00282     for(uint i=0; i < algoSelected_.size(); i++) earliestAlgo_.push_back(9);
00283 
00284     //GMT information
00285     edm::Handle<L1MuGMTReadoutCollection> gmtCollection;
00286     e.getByLabel(l1GtDataSource_,gmtCollection);
00287     
00288 
00289     if (!gmtCollection.isValid()) {
00290       edm::LogInfo("DataNotFound") << "can't find L1MuGMTReadoutCollection with label "
00291                                    << l1GtDataSource_.label() ;
00292     }
00293 
00294     // remember the bx of 1st candidate of each system (9=none)
00295     int bx1st[4] = {9, 9, 9, 9};
00296       
00297     if(tfBitGood){//to avoid single BSC hits
00298 
00299       for(int iebx=0; iebx<=4; iebx++) {
00300         TechnicalTriggerWord tw = gtRecord->technicalTriggerWord(iebx-2);
00301         DecisionWord gtDecisionWord = gtRecord->decisionWord(iebx-2);
00302 
00303         bool denomBitGood = false;
00304 
00305         //check if reference bit is valid
00306         if(denomIsTech_){
00307           if ( ! tw.empty() ) {
00308             if( denomBit_ < tw.size() ){
00309               denomBitGood = true;
00310               if( tw[denomBit_] && earliestDenom_==9 ) earliestDenom_ = iebx; 
00311             }
00312           }
00313         }
00314         else{
00315           if ( ! gtDecisionWord.empty() ) { 
00316             if( denomBit_ < gtDecisionWord.size() ){
00317               denomBitGood = true;
00318               if( gtDecisionWord[denomBit_] && earliestDenom_==9 ) earliestDenom_ = iebx; 
00319             }
00320           }
00321         }
00322 
00323         if(denomBitGood){
00324 
00325           //get earliest tech bx's
00326           if ( ! tw.empty() ) {
00327             for(uint ibit = 0; ibit < techSelected_.size(); ibit++){      
00328               if(techSelected_[ibit] < tw.size()){
00329                 if(tw[techSelected_[ibit]] && earliestTech_[ibit] == 9) earliestTech_[ibit] = iebx;
00330               }
00331             }
00332           }
00333 
00334           //get earliest algo bx's
00335           if(!gtDecisionWord.empty()){      
00336             for(uint ibit = 0; ibit < algoSelected_.size(); ibit++){      
00337               if(algoSelected_[ibit] < gtDecisionWord.size()){
00338                 if(gtDecisionWord[algoSelected_[ibit]] && earliestAlgo_[ibit] == 9) earliestAlgo_[ibit] = iebx;
00339               }
00340             }
00341           }
00342 
00343         }//denomBitGood
00344 
00345       }//bx
00346 
00347 
00348       //get earliest single muon trigger system bx's
00349       if (gmtCollection.isValid()) {
00350 
00351         // get GMT readout collection
00352         L1MuGMTReadoutCollection const* gmtrc = gmtCollection.product();
00353         // get record vector
00354         std::vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords();
00355         // loop over records of individual bx's
00356         std::vector<L1MuGMTReadoutRecord>::const_iterator RRItr;
00357 
00358         for( RRItr = gmt_records.begin(); RRItr != gmt_records.end(); RRItr++ ) {//loop from BX=-2 to BX=2
00359           std::vector<L1MuRegionalCand> INPCands[4] = {
00360             RRItr->getDTBXCands(),
00361             RRItr->getBrlRPCCands(),
00362             RRItr->getCSCCands(),
00363             RRItr->getFwdRPCCands()
00364           };
00365           std::vector<L1MuRegionalCand>::const_iterator INPItr;
00366           int BxInEvent = RRItr->getBxInEvent();
00367           
00368           // find the first non-empty candidate in this bx
00369           for(int i=0; i<4; i++) {//for each single muon trigger system
00370             for( INPItr = INPCands[i].begin(); INPItr != INPCands[i].end(); ++INPItr ) {
00371               if(!INPItr->empty()) {
00372                 if(bx1st[i]==9) bx1st[i]=BxInEvent+2;//must go from 0 to 4 (consistent with above)
00373               }
00374             }      
00375           }
00376           //for(int i=0; i<4; i++) 
00377           //    std::cout << "bx1st[" << i << "] = " << bx1st[i];
00378           //std::cout << std::endl;
00379         }
00380 
00381       }//gmtCollection.isValid
00382 
00383 
00384       //calculated bx difference
00385       if(earliestDenom_ != 9){
00386         for(uint ibit = 0; ibit < techSelected_.size(); ibit++){          
00387           if(earliestTech_[ibit] != 9){
00388             int diff = earliestTech_[ibit] - earliestDenom_ ;
00389             techBxDiff_[ibit]->Fill(diff);
00390             techBxDiffLumi_[ibit]->Fill(e.luminosityBlock(), diff);
00391           }
00392         }
00393         for(uint ibit = 0; ibit < algoSelected_.size(); ibit++){          
00394           if(earliestAlgo_[ibit] != 9){
00395             int diff = earliestAlgo_[ibit] - earliestDenom_ ;
00396             algoBxDiff_[ibit]->Fill(diff);
00397             algoBxDiffLumi_[ibit]->Fill(e.luminosityBlock(), diff);
00398           }
00399         }
00400 
00401         if(bx1st[0] != 9){
00402           int diff = bx1st[0] - earliestDenom_;
00403           dtBxDiff_->Fill(diff);
00404           dtBxDiffLumi_->Fill(e.luminosityBlock(), diff);
00405         }
00406         if(bx1st[1] != 9){
00407           int diff = bx1st[1] - earliestDenom_;
00408           rpcbBxDiff_->Fill(diff);
00409           rpcbBxDiffLumi_->Fill(e.luminosityBlock(), diff);
00410         }
00411         if(bx1st[2] != 9){
00412           int diff = bx1st[2] - earliestDenom_;
00413           cscBxDiff_->Fill(diff);
00414           cscBxDiffLumi_->Fill(e.luminosityBlock(), diff);
00415         }
00416         if(bx1st[3] != 9){
00417           int diff = bx1st[3] - earliestDenom_;
00418           rpcfBxDiff_->Fill(diff);
00419           rpcfBxDiffLumi_->Fill(e.luminosityBlock(), diff);
00420         }
00421 
00422       }
00423 
00424     }//tt41Good
00425 
00426   } // getbylabel succeeded
00427   
00428 
00429 //   // HACK
00430 //   // getting very basic uncalRH
00431 //   edm::Handle<FEDRawDataCollection> theRaw;
00432 //   bool getFed = e.getByLabel(fedRawCollection_, theRaw);
00433 //   if ( ! getFed ) {
00434 //     edm::LogInfo("FEDSizeFilter") << fedRawCollection_ << " not available";
00435 //   }
00436 //   else { // got the fed raw data
00437 //     unsigned int totalFEDsize = 0 ; 
00438 //     for (unsigned int i=fedStart_; i<=fedStop_; ++i) {
00439 //       LogDebug("Parameter") << "Examining fed " << i << " with size "
00440 //                          << theRaw->FEDData(i).size() ;
00441 //       totalFEDsize += theRaw->FEDData(i).size() ; 
00442 //     }
00443 //     pixFedSize_->Fill(totalFEDsize);
00444 //     if( (myGTFEbx!=-1) ) pixFedSizeBx_->Fill(myGTFEbx,totalFEDsize);
00445 
00446 //     LogDebug("Parameter") << "Total FED size: " << totalFEDsize;
00447 //   }      
00448 
00449 //   // HF - stolen from HLTrigger/special
00450 //   // getting very basic uncalRH
00451 //   double maxHFenergy = -1;
00452 //   edm::Handle<HFRecHitCollection> crudeHits;
00453 //   bool getHF = e.getByLabel(HcalRecHitCollection_, crudeHits);
00454 //   if ( ! getHF ) {
00455 //     LogDebug("Status") << HcalRecHitCollection_ << " not available";
00456 //   }
00457 //   else {
00458 
00459 //     LogDebug("Status") << "Filtering, with " << crudeHits->size() 
00460 //                     << " recHits to consider" ;
00461 //     for ( HFRecHitCollection::const_iterator hitItr = crudeHits->begin(); 
00462 //        hitItr != crudeHits->end(); ++hitItr ) {     
00463 //       HFRecHit hit = (*hitItr);
00464      
00465 //       // masking noisy channels
00466 //       std::vector<int>::iterator result;
00467 //       result = std::find( maskedList_.begin(), maskedList_.end(), 
00468 //                        HcalDetId(hit.id()).hashed_index() );    
00469 //       if  (result != maskedList_.end()) 
00470 //      continue; 
00471 //       hfEnergy_->Fill(hit.energy());
00472 //       if( (hit.energy()>maxHFenergy) ) maxHFenergy = hit.energy();
00473 //     }
00474 //   }
00475 
00476 //   if( (maxHFenergy!=-1 && myGTFEbx!=-1) ) hfEnergyMaxTowerBx_->Fill(myGTFEbx,maxHFenergy);
00477 //   // END HACK
00478 
00479   return;
00480  
00481 }
00482 
00483 void L1Scalers::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg, 
00484                                     const edm::EventSetup& iSetup)
00485 {
00486   nLumiBlock_->Fill(lumiSeg.id().luminosityBlock());
00487 
00488 }
00489 
00490 
00492 void L1Scalers::beginRun(const edm::Run& run, const edm::EventSetup& iSetup)
00493 {
00494 }
00495 
00497 void L1Scalers::endRun(const edm::Run& run, const edm::EventSetup& iSetup)
00498 {
00499 }
00500 
00501