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__)
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/AliITSvPPRasymmFMD.h"
32 #include "TPC/AliTPCv2.h"
33 #include "TOF/AliTOFv6T0.h"
34 #include "HMPID/AliHMPIDv2.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"
50 kDefaultPPTrig, kDefaultPbPbTrig
53 const char * pprTrigConfName[] = {
57 Float_t EtaToTheta(Float_t arg);
59 static PprTrigConf_t strig = kDefaultPPTrig;// default PP trigger configuration
63 // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
64 // Theta range given through pseudorapidity limits 22/6/2001
66 // Set Random Number seed
67 gRandom->SetSeed(123456); // Set 0 to use the currecnt time
70 // libraries required by geant321
72 gSystem->Load("libgeant321");
75 new TGeant3TGeo("C++ Interface to Geant3");
80 rl = AliRunLoader::Open("galice.root",
81 AliConfig::GetDefaultEventFolderName(),
85 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
88 rl->SetCompressionLevel(2);
89 rl->SetNumberOfEventsPerFile(3);
90 gAlice->SetRunLoader(rl);
92 // Set the trigger configuration
93 gAlice->SetTriggerDescriptor(pprTrigConfName[strig]);
94 cout<<"Trigger configuration is set to "<<pprTrigConfName[strig]<<endl;
97 // Set External decayer
98 TVirtualMCDecayer *decayer = new AliDecayerPythia();
100 decayer->SetForceDecay(kAll);
102 gMC->SetExternalDecayer(decayer);
103 //=======================================================================
104 // ************* STEERING parameters FOR ALICE SIMULATION **************
105 // --- Specify event type to be tracked through the ALICE setup
106 // --- All positions are in cm, angles in degrees, and P and E in GeV
109 gMC->SetProcess("DCAY",1);
110 gMC->SetProcess("PAIR",1);
111 gMC->SetProcess("COMP",1);
112 gMC->SetProcess("PHOT",1);
113 gMC->SetProcess("PFIS",0);
114 gMC->SetProcess("DRAY",0);
115 gMC->SetProcess("ANNI",1);
116 gMC->SetProcess("BREM",1);
117 gMC->SetProcess("MUNU",1);
118 gMC->SetProcess("CKOV",1);
119 gMC->SetProcess("HADR",1);
120 gMC->SetProcess("LOSS",2);
121 gMC->SetProcess("MULS",1);
122 gMC->SetProcess("RAYL",1);
124 Float_t cut = 1.e-3; // 1MeV cut by default
125 Float_t tofmax = 1.e10;
127 gMC->SetCut("CUTGAM", cut);
128 gMC->SetCut("CUTELE", cut);
129 gMC->SetCut("CUTNEU", cut);
130 gMC->SetCut("CUTHAD", cut);
131 gMC->SetCut("CUTMUO", cut);
132 gMC->SetCut("BCUTE", cut);
133 gMC->SetCut("BCUTM", cut);
134 gMC->SetCut("DCUTE", cut);
135 gMC->SetCut("DCUTM", cut);
136 gMC->SetCut("PPCUTM", cut);
137 gMC->SetCut("TOFMAX", tofmax);
139 // Special generation for Valgrind tests
140 // Each detector is fired by few particles selected
141 // to cover specific cases
144 // The cocktail iitself
146 AliGenCocktail *gener = new AliGenCocktail();
147 gener->SetPhiRange(0, 360);
148 // Set pseudorapidity range from -8 to 8.
149 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
150 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
151 gener->SetThetaRange(thmin,thmax);
152 gener->SetOrigin(0, 0, 0); //vertex position
153 gener->SetSigma(0, 0, 0); //Sigma in (X,Y,Z) (cm) on IP position
156 // Particle guns for the barrel part (taken from RichConfig)
158 AliGenFixed *pG1=new AliGenFixed(1);
159 pG1->SetPart(kProton);
160 pG1->SetMomentum(2.5);
161 pG1->SetTheta(109.5-3);
163 gener->AddGenerator(pG1,"g1",1);
165 AliGenFixed *pG2=new AliGenFixed(1);
166 pG2->SetPart(kPiPlus);
167 pG2->SetMomentum(1.0);
168 pG2->SetTheta( 90.0-3);
170 gener->AddGenerator(pG2,"g2",1);
172 AliGenFixed *pG3=new AliGenFixed(1);
173 pG3->SetPart(kPiMinus);
174 pG3->SetMomentum(1.5);
175 pG3->SetTheta(109.5-3);
177 gener->AddGenerator(pG3,"g3",1);
179 AliGenFixed *pG4=new AliGenFixed(1);
180 pG4->SetPart(kKPlus);
181 pG4->SetMomentum(0.7);
182 pG4->SetTheta( 90.0-3);
184 gener->AddGenerator(pG4,"g4",1);
186 AliGenFixed *pG5=new AliGenFixed(1);
187 pG5->SetPart(kKMinus);
188 pG5->SetMomentum(1.0);
189 pG5->SetTheta( 70.0-3);
191 gener->AddGenerator(pG5,"g5",1);
193 AliGenFixed *pG6=new AliGenFixed(1);
194 pG6->SetPart(kProtonBar);
195 pG6->SetMomentum(2.5);
196 pG6->SetTheta( 90.0-3);
198 gener->AddGenerator(pG6,"g6",1);
200 AliGenFixed *pG7=new AliGenFixed(1);
201 pG7->SetPart(kPiMinus);
202 pG7->SetMomentum(0.7);
203 pG7->SetTheta( 70.0-3);
205 gener->AddGenerator(pG7,"g7",1);
209 AliGenFixed *pG8=new AliGenFixed(1);
210 pG8->SetPart(kElectron);
211 pG8->SetMomentum(1.2);
212 pG8->SetTheta( 95.0);
214 gener->AddGenerator(pG8,"g8",1);
216 AliGenFixed *pG9=new AliGenFixed(1);
217 pG9->SetPart(kPositron);
218 pG9->SetMomentum(1.2);
219 pG9->SetTheta( 85.0);
221 gener->AddGenerator(pG9,"g9",1);
225 AliGenBox *gphos = new AliGenBox(1);
226 gphos->SetMomentumRange(10,11.);
227 gphos->SetPhiRange(270.5,270.7);
228 gphos->SetThetaRange(90.5,90.7);
229 gphos->SetPart(kGamma);
230 gener->AddGenerator(gphos,"GENBOX GAMMA for PHOS",1);
234 AliGenBox *gemcal = new AliGenBox(1);
235 gemcal->SetMomentumRange(10,11.);
236 gemcal->SetPhiRange(90.5,199.5);
237 gemcal->SetThetaRange(90.5,90.7);
238 gemcal->SetPart(kGamma);
239 gener->AddGenerator(gemcal,"GENBOX GAMMA for EMCAL",1);
242 AliGenBox * gmuon1 = new AliGenBox(1);
243 gmuon1->SetMomentumRange(20.,20.1);
244 gmuon1->SetPhiRange(0., 360.);
245 gmuon1->SetThetaRange(171.000,178.001);
246 gmuon1->SetPart(kMuonMinus); // Muons
247 gener->AddGenerator(gmuon1,"GENBOX MUON1",1);
249 AliGenBox * gmuon2 = new AliGenBox(1);
250 gmuon2->SetMomentumRange(20.,20.1);
251 gmuon2->SetPhiRange(0., 360.);
252 gmuon2->SetThetaRange(171.000,178.001);
253 gmuon2->SetPart(kMuonPlus); // Muons
254 gener->AddGenerator(gmuon2,"GENBOX MUON1",1);
257 AliGenFixed *gtof=new AliGenFixed(1);
258 gtof->SetPart(kProton);
259 gtof->SetMomentum(2.5);
262 gener->AddGenerator(gtof,"Proton for TOF",1);
268 // Activate this line if you want the vertex smearing to happen
271 //gener->SetVertexSmear(perTrack);
273 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 2);
274 field->SetL3ConstField(1); //Using const. field in the barrel if 0
275 gAlice->SetField(field);
300 //=================== Alice BODY parameters =============================
301 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
305 //=================== MAG parameters ============================
306 // --- Start with Magnet since detector layouts may be depending ---
307 // --- on the selected Magnet dimensions ---
308 AliMAG *MAG = new AliMAG("MAG", "Magnet");
314 //=================== ABSO parameters ============================
315 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
320 //=================== DIPO parameters ============================
322 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
327 //=================== HALL parameters ============================
329 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
335 //=================== FRAME parameters ============================
337 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
342 //=================== SHIL parameters ============================
344 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
350 //=================== PIPE parameters ============================
352 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
357 //=================== ITS parameters ============================
359 AliITSvPPRasymmFMD *ITS = new AliITSvPPRasymmFMD("ITS","ITS PPR detailed version with asymmetric services");
364 //============================ TPC parameters ===================
365 AliTPC *TPC = new AliTPCv2("TPC", "Default");
370 //=================== TOF parameters ============================
371 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
372 // Partial geometry: modules at 2,3,4,6,7,11,12,14,15,16
373 // starting at 6h in positive direction
374 Int_t TOFSectors[18]={-1,-1,0,0,0,-1,0,0,-1,-1,-1,0,0,-1,0,0,0,0};
375 TOF->SetTOFSectors(TOFSectors);
381 //=================== HMPID parameters ===========================
382 AliHMPID *HMPID = new AliHMPIDv2("HMPID", "normal HMPID");
389 //=================== ZDC parameters ============================
391 AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
396 //=================== TRD parameters ============================
398 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
399 AliTRDgeometry *geoTRD = TRD->GetGeometry();
400 // Partial geometry: modules at 2,3,4,6,11,12,14,15
401 // starting at 6h in positive direction
402 geoTRD->SetSMstatus(0,0);
403 geoTRD->SetSMstatus(1,0);
404 geoTRD->SetSMstatus(5,0);
405 geoTRD->SetSMstatus(7,0);
406 geoTRD->SetSMstatus(8,0);
407 geoTRD->SetSMstatus(9,0);
408 geoTRD->SetSMstatus(10,0);
409 geoTRD->SetSMstatus(13,0);
410 geoTRD->SetSMstatus(16,0);
411 geoTRD->SetSMstatus(17,0);
416 //=================== FMD parameters ============================
417 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
422 //=================== MUON parameters ===========================
423 // New MUONv1 version (geometry defined via builders)
424 AliMUON *MUON = new AliMUONv1("MUON","default");
426 //=================== PHOS parameters ===========================
430 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
436 //=================== PMD parameters ============================
437 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
442 //=================== T0 parameters ============================
443 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
448 //=================== EMCAL parameters ============================
449 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "SHISH_77_TRD1_2X2_FINAL_110DEG");
454 //=================== ACORDE parameters ============================
455 AliACORDE *ACORDE = new AliACORDEv0("ACORDE", "normal ACORDE");
460 //=================== VZERO parameters ============================
461 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
467 Float_t EtaToTheta(Float_t arg){
468 return (180./TMath::Pi())*2.*atan(exp(-arg));