CMS 3D CMS Logo

BxTiming Class Reference

#include <DQM/L1TMonitor/interface/BxTiming.h>

Inheritance diagram for BxTiming:

edm::EDAnalyzer

List of all members.

Public Member Functions

 BxTiming (const edm::ParameterSet &)
 ~BxTiming ()

Protected Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (const edm::EventSetup &)
virtual void endJob ()

Private Types

enum  nsys { NSYS = 9 }
enum  syslist {
  ETP = 0, HTP, GCT, CTP,
  CTF, DTP, DTF, RPC,
  GLT
}

Private Member Functions

int calcBxDiff (int bx1, int bx2)
 calculates the difference (closest distance) between two bunch crossing numbers.
int verbose ()

Private Attributes

DQMStoredbe
std::pair< int, intfedRange_ [NSYS]
int fedRef_
edm::InputTag fedSource_
edm::InputTag gtSource_
MonitorElementhBxDiffAllFed
 histograms
MonitorElementhBxDiffAllFedSpread [nspr_]
MonitorElementhBxDiffSysFed [NSYS]
MonitorElementhBxOccyAllFed
MonitorElementhBxOccyAllFedSpread [nspr_]
MonitorElementhBxOccyGtTrigType [nttype_]
MonitorElement ** hBxOccyOneFed
MonitorElement ** hBxOccyTrigBit [NSYS]
std::string histFile_
std::string histFolder_
std::vector< intlistGtBits_
int nBxDiff [1500][nspr_]
int nBxOccy [1500][nspr_]
int nEvt_
int nfed_
bool runInFF_
int verbose_

Static Private Attributes

static const int half_norb_ = norb_ / 2
static const int nbig_ = 10000
static const int norb_ = 3564
static const int nspr_ = 3
static const int nttype_ = 6


Detailed Description

Definition at line 25 of file BxTiming.h.


Member Enumeration Documentation

enum BxTiming::nsys [private]

Enumerator:
NSYS 

Definition at line 80 of file BxTiming.h.

00080 {NSYS=9};

enum BxTiming::syslist [private]

Enumerator:
ETP 
HTP 
GCT 
CTP 
CTF 
DTP 
DTF 
RPC 
GLT 

Definition at line 81 of file BxTiming.h.

00081 {ETP=0, HTP, GCT, CTP, CTF, DTP, DTF, RPC, GLT};


Constructor & Destructor Documentation

BxTiming::BxTiming ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 10 of file BxTiming.cc.

References GenMuonPlsPt100GeV_cfg::cout, dbe, fedRef_, fedSource_, flush(), edm::ParameterSet::getUntrackedParameter(), gtSource_, histFile_, histFolder_, i, FEDNumbering::lastFEDId(), listGtBits_, nEvt_, nfed_, NULL, runInFF_, DQMStore::setCurrentFolder(), DQMStore::setVerbose(), verbose(), and verbose_.

00010                                                  {
00011 
00012   verbose_ = iConfig.getUntrackedParameter<int>("VerboseFlag",0);
00013   if(verbose())
00014     std::cout << "BxTiming::BxTiming()...\n" << std::flush;
00015 
00016   fedRef_ = iConfig.getUntrackedParameter<int>("ReferenceFedId",813);
00017   fedSource_ = iConfig.getUntrackedParameter<edm::InputTag>
00018     ("FedSource",edm::InputTag("source"));
00019   gtSource_ = iConfig.getUntrackedParameter<edm::InputTag>
00020     ("GtSource",edm::InputTag("gtUnpack"));
00021   histFile_ = iConfig.getUntrackedParameter<std::string>
00022     ("HistFile","");
00023   histFolder_ = iConfig.getUntrackedParameter<std::string>
00024     ("HistFolder", "L1T/BXSynch/");
00025 
00026   runInFF_ = iConfig.getUntrackedParameter<bool> ("RunInFilterFarm", true);
00027  if(runInFF_) histFolder_ = "L1T/BXSynch_EvF/";
00028   if(verbose())
00029     std::cout << "Filter farm run setting?" << runInFF_
00030               << "\n" << std::flush;
00031 
00032   listGtBits_ = iConfig.getUntrackedParameter<std::vector<int> > ("GtBitList", std::vector<int>(1,0));
00033   if(listGtBits_.size()==1 && listGtBits_.at(0)==-1) {
00034     int ngtbits = 128;
00035     listGtBits_.reserve(ngtbits); 
00036     for(int i=0; i<ngtbits; i++) 
00037       listGtBits_[i]=i;
00038   }
00039 
00040   if(verbose()) {
00041     std::cout << "BxTiming: gt bits set for timing dqm:";
00042     std::cout << "nbits:" << listGtBits_.size() << " list: " ;
00043     for(size_t i=0; i!=listGtBits_.size(); i++) 
00044       std::cout << listGtBits_.at(i) << " " ;
00045     std::cout << "\n" << std::flush;
00046   }
00047 
00048   nfed_ = FEDNumbering::lastFEDId()+1;
00049 
00050   dbe = NULL;
00051   if (iConfig.getUntrackedParameter<bool>("DQMStore", false)) { 
00052     dbe = edm::Service<DQMStore>().operator->();
00053     dbe->setVerbose(0);
00054   }
00055   
00056   if(dbe!=NULL)
00057     dbe->setCurrentFolder(histFolder_);
00058   
00059   nEvt_ = 0;
00060   
00061   if(verbose())
00062     std::cout << "BxTiming::BxTiming constructor...done.\n" << std::flush;
00063 }

