]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/ConfigPPR.C
Change back to old sector numbering
[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/AliITSvPPRasymmFMD.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, kPythia6, kPythia6Jets
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 = kPythia6;
74 static PprGeo_t sgeo = kHoles;
75 static PprRad_t srad = kGluonRadiation;
76 static PprMag_t smag = k5kG;
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     field->SetL3ConstField(0); //Using const. field in the barrel
206     rl->CdGAFile();
207     gAlice->SetField(field);    
208 //
209     Int_t   iABSO   = 1;
210     Int_t   iDIPO   = 1;
211     Int_t   iFMD    = 1;
212     Int_t   iFRAME  = 1;
213     Int_t   iHALL   = 1;
214     Int_t   iITS    = 1;
215     Int_t   iMAG    = 1;
216     Int_t   iMUON   = 1;
217     Int_t   iPHOS   = 1;
218     Int_t   iPIPE   = 1;
219     Int_t   iPMD    = 1;
220     Int_t   iRICH   = 1;
221     Int_t   iSHIL   = 1;
222     Int_t   iSTART  = 1;
223     Int_t   iTOF    = 1;
224     Int_t   iTPC    = 1;
225     Int_t   iTRD    = 1;
226     Int_t   iZDC    = 1;
227     Int_t   iEMCAL  = 1;
228     Int_t   iVZERO  = 1;
229     Int_t   iCRT    = 0;
230
231     //=================== Alice BODY parameters =============================
232     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
233
234
235     if (iMAG)
236     {
237         //=================== MAG parameters ============================
238         // --- Start with Magnet since detector layouts may be depending ---
239         // --- on the selected Magnet dimensions ---
240         AliMAG *MAG = new AliMAG("MAG", "Magnet");
241     }
242
243
244     if (iABSO)
245     {
246         //=================== ABSO parameters ============================
247         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
248     }
249
250     if (iDIPO)
251     {
252         //=================== DIPO parameters ============================
253
254         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
255     }
256
257     if (iHALL)
258     {
259         //=================== HALL parameters ============================
260
261         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
262     }
263
264
265     if (iFRAME)
266     {
267         //=================== FRAME parameters ============================
268
269         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
270         if (sgeo == kHoles) {
271             FRAME->SetHoles(1);
272         } else {
273             FRAME->SetHoles(0);
274         }
275     }
276
277     if (iSHIL)
278     {
279         //=================== SHIL parameters ============================
280
281         AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
282     }
283
284
285     if (iPIPE)
286     {
287         //=================== PIPE parameters ============================
288
289         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
290     }
291  
292     if(iITS) {
293
294     //=================== ITS parameters ============================
295     //
296     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
297     // almost all other detectors. This involves the fact that the ITS geometry
298     // still has several options to be followed in parallel in order to determine
299     // the best set-up which minimizes the induced background. All the geometries
300     // available to date are described in the following. Read carefully the comments
301     // and use the default version (the only one uncommented) unless you are making
302     // comparisons and you know what you are doing. In this case just uncomment the
303     // ITS geometry you want to use and run Aliroot.
304     //
305     // Detailed geometries:         
306     //
307     //
308     //AliITS *ITS  = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
309     //
310     //AliITS *ITS  = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
311     //
312         AliITSvPPRasymmFMD *ITS  = new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version with asymmetric services");
313         ITS->SetMinorVersion(2);  // don't touch this parameter if you're not an ITS developer
314         ITS->SetReadDet(kTRUE);   // don't touch this parameter if you're not an ITS developer
315     //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
316         ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
317         ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
318         ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
319         ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
320         ITS->SetRails(0);            // 1 --> rails in ; 0 --> rails out
321         ITS->SetCoolingFluid(1);   // 1 --> water ; 0 --> freon
322
323     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
324     // for reconstruction !):
325     //                                                     
326     //
327     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
328     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
329     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
330     //
331     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
332     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
333     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
334     //                      
335     //
336     //
337     // Geant3 <-> EUCLID conversion
338     // ============================
339     //
340     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
341     // media to two ASCII files (called by default ITSgeometry.euc and
342     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
343     // The default (=0) means that you dont want to use this facility.
344     //
345         ITS->SetEUCLID(0);  
346     }
347
348     if (iTPC)
349     {
350         //============================ TPC parameters ================================
351         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
352         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
353         // --- sectors are specified, any value other than that requires at least one 
354         // --- sector (lower or upper)to be specified!
355         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
356         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
357         // --- SecLows - number of lower sectors specified (up to 6)
358         // --- SecUps - number of upper sectors specified (up to 12)
359         // --- Sens - sensitive strips for the Slow Simulator !!!
360         // --- This does NOT work if all S or L-sectors are specified, i.e.
361         // --- if SecAL or SecAU < 0
362         //
363         //
364         //-----------------------------------------------------------------------------
365
366         //  gROOT->LoadMacro("SetTPCParam.C");
367         //  AliTPCParam *param = SetTPCParam();
368         AliTPC *TPC = new AliTPCv2("TPC", "Default");
369
370         // All sectors included 
371         TPC->SetSecAL(-1);
372         TPC->SetSecAU(-1);
373
374     }
375
376
377     if (iTOF) {
378         if (sgeo == kHoles) {
379         //=================== TOF parameters ============================
380             AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
381         } else {
382             AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
383         }
384     }
385
386
387     if (iRICH)
388     {
389         //=================== RICH parameters ===========================
390         AliRICH *RICH = new AliRICHv3("RICH", "normal RICH");
391
392     }
393
394
395     if (iZDC)
396     {
397         //=================== ZDC parameters ============================
398
399         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
400     }
401
402     if (iTRD)
403     {
404         //=================== TRD parameters ============================
405
406         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
407
408         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
409         TRD->SetGasMix(1);
410         if (sgeo == kHoles) {
411             // With hole in front of PHOS
412             TRD->SetPHOShole();
413             // With hole in front of RICH
414             TRD->SetRICHhole();
415         }
416             // Switch on TR
417             AliTRDsim *TRDsim = TRD->CreateTR();
418     }
419
420     if (iFMD)
421     {
422         //=================== FMD parameters ============================
423         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
424    }
425
426     if (iMUON)
427     {
428         //=================== MUON parameters ===========================
429
430         AliMUON *MUON = new AliMUONv1("MUON", "default");
431     }
432     //=================== PHOS parameters ===========================
433
434     if (iPHOS)
435     {
436         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
437     }
438
439
440     if (iPMD)
441     {
442         //=================== PMD parameters ============================
443         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
444     }
445
446     if (iSTART)
447     {
448         //=================== START parameters ============================
449         AliSTART *START = new AliSTARTv1("START", "START Detector");
450     }
451
452     if (iEMCAL)
453     {
454         //=================== EMCAL parameters ============================
455         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "G56_2_55_19");
456     }
457
458      if (iCRT)
459     {
460         //=================== CRT parameters ============================
461         AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
462     }
463
464      if (iVZERO)
465     {
466         //=================== CRT parameters ============================
467         AliVZERO *VZERO = new AliVZEROv2("VZERO", "normal VZERO");
468     }
469  
470              
471 }
472
473 Float_t EtaToTheta(Float_t arg){
474   return (180./TMath::Pi())*2.*atan(exp(-arg));
475 }
476
477
478
479 AliGenerator* GeneratorFactory(PprRun_t srun) {
480     Int_t isw = 3;
481     if (srad == kNoGluonRadiation) isw = 0;
482     
483
484     AliGenerator * gGener = 0x0;
485     switch (srun) {
486     case test50:
487       {
488         comment = comment.Append(":HIJINGparam test 50 particles");
489         AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
490         gener->SetMomentumRange(0, 999999.);
491         gener->SetPhiRange(0., 360.);
492         // Set pseudorapidity range from -8 to 8.
493         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
494         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
495         gener->SetThetaRange(thmin,thmax);
496         gGener=gener;
497       }
498         break;
499     case kParam_8000:
500       {
501         comment = comment.Append(":HIJINGparam N=8000");
502         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
503         gener->SetMomentumRange(0, 999999.);
504         gener->SetPhiRange(0., 360.);
505         // Set pseudorapidity range from -8 to 8.
506         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
507         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
508         gener->SetThetaRange(thmin,thmax);
509         gGener=gener;
510       }
511         break;
512     case kParam_4000:
513       {
514         comment = comment.Append("HIJINGparam N=4000");
515         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
516         gener->SetMomentumRange(0, 999999.);
517         gener->SetPhiRange(0., 360.);
518         // Set pseudorapidity range from -8 to 8.
519         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
520         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
521         gener->SetThetaRange(thmin,thmax);
522         gGener=gener;
523       }
524         break;
525     case kParam_2000:
526       {
527         comment = comment.Append("HIJINGparam N=2000");
528         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
529         gener->SetMomentumRange(0, 999999.);
530         gener->SetPhiRange(0., 360.);
531         // Set pseudorapidity range from -8 to 8.
532         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
533         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
534         gener->SetThetaRange(thmin,thmax);
535         gGener=gener;
536       }
537         break;
538 //
539 //  Hijing Central
540 //
541     case kHijing_cent1:
542       {
543         comment = comment.Append("HIJING cent1");
544         AliGenHijing *gener = HijingStandard();
545 // impact parameter range
546         gener->SetImpactParameterRange(0., 5.);
547         gGener=gener;
548       }
549         break;
550     case kHijing_cent2:
551       {
552         comment = comment.Append("HIJING cent2");
553         AliGenHijing *gener = HijingStandard();
554 // impact parameter range
555         gener->SetImpactParameterRange(0., 2.);
556         gGener=gener;
557         break;
558       }
559 //
560 // Hijing Peripheral 
561 //
562     case kHijing_per1:
563       {
564         comment = comment.Append("HIJING per1");
565         AliGenHijing *gener = HijingStandard();
566 // impact parameter range
567         gener->SetImpactParameterRange(5., 8.6);
568         gGener=gener;
569       }
570         break;
571     case kHijing_per2:
572       {
573         comment = comment.Append("HIJING per2");
574         AliGenHijing *gener = HijingStandard();
575 // impact parameter range
576         gener->SetImpactParameterRange(8.6, 11.2);
577         gGener=gener;
578       }
579         break;
580     case kHijing_per3:
581       {
582         comment = comment.Append("HIJING per3");
583         AliGenHijing *gener = HijingStandard();
584 // impact parameter range
585         gener->SetImpactParameterRange(11.2, 13.2);
586         gGener=gener;
587       }
588         break;
589     case kHijing_per4:
590       {
591         comment = comment.Append("HIJING per4");
592         AliGenHijing *gener = HijingStandard();
593 // impact parameter range
594         gener->SetImpactParameterRange(13.2, 15.);
595         gGener=gener;
596       }
597         break;
598     case kHijing_per5:
599       {
600         comment = comment.Append("HIJING per5");
601         AliGenHijing *gener = HijingStandard();
602 // impact parameter range
603         gener->SetImpactParameterRange(15., 100.);
604         gGener=gener;
605       }
606         break;
607 //
608 //  Jet-Jet
609 //
610     case kHijing_jj25:
611       {
612         comment = comment.Append("HIJING Jet 25 GeV");
613         AliGenHijing *gener = HijingStandard();
614 // impact parameter range
615         gener->SetImpactParameterRange(0., 5.);
616         // trigger
617         gener->SetTrigger(1);
618         gener->SetPtJet(25.);
619         gener->SetRadiation(isw);
620         gener->SetSimpleJets(!isw);
621         gener->SetJetEtaRange(-0.3,0.3);
622         gener->SetJetPhiRange(75., 165.);   
623         gGener=gener;
624       }
625         break;
626
627     case kHijing_jj50:
628       {
629         comment = comment.Append("HIJING Jet 50 GeV");
630         AliGenHijing *gener = HijingStandard();
631 // impact parameter range
632         gener->SetImpactParameterRange(0., 5.);
633         // trigger
634         gener->SetTrigger(1);
635         gener->SetPtJet(50.);
636         gener->SetRadiation(isw);
637         gener->SetSimpleJets(!isw);
638         gener->SetJetEtaRange(-0.3,0.3);
639         gener->SetJetPhiRange(75., 165.);   
640         gGener=gener;
641       }
642         break;
643
644     case kHijing_jj75:
645       {
646         comment = comment.Append("HIJING Jet 75 GeV");
647         AliGenHijing *gener = HijingStandard();
648 // impact parameter range
649         gener->SetImpactParameterRange(0., 5.);
650         // trigger
651         gener->SetTrigger(1);
652         gener->SetPtJet(75.);
653         gener->SetRadiation(isw);
654         gener->SetSimpleJets(!isw);
655         gener->SetJetEtaRange(-0.3,0.3);
656         gener->SetJetPhiRange(75., 165.);   
657         gGener=gener;
658       }
659         break;
660
661     case kHijing_jj100:
662       {
663         comment = comment.Append("HIJING Jet 100 GeV");
664         AliGenHijing *gener = HijingStandard();
665 // impact parameter range
666         gener->SetImpactParameterRange(0., 5.);
667         // trigger
668         gener->SetTrigger(1);
669         gener->SetPtJet(100.);
670         gener->SetRadiation(isw);
671         gener->SetSimpleJets(!isw);
672         gener->SetJetEtaRange(-0.3,0.3);
673         gener->SetJetPhiRange(75., 165.);   
674         gGener=gener;
675       }
676         break;
677
678     case kHijing_jj200:
679       {
680         comment = comment.Append("HIJING Jet 200 GeV");
681         AliGenHijing *gener = HijingStandard();
682 // impact parameter range
683         gener->SetImpactParameterRange(0., 5.);
684         // trigger
685         gener->SetTrigger(1);
686         gener->SetPtJet(200.);
687         gener->SetRadiation(isw);
688         gener->SetSimpleJets(!isw);
689         gener->SetJetEtaRange(-0.3,0.3);
690         gener->SetJetPhiRange(75., 165.);   
691         gGener=gener;
692       }
693         break;
694 //
695 // Gamma-Jet
696 //
697     case kHijing_gj25:
698       {
699         comment = comment.Append("HIJING Gamma 25 GeV");
700         AliGenHijing *gener = HijingStandard();
701 // impact parameter range
702         gener->SetImpactParameterRange(0., 5.);
703         // trigger
704         gener->SetTrigger(2);
705         gener->SetPtJet(25.);
706         gener->SetRadiation(isw);
707         gener->SetSimpleJets(!isw);
708         gener->SetJetEtaRange(-0.12, 0.12);
709         gener->SetJetPhiRange(220., 320.);
710         gGener=gener;
711       }
712         break;
713
714     case kHijing_gj50:
715       {
716         comment = comment.Append("HIJING Gamma 50 GeV");
717         AliGenHijing *gener = HijingStandard();
718 // impact parameter range
719         gener->SetImpactParameterRange(0., 5.);
720         // trigger
721         gener->SetTrigger(2);
722         gener->SetPtJet(50.);
723         gener->SetRadiation(isw);
724         gener->SetSimpleJets(!isw);
725         gener->SetJetEtaRange(-0.12, 0.12);
726         gener->SetJetPhiRange(220., 320.);
727         gGener=gener;
728       }
729         break;
730
731     case kHijing_gj75:
732       {
733         comment = comment.Append("HIJING Gamma 75 GeV");
734         AliGenHijing *gener = HijingStandard();
735 // impact parameter range
736         gener->SetImpactParameterRange(0., 5.);
737         // trigger
738         gener->SetTrigger(2);
739         gener->SetPtJet(75.);
740         gener->SetRadiation(isw);
741         gener->SetSimpleJets(!isw);
742         gener->SetJetEtaRange(-0.12, 0.12);
743         gener->SetJetPhiRange(220., 320.);
744         gGener=gener;
745       }
746         break;
747
748     case kHijing_gj100:
749       {
750         comment = comment.Append("HIJING Gamma 100 GeV");
751         AliGenHijing *gener = HijingStandard();
752 // impact parameter range
753         gener->SetImpactParameterRange(0., 5.);
754         // trigger
755         gener->SetTrigger(2);
756         gener->SetPtJet(100.);
757         gener->SetRadiation(isw);
758         gener->SetSimpleJets(!isw);
759         gener->SetJetEtaRange(-0.12, 0.12);
760         gener->SetJetPhiRange(220., 320.);
761         gGener=gener;
762       }
763         break;
764
765     case kHijing_gj200:
766       {
767         comment = comment.Append("HIJING Gamma 200 GeV");
768         AliGenHijing *gener = HijingStandard();
769 // impact parameter range
770         gener->SetImpactParameterRange(0., 5.);
771         // trigger
772         gener->SetTrigger(2);
773         gener->SetPtJet(200.);
774         gener->SetRadiation(isw);
775         gener->SetSimpleJets(!isw);
776         gener->SetJetEtaRange(-0.12, 0.12);
777         gener->SetJetPhiRange(220., 320.);
778         gGener=gener;
779       }
780         break;
781     case kHijing_pA:
782       {
783         comment = comment.Append("HIJING pA");
784
785         AliGenCocktail *gener  = new AliGenCocktail();
786 //      gener->SetTrackingFlag(0);
787         gener->SetOrigin(0, 0, 0);    // vertex position
788         gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
789         gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
790         gener->SetVertexSmear(kPerEvent); 
791         AliGenHijing   *hijing = new AliGenHijing(-1);
792 // centre of mass energy 
793         hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
794 // impact parameter range
795         hijing->SetImpactParameterRange(0., 15.);
796 // reference frame
797         hijing->SetReferenceFrame("CMS");
798         hijing->SetBoostLHC(1);
799 // projectile
800         hijing->SetProjectile("P", 1, 1);
801         hijing->SetTarget    ("A", 208, 82);
802 // tell hijing to keep the full parent child chain
803         hijing->KeepFullEvent();
804 // enable jet quenching
805         hijing->SetJetQuenching(0);
806 // enable shadowing
807         hijing->SetShadowing(1);
808 // Don't track spectators
809         hijing->SetSpectators(0);
810 // kinematic selection
811         hijing->SetSelectAll(0);
812 //
813         AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
814         AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
815         gray->SetSlowNucleonModel(model);
816         gray->SetDebug(1);
817         gener->AddGenerator(hijing,"Hijing pPb", 1);
818         gener->AddGenerator(gray,  "Gray Particles",1);
819         gGener=gener;
820       }
821         break;
822     case kPythia6:
823       {
824         comment = comment.Append(":Pythia p-p @ 14 TeV");
825         AliGenPythia *gener = new AliGenPythia(-1); 
826         gener->SetMomentumRange(0,999999);
827         gener->SetPhiRange(-180,180);
828         gener->SetThetaRange(0., 180.);
829         gener->SetYRange(-12,12);
830         gener->SetPtRange(0,1000);
831         gener->SetStrucFunc(kCTEQ4L);   
832         gener->SetProcess(kPyMb);
833         gener->SetEnergyCMS(14000.);
834         gGener=gener;
835       }
836     case kPythia6Jets:
837       {
838         gener = new AliGenPythia(-1);
839         gener->SetEnergyCMS(5500.);//        Centre of mass energy
840         gener->SetProcess(kPyJets);//        Process type
841         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
842         gener->SetJetPhiRange(0., 360.);
843         gener->SetJetEtRange(10., 1000.);
844         gener->SetGluonRadiation(1,1);
845         //    gener->SetPtKick(0.);
846         //   Structure function
847         gener->SetStrucFunc(kCTEQ4L);
848         gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
849         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
850         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
851         gGener=gener;
852       }
853     break;
854     }
855     return gGener;
856 }
857
858 AliGenHijing* HijingStandard()
859 {
860     AliGenHijing *gener = new AliGenHijing(-1);
861 // centre of mass energy 
862     gener->SetEnergyCMS(5500.);
863 // reference frame
864     gener->SetReferenceFrame("CMS");
865 // projectile
866      gener->SetProjectile("A", 208, 82);
867      gener->SetTarget    ("A", 208, 82);
868 // tell hijing to keep the full parent child chain
869      gener->KeepFullEvent();
870 // enable jet quenching
871      gener->SetJetQuenching(4);
872 // enable shadowing
873      gener->SetShadowing(1);
874 // neutral pion and heavy particle decays switched off
875      gener->SetDecaysOff(1);
876 // Don't track spectators
877      gener->SetSpectators(0);
878 // kinematic selection
879      gener->SetSelectAll(0);
880      return gener;
881 }
882
883
884