Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 /*
0002  * \file DTDigiTask.cc
0003  *
0004  * \author M. Zanetti - INFN Padova
0005  *
0006  */
0007 
0008 #include "DQM/DTMonitorModule/interface/DTDigiTask.h"
0009 
0010 // Framework
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 
0013 // Digis
0014 #include "DataFormats/MuonDetId/interface/DTLayerId.h"
0015 #include "DataFormats/MuonDetId/interface/DTChamberId.h"
0016 
0017 // Geometry
0018 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0019 #include "Geometry/DTGeometry/interface/DTLayer.h"
0020 #include "Geometry/DTGeometry/interface/DTTopology.h"
0021 
0022 // T0s
0023 #include "CondFormats/DTObjects/interface/DTT0.h"
0024 #include "CondFormats/DTObjects/interface/DTTtrig.h"
0025 
0026 #include <sstream>
0027 #include <cmath>
0028 
0029 using namespace edm;
0030 using namespace std;
0031 
0032 // Contructor
0033 DTDigiTask::DTDigiTask(const edm::ParameterSet& ps)
0034     : muonGeomToken_(esConsumes<edm::Transition::BeginRun>()),
0035       readOutMapToken_(esConsumes<edm::Transition::BeginRun>()),
0036       TtrigToken_(esConsumes<edm::Transition::BeginRun>()),
0037       T0Token_(esConsumes<edm::Transition::BeginRun>()),
0038       statusMapToken_(esConsumes()) {
0039   // switch for the verbosity
0040   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "[DTDigiTask]: Constructor" << endl;
0041 
0042   // The label to retrieve the digis
0043   dtDigiToken_ = consumes<DTDigiCollection>(ps.getUntrackedParameter<InputTag>("dtDigiLabel"));
0044   // Read the configuration parameters
0045   maxTDCHits = ps.getUntrackedParameter<int>("maxTDCHitsPerChamber", 30000);
0046   // Set to true to read the ttrig from DB (useful to determine in-time and out-of-time hits)
0047   readTTrigDB = ps.getUntrackedParameter<bool>("readDB", false);
0048   // Set to true to subtract t0 from test pulses
0049   subtractT0 = ps.getParameter<bool>("performPerWireT0Calibration");
0050   // Tmax value (TDC counts)
0051   defaultTmax = ps.getParameter<int>("defaultTmax");
0052   // Switch from static to dinamic histo booking
0053   doStaticBooking = ps.getUntrackedParameter<bool>("staticBooking", true);
0054 
0055   // Switch for local/global runs
0056   isLocalRun = ps.getUntrackedParameter<bool>("localrun", true);
0057   if (!isLocalRun) {
0058     ltcDigiCollectionToken_ =
0059         consumes<LTCDigiCollection>(ps.getUntrackedParameter<edm::InputTag>("ltcDigiCollectionTag"));
0060   }
0061 
0062   // Setting for the reset of the ME after n (= ResetCycle) luminosity sections
0063   resetCycle = ps.getUntrackedParameter<int>("ResetCycle", 3);
0064   // Check the DB of noisy channels
0065   checkNoisyChannels = ps.getUntrackedParameter<bool>("checkNoisyChannels", false);
0066   // Default TTrig to be used when not reading the TTrig DB
0067   defaultTTrig = ps.getParameter<int>("defaultTtrig");
0068   inTimeHitsLowerBound = ps.getParameter<int>("inTimeHitsLowerBound");
0069   inTimeHitsUpperBound = ps.getParameter<int>("inTimeHitsUpperBound");
0070   timeBoxGranularity = ps.getUntrackedParameter<int>("timeBoxGranularity", 4);
0071   maxTTMounts = ps.getUntrackedParameter<int>("maxTTMounts", 6400);
0072 
0073   doAllHitsOccupancies = ps.getUntrackedParameter<bool>("doAllHitsOccupancies", true);
0074   doNoiseOccupancies = ps.getUntrackedParameter<bool>("doNoiseOccupancies", false);
0075   doInTimeOccupancies = ps.getUntrackedParameter<bool>("doInTimeOccupancies", false);
0076 
0077   // switch on the mode for running on test pulses (different top folder)
0078   tpMode = ps.getUntrackedParameter<bool>("testPulseMode", false);
0079   // switch on/off the filtering of synchronous noise events (cutting on the # of digis)
0080   // time-boxes and occupancy plots are not filled and summary plots are created to report the problem
0081   filterSyncNoise = ps.getUntrackedParameter<bool>("filterSyncNoise", false);
0082   // look for synch noisy events, produce histograms but do not filter them
0083   lookForSyncNoise = ps.getUntrackedParameter<bool>("lookForSyncNoise", false);
0084 
0085   // switch on the mode for running on slice test (different top folder and other customizations)
0086   sliceTestMode = ps.getUntrackedParameter<bool>("sliceTestMode", false);
0087   // time pedestal used to set the minimum in the time box plot
0088   tdcPedestal = ps.getUntrackedParameter<int>("tdcPedestal", 0);
0089 
0090   // switch on production of time-boxes with layer granularity
0091   doLayerTimeBoxes = ps.getUntrackedParameter<bool>("doLayerTimeBoxes", false);
0092 
0093   syncNumTot = 0;
0094   syncNum = 0;
0095 }
0096 
0097 // destructor
0098 DTDigiTask::~DTDigiTask() {
0099   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "DTDigiTask: analyzed " << nevents << " events" << endl;
0100 }
0101 
0102 void DTDigiTask::dqmBeginRun(const edm::Run& run, const edm::EventSetup& context) {
0103   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "[DTDigiTask]: begin run" << endl;
0104   nevents = 0;
0105 
0106   // Get the geometry
0107   muonGeom = &context.getData(muonGeomToken_);
0108 
0109   // map of the channels
0110   mapping = &context.getData(readOutMapToken_);
0111 
0112   // tTrig
0113   if (readTTrigDB)
0114     tTrigMap = &context.getData(TtrigToken_);
0115   // t0s
0116   if (subtractT0)
0117     t0Map = &context.getData(T0Token_);
0118   // FIXME: tMax (not yet from the DB)
0119   tMax = defaultTmax;
0120 
0121   // ----------------------------------------------------------------------
0122 }
0123 
0124 void DTDigiTask::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& run, edm::EventSetup const& context) {
0125   if (doStaticBooking) {  // Static histo booking
0126     // book the event counter
0127     ibooker.setCurrentFolder("DT/EventInfo/Counters");
0128     nEventMonitor = ibooker.bookFloat(tpMode ? "nProcessedEventsDigiTP" : "nProcessedEventsDigi");
0129     ibooker.setCurrentFolder(topFolder());
0130     for (int wh = -2; wh <= 2; ++wh) {  // loop over wheels
0131       if (sliceTestMode && wh != 2)
0132         continue;
0133 
0134       if (doAllHitsOccupancies)
0135         bookHistos(ibooker, wh, string("Occupancies"), "OccupancyAllHits");
0136 
0137       if (doNoiseOccupancies)
0138         bookHistos(ibooker, wh, string("Occupancies"), "OccupancyNoise");
0139       if (doInTimeOccupancies)
0140         bookHistos(ibooker, wh, string("Occupancies"), "OccupancyInTimeHits");
0141 
0142       if (lookForSyncNoise || filterSyncNoise) {
0143         bookHistos(ibooker, wh, string("SynchNoise"), "SyncNoiseEvents");
0144         bookHistos(ibooker, wh, string("SynchNoise"), "SyncNoiseChambs");
0145       }
0146 
0147       for (int st = 1; st <= 4; ++st) {           // loop over stations
0148         for (int sect = 1; sect <= 14; ++sect) {  // loop over sectors
0149           if ((sect == 13 || sect == 14) && st != 4)
0150             continue;
0151 
0152           if (sliceTestMode && (sect != 12 || wh != 2))
0153             continue;
0154 
0155           // Get the chamber ID
0156           const DTChamberId dtChId(wh, st, sect);
0157 
0158           // Occupancies
0159           if (doAllHitsOccupancies) {
0160             bookHistos(ibooker, dtChId, string("Occupancies"), "OccupancyAllHits_perCh");
0161             // set channel mapping
0162             channelsMap(dtChId, "OccupancyAllHits_perCh");
0163           }
0164           if (doNoiseOccupancies)
0165             bookHistos(ibooker, dtChId, string("Occupancies"), "OccupancyNoise_perCh");
0166           if (doInTimeOccupancies)
0167             bookHistos(ibooker, dtChId, string("Occupancies"), "OccupancyInTimeHits_perCh");
0168 
0169           for (int sl = 1; sl <= 3; ++sl) {  // Loop over SLs
0170             if (st == 4 && sl == 2)
0171               continue;
0172             const DTSuperLayerId dtSLId(wh, st, sect, sl);
0173             if (isLocalRun) {
0174               bookHistos(ibooker, dtSLId, string("TimeBoxes"), "TimeBox");
0175             } else {
0176               // TimeBoxes for different triggers
0177               bookHistos(ibooker, dtSLId, string("TimeBoxes"), "TimeBoxDTonly");
0178               bookHistos(ibooker, dtSLId, string("TimeBoxes"), "TimeBoxNoDT");
0179               bookHistos(ibooker, dtSLId, string("TimeBoxes"), "TimeBoxDTalso");
0180             }
0181           }
0182         }
0183       }
0184     }
0185   }
0186 }
0187 
0188 void DTDigiTask::beginLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& context) {
0189   LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "[DTDigiTask]: Begin of LS transition" << endl;
0190 
0191   // Reset the MonitorElements every n (= ResetCycle) Lumi Blocks
0192   int lumiBlock = lumiSeg.id().luminosityBlock();
0193   if (lumiBlock % resetCycle == 0) {
0194     LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask")
0195         << "[DTDigiTask]: Reset at the LS transition : " << lumiBlock << endl;
0196     // Loop over all ME
0197     map<string, map<uint32_t, MonitorElement*> >::const_iterator histosIt = digiHistos.begin();
0198     map<string, map<uint32_t, MonitorElement*> >::const_iterator histosEnd = digiHistos.end();
0199     for (; histosIt != histosEnd; ++histosIt) {
0200       map<uint32_t, MonitorElement*>::const_iterator histoIt = (*histosIt).second.begin();
0201       map<uint32_t, MonitorElement*>::const_iterator histoEnd = (*histosIt).second.end();
0202       for (; histoIt != histoEnd; ++histoIt) {
0203         (*histoIt).second->Reset();
0204       }
0205     }
0206 
0207     // re-set mapping for not real channels in the occupancyHits per chamber
0208     for (int wh = -2; wh <= 2; wh++) {
0209       for (int sect = 1; sect <= 14; sect++) {
0210         for (int st = 1; st <= 4; st++) {
0211           if (sliceTestMode && (sect != 12 || wh != 2)) {
0212             continue;
0213           }
0214           if ((sect == 13 || sect == 14) && st != 4) {
0215             continue;
0216           }
0217           const DTChamberId dtChId(wh, st, sect);
0218           channelsMap(dtChId, "OccupancyAllHits_perCh");
0219         }
0220       }
0221     }
0222 
0223     // loop over wheel summaries
0224     map<string, map<int, MonitorElement*> >::const_iterator whHistosIt = wheelHistos.begin();
0225     map<string, map<int, MonitorElement*> >::const_iterator whHistosEnd = wheelHistos.end();
0226     for (; whHistosIt != whHistosEnd; ++whHistosIt) {
0227       if ((*whHistosIt).first.find("Sync") == string::npos) {  // FIXME skips synch noise plots
0228         map<int, MonitorElement*>::const_iterator histoIt = (*whHistosIt).second.begin();
0229         map<int, MonitorElement*>::const_iterator histoEnd = (*whHistosIt).second.end();
0230         for (; histoIt != histoEnd; ++histoIt) {
0231           (*histoIt).second->Reset();
0232         }
0233       }
0234     }
0235   }
0236 }
0237 
0238 void DTDigiTask::bookHistos(DQMStore::IBooker& ibooker, const DTSuperLayerId& dtSL, string folder, string histoTag) {
0239   // set the folder
0240   stringstream wheel;
0241   wheel << dtSL.wheel();
0242   stringstream station;
0243   station << dtSL.station();
0244   stringstream sector;
0245   sector << dtSL.sector();
0246   stringstream superLayer;
0247   superLayer << dtSL.superlayer();
0248   ibooker.setCurrentFolder(topFolder() + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str());
0249 
0250   // Build the histo name
0251   string histoName =
0252       histoTag + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() + "_SL" + superLayer.str();
0253 
0254   LogTrace("DTDQM|DTMonitorModule|DTDigiTask")
0255       << "[DTDigiTask]: booking SL histo:" << histoName << " (tag: " << histoTag << ") folder: "
0256       << topFolder() + "Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" + sector.str() + "/" + folder
0257       << endl;
0258 
0259   // ttrig and rms are TDC counts
0260   if (readTTrigDB)
0261     tTrigMap->get(dtSL, tTrig, tTrigRMS, kFactor, DTTimeUnits::counts);
0262   else
0263     tTrig = defaultTTrig;
0264 
0265   if (folder == "TimeBoxes") {
0266     string histoTitle = histoName + " (TDC Counts)";
0267 
0268     if (!readTTrigDB) {
0269       (digiHistos[histoTag])[dtSL.rawId()] = ibooker.book1D(
0270           histoName, histoTitle, maxTTMounts / timeBoxGranularity, tdcPedestal, maxTTMounts + tdcPedestal);
0271       if (doLayerTimeBoxes) {  // Book TimeBoxes per layer
0272         for (int layer = 1; layer != 5; ++layer) {
0273           DTLayerId layerId(dtSL, layer);
0274           stringstream layerHistoName;
0275           layerHistoName << histoName << "_L" << layer;
0276           (digiHistos[histoTag])[layerId.rawId()] = ibooker.book1D(layerHistoName.str(),
0277                                                                    layerHistoName.str(),
0278                                                                    maxTTMounts / timeBoxGranularity,
0279                                                                    tdcPedestal,
0280                                                                    maxTTMounts + tdcPedestal);
0281         }
0282       }
0283     } else {
0284       (digiHistos[histoTag])[dtSL.rawId()] =
0285           ibooker.book1D(histoName, histoTitle, 3 * tMax / timeBoxGranularity, tTrig - tMax, tTrig + 2 * tMax);
0286       if (doLayerTimeBoxes) {
0287         // Book TimeBoxes per layer
0288         for (int layer = 1; layer != 5; ++layer) {
0289           DTLayerId layerId(dtSL, layer);
0290           stringstream layerHistoName;
0291           layerHistoName << histoName << "_L" << layer;
0292           (digiHistos[histoTag])[layerId.rawId()] = ibooker.book1D(
0293               layerHistoName.str(), layerHistoName.str(), 3 * tMax / timeBoxGranularity, tTrig - tMax, tTrig + 2 * tMax);
0294         }
0295       }
0296     }
0297   }
0298 
0299   if (folder == "CathodPhotoPeaks") {
0300     ibooker.setCurrentFolder(topFolder() + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" +
0301                              station.str() + "/" + folder);
0302     (digiHistos[histoTag])[dtSL.rawId()] = ibooker.book1D(histoName, histoName, 500, 0, 1000);
0303   }
0304 }
0305 
0306 void DTDigiTask::bookHistos(DQMStore::IBooker& ibooker, const DTChamberId& dtCh, string folder, string histoTag) {
0307   // set the current folder
0308   stringstream wheel;
0309   wheel << dtCh.wheel();
0310   stringstream station;
0311   station << dtCh.station();
0312   stringstream sector;
0313   sector << dtCh.sector();
0314   ibooker.setCurrentFolder(topFolder() + "Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str());
0315 
0316   // build the histo name
0317   string histoName = histoTag + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
0318 
0319   LogTrace("DTDQM|DTMonitorModule|DTDigiTask")
0320       << "[DTDigiTask]: booking chamber histo:"
0321       << " (tag: " << histoTag
0322       << ") folder: " << topFolder() + "Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" + sector.str()
0323       << endl;
0324 
0325   if (folder == "Occupancies") {
0326     const DTChamber* dtchamber = muonGeom->chamber(dtCh);
0327     const std::vector<const DTSuperLayer*>& dtSupLylist = dtchamber->superLayers();
0328     std::vector<const DTSuperLayer*>::const_iterator suly = dtSupLylist.begin();
0329     std::vector<const DTSuperLayer*>::const_iterator sulyend = dtSupLylist.end();
0330 
0331     int nWires = 0;
0332     int firstWire = 0;
0333     int nWires_max = 0;
0334 
0335     while (suly != sulyend) {
0336       const std::vector<const DTLayer*> dtLyList = (*suly)->layers();
0337       std::vector<const DTLayer*>::const_iterator ly = dtLyList.begin();
0338       std::vector<const DTLayer*>::const_iterator lyend = dtLyList.end();
0339       stringstream superLayer;
0340       superLayer << (*suly)->id().superlayer();
0341 
0342       while (ly != lyend) {
0343         nWires = muonGeom->layer((*ly)->id())->specificTopology().channels();
0344         firstWire = muonGeom->layer((*ly)->id())->specificTopology().firstChannel();
0345         stringstream layer;
0346         layer << (*ly)->id().layer();
0347         string histoName_layer = histoName + "_SL" + superLayer.str() + "_L" + layer.str();
0348         if (histoTag == "OccupancyAllHits_perL" || histoTag == "OccupancyNoise_perL" ||
0349             histoTag == "OccupancyInTimeHits_perL")
0350           (digiHistos[histoTag])[(*ly)->id().rawId()] =
0351               ibooker.book1D(histoName_layer, histoName_layer, nWires, firstWire, nWires + firstWire);
0352         ++ly;
0353         if ((nWires + firstWire) > nWires_max)
0354           nWires_max = (nWires + firstWire);
0355       }
0356       ++suly;
0357     }
0358 
0359     if (histoTag != "OccupancyAllHits_perL" && histoTag != "OccupancyNoise_perL" &&
0360         histoTag != "OccupancyInTimeHits_perL") {
0361       // Set the title to show the time interval used (only if unique == not from DB)
0362       string histoTitle = histoName;
0363       if (!readTTrigDB && histoTag == "OccupancyInTimeHits_perCh") {
0364         stringstream title;
0365         int inTimeHitsLowerBoundCorr = int(round(defaultTTrig)) - inTimeHitsLowerBound;
0366         int inTimeHitsUpperBoundCorr = int(round(defaultTTrig)) + defaultTmax + inTimeHitsUpperBound;
0367         title << "Occ. digis in time [" << inTimeHitsLowerBoundCorr << ", " << inTimeHitsUpperBoundCorr
0368               << "] (TDC counts)";
0369         histoTitle = title.str();
0370       }
0371       (digiHistos[histoTag])[dtCh.rawId()] =
0372           ibooker.book2D(histoName, histoTitle, nWires_max, 1, nWires_max + 1, 12, 0, 12);
0373 
0374       for (int i = 1; i <= 12; i++) {
0375         if (i < 5) {
0376           stringstream layer;
0377           string layer_name;
0378           layer << i;
0379           layer >> layer_name;
0380           string label = "SL1: L" + layer_name;
0381           (digiHistos[histoTag])[dtCh.rawId()]->setBinLabel(i, label, 2);
0382         } else if (i > 4 && i < 9) {
0383           stringstream layer;
0384           string layer_name;
0385           layer << (i - 4);
0386           layer >> layer_name;
0387           string label = "SL2: L" + layer_name;
0388           (digiHistos[histoTag])[dtCh.rawId()]->setBinLabel(i, label, 2);
0389         } else if (i > 8 && i < 13) {
0390           stringstream layer;
0391           string layer_name;
0392           layer << (i - 8);
0393           layer >> layer_name;
0394           string label = "SL3: L" + layer_name;
0395           (digiHistos[histoTag])[dtCh.rawId()]->setBinLabel(i, label, 2);
0396         }
0397       }
0398     }
0399   }
0400 }
0401 
0402 void DTDigiTask::bookHistos(DQMStore::IBooker& ibooker, const int wheelId, string folder, string histoTag) {
0403   // Set the current folder
0404   stringstream wheel;
0405   wheel << wheelId;
0406 
0407   // build the histo name
0408   string histoName = histoTag + "_W" + wheel.str();
0409 
0410   LogTrace("DTDQM|DTMonitorModule|DTDigiTask")
0411       << "[DTDigiTask]: booking wheel histo:" << histoName << " (tag: " << histoTag
0412       << ") folder: " << topFolder() + "Wheel" + wheel.str() + "/" << endl;
0413 
0414   if (folder == "Occupancies") {
0415     ibooker.setCurrentFolder(topFolder() + "Wheel" + wheel.str());
0416     string histoTitle = "# of digis per chamber WHEEL: " + wheel.str();
0417     (wheelHistos[histoTag])[wheelId] = ibooker.book2D(histoName, histoTitle, 12, 1, 13, 4, 1, 5);
0418     (wheelHistos[histoTag])[wheelId]->setBinLabel(1, "MB1", 2);
0419     (wheelHistos[histoTag])[wheelId]->setBinLabel(2, "MB2", 2);
0420     (wheelHistos[histoTag])[wheelId]->setBinLabel(3, "MB3", 2);
0421     (wheelHistos[histoTag])[wheelId]->setBinLabel(4, "MB4", 2);
0422     (wheelHistos[histoTag])[wheelId]->setAxisTitle("sector", 1);
0423   } else if (folder == "SynchNoise") {
0424     ibooker.setCurrentFolder("DT/05-Noise/SynchNoise");
0425     if (histoTag == "SyncNoiseEvents") {
0426       string histoTitle = "# of Syncronous-noise events WHEEL: " + wheel.str();
0427       (wheelHistos[histoTag])[wheelId] = ibooker.book2D(histoName, histoTitle, 12, 1, 13, 4, 1, 5);
0428       (wheelHistos[histoTag])[wheelId]->setBinLabel(1, "MB1", 2);
0429       (wheelHistos[histoTag])[wheelId]->setBinLabel(2, "MB2", 2);
0430       (wheelHistos[histoTag])[wheelId]->setBinLabel(3, "MB3", 2);
0431       (wheelHistos[histoTag])[wheelId]->setBinLabel(4, "MB4", 2);
0432       (wheelHistos[histoTag])[wheelId]->setAxisTitle("sector", 1);
0433     } else if (histoTag == "SyncNoiseChambs") {
0434       string histoTitle = "# of Synchornous-noise chamb per evt. WHEEL: " + wheel.str();
0435       (wheelHistos[histoTag])[wheelId] = ibooker.book1D(histoName, histoTitle, 50, 0.5, 50.5);
0436       (wheelHistos[histoTag])[wheelId]->setAxisTitle("# of noisy chambs.", 1);
0437       (wheelHistos[histoTag])[wheelId]->setAxisTitle("# of evts.", 2);
0438     }
0439   }
0440 }
0441 // does the real job
0442 void DTDigiTask::analyze(const edm::Event& event, const edm::EventSetup& c) {
0443   nevents++;
0444   nEventMonitor->Fill(nevents);
0445   if (nevents % 1000 == 0) {
0446     LogTrace("DTDQM|DTMonitorModule|DTDigiTask")
0447         << "[DTDigiTask] Analyze #Run: " << event.id().run() << " #Event: " << event.id().event() << endl;
0448   }
0449 
0450   // Get the ingredients from the event
0451 
0452   // Digi collection
0453   edm::Handle<DTDigiCollection> dtdigis;
0454   event.getByToken(dtDigiToken_, dtdigis);
0455 
0456   // LTC digis
0457   if (!isLocalRun)
0458     event.getByToken(ltcDigiCollectionToken_, ltcdigis);
0459 
0460   // Status map (for noisy channels)
0461   if (checkNoisyChannels) {
0462     // Get the map of noisy channels
0463     statusMap = &c.getData(statusMapToken_);
0464   }
0465 
0466   string histoTag;
0467 
0468   // Check if the digi container is empty
0469   if (dtdigis->begin() == dtdigis->end()) {
0470     LogTrace("DTDQM|DTMonitorModule|DTDigiTask") << "Event " << nevents << " empty." << endl;
0471   }
0472 
0473   if (lookForSyncNoise || filterSyncNoise) {  // dosync
0474     // Count the # of digis per chamber
0475     DTDigiCollection::DigiRangeIterator dtLayerId_It;
0476     for (dtLayerId_It = dtdigis->begin(); dtLayerId_It != dtdigis->end(); dtLayerId_It++) {
0477       DTChamberId chId = ((*dtLayerId_It).first).chamberId();
0478       if (hitMap.find(chId) == hitMap.end()) {  // new chamber
0479         hitMap[chId] = 0;
0480       }
0481       hitMap[chId] += (((*dtLayerId_It).second).second - ((*dtLayerId_It).second).first);
0482     }
0483 
0484     // check chamber with # of digis above threshold and flag them as noisy
0485     map<DTChamberId, int>::const_iterator hitMapIt = hitMap.begin();
0486     map<DTChamberId, int>::const_iterator hitMapEnd = hitMap.end();
0487 
0488     map<int, int> chMap;
0489 
0490     for (; hitMapIt != hitMapEnd; ++hitMapIt) {
0491       if ((hitMapIt->second) > maxTDCHits) {
0492         DTChamberId chId = hitMapIt->first;
0493         int wh = chId.wheel();
0494 
0495         LogTrace("DTDQM|DTMonitorModule|DTDigiTask")
0496             << "[DTDigiTask] Synch noise in chamber: " << chId << " with # digis: " << hitMapIt->second << endl;
0497 
0498         if (chMap.find(wh) == chMap.end()) {
0499           chMap[wh] = 0;
0500         }
0501         chMap[wh]++;
0502 
0503         syncNoisyChambers.insert(chId);
0504 
0505         wheelHistos["SyncNoiseEvents"][wh]->Fill(chId.sector(), chId.station());
0506       }
0507     }
0508 
0509     // fill # of noisy ch per wheel plot
0510     map<int, int>::const_iterator chMapIt = chMap.begin();
0511     map<int, int>::const_iterator chMapEnd = chMap.end();
0512     for (; chMapIt != chMapEnd; ++chMapIt) {
0513       wheelHistos["SyncNoiseChambs"][(*chMapIt).first]->Fill((*chMapIt).second);
0514     }
0515 
0516     // clear the map of # of digis per chamber: not needed anymore
0517     hitMap.clear();
0518 
0519     if (!syncNoisyChambers.empty()) {
0520       LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask") << "[DTDigiTask] Synch Noise in event: " << nevents;
0521       if (filterSyncNoise)
0522         LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask")
0523             << "\tnoisy time-boxes and occupancy will not be filled!" << endl;
0524       syncNumTot++;
0525       syncNum++;
0526     }
0527 
0528     // Logging of "large" synch Noisy events in private DQM
0529     if (syncNoisyChambers.size() > 3) {
0530       time_t eventTime = time_t(event.time().value() >> 32);
0531 
0532       LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask|DTSynchNoise")
0533           << "[DTDigiTask] At least 4 Synch Noisy chambers in Run : " << event.id().run()
0534           << " Lumi : " << event.id().luminosityBlock() << " Event : " << event.id().event()
0535           << " at time : " << ctime(&eventTime) << endl;
0536 
0537       set<DTChamberId>::const_iterator chIt = syncNoisyChambers.begin();
0538       set<DTChamberId>::const_iterator chEnd = syncNoisyChambers.end();
0539 
0540       stringstream synchNoisyCh;
0541       for (; chIt != chEnd; ++chIt) {
0542         synchNoisyCh << " " << (*chIt);
0543       }
0544       LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask|DTSynchNoise")
0545           << "[DTDigiTask] Chamber List :" << synchNoisyCh.str() << endl;
0546     }
0547 
0548     if (nevents % 1000 == 0) {
0549       LogVerbatim("DTDQM|DTMonitorModule|DTDigiTask")
0550           << (syncNumTot * 100. / nevents) << "% sync noise events since the beginning \n"
0551           << (syncNum * 0.1) << "% sync noise events in the last 1000 events " << endl;
0552       syncNum = 0;
0553     }
0554   }
0555 
0556   bool isSyncNoisy = false;
0557 
0558   DTDigiCollection::DigiRangeIterator dtLayerId_It;
0559   for (dtLayerId_It = dtdigis->begin(); dtLayerId_It != dtdigis->end(); ++dtLayerId_It) {  // Loop over layers
0560     isSyncNoisy = false;
0561     // check if chamber labeled as synch noisy
0562     if (filterSyncNoise) {
0563       DTChamberId chId = ((*dtLayerId_It).first).chamberId();
0564       if (syncNoisyChambers.find(chId) != syncNoisyChambers.end()) {
0565         isSyncNoisy = true;
0566       }
0567     }
0568 
0569     for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
0570          digiIt != ((*dtLayerId_It).second).second;
0571          ++digiIt) {  // Loop over all digis
0572 
0573       bool isNoisy = false;
0574       bool isFEMasked = false;
0575       bool isTDCMasked = false;
0576       bool isTrigMask = false;
0577       bool isDead = false;
0578       bool isNohv = false;
0579       if (checkNoisyChannels) {
0580         const DTWireId wireId(((*dtLayerId_It).first), (*digiIt).wire());
0581         statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
0582       }
0583 
0584       // Get the useful IDs
0585       const DTSuperLayerId dtSLId = ((*dtLayerId_It).first).superlayerId();
0586       uint32_t indexSL = dtSLId.rawId();
0587       const DTChamberId dtChId = dtSLId.chamberId();
0588       uint32_t indexCh = dtChId.rawId();
0589       int layer_number = ((*dtLayerId_It).first).layer();
0590       int superlayer_number = dtSLId.superlayer();
0591 
0592       // Read the ttrig DB or set a rough value from config
0593       // ttrig and rms are TDC counts
0594       if (readTTrigDB)
0595         tTrigMap->get(((*dtLayerId_It).first).superlayerId(), tTrig, tTrigRMS, kFactor, DTTimeUnits::counts);
0596       else
0597         tTrig = defaultTTrig;
0598 
0599       int inTimeHitsLowerBoundCorr = int(round(tTrig)) - inTimeHitsLowerBound;
0600       int inTimeHitsUpperBoundCorr = int(round(tTrig)) + tMax + inTimeHitsUpperBound;
0601 
0602       float t0;
0603       float t0RMS;
0604       int tdcTime = (*digiIt).countsTDC();
0605 
0606       if (subtractT0) {
0607         const DTWireId dtWireId(((*dtLayerId_It).first), (*digiIt).wire());
0608         // t0s and rms are TDC counts
0609         t0Map->get(dtWireId, t0, t0RMS, DTTimeUnits::counts);
0610         tdcTime += int(round(t0));
0611       }
0612 
0613       if (sliceTestMode) {
0614         tdcTime = std::max(tdcPedestal + 1, std::min(tdcPedestal + maxTTMounts - 1, tdcTime));
0615         // std::cout << tdcTime << std::endl;
0616       }
0617 
0618       // Fill Time-Boxes
0619       // NOTE: avoid to fill TB and PhotoPeak with noise. Occupancy are filled anyway
0620       if ((!isNoisy) && (!isSyncNoisy)) {  // Discard noisy channels
0621         // TimeBoxes per SL
0622         histoTag = "TimeBox" + triggerSource();
0623 
0624         (digiHistos.find(histoTag)->second).find(indexSL)->second->Fill(tdcTime);
0625         if (doLayerTimeBoxes)
0626           (digiHistos.find(histoTag)->second).find((*dtLayerId_It).first.rawId())->second->Fill(tdcTime);
0627       }
0628 
0629       // Fill Occupancies
0630       if (!isSyncNoisy) {  // Discard synch noisy channels
0631 
0632         if (doAllHitsOccupancies) {  // fill occupancies for all hits
0633           //Occupancies per chamber & layer
0634           histoTag = "OccupancyAllHits_perCh";
0635           map<uint32_t, MonitorElement*>::const_iterator mappedHisto = digiHistos[histoTag].find(indexCh);
0636 
0637           //FR comment the following cannot pass ibooker to analyze method!
0638           /*
0639           if (mappedHisto == digiHistos[histoTag].end()) { // dynamic booking
0640             bookHistos(ibooker, dtChId, string("Occupancies"), histoTag);
0641             mappedHisto = digiHistos[histoTag].find(indexCh);
0642           }
0643       */
0644           mappedHisto->second->Fill((*digiIt).wire(), (layer_number + (superlayer_number - 1) * 4) - 1);
0645 
0646           // Fill the chamber occupancy
0647           histoTag = "OccupancyAllHits";
0648           map<int, MonitorElement*>::const_iterator histoPerWheel = wheelHistos[histoTag].find(dtChId.wheel());
0649 
0650           histoPerWheel->second->Fill(dtChId.sector(), dtChId.station());  // FIXME: normalize to # of layers
0651         }
0652 
0653         if (doNoiseOccupancies) {  // fill occupancies for hits before the ttrig
0654           if (tdcTime < inTimeHitsLowerBoundCorr) {
0655             // FIXME: what about tdcTime > inTimeHitsUpperBoundCorr ???
0656 
0657             // Noise: Before tTrig
0658             //Occupancies Noise per chamber & layer
0659             histoTag = "OccupancyNoise_perCh";
0660             map<uint32_t, MonitorElement*>::const_iterator mappedHisto = digiHistos[histoTag].find(indexCh);
0661 
0662             mappedHisto->second->Fill((*digiIt).wire(), (layer_number + (superlayer_number - 1) * 4) - 1);
0663 
0664             // Fill the chamber occupancy
0665 
0666             histoTag = "OccupancyNoise";
0667             map<int, MonitorElement*>::const_iterator histoPerWheel = wheelHistos[histoTag].find(dtChId.wheel());
0668 
0669             histoPerWheel->second->Fill(dtChId.sector(), dtChId.station());  // FIXME: normalize to # of layers
0670           }
0671         }
0672 
0673         if (doInTimeOccupancies) {  // fill occpunacies for in-time hits only
0674           if (tdcTime > inTimeHitsLowerBoundCorr && tdcTime < inTimeHitsUpperBoundCorr) {
0675             // Physical hits: within the time window
0676 
0677             //Occupancies Signal per chamber & layer
0678             histoTag = "OccupancyInTimeHits_perCh";
0679             map<uint32_t, MonitorElement*>::const_iterator mappedHisto = digiHistos[histoTag].find(indexCh);
0680 
0681             mappedHisto->second->Fill((*digiIt).wire(), (layer_number + (superlayer_number - 1) * 4) - 1);
0682 
0683             // Fill the chamber occupancy
0684             histoTag = "OccupancyInTimeHits";
0685             map<int, MonitorElement*>::const_iterator histoPerWheel = wheelHistos[histoTag].find(dtChId.wheel());
0686 
0687             histoPerWheel->second->Fill(dtChId.sector(), dtChId.station());  // FIXME: normalize to # of layers
0688           }
0689         }
0690       }
0691     }
0692   }
0693 
0694   syncNoisyChambers.clear();
0695 }
0696 
0697 string DTDigiTask::triggerSource() {
0698   string l1ASource;
0699   if (isLocalRun)
0700     return l1ASource;
0701 
0702   for (std::vector<LTCDigi>::const_iterator ltc_it = ltcdigis->begin(); ltc_it != ltcdigis->end(); ltc_it++) {
0703     size_t otherTriggerSum = 0;
0704     for (size_t i = 1; i < 6; i++)
0705       otherTriggerSum += size_t((*ltc_it).HasTriggered(i));
0706 
0707     if ((*ltc_it).HasTriggered(0) && otherTriggerSum == 0)
0708       l1ASource = "DTonly";
0709     else if (!(*ltc_it).HasTriggered(0))
0710       l1ASource = "NoDT";
0711     else if ((*ltc_it).HasTriggered(0) && otherTriggerSum > 0)
0712       l1ASource = "DTalso";
0713   }
0714 
0715   return l1ASource;
0716 }
0717 
0718 string DTDigiTask::topFolder() const {
0719   if (tpMode)
0720     return string("DT/10-TestPulses/");
0721   else if (sliceTestMode)
0722     return string("DT/01-SliceTestDigi/");
0723   return string("DT/01-Digi/");
0724 }
0725 
0726 void DTDigiTask::channelsMap(const DTChamberId& dtCh, string histoTag) {
0727   // n max channels
0728   int nWires_max = (digiHistos[histoTag])[dtCh.rawId()]->getNbinsX();
0729 
0730   // set bin content = -1 for each not real channel. For visualization purposes
0731   for (int sl = 1; sl <= 3; sl++) {
0732     for (int ly = 1; ly <= 4; ly++) {
0733       for (int ch = 1; ch <= nWires_max; ch++) {
0734         int dduId = -1, rosId = -1, robId = -1, tdcId = -1, channelId = -1;
0735         int realCh = mapping->geometryToReadOut(
0736             dtCh.wheel(), dtCh.station(), dtCh.sector(), sl, ly, ch, dduId, rosId, robId, tdcId, channelId);
0737 
0738         // realCh = 0 if the channel exists, while realCh = 1 if it does not exist
0739         if (realCh) {
0740           int lybin = (4 * sl - 4) + ly;
0741           (digiHistos[histoTag])[dtCh.rawId()]->setBinContent(ch, lybin, -1.);
0742         }
0743       }
0744     }
0745   }
0746 }
0747 
0748 // Local Variables:
0749 // show-trailing-whitespace: t
0750 // truncate-lines: t
0751 // End: