CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes

L1TGMT Class Reference

#include <L1TGMT.h>

Inheritance diagram for L1TGMT:
edm::EDAnalyzer

List of all members.

Public Member Functions

 L1TGMT (const edm::ParameterSet &ps)
virtual ~L1TGMT ()

Protected Member Functions

void analyze (const edm::Event &e, const edm::EventSetup &c)
void beginJob (void)
void beginRun (const edm::Run &r, const edm::EventSetup &c)
void endJob (void)

Private Types

enum  ensubs {
  DTTF = 0, RPCb, CSCTF, RPCf,
  GMT
}

Private Member Functions

void book_ (const edm::EventSetup &c)
double phiconv_ (float phi)

Private Attributes

MonitorElementbx_csc_rpc
MonitorElementbx_dt_csc
MonitorElementbx_dt_rpc
MonitorElementbx_number
int bxnum_old_
DQMStoredbe
MonitorElementdbx_chip
MonitorElementdist_eta_csc_rpc
MonitorElementdist_eta_dt_csc
MonitorElementdist_eta_dt_rpc
MonitorElementdist_phi_csc_rpc
MonitorElementdist_phi_dt_csc
MonitorElementdist_phi_dt_rpc
MonitorElementeta_dtcsc_and_rpc
MonitorElementeta_dtcsc_only
MonitorElementeta_rpc_only
MonitorElementetaphi_dtcsc_and_rpc
MonitorElementetaphi_dtcsc_only
MonitorElementetaphi_rpc_only
int evnum_old_
edm::InputTag gmtSource_
ofstream logFile_
bool monitorDaemon_
MonitorElementn_csctf_vs_dttf
MonitorElementn_rpcb_vs_dttf
MonitorElementn_rpcf_vs_csctf
int nev_
int obnum_old_
std::string outputFile_
MonitorElementphi_dtcsc_and_rpc
MonitorElementphi_dtcsc_only
MonitorElementphi_rpc_only
MonitorElementregional_triggers
MonitorElementsubs_bits [5]
MonitorElementsubs_dbx [4]
MonitorElementsubs_eta [5]
MonitorElementsubs_etaphi [5]
MonitorElementsubs_etaqty [5]
MonitorElementsubs_nbx [5]
MonitorElementsubs_phi [5]
MonitorElementsubs_pt [5]
MonitorElementsubs_qty [5]
int trsrc_old_
bool verbose_

Static Private Attributes

static const double piconv_ = 180. / acos(-1.)

Detailed Description

Definition at line 44 of file L1TGMT.h.


Member Enumeration Documentation

enum L1TGMT::ensubs [private]
Enumerator:
DTTF 
RPCb 
CSCTF 
RPCf 
GMT 

Definition at line 71 of file L1TGMT.h.

{ DTTF=0, RPCb, CSCTF, RPCf, GMT };

Constructor & Destructor Documentation

L1TGMT::L1TGMT ( const edm::ParameterSet ps)

Definition at line 25 of file L1TGMT.cc.

References gather_cfg::cout, dbe, edm::ParameterSet::getUntrackedParameter(), NULL, cmsCodeRules::cppFunctionSkipper::operator, outputFile_, DQMStore::setCurrentFolder(), DQMStore::setVerbose(), and verbose_.

  : gmtSource_( ps.getParameter< InputTag >("gmtSource") )
 {

  // verbosity switch
  verbose_ = ps.getUntrackedParameter<bool>("verbose", false);

  if(verbose_) cout << "L1TGMT: constructor...." << endl;


  dbe = NULL;
  if ( ps.getUntrackedParameter<bool>("DQMStore", false) ) 
  {
    dbe = Service<DQMStore>().operator->();
    dbe->setVerbose(0);
  }

  outputFile_ = ps.getUntrackedParameter<string>("outputFile", "");
  if ( outputFile_.size() != 0 ) {
    cout << "L1T Monitoring histograms will be saved to " << outputFile_.c_str() << endl;
  }

  bool disable = ps.getUntrackedParameter<bool>("disableROOToutput", false);
  if(disable){
    outputFile_="";
  }


  if ( dbe !=NULL ) {
    dbe->setCurrentFolder("L1T/L1TGMT");
  }

}
L1TGMT::~L1TGMT ( ) [virtual]

Definition at line 59 of file L1TGMT.cc.

{
}

Member Function Documentation

void L1TGMT::analyze ( const edm::Event e,
const edm::EventSetup c 
) [protected, virtual]

Implements edm::EDAnalyzer.

Definition at line 92 of file L1TGMT.cc.

