28b75255 |
1 | // |
2 | // Configuration for the first physics production 2008 |
3 | // |
4 | |
5 | // One can use the configuration macro in compiled mode by |
6 | // root [0] gSystem->Load("libgeant321"); |
7 | // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\ |
8 | // -I$ALICE_ROOT -I$ALICE/geant3/TGeant3"); |
9 | // root [0] .x grun.C(1,"Config.C++") |
10 | |
11 | #if !defined(__CINT__) || defined(__MAKECINT__) |
12 | #include <Riostream.h> |
13 | #include <TRandom.h> |
14 | #include <TDatime.h> |
15 | #include <TSystem.h> |
16 | #include <TVirtualMC.h> |
17 | #include <TGeant3TGeo.h> |
18 | #include "STEER/AliRunLoader.h" |
19 | #include "STEER/AliRun.h" |
20 | #include "STEER/AliConfig.h" |
21 | #include "PYTHIA6/AliDecayerPythia.h" |
22 | #include "PYTHIA6/AliGenPythia.h" |
23 | #include "TDPMjet/AliGenDPMjet.h" |
24 | #include "STEER/AliMagWrapCheb.h" |
25 | #include "STRUCT/AliBODY.h" |
26 | #include "STRUCT/AliMAG.h" |
27 | #include "STRUCT/AliABSOv3.h" |
28 | #include "STRUCT/AliDIPOv3.h" |
29 | #include "STRUCT/AliHALLv3.h" |
30 | #include "STRUCT/AliFRAMEv2.h" |
31 | #include "STRUCT/AliSHILv3.h" |
32 | #include "STRUCT/AliPIPEv3.h" |
33 | #include "ITS/AliITSv11Hybrid.h" |
34 | #include "TPC/AliTPCv2.h" |
35 | #include "TOF/AliTOFv6T0.h" |
36 | #include "HMPID/AliHMPIDv3.h" |
37 | #include "ZDC/AliZDCv3.h" |
38 | #include "TRD/AliTRDv1.h" |
39 | #include "TRD/AliTRDgeometry.h" |
40 | #include "FMD/AliFMDv1.h" |
41 | #include "MUON/AliMUONv1.h" |
42 | #include "PHOS/AliPHOSv1.h" |
aad1175b |
43 | #include "PHOS/AliPHOSSimParam.h" |
28b75255 |
44 | #include "PMD/AliPMDv1.h" |
45 | #include "T0/AliT0v1.h" |
46 | #include "EMCAL/AliEMCALv2.h" |
47 | #include "ACORDE/AliACORDEv1.h" |
48 | #include "VZERO/AliVZEROv7.h" |
49 | #endif |
50 | |
51 | |
52 | enum PDC06Proc_t |
53 | { |
54 | kPythia6, kPhojet, kRunMax |
55 | }; |
56 | |
57 | const char * pprRunName[] = { |
58 | "kPythia6", "kPhojet" |
59 | }; |
60 | |
61 | //--- Magnetic Field --- |
62 | enum Mag_t |
63 | { |
64 | kNoField, k5kG, kFieldMax |
65 | }; |
66 | |
67 | const char * pprField[] = { |
68 | "kNoField", "k5kG" |
69 | }; |
70 | |
71 | //--- Functions --- |
72 | class AliGenPythia; |
73 | AliGenerator *MbPythia(); |
74 | AliGenerator *MbPhojet(); |
75 | void ProcessEnvironmentVars(); |
76 | |
77 | // Geterator, field, beam energy |
78 | static PDC06Proc_t proc = kPhojet; |
79 | static Mag_t mag = k5kG; |
80 | static Float_t energy = 10000; // energy in CMS |
81 | //========================// |
82 | // Set Random Number seed // |
83 | //========================// |
84 | TDatime dt; |
85 | static UInt_t seed = dt.Get(); |
86 | |
87 | // Comment line |
88 | static TString comment; |
89 | |
90 | void Config() |
91 | { |
92 | |
93 | |
94 | // Get settings from environment variables |
95 | ProcessEnvironmentVars(); |
96 | |
97 | gRandom->SetSeed(seed); |
98 | cerr<<"Seed for random number generation= "<<seed<<endl; |
99 | |
100 | // Libraries required by geant321 |
101 | #if defined(__CINT__) |
102 | gSystem->Load("liblhapdf"); // Parton density functions |
103 | gSystem->Load("libEGPythia6"); // TGenerator interface |
104 | gSystem->Load("libpythia6"); // Pythia |
105 | gSystem->Load("libAliPythia6"); // ALICE specific implementations |
106 | gSystem->Load("libgeant321"); |
107 | #endif |
108 | |
109 | new TGeant3TGeo("C++ Interface to Geant3"); |
110 | |
111 | //======================================================================= |
112 | // Create the output file |
113 | |
114 | |
115 | AliRunLoader* rl=0x0; |
116 | |
117 | cout<<"Config.C: Creating Run Loader ..."<<endl; |
118 | rl = AliRunLoader::Open("galice.root", |
119 | AliConfig::GetDefaultEventFolderName(), |
120 | "recreate"); |
121 | if (rl == 0x0) |
122 | { |
123 | gAlice->Fatal("Config.C","Can not instatiate the Run Loader"); |
124 | return; |
125 | } |
126 | rl->SetCompressionLevel(2); |
127 | rl->SetNumberOfEventsPerFile(1000); |
128 | gAlice->SetRunLoader(rl); |
129 | // gAlice->SetGeometryFromFile("geometry.root"); |
130 | // gAlice->SetGeometryFromCDB(); |
131 | |
132 | // Set the trigger configuration: proton-proton |
133 | gAlice->SetTriggerDescriptor("p-p"); |
134 | |
135 | // |
136 | //======================================================================= |
137 | // ************* STEERING parameters FOR ALICE SIMULATION ************** |
138 | // --- Specify event type to be tracked through the ALICE setup |
139 | // --- All positions are in cm, angles in degrees, and P and E in GeV |
140 | |
141 | |
142 | gMC->SetProcess("DCAY",1); |
143 | gMC->SetProcess("PAIR",1); |
144 | gMC->SetProcess("COMP",1); |
145 | gMC->SetProcess("PHOT",1); |
146 | gMC->SetProcess("PFIS",0); |
147 | gMC->SetProcess("DRAY",0); |
148 | gMC->SetProcess("ANNI",1); |
149 | gMC->SetProcess("BREM",1); |
150 | gMC->SetProcess("MUNU",1); |
151 | gMC->SetProcess("CKOV",1); |
152 | gMC->SetProcess("HADR",1); |
153 | gMC->SetProcess("LOSS",2); |
154 | gMC->SetProcess("MULS",1); |
155 | gMC->SetProcess("RAYL",1); |
156 | |
157 | Float_t cut = 1.e-3; // 1MeV cut by default |
158 | Float_t tofmax = 1.e10; |
159 | |
160 | gMC->SetCut("CUTGAM", cut); |
161 | gMC->SetCut("CUTELE", cut); |
162 | gMC->SetCut("CUTNEU", cut); |
163 | gMC->SetCut("CUTHAD", cut); |
164 | gMC->SetCut("CUTMUO", cut); |
165 | gMC->SetCut("BCUTE", cut); |
166 | gMC->SetCut("BCUTM", cut); |
167 | gMC->SetCut("DCUTE", cut); |
168 | gMC->SetCut("DCUTM", cut); |
169 | gMC->SetCut("PPCUTM", cut); |
170 | gMC->SetCut("TOFMAX", tofmax); |
171 | |
172 | |
173 | |
174 | |
175 | //======================// |
176 | // Set External decayer // |
177 | //======================// |
178 | TVirtualMCDecayer* decayer = new AliDecayerPythia(); |
179 | decayer->SetForceDecay(kAll); |
180 | decayer->Init(); |
181 | gMC->SetExternalDecayer(decayer); |
182 | |
183 | //=========================// |
184 | // Generator Configuration // |
185 | //=========================// |
186 | AliGenerator* gener = 0x0; |
187 | |
188 | if (proc == kPythia6) { |
189 | gener = MbPythia(); |
190 | } else if (proc == kPhojet) { |
191 | gener = MbPhojet(); |
192 | } |
193 | |
194 | |
195 | |
196 | // PRIMARY VERTEX |
197 | // |
198 | gener->SetOrigin(0., 0., 0.); // vertex position |
199 | // |
200 | // |
201 | // Size of the interaction diamond |
202 | // Longitudinal |
aad1175b |
203 | Float_t sigmaz = 5.4 / TMath::Sqrt(2.); // [cm] |
204 | if (energy == 900) |
205 | sigmaz = 10.5 / TMath::Sqrt(2.); // [cm] |
28b75255 |
206 | // |
207 | // Transverse |
208 | Float_t betast = 10; // beta* [m] |
209 | Float_t eps = 3.75e-6; // emittance [m] |
210 | Float_t gamma = energy / 2.0 / 0.938272; // relativistic gamma [1] |
211 | Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.; // [cm] |
212 | printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz); |
213 | |
214 | gener->SetSigma(sigmaxy, sigmaxy, sigmaz); // Sigma in (X,Y,Z) (cm) on IP position |
215 | gener->SetCutVertexZ(3.); // Truncate at 3 sigma |
216 | gener->SetVertexSmear(kPerEvent); |
217 | |
218 | gener->Init(); |
219 | |
220 | // FIELD |
221 | // |
222 | AliMagWrapCheb* field = 0x0; |
223 | if (mag == kNoField) { |
224 | comment = comment.Append(" | L3 field 0.0 T"); |
225 | field = new AliMagWrapCheb("Maps","Maps", 2, 0., 10., AliMagWrapCheb::k2kG); |
226 | } else if (mag == k5kG) { |
227 | comment = comment.Append(" | L3 field 0.5 T"); |
228 | field = new AliMagWrapCheb("Maps","Maps", 2, 1., 10., AliMagWrapCheb::k5kG); |
229 | } |
230 | printf("\n \n Comment: %s \n \n", comment.Data()); |
231 | |
232 | rl->CdGAFile(); |
233 | gAlice->SetField(field); |
234 | |
235 | |
236 | |
237 | Int_t iABSO = 1; |
238 | Int_t iACORDE= 0; |
239 | Int_t iDIPO = 1; |
240 | Int_t iEMCAL = 0; |
241 | Int_t iFMD = 1; |
242 | Int_t iFRAME = 1; |
243 | Int_t iHALL = 1; |
244 | Int_t iITS = 1; |
245 | Int_t iMAG = 1; |
246 | Int_t iMUON = 1; |
247 | Int_t iPHOS = 1; |
248 | Int_t iPIPE = 1; |
249 | Int_t iPMD = 0; |
250 | Int_t iHMPID = 1; |
251 | Int_t iSHIL = 1; |
252 | Int_t iT0 = 1; |
253 | Int_t iTOF = 1; |
254 | Int_t iTPC = 1; |
255 | Int_t iTRD = 1; |
256 | Int_t iVZERO = 1; |
257 | Int_t iZDC = 1; |
258 | |
259 | |
260 | //=================== Alice BODY parameters ============================= |
261 | AliBODY *BODY = new AliBODY("BODY", "Alice envelop"); |
262 | |
263 | |
264 | if (iMAG) |
265 | { |
266 | //=================== MAG parameters ============================ |
267 | // --- Start with Magnet since detector layouts may be depending --- |
268 | // --- on the selected Magnet dimensions --- |
269 | AliMAG *MAG = new AliMAG("MAG", "Magnet"); |
270 | } |
271 | |
272 | |
273 | if (iABSO) |
274 | { |
275 | //=================== ABSO parameters ============================ |
276 | AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber"); |
277 | } |
278 | |
279 | if (iDIPO) |
280 | { |
281 | //=================== DIPO parameters ============================ |
282 | |
283 | AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3"); |
284 | } |
285 | |
286 | if (iHALL) |
287 | { |
288 | //=================== HALL parameters ============================ |
289 | |
290 | AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall"); |
291 | } |
292 | |
293 | |
294 | if (iFRAME) |
295 | { |
296 | //=================== FRAME parameters ============================ |
297 | |
298 | AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame"); |
299 | FRAME->SetHoles(1); |
300 | } |
301 | |
302 | if (iSHIL) |
303 | { |
304 | //=================== SHIL parameters ============================ |
305 | |
306 | AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3"); |
307 | } |
308 | |
309 | |
310 | if (iPIPE) |
311 | { |
312 | //=================== PIPE parameters ============================ |
313 | |
314 | AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe"); |
315 | } |
316 | |
317 | if (iITS) |
318 | { |
319 | //=================== ITS parameters ============================ |
320 | |
321 | AliITS *ITS = new AliITSv11Hybrid("ITS","ITS v11Hybrid"); |
322 | } |
323 | |
324 | if (iTPC) |
325 | { |
326 | //============================ TPC parameters ===================== |
327 | |
328 | AliTPC *TPC = new AliTPCv2("TPC", "Default"); |
329 | } |
330 | |
331 | |
332 | if (iTOF) { |
333 | //=================== TOF parameters ============================ |
334 | |
335 | AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF"); |
336 | } |
337 | |
338 | |
339 | if (iHMPID) |
340 | { |
341 | //=================== HMPID parameters =========================== |
342 | |
343 | AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID"); |
344 | |
345 | } |
346 | |
347 | |
348 | if (iZDC) |
349 | { |
350 | //=================== ZDC parameters ============================ |
351 | |
352 | AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC"); |
353 | } |
354 | |
355 | if (iTRD) |
356 | { |
357 | //=================== TRD parameters ============================ |
358 | |
359 | AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator"); |
360 | AliTRDgeometry *geoTRD = TRD->GetGeometry(); |
361 | // Partial geometry: modules at 0,8,9,17 |
362 | // starting at 3h in positive direction |
363 | geoTRD->SetSMstatus(1,0); |
364 | geoTRD->SetSMstatus(2,0); |
365 | geoTRD->SetSMstatus(3,0); |
366 | geoTRD->SetSMstatus(4,0); |
367 | geoTRD->SetSMstatus(5,0); |
368 | geoTRD->SetSMstatus(6,0); |
369 | geoTRD->SetSMstatus(7,0); |
370 | geoTRD->SetSMstatus(10,0); |
371 | geoTRD->SetSMstatus(11,0); |
372 | geoTRD->SetSMstatus(12,0); |
373 | geoTRD->SetSMstatus(13,0); |
374 | geoTRD->SetSMstatus(14,0); |
375 | geoTRD->SetSMstatus(15,0); |
376 | geoTRD->SetSMstatus(16,0); |
377 | } |
378 | |
379 | if (iFMD) |
380 | { |
381 | //=================== FMD parameters ============================ |
382 | |
383 | AliFMD *FMD = new AliFMDv1("FMD", "normal FMD"); |
384 | } |
385 | |
386 | if (iMUON) |
387 | { |
388 | //=================== MUON parameters =========================== |
389 | // New MUONv1 version (geometry defined via builders) |
390 | |
391 | AliMUON *MUON = new AliMUONv1("MUON", "default"); |
392 | } |
393 | |
394 | if (iPHOS) |
395 | { |
396 | //=================== PHOS parameters =========================== |
397 | |
398 | AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP"); |
aad1175b |
399 | //Set simulation parameters different from the default ones. |
400 | AliPHOSSimParam* simEmc = AliPHOSSimParam::GetInstance() ; |
401 | |
402 | // APD noise of warm (+20C) PHOS: |
403 | // a2 = a1*(Y1/Y2)*(M1/M2), where a1 = 0.012 is APD noise at -25C, |
404 | // Y1 = 4.3 photo-electrons/MeV, Y2 = 1.7 p.e/MeV - light yields at -25C and +20C, |
405 | // M1 = 50, M2 = 50 - APD gain factors chosen for t1 = -25C and t2 = +20C, |
406 | // Y = MeanLightYield*APDEfficiency. |
407 | |
408 | Float_t apdNoise = 0.012*2.5; |
409 | simEmc->SetAPDNoise(apdNoise); |
410 | |
411 | //Raw Light Yield at +20C |
412 | simEmc->SetMeanLightYield(18800); |
413 | |
414 | //ADC channel width at +18C. |
415 | simEmc->SetADCchannelW(0.0125); |
28b75255 |
416 | } |
417 | |
418 | |
419 | if (iPMD) |
420 | { |
421 | //=================== PMD parameters ============================ |
422 | |
423 | AliPMD *PMD = new AliPMDv1("PMD", "normal PMD"); |
424 | } |
425 | |
426 | if (iT0) |
427 | { |
428 | //=================== T0 parameters ============================ |
429 | AliT0 *T0 = new AliT0v1("T0", "T0 Detector"); |
430 | } |
431 | |
432 | if (iEMCAL) |
433 | { |
434 | //=================== EMCAL parameters ============================ |
435 | |
8224b11d |
436 | AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE"); |
28b75255 |
437 | } |
438 | |
439 | if (iACORDE) |
440 | { |
441 | //=================== ACORDE parameters ============================ |
442 | |
443 | AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE"); |
444 | } |
445 | |
446 | if (iVZERO) |
447 | { |
448 | //=================== ACORDE parameters ============================ |
449 | |
450 | AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO"); |
451 | } |
452 | } |
453 | // |
454 | // PYTHIA |
455 | // |
456 | |
457 | AliGenerator* MbPythia() |
458 | { |
459 | comment = comment.Append(" pp at 14 TeV: Pythia low-pt"); |
460 | // |
461 | // Pythia |
462 | AliGenPythia* pythia = new AliGenPythia(-1); |
463 | pythia->SetMomentumRange(0, 999999.); |
464 | pythia->SetThetaRange(0., 180.); |
465 | pythia->SetYRange(-12.,12.); |
466 | pythia->SetPtRange(0,1000.); |
467 | pythia->SetProcess(kPyMb); |
468 | pythia->SetEnergyCMS(energy); |
469 | |
470 | return pythia; |
471 | } |
472 | |
473 | AliGenerator* MbPhojet() |
474 | { |
475 | comment = comment.Append(" pp at 14 TeV: Pythia low-pt"); |
476 | // |
477 | // DPMJET |
478 | #if defined(__CINT__) |
479 | gSystem->Load("libdpmjet"); // Parton density functions |
480 | gSystem->Load("libTDPMjet"); // Parton density functions |
481 | #endif |
482 | AliGenDPMjet* dpmjet = new AliGenDPMjet(-1); |
483 | dpmjet->SetMomentumRange(0, 999999.); |
484 | dpmjet->SetThetaRange(0., 180.); |
485 | dpmjet->SetYRange(-12.,12.); |
486 | dpmjet->SetPtRange(0,1000.); |
487 | dpmjet->SetProcess(kDpmMb); |
488 | dpmjet->SetEnergyCMS(energy); |
489 | |
490 | return dpmjet; |
491 | } |
492 | |
493 | void ProcessEnvironmentVars() |
494 | { |
495 | // Run type |
496 | if (gSystem->Getenv("CONFIG_RUN_TYPE")) { |
497 | for (Int_t iRun = 0; iRun < kRunMax; iRun++) { |
498 | if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) { |
499 | proc = (PDC06Proc_t)iRun; |
500 | cout<<"Run type set to "<<pprRunName[iRun]<<endl; |
501 | } |
502 | } |
503 | } |
504 | |
505 | // Field |
506 | if (gSystem->Getenv("CONFIG_FIELD")) { |
507 | for (Int_t iField = 0; iField < kFieldMax; iField++) { |
508 | if (strcmp(gSystem->Getenv("CONFIG_FIELD"), pprField[iField])==0) { |
509 | mag = (Mag_t)iField; |
510 | cout<<"Field set to "<<pprField[iField]<<endl; |
511 | } |
512 | } |
513 | } |
514 | |
515 | // Energy |
516 | if (gSystem->Getenv("CONFIG_ENERGY")) { |
517 | energy = atoi(gSystem->Getenv("CONFIG_ENERGY")); |
aad1175b |
518 | cout<<"Energy set to "<<energy<<" GeV"<<endl; |
28b75255 |
519 | } |
520 | |
521 | // Random Number seed |
522 | if (gSystem->Getenv("CONFIG_SEED")) { |
523 | seed = atoi(gSystem->Getenv("CONFIG_SEED")); |
524 | } |
525 | } |