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