]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/Config.C
Bug fix
[u/mrichter/AliRoot.git] / FMD / Config.C
1 //____________________________________________________________________
2 //
3 // $Id$
4 //
5 // One can use the configuration macro in compiled mode by
6 // root [0] gSystem->Load("libgeant321");
7 // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
8 //                   -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
9 // root [0] .x grun.C(1,"ConfigPPR.C++")
10 //
11 /** @file    Config.C
12     @author  Christian Holm Christensen <cholm@nbi.dk>
13     @date    Mon Mar 27 12:50:29 2006
14     @brief   Simulation configuration script
15 */
16
17 //____________________________________________________________________
18 // 
19 // Generator types 
20 //
21 enum EG_t {
22   test50,
23   kParam_8000,                  //
24   kParam_4000,                  //
25   kParam_2000,                  //
26   kParam_fmd,                   //
27   kHijing_cent1,                //
28   kHijing_cent2,                //
29   kHijing_per1,                 //
30   kHijing_per2,                 //
31   kHijing_per3,                 //
32   kHijing_per4,                 //
33   kHijing_per5,                 //
34   kHijing_jj25,                 //
35   kHijing_jj50,                 //
36   kHijing_jj75,                 //
37   kHijing_jj100,                //
38   kHijing_jj200,                //
39   kHijing_gj25,                 //
40   kHijing_gj50,                 //
41   kHijing_gj75,                 //
42   kHijing_gj100,                //
43   kHijing_gj200,                //
44   kHijing_pA,                   //
45   kPythia6,                     //
46   kPythia6Jets20_24,            //
47   kPythia6Jets24_29,            //
48   kPythia6Jets29_35,            //
49   kPythia6Jets35_42,            //
50   kPythia6Jets42_50,            //
51   kPythia6Jets50_60,            //
52   kPythia6Jets60_72,            //
53   kPythia6Jets72_86,            //
54   kPythia6Jets86_104,           //
55   kPythia6Jets104_125,          //
56   kPythia6Jets125_150,          //
57   kPythia6Jets150_180,          //
58   kD0PbPb5500,                  //
59   kCharmSemiElPbPb5500,         //
60   kBeautySemiElPbPb5500,        //
61   kCocktailTRD,                 //
62   kPyJJ,                        //
63   kPyGJ,                        //
64   kMuonCocktailCent1,           //
65   kMuonCocktailPer1,            //
66   kMuonCocktailPer4,            //
67   kMuonCocktailCent1HighPt,     //
68   kMuonCocktailPer1HighPt,      //
69   kMuonCocktailPer4HighPt,      //
70   kMuonCocktailCent1Single,     //
71   kMuonCocktailPer1Single,      //
72   kMuonCocktailPer4Single,
73   kFMD1Flat, 
74   kFMD2Flat, 
75   kFMD3Flat,
76   kFMDFlat,
77   kEgMax
78 };
79
80 //____________________________________________________________________
81 // 
82 // Generator types names
83 //
84 const char* egName[kEgMax] = {
85   "test50",
86   "kParam_8000",                //
87   "kParam_4000",                //
88   "kParam_2000",                //
89   "kParam_fmd",                 //
90   "kHijing_cent1",              //
91   "kHijing_cent2",              //
92   "kHijing_per1",               //
93   "kHijing_per2",               //
94   "kHijing_per3",               //
95   "kHijing_per4",               //
96   "kHijing_per5",               //
97   "kHijing_jj25",               //
98   "kHijing_jj50",               //
99   "kHijing_jj75",               //
100   "kHijing_jj100",              //
101   "kHijing_jj200",              //
102   "kHijing_gj25",               //
103   "kHijing_gj50",               //
104   "kHijing_gj75",               //
105   "kHijing_gj100",              //
106   "kHijing_gj200",              //
107   "kHijing_pA",                 //
108   "kPythia6",                   //
109   "kPythia6Jets20_24",          //
110   "kPythia6Jets24_29",          //
111   "kPythia6Jets29_35",          //
112   "kPythia6Jets35_42",          //
113   "kPythia6Jets42_50",          //
114   "kPythia6Jets50_60",          //
115   "kPythia6Jets60_72",          //
116   "kPythia6Jets72_86",          //
117   "kPythia6Jets86_104",         //
118   "kPythia6Jets104_125",        //
119   "kPythia6Jets125_150",        //
120   "kPythia6Jets150_180",        //
121   "kD0PbPb5500",                //
122   "kCharmSemiElPbPb5500",       //
123   "kBeautySemiElPbPb5500",      //
124   "kCocktailTRD",               //
125   "kPyJJ",                      //
126   "kPyGJ",                      //
127   "kMuonCocktailCent1",         //
128   "kMuonCocktailPer1",          //
129   "kMuonCocktailPer4",          //
130   "kMuonCocktailCent1HighPt",   //
131   "kMuonCocktailPer1HighPt",    //
132   "kMuonCocktailPer4HighPt",    //
133   "kMuonCocktailCent1Single",   //
134   "kMuonCocktailPer1Single",    //
135   "kMuonCocktailPer4Single",
136   "kFMD1Flat",
137   "kFMD2Flat",
138   "kFMD3Flat",
139   "kFMDFlat"
140 };
141
142 //____________________________________________________________________
143 enum Geo_t {
144   kHoles,                       //
145   kNoHoles                      //
146 };
147
148 //____________________________________________________________________
149 enum Rad_t {
150   kGluonRadiation,              //
151   kNoGluonRadiation             //
152 };
153
154 //____________________________________________________________________
155 enum MC_t {
156   kFLUKA, 
157   kGEANT3, 
158   kGEANT4, 
159   kGEANT3TGEO,
160 };
161
162 //____________________________________________________________________
163 // Functions
164 Float_t       EtaToTheta(Float_t eta);
165 EG_t          LookupEG(const Char_t* name);
166 AliGenerator* GeneratorFactory(EG_t eg, Rad_t rad);
167 AliGenHijing* HijingStandard();
168 void          ProcessEnvironmentVars(EG_t& eg, Int_t& seed);
169
170 //____________________________________________________________________
171 void 
172 Config()
173 {
174   //____________________________________________________________________
175   // This part for configuration    
176   EG_t  eg   = kPythia6;
177   Geo_t geo  = kNoHoles;
178   Rad_t rad  = kGluonRadiation;
179   AliMagF::BMap_t mag  = AliMagF::k5kG;
180   Int_t seed = 12345; //Set 0 to use the current time
181   MC_t  mc   = kGEANT3TGEO;
182   
183   //____________________________________________________________________
184   // Get settings from environment variables
185   ProcessEnvironmentVars(eg, seed);
186
187   //____________________________________________________________________
188   // Set Random Number seed
189   gRandom->SetSeed(seed);
190   cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
191
192
193   //__________________________________________________________________
194   switch (mc) {
195   case kFLUKA: 
196     // 
197     // libraries required by fluka21
198     // 
199     gSystem->Load("libGeom");
200     cout << "\t* Loading TFluka..." << endl;  
201     gSystem->Load("libTFluka");    
202     gSystem->MakeDirectory("peg");
203     // 
204     // FLUKA MC
205     //
206     cout << "\t* Instantiating TFluka..." << endl;
207     new TFluka("C++ Interface to Fluka", 0/*verbosity*/);
208     break;
209   case kGEANT3: 
210     {
211       //
212       // Libraries needed by GEANT 3.21 
213       //
214       gSystem->Load("$ALICE_ROOT/lib/tgt_$ALICE_TARGET/libpythia6.so");
215       gSystem->Load("libgeant321");
216       
217       // 
218       // GEANT 3.21 MC 
219       // 
220       TGeant3* gmc = new TGeant3("C++ Interface to Geant3");
221       gmc->SetSWIT(4, 1000);
222     }
223     break;
224   case kGEANT3TGEO:
225     {
226       //
227       // Libraries needed by GEANT 3.21 
228       //
229       gSystem->Load("$ALICE_ROOT/lib/tgt_$ALICE_TARGET/liblhapdf.so");
230       gSystem->Load("$ALICE_ROOT/lib/tgt_$ALICE_TARGET/libpythia6.so");
231       gSystem->Load("libEGPythia6.so"); //<- For non-debian (sigh!)
232       // gSystem->Load("EGPythia6.so");
233       gSystem->Load("libgeant321");
234     
235       // 
236       // GEANT 3.21 MC 
237       // 
238       TGeant3TGeo* gmc  = new TGeant3TGeo("C++ Interface to Geant3");
239       gmc->SetSWIT(4, 1000);
240       Printf("Making a TGeant3TGeo objet");
241     }
242     break;
243   default:
244     gAlice->Fatal("Config.C", "No MC type chosen");
245     return;
246   }
247
248   //__________________________________________________________________
249   AliRunLoader* rl = 0;
250
251   cout<<"Config.C: Creating Run Loader ..."<<endl;
252   rl = AliRunLoader::Open("galice.root",
253                           AliConfig::GetDefaultEventFolderName(),
254                           "recreate");
255   if (!rl) {
256     gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
257     return;
258   }
259   rl->SetCompressionLevel(2);
260   rl->SetNumberOfEventsPerFile(3);
261   gAlice->SetRunLoader(rl);
262
263   //__________________________________________________________________
264   // For FLUKA 
265   switch (mc) {
266   case kFLUKA: 
267     {
268       //
269       // Use kTRUE as argument to generate alice.pemf first
270       //
271       TString alice_pemf(gSystem->Which(".", "peg/mat17.pemf"));
272       if (!alice_pemf.IsNull()) 
273         ((TFluka*)gMC)->SetGeneratePemf(kFALSE);
274       else
275         ((TFluka*)gMC)->SetGeneratePemf(kTRUE);
276       TString flupro(gSystem->Getenv("FLUPRO"));
277       if (flupro.IsNull()) 
278         Fatal("Config.C", "Environment variable FLUPRO not set");
279 #if 0
280       char* files[] = { "brems_fin.bin", "cohff.bin", "elasct.bin", 
281                         "gxsect.bin", "nuclear.bin", "sigmapi.bin", 
282                         0 };
283       char* file = files[0];
284       while (file) {
285         TString which(gSystem->Which(".", file));
286         if (which.IsNull()) {
287           if (gSystem->Symlink(Form("%s/%s", flupro.Data(), file), file)!=0) 
288             Fatal("Config.C", "Couldn't link $(FLUPRO)/%s -> .", file);
289         }
290         file++;
291       }
292 #endif
293       TString neuxsc(gSystem->Which(".", "neuxsc.bin"));
294       if (neuxsc.IsNull()) 
295         gSystem->Symlink(Form("%s/neuxsc_72.bin", flupro.Data()), 
296                          "neuxsc.bin"); 
297       gSystem->CopyFile("$(FLUPRO)/random.dat", "old.seed", kTRUE);
298     }
299     break;
300   }
301
302   //__________________________________________________________________
303   //
304   // Set External decayer
305 #if 0
306   AliDecayer *decayer = new AliDecayerPythia();
307   switch (eg) {
308   case kD0PbPb5500:           decayer->SetForceDecay(kHadronicD);      break;
309   case kCharmSemiElPbPb5500:  decayer->SetForceDecay(kSemiElectronic); break;
310   case kBeautySemiElPbPb5500: decayer->SetForceDecay(kSemiElectronic); break;
311   default:                    decayer->SetForceDecay(kAll);            break;
312   }
313   decayer->Init();
314   gMC->SetExternalDecayer(decayer);
315 #endif
316
317   //__________________________________________________________________
318   // *********** STEERING parameters FOR ALICE SIMULATION ************
319   // - Specify event type to be tracked through the ALICE setup
320   // - All positions are in cm, angles in degrees, and P and E in GeV 
321   gMC->SetProcess("DCAY",1);
322   gMC->SetProcess("PAIR",1);
323   gMC->SetProcess("COMP",1);
324   gMC->SetProcess("PHOT",1);
325   gMC->SetProcess("PFIS",0);
326   gMC->SetProcess("DRAY",0);
327   gMC->SetProcess("ANNI",1);
328   gMC->SetProcess("BREM",1);
329   gMC->SetProcess("MUNU",1);
330   gMC->SetProcess("CKOV",1);
331   gMC->SetProcess("HADR",1);
332   gMC->SetProcess("LOSS",1); // 0:none 1,3:dray 2:nodray 4:nofluct (def:2)
333   gMC->SetProcess("MULS",1);
334   gMC->SetProcess("RAYL",1);
335
336   Float_t cut = 1.e-3;        // 1MeV cut by default
337   Float_t tofmax = 1.e10;
338
339   gMC->SetCut("CUTGAM", cut);
340   gMC->SetCut("CUTELE", cut);
341   gMC->SetCut("CUTNEU", cut);
342   gMC->SetCut("CUTHAD", cut);
343   gMC->SetCut("CUTMUO", cut);
344   gMC->SetCut("BCUTE",  cut); 
345   gMC->SetCut("BCUTM",  cut); 
346   gMC->SetCut("DCUTE",  cut); 
347   gMC->SetCut("DCUTM",  cut); 
348   gMC->SetCut("PPCUTM", cut);
349   gMC->SetCut("TOFMAX", tofmax); 
350
351   
352   //__________________________________________________________________
353   // Generator Configuration
354   AliGenerator* gener = GeneratorFactory(eg, rad);
355   gener->SetOrigin(0, 0, 0);    // vertex position
356   gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
357   gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
358   gener->SetVertexSmear(kPerEvent); 
359   gener->SetTrackingFlag(1);
360   gener->Init();
361     
362   //__________________________________________________________________
363   // Field (L3 0.4 T)
364   AliMagF* field = new AliMagF("Maps","Maps",-1., -1., mag);
365   TGeoGlobalMagField::Instance()->SetField(field);
366
367   rl->CdGAFile();
368
369   TFile* magF = TFile::Open("mag.root", "RECREATE");
370   field->Write("mag");
371   magF->Close();
372
373   //__________________________________________________________________
374   // 
375   // Used detectors 
376   // 
377   Bool_t useABSO  = kFALSE; 
378   Bool_t useACORDE= kFALSE; 
379   Bool_t useDIPO  = kFALSE; 
380   Bool_t useFMD   = kTRUE; 
381   Bool_t useFRAME = kFALSE; 
382   Bool_t useHALL  = kFALSE; 
383   Bool_t useITS   = kTRUE;
384   Bool_t useMAG   = kFALSE; 
385   Bool_t useMUON  = kFALSE; 
386   Bool_t usePHOS  = kFALSE; 
387   Bool_t usePIPE  = kTRUE; 
388   Bool_t usePMD   = kFALSE; 
389   Bool_t useHMPID = kFALSE; 
390   Bool_t useSHIL  = kFALSE; 
391   Bool_t useT0    = kTRUE; 
392   Bool_t useTOF   = kFALSE; 
393   Bool_t useTPC   = kFALSE;
394   Bool_t useTRD   = kFALSE; 
395   Bool_t useZDC   = kFALSE; 
396   Bool_t useEMCAL = kFALSE; 
397   Bool_t useVZERO = kTRUE;
398
399   cout << "\t* Creating the detectors ..." << endl;
400   // ================= Alice BODY parameters =========================
401   AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
402   
403   
404   if (useMAG)   AliMAG  *MAG  = new AliMAG("MAG", "Magnet");
405   if (useABSO)  AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
406   if (useDIPO)  AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
407   if (useHALL)  AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
408   if (useFRAME) {
409     AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
410     switch (geo) {
411     case kHoles: FRAME->SetHoles(1); break;
412     default:     FRAME->SetHoles(0); break;
413     }
414   }
415   if (useSHIL)   AliSHIL   *SHIL   = new AliSHILv3("SHIL", "Shielding v3");
416   if (usePIPE)   AliPIPE   *PIPE   = new AliPIPEv3("PIPE", "Beam Pipe");
417   if (useITS)    AliITS    *ITS    = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
418   if (useTPC)    AliTPC    *TPC    = new AliTPCv2("TPC", "Default");
419   if (useTOF)    AliTOF    *TOF    = new AliTOFv6T0("TOF", "normal TOF");
420   if (useHMPID)  AliHMPID  *HMPID  = new AliHMPIDv1("HMPID", "normal HMPID");
421   if (useZDC)    AliZDC    *ZDC    = new AliZDCv3("ZDC", "normal ZDC");
422   if (useTRD)    AliTRD    *TRD    = new AliTRDv1("TRD", "TRD slow simulator");
423   if (useFMD)    AliFMD    *FMD    = new AliFMDv1("FMD", "normal FMD");
424   if (useMUON)   AliMUON   *MUON   = new AliMUONv1("MUON", "default");
425   if (usePHOS)   AliPHOS   *PHOS   = new AliPHOSv1("PHOS", "IHEP");
426   if (usePMD)    AliPMD    *PMD    = new AliPMDv1("PMD", "normal PMD");
427   if (useT0)     AliT0     *T0     = new AliT0v1("T0", "T0 Detector");
428   if (useEMCAL)  AliEMCAL  *EMCAL  = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
429   if (useACORDE) AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
430   if (useVZERO)  AliVZERO  *VZERO  = new AliVZEROv7("VZERO", "normal VZERO");
431 }
432
433 //____________________________________________________________________
434 Float_t EtaToTheta(Float_t arg)
435 {
436   return (180./TMath::Pi())*2.*TMath::ATan(TMath::Exp(-arg));
437 }
438
439 //____________________________________________________________________
440 Int_t 
441 LookupEG(const Char_t* name) 
442 {
443   TString n(name);
444   for (Int_t i = 0; i < kEgMax; i++) {
445     if (n == egName[i]) return i;
446   }
447   return -1;
448 }
449
450 //____________________________________________________________________  
451 AliGenerator* 
452 GeneratorFactory(EG_t eg, Rad_t rad)  
453 {
454   Int_t isw = 3;
455   if (rad == kNoGluonRadiation) isw = 0;
456   
457   // Possibly load libAliPythia6
458   switch (eg) {
459   case kPythia6:
460   case kPythia6Jets20_24:
461   case kPythia6Jets24_29:
462   case kPythia6Jets29_35:
463   case kPythia6Jets35_42:
464   case kPythia6Jets42_50:
465   case kPythia6Jets50_60:
466   case kPythia6Jets60_72:
467   case kPythia6Jets72_86:
468   case kPythia6Jets86_104:
469   case kPythia6Jets104_125:
470   case kPythia6Jets125_150:
471   case kPythia6Jets150_180:
472     gSystem->Load("liblhapdf.so");
473     // gSystem->Load("/usr/lib/libpythia.so");
474     // gSystem->ListLibraries();
475     gSystem->Load("EGPythia6.so");
476     gSystem->Load("libAliPythia6");
477     break;
478   }
479   
480   AliGenerator * gGener = 0;
481   switch (eg) {
482   case test50:
483     {
484       AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
485       gener->SetMomentumRange(0, 999999.);
486       gener->SetPhiRange(0., 360.);
487       // Set pseudorapidity range from -8 to 8.
488       Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
489       Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
490       gener->SetThetaRange(thmin,thmax);
491       gGener=gener;
492     }
493     break;
494   case kParam_8000:
495     {
496       AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
497       gener->SetMomentumRange(0, 999999.);
498       gener->SetPhiRange(0., 360.);
499       // Set pseudorapidity range from -8 to 8.
500       Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
501       Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
502       gener->SetThetaRange(thmin,thmax);
503       gGener=gener;
504     }
505     break;
506   case kParam_4000:
507     {
508       AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
509       gener->SetMomentumRange(0, 999999.);
510       gener->SetPhiRange(0., 360.);
511       // Set pseudorapidity range from -8 to 8.
512       Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
513       Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
514       gener->SetThetaRange(thmin,thmax);
515       gGener=gener;
516     }
517     break;
518   case kParam_2000:
519     {
520       AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
521       gener->SetMomentumRange(0, 999999.);
522       gener->SetPhiRange(0., 360.);
523       // Set pseudorapidity range from -8 to 8.
524       Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
525       Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
526       gener->SetThetaRange(thmin,thmax);
527       gGener=gener;
528     }
529     break;
530   case kParam_fmd:
531     {
532       AliGenHIJINGpara *gener = new AliGenHIJINGpara(500);
533       gener->SetMomentumRange(0, 999999.);
534       gener->SetPhiRange(0., 360.);
535       // Set pseudorapidity range from -8 to 8.
536       Float_t thmin = EtaToTheta(6);   // theta min. <---> eta max
537       Float_t thmax = EtaToTheta(2);  // theta max. <---> eta min 
538       gener->SetThetaRange(thmin,thmax);
539       gGener=gener;
540     }
541     break;
542     //
543     //  Hijing Central
544     //
545   case kHijing_cent1:
546     {
547       AliGenHijing *gener = HijingStandard();
548       // impact parameter range
549       gener->SetImpactParameterRange(0., 5.);
550       gGener=gener;
551     }
552     break;
553   case kHijing_cent2:
554     {
555       AliGenHijing *gener = HijingStandard();
556       // impact parameter range
557       gener->SetImpactParameterRange(0., 2.);
558       gGener=gener;
559     }
560     break;
561     //
562     // Hijing Peripheral 
563     //
564   case kHijing_per1:
565     {
566       AliGenHijing *gener = HijingStandard();
567       // impact parameter range
568       gener->SetImpactParameterRange(5., 8.6);
569       gGener=gener;
570     }
571     break;
572   case kHijing_per2:
573     {
574       AliGenHijing *gener = HijingStandard();
575       // impact parameter range
576       gener->SetImpactParameterRange(8.6, 11.2);
577       gGener=gener;
578     }
579     break;
580   case kHijing_per3:
581     {
582       AliGenHijing *gener = HijingStandard();
583       // impact parameter range
584       gener->SetImpactParameterRange(11.2, 13.2);
585       gGener=gener;
586     }
587     break;
588   case kHijing_per4:
589     {
590       AliGenHijing *gener = HijingStandard();
591       // impact parameter range
592       gener->SetImpactParameterRange(13.2, 15.);
593       gGener=gener;
594     }
595     break;
596   case kHijing_per5:
597     {
598       AliGenHijing *gener = HijingStandard();
599       // impact parameter range
600       gener->SetImpactParameterRange(15., 100.);
601       gGener=gener;
602     }
603     break;
604     //
605     //  Jet-Jet
606     //
607   case kHijing_jj25:
608     {
609       AliGenHijing *gener = HijingStandard();
610       // impact parameter range
611       gener->SetImpactParameterRange(0., 5.);
612       // trigger
613       gener->SetTrigger(1);
614       gener->SetPtJet(25.);
615       gener->SetRadiation(isw);
616       gener->SetSimpleJets(!isw);
617       gener->SetJetEtaRange(-0.3,0.3);
618       gener->SetJetPhiRange(75., 165.);   
619       gGener=gener;
620     }
621     break;
622
623   case kHijing_jj50:
624     {
625       AliGenHijing *gener = HijingStandard();
626       // impact parameter range
627       gener->SetImpactParameterRange(0., 5.);
628       // trigger
629       gener->SetTrigger(1);
630       gener->SetPtJet(50.);
631       gener->SetRadiation(isw);
632       gener->SetSimpleJets(!isw);
633       gener->SetJetEtaRange(-0.3,0.3);
634       gener->SetJetPhiRange(75., 165.);   
635       gGener=gener;
636     }
637     break;
638
639   case kHijing_jj75:
640     {
641       AliGenHijing *gener = HijingStandard();
642       // impact parameter range
643       gener->SetImpactParameterRange(0., 5.);
644       // trigger
645       gener->SetTrigger(1);
646       gener->SetPtJet(75.);
647       gener->SetRadiation(isw);
648       gener->SetSimpleJets(!isw);
649       gener->SetJetEtaRange(-0.3,0.3);
650       gener->SetJetPhiRange(75., 165.);   
651       gGener=gener;
652     }
653     break;
654
655   case kHijing_jj100:
656     {
657       AliGenHijing *gener = HijingStandard();
658       // impact parameter range
659       gener->SetImpactParameterRange(0., 5.);
660       // trigger
661       gener->SetTrigger(1);
662       gener->SetPtJet(100.);
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   case kHijing_jj200:
672     {
673       AliGenHijing *gener = HijingStandard();
674       // impact parameter range
675       gener->SetImpactParameterRange(0., 5.);
676       // trigger
677       gener->SetTrigger(1);
678       gener->SetPtJet(200.);
679       gener->SetRadiation(isw);
680       gener->SetSimpleJets(!isw);
681       gener->SetJetEtaRange(-0.3,0.3);
682       gener->SetJetPhiRange(75., 165.);   
683       gGener=gener;
684     }
685     break;
686     //
687     // Gamma-Jet
688     //
689   case kHijing_gj25:
690     {
691       AliGenHijing *gener = HijingStandard();
692       // impact parameter range
693       gener->SetImpactParameterRange(0., 5.);
694       // trigger
695       gener->SetTrigger(2);
696       gener->SetPtJet(25.);
697       gener->SetRadiation(isw);
698       gener->SetSimpleJets(!isw);
699       gener->SetJetEtaRange(-0.12, 0.12);
700       gener->SetJetPhiRange(220., 320.);
701       gGener=gener;
702     }
703     break;
704
705   case kHijing_gj50:
706     {
707       AliGenHijing *gener = HijingStandard();
708       // impact parameter range
709       gener->SetImpactParameterRange(0., 5.);
710       // trigger
711       gener->SetTrigger(2);
712       gener->SetPtJet(50.);
713       gener->SetRadiation(isw);
714       gener->SetSimpleJets(!isw);
715       gener->SetJetEtaRange(-0.12, 0.12);
716       gener->SetJetPhiRange(220., 320.);
717       gGener=gener;
718     }
719     break;
720
721   case kHijing_gj75:
722     {
723       AliGenHijing *gener = HijingStandard();
724       // impact parameter range
725       gener->SetImpactParameterRange(0., 5.);
726       // trigger
727       gener->SetTrigger(2);
728       gener->SetPtJet(75.);
729       gener->SetRadiation(isw);
730       gener->SetSimpleJets(!isw);
731       gener->SetJetEtaRange(-0.12, 0.12);
732       gener->SetJetPhiRange(220., 320.);
733       gGener=gener;
734     }
735     break;
736
737   case kHijing_gj100:
738     {
739       AliGenHijing *gener = HijingStandard();
740       // impact parameter range
741       gener->SetImpactParameterRange(0., 5.);
742       // trigger
743       gener->SetTrigger(2);
744       gener->SetPtJet(100.);
745       gener->SetRadiation(isw);
746       gener->SetSimpleJets(!isw);
747       gener->SetJetEtaRange(-0.12, 0.12);
748       gener->SetJetPhiRange(220., 320.);
749       gGener=gener;
750     }
751     break;
752
753   case kHijing_gj200:
754     {
755       AliGenHijing *gener = HijingStandard();
756       // impact parameter range
757       gener->SetImpactParameterRange(0., 5.);
758       // trigger
759       gener->SetTrigger(2);
760       gener->SetPtJet(200.);
761       gener->SetRadiation(isw);
762       gener->SetSimpleJets(!isw);
763       gener->SetJetEtaRange(-0.12, 0.12);
764       gener->SetJetPhiRange(220., 320.);
765       gGener=gener;
766     }
767     break;
768   case kHijing_pA:
769     {
770
771       AliGenCocktail *gener  = new AliGenCocktail();
772
773       AliGenHijing   *hijing = new AliGenHijing(-1);
774       // centre of mass energy 
775       hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
776       // impact parameter range
777       hijing->SetImpactParameterRange(0., 15.);
778       // reference frame
779       hijing->SetReferenceFrame("CMS");
780       hijing->SetBoostLHC(1);
781       // projectile
782       hijing->SetProjectile("P", 1, 1);
783       hijing->SetTarget    ("A", 208, 82);
784       // tell hijing to keep the full parent child chain
785       hijing->KeepFullEvent();
786       // enable jet quenching
787       hijing->SetJetQuenching(0);
788       // enable shadowing
789       hijing->SetShadowing(1);
790       // Don't track spectators
791       hijing->SetSpectators(0);
792       // kinematic selection
793       hijing->SetSelectAll(0);
794       //
795       AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
796       AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
797       gray->SetSlowNucleonModel(model);
798       gray->SetDebug(1);
799       gener->AddGenerator(hijing,"Hijing pPb", 1);
800       gener->AddGenerator(gray,  "Gray Particles",1);
801       gGener=gener;
802     }
803     break;
804   case kPythia6:
805     {
806       AliGenPythia *gener = new AliGenPythia(-1); 
807       gener->SetMomentumRange(0,999999);
808       gener->SetThetaRange(0., 180.);
809       gener->SetYRange(-12,12);
810       gener->SetPtRange(0,1000);
811       gener->SetProcess(kPyMb);
812       gener->SetEnergyCMS(14000.);
813       gGener=gener;
814     }
815     break;
816   case kPythia6Jets20_24:
817     {
818       AliGenPythia * gener = new AliGenPythia(-1);
819       gener->SetEnergyCMS(5500.);//        Centre of mass energy
820       gener->SetProcess(kPyJets);//        Process type
821       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
822       gener->SetJetPhiRange(0., 360.);
823       gener->SetJetEtRange(10., 1000.);
824       gener->SetGluonRadiation(1,1);
825       //    gener->SetPtKick(0.);
826       //   Structure function
827       gener->SetStrucFunc(kCTEQ4L);
828       gener->SetPtHard(20., 24.);// Pt transfer of the hard scattering
829       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
830       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
831       gGener=gener;
832     }
833     break;
834   case kPythia6Jets24_29:
835     {
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(24., 29.);// 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       gGener=gener;
850     }
851     break;
852   case kPythia6Jets29_35:
853     {
854       AliGenPythia * gener = new AliGenPythia(-1);
855       gener->SetEnergyCMS(5500.);//        Centre of mass energy
856       gener->SetProcess(kPyJets);//        Process type
857       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
858       gener->SetJetPhiRange(0., 360.);
859       gener->SetJetEtRange(10., 1000.);
860       gener->SetGluonRadiation(1,1);
861       //    gener->SetPtKick(0.);
862       //   Structure function
863       gener->SetStrucFunc(kCTEQ4L);
864       gener->SetPtHard(29., 35.);// Pt transfer of the hard scattering
865       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
866       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
867       gGener=gener;
868     }
869     break;
870   case kPythia6Jets35_42:
871     {
872       AliGenPythia * gener = new AliGenPythia(-1);
873       gener->SetEnergyCMS(5500.);//        Centre of mass energy
874       gener->SetProcess(kPyJets);//        Process type
875       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
876       gener->SetJetPhiRange(0., 360.);
877       gener->SetJetEtRange(10., 1000.);
878       gener->SetGluonRadiation(1,1);
879       //    gener->SetPtKick(0.);
880       //   Structure function
881       gener->SetStrucFunc(kCTEQ4L);
882       gener->SetPtHard(35., 42.);// Pt transfer of the hard scattering
883       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
884       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
885       gGener=gener;
886     }
887     break;
888   case kPythia6Jets42_50:
889     {
890       AliGenPythia * gener = new AliGenPythia(-1);
891       gener->SetEnergyCMS(5500.);//        Centre of mass energy
892       gener->SetProcess(kPyJets);//        Process type
893       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
894       gener->SetJetPhiRange(0., 360.);
895       gener->SetJetEtRange(10., 1000.);
896       gener->SetGluonRadiation(1,1);
897       //    gener->SetPtKick(0.);
898       //   Structure function
899       gener->SetStrucFunc(kCTEQ4L);
900       gener->SetPtHard(42., 50.);// Pt transfer of the hard scattering
901       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
902       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
903       gGener=gener;
904     }
905     break;
906   case kPythia6Jets50_60:
907     {
908       AliGenPythia * gener = new AliGenPythia(-1);
909       gener->SetEnergyCMS(5500.);//        Centre of mass energy
910       gener->SetProcess(kPyJets);//        Process type
911       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
912       gener->SetJetPhiRange(0., 360.);
913       gener->SetJetEtRange(10., 1000.);
914       gener->SetGluonRadiation(1,1);
915       //    gener->SetPtKick(0.);
916       //   Structure function
917       gener->SetStrucFunc(kCTEQ4L);
918       gener->SetPtHard(50., 60.);// Pt transfer of the hard scattering
919       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
920       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
921       gGener=gener;
922     }
923     break;
924   case kPythia6Jets60_72:
925     {
926       AliGenPythia * gener = new AliGenPythia(-1);
927       gener->SetEnergyCMS(5500.);//        Centre of mass energy
928       gener->SetProcess(kPyJets);//        Process type
929       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
930       gener->SetJetPhiRange(0., 360.);
931       gener->SetJetEtRange(10., 1000.);
932       gener->SetGluonRadiation(1,1);
933       //    gener->SetPtKick(0.);
934       //   Structure function
935       gener->SetStrucFunc(kCTEQ4L);
936       gener->SetPtHard(60., 72.);// Pt transfer of the hard scattering
937       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
938       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
939       gGener=gener;
940     }
941     break;
942   case kPythia6Jets72_86:
943     {
944       AliGenPythia * gener = new AliGenPythia(-1);
945       gener->SetEnergyCMS(5500.);//        Centre of mass energy
946       gener->SetProcess(kPyJets);//        Process type
947       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
948       gener->SetJetPhiRange(0., 360.);
949       gener->SetJetEtRange(10., 1000.);
950       gener->SetGluonRadiation(1,1);
951       //    gener->SetPtKick(0.);
952       //   Structure function
953       gener->SetStrucFunc(kCTEQ4L);
954       gener->SetPtHard(72., 86.);// Pt transfer of the hard scattering
955       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
956       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
957       gGener=gener;
958     }
959     break;
960   case kPythia6Jets86_104:
961     {
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       AliGenPythia * gener = new AliGenPythia(-1);
981       gener->SetEnergyCMS(5500.);//        Centre of mass energy
982       gener->SetProcess(kPyJets);//        Process type
983       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
984       gener->SetJetPhiRange(0., 360.);
985       gener->SetJetEtRange(10., 1000.);
986       gener->SetGluonRadiation(1,1);
987       //    gener->SetPtKick(0.);
988       //   Structure function
989       gener->SetStrucFunc(kCTEQ4L);
990       gener->SetPtHard(104., 125.);// Pt transfer of the hard scattering
991       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
992       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
993       gGener=gener;
994     }
995     break;
996   case kPythia6Jets125_150:
997     {
998       AliGenPythia * gener = new AliGenPythia(-1);
999       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1000       gener->SetProcess(kPyJets);//        Process type
1001       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1002       gener->SetJetPhiRange(0., 360.);
1003       gener->SetJetEtRange(10., 1000.);
1004       gener->SetGluonRadiation(1,1);
1005       //    gener->SetPtKick(0.);
1006       //   Structure function
1007       gener->SetStrucFunc(kCTEQ4L);
1008       gener->SetPtHard(125., 150.);// Pt transfer of the hard scattering
1009       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1010       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1011       gGener=gener;
1012     }
1013     break;
1014   case kPythia6Jets150_180:
1015     {
1016       AliGenPythia * gener = new AliGenPythia(-1);
1017       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1018       gener->SetProcess(kPyJets);//        Process type
1019       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1020       gener->SetJetPhiRange(0., 360.);
1021       gener->SetJetEtRange(10., 1000.);
1022       gener->SetGluonRadiation(1,1);
1023       //    gener->SetPtKick(0.);
1024       //   Structure function
1025       gener->SetStrucFunc(kCTEQ4L);
1026       gener->SetPtHard(150., 180.);// Pt transfer of the hard scattering
1027       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1028       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1029       gGener=gener;
1030     }
1031     break;
1032   case kD0PbPb5500:
1033     {
1034       AliGenPythia * gener = new AliGenPythia(10);
1035       gener->SetProcess(kPyD0PbPbMNR);
1036       gener->SetStrucFunc(kCTEQ4L);
1037       gener->SetPtHard(2.1,-1.0);
1038       gener->SetEnergyCMS(5500.);
1039       gener->SetNuclei(208,208);
1040       gener->SetForceDecay(kHadronicD);
1041       gener->SetYRange(-2,2);
1042       gener->SetFeedDownHigherFamily(kFALSE);
1043       gener->SetStackFillOpt(AliGenPythia::kParentSelection);
1044       gener->SetCountMode(AliGenPythia::kCountParents);
1045       gGener=gener;
1046     }
1047     break;
1048   case kCharmSemiElPbPb5500:
1049     {
1050       AliGenPythia * gener = new AliGenPythia(10);
1051       gener->SetProcess(kPyCharmPbPbMNR);
1052       gener->SetStrucFunc(kCTEQ4L);
1053       gener->SetPtHard(2.1,-1.0);
1054       gener->SetEnergyCMS(5500.);
1055       gener->SetNuclei(208,208);
1056       gener->SetForceDecay(kSemiElectronic);
1057       gener->SetYRange(-2,2);
1058       gener->SetFeedDownHigherFamily(kFALSE);
1059       gener->SetCountMode(AliGenPythia::kCountParents);
1060       gGener=gener;
1061     }
1062     break;
1063   case kBeautySemiElPbPb5500:
1064     {
1065       AliGenPythia *gener = new AliGenPythia(10);
1066       gener->SetProcess(kPyBeautyPbPbMNR);
1067       gener->SetStrucFunc(kCTEQ4L);
1068       gener->SetPtHard(2.75,-1.0);
1069       gener->SetEnergyCMS(5500.);
1070       gener->SetNuclei(208,208);
1071       gener->SetForceDecay(kSemiElectronic);
1072       gener->SetYRange(-2,2);
1073       gener->SetFeedDownHigherFamily(kFALSE);
1074       gener->SetCountMode(AliGenPythia::kCountParents);
1075       gGener=gener;
1076     }
1077     break;
1078   case kCocktailTRD:
1079     {
1080       AliGenCocktail *gener  = new AliGenCocktail();
1081
1082       AliGenParam *jpsi = new AliGenParam(10,
1083                                           new AliGenMUONlib(),
1084                                           AliGenMUONlib::kJpsiFamily,
1085                                           "Vogt PbPb");
1086
1087       jpsi->SetPtRange(0, 100);
1088       jpsi->SetYRange(-1., +1.);
1089       jpsi->SetForceDecay(kDiElectron);
1090
1091       AliGenParam *ups = new AliGenParam(10,
1092                                          new AliGenMUONlib(),
1093                                          AliGenMUONlib::kUpsilonFamily,
1094                                          "Vogt PbPb");
1095       ups->SetPtRange(0, 100);
1096       ups->SetYRange(-1., +1.);
1097       ups->SetForceDecay(kDiElectron);
1098         
1099       AliGenParam *charm = new AliGenParam(10,
1100                                            new AliGenMUONlib(), 
1101                                            AliGenMUONlib::kCharm,
1102                                            "central");
1103       charm->SetPtRange(0, 100);
1104       charm->SetYRange(-1.5, +1.5);
1105       charm->SetForceDecay(kSemiElectronic);
1106         
1107         
1108       AliGenParam *beauty = new AliGenParam(10,
1109                                             new AliGenMUONlib(), 
1110                                             AliGenMUONlib::kBeauty,
1111                                             "central");
1112       beauty->SetPtRange(0, 100);
1113       beauty->SetYRange(-1.5, +1.5);
1114       beauty->SetForceDecay(kSemiElectronic);
1115
1116       gener->AddGenerator(jpsi,"J/psi",1);
1117       gener->AddGenerator(ups,"Upsilon",1);
1118       gener->AddGenerator(charm,"Charm",1);
1119       gener->AddGenerator(beauty,"Beauty",1);
1120       gGener=gener;
1121     }
1122     break;
1123   case kPyJJ:
1124     {
1125       AliGenPythia *gener = new AliGenPythia(-1);
1126       gener->SetEnergyCMS(5500.);
1127       gener->SetProcess(kPyJets);
1128       Double_t ptHardMin=10.0, ptHardMax=-1.0;
1129       gener->SetPtHard(ptHardMin,ptHardMax);
1130       gener->SetYHard(-0.7,0.7);
1131       gener->SetJetEtaRange(-0.2,0.2);
1132       gener->SetEventListRange(0,1);
1133       gGener=gener;
1134     }
1135     break;
1136   case kPyGJ:
1137     {
1138       AliGenPythia *gener = new AliGenPythia(-1);
1139       gener->SetEnergyCMS(5500.);
1140       gener->SetProcess(kPyDirectGamma);
1141       Double_t ptHardMin=10.0, ptHardMax=-1.0;
1142       gener->SetPtHard(ptHardMin,ptHardMax);
1143       gener->SetYHard(-1.0,1.0);
1144       gener->SetGammaEtaRange(-0.13,0.13);
1145       gener->SetGammaPhiRange(210.,330.);
1146       gener->SetEventListRange(0,1);
1147       gGener=gener;
1148     }
1149     break;
1150   case kMuonCocktailCent1:
1151     {
1152       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1153       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1154       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1155       gener->SetYRange(-4.0,-2.4);
1156       gener->SetMuonPtCut(0.8);
1157       gener->SetMuonThetaCut(171.,178.);
1158       gener->SetMuonMultiplicity(2);
1159       gener->SetNumberOfCollisions(1626.);  //Centrality class Cent1 for PDC04
1160       gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1161       gGener=gener;
1162     }
1163     break;
1164   case kMuonCocktailPer1:
1165     {
1166       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1167       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1168       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1169       gener->SetYRange(-4.0,-2.4);
1170       gener->SetMuonPtCut(0.8);
1171       gener->SetMuonThetaCut(171.,178.);
1172       gener->SetMuonMultiplicity(2);
1173       gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1174       gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1175       gGener=gener;
1176     }
1177     break;
1178   case kMuonCocktailPer4:
1179     {
1180       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1181       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1182       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1183       gener->SetYRange(-4.0,-2.4);
1184       gener->SetMuonPtCut(0.8);
1185       gener->SetMuonThetaCut(171.,178.);
1186       gener->SetMuonMultiplicity(2);
1187       gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1188       gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1189       gGener=gener;
1190     }
1191     break;
1192   case kMuonCocktailCent1HighPt:
1193     {
1194       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1195       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1196       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1197       gener->SetYRange(-4.0,-2.4);
1198       gener->SetMuonPtCut(2.5);
1199       gener->SetMuonThetaCut(171.,178.);
1200       gener->SetMuonMultiplicity(2);
1201       gener->SetNumberOfCollisions(1626.);  //Centrality class Cent1 for PDC04
1202       gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1203       gGener=gener;
1204     }
1205     break;
1206   case kMuonCocktailPer1HighPt :
1207     {
1208       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1209       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1210       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1211       gener->SetYRange(-4.0,-2.4);
1212       gener->SetMuonPtCut(2.5);
1213       gener->SetMuonThetaCut(171.,178.);
1214       gener->SetMuonMultiplicity(2);
1215       gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1216       gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1217       gGener=gener;
1218     }
1219     break;
1220   case kMuonCocktailPer4HighPt:
1221     {
1222       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1223       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1224       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1225       gener->SetYRange(-4.0,-2.4);
1226       gener->SetMuonPtCut(2.5);
1227       gener->SetMuonThetaCut(171.,178.);
1228       gener->SetMuonMultiplicity(2);
1229       gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1230       gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1231       gGener=gener;
1232     }
1233     break;
1234   case kMuonCocktailCent1Single:
1235     {
1236       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1237       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1238       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1239       gener->SetYRange(-4.0,-2.4);
1240       gener->SetMuonPtCut(0.8);
1241       gener->SetMuonThetaCut(171.,178.);
1242       gener->SetMuonMultiplicity(1);
1243       gener->SetNumberOfCollisions(1626.);  //Centrality class Cent1 for PDC04
1244       gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1245       gGener=gener;
1246     }
1247     break;
1248   case kMuonCocktailPer1Single :
1249     {
1250       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1251       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1252       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1253       gener->SetYRange(-4.0,-2.4);
1254       gener->SetMuonPtCut(0.8);
1255       gener->SetMuonThetaCut(171.,178.);
1256       gener->SetMuonMultiplicity(1);
1257       gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1258       gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1259       gGener=gener;
1260     }
1261     break;
1262   case kMuonCocktailPer4Single:
1263     {
1264       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1265       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1266       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1267       gener->SetYRange(-4.0,-2.4);
1268       gener->SetMuonPtCut(0.8);
1269       gener->SetMuonThetaCut(171.,178.);
1270       gener->SetMuonMultiplicity(1);
1271       gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1272       gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1273       gGener=gener;
1274     }
1275     break;
1276   case kFMD1Flat: 
1277     {
1278       AliGenBox* gener = new AliGenBox(2000);
1279       gener->SetPart(kPiPlus);
1280       gener->SetMomentumRange(3,4);
1281       gener->SetPhiRange(0, 360);
1282       gener->SetThetaRange(0.77, 3.08);
1283       gGener = gener;
1284     }
1285     break;
1286   case kFMD2Flat: 
1287     {
1288       AliGenBox* gener = new AliGenBox(100);
1289       gener->SetPart(kPiPlus);
1290       gener->SetMomentumRange(3,4);
1291       gener->SetPhiRange(0, 360);
1292       gener->SetThetaRange(2.95, 20.42);
1293       gGener = gener;
1294     }
1295     break;
1296   case kFMD3Flat: 
1297     {
1298       AliGenBox* gener = new AliGenBox(2000);
1299       gener->SetPart(kPiPlus);
1300       gener->SetMomentumRange(3,4);
1301       gener->SetPhiRange(0, 360);
1302       gener->SetThetaRange(155.97, 176.73);
1303       gGener = gener;
1304     }
1305     break;
1306   case kFMDFlat:
1307     {
1308       AliGenCocktail* gener = new AliGenCocktail();
1309       gener->SetMomentumRange(3,4);
1310       gener->SetPhiRange(0, 360);
1311       AliGenBox* gener3 = new AliGenBox(2000);
1312       gener3->SetThetaRange(155.97, 176.73);
1313       gener3->SetPart(kPiPlus);
1314       gener->AddGenerator(gener3, "FMD3", .33);
1315       AliGenBox* gener2 = new AliGenBox(2000);
1316       gener2->SetThetaRange(2.95, 20.42);
1317       gener2->SetPart(kPiPlus);
1318       gener->AddGenerator(gener2, "FMD2", .33);
1319       AliGenBox* gener1 = new AliGenBox(2000);
1320       gener1->SetThetaRange(0.77, 3.08);
1321       gener1->SetPart(kPiPlus);
1322       gener->AddGenerator(gener1, "FMD1", .34);
1323       gGener = gener;
1324     }
1325     break;
1326     
1327   default: break;
1328   }
1329   return gGener;
1330 }
1331
1332 //____________________________________________________________________
1333 AliGenHijing* 
1334 HijingStandard()
1335 {
1336   AliGenHijing *gener = new AliGenHijing(-1);
1337   // centre of mass energy 
1338   gener->SetEnergyCMS(5500.);
1339   // reference frame
1340   gener->SetReferenceFrame("CMS");
1341   // projectile
1342   gener->SetProjectile("A", 208, 82);
1343   gener->SetTarget    ("A", 208, 82);
1344   // tell hijing to keep the full parent child chain
1345   gener->KeepFullEvent();
1346   // enable jet quenching
1347   gener->SetJetQuenching(1);
1348   // enable shadowing
1349   gener->SetShadowing(1);
1350   // neutral pion and heavy particle decays switched off
1351   gener->SetDecaysOff(1);
1352   // Don't track spectators
1353   gener->SetSpectators(0);
1354   // kinematic selection
1355   gener->SetSelectAll(0);
1356   return gener;
1357 }
1358
1359
1360 //____________________________________________________________________
1361 void 
1362 ProcessEnvironmentVars(EG_t& eg, Int_t& seed)
1363 {
1364   // Run type
1365   if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
1366     Int_t eg1 = LookupEG(gSystem->Getenv("CONFIG_RUN_TYPE"));
1367     if  (eg1 >= 0) eg = EG_t(eg1);
1368   }
1369   // Random Number seed
1370   if (gSystem->Getenv("CONFIG_SEED")) {
1371     seed = atoi(gSystem->Getenv("CONFIG_SEED"));
1372   }
1373 }
1374
1375 //____________________________________________________________________
1376 //
1377 // EOF
1378 //