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