]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/QA/Config.C
Removing obsolete macros, changes to work with the trunk, clean-up
[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/AliMagF.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"
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 //--- Trigger config ---
103 enum TrigConf_t
104 {
105     kDefaultPPTrig, kDefaultPbPbTrig
106 };
107
108 const char * TrigConfName[] = {
109     "p-p","Pb-Pb"
110 };
111
112 Float_t eCMS=5500;
113 PDC07Proc_t proc = kPyJetJetEMCAL;
114
115 enum  PprGeo_t
116   {
117     kHoles, kNoHoles
118   };
119
120 static PprGeo_t geo = kHoles;
121 class AliGenPythia ;
122 //--- Functions ---
123 AliGenPythia *PythiaHard(PDC07Proc_t );
124 AliGenPythia *PythiaHVQ(PDC07Proc_t );
125 AliGenerator *MbCocktail();
126 AliGenerator *PyMbTriggered(Int_t pdg);
127 void ProcessEnvironmentVars();
128
129 // This part for configuration
130 static DecayHvFl_t   decHvFl  = kNature;
131 static AliMagF::BMap_t mag    = AliMagF::k5kG;
132 static TrigConf_t    trig     = kDefaultPPTrig; // default pp trigger configuration
133 static Int_t         runNumber= 0;
134 static Int_t         eventNumber= 0;
135 //========================//
136 // Set Random Number seed //
137 //========================//
138 TDatime dt;
139 static UInt_t seed    = dt.Get();
140
141 // nEvts = -1  : you get 1 QQbar pair and all the fragmentation and
142 //               decay chain
143 // nEvts = N>0 : you get N charm / beauty Hadrons
144 Int_t nEvts = -1;
145 // stars = kTRUE : all heavy resonances and their decay stored
146 //       = kFALSE: only final heavy hadrons and their decays stored
147 Bool_t stars = kTRUE;
148
149 // To be used only with kCharmpp14000wmi and kBeautypp14000wmi
150 // To get a "reasonable" agreement with MNR results, events have to be
151 // generated with the minimum ptHard set to 2.76 GeV.
152 // To get a "perfect" agreement with MNR results, events have to be
153 // generated in four ptHard bins with the following relative
154 // normalizations:
155 //  CHARM
156 // 2.76-3 GeV: 25%
157 //    3-4 GeV: 40%
158 //    4-8 GeV: 29%
159 //     >8 GeV:  6%
160 //  BEAUTY
161 // 2.76-4 GeV:  5%
162 //    4-6 GeV: 31%
163 //    6-8 GeV: 28%
164 //     >8 GeV: 36%
165
166 Float_t ptHardMin =  10.;
167 Float_t ptHardMax =  20.;
168 Float_t ptGammaPi0Min = 1.;
169 Int_t iquenching = 0;
170 Float_t qhat = 20.;
171
172 // Comment line
173 static TString comment;
174
175 void Config()
176 {
177         // Get settings from environment variables
178   ProcessEnvironmentVars();
179
180   // Random Number seed
181   if (gSystem->Getenv("CONFIG_SEED")) {
182     seed = atoi(gSystem->Getenv("CONFIG_SEED"));
183   }
184   else if(gSystem->Getenv("DC_EVENT") && gSystem->Getenv("DC_RUN")){
185     seed = runNumber * 100000 + eventNumber;
186   }
187   gRandom->SetSeed(seed);
188   cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
189
190   // libraries required by geant321
191 #if defined(__CINT__)
192         gSystem->Load("liblhapdf");      // Parton density functions
193         gSystem->Load("libEGPythia6");   // TGenerator interface
194         gSystem->Load("libpythia6");     // Pythia
195         gSystem->Load("libAliPythia6");  // ALICE specific implementations
196         gSystem->Load("libgeant321");
197 #endif
198   new TGeant3TGeo("C++ Interface to Geant3");
199
200   // Output every 100 tracks
201   ((TGeant3*)gMC)->SetSWIT(4,100);
202
203   //=======================================================================
204   // Run loader
205   AliRunLoader* rl=0x0;
206
207   AliLog::Message(AliLog::kInfo, "Creating Run Loader", "", "", "Config()"," ConfigPPR.C", __LINE__);
208
209   rl = AliRunLoader::Open("galice.root",
210                           AliConfig::GetDefaultEventFolderName(),
211                           "recreate");
212   if (rl == 0x0)
213     {
214       gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
215       return;
216     }
217   rl->SetCompressionLevel(2);
218   rl->SetNumberOfEventsPerFile(100);
219   gAlice->SetRunLoader(rl);
220
221   // Set the trigger configuration
222   AliSimulation::Instance()->SetTriggerConfig(TrigConfName[trig]);
223   cout<<"Trigger configuration is set to  "<<TrigConfName[trig]<<endl;
224
225   //======================//
226   // Set External decayer //
227   //======================//
228   TVirtualMCDecayer* decayer = new AliDecayerPythia();
229   // DECAYS
230   //
231   switch(decHvFl) {
232     case kNature:
233       decayer->SetForceDecay(kAll);
234       break;
235     case kHadr:
236       decayer->SetForceDecay(kHadronicD);
237       break;
238     case kSemiEl:
239       decayer->SetForceDecay(kSemiElectronic);
240       break;
241     case kSemiMu:
242       decayer->SetForceDecay(kSemiMuonic);
243       break;
244   }
245   decayer->Init();
246   gMC->SetExternalDecayer(decayer);
247   if(proc == kPyJetJetPHOSv2) // in this case we decay the pi0
248     decayer->SetForceDecay(kNeutralPion);
249
250   //
251   //=======================================================================
252   // ************* STEERING parameters FOR ALICE SIMULATION **************
253   // --- Specify event type to be tracked through the ALICE setup
254   // --- All positions are in cm, angles in degrees, and P and E in GeV
255
256
257     gMC->SetProcess("DCAY",1);
258     gMC->SetProcess("PAIR",1);
259     gMC->SetProcess("COMP",1);
260     gMC->SetProcess("PHOT",1);
261     gMC->SetProcess("PFIS",0);
262     gMC->SetProcess("DRAY",0);
263     gMC->SetProcess("ANNI",1);
264     gMC->SetProcess("BREM",1);
265     gMC->SetProcess("MUNU",1);
266     gMC->SetProcess("CKOV",1);
267     gMC->SetProcess("HADR",1);
268     gMC->SetProcess("LOSS",2);
269     gMC->SetProcess("MULS",1);
270     gMC->SetProcess("RAYL",1);
271
272     Float_t cut = 1.e-3;        // 1MeV cut by default
273     Float_t tofmax = 1.e10;
274
275     gMC->SetCut("CUTGAM", cut);
276     gMC->SetCut("CUTELE", cut);
277     gMC->SetCut("CUTNEU", cut);
278     gMC->SetCut("CUTHAD", cut);
279     gMC->SetCut("CUTMUO", cut);
280     gMC->SetCut("BCUTE",  cut);
281     gMC->SetCut("BCUTM",  cut);
282     gMC->SetCut("DCUTE",  cut);
283     gMC->SetCut("DCUTM",  cut);
284     gMC->SetCut("PPCUTM", cut);
285     gMC->SetCut("TOFMAX", tofmax);
286
287
288   //=========================//
289   // Generator Configuration //
290   //=========================//
291   AliGenPythia* gener = 0x0;
292
293   if (proc <=   kBeautypp14000wmi) {
294     AliGenPythia *pythia = PythiaHVQ(proc);
295     // FeedDown option
296     pythia->SetFeedDownHigherFamily(kFALSE);
297     // Stack filling option
298     if(!stars) pythia->SetStackFillOpt(AliGenPythia::kParentSelection);
299     // Set Count mode
300     if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
301     //
302     // DECAYS
303     //
304     switch(decHvFl) {
305     case kNature:
306       pythia->SetForceDecay(kAll);
307       break;
308     case kHadr:
309       pythia->SetForceDecay(kHadronicD);
310       break;
311     case kSemiEl:
312       pythia->SetForceDecay(kSemiElectronic);
313       break;
314     case kSemiMu:
315       pythia->SetForceDecay(kSemiMuonic);
316       break;
317     }
318     //
319     // GEOM & KINE CUTS
320     //
321     pythia->SetMomentumRange(0,99999999);
322     pythia->SetPhiRange(0., 360.);
323     pythia->SetThetaRange(0,180);
324     switch(ycut) {
325     case kFull:
326       pythia->SetYRange(-999,999);
327       break;
328     case kBarrel:
329       pythia->SetYRange(-2,2);
330       break;
331     case kMuonArm:
332       pythia->SetYRange(1,6);
333       break;
334     }
335     gener = pythia;
336   } else if (proc == kPyMbNoHvq) {
337     gener = MbCocktail();
338   } else if (proc == kPyOmegaMinus) {
339     gener = PyMbTriggered(3334);
340   } else if (proc == kPyOmegaPlus) {
341     gener = PyMbTriggered(-3334);
342   } else if (proc <= kPyGammaBremsEMCAL) {
343     AliGenPythia *pythia = PythiaHard(proc);
344     // FeedDown option
345     pythia->SetFeedDownHigherFamily(kFALSE);
346     // Set Count mode
347     if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
348
349       //
350       // GEOM & KINE CUTS
351       //
352     pythia->SetMomentumRange(0,99999999);
353     //       pythia->SetJetEtaRange(-1.5, 1.5);//  Final state kinematic cuts
354     //       pythia->SetJetPhiRange(0., 360.);
355     //       pythia->SetThetaRange(45,135);
356
357    if(proc == kPyJetJetPHOSv2)
358                  pythia->SetForceDecay(kNeutralPion);
359          else
360                  pythia->SetForceDecay(kAll);
361                 
362                 pythia->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
363     pythia->SetPtKick(5); // set the intrinsic kt to 5 GeV/c
364     gener = pythia;
365         }
366
367
368   // PRIMARY VERTEX
369
370   gener->SetOrigin(0., 0., 0.);    // vertex position
371
372   // Size of the interaction diamond
373   // Longitudinal
374   Float_t sigmaz  = 7.55 / TMath::Sqrt(2.); // [cm]
375
376   // Transverse
377   Float_t betast  = 10;                 // beta* [m]
378   Float_t eps     = 3.75e-6;            // emittance [m]
379   Float_t gamma   = 7000. / 0.938272;   // relativistic gamma [1]
380   Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.;  // [cm]
381   printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
382
383   gener->SetSigma(sigmaxy, sigmaxy, sigmaz);      // Sigma in (X,Y,Z) (cm) on IP position
384   gener->SetCutVertexZ(3.);        // Truncate at 3 sigma
385   gener->SetVertexSmear(kPerEvent);
386
387   gener->Init();
388
389   //Quenching
390   gener->SetQuench(iquenching);
391   if(iquenching == 1){
392     Float_t k = 6e5*(qhat/1.7) ; //qhat=1.7, k = 6e5, default  value
393     AliPythia::Instance()->InitQuenching(0.,0.1,k,0,0.95,6);
394
395  }
396   // FIELD
397
398   if (mag == AliMagF::k2kG) {
399     comment = comment.Append(" | L3 field 0.2 T");
400   } else if (mag == AliMagF::k5kG) {
401     comment = comment.Append(" | L3 field 0.5 T");
402   }
403   printf("\n \n Comment: %s \n \n", comment.Data());
404   // to use constant field in the barrel use:
405   //  AliMagF* field = new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kGUniform);
406   AliMagF* field = new AliMagF("Maps","Maps", -1., -1., mag);
407   TGeoGlobalMagField::Instance()->SetField(field);
408   rl->CdGAFile();
409
410   Int_t iABSO  = 1;
411   Int_t iACORDE= 1;
412   Int_t iDIPO  = 1;
413   Int_t iEMCAL = 1;
414   Int_t iFMD   = 1;
415   Int_t iFRAME = 1;
416   Int_t iHALL  = 1;
417   Int_t iITS   = 1;
418   Int_t iMAG   = 1;
419   Int_t iMUON  = 1;
420   Int_t iPHOS  = 1;
421   Int_t iPIPE  = 1;
422   Int_t iPMD   = 1;
423   Int_t iHMPID = 1;
424   Int_t iSHIL  = 1;
425   Int_t iT0    = 1;
426   Int_t iTOF   = 1;
427   Int_t iTPC   = 1;
428   Int_t iTRD   = 1;
429   Int_t iVZERO = 1;
430   Int_t iZDC   = 1;
431
432
433     //=================== Alice BODY parameters =============================
434     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
435
436
437     if (iMAG)
438     {
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");
443     }
444 //
445 //
446     if (iABSO)
447     {
448         //=================== ABSO parameters ============================
449         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
450     }
451
452     if (iDIPO)
453     {
454         //=================== DIPO parameters ============================
455
456         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
457     }
458
459     if (iHALL)
460     {
461         //=================== HALL parameters ============================
462
463         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
464     }
465
466
467     if (iFRAME)
468     {
469         //=================== FRAME parameters ============================
470
471         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
472         if (geo == kHoles) FRAME->SetHoles(1);
473         else FRAME->SetHoles(0);
474
475     }
476
477     if (iSHIL)
478     {
479         //=================== SHIL parameters ============================
480
481         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
482     }
483
484
485     if (iPIPE)
486     {
487         //=================== PIPE parameters ============================
488
489         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
490     }
491
492     if(iITS) {
493
494     //=================== ITS parameters ============================
495     //
496
497       AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
498     }
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 AliHMPIDv3("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     }
536
537     if (iFMD)
538     {
539         //=================== FMD parameters ============================
540         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
541    }
542
543     if (iMUON)
544     {
545         //=================== MUON parameters ===========================
546         // New MUONv1 version (geometry defined via builders)
547         AliMUON *MUON = new AliMUONv1("MUON", "default");
548     }
549     //=================== PHOS parameters ===========================
550
551     if (iPHOS)
552     {
553         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
554     }
555
556
557     if (iPMD)
558     {
559         //=================== PMD parameters ============================
560         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
561     }
562
563     if (iT0)
564     {
565         //=================== T0 parameters ============================
566         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
567     }
568
569     if (iEMCAL)
570     {
571         //=================== EMCAL parameters ============================
572         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE scTh=0.176 pbTh=0.144");
573     }
574
575      if (iACORDE)
576     {
577         //=================== CRT parameters ============================
578         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
579     }
580
581      if (iVZERO)
582     {
583         //=================== CRT parameters ============================
584         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
585     }
586 }
587
588 //           PYTHIA
589
590 AliGenPythia *PythiaHVQ(PDC07Proc_t proc) {
591 //*******************************************************************//
592 // Configuration file for charm / beauty generation with PYTHIA      //
593 //                                                                   //
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.                                       //
597 //                                                                   //
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;
605
606   switch(proc) {
607   case kCharmPbPb5500:
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);
615       break;
616   case kCharmpPb8800:
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);
625       break;
626   case kCharmpp14000:
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.);
633       break;
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.);
641       break;
642   case kD0PbPb5500:
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);
650       break;
651   case kD0pPb8800:
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);
660       break;
661   case kD0pp14000:
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.);
668       break;
669   case kDPlusPbPb5500:
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);
677       break;
678   case kDPluspPb8800:
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);
687       break;
688   case kDPluspp14000:
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.);
695       break;
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);
704       break;
705   case kBeautypPb8800:
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);
714       break;
715   case kBeautypp14000:
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.);
722       break;
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.);
730       break;
731   }
732
733   return gener;
734 }
735
736 AliGenPythia *PythiaHard(PDC07Proc_t proc) {
737 //*******************************************************************//
738 // Configuration file for hard QCD processes generation with PYTHIA  //
739 //                                                                   //
740 //*******************************************************************//
741   AliGenPythia * gener = 0x0;
742
743   switch(proc) {
744
745   case kPyJetJet:
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);
755
756     return gener;
757
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
768       break;
769   case kPyJetJetPHOS:
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);
781
782       printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
783       break;
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);
797           break;
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);
812         break;
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
823     break;
824   case kPyJetJetEMCAL:
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);
837     break;
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);
850
851       printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
852       break;
853
854   }
855
856   return gener;
857 }
858
859 AliGenerator* MbCocktail()
860 {
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();
864
865 //    Pythia
866
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();
875
876 //   J/Psi parameterisation
877
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);
883
884       gener->AddGenerator(pythia, "Pythia", 1.);
885       gener->AddGenerator(jpsi,   "J/Psi", 8.e-4);
886
887       return gener;
888 }
889
890 AliGenerator* PyMbTriggered(Int_t pdg)
891 {
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);
900     return pythia;
901 }
902
903 void ProcessEnvironmentVars()
904 {
905   cout << "######################################" << endl;
906   cout << "## Processing environment variables ##" << endl;
907   cout << "######################################" << endl;
908
909     // Run Number
910     if (gSystem->Getenv("DC_RUN")) {
911       runNumber = atoi(gSystem->Getenv("DC_RUN"));
912     }
913     //cout<<"Run number "<<runNumber<<endl;
914
915     // Event Number
916     if (gSystem->Getenv("DC_EVENT")) {
917       eventNumber = atoi(gSystem->Getenv("DC_EVENT"));
918     }
919     //cout<<"Event number "<<eventNumber<<endl;
920
921
922
923     cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
924     cout<<"Seed for random number generation= "<< seed<<"  "<< gRandom->GetSeed()<<endl;
925     cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
926
927     // Run type
928     if (gSystem->Getenv("DC_RUN_TYPE")) {
929       cout<<"run type "<<gSystem->Getenv("DC_RUN_TYPE")<<endl;
930       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
931         if (strcmp(gSystem->Getenv("DC_RUN_TYPE"), pprRunName[iRun])==0) {
932           proc = (PDC07Proc_t)iRun;
933         }
934       }
935     }
936     else
937       cout << "Environment variable DC_RUN_TYPE is not defined" << endl;
938
939     if (gSystem->Getenv("ECMS"))
940       eCMS = atof(gSystem->Getenv("ECMS"));
941     if (gSystem->Getenv("PTHARDMIN"))
942       ptHardMin = atof(gSystem->Getenv("PTHARDMIN"));
943     if (gSystem->Getenv("PTHARDMAX"))
944       ptHardMax = atof(gSystem->Getenv("PTHARDMAX"));
945     if (gSystem->Getenv("PTGAMMAPI0MIN"))
946       ptGammaPi0Min = atof(gSystem->Getenv("PTGAMMAPI0MIN"));
947    if (gSystem->Getenv("QUENCHING"))
948       iquenching = atof(gSystem->Getenv("QUENCHING"));
949    if (gSystem->Getenv("QHAT"))
950       qhat = atof(gSystem->Getenv("QHAT"));
951
952     cout<<">> Run type set to "<<pprRunName[proc]<<endl;
953     cout<<">> Collision energy set to "<<eCMS <<endl;
954     cout<<">> ptHard limits: "<<ptHardMin<<" to " <<ptHardMax<<" GeV"<<endl;
955     cout<<">> pt gamma/pi0 threshold "<< ptGammaPi0Min<<" GeV "<<endl;
956     cout<<">> quenching on? "<< iquenching<<" qhat "<<qhat<<endl;
957
958     cout << "######################################" << endl;
959     cout << "######################################" << endl;
960 }