CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
DTSegmentAnalysisTask Class Reference

#include <DTSegmentAnalysisTask.h>

Inheritance diagram for DTSegmentAnalysisTask:
edm::EDAnalyzer

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &setup)
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &eSetup)
 Summary. More...
 
void beginRun (const edm::Run &, const edm::EventSetup &)
 BeginRun. More...
 
 DTSegmentAnalysisTask (const edm::ParameterSet &pset)
 Constructor. More...
 
void endJob ()
 Endjob. More...
 
void endLuminosityBlock (edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &eSetup)
 
virtual ~DTSegmentAnalysisTask ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Member Functions

void bookHistos (DTChamberId chamberId)
 
void fillHistos (DTChamberId chamberId, int nHits, float chi2)
 

Private Attributes

bool checkNoisyChannels
 
bool detailedAnalysis
 
edm::ESHandle< DTGeometrydtGeom
 
std::map< DTChamberId,
std::vector< MonitorElement * > > 
histosPerCh
 
std::map< int, std::map< int,
DTTimeEvolutionHisto * > > 
histoTimeEvol
 
bool hltDQMMode
 
DTTimeEvolutionHistohNevtPerLS
 
int nEventsInLS
 
int nLSTimeBin
 
int nTimeBins
 
edm::ParameterSet parameters
 
bool slideTimeBins
 
std::map< int, MonitorElement * > summaryHistos
 
DQMStoretheDbe
 
std::string theRecHits4DLabel
 
std::string topHistoFolder
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

DQM Analysis of 4D DT segments, it produces plots about:

Date:
2009/07/16 08:38:39
Revision:
1.10
Author
G. Cerminara - INFN Torino

Definition at line 36 of file DTSegmentAnalysisTask.h.

Constructor & Destructor Documentation

DTSegmentAnalysisTask::DTSegmentAnalysisTask ( const edm::ParameterSet pset)

Constructor.

Definition at line 43 of file DTSegmentAnalysisTask.cc.

References checkNoisyChannels, detailedAnalysis, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), hltDQMMode, nLSTimeBin, nTimeBins, cmsCodeRules.cppFunctionSkipper::operator, slideTimeBins, theDbe, theRecHits4DLabel, and topHistoFolder.

43  : nEventsInLS(0), hNevtPerLS(0) {
44 
45  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Constructor called!";
46 
47  // switch for detailed analysis
48  detailedAnalysis = pset.getUntrackedParameter<bool>("detailedAnalysis","false");
49  // the name of the 4D rec hits collection
50  theRecHits4DLabel = pset.getParameter<string>("recHits4DLabel");
51  // Get the map of noisy channels
52  checkNoisyChannels = pset.getUntrackedParameter<bool>("checkNoisyChannels","false");
53  // # of bins in the time histos
54  nTimeBins = pset.getUntrackedParameter<int>("nTimeBins",100);
55  // # of LS per bin in the time histos
56  nLSTimeBin = pset.getUntrackedParameter<int>("nLSTimeBin",2);
57  // switch on/off sliding bins in time histos
58  slideTimeBins = pset.getUntrackedParameter<bool>("slideTimeBins",true);
59 
60  // Get the DQM needed services
62 
63  // top folder for the histograms in DQMStore
64  topHistoFolder = pset.getUntrackedParameter<string>("topHistoFolder","DT/02-Segments");
65  // hlt DQM mode
66  hltDQMMode = pset.getUntrackedParameter<bool>("hltDQMMode",false);
67 
68  }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTTimeEvolutionHisto * hNevtPerLS
DTSegmentAnalysisTask::~DTSegmentAnalysisTask ( )
virtual

Destructor.

Definition at line 71 of file DTSegmentAnalysisTask.cc.

71  {
72  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Destructor called!";
73 }

Member Function Documentation

