]> git.uio.no Git - u/mrichter/AliRoot.git/blob - EMCAL/Config.C
cuts on Q out, side, long added
[u/mrichter/AliRoot.git] / EMCAL / Config.C
1 enum PprRun_t 
2 {
3     kTest50,
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_jj125, 
8     kHijing_gj25,  kHijing_gj50, kHijing_gj75, kHijing_gj100, kHijing_gj125, 
9     kJetPlusBg,    kGammaPlusBg, 
10     kParam_8000_Ecal, kParam_4000_Ecal, 
11     kJets_100,
12     kGammaGun, kGammaBox
13 };
14
15 enum PprGeo_t 
16 {
17     kHoles, kNoHoles
18 };
19
20 enum PprRad_t
21 {
22     kGluonRadiation, kNoGluonRadiation
23 };
24
25 // This part for configuration    
26 //static PprRun_t run = kParam_8000_Ecal;
27 static PprRun_t run = kJets_100;
28 static PprGeo_t geo = kHoles;
29 static PprRad_t rad = kNoGluonRadiation;
30 // Comment line 
31 static TString  comment;
32
33 void Config()
34 {
35
36     // 7-DEC-2000 09:00
37     // Switch on Transition adiation simulation. 6/12/00 18:00
38     // iZDC=1  7/12/00 09:00
39     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
40     // Theta range given through pseudorapidity limits 22/6/2001
41
42     // Set Random Number seed
43     // gRandom->SetSeed(12345);
44
45     new AliGeant3("C++ Interface to Geant3");
46
47     if (!gSystem->Getenv("CONFIG_FILE"))
48     {
49         TFile  *rootfile = new TFile("galice.root", "recreate");
50
51         rootfile->SetCompressionLevel(2);
52     }
53
54     TGeant3 *geant3 = (TGeant3 *) gMC;
55
56     //
57     // Set External decayer
58     AliDecayer *decayer = new AliDecayerPythia();
59
60     decayer->SetForceDecay(kAll);
61     decayer->Init();
62     gMC->SetExternalDecayer(decayer);
63     //
64     //
65     //=======================================================================
66     // ******* GEANT STEERING parameters FOR ALICE SIMULATION *******
67     geant3->SetTRIG(1);         //Number of events to be processed 
68     geant3->SetSWIT(4, 10);
69     geant3->SetDEBU(0, 0, 1);
70     geant3->SetDCAY(1);
71     geant3->SetPAIR(1);
72     geant3->SetCOMP(1);
73     geant3->SetPHOT(1);
74     geant3->SetPFIS(0);
75     geant3->SetDRAY(0);
76     geant3->SetANNI(1);
77     geant3->SetBREM(1);
78     geant3->SetMUNU(1);
79     geant3->SetCKOV(1);
80     geant3->SetHADR(1);         //Select pure GEANH (HADR 1) or GEANH/NUCRIN (HADR 3)
81     geant3->SetLOSS(2);
82     geant3->SetMULS(1);
83     geant3->SetRAYL(1);
84     geant3->SetAUTO(1);         //Select automatic STMIN etc... calc. (AUTO 1) or manual (AUTO 0)
85     geant3->SetABAN(0);         //Restore 3.16 behaviour for abandoned tracks
86     geant3->SetOPTI(2);         //Select optimisation level for GEANT geometry searches (0,1,2)
87     geant3->SetERAN(5.e-7);
88
89     Float_t cut    = 1.e-3;        // 1MeV cut by default
90     Float_t tofmax = 1.e10;
91
92     //             GAM ELEC NHAD CHAD MUON EBREM MUHAB EDEL MUDEL MUPA TOFMAX
93     geant3->SetCUTS(cut, cut, cut, cut, cut, cut, cut, cut, cut, cut, tofmax);
94     gAlice->TrackingLimits(900, 900);
95     
96     //
97     //=======================================================================
98     // ************* STEERING parameters FOR ALICE SIMULATION **************
99     // --- Specify event type to be tracked through the ALICE setup
100     // --- All positions are in cm, angles in degrees, and P and E in GeV
101
102
103 // Generator Configuration
104     gAlice->SetDebug(1);
105
106     AliGenerator*  gener = GeneratorFactory(run);
107     gener->Init();
108     gener->SetTrackingFlag(1);
109     
110
111
112     if (rad == kGluonRadiation)
113     {
114         //coment= comment.Append(" | Gluon Radiation On");
115         
116     } else {
117         //coment= comment.Append(" | Gluon Radiation Off");
118     }
119
120     if (geo == kHoles)
121     {
122         //coment= comment.Append(" | Holes for PHOS/RICH");
123         
124     } else {
125         //coment= comment.Append(" | No holes for PHOS/RICH");
126     }
127
128     printf("\n \n Comment: %s \n \n", (char*) comment);
129     
130     
131 // Field (L3 0.4 T)
132
133     AliMagFCM* field = new AliMagFCM(
134         "Map2","$(ALICE_ROOT)/data/field01.dat", 2, 1., 10.);
135     field->SetSolenoidField(4.);
136     gAlice->SetField(field);    
137
138     Int_t iMAG    = 0;
139     Int_t iITS    = 0;
140     Int_t iTPC    = 0;
141     Int_t iTOF    = 0;
142     Int_t iRICH   = 0;
143     Int_t iZDC    = 0;
144     Int_t iCASTOR = 0;
145     Int_t iTRD    = 0;
146     Int_t iABSO   = 0;
147     Int_t iDIPO   = 0;
148     Int_t iHALL   = 0;
149     Int_t iFRAME  = 0;
150     Int_t iSHIL   = 0;
151     Int_t iPIPE   = 0;
152     Int_t iFMD    = 0;
153     Int_t iMUON   = 0;
154     Int_t iPHOS   = 0;
155     Int_t iPMD    = 0;
156     Int_t iSTART  = 0;
157     Int_t iVZERO  = 0;
158     Int_t iEMCAL  = 1;    
159 // 
160 //  enable/disable StepManager()
161 //
162
163     Int_t   enableABSO   = 0;
164     Int_t   enableCASTOR = 0;
165     Int_t   enableDIPO   = 0;
166     Int_t   enableFMD    = 0;
167     Int_t   enableFRAME  = 0;
168     Int_t   enableHALL   = 0;
169     Int_t   enableITS    = 0;
170     Int_t   enableMAG    = 0;
171     Int_t   enableMUON   = 0;
172     Int_t   enablePHOS   = 1;
173     Int_t   enablePIPE   = 0;
174     Int_t   enablePMD    = 0;
175     Int_t   enableRICH   = 0;
176     Int_t   enableSHIL   = 0;
177     Int_t   enableSTART  = 0;
178     Int_t   enableTOF    = 0;
179     Int_t   enableTPC    = 0;
180     Int_t   enableTRD    = 0;
181     Int_t   enableZDC    = 0;
182     Int_t   enableEMCAL  = 1;
183
184     //=================== Alice BODY parameters =============================
185     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
186
187
188     if (iMAG)
189     {
190         //=================== MAG parameters ============================
191         // --- Start with Magnet since detector layouts may be depending ---
192         // --- on the selected Magnet dimensions ---
193         AliMAG *MAG = new AliMAG("MAG", "Magnet");
194     }
195
196
197     if (iABSO)
198     {
199         //=================== ABSO parameters ============================
200         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
201     }
202
203     if (iDIPO)
204     {
205         //=================== DIPO parameters ============================
206
207         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
208     }
209
210     if (iHALL)
211     {
212         //=================== HALL parameters ============================
213
214         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
215     }
216
217
218     if (iFRAME)
219     {
220         //=================== FRAME parameters ============================
221
222         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
223         if (geo == kHoles) {
224             FRAME->SetHoles(1);
225         } else {
226             FRAME->SetHoles(0);
227         }
228     }
229
230     if (iSHIL)
231     {
232         //=================== SHIL parameters ============================
233
234         AliSHIL *SHIL = new AliSHILv0("SHIL", "Shielding");
235     }
236
237
238     if (iPIPE)
239     {
240         //=================== PIPE parameters ============================
241
242         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
243     }
244  
245   if(iITS) {
246
247 //=================== ITS parameters ============================
248     //
249     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
250     // almost all other detectors. This involves the fact that the ITS geometry
251     // still has several options to be followed in parallel in order to determine
252     // the best set-up which minimizes the induced background. All the geometries
253     // available to date are described in the following. Read carefully the comments
254     // and use the default version (the only one uncommented) unless you are making
255     // comparisons and you know what you are doing. In this case just uncomment the
256     // ITS geometry you want to use and run Aliroot.
257     //
258     // Detailed geometries:         
259     //
260     //
261     //AliITS *ITS  = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
262     //
263     //AliITS *ITS  = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
264     //
265     AliITSvPPRasymm *ITS  = new AliITSvPPRasymm("ITS","New ITS PPR detailed version with asymmetric services");
266     ITS->SetMinorVersion(2);                                     // don't touch this parameter if you're not an ITS developer
267     ITS->SetReadDet(kFALSE);                                     // don't touch this parameter if you're not an ITS developer
268     //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
269     ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
270     ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
271     ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
272     ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
273     ITS->SetRails(1);        // 1 --> rails in ; 0 --> rails out
274     ITS->SetCoolingFluid(1);   // 1 --> water ; 0 --> freon
275     //
276     //AliITSvPPRsymm *ITS  = new AliITSvPPRsymm("ITS","New ITS PPR detailed version with symmetric services");
277     //ITS->SetMinorVersion(2);                                       // don't touch this parameter if you're not an ITS developer
278     //ITS->SetReadDet(kFALSE);                                       // don't touch this parameter if you're not an ITS developer
279     //ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRsymm2.det"); // don't touch this parameter if you're not an ITS developer
280     //ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
281     //ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
282     //ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
283     //ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
284     //ITS->SetRails(1);              // 1 --> rails in ; 0 --> rails out
285     //ITS->SetCoolingFluid(1);       // 1 --> water ; 0 --> freon
286     //
287     //
288     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
289     // for reconstruction !):
290     //                                                     
291     //
292     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
293     //ITS->SetRails(1);                // 1 --> rails in ; 0 --> rails out
294     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
295     //
296     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
297     //ITS->SetRails(1);                // 1 --> rails in ; 0 --> rails out
298     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
299     //                      
300     //
301     //
302     // Geant3 <-> EUCLID conversion
303     // ============================
304     //
305     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
306     // media to two ASCII files (called by default ITSgeometry.euc and
307     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
308     // The default (=0) means that you dont want to use this facility.
309     //
310     ITS->SetEUCLID(0);  
311     if (!enableITS) ITS->DisableStepManager();
312   }
313   
314
315     if (iTPC)
316     {
317         //============================ TPC parameters ================================
318         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
319         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
320         // --- sectors are specified, any value other than that requires at least one 
321         // --- sector (lower or upper)to be specified!
322         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
323         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
324         // --- SecLows - number of lower sectors specified (up to 6)
325         // --- SecUps - number of upper sectors specified (up to 12)
326         // --- Sens - sensitive strips for the Slow Simulator !!!
327         // --- This does NOT work if all S or L-sectors are specified, i.e.
328         // --- if SecAL or SecAU < 0
329         //
330         //
331         //-----------------------------------------------------------------------------
332
333         //  gROOT->LoadMacro("SetTPCParam.C");
334         //  AliTPCParam *param = SetTPCParam();
335         AliTPC *TPC = new AliTPCv2("TPC", "Default");
336
337         // All sectors included 
338         TPC->SetSecAL(-1);
339         TPC->SetSecAU(-1);
340         if (!enableTPC) TPC->DisableStepManager();
341     }
342
343
344     if (iRICH)
345     {
346         //=================== RICH parameters ===========================
347         AliRICH *RICH = new AliRICHv1("RICH", "normal RICH");
348         if (!enableRICH) RICH->DisableStepManager();
349     }
350
351
352     if (iZDC)
353     {
354         //=================== ZDC parameters ============================
355
356         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
357         if (!enableZDC) ZDC->DisableStepManager();
358     }
359
360     if (iCASTOR)
361     {
362         //=================== CASTOR parameters ============================
363
364         AliCASTOR *CASTOR = new AliCASTORv1("CASTOR", "normal CASTOR");
365         if (!enableCASTOR) CASTOR->DisableStepManager();
366     }
367
368     if (iTRD)
369     {
370         //=================== TRD parameters ============================
371
372         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
373
374         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
375         TRD->SetGasMix(1);
376         if (geo == kHoles) {
377             // With hole in front of PHOS
378             TRD->SetPHOShole();
379             // With hole in front of RICH
380             TRD->SetRICHhole();
381         }
382             // Switch on TR
383             AliTRDsim *TRDsim = TRD->CreateTR();
384             if (!enableTRD) TRD->DisableStepManager();
385     }
386
387     if (iFMD)
388     {
389         //=================== FMD parameters ============================
390
391         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
392         FMD->SetRingsSi1(256);
393         FMD->SetRingsSi2(64);
394         FMD->SetSectorsSi1(20);
395         FMD->SetSectorsSi2(24);
396         if (!enableFMD) FMD->DisableStepManager();
397    }
398
399     if (iMUON)
400     {
401         //=================== MUON parameters ===========================
402
403         AliMUON *MUON = new AliMUONv1("MUON", "default");
404         if (!enableMUON) MUON->DisableStepManager();
405     }
406     //=================== PHOS parameters ===========================
407
408     if (iPHOS)
409     {
410         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
411         if (!enablePHOS) PHOS->DisableStepManager();
412     }
413
414
415     if (iPMD)
416     {
417         //=================== PMD parameters ============================
418         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
419         PMD->SetPAR(1., 1., 0.8, 0.02);
420         PMD->SetIN(6., 18., -580., 27., 27.);
421         PMD->SetGEO(0.0, 0.2, 4.);
422         PMD->SetPadSize(0.8, 1.0, 1.0, 1.5);
423         if (!enablePMD) PMD->DisableStepManager();
424     }
425
426     if (iTOF) {
427         if (geo == kHoles) {
428         //=================== TOF parameters ============================
429             AliTOF *TOF = new AliTOFv2FHoles("TOF", "TOF with Holes");
430         } else {
431             AliTOF *TOF = new AliTOFv4("TOF", "normal TOF");
432         }
433         if (!enableTOF) TOF->DisableStepManager();
434     }
435
436     if (iSTART)
437     {
438         //=================== START parameters ============================
439         AliSTART *START = new AliSTARTv1("START", "START Detector");
440         if (!enableSTART) START->DisableStepManager();
441     }
442
443     if (iEMCAL && !iRICH)
444     {
445         //=================== EMCAL parameters ============================
446         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCALArch1b");
447         if (!enableEMCAL) EMCAL->DisableStepManager();
448     }
449 }
450
451 Float_t EtaToTheta(Float_t arg){
452   return (180./TMath::Pi())*2.*atan(exp(-arg));
453 }
454
455
456
457 AliGenerator* GeneratorFactory(PprRun_t run) {
458     Int_t isw = 3;
459     if (rad == kNoGluonRadiation) isw = 0;
460     
461
462     switch (run) {
463     case kGammaGun:
464         gener = new AliGenFixed(1);
465         gener->SetMomentum(100.);
466         gener->SetPhi(240.);
467         gener->SetTheta(91.);
468         gener->SetPart(kElectron);
469      break;
470
471     case kGammaBox:
472         gener = new AliGenBox(100);
473         gener->SetMomentumRange(100., 100.1);
474         gener->SetPhiRange(0,360);
475         gener->SetThetaRange(45., 135.);
476         gener->SetPart(kElectron);
477         break;
478         
479     case kTest50:
480         //coment= comment.Append(":HIJINGparam test 50 particles");
481         gener = new AliGenHIJINGparaBa(50);
482         gener->SetMomentumRange(0, 999999.);
483         gener->SetPhiRange(-180., 180.);
484         // Set pseudorapidity range from -8 to 8.
485         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
486         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
487         gener->SetThetaRange(thmin,thmax);
488         break;
489
490     case kParam_8000:
491         //coment= comment.Append(":HIJINGparam N=8000");
492         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
493         gener->SetMomentumRange(0, 999999.);
494         gener->SetPhiRange(-180., 180.);
495         // Set pseudorapidity range from -8 to 8.
496         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
497         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
498         gener->SetThetaRange(thmin,thmax);
499         break;
500     case kParam_4000:
501         //coment= comment.Append("HIJINGparam N=4000");
502         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
503         gener->SetMomentumRange(0, 999999.);
504         gener->SetPhiRange(-180., 180.);
505         // Set pseudorapidity range from -8 to 8.
506         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
507         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
508         gener->SetThetaRange(thmin,thmax);
509         break;
510     case kParam_2000:
511         //coment= comment.Append("HIJINGparam N=2000");
512         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
513         gener->SetMomentumRange(0, 999999.);
514         gener->SetPhiRange(-180., 180.);
515         // Set pseudorapidity range from -8 to 8.
516         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
517         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
518         gener->SetThetaRange(thmin,thmax);
519         break;
520
521     case kParam_8000_Ecal:
522         //coment= comment.Append(":HIJINGparam N=8000 baryons");
523         AliGenHIJINGparaBa *gener = new AliGenHIJINGparaBa(82534);
524         gener->SetMomentumRange(0, 999999.);
525         gener->SetPhiRange(-180., 180.);
526         // Set pseudorapidity range from -8 to 8.
527         Float_t thmin = EtaToTheta( 5);   // theta min. <---> eta max
528         Float_t thmax = EtaToTheta(-5);  // theta max. <---> eta min 
529         gener->SetThetaRange(thmin,thmax);
530         break;
531
532     case kParam_4000_Ecal:
533         //coment= comment.Append(":HIJINGparam N=8000 baryons");
534         AliGenHIJINGparaBa *gener = new AliGenHIJINGparaBa(82534/2.);
535         gener->SetMomentumRange(0, 999999.);
536         gener->SetPhiRange(-180., 180.);
537         // Set pseudorapidity range from -8 to 8.
538         Float_t thmin = EtaToTheta( 5);   // theta min. <---> eta max
539         Float_t thmax = EtaToTheta(-5);  // theta max. <---> eta min 
540         gener->SetThetaRange(thmin,thmax);
541         break;
542 //
543 //  Hijing Central
544 //
545     case kHijing_cent1:
546         //coment= comment.Append("HIJING cent1");
547         AliGenHijing *gener = HijingStandard();
548 // impact parameter range
549         gener->SetImpactParameterRange(0., 5.);
550         break;
551     case kHijing_cent2:
552         //coment= comment.Append("HIJING cent2");
553         AliGenHijing *gener = HijingStandard();
554 // impact parameter range
555         gener->SetImpactParameterRange(0., 2.);
556         break;
557 //
558 // Hijing Peripheral 
559 //
560     case kHijing_per1:
561         //coment= comment.Append("HIJING per1");
562         AliGenHijing *gener = HijingStandard();
563 // impact parameter range
564         gener->SetImpactParameterRange(5., 8.6);
565         break;
566     case kHijing_per2:
567         //coment= comment.Append("HIJING per2");
568         AliGenHijing *gener = HijingStandard();
569 // impact parameter range
570         gener->SetImpactParameterRange(8.6, 11.2);
571         break;
572     case kHijing_per3:
573         //coment= comment.Append("HIJING per3");
574         AliGenHijing *gener = HijingStandard();
575 // impact parameter range
576         gener->SetImpactParameterRange(11.2, 13.2);
577         break;
578     case kHijing_per4:
579         //coment= comment.Append("HIJING per4");
580         AliGenHijing *gener = HijingStandard();
581 // impact parameter range
582         gener->SetImpactParameterRange(13.2, 15.);
583         break;
584     case kHijing_per5:
585         //coment= comment.Append("HIJING per5");
586         AliGenHijing *gener = HijingStandard();
587 // impact parameter range
588         gener->SetImpactParameterRange(15., 100.);
589         break;
590 //
591 //  Jet-Jet
592 //
593     case kHijing_jj25:
594         //coment= comment.Append("HIJING Jet 25 GeV");
595         AliGenHijing *gener = HijingStandard();
596 // impact parameter range
597         gener->SetImpactParameterRange(0., 5.);
598         // trigger
599         gener->SetTrigger(1);
600         gener->SetPtMinJet(25.);
601         gener->SetRadiation(isw);
602         gener->SetJetEtaRange(-0.3,0.3);
603         gener->SetJetPhiRange(15.,105.);   
604         break;
605
606     case kHijing_jj50:
607         //coment= comment.Append("HIJING Jet 50 GeV");
608         AliGenHijing *gener = HijingStandard();
609 // impact parameter range
610         gener->SetImpactParameterRange(0., 5.);
611         // trigger
612         gener->SetTrigger(1);
613         gener->SetPtMinJet(50.);
614         gener->SetRadiation(isw);
615         gener->SetJetEtaRange(-0.3,0.3);
616         gener->SetJetPhiRange(15.,105.);   
617         break;
618
619     case kHijing_jj75:
620         //coment= comment.Append("HIJING Jet 75 GeV");
621         AliGenHijing *gener = HijingStandard();
622 // impact parameter range
623         gener->SetImpactParameterRange(0., 5.);
624         // trigger
625         gener->SetTrigger(1);
626         gener->SetPtMinJet(75.);
627         gener->SetRadiation(isw);
628         gener->SetJetEtaRange(-0.3,0.3);
629         gener->SetJetPhiRange(15.,105.);   
630         break;
631
632     case kHijing_jj100:
633         //coment= comment.Append("HIJING Jet 100 GeV");
634         AliGenHijing *gener = HijingStandard();
635 // impact parameter range
636         gener->SetImpactParameterRange(10., 10.1);
637         // trigger
638         gener->SetTrigger(1);
639         gener->SetPtMinJet(100.);
640         gener->SetRadiation(isw);
641         gener->SetJetEtaRange(-0.3,0.3);
642         gener->SetJetPhiRange(15.,105.);   
643         break;
644
645     case kHijing_jj125:
646         //coment= comment.Append("HIJING Jet 125 GeV");
647         AliGenHijing *gener = HijingStandard();
648 // impact parameter range
649         gener->SetImpactParameterRange(0., 5.);
650         // trigger
651         gener->SetTrigger(1);
652         gener->SetPtMinJet(125.);
653         gener->SetRadiation(isw);
654         gener->SetJetEtaRange(-0.3,0.3);
655         gener->SetJetPhiRange(15.,105.);   
656         break;
657 //
658 // Gamma-Jet
659 //
660     case kHijing_gj25:
661         //coment= comment.Append("HIJING Gamma 25 GeV");
662         AliGenHijing *gener = HijingStandard();
663 // impact parameter range
664         gener->SetImpactParameterRange(0., 5.);
665         // trigger
666         gener->SetTrigger(2);
667         gener->SetPtMinJet(25.);
668         gener->SetRadiation(isw);
669         gener->SetJetEtaRange(-0.3,0.3);
670         gener->SetJetPhiRange(15.,105.);   
671         break;
672
673     case kHijing_gj50:
674         //coment= comment.Append("HIJING Gamma 50 GeV");
675         AliGenHijing *gener = HijingStandard();
676 // impact parameter range
677         gener->SetImpactParameterRange(0., 5.);
678         // trigger
679         gener->SetTrigger(2);
680         gener->SetPtMinJet(50.);
681         gener->SetRadiation(isw);
682         gener->SetJetEtaRange(-0.3,0.3);
683         gener->SetJetPhiRange(15.,105.);   
684         break;
685
686     case kHijing_gj75:
687         //coment= comment.Append("HIJING Gamma 75 GeV");
688         AliGenHijing *gener = HijingStandard();
689 // impact parameter range
690         gener->SetImpactParameterRange(0., 5.);
691         // trigger
692         gener->SetTrigger(2);
693         gener->SetPtMinJet(75.);
694         gener->SetRadiation(isw);
695         gener->SetJetEtaRange(-0.3,0.3);
696         gener->SetJetPhiRange(15.,105.);   
697         break;
698
699     case kHijing_gj100:
700         //coment= comment.Append("HIJING Gamma 100 GeV");
701         AliGenHijing *gener = HijingStandard();
702 // impact parameter range
703         gener->SetImpactParameterRange(0., 5.);
704         // trigger
705         gener->SetTrigger(2);
706         gener->SetPtMinJet(100.);
707         gener->SetRadiation(isw);
708         gener->SetJetEtaRange(-0.3,0.3);
709         gener->SetJetPhiRange(15.,105.);   
710         break;
711
712     case kHijing_gj125:
713         //coment= comment.Append("HIJING Gamma 125 GeV");
714         AliGenHijing *gener = HijingStandard();
715 // impact parameter range
716         gener->SetImpactParameterRange(0., 5.);
717         // trigger
718         gener->SetTrigger(2);
719         gener->SetPtMinJet(125.);
720         gener->SetRadiation(isw);
721         gener->SetJetEtaRange(-0.3,0.3);
722         gener->SetJetPhiRange(15.,105.);   
723         break;
724     case kJetPlusBg:
725         AliGenCocktail *gener = new AliGenCocktail();
726         gener->SetMomentumRange(0, 999999.);
727         gener->SetPhiRange(-180., 180.);
728         // Set pseudorapidity range from -8 to 8.
729         Float_t thmin = EtaToTheta( 5.);   // theta min. <---> eta max
730         Float_t thmax = EtaToTheta(-5.);  // theta max. <---> eta min 
731         gener->SetThetaRange(thmin,thmax);
732
733 //
734 //      Underlying Event
735 //
736 //      AliGenHIJINGparaBa *bg = new AliGenHIJINGparaBa(82534);
737         AliGenHIJINGparaBa *bg = new AliGenHIJINGparaBa(10);
738 //
739 //      Jets from Pythia
740 //
741         AliGenPythia *jets = new AliGenPythia(-1);
742 //   Centre of mass energy 
743         jets->SetEnergyCMS(5500.);
744 //   Process type
745         jets->SetProcess(kPyJets);
746 //   final state kinematic cuts
747         jets->SetJetEtaRange(-0.3, 0.3);
748         jets->SetJetPhiRange(15., 105.);
749 //   Structure function
750         jets->SetStrucFunc(kGRV_LO_98);
751 //   
752 //   Pt transfer of the hard scattering
753         jets->SetPtHard(100.,100.1);
754 //   Decay type (semielectronic, semimuonic, nodecay)
755         jets->SetForceDecay(kAll);
756 //
757 //      Add all to cockail ...    
758 //
759         gener->AddGenerator(bg,"Underlying Event", 1);
760 //      gener->AddGenerator(jets,"Jets",1);
761      break;     
762     case kGammaPlusBg:
763         AliGenCocktail *gener = new AliGenCocktail();
764         gener->SetMomentumRange(0, 999999.);
765         gener->SetPhiRange(-180., 180.);
766         // Set pseudorapidity range from -8 to 8.
767         Float_t thmin = EtaToTheta( 5.);   // theta min. <---> eta max
768         Float_t thmax = EtaToTheta(-5.);  // theta max. <---> eta min 
769         gener->SetThetaRange(thmin,thmax);
770
771 //
772 //      Underlying Event
773 //
774         AliGenHIJINGparaBa *bg = new AliGenHIJINGparaBa(82534);
775 //
776 //      Jets from Pythia
777 //
778         AliGenPythia *jets = new AliGenPythia(-1);
779 //   Centre of mass energy 
780         jets->SetEnergyCMS(5500.);
781 //   Process type
782         jets->SetProcess(kPyDirectGamma);
783 //   final state kinematic cuts
784         jets->SetJetEtaRange(-0.3, 0.3);
785         jets->SetJetPhiRange(15., 105.);
786         jets->SetGammaEtaRange(-0.12, 0.12);
787         jets->SetGammaPhiRange(220., 320.);
788 //   Structure function
789         jets->SetStrucFunc(kGRV_LO_98);
790 //   
791 //   Pt transfer of the hard scattering
792         jets->SetPtHard(100.,100.1);
793 //   Decay type (semielectronic, semimuonic, nodecay)
794         jets->SetForceDecay(kAll);
795 //
796 //      Add all to cockail ...    
797 //
798         gener->AddGenerator(bg,"Underlying Event", 1);
799         gener->AddGenerator(jets,"Jets",1);
800         break;
801     case kJets_100:
802 //  100 GeV Jets  
803         AliGenPythia *gener = PythiaJets(100.);
804         break;
805     }
806
807     return gener;
808 }
809
810 AliGenHijing* HijingStandard()
811 {
812     AliGenHijing *gener = new AliGenHijing(-1);
813 // centre of mass energy 
814     gener->SetEnergyCMS(5500.);
815 // reference frame
816     gener->SetReferenceFrame("CMS");
817 // projectile
818     gener->SetProjectile("A", 208, 82);
819     gener->SetTarget    ("A", 208, 82);
820 // tell hijing to keep the full parent child chain
821     gener->KeepFullEvent();
822 // enable jet quenching
823     gener->SetJetQuenching(1);
824 // enable shadowing
825     gener->SetShadowing(1);
826 // neutral pion and heavy particle decays switched off
827     gener->SetDecaysOff(1);
828 // Don't track spectators
829     gener->SetSpectators(0);
830 // kinematic selection
831     gener->SetSelectAll(0);
832     return gener;
833 }
834
835 AliGenPythia* PythiaJets(Float_t energy)
836 {
837     AliGenPythia *gener = new AliGenPythia(-1);
838 //   Centre of mass energy 
839     gener->SetEnergyCMS(5500.);
840 //   Process type
841     gener->SetProcess(kPyJets);
842 //   final state kinematic cuts
843     gener->SetJetEtaRange(-0.3, 0.3);
844     gener->SetJetPhiRange(15., 105.);
845 //   Structure function
846     gener->SetStrucFunc(kGRV_LO_98);
847 //   
848 //   Pt transfer of the hard scattering
849     gener->SetPtHard(energy, energy+0.1);
850 //   Decay type (semielectronic, semimuonic, nodecay)
851     gener->SetForceDecay(kAll);
852 //
853     return gener;
854 }
855
856
857 AliGenPythia* PythiaGamma(Float_t energy)
858 {
859     AliGenPythia *gener = new AliGenPythia(-1);
860 //   Centre of mass energy 
861     gener->SetEnergyCMS(5500.);
862 //   Process type
863     gener->SetProcess(kPyDirectGamma);
864 //   final state kinematic cuts
865     gener->SetJetEtaRange(-0.3, 0.3);
866     gener->SetJetPhiRange(15., 105.);
867     gener->SetGammaEtaRange(-0.12, 0.12);
868     gener->SetGammaPhiRange(220., 320.);
869 //   Structure function
870     gener->SetStrucFunc(kGRV_LO_98);
871 //   
872 //   Pt transfer of the hard scattering
873     gener->SetPtHard(energy, energy+0.1);
874 //   Decay type (semielectronic, semimuonic, nodecay)
875     gener->SetForceDecay(kAll);
876 //
877     return gener;
878 }
879
880
881
882
883
884
885
886
887
888
889