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