Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:07:52

0001 /*
0002  * \file L1TGMT.cc
0003  *
0004  * \author J. Berryhill, I. Mikulec
0005  *
0006  */
0007 
0008 #include "DQM/L1TMonitor/interface/L1TGMT.h"
0009 #include "DQMServices/Core/interface/DQMStore.h"
0010 
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/ESHandle.h"
0013 #include "CondFormats/L1TObjects/interface/L1MuTriggerScales.h"
0014 #include "CondFormats/DataRecord/interface/L1MuTriggerScalesRcd.h"
0015 #include "CondFormats/L1TObjects/interface/L1MuTriggerPtScale.h"
0016 #include "CondFormats/DataRecord/interface/L1MuTriggerPtScaleRcd.h"
0017 
0018 using namespace std;
0019 using namespace edm;
0020 
0021 const double L1TGMT::piconv_ = 180. / acos(-1.);
0022 
0023 L1TGMT::L1TGMT(const ParameterSet& ps)
0024     : verbose_(ps.getUntrackedParameter<bool>("verbose", false))  // verbosity switch
0025       ,
0026       gmtSource_(consumes<L1MuGMTReadoutCollection>(ps.getParameter<InputTag>("gmtSource"))),
0027       bxnum_old_(0),
0028       obnum_old_(0),
0029       trsrc_old_(0) {
0030   if (verbose_)
0031     cout << "L1TGMT: constructor...." << endl;
0032   l1muTrigscaleToken_ = esConsumes<edm::Transition::BeginRun>();
0033   l1TrigptscaleToken_ = esConsumes<edm::Transition::BeginRun>();
0034 }
0035 
0036 L1TGMT::~L1TGMT() {}
0037 
0038 void L1TGMT::analyze(const Event& e, const EventSetup& c) {
0039   if (verbose_)
0040     cout << "L1TGMT: analyze...." << endl;
0041 
0042   edm::Handle<L1MuGMTReadoutCollection> pCollection;
0043   e.getByToken(gmtSource_, pCollection);
0044 
0045   if (!pCollection.isValid()) {
0046     edm::LogInfo("DataNotFound") << "can't find L1MuGMTReadoutCollection";
0047     return;
0048   }
0049 
0050   // remember the bx of 1st candidate of each system (9=none)
0051   int bx1st[4] = {9, 9, 9, 9};
0052 
0053   // get GMT readout collection
0054   L1MuGMTReadoutCollection const* gmtrc = pCollection.product();
0055   // get record vector
0056   vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords();
0057   // loop over records of individual bx's
0058   vector<L1MuGMTReadoutRecord>::const_iterator RRItr;
0059 
0060   for (RRItr = gmt_records.begin(); RRItr != gmt_records.end(); RRItr++) {
0061     vector<L1MuRegionalCand> INPCands[4] = {
0062         RRItr->getDTBXCands(), RRItr->getBrlRPCCands(), RRItr->getCSCCands(), RRItr->getFwdRPCCands()};
0063     vector<L1MuGMTExtendedCand> GMTCands = RRItr->getGMTCands();
0064 
0065     vector<L1MuRegionalCand>::const_iterator INPItr;
0066     vector<L1MuGMTExtendedCand>::const_iterator GMTItr;
0067     vector<L1MuGMTExtendedCand>::const_iterator GMTItr2;
0068 
0069     int BxInEvent = RRItr->getBxInEvent();
0070 
0071     // count non-empty candidates in this bx
0072     int nSUBS[5] = {0, 0, 0, 0, 0};
0073     for (int i = 0; i < 4; i++) {
0074       for (INPItr = INPCands[i].begin(); INPItr != INPCands[i].end(); ++INPItr) {
0075         if (!INPItr->empty()) {
0076           nSUBS[i]++;
0077           if (bx1st[i] == 9)
0078             bx1st[i] = BxInEvent;
0079         }
0080       }
0081       subs_nbx[i]->Fill(float(nSUBS[i]), float(BxInEvent));
0082     }
0083 
0084     for (GMTItr = GMTCands.begin(); GMTItr != GMTCands.end(); ++GMTItr) {
0085       if (!GMTItr->empty())
0086         nSUBS[GMT]++;
0087     }
0088     subs_nbx[GMT]->Fill(float(nSUBS[GMT]), float(BxInEvent));
0089 
0090     ////////////////////////////////////////////////////////////////////////////////////////////
0091     // from here care only about the L1A bunch crossing
0092     if (BxInEvent != 0)
0093       continue;
0094 
0095     // get the absolute bx number of the L1A
0096     int Bx = RRItr->getBxNr();
0097 
0098     bx_number->Fill(double(Bx));
0099 
0100     for (int i = 0; i < 4; i++) {
0101       for (INPItr = INPCands[i].begin(); INPItr != INPCands[i].end(); ++INPItr) {
0102         if (INPItr->empty())
0103           continue;
0104         subs_eta[i]->Fill(INPItr->etaValue());
0105         subs_phi[i]->Fill(phiconv_(INPItr->phiValue()));
0106         subs_pt[i]->Fill(INPItr->ptValue());
0107         subs_qty[i]->Fill(INPItr->quality());
0108         subs_etaphi[i]->Fill(INPItr->etaValue(), phiconv_(INPItr->phiValue()));
0109         subs_etaqty[i]->Fill(INPItr->etaValue(), INPItr->quality());
0110         int word = INPItr->getDataWord();
0111         for (int j = 0; j < 32; j++) {
0112           if (word & (1 << j))
0113             subs_bits[i]->Fill(float(j));
0114         }
0115       }
0116     }
0117 
0118     for (GMTItr = GMTCands.begin(); GMTItr != GMTCands.end(); ++GMTItr) {
0119       if (GMTItr->empty())
0120         continue;
0121       subs_eta[GMT]->Fill(GMTItr->etaValue());
0122       subs_phi[GMT]->Fill(phiconv_(GMTItr->phiValue()));
0123       subs_pt[GMT]->Fill(GMTItr->ptValue());
0124       subs_qty[GMT]->Fill(GMTItr->quality());
0125       subs_etaphi[GMT]->Fill(GMTItr->etaValue(), phiconv_(GMTItr->phiValue()));
0126       subs_etaqty[GMT]->Fill(GMTItr->etaValue(), GMTItr->quality());
0127       int word = GMTItr->getDataWord();
0128       for (int j = 0; j < 32; j++) {
0129         if (word & (1 << j))
0130           subs_bits[GMT]->Fill(float(j));
0131       }
0132 
0133       if (GMTItr->isMatchedCand()) {
0134         if (GMTItr->quality() > 3) {
0135           eta_dtcsc_and_rpc->Fill(GMTItr->etaValue());
0136           phi_dtcsc_and_rpc->Fill(phiconv_(GMTItr->phiValue()));
0137           etaphi_dtcsc_and_rpc->Fill(GMTItr->etaValue(), phiconv_(GMTItr->phiValue()));
0138         }
0139       } else if (GMTItr->isRPC()) {
0140         if (GMTItr->quality() > 3) {
0141           eta_rpc_only->Fill(GMTItr->etaValue());
0142           phi_rpc_only->Fill(phiconv_(GMTItr->phiValue()));
0143           etaphi_rpc_only->Fill(GMTItr->etaValue(), phiconv_(GMTItr->phiValue()));
0144         }
0145       } else {
0146         if (GMTItr->quality() > 3) {
0147           eta_dtcsc_only->Fill(GMTItr->etaValue());
0148           phi_dtcsc_only->Fill(phiconv_(GMTItr->phiValue()));
0149           etaphi_dtcsc_only->Fill(GMTItr->etaValue(), phiconv_(GMTItr->phiValue()));
0150         }
0151 
0152         if (GMTItr != GMTCands.end()) {
0153           for (GMTItr2 = GMTCands.begin(); GMTItr2 != GMTCands.end(); ++GMTItr2) {
0154             if (GMTItr2 == GMTItr)
0155               continue;
0156             if (GMTItr2->empty())
0157               continue;
0158             if (GMTItr2->isRPC()) {
0159               if (GMTItr->isFwd()) {
0160                 dist_eta_csc_rpc->Fill(GMTItr->etaValue() - GMTItr2->etaValue());
0161                 dist_phi_csc_rpc->Fill(phiconv_(GMTItr->phiValue()) - phiconv_(GMTItr2->phiValue()));
0162               } else {
0163                 dist_eta_dt_rpc->Fill(GMTItr->etaValue() - GMTItr2->etaValue());
0164                 dist_phi_dt_rpc->Fill(phiconv_(GMTItr->phiValue()) - phiconv_(GMTItr2->phiValue()));
0165               }
0166             } else {
0167               if (!(GMTItr->isFwd()) && GMTItr2->isFwd()) {
0168                 dist_eta_dt_csc->Fill(GMTItr->etaValue() - GMTItr2->etaValue());
0169                 dist_phi_dt_csc->Fill(phiconv_(GMTItr->phiValue()) - phiconv_(GMTItr2->phiValue()));
0170               } else if (GMTItr->isFwd() && !(GMTItr2->isFwd())) {
0171                 dist_eta_dt_csc->Fill(GMTItr2->etaValue() - GMTItr->etaValue());
0172                 dist_phi_dt_csc->Fill(phiconv_(GMTItr->phiValue()) - phiconv_(GMTItr2->phiValue()));
0173               }
0174             }
0175           }
0176         }
0177       }
0178     }
0179 
0180     n_rpcb_vs_dttf->Fill(float(nSUBS[DTTF]), float(nSUBS[RPCb]));
0181     n_rpcf_vs_csctf->Fill(float(nSUBS[CSCTF]), float(nSUBS[RPCf]));
0182     n_csctf_vs_dttf->Fill(float(nSUBS[DTTF]), float(nSUBS[CSCTF]));
0183 
0184     regional_triggers->Fill(-1.);  // fill underflow for normalization
0185     if (nSUBS[GMT])
0186       regional_triggers->Fill(0.);  // fill all muon bin
0187     int ioff = 1;
0188     for (int i = 0; i < 4; i++) {
0189       if (nSUBS[i])
0190         regional_triggers->Fill(float(5 * i + nSUBS[i] + ioff));
0191     }
0192     if (nSUBS[DTTF] && (nSUBS[RPCb] || nSUBS[RPCf]))
0193       regional_triggers->Fill(22.);
0194     if (nSUBS[DTTF] && nSUBS[CSCTF])
0195       regional_triggers->Fill(23.);
0196     if (nSUBS[CSCTF] && (nSUBS[RPCb] || nSUBS[RPCf]))
0197       regional_triggers->Fill(24.);
0198     if (nSUBS[DTTF] && nSUBS[CSCTF] && (nSUBS[RPCb] || nSUBS[RPCf]))
0199       regional_triggers->Fill(25.);
0200 
0201     // fill only if previous event corresponds to previous trigger
0202     //    if( (Ev - evnum_old_) == 1 && bxnum_old_ > -1 ) {
0203     // assume getting all events in a sequence (usefull only from reco data)
0204     if (bxnum_old_ > -1) {
0205       float dBx = Bx - bxnum_old_ + 3564.0 * (e.orbitNumber() - obnum_old_);
0206       for (int id = 0; id < 4; id++) {
0207         if (trsrc_old_ & (1 << id)) {
0208           for (int i = 0; i < 4; i++) {
0209             if (nSUBS[i])
0210               subs_dbx[i]->Fill(dBx, float(id));
0211           }
0212         }
0213       }
0214     }
0215 
0216     // save quantities for the next event
0217     bxnum_old_ = Bx;
0218     obnum_old_ = e.orbitNumber();
0219     trsrc_old_ = 0;
0220     for (int i = 0; i < 4; i++) {
0221       if (nSUBS[i])
0222         trsrc_old_ |= (1 << i);
0223     }
0224   }
0225 
0226   if (bx1st[DTTF] < 9 && bx1st[RPCb] < 9)
0227     bx_dt_rpc->Fill(bx1st[DTTF], bx1st[RPCb]);
0228   if (bx1st[CSCTF] < 9 && bx1st[RPCf] < 9)
0229     bx_csc_rpc->Fill(bx1st[CSCTF], bx1st[RPCf]);
0230   if (bx1st[DTTF] < 9 && bx1st[CSCTF] < 9)
0231     bx_dt_csc->Fill(bx1st[DTTF], bx1st[CSCTF]);
0232 }
0233 
0234 double L1TGMT::phiconv_(float phi) {
0235   double phiout = double(phi);
0236   phiout *= piconv_;
0237   phiout += 0.001;  // add a small value to get off the bin edge
0238   return phiout;
0239 }
0240 
0241 void L1TGMT::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const&, edm::EventSetup const& c) {
0242   std::string subs[5] = {"DTTF", "RPCb", "CSCTF", "RPCf", "GMT"};
0243 
0244   const L1MuTriggerScales* scales = &c.getData(l1muTrigscaleToken_);
0245   const L1MuTriggerPtScale* scalept = &c.getData(l1TrigptscaleToken_);
0246 
0247   ibooker.setCurrentFolder("L1T/L1TGMT");
0248 
0249   int nqty = 8;
0250   double qtymin = -0.5;
0251   double qtymax = 7.5;
0252 
0253   float phiscale[145];
0254   int nphiscale;
0255   {
0256     int nbins = scales->getPhiScale()->getNBins();
0257     if (nbins > 144)
0258       nbins = 144;
0259     for (int j = 0; j <= nbins; j++) {
0260       phiscale[j] = piconv_ * scales->getPhiScale()->getValue(j);
0261     }
0262     nphiscale = nbins;
0263   }
0264 
0265   float qscale[9];
0266   {
0267     for (int j = 0; j < 9; j++) {
0268       qscale[j] = -0.5 + j;
0269     }
0270   }
0271 
0272   // pt scale first bin reserved for empty muon
0273   float ptscale[32];
0274   int nptscale;
0275   {
0276     int nbins = scalept->getPtScale()->getNBins() - 1;
0277     if (nbins > 31)
0278       nbins = 31;
0279     for (int j = 1; j <= nbins; j++) {
0280       ptscale[j - 1] = scalept->getPtScale()->getValue(j);
0281     }
0282     ptscale[nbins] = ptscale[nbins - 1] + 10.;  // make reasonable size last bin
0283     nptscale = nbins;
0284   }
0285 
0286   float etascale[5][66];
0287   int netascale[5];
0288   // DTTF eta scale
0289   {
0290     int nbins = scales->getRegionalEtaScale(DTTF)->getNBins();
0291     if (nbins > 65)
0292       nbins = 65;
0293     for (int j = 0; j <= nbins; j++) {
0294       etascale[DTTF][j] = scales->getRegionalEtaScale(DTTF)->getValue(j);
0295     }
0296     netascale[DTTF] = nbins;
0297   }
0298   // RPCb etascale
0299   {
0300     int nbins = scales->getRegionalEtaScale(RPCb)->getNBins();
0301     if (nbins > 65)
0302       nbins = 65;
0303     for (int j = 0; j <= nbins; j++) {
0304       etascale[RPCb][j] = scales->getRegionalEtaScale(RPCb)->getValue(j);
0305     }
0306     netascale[RPCb] = nbins;
0307   }
0308   // CSCTF etascale
0309   // special case - need to mirror 2*32 bins
0310   {
0311     int nbins = scales->getRegionalEtaScale(CSCTF)->getNBins();
0312     if (nbins > 32)
0313       nbins = 32;
0314 
0315     int i = 0;
0316     for (int j = nbins; j >= 0; j--, i++) {
0317       etascale[CSCTF][i] = (-1) * scales->getRegionalEtaScale(CSCTF)->getValue(j);
0318     }
0319     for (int j = 0; j <= nbins; j++, i++) {
0320       etascale[CSCTF][i] = scales->getRegionalEtaScale(CSCTF)->getValue(j);
0321     }
0322     netascale[CSCTF] = i - 1;
0323   }
0324   // RPCf etascale
0325   {
0326     int nbins = scales->getRegionalEtaScale(RPCf)->getNBins();
0327     if (nbins > 65)
0328       nbins = 65;
0329     for (int j = 0; j <= nbins; j++) {
0330       etascale[RPCf][j] = scales->getRegionalEtaScale(RPCf)->getValue(j);
0331     }
0332     netascale[RPCf] = nbins;
0333   }
0334   // GMT etascale
0335   {
0336     int nbins = scales->getGMTEtaScale()->getNBins();
0337     if (nbins > 32)
0338       nbins = 32;
0339 
0340     int i = 0;
0341     for (int j = nbins; j > 0; j--, i++) {
0342       etascale[GMT][i] = (-1) * scales->getGMTEtaScale()->getValue(j);
0343     }
0344     for (int j = 0; j <= nbins; j++, i++) {
0345       etascale[GMT][i] = scales->getGMTEtaScale()->getValue(j);
0346     }
0347     netascale[GMT] = i - 1;
0348   }
0349 
0350   std::string hname("");
0351   std::string htitle("");
0352 
0353   for (int i = 0; i < 5; i++) {
0354     hname = subs[i] + "_nbx";
0355     htitle = subs[i] + " multiplicity in bx";
0356     subs_nbx[i] = ibooker.book2D(hname.data(), htitle.data(), 4, 1., 5., 5, -2.5, 2.5);
0357     subs_nbx[i]->setAxisTitle(subs[i] + " candidates", 1);
0358     subs_nbx[i]->setAxisTitle("bx wrt L1A", 2);
0359 
0360     hname = subs[i] + "_eta";
0361     htitle = subs[i] + " eta value";
0362     subs_eta[i] = ibooker.book1D(hname.data(), htitle.data(), netascale[i], etascale[i]);
0363     subs_eta[i]->setAxisTitle("eta", 1);
0364 
0365     hname = subs[i] + "_phi";
0366     htitle = subs[i] + " phi value";
0367     subs_phi[i] = ibooker.book1D(hname.data(), htitle.data(), nphiscale, phiscale);
0368     subs_phi[i]->setAxisTitle("phi (deg)", 1);
0369 
0370     hname = subs[i] + "_pt";
0371     htitle = subs[i] + " pt value";
0372     subs_pt[i] = ibooker.book1D(hname.data(), htitle.data(), nptscale, ptscale);
0373     subs_pt[i]->setAxisTitle("L1 pT (GeV)", 1);
0374 
0375     hname = subs[i] + "_qty";
0376     htitle = subs[i] + " qty value";
0377     subs_qty[i] = ibooker.book1D(hname.data(), htitle.data(), nqty, qtymin, qtymax);
0378     subs_qty[i]->setAxisTitle(subs[i] + " quality", 1);
0379 
0380     hname = subs[i] + "_etaphi";
0381     htitle = subs[i] + " phi vs eta";
0382     subs_etaphi[i] = ibooker.book2D(hname.data(), htitle.data(), netascale[i], etascale[i], nphiscale, phiscale);
0383     subs_etaphi[i]->setAxisTitle("eta", 1);
0384     subs_etaphi[i]->setAxisTitle("phi (deg)", 2);
0385 
0386     hname = subs[i] + "_etaqty";
0387     htitle = subs[i] + " qty vs eta";
0388     subs_etaqty[i] = ibooker.book2D(hname.data(), htitle.data(), netascale[i], etascale[i], nqty, qscale);
0389     subs_etaqty[i]->setAxisTitle("eta", 1);
0390     subs_etaqty[i]->setAxisTitle(subs[i] + " quality", 2);
0391 
0392     hname = subs[i] + "_bits";
0393     htitle = subs[i] + " bit population";
0394     subs_bits[i] = ibooker.book1D(hname.data(), htitle.data(), 32, -0.5, 31.5);
0395     subs_bits[i]->setAxisTitle("bit number", 1);
0396   }
0397 
0398   regional_triggers = ibooker.book1D("Regional_trigger", "Muon trigger contribution", 27, 0., 27.);
0399   regional_triggers->setAxisTitle("regional trigger", 1);
0400   int ib = 1;
0401   regional_triggers->setBinLabel(ib++, "All muons", 1);
0402   ib++;
0403   regional_triggers->setBinLabel(ib++, "DT 1mu", 1);
0404   regional_triggers->setBinLabel(ib++, "DT 2mu", 1);
0405   regional_triggers->setBinLabel(ib++, "DT 3mu", 1);
0406   regional_triggers->setBinLabel(ib++, "DT 4mu", 1);
0407   ib++;
0408   regional_triggers->setBinLabel(ib++, "RPCb 1mu", 1);
0409   regional_triggers->setBinLabel(ib++, "RPCb 2mu", 1);
0410   regional_triggers->setBinLabel(ib++, "RPCb 3mu", 1);
0411   regional_triggers->setBinLabel(ib++, "RPCb 4mu", 1);
0412   ib++;
0413   regional_triggers->setBinLabel(ib++, "CSC 1mu", 1);
0414   regional_triggers->setBinLabel(ib++, "CSC 2mu", 1);
0415   regional_triggers->setBinLabel(ib++, "CSC 3mu", 1);
0416   regional_triggers->setBinLabel(ib++, "CSC 4mu", 1);
0417   ib++;
0418   regional_triggers->setBinLabel(ib++, "RPCf 1mu", 1);
0419   regional_triggers->setBinLabel(ib++, "RPCf 2mu", 1);
0420   regional_triggers->setBinLabel(ib++, "RPCf 3mu", 1);
0421   regional_triggers->setBinLabel(ib++, "RPCf 4mu", 1);
0422   ib++;
0423   regional_triggers->setBinLabel(ib++, "DT & RPC", 1);
0424   regional_triggers->setBinLabel(ib++, "DT & CSC", 1);
0425   regional_triggers->setBinLabel(ib++, "CSC & RPC", 1);
0426   regional_triggers->setBinLabel(ib++, "DT & CSC & RPC", 1);
0427 
0428   bx_number = ibooker.book1D("Bx_Number", "Bx number ROP chip", 3564, 0., 3564.);
0429   bx_number->setAxisTitle("bx number", 1);
0430 
0431   dbx_chip = ibooker.bookProfile("dbx_Chip", "bx count difference wrt ROP chip", 5, 0., 5., 100, -4000., 4000., "i");
0432   dbx_chip->setAxisTitle("chip name", 1);
0433   dbx_chip->setAxisTitle("delta bx", 2);
0434   dbx_chip->setBinLabel(1, "IND", 1);
0435   dbx_chip->setBinLabel(2, "INB", 1);
0436   dbx_chip->setBinLabel(3, "INC", 1);
0437   dbx_chip->setBinLabel(4, "INF", 1);
0438   dbx_chip->setBinLabel(5, "SRT", 1);
0439 
0440   eta_dtcsc_and_rpc =
0441       ibooker.book1D("eta_DTCSC_and_RPC", "eta of confirmed GMT candidates", netascale[GMT], etascale[GMT]);
0442   eta_dtcsc_and_rpc->setAxisTitle("eta", 1);
0443 
0444   eta_dtcsc_only =
0445       ibooker.book1D("eta_DTCSC_only", "eta of unconfirmed DT/CSC candidates", netascale[GMT], etascale[GMT]);
0446   eta_dtcsc_only->setAxisTitle("eta", 1);
0447 
0448   eta_rpc_only = ibooker.book1D("eta_RPC_only", "eta of unconfirmed RPC candidates", netascale[GMT], etascale[GMT]);
0449   eta_rpc_only->setAxisTitle("eta", 1);
0450 
0451   phi_dtcsc_and_rpc = ibooker.book1D("phi_DTCSC_and_RPC", "phi of confirmed GMT candidates", nphiscale, phiscale);
0452   phi_dtcsc_and_rpc->setAxisTitle("phi (deg)", 1);
0453 
0454   phi_dtcsc_only = ibooker.book1D("phi_DTCSC_only", "phi of unconfirmed DT/CSC candidates", nphiscale, phiscale);
0455   phi_dtcsc_only->setAxisTitle("phi (deg)", 1);
0456 
0457   phi_rpc_only = ibooker.book1D("phi_RPC_only", "phi of unconfirmed RPC candidates", nphiscale, phiscale);
0458   phi_rpc_only->setAxisTitle("phi (deg)", 1);
0459 
0460   etaphi_dtcsc_and_rpc = ibooker.book2D("etaphi_DTCSC_and_RPC",
0461                                         "eta vs phi map of confirmed GMT candidates",
0462                                         netascale[GMT],
0463                                         etascale[GMT],
0464                                         nphiscale,
0465                                         phiscale);
0466   etaphi_dtcsc_and_rpc->setAxisTitle("eta", 1);
0467   etaphi_dtcsc_and_rpc->setAxisTitle("phi (deg)", 2);
0468 
0469   etaphi_dtcsc_only = ibooker.book2D("etaphi_DTCSC_only",
0470                                      "eta vs phi map of unconfirmed DT/CSC candidates",
0471                                      netascale[GMT],
0472                                      etascale[GMT],
0473                                      nphiscale,
0474                                      phiscale);
0475   etaphi_dtcsc_only->setAxisTitle("eta", 1);
0476   etaphi_dtcsc_only->setAxisTitle("phi (deg)", 2);
0477 
0478   etaphi_rpc_only = ibooker.book2D("etaphi_RPC_only",
0479                                    "eta vs phi map of unconfirmed RPC candidates",
0480                                    netascale[GMT],
0481                                    etascale[GMT],
0482                                    nphiscale,
0483                                    phiscale);
0484   etaphi_rpc_only->setAxisTitle("eta", 1);
0485   etaphi_rpc_only->setAxisTitle("phi (deg)", 2);
0486 
0487   dist_phi_dt_rpc = ibooker.book1D("dist_phi_DT_RPC", "Dphi between DT and RPC candidates", 100, -125., 125.);
0488   dist_phi_dt_rpc->setAxisTitle("delta phi (deg)", 1);
0489 
0490   dist_phi_csc_rpc = ibooker.book1D("dist_phi_CSC_RPC", "Dphi between CSC and RPC candidates", 100, -125., 125.);
0491   dist_phi_csc_rpc->setAxisTitle("delta phi (deg)", 1);
0492 
0493   dist_phi_dt_csc = ibooker.book1D("dist_phi_DT_CSC", "Dphi between DT and CSC candidates", 100, -125., 125.);
0494   dist_phi_dt_csc->setAxisTitle("delta phi (deg)", 1);
0495 
0496   dist_eta_dt_rpc = ibooker.book1D("dist_eta_DT_RPC", "Deta between DT and RPC candidates", 40, -1., 1.);
0497   dist_eta_dt_rpc->setAxisTitle("delta eta", 1);
0498 
0499   dist_eta_csc_rpc = ibooker.book1D("dist_eta_CSC_RPC", "Deta between CSC and RPC candidates", 40, -1., 1.);
0500   dist_eta_csc_rpc->setAxisTitle("delta eta", 1);
0501 
0502   dist_eta_dt_csc = ibooker.book1D("dist_eta_DT_CSC", "Deta between DT and CSC candidates", 40, -1., 1.);
0503   dist_eta_dt_csc->setAxisTitle("delta eta", 1);
0504 
0505   n_rpcb_vs_dttf = ibooker.book2D("n_RPCb_vs_DTTF", "n cands RPCb vs DTTF", 5, -0.5, 4.5, 5, -0.5, 4.5);
0506   n_rpcb_vs_dttf->setAxisTitle("DTTF candidates", 1);
0507   n_rpcb_vs_dttf->setAxisTitle("barrel RPC candidates", 2);
0508 
0509   n_rpcf_vs_csctf = ibooker.book2D("n_RPCf_vs_CSCTF", "n cands RPCf vs CSCTF", 5, -0.5, 4.5, 5, -0.5, 4.5);
0510   n_rpcf_vs_csctf->setAxisTitle("CSCTF candidates", 1);
0511   n_rpcf_vs_csctf->setAxisTitle("endcap RPC candidates", 2);
0512 
0513   n_csctf_vs_dttf = ibooker.book2D("n_CSCTF_vs_DTTF", "n cands CSCTF vs DTTF", 5, -0.5, 4.5, 5, -0.5, 4.5);
0514   n_csctf_vs_dttf->setAxisTitle("DTTF candidates", 1);
0515   n_csctf_vs_dttf->setAxisTitle("CSCTF candidates", 2);
0516 
0517   bx_dt_rpc = ibooker.book2D("bx_DT_vs_RPC", "1st bx DT vs. RPC", 5, -2.5, 2.5, 5, -2.5, 2.5);
0518   bx_dt_rpc->setAxisTitle("bx of 1st DTTF candidate", 1);
0519   bx_dt_rpc->setAxisTitle("bx of 1st RPCb candidate", 2);
0520 
0521   bx_csc_rpc = ibooker.book2D("bx_CSC_vs_RPC", "1st bx CSC vs. RPC", 5, -2.5, 2.5, 5, -2.5, 2.5);
0522   bx_csc_rpc->setAxisTitle("bx of 1st CSCTF candidate", 1);
0523   bx_csc_rpc->setAxisTitle("bx of 1st RPCf candidate", 2);
0524 
0525   bx_dt_csc = ibooker.book2D("bx_DT_vs_CSC", "1st bx DT vs. CSC", 5, -2.5, 2.5, 5, -2.5, 2.5);
0526   bx_dt_csc->setAxisTitle("bx of 1st DTTF candidate", 1);
0527   bx_dt_csc->setAxisTitle("bx of 1st CSCTF candidate", 2);
0528 
0529   for (int i = 0; i < 4; i++) {
0530     hname = subs[i] + "_dbx";
0531     htitle = "dBx " + subs[i] + " to previous event";
0532     subs_dbx[i] = ibooker.book2D(hname.data(), htitle.data(), 1000, 0., 1000., 4, 0., 4.);
0533     for (int j = 0; j < 4; j++) {
0534       subs_dbx[i]->setBinLabel((j + 1), subs[j], 2);
0535     }
0536   }
0537 }