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