1 //____________________________________________________________________
5 // One can use the configuration macro in compiled mode by
6 // root [0] gSystem->Load("libgeant321");
7 // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
8 // -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
9 // root [0] .x grun.C(1,"ConfigPPR.C++")
12 @author Christian Holm Christensen <cholm@nbi.dk>
13 @date Mon Mar 27 12:50:29 2006
14 @brief Simulation configuration script
16 #if !defined(__CINT__) || defined(__MAKECINT__)
17 #include <Riostream.h>
20 #include <TVirtualMC.h>
22 #include "STEER/AliRunLoader.h"
23 #include "STEER/AliRun.h"
24 #include "STEER/AliConfig.h"
25 #include "STEER/AliGenerator.h"
26 #include "PYTHIA6/AliDecayerPythia.h"
27 #include "EVGEN/AliGenHIJINGpara.h"
28 #include "THijing/AliGenHijing.h"
29 #include "EVGEN/AliGenCocktail.h"
30 #include "EVGEN/AliGenSlowNucleons.h"
31 #include "EVGEN/AliSlowNucleonModelExp.h"
32 #include "EVGEN/AliGenParam.h"
33 #include "EVGEN/AliGenMUONlib.h"
34 #include "EVGEN/AliGenMUONCocktail.h"
35 #include "PYTHIA6/AliGenPythia.h"
36 #include "STEER/AliMagFMaps.h"
37 #include "STRUCT/AliBODY.h"
38 #include "STRUCT/AliMAG.h"
39 #include "STRUCT/AliABSOv0.h"
40 #include "STRUCT/AliDIPOv2.h"
41 #include "STRUCT/AliHALL.h"
42 #include "STRUCT/AliFRAMEv2.h"
43 #include "STRUCT/AliSHILv2.h"
44 #include "STRUCT/AliPIPEv0.h"
45 #include "ITS/AliITSvPPRasymmFMD.h"
46 #include "TPC/AliTPCv2.h"
47 #include "TOF/AliTOFv4T0.h"
48 #include "RICH/AliRICHv1.h"
49 #include "ZDC/AliZDCv2.h"
50 #include "TRD/AliTRDv1.h"
51 #include "FMD/AliFMDv1.h"
52 #include "MUON/AliMUONv1.h"
53 #include "MUON/AliMUONSt1GeometryBuilderV2.h"
54 #include "MUON/AliMUONSt2GeometryBuilder.h"
55 #include "MUON/AliMUONSlatGeometryBuilder.h"
56 #include "MUON/AliMUONTriggerGeometryBuilder.h"
57 #include "PHOS/AliPHOSv1.h"
58 #include "PMD/AliPMDv1.h"
59 #include "START/AliSTARTv1.h"
60 #include "EMCAL/AliEMCALv1.h"
61 #include "CRT/AliCRTv0.h"
62 #include "VZERO/AliVZEROv2.h"
65 //____________________________________________________________________
100 kPythia6Jets60_72, //
101 kPythia6Jets72_86, //
102 kPythia6Jets86_104, //
103 kPythia6Jets104_125, //
104 kPythia6Jets125_150, //
105 kPythia6Jets150_180, //
107 kCharmSemiElPbPb5500, //
108 kBeautySemiElPbPb5500, //
112 kMuonCocktailCent1, //
113 kMuonCocktailPer1, //
114 kMuonCocktailPer4, //
115 kMuonCocktailCent1HighPt, //
116 kMuonCocktailPer1HighPt, //
117 kMuonCocktailPer4HighPt, //
118 kMuonCocktailCent1Single, //
119 kMuonCocktailPer1Single, //
120 kMuonCocktailPer4Single,
128 //____________________________________________________________________
130 // Generator types names
132 const char* egName[kEgMax] = {
157 "kPythia6Jets20_24", //
158 "kPythia6Jets24_29", //
159 "kPythia6Jets29_35", //
160 "kPythia6Jets35_42", //
161 "kPythia6Jets42_50", //
162 "kPythia6Jets50_60", //
163 "kPythia6Jets60_72", //
164 "kPythia6Jets72_86", //
165 "kPythia6Jets86_104", //
166 "kPythia6Jets104_125", //
167 "kPythia6Jets125_150", //
168 "kPythia6Jets150_180", //
170 "kCharmSemiElPbPb5500", //
171 "kBeautySemiElPbPb5500", //
175 "kMuonCocktailCent1", //
176 "kMuonCocktailPer1", //
177 "kMuonCocktailPer4", //
178 "kMuonCocktailCent1HighPt", //
179 "kMuonCocktailPer1HighPt", //
180 "kMuonCocktailPer4HighPt", //
181 "kMuonCocktailCent1Single", //
182 "kMuonCocktailPer1Single", //
183 "kMuonCocktailPer4Single",
190 //____________________________________________________________________
196 //____________________________________________________________________
202 //____________________________________________________________________
209 //____________________________________________________________________
217 //____________________________________________________________________
219 Float_t EtaToTheta(Float_t eta);
220 Eg_t LookupEG(const Char_t* name);
221 AliGenerator* GeneratorFactory(EG_t eg, Rad_t rad, TString& comment);
222 AliGenHijing* HijingStandard();
223 void ProcessEnvironmentVars(EG_t& eg, Int_t& seed);
225 //____________________________________________________________________
229 //____________________________________________________________________
230 // This part for configuration
231 //static EG_t eg = test50;
232 //EG_t eg = kParam_fmd;
233 EG_t eg = kParam_2000; // kPythia;
234 Geo_t geo = kNoHoles;
235 Rad_t rad = kGluonRadiation;
237 Int_t seed = 12345; //Set 0 to use the current time
238 MC_t mc = kGEANT3TGEO;
240 //____________________________________________________________________
242 static TString comment;
244 //____________________________________________________________________
245 // Get settings from environment variables
246 ProcessEnvironmentVars(eg, seed);
248 //____________________________________________________________________
249 // Set Random Number seed
250 gRandom->SetSeed(seed);
251 cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl;
254 //__________________________________________________________________
258 // libraries required by fluka21
260 gSystem->Load("libGeom");
261 cout << "\t* Loading TFluka..." << endl;
262 gSystem->Load("libTFluka");
263 gSystem->MakeDirectory("peg");
267 cout << "\t* Instantiating TFluka..." << endl;
268 new TFluka("C++ Interface to Fluka", 0/*verbosity*/);
273 // Libraries needed by GEANT 3.21
275 gSystem->Load("libgeant321");
280 TGeant3* gmc = new TGeant3("C++ Interface to Geant3");
281 gmc->SetSWIT(4, 1000);
287 // Libraries needed by GEANT 3.21
289 gSystem->Load("libgeant321");
294 TGeant3TGeo* gmc = new TGeant3TGeo("C++ Interface to Geant3");
295 gmc->SetSWIT(4, 1000);
296 Printf("Making a TGeant3TGeo objet");
300 gAlice->Fatal("Config.C", "No MC type chosen");
304 //__________________________________________________________________
305 AliRunLoader* rl = 0;
307 cout<<"Config.C: Creating Run Loader ..."<<endl;
308 rl = AliRunLoader::Open("galice.root",
309 AliConfig::GetDefaultEventFolderName(),
312 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
315 rl->SetCompressionLevel(2);
316 rl->SetNumberOfEventsPerFile(3);
317 gAlice->SetRunLoader(rl);
319 //__________________________________________________________________
325 // Use kTRUE as argument to generate alice.pemf first
327 TString alice_pemf(gSystem->Which(".", "peg/mat17.pemf"));
328 if (!alice_pemf.IsNull())
329 ((TFluka*)gMC)->SetGeneratePemf(kFALSE);
331 ((TFluka*)gMC)->SetGeneratePemf(kTRUE);
332 TString flupro(gSystem->Getenv("FLUPRO"));
334 Fatal("Config.C", "Environment variable FLUPRO not set");
336 char* files[] = { "brems_fin.bin", "cohff.bin", "elasct.bin",
337 "gxsect.bin", "nuclear.bin", "sigmapi.bin",
339 char* file = files[0];
341 TString which(gSystem->Which(".", file));
342 if (which.IsNull()) {
343 if (gSystem->Symlink(Form("%s/%s", flupro.Data(), file), file)!=0)
344 Fatal("Config.C", "Couldn't link $(FLUPRO)/%s -> .", file);
349 TString neuxsc(gSystem->Which(".", "neuxsc.bin"));
351 gSystem->Symlink(Form("%s/neuxsc_72.bin", flupro.Data()),
353 gSystem->CopyFile("$(FLUPRO)/random.dat", "old.seed", kTRUE);
358 //__________________________________________________________________
360 // Set External decayer
362 AliDecayer *decayer = new AliDecayerPythia();
364 case kD0PbPb5500: decayer->SetForceDecay(kHadronicD); break;
365 case kCharmSemiElPbPb5500: decayer->SetForceDecay(kSemiElectronic); break;
366 case kBeautySemiElPbPb5500: decayer->SetForceDecay(kSemiElectronic); break;
367 default: decayer->SetForceDecay(kAll); break;
370 gMC->SetExternalDecayer(decayer);
373 //__________________________________________________________________
374 // *********** STEERING parameters FOR ALICE SIMULATION ************
375 // - Specify event type to be tracked through the ALICE setup
376 // - All positions are in cm, angles in degrees, and P and E in GeV
377 gMC->SetProcess("DCAY",1);
378 gMC->SetProcess("PAIR",1);
379 gMC->SetProcess("COMP",1);
380 gMC->SetProcess("PHOT",1);
381 gMC->SetProcess("PFIS",0);
382 gMC->SetProcess("DRAY",0);
383 gMC->SetProcess("ANNI",1);
384 gMC->SetProcess("BREM",1);
385 gMC->SetProcess("MUNU",1);
386 gMC->SetProcess("CKOV",1);
387 gMC->SetProcess("HADR",1);
388 gMC->SetProcess("LOSS",2);
389 gMC->SetProcess("MULS",1);
390 gMC->SetProcess("RAYL",1);
392 Float_t cut = 1.e-3; // 1MeV cut by default
393 Float_t tofmax = 1.e10;
395 gMC->SetCut("CUTGAM", cut);
396 gMC->SetCut("CUTELE", cut);
397 gMC->SetCut("CUTNEU", cut);
398 gMC->SetCut("CUTHAD", cut);
399 gMC->SetCut("CUTMUO", cut);
400 gMC->SetCut("BCUTE", cut);
401 gMC->SetCut("BCUTM", cut);
402 gMC->SetCut("DCUTE", cut);
403 gMC->SetCut("DCUTM", cut);
404 gMC->SetCut("PPCUTM", cut);
405 gMC->SetCut("TOFMAX", tofmax);
408 //__________________________________________________________________
409 // Generator Configuration
410 AliGenerator* gener = GeneratorFactory(eg, rad, comment);
411 gener->SetOrigin(0, 0, 0); // vertex position
412 gener->SetSigma(0, 0, 5.3); // Sigma in (X,Y,Z) (cm) on IP position
413 gener->SetCutVertexZ(1.); // Truncate at 1 sigma
414 gener->SetVertexSmear(kPerEvent);
415 gener->SetTrackingFlag(1);
418 //__________________________________________________________________
423 case k2kG: comment = comment.Append(" | L3 field 0.2 T"); break;
424 case k4kG: comment = comment.Append(" | L3 field 0.4 T"); break;
425 case k5kG: comment = comment.Append(" | L3 field 0.5 T"); break;
429 case kGluonRadiation:
430 comment = comment.Append(" | Gluon Radiation On"); break;
432 comment = comment.Append(" | Gluon Radiation Off"); break;
436 case kHoles: comment = comment.Append(" | Holes for PHOS/RICH"); break;
437 default: comment = comment.Append(" | No holes for PHOS/RICH"); break;
440 std::cout << "\n\n Comment: " << comment << "\n" << std::endl;
442 //__________________________________________________________________
444 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag);
445 field->SetL3ConstField(0); //Using const. field in the barrel
447 gAlice->SetField(field);
449 //__________________________________________________________________
453 Bool_t useABSO = kFALSE;
454 Bool_t useCRT = kFALSE;
455 Bool_t useDIPO = kFALSE;
456 Bool_t useFMD = kTRUE;
457 Bool_t useFRAME = kFALSE;
458 Bool_t useHALL = kFALSE;
459 Bool_t useITS = kFALSE;
460 Bool_t useMAG = kFALSE;
461 Bool_t useMUON = kFALSE;
462 Bool_t usePHOS = kFALSE;
463 Bool_t usePIPE = kFALSE;
464 Bool_t usePMD = kFALSE;
465 Bool_t useRICH = kFALSE;
466 Bool_t useSHIL = kFALSE;
467 Bool_t useSTART = kFALSE;
468 Bool_t useTOF = kFALSE;
469 Bool_t useTPC = kFALSE;
470 Bool_t useTRD = kFALSE;
471 Bool_t useZDC = kFALSE;
472 Bool_t useEMCAL = kFALSE;
473 Bool_t useVZERO = kFALSE;
475 cout << "\t* Creating the detectors ..." << endl;
476 // ================= Alice BODY parameters =========================
477 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
481 // =================== MAG parameters ============================
482 // Start with Magnet since detector layouts may be depending on
483 // the selected Magnet dimensions
484 AliMAG *MAG = new AliMAG("MAG", "Magnet");
488 // =================== ABSO parameters ===========================
489 AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
493 // =================== DIPO parameters ===========================
494 AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
498 // =================== HALL parameters ===========================
499 AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
504 // ================== FRAME parameters ===========================
505 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
507 case kHoles: FRAME->SetHoles(1); break;
508 default: FRAME->SetHoles(0); break;
513 // ================== SHIL parameters ============================
514 AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
519 // ================== PIPE parameters ============================
520 AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
524 // =================== ITS parameters ============================
526 // As the innermost detector in ALICE, the Inner Tracking System
527 // "impacts" on almost all other detectors. This involves the fact
528 // that the ITS geometry still has several options to be followed
529 // in parallel in order to determine the best set-up which
530 // minimizes the induced background. All the geometries available
531 // to date are described in the following. Read carefully the
532 // comments and use the default version (the only one uncommented)
533 // unless you are making comparisons and you know what you are
534 // doing. In this case just uncomment the ITS geometry you want to
535 // use and run Aliroot.
537 // Detailed geometries:
541 // new AliITSv5symm("ITS", "Updated ITS TDR detailed version "
542 // "with symmetric services");
544 // new AliITSv5asymm("ITS","Updates ITS TDR detailed version "
545 // "with asymmetric services");
547 AliITSvPPRasymmFMD *ITS =
548 new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version "
549 "with asymmetric services");
550 // don't touch this parameter if you're not an ITS developer
551 ITS->SetMinorVersion(2);
552 // don't touch this parameter if you're not an ITS developer
553 ITS->SetReadDet(kTRUE);
554 // don't touch this parameter if you're not an ITS developer
555 // ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");
556 // detector thickness on layer 1 must be in the range [100,300]
557 ITS->SetThicknessDet1(200.);
558 // detector thickness on layer 2 must be in the range [100,300]
559 ITS->SetThicknessDet2(200.);
560 // chip thickness on layer 1 must be in the range [150,300]
561 ITS->SetThicknessChip1(200.);
562 // chip thickness on layer 2 must be in the range [150,300]
563 ITS->SetThicknessChip2(200.);
564 // 1 --> rails in ; 0 --> rails out
566 // 1 --> water ; 0 --> freon
567 ITS->SetCoolingFluid(1);
569 // Coarse geometries (warning: no hits are produced with these
570 // coarse geometries and they unuseful for reconstruction !):
573 // AliITSvPPRcoarseasymm *ITS =
574 // new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version "
575 // "with asymmetric services");
576 // 1 --> rails in ; 0 --> rails out
578 // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
579 // ITS->SetSupportMaterial(0);
582 // new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version "
583 // "with symmetric services");
584 // 1 --> rails in ; 0 --> rails out
586 // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
587 // ITS->SetSupportMaterial(0);
589 // Geant3 <-> EUCLID conversion
590 // ============================
592 // SetEUCLID is a flag to output (=1) or not to output (=0) both
593 // geometry and media to two ASCII files (called by default
594 // ITSgeometry.euc and ITSgeometry.tme) in a format understandable
595 // to the CAD system EUCLID. The default (=0) means that you dont
596 // want to use this facility.
602 // =================== TPC parameters ============================
604 // This allows the user to specify sectors for the SLOW (TPC
605 // geometry 2) Simulator. SecAL (SecAU) <0 means that ALL lower
606 // (upper) sectors are specified, any value other than that
607 // requires at least one sector (lower or upper)to be specified!
610 // sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
611 // sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
613 // SecLows - number of lower sectors specified (up to 6)
614 // SecUps - number of upper sectors specified (up to 12)
615 // Sens - sensitive strips for the Slow Simulator !!!
617 // This does NOT work if all S or L-sectors are specified, i.e.
618 // if SecAL or SecAU < 0
621 //----------------------------------------------------------------
622 // gROOT->LoadMacro("SetTPCParam.C");
623 // AliTPCParam *param = SetTPCParam();
624 AliTPC *TPC = new AliTPCv2("TPC", "Default");
628 // ================== TOF parameters =============================
629 AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
633 // ================== RICH parameters ============================
634 AliRICH *RICH = new AliRICHv1("RICH", "normal RICH");
639 // ================== ZDC parameters =============================
640 AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
644 // ================== TRD parameters =============================
645 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
647 // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
650 // With hole in front of PHOS
652 // With hole in front of RICH
656 AliTRDsim *TRDsim = TRD->CreateTR();
660 // =================== FMD parameters ============================
661 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
662 // FMD->UseDetailed(kFALSE);
663 // FMD->UseAssembly();
668 // =================== MUON parameters ===========================
669 AliMUON *MUON = new AliMUONv1("MUON", "default");
670 // MUON->AddGeometryBuilder(new AliMUONSt1GeometryBuilder(MUON));
671 // MUON->AddGeometryBuilder(new AliMUONSt2GeometryBuilder(MUON));
672 // MUON->AddGeometryBuilder(new AliMUONSlatGeometryBuilder(MUON));
673 // MUON->AddGeometryBuilder(new AliMUONTriggerGeometryBuilder(MUON));
677 // =================== PHOS parameters ===========================
678 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
682 // =================== PMD parameters ============================
683 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
687 // =================== START parameters ==========================
688 AliSTART *START = new AliSTARTv1("START", "START Detector");
692 // =================== EMCAL parameters ==========================
693 AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCAL_55_25");
697 // =================== CRT parameters ============================
698 AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
702 // =================== V0 parameters =============================
703 AliVZERO *VZERO = new AliVZEROv3("VZERO", "normal VZERO");
707 //____________________________________________________________________
708 Float_t EtaToTheta(Float_t arg)
710 return (180./TMath::Pi())*2.*TMath::ATan(TMath::Exp(-arg));
713 //____________________________________________________________________
715 LookupEG(const Char_t* name)
718 for (Int_t i = 0; i < kEgMax; i++) {
719 if (n == egName[i]) return i;
724 //____________________________________________________________________
726 GeneratorFactory(EG_t eg, Rad_t rad, TString& comment)
729 if (rad == kNoGluonRadiation) isw = 0;
732 AliGenerator * gGener = 0;
736 comment = comment.Append(":HIJINGparam test 50 particles");
737 AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
738 gener->SetMomentumRange(0, 999999.);
739 gener->SetPhiRange(0., 360.);
740 // Set pseudorapidity range from -8 to 8.
741 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
742 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
743 gener->SetThetaRange(thmin,thmax);
749 comment = comment.Append(":HIJINGparam N=8000");
750 AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
751 gener->SetMomentumRange(0, 999999.);
752 gener->SetPhiRange(0., 360.);
753 // Set pseudorapidity range from -8 to 8.
754 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
755 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
756 gener->SetThetaRange(thmin,thmax);
762 comment = comment.Append("HIJINGparam N=4000");
763 AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
764 gener->SetMomentumRange(0, 999999.);
765 gener->SetPhiRange(0., 360.);
766 // Set pseudorapidity range from -8 to 8.
767 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
768 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
769 gener->SetThetaRange(thmin,thmax);
775 comment = comment.Append("HIJINGparam N=2000");
776 AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
777 gener->SetMomentumRange(0, 999999.);
778 gener->SetPhiRange(0., 360.);
779 // Set pseudorapidity range from -8 to 8.
780 Float_t thmin = EtaToTheta(8); // theta min. <---> eta max
781 Float_t thmax = EtaToTheta(-8); // theta max. <---> eta min
782 gener->SetThetaRange(thmin,thmax);
788 comment = comment.Append("HIJINGparam N=100");
789 AliGenHIJINGpara *gener = new AliGenHIJINGpara(500);
790 gener->SetMomentumRange(0, 999999.);
791 gener->SetPhiRange(0., 360.);
792 // Set pseudorapidity range from -8 to 8.
793 Float_t thmin = EtaToTheta(6); // theta min. <---> eta max
794 Float_t thmax = EtaToTheta(2); // theta max. <---> eta min
795 gener->SetThetaRange(thmin,thmax);
804 comment = comment.Append("HIJING cent1");
805 AliGenHijing *gener = HijingStandard();
806 // impact parameter range
807 gener->SetImpactParameterRange(0., 5.);
813 comment = comment.Append("HIJING cent2");
814 AliGenHijing *gener = HijingStandard();
815 // impact parameter range
816 gener->SetImpactParameterRange(0., 2.);
825 comment = comment.Append("HIJING per1");
826 AliGenHijing *gener = HijingStandard();
827 // impact parameter range
828 gener->SetImpactParameterRange(5., 8.6);
834 comment = comment.Append("HIJING per2");
835 AliGenHijing *gener = HijingStandard();
836 // impact parameter range
837 gener->SetImpactParameterRange(8.6, 11.2);
843 comment = comment.Append("HIJING per3");
844 AliGenHijing *gener = HijingStandard();
845 // impact parameter range
846 gener->SetImpactParameterRange(11.2, 13.2);
852 comment = comment.Append("HIJING per4");
853 AliGenHijing *gener = HijingStandard();
854 // impact parameter range
855 gener->SetImpactParameterRange(13.2, 15.);
861 comment = comment.Append("HIJING per5");
862 AliGenHijing *gener = HijingStandard();
863 // impact parameter range
864 gener->SetImpactParameterRange(15., 100.);
873 comment = comment.Append("HIJING Jet 25 GeV");
874 AliGenHijing *gener = HijingStandard();
875 // impact parameter range
876 gener->SetImpactParameterRange(0., 5.);
878 gener->SetTrigger(1);
879 gener->SetPtJet(25.);
880 gener->SetRadiation(isw);
881 gener->SetSimpleJets(!isw);
882 gener->SetJetEtaRange(-0.3,0.3);
883 gener->SetJetPhiRange(75., 165.);
890 comment = comment.Append("HIJING Jet 50 GeV");
891 AliGenHijing *gener = HijingStandard();
892 // impact parameter range
893 gener->SetImpactParameterRange(0., 5.);
895 gener->SetTrigger(1);
896 gener->SetPtJet(50.);
897 gener->SetRadiation(isw);
898 gener->SetSimpleJets(!isw);
899 gener->SetJetEtaRange(-0.3,0.3);
900 gener->SetJetPhiRange(75., 165.);
907 comment = comment.Append("HIJING Jet 75 GeV");
908 AliGenHijing *gener = HijingStandard();
909 // impact parameter range
910 gener->SetImpactParameterRange(0., 5.);
912 gener->SetTrigger(1);
913 gener->SetPtJet(75.);
914 gener->SetRadiation(isw);
915 gener->SetSimpleJets(!isw);
916 gener->SetJetEtaRange(-0.3,0.3);
917 gener->SetJetPhiRange(75., 165.);
924 comment = comment.Append("HIJING Jet 100 GeV");
925 AliGenHijing *gener = HijingStandard();
926 // impact parameter range
927 gener->SetImpactParameterRange(0., 5.);
929 gener->SetTrigger(1);
930 gener->SetPtJet(100.);
931 gener->SetRadiation(isw);
932 gener->SetSimpleJets(!isw);
933 gener->SetJetEtaRange(-0.3,0.3);
934 gener->SetJetPhiRange(75., 165.);
941 comment = comment.Append("HIJING Jet 200 GeV");
942 AliGenHijing *gener = HijingStandard();
943 // impact parameter range
944 gener->SetImpactParameterRange(0., 5.);
946 gener->SetTrigger(1);
947 gener->SetPtJet(200.);
948 gener->SetRadiation(isw);
949 gener->SetSimpleJets(!isw);
950 gener->SetJetEtaRange(-0.3,0.3);
951 gener->SetJetPhiRange(75., 165.);
960 comment = comment.Append("HIJING Gamma 25 GeV");
961 AliGenHijing *gener = HijingStandard();
962 // impact parameter range
963 gener->SetImpactParameterRange(0., 5.);
965 gener->SetTrigger(2);
966 gener->SetPtJet(25.);
967 gener->SetRadiation(isw);
968 gener->SetSimpleJets(!isw);
969 gener->SetJetEtaRange(-0.12, 0.12);
970 gener->SetJetPhiRange(220., 320.);
977 comment = comment.Append("HIJING Gamma 50 GeV");
978 AliGenHijing *gener = HijingStandard();
979 // impact parameter range
980 gener->SetImpactParameterRange(0., 5.);
982 gener->SetTrigger(2);
983 gener->SetPtJet(50.);
984 gener->SetRadiation(isw);
985 gener->SetSimpleJets(!isw);
986 gener->SetJetEtaRange(-0.12, 0.12);
987 gener->SetJetPhiRange(220., 320.);
994 comment = comment.Append("HIJING Gamma 75 GeV");
995 AliGenHijing *gener = HijingStandard();
996 // impact parameter range
997 gener->SetImpactParameterRange(0., 5.);
999 gener->SetTrigger(2);
1000 gener->SetPtJet(75.);
1001 gener->SetRadiation(isw);
1002 gener->SetSimpleJets(!isw);
1003 gener->SetJetEtaRange(-0.12, 0.12);
1004 gener->SetJetPhiRange(220., 320.);
1011 comment = comment.Append("HIJING Gamma 100 GeV");
1012 AliGenHijing *gener = HijingStandard();
1013 // impact parameter range
1014 gener->SetImpactParameterRange(0., 5.);
1016 gener->SetTrigger(2);
1017 gener->SetPtJet(100.);
1018 gener->SetRadiation(isw);
1019 gener->SetSimpleJets(!isw);
1020 gener->SetJetEtaRange(-0.12, 0.12);
1021 gener->SetJetPhiRange(220., 320.);
1028 comment = comment.Append("HIJING Gamma 200 GeV");
1029 AliGenHijing *gener = HijingStandard();
1030 // impact parameter range
1031 gener->SetImpactParameterRange(0., 5.);
1033 gener->SetTrigger(2);
1034 gener->SetPtJet(200.);
1035 gener->SetRadiation(isw);
1036 gener->SetSimpleJets(!isw);
1037 gener->SetJetEtaRange(-0.12, 0.12);
1038 gener->SetJetPhiRange(220., 320.);
1044 comment = comment.Append("HIJING pA");
1046 AliGenCocktail *gener = new AliGenCocktail();
1048 AliGenHijing *hijing = new AliGenHijing(-1);
1049 // centre of mass energy
1050 hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
1051 // impact parameter range
1052 hijing->SetImpactParameterRange(0., 15.);
1054 hijing->SetReferenceFrame("CMS");
1055 hijing->SetBoostLHC(1);
1057 hijing->SetProjectile("P", 1, 1);
1058 hijing->SetTarget ("A", 208, 82);
1059 // tell hijing to keep the full parent child chain
1060 hijing->KeepFullEvent();
1061 // enable jet quenching
1062 hijing->SetJetQuenching(0);
1064 hijing->SetShadowing(1);
1065 // Don't track spectators
1066 hijing->SetSpectators(0);
1067 // kinematic selection
1068 hijing->SetSelectAll(0);
1070 AliGenSlowNucleons* gray = new AliGenSlowNucleons(1);
1071 AliSlowNucleonModel* model = new AliSlowNucleonModelExp();
1072 gray->SetSlowNucleonModel(model);
1074 gener->AddGenerator(hijing,"Hijing pPb", 1);
1075 gener->AddGenerator(gray, "Gray Particles",1);
1081 comment = comment.Append(":Pythia p-p @ 14 TeV");
1082 AliGenPythia *gener = new AliGenPythia(-1);
1083 gener->SetMomentumRange(0,999999);
1084 gener->SetThetaRange(0., 180.);
1085 gener->SetYRange(-12,12);
1086 gener->SetPtRange(0,1000);
1087 gener->SetProcess(kPyMb);
1088 gener->SetEnergyCMS(14000.);
1092 case kPythia6Jets20_24:
1094 comment = comment.Append(":Pythia jets 20-24 GeV @ 5.5 TeV");
1095 AliGenPythia * gener = new AliGenPythia(-1);
1096 gener->SetEnergyCMS(5500.);// Centre of mass energy
1097 gener->SetProcess(kPyJets);// Process type
1098 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1099 gener->SetJetPhiRange(0., 360.);
1100 gener->SetJetEtRange(10., 1000.);
1101 gener->SetGluonRadiation(1,1);
1102 // gener->SetPtKick(0.);
1103 // Structure function
1104 gener->SetStrucFunc(kCTEQ4L);
1105 gener->SetPtHard(20., 24.);// Pt transfer of the hard scattering
1106 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1107 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1111 case kPythia6Jets24_29:
1113 comment = comment.Append(":Pythia jets 24-29 GeV @ 5.5 TeV");
1114 AliGenPythia * gener = new AliGenPythia(-1);
1115 gener->SetEnergyCMS(5500.);// Centre of mass energy
1116 gener->SetProcess(kPyJets);// Process type
1117 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1118 gener->SetJetPhiRange(0., 360.);
1119 gener->SetJetEtRange(10., 1000.);
1120 gener->SetGluonRadiation(1,1);
1121 // gener->SetPtKick(0.);
1122 // Structure function
1123 gener->SetStrucFunc(kCTEQ4L);
1124 gener->SetPtHard(24., 29.);// Pt transfer of the hard scattering
1125 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1126 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1130 case kPythia6Jets29_35:
1132 comment = comment.Append(":Pythia jets 29-35 GeV @ 5.5 TeV");
1133 AliGenPythia * gener = new AliGenPythia(-1);
1134 gener->SetEnergyCMS(5500.);// Centre of mass energy
1135 gener->SetProcess(kPyJets);// Process type
1136 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1137 gener->SetJetPhiRange(0., 360.);
1138 gener->SetJetEtRange(10., 1000.);
1139 gener->SetGluonRadiation(1,1);
1140 // gener->SetPtKick(0.);
1141 // Structure function
1142 gener->SetStrucFunc(kCTEQ4L);
1143 gener->SetPtHard(29., 35.);// Pt transfer of the hard scattering
1144 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1145 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1149 case kPythia6Jets35_42:
1151 comment = comment.Append(":Pythia jets 35-42 GeV @ 5.5 TeV");
1152 AliGenPythia * gener = new AliGenPythia(-1);
1153 gener->SetEnergyCMS(5500.);// Centre of mass energy
1154 gener->SetProcess(kPyJets);// Process type
1155 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1156 gener->SetJetPhiRange(0., 360.);
1157 gener->SetJetEtRange(10., 1000.);
1158 gener->SetGluonRadiation(1,1);
1159 // gener->SetPtKick(0.);
1160 // Structure function
1161 gener->SetStrucFunc(kCTEQ4L);
1162 gener->SetPtHard(35., 42.);// Pt transfer of the hard scattering
1163 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1164 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1168 case kPythia6Jets42_50:
1170 comment = comment.Append(":Pythia jets 42-50 GeV @ 5.5 TeV");
1171 AliGenPythia * gener = new AliGenPythia(-1);
1172 gener->SetEnergyCMS(5500.);// Centre of mass energy
1173 gener->SetProcess(kPyJets);// Process type
1174 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1175 gener->SetJetPhiRange(0., 360.);
1176 gener->SetJetEtRange(10., 1000.);
1177 gener->SetGluonRadiation(1,1);
1178 // gener->SetPtKick(0.);
1179 // Structure function
1180 gener->SetStrucFunc(kCTEQ4L);
1181 gener->SetPtHard(42., 50.);// Pt transfer of the hard scattering
1182 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1183 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1187 case kPythia6Jets50_60:
1189 comment = comment.Append(":Pythia jets 50-60 GeV @ 5.5 TeV");
1190 AliGenPythia * gener = new AliGenPythia(-1);
1191 gener->SetEnergyCMS(5500.);// Centre of mass energy
1192 gener->SetProcess(kPyJets);// Process type
1193 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1194 gener->SetJetPhiRange(0., 360.);
1195 gener->SetJetEtRange(10., 1000.);
1196 gener->SetGluonRadiation(1,1);
1197 // gener->SetPtKick(0.);
1198 // Structure function
1199 gener->SetStrucFunc(kCTEQ4L);
1200 gener->SetPtHard(50., 60.);// Pt transfer of the hard scattering
1201 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1202 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1206 case kPythia6Jets60_72:
1208 comment = comment.Append(":Pythia jets 60-72 GeV @ 5.5 TeV");
1209 AliGenPythia * gener = new AliGenPythia(-1);
1210 gener->SetEnergyCMS(5500.);// Centre of mass energy
1211 gener->SetProcess(kPyJets);// Process type
1212 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1213 gener->SetJetPhiRange(0., 360.);
1214 gener->SetJetEtRange(10., 1000.);
1215 gener->SetGluonRadiation(1,1);
1216 // gener->SetPtKick(0.);
1217 // Structure function
1218 gener->SetStrucFunc(kCTEQ4L);
1219 gener->SetPtHard(60., 72.);// Pt transfer of the hard scattering
1220 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1221 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1225 case kPythia6Jets72_86:
1227 comment = comment.Append(":Pythia jets 72-86 GeV @ 5.5 TeV");
1228 AliGenPythia * gener = new AliGenPythia(-1);
1229 gener->SetEnergyCMS(5500.);// Centre of mass energy
1230 gener->SetProcess(kPyJets);// Process type
1231 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1232 gener->SetJetPhiRange(0., 360.);
1233 gener->SetJetEtRange(10., 1000.);
1234 gener->SetGluonRadiation(1,1);
1235 // gener->SetPtKick(0.);
1236 // Structure function
1237 gener->SetStrucFunc(kCTEQ4L);
1238 gener->SetPtHard(72., 86.);// Pt transfer of the hard scattering
1239 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1240 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1244 case kPythia6Jets86_104:
1246 comment = comment.Append(":Pythia jets 86-104 GeV @ 5.5 TeV");
1247 AliGenPythia * gener = new AliGenPythia(-1);
1248 gener->SetEnergyCMS(5500.);// Centre of mass energy
1249 gener->SetProcess(kPyJets);// Process type
1250 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1251 gener->SetJetPhiRange(0., 360.);
1252 gener->SetJetEtRange(10., 1000.);
1253 gener->SetGluonRadiation(1,1);
1254 // gener->SetPtKick(0.);
1255 // Structure function
1256 gener->SetStrucFunc(kCTEQ4L);
1257 gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
1258 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1259 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1263 case kPythia6Jets104_125:
1265 comment = comment.Append(":Pythia jets 105-125 GeV @ 5.5 TeV");
1266 AliGenPythia * gener = new AliGenPythia(-1);
1267 gener->SetEnergyCMS(5500.);// Centre of mass energy
1268 gener->SetProcess(kPyJets);// Process type
1269 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1270 gener->SetJetPhiRange(0., 360.);
1271 gener->SetJetEtRange(10., 1000.);
1272 gener->SetGluonRadiation(1,1);
1273 // gener->SetPtKick(0.);
1274 // Structure function
1275 gener->SetStrucFunc(kCTEQ4L);
1276 gener->SetPtHard(104., 125.);// Pt transfer of the hard scattering
1277 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1278 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1282 case kPythia6Jets125_150:
1284 comment = comment.Append(":Pythia jets 125-150 GeV @ 5.5 TeV");
1285 AliGenPythia * gener = new AliGenPythia(-1);
1286 gener->SetEnergyCMS(5500.);// Centre of mass energy
1287 gener->SetProcess(kPyJets);// Process type
1288 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1289 gener->SetJetPhiRange(0., 360.);
1290 gener->SetJetEtRange(10., 1000.);
1291 gener->SetGluonRadiation(1,1);
1292 // gener->SetPtKick(0.);
1293 // Structure function
1294 gener->SetStrucFunc(kCTEQ4L);
1295 gener->SetPtHard(125., 150.);// Pt transfer of the hard scattering
1296 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1297 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1301 case kPythia6Jets150_180:
1303 comment = comment.Append(":Pythia jets 150-180 GeV @ 5.5 TeV");
1304 AliGenPythia * gener = new AliGenPythia(-1);
1305 gener->SetEnergyCMS(5500.);// Centre of mass energy
1306 gener->SetProcess(kPyJets);// Process type
1307 gener->SetJetEtaRange(-0.5, 0.5);// Final state kinematic cuts
1308 gener->SetJetPhiRange(0., 360.);
1309 gener->SetJetEtRange(10., 1000.);
1310 gener->SetGluonRadiation(1,1);
1311 // gener->SetPtKick(0.);
1312 // Structure function
1313 gener->SetStrucFunc(kCTEQ4L);
1314 gener->SetPtHard(150., 180.);// Pt transfer of the hard scattering
1315 gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
1316 gener->SetForceDecay(kAll);// Decay type (semielectronic, etc.)
1322 comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
1323 AliGenPythia * gener = new AliGenPythia(10);
1324 gener->SetProcess(kPyD0PbPbMNR);
1325 gener->SetStrucFunc(kCTEQ4L);
1326 gener->SetPtHard(2.1,-1.0);
1327 gener->SetEnergyCMS(5500.);
1328 gener->SetNuclei(208,208);
1329 gener->SetForceDecay(kHadronicD);
1330 gener->SetYRange(-2,2);
1331 gener->SetFeedDownHigherFamily(kFALSE);
1332 gener->SetStackFillOpt(AliGenPythia::kParentSelection);
1333 gener->SetCountMode(AliGenPythia::kCountParents);
1337 case kCharmSemiElPbPb5500:
1339 comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
1340 AliGenPythia * gener = new AliGenPythia(10);
1341 gener->SetProcess(kPyCharmPbPbMNR);
1342 gener->SetStrucFunc(kCTEQ4L);
1343 gener->SetPtHard(2.1,-1.0);
1344 gener->SetEnergyCMS(5500.);
1345 gener->SetNuclei(208,208);
1346 gener->SetForceDecay(kSemiElectronic);
1347 gener->SetYRange(-2,2);
1348 gener->SetFeedDownHigherFamily(kFALSE);
1349 gener->SetCountMode(AliGenPythia::kCountParents);
1353 case kBeautySemiElPbPb5500:
1355 comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
1356 AliGenPythia *gener = new AliGenPythia(10);
1357 gener->SetProcess(kPyBeautyPbPbMNR);
1358 gener->SetStrucFunc(kCTEQ4L);
1359 gener->SetPtHard(2.75,-1.0);
1360 gener->SetEnergyCMS(5500.);
1361 gener->SetNuclei(208,208);
1362 gener->SetForceDecay(kSemiElectronic);
1363 gener->SetYRange(-2,2);
1364 gener->SetFeedDownHigherFamily(kFALSE);
1365 gener->SetCountMode(AliGenPythia::kCountParents);
1371 comment = comment.Append(" Cocktail for TRD at 5.5 TeV");
1372 AliGenCocktail *gener = new AliGenCocktail();
1374 AliGenParam *jpsi = new AliGenParam(10,
1375 new AliGenMUONlib(),
1376 AliGenMUONlib::kJpsiFamily,
1379 jpsi->SetPtRange(0, 100);
1380 jpsi->SetYRange(-1., +1.);
1381 jpsi->SetForceDecay(kDiElectron);
1383 AliGenParam *ups = new AliGenParam(10,
1384 new AliGenMUONlib(),
1385 AliGenMUONlib::kUpsilonFamily,
1387 ups->SetPtRange(0, 100);
1388 ups->SetYRange(-1., +1.);
1389 ups->SetForceDecay(kDiElectron);
1391 AliGenParam *charm = new AliGenParam(10,
1392 new AliGenMUONlib(),
1393 AliGenMUONlib::kCharm,
1395 charm->SetPtRange(0, 100);
1396 charm->SetYRange(-1.5, +1.5);
1397 charm->SetForceDecay(kSemiElectronic);
1400 AliGenParam *beauty = new AliGenParam(10,
1401 new AliGenMUONlib(),
1402 AliGenMUONlib::kBeauty,
1404 beauty->SetPtRange(0, 100);
1405 beauty->SetYRange(-1.5, +1.5);
1406 beauty->SetForceDecay(kSemiElectronic);
1408 gener->AddGenerator(jpsi,"J/psi",1);
1409 gener->AddGenerator(ups,"Upsilon",1);
1410 gener->AddGenerator(charm,"Charm",1);
1411 gener->AddGenerator(beauty,"Beauty",1);
1417 comment = comment.Append(" Jet-jet at 5.5 TeV");
1418 AliGenPythia *gener = new AliGenPythia(-1);
1419 gener->SetEnergyCMS(5500.);
1420 gener->SetProcess(kPyJets);
1421 Double_t ptHardMin=10.0, ptHardMax=-1.0;
1422 gener->SetPtHard(ptHardMin,ptHardMax);
1423 gener->SetYHard(-0.7,0.7);
1424 gener->SetJetEtaRange(-0.2,0.2);
1425 gener->SetEventListRange(0,1);
1431 comment = comment.Append(" Gamma-jet at 5.5 TeV");
1432 AliGenPythia *gener = new AliGenPythia(-1);
1433 gener->SetEnergyCMS(5500.);
1434 gener->SetProcess(kPyDirectGamma);
1435 Double_t ptHardMin=10.0, ptHardMax=-1.0;
1436 gener->SetPtHard(ptHardMin,ptHardMax);
1437 gener->SetYHard(-1.0,1.0);
1438 gener->SetGammaEtaRange(-0.13,0.13);
1439 gener->SetGammaPhiRange(210.,330.);
1440 gener->SetEventListRange(0,1);
1444 case kMuonCocktailCent1:
1446 comment = comment.Append(" Muon Cocktail Cent1");
1447 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1448 gener->SetPtRange(1.0,100.); // Transverse momentum range
1449 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1450 gener->SetYRange(-4.0,-2.4);
1451 gener->SetMuonPtCut(0.8);
1452 gener->SetMuonThetaCut(171.,178.);
1453 gener->SetMuonMultiplicity(2);
1454 gener->SetNumberOfCollisions(1626.); //Centrality class Cent1 for PDC04
1455 gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1459 case kMuonCocktailPer1:
1461 comment = comment.Append(" Muon Cocktail Per1");
1462 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1463 gener->SetPtRange(1.0,100.); // Transverse momentum range
1464 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1465 gener->SetYRange(-4.0,-2.4);
1466 gener->SetMuonPtCut(0.8);
1467 gener->SetMuonThetaCut(171.,178.);
1468 gener->SetMuonMultiplicity(2);
1469 gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1470 gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1474 case kMuonCocktailPer4:
1476 comment = comment.Append(" Muon Cocktail Per4");
1477 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1478 gener->SetPtRange(1.0,100.); // Transverse momentum range
1479 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1480 gener->SetYRange(-4.0,-2.4);
1481 gener->SetMuonPtCut(0.8);
1482 gener->SetMuonThetaCut(171.,178.);
1483 gener->SetMuonMultiplicity(2);
1484 gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1485 gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1489 case kMuonCocktailCent1HighPt:
1491 comment = comment.Append(" Muon Cocktail HighPt Cent1");
1492 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1493 gener->SetPtRange(1.0,100.); // Transverse momentum range
1494 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1495 gener->SetYRange(-4.0,-2.4);
1496 gener->SetMuonPtCut(2.5);
1497 gener->SetMuonThetaCut(171.,178.);
1498 gener->SetMuonMultiplicity(2);
1499 gener->SetNumberOfCollisions(1626.); //Centrality class Cent1 for PDC04
1500 gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1504 case kMuonCocktailPer1HighPt :
1506 comment = comment.Append(" Muon Cocktail HighPt Per1");
1507 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1508 gener->SetPtRange(1.0,100.); // Transverse momentum range
1509 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1510 gener->SetYRange(-4.0,-2.4);
1511 gener->SetMuonPtCut(2.5);
1512 gener->SetMuonThetaCut(171.,178.);
1513 gener->SetMuonMultiplicity(2);
1514 gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1515 gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1519 case kMuonCocktailPer4HighPt:
1521 comment = comment.Append(" Muon Cocktail HighPt Per4");
1522 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1523 gener->SetPtRange(1.0,100.); // Transverse momentum range
1524 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1525 gener->SetYRange(-4.0,-2.4);
1526 gener->SetMuonPtCut(2.5);
1527 gener->SetMuonThetaCut(171.,178.);
1528 gener->SetMuonMultiplicity(2);
1529 gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1530 gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1534 case kMuonCocktailCent1Single:
1536 comment = comment.Append(" Muon Cocktail Single Cent1");
1537 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1538 gener->SetPtRange(1.0,100.); // Transverse momentum range
1539 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1540 gener->SetYRange(-4.0,-2.4);
1541 gener->SetMuonPtCut(0.8);
1542 gener->SetMuonThetaCut(171.,178.);
1543 gener->SetMuonMultiplicity(1);
1544 gener->SetNumberOfCollisions(1626.); //Centrality class Cent1 for PDC04
1545 gener->SetNumberOfParticipants(359.4);//Centrality class Cent1 for PDC04
1549 case kMuonCocktailPer1Single :
1551 comment = comment.Append(" Muon Cocktail Single Per1");
1552 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1553 gener->SetPtRange(1.0,100.); // Transverse momentum range
1554 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1555 gener->SetYRange(-4.0,-2.4);
1556 gener->SetMuonPtCut(0.8);
1557 gener->SetMuonThetaCut(171.,178.);
1558 gener->SetMuonMultiplicity(1);
1559 gener->SetNumberOfCollisions(820.0);//Centrality class Per1 for PDC04
1560 gener->SetNumberOfParticipants(229.3);//Centrality class Per1 for PDC04
1564 case kMuonCocktailPer4Single:
1566 comment = comment.Append(" Muon Cocktail Single Per4");
1567 AliGenMUONCocktail * gener = new AliGenMUONCocktail();
1568 gener->SetPtRange(1.0,100.); // Transverse momentum range
1569 gener->SetPhiRange(0.,360.); // Azimuthal angle range
1570 gener->SetYRange(-4.0,-2.4);
1571 gener->SetMuonPtCut(0.8);
1572 gener->SetMuonThetaCut(171.,178.);
1573 gener->SetMuonMultiplicity(1);
1574 gener->SetNumberOfCollisions(13.6);//Centrality class Per4 for PDC04
1575 gener->SetNumberOfParticipants(13.3);//Centrality class Per4 for PDC04
1581 comment = comment.Append(" Flat in FMD1 range");
1582 AliGenBox* gener = new AliGenBox(2000);
1583 gener->SetPart(211);
1584 gener->SetMomentumRange(3,4);
1585 gener->SetPhiRange(0, 360);
1586 gener->SetThetaRange(0.77, 3.08);
1592 comment = comment.Append(" Flat in FMD2 range");
1593 AliGenBox* gener = new AliGenBox(2000);
1594 gener->SetPart(211);
1595 gener->SetMomentumRange(3,4);
1596 gener->SetPhiRange(0, 360);
1597 gener->SetThetaRange(2.95, 20.42);
1603 comment = comment.Append(" Flat in FMD3 range");
1604 AliGenBox* gener = new AliGenBox(2000);
1605 gener->SetPart(211);
1606 gener->SetMomentumRange(3,4);
1607 gener->SetPhiRange(0, 360);
1608 gener->SetThetaRange(155.97, 176.73);
1614 comment = comment.Append(" Flat in FMD range");
1615 AliGenCocktail* gener = AliGenCocktail("FMD cocktail");
1616 gener->SetPart(211);
1617 gener->SetMomentumRange(3,4);
1618 gener->SetPhiRange(0, 360);
1619 AliGenBox* gener3 = new AliGenBox(2000);
1620 gener3->SetThetaRange(155.97, 176.73);
1621 gener->AddGenerator(gener3, "FMD3", .33);
1622 AliGenBox* gener2 = new AliGenBox(2000);
1623 gener2->SetThetaRange(2.95, 20.42);
1624 gener->AddGenerator(gener2, "FMD2", .33);
1625 AliGenBox* gener1 = new AliGenBox(2000);
1626 gener1->SetThetaRange(0.77, 3.08);
1627 gener->AddGenerator(gener1, "FMD1", .34);
1637 //____________________________________________________________________
1641 AliGenHijing *gener = new AliGenHijing(-1);
1642 // centre of mass energy
1643 gener->SetEnergyCMS(5500.);
1645 gener->SetReferenceFrame("CMS");
1647 gener->SetProjectile("A", 208, 82);
1648 gener->SetTarget ("A", 208, 82);
1649 // tell hijing to keep the full parent child chain
1650 gener->KeepFullEvent();
1651 // enable jet quenching
1652 gener->SetJetQuenching(1);
1654 gener->SetShadowing(1);
1655 // neutral pion and heavy particle decays switched off
1656 gener->SetDecaysOff(1);
1657 // Don't track spectators
1658 gener->SetSpectators(0);
1659 // kinematic selection
1660 gener->SetSelectAll(0);
1665 //____________________________________________________________________
1667 ProcessEnvironmentVars(EG_t& eg, Int_t& seed)
1670 if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
1671 Int_t eg1 = LookupEG(gSystem->Getenv("CONFIG_RUN_TYPE"));
1672 if (eg1 >= 0) eg = EG_t(eg1);
1674 // Random Number seed
1675 if (gSystem->Getenv("CONFIG_SEED")) {
1676 seed = atoi(gSystem->Getenv("CONFIG_SEED"));
1680 //____________________________________________________________________