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