]> git.uio.no Git - u/mrichter/AliRoot.git/blob - MUON/gen.C
Update From Debojit
[u/mrichter/AliRoot.git] / MUON / gen.C
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 /* $Id$ */
17
18 /// \ingroup macros
19 /// \file gen.C
20 /// \brief Configuration macro for generation of kinematics tree which 
21 /// can be then used as an external event generator.
22 ///
23 /// According to: $ALICE_ROOT/test/genkine/gen/fastGen.C
24
25 #if !defined(__CINT__) || defined(__MAKECINT__)
26 #include <Riostream.h>
27 #include <TH1F.h>
28 #include <TStopwatch.h>
29 #include <TDatime.h>
30 #include <TRandom.h>
31 #include <TDatabasePDG.h>
32 #include <TParticle.h>
33 #include <TArrayI.h>
34
35 #include "AliGenerator.h"
36 #include "AliPDG.h"
37 #include "AliRunLoader.h"
38 #include "AliRun.h"
39 #include "AliStack.h"
40 #include "AliHeader.h"
41 #include "PYTHIA6/AliGenPythia.h"
42 #include "PYTHIA6/AliPythia.h"
43 #endif
44
45 void gen(Int_t nev = 1, 
46          const char* genConfig = "$ALICE_ROOT/MUON/genTestConfig.C")
47 {
48   // Load libraries
49   // gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include -I$ALICE_ROOT");
50   gSystem->Load("liblhapdf.so");      // Parton density functions
51   gSystem->Load("libEGPythia6.so");   // TGenerator interface
52   gSystem->Load("libpythia6.so");     // Pythia
53   gSystem->Load("libAliPythia6.so");  // ALICE specific implementations
54
55   AliPDG::AddParticlesToPdgDataBase();
56   TDatabasePDG::Instance();
57
58   // Run loader
59   AliRunLoader* rl = AliRunLoader::Open("galice.root","FASTRUN","recreate");
60   
61   rl->SetCompressionLevel(2);
62   rl->SetNumberOfEventsPerFile(nev);
63   rl->LoadKinematics("RECREATE");
64   rl->MakeTree("E");
65   gAlice->SetRunLoader(rl);
66   
67   //  Create stack
68   rl->MakeStack();
69   AliStack* stack = rl->Stack();
70   
71   //  Header
72   AliHeader* header = rl->GetHeader();
73   
74   //  Create and Initialize Generator
75   gROOT->LoadMacro(genConfig);
76   AliGenerator* gener = genConfig();
77
78   // Go to galice.root
79   rl->CdGAFile();
80
81   // Forbid some decays. Do it after gener->Init(0, because
82   // the initialization of the generator includes reading of the decay table.
83   // ...
84
85   //
86   // Event Loop
87   //
88   
89   TStopwatch timer;
90   timer.Start();
91   for (Int_t iev = 0; iev < nev; iev++) {
92     
93     cout <<"Event number "<< iev << endl;
94     
95     // Initialize event
96     header->Reset(0,iev);
97     rl->SetEventNumber(iev);
98     stack->Reset();
99     rl->MakeTree("K");
100     
101     // Generate event
102     stack->Reset();
103     stack->ConnectTree(rl->TreeK());
104     gener->Generate();
105     cout << "Number of particles " << stack->GetNprimary() << endl;
106     
107     // Finish event
108     header->SetNprimary(stack->GetNprimary());
109     header->SetNtrack(stack->GetNtrack());  
110     
111     // I/O
112     stack->FinishEvent();
113     header->SetStack(stack);
114     rl->TreeE()->Fill();
115     rl->WriteKinematics("OVERWRITE");
116     
117   } // event loop
118   timer.Stop();
119   timer.Print();
120   
121   //                         Termination
122   //  Generator
123   gener->FinishRun();
124   //  Write file
125   rl->WriteHeader("OVERWRITE");
126   gener->Write();
127   rl->Write();
128 }