]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/gun/Config.C
New default geometry descriptions: HMPIDv2 and ITSv11Hybrid
[u/mrichter/AliRoot.git] / test / gun / 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/AliMagFMaps.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/AliACORDEv0.h"
45 #include "VZERO/AliVZEROv7.h"
46 #endif
47
48 enum PprTrigConf_t
49 {
50     kDefaultPPTrig, kDefaultPbPbTrig
51 };
52
53 const char * pprTrigConfName[] = {
54     "p-p","Pb-Pb"
55 };
56
57 enum PprGeo_t
58   {
59     kHoles, kNoHoles
60   };
61
62 Float_t EtaToTheta(Float_t arg);
63
64 static PprTrigConf_t strig = kDefaultPPTrig;// default PP trigger configuration
65 static PprGeo_t geo = kHoles;
66
67 void Config()
68 {
69     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
70     // Theta range given through pseudorapidity limits 22/6/2001
71
72     // Set Random Number seed
73   gRandom->SetSeed(123456); // Set 0 to use the currecnt time
74
75
76    // libraries required by geant321
77 #if defined(__CINT__)
78     gSystem->Load("liblhapdf");
79     gSystem->Load("libEGPythia6");
80     gSystem->Load("libpythia6");
81     gSystem->Load("libAliPythia6");
82     gSystem->Load("libgeant321");
83 #endif
84
85     new     TGeant3TGeo("C++ Interface to Geant3");
86
87     AliRunLoader* rl=0x0;
88
89
90     rl = AliRunLoader::Open("galice.root",
91                             AliConfig::GetDefaultEventFolderName(),
92                             "recreate");
93     if (rl == 0x0)
94       {
95         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
96         return;
97       }
98     rl->SetCompressionLevel(2);
99     rl->SetNumberOfEventsPerFile(3);
100     gAlice->SetRunLoader(rl);
101
102     // Set the trigger configuration
103     gAlice->SetTriggerDescriptor(pprTrigConfName[strig]);
104     cout<<"Trigger configuration is set to  "<<pprTrigConfName[strig]<<endl;
105
106     //
107     // Set External decayer
108     TVirtualMCDecayer *decayer = new AliDecayerPythia();
109
110     decayer->SetForceDecay(kAll);
111     decayer->Init();
112     gMC->SetExternalDecayer(decayer);
113     //=======================================================================
114     // ************* STEERING parameters FOR ALICE SIMULATION **************
115     // --- Specify event type to be tracked through the ALICE setup
116     // --- All positions are in cm, angles in degrees, and P and E in GeV
117
118
119     gMC->SetProcess("DCAY",1);
120     gMC->SetProcess("PAIR",1);
121     gMC->SetProcess("COMP",1);
122     gMC->SetProcess("PHOT",1);
123     gMC->SetProcess("PFIS",0);
124     gMC->SetProcess("DRAY",0);
125     gMC->SetProcess("ANNI",1);
126     gMC->SetProcess("BREM",1);
127     gMC->SetProcess("MUNU",1);
128     gMC->SetProcess("CKOV",1);
129     gMC->SetProcess("HADR",1);
130     gMC->SetProcess("LOSS",2);
131     gMC->SetProcess("MULS",1);
132     gMC->SetProcess("RAYL",1);
133
134     Float_t cut = 1.e-3;        // 1MeV cut by default
135     Float_t tofmax = 1.e10;
136
137     gMC->SetCut("CUTGAM", cut);
138     gMC->SetCut("CUTELE", cut);
139     gMC->SetCut("CUTNEU", cut);
140     gMC->SetCut("CUTHAD", cut);
141     gMC->SetCut("CUTMUO", cut);
142     gMC->SetCut("BCUTE",  cut); 
143     gMC->SetCut("BCUTM",  cut); 
144     gMC->SetCut("DCUTE",  cut); 
145     gMC->SetCut("DCUTM",  cut); 
146     gMC->SetCut("PPCUTM", cut);
147     gMC->SetCut("TOFMAX", tofmax); 
148
149     // Special generation for Valgrind tests
150     // Each detector is fired by few particles selected 
151     // to cover specific cases
152
153
154     // The cocktail iitself
155
156     AliGenCocktail *gener = new AliGenCocktail();
157     gener->SetPhiRange(0, 360);
158     // Set pseudorapidity range from -8 to 8.
159     Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
160     Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
161     gener->SetThetaRange(thmin,thmax);
162     gener->SetOrigin(0, 0, 0);  //vertex position
163     gener->SetSigma(0, 0, 0);   //Sigma in (X,Y,Z) (cm) on IP position
164
165
166     // Particle guns for the barrel part (taken from RichConfig)
167
168     AliGenFixed *pG1=new AliGenFixed(1);
169     pG1->SetPart(kProton);
170     pG1->SetMomentum(2.5);
171     pG1->SetTheta(109.5-3);
172     pG1->SetPhi(10);
173     gener->AddGenerator(pG1,"g1",1);
174     
175     AliGenFixed *pG2=new AliGenFixed(1);
176     pG2->SetPart(kPiPlus);
177     pG2->SetMomentum(1.0);
178     pG2->SetTheta( 90.0-3);
179     pG2->SetPhi(10);
180     gener->AddGenerator(pG2,"g2",1);
181
182     AliGenFixed *pG3=new AliGenFixed(1);
183     pG3->SetPart(kPiMinus);
184     pG3->SetMomentum(1.5);
185     pG3->SetTheta(109.5-3);
186     pG3->SetPhi(30);
187     gener->AddGenerator(pG3,"g3",1);
188     
189     AliGenFixed *pG4=new AliGenFixed(1);
190     pG4->SetPart(kKPlus);
191     pG4->SetMomentum(0.7);
192     pG4->SetTheta( 90.0-3);
193     pG4->SetPhi(30);
194     gener->AddGenerator(pG4,"g4",1);
195     
196     AliGenFixed *pG5=new AliGenFixed(1);
197     pG5->SetPart(kKMinus);
198     pG5->SetMomentum(1.0);
199     pG5->SetTheta( 70.0-3);
200     pG5->SetPhi(30);
201     gener->AddGenerator(pG5,"g5",1);
202     
203     AliGenFixed *pG6=new AliGenFixed(1);
204     pG6->SetPart(kProtonBar);
205     pG6->SetMomentum(2.5);
206     pG6->SetTheta( 90.0-3);
207     pG6->SetPhi(50);
208     gener->AddGenerator(pG6,"g6",1);
209     
210     AliGenFixed *pG7=new AliGenFixed(1);
211     pG7->SetPart(kPiMinus);
212     pG7->SetMomentum(0.7);
213     pG7->SetTheta( 70.0-3);
214     pG7->SetPhi(50);
215     gener->AddGenerator(pG7,"g7",1);
216
217     // Electrons for TRD
218
219     AliGenFixed *pG8=new AliGenFixed(1);
220     pG8->SetPart(kElectron);
221     pG8->SetMomentum(1.2);
222     pG8->SetTheta( 95.0);
223     pG8->SetPhi(190);
224     gener->AddGenerator(pG8,"g8",1);
225
226     AliGenFixed *pG9=new AliGenFixed(1);
227     pG9->SetPart(kPositron);
228     pG9->SetMomentum(1.2);
229     pG9->SetTheta( 85.0);
230     pG9->SetPhi(190);
231     gener->AddGenerator(pG9,"g9",1);
232
233     // PHOS
234
235     AliGenBox *gphos = new AliGenBox(1);
236     gphos->SetMomentumRange(10,11.);
237     gphos->SetPhiRange(270.5,270.7);
238     gphos->SetThetaRange(90.5,90.7);
239     gphos->SetPart(kGamma);
240     gener->AddGenerator(gphos,"GENBOX GAMMA for PHOS",1);
241
242     // EMCAL
243
244     AliGenBox *gemcal = new AliGenBox(1);
245     gemcal->SetMomentumRange(10,11.);
246     gemcal->SetPhiRange(90.5,199.5);
247     gemcal->SetThetaRange(90.5,90.7);
248     gemcal->SetPart(kGamma);
249     gener->AddGenerator(gemcal,"GENBOX GAMMA for EMCAL",1);
250
251     // MUON
252     AliGenBox * gmuon1 = new AliGenBox(1);
253     gmuon1->SetMomentumRange(20.,20.1);
254     gmuon1->SetPhiRange(0., 360.);         
255     gmuon1->SetThetaRange(171.000,178.001);
256     gmuon1->SetPart(kMuonMinus);           // Muons
257     gener->AddGenerator(gmuon1,"GENBOX MUON1",1);
258
259     AliGenBox * gmuon2 = new AliGenBox(1);
260     gmuon2->SetMomentumRange(20.,20.1);
261     gmuon2->SetPhiRange(0., 360.);         
262     gmuon2->SetThetaRange(171.000,178.001);
263     gmuon2->SetPart(kMuonPlus);           // Muons
264     gener->AddGenerator(gmuon2,"GENBOX MUON1",1);
265
266     //TOF
267     AliGenFixed *gtof=new AliGenFixed(1);
268     gtof->SetPart(kProton);
269     gtof->SetMomentum(2.5);
270     gtof->SetTheta(95);
271     gtof->SetPhi(340);
272     gener->AddGenerator(gtof,"Proton for TOF",1);
273
274     //FMD1
275     AliGenFixed *gfmd1=new AliGenFixed(1);
276     gfmd1->SetPart(kGamma);
277     gfmd1->SetMomentum(25);
278     gfmd1->SetTheta(1.8);
279     gfmd1->SetPhi(10);
280     gener->AddGenerator(gfmd1,"Gamma for FMD1",1);
281     
282     //FMD2i
283     AliGenFixed *gfmd2i=new AliGenFixed(1);
284     gfmd2i->SetPart(kPiPlus);
285     gfmd2i->SetMomentum(1.5);
286     gfmd2i->SetTheta(7.3);
287     gfmd2i->SetPhi(20);
288     gener->AddGenerator(gfmd2i,"Pi+ for FMD2i",1);
289     
290     //FMD2o
291     AliGenFixed *gfmd2o=new AliGenFixed(1);
292     gfmd2o->SetPart(kPiMinus);
293     gfmd2o->SetMomentum(1.5);
294     gfmd2o->SetTheta(16.1);
295     gfmd2o->SetPhi(30);
296     gener->AddGenerator(gfmd2o,"Pi- for FMD2o",1);
297     
298     //FMD3o
299     AliGenFixed *gfmd3o=new AliGenFixed(1);
300     gfmd3o->SetPart(kPiPlus);
301     gfmd3o->SetMomentum(1.5);
302     gfmd3o->SetTheta(163.9);
303     gfmd3o->SetPhi(40);
304     gener->AddGenerator(gfmd3o,"Pi+ for FMD3o",1);
305     
306     //FMD3i
307     AliGenFixed *gfmd3i=new AliGenFixed(1);
308     gfmd3i->SetPart(kPiMinus);
309     gfmd3i->SetMomentum(1.5);
310     gfmd3i->SetTheta(170.5);
311     gfmd3i->SetPhi(50);
312     gener->AddGenerator(gfmd3i,"Pi- for FMD3i",1);
313     
314     gener->Init();
315
316
317     // 
318     // Activate this line if you want the vertex smearing to happen
319     // track by track
320     //
321     //gener->SetVertexSmear(perTrack); 
322     // Field (L3 0.5 T)
323     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 2);
324     field->SetL3ConstField(1); //Using const. field in the barrel if 0
325     gAlice->SetField(field);    
326
327
328     Int_t   iABSO  =  1;
329     Int_t   iDIPO  =  1;
330     Int_t   iFMD   =  1;
331     Int_t   iFRAME =  1;
332     Int_t   iHALL  =  1;
333     Int_t   iITS   =  1;
334     Int_t   iMAG   =  1;
335     Int_t   iMUON  =  1;
336     Int_t   iPHOS  =  1;
337     Int_t   iPIPE  =  1;
338     Int_t   iPMD   =  1;
339     Int_t   iHMPID  =  1;
340     Int_t   iSHIL  =  1;
341     Int_t   iT0 =  1;
342     Int_t   iTOF   =  1;
343     Int_t   iTPC   =  1;
344     Int_t   iTRD   =  1;
345     Int_t   iZDC   =  1;
346     Int_t   iEMCAL =  1;
347     Int_t   iACORDE   =  0;
348     Int_t   iVZERO =  1;
349     rl->CdGAFile();
350     //=================== Alice BODY parameters =============================
351     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
352
353     if (iMAG)
354     {
355         //=================== MAG parameters ============================
356         // --- Start with Magnet since detector layouts may be depending ---
357         // --- on the selected Magnet dimensions ---
358         AliMAG *MAG = new AliMAG("MAG", "Magnet");
359     }
360
361
362     if (iABSO)
363     {
364         //=================== ABSO parameters ============================
365         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
366     }
367
368     if (iDIPO)
369     {
370         //=================== DIPO parameters ============================
371
372         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
373     }
374
375     if (iHALL)
376     {
377         //=================== HALL parameters ============================
378
379         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
380     }
381
382
383     if (iFRAME)
384     {
385         //=================== FRAME parameters ============================
386
387         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
388         if (geo == kHoles) {
389           FRAME->SetHoles(1);
390         } else {
391           FRAME->SetHoles(0);
392         }
393     }
394
395     if (iSHIL)
396     {
397         //=================== SHIL parameters ============================
398
399         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
400     }
401
402
403     if (iPIPE)
404     {
405         //=================== PIPE parameters ============================
406
407         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
408     }
409  
410     if (iITS)
411     {
412         //=================== ITS parameters ============================
413
414         AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
415     }
416
417     if (iTPC)
418     {
419         //============================ TPC parameters ===================
420         AliTPC *TPC = new AliTPCv2("TPC", "Default");
421     }
422
423
424     if (iTOF) {
425         //=================== TOF parameters ============================
426         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
427         // Partial geometry: modules at 2,3,4,6,7,11,12,14,15,16
428         // starting at 6h in positive direction
429         Int_t TOFSectors[18]={-1,-1,0,0,0,-1,0,0,-1,-1,-1,0,0,-1,0,0,0,0};
430         TOF->SetTOFSectors(TOFSectors);
431      }
432
433
434     if (iHMPID)
435     {
436         //=================== HMPID parameters ===========================
437         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
438
439     }
440
441
442     if (iZDC)
443     {
444         //=================== ZDC parameters ============================
445
446         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
447     }
448
449     if (iTRD)
450     {
451         //=================== TRD parameters ============================
452
453         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
454         AliTRDgeometry *geoTRD = TRD->GetGeometry();
455         // Partial geometry: modules at 2,3,4,6,11,12,14,15
456         // starting at 6h in positive direction
457         geoTRD->SetSMstatus(0,0);
458         geoTRD->SetSMstatus(1,0);
459         geoTRD->SetSMstatus(5,0);
460         geoTRD->SetSMstatus(7,0);
461         geoTRD->SetSMstatus(8,0);
462         geoTRD->SetSMstatus(9,0);
463         geoTRD->SetSMstatus(10,0);
464         geoTRD->SetSMstatus(13,0);
465         geoTRD->SetSMstatus(16,0);
466         geoTRD->SetSMstatus(17,0);
467     }
468
469     if (iFMD)
470     {
471         //=================== FMD parameters ============================
472         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
473    }
474
475     if (iMUON)
476     {
477         //=================== MUON parameters ===========================
478         // New MUONv1 version (geometry defined via builders)
479         AliMUON *MUON = new AliMUONv1("MUON","default");
480     }
481     //=================== PHOS parameters ===========================
482
483     if (iPHOS)
484     {
485         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
486     }
487
488
489     if (iPMD)
490     {
491         //=================== PMD parameters ============================
492         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
493     }
494
495     if (iT0)
496     {
497         //=================== T0 parameters ============================
498         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
499     }
500
501     if (iEMCAL)
502     {
503         //=================== EMCAL parameters ============================
504         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "SHISH_77_TRD1_2X2_FINAL_110DEG");
505     }
506
507      if (iACORDE)
508     {
509         //=================== ACORDE parameters ============================
510         AliACORDE *ACORDE = new AliACORDEv0("ACORDE", "normal ACORDE");
511     }
512
513      if (iVZERO)
514     {
515         //=================== VZERO parameters ============================
516         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
517     }
518
519
520 }
521
522 Float_t EtaToTheta(Float_t arg){
523   return (180./TMath::Pi())*2.*atan(exp(-arg));
524 }