]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/Config.C
Do not get main runloader from gAlice
[u/mrichter/AliRoot.git] / macros / Config.C
1 // One can use the configuration macro in compiled mode by
2 // root [0] gSystem->Load("libgeant321");
3 // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
4 //                   -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
5 // root [0] .x grun.C(1,"Config.C++")
6
7 #if !defined(__CINT__) || defined(__MAKECINT__)
8 #include <Riostream.h>
9 #include <TPDGCode.h>
10 #include <TRandom.h>
11 #include <TSystem.h>
12 #include <TVirtualMC.h>
13 #include <TGeant3TGeo.h>
14 #include "STEER/AliRunLoader.h"
15 #include "STEER/AliRun.h"
16 #include "STEER/AliConfig.h"
17 #include "PYTHIA6/AliDecayerPythia.h"
18 #include "EVGEN/AliGenCocktail.h"
19 #include "EVGEN/AliGenHIJINGpara.h"
20 #include "STEER/AliMagFMaps.h"
21 #include "STRUCT/AliBODY.h"
22 #include "STRUCT/AliMAG.h"
23 #include "STRUCT/AliABSOv3.h"
24 #include "STRUCT/AliDIPOv3.h"
25 #include "STRUCT/AliHALLv3.h"
26 #include "STRUCT/AliFRAMEv2.h"
27 #include "STRUCT/AliSHILv3.h"
28 #include "STRUCT/AliPIPEv3.h"
29 #include "ITS/AliITSv11Hybrid.h"
30 #include "TPC/AliTPCv2.h"
31 #include "TOF/AliTOFv6T0.h"
32 #include "HMPID/AliHMPIDv3.h"
33 #include "ZDC/AliZDCv3.h"
34 #include "TRD/AliTRDv1.h"
35 #include "FMD/AliFMDv1.h"
36 #include "MUON/AliMUONv1.h"
37 #include "PHOS/AliPHOSv1.h"
38 #include "PMD/AliPMDv1.h"
39 #include "T0/AliT0v1.h"
40 #include "EMCAL/AliEMCALv2.h"
41 #include "ACORDE/AliACORDEv1.h"
42 #include "VZERO/AliVZEROv7.h"
43 #endif
44
45 Float_t EtaToTheta(Float_t arg);
46 void    LoadPythia();
47
48
49 void Config()
50 {
51     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
52     // Theta range given through pseudorapidity limits 22/6/2001
53
54     // Set Random Number seed
55   //gRandom->SetSeed(123456); // Set 0 to use the current time
56   
57   AliLog::Message(AliLog::kInfo, Form("Seed for random number generation = %d",gRandom->GetSeed()), "Config.C", "Config.C", "Config()","Config.C", __LINE__);
58
59   // Load Pythia libraries
60   LoadPythia();
61   // Libraries required by geant321
62 #if defined(__CINT__)
63     gSystem->Load("libgeant321");
64 #endif
65
66     new     TGeant3TGeo("C++ Interface to Geant3");
67
68     AliRunLoader* rl=0x0;
69
70     AliLog::Message(AliLog::kInfo, "Creating Run Loader", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
71
72     rl = AliRunLoader::Open("galice.root",
73                             AliConfig::GetDefaultEventFolderName(),
74                             "recreate");
75     if (rl == 0x0)
76       {
77         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
78         return;
79       }
80     rl->SetCompressionLevel(2);
81     rl->SetNumberOfEventsPerFile(3);
82     gAlice->SetRunLoader(rl);
83     
84     // gAlice->SetGeometryFromFile("geometry.root");
85
86     // Uncomment if you want to load geometry from OCDB!   >>>>
87 /*    
88     if(!AliCDBManager::Instance()->IsDefaultStorageSet()){
89          cout << "#####################################################" << endl;
90          cout << "#                                                   #" << endl;
91          cout << "#     WARNING: CDB DEFAULT STORAGE NOT SET !!!      #" << endl;
92          cout << "#     SETTING IT TO local://$ALICE_ROOT !!!         #" << endl;
93          cout << "#                                                   #" << endl;
94          cout << "#####################################################" << endl;
95           
96          AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT");
97     }
98     
99     if(AliCDBManager::Instance()->GetRun() < 0){
100          cout << "#####################################################" << endl;
101          cout << "#                                                   #" << endl;
102          cout << "#     WARNING: RUN NUMBER NOT SET !!!               #" << endl;
103          cout << "#     SETTING IT TO 0 !!!                           #" << endl;
104          cout << "#                                                   #" << endl;
105          cout << "#####################################################" << endl;
106           
107          AliCDBManager::Instance()->SetRun(0);
108     }
109     gAlice->SetGeometryFromCDB();
110 */
111     // Uncomment if you want to load geometry from OCDB!   <<<<
112
113     // Set the trigger configuration
114     gAlice->SetTriggerDescriptor("Pb-Pb");
115     cout<<"Trigger configuration is set to  Pb-Pb"<<endl;
116
117     //
118     // Set External decayer
119     TVirtualMCDecayer *decayer = new AliDecayerPythia();
120
121     decayer->SetForceDecay(kAll);
122     decayer->Init();
123     gMC->SetExternalDecayer(decayer);
124     //=======================================================================
125     // ************* STEERING parameters FOR ALICE SIMULATION **************
126     // --- Specify event type to be tracked through the ALICE setup
127     // --- All positions are in cm, angles in degrees, and P and E in GeV
128
129
130     gMC->SetProcess("DCAY",1);
131     gMC->SetProcess("PAIR",1);
132     gMC->SetProcess("COMP",1);
133     gMC->SetProcess("PHOT",1);
134     gMC->SetProcess("PFIS",0);
135     gMC->SetProcess("DRAY",0);
136     gMC->SetProcess("ANNI",1);
137     gMC->SetProcess("BREM",1);
138     gMC->SetProcess("MUNU",1);
139     gMC->SetProcess("CKOV",1);
140     gMC->SetProcess("HADR",1);
141     gMC->SetProcess("LOSS",2);
142     gMC->SetProcess("MULS",1);
143     gMC->SetProcess("RAYL",1);
144
145     Float_t cut = 1.e-3;        // 1MeV cut by default
146     Float_t tofmax = 1.e10;
147
148     gMC->SetCut("CUTGAM", cut);
149     gMC->SetCut("CUTELE", cut);
150     gMC->SetCut("CUTNEU", cut);
151     gMC->SetCut("CUTHAD", cut);
152     gMC->SetCut("CUTMUO", cut);
153     gMC->SetCut("BCUTE",  cut); 
154     gMC->SetCut("BCUTM",  cut); 
155     gMC->SetCut("DCUTE",  cut); 
156     gMC->SetCut("DCUTM",  cut); 
157     gMC->SetCut("PPCUTM", cut);
158     gMC->SetCut("TOFMAX", tofmax); 
159
160
161     int     nParticles = 100;
162     if (gSystem->Getenv("CONFIG_NPARTICLES"))
163     {
164         nParticles = atoi(gSystem->Getenv("CONFIG_NPARTICLES"));
165     }
166
167     AliGenCocktail *gener = new AliGenCocktail();
168     gener->SetPhiRange(0, 360);
169     // Set pseudorapidity range from -8 to 8.
170     Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
171     Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
172     gener->SetThetaRange(thmin,thmax);
173     gener->SetOrigin(0, 0, 0);  //vertex position
174     gener->SetSigma(0, 0, 0);   //Sigma in (X,Y,Z) (cm) on IP position
175
176     AliGenHIJINGpara *hijingparam = new AliGenHIJINGpara(nParticles);
177     hijingparam->SetMomentumRange(0.2, 999);
178     gener->AddGenerator(hijingparam,"HIJING PARAM",1);
179
180 //    AliGenBox *genbox = new AliGenBox(nParticles);
181 //    genbox->SetPart(kGamma);
182 //    genbox->SetPtRange(0.3, 10.00);
183 //    gener->AddGenerator(genbox,"GENBOX GAMMA for PHOS",1);
184     gener->Init();
185
186
187     // 
188     // Activate this line if you want the vertex smearing to happen
189     // track by track
190     //
191     //gener->SetVertexSmear(perTrack); 
192     // Field (L3 0.4 T)
193     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
194     gAlice->SetField(field);    
195
196
197     Int_t   iABSO   = 1;
198     Int_t   iDIPO   = 1;
199     Int_t   iFMD    = 1;
200     Int_t   iFRAME  = 1;
201     Int_t   iHALL   = 1;
202     Int_t   iITS    = 1;
203     Int_t   iMAG    = 1;
204     Int_t   iMUON   = 1;
205     Int_t   iPHOS   = 1;
206     Int_t   iPIPE   = 1;
207     Int_t   iPMD    = 1;
208     Int_t   iHMPID  = 1;
209     Int_t   iSHIL   = 1;
210     Int_t   iT0     = 1;
211     Int_t   iTOF    = 1;
212     Int_t   iTPC    = 1;
213     Int_t   iTRD    = 1;
214     Int_t   iZDC    = 1;
215     Int_t   iEMCAL  = 1;
216     Int_t   iACORDE = 1;
217     Int_t   iVZERO  = 1;
218     rl->CdGAFile();
219     //=================== Alice BODY parameters =============================
220     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
221
222     if (iMAG)
223     {
224         //=================== MAG parameters ============================
225         // --- Start with Magnet since detector layouts may be depending ---
226         // --- on the selected Magnet dimensions ---
227         AliMAG *MAG = new AliMAG("MAG", "Magnet");
228     }
229
230
231     if (iABSO)
232     {
233         //=================== ABSO parameters ============================
234         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
235     }
236
237     if (iDIPO)
238     {
239         //=================== DIPO parameters ============================
240
241         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
242     }
243
244     if (iHALL)
245     {
246         //=================== HALL parameters ============================
247
248         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
249     }
250
251
252     if (iFRAME)
253     {
254         //=================== FRAME parameters ============================
255
256         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
257         FRAME->SetHoles(1);
258     }
259
260     if (iSHIL)
261     {
262         //=================== SHIL parameters ============================
263
264         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
265     }
266
267
268     if (iPIPE)
269     {
270         //=================== PIPE parameters ============================
271
272         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
273     }
274  
275     if (iITS)
276     {
277         //=================== ITS parameters ============================
278
279         AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
280     }
281
282     if (iTPC)
283     {
284         //============================ TPC parameters ===================
285         AliTPC *TPC = new AliTPCv2("TPC", "Default");
286     }
287
288
289     if (iTOF) {
290         //=================== TOF parameters ============================
291         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
292     }
293
294
295     if (iHMPID)
296     {
297         //=================== HMPID parameters ===========================
298         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
299
300     }
301
302
303     if (iZDC)
304     {
305         //=================== ZDC parameters ============================
306
307         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
308     }
309
310     if (iTRD)
311     {
312         //=================== TRD parameters ============================
313
314         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
315     }
316
317     if (iFMD)
318     {
319         //=================== FMD parameters ============================
320         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
321    }
322
323     if (iMUON)
324     {
325         //=================== MUON parameters ===========================
326         // New MUONv1 version (geometry defined via builders)
327         AliMUON *MUON = new AliMUONv1("MUON", "default");
328     }
329     //=================== PHOS parameters ===========================
330
331     if (iPHOS)
332     {
333         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
334     }
335
336
337     if (iPMD)
338     {
339         //=================== PMD parameters ============================
340         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
341     }
342
343     if (iT0)
344     {
345         //=================== T0 parameters ============================
346         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
347     }
348
349     if (iEMCAL)
350     {
351         //=================== EMCAL parameters ============================
352         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
353     }
354
355      if (iACORDE)
356     {
357         //=================== ACORDE parameters ============================
358         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
359     }
360
361      if (iVZERO)
362     {
363         //=================== VZERO parameters ============================
364         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
365     }
366
367      AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
368
369 }
370
371 Float_t EtaToTheta(Float_t arg){
372   return (180./TMath::Pi())*2.*atan(exp(-arg));
373 }
374
375
376 void LoadPythia()
377 {
378     // Load Pythia related libraries
379     gSystem->Load("liblhapdf.so");      // Parton density functions
380     gSystem->Load("libEGPythia6.so");   // TGenerator interface
381     gSystem->Load("libpythia6.so");     // Pythia
382     gSystem->Load("libAliPythia6.so");  // ALICE specific implementations
383 }