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