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