00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051
00052
00053
00054
00055
00056
00057
00058
00059
00060
00061
00062
00063
00064
00065
00066
00067
00068
00069
00070
00071
00072
00073
00074
00075 #include "TPythia8.h"
00076
00077 #include "TClonesArray.h"
00078 #include "TParticle.h"
00079 #include "TDatabasePDG.h"
00080 #include "TLorentzVector.h"
00081
00082 ClassImp(TPythia8)
00083
00084 TPythia8* TPythia8::fgInstance = 0;
00085
00086
00087 TPythia8::TPythia8():
00088 TGenerator("TPythia8", "TPythia8"),
00089 fPythia(0),
00090 fNumberOfParticles(0)
00091 {
00092
00093 if (fgInstance)
00094 Fatal("TPythia8", "There's already an instance of TPythia8");
00095
00096 delete fParticles;
00097
00098 fParticles = new TClonesArray("TParticle",50);
00099 fPythia = new Pythia8::Pythia();
00100 }
00101
00102
00103 TPythia8::TPythia8(const char *xmlDir):
00104 TGenerator("TPythia8", "TPythia8"),
00105 fPythia(0),
00106 fNumberOfParticles(0)
00107 {
00108
00109 if (fgInstance)
00110 Fatal("TPythia8", "There's already an instance of TPythia8");
00111
00112 delete fParticles;
00113
00114 fParticles = new TClonesArray("TParticle",50);
00115 fPythia = new Pythia8::Pythia(xmlDir);
00116 }
00117
00118
00119 TPythia8::~TPythia8()
00120 {
00121
00122 if (fParticles) {
00123 fParticles->Delete();
00124 delete fParticles;
00125 fParticles = 0;
00126 }
00127 delete fPythia;
00128 }
00129
00130
00131 TPythia8* TPythia8::Instance()
00132 {
00133
00134 return fgInstance ? fgInstance : (fgInstance = new TPythia8()) ;
00135 }
00136
00137
00138 Bool_t TPythia8::Initialize(Int_t idAin, Int_t idBin, Double_t ecms)
00139 {
00140
00141 AddParticlesToPdgDataBase();
00142 return fPythia->init(idAin, idBin, ecms);
00143 }
00144
00145
00146 void TPythia8::GenerateEvent()
00147 {
00148
00149 fPythia->next();
00150 fNumberOfParticles = fPythia->event.size() - 1;
00151 ImportParticles();
00152 }
00153
00154
00155 Int_t TPythia8::ImportParticles(TClonesArray *particles, Option_t *option)
00156 {
00157
00158 if (particles == 0) return 0;
00159 TClonesArray &clonesParticles = *particles;
00160 clonesParticles.Clear();
00161 Int_t nparts=0;
00162 Int_t i;
00163 fNumberOfParticles = fPythia->event.size() - 1;
00164
00165 if (!strcmp(option,"") || !strcmp(option,"Final")) {
00166 for (i = 0; i <= fNumberOfParticles; i++) {
00167 if (fPythia->event[i].id() == 90) continue;
00168 if (fPythia->event[i].isFinal()) {
00169 new(clonesParticles[nparts]) TParticle(
00170 fPythia->event[i].id(),
00171 fPythia->event[i].status(),
00172 fPythia->event[i].mother1() - 1,
00173 fPythia->event[i].mother2() - 1,
00174 fPythia->event[i].daughter1() - 1,
00175 fPythia->event[i].daughter2() - 1,
00176 fPythia->event[i].px(),
00177 fPythia->event[i].py(),
00178 fPythia->event[i].pz(),
00179 fPythia->event[i].e(),
00180 fPythia->event[i].xProd(),
00181 fPythia->event[i].yProd(),
00182 fPythia->event[i].zProd(),
00183 fPythia->event[i].tProd());
00184 nparts++;
00185 }
00186 }
00187 } else if (!strcmp(option,"All")) {
00188 for (i = 0; i <= fNumberOfParticles; i++) {
00189 if (fPythia->event[i].id() == 90) continue;
00190 new(clonesParticles[nparts]) TParticle(
00191 fPythia->event[i].id(),
00192 fPythia->event[i].status(),
00193 fPythia->event[i].mother1() - 1,
00194 fPythia->event[i].mother2() - 1,
00195 fPythia->event[i].daughter1() - 1,
00196 fPythia->event[i].daughter2() - 1,
00197 fPythia->event[i].px(),
00198 fPythia->event[i].py(),
00199 fPythia->event[i].pz(),
00200 fPythia->event[i].e(),
00201 fPythia->event[i].xProd(),
00202 fPythia->event[i].yProd(),
00203 fPythia->event[i].zProd(),
00204 fPythia->event[i].tProd());
00205 nparts++;
00206 }
00207 }
00208 return nparts;
00209 }
00210
00211
00212 TObjArray* TPythia8::ImportParticles(Option_t* )
00213 {
00214
00215 fParticles->Clear();
00216 Int_t numpart = fPythia->event.size() - 1;
00217 TClonesArray &a = *((TClonesArray*)fParticles);
00218 for (Int_t i = 1; i <= numpart; i++) {
00219 new(a[i]) TParticle(
00220 fPythia->event[i].id(),
00221 fPythia->event[i].status(),
00222 fPythia->event[i].mother1() - 1,
00223 fPythia->event[i].mother2() - 1,
00224 fPythia->event[i].daughter1() - 1,
00225 fPythia->event[i].daughter2() - 1,
00226 fPythia->event[i].px(),
00227 fPythia->event[i].py(),
00228 fPythia->event[i].pz(),
00229 fPythia->event[i].e(),
00230 fPythia->event[i].xProd(),
00231 fPythia->event[i].yProd(),
00232 fPythia->event[i].zProd(),
00233 fPythia->event[i].tProd());
00234 }
00235 return fParticles;
00236 }
00237
00238
00239 Int_t TPythia8::GetN() const
00240 {
00241
00242 return (fPythia->event.size() - 1);
00243 }
00244
00245
00246 void TPythia8::ReadString(const char* string) const
00247 {
00248
00249 fPythia->readString(string);
00250 }
00251
00252
00253 void TPythia8::ReadConfigFile(const char* string) const
00254 {
00255
00256 fPythia->readFile(string);
00257 }
00258
00259
00260 void TPythia8::PrintStatistics() const
00261 {
00262
00263 fPythia->statistics();
00264 }
00265
00266
00267 void TPythia8::EventListing() const
00268 {
00269
00270 fPythia->event.list();
00271 }
00272
00273
00274 void TPythia8::AddParticlesToPdgDataBase()
00275 {
00276
00277
00278 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
00279 pdgDB->AddParticle("string","string", 0, kTRUE,
00280 0, 0, "QCD string", 90);
00281 pdgDB->AddParticle("rho_diff0", "rho_diff0", 0, kTRUE,
00282 0, 0, "QCD diffr. state", 9900110);
00283 pdgDB->AddParticle("pi_diffr+", "pi_diffr+", 0, kTRUE,
00284 0, 1, "QCD diffr. state", 9900210);
00285 pdgDB->AddParticle("omega_di", "omega_di", 0, kTRUE,
00286 0, 0, "QCD diffr. state", 9900220);
00287 pdgDB->AddParticle("phi_diff","phi_diff", 0, kTRUE,
00288 0, 0, "QCD diffr. state", 9900330);
00289 pdgDB->AddParticle("J/psi_di", "J/psi_di", 0, kTRUE,
00290 0, 0, "QCD diffr. state", 9900440);
00291 pdgDB->AddParticle("n_diffr0","n_diffr0",0,kTRUE,
00292 0, 0, "QCD diffr. state", 9902110);
00293 pdgDB->AddParticle("p_diffr+","p_diffr+", 0, kTRUE,
00294 0, 1, "QCD diffr. state", 9902210);
00295 }
00296