]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/QA/Config.C
Introducing event specie in QA (Yves)
[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/AliMagFCheb.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 //--- 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         // Get settings from environment variables
184   ProcessEnvironmentVars();
185
186   // Random Number seed
187   if (gSystem->Getenv("CONFIG_SEED")) {
188     seed = atoi(gSystem->Getenv("CONFIG_SEED"));
189   }
190   else if(gSystem->Getenv("DC_EVENT") && gSystem->Getenv("DC_RUN")){
191     seed = runNumber * 100000 + eventNumber;
192   }
193   gRandom->SetSeed(seed);
194   cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
195
196   // libraries required by geant321
197 #if defined(__CINT__)
198         gSystem->Load("liblhapdf");      // Parton density functions
199         gSystem->Load("libEGPythia6");   // TGenerator interface
200         gSystem->Load("libpythia6");     // Pythia
201         gSystem->Load("libAliPythia6");  // ALICE specific implementations
202         gSystem->Load("libgeant321");
203 #endif
204   new TGeant3TGeo("C++ Interface to Geant3");
205
206   // Output every 100 tracks
207   ((TGeant3*)gMC)->SetSWIT(4,100);
208
209   //=======================================================================
210   // Run loader
211   AliRunLoader* rl=0x0;
212
213   AliLog::Message(AliLog::kInfo, "Creating Run Loader", "", "", "Config()"," ConfigPPR.C", __LINE__);
214
215   rl = AliRunLoader::Open("galice.root",
216                           AliConfig::GetDefaultEventFolderName(),
217                           "recreate");
218   if (rl == 0x0)
219     {
220       gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
221       return;
222     }
223   rl->SetCompressionLevel(2);
224   rl->SetNumberOfEventsPerFile(100);
225   gAlice->SetRunLoader(rl);
226
227   // Set the trigger configuration
228   gAlice->SetTriggerDescriptor(TrigConfName[trig]);
229   cout<<"Trigger configuration is set to  "<<TrigConfName[trig]<<endl;
230
231   //======================//
232   // Set External decayer //
233   //======================//
234   TVirtualMCDecayer* decayer = new AliDecayerPythia();
235   // DECAYS
236   //
237   switch(decHvFl) {
238     case kNature:
239       decayer->SetForceDecay(kAll);
240       break;
241     case kHadr:
242       decayer->SetForceDecay(kHadronicD);
243       break;
244     case kSemiEl:
245       decayer->SetForceDecay(kSemiElectronic);
246       break;
247     case kSemiMu:
248       decayer->SetForceDecay(kSemiMuonic);
249       break;
250   }
251   decayer->Init();
252   gMC->SetExternalDecayer(decayer);
253   if(proc == kPyJetJetPHOSv2) // in this case we decay the pi0
254     decayer->SetForceDecay(kNeutralPion);
255
256   //
257   //=======================================================================
258   // ************* STEERING parameters FOR ALICE SIMULATION **************
259   // --- Specify event type to be tracked through the ALICE setup
260   // --- All positions are in cm, angles in degrees, and P and E in GeV
261
262
263     gMC->SetProcess("DCAY",1);
264     gMC->SetProcess("PAIR",1);
265     gMC->SetProcess("COMP",1);
266     gMC->SetProcess("PHOT",1);
267     gMC->SetProcess("PFIS",0);
268     gMC->SetProcess("DRAY",0);
269     gMC->SetProcess("ANNI",1);
270     gMC->SetProcess("BREM",1);
271     gMC->SetProcess("MUNU",1);
272     gMC->SetProcess("CKOV",1);
273     gMC->SetProcess("HADR",1);
274     gMC->SetProcess("LOSS",2);
275     gMC->SetProcess("MULS",1);
276     gMC->SetProcess("RAYL",1);
277
278     Float_t cut = 1.e-3;        // 1MeV cut by default
279     Float_t tofmax = 1.e10;
280
281     gMC->SetCut("CUTGAM", cut);
282     gMC->SetCut("CUTELE", cut);
283     gMC->SetCut("CUTNEU", cut);
284     gMC->SetCut("CUTHAD", cut);
285     gMC->SetCut("CUTMUO", cut);
286     gMC->SetCut("BCUTE",  cut);
287     gMC->SetCut("BCUTM",  cut);
288     gMC->SetCut("DCUTE",  cut);
289     gMC->SetCut("DCUTM",  cut);
290     gMC->SetCut("PPCUTM", cut);
291     gMC->SetCut("TOFMAX", tofmax);
292
293
294   //=========================//
295   // Generator Configuration //
296   //=========================//
297   AliGenPythia* gener = 0x0;
298
299   if (proc <=   kBeautypp14000wmi) {
300     AliGenPythia *pythia = PythiaHVQ(proc);
301     // FeedDown option
302     pythia->SetFeedDownHigherFamily(kFALSE);
303     // Stack filling option
304     if(!stars) pythia->SetStackFillOpt(AliGenPythia::kParentSelection);
305     // Set Count mode
306     if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
307     //
308     // DECAYS
309     //
310     switch(decHvFl) {
311     case kNature:
312       pythia->SetForceDecay(kAll);
313       break;
314     case kHadr:
315       pythia->SetForceDecay(kHadronicD);
316       break;
317     case kSemiEl:
318       pythia->SetForceDecay(kSemiElectronic);
319       break;
320     case kSemiMu:
321       pythia->SetForceDecay(kSemiMuonic);
322       break;
323     }
324     //
325     // GEOM & KINE CUTS
326     //
327     pythia->SetMomentumRange(0,99999999);
328     pythia->SetPhiRange(0., 360.);
329     pythia->SetThetaRange(0,180);
330     switch(ycut) {
331     case kFull:
332       pythia->SetYRange(-999,999);
333       break;
334     case kBarrel:
335       pythia->SetYRange(-2,2);
336       break;
337     case kMuonArm:
338       pythia->SetYRange(1,6);
339       break;
340     }
341     gener = pythia;
342   } else if (proc == kPyMbNoHvq) {
343     gener = MbCocktail();
344   } else if (proc == kPyOmegaMinus) {
345     gener = PyMbTriggered(3334);
346   } else if (proc == kPyOmegaPlus) {
347     gener = PyMbTriggered(-3334);
348   } else if (proc <= kPyGammaBremsEMCAL) {
349     AliGenPythia *pythia = PythiaHard(proc);
350     // FeedDown option
351     pythia->SetFeedDownHigherFamily(kFALSE);
352     // Set Count mode
353     if(nEvts>0) pythia->SetCountMode(AliGenPythia::kCountParents);
354
355       //
356       // GEOM & KINE CUTS
357       //
358     pythia->SetMomentumRange(0,99999999);
359     //       pythia->SetJetEtaRange(-1.5, 1.5);//  Final state kinematic cuts
360     //       pythia->SetJetPhiRange(0., 360.);
361     //       pythia->SetThetaRange(45,135);
362
363    if(proc == kPyJetJetPHOSv2)
364                  pythia->SetForceDecay(kNeutralPion);
365          else
366                  pythia->SetForceDecay(kAll);
367                 
368                 pythia->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
369     pythia->SetPtKick(5); // set the intrinsic kt to 5 GeV/c
370     gener = pythia;
371         }
372
373
374   // PRIMARY VERTEX
375
376   gener->SetOrigin(0., 0., 0.);    // vertex position
377
378   // Size of the interaction diamond
379   // Longitudinal
380   Float_t sigmaz  = 7.55 / TMath::Sqrt(2.); // [cm]
381
382   // Transverse
383   Float_t betast  = 10;                 // beta* [m]
384   Float_t eps     = 3.75e-6;            // emittance [m]
385   Float_t gamma   = 7000. / 0.938272;   // relativistic gamma [1]
386   Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.;  // [cm]
387   printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
388
389   gener->SetSigma(sigmaxy, sigmaxy, sigmaz);      // Sigma in (X,Y,Z) (cm) on IP position
390   gener->SetCutVertexZ(3.);        // Truncate at 3 sigma
391   gener->SetVertexSmear(kPerEvent);
392
393   gener->Init();
394
395   //Quenching
396   gener->SetQuench(iquenching);
397   if(iquenching == 1){
398     Float_t k = 6e5*(qhat/1.7) ; //qhat=1.7, k = 6e5, default  value
399     AliPythia::Instance()->InitQuenching(0.,0.1,k,0,0.95,6);
400
401  }
402   // FIELD
403
404   if (mag == k2kG) {
405     comment = comment.Append(" | L3 field 0.2 T");
406   } else if (mag == k4kG) {
407     comment = comment.Append(" | L3 field 0.4 T");
408   } else if (mag == k5kG) {
409     comment = comment.Append(" | L3 field 0.5 T");
410   }
411   printf("\n \n Comment: %s \n \n", comment.Data());
412
413   AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., mag);
414   field->SetL3ConstField(0); //Using const. field in the barrel
415   rl->CdGAFile();
416   gAlice->SetField(field);    
417
418
419   Int_t iABSO  = 1;
420   Int_t iACORDE= 1;
421   Int_t iDIPO  = 1;
422   Int_t iEMCAL = 1;
423   Int_t iFMD   = 1;
424   Int_t iFRAME = 1;
425   Int_t iHALL  = 1;
426   Int_t iITS   = 1;
427   Int_t iMAG   = 1;
428   Int_t iMUON  = 1;
429   Int_t iPHOS  = 1;
430   Int_t iPIPE  = 1;
431   Int_t iPMD   = 1;
432   Int_t iHMPID = 1;
433   Int_t iSHIL  = 1;
434   Int_t iT0    = 1;
435   Int_t iTOF   = 1;
436   Int_t iTPC   = 1;
437   Int_t iTRD   = 1;
438   Int_t iVZERO = 1;
439   Int_t iZDC   = 1;
440
441
442     //=================== Alice BODY parameters =============================
443     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
444
445
446     if (iMAG)
447     {
448         //=================== MAG parameters ============================
449         // --- Start with Magnet since detector layouts may be depending ---
450         // --- on the selected Magnet dimensions ---
451         AliMAG *MAG = new AliMAG("MAG", "Magnet");
452     }
453 //
454 //
455     if (iABSO)
456     {
457         //=================== ABSO parameters ============================
458         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
459     }
460
461     if (iDIPO)
462     {
463         //=================== DIPO parameters ============================
464
465         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
466     }
467
468     if (iHALL)
469     {
470         //=================== HALL parameters ============================
471
472         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
473     }
474
475
476     if (iFRAME)
477     {
478         //=================== FRAME parameters ============================
479
480         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
481         if (geo == kHoles) FRAME->SetHoles(1);
482         else FRAME->SetHoles(0);
483
484     }
485
486     if (iSHIL)
487     {
488         //=================== SHIL parameters ============================
489
490         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
491     }
492
493
494     if (iPIPE)
495     {
496         //=================== PIPE parameters ============================
497
498         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
499     }
500
501     if(iITS) {
502
503     //=================== ITS parameters ============================
504     //
505
506       AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
507     }
508
509
510     if (iTPC)
511     {
512       //============================ TPC parameters =====================
513         AliTPC *TPC = new AliTPCv2("TPC", "Default");
514     }
515
516
517     if (iTOF) {
518         //=================== TOF parameters ============================
519         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
520    }
521
522
523     if (iHMPID)
524     {
525         //=================== HMPID parameters ===========================
526         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
527
528     }
529
530
531     if (iZDC)
532     {
533         //=================== ZDC parameters ============================
534
535         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
536     }
537
538    if (iTRD)
539     {
540         //=================== TRD parameters ============================
541  
542         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
543         AliTRDgeometry *geoTRD = TRD->GetGeometry();
544     }
545
546     if (iFMD)
547     {
548         //=================== FMD parameters ============================
549         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
550    }
551
552     if (iMUON)
553     {
554         //=================== MUON parameters ===========================
555         // New MUONv1 version (geometry defined via builders)
556         AliMUON *MUON = new AliMUONv1("MUON", "default");
557     }
558     //=================== PHOS parameters ===========================
559
560     if (iPHOS)
561     {
562         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
563     }
564
565
566     if (iPMD)
567     {
568         //=================== PMD parameters ============================
569         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
570     }
571
572     if (iT0)
573     {
574         //=================== T0 parameters ============================
575         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
576     }
577
578     if (iEMCAL)
579     {
580         //=================== EMCAL parameters ============================
581         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE scTh=0.176 pbTh=0.144");
582     }
583
584      if (iACORDE)
585     {
586         //=================== CRT parameters ============================
587         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
588     }
589
590      if (iVZERO)
591     {
592         //=================== CRT parameters ============================
593         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
594     }
595 }
596
597 //           PYTHIA
598
599 AliGenPythia *PythiaHVQ(PDC07Proc_t proc) {
600 //*******************************************************************//
601 // Configuration file for charm / beauty generation with PYTHIA      //
602 //                                                                   //
603 // The parameters have been tuned in order to reproduce the inclusive//
604 // heavy quark pt distribution given by the NLO pQCD calculation by  //
605 // Mangano, Nason and Ridolfi.                                       //
606 //                                                                   //
607 // For details and for the NORMALIZATION of the yields see:          //
608 //   N.Carrer and A.Dainese,                                         //
609 //   "Charm and beauty production at the LHC",                       //
610 //   ALICE-INT-2003-019, [arXiv:hep-ph/0311225];                     //
611 //   PPR Chapter 6.6, CERN/LHCC 2005-030 (2005).                     //
612 //*******************************************************************//
613   AliGenPythia * gener = 0x0;
614
615   switch(proc) {
616   case kCharmPbPb5500:
617       comment = comment.Append(" Charm in Pb-Pb at 5.5 TeV");
618       gener = new AliGenPythia(nEvts);
619       gener->SetProcess(kPyCharmPbPbMNR);
620       gener->SetStrucFunc(kCTEQ4L);
621       gener->SetPtHard(2.1,-1.0);
622       gener->SetEnergyCMS(5500.);
623       gener->SetNuclei(208,208);
624       break;
625   case kCharmpPb8800:
626       comment = comment.Append(" Charm in p-Pb at 8.8 TeV");
627       gener = new AliGenPythia(nEvts);
628       gener->SetProcess(kPyCharmpPbMNR);
629       gener->SetStrucFunc(kCTEQ4L);
630       gener->SetPtHard(2.1,-1.0);
631       gener->SetEnergyCMS(8800.);
632       gener->SetProjectile("P",1,1);
633       gener->SetTarget("Pb",208,82);
634       break;
635   case kCharmpp14000:
636       comment = comment.Append(" Charm in pp at 14 TeV");
637       gener = new AliGenPythia(nEvts);
638       gener->SetProcess(kPyCharmppMNR);
639       gener->SetStrucFunc(kCTEQ4L);
640       gener->SetPtHard(2.1,-1.0);
641       gener->SetEnergyCMS(14000.);
642       break;
643   case kCharmpp14000wmi:
644       comment = comment.Append(" Charm in pp at 14 TeV with mult. interactions");
645       gener = new AliGenPythia(-1);
646       gener->SetProcess(kPyCharmppMNRwmi);
647       gener->SetStrucFunc(kCTEQ5L);
648       gener->SetPtHard(ptHardMin,ptHardMax);
649       gener->SetEnergyCMS(14000.);
650       break;
651   case kD0PbPb5500:
652       comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
653       gener = new AliGenPythia(nEvts);
654       gener->SetProcess(kPyD0PbPbMNR);
655       gener->SetStrucFunc(kCTEQ4L);
656       gener->SetPtHard(2.1,-1.0);
657       gener->SetEnergyCMS(5500.);
658       gener->SetNuclei(208,208);
659       break;
660   case kD0pPb8800:
661       comment = comment.Append(" D0 in p-Pb at 8.8 TeV");
662       gener = new AliGenPythia(nEvts);
663       gener->SetProcess(kPyD0pPbMNR);
664       gener->SetStrucFunc(kCTEQ4L);
665       gener->SetPtHard(2.1,-1.0);
666       gener->SetEnergyCMS(8800.);
667       gener->SetProjectile("P",1,1);
668       gener->SetTarget("Pb",208,82);
669       break;
670   case kD0pp14000:
671       comment = comment.Append(" D0 in pp at 14 TeV");
672       gener = new AliGenPythia(nEvts);
673       gener->SetProcess(kPyD0ppMNR);
674       gener->SetStrucFunc(kCTEQ4L);
675       gener->SetPtHard(2.1,-1.0);
676       gener->SetEnergyCMS(14000.);
677       break;
678   case kDPlusPbPb5500:
679       comment = comment.Append(" DPlus in Pb-Pb at 5.5 TeV");
680       gener = new AliGenPythia(nEvts);
681       gener->SetProcess(kPyDPlusPbPbMNR);
682       gener->SetStrucFunc(kCTEQ4L);
683       gener->SetPtHard(2.1,-1.0);
684       gener->SetEnergyCMS(5500.);
685       gener->SetNuclei(208,208);
686       break;
687   case kDPluspPb8800:
688       comment = comment.Append(" DPlus in p-Pb at 8.8 TeV");
689       gener = new AliGenPythia(nEvts);
690       gener->SetProcess(kPyDPluspPbMNR);
691       gener->SetStrucFunc(kCTEQ4L);
692       gener->SetPtHard(2.1,-1.0);
693       gener->SetEnergyCMS(8800.);
694       gener->SetProjectile("P",1,1);
695       gener->SetTarget("Pb",208,82);
696       break;
697   case kDPluspp14000:
698       comment = comment.Append(" DPlus in pp at 14 TeV");
699       gener = new AliGenPythia(nEvts);
700       gener->SetProcess(kPyDPlusppMNR);
701       gener->SetStrucFunc(kCTEQ4L);
702       gener->SetPtHard(2.1,-1.0);
703       gener->SetEnergyCMS(14000.);
704       break;
705   case kBeautyPbPb5500:
706       comment = comment.Append(" Beauty in Pb-Pb at 5.5 TeV");
707       gener = new AliGenPythia(nEvts);
708       gener->SetProcess(kPyBeautyPbPbMNR);
709       gener->SetStrucFunc(kCTEQ4L);
710       gener->SetPtHard(2.75,-1.0);
711       gener->SetEnergyCMS(5500.);
712       gener->SetNuclei(208,208);
713       break;
714   case kBeautypPb8800:
715       comment = comment.Append(" Beauty in p-Pb at 8.8 TeV");
716       gener = new AliGenPythia(nEvts);
717       gener->SetProcess(kPyBeautypPbMNR);
718       gener->SetStrucFunc(kCTEQ4L);
719       gener->SetPtHard(2.75,-1.0);
720       gener->SetEnergyCMS(8800.);
721       gener->SetProjectile("P",1,1);
722       gener->SetTarget("Pb",208,82);
723       break;
724   case kBeautypp14000:
725       comment = comment.Append(" Beauty in pp at 14 TeV");
726       gener = new AliGenPythia(nEvts);
727       gener->SetProcess(kPyBeautyppMNR);
728       gener->SetStrucFunc(kCTEQ4L);
729       gener->SetPtHard(2.75,-1.0);
730       gener->SetEnergyCMS(14000.);
731       break;
732   case kBeautypp14000wmi:
733       comment = comment.Append(" Beauty in pp at 14 TeV with mult. interactions");
734       gener = new AliGenPythia(-1);
735       gener->SetProcess(kPyBeautyppMNRwmi);
736       gener->SetStrucFunc(kCTEQ5L);
737       gener->SetPtHard(ptHardMin,ptHardMax);
738       gener->SetEnergyCMS(14000.);
739       break;
740   }
741
742   return gener;
743 }
744
745 AliGenPythia *PythiaHard(PDC07Proc_t proc) {
746 //*******************************************************************//
747 // Configuration file for hard QCD processes generation with PYTHIA  //
748 //                                                                   //
749 //*******************************************************************//
750   AliGenPythia * gener = 0x0;
751
752   switch(proc) {
753
754   case kPyJetJet:
755     comment = comment.Append(" pp->jet + jet over at 14 TeV, no restriction");
756     AliGenPythia * gener = new AliGenPythia(nEvts);
757     gener->SetEnergyCMS(eCMS);//        Centre of mass energy
758     gener->SetProcess(kPyJets);//        Process type
759     gener->SetJetEtaRange(-1.5, 1.5);//  Final state kinematic cuts
760     gener->SetJetPhiRange(0., 360.);
761     gener->SetJetEtRange(10., 1000.);
762     gener->SetPtHard(ptHardMin, ptHardMax);// Pt transfer of the hard scattering
763     gener->SetStrucFunc(kCTEQ4L);
764
765     return gener;
766
767   case kPyGammaJetPHOS:
768       comment = comment.Append(" pp->jet + gamma over PHOS");
769       gener = new AliGenPythia(nEvts);
770       gener->SetEnergyCMS(eCMS);
771       gener->SetProcess(kPyDirectGamma);
772       gener->SetStrucFunc(kCTEQ4L);
773       gener->SetPtHard(ptHardMin,ptHardMax);
774       //gener->SetYHard(-1.0,1.0);
775       gener->SetGammaEtaRange(-0.13,0.13);
776       gener->SetGammaPhiRange(218.,322.);//Over 5 modules +-2 deg
777       break;
778   case kPyJetJetPHOS:
779       comment = comment.Append(" pp->jet + jet over PHOS");
780       gener = new AliGenPythia(nEvts);
781       gener->SetEnergyCMS(eCMS);
782       gener->SetProcess(kPyJets);
783       gener->SetStrucFunc(kCTEQ4L);
784       gener->SetPtHard(ptHardMin,ptHardMax);
785       //gener->SetYHard(-1.0,1.0);
786       gener->SetJetEtaRange(-1.,1.);
787       gener->SetJetPhiRange(200.,340.);
788       gener->SetPi0InPHOS(kTRUE);
789       gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
790
791       printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
792       break;
793   case kPyGammaBremsPHOS:
794       comment = comment.Append(" pp->jet + jet+bremsphoton over PHOS at 14 TeV");
795       gener = new AliGenPythia(nEvts);
796       gener->SetEnergyCMS(eCMS);
797       gener->SetProcess(kPyJets);
798       gener->SetStrucFunc(kCTEQ4L);
799       gener->SetPtHard(ptHardMin,ptHardMax);
800       //gener->SetYHard(-1.0,1.0);
801       gener->SetJetEtaRange(-1.,1.);
802       gener->SetJetPhiRange(200.,340.);
803       gener->SetFragPhotonInPHOS(kTRUE);
804       gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
805           printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
806           break;
807         case kPyJetJetPHOSv2:
808         comment = comment.Append(" pp->jet + jet over PHOS version2 ");
809         gener = new AliGenPythia(nEvts);
810         gener->SetEnergyCMS(eCMS);
811         gener->SetProcess(kPyJets);
812         gener->SetStrucFunc(kCTEQ4L);
813         gener->SetPtHard(ptHardMin,ptHardMax);
814         //gener->SetYHard(-1.0,1.0);
815         gener->SetJetEtaRange(-1.,1.);
816         gener->SetJetPhiRange(200.,340.);
817         //gener->SetPi0InPHOS(kTRUE);
818         gener->SetPhotonInPHOSeta(kTRUE);
819         gener->SetPhotonMinPt(ptGammaPi0Min);
820         gener->SetForceDecay(kAll);
821         break;
822   case kPyGammaJetEMCAL:
823     comment = comment.Append(" pp->jet + gamma over EMCAL at 14 TeV");
824     gener = new AliGenPythia(nEvts);
825     gener->SetEnergyCMS(eCMS);
826     gener->SetProcess(kPyDirectGamma);
827     gener->SetStrucFunc(kCTEQ4L);
828     gener->SetPtHard(ptHardMin,ptHardMax);
829     //gener->SetYHard(-1.0,1.0);
830     gener->SetGammaEtaRange(-0.71,0.71);
831     gener->SetGammaPhiRange(78.,192.);//Over 6 supermodules +-2 deg
832     break;
833   case kPyJetJetEMCAL:
834     comment = comment.Append(" pp->jet + jet over EMCAL at 14 TeV");
835     gener = new AliGenPythia(nEvts);
836     gener->SetEnergyCMS(eCMS);
837     gener->SetProcess(kPyJets);
838     gener->SetStrucFunc(kCTEQ4L);
839     gener->SetPtHard(ptHardMin,ptHardMax);
840     //gener->SetYHard(-1.0,1.0);
841     gener->SetJetEtaRange(-1,1);
842     gener->SetJetPhiRange(60.,210.);
843     gener->SetPi0InEMCAL(kTRUE);
844     gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
845     printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
846     break;
847   case kPyGammaBremsEMCAL:
848       comment = comment.Append(" pp->jet + jet+bremsphoton over EMCAL at 14 TeV");
849       gener = new AliGenPythia(nEvts);
850       gener->SetEnergyCMS(eCMS);
851       gener->SetProcess(kPyJets);
852       gener->SetStrucFunc(kCTEQ4L);
853       gener->SetPtHard(ptHardMin,ptHardMax);
854       //gener->SetYHard(-1.0,1.0);
855       gener->SetJetEtaRange(-1,1);
856       gener->SetJetPhiRange(60.,210.); //Over 2 uncovered PHOS modules
857       gener->SetFragPhotonInEMCAL(kTRUE);
858       gener->SetFragPhotonOrPi0MinPt(ptGammaPi0Min);
859
860       printf("\n \n Event generator: Minimum pT of particle in calorimeter %f \n \n", ptGammaPi0Min);
861       break;
862
863   }
864
865   return gener;
866 }
867
868 AliGenerator* MbCocktail()
869 {
870       comment = comment.Append(" pp at 14 TeV: Pythia low-pt, no heavy quarks + J/Psi from parameterisation");
871       AliGenCocktail * gener = new AliGenCocktail();
872       gener->UsePerEventRates();
873
874 //    Pythia
875
876       AliGenPythia* pythia = new AliGenPythia(-1);
877       pythia->SetMomentumRange(0, 999999.);
878       pythia->SetThetaRange(0., 180.);
879       pythia->SetYRange(-12.,12.);
880       pythia->SetPtRange(0,1000.);
881       pythia->SetProcess(kPyMb);
882       pythia->SetEnergyCMS(14000.);
883       pythia->SwitchHFOff();
884
885 //   J/Psi parameterisation
886
887       AliGenParam* jpsi = new AliGenParam(1, AliGenMUONlib::kJpsi, "CDF scaled", "Jpsi");
888       jpsi->SetPtRange(0.,100.);
889       jpsi->SetYRange(-8., 8.);
890       jpsi->SetPhiRange(0., 360.);
891       jpsi->SetForceDecay(kAll);
892
893       gener->AddGenerator(pythia, "Pythia", 1.);
894       gener->AddGenerator(jpsi,   "J/Psi", 8.e-4);
895
896       return gener;
897 }
898
899 AliGenerator* PyMbTriggered(Int_t pdg)
900 {
901     AliGenPythia* pythia = new AliGenPythia(-1);
902     pythia->SetMomentumRange(0, 999999.);
903     pythia->SetThetaRange(0., 180.);
904     pythia->SetYRange(-12.,12.);
905     pythia->SetPtRange(0,1000.);
906     pythia->SetProcess(kPyMb);
907     pythia->SetEnergyCMS(14000.);
908     pythia->SetTriggerParticle(pdg, 0.9);
909     return pythia;
910 }
911
912 void ProcessEnvironmentVars()
913 {
914   cout << "######################################" << endl;
915   cout << "## Processing environment variables ##" << endl;
916   cout << "######################################" << endl;
917
918     // Run Number
919     if (gSystem->Getenv("DC_RUN")) {
920       runNumber = atoi(gSystem->Getenv("DC_RUN"));
921     }
922     //cout<<"Run number "<<runNumber<<endl;
923
924     // Event Number
925     if (gSystem->Getenv("DC_EVENT")) {
926       eventNumber = atoi(gSystem->Getenv("DC_EVENT"));
927     }
928     //cout<<"Event number "<<eventNumber<<endl;
929
930
931
932     cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
933     cout<<"Seed for random number generation= "<< seed<<"  "<< gRandom->GetSeed()<<endl;
934     cout<<"////////////////////////////////////////////////////////////////////////////////////"<<endl;
935
936     // Run type
937     if (gSystem->Getenv("DC_RUN_TYPE")) {
938       cout<<"run type "<<gSystem->Getenv("DC_RUN_TYPE")<<endl;
939       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
940         if (strcmp(gSystem->Getenv("DC_RUN_TYPE"), pprRunName[iRun])==0) {
941           proc = (PDC07Proc_t)iRun;
942         }
943       }
944     }
945     else
946       cout << "Environment variable DC_RUN_TYPE is not defined" << endl;
947
948     if (gSystem->Getenv("ECMS"))
949       eCMS = atof(gSystem->Getenv("ECMS"));
950     if (gSystem->Getenv("PTHARDMIN"))
951       ptHardMin = atof(gSystem->Getenv("PTHARDMIN"));
952     if (gSystem->Getenv("PTHARDMAX"))
953       ptHardMax = atof(gSystem->Getenv("PTHARDMAX"));
954     if (gSystem->Getenv("PTGAMMAPI0MIN"))
955       ptGammaPi0Min = atof(gSystem->Getenv("PTGAMMAPI0MIN"));
956    if (gSystem->Getenv("QUENCHING"))
957       iquenching = atof(gSystem->Getenv("QUENCHING"));
958    if (gSystem->Getenv("QHAT"))
959       qhat = atof(gSystem->Getenv("QHAT"));
960
961     cout<<">> Run type set to "<<pprRunName[proc]<<endl;
962     cout<<">> Collision energy set to "<<eCMS <<endl;
963     cout<<">> ptHard limits: "<<ptHardMin<<" to " <<ptHardMax<<" GeV"<<endl;
964     cout<<">> pt gamma/pi0 threshold "<< ptGammaPi0Min<<" GeV "<<endl;
965     cout<<">> quenching on? "<< iquenching<<" qhat "<<qhat<<endl;
966
967     cout << "######################################" << endl;
968     cout << "######################################" << endl;
969 }