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