References begin, bx_csc_rpc, bx_dt_csc, bx_dt_rpc, bx_number, bxnum_old_, gather_cfg::cout, CSCTF, dist_eta_csc_rpc, dist_eta_dt_csc, dist_eta_dt_rpc, dist_phi_csc_rpc, dist_phi_dt_csc, dist_phi_dt_rpc, DTTF, eta_dtcsc_and_rpc, eta_dtcsc_only, eta_rpc_only, etaphi_dtcsc_and_rpc, etaphi_dtcsc_only, etaphi_rpc_only, evnum_old_, MonitorElement::Fill(), edm::Event::getByLabel(), L1MuGMTReadoutCollection::getRecords(), GMT, gmtSource_, i, edm::HandleBase::isValid(), j, edm::InputTag::label(), n_csctf_vs_dttf, n_rpcb_vs_dttf, n_rpcf_vs_csctf, nev_, obnum_old_, edm::EventBase::orbitNumber(), phi_dtcsc_and_rpc, phi_dtcsc_only, phi_rpc_only, phiconv_(), edm::Handle< T >::product(), regional_triggers, RPCb, RPCf, subs_bits, subs_dbx, subs_eta, subs_etaphi, subs_etaqty, subs_nbx, subs_phi, subs_pt, subs_qty, trsrc_old_, and verbose_.