void DTSegmentAnalysisTask::analyze ( const edm::Event event,
const edm::EventSetup setup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 126 of file DTSegmentAnalysisTask.cc.

References checkNoisyChannels, edm::EventID::event(), fillHistos(), edm::EventSetup::get(), edm::EventBase::id(), edm::HandleBase::isValid(), nEventsInLS, findQualityFiles::size, DTRecSegment2D::specificRecHits(), crabStatusFromReport::statusMap, and theRecHits4DLabel.

126  {
127  nEventsInLS++;
128  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Analyze #Run: " << event.id().run()
129  << " #Event: " << event.id().event();
130  if(!(event.id().event()%1000))
131  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "[DTSegmentAnalysisTask] Analyze #Run: " << event.id().run()
132  << " #Event: " << event.id().event();
133 
135  if(checkNoisyChannels) {
136  setup.get<DTStatusFlagRcd>().get(statusMap);
137  }
138 
139 
140  // -- 4D segment analysis -----------------------------------------------------
141 
142  // Get the 4D segment collection from the event
144  event.getByLabel(theRecHits4DLabel, all4DSegments);
145 
146  if(!all4DSegments.isValid()) return;
147 
148  // Loop over all chambers containing a segment
150  for (chamberId = all4DSegments->id_begin();
151  chamberId != all4DSegments->id_end();
152  ++chamberId){
153  // Get the range for the corresponding ChamerId
154  DTRecSegment4DCollection::range range = all4DSegments->get(*chamberId);
155 
156  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << " Chamber: " << *chamberId << " has " << distance(range.first, range.second)
157  << " 4D segments";
158 
159  // Loop over the rechits of this ChamerId
160  for (DTRecSegment4DCollection::const_iterator segment4D = range.first;
161  segment4D!=range.second;
162  ++segment4D){
163 
164  //FOR NOISY CHANNELS////////////////////////////////
165  bool segmNoisy = false;
166  if(checkNoisyChannels) {
167 
168  if((*segment4D).hasPhi()){
169  const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
170  vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
171  map<DTSuperLayerId,vector<DTRecHit1D> > hitsBySLMap;
172  for(vector<DTRecHit1D>::const_iterator hit = phiHits.begin();
173  hit != phiHits.end(); ++hit) {
174  DTWireId wireId = (*hit).wireId();
175 
176  // Check for noisy channels to skip them
177  bool isNoisy = false;
178  bool isFEMasked = false;
179  bool isTDCMasked = false;
180  bool isTrigMask = false;
181  bool isDead = false;
182  bool isNohv = false;
183  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
184  if(isNoisy) {
185  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "Wire: " << wireId << " is noisy, skipping!";
186  segmNoisy = true;
187  }
188  }
189  }
190 
191  if((*segment4D).hasZed()) {
192  const DTSLRecSegment2D* zSeg = (*segment4D).zSegment(); // zSeg lives in the SL RF
193  // Check for noisy channels to skip them
194  vector<DTRecHit1D> zHits = zSeg->specificRecHits();
195  for(vector<DTRecHit1D>::const_iterator hit = zHits.begin();
196  hit != zHits.end(); ++hit) {
197  DTWireId wireId = (*hit).wireId();
198  bool isNoisy = false;
199  bool isFEMasked = false;
200  bool isTDCMasked = false;
201  bool isTrigMask = false;
202  bool isDead = false;
203  bool isNohv = false;
204  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
205  if(isNoisy) {
206  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << "Wire: " << wireId << " is noisy, skipping!";
207  segmNoisy = true;
208  }
209  }
210  }
211 
212  } // end of switch on noisy channels
213  if (segmNoisy) {
214  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask")<<"skipping the segment: it contains noisy cells";
215  continue;
216  }
217  //END FOR NOISY CHANNELS////////////////////////////////
218 
219  int nHits=0;
220  LocalPoint segment4DLocalPos = (*segment4D).localPosition();
221  LocalVector segment4DLocalDirection = (*segment4D).localDirection();
222  if((*segment4D).hasPhi())
223  nHits = (((*segment4D).phiSegment())->specificRecHits()).size();
224  if((*segment4D).hasZed())
225  nHits = nHits + ((((*segment4D).zSegment())->specificRecHits()).size());
226 
227  fillHistos(*chamberId,
228  nHits,
229  (*segment4D).chi2()/(*segment4D).degreesOfFreedom());
230  }
231  }
232 
233  // -----------------------------------------------------------------------------
234 }
EventNumber_t event() const
Definition: EventID.h:44
void fillHistos(DTChamberId chamberId, int nHits, float chi2)
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:53
identifier iterator
Definition: RangeMap.h:139
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:46
bool isValid() const
Definition: HandleBase.h:76
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
const T & get() const
Definition: EventSetup.h:55
edm::EventID id() const
Definition: EventBase.h:56
tuple size
Write out results.
void DTSegmentAnalysisTask::beginLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
edm::EventSetup const &  eSetup 
)
virtual

