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