CMS 3D CMS Logo

FastTrackerRecHitMatcher.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 
4 // framework stuff
11 
12 // fast tracker rechits
17 
18 // geometry stuff
25 
26 
27 
28 // sim stuff
31 
33 
34  public:
35 
38  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
39 
40  private:
41 
42  void produce(edm::Event&, const edm::EventSetup&) override;
43 
44  // ---------- typedefs -----------------------------
45  typedef std::pair<LocalPoint,LocalPoint> StripPosition;
46 
47  // ----------internal functions ---------------------------
48 
49  // create projected hit
50  std::unique_ptr<FastTrackerRecHit> projectOnly( const FastSingleTrackerRecHit *originalRH,
51  const GeomDet * monoDet,
52  const GluedGeomDet* gluedDet,
53  LocalVector& ldir) const;
54 
55  // create matched hit
56  std::unique_ptr<FastTrackerRecHit> match( const FastSingleTrackerRecHit *monoRH,
57  const FastSingleTrackerRecHit *stereoRH,
58  const GluedGeomDet* gluedDet,
59  LocalVector& trackdirection,
60  bool stereLayerFirst) const;
61 
62  StripPosition project(const GeomDetUnit *det,
63  const GluedGeomDet* glueddet,
64  const StripPosition& strip,
65  const LocalVector& trackdirection)const;
66 
68  if(!recHit->isSingle()){
69  throw cms::Exception("FastTrackerRecHitMatcher") << "all rechits in simHit2RecHitMap must be instances of FastSingleTrackerRecHit. recHit's rtti: " << recHit->rtti() << std::endl;
70  }
71  return dynamic_cast<const FastSingleTrackerRecHit *>(recHit);
72  }
73 
74  // ----------member data ---------------------------
77 
78 };
79 
81 
82 {
83  simHitsToken = consumes<edm::PSimHitContainer>(iConfig.getParameter<edm::InputTag>("simHits"));
84  simHit2RecHitMapToken = consumes<FastTrackerRecHitRefCollection>(iConfig.getParameter<edm::InputTag>("simHit2RecHitMap"));
85 
86  produces<FastTrackerRecHitCollection>();
87  produces<FastTrackerRecHitRefCollection>("simHit2RecHitMap");
88 }
89 
90 
92 {
93  // services
95  iSetup.get<TrackerDigiGeometryRecord> ().get (geometry);
96 
97  // input
99  iEvent.getByToken(simHitsToken,simHits);
100 
102  iEvent.getByToken(simHit2RecHitMapToken,simHit2RecHitMap);
103 
104  // output
105  auto output_recHits = std::make_unique<FastTrackerRecHitCollection>();
106  auto output_simHit2RecHitMap = std::make_unique<FastTrackerRecHitRefCollection>(simHit2RecHitMap->size(),FastTrackerRecHitRef());
108 
109  bool skipNext = false;
110  for(unsigned simHitCounter = 0;simHitCounter < simHits->size();++simHitCounter){
111 
112  // skip hit in case it was matched to previous one
113  if(skipNext){
114  skipNext = false;
115  continue;
116  }
117  skipNext = false;
118 
119  // get simHit and associated recHit
120  const PSimHit & simHit = (*simHits)[simHitCounter];
121  const FastTrackerRecHitRef & recHitRef = (*simHit2RecHitMap)[simHitCounter];
122 
123  // skip simHits w/o associated recHit
124  if(recHitRef.isNull())
125  continue;
126 
127  // cast
128  const FastSingleTrackerRecHit * recHit = _cast2Single(recHitRef.get());
129 
130  // get subdetector id
131  DetId detid = recHit->geographicalId();
132  unsigned int subdet = detid.subdetId();
133 
134  // treat pixel hits
135  if(subdet <= 2){
136  (*output_simHit2RecHitMap)[simHitCounter] = recHitRef;
137  }
138 
139  // treat strip hits
140  else{
141 
142  StripSubdetector stripSubDetId(detid);
143 
144  // treat regular regular strip hits
145  if(!stripSubDetId.glued()){
146  (*output_simHit2RecHitMap)[simHitCounter] = recHitRef;
147  }
148 
149  // treat strip hits on glued layers
150  else{
151 
152  // Obtain direction of simtrack at simhit in local coordinates of glued module
153  // - direction of simtrack at simhit, in coordindates of the single module
154  LocalVector localSimTrackDir = simHit.localDirection();
155  // - transform to global coordinates
156  GlobalVector globalSimTrackDir= recHit->det()->surface().toGlobal(localSimTrackDir);
157  // - transform to local coordinates of glued module
158  const GluedGeomDet* gluedDet = (const GluedGeomDet*)geometry->idToDet(DetId(stripSubDetId.glued()));
159  LocalVector gluedLocalSimTrackDir = gluedDet->surface().toLocal(globalSimTrackDir);
160 
161  // check whether next hit is partner
162  const FastSingleTrackerRecHit * partnerRecHit = nullptr;
163  // - there must be a next hit
164  if(simHitCounter + 1 < simHits->size()){
165  const FastTrackerRecHitRef & nextRecHitRef = (*simHit2RecHitMap)[simHitCounter + 1];
166  const PSimHit & nextSimHit = (*simHits)[simHitCounter + 1];
167  // - partner hit must not be null
168  // - simHit and partner simHit must belong to same simTrack
169  // - partner hit must be on the module glued to the module of the hit
170  if( (!nextRecHitRef.isNull())
171  && simHit.trackId() == nextSimHit.trackId()
172  && StripSubdetector( nextRecHitRef->geographicalId() ).partnerDetId() == detid.rawId() ) {
173  partnerRecHit = _cast2Single(nextRecHitRef.get());
174  skipNext = true;
175  }
176  }
177 
178  std::unique_ptr<FastTrackerRecHit> newRecHit(nullptr);
179 
180  // if partner found: create a matched hit
181  if( partnerRecHit ){
182  newRecHit = match( stripSubDetId.stereo() ? partnerRecHit : recHit,
183  stripSubDetId.stereo() ? recHit : partnerRecHit,
184  gluedDet , gluedLocalSimTrackDir,
185  stripSubDetId.stereo());
186  }
187  // else: create projected hit
188  else{
189  newRecHit = projectOnly( recHit , geometry->idToDet(detid),gluedDet, gluedLocalSimTrackDir );
190  }
191  output_recHits->push_back(std::move(newRecHit));
192  (*output_simHit2RecHitMap)[simHitCounter] = FastTrackerRecHitRef(output_recHits_refProd,output_recHits->size()-1);
193  }
194  }
195  }
196 
197  iEvent.put(std::move(output_recHits));
198  iEvent.put(std::move(output_simHit2RecHitMap),"simHit2RecHitMap");
199 
200 }
201 
202 std::unique_ptr<FastTrackerRecHit> FastTrackerRecHitMatcher::match(const FastSingleTrackerRecHit *monoRH,
203  const FastSingleTrackerRecHit *stereoRH,
204  const GluedGeomDet* gluedDet,
205  LocalVector& trackdirection,
206  bool stereoHitFirst) const
207 {
208 
209  // stripdet = mono
210  // partnerstripdet = stereo
211  const GeomDetUnit* stripdet = gluedDet->monoDet();
212  const GeomDetUnit* partnerstripdet = gluedDet->stereoDet();
213  const StripTopology& topol=(const StripTopology&)stripdet->topology();
214 
216 
217  // position of the initial and final point of the strip (RPHI cluster) in local strip coordinates
218  MeasurementPoint RPHIpoint=topol.measurementPosition(monoRH->localPosition());
219  MeasurementPoint RPHIpointini=MeasurementPoint(RPHIpoint.x(),-0.5);
220  MeasurementPoint RPHIpointend=MeasurementPoint(RPHIpoint.x(),0.5);
221 
222  // position of the initial and final point of the strip in local coordinates (mono det)
223  StripPosition stripmono=StripPosition(topol.localPosition(RPHIpointini),topol.localPosition(RPHIpointend));
224 
225  if(trackdirection.mag2()<FLT_MIN){// in case of no track hypothesis assume a track from the origin through the center of the strip
226  LocalPoint lcenterofstrip=monoRH->localPosition();
227  GlobalPoint gcenterofstrip=(stripdet->surface()).toGlobal(lcenterofstrip);
228  GlobalVector gtrackdirection=gcenterofstrip-GlobalPoint(0,0,0);
229  trackdirection=(gluedDet->surface()).toLocal(gtrackdirection);
230  }
231 
232  //project mono hit on glued det
233  StripPosition projectedstripmono=project(stripdet,gluedDet,stripmono,trackdirection);
234  const StripTopology& partnertopol=(const StripTopology&)partnerstripdet->topology();
235 
236  //error calculation (the part that depends on mono RH only)
237  LocalVector RPHIpositiononGluedendvector=projectedstripmono.second-projectedstripmono.first;
238  double c1=sin(RPHIpositiononGluedendvector.phi());
239  double s1=-cos(RPHIpositiononGluedendvector.phi());
240  MeasurementError errormonoRH=topol.measurementError(monoRH->localPosition(),monoRH->localPositionError());
241  double pitch=topol.localPitch(monoRH->localPosition());
242  double sigmap12=errormonoRH.uu()*pitch*pitch;
243 
244  // position of the initial and final point of the strip (STEREO cluster)
245  MeasurementPoint STEREOpoint=partnertopol.measurementPosition(stereoRH->localPosition());
246 
247  MeasurementPoint STEREOpointini=MeasurementPoint(STEREOpoint.x(),-0.5);
248  MeasurementPoint STEREOpointend=MeasurementPoint(STEREOpoint.x(),0.5);
249 
250  // position of the initial and final point of the strip in local coordinates (stereo det)
251  StripPosition stripstereo(partnertopol.localPosition(STEREOpointini),partnertopol.localPosition(STEREOpointend));
252 
253  //project stereo hit on glued det
254  StripPosition projectedstripstereo=project(partnerstripdet,gluedDet,stripstereo,trackdirection);
255 
256  //perform the matching
257  //(x2-x1)(y-y1)=(y2-y1)(x-x1)
259  m(0,0)=-(projectedstripmono.second.y()-projectedstripmono.first.y()); m(0,1)=(projectedstripmono.second.x()-projectedstripmono.first.x());
260  m(1,0)=-(projectedstripstereo.second.y()-projectedstripstereo.first.y()); m(1,1)=(projectedstripstereo.second.x()-projectedstripstereo.first.x());
261  c(0)=m(0,1)*projectedstripmono.first.y()+m(0,0)*projectedstripmono.first.x();
262  c(1)=m(1,1)*projectedstripstereo.first.y()+m(1,0)*projectedstripstereo.first.x();
263  m.Invert(); solution = m * c;
264  position=LocalPoint(solution(0),solution(1));
265 
266 
267  //
268  // temporary fix by tommaso
269  //
270 
271 
272  LocalError tempError (100,0,100);
273 
274  // calculate the error
275  LocalVector stereopositiononGluedendvector=projectedstripstereo.second-projectedstripstereo.first;
276  double c2=sin(stereopositiononGluedendvector.phi()); double s2=-cos(stereopositiononGluedendvector.phi());
277  MeasurementError errorstereoRH=partnertopol.measurementError(stereoRH->localPosition(), stereoRH->localPositionError());
278  pitch=partnertopol.localPitch(stereoRH->localPosition());
279  double sigmap22=errorstereoRH.uu()*pitch*pitch;
280  double diff=(c1*s2-c2*s1);
281  double invdet2=1/(diff*diff);
282  float xx=invdet2*(sigmap12*s2*s2+sigmap22*s1*s1);
283  float xy=-invdet2*(sigmap12*c2*s2+sigmap22*c1*s1);
284  float yy=invdet2*(sigmap12*c2*c2+sigmap22*c1*c1);
285  LocalError error=LocalError(xx,xy,yy);
286 
287  //Added by DAO to make sure y positions are zero.
288  DetId det(monoRH->geographicalId());
289  if(det.subdetId() > 2) {
290  return std::make_unique<FastMatchedTrackerRecHit>(position, error, *gluedDet, *monoRH, *stereoRH,stereoHitFirst);
291 
292  }
293 
294  else {
295  throw cms::Exception("FastTrackerRecHitMatcher") << "Matched Pixel!?";
296  }
297 }
298 
299 
302  const GluedGeomDet* glueddet,
303  const StripPosition& strip,
304  const LocalVector& trackdirection)const
305 {
306 
307  GlobalPoint globalpointini=(det->surface()).toGlobal(strip.first);
308  GlobalPoint globalpointend=(det->surface()).toGlobal(strip.second);
309 
310  // position of the initial and final point of the strip in glued local coordinates
311  LocalPoint positiononGluedini=(glueddet->surface()).toLocal(globalpointini);
312  LocalPoint positiononGluedend=(glueddet->surface()).toLocal(globalpointend);
313 
314  //correct the position with the track direction
315 
316  float scale=-positiononGluedini.z()/trackdirection.z();
317 
318  LocalPoint projpositiononGluedini= positiononGluedini + scale*trackdirection;
319  LocalPoint projpositiononGluedend= positiononGluedend + scale*trackdirection;
320 
321  return StripPosition(projpositiononGluedini,projpositiononGluedend);
322 }
323 
324 
325 
326 std::unique_ptr<FastTrackerRecHit> FastTrackerRecHitMatcher::projectOnly( const FastSingleTrackerRecHit *originalRH,
327  const GeomDet * monoDet,
328  const GluedGeomDet* gluedDet,
329  LocalVector& ldir) const
330 {
331  LocalPoint position(originalRH->localPosition().x(), 0.,0.);
332  const BoundPlane& gluedPlane = gluedDet->surface();
333  const BoundPlane& hitPlane = monoDet->surface();
334 
335  double delta = gluedPlane.localZ( hitPlane.position());
336 
337  LocalPoint lhitPos = gluedPlane.toLocal( monoDet->surface().toGlobal(position ) );
338  LocalPoint projectedHitPos = lhitPos - ldir * delta/ldir.z();
339 
340  LocalVector hitXAxis = gluedPlane.toLocal( hitPlane.toGlobal( LocalVector(1,0,0)));
341  LocalError hitErr = originalRH->localPositionError();
342 
343  if (gluedPlane.normalVector().dot( hitPlane.normalVector()) < 0) {
344  // the two planes are inverted, and the correlation element must change sign
345  hitErr = LocalError( hitErr.xx(), -hitErr.xy(), hitErr.yy());
346  }
347  LocalError rotatedError = hitErr.rotate( hitXAxis.x(), hitXAxis.y());
348 
349 
350  const GeomDetUnit *gluedMonoDet = gluedDet->monoDet();
351  const GeomDetUnit *gluedStereoDet = gluedDet->stereoDet();
352  int isMono = 0;
353  int isStereo = 0;
354 
355  if(monoDet->geographicalId()==gluedMonoDet->geographicalId()) isMono = 1;
356  if(monoDet->geographicalId()==gluedStereoDet->geographicalId()) isStereo = 1;
357  //Added by DAO to make sure y positions are zero and correct Mono or stereo Det is filled.
358 
359  if ((isMono && isStereo)||(!isMono&&!isStereo)) throw cms::Exception("FastTrackerRecHitMatcher") << "Something wrong with DetIds.";
360  return std::make_unique<FastProjectedTrackerRecHit>(projectedHitPos, rotatedError, *gluedDet, *originalRH);
361 }
362 
363 
364 
365 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
367  //The following says we do not know what parameters are allowed so do no validation
368  // Please change this to state exactly what you do use, even if it is no parameters
370  desc.setUnknown();
371  descriptions.addDefault(desc);
372 }
373 
374 //define this as a plug-in
dbl * delta
Definition: mlp_gen.cc:36
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:106
T getParameter(std::string const &) const
float xx() const
Definition: LocalError.h:24
T mag2() const
Definition: PV3DBase.h:66
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:32
LocalError localPositionError() const final
const GeomDetUnit * monoDet() const
Definition: GluedGeomDet.h:20
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
virtual const Topology & topology() const
Definition: GeomDet.cc:81
std::unique_ptr< FastTrackerRecHit > projectOnly(const FastSingleTrackerRecHit *originalRH, const GeomDet *monoDet, const GluedGeomDet *gluedDet, LocalVector &ldir) const
void produce(edm::Event &, const edm::EventSetup &) override
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
T y() const
Definition: PV3DBase.h:63
std::unique_ptr< FastTrackerRecHit > match(const FastSingleTrackerRecHit *monoRH, const FastSingleTrackerRecHit *stereoRH, const GluedGeomDet *gluedDet, LocalVector &trackdirection, bool stereLayerFirst) const
LocalVector toLocal(const reco::Track::Vector &v, const Surface &s)
const Plane & surface() const
The nominal surface of the GeomDet.
Definition: GeomDet.h:42
edm::Ref< FastTrackerRecHitCollection > FastTrackerRecHitRef
FastTrackerRecHitMatcher(const edm::ParameterSet &)
ROOT::Math::SMatrix< double, 2, 2, ROOT::Math::MatRepStd< double, 2, 2 > > AlgebraicMatrix22
float xy() const
Definition: LocalError.h:25
Measurement2DPoint MeasurementPoint
Measurement points are two-dimensional by default.
int iEvent
Definition: GenABIO.cc:224
unsigned int glued() const
glued
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void addDefault(ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< FastTrackerRecHitRefCollection > simHit2RecHitMapToken
float yy() const
Definition: LocalError.h:26
virtual LocalPoint localPosition(float strip) const =0
const GeomDet * det() const
LocalPoint toLocal(const GlobalPoint &gp) const
T z() const
Definition: PV3DBase.h:64
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
virtual MeasurementPoint measurementPosition(const LocalPoint &) const =0
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
edm::EDGetTokenT< edm::PSimHitContainer > simHitsToken
DetId geographicalId() const
The label of this GeomDet.
Definition: GeomDet.h:79
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:243
StripPosition project(const GeomDetUnit *det, const GluedGeomDet *glueddet, const StripPosition &strip, const LocalVector &trackdirection) const
bool isNull() const
Checks for null.
Definition: Ref.h:248
RefProd< PROD > getRefBeforePut()
Definition: Event.h:150
LocalVector localDirection() const
Obsolete. Same as momentumAtEntry().unit(), for backward compatibility.
Definition: PSimHit.h:58
Definition: DetId.h:18
unsigned int stereo() const
stereo
const FastSingleTrackerRecHit * _cast2Single(const FastTrackerRecHit *recHit) const
virtual float localPitch(const LocalPoint &) const =0
std::pair< LocalPoint, LocalPoint > StripPosition
ESHandle< TrackerGeometry > geometry
static int position[264][3]
Definition: ReadPGInfo.cc:509
T get() const
Definition: EventSetup.h:71
const TrackerGeomDet * idToDet(DetId) const override
unsigned int trackId() const
Definition: PSimHit.h:106
trackerHitRTTI::RTTI rtti() const
LocalPoint localPosition() const final
DetId geographicalId() const
virtual MeasurementError measurementError(const LocalPoint &, const LocalError &) const =0
LocalError rotate(float x, float y) const
Return a new LocalError, rotated by an angle defined by the direction (x,y)
Definition: LocalError.h:39
T x() const
Definition: PV2DBase.h:45
T x() const
Definition: PV3DBase.h:62
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ROOT::Math::SVector< double, 2 > AlgebraicVector2
def move(src, dest)
Definition: eostools.py:511
const GeomDetUnit * stereoDet() const
Definition: GluedGeomDet.h:21
bool isSingle() const