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