]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - test/PbPbbench/Config.C
specified projecile and target to the generator
[u/mrichter/AliRoot.git] / test / PbPbbench / Config.C
... / ...
CommitLineData
1// One can use the configuration macro in compiled mode by
2// root [0] gSystem->Load("libgeant321");
3// root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
4// -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
5// root [0] .x grun.C(1,"ConfigPPR.C++")
6
7#if !defined(__CINT__) || defined(__MAKECINT__)
8#include <Riostream.h>
9#include <TRandom.h>
10#include <TSystem.h>
11#include <TVirtualMC.h>
12#include <TGeant3TGeo.h>
13#include <TPDGCode.h>
14#include <TF1.h>
15#include "STEER/AliRunLoader.h"
16#include "STEER/AliRun.h"
17#include "STEER/AliConfig.h"
18#include "STEER/AliGenerator.h"
19#include "STEER/AliLog.h"
20#include "PYTHIA6/AliDecayerPythia.h"
21#include "EVGEN/AliGenHIJINGpara.h"
22#include "THijing/AliGenHijing.h"
23#include "EVGEN/AliGenCocktail.h"
24#include "EVGEN/AliGenSlowNucleons.h"
25#include "EVGEN/AliSlowNucleonModelExp.h"
26#include "EVGEN/AliGenParam.h"
27#include "EVGEN/AliGenMUONlib.h"
28#include "EVGEN/AliGenSTRANGElib.h"
29#include "EVGEN/AliGenMUONCocktail.h"
30#include "EVGEN/AliGenCocktail.h"
31#include "EVGEN/AliGenGeVSim.h"
32#include "EVGEN/AliGeVSimParticle.h"
33#include "PYTHIA6/AliGenPythia.h"
34#include "STEER/AliMagF.h"
35#include "STRUCT/AliBODY.h"
36#include "STRUCT/AliMAG.h"
37#include "STRUCT/AliABSOv3.h"
38#include "STRUCT/AliDIPOv3.h"
39#include "STRUCT/AliHALLv3.h"
40#include "STRUCT/AliFRAMEv2.h"
41#include "STRUCT/AliSHILv3.h"
42#include "STRUCT/AliPIPEv3.h"
43#include "ITS/AliITSv11Hybrid.h"
44#include "TPC/AliTPCv2.h"
45#include "TOF/AliTOFv6T0.h"
46#include "HMPID/AliHMPIDv3.h"
47#include "ZDC/AliZDCv3.h"
48#include "TRD/AliTRDv1.h"
49#include "FMD/AliFMDv1.h"
50#include "MUON/AliMUONv1.h"
51#include "PHOS/AliPHOSv1.h"
52#include "PMD/AliPMDv1.h"
53#include "T0/AliT0v1.h"
54#include "EMCAL/AliEMCALv2.h"
55#include "ACORDE/AliACORDEv1.h"
56#include "VZERO/AliVZEROv7.h"
57#endif
58
59enum PprRun_t
60{
61 test50,
62 kParam_8000, kParam_4000, kParam_2000,
63 kHijing_cent1, kHijing_cent2,
64 kHijing_per1, kHijing_per2, kHijing_per3, kHijing_per4, kHijing_per5,
65 kHijing_jj25, kHijing_jj50, kHijing_jj75, kHijing_jj100, kHijing_jj200,
66 kHijing_gj25, kHijing_gj50, kHijing_gj75, kHijing_gj100, kHijing_gj200,
67 kHijing_pA, kPythia6,
68 kPythia6Jets20_24, kPythia6Jets24_29, kPythia6Jets29_35,
69 kPythia6Jets35_42, kPythia6Jets42_50, kPythia6Jets50_60,
70 kPythia6Jets60_72, kPythia6Jets72_86, kPythia6Jets86_104,
71 kPythia6Jets104_125, kPythia6Jets125_150, kPythia6Jets150_180,
72 kD0PbPb5500, kCharmSemiElPbPb5500, kBeautySemiElPbPb5500,
73 kCocktailTRD, kPyJJ, kPyGJ,
74 kMuonCocktailCent1, kMuonCocktailPer1, kMuonCocktailPer4,
75 kMuonCocktailCent1HighPt, kMuonCocktailPer1HighPt, kMuonCocktailPer4HighPt,
76 kMuonCocktailCent1Single, kMuonCocktailPer1Single, kMuonCocktailPer4Single,
77 kFlow_2_2000, kFlow_10_2000, kFlow_6_2000, kFlow_6_5000,
78 kHIJINGplus, kRunMax
79};
80
81const char* pprRunName[] = {
82 "test50",
83 "kParam_8000", "kParam_4000", "kParam_2000",
84 "kHijing_cent1", "kHijing_cent2",
85 "kHijing_per1", "kHijing_per2", "kHijing_per3", "kHijing_per4",
86 "kHijing_per5",
87 "kHijing_jj25", "kHijing_jj50", "kHijing_jj75", "kHijing_jj100",
88 "kHijing_jj200",
89 "kHijing_gj25", "kHijing_gj50", "kHijing_gj75", "kHijing_gj100",
90 "kHijing_gj200", "kHijing_pA", "kPythia6",
91 "kPythia6Jets20_24", "kPythia6Jets24_29", "kPythia6Jets29_35",
92 "kPythia6Jets35_42", "kPythia6Jets42_50", "kPythia6Jets50_60",
93 "kPythia6Jets60_72", "kPythia6Jets72_86", "kPythia6Jets86_104",
94 "kPythia6Jets104_125", "kPythia6Jets125_150", "kPythia6Jets150_180",
95 "kD0PbPb5500", "kCharmSemiElPbPb5500", "kBeautySemiElPbPb5500",
96 "kCocktailTRD", "kPyJJ", "kPyGJ",
97 "kMuonCocktailCent1", "kMuonCocktailPer1", "kMuonCocktailPer4",
98 "kMuonCocktailCent1HighPt", "kMuonCocktailPer1HighPt", "kMuonCocktailPer4HighPt",
99 "kMuonCocktailCent1Single", "kMuonCocktailPer1Single", "kMuonCocktailPer4Single",
100 "kFlow_2_2000", "kFlow_10_2000", "kFlow_6_2000", "kFlow_6_5000", "kHIJINGplus"
101};
102
103enum PprRad_t
104{
105 kGluonRadiation, kNoGluonRadiation
106};
107
108enum PprTrigConf_t
109{
110 kDefaultPPTrig, kDefaultPbPbTrig
111};
112
113const char * pprTrigConfName[] = {
114 "p-p","Pb-Pb"
115};
116
117// This part for configuration
118
119static PprRun_t srun = kHijing_per2;
120static PprRad_t srad = kGluonRadiation;
121static AliMagF::BMap_t smag = AliMagF::k5kG;
122static Int_t sseed = 12345; //Set 0 to use the current time
123static PprTrigConf_t strig = kDefaultPbPbTrig; // default pp trigger configuration
124
125// Comment line
126static TString comment;
127
128// Functions
129Float_t EtaToTheta(Float_t arg);
130AliGenerator* GeneratorFactory(PprRun_t srun);
131AliGenHijing* HijingStandard();
132AliGenGeVSim* GeVSimStandard(Float_t, Float_t);
133void ProcessEnvironmentVars();
134
135void Config()
136{
137 // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
138 // Theta range given through pseudorapidity limits 22/6/2001
139
140 // Get settings from environment variables
141 ProcessEnvironmentVars();
142
143 // Set Random Number seed
144 gRandom->SetSeed(sseed);
145 cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl;
146
147
148 // libraries required by geant321
149#if defined(__CINT__)
150 gSystem->Load("liblhapdf");
151 gSystem->Load("libEGPythia6");
152 gSystem->Load("libpythia6");
153 gSystem->Load("libAliPythia6");
154 gSystem->Load("libgeant321");
155#endif
156
157 new TGeant3TGeo("C++ Interface to Geant3");
158
159 // Output every 100 tracks
160 ((TGeant3*)gMC)->SetSWIT(4,100);
161
162 AliRunLoader* rl=0x0;
163
164 AliLog::Message(AliLog::kInfo, "Creating Run Loader", "", "", "Config()"," ConfigPPR.C", __LINE__);
165
166 rl = AliRunLoader::Open("galice.root",
167 AliConfig::GetDefaultEventFolderName(),
168 "recreate");
169 if (rl == 0x0)
170 {
171 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
172 return;
173 }
174 rl->SetCompressionLevel(2);
175 rl->SetNumberOfEventsPerFile(100);
176 gAlice->SetRunLoader(rl);
177
178 // Set the trigger configuration
179 gAlice->SetTriggerDescriptor(pprTrigConfName[strig]);
180 cout<<"Trigger configuration is set to "<<pprTrigConfName[strig]<<endl;
181
182 //
183 // Set External decayer
184 AliDecayer *decayer = new AliDecayerPythia();
185
186
187 switch (srun) {
188 case kD0PbPb5500:
189 decayer->SetForceDecay(kHadronicD);
190 break;
191 case kCharmSemiElPbPb5500:
192 decayer->SetForceDecay(kSemiElectronic);
193 break;
194 case kBeautySemiElPbPb5500:
195 decayer->SetForceDecay(kSemiElectronic);
196 break;
197 default:
198 decayer->SetForceDecay(kAll);
199 break;
200 }
201 decayer->Init();
202 gMC->SetExternalDecayer(decayer);
203 //
204 //
205 //=======================================================================
206 //
207 //=======================================================================
208 // ************* STEERING parameters FOR ALICE SIMULATION **************
209 // --- Specify event type to be tracked through the ALICE setup
210 // --- All positions are in cm, angles in degrees, and P and E in GeV
211
212 gMC->SetProcess("DCAY",1);
213 gMC->SetProcess("PAIR",1);
214 gMC->SetProcess("COMP",1);
215 gMC->SetProcess("PHOT",1);
216 gMC->SetProcess("PFIS",0);
217 gMC->SetProcess("DRAY",0);
218 gMC->SetProcess("ANNI",1);
219 gMC->SetProcess("BREM",1);
220 gMC->SetProcess("MUNU",1);
221 gMC->SetProcess("CKOV",1);
222 gMC->SetProcess("HADR",1);
223 gMC->SetProcess("LOSS",2);
224 gMC->SetProcess("MULS",1);
225 gMC->SetProcess("RAYL",1);
226
227 Float_t cut = 1.e-3; // 1MeV cut by default
228 Float_t tofmax = 1.e10;
229
230 gMC->SetCut("CUTGAM", cut);
231 gMC->SetCut("CUTELE", cut);
232 gMC->SetCut("CUTNEU", cut);
233 gMC->SetCut("CUTHAD", cut);
234 gMC->SetCut("CUTMUO", cut);
235 gMC->SetCut("BCUTE", cut);
236 gMC->SetCut("BCUTM", cut);
237 gMC->SetCut("DCUTE", cut);
238 gMC->SetCut("DCUTM", cut);
239 gMC->SetCut("PPCUTM", cut);
240 gMC->SetCut("TOFMAX", tofmax);
241
242 // Generator Configuration
243 AliGenerator* gener = GeneratorFactory(srun);
244 gener->SetOrigin(0, 0, 0); // vertex position
245 gener->SetSigma(0, 0, 5.3); // Sigma in (X,Y,Z) (cm) on IP position
246 gener->SetCutVertexZ(1.); // Truncate at 1 sigma
247 gener->SetVertexSmear(kPerEvent);
248 gener->SetTrackingFlag(1);
249 gener->Init();
250
251 if (smag == AliMagF::k2kG) {
252 comment = comment.Append(" | L3 field 0.2 T");
253 } else if (smag == AliMagF::k5kG) {
254 comment = comment.Append(" | L3 field 0.5 T");
255 }
256
257
258 if (srad == kGluonRadiation)
259 {
260 comment = comment.Append(" | Gluon Radiation On");
261
262 } else {
263 comment = comment.Append(" | Gluon Radiation Off");
264 }
265
266 printf("\n \n Comment: %s \n \n", comment.Data());
267
268
269// Field
270 TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 2, 1., 1., 10., smag));
271 rl->CdGAFile();
272//
273 Int_t iABSO = 1;
274 Int_t iDIPO = 1;
275 Int_t iFMD = 1;
276 Int_t iFRAME = 1;
277 Int_t iHALL = 1;
278 Int_t iITS = 1;
279 Int_t iMAG = 1;
280 Int_t iMUON = 1;
281 Int_t iPHOS = 1;
282 Int_t iPIPE = 1;
283 Int_t iPMD = 1;
284 Int_t iHMPID = 1;
285 Int_t iSHIL = 1;
286 Int_t iT0 = 1;
287 Int_t iTOF = 1;
288 Int_t iTPC = 1;
289 Int_t iTRD = 1;
290 Int_t iZDC = 1;
291 Int_t iEMCAL = 1;
292 Int_t iVZERO = 1;
293 Int_t iACORDE = 0;
294
295 //=================== Alice BODY parameters =============================
296 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
297
298
299 if (iMAG)
300 {
301 //=================== MAG parameters ============================
302 // --- Start with Magnet since detector layouts may be depending ---
303 // --- on the selected Magnet dimensions ---
304 AliMAG *MAG = new AliMAG("MAG", "Magnet");
305 }
306
307
308 if (iABSO)
309 {
310 //=================== ABSO parameters ============================
311 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
312 }
313
314 if (iDIPO)
315 {
316 //=================== DIPO parameters ============================
317
318 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
319 }
320
321 if (iHALL)
322 {
323 //=================== HALL parameters ============================
324
325 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
326 }
327
328
329 if (iFRAME)
330 {
331 //=================== FRAME parameters ============================
332
333 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
334 FRAME->SetHoles(1);
335 }
336
337 if (iSHIL)
338 {
339 //=================== SHIL parameters ============================
340
341 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
342 }
343
344
345 if (iPIPE)
346 {
347 //=================== PIPE parameters ============================
348
349 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
350 }
351
352 if (iITS)
353 {
354 //=================== ITS parameters ============================
355
356 AliITS *ITS = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
357 }
358
359 if (iTPC)
360 {
361 //============================ TPC parameters =====================
362 AliTPC *TPC = new AliTPCv2("TPC", "Default");
363 }
364
365
366 if (iTOF) {
367 //=================== TOF parameters ============================
368 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
369 }
370
371
372 if (iHMPID)
373 {
374 //=================== HMPID parameters ===========================
375 AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
376
377 }
378
379
380 if (iZDC)
381 {
382 //=================== ZDC parameters ============================
383
384 AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
385 }
386
387 if (iTRD)
388 {
389 //=================== TRD parameters ============================
390
391 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
392 }
393
394 if (iFMD)
395 {
396 //=================== FMD parameters ============================
397 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
398 }
399
400 if (iMUON)
401 {
402 //=================== MUON parameters ===========================
403 // New MUONv1 version (geometry defined via builders)
404 AliMUON *MUON = new AliMUONv1("MUON", "default");
405 }
406 //=================== PHOS parameters ===========================
407
408 if (iPHOS)
409 {
410 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
411 }
412
413
414 if (iPMD)
415 {
416 //=================== PMD parameters ============================
417 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
418 }
419
420 if (iT0)
421 {
422 //=================== T0 parameters ============================
423 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
424 }
425
426 if (iEMCAL)
427 {
428 //=================== EMCAL parameters ============================
429 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
430 }
431
432 if (iACORDE)
433 {
434 //=================== ACORDE parameters ============================
435 AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
436 }
437
438 if (iVZERO)
439 {
440 //=================== VZERO parameters ============================
441 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
442 }
443
444
445}
446
447Float_t EtaToTheta(Float_t arg){
448 return (180./TMath::Pi())*2.*atan(exp(-arg));
449}
450
451
452
453AliGenerator* GeneratorFactory(PprRun_t srun) {
454 Int_t isw = 3;
455 if (srad == kNoGluonRadiation) isw = 0;
456
457
458 AliGenerator * gGener = 0x0;
459 switch (srun) {
460 case test50:
461 {
462 comment = comment.Append(":HIJINGparam test 50 particles");
463 AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
464 gener->SetMomentumRange(0, 999999.);
465 gener->SetPhiRange(0., 360.);
466 // Set pseudorapidity range from -8 to 8.
467 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
468 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
469 gener->SetThetaRange(thmin,thmax);
470 gGener=gener;
471 }
472 break;
473 case kParam_8000:
474 {
475 comment = comment.Append(":HIJINGparam N=8000");
476 AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
477 gener->SetMomentumRange(0, 999999.);
478 gener->SetPhiRange(0., 360.);
479 // Set pseudorapidity range from -8 to 8.
480 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
481 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
482 gener->SetThetaRange(thmin,thmax);
483 gGener=gener;
484 }
485 break;
486 case kParam_4000:
487 {
488 comment = comment.Append("HIJINGparam N=4000");
489 AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
490 gener->SetMomentumRange(0, 999999.);
491 gener->SetPhiRange(0., 360.);
492 // Set pseudorapidity range from -8 to 8.
493 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
494 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
495 gener->SetThetaRange(thmin,thmax);
496 gGener=gener;
497 }
498 break;
499 case kParam_2000:
500 {
501 comment = comment.Append("HIJINGparam N=2000");
502 AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
503 gener->SetMomentumRange(0, 999999.);
504 gener->SetPhiRange(0., 360.);
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 gGener=gener;
510 }
511 break;
512 //
513 // Hijing Central
514 //
515 case kHijing_cent1:
516 {
517 comment = comment.Append("HIJING cent1");
518 AliGenHijing *gener = HijingStandard();
519 // impact parameter range
520 gener->SetImpactParameterRange(0., 5.);
521 gGener=gener;
522 }
523 break;
524 case kHijing_cent2:
525 {
526 comment = comment.Append("HIJING cent2");
527 AliGenHijing *gener = HijingStandard();
528 // impact parameter range
529 gener->SetImpactParameterRange(0., 2.);
530 gGener=gener;
531 }
532 break;
533 //
534 // Hijing Peripheral
535 //
536 case kHijing_per1:
537 {
538 comment = comment.Append("HIJING per1");
539 AliGenHijing *gener = HijingStandard();
540 // impact parameter range
541 gener->SetImpactParameterRange(5., 8.6);
542 gGener=gener;
543 }
544 break;
545 case kHijing_per2:
546 {
547 comment = comment.Append("HIJING per2");
548 AliGenHijing *gener = HijingStandard();
549 // impact parameter range
550 gener->SetImpactParameterRange(8.6, 11.2);
551 gGener=gener;
552 }
553 break;
554 case kHijing_per3:
555 {
556 comment = comment.Append("HIJING per3");
557 AliGenHijing *gener = HijingStandard();
558 // impact parameter range
559 gener->SetImpactParameterRange(11.2, 13.2);
560 gGener=gener;
561 }
562 break;
563 case kHijing_per4:
564 {
565 comment = comment.Append("HIJING per4");
566 AliGenHijing *gener = HijingStandard();
567 // impact parameter range
568 gener->SetImpactParameterRange(13.2, 15.);
569 gGener=gener;
570 }
571 break;
572 case kHijing_per5:
573 {
574 comment = comment.Append("HIJING per5");
575 AliGenHijing *gener = HijingStandard();
576 // impact parameter range
577 gener->SetImpactParameterRange(15., 100.);
578 gGener=gener;
579 }
580 break;
581 //
582 // Jet-Jet
583 //
584 case kHijing_jj25:
585 {
586 comment = comment.Append("HIJING Jet 25 GeV");
587 AliGenHijing *gener = HijingStandard();
588 // impact parameter range
589 gener->SetImpactParameterRange(0., 5.);
590 // trigger
591 gener->SetTrigger(1);
592 gener->SetPtJet(25.);
593 gener->SetRadiation(isw);
594 gener->SetSimpleJets(!isw);
595 gener->SetJetEtaRange(-0.3,0.3);
596 gener->SetJetPhiRange(75., 165.);
597 gGener=gener;
598 }
599 break;
600
601 case kHijing_jj50:
602 {
603 comment = comment.Append("HIJING Jet 50 GeV");
604 AliGenHijing *gener = HijingStandard();
605 // impact parameter range
606 gener->SetImpactParameterRange(0., 5.);
607 // trigger
608 gener->SetTrigger(1);
609 gener->SetPtJet(50.);
610 gener->SetRadiation(isw);
611 gener->SetSimpleJets(!isw);
612 gener->SetJetEtaRange(-0.3,0.3);
613 gener->SetJetPhiRange(75., 165.);
614 gGener=gener;
615 }
616 break;
617
618 case kHijing_jj75:
619 {
620 comment = 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->SetPtJet(75.);
627 gener->SetRadiation(isw);
628 gener->SetSimpleJets(!isw);
629 gener->SetJetEtaRange(-0.3,0.3);
630 gener->SetJetPhiRange(75., 165.);
631 gGener=gener;
632 }
633 break;
634
635 case kHijing_jj100:
636 {
637 comment = comment.Append("HIJING Jet 100 GeV");
638 AliGenHijing *gener = HijingStandard();
639 // impact parameter range
640 gener->SetImpactParameterRange(0., 5.);
641 // trigger
642 gener->SetTrigger(1);
643 gener->SetPtJet(100.);
644 gener->SetRadiation(isw);
645 gener->SetSimpleJets(!isw);
646 gener->SetJetEtaRange(-0.3,0.3);
647 gener->SetJetPhiRange(75., 165.);
648 gGener=gener;
649 }
650 break;
651
652 case kHijing_jj200:
653 {
654 comment = comment.Append("HIJING Jet 200 GeV");
655 AliGenHijing *gener = HijingStandard();
656 // impact parameter range
657 gener->SetImpactParameterRange(0., 5.);
658 // trigger
659 gener->SetTrigger(1);
660 gener->SetPtJet(200.);
661 gener->SetRadiation(isw);
662 gener->SetSimpleJets(!isw);
663 gener->SetJetEtaRange(-0.3,0.3);
664 gener->SetJetPhiRange(75., 165.);
665 gGener=gener;
666 }
667 break;
668 //
669 // Gamma-Jet
670 //
671 case kHijing_gj25:
672 {
673 comment = comment.Append("HIJING Gamma 25 GeV");
674 AliGenHijing *gener = HijingStandard();
675 // impact parameter range
676 gener->SetImpactParameterRange(0., 5.);
677 // trigger
678 gener->SetTrigger(2);
679 gener->SetPtJet(25.);
680 gener->SetRadiation(isw);
681 gener->SetSimpleJets(!isw);
682 gener->SetJetEtaRange(-0.12, 0.12);
683 gener->SetJetPhiRange(220., 320.);
684 gGener=gener;
685 }
686 break;
687
688 case kHijing_gj50:
689 {
690 comment = comment.Append("HIJING Gamma 50 GeV");
691 AliGenHijing *gener = HijingStandard();
692 // impact parameter range
693 gener->SetImpactParameterRange(0., 5.);
694 // trigger
695 gener->SetTrigger(2);
696 gener->SetPtJet(50.);
697 gener->SetRadiation(isw);
698 gener->SetSimpleJets(!isw);
699 gener->SetJetEtaRange(-0.12, 0.12);
700 gener->SetJetPhiRange(220., 320.);
701 gGener=gener;
702 }
703 break;
704
705 case kHijing_gj75:
706 {
707 comment = comment.Append("HIJING Gamma 75 GeV");
708 AliGenHijing *gener = HijingStandard();
709 // impact parameter range
710 gener->SetImpactParameterRange(0., 5.);
711 // trigger
712 gener->SetTrigger(2);
713 gener->SetPtJet(75.);
714 gener->SetRadiation(isw);
715 gener->SetSimpleJets(!isw);
716 gener->SetJetEtaRange(-0.12, 0.12);
717 gener->SetJetPhiRange(220., 320.);
718 gGener=gener;
719 }
720 break;
721
722 case kHijing_gj100:
723 {
724 comment = comment.Append("HIJING Gamma 100 GeV");
725 AliGenHijing *gener = HijingStandard();
726 // impact parameter range
727 gener->SetImpactParameterRange(0., 5.);
728 // trigger
729 gener->SetTrigger(2);
730 gener->SetPtJet(100.);
731 gener->SetRadiation(isw);
732 gener->SetSimpleJets(!isw);
733 gener->SetJetEtaRange(-0.12, 0.12);
734 gener->SetJetPhiRange(220., 320.);
735 gGener=gener;
736 }
737 break;
738
739 case kHijing_gj200:
740 {
741 comment = comment.Append("HIJING Gamma 200 GeV");
742 AliGenHijing *gener = HijingStandard();
743 // impact parameter range
744 gener->SetImpactParameterRange(0., 5.);
745 // trigger
746 gener->SetTrigger(2);
747 gener->SetPtJet(200.);
748 gener->SetRadiation(isw);
749 gener->SetSimpleJets(!isw);
750 gener->SetJetEtaRange(-0.12, 0.12);
751 gener->SetJetPhiRange(220., 320.);
752 gGener=gener;
753 }
754 break;
755 case kHijing_pA:
756 {
757 comment = comment.Append("HIJING pA");
758
759 AliGenCocktail *gener = new AliGenCocktail();
760
761 AliGenHijing *hijing = new AliGenHijing(-1);
762 // centre of mass energy
763 hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
764 // impact parameter range
765 hijing->SetImpactParameterRange(0., 15.);
766 // reference frame
767 hijing->SetReferenceFrame("CMS");
768 hijing->SetBoostLHC(1);
769 // projectile
770 hijing->SetProjectile("P", 1, 1);
771 hijing->SetTarget ("A", 208, 82);
772 // tell hijing to keep the full parent child chain
773 hijing->KeepFullEvent();
774 // enable jet quenching
775 hijing->SetJetQuenching(0);
776 // enable shadowing
777 hijing->SetShadowing(1);
778 // Don't track spectators
779 hijing->SetSpectators(0);
780 // kinematic selection
781 hijing->SetSelectAll(0);
782 //
783 AliGenSlowNucleons* gray = new AliGenSlowNucleons(1);
784 AliSlowNucleonModel* model = new AliSlowNucleonModelExp();
785 gray->SetSlowNucleonModel(model);
786 gray->SetDebug(1);
787 gener->AddGenerator(hijing,"Hijing pPb", 1);
788 gener->AddGenerator(gray, "Gray Particles",1);
789 gGener=gener;
790 }
791 break;
792 case kPythia6:
793 {
794 comment = comment.Append(":Pythia p-p @ 14 TeV");
795 AliGenPythia *gener = new AliGenPythia(-1);
796 gener->SetMomentumRange(0,999999);
797 gener->SetThetaRange(0., 180.);
798 gener->SetYRange(-12,12);
799 gener->SetPtRange(0,1000);
800 gener->SetProcess(kPyMb);
801 gener->SetEnergyCMS(14000.);
802 gener->SetProjectile("p", 1, 1) ;
803 gener->SetTarget("p", 1, 1) ;
804 gGener=gener;
805 }
806 break;
807 case kPythia6Jets20_24:
808 {
809 comment = comment.Append(":Pythia jets 20-24 GeV @ 5.5 TeV");
810 AliGenPythia * gener = new AliGenPythia(-1);
811 gener->SetEnergyCMS(5500.);// Centre of mass energy
812 gener->SetProcess(kPyJets);// Process type
813 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
814 gener->SetJetPhiRange(0., 360.);
815 gener->SetJetEtRange(10., 1000.);
816 gener->SetGluonRadiation(1,1);
817 // gener->SetPtKick(0.);
818 // Structure function
819 gener->SetStrucFunc(kCTEQ4L);
820 gener->SetPtHard(20., 24.);// Pt transfer of the hard scattering
821 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
822 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
823 gener->SetProjectile("p", 1, 1) ;
824 gener->SetTarget("p", 1, 1) ;
825 gGener=gener;
826 }
827 break;
828 case kPythia6Jets24_29:
829 {
830 comment = comment.Append(":Pythia jets 24-29 GeV @ 5.5 TeV");
831 AliGenPythia * gener = new AliGenPythia(-1);
832 gener->SetEnergyCMS(5500.);// Centre of mass energy
833 gener->SetProcess(kPyJets);// Process type
834 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
835 gener->SetJetPhiRange(0., 360.);
836 gener->SetJetEtRange(10., 1000.);
837 gener->SetGluonRadiation(1,1);
838 // gener->SetPtKick(0.);
839 // Structure function
840 gener->SetStrucFunc(kCTEQ4L);
841 gener->SetPtHard(24., 29.);// Pt transfer of the hard scattering
842 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
843 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
844 gener->SetProjectile("p", 1, 1) ;
845 gener->SetTarget("p", 1, 1) ;
846 gGener=gener;
847 }
848 break;
849 case kPythia6Jets29_35:
850 {
851 comment = comment.Append(":Pythia jets 29-35 GeV @ 5.5 TeV");
852 AliGenPythia * gener = new AliGenPythia(-1);
853 gener->SetEnergyCMS(5500.);// Centre of mass energy
854 gener->SetProcess(kPyJets);// Process type
855 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
856 gener->SetJetPhiRange(0., 360.);
857 gener->SetJetEtRange(10., 1000.);
858 gener->SetGluonRadiation(1,1);
859 // gener->SetPtKick(0.);
860 // Structure function
861 gener->SetStrucFunc(kCTEQ4L);
862 gener->SetPtHard(29., 35.);// Pt transfer of the hard scattering
863 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
864 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
865 gener->SetProjectile("p", 1, 1) ;
866 gener->SetTarget("p", 1, 1) ;
867 gGener=gener;
868 }
869 break;
870 case kPythia6Jets35_42:
871 {
872 comment = comment.Append(":Pythia jets 35-42 GeV @ 5.5 TeV");
873 AliGenPythia * gener = new AliGenPythia(-1);
874 gener->SetEnergyCMS(5500.);// Centre of mass energy
875 gener->SetProcess(kPyJets);// Process type
876 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
877 gener->SetJetPhiRange(0., 360.);
878 gener->SetJetEtRange(10., 1000.);
879 gener->SetGluonRadiation(1,1);
880 // gener->SetPtKick(0.);
881 // Structure function
882 gener->SetStrucFunc(kCTEQ4L);
883 gener->SetPtHard(35., 42.);// Pt transfer of the hard scattering
884 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
885 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
886 gener->SetProjectile("p", 1, 1) ;
887 gener->SetTarget("p", 1, 1) ;
888 gGener=gener;
889 }
890 break;
891 case kPythia6Jets42_50:
892 {
893 comment = comment.Append(":Pythia jets 42-50 GeV @ 5.5 TeV");
894 AliGenPythia * gener = new AliGenPythia(-1);
895 gener->SetEnergyCMS(5500.);// Centre of mass energy
896 gener->SetProcess(kPyJets);// Process type
897 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
898 gener->SetJetPhiRange(0., 360.);
899 gener->SetJetEtRange(10., 1000.);
900 gener->SetGluonRadiation(1,1);
901 // gener->SetPtKick(0.);
902 // Structure function
903 gener->SetStrucFunc(kCTEQ4L);
904 gener->SetPtHard(42., 50.);// Pt transfer of the hard scattering
905 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
906 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
907 gener->SetProjectile("p", 1, 1) ;
908 gener->SetTarget("p", 1, 1) ;
909 gGener=gener;
910 }
911 break;
912 case kPythia6Jets50_60:
913 {
914 comment = comment.Append(":Pythia jets 50-60 GeV @ 5.5 TeV");
915 AliGenPythia * gener = new AliGenPythia(-1);
916 gener->SetEnergyCMS(5500.);// Centre of mass energy
917 gener->SetProcess(kPyJets);// Process type
918 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
919 gener->SetJetPhiRange(0., 360.);
920 gener->SetJetEtRange(10., 1000.);
921 gener->SetGluonRadiation(1,1);
922 // gener->SetPtKick(0.);
923 // Structure function
924 gener->SetStrucFunc(kCTEQ4L);
925 gener->SetPtHard(50., 60.);// Pt transfer of the hard scattering
926 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
927 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
928 gGener=gener;
929 }
930 break;
931 case kPythia6Jets60_72:
932 {
933 comment = comment.Append(":Pythia jets 60-72 GeV @ 5.5 TeV");
934 AliGenPythia * gener = new AliGenPythia(-1);
935 gener->SetEnergyCMS(5500.);// Centre of mass energy
936 gener->SetProcess(kPyJets);// Process type
937 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
938 gener->SetJetPhiRange(0., 360.);
939 gener->SetJetEtRange(10., 1000.);
940 gener->SetGluonRadiation(1,1);
941 // gener->SetPtKick(0.);
942 // Structure function
943 gener->SetStrucFunc(kCTEQ4L);
944 gener->SetPtHard(60., 72.);// Pt transfer of the hard scattering
945 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
946 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
947 gener->SetProjectile("p", 1, 1) ;
948 gener->SetTarget("p", 1, 1) ;
949 gGener=gener;
950 }
951 break;
952 case kPythia6Jets72_86:
953 {
954 comment = comment.Append(":Pythia jets 72-86 GeV @ 5.5 TeV");
955 AliGenPythia * gener = new AliGenPythia(-1);
956 gener->SetEnergyCMS(5500.);// Centre of mass energy
957 gener->SetProcess(kPyJets);// Process type
958 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
959 gener->SetJetPhiRange(0., 360.);
960 gener->SetJetEtRange(10., 1000.);
961 gener->SetGluonRadiation(1,1);
962 // gener->SetPtKick(0.);
963 // Structure function
964 gener->SetStrucFunc(kCTEQ4L);
965 gener->SetPtHard(72., 86.);// Pt transfer of the hard scattering
966 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
967 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
968 gener->SetProjectile("p", 1, 1) ;
969 gener->SetTarget("p", 1, 1) ;
970 gGener=gener;
971 }
972 break;
973 case kPythia6Jets86_104:
974 {
975 comment = comment.Append(":Pythia jets 86-104 GeV @ 5.5 TeV");
976 AliGenPythia * gener = new AliGenPythia(-1);
977 gener->SetEnergyCMS(5500.);// Centre of mass energy
978 gener->SetProcess(kPyJets);// Process type
979 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
980 gener->SetJetPhiRange(0., 360.);
981 gener->SetJetEtRange(10., 1000.);
982 gener->SetGluonRadiation(1,1);
983 // gener->SetPtKick(0.);
984 // Structure function
985 gener->SetStrucFunc(kCTEQ4L);
986 gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
987 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
988 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
989 gener->SetProjectile("p", 1, 1) ;
990 gener->SetTarget("p", 1, 1) ;
991 gGener=gener;
992 }
993 break;
994 case kPythia6Jets104_125:
995 {
996 comment = comment.Append(":Pythia jets 105-125 GeV @ 5.5 TeV");
997 AliGenPythia * gener = new AliGenPythia(-1);
998 gener->SetEnergyCMS(5500.);// Centre of mass energy
999 gener->SetProcess(kPyJets);// Process type
1000 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1001 gener->SetJetPhiRange(0., 360.);
1002 gener->SetJetEtRange(10., 1000.);
1003 gener->SetGluonRadiation(1,1);
1004 // gener->SetPtKick(0.);
1005 // Structure function
1006 gener->SetStrucFunc(kCTEQ4L);
1007 gener->SetPtHard(104., 125.);// Pt transfer of the hard scattering
1008 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1009 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1010 gener->SetProjectile("p", 1, 1) ;
1011 gener->SetTarget("p", 1, 1) ;
1012 gGener=gener;
1013 }
1014 break;
1015 case kPythia6Jets125_150:
1016 {
1017 comment = comment.Append(":Pythia jets 125-150 GeV @ 5.5 TeV");
1018 AliGenPythia * gener = new AliGenPythia(-1);
1019 gener->SetEnergyCMS(5500.);// Centre of mass energy
1020 gener->SetProcess(kPyJets);// Process type
1021 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1022 gener->SetJetPhiRange(0., 360.);
1023 gener->SetJetEtRange(10., 1000.);
1024 gener->SetGluonRadiation(1,1);
1025 // gener->SetPtKick(0.);
1026 // Structure function
1027 gener->SetStrucFunc(kCTEQ4L);
1028 gener->SetPtHard(125., 150.);// Pt transfer of the hard scattering
1029 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1030 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1031 gGener=gener;
1032 }
1033 break;
1034 case kPythia6Jets150_180:
1035 {
1036 comment = comment.Append(":Pythia jets 150-180 GeV @ 5.5 TeV");
1037 AliGenPythia * gener = new AliGenPythia(-1);
1038 gener->SetEnergyCMS(5500.);// Centre of mass energy
1039 gener->SetProcess(kPyJets);// Process type
1040 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1041 gener->SetJetPhiRange(0., 360.);
1042 gener->SetJetEtRange(10., 1000.);
1043 gener->SetGluonRadiation(1,1);
1044 // gener->SetPtKick(0.);
1045 // Structure function
1046 gener->SetStrucFunc(kCTEQ4L);
1047 gener->SetPtHard(150., 180.);// Pt transfer of the hard scattering
1048 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1049 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1050 gener->SetProjectile("p", 1, 1) ;
1051 gener->SetTarget("p", 1, 1) ;
1052 gGener=gener;
1053 }
1054 break;
1055 case kD0PbPb5500:
1056 {
1057 comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
1058 AliGenPythia * gener = new AliGenPythia(10);
1059 gener->SetProcess(kPyD0PbPbMNR);
1060 gener->SetStrucFunc(kCTEQ4L);
1061 gener->SetPtHard(2.1,-1.0);
1062 gener->SetEnergyCMS(5500.);
1063 gener->SetNuclei(208,208);
1064 gener->SetForceDecay(kHadronicD);
1065 gener->SetYRange(-2,2);
1066 gener->SetFeedDownHigherFamily(kFALSE);
1067 gener->SetStackFillOpt(AliGenPythia::kParentSelection);
1068 gener->SetCountMode(AliGenPythia::kCountParents);
1069 gener->SetProjectile("A", 208, 82) ;
1070 gener->SetTarget("A", 208, 82) ;
1071 gGener=gener;
1072 }
1073 break;
1074 case kCharmSemiElPbPb5500:
1075 {
1076 comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
1077 AliGenPythia * gener = new AliGenPythia(10);
1078 gener->SetProcess(kPyCharmPbPbMNR);
1079 gener->SetStrucFunc(kCTEQ4L);
1080 gener->SetPtHard(2.1,-1.0);
1081 gener->SetEnergyCMS(5500.);
1082 gener->SetNuclei(208,208);
1083 gener->SetForceDecay(kSemiElectronic);
1084 gener->SetYRange(-2,2);
1085 gener->SetFeedDownHigherFamily(kFALSE);
1086 gener->SetCountMode(AliGenPythia::kCountParents);
1087 gener->SetProjectile("A", 208, 82) ;
1088 gener->SetTarget("A", 208, 82) ;
1089 gGener=gener;
1090 }
1091 break;
1092 case kBeautySemiElPbPb5500:
1093 {
1094 comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
1095 AliGenPythia *gener = new AliGenPythia(10);
1096 gener->SetProcess(kPyBeautyPbPbMNR);
1097 gener->SetStrucFunc(kCTEQ4L);
1098 gener->SetPtHard(2.75,-1.0);
1099 gener->SetEnergyCMS(5500.);
1100 gener->SetNuclei(208,208);
1101 gener->SetForceDecay(kSemiElectronic);
1102 gener->SetYRange(-2,2);
1103 gener->SetFeedDownHigherFamily(kFALSE);
1104 gener->SetCountMode(AliGenPythia::kCountParents);
1105 gener->SetProjectile("A", 208, 82) ;
1106 gener->SetTarget("A", 208, 82) ;
1107 gGener=gener;
1108 }
1109 break;
1110 case kCocktailTRD:
1111 {
1112 comment = comment.Append(" Cocktail for TRD at 5.5 TeV");
1113 AliGenCocktail *gener = new AliGenCocktail();
1114
1115 AliGenParam *phi = new AliGenParam(10,
1116 new AliGenMUONlib(),
1117 AliGenMUONlib::kPhi,
1118 "Vogt PbPb");
1119
1120 phi->SetPtRange(0, 100);
1121 phi->SetYRange(-1., +1.);
1122 phi->SetForceDecay(kDiElectron);
1123
1124 AliGenParam *omega = new AliGenParam(10,
1125 new AliGenMUONlib(),
1126 AliGenMUONlib::kOmega,
1127 "Vogt PbPb");
1128
1129 omega->SetPtRange(0, 100);
1130 omega->SetYRange(-1., +1.);
1131 omega->SetForceDecay(kDiElectron);
1132
1133 AliGenParam *jpsi = new AliGenParam(10,
1134 new AliGenMUONlib(),
1135 AliGenMUONlib::kJpsiFamily,
1136 "Vogt PbPb");
1137
1138 jpsi->SetPtRange(0, 100);
1139 jpsi->SetYRange(-1., +1.);
1140 jpsi->SetForceDecay(kDiElectron);
1141
1142 AliGenParam *ups = new AliGenParam(10,
1143 new AliGenMUONlib(),
1144 AliGenMUONlib::kUpsilonFamily,
1145 "Vogt PbPb");
1146 ups->SetPtRange(0, 100);
1147 ups->SetYRange(-1., +1.);
1148 ups->SetForceDecay(kDiElectron);
1149
1150 AliGenParam *charm = new AliGenParam(10,
1151 new AliGenMUONlib(),
1152 AliGenMUONlib::kCharm,
1153 "central");
1154 charm->SetPtRange(0, 100);
1155 charm->SetYRange(-1.5, +1.5);
1156 charm->SetForceDecay(kSemiElectronic);
1157
1158
1159 AliGenParam *beauty = new AliGenParam(10,
1160 new AliGenMUONlib(),
1161 AliGenMUONlib::kBeauty,
1162 "central");
1163 beauty->SetPtRange(0, 100);
1164 beauty->SetYRange(-1.5, +1.5);
1165 beauty->SetForceDecay(kSemiElectronic);
1166
1167 AliGenParam *beautyJ = new AliGenParam(10,
1168 new AliGenMUONlib(),
1169 AliGenMUONlib::kBeauty,
1170 "central");
1171 beautyJ->SetPtRange(0, 100);
1172 beautyJ->SetYRange(-1.5, +1.5);
1173 beautyJ->SetForceDecay(kBJpsiDiElectron);
1174
1175 gener->AddGenerator(phi,"Phi",1);
1176 gener->AddGenerator(omega,"Omega",1);
1177 gener->AddGenerator(jpsi,"J/psi",1);
1178 gener->AddGenerator(ups,"Upsilon",1);
1179 gener->AddGenerator(charm,"Charm",1);
1180 gener->AddGenerator(beauty,"Beauty",1);
1181 gener->AddGenerator(beautyJ,"J/Psi from Beauty",1);
1182 gener->SetProjectile("A", 208, 82) ;
1183 gener->SetTarget("A", 208, 82) ;
1184 gGener=gener;
1185 }
1186 break;
1187 case kPyJJ:
1188 {
1189 comment = comment.Append(" Jet-jet at 5.5 TeV");
1190 AliGenPythia *gener = new AliGenPythia(-1);
1191 gener->SetEnergyCMS(5500.);
1192 gener->SetProcess(kPyJets);
1193 Double_t ptHardMin=10.0, ptHardMax=-1.0;
1194 gener->SetPtHard(ptHardMin,ptHardMax);
1195 gener->SetYHard(-0.7,0.7);
1196 gener->SetJetEtaRange(-0.2,0.2);
1197 gener->SetEventListRange(0,1);
1198 gener->SetProjectile("p", 1, 1) ;
1199 gener->SetTarget("p", 1, 1) ;
1200 gGener=gener;
1201 }
1202 break;
1203 case kPyGJ:
1204 {
1205 comment = comment.Append(" Gamma-jet at 5.5 TeV");
1206 AliGenPythia *gener = new AliGenPythia(-1);
1207 gener->SetEnergyCMS(5500.);
1208 gener->SetProcess(kPyDirectGamma);
1209 Double_t ptHardMin=10.0, ptHardMax=-1.0;
1210 gener->SetPtHard(ptHardMin,ptHardMax);
1211 gener->SetYHard(-1.0,1.0);
1212 gener->SetGammaEtaRange(-0.13,0.13);
1213 gener->SetGammaPhiRange(210.,330.);
1214 gener->SetEventListRange(0,1);
1215 gener->SetProjectile("p", 1, 1) ;
1216 gener->SetTarget("p", 1, 1) ;
1217 gGener=gener;
1218 }
1219 break;
1220 case kMuonCocktailCent1:
1221 {
1222 comment = comment.Append(" Muon Cocktail Cent1");
1223 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1224 gener->SetPtRange(0.4,100.); // Transverse momentum range
1225 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1226 gener->SetYRange(-4.0,-2.4);
1227 gener->SetMuonPtCut(0.8);
1228 gener->SetMuonThetaCut(171.,178.);
1229 gener->SetMuonMultiplicity(2);
1230 gener->SetImpactParameterRange(0.,5.); //Centrality class Cent1 for PDC04
1231 gGener=gener;
1232 gener->SetProjectile("A", 208, 82) ;
1233 gener->SetTarget("A", 208, 82) ;
1234 }
1235 break;
1236 case kMuonCocktailPer1:
1237 {
1238 comment = comment.Append(" Muon Cocktail Per1");
1239 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1240 gener->SetPtRange(0.0,100.); // Transverse momentum range
1241 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1242 gener->SetYRange(-4.0,-2.4);
1243 gener->SetMuonPtCut(0.8);
1244 gener->SetMuonThetaCut(171.,178.);
1245 gener->SetMuonMultiplicity(2);
1246 gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1247 gener->SetProjectile("A", 208, 82) ;
1248 gener->SetTarget("A", 208, 82) ;
1249 gGener=gener;
1250 }
1251 break;
1252 case kMuonCocktailPer4:
1253 {
1254 comment = comment.Append(" Muon Cocktail Per4");
1255 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1256 gener->SetPtRange(0.0,100.); // Transverse momentum range
1257 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1258 gener->SetYRange(-4.0,-2.4);
1259 gener->SetMuonPtCut(0.8);
1260 gener->SetMuonThetaCut(171.,178.);
1261 gener->SetMuonMultiplicity(2);
1262 gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1263 gener->SetProjectile("A", 208, 82) ;
1264 gener->SetTarget("A", 208, 82) ;
1265 gGener=gener;
1266 }
1267 break;
1268 case kMuonCocktailCent1HighPt:
1269 {
1270 comment = comment.Append(" Muon Cocktail HighPt Cent1");
1271 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1272 gener->SetPtRange(0.0,100.); // Transverse momentum range
1273 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1274 gener->SetYRange(-4.0,-2.4);
1275 gener->SetMuonPtCut(2.5);
1276 gener->SetMuonThetaCut(171.,178.);
1277 gener->SetMuonMultiplicity(2);
1278 gener->SetImpactParameterRange(0.,5.); //Centrality class Cent1 for PDC04
1279 gener->SetProjectile("A", 208, 82) ;
1280 gener->SetTarget("A", 208, 82) ;
1281 gGener=gener;
1282 }
1283 break;
1284 case kMuonCocktailPer1HighPt :
1285 {
1286 comment = comment.Append(" Muon Cocktail HighPt Per1");
1287 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1288 gener->SetPtRange(0.0,100.); // Transverse momentum range
1289 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1290 gener->SetYRange(-4.0,-2.4);
1291 gener->SetMuonPtCut(2.5);
1292 gener->SetMuonThetaCut(171.,178.);
1293 gener->SetMuonMultiplicity(2);
1294 gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1295 gener->SetProjectile("A", 208, 82) ;
1296 gener->SetTarget("A", 208, 82) ;
1297 gGener=gener;
1298 }
1299 break;
1300 case kMuonCocktailPer4HighPt:
1301 {
1302 comment = comment.Append(" Muon Cocktail HighPt Per4");
1303 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1304 gener->SetPtRange(0.0,100.); // Transverse momentum range
1305 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1306 gener->SetYRange(-4.0,-2.4);
1307 gener->SetMuonPtCut(2.5);
1308 gener->SetMuonThetaCut(171.,178.);
1309 gener->SetMuonMultiplicity(2);
1310 gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1311 gener->SetProjectile("A", 208, 82) ;
1312 gener->SetTarget("A", 208, 82) ;
1313 gGener=gener;
1314 }
1315 break;
1316 case kMuonCocktailCent1Single:
1317 {
1318 comment = comment.Append(" Muon Cocktail Single Cent1");
1319 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1320 gener->SetPtRange(0.0,100.); // Transverse momentum range
1321 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1322 gener->SetYRange(-4.0,-2.4);
1323 gener->SetMuonPtCut(0.8);
1324 gener->SetMuonThetaCut(171.,178.);
1325 gener->SetMuonMultiplicity(1);
1326 gener->SetImpactParameterRange(0.,5.); //Centrality class Cent1 for PDC04
1327 gener->SetProjectile("A", 208, 82) ;
1328 gener->SetTarget("A", 208, 82) ;
1329 gGener=gener;
1330 }
1331 break;
1332 case kMuonCocktailPer1Single :
1333 {
1334 comment = comment.Append(" Muon Cocktail Single Per1");
1335 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1336 gener->SetPtRange(0.0,100.); // Transverse momentum range
1337 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1338 gener->SetYRange(-4.0,-2.4);
1339 gener->SetMuonPtCut(0.8);
1340 gener->SetMuonThetaCut(171.,178.);
1341 gener->SetMuonMultiplicity(1);
1342 gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1343 gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1344 gener->SetProjectile("A", 208, 82) ;
1345 gener->SetTarget("A", 208, 82) ;
1346 gGener=gener;
1347 }
1348 break;
1349 case kMuonCocktailPer4Single:
1350 {
1351 comment = comment.Append(" Muon Cocktail Single Per4");
1352 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1353 gener->SetPtRange(0.0,100.); // Transverse momentum range
1354 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1355 gener->SetYRange(-4.0,-2.4);
1356 gener->SetMuonPtCut(0.8);
1357 gener->SetMuonThetaCut(171.,178.);
1358 gener->SetMuonMultiplicity(1);
1359 gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1360 gener->SetProjectile("A", 208, 82) ;
1361 gener->SetTarget("A", 208, 82) ;
1362 gGener=gener;
1363 }
1364 break;
1365 case kFlow_2_2000:
1366 {
1367 comment = comment.Append(" Flow with dN/deta = 2000, vn = 2%");
1368 gGener = GeVSimStandard(2000., 2.);
1369 }
1370 break;
1371
1372 case kFlow_10_2000:
1373 {
1374 comment = comment.Append(" Flow with dN/deta = 2000, vn = 10%");
1375 gGener = GeVSimStandard(2000., 10.);
1376 }
1377 break;
1378
1379 case kFlow_6_2000:
1380 {
1381 comment = comment.Append(" Flow with dN/deta = 2000, vn = 6%");
1382 gGener = GeVSimStandard(2000., 6.);
1383 }
1384 break;
1385
1386 case kFlow_6_5000:
1387 {
1388 comment = comment.Append(" Flow with dN/deta = 5000, vn = 6%");
1389 gGener = GeVSimStandard(5000., 6.);
1390 }
1391 break;
1392 case kHIJINGplus:
1393 {
1394 //
1395 // The cocktail
1396 AliGenCocktail *gener = new AliGenCocktail();
1397
1398 //
1399 // Charm production by Pythia
1400 AliGenPythia * genpyc = new AliGenPythia(230);
1401 genpyc->SetProcess(kPyCharmPbPbMNR);
1402 genpyc->SetStrucFunc(kCTEQ4L);
1403 genpyc->SetPtHard(2.1,-1.0);
1404 genpyc->SetEnergyCMS(5500.);
1405 genpyc->SetNuclei(208,208);
1406 genpyc->SetYRange(-999,999);
1407 genpyc->SetForceDecay(kAll);
1408 genpyc->SetFeedDownHigherFamily(kFALSE);
1409 genpyc->SetCountMode(AliGenPythia::kCountParents);
1410 //
1411 // Beauty production by Pythia
1412 AliGenPythia * genpyb = new AliGenPythia(9);
1413 genpyb->SetProcess(kPyBeautyPbPbMNR);
1414 genpyb->SetStrucFunc(kCTEQ4L);
1415 genpyb->SetPtHard(2.75,-1.0);
1416 genpyb->SetEnergyCMS(5500.);
1417 genpyb->SetNuclei(208,208);
1418 genpyb->SetYRange(-999,999);
1419 genpyb->SetForceDecay(kAll);
1420 genpyb->SetFeedDownHigherFamily(kFALSE);
1421 genpyb->SetCountMode(AliGenPythia::kCountParents);
1422 //
1423 // Hyperons
1424 //
1425 AliGenSTRANGElib *lib = new AliGenSTRANGElib();
1426 Int_t particle;
1427 // Xi
1428 particle = kXiMinus;
1429 AliGenParam *genXi = new AliGenParam(16,particle,lib->GetPt(particle),lib->GetY(particle),lib->GetIp(particle));
1430 genXi->SetPtRange(0., 12.);
1431 genXi->SetYRange(-1.1, 1.1);
1432 genXi->SetForceDecay(kNoDecay);
1433
1434 //
1435 // Omega
1436 particle = kOmegaMinus;
1437 AliGenParam *genOmega = new AliGenParam(10,particle,lib->GetPt(particle),lib->GetY(particle),lib->GetIp(particle));
1438 genOmega->SetPtRange(0, 12.);
1439 genOmega->SetYRange(-1.1, 1.1);
1440 genOmega->SetForceDecay(kNoDecay);
1441
1442 //
1443 // Central Hijing
1444 AliGenHijing *genHi = HijingStandard();
1445 genHi->SwitchOffHeavyQuarks(kTRUE);
1446 genHi->SetImpactParameterRange(0.,5.);
1447 //
1448 // Add everything to the cocktail and shake ...
1449 gener->AddGenerator(genHi, "Hijing cent1", 1);
1450 gener->AddGenerator(genpyc, "Extra charm", 1);
1451 gener->AddGenerator(genpyb, "Extra beauty", 1);
1452 gener->AddGenerator(genXi, "Xi" , 1);
1453 gener->AddGenerator(genOmega, "Omega", 1);
1454 gener->SetProjectile("A", 208, 82) ;
1455 gener->SetTarget("A", 208, 82) ;
1456 gGener = gener;
1457 }
1458 break;
1459 default: break;
1460 }
1461
1462 return gGener;
1463}
1464
1465AliGenHijing* HijingStandard()
1466{
1467 AliGenHijing *gener = new AliGenHijing(-1);
1468// centre of mass energy
1469 gener->SetEnergyCMS(5500.);
1470// reference frame
1471 gener->SetReferenceFrame("CMS");
1472// projectile
1473 gener->SetProjectile("A", 208, 82);
1474 gener->SetTarget ("A", 208, 82);
1475// tell hijing to keep the full parent child chain
1476 gener->KeepFullEvent();
1477// enable jet quenching
1478 gener->SetJetQuenching(1);
1479// enable shadowing
1480 gener->SetShadowing(1);
1481// neutral pion and heavy particle decays switched off
1482 gener->SetDecaysOff(1);
1483// Don't track spectators
1484 gener->SetSpectators(0);
1485// kinematic selection
1486 gener->SetSelectAll(0);
1487 return gener;
1488}
1489
1490AliGenGeVSim* GeVSimStandard(Float_t mult, Float_t vn)
1491{
1492 AliGenGeVSim* gener = new AliGenGeVSim(0);
1493//
1494// Mult is the number of charged particles in |eta| < 0.5
1495// Vn is in (%)
1496//
1497// Sigma of the Gaussian dN/deta
1498 Float_t sigma_eta = 2.75;
1499//
1500// Maximum eta
1501 Float_t etamax = 7.00;
1502//
1503//
1504// Scale from multiplicity in |eta| < 0.5 to |eta| < |etamax|
1505 Float_t mm = mult * (TMath::Erf(etamax/sigma_eta/sqrt(2.)) / TMath::Erf(0.5/sigma_eta/sqrt(2.)));
1506//
1507// Scale from charged to total multiplicity
1508//
1509 mm *= 1.587;
1510//
1511// Vn
1512 vn /= 100.;
1513//
1514// Define particles
1515//
1516//
1517// 78% Pions (26% pi+, 26% pi-, 26% p0) T = 250 MeV
1518 AliGeVSimParticle *pp = new AliGeVSimParticle(kPiPlus, 1, 0.26 * mm, 0.25, sigma_eta) ;
1519 AliGeVSimParticle *pm = new AliGeVSimParticle(kPiMinus, 1, 0.26 * mm, 0.25, sigma_eta) ;
1520 AliGeVSimParticle *p0 = new AliGeVSimParticle(kPi0, 1, 0.26 * mm, 0.25, sigma_eta) ;
1521//
1522// 12% Kaons (3% K0short, 3% K0long, 3% K+, 3% K-) T = 300 MeV
1523 AliGeVSimParticle *ks = new AliGeVSimParticle(kK0Short, 1, 0.03 * mm, 0.30, sigma_eta) ;
1524 AliGeVSimParticle *kl = new AliGeVSimParticle(kK0Long, 1, 0.03 * mm, 0.30, sigma_eta) ;
1525 AliGeVSimParticle *kp = new AliGeVSimParticle(kKPlus, 1, 0.03 * mm, 0.30, sigma_eta) ;
1526 AliGeVSimParticle *km = new AliGeVSimParticle(kKMinus, 1, 0.03 * mm, 0.30, sigma_eta) ;
1527//
1528// 10% Protons / Neutrons (5% Protons, 5% Neutrons) T = 250 MeV
1529 AliGeVSimParticle *pr = new AliGeVSimParticle(kProton, 1, 0.05 * mm, 0.25, sigma_eta) ;
1530 AliGeVSimParticle *ne = new AliGeVSimParticle(kNeutron, 1, 0.05 * mm, 0.25, sigma_eta) ;
1531//
1532// Set Elliptic Flow properties
1533
1534 Float_t pTsaturation = 2. ;
1535
1536 pp->SetEllipticParam(vn,pTsaturation,0.) ;
1537 pm->SetEllipticParam(vn,pTsaturation,0.) ;
1538 p0->SetEllipticParam(vn,pTsaturation,0.) ;
1539 pr->SetEllipticParam(vn,pTsaturation,0.) ;
1540 ne->SetEllipticParam(vn,pTsaturation,0.) ;
1541 ks->SetEllipticParam(vn,pTsaturation,0.) ;
1542 kl->SetEllipticParam(vn,pTsaturation,0.) ;
1543 kp->SetEllipticParam(vn,pTsaturation,0.) ;
1544 km->SetEllipticParam(vn,pTsaturation,0.) ;
1545//
1546// Set Direct Flow properties
1547 pp->SetDirectedParam(vn,1.0,0.) ;
1548 pm->SetDirectedParam(vn,1.0,0.) ;
1549 p0->SetDirectedParam(vn,1.0,0.) ;
1550 pr->SetDirectedParam(vn,1.0,0.) ;
1551 ne->SetDirectedParam(vn,1.0,0.) ;
1552 ks->SetDirectedParam(vn,1.0,0.) ;
1553 kl->SetDirectedParam(vn,1.0,0.) ;
1554 kp->SetDirectedParam(vn,1.0,0.) ;
1555 km->SetDirectedParam(vn,1.0,0.) ;
1556//
1557// Add particles to the list
1558 gener->AddParticleType(pp) ;
1559 gener->AddParticleType(pm) ;
1560 gener->AddParticleType(p0) ;
1561 gener->AddParticleType(pr) ;
1562 gener->AddParticleType(ne) ;
1563 gener->AddParticleType(ks) ;
1564 gener->AddParticleType(kl) ;
1565 gener->AddParticleType(kp) ;
1566 gener->AddParticleType(km) ;
1567//
1568// Random Ev.Plane ----------------------------------
1569 TF1 *rpa = new TF1("gevsimPsiRndm","1", 0, 360);
1570// --------------------------------------------------
1571 gener->SetPtRange(0., 9.) ; // Use a resonable range! (used for bin size in numerical integration)
1572 gener->SetPhiRange(0, 360);
1573 //
1574 // Set pseudorapidity range
1575 Float_t thmin = EtaToTheta(+etamax);
1576 Float_t thmax = EtaToTheta(-etamax);
1577 gener->SetThetaRange(thmin,thmax);
1578 return gener;
1579}
1580
1581
1582
1583void ProcessEnvironmentVars()
1584{
1585 // Run type
1586 if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
1587 for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
1588 if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
1589 srun = (PprRun_t)iRun;
1590 cout<<"Run type set to "<<pprRunName[iRun]<<endl;
1591 }
1592 }
1593 }
1594
1595 // Random Number seed
1596 if (gSystem->Getenv("CONFIG_SEED")) {
1597 sseed = atoi(gSystem->Getenv("CONFIG_SEED"));
1598 }
1599}