]> git.uio.no Git - u/mrichter/AliRoot.git/blob - prod/LHC09d5/Config.C
b9cf8e85a9b886dc6a0d720a70709d0cf77e67ca
[u/mrichter/AliRoot.git] / prod / LHC09d5 / 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, kPhojet, kRunMax
55 };
56
57 const char * pprRunName[] = {
58     "kPythia6", "kPythia6D6T", "kPhojet" 
59 };
60
61 enum Mag_t
62 {
63   kNoField, k5kG, kFieldMax
64 };
65
66 const char * pprField[] = {
67   "kNoField", "k5kG"
68 };
69
70 //--- Functions ---
71 class AliGenPythia;
72 AliGenerator *MbPythia();
73 AliGenerator *MbPythiaTuneD6T();
74 AliGenerator *MbPhojet();
75 void ProcessEnvironmentVars();
76
77 // Geterator, field, beam energy
78 static PDC06Proc_t   proc     = kPhojet;
79 static Mag_t         mag      = k5kG;
80 static Float_t       energy   = 10000; // energy in CMS
81 //========================//
82 // Set Random Number seed //
83 //========================//
84 TDatime dt;
85 static UInt_t seed    = dt.Get();
86
87 // Comment line
88 static TString comment;
89
90 void Config()
91 {
92     
93
94   // Get settings from environment variables
95   ProcessEnvironmentVars();
96
97   gRandom->SetSeed(seed);
98   cerr<<"Seed for random number generation= "<<seed<<endl; 
99
100   // Libraries required by geant321
101 #if defined(__CINT__)
102   gSystem->Load("liblhapdf");      // Parton density functions
103   gSystem->Load("libEGPythia6");   // TGenerator interface
104   if (proc != kPythia6D6T) {
105       gSystem->Load("libpythia6");     // Pythia 6.2
106   } else {
107       gSystem->Load("libqpythia");     // Pythia 6.4
108   }
109   gSystem->Load("libAliPythia6");  // ALICE specific implementations
110   gSystem->Load("libgeant321");
111 #endif
112
113   new TGeant3TGeo("C++ Interface to Geant3");
114
115   //=======================================================================
116   //  Create the output file
117
118    
119   AliRunLoader* rl=0x0;
120
121   cout<<"Config.C: Creating Run Loader ..."<<endl;
122   rl = AliRunLoader::Open("galice.root",
123                           AliConfig::GetDefaultEventFolderName(),
124                           "recreate");
125   if (rl == 0x0)
126     {
127       gAlice->Fatal("Config.C","Can not instatiate the Run Loader");
128       return;
129     }
130   rl->SetCompressionLevel(2);
131   rl->SetNumberOfEventsPerFile(1000);
132   gAlice->SetRunLoader(rl);
133   // gAlice->SetGeometryFromFile("geometry.root");
134   // gAlice->SetGeometryFromCDB();
135   
136   // Set the trigger configuration: proton-proton
137   gAlice->SetTriggerDescriptor("p-p");
138
139   //
140   //=======================================================================
141   // ************* STEERING parameters FOR ALICE SIMULATION **************
142   // --- Specify event type to be tracked through the ALICE setup
143   // --- All positions are in cm, angles in degrees, and P and E in GeV
144
145
146     gMC->SetProcess("DCAY",1);
147     gMC->SetProcess("PAIR",1);
148     gMC->SetProcess("COMP",1);
149     gMC->SetProcess("PHOT",1);
150     gMC->SetProcess("PFIS",0);
151     gMC->SetProcess("DRAY",0);
152     gMC->SetProcess("ANNI",1);
153     gMC->SetProcess("BREM",1);
154     gMC->SetProcess("MUNU",1);
155     gMC->SetProcess("CKOV",1);
156     gMC->SetProcess("HADR",1);
157     gMC->SetProcess("LOSS",2);
158     gMC->SetProcess("MULS",1);
159     gMC->SetProcess("RAYL",1);
160
161     Float_t cut = 1.e-3;        // 1MeV cut by default
162     Float_t tofmax = 1.e10;
163
164     gMC->SetCut("CUTGAM", cut);
165     gMC->SetCut("CUTELE", cut);
166     gMC->SetCut("CUTNEU", cut);
167     gMC->SetCut("CUTHAD", cut);
168     gMC->SetCut("CUTMUO", cut);
169     gMC->SetCut("BCUTE",  cut); 
170     gMC->SetCut("BCUTM",  cut); 
171     gMC->SetCut("DCUTE",  cut); 
172     gMC->SetCut("DCUTM",  cut); 
173     gMC->SetCut("PPCUTM", cut);
174     gMC->SetCut("TOFMAX", tofmax); 
175
176
177
178
179   //======================//
180   // Set External decayer //
181   //======================//
182   TVirtualMCDecayer* decayer = new AliDecayerPythia();
183   decayer->SetForceDecay(kAll);
184   decayer->Init();
185   gMC->SetExternalDecayer(decayer);
186
187   //=========================//
188   // Generator Configuration //
189   //=========================//
190   AliGenerator* gener = 0x0;
191   
192   if (proc == kPythia6) {
193       gener = MbPythia();
194   } else if (proc == kPythia6D6T) {
195       gener = MbPythiaTuneD6T();
196   } else if (proc == kPhojet) {
197       gener = MbPhojet();
198   }
199   
200   
201
202   // PRIMARY VERTEX
203   //
204   gener->SetOrigin(-0.07, 0.25, -1.5);    // vertex position
205   //
206   //
207   // Size of the interaction diamond
208   // Longitudinal
209   Float_t sigmaz  = 5.4 / TMath::Sqrt(2.); // [cm]
210   if (energy == 900)
211     sigmaz  = 10.5 / TMath::Sqrt(2.); // [cm]
212
213   if (energy == 7000)
214     sigmaz  = 6.3 / TMath::Sqrt(2.); // [cm]
215
216   //
217   // Transverse
218   Float_t betast  = 10;                 // beta* [m]
219   Float_t eps     = 3.75e-6;            // emittance [m]
220   Float_t gamma   = energy / 2.0 / 0.938272;  // relativistic gamma [1]
221   Float_t sigmaxy = TMath::Sqrt(eps * betast / gamma) / TMath::Sqrt(2.) * 100.;  // [cm]
222   printf("\n \n Diamond size x-y: %10.3e z: %10.3e\n \n", sigmaxy, sigmaz);
223     
224   gener->SetSigma(sigmaxy, sigmaxy, sigmaz);      // Sigma in (X,Y,Z) (cm) on IP position
225   gener->SetCutVertexZ(3.);        // Truncate at 3 sigma
226   gener->SetVertexSmear(kPerEvent);
227
228   gener->Init();
229   if (proc == kPythia6D6T) {
230         // F I X 
231     AliPythia::Instance()->Pytune(109); 
232     AliPythia::Instance()->Initialize("CMS","p","p", energy);
233         // F I X
234   }
235   // FIELD
236   //
237   AliMagF* field = 0x0;
238   if (mag == kNoField) {
239     comment = comment.Append(" | L3 field 0.0 T");
240     field = new AliMagF("Maps","Maps", 0., 0., AliMagF::k5kGUniform,AliMagF::kBeamTypepp, energy/2.0);
241   } else if (mag == k5kG) {
242     comment = comment.Append(" | L3 field 0.5 T");
243     field = new AliMagF("Maps","Maps", -1., -1., AliMagF::k5kG, AliMagF::kBeamTypepp, energy/2.0);
244   }
245
246   printf("\n \n Comment: %s \n \n", comment.Data());
247
248   TGeoGlobalMagField::Instance()->SetField(field);
249     
250   rl->CdGAFile();
251   
252   Int_t iABSO  = 1;
253   Int_t iACORDE= 0;
254   Int_t iDIPO  = 1;
255   Int_t iEMCAL = 1;
256   Int_t iFMD   = 1;
257   Int_t iFRAME = 1;
258   Int_t iHALL  = 1;
259   Int_t iITS   = 1;
260   Int_t iMAG   = 1;
261   Int_t iMUON  = 1;
262   Int_t iPHOS  = 1;
263   Int_t iPIPE  = 1;
264   Int_t iPMD   = 1;
265   Int_t iHMPID = 1;
266   Int_t iSHIL  = 1;
267   Int_t iT0    = 1;
268   Int_t iTOF   = 1;
269   Int_t iTPC   = 1;
270   Int_t iTRD   = 1;
271   Int_t iVZERO = 1;
272   Int_t iZDC   = 1;
273   
274
275     //=================== Alice BODY parameters =============================
276     AliBODY *BODY = new AliBODY("BODY", "Alice envelop");
277
278
279     if (iMAG)
280     {
281         //=================== MAG parameters ============================
282         // --- Start with Magnet since detector layouts may be depending ---
283         // --- on the selected Magnet dimensions ---
284         AliMAG *MAG = new AliMAG("MAG", "Magnet");
285     }
286
287
288     if (iABSO)
289     {
290         //=================== ABSO parameters ============================
291         AliABSO *ABSO = new AliABSOv3("ABSO", "Muon Absorber");
292     }
293
294     if (iDIPO)
295     {
296         //=================== DIPO parameters ============================
297
298         AliDIPO *DIPO = new AliDIPOv3("DIPO", "Dipole version 3");
299     }
300
301     if (iHALL)
302     {
303         //=================== HALL parameters ============================
304
305         AliHALL *HALL = new AliHALLv3("HALL", "Alice Hall");
306     }
307
308
309     if (iFRAME)
310     {
311         //=================== FRAME parameters ============================
312
313         AliFRAMEv2 *FRAME = new AliFRAMEv2("FRAME", "Space Frame");
314         FRAME->SetHoles(1);
315     }
316
317     if (iSHIL)
318     {
319         //=================== SHIL parameters ============================
320
321         AliSHIL *SHIL = new AliSHILv3("SHIL", "Shielding Version 3");
322     }
323
324
325     if (iPIPE)
326     {
327         //=================== PIPE parameters ============================
328
329         AliPIPE *PIPE = new AliPIPEv3("PIPE", "Beam Pipe");
330     }
331  
332     if (iITS)
333     {
334         //=================== ITS parameters ============================
335
336         AliITS *ITS  = new AliITSv11Hybrid("ITS","ITS v11Hybrid");
337     }
338
339     if (iTPC)
340     {
341       //============================ TPC parameters =====================
342
343         AliTPC *TPC = new AliTPCv2("TPC", "Default");
344     }
345
346
347     if (iTOF) {
348         //=================== TOF parameters ============================
349
350         AliTOF *TOF = new AliTOFv6T0("TOF", "normal TOF");
351     }
352
353
354     if (iHMPID)
355     {
356         //=================== HMPID parameters ===========================
357
358         AliHMPID *HMPID = new AliHMPIDv3("HMPID", "normal HMPID");
359
360     }
361
362
363     if (iZDC)
364     {
365         //=================== ZDC parameters ============================
366
367         AliZDC *ZDC = new AliZDCv3("ZDC", "normal ZDC");
368     }
369
370     if (iTRD)
371     {
372         //=================== TRD parameters ============================
373
374         AliTRD *TRD = new AliTRDv1("TRD", "TRD slow simulator");
375         AliTRDgeometry *geoTRD = TRD->GetGeometry();
376         // Partial geometry: modules at 0,1,7,8,9,16,17
377         // starting at 3h in positive direction
378         geoTRD->SetSMstatus(2,0);
379         geoTRD->SetSMstatus(3,0);
380         geoTRD->SetSMstatus(4,0);
381         geoTRD->SetSMstatus(5,0);
382         geoTRD->SetSMstatus(6,0);
383         geoTRD->SetSMstatus(10,0);
384         geoTRD->SetSMstatus(11,0);
385         geoTRD->SetSMstatus(12,0);
386         geoTRD->SetSMstatus(13,0);
387         geoTRD->SetSMstatus(14,0);
388         geoTRD->SetSMstatus(15,0);
389     }
390
391     if (iFMD)
392     {
393         //=================== FMD parameters ============================
394
395         AliFMD *FMD = new AliFMDv1("FMD", "normal FMD");
396    }
397
398     if (iMUON)
399     {
400         //=================== MUON parameters ===========================
401         // New MUONv1 version (geometry defined via builders)
402
403         AliMUON *MUON = new AliMUONv1("MUON", "default");
404     }
405
406     if (iPHOS)
407     {
408         //=================== PHOS parameters ===========================
409
410         AliPHOS *PHOS = new AliPHOSv1("PHOS", "IHEP");
411         //Set simulation parameters different from the default ones.
412         AliPHOSSimParam* simEmc = AliPHOSSimParam::GetInstance() ;
413   
414         // APD noise of warm (+20C) PHOS:
415         // a2 = a1*(Y1/Y2)*(M1/M2), where a1 = 0.012 is APD noise at -25C,
416         // Y1 = 4.3 photo-electrons/MeV, Y2 = 1.7 p.e/MeV - light yields at -25C and +20C,
417         // M1 = 50, M2 = 50 - APD gain factors chosen for t1 = -25C and t2 = +20C,
418         // Y = MeanLightYield*APDEfficiency.
419
420         Float_t apdNoise = 0.012*2.5; 
421         simEmc->SetAPDNoise(apdNoise);
422
423         //Raw Light Yield at +20C
424         simEmc->SetMeanLightYield(18800);
425
426         //ADC channel width at +18C.
427         simEmc->SetADCchannelW(0.0125);
428     }
429
430
431     if (iPMD)
432     {
433         //=================== PMD parameters ============================
434
435         AliPMD *PMD = new AliPMDv1("PMD", "normal PMD");
436     }
437
438     if (iT0)
439     {
440         //=================== T0 parameters ============================
441         AliT0 *T0 = new AliT0v1("T0", "T0 Detector");
442     }
443
444     if (iEMCAL)
445     {
446         //=================== EMCAL parameters ============================
447
448         AliEMCAL *EMCAL = new AliEMCALv2("EMCAL", "EMCAL_COMPLETE");
449     }
450
451      if (iACORDE)
452     {
453         //=================== ACORDE parameters ============================
454
455         AliACORDE *ACORDE = new AliACORDEv1("ACORDE", "normal ACORDE");
456     }
457
458      if (iVZERO)
459     {
460         //=================== ACORDE parameters ============================
461
462         AliVZERO *VZERO = new AliVZEROv7("VZERO", "normal VZERO");
463     }
464 }
465 //
466 //           PYTHIA
467 //
468
469 AliGenerator* MbPythia()
470 {
471       comment = comment.Append(" pp at 14 TeV: Pythia low-pt");
472 //
473 //    Pythia
474       AliGenPythia* pythia = new AliGenPythia(-1); 
475       pythia->SetMomentumRange(0, 999999.);
476       pythia->SetThetaRange(0., 180.);
477       pythia->SetYRange(-12.,12.);
478       pythia->SetPtRange(0,1000.);
479       pythia->SetProcess(kPyMb);
480       pythia->SetEnergyCMS(energy);
481       
482       return pythia;
483 }
484
485 AliGenerator* MbPythiaTuneD6T()
486 {
487       comment = comment.Append(" pp at 14 TeV: Pythia low-pt");
488 //
489 //    Pythia
490       AliGenPythia* pythia = new AliGenPythia(-1); 
491       pythia->SetMomentumRange(0, 999999.);
492       pythia->SetThetaRange(0., 180.);
493       pythia->SetYRange(-12.,12.);
494       pythia->SetPtRange(0,1000.);
495       pythia->SetProcess(kPyMb);
496       pythia->SetEnergyCMS(energy);
497 //    Tune
498 //    109     D6T : Rick Field's CDF Tune D6T (NB: needs CTEQ6L pdfs externally)
499 //      pythia->SetTune(109); // F I X 
500       pythia->SetStrucFunc(kCTEQ6l);
501 //
502       return pythia;
503 }
504
505 AliGenerator* MbPhojet()
506 {
507       comment = comment.Append(" pp at 14 TeV: Pythia low-pt");
508 //
509 //    DPMJET
510 #if defined(__CINT__)
511   gSystem->Load("libdpmjet");      // Parton density functions
512   gSystem->Load("libTDPMjet");      // Parton density functions
513 #endif
514       AliGenDPMjet* dpmjet = new AliGenDPMjet(-1); 
515       dpmjet->SetMomentumRange(0, 999999.);
516       dpmjet->SetThetaRange(0., 180.);
517       dpmjet->SetYRange(-12.,12.);
518       dpmjet->SetPtRange(0,1000.);
519       dpmjet->SetProcess(kDpmMb);
520       dpmjet->SetEnergyCMS(energy);
521
522       return dpmjet;
523 }
524
525 void ProcessEnvironmentVars()
526 {
527     // Run type
528     if (gSystem->Getenv("CONFIG_RUN_TYPE")) {
529       for (Int_t iRun = 0; iRun < kRunMax; iRun++) {
530         if (strcmp(gSystem->Getenv("CONFIG_RUN_TYPE"), pprRunName[iRun])==0) {
531           proc = (PDC06Proc_t)iRun;
532           cout<<"Run type set to "<<pprRunName[iRun]<<endl;
533         }
534       }
535     }
536
537     // Field
538     if (gSystem->Getenv("CONFIG_FIELD")) {
539       for (Int_t iField = 0; iField < kFieldMax; iField++) {
540         if (strcmp(gSystem->Getenv("CONFIG_FIELD"), pprField[iField])==0) {
541           mag = (Mag_t)iField;
542           cout<<"Field set to "<<pprField[iField]<<endl;
543         }
544       }
545     }
546
547     // Energy
548     if (gSystem->Getenv("CONFIG_ENERGY")) {
549       energy = atoi(gSystem->Getenv("CONFIG_ENERGY"));
550       cout<<"Energy set to "<<energy<<" GeV"<<endl;
551     }
552
553     // Random Number seed
554     if (gSystem->Getenv("CONFIG_SEED")) {
555       seed = atoi(gSystem->Getenv("CONFIG_SEED"));
556     }
557 }