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