CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/FastSimulation/TrajectoryManager/src/TrajectoryManager.cc

Go to the documentation of this file.
00001 //Framework Headers
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 
00004 //CMSSW Headers
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 // Tracker reco geometry headers 
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 //#include "Geometry/TrackerGeometryBuilder/interface/TrackerGeometry.h"
00020 
00021 //FAMOS Headers
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 //#include "FastSimulation/Utilities/interface/Histos.h"
00034 //#include "FastSimulation/Utilities/interface/FamosLooses.h"
00035 // Numbering scheme
00036 
00037 //#define FAMOS_DEBUG
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)), // reserve space for layers here
00064   theNegLayerOffset(27),
00065   //  myHistos(0),
00066   random(engine)
00067 
00068 {
00069   
00070   // Initialize Bthe stable particle decay engine 
00071   if ( decays.getParameter<bool>("ActivateDecays") ) { 
00072 //    int seed = (int) ( 900000000. * random->flatShoot() );
00073 //    double comE = decays.getParameter<double>("comEnergy");
00074 //    myDecayEngine = new Pythia6Decays(seed,comE);
00075     myDecayEngine = new Pythia6Decays();
00076     distCut = decays.getParameter<double>("DistCut");
00077   }
00078    // new improvement: Muon brem effects 27-Fev-2011-S.Fonseca UERJ/Brazil
00079   // Initialize the Material Effects updator, if needed
00080   if ( matEff.getParameter<bool>("PairProduction") || 
00081        matEff.getParameter<bool>("Bremsstrahlung") ||
00082        matEff.getParameter<bool>("MuonBremsstrahlung") ||
00083        matEff.getParameter<bool>("EnergyLoss") || 
00084        matEff.getParameter<bool>("MultipleScattering") || 
00085        matEff.getParameter<bool>("NuclearInteraction")
00086        )
00087        theMaterialEffects = new MaterialEffects(matEff,random);
00088 
00089   // Save SimHits according to Optiom
00090   // Only the hits from first half loop is saved
00091   firstLoop = simHits.getUntrackedParameter<bool>("firstLoop",true);
00092   // Only if pT>pTmin are the hits saved
00093   pTmin = simHits.getUntrackedParameter<double>("pTmin",0.5);
00094 
00095   // Get the Famos Histos pointer
00096   //  myHistos = Histos::instance();
00097 
00098   // Initialize a few histograms
00099   /* 
00100   myHistos->book("h300",1210,-121.,121.,1210,-121.,121.);
00101   myHistos->book("h301",1200,-300.,300.,1210,-121.,121.);
00102   */
00103 
00104   
00105 }
00106 
00107 void 
00108 TrajectoryManager::initializeRecoGeometry(const GeometricSearchTracker* geomSearchTracker,
00109                                           const TrackerInteractionGeometry* interactionGeometry,
00110                                           const MagneticFieldMap* aFieldMap)
00111 {
00112   
00113   // Initialize the reco tracker geometry
00114   theGeomSearchTracker = geomSearchTracker;
00115   
00116   // Initialize the simplified tracker geometry
00117   _theGeometry = interactionGeometry;
00118 
00119   initializeLayerMap();
00120 
00121   // Initialize the magnetic field
00122   _theFieldMap = aFieldMap;
00123 
00124 }
00125 
00126 void 
00127 TrajectoryManager::initializeTrackerGeometry(const TrackerGeometry* geomTracker) { 
00128   
00129   theGeomTracker = geomTracker;
00130 
00131 }
00132 
00133 const TrackerInteractionGeometry*
00134 TrajectoryManager::theGeometry() {
00135   return _theGeometry;
00136 }
00137 
00138 TrajectoryManager::~TrajectoryManager() {
00139 
00140   if ( myDecayEngine ) delete myDecayEngine;
00141   if ( theMaterialEffects ) delete theMaterialEffects;
00142 
00143   //Write the histograms
00144   //myHistos->put("histos.root");
00145   //  if ( myHistos ) delete myHistos;
00146 
00147 }
00148 
00149 void
00150 TrajectoryManager::reconstruct()
00151 {
00152 
00153   // Clear the hits of the previous event
00154   //  thePSimHits->clear();
00155   thePSimHits.clear();
00156 
00157   // The new event
00158   XYZTLorentzVector myBeamPipe = XYZTLorentzVector(0.,2.5, 9999999.,0.);
00159 
00160   std::list<TrackerLayer>::const_iterator cyliter;
00161 
00162   // bool debug = mySimEvent->id().event() == 8;
00163 
00164   // Loop over the particles (watch out: increasing upper limit!)
00165   for( int fsimi=0; fsimi < (int) mySimEvent->nTracks(); ++fsimi) {
00166 
00167     // If the particle has decayed inside the beampipe, or decays 
00168     // immediately, there is nothing to do
00169     //if ( debug ) std::cout << mySimEvent->track(fsimi) << std::endl;
00170     //if ( debug ) std::cout << "Not yet at end vertex ? " << mySimEvent->track(fsimi).notYetToEndVertex(myBeamPipe) << std::endl;
00171     if( !mySimEvent->track(fsimi).notYetToEndVertex(myBeamPipe) ) continue;
00172     mySimEvent->track(fsimi).setPropagate();
00173 
00174     // Get the geometry elements 
00175     cyliter = _theGeometry->cylinderBegin();
00176     // Prepare the propagation  
00177     ParticlePropagator PP(mySimEvent->track(fsimi),_theFieldMap,random);
00178     //The real work starts here
00179     int success = 1;
00180     int sign = +1;
00181     int loop = 0;
00182     int cyl = 0;
00183 
00184     // Find the initial cylinder to propagate to.      
00185     for ( ; cyliter != _theGeometry->cylinderEnd() ; ++cyliter ) {
00186       
00187       PP.setPropagationConditions(*cyliter);
00188       if ( PP.inside() && !PP.onSurface() ) break;
00189       ++cyl;
00190 
00191     }
00192 
00193     // The particle has a pseudo-rapidity (position or momentum direction) 
00194     // in excess of 3.0. Just simply go to the last tracker layer
00195     // without bothering with all the details of the propagation and 
00196     // material effects.
00197     // 08/02/06 - pv: increase protection from 0.99 (eta=2.9932) to 0.9998 (eta=4.9517)
00198     //                to simulate material effects at large eta 
00199     // if above 0.99: propagate to the last tracker cylinder where the material is concentrated!
00200     double ppcos2T =  PP.cos2Theta();
00201     double ppcos2V =  PP.cos2ThetaV();
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     // if above 0.9998: don't propagate at all (only to the calorimeters directly)
00208     } else if ( ppcos2T > 0.9998 && ( cyl == 0 || ppcos2V > 0.9998 ) ) { 
00209       cyliter = _theGeometry->cylinderEnd();
00210     }
00211         
00212     // Loop over the cylinders
00213     while ( cyliter != _theGeometry->cylinderEnd() &&
00214             loop<100 &&                            // No more than 100 loops
00215             mySimEvent->track(fsimi).notYetToEndVertex(PP.vertex())) { // The particle decayed
00216 
00217       // Skip layers with no material (kept just for historical reasons)
00218       if ( cyliter->surface().mediumProperties()->radLen() < 1E-10 ) { 
00219         ++cyliter; ++cyl;
00220         continue;
00221       }
00222       
00223       // Pathological cases:
00224       // To prevent from interacting twice in a row with the same layer
00225       //      bool escapeBarrel    = (PP.getSuccess() == -1 && success == 1);
00226       bool escapeBarrel    = PP.getSuccess() == -1;
00227       bool escapeEndcap    = (PP.getSuccess() == -2 && success == 1);
00228       // To break the loop
00229       bool fullPropagation = 
00230         (PP.getSuccess() <= 0 && success==0) || escapeEndcap;
00231 
00232       if ( escapeBarrel ) {
00233         ++cyliter; ++cyl;
00234         while (cyliter != _theGeometry->cylinderEnd() && cyliter->forward() ) {
00235           sign=1; ++cyliter; ++cyl;
00236         }
00237 
00238         if ( cyliter == _theGeometry->cylinderEnd()  ) {
00239           --cyliter; --cyl; fullPropagation=true; 
00240         }
00241 
00242       }
00243 
00244       // Define the propagation conditions
00245       PP.setPropagationConditions(*cyliter,!fullPropagation);
00246       if ( escapeEndcap ) PP.increaseRCyl(0.0005);
00247 
00248       // Remember last propagation outcome
00249       success = PP.getSuccess();
00250 
00251       // Propagation was not successful :
00252       // Change the sign of the cylinder increment and count the loops
00253       if ( !PP.propagateToBoundSurface(*cyliter) || 
00254            PP.getSuccess()<=0) {
00255         sign = -sign;
00256         ++loop;
00257       }
00258 
00259       // The particle may have decayed on its way... in which the daughters
00260       // have to be added to the event record
00261       if ( PP.hasDecayed() || (!mySimEvent->track(fsimi).nDaughters() && PP.PDGcTau()<1E-3 ) ) { 
00262         updateWithDaughters(PP,fsimi);
00263         break;
00264       }
00265 
00266       // Exit by the endcaps or innermost cylinder :
00267       // Positive cylinder increment
00268       if ( PP.getSuccess()==2 || cyliter==_theGeometry->cylinderBegin() ) 
00269         sign = +1; 
00270           
00271       // Successful propagation to a cylinder, with some Material :
00272       if( PP.getSuccess() > 0 && PP.onFiducial() ) {
00273 
00274         bool saveHit = 
00275           ( (loop==0 && sign>0) || !firstLoop ) &&   // Save only first half loop
00276           PP.charge()!=0. &&                         // Consider only charged particles
00277           cyliter->sensitive() &&                    // Consider only sensitive layers
00278           PP.Perp2()>pTmin*pTmin;                    // Consider only pT > pTmin
00279         
00280         // Material effects are simulated there
00281         if ( theMaterialEffects ) 
00282           theMaterialEffects->interact(*mySimEvent,*cyliter,PP,fsimi); 
00283 
00284         // There is a PP.setXYZT=(0,0,0,0) if bremss fails
00285         saveHit &= PP.E()>1E-6;
00286 
00287         if ( saveHit ) { 
00288 
00289           // Consider only active layers
00290           if ( cyliter->sensitive() ) {
00291             // Add information to the FSimTrack (not yet available)
00292             //      myTrack.addSimHit(PP,layer);
00293 
00294             // Return one or two (for overlap regions) PSimHits in the full 
00295             // tracker geometry
00296             if ( theGeomTracker ) 
00297               createPSimHits(*cyliter, PP, thePSimHits[fsimi], fsimi,mySimEvent->track(fsimi).type());
00298 
00299           }
00300         }
00301 
00302         // Fill Histos (~poor man event display)
00303         /* 
00304         myHistos->fill("h300",PP.x(),PP.y());
00305         if ( sin(PP.vertex().phi()) > 0. ) 
00306           myHistos->fill("h301",PP.z(),PP.vertex().perp());
00307         else
00308           myHistos->fill("h301",PP.z(),-PP.vertex().perp());
00309         */
00310 
00311         //The particle may have lost its energy in the material
00312         if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.vertex()) && 
00313              !mySimEvent->filter().accept(PP)  ) 
00314           mySimEvent->addSimVertex(PP.vertex(),fsimi, FSimVertexType::END_VERTEX);
00315           
00316       }
00317 
00318       // Stop here if the particle has reached an end
00319       if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.vertex()) ) {
00320 
00321         // Otherwise increment the cylinder iterator
00322         //      do { 
00323         if (sign==1) {++cyliter;++cyl;}
00324         else         {--cyliter;--cyl;}
00325 
00326         // Check if the last surface has been reached 
00327         if( cyliter==_theGeometry->cylinderEnd()) {
00328 
00329           // Try to propagate to the ECAL in half a loop
00330           // Note: Layer1 = ECAL Barrel entrance, or Preshower
00331           // entrance, or ECAL Endcap entrance (in the corner)
00332           PP.propagateToEcal();
00333           // PP.propagateToPreshowerLayer1();
00334 
00335           // If it is not possible, try go back to the last cylinder
00336           if(PP.getSuccess()==0) {
00337             --cyliter; --cyl; sign = -sign;
00338             PP.setPropagationConditions(*cyliter);
00339             PP.propagateToBoundSurface(*cyliter);
00340 
00341             // If there is definitely no way, leave it here.
00342             if(PP.getSuccess()<0) {++cyliter; ++cyl;}
00343 
00344           }
00345 
00346           // Check if the particle has decayed on the way to ECAL
00347           if ( PP.hasDecayed() )
00348             updateWithDaughters(PP,fsimi);
00349 
00350         }
00351       }
00352 
00353     }
00354 
00355     // Propagate all particles without a end vertex to the Preshower, 
00356     // theECAL and the HCAL.
00357     if ( mySimEvent->track(fsimi).notYetToEndVertex(PP.vertex()) )
00358       propagateToCalorimeters(PP,fsimi);
00359 
00360   }
00361 
00362   // Save the information from Nuclear Interaction Generation
00363   if ( theMaterialEffects ) theMaterialEffects->save();
00364 
00365 }
00366 
00367 void 
00368 TrajectoryManager::propagateToCalorimeters(ParticlePropagator& PP, int fsimi) {
00369 
00370   FSimTrack& myTrack = mySimEvent->track(fsimi);
00371 
00372   // Set the position and momentum at the end of the tracker volume
00373   myTrack.setTkPosition(PP.vertex().Vect());
00374   myTrack.setTkMomentum(PP.momentum());
00375 
00376   // Propagate to Preshower Layer 1 
00377   PP.propagateToPreshowerLayer1(false);
00378   if ( PP.hasDecayed() ) {
00379     updateWithDaughters(PP,fsimi);
00380     return;
00381   }
00382   if ( myTrack.notYetToEndVertex(PP.vertex()) && PP.getSuccess() > 0 )
00383     myTrack.setLayer1(PP,PP.getSuccess());
00384   
00385   // Propagate to Preshower Layer 2 
00386   PP.propagateToPreshowerLayer2(false);
00387   if ( PP.hasDecayed() ) { 
00388     updateWithDaughters(PP,fsimi);
00389     return;
00390   }
00391   if ( myTrack.notYetToEndVertex(PP.vertex()) && PP.getSuccess() > 0 )
00392     myTrack.setLayer2(PP,PP.getSuccess());
00393 
00394   // Propagate to Ecal Endcap
00395   PP.propagateToEcalEntrance(false);
00396   if ( PP.hasDecayed() ) { 
00397     updateWithDaughters(PP,fsimi);
00398     return;
00399   }
00400   if ( myTrack.notYetToEndVertex(PP.vertex()) )
00401     myTrack.setEcal(PP,PP.getSuccess());
00402 
00403   // Propagate to HCAL entrance
00404   PP.propagateToHcalEntrance(false);
00405   if ( PP.hasDecayed() ) { 
00406     updateWithDaughters(PP,fsimi);
00407     return;
00408   }
00409   if ( myTrack.notYetToEndVertex(PP.vertex()) )
00410     myTrack.setHcal(PP,PP.getSuccess());
00411 
00412   // Propagate to VFCAL entrance
00413   PP.propagateToVFcalEntrance(false);
00414   if ( PP.hasDecayed() ) { 
00415     updateWithDaughters(PP,fsimi);
00416     return;
00417   }
00418   if ( myTrack.notYetToEndVertex(PP.vertex()) )
00419     myTrack.setVFcal(PP,PP.getSuccess());
00420     
00421 }
00422 
00423 bool
00424 TrajectoryManager::propagateToLayer(ParticlePropagator& PP, unsigned layer) {
00425 
00426   std::list<TrackerLayer>::const_iterator cyliter;
00427   bool done = false;
00428 
00429   // Get the geometry elements 
00430   cyliter = _theGeometry->cylinderBegin();
00431 
00432   // Find the layer to propagate to.      
00433   for ( ; cyliter != _theGeometry->cylinderEnd() ; ++cyliter ) {
00434 
00435     if ( layer != cyliter->layerNumber() ) continue;
00436       
00437     PP.setPropagationConditions(*cyliter);
00438 
00439     done =  
00440       PP.propagateToBoundSurface(*cyliter) &&
00441       PP.getSuccess() > 0 && 
00442       PP.onFiducial();
00443 
00444     break;
00445     
00446   }
00447 
00448   return done;
00449 
00450 }
00451 
00452 void
00453 TrajectoryManager::updateWithDaughters(ParticlePropagator& PP, int fsimi) {
00454 
00455 
00456   // The particle was already decayed in the GenEvent, but still the particle was 
00457   // allowed to propagate (for magnetic field bending, for material effects, etc...)
00458   // Just modify the momentum of the daughters in that case 
00459   unsigned nDaugh = mySimEvent->track(fsimi).nDaughters();
00460   if ( nDaugh ) {
00461 
00462     // Move the vertex
00463     unsigned vertexId = mySimEvent->track(fsimi).endVertex().id();
00464     mySimEvent->vertex(vertexId).setPosition(PP.vertex());
00465 
00466     // Before-propagation and after-propagation momentum and vertex position
00467     XYZTLorentzVector momentumBefore = mySimEvent->track(fsimi).momentum();
00468     XYZTLorentzVector momentumAfter = PP.momentum();
00469     double magBefore = std::sqrt(momentumBefore.Vect().mag2());
00470     double magAfter = std::sqrt(momentumAfter.Vect().mag2());
00471     // Rotation to be applied
00472     XYZVector axis = momentumBefore.Vect().Cross(momentumAfter.Vect());
00473     double angle = std::acos(momentumBefore.Vect().Dot(momentumAfter.Vect())/(magAfter*magBefore));
00474     Rotation r(axis,angle);
00475     // Rescaling to be applied
00476     double rescale = magAfter/magBefore;
00477 
00478     // Move, rescale and rotate daugthers, grand-daughters, etc. 
00479     moveAllDaughters(fsimi,r,rescale);
00480 
00481   // The particle is not decayed in the GenEvent, decay it with PYTHIA 
00482   } else { 
00483 
00484     // Decays are not activated : do nothing
00485     if ( !myDecayEngine ) return;
00486     
00487     // Invoke PYDECY to decay the particle and get the daughters
00488     const DaughterParticleList& daughters = myDecayEngine->particleDaughters(PP);
00489     
00490     // Update the FSimEvent with an end vertex and with the daughters
00491     if ( daughters.size() ) { 
00492       double distMin = 1E99;
00493       int theClosestChargedDaughterId = -1;
00494       DaughterParticleIterator daughter = daughters.begin();
00495       
00496       int ivertex = mySimEvent->addSimVertex(daughter->vertex(),fsimi, 
00497                                              FSimVertexType::DECAY_VERTEX);
00498       
00499       if ( ivertex != -1 ) {
00500         for ( ; daughter != daughters.end(); ++daughter) {
00501           int theDaughterId = mySimEvent->addSimTrack(&(*daughter), ivertex);
00502           // Find the closest charged daughter (if charged mother)
00503           if ( PP.charge() * daughter->charge() > 1E-10 ) {
00504             double dist = (daughter->Vect().Unit().Cross(PP.Vect().Unit())).R();
00505             if ( dist < distCut && dist < distMin ) { 
00506               distMin = dist;
00507               theClosestChargedDaughterId = theDaughterId;
00508             }
00509           }
00510         }
00511       }
00512       // Attach mother and closest daughter sp as to cheat tracking ;-)
00513       if ( theClosestChargedDaughterId >=0 ) 
00514         mySimEvent->track(fsimi).setClosestDaughterId(theClosestChargedDaughterId);
00515     }
00516 
00517   }
00518 
00519 }
00520 
00521 
00522 void
00523 TrajectoryManager::moveAllDaughters(int fsimi, const Rotation& r, double rescale) { 
00524 
00525   //
00526   for ( unsigned idaugh=0; idaugh < (unsigned)(mySimEvent->track(fsimi).nDaughters()); ++idaugh) { 
00527     // Initial momentum of the daughter
00528     XYZTLorentzVector daughMomentum (mySimEvent->track(fsimi).daughter(idaugh).momentum()); 
00529     // Rotate and rescale
00530     XYZVector newMomentum (r * daughMomentum.Vect()); 
00531     newMomentum *= rescale;
00532     double newEnergy = std::sqrt(newMomentum.mag2() + daughMomentum.mag2());
00533     // Set the new momentum
00534     mySimEvent->track(fsimi).setMomentum(XYZTLorentzVector(newMomentum.X(),newMomentum.Y(),newMomentum.Z(),newEnergy));
00535     // Watch out : recursive call to get all grand-daughters
00536     int fsimDaug = mySimEvent->track(fsimi).daughter(idaugh).id();
00537     moveAllDaughters(fsimDaug,r,rescale);
00538   }
00539 }
00540 
00541 void
00542 TrajectoryManager::createPSimHits(const TrackerLayer& layer,
00543                                   const ParticlePropagator& PP,
00544                                   std::map<double,PSimHit>& theHitMap,
00545                                   int trackID, int partID) {
00546 
00547   // Propagate the particle coordinates to the closest tracker detector(s) 
00548   // in this layer and create the PSimHit(s)
00549 
00550   //  const MagneticField& mf = MagneticFieldMap::instance()->magneticField();
00551   // This solution is actually much faster !
00552   LocalMagneticField mf(PP.getMagneticField());
00553   AnalyticalPropagator alongProp(&mf, anyDirection);
00554   InsideBoundsMeasurementEstimator est;
00555 
00556   typedef GeometricSearchDet::DetWithState   DetWithState;
00557   const DetLayer* tkLayer = detLayer(layer,PP.Z());
00558 
00559   TrajectoryStateOnSurface trajState = makeTrajectoryState( tkLayer, PP, &mf);
00560   float thickness = theMaterialEffects ? theMaterialEffects->thickness() : 0.;
00561   float eloss = theMaterialEffects ? theMaterialEffects->energyLoss() : 0.;
00562 
00563   // Find, in the corresponding layers, the detectors compatible 
00564   // with the current track 
00565   std::vector<DetWithState> compat 
00566     = tkLayer->compatibleDets( trajState, alongProp, est);
00567 
00568   // And create the corresponding PSimHits
00569   std::map<double,PSimHit> theTrackHits;
00570   for (std::vector<DetWithState>::const_iterator i=compat.begin(); i!=compat.end(); i++) {
00571     // Correct Eloss for last 3 rings of TEC (thick sensors, 0.05 cm)
00572     // Disgusting fudge factor ! 
00573       makePSimHits( i->first, i->second, theHitMap, trackID, eloss, thickness, partID);
00574   }
00575 }
00576 
00577 TrajectoryStateOnSurface 
00578 TrajectoryManager::makeTrajectoryState( const DetLayer* layer, 
00579                                         const ParticlePropagator& pp,
00580                                         const MagneticField* field) const
00581 {
00582   GlobalPoint  pos( pp.X(), pp.Y(), pp.Z());
00583   GlobalVector mom( pp.Px(), pp.Py(), pp.Pz());
00584   ReferenceCountingPointer<TangentPlane> plane = layer->surface().tangentPlane(pos);
00585   return TrajectoryStateOnSurface
00586     (GlobalTrajectoryParameters( pos, mom, TrackCharge( pp.charge()), field), *plane);
00587 }
00588 
00589 void 
00590 TrajectoryManager::makePSimHits( const GeomDet* det, 
00591                                  const TrajectoryStateOnSurface& ts,
00592                                  std::map<double,PSimHit>& theHitMap,
00593                                  int tkID, float el, float thick, int pID ) 
00594 {
00595 
00596   std::vector< const GeomDet*> comp = det->components();
00597   if (!comp.empty()) {
00598     for (std::vector< const GeomDet*>::const_iterator i = comp.begin();
00599          i != comp.end(); i++) {
00600       const GeomDetUnit* du = dynamic_cast<const GeomDetUnit*>(*i);
00601       if (du != 0)
00602         theHitMap.insert(theHitMap.end(),makeSinglePSimHit( *du, ts, tkID, el, thick, pID));
00603     }
00604   }
00605   else {
00606     const GeomDetUnit* du = dynamic_cast<const GeomDetUnit*>(det);
00607     if (du != 0)
00608       theHitMap.insert(theHitMap.end(),makeSinglePSimHit( *du, ts, tkID, el, thick, pID));
00609   }
00610 
00611 
00612 }
00613 
00614 std::pair<double,PSimHit> 
00615 TrajectoryManager::makeSinglePSimHit( const GeomDetUnit& det,
00616                                       const TrajectoryStateOnSurface& ts, 
00617                                       int tkID, float el, float thick, int pID) const
00618 {
00619 
00620   const float onSurfaceTolarance = 0.01; // 10 microns
00621 
00622   LocalPoint lpos;
00623   LocalVector lmom;
00624   if ( fabs( det.toLocal(ts.globalPosition()).z()) < onSurfaceTolarance) {
00625     lpos = ts.localPosition();
00626     lmom = ts.localMomentum();
00627   }
00628   else {
00629     HelixArbitraryPlaneCrossing crossing( ts.globalPosition().basicVector(),
00630                                           ts.globalMomentum().basicVector(),
00631                                           ts.transverseCurvature(),
00632                                           anyDirection);
00633     std::pair<bool,double> path = crossing.pathLength(det.surface());
00634     if (!path.first) {
00635       // edm::LogWarning("FastTracking") << "TrajectoryManager ERROR: crossing with det failed, skipping PSimHit";
00636       return  std::pair<double,PSimHit>(0.,PSimHit());
00637     }
00638     lpos = det.toLocal( GlobalPoint( crossing.position(path.second)));
00639     lmom = det.toLocal( GlobalVector( crossing.direction(path.second)));
00640     lmom = lmom.unit() * ts.localMomentum().mag();
00641   }
00642 
00643   // The module (half) thickness 
00644   const BoundPlane& theDetPlane = det.surface();
00645   float halfThick = 0.5*theDetPlane.bounds().thickness();
00646   // The Energy loss rescaled to the module thickness
00647   float eloss = el;
00648   if ( thick > 0. ) {
00649     // Total thickness is in radiation lengths, 1 radlen = 9.36 cm
00650     // Sensitive module thickness is about 30 microns larger than 
00651     // the module thickness itself
00652     eloss *= (2.* halfThick - 0.003) / (9.36 * thick);
00653   }
00654   // The entry and exit points, and the time of flight
00655   float pZ = lmom.z();
00656   LocalPoint entry = lpos + (-halfThick/pZ) * lmom;
00657   LocalPoint exit = lpos + halfThick/pZ * lmom;
00658   float tof = ts.globalPosition().mag() / 30. ; // in nanoseconds, FIXME: very approximate
00659 
00660   // If a hadron suffered a nuclear interaction, just assign the hits of the closest 
00661   // daughter to the mother's track. The same applies to a charged particle decay into
00662   // another charged particle.
00663   int localTkID = tkID;
00664   if ( mySimEvent->track(tkID).mother().closestDaughterId() == tkID )
00665     localTkID = mySimEvent->track(tkID).mother().id();
00666 
00667   // FIXME: fix the track ID and the particle ID
00668   PSimHit hit( entry, exit, lmom.mag(), tof, eloss, pID,
00669                   det.geographicalId().rawId(), localTkID,
00670                   lmom.theta(),
00671                   lmom.phi());
00672 
00673   // Check that the PSimHit is physically on the module!
00674   unsigned subdet = DetId(hit.detUnitId()).subdetId(); 
00675   double boundX = theDetPlane.bounds().width()/2.;
00676   double boundY = theDetPlane.bounds().length()/2.;
00677 
00678   // Special treatment for TID and TEC trapeziodal modules
00679   if ( subdet == 4 || subdet == 6 ) 
00680     boundX *=  1. - hit.localPosition().y()/theDetPlane.position().perp();
00681 
00682 #ifdef FAMOS_DEBUG
00683   unsigned detid  = DetId(hit.detUnitId()).rawId();
00684   unsigned stereo = 0;
00685   unsigned theLayer = 0;
00686   unsigned theRing = 0;
00687   switch (subdet) { 
00688   case 1: 
00689     {
00690       PXBDetId module(detid);
00691       theLayer = module.layer();
00692       std::cout << "\tPixel Barrel Layer " << theLayer << std::endl;
00693       stereo = 1;
00694       break;
00695     }
00696   case 2: 
00697     {
00698       PXFDetId module(detid);
00699       theLayer = module.disk();
00700       std::cout << "\tPixel Forward Disk " << theLayer << std::endl;
00701       stereo = 1;
00702       break;
00703     }
00704   case 3:
00705     {
00706       TIBDetId module(detid);
00707       theLayer  = module.layer();
00708       std::cout << "\tTIB Layer " << theLayer << std::endl;
00709       stereo = module.stereo();
00710       break;
00711     }
00712   case 4:
00713     {
00714       TIDDetId module(detid);
00715       theLayer = module.wheel();
00716       theRing  = module.ring();
00717       unsigned int theSide = module.side();
00718       if ( theSide == 1 ) 
00719         std::cout << "\tTID Petal Back " << std::endl; 
00720       else
00721         std::cout << "\tTID Petal Front" << std::endl; 
00722       std::cout << "\tTID Layer " << theLayer << std::endl;
00723       std::cout << "\tTID Ring " << theRing << std::endl;
00724       stereo = module.stereo();
00725       break;
00726     }
00727   case 5:
00728     {
00729       TOBDetId module(detid);
00730       theLayer  = module.layer();
00731       stereo = module.stereo();
00732       std::cout << "\tTOB Layer " << theLayer << std::endl;
00733       break;
00734     }
00735   case 6:
00736     {
00737       TECDetId module(detid);
00738       theLayer = module.wheel();
00739       theRing  = module.ring();
00740       unsigned int theSide = module.petal()[0];
00741       if ( theSide == 1 ) 
00742         std::cout << "\tTEC Petal Back " << std::endl; 
00743       else
00744         std::cout << "\tTEC Petal Front" << std::endl; 
00745       std::cout << "\tTEC Layer " << theLayer << std::endl;
00746       std::cout << "\tTEC Ring " << theRing << std::endl;
00747       stereo = module.stereo();
00748       break;
00749     }
00750   default:
00751     {
00752       stereo = 0;
00753       break;
00754     }
00755   }
00756   
00757   std::cout << "Thickness = " << 2.*halfThick-0.003 << "; " << thick * 9.36 << std::endl
00758             << "Length    = " << det.surface().bounds().length() << std::endl
00759             << "Width     = " << det.surface().bounds().width() << std::endl;
00760     
00761   std::cout << "Hit position = " 
00762             << hit.localPosition().x() << " " 
00763             << hit.localPosition().y() << " " 
00764             << hit.localPosition().z() << std::endl;
00765 #endif
00766 
00767   // Check if the hit is on the physical volume of the module
00768   // (It happens that it is not, in the case of double sided modules,
00769   //  because the envelope of the gluedDet is larger than each of 
00770   //  the mono and the stereo modules)
00771 
00772   double dist = 0.;
00773   GlobalPoint IP (mySimEvent->track(localTkID).vertex().position().x(),
00774                   mySimEvent->track(localTkID).vertex().position().y(),
00775                   mySimEvent->track(localTkID).vertex().position().z());
00776 
00777   dist = ( fabs(hit.localPosition().x()) > boundX  || 
00778            fabs(hit.localPosition().y()) > boundY ) ?  
00779     // Will be used later as a flag to reject the PSimHit!
00780     -( det.surface().toGlobal(hit.localPosition()) - IP ).mag2() 
00781     : 
00782     // These hits are kept!
00783      ( det.surface().toGlobal(hit.localPosition()) - IP ).mag2();
00784 
00785   // Fill Histos (~poor man event display)
00786   /* 
00787      GlobalPoint gpos( det.toGlobal(hit.localPosition()));
00788      myHistos->fill("h300",gpos.x(),gpos.y());
00789      if ( sin(gpos.phi()) > 0. ) 
00790      myHistos->fill("h301",gpos.z(),gpos.perp());
00791      else
00792      myHistos->fill("h301",gpos.z(),-gpos.perp());
00793   */
00794   
00795   return std::pair<double,PSimHit>(dist,hit);
00796 
00797 }
00798 
00799 void 
00800 TrajectoryManager::initializeLayerMap()
00801 {
00802 
00803 // These are the BoundSurface&, the BoundDisk* and the BoundCylinder* for that layer
00804 //   const BoundSurface& theSurface = layer.surface();
00805 //   BoundDisk* theDisk = layer.disk();  // non zero for endcaps
00806 //   BoundCylinder* theCylinder = layer.cylinder(); // non zero for barrel
00807 //   int theLayer = layer.layerNumber(); // 1->3 PixB, 4->5 PixD, 
00808 //                                       // 6->9 TIB, 10->12 TID, 
00809 //                                       // 13->18 TOB, 19->27 TEC
00810 
00813 
00814   std::vector< BarrelDetLayer*>   barrelLayers = 
00815     theGeomSearchTracker->barrelLayers();
00816   LogDebug("FastTracking") << "Barrel DetLayer dump: ";
00817   for (std::vector< BarrelDetLayer*>::const_iterator bl=barrelLayers.begin();
00818        bl != barrelLayers.end(); ++bl) {
00819     LogDebug("FastTracking")<< "radius " << (**bl).specificSurface().radius(); 
00820   }
00821 
00822   std::vector< ForwardDetLayer*>  posForwardLayers = 
00823     theGeomSearchTracker->posForwardLayers();
00824   LogDebug("FastTracking") << "Positive Forward DetLayer dump: ";
00825   for (std::vector< ForwardDetLayer*>::const_iterator fl=posForwardLayers.begin();
00826        fl != posForwardLayers.end(); ++fl) {
00827     LogDebug("FastTracking") << "Z pos "
00828                             << (**fl).surface().position().z()
00829                             << " radii " 
00830                             << (**fl).specificSurface().innerRadius() 
00831                             << ", " 
00832                             << (**fl).specificSurface().outerRadius(); 
00833   }
00834 
00835   const float rTolerance = 1.5;
00836   const float zTolerance = 3.;
00837 
00838   LogDebug("FastTracking")<< "Dump of TrackerInteractionGeometry cylinders:";
00839   for( std::list<TrackerLayer>::const_iterator i=_theGeometry->cylinderBegin();
00840        i!=_theGeometry->cylinderEnd(); ++i) {
00841     const BoundCylinder* cyl = i->cylinder();
00842     const BoundDisk* disk = i->disk();
00843 
00844     LogDebug("FastTracking") << "Famos Layer no " << i->layerNumber()
00845                             << " is sensitive? " << i->sensitive()
00846                             << " pos " << i->surface().position();
00847     if (!i->sensitive()) continue;
00848 
00849     if (cyl != 0) {
00850       LogDebug("FastTracking") << " cylinder radius " << cyl->radius();
00851       bool found = false;
00852       for (std::vector< BarrelDetLayer*>::const_iterator 
00853              bl=barrelLayers.begin(); bl != barrelLayers.end(); ++bl) {
00854         if (fabs( cyl->radius() - (**bl).specificSurface().radius()) < rTolerance) {
00855           theLayerMap[i->layerNumber()] = *bl;
00856           found = true;
00857           LogDebug("FastTracking")<< "Corresponding DetLayer found with radius "
00858                                  << (**bl).specificSurface().radius();
00859           break;
00860         }
00861       }
00862       if (!found) {
00863         edm::LogWarning("FastTracking") << " Trajectory manager FAILED to find a corresponding DetLayer!";
00864       }
00865     }
00866     else {
00867       LogDebug("FastTracking") << " disk radii " << disk->innerRadius() 
00868                  << ", " << disk->outerRadius();
00869       bool found = false;
00870       for (std::vector< ForwardDetLayer*>::const_iterator fl=posForwardLayers.begin();
00871            fl != posForwardLayers.end(); ++fl) {
00872         if (fabs( disk->position().z() - (**fl).surface().position().z()) < zTolerance) {
00873           theLayerMap[i->layerNumber()] = *fl;
00874           found = true;
00875           LogDebug("FastTracking") << "Corresponding DetLayer found with Z pos "
00876                                   << (**fl).surface().position().z()
00877                                   << " and radii " 
00878                                   << (**fl).specificSurface().innerRadius() 
00879                                   << ", " 
00880                                   << (**fl).specificSurface().outerRadius(); 
00881           break;
00882         }
00883       }
00884       if (!found) {
00885         edm::LogWarning("FastTracking") << "FAILED to find a corresponding DetLayer!";
00886       }
00887     }
00888   }
00889 
00890   // Put the negative layers in the same map but with an offset
00891   std::vector< ForwardDetLayer*>  negForwardLayers = theGeomSearchTracker->negForwardLayers();
00892   for (std::vector< ForwardDetLayer*>::const_iterator nl=negForwardLayers.begin();
00893        nl != negForwardLayers.end(); ++nl) {
00894     for (int i=0; i<=theNegLayerOffset; i++) {
00895       if (theLayerMap[i] == 0) continue;
00896       if ( fabs( (**nl).surface().position().z() +theLayerMap[i]-> surface().position().z()) < zTolerance) {
00897         theLayerMap[i+theNegLayerOffset] = *nl;
00898         break;
00899       }
00900     }
00901   }  
00902 
00903 }
00904 
00905 const DetLayer*  
00906 TrajectoryManager::detLayer( const TrackerLayer& layer, float zpos) const
00907 {
00908   if (zpos > 0 || !layer.forward() ) return theLayerMap[layer.layerNumber()];
00909   else return theLayerMap[layer.layerNumber()+theNegLayerOffset];
00910 }
00911 
00912 void 
00913 TrajectoryManager::loadSimHits(edm::PSimHitContainer & c) const
00914 {
00915 
00916   std::map<unsigned,std::map<double,PSimHit> >::const_iterator itrack = thePSimHits.begin();
00917   std::map<unsigned,std::map<double,PSimHit> >::const_iterator itrackEnd = thePSimHits.end();
00918   for ( ; itrack != itrackEnd; ++itrack ) {
00919     std::map<double,PSimHit>::const_iterator it = (itrack->second).begin();
00920     std::map<double,PSimHit>::const_iterator itEnd = (itrack->second).end();
00921     for( ; it!= itEnd; ++it) { 
00922       /*
00923       DetId theDetUnitId((it->second).detUnitId());
00924       const GeomDet* theDet = theGeomTracker->idToDet(theDetUnitId);
00925       std::cout << "Track/z/r after : "
00926                 << (it->second).trackId() << " " 
00927                 << theDet->surface().toGlobal((it->second).localPosition()).z() << " " 
00928                 << theDet->surface().toGlobal((it->second).localPosition()).perp() << std::endl;
00929       */
00930       // Keep only those hits that are on the physical volume of a module
00931       // (The other hits have been assigned a negative <double> value. 
00932       if ( it->first > 0. ) c.push_back(it->second); 
00933     }
00934   }
00935 
00936 }