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