00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #include "DTSegmentAnalysisTask.h"
00012
00013
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/Framework/interface/EventSetup.h"
00016 #include "FWCore/Framework/interface/LuminosityBlock.h"
00017 #include "FWCore/Framework/interface/ESHandle.h"
00018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00019 #include "FWCore/ServiceRegistry/interface/Service.h"
00020
00021 #include "DQMServices/Core/interface/DQMStore.h"
00022 #include "DQMServices/Core/interface/MonitorElement.h"
00023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00024
00025
00026 #include "DataFormats/GeometryVector/interface/Pi.h"
00027 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00028 #include "Geometry/DTGeometry/interface/DTTopology.h"
00029 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00030
00031 #include "DataFormats/DTRecHit/interface/DTRecSegment4DCollection.h"
00032
00033 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00034 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00035
00036 #include "DQM/DTMonitorModule/interface/DTTimeEvolutionHisto.h"
00037
00038 #include <iterator>
00039
00040 using namespace edm;
00041 using namespace std;
00042
00043 DTSegmentAnalysisTask::DTSegmentAnalysisTask(const edm::ParameterSet& pset) : nevents(0) , nEventsInLS(0), hNevtPerLS(0) {
00044
00045 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Constructor called!";
00046
00047
00048 detailedAnalysis = pset.getUntrackedParameter<bool>("detailedAnalysis",false);
00049
00050 theRecHits4DLabel = pset.getParameter<string>("recHits4DLabel");
00051
00052 checkNoisyChannels = pset.getUntrackedParameter<bool>("checkNoisyChannels",false);
00053
00054 nTimeBins = pset.getUntrackedParameter<int>("nTimeBins",100);
00055
00056 nLSTimeBin = pset.getUntrackedParameter<int>("nLSTimeBin",2);
00057
00058 slideTimeBins = pset.getUntrackedParameter<bool>("slideTimeBins",true);
00059
00060
00061 theDbe = edm::Service<DQMStore>().operator->();
00062
00063
00064 topHistoFolder = pset.getUntrackedParameter<string>("topHistoFolder","DT/02-Segments");
00065
00066 hltDQMMode = pset.getUntrackedParameter<bool>("hltDQMMode",false);
00067
00068 }
00069
00070
00071 DTSegmentAnalysisTask::~DTSegmentAnalysisTask(){
00072 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Destructor called!";
00073 }
00074
00075
00076 void DTSegmentAnalysisTask::beginRun(const Run& run, const edm::EventSetup& context){
00077
00078
00079 context.get<MuonGeometryRecord>().get(dtGeom);
00080
00081 theDbe->setCurrentFolder("DT/EventInfo/Counters");
00082 nEventMonitor = theDbe->bookFloat("nProcessedEventsSegment");
00083
00084
00085 vector<DTChamber*> chambers = dtGeom->chambers();
00086 vector<DTChamber*>::const_iterator ch_it = chambers.begin();
00087 vector<DTChamber*>::const_iterator ch_end = chambers.end();
00088 for (; ch_it != ch_end; ++ch_it) {
00089 bookHistos((*ch_it)->id());
00090 }
00091
00092
00093 int modeTimeHisto = 0;
00094 if(!slideTimeBins) modeTimeHisto = 1;
00095 for(int wheel = -2; wheel != 3; ++wheel) {
00096 for(int sector = 1; sector <= 12; ++sector) {
00097 stringstream wheelstr; wheelstr << wheel;
00098 stringstream sectorstr; sectorstr << sector;
00099 string sectorHistoName = "NSegmPerEvent_W" + wheelstr.str()
00100 + "_Sec" + sectorstr.str();
00101 string sectorHistoTitle = "# segm. W" + wheelstr.str() + " Sect." + sectorstr.str();
00102
00103 theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheelstr.str() +
00104 "/Sector" + sectorstr.str());
00105 histoTimeEvol[wheel][sector] = new DTTimeEvolutionHisto(&(*theDbe),sectorHistoName,sectorHistoTitle,
00106 nTimeBins,nLSTimeBin,slideTimeBins,modeTimeHisto);
00107 }
00108 }
00109
00110 if(hltDQMMode) theDbe->setCurrentFolder(topHistoFolder);
00111 else theDbe->setCurrentFolder("DT/EventInfo/");
00112
00113 hNevtPerLS = new DTTimeEvolutionHisto(&(*theDbe),"NevtPerLS","# evt.",nTimeBins,nLSTimeBin,slideTimeBins,2);
00114
00115 }
00116
00117
00118 void DTSegmentAnalysisTask::endJob() {
00119
00120 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") <<"[DTSegmentAnalysisTask] endjob called!";
00121
00122 delete hNevtPerLS;
00123
00124
00125 }
00126
00127
00128
00129 void DTSegmentAnalysisTask::analyze(const edm::Event& event, const edm::EventSetup& setup) {
00130
00131 nevents++;
00132 nEventMonitor->Fill(nevents);
00133
00134 nEventsInLS++;
00135 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Analyze #Run: " << event.id().run()
00136 << " #Event: " << event.id().event();
00137 if(!(event.id().event()%1000))
00138 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Analyze #Run: " << event.id().run()
00139 << " #Event: " << event.id().event();
00140
00141 ESHandle<DTStatusFlag> statusMap;
00142 if(checkNoisyChannels) {
00143 setup.get<DTStatusFlagRcd>().get(statusMap);
00144 }
00145
00146
00147
00148
00149
00150 edm::Handle<DTRecSegment4DCollection> all4DSegments;
00151 event.getByLabel(theRecHits4DLabel, all4DSegments);
00152
00153 if(!all4DSegments.isValid()) return;
00154
00155
00156 DTRecSegment4DCollection::id_iterator chamberId;
00157 for (chamberId = all4DSegments->id_begin();
00158 chamberId != all4DSegments->id_end();
00159 ++chamberId){
00160
00161 DTRecSegment4DCollection::range range = all4DSegments->get(*chamberId);
00162
00163 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << " Chamber: " << *chamberId << " has " << distance(range.first, range.second)
00164 << " 4D segments";
00165
00166
00167 for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
00168 segment4D!=range.second;
00169 ++segment4D){
00170
00171
00172 bool segmNoisy = false;
00173 if(checkNoisyChannels) {
00174
00175 if((*segment4D).hasPhi()){
00176 const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
00177 vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
00178 map<DTSuperLayerId,vector<DTRecHit1D> > hitsBySLMap;
00179 for(vector<DTRecHit1D>::const_iterator hit = phiHits.begin();
00180 hit != phiHits.end(); ++hit) {
00181 DTWireId wireId = (*hit).wireId();
00182
00183
00184 bool isNoisy = false;
00185 bool isFEMasked = false;
00186 bool isTDCMasked = false;
00187 bool isTrigMask = false;
00188 bool isDead = false;
00189 bool isNohv = false;
00190 statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
00191 if(isNoisy) {
00192 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "Wire: " << wireId << " is noisy, skipping!";
00193 segmNoisy = true;
00194 }
00195 }
00196 }
00197
00198 if((*segment4D).hasZed()) {
00199 const DTSLRecSegment2D* zSeg = (*segment4D).zSegment();
00200
00201 vector<DTRecHit1D> zHits = zSeg->specificRecHits();
00202 for(vector<DTRecHit1D>::const_iterator hit = zHits.begin();
00203 hit != zHits.end(); ++hit) {
00204 DTWireId wireId = (*hit).wireId();
00205 bool isNoisy = false;
00206 bool isFEMasked = false;
00207 bool isTDCMasked = false;
00208 bool isTrigMask = false;
00209 bool isDead = false;
00210 bool isNohv = false;
00211 statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
00212 if(isNoisy) {
00213 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "Wire: " << wireId << " is noisy, skipping!";
00214 segmNoisy = true;
00215 }
00216 }
00217 }
00218
00219 }
00220 if (segmNoisy) {
00221 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")<<"skipping the segment: it contains noisy cells";
00222 continue;
00223 }
00224
00225
00226 int nHits=0;
00227 LocalPoint segment4DLocalPos = (*segment4D).localPosition();
00228 LocalVector segment4DLocalDirection = (*segment4D).localDirection();
00229 if((*segment4D).hasPhi())
00230 nHits = (((*segment4D).phiSegment())->specificRecHits()).size();
00231 if((*segment4D).hasZed())
00232 nHits = nHits + ((((*segment4D).zSegment())->specificRecHits()).size());
00233
00234 fillHistos(*chamberId,
00235 nHits,
00236 (*segment4D).chi2()/(*segment4D).degreesOfFreedom());
00237 }
00238 }
00239
00240
00241 }
00242
00243
00244
00245 void DTSegmentAnalysisTask::bookHistos(DTChamberId chamberId) {
00246
00247 edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << " Booking histos for chamber: " << chamberId;
00248
00249
00250
00251 stringstream wheel; wheel << chamberId.wheel();
00252 stringstream station; station << chamberId.station();
00253 stringstream sector; sector << chamberId.sector();
00254
00255 string chamberHistoName =
00256 "_W" + wheel.str() +
00257 "_St" + station.str() +
00258 "_Sec" + sector.str();
00259
00260
00261 for(int wh=-2; wh<=2; wh++){
00262 stringstream wheel; wheel << wh;
00263 theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str());
00264 string histoName = "numberOfSegments_W" + wheel.str();
00265 summaryHistos[wh] = theDbe->book2D(histoName.c_str(),histoName.c_str(),12,1,13,4,1,5);
00266 summaryHistos[wh]->setAxisTitle("Sector",1);
00267 summaryHistos[wh]->setBinLabel(1,"1",1);
00268 summaryHistos[wh]->setBinLabel(2,"2",1);
00269 summaryHistos[wh]->setBinLabel(3,"3",1);
00270 summaryHistos[wh]->setBinLabel(4,"4",1);
00271 summaryHistos[wh]->setBinLabel(5,"5",1);
00272 summaryHistos[wh]->setBinLabel(6,"6",1);
00273 summaryHistos[wh]->setBinLabel(7,"7",1);
00274 summaryHistos[wh]->setBinLabel(8,"8",1);
00275 summaryHistos[wh]->setBinLabel(9,"9",1);
00276 summaryHistos[wh]->setBinLabel(10,"10",1);
00277 summaryHistos[wh]->setBinLabel(11,"11",1);
00278 summaryHistos[wh]->setBinLabel(12,"12",1);
00279 summaryHistos[wh]->setBinLabel(1,"MB1",2);
00280 summaryHistos[wh]->setBinLabel(2,"MB2",2);
00281 summaryHistos[wh]->setBinLabel(3,"MB3",2);
00282 summaryHistos[wh]->setBinLabel(4,"MB4",2);
00283 }
00284
00285
00286 theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str() +
00287 "/Sector" + sector.str() +
00288 "/Station" + station.str());
00289
00290
00291 vector<MonitorElement *> histos;
00292 histos.push_back(theDbe->book1D("h4DSegmNHits"+chamberHistoName,
00293 "# of hits per segment",
00294 16, 0.5, 16.5));
00295 if(detailedAnalysis){
00296 histos.push_back(theDbe->book1D("h4DChi2"+chamberHistoName,
00297 "4D Segment reduced Chi2",
00298 20, 0, 20));
00299 }
00300 histosPerCh[chamberId] = histos;
00301 }
00302
00303
00304
00305 void DTSegmentAnalysisTask::fillHistos(DTChamberId chamberId,
00306 int nHits,
00307 float chi2) {
00308 int sector = chamberId.sector();
00309 if(chamberId.sector()==13) {
00310 sector = 4;
00311 } else if(chamberId.sector()==14) {
00312 sector = 10;
00313 }
00314
00315 summaryHistos[chamberId.wheel()]->Fill(sector,chamberId.station());
00316 histoTimeEvol[chamberId.wheel()][sector]->accumulateValueTimeSlot(1);
00317
00318 vector<MonitorElement *> histos = histosPerCh[chamberId];
00319 histos[0]->Fill(nHits);
00320 if(detailedAnalysis){
00321 histos[1]->Fill(chi2);
00322 }
00323
00324 }
00325
00326
00327 void DTSegmentAnalysisTask::endLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& eSetup) {
00328
00329 hNevtPerLS->updateTimeSlot(lumiSeg.luminosityBlock(), nEventsInLS);
00330
00331 for(int wheel = -2; wheel != 3; ++wheel) {
00332 for(int sector = 1; sector <= 12; ++sector) {
00333 histoTimeEvol[wheel][sector]->updateTimeSlot(lumiSeg.luminosityBlock(), nEventsInLS);
00334 }
00335 }
00336 }
00337
00338
00339 void DTSegmentAnalysisTask::beginLuminosityBlock(LuminosityBlock const& lumiSeg, EventSetup const& eSetup) {
00340 nEventsInLS = 0;
00341 }
00342
00343