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