Summary.

Reimplemented from edm::EDAnalyzer.

Definition at line 332 of file DTSegmentAnalysisTask.cc.

References nEventsInLS.

332  {
333  nEventsInLS = 0;
334 }
void DTSegmentAnalysisTask::beginRun ( const edm::Run run,
const edm::EventSetup context 
)
virtual

BeginRun.

Reimplemented from edm::EDAnalyzer.

Definition at line 76 of file DTSegmentAnalysisTask.cc.

References bookHistos(), chambers, dtGeom, edm::EventSetup::get(), histoTimeEvol, hltDQMMode, hNevtPerLS, nLSTimeBin, nTimeBins, DQMStore::setCurrentFolder(), slideTimeBins, theDbe, and topHistoFolder.

76  {
77 
78  // Get the DT Geometry
79  context.get<MuonGeometryRecord>().get(dtGeom);
80 
81  // loop over all the DT chambers & book the histos
82  vector<DTChamber*> chambers = dtGeom->chambers();
83  vector<DTChamber*>::const_iterator ch_it = chambers.begin();
84  vector<DTChamber*>::const_iterator ch_end = chambers.end();
85  for (; ch_it != ch_end; ++ch_it) {
86  bookHistos((*ch_it)->id());
87  }
88 
89  // book sector time-evolution histos
90  int modeTimeHisto = 0;
91  if(!slideTimeBins) modeTimeHisto = 1;
92  for(int wheel = -2; wheel != 3; ++wheel) { // loop over wheels
93  for(int sector = 1; sector <= 12; ++sector) { // loop over sectors
94  stringstream wheelstr; wheelstr << wheel;
95  stringstream sectorstr; sectorstr << sector;
96  string sectorHistoName = "NSegmPerEvent_W" + wheelstr.str()
97  + "_Sec" + sectorstr.str();
98  string sectorHistoTitle = "# segm. W" + wheelstr.str() + " Sect." + sectorstr.str();
99 
100  theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheelstr.str() +
101  "/Sector" + sectorstr.str());
102  histoTimeEvol[wheel][sector] = new DTTimeEvolutionHisto(&(*theDbe),sectorHistoName,sectorHistoTitle,
103  nTimeBins,nLSTimeBin,slideTimeBins,modeTimeHisto);
104  }
105  }
106 
108  else theDbe->setCurrentFolder("DT/EventInfo/");
109 
110  hNevtPerLS = new DTTimeEvolutionHisto(&(*theDbe),"NevtPerLS","# evt.",nTimeBins,nLSTimeBin,slideTimeBins,2);
111 
112 }
DTTimeEvolutionHisto * hNevtPerLS
static char chambers[TOTALCHAMBERS][20]
Definition: ReadPGInfo.cc:243
edm::ESHandle< DTGeometry > dtGeom
const T & get() const
Definition: EventSetup.h:55
std::map< int, std::map< int, DTTimeEvolutionHisto * > > histoTimeEvol
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:232
void bookHistos(DTChamberId chamberId)
void DTSegmentAnalysisTask::bookHistos ( DTChamberId  chamberId)
private

Definition at line 238 of file DTSegmentAnalysisTask.cc.

References DQMStore::book1D(), DQMStore::book2D(), detailedAnalysis, mergeVDriftHistosByStation::histos, histosPerCh, DTChamberId::sector(), DQMStore::setCurrentFolder(), DTChamberId::station(), relativeConstraints::station, summaryHistos, theDbe, topHistoFolder, and DTChamberId::wheel().

Referenced by beginRun().

