]> git.uio.no Git - u/mrichter/AliRoot.git/blame - macros/ConfigPPR.C
Adding HERWIG
[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;
36d4cb82 136 Int_t iDIPO = 1;
137 Int_t iFMD = 1;
138 Int_t iFRAME = 1;
139 Int_t iHALL = 1;
140 Int_t iITS = 1;
141 Int_t iMAG = 1;
142 Int_t iMUON = 1;
143 Int_t iPHOS = 1;
144 Int_t iPIPE = 1;
145 Int_t iPMD = 1;
146 Int_t iRICH = 1;
147 Int_t iSHIL = 1;
148 Int_t iSTART = 1;
149 Int_t iTOF = 1;
150 Int_t iTPC = 1;
151 Int_t iTRD = 1;
152 Int_t iZDC = 1;
153 Int_t iEMCAL = 0;
154
155 //=================== Alice BODY parameters =============================
156 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
157
158
159 if (iMAG)
160 {
161 //=================== MAG parameters ============================
162 // --- Start with Magnet since detector layouts may be depending ---
163 // --- on the selected Magnet dimensions ---
164 AliMAG *MAG = new AliMAG("MAG", "Magnet");
165 }
166
167
168 if (iABSO)
169 {
170 //=================== ABSO parameters ============================
171 AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
172 }
173
174 if (iDIPO)
175 {
176 //=================== DIPO parameters ============================
177
178 AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
179 }
180
181 if (iHALL)
182 {
183 //=================== HALL parameters ============================
184
185 AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
186 }
187
188
189 if (iFRAME)
190 {
191 //=================== FRAME parameters ============================
192
193 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
194 if (geo == kHoles) {
195 FRAME->SetHoles(1);
196 } else {
197 FRAME->SetHoles(0);
198 }
199 }
200
201 if (iSHIL)
202 {
203 //=================== SHIL parameters ============================
204
205 AliSHIL *SHIL = new AliSHILv0("SHIL", "Shielding");
206 }
207
208
209 if (iPIPE)
210 {
211 //=================== PIPE parameters ============================
212
213 AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
214 }
215
216 if(iITS) {
217
218//=================== ITS parameters ============================
219 //
220 // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
221 // almost all other detectors. This involves the fact that the ITS geometry
222 // still has several options to be followed in parallel in order to determine
223 // the best set-up which minimizes the induced background. All the geometries
224 // available to date are described in the following. Read carefully the comments
225 // and use the default version (the only one uncommented) unless you are making
226 // comparisons and you know what you are doing. In this case just uncomment the
227 // ITS geometry you want to use and run Aliroot.
228 //
229 // Detailed geometries:
230 //
231 //
232 //AliITS *ITS = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
233 //
234 //AliITS *ITS = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
235 //
236 AliITSvPPRasymm *ITS = new AliITSvPPRasymm("ITS","New ITS PPR detailed version with asymmetric services");
237 ITS->SetMinorVersion(2); // don't touch this parameter if you're not an ITS developer
238 ITS->SetReadDet(kFALSE); // don't touch this parameter if you're not an ITS developer
239 // ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det"); // don't touch this parameter if you're not an ITS developer
240 ITS->SetThicknessDet1(200.); // detector thickness on layer 1 must be in the range [100,300]
241 ITS->SetThicknessDet2(200.); // detector thickness on layer 2 must be in the range [100,300]
242 ITS->SetThicknessChip1(200.); // chip thickness on layer 1 must be in the range [150,300]
243 ITS->SetThicknessChip2(200.); // chip thickness on layer 2 must be in the range [150,300]
244 ITS->SetRails(1); // 1 --> rails in ; 0 --> rails out
245 ITS->SetCoolingFluid(1); // 1 --> water ; 0 --> freon
246 //
247 //AliITSvPPRsymm *ITS = new AliITSvPPRsymm("ITS","New ITS PPR detailed version with symmetric services");
248 //ITS->SetMinorVersion(2); // don't touch this parameter if you're not an ITS developer
249 //ITS->SetReadDet(kFALSE); // don't touch this parameter if you're not an ITS developer
250 //ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRsymm2.det"); // don't touch this parameter if you're not an ITS developer
251 //ITS->SetThicknessDet1(200.); // detector thickness on layer 1 must be in the range [100,300]
252 //ITS->SetThicknessDet2(200.); // detector thickness on layer 2 must be in the range [100,300]
253 //ITS->SetThicknessChip1(200.); // chip thickness on layer 1 must be in the range [150,300]
254 //ITS->SetThicknessChip2(200.); // chip thickness on layer 2 must be in the range [150,300]
255 //ITS->SetRails(1); // 1 --> rails in ; 0 --> rails out
256 //ITS->SetCoolingFluid(1); // 1 --> water ; 0 --> freon
257 //
258 //
259 // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful
260 // for reconstruction !):
261 //
262 //
263 //AliITSvPPRcoarseasymm *ITS = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
264 //ITS->SetRails(1); // 1 --> rails in ; 0 --> rails out
265 //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
266 //
267 //AliITS *ITS = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
268 //ITS->SetRails(1); // 1 --> rails in ; 0 --> rails out
269 //ITS->SetSupportMaterial(0); // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
270 //
271 //
272 //
273 // Geant3 <-> EUCLID conversion
274 // ============================
275 //
276 // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
277 // media to two ASCII files (called by default ITSgeometry.euc and
278 // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
279 // The default (=0) means that you dont want to use this facility.
280 //
281 ITS->SetEUCLID(0);
282 }
283
284
285 if (iTPC)
286 {
287 //============================ TPC parameters ================================
288 // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
289 // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
290 // --- sectors are specified, any value other than that requires at least one
291 // --- sector (lower or upper)to be specified!
292 // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
293 // --- sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
294 // --- SecLows - number of lower sectors specified (up to 6)
295 // --- SecUps - number of upper sectors specified (up to 12)
296 // --- Sens - sensitive strips for the Slow Simulator !!!
297 // --- This does NOT work if all S or L-sectors are specified, i.e.
298 // --- if SecAL or SecAU < 0
299 //
300 //
301 //-----------------------------------------------------------------------------
302
303 // gROOT->LoadMacro("SetTPCParam.C");
304 // AliTPCParam *param = SetTPCParam();
305 AliTPC *TPC = new AliTPCv2("TPC", "Default");
306
307 // All sectors included
308 TPC->SetSecAL(-1);
309 TPC->SetSecAU(-1);
310
311 }
312
313 if (iTOF) {
314 if (geo == kHoles) {
315 //=================== TOF parameters ============================
316 AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
317 } else {
318 AliTOF *TOF = new AliTOFv4("TOF", "normal TOF");
319 }
320 }
321
322 if (iRICH)
323 {
324 //=================== RICH parameters ===========================
325 AliRICH *RICH = new AliRICHv3("RICH", "normal RICH");
326
327 }
328
329
330 if (iZDC)
331 {
332 //=================== ZDC parameters ============================
333
334 AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
335 }
336
36d4cb82 337 if (iTRD)
338 {
339 //=================== TRD parameters ============================
340
341 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
342
343 // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
344 TRD->SetGasMix(1);
345 if (geo == kHoles) {
346 // With hole in front of PHOS
347 TRD->SetPHOShole();
348 // With hole in front of RICH
349 TRD->SetRICHhole();
350 }
351 // Switch on TR
352 AliTRDsim *TRDsim = TRD->CreateTR();
353 }
354
355 if (iFMD)
356 {
357 //=================== FMD parameters ============================
358
359 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
360 FMD->SetRingsSi1(256);
361 FMD->SetRingsSi2(64);
362 FMD->SetSectorsSi1(20);
363 FMD->SetSectorsSi2(24);
364 }
365
366 if (iMUON)
367 {
368 //=================== MUON parameters ===========================
369
370 AliMUON *MUON = new AliMUONv1("MUON", "default");
371 }
372 //=================== PHOS parameters ===========================
373
374 if (iPHOS)
375 {
376 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
377 }
378
379
380 if (iPMD)
381 {
382 //=================== PMD parameters ============================
383 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
384 PMD->SetPAR(1., 1., 0.8, 0.02);
385 PMD->SetIN(6., 18., -580., 27., 27.);
386 PMD->SetGEO(0.0, 0.2, 4.);
387 PMD->SetPadSize(0.8, 1.0, 1.0, 1.5);
388 }
389
390 if (iSTART)
391 {
392 //=================== START parameters ============================
393 AliSTART *START = new AliSTARTv1("START", "START Detector");
394 }
395
396 if (iEMCAL && !iRICH)
397 {
398 //=================== EMCAL parameters ============================
399 AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCALArch1a");
400 }
401}
402
403Float_t EtaToTheta(Float_t arg){
404 return (180./TMath::Pi())*2.*atan(exp(-arg));
405}
406
407
408
409AliGenerator* GeneratorFactory(PprRun_t run) {
410 Int_t isw = 3;
411 if (rad == kNoGluonRadiation) isw = 0;
412
413
414 switch (run) {
415 case test50:
416 comment = comment.Append(":HIJINGparam test 50 particles");
417 AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
418 gener->SetMomentumRange(0, 999999.);
419 gener->SetPhiRange(-180., 180.);
420 // Set pseudorapidity range from -8 to 8.
421 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
422 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
423 gener->SetThetaRange(thmin,thmax);
424 break;
425 case kParam_8000:
426 comment = comment.Append(":HIJINGparam N=8000");
427 AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
428 gener->SetMomentumRange(0, 999999.);
429 gener->SetPhiRange(-180., 180.);
430 // Set pseudorapidity range from -8 to 8.
431 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
432 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
433 gener->SetThetaRange(thmin,thmax);
434 break;
435 case kParam_4000:
436 comment = comment.Append("HIJINGparam N=4000");
437 AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
438 gener->SetMomentumRange(0, 999999.);
439 gener->SetPhiRange(-180., 180.);
440 // Set pseudorapidity range from -8 to 8.
441 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
442 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
443 gener->SetThetaRange(thmin,thmax);
444 break;
445 case kParam_2000:
446 comment = comment.Append("HIJINGparam N=2000");
447 AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
448 gener->SetMomentumRange(0, 999999.);
449 gener->SetPhiRange(-180., 180.);
450 // Set pseudorapidity range from -8 to 8.
451 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
452 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
453 gener->SetThetaRange(thmin,thmax);
454 break;
455//
456// Hijing Central
457//
458 case kHijing_cent1:
459 comment = comment.Append("HIJING cent1");
460 AliGenHijing *gener = HijingStandard();
461// impact parameter range
462 gener->SetImpactParameterRange(0., 5.);
463 break;
464 case kHijing_cent2:
465 comment = comment.Append("HIJING cent2");
466 AliGenHijing *gener = HijingStandard();
467// impact parameter range
468 gener->SetImpactParameterRange(0., 2.);
469 break;
470//
471// Hijing Peripheral
472//
473 case kHijing_per1:
474 comment = comment.Append("HIJING per1");
475 AliGenHijing *gener = HijingStandard();
476// impact parameter range
477 gener->SetImpactParameterRange(5., 8.6);
478 break;
479 case kHijing_per2:
480 comment = comment.Append("HIJING per2");
481 AliGenHijing *gener = HijingStandard();
482// impact parameter range
483 gener->SetImpactParameterRange(8.6, 11.2);
484 break;
485 case kHijing_per3:
486 comment = comment.Append("HIJING per3");
487 AliGenHijing *gener = HijingStandard();
488// impact parameter range
489 gener->SetImpactParameterRange(11.2, 13.2);
490 break;
491 case kHijing_per4:
492 comment = comment.Append("HIJING per4");
493 AliGenHijing *gener = HijingStandard();
494// impact parameter range
495 gener->SetImpactParameterRange(13.2, 15.);
496 break;
497 case kHijing_per5:
498 comment = comment.Append("HIJING per5");
499 AliGenHijing *gener = HijingStandard();
500// impact parameter range
501 gener->SetImpactParameterRange(15., 100.);
502 break;
503//
504// Jet-Jet
505//
506 case kHijing_jj25:
507 comment = comment.Append("HIJING Jet 25 GeV");
508 AliGenHijing *gener = HijingStandard();
509// impact parameter range
510 gener->SetImpactParameterRange(0., 5.);
511 // trigger
512 gener->SetTrigger(1);
513 gener->SetPtMinJet(25.);
514 gener->SetRadiation(isw);
515 gener->JetEtaRange(-0.3,0.3);
516 gener->JetPhiRange(15.,105.);
517 break;
518
519 case kHijing_jj50:
520 comment = comment.Append("HIJING Jet 50 GeV");
521 AliGenHijing *gener = HijingStandard();
522// impact parameter range
523 gener->SetImpactParameterRange(0., 5.);
524 // trigger
525 gener->SetTrigger(1);
526 gener->SetPtMinJet(50.);
527 gener->SetRadiation(isw);
528 gener->JetEtaRange(-0.3,0.3);
529 gener->JetPhiRange(15.,105.);
530 break;
531
532 case kHijing_jj75:
533 comment = comment.Append("HIJING Jet 75 GeV");
534 AliGenHijing *gener = HijingStandard();
535// impact parameter range
536 gener->SetImpactParameterRange(0., 5.);
537 // trigger
538 gener->SetTrigger(1);
539 gener->SetPtMinJet(75.);
540 gener->SetRadiation(isw);
541 gener->JetEtaRange(-0.3,0.3);
542 gener->JetPhiRange(15.,105.);
543 break;
544
545 case kHijing_jj100:
546 comment = comment.Append("HIJING Jet 100 GeV");
547 AliGenHijing *gener = HijingStandard();
548// impact parameter range
549 gener->SetImpactParameterRange(0., 5.);
550 // trigger
551 gener->SetTrigger(1);
552 gener->SetPtMinJet(100.);
553 gener->SetRadiation(isw);
554 gener->JetEtaRange(-0.3,0.3);
555 gener->JetPhiRange(15.,105.);
556 break;
557
558 case kHijing_jj125:
559 comment = comment.Append("HIJING Jet 125 GeV");
560 AliGenHijing *gener = HijingStandard();
561// impact parameter range
562 gener->SetImpactParameterRange(0., 5.);
563 // trigger
564 gener->SetTrigger(1);
565 gener->SetPtMinJet(125.);
566 gener->SetRadiation(isw);
567 gener->JetEtaRange(-0.3,0.3);
568 gener->JetPhiRange(15.,105.);
569 break;
570//
571// Gamma-Jet
572//
573 case kHijing_gj25:
574 comment = comment.Append("HIJING Gamma 25 GeV");
575 AliGenHijing *gener = HijingStandard();
576// impact parameter range
577 gener->SetImpactParameterRange(0., 5.);
578 // trigger
579 gener->SetTrigger(2);
580 gener->SetPtMinJet(25.);
581 gener->SetRadiation(isw);
582 gener->JetEtaRange(-0.3,0.3);
583 gener->JetPhiRange(15.,105.);
584 break;
585
586 case kHijing_gj50:
587 comment = comment.Append("HIJING Gamma 50 GeV");
588 AliGenHijing *gener = HijingStandard();
589// impact parameter range
590 gener->SetImpactParameterRange(0., 5.);
591 // trigger
592 gener->SetTrigger(2);
593 gener->SetPtMinJet(50.);
594 gener->SetRadiation(isw);
595 gener->JetEtaRange(-0.3,0.3);
596 gener->JetPhiRange(15.,105.);
597 break;
598
599 case kHijing_gj75:
600 comment = comment.Append("HIJING Gamma 75 GeV");
601 AliGenHijing *gener = HijingStandard();
602// impact parameter range
603 gener->SetImpactParameterRange(0., 5.);
604 // trigger
605 gener->SetTrigger(2);
606 gener->SetPtMinJet(75.);
607 gener->SetRadiation(isw);
608 gener->JetEtaRange(-0.3,0.3);
609 gener->JetPhiRange(15.,105.);
610 break;
611
612 case kHijing_gj100:
613 comment = comment.Append("HIJING Gamma 100 GeV");
614 AliGenHijing *gener = HijingStandard();
615// impact parameter range
616 gener->SetImpactParameterRange(0., 5.);
617 // trigger
618 gener->SetTrigger(2);
619 gener->SetPtMinJet(100.);
620 gener->SetRadiation(isw);
621 gener->JetEtaRange(-0.3,0.3);
622 gener->JetPhiRange(15.,105.);
623 break;
624
625 case kHijing_gj125:
626 comment = comment.Append("HIJING Gamma 125 GeV");
627 AliGenHijing *gener = HijingStandard();
628// impact parameter range
629 gener->SetImpactParameterRange(0., 5.);
630 // trigger
631 gener->SetTrigger(2);
632 gener->SetPtMinJet(125.);
633 gener->SetRadiation(isw);
634 gener->JetEtaRange(-0.3,0.3);
635 gener->JetPhiRange(15.,105.);
636 break;
637 }
638 return gener;
639}
640
641AliGenHijing* HijingStandard()
642{
643 AliGenHijing *gener = new AliGenHijing(-1);
644// centre of mass energy
645 gener->SetEnergyCMS(5500.);
646// reference frame
647 gener->SetReferenceFrame("CMS");
648// projectile
649 gener->SetProjectile("A", 208, 82);
650 gener->SetTarget ("A", 208, 82);
651// tell hijing to keep the full parent child chain
652 gener->KeepFullEvent();
653// enable jet quenching
654 gener->SetJetQuenching(1);
655// enable shadowing
656 gener->SetShadowing(1);
657// neutral pion and heavy particle decays switched off
658 gener->SetDecaysOff(1);
659// Don't track spectators
660 gener->SetSpectators(0);
661// kinematic selection
662 gener->SetSelectAll(0);
663 return gener;
664}
665
666
667