CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_5_3_14/src/FastSimulation/EgammaElectronAlgos/src/FastPixelHitMatcher.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:    EgammaElectronAlgos
00004 // Class:      FastPixelHitMatcher
00005 // 
00013 //
00014 // Original Author: Patrick Janot
00015 //
00016 //
00017 
00018 #include "FastSimulation/EgammaElectronAlgos/interface/FastPixelHitMatcher.h"
00019 #include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00020 #include "FastSimulation/TrackerSetup/interface/TrackerInteractionGeometry.h"
00021 #include "FastSimulation/TrackerSetup/interface/TrackerLayer.h"
00022 #include "FastSimulation/Tracking/interface/TrackerRecHit.h"
00023 #include "FastSimulation/ParticlePropagator/interface/ParticlePropagator.h"
00024 #include "DataFormats/SiPixelDetId/interface/PXBDetId.h"
00025 #include "DataFormats/SiPixelDetId/interface/PXFDetId.h"
00026 
00027 #include "DataFormats/DetId/interface/DetId.h"
00028 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00029 
00030 //#define FAMOS_DEBUG
00031 
00032 FastPixelHitMatcher::FastPixelHitMatcher(float ephi1min, float ephi1max, 
00033                                      float pphi1min, float pphi1max, 
00034                                      float phi2min, float phi2max, 
00035                                      float z2minB, float z2maxB,
00036                                      float r2minF, float r2maxF,
00037                                      float rMinI, float rMaxI, 
00038                                      bool searchInTIDTEC) :
00039   ephi1min(ephi1min), ephi1max(ephi1max), 
00040   pphi1min(pphi1min), pphi1max(pphi1max), 
00041   phi2min(phi2min), phi2max(phi2max), 
00042   z2minB(z2minB), z2maxB(z2maxB), 
00043   r2minF(r2minF), r2maxF(r2maxF),
00044   rMinI(rMinI), rMaxI(rMaxI), 
00045   searchInTIDTEC(searchInTIDTEC),
00046   theTrackerGeometry(0),
00047   theMagneticField(0),
00048   theGeomSearchTracker(0),
00049   _theGeometry(0),
00050   thePixelLayers(50,static_cast<TrackerLayer*>(0)),
00051   vertex(0.) {}
00052 
00053 FastPixelHitMatcher::~FastPixelHitMatcher() { }
00054 
00055 void 
00056 FastPixelHitMatcher::setES(const MagneticFieldMap* aFieldMap, 
00057                          const TrackerGeometry* aTrackerGeometry, 
00058                          const GeometricSearchTracker* geomSearchTracker,
00059                          const TrackerInteractionGeometry* interactionGeometry) {
00060   
00061   // initialize the tracker geometry and the magnetic field map
00062   theTrackerGeometry = aTrackerGeometry; 
00063   //theMagneticField = aMagField;
00064   theGeomSearchTracker = geomSearchTracker;
00065   _theGeometry = interactionGeometry;
00066   theFieldMap = aFieldMap;
00067   
00068   // Initialize (if not already done) the simplified magnetic field geometry
00069   // MagneticFieldMap::instance( theMagneticField, _theGeometry );
00070   
00071   // The pixel layers in the simplified geometry 
00072   unsigned layer = 1;
00073   std::list<TrackerLayer>::const_iterator cyliter = _theGeometry->cylinderBegin();
00074   for ( ; cyliter != _theGeometry->cylinderEnd() ; ++cyliter ) {
00075     if ( layer != cyliter->layerNumber() ) continue;
00076     thePixelLayers[layer++] = &(*cyliter);
00077   }
00078   
00079 }
00080 
00081 std::vector< std::pair<FastPixelHitMatcher::ConstRecHitPointer, 
00082                        FastPixelHitMatcher::ConstRecHitPointer> > 
00083 FastPixelHitMatcher::compatibleHits(const GlobalPoint& thePos,
00084                                   const GlobalPoint& theVertex,
00085                                   float energy,
00086                                   std::vector<TrackerRecHit>& theHits) { 
00087   
00088   std::vector< std::pair<FastPixelHitMatcher::ConstRecHitPointer, 
00089     FastPixelHitMatcher::ConstRecHitPointer> > result;
00090 #ifdef FAMOS_DEBUG
00091   std::cout << "[FastPixelHitMatcher::compatibleHits] entering .. " << std::endl;
00092 #endif
00093   
00094   double zCluster = thePos.z();
00095   double rCluster = thePos.perp();
00096   
00097   // The cluster inferred energy-momentum
00098   double theLength = thePos.mag();
00099   XYZTLorentzVector theMom(thePos.x(), thePos.y(), thePos.z(), theLength);
00100   theMom *= energy / theLength;
00101   XYZTLorentzVector theVert(thePos.x(),thePos.y(),thePos.z(),0.);
00102   XYZTLorentzVector theNominalVertex(theVertex.x(), theVertex.y(), theVertex.z(), 0.);
00103   
00104   // The corresponding RawParticles (to be propagated for e- and e+
00105   ParticlePropagator myElec(theMom,theVert,-1.,theFieldMap);
00106   ParticlePropagator myPosi(theMom,theVert,+1.,theFieldMap); 
00107 #ifdef FAMOS_DEBUG
00108   std::cout << "elec/posi before propagation " << std::endl << myElec << std::endl << myPosi << std::endl;
00109 #endif
00110   
00111   // Propagate the e- and the e+ hypothesis to the nominal vertex
00112   // by modifying the pT direction in an appropriate manner.
00113   myElec.propagateToNominalVertex(theNominalVertex);
00114   myPosi.propagateToNominalVertex(theNominalVertex);
00115 #ifdef FAMOS_DEBUG
00116   std::cout << "elec/posi after propagation " << std::endl << myElec << std::endl << myPosi << std::endl;
00117 #endif
00118   
00119   // Look for an appropriate see in the pixel detector
00120   bool thereIsASeed = false;
00121   unsigned nHits = theHits.size();
00122   
00123   for ( unsigned firstHit=0; firstHit<nHits-1; ++firstHit ) { 
00124     for ( unsigned secondHit=firstHit+1; secondHit<nHits; ++secondHit ) {      
00125 
00126       // Is there a seed associated to this pair of Pixel hits?
00127       thereIsASeed = isASeed(myElec,myPosi,theVertex,
00128                              rCluster,zCluster,
00129                              theHits[firstHit],
00130                              theHits[secondHit]);
00131 
00132 #ifdef FAMOS_DEBUG
00133       std::cout << "Is there a seed with hits " << firstHit << " & "<< secondHit << "? " << thereIsASeed << std::endl;
00134 #endif
00135       if ( !thereIsASeed ) continue;
00136       
00137       ConstRecHitPointer theFirstHit = 
00138         GenericTransientTrackingRecHit::build(theHits[firstHit].geomDet(),
00139                                               theHits[firstHit].hit());
00140       ConstRecHitPointer theSecondHit = 
00141         GenericTransientTrackingRecHit::build(theHits[secondHit].geomDet(),
00142                                               theHits[secondHit].hit());
00143       result.push_back(std::pair<
00144                        FastPixelHitMatcher::ConstRecHitPointer,
00145                        FastPixelHitMatcher::ConstRecHitPointer>(theFirstHit,theSecondHit));
00146       
00147     }
00148   }
00149   
00150   return result;
00151 }
00152 
00153 bool FastPixelHitMatcher::isASeed(const ParticlePropagator& myElec,
00154                                 const ParticlePropagator& myPosi,
00155                                 const GlobalPoint& theVertex,
00156                                 double rCluster,
00157                                 double zCluster,
00158                                 const TrackerRecHit& hit1,
00159                                 const TrackerRecHit& hit2) {
00160   
00161   // Check that the two hits are not on the same layer
00162   if ( hit2.isOnTheSameLayer(hit1) ) return false;
00163 #ifdef FAMOS_DEBUG
00164   std::cout << "isASeed: The two hits are not on the same layer - OK " << std::endl;
00165 #endif
00166 
00167   // Check that the first hit is on PXB or PXD
00168   if ( hit1.subDetId() > 2 ) return false;
00169 #ifdef FAMOS_DEBUG
00170   std::cout << "isASeed: The first hits is on the pixel detector " << std::endl;
00171 #endif
00172 
00173   // Impose the track to originate from zVertex = 0. and check the 
00174   // compatibility with the first hit (beam spot constraint)
00175   GlobalPoint firstHit = hit1.globalPosition();
00176   bool rzok = checkRZCompatibility(zCluster, rCluster, 0., z1min, z1max, firstHit, hit1.subDetId()>1);
00177 #ifdef FAMOS_DEBUG
00178   std::cout << "isASeed: rzok (1) = " << rzok << std::endl;
00179 #endif
00180   if ( !rzok ) return false;
00181   
00182   // Refine the Z vertex by imposing the track to pass through the first RecHit, 
00183   // and check the compatibility with the second rechit 
00184   GlobalPoint secondHit = hit2.globalPosition(); 
00185   rzok = false;
00186 
00187   // The orgin Z vertex for thet track passing through the first rechit
00188   vertex = zVertex(zCluster, rCluster, firstHit);
00189 
00190   // Compute R (forward) or Z (barrel) predicted for the second hit and check compatibility
00191   if ( hit2.subDetId() == 1 ) {
00192     rzok = checkRZCompatibility(zCluster, rCluster, vertex, z2minB, z2maxB, secondHit, false);
00193   } else if ( hit2.subDetId() == 2 ) {  
00194     rzok = checkRZCompatibility(zCluster, rCluster, vertex, r2minF, r2maxF, secondHit, true);
00195   } else { 
00196     rzok = checkRZCompatibility(zCluster, rCluster, vertex, rMinI, rMaxI, secondHit, true);
00197   }
00198 #ifdef FAMOS_DEBUG
00199   std::cout << "isASeed: rzok (2) = " << rzok << std::endl;
00200 #endif
00201   if ( !rzok ) return false;
00202 
00203   // Propagate the inferred electron (positron) to the first layer,
00204   // check the compatibility with the first hit, and propagate back
00205   // to the nominal vertex with the hit constraint
00206   ParticlePropagator elec(myElec);
00207   ParticlePropagator posi(myPosi);
00208 #ifdef FAMOS_DEBUG
00209   std::cout << "isASeed: elec1 to be propagated to first layer" << std::endl;
00210 #endif
00211   bool elec1 = propagateToLayer(elec,theVertex,firstHit,
00212                                 ephi1min,ephi1max,hit1.cylinderNumber());
00213 #ifdef FAMOS_DEBUG
00214   std::cout << "isASeed: posi1 to be propagated to first layer" << std::endl;
00215 #endif
00216   bool posi1 = propagateToLayer(posi,theVertex,firstHit,
00217                                 pphi1min,pphi1max,hit1.cylinderNumber());
00218 
00219 #ifdef FAMOS_DEBUG
00220   std::cout << "isASeed: elec1 / posi1 " << elec1 << " " << posi1 << std::endl;
00221 #endif
00222   // Neither the electron not the positron hypothesis work...
00223   if ( !elec1 && !posi1 ) return false;
00224   
00225   // Otherwise, propagate to the second layer, check the compatibility
00226   // with the second hit and propagate back to the nominal vertex with 
00227   // the hit constraint
00228 #ifdef FAMOS_DEBUG
00229   std::cout << "isASeed: elec2 to be propagated to second layer" << std::endl;
00230 #endif
00231   bool elec2 = elec1 && propagateToLayer(elec,theVertex,secondHit,
00232                                          phi2min,phi2max,hit2.cylinderNumber());
00233   
00234 #ifdef FAMOS_DEBUG
00235   std::cout << "isASeed: posi2 to be propagated to second layer" << std::endl;
00236 #endif
00237   bool posi2 = posi1 && propagateToLayer(posi,theVertex,secondHit,
00238                                          phi2min,phi2max,hit2.cylinderNumber());
00239   
00240 #ifdef FAMOS_DEBUG
00241   std::cout << "isASeed: elec2 / posi2 " << elec2 << " " << posi2 << std::endl;
00242 #endif
00243   if ( !elec2 && !posi2 ) return false;
00244 
00245   return true;
00246 
00247 }
00248   
00249 bool 
00250 FastPixelHitMatcher::propagateToLayer(ParticlePropagator& myPart,
00251                                     const GlobalPoint& theVertex,
00252                                     GlobalPoint& theHit,
00253                                     double phimin, double phimax,
00254                                     unsigned layer) {
00255 
00256   // Set the z position of the particle to the predicted one
00257   XYZTLorentzVector theNominalVertex(theVertex.x(), theVertex.y(), vertex, 0.);
00258   myPart.setVertex(theNominalVertex);
00259 #ifdef FAMOS_DEBUG
00260   std::cout << "propagateToLayer: propagateToLayer: Before propagation (0) " << myPart << std::endl;
00261 #endif
00262 
00263   // Propagate the inferred electron (positron) to the first layer
00264   // Use the radius (barrel) or the z (forward) of the hit instead 
00265   // of the inaccurate layer radius and z.
00266   double rCyl = thePixelLayers[layer]->forward() ? 999. : theHit.perp();
00267   double zCyl = thePixelLayers[layer]->forward() ? fabs(theHit.z()) : 999.;
00268   BaseParticlePropagator* myBasePart = (BaseParticlePropagator*)(&myPart);
00269   myBasePart->setPropagationConditions(rCyl,zCyl);
00270   // myPart.setPropagationConditions(*(thePixelLayers[layer]));
00271 
00272   // bool success = myPart.propagateToBoundSurface(*(thePixelLayers[layer]));
00273   bool success = myPart.propagate();
00274 #ifdef FAMOS_DEBUG
00275   std::cout << "propagateToLayer: Success ? " << success << std::endl;
00276   std::cout << "propagateToLayer: After  propagation (1) " << myPart << std::endl;
00277   std::cout << "propagateToLayer: The hit               " << theHit << std::endl; 
00278 #endif
00279       
00280   // Check that propagated particle is within the proper phi range.
00281   if ( success ) {
00282     double dphi = myPart.vertex().phi() - theHit.phi();
00283     if ( dphi < -M_PI ) 
00284       dphi = dphi + 2.*M_PI;
00285     else if ( dphi > M_PI ) 
00286       dphi = dphi - 2.*M_PI;
00287 #ifdef FAMOS_DEBUG
00288     std::cout << "propagateToLayer: Phi range ? " << phimin << " < " << dphi << " < " << phimax << std::endl; 
00289 #endif
00290     if ( dphi < phimin || dphi > phimax ) success = false;
00291   }
00292       
00293   // Impose the track to go through the hit and propagate back to 
00294   // the nominal vertex
00295   if ( success ) {
00296     myPart.setVertex( XYZTLorentzVector(theHit.x(), theHit.y(), theHit.z(), 0.) );
00297     myPart.propagateToNominalVertex(theNominalVertex);
00298 #ifdef FAMOS_DEBUG
00299     std::cout << "propagateToLayer: After  propagation (2) " << myPart << std::endl;
00300 #endif
00301   }
00302 
00303   return success;
00304 
00305 }
00306                                          
00307 
00308 double 
00309 FastPixelHitMatcher::zVertex(double zCluster,
00310                            double rCluster,
00311                            GlobalPoint& theHit)
00312 {
00313 
00314   // Refine the Z vertex by imposing the track to pass through the RecHit
00315   double pxHitz = theHit.z();
00316   double pxHitr = theHit.perp();
00317   return pxHitz - pxHitr*(zCluster-pxHitz)/(rCluster-pxHitr);
00318 
00319 }
00320 
00321 
00322 bool
00323 FastPixelHitMatcher::checkRZCompatibility(double zCluster,double rCluster, 
00324                                         double zVertex, 
00325                                         float rzMin, float rzMax,
00326                                         GlobalPoint& theHit, 
00327                                         bool forward) 
00328 {
00329 
00330   // The hit position
00331   double zHit = theHit.z();
00332   double rHit = theHit.perp();
00333 
00334   // Compute the intersection of a line joining the cluster position 
00335   // and the predicted z origin (zVertex) with the layer hit 
00336   // (returns R for a forward layer, and Z for a barrel layer)
00337   double checkRZ = forward ?
00338     (zHit-zVertex)/(zCluster-zVertex) * rCluster
00339     :
00340     zVertex + rHit * (zCluster-zVertex)/rCluster;
00341 
00342   // This value is then compared to check with the actual hit position 
00343   // (in R for a forward layer, in Z for a barrel layer)
00344 
00345   return forward ?
00346     checkRZ+rzMin < rHit && rHit < checkRZ+rzMax 
00347     :
00348     checkRZ+rzMin < zHit && zHit < checkRZ+rzMax;
00349     
00350 }
00351 
00352 
00353 
00354 
00355 
00356 
00357