]> git.uio.no Git - u/mrichter/AliRoot.git/blob - test/vmctest/production/Config.C
Added macros for production using Geant4. 3 different physics lists can be used:...
[u/mrichter/AliRoot.git] / test / vmctest / production / Config.C
1 //
2 // Configuration for the first physics production 2008
3 //
4
5 // One can use the configuration macro in compiled mode by
6 // root [0] gSystem->Load("libgeant321");
7 // root [0] gSystem->SetIncludePath("-I$ROOTSYS/include -I$ALICE_ROOT/include\
8 //                   -I$ALICE_ROOT -I$ALICE/geant3/TGeant3");
9 // root [0] .x grun.C(1,"Config.C++")
10
11 #if !defined(__CINT__) || defined(__MAKECINT__)
12 #include <Riostream.h>
13 #include <TRandom.h>
14 #include <TDatime.h>
15 #include <TSystem.h>
16 #include <TVirtualMC.h>
17 #include <TGeant3TGeo.h>
18 #include "STEER/AliRunLoader.h"
19 #include "STEER/AliRun.h"
20 #include "STEER/AliConfig.h"
21 #include "PYTHIA6/AliDecayerPythia.h"
22 #include "PYTHIA6/AliGenPythia.h"
23 #include "TDPMjet/AliGenDPMjet.h"
24 #include "STEER/AliMagFCheb.h"
25 #include "STRUCT/AliBODY.h"
26 #include "STRUCT/AliMAG.h"
27 #include "STRUCT/AliABSOv3.h"
28 #include "STRUCT/AliDIPOv3.h"
29 #include "STRUCT/AliHALLv3.h"
30 #include "STRUCT/AliFRAMEv2.h"
31 #include "STRUCT/AliSHILv3.h"
32 #include "STRUCT/AliPIPEv3.h"
33 #include "ITS/AliITSv11Hybrid.h"
34 #include "TPC/AliTPCv2.h"
35 #include "TOF/AliTOFv6T0.h"
36 #include "HMPID/AliHMPIDv3.h"
37 #include "ZDC/AliZDCv3.h"
38 #include "TRD/AliTRDv1.h"
39 #include "TRD/AliTRDgeometry.h"
40 #include "FMD/AliFMDv1.h"
41 #include "MUON/AliMUONv1.h"
42 #include "PHOS/AliPHOSv1.h"
43 #include "PHOS/AliPHOSSimParam.h"
44 #include "PMD/AliPMDv1.h"
45 #include "T0/AliT0v1.h"
46 #include "EMCAL/AliEMCALv2.h"
47 #include "ACORDE/AliACORDEv1.h"
48 #include "VZERO/AliVZEROv7.h"
49 #endif
50
51
52 enum PDC06Proc_t 
53 {
54   kPythia6, kPythia6D6T, kPythia6D6T_Flat, kPythia6ATLAS, kPythia6ATLAS_Flat, kPythia6_Perugia0, kPhojet, kRunMax
55 };
56
57 const char * pprRunName[] = {
58     "kPythia6", "kPythia6D6T", "kPythia6D6T_Flat", "kPythia6ATLAS", "kPythia6ATLAS_Flat", "kPythia6_Perugia0", "kPhojet" 
59 };
60
61 enum Mag_t
62 {
63   kNoField, k5kG, kFieldMax
64 };
65
66 const char * pprField[] = {
67   "kNoField", "k5kG"
68 };
69
70 enum PhysicsList_t 
71 {
72   QGSP_BERT_EMV,CHIPS,QGSP_BERT_CHIPS,QGSP_BERT_EMV_OPTICAL, CHIPS_OPTICAL, QGSP_BERT_CHIPS_OPTICAL, kListMax
73 };
74
75 const char * physicsListName[] = {
76   "QGSP_BERT_EMV",         "CHIPS",         "QGSP_BERT_CHIPS", 
77   "QGSP_BERT_EMV_OPTICAL", "CHIPS_OPTICAL", "QGSP_BERT_CHIPS_OPTICAL"
78 };
79
80 //--- Functions ---
81 class AliGenPythia;
82 AliGenerator *MbPythia();
83 AliGenerator *MbPythiaTuneD6T();
84 AliGenerator *MbPhojet();
85 void ProcessEnvironmentVars();
86
87 // Geterator, field, beam energy
88 static PDC06Proc_t   proc         = kPhojet;
89 static Mag_t         mag          = k5kG;
90 static Float_t       energy       = 10000; // energy in CMS
91 static PhysicsList_t physicslist  = QGSP_BERT_EMV;
92 //========================//
93 // Set Random Number seed //
94 //========================//
95 TDatime dt;
96 static UInt_t seed    = dt.Get();
97
98 // Comment line
99 static TString comment;
100
101 void Config()
102 {
103     
104
105   // Get settings from environment variables
106   ProcessEnvironmentVars();
107
108   gRandom->SetSeed(seed);
109   cerr<<"Seed for random number generation= "<<seed<<endl; 
110
111   // Libraries required by geant321
112 #if defined(__CINT__)
113   gSystem->Load("liblhapdf");       // Parton density functions
114   gSystem->Load("libEGPythia6");    // TGenerator interface
115
116   if (proc == kPythia6 || proc == kPhojet) {
117       gSystem->Load("libpythia6");        // Pythia 6.2 
118   } else {
119       gSystem->Load("libpythia6.4.21");   // Pythia 6.4
120   } 
121
122   gSystem->Load("libAliPythia6");   // ALICE specific implementations
123   //gSystem->Load("libgeant321");
124 #endif
125
126   //  new TGeant3TGeo("C++ Interface to Geant3");
127
128   //=======================================================================
129   //  Create the output file
130
131    
132   AliRunLoader* rl=0x0;
133
134   cout<<"Config.C: Creating Run Loader ..."<<endl;
135   rl = AliRunLoader::Open("galice.root",
136                           AliConfig::GetDefaultEventFolderName(),
137                           "recreate");
138   if (rl == 0x0)
139     {
140       gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
141       return;
142     }
143   rl->SetCompressionLevel(2);
144   rl->SetNumberOfEventsPerFile(1000);
145   gAlice->SetRunLoader(rl);
146   // gAlice->SetGeometryFromFile("geometry.root");
147   // gAlice->SetGeometryFromCDB();
148   
149   // Set the trigger configuration: proton-proton
150   gAlice->SetTriggerDescriptor("p-p");
151   //  AliSimulation::Instance()->SetTriggerConfig("p-p");
152
153  
154   //
155   // FIELD
156   //
157   AliMagF* field = 0x0;
158   if (mag == kNoField) {
159     comment = comment.Append(" | L3 field 0.0 T");
160     field = new AliMagF("Maps","Maps", 0., 0., AliMagF::k5kGUniform,AliMagF::kBeamTypepp, energy/2.0);
161   } else if (mag == k5kG) {
162     comment = comment.Append(" | L3 field 0.5 T");
163     field = new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG, AliMagF::kBeamTypepp, energy/2.0);
164   }
165
166   printf("\n \n Comment: %s \n \n", comment.Data());
167
168   TGeoGlobalMagField::Instance()->SetField(field);
169     
170   rl->CdGAFile();
171   
172   Int_t iABSO  = 1;
173   Int_t iACORDE= 0;
174   Int_t iDIPO  = 1;
175   Int_t iEMCAL = 1;
176   Int_t iFMD   = 1;
177   Int_t iFRAME = 1;
178   Int_t iHALL  = 1;
179   Int_t iITS   = 1;
180   Int_t iMAG   = 1;
181   Int_t iMUON  = 1;
182   Int_t iPHOS  = 1;
183   Int_t iPIPE  = 1;
184   Int_t iPMD   = 1;
185   Int_t iHMPID = 1;
186   Int_t iSHIL  = 1;
187   Int_t iT0    = 1;
188   Int_t iTOF   = 1;
189   Int_t iTPC   = 1;
190   Int_t iTRD   = 1;
191   Int_t iVZERO = 1;
192   Int_t iZDC   = 1;
193   
194
195     //=================== Alice BODY parameters =============================
196     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
197
198
199     if (iMAG)
200     {
201         //=================== MAG parameters ============================
202         // --- Start with Magnet since detector layouts may be depending ---
203         // --- on the selected Magnet dimensions ---
204         AliMAG *MAG = new AliMAG("MAG", "Magnet");
205     }
206
207
208     if (iABSO)
209     {
210         //=================== ABSO parameters ============================
211         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
212     }
213
214     if (iDIPO)
215     {
216         //=================== DIPO parameters ============================
217
218         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
219     }
220
221     if (iHALL)
222     {
223         //=================== HALL parameters ============================
224
225         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
226     }
227
228
229     if (iFRAME)
230     {
231         //=================== FRAME parameters ============================
232
233         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
234         FRAME->SetHoles(1);
235     }
236
237     if (iSHIL)
238     {
239         //=================== SHIL parameters ============================
240
241         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
242     }
243
244
245     if (iPIPE)
246     {
247         //=================== PIPE parameters ============================
248
249         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
250     }
251  
252     if (iITS)
253     {
254         //=================== ITS parameters ============================
255
256         AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
257     }
258
259     if (iTPC)
260     {
261       //============================ TPC parameters =====================
262
263         AliTPC *TPC = new AliTPCv2("TPC", "Default");
264         TPC->SetPrimaryIonisation(); // not used with Geant3
265     }
266
267
268     if (iTOF) {
269         //=================== TOF parameters ============================
270
271         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
272     }
273
274
275     if (iHMPID)
276     {
277         //=================== HMPID parameters ===========================
278
279         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
280
281     }
282
283
284     if (iZDC)
285     {
286         //=================== ZDC parameters ============================
287
288         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
289     }
290
291     if (iTRD)
292     {
293         //=================== TRD parameters ============================
294
295         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
296         AliTRDgeometry *geoTRD = TRD->GetGeometry();
297         // Partial geometry: modules at 0,1,7,8,9,16,17
298         // starting at 3h in positive direction
299         geoTRD->SetSMstatus(2,0);
300         geoTRD->SetSMstatus(3,0);
301         geoTRD->SetSMstatus(4,0);
302         geoTRD->SetSMstatus(5,0);
303         geoTRD->SetSMstatus(6,0);
304         geoTRD->SetSMstatus(11,0);
305         geoTRD->SetSMstatus(12,0);
306         geoTRD->SetSMstatus(13,0);
307         geoTRD->SetSMstatus(14,0);
308         geoTRD->SetSMstatus(15,0);
309         geoTRD->SetSMstatus(16,0);
310     }
311
312     if (iFMD)
313     {
314         //=================== FMD parameters ============================
315
316         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
317    }
318
319     if (iMUON)
320     {
321         //=================== MUON parameters ===========================
322         // New MUONv1 version (geometry defined via builders)
323
324         AliMUON *MUON = new AliMUONv1("MUON", "default");
325     }
326
327     if (iPHOS)
328     {
329         //=================== PHOS parameters ===========================
330
331         AliPHOS *PHOS = new AliPHOSv1("PHOS", "noCPV");
332         //Set simulation parameters different from the default ones.
333         AliPHOSSimParam* simEmc = AliPHOSSimParam::GetInstance() ;
334   
335         // APD noise of warm (+20C) PHOS:
336         // a2 = a1*(Y1/Y2)*(M1/M2), where a1 = 0.012 is APD noise at -25C,
337         // Y1 = 4.3 photo-electrons/MeV, Y2 = 1.7 p.e/MeV - light yields at -25C and +20C,
338         // M1 = 50, M2 = 50 - APD gain factors chosen for t1 = -25C and t2 = +20C,
339         // Y = MeanLightYield*APDEfficiency.
340
341         Float_t apdNoise = 0.012*2.5; 
342         simEmc->SetAPDNoise(apdNoise);
343
344         //Raw Light Yield at +20C
345         simEmc->SetMeanLightYield(18800);
346
347         //ADC channel width at +18C.
348         simEmc->SetADCchannelW(0.0125);
349     }
350
351
352     if (iPMD)
353     {
354         //=================== PMD parameters ============================
355
356         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
357     }
358
359     if (iT0)
360     {
361         //=================== T0 parameters ============================
362         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
363     }
364
365     if (iEMCAL)
366     {
367         //=================== EMCAL parameters ============================
368
369         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_FIRSTYEAR");
370     }
371
372      if (iACORDE)
373     {
374         //=================== ACORDE parameters ============================
375
376         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
377     }
378
379      if (iVZERO)
380     {
381         //=================== ACORDE parameters ============================
382
383         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
384     }
385
386
387
388   // Load Geant4 + Geant4 VMC libraries
389   //
390   if (gClassTable->GetID("TGeant4") == -1) {
391     // Load Geant4 libraries 
392     if (!gInterpreter->IsLoaded("$G4INSTALL/macro/g4libs.C")) {
393       gROOT->LoadMacro("$G4INSTALL/macro/g4libs.C");
394       gInterpreter->ProcessLine("g4libs()");
395     }
396   }    
397
398
399   // Create Geant4 VMC
400   //  
401   TGeant4 *geant4 = 0;
402   if ( ! gMC ) {
403     TG4RunConfiguration* runConfiguration=0x0;
404     for (Int_t iList = 0; iList < kListMax; iList++) {
405       if(iList<kListMax/2){
406         if(physicslist == iList){
407           runConfiguration = 
408             new TG4RunConfiguration("geomRoot", 
409                                     physicsListName[iList], 
410                                     "specialCuts+stackPopper+stepLimiter",
411                                     true);
412         }
413       }
414       else if(iList>=kListMax/2){//add "optical" PL to HadronPhysicsList
415         if(physicslist == iList){
416           runConfiguration = 
417             new TG4RunConfiguration("geomRoot", 
418                                     Form("%s+optical",physicsListName[iList-kListMax/2]), 
419                                     "specialCuts+stackPopper+stepLimiter",
420                                     true);
421         }
422       }
423     }
424     geant4 = new TGeant4("TGeant4", "The Geant4 Monte Carlo", runConfiguration);
425     cout << "Geant4 has been created." << endl;
426   } 
427   else {
428     cout << "Monte Carlo has been already created." << endl;
429   }  
430   
431   
432   
433   // Customization of Geant4 VMC
434   //
435   geant4->ProcessGeantCommand("/mcVerbose/all 1");  
436   geant4->ProcessGeantCommand("/mcVerbose/geometryManager 1");  
437   geant4->ProcessGeantCommand("/mcVerbose/opGeometryManager 1");  
438   geant4->ProcessGeantCommand("/mcTracking/loopVerbose 0");     
439   geant4->ProcessGeantCommand("/mcPhysics/rangeCuts 0.01 mm"); 
440   geant4->ProcessGeantCommand("/mcPhysics/selectOpProcess Scintillation");
441   geant4->ProcessGeantCommand("/mcPhysics/setOpProcessActivation false");
442   geant4->ProcessGeantCommand("/mcVerbose/composedPhysicsList 2");  
443
444
445
446
447
448
449  //
450   //=======================================================================
451   // ************* STEERING parameters FOR ALICE SIMULATION **************
452   // --- Specify event type to be tracked through the ALICE setup
453   // --- All positions are in cm, angles in degrees, and P and E in GeV
454
455
456     gMC->SetProcess("DCAY",1);
457     gMC->SetProcess("PAIR",1);
458     gMC->SetProcess("COMP",1);
459     gMC->SetProcess("PHOT",1);
460     gMC->SetProcess("PFIS",0);
461     gMC->SetProcess("DRAY",0);
462     gMC->SetProcess("ANNI",1);
463     gMC->SetProcess("BREM",1);
464     gMC->SetProcess("MUNU",1);
465     gMC->SetProcess("CKOV",1);
466     gMC->SetProcess("HADR",1);
467     gMC->SetProcess("LOSS",2);
468     gMC->SetProcess("MULS",1);
469     gMC->SetProcess("RAYL",1);
470
471     Float_t cut = 1.e-3;        // 1MeV cut by default
472     Float_t tofmax = 1.e10;
473
474     gMC->SetCut("CUTGAM", cut);
475     gMC->SetCut("CUTELE", cut);
476     gMC->SetCut("CUTNEU", cut);
477     gMC->SetCut("CUTHAD", cut);
478     gMC->SetCut("CUTMUO", cut);
479     gMC->SetCut("BCUTE",  cut); 
480     gMC->SetCut("BCUTM",  cut); 
481     gMC->SetCut("DCUTE",  cut); 
482     gMC->SetCut("DCUTM",  cut); 
483     gMC->SetCut("PPCUTM", cut);
484     gMC->SetCut("TOFMAX", tofmax); 
485
486
487
488
489   //======================//
490   // Set External decayer //
491   //======================//
492   TVirtualMCDecayer* decayer = new AliDecayerPythia();
493   decayer->SetForceDecay(kAll);
494   decayer->Init();
495   gMC->SetExternalDecayer(decayer);
496
497   //=========================//
498   // Generator Configuration //
499   //=========================//
500   AliGenerator* gener = 0x0;
501   
502   if (proc == kPythia6) {
503       gener = MbPythia();
504   } else if (proc == kPythia6D6T) {
505       gener = MbPythiaTuneD6T();
506   } else if (proc == kPythia6D6T_Flat) {
507       gener = MbPythiaTuneD6T_Flat();
508   } else if (proc == kPythia6ATLAS) {
509       gener = MbPythiaTuneATLAS();
510   } else if (proc == kPythia6ATLAS_Flat) {
511       gener = MbPythiaTuneATLAS_Flat();
512   } else if (proc == kPhojet) {
513       gener = MbPhojet();
514   } else if (proc == kPythia6_Perugia0) {
515       gener = MbPythiaTunePerugia0();
516   }
517
518   //
519   // PRIMARY VERTEX
520   //
521    Double_t xv = 0;
522    Double_t yv = 0;
523    Double_t zv = 0;
524
525    Double_t vtxPos[3];
526    Double_t vtxErr[3];
527    for (Int_t i = 0; i < 3; i++) {
528        vtxPos[i] = 0.;
529        vtxErr[i] = 0.;
530    }
531
532    if(AliCDBManager::Instance()->IsDefaultStorageSet()){
533       AliCDBEntry* entry = AliCDBManager::Instance()->Get("GRP/Calib/MeanVertexSPD");
534       AliESDVertex* vertex = dynamic_cast<AliESDVertex*> (entry->GetObject());
535       vertex->GetXYZ(vtxPos);
536       vertex->GetSigmaXYZ(vtxErr);
537     }
538
539     printf("Vertex position from OCDB entry: x = %13.3f, y = %13.3f, z = %13.3f (sigma = %13.3f)\n",
540     vtxPos[0], vtxPos[1], vtxPos[2], vtxErr[2]);
541
542      gener->SetOrigin(vtxPos[0], vtxPos[1], vtxPos[2]);   // vertex position
543   //
544   //
545   // Size of the interaction diamond
546   // Longitudinal
547      sigmaz = vtxErr[2];
548 //
549   // Transverse
550   Float_t betast  = 10;                 // beta* [m]
551   Float_t eps     = 2.5e-6;            // emittance [m]
552   Float_t gamma   = energy / 2.0 / 0.938272;  // relativistic gamma [1]
553   Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.;  // [cm]
554  
555   printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
556
557   gener->SetSigma(sigmaxy, sigmaxy, sigmaz);      // Sigma in (X,Y,Z) (cm) on IP position
558   gener->SetVertexSmear(kPerEvent);
559   gener->Init();
560
561
562
563
564
565
566
567
568
569
570
571
572
573 }
574 //
575 //           PYTHIA
576 //
577
578 AliGenerator* MbPythia()
579 {
580       comment = comment.Append(" pp: Pythia low-pt");
581 //
582 //    Pythia
583       AliGenPythia* pythia = new AliGenPythia(-1); 
584       pythia->SetMomentumRange(0, 999999.);
585       pythia->SetThetaRange(0., 180.);
586       pythia->SetYRange(-12.,12.);
587       pythia->SetPtRange(0,1000.);
588       pythia->SetProcess(kPyMb);
589       pythia->SetEnergyCMS(energy);
590       
591       return pythia;
592 }
593
594 AliGenerator* MbPythiaTuneD6T()
595 {
596       comment = comment.Append(" pp: Pythia low-pt");
597 //
598 //    Pythia
599       AliGenPythia* pythia = new AliGenPythia(-1); 
600       pythia->SetMomentumRange(0, 999999.);
601       pythia->SetThetaRange(0., 180.);
602       pythia->SetYRange(-12.,12.);
603       pythia->SetPtRange(0,1000.);
604       pythia->SetProcess(kPyMb);
605       pythia->SetEnergyCMS(energy);
606 //    Tune
607 //    109     D6T : Rick Field's CDF Tune D6T (NB: needs CTEQ6L pdfs externally)
608       pythia->SetTune(109); // F I X 
609       pythia->SetStrucFunc(kCTEQ6l);
610 //
611       return pythia;
612 }
613
614 AliGenerator* MbPythiaTunePerugia0()
615 {
616       comment = comment.Append(" pp: Pythia Minimum Bias Tune Perugia0");
617 //
618 //    Pythia
619       AliGenPythia* pythia = new AliGenPythia(-1); 
620       pythia->SetMomentumRange(0, 999999.);
621       pythia->SetThetaRange(0., 180.);
622       pythia->SetYRange(-12.,12.);
623       pythia->SetPtRange(0, 1000.);
624       pythia->SetProcess(kPyMb);
625       pythia->SetEnergyCMS(energy);
626 //    Tune
627       pythia->SetTune(320);
628 //
629       return pythia;
630 }
631
632 AliGenerator* MbPythiaTuneATLAS()
633 {
634       comment = comment.Append(" pp: Pythia low-pt");
635 //
636 //    Pythia
637       AliGenPythia* pythia = new AliGenPythia(-1); 
638       pythia->SetMomentumRange(0, 999999.);
639       pythia->SetThetaRange(0., 180.);
640       pythia->SetYRange(-12.,12.);
641       pythia->SetPtRange(0,1000.);
642       pythia->SetProcess(kPyMb);
643       pythia->SetEnergyCMS(energy);
644 //    Tune
645 //    C   306 ATLAS-CSC: Arthur Moraes' (new) ATLAS tune (needs CTEQ6L externally)
646       pythia->SetTune(306);
647       pythia->SetStrucFunc(kCTEQ6l);
648 //
649       return pythia;
650 }
651
652 AliGenerator* MbPythiaTuneATLAS_Flat()
653 {
654       AliGenPythia* pythia = MbPythiaTuneATLAS();
655       
656       comment = comment.Append("; flat multiplicity distribution");
657       
658       // set high multiplicity trigger
659       // this weight achieves a flat multiplicity distribution
660       TH1 *weight = new TH1D("weight","weight",201,-0.5,200.5);
661       weight->SetBinContent(1,5.49443);
662       weight->SetBinContent(2,8.770816);
663       weight->SetBinContent(6,0.4568624);
664       weight->SetBinContent(7,0.2919915);
665       weight->SetBinContent(8,0.6674189);
666       weight->SetBinContent(9,0.364737);
667       weight->SetBinContent(10,0.8818444);
668       weight->SetBinContent(11,0.531885);
669       weight->SetBinContent(12,1.035197);
670       weight->SetBinContent(13,0.9394057);
671       weight->SetBinContent(14,0.9643193);
672       weight->SetBinContent(15,0.94543);
673       weight->SetBinContent(16,0.9426507);
674       weight->SetBinContent(17,0.9423649);
675       weight->SetBinContent(18,0.789456);
676       weight->SetBinContent(19,1.149026);
677       weight->SetBinContent(20,1.100491);
678       weight->SetBinContent(21,0.6350525);
679       weight->SetBinContent(22,1.351941);
680       weight->SetBinContent(23,0.03233504);
681       weight->SetBinContent(24,0.9574557);
682       weight->SetBinContent(25,0.868133);
683       weight->SetBinContent(26,1.030998);
684       weight->SetBinContent(27,1.08897);
685       weight->SetBinContent(28,1.251382);
686       weight->SetBinContent(29,0.1391099);
687       weight->SetBinContent(30,1.192876);
688       weight->SetBinContent(31,0.448944);
689       weight->SetBinContent(32,1);
690       weight->SetBinContent(33,1);
691       weight->SetBinContent(34,1);
692       weight->SetBinContent(35,1);
693       weight->SetBinContent(36,0.9999997);
694       weight->SetBinContent(37,0.9999997);
695       weight->SetBinContent(38,0.9999996);
696       weight->SetBinContent(39,0.9999996);
697       weight->SetBinContent(40,0.9999995);
698       weight->SetBinContent(41,0.9999993);
699       weight->SetBinContent(42,1);
700       weight->SetBinContent(43,1);
701       weight->SetBinContent(44,1);
702       weight->SetBinContent(45,1);
703       weight->SetBinContent(46,1);
704       weight->SetBinContent(47,0.9999999);
705       weight->SetBinContent(48,0.9999998);
706       weight->SetBinContent(49,0.9999998);
707       weight->SetBinContent(50,0.9999999);
708       weight->SetBinContent(51,0.9999999);
709       weight->SetBinContent(52,0.9999999);
710       weight->SetBinContent(53,0.9999999);
711       weight->SetBinContent(54,0.9999998);
712       weight->SetBinContent(55,0.9999998);
713       weight->SetBinContent(56,0.9999998);
714       weight->SetBinContent(57,0.9999997);
715       weight->SetBinContent(58,0.9999996);
716       weight->SetBinContent(59,0.9999995);
717       weight->SetBinContent(60,1);
718       weight->SetBinContent(61,1);
719       weight->SetBinContent(62,1);
720       weight->SetBinContent(63,1);
721       weight->SetBinContent(64,1);
722       weight->SetBinContent(65,0.9999999);
723       weight->SetBinContent(66,0.9999998);
724       weight->SetBinContent(67,0.9999998);
725       weight->SetBinContent(68,0.9999999);
726       weight->SetBinContent(69,1);
727       weight->SetBinContent(70,1);
728       weight->SetBinContent(71,0.9999997);
729       weight->SetBinContent(72,0.9999995);
730       weight->SetBinContent(73,0.9999994);
731       weight->SetBinContent(74,1);
732       weight->SetBinContent(75,1);
733       weight->SetBinContent(76,1);
734       weight->SetBinContent(77,1);
735       weight->SetBinContent(78,0.9999999);
736       weight->SetBinContent(79,1);
737       weight->SetBinContent(80,1);
738       weight->SetEntries(526);
739         
740       Int_t limit = weight->GetRandom();
741       pythia->SetTriggerChargedMultiplicity(limit, 1.4);
742       
743       comment = comment.Append(Form("; multiplicity threshold set to %d in |eta| < 1.4", limit));
744
745       return pythia;
746 }
747 AliGenerator* MbPythiaTuneD6T_Flat()
748 {
749       AliGenPythia* pythia = MbPythiaTuneD6T();
750       
751       comment = comment.Append("; flat multiplicity distribution");
752       
753       // set high multiplicity trigger
754       // this weight achieves a flat multiplicity distribution
755       TH1 *weight = new TH1D("weight","weight",201,-0.5,200.5);
756       weight->SetBinContent(1,5.49443);
757       weight->SetBinContent(2,8.770816);
758       weight->SetBinContent(6,0.4568624);
759       weight->SetBinContent(7,0.2919915);
760       weight->SetBinContent(8,0.6674189);
761       weight->SetBinContent(9,0.364737);
762       weight->SetBinContent(10,0.8818444);
763       weight->SetBinContent(11,0.531885);
764       weight->SetBinContent(12,1.035197);
765       weight->SetBinContent(13,0.9394057);
766       weight->SetBinContent(14,0.9643193);
767       weight->SetBinContent(15,0.94543);
768       weight->SetBinContent(16,0.9426507);
769       weight->SetBinContent(17,0.9423649);
770       weight->SetBinContent(18,0.789456);
771       weight->SetBinContent(19,1.149026);
772       weight->SetBinContent(20,1.100491);
773       weight->SetBinContent(21,0.6350525);
774       weight->SetBinContent(22,1.351941);
775       weight->SetBinContent(23,0.03233504);
776       weight->SetBinContent(24,0.9574557);
777       weight->SetBinContent(25,0.868133);
778       weight->SetBinContent(26,1.030998);
779       weight->SetBinContent(27,1.08897);
780       weight->SetBinContent(28,1.251382);
781       weight->SetBinContent(29,0.1391099);
782       weight->SetBinContent(30,1.192876);
783       weight->SetBinContent(31,0.448944);
784       weight->SetBinContent(32,1);
785       weight->SetBinContent(33,1);
786       weight->SetBinContent(34,1);
787       weight->SetBinContent(35,1);
788       weight->SetBinContent(36,0.9999997);
789       weight->SetBinContent(37,0.9999997);
790       weight->SetBinContent(38,0.9999996);
791       weight->SetBinContent(39,0.9999996);
792       weight->SetBinContent(40,0.9999995);
793       weight->SetBinContent(41,0.9999993);
794       weight->SetBinContent(42,1);
795       weight->SetBinContent(43,1);
796       weight->SetBinContent(44,1);
797       weight->SetBinContent(45,1);
798       weight->SetBinContent(46,1);
799       weight->SetBinContent(47,0.9999999);
800       weight->SetBinContent(48,0.9999998);
801       weight->SetBinContent(49,0.9999998);
802       weight->SetBinContent(50,0.9999999);
803       weight->SetBinContent(51,0.9999999);
804       weight->SetBinContent(52,0.9999999);
805       weight->SetBinContent(53,0.9999999);
806       weight->SetBinContent(54,0.9999998);
807       weight->SetBinContent(55,0.9999998);
808       weight->SetBinContent(56,0.9999998);
809       weight->SetBinContent(57,0.9999997);
810       weight->SetBinContent(58,0.9999996);
811       weight->SetBinContent(59,0.9999995);
812       weight->SetBinContent(60,1);
813       weight->SetBinContent(61,1);
814       weight->SetBinContent(62,1);
815       weight->SetBinContent(63,1);
816       weight->SetBinContent(64,1);
817       weight->SetBinContent(65,0.9999999);
818       weight->SetBinContent(66,0.9999998);
819       weight->SetBinContent(67,0.9999998);
820       weight->SetBinContent(68,0.9999999);
821       weight->SetBinContent(69,1);
822       weight->SetBinContent(70,1);
823       weight->SetBinContent(71,0.9999997);
824       weight->SetBinContent(72,0.9999995);
825       weight->SetBinContent(73,0.9999994);
826       weight->SetBinContent(74,1);
827       weight->SetBinContent(75,1);
828       weight->SetBinContent(76,1);
829       weight->SetBinContent(77,1);
830       weight->SetBinContent(78,0.9999999);
831       weight->SetBinContent(79,1);
832       weight->SetBinContent(80,1);
833       weight->SetEntries(526);
834         
835       Int_t limit = weight->GetRandom();
836       pythia->SetTriggerChargedMultiplicity(limit, 1.4);
837       
838       comment = comment.Append(Form("; multiplicity threshold set to %d in |eta| < 1.4", limit));
839
840       return pythia;
841 }
842
843 AliGenerator* MbPhojet()
844 {
845       comment = comment.Append(" pp: Pythia low-pt");
846 //
847 //    DPMJET
848 #if defined(__CINT__)
849   gSystem->Load("libdpmjet");      // Parton density functions
850   gSystem->Load("libTDPMjet");      // Parton density functions
851 #endif
852       AliGenDPMjet* dpmjet = new AliGenDPMjet(-1); 
853       dpmjet->SetMomentumRange(0, 999999.);
854       dpmjet->SetThetaRange(0., 180.);
855       dpmjet->SetYRange(-12.,12.);
856       dpmjet->SetPtRange(0,1000.);
857       dpmjet->SetProcess(kDpmMb);
858       dpmjet->SetEnergyCMS(energy);
859
860       return dpmjet;
861 }
862
863 void ProcessEnvironmentVars()
864 {
865     // Run type
866     if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
867       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
868         if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
869           proc = (PDC06Proc_t)iRun;
870           cout<<"Run type set to "<<pprRunName[iRun]<<endl;
871         }
872       }
873     }
874
875     // Field
876     if (gSystem->Getenv("CONFIG_FIELD")) {
877       for (Int_t iField = 0; iField < kFieldMax; iField++) {
878         if (strcmp(gSystem->Getenv("CONFIG_FIELD"), pprField[iField])==0) {
879           mag = (Mag_t)iField;
880           cout<<"Field set to "<<pprField[iField]<<endl;
881         }
882       }
883     }
884
885     // Energy
886     if (gSystem->Getenv("CONFIG_ENERGY")) {
887       energy = atoi(gSystem->Getenv("CONFIG_ENERGY"));
888       cout<<"Energy set to "<<energy<<" GeV"<<endl;
889     }
890
891     // Random Number seed
892     if (gSystem->Getenv("CONFIG_SEED")) {
893       seed = atoi(gSystem->Getenv("CONFIG_SEED"));
894     }
895
896     // Physics lists
897     if (gSystem->Getenv("CONFIG_PHYSICSLIST")) {
898       for (Int_t iList = 0; iList < kListMax; iList++) {
899         if (strcmp(gSystem->Getenv("CONFIG_PHYSICSLIST"), physicsListName[iList])==0) {
900           physicslist = (PhysicsList_t)iList;
901           cout<<"Physics list set to "<< physicsListName[iList]<<endl;
902         }
903       }
904     }
905
906 }