]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - macros/Config.C
Updated detector and mag. field configuration in sync with ConfigPPR.C
[u/mrichter/AliRoot.git] / macros / Config.C
... / ...
CommitLineData
1static Int_t eventsPerRun = 100;
2void Config()
3{
4 // 7-DEC-2000 09:00
5 // Switch on Transition Radiation simulation. 6/12/00 18:00
6 // iZDC=1 7/12/00 09:00
7 // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
8 // Theta range given through pseudorapidity limits 22/6/2001
9
10 // Set Random Number seed
11 // gRandom->SetSeed(12345);
12
13 new AliGeant3("C++ Interface to Geant3");
14
15 if (!gSystem->Getenv("CONFIG_FILE"))
16 {
17 TFile *rootfile = new TFile("galice.root", "recreate");
18
19 rootfile->SetCompressionLevel(2);
20 }
21
22 TGeant3 *geant3 = (TGeant3 *) gMC;
23
24 //
25 // Set External decayer
26 AliDecayer *decayer = new AliDecayerPythia();
27
28 decayer->SetForceDecay(kAll);
29 decayer->Init();
30 gMC->SetExternalDecayer(decayer);
31 //
32 //
33 //=======================================================================
34 // ******* GEANT STEERING parameters FOR ALICE SIMULATION *******
35 geant3->SetTRIG(1); //Number of events to be processed
36 geant3->SetSWIT(4, 10);
37 geant3->SetDEBU(0, 0, 1);
38 //geant3->SetSWIT(2,2);
39 geant3->SetDCAY(1);
40 geant3->SetPAIR(1);
41 geant3->SetCOMP(1);
42 geant3->SetPHOT(1);
43 geant3->SetPFIS(0);
44 geant3->SetDRAY(0);
45 geant3->SetANNI(1);
46 geant3->SetBREM(1);
47 geant3->SetMUNU(1);
48 geant3->SetCKOV(1);
49 geant3->SetHADR(1); //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3)
50 geant3->SetLOSS(2);
51 geant3->SetMULS(1);
52 geant3->SetRAYL(1);
53 geant3->SetAUTO(1); //Select automatic STMIN etc... calc. (AUTO 1) or manual (AUTO 0)
54 geant3->SetABAN(0); //Restore 3.16 behaviour for abandoned tracks
55 geant3->SetOPTI(2); //Select optimisation level for GEANT geometry searches (0,1,2)
56 geant3->SetERAN(5.e-7);
57
58 Float_t cut = 1.e-3; // 1MeV cut by default
59 Float_t tofmax = 1.e10;
60
61 // GAM ELEC NHAD CHAD MUON EBREM MUHAB EDEL MUDEL MUPA TOFMAX
62 geant3->SetCUTS(cut, cut, cut, cut, cut, cut, cut, cut, cut, cut,
63 tofmax);
64 //
65 //=======================================================================
66 // ************* STEERING parameters FOR ALICE SIMULATION **************
67 // --- Specify event type to be tracked through the ALICE setup
68 // --- All positions are in cm, angles in degrees, and P and E in GeV
69 if (gSystem->Getenv("CONFIG_NPARTICLES"))
70 {
71 int nParticles = atoi(gSystem->Getenv("CONFIG_NPARTICLES"));
72 } else
73 {
74 int nParticles = 50;
75 }
76 AliGenHIJINGpara *gener = new AliGenHIJINGpara(nParticles);
77
78 gener->SetMomentumRange(0, 999);
79 gener->SetPhiRange(0, 360);
80 // Set pseudorapidity range from -8 to 8.
81 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
82 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
83 gener->SetThetaRange(thmin,thmax);
84 gener->SetOrigin(0, 0, 0); //vertex position
85 gener->SetSigma(0, 0, 0); //Sigma in (X,Y,Z) (cm) on IP position
86 gener->Init();
87 //
88 // Activate this line if you want the vertex smearing to happen
89 // track by track
90 //
91 //gener->SetVertexSmear(perTrack);
92 // Field (L3 0.4 T)
93 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
94 rootfile->cd();
95 gAlice->SetField(field);
96
97
98 Int_t iABSO = 1;
99 Int_t iDIPO = 1;
100 Int_t iFMD = 1;
101 Int_t iFRAME = 1;
102 Int_t iHALL = 1;
103 Int_t iITS = 1;
104 Int_t iMAG = 1;
105 Int_t iMUON = 1;
106 Int_t iPHOS = 1;
107 Int_t iPIPE = 1;
108 Int_t iPMD = 1;
109 Int_t iRICH = 1;
110 Int_t iSHIL = 1;
111 Int_t iSTART = 1;
112 Int_t iTOF = 1;
113 Int_t iTPC = 1;
114 Int_t iTRD = 1;
115 Int_t iZDC = 1;
116 Int_t iEMCAL = 1;
117 Int_t iCRT = 0;
118
119 //=================== Alice BODY parameters =============================
120 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
121
122 if (iMAG)
123 {
124 //=================== MAG parameters ============================
125 // --- Start with Magnet since detector layouts may be depending ---
126 // --- on the selected Magnet dimensions ---
127 AliMAG *MAG = new AliMAG("MAG", "Magnet");
128 }
129
130
131 if (iABSO)
132 {
133 //=================== ABSO parameters ============================
134 AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
135 }
136
137 if (iDIPO)
138 {
139 //=================== DIPO parameters ============================
140
141 AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
142 }
143
144 if (iHALL)
145 {
146 //=================== HALL parameters ============================
147
148 AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
149 }
150
151
152 if (iFRAME)
153 {
154 //=================== FRAME parameters ============================
155
156 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
157 if (geo == kHoles) {
158 FRAME->SetHoles(1);
159 } else {
160 FRAME->SetHoles(0);
161 }
162 }
163
164 if (iSHIL)
165 {
166 //=================== SHIL parameters ============================
167
168 AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
169 }
170
171
172 if (iPIPE)
173 {
174 //=================== PIPE parameters ============================
175
176 AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
177 }
178
179 if(iITS) {
180
181 //=================== ITS parameters ============================
182 //
183 // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
184 // almost all other detectors. This involves the fact that the ITS geometry
185 // still has several options to be followed in parallel in order to determine
186 // the best set-up which minimizes the induced background. All the geometries
187 // available to date are described in the following. Read carefully the comments
188 // and use the default version (the only one uncommented) unless you are making
189 // comparisons and you know what you are doing. In this case just uncomment the
190 // ITS geometry you want to use and run Aliroot.
191 //
192 // Detailed geometries:
193 //
194 //
195 //AliITS *ITS = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
196 //
197 //AliITS *ITS = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
198 //
199 AliITSvPPRasymm *ITS = new AliITSvPPRasymm("ITS","New ITS PPR detailed version with asymmetric services");
200 ITS->SetMinorVersion(2); // don't touch this parameter if you're not an ITS developer
201 ITS->SetReadDet(kFALSE); // don't touch this parameter if you're not an ITS developer
202 // ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det"); // don't touch this parameter if you're not an ITS developer
203 ITS->SetThicknessDet1(200.); // detector thickness on layer 1 must be in the range [100,300]
204 ITS->SetThicknessDet2(200.); // detector thickness on layer 2 must be in the range [100,300]
205 ITS->SetThicknessChip1(200.); // chip thickness on layer 1 must be in the range [150,300]
206 ITS->SetThicknessChip2(200.); // chip thickness on layer 2 must be in the range [150,300]
207 ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out
208 ITS->SetCoolingFluid(1); // 1 --> water ; 0 --> freon
209 //
210 //AliITSvPPRsymm *ITS = new AliITSvPPRsymm("ITS","New ITS PPR detailed version with symmetric services");
211 //ITS->SetMinorVersion(2); // don't touch this parameter if you're not an ITS developer
212 //ITS->SetReadDet(kFALSE); // don't touch this parameter if you're not an ITS developer
213 //ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRsymm2.det"); // don't touch this parameter if you're not an ITS developer
214 //ITS->SetThicknessDet1(200.); // detector thickness on layer 1 must be in the range [100,300]
215 //ITS->SetThicknessDet2(200.); // detector thickness on layer 2 must be in the range [100,300]
216 //ITS->SetThicknessChip1(200.); // chip thickness on layer 1 must be in the range [150,300]
217 //ITS->SetThicknessChip2(200.); // chip thickness on layer 2 must be in the range [150,300]
218 //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out
219 //ITS->SetCoolingFluid(1); // 1 --> water ; 0 --> freon
220 //
221 //
222 // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful
223 // for reconstruction !):
224 //
225 //
226 //AliITSvPPRcoarseasymm *ITS = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
227 //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out
228 //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
229 //
230 //AliITS *ITS = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
231 //ITS->SetRails(0); // 1 --> rails in ; 0 --> rails out
232 //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
233 //
234 //
235 //
236 // Geant3 <-> EUCLID conversion
237 // ============================
238 //
239 // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
240 // media to two ASCII files (called by default ITSgeometry.euc and
241 // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
242 // The default (=0) means that you dont want to use this facility.
243 //
244 ITS->SetEUCLID(0);
245 }
246
247 if (iTPC)
248 {
249 //============================ TPC parameters ================================
250 // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
251 // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
252 // --- sectors are specified, any value other than that requires at least one
253 // --- sector (lower or upper)to be specified!
254 // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
255 // --- sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
256 // --- SecLows - number of lower sectors specified (up to 6)
257 // --- SecUps - number of upper sectors specified (up to 12)
258 // --- Sens - sensitive strips for the Slow Simulator !!!
259 // --- This does NOT work if all S or L-sectors are specified, i.e.
260 // --- if SecAL or SecAU < 0
261 //
262 //
263 //-----------------------------------------------------------------------------
264
265 // gROOT->LoadMacro("SetTPCParam.C");
266 // AliTPCParam *param = SetTPCParam();
267 AliTPC *TPC = new AliTPCv2("TPC", "Default");
268
269 // All sectors included
270 TPC->SetSecAL(-1);
271 TPC->SetSecAU(-1);
272
273 }
274
275
276 if (iTOF) {
277 if (geo == kHoles) {
278 //=================== TOF parameters ============================
279 AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
280 } else {
281 AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
282 }
283 }
284
285
286 if (iRICH)
287 {
288 //=================== RICH parameters ===========================
289 AliRICH *RICH = new AliRICHv3("RICH", "normal RICH");
290
291 }
292
293
294 if (iZDC)
295 {
296 //=================== ZDC parameters ============================
297
298 AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
299 }
300
301 if (iTRD)
302 {
303 //=================== TRD parameters ============================
304
305 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
306
307 // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
308 TRD->SetGasMix(1);
309 if (geo == kHoles) {
310 // With hole in front of PHOS
311 TRD->SetPHOShole();
312 // With hole in front of RICH
313 TRD->SetRICHhole();
314 }
315 // Switch on TR
316 AliTRDsim *TRDsim = TRD->CreateTR();
317 }
318
319 if (iFMD)
320 {
321 //=================== FMD parameters ============================
322 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
323 FMD->SetRingsSi1(256);
324 FMD->SetRingsSi2(128);
325 FMD->SetSectorsSi1(20);
326 FMD->SetSectorsSi2(40);
327 }
328
329 if (iMUON)
330 {
331 //=================== MUON parameters ===========================
332
333 AliMUON *MUON = new AliMUONv1("MUON", "default");
334 }
335 //=================== PHOS parameters ===========================
336
337 if (iPHOS)
338 {
339 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
340 }
341
342
343 if (iPMD)
344 {
345 //=================== PMD parameters ============================
346 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
347 }
348
349 if (iSTART)
350 {
351 //=================== START parameters ============================
352 AliSTART *START = new AliSTARTv1("START", "START Detector");
353 }
354
355 if (iEMCAL)
356 {
357 //=================== EMCAL parameters ============================
358 AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCALArch1a");
359 }
360
361 if (iCRT)
362 {
363 //=================== CRT parameters ============================
364 AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
365 }
366
367 if (iVZERO)
368 {
369 //=================== CRT parameters ============================
370 AliVZERO *VZERO = new AliVZEROv2("VZERO", "normal VZERO");
371 }
372
373}
374
375Float_t EtaToTheta(Float_t arg){
376 return (180./TMath::Pi())*2.*atan(exp(-arg));
377}