#include "TPythia8.h"
#include "TClonesArray.h"
#include "TParticle.h"
#include "TLorentzVector.h"
ClassImp(TPythia8)
TPythia8*  TPythia8::fgInstance = 0;
TPythia8::TPythia8():
    TGenerator("TPythia8", "TPythia8"),
    fPythia(0),
    fNumberOfParticles(0)
{
   
   if (fgInstance) 
      Fatal("TPythia8", "There's already an instance of TPythia8");
  
   delete fParticles; 
    
   fParticles = new TClonesArray("TParticle",50);
   fPythia    = new Pythia8::Pythia();
}
TPythia8::TPythia8(const char *xmlDir):
    TGenerator("TPythia8", "TPythia8"),
    fPythia(0),
    fNumberOfParticles(0)
{
   
   if (fgInstance) 
      Fatal("TPythia8", "There's already an instance of TPythia8");
  
   delete fParticles; 
    
   fParticles = new TClonesArray("TParticle",50);
   fPythia    = new Pythia8::Pythia(xmlDir);
}
TPythia8::~TPythia8()
{
   
   if (fParticles) {
      fParticles->Delete();
      delete fParticles;
      fParticles = 0;
   }
   delete fPythia;
}
TPythia8* TPythia8::Instance()
{
   
   return fgInstance ? fgInstance : (fgInstance = new TPythia8()) ;
}
Bool_t TPythia8::Initialize(Int_t idAin, Int_t idBin, Double_t ecms)
{
   
   AddParticlesToPdgDataBase();
   return fPythia->init(idAin, idBin, ecms);
}
void TPythia8::GenerateEvent()
{
   
   fPythia->next();
   fNumberOfParticles  = fPythia->event.size() - 1;
   ImportParticles();
}
Int_t TPythia8::ImportParticles(TClonesArray *particles, Option_t *option)
{
   
   if (particles == 0) return 0;
   TClonesArray &clonesParticles = *particles;
   clonesParticles.Clear();
   Int_t nparts=0;
   Int_t i;
   fNumberOfParticles  = fPythia->event.size() - 1;
   if (!strcmp(option,"") || !strcmp(option,"Final")) {
      for (i = 0; i <= fNumberOfParticles; i++) {
	if (fPythia->event[i].id() == 90) continue;
         if (fPythia->event[i].isFinal()) {
            new(clonesParticles[nparts]) TParticle(
                fPythia->event[i].id(),
                fPythia->event[i].isFinal(),
                fPythia->event[i].mother1() - 1,
                fPythia->event[i].mother2() - 1,
                fPythia->event[i].daughter1() - 1, 
                fPythia->event[i].daughter2() - 1,
                fPythia->event[i].px(),     
                fPythia->event[i].py(),     
                fPythia->event[i].pz(),     
                fPythia->event[i].e(),      
                fPythia->event[i].xProd(),  
                fPythia->event[i].yProd(),  
                fPythia->event[i].zProd(),  
                fPythia->event[i].tProd()); 
		nparts++;
	    } 
	} 
    } else if (!strcmp(option,"All")) {
	for (i = 0; i <= fNumberOfParticles; i++) {
	  if (fPythia->event[i].id() == 90) continue;
	    new(clonesParticles[nparts]) TParticle(
		fPythia->event[i].id(),
		fPythia->event[i].isFinal(),
		fPythia->event[i].mother1() - 1,
		fPythia->event[i].mother2() - 1,
		fPythia->event[i].daughter1() - 1,
		fPythia->event[i].daughter2() - 1,
		fPythia->event[i].px(),       
		fPythia->event[i].py(),       
		fPythia->event[i].pz(),       
		fPythia->event[i].e(),        
		fPythia->event[i].xProd(),    
		fPythia->event[i].yProd(),    
		fPythia->event[i].zProd(),    
		fPythia->event[i].tProd());   
	    nparts++;
	} 
    }
    return nparts;
}
TObjArray* TPythia8::ImportParticles(Option_t* )
{
   
   fParticles->Clear();
   Int_t numpart   = fPythia->event.size() - 1;
   TClonesArray &a = *((TClonesArray*)fParticles);
   for (Int_t i = 1; i <= numpart; i++) {
      new(a[i]) TParticle(
         fPythia->event[i].id(),
         fPythia->event[i].isFinal(),
         fPythia->event[i].mother1()  - 1,
         fPythia->event[i].mother2()  - 1,
         fPythia->event[i].daughter1() - 1,
         fPythia->event[i].daughter2() - 1,
         fPythia->event[i].px(),       
         fPythia->event[i].py(),       
         fPythia->event[i].pz(),       
         fPythia->event[i].e(),        
         fPythia->event[i].xProd(),    
         fPythia->event[i].yProd(),    
         fPythia->event[i].zProd(),    
         fPythia->event[i].tProd());   
   }
   return fParticles;
}
Int_t TPythia8::GetN() const
{
   
   return (fPythia->event.size() - 1);
}
void TPythia8::ReadString(const char* string) const
{
   
   fPythia->readString(string);
}
void  TPythia8::ReadConfigFile(const char* string) const
{
  
  fPythia->readFile(string);
}
void TPythia8::PrintStatistics() const
{
   
   fPythia->statistics();
}
void TPythia8::EventListing() const
{
   
   fPythia->event.list();
}
void TPythia8::AddParticlesToPdgDataBase()
{
   
   TDatabasePDG *pdgDB = TDatabasePDG::Instance();
   pdgDB->AddParticle("string","string", 0, kTRUE,
                      0, 0, "QCD string", 90);
   pdgDB->AddParticle("rho_diff0", "rho_diff0", 0, kTRUE,
                      0, 0, "QCD diffr. state", 9900110);
   pdgDB->AddParticle("pi_diffr+", "pi_diffr+", 0, kTRUE,
                      0, 1, "QCD diffr. state", 9900210);
   pdgDB->AddParticle("omega_di", "omega_di", 0, kTRUE,
                      0, 0, "QCD diffr. state", 9900220);
   pdgDB->AddParticle("phi_diff","phi_diff", 0, kTRUE,
                      0, 0, "QCD diffr. state", 9900330);
   pdgDB->AddParticle("J/psi_di", "J/psi_di", 0, kTRUE,
                      0, 0, "QCD diffr. state", 9900440);
   pdgDB->AddParticle("n_diffr0","n_diffr0",0,kTRUE,
                      0, 0, "QCD diffr. state", 9902110);
   pdgDB->AddParticle("p_diffr+","p_diffr+", 0, kTRUE,
                      0, 1, "QCD diffr. state", 9902210);
}
Last change: Mon Sep 22 15:42:19 2008
Last generated: 2008-09-22 15:42
This page has been automatically generated. If you have any comments or suggestions about the page layout send a mail to ROOT support, or contact the developers with any questions or problems regarding ROOT.