Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-11 04:32:28

0001 /*
0002  * \file DTLocalTriggerBaseTask.cc
0003  *
0004  * \author C. Battilana - CIEMAT
0005  *
0006 */
0007 
0008 #include "DQM/DTMonitorModule/src/DTLocalTriggerBaseTask.h"
0009 
0010 // Framework
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 
0013 // DT DQM
0014 #include "DQM/DTMonitorModule/interface/DTTimeEvolutionHisto.h"
0015 
0016 // DT trigger
0017 #include "DQM/DTMonitorModule/interface/DTTrigGeomUtils.h"
0018 
0019 // Geometry
0020 #include "DataFormats/GeometryVector/interface/Pi.h"
0021 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0022 #include "Geometry/DTGeometry/interface/DTLayer.h"
0023 #include "Geometry/DTGeometry/interface/DTTopology.h"
0024 
0025 //Root
0026 #include "TH1.h"
0027 #include "TAxis.h"
0028 
0029 #include <sstream>
0030 #include <iostream>
0031 #include <fstream>
0032 
0033 using namespace edm;
0034 using namespace std;
0035 
0036 class DTTPGCompareUnit {
0037 public:
0038   DTTPGCompareUnit() { m_qual = -1; }
0039   ~DTTPGCompareUnit() {}
0040 
0041   void setTM(int qual, int bx) {
0042     m_qual = qual;
0043     m_BX = bx;
0044   }
0045 
0046   int qualTM() const { return m_qual; }
0047 
0048 private:
0049   int m_qual;
0050   int m_BX;
0051 };
0052 
0053 DTLocalTriggerBaseTask::DTLocalTriggerBaseTask(const edm::ParameterSet& ps)
0054     : m_nEvents(0), m_nLumis(0), m_trigGeomUtils(nullptr), muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
0055   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerBaseTask") << "[DTLocalTriggerBaseTask]: Constructor" << endl;
0056 
0057   m_tpMode = ps.getUntrackedParameter<bool>("testPulseMode");
0058   m_detailedAnalysis = ps.getUntrackedParameter<bool>("detailedAnalysis");
0059 
0060   m_targetBXTM = ps.getUntrackedParameter<int>("targetBXTM");
0061   m_bestAccRange = ps.getUntrackedParameter<int>("bestTrigAccRange");
0062 
0063   m_processTM = ps.getUntrackedParameter<bool>("processTM");
0064   m_processAB7 = ps.getUntrackedParameter<bool>("processAB7");
0065 
0066   m_tm_phiIn_Token = consumes<L1MuDTChambPhContainer>(ps.getUntrackedParameter<InputTag>("inputTagTMphIn"));
0067   m_tm_phiOut_Token = consumes<L1MuDTChambPhContainer>(ps.getUntrackedParameter<InputTag>("inputTagTMphOut"));
0068   m_tm_theta_Token = consumes<L1MuDTChambThContainer>(ps.getUntrackedParameter<InputTag>("inputTagTMth"));
0069   m_ab7_phi_Token = consumes<L1Phase2MuDTPhContainer>(ps.getUntrackedParameter<InputTag>("inputTagAB7"));
0070 
0071   if (m_processTM)
0072     m_types.push_back("TM");
0073 
0074   if (m_processAB7)
0075     m_types.push_back("AB7");
0076 
0077   if (m_tpMode) {
0078     topFolder("TM") = "DT/11-LocalTriggerTP-TM/";
0079     topFolder("AB7") = "DT/12-LocalTriggerTP-SliceTest/";
0080   } else {
0081     topFolder("TM") = "DT/03-LocalTrigger-TM/";
0082     topFolder("AB7") = "DT/04-LocalTrigger-SliceTest/";
0083   }
0084 
0085   m_params = ps;
0086 }
0087 
0088 DTLocalTriggerBaseTask::~DTLocalTriggerBaseTask() {
0089   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerBaseTask")
0090       << "[DTLocalTriggerBaseTask]: analyzed " << m_nEvents << " events" << endl;
0091   if (m_trigGeomUtils)
0092     delete m_trigGeomUtils;
0093 }
0094 
0095 void DTLocalTriggerBaseTask::bookHistograms(DQMStore::IBooker& ibooker,
0096                                             edm::Run const& iRun,
0097                                             edm::EventSetup const& context) {
0098   ibooker.setCurrentFolder("DT/EventInfo/Counters");
0099   m_nEventMonitor = ibooker.bookFloat("nProcessedEventsTrigger");
0100   for (int wh = -2; wh < 3; ++wh) {
0101     for (int stat = 1; stat < 5; ++stat) {
0102       for (int sect = 1; sect < 13; ++sect) {
0103         bookHistos(ibooker, DTChamberId(wh, stat, sect));
0104       }
0105     }
0106   }
0107 }
0108 
0109 void DTLocalTriggerBaseTask::beginLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
0110   m_nEventsInLS = 0;
0111   m_nLumis++;
0112   int resetCycle = m_params.getUntrackedParameter<int>("ResetCycle");
0113 
0114   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerBaseTask")
0115       << "[DTLocalTriggerBaseTask]: Begin of LS transition" << endl;
0116 
0117   if (m_nLumis % resetCycle == 0)
0118     for (auto& histosInChamb : m_chamberHistos)
0119       for (auto& histo : histosInChamb.second)
0120         histo.second->Reset();
0121 }
0122 
0123 void DTLocalTriggerBaseTask::endLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& context) {
0124   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerBaseTask") << "[DTLocalTriggerBaseTask]: End of LS transition" << endl;
0125 
0126   for (auto& trendHisto : m_trendHistos)
0127     trendHisto.second->updateTimeSlot(lumiSeg.luminosityBlock(), m_nEventsInLS);
0128 }
0129 
0130 void DTLocalTriggerBaseTask::dqmBeginRun(const edm::Run& run, const edm::EventSetup& context) {
0131   LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerBaseTask") << "[DTLocalTriggerBaseTask]: BeginRun" << endl;
0132 
0133   geom = &context.getData(muonGeomToken_);
0134   m_trigGeomUtils = new DTTrigGeomUtils(geom);
0135 }
0136 
0137 void DTLocalTriggerBaseTask::analyze(const edm::Event& e, const edm::EventSetup& c) {
0138   m_nEvents++;
0139   m_nEventsInLS++;
0140   m_nEventMonitor->Fill(m_nEvents);
0141 
0142   m_compMapIn.clear();
0143   m_compMapOut.clear();
0144 
0145   Handle<L1MuDTChambPhContainer> phiInTrigsTM;
0146   Handle<L1MuDTChambPhContainer> phiOutTrigsTM;
0147   Handle<L1MuDTChambThContainer> thetaTrigsTM;
0148   Handle<L1Phase2MuDTPhContainer> phiTrigsAB7;
0149 
0150   if (m_processTM) {
0151     e.getByToken(m_tm_phiIn_Token, phiInTrigsTM);
0152     e.getByToken(m_tm_phiOut_Token, phiOutTrigsTM);
0153     e.getByToken(m_tm_theta_Token, thetaTrigsTM);
0154 
0155     if (phiInTrigsTM.isValid() && phiOutTrigsTM.isValid() && thetaTrigsTM.isValid()) {
0156       runTMAnalysis(phiInTrigsTM->getContainer(), phiOutTrigsTM->getContainer(), thetaTrigsTM->getContainer());
0157     } else {
0158       LogVerbatim("DTDQM|DTMonitorModule|DTLocalTriggerBaseTask")
0159           << "[DTLocalTriggerBaseTask]: one or more TM tokens not valid!" << endl;
0160       return;
0161     }
0162   }
0163 
0164   if (m_processAB7) {
0165     e.getByToken(m_ab7_phi_Token, phiTrigsAB7);
0166 
0167     if (phiTrigsAB7.isValid()) {
0168       runAB7Analysis(phiTrigsAB7->getContainer());
0169     } else {
0170       LogVerbatim("DTDQM|DTMonitorModule|DTLocalTriggerBaseTask")
0171           << "[DTLocalTriggerBaseTask]: AB7 token not valid!" << endl;
0172     }
0173   }
0174 }
0175 
0176 void DTLocalTriggerBaseTask::bookHistos(DQMStore::IBooker& ibooker, const DTChamberId& dtCh) {
0177   uint32_t rawId = dtCh.rawId();
0178 
0179   stringstream wheel;
0180   wheel << dtCh.wheel();
0181   stringstream station;
0182   station << dtCh.station();
0183   stringstream sector;
0184   sector << dtCh.sector();
0185 
0186   map<string, int> minBX;
0187   map<string, int> maxBX;
0188 
0189   minBX["TM"] = m_params.getUntrackedParameter<int>("minBXTM");
0190   maxBX["TM"] = m_params.getUntrackedParameter<int>("maxBXTM");
0191   minBX["AB7"] = m_params.getUntrackedParameter<int>("minBXAB7");
0192   maxBX["AB7"] = m_params.getUntrackedParameter<int>("maxBXAB7");
0193 
0194   string chTag = "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
0195   string labelInOut = "";
0196 
0197   for (const auto& type : m_types) {
0198     LogTrace("DTDQM|DTMonitorModule|DTLocalTriggerBaseTask")
0199         << "[DTLocalTriggerBaseTask]: booking histos for " << topFolder(type) << "Wheel" << wheel.str() << "/Sector"
0200         << sector.str() << "/Station" << station.str() << endl;
0201 
0202     if (type == "AB7" && (dtCh.wheel() != 2 || dtCh.sector() != 12))
0203       continue;
0204 
0205     vector<string> plotLabels;
0206     vector<string> folderLabels;
0207 
0208     if (type == "TM") {
0209       plotLabels.push_back("_In");
0210       plotLabels.push_back("_Out");
0211       folderLabels.push_back("/LocalTriggerPhiIn");
0212       folderLabels.push_back("/LocalTriggerPhiOut");
0213     }
0214     if (type == "AB7") {
0215       plotLabels.push_back("");
0216       folderLabels.push_back("/LocalTriggerPhi");
0217     }
0218 
0219     for (size_t iLabel = 0; iLabel < plotLabels.size(); ++iLabel) {
0220       // Book Phi View Related Plots
0221 
0222       auto plotLabel = plotLabels.at(iLabel);
0223       ibooker.setCurrentFolder(topFolder(type) + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" +
0224                                station.str() + folderLabels.at(iLabel));
0225 
0226       int nQualities = type == "AB7" ? 11 : 7;
0227 
0228       string histoTag = type + "_BXvsQual" + plotLabel;
0229       m_chamberHistos[rawId][histoTag] = ibooker.book2D(histoTag + chTag,
0230                                                         "BX vs trigger quality",
0231                                                         nQualities,
0232                                                         -0.5,
0233                                                         nQualities - 0.5,
0234                                                         (int)(maxBX[type] - minBX[type] + 1),
0235                                                         minBX[type] - .5,
0236                                                         maxBX[type] + .5);
0237       if (type == "AB7")
0238         setQLabelsPh2((m_chamberHistos[rawId])[histoTag], 1);
0239       else
0240         setQLabels((m_chamberHistos[rawId])[histoTag], 1);
0241 
0242       if (!m_tpMode && !(type == "AB7")) {
0243         histoTag = type + "_BestQual" + plotLabel;
0244         m_chamberHistos[rawId][histoTag] =
0245             ibooker.book1D(histoTag + chTag, "Trigger quality of best primitives", 7, -0.5, 6.5);
0246         setQLabels(m_chamberHistos[rawId][histoTag], 1);
0247 
0248         histoTag = type + "_Flag1stvsQual" + plotLabel;
0249         m_chamberHistos[dtCh.rawId()][histoTag] =
0250             ibooker.book2D(histoTag + chTag, "1st/2nd trig flag vs quality", 7, -0.5, 6.5, 2, -0.5, 1.5);
0251         setQLabels(m_chamberHistos[rawId][histoTag], 1);
0252 
0253         histoTag = type + "_FlagUpDownvsQual" + plotLabel;
0254         m_chamberHistos[dtCh.rawId()][histoTag] =
0255             ibooker.book2D(histoTag + chTag, "Up/Down trig flag vs quality", 7, -0.5, 6.5, 2, -0.5, 1.5);
0256         setQLabels(m_chamberHistos[rawId][histoTag], 1);
0257       }
0258 
0259       if (type == "TM") {
0260         float minPh, maxPh;
0261         int nBinsPh;
0262         m_trigGeomUtils->phiRange(dtCh, minPh, maxPh, nBinsPh);
0263 
0264         histoTag = type + "_QualvsPhirad" + plotLabel;
0265         m_chamberHistos[rawId][histoTag] =
0266             ibooker.book2D(histoTag + chTag, "Trigger quality vs local position", nBinsPh, minPh, maxPh, 7, -0.5, 6.5);
0267         setQLabels(m_chamberHistos[rawId][histoTag], 2);
0268 
0269         if (plotLabel == "_Out" && !m_tpMode) {
0270           histoTag = type + "_RPCBitvsQual" + plotLabel;
0271           m_chamberHistos[rawId][histoTag] =
0272               ibooker.book2D(histoTag + chTag, "RPC bit vs DT trigger quality", 9, -1.5, 7.5, 3, -0.5, 2.5);
0273           //setQLabels((m_chamberHistos[dtCh.rawId()])[histoTag], 2);
0274         }
0275 
0276         if (m_detailedAnalysis && !m_tpMode) {
0277           histoTag = type + "_QualvsPhibend" + plotLabel;
0278           m_chamberHistos[rawId][histoTag] =
0279               ibooker.book2D(histoTag + chTag, "Trigger quality vs local direction", 200, -40., 40., 7, -0.5, 6.5);
0280           setQLabels((m_chamberHistos[dtCh.rawId()])[histoTag], 2);
0281         }
0282       }
0283     }
0284 
0285     // Book Theta View Related Plots
0286     ibooker.setCurrentFolder(topFolder(type) + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" +
0287                              station.str() + "/LocalTriggerTheta");
0288 
0289     string histoTag = "";
0290     if (type == "TM" && dtCh.station() != 4) {
0291       histoTag = type + "_PositionvsBX";
0292       m_chamberHistos[rawId][histoTag] = ibooker.book2D(histoTag + chTag,
0293                                                         "Theta trigger position vs BX",
0294                                                         (int)(maxBX[type] - minBX[type] + 1),
0295                                                         minBX[type] - .5,
0296                                                         maxBX[type] + .5,
0297                                                         7,
0298                                                         -0.5,
0299                                                         6.5);
0300       histoTag = type + "_PositionvsQual";
0301       m_chamberHistos[rawId][histoTag] =
0302           ibooker.book2D(histoTag + chTag, "Theta trigger position vs quality", 3, 0.5, 3.5, 7, -0.5, 6.5);
0303       setQLabelsTheta(m_chamberHistos[rawId][histoTag], 1);
0304       histoTag = type + "_ThetaBXvsQual";
0305       m_chamberHistos[rawId][histoTag] = ibooker.book2D(histoTag + chTag,
0306                                                         "BX vs trigger quality",
0307                                                         3,
0308                                                         0.5,
0309                                                         3.5,
0310                                                         (int)(maxBX[type] - minBX[type] + 1),
0311                                                         minBX[type] - .5,
0312                                                         maxBX[type] + .5);
0313       setQLabelsTheta(m_chamberHistos[rawId][histoTag], 1);
0314     }
0315   }
0316 }
0317 
0318 void DTLocalTriggerBaseTask::runTMAnalysis(std::vector<L1MuDTChambPhDigi> const* phInTrigs,
0319                                            std::vector<L1MuDTChambPhDigi> const* phOutTrigs,
0320                                            std::vector<L1MuDTChambThDigi> const* thTrigs) {
0321   vector<L1MuDTChambPhDigi>::const_iterator iph = phInTrigs->begin();
0322   vector<L1MuDTChambPhDigi>::const_iterator iphe = phInTrigs->end();
0323 
0324   for (; iph != iphe; ++iph) {
0325     int wh = iph->whNum();
0326     int sec = iph->scNum() + 1;  // B(O)MTF->DT Convention
0327     int st = iph->stNum();
0328     int qual = iph->code();
0329     int is1st = iph->Ts2Tag() ? 1 : 0;
0330     int bx = iph->bxNum() - is1st;
0331     int updown = iph->UpDownTag();
0332 
0333     if (qual < 0 || qual > 6)
0334       continue;  // Check that quality is in a valid range
0335 
0336     DTChamberId dtChId(wh, st, sec);
0337     uint32_t rawId = dtChId.rawId();
0338 
0339     float pos = m_trigGeomUtils->trigPos(&(*iph));
0340     float dir = m_trigGeomUtils->trigDir(&(*iph));
0341 
0342     if (abs(bx - m_targetBXTM) <= m_bestAccRange && m_compMapIn[rawId].qualTM() <= qual)
0343       m_compMapIn[rawId].setTM(qual, bx);
0344 
0345     map<string, MonitorElement*>& innerME = m_chamberHistos[rawId];
0346     if (m_tpMode) {
0347       innerME["TM_BXvsQual_In"]->Fill(qual, bx);       // SM BX vs Qual Phi view (1st tracks)
0348       innerME["TM_QualvsPhirad_In"]->Fill(pos, qual);  // SM Qual vs radial angle Phi view
0349     } else {
0350       innerME["TM_BXvsQual_In"]->Fill(qual, bx);              // SM BX vs Qual Phi view (1st tracks)
0351       innerME["TM_Flag1stvsQual_In"]->Fill(qual, is1st);      // SM Qual 1st/2nd track flag Phi view
0352       innerME["TM_FlagUpDownvsQual_In"]->Fill(qual, updown);  // SM Qual Up/Down track flag Phi view
0353       if (!is1st)
0354         innerME["TM_QualvsPhirad_In"]->Fill(pos, qual);  // SM Qual vs radial angle Phi view ONLY for 1st tracks
0355       if (m_detailedAnalysis) {
0356         innerME["TM_QualvsPhibend_In"]->Fill(dir, qual);  // SM Qual vs bending Phi view
0357       }
0358     }
0359   }
0360 
0361   iph = phOutTrigs->begin();
0362   iphe = phOutTrigs->end();
0363 
0364   for (; iph != iphe; ++iph) {
0365     int wh = iph->whNum();
0366     int sec = iph->scNum() + 1;  // B(O)MTF->DT Convention
0367     int st = iph->stNum();
0368     int qual = iph->code();
0369     int is1st = iph->Ts2Tag() ? 1 : 0;
0370     int rpcBit = iph->RpcBit();
0371     int bx = iph->bxNum() - is1st;
0372     int updown = iph->UpDownTag();
0373     if (qual < 0 || qual > 6)
0374       continue;  // Check that quality is in a valid range
0375 
0376     DTChamberId dtChId(wh, st, sec);
0377     uint32_t rawId = dtChId.rawId();
0378 
0379     float pos = m_trigGeomUtils->trigPos(&(*iph));
0380     float dir = m_trigGeomUtils->trigDir(&(*iph));
0381 
0382     if (abs(bx - m_targetBXTM) <= m_bestAccRange && m_compMapOut[rawId].qualTM() <= qual)
0383       m_compMapOut[rawId].setTM(qual, bx);
0384 
0385     map<string, MonitorElement*>& innerME = m_chamberHistos[rawId];
0386     if (m_tpMode) {
0387       innerME["TM_BXvsQual_Out"]->Fill(qual, bx);       // SM BX vs Qual Phi view (1st tracks)
0388       innerME["TM_QualvsPhirad_Out"]->Fill(pos, qual);  // SM Qual vs radial angle Phi view
0389     } else {
0390       innerME["TM_BXvsQual_Out"]->Fill(qual, bx);              // SM BX vs Qual Phi view (1st tracks)
0391       innerME["TM_RPCBitvsQual_Out"]->Fill(qual, rpcBit);      // SM RPC bitvs Qual Phi view
0392       innerME["TM_Flag1stvsQual_Out"]->Fill(qual, is1st);      // SM Qual 1st/2nd track flag Phi view
0393       innerME["TM_FlagUpDownvsQual_Out"]->Fill(qual, updown);  // SM Qual Up/Down track flag Phi view
0394 
0395       if (!is1st)
0396         innerME["TM_QualvsPhirad_Out"]->Fill(pos, qual);  // SM Qual vs radial angle Phi view ONLY for 1st tracks
0397       if (m_detailedAnalysis) {
0398         innerME["TM_QualvsPhibend_Out"]->Fill(dir, qual);  // SM Qual vs bending Phi view
0399       }
0400     }
0401   }
0402 
0403   vector<L1MuDTChambThDigi>::const_iterator ith = thTrigs->begin();
0404   vector<L1MuDTChambThDigi>::const_iterator ithe = thTrigs->end();
0405 
0406   for (; ith != ithe; ++ith) {
0407     int wh = ith->whNum();
0408     int sec = ith->scNum() + 1;  // B(O)MTF -> DT Convention
0409     int st = ith->stNum();
0410     int bx = ith->bxNum();
0411 
0412     int thcode[7];
0413 
0414     for (int pos = 0; pos < 7; pos++) {
0415       thcode[pos] = ith->code(pos);
0416       if (ith->position(pos) == 0 && ith->quality(pos) == 1)
0417         thcode[pos] = 3;
0418     }
0419 
0420     DTChamberId dtChId(wh, st, sec);
0421     uint32_t rawId = dtChId.rawId();
0422 
0423     map<string, MonitorElement*>& innerME = m_chamberHistos[rawId];
0424 
0425     for (int pos = 0; pos < 7; pos++)
0426       if (thcode[pos] > 0) {                                   //Fired
0427         innerME["TM_PositionvsBX"]->Fill(bx, pos);             // SM BX vs Position Theta view
0428         innerME["TM_PositionvsQual"]->Fill(thcode[pos], pos);  //code = pos + qual; so 0, 1, 2 for 0, L, H resp.
0429         innerME["TM_ThetaBXvsQual"]->Fill(thcode[pos], bx);    //code = pos + qual; so 0, 1, 2 for 0, L, H resp.
0430       }
0431   }
0432   // Fill Quality plots with best TM triggers (phi view In)
0433   if (!m_tpMode) {
0434     for (auto& comp : m_compMapIn) {
0435       int bestQual = comp.second.qualTM();
0436       if (bestQual > -1)
0437         m_chamberHistos[comp.first]["TM_BestQual_In"]->Fill(bestQual);  // SM Best Qual Trigger Phi view
0438     }
0439   }
0440 
0441   // Fill Quality plots with best TM triggers (phi view Out)
0442   if (!m_tpMode) {
0443     for (auto& comp : m_compMapOut) {
0444       int bestQual = comp.second.qualTM();
0445       if (bestQual > -1)
0446         m_chamberHistos[comp.first]["TM_BestQual_Out"]->Fill(bestQual);  // SM Best Qual Trigger Phi view
0447     }
0448   }
0449 }
0450 
0451 void DTLocalTriggerBaseTask::runAB7Analysis(std::vector<L1Phase2MuDTPhDigi> const* phTrigs) {
0452   vector<L1Phase2MuDTPhDigi>::const_iterator iph = phTrigs->begin();
0453   vector<L1Phase2MuDTPhDigi>::const_iterator iphe = phTrigs->end();
0454   for (; iph != iphe; ++iph) {
0455     int wh = iph->whNum();
0456     int sec = iph->scNum() + 1;  // B(O)MTF->DT Convention
0457     int st = iph->stNum();
0458     int qual = iph->quality();
0459     int bx = iph->bxNum();
0460 
0461     DTChamberId dtChId(wh, st, sec);
0462     uint32_t rawId = dtChId.rawId();
0463 
0464     map<string, MonitorElement*>& innerME = m_chamberHistos[rawId];
0465     innerME["AB7_BXvsQual"]->Fill(qual, bx);
0466   }
0467 }
0468 
0469 void DTLocalTriggerBaseTask::setQLabels(MonitorElement* me, short int iaxis) {
0470   TH1* histo = me->getTH1();
0471   if (!histo)
0472     return;
0473 
0474   TAxis* axis = nullptr;
0475   if (iaxis == 1) {
0476     axis = histo->GetXaxis();
0477   } else if (iaxis == 2) {
0478     axis = histo->GetYaxis();
0479   }
0480   if (!axis)
0481     return;
0482 
0483   string labels[7] = {"LI", "LO", "HI", "HO", "LL", "HL", "HH"};
0484   int istart = axis->GetXmin() < -1 ? 2 : 1;
0485   for (int i = 0; i < 7; i++) {
0486     axis->SetBinLabel(i + istart, labels[i].c_str());
0487   }
0488 }
0489 
0490 void DTLocalTriggerBaseTask::setQLabelsTheta(MonitorElement* me, short int iaxis) {
0491   TH1* histo = me->getTH1();
0492   if (!histo)
0493     return;
0494 
0495   TAxis* axis = nullptr;
0496   if (iaxis == 1) {
0497     axis = histo->GetXaxis();
0498   } else if (iaxis == 2) {
0499     axis = histo->GetYaxis();
0500   }
0501   if (!axis)
0502     return;
0503 
0504   string labels[3] = {"L", "H", "err"};
0505   int istart = axis->GetXmin() < -1 ? 2 : 1;
0506   for (int i = 0; i < 3; i++) {
0507     axis->SetBinLabel(i + istart, labels[i].c_str());
0508   }
0509 }
0510 
0511 void DTLocalTriggerBaseTask::setQLabelsPh2(MonitorElement* me, short int iaxis) {
0512   TH1* histo = me->getTH1();
0513   if (!histo)
0514     return;
0515 
0516   TAxis* axis = nullptr;
0517   if (iaxis == 1) {
0518     axis = histo->GetXaxis();
0519   } else if (iaxis == 2) {
0520     axis = histo->GetYaxis();
0521   }
0522   if (!axis)
0523     return;
0524 
0525   string labels[11] = {"", "L only", "L multiple", "H only", "H multiple", "3+2", "LL", "4+2", "HL", "HH", ""};
0526   int istart = axis->GetXmin() < -1 ? 2 : 1;
0527   for (int i = 0; i < 11; i++) {
0528     axis->SetBinLabel(i + istart, labels[i].c_str());
0529   }
0530 }
0531 
0532 // Local Variables:
0533 // show-trailing-whitespace: t
0534 // truncate-lines: t
0535 // End: