bdac0d367c8bc2be766b1ae8cabd3d54f09f28dc
[u/mrichter/AliRoot.git] / macros / ConfigPPR.C
1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <Riostream.h>
3 #include <TRandom.h>
4 #include <TSystem.h>
5 #include <TVirtualMC.h>
6 #include <TGeant3.h>
7 #include "STEER/AliRunLoader.h"
8 #include "STEER/AliRun.h"
9 #include "STEER/AliConfig.h"
10 #include "STEER/AliGenerator.h"
11 #include "PYTHIA6/AliDecayerPythia.h"
12 #include "EVGEN/AliGenHIJINGpara.h"
13 #include "THijing/AliGenHijing.h"
14 #include "EVGEN/AliGenCocktail.h"
15 #include "EVGEN/AliGenSlowNucleons.h"
16 #include "EVGEN/AliSlowNucleonModelExp.h"
17 #include "EVGEN/AliGenParam.h"
18 #include "EVGEN/AliGenGSIlib.h"
19 #include "PYTHIA6/AliGenPythia.h"
20 #include "STEER/AliMagFMaps.h"
21 #include "STRUCT/AliBODY.h"
22 #include "STRUCT/AliMAG.h"
23 #include "STRUCT/AliABSOv0.h"
24 #include "STRUCT/AliDIPOv2.h"
25 #include "STRUCT/AliHALL.h"
26 #include "STRUCT/AliFRAMEv2.h"
27 #include "STRUCT/AliSHILv2.h"
28 #include "STRUCT/AliPIPEv0.h"
29 #include "ITS/AliITSvPPRasymmFMD.h"
30 #include "TPC/AliTPCv2.h"
31 #include "TOF/AliTOFv4T0.h"
32 #include "RICH/AliRICHv1.h"
33 #include "ZDC/AliZDCv2.h"
34 #include "TRD/AliTRDv1.h"
35 #include "FMD/AliFMDv1.h"
36 #include "MUON/AliMUONv1.h"
37 #include "MUON/AliMUONSt1GeometryBuilder.h"
38 #include "MUON/AliMUONSt2GeometryBuilder.h"
39 #include "MUON/AliMUONSlatGeometryBuilder.h"
40 #include "MUON/AliMUONTriggerGeometryBuilder.h"
41 #include "PHOS/AliPHOSv1.h"
42 #include "PMD/AliPMDv1.h"
43 #include "START/AliSTARTv1.h"
44 #include "EMCAL/AliEMCALv1.h"
45 #include "CRT/AliCRTv0.h"
46 #include "VZERO/AliVZEROv3.h"
47 #endif
48
49 enum PprRun_t 
50 {
51     test50,
52     kParam_8000,   kParam_4000,  kParam_2000, 
53     kHijing_cent1, kHijing_cent2, 
54     kHijing_per1,  kHijing_per2, kHijing_per3, kHijing_per4,  kHijing_per5,
55     kHijing_jj25,  kHijing_jj50, kHijing_jj75, kHijing_jj100, kHijing_jj200, 
56     kHijing_gj25,  kHijing_gj50, kHijing_gj75, kHijing_gj100, kHijing_gj200,
57     kHijing_pA, kPythia6, kPythia6Jets, kD0PbPb5500, kD_TRD, kB_TRD, kJpsi_TRD,
58     kU_TRD, kPyJJ, kPyGJ, kRunMax
59 };
60
61 const char* pprRunName[kRunMax] = {
62     "test50",
63     "kParam_8000",   "kParam_4000",  "kParam_2000", 
64     "kHijing_cent1", "kHijing_cent2", 
65     "kHijing_per1",  "kHijing_per2", "kHijing_per3", "kHijing_per4",  
66     "kHijing_per5",
67     "kHijing_jj25",  "kHijing_jj50", "kHijing_jj75", "kHijing_jj100", 
68     "kHijing_jj200", 
69     "kHijing_gj25",  "kHijing_gj50", "kHijing_gj75", "kHijing_gj100", 
70     "kHijing_gj200",
71     "kHijing_pA", "kPythia6", "kPythia6Jets", "kD0PbPb5500", "kD_TRD", 
72     "kB_TRD", "kJpsi_TRD",
73     "kU_TRD", "kPyJJ", "kPyGJ"
74 };
75
76 enum PprGeo_t 
77 {
78     kHoles, kNoHoles
79 };
80
81 enum PprRad_t
82 {
83     kGluonRadiation, kNoGluonRadiation
84 };
85
86 enum PprMag_t
87 {
88     k2kG, k4kG, k5kG
89 };
90
91
92 // This part for configuration    
93 //static PprRun_t srun = test50;
94 static PprRun_t srun = kPythia6;
95 static PprGeo_t sgeo = kHoles;
96 static PprRad_t srad = kGluonRadiation;
97 static PprMag_t smag = k5kG;
98 static Int_t    sseed = 12345; //Set 0 to use the current time
99
100 // Comment line 
101 static TString  comment;
102
103 // Functions
104 Float_t EtaToTheta(Float_t arg);
105 AliGenerator* GeneratorFactory(PprRun_t srun);
106 AliGenHijing* HijingStandard();
107 void ProcessEnvironmentVars();
108
109 void Config()
110 {
111     // ThetaRange is (0., 180.). It was (0.28,179.72) 7/12/00 09:00
112     // Theta range given through pseudorapidity limits 22/6/2001
113
114     // Get settings from environment variables
115     ProcessEnvironmentVars();
116
117     // Set Random Number seed
118     gRandom->SetSeed(sseed);
119     cout<<"Seed for random number generation= "<<gRandom->GetSeed()<<endl; 
120
121
122    // libraries required by geant321
123 #if defined(__CINT__)
124     gSystem->Load("libgeant321");
125 #endif
126
127     new     TGeant3("C++ Interface to Geant3");
128
129     AliRunLoader* rl=0x0;
130
131     cout<<"Config.C: Creating Run Loader ..."<<endl;
132     rl = AliRunLoader::Open("galice.root",
133                             AliConfig::fgkDefaultEventFolderName,
134                             "recreate");
135     if (rl == 0x0)
136       {
137         gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
138         return;
139       }
140     rl->SetCompressionLevel(2);
141     rl->SetNumberOfEventsPerFile(3);
142     gAlice->SetRunLoader(rl);
143
144     //
145     // Set External decayer
146     AliDecayer *decayer = new AliDecayerPythia();
147
148     decayer->SetForceDecay(kAll);
149     decayer->Init();
150     gMC->SetExternalDecayer(decayer);
151     //
152     //
153     //=======================================================================
154     //
155     //=======================================================================
156     // ************* STEERING parameters FOR ALICE SIMULATION **************
157     // --- Specify event type to be tracked through the ALICE setup
158     // --- All positions are in cm, angles in degrees, and P and E in GeV
159
160     gMC->SetProcess("DCAY",1);
161     gMC->SetProcess("PAIR",1);
162     gMC->SetProcess("COMP",1);
163     gMC->SetProcess("PHOT",1);
164     gMC->SetProcess("PFIS",0);
165     gMC->SetProcess("DRAY",0);
166     gMC->SetProcess("ANNI",1);
167     gMC->SetProcess("BREM",1);
168     gMC->SetProcess("MUNU",1);
169     gMC->SetProcess("CKOV",1);
170     gMC->SetProcess("HADR",1);
171     gMC->SetProcess("LOSS",2);
172     gMC->SetProcess("MULS",1);
173     gMC->SetProcess("RAYL",1);
174
175     Float_t cut = 1.e-3;        // 1MeV cut by default
176     Float_t tofmax = 1.e10;
177
178     gMC->SetCut("CUTGAM", cut);
179     gMC->SetCut("CUTELE", cut);
180     gMC->SetCut("CUTNEU", cut);
181     gMC->SetCut("CUTHAD", cut);
182     gMC->SetCut("CUTMUO", cut);
183     gMC->SetCut("BCUTE",  cut); 
184     gMC->SetCut("BCUTM",  cut); 
185     gMC->SetCut("DCUTE",  cut); 
186     gMC->SetCut("DCUTM",  cut); 
187     gMC->SetCut("PPCUTM", cut);
188     gMC->SetCut("TOFMAX", tofmax); 
189
190
191 // Generator Configuration
192     gAlice->SetDebug(0);
193     AliGenerator* gener = GeneratorFactory(srun);
194     gener->SetOrigin(0, 0, 0);    // vertex position
195     gener->SetSigma(0, 0, 5.3);   // Sigma in (X,Y,Z) (cm) on IP position
196     gener->SetCutVertexZ(1.);     // Truncate at 1 sigma
197     gener->SetVertexSmear(kPerEvent); 
198     gener->SetTrackingFlag(1);
199     gener->Init();
200     
201     if (smag == k2kG) {
202         comment = comment.Append(" | L3 field 0.2 T");
203     } else if (smag == k4kG) {
204         comment = comment.Append(" | L3 field 0.4 T");
205     } else if (smag == k5kG) {
206         comment = comment.Append(" | L3 field 0.5 T");
207     }
208     
209     
210     if (srad == kGluonRadiation)
211     {
212         comment = comment.Append(" | Gluon Radiation On");
213         
214     } else {
215         comment = comment.Append(" | Gluon Radiation Off");
216     }
217
218     if (sgeo == kHoles)
219     {
220         comment = comment.Append(" | Holes for PHOS/RICH");
221         
222     } else {
223         comment = comment.Append(" | No holes for PHOS/RICH");
224     }
225
226     printf("\n \n Comment: %s \n \n", comment.Data());
227     
228     
229 // Field (L3 0.4 T)
230     AliMagFMaps* field = new AliMagFMaps("Maps","Maps", 2, 1., 10., smag);
231     field->SetL3ConstField(0); //Using const. field in the barrel
232     rl->CdGAFile();
233     gAlice->SetField(field);    
234 //
235     Int_t   iABSO   = 1;
236     Int_t   iDIPO   = 1;
237     Int_t   iFMD    = 1;
238     Int_t   iFRAME  = 1;
239     Int_t   iHALL   = 1;
240     Int_t   iITS    = 1;
241     Int_t   iMAG    = 1;
242     Int_t   iMUON   = 1;
243     Int_t   iPHOS   = 1;
244     Int_t   iPIPE   = 1;
245     Int_t   iPMD    = 1;
246     Int_t   iRICH   = 1;
247     Int_t   iSHIL   = 1;
248     Int_t   iSTART  = 1;
249     Int_t   iTOF    = 1;
250     Int_t   iTPC    = 1;
251     Int_t   iTRD    = 1;
252     Int_t   iZDC    = 1;
253     Int_t   iEMCAL  = 1;
254     Int_t   iVZERO  = 1;
255     Int_t   iCRT    = 0;
256
257     //=================== Alice BODY parameters =============================
258     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
259
260
261     if (iMAG)
262     {
263         //=================== MAG parameters ============================
264         // --- Start with Magnet since detector layouts may be depending ---
265         // --- on the selected Magnet dimensions ---
266         AliMAG *MAG = new AliMAG("MAG", "Magnet");
267     }
268
269
270     if (iABSO)
271     {
272         //=================== ABSO parameters ============================
273         AliABSO *ABSO = new AliABSOv0("ABSO", "Muon Absorber");
274     }
275
276     if (iDIPO)
277     {
278         //=================== DIPO parameters ============================
279
280         AliDIPO *DIPO = new AliDIPOv2("DIPO", "Dipole version 2");
281     }
282
283     if (iHALL)
284     {
285         //=================== HALL parameters ============================
286
287         AliHALL *HALL = new AliHALL("HALL", "Alice Hall");
288     }
289
290
291     if (iFRAME)
292     {
293         //=================== FRAME parameters ============================
294
295         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
296         if (sgeo == kHoles) {
297             FRAME->SetHoles(1);
298         } else {
299             FRAME->SetHoles(0);
300         }
301     }
302
303     if (iSHIL)
304     {
305         //=================== SHIL parameters ============================
306
307         AliSHIL *SHIL = new AliSHILv2("SHIL", "Shielding Version 2");
308     }
309
310
311     if (iPIPE)
312     {
313         //=================== PIPE parameters ============================
314
315         AliPIPE *PIPE = new AliPIPEv0("PIPE", "Beam Pipe");
316     }
317  
318     if(iITS) {
319
320     //=================== ITS parameters ============================
321     //
322     // As the innermost detector in ALICE, the Inner Tracking System "impacts" on
323     // almost all other detectors. This involves the fact that the ITS geometry
324     // still has several options to be followed in parallel in order to determine
325     // the best set-up which minimizes the induced background. All the geometries
326     // available to date are described in the following. Read carefully the comments
327     // and use the default version (the only one uncommented) unless you are making
328     // comparisons and you know what you are doing. In this case just uncomment the
329     // ITS geometry you want to use and run Aliroot.
330     //
331     // Detailed geometries:         
332     //
333     //
334     //AliITS *ITS  = new AliITSv5symm("ITS","Updated ITS TDR detailed version with symmetric services");
335     //
336     //AliITS *ITS  = new AliITSv5asymm("ITS","Updates ITS TDR detailed version with asymmetric services");
337     //
338         AliITSvPPRasymmFMD *ITS  = new AliITSvPPRasymmFMD("ITS","New ITS PPR detailed version with asymmetric services");
339         ITS->SetMinorVersion(2);  // don't touch this parameter if you're not an ITS developer
340         ITS->SetReadDet(kTRUE);   // don't touch this parameter if you're not an ITS developer
341     //    ITS->SetWriteDet("$ALICE_ROOT/ITS/ITSgeometry_vPPRasymm2.det");  // don't touch this parameter if you're not an ITS developer
342         ITS->SetThicknessDet1(200.);   // detector thickness on layer 1 must be in the range [100,300]
343         ITS->SetThicknessDet2(200.);   // detector thickness on layer 2 must be in the range [100,300]
344         ITS->SetThicknessChip1(200.);  // chip thickness on layer 1 must be in the range [150,300]
345         ITS->SetThicknessChip2(200.);  // chip thickness on layer 2 must be in the range [150,300]
346         ITS->SetRails(0);            // 1 --> rails in ; 0 --> rails out
347         ITS->SetCoolingFluid(1);   // 1 --> water ; 0 --> freon
348
349     // Coarse geometries (warning: no hits are produced with these coarse geometries and they unuseful 
350     // for reconstruction !):
351     //                                                     
352     //
353     //AliITSvPPRcoarseasymm *ITS  = new AliITSvPPRcoarseasymm("ITS","New ITS PPR coarse version with asymmetric services");
354     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
355     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
356     //
357     //AliITS *ITS  = new AliITSvPPRcoarsesymm("ITS","New ITS PPR coarse version with symmetric services");
358     //ITS->SetRails(0);                // 1 --> rails in ; 0 --> rails out
359     //ITS->SetSupportMaterial(0);      // 0 --> Copper ; 1 --> Aluminum ; 2 --> Carbon
360     //                      
361     //
362     //
363     // Geant3 <-> EUCLID conversion
364     // ============================
365     //
366     // SetEUCLID is a flag to output (=1) or not to output (=0) both geometry and
367     // media to two ASCII files (called by default ITSgeometry.euc and
368     // ITSgeometry.tme) in a format understandable to the CAD system EUCLID.
369     // The default (=0) means that you dont want to use this facility.
370     //
371         ITS->SetEUCLID(0);  
372     }
373
374     if (iTPC)
375     {
376         //============================ TPC parameters ================================
377         // --- This allows the user to specify sectors for the SLOW (TPC geometry 2)
378         // --- Simulator. SecAL (SecAU) <0 means that ALL lower (upper)
379         // --- sectors are specified, any value other than that requires at least one 
380         // --- sector (lower or upper)to be specified!
381         // --- Reminder: sectors 1-24 are lower sectors (1-12 -> z>0, 13-24 -> z<0)
382         // ---           sectors 25-72 are the upper ones (25-48 -> z>0, 49-72 -> z<0)
383         // --- SecLows - number of lower sectors specified (up to 6)
384         // --- SecUps - number of upper sectors specified (up to 12)
385         // --- Sens - sensitive strips for the Slow Simulator !!!
386         // --- This does NOT work if all S or L-sectors are specified, i.e.
387         // --- if SecAL or SecAU < 0
388         //
389         //
390         //-----------------------------------------------------------------------------
391
392         //  gROOT->LoadMacro("SetTPCParam.C");
393         //  AliTPCParam *param = SetTPCParam();
394         AliTPC *TPC = new AliTPCv2("TPC", "Default");
395
396         // All sectors included 
397         TPC->SetSecAL(-1);
398         TPC->SetSecAU(-1);
399
400     }
401
402
403     if (iTOF) {
404         //=================== TOF parameters ============================
405         AliTOF *TOF = new AliTOFv4T0("TOF", "normal TOF");
406     }
407
408
409     if (iRICH)
410     {
411         //=================== RICH parameters ===========================
412         AliRICH *RICH = new AliRICHv1("RICH", "normal RICH");
413
414     }
415
416
417     if (iZDC)
418     {
419         //=================== ZDC parameters ============================
420
421         AliZDC *ZDC = new AliZDCv2("ZDC", "normal ZDC");
422     }
423
424     if (iTRD)
425     {
426         //=================== TRD parameters ============================
427
428         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
429
430         // Select the gas mixture (0: 97% Xe + 3% isobutane, 1: 90% Xe + 10% CO2)
431         TRD->SetGasMix(1);
432         if (sgeo == kHoles) {
433             // With hole in front of PHOS
434             TRD->SetPHOShole();
435             // With hole in front of RICH
436             TRD->SetRICHhole();
437         }
438             // Switch on TR
439             AliTRDsim *TRDsim = TRD->CreateTR();
440     }
441
442     if (iFMD)
443     {
444         //=================== FMD parameters ============================
445         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
446    }
447
448     if (iMUON)
449     {
450         //=================== MUON parameters ===========================
451
452         AliMUON *MUON = new AliMUONv1("MUON", "default");
453         MUON->AddGeometryBuilder(new AliMUONSt1GeometryBuilder(MUON));
454         MUON->AddGeometryBuilder(new AliMUONSt2GeometryBuilder(MUON));
455         MUON->AddGeometryBuilder(new AliMUONSlatGeometryBuilder(MUON));
456         MUON->AddGeometryBuilder(new AliMUONTriggerGeometryBuilder(MUON));
457     }
458     //=================== PHOS parameters ===========================
459
460     if (iPHOS)
461     {
462         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
463     }
464
465
466     if (iPMD)
467     {
468         //=================== PMD parameters ============================
469         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
470     }
471
472     if (iSTART)
473     {
474         //=================== START parameters ============================
475         AliSTART *START = new AliSTARTv1("START", "START Detector");
476     }
477
478     if (iEMCAL)
479     {
480         //=================== EMCAL parameters ============================
481         AliEMCAL *EMCAL = new AliEMCALv1("EMCAL", "EMCAL_55_25");
482     }
483
484      if (iCRT)
485     {
486         //=================== CRT parameters ============================
487         AliCRT *CRT = new AliCRTv0("CRT", "normal ACORDE");
488     }
489
490      if (iVZERO)
491     {
492         //=================== CRT parameters ============================
493         AliVZERO *VZERO = new AliVZEROv3("VZERO", "normal VZERO");
494     }
495  
496              
497 }
498
499 Float_t EtaToTheta(Float_t arg){
500   return (180./TMath::Pi())*2.*atan(exp(-arg));
501 }
502
503
504
505 AliGenerator* GeneratorFactory(PprRun_t srun) {
506     Int_t isw = 3;
507     if (srad == kNoGluonRadiation) isw = 0;
508     
509
510     AliGenerator * gGener = 0x0;
511     switch (srun) {
512     case test50:
513       {
514         comment = comment.Append(":HIJINGparam test 50 particles");
515         AliGenHIJINGpara *gener = new AliGenHIJINGpara(50);
516         gener->SetMomentumRange(0, 999999.);
517         gener->SetPhiRange(0., 360.);
518         // Set pseudorapidity range from -8 to 8.
519         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
520         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
521         gener->SetThetaRange(thmin,thmax);
522         gGener=gener;
523       }
524       break;
525     case kParam_8000:
526       {
527         comment = comment.Append(":HIJINGparam N=8000");
528         AliGenHIJINGpara *gener = new AliGenHIJINGpara(86030);
529         gener->SetMomentumRange(0, 999999.);
530         gener->SetPhiRange(0., 360.);
531         // Set pseudorapidity range from -8 to 8.
532         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
533         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
534         gener->SetThetaRange(thmin,thmax);
535         gGener=gener;
536       }
537       break;
538     case kParam_4000:
539       {
540         comment = comment.Append("HIJINGparam N=4000");
541         AliGenHIJINGpara *gener = new AliGenHIJINGpara(43015);
542         gener->SetMomentumRange(0, 999999.);
543         gener->SetPhiRange(0., 360.);
544         // Set pseudorapidity range from -8 to 8.
545         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
546         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
547         gener->SetThetaRange(thmin,thmax);
548         gGener=gener;
549       }
550         break;
551     case kParam_2000:
552       {
553         comment = comment.Append("HIJINGparam N=2000");
554         AliGenHIJINGpara *gener = new AliGenHIJINGpara(21507);
555         gener->SetMomentumRange(0, 999999.);
556         gener->SetPhiRange(0., 360.);
557         // Set pseudorapidity range from -8 to 8.
558         Float_t thmin = EtaToTheta(8);   // theta min. <---> eta max
559         Float_t thmax = EtaToTheta(-8);  // theta max. <---> eta min 
560         gener->SetThetaRange(thmin,thmax);
561         gGener=gener;
562       }
563       break;
564 //
565 //  Hijing Central
566 //
567     case kHijing_cent1:
568       {
569         comment = comment.Append("HIJING cent1");
570         AliGenHijing *gener = HijingStandard();
571 // impact parameter range
572         gener->SetImpactParameterRange(0., 5.);
573         gGener=gener;
574       }
575       break;
576     case kHijing_cent2:
577       {
578         comment = comment.Append("HIJING cent2");
579         AliGenHijing *gener = HijingStandard();
580 // impact parameter range
581         gener->SetImpactParameterRange(0., 2.);
582         gGener=gener;
583       }
584       break;
585 //
586 // Hijing Peripheral 
587 //
588     case kHijing_per1:
589       {
590         comment = comment.Append("HIJING per1");
591         AliGenHijing *gener = HijingStandard();
592 // impact parameter range
593         gener->SetImpactParameterRange(5., 8.6);
594         gGener=gener;
595       }
596       break;
597     case kHijing_per2:
598       {
599         comment = comment.Append("HIJING per2");
600         AliGenHijing *gener = HijingStandard();
601 // impact parameter range
602         gener->SetImpactParameterRange(8.6, 11.2);
603         gGener=gener;
604       }
605       break;
606     case kHijing_per3:
607       {
608         comment = comment.Append("HIJING per3");
609         AliGenHijing *gener = HijingStandard();
610 // impact parameter range
611         gener->SetImpactParameterRange(11.2, 13.2);
612         gGener=gener;
613       }
614       break;
615     case kHijing_per4:
616       {
617         comment = comment.Append("HIJING per4");
618         AliGenHijing *gener = HijingStandard();
619 // impact parameter range
620         gener->SetImpactParameterRange(13.2, 15.);
621         gGener=gener;
622       }
623       break;
624     case kHijing_per5:
625       {
626         comment = comment.Append("HIJING per5");
627         AliGenHijing *gener = HijingStandard();
628 // impact parameter range
629         gener->SetImpactParameterRange(15., 100.);
630         gGener=gener;
631       }
632       break;
633 //
634 //  Jet-Jet
635 //
636     case kHijing_jj25:
637       {
638         comment = comment.Append("HIJING Jet 25 GeV");
639         AliGenHijing *gener = HijingStandard();
640 // impact parameter range
641         gener->SetImpactParameterRange(0., 5.);
642         // trigger
643         gener->SetTrigger(1);
644         gener->SetPtJet(25.);
645         gener->SetRadiation(isw);
646         gener->SetSimpleJets(!isw);
647         gener->SetJetEtaRange(-0.3,0.3);
648         gener->SetJetPhiRange(75., 165.);   
649         gGener=gener;
650       }
651       break;
652
653     case kHijing_jj50:
654       {
655         comment = comment.Append("HIJING Jet 50 GeV");
656         AliGenHijing *gener = HijingStandard();
657 // impact parameter range
658         gener->SetImpactParameterRange(0., 5.);
659         // trigger
660         gener->SetTrigger(1);
661         gener->SetPtJet(50.);
662         gener->SetRadiation(isw);
663         gener->SetSimpleJets(!isw);
664         gener->SetJetEtaRange(-0.3,0.3);
665         gener->SetJetPhiRange(75., 165.);   
666         gGener=gener;
667       }
668         break;
669
670     case kHijing_jj75:
671       {
672         comment = comment.Append("HIJING Jet 75 GeV");
673         AliGenHijing *gener = HijingStandard();
674 // impact parameter range
675         gener->SetImpactParameterRange(0., 5.);
676         // trigger
677         gener->SetTrigger(1);
678         gener->SetPtJet(75.);
679         gener->SetRadiation(isw);
680         gener->SetSimpleJets(!isw);
681         gener->SetJetEtaRange(-0.3,0.3);
682         gener->SetJetPhiRange(75., 165.);   
683         gGener=gener;
684       }
685       break;
686
687     case kHijing_jj100:
688       {
689         comment = comment.Append("HIJING Jet 100 GeV");
690         AliGenHijing *gener = HijingStandard();
691 // impact parameter range
692         gener->SetImpactParameterRange(0., 5.);
693         // trigger
694         gener->SetTrigger(1);
695         gener->SetPtJet(100.);
696         gener->SetRadiation(isw);
697         gener->SetSimpleJets(!isw);
698         gener->SetJetEtaRange(-0.3,0.3);
699         gener->SetJetPhiRange(75., 165.);   
700         gGener=gener;
701       }
702       break;
703
704     case kHijing_jj200:
705       {
706         comment = comment.Append("HIJING Jet 200 GeV");
707         AliGenHijing *gener = HijingStandard();
708 // impact parameter range
709         gener->SetImpactParameterRange(0., 5.);
710         // trigger
711         gener->SetTrigger(1);
712         gener->SetPtJet(200.);
713         gener->SetRadiation(isw);
714         gener->SetSimpleJets(!isw);
715         gener->SetJetEtaRange(-0.3,0.3);
716         gener->SetJetPhiRange(75., 165.);   
717         gGener=gener;
718       }
719       break;
720 //
721 // Gamma-Jet
722 //
723     case kHijing_gj25:
724       {
725         comment = comment.Append("HIJING Gamma 25 GeV");
726         AliGenHijing *gener = HijingStandard();
727 // impact parameter range
728         gener->SetImpactParameterRange(0., 5.);
729         // trigger
730         gener->SetTrigger(2);
731         gener->SetPtJet(25.);
732         gener->SetRadiation(isw);
733         gener->SetSimpleJets(!isw);
734         gener->SetJetEtaRange(-0.12, 0.12);
735         gener->SetJetPhiRange(220., 320.);
736         gGener=gener;
737       }
738       break;
739
740     case kHijing_gj50:
741       {
742         comment = comment.Append("HIJING Gamma 50 GeV");
743         AliGenHijing *gener = HijingStandard();
744 // impact parameter range
745         gener->SetImpactParameterRange(0., 5.);
746         // trigger
747         gener->SetTrigger(2);
748         gener->SetPtJet(50.);
749         gener->SetRadiation(isw);
750         gener->SetSimpleJets(!isw);
751         gener->SetJetEtaRange(-0.12, 0.12);
752         gener->SetJetPhiRange(220., 320.);
753         gGener=gener;
754       }
755       break;
756
757     case kHijing_gj75:
758       {
759         comment = comment.Append("HIJING Gamma 75 GeV");
760         AliGenHijing *gener = HijingStandard();
761 // impact parameter range
762         gener->SetImpactParameterRange(0., 5.);
763         // trigger
764         gener->SetTrigger(2);
765         gener->SetPtJet(75.);
766         gener->SetRadiation(isw);
767         gener->SetSimpleJets(!isw);
768         gener->SetJetEtaRange(-0.12, 0.12);
769         gener->SetJetPhiRange(220., 320.);
770         gGener=gener;
771       }
772       break;
773
774     case kHijing_gj100:
775       {
776         comment = comment.Append("HIJING Gamma 100 GeV");
777         AliGenHijing *gener = HijingStandard();
778 // impact parameter range
779         gener->SetImpactParameterRange(0., 5.);
780         // trigger
781         gener->SetTrigger(2);
782         gener->SetPtJet(100.);
783         gener->SetRadiation(isw);
784         gener->SetSimpleJets(!isw);
785         gener->SetJetEtaRange(-0.12, 0.12);
786         gener->SetJetPhiRange(220., 320.);
787         gGener=gener;
788       }
789       break;
790
791     case kHijing_gj200:
792       {
793         comment = comment.Append("HIJING Gamma 200 GeV");
794         AliGenHijing *gener = HijingStandard();
795 // impact parameter range
796         gener->SetImpactParameterRange(0., 5.);
797         // trigger
798         gener->SetTrigger(2);
799         gener->SetPtJet(200.);
800         gener->SetRadiation(isw);
801         gener->SetSimpleJets(!isw);
802         gener->SetJetEtaRange(-0.12, 0.12);
803         gener->SetJetPhiRange(220., 320.);
804         gGener=gener;
805       }
806       break;
807     case kHijing_pA:
808       {
809         comment = comment.Append("HIJING pA");
810
811         AliGenCocktail *gener  = new AliGenCocktail();
812
813         AliGenHijing   *hijing = new AliGenHijing(-1);
814 // centre of mass energy 
815         hijing->SetEnergyCMS(TMath::Sqrt(82./208.) * 14000.);
816 // impact parameter range
817         hijing->SetImpactParameterRange(0., 15.);
818 // reference frame
819         hijing->SetReferenceFrame("CMS");
820         hijing->SetBoostLHC(1);
821 // projectile
822         hijing->SetProjectile("P", 1, 1);
823         hijing->SetTarget    ("A", 208, 82);
824 // tell hijing to keep the full parent child chain
825         hijing->KeepFullEvent();
826 // enable jet quenching
827         hijing->SetJetQuenching(0);
828 // enable shadowing
829         hijing->SetShadowing(1);
830 // Don't track spectators
831         hijing->SetSpectators(0);
832 // kinematic selection
833         hijing->SetSelectAll(0);
834 //
835         AliGenSlowNucleons*  gray    = new AliGenSlowNucleons(1);
836         AliSlowNucleonModel* model   = new AliSlowNucleonModelExp();
837         gray->SetSlowNucleonModel(model);
838         gray->SetDebug(1);
839         gener->AddGenerator(hijing,"Hijing pPb", 1);
840         gener->AddGenerator(gray,  "Gray Particles",1);
841         gGener=gener;
842       }
843       break;
844     case kPythia6:
845       {
846         comment = comment.Append(":Pythia p-p @ 14 TeV");
847         AliGenPythia *gener = new AliGenPythia(-1); 
848         gener->SetMomentumRange(0,999999);
849         gener->SetPhiRange(-180,180);
850         gener->SetThetaRange(0., 180.);
851         gener->SetYRange(-12,12);
852         gener->SetPtRange(0,1000);
853         gener->SetStrucFunc(kCTEQ4L);   
854         gener->SetProcess(kPyMb);
855         gener->SetEnergyCMS(14000.);
856         gGener=gener;
857       }
858       break;
859     case kPythia6Jets:
860       {
861         comment = comment.Append(":Pythia jets @ 5.5 TeV");
862         AliGenPythia * gener = new AliGenPythia(-1);
863         gener->SetEnergyCMS(5500.);//        Centre of mass energy
864         gener->SetProcess(kPyJets);//        Process type
865         gener->SetJetEtaRange(-0.5, 0.5);//  Final state kinematic cuts
866         gener->SetJetPhiRange(0., 360.);
867         gener->SetJetEtRange(10., 1000.);
868         gener->SetGluonRadiation(1,1);
869         //    gener->SetPtKick(0.);
870         //   Structure function
871         gener->SetStrucFunc(kCTEQ4L);
872         gener->SetPtHard(86., 104.);// Pt transfer of the hard scattering
873         gener->SetPycellParameters(2., 274, 432, 0., 4., 5., 1.0);
874         gener->SetForceDecay(kAll);//  Decay type (semielectronic, etc.)
875         gGener=gener;
876       }
877       break;
878     case kD0PbPb5500:
879       {
880         comment = comment.Append(" D0 in Pb-Pb at 5.5 TeV");
881         AliGenPythia * gener = new AliGenPythia(10);
882         gener->SetProcess(kPyD0PbPbMNR);
883         gener->SetStrucFunc(kCTEQ4L);
884         gener->SetPtHard(2.1,-1.0);
885         gener->SetEnergyCMS(5500.);
886         gener->SetNuclei(208,208);
887         gGener=gener;
888       }
889       break;
890     case kD_TRD:
891       {
892         comment = comment.Append(" D0 for TRD at 5.5 TeV");
893         AliGenPythia *gener = new AliGenPythia(1);
894         gener->SetCutOnChild(0);
895         gener->SetStrucFunc(kCTEQ4L);
896         gener->SetProcess(kPyCharm);
897         gener->SetPtHard(0.,-1);
898         gener->SetEnergyCMS(5500.);
899         gGener=gener;
900       }
901       break;
902     case kB_TRD:
903       {
904         comment = comment.Append(" B for TRD at 5.5 TeV");
905         AliGenPythia *gener = new AliGenPythia(1);
906         gener->SetCutOnChild(0);
907         gener->SetStrucFunc(kCTEQ4L);
908         gener->SetProcess(kPyBeauty);
909         gener->SetPtHard(0.,-1);
910         gener->SetEnergyCMS(5500.);
911         gGener=gener;
912       }
913       break;
914     case kJpsi_TRD:
915       {
916         comment = comment.Append(" J/psi for TRD at 5.5 TeV");
917         AliGenParam *gener = new AliGenParam(1,new AliGenGSIlib(),
918                                              AliGenGSIlib::kJPsi,"MUON");
919         gener->SetMomentumRange(0,999);
920         gener->SetPtRange(0,30.);
921         gener->SetPhiRange(0., 360.);
922         gener->SetYRange(-0.9,+0.9);
923         gener->SetForceDecay(kDiElectron);
924         gGener=gener;
925       }
926       break;
927     case kU_TRD:
928       {
929         comment = comment.Append(" Upsilon for TRD at 5.5 TeV");
930         AliGenParam *gener = new AliGenParam(1,new AliGenGSIlib(),
931                                              AliGenGSIlib::kUpsilon,"RITMAN");
932         gener->SetMomentumRange(0,999);
933         gener->SetPtRange(0,30.);
934         gener->SetPhiRange(0., 360.);
935         gener->SetYRange(-0.9,0.9);
936         gener->SetForceDecay(kDiElectron);
937         gGener=gener;
938       }
939       break;
940     case kPyJJ:
941       {
942         comment = comment.Append(" Jet-jet at 5.5 TeV");
943         AliGenPythia *gener = new AliGenPythia(-1);
944         gener->SetEnergyCMS(5500.);
945         gener->SetProcess(kPyJets);
946         Double_t ptHardMin=10.0, ptHardMax=-1.0;
947         gener->SetPtHard(ptHardMin,ptHardMax);
948         gener->SetYHard(-0.7,0.7);
949         gener->SetJetEtaRange(-0.2,0.2);
950         gener->SetEventListRange(0,1);
951         gGener=gener;
952       }
953       break;
954     case kPyGJ:
955       {
956         comment = comment.Append(" Gamma-jet at 5.5 TeV");
957         AliGenPythia *gener = new AliGenPythia(-1);
958         gener->SetEnergyCMS(5500.);
959         gener->SetProcess(kPyDirectGamma);
960         Double_t ptHardMin=10.0, ptHardMax=-1.0;
961         gener->SetPtHard(ptHardMin,ptHardMax);
962         gener->SetYHard(-1.0,1.0);
963         gener->SetGammaEtaRange(-0.13,0.13);
964         gener->SetGammaPhiRange(210.,330.);
965         gener->SetEventListRange(0,1);
966         gGener=gener;
967       }
968       break;
969     default: break;
970     }
971     return gGener;
972 }
973
974 AliGenHijing* HijingStandard()
975 {
976     AliGenHijing *gener = new AliGenHijing(-1);
977 // centre of mass energy 
978     gener->SetEnergyCMS(5500.);
979 // reference frame
980     gener->SetReferenceFrame("CMS");
981 // projectile
982      gener->SetProjectile("A", 208, 82);
983      gener->SetTarget    ("A", 208, 82);
984 // tell hijing to keep the full parent child chain
985      gener->KeepFullEvent();
986 // enable jet quenching
987      gener->SetJetQuenching(1);
988 // enable shadowing
989      gener->SetShadowing(1);
990 // neutral pion and heavy particle decays switched off
991      gener->SetDecaysOff(1);
992 // Don't track spectators
993      gener->SetSpectators(0);
994 // kinematic selection
995      gener->SetSelectAll(0);
996      return gener;
997 }
998
999
1000 void ProcessEnvironmentVars()
1001 {
1002     // Run type
1003     if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
1004       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
1005         if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
1006           srun = (PprRun_t)iRun;
1007           cout<<"Run type set to "<<pprRunName[iRun]<<endl;
1008         }
1009       }
1010     }
1011
1012     // Random Number seed
1013     if (gSystem->Getenv("CONFIG_SEED")) {
1014       sseed = atoi(gSystem->Getenv("CONFIG_SEED"));
1015     }
1016 }