CMS 3D CMS Logo

CSCSegtoRPC.cc
Go to the documentation of this file.
11 #include "CSCSegtoRPC.h"
12 #include "CSCStationIndex.h"
13 #include "CSCObjectMap.h"
15 
17  : rpcGeoToken_(iC.esConsumes()), cscGeoToken_(iC.esConsumes()), cscMapToken_(iC.esConsumes()) {
18  minBX = iConfig.getParameter<int>("minBX");
19  maxBX = iConfig.getParameter<int>("maxBX");
20 }
21 
22 std::unique_ptr<RPCRecHitCollection> CSCSegtoRPC::thePoints(const CSCSegmentCollection* allCSCSegments,
23  const edm::EventSetup& iSetup,
24  bool debug,
25  double eyr) {
29 
30  double MaxD = 80.;
31 
32  if (debug)
33  LogDebug("CSCSegtoRPC") << "CSC \t Number of CSC Segments in this event = " << allCSCSegments->size() << std::endl;
34 
35  auto _ThePoints = std::make_unique<RPCRecHitCollection>();
36  edm::OwnVector<RPCRecHit> RPCPointVector;
37 
38  if (allCSCSegments->size() == 0) {
39  if (debug)
40  LogDebug("CSCSegtoRPC") << "CSC 0 segments skiping event" << std::endl;
41  } else {
42  std::map<CSCDetId, int> CSCSegmentsCounter;
44 
45  int segmentsInThisEventInTheEndcap = 0;
46 
47  for (segment = allCSCSegments->begin(); segment != allCSCSegments->end(); ++segment) {
48  CSCSegmentsCounter[segment->cscDetId()]++;
49  segmentsInThisEventInTheEndcap++;
50  }
51 
52  float myTime = -9999.;
53  float myTimeErr = -9999.;
54  int myBx = -99;
55  const float bunchCrossTimeDiff = 25.; // time between bunch crossings
56 
57  if (debug)
58  LogDebug("CSCSegtoRPC") << "CSC \t loop over all the CSCSegments " << std::endl;
59  for (segment = allCSCSegments->begin(); segment != allCSCSegments->end(); ++segment) {
60  CSCDetId CSCId = segment->cscDetId();
61 
62  myTime = segment->time();
63  myBx = round(myTime / bunchCrossTimeDiff);
64  if (!(myBx <= maxBX && myBx >= minBX)) // rpc read in bx in [-2, 2]
65  continue;
66 
67  if (debug)
68  LogDebug("CSCSegtoRPC") << "CSC \t \t This Segment is in Chamber id: " << CSCId << std::endl;
69  if (debug)
70  LogDebug("CSCSegtoRPC") << "CSC \t \t Number of segments in this CSC = " << CSCSegmentsCounter[CSCId]
71  << std::endl;
72  if (debug)
73  LogDebug("CSCSegtoRPC")
74  << "CSC \t \t Is the only one in this CSC? is not ind the ring 1? Are there more than 2 "
75  "segments in the event?"
76  << std::endl;
77 
78  if (CSCSegmentsCounter[CSCId] == 1 && CSCId.ring() != 1 && allCSCSegments->size() >= 2) {
79  if (debug)
80  LogDebug("CSCSegtoRPC") << "CSC \t \t yes" << std::endl;
81  int cscEndCap = CSCId.endcap();
82  int cscStation = CSCId.station();
83  int cscRing = CSCId.ring();
84  int rpcRegion = 1;
85  if (cscEndCap == 2)
86  rpcRegion = -1; //Relacion entre las endcaps
87  int rpcRing = cscRing;
88  if (cscRing == 4)
89  rpcRing = 1;
90  int rpcStation = cscStation;
91  int rpcSegment = CSCId.chamber();
92 
93  LocalPoint segmentPosition = segment->localPosition();
94  LocalVector segmentDirection = segment->localDirection();
95  float dz = segmentDirection.z();
96 
97  if (debug)
98  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Information about the segment"
99  << "RecHits =" << segment->nRecHits() << "Angle =" << acos(dz) * 180 / 3.1415926
100  << std::endl;
101 
102  if (debug)
103  LogDebug("CSCSegtoRPC")
104  << "CSC \t \t Is a good Segment? dim = 4, 4 <= nRecHits <= 10 Incident angle int range 45 < "
105  << acos(dz) * 180 / 3.1415926 << " < 135? " << std::endl;
106 
107  if ((segment->dimension() == 4) && (segment->nRecHits() <= 10 && segment->nRecHits() >= 4)) {
108  if (debug)
109  LogDebug("CSCSegtoRPC") << "CSC \t \t yes" << std::endl;
110  if (debug)
111  LogDebug("CSCSegtoRPC") << "CSC \t \t CSC Segment Dimension " << segment->dimension() << std::endl;
112 
113  float Xo = segmentPosition.x();
114  float Yo = segmentPosition.y();
115  float Zo = segmentPosition.z();
116  float dx = segmentDirection.x();
117  float dy = segmentDirection.y();
118  float dz = segmentDirection.z();
119 
120  if (debug)
121  LogDebug("CSCSegtoRPC") << "Creating the CSCIndex" << std::endl;
122  CSCStationIndex theindex(rpcRegion, rpcStation, rpcRing, rpcSegment);
123  if (debug)
124  LogDebug("CSCSegtoRPC") << "Getting the Rolls for the given index" << std::endl;
125  std::set<RPCDetId> rollsForThisCSC = cscMap->getRolls(theindex);
126  if (debug)
127  LogDebug("CSCSegtoRPC") << "CSC \t \t Getting chamber from Geometry" << std::endl;
128  const CSCChamber* TheChamber = cscGeo->chamber(CSCId);
129  if (debug)
130  LogDebug("CSCSegtoRPC") << "CSC \t \t Getting ID from Chamber" << std::endl;
131  const CSCDetId TheId = TheChamber->id();
132 
133  if (debug)
134  LogDebug("CSCSegtoRPC") << "CSC \t \t Number of rolls for this CSC = " << rollsForThisCSC.size()
135  << std::endl;
136 
137  if (debug)
138  LogDebug("CSCSegtoRPC") << "CSC \t \t Printing The Id" << TheId << std::endl;
139 
140  if (rpcRing != 1) { //They don't exist in Run3!
141 
142  assert(!rollsForThisCSC.empty());
143 
144  if (debug)
145  LogDebug("CSCSegtoRPC") << "CSC \t \t Loop over all the rolls asociated to this CSC" << std::endl;
146  for (std::set<RPCDetId>::iterator iteraRoll = rollsForThisCSC.begin(); iteraRoll != rollsForThisCSC.end();
147  iteraRoll++) {
148  const RPCRoll* rollasociated = rpcGeo->roll(*iteraRoll);
149  RPCDetId rpcId = rollasociated->id();
150 
151  if (debug)
152  LogDebug("CSCSegtoRPC") << "CSC \t \t \t We are in the roll getting the surface" << rpcId << std::endl;
153  const BoundPlane& RPCSurface = rollasociated->surface();
154 
155  if (debug)
156  LogDebug("CSCSegtoRPC") << "CSC \t \t \t RollID: " << rpcId << std::endl;
157 
158  if (debug)
159  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Doing the extrapolation to this roll" << std::endl;
160  if (debug)
161  LogDebug("CSCSegtoRPC") << "CSC \t \t \t CSC Segment Direction in CSCLocal " << segmentDirection
162  << std::endl;
163  if (debug)
164  LogDebug("CSCSegtoRPC") << "CSC \t \t \t CSC Segment Point in CSCLocal " << segmentPosition
165  << std::endl;
166 
167  GlobalPoint CenterPointRollGlobal = RPCSurface.toGlobal(LocalPoint(0, 0, 0));
168  if (debug)
169  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Center (0,0,0) of the Roll in Global" << CenterPointRollGlobal
170  << std::endl;
171  GlobalPoint CenterPointCSCGlobal = TheChamber->toGlobal(LocalPoint(0, 0, 0));
172  if (debug)
173  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Center (0,0,0) of the CSC in Global" << CenterPointCSCGlobal
174  << std::endl;
175  GlobalPoint segmentPositionInGlobal =
176  TheChamber->toGlobal(segmentPosition); //new way to convert to global
177  if (debug)
178  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Segment Position in Global" << segmentPositionInGlobal
179  << std::endl;
180  LocalPoint CenterRollinCSCFrame = TheChamber->toLocal(CenterPointRollGlobal);
181 
182  if (debug) { //to check CSC RPC phi relation!
183  float rpcphi = 0;
184  float cscphi = 0;
185 
186  (CenterPointRollGlobal.barePhi() < 0) ? rpcphi = 2 * 3.141592 + CenterPointRollGlobal.barePhi()
187  : rpcphi = CenterPointRollGlobal.barePhi();
188 
189  (CenterPointCSCGlobal.barePhi() < 0) ? cscphi = 2 * 3.1415926536 + CenterPointCSCGlobal.barePhi()
190  : cscphi = CenterPointCSCGlobal.barePhi();
191 
192  float df = fabs(cscphi - rpcphi);
193  float dr = fabs(CenterPointRollGlobal.perp() - CenterPointCSCGlobal.perp());
194  float diffz = CenterPointRollGlobal.z() - CenterPointCSCGlobal.z();
195  float dfg = df * 180. / 3.14159265;
196 
197  if (debug)
198  LogDebug("CSCSegtoRPC") << "CSC \t \t \t z of RPC=" << CenterPointRollGlobal.z() << "z of CSC"
199  << CenterPointCSCGlobal.z() << " dfg=" << dfg << std::endl;
200 
201  RPCGeomServ rpcsrv(rpcId);
202 
203  if (dr > 200. || fabs(dz) > 55. || dfg > 1.) {
204  if (debug)
205  LogDebug("CSCSegtoRPC")
206  << "\t \t \t CSC Station= " << CSCId.station() << " Ring= " << CSCId.ring()
207  << " Chamber= " << CSCId.chamber() << " cscphi=" << cscphi * 180 / 3.14159265
208  << "\t RPC Station= " << rpcId.station() << " ring= " << rpcId.ring() << " segment =-> "
209  << rpcsrv.segment() << " rollphi=" << rpcphi * 180 / 3.14159265 << "\t dfg=" << dfg
210  << " dz=" << diffz << " dr=" << dr << std::endl;
211  }
212  }
213 
214  float D = CenterRollinCSCFrame.z();
215 
216  float X = Xo + dx * D / dz;
217  float Y = Yo + dy * D / dz;
218  float Z = D;
219 
220  const TrapezoidalStripTopology* top_ =
221  dynamic_cast<const TrapezoidalStripTopology*>(&(rollasociated->topology()));
222  LocalPoint xmin = top_->localPosition(0.);
223  if (debug)
224  LogDebug("CSCSegtoRPC") << "CSC \t \t \t xmin of this Roll " << xmin << "cm" << std::endl;
225  LocalPoint xmax = top_->localPosition((float)rollasociated->nstrips());
226  if (debug)
227  LogDebug("CSCSegtoRPC") << "CSC \t \t \t xmax of this Roll " << xmax << "cm" << std::endl;
228  float rsize = fabs(xmax.x() - xmin.x());
229  if (debug)
230  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Roll Size " << rsize << "cm" << std::endl;
231  float stripl = top_->stripLength();
232  float stripw = top_->pitch();
233 
234  if (debug)
235  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Strip Lenght " << stripl << "cm" << std::endl;
236  if (debug)
237  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Strip Width " << stripw << "cm" << std::endl;
238 
239  if (debug)
240  LogDebug("CSCSegtoRPC") << "CSC \t \t \t X Predicted in CSCLocal= " << X << "cm" << std::endl;
241  if (debug)
242  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Y Predicted in CSCLocal= " << Y << "cm" << std::endl;
243  if (debug)
244  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Z Predicted in CSCLocal= " << Z << "cm" << std::endl;
245 
246  float extrapolatedDistance = sqrt((X - Xo) * (X - Xo) + (Y - Yo) * (Y - Yo) + (Z - Zo) * (Z - Zo));
247 
248  if (debug)
249  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Is the distance of extrapolation less than MaxD? ="
250  << extrapolatedDistance << "cm"
251  << " MaxD=" << MaxD << "cm" << std::endl;
252 
253  if (extrapolatedDistance <= MaxD) {
254  if (debug)
255  LogDebug("CSCSegtoRPC") << "CSC \t \t \t yes" << std::endl;
256 
257  GlobalPoint GlobalPointExtrapolated = TheChamber->toGlobal(LocalPoint(X, Y, Z));
258  if (debug)
259  LogDebug("CSCSegtoRPC")
260  << "CSC \t \t \t Point ExtraPolated in Global" << GlobalPointExtrapolated << std::endl;
261 
262  LocalPoint PointExtrapolatedRPCFrame = RPCSurface.toLocal(GlobalPointExtrapolated);
263 
264  if (debug)
265  LogDebug("CSCSegtoRPC")
266  << "CSC \t \t \t Point Extrapolated in RPCLocal" << PointExtrapolatedRPCFrame << std::endl;
267  if (debug)
268  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Corner of the Roll = (" << rsize * eyr << "," << stripl * eyr
269  << ")" << std::endl;
270  if (debug)
271  LogDebug("CSCSegtoRPC")
272  << "CSC \t \t \t Info About the Point Extrapolated in X Abs ("
273  << fabs(PointExtrapolatedRPCFrame.x()) << "," << fabs(PointExtrapolatedRPCFrame.y()) << ","
274  << fabs(PointExtrapolatedRPCFrame.z()) << ")" << std::endl;
275  if (debug)
276  LogDebug("CSCSegtoRPC") << "CSC \t \t \t dz=" << fabs(PointExtrapolatedRPCFrame.z())
277  << " dx=" << fabs(PointExtrapolatedRPCFrame.x())
278  << " dy=" << fabs(PointExtrapolatedRPCFrame.y()) << std::endl;
279 
280  if (debug)
281  LogDebug("CSCSegtoRPC") << "CSC \t \t \t Does the extrapolation go inside this roll????" << std::endl;
282 
283  if (fabs(PointExtrapolatedRPCFrame.z()) < 1. && fabs(PointExtrapolatedRPCFrame.x()) < rsize * eyr &&
284  fabs(PointExtrapolatedRPCFrame.y()) < stripl * eyr) {
285  if (debug)
286  LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t yes" << std::endl;
287  if (debug)
288  LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t Creating the RecHit" << std::endl;
289 
290  RPCRecHit RPCPoint(rpcId, myBx, PointExtrapolatedRPCFrame);
291  RPCPoint.setTimeAndError(myTime, myTimeErr);
292 
293  if (debug)
294  LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t Clearing the vector" << std::endl;
295  RPCPointVector.clear();
296  if (debug)
297  LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t Pushing back" << std::endl;
298  RPCPointVector.push_back(RPCPoint);
299  if (debug)
300  LogDebug("CSCSegtoRPC") << "CSC \t \t \t \t Putting the vector" << std::endl;
301  _ThePoints->put(rpcId, RPCPointVector.begin(), RPCPointVector.end());
302  }
303  }
304  }
305  }
306  }
307  }
308  }
309  }
310  return _ThePoints;
311 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
T perp() const
Definition: PV3DBase.h:69
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
const CSCChamber * chamber(CSCDetId id) const
Return the chamber corresponding to given DetId.
Definition: CSCGeometry.cc:100
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:58
T z() const
Definition: PV3DBase.h:61
#define X(str)
Definition: MuonsGrabber.cc:38
int ring() const
Definition: RPCDetId.h:59
edm::ESGetToken< CSCGeometry, MuonGeometryRecord > cscGeoToken_
Definition: CSCSegtoRPC.h:28
const RPCRoll * roll(RPCDetId id) const
Return a roll given its id.
Definition: RPCGeometry.cc:50
int nstrips() const
Definition: RPCRoll.cc:24
assert(be >=bs)
T barePhi() const
Definition: PV3DBase.h:65
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
T sqrt(T t)
Definition: SSEVec.h:19
int chamber() const
Definition: CSCDetId.h:62
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:130
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
void setTimeAndError(float time, float err)
Set the time and its error.
Definition: RPCRecHit.h:65
virtual int segment()
Definition: RPCGeomServ.cc:361
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:37
#define debug
Definition: HDRShower.cc:19
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:141
int station() const
Definition: CSCDetId.h:79
int endcap() const
Definition: CSCDetId.h:85
int station() const
Definition: RPCDetId.h:78
const Topology & topology() const override
Definition: RPCRoll.cc:18
RPCDetId id() const
Definition: RPCRoll.cc:16
std::unique_ptr< RPCRecHitCollection > thePoints(CSCSegmentCollection const *allCSCSegments, edm::EventSetup const &iSetup, bool debug, double eyr)
Definition: CSCSegtoRPC.cc:22
int ring() const
Definition: CSCDetId.h:68
float stripLength() const override
det heigth (strip length in the middle)
edm::ESGetToken< CSCObjectMap, MuonGeometryRecord > cscMapToken_
Definition: CSCSegtoRPC.h:29
CSCSegtoRPC(edm::ConsumesCollector iC, const edm::ParameterSet &)
Definition: CSCSegtoRPC.cc:16
edm::ESGetToken< RPCGeometry, MuonGeometryRecord > rpcGeoToken_
Definition: CSCSegtoRPC.h:27
std::set< RPCDetId > const & getRolls(CSCStationIndex index) const
Definition: CSCObjectMap.cc:42
LocalPoint localPosition(float strip) const override
#define LogDebug(id)
CSCDetId id() const
Get the (concrete) DetId.
Definition: CSCChamber.h:34