CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
DTSegtoRPC Class Reference

#include <DTSegtoRPC.h>

Public Member Functions

 DTSegtoRPC (edm::ConsumesCollector iC, const edm::ParameterSet &)
 
std::unique_ptr< RPCRecHitCollectionthePoints (DTRecSegment4DCollection const *all4DSegments, edm::EventSetup const &iSetup, bool debug, double eyr)
 

Private Attributes

edm::ESGetToken< DTGeometry, MuonGeometryRecorddtGeoToken_
 
edm::ESGetToken< DTObjectMap, MuonGeometryRecorddtMapToken_
 
bool incldt
 
bool incldtMB4
 
double MaxD
 
double MaxDistanceBetweenSegments
 
double MaxDrb4
 
int maxPhiBX
 
double MinCosAng
 
int minPhiBX
 
edm::ESGetToken< RPCGeometry, MuonGeometryRecordrpcGeoToken_
 

Detailed Description

Definition at line 18 of file DTSegtoRPC.h.

Constructor & Destructor Documentation

◆ DTSegtoRPC()

DTSegtoRPC::DTSegtoRPC ( edm::ConsumesCollector  iC,
const edm::ParameterSet iConfig 
)
explicit

Definition at line 40 of file DTSegtoRPC.cc.

References edm::ParameterSet::getParameter(), incldt, incldtMB4, MaxD, MaxDistanceBetweenSegments, MaxDrb4, maxPhiBX, MinCosAng, and minPhiBX.

42  minPhiBX = iConfig.getParameter<int>("minBX");
43  maxPhiBX = iConfig.getParameter<int>("maxBX");
44  incldt = true;
45  incldtMB4 = true;
46 
47  //By now hard coded parameters
48  MinCosAng = 0.85;
49  MaxD = 80.;
50  MaxDrb4 = 150.;
52 }
double MaxDistanceBetweenSegments
Definition: DTSegtoRPC.h:36
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
int minPhiBX
Definition: DTSegtoRPC.h:37
edm::ESGetToken< DTObjectMap, MuonGeometryRecord > dtMapToken_
Definition: DTSegtoRPC.h:29
double MaxD
Definition: DTSegtoRPC.h:34
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeoToken_
Definition: DTSegtoRPC.h:28
double MaxDrb4
Definition: DTSegtoRPC.h:35
bool incldtMB4
Definition: DTSegtoRPC.h:32
int maxPhiBX
Definition: DTSegtoRPC.h:38
bool incldt
Definition: DTSegtoRPC.h:31
edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeoToken_
Definition: DTSegtoRPC.h:27
double MinCosAng
Definition: DTSegtoRPC.h:33

Member Function Documentation

◆ thePoints()

std::unique_ptr< RPCRecHitCollection > DTSegtoRPC::thePoints ( DTRecSegment4DCollection const *  all4DSegments,
edm::EventSetup const &  iSetup,
bool  debug,
double  eyr 
)

Definition at line 54 of file DTSegtoRPC.cc.

References cms::cuda::assert(), debug, distsector(), distwheel(), dtGeoToken_, dtMapToken_, PVValHelper::dx, PVValHelper::dy, PixelTestBeamValidation_cfi::Dy, PVValHelper::dz, spr::find(), edm::EventSetup::getHandle(), DTObjectMap::getRolls(), RPCRoll::id(), DTGeometry::idToDet(), incldt, incldtMB4, RectangularStripTopology::localPosition(), LogDebug, visualization-live-secondInstance_cfg::m, mag(), MaxD, MaxDistanceBetweenSegments, MaxDrb4, MinCosAng, minPhiBX, RPCGeomServ::name(), RPCRoll::nstrips(), RectangularStripTopology::pitch(), DetId::rawId(), RPCGeometry::roll(), rpcGeoToken_, DTChamberId::sector(), RPCRecHit::setTimeAndError(), mathSSE::sqrt(), DTChamberId::station(), AlCaHLTBitMon_QueryRunRegistry::string, RectangularStripTopology::stripLength(), GeomDet::surface(), RPCRoll::topology(), DTChamberId::wheel(), X, PV3DBase< T, PVType, FrameType >::x(), TrackerOfflineValidation_Dqm_cff::xmax, TrackerOfflineValidation_Dqm_cff::xmin, beamSpotPI::Y, PV3DBase< T, PVType, FrameType >::y(), beamSpotPI::Z, and PV3DBase< T, PVType, FrameType >::z().

