]> git.uio.no Git - u/mrichter/AliRoot.git/blob - ITS/UPGRADE/Config.C
Using EMCAL_COMPLETEV1 (Gustavo)
[u/mrichter/AliRoot.git] / ITS / UPGRADE / 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 "EVGEN/AliGenFixed.h"
21 #include "EVGEN/AliGenBox.h"
22 #include "STEER/AliMagWrapCheb.h"
23 #include "STRUCT/AliBODY.h"
24 #include "STRUCT/AliMAG.h"
25 #include "STRUCT/AliABSOv3.h"
26 #include "STRUCT/AliDIPOv3.h"
27 #include "STRUCT/AliHALLv3.h"
28 #include "STRUCT/AliFRAMEv2.h"
29 #include "STRUCT/AliSHILv3.h"
30 #include "STRUCT/AliPIPEv3.h"
31 #include "ITS/AliITSv11Hybrid.h"
32 #include "TPC/AliTPCv2.h"
33 #include "TOF/AliTOFv6T0.h"
34 #include "HMPID/AliHMPIDv3.h"
35 #include "ZDC/AliZDCv3.h"
36 #include "TRD/AliTRDv1.h"
37 #include "TRD/AliTRDgeometry.h"
38 #include "FMD/AliFMDv1.h"
39 #include "MUON/AliMUONv1.h"
40 #include "PHOS/AliPHOSv1.h"
41 #include "PMD/AliPMDv1.h"
42 #include "T0/AliT0v1.h"
43 #include "EMCAL/AliEMCALv2.h"
44 #include "ACORDE/AliACORDEv1.h"
45 #include "VZERO/AliVZEROv7.h"
46 #include <TVirtualMagField.h>
47 #endif
48
49 /* $Id$ */
50 enum PprTrigConf_t
51 {
52   kDefaultPPTrig, kDefaultPbPbTrig
53 };
54
55 const char * pprTrigConfName[] = {
56   "p-p","Pb-Pb"
57 };
58
59 Float_t EtaToTheta(Float_t arg);
60
61 static PprTrigConf_t strig = kDefaultPPTrig;// default PP trigger configuration
62
63 void Config()
64 {
65   // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
66   // Theta range given through pseudorapidity limits 22/6/2001
67
68   // Set Random Number seed
69   gRandom->SetSeed(0); // Set 0 to use the currecnt time
70
71
72   // libraries required by geant321
73 #if defined(__CINT__)
74   gSystem->Load("liblhapdf");
75   gSystem->Load("libEGPythia6");
76   gSystem->Load("libpythia6");
77   gSystem->Load("libAliPythia6");
78   gSystem->Load("libgeant321");
79 #endif
80
81   new     TGeant3TGeo("C++ Interface to Geant3");
82
83   AliRunLoader* rl=0x0;
84
85
86   rl = AliRunLoader::Open("galice.root",
87                           AliConfig::GetDefaultEventFolderName(),
88                           "recreate");
89   if (rl == 0x0)
90     {
91       gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
92       return;
93     }
94   rl->SetCompressionLevel(2);
95   rl->SetNumberOfEventsPerFile(2000);
96   gAlice->SetRunLoader(rl);
97
98   // Set the trigger configuration
99   // gAlice->SetTriggerDescriptor(pprTrigConfName[strig]);
100   //cout<<"Trigger configuration is set to  "<<pprTrigConfName[strig]<<endl;
101   AliSimulation::Instance()->SetTriggerConfig(pprTrigConfName[strig]);
102   cout<<"Trigger configuration is set to  pprTrigConfName[strig] "<<endl;
103
104   //
105   // Set External decayer
106   TVirtualMCDecayer *decayer = new AliDecayerPythia();
107
108   decayer->SetForceDecay(kAll);
109   decayer->Init();
110   gMC->SetExternalDecayer(decayer);
111   //=======================================================================
112   // ************* STEERING parameters FOR ALICE SIMULATION **************
113   // --- Specify event type to be tracked through the ALICE setup
114   // --- All positions are in cm, angles in degrees, and P and E in GeV
115
116
117   gMC->SetProcess("DCAY",1);
118   gMC->SetProcess("PAIR",1);
119   gMC->SetProcess("COMP",1);
120   gMC->SetProcess("PHOT",1);
121   gMC->SetProcess("PFIS",0);
122   gMC->SetProcess("DRAY",0);
123   gMC->SetProcess("ANNI",1);
124   gMC->SetProcess("BREM",1);
125   gMC->SetProcess("MUNU",1);
126   gMC->SetProcess("CKOV",1);
127   gMC->SetProcess("HADR",0);
128   gMC->SetProcess("LOSS",2);
129   gMC->SetProcess("MULS",1);
130   gMC->SetProcess("RAYL",1);
131
132   Float_t cut = 1.e-3;        // 1MeV cut by default
133   Float_t tofmax = 1.e10;
134
135   gMC->SetCut("CUTGAM", cut);
136   gMC->SetCut("CUTELE", cut);
137   gMC->SetCut("CUTNEU", cut);
138   gMC->SetCut("CUTHAD", cut);
139   gMC->SetCut("CUTMUO", cut);
140   gMC->SetCut("BCUTE",  cut); 
141   gMC->SetCut("BCUTM",  cut); 
142   gMC->SetCut("DCUTE",  cut); 
143   gMC->SetCut("DCUTM",  cut); 
144   gMC->SetCut("PPCUTM", cut);
145   gMC->SetCut("TOFMAX", tofmax); 
146
147   // Special generation for Valgrind tests
148   // Each detector is fired by few particles selected 
149   // to cover specific cases
150
151
152   // The cocktail itself
153
154   AliGenCocktail *gener = new AliGenCocktail();
155   gener->SetPhiRange(0, 360);
156   // Set pseudorapidity range from -8 to 8.
157   Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
158   Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
159   gener->SetThetaRange(thmin,thmax);
160   gener->SetOrigin(0., 0., 0);  //vertex position
161   gener->SetSigma(0., 0., 0);   //Sigma in (X,Y,Z) (cm) on IP position
162
163   AliGenBox *gbox2 = new AliGenBox(1);
164   gbox2->SetPtRange(0.99,1.0001);  //
165   gbox2->SetPhiRange(0,360.);
166   gbox2->SetThetaRange(40.,140.);
167   gbox2->SetPart(kPiMinus);
168   gener->AddGenerator(gbox2,"GENBOX PIONS for ITS",1);
169
170     
171     
172   gener->Init();
173
174
175   // 
176   // Activate this line if you want the vertex smearing to happen
177   // track by track
178   //
179   //VertexSmear_t perTrack;
180   //gener->SetVertexSmear(perTrack); 
181   // Field (L3 0.5 T)
182   //AliMagF* field = new AliMagF("map","map",2, -1.,1., 15, AliMagF::k5kGUniform);
183   //TGeoGlobalMagField::Instance()->SetField(field);
184   TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG));
185
186   Int_t   iABSO  =  0;
187   Int_t   iDIPO  =  0;
188   Int_t   iFMD   =  0;
189   Int_t   iFRAME =  0;
190   Int_t   iHALL  =  0;
191   Int_t   iITS   =  1;
192   Int_t   iMAG   =  0;
193   Int_t   iMUON  =  0;
194   Int_t   iPHOS  =  0;
195   Int_t   iPIPE  =  1;
196   Int_t   iPMD   =  0;
197   Int_t   iHMPID  =  0;
198   Int_t   iSHIL  =  0;
199   Int_t   iT0 =  0;
200   Int_t   iTOF   =  0;
201   Int_t   iTPC   =  0;
202   Int_t   iTRD   =  0;
203   Int_t   iZDC   =  0;
204   Int_t   iEMCAL =  0;
205   Int_t   iACORDE   =  0;
206   Int_t   iVZERO =  0;
207   rl->CdGAFile();
208   //=================== Alice BODY parameters =============================
209   AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
210
211   if (iMAG)
212     {
213       //=================== MAG parameters ============================
214       // --- Start with Magnet since detector layouts may be depending ---
215       // --- on the selected Magnet dimensions ---
216       AliMAG *MAG = new AliMAG("MAG", "Magnet");
217     }
218
219
220   if (iABSO)
221     {
222       //=================== ABSO parameters ============================
223       AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
224     }
225
226   if (iDIPO)
227     {
228       //=================== DIPO parameters ============================
229
230       AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
231     }
232
233   if (iHALL)
234     {
235       //=================== HALL parameters ============================
236
237       AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
238     }
239
240
241   if (iFRAME)
242     {
243       //=================== FRAME parameters ============================
244
245       AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
246       FRAME->SetHoles(1);
247     }
248
249   if (iSHIL)
250     {
251       //=================== SHIL parameters ============================
252
253       AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
254     }
255
256
257   if (iPIPE)
258     {
259       //=================== PIPE parameters ============================
260
261       AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
262     }
263  
264   if (iITS)
265     {
266       //=================== ITS parameters ============================
267
268       Int_t  nlayers =6;
269       TArrayD xsize(nlayers);
270       TArrayD zsize(nlayers);
271       TArrayD spess(nlayers);
272       TArrayD radii(nlayers);
273        
274       TArrayD spessCu(nlayers);
275       TArrayD radiiCu(nlayers);
276       TArrayS copper(nlayers);
277       TArrayD halfLengths(nlayers);
278
279       Bool_t bp=kFALSE;
280       Double_t radiusBP = 2.9;
281       Double_t widthBP = 0.08;
282       Double_t halflengthBP = 80.;
283
284       for(Int_t i=0; i<nlayers; i++){
285         Double_t xsz[6]={31.18*1e-04,41.6*1e-04,131.6*1e-04,131.6*1e-04,69.3*1e-04,69.3*1e-04};
286         Double_t zsz[6]={416*1e-04,416*1e-04,97*1e-04,97*1e-04,2875*1e-04,2875*1e-04};
287         Double_t Lhalf[6]={80.,80.,80.,80.,80.,80.};
288         Int_t npixHalf[6];
289         npixHalf[i]=(Int_t)(Lhalf[i]/zsz[i]); 
290         Double_t HalfL[6];
291         HalfL[i]=npixHalf[i]*zsz[i];
292         Double_t r[6]={4.,7.6,14.9,23.8,39.1,43.6};
293         radii.AddAt(r[i],i);
294         Int_t npixR[6];
295         npixR[i] = (Int_t)(2*TMath::Pi()*r[i]/xsz[i]);
296         Double_t xszInt[6];
297         xszInt[i]= 2*TMath::Pi()*r[i]/npixR[i];
298         xsize.AddAt(xszInt[i],i);
299         zsize.AddAt(zsz[i],i);
300         Double_t thick[6]={75.*1e-04,150.*1e-04,150.*1e-04,150.*1e-04,150.*1e-04,150.*1e-04};
301
302         radiiCu.AddAt(r[i]+thick[i],i);
303         halfLengths.AddAt(HalfL[i],i);
304         spess.AddAt(thick[i],i);
305         spessCu.AddAt(0.015,i);//micron
306         Int_t c[6]={1,1,1,1,1,1};
307         copper.AddAt(c[i],i);
308       }
309       AliITSupgrade *ITS  = new AliITSupgrade("ITS","ITS upgrade",spess,radii,halfLengths,radiiCu,spessCu,copper,bp,radiusBP, widthBP, halflengthBP );
310      
311       ITS->SetFullSegmentation(xsize,zsize);
312     }
313
314
315   if (iTPC)
316     {
317       //============================ TPC parameters ===================
318       AliTPC *TPC = new AliTPCv2("TPC", "Default");
319     }
320
321
322   if (iTOF) {
323     //=================== TOF parameters ============================
324     AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
325   }
326
327
328   if (iHMPID)
329     {
330       //=================== HMPID parameters ===========================
331       AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
332
333     }
334
335
336   if (iZDC)
337     {
338       //=================== ZDC parameters ============================
339
340       AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
341     }
342
343   if (iTRD)
344     {
345       //=================== TRD parameters ============================
346
347       AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
348       AliTRDgeometry *geoTRD = TRD->GetGeometry();
349       // Partial geometry: modules at 2,3,4,6,11,12,14,15
350       // starting at 6h in positive direction
351       geoTRD->SetSMstatus(0,0);
352       geoTRD->SetSMstatus(1,0);
353       geoTRD->SetSMstatus(5,0);
354       geoTRD->SetSMstatus(7,0);
355       geoTRD->SetSMstatus(8,0);
356       geoTRD->SetSMstatus(9,0);
357       geoTRD->SetSMstatus(10,0);
358       geoTRD->SetSMstatus(13,0);
359       geoTRD->SetSMstatus(16,0);
360       geoTRD->SetSMstatus(17,0);
361     }
362
363   if (iFMD)
364     {
365       //=================== FMD parameters ============================
366       AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
367     }
368
369   if (iMUON)
370     {
371       //=================== MUON parameters ===========================
372       // New MUONv1 version (geometry defined via builders)
373       AliMUON *MUON = new AliMUONv1("MUON","default");
374     }
375   //=================== PHOS parameters ===========================
376
377   if (iPHOS)
378     {
379       AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
380     }
381
382
383   if (iPMD)
384     {
385       //=================== PMD parameters ============================
386       AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
387     }
388
389   if (iT0)
390     {
391       //=================== T0 parameters ============================
392       AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
393     }
394
395   if (iEMCAL)
396     {
397       //=================== EMCAL parameters ============================
398       AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETEV1");
399     }
400
401   if (iACORDE)
402     {
403       //=================== ACORDE parameters ============================
404       AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
405     }
406
407   if (iVZERO)
408     {
409       //=================== VZERO parameters ============================
410       AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
411     }
412
413
414 }
415
416 Float_t EtaToTheta(Float_t arg){
417   return (180./TMath::Pi())*2.*atan(exp(-arg));
418 }