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