CMS 3D CMS Logo

MuonDTDigis.cc
Go to the documentation of this file.
1 
7 #include "MuonDTDigis.h"
8 
9 #include <iostream>
10 #include <string>
11 
12 #include "TFile.h"
13 
14 #include "SimMuon/DTDigitizer/test/Histograms.h"
15 
16 using namespace edm;
17 using namespace std;
18 
20 
21  // ----------------------
22  // Get the debug parameter for verbose output
23  verbose_ = pset.getUntrackedParameter<bool>("verbose",false);
24 
25  // the name of the Digi collection
26  SimHitToken_ = consumes<edm::PSimHitContainer>(pset.getParameter<edm::InputTag>("SimHitLabel"));
27 
28  // the name of the Digi collection
29  DigiToken_ = consumes<DTDigiCollection> (pset.getParameter<edm::InputTag>(("DigiLabel")));
30 
31  hDigis_global = std::make_unique<hDigis>("Global");
32  hDigis_W0 = std::make_unique<hDigis>("Wheel0");
33  hDigis_W1 = std::make_unique<hDigis>("Wheel1");
34  hDigis_W2 = std::make_unique<hDigis>("Wheel2");
35  hAllHits = std::make_unique<hHits>("AllHits");
36 }
37 
39  if(verbose_)
40  cout << "[MuonDTDigis] Destructor called" << endl;
41 }
42 
43 void MuonDTDigis::bookHistograms(DQMStore::IBooker & iBooker, edm::Run const & iRun, edm::EventSetup const & /* iSetup */)
44 {
45  meDigiTimeBox_ = nullptr;
46  meDigiTimeBox_wheel2m_ = nullptr;
47  meDigiTimeBox_wheel1m_ = nullptr;
48  meDigiTimeBox_wheel0_ = nullptr;
49  meDigiTimeBox_wheel1p_ = nullptr;
50  meDigiTimeBox_wheel2p_ = nullptr;
51  meDigiEfficiency_ = nullptr;
52  meDigiEfficiencyMu_ = nullptr;
53  meDoubleDigi_ = nullptr;
54  meSimvsDigi_ = nullptr;
55  meWire_DoubleDigi_ = nullptr;
56 
57  meMB1_sim_occup_ = nullptr;
58  meMB1_digi_occup_ = nullptr;
59  meMB2_sim_occup_ = nullptr;
60  meMB2_digi_occup_ = nullptr;
61  meMB3_sim_occup_ = nullptr;
62  meMB3_digi_occup_ = nullptr;
63  meMB4_sim_occup_ = nullptr;
64  meMB4_digi_occup_ = nullptr;
65 
66  meDigiHisto_ = nullptr;
67 
68  // ----------------------
69  // We go
70  // ----------------------
71 
72  Char_t histo_n[100];
73  Char_t histo_t[100];
74 
75  iBooker.setCurrentFolder("MuonDTDigisV/DTDigiValidationTask");
76 
77  sprintf (histo_n, "DigiTimeBox" );
78  sprintf (histo_t, "Digi Time Box" );
79  meDigiTimeBox_ = iBooker.book1D(histo_n, histo_t, 1536,0,1200);
80 
81  sprintf (histo_n, "DigiTimeBox_wheel2m" );
82  sprintf (histo_t, "Digi Time Box wheel -2" );
83  meDigiTimeBox_wheel2m_ = iBooker.book1D(histo_n, histo_t, 384,0,1200);
84 
85  sprintf (histo_n, "DigiTimeBox_wheel1m" );
86  sprintf (histo_t, "Digi Time Box wheel -1" );
87  meDigiTimeBox_wheel1m_ = iBooker.book1D(histo_n, histo_t, 384,0,1200);
88 
89  sprintf (histo_n, "DigiTimeBox_wheel0" );
90  sprintf (histo_t, "Digi Time Box wheel 0" );
91  meDigiTimeBox_wheel0_ = iBooker.book1D(histo_n, histo_t, 384,0,1200);
92 
93  sprintf (histo_n, "DigiTimeBox_wheel1p" );
94  sprintf (histo_t, "Digi Time Box wheel 1" );
95  meDigiTimeBox_wheel1p_ = iBooker.book1D(histo_n, histo_t, 384,0,1200);
96 
97  sprintf (histo_n, "DigiTimeBox_wheel2p" );
98  sprintf (histo_t, "Digi Time Box wheel 2" );
99  meDigiTimeBox_wheel2p_ = iBooker.book1D(histo_n, histo_t, 384,0,1200);
100 
101  sprintf (histo_n, "DigiEfficiencyMu" );
102  sprintf (histo_t, "Ratio (#Digis Mu)/(#SimHits Mu)" );
103  meDigiEfficiencyMu_ = iBooker.book1D(histo_n, histo_t, 100, 0., 5.);
104 
105  sprintf (histo_n, "DigiEfficiency" );
106  sprintf (histo_t, "Ratio (#Digis)/(#SimHits)" );
107  meDigiEfficiency_ = iBooker.book1D(histo_n, histo_t, 100, 0., 5.);
108 
109  sprintf (histo_n, "Number_Digi_per_layer" );
110  sprintf (histo_t, "Number_Digi_per_layer" );
111  meDoubleDigi_ = iBooker.book1D(histo_n, histo_t, 10,0.,10.);
112 
113  sprintf (histo_n, "Number_simhit_vs_digi" );
114  sprintf (histo_t, "Number_simhit_vs_digi" );
115  meSimvsDigi_ = iBooker.book2D(histo_n, histo_t, 100, 0., 140., 100, 0., 140.);
116 
117  sprintf (histo_n, "Wire_Number_with_double_Digi" );
118  sprintf (histo_t, "Wire_Number_with_double_Digi" );
119  meWire_DoubleDigi_ = iBooker.book1D(histo_n, histo_t, 100,0.,100.);
120 
121  sprintf (histo_n, "Simhit_occupancy_MB1" );
122  sprintf (histo_t, "Simhit_occupancy_MB1" );
123  meMB1_sim_occup_ = iBooker.book1D(histo_n, histo_t, 55, 0., 55. );
124 
125  sprintf (histo_n, "Digi_occupancy_MB1" );
126  sprintf (histo_t, "Digi_occupancy_MB1" );
127  meMB1_digi_occup_ = iBooker.book1D(histo_n, histo_t, 55, 0., 55. );
128 
129  sprintf (histo_n, "Simhit_occupancy_MB2" );
130  sprintf (histo_t, "Simhit_occupancy_MB2" );
131  meMB2_sim_occup_ = iBooker.book1D(histo_n, histo_t, 63, 0., 63. );
132 
133  sprintf (histo_n, "Digi_occupancy_MB2" );
134  sprintf (histo_t, "Digi_occupancy_MB2" );
135  meMB2_digi_occup_ = iBooker.book1D(histo_n, histo_t, 63, 0., 63. );
136 
137  sprintf (histo_n, "Simhit_occupancy_MB3" );
138  sprintf (histo_t, "Simhit_occupancy_MB3" );
139  meMB3_sim_occup_ = iBooker.book1D(histo_n, histo_t, 75, 0., 75. );
140 
141  sprintf (histo_n, "Digi_occupancy_MB3" );
142  sprintf (histo_t, "Digi_occupancy_MB3" );
143  meMB3_digi_occup_ = iBooker.book1D(histo_n, histo_t, 75, 0., 75. );
144 
145  sprintf (histo_n, "Simhit_occupancy_MB4" );
146  sprintf (histo_t, "Simhit_occupancy_MB4" );
147  meMB4_sim_occup_ = iBooker.book1D(histo_n, histo_t, 99, 0., 99. );
148 
149  sprintf (histo_n, "Digi_occupancy_MB4" );
150  sprintf (histo_t, "Digi_occupancy_MB4" );
151  meMB4_digi_occup_ = iBooker.book1D(histo_n, histo_t, 99, 0., 99. );
152 
153  // Begona's option
154  char stringcham[40];
155  for ( int slnum = 1; slnum < 62; ++slnum ) {
156  sprintf(stringcham, "DigiTimeBox_slid_%d", slnum) ;
157  meDigiHisto_ = iBooker.book1D(stringcham, stringcham, 100,0,1200);
158  meDigiTimeBox_SL_.push_back(meDigiHisto_);
159  }
160 }
161 
162 void MuonDTDigis::analyze(const Event & event, const EventSetup& eventSetup){
163 
164  if(verbose_)
165  cout << "--- [MuonDTDigis] Analysing Event: #Run: " << event.id().run()
166  << " #Event: " << event.id().event() << endl;
167 
168  // Get the DT Geometry
169  ESHandle<DTGeometry> muonGeom;
170  eventSetup.get<MuonGeometryRecord>().get(muonGeom);
171 
172  // Get the Digi collection from the event
173  Handle<DTDigiCollection> dtDigis;
174  event.getByToken(DigiToken_, dtDigis);
175 
176  // Get simhits
178  event.getByToken(SimHitToken_, simHits);
179 
180 
181  int num_mudigis;
182  int num_musimhits;
183  int num_digis;
184  int num_digis_layer;
185  int cham_num ;
186  int wire_touched;
187  num_digis = 0;
188  num_mudigis = 0;
189  num_musimhits = 0;
190  DTWireIdMap wireMap;
191 
192  for(vector<PSimHit>::const_iterator hit = simHits->begin();
193  hit != simHits->end(); hit++){
194  // Create the id of the wire, the simHits in the DT known also the wireId
195  DTWireId wireId(hit->detUnitId());
196  // cout << " PSimHits wire id " << wireId << " part type " << hit->particleType() << endl;
197 
198  // Fill the map
199  wireMap[wireId].push_back(&(*hit));
200 
201  LocalPoint entryP = hit->entryPoint();
202  LocalPoint exitP = hit->exitPoint();
203  int partType = hit->particleType();
204  if ( abs(partType) == 13 ) num_musimhits++;
205 
206  if ( wireId.station() == 1 && abs(partType) == 13 ) meMB1_sim_occup_->Fill(wireId.wire());
207  if ( wireId.station() == 2 && abs(partType) == 13 ) meMB2_sim_occup_->Fill(wireId.wire());
208  if ( wireId.station() == 3 && abs(partType) == 13 ) meMB3_sim_occup_->Fill(wireId.wire());
209  if ( wireId.station() == 4 && abs(partType) == 13 ) meMB4_sim_occup_->Fill(wireId.wire());
210 
211 
212  float path = (exitP-entryP).mag();
213  float path_x = fabs((exitP-entryP).x());
214 
215  hAllHits->Fill(entryP.x(),exitP.x(),
216  entryP.y(),exitP.y(),
217  entryP.z(),exitP.z(),
218  path , path_x,
219  partType, hit->processType(),
220  hit->pabs());
221 
222  }
223 
224  // cout << "num muon simhits " << num_musimhits << endl;
225 
227  for (detUnitIt=dtDigis->begin();
228  detUnitIt!=dtDigis->end();
229  ++detUnitIt){
230 
231  const DTLayerId& id = (*detUnitIt).first;
232  const DTDigiCollection::Range& range = (*detUnitIt).second;
233 
234  num_digis_layer = 0 ;
235  cham_num = 0 ;
236  wire_touched = 0;
237 
238  // Loop over the digis of this DetUnit
239  for (DTDigiCollection::const_iterator digiIt = range.first;
240  digiIt!=range.second;
241  ++digiIt){
242  // cout<<" Wire: "<<(*digiIt).wire()<<endl
243  // <<" digi time (ns): "<<(*digiIt).time()<<endl;
244 
245  num_digis++;
246  num_digis_layer++;
247  if (num_digis_layer > 1 )
248  {
249  if ( (*digiIt).wire() == wire_touched )
250  {
251  meWire_DoubleDigi_->Fill((*digiIt).wire());
252  // cout << "old & new wire " << wire_touched << " " << (*digiIt).wire() << endl;
253  }
254  }
255  wire_touched = (*digiIt).wire();
256 
257  meDigiTimeBox_->Fill((*digiIt).time());
258  if (id.wheel() == -2 ) meDigiTimeBox_wheel2m_->Fill((*digiIt).time());
259  if (id.wheel() == -1 ) meDigiTimeBox_wheel1m_->Fill((*digiIt).time());
260  if (id.wheel() == 0 ) meDigiTimeBox_wheel0_ ->Fill((*digiIt).time());
261  if (id.wheel() == 1 ) meDigiTimeBox_wheel1p_->Fill((*digiIt).time());
262  if (id.wheel() == 2 ) meDigiTimeBox_wheel2p_->Fill((*digiIt).time());
263 
264  // Superlayer number and fill histo with digi timebox
265 
266  cham_num = (id.wheel() +2)*12 + (id.station() -1)*3 + id.superlayer();
267  // cout << " Histo number " << cham_num << endl;
268 
269  meDigiTimeBox_SL_[cham_num]->Fill((*digiIt).time());
270 
271  // cout << " size de digis " << (*digiIt).size() << endl;
272 
273  DTWireId wireId(id,(*digiIt).wire());
274  if (wireId.station() == 1 ) meMB1_digi_occup_->Fill((*digiIt).wire());
275  if (wireId.station() == 2 ) meMB2_digi_occup_->Fill((*digiIt).wire());
276  if (wireId.station() == 3 ) meMB3_digi_occup_->Fill((*digiIt).wire());
277  if (wireId.station() == 4 ) meMB4_digi_occup_->Fill((*digiIt).wire());
278 
279  int mu=0;
280  float theta = 0;
281 
282  for(vector<const PSimHit*>::iterator hit = wireMap[wireId].begin();
283  hit != wireMap[wireId].end(); hit++)
284  if( abs((*hit)->particleType()) == 13){
285  theta = atan( (*hit)->momentumAtEntry().x()/ (-(*hit)->momentumAtEntry().z()) )*180/M_PI;
286  // cout<<"momentum x: "<<(*hit)->momentumAtEntry().x()<<endl
287  // <<"momentum z: "<<(*hit)->momentumAtEntry().z()<<endl
288  // <<"atan: "<<theta<<endl;
289  mu++;
290  }
291 
292  if( mu ) num_mudigis++;
293 
294  if(mu && theta){
295  hDigis_global->Fill((*digiIt).time(),theta,id.superlayer());
296  //filling digi histos for wheel and for RZ and RPhi
297  WheelHistos(id.wheel())->Fill((*digiIt).time(),theta,id.superlayer());
298  }
299 
300  }// for digis in layer
301 
302  meDoubleDigi_->Fill( (float)num_digis_layer );
303 
304  }// for layers
305 
306  //cout << "num_digis " << num_digis << "mu digis " << num_mudigis << endl;
307 
308  if (num_musimhits != 0) {
309  meDigiEfficiencyMu_->Fill( (float)num_mudigis/(float)num_musimhits );
310  meDigiEfficiency_->Fill( (float)num_digis/(float)num_musimhits );
311  }
312 
313  meSimvsDigi_->Fill( (float)num_musimhits, (float)num_digis ) ;
314  // cout<<"--------------"<<endl;
315 
316 }
317 
319  switch(abs(wheel)){
320 
321  case 0: return hDigis_W0.get();
322 
323  case 1: return hDigis_W1.get();
324 
325  case 2: return hDigis_W2.get();
326 
327  default: return nullptr;
328  }
329 }
333 
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::map< DTWireId, std::vector< const PSimHit * > > DTWireIdMap
Definition: MuonDTDigis.h:67
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Geom::Theta< T > theta() const
T y() const
Definition: PV3DBase.h:63
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: MuonDTDigis.cc:43
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
T z() const
Definition: PV3DBase.h:64
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const int mu
Definition: Constants.h:22
hDigis * WheelHistos(int wheel)
Definition: MuonDTDigis.cc:318
#define M_PI
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:136
~MuonDTDigis() override
Definition: MuonDTDigis.cc:38
std::vector< DTDigi >::const_iterator const_iterator
#define begin
Definition: vmac.h:32
HLT enums.
T get() const
Definition: EventSetup.h:63
std::pair< const_iterator, const_iterator > Range
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
Definition: MuonDTDigis.cc:162
T x() const
Definition: PV3DBase.h:62
Definition: event.py:1
Definition: Run.h:44
MuonDTDigis(const edm::ParameterSet &pset)
Definition: MuonDTDigis.cc:19