2 // Configuration for the Physics Data Challenge 2007
5 // Assuming inel = 70 mb (PPR v.1, p.64)
7 // 84.98% of MSEL=0 events (including diffractive) with
8 // QQbar switched off (these events will include injected J/psi) => kPyMbNoHvq
10 // 14.14% of MSEL=1 events with ccbar (in 4 subsamples) => kCharmpp14000wmi
11 // bin 1 25% (3.535%): 2.76 < pthard < 3 GeV/c
12 // bin 2 40% (5.656%): 3 < pthard < 4 GeV/c
13 // bin 3 29% (4.101%): 4 < pthard < 8 GeV/c
14 // bin 4 6% (0.848%): pthard > 8 GeV/c
16 // 0.73% of MSEL=1 events with bbbar (in 4 subsamples) => kBeautypp14000wmi
17 // bin 1 5% (0.037%): 2.76 < pthard < 4 GeV/c
18 // bin 2 31% (0.226%): 4 < pthard < 6 GeV/c
19 // bin 3 28% (0.204%): 6 < pthard < 8 GeV/c
20 // bin 4 36% (0.263%): pthard >8 GeV/c
22 // 0.075% of MSEL=0 events with QQbar switched off and 1 Omega- => kPyOmegaMinus
23 // 0.075% of MSEL=0 events with QQbar switched off and 1 OmegaBar+ => kPyOmegaPlus
25 // One can use the configuration macro in compiled mode by
26 // root [0] gSystem->Load("libgeant321");
27 // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
28 // -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
29 // root [0] .x grun.C(1,"Config_PDC07.C++")
31 #if !defined(__CINT__) || defined(__MAKECINT__)
32 #include <Riostream.h>
36 #include <TVirtualMC.h>
37 #include <TGeant3TGeo.h>
38 #include "EVGEN/AliGenCocktail.h"
39 #include "EVGEN/AliGenParam.h"
40 #include "EVGEN/AliGenMUONlib.h"
41 #include "STEER/AliRunLoader.h"
42 #include "STEER/AliRun.h"
43 #include "STEER/AliConfig.h"
44 #include "PYTHIA6/AliDecayerPythia.h"
45 #include "PYTHIA6/AliGenPythia.h"
46 #include "STEER/AliMagFMaps.h"
47 #include "STRUCT/AliBODY.h"
48 #include "STRUCT/AliMAG.h"
49 #include "STRUCT/AliABSOv3.h"
50 #include "STRUCT/AliDIPOv3.h"
51 #include "STRUCT/AliHALLv3.h"
52 #include "STRUCT/AliFRAMEv2.h"
53 #include "STRUCT/AliSHILv3.h"
54 #include "STRUCT/AliPIPEv3.h"
55 #include "ITS/AliITSgeom.h"
56 #include "ITS/AliITSvPPRasymmFMD.h"
57 #include "TPC/AliTPCv2.h"
58 #include "TOF/AliTOFv6T0.h"
59 #include "HMPID/AliHMPIDv2.h"
60 #include "ZDC/AliZDCv3.h"
61 #include "TRD/AliTRDv1.h"
62 #include "FMD/AliFMDv1.h"
63 #include "MUON/AliMUONv1.h"
64 #include "PHOS/AliPHOSv1.h"
65 #include "PMD/AliPMDv1.h"
66 #include "T0/AliT0v1.h"
67 #include "EMCAL/AliEMCALv2.h"
68 #include "ACORDE/AliACORDEv0.h"
69 #include "VZERO/AliVZEROv7.h"
75 //--- Heavy Flavour Production ---
76 kCharmPbPb5500, kCharmpPb8800, kCharmpp14000, kCharmpp14000wmi,
77 kD0PbPb5500, kD0pPb8800, kD0pp14000,
78 kDPlusPbPb5500, kDPluspPb8800, kDPluspp14000,
79 kBeautyPbPb5500, kBeautypPb8800, kBeautypp14000, kBeautypp14000wmi,
81 kPyMbNoHvq, kPyOmegaPlus, kPyOmegaMinus, kPyJetJet,
82 kPyGammaJetPHOS, kPyJetJetPHOS, kPyJetJetPHOSv2, kPyGammaBremsPHOS,
83 kPyGammaJetEMCAL, kPyJetJetEMCAL, kPyGammaBremsEMCAL, kRunMax
86 const char * pprRunName[] = {
87 "kCharmPbPb5500", "kCharmpPb8800", "kCharmpp14000", "kCharmpp14000wmi",
88 "kD0PbPb5500", "kD0pPb8800", "kD0pp14000",
89 "kDPlusPbPb5500", "kDPluspPb8800", "kDPluspp14000",
90 "kBeautyPbPb5500", "kBeautypPb8800", "kBeautypp14000", "kBeautypp14000wmi",
91 "kPyMbNoHvq", "kPyOmegaPlus", "kPyOmegaMinus", "kPyJetJet",
92 "kPyGammaJetPHOS", "kPyJetJetPHOS", "kPyJetJetPHOSv2", "kPyGammaBremsPHOS",
93 "kPyGammaJetEMCAL", "kPyJetJetEMCAL", "kPyGammaBremsEMCAL"
100 kNature, kHadr, kSemiEl, kSemiMu
102 //--- Magnetic Field ---
108 //--- Trigger config ---
111 kDefaultPPTrig, kDefaultPbPbTrig
114 const char * TrigConfName[] = {
119 PDC07Proc_t proc = kPyJetJetEMCAL;
126 static PprGeo_t geo = kHoles;
129 AliGenPythia *PythiaHard(PDC07Proc_t );
130 AliGenPythia *PythiaHVQ(PDC07Proc_t );
131 AliGenerator *MbCocktail();
132 AliGenerator *PyMbTriggered(Int_t pdg);
133 void ProcessEnvironmentVars();
135 // This part for configuration
136 static DecayHvFl_t decHvFl = kNature;
137 static Mag_t mag = k5kG;
138 static TrigConf_t trig = kDefaultPPTrig; // default pp trigger configuration
139 static Int_t runNumber= 0;
140 static Int_t eventNumber= 0;
141 //========================//
142 // Set Random Number seed //
143 //========================//
145 static UInt_t seed = dt.Get();
147 // nEvts = -1 : you get 1 QQbar pair and all the fragmentation and
149 // nEvts = N>0 : you get N charm / beauty Hadrons
151 // stars = kTRUE : all heavy resonances and their decay stored
152 // = kFALSE: only final heavy hadrons and their decays stored
153 Bool_t stars = kTRUE;
155 // To be used only with kCharmpp14000wmi and kBeautypp14000wmi
156 // To get a "reasonable" agreement with MNR results, events have to be
157 // generated with the minimum ptHard set to 2.76 GeV.
158 // To get a "perfect" agreement with MNR results, events have to be
159 // generated in four ptHard bins with the following relative
172 Float_t ptHardMin = 10.;
173 Float_t ptHardMax = 20.;
174 Float_t ptGammaPi0Min = 1.;
175 Int_t iquenching = 0;
179 static TString comment;
185 // Get settings from environment variables
186 ProcessEnvironmentVars();
188 // libraries required by geant321
189 #if defined(__CINT__)
190 gSystem->Load("liblhapdf.so"); // Parton density functions
191 gSystem->Load("libpythia6.so"); // Pythia
192 gSystem->Load("libgeant321");
193 gSystem->Load("libEG");
194 gSystem->Load("libEGPythia6.so"); // TGenerator interface
195 gSystem->Load("libAliPythia6.so"); // ALICE specific implementations
197 new TGeant3TGeo("C++ Interface to Geant3");
199 // Output every 100 tracks
200 ((TGeant3*)gMC)->SetSWIT(4,100);
202 //=======================================================================
204 AliRunLoader* rl=0x0;
205 rl = AliRunLoader::Open("galice.root",
206 AliConfig::GetDefaultEventFolderName(),
210 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
213 rl->SetCompressionLevel(2);
214 rl->SetNumberOfEventsPerFile(1000);
215 gAlice->SetRunLoader(rl);
218 //gAlice->SetRunNumber(runNumber);
220 // Set the trigger configuration
221 gAlice->SetTriggerDescriptor(TrigConfName[trig]);
222 cout<<"Trigger configuration is set to "<<TrigConfName[trig]<<endl;
225 //=======================================================================
226 // ************* STEERING parameters FOR ALICE SIMULATION **************
227 // --- Specify event type to be tracked through the ALICE setup
228 // --- All positions are in cm, angles in degrees, and P and E in GeV
231 gMC->SetProcess("DCAY",1);
232 gMC->SetProcess("PAIR",1);
233 gMC->SetProcess("COMP",1);
234 gMC->SetProcess("PHOT",1);
235 gMC->SetProcess("PFIS",0);
236 gMC->SetProcess("DRAY",0);
237 gMC->SetProcess("ANNI",1);
238 gMC->SetProcess("BREM",1);
239 gMC->SetProcess("MUNU",1);
240 gMC->SetProcess("CKOV",1);
241 gMC->SetProcess("HADR",1);
242 gMC->SetProcess("LOSS",2);
243 gMC->SetProcess("MULS",1);
244 gMC->SetProcess("RAYL",1);
246 Float_t cut = 1.e-3; // 1MeV cut by default
247 Float_t tofmax = 1.e10;
249 gMC->SetCut("CUTGAM", cut);
250 gMC->SetCut("CUTELE", cut);
251 gMC->SetCut("CUTNEU", cut);
252 gMC->SetCut("CUTHAD", cut);
253 gMC->SetCut("CUTMUO", cut);
254 gMC->SetCut("BCUTE", cut);
255 gMC->SetCut("BCUTM", cut);
256 gMC->SetCut("DCUTE", cut);
257 gMC->SetCut("DCUTM", cut);
258 gMC->SetCut("PPCUTM", cut);
259 gMC->SetCut("TOFMAX", tofmax);
261 //======================//
262 // Set External decayer //
263 //======================//
264 TVirtualMCDecayer* decayer = new AliDecayerPythia();
269 decayer->SetForceDecay(kAll);
272 decayer->SetForceDecay(kHadronicD);
275 decayer->SetForceDecay(kSemiElectronic);
278 decayer->SetForceDecay(kSemiMuonic);
282 gMC->SetExternalDecayer(decayer);
283 if(proc == kPyJetJetPHOSv2) // in this case we decay the pi0
284 decayer->SetForceDecay(kNeutralPion);
286 //=========================//
287 // Generator Configuration //
288 //=========================//
289 AliGenPythia* gener = 0x0;
291 if (proc <= kBeautypp14000wmi) {
292 AliGenPythia *pythia = PythiaHVQ(proc);
294 pythia->SetFeedDownHigherFamily(kFALSE);
295 // Stack filling option
296 if(!stars) pythia->SetStackFillOpt(AliGenPythia::kParentSelection);
298 if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
304 pythia->SetForceDecay(kAll);
307 pythia->SetForceDecay(kHadronicD);
310 pythia->SetForceDecay(kSemiElectronic);
313 pythia->SetForceDecay(kSemiMuonic);
319 pythia->SetMomentumRange(0,99999999);
320 pythia->SetPhiRange(0., 360.);
321 pythia->SetThetaRange(0,180);
324 pythia->SetYRange(-999,999);
327 pythia->SetYRange(-2,2);
330 pythia->SetYRange(1,6);
334 } else if (proc == kPyMbNoHvq) {
335 gener = MbCocktail();
336 } else if (proc == kPyOmegaMinus) {
337 gener = PyMbTriggered(3334);
338 } else if (proc == kPyOmegaPlus) {
339 gener = PyMbTriggered(-3334);
340 } else if (proc <= kPyGammaBremsEMCAL) {
341 AliGenPythia *pythia = PythiaHard(proc);
343 pythia->SetFeedDownHigherFamily(kFALSE);
345 if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
350 pythia->SetMomentumRange(0,99999999);
351 // pythia->SetJetEtaRange(-1.5, 1.5);// Final state kinematic cuts
352 // pythia->SetJetPhiRange(0., 360.);
353 // pythia->SetThetaRange(45,135);
355 if(proc == kPyJetJetPHOSv2)
356 pythia->SetForceDecay(kNeutralPion);
358 pythia->SetForceDecay(kAll);
360 pythia->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
361 pythia->SetPtKick(5); // set the intrinsic kt to 5 GeV/c
368 gener->SetOrigin(0., 0., 0.); // vertex position
370 // Size of the interaction diamond
372 Float_t sigmaz = 7.55 / TMath::Sqrt(2.); // [cm]
375 Float_t betast = 10; // beta* [m]
376 Float_t eps = 3.75e-6; // emittance [m]
377 Float_t gamma = 7000. / 0.938272; // relativistic gamma [1]
378 Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.; // [cm]
379 printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
381 gener->SetSigma(sigmaxy, sigmaxy, sigmaz); // Sigma in (X,Y,Z) (cm) on IP position
382 gener->SetCutVertexZ(3.); // Truncate at 3 sigma
383 gener->SetVertexSmear(kPerEvent);
388 gener->SetQuench(iquenching);
390 Float_t k = 6e5*(qhat/1.7) ; //qhat=1.7, k = 6e5, default value
391 AliPythia::Instance()->InitQuenching(0.,0.1,k,0,0.95,6);
397 comment = comment.Append(" | L3 field 0.2 T");
398 } else if (mag == k4kG) {
399 comment = comment.Append(" | L3 field 0.4 T");
400 } else if (mag == k5kG) {
401 comment = comment.Append(" | L3 field 0.5 T");
403 printf("\n \n Comment: %s \n \n", comment.Data());
405 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag);
406 field->SetL3ConstField(0); //Using const. field in the barrel
408 gAlice->SetField(field);
435 //=================== Alice BODY parameters =============================
436 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
441 //=================== MAG parameters ============================
442 // --- Start with Magnet since detector layouts may be depending ---
443 // --- on the selected Magnet dimensions ---
444 AliMAG *MAG = new AliMAG("MAG", "Magnet");
450 //=================== ABSO parameters ============================
451 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
456 //=================== DIPO parameters ============================
458 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
463 //=================== HALL parameters ============================
465 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
471 //=================== FRAME parameters ============================
473 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
474 if (geo == kHoles) FRAME->SetHoles(1);
475 else FRAME->SetHoles(0);
481 //=================== SHIL parameters ============================
483 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
489 //=================== PIPE parameters ============================
491 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
496 //=================== ITS parameters ============================
499 AliITS *ITS = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
505 //============================ TPC parameters =====================
506 AliTPC *TPC = new AliTPCv2("TPC", "Default");
511 //=================== TOF parameters ============================
512 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
518 //=================== HMPID parameters ===========================
519 AliHMPID *HMPID = new AliHMPIDv2("HMPID", "normal HMPID");
526 //=================== ZDC parameters ============================
528 AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
533 //=================== TRD parameters ============================
535 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
536 AliTRDgeometry *geoTRD = TRD->GetGeometry();
541 //=================== FMD parameters ============================
542 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
547 //=================== MUON parameters ===========================
548 // New MUONv1 version (geometry defined via builders)
549 AliMUON *MUON = new AliMUONv1("MUON", "default");
551 //=================== PHOS parameters ===========================
555 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
561 //=================== PMD parameters ============================
562 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
567 //=================== T0 parameters ============================
568 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
573 //=================== EMCAL parameters ============================
574 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "ESHISH_77_TRD1_2X2_FINAL_110DEG scTh=0.176 pbTh=0.144");
579 //=================== CRT parameters ============================
580 AliACORDE *ACORDE = new AliACORDEv0("CRT", "normal ACORDE");
585 //=================== CRT parameters ============================
586 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
592 AliGenPythia *PythiaHVQ(PDC07Proc_t proc) {
593 //*******************************************************************//
594 // Configuration file for charm / beauty generation with PYTHIA //
596 // The parameters have been tuned in order to reproduce the inclusive//
597 // heavy quark pt distribution given by the NLO pQCD calculation by //
598 // Mangano, Nason and Ridolfi. //
600 // For details and for the NORMALIZATION of the yields see: //
601 // N.Carrer and A.Dainese, //
602 // "Charm and beauty production at the LHC", //
603 // ALICE-INT-2003-019, [arXiv:hep-ph/0311225]; //
604 // PPR Chapter 6.6, CERN/LHCC 2005-030 (2005). //
605 //*******************************************************************//
606 AliGenPythia * gener = 0x0;
610 comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
611 gener = new AliGenPythia(nEvts);
612 gener->SetProcess(kPyCharmPbPbMNR);
613 gener->SetStrucFunc(kCTEQ4L);
614 gener->SetPtHard(2.1,-1.0);
615 gener->SetEnergyCMS(5500.);
616 gener->SetNuclei(208,208);
619 comment = comment.Append(" Charm in p-Pb at 8.8 TeV");
620 gener = new AliGenPythia(nEvts);
621 gener->SetProcess(kPyCharmpPbMNR);
622 gener->SetStrucFunc(kCTEQ4L);
623 gener->SetPtHard(2.1,-1.0);
624 gener->SetEnergyCMS(8800.);
625 gener->SetProjectile("P",1,1);
626 gener->SetTarget("Pb",208,82);
629 comment = comment.Append(" Charm in pp at 14 TeV");
630 gener = new AliGenPythia(nEvts);
631 gener->SetProcess(kPyCharmppMNR);
632 gener->SetStrucFunc(kCTEQ4L);
633 gener->SetPtHard(2.1,-1.0);
634 gener->SetEnergyCMS(14000.);
636 case kCharmpp14000wmi:
637 comment = comment.Append(" Charm in pp at 14 TeV with mult. interactions");
638 gener = new AliGenPythia(-1);
639 gener->SetProcess(kPyCharmppMNRwmi);
640 gener->SetStrucFunc(kCTEQ5L);
641 gener->SetPtHard(ptHardMin,ptHardMax);
642 gener->SetEnergyCMS(14000.);
645 comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
646 gener = new AliGenPythia(nEvts);
647 gener->SetProcess(kPyD0PbPbMNR);
648 gener->SetStrucFunc(kCTEQ4L);
649 gener->SetPtHard(2.1,-1.0);
650 gener->SetEnergyCMS(5500.);
651 gener->SetNuclei(208,208);
654 comment = comment.Append(" D0 in p-Pb at 8.8 TeV");
655 gener = new AliGenPythia(nEvts);
656 gener->SetProcess(kPyD0pPbMNR);
657 gener->SetStrucFunc(kCTEQ4L);
658 gener->SetPtHard(2.1,-1.0);
659 gener->SetEnergyCMS(8800.);
660 gener->SetProjectile("P",1,1);
661 gener->SetTarget("Pb",208,82);
664 comment = comment.Append(" D0 in pp at 14 TeV");
665 gener = new AliGenPythia(nEvts);
666 gener->SetProcess(kPyD0ppMNR);
667 gener->SetStrucFunc(kCTEQ4L);
668 gener->SetPtHard(2.1,-1.0);
669 gener->SetEnergyCMS(14000.);
672 comment = comment.Append(" DPlus in Pb-Pb at 5.5 TeV");
673 gener = new AliGenPythia(nEvts);
674 gener->SetProcess(kPyDPlusPbPbMNR);
675 gener->SetStrucFunc(kCTEQ4L);
676 gener->SetPtHard(2.1,-1.0);
677 gener->SetEnergyCMS(5500.);
678 gener->SetNuclei(208,208);
681 comment = comment.Append(" DPlus in p-Pb at 8.8 TeV");
682 gener = new AliGenPythia(nEvts);
683 gener->SetProcess(kPyDPluspPbMNR);
684 gener->SetStrucFunc(kCTEQ4L);
685 gener->SetPtHard(2.1,-1.0);
686 gener->SetEnergyCMS(8800.);
687 gener->SetProjectile("P",1,1);
688 gener->SetTarget("Pb",208,82);
691 comment = comment.Append(" DPlus in pp at 14 TeV");
692 gener = new AliGenPythia(nEvts);
693 gener->SetProcess(kPyDPlusppMNR);
694 gener->SetStrucFunc(kCTEQ4L);
695 gener->SetPtHard(2.1,-1.0);
696 gener->SetEnergyCMS(14000.);
698 case kBeautyPbPb5500:
699 comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
700 gener = new AliGenPythia(nEvts);
701 gener->SetProcess(kPyBeautyPbPbMNR);
702 gener->SetStrucFunc(kCTEQ4L);
703 gener->SetPtHard(2.75,-1.0);
704 gener->SetEnergyCMS(5500.);
705 gener->SetNuclei(208,208);
708 comment = comment.Append(" Beauty in p-Pb at 8.8 TeV");
709 gener = new AliGenPythia(nEvts);
710 gener->SetProcess(kPyBeautypPbMNR);
711 gener->SetStrucFunc(kCTEQ4L);
712 gener->SetPtHard(2.75,-1.0);
713 gener->SetEnergyCMS(8800.);
714 gener->SetProjectile("P",1,1);
715 gener->SetTarget("Pb",208,82);
718 comment = comment.Append(" Beauty in pp at 14 TeV");
719 gener = new AliGenPythia(nEvts);
720 gener->SetProcess(kPyBeautyppMNR);
721 gener->SetStrucFunc(kCTEQ4L);
722 gener->SetPtHard(2.75,-1.0);
723 gener->SetEnergyCMS(14000.);
725 case kBeautypp14000wmi:
726 comment = comment.Append(" Beauty in pp at 14 TeV with mult. interactions");
727 gener = new AliGenPythia(-1);
728 gener->SetProcess(kPyBeautyppMNRwmi);
729 gener->SetStrucFunc(kCTEQ5L);
730 gener->SetPtHard(ptHardMin,ptHardMax);
731 gener->SetEnergyCMS(14000.);
738 AliGenPythia *PythiaHard(PDC07Proc_t proc) {
739 //*******************************************************************//
740 // Configuration file for hard QCD processes generation with PYTHIA //
742 //*******************************************************************//
743 AliGenPythia * gener = 0x0;
748 comment = comment.Append(" pp->jet + jet over at 14 TeV, no restriction");
749 AliGenPythia * gener = new AliGenPythia(nEvts);
750 gener->SetEnergyCMS(eCMS);// Centre of mass energy
751 gener->SetProcess(kPyJets);// Process type
752 gener->SetJetEtaRange(-1.5, 1.5);// Final state kinematic cuts
753 gener->SetJetPhiRange(0., 360.);
754 gener->SetJetEtRange(10., 1000.);
755 gener->SetPtHard(ptHardMin, ptHardMax);// Pt transfer of the hard scattering
756 gener->SetStrucFunc(kCTEQ4L);
760 case kPyGammaJetPHOS:
761 comment = comment.Append(" pp->jet + gamma over PHOS");
762 gener = new AliGenPythia(nEvts);
763 gener->SetEnergyCMS(eCMS);
764 gener->SetProcess(kPyDirectGamma);
765 gener->SetStrucFunc(kCTEQ4L);
766 gener->SetPtHard(ptHardMin,ptHardMax);
767 //gener->SetYHard(-1.0,1.0);
768 gener->SetGammaEtaRange(-0.13,0.13);
769 gener->SetGammaPhiRange(218.,322.);//Over 5 modules +-2 deg
772 comment = comment.Append(" pp->jet + jet over PHOS");
773 gener = new AliGenPythia(nEvts);
774 gener->SetEnergyCMS(eCMS);
775 gener->SetProcess(kPyJets);
776 gener->SetStrucFunc(kCTEQ4L);
777 gener->SetPtHard(ptHardMin,ptHardMax);
778 //gener->SetYHard(-1.0,1.0);
779 gener->SetJetEtaRange(-1.,1.);
780 gener->SetJetPhiRange(200.,340.);
781 gener->SetPi0InPHOS(kTRUE);
782 gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
784 printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
786 case kPyGammaBremsPHOS:
787 comment = comment.Append(" pp->jet + jet+bremsphoton over PHOS at 14 TeV");
788 gener = new AliGenPythia(nEvts);
789 gener->SetEnergyCMS(eCMS);
790 gener->SetProcess(kPyJets);
791 gener->SetStrucFunc(kCTEQ4L);
792 gener->SetPtHard(ptHardMin,ptHardMax);
793 //gener->SetYHard(-1.0,1.0);
794 gener->SetJetEtaRange(-1.,1.);
795 gener->SetJetPhiRange(200.,340.);
796 gener->SetFragPhotonInPHOS(kTRUE);
797 gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
798 printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
800 case kPyJetJetPHOSv2:
801 comment = comment.Append(" pp->jet + jet over PHOS version2 ");
802 gener = new AliGenPythia(nEvts);
803 gener->SetEnergyCMS(eCMS);
804 gener->SetProcess(kPyJets);
805 gener->SetStrucFunc(kCTEQ4L);
806 gener->SetPtHard(ptHardMin,ptHardMax);
807 //gener->SetYHard(-1.0,1.0);
808 gener->SetJetEtaRange(-1.,1.);
809 gener->SetJetPhiRange(200.,340.);
810 //gener->SetPi0InPHOS(kTRUE);
811 gener->SetPhotonInPHOSeta(kTRUE);
812 gener->SetPhotonMinPt(ptGammaPi0Min);
813 gener->SetForceDecay(kAll);
815 case kPyGammaJetEMCAL:
816 comment = comment.Append(" pp->jet + gamma over EMCAL at 14 TeV");
817 gener = new AliGenPythia(nEvts);
818 gener->SetEnergyCMS(eCMS);
819 gener->SetProcess(kPyDirectGamma);
820 gener->SetStrucFunc(kCTEQ4L);
821 gener->SetPtHard(ptHardMin,ptHardMax);
822 //gener->SetYHard(-1.0,1.0);
823 gener->SetGammaEtaRange(-0.71,0.71);
824 gener->SetGammaPhiRange(78.,192.);//Over 6 supermodules +-2 deg
827 comment = comment.Append(" pp->jet + jet over EMCAL at 14 TeV");
828 gener = new AliGenPythia(nEvts);
829 gener->SetEnergyCMS(eCMS);
830 gener->SetProcess(kPyJets);
831 gener->SetStrucFunc(kCTEQ4L);
832 gener->SetPtHard(ptHardMin,ptHardMax);
833 //gener->SetYHard(-1.0,1.0);
834 gener->SetJetEtaRange(-1,1);
835 gener->SetJetPhiRange(60.,210.);
836 gener->SetPi0InEMCAL(kTRUE);
837 gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
838 printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
840 case kPyGammaBremsEMCAL:
841 comment = comment.Append(" pp->jet + jet+bremsphoton over EMCAL at 14 TeV");
842 gener = new AliGenPythia(nEvts);
843 gener->SetEnergyCMS(eCMS);
844 gener->SetProcess(kPyJets);
845 gener->SetStrucFunc(kCTEQ4L);
846 gener->SetPtHard(ptHardMin,ptHardMax);
847 //gener->SetYHard(-1.0,1.0);
848 gener->SetJetEtaRange(-1,1);
849 gener->SetJetPhiRange(60.,210.); //Over 2 uncovered PHOS modules
850 gener->SetFragPhotonInEMCAL(kTRUE);
851 gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
853 printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
861 AliGenerator* MbCocktail()
863 comment = comment.Append(" pp at 14 TeV: Pythia low-pt, no heavy quarks + J/Psi from parameterisation");
864 AliGenCocktail * gener = new AliGenCocktail();
865 gener->UsePerEventRates();
869 AliGenPythia* pythia = new AliGenPythia(-1);
870 pythia->SetMomentumRange(0, 999999.);
871 pythia->SetThetaRange(0., 180.);
872 pythia->SetYRange(-12.,12.);
873 pythia->SetPtRange(0,1000.);
874 pythia->SetProcess(kPyMb);
875 pythia->SetEnergyCMS(14000.);
876 pythia->SwitchHFOff();
878 // J/Psi parameterisation
880 AliGenParam* jpsi = new AliGenParam(1, AliGenMUONlib::kJpsi, "CDF scaled", "Jpsi");
881 jpsi->SetPtRange(0.,100.);
882 jpsi->SetYRange(-8., 8.);
883 jpsi->SetPhiRange(0., 360.);
884 jpsi->SetForceDecay(kAll);
886 gener->AddGenerator(pythia, "Pythia", 1.);
887 gener->AddGenerator(jpsi, "J/Psi", 8.e-4);
892 AliGenerator* PyMbTriggered(Int_t pdg)
894 AliGenPythia* pythia = new AliGenPythia(-1);
895 pythia->SetMomentumRange(0, 999999.);
896 pythia->SetThetaRange(0., 180.);
897 pythia->SetYRange(-12.,12.);
898 pythia->SetPtRange(0,1000.);
899 pythia->SetProcess(kPyMb);
900 pythia->SetEnergyCMS(14000.);
901 pythia->SetTriggerParticle(pdg, 0.9);
905 void ProcessEnvironmentVars()
907 cout << "######################################" << endl;
908 cout << "## Processing environment variables ##" << endl;
909 cout << "######################################" << endl;
912 if (gSystem->Getenv("DC_RUN")) {
913 runNumber = atoi(gSystem->Getenv("DC_RUN"));
915 //cout<<"Run number "<<runNumber<<endl;
918 if (gSystem->Getenv("DC_EVENT")) {
919 eventNumber = atoi(gSystem->Getenv("DC_EVENT"));
921 //cout<<"Event number "<<eventNumber<<endl;
923 // Random Number seed
924 if (gSystem->Getenv("CONFIG_SEED")) {
925 seed = atoi(gSystem->Getenv("CONFIG_SEED"));
927 else if(gSystem->Getenv("DC_EVENT") && gSystem->Getenv("DC_RUN")){
928 seed = runNumber * 100000 + eventNumber;
931 gRandom->SetSeed(seed);
933 cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
934 cout<<"Seed for random number generation= "<< seed<<" "<< gRandom->GetSeed()<<endl;
935 cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
938 if (gSystem->Getenv("DC_RUN_TYPE")) {
939 cout<<"run type "<<gSystem->Getenv("DC_RUN_TYPE")<<endl;
940 for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
941 if (strcmp(gSystem->Getenv("DC_RUN_TYPE"), pprRunName[iRun])==0) {
942 proc = (PDC07Proc_t)iRun;
947 cout << "Environment variable DC_RUN_TYPE is not defined" << endl;
949 if (gSystem->Getenv("ECMS"))
950 eCMS = atof(gSystem->Getenv("ECMS"));
951 if (gSystem->Getenv("PTHARDMIN"))
952 ptHardMin = atof(gSystem->Getenv("PTHARDMIN"));
953 if (gSystem->Getenv("PTHARDMAX"))
954 ptHardMax = atof(gSystem->Getenv("PTHARDMAX"));
955 if (gSystem->Getenv("PTGAMMAPI0MIN"))
956 ptGammaPi0Min = atof(gSystem->Getenv("PTGAMMAPI0MIN"));
957 if (gSystem->Getenv("QUENCHING"))
958 iquenching = atof(gSystem->Getenv("QUENCHING"));
959 if (gSystem->Getenv("QHAT"))
960 qhat = atof(gSystem->Getenv("QHAT"));
962 cout<<">> Run type set to "<<pprRunName[proc]<<endl;
963 cout<<">> Collision energy set to "<<eCMS <<endl;
964 cout<<">> ptHard limits: "<<ptHardMin<<" to " <<ptHardMax<<" GeV"<<endl;
965 cout<<">> pt gamma/pi0 threshold "<< ptGammaPi0Min<<" GeV "<<endl;
966 cout<<">> quenching on? "<< iquenching<<" qhat "<<qhat<<endl;
968 cout << "######################################" << endl;
969 cout << "######################################" << endl;