]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FASTSIM/fastGen.C
e68598298892fa4a98a17507d3a46bc79a329437
[u/mrichter/AliRoot.git] / FASTSIM / fastGen.C
1 AliGenerator*  CreateGenerator();
2
3 void fastGen(Int_t nev = 1, char* filename = "galice.root")
4 {
5 //
6 //                        Construction
7 //
8 //  Output file
9     TFile*  file         = new TFile(filename, "recreate");
10 //  Create stack
11     AliStack* stack      = new AliStack(10000);
12     stack->MakeTree(0, filename);
13
14 //  Create Header
15     AliHeader* header    = new AliHeader();
16 //  Create Header Tree
17     TTree* treeE         = new TTree("TE","Headers");
18     treeE->Branch("Header", "AliHeader", &header, 4000, 0);
19     treeE->Write();
20 //
21 //  Create and Initialize Generator
22     AliGenerator *gener = CreateGenerator();
23     gener->Init();
24     gener->SetStack(stack);
25     
26 //
27 //                        Event Loop
28 //
29     Int_t iev;
30      
31     for (iev = 0; iev < nev; iev++) {
32
33         printf("\n \n Event number %d \n \n", iev);
34         
35 //  Initialize event
36         header->Reset(0,iev);
37         stack->Reset();
38         stack->BeginEvent(iev);
39
40 //  Generate event
41         gener->Generate();
42 //  Analysis
43         Int_t npart = stack->GetNprimary();
44         printf("Analyse %d Particles\n", npart);
45         for (Int_t part=0; part<npart; part++) {
46             TParticle *MPart = stack->Particle(part);
47             Int_t mpart  = MPart->GetPdgCode();
48             printf("Particle %d\n", mpart);
49         }
50         
51 //  Finish event
52         header->SetNprimary(stack->GetNprimary());
53         header->SetNtrack(stack->GetNtrack());  
54 //      I/O
55 //      
56         stack->FinishEvent();
57         header->SetStack(stack);
58         treeE->Fill();
59         (stack->TreeK())->Write(0,TObject::kOverwrite);
60     } // event loop
61 //
62 //                         Termination
63 //  Generator
64     gener->FinishRun();
65 //  Header
66     treeE->Write(0,TObject::kOverwrite);
67     delete treeE;   treeE = 0;
68 //  Stack
69     stack->FinishRun();
70 //  Write file
71     gener->Write();
72     file->Write();
73 }
74
75
76 AliGenerator*  CreateGenerator()
77 {
78     gener = new AliGenPythia(1);
79 //
80 //
81 //   vertex position and smearing 
82     gener->SetVertexSmear(kPerEvent);
83 //   structure function
84     gener->SetStrucFunc(kGRVHO);
85 //   charm, beauty, charm_unforced, beauty_unforced, jpsi, jpsi_chi, mb
86     gener->SetProcess(kPyJets);
87 //   Centre of mass energy 
88     gener->SetEnergyCMS(5500.);
89 //   Pt transfer of the hard scattering
90     gener->SetPtHard(50.,50.2);
91 //   Initialize generator    
92     return gener;
93 }
94
95
96
97
98
99
100
101
102
103
104
105