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