{
  
  nev_++; 
  if(verbose_) cout << "L1TGMT: analyze...." << endl;


  edm::Handle<L1MuGMTReadoutCollection> pCollection;
  e.getByLabel(gmtSource_,pCollection);
  
  if (!pCollection.isValid()) {
    edm::LogInfo("DataNotFound") << "can't find L1MuGMTReadoutCollection with label "
    << gmtSource_.label() ;
    return;
  }

  // remember the bx of 1st candidate of each system (9=none)
  int bx1st[4] = {9, 9, 9, 9};

  // get GMT readout collection
  L1MuGMTReadoutCollection const* gmtrc = pCollection.product();
  // get record vector
  vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords();
  // loop over records of individual bx's
  vector<L1MuGMTReadoutRecord>::const_iterator RRItr;
  
  for( RRItr = gmt_records.begin(); RRItr != gmt_records.end(); RRItr++ ) 
  {
    
    vector<L1MuRegionalCand> INPCands[4] = {
        RRItr->getDTBXCands(),
        RRItr->getBrlRPCCands(),
        RRItr->getCSCCands(),
        RRItr->getFwdRPCCands()
    };
    vector<L1MuGMTExtendedCand> GMTCands   = RRItr->getGMTCands();
    
    vector<L1MuRegionalCand>::const_iterator INPItr;
    vector<L1MuGMTExtendedCand>::const_iterator GMTItr;
    vector<L1MuGMTExtendedCand>::const_iterator GMTItr2;
    
    int BxInEvent = RRItr->getBxInEvent();
    
    // count non-empty candidates in this bx
    int nSUBS[5] = {0, 0, 0, 0, 0};
    for(int i=0; i<4; i++) {
      for( INPItr = INPCands[i].begin(); INPItr != INPCands[i].end(); ++INPItr ) {
        if(!INPItr->empty()) {
          nSUBS[i]++;
          if(bx1st[i]==9) bx1st[i]=BxInEvent;
        }
      }      
      subs_nbx[i]->Fill(float(nSUBS[i]),float(BxInEvent));
    }

    for( GMTItr = GMTCands.begin(); GMTItr != GMTCands.end(); ++GMTItr ) {
      if(!GMTItr->empty()) nSUBS[GMT]++;
    }
    subs_nbx[GMT]->Fill(float(nSUBS[GMT]),float(BxInEvent));
    
    // from here care only about the L1A bunch crossing
    if(BxInEvent!=0) continue;
    
    // get the absolute bx number of the L1A
    int Bx = RRItr->getBxNr();
    int Ev = RRItr->getEvNr();
    
    bx_number->Fill(double(Bx));
 
    for(int i=0; i<4; i++) {
      for( INPItr = INPCands[i].begin(); INPItr != INPCands[i].end(); ++INPItr ) {
        if(INPItr->empty()) continue;
        subs_eta[i]->Fill(INPItr->etaValue());
        subs_phi[i]->Fill(phiconv_(INPItr->phiValue()));
        subs_pt[i]->Fill(INPItr->ptValue());
        subs_qty[i]->Fill(INPItr->quality());
        subs_etaphi[i]->Fill(INPItr->etaValue(),phiconv_(INPItr->phiValue()));
        subs_etaqty[i]->Fill(INPItr->etaValue(),INPItr->quality());
        int word = INPItr->getDataWord();
        for( int j=0; j<32; j++ ) {
          if( word&(1<<j) ) subs_bits[i]->Fill(float(j));
        }
      }
    }
        
    for( GMTItr = GMTCands.begin(); GMTItr != GMTCands.end(); ++GMTItr ) {
      if(GMTItr->empty()) continue;
      subs_eta[GMT]->Fill(GMTItr->etaValue());
      subs_phi[GMT]->Fill(phiconv_(GMTItr->phiValue()));
      subs_pt[GMT]->Fill(GMTItr->ptValue());
      subs_qty[GMT]->Fill(GMTItr->quality());
      subs_etaphi[GMT]->Fill(GMTItr->etaValue(),phiconv_(GMTItr->phiValue()));
      subs_etaqty[GMT]->Fill(GMTItr->etaValue(),GMTItr->quality());
      int word = GMTItr->getDataWord();
      for( int j=0; j<32; j++ ) {
        if( word&(1<<j) ) subs_bits[GMT]->Fill(float(j));
      }
      
      if(GMTItr->isMatchedCand()) {
        if(GMTItr->quality()>3) {
          eta_dtcsc_and_rpc->Fill(GMTItr->etaValue());
          phi_dtcsc_and_rpc->Fill(phiconv_(GMTItr->phiValue()));
          etaphi_dtcsc_and_rpc->Fill(GMTItr->etaValue(),phiconv_(GMTItr->phiValue()));
        }
      } else if(GMTItr->isRPC()) {
        if(GMTItr->quality()>3) {
          eta_rpc_only->Fill(GMTItr->etaValue());
          phi_rpc_only->Fill(phiconv_(GMTItr->phiValue()));
          etaphi_rpc_only->Fill(GMTItr->etaValue(),phiconv_(GMTItr->phiValue()));        
        }
      } else {
        if(GMTItr->quality()>3) {
          eta_dtcsc_only->Fill(GMTItr->etaValue());
          phi_dtcsc_only->Fill(phiconv_(GMTItr->phiValue()));
          etaphi_dtcsc_only->Fill(GMTItr->etaValue(),phiconv_(GMTItr->phiValue()));
        }
        
        if(GMTItr != GMTCands.end()){
          for( GMTItr2 = GMTCands.begin(); GMTItr2 != GMTCands.end(); ++GMTItr2 ) {
            if(GMTItr2==GMTItr) continue;
            if(GMTItr2->empty()) continue;
            if(GMTItr2->isRPC()) {
              if(GMTItr->isFwd()) {
                dist_eta_csc_rpc->Fill( GMTItr->etaValue() - GMTItr2->etaValue() );
                dist_phi_csc_rpc->Fill( phiconv_(GMTItr->phiValue()) - phiconv_(GMTItr2->phiValue()) );
              } else {
                dist_eta_dt_rpc->Fill( GMTItr->etaValue() - GMTItr2->etaValue() );
                dist_phi_dt_rpc->Fill( phiconv_(GMTItr->phiValue()) - phiconv_(GMTItr2->phiValue()) );                
              }
            } else {
              if(!(GMTItr->isFwd()) && GMTItr2->isFwd()) {
                dist_eta_dt_csc->Fill( GMTItr->etaValue() - GMTItr2->etaValue() );
                dist_phi_dt_csc->Fill( phiconv_(GMTItr->phiValue()) - phiconv_(GMTItr2->phiValue()) );
              } else if(GMTItr->isFwd() && !(GMTItr2->isFwd())){
                dist_eta_dt_csc->Fill( GMTItr2->etaValue() - GMTItr->etaValue() );
                dist_phi_dt_csc->Fill( phiconv_(GMTItr->phiValue()) - phiconv_(GMTItr2->phiValue()) );                
              }
            }
          }     
        }
        
      }
      
    }
    
    n_rpcb_vs_dttf ->Fill(float(nSUBS[DTTF]),float(nSUBS[RPCb]));
    n_rpcf_vs_csctf->Fill(float(nSUBS[CSCTF]),float(nSUBS[RPCf]));
    n_csctf_vs_dttf->Fill(float(nSUBS[DTTF]),float(nSUBS[CSCTF]));
    
    regional_triggers->Fill(-1.); // fill underflow for normalization
    if(nSUBS[GMT]) regional_triggers->Fill(0.); // fill all muon bin
    int ioff=1;
    for(int i=0; i<4; i++) {
      if(nSUBS[i]) regional_triggers->Fill(float(5*i+nSUBS[i]+ioff));
    }
    if(nSUBS[DTTF] && (nSUBS[RPCb] || nSUBS[RPCf])) regional_triggers->Fill(22.);
    if(nSUBS[DTTF] && nSUBS[CSCTF]) regional_triggers->Fill(23.);
    if(nSUBS[CSCTF] && (nSUBS[RPCb] || nSUBS[RPCf])) regional_triggers->Fill(24.);
    if(nSUBS[DTTF] && nSUBS[CSCTF] && (nSUBS[RPCb] || nSUBS[RPCf])) regional_triggers->Fill(25.);
        
    // fill only if previous event corresponds to previous trigger
//    if( (Ev - evnum_old_) == 1 && bxnum_old_ > -1 ) {
    // assume getting all events in a sequence (usefull only from reco data)
      if( bxnum_old_ > -1 ) {
      int dBx = Bx - bxnum_old_ + 3564*(e.orbitNumber() - obnum_old_);
      for(int id = 0; id<4; id++) {
        if( trsrc_old_&(1<<id) ) {
          for(int i=0; i<4; i++) {
            if(nSUBS[i]) subs_dbx[i]->Fill(float(dBx),float(id));
          }
        }
      }
      
    }
    
    // save quantities for the next event
    evnum_old_ = Ev;
    bxnum_old_ = Bx;
    obnum_old_ = e.orbitNumber();
    trsrc_old_ = 0;
    for(int i=0; i<4; i++) {
      if(nSUBS[i])  trsrc_old_ |= (1<<i);
    }
  }
  
  if(bx1st[DTTF]<9  && bx1st[RPCb]<9)  bx_dt_rpc->Fill(bx1st[DTTF], bx1st[RPCb]);
  if(bx1st[CSCTF]<9 && bx1st[RPCf]<9) bx_csc_rpc->Fill(bx1st[CSCTF],bx1st[RPCf]);
  if(bx1st[DTTF]<9 && bx1st[CSCTF]<9)  bx_dt_csc->Fill(bx1st[DTTF], bx1st[CSCTF]);

}
void L1TGMT::beginJob ( void  ) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 63 of file L1TGMT.cc.

