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++")
7 #if !defined(__CINT__) || defined(__MAKECINT__)
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/AliABSOv3.h"
25 #include "STRUCT/AliDIPOv3.h"
26 #include "STRUCT/AliHALLv3.h"
27 #include "STRUCT/AliFRAMEv2.h"
28 #include "STRUCT/AliSHILv3.h"
29 #include "STRUCT/AliPIPEv3.h"
30 #include "ITS/AliITSvPPRasymmFMD.h"
31 #include "TPC/AliTPCv2.h"
32 #include "TOF/AliTOFv6T0.h"
33 #include "HMPID/AliHMPIDv2.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"
49 kDefaultPPTrig, kDefaultPbPbTrig
52 const char * pprTrigConfName[] = {
56 Float_t EtaToTheta(Float_t arg);
58 static PprTrigConf_t strig = kDefaultPPTrig;// default PP trigger configuration
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
65 // Set Random Number seed
66 gRandom->SetSeed(123456); // Set 0 to use the currecnt time
69 // libraries required by geant321
71 gSystem->Load("libgeant321");
74 new TGeant3TGeo("C++ Interface to Geant3");
79 rl = AliRunLoader::Open("galice.root",
80 AliConfig::GetDefaultEventFolderName(),
84 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
87 rl->SetCompressionLevel(2);
88 rl->SetNumberOfEventsPerFile(3);
89 gAlice->SetRunLoader(rl);
91 // Set the trigger configuration
92 gAlice->SetTriggerDescriptor(pprTrigConfName[strig]);
93 cout<<"Trigger configuration is set to "<<pprTrigConfName[strig]<<endl;
96 // Set External decayer
97 TVirtualMCDecayer *decayer = new AliDecayerPythia();
99 decayer->SetForceDecay(kAll);
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
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);
123 Float_t cut = 1.e-3; // 1MeV cut by default
124 Float_t tofmax = 1.e10;
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);
138 // Special generation for Valgrind tests
139 // Each detector is fired by few particles selected
140 // to cover specific cases
143 // The cocktail iitself
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
155 // Particle guns for the barrel part (taken from RichConfig)
157 AliGenFixed *pG1=new AliGenFixed(1);
159 pG1->SetMomentum(2.5);
160 pG1->SetTheta(109.5-3);
162 gener->AddGenerator(pG1,"g1",1);
164 AliGenFixed *pG2=new AliGenFixed(1);
166 pG2->SetMomentum(1.0);
167 pG2->SetTheta( 90.0-3);
169 gener->AddGenerator(pG2,"g2",1);
171 AliGenFixed *pG3=new AliGenFixed(1);
173 pG3->SetMomentum(1.5);
174 pG3->SetTheta(109.5-3);
176 gener->AddGenerator(pG3,"g3",1);
178 AliGenFixed *pG4=new AliGenFixed(1);
180 pG4->SetMomentum(0.7);
181 pG4->SetTheta( 90.0-3);
183 gener->AddGenerator(pG4,"g4",1);
185 AliGenFixed *pG5=new AliGenFixed(1);
187 pG5->SetMomentum(1.0);
188 pG5->SetTheta( 70.0-3);
190 gener->AddGenerator(pG5,"g5",1);
192 AliGenFixed *pG6=new AliGenFixed(1);
194 pG6->SetMomentum(2.5);
195 pG6->SetTheta( 90.0-3);
197 gener->AddGenerator(pG6,"g6",1);
199 AliGenFixed *pG7=new AliGenFixed(1);
201 pG7->SetMomentum(0.7);
202 pG7->SetTheta( 70.0-3);
204 gener->AddGenerator(pG7,"g7",1);
208 AliGenFixed *pG8=new AliGenFixed(1);
210 pG8->SetMomentum(1.2);
211 pG8->SetTheta( 95.0);
213 gener->AddGenerator(pG8,"g8",1);
215 AliGenFixed *pG9=new AliGenFixed(1);
217 pG9->SetMomentum(1.2);
218 pG9->SetTheta( 85.0);
220 gener->AddGenerator(pG9,"g9",1);
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);
229 gener->AddGenerator(gphos,"GENBOX GAMMA for PHOS",1);
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);
238 gener->AddGenerator(gemcal,"GENBOX GAMMA for EMCAL",1);
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);
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);
256 AliGenFixed *gtof=new AliGenFixed(1);
258 gtof->SetMomentum(2.5);
261 gener->AddGenerator(gtof,"Proton for TOF",1);
267 // Activate this line if you want the vertex smearing to happen
270 //gener->SetVertexSmear(perTrack);
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);
299 //=================== Alice BODY parameters =============================
300 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
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");
313 //=================== ABSO parameters ============================
314 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
319 //=================== DIPO parameters ============================
321 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
326 //=================== HALL parameters ============================
328 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
334 //=================== FRAME parameters ============================
336 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
341 //=================== SHIL parameters ============================
343 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
349 //=================== PIPE parameters ============================
351 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
356 //=================== ITS parameters ============================
358 AliITSvPPRasymmFMD *ITS = new AliITSvPPRasymmFMD("ITS","ITS PPR detailed version with asymmetric services");
363 //============================ TPC parameters ===================
364 AliTPC *TPC = new AliTPCv2("TPC", "Default");
369 //=================== TOF parameters ============================
370 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
371 // Partial geometry: modules at 2,3,4,6,7,11,12,14,15,16
372 // starting at 6h in positive direction
373 Int_t TOFSectors[18]={-1,-1,0,0,0,-1,0,0,-1,-1,-1,0,0,-1,0,0,0,0};
374 TOF->SetTOFSectors(TOFSectors);
380 //=================== HMPID parameters ===========================
381 AliHMPID *HMPID = new AliHMPIDv2("HMPID", "normal HMPID");
388 //=================== ZDC parameters ============================
390 AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
395 //=================== TRD parameters ============================
397 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
398 AliTRDgeometry *geoTRD = TRD->GetGeometry();
399 // Partial geometry: modules at 2,3,4,6,11,12,14,15
400 // starting at 6h in positive direction
401 geoTRD->SetSMstatus(0,0);
402 geoTRD->SetSMstatus(1,0);
403 geoTRD->SetSMstatus(5,0);
404 geoTRD->SetSMstatus(7,0);
405 geoTRD->SetSMstatus(8,0);
406 geoTRD->SetSMstatus(9,0);
407 geoTRD->SetSMstatus(10,0);
408 geoTRD->SetSMstatus(13,0);
409 geoTRD->SetSMstatus(16,0);
410 geoTRD->SetSMstatus(17,0);
415 //=================== FMD parameters ============================
416 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
421 //=================== MUON parameters ===========================
422 // New MUONv1 version (geometry defined via builders)
423 AliMUON *MUON = new AliMUONv1("MUON","default");
425 //=================== PHOS parameters ===========================
429 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
435 //=================== PMD parameters ============================
436 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
441 //=================== T0 parameters ============================
442 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
447 //=================== EMCAL parameters ============================
448 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "SHISH_77_TRD1_2X2_FINAL_110DEG");
453 //=================== ACORDE parameters ============================
454 AliACORDE *ACORDE = new AliACORDEv0("ACORDE", "normal ACORDE");
459 //=================== VZERO parameters ============================
460 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
466 Float_t EtaToTheta(Float_t arg){
467 return (180./TMath::Pi())*2.*atan(exp(-arg));