]> git.uio.no Git - u/mrichter/AliRoot.git/blob - macros/ConfigPPR.C
NewIO Implementation
[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/AliITSvPPRasymm.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         AliITSvPPRasymm *ITS  = new AliITSvPPRasymm("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(kFALSE);                                         // 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     //AliITSvPPRsymm *ITS  = new AliITSvPPRsymm("ITS","New ITS PPR detailed version with symmetric services");
323     //ITS->SetMinorVersion(2);                                       // don't touch this parameter if you're not an ITS developer
324     //ITS->SetReadDet(kFALSE);                                       // don't touch this parameter if you're not an ITS developer
325     //ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRsymm2.det"); // don't touch this parameter if you're not an ITS developer
326     //ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
327     //ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
328     //ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
329     //ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
330     //ITS->SetRails(0);              // 1 --> rails in ; 0 --> rails out
331     //ITS->SetCoolingFluid(1);       // 1 --> water ; 0 --> freon
332     //
333     //
334     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
335     // for reconstruction !):
336     //                                                     
337     //
338     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
339     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
340     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
341     //
342     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
343     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
344     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
345     //                      
346     //
347     //
348     // Geant3 <-> EUCLID conversion
349     // ============================
350     //
351     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
352     // media to two ASCII files (called by default ITSgeometry.euc and
353     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
354     // The default (=0) means that you dont want to use this facility.
355     //
356         ITS->SetEUCLID(0);  
357     }
358
359     if (iTPC)
360     {
361         //============================ TPC parameters ================================
362         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
363         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
364         // --- sectors are specified, any value other than that requires at least one 
365         // --- sector (lower or upper)to be specified!
366         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
367         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
368         // --- SecLows - number of lower sectors specified (up to 6)
369         // --- SecUps - number of upper sectors specified (up to 12)
370         // --- Sens - sensitive strips for the Slow Simulator !!!
371         // --- This does NOT work if all S or L-sectors are specified, i.e.
372         // --- if SecAL or SecAU < 0
373         //
374         //
375         //-----------------------------------------------------------------------------
376
377         //  gROOT->LoadMacro("SetTPCParam.C");
378         //  AliTPCParam *param = SetTPCParam();
379         AliTPC *TPC = new AliTPCv2("TPC", "Default");
380
381         // All sectors included 
382         TPC->SetSecAL(-1);
383         TPC->SetSecAU(-1);
384
385     }
386
387
388     if (iTOF) {
389         if (sgeo == kHoles) {
390         //=================== TOF parameters ============================
391             AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
392         } else {
393             AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
394         }
395     }
396
397
398     if (iRICH)
399     {
400         //=================== RICH parameters ===========================
401         AliRICH *RICH = new AliRICHv3("RICH", "normal RICH");
402
403     }
404
405
406     if (iZDC)
407     {
408         //=================== ZDC parameters ============================
409
410         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
411     }
412
413     if (iTRD)
414     {
415         //=================== TRD parameters ============================
416
417         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
418
419         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
420         TRD->SetGasMix(1);
421         if (sgeo == kHoles) {
422             // With hole in front of PHOS
423             TRD->SetPHOShole();
424             // With hole in front of RICH
425             TRD->SetRICHhole();
426         }
427             // Switch on TR
428             AliTRDsim *TRDsim = TRD->CreateTR();
429     }
430
431     if (iFMD)
432     {
433         //=================== FMD parameters ============================
434         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
435         FMD->SetRingsSi1(256);
436         FMD->SetRingsSi2(128);
437         FMD->SetSectorsSi1(20);
438         FMD->SetSectorsSi2(40);      
439    }
440
441     if (iMUON)
442     {
443         //=================== MUON parameters ===========================
444
445         AliMUON *MUON = new AliMUONv1("MUON", "default");
446     }
447     //=================== PHOS parameters ===========================
448
449     if (iPHOS)
450     {
451         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
452     }
453
454
455     if (iPMD)
456     {
457         //=================== PMD parameters ============================
458         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
459     }
460
461     if (iSTART)
462     {
463         //=================== START parameters ============================
464         AliSTART *START = new AliSTARTv1("START", "START Detector");
465     }
466
467     if (iEMCAL)
468     {
469         //=================== EMCAL parameters ============================
470         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "G56_2_55_19");
471     }
472
473      if (iCRT)
474     {
475         //=================== CRT parameters ============================
476         AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
477     }
478
479      if (iVZERO)
480     {
481         //=================== CRT parameters ============================
482         AliVZERO *VZERO = new AliVZEROv2("VZERO", "normal VZERO");
483     }
484  
485              
486 }
487
488 Float_t EtaToTheta(Float_t arg){
489   return (180./TMath::Pi())*2.*atan(exp(-arg));
490 }
491
492
493
494 AliGenerator* GeneratorFactory(PprRun_t srun) {
495     Int_t isw = 3;
496     if (srad == kNoGluonRadiation) isw = 0;
497     
498
499     AliGenerator * gGener = 0x0;
500     switch (srun) {
501     case test50:
502       {
503         comment = comment.Append(":HIJINGparam test 50 particles");
504         AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
505         gener->SetMomentumRange(0, 999999.);
506         gener->SetPhiRange(0., 360.);
507         // Set pseudorapidity range from -8 to 8.
508         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
509         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
510         gener->SetThetaRange(thmin,thmax);
511         gGener=gener;
512       }
513         break;
514     case kParam_8000:
515       {
516         comment = comment.Append(":HIJINGparam N=8000");
517         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
518         gener->SetMomentumRange(0, 999999.);
519         gener->SetPhiRange(0., 360.);
520         // Set pseudorapidity range from -8 to 8.
521         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
522         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
523         gener->SetThetaRange(thmin,thmax);
524         gGener=gener;
525       }
526         break;
527     case kParam_4000:
528       {
529         comment = comment.Append("HIJINGparam N=4000");
530         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
531         gener->SetMomentumRange(0, 999999.);
532         gener->SetPhiRange(0., 360.);
533         // Set pseudorapidity range from -8 to 8.
534         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
535         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
536         gener->SetThetaRange(thmin,thmax);
537         gGener=gener;
538       }
539         break;
540     case kParam_2000:
541       {
542         comment = comment.Append("HIJINGparam N=2000");
543         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
544         gener->SetMomentumRange(0, 999999.);
545         gener->SetPhiRange(0., 360.);
546         // Set pseudorapidity range from -8 to 8.
547         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
548         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
549         gener->SetThetaRange(thmin,thmax);
550         gGener=gener;
551       }
552         break;
553 //
554 //  Hijing Central
555 //
556     case kHijing_cent1:
557       {
558         comment = comment.Append("HIJING cent1");
559         AliGenHijing *gener = HijingStandard();
560 // impact parameter range
561         gener->SetImpactParameterRange(0., 5.);
562         gGener=gener;
563       }
564         break;
565     case kHijing_cent2:
566       {
567         comment = comment.Append("HIJING cent2");
568         AliGenHijing *gener = HijingStandard();
569 // impact parameter range
570         gener->SetImpactParameterRange(0., 2.);
571         gGener=gener;
572         break;
573       }
574 //
575 // Hijing Peripheral 
576 //
577     case kHijing_per1:
578       {
579         comment = comment.Append("HIJING per1");
580         AliGenHijing *gener = HijingStandard();
581 // impact parameter range
582         gener->SetImpactParameterRange(5., 8.6);
583         gGener=gener;
584       }
585         break;
586     case kHijing_per2:
587       {
588         comment = comment.Append("HIJING per2");
589         AliGenHijing *gener = HijingStandard();
590 // impact parameter range
591         gener->SetImpactParameterRange(8.6, 11.2);
592         gGener=gener;
593       }
594         break;
595     case kHijing_per3:
596       {
597         comment = comment.Append("HIJING per3");
598         AliGenHijing *gener = HijingStandard();
599 // impact parameter range
600         gener->SetImpactParameterRange(11.2, 13.2);
601         gGener=gener;
602       }
603         break;
604     case kHijing_per4:
605       {
606         comment = comment.Append("HIJING per4");
607         AliGenHijing *gener = HijingStandard();
608 // impact parameter range
609         gener->SetImpactParameterRange(13.2, 15.);
610         gGener=gener;
611       }
612         break;
613     case kHijing_per5:
614       {
615         comment = comment.Append("HIJING per5");
616         AliGenHijing *gener = HijingStandard();
617 // impact parameter range
618         gener->SetImpactParameterRange(15., 100.);
619         gGener=gener;
620       }
621         break;
622 //
623 //  Jet-Jet
624 //
625     case kHijing_jj25:
626       {
627         comment = comment.Append("HIJING Jet 25 GeV");
628         AliGenHijing *gener = HijingStandard();
629 // impact parameter range
630         gener->SetImpactParameterRange(0., 5.);
631         // trigger
632         gener->SetTrigger(1);
633         gener->SetPtJet(25.);
634         gener->SetRadiation(isw);
635         gener->SetSimpleJets(!isw);
636         gener->SetJetEtaRange(-0.3,0.3);
637         gener->SetJetPhiRange(75., 165.);   
638         gGener=gener;
639       }
640         break;
641
642     case kHijing_jj50:
643       {
644         comment = comment.Append("HIJING Jet 50 GeV");
645         AliGenHijing *gener = HijingStandard();
646 // impact parameter range
647         gener->SetImpactParameterRange(0., 5.);
648         // trigger
649         gener->SetTrigger(1);
650         gener->SetPtJet(50.);
651         gener->SetRadiation(isw);
652         gener->SetSimpleJets(!isw);
653         gener->SetJetEtaRange(-0.3,0.3);
654         gener->SetJetPhiRange(75., 165.);   
655         gGener=gener;
656       }
657         break;
658
659     case kHijing_jj75:
660       {
661         comment = comment.Append("HIJING Jet 75 GeV");
662         AliGenHijing *gener = HijingStandard();
663 // impact parameter range
664         gener->SetImpactParameterRange(0., 5.);
665         // trigger
666         gener->SetTrigger(1);
667         gener->SetPtJet(75.);
668         gener->SetRadiation(isw);
669         gener->SetSimpleJets(!isw);
670         gener->SetJetEtaRange(-0.3,0.3);
671         gener->SetJetPhiRange(75., 165.);   
672         gGener=gener;
673       }
674         break;
675
676     case kHijing_jj100:
677       {
678         comment = comment.Append("HIJING Jet 100 GeV");
679         AliGenHijing *gener = HijingStandard();
680 // impact parameter range
681         gener->SetImpactParameterRange(0., 5.);
682         // trigger
683         gener->SetTrigger(1);
684         gener->SetPtJet(100.);
685         gener->SetRadiation(isw);
686         gener->SetSimpleJets(!isw);
687         gener->SetJetEtaRange(-0.3,0.3);
688         gener->SetJetPhiRange(75., 165.);   
689         gGener=gener;
690       }
691         break;
692
693     case kHijing_jj200:
694       {
695         comment = comment.Append("HIJING Jet 200 GeV");
696         AliGenHijing *gener = HijingStandard();
697 // impact parameter range
698         gener->SetImpactParameterRange(0., 5.);
699         // trigger
700         gener->SetTrigger(1);
701         gener->SetPtJet(200.);
702         gener->SetRadiation(isw);
703         gener->SetSimpleJets(!isw);
704         gener->SetJetEtaRange(-0.3,0.3);
705         gener->SetJetPhiRange(75., 165.);   
706         gGener=gener;
707       }
708         break;
709 //
710 // Gamma-Jet
711 //
712     case kHijing_gj25:
713       {
714         comment = comment.Append("HIJING Gamma 25 GeV");
715         AliGenHijing *gener = HijingStandard();
716 // impact parameter range
717         gener->SetImpactParameterRange(0., 5.);
718         // trigger
719         gener->SetTrigger(2);
720         gener->SetPtJet(25.);
721         gener->SetRadiation(isw);
722         gener->SetSimpleJets(!isw);
723         gener->SetJetEtaRange(-0.12, 0.12);
724         gener->SetJetPhiRange(220., 320.);
725         gGener=gener;
726       }
727         break;
728
729     case kHijing_gj50:
730       {
731         comment = comment.Append("HIJING Gamma 50 GeV");
732         AliGenHijing *gener = HijingStandard();
733 // impact parameter range
734         gener->SetImpactParameterRange(0., 5.);
735         // trigger
736         gener->SetTrigger(2);
737         gener->SetPtJet(50.);
738         gener->SetRadiation(isw);
739         gener->SetSimpleJets(!isw);
740         gener->SetJetEtaRange(-0.12, 0.12);
741         gener->SetJetPhiRange(220., 320.);
742         gGener=gener;
743       }
744         break;
745
746     case kHijing_gj75:
747       {
748         comment = comment.Append("HIJING Gamma 75 GeV");
749         AliGenHijing *gener = HijingStandard();
750 // impact parameter range
751         gener->SetImpactParameterRange(0., 5.);
752         // trigger
753         gener->SetTrigger(2);
754         gener->SetPtJet(75.);
755         gener->SetRadiation(isw);
756         gener->SetSimpleJets(!isw);
757         gener->SetJetEtaRange(-0.12, 0.12);
758         gener->SetJetPhiRange(220., 320.);
759         gGener=gener;
760       }
761         break;
762
763     case kHijing_gj100:
764       {
765         comment = comment.Append("HIJING Gamma 100 GeV");
766         AliGenHijing *gener = HijingStandard();
767 // impact parameter range
768         gener->SetImpactParameterRange(0., 5.);
769         // trigger
770         gener->SetTrigger(2);
771         gener->SetPtJet(100.);
772         gener->SetRadiation(isw);
773         gener->SetSimpleJets(!isw);
774         gener->SetJetEtaRange(-0.12, 0.12);
775         gener->SetJetPhiRange(220., 320.);
776         gGener=gener;
777       }
778         break;
779
780     case kHijing_gj200:
781       {
782         comment = comment.Append("HIJING Gamma 200 GeV");
783         AliGenHijing *gener = HijingStandard();
784 // impact parameter range
785         gener->SetImpactParameterRange(0., 5.);
786         // trigger
787         gener->SetTrigger(2);
788         gener->SetPtJet(200.);
789         gener->SetRadiation(isw);
790         gener->SetSimpleJets(!isw);
791         gener->SetJetEtaRange(-0.12, 0.12);
792         gener->SetJetPhiRange(220., 320.);
793         gGener=gener;
794       }
795         break;
796     case kHijing_pA:
797       {
798         comment = comment.Append("HIJING pA");
799
800         AliGenCocktail *gener  = new AliGenCocktail();
801 //      gener->SetTrackingFlag(0);
802         gener->SetOrigin(0, 0, 0);    // vertex position
803         gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
804         gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
805         gener->SetVertexSmear(kPerEvent); 
806         AliGenHijing   *hijing = new AliGenHijing(-1);
807 // centre of mass energy 
808         hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
809 // impact parameter range
810         hijing->SetImpactParameterRange(0., 15.);
811 // reference frame
812         hijing->SetReferenceFrame("CMS");
813         hijing->SetBoostLHC(1);
814 // projectile
815         hijing->SetProjectile("P", 1, 1);
816         hijing->SetTarget    ("A", 208, 82);
817 // tell hijing to keep the full parent child chain
818         hijing->KeepFullEvent();
819 // enable jet quenching
820         hijing->SetJetQuenching(0);
821 // enable shadowing
822         hijing->SetShadowing(1);
823 // Don't track spectators
824         hijing->SetSpectators(0);
825 // kinematic selection
826         hijing->SetSelectAll(0);
827 //
828         AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
829         AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
830         gray->SetSlowNucleonModel(model);
831         gray->SetDebug(1);
832         gener->AddGenerator(hijing,"Hijing pPb", 1);
833         gener->AddGenerator(gray,  "Gray Particles",1);
834         gGener=gener;
835       }
836         break;
837     case kPythia:
838       {
839         comment = comment.Append(":Pythia p-p @ 14 TeV");
840         AliGenPythia *gener = new AliGenPythia(-1); 
841         gener->SetMomentumRange(0,999999);
842         gener->SetPhiRange(-180,180);
843         gener->SetThetaRange(0., 180.);
844         gener->SetYRange(-12,12);
845         gener->SetPtRange(0,1000);
846         gener->SetStrucFunc(kCTEQ4L);   
847         gener->SetProcess(kPyMb);
848         gener->SetEnergyCMS(14000.);
849         gGener=gener;
850       }
851     break;
852     }
853     return gGener;
854 }
855
856 AliGenHijing* HijingStandard()
857 {
858     AliGenHijing *gener = new AliGenHijing(-1);
859 // centre of mass energy 
860     gener->SetEnergyCMS(5500.);
861 // reference frame
862     gener->SetReferenceFrame("CMS");
863 // projectile
864      gener->SetProjectile("A", 208, 82);
865      gener->SetTarget    ("A", 208, 82);
866 // tell hijing to keep the full parent child chain
867      gener->KeepFullEvent();
868 // enable jet quenching
869      gener->SetJetQuenching(4);
870 // enable shadowing
871      gener->SetShadowing(1);
872 // neutral pion and heavy particle decays switched off
873      gener->SetDecaysOff(1);
874 // Don't track spectators
875      gener->SetSpectators(0);
876 // kinematic selection
877      gener->SetSelectAll(0);
878      return gener;
879 }
880
881
882