]> git.uio.no Git - u/mrichter/AliRoot.git/blame - test/ppbench/Config.C
added a setter of fEventSpecie
[u/mrichter/AliRoot.git] / test / ppbench / Config.C
CommitLineData
27b7c1fb 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"
f7a1cc68 34#include "STEER/AliMagF.h"
27b7c1fb 35#include "STRUCT/AliBODY.h"
36#include "STRUCT/AliMAG.h"
96fb3353 37#include "STRUCT/AliABSOv3.h"
38#include "STRUCT/AliDIPOv3.h"
39#include "STRUCT/AliHALLv3.h"
27b7c1fb 40#include "STRUCT/AliFRAMEv2.h"
96fb3353 41#include "STRUCT/AliSHILv3.h"
42#include "STRUCT/AliPIPEv3.h"
99887042 43#include "ITS/AliITSv11Hybrid.h"
27b7c1fb 44#include "TPC/AliTPCv2.h"
0d7ef405 45#include "TOF/AliTOFv6T0.h"
99887042 46#include "HMPID/AliHMPIDv3.h"
cb14d625 47#include "ZDC/AliZDCv3.h"
27b7c1fb 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"
f7882672 55#include "ACORDE/AliACORDEv1.h"
27b7c1fb 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
27b7c1fb 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 = kPythia6;
120static PprRad_t srad = kGluonRadiation;
f7a1cc68 121static AliMagF::BMap_t smag = AliMagF::k5kG;
27b7c1fb 122static Int_t sseed = 12345; //Set 0 to use the current time
123static PprTrigConf_t strig = kDefaultPPTrig; // 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__)
f1e07bdc 150 gSystem->Load("liblhapdf");
151 gSystem->Load("libEGPythia6");
152 gSystem->Load("libpythia6");
153 gSystem->Load("libAliPythia6");
27b7c1fb 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
f7a1cc68 251 if (smag == AliMagF::k2kG) {
27b7c1fb 252 comment = comment.Append(" | L3 field 0.2 T");
f7a1cc68 253 } else if (smag == AliMagF::k5kG) {
27b7c1fb 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
f7a1cc68 270 TGeoGlobalMagField::Instance()->SetField(new AliMagF("Maps","Maps", 2, 1., 1., 10., AliMagF::k5kG));
271
27b7c1fb 272 rl->CdGAFile();
27b7c1fb 273//
274 Int_t iABSO = 1;
275 Int_t iDIPO = 1;
276 Int_t iFMD = 1;
277 Int_t iFRAME = 1;
278 Int_t iHALL = 1;
279 Int_t iITS = 1;
280 Int_t iMAG = 1;
281 Int_t iMUON = 1;
282 Int_t iPHOS = 1;
283 Int_t iPIPE = 1;
284 Int_t iPMD = 1;
285 Int_t iHMPID = 1;
286 Int_t iSHIL = 1;
287 Int_t iT0 = 1;
288 Int_t iTOF = 1;
289 Int_t iTPC = 1;
290 Int_t iTRD = 1;
291 Int_t iZDC = 1;
292 Int_t iEMCAL = 1;
293 Int_t iVZERO = 1;
294 Int_t iACORDE = 0;
295
296 //=================== Alice BODY parameters =============================
297 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
298
299
300 if (iMAG)
301 {
302 //=================== MAG parameters ============================
303 // --- Start with Magnet since detector layouts may be depending ---
304 // --- on the selected Magnet dimensions ---
305 AliMAG *MAG = new AliMAG("MAG", "Magnet");
306 }
307
308
309 if (iABSO)
310 {
311 //=================== ABSO parameters ============================
96fb3353 312 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
27b7c1fb 313 }
314
315 if (iDIPO)
316 {
317 //=================== DIPO parameters ============================
318
96fb3353 319 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
27b7c1fb 320 }
321
322 if (iHALL)
323 {
324 //=================== HALL parameters ============================
325
96fb3353 326 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
27b7c1fb 327 }
328
329
330 if (iFRAME)
331 {
332 //=================== FRAME parameters ============================
333
334 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
bab433b3 335 FRAME->SetHoles(1);
27b7c1fb 336 }
337
338 if (iSHIL)
339 {
340 //=================== SHIL parameters ============================
341
96fb3353 342 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
27b7c1fb 343 }
344
345
346 if (iPIPE)
347 {
348 //=================== PIPE parameters ============================
349
96fb3353 350 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
27b7c1fb 351 }
352
e89fa83b 353 if (iITS)
354 {
355 //=================== ITS parameters ============================
27b7c1fb 356
99887042 357 AliITS *ITS = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
27b7c1fb 358 }
359
360 if (iTPC)
361 {
362 //============================ TPC parameters =====================
363 AliTPC *TPC = new AliTPCv2("TPC", "Default");
364 }
365
366
367 if (iTOF) {
368 //=================== TOF parameters ============================
0d7ef405 369 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
27b7c1fb 370 }
371
372
373 if (iHMPID)
374 {
375 //=================== HMPID parameters ===========================
99887042 376 AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
27b7c1fb 377
378 }
379
380
381 if (iZDC)
382 {
383 //=================== ZDC parameters ============================
384
cb14d625 385 AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
27b7c1fb 386 }
387
388 if (iTRD)
389 {
390 //=================== TRD parameters ============================
391
392 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
393 }
394
395 if (iFMD)
396 {
397 //=================== FMD parameters ============================
398 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
399 }
400
401 if (iMUON)
402 {
403 //=================== MUON parameters ===========================
404 // New MUONv1 version (geometry defined via builders)
405 AliMUON *MUON = new AliMUONv1("MUON", "default");
406 }
407 //=================== PHOS parameters ===========================
408
409 if (iPHOS)
410 {
411 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
412 }
413
414
415 if (iPMD)
416 {
417 //=================== PMD parameters ============================
418 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
419 }
420
421 if (iT0)
422 {
423 //=================== T0 parameters ============================
424 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
425 }
426
427 if (iEMCAL)
428 {
429 //=================== EMCAL parameters ============================
8224b11d 430 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
27b7c1fb 431 }
432
433 if (iACORDE)
434 {
435 //=================== ACORDE parameters ============================
f7882672 436 AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
27b7c1fb 437 }
438
439 if (iVZERO)
440 {
441 //=================== VZERO parameters ============================
442 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
443 }
444
445
446}
447
448Float_t EtaToTheta(Float_t arg){
449 return (180./TMath::Pi())*2.*atan(exp(-arg));
450}
451
452
453
454AliGenerator* GeneratorFactory(PprRun_t srun) {
455 Int_t isw = 3;
456 if (srad == kNoGluonRadiation) isw = 0;
457
458
459 AliGenerator * gGener = 0x0;
460 switch (srun) {
461 case test50:
462 {
463 comment = comment.Append(":HIJINGparam test 50 particles");
464 AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
465 gener->SetMomentumRange(0, 999999.);
466 gener->SetPhiRange(0., 360.);
467 // Set pseudorapidity range from -8 to 8.
468 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
469 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
470 gener->SetThetaRange(thmin,thmax);
471 gGener=gener;
472 }
473 break;
474 case kParam_8000:
475 {
476 comment = comment.Append(":HIJINGparam N=8000");
477 AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
478 gener->SetMomentumRange(0, 999999.);
479 gener->SetPhiRange(0., 360.);
480 // Set pseudorapidity range from -8 to 8.
481 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
482 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
483 gener->SetThetaRange(thmin,thmax);
484 gGener=gener;
485 }
486 break;
487 case kParam_4000:
488 {
489 comment = comment.Append("HIJINGparam N=4000");
490 AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
491 gener->SetMomentumRange(0, 999999.);
492 gener->SetPhiRange(0., 360.);
493 // Set pseudorapidity range from -8 to 8.
494 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
495 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
496 gener->SetThetaRange(thmin,thmax);
497 gGener=gener;
498 }
499 break;
500 case kParam_2000:
501 {
502 comment = comment.Append("HIJINGparam N=2000");
503 AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
504 gener->SetMomentumRange(0, 999999.);
505 gener->SetPhiRange(0., 360.);
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 gGener=gener;
511 }
512 break;
513//
514// Hijing Central
515//
516 case kHijing_cent1:
517 {
518 comment = comment.Append("HIJING cent1");
519 AliGenHijing *gener = HijingStandard();
520// impact parameter range
521 gener->SetImpactParameterRange(0., 5.);
522 gGener=gener;
523 }
524 break;
525 case kHijing_cent2:
526 {
527 comment = comment.Append("HIJING cent2");
528 AliGenHijing *gener = HijingStandard();
529// impact parameter range
530 gener->SetImpactParameterRange(0., 2.);
531 gGener=gener;
532 }
533 break;
534//
535// Hijing Peripheral
536//
537 case kHijing_per1:
538 {
539 comment = comment.Append("HIJING per1");
540 AliGenHijing *gener = HijingStandard();
541// impact parameter range
542 gener->SetImpactParameterRange(5., 8.6);
543 gGener=gener;
544 }
545 break;
546 case kHijing_per2:
547 {
548 comment = comment.Append("HIJING per2");
549 AliGenHijing *gener = HijingStandard();
550// impact parameter range
551 gener->SetImpactParameterRange(8.6, 11.2);
552 gGener=gener;
553 }
554 break;
555 case kHijing_per3:
556 {
557 comment = comment.Append("HIJING per3");
558 AliGenHijing *gener = HijingStandard();
559// impact parameter range
560 gener->SetImpactParameterRange(11.2, 13.2);
561 gGener=gener;
562 }
563 break;
564 case kHijing_per4:
565 {
566 comment = comment.Append("HIJING per4");
567 AliGenHijing *gener = HijingStandard();
568// impact parameter range
569 gener->SetImpactParameterRange(13.2, 15.);
570 gGener=gener;
571 }
572 break;
573 case kHijing_per5:
574 {
575 comment = comment.Append("HIJING per5");
576 AliGenHijing *gener = HijingStandard();
577// impact parameter range
578 gener->SetImpactParameterRange(15., 100.);
579 gGener=gener;
580 }
581 break;
582//
583// Jet-Jet
584//
585 case kHijing_jj25:
586 {
587 comment = comment.Append("HIJING Jet 25 GeV");
588 AliGenHijing *gener = HijingStandard();
589// impact parameter range
590 gener->SetImpactParameterRange(0., 5.);
591 // trigger
592 gener->SetTrigger(1);
593 gener->SetPtJet(25.);
594 gener->SetRadiation(isw);
595 gener->SetSimpleJets(!isw);
596 gener->SetJetEtaRange(-0.3,0.3);
597 gener->SetJetPhiRange(75., 165.);
598 gGener=gener;
599 }
600 break;
601
602 case kHijing_jj50:
603 {
604 comment = comment.Append("HIJING Jet 50 GeV");
605 AliGenHijing *gener = HijingStandard();
606// impact parameter range
607 gener->SetImpactParameterRange(0., 5.);
608 // trigger
609 gener->SetTrigger(1);
610 gener->SetPtJet(50.);
611 gener->SetRadiation(isw);
612 gener->SetSimpleJets(!isw);
613 gener->SetJetEtaRange(-0.3,0.3);
614 gener->SetJetPhiRange(75., 165.);
615 gGener=gener;
616 }
617 break;
618
619 case kHijing_jj75:
620 {
621 comment = comment.Append("HIJING Jet 75 GeV");
622 AliGenHijing *gener = HijingStandard();
623// impact parameter range
624 gener->SetImpactParameterRange(0., 5.);
625 // trigger
626 gener->SetTrigger(1);
627 gener->SetPtJet(75.);
628 gener->SetRadiation(isw);
629 gener->SetSimpleJets(!isw);
630 gener->SetJetEtaRange(-0.3,0.3);
631 gener->SetJetPhiRange(75., 165.);
632 gGener=gener;
633 }
634 break;
635
636 case kHijing_jj100:
637 {
638 comment = comment.Append("HIJING Jet 100 GeV");
639 AliGenHijing *gener = HijingStandard();
640// impact parameter range
641 gener->SetImpactParameterRange(0., 5.);
642 // trigger
643 gener->SetTrigger(1);
644 gener->SetPtJet(100.);
645 gener->SetRadiation(isw);
646 gener->SetSimpleJets(!isw);
647 gener->SetJetEtaRange(-0.3,0.3);
648 gener->SetJetPhiRange(75., 165.);
649 gGener=gener;
650 }
651 break;
652
653 case kHijing_jj200:
654 {
655 comment = comment.Append("HIJING Jet 200 GeV");
656 AliGenHijing *gener = HijingStandard();
657// impact parameter range
658 gener->SetImpactParameterRange(0., 5.);
659 // trigger
660 gener->SetTrigger(1);
661 gener->SetPtJet(200.);
662 gener->SetRadiation(isw);
663 gener->SetSimpleJets(!isw);
664 gener->SetJetEtaRange(-0.3,0.3);
665 gener->SetJetPhiRange(75., 165.);
666 gGener=gener;
667 }
668 break;
669//
670// Gamma-Jet
671//
672 case kHijing_gj25:
673 {
674 comment = comment.Append("HIJING Gamma 25 GeV");
675 AliGenHijing *gener = HijingStandard();
676// impact parameter range
677 gener->SetImpactParameterRange(0., 5.);
678 // trigger
679 gener->SetTrigger(2);
680 gener->SetPtJet(25.);
681 gener->SetRadiation(isw);
682 gener->SetSimpleJets(!isw);
683 gener->SetJetEtaRange(-0.12, 0.12);
684 gener->SetJetPhiRange(220., 320.);
685 gGener=gener;
686 }
687 break;
688
689 case kHijing_gj50:
690 {
691 comment = 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->SetRadiation(isw);
699 gener->SetSimpleJets(!isw);
700 gener->SetJetEtaRange(-0.12, 0.12);
701 gener->SetJetPhiRange(220., 320.);
702 gGener=gener;
703 }
704 break;
705
706 case kHijing_gj75:
707 {
708 comment = comment.Append("HIJING Gamma 75 GeV");
709 AliGenHijing *gener = HijingStandard();
710// impact parameter range
711 gener->SetImpactParameterRange(0., 5.);
712 // trigger
713 gener->SetTrigger(2);
714 gener->SetPtJet(75.);
715 gener->SetRadiation(isw);
716 gener->SetSimpleJets(!isw);
717 gener->SetJetEtaRange(-0.12, 0.12);
718 gener->SetJetPhiRange(220., 320.);
719 gGener=gener;
720 }
721 break;
722
723 case kHijing_gj100:
724 {
725 comment = comment.Append("HIJING Gamma 100 GeV");
726 AliGenHijing *gener = HijingStandard();
727// impact parameter range
728 gener->SetImpactParameterRange(0., 5.);
729 // trigger
730 gener->SetTrigger(2);
731 gener->SetPtJet(100.);
732 gener->SetRadiation(isw);
733 gener->SetSimpleJets(!isw);
734 gener->SetJetEtaRange(-0.12, 0.12);
735 gener->SetJetPhiRange(220., 320.);
736 gGener=gener;
737 }
738 break;
739
740 case kHijing_gj200:
741 {
742 comment = comment.Append("HIJING Gamma 200 GeV");
743 AliGenHijing *gener = HijingStandard();
744// impact parameter range
745 gener->SetImpactParameterRange(0., 5.);
746 // trigger
747 gener->SetTrigger(2);
748 gener->SetPtJet(200.);
749 gener->SetRadiation(isw);
750 gener->SetSimpleJets(!isw);
751 gener->SetJetEtaRange(-0.12, 0.12);
752 gener->SetJetPhiRange(220., 320.);
753 gGener=gener;
754 }
755 break;
756 case kHijing_pA:
757 {
758 comment = comment.Append("HIJING pA");
759
760 AliGenCocktail *gener = new AliGenCocktail();
761
762 AliGenHijing *hijing = new AliGenHijing(-1);
763// centre of mass energy
764 hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
765// impact parameter range
766 hijing->SetImpactParameterRange(0., 15.);
767// reference frame
768 hijing->SetReferenceFrame("CMS");
769 hijing->SetBoostLHC(1);
770// projectile
771 hijing->SetProjectile("P", 1, 1);
772 hijing->SetTarget ("A", 208, 82);
773// tell hijing to keep the full parent child chain
774 hijing->KeepFullEvent();
775// enable jet quenching
776 hijing->SetJetQuenching(0);
777// enable shadowing
778 hijing->SetShadowing(1);
779// Don't track spectators
780 hijing->SetSpectators(0);
781// kinematic selection
782 hijing->SetSelectAll(0);
783//
784 AliGenSlowNucleons* gray = new AliGenSlowNucleons(1);
785 AliSlowNucleonModel* model = new AliSlowNucleonModelExp();
786 gray->SetSlowNucleonModel(model);
787 gray->SetDebug(1);
788 gener->AddGenerator(hijing,"Hijing pPb", 1);
789 gener->AddGenerator(gray, "Gray Particles",1);
790 gGener=gener;
791 }
792 break;
793 case kPythia6:
794 {
795 comment = comment.Append(":Pythia p-p @ 14 TeV");
796 AliGenPythia *gener = new AliGenPythia(-1);
797 gener->SetMomentumRange(0,999999);
798 gener->SetThetaRange(0., 180.);
799 gener->SetYRange(-12,12);
800 gener->SetPtRange(0,1000);
801 gener->SetProcess(kPyMb);
802 gener->SetEnergyCMS(14000.);
803 gGener=gener;
804 }
805 break;
806 case kPythia6Jets20_24:
807 {
808 comment = comment.Append(":Pythia jets 20-24 GeV @ 5.5 TeV");
809 AliGenPythia * gener = new AliGenPythia(-1);
810 gener->SetEnergyCMS(5500.);// Centre of mass energy
811 gener->SetProcess(kPyJets);// Process type
812 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
813 gener->SetJetPhiRange(0., 360.);
814 gener->SetJetEtRange(10., 1000.);
815 gener->SetGluonRadiation(1,1);
816 // gener->SetPtKick(0.);
817 // Structure function
818 gener->SetStrucFunc(kCTEQ4L);
819 gener->SetPtHard(20., 24.);// Pt transfer of the hard scattering
820 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
821 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
822 gGener=gener;
823 }
824 break;
825 case kPythia6Jets24_29:
826 {
827 comment = comment.Append(":Pythia jets 24-29 GeV @ 5.5 TeV");
828 AliGenPythia * gener = new AliGenPythia(-1);
829 gener->SetEnergyCMS(5500.);// Centre of mass energy
830 gener->SetProcess(kPyJets);// Process type
831 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
832 gener->SetJetPhiRange(0., 360.);
833 gener->SetJetEtRange(10., 1000.);
834 gener->SetGluonRadiation(1,1);
835 // gener->SetPtKick(0.);
836 // Structure function
837 gener->SetStrucFunc(kCTEQ4L);
838 gener->SetPtHard(24., 29.);// Pt transfer of the hard scattering
839 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
840 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
841 gGener=gener;
842 }
843 break;
844 case kPythia6Jets29_35:
845 {
846 comment = comment.Append(":Pythia jets 29-35 GeV @ 5.5 TeV");
847 AliGenPythia * gener = new AliGenPythia(-1);
848 gener->SetEnergyCMS(5500.);// Centre of mass energy
849 gener->SetProcess(kPyJets);// Process type
850 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
851 gener->SetJetPhiRange(0., 360.);
852 gener->SetJetEtRange(10., 1000.);
853 gener->SetGluonRadiation(1,1);
854 // gener->SetPtKick(0.);
855 // Structure function
856 gener->SetStrucFunc(kCTEQ4L);
857 gener->SetPtHard(29., 35.);// Pt transfer of the hard scattering
858 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
859 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
860 gGener=gener;
861 }
862 break;
863 case kPythia6Jets35_42:
864 {
865 comment = comment.Append(":Pythia jets 35-42 GeV @ 5.5 TeV");
866 AliGenPythia * gener = new AliGenPythia(-1);
867 gener->SetEnergyCMS(5500.);// Centre of mass energy
868 gener->SetProcess(kPyJets);// Process type
869 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
870 gener->SetJetPhiRange(0., 360.);
871 gener->SetJetEtRange(10., 1000.);
872 gener->SetGluonRadiation(1,1);
873 // gener->SetPtKick(0.);
874 // Structure function
875 gener->SetStrucFunc(kCTEQ4L);
876 gener->SetPtHard(35., 42.);// Pt transfer of the hard scattering
877 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
878 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
879 gGener=gener;
880 }
881 break;
882 case kPythia6Jets42_50:
883 {
884 comment = comment.Append(":Pythia jets 42-50 GeV @ 5.5 TeV");
885 AliGenPythia * gener = new AliGenPythia(-1);
886 gener->SetEnergyCMS(5500.);// Centre of mass energy
887 gener->SetProcess(kPyJets);// Process type
888 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
889 gener->SetJetPhiRange(0., 360.);
890 gener->SetJetEtRange(10., 1000.);
891 gener->SetGluonRadiation(1,1);
892 // gener->SetPtKick(0.);
893 // Structure function
894 gener->SetStrucFunc(kCTEQ4L);
895 gener->SetPtHard(42., 50.);// Pt transfer of the hard scattering
896 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
897 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
898 gGener=gener;
899 }
900 break;
901 case kPythia6Jets50_60:
902 {
903 comment = comment.Append(":Pythia jets 50-60 GeV @ 5.5 TeV");
904 AliGenPythia * gener = new AliGenPythia(-1);
905 gener->SetEnergyCMS(5500.);// Centre of mass energy
906 gener->SetProcess(kPyJets);// Process type
907 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
908 gener->SetJetPhiRange(0., 360.);
909 gener->SetJetEtRange(10., 1000.);
910 gener->SetGluonRadiation(1,1);
911 // gener->SetPtKick(0.);
912 // Structure function
913 gener->SetStrucFunc(kCTEQ4L);
914 gener->SetPtHard(50., 60.);// Pt transfer of the hard scattering
915 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
916 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
917 gGener=gener;
918 }
919 break;
920 case kPythia6Jets60_72:
921 {
922 comment = comment.Append(":Pythia jets 60-72 GeV @ 5.5 TeV");
923 AliGenPythia * gener = new AliGenPythia(-1);
924 gener->SetEnergyCMS(5500.);// Centre of mass energy
925 gener->SetProcess(kPyJets);// Process type
926 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
927 gener->SetJetPhiRange(0., 360.);
928 gener->SetJetEtRange(10., 1000.);
929 gener->SetGluonRadiation(1,1);
930 // gener->SetPtKick(0.);
931 // Structure function
932 gener->SetStrucFunc(kCTEQ4L);
933 gener->SetPtHard(60., 72.);// Pt transfer of the hard scattering
934 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
935 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
936 gGener=gener;
937 }
938 break;
939 case kPythia6Jets72_86:
940 {
941 comment = comment.Append(":Pythia jets 72-86 GeV @ 5.5 TeV");
942 AliGenPythia * gener = new AliGenPythia(-1);
943 gener->SetEnergyCMS(5500.);// Centre of mass energy
944 gener->SetProcess(kPyJets);// Process type
945 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
946 gener->SetJetPhiRange(0., 360.);
947 gener->SetJetEtRange(10., 1000.);
948 gener->SetGluonRadiation(1,1);
949 // gener->SetPtKick(0.);
950 // Structure function
951 gener->SetStrucFunc(kCTEQ4L);
952 gener->SetPtHard(72., 86.);// Pt transfer of the hard scattering
953 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
954 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
955 gGener=gener;
956 }
957 break;
958 case kPythia6Jets86_104:
959 {
960 comment = comment.Append(":Pythia jets 86-104 GeV @ 5.5 TeV");
961 AliGenPythia * gener = new AliGenPythia(-1);
962 gener->SetEnergyCMS(5500.);// Centre of mass energy
963 gener->SetProcess(kPyJets);// Process type
964 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
965 gener->SetJetPhiRange(0., 360.);
966 gener->SetJetEtRange(10., 1000.);
967 gener->SetGluonRadiation(1,1);
968 // gener->SetPtKick(0.);
969 // Structure function
970 gener->SetStrucFunc(kCTEQ4L);
971 gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
972 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
973 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
974 gGener=gener;
975 }
976 break;
977 case kPythia6Jets104_125:
978 {
979 comment = comment.Append(":Pythia jets 105-125 GeV @ 5.5 TeV");
980 AliGenPythia * gener = new AliGenPythia(-1);
981 gener->SetEnergyCMS(5500.);// Centre of mass energy
982 gener->SetProcess(kPyJets);// Process type
983 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
984 gener->SetJetPhiRange(0., 360.);
985 gener->SetJetEtRange(10., 1000.);
986 gener->SetGluonRadiation(1,1);
987 // gener->SetPtKick(0.);
988 // Structure function
989 gener->SetStrucFunc(kCTEQ4L);
990 gener->SetPtHard(104., 125.);// Pt transfer of the hard scattering
991 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
992 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
993 gGener=gener;
994 }
995 break;
996 case kPythia6Jets125_150:
997 {
998 comment = comment.Append(":Pythia jets 125-150 GeV @ 5.5 TeV");
999 AliGenPythia * gener = new AliGenPythia(-1);
1000 gener->SetEnergyCMS(5500.);// Centre of mass energy
1001 gener->SetProcess(kPyJets);// Process type
1002 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1003 gener->SetJetPhiRange(0., 360.);
1004 gener->SetJetEtRange(10., 1000.);
1005 gener->SetGluonRadiation(1,1);
1006 // gener->SetPtKick(0.);
1007 // Structure function
1008 gener->SetStrucFunc(kCTEQ4L);
1009 gener->SetPtHard(125., 150.);// Pt transfer of the hard scattering
1010 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1011 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1012 gGener=gener;
1013 }
1014 break;
1015 case kPythia6Jets150_180:
1016 {
1017 comment = comment.Append(":Pythia jets 150-180 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(150., 180.);// 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 kD0PbPb5500:
1035 {
1036 comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
1037 AliGenPythia * gener = new AliGenPythia(10);
1038 gener->SetProcess(kPyD0PbPbMNR);
1039 gener->SetStrucFunc(kCTEQ4L);
1040 gener->SetPtHard(2.1,-1.0);
1041 gener->SetEnergyCMS(5500.);
1042 gener->SetNuclei(208,208);
1043 gener->SetForceDecay(kHadronicD);
1044 gener->SetYRange(-2,2);
1045 gener->SetFeedDownHigherFamily(kFALSE);
1046 gener->SetStackFillOpt(AliGenPythia::kParentSelection);
1047 gener->SetCountMode(AliGenPythia::kCountParents);
1048 gGener=gener;
1049 }
1050 break;
1051 case kCharmSemiElPbPb5500:
1052 {
1053 comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
1054 AliGenPythia * gener = new AliGenPythia(10);
1055 gener->SetProcess(kPyCharmPbPbMNR);
1056 gener->SetStrucFunc(kCTEQ4L);
1057 gener->SetPtHard(2.1,-1.0);
1058 gener->SetEnergyCMS(5500.);
1059 gener->SetNuclei(208,208);
1060 gener->SetForceDecay(kSemiElectronic);
1061 gener->SetYRange(-2,2);
1062 gener->SetFeedDownHigherFamily(kFALSE);
1063 gener->SetCountMode(AliGenPythia::kCountParents);
1064 gGener=gener;
1065 }
1066 break;
1067 case kBeautySemiElPbPb5500:
1068 {
1069 comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
1070 AliGenPythia *gener = new AliGenPythia(10);
1071 gener->SetProcess(kPyBeautyPbPbMNR);
1072 gener->SetStrucFunc(kCTEQ4L);
1073 gener->SetPtHard(2.75,-1.0);
1074 gener->SetEnergyCMS(5500.);
1075 gener->SetNuclei(208,208);
1076 gener->SetForceDecay(kSemiElectronic);
1077 gener->SetYRange(-2,2);
1078 gener->SetFeedDownHigherFamily(kFALSE);
1079 gener->SetCountMode(AliGenPythia::kCountParents);
1080 gGener=gener;
1081 }
1082 break;
1083 case kCocktailTRD:
1084 {
1085 comment = comment.Append(" Cocktail for TRD at 5.5 TeV");
1086 AliGenCocktail *gener = new AliGenCocktail();
1087
1088 AliGenParam *phi = new AliGenParam(10,
1089 new AliGenMUONlib(),
1090 AliGenMUONlib::kPhi,
1091 "Vogt PbPb");
1092
1093 phi->SetPtRange(0, 100);
1094 phi->SetYRange(-1., +1.);
1095 phi->SetForceDecay(kDiElectron);
1096
1097 AliGenParam *omega = new AliGenParam(10,
1098 new AliGenMUONlib(),
1099 AliGenMUONlib::kOmega,
1100 "Vogt PbPb");
1101
1102 omega->SetPtRange(0, 100);
1103 omega->SetYRange(-1., +1.);
1104 omega->SetForceDecay(kDiElectron);
1105
1106 AliGenParam *jpsi = new AliGenParam(10,
1107 new AliGenMUONlib(),
1108 AliGenMUONlib::kJpsiFamily,
1109 "Vogt PbPb");
1110
1111 jpsi->SetPtRange(0, 100);
1112 jpsi->SetYRange(-1., +1.);
1113 jpsi->SetForceDecay(kDiElectron);
1114
1115 AliGenParam *ups = new AliGenParam(10,
1116 new AliGenMUONlib(),
1117 AliGenMUONlib::kUpsilonFamily,
1118 "Vogt PbPb");
1119 ups->SetPtRange(0, 100);
1120 ups->SetYRange(-1., +1.);
1121 ups->SetForceDecay(kDiElectron);
1122
1123 AliGenParam *charm = new AliGenParam(10,
1124 new AliGenMUONlib(),
1125 AliGenMUONlib::kCharm,
1126 "central");
1127 charm->SetPtRange(0, 100);
1128 charm->SetYRange(-1.5, +1.5);
1129 charm->SetForceDecay(kSemiElectronic);
1130
1131
1132 AliGenParam *beauty = new AliGenParam(10,
1133 new AliGenMUONlib(),
1134 AliGenMUONlib::kBeauty,
1135 "central");
1136 beauty->SetPtRange(0, 100);
1137 beauty->SetYRange(-1.5, +1.5);
1138 beauty->SetForceDecay(kSemiElectronic);
1139
1140 AliGenParam *beautyJ = new AliGenParam(10,
1141 new AliGenMUONlib(),
1142 AliGenMUONlib::kBeauty,
1143 "central");
1144 beautyJ->SetPtRange(0, 100);
1145 beautyJ->SetYRange(-1.5, +1.5);
1146 beautyJ->SetForceDecay(kBJpsiDiElectron);
1147
1148 gener->AddGenerator(phi,"Phi",1);
1149 gener->AddGenerator(omega,"Omega",1);
1150 gener->AddGenerator(jpsi,"J/psi",1);
1151 gener->AddGenerator(ups,"Upsilon",1);
1152 gener->AddGenerator(charm,"Charm",1);
1153 gener->AddGenerator(beauty,"Beauty",1);
1154 gener->AddGenerator(beautyJ,"J/Psi from Beauty",1);
1155 gGener=gener;
1156 }
1157 break;
1158 case kPyJJ:
1159 {
1160 comment = comment.Append(" Jet-jet at 5.5 TeV");
1161 AliGenPythia *gener = new AliGenPythia(-1);
1162 gener->SetEnergyCMS(5500.);
1163 gener->SetProcess(kPyJets);
1164 Double_t ptHardMin=10.0, ptHardMax=-1.0;
1165 gener->SetPtHard(ptHardMin,ptHardMax);
1166 gener->SetYHard(-0.7,0.7);
1167 gener->SetJetEtaRange(-0.2,0.2);
1168 gener->SetEventListRange(0,1);
1169 gGener=gener;
1170 }
1171 break;
1172 case kPyGJ:
1173 {
1174 comment = comment.Append(" Gamma-jet at 5.5 TeV");
1175 AliGenPythia *gener = new AliGenPythia(-1);
1176 gener->SetEnergyCMS(5500.);
1177 gener->SetProcess(kPyDirectGamma);
1178 Double_t ptHardMin=10.0, ptHardMax=-1.0;
1179 gener->SetPtHard(ptHardMin,ptHardMax);
1180 gener->SetYHard(-1.0,1.0);
1181 gener->SetGammaEtaRange(-0.13,0.13);
1182 gener->SetGammaPhiRange(210.,330.);
1183 gener->SetEventListRange(0,1);
1184 gGener=gener;
1185 }
1186 break;
1187 case kMuonCocktailCent1:
1188 {
1189 comment = comment.Append(" Muon Cocktail Cent1");
1190 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1191 gener->SetPtRange(0.4,100.); // Transverse momentum range
1192 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1193 gener->SetYRange(-4.0,-2.4);
1194 gener->SetMuonPtCut(0.8);
1195 gener->SetMuonThetaCut(171.,178.);
1196 gener->SetMuonMultiplicity(2);
1197 gener->SetImpactParameterRange(0.,5.); //Centrality class Cent1 for PDC04
1198 gGener=gener;
1199 }
1200 break;
1201 case kMuonCocktailPer1:
1202 {
1203 comment = comment.Append(" Muon Cocktail Per1");
1204 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1205 gener->SetPtRange(0.0,100.); // Transverse momentum range
1206 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1207 gener->SetYRange(-4.0,-2.4);
1208 gener->SetMuonPtCut(0.8);
1209 gener->SetMuonThetaCut(171.,178.);
1210 gener->SetMuonMultiplicity(2);
1211 gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1212 gGener=gener;
1213 }
1214 break;
1215 case kMuonCocktailPer4:
1216 {
1217 comment = comment.Append(" Muon Cocktail Per4");
1218 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1219 gener->SetPtRange(0.0,100.); // Transverse momentum range
1220 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1221 gener->SetYRange(-4.0,-2.4);
1222 gener->SetMuonPtCut(0.8);
1223 gener->SetMuonThetaCut(171.,178.);
1224 gener->SetMuonMultiplicity(2);
1225 gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1226 gGener=gener;
1227 }
1228 break;
1229 case kMuonCocktailCent1HighPt:
1230 {
1231 comment = comment.Append(" Muon Cocktail HighPt Cent1");
1232 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1233 gener->SetPtRange(0.0,100.); // Transverse momentum range
1234 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1235 gener->SetYRange(-4.0,-2.4);
1236 gener->SetMuonPtCut(2.5);
1237 gener->SetMuonThetaCut(171.,178.);
1238 gener->SetMuonMultiplicity(2);
1239 gener->SetImpactParameterRange(0.,5.); //Centrality class Cent1 for PDC04
1240 gGener=gener;
1241 }
1242 break;
1243 case kMuonCocktailPer1HighPt :
1244 {
1245 comment = comment.Append(" Muon Cocktail HighPt Per1");
1246 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1247 gener->SetPtRange(0.0,100.); // Transverse momentum range
1248 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1249 gener->SetYRange(-4.0,-2.4);
1250 gener->SetMuonPtCut(2.5);
1251 gener->SetMuonThetaCut(171.,178.);
1252 gener->SetMuonMultiplicity(2);
1253 gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1254 gGener=gener;
1255 }
1256 break;
1257 case kMuonCocktailPer4HighPt:
1258 {
1259 comment = comment.Append(" Muon Cocktail HighPt Per4");
1260 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1261 gener->SetPtRange(0.0,100.); // Transverse momentum range
1262 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1263 gener->SetYRange(-4.0,-2.4);
1264 gener->SetMuonPtCut(2.5);
1265 gener->SetMuonThetaCut(171.,178.);
1266 gener->SetMuonMultiplicity(2);
1267 gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1268 gGener=gener;
1269 }
1270 break;
1271 case kMuonCocktailCent1Single:
1272 {
1273 comment = comment.Append(" Muon Cocktail Single Cent1");
1274 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1275 gener->SetPtRange(0.0,100.); // Transverse momentum range
1276 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1277 gener->SetYRange(-4.0,-2.4);
1278 gener->SetMuonPtCut(0.8);
1279 gener->SetMuonThetaCut(171.,178.);
1280 gener->SetMuonMultiplicity(1);
1281 gener->SetImpactParameterRange(0.,5.); //Centrality class Cent1 for PDC04
1282 gGener=gener;
1283 }
1284 break;
1285 case kMuonCocktailPer1Single :
1286 {
1287 comment = comment.Append(" Muon Cocktail Single Per1");
1288 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1289 gener->SetPtRange(0.0,100.); // Transverse momentum range
1290 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1291 gener->SetYRange(-4.0,-2.4);
1292 gener->SetMuonPtCut(0.8);
1293 gener->SetMuonThetaCut(171.,178.);
1294 gener->SetMuonMultiplicity(1);
1295 gener->SetImpactParameterRange(5.,8.6);//Centrality class Per1 for PDC04
1296 gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1297 gGener=gener;
1298 }
1299 break;
1300 case kMuonCocktailPer4Single:
1301 {
1302 comment = comment.Append(" Muon Cocktail Single 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(0.8);
1308 gener->SetMuonThetaCut(171.,178.);
1309 gener->SetMuonMultiplicity(1);
1310 gener->SetImpactParameterRange(13.2,15.0);//Centrality class Per4 for PDC04
1311 gGener=gener;
1312 }
1313 break;
1314 case kFlow_2_2000:
1315 {
1316 comment = comment.Append(" Flow with dN/deta = 2000, vn = 2%");
1317 gGener = GeVSimStandard(2000., 2.);
1318 }
1319 break;
1320
1321 case kFlow_10_2000:
1322 {
1323 comment = comment.Append(" Flow with dN/deta = 2000, vn = 10%");
1324 gGener = GeVSimStandard(2000., 10.);
1325 }
1326 break;
1327
1328 case kFlow_6_2000:
1329 {
1330 comment = comment.Append(" Flow with dN/deta = 2000, vn = 6%");
1331 gGener = GeVSimStandard(2000., 6.);
1332 }
1333 break;
1334
1335 case kFlow_6_5000:
1336 {
1337 comment = comment.Append(" Flow with dN/deta = 5000, vn = 6%");
1338 gGener = GeVSimStandard(5000., 6.);
1339 }
1340 break;
1341 case kHIJINGplus:
1342 {
1343 //
1344 // The cocktail
1345 AliGenCocktail *gener = new AliGenCocktail();
1346
1347 //
1348 // Charm production by Pythia
1349 AliGenPythia * genpyc = new AliGenPythia(230);
1350 genpyc->SetProcess(kPyCharmPbPbMNR);
1351 genpyc->SetStrucFunc(kCTEQ4L);
1352 genpyc->SetPtHard(2.1,-1.0);
1353 genpyc->SetEnergyCMS(5500.);
1354 genpyc->SetNuclei(208,208);
1355 genpyc->SetYRange(-999,999);
1356 genpyc->SetForceDecay(kAll);
1357 genpyc->SetFeedDownHigherFamily(kFALSE);
1358 genpyc->SetCountMode(AliGenPythia::kCountParents);
1359 //
1360 // Beauty production by Pythia
1361 AliGenPythia * genpyb = new AliGenPythia(9);
1362 genpyb->SetProcess(kPyBeautyPbPbMNR);
1363 genpyb->SetStrucFunc(kCTEQ4L);
1364 genpyb->SetPtHard(2.75,-1.0);
1365 genpyb->SetEnergyCMS(5500.);
1366 genpyb->SetNuclei(208,208);
1367 genpyb->SetYRange(-999,999);
1368 genpyb->SetForceDecay(kAll);
1369 genpyb->SetFeedDownHigherFamily(kFALSE);
1370 genpyb->SetCountMode(AliGenPythia::kCountParents);
1371 //
1372 // Hyperons
1373 //
1374 AliGenSTRANGElib *lib = new AliGenSTRANGElib();
1375 Int_t particle;
1376 // Xi
7d245e03 1377 particle = kXiMinus;
27b7c1fb 1378 AliGenParam *genXi = new AliGenParam(16,particle,lib->GetPt(particle),lib->GetY(particle),lib->GetIp(particle));
1379 genXi->SetPtRange(0., 12.);
1380 genXi->SetYRange(-1.1, 1.1);
1381 genXi->SetForceDecay(kNoDecay);
1382
1383 //
1384 // Omega
7d245e03 1385 particle = kOmegaMinus;
27b7c1fb 1386 AliGenParam *genOmega = new AliGenParam(10,particle,lib->GetPt(particle),lib->GetY(particle),lib->GetIp(particle));
1387 genOmega->SetPtRange(0, 12.);
1388 genOmega->SetYRange(-1.1, 1.1);
1389 genOmega->SetForceDecay(kNoDecay);
1390
1391 //
1392 // Central Hijing
1393 AliGenHijing *genHi = HijingStandard();
1394 genHi->SwitchOffHeavyQuarks(kTRUE);
1395 genHi->SetImpactParameterRange(0.,5.);
1396 //
1397 // Add everything to the cocktail and shake ...
1398 gener->AddGenerator(genHi, "Hijing cent1", 1);
1399 gener->AddGenerator(genpyc, "Extra charm", 1);
1400 gener->AddGenerator(genpyb, "Extra beauty", 1);
1401 gener->AddGenerator(genXi, "Xi" , 1);
1402 gener->AddGenerator(genOmega, "Omega", 1);
1403 gGener = gener;
1404 }
1405 break;
1406 default: break;
1407 }
1408
1409 return gGener;
1410}
1411
1412AliGenHijing* HijingStandard()
1413{
1414 AliGenHijing *gener = new AliGenHijing(-1);
1415// centre of mass energy
1416 gener->SetEnergyCMS(5500.);
1417// reference frame
1418 gener->SetReferenceFrame("CMS");
1419// projectile
1420 gener->SetProjectile("A", 208, 82);
1421 gener->SetTarget ("A", 208, 82);
1422// tell hijing to keep the full parent child chain
1423 gener->KeepFullEvent();
1424// enable jet quenching
1425 gener->SetJetQuenching(1);
1426// enable shadowing
1427 gener->SetShadowing(1);
1428// neutral pion and heavy particle decays switched off
1429 gener->SetDecaysOff(1);
1430// Don't track spectators
1431 gener->SetSpectators(0);
1432// kinematic selection
1433 gener->SetSelectAll(0);
1434 return gener;
1435}
1436
1437AliGenGeVSim* GeVSimStandard(Float_t mult, Float_t vn)
1438{
1439 AliGenGeVSim* gener = new AliGenGeVSim(0);
1440//
1441// Mult is the number of charged particles in |eta| < 0.5
1442// Vn is in (%)
1443//
1444// Sigma of the Gaussian dN/deta
1445 Float_t sigma_eta = 2.75;
1446//
1447// Maximum eta
1448 Float_t etamax = 7.00;
1449//
1450//
1451// Scale from multiplicity in |eta| < 0.5 to |eta| < |etamax|
1452 Float_t mm = mult * (TMath::Erf(etamax/sigma_eta/sqrt(2.)) / TMath::Erf(0.5/sigma_eta/sqrt(2.)));
1453//
1454// Scale from charged to total multiplicity
1455//
1456 mm *= 1.587;
1457//
1458// Vn
1459 vn /= 100.;
1460//
1461// Define particles
1462//
1463//
1464// 78% Pions (26% pi+, 26% pi-, 26% p0) T = 250 MeV
1465 AliGeVSimParticle *pp = new AliGeVSimParticle(kPiPlus, 1, 0.26 * mm, 0.25, sigma_eta) ;
1466 AliGeVSimParticle *pm = new AliGeVSimParticle(kPiMinus, 1, 0.26 * mm, 0.25, sigma_eta) ;
1467 AliGeVSimParticle *p0 = new AliGeVSimParticle(kPi0, 1, 0.26 * mm, 0.25, sigma_eta) ;
1468//
1469// 12% Kaons (3% K0short, 3% K0long, 3% K+, 3% K-) T = 300 MeV
1470 AliGeVSimParticle *ks = new AliGeVSimParticle(kK0Short, 1, 0.03 * mm, 0.30, sigma_eta) ;
1471 AliGeVSimParticle *kl = new AliGeVSimParticle(kK0Long, 1, 0.03 * mm, 0.30, sigma_eta) ;
1472 AliGeVSimParticle *kp = new AliGeVSimParticle(kKPlus, 1, 0.03 * mm, 0.30, sigma_eta) ;
1473 AliGeVSimParticle *km = new AliGeVSimParticle(kKMinus, 1, 0.03 * mm, 0.30, sigma_eta) ;
1474//
1475// 10% Protons / Neutrons (5% Protons, 5% Neutrons) T = 250 MeV
1476 AliGeVSimParticle *pr = new AliGeVSimParticle(kProton, 1, 0.05 * mm, 0.25, sigma_eta) ;
1477 AliGeVSimParticle *ne = new AliGeVSimParticle(kNeutron, 1, 0.05 * mm, 0.25, sigma_eta) ;
1478//
1479// Set Elliptic Flow properties
1480
1481 Float_t pTsaturation = 2. ;
1482
1483 pp->SetEllipticParam(vn,pTsaturation,0.) ;
1484 pm->SetEllipticParam(vn,pTsaturation,0.) ;
1485 p0->SetEllipticParam(vn,pTsaturation,0.) ;
1486 pr->SetEllipticParam(vn,pTsaturation,0.) ;
1487 ne->SetEllipticParam(vn,pTsaturation,0.) ;
1488 ks->SetEllipticParam(vn,pTsaturation,0.) ;
1489 kl->SetEllipticParam(vn,pTsaturation,0.) ;
1490 kp->SetEllipticParam(vn,pTsaturation,0.) ;
1491 km->SetEllipticParam(vn,pTsaturation,0.) ;
1492//
1493// Set Direct Flow properties
1494 pp->SetDirectedParam(vn,1.0,0.) ;
1495 pm->SetDirectedParam(vn,1.0,0.) ;
1496 p0->SetDirectedParam(vn,1.0,0.) ;
1497 pr->SetDirectedParam(vn,1.0,0.) ;
1498 ne->SetDirectedParam(vn,1.0,0.) ;
1499 ks->SetDirectedParam(vn,1.0,0.) ;
1500 kl->SetDirectedParam(vn,1.0,0.) ;
1501 kp->SetDirectedParam(vn,1.0,0.) ;
1502 km->SetDirectedParam(vn,1.0,0.) ;
1503//
1504// Add particles to the list
1505 gener->AddParticleType(pp) ;
1506 gener->AddParticleType(pm) ;
1507 gener->AddParticleType(p0) ;
1508 gener->AddParticleType(pr) ;
1509 gener->AddParticleType(ne) ;
1510 gener->AddParticleType(ks) ;
1511 gener->AddParticleType(kl) ;
1512 gener->AddParticleType(kp) ;
1513 gener->AddParticleType(km) ;
1514//
1515// Random Ev.Plane ----------------------------------
1516 TF1 *rpa = new TF1("gevsimPsiRndm","1", 0, 360);
1517// --------------------------------------------------
1518 gener->SetPtRange(0., 9.) ; // Use a resonable range! (used for bin size in numerical integration)
1519 gener->SetPhiRange(0, 360);
1520 //
1521 // Set pseudorapidity range
1522 Float_t thmin = EtaToTheta(+etamax);
1523 Float_t thmax = EtaToTheta(-etamax);
1524 gener->SetThetaRange(thmin,thmax);
1525 return gener;
1526}
1527
1528
1529
1530void ProcessEnvironmentVars()
1531{
1532 // Run type
1533 if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
1534 for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
1535 if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
1536 srun = (PprRun_t)iRun;
1537 cout<<"Run type set to "<<pprRunName[iRun]<<endl;
1538 }
1539 }
1540 }
1541
1542 // Random Number seed
1543 if (gSystem->Getenv("CONFIG_SEED")) {
1544 sseed = atoi(gSystem->Getenv("CONFIG_SEED"));
1545 }
1546}