]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFconfig.C
Remove a run from LHC10e wo AODs; added run list for the new PbPb MC LHC11a10b_bis
[u/mrichter/AliRoot.git] / TOF / AliTOFconfig.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <Riostream.h>
3 #include <TRandom.h>
4 #include <TSystem.h>
5 #include <TVirtualMC.h>
6 #include <TGeant3.h>
7 #include "STEER/AliRunLoader.h"
8 #include "STEER/AliRun.h"
9 #include "STEER/AliConfig.h"
10 #include "STEER/AliGenerator.h"
11 #include "PYTHIA6/AliDecayerPythia.h"
12 #include "EVGEN/AliGenHIJINGpara.h"
13 #include "THijing/AliGenHijing.h"
14 #include "EVGEN/AliGenCocktail.h"
15 #include "EVGEN/AliGenSlowNucleons.h"
16 #include "EVGEN/AliSlowNucleonModelExp.h"
17 #include "PYTHIA6/AliGenPythia.h"
18 #include "STEER/AliMagF.h"
19 #include "STRUCT/AliBODY.h"
20 #include "STRUCT/AliMAG.h"
21 #include "STRUCT/AliABSOv0.h"
22 #include "STRUCT/AliDIPOv2.h"
23 #include "STRUCT/AliHALL.h"
24 #include "STRUCT/AliFRAMEv2.h"
25 #include "STRUCT/AliSHILv2.h"
26 #include "STRUCT/AliPIPEv0.h"
27 #include "ITS/AliITSvPPRasymmFMD.h"
28 #include "TPC/AliTPCv2.h"
29 #include "TOF/AliTOFv4T0.h"
30 #include "HMPID/AliHMPIDv1.h"
31 #include "ZDC/AliZDCv2.h"
32 #include "TRD/AliTRDv1.h"
33 #include "FMD/AliFMDv1.h"
34 #include "MUON/AliMUONv1.h"
35 #include "PHOS/AliPHOSv1.h"
36 #include "PMD/AliPMDv1.h"
37 #include "T0/AliT0v1.h"
38 #include "EMCAL/AliEMCALv1.h"
39 #include "ACORDE/AliACORDEv1.h"
40 #include "VZERO/AliVZEROv3.h"
41 #endif
42
43 enum PprRun_t 
44 {
45     test50,
46     kParam_8000,   kParam_4000,  kParam_2000, 
47     kHijing_cent1, kHijing_cent2, 
48     kHijing_per1,  kHijing_per2, kHijing_per3, kHijing_per4,  kHijing_per5,
49     kHijing_jj25,  kHijing_jj50, kHijing_jj75, kHijing_jj100, kHijing_jj200, 
50     kHijing_gj25,  kHijing_gj50, kHijing_gj75, kHijing_gj100, kHijing_gj200,
51     kHijing_pA, kPythia6, kPythia6Jets, kD0PbPb5500
52 };
53
54 enum PprGeo_t 
55 {
56     kHoles, kNoHoles
57 };
58
59 enum PprRad_t
60 {
61     kGluonRadiation, kNoGluonRadiation
62 };
63
64 // This part for configuration    
65 static PprRun_t srun = test50;
66 //static PprRun_t srun = kPythia6;
67 static PprGeo_t sgeo = kNoHoles;
68 static PprRad_t srad = kGluonRadiation;
69 //static AliMagF::BMap_t smag = AliMagF::k5kG;
70
71 // Comment line 
72 static TString  comment;
73
74 // Functions
75 Float_t EtaToTheta(Float_t arg);
76 AliGenerator* GeneratorFactory(PprRun_t srun);
77 AliGenHijing* HijingStandard();
78
79 void Config()
80 {
81     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
82     // Theta range given through pseudorapidity limits 22/6/2001
83
84     // Set Random Number seed
85   gRandom->SetSeed(12345); //Set 0 to use the current time
86     cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
87
88
89    // libraries required by geant321
90 #if defined(__CINT__)
91     gSystem->Load("libgeant321");
92 #endif
93
94     new     TGeant3("C++ Interface to Geant3");
95
96     AliRunLoader* rl=0x0;
97
98     cout<<"Config.C: Creating Run Loader ..."<<endl;
99     rl = AliRunLoader::Open("galice.root",
100                             AliConfig::GetDefaultEventFolderName(),
101                             "recreate");
102     if (rl == 0x0)
103       {
104         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
105         return;
106       }
107     rl->SetCompressionLevel(2);
108     rl->SetNumberOfEventsPerFile(3);
109     gAlice->SetRunLoader(rl);
110
111     //
112     // Set External decayer
113     AliDecayer *decayer = new AliDecayerPythia();
114
115     decayer->SetForceDecay(kAll);
116     decayer->Init();
117     gMC->SetExternalDecayer(decayer);
118     //
119     //
120     //=======================================================================
121     //
122     //=======================================================================
123     // ************* STEERING parameters FOR ALICE SIMULATION **************
124     // --- Specify event type to be tracked through the ALICE setup
125     // --- All positions are in cm, angles in degrees, and P and E in GeV
126
127     gMC->SetProcess("DCAY",1);
128     gMC->SetProcess("PAIR",1);
129     gMC->SetProcess("COMP",1);
130     gMC->SetProcess("PHOT",1);
131     gMC->SetProcess("PFIS",0);
132     gMC->SetProcess("DRAY",0);
133     gMC->SetProcess("ANNI",1);
134     gMC->SetProcess("BREM",1);
135     gMC->SetProcess("MUNU",1);
136     gMC->SetProcess("CKOV",1);
137     gMC->SetProcess("HADR",1);
138     gMC->SetProcess("LOSS",2);
139     gMC->SetProcess("MULS",1);
140     gMC->SetProcess("RAYL",1);
141
142     Float_t cut = 1.e-3;        // 1MeV cut by default
143     Float_t tofmax = 1.e10;
144
145     gMC->SetCut("CUTGAM", cut);
146     gMC->SetCut("CUTELE", cut);
147     gMC->SetCut("CUTNEU", cut);
148     gMC->SetCut("CUTHAD", cut);
149     gMC->SetCut("CUTMUO", cut);
150     gMC->SetCut("BCUTE",  cut); 
151     gMC->SetCut("BCUTM",  cut); 
152     gMC->SetCut("DCUTE",  cut); 
153     gMC->SetCut("DCUTM",  cut); 
154     gMC->SetCut("PPCUTM", cut);
155     gMC->SetCut("TOFMAX", tofmax); 
156
157
158 // Generator Configuration
159     gAlice->SetDebug(1);
160     AliGenerator* gener = GeneratorFactory(srun);
161     gener->SetOrigin(0, 0, 0);    // vertex position
162     gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
163     gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
164     gener->SetVertexSmear(kPerEvent); 
165     gener->SetTrackingFlag(1);
166     gener->Init();
167
168     /*    
169     if (smag == AliMagF::k2kG) {
170         comment = comment.Append(" | L3 field 0.2 T");
171     } else if (smag == AliMagF::k5kG) {
172         comment = comment.Append(" | L3 field 0.5 T");
173     }
174     */
175     
176     if (srad == kGluonRadiation)
177     {
178         comment = comment.Append(" | Gluon Radiation On");
179         
180     } else {
181         comment = comment.Append(" | Gluon Radiation Off");
182     }
183
184     if (sgeo == kHoles)
185     {
186         comment = comment.Append(" | Holes for PHOS/HMPID");
187         
188     } else {
189         comment = comment.Append(" | No holes for PHOS/HMPID");
190     }
191
192     printf("\n \n Comment: %s \n \n", comment.Data());
193     
194     
195 // Field (L3 +-0.2/0.5 T)
196     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps",-1,-1, AliMagF::k5kG)); // +0.5/-0.7
197     /*
198     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 1, 1, AliMagF::k5kG)); // -0.5/-0.7
199     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps",-1,-1, AliMagF::k2kG)); // +0.2/-0.7
200     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 1, 1, AliMagF::k2kG)); // -0.2/+0.7
201     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 1, 0, AliMagF::k5kG)); // +-0.5/0
202     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 1, 0, AliMagF::k2kG)); // +-0.2/0
203     TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 0, 0, AliMagF::k5kGUniform)); // 0/0
204     */
205
206 //
207     Int_t   iABSO   = 0;
208     Int_t   iDIPO   = 0;
209     Int_t   iFMD    = 0;
210     Int_t   iFRAME  = 1;
211     Int_t   iHALL   = 0;
212     Int_t   iITS    = 1;
213     Int_t   iMAG    = 1;
214     Int_t   iMUON   = 0;
215     Int_t   iPHOS   = 0;
216     Int_t   iPIPE   = 1;
217     Int_t   iPMD    = 0;
218     Int_t   iHMPID   = 0;
219     Int_t   iSHIL   = 0;
220     Int_t   iT0  = 0;
221     Int_t   iTOF    = 1;
222     Int_t   iTPC    = 1;
223     Int_t   iTRD    = 1;
224     Int_t   iZDC    = 0;
225     Int_t   iEMCAL  = 0;
226     Int_t   iVZERO  = 0;
227     Int_t   iACORDE    = 0;
228
229     //=================== Alice BODY parameters =============================
230     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
231
232
233     if (iMAG)
234     {
235         //=================== MAG parameters ============================
236         // --- Start with Magnet since detector layouts may be depending ---
237         // --- on the selected Magnet dimensions ---
238         AliMAG *MAG = new AliMAG("MAG", "Magnet");
239     }
240
241
242     if (iABSO)
243     {
244         //=================== ABSO parameters ============================
245         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
246     }
247
248     if (iDIPO)
249     {
250         //=================== DIPO parameters ============================
251
252         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
253     }
254
255     if (iHALL)
256     {
257         //=================== HALL parameters ============================
258
259         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
260     }
261
262
263     if (iFRAME)
264     {
265         //=================== FRAME parameters ============================
266
267         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
268         if (sgeo == kHoles) {
269             FRAME->SetHoles(1);
270         } else {
271             FRAME->SetHoles(0);
272         }
273     }
274
275     if (iSHIL)
276     {
277         //=================== SHIL parameters ============================
278
279         AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
280     }
281
282
283     if (iPIPE)
284     {
285         //=================== PIPE parameters ============================
286
287         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
288     }
289  
290     if(iITS) {
291
292     //=================== ITS parameters ============================
293     //
294     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
295     // almost all other detectors. This involves the fact that the ITS geometry
296     // still has several options to be followed in parallel in order to determine
297     // the best set-up which minimizes the induced background. All the geometries
298     // available to date are described in the following. Read carefully the comments
299     // and use the default version (the only one uncommented) unless you are making
300     // comparisons and you know what you are doing. In this case just uncomment the
301     // ITS geometry you want to use and run Aliroot.
302     //
303     // Detailed geometries:         
304     //
305     //
306     //AliITS *ITS  = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
307     //
308     //AliITS *ITS  = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
309     //
310         AliITSvPPRasymmFMD *ITS  = new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version with asymmetric services");
311         ITS->SetMinorVersion(2);  // don't touch this parameter if you're not an ITS developer
312         ITS->SetReadDet(kTRUE);   // don't touch this parameter if you're not an ITS developer
313     //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
314         ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
315         ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
316         ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
317         ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
318         ITS->SetRails(0);            // 1 --> rails in ; 0 --> rails out
319         ITS->SetCoolingFluid(1);   // 1 --> water ; 0 --> freon
320
321     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
322     // for reconstruction !):
323     //                                                     
324     //
325     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
326     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
327     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
328     //
329     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
330     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
331     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
332     //                      
333     //
334     //
335     // Geant3 <-> EUCLID conversion
336     // ============================
337     //
338     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
339     // media to two ASCII files (called by default ITSgeometry.euc and
340     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
341     // The default (=0) means that you dont want to use this facility.
342     //
343         ITS->SetEUCLID(0);  
344     }
345
346     if (iTPC)
347     {
348         //============================ TPC parameters ================================
349         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
350         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
351         // --- sectors are specified, any value other than that requires at least one 
352         // --- sector (lower or upper)to be specified!
353         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
354         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
355         // --- SecLows - number of lower sectors specified (up to 6)
356         // --- SecUps - number of upper sectors specified (up to 12)
357         // --- Sens - sensitive strips for the Slow Simulator !!!
358         // --- This does NOT work if all S or L-sectors are specified, i.e.
359         // --- if SecAL or SecAU < 0
360         //
361         //
362         //-----------------------------------------------------------------------------
363
364         //  gROOT->LoadMacro("SetTPCParam.C");
365         //  AliTPCParam *param = SetTPCParam();
366         AliTPC *TPC = new AliTPCv2("TPC", "Default");
367
368         // All sectors included 
369         TPC->SetSecAL(-1);
370         TPC->SetSecAU(-1);
371
372     }
373
374
375     if (iTOF) {
376         //=================== TOF parameters ============================
377         AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
378     }
379
380
381     if (iHMPID)
382     {
383         //=================== HMPID parameters ===========================
384         AliHMPID *HMPID = new AliHMPIDv1("HMPID", "normal HMPID");
385
386     }
387
388
389     if (iZDC)
390     {
391         //=================== ZDC parameters ============================
392
393         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
394     }
395
396     if (iTRD)
397     {
398         //=================== TRD parameters ============================
399
400         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
401
402         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
403         TRD->SetGasMix(1);
404         if (sgeo == kHoles) {
405             // With hole in front of PHOS
406             TRD->SetPHOShole();
407             // With hole in front of HMPID
408             TRD->SetHMPIDhole();
409         }
410             // Switch on TR
411             AliTRDsim *TRDsim = TRD->CreateTR();
412     }
413
414     if (iFMD)
415     {
416         //=================== FMD parameters ============================
417         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
418    }
419
420     if (iMUON)
421     {
422         //=================== MUON parameters ===========================
423
424         AliMUON *MUON = new AliMUONv1("MUON", "default");
425     }
426     //=================== PHOS parameters ===========================
427
428     if (iPHOS)
429     {
430         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
431     }
432
433
434     if (iPMD)
435     {
436         //=================== PMD parameters ============================
437         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
438     }
439
440     if (iT0)
441     {
442         //=================== T0 parameters ============================
443         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
444     }
445
446     if (iEMCAL)
447     {
448         //=================== EMCAL parameters ============================
449         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "G56_2_55_19");
450     }
451
452      if (iACORDE)
453     {
454         //=================== ACORDE parameters ============================
455         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
456     }
457
458      if (iVZERO)
459     {
460         //=================== ACORDE parameters ============================
461         AliVZERO *VZERO = new AliVZEROv3("VZERO", "normal VZERO");
462     }
463  
464              
465 }
466
467 Float_t EtaToTheta(Float_t arg){
468   return (180./TMath::Pi())*2.*atan(exp(-arg));
469 }
470
471
472
473 AliGenerator* GeneratorFactory(PprRun_t srun) {
474     Int_t isw = 3;
475     if (srad == kNoGluonRadiation) isw = 0;
476     
477
478     AliGenerator * gGener = 0x0;
479     switch (srun) {
480     case test50:
481       {
482         comment = comment.Append(":HIJINGparam test 50 particles");
483         AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
484         gener->SetMomentumRange(0, 999999.);
485         gener->SetPhiRange(0., 360.);
486         // Set pseudorapidity range from -8 to 8.
487         Float_t thmin = EtaToTheta(1);   // theta min. <---> eta max
488         Float_t thmax = EtaToTheta(-1);  // theta max. <---> eta min 
489         gener->SetThetaRange(thmin,thmax);
490         gGener=gener;
491       }
492       break;
493     case kParam_8000:
494       {
495         comment = comment.Append(":HIJINGparam N=8000");
496         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
497         gener->SetMomentumRange(0, 999999.);
498         gener->SetPhiRange(0., 360.);
499         // Set pseudorapidity range from -8 to 8.
500         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
501         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
502         gener->SetThetaRange(thmin,thmax);
503         gGener=gener;
504       }
505       break;
506     case kParam_4000:
507       {
508         comment = comment.Append("HIJINGparam N=4000");
509         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
510         gener->SetMomentumRange(0, 999999.);
511         gener->SetPhiRange(0., 360.);
512         // Set pseudorapidity range from -8 to 8.
513         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
514         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
515         gener->SetThetaRange(thmin,thmax);
516         gGener=gener;
517       }
518         break;
519     case kParam_2000:
520       {
521         comment = comment.Append("HIJINGparam N=2000");
522         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
523         gener->SetMomentumRange(0, 999999.);
524         gener->SetPhiRange(0., 360.);
525         // Set pseudorapidity range from -8 to 8.
526         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
527         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
528         gener->SetThetaRange(thmin,thmax);
529         gGener=gener;
530       }
531       break;
532 //
533 //  Hijing Central
534 //
535     case kHijing_cent1:
536       {
537         comment = comment.Append("HIJING cent1");
538         AliGenHijing *gener = HijingStandard();
539 // impact parameter range
540         gener->SetImpactParameterRange(0., 5.);
541         gGener=gener;
542       }
543       break;
544     case kHijing_cent2:
545       {
546         comment = comment.Append("HIJING cent2");
547         AliGenHijing *gener = HijingStandard();
548 // impact parameter range
549         gener->SetImpactParameterRange(0., 2.);
550         gGener=gener;
551       }
552       break;
553 //
554 // Hijing Peripheral 
555 //
556     case kHijing_per1:
557       {
558         comment = comment.Append("HIJING per1");
559         AliGenHijing *gener = HijingStandard();
560 // impact parameter range
561         gener->SetImpactParameterRange(5., 8.6);
562         gGener=gener;
563       }
564       break;
565     case kHijing_per2:
566       {
567         comment = comment.Append("HIJING per2");
568         AliGenHijing *gener = HijingStandard();
569 // impact parameter range
570         gener->SetImpactParameterRange(8.6, 11.2);
571         gGener=gener;
572       }
573       break;
574     case kHijing_per3:
575       {
576         comment = comment.Append("HIJING per3");
577         AliGenHijing *gener = HijingStandard();
578 // impact parameter range
579         gener->SetImpactParameterRange(11.2, 13.2);
580         gGener=gener;
581       }
582       break;
583     case kHijing_per4:
584       {
585         comment = comment.Append("HIJING per4");
586         AliGenHijing *gener = HijingStandard();
587 // impact parameter range
588         gener->SetImpactParameterRange(13.2, 15.);
589         gGener=gener;
590       }
591       break;
592     case kHijing_per5:
593       {
594         comment = comment.Append("HIJING per5");
595         AliGenHijing *gener = HijingStandard();
596 // impact parameter range
597         gener->SetImpactParameterRange(15., 100.);
598         gGener=gener;
599       }
600       break;
601 //
602 //  Jet-Jet
603 //
604     case kHijing_jj25:
605       {
606         comment = comment.Append("HIJING Jet 25 GeV");
607         AliGenHijing *gener = HijingStandard();
608 // impact parameter range
609         gener->SetImpactParameterRange(0., 5.);
610         // trigger
611         gener->SetTrigger(1);
612         gener->SetPtJet(25.);
613         gener->SetRadiation(isw);
614         gener->SetSimpleJets(!isw);
615         gener->SetJetEtaRange(-0.3,0.3);
616         gener->SetJetPhiRange(75., 165.);   
617         gGener=gener;
618       }
619       break;
620
621     case kHijing_jj50:
622       {
623         comment = comment.Append("HIJING Jet 50 GeV");
624         AliGenHijing *gener = HijingStandard();
625 // impact parameter range
626         gener->SetImpactParameterRange(0., 5.);
627         // trigger
628         gener->SetTrigger(1);
629         gener->SetPtJet(50.);
630         gener->SetRadiation(isw);
631         gener->SetSimpleJets(!isw);
632         gener->SetJetEtaRange(-0.3,0.3);
633         gener->SetJetPhiRange(75., 165.);   
634         gGener=gener;
635       }
636         break;
637
638     case kHijing_jj75:
639       {
640         comment = comment.Append("HIJING Jet 75 GeV");
641         AliGenHijing *gener = HijingStandard();
642 // impact parameter range
643         gener->SetImpactParameterRange(0., 5.);
644         // trigger
645         gener->SetTrigger(1);
646         gener->SetPtJet(75.);
647         gener->SetRadiation(isw);
648         gener->SetSimpleJets(!isw);
649         gener->SetJetEtaRange(-0.3,0.3);
650         gener->SetJetPhiRange(75., 165.);   
651         gGener=gener;
652       }
653       break;
654
655     case kHijing_jj100:
656       {
657         comment = comment.Append("HIJING Jet 100 GeV");
658         AliGenHijing *gener = HijingStandard();
659 // impact parameter range
660         gener->SetImpactParameterRange(0., 5.);
661         // trigger
662         gener->SetTrigger(1);
663         gener->SetPtJet(100.);
664         gener->SetRadiation(isw);
665         gener->SetSimpleJets(!isw);
666         gener->SetJetEtaRange(-0.3,0.3);
667         gener->SetJetPhiRange(75., 165.);   
668         gGener=gener;
669       }
670       break;
671
672     case kHijing_jj200:
673       {
674         comment = comment.Append("HIJING Jet 200 GeV");
675         AliGenHijing *gener = HijingStandard();
676 // impact parameter range
677         gener->SetImpactParameterRange(0., 5.);
678         // trigger
679         gener->SetTrigger(1);
680         gener->SetPtJet(200.);
681         gener->SetRadiation(isw);
682         gener->SetSimpleJets(!isw);
683         gener->SetJetEtaRange(-0.3,0.3);
684         gener->SetJetPhiRange(75., 165.);   
685         gGener=gener;
686       }
687       break;
688 //
689 // Gamma-Jet
690 //
691     case kHijing_gj25:
692       {
693         comment = comment.Append("HIJING Gamma 25 GeV");
694         AliGenHijing *gener = HijingStandard();
695 // impact parameter range
696         gener->SetImpactParameterRange(0., 5.);
697         // trigger
698         gener->SetTrigger(2);
699         gener->SetPtJet(25.);
700         gener->SetRadiation(isw);
701         gener->SetSimpleJets(!isw);
702         gener->SetJetEtaRange(-0.12, 0.12);
703         gener->SetJetPhiRange(220., 320.);
704         gGener=gener;
705       }
706       break;
707
708     case kHijing_gj50:
709       {
710         comment = comment.Append("HIJING Gamma 50 GeV");
711         AliGenHijing *gener = HijingStandard();
712 // impact parameter range
713         gener->SetImpactParameterRange(0., 5.);
714         // trigger
715         gener->SetTrigger(2);
716         gener->SetPtJet(50.);
717         gener->SetRadiation(isw);
718         gener->SetSimpleJets(!isw);
719         gener->SetJetEtaRange(-0.12, 0.12);
720         gener->SetJetPhiRange(220., 320.);
721         gGener=gener;
722       }
723       break;
724
725     case kHijing_gj75:
726       {
727         comment = comment.Append("HIJING Gamma 75 GeV");
728         AliGenHijing *gener = HijingStandard();
729 // impact parameter range
730         gener->SetImpactParameterRange(0., 5.);
731         // trigger
732         gener->SetTrigger(2);
733         gener->SetPtJet(75.);
734         gener->SetRadiation(isw);
735         gener->SetSimpleJets(!isw);
736         gener->SetJetEtaRange(-0.12, 0.12);
737         gener->SetJetPhiRange(220., 320.);
738         gGener=gener;
739       }
740       break;
741
742     case kHijing_gj100:
743       {
744         comment = comment.Append("HIJING Gamma 100 GeV");
745         AliGenHijing *gener = HijingStandard();
746 // impact parameter range
747         gener->SetImpactParameterRange(0., 5.);
748         // trigger
749         gener->SetTrigger(2);
750         gener->SetPtJet(100.);
751         gener->SetRadiation(isw);
752         gener->SetSimpleJets(!isw);
753         gener->SetJetEtaRange(-0.12, 0.12);
754         gener->SetJetPhiRange(220., 320.);
755         gGener=gener;
756       }
757       break;
758
759     case kHijing_gj200:
760       {
761         comment = comment.Append("HIJING Gamma 200 GeV");
762         AliGenHijing *gener = HijingStandard();
763 // impact parameter range
764         gener->SetImpactParameterRange(0., 5.);
765         // trigger
766         gener->SetTrigger(2);
767         gener->SetPtJet(200.);
768         gener->SetRadiation(isw);
769         gener->SetSimpleJets(!isw);
770         gener->SetJetEtaRange(-0.12, 0.12);
771         gener->SetJetPhiRange(220., 320.);
772         gGener=gener;
773       }
774       break;
775     case kHijing_pA:
776       {
777         comment = comment.Append("HIJING pA");
778
779         AliGenCocktail *gener  = new AliGenCocktail();
780 //      gener->SetTrackingFlag(0);
781         gener->SetOrigin(0, 0, 0);    // vertex position
782         gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
783         gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
784         gener->SetVertexSmear(kPerEvent); 
785         AliGenHijing   *hijing = new AliGenHijing(-1);
786 // centre of mass energy 
787         hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
788 // impact parameter range
789         hijing->SetImpactParameterRange(0., 15.);
790 // reference frame
791         hijing->SetReferenceFrame("CMS");
792         hijing->SetBoostLHC(1);
793 // projectile
794         hijing->SetProjectile("P", 1, 1);
795         hijing->SetTarget    ("A", 208, 82);
796 // tell hijing to keep the full parent child chain
797         hijing->KeepFullEvent();
798 // enable jet quenching
799         hijing->SetJetQuenching(0);
800 // enable shadowing
801         hijing->SetShadowing(1);
802 // Don't track spectators
803         hijing->SetSpectators(0);
804 // kinematic selection
805         hijing->SetSelectAll(0);
806 //
807         AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
808         AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
809         gray->SetSlowNucleonModel(model);
810         gray->SetDebug(1);
811         gener->AddGenerator(hijing,"Hijing pPb", 1);
812         gener->AddGenerator(gray,  "Gray Particles",1);
813         gGener=gener;
814       }
815       break;
816     case kPythia6:
817       {
818         comment = comment.Append(":Pythia p-p @ 14 TeV");
819         AliGenPythia *gener = new AliGenPythia(-1); 
820         gener->SetMomentumRange(0,999999);
821         gener->SetPhiRange(-180,180);
822         gener->SetThetaRange(0., 180.);
823         gener->SetYRange(-12,12);
824         gener->SetPtRange(0,1000);
825         gener->SetStrucFunc(kCTEQ4L);   
826         gener->SetProcess(kPyMb);
827         gener->SetEnergyCMS(14000.);
828         gGener=gener;
829       }
830       break;
831     case kPythia6Jets:
832       {
833         comment = comment.Append(":Pythia jets @ 5.5 TeV");
834         AliGenPythia * gener = new AliGenPythia(-1);
835         gener->SetEnergyCMS(5500.);//        Centre of mass energy
836         gener->SetProcess(kPyJets);//        Process type
837         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
838         gener->SetJetPhiRange(0., 360.);
839         gener->SetJetEtRange(10., 1000.);
840         gener->SetGluonRadiation(1,1);
841         //    gener->SetPtKick(0.);
842         //   Structure function
843         gener->SetStrucFunc(kCTEQ4L);
844         gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
845         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
846         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
847         gGener=gener;
848       }
849       break;
850     case kD0PbPb5500:
851       {
852         comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
853         AliGenPythia * gener = new AliGenPythia(10);
854         gener->SetProcess(kPyD0PbPbMNR);
855         gener->SetStrucFunc(kCTEQ4L);
856         gener->SetPtHard(2.1,-1.0);
857         gener->SetEnergyCMS(5500.);
858         gener->SetNuclei(208,208);
859         gGener=gener;
860       }
861       break;
862     }
863     return gGener;
864 }
865
866 AliGenHijing* HijingStandard()
867 {
868     AliGenHijing *gener = new AliGenHijing(-1);
869 // centre of mass energy 
870     gener->SetEnergyCMS(5500.);
871 // reference frame
872     gener->SetReferenceFrame("CMS");
873 // projectile
874      gener->SetProjectile("A", 208, 82);
875      gener->SetTarget    ("A", 208, 82);
876 // tell hijing to keep the full parent child chain
877      gener->KeepFullEvent();
878 // enable jet quenching
879      gener->SetJetQuenching(4);
880 // enable shadowing
881      gener->SetShadowing(1);
882 // neutral pion and heavy particle decays switched off
883      gener->SetDecaysOff(1);
884 // Don't track spectators
885      gener->SetSpectators(0);
886 // kinematic selection
887      gener->SetSelectAll(0);
888      return gener;
889 }
890
891
892