First commit.
[u/mrichter/AliRoot.git] / macros / ConfigPPR.C
CommitLineData
36d4cb82 1enum PprRun_t
2{
3 test50,
4 kParam_8000, kParam_4000, kParam_2000,
5 kHijing_cent1, kHijing_cent2,
6 kHijing_per1, kHijing_per2, kHijing_per3, kHijing_per4, kHijing_per5,
7 kHijing_jj25, kHijing_jj50, kHijing_jj75, kHijing_jj100, kHijing_jj125,
8 kHijing_gj25, kHijing_gj50, kHijing_gj75, kHijing_gj100, kHijing_gj125
9};
10
11enum PprGeo_t
12{
13 kHoles, kNoHoles
14};
15
16enum PprRad_t
17{
18 kGluonRadiation, kNoGluonRadiation
19};
20
21// This part for configuration
22static PprRun_t run = test50;
23static PprGeo_t geo = kHoles;
24static PprRad_t rad = kGluonRadiation;
25// Comment line
26static TString comment;
27
28
29
30void Config()
31{
32
33 // 7-DEC-2000 09:00
34 // Switch on Transition adiation simulation. 6/12/00 18:00
35 // iZDC=1 7/12/00 09:00
36 // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
37 // Theta range given through pseudorapidity limits 22/6/2001
38
39 // Set Random Number seed
40 // gRandom->SetSeed(12345);
41
42 new AliGeant3("C++ Interface to Geant3");
43
44 if (!gSystem->Getenv("CONFIG_FILE"))
45 {
46 TFile *rootfile = new TFile("galice.root", "recreate");
47
48 rootfile->SetCompressionLevel(2);
49 }
50
51 TGeant3 *geant3 = (TGeant3 *) gMC;
52
53 //
54 // Set External decayer
55 AliDecayer *decayer = new AliDecayerPythia();
56
57 decayer->SetForceDecay(kAll);
58 decayer->Init();
59 gMC->SetExternalDecayer(decayer);
60 //
61 //
62 //=======================================================================
63 // ******* GEANT STEERING parameters FOR ALICE SIMULATION *******
64 geant3->SetTRIG(1); //Number of events to be processed
65 geant3->SetSWIT(4, 10);
66 geant3->SetDEBU(0, 0, 1);
67 //geant3->SetSWIT(2,2);
68 geant3->SetDCAY(1);
69 geant3->SetPAIR(1);
70 geant3->SetCOMP(1);
71 geant3->SetPHOT(1);
72 geant3->SetPFIS(0);
73 geant3->SetDRAY(0);
74 geant3->SetANNI(1);
75 geant3->SetBREM(1);
76 geant3->SetMUNU(1);
77 geant3->SetCKOV(1);
78 geant3->SetHADR(1); //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3)
79 geant3->SetLOSS(2);
80 geant3->SetMULS(1);
81 geant3->SetRAYL(1);
82 geant3->SetAUTO(1); //Select automatic STMIN etc... calc. (AUTO 1) or manual (AUTO 0)
83 geant3->SetABAN(0); //Restore 3.16 behaviour for abandoned tracks
84 geant3->SetOPTI(2); //Select optimisation level for GEANT geometry searches (0,1,2)
85 geant3->SetERAN(5.e-7);
86
87 Float_t cut = 1.e-3; // 1MeV cut by default
88 Float_t tofmax = 1.e10;
89
90 // GAM ELEC NHAD CHAD MUON EBREM MUHAB EDEL MUDEL MUPA TOFMAX
91 geant3->SetCUTS(cut, cut, cut, cut, cut, cut, cut, cut, cut, cut,
92 tofmax);
93 //
94 //=======================================================================
95 // ************* STEERING parameters FOR ALICE SIMULATION **************
96 // --- Specify event type to be tracked through the ALICE setup
97 // --- All positions are in cm, angles in degrees, and P and E in GeV
98
99
100// Generator Configuration
101 gAlice->SetDebug(1);
102 AliGenerator* gener = GeneratorFactory(run);
103 gener->SetOrigin(0, 0, 0); // vertex position
104 gener->SetSigma(0, 0, 5.3); // Sigma in (X,Y,Z) (cm) on IP position
105 gener->SetCutVertexZ(1.); // Truncate at 1 sigma
106 gener->SetVertexSmear(kPerEvent);
107 gener->SetTrackingFlag(1);
108 gener->Init();
109
110 if (rad == kGluonRadiation)
111 {
112 comment = comment.Append(" | Gluon Radiation On");
113
114 } else {
115 comment = comment.Append(" | Gluon Radiation Off");
116 }
117
118 if (geo == kHoles)
119 {
120 comment = comment.Append(" | Holes for PHOS/RICH");
121
122 } else {
123 comment = comment.Append(" | No holes for PHOS/RICH");
124 }
125
126 printf("\n \n Comment: %s \n \n", (char*) comment);
127
128
129// Field (L3 0.4 T)
130 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., 1);
131 rootfile->cd();
132 gAlice->SetField(field);
133
134//
135 Int_t iABSO = 1;
136 Int_t iCASTOR = 1;
137 Int_t iDIPO = 1;
138 Int_t iFMD = 1;
139 Int_t iFRAME = 1;
140 Int_t iHALL = 1;
141 Int_t iITS = 1;
142 Int_t iMAG = 1;
143 Int_t iMUON = 1;
144 Int_t iPHOS = 1;
145 Int_t iPIPE = 1;
146 Int_t iPMD = 1;
147 Int_t iRICH = 1;
148 Int_t iSHIL = 1;
149 Int_t iSTART = 1;
150 Int_t iTOF = 1;
151 Int_t iTPC = 1;
152 Int_t iTRD = 1;
153 Int_t iZDC = 1;
154 Int_t iEMCAL = 0;
155
156 //=================== Alice BODY parameters =============================
157 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
158
159
160 if (iMAG)
161 {
162 //=================== MAG parameters ============================
163 // --- Start with Magnet since detector layouts may be depending ---
164 // --- on the selected Magnet dimensions ---
165 AliMAG *MAG = new AliMAG("MAG", "Magnet");
166 }
167
168
169 if (iABSO)
170 {
171 //=================== ABSO parameters ============================
172 AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
173 }
174
175 if (iDIPO)
176 {
177 //=================== DIPO parameters ============================
178
179 AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
180 }
181
182 if (iHALL)
183 {
184 //=================== HALL parameters ============================
185
186 AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
187 }
188
189
190 if (iFRAME)
191 {
192 //=================== FRAME parameters ============================
193
194 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
195 if (geo == kHoles) {
196 FRAME->SetHoles(1);
197 } else {
198 FRAME->SetHoles(0);
199 }
200 }
201
202 if (iSHIL)
203 {
204 //=================== SHIL parameters ============================
205
206 AliSHIL *SHIL = new AliSHILv0("SHIL", "Shielding");
207 }
208
209
210 if (iPIPE)
211 {
212 //=================== PIPE parameters ============================
213
214 AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
215 }
216
217 if(iITS) {
218
219//=================== ITS parameters ============================
220 //
221 // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
222 // almost all other detectors. This involves the fact that the ITS geometry
223 // still has several options to be followed in parallel in order to determine
224 // the best set-up which minimizes the induced background. All the geometries
225 // available to date are described in the following. Read carefully the comments
226 // and use the default version (the only one uncommented) unless you are making
227 // comparisons and you know what you are doing. In this case just uncomment the
228 // ITS geometry you want to use and run Aliroot.
229 //
230 // Detailed geometries:
231 //
232 //
233 //AliITS *ITS = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
234 //
235 //AliITS *ITS = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
236 //
237 AliITSvPPRasymm *ITS = new AliITSvPPRasymm("ITS","New ITS PPR detailed version with asymmetric services");
238 ITS->SetMinorVersion(2); // don't touch this parameter if you're not an ITS developer
239 ITS->SetReadDet(kFALSE); // don't touch this parameter if you're not an ITS developer
240 // ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det"); // don't touch this parameter if you're not an ITS developer
241 ITS->SetThicknessDet1(200.); // detector thickness on layer 1 must be in the range [100,300]
242 ITS->SetThicknessDet2(200.); // detector thickness on layer 2 must be in the range [100,300]
243 ITS->SetThicknessChip1(200.); // chip thickness on layer 1 must be in the range [150,300]
244 ITS->SetThicknessChip2(200.); // chip thickness on layer 2 must be in the range [150,300]
245 ITS->SetRails(1); // 1 --> rails in ; 0 --> rails out
246 ITS->SetCoolingFluid(1); // 1 --> water ; 0 --> freon
247 //
248 //AliITSvPPRsymm *ITS = new AliITSvPPRsymm("ITS","New ITS PPR detailed version with symmetric services");
249 //ITS->SetMinorVersion(2); // don't touch this parameter if you're not an ITS developer
250 //ITS->SetReadDet(kFALSE); // don't touch this parameter if you're not an ITS developer
251 //ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRsymm2.det"); // don't touch this parameter if you're not an ITS developer
252 //ITS->SetThicknessDet1(200.); // detector thickness on layer 1 must be in the range [100,300]
253 //ITS->SetThicknessDet2(200.); // detector thickness on layer 2 must be in the range [100,300]
254 //ITS->SetThicknessChip1(200.); // chip thickness on layer 1 must be in the range [150,300]
255 //ITS->SetThicknessChip2(200.); // chip thickness on layer 2 must be in the range [150,300]
256 //ITS->SetRails(1); // 1 --> rails in ; 0 --> rails out
257 //ITS->SetCoolingFluid(1); // 1 --> water ; 0 --> freon
258 //
259 //
260 // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful
261 // for reconstruction !):
262 //
263 //
264 //AliITSvPPRcoarseasymm *ITS = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
265 //ITS->SetRails(1); // 1 --> rails in ; 0 --> rails out
266 //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
267 //
268 //AliITS *ITS = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
269 //ITS->SetRails(1); // 1 --> rails in ; 0 --> rails out
270 //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
271 //
272 //
273 //
274 // Geant3 <-> EUCLID conversion
275 // ============================
276 //
277 // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
278 // media to two ASCII files (called by default ITSgeometry.euc and
279 // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
280 // The default (=0) means that you dont want to use this facility.
281 //
282 ITS->SetEUCLID(0);
283 }
284
285
286 if (iTPC)
287 {
288 //============================ TPC parameters ================================
289 // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
290 // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
291 // --- sectors are specified, any value other than that requires at least one
292 // --- sector (lower or upper)to be specified!
293 // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
294 // --- sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
295 // --- SecLows - number of lower sectors specified (up to 6)
296 // --- SecUps - number of upper sectors specified (up to 12)
297 // --- Sens - sensitive strips for the Slow Simulator !!!
298 // --- This does NOT work if all S or L-sectors are specified, i.e.
299 // --- if SecAL or SecAU < 0
300 //
301 //
302 //-----------------------------------------------------------------------------
303
304 // gROOT->LoadMacro("SetTPCParam.C");
305 // AliTPCParam *param = SetTPCParam();
306 AliTPC *TPC = new AliTPCv2("TPC", "Default");
307
308 // All sectors included
309 TPC->SetSecAL(-1);
310 TPC->SetSecAU(-1);
311
312 }
313
314 if (iTOF) {
315 if (geo == kHoles) {
316 //=================== TOF parameters ============================
317 AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
318 } else {
319 AliTOF *TOF = new AliTOFv4("TOF", "normal TOF");
320 }
321 }
322
323 if (iRICH)
324 {
325 //=================== RICH parameters ===========================
326 AliRICH *RICH = new AliRICHv3("RICH", "normal RICH");
327
328 }
329
330
331 if (iZDC)
332 {
333 //=================== ZDC parameters ============================
334
335 AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
336 }
337
338 if (iCASTOR)
339 {
340 //=================== CASTOR parameters ============================
341
342 AliCASTOR *CASTOR = new AliCASTORv1("CASTOR", "normal CASTOR");
343 }
344
345 if (iTRD)
346 {
347 //=================== TRD parameters ============================
348
349 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
350
351 // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
352 TRD->SetGasMix(1);
353 if (geo == kHoles) {
354 // With hole in front of PHOS
355 TRD->SetPHOShole();
356 // With hole in front of RICH
357 TRD->SetRICHhole();
358 }
359 // Switch on TR
360 AliTRDsim *TRDsim = TRD->CreateTR();
361 }
362
363 if (iFMD)
364 {
365 //=================== FMD parameters ============================
366
367 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
368 FMD->SetRingsSi1(256);
369 FMD->SetRingsSi2(64);
370 FMD->SetSectorsSi1(20);
371 FMD->SetSectorsSi2(24);
372 }
373
374 if (iMUON)
375 {
376 //=================== MUON parameters ===========================
377
378 AliMUON *MUON = new AliMUONv1("MUON", "default");
379 }
380 //=================== PHOS parameters ===========================
381
382 if (iPHOS)
383 {
384 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
385 }
386
387
388 if (iPMD)
389 {
390 //=================== PMD parameters ============================
391 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
392 PMD->SetPAR(1., 1., 0.8, 0.02);
393 PMD->SetIN(6., 18., -580., 27., 27.);
394 PMD->SetGEO(0.0, 0.2, 4.);
395 PMD->SetPadSize(0.8, 1.0, 1.0, 1.5);
396 }
397
398 if (iSTART)
399 {
400 //=================== START parameters ============================
401 AliSTART *START = new AliSTARTv1("START", "START Detector");
402 }
403
404 if (iEMCAL && !iRICH)
405 {
406 //=================== EMCAL parameters ============================
407 AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCALArch1a");
408 }
409}
410
411Float_t EtaToTheta(Float_t arg){
412 return (180./TMath::Pi())*2.*atan(exp(-arg));
413}
414
415
416
417AliGenerator* GeneratorFactory(PprRun_t run) {
418 Int_t isw = 3;
419 if (rad == kNoGluonRadiation) isw = 0;
420
421
422 switch (run) {
423 case test50:
424 comment = comment.Append(":HIJINGparam test 50 particles");
425 AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
426 gener->SetMomentumRange(0, 999999.);
427 gener->SetPhiRange(-180., 180.);
428 // Set pseudorapidity range from -8 to 8.
429 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
430 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
431 gener->SetThetaRange(thmin,thmax);
432 break;
433 case kParam_8000:
434 comment = comment.Append(":HIJINGparam N=8000");
435 AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
436 gener->SetMomentumRange(0, 999999.);
437 gener->SetPhiRange(-180., 180.);
438 // Set pseudorapidity range from -8 to 8.
439 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
440 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
441 gener->SetThetaRange(thmin,thmax);
442 break;
443 case kParam_4000:
444 comment = comment.Append("HIJINGparam N=4000");
445 AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
446 gener->SetMomentumRange(0, 999999.);
447 gener->SetPhiRange(-180., 180.);
448 // Set pseudorapidity range from -8 to 8.
449 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
450 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
451 gener->SetThetaRange(thmin,thmax);
452 break;
453 case kParam_2000:
454 comment = comment.Append("HIJINGparam N=2000");
455 AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
456 gener->SetMomentumRange(0, 999999.);
457 gener->SetPhiRange(-180., 180.);
458 // Set pseudorapidity range from -8 to 8.
459 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
460 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
461 gener->SetThetaRange(thmin,thmax);
462 break;
463//
464// Hijing Central
465//
466 case kHijing_cent1:
467 comment = comment.Append("HIJING cent1");
468 AliGenHijing *gener = HijingStandard();
469// impact parameter range
470 gener->SetImpactParameterRange(0., 5.);
471 break;
472 case kHijing_cent2:
473 comment = comment.Append("HIJING cent2");
474 AliGenHijing *gener = HijingStandard();
475// impact parameter range
476 gener->SetImpactParameterRange(0., 2.);
477 break;
478//
479// Hijing Peripheral
480//
481 case kHijing_per1:
482 comment = comment.Append("HIJING per1");
483 AliGenHijing *gener = HijingStandard();
484// impact parameter range
485 gener->SetImpactParameterRange(5., 8.6);
486 break;
487 case kHijing_per2:
488 comment = comment.Append("HIJING per2");
489 AliGenHijing *gener = HijingStandard();
490// impact parameter range
491 gener->SetImpactParameterRange(8.6, 11.2);
492 break;
493 case kHijing_per3:
494 comment = comment.Append("HIJING per3");
495 AliGenHijing *gener = HijingStandard();
496// impact parameter range
497 gener->SetImpactParameterRange(11.2, 13.2);
498 break;
499 case kHijing_per4:
500 comment = comment.Append("HIJING per4");
501 AliGenHijing *gener = HijingStandard();
502// impact parameter range
503 gener->SetImpactParameterRange(13.2, 15.);
504 break;
505 case kHijing_per5:
506 comment = comment.Append("HIJING per5");
507 AliGenHijing *gener = HijingStandard();
508// impact parameter range
509 gener->SetImpactParameterRange(15., 100.);
510 break;
511//
512// Jet-Jet
513//
514 case kHijing_jj25:
515 comment = comment.Append("HIJING Jet 25 GeV");
516 AliGenHijing *gener = HijingStandard();
517// impact parameter range
518 gener->SetImpactParameterRange(0., 5.);
519 // trigger
520 gener->SetTrigger(1);
521 gener->SetPtMinJet(25.);
522 gener->SetRadiation(isw);
523 gener->JetEtaRange(-0.3,0.3);
524 gener->JetPhiRange(15.,105.);
525 break;
526
527 case kHijing_jj50:
528 comment = comment.Append("HIJING Jet 50 GeV");
529 AliGenHijing *gener = HijingStandard();
530// impact parameter range
531 gener->SetImpactParameterRange(0., 5.);
532 // trigger
533 gener->SetTrigger(1);
534 gener->SetPtMinJet(50.);
535 gener->SetRadiation(isw);
536 gener->JetEtaRange(-0.3,0.3);
537 gener->JetPhiRange(15.,105.);
538 break;
539
540 case kHijing_jj75:
541 comment = comment.Append("HIJING Jet 75 GeV");
542 AliGenHijing *gener = HijingStandard();
543// impact parameter range
544 gener->SetImpactParameterRange(0., 5.);
545 // trigger
546 gener->SetTrigger(1);
547 gener->SetPtMinJet(75.);
548 gener->SetRadiation(isw);
549 gener->JetEtaRange(-0.3,0.3);
550 gener->JetPhiRange(15.,105.);
551 break;
552
553 case kHijing_jj100:
554 comment = comment.Append("HIJING Jet 100 GeV");
555 AliGenHijing *gener = HijingStandard();
556// impact parameter range
557 gener->SetImpactParameterRange(0., 5.);
558 // trigger
559 gener->SetTrigger(1);
560 gener->SetPtMinJet(100.);
561 gener->SetRadiation(isw);
562 gener->JetEtaRange(-0.3,0.3);
563 gener->JetPhiRange(15.,105.);
564 break;
565
566 case kHijing_jj125:
567 comment = comment.Append("HIJING Jet 125 GeV");
568 AliGenHijing *gener = HijingStandard();
569// impact parameter range
570 gener->SetImpactParameterRange(0., 5.);
571 // trigger
572 gener->SetTrigger(1);
573 gener->SetPtMinJet(125.);
574 gener->SetRadiation(isw);
575 gener->JetEtaRange(-0.3,0.3);
576 gener->JetPhiRange(15.,105.);
577 break;
578//
579// Gamma-Jet
580//
581 case kHijing_gj25:
582 comment = comment.Append("HIJING Gamma 25 GeV");
583 AliGenHijing *gener = HijingStandard();
584// impact parameter range
585 gener->SetImpactParameterRange(0., 5.);
586 // trigger
587 gener->SetTrigger(2);
588 gener->SetPtMinJet(25.);
589 gener->SetRadiation(isw);
590 gener->JetEtaRange(-0.3,0.3);
591 gener->JetPhiRange(15.,105.);
592 break;
593
594 case kHijing_gj50:
595 comment = comment.Append("HIJING Gamma 50 GeV");
596 AliGenHijing *gener = HijingStandard();
597// impact parameter range
598 gener->SetImpactParameterRange(0., 5.);
599 // trigger
600 gener->SetTrigger(2);
601 gener->SetPtMinJet(50.);
602 gener->SetRadiation(isw);
603 gener->JetEtaRange(-0.3,0.3);
604 gener->JetPhiRange(15.,105.);
605 break;
606
607 case kHijing_gj75:
608 comment = comment.Append("HIJING Gamma 75 GeV");
609 AliGenHijing *gener = HijingStandard();
610// impact parameter range
611 gener->SetImpactParameterRange(0., 5.);
612 // trigger
613 gener->SetTrigger(2);
614 gener->SetPtMinJet(75.);
615 gener->SetRadiation(isw);
616 gener->JetEtaRange(-0.3,0.3);
617 gener->JetPhiRange(15.,105.);
618 break;
619
620 case kHijing_gj100:
621 comment = comment.Append("HIJING Gamma 100 GeV");
622 AliGenHijing *gener = HijingStandard();
623// impact parameter range
624 gener->SetImpactParameterRange(0., 5.);
625 // trigger
626 gener->SetTrigger(2);
627 gener->SetPtMinJet(100.);
628 gener->SetRadiation(isw);
629 gener->JetEtaRange(-0.3,0.3);
630 gener->JetPhiRange(15.,105.);
631 break;
632
633 case kHijing_gj125:
634 comment = comment.Append("HIJING Gamma 125 GeV");
635 AliGenHijing *gener = HijingStandard();
636// impact parameter range
637 gener->SetImpactParameterRange(0., 5.);
638 // trigger
639 gener->SetTrigger(2);
640 gener->SetPtMinJet(125.);
641 gener->SetRadiation(isw);
642 gener->JetEtaRange(-0.3,0.3);
643 gener->JetPhiRange(15.,105.);
644 break;
645 }
646 return gener;
647}
648
649AliGenHijing* HijingStandard()
650{
651 AliGenHijing *gener = new AliGenHijing(-1);
652// centre of mass energy
653 gener->SetEnergyCMS(5500.);
654// reference frame
655 gener->SetReferenceFrame("CMS");
656// projectile
657 gener->SetProjectile("A", 208, 82);
658 gener->SetTarget ("A", 208, 82);
659// tell hijing to keep the full parent child chain
660 gener->KeepFullEvent();
661// enable jet quenching
662 gener->SetJetQuenching(1);
663// enable shadowing
664 gener->SetShadowing(1);
665// neutral pion and heavy particle decays switched off
666 gener->SetDecaysOff(1);
667// Don't track spectators
668 gener->SetSpectators(0);
669// kinematic selection
670 gener->SetSelectAll(0);
671 return gener;
672}
673
674
675