57  {
58  auto _ThePoints = std::make_unique<RPCRecHitCollection>();
59  edm::OwnVector<RPCRecHit> RPCPointVector;
60  std::vector<uint32_t> extrapolatedRolls;
61 
62  if (all4DSegments->size() > 8) {
63  if (debug)
64  LogDebug("DTSegtoRPC") << "Too many segments in this event we are not doing the extrapolation" << std::endl;
65  } else {
66  edm::ESHandle<RPCGeometry> rpcGeo = iSetup.getHandle(rpcGeoToken_);
67  edm::ESHandle<DTGeometry> dtGeo = iSetup.getHandle(dtGeoToken_);
68  edm::ESHandle<DTObjectMap> dtMap = iSetup.getHandle(dtMapToken_);
69 
70  std::map<DTChamberId, int> DTSegmentCounter;
72 
73  for (segment = all4DSegments->begin(); segment != all4DSegments->end(); ++segment) {
74  DTSegmentCounter[segment->chamberId()]++;
75  }
76 
77  const float bunchCrossTimeDiff = 25.; // time between bunch crossings
78 
79  if (incldt) {
80  for (segment = all4DSegments->begin(); segment != all4DSegments->end(); ++segment) {
81  DTChamberId DTId = segment->chamberId();
82 
83  if (debug)
84  LogDebug("DTSegtoRPC") << "DT \t \t This Segment is in Chamber id: " << DTId << std::endl;
85  if (debug)
86  LogDebug("DTSegtoRPC") << "DT \t \t Number of segments in this DT = " << DTSegmentCounter[DTId] << std::endl;
87  if (debug)
88  LogDebug("DTSegtoRPC") << "DT \t \t Is the only one in this DT? and is not in the 4th Station?" << std::endl;
89 
90  if (DTSegmentCounter[DTId] != 1 || DTId.station() == 4) {
91  if (debug)
92  LogDebug("DTSegtoRPC") << "DT \t \t More than one segment in this chamber, or we are in Station 4"
93  << std::endl;
94  continue;
95  }
96 
97  int dtWheel = DTId.wheel();
98  int dtStation = DTId.station();
99  int dtSector = DTId.sector();
100 
101  LocalPoint segmentPosition = segment->localPosition();
102  LocalVector segmentDirection = segment->localDirection();
103 
104  const GeomDet* gdet = dtGeo->idToDet(segment->geographicalId());
105  const BoundPlane& DTSurface = gdet->surface();
106 
107  //check if the dimension of the segment is 4
108 
109  if (debug)
110  LogDebug("DTSegtoRPC") << "DT \t \t Is the segment 4D?" << std::endl;
111 
112  if (segment->dimension() != 4) {
113  if (debug)
114  LogDebug("DTSegtoRPC") << "DT \t \t no" << std::endl;
115  continue;
116  }
117 
118  if (debug)
119  LogDebug("DTSegtoRPC") << "DT \t \t yes" << std::endl;
120  if (debug)
121  LogDebug("DTSegtoRPC") << "DT \t \t DT Segment Dimension " << segment->dimension() << std::endl;
122 
123  float Xo = segmentPosition.x();
124  float Yo = segmentPosition.y();
125  float Zo = segmentPosition.z();
126  float dx = segmentDirection.x();
127  float dy = segmentDirection.y();
128  float dz = segmentDirection.z();
129 
130  float myPhiTime = -9999.;
131  float myPhiTimeErr = -9999.;
132  int myPhiBx = -99;
133 
134  if (segment->hasPhi()) {
135  myPhiTime = segment->phiSegment()->t0();
136  myPhiBx = round(segment->phiSegment()->t0() / bunchCrossTimeDiff);
137  if (debug)
138  LogDebug("DTSegtoRPC") << "segment t0 = " << myPhiTime << "\tround = " << myPhiBx << std::endl;
139  }
140  if (!(myPhiBx <= maxPhiBX && myPhiBx >= minPhiBX)) //limit only to the RPC readout range
141  continue;
142 
143  if (debug)
144  LogDebug("DTSegtoRPC") << "Creating the DTIndex" << std::endl;
145  DTStationIndex theindex(0, dtWheel, dtSector, dtStation);
146  if (debug)
147  LogDebug("DTSegtoRPC") << "Getting the Rolls for the given index" << std::endl;
148  std::set<RPCDetId> rollsForThisDT = dtMap->getRolls(theindex);
149 
150  if (debug)
151  LogDebug("DTSegtoRPC") << "DT \t \t Number of rolls for this DT = " << rollsForThisDT.size() << std::endl;
152 
153  assert(!rollsForThisDT.empty());
154 
155  if (debug)
156  LogDebug("DTSegtoRPC") << "DT \t \t Loop over all the rolls asociated to this DT" << std::endl;
157  for (std::set<RPCDetId>::iterator iteraRoll = rollsForThisDT.begin(); iteraRoll != rollsForThisDT.end();
158  iteraRoll++) {
159  const RPCRoll* rollasociated = rpcGeo->roll(*iteraRoll);
160  RPCDetId rpcId = rollasociated->id();
161  const BoundPlane& RPCSurface = rollasociated->surface();
162 
163  RPCGeomServ rpcsrv(rpcId);
164  std::string nameRoll = rpcsrv.name();
165 
166  if (debug)
167  LogDebug("DTSegtoRPC") << "DT \t \t \t RollName: " << nameRoll << std::endl;
168  if (debug)
169  LogDebug("DTSegtoRPC") << "DT \t \t \t Doing the extrapolation to this roll" << std::endl;
170  if (debug)
171  LogDebug("DTSegtoRPC") << "DT \t \t \t DT Segment Direction in DTLocal " << segmentDirection << std::endl;
172  if (debug)
173  LogDebug("DTSegtoRPC") << "DT \t \t \t DT Segment Point in DTLocal " << segmentPosition << std::endl;
174 
175  GlobalPoint CenterPointRollGlobal = RPCSurface.toGlobal(LocalPoint(0, 0, 0));
176 
177  LocalPoint CenterRollinDTFrame = DTSurface.toLocal(CenterPointRollGlobal);
178 
179  if (debug)
180  LogDebug("DTSegtoRPC") << "DT \t \t \t Center (0,0,0) Roll In DTLocal" << CenterRollinDTFrame << std::endl;
181  if (debug)
182  LogDebug("DTSegtoRPC") << "DT \t \t \t Center (0,0,0) of the Roll in Global" << CenterPointRollGlobal
183  << std::endl;
184 
185  float D = CenterRollinDTFrame.z();
186 
187  float X = Xo + dx * D / dz;
188  float Y = Yo + dy * D / dz;
189  float Z = D;
190 
191  const RectangularStripTopology* top_ =
192  dynamic_cast<const RectangularStripTopology*>(&(rollasociated->topology()));
193  LocalPoint xmin = top_->localPosition(0.);
194  if (debug)
195  LogDebug("DTSegtoRPC") << "DT \t \t \t xmin of this Roll " << xmin << "cm" << std::endl;
196  LocalPoint xmax = top_->localPosition((float)rollasociated->nstrips());
197  if (debug)
198  LogDebug("DTSegtoRPC") << "DT \t \t \t xmax of this Roll " << xmax << "cm" << std::endl;
199  float rsize = fabs(xmax.x() - xmin.x());
200  if (debug)
201  LogDebug("DTSegtoRPC") << "DT \t \t \t Roll Size " << rsize << "cm" << std::endl;
202  float stripl = top_->stripLength();
203 
204  float stripw = top_->pitch();
205 
206  if (debug)
207  LogDebug("DTSegtoRPC") << "DT \t \t \t Strip Lenght " << stripl << "cm" << std::endl;
208  if (debug)
209  LogDebug("DTSegtoRPC") << "DT \t \t \t Strip Width " << stripw << "cm" << std::endl;
210  if (debug)
211  LogDebug("DTSegtoRPC") << "DT \t \t \t X Predicted in DTLocal= " << X << "cm" << std::endl;
212  if (debug)
213  LogDebug("DTSegtoRPC") << "DT \t \t \t Y Predicted in DTLocal= " << Y << "cm" << std::endl;
214  if (debug)
215  LogDebug("DTSegtoRPC") << "DT \t \t \t Z Predicted in DTLocal= " << Z << "cm" << std::endl;
216 
217  float extrapolatedDistance = sqrt((X - Xo) * (X - Xo) + (Y - Yo) * (Y - Yo) + (Z - Zo) * (Z - Zo));
218 
219  if (debug)
220  LogDebug("DTSegtoRPC") << "DT \t \t \t Is the distance of extrapolation less than MaxD? ="
221  << extrapolatedDistance << "cm"
222  << "MaxD=" << MaxD << "cm" << std::endl;
223 
224  if (extrapolatedDistance <= MaxD) {
225  if (debug)
226  LogDebug("DTSegtoRPC") << "DT \t \t \t yes" << std::endl;
227  GlobalPoint GlobalPointExtrapolated = DTSurface.toGlobal(LocalPoint(X, Y, Z));
228  if (debug)
229  LogDebug("DTSegtoRPC") << "DT \t \t \t Point ExtraPolated in Global" << GlobalPointExtrapolated
230  << std::endl;
231  LocalPoint PointExtrapolatedRPCFrame = RPCSurface.toLocal(GlobalPointExtrapolated);
232 
233  if (debug)
234  LogDebug("DTSegtoRPC") << "DT \t \t \t Point Extrapolated in RPCLocal" << PointExtrapolatedRPCFrame
235  << std::endl;
236  if (debug)
237  LogDebug("DTSegtoRPC") << "DT \t \t \t Corner of the Roll = (" << rsize * eyr << "," << stripl * eyr
238  << ")" << std::endl;
239  if (debug)
240  LogDebug("DTSegtoRPC") << "DT \t \t \t Info About the Point Extrapolated in X Abs ("
241  << fabs(PointExtrapolatedRPCFrame.x()) << ","
242  << fabs(PointExtrapolatedRPCFrame.y()) << ","
243  << fabs(PointExtrapolatedRPCFrame.z()) << ")" << std::endl;
244  if (debug)
245  LogDebug("DTSegtoRPC") << "DT \t \t \t Does the extrapolation go inside this roll?" << std::endl;
246 
247  if (fabs(PointExtrapolatedRPCFrame.z()) < 1. && fabs(PointExtrapolatedRPCFrame.x()) < rsize * eyr &&
248  fabs(PointExtrapolatedRPCFrame.y()) < stripl * eyr) {
249  if (debug)
250  LogDebug("DTSegtoRPC") << "DT \t \t \t \t yes" << std::endl;
251  if (debug)
252  LogDebug("DTSegtoRPC") << "DT \t \t \t \t Creating the RecHit" << std::endl;
253 
254  RPCRecHit RPCPoint(rpcId, myPhiBx, PointExtrapolatedRPCFrame);
255  RPCPoint.setTimeAndError(myPhiTime, myPhiTimeErr);
256 
257  if (debug)
258  LogDebug("DTSegtoRPC") << "DT \t \t \t \t Clearing the vector" << std::endl;
259  RPCPointVector.clear();
260  if (debug)
261  LogDebug("DTSegtoRPC") << "DT \t \t \t \t Pushing back" << std::endl;
262  RPCPointVector.push_back(RPCPoint);
263  if (debug)
264  LogDebug("DTSegtoRPC") << "DT \t \t \t \t Putting the vector" << std::endl;
265  _ThePoints->put(rpcId, RPCPointVector.begin(), RPCPointVector.end());
266 
267  if (debug)
268  LogDebug("DTSegtoRPC") << "DT \t \t \t \t Filling container with " << nameRoll
269  << " Point.x=" << PointExtrapolatedRPCFrame.x()
270  << " Point.y=" << PointExtrapolatedRPCFrame.y()
271  << " size=" << RPCPointVector.size() << std::endl;
272 
273  } else {
274  if (debug)
275  LogDebug("DTSegtoRPC") << "DT \t \t \t \t No the prediction is outside of this roll" << std::endl;
276  } //Condition for the right match
277  } else {
278  if (debug)
279  LogDebug("DTSegtoRPC") << "DT \t \t \t No, Exrtrapolation too long!, canceled" << std::endl;
280  } //D so big
281  } //loop over all the rolls asociated
282  }
283  }
284 
285  if (incldtMB4) {
286  if (all4DSegments->size() > 0) {
287  if (debug)
288  LogDebug("DTSegtoRPC") << "MB4 \t \t Loop Over all4DSegments " << all4DSegments->size() << std::endl;
289  extrapolatedRolls.clear();
290  for (segment = all4DSegments->begin(); segment != all4DSegments->end(); ++segment) {
291  DTChamberId DTId = segment->chamberId();
292 
293  float myPhiTime = -9999.;
294  float myPhiTimeErr = -9999.;
295  int myPhiBx = -99;
296 
297  if (segment->hasPhi()) {
298  myPhiTime = segment->phiSegment()->t0();
299  myPhiBx = round(segment->phiSegment()->t0() / bunchCrossTimeDiff);
300  if (debug)
301  LogDebug("DTSegtoRPC") << "segment t0 = " << myPhiTime << "\tround = " << myPhiBx << std::endl;
302  }
303  if (!(myPhiBx <= maxPhiBX && myPhiBx >= minPhiBX)) //limit only to the RPC readout range
304  continue;
305 
306  if (debug)
307  LogDebug("DTSegtoRPC") << "MB4 \t \t This Segment is in Chamber id: " << DTId << std::endl;
308  if (debug)
309  LogDebug("DTSegtoRPC") << "MB4 \t \t Number of segments in this DT = " << DTSegmentCounter[DTId]
310  << std::endl;
311  if (debug)
312  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Is the only one in this DT? and is in the Station 4?" << std::endl;
313 
314  if (DTSegmentCounter[DTId] == 1 && DTId.station() == 4) {
315  if (debug)
316  LogDebug("DTSegtoRPC") << "MB4 \t \t \t yes" << std::endl;
317  int dtWheel = DTId.wheel();
318  int dtStation = DTId.station();
319  int dtSector = DTId.sector();
320 
321  LocalPoint segmentPosition = segment->localPosition();
322  LocalVector segmentDirection = segment->localDirection();
323 
324  if (debug)
325  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t The Segment in MB4 is 2D?" << std::endl;
326  if (segment->dimension() == 2) {
327  if (debug)
328  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t yes" << std::endl;
329  const LocalVector& segmentDirectionMB4 = segmentDirection;
330  const LocalPoint& segmentPositionMB4 = segmentPosition;
331 
332  const BoundPlane& DTSurface4 = dtGeo->idToDet(DTId)->surface();
333 
335 
336  if (debug)
337  LogDebug("DTSegtoRPC")
338  << "MB4 \t \t \t \t Loop on segments in =sector && MB3 && adjacent sectors && y dim=4" << std::endl;
339  for (segMB3 = all4DSegments->begin(); segMB3 != all4DSegments->end(); ++segMB3) {
340  DTChamberId dtid3 = segMB3->chamberId();
341 
342  if (debug)
343  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t Segment in Chamber =" << dtid3 << std::endl;
344 
345  if (distsector(dtid3.sector(), DTId.sector()) <=
346  1 //The DT sector could be 13 or 14 and because is corrected in the calculation of the distance.
347  && distwheel(dtid3.wheel(), DTId.wheel()) <=
348  1 //The we could have segments in neighbohr wheels in pp collisions
349  && dtid3.station() == 3 && DTSegmentCounter[dtid3] == 1 && segMB3->dimension() == 4) {
350  if (debug)
351  LogDebug("DTSegtoRPC")
352  << "MB4 \t \t \t \t distsector =" << distsector(dtid3.sector(), DTId.sector()) << std::endl;
353  if (debug)
354  LogDebug("DTSegtoRPC")
355  << "MB4 \t \t \t \t distwheel =" << distwheel(dtid3.wheel(), DTId.wheel()) << std::endl;
356 
357  const GeomDet* gdet3 = dtGeo->idToDet(segMB3->geographicalId());
358  const BoundPlane& DTSurface3 = gdet3->surface();
359 
360  LocalVector segmentDirectionMB3 = segMB3->localDirection();
361  GlobalPoint segmentPositionMB3inGlobal = DTSurface3.toGlobal(segMB3->localPosition());
362  GlobalPoint segmentPositionMB4inGlobal = DTSurface4.toGlobal(segmentPosition);
363 
364  LocalVector segDirMB3inMB4Frame = DTSurface4.toLocal(DTSurface3.toGlobal(segmentDirectionMB3));
365 
366  GlobalVector segDirMB4inGlobalFrame = DTSurface4.toGlobal(segmentDirectionMB4);
367  GlobalVector segDirMB3inGlobalFrame = DTSurface3.toGlobal(segmentDirectionMB3);
368 
369  float dx = segDirMB4inGlobalFrame.x();
370  float dy = segDirMB4inGlobalFrame.y();
371 
372  float dx3 = segDirMB3inGlobalFrame.x();
373  float dy3 = segDirMB3inGlobalFrame.y();
374 
375  double cosAng = fabs(dx * dx3 + dy * dy3 / sqrt((dx3 * dx3 + dy3 * dy3) * (dx * dx + dy * dy)));
376 
377  if (debug)
378  LogDebug("DTSegtoRPC")
379  << "MB4 \t \t \t \t cosAng" << cosAng << "Beetween " << dtid3 << " and " << DTId << std::endl;
380 
381  if (debug) {
382  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t dx=" << dx << " dy=" << dy << std::endl;
383  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t dx3=" << dx3 << " dy3=" << dy << std::endl;
384  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t cosAng=" << cosAng << std::endl;
385  }
386 
387  float DistanceBetweenSegments = ((segmentPositionMB3inGlobal) - (segmentPositionMB4inGlobal)).mag();
388 
389  if (cosAng > MinCosAng && DistanceBetweenSegments < MaxDistanceBetweenSegments) {
390  if (debug)
391  LogDebug("DTSegtoRPC")
392  << "MB4 \t \t \t \t Distance between segments=" << DistanceBetweenSegments << std::endl;
393 
394  if (debug)
395  LogDebug("DTSegtoRPC")
396  << "MB4 \t \t We found compatible Segments (similar direction and close enough) in " << dtid3
397  << " and " << DTId << std::endl;
398 
399  if (dtSector == 13) {
400  dtSector = 4;
401  }
402  if (dtSector == 14) {
403  dtSector = 10;
404  }
405 
406  if (debug)
407  LogDebug("DTSegtoRPC") << "Creating the DTIndex" << std::endl;
408  DTStationIndex theindex(0, dtWheel, dtSector, dtStation);
409  if (debug)
410  LogDebug("DTSegtoRPC") << "Getting the Rolls for the given index" << std::endl;
411  std::set<RPCDetId> rollsForThisDT = dtMap->getRolls(theindex);
412 
413  if (debug)
414  LogDebug("DTSegtoRPC")
415  << "MB4 \t \t Number of rolls for this DT = " << rollsForThisDT.size() << std::endl;
416 
417  assert(!rollsForThisDT.empty());
418 
419  if (debug)
420  LogDebug("DTSegtoRPC") << "MB4 \t \t Loop over all the rolls asociated to this DT" << std::endl;
421  for (std::set<RPCDetId>::iterator iteraRoll = rollsForThisDT.begin();
422  iteraRoll != rollsForThisDT.end();
423  iteraRoll++) {
424  const RPCRoll* rollasociated = rpcGeo->roll(*iteraRoll); //roll asociado a MB4
425  RPCDetId rpcId = rollasociated->id();
426  const BoundPlane& RPCSurfaceRB4 = rollasociated->surface(); //surface MB4
427 
428  RPCGeomServ rpcsrv(rpcId);
429  std::string nameRoll = rpcsrv.name();
430 
431  if (debug)
432  LogDebug("DTSegtoRPC") << "MB4 \t \t \t RollName: " << nameRoll << std::endl;
433  if (debug)
434  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Doing the extrapolation to this roll" << std::endl;
435 
436  GlobalPoint CenterPointRollGlobal = RPCSurfaceRB4.toGlobal(LocalPoint(0, 0, 0));
437  LocalPoint CenterRollinMB4Frame = DTSurface4.toLocal(CenterPointRollGlobal); //In MB4
438  LocalPoint segmentPositionMB3inMB4Frame =
439  DTSurface4.toLocal(segmentPositionMB3inGlobal); //In MB4
440 
441  LocalVector segmentDirectionMB3inMB4Frame = DTSurface4.toLocal(segDirMB3inGlobalFrame); //In MB4
442 
443  //The exptrapolation is done in MB4 frame. for local x and z is done from MB4,
444  float Dxz = CenterRollinMB4Frame.z();
445  float Xo4 = segmentPositionMB4.x();
446  float dxl = segmentDirectionMB4.x(); //dx local for MB4 segment in MB4 Frame
447  float dzl = segmentDirectionMB4.z(); //dx local for MB4 segment in MB4 Frame
448 
449  float X = Xo4 + dxl * Dxz / dzl; //In MB4 frame
450  float Z = Dxz; //In MB4 frame
451 
452  //for local y is done from MB3
453  float Yo34 = segmentPositionMB3inMB4Frame.y();
454  float dy34 = segmentDirectionMB3inMB4Frame.y();
455  float dz34 = segmentDirectionMB3inMB4Frame.z();
456  float Dy =
457  Dxz -
458  (segmentPositionMB3inMB4Frame.z()); //Distance beetween the segment in MB3 and the RB4 surface
459 
460  if (debug)
461  LogDebug("DTSegtoRPC")
462  << "MB4 \t \t \t The distance to extrapolate in Y from MB3 is " << Dy << "cm" << std::endl;
463 
464  float Y = Yo34 + dy34 * Dy / dz34; //In MB4 Frame
465 
466  const RectangularStripTopology* top_ = dynamic_cast<const RectangularStripTopology*>(
467  &(rollasociated->topology())); //Topology roll asociated MB4
468  LocalPoint xmin = top_->localPosition(0.);
469  LocalPoint xmax = top_->localPosition((float)rollasociated->nstrips());
470  float rsize = fabs(xmax.x() - xmin.x());
471  float stripl = top_->stripLength();
472  float stripw = top_->pitch();
473 
474  if (debug)
475  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Strip Lenght " << stripl << "cm" << std::endl;
476  if (debug)
477  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Strip Width " << stripw << "cm" << std::endl;
478 
479  if (debug)
480  LogDebug("DTSegtoRPC") << "MB4 \t \t \t X Predicted in MB4DTLocal= " << X << "cm" << std::endl;
481  if (debug)
482  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Y Predicted in MB4DTLocal= " << Y << "cm" << std::endl;
483  if (debug)
484  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Z Predicted in MB4DTLocal= " << Z << "cm" << std::endl;
485 
486  float extrapolatedDistance = sqrt((Y - Yo34) * (Y - Yo34) + Dy * Dy);
487 
488  if (debug)
489  LogDebug("DTSegtoRPC")
490  << "MB4 \t \t \t segmentPositionMB3inMB4Frame" << segmentPositionMB3inMB4Frame << std::endl;
491  if (debug)
492  LogDebug("DTSegtoRPC")
493  << "MB4 \t \t \t segmentPositionMB4inMB4Frame" << segmentPosition << std::endl;
494 
495  if (debug)
496  LogDebug("DTSegtoRPC")
497  << "MB4 \t \t \t segmentDirMB3inMB4Frame" << segDirMB3inMB4Frame << std::endl;
498  if (debug)
499  LogDebug("DTSegtoRPC")
500  << "MB4 \t \t \t segmentDirMB4inMB4Frame" << segmentDirectionMB4 << std::endl;
501 
502  if (debug)
503  LogDebug("DTSegtoRPC")
504  << "MB4 \t \t \t CenterRB4PositioninMB4Frame" << CenterRollinMB4Frame << std::endl;
505 
506  if (debug)
507  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Is the extrapolation distance =" << extrapolatedDistance
508  << "less than " << MaxDrb4 << std::endl;
509 
510  if (extrapolatedDistance <= MaxDrb4) {
511  if (debug)
512  LogDebug("DTSegtoRPC") << "MB4 \t \t \t yes" << std::endl;
513 
514  GlobalPoint GlobalPointExtrapolated = DTSurface4.toGlobal(LocalPoint(X, Y, Z));
515 
516  if (debug)
517  LogDebug("DTSegtoRPC")
518  << "MB4 \t \t \t Point ExtraPolated in Global" << GlobalPointExtrapolated << std::endl;
519 
520  LocalPoint PointExtrapolatedRPCFrame = RPCSurfaceRB4.toLocal(GlobalPointExtrapolated);
521 
522  if (debug)
523  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Point Extrapolated in RPCLocal"
524  << PointExtrapolatedRPCFrame << std::endl;
525  if (debug)
526  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Corner of the Roll = (" << rsize * eyr << ","
527  << stripl * eyr << ")" << std::endl;
528  if (debug)
529  LogDebug("DTSegtoRPC")
530  << "MB4 \t \t \t Info About the Point Extrapolated in X Abs ("
531  << fabs(PointExtrapolatedRPCFrame.x()) << "," << fabs(PointExtrapolatedRPCFrame.y())
532  << "," << fabs(PointExtrapolatedRPCFrame.z()) << ")" << std::endl;
533 
534  if (debug)
535  LogDebug("DTSegtoRPC")
536  << "MB4 \t \t \t Does the extrapolation go inside this roll?" << std::endl;
537 
538  if (fabs(PointExtrapolatedRPCFrame.z()) < 5. &&
539  fabs(PointExtrapolatedRPCFrame.x()) < rsize * eyr &&
540  fabs(PointExtrapolatedRPCFrame.y()) < stripl * eyr) {
541  if (debug)
542  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t yes" << std::endl;
543  if (debug)
544  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t Creating the RecHit" << std::endl;
545 
546  RPCRecHit RPCPointMB4(rpcId, myPhiBx, PointExtrapolatedRPCFrame);
547  RPCPointMB4.setTimeAndError(myPhiTime, myPhiTimeErr);
548 
549  if (debug)
550  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t Clearing the RPCPointVector" << std::endl;
551  RPCPointVector.clear();
552  if (debug)
553  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t Pushing Back" << std::endl;
554  RPCPointVector.push_back(RPCPointMB4);
555  if (debug)
556  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t Putting for " << rpcId << std::endl;
557  if (debug)
558  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t Filling container with " << nameRoll
559  << " Point.x=" << PointExtrapolatedRPCFrame.x()
560  << " Point.y=" << PointExtrapolatedRPCFrame.y()
561  << " size=" << RPCPointVector.size() << std::endl;
562  if (debug)
563  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t Number of rolls already extrapolated in RB4 = "
564  << extrapolatedRolls.size() << std::endl;
565  if (find(extrapolatedRolls.begin(), extrapolatedRolls.end(), rpcId.rawId()) ==
566  extrapolatedRolls.end()) {
567  extrapolatedRolls.push_back(rpcId.rawId());
568  _ThePoints->put(rpcId, RPCPointVector.begin(), RPCPointVector.end());
569  } else {
570  if (debug)
571  LogDebug("DTSegtoRPC")
572  << "MB4 \t \t \t \t roll already extrapolated " << rpcId << std::endl;
573  }
574  if (debug)
575  LogDebug("DTSegtoRPC")
576  << "MB4 \t \t \t \t Extrapolations done after this point = " << extrapolatedRolls.size()
577  << std::endl;
578  if (debug)
579  for (uint32_t m = 0; m < extrapolatedRolls.size(); m++)
580  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t" << extrapolatedRolls.at(m) << std::endl;
581  } else {
582  if (debug)
583  LogDebug("DTSegtoRPC")
584  << "MB4 \t \t \t \t No the prediction is outside of this roll" << std::endl;
585  }
586  } //Condition for the right match
587  else {
588  if (debug)
589  LogDebug("DTSegtoRPC") << "MB4 \t \t \t No, Exrtrapolation too long!, canceled" << std::endl;
590  }
591  } //loop over all the rollsasociated
592  } else {
593  if (debug)
594  LogDebug("DTSegtoRPC")
595  << "MB4 \t \t \t \t I found segments in MB4 and MB3 adjacent wheel and/or sector but "
596  "not compatibles, Diferent Directions"
597  << std::endl;
598  }
599  } else {
600  if (debug)
601  LogDebug("DTSegtoRPC")
602  << "MB4 \t \t \t No the same station or same wheel or segment dim in mb3 not 4D" << std::endl;
603  }
604  } //loop over all the segments looking for one in MB3
605  } else {
606  if (debug)
607  LogDebug("DTSegtoRPC") << "MB4 \t \t \t Is NOT a 2D Segment" << std::endl;
608  }
609  } else {
610  if (debug)
611  LogDebug("DTSegtoRPC") << "MB4 \t \t \t \t There is not just one segment or is not in station 4"
612  << std::endl;
613  } //De aca para abajo esta en dtpart.inl
614  }
615  } else {
616  if (debug)
617  LogDebug("DTSegtoRPC") << "MB4 \t This event doesn't have 4D Segment" << std::endl;
618  }
619  }
620  }
621 
622  return _ThePoints;
623 }
int station() const
Return the station number.
Definition: DTChamberId.h:45
double MaxDistanceBetweenSegments
Definition: DTSegtoRPC.h:36
LocalPoint localPosition(float strip) const override
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
int minPhiBX
Definition: DTSegtoRPC.h:37
T z() const
Definition: PV3DBase.h:61
int distsector(int sector1, int sector2)
Definition: DTSegtoRPC.cc:18
#define X(str)
Definition: MuonsGrabber.cc:38
float pitch() const override
edm::ESGetToken< DTObjectMap, MuonGeometryRecord > dtMapToken_
Definition: DTSegtoRPC.h:29
const RPCRoll * roll(RPCDetId id) const
Return a roll given its id.
Definition: RPCGeometry.cc:50
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
int nstrips() const
Definition: RPCRoll.cc:24
assert(be >=bs)
std::set< RPCDetId > const & getRolls(DTStationIndex index) const
Definition: DTObjectMap.cc:36
float stripLength() const override
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
T x() const
Definition: PV3DBase.h:59
T y() const
Definition: PV3DBase.h:60
double MaxD
Definition: DTSegtoRPC.h:34
T sqrt(T t)
Definition: SSEVec.h:19
edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeoToken_
Definition: DTSegtoRPC.h:28
const GeomDet * idToDet(DetId) const override
Definition: DTGeometry.cc:77
double MaxDrb4
Definition: DTSegtoRPC.h:35
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
#define debug
Definition: HDRShower.cc:19
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
bool incldtMB4
Definition: DTSegtoRPC.h:32
bool incldt
Definition: DTSegtoRPC.h:31
int distwheel(int wheel1, int wheel2)
Definition: DTSegtoRPC.cc:35
edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeoToken_
Definition: DTSegtoRPC.h:27
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:42
int sector() const
Definition: DTChamberId.h:52
const Topology & topology() const override
Definition: RPCRoll.cc:18
RPCDetId id() const
Definition: RPCRoll.cc:16
double MinCosAng
Definition: DTSegtoRPC.h:33
#define LogDebug(id)

Member Data Documentation

◆ dtGeoToken_

edm::ESGetToken<DTGeometry, MuonGeometryRecord> DTSegtoRPC::dtGeoToken_
private

Definition at line 28 of file DTSegtoRPC.h.

Referenced by thePoints().

◆ dtMapToken_

edm::ESGetToken<DTObjectMap, MuonGeometryRecord> DTSegtoRPC::dtMapToken_
private

Definition at line 29 of file DTSegtoRPC.h.

Referenced by thePoints().

◆ incldt

bool DTSegtoRPC::incldt
private

Definition at line 31 of file DTSegtoRPC.h.

Referenced by DTSegtoRPC(), and thePoints().

◆ incldtMB4

bool DTSegtoRPC::incldtMB4
private

Definition at line 32 of file DTSegtoRPC.h.

Referenced by DTSegtoRPC(), and thePoints().

◆ MaxD

double DTSegtoRPC::MaxD
private

Definition at line 34 of file DTSegtoRPC.h.

Referenced by DTSegtoRPC(), and thePoints().

◆ MaxDistanceBetweenSegments

double DTSegtoRPC::MaxDistanceBetweenSegments
private

Definition at line 36 of file DTSegtoRPC.h.

Referenced by DTSegtoRPC(), and thePoints().

◆ MaxDrb4

double DTSegtoRPC::MaxDrb4
private

Definition at line 35 of file DTSegtoRPC.h.

Referenced by DTSegtoRPC(), and thePoints().

◆ maxPhiBX

int DTSegtoRPC::maxPhiBX
private

Definition at line 38 of file DTSegtoRPC.h.

Referenced by DTSegtoRPC().

◆ MinCosAng

double DTSegtoRPC::MinCosAng
private

Definition at line 33 of file DTSegtoRPC.h.

Referenced by DTSegtoRPC(), and thePoints().

◆ minPhiBX

int DTSegtoRPC::minPhiBX
private

Definition at line 37 of file DTSegtoRPC.h.

Referenced by DTSegtoRPC(), and thePoints().

◆ rpcGeoToken_

edm::ESGetToken<RPCGeometry, MuonGeometryRecord> DTSegtoRPC::rpcGeoToken_
private

Definition at line 27 of file DTSegtoRPC.h.

Referenced by thePoints().