]> git.uio.no Git - u/mrichter/AliRoot.git/blob - FMD/Config.C
Simple scripts to run various steps.
[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/AliACORDEv1.h"
62 #include "VZERO/AliVZEROv2.h"
63 #endif
64
65 //____________________________________________________________________
66 // 
67 // Generator types 
68 //
69 enum EG_t {
70   test50,
71   kParam_8000,                  //
72   kParam_4000,                  //
73   kParam_2000,                  //
74   kParam_fmd,                   //
75   kHijing_cent1,                //
76   kHijing_cent2,                //
77   kHijing_per1,                 //
78   kHijing_per2,                 //
79   kHijing_per3,                 //
80   kHijing_per4,                 //
81   kHijing_per5,                 //
82   kHijing_jj25,                 //
83   kHijing_jj50,                 //
84   kHijing_jj75,                 //
85   kHijing_jj100,                //
86   kHijing_jj200,                //
87   kHijing_gj25,                 //
88   kHijing_gj50,                 //
89   kHijing_gj75,                 //
90   kHijing_gj100,                //
91   kHijing_gj200,                //
92   kHijing_pA,                   //
93   kPythia6,                     //
94   kPythia6Jets20_24,            //
95   kPythia6Jets24_29,            //
96   kPythia6Jets29_35,            //
97   kPythia6Jets35_42,            //
98   kPythia6Jets42_50,            //
99   kPythia6Jets50_60,            //
100   kPythia6Jets60_72,            //
101   kPythia6Jets72_86,            //
102   kPythia6Jets86_104,           //
103   kPythia6Jets104_125,          //
104   kPythia6Jets125_150,          //
105   kPythia6Jets150_180,          //
106   kD0PbPb5500,                  //
107   kCharmSemiElPbPb5500,         //
108   kBeautySemiElPbPb5500,        //
109   kCocktailTRD,                 //
110   kPyJJ,                        //
111   kPyGJ,                        //
112   kMuonCocktailCent1,           //
113   kMuonCocktailPer1,            //
114   kMuonCocktailPer4,            //
115   kMuonCocktailCent1HighPt,     //
116   kMuonCocktailPer1HighPt,      //
117   kMuonCocktailPer4HighPt,      //
118   kMuonCocktailCent1Single,     //
119   kMuonCocktailPer1Single,      //
120   kMuonCocktailPer4Single,
121   kFMD1Flat, 
122   kFMD2Flat, 
123   kFMD3Flat,
124   kFMDFlat,
125   kEgMax
126 };
127
128 //____________________________________________________________________
129 // 
130 // Generator types names
131 //
132 const char* egName[kEgMax] = {
133   "test50",
134   "kParam_8000",                //
135   "kParam_4000",                //
136   "kParam_2000",                //
137   "kParam_fmd",                 //
138   "kHijing_cent1",              //
139   "kHijing_cent2",              //
140   "kHijing_per1",               //
141   "kHijing_per2",               //
142   "kHijing_per3",               //
143   "kHijing_per4",               //
144   "kHijing_per5",               //
145   "kHijing_jj25",               //
146   "kHijing_jj50",               //
147   "kHijing_jj75",               //
148   "kHijing_jj100",              //
149   "kHijing_jj200",              //
150   "kHijing_gj25",               //
151   "kHijing_gj50",               //
152   "kHijing_gj75",               //
153   "kHijing_gj100",              //
154   "kHijing_gj200",              //
155   "kHijing_pA",                 //
156   "kPythia6",                   //
157   "kPythia6Jets20_24",          //
158   "kPythia6Jets24_29",          //
159   "kPythia6Jets29_35",          //
160   "kPythia6Jets35_42",          //
161   "kPythia6Jets42_50",          //
162   "kPythia6Jets50_60",          //
163   "kPythia6Jets60_72",          //
164   "kPythia6Jets72_86",          //
165   "kPythia6Jets86_104",         //
166   "kPythia6Jets104_125",        //
167   "kPythia6Jets125_150",        //
168   "kPythia6Jets150_180",        //
169   "kD0PbPb5500",                //
170   "kCharmSemiElPbPb5500",       //
171   "kBeautySemiElPbPb5500",      //
172   "kCocktailTRD",               //
173   "kPyJJ",                      //
174   "kPyGJ",                      //
175   "kMuonCocktailCent1",         //
176   "kMuonCocktailPer1",          //
177   "kMuonCocktailPer4",          //
178   "kMuonCocktailCent1HighPt",   //
179   "kMuonCocktailPer1HighPt",    //
180   "kMuonCocktailPer4HighPt",    //
181   "kMuonCocktailCent1Single",   //
182   "kMuonCocktailPer1Single",    //
183   "kMuonCocktailPer4Single",
184   "kFMD1Flat",
185   "kFMD2Flat",
186   "kFMD3Flat",
187   "kFMDFlat"
188 };
189
190 //____________________________________________________________________
191 enum Geo_t {
192   kHoles,                       //
193   kNoHoles                      //
194 };
195
196 //____________________________________________________________________
197 enum Rad_t {
198   kGluonRadiation,              //
199   kNoGluonRadiation             //
200 };
201
202 //____________________________________________________________________
203 enum Mag_t {
204   k2kG,                         //
205   k4kG,                         //
206   k5kG                          //
207 };
208
209 //____________________________________________________________________
210 enum MC_t {
211   kFLUKA, 
212   kGEANT3, 
213   kGEANT4, 
214   kGEANT3TGEO,
215 };
216
217 //____________________________________________________________________
218 // Functions
219 Float_t       EtaToTheta(Float_t eta);
220 Eg_t          LookupEG(const Char_t* name);
221 AliGenerator* GeneratorFactory(EG_t eg, Rad_t rad, TString& comment);
222 AliGenHijing* HijingStandard();
223 void          ProcessEnvironmentVars(EG_t& eg, Int_t& seed);
224
225 //____________________________________________________________________
226 void 
227 Config()
228 {
229   //____________________________________________________________________
230   // This part for configuration    
231   // EG_t  eg   = test50;
232   // EG_t  eg   = kParam_fmd;
233   EG_t  eg   = kFMDFlat; // kParam_2000; // kPythia;
234   // EG_t  eg   = kFMDFlat;
235   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   TFile* magF = TFile::Open("mag.root", "RECREATE");
451   field->Write("mag");
452   magF->Close();
453
454   //__________________________________________________________________
455   // 
456   // Used detectors 
457   // 
458   Bool_t useABSO  = kFALSE; 
459   Bool_t useACORDE= kFALSE; 
460   Bool_t useDIPO  = kFALSE; 
461   Bool_t useFMD   = kTRUE; 
462   Bool_t useFRAME = kFALSE; 
463   Bool_t useHALL  = kFALSE; 
464   Bool_t useITS   = kFALSE;
465   Bool_t useMAG   = kFALSE; 
466   Bool_t useMUON  = kFALSE; 
467   Bool_t usePHOS  = kFALSE; 
468   Bool_t usePIPE  = kFALSE; 
469   Bool_t usePMD   = kFALSE; 
470   Bool_t useHMPID = kFALSE; 
471   Bool_t useSHIL  = kFALSE; 
472   Bool_t useT0    = kFALSE; 
473   Bool_t useTOF   = kFALSE; 
474   Bool_t useTPC   = kFALSE;
475   Bool_t useTRD   = kFALSE; 
476   Bool_t useZDC   = kFALSE; 
477   Bool_t useEMCAL = kFALSE; 
478   Bool_t useVZERO = kFALSE;
479
480   cout << "\t* Creating the detectors ..." << endl;
481   // ================= Alice BODY parameters =========================
482   AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
483   
484   
485   if (useMAG) {
486     // =================== MAG parameters ============================
487     // Start with Magnet since detector layouts may be depending on
488     // the selected Magnet dimensions 
489     AliMAG *MAG = new AliMAG("MAG", "Magnet");
490   }
491
492   if (useABSO) {
493     // =================== ABSO parameters ===========================
494     AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
495   }
496
497   if (useDIPO) {
498     // =================== DIPO parameters ===========================
499     AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
500   }
501
502   if (useHALL) {
503     // =================== HALL parameters ===========================
504     AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
505   }
506
507
508   if (useFRAME) {
509     // ================== FRAME parameters ===========================
510     AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
511     switch (geo) {
512     case kHoles: FRAME->SetHoles(1); break;
513     default:     FRAME->SetHoles(0); break;
514     }
515   }
516
517   if (useSHIL) {
518     // ================== SHIL parameters ============================
519     AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
520   }
521
522
523   if (usePIPE) {
524     // ================== PIPE parameters ============================
525     AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
526   }
527   
528   if (useITS) {
529      AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
530 #if 0
531     // =================== ITS parameters ============================
532     //
533     // As the innermost detector in ALICE, the Inner Tracking System
534     // "impacts" on almost all other detectors. This involves the fact
535     // that the ITS geometry still has several options to be followed
536     // in parallel in order to determine the best set-up which
537     // minimizes the induced background. All the geometries available
538     // to date are described in the following. Read carefully the
539     // comments and use the default version (the only one uncommented)
540     // unless you are making comparisons and you know what you are
541     // doing. In this case just uncomment the ITS geometry you want to
542     // use and run Aliroot.
543     //
544     // Detailed geometries:
545     //
546     //
547     // AliITS *ITS = 
548     //   new AliITSv5symm("ITS", "Updated ITS TDR detailed version "
549     //                    "with symmetric services");
550     // AliITS *ITS  = 
551     //   new AliITSv5asymm("ITS","Updates ITS TDR detailed version "
552     //                     "with asymmetric services");
553     //
554     AliITSvPPRasymmFMD *ITS  = 
555       new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version "
556                              "with asymmetric services");
557      // don't touch this parameter if you're not an ITS developer
558     ITS->SetMinorVersion(2); 
559     // don't touch this parameter if you're not an ITS developer
560     ITS->SetReadDet(kTRUE);
561     // don't touch this parameter if you're not an ITS developer
562     // ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  
563     // detector thickness on layer 1 must be in the range [100,300]
564     ITS->SetThicknessDet1(200.);   
565     // detector thickness on layer 2 must be in the range [100,300]
566     ITS->SetThicknessDet2(200.);   
567     // chip thickness on layer 1 must be in the range [150,300]
568     ITS->SetThicknessChip1(200.);  
569     // chip thickness on layer 2 must be in the range [150,300]
570     ITS->SetThicknessChip2(200.);
571     // 1 --> rails in ; 0 --> rails out
572     ITS->SetRails(0);          
573     // 1 --> water ; 0 --> freon
574     ITS->SetCoolingFluid(1);   
575
576     // Coarse geometries (warning: no hits are produced with these
577     // coarse geometries and they unuseful for reconstruction !):
578     //
579     //
580     // AliITSvPPRcoarseasymm *ITS  = 
581     //   new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version "
582     //                             "with asymmetric services");
583     // 1 --> rails in ; 0 --> rails out
584     // ITS->SetRails(0);
585     // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
586     // ITS->SetSupportMaterial(0);      
587     //
588     // AliITS *ITS  = 
589     //  new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version "
590     //                           "with symmetric services");
591     // 1 --> rails in ; 0 --> rails out
592     // ITS->SetRails(0);                
593     // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
594     // ITS->SetSupportMaterial(0);      
595     //
596     // Geant3 <-> EUCLID conversion
597     // ============================
598     //
599     // SetEUCLID is a flag to output (=1) or not to output (=0) both
600     // geometry and media to two ASCII files (called by default
601     // ITSgeometry.euc and ITSgeometry.tme) in a format understandable
602     // to the CAD system EUCLID.  The default (=0) means that you dont
603     // want to use this facility.
604     //
605     ITS->SetEUCLID(0);
606 #endif
607   }
608
609   if (useTPC) {
610     // =================== TPC parameters ============================
611     //
612     // This allows the user to specify sectors for the SLOW (TPC
613     // geometry 2) Simulator. SecAL (SecAU) <0 means that ALL lower
614     // (upper) sectors are specified, any value other than that
615     // requires at least one sector (lower or upper)to be specified!
616     //
617     // Reminder: 
618     //   sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
619     //   sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
620     //
621     //   SecLows - number of lower sectors specified (up to 6)
622     //   SecUps  - number of upper sectors specified (up to 12)
623     //   Sens    - sensitive strips for the Slow Simulator !!!
624     //
625     // This does NOT work if all S or L-sectors are specified, i.e.
626     // if SecAL or SecAU < 0
627     //
628     //
629     //----------------------------------------------------------------
630     //  gROOT->LoadMacro("SetTPCParam.C");
631     //  AliTPCParam *param = SetTPCParam();
632     AliTPC *TPC = new AliTPCv2("TPC", "Default");
633   }
634
635   if (useTOF) {
636     // ================== TOF parameters =============================
637     AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
638   }
639
640   if (useHMPID) {
641     // ================== HMPID parameters ============================
642     AliHMPID *HMPID = new AliHMPIDv1("HMPID", "normal HMPID");
643
644   }
645
646   if (useZDC) {
647     // ================== ZDC parameters =============================
648     AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
649   }
650
651   if (useTRD) {
652     // ================== TRD parameters =============================
653     AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
654
655     // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
656     TRD->SetGasMix(1);
657     if (geo == kHoles) {
658       // With hole in front of PHOS
659       TRD->SetPHOShole();
660       // With hole in front of HMPID
661       TRD->SetHMPIDhole();
662     }
663     // Switch on TR
664     AliTRDsim *TRDsim = TRD->CreateTR();
665   }
666
667   if (useFMD) {
668     // =================== FMD parameters ============================
669     AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
670     // FMD->UseDetailed(kFALSE);
671     // FMD->UseAssembly();
672     // FMD->UseOld();
673   }
674
675   if (useMUON) {
676     // =================== MUON parameters ===========================
677     AliMUON *MUON = new AliMUONv1("MUON", "default");
678     // MUON->AddGeometryBuilder(new AliMUONSt1GeometryBuilder(MUON));
679     // MUON->AddGeometryBuilder(new AliMUONSt2GeometryBuilder(MUON));
680     // MUON->AddGeometryBuilder(new AliMUONSlatGeometryBuilder(MUON));
681     // MUON->AddGeometryBuilder(new AliMUONTriggerGeometryBuilder(MUON));
682   }
683
684   if (usePHOS) {
685     // =================== PHOS parameters ===========================
686     AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
687   }
688
689   if (usePMD) {
690     // =================== PMD parameters ============================
691     AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
692   }
693
694   if (useT0) {
695     // =================== T0 parameters ==========================
696     AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
697   }
698
699   if (useEMCAL) {
700     // =================== EMCAL parameters ==========================
701     AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
702   }
703
704   if (useACORDE) {
705     // =================== ACORDE parameters ============================
706     AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
707   }
708
709   if (useVZERO) {
710     // =================== V0 parameters =============================
711     AliVZERO *VZERO = new AliVZEROv3("VZERO", "normal VZERO");
712   }
713 }
714
715 //____________________________________________________________________
716 Float_t EtaToTheta(Float_t arg)
717 {
718   return (180./TMath::Pi())*2.*TMath::ATan(TMath::Exp(-arg));
719 }
720
721 //____________________________________________________________________
722 Int_t 
723 LookupEG(const Char_t* name) 
724 {
725   TString n(name);
726   for (Int_t i = 0; i < kEgMax; i++) {
727     if (n == egName[i]) return i;
728   }
729   return -1;
730 }
731
732 //____________________________________________________________________  
733 AliGenerator* 
734 GeneratorFactory(EG_t eg, Rad_t rad, TString& comment)  
735 {
736   Int_t isw = 3;
737   if (rad == kNoGluonRadiation) isw = 0;
738   
739   
740   AliGenerator * gGener = 0;
741   switch (eg) {
742   case test50:
743     {
744       comment = comment.Append(":HIJINGparam test 50 particles");
745       AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
746       gener->SetMomentumRange(0, 999999.);
747       gener->SetPhiRange(0., 360.);
748       // Set pseudorapidity range from -8 to 8.
749       Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
750       Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
751       gener->SetThetaRange(thmin,thmax);
752       gGener=gener;
753     }
754     break;
755   case kParam_8000:
756     {
757       comment = comment.Append(":HIJINGparam N=8000");
758       AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
759       gener->SetMomentumRange(0, 999999.);
760       gener->SetPhiRange(0., 360.);
761       // Set pseudorapidity range from -8 to 8.
762       Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
763       Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
764       gener->SetThetaRange(thmin,thmax);
765       gGener=gener;
766     }
767     break;
768   case kParam_4000:
769     {
770       comment = comment.Append("HIJINGparam N=4000");
771       AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
772       gener->SetMomentumRange(0, 999999.);
773       gener->SetPhiRange(0., 360.);
774       // Set pseudorapidity range from -8 to 8.
775       Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
776       Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
777       gener->SetThetaRange(thmin,thmax);
778       gGener=gener;
779     }
780     break;
781   case kParam_2000:
782     {
783       comment = comment.Append("HIJINGparam N=2000");
784       AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
785       gener->SetMomentumRange(0, 999999.);
786       gener->SetPhiRange(0., 360.);
787       // Set pseudorapidity range from -8 to 8.
788       Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
789       Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
790       gener->SetThetaRange(thmin,thmax);
791       gGener=gener;
792     }
793     break;
794   case kParam_fmd:
795     {
796       comment = comment.Append("HIJINGparam N=100");
797       AliGenHIJINGpara *gener = new AliGenHIJINGpara(500);
798       gener->SetMomentumRange(0, 999999.);
799       gener->SetPhiRange(0., 360.);
800       // Set pseudorapidity range from -8 to 8.
801       Float_t thmin = EtaToTheta(6);   // theta min. <---> eta max
802       Float_t thmax = EtaToTheta(2);  // theta max. <---> eta min 
803       gener->SetThetaRange(thmin,thmax);
804       gGener=gener;
805     }
806     break;
807     //
808     //  Hijing Central
809     //
810   case kHijing_cent1:
811     {
812       comment = comment.Append("HIJING cent1");
813       AliGenHijing *gener = HijingStandard();
814       // impact parameter range
815       gener->SetImpactParameterRange(0., 5.);
816       gGener=gener;
817     }
818     break;
819   case kHijing_cent2:
820     {
821       comment = comment.Append("HIJING cent2");
822       AliGenHijing *gener = HijingStandard();
823       // impact parameter range
824       gener->SetImpactParameterRange(0., 2.);
825       gGener=gener;
826     }
827     break;
828     //
829     // Hijing Peripheral 
830     //
831   case kHijing_per1:
832     {
833       comment = comment.Append("HIJING per1");
834       AliGenHijing *gener = HijingStandard();
835       // impact parameter range
836       gener->SetImpactParameterRange(5., 8.6);
837       gGener=gener;
838     }
839     break;
840   case kHijing_per2:
841     {
842       comment = comment.Append("HIJING per2");
843       AliGenHijing *gener = HijingStandard();
844       // impact parameter range
845       gener->SetImpactParameterRange(8.6, 11.2);
846       gGener=gener;
847     }
848     break;
849   case kHijing_per3:
850     {
851       comment = comment.Append("HIJING per3");
852       AliGenHijing *gener = HijingStandard();
853       // impact parameter range
854       gener->SetImpactParameterRange(11.2, 13.2);
855       gGener=gener;
856     }
857     break;
858   case kHijing_per4:
859     {
860       comment = comment.Append("HIJING per4");
861       AliGenHijing *gener = HijingStandard();
862       // impact parameter range
863       gener->SetImpactParameterRange(13.2, 15.);
864       gGener=gener;
865     }
866     break;
867   case kHijing_per5:
868     {
869       comment = comment.Append("HIJING per5");
870       AliGenHijing *gener = HijingStandard();
871       // impact parameter range
872       gener->SetImpactParameterRange(15., 100.);
873       gGener=gener;
874     }
875     break;
876     //
877     //  Jet-Jet
878     //
879   case kHijing_jj25:
880     {
881       comment = comment.Append("HIJING Jet 25 GeV");
882       AliGenHijing *gener = HijingStandard();
883       // impact parameter range
884       gener->SetImpactParameterRange(0., 5.);
885       // trigger
886       gener->SetTrigger(1);
887       gener->SetPtJet(25.);
888       gener->SetRadiation(isw);
889       gener->SetSimpleJets(!isw);
890       gener->SetJetEtaRange(-0.3,0.3);
891       gener->SetJetPhiRange(75., 165.);   
892       gGener=gener;
893     }
894     break;
895
896   case kHijing_jj50:
897     {
898       comment = comment.Append("HIJING Jet 50 GeV");
899       AliGenHijing *gener = HijingStandard();
900       // impact parameter range
901       gener->SetImpactParameterRange(0., 5.);
902       // trigger
903       gener->SetTrigger(1);
904       gener->SetPtJet(50.);
905       gener->SetRadiation(isw);
906       gener->SetSimpleJets(!isw);
907       gener->SetJetEtaRange(-0.3,0.3);
908       gener->SetJetPhiRange(75., 165.);   
909       gGener=gener;
910     }
911     break;
912
913   case kHijing_jj75:
914     {
915       comment = comment.Append("HIJING Jet 75 GeV");
916       AliGenHijing *gener = HijingStandard();
917       // impact parameter range
918       gener->SetImpactParameterRange(0., 5.);
919       // trigger
920       gener->SetTrigger(1);
921       gener->SetPtJet(75.);
922       gener->SetRadiation(isw);
923       gener->SetSimpleJets(!isw);
924       gener->SetJetEtaRange(-0.3,0.3);
925       gener->SetJetPhiRange(75., 165.);   
926       gGener=gener;
927     }
928     break;
929
930   case kHijing_jj100:
931     {
932       comment = comment.Append("HIJING Jet 100 GeV");
933       AliGenHijing *gener = HijingStandard();
934       // impact parameter range
935       gener->SetImpactParameterRange(0., 5.);
936       // trigger
937       gener->SetTrigger(1);
938       gener->SetPtJet(100.);
939       gener->SetRadiation(isw);
940       gener->SetSimpleJets(!isw);
941       gener->SetJetEtaRange(-0.3,0.3);
942       gener->SetJetPhiRange(75., 165.);   
943       gGener=gener;
944     }
945     break;
946
947   case kHijing_jj200:
948     {
949       comment = comment.Append("HIJING Jet 200 GeV");
950       AliGenHijing *gener = HijingStandard();
951       // impact parameter range
952       gener->SetImpactParameterRange(0., 5.);
953       // trigger
954       gener->SetTrigger(1);
955       gener->SetPtJet(200.);
956       gener->SetRadiation(isw);
957       gener->SetSimpleJets(!isw);
958       gener->SetJetEtaRange(-0.3,0.3);
959       gener->SetJetPhiRange(75., 165.);   
960       gGener=gener;
961     }
962     break;
963     //
964     // Gamma-Jet
965     //
966   case kHijing_gj25:
967     {
968       comment = comment.Append("HIJING Gamma 25 GeV");
969       AliGenHijing *gener = HijingStandard();
970       // impact parameter range
971       gener->SetImpactParameterRange(0., 5.);
972       // trigger
973       gener->SetTrigger(2);
974       gener->SetPtJet(25.);
975       gener->SetRadiation(isw);
976       gener->SetSimpleJets(!isw);
977       gener->SetJetEtaRange(-0.12, 0.12);
978       gener->SetJetPhiRange(220., 320.);
979       gGener=gener;
980     }
981     break;
982
983   case kHijing_gj50:
984     {
985       comment = comment.Append("HIJING Gamma 50 GeV");
986       AliGenHijing *gener = HijingStandard();
987       // impact parameter range
988       gener->SetImpactParameterRange(0., 5.);
989       // trigger
990       gener->SetTrigger(2);
991       gener->SetPtJet(50.);
992       gener->SetRadiation(isw);
993       gener->SetSimpleJets(!isw);
994       gener->SetJetEtaRange(-0.12, 0.12);
995       gener->SetJetPhiRange(220., 320.);
996       gGener=gener;
997     }
998     break;
999
1000   case kHijing_gj75:
1001     {
1002       comment = comment.Append("HIJING Gamma 75 GeV");
1003       AliGenHijing *gener = HijingStandard();
1004       // impact parameter range
1005       gener->SetImpactParameterRange(0., 5.);
1006       // trigger
1007       gener->SetTrigger(2);
1008       gener->SetPtJet(75.);
1009       gener->SetRadiation(isw);
1010       gener->SetSimpleJets(!isw);
1011       gener->SetJetEtaRange(-0.12, 0.12);
1012       gener->SetJetPhiRange(220., 320.);
1013       gGener=gener;
1014     }
1015     break;
1016
1017   case kHijing_gj100:
1018     {
1019       comment = comment.Append("HIJING Gamma 100 GeV");
1020       AliGenHijing *gener = HijingStandard();
1021       // impact parameter range
1022       gener->SetImpactParameterRange(0., 5.);
1023       // trigger
1024       gener->SetTrigger(2);
1025       gener->SetPtJet(100.);
1026       gener->SetRadiation(isw);
1027       gener->SetSimpleJets(!isw);
1028       gener->SetJetEtaRange(-0.12, 0.12);
1029       gener->SetJetPhiRange(220., 320.);
1030       gGener=gener;
1031     }
1032     break;
1033
1034   case kHijing_gj200:
1035     {
1036       comment = comment.Append("HIJING Gamma 200 GeV");
1037       AliGenHijing *gener = HijingStandard();
1038       // impact parameter range
1039       gener->SetImpactParameterRange(0., 5.);
1040       // trigger
1041       gener->SetTrigger(2);
1042       gener->SetPtJet(200.);
1043       gener->SetRadiation(isw);
1044       gener->SetSimpleJets(!isw);
1045       gener->SetJetEtaRange(-0.12, 0.12);
1046       gener->SetJetPhiRange(220., 320.);
1047       gGener=gener;
1048     }
1049     break;
1050   case kHijing_pA:
1051     {
1052       comment = comment.Append("HIJING pA");
1053
1054       AliGenCocktail *gener  = new AliGenCocktail();
1055
1056       AliGenHijing   *hijing = new AliGenHijing(-1);
1057       // centre of mass energy 
1058       hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
1059       // impact parameter range
1060       hijing->SetImpactParameterRange(0., 15.);
1061       // reference frame
1062       hijing->SetReferenceFrame("CMS");
1063       hijing->SetBoostLHC(1);
1064       // projectile
1065       hijing->SetProjectile("P", 1, 1);
1066       hijing->SetTarget    ("A", 208, 82);
1067       // tell hijing to keep the full parent child chain
1068       hijing->KeepFullEvent();
1069       // enable jet quenching
1070       hijing->SetJetQuenching(0);
1071       // enable shadowing
1072       hijing->SetShadowing(1);
1073       // Don't track spectators
1074       hijing->SetSpectators(0);
1075       // kinematic selection
1076       hijing->SetSelectAll(0);
1077       //
1078       AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
1079       AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
1080       gray->SetSlowNucleonModel(model);
1081       gray->SetDebug(1);
1082       gener->AddGenerator(hijing,"Hijing pPb", 1);
1083       gener->AddGenerator(gray,  "Gray Particles",1);
1084       gGener=gener;
1085     }
1086     break;
1087   case kPythia6:
1088     {
1089       comment = comment.Append(":Pythia p-p @ 14 TeV");
1090       AliGenPythia *gener = new AliGenPythia(-1); 
1091       gener->SetMomentumRange(0,999999);
1092       gener->SetThetaRange(0., 180.);
1093       gener->SetYRange(-12,12);
1094       gener->SetPtRange(0,1000);
1095       gener->SetProcess(kPyMb);
1096       gener->SetEnergyCMS(14000.);
1097       gGener=gener;
1098     }
1099     break;
1100   case kPythia6Jets20_24:
1101     {
1102       comment = comment.Append(":Pythia jets 20-24 GeV @ 5.5 TeV");
1103       AliGenPythia * gener = new AliGenPythia(-1);
1104       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1105       gener->SetProcess(kPyJets);//        Process type
1106       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1107       gener->SetJetPhiRange(0., 360.);
1108       gener->SetJetEtRange(10., 1000.);
1109       gener->SetGluonRadiation(1,1);
1110       //    gener->SetPtKick(0.);
1111       //   Structure function
1112       gener->SetStrucFunc(kCTEQ4L);
1113       gener->SetPtHard(20., 24.);// Pt transfer of the hard scattering
1114       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1115       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1116       gGener=gener;
1117     }
1118     break;
1119   case kPythia6Jets24_29:
1120     {
1121       comment = comment.Append(":Pythia jets 24-29 GeV @ 5.5 TeV");
1122       AliGenPythia * gener = new AliGenPythia(-1);
1123       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1124       gener->SetProcess(kPyJets);//        Process type
1125       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1126       gener->SetJetPhiRange(0., 360.);
1127       gener->SetJetEtRange(10., 1000.);
1128       gener->SetGluonRadiation(1,1);
1129       //    gener->SetPtKick(0.);
1130       //   Structure function
1131       gener->SetStrucFunc(kCTEQ4L);
1132       gener->SetPtHard(24., 29.);// Pt transfer of the hard scattering
1133       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1134       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1135       gGener=gener;
1136     }
1137     break;
1138   case kPythia6Jets29_35:
1139     {
1140       comment = comment.Append(":Pythia jets 29-35 GeV @ 5.5 TeV");
1141       AliGenPythia * gener = new AliGenPythia(-1);
1142       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1143       gener->SetProcess(kPyJets);//        Process type
1144       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1145       gener->SetJetPhiRange(0., 360.);
1146       gener->SetJetEtRange(10., 1000.);
1147       gener->SetGluonRadiation(1,1);
1148       //    gener->SetPtKick(0.);
1149       //   Structure function
1150       gener->SetStrucFunc(kCTEQ4L);
1151       gener->SetPtHard(29., 35.);// Pt transfer of the hard scattering
1152       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1153       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1154       gGener=gener;
1155     }
1156     break;
1157   case kPythia6Jets35_42:
1158     {
1159       comment = comment.Append(":Pythia jets 35-42 GeV @ 5.5 TeV");
1160       AliGenPythia * gener = new AliGenPythia(-1);
1161       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1162       gener->SetProcess(kPyJets);//        Process type
1163       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1164       gener->SetJetPhiRange(0., 360.);
1165       gener->SetJetEtRange(10., 1000.);
1166       gener->SetGluonRadiation(1,1);
1167       //    gener->SetPtKick(0.);
1168       //   Structure function
1169       gener->SetStrucFunc(kCTEQ4L);
1170       gener->SetPtHard(35., 42.);// Pt transfer of the hard scattering
1171       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1172       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1173       gGener=gener;
1174     }
1175     break;
1176   case kPythia6Jets42_50:
1177     {
1178       comment = comment.Append(":Pythia jets 42-50 GeV @ 5.5 TeV");
1179       AliGenPythia * gener = new AliGenPythia(-1);
1180       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1181       gener->SetProcess(kPyJets);//        Process type
1182       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1183       gener->SetJetPhiRange(0., 360.);
1184       gener->SetJetEtRange(10., 1000.);
1185       gener->SetGluonRadiation(1,1);
1186       //    gener->SetPtKick(0.);
1187       //   Structure function
1188       gener->SetStrucFunc(kCTEQ4L);
1189       gener->SetPtHard(42., 50.);// Pt transfer of the hard scattering
1190       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1191       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1192       gGener=gener;
1193     }
1194     break;
1195   case kPythia6Jets50_60:
1196     {
1197       comment = comment.Append(":Pythia jets 50-60 GeV @ 5.5 TeV");
1198       AliGenPythia * gener = new AliGenPythia(-1);
1199       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1200       gener->SetProcess(kPyJets);//        Process type
1201       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1202       gener->SetJetPhiRange(0., 360.);
1203       gener->SetJetEtRange(10., 1000.);
1204       gener->SetGluonRadiation(1,1);
1205       //    gener->SetPtKick(0.);
1206       //   Structure function
1207       gener->SetStrucFunc(kCTEQ4L);
1208       gener->SetPtHard(50., 60.);// Pt transfer of the hard scattering
1209       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1210       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1211       gGener=gener;
1212     }
1213     break;
1214   case kPythia6Jets60_72:
1215     {
1216       comment = comment.Append(":Pythia jets 60-72 GeV @ 5.5 TeV");
1217       AliGenPythia * gener = new AliGenPythia(-1);
1218       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1219       gener->SetProcess(kPyJets);//        Process type
1220       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1221       gener->SetJetPhiRange(0., 360.);
1222       gener->SetJetEtRange(10., 1000.);
1223       gener->SetGluonRadiation(1,1);
1224       //    gener->SetPtKick(0.);
1225       //   Structure function
1226       gener->SetStrucFunc(kCTEQ4L);
1227       gener->SetPtHard(60., 72.);// Pt transfer of the hard scattering
1228       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1229       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1230       gGener=gener;
1231     }
1232     break;
1233   case kPythia6Jets72_86:
1234     {
1235       comment = comment.Append(":Pythia jets 72-86 GeV @ 5.5 TeV");
1236       AliGenPythia * gener = new AliGenPythia(-1);
1237       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1238       gener->SetProcess(kPyJets);//        Process type
1239       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1240       gener->SetJetPhiRange(0., 360.);
1241       gener->SetJetEtRange(10., 1000.);
1242       gener->SetGluonRadiation(1,1);
1243       //    gener->SetPtKick(0.);
1244       //   Structure function
1245       gener->SetStrucFunc(kCTEQ4L);
1246       gener->SetPtHard(72., 86.);// Pt transfer of the hard scattering
1247       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1248       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1249       gGener=gener;
1250     }
1251     break;
1252   case kPythia6Jets86_104:
1253     {
1254       comment = comment.Append(":Pythia jets 86-104 GeV @ 5.5 TeV");
1255       AliGenPythia * gener = new AliGenPythia(-1);
1256       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1257       gener->SetProcess(kPyJets);//        Process type
1258       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1259       gener->SetJetPhiRange(0., 360.);
1260       gener->SetJetEtRange(10., 1000.);
1261       gener->SetGluonRadiation(1,1);
1262       //    gener->SetPtKick(0.);
1263       //   Structure function
1264       gener->SetStrucFunc(kCTEQ4L);
1265       gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
1266       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1267       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1268       gGener=gener;
1269     }
1270     break;
1271   case kPythia6Jets104_125:
1272     {
1273       comment = comment.Append(":Pythia jets 105-125 GeV @ 5.5 TeV");
1274       AliGenPythia * gener = new AliGenPythia(-1);
1275       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1276       gener->SetProcess(kPyJets);//        Process type
1277       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1278       gener->SetJetPhiRange(0., 360.);
1279       gener->SetJetEtRange(10., 1000.);
1280       gener->SetGluonRadiation(1,1);
1281       //    gener->SetPtKick(0.);
1282       //   Structure function
1283       gener->SetStrucFunc(kCTEQ4L);
1284       gener->SetPtHard(104., 125.);// Pt transfer of the hard scattering
1285       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1286       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1287       gGener=gener;
1288     }
1289     break;
1290   case kPythia6Jets125_150:
1291     {
1292       comment = comment.Append(":Pythia jets 125-150 GeV @ 5.5 TeV");
1293       AliGenPythia * gener = new AliGenPythia(-1);
1294       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1295       gener->SetProcess(kPyJets);//        Process type
1296       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1297       gener->SetJetPhiRange(0., 360.);
1298       gener->SetJetEtRange(10., 1000.);
1299       gener->SetGluonRadiation(1,1);
1300       //    gener->SetPtKick(0.);
1301       //   Structure function
1302       gener->SetStrucFunc(kCTEQ4L);
1303       gener->SetPtHard(125., 150.);// Pt transfer of the hard scattering
1304       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1305       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1306       gGener=gener;
1307     }
1308     break;
1309   case kPythia6Jets150_180:
1310     {
1311       comment = comment.Append(":Pythia jets 150-180 GeV @ 5.5 TeV");
1312       AliGenPythia * gener = new AliGenPythia(-1);
1313       gener->SetEnergyCMS(5500.);//        Centre of mass energy
1314       gener->SetProcess(kPyJets);//        Process type
1315       gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
1316       gener->SetJetPhiRange(0., 360.);
1317       gener->SetJetEtRange(10., 1000.);
1318       gener->SetGluonRadiation(1,1);
1319       //    gener->SetPtKick(0.);
1320       //   Structure function
1321       gener->SetStrucFunc(kCTEQ4L);
1322       gener->SetPtHard(150., 180.);// Pt transfer of the hard scattering
1323       gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1324       gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
1325       gGener=gener;
1326     }
1327     break;
1328   case kD0PbPb5500:
1329     {
1330       comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
1331       AliGenPythia * gener = new AliGenPythia(10);
1332       gener->SetProcess(kPyD0PbPbMNR);
1333       gener->SetStrucFunc(kCTEQ4L);
1334       gener->SetPtHard(2.1,-1.0);
1335       gener->SetEnergyCMS(5500.);
1336       gener->SetNuclei(208,208);
1337       gener->SetForceDecay(kHadronicD);
1338       gener->SetYRange(-2,2);
1339       gener->SetFeedDownHigherFamily(kFALSE);
1340       gener->SetStackFillOpt(AliGenPythia::kParentSelection);
1341       gener->SetCountMode(AliGenPythia::kCountParents);
1342       gGener=gener;
1343     }
1344     break;
1345   case kCharmSemiElPbPb5500:
1346     {
1347       comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
1348       AliGenPythia * gener = new AliGenPythia(10);
1349       gener->SetProcess(kPyCharmPbPbMNR);
1350       gener->SetStrucFunc(kCTEQ4L);
1351       gener->SetPtHard(2.1,-1.0);
1352       gener->SetEnergyCMS(5500.);
1353       gener->SetNuclei(208,208);
1354       gener->SetForceDecay(kSemiElectronic);
1355       gener->SetYRange(-2,2);
1356       gener->SetFeedDownHigherFamily(kFALSE);
1357       gener->SetCountMode(AliGenPythia::kCountParents);
1358       gGener=gener;
1359     }
1360     break;
1361   case kBeautySemiElPbPb5500:
1362     {
1363       comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
1364       AliGenPythia *gener = new AliGenPythia(10);
1365       gener->SetProcess(kPyBeautyPbPbMNR);
1366       gener->SetStrucFunc(kCTEQ4L);
1367       gener->SetPtHard(2.75,-1.0);
1368       gener->SetEnergyCMS(5500.);
1369       gener->SetNuclei(208,208);
1370       gener->SetForceDecay(kSemiElectronic);
1371       gener->SetYRange(-2,2);
1372       gener->SetFeedDownHigherFamily(kFALSE);
1373       gener->SetCountMode(AliGenPythia::kCountParents);
1374       gGener=gener;
1375     }
1376     break;
1377   case kCocktailTRD:
1378     {
1379       comment = comment.Append(" Cocktail for TRD at 5.5 TeV");
1380       AliGenCocktail *gener  = new AliGenCocktail();
1381
1382       AliGenParam *jpsi = new AliGenParam(10,
1383                                           new AliGenMUONlib(),
1384                                           AliGenMUONlib::kJpsiFamily,
1385                                           "Vogt PbPb");
1386
1387       jpsi->SetPtRange(0, 100);
1388       jpsi->SetYRange(-1., +1.);
1389       jpsi->SetForceDecay(kDiElectron);
1390
1391       AliGenParam *ups = new AliGenParam(10,
1392                                          new AliGenMUONlib(),
1393                                          AliGenMUONlib::kUpsilonFamily,
1394                                          "Vogt PbPb");
1395       ups->SetPtRange(0, 100);
1396       ups->SetYRange(-1., +1.);
1397       ups->SetForceDecay(kDiElectron);
1398         
1399       AliGenParam *charm = new AliGenParam(10,
1400                                            new AliGenMUONlib(), 
1401                                            AliGenMUONlib::kCharm,
1402                                            "central");
1403       charm->SetPtRange(0, 100);
1404       charm->SetYRange(-1.5, +1.5);
1405       charm->SetForceDecay(kSemiElectronic);
1406         
1407         
1408       AliGenParam *beauty = new AliGenParam(10,
1409                                             new AliGenMUONlib(), 
1410                                             AliGenMUONlib::kBeauty,
1411                                             "central");
1412       beauty->SetPtRange(0, 100);
1413       beauty->SetYRange(-1.5, +1.5);
1414       beauty->SetForceDecay(kSemiElectronic);
1415
1416       gener->AddGenerator(jpsi,"J/psi",1);
1417       gener->AddGenerator(ups,"Upsilon",1);
1418       gener->AddGenerator(charm,"Charm",1);
1419       gener->AddGenerator(beauty,"Beauty",1);
1420       gGener=gener;
1421     }
1422     break;
1423   case kPyJJ:
1424     {
1425       comment = comment.Append(" Jet-jet at 5.5 TeV");
1426       AliGenPythia *gener = new AliGenPythia(-1);
1427       gener->SetEnergyCMS(5500.);
1428       gener->SetProcess(kPyJets);
1429       Double_t ptHardMin=10.0, ptHardMax=-1.0;
1430       gener->SetPtHard(ptHardMin,ptHardMax);
1431       gener->SetYHard(-0.7,0.7);
1432       gener->SetJetEtaRange(-0.2,0.2);
1433       gener->SetEventListRange(0,1);
1434       gGener=gener;
1435     }
1436     break;
1437   case kPyGJ:
1438     {
1439       comment = comment.Append(" Gamma-jet at 5.5 TeV");
1440       AliGenPythia *gener = new AliGenPythia(-1);
1441       gener->SetEnergyCMS(5500.);
1442       gener->SetProcess(kPyDirectGamma);
1443       Double_t ptHardMin=10.0, ptHardMax=-1.0;
1444       gener->SetPtHard(ptHardMin,ptHardMax);
1445       gener->SetYHard(-1.0,1.0);
1446       gener->SetGammaEtaRange(-0.13,0.13);
1447       gener->SetGammaPhiRange(210.,330.);
1448       gener->SetEventListRange(0,1);
1449       gGener=gener;
1450     }
1451     break;
1452   case kMuonCocktailCent1:
1453     {
1454       comment = comment.Append(" Muon Cocktail Cent1");
1455       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1456       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1457       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1458       gener->SetYRange(-4.0,-2.4);
1459       gener->SetMuonPtCut(0.8);
1460       gener->SetMuonThetaCut(171.,178.);
1461       gener->SetMuonMultiplicity(2);
1462       gener->SetNumberOfCollisions(1626.);  //Centrality class Cent1 for PDC04
1463       gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1464       gGener=gener;
1465     }
1466     break;
1467   case kMuonCocktailPer1:
1468     {
1469       comment = comment.Append(" Muon Cocktail Per1");
1470       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1471       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1472       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1473       gener->SetYRange(-4.0,-2.4);
1474       gener->SetMuonPtCut(0.8);
1475       gener->SetMuonThetaCut(171.,178.);
1476       gener->SetMuonMultiplicity(2);
1477       gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1478       gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1479       gGener=gener;
1480     }
1481     break;
1482   case kMuonCocktailPer4:
1483     {
1484       comment = comment.Append(" Muon Cocktail Per4");
1485       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1486       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1487       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1488       gener->SetYRange(-4.0,-2.4);
1489       gener->SetMuonPtCut(0.8);
1490       gener->SetMuonThetaCut(171.,178.);
1491       gener->SetMuonMultiplicity(2);
1492       gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1493       gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1494       gGener=gener;
1495     }
1496     break;
1497   case kMuonCocktailCent1HighPt:
1498     {
1499       comment = comment.Append(" Muon Cocktail HighPt Cent1");
1500       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1501       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1502       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1503       gener->SetYRange(-4.0,-2.4);
1504       gener->SetMuonPtCut(2.5);
1505       gener->SetMuonThetaCut(171.,178.);
1506       gener->SetMuonMultiplicity(2);
1507       gener->SetNumberOfCollisions(1626.);  //Centrality class Cent1 for PDC04
1508       gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1509       gGener=gener;
1510     }
1511     break;
1512   case kMuonCocktailPer1HighPt :
1513     {
1514       comment = comment.Append(" Muon Cocktail HighPt Per1");
1515       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1516       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1517       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1518       gener->SetYRange(-4.0,-2.4);
1519       gener->SetMuonPtCut(2.5);
1520       gener->SetMuonThetaCut(171.,178.);
1521       gener->SetMuonMultiplicity(2);
1522       gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1523       gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1524       gGener=gener;
1525     }
1526     break;
1527   case kMuonCocktailPer4HighPt:
1528     {
1529       comment = comment.Append(" Muon Cocktail HighPt Per4");
1530       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1531       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1532       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1533       gener->SetYRange(-4.0,-2.4);
1534       gener->SetMuonPtCut(2.5);
1535       gener->SetMuonThetaCut(171.,178.);
1536       gener->SetMuonMultiplicity(2);
1537       gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1538       gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1539       gGener=gener;
1540     }
1541     break;
1542   case kMuonCocktailCent1Single:
1543     {
1544       comment = comment.Append(" Muon Cocktail Single Cent1");
1545       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1546       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1547       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1548       gener->SetYRange(-4.0,-2.4);
1549       gener->SetMuonPtCut(0.8);
1550       gener->SetMuonThetaCut(171.,178.);
1551       gener->SetMuonMultiplicity(1);
1552       gener->SetNumberOfCollisions(1626.);  //Centrality class Cent1 for PDC04
1553       gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1554       gGener=gener;
1555     }
1556     break;
1557   case kMuonCocktailPer1Single :
1558     {
1559       comment = comment.Append(" Muon Cocktail Single Per1");
1560       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1561       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1562       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1563       gener->SetYRange(-4.0,-2.4);
1564       gener->SetMuonPtCut(0.8);
1565       gener->SetMuonThetaCut(171.,178.);
1566       gener->SetMuonMultiplicity(1);
1567       gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1568       gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1569       gGener=gener;
1570     }
1571     break;
1572   case kMuonCocktailPer4Single:
1573     {
1574       comment = comment.Append(" Muon Cocktail Single Per4");
1575       AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1576       gener->SetPtRange(1.0,100.);       // Transverse momentum range   
1577       gener->SetPhiRange(0.,360.);    // Azimuthal angle range  
1578       gener->SetYRange(-4.0,-2.4);
1579       gener->SetMuonPtCut(0.8);
1580       gener->SetMuonThetaCut(171.,178.);
1581       gener->SetMuonMultiplicity(1);
1582       gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1583       gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1584       gGener=gener;
1585     }
1586     break;
1587   case kFMD1Flat: 
1588     {
1589       comment = comment.Append(" Flat in FMD1 range");
1590       AliGenBox* gener = new AliGenBox(2000);
1591       gener->SetPart(kPiPlus);
1592       gener->SetMomentumRange(3,4);
1593       gener->SetPhiRange(0, 360);
1594       gener->SetThetaRange(0.77, 3.08);
1595       gGener = gener;
1596     }
1597     break;
1598   case kFMD2Flat: 
1599     {
1600       comment = comment.Append(" Flat in FMD2 range");
1601       AliGenBox* gener = new AliGenBox(2000);
1602       gener->SetPart(kPiPlus);
1603       gener->SetMomentumRange(3,4);
1604       gener->SetPhiRange(0, 360);
1605       gener->SetThetaRange(2.95, 20.42);
1606       gGener = gener;
1607     }
1608     break;
1609   case kFMD3Flat: 
1610     {
1611       comment = comment.Append(" Flat in FMD3 range");
1612       AliGenBox* gener = new AliGenBox(2000);
1613       gener->SetPart(kPiPlus);
1614       gener->SetMomentumRange(3,4);
1615       gener->SetPhiRange(0, 360);
1616       gener->SetThetaRange(155.97, 176.73);
1617       gGener = gener;
1618     }
1619     break;
1620   case kFMDFlat:
1621     {
1622       comment = comment.Append(" Flat in FMD range");
1623       AliGenCocktail* gener = new AliGenCocktail();
1624       gener->SetMomentumRange(3,4);
1625       gener->SetPhiRange(0, 360);
1626       AliGenBox* gener3 = new AliGenBox(2000);
1627       gener3->SetThetaRange(155.97, 176.73);
1628       gener3->SetPart(kPiPlus);
1629       gener->AddGenerator(gener3, "FMD3", .33);
1630       AliGenBox* gener2 = new AliGenBox(2000);
1631       gener2->SetThetaRange(2.95, 20.42);
1632       gener2->SetPart(kPiPlus);
1633       gener->AddGenerator(gener2, "FMD2", .33);
1634       AliGenBox* gener1 = new AliGenBox(2000);
1635       gener1->SetThetaRange(0.77, 3.08);
1636       gener1->SetPart(kPiPlus);
1637       gener->AddGenerator(gener1, "FMD1", .34);
1638       gGener = gener;
1639     }
1640     break;
1641     
1642   default: break;
1643   }
1644   return gGener;
1645 }
1646
1647 //____________________________________________________________________
1648 AliGenHijing* 
1649 HijingStandard()
1650 {
1651   AliGenHijing *gener = new AliGenHijing(-1);
1652   // centre of mass energy 
1653   gener->SetEnergyCMS(5500.);
1654   // reference frame
1655   gener->SetReferenceFrame("CMS");
1656   // projectile
1657   gener->SetProjectile("A", 208, 82);
1658   gener->SetTarget    ("A", 208, 82);
1659   // tell hijing to keep the full parent child chain
1660   gener->KeepFullEvent();
1661   // enable jet quenching
1662   gener->SetJetQuenching(1);
1663   // enable shadowing
1664   gener->SetShadowing(1);
1665   // neutral pion and heavy particle decays switched off
1666   gener->SetDecaysOff(1);
1667   // Don't track spectators
1668   gener->SetSpectators(0);
1669   // kinematic selection
1670   gener->SetSelectAll(0);
1671   return gener;
1672 }
1673
1674
1675 //____________________________________________________________________
1676 void 
1677 ProcessEnvironmentVars(EG_t& eg, Int_t& seed)
1678 {
1679   // Run type
1680   if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
1681     Int_t eg1 = LookupEG(gSystem->Getenv("CONFIG_RUN_TYPE"));
1682     if  (eg1 >= 0) eg = EG_t(eg1);
1683   }
1684   // Random Number seed
1685   if (gSystem->Getenv("CONFIG_SEED")) {
1686     seed = atoi(gSystem->Getenv("CONFIG_SEED"));
1687   }
1688 }
1689
1690 //____________________________________________________________________
1691 //
1692 // EOF
1693 //