]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/ConfigPPR.C
ntracks declared Int_t
[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_jj200, 
8     kHijing_gj25,  kHijing_gj50, kHijing_gj75, kHijing_gj100, kHijing_gj200
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 TGeant3("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(0);            // 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(0);              // 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(0);                // 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(0);                // 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     if (iTPC)
304     {
305         //============================ TPC parameters ================================
306         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
307         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
308         // --- sectors are specified, any value other than that requires at least one 
309         // --- sector (lower or upper)to be specified!
310         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
311         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
312         // --- SecLows - number of lower sectors specified (up to 6)
313         // --- SecUps - number of upper sectors specified (up to 12)
314         // --- Sens - sensitive strips for the Slow Simulator !!!
315         // --- This does NOT work if all S or L-sectors are specified, i.e.
316         // --- if SecAL or SecAU < 0
317         //
318         //
319         //-----------------------------------------------------------------------------
320
321         //  gROOT->LoadMacro("SetTPCParam.C");
322         //  AliTPCParam *param = SetTPCParam();
323         AliTPC *TPC = new AliTPCv2("TPC", "Default");
324
325         // All sectors included 
326         TPC->SetSecAL(-1);
327         TPC->SetSecAU(-1);
328
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 AliTOFv4T0("TOF", "normal TOF");
338         }
339     }
340
341
342     if (iRICH)
343     {
344         //=================== RICH parameters ===========================
345         AliRICH *RICH = new AliRICHv3("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 (iTRD)
358     {
359         //=================== TRD parameters ============================
360
361         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
362
363         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
364         TRD->SetGasMix(1);
365         if (geo == kHoles) {
366             // With hole in front of PHOS
367             TRD->SetPHOShole();
368             // With hole in front of RICH
369             TRD->SetRICHhole();
370         }
371             // Switch on TR
372             AliTRDsim *TRDsim = TRD->CreateTR();
373     }
374
375     if (iFMD)
376     {
377         //=================== FMD parameters ============================
378         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
379         FMD->SetRingsSi1(256);
380         FMD->SetRingsSi2(128);
381         FMD->SetSectorsSi1(20);
382         FMD->SetSectorsSi2(40);      
383    }
384
385     if (iMUON)
386     {
387         //=================== MUON parameters ===========================
388
389         AliMUON *MUON = new AliMUONv1("MUON", "default");
390     }
391     //=================== PHOS parameters ===========================
392
393     if (iPHOS)
394     {
395         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
396     }
397
398
399     if (iPMD)
400     {
401         //=================== PMD parameters ============================
402         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
403     }
404
405     if (iSTART)
406     {
407         //=================== START parameters ============================
408         AliSTART *START = new AliSTARTv1("START", "START Detector");
409     }
410
411     if (iEMCAL)
412     {
413         //=================== EMCAL parameters ============================
414         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCALArch1a");
415     }
416
417      if (iCRT)
418     {
419         //=================== CRT parameters ============================
420         AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
421     }
422
423      if (iVZERO)
424     {
425         //=================== CRT parameters ============================
426         AliVZERO *VZERO = new AliVZEROv2("VZERO", "normal VZERO");
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 kParam_8000:
455         comment = comment.Append(":HIJINGparam N=8000");
456         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
457         gener->SetMomentumRange(0, 999999.);
458         gener->SetPhiRange(-180., 180.);
459         // Set pseudorapidity range from -8 to 8.
460         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
461         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
462         gener->SetThetaRange(thmin,thmax);
463         break;
464     case kParam_4000:
465         comment = comment.Append("HIJINGparam N=4000");
466         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
467         gener->SetMomentumRange(0, 999999.);
468         gener->SetPhiRange(-180., 180.);
469         // Set pseudorapidity range from -8 to 8.
470         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
471         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
472         gener->SetThetaRange(thmin,thmax);
473         break;
474     case kParam_2000:
475         comment = comment.Append("HIJINGparam N=2000");
476         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
477         gener->SetMomentumRange(0, 999999.);
478         gener->SetPhiRange(-180., 180.);
479         // Set pseudorapidity range from -8 to 8.
480         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
481         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
482         gener->SetThetaRange(thmin,thmax);
483         break;
484 //
485 //  Hijing Central
486 //
487     case kHijing_cent1:
488         comment = comment.Append("HIJING cent1");
489         AliGenHijing *gener = HijingStandard();
490 // impact parameter range
491         gener->SetImpactParameterRange(0., 5.);
492         break;
493     case kHijing_cent2:
494         comment = comment.Append("HIJING cent2");
495         AliGenHijing *gener = HijingStandard();
496 // impact parameter range
497         gener->SetImpactParameterRange(0., 2.);
498         break;
499 //
500 // Hijing Peripheral 
501 //
502     case kHijing_per1:
503         comment = comment.Append("HIJING per1");
504         AliGenHijing *gener = HijingStandard();
505 // impact parameter range
506         gener->SetImpactParameterRange(5., 8.6);
507         break;
508     case kHijing_per2:
509         comment = comment.Append("HIJING per2");
510         AliGenHijing *gener = HijingStandard();
511 // impact parameter range
512         gener->SetImpactParameterRange(8.6, 11.2);
513         break;
514     case kHijing_per3:
515         comment = comment.Append("HIJING per3");
516         AliGenHijing *gener = HijingStandard();
517 // impact parameter range
518         gener->SetImpactParameterRange(11.2, 13.2);
519         break;
520     case kHijing_per4:
521         comment = comment.Append("HIJING per4");
522         AliGenHijing *gener = HijingStandard();
523 // impact parameter range
524         gener->SetImpactParameterRange(13.2, 15.);
525         break;
526     case kHijing_per5:
527         comment = comment.Append("HIJING per5");
528         AliGenHijing *gener = HijingStandard();
529 // impact parameter range
530         gener->SetImpactParameterRange(15., 100.);
531         break;
532 //
533 //  Jet-Jet
534 //
535     case kHijing_jj25:
536         comment = comment.Append("HIJING Jet 25 GeV");
537         AliGenHijing *gener = HijingStandard();
538 // impact parameter range
539         gener->SetImpactParameterRange(0., 5.);
540         // trigger
541         gener->SetTrigger(1);
542         gener->SetPtJet(25.);
543         gener->SetRadiation(isw);
544         gener->SetSimpleJets(!isw);
545         gener->SetJetEtaRange(-0.3,0.3);
546         gener->SetJetPhiRange(75., 165.);   
547         break;
548
549     case kHijing_jj50:
550         comment = comment.Append("HIJING Jet 50 GeV");
551         AliGenHijing *gener = HijingStandard();
552 // impact parameter range
553         gener->SetImpactParameterRange(0., 5.);
554         // trigger
555         gener->SetTrigger(1);
556         gener->SetPtJet(50.);
557         gener->SetRadiation(isw);
558         gener->SetSimpleJets(!isw);
559         gener->SetJetEtaRange(-0.3,0.3);
560         gener->SetJetPhiRange(75., 165.);   
561         break;
562
563     case kHijing_jj75:
564         comment = comment.Append("HIJING Jet 75 GeV");
565         AliGenHijing *gener = HijingStandard();
566 // impact parameter range
567         gener->SetImpactParameterRange(0., 5.);
568         // trigger
569         gener->SetTrigger(1);
570         gener->SetPtJet(75.);
571         gener->SetRadiation(isw);
572         gener->SetSimpleJets(!isw);
573         gener->SetJetEtaRange(-0.3,0.3);
574         gener->SetJetPhiRange(75., 165.);   
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->SetPtJet(100.);
585         gener->SetRadiation(isw);
586         gener->SetSimpleJets(!isw);
587         gener->SetJetEtaRange(-0.3,0.3);
588         gener->SetJetPhiRange(75., 165.);   
589         break;
590
591     case kHijing_jj200:
592         comment = comment.Append("HIJING Jet 200 GeV");
593         AliGenHijing *gener = HijingStandard();
594 // impact parameter range
595         gener->SetImpactParameterRange(0., 5.);
596         // trigger
597         gener->SetTrigger(1);
598         gener->SetPtJet(200.);
599         gener->SetRadiation(isw);
600         gener->SetSimpleJets(!isw);
601         gener->SetJetEtaRange(-0.3,0.3);
602         gener->SetJetPhiRange(75., 165.);   
603         break;
604 //
605 // Gamma-Jet
606 //
607     case kHijing_gj25:
608         comment = comment.Append("HIJING Gamma 25 GeV");
609         AliGenHijing *gener = HijingStandard();
610 // impact parameter range
611         gener->SetImpactParameterRange(0., 5.);
612         // trigger
613         gener->SetTrigger(2);
614         gener->SetPtJet(25.);
615         gener->SetRadiation(isw);
616         gener->SetSimpleJets(!isw);
617         gener->SetJetEtaRange(-0.12, 0.12);
618         gener->SetJetPhiRange(220., 320.);
619         break;
620
621     case kHijing_gj50:
622         comment = comment.Append("HIJING Gamma 50 GeV");
623         AliGenHijing *gener = HijingStandard();
624 // impact parameter range
625         gener->SetImpactParameterRange(0., 5.);
626         // trigger
627         gener->SetTrigger(2);
628         gener->SetPtJet(50.);
629         gener->SetRadiation(isw);
630         gener->SetSimpleJets(!isw);
631         gener->SetJetEtaRange(-0.12, 0.12);
632         gener->SetJetPhiRange(220., 320.);
633         break;
634
635     case kHijing_gj75:
636         comment = comment.Append("HIJING Gamma 75 GeV");
637         AliGenHijing *gener = HijingStandard();
638 // impact parameter range
639         gener->SetImpactParameterRange(0., 5.);
640         // trigger
641         gener->SetTrigger(2);
642         gener->SetPtJet(75.);
643         gener->SetRadiation(isw);
644         gener->SetSimpleJets(!isw);
645         gener->SetJetEtaRange(-0.12, 0.12);
646         gener->SetJetPhiRange(220., 320.);
647         break;
648
649     case kHijing_gj100:
650         comment = comment.Append("HIJING Gamma 100 GeV");
651         AliGenHijing *gener = HijingStandard();
652 // impact parameter range
653         gener->SetImpactParameterRange(0., 5.);
654         // trigger
655         gener->SetTrigger(2);
656         gener->SetPtJet(100.);
657         gener->SetRadiation(isw);
658         gener->SetSimpleJets(!isw);
659         gener->SetJetEtaRange(-0.12, 0.12);
660         gener->SetJetPhiRange(220., 320.);
661         break;
662
663     case kHijing_gj200:
664         comment = comment.Append("HIJING Gamma 200 GeV");
665         AliGenHijing *gener = HijingStandard();
666 // impact parameter range
667         gener->SetImpactParameterRange(0., 5.);
668         // trigger
669         gener->SetTrigger(2);
670         gener->SetPtJet(200.);
671         gener->SetRadiation(isw);
672         gener->SetSimpleJets(!isw);
673         gener->SetJetEtaRange(-0.12, 0.12);
674         gener->SetJetPhiRange(220., 320.);
675         break;
676     }
677     return gener;
678 }
679
680 AliGenHijing* HijingStandard()
681 {
682     AliGenHijing *gener = new AliGenHijing(-1);
683 // centre of mass energy 
684     gener->SetEnergyCMS(5500.);
685 // reference frame
686     gener->SetReferenceFrame("CMS");
687 // projectile
688      gener->SetProjectile("A", 208, 82);
689      gener->SetTarget    ("A", 208, 82);
690 // tell hijing to keep the full parent child chain
691      gener->KeepFullEvent();
692 // enable jet quenching
693      gener->SetJetQuenching(4);
694 // enable shadowing
695      gener->SetShadowing(1);
696 // neutral pion and heavy particle decays switched off
697      gener->SetDecaysOff(1);
698 // Don't track spectators
699      gener->SetSpectators(0);
700 // kinematic selection
701      gener->SetSelectAll(0);
702      return gener;
703 }
704
705
706