CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTDigiForNoiseTask.cc
Go to the documentation of this file.
1  /*
2  * \file DTDigiForNoiseTask.cc
3  *
4  * \author G. Mila - INFN Torino
5  *
6  */
7 
9 
10 // Framework
12 
13 // Digis
16 
17 // Geometry
22 
25 
26 #include <stdio.h>
27 #include <sstream>
28 #include <math.h>
29 
30 using namespace edm;
31 using namespace std;
32 
33 
35 
36  debug = ps.getUntrackedParameter<bool>("debug", false);
37  dtDigisToken_ = consumes<DTDigiCollection>(
38  edm::InputTag(ps.getUntrackedParameter<std::string>("diDigisLabel", "dtunpacker")));
39 
40  if(debug)
41  cout<<"[DTDigiForNoiseTask]: Constructor"<<endl;
42 
43  parameters = ps;
44 
45  dbe = edm::Service<DQMStore>().operator->();
46 
47 }
48 
49 
51 
52  if(debug)
53  cout << "DTDigiForNoiseTask: analyzed " << nevents << " events" << endl;
54 
55 }
56 
57 
59 
60  if(debug)
61  cout<<"[DTDigiForNoiseTask] endjob called!"<<endl;
62 
63  dbe->rmdir("DT/DTDigiForNoiseTask");
64 
65 }
66 
67 
69 
70  if(debug)
71  cout<<"[DTDigiForNoiseTask]: BeginJob"<<endl;
72 
73  nevents = 0;
74 
75 }
76 
78 
79  if(debug)
80  cout<<"[DTDigiForNoiseTask]: BeginRun"<<endl;
81 
82  // Get the geometry
83  context.get<MuonGeometryRecord>().get(muonGeom);
84 
85 }
86 
88 
89  if(debug)
90  cout<<"[DTDigiForNoiseTask]: Begin of LS transition"<<endl;
91 
92  if(lumiSeg.id().luminosityBlock()%parameters.getUntrackedParameter<int>("ResetCycle", 3) == 0) {
93  for(map< DTLayerId, MonitorElement* > ::const_iterator histo = digiHistos.begin();
94  histo != digiHistos.end();
95  histo++) {
96  (*histo).second->Reset();
97  }
98  }
99 
100 }
101 
102 
104 
105  if (debug) cout<<"[DTDigiForNoiseTask]: booking"<<endl;
106 
107  const DTSuperLayerId dtSLId = lId.superlayerId();
108  const DTChamberId dtChId = dtSLId.chamberId();
109  stringstream layer; layer << lId.layer();
110  stringstream superLayer; superLayer << dtSLId.superlayer();
111  stringstream wheel; wheel << dtChId.wheel();
112  stringstream station; station << dtChId.station();
113  stringstream sector; sector << dtChId.sector();
114 
115  dbe->setCurrentFolder("DT/DTDigiForNoiseTask/Wheel" + wheel.str() +
116  "/Station" + station.str() +
117  "/Sector" + sector.str() + "/DigiPerEvent");
118 
119  if (debug){
120  cout<<"[DTDigiForNoiseTask]: folder "<< "DT/DTDigiTask/Wheel" + wheel.str() +
121  "/Station" + station.str() +
122  "/Sector" + sector.str() + "/DigiPerEvent"<<endl;
123  }
124 
125  string histoName =
126  "DigiPerEvent_W" + wheel.str()
127  + "_St" + station.str()
128  + "_Sec" + sector.str()
129  + "_SL" + superLayer.str()
130  + "_L" + layer.str();
131 
132  if (debug) cout<<"[DTDigiTask]: histoName "<<histoName<<endl;
133 
134  const DTTopology& dtTopo = muonGeom->layer(lId)->specificTopology();
135  const int firstWire = dtTopo.firstChannel();
136  const int lastWire = dtTopo.lastChannel();
137  int nWires = lastWire-firstWire+1;
138 
139  digiHistos[lId] = dbe->book2D(histoName,histoName,nWires,firstWire,lastWire,10,-0.5,9.5);
140 
141 }
142 
143 
145 
146  nevents++;
147  // cout << "events: " << nevents << endl;
148  if (nevents%1000 == 0 && debug) {}
149 
151  e.getByToken(dtDigisToken_, dtdigis);
152 
153  std::map< int,int > DigiPerWirePerEvent;
154 
155  // Loop over all the chambers
156  vector<DTChamber*>::const_iterator ch_it = muonGeom->chambers().begin();
157  vector<DTChamber*>::const_iterator ch_end = muonGeom->chambers().end();
158  // Loop over the SLs
159  for (; ch_it != ch_end; ++ch_it) {
160  // DTChamberId ch = (*ch_it)->id();
161  vector<const DTSuperLayer*>::const_iterator sl_it = (*ch_it)->superLayers().begin();
162  vector<const DTSuperLayer*>::const_iterator sl_end = (*ch_it)->superLayers().end();
163  // Loop over the SLs
164  for(; sl_it != sl_end; ++sl_it) {
165  vector<const DTLayer*>::const_iterator l_it = (*sl_it)->layers().begin();
166  vector<const DTLayer*>::const_iterator l_end = (*sl_it)->layers().end();
167  // Loop over the Ls
168  for(; l_it != l_end; ++l_it) {
169  DTLayerId layerId = (*l_it)->id();
170 
171  DTDigiCollection::Range layerDigi= dtdigis->get(layerId);
172  if(layerDigi.first != layerDigi.second){
173 
174  const DTTopology& dtTopo = muonGeom->layer(layerId)->specificTopology();
175  const int firstWire = dtTopo.firstChannel();
176  const int lastWire = dtTopo.lastChannel();
177 
178  if (digiHistos.find(layerId) == digiHistos.end())
179  bookHistos(layerId);
180 
181  if (digiHistos.find(layerId) != digiHistos.end()){
182  for (int wire=firstWire; wire-lastWire <= 0; wire++) {
183  DigiPerWirePerEvent[wire]= 0;
184  }
185 
186  for (DTDigiCollection::const_iterator digi = layerDigi.first;
187  digi!=layerDigi.second;
188  ++digi){
189  DigiPerWirePerEvent[(*digi).wire()]+=1;
190  }
191 
192  for (int wire=firstWire; wire-lastWire<=0; wire++) {
193  digiHistos.find(layerId)->second->Fill(wire,DigiPerWirePerEvent[wire]);
194  }
195  }
196  }
197 
198  } //Loop Ls
199  } //Loop SLs
200  } //Loop over chambers
201 
202 }
203 
204 // Local Variables:
205 // show-trailing-whitespace: t
206 // truncate-lines: t
207 // End:
LuminosityBlockID id() const
T getUntrackedParameter(std::string const &, T const &) const
dictionary parameters
Definition: Parameters.py:2
void beginJob()
BeginJob.
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
DTChamberId chamberId() const
Return the corresponding ChamberId.
void beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, edm::EventSetup const &context)
To reset the MEs.
int layer() const
Return the layer number.
Definition: DTLayerId.h:53
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:59
int firstChannel() const
Returns the wire number of the first wire.
Definition: DTTopology.h:78
void bookHistos()
Definition: Histogram.h:33
int lastChannel() const
Returns the wire number of the last wire.
Definition: DTTopology.h:80
void bookHistos(const DTLayerId &dtSL)
Book the ME.
int nevents
void analyze(const edm::Event &e, const edm::EventSetup &c)
Analyze.
virtual ~DTDigiForNoiseTask()
Destructor.
int superlayer() const
Return the superlayer number (deprecated method name)
#define debug
Definition: HDRShower.cc:19
const T & get() const
Definition: EventSetup.h:55
LuminosityBlockNumber_t luminosityBlock() const
std::vector< DTDigi >::const_iterator const_iterator
int sector() const
Definition: DTChamberId.h:61
void beginRun(const edm::Run &, const edm::EventSetup &)
BeginRun.
std::pair< const_iterator, const_iterator > Range
tuple cout
Definition: gather_cfg.py:121
DTDigiForNoiseTask(const edm::ParameterSet &ps)
Constructor.
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
Definition: Run.h:41