00001
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003
00004
00005 #include "DataFormats/GeometrySurface/interface/BoundDisk.h"
00006 #include "DataFormats/GeometrySurface/interface/BoundCylinder.h"
00007 #include "DataFormats/GeometrySurface/interface/Surface.h"
00008 #include "DataFormats/GeometrySurface/interface/TangentPlane.h"
00009
00010
00011 #include "TrackingTools/DetLayers/interface/DetLayer.h"
00012 #include "TrackingTools/DetLayers/interface/BarrelDetLayer.h"
00013 #include "TrackingTools/DetLayers/interface/ForwardDetLayer.h"
00014 #include "TrackingTools/GeomPropagators/interface/AnalyticalPropagator.h"
00015 #include "FastSimulation/TrajectoryManager/interface/InsideBoundsMeasurementEstimator.h"
00016 #include "Geometry/CommonDetUnit/interface/GeomDetUnit.h"
00017 #include "TrackingTools/GeomPropagators/interface/HelixArbitraryPlaneCrossing.h"
00018 #include "RecoTracker/TkDetLayers/interface/GeometricSearchTracker.h"
00019
00020
00021
00022 #include "FastSimulation/TrajectoryManager/interface/TrajectoryManager.h"
00023 #include "FastSimulation/TrajectoryManager/interface/LocalMagneticField.h"
00024 #include "FastSimulation/ParticlePropagator/interface/ParticlePropagator.h"
00025 #include "FastSimulation/TrackerSetup/interface/TrackerInteractionGeometry.h"
00026 #include "FastSimulation/ParticleDecay/interface/PythiaDecays.h"
00027 #include "FastSimulation/Event/interface/FSimEvent.h"
00028 #include "FastSimulation/Event/interface/FSimVertex.h"
00029 #include "FastSimulation/Event/interface/KineParticleFilter.h"
00030
00031 #include "FastSimulation/Utilities/interface/RandomEngine.h"
00032
00033
00034
00035
00036
00037
00038 #ifdef FAMOS_DEBUG
00039 #include "DataFormats/TrackerCommon/interface/TrackerTopology.h"
00040 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00041 #endif
00042
00043 #include <list>
00044
00045 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00046
00047 TrajectoryManager::TrajectoryManager(FSimEvent* aSimEvent,
00048 const edm::ParameterSet& matEff,
00049 const edm::ParameterSet& simHits,
00050 const edm::ParameterSet& decays,
00051 const RandomEngine* engine) :
00052 mySimEvent(aSimEvent),
00053 _theGeometry(0),
00054 _theFieldMap(0),
00055 theMaterialEffects(0),
00056 myDecayEngine(0),
00057 theGeomTracker(0),
00058 theGeomSearchTracker(0),
00059 theLayerMap(56, static_cast<const DetLayer*>(0)),
00060 theNegLayerOffset(27),
00061
00062 random(engine),
00063 use_hardcoded(1)
00064
00065 {
00066
00067 use_hardcoded = matEff.getParameter<bool>("use_hardcoded_geometry");
00068
00069
00070 if ( decays.getParameter<bool>("ActivateDecays") && ( decays.getParameter<std::string>("Decayer") == "pythia6" || decays.getParameter<std::string>("Decayer") == "pythia8" ) ) {
00071 decayer = decays.getParameter<std::string>("Decayer");
00072 myDecayEngine = new PythiaDecays(decayer);
00073 distCut = decays.getParameter<double>("DistCut");
00074 } else if (! ( decays.getParameter<std::string>("Decayer") == "pythia6" || decays.getParameter<std::string>("Decayer") == "pythia8" ) )
00075 std::cout << "No valid decayer has been selected! No decay performed..." << std::endl;
00076
00077 if ( matEff.getParameter<bool>("PairProduction") ||
00078 matEff.getParameter<bool>("Bremsstrahlung") ||
00079 matEff.getParameter<bool>("MuonBremsstrahlung") ||
00080 matEff.getParameter<bool>("EnergyLoss") ||
00081 matEff.getParameter<bool>("MultipleScattering") ||
00082 matEff.getParameter<bool>("NuclearInteraction")
00083 )
00084 theMaterialEffects = new MaterialEffects(matEff,random);
00085
00086
00087
00088 firstLoop = simHits.getUntrackedParameter<bool>("firstLoop",true);
00089
00090 pTmin = simHits.getUntrackedParameter<double>("pTmin",0.5);
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103 }
00104
00105 void
00106 TrajectoryManager::initializeRecoGeometry(const GeometricSearchTracker* geomSearchTracker,
00107 const TrackerInteractionGeometry* interactionGeometry,
00108 const MagneticFieldMap* aFieldMap)
00109 {
00110
00111
00112 theGeomSearchTracker = geomSearchTracker;
00113
00114
00115 _theGeometry = interactionGeometry;
00116
00117 initializeLayerMap();
00118
00119
00120 _theFieldMap = aFieldMap;
00121
00122 }
00123
00124 void
00125 TrajectoryManager::initializeTrackerGeometry(const TrackerGeometry* geomTracker) {
00126
00127 theGeomTracker = geomTracker;
00128
00129 }
00130
00131 const TrackerInteractionGeometry*
00132 TrajectoryManager::theGeometry() {
00133 return _theGeometry;
00134 }
00135
00136 TrajectoryManager::~TrajectoryManager() {
00137
00138 if ( myDecayEngine ) delete myDecayEngine;
00139 if ( theMaterialEffects ) delete theMaterialEffects;
00140
00141
00142
00143
00144
00145
00146 }
00147
00148 void
00149 TrajectoryManager::reconstruct(const TrackerTopology *tTopo)
00150 {
00151
00152
00153
00154 thePSimHits.clear();
00155
00156
00157 XYZTLorentzVector myBeamPipe = XYZTLorentzVector(0.,2.5, 9999999.,0.);
00158
00159 std::list<TrackerLayer>::const_iterator cyliter;
00160
00161
00162
00163
00164 for( int fsimi=0; fsimi < (int) mySimEvent->nTracks(); ++fsimi) {
00165
00166
00167
00168
00169 if( !mySimEvent->track(fsimi).notYetToEndVertex(myBeamPipe) ) continue;
00170 mySimEvent->track(fsimi).setPropagate();
00171
00172
00173 cyliter = _theGeometry->cylinderBegin();
00174
00175 ParticlePropagator PP(mySimEvent->track(fsimi),_theFieldMap,random);
00176
00177 int success = 1;
00178 int sign = +1;
00179 int loop = 0;
00180 int cyl = 0;
00181
00182
00183 for ( ; cyliter != _theGeometry->cylinderEnd() ; ++cyliter ) {
00184
00185 PP.setPropagationConditions(*cyliter);
00186 if ( PP.inside() && !PP.onSurface() ) break;
00187 ++cyl;
00188
00189 }
00190
00191
00192
00193
00194
00195
00196
00197
00198 double ppcos2T = PP.cos2Theta();
00199 double ppcos2V = PP.cos2ThetaV();
00200
00201 if(use_hardcoded){
00202 if ( ( ppcos2T > 0.99 && ppcos2T < 0.9998 ) && ( cyl == 0 || ( ppcos2V > 0.99 && ppcos2V < 0.9998 ) ) ){
00203 if ( cyliter != _theGeometry->cylinderEnd() ) {
00204 cyliter = _theGeometry->cylinderEnd();
00205 --cyliter;
00206 }
00207
00208 } else if ( ppcos2T > 0.9998 && ( cyl == 0 || ppcos2V > 0.9998 ) ) {
00209 cyliter = _theGeometry->cylinderEnd();
00210 }
00211 }
00212 else {
00213 if ( ppcos2T > 0.9998 && ( cyl == 0 || ppcos2V > 0.9998 ) ) {
00214 cyliter = _theGeometry->cylinderEnd();
00215 }
00216 }
00217
00218
00219 while ( cyliter != _theGeometry->cylinderEnd() &&
00220 loop<100 &&
00221 mySimEvent->track(fsimi).notYetToEndVertex(PP.vertex())) {
00222
00223
00224 if ( cyliter->surface().mediumProperties().radLen() < 1E-10 ) {
00225 ++cyliter; ++cyl;
00226 continue;
00227 }
00228
00229
00230
00231
00232 bool escapeBarrel = PP.getSuccess() == -1;
00233 bool escapeEndcap = (PP.getSuccess() == -2 && success == 1);
00234
00235 bool fullPropagation =
00236 (PP.getSuccess() <= 0 && success==0) || escapeEndcap;
00237
00238 if ( escapeBarrel ) {
00239 ++cyliter; ++cyl;
00240 while (cyliter != _theGeometry->cylinderEnd() && cyliter->forward() ) {
00241 sign=1; ++cyliter; ++cyl;
00242 }
00243
00244 if ( cyliter == _theGeometry->cylinderEnd() ) {
00245 --cyliter; --cyl; fullPropagation=true;
00246 }
00247
00248 }
00249
00250
00251 PP.setPropagationConditions(*cyliter,!fullPropagation);
00252 if ( escapeEndcap ) PP.increaseRCyl(0.0005);
00253
00254
00255 success = PP.getSuccess();
00256
00257
00258
00259 if ( !PP.propagateToBoundSurface(*cyliter) ||
00260 PP.getSuccess()<=0) {
00261 sign = -sign;
00262 ++loop;
00263 }
00264
00265
00266
00267 if ( PP.hasDecayed() || (!mySimEvent->track(fsimi).nDaughters() && PP.PDGcTau()<1E-3 ) ) {
00268 updateWithDaughters(PP,fsimi);
00269 break;
00270 }
00271
00272
00273
00274 if ( PP.getSuccess()==2 || cyliter==_theGeometry->cylinderBegin() )
00275 sign = +1;
00276
00277
00278 if( PP.getSuccess() > 0 && PP.onFiducial() ) {
00279
00280 bool saveHit =
00281 ( (loop==0 && sign>0) || !firstLoop ) &&
00282 PP.charge()!=0. &&
00283 cyliter->sensitive() &&
00284 PP.Perp2()>pTmin*pTmin;
00285
00286
00287 if ( theMaterialEffects )
00288 theMaterialEffects->interact(*mySimEvent,*cyliter,PP,fsimi);
00289
00290
00291 saveHit &= PP.E()>1E-6;
00292
00293 if ( saveHit ) {
00294
00295 if ( cyliter->sensitive() ) {
00296
00297
00298
00299
00300
00301 if ( theGeomTracker )
00302 createPSimHits(*cyliter, PP, thePSimHits[fsimi], fsimi,mySimEvent->track(fsimi).type(), tTopo);
00303
00304
00305
00306
00307
00308
00309
00310
00311
00312 }
00313 }
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324
00325 if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.vertex()) &&
00326 !mySimEvent->filter().accept(PP) )
00327 mySimEvent->addSimVertex(PP.vertex(),fsimi, FSimVertexType::END_VERTEX);
00328
00329 }
00330
00331
00332 if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.vertex()) ) {
00333
00334
00335
00336 if (sign==1) {++cyliter;++cyl;}
00337 else {--cyliter;--cyl;}
00338
00339
00340 if( cyliter==_theGeometry->cylinderEnd()) {
00341
00342
00343
00344
00345 PP.propagateToEcal();
00346
00347
00348
00349 if(PP.getSuccess()==0) {
00350 --cyliter; --cyl; sign = -sign;
00351 PP.setPropagationConditions(*cyliter);
00352 PP.propagateToBoundSurface(*cyliter);
00353
00354
00355 if(PP.getSuccess()<0) {++cyliter; ++cyl;}
00356
00357 }
00358
00359
00360 if ( PP.hasDecayed() )
00361 updateWithDaughters(PP,fsimi);
00362
00363 }
00364 }
00365
00366 }
00367
00368
00369
00370 if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.vertex()) )
00371 propagateToCalorimeters(PP,fsimi);
00372
00373 }
00374
00375
00376 if ( theMaterialEffects ) theMaterialEffects->save();
00377
00378 }
00379
00380 void
00381 TrajectoryManager::propagateToCalorimeters(ParticlePropagator& PP, int fsimi) {
00382
00383 FSimTrack& myTrack = mySimEvent->track(fsimi);
00384
00385
00386 myTrack.setTkPosition(PP.vertex().Vect());
00387 myTrack.setTkMomentum(PP.momentum());
00388
00389
00390 PP.propagateToPreshowerLayer1(false);
00391 if ( PP.hasDecayed() ) {
00392 updateWithDaughters(PP,fsimi);
00393 return;
00394 }
00395 if ( myTrack.notYetToEndVertex(PP.vertex()) && PP.getSuccess() > 0 )
00396 myTrack.setLayer1(PP,PP.getSuccess());
00397
00398
00399 PP.propagateToPreshowerLayer2(false);
00400 if ( PP.hasDecayed() ) {
00401 updateWithDaughters(PP,fsimi);
00402 return;
00403 }
00404 if ( myTrack.notYetToEndVertex(PP.vertex()) && PP.getSuccess() > 0 )
00405 myTrack.setLayer2(PP,PP.getSuccess());
00406
00407
00408 PP.propagateToEcalEntrance(false);
00409 if ( PP.hasDecayed() ) {
00410 updateWithDaughters(PP,fsimi);
00411 return;
00412 }
00413 if ( myTrack.notYetToEndVertex(PP.vertex()) )
00414 myTrack.setEcal(PP,PP.getSuccess());
00415
00416
00417 PP.propagateToHcalEntrance(false);
00418 if ( PP.hasDecayed() ) {
00419 updateWithDaughters(PP,fsimi);
00420 return;
00421 }
00422 if ( myTrack.notYetToEndVertex(PP.vertex()) )
00423 myTrack.setHcal(PP,PP.getSuccess());
00424
00425
00426 PP.propagateToVFcalEntrance(false);
00427 if ( PP.hasDecayed() ) {
00428 updateWithDaughters(PP,fsimi);
00429 return;
00430 }
00431 if ( myTrack.notYetToEndVertex(PP.vertex()) )
00432 myTrack.setVFcal(PP,PP.getSuccess());
00433
00434 }
00435
00436 bool
00437 TrajectoryManager::propagateToLayer(ParticlePropagator& PP, unsigned layer) {
00438
00439 std::list<TrackerLayer>::const_iterator cyliter;
00440 bool done = false;
00441
00442
00443 cyliter = _theGeometry->cylinderBegin();
00444
00445
00446 for ( ; cyliter != _theGeometry->cylinderEnd() ; ++cyliter ) {
00447
00448 if ( layer != cyliter->layerNumber() ) continue;
00449
00450 PP.setPropagationConditions(*cyliter);
00451
00452 done =
00453 PP.propagateToBoundSurface(*cyliter) &&
00454 PP.getSuccess() > 0 &&
00455 PP.onFiducial();
00456
00457 break;
00458
00459 }
00460
00461 return done;
00462
00463 }
00464
00465 void
00466 TrajectoryManager::updateWithDaughters(ParticlePropagator& PP, int fsimi) {
00467
00468
00469
00470
00471 unsigned nDaugh = mySimEvent->track(fsimi).nDaughters();
00472 if ( nDaugh ) {
00473
00474
00475 unsigned vertexId = mySimEvent->track(fsimi).endVertex().id();
00476 mySimEvent->vertex(vertexId).setPosition(PP.vertex());
00477
00478
00479 XYZTLorentzVector momentumBefore = mySimEvent->track(fsimi).momentum();
00480 XYZTLorentzVector momentumAfter = PP.momentum();
00481 double magBefore = std::sqrt(momentumBefore.Vect().mag2());
00482 double magAfter = std::sqrt(momentumAfter.Vect().mag2());
00483
00484 XYZVector axis = momentumBefore.Vect().Cross(momentumAfter.Vect());
00485 double angle = std::acos(momentumBefore.Vect().Dot(momentumAfter.Vect())/(magAfter*magBefore));
00486 Rotation r(axis,angle);
00487
00488 double rescale = magAfter/magBefore;
00489
00490
00491 moveAllDaughters(fsimi,r,rescale);
00492
00493
00494 } else {
00495
00496
00497 if ( !myDecayEngine ) return;
00498
00499
00500 const DaughterParticleList& daughters = (decayer == "pythia6") ? myDecayEngine->particleDaughtersPy6(PP) : myDecayEngine->particleDaughtersPy8(PP);
00501
00502
00503 if ( daughters.size() ) {
00504 double distMin = 1E99;
00505 int theClosestChargedDaughterId = -1;
00506 DaughterParticleIterator daughter = daughters.begin();
00507
00508 int ivertex = mySimEvent->addSimVertex(daughter->vertex(),fsimi,
00509 FSimVertexType::DECAY_VERTEX);
00510
00511 if ( ivertex != -1 ) {
00512 for ( ; daughter != daughters.end(); ++daughter) {
00513 int theDaughterId = mySimEvent->addSimTrack(&(*daughter), ivertex);
00514
00515 if ( PP.charge() * daughter->charge() > 1E-10 ) {
00516 double dist = (daughter->Vect().Unit().Cross(PP.Vect().Unit())).R();
00517 if ( dist < distCut && dist < distMin ) {
00518 distMin = dist;
00519 theClosestChargedDaughterId = theDaughterId;
00520 }
00521 }
00522 }
00523 }
00524
00525 if ( theClosestChargedDaughterId >=0 )
00526 mySimEvent->track(fsimi).setClosestDaughterId(theClosestChargedDaughterId);
00527 }
00528
00529 }
00530
00531 }
00532
00533
00534 void
00535 TrajectoryManager::moveAllDaughters(int fsimi, const Rotation& r, double rescale) {
00536
00537
00538 for ( unsigned idaugh=0; idaugh < (unsigned)(mySimEvent->track(fsimi).nDaughters()); ++idaugh) {
00539
00540 XYZTLorentzVector daughMomentum (mySimEvent->track(fsimi).daughter(idaugh).momentum());
00541
00542 XYZVector newMomentum (r * daughMomentum.Vect());
00543 newMomentum *= rescale;
00544 double newEnergy = std::sqrt(newMomentum.mag2() + daughMomentum.mag2());
00545
00546 mySimEvent->track(fsimi).setMomentum(XYZTLorentzVector(newMomentum.X(),newMomentum.Y(),newMomentum.Z(),newEnergy));
00547
00548 int fsimDaug = mySimEvent->track(fsimi).daughter(idaugh).id();
00549 moveAllDaughters(fsimDaug,r,rescale);
00550 }
00551 }
00552
00553 void
00554 TrajectoryManager::createPSimHits(const TrackerLayer& layer,
00555 const ParticlePropagator& PP,
00556 std::map<double,PSimHit>& theHitMap,
00557 int trackID, int partID, const TrackerTopology *tTopo) {
00558
00559
00560
00561
00562
00563
00564 LocalMagneticField mf(PP.getMagneticField());
00565 AnalyticalPropagator alongProp(&mf, anyDirection);
00566 InsideBoundsMeasurementEstimator est;
00567
00568
00569
00570
00571
00572 typedef GeometricSearchDet::DetWithState DetWithState;
00573 const DetLayer* tkLayer = detLayer(layer,PP.Z());
00574
00575 TrajectoryStateOnSurface trajState = makeTrajectoryState( tkLayer, PP, &mf);
00576 float thickness = theMaterialEffects ? theMaterialEffects->thickness() : 0.;
00577 float eloss = theMaterialEffects ? theMaterialEffects->energyLoss() : 0.;
00578
00579
00580
00581 std::vector<DetWithState> compat
00582 = tkLayer->compatibleDets( trajState, alongProp, est);
00583
00584
00585 std::map<double,PSimHit> theTrackHits;
00586 for (std::vector<DetWithState>::const_iterator i=compat.begin(); i!=compat.end(); i++) {
00587
00588
00589 makePSimHits( i->first, i->second, theHitMap, trackID, eloss, thickness, partID,tTopo);
00590 }
00591
00592 }
00593
00594 TrajectoryStateOnSurface
00595 TrajectoryManager::makeTrajectoryState( const DetLayer* layer,
00596 const ParticlePropagator& pp,
00597 const MagneticField* field) const
00598 {
00599 GlobalPoint pos( pp.X(), pp.Y(), pp.Z());
00600 GlobalVector mom( pp.Px(), pp.Py(), pp.Pz());
00601 ReferenceCountingPointer<TangentPlane> plane = layer->surface().tangentPlane(pos);
00602 return TrajectoryStateOnSurface
00603 (GlobalTrajectoryParameters( pos, mom, TrackCharge( pp.charge()), field), *plane);
00604 }
00605
00606 void
00607 TrajectoryManager::makePSimHits( const GeomDet* det,
00608 const TrajectoryStateOnSurface& ts,
00609 std::map<double,PSimHit>& theHitMap,
00610 int tkID, float el, float thick, int pID,
00611 const TrackerTopology *tTopo)
00612 {
00613
00614 std::vector< const GeomDet*> comp = det->components();
00615 if (!comp.empty()) {
00616 for (std::vector< const GeomDet*>::const_iterator i = comp.begin();
00617 i != comp.end(); i++) {
00618 const GeomDetUnit* du = dynamic_cast<const GeomDetUnit*>(*i);
00619 if (du != 0)
00620 theHitMap.insert(theHitMap.end(),makeSinglePSimHit( *du, ts, tkID, el, thick, pID,tTopo));
00621 }
00622 }
00623 else {
00624 const GeomDetUnit* du = dynamic_cast<const GeomDetUnit*>(det);
00625 if (du != 0)
00626 theHitMap.insert(theHitMap.end(),makeSinglePSimHit( *du, ts, tkID, el, thick, pID,tTopo));
00627 }
00628
00629 }
00630
00631 std::pair<double,PSimHit>
00632 TrajectoryManager::makeSinglePSimHit( const GeomDetUnit& det,
00633 const TrajectoryStateOnSurface& ts,
00634 int tkID, float el, float thick, int pID,
00635 const TrackerTopology *tTopo) const
00636 {
00637
00638 const float onSurfaceTolarance = 0.01;
00639
00640 LocalPoint lpos;
00641 LocalVector lmom;
00642 if ( fabs( det.toLocal(ts.globalPosition()).z()) < onSurfaceTolarance) {
00643 lpos = ts.localPosition();
00644 lmom = ts.localMomentum();
00645 }
00646 else {
00647 HelixArbitraryPlaneCrossing crossing( ts.globalPosition().basicVector(),
00648 ts.globalMomentum().basicVector(),
00649 ts.transverseCurvature(),
00650 anyDirection);
00651 std::pair<bool,double> path = crossing.pathLength(det.surface());
00652 if (!path.first) {
00653
00654 return std::pair<double,PSimHit>(0.,PSimHit());
00655 }
00656 lpos = det.toLocal( GlobalPoint( crossing.position(path.second)));
00657 lmom = det.toLocal( GlobalVector( crossing.direction(path.second)));
00658 lmom = lmom.unit() * ts.localMomentum().mag();
00659 }
00660
00661
00662 const BoundPlane& theDetPlane = det.surface();
00663 float halfThick = 0.5*theDetPlane.bounds().thickness();
00664
00665 float eloss = el;
00666 if ( thick > 0. ) {
00667
00668
00669
00670 eloss *= (2.* halfThick - 0.003) / (9.36 * thick);
00671 }
00672
00673 float pZ = lmom.z();
00674 LocalPoint entry = lpos + (-halfThick/pZ) * lmom;
00675 LocalPoint exit = lpos + halfThick/pZ * lmom;
00676 float tof = ts.globalPosition().mag() / 30. ;
00677
00678
00679
00680
00681 int localTkID = tkID;
00682 if ( mySimEvent->track(tkID).mother().closestDaughterId() == tkID )
00683 localTkID = mySimEvent->track(tkID).mother().id();
00684
00685
00686 PSimHit hit( entry, exit, lmom.mag(), tof, eloss, pID,
00687 det.geographicalId().rawId(), localTkID,
00688 lmom.theta(),
00689 lmom.phi());
00690
00691
00692 unsigned subdet = DetId(hit.detUnitId()).subdetId();
00693 double boundX = theDetPlane.bounds().width()/2.;
00694 double boundY = theDetPlane.bounds().length()/2.;
00695
00696
00697 if ( subdet == 4 || subdet == 6 )
00698 boundX *= 1. - hit.localPosition().y()/theDetPlane.position().perp();
00699
00700 #ifdef FAMOS_DEBUG
00701 unsigned detid = DetId(hit.detUnitId()).rawId();
00702 unsigned stereo = 0;
00703 unsigned theLayer = 0;
00704 unsigned theRing = 0;
00705 switch (subdet) {
00706 case 1:
00707 {
00708
00709 theLayer = tTopo->pxbLayer(detid);
00710 std::cout << "\tPixel Barrel Layer " << theLayer << std::endl;
00711 stereo = 1;
00712 break;
00713 }
00714 case 2:
00715 {
00716
00717 theLayer = tTopo->pxfDisk(detid);
00718 std::cout << "\tPixel Forward Disk " << theLayer << std::endl;
00719 stereo = 1;
00720 break;
00721 }
00722 case 3:
00723 {
00724
00725 theLayer = tTopo->tibLayer(detid);
00726 std::cout << "\tTIB Layer " << theLayer << std::endl;
00727 stereo = module.stereo();
00728 break;
00729 }
00730 case 4:
00731 {
00732
00733 theLayer = tTopo->tidWheel(detid);
00734 theRing = tTopo->tidRing(detid);
00735 unsigned int theSide = module.side();
00736 if ( theSide == 1 )
00737 std::cout << "\tTID Petal Back " << std::endl;
00738 else
00739 std::cout << "\tTID Petal Front" << std::endl;
00740 std::cout << "\tTID Layer " << theLayer << std::endl;
00741 std::cout << "\tTID Ring " << theRing << std::endl;
00742 stereo = module.stereo();
00743 break;
00744 }
00745 case 5:
00746 {
00747
00748 theLayer = tTopo->tobLayer(detid);
00749 stereo = tTopo->tobStereo(detid);
00750 std::cout << "\tTOB Layer " << theLayer << std::endl;
00751 break;
00752 }
00753 case 6:
00754 {
00755
00756 theLayer = tTopo->tecWheel(detid);
00757 theRing = tTopo->tecRing(detid);
00758 unsigned int theSide = module.petal()[0];
00759 if ( theSide == 1 )
00760 std::cout << "\tTEC Petal Back " << std::endl;
00761 else
00762 std::cout << "\tTEC Petal Front" << std::endl;
00763 std::cout << "\tTEC Layer " << theLayer << std::endl;
00764 std::cout << "\tTEC Ring " << theRing << std::endl;
00765 stereo = module.stereo();
00766 break;
00767 }
00768 default:
00769 {
00770 stereo = 0;
00771 break;
00772 }
00773 }
00774
00775 std::cout << "Thickness = " << 2.*halfThick-0.003 << "; " << thick * 9.36 << std::endl
00776 << "Length = " << det.surface().bounds().length() << std::endl
00777 << "Width = " << det.surface().bounds().width() << std::endl;
00778
00779 std::cout << "Hit position = "
00780 << hit.localPosition().x() << " "
00781 << hit.localPosition().y() << " "
00782 << hit.localPosition().z() << std::endl;
00783 #endif
00784
00785
00786
00787
00788
00789
00790 double dist = 0.;
00791 GlobalPoint IP (mySimEvent->track(localTkID).vertex().position().x(),
00792 mySimEvent->track(localTkID).vertex().position().y(),
00793 mySimEvent->track(localTkID).vertex().position().z());
00794
00795 dist = ( fabs(hit.localPosition().x()) > boundX ||
00796 fabs(hit.localPosition().y()) > boundY ) ?
00797
00798 -( det.surface().toGlobal(hit.localPosition()) - IP ).mag2()
00799 :
00800
00801 ( det.surface().toGlobal(hit.localPosition()) - IP ).mag2();
00802
00803
00804
00805
00806
00807
00808
00809
00810
00811
00812
00813
00814
00815 return std::pair<double,PSimHit>(dist,hit);
00816
00817 }
00818
00819 void
00820 TrajectoryManager::initializeLayerMap()
00821 {
00822
00823
00824
00825
00826
00827
00828
00829
00830
00833
00834 std::vector< BarrelDetLayer*> barrelLayers =
00835 theGeomSearchTracker->barrelLayers();
00836 LogDebug("FastTracking") << "Barrel DetLayer dump: ";
00837 for (std::vector< BarrelDetLayer*>::const_iterator bl=barrelLayers.begin();
00838 bl != barrelLayers.end(); ++bl) {
00839 LogDebug("FastTracking")<< "radius " << (**bl).specificSurface().radius();
00840 }
00841
00842 std::vector< ForwardDetLayer*> posForwardLayers =
00843 theGeomSearchTracker->posForwardLayers();
00844 LogDebug("FastTracking") << "Positive Forward DetLayer dump: ";
00845 for (std::vector< ForwardDetLayer*>::const_iterator fl=posForwardLayers.begin();
00846 fl != posForwardLayers.end(); ++fl) {
00847 LogDebug("FastTracking") << "Z pos "
00848 << (**fl).surface().position().z()
00849 << " radii "
00850 << (**fl).specificSurface().innerRadius()
00851 << ", "
00852 << (**fl).specificSurface().outerRadius();
00853 }
00854
00855 const float rTolerance = 1.5;
00856 const float zTolerance = 3.;
00857
00858 LogDebug("FastTracking")<< "Dump of TrackerInteractionGeometry cylinders:";
00859 for( std::list<TrackerLayer>::const_iterator i=_theGeometry->cylinderBegin();
00860 i!=_theGeometry->cylinderEnd(); ++i) {
00861 const BoundCylinder* cyl = i->cylinder();
00862 const BoundDisk* disk = i->disk();
00863
00864 LogDebug("FastTracking") << "Famos Layer no " << i->layerNumber()
00865 << " is sensitive? " << i->sensitive()
00866 << " pos " << i->surface().position();
00867 if (!i->sensitive()) continue;
00868
00869 if (cyl != 0) {
00870 LogDebug("FastTracking") << " cylinder radius " << cyl->radius();
00871 bool found = false;
00872 for (std::vector< BarrelDetLayer*>::const_iterator
00873 bl=barrelLayers.begin(); bl != barrelLayers.end(); ++bl) {
00874
00875 if (fabs( cyl->radius() - (**bl).specificSurface().radius()) < rTolerance) {
00876 theLayerMap[i->layerNumber()] = *bl;
00877 found = true;
00878 LogDebug("FastTracking")<< "Corresponding DetLayer found with radius "
00879 << (**bl).specificSurface().radius();
00880 break;
00881 }
00882 }
00883 if (!found) {
00884 edm::LogWarning("FastTracking") << " Trajectory manager FAILED to find a corresponding DetLayer!";
00885 }
00886 }
00887 else {
00888 LogDebug("FastTracking") << " disk radii " << disk->innerRadius()
00889 << ", " << disk->outerRadius();
00890 bool found = false;
00891 for (std::vector< ForwardDetLayer*>::const_iterator fl=posForwardLayers.begin();
00892 fl != posForwardLayers.end(); ++fl) {
00893
00894 if (fabs( disk->position().z() - (**fl).surface().position().z()) < zTolerance) {
00895 theLayerMap[i->layerNumber()] = *fl;
00896 found = true;
00897 LogDebug("FastTracking") << "Corresponding DetLayer found with Z pos "
00898 << (**fl).surface().position().z()
00899 << " and radii "
00900 << (**fl).specificSurface().innerRadius()
00901 << ", "
00902 << (**fl).specificSurface().outerRadius();
00903 break;
00904 }
00905 }
00906 if (!found) {
00907 edm::LogWarning("FastTracking") << "FAILED to find a corresponding DetLayer!";
00908 }
00909 }
00910 }
00911
00912
00913 std::vector< ForwardDetLayer*> negForwardLayers = theGeomSearchTracker->negForwardLayers();
00914 for (std::vector< ForwardDetLayer*>::const_iterator nl=negForwardLayers.begin();
00915 nl != negForwardLayers.end(); ++nl) {
00916 for (int i=0; i<=theNegLayerOffset; i++) {
00917 if (theLayerMap[i] == 0) continue;
00918 if ( fabs( (**nl).surface().position().z() +theLayerMap[i]-> surface().position().z()) < zTolerance) {
00919 theLayerMap[i+theNegLayerOffset] = *nl;
00920 break;
00921 }
00922 }
00923 }
00924
00925 }
00926
00927 const DetLayer*
00928 TrajectoryManager::detLayer( const TrackerLayer& layer, float zpos) const
00929 {
00930 if (zpos > 0 || !layer.forward() ) return theLayerMap[layer.layerNumber()];
00931 else return theLayerMap[layer.layerNumber()+theNegLayerOffset];
00932 }
00933
00934 void
00935 TrajectoryManager::loadSimHits(edm::PSimHitContainer & c) const
00936 {
00937
00938 std::map<unsigned,std::map<double,PSimHit> >::const_iterator itrack = thePSimHits.begin();
00939 std::map<unsigned,std::map<double,PSimHit> >::const_iterator itrackEnd = thePSimHits.end();
00940 for ( ; itrack != itrackEnd; ++itrack ) {
00941 std::map<double,PSimHit>::const_iterator it = (itrack->second).begin();
00942 std::map<double,PSimHit>::const_iterator itEnd = (itrack->second).end();
00943 for( ; it!= itEnd; ++it) {
00944
00945
00946
00947
00948
00949
00950
00951
00952
00953
00954 if ( it->first > 0. ) c.push_back(it->second);
00955 }
00956 }
00957
00958 }