BxTiming::~BxTiming (  ) 

Definition at line 65 of file BxTiming.cc.

00065 {}


Member Function Documentation

void BxTiming::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [protected, virtual]

get the raw data

Implements edm::EDAnalyzer.

Definition at line 254 of file BxTiming.cc.

References calcBxDiff(), GenMuonPlsPt100GeV_cfg::cout, FEDRawData::data(), data, lat::endl(), fedRange_, fedRef_, fedSource_, MonitorElement::Fill(), first, flush(), edm::Event::getByLabel(), GLT, gtSource_, hBxDiffAllFed, hBxDiffAllFedSpread, hBxDiffSysFed, hBxOccyAllFed, hBxOccyAllFedSpread, hBxOccyGtTrigType, hBxOccyOneFed, hBxOccyTrigBit, header, i, int, edm::Handle< T >::isValid(), j, k, FEDNumbering::lastFEDId(), listGtBits_, nBxDiff, nBxOccy, nEvt_, nfed_, nspr_, NSYS, nttype_, runInFF_, edm::second(), MonitorElement::setBinContent(), FEDRawData::size(), size, and verbose().

00254                                                                      {
00255   
00256   if(verbose())
00257     std::cout << "BxTiming::analyze()  start\n" << std::flush;
00258 
00259   nEvt_++;
00260 
00262   edm::Handle<FEDRawDataCollection> rawdata;
00263   iEvent.getByLabel(fedSource_, rawdata);
00264   //iEvent.getByType(rawdata);
00265 
00266   // get the GT bits
00267   edm::Handle<L1GlobalTriggerReadoutRecord> gtdata;
00268   iEvent.getByLabel(gtSource_, gtdata);
00269   std::vector<bool> gtbits;
00270   int ngtbits = 128;
00271   gtbits.reserve(ngtbits); for(int i=0; i<ngtbits; i++) gtbits[i]=false;
00272   if(gtdata.isValid())
00273     gtbits = gtdata->decisionWord();
00274   
00275   if(gtbits.size()==0) {
00276     gtbits.push_back(true); // gtdata->decision();
00277     if(verbose())
00278       std::cout << "BxTiming::analyze() | unexpected empty decision bits!";
00279   }
00280 
00281   if(verbose()) {
00282     std::cout << "BxTiming::analyze()  gt data valid:" << (int)(gtdata.isValid()?0:1)
00283               << " decision word size:" << (int)(gtbits.size()) << "  bits: ";
00284     for(size_t i=0; i!=gtbits.size(); i++) {
00285       int ii = gtbits.at(i)? 1:0;
00286       std::cout << ii;
00287     }
00288     std::cout << ".\n" << std::flush;
00289   }
00290 
00291 
00292   // get reference bx
00293   int bxRef = FEDHeader(rawdata->FEDData(fedRef_).data()).bxID();
00294 
00295   // triggerType
00296   //trigger types: physics (1), calibration (2), random (3), traced physics (5),  test (6) 
00297   int ttype = FEDHeader(rawdata->FEDData(812).data()).triggerType();
00298 
00299   // loop over feds
00300   for (int i = 0; i<FEDNumbering::lastFEDId(); i++){
00301     const FEDRawData& data = rawdata->FEDData(i);
00302     size_t size=data.size();
00303     if(!size) continue;
00304     FEDHeader header(data.data());
00305     //int lvl1id = header.lvl1ID();//Level-1 event number generated by the TTC system
00306     int bx = header.bxID(); // The bunch crossing number
00307 
00308     int bxDiff = calcBxDiff(bx,bxRef); // deviation from reference bx
00309 
00310     //min
00311     if(nBxDiff[i][1]>bxDiff) nBxDiff[i][1] = bxDiff;
00312     if(nBxOccy[i][1]>bx    ) nBxOccy[i][1] = bx;
00313     //max
00314     if(nBxDiff[i][2]<bxDiff) nBxDiff[i][2] = bxDiff;
00315     if(nBxOccy[i][2]<bx    ) nBxOccy[i][2] = bx;
00316 
00317     if(verbose())
00318       std::cout << " fed:" <<  i 
00319                 << " bx:" << bx 
00320                 << " bxRef:" << bxRef
00321                 << " diff:" << bxDiff 
00322                 << " nBxDiff"<<" del:"<<nBxDiff[i][0]<<" min:"<<nBxDiff[i][1]<<" max:"<<nBxDiff[i][2]
00323                 << " nBxOccy"<<" del:"<<nBxOccy[i][0]<<" min:"<<nBxOccy[i][1]<<" max:"<<nBxOccy[i][2]
00324                 << "\n" << std::flush;
00325 
00326     hBxDiffAllFed->Fill(i,bxDiff);
00327     
00328     //if(ttype==1) //skip if not a physics trigger
00329     hBxOccyAllFed->Fill(bx);
00330 
00331 
00332     // done if running in filter farm
00333     if(runInFF_)
00334       continue;
00335 
00336     for(int j=0; j<NSYS; j++)
00337       if(i>=fedRange_[j].first && i<=fedRange_[j].second)
00338           hBxDiffSysFed[j]->Fill(i,bxDiff);
00339           
00340     for(size_t k=0; k!=listGtBits_.size(); k++) {
00341       if((int)gtbits.size() <= listGtBits_.at(k)) {
00342         if(verbose()) 
00343           std::cout << "BxTiming analyze | problem with vector size!\n" << std::endl;
00344         continue;
00345       }
00346       else if(!gtbits.at(listGtBits_.at(k))) 
00347         continue;
00348       for(int j=0; j<NSYS; j++) {
00349         if(i>=fedRange_[j].first && i<=fedRange_[j].second) {
00350           hBxOccyTrigBit[j][k]->Fill(bx);
00351         }
00352       }
00353     }
00354 
00355     if(i>=fedRange_[GLT].first && i<=fedRange_[GLT].second) //GT fed
00356       if(ttype<nttype_)
00357         hBxOccyGtTrigType[ttype-1]->Fill(bx);
00358 
00359     if(ttype!=1) continue; //skip if not a physics trigger
00360     //hBxOccyAllFed->Fill(bx);
00361     hBxOccyOneFed[i]->Fill(bx);
00362 
00363   }
00364 
00365   for(int i=0; i<nfed_;i++) {
00366     nBxDiff[i][0]=nBxDiff[i][2]-nBxDiff[i][1]; 
00367     nBxOccy[i][0]=nBxOccy[i][2]-nBxOccy[i][1];
00368     if(nBxDiff[i][0]<0 || nBxOccy[i][0]<0) continue;
00369     for(int j=0; j<nspr_; j++) {
00370       hBxDiffAllFedSpread[j]->setBinContent(i,nBxDiff[i][j]);      
00371       hBxOccyAllFedSpread[j]->setBinContent(i,nBxOccy[i][j]);      
00372     }
00373     if(verbose())
00374       std::cout << "BxTiming fed:" << i 
00375                 << " Bx-Bx(" << fedRef_ << ")::" 
00376                 << " del:" << nBxDiff[i][0]
00377                 << " min:" << nBxDiff[i][1]
00378                 << " max:" << nBxDiff[i][2]
00379                 << " Occy: "
00380                 << " del:" << nBxOccy[i][0]
00381                 << " min:" << nBxOccy[i][1]
00382                 << " max:" << nBxOccy[i][2]
00383                 <<"\n" << std::flush;
00384 
00385   }
00386 
00387 
00388   if(verbose())
00389     std::cout << "BxTiming::analyze() end.\n" << std::flush;
00390 }

