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