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