References bxnum_old_, evnum_old_, and nev_.

{

  nev_ = 0;
  evnum_old_ = -1;
  bxnum_old_ = -1;

}
void L1TGMT::beginRun ( const edm::Run r,
const edm::EventSetup c 
) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 72 of file L1TGMT.cc.

References book_(), and nev_.

{
  
  if(nev_==0) {
      book_(c);
    }

}
void L1TGMT::book_ ( const edm::EventSetup c) [private]

Definition at line 291 of file L1TGMT.cc.

References DQMStore::book1D(), DQMStore::book2D(), DQMStore::bookProfile(), bx_csc_rpc, bx_dt_csc, bx_dt_rpc, bx_number, CSCTF, AlCaHLTBitMon_QueryRunRegistry::data, dbe, dbx_chip, dist_eta_csc_rpc, dist_eta_dt_csc, dist_eta_dt_rpc, dist_phi_csc_rpc, dist_phi_dt_csc, dist_phi_dt_rpc, DTTF, eta_dtcsc_and_rpc, eta_dtcsc_only, eta_rpc_only, etaphi_dtcsc_and_rpc, etaphi_dtcsc_only, etaphi_rpc_only, edm::EventSetup::get(), L1MuTriggerScales::getGMTEtaScale(), L1MuScale::getNBins(), L1MuTriggerScales::getPhiScale(), L1MuTriggerPtScale::getPtScale(), L1MuTriggerScales::getRegionalEtaScale(), L1MuScale::getValue(), GMT, i, j, n_csctf_vs_dttf, n_rpcb_vs_dttf, n_rpcf_vs_csctf, pileupCalc::nbins, cmsCodeRules::cppFunctionSkipper::operator, phi_dtcsc_and_rpc, phi_dtcsc_only, phi_rpc_only, piconv_, edm::ESHandle< T >::product(), regional_triggers, DQMStore::rmdir(), RPCb, RPCf, MonitorElement::setAxisTitle(), MonitorElement::setBinLabel(), DQMStore::setCurrentFolder(), subs_bits, subs_dbx, subs_eta, subs_etaphi, subs_etaqty, subs_nbx, subs_phi, subs_pt, and subs_qty.

Referenced by beginRun().

