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