CMS 3D CMS Logo

RPCDigiValid.cc
Go to the documentation of this file.
2 
5 
11 
12 #include <cmath>
13 #include <fmt/format.h>
14 
15 using namespace std;
16 using namespace edm;
17 
19  // Init the tokens for run data retrieval - stanislav
20  // ps.getUntackedParameter<InputTag> retrieves a InputTag from the
21  // configuration. The second param is default value module, instance and
22  // process labels may be passed in a single string if separated by colon ':'
23  // (@see the edm::InputTag constructor documentation)
24  simHitToken_ = consumes<PSimHitContainer>(
25  ps.getUntrackedParameter<edm::InputTag>("simHitTag", edm::InputTag("g4SimHits:MuonRPCHits")));
26  rpcDigiToken_ = consumes<RPCDigiCollection>(
27  ps.getUntrackedParameter<edm::InputTag>("rpcDigiTag", edm::InputTag("simMuonRPCDigis")));
28 
29  isDigiTimeAvailable_ = ps.getUntrackedParameter<bool>("digiTime", false);
30 
31  rpcGeomToken_ = esConsumes();
32 }
33 
35  // Get the RPC Geometry
36  auto rpcGeom = eventSetup.getHandle(rpcGeomToken_);
37 
38  edm::Handle<PSimHitContainer> simHitHandle;
39  edm::Handle<RPCDigiCollection> rpcDigisHandle;
40  event.getByToken(simHitToken_, simHitHandle);
41  event.getByToken(rpcDigiToken_, rpcDigisHandle);
42 
43  // loop over Simhit
44  std::map<const RPCRoll *, std::vector<double>> detToSimHitXsMap;
45  for (auto simIt = simHitHandle->begin(); simIt != simHitHandle->end(); ++simIt) {
46  const RPCDetId rsid = simIt->detUnitId();
47  const RPCRoll *roll = dynamic_cast<const RPCRoll *>(rpcGeom->roll(rsid));
48  if (!roll)
49  continue;
50 
51  if (detToSimHitXsMap.find(roll) == detToSimHitXsMap.end())
52  detToSimHitXsMap[roll] = std::vector<double>();
53  detToSimHitXsMap[roll].push_back(simIt->localPosition().x());
54 
55  const int region = rsid.region();
56  const GlobalPoint gp = roll->toGlobal(simIt->localPosition());
57 
58  hRZ_->Fill(gp.z(), gp.perp());
59 
60  if (region == 0) {
61  // Barrel
62  hXY_Barrel_->Fill(gp.x(), gp.y());
63 
64  const int station = rsid.station();
65  const int layer = rsid.layer();
66  const int stla = (station <= 2) ? (2 * (station - 1) + layer) : (station + 2);
67 
68  auto match = hZPhi_.find(stla);
69  if (match != hZPhi_.end()) {
70  const double phiInDeg = 180. * gp.barePhi() / TMath::Pi();
71  match->second->Fill(gp.z(), phiInDeg);
72  }
73  } else {
74  // Endcap
75  const int disk = region * rsid.station();
76  auto match = hXY_Endcap_.find(disk);
77  if (match != hXY_Endcap_.end())
78  match->second->Fill(gp.x(), gp.y());
79  }
80  }
81  for (auto detToSimHitXs : detToSimHitXsMap) {
82  hNSimHitPerRoll_->Fill(detToSimHitXs.second.size());
83  }
84 
85  // loop over Digis
86  std::map<const RPCRoll *, std::vector<double>> detToDigiXsMap;
87  for (auto detUnitIt = rpcDigisHandle->begin(); detUnitIt != rpcDigisHandle->end(); ++detUnitIt) {
88  const RPCDetId rsid = (*detUnitIt).first;
89  const RPCRoll *roll = dynamic_cast<const RPCRoll *>(rpcGeom->roll(rsid));
90  if (!roll)
91  continue;
92  const int region = rsid.region();
93 
94  const RPCDigiCollection::Range &range = (*detUnitIt).second;
95 
96  for (auto digiIt = range.first; digiIt != range.second; ++digiIt) {
97  // Strip profile
98  const int strip = digiIt->strip();
99  hStripProf_->Fill(strip);
100 
101  if (region == 0) {
102  const int station = rsid.station();
103  if (station == 1 or station == 2)
104  hStripProf_RB12_->Fill(strip);
105  else if (station == 3 or station == 4)
106  hStripProf_RB34_->Fill(strip);
107  } else {
108  if (roll->isIRPC())
109  hStripProf_IRPC_->Fill(strip);
110  else
111  hStripProf_Endcap_->Fill(strip);
112  }
113 
114  // Bunch crossing
115  const int bx = digiIt->bx();
116  hBxDist_->Fill(bx);
117  // bx for 4 endcaps
118  if (rsid.station() == 4) {
119  if (region == 1) {
120  hBxDisc_4Plus_->Fill(bx);
121  } else if (region == -1) {
122  hBxDisc_4Min_->Fill(bx);
123  }
124  }
125 
126  // Fill timing information
127  if (isDigiTimeAvailable_) {
128  const double digiTime = digiIt->hasTime() ? digiIt->time() : digiIt->bx() * 25;
129  hDigiTimeAll_->Fill(digiTime);
130  if (digiIt->hasTime()) {
131  hDigiTime_->Fill(digiTime);
132  if (roll->isIRPC())
133  hDigiTimeIRPC_->Fill(digiTime);
134  else
135  hDigiTimeNoIRPC_->Fill(digiTime);
136  }
137  }
138 
139  // Keep digi position
140  const double digiX = roll->centreOfStrip(digiIt->strip()).x();
141  if (detToDigiXsMap.find(roll) == detToDigiXsMap.end())
142  detToDigiXsMap[roll] = std::vector<double>();
143  detToDigiXsMap[roll].push_back(digiX);
144  }
145  }
146  for (auto detToDigiXs : detToDigiXsMap) {
147  const auto digiXs = detToDigiXs.second;
148  const int nDigi = digiXs.size();
149  hNDigiPerRoll_->Fill(nDigi);
150 
151  // Fill residual plots, only for nDigi==1 and nSimHit==1
152  const auto roll = detToDigiXs.first;
153  const auto detId = roll->id();
154  if (nDigi != 1)
155  continue;
156  if (detToSimHitXsMap.find(roll) == detToSimHitXsMap.end())
157  continue;
158 
159  const auto simHitXs = detToSimHitXsMap[roll];
160  const int nSimHit = simHitXs.size();
161  if (nSimHit != 1)
162  continue;
163 
164  const double dx = digiXs[0] - simHitXs[0];
165  hRes_->Fill(dx);
166  if (roll->isBarrel()) {
167  const int wheel = detId.ring(); // ring() is wheel number for Barrel
168  const int station = detId.station();
169  const int layer = detId.layer();
170  const int stla = (station <= 2) ? (2 * (station - 1) + layer) : (station + 2);
171 
172  auto matchLayer = hResBarrelLayers_.find(stla);
173  if (matchLayer != hResBarrelLayers_.end())
174  matchLayer->second->Fill(dx);
175 
176  auto matchWheel = hResBarrelWheels_.find(wheel);
177  if (matchWheel != hResBarrelWheels_.end())
178  matchWheel->second->Fill(dx);
179  } else {
180  const int disk = detId.region() * detId.station();
181  auto matchDisk = hResEndcapDisks_.find(disk);
182  if (matchDisk != hResEndcapDisks_.end())
183  matchDisk->second->Fill(dx);
184 
185  auto matchRing = hResEndcapRings_.find(detId.ring());
186  if (matchRing != hResEndcapRings_.end())
187  matchRing->second->Fill(dx);
188  }
189  }
190 }
191 
193  booker.setCurrentFolder("RPCDigisV/RPCDigis");
194 
195  // Define binnings of 2D-histograms
196  const double maxZ = 1100;
197  const int nbinsZ = 220; // bin width: 10cm
198  const double maxXY = 800;
199  const int nbinsXY = 160; // bin width: 10cm
200  const double minR = 100, maxR = 800;
201  const int nbinsR = 70; // bin width: 10cm
202  const int nbinsPhi = 72; // bin width: 5 degree
203  const double maxBarrelZ = 700;
204  const int nbinsBarrelZ = 140; // bin width: 10cm
205 
206  // RZ plot
207  hRZ_ = booker.book2D("RZ", "R-Z view;Z (cm);R (cm)", nbinsZ, -maxZ, maxZ, nbinsR, minR, maxR);
208  hRZ_->setOption("colz");
209 
210  // XY plots
211  hXY_Barrel_ = booker.book2D("XY_Barrel", "X-Y view of Barrel", nbinsXY, -maxXY, maxXY, nbinsXY, -maxXY, maxXY);
212  hXY_Barrel_->setOption("colz");
213  for (int disk = 1; disk <= 4; ++disk) {
214  const std::string meNameP = fmt::format("XY_Endcap_p{:1d}", disk);
215  const std::string meNameN = fmt::format("XY_Endcap_m{:1d}", disk);
216  const std::string meTitleP = fmt::format("X-Y view of Endcap{:+1d};X (cm);Y (cm)", disk);
217  const std::string meTitleN = fmt::format("X-Y view of Endcap{:+1d};X (cm);Y (cm)", -disk);
218  hXY_Endcap_[disk] = booker.book2D(meNameP, meTitleP, nbinsXY, -maxXY, maxXY, nbinsXY, -maxXY, maxXY);
219  hXY_Endcap_[-disk] = booker.book2D(meNameN, meTitleN, nbinsXY, -maxXY, maxXY, nbinsXY, -maxXY, maxXY);
220  hXY_Endcap_[disk]->setOption("colz");
221  hXY_Endcap_[-disk]->setOption("colz");
222  }
223 
224  // Z-phi plots
225  for (int layer = 1; layer <= 6; ++layer) {
226  const std::string meName = fmt::format("ZPhi_Layer{:1d}", layer);
227  const std::string meTitle = fmt::format("Z-#phi view of Layer{:1d};Z (cm);#phi (degree)", layer);
228  hZPhi_[layer] = booker.book2D(meName, meTitle, nbinsBarrelZ, -maxBarrelZ, maxBarrelZ, nbinsPhi, -180, 180);
229  hZPhi_[layer]->setOption("colz");
230  }
231 
232  // Strip profile
233  hStripProf_ = booker.book1D("Strip_Profile", "Strip_Profile;Strip Number", 100, 0, 100);
234  hStripProf_RB12_ = booker.book1D("Strip_Profile_RB12", "Strip Profile RB1 and RB2;Strip Number", 92, 0, 92);
235  hStripProf_RB34_ = booker.book1D("Strip_Profile_RB34", "Strip Profile RB3 and RB4;Strip Number", 62, 0, 62);
236  hStripProf_Endcap_ = booker.book1D("Strip_Profile_Endcap", "Strip Profile Endcap;Strip Number", 40, 0, 40);
237  hStripProf_IRPC_ = booker.book1D("Strip_Profile_IRPC", "Strip Profile IRPC;Strip Number", 100, 0, 100);
238 
239  // Bunch crossing
240  hBxDist_ = booker.book1D("Bunch_Crossing", "Bunch Crossing;Bunch crossing", 20, -10., 10.);
241  hBxDisc_4Plus_ = booker.book1D("BxDisc_4Plus", "BxDisc_4Plus", 20, -10., 10.);
242  hBxDisc_4Min_ = booker.book1D("BxDisc_4Min", "BxDisc_4Min", 20, -10., 10.);
243 
244  // Timing informations
245  if (isDigiTimeAvailable_) {
246  hDigiTimeAll_ =
247  booker.book1D("DigiTimeAll", "Digi time including present electronics;Digi time (ns)", 100, -12.5, 12.5);
248  hDigiTime_ = booker.book1D("DigiTime", "Digi time only with timing information;Digi time (ns)", 100, -12.5, 12.5);
249  hDigiTimeIRPC_ = booker.book1D("DigiTimeIRPC", "IRPC Digi time;Digi time (ns)", 100, -12.5, 12.5);
250  hDigiTimeNoIRPC_ = booker.book1D("DigiTimeNoIRPC", "non-IRPC Digi time;Digi time (ns)", 100, -12.5, 12.5);
251  }
252 
253  // SimHit and Digi multiplicity per roll
254  hNSimHitPerRoll_ = booker.book1D("NSimHitPerRoll", "SimHit multiplicity per Roll;Multiplicity", 10, 0, 10);
255  hNDigiPerRoll_ = booker.book1D("NDigiPerRoll", "Digi multiplicity per Roll;Multiplicity", 10, 0, 10);
256 
257  // Residual of SimHit-Digi x-position
258  hRes_ = booker.book1D("Digi_SimHit_Difference", "Digi-SimHit difference;dx (cm)", 100, -8, 8);
259 
260  for (int layer = 1; layer <= 6; ++layer) {
261  const std::string meName = fmt::format("Residual_Barrel_Layer{:1d}", layer);
262  const std::string meTitle = fmt::format("Residual of Barrel Layer{:1d};dx (cm)", layer);
263  hResBarrelLayers_[layer] = booker.book1D(meName, meTitle, 100, -8, 8);
264  }
265 
266  hResBarrelWheels_[-2] = booker.book1D("Residual_Barrel_Wheel_m2", "Residual of Barrel Wheel-2;dx (cm)", 100, -8, 8);
267  hResBarrelWheels_[-1] = booker.book1D("Residual_Barrel_Wheel_m1", "Residual of Barrel Wheel-1;dx (cm)", 100, -8, 8);
268  hResBarrelWheels_[+0] = booker.book1D("Residual_Barrel_Wheel_00", "Residual of Barrel Wheel 0;dx (cm)", 100, -8, 8);
269  hResBarrelWheels_[+1] = booker.book1D("Residual_Barrel_Wheel_p1", "Residual of Barrel Wheel+1;dx (cm)", 100, -8, 8);
270  hResBarrelWheels_[+2] = booker.book1D("Residual_Barrel_Wheel_p2", "Residual of Barrel Wheel+2;dx (cm)", 100, -8, 8);
271 
272  for (int disk = 1; disk <= 4; ++disk) {
273  const std::string meNameP = fmt::format("Residual_Endcap_Disk_p{:1d}", disk);
274  const std::string meNameN = fmt::format("Residual_Endcap_Disk_m{:1d}", disk);
275  const std::string meTitleP = fmt::format("Residual of Endcap Disk{:+1d};dx (cm)", disk);
276  const std::string meTitleN = fmt::format("Residual of Endcap Disk{:+1d};dx (cm)", -disk);
277  hResEndcapDisks_[+disk] = booker.book1D(meNameP, meTitleP, 100, -8, 8);
278  hResEndcapDisks_[-disk] = booker.book1D(meNameN, meTitleN, 100, -8, 8);
279  }
280 
281  hResEndcapRings_[1] = booker.book1D("Residual_Endcap_Ring1", "Residual of Endcap Ring1;dx (cm)", 100, -12, 12);
282  hResEndcapRings_[2] = booker.book1D("Residual_Endcap_Ring2", "Residual of Endcap Ring2;dx (cm)", 100, -8, 8);
283  hResEndcapRings_[3] = booker.book1D("Residual_Endcap_Ring3", "Residual of Endcap Ring3;dx (cm)", 100, -8, 8);
284 }
const double Pi
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:36
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
virtual void setOption(const char *option)
RPCDigiValid(const edm::ParameterSet &ps)
Definition: RPCDigiValid.cc:18
T getUntrackedParameter(std::string const &, T const &) const
void analyze(const edm::Event &e, const edm::EventSetup &c) override
Definition: RPCDigiValid.cc:34
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
std::pair< const_iterator, const_iterator > Range
MonitorElement * book2D(TString const &name, TString const &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, FUNC onbooking=NOOP())
Definition: DQMStore.h:221
int station() const
Definition: RPCDetId.h:78
int region() const
Region id: 0 for Barrel, +/-1 For +/- Endcap.
Definition: RPCDetId.h:53
HLT enums.
float x
int layer() const
Definition: RPCDetId.h:85
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
Definition: event.py:1
Definition: Run.h:45