{

  std::string subs[5] = { "DTTF", "RPCb", "CSCTF", "RPCf", "GMT" };

  edm::ESHandle< L1MuTriggerScales > trigscales_h;
  c.get< L1MuTriggerScalesRcd >().get( trigscales_h );
  const L1MuTriggerScales* scales = trigscales_h.product();

  edm::ESHandle< L1MuTriggerPtScale > trigptscale_h;
  c.get< L1MuTriggerPtScaleRcd >().get( trigptscale_h );
  const L1MuTriggerPtScale* scalept = trigptscale_h.product();  

  // get hold of back-end interface
  DQMStore* dbe = 0;
  dbe = Service<DQMStore>().operator->();

  if ( dbe ) {
    dbe->setCurrentFolder("L1T/L1TGMT");
    dbe->rmdir("L1T/L1TGMT");
  }


  if ( dbe ) 
  {
    dbe->setCurrentFolder("L1T/L1TGMT");
    
    int nqty=8; double qtymin=-0.5; double qtymax=7.5;

    float phiscale[145];
    int nphiscale;
    {
      int nbins = scales->getPhiScale()->getNBins();
      if(nbins>144) nbins=144;
      for(int j=0; j<=nbins; j++) {
        phiscale[j] = piconv_ * scales->getPhiScale()->getValue(j);
      }
      nphiscale = nbins;
    }

    float qscale[9];
    {
      for(int j=0; j<9; j++) {
        qscale[j] = -0.5 + j;
      }
    } 
    
    // pt scale first bin reserved for empty muon
    float ptscale[32];
    int nptscale;
    {
      int nbins = scalept->getPtScale()->getNBins() - 1;
      if(nbins>31) nbins=31;
      for(int j=1; j<=nbins; j++) {
        ptscale[j-1] = scalept->getPtScale()->getValue(j);
      }
      ptscale[nbins]=ptscale[nbins-1]+10.; // make reasonable size last bin
      nptscale = nbins;
    }
      
    float etascale[5][66];
    int netascale[5];
    // DTTF eta scale
    {
      int nbins = scales->getRegionalEtaScale(DTTF)->getNBins();
      if(nbins>65) nbins = 65;
      for(int j=0; j<=nbins; j++) {
        etascale[DTTF][j] = scales->getRegionalEtaScale(DTTF)->getValue(j);
      }
      netascale[DTTF]=nbins;
    }
    // RPCb etascale
    {
      int nbins = scales->getRegionalEtaScale(RPCb)->getNBins();
      if(nbins>65) nbins = 65;
      for(int j=0; j<=nbins; j++) {
        etascale[RPCb][j] = scales->getRegionalEtaScale(RPCb)->getValue(j);
      }
      netascale[RPCb]=nbins;
    }
    // CSCTF etascale
    // special case - need to mirror 2*32 bins
    {
      int nbins = scales->getRegionalEtaScale(CSCTF)->getNBins();
      if(nbins>32) nbins = 32;
      
      int i=0;
      for(int j=nbins; j>=0; j--,i++) {
        etascale[CSCTF][i] = (-1) * scales->getRegionalEtaScale(CSCTF)->getValue(j);
      }
      for(int j=0; j<=nbins; j++,i++) {
        etascale[CSCTF][i] = scales->getRegionalEtaScale(CSCTF)->getValue(j);
      }
      netascale[CSCTF]=i-1;
    }
    // RPCf etascale
    {
      int nbins = scales->getRegionalEtaScale(RPCf)->getNBins();
      if(nbins>65) nbins = 65;
      for(int j=0; j<=nbins; j++) {
        etascale[RPCf][j] = scales->getRegionalEtaScale(RPCf)->getValue(j);
      }
      netascale[RPCf]=nbins;
    }
    // GMT etascale
    {
      int nbins = scales->getGMTEtaScale()->getNBins();
      if(nbins>32) nbins = 32;
      
      int i=0;
      for(int j=nbins; j>0; j--,i++) {
        etascale[GMT][i] = (-1) * scales->getGMTEtaScale()->getValue(j);
      }
      for(int j=0; j<=nbins; j++,i++) {
        etascale[GMT][i] = scales->getGMTEtaScale()->getValue(j);
      }
      netascale[GMT]=i-1;
    }
    
    
    std::string hname("");
    std::string htitle("");
    
    for(int i=0; i<5; i++) {
      
      hname = subs[i] + "_nbx"; htitle = subs[i] + " multiplicity in bx";
      subs_nbx[i] = dbe->book2D(hname.data(),htitle.data(), 4, 1., 5., 5, -2.5, 2.5);
      subs_nbx[i]->setAxisTitle(subs[i] + " candidates",1);
      subs_nbx[i]->setAxisTitle("bx wrt L1A",2);
      
      hname = subs[i] + "_eta"; htitle = subs[i] + " eta value";
      subs_eta[i] = dbe->book1D(hname.data(),htitle.data(), netascale[i], etascale[i]);
      subs_eta[i]->setAxisTitle("eta",1);
      
      hname = subs[i] + "_phi"; htitle = subs[i] + " phi value";
      subs_phi[i] = dbe->book1D(hname.data(),htitle.data(), nphiscale, phiscale);
      subs_phi[i]->setAxisTitle("phi (deg)",1);
      
      hname = subs[i] + "_pt"; htitle = subs[i] + " pt value";
      subs_pt[i]  = dbe->book1D(hname.data(),htitle.data(), nptscale, ptscale);
      subs_pt[i]->setAxisTitle("L1 pT (GeV)",1);
      
      hname = subs[i] + "_qty"; htitle = subs[i] + " qty value";
      subs_qty[i] = dbe->book1D(hname.data(),htitle.data(), nqty, qtymin, qtymax);
      subs_qty[i]->setAxisTitle(subs[i] + " quality",1);
      
      hname = subs[i] + "_etaphi"; htitle = subs[i] + " phi vs eta";
      subs_etaphi[i] = dbe->book2D(hname.data(),htitle.data(), netascale[i], etascale[i], nphiscale, phiscale);
      subs_etaphi[i]->setAxisTitle("eta",1);
      subs_etaphi[i]->setAxisTitle("phi (deg)",2);
      
      hname = subs[i] + "_etaqty"; htitle = subs[i] + " qty vs eta";
      subs_etaqty[i] = dbe->book2D(hname.data(),htitle.data(), netascale[i], etascale[i], nqty, qscale);
      subs_etaqty[i]->setAxisTitle("eta",1);
      subs_etaqty[i]->setAxisTitle(subs[i] + " quality",2);
      
      hname = subs[i] + "_bits"; htitle = subs[i] + " bit population";
      subs_bits[i] = dbe->book1D(hname.data(),htitle.data(), 32, -0.5, 31.5);
      subs_bits[i]->setAxisTitle("bit number",1);
    }
    
    regional_triggers = dbe->book1D("Regional_trigger","Muon trigger contribution", 27, 0., 27.);
    regional_triggers->setAxisTitle("regional trigger",1);
    int ib=1;
    regional_triggers->setBinLabel(ib++,"All muons",1);
    ib++;
    regional_triggers->setBinLabel(ib++,"DT 1mu",1);
    regional_triggers->setBinLabel(ib++,"DT 2mu",1);
    regional_triggers->setBinLabel(ib++,"DT 3mu",1);
    regional_triggers->setBinLabel(ib++,"DT 4mu",1);
    ib++;
    regional_triggers->setBinLabel(ib++,"RPCb 1mu",1);
    regional_triggers->setBinLabel(ib++,"RPCb 2mu",1);
    regional_triggers->setBinLabel(ib++,"RPCb 3mu",1);
    regional_triggers->setBinLabel(ib++,"RPCb 4mu",1);
    ib++;
    regional_triggers->setBinLabel(ib++,"CSC 1mu",1);
    regional_triggers->setBinLabel(ib++,"CSC 2mu",1);
    regional_triggers->setBinLabel(ib++,"CSC 3mu",1);
    regional_triggers->setBinLabel(ib++,"CSC 4mu",1);
    ib++;
    regional_triggers->setBinLabel(ib++,"RPCf 1mu",1);
    regional_triggers->setBinLabel(ib++,"RPCf 2mu",1);
    regional_triggers->setBinLabel(ib++,"RPCf 3mu",1);
    regional_triggers->setBinLabel(ib++,"RPCf 4mu",1);
    ib++;
    regional_triggers->setBinLabel(ib++,"DT & RPC",1);
    regional_triggers->setBinLabel(ib++,"DT & CSC",1);
    regional_triggers->setBinLabel(ib++,"CSC & RPC",1);
    regional_triggers->setBinLabel(ib++,"DT & CSC & RPC",1);

    
    bx_number = dbe->book1D("Bx_Number","Bx number ROP chip", 3564, 0., 3564.);
    bx_number->setAxisTitle("bx number",1);
    
    dbx_chip = dbe->bookProfile("dbx_Chip","bx count difference wrt ROP chip", 5, 0., 5.,100,-4000.,4000.,"i");
    dbx_chip->setAxisTitle("chip name",1);
    dbx_chip->setAxisTitle("delta bx",2);
    dbx_chip->setBinLabel(1,"IND",1);
    dbx_chip->setBinLabel(2,"INB",1);
    dbx_chip->setBinLabel(3,"INC",1);
    dbx_chip->setBinLabel(4,"INF",1);
    dbx_chip->setBinLabel(5,"SRT",1);
    
    eta_dtcsc_and_rpc = dbe->book1D("eta_DTCSC_and_RPC","eta of confirmed GMT candidates",
        netascale[GMT], etascale[GMT]);
    eta_dtcsc_and_rpc->setAxisTitle("eta",1);
    
    eta_dtcsc_only = dbe->book1D("eta_DTCSC_only","eta of unconfirmed DT/CSC candidates",
        netascale[GMT], etascale[GMT]);
    eta_dtcsc_only->setAxisTitle("eta",1);
    
    eta_rpc_only = dbe->book1D("eta_RPC_only","eta of unconfirmed RPC candidates",
        netascale[GMT], etascale[GMT]);
    eta_rpc_only->setAxisTitle("eta",1);
    
    phi_dtcsc_and_rpc = dbe->book1D("phi_DTCSC_and_RPC","phi of confirmed GMT candidates",
        nphiscale, phiscale);
    phi_dtcsc_and_rpc->setAxisTitle("phi (deg)",1);
    
    phi_dtcsc_only = dbe->book1D("phi_DTCSC_only","phi of unconfirmed DT/CSC candidates",
        nphiscale, phiscale);
    phi_dtcsc_only->setAxisTitle("phi (deg)",1);
    
    phi_rpc_only = dbe->book1D("phi_RPC_only","phi of unconfirmed RPC candidates",
        nphiscale, phiscale);
    phi_rpc_only->setAxisTitle("phi (deg)",1);
    
    etaphi_dtcsc_and_rpc = dbe->book2D("etaphi_DTCSC_and_RPC","eta vs phi map of confirmed GMT candidates",
        netascale[GMT], etascale[GMT], nphiscale, phiscale);
    etaphi_dtcsc_and_rpc->setAxisTitle("eta",1);
    etaphi_dtcsc_and_rpc->setAxisTitle("phi (deg)",2);
    
    etaphi_dtcsc_only = dbe->book2D("etaphi_DTCSC_only","eta vs phi map of unconfirmed DT/CSC candidates",
        netascale[GMT], etascale[GMT], nphiscale, phiscale);
    etaphi_dtcsc_only->setAxisTitle("eta",1);
    etaphi_dtcsc_only->setAxisTitle("phi (deg)",2);
    
    etaphi_rpc_only = dbe->book2D("etaphi_RPC_only","eta vs phi map of unconfirmed RPC candidates",
        netascale[GMT], etascale[GMT], nphiscale, phiscale);
    etaphi_rpc_only->setAxisTitle("eta",1);
    etaphi_rpc_only->setAxisTitle("phi (deg)",2);
    
    
    dist_phi_dt_rpc = dbe->book1D("dist_phi_DT_RPC","Dphi between DT and RPC candidates", 100, -125., 125.);
    dist_phi_dt_rpc->setAxisTitle("delta phi (deg)",1);

    dist_phi_csc_rpc = dbe->book1D("dist_phi_CSC_RPC","Dphi between CSC and RPC candidates", 100, -125., 125.);
    dist_phi_csc_rpc->setAxisTitle("delta phi (deg)",1);

    dist_phi_dt_csc = dbe->book1D("dist_phi_DT_CSC","Dphi between DT and CSC candidates", 100, -125., 125.);
    dist_phi_dt_csc->setAxisTitle("delta phi (deg)",1);
       

    dist_eta_dt_rpc = dbe->book1D("dist_eta_DT_RPC","Deta between DT and RPC candidates", 40, -1., 1.);
    dist_eta_dt_rpc->setAxisTitle("delta eta",1);

    dist_eta_csc_rpc = dbe->book1D("dist_eta_CSC_RPC","Deta between CSC and RPC candidates", 40, -1., 1.);
    dist_eta_csc_rpc->setAxisTitle("delta eta",1);

    dist_eta_dt_csc = dbe->book1D("dist_eta_DT_CSC","Deta between DT and CSC candidates", 40, -1., 1.);
    dist_eta_dt_csc->setAxisTitle("delta eta",1);

       
    n_rpcb_vs_dttf  = dbe->book2D("n_RPCb_vs_DTTF",  "n cands RPCb vs DTTF",  5, -0.5, 4.5, 5, -0.5, 4.5);
    n_rpcb_vs_dttf->setAxisTitle("DTTF candidates",1);
    n_rpcb_vs_dttf->setAxisTitle("barrel RPC candidates",2);
    
    n_rpcf_vs_csctf = dbe->book2D("n_RPCf_vs_CSCTF", "n cands RPCf vs CSCTF", 5, -0.5, 4.5, 5, -0.5, 4.5);
    n_rpcf_vs_csctf->setAxisTitle("CSCTF candidates",1);
    n_rpcf_vs_csctf->setAxisTitle("endcap RPC candidates",2);
    
    n_csctf_vs_dttf = dbe->book2D("n_CSCTF_vs_DTTF", "n cands CSCTF vs DTTF", 5, -0.5, 4.5, 5, -0.5, 4.5);
    n_csctf_vs_dttf->setAxisTitle("DTTF candidates",1);
    n_csctf_vs_dttf->setAxisTitle("CSCTF candidates",2);
    
    bx_dt_rpc  = dbe->book2D("bx_DT_vs_RPC",  "1st bx DT vs. RPC",  5, -2.5, 2.5, 5, -2.5, 2.5);
    bx_dt_rpc->setAxisTitle("bx of 1st DTTF candidate",1);
    bx_dt_rpc->setAxisTitle("bx of 1st RPCb candidate",2);
    
    bx_csc_rpc  = dbe->book2D("bx_CSC_vs_RPC",  "1st bx CSC vs. RPC",  5, -2.5, 2.5, 5, -2.5, 2.5);
    bx_csc_rpc->setAxisTitle("bx of 1st CSCTF candidate",1);
    bx_csc_rpc->setAxisTitle("bx of 1st RPCf candidate",2);
    
    bx_dt_csc  = dbe->book2D("bx_DT_vs_CSC",  "1st bx DT vs. CSC",  5, -2.5, 2.5, 5, -2.5, 2.5);
    bx_dt_csc->setAxisTitle("bx of 1st DTTF candidate",1);
    bx_dt_csc->setAxisTitle("bx of 1st CSCTF candidate",2);
    
    
    for(int i=0; i<4; i++) {
      hname = subs[i] + "_dbx"; htitle = "dBx " + subs[i] + " to previous event";
      subs_dbx[i] = dbe->book2D(hname.data(),htitle.data(), 1000, 0., 1000., 4, 0., 4.);
      for(int j=0; j<4; j++) {
        subs_dbx[i]->setBinLabel((j+1),subs[j].data(),2);
      }
    }        
  }  
}
void L1TGMT::endJob ( void  ) [protected, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 82 of file L1TGMT.cc.

References gather_cfg::cout, dbe, nev_, outputFile_, DQMStore::save(), and verbose_.

{
  if(verbose_) cout << "L1TGMT: end job...." << endl;
  LogInfo("EndJob") << "analyzed " << nev_ << " events"; 

 if ( outputFile_.size() != 0  && dbe ) dbe->save(outputFile_);

 return;
}
double L1TGMT::phiconv_ ( float  phi) [private]

Definition at line 284 of file L1TGMT.cc.

References piconv_.

Referenced by analyze().

                                 {
  double phiout = double(phi);
  phiout *= piconv_;
  phiout += 0.001; // add a small value to get off the bin edge
  return phiout;
}

Member Data Documentation

Definition at line 102 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 103 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 101 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 84 of file L1TGMT.h.

Referenced by analyze(), and book_().

int L1TGMT::bxnum_old_ [private]

Definition at line 119 of file L1TGMT.h.

Referenced by analyze(), and beginJob().

DQMStore* L1TGMT::dbe [private]

Definition at line 69 of file L1TGMT.h.

Referenced by book_(), endJob(), and L1TGMT().

Definition at line 85 of file L1TGMT.h.

Referenced by book_().

Definition at line 99 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 100 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 98 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 96 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 97 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 95 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 86 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 87 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 88 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 92 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 93 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 94 of file L1TGMT.h.

Referenced by analyze(), and book_().

int L1TGMT::evnum_old_ [private]

Definition at line 118 of file L1TGMT.h.

Referenced by analyze(), and beginJob().

Definition at line 116 of file L1TGMT.h.

Referenced by analyze().

ofstream L1TGMT::logFile_ [private]

Definition at line 115 of file L1TGMT.h.

bool L1TGMT::monitorDaemon_ [private]

Definition at line 114 of file L1TGMT.h.

Definition at line 107 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 105 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 106 of file L1TGMT.h.

Referenced by analyze(), and book_().

int L1TGMT::nev_ [private]

Definition at line 111 of file L1TGMT.h.

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

int L1TGMT::obnum_old_ [private]

Definition at line 120 of file L1TGMT.h.

Referenced by analyze().

std::string L1TGMT::outputFile_ [private]

Definition at line 112 of file L1TGMT.h.

Referenced by endJob(), and L1TGMT().

Definition at line 89 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 90 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 91 of file L1TGMT.h.

Referenced by analyze(), and book_().

const double L1TGMT::piconv_ = 180. / acos(-1.) [static, private]

Definition at line 123 of file L1TGMT.h.

Referenced by book_(), and phiconv_().

Definition at line 82 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 80 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 109 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 74 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 78 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 79 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 73 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 75 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 76 of file L1TGMT.h.

Referenced by analyze(), and book_().

Definition at line 77 of file L1TGMT.h.

Referenced by analyze(), and book_().

int L1TGMT::trsrc_old_ [private]

Definition at line 121 of file L1TGMT.h.

Referenced by analyze().

bool L1TGMT::verbose_ [private]

Definition at line 113 of file L1TGMT.h.

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