]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/macros/Config.C
Extra header added to the list
[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 AliGenerator *GenParamCalo(Int_t nPart, Int_t type, TString calo);
48
49 Int_t  year = 2012;
50 Bool_t checkGeoAndRun = kFALSE;
51
52 void Config()
53 {
54   //AliLog::SetGlobalDebugLevel(2);
55   
56   // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
57   // Theta range given through pseudorapidity limits 22/6/2001
58   
59   // Set Random Number seed
60   //gRandom->SetSeed(123456); // Set 0 to use the currecnt time
61   AliLog::Message(AliLog::kInfo, Form("Seed for random number generation = %d",gRandom->GetSeed()), "Config.C", "Config.C", "Config()","Config.C", __LINE__);
62   
63   // Load Pythia libraries                                        
64   LoadPythia();
65   
66   // libraries required by geant321
67 #if defined(__CINT__)
68   gSystem->Load("libgeant321");
69 #endif
70   
71   new     TGeant3TGeo("C++ Interface to Geant3");
72   
73   AliRunLoader* rl=0x0;
74   
75   AliLog::Message(AliLog::kInfo, "Creating Run Loader", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
76   
77   rl = AliRunLoader::Open("galice.root",
78                           AliConfig::GetDefaultEventFolderName(),
79                           "recreate");
80   if (rl == 0x0)
81   {
82     gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
83     return;
84   }
85   rl->SetCompressionLevel(2);
86   rl->SetNumberOfEventsPerFile(10000);
87   gAlice->SetRunLoader(rl);
88   
89   // Set the trigger configuration
90   AliSimulation::Instance()->SetTriggerConfig("Pb-Pb");
91   cout<<"Trigger configuration is set to  Pb-Pb"<<endl;
92   
93   //
94   // Set External decayer
95   TVirtualMCDecayer *decayer = new AliDecayerPythia();
96   
97   decayer->SetForceDecay(kAll);
98   decayer->Init();
99   gMC->SetExternalDecayer(decayer);
100   //=======================================================================
101   // ************* STEERING parameters FOR ALICE SIMULATION **************
102   // --- Specify event type to be tracked through the ALICE setup
103   // --- All positions are in cm, angles in degrees, and P and E in GeV
104   
105   
106   gMC->SetProcess("DCAY",1);
107   gMC->SetProcess("PAIR",1);
108   gMC->SetProcess("COMP",1);
109   gMC->SetProcess("PHOT",1);
110   gMC->SetProcess("PFIS",0);
111   gMC->SetProcess("DRAY",0);
112   gMC->SetProcess("ANNI",1);
113   gMC->SetProcess("BREM",1);
114   gMC->SetProcess("MUNU",1);
115   gMC->SetProcess("CKOV",1);
116   gMC->SetProcess("HADR",1);
117   gMC->SetProcess("LOSS",2);
118   gMC->SetProcess("MULS",1);
119   gMC->SetProcess("RAYL",1);
120   
121   Float_t cut    = 1.e-3;  // 1MeV cut by default
122   Float_t tofmax = 1.e10;
123   
124   gMC->SetCut("CUTGAM", cut);
125   gMC->SetCut("CUTELE", cut);
126   gMC->SetCut("CUTNEU", cut);
127   gMC->SetCut("CUTHAD", cut);
128   gMC->SetCut("CUTMUO", cut);
129   gMC->SetCut("BCUTE",  cut); 
130   gMC->SetCut("BCUTM",  cut); 
131   gMC->SetCut("DCUTE",  cut); 
132   gMC->SetCut("DCUTM",  cut); 
133   gMC->SetCut("PPCUTM", cut);
134   gMC->SetCut("TOFMAX", tofmax); 
135   
136   
137   int     nParticles = 1;
138   if (gSystem->Getenv("CONFIG_NPARTICLES"))
139   {
140     nParticles = atoi(gSystem->Getenv("CONFIG_NPARTICLES"));
141   }
142   
143   if (gSystem->Getenv("CONFIG_YEAR"))
144   {
145     year = atoi(gSystem->Getenv("CONFIG_YEAR"));
146   }
147   
148   AliGenBox *gener = new AliGenBox(nParticles);
149   gener->SetMomentumRange(10.,10.);
150   
151   gener->SetPart(22);
152   
153   if     (year == 2010)
154     gener->SetPhiRange(80.0,120.0);
155   else if(year == 2011)
156     gener->SetPhiRange(80.0,180.0);
157   else if(year == 2012 || year == 2013)
158     gener->SetPhiRange(80.0,190.0);
159   else
160     gener->SetPhiRange(80.0,330.0); // Include DCal
161   
162   gener->SetThetaRange(EtaToTheta(0.7), EtaToTheta(-0.7));
163   
164 //  AliGenLib* lib   = new AliGenPHOSlib();
165 //  Int_t      type  = AliGenPHOSlib::kEtaFlat;
166 //  AliGenParam *gener = new AliGenParam(1,lib,type,"");
167 //      gener->SetMomentumRange(0,999);
168 //      gener->SetPtRange(1,30);
169 //      gener->SetPhiRange(80, 200.);
170 //      gener->SetYRange(-2,2);
171 //      gener->SetThetaRange(EtaToTheta(0.7),EtaToTheta(-0.7));
172 //      gener->SetCutOnChild(1);
173 //  gener->SetChildPtRange(0.1,30);
174 //      gener->SetChildThetaRange(EtaToTheta(0.7),EtaToTheta(-0.7));
175 //  gener->SetChildPhiRange(80, 180.);
176 //      gener->SetOrigin(0,0,0);        //vertex position
177 //      gener->SetSigma(0,0,5.3);       //Sigma in (X,Y,Z) (cm) on IP position
178 //      gener->SetForceDecay(kGammaEM);
179 //
180 //  //gener->SetTrackingFlag(0);
181   
182 //  AliGenCocktail *gener = new AliGenCocktail();
183 //  gener->SetProjectile("A", 208, 82);
184 //  gener->SetTarget    ("A", 208, 82);
185 //  
186 //  // 1 Pi0 in EMCAL, 2010 configuration, 4 SM
187 //  AliGenParam *gEMCPi0 = GenParamCalo(1, AliGenPHOSlib::kPi0Flat, "EMCAL");
188 //  gener->AddGenerator(gEMCPi0,"pi0EMC", 1);
189 //  
190 //  // 1 Eta in EMCAL, 2010 configuration, 4 SM
191 //  AliGenParam *gEMCEta = GenParamCalo(1, AliGenPHOSlib::kEtaFlat, "EMCAL");
192 //  gener->AddGenerator(gEMCEta,"etaEMC", 1);
193 //  
194  
195   gener->Init();
196   
197   // 
198   // Activate this line if you want the vertex smearing to happen
199   // track by track
200   //
201   //gener->SetVertexSmear(perTrack); 
202   // Field (L3 0.5 T)
203   TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG));
204   
205   Int_t   iABSO  =  0;
206   Int_t   iDIPO  =  0;
207   Int_t   iFMD   =  0;
208   Int_t   iFRAME =  0;
209   Int_t   iHALL  =  0;
210   Int_t   iITS   =  0;
211   Int_t   iMAG   =  0;
212   Int_t   iMUON  =  0;
213   Int_t   iPHOS  =  0;
214   Int_t   iPIPE  =  0;
215   Int_t   iPMD   =  0;
216   Int_t   iHMPID =  0;
217   Int_t   iSHIL  =  0;
218   Int_t   iT0    =  0;
219   Int_t   iTOF   =  0;
220   Int_t   iTPC   =  0;
221   Int_t   iTRD   =  0;
222   Int_t   iZDC   =  0;
223   Int_t   iEMCAL =  1;
224   Int_t   iACORDE=  0;
225   Int_t   iVZERO =  0;
226   rl->CdGAFile();
227   //=================== Alice BODY parameters =============================
228   AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
229   
230   if (iMAG)
231   {
232     //=================== MAG parameters ============================
233     // --- Start with Magnet since detector layouts may be depending ---
234     // --- on the selected Magnet dimensions ---
235     AliMAG *MAG = new AliMAG("MAG", "Magnet");
236   }
237   
238   
239   if (iABSO)
240   {
241     //=================== ABSO parameters ============================
242     AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
243   }
244   
245   if (iDIPO)
246   {
247     //=================== DIPO parameters ============================
248     
249     AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
250   }
251   
252   if (iHALL)
253   {
254     //=================== HALL parameters ============================
255     
256     AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
257   }
258   
259   
260   if (iFRAME)
261   {
262     //=================== FRAME parameters ============================
263     
264     AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
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     AliITS *ITS  = new AliITSv11("ITS","ITS v11");
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 AliHMPIDv3("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     AliTRDgeometry *geoTRD = TRD->GetGeometry();
323     // starting at 3h in positive direction
324     if(year==2011 || year == 2010)
325     { // not sure if good for 2010
326       // Partial geometry: modules at 0,1,7,8,9,10,11,15,16,17
327       printf("*** TRD configuration for 2011\n");
328       geoTRD->SetSMstatus(2,0);
329       geoTRD->SetSMstatus(3,0);
330       geoTRD->SetSMstatus(4,0);
331       geoTRD->SetSMstatus(5,0);
332       geoTRD->SetSMstatus(6,0);
333       geoTRD->SetSMstatus(12,0);
334       geoTRD->SetSMstatus(13,0);
335       geoTRD->SetSMstatus(14,0);
336     }
337     else if(year==2012)
338     {
339       printf("*** TRD configuration for 2012\n");
340       geoTRD->SetSMstatus(4,0);
341       geoTRD->SetSMstatus(5,0);
342       geoTRD->SetSMstatus(12,0);
343       geoTRD->SetSMstatus(13,0);
344       geoTRD->SetSMstatus(14,0); 
345     }
346   }
347   
348   if (iFMD)
349   {
350     //=================== FMD parameters ============================
351     AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
352   }
353   
354   if (iMUON)
355   {
356     //=================== MUON parameters ===========================
357     // New MUONv1 version (geometry defined via builders)
358     AliMUON *MUON = new AliMUONv1("MUON", "default");
359   }
360   //=================== PHOS parameters ===========================
361   
362   if (iPHOS)
363   {
364     AliPHOS *PHOS = new AliPHOSv1("PHOS", "noCPV_Modules123");
365   }
366   
367   
368   if (iPMD)
369   {
370     //=================== PMD parameters ============================
371     AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
372   }
373   
374   if (iT0)
375   {
376     //=================== T0 parameters ============================
377     AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
378   }
379   
380   if (iEMCAL)
381   {
382     //=================== EMCAL parameters ============================
383     AliEMCAL *EMCAL = 0;
384     if      (year == 2010)  // d phi = 40 degrees
385       EMCAL = new AliEMCALv2("EMCAL", "EMCAL_FIRSTYEARV1",    checkGeoAndRun);
386     else if (year == 2011)  // d phi = 100 degrees
387       EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETEV1",     checkGeoAndRun);
388     else if (year == 2012 || year == 2013)   // d phi = 107 degrees
389       EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE12SMV1", checkGeoAndRun);
390     else
391       EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE12SMV1_DCAL_8SM", checkGeoAndRun); // EMCAL+DCAL dphi = 107 (EMCAL) + 33 (gap) + 67 (DCAL)
392   }
393   
394   if (iACORDE)
395   {
396     //=================== ACORDE parameters ============================
397     AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
398   }
399   
400   if (iVZERO)
401   {
402     //=================== ACORDE parameters ============================
403     AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
404   }
405   
406   AliLog::Message(AliLog::kInfo, "End of Config", "Config.C", "Config.C", "Config()"," Config.C", __LINE__);
407   
408 }
409
410 Float_t EtaToTheta(Float_t arg){
411   return (180./TMath::Pi())*2.*atan(exp(-arg));
412 }
413
414 void LoadPythia()
415 {
416   // Load Pythia related libraries                                                                
417   gSystem->Load("liblhapdf");      // Parton density functions                                 
418   gSystem->Load("libEGPythia6");   // TGenerator interface                                     
419   gSystem->Load("libpythia6");     // Pythia                                                   
420   gSystem->Load("libAliPythia6");  // ALICE specific
421   // implementations                           
422 }
423
424
425 AliGenerator * GenParamCalo(Int_t nPart, Int_t type, TString calo)
426 {
427   // nPart of type (Pi0, Eta, Pi0Flat, EtaFlat, ...) in EMCAL or PHOS
428   // CAREFUL EMCAL year 2010 configuration
429   AliGenParam *gener = new AliGenParam(nPart,new AliGenPHOSlib(),type,"");
430   
431   // meson cuts
432   gener->SetMomentumRange(0,999);
433   gener->SetYRange(-2,2);
434   gener->SetPtRange(1,30);
435   // photon cuts
436   gener->SetForceDecay(kGammaEM); // Ensure the decays are photons
437   gener->SetCutOnChild(1);
438   gener->SetChildPtRange(0.,30);
439   
440   if(calo=="EMCAL")
441   {
442     //meson acceptance
443     gener->SetPhiRange(80., 100.); // year 2010
444     //gener->SetPhiRange(80., 180.); // year 2011
445     gener->SetThetaRange(EtaToTheta(0.7),EtaToTheta(-0.7));
446     //decay acceptance
447     gener->SetChildThetaRange(EtaToTheta(0.7),EtaToTheta(-0.7));
448     gener->SetChildPhiRange(80., 100.); // year 2010
449     //gener->SetChildPhiRange(80., 180.); // year 2011
450   }
451   else if(calo=="PHOS")
452   {
453     //meson acceptance
454     gener->SetPhiRange(260., 320.);
455     gener->SetThetaRange(EtaToTheta(0.13),EtaToTheta(-0.13));
456     //decay acceptance
457     gener->SetChildThetaRange(EtaToTheta(0.13),EtaToTheta(-0.13));
458     gener->SetChildPhiRange(260., 320.);
459   }
460   
461   return gener;
462   
463 }
464