]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFconfig.C
Adding the track fit residuals as a consequence of the ExB distortions (Marian)
[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/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/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/AliACORDEv0.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 enum PprMag_t
65 {
66     k2kG, k4kG, k5kG
67 };
68
69
70 // This part for configuration    
71 static PprRun_t srun = test50;
72 //static PprRun_t srun = kPythia6;
73 static PprGeo_t sgeo = kNoHoles;
74 static PprRad_t srad = kGluonRadiation;
75 static PprMag_t smag = k5kG;
76
77 // Comment line 
78 static TString  comment;
79
80 // Functions
81 Float_t EtaToTheta(Float_t arg);
82 AliGenerator* GeneratorFactory(PprRun_t srun);
83 AliGenHijing* HijingStandard();
84
85 void Config()
86 {
87     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
88     // Theta range given through pseudorapidity limits 22/6/2001
89
90     // Set Random Number seed
91   gRandom->SetSeed(12345); //Set 0 to use the current time
92     cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
93
94
95    // libraries required by geant321
96 #if defined(__CINT__)
97     gSystem->Load("libgeant321");
98 #endif
99
100     new     TGeant3("C++ Interface to Geant3");
101
102     AliRunLoader* rl=0x0;
103
104     cout<<"Config.C: Creating Run Loader ..."<<endl;
105     rl = AliRunLoader::Open("galice.root",
106                             AliConfig::GetDefaultEventFolderName(),
107                             "recreate");
108     if (rl == 0x0)
109       {
110         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
111         return;
112       }
113     rl->SetCompressionLevel(2);
114     rl->SetNumberOfEventsPerFile(3);
115     gAlice->SetRunLoader(rl);
116
117     //
118     // Set External decayer
119     AliDecayer *decayer = new AliDecayerPythia();
120
121     decayer->SetForceDecay(kAll);
122     decayer->Init();
123     gMC->SetExternalDecayer(decayer);
124     //
125     //
126     //=======================================================================
127     //
128     //=======================================================================
129     // ************* STEERING parameters FOR ALICE SIMULATION **************
130     // --- Specify event type to be tracked through the ALICE setup
131     // --- All positions are in cm, angles in degrees, and P and E in GeV
132
133     gMC->SetProcess("DCAY",1);
134     gMC->SetProcess("PAIR",1);
135     gMC->SetProcess("COMP",1);
136     gMC->SetProcess("PHOT",1);
137     gMC->SetProcess("PFIS",0);
138     gMC->SetProcess("DRAY",0);
139     gMC->SetProcess("ANNI",1);
140     gMC->SetProcess("BREM",1);
141     gMC->SetProcess("MUNU",1);
142     gMC->SetProcess("CKOV",1);
143     gMC->SetProcess("HADR",1);
144     gMC->SetProcess("LOSS",2);
145     gMC->SetProcess("MULS",1);
146     gMC->SetProcess("RAYL",1);
147
148     Float_t cut = 1.e-3;        // 1MeV cut by default
149     Float_t tofmax = 1.e10;
150
151     gMC->SetCut("CUTGAM", cut);
152     gMC->SetCut("CUTELE", cut);
153     gMC->SetCut("CUTNEU", cut);
154     gMC->SetCut("CUTHAD", cut);
155     gMC->SetCut("CUTMUO", cut);
156     gMC->SetCut("BCUTE",  cut); 
157     gMC->SetCut("BCUTM",  cut); 
158     gMC->SetCut("DCUTE",  cut); 
159     gMC->SetCut("DCUTM",  cut); 
160     gMC->SetCut("PPCUTM", cut);
161     gMC->SetCut("TOFMAX", tofmax); 
162
163
164 // Generator Configuration
165     gAlice->SetDebug(1);
166     AliGenerator* gener = GeneratorFactory(srun);
167     gener->SetOrigin(0, 0, 0);    // vertex position
168     gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
169     gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
170     gener->SetVertexSmear(kPerEvent); 
171     gener->SetTrackingFlag(1);
172     gener->Init();
173     
174     if (smag == k2kG) {
175         comment = comment.Append(" | L3 field 0.2 T");
176     } else if (smag == k4kG) {
177         comment = comment.Append(" | L3 field 0.4 T");
178     } else if (smag == k5kG) {
179         comment = comment.Append(" | L3 field 0.5 T");
180     }
181     
182     
183     if (srad == kGluonRadiation)
184     {
185         comment = comment.Append(" | Gluon Radiation On");
186         
187     } else {
188         comment = comment.Append(" | Gluon Radiation Off");
189     }
190
191     if (sgeo == kHoles)
192     {
193         comment = comment.Append(" | Holes for PHOS/HMPID");
194         
195     } else {
196         comment = comment.Append(" | No holes for PHOS/HMPID");
197     }
198
199     printf("\n \n Comment: %s \n \n", comment.Data());
200     
201     
202 // Field (L3 0.4 T)
203     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., smag);
204     field->SetL3ConstField(0); //Using const. field in the barrel
205     rl->CdGAFile();
206     gAlice->SetField(field);    
207 //
208     Int_t   iABSO   = 0;
209     Int_t   iDIPO   = 0;
210     Int_t   iFMD    = 0;
211     Int_t   iFRAME  = 1;
212     Int_t   iHALL   = 0;
213     Int_t   iITS    = 1;
214     Int_t   iMAG    = 1;
215     Int_t   iMUON   = 0;
216     Int_t   iPHOS   = 0;
217     Int_t   iPIPE   = 1;
218     Int_t   iPMD    = 0;
219     Int_t   iHMPID   = 0;
220     Int_t   iSHIL   = 0;
221     Int_t   iT0  = 0;
222     Int_t   iTOF    = 1;
223     Int_t   iTPC    = 1;
224     Int_t   iTRD    = 1;
225     Int_t   iZDC    = 0;
226     Int_t   iEMCAL  = 0;
227     Int_t   iVZERO  = 0;
228     Int_t   iACORDE    = 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         //=================== TOF parameters ============================
378         AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
379     }
380
381
382     if (iHMPID)
383     {
384         //=================== HMPID parameters ===========================
385         AliHMPID *HMPID = new AliHMPIDv1("HMPID", "normal HMPID");
386
387     }
388
389
390     if (iZDC)
391     {
392         //=================== ZDC parameters ============================
393
394         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
395     }
396
397     if (iTRD)
398     {
399         //=================== TRD parameters ============================
400
401         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
402
403         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
404         TRD->SetGasMix(1);
405         if (sgeo == kHoles) {
406             // With hole in front of PHOS
407             TRD->SetPHOShole();
408             // With hole in front of HMPID
409             TRD->SetHMPIDhole();
410         }
411             // Switch on TR
412             AliTRDsim *TRDsim = TRD->CreateTR();
413     }
414
415     if (iFMD)
416     {
417         //=================== FMD parameters ============================
418         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
419    }
420
421     if (iMUON)
422     {
423         //=================== MUON parameters ===========================
424
425         AliMUON *MUON = new AliMUONv1("MUON", "default");
426     }
427     //=================== PHOS parameters ===========================
428
429     if (iPHOS)
430     {
431         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
432     }
433
434
435     if (iPMD)
436     {
437         //=================== PMD parameters ============================
438         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
439     }
440
441     if (iT0)
442     {
443         //=================== T0 parameters ============================
444         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
445     }
446
447     if (iEMCAL)
448     {
449         //=================== EMCAL parameters ============================
450         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "G56_2_55_19");
451     }
452
453      if (iACORDE)
454     {
455         //=================== ACORDE parameters ============================
456         AliACORDE *ACORDE = new AliACORDEv0("ACORDE", "normal ACORDE");
457     }
458
459      if (iVZERO)
460     {
461         //=================== ACORDE parameters ============================
462         AliVZERO *VZERO = new AliVZEROv3("VZERO", "normal VZERO");
463     }
464  
465              
466 }
467
468 Float_t EtaToTheta(Float_t arg){
469   return (180./TMath::Pi())*2.*atan(exp(-arg));
470 }
471
472
473
474 AliGenerator* GeneratorFactory(PprRun_t srun) {
475     Int_t isw = 3;
476     if (srad == kNoGluonRadiation) isw = 0;
477     
478
479     AliGenerator * gGener = 0x0;
480     switch (srun) {
481     case test50:
482       {
483         comment = comment.Append(":HIJINGparam test 50 particles");
484         AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
485         gener->SetMomentumRange(0, 999999.);
486         gener->SetPhiRange(0., 360.);
487         // Set pseudorapidity range from -8 to 8.
488         Float_t thmin = EtaToTheta(1);   // theta min. <---> eta max
489         Float_t thmax = EtaToTheta(-1);  // theta max. <---> eta min 
490         gener->SetThetaRange(thmin,thmax);
491         gGener=gener;
492       }
493       break;
494     case kParam_8000:
495       {
496         comment = comment.Append(":HIJINGparam N=8000");
497         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
498         gener->SetMomentumRange(0, 999999.);
499         gener->SetPhiRange(0., 360.);
500         // Set pseudorapidity range from -8 to 8.
501         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
502         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
503         gener->SetThetaRange(thmin,thmax);
504         gGener=gener;
505       }
506       break;
507     case kParam_4000:
508       {
509         comment = comment.Append("HIJINGparam N=4000");
510         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
511         gener->SetMomentumRange(0, 999999.);
512         gener->SetPhiRange(0., 360.);
513         // Set pseudorapidity range from -8 to 8.
514         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
515         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
516         gener->SetThetaRange(thmin,thmax);
517         gGener=gener;
518       }
519         break;
520     case kParam_2000:
521       {
522         comment = comment.Append("HIJINGparam N=2000");
523         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
524         gener->SetMomentumRange(0, 999999.);
525         gener->SetPhiRange(0., 360.);
526         // Set pseudorapidity range from -8 to 8.
527         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
528         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
529         gener->SetThetaRange(thmin,thmax);
530         gGener=gener;
531       }
532       break;
533 //
534 //  Hijing Central
535 //
536     case kHijing_cent1:
537       {
538         comment = comment.Append("HIJING cent1");
539         AliGenHijing *gener = HijingStandard();
540 // impact parameter range
541         gener->SetImpactParameterRange(0., 5.);
542         gGener=gener;
543       }
544       break;
545     case kHijing_cent2:
546       {
547         comment = comment.Append("HIJING cent2");
548         AliGenHijing *gener = HijingStandard();
549 // impact parameter range
550         gener->SetImpactParameterRange(0., 2.);
551         gGener=gener;
552       }
553       break;
554 //
555 // Hijing Peripheral 
556 //
557     case kHijing_per1:
558       {
559         comment = comment.Append("HIJING per1");
560         AliGenHijing *gener = HijingStandard();
561 // impact parameter range
562         gener->SetImpactParameterRange(5., 8.6);
563         gGener=gener;
564       }
565       break;
566     case kHijing_per2:
567       {
568         comment = comment.Append("HIJING per2");
569         AliGenHijing *gener = HijingStandard();
570 // impact parameter range
571         gener->SetImpactParameterRange(8.6, 11.2);
572         gGener=gener;
573       }
574       break;
575     case kHijing_per3:
576       {
577         comment = comment.Append("HIJING per3");
578         AliGenHijing *gener = HijingStandard();
579 // impact parameter range
580         gener->SetImpactParameterRange(11.2, 13.2);
581         gGener=gener;
582       }
583       break;
584     case kHijing_per4:
585       {
586         comment = comment.Append("HIJING per4");
587         AliGenHijing *gener = HijingStandard();
588 // impact parameter range
589         gener->SetImpactParameterRange(13.2, 15.);
590         gGener=gener;
591       }
592       break;
593     case kHijing_per5:
594       {
595         comment = comment.Append("HIJING per5");
596         AliGenHijing *gener = HijingStandard();
597 // impact parameter range
598         gener->SetImpactParameterRange(15., 100.);
599         gGener=gener;
600       }
601       break;
602 //
603 //  Jet-Jet
604 //
605     case kHijing_jj25:
606       {
607         comment = comment.Append("HIJING Jet 25 GeV");
608         AliGenHijing *gener = HijingStandard();
609 // impact parameter range
610         gener->SetImpactParameterRange(0., 5.);
611         // trigger
612         gener->SetTrigger(1);
613         gener->SetPtJet(25.);
614         gener->SetRadiation(isw);
615         gener->SetSimpleJets(!isw);
616         gener->SetJetEtaRange(-0.3,0.3);
617         gener->SetJetPhiRange(75., 165.);   
618         gGener=gener;
619       }
620       break;
621
622     case kHijing_jj50:
623       {
624         comment = comment.Append("HIJING Jet 50 GeV");
625         AliGenHijing *gener = HijingStandard();
626 // impact parameter range
627         gener->SetImpactParameterRange(0., 5.);
628         // trigger
629         gener->SetTrigger(1);
630         gener->SetPtJet(50.);
631         gener->SetRadiation(isw);
632         gener->SetSimpleJets(!isw);
633         gener->SetJetEtaRange(-0.3,0.3);
634         gener->SetJetPhiRange(75., 165.);   
635         gGener=gener;
636       }
637         break;
638
639     case kHijing_jj75:
640       {
641         comment = comment.Append("HIJING Jet 75 GeV");
642         AliGenHijing *gener = HijingStandard();
643 // impact parameter range
644         gener->SetImpactParameterRange(0., 5.);
645         // trigger
646         gener->SetTrigger(1);
647         gener->SetPtJet(75.);
648         gener->SetRadiation(isw);
649         gener->SetSimpleJets(!isw);
650         gener->SetJetEtaRange(-0.3,0.3);
651         gener->SetJetPhiRange(75., 165.);   
652         gGener=gener;
653       }
654       break;
655
656     case kHijing_jj100:
657       {
658         comment = comment.Append("HIJING Jet 100 GeV");
659         AliGenHijing *gener = HijingStandard();
660 // impact parameter range
661         gener->SetImpactParameterRange(0., 5.);
662         // trigger
663         gener->SetTrigger(1);
664         gener->SetPtJet(100.);
665         gener->SetRadiation(isw);
666         gener->SetSimpleJets(!isw);
667         gener->SetJetEtaRange(-0.3,0.3);
668         gener->SetJetPhiRange(75., 165.);   
669         gGener=gener;
670       }
671       break;
672
673     case kHijing_jj200:
674       {
675         comment = comment.Append("HIJING Jet 200 GeV");
676         AliGenHijing *gener = HijingStandard();
677 // impact parameter range
678         gener->SetImpactParameterRange(0., 5.);
679         // trigger
680         gener->SetTrigger(1);
681         gener->SetPtJet(200.);
682         gener->SetRadiation(isw);
683         gener->SetSimpleJets(!isw);
684         gener->SetJetEtaRange(-0.3,0.3);
685         gener->SetJetPhiRange(75., 165.);   
686         gGener=gener;
687       }
688       break;
689 //
690 // Gamma-Jet
691 //
692     case kHijing_gj25:
693       {
694         comment = comment.Append("HIJING Gamma 25 GeV");
695         AliGenHijing *gener = HijingStandard();
696 // impact parameter range
697         gener->SetImpactParameterRange(0., 5.);
698         // trigger
699         gener->SetTrigger(2);
700         gener->SetPtJet(25.);
701         gener->SetRadiation(isw);
702         gener->SetSimpleJets(!isw);
703         gener->SetJetEtaRange(-0.12, 0.12);
704         gener->SetJetPhiRange(220., 320.);
705         gGener=gener;
706       }
707       break;
708
709     case kHijing_gj50:
710       {
711         comment = comment.Append("HIJING Gamma 50 GeV");
712         AliGenHijing *gener = HijingStandard();
713 // impact parameter range
714         gener->SetImpactParameterRange(0., 5.);
715         // trigger
716         gener->SetTrigger(2);
717         gener->SetPtJet(50.);
718         gener->SetRadiation(isw);
719         gener->SetSimpleJets(!isw);
720         gener->SetJetEtaRange(-0.12, 0.12);
721         gener->SetJetPhiRange(220., 320.);
722         gGener=gener;
723       }
724       break;
725
726     case kHijing_gj75:
727       {
728         comment = comment.Append("HIJING Gamma 75 GeV");
729         AliGenHijing *gener = HijingStandard();
730 // impact parameter range
731         gener->SetImpactParameterRange(0., 5.);
732         // trigger
733         gener->SetTrigger(2);
734         gener->SetPtJet(75.);
735         gener->SetRadiation(isw);
736         gener->SetSimpleJets(!isw);
737         gener->SetJetEtaRange(-0.12, 0.12);
738         gener->SetJetPhiRange(220., 320.);
739         gGener=gener;
740       }
741       break;
742
743     case kHijing_gj100:
744       {
745         comment = comment.Append("HIJING Gamma 100 GeV");
746         AliGenHijing *gener = HijingStandard();
747 // impact parameter range
748         gener->SetImpactParameterRange(0., 5.);
749         // trigger
750         gener->SetTrigger(2);
751         gener->SetPtJet(100.);
752         gener->SetRadiation(isw);
753         gener->SetSimpleJets(!isw);
754         gener->SetJetEtaRange(-0.12, 0.12);
755         gener->SetJetPhiRange(220., 320.);
756         gGener=gener;
757       }
758       break;
759
760     case kHijing_gj200:
761       {
762         comment = comment.Append("HIJING Gamma 200 GeV");
763         AliGenHijing *gener = HijingStandard();
764 // impact parameter range
765         gener->SetImpactParameterRange(0., 5.);
766         // trigger
767         gener->SetTrigger(2);
768         gener->SetPtJet(200.);
769         gener->SetRadiation(isw);
770         gener->SetSimpleJets(!isw);
771         gener->SetJetEtaRange(-0.12, 0.12);
772         gener->SetJetPhiRange(220., 320.);
773         gGener=gener;
774       }
775       break;
776     case kHijing_pA:
777       {
778         comment = comment.Append("HIJING pA");
779
780         AliGenCocktail *gener  = new AliGenCocktail();
781 //      gener->SetTrackingFlag(0);
782         gener->SetOrigin(0, 0, 0);    // vertex position
783         gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
784         gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
785         gener->SetVertexSmear(kPerEvent); 
786         AliGenHijing   *hijing = new AliGenHijing(-1);
787 // centre of mass energy 
788         hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
789 // impact parameter range
790         hijing->SetImpactParameterRange(0., 15.);
791 // reference frame
792         hijing->SetReferenceFrame("CMS");
793         hijing->SetBoostLHC(1);
794 // projectile
795         hijing->SetProjectile("P", 1, 1);
796         hijing->SetTarget    ("A", 208, 82);
797 // tell hijing to keep the full parent child chain
798         hijing->KeepFullEvent();
799 // enable jet quenching
800         hijing->SetJetQuenching(0);
801 // enable shadowing
802         hijing->SetShadowing(1);
803 // Don't track spectators
804         hijing->SetSpectators(0);
805 // kinematic selection
806         hijing->SetSelectAll(0);
807 //
808         AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
809         AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
810         gray->SetSlowNucleonModel(model);
811         gray->SetDebug(1);
812         gener->AddGenerator(hijing,"Hijing pPb", 1);
813         gener->AddGenerator(gray,  "Gray Particles",1);
814         gGener=gener;
815       }
816       break;
817     case kPythia6:
818       {
819         comment = comment.Append(":Pythia p-p @ 14 TeV");
820         AliGenPythia *gener = new AliGenPythia(-1); 
821         gener->SetMomentumRange(0,999999);
822         gener->SetPhiRange(-180,180);
823         gener->SetThetaRange(0., 180.);
824         gener->SetYRange(-12,12);
825         gener->SetPtRange(0,1000);
826         gener->SetStrucFunc(kCTEQ4L);   
827         gener->SetProcess(kPyMb);
828         gener->SetEnergyCMS(14000.);
829         gGener=gener;
830       }
831       break;
832     case kPythia6Jets:
833       {
834         comment = comment.Append(":Pythia jets @ 5.5 TeV");
835         AliGenPythia * gener = new AliGenPythia(-1);
836         gener->SetEnergyCMS(5500.);//        Centre of mass energy
837         gener->SetProcess(kPyJets);//        Process type
838         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
839         gener->SetJetPhiRange(0., 360.);
840         gener->SetJetEtRange(10., 1000.);
841         gener->SetGluonRadiation(1,1);
842         //    gener->SetPtKick(0.);
843         //   Structure function
844         gener->SetStrucFunc(kCTEQ4L);
845         gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
846         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
847         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
848         gGener=gener;
849       }
850       break;
851     case kD0PbPb5500:
852       {
853         comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
854         AliGenPythia * gener = new AliGenPythia(10);
855         gener->SetProcess(kPyD0PbPbMNR);
856         gener->SetStrucFunc(kCTEQ4L);
857         gener->SetPtHard(2.1,-1.0);
858         gener->SetEnergyCMS(5500.);
859         gener->SetNuclei(208,208);
860         gGener=gener;
861       }
862       break;
863     }
864     return gGener;
865 }
866
867 AliGenHijing* HijingStandard()
868 {
869     AliGenHijing *gener = new AliGenHijing(-1);
870 // centre of mass energy 
871     gener->SetEnergyCMS(5500.);
872 // reference frame
873     gener->SetReferenceFrame("CMS");
874 // projectile
875      gener->SetProjectile("A", 208, 82);
876      gener->SetTarget    ("A", 208, 82);
877 // tell hijing to keep the full parent child chain
878      gener->KeepFullEvent();
879 // enable jet quenching
880      gener->SetJetQuenching(4);
881 // enable shadowing
882      gener->SetShadowing(1);
883 // neutral pion and heavy particle decays switched off
884      gener->SetDecaysOff(1);
885 // Don't track spectators
886      gener->SetSpectators(0);
887 // kinematic selection
888      gener->SetSelectAll(0);
889      return gener;
890 }
891
892
893