]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/QA/Config.C
typo
[u/mrichter/AliRoot.git] / test / QA / Config.C
1 //
2 // Configuration for the Physics Data Challenge 2007
3 //
4 //
5 // Assuming inel = 70 mb (PPR v.1, p.64)
6 //
7 // 84.98% of MSEL=0 events (including diffractive) with
8 // QQbar switched off (these events will include injected J/psi) => kPyMbNoHvq
9 //
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
15 //
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
21 //
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
24
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++")
30
31 #if !defined(__CINT__) || defined(__MAKECINT__)
32 #include <Riostream.h>
33 #include <TRandom.h>
34 #include <TDatime.h>
35 #include <TSystem.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"
70 #endif
71
72
73 enum PDC07Proc_t
74 {
75 //--- Heavy Flavour Production ---
76   kCharmPbPb5500,  kCharmpPb8800,  kCharmpp14000,  kCharmpp14000wmi,
77   kD0PbPb5500,     kD0pPb8800,     kD0pp14000,
78   kDPlusPbPb5500,  kDPluspPb8800,  kDPluspp14000,
79   kBeautyPbPb5500, kBeautypPb8800, kBeautypp14000, kBeautypp14000wmi,
80 // -- Pythia Mb
81   kPyMbNoHvq, kPyOmegaPlus, kPyOmegaMinus, kPyJetJet,
82   kPyGammaJetPHOS, kPyJetJetPHOS, kPyJetJetPHOSv2, kPyGammaBremsPHOS,
83  kPyGammaJetEMCAL, kPyJetJetEMCAL, kPyGammaBremsEMCAL, kRunMax
84 };
85
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"
94 };
95
96
97 //--- Decay Mode ---
98 enum DecayHvFl_t
99 {
100   kNature,  kHadr, kSemiEl, kSemiMu
101 };
102 //--- Magnetic Field ---
103 enum Mag_t
104 {
105     k2kG, k4kG, k5kG
106 };
107
108 //--- Trigger config ---
109 enum TrigConf_t
110 {
111     kDefaultPPTrig, kDefaultPbPbTrig
112 };
113
114 const char * TrigConfName[] = {
115     "p-p","Pb-Pb"
116 };
117
118 Float_t eCMS=5500;
119 PDC07Proc_t proc = kPyJetJetEMCAL;
120
121 enum  PprGeo_t
122   {
123     kHoles, kNoHoles
124   };
125
126 static PprGeo_t geo = kHoles;
127 class AliGenPythia ;
128 //--- Functions ---
129 AliGenPythia *PythiaHard(PDC07Proc_t );
130 AliGenPythia *PythiaHVQ(PDC07Proc_t );
131 AliGenerator *MbCocktail();
132 AliGenerator *PyMbTriggered(Int_t pdg);
133 void ProcessEnvironmentVars();
134
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 //========================//
144 TDatime dt;
145 static UInt_t seed    = dt.Get();
146
147 // nEvts = -1  : you get 1 QQbar pair and all the fragmentation and
148 //               decay chain
149 // nEvts = N>0 : you get N charm / beauty Hadrons
150 Int_t nEvts = -1;
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;
154
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
160 // normalizations:
161 //  CHARM
162 // 2.76-3 GeV: 25%
163 //    3-4 GeV: 40%
164 //    4-8 GeV: 29%
165 //     >8 GeV:  6%
166 //  BEAUTY
167 // 2.76-4 GeV:  5%
168 //    4-6 GeV: 31%
169 //    6-8 GeV: 28%
170 //     >8 GeV: 36%
171
172 Float_t ptHardMin =  10.;
173 Float_t ptHardMax =  20.;
174 Float_t ptGammaPi0Min = 1.;
175 Int_t iquenching = 0;
176 Float_t qhat = 20.;
177
178 // Comment line
179 static TString comment;
180
181 void Config()
182 {
183
184
185   // Get settings from environment variables
186   ProcessEnvironmentVars();
187
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
196 #endif
197   new TGeant3TGeo("C++ Interface to Geant3");
198
199   // Output every 100 tracks
200   ((TGeant3*)gMC)->SetSWIT(4,100);
201
202   //=======================================================================
203   // Run loader
204   AliRunLoader* rl=0x0;
205   rl = AliRunLoader::Open("galice.root",
206                           AliConfig::GetDefaultEventFolderName(),
207                           "recreate");
208   if (rl == 0x0)
209     {
210       gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
211       return;
212     }
213   rl->SetCompressionLevel(2);
214   rl->SetNumberOfEventsPerFile(1000);
215   gAlice->SetRunLoader(rl);
216
217   // Run number
218   //gAlice->SetRunNumber(runNumber);
219
220   // Set the trigger configuration
221   gAlice->SetTriggerDescriptor(TrigConfName[trig]);
222   cout<<"Trigger configuration is set to  "<<TrigConfName[trig]<<endl;
223
224   //
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
229
230
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);
245
246     Float_t cut = 1.e-3;        // 1MeV cut by default
247     Float_t tofmax = 1.e10;
248
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);
260
261   //======================//
262   // Set External decayer //
263   //======================//
264   TVirtualMCDecayer* decayer = new AliDecayerPythia();
265   // DECAYS
266   //
267   switch(decHvFl) {
268   case kNature:
269     decayer->SetForceDecay(kAll);
270     break;
271   case kHadr:
272     decayer->SetForceDecay(kHadronicD);
273     break;
274   case kSemiEl:
275     decayer->SetForceDecay(kSemiElectronic);
276     break;
277   case kSemiMu:
278     decayer->SetForceDecay(kSemiMuonic);
279     break;
280   }
281   decayer->Init();
282   gMC->SetExternalDecayer(decayer);
283   if(proc == kPyJetJetPHOSv2) // in this case we decay the pi0
284     decayer->SetForceDecay(kNeutralPion);
285
286   //=========================//
287   // Generator Configuration //
288   //=========================//
289   AliGenPythia* gener = 0x0;
290
291   if (proc <=   kBeautypp14000wmi) {
292     AliGenPythia *pythia = PythiaHVQ(proc);
293     // FeedDown option
294     pythia->SetFeedDownHigherFamily(kFALSE);
295     // Stack filling option
296     if(!stars) pythia->SetStackFillOpt(AliGenPythia::kParentSelection);
297     // Set Count mode
298     if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
299     //
300     // DECAYS
301     //
302     switch(decHvFl) {
303     case kNature:
304       pythia->SetForceDecay(kAll);
305       break;
306     case kHadr:
307       pythia->SetForceDecay(kHadronicD);
308       break;
309     case kSemiEl:
310       pythia->SetForceDecay(kSemiElectronic);
311       break;
312     case kSemiMu:
313       pythia->SetForceDecay(kSemiMuonic);
314       break;
315     }
316     //
317     // GEOM & KINE CUTS
318     //
319     pythia->SetMomentumRange(0,99999999);
320     pythia->SetPhiRange(0., 360.);
321     pythia->SetThetaRange(0,180);
322     switch(ycut) {
323     case kFull:
324       pythia->SetYRange(-999,999);
325       break;
326     case kBarrel:
327       pythia->SetYRange(-2,2);
328       break;
329     case kMuonArm:
330       pythia->SetYRange(1,6);
331       break;
332     }
333     gener = pythia;
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);
342     // FeedDown option
343     pythia->SetFeedDownHigherFamily(kFALSE);
344     // Set Count mode
345     if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
346
347       //
348       // GEOM & KINE CUTS
349       //
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);
354
355    if(proc == kPyJetJetPHOSv2)
356          pythia->SetForceDecay(kNeutralPion);
357    else
358          pythia->SetForceDecay(kAll);
359
360     pythia->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
361     pythia->SetPtKick(5); // set the intrinsic kt to 5 GeV/c
362     gener = pythia;
363   }
364
365
366   // PRIMARY VERTEX
367
368   gener->SetOrigin(0., 0., 0.);    // vertex position
369
370   // Size of the interaction diamond
371   // Longitudinal
372   Float_t sigmaz  = 7.55 / TMath::Sqrt(2.); // [cm]
373
374   // Transverse
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);
380
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);
384
385   gener->Init();
386
387   //Quenching
388   gener->SetQuench(iquenching);
389   if(iquenching == 1){
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);
392
393   }
394   // FIELD
395
396   if (mag == k2kG) {
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");
402   }
403   printf("\n \n Comment: %s \n \n", comment.Data());
404
405   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag);
406   field->SetL3ConstField(0); //Using const. field in the barrel
407   rl->CdGAFile();
408   gAlice->SetField(field);
409
410
411
412   Int_t iABSO  = 1;
413   Int_t iACORDE   = 0;
414   Int_t iDIPO  = 1;
415   Int_t iEMCAL = 1;
416   Int_t iFMD   = 1;
417   Int_t iFRAME = 1;
418   Int_t iHALL  = 1;
419   Int_t iITS   = 1;
420   Int_t iMAG   = 1;
421   Int_t iMUON  = 1;
422   Int_t iPHOS  = 1;
423   Int_t iPIPE  = 1;
424   Int_t iPMD   = 1;
425   Int_t iHMPID = 1;
426   Int_t iSHIL  = 1;
427   Int_t iT0    = 1;
428   Int_t iTOF   = 1;
429   Int_t iTPC   = 1;
430   Int_t iTRD   = 1;
431   Int_t iVZERO = 1;
432   Int_t iZDC   = 1;
433
434
435     //=================== Alice BODY parameters =============================
436     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
437
438
439     if (iMAG)
440     {
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");
445     }
446
447
448     if (iABSO)
449     {
450         //=================== ABSO parameters ============================
451         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
452     }
453
454     if (iDIPO)
455     {
456         //=================== DIPO parameters ============================
457
458         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
459     }
460
461     if (iHALL)
462     {
463         //=================== HALL parameters ============================
464
465         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
466     }
467
468
469     if (iFRAME)
470     {
471         //=================== FRAME parameters ============================
472
473         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
474         if (geo == kHoles) FRAME->SetHoles(1);
475         else FRAME->SetHoles(0);
476
477     }
478
479     if (iSHIL)
480     {
481         //=================== SHIL parameters ============================
482
483         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
484     }
485
486
487     if (iPIPE)
488     {
489         //=================== PIPE parameters ============================
490
491         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
492     }
493
494     if(iITS) {
495
496     //=================== ITS parameters ============================
497     //
498         AliITSvPPRasymmFMD *ITS  = new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version with asymmetric services");
499     }
500
501     if (iTPC)
502     {
503       //============================ TPC parameters =====================
504         AliTPC *TPC = new AliTPCv2("TPC", "Default");
505     }
506
507
508     if (iTOF) {
509         //=================== TOF parameters ============================
510         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
511    }
512
513
514     if (iHMPID)
515     {
516         //=================== HMPID parameters ===========================
517         AliHMPID *HMPID = new AliHMPIDv2("HMPID", "normal HMPID");
518
519     }
520
521
522     if (iZDC)
523     {
524         //=================== ZDC parameters ============================
525
526         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
527     }
528
529     if (iTRD)
530     {
531         //=================== TRD parameters ============================
532
533         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
534         AliTRDgeometry *geoTRD = TRD->GetGeometry();
535         // Partial geometry: modules at 2,3,4,6,11,12,14,15
536         // starting at 6h in positive direction
537         //Hole on SM 13 and 14 for PHOS
538         geoTRD->SetSMstatus(13,0);
539         geoTRD->SetSMstatus(14,0);
540         geoTRD->SetSMstatus(15,0);
541     }
542
543     if (iFMD)
544     {
545         //=================== FMD parameters ============================
546         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
547    }
548
549     if (iMUON)
550     {
551         //=================== MUON parameters ===========================
552         // New MUONv1 version (geometry defined via builders)
553         AliMUON *MUON = new AliMUONv1("MUON", "default");
554     }
555     //=================== PHOS parameters ===========================
556
557     if (iPHOS)
558     {
559         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
560     }
561
562
563     if (iPMD)
564     {
565         //=================== PMD parameters ============================
566         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
567     }
568
569     if (iT0)
570     {
571         //=================== T0 parameters ============================
572         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
573     }
574
575     if (iEMCAL)
576     {
577         //=================== EMCAL parameters ============================
578         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "ESHISH_77_TRD1_2X2_FINAL_110DEG scTh=0.176 pbTh=0.144");
579     }
580
581      if (iACORDE)
582     {
583         //=================== CRT parameters ============================
584         AliACORDE *ACORDE = new AliACORDEv0("CRT", "normal ACORDE");
585     }
586
587      if (iVZERO)
588     {
589         //=================== CRT parameters ============================
590         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
591     }
592 }
593
594 //           PYTHIA
595
596 AliGenPythia *PythiaHVQ(PDC07Proc_t proc) {
597 //*******************************************************************//
598 // Configuration file for charm / beauty generation with PYTHIA      //
599 //                                                                   //
600 // The parameters have been tuned in order to reproduce the inclusive//
601 // heavy quark pt distribution given by the NLO pQCD calculation by  //
602 // Mangano, Nason and Ridolfi.                                       //
603 //                                                                   //
604 // For details and for the NORMALIZATION of the yields see:          //
605 //   N.Carrer and A.Dainese,                                         //
606 //   "Charm and beauty production at the LHC",                       //
607 //   ALICE-INT-2003-019, [arXiv:hep-ph/0311225];                     //
608 //   PPR Chapter 6.6, CERN/LHCC 2005-030 (2005).                     //
609 //*******************************************************************//
610   AliGenPythia * gener = 0x0;
611
612   switch(proc) {
613   case kCharmPbPb5500:
614       comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
615       gener = new AliGenPythia(nEvts);
616       gener->SetProcess(kPyCharmPbPbMNR);
617       gener->SetStrucFunc(kCTEQ4L);
618       gener->SetPtHard(2.1,-1.0);
619       gener->SetEnergyCMS(5500.);
620       gener->SetNuclei(208,208);
621       break;
622   case kCharmpPb8800:
623       comment = comment.Append(" Charm in p-Pb at 8.8 TeV");
624       gener = new AliGenPythia(nEvts);
625       gener->SetProcess(kPyCharmpPbMNR);
626       gener->SetStrucFunc(kCTEQ4L);
627       gener->SetPtHard(2.1,-1.0);
628       gener->SetEnergyCMS(8800.);
629       gener->SetProjectile("P",1,1);
630       gener->SetTarget("Pb",208,82);
631       break;
632   case kCharmpp14000:
633       comment = comment.Append(" Charm in pp at 14 TeV");
634       gener = new AliGenPythia(nEvts);
635       gener->SetProcess(kPyCharmppMNR);
636       gener->SetStrucFunc(kCTEQ4L);
637       gener->SetPtHard(2.1,-1.0);
638       gener->SetEnergyCMS(14000.);
639       break;
640   case kCharmpp14000wmi:
641       comment = comment.Append(" Charm in pp at 14 TeV with mult. interactions");
642       gener = new AliGenPythia(-1);
643       gener->SetProcess(kPyCharmppMNRwmi);
644       gener->SetStrucFunc(kCTEQ5L);
645       gener->SetPtHard(ptHardMin,ptHardMax);
646       gener->SetEnergyCMS(14000.);
647       break;
648   case kD0PbPb5500:
649       comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
650       gener = new AliGenPythia(nEvts);
651       gener->SetProcess(kPyD0PbPbMNR);
652       gener->SetStrucFunc(kCTEQ4L);
653       gener->SetPtHard(2.1,-1.0);
654       gener->SetEnergyCMS(5500.);
655       gener->SetNuclei(208,208);
656       break;
657   case kD0pPb8800:
658       comment = comment.Append(" D0 in p-Pb at 8.8 TeV");
659       gener = new AliGenPythia(nEvts);
660       gener->SetProcess(kPyD0pPbMNR);
661       gener->SetStrucFunc(kCTEQ4L);
662       gener->SetPtHard(2.1,-1.0);
663       gener->SetEnergyCMS(8800.);
664       gener->SetProjectile("P",1,1);
665       gener->SetTarget("Pb",208,82);
666       break;
667   case kD0pp14000:
668       comment = comment.Append(" D0 in pp at 14 TeV");
669       gener = new AliGenPythia(nEvts);
670       gener->SetProcess(kPyD0ppMNR);
671       gener->SetStrucFunc(kCTEQ4L);
672       gener->SetPtHard(2.1,-1.0);
673       gener->SetEnergyCMS(14000.);
674       break;
675   case kDPlusPbPb5500:
676       comment = comment.Append(" DPlus in Pb-Pb at 5.5 TeV");
677       gener = new AliGenPythia(nEvts);
678       gener->SetProcess(kPyDPlusPbPbMNR);
679       gener->SetStrucFunc(kCTEQ4L);
680       gener->SetPtHard(2.1,-1.0);
681       gener->SetEnergyCMS(5500.);
682       gener->SetNuclei(208,208);
683       break;
684   case kDPluspPb8800:
685       comment = comment.Append(" DPlus in p-Pb at 8.8 TeV");
686       gener = new AliGenPythia(nEvts);
687       gener->SetProcess(kPyDPluspPbMNR);
688       gener->SetStrucFunc(kCTEQ4L);
689       gener->SetPtHard(2.1,-1.0);
690       gener->SetEnergyCMS(8800.);
691       gener->SetProjectile("P",1,1);
692       gener->SetTarget("Pb",208,82);
693       break;
694   case kDPluspp14000:
695       comment = comment.Append(" DPlus in pp at 14 TeV");
696       gener = new AliGenPythia(nEvts);
697       gener->SetProcess(kPyDPlusppMNR);
698       gener->SetStrucFunc(kCTEQ4L);
699       gener->SetPtHard(2.1,-1.0);
700       gener->SetEnergyCMS(14000.);
701       break;
702   case kBeautyPbPb5500:
703       comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
704       gener = new AliGenPythia(nEvts);
705       gener->SetProcess(kPyBeautyPbPbMNR);
706       gener->SetStrucFunc(kCTEQ4L);
707       gener->SetPtHard(2.75,-1.0);
708       gener->SetEnergyCMS(5500.);
709       gener->SetNuclei(208,208);
710       break;
711   case kBeautypPb8800:
712       comment = comment.Append(" Beauty in p-Pb at 8.8 TeV");
713       gener = new AliGenPythia(nEvts);
714       gener->SetProcess(kPyBeautypPbMNR);
715       gener->SetStrucFunc(kCTEQ4L);
716       gener->SetPtHard(2.75,-1.0);
717       gener->SetEnergyCMS(8800.);
718       gener->SetProjectile("P",1,1);
719       gener->SetTarget("Pb",208,82);
720       break;
721   case kBeautypp14000:
722       comment = comment.Append(" Beauty in pp at 14 TeV");
723       gener = new AliGenPythia(nEvts);
724       gener->SetProcess(kPyBeautyppMNR);
725       gener->SetStrucFunc(kCTEQ4L);
726       gener->SetPtHard(2.75,-1.0);
727       gener->SetEnergyCMS(14000.);
728       break;
729   case kBeautypp14000wmi:
730       comment = comment.Append(" Beauty in pp at 14 TeV with mult. interactions");
731       gener = new AliGenPythia(-1);
732       gener->SetProcess(kPyBeautyppMNRwmi);
733       gener->SetStrucFunc(kCTEQ5L);
734       gener->SetPtHard(ptHardMin,ptHardMax);
735       gener->SetEnergyCMS(14000.);
736       break;
737   }
738
739   return gener;
740 }
741
742 AliGenPythia *PythiaHard(PDC07Proc_t proc) {
743 //*******************************************************************//
744 // Configuration file for hard QCD processes generation with PYTHIA  //
745 //                                                                   //
746 //*******************************************************************//
747   AliGenPythia * gener = 0x0;
748
749   switch(proc) {
750
751   case kPyJetJet:
752     comment = comment.Append(" pp->jet + jet over at 14 TeV, no restriction");
753     AliGenPythia * gener = new AliGenPythia(nEvts);
754     gener->SetEnergyCMS(eCMS);//        Centre of mass energy
755     gener->SetProcess(kPyJets);//        Process type
756     gener->SetJetEtaRange(-1.5, 1.5);//  Final state kinematic cuts
757     gener->SetJetPhiRange(0., 360.);
758     gener->SetJetEtRange(10., 1000.);
759     gener->SetPtHard(ptHardMin, ptHardMax);// Pt transfer of the hard scattering
760     gener->SetStrucFunc(kCTEQ4L);
761
762     return gener;
763
764   case kPyGammaJetPHOS:
765       comment = comment.Append(" pp->jet + gamma over PHOS");
766       gener = new AliGenPythia(nEvts);
767       gener->SetEnergyCMS(eCMS);
768       gener->SetProcess(kPyDirectGamma);
769       gener->SetStrucFunc(kCTEQ4L);
770       gener->SetPtHard(ptHardMin,ptHardMax);
771       //gener->SetYHard(-1.0,1.0);
772       gener->SetGammaEtaRange(-0.13,0.13);
773       gener->SetGammaPhiRange(218.,322.);//Over 5 modules +-2 deg
774       break;
775   case kPyJetJetPHOS:
776       comment = comment.Append(" pp->jet + jet over PHOS");
777       gener = new AliGenPythia(nEvts);
778       gener->SetEnergyCMS(eCMS);
779       gener->SetProcess(kPyJets);
780       gener->SetStrucFunc(kCTEQ4L);
781       gener->SetPtHard(ptHardMin,ptHardMax);
782       //gener->SetYHard(-1.0,1.0);
783       gener->SetJetEtaRange(-1.,1.);
784       gener->SetJetPhiRange(200.,340.);
785       gener->SetPi0InPHOS(kTRUE);
786       gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
787
788       printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
789       break;
790   case kPyGammaBremsPHOS:
791       comment = comment.Append(" pp->jet + jet+bremsphoton over PHOS at 14 TeV");
792       gener = new AliGenPythia(nEvts);
793       gener->SetEnergyCMS(eCMS);
794       gener->SetProcess(kPyJets);
795       gener->SetStrucFunc(kCTEQ4L);
796       gener->SetPtHard(ptHardMin,ptHardMax);
797       //gener->SetYHard(-1.0,1.0);
798       gener->SetJetEtaRange(-1.,1.);
799       gener->SetJetPhiRange(200.,340.);
800       gener->SetFragPhotonInPHOS(kTRUE);
801       gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
802           printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
803           break;
804         case kPyJetJetPHOSv2:
805         comment = comment.Append(" pp->jet + jet over PHOS version2 ");
806         gener = new AliGenPythia(nEvts);
807         gener->SetEnergyCMS(eCMS);
808         gener->SetProcess(kPyJets);
809         gener->SetStrucFunc(kCTEQ4L);
810         gener->SetPtHard(ptHardMin,ptHardMax);
811         //gener->SetYHard(-1.0,1.0);
812         gener->SetJetEtaRange(-1.,1.);
813         gener->SetJetPhiRange(200.,340.);
814         //gener->SetPi0InPHOS(kTRUE);
815         gener->SetPhotonInPHOSeta(kTRUE);
816         gener->SetPhotonMinPt(ptGammaPi0Min);
817         gener->SetForceDecay(kAll);
818         break;
819   case kPyGammaJetEMCAL:
820     comment = comment.Append(" pp->jet + gamma over EMCAL at 14 TeV");
821     gener = new AliGenPythia(nEvts);
822     gener->SetEnergyCMS(eCMS);
823     gener->SetProcess(kPyDirectGamma);
824     gener->SetStrucFunc(kCTEQ4L);
825     gener->SetPtHard(ptHardMin,ptHardMax);
826     //gener->SetYHard(-1.0,1.0);
827     gener->SetGammaEtaRange(-0.71,0.71);
828     gener->SetGammaPhiRange(78.,192.);//Over 6 supermodules +-2 deg
829     break;
830   case kPyJetJetEMCAL:
831     comment = comment.Append(" pp->jet + jet over EMCAL at 14 TeV");
832     gener = new AliGenPythia(nEvts);
833     gener->SetEnergyCMS(eCMS);
834     gener->SetProcess(kPyJets);
835     gener->SetStrucFunc(kCTEQ4L);
836     gener->SetPtHard(ptHardMin,ptHardMax);
837     //gener->SetYHard(-1.0,1.0);
838     gener->SetJetEtaRange(-1,1);
839     gener->SetJetPhiRange(60.,210.);
840     gener->SetPi0InEMCAL(kTRUE);
841     gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
842     printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
843     break;
844   case kPyGammaBremsEMCAL:
845       comment = comment.Append(" pp->jet + jet+bremsphoton over EMCAL at 14 TeV");
846       gener = new AliGenPythia(nEvts);
847       gener->SetEnergyCMS(eCMS);
848       gener->SetProcess(kPyJets);
849       gener->SetStrucFunc(kCTEQ4L);
850       gener->SetPtHard(ptHardMin,ptHardMax);
851       //gener->SetYHard(-1.0,1.0);
852       gener->SetJetEtaRange(-1,1);
853       gener->SetJetPhiRange(60.,210.); //Over 2 uncovered PHOS modules
854       gener->SetFragPhotonInEMCAL(kTRUE);
855       gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
856
857       printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
858       break;
859
860   }
861
862   return gener;
863 }
864
865 AliGenerator* MbCocktail()
866 {
867       comment = comment.Append(" pp at 14 TeV: Pythia low-pt, no heavy quarks + J/Psi from parameterisation");
868       AliGenCocktail * gener = new AliGenCocktail();
869       gener->UsePerEventRates();
870
871 //    Pythia
872
873       AliGenPythia* pythia = new AliGenPythia(-1);
874       pythia->SetMomentumRange(0, 999999.);
875       pythia->SetThetaRange(0., 180.);
876       pythia->SetYRange(-12.,12.);
877       pythia->SetPtRange(0,1000.);
878       pythia->SetProcess(kPyMb);
879       pythia->SetEnergyCMS(14000.);
880       pythia->SwitchHFOff();
881
882 //   J/Psi parameterisation
883
884       AliGenParam* jpsi = new AliGenParam(1, AliGenMUONlib::kJpsi, "CDF scaled", "Jpsi");
885       jpsi->SetPtRange(0.,100.);
886       jpsi->SetYRange(-8., 8.);
887       jpsi->SetPhiRange(0., 360.);
888       jpsi->SetForceDecay(kAll);
889
890       gener->AddGenerator(pythia, "Pythia", 1.);
891       gener->AddGenerator(jpsi,   "J/Psi", 8.e-4);
892
893       return gener;
894 }
895
896 AliGenerator* PyMbTriggered(Int_t pdg)
897 {
898     AliGenPythia* pythia = new AliGenPythia(-1);
899     pythia->SetMomentumRange(0, 999999.);
900     pythia->SetThetaRange(0., 180.);
901     pythia->SetYRange(-12.,12.);
902     pythia->SetPtRange(0,1000.);
903     pythia->SetProcess(kPyMb);
904     pythia->SetEnergyCMS(14000.);
905     pythia->SetTriggerParticle(pdg, 0.9);
906     return pythia;
907 }
908
909 void ProcessEnvironmentVars()
910 {
911   cout << "######################################" << endl;
912   cout << "## Processing environment variables ##" << endl;
913   cout << "######################################" << endl;
914
915     // Run Number
916     if (gSystem->Getenv("DC_RUN")) {
917       runNumber = atoi(gSystem->Getenv("DC_RUN"));
918     }
919     //cout<<"Run number "<<runNumber<<endl;
920
921     // Event Number
922     if (gSystem->Getenv("DC_EVENT")) {
923       eventNumber = atoi(gSystem->Getenv("DC_EVENT"));
924     }
925     //cout<<"Event number "<<eventNumber<<endl;
926
927    // Random Number seed
928     if (gSystem->Getenv("CONFIG_SEED")) {
929       seed = atoi(gSystem->Getenv("CONFIG_SEED"));
930     }
931     else if(gSystem->Getenv("DC_EVENT") && gSystem->Getenv("DC_RUN")){
932       seed = runNumber * 100000 + eventNumber;
933     }
934
935     gRandom->SetSeed(seed);
936
937     cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
938     cout<<"Seed for random number generation= "<< seed<<"  "<< gRandom->GetSeed()<<endl;
939     cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
940
941     // Run type
942     if (gSystem->Getenv("DC_RUN_TYPE")) {
943       cout<<"run type "<<gSystem->Getenv("DC_RUN_TYPE")<<endl;
944       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
945         if (strcmp(gSystem->Getenv("DC_RUN_TYPE"), pprRunName[iRun])==0) {
946           proc = (PDC07Proc_t)iRun;
947         }
948       }
949     }
950     else
951       cout << "Environment variable DC_RUN_TYPE is not defined" << endl;
952
953     if (gSystem->Getenv("ECMS"))
954       eCMS = atof(gSystem->Getenv("ECMS"));
955     if (gSystem->Getenv("PTHARDMIN"))
956       ptHardMin = atof(gSystem->Getenv("PTHARDMIN"));
957     if (gSystem->Getenv("PTHARDMAX"))
958       ptHardMax = atof(gSystem->Getenv("PTHARDMAX"));
959     if (gSystem->Getenv("PTGAMMAPI0MIN"))
960       ptGammaPi0Min = atof(gSystem->Getenv("PTGAMMAPI0MIN"));
961    if (gSystem->Getenv("QUENCHING"))
962       iquenching = atof(gSystem->Getenv("QUENCHING"));
963    if (gSystem->Getenv("QHAT"))
964       qhat = atof(gSystem->Getenv("QHAT"));
965
966     cout<<">> Run type set to "<<pprRunName[proc]<<endl;
967     cout<<">> Collision energy set to "<<eCMS <<endl;
968     cout<<">> ptHard limits: "<<ptHardMin<<" to " <<ptHardMax<<" GeV"<<endl;
969     cout<<">> pt gamma/pi0 threshold "<< ptGammaPi0Min<<" GeV "<<endl;
970     cout<<">> quenching on? "<< iquenching<<" qhat "<<qhat<<endl;
971
972     cout << "######################################" << endl;
973     cout << "######################################" << endl;
974 }