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