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