]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/ConfigPPR.C
Tracking in non-uniform nmagnetic field (Yu.Belikov)
[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->SetSegmentationType(2);// default wise to old (1), new (2)
503         MUON->AddGeometryBuilder(new AliMUONSt1GeometryBuilderV2(MUON));
504         MUON->AddGeometryBuilder(new AliMUONSt2GeometryBuilderV2(MUON));
505         MUON->AddGeometryBuilder(new AliMUONSlatGeometryBuilder(MUON));
506         MUON->AddGeometryBuilder(new AliMUONTriggerGeometryBuilder(MUON));
507     }
508     //=================== PHOS parameters ===========================
509
510     if (iPHOS)
511     {
512         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
513     }
514
515
516     if (iPMD)
517     {
518         //=================== PMD parameters ============================
519         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
520     }
521
522     if (iSTART)
523     {
524         //=================== START parameters ============================
525         AliSTART *START = new AliSTARTv1("START", "START Detector");
526     }
527
528     if (iEMCAL)
529     {
530         //=================== EMCAL parameters ============================
531         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCAL_55_25");
532     }
533
534      if (iCRT)
535     {
536         //=================== CRT parameters ============================
537         AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
538     }
539
540      if (iVZERO)
541     {
542         //=================== CRT parameters ============================
543         AliVZERO *VZERO = new AliVZEROv5("VZERO", "normal VZERO");
544     }
545  
546              
547 }
548
549 Float_t EtaToTheta(Float_t arg){
550   return (180./TMath::Pi())*2.*atan(exp(-arg));
551 }
552
553
554
555 AliGenerator* GeneratorFactory(PprRun_t srun) {
556     Int_t isw = 3;
557     if (srad == kNoGluonRadiation) isw = 0;
558     
559
560     AliGenerator * gGener = 0x0;
561     switch (srun) {
562     case test50:
563       {
564         comment = comment.Append(":HIJINGparam test 50 particles");
565         AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
566         gener->SetMomentumRange(0, 999999.);
567         gener->SetPhiRange(0., 360.);
568         // Set pseudorapidity range from -8 to 8.
569         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
570         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
571         gener->SetThetaRange(thmin,thmax);
572         gGener=gener;
573       }
574       break;
575     case kParam_8000:
576       {
577         comment = comment.Append(":HIJINGparam N=8000");
578         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
579         gener->SetMomentumRange(0, 999999.);
580         gener->SetPhiRange(0., 360.);
581         // Set pseudorapidity range from -8 to 8.
582         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
583         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
584         gener->SetThetaRange(thmin,thmax);
585         gGener=gener;
586       }
587       break;
588     case kParam_4000:
589       {
590         comment = comment.Append("HIJINGparam N=4000");
591         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
592         gener->SetMomentumRange(0, 999999.);
593         gener->SetPhiRange(0., 360.);
594         // Set pseudorapidity range from -8 to 8.
595         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
596         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
597         gener->SetThetaRange(thmin,thmax);
598         gGener=gener;
599       }
600         break;
601     case kParam_2000:
602       {
603         comment = comment.Append("HIJINGparam N=2000");
604         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
605         gener->SetMomentumRange(0, 999999.);
606         gener->SetPhiRange(0., 360.);
607         // Set pseudorapidity range from -8 to 8.
608         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
609         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
610         gener->SetThetaRange(thmin,thmax);
611         gGener=gener;
612       }
613       break;
614 //
615 //  Hijing Central
616 //
617     case kHijing_cent1:
618       {
619         comment = comment.Append("HIJING cent1");
620         AliGenHijing *gener = HijingStandard();
621 // impact parameter range
622         gener->SetImpactParameterRange(0., 5.);
623         gGener=gener;
624       }
625       break;
626     case kHijing_cent2:
627       {
628         comment = comment.Append("HIJING cent2");
629         AliGenHijing *gener = HijingStandard();
630 // impact parameter range
631         gener->SetImpactParameterRange(0., 2.);
632         gGener=gener;
633       }
634       break;
635 //
636 // Hijing Peripheral 
637 //
638     case kHijing_per1:
639       {
640         comment = comment.Append("HIJING per1");
641         AliGenHijing *gener = HijingStandard();
642 // impact parameter range
643         gener->SetImpactParameterRange(5., 8.6);
644         gGener=gener;
645       }
646       break;
647     case kHijing_per2:
648       {
649         comment = comment.Append("HIJING per2");
650         AliGenHijing *gener = HijingStandard();
651 // impact parameter range
652         gener->SetImpactParameterRange(8.6, 11.2);
653         gGener=gener;
654       }
655       break;
656     case kHijing_per3:
657       {
658         comment = comment.Append("HIJING per3");
659         AliGenHijing *gener = HijingStandard();
660 // impact parameter range
661         gener->SetImpactParameterRange(11.2, 13.2);
662         gGener=gener;
663       }
664       break;
665     case kHijing_per4:
666       {
667         comment = comment.Append("HIJING per4");
668         AliGenHijing *gener = HijingStandard();
669 // impact parameter range
670         gener->SetImpactParameterRange(13.2, 15.);
671         gGener=gener;
672       }
673       break;
674     case kHijing_per5:
675       {
676         comment = comment.Append("HIJING per5");
677         AliGenHijing *gener = HijingStandard();
678 // impact parameter range
679         gener->SetImpactParameterRange(15., 100.);
680         gGener=gener;
681       }
682       break;
683 //
684 //  Jet-Jet
685 //
686     case kHijing_jj25:
687       {
688         comment = comment.Append("HIJING Jet 25 GeV");
689         AliGenHijing *gener = HijingStandard();
690 // impact parameter range
691         gener->SetImpactParameterRange(0., 5.);
692         // trigger
693         gener->SetTrigger(1);
694         gener->SetPtJet(25.);
695         gener->SetRadiation(isw);
696         gener->SetSimpleJets(!isw);
697         gener->SetJetEtaRange(-0.3,0.3);
698         gener->SetJetPhiRange(75., 165.);   
699         gGener=gener;
700       }
701       break;
702
703     case kHijing_jj50:
704       {
705         comment = comment.Append("HIJING Jet 50 GeV");
706         AliGenHijing *gener = HijingStandard();
707 // impact parameter range
708         gener->SetImpactParameterRange(0., 5.);
709         // trigger
710         gener->SetTrigger(1);
711         gener->SetPtJet(50.);
712         gener->SetRadiation(isw);
713         gener->SetSimpleJets(!isw);
714         gener->SetJetEtaRange(-0.3,0.3);
715         gener->SetJetPhiRange(75., 165.);   
716         gGener=gener;
717       }
718         break;
719
720     case kHijing_jj75:
721       {
722         comment = comment.Append("HIJING Jet 75 GeV");
723         AliGenHijing *gener = HijingStandard();
724 // impact parameter range
725         gener->SetImpactParameterRange(0., 5.);
726         // trigger
727         gener->SetTrigger(1);
728         gener->SetPtJet(75.);
729         gener->SetRadiation(isw);
730         gener->SetSimpleJets(!isw);
731         gener->SetJetEtaRange(-0.3,0.3);
732         gener->SetJetPhiRange(75., 165.);   
733         gGener=gener;
734       }
735       break;
736
737     case kHijing_jj100:
738       {
739         comment = comment.Append("HIJING Jet 100 GeV");
740         AliGenHijing *gener = HijingStandard();
741 // impact parameter range
742         gener->SetImpactParameterRange(0., 5.);
743         // trigger
744         gener->SetTrigger(1);
745         gener->SetPtJet(100.);
746         gener->SetRadiation(isw);
747         gener->SetSimpleJets(!isw);
748         gener->SetJetEtaRange(-0.3,0.3);
749         gener->SetJetPhiRange(75., 165.);   
750         gGener=gener;
751       }
752       break;
753
754     case kHijing_jj200:
755       {
756         comment = comment.Append("HIJING Jet 200 GeV");
757         AliGenHijing *gener = HijingStandard();
758 // impact parameter range
759         gener->SetImpactParameterRange(0., 5.);
760         // trigger
761         gener->SetTrigger(1);
762         gener->SetPtJet(200.);
763         gener->SetRadiation(isw);
764         gener->SetSimpleJets(!isw);
765         gener->SetJetEtaRange(-0.3,0.3);
766         gener->SetJetPhiRange(75., 165.);   
767         gGener=gener;
768       }
769       break;
770 //
771 // Gamma-Jet
772 //
773     case kHijing_gj25:
774       {
775         comment = comment.Append("HIJING Gamma 25 GeV");
776         AliGenHijing *gener = HijingStandard();
777 // impact parameter range
778         gener->SetImpactParameterRange(0., 5.);
779         // trigger
780         gener->SetTrigger(2);
781         gener->SetPtJet(25.);
782         gener->SetRadiation(isw);
783         gener->SetSimpleJets(!isw);
784         gener->SetJetEtaRange(-0.12, 0.12);
785         gener->SetJetPhiRange(220., 320.);
786         gGener=gener;
787       }
788       break;
789
790     case kHijing_gj50:
791       {
792         comment = comment.Append("HIJING Gamma 50 GeV");
793         AliGenHijing *gener = HijingStandard();
794 // impact parameter range
795         gener->SetImpactParameterRange(0., 5.);
796         // trigger
797         gener->SetTrigger(2);
798         gener->SetPtJet(50.);
799         gener->SetRadiation(isw);
800         gener->SetSimpleJets(!isw);
801         gener->SetJetEtaRange(-0.12, 0.12);
802         gener->SetJetPhiRange(220., 320.);
803         gGener=gener;
804       }
805       break;
806
807     case kHijing_gj75:
808       {
809         comment = comment.Append("HIJING Gamma 75 GeV");
810         AliGenHijing *gener = HijingStandard();
811 // impact parameter range
812         gener->SetImpactParameterRange(0., 5.);
813         // trigger
814         gener->SetTrigger(2);
815         gener->SetPtJet(75.);
816         gener->SetRadiation(isw);
817         gener->SetSimpleJets(!isw);
818         gener->SetJetEtaRange(-0.12, 0.12);
819         gener->SetJetPhiRange(220., 320.);
820         gGener=gener;
821       }
822       break;
823
824     case kHijing_gj100:
825       {
826         comment = comment.Append("HIJING Gamma 100 GeV");
827         AliGenHijing *gener = HijingStandard();
828 // impact parameter range
829         gener->SetImpactParameterRange(0., 5.);
830         // trigger
831         gener->SetTrigger(2);
832         gener->SetPtJet(100.);
833         gener->SetRadiation(isw);
834         gener->SetSimpleJets(!isw);
835         gener->SetJetEtaRange(-0.12, 0.12);
836         gener->SetJetPhiRange(220., 320.);
837         gGener=gener;
838       }
839       break;
840
841     case kHijing_gj200:
842       {
843         comment = comment.Append("HIJING Gamma 200 GeV");
844         AliGenHijing *gener = HijingStandard();
845 // impact parameter range
846         gener->SetImpactParameterRange(0., 5.);
847         // trigger
848         gener->SetTrigger(2);
849         gener->SetPtJet(200.);
850         gener->SetRadiation(isw);
851         gener->SetSimpleJets(!isw);
852         gener->SetJetEtaRange(-0.12, 0.12);
853         gener->SetJetPhiRange(220., 320.);
854         gGener=gener;
855       }
856       break;
857     case kHijing_pA:
858       {
859         comment = comment.Append("HIJING pA");
860
861         AliGenCocktail *gener  = new AliGenCocktail();
862
863         AliGenHijing   *hijing = new AliGenHijing(-1);
864 // centre of mass energy 
865         hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
866 // impact parameter range
867         hijing->SetImpactParameterRange(0., 15.);
868 // reference frame
869         hijing->SetReferenceFrame("CMS");
870         hijing->SetBoostLHC(1);
871 // projectile
872         hijing->SetProjectile("P", 1, 1);
873         hijing->SetTarget    ("A", 208, 82);
874 // tell hijing to keep the full parent child chain
875         hijing->KeepFullEvent();
876 // enable jet quenching
877         hijing->SetJetQuenching(0);
878 // enable shadowing
879         hijing->SetShadowing(1);
880 // Don't track spectators
881         hijing->SetSpectators(0);
882 // kinematic selection
883         hijing->SetSelectAll(0);
884 //
885         AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
886         AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
887         gray->SetSlowNucleonModel(model);
888         gray->SetDebug(1);
889         gener->AddGenerator(hijing,"Hijing pPb", 1);
890         gener->AddGenerator(gray,  "Gray Particles",1);
891         gGener=gener;
892       }
893       break;
894     case kPythia6:
895       {
896         comment = comment.Append(":Pythia p-p @ 14 TeV");
897         AliGenPythia *gener = new AliGenPythia(-1); 
898         gener->SetMomentumRange(0,999999);
899         gener->SetThetaRange(0., 180.);
900         gener->SetYRange(-12,12);
901         gener->SetPtRange(0,1000);
902         gener->SetProcess(kPyMb);
903         gener->SetEnergyCMS(14000.);
904         gGener=gener;
905       }
906       break;
907     case kPythia6Jets20_24:
908       {
909         comment = comment.Append(":Pythia jets 20-24 GeV @ 5.5 TeV");
910         AliGenPythia * gener = new AliGenPythia(-1);
911         gener->SetEnergyCMS(5500.);//        Centre of mass energy
912         gener->SetProcess(kPyJets);//        Process type
913         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
914         gener->SetJetPhiRange(0., 360.);
915         gener->SetJetEtRange(10., 1000.);
916         gener->SetGluonRadiation(1,1);
917         //    gener->SetPtKick(0.);
918         //   Structure function
919         gener->SetStrucFunc(kCTEQ4L);
920         gener->SetPtHard(20., 24.);// Pt transfer of the hard scattering
921         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
922         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
923         gGener=gener;
924       }
925       break;
926     case kPythia6Jets24_29:
927       {
928         comment = comment.Append(":Pythia jets 24-29 GeV @ 5.5 TeV");
929         AliGenPythia * gener = new AliGenPythia(-1);
930         gener->SetEnergyCMS(5500.);//        Centre of mass energy
931         gener->SetProcess(kPyJets);//        Process type
932         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
933         gener->SetJetPhiRange(0., 360.);
934         gener->SetJetEtRange(10., 1000.);
935         gener->SetGluonRadiation(1,1);
936         //    gener->SetPtKick(0.);
937         //   Structure function
938         gener->SetStrucFunc(kCTEQ4L);
939         gener->SetPtHard(24., 29.);// Pt transfer of the hard scattering
940         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
941         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
942         gGener=gener;
943       }
944       break;
945     case kPythia6Jets29_35:
946       {
947         comment = comment.Append(":Pythia jets 29-35 GeV @ 5.5 TeV");
948         AliGenPythia * gener = new AliGenPythia(-1);
949         gener->SetEnergyCMS(5500.);//        Centre of mass energy
950         gener->SetProcess(kPyJets);//        Process type
951         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
952         gener->SetJetPhiRange(0., 360.);
953         gener->SetJetEtRange(10., 1000.);
954         gener->SetGluonRadiation(1,1);
955         //    gener->SetPtKick(0.);
956         //   Structure function
957         gener->SetStrucFunc(kCTEQ4L);
958         gener->SetPtHard(29., 35.);// Pt transfer of the hard scattering
959         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
960         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
961         gGener=gener;
962       }
963       break;
964     case kPythia6Jets35_42:
965       {
966         comment = comment.Append(":Pythia jets 35-42 GeV @ 5.5 TeV");
967         AliGenPythia * gener = new AliGenPythia(-1);
968         gener->SetEnergyCMS(5500.);//        Centre of mass energy
969         gener->SetProcess(kPyJets);//        Process type
970         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
971         gener->SetJetPhiRange(0., 360.);
972         gener->SetJetEtRange(10., 1000.);
973         gener->SetGluonRadiation(1,1);
974         //    gener->SetPtKick(0.);
975         //   Structure function
976         gener->SetStrucFunc(kCTEQ4L);
977         gener->SetPtHard(35., 42.);// Pt transfer of the hard scattering
978         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
979         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
980         gGener=gener;
981       }
982       break;
983     case kPythia6Jets42_50:
984       {
985         comment = comment.Append(":Pythia jets 42-50 GeV @ 5.5 TeV");
986         AliGenPythia * gener = new AliGenPythia(-1);
987         gener->SetEnergyCMS(5500.);//        Centre of mass energy
988         gener->SetProcess(kPyJets);//        Process type
989         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
990         gener->SetJetPhiRange(0., 360.);
991         gener->SetJetEtRange(10., 1000.);
992         gener->SetGluonRadiation(1,1);
993         //    gener->SetPtKick(0.);
994         //   Structure function
995         gener->SetStrucFunc(kCTEQ4L);
996         gener->SetPtHard(42., 50.);// Pt transfer of the hard scattering
997         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
998         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
999         gGener=gener;
1000       }
1001       break;
1002     case kPythia6Jets50_60:
1003       {
1004         comment = comment.Append(":Pythia jets 50-60 GeV @ 5.5 TeV");
1005         AliGenPythia * gener = new AliGenPythia(-1);
1006         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1007         gener->SetProcess(kPyJets);//        Process type
1008         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1009         gener->SetJetPhiRange(0., 360.);
1010         gener->SetJetEtRange(10., 1000.);
1011         gener->SetGluonRadiation(1,1);
1012         //    gener->SetPtKick(0.);
1013         //   Structure function
1014         gener->SetStrucFunc(kCTEQ4L);
1015         gener->SetPtHard(50., 60.);// Pt transfer of the hard scattering
1016         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1017         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1018         gGener=gener;
1019       }
1020       break;
1021     case kPythia6Jets60_72:
1022       {
1023         comment = comment.Append(":Pythia jets 60-72 GeV @ 5.5 TeV");
1024         AliGenPythia * gener = new AliGenPythia(-1);
1025         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1026         gener->SetProcess(kPyJets);//        Process type
1027         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1028         gener->SetJetPhiRange(0., 360.);
1029         gener->SetJetEtRange(10., 1000.);
1030         gener->SetGluonRadiation(1,1);
1031         //    gener->SetPtKick(0.);
1032         //   Structure function
1033         gener->SetStrucFunc(kCTEQ4L);
1034         gener->SetPtHard(60., 72.);// Pt transfer of the hard scattering
1035         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1036         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1037         gGener=gener;
1038       }
1039       break;
1040     case kPythia6Jets72_86:
1041       {
1042         comment = comment.Append(":Pythia jets 72-86 GeV @ 5.5 TeV");
1043         AliGenPythia * gener = new AliGenPythia(-1);
1044         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1045         gener->SetProcess(kPyJets);//        Process type
1046         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1047         gener->SetJetPhiRange(0., 360.);
1048         gener->SetJetEtRange(10., 1000.);
1049         gener->SetGluonRadiation(1,1);
1050         //    gener->SetPtKick(0.);
1051         //   Structure function
1052         gener->SetStrucFunc(kCTEQ4L);
1053         gener->SetPtHard(72., 86.);// Pt transfer of the hard scattering
1054         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1055         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1056         gGener=gener;
1057       }
1058       break;
1059     case kPythia6Jets86_104:
1060       {
1061         comment = comment.Append(":Pythia jets 86-104 GeV @ 5.5 TeV");
1062         AliGenPythia * gener = new AliGenPythia(-1);
1063         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1064         gener->SetProcess(kPyJets);//        Process type
1065         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1066         gener->SetJetPhiRange(0., 360.);
1067         gener->SetJetEtRange(10., 1000.);
1068         gener->SetGluonRadiation(1,1);
1069         //    gener->SetPtKick(0.);
1070         //   Structure function
1071         gener->SetStrucFunc(kCTEQ4L);
1072         gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
1073         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1074         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1075         gGener=gener;
1076       }
1077       break;
1078     case kPythia6Jets104_125:
1079       {
1080         comment = comment.Append(":Pythia jets 105-125 GeV @ 5.5 TeV");
1081         AliGenPythia * gener = new AliGenPythia(-1);
1082         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1083         gener->SetProcess(kPyJets);//        Process type
1084         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1085         gener->SetJetPhiRange(0., 360.);
1086         gener->SetJetEtRange(10., 1000.);
1087         gener->SetGluonRadiation(1,1);
1088         //    gener->SetPtKick(0.);
1089         //   Structure function
1090         gener->SetStrucFunc(kCTEQ4L);
1091         gener->SetPtHard(104., 125.);// Pt transfer of the hard scattering
1092         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1093         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1094         gGener=gener;
1095       }
1096       break;
1097     case kPythia6Jets125_150:
1098       {
1099         comment = comment.Append(":Pythia jets 125-150 GeV @ 5.5 TeV");
1100         AliGenPythia * gener = new AliGenPythia(-1);
1101         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1102         gener->SetProcess(kPyJets);//        Process type
1103         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1104         gener->SetJetPhiRange(0., 360.);
1105         gener->SetJetEtRange(10., 1000.);
1106         gener->SetGluonRadiation(1,1);
1107         //    gener->SetPtKick(0.);
1108         //   Structure function
1109         gener->SetStrucFunc(kCTEQ4L);
1110         gener->SetPtHard(125., 150.);// Pt transfer of the hard scattering
1111         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1112         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1113         gGener=gener;
1114       }
1115       break;
1116     case kPythia6Jets150_180:
1117       {
1118         comment = comment.Append(":Pythia jets 150-180 GeV @ 5.5 TeV");
1119         AliGenPythia * gener = new AliGenPythia(-1);
1120         gener->SetEnergyCMS(5500.);//        Centre of mass energy
1121         gener->SetProcess(kPyJets);//        Process type
1122         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1123         gener->SetJetPhiRange(0., 360.);
1124         gener->SetJetEtRange(10., 1000.);
1125         gener->SetGluonRadiation(1,1);
1126         //    gener->SetPtKick(0.);
1127         //   Structure function
1128         gener->SetStrucFunc(kCTEQ4L);
1129         gener->SetPtHard(150., 180.);// Pt transfer of the hard scattering
1130         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1131         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1132         gGener=gener;
1133       }
1134       break;
1135     case kD0PbPb5500:
1136       {
1137         comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
1138         AliGenPythia * gener = new AliGenPythia(10);
1139         gener->SetProcess(kPyD0PbPbMNR);
1140         gener->SetStrucFunc(kCTEQ4L);
1141         gener->SetPtHard(2.1,-1.0);
1142         gener->SetEnergyCMS(5500.);
1143         gener->SetNuclei(208,208);
1144         gener->SetForceDecay(kHadronicD);
1145         gener->SetYRange(-2,2);
1146         gener->SetFeedDownHigherFamily(kFALSE);
1147         gener->SetStackFillOpt(AliGenPythia::kParentSelection);
1148         gener->SetCountMode(AliGenPythia::kCountParents);
1149         gGener=gener;
1150       }
1151       break;
1152     case kCharmSemiElPbPb5500:
1153       {
1154         comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
1155         AliGenPythia * gener = new AliGenPythia(10);
1156         gener->SetProcess(kPyCharmPbPbMNR);
1157         gener->SetStrucFunc(kCTEQ4L);
1158         gener->SetPtHard(2.1,-1.0);
1159         gener->SetEnergyCMS(5500.);
1160         gener->SetNuclei(208,208);
1161         gener->SetForceDecay(kSemiElectronic);
1162         gener->SetYRange(-2,2);
1163         gener->SetFeedDownHigherFamily(kFALSE);
1164         gener->SetCountMode(AliGenPythia::kCountParents);
1165         gGener=gener;
1166       }
1167       break;
1168     case kBeautySemiElPbPb5500:
1169       {
1170         comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
1171         AliGenPythia *gener = new AliGenPythia(10);
1172         gener->SetProcess(kPyBeautyPbPbMNR);
1173         gener->SetStrucFunc(kCTEQ4L);
1174         gener->SetPtHard(2.75,-1.0);
1175         gener->SetEnergyCMS(5500.);
1176         gener->SetNuclei(208,208);
1177         gener->SetForceDecay(kSemiElectronic);
1178         gener->SetYRange(-2,2);
1179         gener->SetFeedDownHigherFamily(kFALSE);
1180         gener->SetCountMode(AliGenPythia::kCountParents);
1181         gGener=gener;
1182       }
1183       break;
1184     case kCocktailTRD:
1185       {
1186         comment = comment.Append(" Cocktail for TRD at 5.5 TeV");
1187         AliGenCocktail *gener  = new AliGenCocktail();
1188         
1189         AliGenParam *phi = new AliGenParam(10,
1190                                            new AliGenMUONlib(),
1191                                            AliGenMUONlib::kPhi,
1192                                            "Vogt PbPb");
1193
1194         phi->SetPtRange(0, 100);
1195         phi->SetYRange(-1., +1.);
1196         phi->SetForceDecay(kDiElectron);
1197
1198         AliGenParam *omega = new AliGenParam(10,
1199                                              new AliGenMUONlib(),
1200                                              AliGenMUONlib::kOmega,
1201                                              "Vogt PbPb");
1202
1203         omega->SetPtRange(0, 100);
1204         omega->SetYRange(-1., +1.);
1205         omega->SetForceDecay(kDiElectron);
1206         
1207         AliGenParam *jpsi = new AliGenParam(10,
1208                                             new AliGenMUONlib(),
1209                                             AliGenMUONlib::kJpsiFamily,
1210                                             "Vogt PbPb");
1211
1212         jpsi->SetPtRange(0, 100);
1213         jpsi->SetYRange(-1., +1.);
1214         jpsi->SetForceDecay(kDiElectron);
1215
1216         AliGenParam *ups = new AliGenParam(10,
1217                                            new AliGenMUONlib(),
1218                                            AliGenMUONlib::kUpsilonFamily,
1219                                            "Vogt PbPb");
1220         ups->SetPtRange(0, 100);
1221         ups->SetYRange(-1., +1.);
1222         ups->SetForceDecay(kDiElectron);
1223         
1224         AliGenParam *charm = new AliGenParam(10,
1225                                              new AliGenMUONlib(), 
1226                                              AliGenMUONlib::kCharm,
1227                                              "central");
1228         charm->SetPtRange(0, 100);
1229         charm->SetYRange(-1.5, +1.5);
1230         charm->SetForceDecay(kSemiElectronic);
1231         
1232         
1233         AliGenParam *beauty = new AliGenParam(10,
1234                                               new AliGenMUONlib(), 
1235                                               AliGenMUONlib::kBeauty,
1236                                               "central");
1237         beauty->SetPtRange(0, 100);
1238         beauty->SetYRange(-1.5, +1.5);
1239         beauty->SetForceDecay(kSemiElectronic);
1240
1241         AliGenParam *beautyJ = new AliGenParam(10,
1242                                                new AliGenMUONlib(), 
1243                                                AliGenMUONlib::kBeauty,
1244                                                "central");
1245         beautyJ->SetPtRange(0, 100);
1246         beautyJ->SetYRange(-1.5, +1.5);
1247         beautyJ->SetForceDecay(kBJpsiDiElectron);
1248
1249         gener->AddGenerator(phi,"Phi",1);
1250         gener->AddGenerator(omega,"Omega",1);
1251         gener->AddGenerator(jpsi,"J/psi",1);
1252         gener->AddGenerator(ups,"Upsilon",1);
1253         gener->AddGenerator(charm,"Charm",1);
1254         gener->AddGenerator(beauty,"Beauty",1);
1255         gener->AddGenerator(beautyJ,"J/Psi from Beauty",1);
1256         gGener=gener;
1257       }
1258       break;
1259     case kPyJJ:
1260       {
1261         comment = comment.Append(" Jet-jet at 5.5 TeV");
1262         AliGenPythia *gener = new AliGenPythia(-1);
1263         gener->SetEnergyCMS(5500.);
1264         gener->SetProcess(kPyJets);
1265         Double_t ptHardMin=10.0, ptHardMax=-1.0;
1266         gener->SetPtHard(ptHardMin,ptHardMax);
1267         gener->SetYHard(-0.7,0.7);
1268         gener->SetJetEtaRange(-0.2,0.2);
1269         gener->SetEventListRange(0,1);
1270         gGener=gener;
1271       }
1272       break;
1273     case kPyGJ:
1274       {
1275         comment = comment.Append(" Gamma-jet at 5.5 TeV");
1276         AliGenPythia *gener = new AliGenPythia(-1);
1277         gener->SetEnergyCMS(5500.);
1278         gener->SetProcess(kPyDirectGamma);
1279         Double_t ptHardMin=10.0, ptHardMax=-1.0;
1280         gener->SetPtHard(ptHardMin,ptHardMax);
1281         gener->SetYHard(-1.0,1.0);
1282         gener->SetGammaEtaRange(-0.13,0.13);
1283         gener->SetGammaPhiRange(210.,330.);
1284         gener->SetEventListRange(0,1);
1285         gGener=gener;
1286       }
1287       break;
1288     case kMuonCocktailCent1:
1289       {
1290         comment = comment.Append(" Muon Cocktail Cent1");
1291         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1292         gener->SetPtRange(0.4,100.);       // Transverse momentum range   
1293         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1294         gener->SetYRange(-4.0,-2.4);
1295         gener->SetMuonPtCut(0.8);
1296         gener->SetMuonThetaCut(171.,178.);
1297         gener->SetMuonMultiplicity(2);
1298         gener->SetImpactParameterRange(12.,16.);  //Centrality class Cent1 for PDC04
1299         gGener=gener;
1300       }
1301       break;
1302     case kMuonCocktailPer1:
1303       {
1304         comment = comment.Append(" Muon Cocktail Per1");
1305         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1306         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1307         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1308         gener->SetYRange(-4.0,-2.4);
1309         gener->SetMuonPtCut(0.8);
1310         gener->SetMuonThetaCut(171.,178.);
1311         gener->SetMuonMultiplicity(2);
1312         gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1313         gGener=gener;
1314       }
1315       break;
1316     case kMuonCocktailPer4:
1317       {
1318         comment = comment.Append(" Muon Cocktail Per4");
1319         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1320         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1321         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1322         gener->SetYRange(-4.0,-2.4);
1323         gener->SetMuonPtCut(0.8);
1324         gener->SetMuonThetaCut(171.,178.);
1325         gener->SetMuonMultiplicity(2);
1326         gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1327         gGener=gener;
1328       }
1329       break;
1330     case kMuonCocktailCent1HighPt:
1331       {
1332         comment = comment.Append(" Muon Cocktail HighPt Cent1");
1333         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1334         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1335         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1336         gener->SetYRange(-4.0,-2.4);
1337         gener->SetMuonPtCut(2.5);
1338         gener->SetMuonThetaCut(171.,178.);
1339         gener->SetMuonMultiplicity(2);
1340         gener->SetImpactParameterRange(0.,5.);  //Centrality class Cent1 for PDC04
1341         gGener=gener;
1342       }
1343       break;
1344     case kMuonCocktailPer1HighPt :
1345       {
1346         comment = comment.Append(" Muon Cocktail HighPt Per1");
1347         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1348         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1349         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1350         gener->SetYRange(-4.0,-2.4);
1351         gener->SetMuonPtCut(2.5);
1352         gener->SetMuonThetaCut(171.,178.);
1353         gener->SetMuonMultiplicity(2);
1354         gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1355         gGener=gener;
1356       }
1357       break;
1358     case kMuonCocktailPer4HighPt:
1359       {
1360         comment = comment.Append(" Muon Cocktail HighPt Per4");
1361         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1362         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1363         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1364         gener->SetYRange(-4.0,-2.4);
1365         gener->SetMuonPtCut(2.5);
1366         gener->SetMuonThetaCut(171.,178.);
1367         gener->SetMuonMultiplicity(2);
1368         gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1369         gGener=gener;
1370       }
1371       break;
1372     case kMuonCocktailCent1Single:
1373       {
1374         comment = comment.Append(" Muon Cocktail Single Cent1");
1375         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1376         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1377         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1378         gener->SetYRange(-4.0,-2.4);
1379         gener->SetMuonPtCut(0.8);
1380         gener->SetMuonThetaCut(171.,178.);
1381         gener->SetMuonMultiplicity(1);
1382         gener->SetImpactParameterRange(0.,5.);  //Centrality class Cent1 for PDC04
1383         gGener=gener;
1384       }
1385       break;
1386     case kMuonCocktailPer1Single :
1387       {
1388         comment = comment.Append(" Muon Cocktail Single Per1");
1389         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1390         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1391         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1392         gener->SetYRange(-4.0,-2.4);
1393         gener->SetMuonPtCut(0.8);
1394         gener->SetMuonThetaCut(171.,178.);
1395         gener->SetMuonMultiplicity(1);
1396         gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1397         gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1398         gGener=gener;
1399       }
1400       break;
1401     case kMuonCocktailPer4Single:
1402       {
1403         comment = comment.Append(" Muon Cocktail Single Per4");
1404         AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1405         gener->SetPtRange(0.0,100.);       // Transverse momentum range   
1406         gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1407         gener->SetYRange(-4.0,-2.4);
1408         gener->SetMuonPtCut(0.8);
1409         gener->SetMuonThetaCut(171.,178.);
1410         gener->SetMuonMultiplicity(1);
1411         gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1412         gGener=gener;
1413       }
1414       break;
1415     case kFlow_2_2000:
1416     {
1417         comment = comment.Append(" Flow with dN/deta  = 2000, vn = 2%");
1418         gGener = GeVSimStandard(2000., 2.);
1419     }
1420     break;
1421     
1422     case kFlow_10_2000:
1423     {
1424         comment = comment.Append(" Flow with dN/deta  = 2000, vn = 10%");
1425         gGener = GeVSimStandard(2000., 10.);
1426     }
1427     break;
1428     
1429     case kFlow_6_2000:
1430     {
1431         comment = comment.Append(" Flow with dN/deta  = 2000, vn = 6%");
1432         gGener = GeVSimStandard(2000., 6.);
1433     }
1434     break;
1435     
1436     case kFlow_6_5000:
1437     {
1438         comment = comment.Append(" Flow with dN/deta  = 5000, vn = 6%");
1439         gGener = GeVSimStandard(5000., 6.);
1440     }
1441     break;
1442     case kHIJINGplus:
1443     {
1444         //
1445         // The cocktail
1446         AliGenCocktail *gener  = new AliGenCocktail();
1447
1448         //
1449         // Charm production by Pythia
1450         AliGenPythia * genpyc = new AliGenPythia(230);
1451         genpyc->SetProcess(kPyCharmPbPbMNR);
1452         genpyc->SetStrucFunc(kCTEQ4L);
1453         genpyc->SetPtHard(2.1,-1.0);
1454         genpyc->SetEnergyCMS(5500.);
1455         genpyc->SetNuclei(208,208);
1456         genpyc->SetYRange(-999,999);
1457         genpyc->SetForceDecay(kAll);
1458         genpyc->SetFeedDownHigherFamily(kFALSE);
1459         genpyc->SetCountMode(AliGenPythia::kCountParents);
1460         //
1461         // Beauty production by Pythia
1462         AliGenPythia * genpyb = new AliGenPythia(9);
1463         genpyb->SetProcess(kPyBeautyPbPbMNR);
1464         genpyb->SetStrucFunc(kCTEQ4L);
1465         genpyb->SetPtHard(2.75,-1.0);
1466         genpyb->SetEnergyCMS(5500.);
1467         genpyb->SetNuclei(208,208);
1468         genpyb->SetYRange(-999,999);
1469         genpyb->SetForceDecay(kAll);
1470         genpyb->SetFeedDownHigherFamily(kFALSE);
1471         genpyb->SetCountMode(AliGenPythia::kCountParents);
1472         //
1473         // Hyperons
1474         //
1475         AliGenSTRANGElib *lib = new AliGenSTRANGElib();
1476         Int_t particle;
1477         // Xi
1478         particle = AliGenSTRANGElib::kXiMinus;
1479         AliGenParam *genXi = new AliGenParam(16,particle,lib->GetPt(particle),lib->GetY(particle),lib->GetIp(particle));
1480         genXi->SetPtRange(0., 12.);
1481         genXi->SetYRange(-1.1, 1.1);
1482         genXi->SetForceDecay(kNoDecay); 
1483  
1484         //
1485         // Omega
1486         particle = AliGenSTRANGElib::kOmegaMinus;
1487         AliGenParam *genOmega = new AliGenParam(10,particle,lib->GetPt(particle),lib->GetY(particle),lib->GetIp(particle));     
1488         genOmega->SetPtRange(0, 12.);
1489         genOmega->SetYRange(-1.1, 1.1);
1490         genOmega->SetForceDecay(kNoDecay);
1491         
1492         //
1493         // Central Hijing 
1494         AliGenHijing *genHi = HijingStandard();
1495         genHi->SwitchOffHeavyQuarks(kTRUE);
1496         genHi->SetImpactParameterRange(0.,5.);
1497         //
1498         // Add everything to the cocktail and shake ...
1499         gener->AddGenerator(genHi,    "Hijing cent1", 1);
1500         gener->AddGenerator(genpyc,   "Extra charm",  1);
1501         gener->AddGenerator(genpyb,   "Extra beauty", 1);
1502         gener->AddGenerator(genXi,    "Xi"          , 1);
1503         gener->AddGenerator(genOmega, "Omega",        1);
1504         gGener = gener;
1505     }
1506     break;
1507     default: break;
1508     }
1509     
1510     return gGener;
1511 }
1512
1513 AliGenHijing* HijingStandard()
1514 {
1515     AliGenHijing *gener = new AliGenHijing(-1);
1516 // centre of mass energy 
1517     gener->SetEnergyCMS(5500.);
1518 // reference frame
1519     gener->SetReferenceFrame("CMS");
1520 // projectile
1521      gener->SetProjectile("A", 208, 82);
1522      gener->SetTarget    ("A", 208, 82);
1523 // tell hijing to keep the full parent child chain
1524      gener->KeepFullEvent();
1525 // enable jet quenching
1526      gener->SetJetQuenching(1);
1527 // enable shadowing
1528      gener->SetShadowing(1);
1529 // neutral pion and heavy particle decays switched off
1530      gener->SetDecaysOff(1);
1531 // Don't track spectators
1532      gener->SetSpectators(0);
1533 // kinematic selection
1534      gener->SetSelectAll(0);
1535      return gener;
1536 }
1537
1538 AliGenGeVSim* GeVSimStandard(Float_t mult, Float_t vn)
1539 {
1540     AliGenGeVSim* gener = new AliGenGeVSim(0);
1541 //
1542 // Mult is the number of charged particles in |eta| < 0.5
1543 // Vn is in (%)
1544 //
1545 // Sigma of the Gaussian dN/deta
1546     Float_t sigma_eta  = 2.75;
1547 //
1548 // Maximum eta
1549     Float_t etamax     = 7.00;
1550 //
1551 //
1552 // Scale from multiplicity in |eta| < 0.5 to |eta| < |etamax|   
1553     Float_t mm = mult * (TMath::Erf(etamax/sigma_eta/sqrt(2.)) / TMath::Erf(0.5/sigma_eta/sqrt(2.))); 
1554 //
1555 // Scale from charged to total multiplicity
1556 // 
1557     mm *= 1.587;
1558 //
1559 // Vn 
1560     vn /= 100.;          
1561 //
1562 // Define particles
1563 //
1564 //
1565 // 78% Pions (26% pi+, 26% pi-, 26% p0)              T = 250 MeV
1566     AliGeVSimParticle *pp =  new AliGeVSimParticle(kPiPlus,  1, 0.26 * mm, 0.25, sigma_eta) ;
1567     AliGeVSimParticle *pm =  new AliGeVSimParticle(kPiMinus, 1, 0.26 * mm, 0.25, sigma_eta) ;
1568     AliGeVSimParticle *p0 =  new AliGeVSimParticle(kPi0,     1, 0.26 * mm, 0.25, sigma_eta) ;
1569 //
1570 // 12% Kaons (3% K0short, 3% K0long, 3% K+, 3% K-)   T = 300 MeV
1571     AliGeVSimParticle *ks =  new AliGeVSimParticle(kK0Short, 1, 0.03 * mm, 0.30, sigma_eta) ;
1572     AliGeVSimParticle *kl =  new AliGeVSimParticle(kK0Long,  1, 0.03 * mm, 0.30, sigma_eta) ;
1573     AliGeVSimParticle *kp =  new AliGeVSimParticle(kKPlus,   1, 0.03 * mm, 0.30, sigma_eta) ;
1574     AliGeVSimParticle *km =  new AliGeVSimParticle(kKMinus,  1, 0.03 * mm, 0.30, sigma_eta) ;
1575 //
1576 // 10% Protons / Neutrons (5% Protons, 5% Neutrons)  T = 250 MeV
1577     AliGeVSimParticle *pr =  new AliGeVSimParticle(kProton,  1, 0.05 * mm, 0.25, sigma_eta) ;
1578     AliGeVSimParticle *ne =  new AliGeVSimParticle(kNeutron, 1, 0.05 * mm, 0.25, sigma_eta) ;
1579 //
1580 // Set Elliptic Flow properties         
1581     pp->SetEllipticParam(vn,1.0,0.) ;
1582     pm->SetEllipticParam(vn,1.0,0.) ;
1583     p0->SetEllipticParam(vn,1.0,0.) ;
1584     pr->SetEllipticParam(vn,1.0,0.) ;
1585     ne->SetEllipticParam(vn,1.0,0.) ;
1586     ks->SetEllipticParam(vn,1.0,0.) ;
1587     kl->SetEllipticParam(vn,1.0,0.) ;
1588     kp->SetEllipticParam(vn,1.0,0.) ;
1589     km->SetEllipticParam(vn,1.0,0.) ;
1590 //
1591 // Set Direct Flow properties   
1592     pp->SetDirectedParam(vn,1.0,0.) ;
1593     pm->SetDirectedParam(vn,1.0,0.) ;
1594     p0->SetDirectedParam(vn,1.0,0.) ;
1595     pr->SetDirectedParam(vn,1.0,0.) ;
1596     ne->SetDirectedParam(vn,1.0,0.) ;
1597     ks->SetDirectedParam(vn,1.0,0.) ;
1598     kl->SetDirectedParam(vn,1.0,0.) ;
1599     kp->SetDirectedParam(vn,1.0,0.) ;
1600     km->SetDirectedParam(vn,1.0,0.) ;
1601 //
1602 // Add particles to the list
1603     gener->AddParticleType(pp) ;
1604     gener->AddParticleType(pm) ;
1605     gener->AddParticleType(p0) ;
1606     gener->AddParticleType(pr) ;
1607     gener->AddParticleType(ne) ;
1608     gener->AddParticleType(ks) ;
1609     gener->AddParticleType(kl) ;
1610     gener->AddParticleType(kp) ;
1611     gener->AddParticleType(km) ;
1612 //      
1613 // Random Ev.Plane ----------------------------------
1614     TF1 *rpa = new TF1("gevsimPsiRndm","1", 0, 360);
1615 // --------------------------------------------------
1616     gener->SetPtRange(0., 9.) ; // Use a resonable range! (used for bin size in numerical integration)
1617     gener->SetPhiRange(0, 360);
1618     //
1619     // Set pseudorapidity range 
1620     Float_t thmin = EtaToTheta(+etamax);   
1621     Float_t thmax = EtaToTheta(-etamax);   
1622     gener->SetThetaRange(thmin,thmax);     
1623     return gener;
1624 }
1625
1626
1627
1628 void ProcessEnvironmentVars()
1629 {
1630     // Run type
1631     if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
1632       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
1633         if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
1634           srun = (PprRun_t)iRun;
1635           cout<<"Run type set to "<<pprRunName[iRun]<<endl;
1636         }
1637       }
1638     }
1639
1640     // Random Number seed
1641     if (gSystem->Getenv("CONFIG_SEED")) {
1642       sseed = atoi(gSystem->Getenv("CONFIG_SEED"));
1643     }
1644 }