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