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