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