]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/ConfigPPR.C
Update of input distributions.
[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     // libraries required by geant321
51     gSystem->Load("libgeant321");
52
53     new TGeant3("C++ Interface to Geant3");
54
55     if (!gSystem->Getenv("CONFIG_FILE"))
56     {
57         TFile  *rootfile = new TFile("galice.root", "recreate");
58
59         rootfile->SetCompressionLevel(2);
60     }
61
62     TGeant3 *geant3 = (TGeant3 *) gMC;
63
64     //
65     // Set External decayer
66     AliDecayer *decayer = new AliDecayerPythia();
67
68     decayer->SetForceDecay(kAll);
69     decayer->Init();
70     gMC->SetExternalDecayer(decayer);
71     //
72     //
73     //=======================================================================
74     // ******* GEANT STEERING parameters FOR ALICE SIMULATION *******
75     geant3->SetTRIG(1);         //Number of events to be processed 
76     geant3->SetSWIT(4, 10);
77     geant3->SetDEBU(0, 0, 1);
78     //geant3->SetSWIT(2,2);
79     geant3->SetDCAY(1);
80     geant3->SetPAIR(1);
81     geant3->SetCOMP(1);
82     geant3->SetPHOT(1);
83     geant3->SetPFIS(0);
84     geant3->SetDRAY(0);
85     geant3->SetANNI(1);
86     geant3->SetBREM(1);
87     geant3->SetMUNU(1);
88     geant3->SetCKOV(1);
89     geant3->SetHADR(1);         //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3)
90     geant3->SetLOSS(2);
91     geant3->SetMULS(1);
92     geant3->SetRAYL(1);
93     geant3->SetAUTO(1);         //Select automatic STMIN etc... calc. (AUTO 1) or manual (AUTO 0)
94     geant3->SetABAN(0);         //Restore 3.16 behaviour for abandoned tracks
95     geant3->SetOPTI(2);         //Select optimisation level for GEANT geometry searches (0,1,2)
96     geant3->SetERAN(5.e-7);
97
98     Float_t cut = 1.e-3;        // 1MeV cut by default
99     Float_t tofmax = 1.e10;
100
101     //             GAM ELEC NHAD CHAD MUON EBREM MUHAB EDEL MUDEL MUPA TOFMAX
102     geant3->SetCUTS(cut, cut, cut, cut, cut, cut, cut, cut, cut, cut,
103                     tofmax);
104     //
105     //=======================================================================
106     // ************* STEERING parameters FOR ALICE SIMULATION **************
107     // --- Specify event type to be tracked through the ALICE setup
108     // --- All positions are in cm, angles in degrees, and P and E in GeV
109
110
111 // Generator Configuration
112     gAlice->SetDebug(1);
113     AliGenerator* gener = GeneratorFactory(run);
114     gener->SetOrigin(0, 0, 0);    // vertex position
115     gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
116     gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
117     gener->SetVertexSmear(kPerEvent); 
118     gener->SetTrackingFlag(1);
119     gener->Init();
120     
121     if (mag == k2kG) {
122         comment = comment.Append(" | L3 field 0.2 T");
123     } else if (mag == k4kG) {
124         comment = comment.Append(" | L3 field 0.4 T");
125     } else if (mag == k5kG) {
126         comment = comment.Append(" | L3 field 0.5 T");
127     }
128     
129     
130     if (rad == kGluonRadiation)
131     {
132         comment = comment.Append(" | Gluon Radiation On");
133         
134     } else {
135         comment = comment.Append(" | Gluon Radiation Off");
136     }
137
138     if (geo == kHoles)
139     {
140         comment = comment.Append(" | Holes for PHOS/RICH");
141         
142     } else {
143         comment = comment.Append(" | No holes for PHOS/RICH");
144     }
145
146     printf("\n \n Comment: %s \n \n", (char*) comment);
147     
148     
149 // Field (L3 0.4 T)
150     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag);
151     rootfile->cd();
152     gAlice->SetField(field);    
153     
154 //
155     Int_t   iABSO   = 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  = 1;
174     Int_t   iVZERO  = 1;
175     Int_t   iCRT    = 0;
176
177     //=================== Alice BODY parameters =============================
178     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
179
180
181     if (iMAG)
182     {
183         //=================== MAG parameters ============================
184         // --- Start with Magnet since detector layouts may be depending ---
185         // --- on the selected Magnet dimensions ---
186         AliMAG *MAG = new AliMAG("MAG", "Magnet");
187     }
188
189
190     if (iABSO)
191     {
192         //=================== ABSO parameters ============================
193         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
194     }
195
196     if (iDIPO)
197     {
198         //=================== DIPO parameters ============================
199
200         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
201     }
202
203     if (iHALL)
204     {
205         //=================== HALL parameters ============================
206
207         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
208     }
209
210
211     if (iFRAME)
212     {
213         //=================== FRAME parameters ============================
214
215         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
216         if (geo == kHoles) {
217             FRAME->SetHoles(1);
218         } else {
219             FRAME->SetHoles(0);
220         }
221     }
222
223     if (iSHIL)
224     {
225         //=================== SHIL parameters ============================
226
227         AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
228     }
229
230
231     if (iPIPE)
232     {
233         //=================== PIPE parameters ============================
234
235         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
236     }
237  
238     if(iITS) {
239
240     //=================== ITS parameters ============================
241     //
242     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
243     // almost all other detectors. This involves the fact that the ITS geometry
244     // still has several options to be followed in parallel in order to determine
245     // the best set-up which minimizes the induced background. All the geometries
246     // available to date are described in the following. Read carefully the comments
247     // and use the default version (the only one uncommented) unless you are making
248     // comparisons and you know what you are doing. In this case just uncomment the
249     // ITS geometry you want to use and run Aliroot.
250     //
251     // Detailed geometries:         
252     //
253     //
254     //AliITS *ITS  = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
255     //
256     //AliITS *ITS  = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
257     //
258         AliITSvPPRasymm *ITS  = new AliITSvPPRasymm("ITS","New ITS PPR detailed version with asymmetric services");
259         ITS->SetMinorVersion(2);                                         // don't touch this parameter if you're not an ITS developer
260         ITS->SetReadDet(kFALSE);                                         // don't touch this parameter if you're not an ITS developer
261     //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
262         ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
263         ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
264         ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
265         ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
266         ITS->SetRails(0);            // 1 --> rails in ; 0 --> rails out
267         ITS->SetCoolingFluid(1);   // 1 --> water ; 0 --> freon
268         //
269     //AliITSvPPRsymm *ITS  = new AliITSvPPRsymm("ITS","New ITS PPR detailed version with symmetric services");
270     //ITS->SetMinorVersion(2);                                       // don't touch this parameter if you're not an ITS developer
271     //ITS->SetReadDet(kFALSE);                                       // don't touch this parameter if you're not an ITS developer
272     //ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRsymm2.det"); // don't touch this parameter if you're not an ITS developer
273     //ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
274     //ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
275     //ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
276     //ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
277     //ITS->SetRails(0);              // 1 --> rails in ; 0 --> rails out
278     //ITS->SetCoolingFluid(1);       // 1 --> water ; 0 --> freon
279     //
280     //
281     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
282     // for reconstruction !):
283     //                                                     
284     //
285     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
286     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
287     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
288     //
289     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
290     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
291     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
292     //                      
293     //
294     //
295     // Geant3 <-> EUCLID conversion
296     // ============================
297     //
298     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
299     // media to two ASCII files (called by default ITSgeometry.euc and
300     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
301     // The default (=0) means that you dont want to use this facility.
302     //
303         ITS->SetEUCLID(0);  
304     }
305
306     if (iTPC)
307     {
308         //============================ TPC parameters ================================
309         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
310         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
311         // --- sectors are specified, any value other than that requires at least one 
312         // --- sector (lower or upper)to be specified!
313         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
314         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
315         // --- SecLows - number of lower sectors specified (up to 6)
316         // --- SecUps - number of upper sectors specified (up to 12)
317         // --- Sens - sensitive strips for the Slow Simulator !!!
318         // --- This does NOT work if all S or L-sectors are specified, i.e.
319         // --- if SecAL or SecAU < 0
320         //
321         //
322         //-----------------------------------------------------------------------------
323
324         //  gROOT->LoadMacro("SetTPCParam.C");
325         //  AliTPCParam *param = SetTPCParam();
326         AliTPC *TPC = new AliTPCv2("TPC", "Default");
327
328         // All sectors included 
329         TPC->SetSecAL(-1);
330         TPC->SetSecAU(-1);
331
332     }
333
334
335     if (iTOF) {
336         if (geo == kHoles) {
337         //=================== TOF parameters ============================
338             AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
339         } else {
340             AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
341         }
342     }
343
344
345     if (iRICH)
346     {
347         //=================== RICH parameters ===========================
348         AliRICH *RICH = new AliRICHv3("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 (iTRD)
361     {
362         //=================== TRD parameters ============================
363
364         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
365
366         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
367         TRD->SetGasMix(1);
368         if (geo == kHoles) {
369             // With hole in front of PHOS
370             TRD->SetPHOShole();
371             // With hole in front of RICH
372             TRD->SetRICHhole();
373         }
374             // Switch on TR
375             AliTRDsim *TRDsim = TRD->CreateTR();
376     }
377
378     if (iFMD)
379     {
380         //=================== FMD parameters ============================
381         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
382         FMD->SetRingsSi1(256);
383         FMD->SetRingsSi2(128);
384         FMD->SetSectorsSi1(20);
385         FMD->SetSectorsSi2(40);      
386    }
387
388     if (iMUON)
389     {
390         //=================== MUON parameters ===========================
391
392         AliMUON *MUON = new AliMUONv1("MUON", "default");
393     }
394     //=================== PHOS parameters ===========================
395
396     if (iPHOS)
397     {
398         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
399     }
400
401
402     if (iPMD)
403     {
404         //=================== PMD parameters ============================
405         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
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->SetPtJet(25.);
546         gener->SetRadiation(isw);
547         gener->SetSimpleJets(!isw);
548         gener->SetJetEtaRange(-0.3,0.3);
549         gener->SetJetPhiRange(75., 165.);   
550         break;
551
552     case kHijing_jj50:
553         comment = comment.Append("HIJING Jet 50 GeV");
554         AliGenHijing *gener = HijingStandard();
555 // impact parameter range
556         gener->SetImpactParameterRange(0., 5.);
557         // trigger
558         gener->SetTrigger(1);
559         gener->SetPtJet(50.);
560         gener->SetRadiation(isw);
561         gener->SetSimpleJets(!isw);
562         gener->SetJetEtaRange(-0.3,0.3);
563         gener->SetJetPhiRange(75., 165.);   
564         break;
565
566     case kHijing_jj75:
567         comment = comment.Append("HIJING Jet 75 GeV");
568         AliGenHijing *gener = HijingStandard();
569 // impact parameter range
570         gener->SetImpactParameterRange(0., 5.);
571         // trigger
572         gener->SetTrigger(1);
573         gener->SetPtJet(75.);
574         gener->SetRadiation(isw);
575         gener->SetSimpleJets(!isw);
576         gener->SetJetEtaRange(-0.3,0.3);
577         gener->SetJetPhiRange(75., 165.);   
578         break;
579
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->SetPtJet(100.);
588         gener->SetRadiation(isw);
589         gener->SetSimpleJets(!isw);
590         gener->SetJetEtaRange(-0.3,0.3);
591         gener->SetJetPhiRange(75., 165.);   
592         break;
593
594     case kHijing_jj200:
595         comment = comment.Append("HIJING Jet 200 GeV");
596         AliGenHijing *gener = HijingStandard();
597 // impact parameter range
598         gener->SetImpactParameterRange(0., 5.);
599         // trigger
600         gener->SetTrigger(1);
601         gener->SetPtJet(200.);
602         gener->SetRadiation(isw);
603         gener->SetSimpleJets(!isw);
604         gener->SetJetEtaRange(-0.3,0.3);
605         gener->SetJetPhiRange(75., 165.);   
606         break;
607 //
608 // Gamma-Jet
609 //
610     case kHijing_gj25:
611         comment = comment.Append("HIJING Gamma 25 GeV");
612         AliGenHijing *gener = HijingStandard();
613 // impact parameter range
614         gener->SetImpactParameterRange(0., 5.);
615         // trigger
616         gener->SetTrigger(2);
617         gener->SetPtJet(25.);
618         gener->SetRadiation(isw);
619         gener->SetSimpleJets(!isw);
620         gener->SetJetEtaRange(-0.12, 0.12);
621         gener->SetJetPhiRange(220., 320.);
622         break;
623
624     case kHijing_gj50:
625         comment = comment.Append("HIJING Gamma 50 GeV");
626         AliGenHijing *gener = HijingStandard();
627 // impact parameter range
628         gener->SetImpactParameterRange(0., 5.);
629         // trigger
630         gener->SetTrigger(2);
631         gener->SetPtJet(50.);
632         gener->SetRadiation(isw);
633         gener->SetSimpleJets(!isw);
634         gener->SetJetEtaRange(-0.12, 0.12);
635         gener->SetJetPhiRange(220., 320.);
636         break;
637
638     case kHijing_gj75:
639         comment = comment.Append("HIJING Gamma 75 GeV");
640         AliGenHijing *gener = HijingStandard();
641 // impact parameter range
642         gener->SetImpactParameterRange(0., 5.);
643         // trigger
644         gener->SetTrigger(2);
645         gener->SetPtJet(75.);
646         gener->SetRadiation(isw);
647         gener->SetSimpleJets(!isw);
648         gener->SetJetEtaRange(-0.12, 0.12);
649         gener->SetJetPhiRange(220., 320.);
650         break;
651
652     case kHijing_gj100:
653         comment = comment.Append("HIJING Gamma 100 GeV");
654         AliGenHijing *gener = HijingStandard();
655 // impact parameter range
656         gener->SetImpactParameterRange(0., 5.);
657         // trigger
658         gener->SetTrigger(2);
659         gener->SetPtJet(100.);
660         gener->SetRadiation(isw);
661         gener->SetSimpleJets(!isw);
662         gener->SetJetEtaRange(-0.12, 0.12);
663         gener->SetJetPhiRange(220., 320.);
664         break;
665
666     case kHijing_gj200:
667         comment = comment.Append("HIJING Gamma 200 GeV");
668         AliGenHijing *gener = HijingStandard();
669 // impact parameter range
670         gener->SetImpactParameterRange(0., 5.);
671         // trigger
672         gener->SetTrigger(2);
673         gener->SetPtJet(200.);
674         gener->SetRadiation(isw);
675         gener->SetSimpleJets(!isw);
676         gener->SetJetEtaRange(-0.12, 0.12);
677         gener->SetJetPhiRange(220., 320.);
678         break;
679     }
680     return gener;
681 }
682
683 AliGenHijing* HijingStandard()
684 {
685     AliGenHijing *gener = new AliGenHijing(-1);
686 // centre of mass energy 
687     gener->SetEnergyCMS(5500.);
688 // reference frame
689     gener->SetReferenceFrame("CMS");
690 // projectile
691      gener->SetProjectile("A", 208, 82);
692      gener->SetTarget    ("A", 208, 82);
693 // tell hijing to keep the full parent child chain
694      gener->KeepFullEvent();
695 // enable jet quenching
696      gener->SetJetQuenching(4);
697 // enable shadowing
698      gener->SetShadowing(1);
699 // neutral pion and heavy particle decays switched off
700      gener->SetDecaysOff(1);
701 // Don't track spectators
702      gener->SetSpectators(0);
703 // kinematic selection
704      gener->SetSelectAll(0);
705      return gener;
706 }
707
708
709