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