void BxTiming::beginJob ( const edm::EventSetup  )  [protected, virtual]

initialize counters

book the histograms

labeling (cosmetics added here)

Reimplemented from edm::EDAnalyzer.

Definition at line 68 of file BxTiming.cc.

References DQMStore::book1D(), DQMStore::bookProfile(), GenMuonPlsPt100GeV_cfg::cout, dbe, fedRange_, fedRef_, first, flush(), FEDNumbering::getCSCFEDIds(), FEDNumbering::getCSCTFFEDIds(), FEDNumbering::getDTFEDIds(), FEDNumbering::getDTTFFEDIds(), FEDNumbering::getEcalFEDIds(), FEDNumbering::getHcalFEDIds(), FEDNumbering::getRPCFEDIds(), FEDNumbering::getTriggerGCTFEDIds(), FEDNumbering::getTriggerGTPFEDIds(), hBxDiffAllFed, hBxDiffAllFedSpread, hBxDiffSysFed, hBxOccyAllFed, hBxOccyAllFedSpread, hBxOccyGtTrigType, hBxOccyOneFed, hBxOccyTrigBit, histFolder_, i, j, listGtBits_, nbig_, nBxDiff, nBxOccy, nfed_, norb_, nspr_, NSYS, nttype_, runInFF_, edm::second(), MonitorElement::setAxisTitle(), DQMStore::setCurrentFolder(), and verbose().

