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