CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_10_patch1/src/FastSimulation/Event/src/FBaseSimEvent.cc

Go to the documentation of this file.
00001 //HepMC Headers
00002 #include "HepMC/GenEvent.h"
00003 #include "HepMC/GenVertex.h"
00004 #include "HepMC/GenParticle.h"
00005 
00006 //Framework Headers
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 
00009 //CMSSW Data Formats
00010 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
00011 #include "DataFormats/Candidate/interface/Candidate.h"
00012 
00013 //FAMOS Headers
00014 #include "FastSimulation/Event/interface/FBaseSimEvent.h"
00015 #include "FastSimulation/Event/interface/FSimTrack.h"
00016 #include "FastSimulation/Event/interface/FSimVertex.h"
00017 #include "FastSimulation/Event/interface/KineParticleFilter.h"
00018 #include "FastSimulation/BaseParticlePropagator/interface/BaseParticlePropagator.h"
00019 #include "FastSimulation/Event/interface/BetaFuncPrimaryVertexGenerator.h"
00020 #include "FastSimulation/Event/interface/GaussianPrimaryVertexGenerator.h"
00021 #include "FastSimulation/Event/interface/FlatPrimaryVertexGenerator.h"
00022 #include "FastSimulation/Event/interface/NoPrimaryVertexGenerator.h"
00023 
00024 #include "FastSimDataFormats/NuclearInteractions/interface/FSimVertexType.h"
00025 
00026 using namespace HepPDT;
00027 
00028 // system include
00029 #include <iostream>
00030 #include <iomanip>
00031 #include <map>
00032 #include <string>
00033 
00034 FBaseSimEvent::FBaseSimEvent(const edm::ParameterSet& kine) 
00035   :
00036   nSimTracks(0),
00037   nSimVertices(0),
00038   nGenParticles(0),
00039   nChargedParticleTracks(0),
00040   initialSize(5000),
00041   random(0)
00042 {
00043 
00044   theVertexGenerator = new NoPrimaryVertexGenerator();
00045   theBeamSpot = math::XYZPoint(0.0,0.0,0.0);
00046 
00047   // Initialize the vectors of particles and vertices
00048   theGenParticles = new std::vector<HepMC::GenParticle*>(); 
00049   theSimTracks = new std::vector<FSimTrack>;
00050   theSimVertices = new std::vector<FSimVertex>;
00051   theChargedTracks = new std::vector<unsigned>();
00052   theFSimVerticesType = new FSimVertexTypeCollection();
00053 
00054   // Reserve some size to avoid mutiple copies
00055   /* */
00056   theSimTracks->resize(initialSize);
00057   theSimVertices->resize(initialSize);
00058   theGenParticles->resize(initialSize);
00059   theChargedTracks->resize(initialSize);
00060   theFSimVerticesType->resize(initialSize);
00061   theTrackSize = initialSize;
00062   theVertexSize = initialSize;
00063   theGenSize = initialSize;
00064   theChargedSize = initialSize;
00065   /* */
00066 
00067   // Initialize the Particle filter
00068   myFilter = new KineParticleFilter(kine);
00069 
00070 }
00071 
00072 FBaseSimEvent::FBaseSimEvent(const edm::ParameterSet& vtx,
00073                              const edm::ParameterSet& kine,
00074                              const RandomEngine* engine) 
00075   :
00076   nSimTracks(0),
00077   nSimVertices(0),
00078   nGenParticles(0),
00079   nChargedParticleTracks(0), 
00080   initialSize(5000),
00081   theVertexGenerator(0), 
00082   random(engine)
00083 {
00084 
00085   // Initialize the vertex generator
00086   std::string vtxType = vtx.getParameter<std::string>("type");
00087   if ( vtxType == "Gaussian" ) 
00088     theVertexGenerator = new GaussianPrimaryVertexGenerator(vtx,random);
00089   else if ( vtxType == "Flat" ) 
00090     theVertexGenerator = new FlatPrimaryVertexGenerator(vtx,random);
00091   else if ( vtxType == "BetaFunc" )
00092     theVertexGenerator = new BetaFuncPrimaryVertexGenerator(vtx,random);
00093   else
00094     theVertexGenerator = new NoPrimaryVertexGenerator();
00095   // Initialize the beam spot, if not read from the DataBase
00096   theBeamSpot = math::XYZPoint(0.0,0.0,0.0);
00097 
00098   // Initialize the distance from (0,0,0) after which *generated* particles are 
00099   // no longer considered - because the mother could have interacted before.
00100   // unit : cm x cm
00101   lateVertexPosition = 2.5*2.5;
00102 
00103   // Initialize the vectors of particles and vertices
00104   theGenParticles = new std::vector<HepMC::GenParticle*>(); 
00105   theSimTracks = new std::vector<FSimTrack>;
00106   theSimVertices = new std::vector<FSimVertex>;
00107   theChargedTracks = new std::vector<unsigned>();
00108   theFSimVerticesType = new FSimVertexTypeCollection();
00109 
00110   // Reserve some size to avoid mutiple copies
00111   /* */
00112   theSimTracks->resize(initialSize);
00113   theSimVertices->resize(initialSize);
00114   theGenParticles->resize(initialSize);
00115   theChargedTracks->resize(initialSize);
00116   theFSimVerticesType->resize(initialSize);
00117   theTrackSize = initialSize;
00118   theVertexSize = initialSize;
00119   theGenSize = initialSize;
00120   theChargedSize = initialSize;
00121   /* */
00122 
00123   // Initialize the Particle filter
00124   myFilter = new KineParticleFilter(kine);
00125 
00126 }
00127  
00128 FBaseSimEvent::~FBaseSimEvent(){
00129 
00130   // Clear the vectors
00131   theGenParticles->clear();
00132   theSimTracks->clear();
00133   theSimVertices->clear();
00134   theChargedTracks->clear();
00135   theFSimVerticesType->clear();
00136 
00137   // Delete 
00138   delete theGenParticles;
00139   delete theSimTracks;
00140   delete theSimVertices;
00141   delete theChargedTracks;
00142   delete theFSimVerticesType;
00143   delete myFilter;
00144 
00145 }
00146 
00147 void 
00148 FBaseSimEvent::initializePdt(const HepPDT::ParticleDataTable* aPdt) { 
00149 
00150   pdt = aPdt; 
00151 
00152 }
00153 
00154 /*
00155 const HepPDT::ParticleDataTable*
00156 FBaseSimEvent::theTable() const {
00157   return pdt;
00158 }
00159 */
00160 
00161 void
00162 FBaseSimEvent::fill(const HepMC::GenEvent& myGenEvent) {
00163   
00164   // Clear old vectors
00165   clear();
00166 
00167   // Add the particles in the FSimEvent
00168   addParticles(myGenEvent);
00169 
00170   /*
00171   std::cout << "The MC truth! " << std::endl;
00172   printMCTruth(myGenEvent);
00173 
00174   std::cout << std::endl  << "The FAMOS event! " << std::endl;
00175   print();
00176   */
00177 
00178 }
00179 
00180 void
00181 FBaseSimEvent::fill(const reco::GenParticleCollection& myGenParticles) {
00182   
00183   // Clear old vectors
00184   clear();
00185 
00186   // Add the particles in the FSimEvent
00187   addParticles(myGenParticles);
00188 
00189 }
00190 
00191 void
00192 FBaseSimEvent::fill(const std::vector<SimTrack>& simTracks, 
00193                     const std::vector<SimVertex>& simVertices) {
00194 
00195   // Watch out there ! A SimVertex is in mm (stupid), 
00196   //            while a FSimVertex is in cm (clever).
00197   
00198   clear();
00199 
00200   unsigned nVtx = simVertices.size();
00201   unsigned nTks = simTracks.size();
00202 
00203   // Empty event, do nothin'
00204   if ( nVtx == 0 ) return;
00205 
00206   // Two arrays for internal use.
00207   std::vector<int> myVertices(nVtx,-1);
00208   std::vector<int> myTracks(nTks,-1);
00209 
00210   // create a map associating geant particle id and position in the 
00211   // event SimTrack vector
00212   
00213   std::map<unsigned, unsigned> geantToIndex;
00214   for( unsigned it=0; it<simTracks.size(); ++it ) {
00215     geantToIndex[ simTracks[it].trackId() ] = it;
00216   }  
00217 
00218   // Create also a map associating a SimTrack with its endVertex
00219   /*
00220   std::map<unsigned, unsigned> endVertex;
00221   for ( unsigned iv=0; iv<simVertices.size(); ++iv ) { 
00222     endVertex[ simVertices[iv].parentIndex() ] = iv;
00223   }
00224   */
00225 
00226   // Set the main vertex for the kine particle filter
00227   // SimVertices were in mm until 110_pre2
00228   //  HepLorentzVector primaryVertex = simVertices[0].position()/10.;
00229   // SImVertices are now in cm
00230   // Also : position is copied until SimVertex switches to Mathcore.
00231   //  XYZTLorentzVector primaryVertex = simVertices[0].position();
00232   // The next 5 lines to be then replaced by the previous line
00233   XYZTLorentzVector primaryVertex(simVertices[0].position().x(),
00234                                   simVertices[0].position().y(),
00235                                   simVertices[0].position().z(),
00236                                   simVertices[0].position().t());
00237   //
00238   myFilter->setMainVertex(primaryVertex);
00239   // Add the main vertex to the list.
00240   addSimVertex(myFilter->vertex(), -1, FSimVertexType::PRIMARY_VERTEX);
00241   myVertices[0] = 0;
00242 
00243   for( unsigned trackId=0; trackId<nTks; ++trackId ) {
00244 
00245     // The track
00246     const SimTrack& track = simTracks[trackId];
00247     //    std::cout << std::endl << "SimTrack " << trackId << " " << track << std::endl;
00248 
00249     // The origin vertex
00250     int vertexId = track.vertIndex();
00251     const SimVertex& vertex = simVertices[vertexId];
00252     //std::cout << "Origin vertex " << vertexId << " " << vertex << std::endl;
00253 
00254     // The mother track 
00255     int motherId = -1;
00256     if( !vertex.noParent() ) { // there is a parent to this vertex
00257       // geant id of the mother
00258       unsigned motherGeantId =   vertex.parentIndex(); 
00259       std::map<unsigned, unsigned >::iterator association  
00260         = geantToIndex.find( motherGeantId );
00261       if(association != geantToIndex.end() )
00262         motherId = association->second;
00263     }
00264     int originId = motherId == - 1 ? -1 : myTracks[motherId];
00265     //std::cout << "Origin id " << originId << std::endl;
00266 
00267     /*
00268     if ( endVertex.find(trackId) != endVertex.end() ) 
00269       std::cout << "End vertex id = " << endVertex[trackId] << std::endl;
00270     else
00271       std::cout << "No endVertex !!! " << std::endl;
00272     std::cout << "Tracker surface position " << track.trackerSurfacePosition() << std::endl;
00273     */
00274 
00275     // Add the vertex (if it does not already exist!)
00276     XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
00277                                vertex.position().pz(),vertex.position().e());
00278     if ( myVertices[vertexId] == -1 ) 
00279       // Momentum and position are copied until SimTrack and SimVertex
00280       // switch to Mathcore.
00281       //      myVertices[vertexId] = addSimVertex(vertex.position(),originId); 
00282       // The next line to be then replaced by the previous line
00283       myVertices[vertexId] = addSimVertex(position,originId); 
00284 
00285     // Add the track (with protection for brem'ing electrons and muons)
00286     int motherType = motherId == -1 ? 0 : simTracks[motherId].type();
00287 
00288     bool notBremInDetector =
00289       (abs(motherType) != 11 && abs(motherType) != 13) ||
00290       motherType != track.type() ||
00291       position.Perp2() < lateVertexPosition ;
00292 
00293     if ( notBremInDetector ) {
00294       // Momentum and position are copied until SimTrack and SimVertex
00295       // switch to Mathcore.
00296       //      RawParticle part(track.momentum(), vertex.position());
00297       // The next 3 lines to be then replaced by the previous line
00298       XYZTLorentzVector momentum(track.momentum().px(),track.momentum().py(),
00299                                  track.momentum().pz(),track.momentum().e());
00300       RawParticle part(momentum,position);
00301       //
00302       part.setID(track.type()); 
00303       //std::cout << "Ctau  = " << part.PDGcTau() << std::endl;
00304       // Don't save tracks that have decayed immediately but for which no daughters
00305       // were saved (probably due to cuts on E, pT and eta)
00306       //  if ( part.PDGcTau() > 0.1 || endVertex.find(trackId) != endVertex.end() ) 
00307       myTracks[trackId] = addSimTrack(&part,myVertices[vertexId],track.genpartIndex());
00308       if ( myTracks[trackId] >= 0 ) { 
00309         (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
00310         (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
00311       }
00312     } else {
00313 
00314       myTracks[trackId] = myTracks[motherId];
00315       if ( myTracks[trackId] >= 0 ) { 
00316         (*theSimTracks)[ myTracks[trackId] ].setTkPosition(track.trackerSurfacePosition());
00317         (*theSimTracks)[ myTracks[trackId] ].setTkMomentum(track.trackerSurfaceMomentum());
00318       }
00319     }
00320     
00321   }
00322 
00323   // Now loop over the remaining end vertices !
00324   for( unsigned vertexId=0; vertexId<nVtx; ++vertexId ) {
00325 
00326     // if the vertex is already saved, just ignore.
00327     if ( myVertices[vertexId] != -1 ) continue;
00328 
00329     // The yet unused vertex
00330     const SimVertex& vertex = simVertices[vertexId];
00331 
00332     // The mother track 
00333     int motherId = -1;
00334     if( !vertex.noParent() ) { // there is a parent to this vertex
00335 
00336       // geant id of the mother
00337       unsigned motherGeantId =   vertex.parentIndex(); 
00338       std::map<unsigned, unsigned >::iterator association  
00339         = geantToIndex.find( motherGeantId );
00340       if(association != geantToIndex.end() )
00341         motherId = association->second;
00342     }
00343     int originId = motherId == - 1 ? -1 : myTracks[motherId];
00344 
00345     // Add the vertex
00346     // Momentum and position are copied until SimTrack and SimVertex
00347     // switch to Mathcore.
00348     //    myVertices[vertexId] = addSimVertex(vertex.position(),originId);
00349     // The next 3 lines to be then replaced by the previous line
00350     XYZTLorentzVector position(vertex.position().px(),vertex.position().py(),
00351                                vertex.position().pz(),vertex.position().e());
00352     myVertices[vertexId] = addSimVertex(position,originId); 
00353   }
00354 
00355   // Finally, propagate all particles to the calorimeters
00356   BaseParticlePropagator myPart;
00357   XYZTLorentzVector mom;
00358   XYZTLorentzVector pos;
00359 
00360 
00361   // Loop over the tracks
00362   for( int fsimi=0; fsimi < (int)nTracks() ; ++fsimi) {
00363 
00364     
00365     FSimTrack& myTrack = track(fsimi);
00366     double trackerSurfaceTime = myTrack.vertex().position().t() 
00367                               + myTrack.momentum().e()/myTrack.momentum().pz()
00368                               * ( myTrack.trackerSurfacePosition().z()
00369                                 - myTrack.vertex().position().z() );
00370     pos = XYZTLorentzVector(myTrack.trackerSurfacePosition().x(),
00371                             myTrack.trackerSurfacePosition().y(),
00372                             myTrack.trackerSurfacePosition().z(),
00373                                     trackerSurfaceTime);
00374     mom = XYZTLorentzVector(myTrack.trackerSurfaceMomentum().x(),
00375                             myTrack.trackerSurfaceMomentum().y(),
00376                             myTrack.trackerSurfaceMomentum().z(),
00377                             myTrack.trackerSurfaceMomentum().t());
00378 
00379     if ( mom.T() >  0. ) {  
00380       // The particle to be propagated
00381       myPart = BaseParticlePropagator(RawParticle(mom,pos),0.,0.,4.);
00382       myPart.setCharge(myTrack.charge());
00383       
00384       // Propagate to Preshower layer 1
00385       myPart.propagateToPreshowerLayer1(false);
00386       if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
00387         myTrack.setLayer1(myPart,myPart.getSuccess());
00388       
00389       // Propagate to Preshower Layer 2 
00390       myPart.propagateToPreshowerLayer2(false);
00391       if ( myTrack.notYetToEndVertex(myPart.vertex()) && myPart.getSuccess()>0 )
00392         myTrack.setLayer2(myPart,myPart.getSuccess());
00393       
00394       // Propagate to Ecal Endcap
00395       myPart.propagateToEcalEntrance(false);
00396       if ( myTrack.notYetToEndVertex(myPart.vertex()) )
00397         myTrack.setEcal(myPart,myPart.getSuccess());
00398       
00399       // Propagate to HCAL entrance
00400       myPart.propagateToHcalEntrance(false);
00401       if ( myTrack.notYetToEndVertex(myPart.vertex()) )
00402         myTrack.setHcal(myPart,myPart.getSuccess());
00403       
00404       // Attempt propagation to HF for low pt and high eta 
00405       if ( myPart.cos2ThetaV()>0.8 || mom.T() < 3. ) {
00406         // Propagate to VFCAL entrance
00407         myPart.propagateToVFcalEntrance(false);
00408         if ( myTrack.notYetToEndVertex(myPart.vertex()) )
00409         myTrack.setVFcal(myPart,myPart.getSuccess());
00410         
00411         // Otherwise propagate to the HCAL exit and HO.
00412       } else { 
00413         // Propagate to HCAL exit
00414         myPart.propagateToHcalExit(false);
00415         if ( myTrack.notYetToEndVertex(myPart.vertex()) )
00416           myTrack.setHcalExit(myPart,myPart.getSuccess());     
00417         // Propagate to HOLayer entrance
00418         myPart.setMagneticField(0);
00419         myPart.propagateToHOLayer(false);
00420         if ( myTrack.notYetToEndVertex(myPart.vertex()) )
00421           myTrack.setHO(myPart,myPart.getSuccess());
00422       } 
00423     }
00424   }
00425 }
00426 
00427 
00428 void
00429 FBaseSimEvent::addParticles(const HepMC::GenEvent& myGenEvent) {
00430 
00432   int genEventSize = myGenEvent.particles_size();
00433   std::vector<int> myGenVertices(genEventSize, static_cast<int>(0));
00434 
00435   // If no particles, no work to be done !
00436   if ( myGenEvent.particles_empty() ) return;
00437 
00438   // Are there particles in the FSimEvent already ? 
00439   int offset = nGenParts();
00440 
00441   // Primary vertex (already smeared by the SmearedVtx module)
00442   HepMC::GenVertex* primaryVertex = *(myGenEvent.vertices_begin());
00443 
00444   // Beginning of workaround a bug in pythia particle gun
00445   unsigned primaryMother = primaryVertex->particles_in_size();
00446   if ( primaryMother ) {
00447     unsigned partId = (*(primaryVertex->particles_in_const_begin()))->pdg_id();
00448     if ( abs(partId) == 2212 ) primaryMother = 0;
00449   }
00450   // End of workaround a bug in pythia particle gun
00451 
00452   XYZTLorentzVector primaryVertexPosition(primaryVertex->position().x()/10.,
00453                                           primaryVertex->position().y()/10.,
00454                                           primaryVertex->position().z()/10.,
00455                                           primaryVertex->position().t()/10.);
00456   // Actually this is the true end of the workaround
00457   primaryVertexPosition *= (1-primaryMother);
00458   // THE END.
00459 
00460   // Smear the main vertex if needed
00461   // Now takes the origin from the database
00462   XYZTLorentzVector smearedVertex; 
00463   if ( primaryVertexPosition.Vect().Mag2() < 1E-16 ) {
00464     theVertexGenerator->generate();
00465     smearedVertex = XYZTLorentzVector(
00466       theVertexGenerator->X()-theVertexGenerator->beamSpot().X()+theBeamSpot.X(),
00467       theVertexGenerator->Y()-theVertexGenerator->beamSpot().Y()+theBeamSpot.Y(),
00468       theVertexGenerator->Z()-theVertexGenerator->beamSpot().Z()+theBeamSpot.Z(),
00469       0.);
00470   }
00471 
00472   // Set the main vertex
00473   myFilter->setMainVertex(primaryVertexPosition+smearedVertex);
00474 
00475   // This is the smeared main vertex
00476   int mainVertex = addSimVertex(myFilter->vertex(), -1, FSimVertexType::PRIMARY_VERTEX);
00477 
00478   HepMC::GenEvent::particle_const_iterator piter;
00479   HepMC::GenEvent::particle_const_iterator pbegin = myGenEvent.particles_begin();
00480   HepMC::GenEvent::particle_const_iterator pend = myGenEvent.particles_end();
00481 
00482   int initialBarcode = 0; 
00483   if ( pbegin != pend ) initialBarcode = (*pbegin)->barcode();
00484   // Loop on the particles of the generated event
00485   for ( piter = pbegin; piter != pend; ++piter ) {
00486 
00487     // This is the generated particle pointer - for the signal event only
00488     HepMC::GenParticle* p = *piter;
00489 
00490     if  ( !offset ) {
00491       (*theGenParticles)[nGenParticles++] = p;
00492       if ( nGenParticles/theGenSize*theGenSize == nGenParticles ) { 
00493         theGenSize *= 2;
00494         theGenParticles->resize(theGenSize);
00495       }
00496 
00497     }
00498 
00499     // Reject particles with late origin vertex (i.e., coming from late decays)
00500     // This should not happen, but one never knows what users may be up to!
00501     // For example exotic particles might decay late - keep the decay products in the case.
00502     XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
00503     HepMC::GenVertex* productionVertex = p->production_vertex();
00504     if ( productionVertex ) { 
00505       unsigned productionMother = productionVertex->particles_in_size();
00506       if ( productionMother ) {
00507         unsigned motherId = (*(productionVertex->particles_in_const_begin()))->pdg_id();
00508         if ( abs(motherId) < 1000000 ) 
00509           productionVertexPosition = 
00510             XYZTLorentzVector(productionVertex->position().x()/10.,
00511                               productionVertex->position().y()/10.,
00512                               productionVertex->position().z()/10.,
00513                               productionVertex->position().t()/10.) + smearedVertex;
00514       }
00515     }
00516     if ( !myFilter->accept(productionVertexPosition) ) continue;
00517 
00518     int abspdgId = abs(p->pdg_id());
00519     HepMC::GenVertex* endVertex = p->end_vertex();
00520 
00521     // Keep only: 
00522     // 1) Stable particles (watch out! New status code = 1001!)
00523     bool testStable = p->status()%1000==1;
00524     // Declare stable standard particles that decay after a macroscopic path length
00525     // (except if exotic)
00526     if ( p->status() == 2 && abspdgId < 1000000) {
00527       if ( endVertex ) { 
00528         XYZTLorentzVector decayPosition = 
00529           XYZTLorentzVector(endVertex->position().x()/10.,
00530                             endVertex->position().y()/10.,
00531                             endVertex->position().z()/10.,
00532                             endVertex->position().t()/10.) + smearedVertex;
00533         // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
00534         if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
00535       }
00536     }      
00537 
00538     // 2) or particles with stable daughters (watch out! New status code = 1001!)
00539     bool testDaugh = false;
00540     if ( !testStable && 
00541          p->status() == 2 &&
00542          endVertex && 
00543          endVertex->particles_out_size() ) { 
00544       HepMC::GenVertex::particles_out_const_iterator firstDaughterIt = 
00545         endVertex->particles_out_const_begin();
00546       HepMC::GenVertex::particles_out_const_iterator lastDaughterIt = 
00547         endVertex->particles_out_const_end();
00548       for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
00549         HepMC::GenParticle* daugh = *firstDaughterIt;
00550         if ( daugh->status()%1000==1 ) {
00551           // Check that it is not a "prompt electron or muon brem":
00552           if (abspdgId == 11 || abspdgId == 13) {
00553             if ( endVertex ) { 
00554               XYZTLorentzVector endVertexPosition = XYZTLorentzVector(endVertex->position().x()/10.,
00555                                                                       endVertex->position().y()/10.,
00556                                                                       endVertex->position().z()/10.,
00557                                                                       endVertex->position().t()/10.);
00558               // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
00559               if ( endVertexPosition.Perp2() < lateVertexPosition ) {
00560                 break;
00561               }
00562             }
00563           }
00564           testDaugh=true;
00565           break;
00566         }
00567       }
00568     }
00569 
00570     // 3) or particles that fly more than one micron.
00571     double dist = 0.;
00572     if ( !testStable && !testDaugh && p->production_vertex() ) {
00573       XYZTLorentzVector 
00574         productionVertexPosition(p->production_vertex()->position().x()/10.,
00575                                  p->production_vertex()->position().y()/10.,
00576                                  p->production_vertex()->position().z()/10.,
00577                                  p->production_vertex()->position().t()/10.);
00578       dist = (primaryVertexPosition-productionVertexPosition).Vect().Mag2();
00579     }
00580     bool testDecay = ( dist > 1e-8 ) ? true : false; 
00581 
00582     // Save the corresponding particle and vertices
00583     if ( testStable || testDaugh || testDecay ) {
00584       
00585       /*
00586       const HepMC::GenParticle* mother = p->production_vertex() ?
00587         *(p->production_vertex()->particles_in_const_begin()) : 0;
00588       */
00589 
00590       int motherBarcode = p->production_vertex() && 
00591         p->production_vertex()->particles_in_const_begin() !=
00592         p->production_vertex()->particles_in_const_end() ?
00593         (*(p->production_vertex()->particles_in_const_begin()))->barcode() : 0;
00594 
00595       int originVertex = 
00596         motherBarcode && myGenVertices[motherBarcode-initialBarcode] ?
00597         myGenVertices[motherBarcode-initialBarcode] : mainVertex;
00598 
00599       XYZTLorentzVector momentum(p->momentum().px(),
00600                                  p->momentum().py(),
00601                                  p->momentum().pz(),
00602                                  p->momentum().e());
00603       RawParticle part(momentum, vertex(originVertex).position());
00604       part.setID(p->pdg_id());
00605 
00606       // Add the particle to the event and to the various lists
00607       
00608       int theTrack = testStable && p->end_vertex() ? 
00609         // The particle is scheduled to decay
00610         addSimTrack(&part,originVertex, nGenParts()-offset,p->end_vertex()) :
00611         // The particle is not scheduled to decay 
00612         addSimTrack(&part,originVertex, nGenParts()-offset);
00613 
00614       if ( 
00615           // This one deals with particles with no end vertex
00616           !p->end_vertex() ||
00617           // This one deals with particles that have a pre-defined
00618           // decay proper time, but have not decayed yet
00619           ( testStable && p->end_vertex() &&  !p->end_vertex()->particles_out_size() ) 
00620           // In both case, just don't add a end vertex in the FSimEvent 
00621           ) continue; 
00622       
00623       // Add the vertex to the event and to the various lists
00624       XYZTLorentzVector decayVertex = 
00625         XYZTLorentzVector(p->end_vertex()->position().x()/10.,
00626                           p->end_vertex()->position().y()/10.,
00627                           p->end_vertex()->position().z()/10.,
00628                           p->end_vertex()->position().t()/10.) +
00629         smearedVertex;
00630       //        vertex(mainVertex).position();
00631       int theVertex = addSimVertex(decayVertex,theTrack, FSimVertexType::DECAY_VERTEX);
00632 
00633       if ( theVertex != -1 ) myGenVertices[p->barcode()-initialBarcode] = theVertex;
00634 
00635       // There we are !
00636     }
00637   }
00638 
00639 }
00640 
00641 void
00642 FBaseSimEvent::addParticles(const reco::GenParticleCollection& myGenParticles) {
00643 
00644   // If no particles, no work to be done !
00645   unsigned int nParticles = myGenParticles.size();
00646   nGenParticles = nParticles;
00647 
00648   if ( !nParticles ) return;
00649 
00651   std::map<const reco::Candidate*,int> myGenVertices;
00652 
00653   // Are there particles in the FSimEvent already ? 
00654   int offset = nTracks();
00655 
00656   // Skip the incoming protons
00657   nGenParticles = 0;
00658   unsigned int ip = 0;
00659   if ( nParticles > 1 && 
00660        myGenParticles[0].pdgId() == 2212 &&
00661        myGenParticles[1].pdgId() == 2212 ) { 
00662     ip = 2;
00663     nGenParticles = 2;
00664   }
00665 
00666   // Primary vertex (already smeared by the SmearedVtx module)
00667   XYZTLorentzVector primaryVertex (myGenParticles[ip].vx(),
00668                                    myGenParticles[ip].vy(),
00669                                    myGenParticles[ip].vz(),
00670                                    0.);
00671 
00672   // Smear the main vertex if needed
00673   XYZTLorentzVector smearedVertex;
00674   if ( primaryVertex.mag() < 1E-8 ) {
00675     theVertexGenerator->generate();
00676     smearedVertex = XYZTLorentzVector(
00677       theVertexGenerator->X()-theVertexGenerator->beamSpot().X()+theBeamSpot.X(),
00678       theVertexGenerator->Y()-theVertexGenerator->beamSpot().Y()+theBeamSpot.Y(),
00679       theVertexGenerator->Z()-theVertexGenerator->beamSpot().Z()+theBeamSpot.Z(),
00680       0.);
00681   }
00682 
00683   // Set the main vertex
00684   myFilter->setMainVertex(primaryVertex+smearedVertex);
00685 
00686   // This is the smeared main vertex
00687   int mainVertex = addSimVertex(myFilter->vertex(), -1, FSimVertexType::PRIMARY_VERTEX);
00688 
00689   // Loop on the particles of the generated event
00690   for ( ; ip<nParticles; ++ip ) { 
00691     
00692     // nGenParticles = ip;
00693     
00694     nGenParticles++;
00695     const reco::GenParticle& p = myGenParticles[ip];
00696 
00697     // Reject particles with late origin vertex (i.e., coming from late decays)
00698     // This should not happen, but one never knows what users may be up to!
00699     // For example exotic particles might decay late - keep the decay products in the case.
00700     XYZTLorentzVector productionVertexPosition(0.,0.,0.,0.);
00701     const reco::Candidate* productionMother = p.numberOfMothers() ? p.mother(0) : 0;
00702     if ( productionMother ) {
00703       unsigned motherId = productionMother->pdgId();
00704       if ( abs(motherId) < 1000000 )
00705         productionVertexPosition = XYZTLorentzVector(p.vx(), p.vy(), p.vz(), 0.) + smearedVertex;
00706     }
00707     if ( !myFilter->accept(productionVertexPosition) ) continue;
00708 
00709     // Keep only: 
00710     // 1) Stable particles
00711     bool testStable = p.status()%1000==1;
00712     // Declare stable standard particles that decay after a macroscopic path length 
00713     // (except if exotic particle)
00714     if ( p.status() == 2 && abs(p.pdgId()) < 1000000 ) {
00715       unsigned int nDaughters = p.numberOfDaughters();
00716       if ( nDaughters ) { 
00717         const reco::Candidate* daughter = p.daughter(0);
00718         XYZTLorentzVector decayPosition = 
00719           XYZTLorentzVector(daughter->vx(), daughter->vy(), daughter->vz(), 0.) + smearedVertex;
00720         // If the particle flew enough to be beyond the beam pipe enveloppe, just declare it stable
00721         if ( decayPosition.Perp2() > lateVertexPosition ) testStable = true;
00722       }
00723     }
00724 
00725     // 2) or particles with stable daughters
00726     bool testDaugh = false;
00727     unsigned int nDaughters = p.numberOfDaughters();
00728     if ( !testStable  && 
00729          //      p.status() == 2 && 
00730          nDaughters ) {  
00731       for ( unsigned iDaughter=0; iDaughter<nDaughters; ++iDaughter ) {
00732         const reco::Candidate* daughter = p.daughter(iDaughter);
00733         if ( daughter->status()%1000==1 ) {
00734           testDaugh=true;
00735           break;
00736         }
00737       }
00738     }
00739 
00740     // 3) or particles that fly more than one micron.
00741     double dist = 0.;
00742     if ( !testStable && !testDaugh ) {
00743       XYZTLorentzVector productionVertex(p.vx(),p.vy(),p.vz(),0.);
00744       dist = (primaryVertex-productionVertex).Vect().Mag2();
00745     }
00746     bool testDecay = ( dist > 1e-8 ) ? true : false; 
00747 
00748     // Save the corresponding particle and vertices
00749     if ( testStable || testDaugh || testDecay ) {
00750       
00751       const reco::Candidate* mother = p.numberOfMothers() ? p.mother(0) : 0;
00752 
00753       int originVertex = 
00754         mother &&  
00755         myGenVertices.find(mother) != myGenVertices.end() ? 
00756         myGenVertices[mother] : mainVertex;
00757       
00758       XYZTLorentzVector momentum(p.px(),p.py(),p.pz(),p.energy());
00759       RawParticle part(momentum, vertex(originVertex).position());
00760       part.setID(p.pdgId());
00761 
00762       // Add the particle to the event and to the various lists
00763       int theTrack = addSimTrack(&part,originVertex, nGenParts()-offset);
00764 
00765       // It there an end vertex ?
00766       if ( !nDaughters ) continue; 
00767       const reco::Candidate* daughter = p.daughter(0);
00768 
00769       // Add the vertex to the event and to the various lists
00770       XYZTLorentzVector decayVertex = 
00771         XYZTLorentzVector(daughter->vx(), daughter->vy(),
00772                           daughter->vz(), 0.) + smearedVertex;
00773       int theVertex = addSimVertex(decayVertex,theTrack, FSimVertexType::DECAY_VERTEX);
00774 
00775       if ( theVertex != -1 ) myGenVertices[&p] = theVertex;
00776 
00777       // There we are !
00778     }
00779   }
00780 
00781   // There is no GenParticle's in that case...
00782   // nGenParticles=0;
00783 
00784 }
00785 
00786 int 
00787 FBaseSimEvent::addSimTrack(const RawParticle* p, int iv, int ig, 
00788                            const HepMC::GenVertex* ev) { 
00789   
00790   // Check that the particle is in the Famos "acceptance"
00791   // Keep all primaries of pile-up events, though
00792   if ( !myFilter->accept(p) && ig >= -1 ) return -1;
00793 
00794   // The new track index
00795   int trackId = nSimTracks++;
00796   if ( nSimTracks/theTrackSize*theTrackSize == nSimTracks ) {
00797     theTrackSize *= 2;
00798     theSimTracks->resize(theTrackSize);
00799   }
00800 
00801   // Attach the particle to the origin vertex, and to the mother
00802   vertex(iv).addDaughter(trackId);
00803   if ( !vertex(iv).noParent() ) {
00804     track(vertex(iv).parent().id()).addDaughter(trackId);
00805 
00806     if ( ig == -1 ) {
00807       int motherId = track(vertex(iv).parent().id()).genpartIndex();
00808       if ( motherId < -1 ) ig = motherId;
00809     }
00810   }
00811     
00812   // Some transient information for FAMOS internal use
00813   (*theSimTracks)[trackId] = ev ? 
00814     // A proper decay time is scheduled
00815     FSimTrack(p,iv,ig,trackId,this,
00816               ev->position().t()/10.
00817               * p->PDGmass()
00818               / std::sqrt(p->momentum().Vect().Mag2())) : 
00819     // No proper decay time is scheduled
00820     FSimTrack(p,iv,ig,trackId,this);
00821 
00822   return trackId;
00823 
00824 }
00825 
00826 int
00827 FBaseSimEvent::addSimVertex(const XYZTLorentzVector& v, int im, FSimVertexType::VertexType type) {
00828   
00829   // Check that the vertex is in the Famos "acceptance"
00830   if ( !myFilter->accept(v) ) return -1;
00831 
00832   // The number of vertices
00833   int vertexId = nSimVertices++;
00834   if ( nSimVertices/theVertexSize*theVertexSize == nSimVertices ) {
00835     theVertexSize *= 2;
00836     theSimVertices->resize(theVertexSize);
00837     theFSimVerticesType->resize(theVertexSize);
00838   }
00839 
00840   // Attach the end vertex to the particle (if accepted)
00841   if ( im !=-1 ) track(im).setEndVertex(vertexId);
00842 
00843   // Some transient information for FAMOS internal use
00844   (*theSimVertices)[vertexId] = FSimVertex(v,im,vertexId,this);
00845 
00846   (*theFSimVerticesType)[vertexId] = FSimVertexType(type);
00847 
00848   return vertexId;
00849 
00850 }
00851 
00852 void
00853 FBaseSimEvent::printMCTruth(const HepMC::GenEvent& myGenEvent) {
00854   
00855   std::cout << "Id  Gen Name       eta    phi     pT     E    Vtx1   " 
00856             << " x      y      z   " 
00857             << "Moth  Vtx2  eta   phi     R      Z   Da1  Da2 Ecal?" << std::endl;
00858 
00859   for ( HepMC::GenEvent::particle_const_iterator 
00860           piter  = myGenEvent.particles_begin();
00861           piter != myGenEvent.particles_end(); 
00862         ++piter ) {
00863     
00864     HepMC::GenParticle* p = *piter;
00865      /* */
00866     int partId = p->pdg_id();
00867     std::string name;
00868 
00869     if ( pdt->particle(ParticleID(partId)) !=0 ) {
00870       name = (pdt->particle(ParticleID(partId)))->name();
00871     } else {
00872       name = "none";
00873     }
00874   
00875     XYZTLorentzVector momentum1(p->momentum().px(),
00876                                 p->momentum().py(),
00877                                 p->momentum().pz(),
00878                                 p->momentum().e());
00879 
00880     int vertexId1 = 0;
00881 
00882     if ( !p->production_vertex() ) continue;
00883 
00884     XYZVector vertex1 (p->production_vertex()->position().x()/10.,
00885                        p->production_vertex()->position().y()/10.,
00886                        p->production_vertex()->position().z()/10.);
00887     vertexId1 = p->production_vertex()->barcode();
00888     
00889     std::cout.setf(std::ios::fixed, std::ios::floatfield);
00890     std::cout.setf(std::ios::right, std::ios::adjustfield);
00891     
00892     std::cout << std::setw(4) << p->barcode() << " " 
00893          << name;
00894     
00895     for(unsigned int k=0;k<11-name.length() && k<12; k++) std::cout << " ";  
00896     
00897     double eta = momentum1.eta();
00898     if ( eta > +10. ) eta = +10.;
00899     if ( eta < -10. ) eta = -10.;
00900     std::cout << std::setw(6) << std::setprecision(2) << eta << " " 
00901               << std::setw(6) << std::setprecision(2) << momentum1.phi() << " " 
00902               << std::setw(7) << std::setprecision(2) << momentum1.pt() << " " 
00903               << std::setw(7) << std::setprecision(2) << momentum1.e() << " " 
00904               << std::setw(4) << vertexId1 << " " 
00905               << std::setw(6) << std::setprecision(1) << vertex1.x() << " " 
00906               << std::setw(6) << std::setprecision(1) << vertex1.y() << " " 
00907               << std::setw(6) << std::setprecision(1) << vertex1.z() << " ";
00908 
00909     const HepMC::GenParticle* mother = 
00910       *(p->production_vertex()->particles_in_const_begin());
00911 
00912     if ( mother )
00913       std::cout << std::setw(4) << mother->barcode() << " ";
00914     else 
00915       std::cout << "     " ;
00916     
00917     if ( p->end_vertex() ) {  
00918       XYZTLorentzVector vertex2(p->end_vertex()->position().x()/10.,
00919                                 p->end_vertex()->position().y()/10.,
00920                                 p->end_vertex()->position().z()/10.,
00921                                 p->end_vertex()->position().t()/10.);
00922       int vertexId2 = p->end_vertex()->barcode();
00923 
00924       std::vector<const HepMC::GenParticle*> children;
00925       HepMC::GenVertex::particles_out_const_iterator firstDaughterIt = 
00926         p->end_vertex()->particles_out_const_begin();
00927       HepMC::GenVertex::particles_out_const_iterator lastDaughterIt = 
00928         p->end_vertex()->particles_out_const_end();
00929       for ( ; firstDaughterIt != lastDaughterIt ; ++firstDaughterIt ) {
00930         children.push_back(*firstDaughterIt);
00931       }      
00932 
00933       std::cout << std::setw(4) << vertexId2 << " "
00934                 << std::setw(6) << std::setprecision(2) << vertex2.eta() << " " 
00935                 << std::setw(6) << std::setprecision(2) << vertex2.phi() << " " 
00936                 << std::setw(5) << std::setprecision(1) << vertex2.pt() << " " 
00937                 << std::setw(6) << std::setprecision(1) << vertex2.z() << " ";
00938       for ( unsigned id=0; id<children.size(); ++id )
00939         std::cout << std::setw(4) << children[id]->barcode() << " ";
00940     }
00941     std::cout << std::endl;
00942 
00943   }
00944 
00945 }
00946 
00947 void
00948 FBaseSimEvent::print() const {
00949 
00950   std::cout << "  Id  Gen Name       eta    phi     pT     E    Vtx1   " 
00951             << " x      y      z   " 
00952             << "Moth  Vtx2  eta   phi     R      Z   Daughters Ecal?" << std::endl;
00953 
00954   for( int i=0; i<(int)nTracks(); i++ ) 
00955     std::cout << track(i) << std::endl;
00956 
00957   for( int i=0; i<(int)nVertices(); i++ ) 
00958     std::cout << "i = " << i << "  " << vertexType(i) << std::endl;
00959 
00960   
00961 
00962 }
00963 
00964 void 
00965 FBaseSimEvent::clear() {
00966 
00967   nSimTracks = 0;
00968   nSimVertices = 0;
00969   nGenParticles = 0;
00970   nChargedParticleTracks = 0;
00971 
00972 }
00973 
00974 void 
00975 FBaseSimEvent::addChargedTrack(int id) { 
00976   (*theChargedTracks)[nChargedParticleTracks++] = id;
00977   if ( nChargedParticleTracks/theChargedSize*theChargedSize 
00978        == nChargedParticleTracks ) {
00979     theChargedSize *= 2;
00980     theChargedTracks->resize(theChargedSize);
00981   }
00982 }
00983 
00984 int
00985 FBaseSimEvent::chargedTrack(int id) const {
00986   if (id>=0 && id<(int)nChargedParticleTracks) 
00987     return (*theChargedTracks)[id]; 
00988   else 
00989     return -1;
00990 }
00991 
00992 /* 
00993 const SimTrack & 
00994 FBaseSimEvent::embdTrack(int i) const {  
00995   return (*theSimTracks)[i].simTrack();
00996 }
00997 
00998 const SimVertex & 
00999 FBaseSimEvent::embdVertex(int i) const { 
01000   return (*theSimVertices)[i].simVertex();
01001 }
01002 */
01003 
01004 const HepMC::GenParticle* 
01005 FBaseSimEvent::embdGenpart(int i) const {
01006   return (*theGenParticles)[i]; 
01007 }
01008 
01009 /*
01010 FSimTrack&  
01011 FBaseSimEvent::track(int id) const { 
01012   return (*theSimTracks)[id];
01013 }
01014 
01015 
01016 FSimVertex&  
01017 FBaseSimEvent::vertex(int id) const { 
01018   return (*theSimVertices)[id];
01019 }
01020 */