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