]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/macros/Config.C
pdf figures to compile with pdflatex
[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/AliITSv11.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 = 2012;
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 = 2;
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   
150   if     (year == 2010)
151     gener->SetPhiRange(80.0,120.0);
152   else if(year == 2011)
153     gener->SetPhiRange(80.0,180.0);
154   else
155     gener->SetPhiRange(80.0,190.0);
156   
157   gener->SetThetaRange(EtaToTheta(0.7), EtaToTheta(-0.7));
158   
159   gener->SetOrigin(0,0,0);        //vertex position
160   gener->SetSigma(0,0,0);         //Sigma in (X,Y,Z) (cm) on IP position
161   gener->SetPart(kGamma);
162   gener->Init();
163   
164   // 
165   // Activate this line if you want the vertex smearing to happen
166   // track by track
167   //
168   //gener->SetVertexSmear(perTrack); 
169   // Field (L3 0.5 T)
170   TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG));
171   
172   Int_t   iABSO  =  1;
173   Int_t   iDIPO  =  1;
174   Int_t   iFMD   =  1;
175   Int_t   iFRAME =  1;
176   Int_t   iHALL  =  1;
177   Int_t   iITS   =  1;
178   Int_t   iMAG   =  1;
179   Int_t   iMUON  =  1;
180   Int_t   iPHOS  =  1;
181   Int_t   iPIPE  =  1;
182   Int_t   iPMD   =  1;
183   Int_t   iHMPID =  1;
184   Int_t   iSHIL  =  1;
185   Int_t   iT0    =  1;
186   Int_t   iTOF   =  1;
187   Int_t   iTPC   =  1;
188   Int_t   iTRD   =  1;
189   Int_t   iZDC   =  1;
190   Int_t   iEMCAL =  1;
191   Int_t   iACORDE=  0;
192   Int_t   iVZERO =  1;
193   rl->CdGAFile();
194   //=================== Alice BODY parameters =============================
195   AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
196   
197   if (iMAG)
198   {
199     //=================== MAG parameters ============================
200     // --- Start with Magnet since detector layouts may be depending ---
201     // --- on the selected Magnet dimensions ---
202     AliMAG *MAG = new AliMAG("MAG", "Magnet");
203   }
204   
205   
206   if (iABSO)
207   {
208     //=================== ABSO parameters ============================
209     AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
210   }
211   
212   if (iDIPO)
213   {
214     //=================== DIPO parameters ============================
215     
216     AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
217   }
218   
219   if (iHALL)
220   {
221     //=================== HALL parameters ============================
222     
223     AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
224   }
225   
226   
227   if (iFRAME)
228   {
229     //=================== FRAME parameters ============================
230     
231     AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
232   }
233   
234   if (iSHIL)
235   {
236     //=================== SHIL parameters ============================
237     
238     AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
239   }
240   
241   
242   if (iPIPE)
243   {
244     //=================== PIPE parameters ============================
245     
246     AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
247   }
248   
249   if(iITS) {
250     
251     //=================== ITS parameters ============================
252     
253     AliITS *ITS  = new AliITSv11("ITS","ITS v11");
254   }
255   
256   if (iTPC)
257   {
258     //============================ TPC parameters ===================
259     AliTPC *TPC = new AliTPCv2("TPC", "Default");
260   }
261   
262   
263   if (iTOF) {
264     //=================== TOF parameters ============================
265     AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
266   }
267   
268   
269   if (iHMPID)
270   {
271     //=================== HMPID parameters ===========================
272     AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
273     
274   }
275   
276   
277   if (iZDC)
278   {
279     //=================== ZDC parameters ============================
280     
281     AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
282   }
283   
284   if (iTRD)
285   {
286     //=================== TRD parameters ============================
287     
288     AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
289     AliTRDgeometry *geoTRD = TRD->GetGeometry();
290     // starting at 3h in positive direction
291     if(year==2011 || year == 2010)
292     { // not sure if good for 2010
293       // Partial geometry: modules at 0,1,7,8,9,10,11,15,16,17
294       printf("*** TRD configuration for 2011\n");
295       geoTRD->SetSMstatus(2,0);
296       geoTRD->SetSMstatus(3,0);
297       geoTRD->SetSMstatus(4,0);
298       geoTRD->SetSMstatus(5,0);
299       geoTRD->SetSMstatus(6,0);
300       geoTRD->SetSMstatus(12,0);
301       geoTRD->SetSMstatus(13,0);
302       geoTRD->SetSMstatus(14,0);
303     }
304     else if(year==2012)
305     {
306       printf("*** TRD configuration for 2012\n");
307       geoTRD->SetSMstatus(4,0);
308       geoTRD->SetSMstatus(5,0);
309       geoTRD->SetSMstatus(12,0);
310       geoTRD->SetSMstatus(13,0);
311       geoTRD->SetSMstatus(14,0); 
312     }
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", "noCPV_Modules123");
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     if      (year == 2010)  // d phi = 40 degrees
351       AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_FIRSTYEARV1",    checkGeoAndRun);
352     else if (year == 2011)  // d phi = 100 degrees
353       AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETEV1",     checkGeoAndRun);
354     else if (year > 2011)   // d phi = 110 degrees
355       AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE12SMV1", checkGeoAndRun);
356     else // Old configuration with 110 degrees but not perfect geometry
357       AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE",       checkGeoAndRun);
358   }
359   
360   if (iACORDE)
361   {
362     //=================== ACORDE parameters ============================
363     AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
364   }
365   
366   if (iVZERO)
367   {
368     //=================== ACORDE parameters ============================
369     AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
370   }
371   
372   AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
373   
374 }
375
376 Float_t EtaToTheta(Float_t arg){
377   return (180./TMath::Pi())*2.*atan(exp(-arg));
378 }
379
380 void LoadPythia()
381 {
382   // Load Pythia related libraries                                                                
383   gSystem->Load("liblhapdf.so");      // Parton density functions                                 
384   gSystem->Load("libEGPythia6.so");   // TGenerator interface                                     
385   gSystem->Load("libpythia6.so");     // Pythia                                                   
386   gSystem->Load("libAliPythia6.so");  // ALICE specific
387   // implementations                           
388 }