Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 /*
0003  *  See header file for a description of this class.
0004  *
0005  *  \authors G. Mila , G. Cerminara - INFN Torino
0006  */
0007 
0008 #include "DTNoiseTask.h"
0009 
0010 // Framework
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/ServiceRegistry/interface/Service.h"
0015 
0016 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0017 
0018 // Geometry
0019 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0020 #include "Geometry/DTGeometry/interface/DTLayer.h"
0021 #include "Geometry/DTGeometry/interface/DTTopology.h"
0022 
0023 // Digi
0024 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
0025 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
0026 #include "CondFormats/DTObjects/interface/DTReadOutMapping.h"
0027 
0028 // Database
0029 #include "CondFormats/DTObjects/interface/DTTtrig.h"
0030 
0031 #include <sstream>
0032 #include <string>
0033 
0034 using namespace edm;
0035 using namespace std;
0036 
0037 DTNoiseTask::DTNoiseTask(const ParameterSet& ps)
0038     : evtNumber(0),
0039       muonGeomToken_(esConsumes<edm::Transition::BeginRun>()),
0040       tTrigMapToken_(esConsumes<edm::Transition::BeginRun>()) {
0041   LogVerbatim("DTNoiseTask") << "[DTNoiseTask]: Constructor" << endl;
0042 
0043   //switch for timeBox booking
0044   doTimeBoxHistos = ps.getUntrackedParameter<bool>("doTbHistos", false);
0045 
0046   // The label to retrieve the digis
0047   dtDigiToken_ = consumes<DTDigiCollection>(ps.getUntrackedParameter<InputTag>("dtDigiLabel"));
0048 
0049   // the name of the 4D rec hits collection
0050   recHits4DToken_ =
0051       consumes<DTRecSegment4DCollection>(edm::InputTag(ps.getUntrackedParameter<string>("recHits4DLabel")));
0052 
0053   // switch for segment veto
0054   doSegmentVeto = ps.getUntrackedParameter<bool>("doSegmentVeto", false);
0055 
0056   // safe margin (ns) between ttrig and beginning of counting area
0057   safeMargin = ps.getUntrackedParameter<double>("safeMargin", 200.);
0058 }
0059 
0060 DTNoiseTask::~DTNoiseTask() {}
0061 
0062 /// Analyze
0063 void DTNoiseTask::analyze(const edm::Event& e, const edm::EventSetup& c) {
0064   evtNumber++;
0065   nEventMonitor->Fill(evtNumber);
0066 
0067   if (evtNumber % 1000 == 0)
0068     LogVerbatim("DTNoiseTask") << "[DTNoiseTask]: Analyzing evt number :" << evtNumber << endl;
0069 
0070   // map of the chambers with at least 1 segment
0071   std::map<DTChamberId, int> segmentsChId;
0072 
0073   // Get the 4D segment collection from the event
0074   edm::Handle<DTRecSegment4DCollection> all4DSegments;
0075   if (doSegmentVeto) {
0076     e.getByToken(recHits4DToken_, all4DSegments);
0077 
0078     // Loop over all chambers containing a segment and look for the number of segments
0079     DTRecSegment4DCollection::id_iterator chamberId;
0080     for (chamberId = all4DSegments->id_begin(); chamberId != all4DSegments->id_end(); ++chamberId) {
0081       segmentsChId[*chamberId] = 1;
0082     }
0083   }
0084 
0085   // Get the digis from the event
0086   edm::Handle<DTDigiCollection> dtdigis;
0087   e.getByToken(dtDigiToken_, dtdigis);
0088 
0089   // LOOP OVER ALL THE DIGIS OF THE EVENT
0090   DTDigiCollection::DigiRangeIterator dtLayerId_It;
0091   for (dtLayerId_It = dtdigis->begin(); dtLayerId_It != dtdigis->end(); ++dtLayerId_It) {
0092     for (DTDigiCollection::const_iterator digiIt = ((*dtLayerId_It).second).first;
0093          digiIt != ((*dtLayerId_It).second).second;
0094          ++digiIt) {
0095       //Check the TDC trigger width
0096       int tdcTime = (*digiIt).countsTDC();
0097       double upperLimit = tTrigStMap[(*dtLayerId_It).first.superlayerId().chamberId()] - safeMargin;
0098       if (doTimeBoxHistos)
0099         tbHistos[(*dtLayerId_It).first.superlayerId()]->Fill(tdcTime);
0100       if (tdcTime > upperLimit)
0101         continue;
0102 
0103       //Check the chamber has no 4D segments (optional)
0104       if (doSegmentVeto && segmentsChId.find((*dtLayerId_It).first.superlayerId().chamberId()) != segmentsChId.end())
0105         continue;
0106 
0107       // fill the occupancy histo
0108       // FIXME: needs to be optimized: no need to rescale the histo for each digi
0109       TH2D* noise_root = noiseHistos[(*dtLayerId_It).first.superlayerId().chamberId()]->getTH2D();
0110       double normalization = 0;
0111       if (mapEvt.find((*dtLayerId_It).first.superlayerId().chamberId()) != mapEvt.end()) {
0112         LogVerbatim("DTNoiseTask") << " Last fill: # of events: "
0113                                    << mapEvt[(*dtLayerId_It).first.superlayerId().chamberId()] << endl;
0114         normalization = 1e-9 * upperLimit * mapEvt[(*dtLayerId_It).first.superlayerId().chamberId()];
0115         // revert back to # of entries
0116         noise_root->Scale(normalization);
0117       }
0118       int yBin = (*dtLayerId_It).first.layer() + (4 * ((*dtLayerId_It).first.superlayerId().superlayer() - 1));
0119       noise_root->Fill((*digiIt).wire(), yBin);
0120       // normalize the occupancy histo
0121       mapEvt[(*dtLayerId_It).first.superlayerId().chamberId()] = evtNumber;
0122       LogVerbatim("DTNoiseTask") << (*dtLayerId_It).first << " wire: " << (*digiIt).wire()
0123                                  << " # counts: " << noise_root->GetBinContent((*digiIt).wire(), yBin)
0124                                  << " Time interval: " << upperLimit << " # of events: " << evtNumber << endl;
0125       ;
0126       normalization = double(1e-9 * upperLimit * mapEvt[(*dtLayerId_It).first.superlayerId().chamberId()]);
0127       // update the rate
0128       noise_root->Scale(1. / normalization);
0129       LogVerbatim("DTNoiseTask") << "    noise rate: " << noise_root->GetBinContent((*digiIt).wire(), yBin) << endl;
0130     }
0131   }
0132 }
0133 
0134 void DTNoiseTask::bookHistos(DQMStore::IBooker& ibooker, DTChamberId chId) {
0135   // set the folder
0136   stringstream wheel;
0137   wheel << chId.wheel();
0138   stringstream station;
0139   station << chId.station();
0140   stringstream sector;
0141   sector << chId.sector();
0142 
0143   ibooker.setCurrentFolder("DT/05-Noise/Wheel" + wheel.str() +
0144                            //           "/Station" + station.str() +
0145                            "/Sector" + sector.str());
0146 
0147   // Build the histo name
0148   string histoName = string("NoiseRate") + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str();
0149 
0150   LogVerbatim("DTNoiseTask") << "[DTNoiseTask]: booking chamber histo:" << endl;
0151   LogVerbatim("DTNoiseTask") << "              folder "
0152                              << "DT/05-Noise/Wheel" + wheel.str() +
0153                                     //     "/Station" + station.str() +
0154                                     "/Sector" + sector.str() + "/"
0155                              << endl;
0156   LogVerbatim("DTNoiseTask") << "              histoName " << histoName << endl;
0157 
0158   // Get the chamber from the geometry
0159   int nWires_max = 0;
0160   const DTChamber* dtchamber = dtGeom->chamber(chId);
0161   const vector<const DTSuperLayer*>& superlayers = dtchamber->superLayers();
0162 
0163   // Loop over layers and find the max # of wires
0164   for (vector<const DTSuperLayer*>::const_iterator sl = superlayers.begin(); sl != superlayers.end();
0165        ++sl) {  // loop over SLs
0166     vector<const DTLayer*> layers = (*sl)->layers();
0167     for (vector<const DTLayer*>::const_iterator lay = layers.begin(); lay != layers.end(); ++lay) {  // loop over layers
0168       int nWires = (*lay)->specificTopology().channels();
0169       if (nWires > nWires_max)
0170         nWires_max = nWires;
0171     }
0172   }
0173 
0174   noiseHistos[chId] =
0175       ibooker.book2DD(histoName, "Noise rate (Hz) per channel", nWires_max, 1, nWires_max + 1, 12, 1, 13);
0176   noiseHistos[chId]->setAxisTitle("wire number", 1);
0177   noiseHistos[chId]->setBinLabel(1, "SL1-L1", 2);
0178   noiseHistos[chId]->setBinLabel(2, "SL1-L2", 2);
0179   noiseHistos[chId]->setBinLabel(3, "SL1-L3", 2);
0180   noiseHistos[chId]->setBinLabel(4, "SL1-L4", 2);
0181   noiseHistos[chId]->setBinLabel(5, "SL2-L1", 2);
0182   noiseHistos[chId]->setBinLabel(6, "SL2-L2", 2);
0183   noiseHistos[chId]->setBinLabel(7, "SL2-L3", 2);
0184   noiseHistos[chId]->setBinLabel(8, "SL2-L4", 2);
0185   noiseHistos[chId]->setBinLabel(9, "SL3-L1", 2);
0186   noiseHistos[chId]->setBinLabel(10, "SL3-L2", 2);
0187   noiseHistos[chId]->setBinLabel(11, "SL3-L3", 2);
0188   noiseHistos[chId]->setBinLabel(12, "SL3-L4", 2);
0189 }
0190 
0191 void DTNoiseTask::bookHistos(DQMStore::IBooker& ibooker, DTSuperLayerId slId) {
0192   // set the folder
0193   stringstream wheel;
0194   wheel << slId.chamberId().wheel();
0195   stringstream station;
0196   station << slId.chamberId().station();
0197   stringstream sector;
0198   sector << slId.chamberId().sector();
0199   stringstream superlayer;
0200   superlayer << slId.superlayer();
0201 
0202   ibooker.setCurrentFolder("DT/05-Noise/Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" + sector.str());
0203 
0204   // Build the histo name
0205   string histoName =
0206       string("TimeBox") + "_W" + wheel.str() + "_St" + station.str() + "_Sec" + sector.str() + "_SL" + superlayer.str();
0207 
0208   LogVerbatim("DTNoiseTask") << "[DTNoiseTask]: booking SL histo:" << endl;
0209   LogVerbatim("DTNoiseTask") << "              folder "
0210                              << "DT/05-Noise/Wheel" + wheel.str() + "/Station" + station.str() + "/Sector" +
0211                                     sector.str() + "/"
0212                              << endl;
0213   LogVerbatim("DTNoiseTask") << "              histoName " << histoName << endl;
0214 
0215   tbHistos[slId] = ibooker.book1D(histoName, "Time Box (TDC counts)", 1000, 0, 6000);
0216 }
0217 
0218 void DTNoiseTask::dqmBeginRun(const Run& run, const EventSetup& setup) {
0219   LogVerbatim("DTNoiseTask") << "[DTNoiseTask]: Begin of run" << endl;
0220 
0221   // tTrig Map
0222   tTrigMap = &setup.getData(tTrigMapToken_);
0223 
0224   // get the geometry
0225   dtGeom = &setup.getData(muonGeomToken_);
0226 }
0227 
0228 void DTNoiseTask::bookHistograms(DQMStore::IBooker& ibooker, edm::Run const& run, edm::EventSetup const& setup) {
0229   ibooker.setCurrentFolder("DT/EventInfo/Counters");
0230   nEventMonitor = ibooker.bookFloat("nProcessedEventsNoise");
0231 
0232   // Loop over all the chambers
0233   vector<const DTChamber*>::const_iterator ch_it = dtGeom->chambers().begin();
0234   vector<const DTChamber*>::const_iterator ch_end = dtGeom->chambers().end();
0235   for (; ch_it != ch_end; ++ch_it) {
0236     DTChamberId chId = (*ch_it)->id();
0237     // histo booking
0238     bookHistos(ibooker, chId);
0239     vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
0240     vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
0241     // Loop over the SLs
0242     for (; sl_it != sl_end; ++sl_it) {
0243       DTSuperLayerId slId = (*sl_it)->id();
0244       if (doTimeBoxHistos)
0245         bookHistos(ibooker, slId);
0246       float tTrig, tTrigRMS, kFactor;
0247       tTrigMap->get(slId, tTrig, tTrigRMS, kFactor, DTTimeUnits::ns);
0248       // tTrig mapping per station
0249       // check that the ttrig is the lowest of the 3 SLs
0250       if (tTrigStMap.find(chId) == tTrigStMap.end() ||
0251           (tTrigStMap.find(chId) != tTrigStMap.end() && tTrig < tTrigStMap[chId]))
0252         tTrigStMap[chId] = tTrig;
0253     }
0254   }
0255 }
0256 
0257 void DTNoiseTask::endLuminosityBlock(const LuminosityBlock& lumiSeg, const EventSetup& setup) {
0258   LogVerbatim("DTNoiseTask") << "[DTNoiseTask]: End LS, update rates in all histos" << endl;
0259 
0260   // update the rate of all histos (usefull for histos with few entries: they are not updated very often
0261   for (map<DTChamberId, MonitorElement*>::const_iterator meAndChamber = noiseHistos.begin();
0262        meAndChamber != noiseHistos.end();
0263        ++meAndChamber) {
0264     DTChamberId chId = (*meAndChamber).first;
0265     TH2D* noise_root = (*meAndChamber).second->getTH2D();
0266     double upperLimit = tTrigStMap[chId] - safeMargin;
0267 
0268     double normalization = 0;
0269     if (mapEvt.find(chId) != mapEvt.end()) {
0270       LogVerbatim("DTNoiseTask") << " Ch: " << chId << " Last fill: # of events: " << mapEvt[chId] << endl;
0271       normalization = 1e-9 * upperLimit * mapEvt[chId];
0272       // revert back to # of entries
0273       noise_root->Scale(normalization);
0274     }
0275     //check that event analyzed != 0 might happen oline
0276     if (evtNumber) {
0277       // set the # of events analyzed until this update
0278       LogVerbatim("DTNoiseTask") << "          Update for events: " << evtNumber << endl;
0279       mapEvt[chId] = evtNumber;
0280       // update the rate
0281       normalization = double(1e-9 * upperLimit * evtNumber);
0282       noise_root->Scale(1. / normalization);
0283     }
0284   }
0285 }
0286 
0287 // Local Variables:
0288 // show-trailing-whitespace: t
0289 // truncate-lines: t
0290 // End: