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