1d6c81d7889489ac9cff02845c6fe754f00343c1
[u/mrichter/AliRoot.git] / macros / ConfigPPR.C
1 enum 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
11 enum PprGeo_t 
12 {
13     kHoles, kNoHoles
14 };
15
16 enum PprRad_t
17 {
18     kGluonRadiation, kNoGluonRadiation
19 };
20
21 enum PprMag_t
22 {
23     k2kG, k4kG, k5kG
24 };
25
26
27 // This part for configuration    
28 static PprRun_t run = test50;
29 static PprGeo_t geo = kHoles;
30 static PprRad_t rad = kGluonRadiation;
31 static PprMag_t mag = k4kG;
32
33 // Comment line 
34 static TString  comment;
35
36
37
38 void Config()
39 {
40
41     // 7-DEC-2000 09:00
42     // Switch on Transition adiation simulation. 6/12/00 18:00
43     // iZDC=1  7/12/00 09:00
44     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
45     // Theta range given through pseudorapidity limits 22/6/2001
46
47     // Set Random Number seed
48     // gRandom->SetSeed(12345);
49
50     new AliGeant3("C++ Interface to Geant3");
51
52     if (!gSystem->Getenv("CONFIG_FILE"))
53     {
54         TFile  *rootfile = new TFile("galice.root", "recreate");
55
56         rootfile->SetCompressionLevel(2);
57     }
58
59     TGeant3 *geant3 = (TGeant3 *) gMC;
60
61     //
62     // Set External decayer
63     AliDecayer *decayer = new AliDecayerPythia();
64
65     decayer->SetForceDecay(kAll);
66     decayer->Init();
67     gMC->SetExternalDecayer(decayer);
68     //
69     //
70     //=======================================================================
71     // ******* GEANT STEERING parameters FOR ALICE SIMULATION *******
72     geant3->SetTRIG(1);         //Number of events to be processed 
73     geant3->SetSWIT(4, 10);
74     geant3->SetDEBU(0, 0, 1);
75     //geant3->SetSWIT(2,2);
76     geant3->SetDCAY(1);
77     geant3->SetPAIR(1);
78     geant3->SetCOMP(1);
79     geant3->SetPHOT(1);
80     geant3->SetPFIS(0);
81     geant3->SetDRAY(0);
82     geant3->SetANNI(1);
83     geant3->SetBREM(1);
84     geant3->SetMUNU(1);
85     geant3->SetCKOV(1);
86     geant3->SetHADR(1);         //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3)
87     geant3->SetLOSS(2);
88     geant3->SetMULS(1);
89     geant3->SetRAYL(1);
90     geant3->SetAUTO(1);         //Select automatic STMIN etc... calc. (AUTO 1) or manual (AUTO 0)
91     geant3->SetABAN(0);         //Restore 3.16 behaviour for abandoned tracks
92     geant3->SetOPTI(2);         //Select optimisation level for GEANT geometry searches (0,1,2)
93     geant3->SetERAN(5.e-7);
94
95     Float_t cut = 1.e-3;        // 1MeV cut by default
96     Float_t tofmax = 1.e10;
97
98     //             GAM ELEC NHAD CHAD MUON EBREM MUHAB EDEL MUDEL MUPA TOFMAX
99     geant3->SetCUTS(cut, cut, cut, cut, cut, cut, cut, cut, cut, cut,
100                     tofmax);
101     //
102     //=======================================================================
103     // ************* STEERING parameters FOR ALICE SIMULATION **************
104     // --- Specify event type to be tracked through the ALICE setup
105     // --- All positions are in cm, angles in degrees, and P and E in GeV
106
107
108 // Generator Configuration
109     gAlice->SetDebug(1);
110     AliGenerator* gener = GeneratorFactory(run);
111     gener->SetOrigin(0, 0, 0);    // vertex position
112     gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
113     gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
114     gener->SetVertexSmear(kPerEvent); 
115     gener->SetTrackingFlag(1);
116     gener->Init();
117     
118     if (mag == k2kG) {
119         comment = comment.Append(" | L3 field 0.2 T");
120     } else if (mag == k4kG) {
121         comment = comment.Append(" | L3 field 0.4 T");
122     } else if (mag == k5kG) {
123         comment = comment.Append(" | L3 field 0.5 T");
124     }
125     
126     
127     if (rad == kGluonRadiation)
128     {
129         comment = comment.Append(" | Gluon Radiation On");
130         
131     } else {
132         comment = comment.Append(" | Gluon Radiation Off");
133     }
134
135     if (geo == kHoles)
136     {
137         comment = comment.Append(" | Holes for PHOS/RICH");
138         
139     } else {
140         comment = comment.Append(" | No holes for PHOS/RICH");
141     }
142
143     printf("\n \n Comment: %s \n \n", (char*) comment);
144     
145     
146 // Field (L3 0.4 T)
147     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag);
148     rootfile->cd();
149     gAlice->SetField(field);    
150     
151 //
152     Int_t   iABSO   = 1;
153     Int_t   iDIPO   = 1;
154     Int_t   iFMD    = 1;
155     Int_t   iFRAME  = 1;
156     Int_t   iHALL   = 1;
157     Int_t   iITS    = 1;
158     Int_t   iMAG    = 1;
159     Int_t   iMUON   = 1;
160     Int_t   iPHOS   = 1;
161     Int_t   iPIPE   = 1;
162     Int_t   iPMD    = 1;
163     Int_t   iRICH   = 1;
164     Int_t   iSHIL   = 1;
165     Int_t   iSTART  = 1;
166     Int_t   iTOF    = 1;
167     Int_t   iTPC    = 1;
168     Int_t   iTRD    = 1;
169     Int_t   iZDC    = 1;
170     Int_t   iEMCAL  = 1;
171     Int_t   iVZERO  = 1;
172     Int_t   iCRT    = 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 Version 2");
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 AliRICHv3("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 (iTRD)
357     {
358         //=================== TRD parameters ============================
359
360         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
361
362         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
363         TRD->SetGasMix(1);
364         if (geo == kHoles) {
365             // With hole in front of PHOS
366             TRD->SetPHOShole();
367             // With hole in front of RICH
368             TRD->SetRICHhole();
369         }
370             // Switch on TR
371             AliTRDsim *TRDsim = TRD->CreateTR();
372     }
373
374     if (iFMD)
375     {
376         //=================== FMD parameters ============================
377         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
378         FMD->SetRingsSi1(256);
379         FMD->SetRingsSi2(128);
380         FMD->SetSectorsSi1(20);
381         FMD->SetSectorsSi2(40);      
382    }
383
384     if (iMUON)
385     {
386         //=================== MUON parameters ===========================
387
388         AliMUON *MUON = new AliMUONv1("MUON", "default");
389     }
390     //=================== PHOS parameters ===========================
391
392     if (iPHOS)
393     {
394         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
395     }
396
397
398     if (iPMD)
399     {
400         //=================== PMD parameters ============================
401         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
402         PMD->SetPAR(1., 1., 0.8, 0.02);
403         PMD->SetIN(6., 18., -580., 27., 27.);
404         PMD->SetGEO(0.0, 0.2, 4.);
405         PMD->SetPadSize(0.8, 1.0, 1.0, 1.5);
406     }
407
408     if (iSTART)
409     {
410         //=================== START parameters ============================
411         AliSTART *START = new AliSTARTv1("START", "START Detector");
412     }
413
414     if (iEMCAL)
415     {
416         //=================== EMCAL parameters ============================
417         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCALArch1a");
418     }
419
420      if (iCRT)
421     {
422         //=================== CRT parameters ============================
423         AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
424     }
425
426      if (iVZERO)
427     {
428         //=================== CRT parameters ============================
429         AliVZERO *VZERO = new AliVZEROv2("VZERO", "normal VZERO");
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 kParam_8000:
458         comment = comment.Append(":HIJINGparam N=8000");
459         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
460         gener->SetMomentumRange(0, 999999.);
461         gener->SetPhiRange(-180., 180.);
462         // Set pseudorapidity range from -8 to 8.
463         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
464         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
465         gener->SetThetaRange(thmin,thmax);
466         break;
467     case kParam_4000:
468         comment = comment.Append("HIJINGparam N=4000");
469         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
470         gener->SetMomentumRange(0, 999999.);
471         gener->SetPhiRange(-180., 180.);
472         // Set pseudorapidity range from -8 to 8.
473         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
474         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
475         gener->SetThetaRange(thmin,thmax);
476         break;
477     case kParam_2000:
478         comment = comment.Append("HIJINGparam N=2000");
479         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
480         gener->SetMomentumRange(0, 999999.);
481         gener->SetPhiRange(-180., 180.);
482         // Set pseudorapidity range from -8 to 8.
483         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
484         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
485         gener->SetThetaRange(thmin,thmax);
486         break;
487 //
488 //  Hijing Central
489 //
490     case kHijing_cent1:
491         comment = comment.Append("HIJING cent1");
492         AliGenHijing *gener = HijingStandard();
493 // impact parameter range
494         gener->SetImpactParameterRange(0., 5.);
495         break;
496     case kHijing_cent2:
497         comment = comment.Append("HIJING cent2");
498         AliGenHijing *gener = HijingStandard();
499 // impact parameter range
500         gener->SetImpactParameterRange(0., 2.);
501         break;
502 //
503 // Hijing Peripheral 
504 //
505     case kHijing_per1:
506         comment = comment.Append("HIJING per1");
507         AliGenHijing *gener = HijingStandard();
508 // impact parameter range
509         gener->SetImpactParameterRange(5., 8.6);
510         break;
511     case kHijing_per2:
512         comment = comment.Append("HIJING per2");
513         AliGenHijing *gener = HijingStandard();
514 // impact parameter range
515         gener->SetImpactParameterRange(8.6, 11.2);
516         break;
517     case kHijing_per3:
518         comment = comment.Append("HIJING per3");
519         AliGenHijing *gener = HijingStandard();
520 // impact parameter range
521         gener->SetImpactParameterRange(11.2, 13.2);
522         break;
523     case kHijing_per4:
524         comment = comment.Append("HIJING per4");
525         AliGenHijing *gener = HijingStandard();
526 // impact parameter range
527         gener->SetImpactParameterRange(13.2, 15.);
528         break;
529     case kHijing_per5:
530         comment = comment.Append("HIJING per5");
531         AliGenHijing *gener = HijingStandard();
532 // impact parameter range
533         gener->SetImpactParameterRange(15., 100.);
534         break;
535 //
536 //  Jet-Jet
537 //
538     case kHijing_jj25:
539         comment = comment.Append("HIJING Jet 25 GeV");
540         AliGenHijing *gener = HijingStandard();
541 // impact parameter range
542         gener->SetImpactParameterRange(0., 5.);
543         // trigger
544         gener->SetTrigger(1);
545         gener->SetPtMinJet(25.);
546         gener->SetRadiation(isw);
547         gener->JetEtaRange(-0.3,0.3);
548         gener->JetPhiRange(15.,105.);   
549         break;
550
551     case kHijing_jj50:
552         comment = comment.Append("HIJING Jet 50 GeV");
553         AliGenHijing *gener = HijingStandard();
554 // impact parameter range
555         gener->SetImpactParameterRange(0., 5.);
556         // trigger
557         gener->SetTrigger(1);
558         gener->SetPtMinJet(50.);
559         gener->SetRadiation(isw);
560         gener->JetEtaRange(-0.3,0.3);
561         gener->JetPhiRange(15.,105.);   
562         break;
563
564     case kHijing_jj75:
565         comment = comment.Append("HIJING Jet 75 GeV");
566         AliGenHijing *gener = HijingStandard();
567 // impact parameter range
568         gener->SetImpactParameterRange(0., 5.);
569         // trigger
570         gener->SetTrigger(1);
571         gener->SetPtMinJet(75.);
572         gener->SetRadiation(isw);
573         gener->JetEtaRange(-0.3,0.3);
574         gener->JetPhiRange(15.,105.);   
575         break;
576
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         gener->JetEtaRange(-0.3,0.3);
587         gener->JetPhiRange(15.,105.);   
588         break;
589
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         gener->JetEtaRange(-0.3,0.3);
600         gener->JetPhiRange(15.,105.);   
601         break;
602 //
603 // Gamma-Jet
604 //
605     case kHijing_gj25:
606         comment = comment.Append("HIJING Gamma 25 GeV");
607         AliGenHijing *gener = HijingStandard();
608 // impact parameter range
609         gener->SetImpactParameterRange(0., 5.);
610         // trigger
611         gener->SetTrigger(2);
612         gener->SetPtMinJet(25.);
613         gener->SetRadiation(isw);
614         gener->JetEtaRange(-0.3,0.3);
615         gener->JetPhiRange(15.,105.);   
616         break;
617
618     case kHijing_gj50:
619         comment = comment.Append("HIJING Gamma 50 GeV");
620         AliGenHijing *gener = HijingStandard();
621 // impact parameter range
622         gener->SetImpactParameterRange(0., 5.);
623         // trigger
624         gener->SetTrigger(2);
625         gener->SetPtMinJet(50.);
626         gener->SetRadiation(isw);
627         gener->JetEtaRange(-0.3,0.3);
628         gener->JetPhiRange(15.,105.);   
629         break;
630
631     case kHijing_gj75:
632         comment = comment.Append("HIJING Gamma 75 GeV");
633         AliGenHijing *gener = HijingStandard();
634 // impact parameter range
635         gener->SetImpactParameterRange(0., 5.);
636         // trigger
637         gener->SetTrigger(2);
638         gener->SetPtMinJet(75.);
639         gener->SetRadiation(isw);
640         gener->JetEtaRange(-0.3,0.3);
641         gener->JetPhiRange(15.,105.);   
642         break;
643
644     case kHijing_gj100:
645         comment = comment.Append("HIJING Gamma 100 GeV");
646         AliGenHijing *gener = HijingStandard();
647 // impact parameter range
648         gener->SetImpactParameterRange(0., 5.);
649         // trigger
650         gener->SetTrigger(2);
651         gener->SetPtMinJet(100.);
652         gener->SetRadiation(isw);
653         gener->JetEtaRange(-0.3,0.3);
654         gener->JetPhiRange(15.,105.);   
655         break;
656
657     case kHijing_gj125:
658         comment = comment.Append("HIJING Gamma 125 GeV");
659         AliGenHijing *gener = HijingStandard();
660 // impact parameter range
661         gener->SetImpactParameterRange(0., 5.);
662         // trigger
663         gener->SetTrigger(2);
664         gener->SetPtMinJet(125.);
665         gener->SetRadiation(isw);
666         gener->JetEtaRange(-0.3,0.3);
667         gener->JetPhiRange(15.,105.);   
668         break;
669     }
670     return gener;
671 }
672
673 AliGenHijing* HijingStandard()
674 {
675     AliGenHijing *gener = new AliGenHijing(-1);
676 // centre of mass energy 
677     gener->SetEnergyCMS(5500.);
678 // reference frame
679     gener->SetReferenceFrame("CMS");
680 // projectile
681      gener->SetProjectile("A", 208, 82);
682      gener->SetTarget    ("A", 208, 82);
683 // tell hijing to keep the full parent child chain
684      gener->KeepFullEvent();
685 // enable jet quenching
686      gener->SetJetQuenching(1);
687 // enable shadowing
688      gener->SetShadowing(1);
689 // neutral pion and heavy particle decays switched off
690      gener->SetDecaysOff(1);
691 // Don't track spectators
692      gener->SetSpectators(0);
693 // kinematic selection
694      gener->SetSelectAll(0);
695      return gener;
696 }
697
698
699