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