]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/ConfigPPR.C
Obsolete FMD parameters removed
[u/mrichter/AliRoot.git] / macros / ConfigPPR.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <Riostream.h>
3 #include <TRandom.h>
4 #include <TSystem.h>
5 #include <TVirtualMC.h>
6 #include <TGeant3.h>
7 #include "STEER/AliRunLoader.h"
8 #include "STEER/AliRun.h"
9 #include "STEER/AliConfig.h"
10 #include "STEER/AliGenerator.h"
11 #include "PYTHIA6/AliDecayerPythia.h"
12 #include "EVGEN/AliGenHIJINGpara.h"
13 #include "THijing/AliGenHijing.h"
14 #include "EVGEN/AliGenCocktail.h"
15 #include "EVGEN/AliGenSlowNucleons.h"
16 #include "EVGEN/AliSlowNucleonModelExp.h"
17 #include "PYTHIA6/AliGenPythia.h"
18 #include "STEER/AliMagFMaps.h"
19 #include "STRUCT/AliBODY.h"
20 #include "STRUCT/AliMAG.h"
21 #include "STRUCT/AliABSOv0.h"
22 #include "STRUCT/AliDIPOv2.h"
23 #include "STRUCT/AliHALL.h"
24 #include "STRUCT/AliFRAMEv2.h"
25 #include "STRUCT/AliSHILv2.h"
26 #include "STRUCT/AliPIPEv0.h"
27 #include "ITS/AliITSvPPRasymm.h"
28 #include "TPC/AliTPCv2.h"
29 #include "TOF/AliTOFv2FHoles.h"
30 #include "TOF/AliTOFv4T0.h"
31 #include "RICH/AliRICHv3.h"
32 #include "ZDC/AliZDCv2.h"
33 #include "TRD/AliTRDv1.h"
34 #include "FMD/AliFMDv1.h"
35 #include "MUON/AliMUONv1.h"
36 #include "PHOS/AliPHOSv1.h"
37 #include "PMD/AliPMDv1.h"
38 #include "START/AliSTARTv1.h"
39 #include "EMCAL/AliEMCALv1.h"
40 #include "CRT/AliCRTv0.h"
41 #include "VZERO/AliVZEROv2.h"
42 #endif
43
44 enum PprRun_t 
45 {
46     test50,
47     kParam_8000,   kParam_4000,  kParam_2000, 
48     kHijing_cent1, kHijing_cent2, 
49     kHijing_per1,  kHijing_per2, kHijing_per3, kHijing_per4,  kHijing_per5,
50     kHijing_jj25,  kHijing_jj50, kHijing_jj75, kHijing_jj100, kHijing_jj200, 
51     kHijing_gj25,  kHijing_gj50, kHijing_gj75, kHijing_gj100, kHijing_gj200,
52     kHijing_pA, kPythia
53 };
54
55 enum PprGeo_t 
56 {
57     kHoles, kNoHoles
58 };
59
60 enum PprRad_t
61 {
62     kGluonRadiation, kNoGluonRadiation
63 };
64
65 enum PprMag_t
66 {
67     k2kG, k4kG, k5kG
68 };
69
70
71 // This part for configuration    
72 //static PprRun_t srun = test50;
73 static PprRun_t srun = kPythia;
74 static PprGeo_t sgeo = kHoles;
75 static PprRad_t srad = kGluonRadiation;
76 static PprMag_t smag = k4kG;
77
78 // Comment line 
79 static TString  comment;
80
81 // Functions
82 Float_t EtaToTheta(Float_t arg);
83 AliGenerator* GeneratorFactory(PprRun_t srun);
84 AliGenHijing* HijingStandard();
85
86 void Config()
87 {
88     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
89     // Theta range given through pseudorapidity limits 22/6/2001
90
91     // Set Random Number seed
92   gRandom->SetSeed(12345); //Set 0 to use the current time
93     cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
94
95
96    // libraries required by geant321
97 #if defined(__CINT__)
98     gSystem->Load("libgeant321");
99 #endif
100
101     new     TGeant3("C++ Interface to Geant3");
102
103     AliRunLoader* rl=0x0;
104
105     cout<<"Config.C: Creating Run Loader ..."<<endl;
106     rl = AliRunLoader::Open("galice.root",
107                             AliConfig::fgkDefaultEventFolderName,
108                             "recreate");
109     if (rl == 0x0)
110       {
111         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
112         return;
113       }
114     rl->SetCompressionLevel(2);
115     rl->SetNumberOfEventsPerFile(3);
116     gAlice->SetRunLoader(rl);
117
118     //
119     // Set External decayer
120     AliDecayer *decayer = new AliDecayerPythia();
121
122     decayer->SetForceDecay(kAll);
123     decayer->Init();
124     gMC->SetExternalDecayer(decayer);
125     //
126     //
127     //=======================================================================
128     //
129     //=======================================================================
130     // ************* STEERING parameters FOR ALICE SIMULATION **************
131     // --- Specify event type to be tracked through the ALICE setup
132     // --- All positions are in cm, angles in degrees, and P and E in GeV
133
134     gMC->SetProcess("DCAY",1);
135     gMC->SetProcess("PAIR",1);
136     gMC->SetProcess("COMP",1);
137     gMC->SetProcess("PHOT",1);
138     gMC->SetProcess("PFIS",0);
139     gMC->SetProcess("DRAY",0);
140     gMC->SetProcess("ANNI",1);
141     gMC->SetProcess("BREM",1);
142     gMC->SetProcess("MUNU",1);
143     gMC->SetProcess("CKOV",1);
144     gMC->SetProcess("HADR",1);
145     gMC->SetProcess("LOSS",2);
146     gMC->SetProcess("MULS",1);
147     gMC->SetProcess("RAYL",1);
148
149     Float_t cut = 1.e-3;        // 1MeV cut by default
150     Float_t tofmax = 1.e10;
151
152     gMC->SetCut("CUTGAM", cut);
153     gMC->SetCut("CUTELE", cut);
154     gMC->SetCut("CUTNEU", cut);
155     gMC->SetCut("CUTHAD", cut);
156     gMC->SetCut("CUTMUO", cut);
157     gMC->SetCut("BCUTE",  cut); 
158     gMC->SetCut("BCUTM",  cut); 
159     gMC->SetCut("DCUTE",  cut); 
160     gMC->SetCut("DCUTM",  cut); 
161     gMC->SetCut("PPCUTM", cut);
162     gMC->SetCut("TOFMAX", tofmax); 
163
164
165 // Generator Configuration
166     gAlice->SetDebug(1);
167     AliGenerator* gener = GeneratorFactory(srun);
168     gener->SetOrigin(0, 0, 0);    // vertex position
169     gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
170     gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
171     gener->SetVertexSmear(kPerEvent); 
172     gener->SetTrackingFlag(1);
173     gener->Init();
174     
175     if (smag == k2kG) {
176         comment = comment.Append(" | L3 field 0.2 T");
177     } else if (smag == k4kG) {
178         comment = comment.Append(" | L3 field 0.4 T");
179     } else if (smag == k5kG) {
180         comment = comment.Append(" | L3 field 0.5 T");
181     }
182     
183     
184     if (srad == kGluonRadiation)
185     {
186         comment = comment.Append(" | Gluon Radiation On");
187         
188     } else {
189         comment = comment.Append(" | Gluon Radiation Off");
190     }
191
192     if (sgeo == kHoles)
193     {
194         comment = comment.Append(" | Holes for PHOS/RICH");
195         
196     } else {
197         comment = comment.Append(" | No holes for PHOS/RICH");
198     }
199
200     printf("\n \n Comment: %s \n \n", comment.Data());
201     
202     
203 // Field (L3 0.4 T)
204     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., smag);
205     rl->CdGAFile();
206     
207 //
208     Int_t   iABSO   = 1;
209     Int_t   iDIPO   = 1;
210     Int_t   iFMD    = 1;
211     Int_t   iFRAME  = 1;
212     Int_t   iHALL   = 1;
213     Int_t   iITS    = 1;
214     Int_t   iMAG    = 1;
215     Int_t   iMUON   = 1;
216     Int_t   iPHOS   = 1;
217     Int_t   iPIPE   = 1;
218     Int_t   iPMD    = 1;
219     Int_t   iRICH   = 1;
220     Int_t   iSHIL   = 1;
221     Int_t   iSTART  = 1;
222     Int_t   iTOF    = 1;
223     Int_t   iTPC    = 1;
224     Int_t   iTRD    = 1;
225     Int_t   iZDC    = 1;
226     Int_t   iEMCAL  = 1;
227     Int_t   iVZERO  = 1;
228     Int_t   iCRT    = 0;
229
230     //=================== Alice BODY parameters =============================
231     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
232
233
234     if (iMAG)
235     {
236         //=================== MAG parameters ============================
237         // --- Start with Magnet since detector layouts may be depending ---
238         // --- on the selected Magnet dimensions ---
239         AliMAG *MAG = new AliMAG("MAG", "Magnet");
240     }
241
242
243     if (iABSO)
244     {
245         //=================== ABSO parameters ============================
246         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
247     }
248
249     if (iDIPO)
250     {
251         //=================== DIPO parameters ============================
252
253         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
254     }
255
256     if (iHALL)
257     {
258         //=================== HALL parameters ============================
259
260         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
261     }
262
263
264     if (iFRAME)
265     {
266         //=================== FRAME parameters ============================
267
268         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
269         if (sgeo == kHoles) {
270             FRAME->SetHoles(1);
271         } else {
272             FRAME->SetHoles(0);
273         }
274     }
275
276     if (iSHIL)
277     {
278         //=================== SHIL parameters ============================
279
280         AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
281     }
282
283
284     if (iPIPE)
285     {
286         //=================== PIPE parameters ============================
287
288         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
289     }
290  
291     if(iITS) {
292
293     //=================== ITS parameters ============================
294     //
295     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
296     // almost all other detectors. This involves the fact that the ITS geometry
297     // still has several options to be followed in parallel in order to determine
298     // the best set-up which minimizes the induced background. All the geometries
299     // available to date are described in the following. Read carefully the comments
300     // and use the default version (the only one uncommented) unless you are making
301     // comparisons and you know what you are doing. In this case just uncomment the
302     // ITS geometry you want to use and run Aliroot.
303     //
304     // Detailed geometries:         
305     //
306     //
307     //AliITS *ITS  = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
308     //
309     //AliITS *ITS  = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
310     //
311         AliITSvPPRasymm *ITS  = new AliITSvPPRasymm("ITS","New ITS PPR detailed version with asymmetric services");
312         ITS->SetMinorVersion(2);                                         // don't touch this parameter if you're not an ITS developer
313         ITS->SetReadDet(kFALSE);                                         // don't touch this parameter if you're not an ITS developer
314     //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
315         ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
316         ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
317         ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
318         ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
319         ITS->SetRails(0);            // 1 --> rails in ; 0 --> rails out
320         ITS->SetCoolingFluid(1);   // 1 --> water ; 0 --> freon
321         //
322     //AliITSvPPRsymm *ITS  = new AliITSvPPRsymm("ITS","New ITS PPR detailed version with symmetric services");
323     //ITS->SetMinorVersion(2);                                       // don't touch this parameter if you're not an ITS developer
324     //ITS->SetReadDet(kFALSE);                                       // don't touch this parameter if you're not an ITS developer
325     //ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRsymm2.det"); // don't touch this parameter if you're not an ITS developer
326     //ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
327     //ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
328     //ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
329     //ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
330     //ITS->SetRails(0);              // 1 --> rails in ; 0 --> rails out
331     //ITS->SetCoolingFluid(1);       // 1 --> water ; 0 --> freon
332     //
333     //
334     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
335     // for reconstruction !):
336     //                                                     
337     //
338     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
339     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
340     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
341     //
342     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
343     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
344     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
345     //                      
346     //
347     //
348     // Geant3 <-> EUCLID conversion
349     // ============================
350     //
351     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
352     // media to two ASCII files (called by default ITSgeometry.euc and
353     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
354     // The default (=0) means that you dont want to use this facility.
355     //
356         ITS->SetEUCLID(0);  
357     }
358
359     if (iTPC)
360     {
361         //============================ TPC parameters ================================
362         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
363         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
364         // --- sectors are specified, any value other than that requires at least one 
365         // --- sector (lower or upper)to be specified!
366         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
367         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
368         // --- SecLows - number of lower sectors specified (up to 6)
369         // --- SecUps - number of upper sectors specified (up to 12)
370         // --- Sens - sensitive strips for the Slow Simulator !!!
371         // --- This does NOT work if all S or L-sectors are specified, i.e.
372         // --- if SecAL or SecAU < 0
373         //
374         //
375         //-----------------------------------------------------------------------------
376
377         //  gROOT->LoadMacro("SetTPCParam.C");
378         //  AliTPCParam *param = SetTPCParam();
379         AliTPC *TPC = new AliTPCv2("TPC", "Default");
380
381         // All sectors included 
382         TPC->SetSecAL(-1);
383         TPC->SetSecAU(-1);
384
385     }
386
387
388     if (iTOF) {
389         if (sgeo == kHoles) {
390         //=================== TOF parameters ============================
391             AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
392         } else {
393             AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
394         }
395     }
396
397
398     if (iRICH)
399     {
400         //=================== RICH parameters ===========================
401         AliRICH *RICH = new AliRICHv3("RICH", "normal RICH");
402
403     }
404
405
406     if (iZDC)
407     {
408         //=================== ZDC parameters ============================
409
410         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
411     }
412
413     if (iTRD)
414     {
415         //=================== TRD parameters ============================
416
417         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
418
419         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
420         TRD->SetGasMix(1);
421         if (sgeo == kHoles) {
422             // With hole in front of PHOS
423             TRD->SetPHOShole();
424             // With hole in front of RICH
425             TRD->SetRICHhole();
426         }
427             // Switch on TR
428             AliTRDsim *TRDsim = TRD->CreateTR();
429     }
430
431     if (iFMD)
432     {
433         //=================== FMD parameters ============================
434         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
435    }
436
437     if (iMUON)
438     {
439         //=================== MUON parameters ===========================
440
441         AliMUON *MUON = new AliMUONv1("MUON", "default");
442     }
443     //=================== PHOS parameters ===========================
444
445     if (iPHOS)
446     {
447         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
448     }
449
450
451     if (iPMD)
452     {
453         //=================== PMD parameters ============================
454         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
455     }
456
457     if (iSTART)
458     {
459         //=================== START parameters ============================
460         AliSTART *START = new AliSTARTv1("START", "START Detector");
461     }
462
463     if (iEMCAL)
464     {
465         //=================== EMCAL parameters ============================
466         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "G56_2_55_19");
467     }
468
469      if (iCRT)
470     {
471         //=================== CRT parameters ============================
472         AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
473     }
474
475      if (iVZERO)
476     {
477         //=================== CRT parameters ============================
478         AliVZERO *VZERO = new AliVZEROv2("VZERO", "normal VZERO");
479     }
480  
481              
482 }
483
484 Float_t EtaToTheta(Float_t arg){
485   return (180./TMath::Pi())*2.*atan(exp(-arg));
486 }
487
488
489
490 AliGenerator* GeneratorFactory(PprRun_t srun) {
491     Int_t isw = 3;
492     if (srad == kNoGluonRadiation) isw = 0;
493     
494
495     AliGenerator * gGener = 0x0;
496     switch (srun) {
497     case test50:
498       {
499         comment = comment.Append(":HIJINGparam test 50 particles");
500         AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
501         gener->SetMomentumRange(0, 999999.);
502         gener->SetPhiRange(0., 360.);
503         // Set pseudorapidity range from -8 to 8.
504         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
505         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
506         gener->SetThetaRange(thmin,thmax);
507         gGener=gener;
508       }
509         break;
510     case kParam_8000:
511       {
512         comment = comment.Append(":HIJINGparam N=8000");
513         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
514         gener->SetMomentumRange(0, 999999.);
515         gener->SetPhiRange(0., 360.);
516         // Set pseudorapidity range from -8 to 8.
517         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
518         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
519         gener->SetThetaRange(thmin,thmax);
520         gGener=gener;
521       }
522         break;
523     case kParam_4000:
524       {
525         comment = comment.Append("HIJINGparam N=4000");
526         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
527         gener->SetMomentumRange(0, 999999.);
528         gener->SetPhiRange(0., 360.);
529         // Set pseudorapidity range from -8 to 8.
530         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
531         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
532         gener->SetThetaRange(thmin,thmax);
533         gGener=gener;
534       }
535         break;
536     case kParam_2000:
537       {
538         comment = comment.Append("HIJINGparam N=2000");
539         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
540         gener->SetMomentumRange(0, 999999.);
541         gener->SetPhiRange(0., 360.);
542         // Set pseudorapidity range from -8 to 8.
543         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
544         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
545         gener->SetThetaRange(thmin,thmax);
546         gGener=gener;
547       }
548         break;
549 //
550 //  Hijing Central
551 //
552     case kHijing_cent1:
553       {
554         comment = comment.Append("HIJING cent1");
555         AliGenHijing *gener = HijingStandard();
556 // impact parameter range
557         gener->SetImpactParameterRange(0., 5.);
558         gGener=gener;
559       }
560         break;
561     case kHijing_cent2:
562       {
563         comment = comment.Append("HIJING cent2");
564         AliGenHijing *gener = HijingStandard();
565 // impact parameter range
566         gener->SetImpactParameterRange(0., 2.);
567         gGener=gener;
568         break;
569       }
570 //
571 // Hijing Peripheral 
572 //
573     case kHijing_per1:
574       {
575         comment = comment.Append("HIJING per1");
576         AliGenHijing *gener = HijingStandard();
577 // impact parameter range
578         gener->SetImpactParameterRange(5., 8.6);
579         gGener=gener;
580       }
581         break;
582     case kHijing_per2:
583       {
584         comment = comment.Append("HIJING per2");
585         AliGenHijing *gener = HijingStandard();
586 // impact parameter range
587         gener->SetImpactParameterRange(8.6, 11.2);
588         gGener=gener;
589       }
590         break;
591     case kHijing_per3:
592       {
593         comment = comment.Append("HIJING per3");
594         AliGenHijing *gener = HijingStandard();
595 // impact parameter range
596         gener->SetImpactParameterRange(11.2, 13.2);
597         gGener=gener;
598       }
599         break;
600     case kHijing_per4:
601       {
602         comment = comment.Append("HIJING per4");
603         AliGenHijing *gener = HijingStandard();
604 // impact parameter range
605         gener->SetImpactParameterRange(13.2, 15.);
606         gGener=gener;
607       }
608         break;
609     case kHijing_per5:
610       {
611         comment = comment.Append("HIJING per5");
612         AliGenHijing *gener = HijingStandard();
613 // impact parameter range
614         gener->SetImpactParameterRange(15., 100.);
615         gGener=gener;
616       }
617         break;
618 //
619 //  Jet-Jet
620 //
621     case kHijing_jj25:
622       {
623         comment = comment.Append("HIJING Jet 25 GeV");
624         AliGenHijing *gener = HijingStandard();
625 // impact parameter range
626         gener->SetImpactParameterRange(0., 5.);
627         // trigger
628         gener->SetTrigger(1);
629         gener->SetPtJet(25.);
630         gener->SetRadiation(isw);
631         gener->SetSimpleJets(!isw);
632         gener->SetJetEtaRange(-0.3,0.3);
633         gener->SetJetPhiRange(75., 165.);   
634         gGener=gener;
635       }
636         break;
637
638     case kHijing_jj50:
639       {
640         comment = comment.Append("HIJING Jet 50 GeV");
641         AliGenHijing *gener = HijingStandard();
642 // impact parameter range
643         gener->SetImpactParameterRange(0., 5.);
644         // trigger
645         gener->SetTrigger(1);
646         gener->SetPtJet(50.);
647         gener->SetRadiation(isw);
648         gener->SetSimpleJets(!isw);
649         gener->SetJetEtaRange(-0.3,0.3);
650         gener->SetJetPhiRange(75., 165.);   
651         gGener=gener;
652       }
653         break;
654
655     case kHijing_jj75:
656       {
657         comment = comment.Append("HIJING Jet 75 GeV");
658         AliGenHijing *gener = HijingStandard();
659 // impact parameter range
660         gener->SetImpactParameterRange(0., 5.);
661         // trigger
662         gener->SetTrigger(1);
663         gener->SetPtJet(75.);
664         gener->SetRadiation(isw);
665         gener->SetSimpleJets(!isw);
666         gener->SetJetEtaRange(-0.3,0.3);
667         gener->SetJetPhiRange(75., 165.);   
668         gGener=gener;
669       }
670         break;
671
672     case kHijing_jj100:
673       {
674         comment = comment.Append("HIJING Jet 100 GeV");
675         AliGenHijing *gener = HijingStandard();
676 // impact parameter range
677         gener->SetImpactParameterRange(0., 5.);
678         // trigger
679         gener->SetTrigger(1);
680         gener->SetPtJet(100.);
681         gener->SetRadiation(isw);
682         gener->SetSimpleJets(!isw);
683         gener->SetJetEtaRange(-0.3,0.3);
684         gener->SetJetPhiRange(75., 165.);   
685         gGener=gener;
686       }
687         break;
688
689     case kHijing_jj200:
690       {
691         comment = comment.Append("HIJING Jet 200 GeV");
692         AliGenHijing *gener = HijingStandard();
693 // impact parameter range
694         gener->SetImpactParameterRange(0., 5.);
695         // trigger
696         gener->SetTrigger(1);
697         gener->SetPtJet(200.);
698         gener->SetRadiation(isw);
699         gener->SetSimpleJets(!isw);
700         gener->SetJetEtaRange(-0.3,0.3);
701         gener->SetJetPhiRange(75., 165.);   
702         gGener=gener;
703       }
704         break;
705 //
706 // Gamma-Jet
707 //
708     case kHijing_gj25:
709       {
710         comment = comment.Append("HIJING Gamma 25 GeV");
711         AliGenHijing *gener = HijingStandard();
712 // impact parameter range
713         gener->SetImpactParameterRange(0., 5.);
714         // trigger
715         gener->SetTrigger(2);
716         gener->SetPtJet(25.);
717         gener->SetRadiation(isw);
718         gener->SetSimpleJets(!isw);
719         gener->SetJetEtaRange(-0.12, 0.12);
720         gener->SetJetPhiRange(220., 320.);
721         gGener=gener;
722       }
723         break;
724
725     case kHijing_gj50:
726       {
727         comment = comment.Append("HIJING Gamma 50 GeV");
728         AliGenHijing *gener = HijingStandard();
729 // impact parameter range
730         gener->SetImpactParameterRange(0., 5.);
731         // trigger
732         gener->SetTrigger(2);
733         gener->SetPtJet(50.);
734         gener->SetRadiation(isw);
735         gener->SetSimpleJets(!isw);
736         gener->SetJetEtaRange(-0.12, 0.12);
737         gener->SetJetPhiRange(220., 320.);
738         gGener=gener;
739       }
740         break;
741
742     case kHijing_gj75:
743       {
744         comment = comment.Append("HIJING Gamma 75 GeV");
745         AliGenHijing *gener = HijingStandard();
746 // impact parameter range
747         gener->SetImpactParameterRange(0., 5.);
748         // trigger
749         gener->SetTrigger(2);
750         gener->SetPtJet(75.);
751         gener->SetRadiation(isw);
752         gener->SetSimpleJets(!isw);
753         gener->SetJetEtaRange(-0.12, 0.12);
754         gener->SetJetPhiRange(220., 320.);
755         gGener=gener;
756       }
757         break;
758
759     case kHijing_gj100:
760       {
761         comment = comment.Append("HIJING Gamma 100 GeV");
762         AliGenHijing *gener = HijingStandard();
763 // impact parameter range
764         gener->SetImpactParameterRange(0., 5.);
765         // trigger
766         gener->SetTrigger(2);
767         gener->SetPtJet(100.);
768         gener->SetRadiation(isw);
769         gener->SetSimpleJets(!isw);
770         gener->SetJetEtaRange(-0.12, 0.12);
771         gener->SetJetPhiRange(220., 320.);
772         gGener=gener;
773       }
774         break;
775
776     case kHijing_gj200:
777       {
778         comment = comment.Append("HIJING Gamma 200 GeV");
779         AliGenHijing *gener = HijingStandard();
780 // impact parameter range
781         gener->SetImpactParameterRange(0., 5.);
782         // trigger
783         gener->SetTrigger(2);
784         gener->SetPtJet(200.);
785         gener->SetRadiation(isw);
786         gener->SetSimpleJets(!isw);
787         gener->SetJetEtaRange(-0.12, 0.12);
788         gener->SetJetPhiRange(220., 320.);
789         gGener=gener;
790       }
791         break;
792     case kHijing_pA:
793       {
794         comment = comment.Append("HIJING pA");
795
796         AliGenCocktail *gener  = new AliGenCocktail();
797 //      gener->SetTrackingFlag(0);
798         gener->SetOrigin(0, 0, 0);    // vertex position
799         gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
800         gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
801         gener->SetVertexSmear(kPerEvent); 
802         AliGenHijing   *hijing = new AliGenHijing(-1);
803 // centre of mass energy 
804         hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
805 // impact parameter range
806         hijing->SetImpactParameterRange(0., 15.);
807 // reference frame
808         hijing->SetReferenceFrame("CMS");
809         hijing->SetBoostLHC(1);
810 // projectile
811         hijing->SetProjectile("P", 1, 1);
812         hijing->SetTarget    ("A", 208, 82);
813 // tell hijing to keep the full parent child chain
814         hijing->KeepFullEvent();
815 // enable jet quenching
816         hijing->SetJetQuenching(0);
817 // enable shadowing
818         hijing->SetShadowing(1);
819 // Don't track spectators
820         hijing->SetSpectators(0);
821 // kinematic selection
822         hijing->SetSelectAll(0);
823 //
824         AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
825         AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
826         gray->SetSlowNucleonModel(model);
827         gray->SetDebug(1);
828         gener->AddGenerator(hijing,"Hijing pPb", 1);
829         gener->AddGenerator(gray,  "Gray Particles",1);
830         gGener=gener;
831       }
832         break;
833     case kPythia:
834       {
835         comment = comment.Append(":Pythia p-p @ 14 TeV");
836         AliGenPythia *gener = new AliGenPythia(-1); 
837         gener->SetMomentumRange(0,999999);
838         gener->SetPhiRange(-180,180);
839         gener->SetThetaRange(0., 180.);
840         gener->SetYRange(-12,12);
841         gener->SetPtRange(0,1000);
842         gener->SetStrucFunc(kCTEQ4L);   
843         gener->SetProcess(kPyMb);
844         gener->SetEnergyCMS(14000.);
845         gGener=gener;
846       }
847     break;
848     }
849     return gGener;
850 }
851
852 AliGenHijing* HijingStandard()
853 {
854     AliGenHijing *gener = new AliGenHijing(-1);
855 // centre of mass energy 
856     gener->SetEnergyCMS(5500.);
857 // reference frame
858     gener->SetReferenceFrame("CMS");
859 // projectile
860      gener->SetProjectile("A", 208, 82);
861      gener->SetTarget    ("A", 208, 82);
862 // tell hijing to keep the full parent child chain
863      gener->KeepFullEvent();
864 // enable jet quenching
865      gener->SetJetQuenching(4);
866 // enable shadowing
867      gener->SetShadowing(1);
868 // neutral pion and heavy particle decays switched off
869      gener->SetDecaysOff(1);
870 // Don't track spectators
871      gener->SetSpectators(0);
872 // kinematic selection
873      gener->SetSelectAll(0);
874      return gener;
875 }
876
877
878