]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/ConfigPPR.C
pT saturation set to 2GeV/c in GeVSim
[u/mrichter/AliRoot.git] / macros / ConfigPPR.C
1 // One can use the configuration macro in compiled mode by
2 // root [0] gSystem->Load("libgeant321");
3 // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
4 //                   -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
5 // root [0] .x grun.C(1,"ConfigPPR.C++")
6
7 #if !defined(__CINT__) || defined(__MAKECINT__)
8 #include <Riostream.h>
9 #include <TRandom.h>
10 #include <TSystem.h>
11 #include <TVirtualMC.h>
12 #include <TGeant3TGeo.h>
13 #include <TPDGCode.h>
14 #include <TF1.h>
15 #include "STEER/AliRunLoader.h"
16 #include "STEER/AliRun.h"
17 #include "STEER/AliConfig.h"
18 #include "STEER/AliGenerator.h"
19 #include "STEER/AliLog.h"
20 #include "PYTHIA6/AliDecayerPythia.h"
21 #include "EVGEN/AliGenHIJINGpara.h"
22 #include "THijing/AliGenHijing.h"
23 #include "EVGEN/AliGenCocktail.h"
24 #include "EVGEN/AliGenSlowNucleons.h"
25 #include "EVGEN/AliSlowNucleonModelExp.h"
26 #include "EVGEN/AliGenParam.h"
27 #include "EVGEN/AliGenMUONlib.h"
28 #include "EVGEN/AliGenSTRANGElib.h"
29 #include "EVGEN/AliGenMUONCocktail.h"
30 #include "EVGEN/AliGenCocktail.h"
31 #include "EVGEN/AliGenGeVSim.h"
32 #include "EVGEN/AliGeVSimParticle.h"
33 #include "PYTHIA6/AliGenPythia.h"
34 #include "STEER/AliMagFMaps.h"
35 #include "STRUCT/AliBODY.h"
36 #include "STRUCT/AliMAG.h"
37 #include "STRUCT/AliABSOv0.h"
38 #include "STRUCT/AliDIPOv2.h"
39 #include "STRUCT/AliHALL.h"
40 #include "STRUCT/AliFRAMEv2.h"
41 #include "STRUCT/AliSHILv2.h"
42 #include "STRUCT/AliPIPEv0.h"
43 #include "ITS/AliITSvPPRasymmFMD.h"
44 #include "TPC/AliTPCv2.h"
45 #include "TOF/AliTOFv4T0.h"
46 #include "RICH/AliRICHv1.h"
47 #include "ZDC/AliZDCv2.h"
48 #include "TRD/AliTRDv1.h"
49 #include "FMD/AliFMDv1.h"
50 #include "MUON/AliMUONv1.h"
51 #include "MUON/AliMUONSt1GeometryBuilderV2.h"
52 #include "MUON/AliMUONSt2GeometryBuilderV2.h"
53 #include "MUON/AliMUONSlatGeometryBuilder.h"
54 #include "MUON/AliMUONTriggerGeometryBuilder.h"
55 #include "PHOS/AliPHOSv1.h"
56 #include "PMD/AliPMDv1.h"
57 #include "START/AliSTARTv1.h"
58 #include "EMCAL/AliEMCALv1.h"
59 #include "CRT/AliCRTv0.h"
60 #include "VZERO/AliVZEROv5.h"
61 #endif
62
63 enum PprRun_t 
64 {
65     test50,
66     kParam_8000,   kParam_4000,  kParam_2000, 
67     kHijing_cent1, kHijing_cent2, 
68     kHijing_per1,  kHijing_per2, kHijing_per3, kHijing_per4,  kHijing_per5,
69     kHijing_jj25,  kHijing_jj50, kHijing_jj75, kHijing_jj100, kHijing_jj200, 
70     kHijing_gj25,  kHijing_gj50, kHijing_gj75, kHijing_gj100, kHijing_gj200,
71     kHijing_pA, kPythia6, 
72     kPythia6Jets20_24,   kPythia6Jets24_29,   kPythia6Jets29_35,
73     kPythia6Jets35_42,   kPythia6Jets42_50,   kPythia6Jets50_60,
74     kPythia6Jets60_72,   kPythia6Jets72_86,   kPythia6Jets86_104,
75     kPythia6Jets104_125, kPythia6Jets125_150, kPythia6Jets150_180,
76     kD0PbPb5500, kCharmSemiElPbPb5500, kBeautySemiElPbPb5500,
77     kCocktailTRD, kPyJJ, kPyGJ, 
78     kMuonCocktailCent1, kMuonCocktailPer1, kMuonCocktailPer4, 
79     kMuonCocktailCent1HighPt, kMuonCocktailPer1HighPt, kMuonCocktailPer4HighPt,
80     kMuonCocktailCent1Single, kMuonCocktailPer1Single, kMuonCocktailPer4Single,
81     kFlow_2_2000, kFlow_10_2000, kFlow_6_2000, kFlow_6_5000,
82     kHIJINGplus, kRunMax
83 };
84
85 const char* pprRunName[] = {
86     "test50",
87     "kParam_8000",   "kParam_4000",  "kParam_2000", 
88     "kHijing_cent1", "kHijing_cent2", 
89     "kHijing_per1",  "kHijing_per2", "kHijing_per3", "kHijing_per4",  
90     "kHijing_per5",
91     "kHijing_jj25",  "kHijing_jj50", "kHijing_jj75", "kHijing_jj100", 
92     "kHijing_jj200", 
93     "kHijing_gj25",  "kHijing_gj50", "kHijing_gj75", "kHijing_gj100", 
94     "kHijing_gj200", "kHijing_pA", "kPythia6", 
95     "kPythia6Jets20_24",   "kPythia6Jets24_29",   "kPythia6Jets29_35",
96     "kPythia6Jets35_42",   "kPythia6Jets42_50",   "kPythia6Jets50_60",
97     "kPythia6Jets60_72",   "kPythia6Jets72_86",   "kPythia6Jets86_104",
98     "kPythia6Jets104_125", "kPythia6Jets125_150", "kPythia6Jets150_180",
99     "kD0PbPb5500", "kCharmSemiElPbPb5500", "kBeautySemiElPbPb5500",
100     "kCocktailTRD", "kPyJJ", "kPyGJ", 
101     "kMuonCocktailCent1", "kMuonCocktailPer1", "kMuonCocktailPer4",  
102     "kMuonCocktailCent1HighPt", "kMuonCocktailPer1HighPt", "kMuonCocktailPer4HighPt",
103     "kMuonCocktailCent1Single", "kMuonCocktailPer1Single", "kMuonCocktailPer4Single",
104     "kFlow_2_2000", "kFlow_10_2000", "kFlow_6_2000", "kFlow_6_5000", "kHIJINGplus"
105 };
106
107 enum PprGeo_t 
108 {
109     kHoles, kNoHoles
110 };
111
112 enum PprRad_t
113 {
114     kGluonRadiation, kNoGluonRadiation
115 };
116
117 enum PprMag_t
118 {
119     k2kG, k4kG, k5kG
120 };
121
122
123 // This part for configuration    
124 //static PprRun_t srun = test50;
125 static PprRun_t srun = kHIJINGplus;
126 static PprGeo_t sgeo = kNoHoles;
127 static PprRad_t srad = kGluonRadiation;
128 static PprMag_t smag = k5kG;
129 static Int_t    sseed = 12345; //Set 0 to use the current time
130
131 // Comment line 
132 static TString  comment;
133
134 // Functions
135 Float_t EtaToTheta(Float_t arg);
136 AliGenerator* GeneratorFactory(PprRun_t srun);
137 AliGenHijing* HijingStandard();
138 AliGenGeVSim* GeVSimStandard(Float_t, Float_t);
139 void ProcessEnvironmentVars();
140
141 void Config()
142 {
143     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
144     // Theta range given through pseudorapidity limits 22/6/2001
145
146     // Get settings from environment variables
147     ProcessEnvironmentVars();
148
149     // Set Random Number seed
150     gRandom->SetSeed(sseed);
151     cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
152
153
154    // libraries required by geant321
155 #if defined(__CINT__)
156     gSystem->Load("libgeant321");
157 #endif
158
159     new     TGeant3TGeo("C++ Interface to Geant3");
160
161     AliRunLoader* rl=0x0;
162
163     cout<<"Config.C: Creating Run Loader ..."<<endl;
164     rl = AliRunLoader::Open("galice.root",
165                             AliConfig::GetDefaultEventFolderName(),
166                             "recreate");
167     if (rl == 0x0)
168       {
169         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
170         return;
171       }
172     rl->SetCompressionLevel(2);
173     rl->SetNumberOfEventsPerFile(3);
174     gAlice->SetRunLoader(rl);
175
176     //
177     // Set External decayer
178     AliDecayer *decayer = new AliDecayerPythia();
179
180
181     switch (srun) {
182     case kD0PbPb5500:
183       decayer->SetForceDecay(kHadronicD);
184       break;
185     case kCharmSemiElPbPb5500:
186       decayer->SetForceDecay(kSemiElectronic);
187       break;
188     case kBeautySemiElPbPb5500:
189       decayer->SetForceDecay(kSemiElectronic);
190       break;
191     default:
192       decayer->SetForceDecay(kAll);
193       break;
194     }
195     decayer->Init();
196     gMC->SetExternalDecayer(decayer);
197     //
198     //
199     //=======================================================================
200     //
201     //=======================================================================
202     // ************* STEERING parameters FOR ALICE SIMULATION **************
203     // --- Specify event type to be tracked through the ALICE setup
204     // --- All positions are in cm, angles in degrees, and P and E in GeV
205
206     gMC->SetProcess("DCAY",1);
207     gMC->SetProcess("PAIR",1);
208     gMC->SetProcess("COMP",1);
209     gMC->SetProcess("PHOT",1);
210     gMC->SetProcess("PFIS",0);
211     gMC->SetProcess("DRAY",0);
212     gMC->SetProcess("ANNI",1);
213     gMC->SetProcess("BREM",1);
214     gMC->SetProcess("MUNU",1);
215     gMC->SetProcess("CKOV",1);
216     gMC->SetProcess("HADR",1);
217     gMC->SetProcess("LOSS",2);
218     gMC->SetProcess("MULS",1);
219     gMC->SetProcess("RAYL",1);
220
221     Float_t cut = 1.e-3;        // 1MeV cut by default
222     Float_t tofmax = 1.e10;
223
224     gMC->SetCut("CUTGAM", cut);
225     gMC->SetCut("CUTELE", cut);
226     gMC->SetCut("CUTNEU", cut);
227     gMC->SetCut("CUTHAD", cut);
228     gMC->SetCut("CUTMUO", cut);
229     gMC->SetCut("BCUTE",  cut); 
230     gMC->SetCut("BCUTM",  cut); 
231     gMC->SetCut("DCUTE",  cut); 
232     gMC->SetCut("DCUTM",  cut); 
233     gMC->SetCut("PPCUTM", cut);
234     gMC->SetCut("TOFMAX", tofmax); 
235
236     // Debug and log level
237     AliLog::SetGlobalDebugLevel(0);
238     AliLog::SetGlobalLogLevel(AliLog::kError);
239
240     // Generator Configuration
241     AliGenerator* gener = GeneratorFactory(srun);
242     gener->SetOrigin(0, 0, 0);    // vertex position
243     gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
244     gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
245     gener->SetVertexSmear(kPerEvent); 
246     gener->SetTrackingFlag(1);
247     gener->Init();
248     
249     if (smag == k2kG) {
250         comment = comment.Append(" | L3 field 0.2 T");
251     } else if (smag == k4kG) {
252         comment = comment.Append(" | L3 field 0.4 T");
253     } else if (smag == k5kG) {
254         comment = comment.Append(" | L3 field 0.5 T");
255     }
256     
257     
258     if (srad == kGluonRadiation)
259     {
260         comment = comment.Append(" | Gluon Radiation On");
261         
262     } else {
263         comment = comment.Append(" | Gluon Radiation Off");
264     }
265
266     if (sgeo == kHoles)
267     {
268         comment = comment.Append(" | Holes for PHOS/RICH");
269         
270     } else {
271         comment = comment.Append(" | No holes for PHOS/RICH");
272     }
273
274     printf("\n \n Comment: %s \n \n", comment.Data());
275     
276     
277 // Field (L3 0.4 T)
278     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., smag);
279     field->SetL3ConstField(0); //Using const. field in the barrel
280     rl->CdGAFile();
281     gAlice->SetField(field);    
282 //
283     Int_t   iABSO   = 1;
284     Int_t   iDIPO   = 1;
285     Int_t   iFMD    = 1;
286     Int_t   iFRAME  = 1;
287     Int_t   iHALL   = 1;
288     Int_t   iITS    = 1;
289     Int_t   iMAG    = 1;
290     Int_t   iMUON   = 1;
291     Int_t   iPHOS   = 1;
292     Int_t   iPIPE   = 1;
293     Int_t   iPMD    = 1;
294     Int_t   iRICH   = 1;
295     Int_t   iSHIL   = 1;
296     Int_t   iSTART  = 1;
297     Int_t   iTOF    = 1;
298     Int_t   iTPC    = 1;
299     Int_t   iTRD    = 1;
300     Int_t   iZDC    = 1;
301     Int_t   iEMCAL  = 1;
302     Int_t   iVZERO  = 1;
303     Int_t   iCRT    = 0;
304
305     //=================== Alice BODY parameters =============================
306     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
307
308
309     if (iMAG)
310     {
311         //=================== MAG parameters ============================
312         // --- Start with Magnet since detector layouts may be depending ---
313         // --- on the selected Magnet dimensions ---
314         AliMAG *MAG = new AliMAG("MAG", "Magnet");
315     }
316
317
318     if (iABSO)
319     {
320         //=================== ABSO parameters ============================
321         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
322     }
323
324     if (iDIPO)
325     {
326         //=================== DIPO parameters ============================
327
328         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
329     }
330
331     if (iHALL)
332     {
333         //=================== HALL parameters ============================
334
335         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
336     }
337
338
339     if (iFRAME)
340     {
341         //=================== FRAME parameters ============================
342
343         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
344         if (sgeo == kHoles) {
345             FRAME->SetHoles(1);
346         } else {
347             FRAME->SetHoles(0);
348         }
349     }
350
351     if (iSHIL)
352     {
353         //=================== SHIL parameters ============================
354
355         AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
356     }
357
358
359     if (iPIPE)
360     {
361         //=================== PIPE parameters ============================
362
363         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
364     }
365  
366     if(iITS) {
367
368     //=================== ITS parameters ============================
369     //
370     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
371     // almost all other detectors. This involves the fact that the ITS geometry
372     // still has several options to be followed in parallel in order to determine
373     // the best set-up which minimizes the induced background. All the geometries
374     // available to date are described in the following. Read carefully the comments
375     // and use the default version (the only one uncommented) unless you are making
376     // comparisons and you know what you are doing. In this case just uncomment the
377     // ITS geometry you want to use and run Aliroot.
378     //
379     // Detailed geometries:         
380     //
381     //
382     //AliITS *ITS  = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
383     //
384     //AliITS *ITS  = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
385     //
386         AliITSvPPRasymmFMD *ITS  = new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version with asymmetric services");
387         ITS->SetMinorVersion(2);  // don't touch this parameter if you're not an ITS developer
388         ITS->SetReadDet(kTRUE);   // don't touch this parameter if you're not an ITS developer
389     //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
390         ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
391         ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
392         ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
393         ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
394         ITS->SetRails(0);            // 1 --> rails in ; 0 --> rails out
395         ITS->SetCoolingFluid(1);   // 1 --> water ; 0 --> freon
396
397     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
398     // for reconstruction !):
399     //                                                     
400     //
401     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
402     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
403     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
404     //
405     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
406     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
407     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
408     //                      
409     //
410     //
411     // Geant3 <-> EUCLID conversion
412     // ============================
413     //
414     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
415     // media to two ASCII files (called by default ITSgeometry.euc and
416     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
417     // The default (=0) means that you dont want to use this facility.
418     //
419         ITS->SetEUCLID(0);  
420     }
421
422     if (iTPC)
423     {
424         //============================ TPC parameters ================================
425         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
426         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
427         // --- sectors are specified, any value other than that requires at least one 
428         // --- sector (lower or upper)to be specified!
429         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
430         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
431         // --- SecLows - number of lower sectors specified (up to 6)
432         // --- SecUps - number of upper sectors specified (up to 12)
433         // --- Sens - sensitive strips for the Slow Simulator !!!
434         // --- This does NOT work if all S or L-sectors are specified, i.e.
435         // --- if SecAL or SecAU < 0
436         //
437         //
438         //-----------------------------------------------------------------------------
439
440         //  gROOT->LoadMacro("SetTPCParam.C");
441         //  AliTPCParam *param = SetTPCParam();
442         AliTPC *TPC = new AliTPCv2("TPC", "Default");
443
444         // All sectors included 
445         TPC->SetSecAL(-1);
446         TPC->SetSecAU(-1);
447
448     }
449
450
451     if (iTOF) {
452         //=================== TOF parameters ============================
453         AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
454     }
455
456
457     if (iRICH)
458     {
459         //=================== RICH parameters ===========================
460         AliRICH *RICH = new AliRICHv1("RICH", "normal RICH");
461
462     }
463
464
465     if (iZDC)
466     {
467         //=================== ZDC parameters ============================
468
469         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
470     }
471
472     if (iTRD)
473     {
474         //=================== TRD parameters ============================
475
476         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
477
478         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
479         TRD->SetGasMix(1);
480         if (sgeo == kHoles) {
481             // With hole in front of PHOS
482             TRD->SetPHOShole();
483             // With hole in front of RICH
484             TRD->SetRICHhole();
485         }
486             // Switch on TR
487             AliTRDsim *TRDsim = TRD->CreateTR();
488     }
489
490     if (iFMD)
491     {
492         //=================== FMD parameters ============================
493         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
494    }
495
496     if (iMUON)
497     {
498         //=================== MUON parameters ===========================
499         // New MUONv1 version (geometry defined via builders)
500         AliMUON *MUON = new AliMUONv1("MUON", "default");
501         ((AliMUONv1*)MUON)->SetStepManagerVersionDE(true);
502         MUON->AddGeometryBuilder(new AliMUONSt1GeometryBuilderV2(MUON));
503         MUON->AddGeometryBuilder(new AliMUONSt2GeometryBuilderV2(MUON));
504         MUON->AddGeometryBuilder(new AliMUONSlatGeometryBuilder(MUON));
505         MUON->AddGeometryBuilder(new AliMUONTriggerGeometryBuilder(MUON));
506     }
507     //=================== PHOS parameters ===========================
508
509     if (iPHOS)
510     {
511         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
512     }
513
514
515     if (iPMD)
516     {
517         //=================== PMD parameters ============================
518         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
519     }
520
521     if (iSTART)
522     {
523         //=================== START parameters ============================
524         AliSTART *START = new AliSTARTv1("START", "START Detector");
525     }
526
527     if (iEMCAL)
528     {
529         //=================== EMCAL parameters ============================
530         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCAL_55_25");
531     }
532
533      if (iCRT)
534     {
535         //=================== CRT parameters ============================
536         AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
537     }
538
539      if (iVZERO)
540     {
541         //=================== CRT parameters ============================
542         AliVZERO *VZERO = new AliVZEROv5("VZERO", "normal VZERO");
543     }
544  
545              
546 }
547
548 Float_t EtaToTheta(Float_t arg){
549   return (180./TMath::Pi())*2.*atan(exp(-arg));
550 }
551
552
553
554 AliGenerator* GeneratorFactory(PprRun_t srun) {
555     Int_t isw = 3;
556     if (srad == kNoGluonRadiation) isw = 0;
557     
558
559     AliGenerator * gGener = 0x0;
560     switch (srun) {
561     case test50:
562       {
563         comment = comment.Append(":HIJINGparam test 50 particles");
564         AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
565         gener->SetMomentumRange(0, 999999.);
566         gener->SetPhiRange(0., 360.);
567         // Set pseudorapidity range from -8 to 8.
568         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
569         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
570         gener->SetThetaRange(thmin,thmax);
571         gGener=gener;
572       }
573       break;
574     case kParam_8000:
575       {
576         comment = comment.Append(":HIJINGparam N=8000");
577         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
578         gener->SetMomentumRange(0, 999999.);
579         gener->SetPhiRange(0., 360.);
580         // Set pseudorapidity range from -8 to 8.
581         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
582         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
583         gener->SetThetaRange(thmin,thmax);
584         gGener=gener;
585       }
586       break;
587     case kParam_4000:
588       {
589         comment = comment.Append("HIJINGparam N=4000");
590         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
591         gener->SetMomentumRange(0, 999999.);
592         gener->SetPhiRange(0., 360.);
593         // Set pseudorapidity range from -8 to 8.
594         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
595         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
596         gener->SetThetaRange(thmin,thmax);
597         gGener=gener;
598       }
599         break;
600     case kParam_2000:
601       {
602         comment = comment.Append("HIJINGparam N=2000");
603         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
604         gener->SetMomentumRange(0, 999999.);
605         gener->SetPhiRange(0., 360.);
606         // Set pseudorapidity range from -8 to 8.
607         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
608         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
609         gener->SetThetaRange(thmin,thmax);
610         gGener=gener;
611       }
612       break;
613 //
614 //  Hijing Central
615 //
616     case kHijing_cent1:
617       {
618         comment = comment.Append("HIJING cent1");
619         AliGenHijing *gener = HijingStandard();
620 // impact parameter range
621         gener->SetImpactParameterRange(0., 5.);
622         gGener=gener;
623       }
624       break;
625     case kHijing_cent2:
626       {
627         comment = comment.Append("HIJING cent2");
628         AliGenHijing *gener = HijingStandard();
629 // impact parameter range
630         gener->SetImpactParameterRange(0., 2.);
631         gGener=gener;
632       }
633       break;
634 //
635 // Hijing Peripheral 
636 //
637     case kHijing_per1:
638       {
639         comment = comment.Append("HIJING per1");
640         AliGenHijing *gener = HijingStandard();
641 // impact parameter range
642         gener->SetImpactParameterRange(5., 8.6);
643         gGener=gener;
644       }
645       break;
646     case kHijing_per2:
647       {
648         comment = comment.Append("HIJING per2");
649         AliGenHijing *gener = HijingStandard();
650 // impact parameter range
651         gener->SetImpactParameterRange(8.6, 11.2);
652         gGener=gener;
653       }
654       break;
655     case kHijing_per3:
656       {
657         comment = comment.Append("HIJING per3");
658         AliGenHijing *gener = HijingStandard();
659 // impact parameter range
660         gener->SetImpactParameterRange(11.2, 13.2);
661         gGener=gener;
662       }
663       break;
664     case kHijing_per4:
665       {
666         comment = comment.Append("HIJING per4");
667         AliGenHijing *gener = HijingStandard();
668 // impact parameter range
669         gener->SetImpactParameterRange(13.2, 15.);
670         gGener=gener;
671       }
672       break;
673     case kHijing_per5:
674       {
675         comment = comment.Append("HIJING per5");
676         AliGenHijing *gener = HijingStandard();
677 // impact parameter range
678         gener->SetImpactParameterRange(15., 100.);
679         gGener=gener;
680       }
681       break;
682 //
683 //  Jet-Jet
684 //
685     case kHijing_jj25:
686       {
687         comment = comment.Append("HIJING Jet 25 GeV");
688         AliGenHijing *gener = HijingStandard();
689 // impact parameter range
690         gener->SetImpactParameterRange(0., 5.);
691         // trigger
692         gener->SetTrigger(1);
693         gener->SetPtJet(25.);
694         gener->SetRadiation(isw);
695         gener->SetSimpleJets(!isw);
696         gener->SetJetEtaRange(-0.3,0.3);
697         gener->SetJetPhiRange(75., 165.);   
698         gGener=gener;
699       }
700       break;
701
702     case kHijing_jj50:
703       {
704         comment = comment.Append("HIJING Jet 50 GeV");
705         AliGenHijing *gener = HijingStandard();
706 // impact parameter range
707         gener->SetImpactParameterRange(0., 5.);
708         // trigger
709         gener->SetTrigger(1);
710         gener->SetPtJet(50.);
711         gener->SetRadiation(isw);
712         gener->SetSimpleJets(!isw);
713         gener->SetJetEtaRange(-0.3,0.3);
714         gener->SetJetPhiRange(75., 165.);   
715         gGener=gener;
716       }
717         break;
718
719     case kHijing_jj75:
720       {
721         comment = comment.Append("HIJING Jet 75 GeV");
722         AliGenHijing *gener = HijingStandard();
723 // impact parameter range
724         gener->SetImpactParameterRange(0., 5.);
725         // trigger
726         gener->SetTrigger(1);
727         gener->SetPtJet(75.);
728         gener->SetRadiation(isw);
729         gener->SetSimpleJets(!isw);
730         gener->SetJetEtaRange(-0.3,0.3);
731         gener->SetJetPhiRange(75., 165.);   
732         gGener=gener;
733       }
734       break;
735
736     case kHijing_jj100:
737       {
738         comment = comment.Append("HIJING Jet 100 GeV");
739         AliGenHijing *gener = HijingStandard();
740 // impact parameter range
741         gener->SetImpactParameterRange(0., 5.);
742         // trigger
743         gener->SetTrigger(1);
744         gener->SetPtJet(100.);
745         gener->SetRadiation(isw);
746         gener->SetSimpleJets(!isw);
747         gener->SetJetEtaRange(-0.3,0.3);
748         gener->SetJetPhiRange(75., 165.);   
749         gGener=gener;
750       }
751       break;
752
753     case kHijing_jj200:
754       {
755         comment = comment.Append("HIJING Jet 200 GeV");
756         AliGenHijing *gener = HijingStandard();
757 // impact parameter range
758         gener->SetImpactParameterRange(0., 5.);
759         // trigger
760         gener->SetTrigger(1);
761         gener->SetPtJet(200.);
762         gener->SetRadiation(isw);
763         gener->SetSimpleJets(!isw);
764         gener->SetJetEtaRange(-0.3,0.3);
765         gener->SetJetPhiRange(75., 165.);   
766         gGener=gener;
767       }
768       break;
769 //
770 // Gamma-Jet
771 //
772     case kHijing_gj25:
773       {
774         comment = comment.Append("HIJING Gamma 25 GeV");
775         AliGenHijing *gener = HijingStandard();
776 // impact parameter range
777         gener->SetImpactParameterRange(0., 5.);
778         // trigger
779         gener->SetTrigger(2);
780         gener->SetPtJet(25.);
781         gener->SetRadiation(isw);
782         gener->SetSimpleJets(!isw);
783         gener->SetJetEtaRange(-0.12, 0.12);
784         gener->SetJetPhiRange(220., 320.);
785         gGener=gener;
786       }
787       break;
788
789     case kHijing_gj50:
790       {
791         comment = comment.Append("HIJING Gamma 50 GeV");
792         AliGenHijing *gener = HijingStandard();
793 // impact parameter range
794         gener->SetImpactParameterRange(0., 5.);
795         // trigger
796         gener->SetTrigger(2);
797         gener->SetPtJet(50.);
798         gener->SetRadiation(isw);
799         gener->SetSimpleJets(!isw);
800         gener->SetJetEtaRange(-0.12, 0.12);
801         gener->SetJetPhiRange(220., 320.);
802         gGener=gener;
803       }
804       break;
805
806     case kHijing_gj75:
807       {
808         comment = comment.Append("HIJING Gamma 75 GeV");
809         AliGenHijing *gener = HijingStandard();
810 // impact parameter range
811         gener->SetImpactParameterRange(0., 5.);
812         // trigger
813         gener->SetTrigger(2);
814         gener->SetPtJet(75.);
815         gener->SetRadiation(isw);
816         gener->SetSimpleJets(!isw);
817         gener->SetJetEtaRange(-0.12, 0.12);
818         gener->SetJetPhiRange(220., 320.);
819         gGener=gener;
820       }
821       break;
822
823     case kHijing_gj100:
824       {
825         comment = comment.Append("HIJING Gamma 100 GeV");
826         AliGenHijing *gener = HijingStandard();
827 // impact parameter range
828         gener->SetImpactParameterRange(0., 5.);
829         // trigger
830         gener->SetTrigger(2);
831         gener->SetPtJet(100.);
832         gener->SetRadiation(isw);
833         gener->SetSimpleJets(!isw);
834         gener->SetJetEtaRange(-0.12, 0.12);
835         gener->SetJetPhiRange(220., 320.);
836         gGener=gener;
837       }
838       break;
839
840     case kHijing_gj200:
841       {
842         comment = comment.Append("HIJING Gamma 200 GeV");
843         AliGenHijing *gener = HijingStandard();
844 // impact parameter range
845         gener->SetImpactParameterRange(0., 5.);
846         // trigger
847         gener->SetTrigger(2);
848         gener->SetPtJet(200.);
849         gener->SetRadiation(isw);
850         gener->SetSimpleJets(!isw);
851         gener->SetJetEtaRange(-0.12, 0.12);
852         gener->SetJetPhiRange(220., 320.);
853         gGener=gener;
854       }
855       break;
856     case kHijing_pA:
857       {
858         comment = comment.Append("HIJING pA");
859
860         AliGenCocktail *gener  = new AliGenCocktail();
861
862         AliGenHijing   *hijing = new AliGenHijing(-1);
863 // centre of mass energy 
864         hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
865 // impact parameter range
866         hijing->SetImpactParameterRange(0., 15.);
867 // reference frame
868         hijing->SetReferenceFrame("CMS");
869         hijing->SetBoostLHC(1);
870 // projectile
871         hijing->SetProjectile("P", 1, 1);
872         hijing->SetTarget    ("A", 208, 82);
873 // tell hijing to keep the full parent child chain
874         hijing->KeepFullEvent();
875 // enable jet quenching
876         hijing->SetJetQuenching(0);
877 // enable shadowing
878         hijing->SetShadowing(1);
879 // Don't track spectators
880         hijing->SetSpectators(0);
881 // kinematic selection
882         hijing->SetSelectAll(0);
883 //
884         AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
885         AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
886         gray->SetSlowNucleonModel(model);
887         gray->SetDebug(1);
888         gener->AddGenerator(hijing,"Hijing pPb", 1);
889         gener->AddGenerator(gray,  "Gray Particles",1);
890         gGener=gener;
891       }
892       break;
893     case kPythia6:
894       {
895         comment = comment.Append(":Pythia p-p @ 14 TeV");
896         AliGenPythia *gener = new AliGenPythia(-1); 
897         gener->SetMomentumRange(0,999999);
898         gener->SetThetaRange(0., 180.);
899         gener->SetYRange(-12,12);
900         gener->SetPtRange(0,1000);
901         gener->SetProcess(kPyMb);
902         gener->SetEnergyCMS(14000.);
903         gGener=gener;
904       }
905       break;
906     case kPythia6Jets20_24:
907       {
908         comment = comment.Append(":Pythia jets 20-24 GeV @ 5.5 TeV");
909         AliGenPythia * gener = new AliGenPythia(-1);
910         gener->SetEnergyCMS(5500.);//        Centre of mass energy
911         gener->SetProcess(kPyJets);//        Process type
912         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
913         gener->SetJetPhiRange(0., 360.);
914         gener->SetJetEtRange(10., 1000.);
915         gener->SetGluonRadiation(1,1);
916         //    gener->SetPtKick(0.);
917         //   Structure function
918         gener->SetStrucFunc(kCTEQ4L);
919         gener->SetPtHard(20., 24.);// Pt transfer of the hard scattering
920         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
921         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
922         gGener=gener;
923       }
924       break;
925     case kPythia6Jets24_29:
926       {
927         comment = comment.Append(":Pythia jets 24-29 GeV @ 5.5 TeV");
928         AliGenPythia * gener = new AliGenPythia(-1);
929         gener->SetEnergyCMS(5500.);//        Centre of mass energy
930         gener->SetProcess(kPyJets);//        Process type
931         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
932         gener->SetJetPhiRange(0., 360.);
933         gener->SetJetEtRange(10., 1000.);
934         gener->SetGluonRadiation(1,1);
935         //    gener->SetPtKick(0.);
936         //   Structure function
937         gener->SetStrucFunc(kCTEQ4L);
938         gener->SetPtHard(24., 29.);// Pt transfer of the hard scattering
939         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
940         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
941         gGener=gener;
942       }
943       break;
944     case kPythia6Jets29_35:
945       {
946         comment = comment.Append(":Pythia jets 29-35 GeV @ 5.5 TeV");
947         AliGenPythia * gener = new AliGenPythia(-1);
948         gener->SetEnergyCMS(5500.);//        Centre of mass energy
949         gener->SetProcess(kPyJets);//        Process type
950         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
951         gener->SetJetPhiRange(0., 360.);
952         gener->SetJetEtRange(10., 1000.);
953         gener->SetGluonRadiation(1,1);
954         //    gener->SetPtKick(0.);
955         //   Structure function
956         gener->SetStrucFunc(kCTEQ4L);
957         gener->SetPtHard(29., 35.);// Pt transfer of the hard scattering
958         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
959         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
960         gGener=gener;
961       }
962       break;
963     case kPythia6Jets35_42:
964       {
965         comment = comment.Append(":Pythia jets 35-42 GeV @ 5.5 TeV");
966         AliGenPythia * gener = new AliGenPythia(-1);
967         gener->SetEnergyCMS(5500.);//        Centre of mass energy
968         gener->SetProcess(kPyJets);//        Process type
969         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
970         gener->SetJetPhiRange(0., 360.);
971         gener->SetJetEtRange(10., 1000.);
972         gener->SetGluonRadiation(1,1);
973         //    gener->SetPtKick(0.);
974         //   Structure function
975         gener->SetStrucFunc(kCTEQ4L);
976         gener->SetPtHard(35., 42.);// Pt transfer of the hard scattering
977         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
978         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
979         gGener=gener;
980       }
981       break;
982     case kPythia6Jets42_50:
983       {
984         comment = comment.Append(":Pythia jets 42-50 GeV @ 5.5 TeV");
985         AliGenPythia * gener = new AliGenPythia(-1);
986         gener->SetEnergyCMS(5500.);//        Centre of mass energy
987         gener->SetProcess(kPyJets);//        Process type
988         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
989         gener->SetJetPhiRange(0., 360.);
990         gener->SetJetEtRange(10., 1000.);
991         gener->SetGluonRadiation(1,1);
992         //    gener->SetPtKick(0.);
993         //   Structure function
994         gener->SetStrucFunc(kCTEQ4L);
995         gener->SetPtHard(42., 50.);// Pt transfer of the hard scattering
996         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
997         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
998         gGener=gener;
999       }
1000       break;
1001     case kPythia6Jets50_60:
1002       {
1003         comment = comment.Append(":Pythia jets 50-60 GeV @ 5.5 TeV");
1004         AliGenPythia * gener = new AliGenPythia(-1);
1005         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1006         gener->SetProcess(kPyJets);//        Process type
1007         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1008         gener->SetJetPhiRange(0., 360.);
1009         gener->SetJetEtRange(10., 1000.);
1010         gener->SetGluonRadiation(1,1);
1011         //    gener->SetPtKick(0.);
1012         //   Structure function
1013         gener->SetStrucFunc(kCTEQ4L);
1014         gener->SetPtHard(50., 60.);// Pt transfer of the hard scattering
1015         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1016         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1017         gGener=gener;
1018       }
1019       break;
1020     case kPythia6Jets60_72:
1021       {
1022         comment = comment.Append(":Pythia jets 60-72 GeV @ 5.5 TeV");
1023         AliGenPythia * gener = new AliGenPythia(-1);
1024         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1025         gener->SetProcess(kPyJets);//        Process type
1026         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1027         gener->SetJetPhiRange(0., 360.);
1028         gener->SetJetEtRange(10., 1000.);
1029         gener->SetGluonRadiation(1,1);
1030         //    gener->SetPtKick(0.);
1031         //   Structure function
1032         gener->SetStrucFunc(kCTEQ4L);
1033         gener->SetPtHard(60., 72.);// Pt transfer of the hard scattering
1034         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1035         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1036         gGener=gener;
1037       }
1038       break;
1039     case kPythia6Jets72_86:
1040       {
1041         comment = comment.Append(":Pythia jets 72-86 GeV @ 5.5 TeV");
1042         AliGenPythia * gener = new AliGenPythia(-1);
1043         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1044         gener->SetProcess(kPyJets);//        Process type
1045         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1046         gener->SetJetPhiRange(0., 360.);
1047         gener->SetJetEtRange(10., 1000.);
1048         gener->SetGluonRadiation(1,1);
1049         //    gener->SetPtKick(0.);
1050         //   Structure function
1051         gener->SetStrucFunc(kCTEQ4L);
1052         gener->SetPtHard(72., 86.);// Pt transfer of the hard scattering
1053         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1054         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1055         gGener=gener;
1056       }
1057       break;
1058     case kPythia6Jets86_104:
1059       {
1060         comment = comment.Append(":Pythia jets 86-104 GeV @ 5.5 TeV");
1061         AliGenPythia * gener = new AliGenPythia(-1);
1062         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1063         gener->SetProcess(kPyJets);//        Process type
1064         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1065         gener->SetJetPhiRange(0., 360.);
1066         gener->SetJetEtRange(10., 1000.);
1067         gener->SetGluonRadiation(1,1);
1068         //    gener->SetPtKick(0.);
1069         //   Structure function
1070         gener->SetStrucFunc(kCTEQ4L);
1071         gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
1072         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1073         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1074         gGener=gener;
1075       }
1076       break;
1077     case kPythia6Jets104_125:
1078       {
1079         comment = comment.Append(":Pythia jets 105-125 GeV @ 5.5 TeV");
1080         AliGenPythia * gener = new AliGenPythia(-1);
1081         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1082         gener->SetProcess(kPyJets);//        Process type
1083         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1084         gener->SetJetPhiRange(0., 360.);
1085         gener->SetJetEtRange(10., 1000.);
1086         gener->SetGluonRadiation(1,1);
1087         //    gener->SetPtKick(0.);
1088         //   Structure function
1089         gener->SetStrucFunc(kCTEQ4L);
1090         gener->SetPtHard(104., 125.);// Pt transfer of the hard scattering
1091         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1092         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1093         gGener=gener;
1094       }
1095       break;
1096     case kPythia6Jets125_150:
1097       {
1098         comment = comment.Append(":Pythia jets 125-150 GeV @ 5.5 TeV");
1099         AliGenPythia * gener = new AliGenPythia(-1);
1100         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1101         gener->SetProcess(kPyJets);//        Process type
1102         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1103         gener->SetJetPhiRange(0., 360.);
1104         gener->SetJetEtRange(10., 1000.);
1105         gener->SetGluonRadiation(1,1);
1106         //    gener->SetPtKick(0.);
1107         //   Structure function
1108         gener->SetStrucFunc(kCTEQ4L);
1109         gener->SetPtHard(125., 150.);// Pt transfer of the hard scattering
1110         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1111         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1112         gGener=gener;
1113       }
1114       break;
1115     case kPythia6Jets150_180:
1116       {
1117         comment = comment.Append(":Pythia jets 150-180 GeV @ 5.5 TeV");
1118         AliGenPythia * gener = new AliGenPythia(-1);
1119         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1120         gener->SetProcess(kPyJets);//        Process type
1121         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1122         gener->SetJetPhiRange(0., 360.);
1123         gener->SetJetEtRange(10., 1000.);
1124         gener->SetGluonRadiation(1,1);
1125         //    gener->SetPtKick(0.);
1126         //   Structure function
1127         gener->SetStrucFunc(kCTEQ4L);
1128         gener->SetPtHard(150., 180.);// Pt transfer of the hard scattering
1129         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1130         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1131         gGener=gener;
1132       }
1133       break;
1134     case kD0PbPb5500:
1135       {
1136         comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
1137         AliGenPythia * gener = new AliGenPythia(10);
1138         gener->SetProcess(kPyD0PbPbMNR);
1139         gener->SetStrucFunc(kCTEQ4L);
1140         gener->SetPtHard(2.1,-1.0);
1141         gener->SetEnergyCMS(5500.);
1142         gener->SetNuclei(208,208);
1143         gener->SetForceDecay(kHadronicD);
1144         gener->SetYRange(-2,2);
1145         gener->SetFeedDownHigherFamily(kFALSE);
1146         gener->SetStackFillOpt(AliGenPythia::kParentSelection);
1147         gener->SetCountMode(AliGenPythia::kCountParents);
1148         gGener=gener;
1149       }
1150       break;
1151     case kCharmSemiElPbPb5500:
1152       {
1153         comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
1154         AliGenPythia * gener = new AliGenPythia(10);
1155         gener->SetProcess(kPyCharmPbPbMNR);
1156         gener->SetStrucFunc(kCTEQ4L);
1157         gener->SetPtHard(2.1,-1.0);
1158         gener->SetEnergyCMS(5500.);
1159         gener->SetNuclei(208,208);
1160         gener->SetForceDecay(kSemiElectronic);
1161         gener->SetYRange(-2,2);
1162         gener->SetFeedDownHigherFamily(kFALSE);
1163         gener->SetCountMode(AliGenPythia::kCountParents);
1164         gGener=gener;
1165       }
1166       break;
1167     case kBeautySemiElPbPb5500:
1168       {
1169         comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
1170         AliGenPythia *gener = new AliGenPythia(10);
1171         gener->SetProcess(kPyBeautyPbPbMNR);
1172         gener->SetStrucFunc(kCTEQ4L);
1173         gener->SetPtHard(2.75,-1.0);
1174         gener->SetEnergyCMS(5500.);
1175         gener->SetNuclei(208,208);
1176         gener->SetForceDecay(kSemiElectronic);
1177         gener->SetYRange(-2,2);
1178         gener->SetFeedDownHigherFamily(kFALSE);
1179         gener->SetCountMode(AliGenPythia::kCountParents);
1180         gGener=gener;
1181       }
1182       break;
1183     case kCocktailTRD:
1184       {
1185         comment = comment.Append(" Cocktail for TRD at 5.5 TeV");
1186         AliGenCocktail *gener  = new AliGenCocktail();
1187         
1188         AliGenParam *phi = new AliGenParam(10,
1189                                            new AliGenMUONlib(),
1190                                            AliGenMUONlib::kPhi,
1191                                            "Vogt PbPb");
1192
1193         phi->SetPtRange(0, 100);
1194         phi->SetYRange(-1., +1.);
1195         phi->SetForceDecay(kDiElectron);
1196
1197         AliGenParam *omega = new AliGenParam(10,
1198                                              new AliGenMUONlib(),
1199                                              AliGenMUONlib::kOmega,
1200                                              "Vogt PbPb");
1201
1202         omega->SetPtRange(0, 100);
1203         omega->SetYRange(-1., +1.);
1204         omega->SetForceDecay(kDiElectron);
1205         
1206         AliGenParam *jpsi = new AliGenParam(10,
1207                                             new AliGenMUONlib(),
1208                                             AliGenMUONlib::kJpsiFamily,
1209                                             "Vogt PbPb");
1210
1211         jpsi->SetPtRange(0, 100);
1212         jpsi->SetYRange(-1., +1.);
1213         jpsi->SetForceDecay(kDiElectron);
1214
1215         AliGenParam *ups = new AliGenParam(10,
1216                                            new AliGenMUONlib(),
1217                                            AliGenMUONlib::kUpsilonFamily,
1218                                            "Vogt PbPb");
1219         ups->SetPtRange(0, 100);
1220         ups->SetYRange(-1., +1.);
1221         ups->SetForceDecay(kDiElectron);
1222         
1223         AliGenParam *charm = new AliGenParam(10,
1224                                              new AliGenMUONlib(), 
1225                                              AliGenMUONlib::kCharm,
1226                                              "central");
1227         charm->SetPtRange(0, 100);
1228         charm->SetYRange(-1.5, +1.5);
1229         charm->SetForceDecay(kSemiElectronic);
1230         
1231         
1232         AliGenParam *beauty = new AliGenParam(10,
1233                                               new AliGenMUONlib(), 
1234                                               AliGenMUONlib::kBeauty,
1235                                               "central");
1236         beauty->SetPtRange(0, 100);
1237         beauty->SetYRange(-1.5, +1.5);
1238         beauty->SetForceDecay(kSemiElectronic);
1239
1240         AliGenParam *beautyJ = new AliGenParam(10,
1241                                                new AliGenMUONlib(), 
1242                                                AliGenMUONlib::kBeauty,
1243                                                "central");
1244         beautyJ->SetPtRange(0, 100);
1245         beautyJ->SetYRange(-1.5, +1.5);
1246         beautyJ->SetForceDecay(kBJpsiDiElectron);
1247
1248         gener->AddGenerator(phi,"Phi",1);
1249         gener->AddGenerator(omega,"Omega",1);
1250         gener->AddGenerator(jpsi,"J/psi",1);
1251         gener->AddGenerator(ups,"Upsilon",1);
1252         gener->AddGenerator(charm,"Charm",1);
1253         gener->AddGenerator(beauty,"Beauty",1);
1254         gener->AddGenerator(beautyJ,"J/Psi from Beauty",1);
1255         gGener=gener;
1256       }
1257       break;
1258     case kPyJJ:
1259       {
1260         comment = comment.Append(" Jet-jet at 5.5 TeV");
1261         AliGenPythia *gener = new AliGenPythia(-1);
1262         gener->SetEnergyCMS(5500.);
1263         gener->SetProcess(kPyJets);
1264         Double_t ptHardMin=10.0, ptHardMax=-1.0;
1265         gener->SetPtHard(ptHardMin,ptHardMax);
1266         gener->SetYHard(-0.7,0.7);
1267         gener->SetJetEtaRange(-0.2,0.2);
1268         gener->SetEventListRange(0,1);
1269         gGener=gener;
1270       }
1271       break;
1272     case kPyGJ:
1273       {
1274         comment = comment.Append(" Gamma-jet at 5.5 TeV");
1275         AliGenPythia *gener = new AliGenPythia(-1);
1276         gener->SetEnergyCMS(5500.);
1277         gener->SetProcess(kPyDirectGamma);
1278         Double_t ptHardMin=10.0, ptHardMax=-1.0;
1279         gener->SetPtHard(ptHardMin,ptHardMax);
1280         gener->SetYHard(-1.0,1.0);
1281         gener->SetGammaEtaRange(-0.13,0.13);
1282         gener->SetGammaPhiRange(210.,330.);
1283         gener->SetEventListRange(0,1);
1284         gGener=gener;
1285       }
1286       break;
1287     case kMuonCocktailCent1:
1288       {
1289         comment = comment.Append(" Muon Cocktail Cent1");
1290         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1291         gener->SetPtRange(0.4,100.);       // Transverse momentum range   
1292         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1293         gener->SetYRange(-4.0,-2.4);
1294         gener->SetMuonPtCut(0.8);
1295         gener->SetMuonThetaCut(171.,178.);
1296         gener->SetMuonMultiplicity(2);
1297         gener->SetImpactParameterRange(12.,16.);  //Centrality class Cent1 for PDC04
1298         gGener=gener;
1299       }
1300       break;
1301     case kMuonCocktailPer1:
1302       {
1303         comment = comment.Append(" Muon Cocktail Per1");
1304         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1305         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1306         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1307         gener->SetYRange(-4.0,-2.4);
1308         gener->SetMuonPtCut(0.8);
1309         gener->SetMuonThetaCut(171.,178.);
1310         gener->SetMuonMultiplicity(2);
1311         gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1312         gGener=gener;
1313       }
1314       break;
1315     case kMuonCocktailPer4:
1316       {
1317         comment = comment.Append(" Muon Cocktail Per4");
1318         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1319         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1320         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1321         gener->SetYRange(-4.0,-2.4);
1322         gener->SetMuonPtCut(0.8);
1323         gener->SetMuonThetaCut(171.,178.);
1324         gener->SetMuonMultiplicity(2);
1325         gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1326         gGener=gener;
1327       }
1328       break;
1329     case kMuonCocktailCent1HighPt:
1330       {
1331         comment = comment.Append(" Muon Cocktail HighPt Cent1");
1332         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1333         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1334         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1335         gener->SetYRange(-4.0,-2.4);
1336         gener->SetMuonPtCut(2.5);
1337         gener->SetMuonThetaCut(171.,178.);
1338         gener->SetMuonMultiplicity(2);
1339         gener->SetImpactParameterRange(0.,5.);  //Centrality class Cent1 for PDC04
1340         gGener=gener;
1341       }
1342       break;
1343     case kMuonCocktailPer1HighPt :
1344       {
1345         comment = comment.Append(" Muon Cocktail HighPt Per1");
1346         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1347         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1348         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1349         gener->SetYRange(-4.0,-2.4);
1350         gener->SetMuonPtCut(2.5);
1351         gener->SetMuonThetaCut(171.,178.);
1352         gener->SetMuonMultiplicity(2);
1353         gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1354         gGener=gener;
1355       }
1356       break;
1357     case kMuonCocktailPer4HighPt:
1358       {
1359         comment = comment.Append(" Muon Cocktail HighPt Per4");
1360         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1361         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1362         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1363         gener->SetYRange(-4.0,-2.4);
1364         gener->SetMuonPtCut(2.5);
1365         gener->SetMuonThetaCut(171.,178.);
1366         gener->SetMuonMultiplicity(2);
1367         gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1368         gGener=gener;
1369       }
1370       break;
1371     case kMuonCocktailCent1Single:
1372       {
1373         comment = comment.Append(" Muon Cocktail Single Cent1");
1374         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1375         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1376         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1377         gener->SetYRange(-4.0,-2.4);
1378         gener->SetMuonPtCut(0.8);
1379         gener->SetMuonThetaCut(171.,178.);
1380         gener->SetMuonMultiplicity(1);
1381         gener->SetImpactParameterRange(0.,5.);  //Centrality class Cent1 for PDC04
1382         gGener=gener;
1383       }
1384       break;
1385     case kMuonCocktailPer1Single :
1386       {
1387         comment = comment.Append(" Muon Cocktail Single Per1");
1388         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1389         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1390         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1391         gener->SetYRange(-4.0,-2.4);
1392         gener->SetMuonPtCut(0.8);
1393         gener->SetMuonThetaCut(171.,178.);
1394         gener->SetMuonMultiplicity(1);
1395         gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1396         gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1397         gGener=gener;
1398       }
1399       break;
1400     case kMuonCocktailPer4Single:
1401       {
1402         comment = comment.Append(" Muon Cocktail Single Per4");
1403         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1404         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1405         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1406         gener->SetYRange(-4.0,-2.4);
1407         gener->SetMuonPtCut(0.8);
1408         gener->SetMuonThetaCut(171.,178.);
1409         gener->SetMuonMultiplicity(1);
1410         gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1411         gGener=gener;
1412       }
1413       break;
1414     case kFlow_2_2000:
1415     {
1416         comment = comment.Append(" Flow with dN/deta  = 2000, vn = 2%");
1417         gGener = GeVSimStandard(2000., 2.);
1418     }
1419     break;
1420     
1421     case kFlow_10_2000:
1422     {
1423         comment = comment.Append(" Flow with dN/deta  = 2000, vn = 10%");
1424         gGener = GeVSimStandard(2000., 10.);
1425     }
1426     break;
1427     
1428     case kFlow_6_2000:
1429     {
1430         comment = comment.Append(" Flow with dN/deta  = 2000, vn = 6%");
1431         gGener = GeVSimStandard(2000., 6.);
1432     }
1433     break;
1434     
1435     case kFlow_6_5000:
1436     {
1437         comment = comment.Append(" Flow with dN/deta  = 5000, vn = 6%");
1438         gGener = GeVSimStandard(5000., 6.);
1439     }
1440     break;
1441     case kHIJINGplus:
1442     {
1443         //
1444         // The cocktail
1445         AliGenCocktail *gener  = new AliGenCocktail();
1446
1447         //
1448         // Charm production by Pythia
1449         AliGenPythia * genpyc = new AliGenPythia(230);
1450         genpyc->SetProcess(kPyCharmPbPbMNR);
1451         genpyc->SetStrucFunc(kCTEQ4L);
1452         genpyc->SetPtHard(2.1,-1.0);
1453         genpyc->SetEnergyCMS(5500.);
1454         genpyc->SetNuclei(208,208);
1455         genpyc->SetYRange(-999,999);
1456         genpyc->SetForceDecay(kAll);
1457         genpyc->SetFeedDownHigherFamily(kFALSE);
1458         genpyc->SetCountMode(AliGenPythia::kCountParents);
1459         //
1460         // Beauty production by Pythia
1461         AliGenPythia * genpyb = new AliGenPythia(9);
1462         genpyb->SetProcess(kPyBeautyPbPbMNR);
1463         genpyb->SetStrucFunc(kCTEQ4L);
1464         genpyb->SetPtHard(2.75,-1.0);
1465         genpyb->SetEnergyCMS(5500.);
1466         genpyb->SetNuclei(208,208);
1467         genpyb->SetYRange(-999,999);
1468         genpyb->SetForceDecay(kAll);
1469         genpyb->SetFeedDownHigherFamily(kFALSE);
1470         genpyb->SetCountMode(AliGenPythia::kCountParents);
1471         //
1472         // Hyperons
1473         //
1474         AliGenSTRANGElib *lib = new AliGenSTRANGElib();
1475         Int_t particle;
1476         // Xi
1477         particle = AliGenSTRANGElib::kXiMinus;
1478         AliGenParam *genXi = new AliGenParam(16,particle,lib->GetPt(particle),lib->GetY(particle),lib->GetIp(particle));
1479         genXi->SetPtRange(0., 12.);
1480         genXi->SetYRange(-1.1, 1.1);
1481         genXi->SetForceDecay(kNoDecay); 
1482  
1483         //
1484         // Omega
1485         particle = AliGenSTRANGElib::kOmegaMinus;
1486         AliGenParam *genOmega = new AliGenParam(10,particle,lib->GetPt(particle),lib->GetY(particle),lib->GetIp(particle));     
1487         genOmega->SetPtRange(0, 12.);
1488         genOmega->SetYRange(-1.1, 1.1);
1489         genOmega->SetForceDecay(kNoDecay);
1490         
1491         //
1492         // Central Hijing 
1493         AliGenHijing *genHi = HijingStandard();
1494         genHi->SwitchOffHeavyQuarks(kTRUE);
1495         genHi->SetImpactParameterRange(0.,5.);
1496         //
1497         // Add everything to the cocktail and shake ...
1498         gener->AddGenerator(genHi,    "Hijing cent1", 1);
1499         gener->AddGenerator(genpyc,   "Extra charm",  1);
1500         gener->AddGenerator(genpyb,   "Extra beauty", 1);
1501         gener->AddGenerator(genXi,    "Xi"          , 1);
1502         gener->AddGenerator(genOmega, "Omega",        1);
1503         gGener = gener;
1504     }
1505     break;
1506     default: break;
1507     }
1508     
1509     return gGener;
1510 }
1511
1512 AliGenHijing* HijingStandard()
1513 {
1514     AliGenHijing *gener = new AliGenHijing(-1);
1515 // centre of mass energy 
1516     gener->SetEnergyCMS(5500.);
1517 // reference frame
1518     gener->SetReferenceFrame("CMS");
1519 // projectile
1520      gener->SetProjectile("A", 208, 82);
1521      gener->SetTarget    ("A", 208, 82);
1522 // tell hijing to keep the full parent child chain
1523      gener->KeepFullEvent();
1524 // enable jet quenching
1525      gener->SetJetQuenching(1);
1526 // enable shadowing
1527      gener->SetShadowing(1);
1528 // neutral pion and heavy particle decays switched off
1529      gener->SetDecaysOff(1);
1530 // Don't track spectators
1531      gener->SetSpectators(0);
1532 // kinematic selection
1533      gener->SetSelectAll(0);
1534      return gener;
1535 }
1536
1537 AliGenGeVSim* GeVSimStandard(Float_t mult, Float_t vn)
1538 {
1539     AliGenGeVSim* gener = new AliGenGeVSim(0);
1540 //
1541 // Mult is the number of charged particles in |eta| < 0.5
1542 // Vn is in (%)
1543 //
1544 // Sigma of the Gaussian dN/deta
1545     Float_t sigma_eta  = 2.75;
1546 //
1547 // Maximum eta
1548     Float_t etamax     = 7.00;
1549 //
1550 //
1551 // Scale from multiplicity in |eta| < 0.5 to |eta| < |etamax|   
1552     Float_t mm = mult * (TMath::Erf(etamax/sigma_eta/sqrt(2.)) / TMath::Erf(0.5/sigma_eta/sqrt(2.))); 
1553 //
1554 // Scale from charged to total multiplicity
1555 // 
1556     mm *= 1.587;
1557 //
1558 // Vn 
1559     vn /= 100.;          
1560 //
1561 // Define particles
1562 //
1563 //
1564 // 78% Pions (26% pi+, 26% pi-, 26% p0)              T = 250 MeV
1565     AliGeVSimParticle *pp =  new AliGeVSimParticle(kPiPlus,  1, 0.26 * mm, 0.25, sigma_eta) ;
1566     AliGeVSimParticle *pm =  new AliGeVSimParticle(kPiMinus, 1, 0.26 * mm, 0.25, sigma_eta) ;
1567     AliGeVSimParticle *p0 =  new AliGeVSimParticle(kPi0,     1, 0.26 * mm, 0.25, sigma_eta) ;
1568 //
1569 // 12% Kaons (3% K0short, 3% K0long, 3% K+, 3% K-)   T = 300 MeV
1570     AliGeVSimParticle *ks =  new AliGeVSimParticle(kK0Short, 1, 0.03 * mm, 0.30, sigma_eta) ;
1571     AliGeVSimParticle *kl =  new AliGeVSimParticle(kK0Long,  1, 0.03 * mm, 0.30, sigma_eta) ;
1572     AliGeVSimParticle *kp =  new AliGeVSimParticle(kKPlus,   1, 0.03 * mm, 0.30, sigma_eta) ;
1573     AliGeVSimParticle *km =  new AliGeVSimParticle(kKMinus,  1, 0.03 * mm, 0.30, sigma_eta) ;
1574 //
1575 // 10% Protons / Neutrons (5% Protons, 5% Neutrons)  T = 250 MeV
1576     AliGeVSimParticle *pr =  new AliGeVSimParticle(kProton,  1, 0.05 * mm, 0.25, sigma_eta) ;
1577     AliGeVSimParticle *ne =  new AliGeVSimParticle(kNeutron, 1, 0.05 * mm, 0.25, sigma_eta) ;
1578 //
1579
1580     Float_t pTsaturation = 2. ;
1581
1582 // Set Elliptic Flow properties         
1583     pp->SetEllipticParam(vn,pTsaturation,0.) ;
1584     pm->SetEllipticParam(vn,pTsaturation,0.) ;
1585     p0->SetEllipticParam(vn,pTsaturation,0.) ;
1586     pr->SetEllipticParam(vn,pTsaturation,0.) ;
1587     ne->SetEllipticParam(vn,pTsaturation,0.) ;
1588     ks->SetEllipticParam(vn,pTsaturation,0.) ;
1589     kl->SetEllipticParam(vn,pTsaturation,0.) ;
1590     kp->SetEllipticParam(vn,pTsaturation,0.) ;
1591     km->SetEllipticParam(vn,pTsaturation,0.) ;
1592 //
1593 // Set Direct Flow properties   
1594     pp->SetDirectedParam(vn,pTsaturation,0.) ;
1595     pm->SetDirectedParam(vn,pTsaturation,0.) ;
1596     p0->SetDirectedParam(vn,pTsaturation,0.) ;
1597     pr->SetDirectedParam(vn,pTsaturation,0.) ;
1598     ne->SetDirectedParam(vn,pTsaturation,0.) ;
1599     ks->SetDirectedParam(vn,pTsaturation,0.) ;
1600     kl->SetDirectedParam(vn,pTsaturation,0.) ;
1601     kp->SetDirectedParam(vn,pTsaturation,0.) ;
1602     km->SetDirectedParam(vn,pTsaturation,0.) ;
1603 //
1604 // Add particles to the list
1605     gener->AddParticleType(pp) ;
1606     gener->AddParticleType(pm) ;
1607     gener->AddParticleType(p0) ;
1608     gener->AddParticleType(pr) ;
1609     gener->AddParticleType(ne) ;
1610     gener->AddParticleType(ks) ;
1611     gener->AddParticleType(kl) ;
1612     gener->AddParticleType(kp) ;
1613     gener->AddParticleType(km) ;
1614 //      
1615 // Random Ev.Plane ----------------------------------
1616     TF1 *rpa = new TF1("gevsimPsiRndm","1", 0, 360);
1617 // --------------------------------------------------
1618     gener->SetPtRange(0., 9.) ; // Use a resonable range! (used for bin size in numerical integration)
1619     gener->SetPhiRange(0, 360);
1620     //
1621     // Set pseudorapidity range 
1622     Float_t thmin = EtaToTheta(+etamax);   
1623     Float_t thmax = EtaToTheta(-etamax);   
1624     gener->SetThetaRange(thmin,thmax);     
1625     return gener;
1626 }
1627
1628
1629
1630 void ProcessEnvironmentVars()
1631 {
1632     // Run type
1633     if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
1634       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
1635         if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
1636           srun = (PprRun_t)iRun;
1637           cout<<"Run type set to "<<pprRunName[iRun]<<endl;
1638         }
1639       }
1640     }
1641
1642     // Random Number seed
1643     if (gSystem->Getenv("CONFIG_SEED")) {
1644       sseed = atoi(gSystem->Getenv("CONFIG_SEED"));
1645     }
1646 }