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