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