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