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