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