238  {
239 
240  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") << " Booking histos for chamber: " << chamberId;
241 
242 
243  // Compose the chamber name
244  stringstream wheel; wheel << chamberId.wheel();
245  stringstream station; station << chamberId.station();
246  stringstream sector; sector << chamberId.sector();
247 
248  string chamberHistoName =
249  "_W" + wheel.str() +
250  "_St" + station.str() +
251  "_Sec" + sector.str();
252 
253 
254  for(int wh=-2; wh<=2; wh++){
255  stringstream wheel; wheel << wh;
256  theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str());
257  string histoName = "numberOfSegments_W" + wheel.str();
258  summaryHistos[wh] = theDbe->book2D(histoName.c_str(),histoName.c_str(),12,1,13,4,1,5);
259  summaryHistos[wh]->setAxisTitle("Sector",1);
260  summaryHistos[wh]->setBinLabel(1,"1",1);
261  summaryHistos[wh]->setBinLabel(2,"2",1);
262  summaryHistos[wh]->setBinLabel(3,"3",1);
263  summaryHistos[wh]->setBinLabel(4,"4",1);
264  summaryHistos[wh]->setBinLabel(5,"5",1);
265  summaryHistos[wh]->setBinLabel(6,"6",1);
266  summaryHistos[wh]->setBinLabel(7,"7",1);
267  summaryHistos[wh]->setBinLabel(8,"8",1);
268  summaryHistos[wh]->setBinLabel(9,"9",1);
269  summaryHistos[wh]->setBinLabel(10,"10",1);
270  summaryHistos[wh]->setBinLabel(11,"11",1);
271  summaryHistos[wh]->setBinLabel(12,"12",1);
272  summaryHistos[wh]->setBinLabel(1,"MB1",2);
273  summaryHistos[wh]->setBinLabel(2,"MB2",2);
274  summaryHistos[wh]->setBinLabel(3,"MB3",2);
275  summaryHistos[wh]->setBinLabel(4,"MB4",2);
276  }
277 
278 
279  theDbe->setCurrentFolder(topHistoFolder + "/Wheel" + wheel.str() +
280  "/Sector" + sector.str() +
281  "/Station" + station.str());
282 
283  // Create the monitor elements
284  vector<MonitorElement *> histos;
285  histos.push_back(theDbe->book1D("h4DSegmNHits"+chamberHistoName,
286  "# of hits per segment",
287  16, 0.5, 16.5));
288  if(detailedAnalysis){
289  histos.push_back(theDbe->book1D("h4DChi2"+chamberHistoName,
290  "4D Segment reduced Chi2",
291  20, 0, 20));
292  }
293  histosPerCh[chamberId] = histos;
294 }
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:514
std::map< int, MonitorElement * > summaryHistos
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
int sector() const
Definition: DTChamberId.h:63
int station() const
Return the station number.
Definition: DTChamberId.h:53
MonitorElement * book2D(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Book 2D histogram.
Definition: DQMStore.cc:642
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:232
void DTSegmentAnalysisTask::endJob ( void  )
virtual

Endjob.

Reimplemented from edm::EDAnalyzer.

Definition at line 115 of file DTSegmentAnalysisTask.cc.

References hNevtPerLS.

115  {
116 
117  edm::LogVerbatim ("DTDQM|DTMonitorModule|DTSegmentAnalysisTask") <<"[DTSegmentAnalysisTask] endjob called!";
118 
119  delete hNevtPerLS;
120  //theDbe->save("testMonitoring.root");
121 
122 }
DTTimeEvolutionHisto * hNevtPerLS
void DTSegmentAnalysisTask::endLuminosityBlock ( edm::LuminosityBlock const &  lumiSeg,
edm::EventSetup const &  eSetup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 320 of file DTSegmentAnalysisTask.cc.

References histoTimeEvol, hNevtPerLS, edm::LuminosityBlockBase::luminosityBlock(), nEventsInLS, and DTTimeEvolutionHisto::updateTimeSlot().

320  {
321 
322  hNevtPerLS->updateTimeSlot(lumiSeg.luminosityBlock(), nEventsInLS);
323  // book sector time-evolution histos
324  for(int wheel = -2; wheel != 3; ++wheel) {
325  for(int sector = 1; sector <= 12; ++sector) {
326  histoTimeEvol[wheel][sector]->updateTimeSlot(lumiSeg.luminosityBlock(), nEventsInLS);
327  }
328  }
329 }
DTTimeEvolutionHisto * hNevtPerLS
void updateTimeSlot(int ls, int nEventsInLS)
std::map< int, std::map< int, DTTimeEvolutionHisto * > > histoTimeEvol
void DTSegmentAnalysisTask::fillHistos ( DTChamberId  chamberId,
int  nHits,
float  chi2 
)
private

Definition at line 298 of file DTSegmentAnalysisTask.cc.

References detailedAnalysis, mergeVDriftHistosByStation::histos, histosPerCh, histoTimeEvol, DTChamberId::sector(), DTChamberId::station(), summaryHistos, and DTChamberId::wheel().

Referenced by analyze().

300  {
301  int sector = chamberId.sector();
302  if(chamberId.sector()==13) {
303  sector = 4;
304  } else if(chamberId.sector()==14) {
305  sector = 10;
306  }
307 
308  summaryHistos[chamberId.wheel()]->Fill(sector,chamberId.station());
309  histoTimeEvol[chamberId.wheel()][sector]->accumulateValueTimeSlot(1);
310 
311  vector<MonitorElement *> histos = histosPerCh[chamberId];
312  histos[0]->Fill(nHits);
313  if(detailedAnalysis){
314  histos[1]->Fill(chi2);
315  }
316 
317 }
std::map< int, MonitorElement * > summaryHistos
std::map< DTChamberId, std::vector< MonitorElement * > > histosPerCh
std::map< int, std::map< int, DTTimeEvolutionHisto * > > histoTimeEvol
int sector() const
Definition: DTChamberId.h:63
int station() const
Return the station number.
Definition: DTChamberId.h:53
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47

Member Data Documentation

bool DTSegmentAnalysisTask::checkNoisyChannels
private

Definition at line 74 of file DTSegmentAnalysisTask.h.

Referenced by analyze(), and DTSegmentAnalysisTask().

bool DTSegmentAnalysisTask::detailedAnalysis
private

Definition at line 65 of file DTSegmentAnalysisTask.h.

Referenced by bookHistos(), DTSegmentAnalysisTask(), and fillHistos().

edm::ESHandle<DTGeometry> DTSegmentAnalysisTask::dtGeom
private

Definition at line 68 of file DTSegmentAnalysisTask.h.

Referenced by beginRun().

std::map<DTChamberId, std::vector<MonitorElement*> > DTSegmentAnalysisTask::histosPerCh
private

Definition at line 86 of file DTSegmentAnalysisTask.h.

Referenced by bookHistos(), and fillHistos().

std::map<int, std::map<int, DTTimeEvolutionHisto*> > DTSegmentAnalysisTask::histoTimeEvol
private

Definition at line 88 of file DTSegmentAnalysisTask.h.

Referenced by beginRun(), endLuminosityBlock(), and fillHistos().

bool DTSegmentAnalysisTask::hltDQMMode
private

Definition at line 102 of file DTSegmentAnalysisTask.h.

Referenced by beginRun(), and DTSegmentAnalysisTask().

DTTimeEvolutionHisto* DTSegmentAnalysisTask::hNevtPerLS
private

Definition at line 91 of file DTSegmentAnalysisTask.h.

Referenced by beginRun(), endJob(), and endLuminosityBlock().

int DTSegmentAnalysisTask::nEventsInLS
private

Definition at line 90 of file DTSegmentAnalysisTask.h.

Referenced by analyze(), beginLuminosityBlock(), and endLuminosityBlock().

int DTSegmentAnalysisTask::nLSTimeBin
private

Definition at line 96 of file DTSegmentAnalysisTask.h.

Referenced by beginRun(), and DTSegmentAnalysisTask().

int DTSegmentAnalysisTask::nTimeBins
private

Definition at line 94 of file DTSegmentAnalysisTask.h.

Referenced by beginRun(), and DTSegmentAnalysisTask().

edm::ParameterSet DTSegmentAnalysisTask::parameters
private
bool DTSegmentAnalysisTask::slideTimeBins
private

Definition at line 98 of file DTSegmentAnalysisTask.h.

Referenced by beginRun(), and DTSegmentAnalysisTask().

std::map< int, MonitorElement* > DTSegmentAnalysisTask::summaryHistos
private

Definition at line 87 of file DTSegmentAnalysisTask.h.

Referenced by bookHistos(), and fillHistos().

DQMStore* DTSegmentAnalysisTask::theDbe
private

Definition at line 62 of file DTSegmentAnalysisTask.h.

Referenced by beginRun(), bookHistos(), and DTSegmentAnalysisTask().

std::string DTSegmentAnalysisTask::theRecHits4DLabel
private

Definition at line 71 of file DTSegmentAnalysisTask.h.

Referenced by analyze(), and DTSegmentAnalysisTask().

std::string DTSegmentAnalysisTask::topHistoFolder
private

Definition at line 100 of file DTSegmentAnalysisTask.h.

Referenced by beginRun(), bookHistos(), and DTSegmentAnalysisTask().