CMS 3D CMS Logo

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