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