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