]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/Config.C
Fix for #86577: AliReconstruction::SlaveTerminate instantiate ALL Reconstructors
[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/AliMagF.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/OCDB !!!         #" << endl;
93          cout << "#                                                   #" << endl;
94          cout << "#####################################################" << endl;
95           
96          AliCDBManager::Instance()->SetDefaultStorage("local://$ALICE_ROOT/OCDB");
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     AliSimulation::Instance()->SetTriggerConfig("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     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG));
194
195     Int_t   iABSO   = 1;
196     Int_t   iDIPO   = 1;
197     Int_t   iFMD    = 1;
198     Int_t   iFRAME  = 1;
199     Int_t   iHALL   = 1;
200     Int_t   iITS    = 1;
201     Int_t   iMAG    = 1;
202     Int_t   iMUON   = 1;
203     Int_t   iPHOS   = 1;
204     Int_t   iPIPE   = 1;
205     Int_t   iPMD    = 1;
206     Int_t   iHMPID  = 1;
207     Int_t   iSHIL   = 1;
208     Int_t   iT0     = 1;
209     Int_t   iTOF    = 1;
210     Int_t   iTPC    = 1;
211     Int_t   iTRD    = 1;
212     Int_t   iZDC    = 1;
213     Int_t   iEMCAL  = 1;
214     Int_t   iACORDE = 1;
215     Int_t   iVZERO  = 1;
216     rl->CdGAFile();
217     //=================== Alice BODY parameters =============================
218     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
219
220     if (iMAG)
221     {
222         //=================== MAG parameters ============================
223         // --- Start with Magnet since detector layouts may be depending ---
224         // --- on the selected Magnet dimensions ---
225         AliMAG *MAG = new AliMAG("MAG", "Magnet");
226     }
227
228
229     if (iABSO)
230     {
231         //=================== ABSO parameters ============================
232         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
233     }
234
235     if (iDIPO)
236     {
237         //=================== DIPO parameters ============================
238
239         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
240     }
241
242     if (iHALL)
243     {
244         //=================== HALL parameters ============================
245
246         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
247     }
248
249
250     if (iFRAME)
251     {
252         //=================== FRAME parameters ============================
253
254         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
255         FRAME->SetHoles(1);
256     }
257
258     if (iSHIL)
259     {
260         //=================== SHIL parameters ============================
261
262         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
263     }
264
265
266     if (iPIPE)
267     {
268         //=================== PIPE parameters ============================
269
270         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
271     }
272  
273     if (iITS)
274     {
275         //=================== ITS parameters ============================
276
277         AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
278     }
279
280     if (iTPC)
281     {
282         //============================ TPC parameters ===================
283         AliTPC *TPC = new AliTPCv2("TPC", "Default");
284     }
285
286
287     if (iTOF) {
288         //=================== TOF parameters ============================
289         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
290     }
291
292
293     if (iHMPID)
294     {
295         //=================== HMPID parameters ===========================
296         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
297
298     }
299
300
301     if (iZDC)
302     {
303         //=================== ZDC parameters ============================
304
305         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
306     }
307
308     if (iTRD)
309     {
310         //=================== TRD parameters ============================
311
312         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
313     }
314
315     if (iFMD)
316     {
317         //=================== FMD parameters ============================
318         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
319    }
320
321     if (iMUON)
322     {
323         //=================== MUON parameters ===========================
324         // New MUONv1 version (geometry defined via builders)
325         AliMUON *MUON = new AliMUONv1("MUON", "default");
326     }
327     //=================== PHOS parameters ===========================
328
329     if (iPHOS)
330     {
331         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
332     }
333
334
335     if (iPMD)
336     {
337         //=================== PMD parameters ============================
338         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
339     }
340
341     if (iT0)
342     {
343         //=================== T0 parameters ============================
344         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
345     }
346
347     if (iEMCAL)
348     {
349         //=================== EMCAL parameters ============================
350         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETEV1");
351     }
352
353      if (iACORDE)
354     {
355         //=================== ACORDE parameters ============================
356         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
357     }
358
359      if (iVZERO)
360     {
361         //=================== VZERO parameters ============================
362         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
363     }
364
365      AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
366
367 }
368
369 Float_t EtaToTheta(Float_t arg){
370   return (180./TMath::Pi())*2.*atan(exp(-arg));
371 }
372
373
374 void LoadPythia()
375 {
376     // Load Pythia related libraries
377     gSystem->Load("liblhapdf.so");      // Parton density functions
378     gSystem->Load("libEGPythia6.so");   // TGenerator interface
379     gSystem->Load("libpythia6.so");     // Pythia
380     gSystem->Load("libAliPythia6.so");  // ALICE specific implementations
381 }