CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/FastSimulation/ParticleDecay/src/PythiaDecays.cc

Go to the documentation of this file.
00001 // HepMC Headers
00002 #include "HepMC/PythiaWrapper6_4.h"
00003 
00004 // Pythia6 framework integration service Headers
00005 #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Service.h"
00006 
00007 // FAMOS Headers
00008 #include "FastSimulation/ParticlePropagator/interface/ParticlePropagator.h"
00009 #include "FastSimulation/ParticleDecay/interface/PythiaDecays.h"
00010 #include "FastSimulation/ParticleDecay/interface/Pythia6jets.h"
00011 #include "FastSimulation/ParticleDecay/interface/RandomP8.h"
00012 
00013 
00014 // Needed for Pythia6 
00015 #define PYTHIA6PYDECY pythia6pydecy_
00016 
00017 extern "C" {
00018   void PYTHIA6PYDECY(int *ip);
00019 }
00020 
00021 PythiaDecays::PythiaDecays(std::string program)
00022 {
00023   program_=program;
00024   if (program_ == "pythia6") {
00026     pyjets = new Pythia6jets();
00027     pyservice = new gen::Pythia6Service();
00028     // The PYTHIA decay tables will be initialized later 
00029   } else if (program_ == "pythia8") {
00031     pythia.reset(new Pythia8::Pythia);
00032     decayer.reset(new Pythia8::Pythia);
00033 
00034     RandomP8* RP8 = new RandomP8();
00035     pythia->setRndmEnginePtr(RP8);
00036     decayer->setRndmEnginePtr(RP8);
00037   } else {
00038     std::cout << "WARNING: you are requesting an option which is not available in PythiaDecays::PythiaDecays " << std::endl;
00039   }
00040 
00041 }
00042 
00043 PythiaDecays::~PythiaDecays() {
00044   if (program_ == "pythia6") {
00045     delete pyjets;
00046     delete pyservice;
00047   }
00048 }
00049 
00050 const DaughterParticleList&
00051 PythiaDecays::particleDaughtersPy8(ParticlePropagator& particle)
00052 {
00053 
00054   theList.clear();
00055 
00056   /*
00057   std::cout << "###" << std::endl;
00058   std::cout << "particle.PDGname() == " << particle.PDGname() << std::endl; 
00059   std::cout << "particle.status() == " << particle.status() << std::endl;
00060   std::cout << "particle.pid() == " << particle.pid() << std::endl;
00061   std::cout << "particle.momentum().x() == " << particle.momentum().x() << std::endl;
00062   std::cout << "particle.Px() == " << particle.Px() << std::endl;
00063   std::cout << "particle.X() == " << particle.X() << std::endl;
00064   std::cout << "particle.mass() == " << particle.mass() << std::endl;
00065   std::cout << "particle.PDGmass() == " << particle.PDGmass() << std::endl;
00066   std::cout << "particle.PDGcTau() == " << particle.PDGcTau() << std::endl;
00067   std::cout << "(decayer->particleData).tau0( particle.pid() ) == " << (decayer->particleData).tau0( particle.pid() ) << std::endl;
00068   */
00069 
00070   // inspired by method Pythia8Hadronizer::residualDecay() in GeneratorInterface/Pythia8Interface/src/Pythia8Hadronizer.cc
00071   decayer->event.reset();
00072   Pythia8::Particle py8part(  particle.pid(), 93, 0, 0, 0, 0, 0, 0,
00073                      particle.momentum().x(), // note: momentum().x() and Px() are the same
00074                      particle.momentum().y(),
00075                      particle.momentum().z(),
00076                      particle.momentum().t(),
00077                      particle.mass() );
00078   py8part.vProd( particle.X(), particle.Y(), 
00079                  particle.Z(), particle.T() );
00080   py8part.tau( (decayer->particleData).tau0( particle.pid() ) );
00081   decayer->event.append( py8part );
00082 
00083   int nentries = decayer->event.size();
00084   if ( !decayer->event[nentries-1].mayDecay() ) return theList;
00085   decayer->next();
00086   int nentries1 = decayer->event.size();
00087   if ( nentries1 <= nentries ) return theList; //same number of particles, no decays...
00088   Pythia8::Particle& py8daughter = decayer->event[nentries]; // the 1st daughter // DO I NEED THIS LINE?
00089 
00090   theList.resize(nentries,RawParticle());
00091 
00092   for ( int ipart=nentries+1; ipart<nentries1; ipart++ )
00093     {
00094       py8daughter = decayer->event[ipart];
00095       theList[ipart-nentries-1].SetXYZT( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
00096       theList[ipart-nentries-1].setVertex( py8daughter.xProd(),
00097                                            py8daughter.yProd(),
00098                                            py8daughter.zProd(),
00099                                            py8daughter.tProd() );
00100       theList[ipart-nentries-1].setID( py8daughter.id() );
00101       theList[ipart-nentries-1].setMass( py8daughter.m() );
00102     }
00103 
00104   return theList;
00105 }
00106 
00107 const DaughterParticleList&
00108 PythiaDecays::particleDaughtersPy6(ParticlePropagator& particle)
00109 {
00110   gen::Pythia6Service::InstanceWrapper guard(pyservice); // grab Py6 context
00111 
00112   //  Pythia6jets pyjets;
00113   int ip;
00114 
00115   pyjets->k(1,1) = 1;
00116   pyjets->k(1,2) = particle.pid();
00117   pyjets->p(1,1) = particle.Px();
00118   pyjets->p(1,2) = particle.Py();
00119   pyjets->p(1,3) = particle.Pz();
00120   pyjets->p(1,4) = std::max(particle.mass(),particle.e());
00121   pyjets->p(1,5) = particle.mass();
00122   pyjets->v(1,1) = particle.X();
00123   pyjets->v(1,2) = particle.Y();
00124   pyjets->v(1,3) = particle.Z();
00125   pyjets->v(1,4) = particle.T();
00126   pyjets->n() = 1;
00127   
00128   ip = 1;
00129   PYTHIA6PYDECY(&ip);
00130 
00131   // Fill the list of daughters
00132   theList.clear();
00133   if ( pyjets->n()==1 ) return theList; 
00134 
00135   theList.resize(pyjets->n()-1,RawParticle());
00136 
00137   for (int i=2;i<=pyjets->n();++i) {
00138     
00139     theList[i-2].SetXYZT(pyjets->p(i,1),pyjets->p(i,2),
00140                          pyjets->p(i,3),pyjets->p(i,4)); 
00141     theList[i-2].setVertex(pyjets->v(i,1),pyjets->v(i,2),
00142                            pyjets->v(i,3),pyjets->v(i,4));
00143     theList[i-2].setID(pyjets->k(i,2));
00144     theList[i-2].setMass(pyjets->p(i,5));
00145 
00146   }
00147 
00148   return theList;
00149   
00150 }