00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "DQM/L1TMonitor/interface/L1TGMT.h"
00011 #include "DQMServices/Core/interface/DQMStore.h"
00012
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/Framework/interface/ESHandle.h"
00015 #include "CondFormats/L1TObjects/interface/L1MuTriggerScales.h"
00016 #include "CondFormats/DataRecord/interface/L1MuTriggerScalesRcd.h"
00017 #include "CondFormats/L1TObjects/interface/L1MuTriggerPtScale.h"
00018 #include "CondFormats/DataRecord/interface/L1MuTriggerPtScaleRcd.h"
00019
00020 using namespace std;
00021 using namespace edm;
00022
00023 const double L1TGMT::piconv_ = 180. / acos(-1.);
00024
00025 L1TGMT::L1TGMT(const ParameterSet& ps)
00026 : gmtSource_( ps.getParameter< InputTag >("gmtSource") )
00027 {
00028
00029
00030 verbose_ = ps.getUntrackedParameter<bool>("verbose", false);
00031
00032 if(verbose_) cout << "L1TGMT: constructor...." << endl;
00033
00034
00035 dbe = NULL;
00036 if ( ps.getUntrackedParameter<bool>("DQMStore", false) )
00037 {
00038 dbe = Service<DQMStore>().operator->();
00039 dbe->setVerbose(0);
00040 }
00041
00042 outputFile_ = ps.getUntrackedParameter<string>("outputFile", "");
00043 if ( outputFile_.size() != 0 ) {
00044 cout << "L1T Monitoring histograms will be saved to " << outputFile_.c_str() << endl;
00045 }
00046
00047 bool disable = ps.getUntrackedParameter<bool>("disableROOToutput", false);
00048 if(disable){
00049 outputFile_="";
00050 }
00051
00052
00053 if ( dbe !=NULL ) {
00054 dbe->setCurrentFolder("L1T/L1TGMT");
00055 }
00056
00057 }
00058
00059 L1TGMT::~L1TGMT()
00060 {
00061 }
00062
00063 void L1TGMT::beginJob(const EventSetup& c)
00064 {
00065
00066 std::string subs[5] = { "DTTF", "RPCb", "CSCTF", "RPCf", "GMT" };
00067
00068 nev_ = 0;
00069 evnum_old_ = -1;
00070 bxnum_old_ = -1;
00071
00072 edm::ESHandle< L1MuTriggerScales > trigscales_h;
00073 c.get< L1MuTriggerScalesRcd >().get( trigscales_h );
00074 const L1MuTriggerScales* scales = trigscales_h.product();
00075
00076 edm::ESHandle< L1MuTriggerPtScale > trigptscale_h;
00077 c.get< L1MuTriggerPtScaleRcd >().get( trigptscale_h );
00078 const L1MuTriggerPtScale* scalept = trigptscale_h.product();
00079
00080
00081 DQMStore* dbe = 0;
00082 dbe = Service<DQMStore>().operator->();
00083
00084 if ( dbe ) {
00085 dbe->setCurrentFolder("L1T/L1TGMT");
00086 dbe->rmdir("L1T/L1TGMT");
00087 }
00088
00089
00090 if ( dbe )
00091 {
00092 dbe->setCurrentFolder("L1T/L1TGMT");
00093
00094 int nphi=144; double phimin= 0.; double phimax=360.;
00095 int nqty= 8; double qtymin=-0.5; double qtymax=7.5;
00096
00097 float phiscale[145];
00098 {
00099 for(int j=0; j<145; j++) {
00100 phiscale[j] = j*2.5 ;
00101 }
00102 }
00103
00104 float qscale[9];
00105 {
00106 for(int j=0; j<9; j++) {
00107 qscale[j] = -0.5 + j;
00108 }
00109 }
00110
00111 float ptscale[32];
00112 {
00113 int i=0;
00114 for(int j=1; j<32; j++,i++) {
00115 ptscale[i] = scalept->getPtScale()->getLowEdge(j);
00116 }
00117 ptscale[31]=ptscale[30]+10.;
00118 }
00119
00120 float etascale[5][66];
00121 int netascale[5];
00122
00123 {
00124 int i=0;
00125 for(int j=32; j<=63; j++,i++) {
00126 etascale[DTTF][i] = scales->getRegionalEtaScale(DTTF)->getLowEdge(j);
00127 }
00128 for(int j=0; j<=31; j++,i++) {
00129 etascale[DTTF][i] = scales->getRegionalEtaScale(DTTF)->getLowEdge(j);
00130 }
00131 etascale[DTTF][64] = scales->getRegionalEtaScale(DTTF)->getScaleMax();
00132 netascale[DTTF]=64;
00133 }
00134
00135 {
00136 int i=0;
00137 for(int j=48; j<=63; j++,i++) {
00138 etascale[RPCb][i] = scales->getRegionalEtaScale(RPCb)->getLowEdge(j);
00139 }
00140 for(int j=0; j<=16; j++,i++) {
00141 etascale[RPCb][i] = scales->getRegionalEtaScale(RPCb)->getLowEdge(j);
00142 }
00143 etascale[RPCb][33] = scales->getRegionalEtaScale(RPCb)->getScaleMax();
00144 netascale[RPCb]=33;
00145 }
00146
00147 {
00148 etascale[CSCTF][0] = (-1) * scales->getRegionalEtaScale(CSCTF)->getScaleMax();
00149 int i=1;
00150 for(int j=31; j>=0; j--,i++) {
00151 etascale[CSCTF][i] = (-1) * scales->getRegionalEtaScale(CSCTF)->getLowEdge(j);
00152 }
00153 for(int j=0; j<=31; j++,i++) {
00154 etascale[CSCTF][i] = scales->getRegionalEtaScale(CSCTF)->getLowEdge(j);
00155 }
00156 etascale[CSCTF][65] = scales->getRegionalEtaScale(CSCTF)->getScaleMax();
00157 netascale[CSCTF]=65;
00158 }
00159
00160 {
00161 int i=0;
00162 for(int j=48; j<=63; j++,i++) {
00163 etascale[RPCf][i] = scales->getRegionalEtaScale(RPCf)->getLowEdge(j);
00164 }
00165 for(int j=0; j<=16; j++,i++) {
00166 etascale[RPCf][i] = scales->getRegionalEtaScale(RPCf)->getLowEdge(j);
00167 }
00168 etascale[RPCf][33] = scales->getRegionalEtaScale(RPCf)->getScaleMax();
00169 netascale[RPCf]=33;
00170 }
00171
00172 {
00173 etascale[GMT][0] = (-1) * scales->getGMTEtaScale()->getScaleMax();
00174 int i=1;
00175 for(int j=30; j>0; j--,i++) {
00176 etascale[GMT][i] = (-1) * scales->getGMTEtaScale()->getLowEdge(j);
00177 }
00178 for(int j=0; j<=30; j++,i++) {
00179 etascale[GMT][i] = scales->getGMTEtaScale()->getLowEdge(j);
00180 }
00181 etascale[GMT][62] = scales->getGMTEtaScale()->getScaleMax();
00182 netascale[GMT]=62;
00183 }
00184
00185
00186 std::string hname("");
00187 std::string htitle("");
00188
00189 for(int i=0; i<5; i++) {
00190
00191 hname = subs[i] + "_nbx"; htitle = subs[i] + " multiplicity in bx";
00192 subs_nbx[i] = dbe->book2D(hname.data(),htitle.data(), 4, 1., 5., 5, -2.5, 2.5);
00193 subs_nbx[i]->setAxisTitle(subs[i] + " candidates",1);
00194 subs_nbx[i]->setAxisTitle("bx wrt L1A",2);
00195
00196 hname = subs[i] + "_eta"; htitle = subs[i] + " eta value";
00197 subs_eta[i] = dbe->book1D(hname.data(),htitle.data(), netascale[i], etascale[i]);
00198 subs_eta[i]->setAxisTitle("eta",1);
00199
00200 hname = subs[i] + "_phi"; htitle = subs[i] + " phi value";
00201 subs_phi[i] = dbe->book1D(hname.data(),htitle.data(), nphi, phimin, phimax);
00202 subs_phi[i]->setAxisTitle("phi (deg)",1);
00203
00204 hname = subs[i] + "_pt"; htitle = subs[i] + " pt value";
00205 subs_pt[i] = dbe->book1D(hname.data(),htitle.data(), 31, ptscale);
00206 subs_pt[i]->setAxisTitle("L1 pT (GeV)",1);
00207
00208 hname = subs[i] + "_qty"; htitle = subs[i] + " qty value";
00209 subs_qty[i] = dbe->book1D(hname.data(),htitle.data(), nqty, qtymin, qtymax);
00210 subs_qty[i]->setAxisTitle(subs[i] + " quality",1);
00211
00212 hname = subs[i] + "_etaphi"; htitle = subs[i] + " phi vs eta";
00213 subs_etaphi[i] = dbe->book2D(hname.data(),htitle.data(), netascale[i], etascale[i], 144, phiscale);
00214 subs_etaphi[i]->setAxisTitle("eta",1);
00215 subs_etaphi[i]->setAxisTitle("phi (deg)",2);
00216
00217 hname = subs[i] + "_etaqty"; htitle = subs[i] + " qty vs eta";
00218 subs_etaqty[i] = dbe->book2D(hname.data(),htitle.data(), netascale[i], etascale[i], 8, qscale);
00219 subs_etaqty[i]->setAxisTitle("eta",1);
00220 subs_etaqty[i]->setAxisTitle(subs[i] + " quality",2);
00221
00222 hname = subs[i] + "_bits"; htitle = subs[i] + " bit population";
00223 subs_bits[i] = dbe->book1D(hname.data(),htitle.data(), 32, -0.5, 31.5);
00224 subs_bits[i]->setAxisTitle("bit number",1);
00225
00226 hname = subs[i] + "_candlumi"; htitle = "number of " + subs[i] + " candidates per lumisegment";
00227 subs_candlumi[i] = dbe->book1D(hname.data(),htitle.data(), 250, 0., 250.);
00228 subs_candlumi[i]->setAxisTitle("luminosity segment number",1);
00229 }
00230
00231 regional_triggers = dbe->book1D("Regional_trigger","Muon trigger contribution", 4, 0., 4.);
00232 regional_triggers->setAxisTitle("regional trigger",1);
00233 regional_triggers->setBinLabel(1,"DTTF",1);
00234 regional_triggers->setBinLabel(2,"RPCb",1);
00235 regional_triggers->setBinLabel(3,"CSCTF",1);
00236 regional_triggers->setBinLabel(4,"RPCf",1);
00237
00238 bx_number = dbe->book1D("Bx_Number","Bx number ROP chip", 3564, 0., 3564.);
00239 bx_number->setAxisTitle("bx number",1);
00240
00241 dbx_chip = dbe->bookProfile("dbx_Chip","bx count difference wrt ROP chip", 5, 0., 5.,100,-4000.,4000.,"i");
00242 dbx_chip->setAxisTitle("chip name",1);
00243 dbx_chip->setAxisTitle("delta bx",2);
00244 dbx_chip->setBinLabel(1,"IND",1);
00245 dbx_chip->setBinLabel(2,"INB",1);
00246 dbx_chip->setBinLabel(3,"INC",1);
00247 dbx_chip->setBinLabel(4,"INF",1);
00248 dbx_chip->setBinLabel(5,"SRT",1);
00249
00250 eta_dtcsc_and_rpc = dbe->book1D("eta_DTCSC_and_RPC","eta of confirmed GMT candidates",
00251 netascale[GMT], etascale[GMT]);
00252 eta_dtcsc_and_rpc->setAxisTitle("eta",1);
00253
00254 eta_dtcsc_only = dbe->book1D("eta_DTCSC_only","eta of unconfirmed DT/CSC candidates",
00255 netascale[GMT], etascale[GMT]);
00256 eta_dtcsc_only->setAxisTitle("eta",1);
00257
00258 eta_rpc_only = dbe->book1D("eta_RPC_only","eta of unconfirmed RPC candidates",
00259 netascale[GMT], etascale[GMT]);
00260 eta_rpc_only->setAxisTitle("eta",1);
00261
00262 phi_dtcsc_and_rpc = dbe->book1D("phi_DTCSC_and_RPC","phi of confirmed GMT candidates",
00263 nphi, phimin, phimax);
00264 phi_dtcsc_and_rpc->setAxisTitle("phi (deg)",1);
00265
00266 phi_dtcsc_only = dbe->book1D("phi_DTCSC_only","phi of unconfirmed DT/CSC candidates",
00267 nphi, phimin, phimax);
00268 phi_dtcsc_only->setAxisTitle("phi (deg)",1);
00269
00270 phi_rpc_only = dbe->book1D("phi_RPC_only","phi of unconfirmed RPC candidates",
00271 nphi, phimin, phimax);
00272 phi_rpc_only->setAxisTitle("phi (deg)",1);
00273
00274 etaphi_dtcsc_and_rpc = dbe->book2D("etaphi_DTCSC_and_RPC","eta vs phi map of confirmed GMT candidates",
00275 100, -2.5, 2.5, nphi, phimin, phimax);
00276 etaphi_dtcsc_and_rpc->setAxisTitle("eta",1);
00277 etaphi_dtcsc_and_rpc->setAxisTitle("phi (deg)",2);
00278
00279 etaphi_dtcsc_only = dbe->book2D("etaphi_DTCSC_only","eta vs phi map of unconfirmed DT/CSC candidates",
00280 100, -2.5, 2.5, nphi, phimin, phimax);
00281 etaphi_dtcsc_only->setAxisTitle("eta",1);
00282 etaphi_dtcsc_only->setAxisTitle("phi (deg)",2);
00283
00284 etaphi_rpc_only = dbe->book2D("etaphi_RPC_only","eta vs phi map of unconfirmed RPC candidates",
00285 100, -2.5, 2.5, nphi, phimin, phimax);
00286 etaphi_rpc_only->setAxisTitle("eta",1);
00287 etaphi_rpc_only->setAxisTitle("phi (deg)",2);
00288
00289
00290 dist_phi_dt_rpc = dbe->book1D("dist_phi_DT_RPC","Dphi between DT and RPC candidates", 100, -125., 125.);
00291 dist_phi_dt_rpc->setAxisTitle("delta phi (deg)",1);
00292
00293 dist_phi_csc_rpc = dbe->book1D("dist_phi_CSC_RPC","Dphi between CSC and RPC candidates", 100, -125., 125.);
00294 dist_phi_csc_rpc->setAxisTitle("delta phi (deg)",1);
00295
00296 dist_phi_dt_csc = dbe->book1D("dist_phi_DT_CSC","Dphi between DT and CSC candidates", 100, -125., 125.);
00297 dist_phi_dt_csc->setAxisTitle("delta phi (deg)",1);
00298
00299
00300 dist_eta_dt_rpc = dbe->book1D("dist_eta_DT_RPC","Deta between DT and RPC candidates", 40, -1., 1.);
00301 dist_eta_dt_rpc->setAxisTitle("delta eta",1);
00302
00303 dist_eta_csc_rpc = dbe->book1D("dist_eta_CSC_RPC","Deta between CSC and RPC candidates", 40, -1., 1.);
00304 dist_eta_csc_rpc->setAxisTitle("delta eta",1);
00305
00306 dist_eta_dt_csc = dbe->book1D("dist_eta_DT_CSC","Deta between DT and CSC candidates", 40, -1., 1.);
00307 dist_eta_dt_csc->setAxisTitle("delta eta",1);
00308
00309
00310 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);
00311 n_rpcb_vs_dttf->setAxisTitle("DTTF candidates",1);
00312 n_rpcb_vs_dttf->setAxisTitle("barrel RPC candidates",2);
00313
00314 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);
00315 n_rpcf_vs_csctf->setAxisTitle("CSCTF candidates",1);
00316 n_rpcf_vs_csctf->setAxisTitle("endcap RPC candidates",2);
00317
00318 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);
00319 n_csctf_vs_dttf->setAxisTitle("DTTF candidates",1);
00320 n_csctf_vs_dttf->setAxisTitle("CSCTF candidates",2);
00321
00322
00323 for(int i=0; i<4; i++) {
00324 hname = subs[i] + "_dbx"; htitle = "dBx " + subs[i] + " to previous event";
00325 subs_dbx[i] = dbe->book2D(hname.data(),htitle.data(), 100, 0., 100., 4, 0., 4.);
00326 for(int j=0; j<4; j++) {
00327 subs_dbx[i]->setBinLabel((j+1),subs[j].data(),2);
00328 }
00329 }
00330 }
00331 }
00332
00333
00334 void L1TGMT::endJob(void)
00335 {
00336 if(verbose_) cout << "L1TGMT: end job...." << endl;
00337 LogInfo("EndJob") << "analyzed " << nev_ << " events";
00338
00339 if ( outputFile_.size() != 0 && dbe ) dbe->save(outputFile_);
00340
00341 return;
00342 }
00343
00344 void L1TGMT::analyze(const Event& e, const EventSetup& c)
00345 {
00346 nev_++;
00347 if(verbose_) cout << "L1TGMT: analyze...." << endl;
00348
00349
00350 edm::Handle<L1MuGMTReadoutCollection> pCollection;
00351 e.getByLabel(gmtSource_,pCollection);
00352
00353 if (!pCollection.isValid()) {
00354 edm::LogInfo("DataNotFound") << "can't find L1MuGMTReadoutCollection with label "
00355 << gmtSource_.label() ;
00356 return;
00357 }
00358
00359
00360 L1MuGMTReadoutCollection const* gmtrc = pCollection.product();
00361
00362 vector<L1MuGMTReadoutRecord> gmt_records = gmtrc->getRecords();
00363
00364 vector<L1MuGMTReadoutRecord>::const_iterator RRItr;
00365 for( RRItr = gmt_records.begin(); RRItr != gmt_records.end(); RRItr++ )
00366 {
00367
00368 vector<L1MuRegionalCand> INPCands[4] = {
00369 RRItr->getDTBXCands(),
00370 RRItr->getBrlRPCCands(),
00371 RRItr->getCSCCands(),
00372 RRItr->getFwdRPCCands()
00373 };
00374 vector<L1MuGMTExtendedCand> GMTCands = RRItr->getGMTCands();
00375
00376 vector<L1MuRegionalCand>::const_iterator INPItr;
00377 vector<L1MuGMTExtendedCand>::const_iterator GMTItr;
00378 vector<L1MuGMTExtendedCand>::const_iterator GMTItr2;
00379
00380 int BxInEvent = RRItr->getBxInEvent();
00381
00382
00383 int nSUBS[5] = {0, 0, 0, 0, 0};
00384 for(int i=0; i<4; i++) {
00385 for( INPItr = INPCands[i].begin(); INPItr != INPCands[i].end(); ++INPItr ) {
00386 if(!INPItr->empty()) nSUBS[i]++;
00387 }
00388 subs_nbx[i]->Fill(float(nSUBS[i]),float(BxInEvent));
00389 }
00390
00391 for( GMTItr = GMTCands.begin(); GMTItr != GMTCands.end(); ++GMTItr ) {
00392 if(!GMTItr->empty()) nSUBS[GMT]++;
00393 }
00394 subs_nbx[GMT]->Fill(float(nSUBS[GMT]),float(BxInEvent));
00395
00397
00398 if(BxInEvent!=0) continue;
00399
00400
00401 int Bx = RRItr->getBxNr();
00402 int Ev = RRItr->getEvNr();
00403
00404 bx_number->Fill(double(Bx));
00405
00406 for(int i=0; i<4; i++) {
00407 for( INPItr = INPCands[i].begin(); INPItr != INPCands[i].end(); ++INPItr ) {
00408 if(INPItr->empty()) continue;
00409 subs_eta[i]->Fill(INPItr->etaValue());
00410 subs_phi[i]->Fill(piconv_*INPItr->phiValue());
00411 subs_pt[i]->Fill(INPItr->ptValue());
00412 subs_qty[i]->Fill(INPItr->quality());
00413 subs_etaphi[i]->Fill(INPItr->etaValue(),piconv_*INPItr->phiValue());
00414 subs_etaqty[i]->Fill(INPItr->etaValue(),INPItr->quality());
00415 int word = INPItr->getDataWord();
00416 for( int j=0; j<32; j++ ) {
00417 if( word&(1<<j) ) subs_bits[i]->Fill(float(j));
00418 }
00419 }
00420 subs_candlumi[i]->Fill(float(e.luminosityBlock()),float(nSUBS[i]));
00421 }
00422
00423 for( GMTItr = GMTCands.begin(); GMTItr != GMTCands.end(); ++GMTItr ) {
00424 if(GMTItr->empty()) continue;
00425 subs_eta[GMT]->Fill(GMTItr->etaValue());
00426 subs_phi[GMT]->Fill(piconv_*GMTItr->phiValue());
00427 subs_pt[GMT]->Fill(GMTItr->ptValue());
00428 subs_qty[GMT]->Fill(GMTItr->quality());
00429 subs_etaphi[GMT]->Fill(GMTItr->etaValue(),piconv_*GMTItr->phiValue());
00430 subs_etaqty[GMT]->Fill(GMTItr->etaValue(),GMTItr->quality());
00431 int word = GMTItr->getDataWord();
00432 for( int j=0; j<32; j++ ) {
00433 if( word&(1<<j) ) subs_bits[GMT]->Fill(float(j));
00434 }
00435
00436 if(GMTItr->isMatchedCand()) {
00437 if(GMTItr->quality()>3) {
00438 eta_dtcsc_and_rpc->Fill(GMTItr->etaValue());
00439 phi_dtcsc_and_rpc->Fill(piconv_*GMTItr->phiValue());
00440 etaphi_dtcsc_and_rpc->Fill(GMTItr->etaValue(),piconv_*GMTItr->phiValue());
00441 }
00442 } else if(GMTItr->isRPC()) {
00443 if(GMTItr->quality()>3) {
00444 eta_rpc_only->Fill(GMTItr->etaValue());
00445 phi_rpc_only->Fill(piconv_*GMTItr->phiValue());
00446 etaphi_rpc_only->Fill(GMTItr->etaValue(),piconv_*GMTItr->phiValue());
00447 }
00448 } else {
00449 if(GMTItr->quality()>3) {
00450 eta_dtcsc_only->Fill(GMTItr->etaValue());
00451 phi_dtcsc_only->Fill(piconv_*GMTItr->phiValue());
00452 etaphi_dtcsc_only->Fill(GMTItr->etaValue(),piconv_*GMTItr->phiValue());
00453 }
00454
00455 if(GMTItr != GMTCands.end()){
00456 for( GMTItr2 = GMTCands.begin(); GMTItr2 != GMTCands.end(); ++GMTItr2 ) {
00457 if(GMTItr2==GMTItr) continue;
00458 if(GMTItr2->empty()) continue;
00459 if(GMTItr2->isRPC()) {
00460 if(GMTItr->isFwd()) {
00461 dist_eta_csc_rpc->Fill( GMTItr->etaValue() - GMTItr2->etaValue() );
00462 dist_phi_csc_rpc->Fill( piconv_*(GMTItr->phiValue() - GMTItr2->phiValue()) );
00463 } else {
00464 dist_eta_dt_rpc->Fill( GMTItr->etaValue() - GMTItr2->etaValue() );
00465 dist_phi_dt_rpc->Fill( piconv_*(GMTItr->phiValue() - GMTItr2->phiValue()) );
00466 }
00467 } else {
00468 if(!(GMTItr->isFwd()) && GMTItr2->isFwd()) {
00469 dist_eta_dt_csc->Fill( GMTItr->etaValue() - GMTItr2->etaValue() );
00470 dist_phi_dt_csc->Fill( piconv_*(GMTItr->phiValue() - GMTItr2->phiValue()) );
00471 } else if(GMTItr->isFwd() && !(GMTItr2->isFwd())){
00472 dist_eta_dt_csc->Fill( GMTItr2->etaValue() - GMTItr->etaValue() );
00473 dist_phi_dt_csc->Fill( piconv_*(GMTItr2->phiValue() - GMTItr->phiValue()) );
00474 }
00475 }
00476 }
00477 }
00478
00479 }
00480
00481 }
00482 subs_candlumi[GMT]->Fill(float(e.luminosityBlock()),float(nSUBS[GMT]));
00483
00484 n_rpcb_vs_dttf ->Fill(float(nSUBS[DTTF]),float(nSUBS[RPCb]));
00485 n_rpcf_vs_csctf->Fill(float(nSUBS[CSCTF]),float(nSUBS[RPCf]));
00486 n_csctf_vs_dttf->Fill(float(nSUBS[DTTF]),float(nSUBS[CSCTF]));
00487
00488 regional_triggers->Fill(-1.);
00489 for(int i=0; i<4; i++) {
00490 if(nSUBS[i]) regional_triggers->Fill(float(i));
00491 }
00492
00493
00494 if( (Ev - evnum_old_) == 1 && bxnum_old_ > -1 ) {
00495 int dBx = Bx - bxnum_old_;
00496 for(int id = 0; id<4; id++) {
00497 if( trsrc_old_&(1<<id) ) {
00498 for(int i=0; i<4; i++) {
00499 if(nSUBS[i]) subs_dbx[i]->Fill(float(dBx),float(id));
00500 }
00501 }
00502 }
00503
00504 }
00505
00506
00507 evnum_old_ = Ev;
00508 bxnum_old_ = Bx;
00509 trsrc_old_ = 0;
00510 for(int i=0; i<4; i++) {
00511 if(nSUBS[i]) trsrc_old_ |= (1<<i);
00512 }
00513
00514 }
00515
00516 }
00517
00518