00068                                        {
00069 
00070   if(verbose())
00071     std::cout << "BxTiming::beginJob()  start\n" << std::flush;
00072 
00073   DQMStore* dbe = 0;
00074   dbe = edm::Service<DQMStore>().operator->();
00075   if(dbe) {
00076     dbe->setCurrentFolder(histFolder_);
00077     //dbe->rmdir(histFolder_);
00078   }
00079 
00081   for(int i=0; i<nfed_;i++) {
00082     nBxDiff[i][0]=0; nBxDiff[i][1]=nbig_; nBxDiff[i][2]=-1*nbig_;
00083     nBxOccy[i][0]=0; nBxOccy[i][1]=nbig_; nBxOccy[i][2]=-1*nbig_;
00084   }
00085 
00086   std::string lbl("");
00087   std::string SysLabel[NSYS] = {
00088     "ECAL", "HCAL", "GCT", "CSCTPG", "CSCTF", "DTTPG", "DTTF", "RPC", "GT"
00089   };
00090   
00091   std::pair<int,int> fedRange[NSYS] = {
00092     FEDNumbering::getEcalFEDIds(),      //600..670
00093     FEDNumbering::getHcalFEDIds(),      //700..731
00094     FEDNumbering::getTriggerGCTFEDIds(),//745..749
00095     FEDNumbering::getCSCFEDIds(),       //750..757
00096     FEDNumbering::getCSCTFFEDIds(),     //760..760
00097     FEDNumbering::getDTFEDIds(),        //770..775
00098     FEDNumbering::getDTTFFEDIds(),      //780..780
00099     FEDNumbering::getRPCFEDIds(),       //790..795
00100     FEDNumbering::getTriggerGTPFEDIds() //812..813
00101   };
00102   for(int i=0; i<NSYS; i++) fedRange_[i]=fedRange[i];
00103 
00104 
00105   int fedRefSys=-1;
00106   for(int i=0; i<NSYS; i++)
00107     if(fedRef_>=fedRange_[i].first && fedRef_<=fedRange_[i].second)
00108       {fedRefSys=i; break;}
00109   std::string refName("");
00110   std::string spreadLabel[nspr_] = {"Spread","Min", "Max"};
00111   if(fedRefSys>=0)
00112     refName+=SysLabel[fedRefSys];
00113   else
00114     refName+=fedRef_;
00115 
00117 
00118   const int dbx = nbig_;
00119 
00120   if(dbe) {
00121 
00122     dbe->setCurrentFolder(histFolder_);
00123 
00124     hBxDiffAllFed = dbe->bookProfile("BxDiffAllFed", "BxDiffAllFed", 
00125                                      nfed_ + 1, -0.5, nfed_+0.5, 
00126                                      2*dbx+1, -1*dbx-0.5,dbx+0.5
00127                                      );
00128 
00129     for(int i=0; i<nspr_; i++) {
00130       lbl.clear();lbl+="BxDiffAllFed";lbl+=spreadLabel[i];
00131       hBxDiffAllFedSpread[i] = dbe->book1D(lbl.data(),lbl.data(), nfed_ + 1, -0.5, nfed_+0.5); 
00132       lbl.clear();lbl+="BxOccyAllFed";lbl+=spreadLabel[i];
00133       hBxOccyAllFedSpread[i] = dbe->book1D(lbl.data(),lbl.data(), nfed_ + 1, -0.5, nfed_+0.5); 
00134     }
00135 
00136     lbl.clear();lbl+="BxOccyAllFed";
00137     hBxOccyAllFed = dbe->book1D(lbl.data(),lbl.data(),norb_+1,-0.5,norb_+0.5);
00138 
00139   }
00140 
00141   // following histos defined only when not runing in the ff
00142   if(dbe && !runInFF_) {
00143     
00144     dbe->setCurrentFolder(histFolder_);
00145     
00146     for(int i=0; i<NSYS; i++) {
00147       lbl.clear();lbl+=SysLabel[i];lbl+="FedBxDiff"; 
00148       int nfeds = fedRange_[i].second - fedRange_[i].first + 1;
00149       nfeds = (nfeds>0)? nfeds:1;
00150       hBxDiffSysFed[i] = dbe->bookProfile(lbl.data(),lbl.data(), nfeds, 
00151                                           fedRange_[i].first-0.5, fedRange_[i].second+0.5,
00152                                           2*dbx+1,-1*dbx-0.5,dbx+0.5);
00153     }
00154 
00155 
00156     lbl.clear();lbl+="BxOccyAllFed";
00157     hBxOccyOneFed = new MonitorElement*[nfed_];
00158     dbe->setCurrentFolder(histFolder_+"SingleFed/");
00159     for(int i=0; i<nfed_; i++) {
00160       lbl.clear(); lbl+="BxOccyOneFed";
00161       char *ii = new char[1000]; std::sprintf(ii,"%d",i);lbl+=ii;
00162       hBxOccyOneFed[i] = dbe->book1D(lbl.data(),lbl.data(),norb_+1,-0.5,norb_+0.5);
00163       delete ii;
00164     }
00165 
00166     dbe->setCurrentFolder(histFolder_);
00167     for(int i=0; i<nttype_; i++) {
00168       lbl.clear();lbl+="BxOccyGtTrigType";
00169       char *ii = new char[10]; std::sprintf(ii,"%d",i+1);lbl+=ii;
00170       hBxOccyGtTrigType[i] = dbe->book1D(lbl.data(),lbl.data(),norb_+1,-0.5,norb_+0.5);
00171     }
00172 
00173     dbe->setCurrentFolder(histFolder_+"SingleBit/");
00174     for(int i=0; i<NSYS; i++) {
00175       hBxOccyTrigBit[i] = new MonitorElement*[listGtBits_.size()];
00176       for(size_t j=0; j<listGtBits_.size(); j++) {
00177         lbl.clear();lbl+=SysLabel[i];lbl+="BxOccyGtBit"; 
00178         char *ii = new char[1000]; std::sprintf(ii,"%d",listGtBits_.at(j)); lbl+=ii;
00179         hBxOccyTrigBit[i][j] = dbe->book1D(lbl.data(),lbl.data(),norb_+1,-0.5,norb_+0.5);
00180       }
00181     }
00182 
00183   }
00184   
00186   hBxDiffAllFed->setAxisTitle("FED ID",1);
00187   lbl.clear(); lbl+="BX(fed)-BX("; lbl+=refName; lbl+=")";
00188   hBxDiffAllFed->setAxisTitle(lbl,2);
00189   for(int i=0; i<nspr_; i++) {
00190     lbl.clear(); lbl+="BX(fed)-BX("; lbl+=refName; lbl+=") "+spreadLabel[i];
00191     hBxDiffAllFedSpread[i]->setAxisTitle("FED ID",1);
00192     hBxDiffAllFedSpread[i]->setAxisTitle(lbl,2);
00193     lbl.clear(); lbl+="Bx FED occupancy"; lbl+=" "; lbl+=spreadLabel[i]; 
00194     hBxOccyAllFedSpread[i]->setAxisTitle("FED ID",1); 
00195     hBxOccyAllFedSpread[i]->setAxisTitle(lbl,2);
00196   }
00197 
00198   hBxOccyAllFed->setAxisTitle("bx",1);
00199   lbl.clear(); lbl+="Combined FED occupancy";
00200   hBxOccyAllFed->setAxisTitle(lbl,2);
00201  
00202   // skip next if running in filter farm
00203   if(runInFF_)
00204     return;
00205 
00206   for(int i=0; i<NSYS; i++) {
00207     lbl.clear(); lbl+=SysLabel[i]; lbl+=" FED ID";
00208     hBxDiffSysFed[i]->setAxisTitle(lbl,1);
00209     lbl.clear(); lbl+="BX("; lbl+=SysLabel[i]; lbl+=")-BX(";lbl+=refName; lbl+=")";
00210     hBxDiffSysFed[i]->setAxisTitle(lbl,2);
00211   }
00212 
00213   for(int i=0; i<nfed_; i++) {
00214     hBxOccyOneFed[i] ->setAxisTitle("bx",1);
00215     lbl.clear(); lbl+=" FED "; char *ii = new char[1000]; std::sprintf(ii,"%d",i);lbl+=ii; lbl+=" occupancy";
00216     hBxOccyOneFed[i] ->setAxisTitle(lbl,2);
00217     delete ii;
00218   }
00219   for(int i=0; i<nttype_; i++) {
00220     hBxOccyGtTrigType[i]->setAxisTitle("bx",1);
00221     lbl.clear(); lbl+="GT occupancy for trigger type "; char *ii = new char[10]; std::sprintf(ii,"%d",i+1);lbl+=ii;
00222     hBxOccyGtTrigType[i]->setAxisTitle(lbl,2);
00223   }
00224   
00225   for(int i=0; i<NSYS; i++) {
00226     for(size_t j=0; j<listGtBits_.size(); j++) {
00227       hBxOccyTrigBit[i][j]->setAxisTitle("bx",1);
00228       lbl.clear();lbl+=SysLabel[i];lbl+=" Bx occupancy for Trigger bit "; 
00229       char *ii = new char[10]; std::sprintf(ii,"%d",listGtBits_.at(j)); lbl+=ii;
00230       hBxOccyTrigBit[i][j]->setAxisTitle(lbl,2);
00231     }
00232   }
00233     
00234   if(verbose())
00235     std::cout << "BxTiming::beginJob()  end.\n" << std::flush;
00236 }

