CMS 3D CMS Logo

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