CMS 3D CMS Logo

CSCSegtoRPC.cc
Go to the documentation of this file.
14 
15 CSCSegtoRPC::CSCSegtoRPC(const CSCSegmentCollection * allCSCSegments, const edm::EventSetup& iSetup, bool debug, double eyr){
16 
20 
21  iSetup.get<MuonGeometryRecord>().get(rpcGeo);
22  iSetup.get<MuonGeometryRecord>().get(cscGeo);
23  iSetup.get<MuonGeometryRecord>().get(cscMap);
24 
25  MaxD=80.;
26 
27  if(debug) std::cout<<"CSC \t Number of CSC Segments in this event = "<<allCSCSegments->size()<<std::endl;
28 
29  _ThePoints = std::make_unique<RPCRecHitCollection>();
30 
31  if(allCSCSegments->size()==0){
32  if(debug) std::cout<<"CSC 0 segments skiping event"<<std::endl;
33  }else {
34  std::map<CSCDetId,int> CSCSegmentsCounter;
36 
37  int segmentsInThisEventInTheEndcap=0;
38 
39  for (segment = allCSCSegments->begin();segment!=allCSCSegments->end(); ++segment){
40  CSCSegmentsCounter[segment->cscDetId()]++;
41  segmentsInThisEventInTheEndcap++;
42  }
43 
44  if(debug) std::cout<<"CSC \t loop over all the CSCSegments "<<std::endl;
45  for (segment = allCSCSegments->begin();segment!=allCSCSegments->end(); ++segment){
46  CSCDetId CSCId = segment->cscDetId();
47 
48  if(debug) std::cout<<"CSC \t \t This Segment is in Chamber id: "<<CSCId<<std::endl;
49  if(debug) std::cout<<"CSC \t \t Number of segments in this CSC = "<<CSCSegmentsCounter[CSCId]<<std::endl;
50  if(debug) std::cout<<"CSC \t \t Is the only one in this CSC? is not ind the ring 1 or station 4? Are there more than 2 segments in the event?"<<std::endl;
51 
52  if(CSCSegmentsCounter[CSCId]==1 && CSCId.station()!=4 && CSCId.ring()!=1 && allCSCSegments->size()>=2){
53  if(debug) std::cout<<"CSC \t \t yes"<<std::endl;
54  int cscEndCap = CSCId.endcap();
55  int cscStation = CSCId.station();
56  int cscRing = CSCId.ring();
57  int rpcRegion = 1; if(cscEndCap==2) rpcRegion= -1;//Relacion entre las endcaps
58  int rpcRing = cscRing;
59  if(cscRing==4)rpcRing =1;
60  int rpcStation = cscStation;
61  int rpcSegment = CSCId.chamber();
62 
63  LocalPoint segmentPosition= segment->localPosition();
64  LocalVector segmentDirection=segment->localDirection();
65  float dz=segmentDirection.z();
66 
67  if(debug) std::cout<<"CSC \t \t \t Information about the segment"
68  <<"RecHits ="<<segment->nRecHits()
69  <<"Angle ="<<acos(dz)*180/3.1415926<<std::endl;
70 
71  if(debug) std::cout<<"CSC \t \t Is a good Segment? dim = 4, 4 <= nRecHits <= 10 Incident angle int range 45 < "<<acos(dz)*180/3.1415926<<" < 135? "<<std::endl;
72 
73  if((segment->dimension()==4) && (segment->nRecHits()<=10 && segment->nRecHits()>=4)){
74  //&& acos(dz)*180/3.1415926 > 45. && acos(dz)*180/3.1415926 < 135.){
75  //&& segment->chi2()< ??)Add 3 segmentes in the endcaps???
76 
77 
78  if(debug) std::cout<<"CSC \t \t yes"<<std::endl;
79  if(debug) std::cout<<"CSC \t \t CSC Segment Dimension "<<segment->dimension()<<std::endl;
80 
81  float Xo=segmentPosition.x();
82  float Yo=segmentPosition.y();
83  float Zo=segmentPosition.z();
84  float dx=segmentDirection.x();
85  float dy=segmentDirection.y();
86  float dz=segmentDirection.z();
87 
88  if(debug) std::cout<<"Creating the CSCIndex"<<std::endl;
89  CSCStationIndex theindex(rpcRegion,rpcStation,rpcRing,rpcSegment);
90  if(debug) std::cout<<"Getting the Rolls for the given index"<<std::endl;
91  std::set<RPCDetId> rollsForThisCSC = cscMap->getRolls(theindex);
92  if(debug) std::cout<<"CSC \t \t Getting chamber from Geometry"<<std::endl;
93  const CSCChamber* TheChamber=cscGeo->chamber(CSCId);
94  if(debug) std::cout<<"CSC \t \t Getting ID from Chamber"<<std::endl;
95  const CSCDetId TheId=TheChamber->id();
96 
97  if(debug) std::cout<<"CSC \t \t Number of rolls for this CSC = "<<rollsForThisCSC.size()<<std::endl;
98 
99  if(debug) std::cout<<"CSC \t \t Printing The Id"<<TheId<<std::endl;
100 
101  if(rpcRing!=1&&rpcStation!=4){//They don't exist!
102 
103  assert(!rollsForThisCSC.empty());
104 
105  if(debug) std::cout<<"CSC \t \t Loop over all the rolls asociated to this CSC"<<std::endl;
106  for (std::set<RPCDetId>::iterator iteraRoll = rollsForThisCSC.begin();iteraRoll != rollsForThisCSC.end(); iteraRoll++){
107  const RPCRoll* rollasociated = rpcGeo->roll(*iteraRoll);
108  RPCDetId rpcId = rollasociated->id();
109 
110  if(debug) std::cout<<"CSC \t \t \t We are in the roll getting the surface"<<rpcId<<std::endl;
111  const BoundPlane & RPCSurface = rollasociated->surface();
112 
113  if(debug) std::cout<<"CSC \t \t \t RollID: "<<rpcId<<std::endl;
114 
115  if(debug) std::cout<<"CSC \t \t \t Doing the extrapolation to this roll"<<std::endl;
116  if(debug) std::cout<<"CSC \t \t \t CSC Segment Direction in CSCLocal "<<segmentDirection<<std::endl;
117  if(debug) std::cout<<"CSC \t \t \t CSC Segment Point in CSCLocal "<<segmentPosition<<std::endl;
118 
119  GlobalPoint CenterPointRollGlobal = RPCSurface.toGlobal(LocalPoint(0,0,0));
120  if(debug) std::cout<<"CSC \t \t \t Center (0,0,0) of the Roll in Global"<<CenterPointRollGlobal<<std::endl;
121  GlobalPoint CenterPointCSCGlobal = TheChamber->toGlobal(LocalPoint(0,0,0));
122  if(debug) std::cout<<"CSC \t \t \t Center (0,0,0) of the CSC in Global"<<CenterPointCSCGlobal<<std::endl;
123  GlobalPoint segmentPositionInGlobal=TheChamber->toGlobal(segmentPosition); //new way to convert to global
124  if(debug) std::cout<<"CSC \t \t \t Segment Position in Global"<<segmentPositionInGlobal<<std::endl;
125  LocalPoint CenterRollinCSCFrame = TheChamber->toLocal(CenterPointRollGlobal);
126 
127  if(debug){//to check CSC RPC phi relation!
128  float rpcphi=0;
129  float cscphi=0;
130 
131  (CenterPointRollGlobal.barePhi()<0)?
132  rpcphi = 2*3.141592+CenterPointRollGlobal.barePhi():rpcphi=CenterPointRollGlobal.barePhi();
133 
134  (CenterPointCSCGlobal.barePhi()<0)?
135  cscphi = 2*3.1415926536+CenterPointCSCGlobal.barePhi():cscphi=CenterPointCSCGlobal.barePhi();
136 
137  float df=fabs(cscphi-rpcphi);
138  float dr=fabs(CenterPointRollGlobal.perp()-CenterPointCSCGlobal.perp());
139  float diffz=CenterPointRollGlobal.z()-CenterPointCSCGlobal.z();
140  float dfg=df*180./3.14159265;
141 
142  if(debug) std::cout<<"CSC \t \t \t z of RPC="<<CenterPointRollGlobal.z()<<"z of CSC"<<CenterPointCSCGlobal.z()<<" dfg="<<dfg<<std::endl;
143 
144  RPCGeomServ rpcsrv(rpcId);
145 
146  if(dr>200.||fabs(dz)>55.||dfg>1.){
147  //if(rpcRegion==1&&dfg>1.&&dr>100.){
148  if (debug) std::cout
149  <<"\t \t \t CSC Station= "<<CSCId.station()
150  <<" Ring= "<<CSCId.ring()
151  <<" Chamber= "<<CSCId.chamber()
152  <<" cscphi="<<cscphi*180/3.14159265
153  <<"\t RPC Station= "<<rpcId.station()
154  <<" ring= "<<rpcId.ring()
155  <<" segment =-> "<<rpcsrv.segment()
156  <<" rollphi="<<rpcphi*180/3.14159265
157  <<"\t dfg="<<dfg
158  <<" dz="<<diffz
159  <<" dr="<<dr
160  <<std::endl;
161 
162  }
163  }
164 
165  float D=CenterRollinCSCFrame.z();
166 
167  float X=Xo+dx*D/dz;
168  float Y=Yo+dy*D/dz;
169  float Z=D;
170 
171  const TrapezoidalStripTopology* top_=dynamic_cast<const TrapezoidalStripTopology*>(&(rollasociated->topology()));
172  LocalPoint xmin = top_->localPosition(0.);
173  if(debug) std::cout<<"CSC \t \t \t xmin of this Roll "<<xmin<<"cm"<<std::endl;
174  LocalPoint xmax = top_->localPosition((float)rollasociated->nstrips());
175  if(debug) std::cout<<"CSC \t \t \t xmax of this Roll "<<xmax<<"cm"<<std::endl;
176  float rsize = fabs( xmax.x()-xmin.x() );
177  if(debug) std::cout<<"CSC \t \t \t Roll Size "<<rsize<<"cm"<<std::endl;
178  float stripl = top_->stripLength();
179  float stripw = top_->pitch();
180 
181  if(debug) std::cout<<"CSC \t \t \t Strip Lenght "<<stripl<<"cm"<<std::endl;
182  if(debug) std::cout<<"CSC \t \t \t Strip Width "<<stripw<<"cm"<<std::endl;
183 
184  if(debug) std::cout<<"CSC \t \t \t X Predicted in CSCLocal= "<<X<<"cm"<<std::endl;
185  if(debug) std::cout<<"CSC \t \t \t Y Predicted in CSCLocal= "<<Y<<"cm"<<std::endl;
186  if(debug) std::cout<<"CSC \t \t \t Z Predicted in CSCLocal= "<<Z<<"cm"<<std::endl;
187 
188  float extrapolatedDistance = sqrt((X-Xo)*(X-Xo)+(Y-Yo)*(Y-Yo)+(Z-Zo)*(Z-Zo));
189 
190  if(debug) std::cout<<"CSC \t \t \t Is the distance of extrapolation less than MaxD? ="<<extrapolatedDistance<<"cm"<<" MaxD="<<MaxD<<"cm"<<std::endl;
191 
192  if(extrapolatedDistance<=MaxD){
193 
194  if(debug) std::cout<<"CSC \t \t \t yes"<<std::endl;
195 
196  GlobalPoint GlobalPointExtrapolated=TheChamber->toGlobal(LocalPoint(X,Y,Z));
197  if(debug) std::cout<<"CSC \t \t \t Point ExtraPolated in Global"<<GlobalPointExtrapolated<< std::endl;
198 
199 
200  LocalPoint PointExtrapolatedRPCFrame = RPCSurface.toLocal(GlobalPointExtrapolated);
201 
202  if(debug) std::cout<<"CSC \t \t \t Point Extrapolated in RPCLocal"<<PointExtrapolatedRPCFrame<< std::endl;
203  if(debug) std::cout<<"CSC \t \t \t Corner of the Roll = ("<<rsize*eyr<<","<<stripl*eyr<<")"<<std::endl;
204  if(debug) std::cout<<"CSC \t \t \t Info About the Point Extrapolated in X Abs ("<<fabs(PointExtrapolatedRPCFrame.x())<<","
205  <<fabs(PointExtrapolatedRPCFrame.y())<<","<<fabs(PointExtrapolatedRPCFrame.z())<<")"<<std::endl;
206  if(debug) std::cout<<"CSC \t \t \t dz="
207  <<fabs(PointExtrapolatedRPCFrame.z())<<" dx="
208  <<fabs(PointExtrapolatedRPCFrame.x())<<" dy="
209  <<fabs(PointExtrapolatedRPCFrame.y())<<std::endl;
210 
211  if(debug) std::cout<<"CSC \t \t \t Does the extrapolation go inside this roll????"<<std::endl;
212 
213  if(fabs(PointExtrapolatedRPCFrame.z()) < 1. &&
214  fabs(PointExtrapolatedRPCFrame.x()) < rsize*eyr &&
215  fabs(PointExtrapolatedRPCFrame.y()) < stripl*eyr){
216  if(debug) std::cout<<"CSC \t \t \t \t yes"<<std::endl;
217  if(debug) std::cout<<"CSC \t \t \t \t Creating the RecHit"<<std::endl;
218  RPCRecHit RPCPoint(rpcId,0,PointExtrapolatedRPCFrame);
219  if(debug) std::cout<<"CSC \t \t \t \t Clearing the vector"<<std::endl;
221  if(debug) std::cout<<"CSC \t \t \t \t Pushing back"<<std::endl;
222  RPCPointVector.push_back(RPCPoint);
223  if(debug) std::cout<<"CSC \t \t \t \t Putting the vector"<<std::endl;
225  }
226  }
227  }
228  }
229  }
230  }
231  }
232  }
233 }
234 
236 
237 }
int chamber() const
Definition: CSCDetId.h:68
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
LocalPoint localPosition(float strip) const override
T perp() const
Definition: PV3DBase.h:72
CSCDetId id() const
Get the (concrete) DetId.
Definition: CSCChamber.h:37
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
int nstrips() const
Definition: RPCRoll.cc:46
T y() const
Definition: PV3DBase.h:63
#define X(str)
Definition: MuonsGrabber.cc:48
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:69
std::set< RPCDetId > const & getRolls(CSCStationIndex index) const
Definition: CSCObjectMap.cc:48
const Topology & topology() const override
Definition: RPCRoll.cc:30
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
double MaxD
Definition: CSCSegtoRPC.h:23
iterator begin()
Definition: OwnVector.h:244
int endcap() const
Definition: CSCDetId.h:93
std::unique_ptr< RPCRecHitCollection > _ThePoints
Definition: CSCSegtoRPC.h:20
T barePhi() const
Definition: PV3DBase.h:68
void push_back(D *&d)
Definition: OwnVector.h:290
RPCDetId id() const
Definition: RPCRoll.cc:24
int ring() const
Definition: RPCDetId.h:72
T sqrt(T t)
Definition: SSEVec.h:18
T z() const
Definition: PV3DBase.h:64
void clear()
Definition: OwnVector.h:445
float stripLength() const override
det heigth (strip length in the middle)
CSCSegtoRPC(CSCSegmentCollection const *allCSCSegments, edm::EventSetup const &iSetup, bool debug, double eyr)
Definition: CSCSegtoRPC.cc:15
int ring() const
Definition: CSCDetId.h:75
iterator end()
Definition: OwnVector.h:249
virtual int segment()
Definition: RPCGeomServ.cc:469
#define debug
Definition: HDRShower.cc:19
DecomposeProduct< arg, typename Div::arg > D
Definition: Factorize.h:151
const T & get() const
Definition: EventSetup.h:55
const CSCChamber * chamber(CSCDetId id) const
Return the chamber corresponding to given DetId.
Definition: CSCGeometry.cc:118
edm::OwnVector< RPCRecHit > RPCPointVector
Definition: CSCSegtoRPC.h:21
int station() const
Definition: CSCDetId.h:86
T x() const
Definition: PV3DBase.h:62
const RPCRoll * roll(RPCDetId id) const
Return a roll given its id.
Definition: RPCGeometry.cc:75
int station() const
Definition: RPCDetId.h:96