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