Back to home page

Project CMSSW displayed by LXR

 
 

    


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

0001 
0002 /*
0003  *  See header file for a description of this class.
0004  *
0005  *  \author G. Mila - INFN Torino
0006  */
0007 
0008 #include "DQMOffline/CalibMuon/interface/DTnoiseDBValidation.h"
0009 
0010 // Framework
0011 #include "FWCore/Framework/interface/ESHandle.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/EventSetup.h"
0014 #include "FWCore/ServiceRegistry/interface/Service.h"
0015 
0016 #include "DQMServices/Core/interface/DQMStore.h"
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 
0019 // Geometry
0020 #include "DataFormats/MuonDetId/interface/DTLayerId.h"
0021 #include "DataFormats/MuonDetId/interface/DTWireId.h"
0022 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0023 #include "Geometry/DTGeometry/interface/DTTopology.h"
0024 
0025 // Noise record
0026 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
0027 
0028 #include "TFile.h"
0029 #include "TH1F.h"
0030 #include <cmath>
0031 #include <cstdio>
0032 #include <sstream>
0033 
0034 using namespace edm;
0035 using namespace std;
0036 
0037 DTnoiseDBValidation::DTnoiseDBValidation(const ParameterSet &pset)
0038     : labelDBRef_(esConsumes(edm::ESInputTag("", pset.getParameter<string>("labelDBRef")))),
0039       labelDB_(esConsumes(edm::ESInputTag("", pset.getParameter<string>("labelDB")))),
0040       muonGeomToken_(esConsumes<edm::Transition::BeginRun>()) {
0041   LogVerbatim("NoiseDBValidation") << "[DTnoiseDBValidation] Constructor called!";
0042 
0043   // Get the DQM needed services
0044   usesResource("DQMStore");
0045   dbe_ = edm::Service<DQMStore>().operator->();
0046   dbe_->setCurrentFolder("DT/DtCalib/NoiseDBValidation");
0047 
0048   diffTestName_ = "noiseDifferenceInRange";
0049   if (pset.exists("diffTestName"))
0050     diffTestName_ = pset.getParameter<string>("diffTestName");
0051 
0052   wheelTestName_ = "noiseWheelOccInRange";
0053   if (pset.exists("wheelTestName"))
0054     wheelTestName_ = pset.getParameter<string>("wheelTestName");
0055 
0056   stationTestName_ = "noiseStationOccInRange";
0057   if (pset.exists("stationTestName"))
0058     stationTestName_ = pset.getParameter<string>("stationTestName");
0059 
0060   sectorTestName_ = "noiseSectorOccInRange";
0061   if (pset.exists("sectorTestName"))
0062     sectorTestName_ = pset.getParameter<string>("sectorTestName");
0063 
0064   layerTestName_ = "noiseLayerOccInRange";
0065   if (pset.exists("layerTestName"))
0066     layerTestName_ = pset.getParameter<string>("layerTestName");
0067 
0068   outputMEsInRootFile_ = false;
0069   if (pset.exists("OutputFileName")) {
0070     outputMEsInRootFile_ = true;
0071     outputFileName_ = pset.getParameter<std::string>("OutputFileName");
0072   }
0073 }
0074 
0075 DTnoiseDBValidation::~DTnoiseDBValidation() {}
0076 
0077 void DTnoiseDBValidation::beginRun(const edm::Run &run, const EventSetup &setup) {
0078   noiseRefMap_ = &setup.getData(labelDBRef_);
0079 
0080   noiseMap_ = &setup.getData(labelDB_);
0081   ;
0082 
0083   // Get the geometry
0084   dtGeom = &setup.getData(muonGeomToken_);
0085 
0086   LogVerbatim("NoiseDBValidation") << "[DTnoiseDBValidation] Parameters initialization";
0087 
0088   noisyCellsRef_ = 0;
0089   noisyCellsValid_ = 0;
0090 
0091   // Histo booking
0092   diffHisto_ =
0093       dbe_->book1D("noisyCellDiff", "percentual (wrt the previous db) total number of noisy cells", 1, 0.5, 1.5);
0094   diffHisto_->setBinLabel(1, "Diff");
0095   wheelHisto_ = dbe_->book1D("wheelOccupancy", "percentual noisy cells occupancy per wheel", 5, -2.5, 2.5);
0096   wheelHisto_->setBinLabel(1, "Wh-2");
0097   wheelHisto_->setBinLabel(2, "Wh-1");
0098   wheelHisto_->setBinLabel(3, "Wh0");
0099   wheelHisto_->setBinLabel(4, "Wh1");
0100   wheelHisto_->setBinLabel(5, "Wh2");
0101   stationHisto_ = dbe_->book1D("stationOccupancy", "percentual noisy cells occupancy per station", 4, 0.5, 4.5);
0102   stationHisto_->setBinLabel(1, "St1");
0103   stationHisto_->setBinLabel(2, "St2");
0104   stationHisto_->setBinLabel(3, "St3");
0105   stationHisto_->setBinLabel(4, "St4");
0106   sectorHisto_ = dbe_->book1D("sectorOccupancy", "percentual noisy cells occupancy per sector", 12, 0.5, 12.5);
0107   sectorHisto_->setBinLabel(1, "Sect1");
0108   sectorHisto_->setBinLabel(2, "Sect2");
0109   sectorHisto_->setBinLabel(3, "Sect3");
0110   sectorHisto_->setBinLabel(4, "Sect4");
0111   sectorHisto_->setBinLabel(5, "Sect5");
0112   sectorHisto_->setBinLabel(6, "Sect6");
0113   sectorHisto_->setBinLabel(7, "Sect7");
0114   sectorHisto_->setBinLabel(8, "Sect8");
0115   sectorHisto_->setBinLabel(9, "Sect9");
0116   sectorHisto_->setBinLabel(10, "Sect10");
0117   sectorHisto_->setBinLabel(11, "Sect11");
0118   sectorHisto_->setBinLabel(12, "Sect12");
0119   layerHisto_ = dbe_->book1D("layerOccupancy", "percentual noisy cells occupancy per layer", 3, 0.5, 3.5);
0120   layerHisto_->setBinLabel(1, "First 10 bins");
0121   layerHisto_->setBinLabel(2, "Middle bins");
0122   layerHisto_->setBinLabel(3, "Last 10 bins");
0123 
0124   // map initialization
0125   map<int, int> whMap;
0126   whMap.clear();
0127   map<int, int> stMap;
0128   stMap.clear();
0129   map<int, int> sectMap;
0130   sectMap.clear();
0131   map<int, int> layerMap;
0132   layerMap.clear();
0133 
0134   // Loop over reference DB entries
0135   for (DTStatusFlag::const_iterator noise = noiseRefMap_->begin(); noise != noiseRefMap_->end(); noise++) {
0136     DTWireId wireId((*noise).first.wheelId,
0137                     (*noise).first.stationId,
0138                     (*noise).first.sectorId,
0139                     (*noise).first.slId,
0140                     (*noise).first.layerId,
0141                     (*noise).first.cellId);
0142     LogVerbatim("NoiseDBValidation") << "Ref. noisy wire: " << wireId;
0143     ++noisyCellsRef_;
0144   }
0145 
0146   // Loop over validation DB entries
0147   for (DTStatusFlag::const_iterator noise = noiseMap_->begin(); noise != noiseMap_->end(); noise++) {
0148     DTWireId wireId((*noise).first.wheelId,
0149                     (*noise).first.stationId,
0150                     (*noise).first.sectorId,
0151                     (*noise).first.slId,
0152                     (*noise).first.layerId,
0153                     (*noise).first.cellId);
0154     LogVerbatim("NoiseDBValidation") << "Valid. noisy wire: " << wireId;
0155     ++noisyCellsValid_;
0156 
0157     whMap[(*noise).first.wheelId]++;
0158     stMap[(*noise).first.stationId]++;
0159     sectMap[(*noise).first.sectorId]++;
0160 
0161     const DTTopology &dtTopo = dtGeom->layer(wireId.layerId())->specificTopology();
0162     const int lastWire = dtTopo.lastChannel();
0163     if ((*noise).first.cellId <= 10)
0164       layerMap[1]++;
0165     if ((*noise).first.cellId > 10 && (*noise).first.cellId < (lastWire - 10))
0166       layerMap[2]++;
0167     if ((*noise).first.cellId >= (lastWire - 10))
0168       layerMap[3]++;
0169 
0170     const DTChamberId chId = wireId.layerId().superlayerId().chamberId();
0171     if (noiseHistoMap_.find(chId) == noiseHistoMap_.end())
0172       bookHisto(chId);
0173     int binNumber = 4 * (wireId.superLayer() - 1) + wireId.layer();
0174     noiseHistoMap_[chId]->Fill(wireId.wire(), binNumber);
0175   }
0176 
0177   // histo filling
0178   double scale = 1 / double(noisyCellsRef_);
0179   diffHisto_->Fill(1, abs(noisyCellsRef_ - noisyCellsValid_) * scale);
0180 
0181   scale = 1 / double(noisyCellsValid_);
0182   for (map<int, int>::const_iterator wheel = whMap.begin(); wheel != whMap.end(); wheel++) {
0183     wheelHisto_->Fill((*wheel).first, ((*wheel).second) * scale);
0184   }
0185   for (map<int, int>::const_iterator station = stMap.begin(); station != stMap.end(); station++) {
0186     stationHisto_->Fill((*station).first, ((*station).second) * scale);
0187   }
0188   for (map<int, int>::const_iterator sector = sectMap.begin(); sector != sectMap.end(); sector++) {
0189     sectorHisto_->Fill((*sector).first, ((*sector).second) * scale);
0190   }
0191   for (map<int, int>::const_iterator layer = layerMap.begin(); layer != layerMap.end(); layer++) {
0192     layerHisto_->Fill((*layer).first, ((*layer).second) * scale);
0193   }
0194 }
0195 
0196 void DTnoiseDBValidation::endRun(edm::Run const &run, edm::EventSetup const &setup) {
0197   // test on difference histo
0198   // string testCriterionName;
0199   // testCriterionName =
0200   // parameters.getUntrackedParameter<string>("diffTestName","noiseDifferenceInRange");
0201   const QReport *theDiffQReport = diffHisto_->getQReport(diffTestName_);
0202   if (theDiffQReport) {
0203     vector<dqm::me_util::Channel> badChannels = theDiffQReport->getBadChannels();
0204     for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0205          channel++) {
0206       LogWarning("NoiseDBValidation") << " Bad partial difference of noisy channels! Contents : "
0207                                       << (*channel).getContents();
0208     }
0209   }
0210   // testCriterionName =
0211   // parameters.getUntrackedParameter<string>("wheelTestName","noiseWheelOccInRange");
0212   const QReport *theDiffQReport2 = wheelHisto_->getQReport(wheelTestName_);
0213   if (theDiffQReport2) {
0214     vector<dqm::me_util::Channel> badChannels = theDiffQReport2->getBadChannels();
0215     for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0216          channel++) {
0217       int wheel = (*channel).getBin() - 3;
0218       LogWarning("NoiseDBValidation") << " Bad percentual occupancy for wheel : " << wheel
0219                                       << "  Contents : " << (*channel).getContents();
0220     }
0221   }
0222   // testCriterionName =
0223   // parameters.getUntrackedParameter<string>("stationTestName","noiseStationOccInRange");
0224   const QReport *theDiffQReport3 = stationHisto_->getQReport(stationTestName_);
0225   if (theDiffQReport3) {
0226     vector<dqm::me_util::Channel> badChannels = theDiffQReport3->getBadChannels();
0227     for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0228          channel++) {
0229       LogWarning("NoiseDBValidation") << " Bad percentual occupancy for station : " << (*channel).getBin()
0230                                       << "  Contents : " << (*channel).getContents();
0231     }
0232   }
0233   // testCriterionName =
0234   // parameters.getUntrackedParameter<string>("sectorTestName","noiseSectorOccInRange");
0235   const QReport *theDiffQReport4 = sectorHisto_->getQReport(sectorTestName_);
0236   if (theDiffQReport4) {
0237     vector<dqm::me_util::Channel> badChannels = theDiffQReport4->getBadChannels();
0238     for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0239          channel++) {
0240       LogWarning("NoiseDBValidation") << " Bad percentual occupancy for sector : " << (*channel).getBin()
0241                                       << "  Contents : " << (*channel).getContents();
0242     }
0243   }
0244   // testCriterionName =
0245   // parameters.getUntrackedParameter<string>("layerTestName","noiseLayerOccInRange");
0246   const QReport *theDiffQReport5 = layerHisto_->getQReport(layerTestName_);
0247   if (theDiffQReport5) {
0248     vector<dqm::me_util::Channel> badChannels = theDiffQReport5->getBadChannels();
0249     for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin(); channel != badChannels.end();
0250          channel++) {
0251       if ((*channel).getBin() == 1)
0252         LogWarning("NoiseDBValidation") << " Bad percentual occupancy for the first 10 wires! Contents : "
0253                                         << (*channel).getContents();
0254       if ((*channel).getBin() == 2)
0255         LogWarning("NoiseDBValidation") << " Bad percentual occupancy for the middle wires! Contents : "
0256                                         << (*channel).getContents();
0257       if ((*channel).getBin() == 3)
0258         LogWarning("NoiseDBValidation") << " Bad percentual occupancy for the last 10 wires! Contents : "
0259                                         << (*channel).getContents();
0260     }
0261   }
0262 }
0263 
0264 void DTnoiseDBValidation::endJob() {
0265   // Write the histos in a ROOT file
0266   if (outputMEsInRootFile_)
0267     dbe_->save(outputFileName_);
0268 }
0269 
0270 void DTnoiseDBValidation::bookHisto(const DTChamberId &chId) {
0271   stringstream histoName;
0272   histoName << "NoiseOccupancy"
0273             << "_W" << chId.wheel() << "_St" << chId.station() << "_Sec" << chId.sector();
0274 
0275   if (noiseHistoMap_.find(chId) == noiseHistoMap_.end()) {  // Redundant check
0276     // Get the chamber from the geometry
0277     int nWiresMax = 0;
0278     const DTChamber *dtchamber = dtGeom->chamber(chId);
0279     const vector<const DTSuperLayer *> &superlayers = dtchamber->superLayers();
0280 
0281     // Loop over layers and find the max # of wires
0282     for (vector<const DTSuperLayer *>::const_iterator sl = superlayers.begin(); sl != superlayers.end();
0283          ++sl) {  // loop over SLs
0284       vector<const DTLayer *> layers = (*sl)->layers();
0285       for (vector<const DTLayer *>::const_iterator lay = layers.begin(); lay != layers.end();
0286            ++lay) {  // loop over layers
0287         int nWires = (*lay)->specificTopology().channels();
0288         if (nWires > nWiresMax)
0289           nWiresMax = nWires;
0290       }
0291     }
0292 
0293     noiseHistoMap_[chId] = dbe_->book2D(histoName.str(), "Noise occupancy", nWiresMax, 1, (nWiresMax + 1), 12, 1, 13);
0294     for (int i_sl = 1; i_sl <= 3; ++i_sl) {
0295       for (int i_lay = 1; i_lay <= 4; ++i_lay) {
0296         int binNumber = 4 * (i_sl - 1) + i_lay;
0297         stringstream label;
0298         label << "SL" << i_sl << ": L" << i_lay;
0299         noiseHistoMap_[chId]->setBinLabel(binNumber, label.str(), 2);
0300       }
0301     }
0302   }
0303 }