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