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