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