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/AliHMPIDv3.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/AliACORDEv1.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;
183 // Get settings from environment variables
184 ProcessEnvironmentVars();
186 // libraries required by geant321
187 #if defined(__CINT__)
188 gSystem->Load("liblhapdf.so"); // Parton density functions
189 gSystem->Load("libpythia6.so"); // Pythia
190 gSystem->Load("libgeant321");
191 gSystem->Load("libEG");
192 gSystem->Load("libEGPythia6.so"); // TGenerator interface
193 gSystem->Load("libAliPythia6.so"); // ALICE specific implementations
195 new TGeant3TGeo("C++ Interface to Geant3");
197 // Output every 100 tracks
198 ((TGeant3*)gMC)->SetSWIT(4,100);
200 //=======================================================================
202 AliRunLoader* rl=0x0;
203 rl = AliRunLoader::Open("galice.root",
204 AliConfig::GetDefaultEventFolderName(),
208 gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
211 rl->SetCompressionLevel(2);
212 rl->SetNumberOfEventsPerFile(1000);
213 gAlice->SetRunLoader(rl);
216 //gAlice->SetRunNumber(runNumber);
218 // Set the trigger configuration
219 gAlice->SetTriggerDescriptor(TrigConfName[trig]);
220 cout<<"Trigger configuration is set to "<<TrigConfName[trig]<<endl;
223 //=======================================================================
224 // ************* STEERING parameters FOR ALICE SIMULATION **************
225 // --- Specify event type to be tracked through the ALICE setup
226 // --- All positions are in cm, angles in degrees, and P and E in GeV
229 gMC->SetProcess("DCAY",1);
230 gMC->SetProcess("PAIR",1);
231 gMC->SetProcess("COMP",1);
232 gMC->SetProcess("PHOT",1);
233 gMC->SetProcess("PFIS",0);
234 gMC->SetProcess("DRAY",0);
235 gMC->SetProcess("ANNI",1);
236 gMC->SetProcess("BREM",1);
237 gMC->SetProcess("MUNU",1);
238 gMC->SetProcess("CKOV",1);
239 gMC->SetProcess("HADR",1);
240 gMC->SetProcess("LOSS",2);
241 gMC->SetProcess("MULS",1);
242 gMC->SetProcess("RAYL",1);
244 Float_t cut = 1.e-3; // 1MeV cut by default
245 Float_t tofmax = 1.e10;
247 gMC->SetCut("CUTGAM", cut);
248 gMC->SetCut("CUTELE", cut);
249 gMC->SetCut("CUTNEU", cut);
250 gMC->SetCut("CUTHAD", cut);
251 gMC->SetCut("CUTMUO", cut);
252 gMC->SetCut("BCUTE", cut);
253 gMC->SetCut("BCUTM", cut);
254 gMC->SetCut("DCUTE", cut);
255 gMC->SetCut("DCUTM", cut);
256 gMC->SetCut("PPCUTM", cut);
257 gMC->SetCut("TOFMAX", tofmax);
259 //======================//
260 // Set External decayer //
261 //======================//
262 TVirtualMCDecayer* decayer = new AliDecayerPythia();
267 decayer->SetForceDecay(kAll);
270 decayer->SetForceDecay(kHadronicD);
273 decayer->SetForceDecay(kSemiElectronic);
276 decayer->SetForceDecay(kSemiMuonic);
280 gMC->SetExternalDecayer(decayer);
281 if(proc == kPyJetJetPHOSv2) // in this case we decay the pi0
282 decayer->SetForceDecay(kNeutralPion);
284 //=========================//
285 // Generator Configuration //
286 //=========================//
287 AliGenPythia* gener = 0x0;
289 if (proc <= kBeautypp14000wmi) {
290 AliGenPythia *pythia = PythiaHVQ(proc);
292 pythia->SetFeedDownHigherFamily(kFALSE);
293 // Stack filling option
294 if(!stars) pythia->SetStackFillOpt(AliGenPythia::kParentSelection);
296 if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
302 pythia->SetForceDecay(kAll);
305 pythia->SetForceDecay(kHadronicD);
308 pythia->SetForceDecay(kSemiElectronic);
311 pythia->SetForceDecay(kSemiMuonic);
317 pythia->SetMomentumRange(0,99999999);
318 pythia->SetPhiRange(0., 360.);
319 pythia->SetThetaRange(0,180);
322 pythia->SetYRange(-999,999);
325 pythia->SetYRange(-2,2);
328 pythia->SetYRange(1,6);
332 } else if (proc == kPyMbNoHvq) {
333 gener = MbCocktail();
334 } else if (proc == kPyOmegaMinus) {
335 gener = PyMbTriggered(3334);
336 } else if (proc == kPyOmegaPlus) {
337 gener = PyMbTriggered(-3334);
338 } else if (proc <= kPyGammaBremsEMCAL) {
339 AliGenPythia *pythia = PythiaHard(proc);
341 pythia->SetFeedDownHigherFamily(kFALSE);
343 if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
348 pythia->SetMomentumRange(0,99999999);
349 // pythia->SetJetEtaRange(-1.5, 1.5);// Final state kinematic cuts
350 // pythia->SetJetPhiRange(0., 360.);
351 // pythia->SetThetaRange(45,135);
353 if(proc == kPyJetJetPHOSv2)
354 pythia->SetForceDecay(kNeutralPion);
356 pythia->SetForceDecay(kAll);
358 pythia->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
359 pythia->SetPtKick(5); // set the intrinsic kt to 5 GeV/c
366 gener->SetOrigin(0., 0., 0.); // vertex position
368 // Size of the interaction diamond
370 Float_t sigmaz = 7.55 / TMath::Sqrt(2.); // [cm]
373 Float_t betast = 10; // beta* [m]
374 Float_t eps = 3.75e-6; // emittance [m]
375 Float_t gamma = 7000. / 0.938272; // relativistic gamma [1]
376 Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.; // [cm]
377 printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
379 gener->SetSigma(sigmaxy, sigmaxy, sigmaz); // Sigma in (X,Y,Z) (cm) on IP position
380 gener->SetCutVertexZ(3.); // Truncate at 3 sigma
381 gener->SetVertexSmear(kPerEvent);
386 gener->SetQuench(iquenching);
388 Float_t k = 6e5*(qhat/1.7) ; //qhat=1.7, k = 6e5, default value
389 AliPythia::Instance()->InitQuenching(0.,0.1,k,0,0.95,6);
395 comment = comment.Append(" | L3 field 0.2 T");
396 } else if (mag == k4kG) {
397 comment = comment.Append(" | L3 field 0.4 T");
398 } else if (mag == k5kG) {
399 comment = comment.Append(" | L3 field 0.5 T");
401 printf("\n \n Comment: %s \n \n", comment.Data());
403 AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag);
404 field->SetL3ConstField(0); //Using const. field in the barrel
406 gAlice->SetField(field);
433 //=================== Alice BODY parameters =============================
434 AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
439 //=================== MAG parameters ============================
440 // --- Start with Magnet since detector layouts may be depending ---
441 // --- on the selected Magnet dimensions ---
442 AliMAG *MAG = new AliMAG("MAG", "Magnet");
448 //=================== ABSO parameters ============================
449 AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
454 //=================== DIPO parameters ============================
456 AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
461 //=================== HALL parameters ============================
463 AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
469 //=================== FRAME parameters ============================
471 AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
472 if (geo == kHoles) FRAME->SetHoles(1);
473 else FRAME->SetHoles(0);
479 //=================== SHIL parameters ============================
481 AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
487 //=================== PIPE parameters ============================
489 AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
494 //=================== ITS parameters ============================
497 AliITS *ITS = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
503 //============================ TPC parameters =====================
504 AliTPC *TPC = new AliTPCv2("TPC", "Default");
509 //=================== TOF parameters ============================
510 AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
516 //=================== HMPID parameters ===========================
517 AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
524 //=================== ZDC parameters ============================
526 AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
531 //=================== TRD parameters ============================
533 AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
534 AliTRDgeometry *geoTRD = TRD->GetGeometry();
539 //=================== FMD parameters ============================
540 AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
545 //=================== MUON parameters ===========================
546 // New MUONv1 version (geometry defined via builders)
547 AliMUON *MUON = new AliMUONv1("MUON", "default");
549 //=================== PHOS parameters ===========================
553 AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
559 //=================== PMD parameters ============================
560 AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
565 //=================== T0 parameters ============================
566 AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
571 //=================== EMCAL parameters ============================
572 AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE scTh=0.176 pbTh=0.144");
577 //=================== CRT parameters ============================
578 AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
583 //=================== CRT parameters ============================
584 AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
590 AliGenPythia *PythiaHVQ(PDC07Proc_t proc) {
591 //*******************************************************************//
592 // Configuration file for charm / beauty generation with PYTHIA //
594 // The parameters have been tuned in order to reproduce the inclusive//
595 // heavy quark pt distribution given by the NLO pQCD calculation by //
596 // Mangano, Nason and Ridolfi. //
598 // For details and for the NORMALIZATION of the yields see: //
599 // N.Carrer and A.Dainese, //
600 // "Charm and beauty production at the LHC", //
601 // ALICE-INT-2003-019, [arXiv:hep-ph/0311225]; //
602 // PPR Chapter 6.6, CERN/LHCC 2005-030 (2005). //
603 //*******************************************************************//
604 AliGenPythia * gener = 0x0;
608 comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
609 gener = new AliGenPythia(nEvts);
610 gener->SetProcess(kPyCharmPbPbMNR);
611 gener->SetStrucFunc(kCTEQ4L);
612 gener->SetPtHard(2.1,-1.0);
613 gener->SetEnergyCMS(5500.);
614 gener->SetNuclei(208,208);
617 comment = comment.Append(" Charm in p-Pb at 8.8 TeV");
618 gener = new AliGenPythia(nEvts);
619 gener->SetProcess(kPyCharmpPbMNR);
620 gener->SetStrucFunc(kCTEQ4L);
621 gener->SetPtHard(2.1,-1.0);
622 gener->SetEnergyCMS(8800.);
623 gener->SetProjectile("P",1,1);
624 gener->SetTarget("Pb",208,82);
627 comment = comment.Append(" Charm in pp at 14 TeV");
628 gener = new AliGenPythia(nEvts);
629 gener->SetProcess(kPyCharmppMNR);
630 gener->SetStrucFunc(kCTEQ4L);
631 gener->SetPtHard(2.1,-1.0);
632 gener->SetEnergyCMS(14000.);
634 case kCharmpp14000wmi:
635 comment = comment.Append(" Charm in pp at 14 TeV with mult. interactions");
636 gener = new AliGenPythia(-1);
637 gener->SetProcess(kPyCharmppMNRwmi);
638 gener->SetStrucFunc(kCTEQ5L);
639 gener->SetPtHard(ptHardMin,ptHardMax);
640 gener->SetEnergyCMS(14000.);
643 comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
644 gener = new AliGenPythia(nEvts);
645 gener->SetProcess(kPyD0PbPbMNR);
646 gener->SetStrucFunc(kCTEQ4L);
647 gener->SetPtHard(2.1,-1.0);
648 gener->SetEnergyCMS(5500.);
649 gener->SetNuclei(208,208);
652 comment = comment.Append(" D0 in p-Pb at 8.8 TeV");
653 gener = new AliGenPythia(nEvts);
654 gener->SetProcess(kPyD0pPbMNR);
655 gener->SetStrucFunc(kCTEQ4L);
656 gener->SetPtHard(2.1,-1.0);
657 gener->SetEnergyCMS(8800.);
658 gener->SetProjectile("P",1,1);
659 gener->SetTarget("Pb",208,82);
662 comment = comment.Append(" D0 in pp at 14 TeV");
663 gener = new AliGenPythia(nEvts);
664 gener->SetProcess(kPyD0ppMNR);
665 gener->SetStrucFunc(kCTEQ4L);
666 gener->SetPtHard(2.1,-1.0);
667 gener->SetEnergyCMS(14000.);
670 comment = comment.Append(" DPlus in Pb-Pb at 5.5 TeV");
671 gener = new AliGenPythia(nEvts);
672 gener->SetProcess(kPyDPlusPbPbMNR);
673 gener->SetStrucFunc(kCTEQ4L);
674 gener->SetPtHard(2.1,-1.0);
675 gener->SetEnergyCMS(5500.);
676 gener->SetNuclei(208,208);
679 comment = comment.Append(" DPlus in p-Pb at 8.8 TeV");
680 gener = new AliGenPythia(nEvts);
681 gener->SetProcess(kPyDPluspPbMNR);
682 gener->SetStrucFunc(kCTEQ4L);
683 gener->SetPtHard(2.1,-1.0);
684 gener->SetEnergyCMS(8800.);
685 gener->SetProjectile("P",1,1);
686 gener->SetTarget("Pb",208,82);
689 comment = comment.Append(" DPlus in pp at 14 TeV");
690 gener = new AliGenPythia(nEvts);
691 gener->SetProcess(kPyDPlusppMNR);
692 gener->SetStrucFunc(kCTEQ4L);
693 gener->SetPtHard(2.1,-1.0);
694 gener->SetEnergyCMS(14000.);
696 case kBeautyPbPb5500:
697 comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
698 gener = new AliGenPythia(nEvts);
699 gener->SetProcess(kPyBeautyPbPbMNR);
700 gener->SetStrucFunc(kCTEQ4L);
701 gener->SetPtHard(2.75,-1.0);
702 gener->SetEnergyCMS(5500.);
703 gener->SetNuclei(208,208);
706 comment = comment.Append(" Beauty in p-Pb at 8.8 TeV");
707 gener = new AliGenPythia(nEvts);
708 gener->SetProcess(kPyBeautypPbMNR);
709 gener->SetStrucFunc(kCTEQ4L);
710 gener->SetPtHard(2.75,-1.0);
711 gener->SetEnergyCMS(8800.);
712 gener->SetProjectile("P",1,1);
713 gener->SetTarget("Pb",208,82);
716 comment = comment.Append(" Beauty in pp at 14 TeV");
717 gener = new AliGenPythia(nEvts);
718 gener->SetProcess(kPyBeautyppMNR);
719 gener->SetStrucFunc(kCTEQ4L);
720 gener->SetPtHard(2.75,-1.0);
721 gener->SetEnergyCMS(14000.);
723 case kBeautypp14000wmi:
724 comment = comment.Append(" Beauty in pp at 14 TeV with mult. interactions");
725 gener = new AliGenPythia(-1);
726 gener->SetProcess(kPyBeautyppMNRwmi);
727 gener->SetStrucFunc(kCTEQ5L);
728 gener->SetPtHard(ptHardMin,ptHardMax);
729 gener->SetEnergyCMS(14000.);
736 AliGenPythia *PythiaHard(PDC07Proc_t proc) {
737 //*******************************************************************//
738 // Configuration file for hard QCD processes generation with PYTHIA //
740 //*******************************************************************//
741 AliGenPythia * gener = 0x0;
746 comment = comment.Append(" pp->jet + jet over at 14 TeV, no restriction");
747 AliGenPythia * gener = new AliGenPythia(nEvts);
748 gener->SetEnergyCMS(eCMS);// Centre of mass energy
749 gener->SetProcess(kPyJets);// Process type
750 gener->SetJetEtaRange(-1.5, 1.5);// Final state kinematic cuts
751 gener->SetJetPhiRange(0., 360.);
752 gener->SetJetEtRange(10., 1000.);
753 gener->SetPtHard(ptHardMin, ptHardMax);// Pt transfer of the hard scattering
754 gener->SetStrucFunc(kCTEQ4L);
758 case kPyGammaJetPHOS:
759 comment = comment.Append(" pp->jet + gamma over PHOS");
760 gener = new AliGenPythia(nEvts);
761 gener->SetEnergyCMS(eCMS);
762 gener->SetProcess(kPyDirectGamma);
763 gener->SetStrucFunc(kCTEQ4L);
764 gener->SetPtHard(ptHardMin,ptHardMax);
765 //gener->SetYHard(-1.0,1.0);
766 gener->SetGammaEtaRange(-0.13,0.13);
767 gener->SetGammaPhiRange(218.,322.);//Over 5 modules +-2 deg
770 comment = comment.Append(" pp->jet + jet over PHOS");
771 gener = new AliGenPythia(nEvts);
772 gener->SetEnergyCMS(eCMS);
773 gener->SetProcess(kPyJets);
774 gener->SetStrucFunc(kCTEQ4L);
775 gener->SetPtHard(ptHardMin,ptHardMax);
776 //gener->SetYHard(-1.0,1.0);
777 gener->SetJetEtaRange(-1.,1.);
778 gener->SetJetPhiRange(200.,340.);
779 gener->SetPi0InPHOS(kTRUE);
780 gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
782 printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
784 case kPyGammaBremsPHOS:
785 comment = comment.Append(" pp->jet + jet+bremsphoton over PHOS at 14 TeV");
786 gener = new AliGenPythia(nEvts);
787 gener->SetEnergyCMS(eCMS);
788 gener->SetProcess(kPyJets);
789 gener->SetStrucFunc(kCTEQ4L);
790 gener->SetPtHard(ptHardMin,ptHardMax);
791 //gener->SetYHard(-1.0,1.0);
792 gener->SetJetEtaRange(-1.,1.);
793 gener->SetJetPhiRange(200.,340.);
794 gener->SetFragPhotonInPHOS(kTRUE);
795 gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
796 printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
798 case kPyJetJetPHOSv2:
799 comment = comment.Append(" pp->jet + jet over PHOS version2 ");
800 gener = new AliGenPythia(nEvts);
801 gener->SetEnergyCMS(eCMS);
802 gener->SetProcess(kPyJets);
803 gener->SetStrucFunc(kCTEQ4L);
804 gener->SetPtHard(ptHardMin,ptHardMax);
805 //gener->SetYHard(-1.0,1.0);
806 gener->SetJetEtaRange(-1.,1.);
807 gener->SetJetPhiRange(200.,340.);
808 //gener->SetPi0InPHOS(kTRUE);
809 gener->SetPhotonInPHOSeta(kTRUE);
810 gener->SetPhotonMinPt(ptGammaPi0Min);
811 gener->SetForceDecay(kAll);
813 case kPyGammaJetEMCAL:
814 comment = comment.Append(" pp->jet + gamma over EMCAL at 14 TeV");
815 gener = new AliGenPythia(nEvts);
816 gener->SetEnergyCMS(eCMS);
817 gener->SetProcess(kPyDirectGamma);
818 gener->SetStrucFunc(kCTEQ4L);
819 gener->SetPtHard(ptHardMin,ptHardMax);
820 //gener->SetYHard(-1.0,1.0);
821 gener->SetGammaEtaRange(-0.71,0.71);
822 gener->SetGammaPhiRange(78.,192.);//Over 6 supermodules +-2 deg
825 comment = comment.Append(" pp->jet + jet over EMCAL at 14 TeV");
826 gener = new AliGenPythia(nEvts);
827 gener->SetEnergyCMS(eCMS);
828 gener->SetProcess(kPyJets);
829 gener->SetStrucFunc(kCTEQ4L);
830 gener->SetPtHard(ptHardMin,ptHardMax);
831 //gener->SetYHard(-1.0,1.0);
832 gener->SetJetEtaRange(-1,1);
833 gener->SetJetPhiRange(60.,210.);
834 gener->SetPi0InEMCAL(kTRUE);
835 gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
836 printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
838 case kPyGammaBremsEMCAL:
839 comment = comment.Append(" pp->jet + jet+bremsphoton over EMCAL at 14 TeV");
840 gener = new AliGenPythia(nEvts);
841 gener->SetEnergyCMS(eCMS);
842 gener->SetProcess(kPyJets);
843 gener->SetStrucFunc(kCTEQ4L);
844 gener->SetPtHard(ptHardMin,ptHardMax);
845 //gener->SetYHard(-1.0,1.0);
846 gener->SetJetEtaRange(-1,1);
847 gener->SetJetPhiRange(60.,210.); //Over 2 uncovered PHOS modules
848 gener->SetFragPhotonInEMCAL(kTRUE);
849 gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
851 printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
859 AliGenerator* MbCocktail()
861 comment = comment.Append(" pp at 14 TeV: Pythia low-pt, no heavy quarks + J/Psi from parameterisation");
862 AliGenCocktail * gener = new AliGenCocktail();
863 gener->UsePerEventRates();
867 AliGenPythia* pythia = new AliGenPythia(-1);
868 pythia->SetMomentumRange(0, 999999.);
869 pythia->SetThetaRange(0., 180.);
870 pythia->SetYRange(-12.,12.);
871 pythia->SetPtRange(0,1000.);
872 pythia->SetProcess(kPyMb);
873 pythia->SetEnergyCMS(14000.);
874 pythia->SwitchHFOff();
876 // J/Psi parameterisation
878 AliGenParam* jpsi = new AliGenParam(1, AliGenMUONlib::kJpsi, "CDF scaled", "Jpsi");
879 jpsi->SetPtRange(0.,100.);
880 jpsi->SetYRange(-8., 8.);
881 jpsi->SetPhiRange(0., 360.);
882 jpsi->SetForceDecay(kAll);
884 gener->AddGenerator(pythia, "Pythia", 1.);
885 gener->AddGenerator(jpsi, "J/Psi", 8.e-4);
890 AliGenerator* PyMbTriggered(Int_t pdg)
892 AliGenPythia* pythia = new AliGenPythia(-1);
893 pythia->SetMomentumRange(0, 999999.);
894 pythia->SetThetaRange(0., 180.);
895 pythia->SetYRange(-12.,12.);
896 pythia->SetPtRange(0,1000.);
897 pythia->SetProcess(kPyMb);
898 pythia->SetEnergyCMS(14000.);
899 pythia->SetTriggerParticle(pdg, 0.9);
903 void ProcessEnvironmentVars()
905 cout << "######################################" << endl;
906 cout << "## Processing environment variables ##" << endl;
907 cout << "######################################" << endl;
910 if (gSystem->Getenv("DC_RUN")) {
911 runNumber = atoi(gSystem->Getenv("DC_RUN"));
913 //cout<<"Run number "<<runNumber<<endl;
916 if (gSystem->Getenv("DC_EVENT")) {
917 eventNumber = atoi(gSystem->Getenv("DC_EVENT"));
919 //cout<<"Event number "<<eventNumber<<endl;
921 // Random Number seed
922 if (gSystem->Getenv("CONFIG_SEED")) {
923 seed = atoi(gSystem->Getenv("CONFIG_SEED"));
925 else if(gSystem->Getenv("DC_EVENT") && gSystem->Getenv("DC_RUN")){
926 seed = runNumber * 100000 + eventNumber;
929 gRandom->SetSeed(seed);
931 cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
932 cout<<"Seed for random number generation= "<< seed<<" "<< gRandom->GetSeed()<<endl;
933 cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
936 if (gSystem->Getenv("DC_RUN_TYPE")) {
937 cout<<"run type "<<gSystem->Getenv("DC_RUN_TYPE")<<endl;
938 for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
939 if (strcmp(gSystem->Getenv("DC_RUN_TYPE"), pprRunName[iRun])==0) {
940 proc = (PDC07Proc_t)iRun;
945 cout << "Environment variable DC_RUN_TYPE is not defined" << endl;
947 if (gSystem->Getenv("ECMS"))
948 eCMS = atof(gSystem->Getenv("ECMS"));
949 if (gSystem->Getenv("PTHARDMIN"))
950 ptHardMin = atof(gSystem->Getenv("PTHARDMIN"));
951 if (gSystem->Getenv("PTHARDMAX"))
952 ptHardMax = atof(gSystem->Getenv("PTHARDMAX"));
953 if (gSystem->Getenv("PTGAMMAPI0MIN"))
954 ptGammaPi0Min = atof(gSystem->Getenv("PTGAMMAPI0MIN"));
955 if (gSystem->Getenv("QUENCHING"))
956 iquenching = atof(gSystem->Getenv("QUENCHING"));
957 if (gSystem->Getenv("QHAT"))
958 qhat = atof(gSystem->Getenv("QHAT"));
960 cout<<">> Run type set to "<<pprRunName[proc]<<endl;
961 cout<<">> Collision energy set to "<<eCMS <<endl;
962 cout<<">> ptHard limits: "<<ptHardMin<<" to " <<ptHardMax<<" GeV"<<endl;
963 cout<<">> pt gamma/pi0 threshold "<< ptGammaPi0Min<<" GeV "<<endl;
964 cout<<">> quenching on? "<< iquenching<<" qhat "<<qhat<<endl;
966 cout << "######################################" << endl;
967 cout << "######################################" << endl;