test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastPixelHitMatcher.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: EgammaElectronAlgos
4 // Class: FastPixelHitMatcher
5 //
13 //
14 // Original Author: Patrick Janot
15 //
16 //
17 
18 #include "FastPixelHitMatcher.h"
23 
26 
27 //#define FAMOS_DEBUG
28 
29 FastPixelHitMatcher::FastPixelHitMatcher(float ephi1min, float ephi1max,
30  float pphi1min, float pphi1max,
31  float phi2min, float phi2max,
32  float z2minB, float z2maxB,
33  float r2minF, float r2maxF,
34  float rMinI, float rMaxI,
35  bool searchInTIDTEC) :
36  ephi1min(ephi1min), ephi1max(ephi1max),
37  pphi1min(pphi1min), pphi1max(pphi1max),
38  phi2min(phi2min), phi2max(phi2max),
39  z2minB(z2minB), z2maxB(z2maxB),
40  r2minF(r2minF), r2maxF(r2maxF),
41  rMinI(rMinI), rMaxI(rMaxI),
42  searchInTIDTEC(searchInTIDTEC),
43  theTrackerGeometry(0),
44  theMagneticField(0),
45  theGeomSearchTracker(0),
46  _theGeometry(0),
47  thePixelLayers(50,static_cast<TrackerLayer*>(0)),
48  vertex(0.) {}
49 
51 
52 void
54  const TrackerGeometry* aTrackerGeometry,
55  const GeometricSearchTracker* geomSearchTracker,
56  const TrackerInteractionGeometry* interactionGeometry) {
57 
58  // initialize the tracker geometry and the magnetic field map
59  theTrackerGeometry = aTrackerGeometry;
60  //theMagneticField = aMagField;
61  theGeomSearchTracker = geomSearchTracker;
62  _theGeometry = interactionGeometry;
63  theFieldMap = aFieldMap;
64 
65  // Initialize (if not already done) the simplified magnetic field geometry
66  // MagneticFieldMap::instance( theMagneticField, _theGeometry );
67 
68  // The pixel layers in the simplified geometry
69  unsigned layer = 1;
70  std::list<TrackerLayer>::const_iterator cyliter = _theGeometry->cylinderBegin();
71  for ( ; cyliter != _theGeometry->cylinderEnd() ; ++cyliter ) {
72  if ( layer != cyliter->layerNumber() ) continue;
73  thePixelLayers[layer++] = &(*cyliter);
74  }
75 
76 }
77 
78 std::vector< std::pair<FastPixelHitMatcher::ConstRecHitPointer,
81  const GlobalPoint& theVertex,
82  float energy,
83  std::vector<TrajectorySeedHitCandidate>& theHits) {
84 
85  std::vector< std::pair<FastPixelHitMatcher::ConstRecHitPointer,
86  FastPixelHitMatcher::ConstRecHitPointer> > result;
87 #ifdef FAMOS_DEBUG
88  std::cout << "[FastPixelHitMatcher::compatibleHits] entering .. " << std::endl;
89 #endif
90 
91  double zCluster = thePos.z();
92  double rCluster = thePos.perp();
93 
94  // The cluster inferred energy-momentum
95  double theLength = thePos.mag();
96  XYZTLorentzVector theMom(thePos.x(), thePos.y(), thePos.z(), theLength);
97  theMom *= energy / theLength;
98  XYZTLorentzVector theVert(thePos.x(),thePos.y(),thePos.z(),0.);
99  XYZTLorentzVector theNominalVertex(theVertex.x(), theVertex.y(), theVertex.z(), 0.);
100 
101  // The corresponding RawParticles (to be propagated for e- and e+
102  ParticlePropagator myElec(theMom,theVert,-1.,theFieldMap);
103  ParticlePropagator myPosi(theMom,theVert,+1.,theFieldMap);
104 #ifdef FAMOS_DEBUG
105  std::cout << "elec/posi before propagation " << std::endl << myElec << std::endl << myPosi << std::endl;
106 #endif
107 
108  // Propagate the e- and the e+ hypothesis to the nominal vertex
109  // by modifying the pT direction in an appropriate manner.
110  myElec.propagateToNominalVertex(theNominalVertex);
111  myPosi.propagateToNominalVertex(theNominalVertex);
112 #ifdef FAMOS_DEBUG
113  std::cout << "elec/posi after propagation " << std::endl << myElec << std::endl << myPosi << std::endl;
114 #endif
115 
116  // Look for an appropriate see in the pixel detector
117  bool thereIsASeed = false;
118  unsigned nHits = theHits.size();
119 
120  for ( unsigned firstHit=0; firstHit<nHits-1; ++firstHit ) {
121  for ( unsigned secondHit=firstHit+1; secondHit<nHits; ++secondHit ) {
122 
123  // Is there a seed associated to this pair of Pixel hits?
124  thereIsASeed = isASeed(myElec,myPosi,theVertex,
125  rCluster,zCluster,
126  theHits[firstHit],
127  theHits[secondHit]);
128 
129 #ifdef FAMOS_DEBUG
130  std::cout << "Is there a seed with hits " << firstHit << " & "<< secondHit << "? " << thereIsASeed << std::endl;
131 #endif
132  if ( !thereIsASeed ) continue;
133 
134  ConstRecHitPointer theFirstHit =
135  GenericTransientTrackingRecHit::build(theHits[firstHit].geomDet(),
136  theHits[firstHit].hit());
137  ConstRecHitPointer theSecondHit =
138  GenericTransientTrackingRecHit::build(theHits[secondHit].geomDet(),
139  theHits[secondHit].hit());
140  result.push_back(std::pair<
141  FastPixelHitMatcher::ConstRecHitPointer,
142  FastPixelHitMatcher::ConstRecHitPointer>(theFirstHit,theSecondHit));
143 
144  }
145  }
146 
147  return result;
148 }
149 
151  const ParticlePropagator& myPosi,
152  const GlobalPoint& theVertex,
153  double rCluster,
154  double zCluster,
155  const TrajectorySeedHitCandidate& hit1,
156  const TrajectorySeedHitCandidate& hit2) {
157 
158  // Check that the two hits are not on the same layer
159  if ( hit2.isOnTheSameLayer(hit1) ) return false;
160 #ifdef FAMOS_DEBUG
161  std::cout << "isASeed: The two hits are not on the same layer - OK " << std::endl;
162 #endif
163 
164  // Check that the first hit is on PXB or PXD
165  if ( hit1.subDetId() > 2 ) return false;
166 #ifdef FAMOS_DEBUG
167  std::cout << "isASeed: The first hits is on the pixel detector " << std::endl;
168 #endif
169 
170  // Impose the track to originate from zVertex = 0. and check the
171  // compatibility with the first hit (beam spot constraint)
172  GlobalPoint firstHit = hit1.globalPosition();
173  bool rzok = checkRZCompatibility(zCluster, rCluster, 0., z1min, z1max, firstHit, hit1.subDetId()>1);
174 #ifdef FAMOS_DEBUG
175  std::cout << "isASeed: rzok (1) = " << rzok << std::endl;
176 #endif
177  if ( !rzok ) return false;
178 
179  // Refine the Z vertex by imposing the track to pass through the first RecHit,
180  // and check the compatibility with the second rechit
181  GlobalPoint secondHit = hit2.globalPosition();
182  rzok = false;
183 
184  // The orgin Z vertex for thet track passing through the first rechit
185  vertex = zVertex(zCluster, rCluster, firstHit);
186 
187  // Compute R (forward) or Z (barrel) predicted for the second hit and check compatibility
188  if ( hit2.subDetId() == 1 ) {
189  rzok = checkRZCompatibility(zCluster, rCluster, vertex, z2minB, z2maxB, secondHit, false);
190  } else if ( hit2.subDetId() == 2 ) {
191  rzok = checkRZCompatibility(zCluster, rCluster, vertex, r2minF, r2maxF, secondHit, true);
192  } else {
193  rzok = checkRZCompatibility(zCluster, rCluster, vertex, rMinI, rMaxI, secondHit, true);
194  }
195 #ifdef FAMOS_DEBUG
196  std::cout << "isASeed: rzok (2) = " << rzok << std::endl;
197 #endif
198  if ( !rzok ) return false;
199 
200  // Propagate the inferred electron (positron) to the first layer,
201  // check the compatibility with the first hit, and propagate back
202  // to the nominal vertex with the hit constraint
203  ParticlePropagator elec(myElec);
204  ParticlePropagator posi(myPosi);
205 #ifdef FAMOS_DEBUG
206  std::cout << "isASeed: elec1 to be propagated to first layer" << std::endl;
207 #endif
208  bool elec1 = propagateToLayer(elec,theVertex,firstHit,
210 #ifdef FAMOS_DEBUG
211  std::cout << "isASeed: posi1 to be propagated to first layer" << std::endl;
212 #endif
213  bool posi1 = propagateToLayer(posi,theVertex,firstHit,
215 
216 #ifdef FAMOS_DEBUG
217  std::cout << "isASeed: elec1 / posi1 " << elec1 << " " << posi1 << std::endl;
218 #endif
219  // Neither the electron not the positron hypothesis work...
220  if ( !elec1 && !posi1 ) return false;
221 
222  // Otherwise, propagate to the second layer, check the compatibility
223  // with the second hit and propagate back to the nominal vertex with
224  // the hit constraint
225 #ifdef FAMOS_DEBUG
226  std::cout << "isASeed: elec2 to be propagated to second layer" << std::endl;
227 #endif
228  bool elec2 = elec1 && propagateToLayer(elec,theVertex,secondHit,
230 
231 #ifdef FAMOS_DEBUG
232  std::cout << "isASeed: posi2 to be propagated to second layer" << std::endl;
233 #endif
234  bool posi2 = posi1 && propagateToLayer(posi,theVertex,secondHit,
236 
237 #ifdef FAMOS_DEBUG
238  std::cout << "isASeed: elec2 / posi2 " << elec2 << " " << posi2 << std::endl;
239 #endif
240  if ( !elec2 && !posi2 ) return false;
241 
242  return true;
243 
244 }
245 
246 bool
248  const GlobalPoint& theVertex,
249  GlobalPoint& theHit,
250  double phimin, double phimax,
251  unsigned layer) {
252 
253  // Set the z position of the particle to the predicted one
254  XYZTLorentzVector theNominalVertex(theVertex.x(), theVertex.y(), vertex, 0.);
255  myPart.setVertex(theNominalVertex);
256 #ifdef FAMOS_DEBUG
257  std::cout << "propagateToLayer: propagateToLayer: Before propagation (0) " << myPart << std::endl;
258 #endif
259 
260  // Propagate the inferred electron (positron) to the first layer
261  // Use the radius (barrel) or the z (forward) of the hit instead
262  // of the inaccurate layer radius and z.
263  double rCyl = thePixelLayers[layer]->forward() ? 999. : theHit.perp();
264  double zCyl = thePixelLayers[layer]->forward() ? fabs(theHit.z()) : 999.;
265  BaseParticlePropagator* myBasePart = (BaseParticlePropagator*)(&myPart);
266  myBasePart->setPropagationConditions(rCyl,zCyl);
267  // myPart.setPropagationConditions(*(thePixelLayers[layer]));
268 
269  // bool success = myPart.propagateToBoundSurface(*(thePixelLayers[layer]));
270  bool success = myPart.propagate();
271 #ifdef FAMOS_DEBUG
272  std::cout << "propagateToLayer: Success ? " << success << std::endl;
273  std::cout << "propagateToLayer: After propagation (1) " << myPart << std::endl;
274  std::cout << "propagateToLayer: The hit " << theHit << std::endl;
275 #endif
276 
277  // Check that propagated particle is within the proper phi range.
278  if ( success ) {
279  double dphi = myPart.vertex().phi() - theHit.phi();
280  if ( dphi < -M_PI )
281  dphi = dphi + 2.*M_PI;
282  else if ( dphi > M_PI )
283  dphi = dphi - 2.*M_PI;
284 #ifdef FAMOS_DEBUG
285  std::cout << "propagateToLayer: Phi range ? " << phimin << " < " << dphi << " < " << phimax << std::endl;
286 #endif
287  if ( dphi < phimin || dphi > phimax ) success = false;
288  }
289 
290  // Impose the track to go through the hit and propagate back to
291  // the nominal vertex
292  if ( success ) {
293  myPart.setVertex( XYZTLorentzVector(theHit.x(), theHit.y(), theHit.z(), 0.) );
294  myPart.propagateToNominalVertex(theNominalVertex);
295 #ifdef FAMOS_DEBUG
296  std::cout << "propagateToLayer: After propagation (2) " << myPart << std::endl;
297 #endif
298  }
299 
300  return success;
301 
302 }
303 
304 
305 double
307  double rCluster,
308  GlobalPoint& theHit)
309 {
310 
311  // Refine the Z vertex by imposing the track to pass through the RecHit
312  double pxHitz = theHit.z();
313  double pxHitr = theHit.perp();
314  return pxHitz - pxHitr*(zCluster-pxHitz)/(rCluster-pxHitr);
315 
316 }
317 
318 
319 bool
320 FastPixelHitMatcher::checkRZCompatibility(double zCluster,double rCluster,
321  double zVertex,
322  float rzMin, float rzMax,
323  GlobalPoint& theHit,
324  bool forward)
325 {
326 
327  // The hit position
328  double zHit = theHit.z();
329  double rHit = theHit.perp();
330 
331  // Compute the intersection of a line joining the cluster position
332  // and the predicted z origin (zVertex) with the layer hit
333  // (returns R for a forward layer, and Z for a barrel layer)
334  double checkRZ = forward ?
335  (zHit-zVertex)/(zCluster-zVertex) * rCluster
336  :
337  zVertex + rHit * (zCluster-zVertex)/rCluster;
338 
339  // This value is then compared to check with the actual hit position
340  // (in R for a forward layer, in Z for a barrel layer)
341 
342  return forward ?
343  checkRZ+rzMin < rHit && rHit < checkRZ+rzMax
344  :
345  checkRZ+rzMin < zHit && zHit < checkRZ+rzMax;
346 
347 }
348 
349 
350 
351 
352 
353 
354 
T perp() const
Definition: PV3DBase.h:72
unsigned int cylinderNumber() const
The global layer number in the nested cylinder geometry.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
T y() const
Definition: PV3DBase.h:63
std::list< TrackerLayer >::const_iterator cylinderEnd() const
Returns the last pointer in the cylinder list.
void setPropagationConditions(double r, double z, bool firstLoop=true)
Set the propagation characteristics (rCyl, zCyl and first loop only)
GlobalPoint globalPosition() const
The global position.
std::vector< const TrackerLayer * > thePixelLayers
double zVertex(double zCluster, double rCluster, GlobalPoint &theHit)
static RecHitPointer build(const GeomDet *geom, const TrackingRecHit *rh)
T mag() const
Definition: PV3DBase.h:67
TransientTrackingRecHit::ConstRecHitPointer ConstRecHitPointer
T z() const
Definition: PV3DBase.h:64
tuple result
Definition: query.py:137
bool isASeed(const ParticlePropagator &myElec, const ParticlePropagator &myPosi, const GlobalPoint &theVertex, double rCluster, double zCluster, const TrajectorySeedHitCandidate &hit1, const TrajectorySeedHitCandidate &hit2)
#define M_PI
bool propagateToLayer(ParticlePropagator &myPart, const GlobalPoint &theVertex, GlobalPoint &theHit, double phimin, double phimax, unsigned layer)
const GeometricSearchTracker * theGeomSearchTracker
const MagneticFieldMap * theFieldMap
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
std::vector< std::pair< ConstRecHitPointer, ConstRecHitPointer > > compatibleHits(const GlobalPoint &xmeas, const GlobalPoint &vprim, float energy, std::vector< TrajectorySeedHitCandidate > &theHits)
FastPixelHitMatcher(float, float, float, float, float, float, float, float, float, float, float, float, bool)
const TrackerGeometry * theTrackerGeometry
std::list< TrackerLayer >::const_iterator cylinderBegin() const
Returns the first pointer in the cylinder list.
const TrackerInteractionGeometry * _theGeometry
bool checkRZCompatibility(double zCluster, double rCluster, double zVertex, float rzMin, float rzMax, GlobalPoint &theHit, bool forward)
bool propagateToNominalVertex(const XYZTLorentzVector &hit2=XYZTLorentzVector(0., 0., 0., 0.))
void setES(const MagneticFieldMap *aFieldMap, const TrackerGeometry *aTrackerGeometry, const GeometricSearchTracker *geomSearchTracker, const TrackerInteractionGeometry *interactionGeometry)
tuple cout
Definition: gather_cfg.py:121
bool isOnTheSameLayer(const TrajectorySeedHitCandidate &other) const
Check if two hits are on the same layer of the same subdetector.
unsigned int subDetId() const
The subdet Id.
T x() const
Definition: PV3DBase.h:62
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:288
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15