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