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 Attributes
DTSegmentsTask Class Reference

#include <DTSegmentsTask.h>

Inheritance diagram for DTSegmentsTask:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &setup)
 
void beginJob (void)
 book the histos More...
 
void beginRun (const edm::Run &, const edm::EventSetup &)
 
 DTSegmentsTask (const edm::ParameterSet &pset)
 Constructor. More...
 
void endJob ()
 Endjob. More...
 
void endRun (const edm::Run &, const edm::EventSetup &)
 
virtual ~DTSegmentsTask ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesDependentUpon (const std::string &iProcessName, std::vector< const char * > &oModuleLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

Private Attributes

bool debug
 
edm::ParameterSet parameters
 
std::vector< MonitorElement * > phiHistos
 
DQMStoretheDbe
 
edm::EDGetTokenT
< DTRecSegment4DCollection
theRecHits4DLabel_
 
std::vector< MonitorElement * > thetaHistos
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

DQM Analysis of 4D DT segments

Author
G. Mila - INFN Torino

Definition at line 25 of file DTSegmentsTask.h.

Constructor & Destructor Documentation

DTSegmentsTask::DTSegmentsTask ( const edm::ParameterSet pset)

Constructor.

Definition at line 28 of file DTSegmentsTask.cc.

References debug, edm::ParameterSet::getUntrackedParameter(), cppFunctionSkipper::operator, Parameters::parameters, and AlCaHLTBitMon_QueryRunRegistry::string.

28  {
29 
30  debug = pset.getUntrackedParameter<bool>("debug",false);
31 
32  // Get the DQM needed services
34 
35  parameters = pset;
36 
37  // the name of the 4D rec hits collection
38  theRecHits4DLabel_ = consumes<DTRecSegment4DCollection>(parameters.getParameter<std::string>("recHits4DLabel"));
39 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< DTRecSegment4DCollection > theRecHits4DLabel_
edm::ParameterSet parameters
DQMStore * theDbe
DTSegmentsTask::~DTSegmentsTask ( )
virtual

Destructor.

Definition at line 42 of file DTSegmentsTask.cc.

42  {
43 }

Member Function Documentation

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

Implements edm::EDAnalyzer.

Definition at line 130 of file DTSegmentsTask.cc.

References gather_cfg::cout, debug, edm::EventSetup::get(), Parameters::parameters, findQualityFiles::size, DTRecSegment2D::specificRecHits(), and crabStatusFromReport::statusMap.

130  {
131 
132  // Get the map of noisy channels
133  bool checkNoisyChannels = parameters.getUntrackedParameter<bool>("checkNoisyChannels",false);
135  if(checkNoisyChannels) {
136  setup.get<DTStatusFlagRcd>().get(statusMap);
137  }
138 
139  // Get the 4D segment collection from the event
141  event.getByToken(theRecHits4DLabel_, all4DSegments);
142 
143  // Loop over all chambers containing a segment
144  DTRecSegment4DCollection::id_iterator chamberId;
145  for (chamberId = all4DSegments->id_begin(); chamberId != all4DSegments->id_end(); ++chamberId){
146  // Get the range for the corresponding ChamerId
147  DTRecSegment4DCollection::range range = all4DSegments->get(*chamberId);
148  int nsegm = distance(range.first, range.second);
149  if(debug)
150  cout << " Chamber: " << *chamberId << " has " << nsegm
151  << " 4D segments" << endl;
152 
153  // Loop over the rechits of this ChamerId
154  for (DTRecSegment4DCollection::const_iterator segment4D = range.first; segment4D!=range.second; ++segment4D){
155 
156  //FOR NOISY CHANNELS////////////////////////////////
157  bool segmNoisy = false;
158  if((*segment4D).hasPhi()){
159  const DTChamberRecSegment2D* phiSeg = (*segment4D).phiSegment();
160  vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
161  map<DTSuperLayerId,vector<DTRecHit1D> > hitsBySLMap;
162  for(vector<DTRecHit1D>::const_iterator hit = phiHits.begin(); hit != phiHits.end(); ++hit) {
163  DTWireId wireId = (*hit).wireId();
164 
165  // Check for noisy channels to skip them
166  if(checkNoisyChannels) {
167  bool isNoisy = false;
168  bool isFEMasked = false;
169  bool isTDCMasked = false;
170  bool isTrigMask = false;
171  bool isDead = false;
172  bool isNohv = false;
173  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
174  if(isNoisy) {
175  if(debug)
176  cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
177  segmNoisy = true;
178  }
179  }
180  }
181  }
182 
183  if((*segment4D).hasZed()) {
184  const DTSLRecSegment2D* zSeg = (*segment4D).zSegment(); // zSeg lives in the SL RF
185  // Check for noisy channels to skip them
186  vector<DTRecHit1D> zHits = zSeg->specificRecHits();
187  for(vector<DTRecHit1D>::const_iterator hit = zHits.begin();
188  hit != zHits.end(); ++hit) {
189  DTWireId wireId = (*hit).wireId();
190  if(checkNoisyChannels) {
191  bool isNoisy = false;
192  bool isFEMasked = false;
193  bool isTDCMasked = false;
194  bool isTrigMask = false;
195  bool isDead = false;
196  bool isNohv = false;
197  //cout<<"wire id "<<wireId<<endl;
198  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
199  if(isNoisy) {
200  if(debug)
201  cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
202  segmNoisy = true;
203  }
204  }
205  }
206  }
207 
208  if (segmNoisy) {
209  if(debug)
210  cout<<"skipping the segment: it contains noisy cells"<<endl;
211  continue;
212  }
213  //END FOR NOISY CHANNELS////////////////////////////////
214 
215  // Fill the histos
216  int nHits=0;
217  if((*segment4D).hasPhi()){
218  nHits = (((*segment4D).phiSegment())->specificRecHits()).size();
219  if(debug)
220  cout<<"Phi segment with number of hits: "<<nHits<<endl;
221  phiHistos[0]->Fill((*chamberId).wheel(), nHits);
222  phiHistos[1]->Fill((*chamberId).sector(), nHits);
223  phiHistos[2]->Fill((*chamberId).station(), nHits);
224  }
225  if((*segment4D).hasZed()) {
226  nHits = (((*segment4D).zSegment())->specificRecHits()).size();
227  if(debug)
228  cout<<"Zed segment with number of hits: "<<nHits<<endl;
229  thetaHistos[0]->Fill((*chamberId).wheel(), nHits);
230  thetaHistos[1]->Fill((*chamberId).sector(), nHits);
231  thetaHistos[2]->Fill((*chamberId).station(), nHits);
232  }
233 
234  } //loop over segments
235  } // loop over chambers
236 }
T getUntrackedParameter(std::string const &, T const &) const
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
edm::EDGetTokenT< DTRecSegment4DCollection > theRecHits4DLabel_
std::vector< MonitorElement * > thetaHistos
std::vector< MonitorElement * > phiHistos
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
const T & get() const
Definition: EventSetup.h:55
edm::ParameterSet parameters
tuple cout
Definition: gather_cfg.py:121
tuple size
Write out results.
void DTSegmentsTask::beginJob ( void  )
virtual

book the histos

Reimplemented from edm::EDAnalyzer.

Definition at line 46 of file DTSegmentsTask.cc.

46  {
47 }
void DTSegmentsTask::beginRun ( const edm::Run irun,
const edm::EventSetup isetup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 49 of file DTSegmentsTask.cc.

49  {
50 
51  theDbe->setVerbose(1);
52  theDbe->cd();
53  theDbe->setCurrentFolder("Muons/DTSegmentsMonitor");
54 
55  // histos for phi segments
56  phiHistos.push_back(theDbe->book2D("phiSegments_numHitsVsWheel", "phiSegments_numHitsVsWheel", 5, -2.5, 2.5, 20, 0, 20));
57  phiHistos[0]->setBinLabel(1,"W-2",1);
58  phiHistos[0]->setBinLabel(2,"W-1",1);
59  phiHistos[0]->setBinLabel(3,"W0",1);
60  phiHistos[0]->setBinLabel(4,"W1",1);
61  phiHistos[0]->setBinLabel(5,"W2",1);
62  phiHistos.push_back(theDbe->book2D("phiSegments_numHitsVsSector", "phiSegments_numHitsVsSector", 14, 0.5, 14.5, 20, 0, 20));
63  phiHistos[1]->setBinLabel(1,"Sec1",1);
64  phiHistos[1]->setBinLabel(2,"Sec2",1);
65  phiHistos[1]->setBinLabel(3,"Sec3",1);
66  phiHistos[1]->setBinLabel(4,"Sec4",1);
67  phiHistos[1]->setBinLabel(5,"Sec5",1);
68  phiHistos[1]->setBinLabel(6,"Sec6",1);
69  phiHistos[1]->setBinLabel(7,"Sec7",1);
70  phiHistos[1]->setBinLabel(8,"Sec8",1);
71  phiHistos[1]->setBinLabel(9,"Sec9",1);
72  phiHistos[1]->setBinLabel(10,"Sec10",1);
73  phiHistos[1]->setBinLabel(11,"Sec11",1);
74  phiHistos[1]->setBinLabel(12,"Sec12",1);
75  phiHistos[1]->setBinLabel(13,"Sec13",1);
76  phiHistos[1]->setBinLabel(14,"Sec14",1);
77  phiHistos.push_back(theDbe->book2D("phiSegments_numHitsVsStation", "phiSegments_numHitsVsStation", 4, 0.5, 4.5, 20, 0, 20));
78  phiHistos[2]->setBinLabel(1,"St1",1);
79  phiHistos[2]->setBinLabel(2,"St2",1);
80  phiHistos[2]->setBinLabel(3,"St3",1);
81  phiHistos[2]->setBinLabel(4,"St4",1);
82 
83  // histos for theta segments
84  thetaHistos.push_back(theDbe->book2D("thetaSegments_numHitsVsWheel", "thetaSegments_numHitsVsWheel", 5, -2.5, 2.5, 20, 0, 20));
85  thetaHistos[0]->setBinLabel(1,"W-2",1);
86  thetaHistos[0]->setBinLabel(2,"W-1",1);
87  thetaHistos[0]->setBinLabel(3,"W0",1);
88  thetaHistos[0]->setBinLabel(4,"W1",1);
89  thetaHistos[0]->setBinLabel(5,"W2",1);
90  thetaHistos.push_back(theDbe->book2D("thetaSegments_numHitsVsSector", "thetaSegments_numHitsVsSector", 14, 0.5, 14.5, 20, 0, 20));
91  thetaHistos[1]->setBinLabel(1,"Sec1",1);
92  thetaHistos[1]->setBinLabel(2,"Sec2",1);
93  thetaHistos[1]->setBinLabel(3,"Sec3",1);
94  thetaHistos[1]->setBinLabel(4,"Sec4",1);
95  thetaHistos[1]->setBinLabel(5,"Sec5",1);
96  thetaHistos[1]->setBinLabel(6,"Sec6",1);
97  thetaHistos[1]->setBinLabel(7,"Sec7",1);
98  thetaHistos[1]->setBinLabel(8,"Sec8",1);
99  thetaHistos[1]->setBinLabel(9,"Sec9",1);
100  thetaHistos[1]->setBinLabel(10,"Sec10",1);
101  thetaHistos[1]->setBinLabel(11,"Sec11",1);
102  thetaHistos[1]->setBinLabel(12,"Sec12",1);
103  thetaHistos[1]->setBinLabel(13,"Sec13",1);
104  thetaHistos[1]->setBinLabel(14,"Sec14",1);
105  thetaHistos.push_back(theDbe->book2D("thetaSegments_numHitsVsStation", "thetaSegments_numHitsVsStation", 4, 0.5, 4.5, 20, 0, 20));
106  thetaHistos[2]->setBinLabel(1,"St1",1);
107  thetaHistos[2]->setBinLabel(2,"St2",1);
108  thetaHistos[2]->setBinLabel(3,"St3",1);
109  thetaHistos[2]->setBinLabel(4,"St4",1);
110 
111 }
void cd(void)
go to top directory (ie. root)
Definition: DQMStore.cc:644
std::vector< MonitorElement * > thetaHistos
void setVerbose(unsigned level)
Definition: DQMStore.cc:631
std::vector< MonitorElement * > phiHistos
DQMStore * theDbe
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:1082
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:667
void DTSegmentsTask::endJob ( void  )
virtual

Endjob.

Reimplemented from edm::EDAnalyzer.

Definition at line 120 of file DTSegmentsTask.cc.

References dumpDBToFile_GT_ttrig_cfg::outputFileName, Parameters::parameters, and AlCaHLTBitMon_QueryRunRegistry::string.

120  {
121  bool outputMEsInRootFile = parameters.getParameter<bool>("OutputMEsInRootFile");
123  if(outputMEsInRootFile){
124  theDbe->save(outputFileName);
125  }
126 
127  theDbe->rmdir("DT/DTSegmentsTask");
128 }
T getParameter(std::string const &) const
void rmdir(const std::string &fullpath)
Definition: DQMStore.cc:3101
edm::ParameterSet parameters
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, const uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE", const bool resetMEsAfterWriting=false)
Definition: DQMStore.cc:2540
DQMStore * theDbe
void DTSegmentsTask::endRun ( const edm::Run irun,
const edm::EventSetup isetup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 113 of file DTSegmentsTask.cc.

References dumpDBToFile_GT_ttrig_cfg::outputFileName, Parameters::parameters, and AlCaHLTBitMon_QueryRunRegistry::string.

113  {
114  bool outputMEsInRootFile = parameters.getParameter<bool>("OutputMEsInRootFile");
116  if(outputMEsInRootFile){
117  theDbe->save(outputFileName);
118  }
119 }
T getParameter(std::string const &) const
edm::ParameterSet parameters
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, const uint32_t lumi=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE", const bool resetMEsAfterWriting=false)
Definition: DQMStore.cc:2540
DQMStore * theDbe

Member Data Documentation

bool DTSegmentsTask::debug
private
edm::ParameterSet DTSegmentsTask::parameters
private
std::vector<MonitorElement*> DTSegmentsTask::phiHistos
private

Definition at line 56 of file DTSegmentsTask.h.

DQMStore* DTSegmentsTask::theDbe
private

Definition at line 49 of file DTSegmentsTask.h.

edm::EDGetTokenT<DTRecSegment4DCollection> DTSegmentsTask::theRecHits4DLabel_
private

Definition at line 60 of file DTSegmentsTask.h.

std::vector<MonitorElement*> DTSegmentsTask::thetaHistos
private

Definition at line 57 of file DTSegmentsTask.h.