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