CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
DTSegtoRPC.cc
Go to the documentation of this file.
14 
15 #include <ctime>
16 
17 int distsector(int sector1,int sector2){
18  if(sector1==13) sector1=4;
19  if(sector1==14) sector1=10;
20 
21  if(sector2==13) sector2=4;
22  if(sector2==14) sector2=10;
23 
24  int distance = std::abs(sector1 - sector2);
25  if(distance>6) distance = 12-distance;
26  return distance;
27 }
28 
29 int distwheel(int wheel1,int wheel2){
30  int distance = std::abs(wheel1 - wheel2);
31  return distance;
32 }
33 
34 DTSegtoRPC::DTSegtoRPC(const DTRecSegment4DCollection * all4DSegments, const edm::EventSetup& iSetup, bool debug,double eyr){
35 
36  /*
37  MinCosAng=iConfig.getUntrackedParameter<double>("MinCosAng",0.95);
38  MaxD=iConfig.getUntrackedParameter<double>("MaxD",80.);
39  MaxDrb4=iConfig.getUntrackedParameter<double>("MaxDrb4",150.);
40  */
41  incldt=true;
42  incldtMB4=true;
43 
44  //By now hard coded parameters
45  MinCosAng=0.85;
46  MaxD=80.;
47  MaxDrb4=150.;
49  /*
50 
51  //These should be always true expect for debuggin porpouses
52  incldt=true;
53  incldtMB4=true;
54 
55 
56  struct timespec start_time, stop_time;
57  time_t fs;
58  time_t fn;
59  time_t ls;
60  time_t ln;
61  clock_gettime(CLOCK_REALTIME, &start_time);
62  */
63 
64  _ThePoints.reset(new RPCRecHitCollection());
65 
66  if(all4DSegments->size()>8){
67  if(debug) std::cout<<"Too many segments in this event we are not doing the extrapolation"<<std::endl;
68  }else{
72 
73  iSetup.get<MuonGeometryRecord>().get(rpcGeo);
74  iSetup.get<MuonGeometryRecord>().get(dtGeo);
75  iSetup.get<MuonGeometryRecord>().get(dtMap);
76 
77  /*
78  clock_gettime(CLOCK_REALTIME, &stop_time);
79  fs=start_time.tv_sec;
80  fn=start_time.tv_nsec;
81  ls=stop_time.tv_sec;
82  ln=stop_time.tv_nsec;
83  std::cout <<" =================|| "<<ls-fs<<" sec "<<ln-fn<<" us"<<std::endl;
84  clock_gettime(CLOCK_REALTIME, &start_time);
85  */
86 
87  std::map<DTChamberId,int> DTSegmentCounter;
89 
90  for (segment = all4DSegments->begin();segment!=all4DSegments->end(); ++segment){
91  DTSegmentCounter[segment->chamberId()]++;
92  }
93 
94 
95  /*
96  clock_gettime(CLOCK_REALTIME, &stop_time);
97  fs=start_time.tv_sec;
98  fn=start_time.tv_nsec;
99  ls=stop_time.tv_sec;
100  ln=stop_time.tv_nsec;
101  if(debug) std::cout <<" =================||| "<<ls-fs<<" sec "<<ln-fn<<" us"<<std::endl;
102  clock_gettime(CLOCK_REALTIME, &start_time);
103  */
104 
105  if(incldt){
106  for (segment = all4DSegments->begin(); segment != all4DSegments->end(); ++segment){
107 
108  DTChamberId DTId = segment->chamberId();
109 
110  if(debug) std::cout<<"DT \t \t This Segment is in Chamber id: "<<DTId<<std::endl;
111  if(debug) std::cout<<"DT \t \t Number of segments in this DT = "<<DTSegmentCounter[DTId]<<std::endl;
112  if(debug) std::cout<<"DT \t \t Is the only one in this DT? and is not in the 4th Station?"<<std::endl;
113 
114  if(DTSegmentCounter[DTId]!=1 || DTId.station()==4){
115  if(debug) std::cout<<"DT \t \t More than one segment in this chamber, or we are in Station 4"<<std::endl;
116  continue;
117  }
118 
119  int dtWheel = DTId.wheel();
120  int dtStation = DTId.station();
121  int dtSector = DTId.sector();
122 
123  LocalPoint segmentPosition= segment->localPosition();
124  LocalVector segmentDirection=segment->localDirection();
125 
126  const GeomDet* gdet=dtGeo->idToDet(segment->geographicalId());
127  const BoundPlane & DTSurface = gdet->surface();
128 
129  //check if the dimension of the segment is 4
130 
131  if(debug) std::cout<<"DT \t \t Is the segment 4D?"<<std::endl;
132 
133  if(segment->dimension()!=4){
134  if(debug) std::cout<<"DT \t \t no"<<std::endl;
135  continue;
136  }
137 
138  if(debug) std::cout<<"DT \t \t yes"<<std::endl;
139  if(debug) std::cout<<"DT \t \t DT Segment Dimension "<<segment->dimension()<<std::endl;
140 
141  float Xo=segmentPosition.x();
142  float Yo=segmentPosition.y();
143  float Zo=segmentPosition.z();
144  float dx=segmentDirection.x();
145  float dy=segmentDirection.y();
146  float dz=segmentDirection.z();
147 
148  if(debug) std::cout<<"Creating the DTIndex"<<std::endl;
149  DTStationIndex theindex(0,dtWheel,dtSector,dtStation);
150  if(debug) std::cout<<"Getting the Rolls for the given index"<<std::endl;
151  std::set<RPCDetId> rollsForThisDT = dtMap->getRolls(theindex);
152 
153  if(debug) std::cout<<"DT \t \t Number of rolls for this DT = "<<rollsForThisDT.size()<<std::endl;
154 
155  assert(rollsForThisDT.size()>=1);
156 
157  if(debug) std::cout<<"DT \t \t Loop over all the rolls asociated to this DT"<<std::endl;
158  for (std::set<RPCDetId>::iterator iteraRoll = rollsForThisDT.begin();iteraRoll != rollsForThisDT.end(); 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) std::cout<<"DT \t \t \t RollName: "<<nameRoll<<std::endl;
167  if(debug) std::cout<<"DT \t \t \t Doing the extrapolation to this roll"<<std::endl;
168  if(debug) std::cout<<"DT \t \t \t DT Segment Direction in DTLocal "<<segmentDirection<<std::endl;
169  if(debug) std::cout<<"DT \t \t \t DT Segment Point in DTLocal "<<segmentPosition<<std::endl;
170 
171  GlobalPoint CenterPointRollGlobal = RPCSurface.toGlobal(LocalPoint(0,0,0));
172 
173  LocalPoint CenterRollinDTFrame = DTSurface.toLocal(CenterPointRollGlobal);
174 
175  if(debug) std::cout<<"DT \t \t \t Center (0,0,0) Roll In DTLocal"<<CenterRollinDTFrame<<std::endl;
176  if(debug) std::cout<<"DT \t \t \t Center (0,0,0) of the Roll in Global"<<CenterPointRollGlobal<<std::endl;
177 
178  float D=CenterRollinDTFrame.z();
179 
180  float X=Xo+dx*D/dz;
181  float Y=Yo+dy*D/dz;
182  float Z=D;
183 
184  const RectangularStripTopology* top_= dynamic_cast<const RectangularStripTopology*> (&(rollasociated->topology()));
185  LocalPoint xmin = top_->localPosition(0.);
186  if(debug) std::cout<<"DT \t \t \t xmin of this Roll "<<xmin<<"cm"<<std::endl;
187  LocalPoint xmax = top_->localPosition((float)rollasociated->nstrips());
188  if(debug) std::cout<<"DT \t \t \t xmax of this Roll "<<xmax<<"cm"<<std::endl;
189  float rsize = fabs( xmax.x()-xmin.x() );
190  if(debug) std::cout<<"DT \t \t \t Roll Size "<<rsize<<"cm"<<std::endl;
191  float stripl = top_->stripLength();
192 
193  float stripw = top_->pitch();
194 
195  if(debug) std::cout<<"DT \t \t \t Strip Lenght "<<stripl<<"cm"<<std::endl;
196  if(debug) std::cout<<"DT \t \t \t Strip Width "<<stripw<<"cm"<<std::endl;
197  if(debug) std::cout<<"DT \t \t \t X Predicted in DTLocal= "<<X<<"cm"<<std::endl;
198  if(debug) std::cout<<"DT \t \t \t Y Predicted in DTLocal= "<<Y<<"cm"<<std::endl;
199  if(debug) std::cout<<"DT \t \t \t Z Predicted in DTLocal= "<<Z<<"cm"<<std::endl;
200 
201  float extrapolatedDistance = sqrt((X-Xo)*(X-Xo)+(Y-Yo)*(Y-Yo)+(Z-Zo)*(Z-Zo));
202 
203  if(debug) std::cout<<"DT \t \t \t Is the distance of extrapolation less than MaxD? ="<<extrapolatedDistance<<"cm"<<"MaxD="<<MaxD<<"cm"<<std::endl;
204 
205  if(extrapolatedDistance<=MaxD){
206  if(debug) std::cout<<"DT \t \t \t yes"<<std::endl;
207  GlobalPoint GlobalPointExtrapolated = DTSurface.toGlobal(LocalPoint(X,Y,Z));
208  if(debug) std::cout<<"DT \t \t \t Point ExtraPolated in Global"<<GlobalPointExtrapolated<< std::endl;
209  LocalPoint PointExtrapolatedRPCFrame = RPCSurface.toLocal(GlobalPointExtrapolated);
210 
211  if(debug) std::cout<<"DT \t \t \t Point Extrapolated in RPCLocal"<<PointExtrapolatedRPCFrame<< std::endl;
212  if(debug) std::cout<<"DT \t \t \t Corner of the Roll = ("<<rsize*eyr<<","<<stripl*eyr<<")"<<std::endl;
213  if(debug) std::cout<<"DT \t \t \t Info About the Point Extrapolated in X Abs ("<<fabs(PointExtrapolatedRPCFrame.x())<<","
214  <<fabs(PointExtrapolatedRPCFrame.y())<<","<<fabs(PointExtrapolatedRPCFrame.z())<<")"<<std::endl;
215  if(debug) std::cout<<"DT \t \t \t Does the extrapolation go inside this roll?"<<std::endl;
216 
217  if(fabs(PointExtrapolatedRPCFrame.z()) < 1. &&
218  fabs(PointExtrapolatedRPCFrame.x()) < rsize*eyr &&
219  fabs(PointExtrapolatedRPCFrame.y()) < stripl*eyr){
220  if(debug) std::cout<<"DT \t \t \t \t yes"<<std::endl;
221  if(debug) std::cout<<"DT \t \t \t \t Creating the RecHit"<<std::endl;
222  RPCRecHit RPCPoint(rpcId,0,PointExtrapolatedRPCFrame);
223  if(debug) std::cout<<"DT \t \t \t \t Clearing the vector"<<std::endl;
225  if(debug) std::cout<<"DT \t \t \t \t Pushing back"<<std::endl;
226  RPCPointVector.push_back(RPCPoint);
227  if(debug) std::cout<<"DT \t \t \t \t Putting the vector"<<std::endl;
229 
230  if(debug) std::cout<<"DT \t \t \t \t Filling container with "<<nameRoll
231  <<" Point.x="<<PointExtrapolatedRPCFrame.x()<<" Point.y="<<PointExtrapolatedRPCFrame.y()<<" size="<<RPCPointVector.size()<<std::endl;
232 
233  }else {
234  if(debug) std::cout<<"DT \t \t \t \t No the prediction is outside of this roll"<<std::endl;
235  }//Condition for the right match
236  }else{
237  if(debug) std::cout<<"DT \t \t \t No, Exrtrapolation too long!, canceled"<<std::endl;
238  }//D so big
239  }//loop over all the rolls asociated
240  }
241  }
242 
243  if(incldtMB4){
244  if(all4DSegments->size()>0){
245  if(debug) std::cout<<"MB4 \t \t Loop Over all4DSegments "<<all4DSegments->size()<<std::endl;
246  extrapolatedRolls.clear();
247  for (segment = all4DSegments->begin(); segment != all4DSegments->end(); ++segment){
248 
249  DTChamberId DTId = segment->chamberId();
250 
251  if(debug) std::cout<<"MB4 \t \t This Segment is in Chamber id: "<<DTId<<std::endl;
252  if(debug) std::cout<<"MB4 \t \t Number of segments in this DT = "<<DTSegmentCounter[DTId]<<std::endl;
253  if(debug) std::cout<<"MB4 \t \t \t Is the only one in this DT? and is in the Station 4?"<<std::endl;
254 
255  if(DTSegmentCounter[DTId] == 1 && DTId.station()==4){
256 
257  if(debug) std::cout<<"MB4 \t \t \t yes"<<std::endl;
258  int dtWheel = DTId.wheel();
259  int dtStation = DTId.station();
260  int dtSector = DTId.sector();
261 
262  LocalPoint segmentPosition= segment->localPosition();
263  LocalVector segmentDirection=segment->localDirection();
264 
265  if(debug) std::cout<<"MB4 \t \t \t \t The Segment in MB4 is 2D?"<<std::endl;
266  if(segment->dimension()==2){
267  if(debug) std::cout<<"MB4 \t \t \t \t yes"<<std::endl;
268  LocalVector segmentDirectionMB4=segmentDirection;
269  LocalPoint segmentPositionMB4=segmentPosition;
270 
271  const BoundPlane& DTSurface4 = dtGeo->idToDet(DTId)->surface();
272 
274 
275  if(debug) std::cout<<"MB4 \t \t \t \t Loop on segments in =sector && MB3 && adjacent sectors && y dim=4"<<std::endl;
276  for(segMB3=all4DSegments->begin();segMB3!=all4DSegments->end();++segMB3){
277 
278  DTChamberId dtid3 = segMB3->chamberId();
279 
280  if(debug) std::cout<<"MB4 \t \t \t \t Segment in Chamber ="<<dtid3<<std::endl;
281 
282 
283  if(distsector(dtid3.sector(),DTId.sector())<=1 //The DT sector could be 13 or 14 and because is corrected in the calculation of the distance.
284  && distwheel(dtid3.wheel(),DTId.wheel())<=1 //The we could have segments in neighbohr wheels in pp collisions
285  && dtid3.station()==3
286  && DTSegmentCounter[dtid3] == 1
287  && segMB3->dimension()==4){
288 
289  if(debug) std::cout<<"MB4 \t \t \t \t distsector ="<<distsector(dtid3.sector(),DTId.sector())<<std::endl;
290  if(debug) std::cout<<"MB4 \t \t \t \t distwheel ="<<distwheel(dtid3.wheel(),DTId.wheel())<<std::endl;
291 
292  const GeomDet* gdet3=dtGeo->idToDet(segMB3->geographicalId());
293  const BoundPlane & DTSurface3 = gdet3->surface();
294 
295  LocalVector segmentDirectionMB3 = segMB3->localDirection();
296  GlobalPoint segmentPositionMB3inGlobal = DTSurface3.toGlobal(segMB3->localPosition());
297  GlobalPoint segmentPositionMB4inGlobal = DTSurface4.toGlobal(segmentPosition);
298 
299  //LocalVector segDirMB4inMB3Frame=DTSurface3.toLocal(DTSurface4.toGlobal(segmentDirectionMB4));
300  LocalVector segDirMB3inMB4Frame=DTSurface4.toLocal(DTSurface3.toGlobal(segmentDirectionMB3));
301 
302  GlobalVector segDirMB4inGlobalFrame=DTSurface4.toGlobal(segmentDirectionMB4);
303  GlobalVector segDirMB3inGlobalFrame=DTSurface3.toGlobal(segmentDirectionMB3);
304 
305  float dx=segDirMB4inGlobalFrame.x();
306  float dy=segDirMB4inGlobalFrame.y();
307 
308  float dx3=segDirMB3inGlobalFrame.x();
309  float dy3=segDirMB3inGlobalFrame.y();
310 
311  double cosAng=fabs(dx*dx3+dy*dy3/sqrt((dx3*dx3+dy3*dy3)*(dx*dx+dy*dy)));
312 
313  if(debug) std::cout<<"MB4 \t \t \t \t cosAng"<<cosAng<<"Beetween "<<dtid3<<" and "<<DTId<<std::endl;
314 
315  if(debug){
316  std::cout<<"MB4 \t \t \t \t dx="<<dx<<" dy="<<dy<<std::endl;
317  std::cout<<"MB4 \t \t \t \t dx3="<<dx3<<" dy3="<<dy<<std::endl;
318  std::cout<<"MB4 \t \t \t \t cosAng="<<cosAng<<std::endl;
319  }
320 
321  float DistanceBetweenSegments = ((segmentPositionMB3inGlobal) - (segmentPositionMB4inGlobal)).mag();
322 
323  if(cosAng>MinCosAng && DistanceBetweenSegments < MaxDistanceBetweenSegments){
324 
325  if(debug) std::cout<<"MB4 \t \t \t \t Distance between segments="<<DistanceBetweenSegments<<std::endl;
326 
327  if(debug) std::cout<<"MB4 \t \t We found compatible Segments (similar direction and close enough) in "<<dtid3<<" and "<<DTId<<std::endl;
328 
329  if(dtSector==13){
330  dtSector=4;
331  }
332  if(dtSector==14){
333  dtSector=10;
334  }
335 
336  if(debug) std::cout<<"Creating the DTIndex"<<std::endl;
337  DTStationIndex theindex(0,dtWheel,dtSector,dtStation);
338  if(debug) std::cout<<"Getting the Rolls for the given index"<<std::endl;
339  std::set<RPCDetId> rollsForThisDT = dtMap->getRolls(theindex);
340 
341  if(debug) std::cout<<"MB4 \t \t Number of rolls for this DT = "<<rollsForThisDT.size()<<std::endl;
342 
343  assert(rollsForThisDT.size()>=1);
344 
345  if(debug) std::cout<<"MB4 \t \t Loop over all the rolls asociated to this DT"<<std::endl;
346  for (std::set<RPCDetId>::iterator iteraRoll=rollsForThisDT.begin();iteraRoll != rollsForThisDT.end(); iteraRoll++){
347  const RPCRoll* rollasociated = rpcGeo->roll(*iteraRoll); //roll asociado a MB4
348  RPCDetId rpcId = rollasociated->id();
349  const BoundPlane & RPCSurfaceRB4 = rollasociated->surface(); //surface MB4
350 
351  RPCGeomServ rpcsrv(rpcId);
352  std::string nameRoll = rpcsrv.name();
353 
354  if(debug) std::cout<<"MB4 \t \t \t RollName: "<<nameRoll<<std::endl;
355  if(debug) std::cout<<"MB4 \t \t \t Doing the extrapolation to this roll"<<std::endl;
356 
357  GlobalPoint CenterPointRollGlobal=RPCSurfaceRB4.toGlobal(LocalPoint(0,0,0));
358  LocalPoint CenterRollinMB4Frame = DTSurface4.toLocal(CenterPointRollGlobal); //In MB4
359  LocalPoint segmentPositionMB3inMB4Frame = DTSurface4.toLocal(segmentPositionMB3inGlobal); //In MB4
360  //LocalPoint segmentPositionMB3inRB4Frame = RPCSurfaceRB4.toLocal(segmentPositionMB3inGlobal); //In MB4
361  LocalVector segmentDirectionMB3inMB4Frame = DTSurface4.toLocal(segDirMB3inGlobalFrame); //In MB4
362 
363  //The exptrapolation is done in MB4 frame. for local x and z is done from MB4,
364  float Dxz=CenterRollinMB4Frame.z();
365  float Xo4=segmentPositionMB4.x();
366  float dxl=segmentDirectionMB4.x(); //dx local for MB4 segment in MB4 Frame
367  float dzl=segmentDirectionMB4.z(); //dx local for MB4 segment in MB4 Frame
368 
369  float X=Xo4+dxl*Dxz/dzl; //In MB4 frame
370  float Z=Dxz;//In MB4 frame
371 
372  //for local y is done from MB3
373  float Yo34 = segmentPositionMB3inMB4Frame.y();
374  float dy34 = segmentDirectionMB3inMB4Frame.y();
375  float dz34 = segmentDirectionMB3inMB4Frame.z();
376  float Dy=Dxz-(segmentPositionMB3inMB4Frame.z()); //Distance beetween the segment in MB3 and the RB4 surface
377 
378  if(debug) std::cout<<"MB4 \t \t \t The distance to extrapolate in Y from MB3 is "<<Dy<<"cm"<<std::endl;
379 
380  float Y=Yo34+dy34*Dy/dz34;//In MB4 Frame
381 
382  const RectangularStripTopology* top_
383  =dynamic_cast<const RectangularStripTopology*>(&(rollasociated->topology())); //Topology roll asociated MB4
384  LocalPoint xmin = top_->localPosition(0.);
385  LocalPoint xmax = top_->localPosition((float)rollasociated->nstrips());
386  float rsize = fabs( xmax.x()-xmin.x() );
387  float stripl = top_->stripLength();
388  float stripw = top_->pitch();
389 
390 
391  if(debug) std::cout<<"MB4 \t \t \t Strip Lenght "<<stripl<<"cm"<<std::endl;
392  if(debug) std::cout<<"MB4 \t \t \t Strip Width "<<stripw<<"cm"<<std::endl;
393 
394  if(debug) std::cout<<"MB4 \t \t \t X Predicted in MB4DTLocal= "<<X<<"cm"<<std::endl;
395  if(debug) std::cout<<"MB4 \t \t \t Y Predicted in MB4DTLocal= "<<Y<<"cm"<<std::endl;
396  if(debug) std::cout<<"MB4 \t \t \t Z Predicted in MB4DTLocal= "<<Z<<"cm"<<std::endl;
397 
398  float extrapolatedDistance = sqrt((Y-Yo34)*(Y-Yo34)+Dy*Dy);
399 
400  if(debug) std::cout<<"MB4 \t \t \t segmentPositionMB3inMB4Frame"<<segmentPositionMB3inMB4Frame<<std::endl;
401  if(debug) std::cout<<"MB4 \t \t \t segmentPositionMB4inMB4Frame"<<segmentPosition<<std::endl;
402 
403  if(debug) std::cout<<"MB4 \t \t \t segmentDirMB3inMB4Frame"<<segDirMB3inMB4Frame<<std::endl;
404  if(debug) std::cout<<"MB4 \t \t \t segmentDirMB4inMB4Frame"<<segmentDirectionMB4<<std::endl;
405 
406  if(debug) std::cout<<"MB4 \t \t \t CenterRB4PositioninMB4Frame"<<CenterRollinMB4Frame<<std::endl;
407 
408  if(debug) std::cout<<"MB4 \t \t \t Is the extrapolation distance ="<<extrapolatedDistance<<"less than "<<MaxDrb4<<std::endl;
409 
410 
411  if(extrapolatedDistance<=MaxDrb4){
412  if(debug) std::cout<<"MB4 \t \t \t yes"<<std::endl;
413 
414  GlobalPoint GlobalPointExtrapolated = DTSurface4.toGlobal(LocalPoint(X,Y,Z));
415 
416  if(debug) std::cout<<"MB4 \t \t \t Point ExtraPolated in Global"<<GlobalPointExtrapolated<< std::endl;
417 
418  LocalPoint PointExtrapolatedRPCFrame = RPCSurfaceRB4.toLocal(GlobalPointExtrapolated);
419 
420  if(debug) std::cout<<"MB4 \t \t \t Point Extrapolated in RPCLocal"<<PointExtrapolatedRPCFrame<< std::endl;
421  if(debug) std::cout<<"MB4 \t \t \t Corner of the Roll = ("<<rsize*eyr<<","<<stripl*eyr<<")"<<std::endl;
422  if(debug) std::cout<<"MB4 \t \t \t Info About the Point Extrapolated in X Abs ("<<fabs(PointExtrapolatedRPCFrame.x())<<","
423  <<fabs(PointExtrapolatedRPCFrame.y())<<","<<fabs(PointExtrapolatedRPCFrame.z())<<")"<<std::endl;
424 
425  if(debug) std::cout<<"MB4 \t \t \t Does the extrapolation go inside this roll?"<<std::endl;
426 
427  if(fabs(PointExtrapolatedRPCFrame.z()) < 5. &&
428  fabs(PointExtrapolatedRPCFrame.x()) < rsize*eyr &&
429  fabs(PointExtrapolatedRPCFrame.y()) < stripl*eyr){
430  if(debug) std::cout<<"MB4 \t \t \t \t yes"<<std::endl;
431  if(debug) std::cout<<"MB4 \t \t \t \t Creating the RecHit"<<std::endl;
432  RPCRecHit RPCPointMB4(rpcId,0,PointExtrapolatedRPCFrame);
433  if(debug) std::cout<<"MB4 \t \t \t \t Clearing the RPCPointVector"<<std::endl;
435  if(debug) std::cout<<"MB4 \t \t \t \t Pushing Back"<<std::endl;
436  RPCPointVector.push_back(RPCPointMB4);
437  if(debug) std::cout<<"MB4 \t \t \t \t Putting for "<<rpcId<<std::endl;
438  if(debug) std::cout<<"MB4 \t \t \t \t Filling container with "<<nameRoll
439  <<" Point.x="<<PointExtrapolatedRPCFrame.x()<<" Point.y="<<PointExtrapolatedRPCFrame.y()<<" size="<<RPCPointVector.size()<<std::endl;
440  if(debug) std::cout<<"MB4 \t \t \t \t Number of rolls already extrapolated in RB4 = "<<extrapolatedRolls.size()<<std::endl;
441  if(find (extrapolatedRolls.begin(),extrapolatedRolls.end(),rpcId.rawId()) == extrapolatedRolls.end()){
442  extrapolatedRolls.push_back(rpcId.rawId());
444  }else{
445  if(debug) std::cout<<"MB4 \t \t \t \t roll already extrapolated "<<rpcId<<std::endl;
446  }
447  if(debug) std::cout<<"MB4 \t \t \t \t Extrapolations done after this point = "<<extrapolatedRolls.size()<<std::endl;
448  if(debug) for(uint32_t m=0;m<extrapolatedRolls.size();m++) std::cout<<"MB4 \t \t \t \t"<< extrapolatedRolls.at(m)<<std::endl;
449  }else{
450  if(debug) std::cout<<"MB4 \t \t \t \t No the prediction is outside of this roll"<<std::endl;
451  }
452  }//Condition for the right match
453  else{
454  if(debug) std::cout<<"MB4 \t \t \t No, Exrtrapolation too long!, canceled"<<std::endl;
455  }
456  }//loop over all the rollsasociated
457  }else{
458  if(debug) std::cout<<"MB4 \t \t \t \t I found segments in MB4 and MB3 adjacent wheel and/or sector but not compatibles, Diferent Directions"<<std::endl;
459  }
460  }else{//if dtid3.station()==3&&dtid3.sector()==DTId.sector()&&dtid3.wheel()==DTId.wheel()&&segMB3->dim()==4
461  if(debug) std::cout<<"MB4 \t \t \t No the same station or same wheel or segment dim in mb3 not 4D"<<std::endl;
462  }
463  }//loop over all the segments looking for one in MB3
464  }else{
465  if(debug) std::cout<<"MB4 \t \t \t Is NOT a 2D Segment"<<std::endl;
466  }
467  }else{
468  if(debug) std::cout<<"MB4 \t \t \t \t There is not just one segment or is not in station 4"<<std::endl;
469  }//De aca para abajo esta en dtpart.inl
470  }
471  }else{
472  if(debug) std::cout<<"MB4 \t This event doesn't have 4D Segment"<<std::endl;
473  }
474  }
475  }
476 
477 
478 
479 
480  /*
481  clock_gettime(CLOCK_REALTIME, &stop_time);
482  fs=start_time.tv_sec;
483  fn=start_time.tv_nsec;
484  ls=stop_time.tv_sec;
485  ln=stop_time.tv_nsec;
486  std::cout <<" =================|||| "<<ls-fs<<" sec "<<ln-fn<<" us"<<std::endl;
487  */
488 }
489 
490 
491 
493 }
const double Z[kNumberCalorimeter]
double MaxDistanceBetweenSegments
Definition: DTSegtoRPC.h:26
const Topology & topology() const
Definition: RPCRoll.cc:30
virtual float stripLength() const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
int distsector(int sector1, int sector2)
Definition: DTSegtoRPC.cc:17
assert(m_qm.get())
int nstrips() const
Definition: RPCRoll.cc:46
size_type size() const
Definition: OwnVector.h:254
T y() const
Definition: PV3DBase.h:63
#define X(str)
Definition: MuonsGrabber.cc:48
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
edm::OwnVector< RPCRecHit > RPCPointVector
Definition: DTSegtoRPC.h:20
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:40
iterator begin()
Definition: OwnVector.h:234
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
void push_back(D *&d)
Definition: OwnVector.h:280
RPCDetId id() const
Definition: RPCRoll.cc:24
virtual std::string name()
Definition: RPCGeomServ.cc:20
double MaxD
Definition: DTSegtoRPC.h:24
T sqrt(T t)
Definition: SSEVec.h:48
T z() const
Definition: PV3DBase.h:64
void clear()
Definition: OwnVector.h:377
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< uint32_t > extrapolatedRolls
Definition: DTSegtoRPC.h:27
std::unique_ptr< RPCRecHitCollection > _ThePoints
Definition: DTSegtoRPC.h:19
double MaxDrb4
Definition: DTSegtoRPC.h:25
iterator end()
Definition: OwnVector.h:239
edm::RangeMap< RPCDetId, edm::OwnVector< RPCRecHit, edm::ClonePolicy< RPCRecHit > >, edm::ClonePolicy< RPCRecHit > > RPCRecHitCollection
#define debug
Definition: HDRShower.cc:19
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:150
const T & get() const
Definition: EventSetup.h:55
bool incldtMB4
Definition: DTSegtoRPC.h:22
bool incldt
Definition: DTSegtoRPC.h:21
int distwheel(int wheel1, int wheel2)
Definition: DTSegtoRPC.cc:29
virtual float pitch() const
int sector() const
Definition: DTChamberId.h:61
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
tuple cout
Definition: gather_cfg.py:121
double MinCosAng
Definition: DTSegtoRPC.h:23
virtual LocalPoint localPosition(float strip) const
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
T x() const
Definition: PV3DBase.h:62
DTSegtoRPC(DTRecSegment4DCollection const *all4DSegments, edm::EventSetup const &iSetup, bool debug, double eyr)
Definition: DTSegtoRPC.cc:34