00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "DQMOffline/CalibMuon/interface/DTnoiseDBValidation.h"
00011
00012
00013 #include "FWCore/Framework/interface/Event.h"
00014 #include "FWCore/Framework/interface/EventSetup.h"
00015 #include "FWCore/Framework/interface/ESHandle.h"
00016 #include "FWCore/ServiceRegistry/interface/Service.h"
00017
00018 #include "DQMServices/Core/interface/DQMStore.h"
00019 #include "DQMServices/Core/interface/MonitorElement.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021
00022
00023 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00024 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00025 #include "Geometry/DTGeometry/interface/DTLayer.h"
00026 #include "Geometry/DTGeometry/interface/DTTopology.h"
00027
00028
00029 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00030 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00031
00032
00033 #include <stdio.h>
00034 #include <sstream>
00035 #include <math.h>
00036 #include "TFile.h"
00037 #include "TH1F.h"
00038
00039 using namespace edm;
00040 using namespace std;
00041
00042
00043
00044
00045 DTnoiseDBValidation::DTnoiseDBValidation(const ParameterSet& pset) {
00046
00047 cout << "[DTnoiseDBValidation] Constructor called!" << endl;
00048
00049
00050 dbe = edm::Service<DQMStore>().operator->();
00051 dbe->setCurrentFolder("DT/noiseDBValidation");
00052
00053
00054 labelDBRef = pset.getUntrackedParameter<string>("labelDBRef");
00055 labelDB = pset.getUntrackedParameter<string>("labelDB");
00056
00057 parameters = pset;
00058 }
00059
00060
00061 DTnoiseDBValidation::~DTnoiseDBValidation(){}
00062
00063 void DTnoiseDBValidation::beginRun(const edm::Run& run, const EventSetup& setup) {
00064 ESHandle<DTStatusFlag> noise_Ref;
00065 setup.get<DTStatusFlagRcd>().get(labelDBRef, noise_Ref);
00066 noiseRefMap = &*noise_Ref;
00067
00068 ESHandle<DTStatusFlag> noise_toTest;
00069 setup.get<DTStatusFlagRcd>().get(labelDB, noise_toTest);
00070 noiseMap = &*noise_toTest;
00071
00072
00073 setup.get<MuonGeometryRecord>().get(dtGeom);
00074 }
00075
00076 void DTnoiseDBValidation::beginJob() {
00077
00078
00079 metname = "noiseDbValidation";
00080 LogTrace(metname)<<"[DTnoiseDBValidation] Parameters initialization";
00081
00082 outputFileName = parameters.getUntrackedParameter<std::string>("OutputFileName");
00083
00084 noisyCells_Ref=0;
00085 noisyCells_toTest=0;
00086
00087
00088 diffHisto = dbe->book1D("noisyCellDiff", "percentual (wrt the previous db) total number of noisy cells",1, 0.5, 1.5);
00089 diffHisto->setBinLabel(1,"diff");
00090 wheelHisto = dbe->book1D("wheelOccupancy", "percentual noisy cells occupancy per wheel",5, -2.5, 2.5);
00091 wheelHisto->setBinLabel(1,"wh-2");
00092 wheelHisto->setBinLabel(2,"wh-1");
00093 wheelHisto->setBinLabel(3,"wh0");
00094 wheelHisto->setBinLabel(4,"wh1");
00095 wheelHisto->setBinLabel(5,"wh2");
00096 stationHisto = dbe->book1D("stationOccupancy", "percentual noisy cells occupancy per station",4, 0.5, 4.5);
00097 stationHisto->setBinLabel(1,"st1");
00098 stationHisto->setBinLabel(2,"st2");
00099 stationHisto->setBinLabel(3,"st3");
00100 stationHisto->setBinLabel(4,"st4");
00101 sectorHisto = dbe->book1D("sectorOccupancy", "percentual noisy cells occupancy per sector",12, 0.5, 12.5);
00102 sectorHisto->setBinLabel(1,"sect1");
00103 sectorHisto->setBinLabel(2,"sect2");
00104 sectorHisto->setBinLabel(3,"sect3");
00105 sectorHisto->setBinLabel(4,"sect4");
00106 sectorHisto->setBinLabel(5,"sect5");
00107 sectorHisto->setBinLabel(6,"sect6");
00108 sectorHisto->setBinLabel(7,"sect7");
00109 sectorHisto->setBinLabel(8,"sect8");
00110 sectorHisto->setBinLabel(9,"sect9");
00111 sectorHisto->setBinLabel(10,"sect10");
00112 sectorHisto->setBinLabel(11,"sect11");
00113 sectorHisto->setBinLabel(12,"sect12");
00114 layerHisto = dbe->book1D("layerOccupancy", "percentual noisy cells occupancy per layer",3, 0.5, 3.5);
00115 layerHisto->setBinLabel(1,"first 10 bins");
00116 layerHisto->setBinLabel(2,"middle bins");
00117 layerHisto->setBinLabel(3,"last 10 bins");
00118
00119
00120 map<int, int> whMap;
00121 whMap.clear();
00122 map<int, int> stMap;
00123 stMap.clear();
00124 map<int, int> sectMap;
00125 sectMap.clear();
00126 map<int, int> layerMap;
00127 layerMap.clear();
00128
00129
00130 for(DTStatusFlag::const_iterator noise = noiseRefMap->begin();
00131 noise != noiseRefMap->end(); noise++) {
00132 DTWireId wireId((*noise).first.wheelId,
00133 (*noise).first.stationId,
00134 (*noise).first.sectorId,
00135 (*noise).first.slId,
00136 (*noise).first.layerId,
00137 (*noise).first.cellId);
00138 cout<< "Ref Wire: " << wireId<<endl;
00139 noisyCells_Ref++;
00140 }
00141
00142
00143 for(DTStatusFlag::const_iterator noise = noiseMap->begin();
00144 noise != noiseMap->end(); noise++) {
00145 DTWireId wireId((*noise).first.wheelId,
00146 (*noise).first.stationId,
00147 (*noise).first.sectorId,
00148 (*noise).first.slId,
00149 (*noise).first.layerId,
00150 (*noise).first.cellId);
00151 cout<< "toTest Wire: " << wireId<<endl;
00152 noisyCells_toTest++;
00153
00154 whMap[(*noise).first.wheelId]++;
00155 stMap[(*noise).first.stationId]++;
00156 sectMap[(*noise).first.sectorId]++;
00157
00158 const DTTopology& dtTopo = dtGeom->layer(wireId.layerId())->specificTopology();
00159 const int lastWire = dtTopo.lastChannel();
00160 if((*noise).first.cellId<=10)
00161 layerMap[1]++;
00162 if((*noise).first.cellId>10 && (*noise).first.cellId<(lastWire-10))
00163 layerMap[2]++;
00164 if((*noise).first.cellId>=(lastWire-10))
00165 layerMap[3]++;
00166
00167 }
00168
00169
00170 double scale = 1/double(noisyCells_Ref);
00171 diffHisto->Fill(1,abs(noisyCells_Ref-noisyCells_toTest)*scale);
00172
00173 scale = 1/double(noisyCells_toTest);
00174 for(map<int, int >::const_iterator wheel = whMap.begin();
00175 wheel != whMap.end();
00176 wheel++) {
00177 wheelHisto->Fill((*wheel).first, ((*wheel).second)*scale);
00178 }
00179 for(map<int, int >::const_iterator station = stMap.begin();
00180 station != stMap.end();
00181 station++) {
00182 stationHisto->Fill((*station).first, ((*station).second)*scale);
00183 }
00184 for(map<int, int >::const_iterator sector = sectMap.begin();
00185 sector != sectMap.end();
00186 sector++) {
00187 sectorHisto->Fill((*sector).first, ((*sector).second)*scale);
00188 }
00189 for(map<int, int >::const_iterator layer = layerMap.begin();
00190 layer != layerMap.end();
00191 layer++) {
00192 layerHisto->Fill((*layer).first, ((*layer).second)*scale);
00193 }
00194
00195 }
00196
00197 void DTnoiseDBValidation::endJob() {
00198
00199
00200 string testCriterionName = parameters.getUntrackedParameter<string>("diffTestName","noiseDifferenceInRange");
00201 const QReport * theDiffQReport = diffHisto->getQReport(testCriterionName);
00202 if(theDiffQReport) {
00203 vector<dqm::me_util::Channel> badChannels = theDiffQReport->getBadChannels();
00204 for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
00205 channel != badChannels.end(); channel++) {
00206 cout << " Bad partial difference of noisy channels! Contents : "<<(*channel).getContents()<<endl;
00207 }
00208 }
00209 testCriterionName = parameters.getUntrackedParameter<string>("wheelTestName","noiseWheelOccInRange");
00210 const QReport * theDiffQReport2 = wheelHisto->getQReport(testCriterionName);
00211 if(theDiffQReport2) {
00212 vector<dqm::me_util::Channel> badChannels = theDiffQReport2->getBadChannels();
00213 for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
00214 channel != badChannels.end(); channel++) {
00215 int wheel = (*channel).getBin()-3;
00216 cout << " Bad percentual occupancy for wheel : "<<wheel<<" Contents : "<<(*channel).getContents()<<endl;
00217 }
00218 }
00219 testCriterionName = parameters.getUntrackedParameter<string>("stationTestName","noiseStationOccInRange");
00220 const QReport * theDiffQReport3 = stationHisto->getQReport(testCriterionName);
00221 if(theDiffQReport3) {
00222 vector<dqm::me_util::Channel> badChannels = theDiffQReport3->getBadChannels();
00223 for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
00224 channel != badChannels.end(); channel++) {
00225 cout << " Bad percentual occupancy for station : "<<(*channel).getBin()<<" Contents : "<<(*channel).getContents()<<endl;
00226 }
00227 }
00228 testCriterionName = parameters.getUntrackedParameter<string>("sectorTestName","noiseSectorOccInRange");
00229 const QReport * theDiffQReport4 = sectorHisto->getQReport(testCriterionName);
00230 if(theDiffQReport4) {
00231 vector<dqm::me_util::Channel> badChannels = theDiffQReport4->getBadChannels();
00232 for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
00233 channel != badChannels.end(); channel++) {
00234 cout << " Bad percentual occupancy for sector : "<<(*channel).getBin()<<" Contents : "<<(*channel).getContents()<<endl;
00235 }
00236 }
00237 testCriterionName = parameters.getUntrackedParameter<string>("layerTestName","noiseLayerOccInRange");
00238 const QReport * theDiffQReport5 = layerHisto->getQReport(testCriterionName);
00239 if(theDiffQReport5) {
00240 vector<dqm::me_util::Channel> badChannels = theDiffQReport5->getBadChannels();
00241 for (vector<dqm::me_util::Channel>::iterator channel = badChannels.begin();
00242 channel != badChannels.end(); channel++) {
00243 if((*channel).getBin()==1)
00244 cout << " Bad percentual occupancy for the first 10 wires! Contents : "<<(*channel).getContents()<<endl;
00245 if((*channel).getBin()==2)
00246 cout << " Bad percentual occupancy for the middle wires! Contents : "<<(*channel).getContents()<<endl;
00247 if((*channel).getBin()==3)
00248 cout << " Bad percentual occupancy for the last 10 wires! Contents : "<<(*channel).getContents()<<endl;
00249 }
00250 }
00251
00252
00253 dbe->save(outputFileName);
00254
00255 }