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