int BxTiming::calcBxDiff ( int  bx1,
int  bx2 
) [private]

calculates the difference (closest distance) between two bunch crossing numbers.

This is similar to calculating delta phi between two vectors.

Calculates bx1 - bx2 but makes sure that the value is in the range -num_bx_per_orbit / 2 .. + num_bx_per_orbit / 2 .

Definition at line 395 of file BxTiming.cc.

References diff, half_norb_, and norb_.

Referenced by analyze().

00396 {
00397   int diff = bx1 - bx2;
00398   
00399   while (diff < -half_norb_)
00400     diff += norb_;
00401   
00402   while (diff > half_norb_)
00403     diff -= norb_;
00404 
00405   return diff;
00406 }

void BxTiming::endJob ( void   )  [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 239 of file BxTiming.cc.

References GenMuonPlsPt100GeV_cfg::cout, dbe, flush(), histFile_, nEvt_, DQMStore::save(), and verbose().

00239                  {
00240 
00241   if(verbose())
00242     std::cout << "BxTiming::endJob Nevents: " << nEvt_ << "\n" << std::flush;
00243 
00244   if(histFile_.size()!=0  && dbe) 
00245     dbe->save(histFile_);
00246   
00247   if(verbose())
00248     std::cout << "BxTiming::endJob()  end.\n" << std::flush;
00249 }

int BxTiming::verbose (  )  [inline, private]

Definition at line 46 of file BxTiming.h.

References verbose_.

Referenced by analyze(), beginJob(), BxTiming(), and endJob().

00046 {return verbose_;}


Member Data Documentation

DQMStore* BxTiming::dbe [private]

Definition at line 66 of file BxTiming.h.

Referenced by beginJob(), BxTiming(), and endJob().

std::pair<int,int> BxTiming::fedRange_[NSYS] [private]

Definition at line 82 of file BxTiming.h.

Referenced by analyze(), and beginJob().

int BxTiming::fedRef_ [private]

Definition at line 84 of file BxTiming.h.

Referenced by analyze(), beginJob(), and BxTiming().

edm::InputTag BxTiming::fedSource_ [private]

Definition at line 41 of file BxTiming.h.

Referenced by analyze(), and BxTiming().

edm::InputTag BxTiming::gtSource_ [private]

Definition at line 42 of file BxTiming.h.

Referenced by analyze(), and BxTiming().

const int BxTiming::half_norb_ = norb_ / 2 [static, private]

Definition at line 73 of file BxTiming.h.

Referenced by calcBxDiff().

MonitorElement* BxTiming::hBxDiffAllFed [private]

histograms

Definition at line 92 of file BxTiming.h.

Referenced by analyze(), and beginJob().

MonitorElement* BxTiming::hBxDiffAllFedSpread[nspr_] [private]

Definition at line 97 of file BxTiming.h.

Referenced by analyze(), and beginJob().

MonitorElement* BxTiming::hBxDiffSysFed[NSYS] [private]

Definition at line 93 of file BxTiming.h.

Referenced by analyze(), and beginJob().

MonitorElement* BxTiming::hBxOccyAllFed [private]

Definition at line 94 of file BxTiming.h.

Referenced by analyze(), and beginJob().

MonitorElement* BxTiming::hBxOccyAllFedSpread[nspr_] [private]

Definition at line 98 of file BxTiming.h.

Referenced by analyze(), and beginJob().

MonitorElement* BxTiming::hBxOccyGtTrigType[nttype_] [private]

Definition at line 100 of file BxTiming.h.

Referenced by analyze(), and beginJob().

MonitorElement** BxTiming::hBxOccyOneFed [private]

Definition at line 95 of file BxTiming.h.

Referenced by analyze(), and beginJob().

MonitorElement** BxTiming::hBxOccyTrigBit[NSYS] [private]

Definition at line 101 of file BxTiming.h.

Referenced by analyze(), and beginJob().

std::string BxTiming::histFile_ [private]

Definition at line 60 of file BxTiming.h.

Referenced by BxTiming(), and endJob().

std::string BxTiming::histFolder_ [private]

Definition at line 63 of file BxTiming.h.

Referenced by beginJob(), and BxTiming().

std::vector<int> BxTiming::listGtBits_ [private]

Definition at line 78 of file BxTiming.h.

Referenced by analyze(), beginJob(), and BxTiming().

const int BxTiming::nbig_ = 10000 [static, private]

Definition at line 75 of file BxTiming.h.

Referenced by beginJob().

int BxTiming::nBxDiff[1500][nspr_] [private]

Definition at line 88 of file BxTiming.h.

Referenced by analyze(), and beginJob().

int BxTiming::nBxOccy[1500][nspr_] [private]

Definition at line 89 of file BxTiming.h.

Referenced by analyze(), and beginJob().

int BxTiming::nEvt_ [private]

Definition at line 57 of file BxTiming.h.

Referenced by analyze(), BxTiming(), and endJob().

int BxTiming::nfed_ [private]

Definition at line 83 of file BxTiming.h.

Referenced by analyze(), beginJob(), and BxTiming().

const int BxTiming::norb_ = 3564 [static, private]

Definition at line 72 of file BxTiming.h.

Referenced by beginJob(), and calcBxDiff().

const int BxTiming::nspr_ = 3 [static, private]

Definition at line 87 of file BxTiming.h.

Referenced by analyze(), and beginJob().

const int BxTiming::nttype_ = 6 [static, private]

Definition at line 76 of file BxTiming.h.

Referenced by analyze(), and beginJob().

bool BxTiming::runInFF_ [private]

Definition at line 69 of file BxTiming.h.

Referenced by analyze(), beginJob(), and BxTiming().

int BxTiming::verbose_ [private]

Definition at line 45 of file BxTiming.h.

Referenced by BxTiming(), and verbose().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:15:36 2009 for CMSSW by  doxygen 1.5.4